Tools

Convex hull construction

class icet.tools.convex_hull.ConvexHull(concentrations, energies)

Convex hull of concentration and energies.

dimensions

int – Number of independent concentrations needed to specify a point in concentration space (1 for binaries, 2 for ternaries etc)

concentrations

NumPy array (N, dimensions) – Concentration of the N structures on the convex hull.

energies

NumPy array – Energy of the N structures on the convex hull.

structures

list of ints – Indices of structures that constitute the convex hull (indices are defined by the order their concentrations and energies are fed when initializing the ConvexHull object).

extract_low_energy_structures(concentrations, energies, energy_tolerance, structures=None)

Extract structures that are sufficiently close in energy to the convex hull.

Parameters:
  • concentrations (list of lists) – Concentrations of candidate structures. If there is one independent concentration, a list of floats is fine. Otherwise, the concentrations are given as a list of lists, such as [[0.1, 0.2], [0.3, 0.1], …]
  • energies (list of floats) – Energies of candidate structures.
  • energy_tolerance (float) – Every structure that has an energy within energy_tolerance from the convex hull at its concentration will be returned.
  • structures (list) – List of candidate ASE Atoms or other objects corresponding to the concentrations and energies. The same list will be returned, but with the objects too far from the convex hull removed. If None, a list of indices is returned instead.
Returns:

The members of structures whose energy is within energy_tolerance from the convex hull. If structures is left empty, a list of indices corresponding to the order in concentrations/energies is returned instead.

Return type:

list

get_energy_at_convex_hull(target_concentrations)

Get energy of convex hull at specified concentrations.

Parameters:target_concentrations (list of lists) – Concentrations at target points. If there is one independent concentration, a list of floats is fine. Otherwise, the concentrations are given as a list of lists, such as [[0.1, 0.2], [0.3, 0.1], …]
Returns:Energies at the specified target_concentrations. If any concentration is outside the allowed range, NaN is returned.
Return type:NumPy array

Mapping structures

icet.tools.structure_mapping.map_structure_to_reference(input_structure, reference_structure, tolerance_mapping, vacancy_type=None, inert_species=None, tolerance_cell=0.05, tolerance_positions=0.01, verbose=False)

Map a relaxed structure onto a reference structure.

Parameters:
  • input_structure (ASE Atoms object) – relaxed input structure
  • reference_structure (ASE Atoms object) – reference structure, which can but need not represent the primitive cell
  • tolerance_mapping (float) – maximum allowed displacement for mapping an atom in the relaxed (but rescaled) structure to the reference supercell Note: A reasonable choice is up to 20-30% of the first nearest neighbor distance (r1). A value above 50% of r1 will most likely lead to atoms being multiply assigned, whereby the mapping fails.
  • vacancy_type (str) – If this parameter is set to a non-zero string unassigned sites in the reference structure will be assigned to this type. Note 1: By default (None) the method will fail if there are any unassigned sites in the reference structure. Note 2: vacancy_type must be a valid element type as enforced by the ASE Atoms class.
  • inert_species (list of str) – List of chemical symbols (e.g., [‘Au’, ‘Pd’]) that are never substituted for a vacancy. Used to make an initial rescale of the cell and thus increases the probability for a successful mapping. Need not be specified if vacancy_type is None.
  • tolerance_cell (float) – tolerance factor applied when computing permutation matrix to generate supercell (h = P h_p)
  • tolerance_positions (float) – tolerance factor applied when scanning for overlapping positions in Angstrom (forwarded to ase.build.cut)
  • verbose (bool) – turn on verbose output
Returns:

  • the ideal supercell most closely matching the input structure
  • the largets deviation of any input coordinate from its ideal coordinate
  • the average deviation of the input coordinates from the ideal coordinates

Return type:

ASE Atoms, float, float

Example

The following code snippet illustrates the general usage. It first creates a primitive FCC cell, which is latter used as reference structure. To emulate a relaxed structure obtained from, e.g., a density functional theory calculation, the code then creates a 4x4x4 conventional FCC supercell, which is populated with two different atom types, has distorted cell vectors, and random displacements to the atoms. Finally, the present function is used to map the structure back the ideal lattice:

from ase.build import bulk
reference = bulk('Au', a=4.09)
atoms = bulk('Au', cubic=True, a=4.09).repeat(4)
atoms.set_chemical_symbols(10 * ['Ag'] + (len(atoms) - 10) * ['Au'])
atoms.set_cell(atoms.cell * 1.02, scale_atoms=True)
atoms.rattle(0.1)
mapped_atoms = map_structure_to_reference(atoms, reference, 1.0)

Structure enumeration

This module has the purpose of enumerating structures. Given a lattice (possibly with as basis) and a number of elements, the code generates all the derivative superstructures having a certain size defined by the user.

The algorithm was developed by Gus L. W Hart and Rodney W. Forcade in

  • Hart, G. L. W. and Forcade, R. W., Phys. Rev. B 77, 224115 (2008)
  • Hart, G. L. W. and Forcade, R. W., Phys. Rev. B 80, 014120 (2009)
icet.tools.structure_enumeration.enumerate_structures(atoms, sizes, subelements, concentration_restrictions=None, niggli_reduce=None)

Generate enumerated structures, i.e. all inequivalent structures up to a certain size.

The algorithm implemented here was developed by Gus L. W. Hart and Rodney W. Forcade in

The function is sensitive to the boundary conditions of the input structure. An enumeration of, for example, a surface can thus be performed by setting atoms.pbc = [True, True, False].

Parameters:
  • atoms (ASE Atoms) – Primitive structure from which derivative superstructures should be generated.
  • sizes (list of ints) – Maximum number of atoms in the returned structures.
  • subelements (list of str) – Elements to decorate the structure, e.g. [‘Au’, ‘Ag’]
  • concentration_restrictions (dict) – Defines allowed concentration for one or more element in subelements, e.g. {‘Au’: (0, 0.2)} will only enumerate structures in which the Au content is between 0 and 20 %. Concentration is here always defined as the number of atoms of the specified kind divided by the number of all atoms.
  • niggli_reduction (bool) – If True perform a Niggli reduction with spglib for each structure. Default is True if atoms has all boundary conditions periodic, otherwise False.
Yields:

ASE Atoms – Enumerated structure, each and every one of which is unique.

icet.tools.structure_enumeration.get_symmetry_operations(atoms, tol=0.001)

Use spglib to calculate the symmetry operations of atoms. The symmetry operations consist of three parts: rotations, translation and “basis_shifts”. The latter define the way that the sublattices shift upon rotation (correponds to d_Nd in [HarFor09]).

Parameters:atoms (ASE Atoms) – Structure for which the symmetry operations are sought.
Returns:Containing rotations, translations and basis_shifts.
Return type:dict of lists