ease4lmp package

Submodules

ease4lmp.bonded_atoms module

Submodule for a class extending ase.Atoms with bonds.

class ease4lmp.bonded_atoms.BondedAtoms(*args, **kwargs)[source]

Bases: ase.atoms.Atoms

Inheriting ase.Atoms class, and has functionalities to handle (covalent) bonds connecting a pair of atoms.

An instance of this class has bond data as a numpy.ndarray as with other atomic properties like positions and velocities. Shape of this array for bond data is (N, M, 4), where N is the number of atoms and M is the maximum number of bonds per atom. The first axis of the array corresponds to atoms, and the second axis corresponds to bonds connected to each atom. Four integers in the third axis describes each bond: the first integer is a relative index for the other atom, the second to fourth integers are relative image flags (used for resolving a periodic boundary condition) for the other atom in the x, y and z direction, respectively.

In addition, an instance of this class has type of atoms as a one-dimensional numpy.ndarray. Length of this array is equal to the number of atoms. Atom’s type is typically required for making Lammps’ data file and molecule file.

Basically, you can use methods of ase.Atoms in the same manner.

__delitem__(idx)[source]

Delete a selected atom.

This method overrides ase.Atoms.__delitem__. Before calling the parent method to remove a selected atom, this method deletes bond data associated with the atom.

Parameters:

idx: int or list/tuple of int
Index for an atom to be deleted. A list or tuple of integers is also acceptable for multiple atoms.

Examples:

>>> from ease4lmp import BondedAtoms
>>> atoms = BondedAtoms('CO', positions=[(0, 0, 0), (0, 0, 1.1)])
>>> atoms.set_bonds([
...   [[1, 0, 0, 0],[0, 0, 0, 0],[0, 0, 0, 0],[0, 0, 0, 0]],
...   [[-1, 0, 0, 0],[0, 0, 0, 0],[0, 0, 0, 0],[0, 0, 0, 0]]
... ])
>>> atoms.get_bonds()
array([[[ 1,  0,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0]],
       [[-1,  0,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0]]])
>>> del atoms[0]
>>> atoms.get_bonds()
array([[[0, 0, 0, 0],
        [0, 0, 0, 0],
        [0, 0, 0, 0],
        [0, 0, 0, 0]]])
__imul__(m)[source]

Create a super cell and returns it.

This method overrides ase.Atoms.__imul__. After calling the parent method to expand the simulation box, this method recalculates image flags over the expanded bond data. Assuming that two atoms connected with a bond are in different images of the periodic simulation box before expansion. If the two atoms are in the same image after expansion, their image falgs for the opponent must be updated to (0, 0, 0).

Parameters:

m: int or list/tuple of int
Three-element tuple or list defining how many times a periodic simulation box is repeated in each direction. The first, second and third element corresponds to the x, y and z direction, respectively. If an integer is given, a tuple (m, m, m) is used.

Examples:

>>> from ease4lmp import BondedAtoms
>>> wire = BondedAtoms(
...   'Au',
...   positions=[[0.0, 5.0, 5.0]],
...   cell=[2.9, 5.0, 5.0],
...   pbc=[1, 0, 0])
>>> wire.add_bond(0, 0, img2=(1, 0, 0))
>>> wire.get_bonds()
array([[[ 0,  1,  0,  0],
        [ 0, -1,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0]]])
>>> wire *= (3,1,1)
>>> wire.get_bonds()
array([[[ 1,  0,  0,  0],
        [ 2, -1,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0]],
       [[ 1,  0,  0,  0],
        [-1,  0,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0]],
       [[-2,  1,  0,  0],
        [-1,  0,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0]]])
__init__(*args, **kwargs)[source]

Parameters:

args:
A variable number of arguments passed to super().__init__().
kwargs:
A variable number of named arguments passed to super().__init__().
add_bond(atom1, atom2, img1=(0, 0, 0), img2=(0, 0, 0))[source]

Add a bond connecting two atoms.

Parameters:

atom1: int
Index for the first atom.
atom2: int
Index for the second atom.
img1: tuple/list of int
Three-element tuple or list of integers indicating which image of a periodic simulation box the first atom is in.
img2: tuple/list of int
Three-element tuple or list of integers indicating which image of a periodic simulation box the second atom is in.

Note

  • Bond data stored in the first atom: [atom2-atom1, img2[0]-img1[0], img2[1]-img1[1], img2[2]-img1[2]].
  • Bond data stored in the second atom: [atom1-atom2, img1[0]-img2[0], img1[1]-img2[1], img1[2]-img2[2]].
adjust_pbc_bonds()[source]

Adjust bonds under periodic boundary conditions.

This method corrects image flags of each bond if distance between two atoms connected with the bond is longer than a half length of the unit cell side.

Note

If the i th atom connects to one of images of the j th atom, it must not connect to another image of the the j th atom. In addition, all the bond lengths must be less than a half length of the unit cell side.

change_max_bonds(n=4)[source]

Change the maximum number of bonds per atom.

Like original properties of ase.Atoms, bond data is also stored as a numpy.ndarray, which is a homogeneous array of fixed-size items. This means that all the atoms have the same capacity for bond data, and the capacity is fixed until manually changed. This method changes the capacity, namely, the maximum number of bonds per atom. Note that existing data stay unchanged after changing the capacity.

Parameters:

n: int
New maximum number of bonds per atom.
extend(other)[source]

Extend this instance by appending the given instance.

This method returns this extended instance, which is the same behavior with ase.Atoms.extend().

Parameters:

other: instance of ase.Atom/ase.Atoms/BondedAtoms
For ase.Atom or ase.Atoms instance, the instance will be down-casted to BondedAtoms before extension.
get_bonded_angles()[source]

Return a two-dimensional numpy.ndarray describing all angles defined by two consecutive bonds.

Shape of the returned array is (N_a, 3), where N_a is the number of angles defined by two consecutive bonds. Each row (the first axis) corresponds to each angle, and each column (the second axis) consists of absolute index of the first, second and third atom in the each angle.

get_bonded_bonds()[source]

Return a two-dimensional numpy.ndarray describing all bonds.

Shape of the returned array is (N_b, 2), where N_b is the number of bonds. Each row (the first axis) corresponds to each bond, and each column (the second axis) consists of absolute index of the first and second atom connected by the each bond.

get_bonded_dihedrals()[source]

Return a two-dimensional numpy.ndarray describing all dihedrals defined by three consecutive bonds.

Shape of the returned array is (N_d, 4), where N_d is the number of dihedrals defined by three consecutive bonds. Each row (the first axis) corresponds to each dihedral, and each column (the second axis) consists of absolute index of the first, second, third and fourth atom in the each dihedral.

get_bonded_impropers()[source]

Return a two-dimensional numpy.ndarray describing all impropers defined by three bonds connected to the same atom.

Shape of the returned array is (N_i, 4), where N_i is the number of impropers defined by three bonds connected to the same atom. Each row (the first axis) corresponds to each improper, and each column (the second axis) consists of absolute index of the first, second, third and fourth atom in the each improper; the first atom is the center one.

get_bonds()[source]

Return an array of integers containing bond data.

An instance of ase.Atoms stores data in self.arrays; for example, atoms’ positions are accessed by key ‘positions’. For BondedAtoms instance, bond data is also stored in self.arrays with key of ‘bonds’. This method returns a copy of self.arrays["bonds"].

get_bonds_per_atom()[source]

Return a nested list of bond data.

The returned list is two-dimensional: each element of the first list corresponds to each atom, and the second list consists of absolute indices for atoms bonded with the each atom.

Note

This method is different from self.get_bonds() in terms of the following two points:

  • No data about image flags.
  • Contains data of existing bonds only.
get_types()[source]

Return an array of integers containing atom’s type.

An instance of ase.Atoms stores data in self.arrays; for example, atoms’ positions are accessed by key ‘positions’. For BondedAtoms instance, atom’s type data is also stored in self.arrays with key of ‘types’. This method returns a copy of self.arrays["types"].

remove_bond(atom1, atom2, img1=(0, 0, 0), img2=(0, 0, 0))[source]

Remove a bond connecting the given two atoms.

Parameters:

atom1: int
Index for the first atom.
atom2: int
Index for the second atom.
img1: tuple/list of int
Three-element tuple or list of integers indicating which image of a periodic simulation box the first atom is in.
img2: tuple/list of int
Three-element tuple or list of integers indicating which image of a periodic simulation box the second atom is in.
set_bonds(bonds)[source]

Store the given bond data.

This method sets bond data to self.array using self.set_array().

Parameters:

bonds: array-like object of int
Array describing bonds. Shape of the array must be (N, M, 4), where N is the number of atoms and M is the maximum number of bonds per atom. Four integers of innermost-axis are as follows: the first integer is a relative index for the other atom, the second to fourth integers are relative image flags for the other atom in the x, y and z direction, respectively.
set_types(types)[source]

Store the given type of atoms.

This method sets type of atoms to self.array using self.set_array().

Parameters:

types: array-like object of int
Per-atom array containing atom’s type. Length of the array must be equal to the number of atoms.
sort_bonds()[source]

Sort bonds in each atom.

Note that the order of atoms stays unchanged.

wrap(center=(0.5, 0.5, 0.5), pbc=None, eps=1e-07)[source]

Wrap positions to unit cell.

Parameters:

center: array-like object of float
The positons in fractional coordinates that the new positions will be nearest possible to.
pbc: bool or list/tuple of bool
For each axis in the unit cell decides whether the positions will be moved along this axis. By default, the boundary conditions of the Atoms object will be used.
eps: float
Small number to prevent slightly negative coordinates from being wrapped.
ease4lmp.bonded_atoms.create_atoms_from_data(path, atom_style, pbc=False, lammps_unit='real')[source]

Create a BondedAtoms instance from Lammps’ data file.

Parameters:

path: str
File path to Lammps’ data file.
atom_style: str
Specifies an atom style used in Lammps.
pbc: array-like object of bool
For each axis in the unit cell decides whether the positions will be moved along this axis.
lammps_unit: str
Specifies a unit used in Lammps.
ease4lmp.bonded_atoms.create_atoms_from_json(atoms, bonds=[], lammps_unit='real')[source]

Create a BondedAtoms instance from JSON object(s).

Parameters:

atoms: list of dict

JSON describing atoms. All the dictionary contained in the list must have the same keys. The following keys are required (some are optional):

  • xu (float) : In Lammps’ unit.
  • yu (float) : In Lammps’ unit.
  • zu (float) : In Lammps’ unit.
  • type (int) : Optional.
  • mass (float) : In Lammps’ unit. Optional.
  • vx (float) : In Lammps’ unit. Optional.
  • vy (float) : In Lammps’ unit. Optional.
  • vz (float) : In Lammps’ unit. Optional.
  • id (int) : Required when bonds is supplied.
bonds: list of dict

JSON describing bonds. All the dictionary contained in the list must have the same keys. The following keys are required:

  • atom1-id (int) : The first atom id.
  • atom2-id (int) : The second atom id.
lammps_unit: str
Specifies a unit used in Lammps.
ease4lmp.bonded_atoms.create_atoms_from_molecule(path, lammps_unit='real')[source]

Create a BondedAtoms instance from Lammps’ molecule file.

Parameters:

path: str
File path to Lammps’ molecule file.
lammps_unit: str
Specifies a unit used in Lammps.

ease4lmp.lammps_dataformats module

Submodule containing format data for each atom style used in Atoms and Velocities section of Lammps’ data file.

Pease see also the following for details:

ease4lmp.lammps_reader module

Submodule for functions to read Lammps’ data file.

ease4lmp.lammps_reader.read_angles(path)[source]

Read angles data from a Lammps’ data (or molecule) file.

Returned value is a JSON object (list of dict).

Parameters:

path: str
File path to Lammps’ data file (or molecule file).
ease4lmp.lammps_reader.read_atoms_from_data(path, atom_style, mass=True, velocity=True)[source]

Read atoms data from a Lammps’ data file.

Returned value is a JSON object (list of dict).

Parameters:

path: str
File path to Lammps’ data file.
atom_style: str
Specifies an atom style used in Lammps.
mass: bool
Whether to include mass data if Masses section extis.
velocity: bool
Whether to include velocity data if Velocities section extis.
ease4lmp.lammps_reader.read_atoms_from_molecule(path)[source]

Read atoms data from a Lammps’ molecule file.

Returned value is a JSON object (list of dict).

Parameters:

path: str
File path to Lammps’ molecule file.
ease4lmp.lammps_reader.read_bonds(path)[source]

Read bonds data from a Lammps’ data (or molecule) file.

Returned value is a JSON object (list of dict).

Parameters:

path: str
File path to Lammps’ data file (or molecule file).
ease4lmp.lammps_reader.read_box(path)[source]

Read side lengths of the simulation box.

Parameters:

path: str
File path to Lammps’ data file (or molecule file).
ease4lmp.lammps_reader.read_dihedrals(path)[source]

Read dihedrals data from a Lammps’ data (or molecule) file.

Returned value is a JSON object (list of dict).

Parameters:

path: str
File path to Lammps’ data file (or molecule file).
ease4lmp.lammps_reader.read_impropers(path)[source]

Read impropers data from a Lammps’ data (or molecule) file.

Returned value is a JSON object (list of dict).

Parameters:

path: str
File path to Lammps’ data file (or molecule file).

ease4lmp.lammps_sectionwriter module

Submodule for classes writing each section in Lammps’ data file.

class ease4lmp.lammps_sectionwriter.LammpsAngles(sequences, datanames)[source]

Bases: ease4lmp.lammps_sectionwriter.LammpsTopology

An interface to write Angles section in Lammps’ data file (or molecule file).

class ease4lmp.lammps_sectionwriter.LammpsAtoms(atom_style, types, positions, velocities=None, masses=None)[source]

Bases: object

An interface to write Atoms (and Velocities) section in Lammps’ data file (or molecule file).

__init__(atom_style, types, positions, velocities=None, masses=None)[source]

Parameters:

atom_style: str
Specify atom style of Lammps.
types: array-like object of int
Something like returned value of BondedAtoms.get_types().
positions: array-like object of float
Something like returned value of ase.Atoms.get_positions().
velocities: None or array-like object of float
Something like returned value of ase.Atoms.get_velocities(). If None, the Velocities section will not be written. Note that the i th atom in velocities must correspond to the i th atom in types.
masses: None or array-like object of float
A list whose length is equal to the number of atoms. Each element is a mass of each atom used in Lammps. Note that the i th atom in masses must correspond to the i th atom in types.
__weakref__

list of weak references to the object (if defined)

get_num()[source]

Return the number of atoms.

get_num_type()[source]

Return the number of atom types.

get_required_datanames(molecule=False)[source]

Return a set of required data names (keys).

Names of data are returned if the data is required to write Lammps’ data file and has not been set yet.

Parameters:

molecule: bool
Whether to return data names required for Lammps’ molecule file. This method returns data names for Lammps’ data file by default.
set_data(**kwargs)[source]

Store given data in self._data.

Parameters:

kwargs:
A variable number of named arguments. Keywords are data names, and arguments are data values.
set_masses(type2mass)[source]

Store masses in self._data.

Though the parameter type2mass is a dictionary mapping atom’s type to its mass, masses will be stored as a per-atom array.

Parameters:

type2mass: dict from int to float
Dictionary mapping atom’s type to its mass.
shift_positions(shift=(0.0, 0.0, 0.0))[source]

Shift atom positions stored in self._data.

Parameters:

shift: tuple/list of number
A tuple or list defining a Cartesian vector by which atom positions shift.
write_lines(path, velocity=False, mass=False, **kwargs)[source]

Write lines of Atoms and Velocities section.

Parameters:

path: str
File path to Lammps’ data file.
velocity: bool
Whether to write Velocities section or not.
mass: bool
Whether to write Masses section or not.
kwargs:
A variable number of named arguments (currently not in use).
write_lines_for_molecule(path, mass=False, **kwargs)[source]

Write lines of Coords, Types, Charges and Masses section to Lammps’ molecule file.

This method creates a new LammpsDataLines instance to write Coords and Types sections (and Charges section) in Lammps’ molecule file.

Parameters:

path: str
File path to Lammps’ molecule file.
mass: bool
Whether to write Masses section or not.
kwargs:
A variable number of named arguments (currently not in use).
class ease4lmp.lammps_sectionwriter.LammpsBonds(sequences, datanames)[source]

Bases: ease4lmp.lammps_sectionwriter.LammpsTopology

An interface to write Bonds section in Lammps’ data file (or molecule file).

get_maximum_per_atom()[source]

Return the maximum number of topology components per atom.

This method overrides LammpsTopology.get_maximum_per_atom() because data for a bond is linked to the first atom in the bond whereas data for the other topology components is linked to the second atom in the component. See also Lammps’ source code of Atom::data_bonds().

class ease4lmp.lammps_sectionwriter.LammpsDihedrals(sequences, datanames)[source]

Bases: ease4lmp.lammps_sectionwriter.LammpsTopology

An interface to write Dihedrals section in Lammps’ data file (or molecule file).

class ease4lmp.lammps_sectionwriter.LammpsImpropers(sequences, datanames, class2=False, **kwargs)[source]

Bases: ease4lmp.lammps_sectionwriter.LammpsTopology

An interface to write Impropers section in Lammps’ data file (or molecule file).

__init__(sequences, datanames, class2=False, **kwargs)[source]

Parameters:

sequences: array-like object of int
Two-dimensional array describing topology components. Note that the first atom of each sequence is a center atom of improper (all the three bonds connect to the atom).
datanames: tuple of str
Tuple of strings specifying data names written in a line.
class2: bool
Whether forcefiled is CLASS2 or not. If forcefield is CLASS2, a center atom of each improper must be the second atom in a sequence. For other forcefields, the first atom is the center one.
kwargs:
A variable number of named arguments (for receiving excessive arguments).
get_sequence_patterns(atom_types)[source]

Return a set of unique sequences of atom types appearing in the topology components.

This method overrides LammpsTopology.get_sequence_patterns() where original and reverse ordered sequence are considered as equivalent; it is not the case with improper.

Parameters:

atom_types: array-like object of int
A one dimensional array-like object whose i th element is a type of the i th atom.
class ease4lmp.lammps_sectionwriter.LammpsSpecialBonds(bonds_per_atom)[source]

Bases: object

An interface to write Special Bond Counts and Special Bonds section in Lammps’ molecule file.

__init__(bonds_per_atom)[source]

Parameters:

bonds_per_atom: array-like object of int
Two-dimensional list containing bond data for each atom: each element of the first list corresponds to each atom, and the second list consists of absolute indices for atoms bonded with the each atom (Returned value of BondedAtoms.get_bonds_per_atom()).
__weakref__

list of weak references to the object (if defined)

get_maximum_per_atom()[source]

Return the maximum number of special bonds per atom (sum over 1-2, 1-3, and 1-4 special bonds).

write_lines(path)[source]

Write lines of Special Bond Counts and Special Bonds section in Lammps’ molecule file.

Parameters:

path: str
File path to Lammps’ molecule file.
class ease4lmp.lammps_sectionwriter.LammpsTopology(sequences, datanames)[source]

Bases: object

An (abstract) interface to write a section of topology component (Bonds, Angles, Dihedrals, Impropers) in Lammps’ data file (or molecule file).

__init__(sequences, datanames)[source]

Parameters:

sequences: array-like object of int
Two-dimensional array describing topology components.
datanames: tuple of str
Tuple of strings specifying data names written in a line.
__weakref__

list of weak references to the object (if defined)

get_maximum_per_atom()[source]

Return the maximum number of topology components per atom.

For more details see also Lammps’ source code of Atom::data_angles(), Atom::data_dihedrals() and Atom::data_impropers().

get_num()[source]

Return the number of topology components.

The number of topology components is considered as 0 until types of the topology components are set.

get_num_type()[source]

Return the number of topology types.

get_sequence_patterns(atom_types)[source]

Return a set of unique sequences of atom types appearing in the topology components.

Parameters:

atom_types: array-like object of int
A one dimensional array-like object whose i th element is a type of the i th atom.
set_types(seq_to_type, atom_types)[source]

Set type of the topology components.

Parameters:

seq_to_type: dict from tuple to int
Mapping from sequences of atom types (tuple) to types of topology components consisting of those atoms (int). Note that the sequences can be arbitrary order; for example, two sequences of angle (1, 2, 3) and (3, 2, 1) have the same meaning.
write_lines(path)[source]

Write lines of a topology section in Lammps’ data file (or molecule file).

Parameters:

path: str
File path to Lammps’ data file (or molecule file).
ease4lmp.lammps_sectionwriter.create_topology(name, sequences, **kwargs)[source]

A factory function for LammpsTopology’s subclass.

Parameters:

name: str
Name of topology component: ‘bond’, ‘angle’, ‘dihedral’ or ‘improper’.
sequences: array-like object of int
Two-dimensional array describing of topology components. Shapes of the array for ‘bond’, ‘angle’, ‘dihedral’ and ‘improper’ are (N, 2), (N, 3), (N, 4) and (N, 4), respectively. Here N is the number of topology components. Elements of the array are zero-based indices for atoms.
kwargs:

A variable number of named arguments.

  • class2 (bool) : Set to True when using CLASS2 forcefield.

ease4lmp.lammps_units module

Submodule containing unit data for Lammps (currently only ‘real’ and ‘si’ are supported).

Pease see also the followings for details:

A redefinition of SI base units is scheduled on 2019/05/20. I wonder if Lammps changes its unit definitions.

  • Avogadro constant becomes exactly 6.02214076*10^23 [/mol]
  • Boltzmann constant becomes exactly 1.380649*10^−23 [J/K]
  • Planck constant becomes exactly 6.626070150*10^-34 [J*s]
  • Elementary charge becomes exactly 1.602176634*10^−19 [C]

ease4lmp.lammps_writer module

Submodule for an interface to write Lammps’ data file.

class ease4lmp.lammps_writer.ExtendedPrism(cell, pbc=(True, True, True), digits=10)[source]

Bases: ase.calculators.lammpsrun.Prism

Extended ase.calculators.lammpsrun.Prism class.

transform_to_lammps(vectors)[source]

Return transposed vectors.

This method is required to convert vectors from ASE’s cartesian coordinate system to Lammps’ skewed coordinate system.

Note that transposing occurs only if the simulation box is skewed.

Parameters:

vectors: numpy.ndarray
Positions or velocities of atoms.
class ease4lmp.lammps_writer.LammpsWriter(atoms, atom_style, lammps_unit='real', **kwargs)[source]

Bases: object

An interface to write Lammps’ data file (or molecule file).

__init__(atoms, atom_style, lammps_unit='real', **kwargs)[source]

Parameters:

atoms: instance of BondedAtoms or ase.Atoms
Atoms to be written to Lammps’ data file.
atom_style: str
Specifies an atom style used in Lammps.
lammps_unit: str
Specifies a unit used in Lammps.
kwargs:

A variable number of named arguments.

  • class2 (bool) : If Lammps’ data file written by this instance is used with CLASS2 forcefield, class2=True must be set.
__weakref__

list of weak references to the object (if defined)

get_angle_patterns()[source]

Return a set of unique sequences of atom types appearing in all the sequences describing angles.

get_bond_patterns()[source]

Return a set of unique sequences of atom types appearing in all the sequences describing bonds.

get_dihedral_patterns()[source]

Return a set of unique sequences of atom types appearing in all the sequences describing dihedrals.

get_improper_patterns()[source]

Return a set of unique sequences of atom types appearing in all the sequences describing impropers.

get_required_datanames()[source]

Return a set of required data names (keys).

Names of data are returned if the data is required to write Lammps’ data file and has not been set yet.

get_required_datanames_for_molecule()[source]

Return a set of required data names (keys).

Names of data are returned if the data is required to write Lammps’ molecule file and has not been set yet.

get_sequence_patterns(key)[source]

Return a set of unique sequences of atom types appearing in the given topology component.

Parameters:

key: str
Name of topology component: ‘bond’, ‘angle’, ‘dihedral’ or ‘improper’.
print_max_angles_per_atom()[source]

Return the maximum number of angles per atom.

print_max_bonds_per_atom()[source]

Return the maximum number of bonds per atom.

print_max_dihedrals_per_atom()[source]

Return the maximum number of dihedrals per atom.

print_max_impropers_per_atom()[source]

Return the maximum number of impropers per atom.

print_max_specials_per_atom()[source]

Return the maximum number of special bonds per atom.

print_maximum_per_atom(key)[source]

Return the maximum number of given topology components per atom.

Parameters:

key: str
Name of topology component: ‘bond’, ‘angle’, ‘dihedral’ or ‘improper’.
set_angle_types(seq_to_type)[source]

Set types of angles.

Parameters:

seq_to_type: dict form tuple to int
Mapping three-element tuples of atom types to corresponding types of angle. Note that second atom of each tuple should be at the center of angle.

Examples:

>>> from ease4lmp import BondedAtoms, LammpsWriter
>>> from ase.build import molecule
>>> methanol = BondedAtoms(molecule("CH3OH"))
>>> methanol.get_atomic_numbers()
array([6, 8, 1, 1, 1, 1])
>>> methanol.get_distance(1, 3)  # O-H bond
0.9700009076665858
>>> methanol.set_types([1, 2, 3, 4, 3, 3])
>>> bond_list = [(0, 1), (0, 2), (0, 4), (0, 5), (1, 3)]
>>> for t in bond_list:
...   methanol.add_bond(*t)
>>> methanol.set_cell([[10., 0., 0.], [0., 10., 0.], [0., 0., 10.]])
>>> methanol.center()
>>> writer = LammpsWriter(methanol, atom_style="molecular", special_bonds=True)
LammpsAtoms: 'id' have been set
LammpsAtoms: 'type' have been set
LammpsAtoms: 'x' have been set
LammpsAtoms: 'y' have been set
LammpsAtoms: 'z' have been set
LammpsAtoms: 'mass' have been set
LammpsBonds: 'id' have been set
LammpsBonds: 'atom1' have been set
LammpsBonds: 'atom2' have been set
LammpsAngles: 'id' have been set
LammpsAngles: 'atom1' have been set
LammpsAngles: 'atom2' have been set
LammpsAngles: 'atom3' have been set
LammpsDihedrals: 'id' have been set
LammpsDihedrals: 'atom1' have been set
LammpsDihedrals: 'atom2' have been set
LammpsDihedrals: 'atom3' have been set
LammpsDihedrals: 'atom4' have been set
LammpsImpropers: 'id' have been set
LammpsImpropers: 'atom1' have been set
LammpsImpropers: 'atom2' have been set
LammpsImpropers: 'atom3' have been set
LammpsImpropers: 'atom4' have been set
>>> writer.get_angle_patterns()
{(1, 2, 4), (2, 1, 3), (3, 1, 3)}
>>> writer.set_angle_types({
...   (1, 2, 4): 1,
...   (2, 1, 3): 2,
...   (3, 1, 3): 3,
... })
LammpsAngles: 'type' have been set
set_atom_data(**kwargs)[source]

Set data for Atoms (and Velocities) section.

Parameters:

kwargs:
A variable number of named arguments. Keywords are data names, and arguments are data values (list-like objects).
set_bond_types(seq_to_type)[source]

Set types of bonds.

Parameters:

seq_to_type: dict form tuple to int
Mapping two-element tuples of atom types to corresponding types of bond.

Examples:

>>> from ease4lmp import BondedAtoms, LammpsWriter
>>> from ase.build import molecule
>>> methanol = BondedAtoms(molecule("CH3OH"))
>>> methanol.get_atomic_numbers()
array([6, 8, 1, 1, 1, 1])
>>> methanol.get_distance(1, 3)  # O-H bond
0.9700009076665858
>>> methanol.set_types([1, 2, 3, 4, 3, 3])
>>> bond_list = [(0, 1), (0, 2), (0, 4), (0, 5), (1, 3)]
>>> for t in bond_list:
...   methanol.add_bond(*t)
>>> methanol.set_cell([[10., 0., 0.], [0., 10., 0.], [0., 0., 10.]])
>>> methanol.center()
>>> writer = LammpsWriter(methanol, atom_style="molecular", special_bonds=True)
LammpsAtoms: 'id' have been set
LammpsAtoms: 'type' have been set
LammpsAtoms: 'x' have been set
LammpsAtoms: 'y' have been set
LammpsAtoms: 'z' have been set
LammpsAtoms: 'mass' have been set
LammpsBonds: 'id' have been set
LammpsBonds: 'atom1' have been set
LammpsBonds: 'atom2' have been set
LammpsAngles: 'id' have been set
LammpsAngles: 'atom1' have been set
LammpsAngles: 'atom2' have been set
LammpsAngles: 'atom3' have been set
LammpsDihedrals: 'id' have been set
LammpsDihedrals: 'atom1' have been set
LammpsDihedrals: 'atom2' have been set
LammpsDihedrals: 'atom3' have been set
LammpsDihedrals: 'atom4' have been set
LammpsImpropers: 'id' have been set
LammpsImpropers: 'atom1' have been set
LammpsImpropers: 'atom2' have been set
LammpsImpropers: 'atom3' have been set
LammpsImpropers: 'atom4' have been set
>>> writer.get_bond_patterns()
{(1, 2), (1, 3), (2, 4)}
>>> writer.set_bond_types({
...   (1, 2): 1,
...   (1, 3): 2,
...   (2, 4): 3,
... })
LammpsBonds: 'type' have been set
set_dihedral_types(seq_to_type)[source]

Set types of dihedrals.

Parameters:

seq_to_type: dict form tuple to int
Mapping four-element tuples of atom types to corresponding types of dihedral. Note that the four atoms should be connected linearly by three bonds in that order.

Examples:

>>> from ease4lmp import BondedAtoms, LammpsWriter
>>> from ase.build import molecule
>>> methanol = BondedAtoms(molecule("CH3OH"))
>>> methanol.get_atomic_numbers()
array([6, 8, 1, 1, 1, 1])
>>> methanol.get_distance(1, 3)  # O-H bond
0.9700009076665858
>>> methanol.set_types([1, 2, 3, 4, 3, 3])
>>> bond_list = [(0, 1), (0, 2), (0, 4), (0, 5), (1, 3)]
>>> for t in bond_list:
...   methanol.add_bond(*t)
>>> methanol.set_cell([[10., 0., 0.], [0., 10., 0.], [0., 0., 10.]])
>>> methanol.center()
>>> writer = LammpsWriter(methanol, atom_style="molecular", special_bonds=True)
LammpsAtoms: 'id' have been set
LammpsAtoms: 'type' have been set
LammpsAtoms: 'x' have been set
LammpsAtoms: 'y' have been set
LammpsAtoms: 'z' have been set
LammpsAtoms: 'mass' have been set
LammpsBonds: 'id' have been set
LammpsBonds: 'atom1' have been set
LammpsBonds: 'atom2' have been set
LammpsAngles: 'id' have been set
LammpsAngles: 'atom1' have been set
LammpsAngles: 'atom2' have been set
LammpsAngles: 'atom3' have been set
LammpsDihedrals: 'id' have been set
LammpsDihedrals: 'atom1' have been set
LammpsDihedrals: 'atom2' have been set
LammpsDihedrals: 'atom3' have been set
LammpsDihedrals: 'atom4' have been set
LammpsImpropers: 'id' have been set
LammpsImpropers: 'atom1' have been set
LammpsImpropers: 'atom2' have been set
LammpsImpropers: 'atom3' have been set
LammpsImpropers: 'atom4' have been set
>>> writer.get_dihedral_patterns()
{(3, 1, 2, 4)}
>>> writer.set_dihedral_types({
...   (3, 1, 2, 4): 1,
... })
LammpsDihedrals: 'type' have been set
set_improper_types(seq_to_type)[source]

Set types of impropers.

Parameters:

seq_to_type: dict form tuple to int
Mapping four-element tuples of atom types to corresponding types of improper. Note that, if CLASS2 forcefield is used, the second atom of each tuple should be at the center of improper; the first atom should be at the center for other forcefields.

Examples:

>>> from ease4lmp import BondedAtoms, LammpsWriter
>>> from ase.build import molecule
>>> methanol = BondedAtoms(molecule("CH3OH"))
>>> methanol.get_atomic_numbers()
array([6, 8, 1, 1, 1, 1])
>>> methanol.get_distance(1, 3)  # O-H bond
0.9700009076665858
>>> methanol.set_types([1, 2, 3, 4, 3, 3])
>>> bond_list = [(0, 1), (0, 2), (0, 4), (0, 5), (1, 3)]
>>> for t in bond_list:
...   methanol.add_bond(*t)
>>> methanol.set_cell([[10., 0., 0.], [0., 10., 0.], [0., 0., 10.]])
>>> methanol.center()
>>> writer = LammpsWriter(methanol, atom_style="molecular", special_bonds=True)
LammpsAtoms: 'id' have been set
LammpsAtoms: 'type' have been set
LammpsAtoms: 'x' have been set
LammpsAtoms: 'y' have been set
LammpsAtoms: 'z' have been set
LammpsAtoms: 'mass' have been set
LammpsBonds: 'id' have been set
LammpsBonds: 'atom1' have been set
LammpsBonds: 'atom2' have been set
LammpsAngles: 'id' have been set
LammpsAngles: 'atom1' have been set
LammpsAngles: 'atom2' have been set
LammpsAngles: 'atom3' have been set
LammpsDihedrals: 'id' have been set
LammpsDihedrals: 'atom1' have been set
LammpsDihedrals: 'atom2' have been set
LammpsDihedrals: 'atom3' have been set
LammpsDihedrals: 'atom4' have been set
LammpsImpropers: 'id' have been set
LammpsImpropers: 'atom1' have been set
LammpsImpropers: 'atom2' have been set
LammpsImpropers: 'atom3' have been set
LammpsImpropers: 'atom4' have been set
>>> writer.get_improper_patterns()
{(1, 2, 3, 3), (1, 3, 3, 3)}
>>> writer.set_improper_types({
...   (1, 2, 3, 3): 1,
...   (1, 3, 3, 3): 2,
... })
LammpsImpropers: 'type' have been set
set_masses(type2mass)[source]

Set data for Masses section.

Parameters:

type2mass: dict
Dictionary mapping atom’s type (int) to its mass (float).
set_topology_types(**kwargs)[source]

Set types of the given topology components.

Parameters:

kwargs:
A variable number of named arguments. Keywords are names of topology component (‘bond’, ‘angle’, ‘dihedral’, ‘improper’), and arguments are dictionary mapping sequences of atom types to corresponding types of the topology component.
write_lammps_data(path, centering=False, **kwargs)[source]

Write Lammps’ data file.

Parameters:

path: str
File path to Lammps’ data file.
centering: bool
Whether to shift a simulation box so that its center position becomes (0, 0, 0).
kwargs:

A variable number of named arguments.

  • velocity (bool) : Whether to write Velocities section or not.
  • mass (bool) : Whether to write Masses section or not.
  • num_atom_type (bool) : Number of atom types; this overwrites the number of atom types stored in this instance. It is useful when extra atoms will be added to the simulation box later.
  • num_{name}_type (bool) : Number of types of topology component (name is one of bond, angle, dihedral and improper); this overwrites the number of the component types stored in this instance. It is useful when extra atoms will be added to the simulation box later.
write_lammps_molecule(path, special_bonds=True, **kwargs)[source]

Write Lammps’ molecule file.

Parameters:

path: str
File path to Lammps’ molecule file.
special_bonds: bool
Whether to write Special Bond Counts and Special Bonds* section.
kwargs:

A variable number of named arguments.

  • mass (bool) : Whether to write Masses section or not.

Module contents

Extension of Atomic Simulation Environment (ASE) for LAMMPS

class ease4lmp.BondedAtoms(*args, **kwargs)[source]

Bases: ase.atoms.Atoms

Inheriting ase.Atoms class, and has functionalities to handle (covalent) bonds connecting a pair of atoms.

An instance of this class has bond data as a numpy.ndarray as with other atomic properties like positions and velocities. Shape of this array for bond data is (N, M, 4), where N is the number of atoms and M is the maximum number of bonds per atom. The first axis of the array corresponds to atoms, and the second axis corresponds to bonds connected to each atom. Four integers in the third axis describes each bond: the first integer is a relative index for the other atom, the second to fourth integers are relative image flags (used for resolving a periodic boundary condition) for the other atom in the x, y and z direction, respectively.

In addition, an instance of this class has type of atoms as a one-dimensional numpy.ndarray. Length of this array is equal to the number of atoms. Atom’s type is typically required for making Lammps’ data file and molecule file.

Basically, you can use methods of ase.Atoms in the same manner.

__delitem__(idx)[source]

Delete a selected atom.

This method overrides ase.Atoms.__delitem__. Before calling the parent method to remove a selected atom, this method deletes bond data associated with the atom.

Parameters:

idx: int or list/tuple of int
Index for an atom to be deleted. A list or tuple of integers is also acceptable for multiple atoms.

Examples:

>>> from ease4lmp import BondedAtoms
>>> atoms = BondedAtoms('CO', positions=[(0, 0, 0), (0, 0, 1.1)])
>>> atoms.set_bonds([
...   [[1, 0, 0, 0],[0, 0, 0, 0],[0, 0, 0, 0],[0, 0, 0, 0]],
...   [[-1, 0, 0, 0],[0, 0, 0, 0],[0, 0, 0, 0],[0, 0, 0, 0]]
... ])
>>> atoms.get_bonds()
array([[[ 1,  0,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0]],
       [[-1,  0,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0]]])
>>> del atoms[0]
>>> atoms.get_bonds()
array([[[0, 0, 0, 0],
        [0, 0, 0, 0],
        [0, 0, 0, 0],
        [0, 0, 0, 0]]])
__imul__(m)[source]

Create a super cell and returns it.

This method overrides ase.Atoms.__imul__. After calling the parent method to expand the simulation box, this method recalculates image flags over the expanded bond data. Assuming that two atoms connected with a bond are in different images of the periodic simulation box before expansion. If the two atoms are in the same image after expansion, their image falgs for the opponent must be updated to (0, 0, 0).

Parameters:

m: int or list/tuple of int
Three-element tuple or list defining how many times a periodic simulation box is repeated in each direction. The first, second and third element corresponds to the x, y and z direction, respectively. If an integer is given, a tuple (m, m, m) is used.

Examples:

>>> from ease4lmp import BondedAtoms
>>> wire = BondedAtoms(
...   'Au',
...   positions=[[0.0, 5.0, 5.0]],
...   cell=[2.9, 5.0, 5.0],
...   pbc=[1, 0, 0])
>>> wire.add_bond(0, 0, img2=(1, 0, 0))
>>> wire.get_bonds()
array([[[ 0,  1,  0,  0],
        [ 0, -1,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0]]])
>>> wire *= (3,1,1)
>>> wire.get_bonds()
array([[[ 1,  0,  0,  0],
        [ 2, -1,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0]],
       [[ 1,  0,  0,  0],
        [-1,  0,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0]],
       [[-2,  1,  0,  0],
        [-1,  0,  0,  0],
        [ 0,  0,  0,  0],
        [ 0,  0,  0,  0]]])
__init__(*args, **kwargs)[source]

Parameters:

args:
A variable number of arguments passed to super().__init__().
kwargs:
A variable number of named arguments passed to super().__init__().
__module__ = 'ease4lmp.bonded_atoms'
add_bond(atom1, atom2, img1=(0, 0, 0), img2=(0, 0, 0))[source]

Add a bond connecting two atoms.

Parameters:

atom1: int
Index for the first atom.
atom2: int
Index for the second atom.
img1: tuple/list of int
Three-element tuple or list of integers indicating which image of a periodic simulation box the first atom is in.
img2: tuple/list of int
Three-element tuple or list of integers indicating which image of a periodic simulation box the second atom is in.

Note

  • Bond data stored in the first atom: [atom2-atom1, img2[0]-img1[0], img2[1]-img1[1], img2[2]-img1[2]].
  • Bond data stored in the second atom: [atom1-atom2, img1[0]-img2[0], img1[1]-img2[1], img1[2]-img2[2]].
adjust_pbc_bonds()[source]

Adjust bonds under periodic boundary conditions.

This method corrects image flags of each bond if distance between two atoms connected with the bond is longer than a half length of the unit cell side.

Note

If the i th atom connects to one of images of the j th atom, it must not connect to another image of the the j th atom. In addition, all the bond lengths must be less than a half length of the unit cell side.

change_max_bonds(n=4)[source]

Change the maximum number of bonds per atom.

Like original properties of ase.Atoms, bond data is also stored as a numpy.ndarray, which is a homogeneous array of fixed-size items. This means that all the atoms have the same capacity for bond data, and the capacity is fixed until manually changed. This method changes the capacity, namely, the maximum number of bonds per atom. Note that existing data stay unchanged after changing the capacity.

Parameters:

n: int
New maximum number of bonds per atom.
extend(other)[source]

Extend this instance by appending the given instance.

This method returns this extended instance, which is the same behavior with ase.Atoms.extend().

Parameters:

other: instance of ase.Atom/ase.Atoms/BondedAtoms
For ase.Atom or ase.Atoms instance, the instance will be down-casted to BondedAtoms before extension.
get_bonded_angles()[source]

Return a two-dimensional numpy.ndarray describing all angles defined by two consecutive bonds.

Shape of the returned array is (N_a, 3), where N_a is the number of angles defined by two consecutive bonds. Each row (the first axis) corresponds to each angle, and each column (the second axis) consists of absolute index of the first, second and third atom in the each angle.

get_bonded_bonds()[source]

Return a two-dimensional numpy.ndarray describing all bonds.

Shape of the returned array is (N_b, 2), where N_b is the number of bonds. Each row (the first axis) corresponds to each bond, and each column (the second axis) consists of absolute index of the first and second atom connected by the each bond.

get_bonded_dihedrals()[source]

Return a two-dimensional numpy.ndarray describing all dihedrals defined by three consecutive bonds.

Shape of the returned array is (N_d, 4), where N_d is the number of dihedrals defined by three consecutive bonds. Each row (the first axis) corresponds to each dihedral, and each column (the second axis) consists of absolute index of the first, second, third and fourth atom in the each dihedral.

get_bonded_impropers()[source]

Return a two-dimensional numpy.ndarray describing all impropers defined by three bonds connected to the same atom.

Shape of the returned array is (N_i, 4), where N_i is the number of impropers defined by three bonds connected to the same atom. Each row (the first axis) corresponds to each improper, and each column (the second axis) consists of absolute index of the first, second, third and fourth atom in the each improper; the first atom is the center one.

get_bonds()[source]

Return an array of integers containing bond data.

An instance of ase.Atoms stores data in self.arrays; for example, atoms’ positions are accessed by key ‘positions’. For BondedAtoms instance, bond data is also stored in self.arrays with key of ‘bonds’. This method returns a copy of self.arrays["bonds"].

get_bonds_per_atom()[source]

Return a nested list of bond data.

The returned list is two-dimensional: each element of the first list corresponds to each atom, and the second list consists of absolute indices for atoms bonded with the each atom.

Note

This method is different from self.get_bonds() in terms of the following two points:

  • No data about image flags.
  • Contains data of existing bonds only.
get_types()[source]

Return an array of integers containing atom’s type.

An instance of ase.Atoms stores data in self.arrays; for example, atoms’ positions are accessed by key ‘positions’. For BondedAtoms instance, atom’s type data is also stored in self.arrays with key of ‘types’. This method returns a copy of self.arrays["types"].

remove_bond(atom1, atom2, img1=(0, 0, 0), img2=(0, 0, 0))[source]

Remove a bond connecting the given two atoms.

Parameters:

atom1: int
Index for the first atom.
atom2: int
Index for the second atom.
img1: tuple/list of int
Three-element tuple or list of integers indicating which image of a periodic simulation box the first atom is in.
img2: tuple/list of int
Three-element tuple or list of integers indicating which image of a periodic simulation box the second atom is in.
set_bonds(bonds)[source]

Store the given bond data.

This method sets bond data to self.array using self.set_array().

Parameters:

bonds: array-like object of int
Array describing bonds. Shape of the array must be (N, M, 4), where N is the number of atoms and M is the maximum number of bonds per atom. Four integers of innermost-axis are as follows: the first integer is a relative index for the other atom, the second to fourth integers are relative image flags for the other atom in the x, y and z direction, respectively.
set_types(types)[source]

Store the given type of atoms.

This method sets type of atoms to self.array using self.set_array().

Parameters:

types: array-like object of int
Per-atom array containing atom’s type. Length of the array must be equal to the number of atoms.
sort_bonds()[source]

Sort bonds in each atom.

Note that the order of atoms stays unchanged.

wrap(center=(0.5, 0.5, 0.5), pbc=None, eps=1e-07)[source]

Wrap positions to unit cell.

Parameters:

center: array-like object of float
The positons in fractional coordinates that the new positions will be nearest possible to.
pbc: bool or list/tuple of bool
For each axis in the unit cell decides whether the positions will be moved along this axis. By default, the boundary conditions of the Atoms object will be used.
eps: float
Small number to prevent slightly negative coordinates from being wrapped.
ease4lmp.create_atoms_from_json(atoms, bonds=[], lammps_unit='real')[source]

Create a BondedAtoms instance from JSON object(s).

Parameters:

atoms: list of dict

JSON describing atoms. All the dictionary contained in the list must have the same keys. The following keys are required (some are optional):

  • xu (float) : In Lammps’ unit.
  • yu (float) : In Lammps’ unit.
  • zu (float) : In Lammps’ unit.
  • type (int) : Optional.
  • mass (float) : In Lammps’ unit. Optional.
  • vx (float) : In Lammps’ unit. Optional.
  • vy (float) : In Lammps’ unit. Optional.
  • vz (float) : In Lammps’ unit. Optional.
  • id (int) : Required when bonds is supplied.
bonds: list of dict

JSON describing bonds. All the dictionary contained in the list must have the same keys. The following keys are required:

  • atom1-id (int) : The first atom id.
  • atom2-id (int) : The second atom id.
lammps_unit: str
Specifies a unit used in Lammps.
ease4lmp.create_atoms_from_data(path, atom_style, pbc=False, lammps_unit='real')[source]

Create a BondedAtoms instance from Lammps’ data file.

Parameters:

path: str
File path to Lammps’ data file.
atom_style: str
Specifies an atom style used in Lammps.
pbc: array-like object of bool
For each axis in the unit cell decides whether the positions will be moved along this axis.
lammps_unit: str
Specifies a unit used in Lammps.
ease4lmp.create_atoms_from_molecule(path, lammps_unit='real')[source]

Create a BondedAtoms instance from Lammps’ molecule file.

Parameters:

path: str
File path to Lammps’ molecule file.
lammps_unit: str
Specifies a unit used in Lammps.
class ease4lmp.LammpsWriter(atoms, atom_style, lammps_unit='real', **kwargs)[source]

Bases: object

An interface to write Lammps’ data file (or molecule file).

__dict__ = mappingproxy({'__module__': 'ease4lmp.lammps_writer', '__doc__': "An interface to write Lammps' data file (or molecule file).", '__init__': <function LammpsWriter.__init__>, 'get_required_datanames': <function LammpsWriter.get_required_datanames>, 'get_required_datanames_for_molecule': <function LammpsWriter.get_required_datanames_for_molecule>, 'get_sequence_patterns': <function LammpsWriter.get_sequence_patterns>, 'get_bond_patterns': <function LammpsWriter.get_bond_patterns>, 'get_angle_patterns': <function LammpsWriter.get_angle_patterns>, 'get_dihedral_patterns': <function LammpsWriter.get_dihedral_patterns>, 'get_improper_patterns': <function LammpsWriter.get_improper_patterns>, 'print_maximum_per_atom': <function LammpsWriter.print_maximum_per_atom>, 'print_max_bonds_per_atom': <function LammpsWriter.print_max_bonds_per_atom>, 'print_max_angles_per_atom': <function LammpsWriter.print_max_angles_per_atom>, 'print_max_dihedrals_per_atom': <function LammpsWriter.print_max_dihedrals_per_atom>, 'print_max_impropers_per_atom': <function LammpsWriter.print_max_impropers_per_atom>, 'print_max_specials_per_atom': <function LammpsWriter.print_max_specials_per_atom>, 'set_atom_data': <function LammpsWriter.set_atom_data>, 'set_masses': <function LammpsWriter.set_masses>, 'set_topology_types': <function LammpsWriter.set_topology_types>, 'set_bond_types': <function LammpsWriter.set_bond_types>, 'set_angle_types': <function LammpsWriter.set_angle_types>, 'set_dihedral_types': <function LammpsWriter.set_dihedral_types>, 'set_improper_types': <function LammpsWriter.set_improper_types>, 'write_lammps_data': <function LammpsWriter.write_lammps_data>, 'write_lammps_molecule': <function LammpsWriter.write_lammps_molecule>, '__dict__': <attribute '__dict__' of 'LammpsWriter' objects>, '__weakref__': <attribute '__weakref__' of 'LammpsWriter' objects>})
__init__(atoms, atom_style, lammps_unit='real', **kwargs)[source]

Parameters:

atoms: instance of BondedAtoms or ase.Atoms
Atoms to be written to Lammps’ data file.
atom_style: str
Specifies an atom style used in Lammps.
lammps_unit: str
Specifies a unit used in Lammps.
kwargs:

A variable number of named arguments.

  • class2 (bool) : If Lammps’ data file written by this instance is used with CLASS2 forcefield, class2=True must be set.
__module__ = 'ease4lmp.lammps_writer'
__weakref__

list of weak references to the object (if defined)

get_angle_patterns()[source]

Return a set of unique sequences of atom types appearing in all the sequences describing angles.

get_bond_patterns()[source]

Return a set of unique sequences of atom types appearing in all the sequences describing bonds.

get_dihedral_patterns()[source]

Return a set of unique sequences of atom types appearing in all the sequences describing dihedrals.

get_improper_patterns()[source]

Return a set of unique sequences of atom types appearing in all the sequences describing impropers.

get_required_datanames()[source]

Return a set of required data names (keys).

Names of data are returned if the data is required to write Lammps’ data file and has not been set yet.

get_required_datanames_for_molecule()[source]

Return a set of required data names (keys).

Names of data are returned if the data is required to write Lammps’ molecule file and has not been set yet.

get_sequence_patterns(key)[source]

Return a set of unique sequences of atom types appearing in the given topology component.

Parameters:

key: str
Name of topology component: ‘bond’, ‘angle’, ‘dihedral’ or ‘improper’.
print_max_angles_per_atom()[source]

Return the maximum number of angles per atom.

print_max_bonds_per_atom()[source]

Return the maximum number of bonds per atom.

print_max_dihedrals_per_atom()[source]

Return the maximum number of dihedrals per atom.

print_max_impropers_per_atom()[source]

Return the maximum number of impropers per atom.

print_max_specials_per_atom()[source]

Return the maximum number of special bonds per atom.

print_maximum_per_atom(key)[source]

Return the maximum number of given topology components per atom.

Parameters:

key: str
Name of topology component: ‘bond’, ‘angle’, ‘dihedral’ or ‘improper’.
set_angle_types(seq_to_type)[source]

Set types of angles.

Parameters:

seq_to_type: dict form tuple to int
Mapping three-element tuples of atom types to corresponding types of angle. Note that second atom of each tuple should be at the center of angle.

Examples:

>>> from ease4lmp import BondedAtoms, LammpsWriter
>>> from ase.build import molecule
>>> methanol = BondedAtoms(molecule("CH3OH"))
>>> methanol.get_atomic_numbers()
array([6, 8, 1, 1, 1, 1])
>>> methanol.get_distance(1, 3)  # O-H bond
0.9700009076665858
>>> methanol.set_types([1, 2, 3, 4, 3, 3])
>>> bond_list = [(0, 1), (0, 2), (0, 4), (0, 5), (1, 3)]
>>> for t in bond_list:
...   methanol.add_bond(*t)
>>> methanol.set_cell([[10., 0., 0.], [0., 10., 0.], [0., 0., 10.]])
>>> methanol.center()
>>> writer = LammpsWriter(methanol, atom_style="molecular", special_bonds=True)
LammpsAtoms: 'id' have been set
LammpsAtoms: 'type' have been set
LammpsAtoms: 'x' have been set
LammpsAtoms: 'y' have been set
LammpsAtoms: 'z' have been set
LammpsAtoms: 'mass' have been set
LammpsBonds: 'id' have been set
LammpsBonds: 'atom1' have been set
LammpsBonds: 'atom2' have been set
LammpsAngles: 'id' have been set
LammpsAngles: 'atom1' have been set
LammpsAngles: 'atom2' have been set
LammpsAngles: 'atom3' have been set
LammpsDihedrals: 'id' have been set
LammpsDihedrals: 'atom1' have been set
LammpsDihedrals: 'atom2' have been set
LammpsDihedrals: 'atom3' have been set
LammpsDihedrals: 'atom4' have been set
LammpsImpropers: 'id' have been set
LammpsImpropers: 'atom1' have been set
LammpsImpropers: 'atom2' have been set
LammpsImpropers: 'atom3' have been set
LammpsImpropers: 'atom4' have been set
>>> writer.get_angle_patterns()
{(1, 2, 4), (2, 1, 3), (3, 1, 3)}
>>> writer.set_angle_types({
...   (1, 2, 4): 1,
...   (2, 1, 3): 2,
...   (3, 1, 3): 3,
... })
LammpsAngles: 'type' have been set
set_atom_data(**kwargs)[source]

Set data for Atoms (and Velocities) section.

Parameters:

kwargs:
A variable number of named arguments. Keywords are data names, and arguments are data values (list-like objects).
set_bond_types(seq_to_type)[source]

Set types of bonds.

Parameters:

seq_to_type: dict form tuple to int
Mapping two-element tuples of atom types to corresponding types of bond.

Examples:

>>> from ease4lmp import BondedAtoms, LammpsWriter
>>> from ase.build import molecule
>>> methanol = BondedAtoms(molecule("CH3OH"))
>>> methanol.get_atomic_numbers()
array([6, 8, 1, 1, 1, 1])
>>> methanol.get_distance(1, 3)  # O-H bond
0.9700009076665858
>>> methanol.set_types([1, 2, 3, 4, 3, 3])
>>> bond_list = [(0, 1), (0, 2), (0, 4), (0, 5), (1, 3)]
>>> for t in bond_list:
...   methanol.add_bond(*t)
>>> methanol.set_cell([[10., 0., 0.], [0., 10., 0.], [0., 0., 10.]])
>>> methanol.center()
>>> writer = LammpsWriter(methanol, atom_style="molecular", special_bonds=True)
LammpsAtoms: 'id' have been set
LammpsAtoms: 'type' have been set
LammpsAtoms: 'x' have been set
LammpsAtoms: 'y' have been set
LammpsAtoms: 'z' have been set
LammpsAtoms: 'mass' have been set
LammpsBonds: 'id' have been set
LammpsBonds: 'atom1' have been set
LammpsBonds: 'atom2' have been set
LammpsAngles: 'id' have been set
LammpsAngles: 'atom1' have been set
LammpsAngles: 'atom2' have been set
LammpsAngles: 'atom3' have been set
LammpsDihedrals: 'id' have been set
LammpsDihedrals: 'atom1' have been set
LammpsDihedrals: 'atom2' have been set
LammpsDihedrals: 'atom3' have been set
LammpsDihedrals: 'atom4' have been set
LammpsImpropers: 'id' have been set
LammpsImpropers: 'atom1' have been set
LammpsImpropers: 'atom2' have been set
LammpsImpropers: 'atom3' have been set
LammpsImpropers: 'atom4' have been set
>>> writer.get_bond_patterns()
{(1, 2), (1, 3), (2, 4)}
>>> writer.set_bond_types({
...   (1, 2): 1,
...   (1, 3): 2,
...   (2, 4): 3,
... })
LammpsBonds: 'type' have been set
set_dihedral_types(seq_to_type)[source]

Set types of dihedrals.

Parameters:

seq_to_type: dict form tuple to int
Mapping four-element tuples of atom types to corresponding types of dihedral. Note that the four atoms should be connected linearly by three bonds in that order.

Examples:

>>> from ease4lmp import BondedAtoms, LammpsWriter
>>> from ase.build import molecule
>>> methanol = BondedAtoms(molecule("CH3OH"))
>>> methanol.get_atomic_numbers()
array([6, 8, 1, 1, 1, 1])
>>> methanol.get_distance(1, 3)  # O-H bond
0.9700009076665858
>>> methanol.set_types([1, 2, 3, 4, 3, 3])
>>> bond_list = [(0, 1), (0, 2), (0, 4), (0, 5), (1, 3)]
>>> for t in bond_list:
...   methanol.add_bond(*t)
>>> methanol.set_cell([[10., 0., 0.], [0., 10., 0.], [0., 0., 10.]])
>>> methanol.center()
>>> writer = LammpsWriter(methanol, atom_style="molecular", special_bonds=True)
LammpsAtoms: 'id' have been set
LammpsAtoms: 'type' have been set
LammpsAtoms: 'x' have been set
LammpsAtoms: 'y' have been set
LammpsAtoms: 'z' have been set
LammpsAtoms: 'mass' have been set
LammpsBonds: 'id' have been set
LammpsBonds: 'atom1' have been set
LammpsBonds: 'atom2' have been set
LammpsAngles: 'id' have been set
LammpsAngles: 'atom1' have been set
LammpsAngles: 'atom2' have been set
LammpsAngles: 'atom3' have been set
LammpsDihedrals: 'id' have been set
LammpsDihedrals: 'atom1' have been set
LammpsDihedrals: 'atom2' have been set
LammpsDihedrals: 'atom3' have been set
LammpsDihedrals: 'atom4' have been set
LammpsImpropers: 'id' have been set
LammpsImpropers: 'atom1' have been set
LammpsImpropers: 'atom2' have been set
LammpsImpropers: 'atom3' have been set
LammpsImpropers: 'atom4' have been set
>>> writer.get_dihedral_patterns()
{(3, 1, 2, 4)}
>>> writer.set_dihedral_types({
...   (3, 1, 2, 4): 1,
... })
LammpsDihedrals: 'type' have been set
set_improper_types(seq_to_type)[source]

Set types of impropers.

Parameters:

seq_to_type: dict form tuple to int
Mapping four-element tuples of atom types to corresponding types of improper. Note that, if CLASS2 forcefield is used, the second atom of each tuple should be at the center of improper; the first atom should be at the center for other forcefields.

Examples:

>>> from ease4lmp import BondedAtoms, LammpsWriter
>>> from ase.build import molecule
>>> methanol = BondedAtoms(molecule("CH3OH"))
>>> methanol.get_atomic_numbers()
array([6, 8, 1, 1, 1, 1])
>>> methanol.get_distance(1, 3)  # O-H bond
0.9700009076665858
>>> methanol.set_types([1, 2, 3, 4, 3, 3])
>>> bond_list = [(0, 1), (0, 2), (0, 4), (0, 5), (1, 3)]
>>> for t in bond_list:
...   methanol.add_bond(*t)
>>> methanol.set_cell([[10., 0., 0.], [0., 10., 0.], [0., 0., 10.]])
>>> methanol.center()
>>> writer = LammpsWriter(methanol, atom_style="molecular", special_bonds=True)
LammpsAtoms: 'id' have been set
LammpsAtoms: 'type' have been set
LammpsAtoms: 'x' have been set
LammpsAtoms: 'y' have been set
LammpsAtoms: 'z' have been set
LammpsAtoms: 'mass' have been set
LammpsBonds: 'id' have been set
LammpsBonds: 'atom1' have been set
LammpsBonds: 'atom2' have been set
LammpsAngles: 'id' have been set
LammpsAngles: 'atom1' have been set
LammpsAngles: 'atom2' have been set
LammpsAngles: 'atom3' have been set
LammpsDihedrals: 'id' have been set
LammpsDihedrals: 'atom1' have been set
LammpsDihedrals: 'atom2' have been set
LammpsDihedrals: 'atom3' have been set
LammpsDihedrals: 'atom4' have been set
LammpsImpropers: 'id' have been set
LammpsImpropers: 'atom1' have been set
LammpsImpropers: 'atom2' have been set
LammpsImpropers: 'atom3' have been set
LammpsImpropers: 'atom4' have been set
>>> writer.get_improper_patterns()
{(1, 2, 3, 3), (1, 3, 3, 3)}
>>> writer.set_improper_types({
...   (1, 2, 3, 3): 1,
...   (1, 3, 3, 3): 2,
... })
LammpsImpropers: 'type' have been set
set_masses(type2mass)[source]

Set data for Masses section.

Parameters:

type2mass: dict
Dictionary mapping atom’s type (int) to its mass (float).
set_topology_types(**kwargs)[source]

Set types of the given topology components.

Parameters:

kwargs:
A variable number of named arguments. Keywords are names of topology component (‘bond’, ‘angle’, ‘dihedral’, ‘improper’), and arguments are dictionary mapping sequences of atom types to corresponding types of the topology component.
write_lammps_data(path, centering=False, **kwargs)[source]

Write Lammps’ data file.

Parameters:

path: str
File path to Lammps’ data file.
centering: bool
Whether to shift a simulation box so that its center position becomes (0, 0, 0).
kwargs:

A variable number of named arguments.

  • velocity (bool) : Whether to write Velocities section or not.
  • mass (bool) : Whether to write Masses section or not.
  • num_atom_type (bool) : Number of atom types; this overwrites the number of atom types stored in this instance. It is useful when extra atoms will be added to the simulation box later.
  • num_{name}_type (bool) : Number of types of topology component (name is one of bond, angle, dihedral and improper); this overwrites the number of the component types stored in this instance. It is useful when extra atoms will be added to the simulation box later.
write_lammps_molecule(path, special_bonds=True, **kwargs)[source]

Write Lammps’ molecule file.

Parameters:

path: str
File path to Lammps’ molecule file.
special_bonds: bool
Whether to write Special Bond Counts and Special Bonds* section.
kwargs:

A variable number of named arguments.

  • mass (bool) : Whether to write Masses section or not.
ease4lmp.read_bonds(path)[source]

Read bonds data from a Lammps’ data (or molecule) file.

Returned value is a JSON object (list of dict).

Parameters:

path: str
File path to Lammps’ data file (or molecule file).
ease4lmp.read_angles(path)[source]

Read angles data from a Lammps’ data (or molecule) file.

Returned value is a JSON object (list of dict).

Parameters:

path: str
File path to Lammps’ data file (or molecule file).
ease4lmp.read_dihedrals(path)[source]

Read dihedrals data from a Lammps’ data (or molecule) file.

Returned value is a JSON object (list of dict).

Parameters:

path: str
File path to Lammps’ data file (or molecule file).
ease4lmp.read_impropers(path)[source]

Read impropers data from a Lammps’ data (or molecule) file.

Returned value is a JSON object (list of dict).

Parameters:

path: str
File path to Lammps’ data file (or molecule file).
ease4lmp.read_atoms_from_data(path, atom_style, mass=True, velocity=True)[source]

Read atoms data from a Lammps’ data file.

Returned value is a JSON object (list of dict).

Parameters:

path: str
File path to Lammps’ data file.
atom_style: str
Specifies an atom style used in Lammps.
mass: bool
Whether to include mass data if Masses section extis.
velocity: bool
Whether to include velocity data if Velocities section extis.
ease4lmp.read_atoms_from_molecule(path)[source]

Read atoms data from a Lammps’ molecule file.

Returned value is a JSON object (list of dict).

Parameters:

path: str
File path to Lammps’ molecule file.