ferrodispcalc package

Submodules

ferrodispcalc.compute module

ferrodispcalc.compute.calculate_averaged_structure(traj: list[Atoms], select: list[int] | slice | None = None) Atoms[source]
ferrodispcalc.compute.calculate_displacement(traj: list[Atoms] | Atoms, nl: ndarray, select: list[int] | slice | None = None) ndarray[source]
ferrodispcalc.compute.calculate_octahedral_tilt(traj: list[Atoms] | Atoms, nl_bo: ndarray, select: list[int] | slice | None = None) ndarray[source]
ferrodispcalc.compute.calculate_polarization(traj, nl_ba: ndarray, nl_bx: ndarray, born_effective_charge: dict, select: list[int] | slice | None = None) ndarray[source]

ferrodispcalc.config module

Used for common configuration, including type map and born effective charges.

Author: Denan Li Last modified: 2024-07-15 Email: lidenan@westlake.edu.cn

class ferrodispcalc.config.TypeMap(UniPero: Tuple[str] = ('Ba', 'Pb', 'Ca', 'Sr', 'Bi', 'K', 'Na', 'Hf', 'Ti', 'Zr', 'Nb', 'Mg', 'In', 'Zn', 'O'), PSTO: Tuple[str] = ('Sr', 'Pb', 'Ti', 'O'))[source]

Bases: object

PSTO: Tuple[str] = ('Sr', 'Pb', 'Ti', 'O')
UniPero: Tuple[str] = ('Ba', 'Pb', 'Ca', 'Sr', 'Bi', 'K', 'Na', 'Hf', 'Ti', 'Zr', 'Nb', 'Mg', 'In', 'Zn', 'O')

ferrodispcalc.io module

ferrodispcalc.io.read_lammps_data(filename: str, type_map: list[str]) Atoms[source]
ferrodispcalc.io.read_lammps_dump(filename: str, type_map: List[str], select: int | slice | List[int] | None = None, cache: bool = True) Atoms | List[Atoms] | None[source]

Place holder for docstring Note: We must notice you that:

THE COORDINATE IS NOT SHIFTED BACK TO ORIGINAL POSITION. THUS A CONSTANT SHIFT MAY EXIST IN THE COORDINATES COMPARED TO OTHER TOOLS LIKE DPDATA. However, the shift is constant, thus not important

ferrodispcalc.io.read_lammps_log(file_name: str) Dict[source]
ferrodispcalc.io.read_xyz(filename: str, select_frames: list = None, cache: bool = True) list[Atoms][source]

Read an XYZ trajectory file and return a list of ASE Atoms objects.

Parameters: filename (str): Path to the XYZ file. select_frames (list, optional): List of frame indices to read. If None, read all frames. Defaults to None. cache (bool, optional): Whether to cache the read frames in memory. Defaults to True.

ferrodispcalc.neighborlist module

ferrodispcalc.neighborlist.build_neighbor_list(atoms: Atoms, center_elements: list[str], neighbor_elements: list[str], cutoff: float, neighbor_num: int, defect: bool = False) ndarray[source]

Build neighbor lists for selected atom types.

Find the neighbor of selected center atoms within a cutoff radius. And return the 1-based neighbor list array.

Parameters:
  • atoms (ase.Atoms) – ASE Atoms object containing atomic positions, cell and element types.

  • center_elements (list[str]) – Element symbols to treat as centers (e.g. [“Ti”, “Ba”]).

  • neighbor_elements (list[str]) – Element symbols considered as neighbors (e.g. [“O”]).

  • cutoff (float) – Cutoff radius used to search for neighbors (Unit: Angstrom).

  • neighbor_num (int) – Number of neighbors to keep per center. If fewer neighbors are found and defect is False a ValueError is raised; if defect is True missing neighbor slots are filled with the center index.

  • defect (bool, optional) – When True, allow centers with fewer than neighbor_num neighbors and fill missing entries with the center index (default: False).

Returns:

Integer array of shape (n_centers, neighbor_num + 1). Each row contains the center atom index followed by the neighbor atom indices. Indices are 1-based.

Return type:

np.ndarray

Raises:

ValueError – If any center atom has fewer than neighbor_num neighbors and defect is False.

Notes

  • Internally the function uses pymatgen.core.Structure.get_neighbor_list().

  • A small-cell check prints a warning if any cell vector length is less than 4.0 Å.

  • Ensure that the provided atoms object has periodic boundary conditions (PBC) and a valid cell when you expect PBC behavior.

Examples

from ase.io import read from ferrodispcalc.neighborlist import build_neighbor_list

atoms = read(“POSCAR”) nl = build_neighbor_list(

atoms, center_elements=[‘Pb’, ‘Sr’], neighbor_elements=[‘O’], cutoff=4.0, neighbor_num=12, defect=False

)

ferrodispcalc.neighborlist.save_neighbor_list(nl: ndarray, file_name: str, zero_based: bool = False) None[source]

Save a neighbor list array to a text file.

Parameters:
  • nl (np.ndarray) – Neighbor list array as returned by build_neighbor_list(). The function expects an integer array where each row contains a center index followed by neighbor indices.

  • file_name (str) – Output file path. The neighbor list will be written as a text file with fixed-width integer columns.

  • zero_based (bool, optional) – If True, convert indices to zero-based before saving. By default the neighbor list uses 1-based indices and zero_based is False.

Examples

>>> from ferrodispcalc.neighborlist import save_neighbor_list, build_neighbor_list
>>> nl = build_neighbor_list(...)
>>> save_neighbor_list(nl, 'nl.dat', zero_based=False)

ferrodispcalc.utils module

Module contents