Ensembles

Canonical ensemble

class mchammer.ensembles.CanonicalEnsemble(structure, calculator, temperature, user_tag=None, boltzmann_constant=8.617330337217213e-05, random_seed=None, dc_filename=None, data_container=None, data_container_write_period=600, ensemble_data_write_interval=None, trajectory_write_interval=None, sublattice_probabilities=None)[source]

Instances of this class allow one to simulate systems in the canonical ensemble (\(N_iVT\)), i.e. at constant temperature (\(T\)), number of atoms of each species (\(N_i\)), and volume (\(V\)).

The probability for a particular state in the canonical ensemble is proportional to the well-known Boltzmann factor,

\[\rho_{\text{C}} \propto \exp [ - E / k_B T ].\]

Since the concentrations or equivalently the number of atoms of each species is held fixed in the canonical ensemble, a trial step must conserve the concentrations. This is accomplished by randomly picking two unlike atoms and swapping their identities. The swap is accepted with probability

\[P = \min \{ 1, \, \exp [ - \Delta E / k_B T ] \},\]

where \(\Delta E\) is the change in potential energy caused by the swap.

The canonical ensemble provides an ideal framework for studying the properties of a system at a specific concentration. Properties such as potential energy or phenomena such as chemical ordering at a specific temperature can conveniently be studied by simulating at that temperature. The canonical ensemble is also a convenient tool for “optimizing” a system, i.e., finding its lowest energy chemical ordering. In practice, this is usually achieved by simulated annealing, i.e. the system is equilibrated at a high temperature, after which the temperature is continuously lowered until the acceptance probability is almost zero. In a well-behaved system, the chemical ordering at that point corresponds to a low-energy structure, possibly the global minimum at that particular concentration.

Parameters:
  • structure (Atoms) – atomic configuration to be used in the Monte Carlo simulation; also defines the initial occupation vector

  • calculator (BaseCalculator) – calculator to be used for calculating the potential changes that enter the evaluation of the Metropolis criterion

  • temperature (float) – temperature \(T\) in appropriate units [commonly Kelvin]

  • boltzmann_constant (float) – Boltzmann constant \(k_B\) in appropriate units, i.e. units that are consistent with the underlying cluster expansion and the temperature units [default: eV/K]

  • user_tag (str) – human-readable tag for ensemble [default: None]

  • random_seed (int) – seed for the random number generator used in the Monte Carlo simulation

  • dc_filename (str) – name of file the data container associated with the ensemble will be written to; if the file exists it will be read, the data container will be appended, and the file will be updated/overwritten

  • data_container_write_period (float) – period in units of seconds at which the data container is written to file; writing periodically to file provides both a way to examine the progress of the simulation and to back up the data [default: 600 s]

  • ensemble_data_write_interval (int) – interval at which data is written to the data container; this includes for example the current value of the calculator (i.e. usually the energy) as well as ensembles specific fields such as temperature or the number of atoms of different species

  • trajectory_write_interval (int) – interval at which the current occupation vector of the atomic configuration is written to the data container.

  • sublattice_probabilities (List[float]) – probability for picking a sublattice when doing a random swap. This should be as long as the number of sublattices and should sum up to 1.

Example

The following snippet illustrate how to carry out a simple Monte Carlo simulation in the canonical ensemble. 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

>>> # 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
>>> structure = prim.repeat(3)
>>> for k in range(5):
>>>     structure[k].symbol = 'Ag'

>>> # set up and run MC simulation
>>> calc = ClusterExpansionCalculator(structure, ce)
>>> mc = CanonicalEnsemble(structure=structure, calculator=calc,
...                        temperature=600,
...                        dc_filename='myrun_canonical.dc')
>>> mc.run(100)  # carry out 100 trial swaps
attach_observer(observer, tag=None)

Attaches an observer to the ensemble.

If the observer does not have an observation interval, then it will be set to the default_interval len(atoms).

Parameters:
  • observer (BaseObserver) – observer instance to attach

  • tag – name used in data container

property boltzmann_constant: float

Boltzmann constant \(k_B\) (see parameters section above)

property calculator: BaseCalculator

calculator attached to the ensemble

property data_container: BaseDataContainer

data container associated with ensemble

do_canonical_swap(sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • sublattice_index (int) – the sublattice the swap will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_sgc_flip(chemical_potentials, sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • chemical_potentials (Dict[int, float]) – chemical potentials used to calculate the potential difference

  • sublattice_index (int) – the sublattice the flip will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_thermodynamic_swap(sublattice_index, lambda_val, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • sublattice_index (int) – the sublattice the swap will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_vcsgc_flip(phis, kappa, sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • phis (Dict[int, float]) – average constraint parameters

  • kappa (float) – parameter that constrains the variance of the concentration

  • sublattice_index (int) – the sublattice the flip will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

property ensemble_parameters: dict

Returns parameters associated with the ensemble.

get_random_sublattice_index(probability_distribution)

Returns a random sublattice index based on the weights of the sublattice.

Parameters:

probability_distribution – probability distributions for the sublattices

Return type:

int

property observer_interval: int

minimum number of steps to run Monte Carlo simulation without interruption for observation

property observers: Dict[str, BaseObserver]
property random_seed: int

seed used to initialize random number generator

run(number_of_trial_steps)

Samples the ensemble for the given number of trial steps.

Parameters:
  • number_of_trial_steps (int) – number of MC trial steps to run in total

  • reset_step – if True the MC trial step counter and the data container will be reset to zero and empty, respectively.

Raises:

TypeError – if number_of_trial_steps is not an int

property step: int

current trial step counter

property structure: Atoms

current configuration (copy)

property sublattices: Sublattices

sublattices for the configuration being sampled

property temperature: float

Current temperature

update_occupations(sites, species)

Updates the occupation vector of the configuration being sampled. This will change the state of the configuration in both the calculator and the configuration manager.

Parameters:
  • sites (List[int]) – indices of sites of the configuration to change

  • species (List[int]) – new occupations (species) by atomic number

Raises:

ValueError – if input lists are not of the same length

property user_tag: str | None

tag used for labeling the ensemble

write_data_container(outfile)

Updates last state of the Monte Carlo simulation and writes data container to file.

Parameters:

outfile (Union[str, bytes]) – file to which to write

Canonical annealing

class mchammer.ensembles.CanonicalAnnealing(structure, calculator, T_start, T_stop, n_steps, cooling_function='exponential', user_tag=None, boltzmann_constant=8.617330337217213e-05, random_seed=None, dc_filename=None, data_container_write_period=600, ensemble_data_write_interval=None, trajectory_write_interval=None, sublattice_probabilities=None)[source]

Instances of this class allow one to carry out simulated annealing in the canonical ensemble, i.e. the temperature is varied in pre-defined fashion while the composition is kept fixed. See mchammer.ensembles.CanonicalEnsemble for more information about the standard canonical ensemble.

The canonical annealing ensemble can be useful, for example, for finding ground states or generating low energy configurations.

The temperature control scheme is selected via the cooling_function keyword argument, while the initial and final temperature are set via the T_start and T_stop arguments. Several pre-defined temperature control schemes are available including linear and exponential. In the latter case the temperature varies logarithmatically as a function of the MC step, emulating the exponential temperature dependence of the atomic exchange rate encountered in many materials. It is also possible to provide a user defined cooling function via the keyword argument. This function must comply with the following function header:

def cooling_function(step, T_start, T_stop, n_steps):
    T = ...  # compute temperature
    return T

Here step refers to the current MC trial step.

Parameters:
  • structure (Atoms) – atomic configuration to be used in the Monte Carlo simulation; also defines the initial occupation vector

  • calculator (BaseCalculator) – calculator to be used for calculating the potential changes that enter the evaluation of the Metropolis criterion

  • T_start (float) – temperature from which the annealing is started

  • T_stop (float) – final temperature for annealing

  • n_steps (int) – number of steps to take in the annealing simulation

  • cooling_function (str/function) – to use the predefined cooling functions provide a string linear or exponential, otherwise provide a function

  • boltzmann_constant (float) – Boltzmann constant \(k_B\) in appropriate units, i.e. units that are consistent with the underlying cluster expansion and the temperature units [default: eV/K]

  • user_tag (str) – human-readable tag for ensemble [default: None]

  • random_seed (int) – seed for the random number generator used in the Monte Carlo simulation

  • dc_filename (str) – name of file the data container associated with the ensemble will be written to; if the file exists it will be read, the data container will be appended, and the file will be updated/overwritten

  • data_container_write_period (float) – period in units of seconds at which the data container is written to file; writing periodically to file provides both a way to examine the progress of the simulation and to back up the data [default: 600 s]

  • ensemble_data_write_interval (int) – interval at which data is written to the data container; this includes for example the current value of the calculator (i.e. usually the energy) as well as ensembles specific fields such as temperature or the number of atoms of different species

  • trajectory_write_interval (int) – interval at which the current occupation vector of the atomic configuration is written to the data container.

  • sublattice_probabilities (List[float]) – probability for picking a sublattice when doing a random swap. This should be as long as the number of sublattices and should sum up to 1.

property T_start: float

Starting temperature

property T_stop: float

Starting temperature

attach_observer(observer, tag=None)

Attaches an observer to the ensemble.

If the observer does not have an observation interval, then it will be set to the default_interval len(atoms).

Parameters:
  • observer (BaseObserver) – observer instance to attach

  • tag – name used in data container

property boltzmann_constant: float

Boltzmann constant \(k_B\) (see parameters section above)

property calculator: BaseCalculator

calculator attached to the ensemble

property data_container: BaseDataContainer

data container associated with ensemble

do_canonical_swap(sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • sublattice_index (int) – the sublattice the swap will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_sgc_flip(chemical_potentials, sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • chemical_potentials (Dict[int, float]) – chemical potentials used to calculate the potential difference

  • sublattice_index (int) – the sublattice the flip will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_thermodynamic_swap(sublattice_index, lambda_val, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • sublattice_index (int) – the sublattice the swap will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_vcsgc_flip(phis, kappa, sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • phis (Dict[int, float]) – average constraint parameters

  • kappa (float) – parameter that constrains the variance of the concentration

  • sublattice_index (int) – the sublattice the flip will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

property ensemble_parameters: dict

Returns parameters associated with the ensemble.

property estimated_ground_state

Structure with lowest observed potential during run

property estimated_ground_state_potential

Lowest observed potential during run

get_random_sublattice_index(probability_distribution)

Returns a random sublattice index based on the weights of the sublattice.

Parameters:

probability_distribution – probability distributions for the sublattices

Return type:

int

property n_steps: int

Number of steps to carry out

property observer_interval: int

minimum number of steps to run Monte Carlo simulation without interruption for observation

property observers: Dict[str, BaseObserver]
property random_seed: int

seed used to initialize random number generator

run()[source]

Runs the annealing.

property step: int

current trial step counter

property structure: Atoms

current configuration (copy)

property sublattices: Sublattices

sublattices for the configuration being sampled

property temperature: float

Current temperature

update_occupations(sites, species)

Updates the occupation vector of the configuration being sampled. This will change the state of the configuration in both the calculator and the configuration manager.

Parameters:
  • sites (List[int]) – indices of sites of the configuration to change

  • species (List[int]) – new occupations (species) by atomic number

Raises:

ValueError – if input lists are not of the same length

property user_tag: str | None

tag used for labeling the ensemble

write_data_container(outfile)

Updates last state of the Monte Carlo simulation and writes data container to file.

Parameters:

outfile (Union[str, bytes]) – file to which to write

Semi-grand canonical ensemble

class mchammer.ensembles.SemiGrandCanonicalEnsemble(structure, calculator, temperature, chemical_potentials, boltzmann_constant=8.617330337217213e-05, user_tag=None, random_seed=None, dc_filename=None, data_container=None, data_container_write_period=600, ensemble_data_write_interval=None, trajectory_write_interval=None, sublattice_probabilities=None)[source]

Instances of this class allow one to simulate systems in the semi-grand canonical (SGC) ensemble (\(N\Delta\mu_i VT\)), i.e. at constant temperature (\(T\)), total number of sites (\(N=\sum_i N_i\)), relative chemical potentials (\(\Delta\mu_i=\mu_i - \mu_1\), where \(i\) denotes the species), and volume (\(V\)).

The probability for a particular state in the SGC ensemble for a \(m\)-component system can be written

\[\rho_{\text{SGC}} \propto \exp\Big[ - \big( E + \sum_{i>1}^m \Delta\mu_i N_i \big) \big / k_B T \Big]\]

with the relative chemical potentials \(\Delta\mu_i = \mu_i - \mu_1\) and species counts \(N_i\). Unlike the canonical ensemble, the number of the respective species (or, equivalently, the concentrations) are allowed to vary in the SGC ensemble. A trial step thus consists of randomly picking an atom and changing its identity with probability

\[P = \min \Big\{ 1, \, \exp \big[ - \big( \Delta E + \sum_i \Delta \mu_i \Delta N_i \big) \big / k_B T \big] \Big\},\]

where \(\Delta E\) is the change in potential energy caused by the swap.

There exists a simple relation between the differences in chemical potential and the canonical free energy \(F\). In a binary system, this relationship reads

\[\Delta \mu = - \frac{1}{N} \frac{\partial F}{\partial c} ( N, V, T, \langle c \rangle).\]

Here \(c\) denotes concentration (\(c=N_i/N\)) and \(\langle c \rangle\) the average concentration observed in the simulation. By recording \(\langle c \rangle\) while gradually changing \(\Delta \mu\), one can thus in principle calculate the difference in canonical free energy between the pure phases (\(c=0\) or \(1\)) and any concentration by integrating \(\Delta \mu\) over that concentration range. In practice this requires that the average recorded concentration \(\langle c \rangle\) varies continuously with \(\Delta \mu\). This is not the case for materials with multiphase regions (such as miscibility gaps), because in such regions \(\Delta \mu\) maps to multiple concentrations. In a Monte Carlo simulation, this is typically manifested by discontinuous jumps in concentration. Such jumps mark the phase boundaries of a multiphase region and can thus be used to construct the phase diagram. To recover the free energy, however, such systems require sampling in other ensembles, such as the variance-constrained semi-grand canonical ensemble.

Parameters:
  • structure (Atoms) – atomic configuration to be used in the Monte Carlo simulation; also defines the initial occupation vector

  • calculator (BaseCalculator) – calculator to be used for calculating the potential changes that enter the evaluation of the Metropolis criterion

  • temperature (float) – temperature \(T\) in appropriate units [commonly Kelvin]

  • chemical_potentials (Dict[str, float]) – chemical potential for each species \(\mu_i\); the key denotes the species, the value specifies the chemical potential in units that are consistent with the underlying cluster expansion

  • boltzmann_constant (float) – Boltzmann constant \(k_B\) in appropriate units, i.e. units that are consistent with the underlying cluster expansion and the temperature units [default: eV/K]

  • user_tag (str) – human-readable tag for ensemble [default: None]

  • random_seed (int) – seed for the random number generator used in the Monte Carlo simulation

  • dc_filename (str) – name of file the data container associated with the ensemble will be written to; if the file exists it will be read, the data container will be appended, and the file will be updated/overwritten

  • data_container_write_period (float) – period in units of seconds at which the data container is written to file; writing periodically to file provides both a way to examine the progress of the simulation and to back up the data [default: 600 s]

  • ensemble_data_write_interval (int) – interval at which data is written to the data container; this includes for example the current value of the calculator (i.e. usually the energy) as well as ensembles specific fields such as temperature or the number of atoms of different species

  • trajectory_write_interval (int) – interval at which the current occupation vector of the atomic configuration is written to the data container.

  • sublattice_probabilities (List[float]) – probability for picking a sublattice when doing a random flip. This should be as long as the number of sublattices and should sum up to 1.

Example

The following snippet illustrate how to carry out a simple Monte Carlo simulation in the semi-canonical ensemble. 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

>>> # 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])

>>> # set up and run MC simulation (T=600 K, delta_mu=0.8 eV/atom)
>>> structure = prim.repeat(3)
>>> calc = ClusterExpansionCalculator(structure, ce)
>>> mc = SemiGrandCanonicalEnsemble(structure=structure, calculator=calc,
...                                temperature=600,
...                                dc_filename='myrun_sgc.dc',
...                                chemical_potentials={'Ag': 0, 'Au': 0.8})
>>> mc.run(100)  # carry out 100 trial swaps
attach_observer(observer, tag=None)

Attaches an observer to the ensemble.

If the observer does not have an observation interval, then it will be set to the default_interval len(atoms).

Parameters:
  • observer (BaseObserver) – observer instance to attach

  • tag – name used in data container

property boltzmann_constant: float

Boltzmann constant \(k_B\) (see parameters section above)

property calculator: BaseCalculator

calculator attached to the ensemble

property chemical_potentials: Dict[int, float]

chemical potentials \(\mu_i\) (see parameters section above)

property data_container: BaseDataContainer

data container associated with ensemble

do_canonical_swap(sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • sublattice_index (int) – the sublattice the swap will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_sgc_flip(chemical_potentials, sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • chemical_potentials (Dict[int, float]) – chemical potentials used to calculate the potential difference

  • sublattice_index (int) – the sublattice the flip will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_thermodynamic_swap(sublattice_index, lambda_val, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • sublattice_index (int) – the sublattice the swap will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_vcsgc_flip(phis, kappa, sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • phis (Dict[int, float]) – average constraint parameters

  • kappa (float) – parameter that constrains the variance of the concentration

  • sublattice_index (int) – the sublattice the flip will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

property ensemble_parameters: dict

Returns parameters associated with the ensemble.

get_random_sublattice_index(probability_distribution)

Returns a random sublattice index based on the weights of the sublattice.

Parameters:

probability_distribution – probability distributions for the sublattices

Return type:

int

property observer_interval: int

minimum number of steps to run Monte Carlo simulation without interruption for observation

property observers: Dict[str, BaseObserver]
property random_seed: int

seed used to initialize random number generator

run(number_of_trial_steps)

Samples the ensemble for the given number of trial steps.

Parameters:
  • number_of_trial_steps (int) – number of MC trial steps to run in total

  • reset_step – if True the MC trial step counter and the data container will be reset to zero and empty, respectively.

Raises:

TypeError – if number_of_trial_steps is not an int

property step: int

current trial step counter

property structure: Atoms

current configuration (copy)

property sublattices: Sublattices

sublattices for the configuration being sampled

property temperature: float

temperature \(T\) (see parameters section above)

update_occupations(sites, species)

Updates the occupation vector of the configuration being sampled. This will change the state of the configuration in both the calculator and the configuration manager.

Parameters:
  • sites (List[int]) – indices of sites of the configuration to change

  • species (List[int]) – new occupations (species) by atomic number

Raises:

ValueError – if input lists are not of the same length

property user_tag: str | None

tag used for labeling the ensemble

write_data_container(outfile)

Updates last state of the Monte Carlo simulation and writes data container to file.

Parameters:

outfile (Union[str, bytes]) – file to which to write

Variance-constrained semi-grand canonical ensemble

class mchammer.ensembles.VCSGCEnsemble(structure, calculator, temperature, phis, kappa, boltzmann_constant=8.617330337217213e-05, user_tag=None, random_seed=None, dc_filename=None, data_container=None, data_container_write_period=600, ensemble_data_write_interval=None, trajectory_write_interval=None, sublattice_probabilities=None)[source]

Instances of this class allow one to simulate systems in the variance-constrained semi-grand canonical (VCSGC) ensemble (\(N\phi\kappa VT\)), i.e. at constant temperature (\(T\)), total number of sites (\(N=\sum_i N_i\)), and two additional dimensionless parameters \(\phi\) and \(\kappa\), which constrain average and variance of the concentration, respectively.

The below examples treat the binary case, but the generalization of to ternaries and higher-order systems is straight-forward. The probability for a particular state in the VCSGC ensemble for a \(2\)-component system can be written

\[\rho_{\text{VCSGC}} \propto \exp\Big[ - E / k_B T + \kappa N ( c_1 + \phi_1 / 2 )^2 \Big],\]

where \(c_1\) represents the concentration of species 1, i.e. \(c_1=N_1/N\). (Please note that the quantities \(\kappa\) and \(\phi\) correspond, respectively, to \(\bar{\kappa}\) and \(\bar{\phi}\) in [SadErh12].). The \(\phi\) may refer to any of the two species. If \(\phi\) is specified for species A, an equivalent simulation can be carried out by specifying \(\phi_B\) as \(-2-\phi_A\). In general, simulations of \(N\)-component systems requires the specification of \(\phi\) for \(N-1\) elements.

Just like the semi-grand canonical ensemble, the VCSGC ensemble allows concentrations to change. A trial step consists of changing the identity of a randomly chosen atom and accepting the change with probability

\[P = \min \{ 1, \, \exp [ - \Delta E / k_B T + \kappa N \Delta c_1 (\phi_1 + \Delta c_1 + 2 c_1 ) ] \}.\]

Note that for a sufficiently large value of \(\kappa\), say 200, the probability density \(\rho_{\text{VCSGC}}\) is sharply peaked around \(c_1=-\phi_1 / 2\). In practice, this means that we can gradually change \(\phi_1\) from (using some margins) \(-2.1\) to \(0.1\) and take the system continuously from \(c_1 = 0\) to \(1\). The parameter \(\kappa\) constrains the fluctuations (or the variance) of the concentration at each value of \(\phi_1\), with higher values of \(\kappa\) meaning less fluctuations. Unlike the semi-grand canonical ensemble, one value of \(\phi_1\) maps to one and only one concentration also in multiphase regions. Since the derivative of the canonical free energy can be expressed in terms of parameters and observables of the VCSGC ensemble,

\[k_B T \kappa ( \phi_1 + 2 \langle c_1 \rangle ) = - \frac{1}{N} \frac{\partial F}{\partial c_1} (N, V, T, \langle c_1 \rangle ),\]

this ensemble allows for thermodynamic integration across multiphase regions. This means that we can construct phase diagrams by directly comparing the free energies of the different phases. This often makes the VCSGC ensemble more convenient than the semi-grand canonical ensemble when simulating materials with multiphase regions, such as alloys with miscibility gaps.

When using the VCSGC ensemble, please cite Sadigh, B. and Erhart, P., Phys. Rev. B 86, 134204 (2012) [SadErh12].

Parameters:
  • structure (Atoms) – atomic configuration to be used in the Monte Carlo simulation; also defines the initial occupation vector

  • calculator (BaseCalculator) – calculator to be used for calculating the potential changes that enter the evaluation of the Metropolis criterion

  • temperature (float) – temperature \(T\) in appropriate units [commonly Kelvin]

  • phis (Dict[str, float]) – average constraint parameters \(\phi_i\); the key denotes the species; for a N-component sublattice, there should be N - 1 different \(\phi_i\) (referred to as \(\bar{\phi}\) in [SadErh12])

  • kappa (float) – parameter that constrains the variance of the concentration (referred to as \(\bar{\kappa}\) in [SadErh12])

  • boltzmann_constant (float) – Boltzmann constant \(k_B\) in appropriate units, i.e. units that are consistent with the underlying cluster expansion and the temperature units [default: eV/K]

  • user_tag (str) – human-readable tag for ensemble [default: None]

  • random_seed (int) – seed for the random number generator used in the Monte Carlo simulation

  • dc_filename (str) – name of file the data container associated with the ensemble will be written to; if the file exists it will be read, the data container will be appended, and the file will be updated/overwritten

  • data_container_write_period (float) – period in units of seconds at which the data container is written to file; writing periodically to file provides both a way to examine the progress of the simulation and to back up the data [default: 600 s]

  • ensemble_data_write_interval (int) – interval at which data is written to the data container; this includes for example the current value of the calculator (i.e. usually the energy) as well as ensembles specific fields such as temperature or the number of atoms of different species

  • trajectory_write_interval (int) – interval at which the current occupation vector of the atomic configuration is written to the data container.

  • sublattice_probabilities (List[float]) – probability for picking a sublattice when doing a random flip. The list should be as long as the number of sublattices and should sum up to 1.

Example

The following snippet illustrate how to carry out a simple Monte Carlo simulation in the variance-constrained semi-canonical ensemble. 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

>>> # 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])

>>> # set up and run MC simulation
>>> structure = prim.repeat(3)
>>> calc = ClusterExpansionCalculator(structure, ce)
>>> phi = 0.6
>>> mc = VCSGCEnsemble(structure=structure, calculator=calc,
...                   temperature=600,
...                   dc_filename='myrun_vcsgc.dc',
...                   phis={'Au': phi},
...                   kappa=200)
>>> mc.run(100)  # carry out 100 trial swaps
attach_observer(observer, tag=None)

Attaches an observer to the ensemble.

If the observer does not have an observation interval, then it will be set to the default_interval len(atoms).

Parameters:
  • observer (BaseObserver) – observer instance to attach

  • tag – name used in data container

property boltzmann_constant: float

Boltzmann constant \(k_B\) (see parameters section above)

property calculator: BaseCalculator

calculator attached to the ensemble

property data_container: BaseDataContainer

data container associated with ensemble

do_canonical_swap(sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • sublattice_index (int) – the sublattice the swap will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_sgc_flip(chemical_potentials, sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • chemical_potentials (Dict[int, float]) – chemical potentials used to calculate the potential difference

  • sublattice_index (int) – the sublattice the flip will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_thermodynamic_swap(sublattice_index, lambda_val, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • sublattice_index (int) – the sublattice the swap will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_vcsgc_flip(phis, kappa, sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • phis (Dict[int, float]) – average constraint parameters

  • kappa (float) – parameter that constrains the variance of the concentration

  • sublattice_index (int) – the sublattice the flip will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

property ensemble_parameters: dict

Returns parameters associated with the ensemble.

get_random_sublattice_index(probability_distribution)

Returns a random sublattice index based on the weights of the sublattice.

Parameters:

probability_distribution – probability distributions for the sublattices

Return type:

int

property kappa: float

kappa \(\bar{\kappa}\) constrain parameter (see parameters section above)

property observer_interval: int

minimum number of steps to run Monte Carlo simulation without interruption for observation

property observers: Dict[str, BaseObserver]
property phis: Dict[int, float]

phis \(\phi_i\) for all species but one (referred to as \(\bar{\phi}\) in [SadErh12])

property random_seed: int

seed used to initialize random number generator

run(number_of_trial_steps)

Samples the ensemble for the given number of trial steps.

Parameters:
  • number_of_trial_steps (int) – number of MC trial steps to run in total

  • reset_step – if True the MC trial step counter and the data container will be reset to zero and empty, respectively.

Raises:

TypeError – if number_of_trial_steps is not an int

property step: int

current trial step counter

property structure: Atoms

current configuration (copy)

property sublattices: Sublattices

sublattices for the configuration being sampled

property temperature: float

temperature \(T\) (see parameters section above)

update_occupations(sites, species)

Updates the occupation vector of the configuration being sampled. This will change the state of the configuration in both the calculator and the configuration manager.

Parameters:
  • sites (List[int]) – indices of sites of the configuration to change

  • species (List[int]) – new occupations (species) by atomic number

Raises:

ValueError – if input lists are not of the same length

property user_tag: str | None

tag used for labeling the ensemble

write_data_container(outfile)

Updates last state of the Monte Carlo simulation and writes data container to file.

Parameters:

outfile (Union[str, bytes]) – file to which to write

Hybrid ensemble

class mchammer.ensembles.HybridEnsemble(structure, calculator, temperature, ensemble_specs, probabilities=None, boltzmann_constant=8.617330337217213e-05, user_tag=None, random_seed=None, dc_filename=None, data_container=None, data_container_write_period=600, ensemble_data_write_interval=None, trajectory_write_interval=None)[source]

Instances of this class allows one to combine multiple ensembles. In particular, a dictionary should be provided for each ensemble, which must include the type (ensemble) as well as the index of the sublattice (sublattice_index). In addition, it is possible to provide a list of allowed symbols (allowed_symbols), which must represent a subset of the elements that can occupy the sites on the specified sublattice. Note that additional arguments are required for the SGC and VCSGC ensembles, namely chemical potentials (chemical_potentials) for the former and constraint parameters (phis and kappa) for the latter. For more detailed information regarding the different ensembles, please see CanonicalEnsemble, SemiGrandCanonicalEnsemble, and VCSGCEnsemble.

This class is particularly useful for effectively sampling complex multi-component systems with several active sublattices, in which case different ensembles can be defined for each of the latter. The fact that it is possible to set the allowed chemical symbols means that one can vary the concentrations of a few selected species, with the help of a VCSGC or semi-grand canonical ensemble, while still allowing swaps between any two sites, using a canonical ensemble (see also the below example). It is advisable to carefully consider how to define the ensemble probabilities. By default the ensembles are weighted by the sizes of the corresponding sublattices, which should give suitable probabilities in most cases. As is shown in the example below, it might be prudent to provide different values if allowed symbols are provided as well as for cases where there are several ensembles that are active on different sublattices.

Parameters:
  • structure (Atoms) – atomic configuration to be used in the Monte Carlo simulation; also defines the initial occupation vector

  • calculator (BaseCalculator) – calculator to be used for calculating the potential changes that enter the evaluation of the Metropolis criterion

  • temperature (float) – temperature \(T\) in appropriate units [commonly Kelvin]

  • ensemble_specs (List[Dict]) –

    A list of dictionaries, which should contain the following items:

    • ’ensemble’, which could be either ‘vcsgc’; ‘semi-grand’ or ‘canonical’, lowercase and uppercase letters or any combination thereof are accepted (required)

    • ’sublattice_index’, index for the sublattice of interest (required)

    • ’allowed_symbols’, list of allowed chemical symbols (default: read from ClusterSpace)

    • ’chemical_potentials’, a dictionary of chemical potentials for each species \(\mu_i\); the key denotes the species, the value :specifies the chemical potential in units that are :consistent with the underlying cluster expansion (only :applicable and required for SGC ensembles)

    • ’phis ‘, dictionary with average constraint parameters ‘\(\phi_i\); the key denotes the species; for a N-component sublattice, there should be N - 1 different phi_i (referred to as \(\bar{\phi}\) in [SadErh12]; only applicable and required for VCSGC ensembles, see also VCSGCEnsemble)

    • ’kappa’, parameter that constrains the variance of the ‘concentration (referred to as \(\bar{\kappa}\) in [SadErh12]; only applicable :and required for VCSGC ensembles)

  • probabilities (List[float]) – list of floats with the probabilities for choosing a particular ensemble with the same length as ensemble specs. If left unspecified the probabilties are weighted by the sizes of the associated sublattices

  • boltzmann_constant (float) – Boltzmann constant \(k_B\) in appropriate units, i.e. units that are consistent with the underlying cluster expansion and the temperature units [default: eV/K]

  • user_tag (str) – human-readable tag for ensemble [default: None]

  • random_seed (int) – seed for the random number generator used in the Monte Carlo simulation

  • dc_filename (str) – name of file the data container associated with the ensemble will be written to; if the file exists it will be read, the data container will be appended, and the file will be updated/overwritten

  • data_container_write_period (float) – period in units of seconds at which the data container is written to file; writing periodically to file provides both a way to examine the progress of the simulation and to back up the data [default: 600 s]

  • ensemble_data_write_interval (int) – interval at which data is written to the data container; this includes for example the current value of the calculator (i.e. usually the energy) as well as ensembles specific fields such as temperature or the number of atoms of different species

  • trajectory_write_interval (int) – interval at which the current occupation vector of the atomic configuration is written to the data container.

Example

The following snippet illustrates how to carry out a simple Monte Carlo simulation using a combination of one canonical and one VCSGC ensemble. Specifically, the concentration of one species (Au) is kept constant while the others (Ag and Pd) are varied, while swaps are still allowed. 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

>>> # 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', 'Pd'])
>>> ce = ClusterExpansion(
...    cs, [0, 0, 0, 0.1, 0.1, 0.1, -0.02, -0.02, -0.02])

# define structure object
>>> structure = prim.repeat(3)
>>> for i, atom in enumerate(structure):
>>>    if i % 2 == 0:
>>>        atom.symbol = 'Ag'
>>>    elif i % 3 == 0:
>>>        atom.symbol = 'Pd'

# the default probabilities for this case would be [0.5, 0.5], but
# since the VCSGC ensemble only performs flips on a subset of all
# sites on the sublattice, namely those originally occupied by Ag
# and Pd atoms, specific values will be provided
>>> weights = [len(structure),
...           len([s for s in structure.get_chemical_symbols() if s != 'Au'])]
>>> norm = sum(weights)
>>> probabilities = [w / norm for w in weights]

# set up and run MC simulation
>>> calc = ClusterExpansionCalculator(structure, ce)
>>> ensemble_specs = [
...    {'ensemble': 'canonical', 'sublattice_index': 0},
...    {'ensemble': 'vcsgc', 'sublattice_index': 0,
...     'phis': {'Ag': -0.2}, 'kappa': 200,
...     'allowed_symbols':['Ag', 'Pd']}]
>>> mc = HybridEnsemble(structure=structure, calculator=calc,
...                    ensemble_specs=ensemble_specs,
...                    temperature=600, probabilities=probabilities,
...                    dc_filename='myrun_hybrid.dc')
>>> mc.run(100)  # carry out 100 trial steps
attach_observer(observer, tag=None)

Attaches an observer to the ensemble.

If the observer does not have an observation interval, then it will be set to the default_interval len(atoms).

Parameters:
  • observer (BaseObserver) – observer instance to attach

  • tag – name used in data container

property boltzmann_constant: float

Boltzmann constant \(k_B\) (see parameters section above)

property calculator: BaseCalculator

calculator attached to the ensemble

property data_container: BaseDataContainer

data container associated with ensemble

do_canonical_swap(sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • sublattice_index (int) – the sublattice the swap will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_sgc_flip(chemical_potentials, sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • chemical_potentials (Dict[int, float]) – chemical potentials used to calculate the potential difference

  • sublattice_index (int) – the sublattice the flip will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_thermodynamic_swap(sublattice_index, lambda_val, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • sublattice_index (int) – the sublattice the swap will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_vcsgc_flip(phis, kappa, sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • phis (Dict[int, float]) – average constraint parameters

  • kappa (float) – parameter that constrains the variance of the concentration

  • sublattice_index (int) – the sublattice the flip will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

property ensemble_parameters: dict

Returns parameters associated with the ensemble.

get_random_sublattice_index(probability_distribution)

Returns a random sublattice index based on the weights of the sublattice.

Parameters:

probability_distribution – probability distributions for the sublattices

Return type:

int

property observer_interval: int

minimum number of steps to run Monte Carlo simulation without interruption for observation

property observers: Dict[str, BaseObserver]
property probabilities: List[float]

Ensemble propabilities

property random_seed: int

seed used to initialize random number generator

run(number_of_trial_steps)

Samples the ensemble for the given number of trial steps.

Parameters:
  • number_of_trial_steps (int) – number of MC trial steps to run in total

  • reset_step – if True the MC trial step counter and the data container will be reset to zero and empty, respectively.

Raises:

TypeError – if number_of_trial_steps is not an int

property step: int

current trial step counter

property structure: Atoms

current configuration (copy)

property sublattices: Sublattices

sublattices for the configuration being sampled

property temperature: float

Current temperature

property trial_steps_per_ensemble: Dict[str, int]

Number of Monte Carlo trial steps for each ensemble

update_occupations(sites, species)

Updates the occupation vector of the configuration being sampled. This will change the state of the configuration in both the calculator and the configuration manager.

Parameters:
  • sites (List[int]) – indices of sites of the configuration to change

  • species (List[int]) – new occupations (species) by atomic number

Raises:

ValueError – if input lists are not of the same length

property user_tag: str | None

tag used for labeling the ensemble

write_data_container(outfile)

Updates last state of the Monte Carlo simulation and writes data container to file.

Parameters:

outfile (Union[str, bytes]) – file to which to write

Wang-Landau ensemble

For analysis functions see here.

class mchammer.ensembles.WangLandauEnsemble(structure, calculator, energy_spacing, energy_limit_left=None, energy_limit_right=None, trial_move='swap', fill_factor_limit=1e-06, flatness_check_interval=None, flatness_limit=0.8, window_search_penalty=2.0, user_tag=None, dc_filename=None, data_container=None, random_seed=None, data_container_write_period=600, ensemble_data_write_interval=None, trajectory_write_interval=None, sublattice_probabilities=None)[source]

Instances of this class allow one to sample a system using the Wang-Landau (WL) algorithm, see Phys. Rev. Lett. 86, 2050 (2001) [WanLan01a]. The WL algorithm enables one to acquire the density of states (DOS) as a function of energy, from which one can readily calculate many thermodynamic observables as a function of temperature. To this end, the WL algorithm accumulates both the microcanonical entropy \(S(E)\) and a histogram \(H(E)\) on an energy grid with a predefined spacing (energy_spacing).

The algorithm is initialized as follows.

  1. Generate an initial configuration.

  2. Initialize counters for the microcanonical entropy \(S(E)\) and the histogram \(H(E)\) to zero.

  3. Set the fill factor \(f=1\).

It then proceeds as follows.

  1. Propose a new configuration (see trial_move).

  2. Accept or reject the new configuration with probability

    \[P = \min \{ 1, \, \exp [ S(E_\mathrm{new}) - S(E_\mathrm{cur}) ] \},\]

    where \(E_\mathrm{cur}\) and \(E_\mathrm{new}\) are the energies of the current and new configurations, respectively.

  3. Update the microcanonical entropy \(S(E)\leftarrow S(E) + f\) and histogram \(H(E) \leftarrow H(E) + 1\) where \(E\) is the energy of the system at the end of the move.

  4. Check the flatness of the histogram \(H(E)\). If \(H(E) > \chi \langle H(E)\rangle\,\forall E\) reset the histogram \(H(E) = 0\) and reduce the fill factor \(f \leftarrow f / 2\). The parameter \(\chi\) is set via flatness_limit.

  5. If \(f\) is smaller than fill_factor_limit terminate the loop, otherwise return to 1.

The microcanonical entropy \(S(E)\) and the histogram along with related information are written to the data container every time \(f\) is updated. Using the density \(\rho(E) = \exp S(E)\) one can then readily compute various thermodynamic quantities, including, e.g., the average energy:

\[\left<E\right> = \frac{\sum_E E \rho(E) \exp(-E / k_B T)}{ \sum_E \rho(E) \exp(-E / k_B T)}\]
Parameters:
  • structure (Atoms) – atomic configuration to be used in the Wang-Landau simulation; also defines the initial occupation vector

  • calculator (BaseCalculator) – calculator to be used for calculating potential changes

  • trial_move (str) – One can choose between two different trial moves for generating new configurations. In a ‘swap’ move two sites are selected and their occupations are swapped; in a ‘flip’ move one site is selected and its occupation is flipped to a different species. While ‘swap’ moves conserve the concentrations of the species in the system, ‘flip’ moves allow one in principle to sample the full composition space.

  • energy_spacing (float) – defines the bin size of the energy grid on which the microcanonical entropy \(S(E)\), and thus the density \(\exp S(E)\), is evaluated; the spacing should be small enough to capture the features of the density of states; too small values will, however, render the convergence very tedious if not impossible

  • energy_limit_left (float) – defines the lower limit of the energy range within which the microcanonical entropy \(S(E)\) will be sampled. By default (None) no limit is imposed. Setting limits can be useful if only a part of the density of states is required.

  • energy_limit_right (float) – defines the upper limit of the energy range within which the microcanonical entropy \(S(E)\) will be sampled. By default (None) no limit is imposed. Setting limits can be useful if only a part of the density of states is required.

  • fill_factor_limit (float) – If the fill_factor \(f\) falls below this value, the WL sampling loop is terminated.

  • flatness_check_interval (int) – For computational efficiency the flatness condition is only evaluated every flatness_check_interval-th trial step. By default (None) flatness_check_interval is set to 1000 times the number of sites in structure, i.e. 1000 Monte Carlo sweeps.

  • flatness_limit (float) – The histogram \(H(E)\) is deemed sufficiently flat if \(H(E) > \chi \left<H(E)\right>\,\forall E\). flatness_limit sets the parameter \(\chi\).

  • window_search_penalty (float) – If energy_limit_left and/or energy_limit_right have been provided, a modified acceptance probability, \(P=\min\{1,\,\exp[C_\mathrm{WSP}(d_\mathrm{new}- d_\mathrm{cur})]\}\), will be used until a configuration is found within the interval of interest. This parameter, specifically, corresponds to \(C_\mathrm{WSP}\), which controls how strongly moves that lead to an increase in the distance, i.e. difference in energy divided by the energy spacing, to the energy window (\(d_\mathrm{new}> d_\mathrm{cur}\)) should be penalized. A higher value leads to a lower acceptance probability for such moves.

  • user_tag (str) – human-readable tag for ensemble [default: None]

  • dc_filename (str) – name of file the data container associated with the ensemble will be written to; if the file exists it will be read, the data container will be appended, and the file will be updated/overwritten

  • random_seed (int) – seed for the random number generator used in the Monte Carlo simulation

  • ensemble_data_write_interval (int) – interval at which data is written to the data container; this includes for example the current value of the calculator (i.e. usually the energy) as well as ensembles specific fields such as temperature or the number of atoms of different species

  • data_container_write_period (float) – period in units of seconds at which the data container is written to file; writing periodically to file provides both a way to examine the progress of the simulation and to back up the data [default: 600 s]

  • trajectory_write_interval (int) – interval at which the current occupation vector of the atomic configuration is written to the data container.

  • sublattice_probabilities (List[float]) – probability for picking a sublattice when doing a random swap. The list must contain as many elements as there are sublattices and it needs to sum up to 1.

Example

The following snippet illustrates how to carry out a Wang-Landau simulation. For the purpose of demonstration, the parameters of the cluster expansion are set to obtain a two-dimensional square Ising model, one of the systems studied in the original work by Wang and Landau:

>>> from ase import Atoms
>>> from icet import ClusterExpansion, ClusterSpace
>>> from mchammer.calculators import ClusterExpansionCalculator
>>> from mchammer.ensembles import WangLandauEnsemble

>>> # prepare cluster expansion
>>> prim = Atoms('Au', positions=[[0, 0, 0]], cell=[1, 1, 10], pbc=True)
>>> cs = ClusterSpace(prim, cutoffs=[1.1], chemical_symbols=['Ag', 'Au'])
>>> ce = ClusterExpansion(cs, [0, 0, 2])

>>> # prepare initial configuration
>>> structure = prim.repeat((4, 4, 1))
>>> for k in range(8):
...     structure[k].symbol = 'Ag'

>>> # set up and run Wang-Landau simulation
>>> calculator = ClusterExpansionCalculator(structure, ce)
>>> mc = WangLandauEnsemble(structure=structure,
...                         calculator=calculator,
...                         energy_spacing=1,
...                         dc_filename='ising_2d_run.dc')
>>> mc.run(number_of_trial_steps=len(structure)*1000)  # in practice one requires more steps
attach_observer(observer, tag=None)

Attaches an observer to the ensemble.

If the observer does not have an observation interval, then it will be set to the default_interval len(atoms).

Parameters:
  • observer (BaseObserver) – observer instance to attach

  • tag – name used in data container

property calculator: BaseCalculator

calculator attached to the ensemble

property converged: bool | None

True if convergence has been achieved

property data_container: BaseDataContainer

data container associated with ensemble

property ensemble_parameters: dict

Returns parameters associated with the ensemble.

property fill_factor: float

current value of fill factor

property fill_factor_history: Dict[int, float]

evolution of the fill factor in the Wang-Landau algorithm (key=MC trial step, value=fill factor)

property fill_factor_limit: float

If the fill_factor \(f\) falls below this value, the Wang-Landau sampling loop is terminated.

property flatness_check_interval: int

number of MC trial steps between checking the flatness condition

property flatness_limit: float

The histogram \(H(E)\) is deemed sufficiently flat if \(H(E) > \chi \left<H(E)\right>\,\forall E\) where flatness_limit sets the parameter \(\chi\).

get_random_sublattice_index(probability_distribution)

Returns a random sublattice index based on the weights of the sublattice.

Parameters:

probability_distribution – probability distributions for the sublattices

Return type:

int

property observer_interval: int

minimum number of steps to run Monte Carlo simulation without interruption for observation

property observers: Dict[str, BaseObserver]
property random_seed: int

seed used to initialize random number generator

run(number_of_trial_steps)[source]

Samples the ensemble for the given number of trial steps.

Parameters:
  • number_of_trial_steps (int) – maximum number of MC trial steps to run in total (the run will terminate earlier if fill_factor_limit is reached)

  • reset_step – if True the MC trial step counter and the data container will be reset to zero and empty, respectively.

Raises:

TypeError – if number_of_trial_steps is not an int

property step: int

current trial step counter

property structure: Atoms

current configuration (copy)

property sublattices: Sublattices

sublattices for the configuration being sampled

update_occupations(sites, species)

Updates the occupation vector of the configuration being sampled. This will change the state of the configuration in both the calculator and the configuration manager.

Parameters:
  • sites (List[int]) – indices of sites of the configuration to change

  • species (List[int]) – new occupations (species) by atomic number

Raises:

ValueError – if input lists are not of the same length

property user_tag: str | None

tag used for labeling the ensemble

write_data_container(outfile)[source]

Updates last state of the Wang-Landau simulation and writes DataContainer to file.

Parameters:

outfile (Union[str, bytes]) – file to which to write

mchammer.ensembles.wang_landau_ensemble.get_bins_for_parallel_simulations(n_bins, energy_spacing, minimum_energy, maximum_energy, overlap=4, bin_size_exponent=1.0)[source]

Generates a list of energy bins (lower and upper bound) suitable for parallel Wang-Landau simulations. For the latter, the energy range is split up into a several bins (n_bins). Each bin is then sampled in a separate Wang-Landau simulation. Once the density of states in the individual bins has been converged the total density of states can be constructed by patching the segments back together. To this end, one requires some over overlap between the segments (overlap).

The function returns a list of tuples. Each tuple provides the lower (energy_limit_left) and upper (energy_limit_right) bound of one bin, which are then to be used to set energy_limit_left and energy_limit_right when initializing a WangLandauEnsemble instance.

N.B.: The left-most/right-most bin has no lower/upper bound (set to None).

Parameters:
  • n_bins (int) – number of bins

  • energy_spacing (float) – defines the bin size of the energy grid used by the Wang-Landau simulation, see WangLandauEnsemble for details

  • minimum_energy (float) – an estimate for the lowest energy to be encountered in this system

  • maximum_energy (float) – an estimate for the highest energy to be encountered in this system

  • overlap (int) – amount of overlap between bins in units of energy_spacing

  • bin_size_exponent (float) – Expert option: This parameter allows one to generate a non-uniform distribution of bin sizes. If bin_size_exponent is smaller than one bins at the lower and upper end of the energy range (specified via minimum_energy and maximum_energy) will be shrunk relative to the bins in the middle of the energy range. In principle this can be used one to achieve a more even distribution of computational load between the individual Wang-Landau simulations.

Return type:

List[Tuple[float, float]]

Thermodynamic-integration ensemble

For analysis functions see here.

class mchammer.ensembles.ThermodynamicIntegrationEnsemble(structure, calculator, temperature, n_steps, forward, user_tag=None, boltzmann_constant=8.617330337217213e-05, random_seed=None, dc_filename=None, data_container=None, data_container_write_period=600, ensemble_data_write_interval=None, trajectory_write_interval=None, sublattice_probabilities=None)[source]

Instances of this class allow one to find the free energy of the system. Todo this we use the canonncal ensemble with a modified Hamiltonian,

\[H(\lambda) = (1 - \lambda) H_{A} + \lambda H_{B}\]

The Hamiltonian is then sampled continuously from \(\lambda=0\) to \(\lambda=1\). \(H_{B}\) is your cluster expansion and \(H_{A}=0\), is a completely disordered system, with free energy given by the ideal mixing entropy.

The free energy, A, of system B is then given by:

\[A_{B} = A_{A} + \int_{0}^{1} \left\langle\frac{\mathrm{d}H(\lambda)} {\mathrm{d}\lambda}\right\rangle_{H} \mathrm{d}\lambda\]

and since \(A_{A}\) is known it is easy to compute \(A_{B}\)

$lambda$ is parametrized as,

\[\lambda(x) = x^5(70x^4 - 315x^3 + 540x^2 - 420x + 126)\]

where $x = step / (n_steps - 1)$

Parameters:
  • structure (Atoms) – atomic configuration to be used in the Monte Carlo simulation; also defines the initial occupation vector

  • calculator (BaseCalculator) – calculator to be used for calculating the potential changes that enter the evaluation of the Metropolis criterion

  • temperature (float) – temperature \(T\) in appropriate units [commonly Kelvin]

  • n_lambdas (int) – the number of \(\lambda\) values to be sampled between 0 and 1

  • forward (bool) – If this is set to True, the simulation runs from \(H_A\) to \(H_B\), otherwise it runs from \(H_B\) to \(H_A\). \(H_B\) is the cluster expansion and \(H_A = 0\), is the fully disordered system.

  • boltzmann_constant (float) – Boltzmann constant \(k_B\) in appropriate units, i.e. units that are consistent with the underlying cluster expansion and the temperature units [default: eV/K]

  • user_tag (str) – human-readable tag for ensemble [default: None]

  • random_seed (int) – seed for the random number generator used in the Monte Carlo simulation

  • dc_filename (str) – name of file the data container associated with the ensemble will be written to; if the file exists it will be read, the data container will be appended, and the file will be updated/overwritten

  • data_container_write_period (float) – period in units of seconds at which the data container is written to file; writing periodically to file provides both a way to examine the progress of the simulation and to back up the data [default: 600 s]

  • ensemble_data_write_interval (int) – interval at which data is written to the data container; this includes for example the current value of the calculator (i.e. usually the energy) as well as ensembles specific fields such as temperature or the number of atoms of different species

  • trajectory_write_interval (int) – interval at which the current occupation vector of the atomic configuration is written to the data container.

  • sublattice_probabilities (List[float]) – probability for picking a sublattice when doing a random swap. This should be as long as the number of sublattices and should sum up to 1.

Example

The following snippet illustrate how to carry out a simple thermodynamic integration. 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

>>> # 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
>>> structure = prim.repeat(3)
>>> for k in range(5):
>>>     structure[k].symbol = 'Ag'

>>> # set up and run MC simulation
>>> calc = ClusterExpansionCalculator(structure, ce)
>>> mc = ThermodynamicIntegrationEnsemble(structure=structure, calculator=calc,
...                                       temperature=600,
...                                       n_steps=100000,
...                                       forward=True,
...                                       dc_filename='myrun_thermodynamic_integration.dc')
>>> mc.run()
attach_observer(observer, tag=None)

Attaches an observer to the ensemble.

If the observer does not have an observation interval, then it will be set to the default_interval len(atoms).

Parameters:
  • observer (BaseObserver) – observer instance to attach

  • tag – name used in data container

property boltzmann_constant: float

Boltzmann constant \(k_B\) (see parameters section above)

property calculator: BaseCalculator

calculator attached to the ensemble

property data_container: BaseDataContainer

data container associated with ensemble

do_canonical_swap(sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • sublattice_index (int) – the sublattice the swap will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_sgc_flip(chemical_potentials, sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • chemical_potentials (Dict[int, float]) – chemical potentials used to calculate the potential difference

  • sublattice_index (int) – the sublattice the flip will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_thermodynamic_swap(sublattice_index, lambda_val, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • sublattice_index (int) – the sublattice the swap will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

do_vcsgc_flip(phis, kappa, sublattice_index, allowed_species=None)

Carries out one Monte Carlo trial step.

Parameters:
  • phis (Dict[int, float]) – average constraint parameters

  • kappa (float) – parameter that constrains the variance of the concentration

  • sublattice_index (int) – the sublattice the flip will act on

  • allowed_species (Optional[List[int]]) – list of atomic numbers for allowed species

Return type:

Returns 1 or 0 depending on if trial move was accepted or rejected

property ensemble_parameters: dict

Returns parameters associated with the ensemble.

get_random_sublattice_index(probability_distribution)

Returns a random sublattice index based on the weights of the sublattice.

Parameters:

probability_distribution – probability distributions for the sublattices

Return type:

int

property n_steps: int
property observer_interval: int

minimum number of steps to run Monte Carlo simulation without interruption for observation

property observers: Dict[str, BaseObserver]
property random_seed: int

seed used to initialize random number generator

run()[source]

Runs the thermodynamic integration

property step: int

current trial step counter

property structure: Atoms

current configuration (copy)

property sublattices: Sublattices

sublattices for the configuration being sampled

property temperature: float

Current temperature

update_occupations(sites, species)

Updates the occupation vector of the configuration being sampled. This will change the state of the configuration in both the calculator and the configuration manager.

Parameters:
  • sites (List[int]) – indices of sites of the configuration to change

  • species (List[int]) – new occupations (species) by atomic number

Raises:

ValueError – if input lists are not of the same length

property user_tag: str | None

tag used for labeling the ensemble

write_data_container(outfile)

Updates last state of the Monte Carlo simulation and writes data container to file.

Parameters:

outfile (Union[str, bytes]) – file to which to write