Tools

Convex hull construction

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

This class provides functionality for extracting the convex hull of the (free) energy of mixing. It is based on the convex hull calculator in SciPy.

Parameters:
  • concentrations (list of floats / list of lists of floats) – concentrations for each structure listed as [[c1, c2], [c1, c2], ...]; for binaries, in which case there is only one independent concentration, the format [c1, c2, c3, ...] works as well.
  • energies (list of floats) – energy (or energy of mixing) for each structure
concentrations

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

energies

NumPy array – energies of the N structures on the convex hull

dimensions

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

structures

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

Examples

A ConvexHull object is easily initialized by providing lists of concentrations and energies:

hull = ConvexHull(data['concentration'], data['mixing_energy'])

after which one can for example plot the data (assuming a matplotlib axis object ax):

ax.plot(hull.concentrations, hull.energies)

or extract structures at or close to the convex hull:

low_energy_structures = hull.extract_low_energy_structures(
    data['concentration'], data['mixing_energy'],
    tolerance=0.005, structures=list_of_structures)

A complete example can be found in the basic tutorial.

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

Returns structures that lie within a certain tolerance of the convex hull.

Parameters:
  • concentrations (Union[List[float], List[List[float]]]) –

    concentrations of candidate structures

    If there is one independent concentration, a list of floats is sufficient. Otherwise, the concentrations must be provided as a list of lists, such as [[0.1, 0.2], [0.3, 0.1], ...].

  • energies (List[float]) – energies of candidate structures
  • energy_tolerance (float) – include structures with an energy that is at most this far from the convex hull
  • structures (Optional[list]) –

    list of candidate structures, e.g., ASE Atoms objects, corresponding to concentrations and energies

    The list will be returned, but with the objects too far from the convex hull removed. If None, a list of indices is returned instead.

get_energy_at_convex_hull(target_concentrations)

Returns the energy of the convex hull at specified concentrations. If any concentration is outside the allowed range, NaN is returned.

Parameters:target_concentrations (Union[List[float], List[List[float]]]) –

concentrations at target points

If there is one independent concentration, a list of floats is sufficient. Otherwise, the concentrations ought to be provided as a list of lists, such as [[0.1, 0.2], [0.3, 0.1], ...].

Return type:ndarray

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)

Maps a relaxed structure onto a reference structure. The function returns a tuple comprising

  • the ideal supercell most closely matching the input structure,
  • the largest deviation of any input coordinate from its ideal coordinate, and
  • the average deviation of the input coordinates from the ideal coordinates.
Parameters:
  • input_structure (Atoms) – relaxed input structure
  • reference_structure (Atoms) – 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 (Optional[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 (Optional[List[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
  • tolerance_positions (float) – tolerance factor applied when scanning for overlapping positions in Angstrom (forwarded to ase.build.cut())

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)
Return type:Tuple[Atoms, float, float]

Structure enumeration

icet.tools.structure_enumeration.enumerate_structures(atoms, sizes, species, concentration_restrictions=None, niggli_reduce=None)

Yields a sequence of enumerated structures. The function generates all inequivalent structures that are permissible given a certain lattice. Using the species and concentration_restrictions keyword arguments it is possible to specify which species are to be included on which site and in which concentration range.

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

The algorithm implemented here was developed by Gus L. W. Hart and Rodney W. Forcade in Phys. Rev. B 77, 224115 (2008) [HarFor08] and Phys. Rev. B 80, 014120 (2009) [HarFor09].

Parameters:
  • atoms (Atoms) – primitive structure from which derivative superstructures should be generated
  • sizes (List[int]) – number of sites (included in the enumeration)
  • species (list) – species with which to decorate the structure, e.g., ['Au', 'Ag']; see below for more examples
  • concentration_restrictions (Optional[dict]) – allowed concentration range for one or more element in species, e.g., {'Au': (0, 0.2)} will only enumerate structures in which the Au content is between 0 and 20 %; here, concentration is always defined as the number of atoms of the specified kind divided by the number of all atoms.
  • niggli_reduction – if True perform a Niggli reduction with spglib for each structure; the default is True if atoms is periodic in all directions, False otherwise.

Examples

The following code snippet illustrates how to enumerate structures with up to 6 atoms in the unit cell for a binary alloy without any constraints:

from ase.build import bulk
prim = bulk('Ag')
enumerate_structures(atoms=prim, sizes=range(1, 7),
                     species=['Ag', 'Au'])

To limit the concentration range to 10 to 40% Au the code should be modified as follows:

enumerate_structures(atoms=prim, sizes=range(1, 7),
                     species=['Ag', 'Au'],
                     concentration_restrictions={'Au': (0.1, 0.4)})

Often one would like to consider mixing on only one sublattice. This can be achieved as illustrated for a Ga(1-x)Al(x)As alloy as follows:

prim = bulk('GaAs', crystalstructure='zincblende', a=5.65)
enumerate_structures(atoms=prim, sizes=range(1, 9),
                     species=[['Ga', 'Al'], ['As']])
Return type:Atoms
icet.tools.structure_enumeration.get_symmetry_operations(atoms, tolerance=0.001)

Returns symmetry operations permissable for a given structure as obtained via spglib. The symmetry operations consist of three parts: rotation, translation and basis shifts. The latter define the way that sublattices shift upon rotation (correponds to d_Nd in [HarFor09]).

Parameters:
  • atoms (Atoms) – structure for which symmetry operations are sought
  • tolerance (float) – numerical tolerance imposed during symmetry analysis
Return type:

Dict[str, list]