Monte Carlo Observers
- class clease.montecarlo.observers.AcceptanceRate[source]
Observer that tracks the fraction of monte carlo steps that is accepted
- get_averages() Dict[str, float] [source]
Return dictionary with the rate such that it is added to thermodynaic quantities
- property rate: float
Acceptance rate
- class clease.montecarlo.observers.ConcentrationObserver(atoms: Atoms, element: str)[source]
Observer that can be attached to a MC run, to track the concenctration of a particular element. This observer has to be executed on every MC step.
Parameters:
- atoms: Atoms object
Atoms object used for MC
- element: str
The element that should be tracked
- calculate_from_scratch(atoms: Atoms) float [source]
Calculate the concentration of the element in the atoms object.
- interval_ok(interval: int) bool [source]
Every step must be observed, as otherwise we’d miss updates, and the concentration becomes incorerct.
- new_concentration(system_changes: Sequence[SystemChange]) float [source]
Calculate the new consentration after the changes.
- class clease.montecarlo.observers.CorrelationFunctionObserver(calc, names=None)[source]
Track the history of the correlation function.
Parameters:
- calc: clease.calculators.Clease
Clease calculator
- names: list
List with correlation functions to track. If None, all correlation functions are tracked.
- class clease.montecarlo.observers.DiffractionObserver(atoms=None, k_vector=(), active_symbols=(), all_symbols=(), name='reflection1')[source]
Trace the reflection intensity.
Parameters:
- atoms: Atoms
Atoms object used in Monte Carlo
- k_vector: list
Fourier reflection to be traced
- active_symbols: list
List of symbols that reflects
- all_symbols: list
List of all symbols in the simulation
- name: str
Name of the DiffractionObserver (users are given the freedom to set names because they can attach multiple DiffractionObserver instances)
Example:
Consider a system where Al, Mg and Si occupy FCC lattice sites. We want to trace the occurence of Mg layers that are separated by a distance 3*a where a is the lattice parameter. We further assume that the y-axis is normal to the planes we want to trace. In that case, we specify the variables as
>>> from ase.build import bulk >>> import numpy as np >>> a = 4.05 >>> atoms = bulk('Al', crystalstructure='fcc', a=a) >>> k_vector = [0, 2.0*np.pi/(3*a), 0] >>> active_elements = ['Mg'] >>> all_symbols = ['Al', 'Mg', 'Si']
If we do not wish to distinguish Mi and Si (we do not distiguish Mg layer, Si layer or a mixture of the two) the active_elements is changed to
>>> active_elements = ['Mg', 'Si']
- class clease.montecarlo.observers.EnergyEvolution(mc, ignore_reset=False)[source]
Trace the evolution of energy.
- class clease.montecarlo.observers.EnergyPlotUpdater(energy_obs=None, graph=None, mean_plot=None)[source]
- class clease.montecarlo.observers.EntropyProductionRate(buffer_length: int = 10000, logfile: str | Path = 'epr.txt')[source]
Tracks entropy production rate (EPR) using a Gallavotti-Cohen functional.
EPR = 1/N sum_{i=0}^N ln P(i -> j)/P(j -> i)
N is the number of steps of a path and P(i -> j) is the probability of going from state i to state j. The expression is exact in the limit N -> infty. However, this class tracks the terms inside the sum and write them to file. To calculate the time evolution of EPR one can use a windowed average of the resulting data.
References:
- [1] Gourgoulias, Konstantinos, Markos A. Katsoulakis, and Luc Rey-Bellet.
“Information criteria for quantifying loss of reversibility in parallelized KMC.” Journal of Computational Physics 328 (2017): 438-454.
- Parameters:
buffer_length – Length of buffer used to temporarily store the terms in the sum in memory. When the buffer is full, it is flushed to a text file.
logfile – Filename of the file used when the buffer is flushed
- class clease.montecarlo.observers.LowestEnergyStructure(atoms: Atoms, track_cf: bool = False, verbose: bool = False)[source]
Track the lowest energy state visited during an MC run.
- atoms: Atoms object
Atoms object used in Monte Carlo
- track_cf: bool
Whether to keep a copy of the correlation functions for the emin structure. If enabled, this will be stored in the
lowest_energy_cf
variable. Defaults to False.- verbose: bool
If True, progress messages will be printed
- property emin_results: dict
The results dictionary of the lowest energy atoms
- property energy
The energy of the current atoms object (not the emin energy)
- class clease.montecarlo.observers.MCObserver[source]
Base class for all MC observers.
Child observers should override the observe_step() method.
- calculate_from_scratch(atoms: Atoms) None [source]
Method for calculating the tracked value from scratch (i.e. without using fast update methods)
- interval_ok(interval: int) bool [source]
Check if the interval specified on attach is ok. Default is that all intervals are OK
- Parameters:
interval – Interval controlling how often a MC observer will be called.
- class clease.montecarlo.observers.MoveObserver(base_atoms: Atoms, only_accept: bool = False)[source]
Store each step from an MC run to reconstruct the individual atoms objects later.
The interval must be set to 1 when attaching this observer, as otherwise steps may be lost and the reconstruction may end up being incorrect.
- Parameters:
base_atoms (ase.Atoms) – The base atoms object which is run in the MC.
only_accept (bool, optional) – Selects whether the only accepted moves or all the attempted moves are saved. If False, every move will be saved. Defaults to False.
- class clease.montecarlo.observers.MultiStateSGCConcObserver(ref_state: SGCState, thermo_states: List[SGCState], calc: Clease)[source]
Observer that tracks the concentration at severl different temperatures and/or chemical potentials. The observer utilizes the following results. Let \(A\) be an observable, \(\beta=\frac{1}{kT}\), \(\mu\) the chemical potential and \(n\) the number of atoms of one of the species in a binary alloy. The average value of the observable is given by
\[\langle A \rangle = \frac{\sum_{conf} A \exp{(\beta\mu n - \beta E)}}{Z(\beta, \mu)}\]where \(Z\) is the partition function. At a different chemical potential \(\mu\) and inverse temperature \(\beta'\),
\[\langle A'\rangle' = \frac{\sum_{conf} A \exp{(\beta'\mu'n-\beta' E)}}{Z(\beta', \mu')}\]After some algabraic manipulation one arrives at
\[\langle A'\rangle' = \frac{\langle A \exp{((\beta'\mu' - \beta\mu)n - (\beta-\beta') E)\rangle}} {\langle \exp{((\beta'\mu' - \beta\mu)n - (\beta-\beta') E)\rangle}}\]where the averages should be taken at inverse temperature \(\beta\) and chemical potential \(\mu\). It should be noted that the predicted value will not be accurate if \(\mu'\) or \(\beta'\) is very different from the reference values \(\mu\) and \(\beta\).
- Parameters:
ref_state – Reference state
thermo_states – List of SGCStates where the concentration should be tracked
calc – Reference to the calculator attached to the atoms object used in the Monte Carlo simulation
- class clease.montecarlo.observers.SGCObserver(calc: Clease, observe_singlets: bool = False)[source]
Observer mainly intended to track additional quantities needed when running SGC Monte Carlo. This observer has to be executed on every MC step.
Parameters:
- calc: clease.calculators.Clease
Clease calculator
- observe_singlets: bool
Whether the singlet values of the calculator are measured during each observation. Measuring singlets is slightly more expensive, so this is disabled by default.
- class clease.montecarlo.observers.SGCState(temp: float, chem_pot: Dict[str, float])[source]
Represent a thermodynamic state in the semi-grand-canonical ensemble.
- Parameters:
temp – Temperature in kelvin
chem_pot – Chemical potentials of the form {c1_0: -0.2, c1_1: 0.3}. The function clease.tools.species_chempot2eci is useful to convert chemical potentials given for each species to chemical potentials for each singlet.
- property prefix: str
Construct a prefix based on the chemical potentials and the temperature
- class clease.montecarlo.observers.SiteOrderParameter(atoms)[source]
Detect phase transitions by monitoring the average number of sites that are occupied by a different element from the initial structure. This observer has to be executed on every MC step.
Parameters:
- atoms: Atoms object
Atoms object use for Monte Carlo
- get_averages()[source]
Get the average and standard deviation of the number of sites that are different from the initial state.
- class clease.montecarlo.observers.Snapshot(atoms: Atoms, fname: str = 'snapshot.traj', mode: str = 'w')[source]
Store a snapshot in a trajectory file.
- Parameters:
atoms – Instance of the atoms objected modofied by the MC object
fname – Name of the trajectory file. Adds extension ‘.traj’ if none is given.
mode – IO mode used by the ASE TrajectoryWriter (must be w or a)