matscipy.calculators.mcfm.calculator

Classes

MultiClusterForceMixingPotential([atoms, ...])

Subclass of ASE.Calculator facilitating a support for multiple QM clusters.

class matscipy.calculators.mcfm.calculator.MultiClusterForceMixingPotential(atoms=None, classical_calculator=None, qm_calculator=None, qm_cluster=None, forced_qm_list=None, change_bonds=True, calculate_errors=False, calculation_always_required=False, buffer_hops=10, verbose=0, enable_check_state=True)

Bases: Calculator

Subclass of ASE.Calculator facilitating a support for multiple QM clusters. It utilizes the given classical_calculator and qm_calculator to initialize an instace of ForceMixingPotential

Extends:

Calculator

Variables:

implemented_properties {list} – [“energy”, “forces”, “potential_energies”, “stress”] all_changes {list} – [‘positions’, ‘numbers’, ‘cell’, ‘pbc’]

Attributes:
directory
label

Methods

attach_hybrid_data([atoms, ...])

Mark atoms that are treated quantm mmechanically for more comprehensive ooutput

band_structure()

Create band-structure object for plotting.

calculate([atoms, properties, system_changes])

Calculate selected properties of the given Atoms object.

calculate_numerical_forces(atoms[, d])

Calculate numerical forces using finite difference.

calculate_numerical_stress(atoms[, d, voigt])

Calculate numerical stress using finite difference.

calculate_properties(atoms, properties)

This method is experimental; currently for internal use.

check_state(atoms[, tol])

Check for any system changes since last calculation.

combine_qm_mm_forces([atoms, forces, ...])

This combines QM and MM forces

compute_stress(atoms, forces)

Compute total stresses using viral theorem.

evaluate_errors([atoms, heavy_only, r_force])

Use the forces and reference forces to get errors on hybrid atom force evaluations

evaluate_qm_cluster_serial([atoms, cluster, ...])

Evaluate forces for a single QM cluster given the buffer hops

get_magnetic_moments([atoms])

Calculate magnetic moments projected onto atoms.

get_property(name[, atoms, allow_calculation])

Get the named property.

get_stresses([atoms])

the calculator should return intensive stresses, i.e., such that stresses.sum(axis=0) == stress

print_message(message[, limit])

Print a message if the calculators verbosity level is above the given threshold

produce_classical_results([atoms])

Call the classical potential ot obtain forces, potential energy and potential energies per atom

produce_qm_clusters(atoms[, potential_energies])

Update qm clusters based on potential energies per atom

read(label)

Read atoms, parameters and calculated properties from output file.

reset()

Clear all information from old calculation.

set(**kwargs)

Set parameters like set(key1=value1, key2=value2, ...).

set_label(label)

Set label and convert label to directory and prefix.

set_qm_atoms(qm_list[, atoms])

Force a certian set of clusters for qmmm evaluation, If forced_qm_list is assigned, the cluster list is not updated throughout the run

calculation_required

export_properties

get_atoms

get_charges

get_default_parameters

get_dipole_moment

get_forces

get_magnetic_moment

get_potential_energies

get_potential_energy

get_stress

read_atoms

todict

implemented_properties: List[str] = ['energy', 'forces', 'potential_energies', 'stress']

Properties calculator can handle (energy, forces, …)

all_changes = ['positions', 'numbers', 'cell', 'pbc']
__init__(atoms=None, classical_calculator=None, qm_calculator=None, qm_cluster=None, forced_qm_list=None, change_bonds=True, calculate_errors=False, calculation_always_required=False, buffer_hops=10, verbose=0, enable_check_state=True)

Initialize a generic ASE potential without any calculating power, This is only to have access to the necessary functions, all the evaluations will be performes in self.mm_pot and self.qm_calculator

Parameters:
  • atoms (ASE.atoms) – atoms object

  • classical_calculator (AE.calculator) – classical calculator

  • qm_calculator (ASE.calculator) – qm calculator

  • qm_cluster (matscipy.calculators.mcfm.qm_cluster) – flagging/cluster carving utility

  • forced_qm_list (list) – add this list to enforce a set of atoms for qm treatment

  • change_bonds (bool) – call the classical potential to update topology

  • calculate_errors (bool) – evaluate errors after each step

  • calculation_always_required (bool) – as name

  • buffer_hops (int) – number of neighbours hop used to construct the core QM region

  • verbose (int) – For now verbose levels are: 0 - nothing is printed 1 - More data is added to Atoms object 10 - Calculate steps are listed 100 - Information about specific QM clusters (the default is 0)

  • enable_check_state (bool) – Save the atoms after each evaluation to enable meth::check_state

calculate(atoms=None, properties=['energy'], system_changes=['positions', 'numbers', 'cell', 'pbc'])

Calculate selected properties of the given Atoms object. Initially, a classical potential is called to evaluate potential energies for each atom and afterwards a qm_cluster object is employed to analyze them. If no atom is flagged for QM treatment, classical forces are returned. In case some atoms are flagged for QM treatment each Qm cluster is independently send to a qmmm potential to evaluate more accurate forces. The results of qmmm evaluations are used to modify the classical forces and the final array is given

produce_classical_results(atoms=None)

Call the classical potential ot obtain forces, potential energy and potential energies per atom

Parameters:

atoms (ASE.atoms) – atoms object

Returns:

  • forces (np.array) – Atomic forces

  • potential_energy (np.array) – Potential energy of the system

  • potential_energies (np.array) – Per atom potential energies

produce_qm_clusters(atoms, potential_energies=None)

Update qm clusters based on potential energies per atom

Parameters:
  • atoms (ASE.atoms) – atoms object

  • potential_energies (np.array) – Per atom potential energies

evaluate_qm_cluster_serial(atoms=None, cluster=None, clusterNumber=0)

Evaluate forces for a single QM cluster given the buffer hops

Parameters:
  • atoms (ASE.atoms) – atoms object

  • cluster (list) – list of core qm atoms

  • clusterNumber (int) – cluster number

Returns:

Cluster

object with forces

qm_atoms mark core qm list

Return type:

cluster_data

combine_qm_mm_forces(atoms=None, forces=None, cluster_data_list=None)

This combines QM and MM forces

Parameters:
  • atoms (ASE.atoms) – atoms object

  • forces (np.array) – atomic forces

  • cluster_data_list (list of matscipy.calculators.mcfm.ClusterData) – information about the clusters

Returns:

forces – atomic forces

Return type:

np.array

attach_hybrid_data(atoms=None, full_qm_atoms_mask=None, cluster_data=None)

Mark atoms that are treated quantm mmechanically for more comprehensive ooutput

Parameters:
  • atoms (ASE.atoms) – atoms object

  • full_qm_atoms_mask (list) – list of all qm atoms

  • cluster_data_list (list of matscipy.calculators.mcfm.ClusterData) – information about the clusters

band_structure()

Create band-structure object for plotting.

calculate_numerical_forces(atoms, d=0.001)

Calculate numerical forces using finite difference.

All atoms will be displaced by +d and -d in all directions.

calculate_numerical_stress(atoms, d=1e-06, voigt=True)

Calculate numerical stress using finite difference.

calculate_properties(atoms, properties)

This method is experimental; currently for internal use.

calculation_required(atoms, properties)
check_state(atoms, tol=1e-15)

Check for any system changes since last calculation.

default_parameters: Dict[str, Any] = {}

Default parameters

property directory: str
discard_results_on_any_change = False

Whether we purge the results following any change in the set() method.

export_properties()
get_atoms()
get_charges(atoms=None)
get_default_parameters()
get_dipole_moment(atoms=None)
get_forces(atoms=None)
get_magnetic_moment(atoms=None)
get_magnetic_moments(atoms=None)

Calculate magnetic moments projected onto atoms.

get_potential_energies(atoms=None)
get_potential_energy(atoms=None, force_consistent=False)
get_property(name, atoms=None, allow_calculation=True)

Get the named property.

get_stress(atoms=None)
get_stresses(atoms=None)

the calculator should return intensive stresses, i.e., such that stresses.sum(axis=0) == stress

ignored_changes: Set[str] = {}

Properties of Atoms which we ignore for the purposes of cache

property label
read(label)

Read atoms, parameters and calculated properties from output file.

Read result from self.label file. Raise ReadError if the file is not there. If the file is corrupted or contains an error message from the calculation, a ReadError should also be raised. In case of succes, these attributes must set:

atoms: Atoms object

The state of the atoms from last calculation.

parameters: Parameters object

The parameter dictionary.

results: dict

Calculated properties like energy and forces.

The FileIOCalculator.read() method will typically read atoms and parameters and get the results dict by calling the read_results() method.

classmethod read_atoms(restart, **kwargs)
reset()

Clear all information from old calculation.

set(**kwargs)

Set parameters like set(key1=value1, key2=value2, …).

A dictionary containing the parameters that have been changed is returned.

Subclasses must implement a set() method that will look at the chaneged parameters and decide if a call to reset() is needed. If the changed parameters are harmless, like a change in verbosity, then there is no need to call reset().

The special keyword ‘parameters’ can be used to read parameters from a file.

set_label(label)

Set label and convert label to directory and prefix.

Examples:

  • label=’abc’: (directory=’.’, prefix=’abc’)

  • label=’dir1/abc’: (directory=’dir1’, prefix=’abc’)

  • label=None: (directory=’.’, prefix=None)

todict(skip_default=True)
evaluate_errors(atoms=None, heavy_only=False, r_force=None)

Use the forces and reference forces to get errors on hybrid atom force evaluations

Parameters:
  • atoms (ASE.atoms) – atoms object

  • heavy_only (bool) – Do not evaluate errors on hydrogens

  • r_force (np.array) – array with reference forces

set_qm_atoms(qm_list, atoms=None)

Force a certian set of clusters for qmmm evaluation, If forced_qm_list is assigned, the cluster list is not updated throughout the run

Parameters:
  • qm_list (list) – list of atoms

  • atoms (ASE.atoms) – atoms object

compute_stress(atoms, forces)

Compute total stresses using viral theorem. WARNING: only works for non-PBC structures

print_message(message, limit=100)

Print a message if the calculators verbosity level is above the given threshold