Available Calculators¶
- class featomic.calculators.CalculatorBase(name: str, parameters: str)¶
This is the base class for all calculators in featomic, providing the
CalculatorBase.compute()
function.One can initialize a
Calculator
in two ways: either directly with the registered name and JSON parameter string (which are documented in the Calculator reference); or through one of the child class documented below.- Parameters:
name – name used to register this calculator
parameters – JSON parameter string for the calculator
- property parameters¶
parameters (formatted as JSON) used to create this calculator
- compute(systems: IntoSystem | List[IntoSystem], *, gradients: List[str] | None = None, use_native_system: bool = True, selected_samples: Labels | TensorMap | None = None, selected_properties: Labels | TensorMap | None = None, selected_keys: Labels | None = None) TensorMap ¶
Runs a calculation with this calculator on the given
systems
.- Parameters:
systems – single system or list of systems on which to run the calculation. The systems will automatically be wrapped into compatible classes (using
featomic.systems.wrap_system()
). Multiple types of systems are supported, see the documentation offeatomic.IntoSystem
to get the full list.use_native_system – If
True
(this is the default), copy data from thesystems
into RustSimpleSystem
. This can be a lot faster than having to cross the FFI boundary often when accessing the neighbor list. Otherwise the Python neighbor list is used.gradients –
List of gradients to compute. If this is
None
or an empty list[]
, no gradients are computed. Gradients are stored inside the different blocks, and can be accessed withdescriptor.block(...).gradient(<parameter>)
, where<parameter>
is a string describing the gradients. The following gradients are available:"positions"
, for gradients of the representation with respect to atomic positions, with fixed cell matrix parameters. Positions gradients are computed as\[\frac{\partial \langle q \vert A_i \rangle} {\partial \mathbf{r_j}}\]where \(\langle q \vert A_i \rangle\) is the representation around atom \(i\) and \(\mathbf{r_j}\) is the position vector of the atom \(j\).
Note: Position gradients of an atom are computed with respect to all other atoms within the representation. To recover the force one has to accumulate all pairs associated with atom \(i\).
"strain"
, for gradients of the representation with respect to strain. These gradients are typically used to compute the virial, and from there the pressure acting on a system. To compute them, we pretend that all the positions \(\mathbf r\) and unit cell \(\mathbf H\) have been scaled by a strain matrix \(\epsilon\):\[\begin{split}\mathbf r &\rightarrow \mathbf r \left(\mathbb{1} + \epsilon \right)\\ \mathbf H &\rightarrow \mathbf H \left(\mathbb{1} + \epsilon \right)\end{split}\]and then take the gradients of the representation with respect to this matrix:
\[\frac{\partial \langle q \vert A_i \rangle} {\partial \mathbf{\epsilon}}\]"cell"
, for gradients of the representation with respect to the system’s cell parameters. These gradients are computed at fixed positions, and often not what you want when computing gradients explicitly (they are mainly used infeatomic.torch
to integrate with backward propagation). If you are trying to compute the virial or the stress, you should use"strain"
gradients instead.\[\left. \frac{\partial \langle q \vert A_i \rangle} {\partial \mathbf{H}} \right |_\mathbf{r}\]
selected_samples –
Set of samples on which to run the calculation. Use
None
to run the calculation on all samples in thesystems
(this is the default).If
selected_samples
is anmetatensor.TensorMap
, then the samples for each key will be used as-is when computing the representation.If
selected_samples
is a set ofmetatensor.Labels
containing the same variables as the default set of samples for this calculator, then only entries from the full set that also appear inselected_samples
will be used.If
selected_samples
is a set ofmetatensor.Labels
containing a subset of the variables of the default set of samples, then only samples from the default set with the same values for these variables as one of the entries inselected_samples
will be used.selected_properties –
Set of properties to compute. Use
None
to run the calculation on all properties (this is the default).If
selected_properties
is anmetatensor.TensorMap
, then the properties for each key will be used as-is when computing the representation.If
selected_properties
is a set ofmetatensor.Labels
containing the same variables as the default set of properties for this calculator, then only entries from the full set that also appear inselected_properties
will be used.If
selected_properties
is a set ofmetatensor.Labels
containing a subset of the variables of the default set of properties, then only properties from the default set with the same values for these variables as one of the entries inselected_properties
will be used.selected_keys – Selection for the keys to include in the output. If this is
None
, the default set of keys (as determined by the calculator) will be used. Note that this default set of keys can depend on which systems we are running the calculation on.
- class featomic.AtomicComposition(*, per_system)¶
Bases:
CalculatorBase
An atomic composition calculator for obtaining the stoichiometric information.
For
per_system=False
calculator has one propertycount
that is1
for all atoms, and has a sample index that indicates the central atom type.For
per_system=True
a sum for each system is performed. The number of atoms per system is saved. The only sample left is namedsystem
.
- class featomic.NeighborList(*, cutoff, full_neighbor_list, self_pairs=False)¶
Bases:
CalculatorBase
This calculator computes the neighbor list for a given spherical cutoff, and returns the list of distance vectors between all pairs of atoms strictly inside the cutoff.
Users can request either a “full” neighbor list (including an entry for both
i - j
pairs andj - i
pairs) or save memory/computational by only working with “half” neighbor list (only including one entry for eachi/j
pair)Pairs between an atom and it’s own periodic copy can appear when the cutoff is larger than the cell under periodic boundary conditions. Self pairs with a distance of 0 (i.e. self pairs inside the original unit cell) are only included when using
self_pairs=True
.This calculator produces a single property (
"distance"
) with three components ("pair_xyz"
) containing the x, y, and z component of the distance vector of the pair.The samples contain the two atoms indexes, as well as the number of cell boundaries crossed to create this pair.
- class featomic.SortedDistances(*, cutoff, max_neighbors, separate_neighbor_types)¶
Bases:
CalculatorBase
Sorted distances vector representation of an atomic environment.
Each atomic center is represented by a vector of distance to its neighbors within the spherical
cutoff
, sorted from smallest to largest. If there are less neighbors thanmax_neighbors
, the remaining entries are filled withcutoff
instead.Separate atomic types for neighbors are represented separately, meaning that the
max_neighbors
parameter only apply to a single atomic type.For a full description of the hyper-parameters, see the corresponding documentation.
- class featomic.SphericalExpansion(*, cutoff=None, density=None, basis=None, **kwargs)¶
Bases:
CalculatorBase
Spherical expansion of Smooth Overlap of Atomic Positions (SOAP).
The spherical expansion is at the core of representations in the SOAP family of descriptors. The spherical expansion represent atomic density as a collection of Gaussian functions centered on each atom, and then represent the local density around each atom (up to a cutoff) on a basis of radial functions and spherical harmonics. This representation is not rotationally invariant, for that you should use the
SoapPowerSpectrum
class.See this review article for more information on the SOAP representation, and this paper for information on how it is implemented in featomic.
For a full description of the hyper-parameters, see the corresponding documentation.
- class featomic.SphericalExpansionByPair(*, cutoff=None, density=None, basis=None, **kwargs)¶
Bases:
CalculatorBase
Pair-by-pair version of the spherical expansion of Smooth Overlap of Atomic Positions (SOAP).
This is usually an intermediary step when computing the full
SphericalExpansion
, which most users do not need to care about. When computing SOAP, the density around the central atom \(i\) is a sum of pair contributions:\[\rho_i(\mathbf{r}) = \sum_j g(\mathbf{r_{ji}} - \mathbf{r}).\]The full spherical expansion is then computed as a sum over all pairs within the cutoff:
\[ \begin{align}\begin{aligned}\int Y^l_m(\mathbf{r})\ R_n(\mathbf{r}) \, \rho_i(\mathbf{r}) \mathrm{d}\mathbf{r} = \sum_j C^{ij}_{nlm}\\C^{ij}_{nlm} = \int Y^l_m(\mathbf{r}) \ R_n(\mathbf{r}) \, g(\mathbf{r_{ji}} - \mathbf{r}) \, \mathrm{d}\mathbf{r}\end{aligned}\end{align} \]This calculators computes the individual \(C^{ij}_{nlm}\) terms, which can then be used to build multi-center representations such as the ones discussed in “Unified theory of atom-centered representations and message-passing machine-learning schemes”.
For a full description of the hyper-parameters, see the corresponding documentation.
- class featomic.SoapRadialSpectrum(*, cutoff=None, density=None, basis=None, **kwargs)¶
Bases:
CalculatorBase
Radial spectrum of Smooth Overlap of Atomic Positions (SOAP).
The SOAP radial spectrum represent each atom by the radial average of the density of its neighbors. It is very similar to a radial distribution function g(r). It is a 2-body representation, only containing information about the distances between atoms.
See this review article for more information on the SOAP representation, and this paper for information on how it is implemented in featomic.
For a full description of the hyper-parameters, see the corresponding documentation.
- class featomic.SoapPowerSpectrum(*, cutoff=None, density=None, basis=None, **kwargs)¶
Bases:
CalculatorBase
Power spectrum of Smooth Overlap of Atomic Positions (SOAP).
The SOAP power spectrum is the main member of the SOAP family of descriptors. The power spectrum is based on the
SphericalExpansion
coefficients, which are combined to create a rotationally invariant three-body descriptor.See this review article for more information on the SOAP representation, and this paper for information on how it is implemented in featomic.
For a full description of the hyper-parameters, see the corresponding documentation.
See also
featomic.utils.PowerSpectrum
is an implementation that allows to compute the power spectrum from different spherical expansions.
- class featomic.LodeSphericalExpansion(*, density=None, basis=None, k_cutoff=None, **kwargs)¶
Bases:
CalculatorBase
Long-Distance Equivariant (LODE).
The spherical expansion is at the core of representations in the LODE family. The LODE spherical expansion represent atomic density as a collection of ‘decorated’ gaussian functions centered on each atom, and then represent the local density around each atom on a basis of radial functions and spherical harmonics. This representation is not rotationally invariant.
See this article for more information on the LODE representation.
For a full description of the hyper-parameters, see the corresponding documentation.