Force constant models#

ForceConstantPotential#

class hiphive.ForceConstantPotential(cs, parameters, metadata=None)[source]#

A finalized force constant model. Can produce force constants for any structure compatible with the structure for which the model was set up.

Parameters:
  • cs (ClusterSpace) – The cluster space the model is based upon.

  • parameters (numpy.ndarray) – The fitted parameters.

  • metadata (dict) – Meta data; will be pickled when object is written to file.

get_force_constants(atoms)[source]#

Return the force constants of a compatible structure.

Parameters:

atoms (Atoms) – Input structure.

Return type:

ForceConstantModel

Returns:

Force constants in format for the input structure.

property metadata: dict[str, Any]#

Metadata associated with force constant potential.

property orbit_data: List[Dict[str, Any]]#

List of dictionaries containing detailed information for each orbit, e.g., cluster radius and force constant.

property primitive_structure: Atoms#

Atomic structure.

print_tables()[source]#

Prints information concerning the underlying cluster space to stdout,including, e.g., the number of cluster, orbits, and parameters by order and number of bodies.

static read(filename)[source]#

Reads a ForceConstantPotential object from file.

Parameters:

filename (str) – Name of input file to load ForceConstantPotential from.

Returns:

The original object as stored in the file.

Return type:

ForceConstantPotential

property symprec#
write(filename)[source]#

Writes a ForceConstantPotential to file.

Parameters:

filename (str) – Name of file to write ForceConstantPotential to.

ForceConstantCalculator#

class hiphive.calculators.ForceConstantCalculator(fcs, max_disp=3.0)[source]#

This class provides an ASE calculator that can be used in conjunction with integrators and optimizers with the atomic simulation environment (ASE). To initialize an object of this class one must provide the ideal atomic configuration along with a compatible force constant model.

Parameters:
  • fcs (ForceConstants) – Force constants instance; must include the atomic configuration.

  • max_disp (float) – Maximum allowed displacement before calculator raises ValueError.

calculate(atoms=None, properties=['energy'], system_changes=['positions', 'numbers', 'cell', 'pbc', 'initial_charges', 'initial_magmoms'])[source]#

Do the calculation.

Return type:

None

properties: list of str

List of what needs to be calculated. Can be any combination of ‘energy’, ‘forces’, ‘stress’, ‘dipole’, ‘charges’, ‘magmom’ and ‘magmoms’.

system_changes: list of str

List of what has changed since last calculation. Can be any combination of these six: ‘positions’, ‘numbers’, ‘cell’, ‘pbc’, ‘initial_charges’ and ‘initial_magmoms’.

Subclasses need to implement this, but can ignore properties and system_changes if they want. Calculated properties should be inserted into results dictionary like shown in this dummy example:

self.results = {'energy': 0.0,
                'forces': np.zeros((len(atoms), 3)),
                'stress': np.zeros(6),
                'dipole': np.zeros(3),
                'charges': np.zeros(len(atoms)),
                'magmom': 0.0,
                'magmoms': np.zeros(len(atoms))}

The subclass implementation should first call this implementation to set the atoms attribute and create any missing directories.

compute_energy_and_forces()[source]#

Computes energy and forces.

Returns:

energy and forces

Return type:

float, list(list(float))

implemented_properties: List[str] = ['energy', 'forces']#

Properties calculator can handle (energy, forces, …)

ForceConstants#

This module provides functionality for storing and handling of force constants.

class hiphive.force_constants.ForceConstants(supercell)[source]#

Base class for force constants.

assert_acoustic_sum_rules(order=None, tol=1e-06)[source]#

Asserts that force constants obey acoustic sum rules.

Parameters:
  • order (int) – Specifies which order to check, if None all are checked.

  • tol (float) – Numeric tolerance for checking sum rules.

Raises:

AssertionError – If the acoustic sum rules are violated.

property clusters: list#

Sorted list of clusters

compute_gamma_frequencies()[source]#

Returns the Gamma frequencies in THz using the second-order force constants.

Return type:

ndarray

classmethod from_arrays(supercell, fc2_array=None, fc3_array=None)[source]#

Constructs FCs from numpy arrays. One or both of fc2_array and fc3_array must not be None.

Parameters:
  • supercell (Atoms) – Aupercell structure.

  • fc2_array (ndarray) – Second-order force constant in phonopy format, i.e., in shape (N, N, 3, 3).

  • fc3_array (ndarray) – Third-order force constant in phonopy format, i.e., in shape (N, N, N, 3, 3, 3).

classmethod from_dense_dict(fc_dict, supercell)[source]#

All permutations of clusters that are not zero must be listed, if label symmetries are fullfilled will return a SortedForceConstants.

Parameters:
  • fc_dict (dict) – Force constants with keys corresponding to clusters and values to the force constants.

  • supercell (Atoms) – Atomic configuration.

classmethod from_sparse_dict(fc_dict, supercell)[source]#

Assumes label symmetries, meaning only one cluster for each permuation should be included.

Parameters:
  • fc_dict (dict) – Force constants with keys corresponding to clusters and values to the force constants.

  • supercell (Atoms) – Atomic configuration.

get_fc_array(order, format='phonopy')[source]#

Returns force constants in array format for specified order.

Parameters:
  • order (int) – Force constants for this order will be returned.

  • format (str) – Specify which format (shape) the NumPy array should have, possible values are 'phonopy' and 'ase'.

Return type:

ndarray

Returns:

Array with shape (N,)*order + (3,)*order where N is the number of atoms.

get_fc_dict(order=None)[source]#

Returns force constant dictionary for one specific order.

The returned dict may be sparse or may be dense depending on the underlying force constants.

Parameters:

order (int) – Force constants returned for this order.

Return type:

dict

Returns:

Dictionary with keys corresponding to clusters and values to respective force constant.

property n_atoms: int#

Number of atoms.

print_force_constant(cluster)[source]#

Prints force constants for a cluster in a nice format.

Parameters:

cluster (Tuple[int]) – Sites belonging to the cluster.

Return type:

None

classmethod read(fname)[source]#

Reads a ForceConstants object from file.

Parameters:

fname (str) – Name of file from which to read.

classmethod read_phono3py(supercell, fname)[source]#

Reads force constants from a phono3py calculation.

Parameters:
  • supercell (Atoms) – Supercell structure (corresponding to SPOSCAR).

  • fname (str) – Name of third-order force constant file.

classmethod read_phonopy(supercell, fname, format=None)[source]#

Reads force constants from a phonopy calculation.

Parameters:
  • supercell (Atoms) – Supercell structure (corresponding to SPOSCAR).

  • fname (str) – Name of second-order force constant file.

  • format (str) – Format for second-order force constants; possible values: 'text', 'hdf5'.

classmethod read_shengBTE(supercell, fname, prim, symprec=1e-05)[source]#

Reads third order force constants from a shengBTE calculation. The shengBTE force constants will be mapped onto a supercell.

Parameters:
  • supercell (Atoms) – Supercell structure.

  • fname (str) – Name of third-order force constant file.

  • prim (Atoms) – Primitive configuration (must be equivalent to structure used in the shengBTE calculation).

  • symprec – Structural symmetry tolerance.

property supercell: Atoms#

Supercell associated with force constants.

write(fname)[source]#

Writes ForceConstants object to file.

Parameters:

fname (str) – Name of file to which to write.

Return type:

None

write_to_phono3py(fname)[source]#

Writes force constants in phono3py format.

Parameters:

fname (str) – Name of file to which to write third-order force constant.

Return type:

None

write_to_phonopy(fname, format=None)[source]#

Writes force constants in phonopy format.

Parameters:
  • fname (str) – Name of file to which to write second-order force constant.

  • format (str) – Format for second-order force constants; possible values: 'text', 'hdf5'.

Return type:

None

write_to_shengBTE(fname, prim, **kwargs)[source]#

Writes third-order force constants in shengBTE format.

Parameters:
  • fname (str) – Name of file to which to write third-order force constant.

  • prim (Atoms) – Primitive configuration (must be equivalent to structure used in the shengBTE calculation).

Return type:

None

class hiphive.force_constants.RawForceConstants(fc_dict, supercell)[source]#

Force constants without label symmetries.

Parameters:
  • fc_dict (dict) – Force constants with keys corresponding to clusters and values to the force constants, should contain all clusters with nonzero force constants.

  • supercell (Atoms) – Supercell structure.

property orders: List[int]#

Orders for which force constants exist.

class hiphive.force_constants.SortedForceConstants(fc_dict, supercell)[source]#

Force constants with label symmetries.

Parameters:
  • fc_dict (dict) – Force constants with keys corresponding to clusters and values to the force constants; should only contain sorted clusters.

  • supercell (Atoms) – Supercell structure.

property orders: List[int]#

Orders for which force constants exist.

write_to_GPUMD(fname_fc, fname_clusters, order, tol=1e-10)[source]#

Writes force constants of the specified order in GPUMD format.

Parameters:
  • fname_fc (str) – Name of file which contains the look-up force constants.

  • fname_clusters (str) – Name of file which contains the clusters and the fc lookup index.

  • order (int) – Force constants for this order will be written to file.

  • tol (float) – If the norm of a force constant is less than tol then it is not written. If two force-constants are within tol; they are considered equal.

hiphive.force_constants.array_to_dense_dict(fc_array, fc_tol=1e-10)[source]#

Constructs a dense dict from a FC array in phonopy format.

Force constants with norm smaller than fc_tol will be considered zero and not included in the output dictionary.

Parameters:
  • fc_array (ndarray) – Force constant array in phonopy format.

  • fc_tol (float) – Force constants below this value are considered zero and omitted from the output dict.

Return type:

dict

hiphive.force_constants.check_label_symmetries(fc_dict)[source]#

Checks label symmetries for dense FC dict.

Parameters:

fc_dict (dict) – Force constants with keys corresponding to clusters and values to the force constants.

Return type:

bool

hiphive.force_constants.dense_dict_to_sparse_dict(fc_dict)[source]#

Converts dense dict to sparse dict.

This does not check if label symmetry is True, but rather will just keep the sorted clusters and their force constants.

Parameters:

fc_dict (dict) – Force constants with keys corresponding to clusters and values to the force constants.

Return type:

dict

hiphive.force_constants.symbolize_force_constant(fc, tol=1e-10)[source]#

Carries out a symbolic symmetrization of a force constant tensor.

Parameters:
  • fc (ndarray) – Force constant tensor.

  • tol (float) – Tolerance used to decide whether two elements are identical.

Return type:

ndarray

Returns:

Symbolic representation of force constant matrix.

ForceConstantModel#

Contains the force constant modell (FCM) which handles cluster and force constant information of a super cell object.

class hiphive.force_constant_model.ForceConstantModel(atoms, cs)[source]#

Transfers a cluster space onto a super structure.

Contains the full description of all clusters and the force constants within a super cell with periodic boundary conditions.

Parameters:
  • atoms (Atoms) – Configuration to which the cluster space is to be applied.

  • cs (ClusterSpace) – A cluster space compatible with the structure of the atoms.

get_energy_fit_vector(fit_matrix, displacements)[source]#

Returns the energy fit vector. The vector represents the linear relation between the parameters and the energy of the structure.

Parameters:
  • fit_matrix (ndarray) – Fit matrix for the given supercell, can be obtained by get_fit_matrix().

  • displacements (ndarray) – Displacements of each atom in the supercell ((N, 3) array).

Return type:

ndarray

get_fcs_sensing(fcs, sparse=False)[source]#

Creates a fit matrix from force constants directly.

If the underlying cluster space can completely span the force constants the correct parameters should be recovered. The method assumes that the force constants obey crystal, lattice and label symmetries and will naively incorporate only one force constant per orbit into the sensing matrix.

The parameters can be extracted using, e.g., least squares from numpy:

parameters = np.linalg.lstsq(*fcm.get_fcs_sensing(fcs))[0]
Parameters:

fcs (ForceConstants) – Force constants that are compatible with the ideal structure of the model.

Return type:

Union[Tuple[ndarray, ndarray], Tuple[coo_matrix, ndarray]]

Returns:

A tuple comprising the sensing matrix and the flattened, irreducible force constants.

get_fit_matrix(displacements)[source]#

Returns the matrix used to fit the parameters. Represents the linear relation between the parameters and the forces.

Parameters:

displacements (ndarray) – Displacements of each atom in the supercell ((N, 3) array).

Return type:

ndarray

get_force_constants()[source]#

Returns the force constants of the super cell.

Return type:

SortedForceConstants

Returns:

Complete set of force constants.

get_forces(displacements)[source]#

Returns the forces in the system given displacements. The parameters of the model must be set to get any result.

Parameters:

displacements (ndarray) – Displacements of each atom in the supercell (N, 3 array).

Return type:

ndarray

property parameters: ndarray#

Parameters.

static read(f)[source]#

Reads a force constant model from file.

Parameters:

f (Union[str, IO]) – name of input file (str) or stream to load from (IO)

Returns:

Force constant model object as stored in the file.

Return type:

ForceConstantModel

write(f)[source]#

Writes a force constant model to file.

Parameters:

f (Union[str, IO]) – name of input file (str) or stream to write to (IO)

Return type:

None

Constraints#

hiphive.enforce_rotational_sum_rules(cs, parameters, sum_rules=None, alpha=1e-06, **ridge_kwargs)[source]

Enforces rotational sum rules by projecting parameters.

Note

The interface to this function might change in future releases.

Parameters:
  • cs (ClusterSpace) – The underlying cluster space.

  • parameters (ndarray) – Parameters to be constrained.

  • sum_rules (List[str]) – Type of sum rules to enforce; possible values: 'Huang', 'Born-Huang'.

  • alpha (float) – Hyperparameter of the ridge regression algorithm; keyword argument passed to the optimizer; larger values specify stronger regularization, i.e., less correction but higher stability [default: 1e-6].

  • ridge_kwargs (dict) – kwargs to be passed to sklearn Ridge.

Return type:

ndarray

Returns:

Constrained parameters.

Examples

The rotational sum rules can be enforced to the parameters before constructing a force constant potential as illustrated by the following snippet:

cs = ClusterSpace(reference_structure, cutoffs)
sc = StructureContainer(cs)
# add structures to structure container
opt = Optimizer(sc.get_fit_data())
opt.train()
new_params = enforce_rotational_sum_rules(cs, opt.parameters,
    sum_rules=['Huang', 'Born-Huang'])
fcp = ForceConstantPotential(cs, new_params)