ForceConstants API Reference

The ForceConstants class is the basis for all of the kALDo calculations. It holds or calculates the interatomic force constants for the system.

Attributes

class kaldo.forceconstants.ForceConstants(atoms, supercell=(1, 1, 1), third_supercell=None, folder='displacement', distance_threshold=None)[source]

Class for constructing the finite difference object to calculate the second/third order force constant matrices after providing the unit cell geometry and calculator information.

Parameters:
  • atoms (Tabulated xyz files or ASE Atoms object) – The atoms to work on.

  • supercell – Size of supercell given by the number of repetitions (l, m, n) of the small unit cell in each direction. Defaults to (1, 1, 1)

  • third_supercell (tuple, optional) – Same as supercell, but for the third order force constant matrix. If not provided, it’s copied from supercell. Defaults to self.supercell

  • folder (str, optional) – Name to be used for the displacement information folder. Defaults to ‘displacement’

  • distance_threshold (float, optional) – If the distance between two atoms exceeds threshold, the interatomic force is ignored. Defaults to None

Methods

from_folder(folder[, supercell, format, ...])

Create a finite difference object from a folder

unfold_third_order([reduced_third, ...])

This method extrapolates a third order force constant matrix from a unit cell into a matrix for a larger supercell.

classmethod from_folder(folder, supercell=(1, 1, 1), format='numpy', third_energy_threshold=0.0, third_supercell=None, is_acoustic_sum=False, only_second=False, distance_threshold=None)[source]

Create a finite difference object from a folder

Parameters:
  • folder (str) – Chosen folder to load in system information.

  • supercell ((int, int, int), optional) – Number of unit cells in each cartesian direction replicated to form the input structure. Default is (1, 1, 1)

  • format ('numpy', 'eskm', 'lammps', 'shengbte', 'shengbte-qe', 'hiphive') – Format of force constant information being loaded into ForceConstants object. Default is ‘numpy’

  • third_energy_threshold (float, optional) –

    When importing sparse third order force constant matrices, energies below the threshold value in magnitude are ignored. Units: ev/A^3

    Default is None

  • distance_threshold (float, optional) – When calculating force constants, contributions from atoms further than the distance threshold will be ignored.

  • third_supercell ((int, int, int), optional) – Takes in the unit cell for the third order force constant matrix. Default is self.supercell

  • is_acoustic_sum (Bool, optional) – If true, the accoustic sum rule is applied to the dynamical matrix. Default is False

  • Inputs

  • ------

  • numpy (replicated_atoms.xyz, second.npy, third.npz) –

  • eskm (CONFIG, replicated_atoms.xyz, Dyn.form, THIRD) –

  • lammps (replicated_atoms.xyz, Dyn.form, THIRD) –

  • shengbte (CONTROL, POSCAR, FORCE_CONSTANTS_2ND/FORCE_CONSTANTS, FORCE_CONSTANTS_3RD) –

  • shengbte-qe (CONTROL, POSCAR, espresso.ifc2, FORCE_CONSTANTS_3RD) –

  • hiphive (atom_prim.xyz, replicated_atoms.xyz, model2.fcs, model3.fcs) –

Return type:

ForceConstants object

unfold_third_order(reduced_third=None, distance_threshold=None)[source]

This method extrapolates a third order force constant matrix from a unit cell into a matrix for a larger supercell.

Parameters:
  • reduced_third (array, optional) – The third order force constant matrix. Default is self.third

  • distance_threshold (float, optional) – When calculating force constants, contributions from atoms further than the distance threshold will be ignored. Default is self.distance_threshold