LAMMPS File Writers

Write Lammps data

LAMMPS data format.

mbuild.formats.lammpsdata.write_lammpsdata(structure, filename, atom_style='full', unit_style='real', mins=None, maxs=None, pair_coeff_label=None, detect_forcefield_style=True, nbfix_in_data_file=True, use_urey_bradleys=False, use_rb_torsions=True, use_dihedrals=False, sigma_conversion_factor=None, epsilon_conversion_factor=None, mass_conversion_factor=None, charge_conversion_factor=True, zero_dihedral_weighting_factor=False, moleculeID_offset=1)[source]

Output a LAMMPS data file.

Outputs a LAMMPS data file in the ‘full’ atom style format. Default units are ‘real’ units. See http://lammps.sandia.gov/doc/atom_style.html for more information on atom styles.

Parameters:
structureparmed.Structure

ParmEd structure object

filenamestr

Path of the output file

atom_style: str, optional, default=’full’

Defines the style of atoms to be saved in a LAMMPS data file. The following atom styles are currently supported: ‘full’, ‘atomic’, ‘charge’, ‘molecular’ see http://lammps.sandia.gov/doc/atom_style.html for more information on atom styles.

unit_style: str, optional, default=’real’

Defines to unit style to be save in a LAMMPS data file. Defaults to ‘real’ units. Current styles are supported: ‘real’, ‘lj’, ‘metal’ see https://lammps.sandia.gov/doc/99/units.html for more information on unit styles

minslist, optional, default=None

Minimum box dimension in x, y, z directions, nm

maxslist, optional, default=None

Maximum box dimension in x, y, z directions, nm

pair_coeff_labelstr, optional, default=None

Provide a custom label to the pair_coeffs section in the lammps data file. A value of None means a suitable default will be chosen.

detect_forcefield_stylebool, optional, default=True

If True, format lammpsdata parameters based on the contents of the parmed Structure

use_urey_bradleysbool, optional, default=False

If True, will treat angles as CHARMM-style angles with urey bradley terms while looking for structure.urey_bradleys

use_rb_torsionsbool, optional, default=True

If True, will treat dihedrals OPLS-style torsions while looking for structure.rb_torsions

use_dihedralsbool, optional, default=False

If True, will treat dihedrals as CHARMM-style dihedrals while looking for structure.dihedrals

zero_dihedral_weighting_factorbool, optional, default=False

If True, will set weighting parameter to zero in CHARMM-style dihedrals. This should be True if the CHARMM dihedral style is used in non-CHARMM forcefields.

sigma_conversion_factorfloat, optional, default=None

If unit_style is set to ‘lj’, then sigma conversion factor is used to non-dimensionalize. Assume to be in units of nm. If None, will take the largest sigma value in the structure.atoms.sigma values.

epsilon_conversion_factorfloat, optional, default=None

If unit_style is set to ‘lj’, then epsilon conversion factor is used to non-dimensionalize. Assume to be in units of kCal/mol. If None, will take the largest epsilon value in the structure.atoms.epsilon values.

mass_conversion_factorfloat, optional, default=None

If unit_style is set to ‘lj’, then mass conversion factor is used to non-dimensionalize. Assume to be in units of amu. If None, will take the largest mass value in the structure.atoms.mass values.

charge_conversion_factorbool, optional, default=True

If unit_style is set to ‘lj’, then charge conversion factor may or may not be used to non-dimensionalize. Assume to be in elementary charge units. If True, the charges are scaled by np.sqrt(4*np.pi()*eps_0*sigma_conversion_factor*epsilon_conversion_factor). If False, the charges are not scaled and the user must be wary to choose the dielectric constant properly, which may be more convenient to implement an implicit solvent.

moleculeID_offsetint , optional, default=1

Since LAMMPS treats the MoleculeID as an additional set of information to identify what molecule an atom belongs to, this currently behaves as a residue id. This value needs to start at 1 to be considered a real molecule.

Notes

See http://lammps.sandia.gov/doc/2001/data_format.html for a full description of the LAMMPS data format. Currently the following sections are supported (in addition to the header): Masses, Nonbond Coeffs, Bond Coeffs, Angle Coeffs, Dihedral Coeffs, Atoms, Bonds, Angles, Dihedrals, Impropers OPLS and CHARMM forcefield styles are supported, AMBER forcefield styles are NOT

Some of this function has beed adopted from mdtraj’s support of the LAMMPSTRJ trajectory format. See https://github.com/mdtraj/mdtraj/blob/master/mdtraj/formats/lammpstrj.py for details.

unique_typesa sorted list of unique atomtypes for all atoms in the structure.
Defined by:

atomtype : atom.type

unique_bond_types: an enumarated OrderedDict of unique bond types for all bonds in the structure.

Defined by bond parameters and component atomtypes, in order: — k : bond.type.k — req : bond.type.req — atomtypes : sorted((bond.atom1.type, bond.atom2.type))

unique_angle_types: an enumerated OrderedDict of unique angle types for all

angles in the structure. Defined by angle parameters and component atomtypes, in order: — k : angle.type.k — theteq : angle.type.theteq — vertex atomtype: angle.atom2.type — atomtypes: sorted((bond.atom1.type, bond.atom3.type))

unique_dihedral_types: an enumerated OrderedDict of unique dihedrals type

for all dihedrals in the structure. Defined by dihedral parameters and component atomtypes, in order: — c0 : dihedral.type.c0 — c1 : dihedral.type.c1 — c2 : dihedral.type.c2 — c3 : dihedral.type.c3 — c4 : dihedral.type.c4 — c5 : dihedral.type.c5 — scee : dihedral.type.scee — scnb : dihedral.type.scnb — atomtype 1 : dihedral.atom1.type — atomtype 2 : dihedral.atom2.type — atomtype 3 : dihedral.atom3.type — atomtype 4 : dihedral.atom4.type