Skip to content

API

asext

Python package extends functions of ASE (Atomic Simulation Environment).

Developed and maintained by C.Thang Nguyen

Modules:

Attributes:

ASEXT_ROOT = Path(__file__).parent module-attribute

__author__ = 'C.Thang Nguyen' module-attribute

__contact__ = 'http://thangckt.github.io/email' module-attribute

cell

Classes:

  • AseCell
  • CellTransform

    Tranform the cell and atom properties from old_cell to new_cell orientations.

Functions:

AseCell(array: np.ndarray)

Bases: Cell

Methods:

lower_triangular_form() -> tuple[Cell, np.ndarray]

Rename original function Cell.standard_form(), see https://gitlab.com/ase/ase/-/blob/master/ase/cell.py?ref_type=heads#L333

upper_triangular_form() -> tuple[Cell, np.ndarray]

Rotate axes such that the unit cell is an upper triangular matrix.

CellTransform(old_cell: np.ndarray, new_cell: np.ndarray, pure_rotation: bool = True)

Tranform the cell and atom properties from old_cell to new_cell orientations.

The idea is compute a linear transformation that maps the old cell to the new cell. A = solve(old_cell, new_cell) = old_cell^(-1) new_cell

Generally, this linear transformation A can include rotation R + shear/reshape U (stretching and shearing), i.e., A = R * U.

Therefore, this transformation can be used in two ways: 1. Directly apply A that includes both rotation and shear/stretch. (should avoid using this, since it is not clear how to transform properties like stress/forces) 2. Extract only the rotation part R from A (using polar decomposition), and use it to rotate vectors/tensors, ignoring shear/reshape change. - Extract the closest pure rotation R from A (using polar decomposition) - Use that R to rotate positions, forces, stress, etc.

Parameters:

  • old_cell (ndarray) –

    3x3 matrix represent the old cell.

  • new_cell (ndarray) –

    3x3 matrix represent the new cell.

  • pure_rotation (bool, default: True ) –

    If True, only use the rotation part of the transformation. Defaults to True.

Note
  • np.linalg.solve(A, B) solves AX = B for X. May fail if A is singular (square matrix with a determinant of zero, det(A)=0).
  • Rotation matrix is derived from QR decomposition of the cell, following Prism class

Methods:

Attributes:

old_cell = np.asarray(old_cell, dtype=float) instance-attribute
new_cell = np.asarray(new_cell, dtype=float) instance-attribute
R = _polar_rotation(A) instance-attribute
vectors_forward(vec: np.ndarray) -> np.ndarray

Rotate vectors from the old_cell's orient to the new_cell's orient.

Parameters:

  • vec (ndarray) –

    Nx3 matrix represent the vector properties. (positions, forces, etc. each row is a vector)

Returns:

  • ndarray

    np.ndarray: Rotated vectors.

vectors_backward(vec: np.ndarray) -> np.ndarray

Rotate vectors back from the new_cell to the old_cell. Same as Prism.vector_to_ase

tensor_forward(tensor: np.ndarray) -> np.ndarray

Rotate the tensor from the old_cell's orient to the new_cell's orient. (T' = Rᵀ T R) rotates the tensor into the rotated coordinate system

Parameters:

  • tensor (ndarray) –

    3x3 matrix represent the tensor properties. (e.g., 3x3 stress tensor)

Returns:

  • ndarray

    np.ndarray: Transformed tensor.

tensor_backward(tensor: np.ndarray) -> np.ndarray

Rotate the tensor back from the new_cell to the old_cell. Same as Prism.tensor_to_ase (T = R T' Rᵀ) rotates the tensor back into the original coordinate system

make_upper_triangular_cell(atoms: Atoms, zero_tol: float = 1e-12) -> Atoms

Atoms with a box is an upper triangular matrix is a requirement to run NPT class in ASE. [[ ax, ay, az ] [ 0, by, bz ] [ 0, 0, cz ]]

make_lower_triangular_cell(atoms: Atoms, zero_tol: float = 1e-12) -> Atoms

Converts the cell matrix of atoms into a lower triangular, to be used in LAMMPS: [[ ax, 0, 0 ] [ bx, by, 0 ] [ cx, cy, cz ]]

make_triangular_cell_extxyz(extxyz_file: str, form: str = 'lower') -> None

Make the cell of atoms in extxyz file to be triangular. Args: extxyz_file (str): Path to the extxyz file. form (str): 'upper' or 'lower'. Defaults to 'lower'.

_polar_rotation(A: np.ndarray) -> np.ndarray

Return the closest proper rotation to matrix A (polar decomposition).

The purpose of this function is to get only the orientation difference, ignoring any shear/stretch.

Remind: Given a linear transformation A=old_cell^(-1) new_cell (carry old cell vectors into the new cell vectors), we can decompose it into a rotation R and a symmetric positive semi-definite matrix U (which represents stretch/shear) such that A = R * U. The rotation matrix R captures the pure rotational component of the transformation, while U captures the deformation (stretching and shearing) component.

rotate_struct_property(struct: Atoms, new_cell: np.ndarray, wrap: bool = False, custom_vector_props: list[str] | None = None, custom_tensor_props: list[str] | None = None) -> Atoms

Rotate atomic structure and its properties to match a new cell orientation.

Parameters:

  • struct (Atoms) –

    Atoms object.

  • new_cell (ndarray) –

    3x3 matrix represent the new cell.

  • wrap (bool, default: False ) –

    If True, wrap atoms into the new cell.

  • custom_vector_props (list, default: None ) –

    List of vector properties to rotate. This allows to set vector properties with custom names.

  • custom_tensor_props (list, default: None ) –

    List of tensor properties to rotate. This allows to set tensor properties with custom names.

Returns:

  • Atoms

    ase.Atoms: Atoms object with rotated properties.

Note
  • Important note: deepcopy(struct) copies the struct.calc object, but struct.copy() does not.

io

Modules:

lmpdata

Functions:

_lammps_data_to_ase_atoms(data, colnames, cell, celldisp, pbc=False, atomsobj=Atoms, order=True, specorder=None, prismobj=None, units='metal')

Extract positions and other per-atom parameters and create Atoms

:param data: per atom data :param colnames: index for data :param cell: cell dimensions :param celldisp: origin shift :param pbc: periodic boundaries :param atomsobj: function to create ase-Atoms object :param order: sort atoms by id. Might be faster to turn off. Disregarded in case id column is not given in file. :param specorder: list of species to map lammps types to ase-species (usually .dump files to not contain type to species mapping) :param prismobj: Coordinate transformation between lammps and ase :type prismobj: Prism :param units: lammps units for unit transformation between lammps and ase :returns: Atoms object :rtype: Atoms

Notes: - The original function in ase.io.lammpsrun.lammps_data_to_ase_atoms can not recover the atom types. This function is modified to save the atom types if type column is given in the LAMMPS dump file.

read_lammps_dump_text(file: str, index=-1, **kwargs)

Process cleartext lammps dumpfiles

:param fileobj: filestream providing the trajectory data :param index: integer or slice object (default: get the last timestep) :returns: list of Atoms objects :rtype: list

Notes: - This function is a modified version of ase.io.lammpsrun.read_lammps_dump_text to allow storing atom types if type column is given in the LAMMPS dump file.

write_lammps_data(file: str, atoms: Atoms, *, specorder: list = None, reduce_cell: bool = False, force_skew: bool = False, prismobj: Prism = None, write_image_flags: bool = False, masses: bool = False, velocities: bool = False, units: str = 'metal', bonds: bool = True, atom_style: str = 'atomic')

Write atomic structure data to a LAMMPS data file.

Parameters

fd : file|str File to which the output will be written. atoms : Atoms Atoms to be written. specorder : list[str], optional Chemical symbols in the order of LAMMPS atom types, by default None force_skew : bool, optional Force to write the cell as a triclinic <https://docs.lammps.org/Howto_triclinic.html>__ box, by default False reduce_cell : bool, optional Whether the cell shape is reduced or not, by default False prismobj : Prism|None, optional Prism, by default None write_image_flags : bool, default False If True, the image flags, i.e., in which images of the periodic simulation box the atoms are, are written. masses : bool, optional Whether the atomic masses are written or not, by default False velocities : bool, optional Whether the atomic velocities are written or not, by default False units : str, optional LAMMPS units <https://docs.lammps.org/units.html>, by default 'metal' bonds : bool, optional Whether the bonds are written or not. Bonds can only be written for atom_style='full', by default True atom_style : {'atomic', 'charge', 'full'}, optional LAMMPS atom style <https://docs.lammps.org/atom_style.html>, by default 'atomic'

Notes: - This function is a modified version of ase.io.lammpsdata.write_lammps_data to allow writing atom types based on atoms.arrays['type'] if it exists. Otherwise, the atom types are assigned based on the order of specorder or sorted chemical symbols.

_get_types(atoms: Atoms, species: list)
_get_symbols_by_types(atoms: Atoms)

readwrite

Functions:

  • read_extxyz

    Read extxyz file. The exited ase.io.read returns a single Atoms object if file contains only one frame. This function will return a list of Atoms object.

  • write_extxyz

    Write a list of Atoms object to an extxyz file. The exited ase.io.write function does not support writing file if the parent directory does not exist. This function will overcome this problem.

  • read_lmpdump

    Shortcut to ase.io.lammpsrun.read_lammps_dump function.

  • write_lmpdata

    Shortcut to ase.io.lammpsdata.write_lammps_data function.

  • extxyz2lmpdata

    Convert extxyz file to LAMMPS data file.

  • lmpdata2extxyz

    Convert LAMMPS data file to extxyz file.

  • lmpdump2extxyz

    Convert LAMMPS dump file to extxyz file.

read_extxyz(extxyz_file: str, index=':') -> list[Atoms]

Read extxyz file. The exited ase.io.read returns a single Atoms object if file contains only one frame. This function will return a list of Atoms object.

Parameters:

  • extxyz_file (str) –

    Path to the output file.

Returns:

  • list ( list[Atoms] ) –

    List of Atoms object.

Note
  • ase.io.read returns a single Atoms object or a list of Atoms object, depending on the index argument.
    • index=":" will always return a list.
    • index=0 or index=-1 will return a single Atoms object.
  • this function will always return a list of Atoms object, even index=0 or index=-1
write_extxyz(outfile: str, structs: list[Atoms]) -> None

Write a list of Atoms object to an extxyz file. The exited ase.io.write function does not support writing file if the parent directory does not exist. This function will overcome this problem.

Parameters:

  • structs (list) –

    List of Atoms object.

  • outfile (str) –

    Path to the output file.

read_lmpdump(lmpdump_file: str, index=-1, units='metal', **kwargs) -> list[Atoms]

Shortcut to ase.io.lammpsrun.read_lammps_dump function.

Parameters:

  • lmpdump_file (str) –

    Path to the LAMMPS dump file.

  • index (int | list[int], default: -1 ) –

    integer or slice object (default: get the last timestep)

  • order (bool) –

    sort atoms by id. Might be faster to turn off. Default: True

  • specorder (list[str]) –

    list of species to map lammps types to ase-species. Default: None

  • units (str, default: 'metal' ) –

    lammps units for unit transformation between lammps and ase

Returns:

  • list ( list[Atoms] ) –

    List of Atoms object.

write_lmpdata(file: str, atoms: Atoms, *, specorder: list = None, reduce_cell: bool = False, force_skew: bool = False, prismobj: Prism = None, write_image_flags: bool = False, masses: bool = True, velocities: bool = False, units: str = 'metal', bonds: bool = True, atom_style: str = 'atomic') -> None

Shortcut to ase.io.lammpsdata.write_lammps_data function.

Parameters:

  • file (str) –

    File to which the output will be written.

  • atoms (Atoms) –

    Atoms to be written.

  • specorder (list[str], default: None ) –

    Chemical symbols in the order of LAMMPS atom types, by default None

  • force_skew (bool, default: False ) –

    Force to write the cell as a triclinic <https://docs.lammps.org/Howto_triclinic.html> box, by default False

  • reduce_cell (bool, default: False ) –

    Whether the cell shape is reduced or not, by default False

  • prismobj (Prism | None, default: None ) –

    Prism, by default None

  • write_image_flags (bool, default: False ) –

    default False. If True, the image flags, i.e., in which images of the periodic simulation box the atoms are, are written.

  • masses (bool, default: True ) –

    Whether the atomic masses are written or not, by default True

  • velocities (bool, default: False ) –

    Whether the atomic velocities are written or not, by default False

  • units (str, default: 'metal' ) –

    LAMMPS units <https://docs.lammps.org/units.html>, by default 'metal'

  • bonds (bool, default: True ) –

    Whether the bonds are written or not. Bonds can only be written for atom_style='full', by default True

  • atom_style

    {'atomic', 'charge', 'full'}, optional. LAMMPS atom style <https://docs.lammps.org/atom_style.html>, by default 'atomic'

Returns:

  • None

    None

extxyz2lmpdata(extxyz_file: str, lmpdata_file: str, masses: bool = True, units: str = 'metal', atom_style: str = 'atomic', **kwargs) -> list[str]

Convert extxyz file to LAMMPS data file. Note: - need to save 'original_cell' to able to revert the original orientation of the crystal. - Use atoms.arrays['type'] to set atom types when convert from extxyz to lammpsdata file.

lmpdata2extxyz(lmpdata_file: str, extxyz_file: str, original_cell_file: str = None)

Convert LAMMPS data file to extxyz file.

lmpdump2extxyz(lmpdump_file: str, extxyz_file: str, index: int | slice = -1, original_cell_file: str = None, stress_file: str = None, lammps_units: str = 'metal')

Convert LAMMPS dump file to extxyz file.

Parameters:

  • lmpdump_file (str) –

    Path to the LAMMPS dump file.

  • extxyz_file (str) –

    Path to the output extxyz file.

  • original_cell_file (str, default: None ) –

    Path to the text file contains original_cell. It should a simple text file that can write/read with numpy. If not provided, try to find in the same directory as lmpdump_file with the extension .original_cell. Defaults to None.

  • stress_file (str, default: None ) –

    Path to the text file contains stress tensor. Defaults to None.

Restriction
  • Current ver: stress is mapped based on frame_index, it requires that frames in text stress file must be in the same "length and order" as in the LAMMPS dump file.
  • struct.info.get("timestep") is a new feature in ASE 3.25 ?

struct

Functions:

  • strain_struct

    Apply engineering strain to an ASE Atoms structure along lattice vectors a, b, c.

  • perturb_struct

    Perturb the atoms by random displacements. This method adds random displacements to the atomic positions. See more

  • slice_struct

    Slice structure into the first subcell by given numbers along a, b, c (cell vector) directions.

  • align_struct_min_pos

    Align min atoms position to the min cell corner (0,0,0)

  • set_vacuum

    This function sets vacuum along cell vectors a, b, c.

  • check_bad_box_extxyz

    Check structure in extxyz file whether it has bad box.

  • check_bad_box

    Check if a simulation box is "bad" based on given criteria.

strain_struct(struct_in: Atoms, strains: list = [0, 0, 0]) -> Atoms

Apply engineering strain to an ASE Atoms structure along lattice vectors a, b, c.

Parameters:

  • struct (Atoms) –

    ASE Atoms object.

  • strains (list[float], default: [0, 0, 0] ) –

    Engineering strains [ε_a, ε_b, ε_c]. New_length = old_length * (1 + ε).

Returns:

  • atoms ( Atoms ) –

    New strained structure with scaled cell and atom positions.

perturb_struct(struct: Atoms, std_disp: float) -> Atoms

Perturb the atoms by random displacements. This method adds random displacements to the atomic positions. See more

slice_struct(struct_in: Atoms, slice_num=(1, 1, 1), tol=1e-05) -> Atoms

Slice structure into the first subcell by given numbers along a, b, c (cell vector) directions.

align_struct_min_pos(struct: Atoms) -> Atoms

Align min atoms position to the min cell corner (0,0,0)

set_vacuum(struct_in: Atoms, distances: list = [0.0, 0.0, 0.0]) -> Atoms

This function sets vacuum along cell vectors a, b, c.

Parameters:

  • struct (Atoms) –

    ASE Atoms object to add vacuum.

  • distances (list, default: [0.0, 0.0, 0.0] ) –

    Distances to add along cell vectors a, b, c (not x, y, z dims in Cartersian axes). Must be a list of 3 floats.

Returns:

  • struct ( Atoms ) –

    A new Atoms object with an expanded cell and centered atoms.

Notes
  • atoms.center() sets vacuum on both sides of the cell along the specified axis. So the total vacuum is twice the input value. This function is different in that, it set total vacuum equal to the input value.

check_bad_box_extxyz(extxyz_file: str, criteria: dict = {'length_ratio': 100, 'wrap_ratio': 0.5, 'tilt_ratio': 0.5}) -> list[int]

Check structure in extxyz file whether it has bad box. Return: a file remarking the bad box frames.

check_bad_box(struct: Atoms, criteria: dict = {'length_ratio': 20, 'wrap_ratio': 0.5, 'tilt_ratio': 0.5}) -> bool

Check if a simulation box is "bad" based on given criteria.

Args:

struct : ase.Atoms Atoms object containing the atomic structure. criteria : dict A dictionary of criteria to check, which contains pairs of {'criteria_name': threshold_value}. Available criteria: - length_ratio: The ratio of the longest to the shortest cell vector. - Formula: max(|a|, |b|, |c|) / min(|a|, |b|, |c|) - Prevents highly elongated simulation boxes. - wrap_ratio: Checks if one cell vector component is excessively wrapped around another. - Formula: [b_x / a_x, c_y / b_y, c_x / a_x] - Prevents excessive skewing. - tilt_ratio: Measures tilting of cell vectors relative to their axes. - Formula: [b_x / b_y, c_y / c_z, c_x / c_z] - Avoids excessive tilting that may disrupt periodic boundaries.

Returns:

is_bad : bool True if the simulation box violates any of the given criteria, otherwise False.

Raises:

RuntimeError If an unknown criterion key is provided.