guess_core_dist

nanoCAT.bde.guess_core_dist

A module for estimating ideal values for ["optional"]["qd"]["bde"]["core_core_dist"] .

Index

guess_core_core_dist(mol[, atom, dr, r_max, ...])

Guess a value for the ["optional"]["qd"]["bde"]["core_core_dist"] parameter in CAT.

API

nanoCAT.bde.guess_core_dist.guess_core_core_dist(mol, atom=None, dr=0.1, r_max=8.0, window_length=21, polyorder=7)[source]

Guess a value for the ["optional"]["qd"]["bde"]["core_core_dist"] parameter in CAT.

The estimation procedure involves finding the first minimum in the radial distribution function (RDF) of mol. After smoothing the RDF wth a Savitzky-Golay filer, the gradient of the RDF is explored (starting from the RDFs’ global maximum) until a stationary point is found with a positive second derivative (i.e. a minimum).

Examples

>>> from scm.plams import Molecule
>>> from nanoCAT.bde.guess_core_dist import guess_core_core_dist

>>> atom1 = 'Cl'  # equivalent to ('Cl', 'Cl')
>>> atom2 = 'Cl', 'Br'

>>> mol = Molecule(...)

>>> guess_core_core_dist(mol, atom1)
>>> guess_core_core_dist(mol, atom2)
Parameters:
  • mol (array-like [float], shape \((n, 3)\)) – A molecule.

  • atom (str or int, optional) – An atomic number or symbol for defining an atom subset within mol. The RDF is constructed for this subset. Providing a 2-tuple will construct the RDF between these 2 atom subsets.

  • dr (float) – The RDF integration step-size in Angstrom, i.e. the distance between concentric spheres.

  • r_max (float) – The maximum to be evaluated interatomic distance in the RDF.

  • window_length (int) – The length of the filter window (i.e. the number of coefficients) for the Savitzky-Golay filter.

  • polyorder (int) – The order of the polynomial used to fit the samples for the Savitzky-Golay filter.

Returns:

The interatomic radius of the first RDF minimum (following the first maximum).

Return type:

float

Raises:
  • MoleculeError – Raised if atom is not in mol.

  • ValueError – Raised if no minimum is found in the smoothed RDF.

See also

savgol_filter()

Apply a Savitzky-Golay filter to an array.