Index

MonteCarloABC(molecule, package_manager, param)

The base MonteCarloABC class.

ARMC(phi[, iter_len, sub_iter_len])

The Addaptive Rate Monte Carlo class (ARMC).

ARMCPT([swapper])

An ARMC subclass implementing a parallel tempering procedure.

API

class FOX.armc.MonteCarloABC(molecule, package_manager, param, keep_files=False, hdf5_file='armc.hdf5', logger=None, pes_post_process=None, **kwargs)[source]

The base MonteCarloABC class.

property molecule

Get value or set value as a tuple of MultiMolecule instances.

property pes_post_process

Get or set post-processing functions.

property logger

Get or set the logger.

keys()[source]

Return a set-like object providing a view of this instance’s keys.

items()[source]

Return a set-like object providing a view of this instance’s key/value pairs.

values()[source]

Return an object providing a view of this instance’s values.

get(key, default=None)[source]

Return the value for key if it’s available; return default otherwise.

add_pes_evaluator(name, func, err_func, args=(), kwargs=mappingproxy({}), validation=False, ref=None)[source]

Add a callable to this instance for constructing PES-descriptors.

Examples

>>> from FOX import MonteCarlo, MultiMolecule

>>> mc = MonteCarlo(...)
>>> mol = MultiMolecule.from_xyz(...)

# Prepare arguments
>>> name = 'rdf'
>>> func = FOX.MultiMolecule.init_rdf
>>> atom_subset = ['Cd', 'Se', 'O']  # Keyword argument for func

# Add the PES-descriptor constructor
>>> mc.add_pes_evaluator(name, func, kwargs={'atom_subset': atom_subset})
Parameters
  • name (str) – The name under which the PES-descriptor will be stored (e.g. "RDF").

  • func (Callable) – The callable for constructing the PES-descriptor. The callable should take an array-like object as input and return a new array-like object as output.

  • err_func (Callable) – The function for computing the auxilary error.

  • args (Sequence) – A sequence of positional arguments.

  • kwargs (dict or Iterable[dict]) – A dictionary or an iterable of dictionaries with keyword arguments. Providing an iterable allows one to use a unique set of keyword arguments for each molecule in MonteCarlo.molecule.

  • validation (bool) – Whether the PES-descriptor is used exclusively for validation or not.

property clear_jobs

Delete all cp2k output files.

run_jobs()[source]

Run a geometry optimization followed by a molecular dynamics (MD) job.

Returns a new MultiMolecule instance constructed from the MD trajectory and the path to the MD results. If no trajectory is available (i.e. the job crashed) return None instead.

  • The MD job is constructed according to the provided settings in self.job.

Returns

A list of MultiMolecule instance(s) constructed from the MD trajectory. Will return None if one of the jobs crashed

Return type

list[FOX.MultiMolecule], optional

move(idx=None)[source]

Update a random parameter in self.param by a random value from self.move.range.

Performs in inplace update of the 'param' column in self.param. By default the move is applied in a multiplicative manner. self.job.md_settings and self.job.preopt_settings are updated to reflect the change in parameters.

Examples

>>> print(armc.param['param'])
charge   Br      -0.731687
         Cs       0.731687
epsilon  Br Br    1.045000
         Cs Br    0.437800
         Cs Cs    0.300000
sigma    Br Br    0.421190
         Cs Br    0.369909
         Cs Cs    0.592590
Name: param, dtype: float64

>>> for _ in range(1000):  # Perform 1000 random moves
>>>     armc.move()

>>> print(armc.param['param'])
charge   Br      -0.597709
         Cs       0.444592
epsilon  Br Br    0.653053
         Cs Br    1.088848
         Cs Cs    1.025769
sigma    Br Br    0.339293
         Cs Br    0.136361
         Cs Cs    0.101097
Name: param, dtype: float64
Parameters

idx (int, optional) – The column key for param_mapping["param"].

Returns

A tuple with the (new) values in the 'param' column of self.param.

Return type

tuple[float, ...]

get_pes_descriptors(get_first_key=False)[source]

Check if a key is already present in history_dict.

If True, return the matching list of PES descriptors; If False, construct and return a new list of PES descriptors.

  • The PES descriptors are constructed by the provided settings in self.pes.

Parameters

get_first_key (bool) – Keep both the files and the job_cache if this is the first ARMC iteration. Usefull for manual inspection in case cp2k hard-crashes at this point.

Returns

A previous value from history_dict or a new value from an MD calculation & a MultiMolecule instance constructed from the MD simulation. Values are set to np.inf if the MD job crashed.

Return type

dict[str, np.ndarray[np.float64]], dict[str, np.ndarray[np.float64]] and list[FOX.MultiMolecule]

class FOX.armc.ARMC(phi, iter_len=50000, sub_iter_len=100, **kwargs)[source]

The Addaptive Rate Monte Carlo class (ARMC).

A subclass of MonteCarloABC.

iter_len

The total number of ARMC iterations \(\kappa \omega\).

Type

int

super_iter_len

The length of each ARMC subiteration \(\kappa\).

Type

int

sub_iter_len

The length of each ARMC subiteration \(\omega\).

Type

int

phi

A PhiUpdater instance.

Type

PhiUpdaterABC

\**kwargs

Keyword arguments for the MonteCarlo superclass.

Type

Any

acceptance()[source]

Create an empty 1D boolean array for holding the acceptance.

to_yaml_dict(*, path='.', folder='MM_MD_workdir', logfile='armc.log', psf=None)[source]

Convert an ARMC instance into a .yaml readable by ARMC.from_yaml.

Returns

A dictionary.

Return type

dict[str, Any]

do_inner(kappa, omega, acceptance, key_old)[source]

Run the inner loop of the ARMC.__call__() method.

Parameters
  • kappa (int) – The super-iteration, \(\kappa\), in ARMC.__call__().

  • omega (int) – The sub-iteration, \(\omega\), in ARMC.__call__().

  • acceptance (np.ndarray[np.bool_]) – An array with the acceptance over the course of the latest super-iteration

  • key_new (tuple[float, ...]) – A tuple with the latest set of forcefield parameters.

Returns

The latest set of parameters.

Return type

tuple[float, ...]

property apply_phi

Apply phi to value.

to_hdf5(mol_list, accept, aux_new, aux_validation, pes_new, pes_validation, kappa, omega)[source]

Construct a dictionary with the hdf5_kwarg and pass it to to_hdf5().

Parameters
Returns

A dictionary with the hdf5_kwarg argument for to_hdf5().

Return type

dict[str, Any]

get_aux_error(pes_dict, validation=False)[source]

Return the auxiliary error \(\Delta \varepsilon_{QM-MM}\).

The auxiliary error is constructed using the PES descriptors in values with respect to self.ref.

The default function is equivalent to:

\[\Delta \varepsilon_{QM-MM} = \frac{ \sum_{i}^{N} |r_{i}^{QM} - r_{i}^{MM}|^2 } {r_{i}^{QM}}\]
Parameters

pes_dict (dict[str, np.ndarray[np.float64]]) – An dictionary with \(m*n\) PES descriptors each.

Returns

An array with \(m*n\) auxilary errors

Return type

np.ndarray[np.float64], shape \((m, n)\)

restart()[source]

Restart a previously started Addaptive Rate Monte Carlo procedure.

class FOX.armc.ARMCPT(swapper=<function swap_random>, **kwargs)[source]

An ARMC subclass implementing a parallel tempering procedure.

acceptance()[source]

Create an empty 2D boolean array for holding the acceptance.

do_inner(kappa, omega, acceptance, key_old)[source]

Run the inner loop of the ARMC.__call__() method.

Parameters
  • kappa (int) – The super-iteration, \(\kappa\), in ARMC.__call__().

  • omega (int) – The sub-iteration, \(\omega\), in ARMC.__call__().

  • acceptance (np.ndarray[np.bool_]) – An array with the acceptance over the course of the latest super-iteration

  • key_new (tuple[float, ...]) – A tuple with the latest set of forcefield parameters.

Returns

The latest set of parameters.

Return type

tuple[float, ...]

to_yaml_dict(*, path='.', folder='MM_MD_workdir', logfile='armc.log', psf=None)[source]

Convert an ARMC instance into a .yaml readable by ARMC.from_yaml.

Returns

A dictionary.

Return type

dict[str, Any]