Core components#

The low-level functions and classes described here are part of the hiPhive core components and not intended to be accessed directly by users during normal operation.

config#

hiphive.core.config.constraint_vectors_compress_mode: str = 'symbolic'#

Possible values are 'symbolic', 'numeric', None.

hiphive.core.config.constraint_vectors_simplify_before_compress: bool = True#

Possible values are True, False.

hiphive.core.config.constraint_vectors_simplify_before_solve: bool = True#

Possible values are True, False.

hiphive.core.config.constraint_vectors_solve_mode: str = 'symbolic'#

Possible values are 'symbolic', 'numeric'.

hiphive.core.config.eigentensor_compress_before_solve: str = None#

If the compress_mode is None the constraint matrix might be compressed right before the nullspace() solver. Possible values are 'numeric', 'symbolic', None.

hiphive.core.config.eigentensor_compress_mode: str = None#

For each symmetry, the constraint matrix can be reduced to square again. This can be done either by 'symbolic', 'numeric' or not at all (None). Default is None since the matrix is often small enough to fit in memory.

hiphive.core.config.eigentensor_simplify_before_compress: bool = False#

If this is True, before every symbolic compression the values will be simplified by sympy, potentially turning them into exact rational or irrational numbers. This can be useful for systems with non-integer rotation matrices in cartesian coordinates, e.g., hcp. The main purpose is to make the rref more stable against repeating rounding errors.

hiphive.core.config.eigentensor_simplify_before_last_compress: bool = False#

If non compress was used during construction but used before solving.

hiphive.core.config.eigentensor_simplify_before_solve: bool = True#

This might make the nullspace() solver more stable.

hiphive.core.config.eigentensor_solve_mode: str = 'symbolic'#

Possible values are 'symbolic', 'numeric'.

hiphive.core.config.integer_tolerance: float = 1e-12#

This is the tolerance used when the rotation matrices and similar are tested if they are integers.

hiphive.core.config.sum_rule_constraint_mode: str = 'symbolic'#

Possible values are 'symbolic', 'numeric'.

hiphive.core.config.sum_rule_constraint_simplify: bool = True#

Possible values are True, False.

cluster_space_builder#

hiphive.core.cluster_space_builder.build_cluster_space(cluster_space, prototype_structure)[source]#

The permutation list is an indexed fast lookup table for permutation vectors.

clusters#

hiphive.core.clusters.get_clusters(atoms, cutoffs, n_prim, multiplicity=True, use_geometrical_order=False)[source]#

Generates a list of all clusters in the atoms object which includes the center atoms with positions within the cell metric. The cutoff determines up to which order and range clusters should be generated.

With multiplicity set to True clusters like [0, 0] and [3, 3, 4] etc will be generated. This is useful when working with force constants but not so much for cluster expansions.

The geometrical order is the total number of different atoms in the cluster. [0, 0, 1] would have geometrical order 2 and [1, 2, 3, 4] would have order 4. If the key word is True the cutoff criteria will be based on the geometrical order of the cluster. This is based on the observation that many body interactions decrease fast with cutoff but anharmonic interactions can be quite long ranged.

Parameters:
  • atoms (Atoms) – Atomic structure; can be a general ASE:class:Atoms object but must have pbc=False.

  • cutoffs (dict[float]) – Cutoffs by order; the keys specify the order while the values specify the cutoff radii.

  • n_prim (int) – Number of atoms in the primitive cell.

  • multiplicity (bool) – Include clusters in which same atom appears more than once.

  • use_geometrical_order (bool) – Specify if the geometrical order should be used as cutoff order, otherwise the normal order of the cluster is used.

Return type:

list[tuple[int]]

Returns:

A list of clusters where each entry is a tuple of indices, which refer to the atoms in the input supercell.

orbits#

Contains the Orbit class which hold onformation about equivalent clusters.

class hiphive.core.orbits.Orbit[source]#

This class serves as a container for storing data pertaining to an orbit.

orientation_families#

Orientation families of the orbit.

Type:

list[OrientationFamily]

eigensymmetries#

Each eigensymmetry corresponds to a pair where the first index is the symmetry and the second is the permutation.

Type:

list[tuple]

eigentensors#

Decomposition of the force constant into symmetry elements.

Type:

list[np.ndarray]

property prototype_index: int#

Index of cluster that serves as prototype for this orbit. In the code the first symmetry is always the identity so the first orientation family should always correspond to the prototype.

static read(f)[source]#

Load an Orbit from file.

Parameters:

f (Union[str, BinaryIO, TextIO]) – Name of input file (str) or stream to load from (file object).

write(f)[source]#

Write an Orbit instance to a file.

Parameters:

f (Union[str, BinaryIO, TextIO]) – Name of input file (str) or stream to write to (file object).

Return type:

None

class hiphive.core.orbits.OrientationFamily(symmetry_index=None)[source]#

A container for storing information for a “family of orientations”.

An orbit contains many clusters. Some of the clusters can be translated onto each other and other must first be rotated. A set of clusters in the orbit which can all be translated onto each other are oriented in the same way and belong to the same orientation family. The family is characterized by the symmetry (rotation) which relates it to the prototype structure of the orbit. Since the clusters are generally stored sorted the permutation information must also be stored.

Parameters:

symmetry_index (int) – The index of the symmetry corresponding to spglibs symmetry.

symmetry_index#

The index of the symmetry corresponding to the symmetry according to spglib.

Type:

int

cluster_indices#

The indices of the clusters belonging to this family.

Type:

list[int]

permutation_indices#

The indices of the permutation vector.

Type:

list[int]

static read(f)[source]#

Load an OrientationFamily object from a pickle file.

Parameters:

f (Union[str, BinaryIO, TextIO]) – Name of input file (str) or stream to load from (file object).

Return type:

OrientationFamily

write(f)[source]#

Write the object to file.

Parameters:

f (Union[str, BinaryIO, TextIO]) – Name of input file (str) or stream to write to (file object).

hiphive.core.orbits.get_geometrical_radius(positions)[source]#

Compute the geometrical size of a 3-dimensional point cloud. The geometrical size is defined as the average distance to the geometric center.

Parameters:

positions (list[ndarray]) – Positions of points in cloud.

Return type:

float

Returns:

Geometric size of point cloud.

hiphive.core.orbits.get_maximum_distance(positions)[source]#

Compute the maximum distance between any two points in a 3-dimensional point cloud. This is equivalent to the “size” criterion used when imposing a certain (pair) cutoff criterion during construction of a set of clusters.

Parameters:

positions (list[ndarray]) – Positions of points in cloud.

Return type:

float

Returns:

Maximum distance betwee any two points.

hiphive.core.orbits.get_orbits(cluster_list, atom_list, rotation_matrices, translation_vectors, permutations, prim, symprec)[source]#

Generate a list of the orbits for the clusters in a supercell configuration.

This method requires as input a list of the clusters in a supercell configuration as well as a set of symmetry operations (rotations and translations). From this information it will generate a list of the orbits, i.e. the set of symmetry inequivalent clusters each associated with its respective set of equivalent clusters.

Parameters:
  • cluster_list (BiMap) – List of clusters.

  • atom_list (BiMap) – List of atoms in a supercell.

  • rotation_matrices (list[ndarray]) – Rotational symmetries as 3x3 matrices to be imposed (e.g., from spglib).

  • translation_vectors (list[ndarray]) – Translational symmetries as 3 vectors to be imposed (e.g., from spglib).

  • permutations (list) – Look-up table for permutations.

  • prim (Atoms) – Primitive structure.

Return type:

list[Orbit]

Returns:

Orbits associated with the list of input clusters.

rotational_constraints#

Functionality for 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 (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)

translational_constraints#

Functionality for enforcing translational sum rules

atoms#

Collection of functions and classes for handling information concerning atoms and structures, including the relationship between primitive cell and supercells that are derived thereof.

class hiphive.core.atoms.Atom(site, offset)[source]#

Unique representation of an atom in a lattice with a basis

Class for storing information about the position of an atom in a supercell relative to the origin of the underlying primitive cell. This class is used for handling the relationship between a primitive cell and supercells derived thereof.

Parameters:
  • site (int) – site index

  • offset (list(float) or numpy.ndarray) – must contain three elements, offset_x, offset_y, offset_z

property offset#

translational offset of the supercell site relative to the origin of the primitive cell in units of primitive lattice vectors

Type:

list(int)

property site#

index of corresponding site in the primitive basis

Type:

int

class hiphive.core.atoms.Atoms(symbols=None, positions=None, numbers=None, tags=None, momenta=None, masses=None, magmoms=None, charges=None, scaled_positions=None, cell=None, pbc=None, celldisp=None, constraint=None, calculator=None, info=None, velocities=None)[source]#

Minimally augmented version of the ASE Atoms class suitable for handling primitive cell information.

Saves and loads by pickle.

property basis#

scaled coordinates of the sites in the primitive basis

Type:

numpy.ndarray

static read(f)[source]#

Load an hiPhive Atoms object from file.

Parameters:

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

Return type:

hiPhive Atoms object

write(f)[source]#

Writes the object to file.

Note: Only the cell, basis and numbers are stored!

Parameters:

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

hiphive.core.atoms.atom_to_spos(atom, basis)[source]#

Helper function for obtaining the position of a supercell atom in scaled coordinates.

Parameters:
  • atom (hiPhive.Atom) – supercell atom

  • basis (list(list(float)) or numpy.ndarray) – positions of sites in the primitive basis

Returns:

scaled coordinates of an atom in a supercell

Return type:

numpy.ndarray

hiphive.core.atoms.spos_to_atom(spos, basis, tol=0.0001)[source]#

Helper function for transforming a supercell position to the primitive basis.

Parameters:
  • spos (list(list(float)) or numpy.ndarray) – scaled coordinates of an atom in a supercell

  • basis (list(list(float)) or numpy.ndarray) – positions of sites in the primitive basis

  • tol (float) – a general tolerance

Returns:

supercell atom

Return type:

hiphive.Atom

structures#

class hiphive.core.structures.Atom(*args, **kwargs)[source]

This class represents a crystal atom in a given structure.

class hiphive.core.structures.BaseAtom(site, offset)[source]

This class represents an atom placed in an infinite crystal.

site

Site index.

Type:

int

offset

Offset in scaled coordinates.

Type:

list[int]

astype(dtype)[source]

Useful arguments: list, tuple, np.int64

property offset: ndarray

Offset in scaled coordinates.

property site: int

Site index.

class hiphive.core.structures.Structure(atoms, symprec=1e-06)[source]

This class essentially wraps the ASE Atoms class but is a bit more careful with respect to periodic boundary conditions and scaled coordinates. Note that it returns hiphive.Atom objects when queried for individual atoms.

Parameters:
  • atoms (Atoms) – Atomic structure.

  • symprec (float) – Tolerance imposed when rounding scaled coordinates.

atom_from_pos(pos, symprec=None)[source]

Returns atom given a position in Cartesian coordinates.

Parameters:
  • pos (ndarray) – Position in Cartesian coordinates.

  • symprec (float) – Tolerance imposed when rounding scaled coordinates.

Return type:

Atom

property cell

Cell metric.

property spos: ndarray

Scaled coordinates.

class hiphive.core.structures.Supercell(supercell, prim, symprec)[source]

This class tries to represent atoms in a supercell as positioned on the primitive lattice.

Parameters:
  • supercell (Atoms) – Supercell structure.

  • prim (Atoms) – Primitive structure.

  • symprec (float) – Tolerance imposed when rounding scaled coordinates.

index(site, offset)[source]

“ Returns index of atom given a site index and an offset.

Parameters:
  • site (int) – Site index.

  • offset (ndarray) – Offset in scaled coordinates (should be integers).

Return type:

int

class hiphive.core.structures.SupercellAtom(*args, **kwargs)[source]

Represents an atom in a supercell but site and offset given by an underlying primitve cell.

hiphive.core.structures.pos_to_site_offset(pos, cell, basis_spos, symprec)[source]

Returns offset given the position in Cartesian coordinates and the cell metric.

Parameters:
  • pos (ndarray) – Position in Cartesian coordinates.

  • cell (ndarray) – Cell metric.

  • basis_spos (list[ndarray]) – Positions of atoms in basis in scaled coordinates.

  • symprec (float) – Tolerance imposed when rounding the offset.

Return type:

ndarray

hiphive.core.structures.pos_to_spos(pos, cell)[source]

Returns the scaled coordinate given the Cartesian coordinate and the cell metric. Inverse of spos_to_pos().

Parameters:
  • pos (ndarray) – Position in Cartesian coordinates.

  • cell (ndarray) – Cell metric.

Return type:

ndarray

hiphive.core.structures.site_offset_to_pos(site, offset, cell, basis_spos)[source]

Returns the position in Cartesian coordinates given a site index, an offset, the cell metric, and the position in the basis.

Parameters:
  • site (int) – Site index.

  • offset (ndarray) – Offset in scaled coordinates (should be integers).

  • cell (ndarray) – Cell metric.

  • basis_spos (list[ndarray]) – Positions of atoms in basis in scaled coordinates.

Return type:

ndarray

hiphive.core.structures.site_offset_to_spos(site, offset, basis_spos)[source]

Returns the scaled position of an atom at the specified site and offset relative to the basis in scaled coordinates.

Parameters:
  • site (int) – Site index.

  • offset (ndarray) – Offset in scaled coordinates (should be integers).

  • basis_spos (list[ndarray]) – Positions of atoms in basis in scaled coordinates.

Return type:

ndarray

hiphive.core.structures.spos_to_pos(spos, cell)[source]

Returns the Cartesian coordinate given the scaled coordinate and cell metric (cell vectors as rows).

Parameters:
  • spos (ndarray) – Position in scaled coordinates.

  • cell (ndarray) – Cell metric.

Return type:

ndarray

hiphive.core.structures.spos_to_site_offset(spos, basis_spos, symprec)[source]

Returns the site and offset of the atom at the specified scaled coordinate given the scaled positions of the basis atoms.

Parameters:
  • spos (ndarray) – Position in scaled coordinates.

  • basis_spos (list[ndarray]) – Positions of atoms in basis in scaled coordinates.

  • symprec (float) – Tolerance imposed when rounding the offset.

Return type:

ndarray

structure_alignment#

hiphive.core.structure_alignment.align_supercell(supercell, prim, symprec=None)[source]#

Rotates and translates a supercell configuration such that it is aligned with the target primitive cell.

Parameters:
  • supercell (Atoms) – Supercell configuration.

  • prim (Atoms) – Target primitive configuration.

  • symprec (float) – Precision parameter forwarded to spglib.

Return type:

tuple[Atoms, ndarray, ndarray]

Returns:

Aligned supercell configuration as well as rotation matrix (3x3 array) and translation vector (3x1 array) that relate the input to the aligned supercell configuration.

hiphive.core.structure_alignment.are_nonpaired_configurations_equal(atoms1, atoms2)[source]#

Checks whether two configurations are identical. To be considered equal the structures must have the same cell metric, elemental occupation, scaled positions (modulo one), and periodic boundary conditions.

Unlike the __eq__ operator of ase.Atoms the order of the atoms does not matter.

Parameters:
Return type:

bool

Returns:

True if atoms are equal, False otherwise

hiphive.core.structure_alignment.get_primitive_cell(atoms, to_primitive=True, no_idealize=True, symprec=1e-05)[source]#

Returns primitive cell obtained using spglib.

Parameters:
  • atoms (Atoms) – Atomic structure.

  • to_primitive (bool) – Passed to spglib.

  • no_idealize (bool) – Passed to spglib.

  • symprec (float) – Numerical tolerance; passed to spglib.

Return type:

Atoms

hiphive.core.structure_alignment.is_rotation(R, cell_metric=None)[source]#

Checks if a rotation matrix is orthonormal. A cell metric can be passed if the rotation matrix is in scaled coordinates

Parameters:
  • R (ndarray) – Rotation matrix (3x3 array).

  • cell_metric (ndarray) – Cell metric if the rotation is in scaled coordinates.

hiphive.core.structure_alignment.relate_structures(reference, target, symprec=1e-05)[source]#

Finds rotation and translation operations that align two structures with periodic boundary conditions. The rotation and translation in Cartesian coordinates will map the reference structure onto the target.

Aligning the reference with the target structure can be achieved via the transformations:

R, T = relate_structures(atoms_ref, atoms_target)
atoms_ref_rotated = rotate_atoms(atoms_ref, R)
atoms_ref_rotated.translate(T)
atoms_ref_rotated.wrap()
atoms_ref_rotated == atoms_target
Parameters:
  • reference (Atoms) – Reference structure to be mapped.

  • target (Atoms) – Target structure.

Return type:

tuple[ndarray, ndarray]

Returns:

A tuple comprising the rotation matrix in Cartesian coordinates (3x3 array) and the translation vector in Cartesian coordinates.

hiphive.core.structure_alignment.rotate_atoms(atoms, rotation)[source]#

Rotates the cell and positions of atoms and returns a copy.

Parameters:
  • atoms (Union[Atoms, Atoms]) – Atomic structure.

  • rotation (ndarray) – Rotation matrix (3x3 array).

Return type:

Atoms

tensors#

Module containing tensor related functions

hiphive.core.tensors.rotate_tensor(T, R, path=None)[source]#

Equivalent to \(T_{abc} \ldots = T_{ijk} \ldots R_{ia} R_{jb} R_{kc} \ldots\).

hiphive.core.tensors.rotation_to_cart_coord(R, cell)[source]#

Returns the rotation matrix in cart coord given a cell metric.

utilities#

class hiphive.core.utilities.BiMap[source]#

Simple list like structure with fast dict-lookup.

The structure can append objects and supports some simple list interfaces. The lookup is fast since an internal dict stores the indices.