Splined radial integrals

The classes presented here can take arbitrary radial basis function and density; and compute the radial integral that enters many density-based representations such as SOAP and LODE. This enables using arbitrary, user-defined basis functions and density with the existing calculators. Both classes require scipy to be installed in order to perform the numerical integrals.

class featomic.splines.SoapSpliner(*, cutoff: Cutoff, density: AtomicDensity, basis: ExpansionBasis, n_spline_points: int | None = None)

Compute an explicit spline of the radial integral for SOAP calculators.

This allows a great deal of customization in the radial basis function used (any child class of RadialBasis) and atomic density (any child class of AtomicDensity). This way, users can define custom densities and/or basis, and use them with any of the SOAP calculators. For more information about the radial integral, you can refer to this document.

This class should be used only in combination with SOAP calculators like featomic.SphericalExpansion or featomic.SoapPowerSpectrum. For k-space spherical expansions use LodeSpliner.

If density is either featomic.density.Delta or featomic.density.Gaussian the radial integral will be partly solved analytically, for faster and more stable evaluation.

Example

First let’s define the hyper parameters for the spherical expansions. It is important to note that only class-based hyper parameters are supported (in opposition to dict based hyper parameters).

>>> import featomic
>>> cutoff = featomic.cutoff.Cutoff(radius=2, smoothing=None)
>>> density = featomic.density.Gaussian(width=1.0)
>>> basis = featomic.basis.TensorProduct(
...     max_angular=4,
...     radial=featomic.basis.Gto(max_radial=5, radius=2),
...     # use a reduced ``accuracy`` of ``1e-3`` (the default is ``1e-8``)
...     # to speed up calculations.
...     spline_accuracy=1e-3,
... )

From here we can initialize the spliner instance

>>> spliner = SoapSpliner(cutoff=cutoff, density=density, basis=basis)

You can then give the result of SoapSpliner.get_hypers() directly the calculator:

>>> calculator = featomic.SphericalExpansion(**spliner.get_hypers())

See also

LodeSpliner for a spliner class that works with featomic.LodeSphericalExpansion

param cutoff:

description of the local atomic environment, defined by a spherical cutoff

param density:

atomic density that should be expanded for all neighbors in the local atomic environment

param basis:

basis function to use to expand the neighbors’ atomic density.

param n_spline_points:

number of spline points to use. If None, points will be added to the spline until the accuracy is at least the requested basis.spline_accuracy.

get_hypers()

Get the SOAP hyper-parameters for the splined basis and density.

class featomic.splines.LodeSpliner(density: AtomicDensity, basis: ExpansionBasis, k_cutoff: float | None = None, n_spline_points: int | None = None)

Compute an explicit spline of the radial integral for LODE/k-space calculators.

This allows a great deal of customization in the radial basis function used (any child class of RadialBasis). This way, users can define custom basis, and use them with the LODE calculators. For more information about the radial integral, you can refer to this document.

This class should be used only in combination with k-space/LODE calculators like featomic.LodeSphericalExpansion. For real space spherical expansions you should use SoapSpliner.

Example

First let’s define the hyper parameters for the LODE spherical expansions. It is important to note that only class-based hyper parameters are supported (in opposition to dict based hyper parameters).

>>> import featomic
>>> density = featomic.density.SmearedPowerLaw(smearing=1.0, exponent=1)
>>> basis = featomic.basis.TensorProduct(
...     max_angular=4,
...     radial=featomic.basis.Gto(max_radial=5, radius=2),
... )

From here we can initialize the spliner instance

>>> spliner = LodeSpliner(density=density, basis=basis)

You can then give the result of LodeSpliner.get_hypers() directly the calculator:

>>> calculator = featomic.LodeSphericalExpansion(**spliner.get_hypers())

See also

SoapSpliner for a spliner class that works with featomic.SphericalExpansion

param density:

atomic density that should be expanded for all neighbors in the local atomic environment

param basis:

basis function to use to expand the neighbors’ atomic density. Currently only TensorProduct expansion basis are supported.

param k_cutoff:

Spherical cutoff in reciprocal space.

param n_spline_points:

number of spline points to use. If None, points will be added to the spline until the accuracy is at least the requested basis.spline_accuracy.

get_hypers()

Get the LODE hyper-parameters for the splined basis and density.