Observers

SiteOccupancyObserver

class mchammer.observers.SiteOccupancyObserver(cluster_space, structure, sites, interval=None)[source]

This class represents a site occupation factor (SOF) observer.

A SOF observer allows to compute the site occupation factors along the trajectory sampled by a Monte Carlo (MC) simulation.

Parameters:
  • cluster_space (icet.ClusterSpace) – cluster space from which the allowed species are extracted

  • structure (ase.Atoms) – supercell consistent with primitive structure in cluster space; used to determine which species are allowed on each site

  • sites (dict(str, list(int))) – dictionary containing lists of sites that are to be considered; the keys will be taken as the names of the sites; the indices refer to the primitive structure associated with the cluster space

  • interval (int) – the observation interval, defaults to None meaning that if the observer is used in a Monte Carlo simulation, then the Ensemble object will set the interval.

tag

name of observer

Type:

str

interval

observation interval

Type:

int

Example

The following snippet illustrate how to use the site occupancy factor (SOF) observer in a Monte Carlo simulation of a surface slab. Here, the SOF observer is used to monitor the concentrations of different species at the surface, the first subsurface layer, and the remaining “bulk”. A minimal cluster expansion is used with slightly modified surface interactions in order to obtain an example that can be run without much ado. In practice, one should of course use a proper cluster expansion:

>>> from ase.build import fcc111
>>> from icet import ClusterExpansion, ClusterSpace
>>> from mchammer.calculators import ClusterExpansionCalculator
>>> from mchammer.ensembles import CanonicalEnsemble
>>> from mchammer.observers import SiteOccupancyObserver

>>> # prepare reference structure
>>> prim = fcc111('Au', size=(1, 1, 10), vacuum=10.0)
>>> prim.translate((0.1, 0.1, 0.0))
>>> prim.wrap()
>>> prim.pbc = True  # icet requires pbc in all directions

>>> # prepare cluster expansion
>>> cs = ClusterSpace(prim, cutoffs=[3.7], chemical_symbols=['Ag', 'Au'])
>>> params = [0] + 5 * [0] + 10 * [0.1]
>>> params[1] = 0.01
>>> params[6] = 0.12
>>> ce = ClusterExpansion(cs, params)
>>> print(ce)

>>> # prepare initial configuration based on a 2x2 supercell
>>> structure = prim.repeat((2, 2, 1))
>>> for k in range(20):
>>>     structure[k].symbol = 'Ag'

>>> # set up MC simulation
>>> calc = ClusterExpansionCalculator(structure, ce)
>>> mc = CanonicalEnsemble(structure=structure, calculator=calc, temperature=600,
...                        dc_filename='myrun_sof.dc')

>>> # set up observer and attach it to the MC simulation
>>> sites = {'surface': [0, 9], 'subsurface': [1, 8],
...          'bulk': list(range(2, 8))}
>>> sof = SiteOccupancyObserver(cs, structure, sites, interval=len(structure))
>>> mc.attach_observer(sof)

>>> # run 1000 trial steps
>>> mc.run(1000)

After having run this snippet one can access the SOFs via the data container:

>>> print(mc.data_container.data)
get_observable(structure)[source]

Returns the site occupation factors for a given atomic configuration.

Parameters:

structure (Atoms) – input atomic structure.

Return type:

Dict[str, List[float]]

property return_type: type

Data type of the observed data.

BinaryShortRangeOrderObserver

class mchammer.observers.BinaryShortRangeOrderObserver(cluster_space, structure, radius, interval=None)[source]

This class represents a short range order (SRO) observer for a binary system.

Parameters:
  • cluster_space (icet.ClusterSpace) – cluster space used for initialization

  • structure (ase.Atoms) – defines the lattice which the observer will work on

  • interval (int) – the observation interval, defaults to None meaning that if the observer is used in a Monte Carlo simulations, then the Ensemble object will set the interval.

  • radius (float) – the maximum radius for the neigbhor shells considered

tag

human readable observer name (BinaryShortRangeOrderObserver)

Type:

str

interval

observation interval

Type:

int

Example

The following snippet illustrate how to use the short-range order (SRO) observer in a Monte Carlo simulation of a bulk supercell. Here, the parameters of the cluster expansion are set to emulate a simple Ising model in order to obtain an example that can be run without modification. In practice, one should of course use a proper cluster expansion:

>>> from ase.build import bulk
>>> from icet import ClusterExpansion, ClusterSpace
>>> from mchammer.calculators import ClusterExpansionCalculator
>>> from mchammer.ensembles import CanonicalEnsemble
>>> from mchammer.observers import BinaryShortRangeOrderObserver

>>> # prepare cluster expansion
>>> # the setup emulates a second nearest-neighbor (NN) Ising model
>>> # (zerolet and singlet ECIs are zero; only first and second neighbor
>>> # pairs are included)
>>> prim = bulk('Au')
>>> cs = ClusterSpace(prim, cutoffs=[4.3], chemical_symbols=['Ag', 'Au'])
>>> ce = ClusterExpansion(cs, [0, 0, 0.1, -0.02])

>>> # prepare initial configuration
>>> nAg = 10
>>> structure = prim.repeat(3)
>>> structure.set_chemical_symbols(nAg * ['Ag'] + (len(structure) - nAg) * ['Au'])

>>> # set up MC simulation
>>> calc = ClusterExpansionCalculator(structure, ce)
>>> mc = CanonicalEnsemble(structure=structure, calculator=calc, temperature=600,
...                        dc_filename='myrun_sro.dc')

# set up observer and attach it to the MC simulation
sro = BinaryShortRangeOrderObserver(cs, structure, interval=len(structure),
                                    radius=4.3)
mc.attach_observer(sro)

# run 1000 trial steps
mc.run(1000)

After having run this snippet one can access the SRO parameters via the data container:

print(mc.data_container.data)
get_observable(structure)[source]

Returns the value of the property from a cluster expansion model for a given atomic configurations.

Parameters:

structure (Atoms) – input atomic structure

Return type:

Dict[str, float]

property return_type: type

Data type of the observed data.

StructureFactorObserver

class mchammer.observers.StructureFactorObserver(structure, q_points, symbol_pairs=None, form_factors=None, interval=None)[source]

This class represents a structure factor observer.

This observer allows one to compute structure factors along the trajectory sampled by a Monte Carlo (MC) simulation. Structure factors are convenient for monitoring long-range order. The structure factor is defined as:

\[S(\vec{q}) = \frac{1}{\sum_{j=1}^N f_j^2} \sum_{j,k}^N e^{-i \vec{q} \cdot (\vec{R}_k - \vec{R}_j)}\]

In addition to this “total” structure factor, this observer calculates pair-specific structure factors, which correspond to parts of the summation defined above, with the summation restricted to pairs of specific types, e.g., Au-Au, Au-Cu and Cu-Cu in the example below.

Parameters:
  • structure (Atoms) – prototype for the structures for which the structure factor will be computed later; the supercell size (but not its decoration) must be identical; this structure is also used to determine the the possible pairs if symbol_pairs=None

  • q_points (List[np.ndarray]) – array of q-points at which to evaluate the structure factor; the q-points need to be compatible with the supercell in order for the structure factor to be real

  • symbol_pairs (Optional[List[Tuple[str, str]]]) – list of symbol pairs for which structure factors will be computed, e.g. [(‘Al’, ‘Cu’), (‘Al’, ‘Al’)]; if None (default) use all pairs possible based on the input structure

  • form_factors (Optional[Dict[str, float]]) – form factors for each atom type; this can be used to (very coarsely) simulate X-ray or neutron spectra; note that in general the form factors are q-dependent, see, e.g., here and here; by default (None) all form factors are set to 1

  • interval (int) – the observation interval, defaults to None meaning that if the observer is used in a Monte Carlo simulation, the Ensemble object will set the interval.

Raises:

ValueError – if q-point is not consistent with metric of input structure

tag

name of observer

Type:

str

interval

observation interval

Type:

int

Example

The following snippet illustrates how to use the structure factor observer in a simulated annealing run of dummy Cu-Au model to observe the emergence of a long-range ordered L1_2 structure:

>>> import numpy as np
>>> from ase.build import bulk
>>> from icet import ClusterSpace, ClusterExpansion
>>> from mchammer.calculators import ClusterExpansionCalculator
>>> from mchammer.ensembles import CanonicalAnnealing
>>> from mchammer.observers import StructureFactorObserver

>>> # parameters
>>> size = 3
>>> alat = 4.0
>>> symbols = ['Cu', 'Au']

>>> # setup
>>> prim = bulk('Cu', a=alat, cubic=True)
>>> cs = ClusterSpace(prim, [0.9*alat], symbols)
>>> ce = ClusterExpansion(cs, [0, 0, 0.2])

>>> # make supercell
>>> supercell = prim.repeat(size)
>>> ns = int(0.25 * len(supercell))
>>> supercell.symbols[0:ns] = 'Au'
>>> np.random.shuffle(supercell.symbols)

>>> # define q-points to sample
>>> q_points = []
>>> q_points.append(2 * np.pi / alat * np.array([1, 0, 0]))
>>> q_points.append(2 * np.pi / alat * np.array([0, 1, 0]))
>>> q_points.append(2 * np.pi / alat * np.array([0, 0, 1]))

>>> # set up structure factor observer
>>> sfo = StructureFactorObserver(supercell, q_points)

>>> # run simulation
>>> calc = ClusterExpansionCalculator(supercell, ce)
>>> mc = CanonicalAnnealing(supercell, calc,
...                         T_start=900, T_stop=500, cooling_function='linear',
...                         n_steps=400*len(supercell))
>>> mc.attach_observer(sfo)
>>> mc.run()

After having run this snippet, one can access the structure factors via the data container:

>>> dc = mc.data_container
>>> print(dc.data)

The emergence of the ordered low-temperature structure can be monitored by following the temperature dependence of any of the pair-specific structure factors.

property form_factors: Dict[str, float]

Form factors used in structure factor calculation

get_observable(structure)[source]

Returns the structure factors for a given atomic configuration.

Parameters:

structure (Atoms) – input atomic structure

Raises:

ValueError – if input structure is incompatible with structure used for initialization

Return type:

Dict[str, float]

property q_points: List[ndarray]

q-points for which structure factor is calculated

property return_type: type

Data type of the observed data.

ClusterCountObserver

class mchammer.observers.ClusterCountObserver(cluster_space, structure, interval=None, orbit_indices=None)[source]

This class represents a cluster count observer.

A cluster count observer enables one to keep track of the occupation of clusters along the trajectory sampled by a Monte Carlo (MC) simulation. For example, using this observer, several canonical MC simulations could be carried out at different temperatures and the temperature dependence of the number of nearest neigbhors of a particular species could be accessed with this observer.

Parameters:
  • cluster_space (icet.ClusterSpace) – cluster space to define the clusters to be counted

  • structure (ase.Atoms) – defines the lattice that the observer will work on

  • interval (int) – observation interval during the Monte Carlo simulation

  • orbit_indices (List[int]) – only include orbits up to the orbit with this index (default is to include all orbits)

tag

human readable observer name

Type:

str

interval

the observation interval, defaults to None meaning that if the observer is used in a Monte Carlo simulation, then the Ensemble object will set the interval.

Type:

int

get_cluster_counts(structure)[source]

Counts the number of times different clusters appear in the structure and returns this information as a pandas dataframe.

Parameters:

structure (Atoms) – input atomic structure.

Return type:

DataFrame

get_observable(structure)[source]

Returns the value of the property from a cluster expansion model for a given atomic configuration.

Parameters:

structure (Atoms) – input atomic structure

Return type:

dict

property return_type: type

Data type of the observed data.

ClusterExpansionObserver

class mchammer.observers.ClusterExpansionObserver(cluster_expansion, interval=None, tag='ClusterExpansionObserver')[source]

This class represents a cluster expansion (CE) observer.

A CE observer allows to compute a property described by a CE along the trajectory sampled by a Monte Carlo (MC) simulation. In general this CE differs from the CE that is used to generate the trajectory. For example in a canonical MC simulation the latter would usually represent an energy (total or mixing energy) whereas the former CE(s) could map lattice constant or band gap.

Parameters:
  • cluster_expansion (icet.ClusterExpansion cluster expansion model) – to be used for observation

  • tag (str) – human readable observer name (default: ClusterExpansionObserver)

  • interval (int) – observation interval during the Monte Carlo simulation

tag

name of observer

Type:

str

interval

the observation interval, defaults to None meaning that if the observer is used in a Monte Carlo simulation, then the Ensemble object will set the interval.

Type:

int

get_observable(structure)[source]

Returns the value of the property from a cluster expansion model for a given atomic configuration.

Parameters:

structure (Atoms) – input atomic structure.

Return type:

float

property return_type: type

Data type of the observed data.

ConstituentStrainObserver

class mchammer.observers.ConstituentStrainObserver(constituent_strain, interval=None)[source]

This class represents a constituent strain observer. It allows observation of constituent strain energy separate from the energy calculated by the cluster expansion.

Parameters:
  • constituent_strain (ConstituentStrain) – ConstituentStrain object

  • interval (Optional[int]) – observation interval during the Monte Carlo simulation

tag

human readable observer name

Type:

str

interval

the observation interval. If None the ensemble object will set the interval (if the observer is used in a Monte Carlo simulation)

Type:

int

get_observable(structure)[source]

Returns the constituent strain energy for a given atomic configuration.

Parameters:

structure (Atoms) – input atomic structure

Return type:

dict

property return_type: type

Data type of the observed data.