Other functions

Self-consistent phonons

hiphive.self_consistent_phonons.self_consistent_harmonic_model(atoms_ideal, calc, cs, T, alpha, n_iterations, n_structures, QM_statistics=False, parameters_start=None, imag_freq_factor=1.0, fit_kwargs={})[source]

Constructs a set of self-consistent second-order force constants that provides the closest match to the potential energy surface at a the specified temperature.

Parameters
  • atoms_ideal (ase.Atoms) – ideal structure

  • calc (ASE calculator object) – calculator to be used as reference potential

  • cs (ClusterSpace) – clusterspace onto which to project the reference potential

  • T (float) – temperature in K

  • alpha (float) – stepsize in optimization algorithm

  • n_iterations (int) – number of iterations in poor mans

  • n_structures (int) – number of structures to use when fitting

  • QM_statistics (bool) – if the amplitude of the quantum harmonic oscillator shoule be used instead of the classical amplitude in the phonon rattler

  • parameters_start (numpy.ndarray) – parameters from which to start the optimization

  • image_freq_factor (float) – If a squared frequency, w2, is negative then it is set to w2 = imag_freq_factor * np.abs(w2)

  • fit_kwargs (dict) – kwargs to be used in the fitting process (via Optimizer)

Returns

sequence of parameter vectors generated while iterating to self-consistency

Return type

list(numpy.ndarray)

MD tools

hiphive.md_tools.spectral_energy_density.compute_sed(traj, ideal, prim, k_points)[source]

Computes spectral energy density for a trajectory.

Parameters
  • traj (list) – trajectory with atoms objects with velocities

  • ideal (ASE atoms object) – ideal atoms object

  • prim (ASE atoms object) – compatible primitive cell. Must be aligned correctly

  • k_points (list) – list of k points in cart coord (2pi must be included)

Enforcing rotational sum rules

hiphive.core.rotational_constraints.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 (Optional[List[str]]) – type of sum rules to enforce; possible values: ‘Huang’, ‘Born-Huang’

  • alpha (float) – hyperparameter to 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

Returns

constrained parameters

Return type

numpy.ndarray

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)

Utilities

This module contains various support/utility functions.

class hiphive.utilities.Shell(types, distance, count=0)[source]

Neighbor Shell class

Parameters
  • types (list or tuple) – atomic types for neighbor shell

  • distance (float) – interatomic distance for neighbor shell

  • count (int) – number of pairs in the neighbor shell

hiphive.utilities.extract_parameters(fcs, cs, sanity_check=True, lstsq_method='numpy')[source]

Extracts parameters from force constants.

This function can be used to extract parameters to create a ForceConstantPotential from a known set of force constants. The return values come from NumPy’s lstsq function or from SciPy’s sparse lsqr function. Using lstsq_method=’scipy’ might be faster and have a smaller memory footprint for large systems, at the expense of some accuracy. This is due to the use of sparse matrices and an iterative solver.

Parameters
  • fcs (ForceConstants) – force constants

  • cs (ClusterSpace) – cluster space

  • sanity_check (bool) – bool whether or not to perform a sanity check by computing the relative error between the input fcs and the output fcs

  • lstsq_method (str) – method to use when making a least squares fit of a ForceConstantModel to the given fcs, allowed values are ‘numpy’ for np.linalg.lstsq or ‘scipy’ for scipy.sparse.linalg.lsqr

Returns

parameters that together with the ClusterSpace generates the best representation of the FCs

Return type

parameters

hiphive.utilities.find_permutation(atoms, atoms_ref)[source]

Returns the best permutation of atoms for mapping one configuration onto another.

Parameters
  • atoms (Atoms) – configuration to be permuted

  • atoms_ref (Atoms) – configuration onto which to map

Examples

After obtaining the permutation via p = find_permutation(atoms1, atoms2) the reordered structure atoms1[p] will give the closest match to atoms2.

Return type

List[int]

hiphive.utilities.get_displacements(atoms, atoms_ideal, cell_tol=0.0001)[source]

Returns the the smallest possible displacements between a displaced configuration relative to an ideal (reference) configuration.

Notes

Parameters
  • atoms (Atoms) – configuration with displaced atoms

  • atoms_ideal (Atoms) – ideal configuration relative to which displacements are computed

  • cell_tol (float) – cell tolerance; if cell missmatch more than tol value error is raised

Return type

ndarray

hiphive.utilities.get_neighbor_shells(atoms, cutoff, dist_tol=1e-05)[source]

Returns a list of neighbor shells.

Distances are grouped into shells via the following algorithm:

  1. Find smallest atomic distance d_min

  2. Find all pair distances in the range d_min + 1 * dist_tol

  3. Construct a shell from these and pop them from distance list

  4. Go to 1.

Parameters
  • atoms (Atoms) – configuration used for finding shells

  • cutoff (float) – exclude neighbor shells which have a distance larger than this value

  • dist_tol (float) – distance tolerance

Return type

List[Shell]

hiphive.utilities.prepare_structure(atoms, atoms_ideal, calc=None, check_permutation=True)[source]

Prepare a structure in the format suitable for a StructureContainer.

Either forces should be attached to input atoms object as an array, or the atoms object should have a SinglePointCalculator attached to it containing forces, or a calculator (calc) should be supplied.

Parameters
  • atoms (Atoms) – input structure

  • atoms_ideal (Atoms) – reference structure relative to which displacements are computed

  • check_permutation (bool) – whether find_permutation should be used or not

  • calc (Optional[SinglePointCalculator]) – ASE calculator used for computing forces

Returns

prepared ASE atoms object with forces and displacements as arrays

Return type

ASE atoms object

hiphive.utilities.prepare_structures(structures, atoms_ideal, calc=None, check_permutation=True)[source]

Prepares a set of structures in the format suitable for adding them to a StructureContainer.

structures should represent a list of supercells with displacements while atoms_ideal should provide the ideal reference structure (without displacements) for the given structures.

The structures that are returned will have their positions reset to the ideal structures. Displacements and forces will be added as arrays to the atoms objects.

If no calculator is provided, then there must be an ASE SinglePointCalculator <ase.calculators.singlepoint> object attached to the structures or the forces should already be attached as arrays to the structures.

If a calculator is provided then it will be used to compute the forces for all structures.

Example

The following example illustrates the use of this function:

db = connect('dft_training_structures.db')
training_structures = [row.toatoms() for row in db.select()]
training_structures = prepare_structures(training_structures, atoms_ideal)
for s in training_structures:
    sc.add_structure(s)
Parameters
  • structures (List[Atoms]) – list of input displaced structures

  • atoms_ideal (Atoms) – reference structure relative to which displacements are computed

  • calc (Optional[SinglePointCalculator]) – ASE calculator used for computing forces

Return type

list of prepared structures with forces and displacements as arrays