Structures

Preparing structures for training

hiphive provides some utility functions for handling training :structures.

For example get_displacements can be used to calculate the displacements of atoms given a reference structure accounting while accounting for periodic boundary conditions.

In cases where the ordering/indexing of atoms in the training structures differs from the reference structure, one can use find_permutation which finds the permutation that re-orderes the training structures to match the reference.

The function prepare_structures combines these two functions together and adds displacements and forces as arrays to the structures such that they can be added directly a StructureContainer.

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

StructureContainer

hiphive organizes the data for training and testing in a so-called :structure container, which provides methods for adding and accessing these :data. The structure container is essentially a collection of :FitStructure objects, each of which comprises an original atomic :configuration along :with its representation in clusters and reference :forces.

class hiphive.StructureContainer(cs, fit_structure_list=None)[source]

This class serves as a container for structures as well as associated fit properties and fit matrices.

Parameters
  • cs (ClusterSpace) – cluster space that is the basis for the container

  • fit_structure_list (list(FitStructure)) – structures to be added to the container

add_structure(atoms, **meta_data)[source]

Add a structure to the container.

Note that custom information about the atoms object may not be stored inside, for example an ASE SinglePointCalculator will not be kept.

Parameters
  • atoms (ase.Atoms) – the structure to be added; the Atoms object must contain supplementary per-atom arrays with displacements and forces

  • meta_data (dict) – dict with meta_data about the atoms

property cluster_space

copy of the cluster space the structure container is based on

Type

ClusterSpace

property data_shape

tuple of integers representing the shape of the fit data matrix

Type

tuple

delete_all_structures()[source]

Remove all current structures in StructureContainer.

get_fit_data(structures=None)[source]

Return fit data for structures. The fit matrices and target forces for the structures are stacked into NumPy arrays.

Parameters

structures (list, tuple) – list of integers corresponding to structure indices. Defaults to None and in that case returns all fit data available.

Returns

stacked fit matrices, stacked target forces for the structures

Return type

numpy.ndarray, numpy.ndarray

static read(fileobj, read_structures=True)[source]

Restore a StructureContainer object from file.

Parameters
  • f (str or file object) – name of input file (str) or stream to load from (file object)

  • read_structures (bool) – if True the structures will be read; if False only the cluster space will be read

write(f)[source]

Write a StructureContainer instance to a file.

Parameters

f (str or file object) – name of input file (str) or stream to write to (file object)

FitStructure

class hiphive.structure_container.FitStructure(atoms, fit_matrix, **meta_data)[source]

This class holds a structure with displacements and forces as well as the fit matrix.

Parameters
  • atoms (ase.Atoms) – supercell structure

  • fit_matrix (numpy.ndarray) – fit matrix, N, M array with N = 3 * len(atoms)

  • meta_data (dict) – any meta data that needs to be stored in the FitStructure

property atoms

supercell structure

Type

ase.Atoms

property displacements

atomic displacements

Type

numpy.ndarray

property fit_matrix

the fit matrix

Type

numpy.ndarray

property forces

forces

Type

numpy.ndarray

static read(fileobj, read_fit_matrix=True)[source]

Read a OrientationFamily instance from a file.

Parameters
  • fileobj (str or file object) – name of input file (str) or stream to read from (file object)

  • read_fit_matrix (bool) – whether or not to read the fit_matrix

Return type

FitStructure instance

write(fileobj)[source]

Write the instance to file.

Parameters

fileobj (str or file object) – name of input file (str) or stream to write to (file object)

Structure generation

This module enables the generation of structures with displacements that can be used to generate reference forces.

hiphive.structure_generation.generate_mc_rattled_structures(atoms, n_structures, rattle_std, d_min, seed=42, **kwargs)[source]

Returns list of Monte Carlo rattled configurations.

Rattling atom i is carried out as a Monte Carlo move that is accepted with a probability determined from the minimum interatomic distance \(d_{ij}\). If \(\min(d_{ij})\) is smaller than \(d_{min}\) the move is only accepted with a low probability.

This process is repeated for each atom a number of times meaning the magnitude of the final displacements is not directly connected to rattle_std.

Warning

Repeatedly calling this function without providing different seeds will yield identical or correlated results. To avoid this behavior it is recommended to specify a different seed for each call to this function.

Notes

The procedure implemented here might not generate a symmetric distribution for the displacements kwargs will be forwarded to mc_rattle (see user guide for a detailed explanation)

The displacements generated will roughly be n_iter**0.5 * rattle_std for small values of n_iter

Parameters
  • atoms (ase.Atoms) – prototype structure

  • n_structures (int) – number of structures to generate

  • rattle_std (float) – rattle amplitude (standard deviation in normal distribution); note this value is not connected to the final average displacement for the structures

  • d_min (float) – interatomic distance used for computing the probability for each rattle move

  • seed (int) – seed for setting up NumPy random state from which random numbers are generated

  • n_iter (int) – number of Monte Carlo cycles (iterations), larger number of iterations will generate larger displacements (defaults to 10)

Returns

generated structures

Return type

list of ase.Atoms

hiphive.structure_generation.generate_phonon_rattled_structures(atoms, fc2, n_structures, temperature, QM_statistics=False, imag_freq_factor=1.0)[source]

Returns list of phonon-rattled configurations.

Configurations are generated by superimposing harmonic phonon eigenmodes with random amplitudes and phase factors consistent with a certain temperature.

Let \(\boldsymbol{X}_{ai}\) be the phonon modes indexed by atom \(a\) and mode \(i\), \(\omega_i\) the phonon frequencies, and let \(0 < Q_i \leq 1\) and \(0 \leq U_i < 1\) be uniformly random numbers. Then

\[\boldsymbol{R}_a = \boldsymbol{R}^0_a + \left<\frac{k_B T}{m_a} \right>^{1/2} \sum_i \frac{1}{\omega_i} \boldsymbol{X}_{ai} \sqrt{-2 \ln Q_i} \cos(\pi \omega_i U_i)\]

See: West and Estreicher, Physical Review Letters 96, 115504 (2006)

or

\[\boldsymbol{R}_a = \boldsymbol{R}^0_a + \sum_i \left(\frac{\hbar(0.5 + n(\omega_i, T))}{m_a\omega_i} \right)^{1/2} \boldsymbol{X}_{ai} \sqrt{-2 \ln Q_i} \cos(\pi \omega_i U_i)\]

where

\[n = \frac{1}{e^{\hbar\omega_i/k_BT}-1}\]
Parameters
  • atoms (ase.Atoms) – prototype structure

  • fc2 (numpy.ndarray) – second order force constant matrix, with shape (3N, 3N) or (N, N, 3, 3). The conversion will be done internally if.

  • n_structures (int) – number of structures to generate

  • temperature (float) – temperature in Kelvin

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

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

Returns

generated structures

Return type

list(ase.Atoms)

hiphive.structure_generation.generate_rattled_structures(atoms, n_structures, rattle_std, seed=42)[source]

Returns list of rattled configurations.

Displacements are drawn from normal distributions for each Cartesian directions for each atom independently.

Warning

Repeatedly calling this function without providing different seeds will yield identical or correlated results. To avoid this behavior it is recommended to specify a different seed for each call to this function.

Parameters
  • atoms (ase.Atoms) – prototype structure

  • n_structures (int) – number of structures to generate

  • rattle_std (float) – rattle amplitude (standard deviation of the normal distribution)

  • seed (int) – seed for setting up NumPy random state from which random numbers are generated

Returns

generated structures

Return type

list of ase.Atoms

Other functions

hiphive.structure_container.are_configurations_equal(atoms1, atoms2, tol=1e-10)[source]

Compare if two configurations are equal within some tolerance. This includes checking all available arrays in the two atoms objects.

Parameters
Returns

True if atoms are equal, False otherwise

Return type

bool