Other functions

Self-consistent phonons

hiphive.self_consistent_phonons.self_consistent_harmonic_model(atoms_ideal, calc, cs, T, alpha, n_iterations, n_structures, parameters_start=None, 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

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

  • 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)

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)[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.

Parameters
  • fcs (ForceConstants) – force constants

  • cs (ClusterSpace) – cluster space

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)[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

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_structures(structures, atoms_ideal, calc)[source]

Prepares a set of structures in the format suitable for a StructureContainer. This includes retrieving the forces using the provided calculator, resetting the positions to the ideal configuration, and adding arrays with forces and displacements.

Notes

Changes the structures in place.

Parameters
  • structures (List[Atoms]) – list of input configurations

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

  • calc (Calculator) – ASE calculator used for computing forces

Return type

None

Enforcing rotational sum rules

hiphive.core.rotational_constraints.enforce_rotational_sum_rules(cs, parameters, sum_rules, **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 (numpy.ndarray) – parameters to be constrained

  • sum_rules (list(str)) – type of sum rules to enforce; possible values: ‘Huang’, ‘Born-Huang’

  • ridge_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]

  • iterations (int) – number of iterations to run the projection since each step projects the solution down to each nullspace in serial; keyword argument passed to the optimizer [default: 10]

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)