Splined radial integrals¶
Classes for generating splines which can be used as tabulated radial integrals in the various SOAP and LODE calculators.
All classes are based on rascaline.utils.RadialIntegralSplinerBase
. We
provides several ways to compute a radial integral: you may chose and initialize a pre
defined atomic density and radial basis and provide them to
rascaline.utils.SoapSpliner
or rascaline.utils.LodeSpliner
. Both
classes require scipy to be installed in order to perform the numerical integrals.
Alternatively, you can also explicitly provide functions for the radial integral and its
derivative and passing them to rascaline.utils.RadialIntegralFromFunction
.
- class rascaline.utils.RadialIntegralSplinerBase(max_radial: int, max_angular: int, spline_cutoff: float, basis: RadialBasisBase | None, accuracy: float)¶
Bases:
ABC
Base class for splining arbitrary radial integrals.
If
RadialIntegralSplinerBase.radial_integral_derivative()
is not implemented in a child class it will computed based on finite differences.- Parameters:
max_angular – number of radial components
max_radial – number of angular components
spline_cutoff – cutoff radius for the spline interpolation. This is also the maximal value that can be interpolated.
basis – Provide a
RadialBasisBase
instance to orthonormalize the radial integral.accuracy – accuracy of the numerical integration and the splining. Accuracy is reached when either the mean absolute error or the mean relative error gets below the
accuracy
threshold.
- compute(n_spline_points: int | None = None) Dict ¶
Compute the spline for rascaline’s tabulated radial integrals.
- Parameters:
n_spline_points – Use fixed number of spline points instead of find the number based on the provided
accuracy
.- Returns dict:
dictionary for the input as the
radial_basis
parameter of a rascaline calculator.
- abstract radial_integral(n: int, ell: int, positions: ndarray) ndarray ¶
evaluate the radial integral
- class rascaline.utils.SoapSpliner(cutoff: float, max_radial: int, max_angular: int, basis: RadialBasisBase, density: AtomicDensityBase, accuracy: float = 1e-08)¶
Bases:
RadialIntegralSplinerBase
Compute radial integral spline points for real space calculators.
Use only in combination with a real space calculators like
rascaline.SphericalExpansion
orrascaline.SoapPowerSpectrum
. For k-space spherical expansions useLodeSpliner
.If
density
is eitherrascaline.utils.DeltaDensity
orrascaline.utils.GaussianDensity
the radial integral will be partly solved analytical. These simpler expressions result in a faster and more stable evaluation.- Parameters:
cutoff – spherical cutoff for the radial basis
max_radial – number of angular components
max_angular – number of radial components
basis – definition of the radial basis
density – definition of the atomic density
accuracy – accuracy of the numerical integration and the splining. Accuracy is reached when either the mean absolute error or the mean relative error gets below the
accuracy
threshold.
- Raises:
ValueError – if scipy is not available
Example¶
First import the necessary classed and define hyper parameters for the spherical expansions.
>>> from rascaline import SphericalExpansion >>> from rascaline.utils import GaussianDensity, GtoBasis
>>> cutoff = 2 >>> max_radial = 6 >>> max_angular = 4 >>> atomic_gaussian_width = 1.0
Next we initialize our radial basis and the density
>>> basis = GtoBasis(cutoff=cutoff, max_radial=max_radial) >>> density = GaussianDensity(atomic_gaussian_width=atomic_gaussian_width)
And finally the actual spliner instance
>>> spliner = SoapSpliner( ... cutoff=cutoff, ... max_radial=max_radial, ... max_angular=max_angular, ... basis=basis, ... density=density, ... accuracy=1e-3, ... )
Above we reduced
accuracy
from the default value of1e-8
to1e-3
to speed up calculations.As for all spliner classes you can use the output
RadialIntegralSplinerBase.compute()
method directly as theradial_basis
parameter.>>> calculator = SphericalExpansion( ... cutoff=cutoff, ... max_radial=max_radial, ... max_angular=max_angular, ... center_atom_weight=1.0, ... atomic_gaussian_width=atomic_gaussian_width, ... radial_basis=spliner.compute(), ... cutoff_function={"Step": {}}, ... )
You can now use
calculator
to obtain the spherical expansion coefficients of your systems. Note that the the spliner based used here will produce the same coefficients as ifradial_basis={"Gto": {}}
would be used.An additional example using a “rectangular” Laplacian eigenstate (LE) basis is provided in the Laplacian eigenstate basis.
See also
LodeSpliner
for a spliner class that works withrascaline.LodeSphericalExpansion
- class rascaline.utils.LodeSpliner(k_cutoff: float, max_radial: int, max_angular: int, basis: RadialBasisBase, density: AtomicDensityBase, accuracy: float = 1e-08)¶
Bases:
RadialIntegralSplinerBase
Compute radial integral spline points for k-space calculators.
Use only in combination with a k-space/Fourier-space calculators like
rascaline.LodeSphericalExpansion
. For real space spherical expansions useSoapSpliner
.- Parameters:
k_cutoff – spherical reciprocal cutoff
max_radial – number of angular components
max_angular – number of radial components
basis – definition of the radial basis
density – definition of the atomic density
accuracy – accuracy of the numerical integration and the splining. Accuracy is reached when either the mean absolute error or the mean relative error gets below the
accuracy
threshold.
- Raises:
ValueError – if scipy is not available
Example¶
First import the necessary classed and define hyper parameters for the spherical expansions.
>>> from rascaline import LodeSphericalExpansion >>> from rascaline.utils import GaussianDensity, GtoBasis
Note that
cutoff
defined below denotes the maximal distance for the projection of the density. In contrast to SOAP, LODE also takes atoms outside of thiscutoff
into account for the density.>>> cutoff = 2 >>> max_radial = 6 >>> max_angular = 4 >>> atomic_gaussian_width = 1.0
\(1.2 \, \pi \, \sigma\) where \(\sigma\) is the
atomic_gaussian_width
which is a reasonable value for most systems.>>> k_cutoff = 1.2 * np.pi / atomic_gaussian_width
Next we initialize our radial basis and the density
>>> basis = GtoBasis(cutoff=cutoff, max_radial=max_radial) >>> density = GaussianDensity(atomic_gaussian_width=atomic_gaussian_width)
And finally the actual spliner instance
>>> spliner = LodeSpliner( ... k_cutoff=k_cutoff, ... max_radial=max_radial, ... max_angular=max_angular, ... basis=basis, ... density=density, ... )
As for all spliner classes you can use the output
RadialIntegralSplinerBase.compute()
method directly as theradial_basis
parameter.>>> calculator = LodeSphericalExpansion( ... cutoff=cutoff, ... max_radial=max_radial, ... max_angular=max_angular, ... center_atom_weight=1.0, ... atomic_gaussian_width=atomic_gaussian_width, ... potential_exponent=1, ... radial_basis=spliner.compute(), ... )
You can now use
calculator
to obtain the spherical expansion coefficients of your systems. Note that the the spliner based used here will produce the same coefficients as ifradial_basis={"Gto": {}}
would be used.See also
SoapSpliner
for a spliner class that works withrascaline.SphericalExpansion
- class rascaline.utils.RadialIntegralFromFunction(radial_integral: Callable[[int, int, ndarray], ndarray], spline_cutoff: float, max_radial: int, max_angular: int, radial_integral_derivative: Callable[[int, int, ndarray], ndarray] | None = None, center_contribution: ndarray | None = None, accuracy: float = 1e-08)¶
Bases:
RadialIntegralSplinerBase
Compute radial integral spline points based on a provided function.
- Parameters:
radial_integral – Function to compute the radial integral. Function must take
n
,l
, andpositions
as inputs, wheren
andl
are integers andpositions
is a numpy 1-D array that contains the spline points at which the radial integral will be evaluated. The function must return a numpy 1-D array containing the values of the radial integral.spline_cutoff – cutoff radius for the spline interpolation. This is also the maximal value that can be interpolated.
max_radial – number of angular components
max_angular – number of radial components
radial_integral_derivative – The derivative of the radial integral taking the same parameters as
radial_integral
. If it isNone
(default), finite differences are used to calculate the derivative of the radial integral. It is recommended to provide this parameter if possible. Derivatives from finite differences can cause problems when evaluating at the edges of the domain (i.e., at0
andspline_cutoff
) because the function might not be defined outside of the domain.accuracy – accuracy of the numerical integration and the splining. Accuracy is reached when either the mean absolute error or the mean relative error gets below the
accuracy
threshold.center_contribution –
Contribution of the central atom required for LODE calculations. The
center_contribution
is defined as\[c_n = \sqrt{4π}\int_0^\infty dr r^2 R_n(r) g(r)\]where \(g(r)\) is the radially symmetric density function, R_n(r) the radial basis function and \(n\) the current radial channel. This should be pre-computed and provided as a separate parameter.
Example¶
First define a
radial_integral
function>>> def radial_integral(n, ell, r): ... return np.sin(r)
and provide this as input to the spline generator
>>> spliner = RadialIntegralFromFunction( ... radial_integral=radial_integral, ... max_radial=12, ... max_angular=9, ... spline_cutoff=8.0, ... )
Finally, we can use the
spliner
directly in theradial_integral
section of a calculator>>> from rascaline import SoapPowerSpectrum >>> calculator = SoapPowerSpectrum( ... cutoff=8.0, ... max_radial=12, ... max_angular=9, ... center_atom_weight=1.0, ... radial_basis=spliner.compute(), ... atomic_gaussian_width=1.0, # ignored ... cutoff_function={"Step": {}}, ... )
The
atomic_gaussian_width
parameter is required by the calculator but will be will be ignored during the feature computation.A more in depth example using a “rectangular” Laplacian eigenstate (LE) basis is provided in the Splined radial integral how-to guide.