matscipy.calculators.mcfm.qm_cluster_tools.qm_flagging_tool

Functions

exponential_moving_average(oldset[, newset, ...])

Apply the exponential moving average to the given array

update_avg_property_per_atom(atoms, ...)

Update the per atom property using running avarages and store it in atoms.properties[property_str]

Classes

QMFlaggingTool([mediator, ...])

This class is responsible for flagging atoms that move out of their equilibrium

class matscipy.calculators.mcfm.qm_cluster_tools.qm_flagging_tool.QMFlaggingTool(mediator=None, qm_flag_potential_energies=None, small_cluster_hops=3, only_heavy=False, ema_parameter=0.1, energy_cap=None, energy_increase=1)

Bases: BaseQMClusterTool

This class is responsible for flagging atoms that move out of their equilibrium

Methods

create_cluster_around_atom(atoms, atom_id[, ...])

Carve a cluster around the atom with atom_id This function operates on sets and returns a set

expand_cluster(special_atoms_list)

Include extra atoms in the cluster.

find_neighbours(atoms, index)

Find the neighbours of atom i using self.neighbour_list returns a list of [heavy_neighbours, hydrogen_neighbours]

get_energized_list(atoms, data_array, ...)

Produce a list of atoms that are ot be flagged as a QM region based on the properties given in the array according to the tolerance given.

hydrogenate_cluster(atoms, cluster)

Add neigoburing hydrogens to a cluster composed of heavy ions The input should be a set representing heavy ions in a cluster This functions operates on sets

join_clusters([verbose])

This function will join the clusters if they overlap Input is an array of sets each representing individual small cluster

update_qm_region(atoms[, potential_energies])

Update the QM region while the simulation is running

__init__(mediator=None, qm_flag_potential_energies=None, small_cluster_hops=3, only_heavy=False, ema_parameter=0.1, energy_cap=None, energy_increase=1)

This class is responsible for flagging atoms that move out of their equilibrium

Parameters:
  • mediator (matscipy.calculators.mcfm.QMCluster) – class responsible for managing the QM clusters in the simulation

  • qm_flag_potential_energies (np.array) –

    threshholds for flagging indivual atoms. The diensions are (nAtoms, 2) where:

    column 1: threshold to enter the QM regios column 2: threshold to stay the QM region

  • small_cluster_hops (int) – Each flagged atom and atoms around it within small_cluster_hops neighbour hops will generate a single cluster, clusters are later joined.

  • only_heavy (bool) – If True, only consider non-hydrogen atoms in cluster expansion. Hydrogens are added later

  • ema_parameter (float) – parameter lambda in the exponential mean average calculation

  • energy_cap (float) – if not None, cap potential energy per atom at this value

  • energy_increase (int) – Multiplier for potential energy per atom, used to scale it for convininece

get_energized_list(atoms, data_array, property_str, hysteretic_tolerance)

Produce a list of atoms that are ot be flagged as a QM region based on the properties given in the array according to the tolerance given.

Parameters:
  • atoms (ase.Atoms) – Whole structure

  • data_array (array) – an array of per atom data providing information

  • property_str (str) – name of th property so that it can be stored in atoms.properties.

  • hysteretic_tolerance (array) –

    Threshholds for flagging indivual atoms. The diensions are (nAtoms, 2) where:

    column 1: threshold to enter the QM regios column 2: threshold to stay the QM region

Returns:

List of flagged atoms

Return type:

list

create_cluster_around_atom(atoms, atom_id, hydrogenate=False)

Carve a cluster around the atom with atom_id This function operates on sets and returns a set

Parameters:
  • atoms (ase.Atoms) – Whole structure

  • atom_id (int) – Atomic index

  • hydrogenate (bool) – If true, hydrogenate the resulting structure

Returns:

atoms in the new cluster

Return type:

list

join_clusters(verbose=False)

This function will join the clusters if they overlap Input is an array of sets each representing individual small cluster

Parameters:

verbose (bool) – Print messages during calculation

expand_cluster(special_atoms_list)

Include extra atoms in the cluster.

If one of the special atoms is included in one of the clusters, add all other special atoms to this cluster

Parameters:

special_atoms_list (list) – list of the special atoms

update_qm_region(atoms, potential_energies=None)

Update the QM region while the simulation is running

Parameters:
  • atoms (ase.Atoms) – whole structure

  • potential_energies (array) – Potential energy per atom

Returns:

list of individual clusters as lists of atoms

Return type:

list of lists of ints

find_neighbours(atoms, index)

Find the neighbours of atom i using self.neighbour_list returns a list of [heavy_neighbours, hydrogen_neighbours]

Parameters:
  • atoms (ase.Atoms object) – structure in which it is necessary to find the neighbours

  • index (int) – atomic index

Returns:

  • list – non-hydrogen neighbours

  • list – hydrogen neighbours

hydrogenate_cluster(atoms, cluster)

Add neigoburing hydrogens to a cluster composed of heavy ions The input should be a set representing heavy ions in a cluster This functions operates on sets

Parameters:
  • atoms (ase.Atoms object) – structure in which it is necessary to find the neighbours

  • cluster (ase.Atoms object) – sub-structure of the larger struct that needs its dangling bonds hydrogenated

Returns:

The original cluster but now hydrogenated

Return type:

ase.Atoms

matscipy.calculators.mcfm.qm_cluster_tools.qm_flagging_tool.exponential_moving_average(oldset, newset=None, ema_parameter=0.1)

Apply the exponential moving average to the given array

Parameters:
  • oldset (array) – old values

  • newset (array) – new data set

  • ema_parameter (float) – parameter lambda

matscipy.calculators.mcfm.qm_cluster_tools.qm_flagging_tool.update_avg_property_per_atom(atoms, data_array, property_str, ema_parameter)

Update the per atom property using running avarages and store it in atoms.properties[property_str]

Parameters:
  • atoms (ase.Atoms) – structure that need updated values

  • data_array (array) – data that need to be attached to atoms

  • property_str (str) – key for structure properties dictionary

  • ema_parameter (float) – Coefficient for the Exponential Moving Average