Tools

Convex hull construction

This module provides the ConvexHull class.

class icet.tools.convex_hull.ConvexHull(concentrations, energies)[source]

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(float) or list(list(float))) – 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(float)) – energy (or energy of mixing) for each structure
concentrations

concentrations of the N structures on the convex hull

Type:np.ndarray
energies

energies of the N structures on the convex hull

Type:np.ndarray
dimensions

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

Type:int
structures

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)

Type:list(int)

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

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.

Return type:

List[int]

get_energy_at_convex_hull(target_concentrations)[source]

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

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 species 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)
structure = bulk('Au', cubic=True, a=4.09).repeat(4)
structure.set_chemical_symbols(10 * ['Ag'] + (len(structure) - 10) * ['Au'])
structure.set_cell(structure.cell * 1.02, scale_atoms=True)
structure.rattle(0.1)
mapped_structure = map_structure_to_reference(structure, reference, 1.0)
Return type:Tuple[Atoms, float, float]

Structure enumeration

icet.tools.structure_enumeration.enumerate_structures(structure, sizes, chemical_symbols, concentration_restrictions=None, niggli_reduce=None)[source]

Yields a sequence of enumerated structures. The function generates all inequivalent structures that are permissible given a certain lattice. Using the chemical_symbols and concentration_restrictions keyword arguments it is possible to specify which chemical_symbols 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 structure.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:
  • structure (Atoms) – primitive structure from which derivative superstructures should be generated
  • sizes (List[int]) – number of sites (included in enumeration)
  • chemical_symbols (list) – chemical 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 chemical_symbols, 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 structure 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(structure=prim, sizes=range(1, 7),
                     chemical_symbols=['Ag', 'Au'])

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

enumerate_structures(structure=prim, sizes=range(1, 7),
                     chemical_symbols=['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(structure=prim, sizes=range(1, 9),
                     chemical_symbols=[['Ga', 'Al'], ['As']])
Return type:Atoms
icet.tools.structure_enumeration.enumerate_supercells(structure, sizes, niggli_reduce=None)[source]

Yields a sequence of enumerated supercells. The function generates all inequivalent supercells that are permissible given a certain lattice. Any supercell can be reduced to one of the supercells generated.

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 structure.pbc = [True, True, False].

The algorithm is based on 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:
  • structure (Atoms) – primitive structure from which supercells should be generated
  • sizes (List[int]) – number of sites (included in enumeration)
  • niggli_reduction – if True perform a Niggli reduction with spglib for each supercell; the default is True if structure is periodic in all directions, False otherwise.

Examples

The following code snippet illustrates how to enumerate supercells with up to 6 atoms in the unit cell:

from ase.build import bulk
prim = bulk('Ag')
enumerate_supercells(structure=prim, sizes=range(1, 7))
Return type:Atoms
icet.tools.structure_enumeration.get_symmetry_operations(structure, tolerance=0.001)[source]

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:
  • structure (Atoms) – structure for which symmetry operations are sought
  • tolerance (float) – numerical tolerance imposed during symmetry analysis
Return type:

Dict[str, list]

Structure generation and SQS

icet.tools.structure_generation.generate_sqs(cluster_space, max_size, target_concentrations, include_smaller_cells=True, T_start=5.0, T_stop=0.001, n_steps=None, optimality_weight=1.0, random_seed=None, tol=1e-05)[source]

Given a cluster_space, generate a special quasirandom structure (SQS), i.e., a structure that for a given supercell size provides the best possible approximation to a random alloy [ZunWeiFer90].

In the present case, this means that the generated structure will have a cluster vector that as closely as possible matches the cluster vector of an infintely large randomly occupated supercell. Internally the function uses a simulated annealing algorithm and the difference between two cluster vectors is calculated with the measure suggested by A. van de Walle et al. in Calphad 42, 13-18 (2013) [WalTiwJon13] (for more information, see mchammer.calculators.TargetVectorCalculator).

Parameters:
  • cluster_space (ClusterSpace) – a cluster space defining the lattice to be occupated
  • max_size (int) – maximum supercell size
  • target_concentrations (dict) –

    concentration of each species in the target structure (for example {'Ag': 0.5, 'Pd': 0.5}

    Concentrations are always expressed with respect to all atoms in the supercell, which implies that the sum of all concentrations should always be 1. In the case of multiple sublattices, a valid specification would thus be {'Au': 0.25, 'Pd': 0.25, 'H': 0.1, 'V': 0.4}.

  • include_smaller_cells (bool) – if True, search among all supercell sizes including max_size, else search only among those exactly matching max_size
  • T_start (float) – artificial temperature at which the simulated annealing starts
  • T_stop (float) – artifical temperature at which the simulated annealing stops
  • n_steps (Optional[float]) – total number of Monte Carlo steps in the simulation
  • optimality_weight (float) – controls weighting \(L\) of perfect correlations, see mchammer.calculators.TargetVectorCalculator
  • random_seed (Optional[int]) – seed for the random number generator used in the Monte Carlo simulation
  • tol (float) – Numerical tolerance
Return type:

Atoms

icet.tools.structure_generation.generate_sqs_by_enumeration(cluster_space, max_size, target_concentrations, include_smaller_cells=True, optimality_weight=1.0, tol=1e-05)[source]

Given a cluster_space, generate a special quasirandom structure (SQS), i.e., a structure that for a given supercell size provides the best possible approximation to a random alloy [ZunWeiFer90].

In the present case, this means that the generated structure will have a cluster vector that as closely as possible matches the cluster vector of an infintely large randomly occupied supercell. Internally the function uses a simulated annealing algorithm and the difference between two cluster vectors is calculated with the measure suggested by A. van de Walle et al. in Calphad 42, 13-18 (2013) [WalTiwJon13] (for more information, see mchammer.calculators.TargetVectorCalculator).

This functions generates SQS cells by exhaustive enumeration, which means that the generated SQS cell is guaranteed to be optimal with regard to the specified measure and cell size.

Parameters:
  • cluster_space (ClusterSpace) – a cluster space defining the lattice to be occupied
  • max_size (int) – maximum supercell size
  • target_concentrations (dict) –

    concentration of each species in the target structure (for example {'Ag': 0.5, 'Pd': 0.5}

    Concentrations are always expressed with respect to all atoms in the supercell, which implies that the sum of all concentrations should always be 1. In the case of multiple sublattices, a valid specification would thus be {'Au': 0.25, 'Pd': 0.25, 'H': 0.1, 'V': 0.4}.

  • include_smaller_cells (bool) – if True, search among all supercell sizes including max_size, else search only among those exactly matching max_size
  • optimality_weight (float) – controls weighting \(L\) of perfect correlations, see mchammer.calculators.TargetVectorCalculator
  • tol (float) – Numerical tolerance
Return type:

Atoms

icet.tools.structure_generation.generate_target_structure(cluster_space, max_size, target_concentrations, target_cluster_vector, include_smaller_cells=True, T_start=5.0, T_stop=0.001, n_steps=None, optimality_weight=1.0, random_seed=None, tol=1e-05)[source]

Given a cluster_space and a target_cluster_vector, generate a structure that as closely as possible matches that cluster vector. The search is performed among all inequivalent supercells shapes up to a certain size.

Internally the function uses a simulated annealing algorithm and the difference between two cluster vectors is calculated with the measure suggested by A. van de Walle et al. in Calphad 42, 13-18 (2013) [WalTiwJon13] (for more information, see mchammer.calculators.TargetVectorCalculator).

Parameters:
  • cluster_space (ClusterSpace) – a cluster space defining the lattice to be occupied
  • max_size (int) – maximum supercell size
  • target_concentrations (dict) –

    concentration of each species in the target structure (for example {'Ag': 0.5, 'Pd': 0.5}

    Concentrations are always expressed with respect to all atoms in the supercell, which implies that the sum of all concentrations should always be 1. In the case of multiple sublattices, a valid specification would thus be {'Au': 0.25, 'Pd': 0.25, 'H': 0.1, 'V': 0.4}.

  • target_cluster_vector (List[float]) – cluster vector that the generated structure should match as closely as possible
  • include_smaller_cells (bool) – if True, search among all supercell sizes including max_size, else search only among those exactly matching max_size
  • T_start (float) – artificial temperature at which the simulated annealing starts
  • T_stop (float) – artifical temperature at which the simulated annealing stops
  • n_steps (Optional[float]) – total number of Monte Carlo steps in the simulation
  • optimality_weight (float) – controls weighting \(L\) of perfect correlations, see mchammer.calculators.TargetVectorCalculator
  • random_seed (Optional[int]) – seed for the random number generator used in the Monte Carlo simulation
  • tol (float) – Numerical tolerance
Return type:

Atoms