Introduction to atomman: LAMMPS data file conversions

Lucas M. Hale, lucas.hale@nist.gov, Materials Science and Engineering Division, NIST.

Disclaimers

1. Introduction

The atom_data format is the atomic data files used by LAMMPS for importing initial configurations. Currently, atomman offers partial support for the atom_data format based on atomman’s core features. See the LAMMPS documentation for more details on the format.

The supported atom_data fields are

  • Header sections:

    • atoms: number of atoms. Required during load, always listed during dump.

    • atom types: number of atom types. Required during load if Masses is included, always listed during dump.

    • xlo xhi: simulation box boundaries in x dimension. Required during load, always listed during dump.

    • ylo yhi: simulation box boundaries in y dimension. Required during load, always listed during dump.

    • zlo zhi: simulation box boundaries in z dimension. Required during load, always listed during dump.

    • xy xz yz: simulation box tilt factors for triclinic system. Optional during load, listed during dump only if the tilt factors have non-zero values.

  • Body sections:

    • Atoms: lists per-atom properties associated with each atom. The mapping of LAMMPS<->atomman representations of the per-atom properties can be found in atomman.load.atom_data.atoms_prop_info and atomman.dump.atom_data.atoms_prop_info. Required during load, always listed during dump.

    • Velocities: list per-atom velocity properties associated with each atom. The mapping of LAMMPS<->atomman representations of the per-atom properties can be found in atomman.load.atom_data.velocities_prop_info and atomman.dump.atom_data.velocities_prop_info. Optional during load, will be listed during dump if ‘velocity’ is set as a per-atom property.

    • Masses: masses per each atom type. Optional during load, where values are saved to system.masses. As for dump, the Masses section is never included because atomman.lammps.Potential generates LAMMPS mass input command lines. See Section 3.4 below for more details.

Notes on ignored fields:

  • The Pair Coeffs and PairIJ Coeffs fields are ignored as potential parameters are handled with the atomman.lammps.Potential class.

  • Support for bonds, angles, dihedrals, impropers, ellipsoids, lines, triangles, and bodies would require defining how to represent these as objects in atomman. Extending support is possible, but only if there is enough interest and help from experts of these data types.

Library Imports

[1]:
# Standard Python libraries
import os
import datetime

# http://www.numpy.org/
import numpy as np

import atomman as am
import atomman.lammps as lmp
import atomman.unitconvert as uc

# Show atomman version
print('atomman version =', am.__version__)

# Show date of Notebook execution
print('Notebook executed on', datetime.date.today())
atomman version = 1.4.10
Notebook executed on 2023-07-28

2. atomman.load(‘atom_data’)

Update version 1.3.0: Support for image flags, Masses, and atom_style values listed in the data file are added. Errors returned by the method are updated to be more meaningful.

Parameters

  • data (str or file-like object) The atom data content to read. Can be str content, path name, or open file-like object.

  • pbc (list of bool) Three boolean values indicating which System directions are periodic. Default value is (True, True, True).

  • symbols (tuple, optional) Allows the list of element symbols to be assigned during loading.

  • atom_style (str, optional) The LAMMPS atom_style option associated with the data file. Optional as the data can list this value in a comment in the Atoms section header. If not given and not found in data, the default value of ‘atomic’ is used.

  • units (str, optional) The LAMMPS units option associated with the data file. Default value is ‘metal’.

Returns

  • (atomman.System) The corresponding system. Note all property values will be automatically converted to atomman.unitconvert’s working units.

Raises

  • (FileNotFoundError) If data is (likely) a file name and no matching file is found.

  • (FileFormatError) If required content fields not found.

  • (ValueError) If atom_style is both given as a parameter and found in data, but are not the same

2.1. Basic example

bcc Fe unit cell without masses and default atom_style=’atomic’ assumed.

[2]:
bccFe = """
2 atoms
1 atom types
0.0000000000000 2.8665000000000 xlo xhi
0.0000000000000 2.8665000000000 ylo yhi
0.0000000000000 2.8665000000000 zlo zhi

Masses

1 55.845

Atoms

1 1 0.0000000000000 0.0000000000000 0.0000000000000
2 1 1.4332500000000 1.4332500000000 1.4332500000000
"""

# Symbols specified here as file format does not include them
system = am.load('atom_data', bccFe, symbols='Fe')
print(system)
avect =  [ 2.866,  0.000,  0.000]
bvect =  [ 0.000,  2.866,  0.000]
cvect =  [ 0.000,  0.000,  2.866]
origin = [ 0.000,  0.000,  0.000]
natoms = 2
natypes = 1
symbols = ('Fe',)
pbc = [ True  True  True]
per-atom properties = ['atype', 'pos']
     id   atype  pos[0]  pos[1]  pos[2]
      0       1   0.000   0.000   0.000
      1       1   1.433   1.433   1.433

Here, Masses were listed in the file and are saved to system.masses. If the Masses section was not included, all system.masses values would be None.

[3]:
print('masses =',system.masses)
masses = (55.845,)

2.2. Different atom_styles

Files that use atom_style formats different from ‘atomic’ can be read if either the atom_style parameter is specified or if the atom_style is listed as a comment after the “Atoms” header line in the data.

[4]:
# Example 1 for atom_style charge
bccFe = """
2 atoms
1 atom types
0.0000000000000 2.8665000000000 xlo xhi
0.0000000000000 2.8665000000000 ylo yhi
0.0000000000000 2.8665000000000 zlo zhi

Masses

1 55.845

Atoms # charge

1 1 0.0000000000000 0.0000000000000 0.0000000000000 0.0000000000000
2 1 0.0000000000000 1.4332500000000 1.4332500000000 1.4332500000000
"""

# atom_style not needed as 'charge' is specified in the data
system = am.load('atom_data', bccFe, symbols='Fe')
print(system)

# Show that charges are set
print()
print('system.atoms.charge =',system.atoms.charge)
avect =  [ 2.866,  0.000,  0.000]
bvect =  [ 0.000,  2.866,  0.000]
cvect =  [ 0.000,  0.000,  2.866]
origin = [ 0.000,  0.000,  0.000]
natoms = 2
natypes = 1
symbols = ('Fe',)
pbc = [ True  True  True]
per-atom properties = ['atype', 'pos', 'charge']
     id   atype  pos[0]  pos[1]  pos[2]
      0       1   0.000   0.000   0.000
      1       1   1.433   1.433   1.433

system.atoms.charge = [0. 0.]
[5]:
# Example 2 for atom_style charge (no comment in data)
bccFe = """
2 atoms
1 atom types
0.0000000000000 2.8665000000000 xlo xhi
0.0000000000000 2.8665000000000 ylo yhi
0.0000000000000 2.8665000000000 zlo zhi

Masses

1 55.845

Atoms

1 1 0.0000000000000 0.0000000000000 0.0000000000000 0.0000000000000
2 1 0.0000000000000 1.4332500000000 1.4332500000000 1.4332500000000
"""

# atom_style='charge' required
system = am.load('atom_data', bccFe, symbols='Fe', atom_style='charge')
print(system)

# Show that charges are set
print()
print('system.atoms.charge =',system.atoms.charge)
avect =  [ 2.866,  0.000,  0.000]
bvect =  [ 0.000,  2.866,  0.000]
cvect =  [ 0.000,  0.000,  2.866]
origin = [ 0.000,  0.000,  0.000]
natoms = 2
natypes = 1
symbols = ('Fe',)
pbc = [ True  True  True]
per-atom properties = ['atype', 'pos', 'charge']
     id   atype  pos[0]  pos[1]  pos[2]
      0       1   0.000   0.000   0.000
      1       1   1.433   1.433   1.433

system.atoms.charge = [0. 0.]

2.3. Image flag handling

Integer image flags can be added to Atoms table to move atoms to periodic replicas

[6]:
bccFe = """
2 atoms
1 atom types
0.0000000000000 2.8665000000000 xlo xhi
0.0000000000000 2.8665000000000 ylo yhi
0.0000000000000 2.8665000000000 zlo zhi

Masses

1 55.845

Atoms

1 1 0.0000000000000 0.0000000000000 0.0000000000000 0 0 0
2 1 1.4332500000000 1.4332500000000 1.4332500000000 1 0 0
"""

# Symbols specified here as file format
system = am.load('atom_data', bccFe, symbols='Fe')
print(system)
avect =  [ 2.866,  0.000,  0.000]
bvect =  [ 0.000,  2.866,  0.000]
cvect =  [ 0.000,  0.000,  2.866]
origin = [ 0.000,  0.000,  0.000]
natoms = 2
natypes = 1
symbols = ('Fe',)
pbc = [ True  True  True]
per-atom properties = ['atype', 'pos']
     id   atype  pos[0]  pos[1]  pos[2]
      0       1   0.000   0.000   0.000
      1       1   4.300   1.433   1.433

3. System.dump(‘atom_data’)

Updated version 1.3.0: Parameters natypes, potential, and return_pair_info are added to provide better integration and compatibility across potentials with different pair_styles.

Updated version 1.4.0: The return_pair_info parameter has been removed and potential commands are included in the returned info simply if potential is given.

Parameters

  • f (str or file-like object, optional) File path or file-like object to write the content to. If not given, then the content is returned as a str.

  • atom_style (str, optional) The LAMMPS atom_style option associated with the data file. If neither atom_style or potential is given, will set atom_style to ‘atomic’.

  • units (str, optional) The LAMMPS units option associated with the data file. If neitherunits or potential is given, will set units ‘metal’.

  • natypes (int, optional) Allows the natypes value to be manually changed. This is needed if natypes needs to be greater than the current number of atypes. If neither natypes or potential is given, will use system.natypes.

  • potential (atomman.lammps.Potential, optional) Potential-specific values of atom_style, units, and natypes can be extracted from a Potential object. If both potential and any of the individual values are given, the individual values will be used.

  • float_format (str, optional) c-style formatting string for floating point values. Default value is ‘%.13f’.

  • return_info (bool, optional) Indicates if the LAMMPS command lines associated with reading in the file are to be returned as a str. If potential is given, then the commands associated with the potential will be included. Default value is True.

  • safecopy (bool, optional) The LAMMPS data format requires all atoms to be inside box bounds, i.e. “wrapped”. If safecopy is True then a copy of the system is made to keep the original unwrapped. Default value is False.

Returns

  • content (str, optional) The data file contents. Returned if f is not given.

  • read_info (str, optional) The LAMMPS input command lines to read the created data file in. Returned if return_info is True. If return_pair_info is also True and potential is given, the LAMMPS input command lines for the potential are also included.

Raises

  • (ValueError) If return_pair_info is True and return_info is False or potential is not given.

3.1. Basic example

Show dump by loading a simple data file, then dumping it to a file.

NOTE: It is highly recommended that the potential be given (Section 3.3.) if you are using the generated LAMMPS command lines. Doing so ensures the proper options and order of command lines across all native LAMMPS and OpenKIM potentials currently in the NIST database.

[7]:
bccFe = """
2 atoms
1 atom types
0.0000000000000 2.8665000000000 xlo xhi
0.0000000000000 2.8665000000000 ylo yhi
0.0000000000000 2.8665000000000 zlo zhi

Atoms

1 1 0.0000000000000 0.0000000000000 0.0000000000000
2 1 1.4332500000000 1.4332500000000 1.4332500000000
"""
system = am.load('atom_data', bccFe, symbols='Fe')

info = system.dump('atom_data', f='bccFe.dat')

# Show created file is identical to original content
with open('bccFe.dat') as f:
    print(f.read())

2 atoms
1 atom types
0.0000000000000 2.8665000000000 xlo xhi
0.0000000000000 2.8665000000000 ylo yhi
0.0000000000000 2.8665000000000 zlo zhi

Atoms # atomic

1 1 0.0000000000000 0.0000000000000 0.0000000000000
2 1 1.4332500000000 1.4332500000000 1.4332500000000

The returned info contains LAMMPS input command lines associated with reading the file in.

[8]:
print(info)
# Script and atom data file prepared using atomman Python package

units None
atom_style None

boundary p p p
read_data bccFe.dat

3.2. Manual modifications

Now, we’ll dump the same system but change natypes, units, and atom_style parameters.

[9]:
# atom_style='charge' requires per-atom property charge to be set
system.atoms.charge = 0

# units='nano' will return pos in nm instead of Angstroms
# natypes=2 only changes the atom types header
info = system.dump('atom_data', f='bccFe.dat', units='nano', natypes=2, atom_style='charge')

# Show created file with specified changes
with open('bccFe.dat') as f:
    print(f.read())

2 atoms
2 atom types
0.0000000000000 0.2866500000000 xlo xhi
0.0000000000000 0.2866500000000 ylo yhi
0.0000000000000 0.2866500000000 zlo zhi

Atoms # charge

1 1 0.0000000000000 0.0000000000000 0.0000000000000 0.0000000000000
2 1 0.0000000000000 0.1433250000000 0.1433250000000 0.1433250000000

Show info with specified changes

[10]:
print(info)
# Script and atom data file prepared using atomman Python package

units None
atom_style None

boundary p p p
read_data bccFe.dat

3.3. Integrated potential-data file handling

The generation of LAMMPS data files and associated input command lines depends on properties of both the atomic configuration and the choice of interatomic potential. To better support this fact, the atom_data dump style supports a more integrated means of generating the content. This is handled by passing dump(‘atom_style’) a Potential object.

Define an example potential with units=’nano’ and atom_style=’charge’

[11]:
potential_json_example_1 = """{
    "potential-LAMMPS": {
        "key": "NOKEY",
        "id": "demo-example-1",
        "potential": {
            "key": "NOKEY",
            "id": "demo-example" },
        "units": "nano",
        "atom_style": "charge",
        "atom": {
            "element": "Fe" },
        "pair_style": {
            "type": "eam/alloy" },
        "pair_coeff": {
            "term": [
                {
                    "file": "FeDemo.eam.alloy" },
                {
                    "symbols": "True" } ] } } }"""

# Load as Potential
potential = lmp.Potential(potential_json_example_1)

Create data file and info by supplying potential

[12]:
info = system.dump('atom_data', f='bccFe.dat', potential=potential)

# Show created file uses settings from potential
with open('bccFe.dat') as f:
    print(f.read())

2 atoms
1 atom types
0.0000000000000 0.2866500000000 xlo xhi
0.0000000000000 0.2866500000000 ylo yhi
0.0000000000000 0.2866500000000 zlo zhi

Atoms # charge

1 1 0.0000000000000 0.0000000000000 0.0000000000000 0.0000000000000
2 1 0.0000000000000 0.1433250000000 0.1433250000000 0.1433250000000

The LAMMPS command lines for the potential are now included in the returned info.

[13]:
print(info)
units nano
atom_style charge

boundary p p p
read_data bccFe.dat

pair_style eam/alloy
pair_coeff * * FeDemo.eam.alloy Fe
mass 1 55.845


3.4. Notes on masses

The Masses section is not included in the generated data file because atomman specifies mass values by creating LAMMPS mass command lines. This allows for default mass values to be defined based on the chosen potential’s symbols/elements. The mass command lines can be obtained using either Potential.pair_style(), or as described in the last section.

The masses that are listed in the generated mass command lines are selected according to the following order of precedence

  1. Any values of System.masses that are not None.

  2. The value of mass specified in a Potential data model for the atom type’s symbol.

  3. The standard reference mass for the element associated with the atom type’s symbol.

In the previous example, the mass listed was taken as the standard reference value of Fe: 55.845. Some potentials define specific mass values to use, which can be listed in the data model on a per model basis.

[14]:
# mass now explicitly defined for Fe
potential_json_example_2 = """{
    "potential-LAMMPS": {
        "key": "NOKEY",
        "id": "demo-example-2",
        "potential": {
            "key": "NOKEY",
            "id": "demo-example" },
        "units": "nano",
        "atom_style": "charge",
        "atom": {
            "element": "Fe",
            "mass": 55.85 },
        "pair_style": {
            "type": "eam/alloy" },
        "pair_coeff": {
            "term": [
                {
                    "file": "FeDemo.eam.alloy" },
                {
                    "symbols": "True" } ] } } }"""

# Load as Potential
potential = lmp.Potential(potential_json_example_2)
[15]:
# Generate data and info
info = system.dump('atom_data', f='bccFe.dat', potential=potential)

# Show that mass in info is now 55.85
print(info)
units nano
atom_style charge

boundary p p p
read_data bccFe.dat

pair_style eam/alloy
pair_coeff * * FeDemo.eam.alloy Fe
mass 1 55.85


Mass values can also be manually overridden by setting system.masses. This is particularly useful for simulations of different isotopes.

[16]:
# Set mass for Fe-58 isotope
system.masses = am.tools.atomic_mass('Fe-58')

# Generate data and info
info = system.dump('atom_data', f='bccFe.dat', potential=potential)

# Show that mass in info is now 57.93327443
print(info)
units nano
atom_style charge

boundary p p p
read_data bccFe.dat

pair_style eam/alloy
pair_coeff * * FeDemo.eam.alloy Fe
mass 1 57.93327443


cleanup

[17]:
os.remove('bccFe.dat')