Radial Basis

Radial basis functions \(R_{nl}(\boldsymbol{r})\) are besides atomic densities \(\rho_i\) the central ingredients to compute spherical expansion coefficients \(\langle anlm\vert\rho_i\rangle\). Radial basis functions, define how which the atomic density is projected. To be more precise, the actual basis functions are of

\[B_{nlm}(\boldsymbol{r}) = R_{nl}(r)Y_{lm}(\hat{r}) \,,\]

where \(Y_{lm}(\hat{r})\) are the real spherical harmonics evaluated at the point \(\hat{r}\), i.e. at the spherical angles \((\theta, \phi)\) that determine the orientation of the unit vector \(\hat{r} = \boldsymbol{r}/r\).

Radial basis are represented as different child class of rascaline.utils.RadialBasisBase: rascaline.utils.GtoBasis, rascaline.utils.MonomialBasis, and rascaline.utils.SphericalBesselBasis are provided, and you can implement your own by defining a new class.

class rascaline.utils.RadialBasisBase(integration_radius: float)

Bases: ABC

Base class to define radial basis and their evaluation.

The class provides methods to evaluate the radial basis \(R_{nl}(r)\) as well as its (numerical) derivative with respect to positions \(r\).

Parameters:

integration_radius – Value up to which the radial integral should be performed. The usual value is \(\infty\).

abstract compute(n: int, ell: int, integrand_positions: float | ndarray) float | ndarray

Compute the n/l radial basis at all given integrand_positions

Parameters:
  • n – radial channel

  • ell – angular channel

  • integrand_positions – positions to evaluate the radial basis

Returns:

evaluated radial basis

compute_derivative(n: int, ell: int, integrand_positions: ndarray) ndarray

Compute the derivative of the n/l radial basis at all given integrand_positions

This is used for radial integrals with delta-like atomic densities. If not defined in a child class, a numerical derivative based on finite differences of integrand_positions will be used instead.

Parameters:
  • n – radial channel

  • ell – angular channel

  • integrand_positions – positions to evaluate the radial basis

Returns:

evaluated derivative of the radial basis

compute_gram_matrix(max_radial: int, max_angular: int) ndarray

Gram matrix of the current basis.

Parameters:
  • max_radial – number of angular components

  • max_angular – number of radial components

Returns:

orthonormalization matrix of shape (max_angular + 1, max_radial, max_radial)

compute_orthonormalization_matrix(max_radial: int, max_angular: int) ndarray

Compute orthonormalization matrix

Parameters:
  • max_radial – number of angular components

  • max_angular – number of radial components

Returns:

orthonormalization matrix of shape (max_angular + 1, max_radial, max_radial)

class rascaline.utils.GtoBasis(cutoff, max_radial)

Bases: RadialBasisBase

Primitive (not normalized nor orthonormalized) GTO radial basis.

It is defined as

\[R_{nl}(r) = R_n(r) = r^n e^{-\frac{r^2}{2\sigma_n^2}},\]

where \(\sigma_n = \sqrt{n} r_\mathrm{cut}/n_\mathrm{max}\) with \(r_\mathrm{cut}\) being the cutoff and \(n_\mathrm{max}\) the maximal number of radial components.

Parameters:
  • cutoff – spherical cutoff for the radial basis

  • max_radial – number of radial components

compute(n: int, ell: int, integrand_positions: float | ndarray) float | ndarray

Compute the n/l radial basis at all given integrand_positions

Parameters:
  • n – radial channel

  • ell – angular channel

  • integrand_positions – positions to evaluate the radial basis

Returns:

evaluated radial basis

compute_derivative(n: int, ell: int, integrand_positions: float | ndarray) float | ndarray

Compute the derivative of the n/l radial basis at all given integrand_positions

This is used for radial integrals with delta-like atomic densities. If not defined in a child class, a numerical derivative based on finite differences of integrand_positions will be used instead.

Parameters:
  • n – radial channel

  • ell – angular channel

  • integrand_positions – positions to evaluate the radial basis

Returns:

evaluated derivative of the radial basis

class rascaline.utils.MonomialBasis(cutoff)

Bases: RadialBasisBase

Monomial basis.

Basis is consisting of functions

\[R_{nl}(r) = r^{l+2n},\]

where \(n\) runs from \(0,1,...,n_\mathrm{max}-1\). These capture precisely the radial dependence if we compute the Taylor expansion of a generic function defined in 3D space.

Parameters:

cutoff – spherical cutoff for the radial basis

compute(n: int, ell: int, integrand_positions: float | ndarray) float | ndarray

Compute the n/l radial basis at all given integrand_positions

Parameters:
  • n – radial channel

  • ell – angular channel

  • integrand_positions – positions to evaluate the radial basis

Returns:

evaluated radial basis

compute_derivative(n: int, ell: int, integrand_positions: float | ndarray) float | ndarray

Compute the derivative of the n/l radial basis at all given integrand_positions

This is used for radial integrals with delta-like atomic densities. If not defined in a child class, a numerical derivative based on finite differences of integrand_positions will be used instead.

Parameters:
  • n – radial channel

  • ell – angular channel

  • integrand_positions – positions to evaluate the radial basis

Returns:

evaluated derivative of the radial basis

class rascaline.utils.SphericalBesselBasis(cutoff, max_radial, max_angular)

Bases: RadialBasisBase

Spherical Bessel functions used in the Laplacian eigenstate (LE) basis.

Parameters:
  • cutoff – spherical cutoff for the radial basis

  • max_radial – number of angular components

  • max_angular – number of radial components

compute(n: int, ell: int, integrand_positions: float | ndarray) float | ndarray

Compute the n/l radial basis at all given integrand_positions

Parameters:
  • n – radial channel

  • ell – angular channel

  • integrand_positions – positions to evaluate the radial basis

Returns:

evaluated radial basis

compute_derivative(n: int, ell: int, integrand_positions: float | ndarray) float | ndarray

Compute the derivative of the n/l radial basis at all given integrand_positions

This is used for radial integrals with delta-like atomic densities. If not defined in a child class, a numerical derivative based on finite differences of integrand_positions will be used instead.

Parameters:
  • n – radial channel

  • ell – angular channel

  • integrand_positions – positions to evaluate the radial basis

Returns:

evaluated derivative of the radial basis

static compute_zeros(max_angular: int, max_radial: int) ndarray

Zeros of spherical bessel functions.

Code is taken from the Scipy Cookbook.

Parameters:
  • max_radial – number of angular components

  • max_angular – number of radial components

Returns:

computed zeros of the spherical bessel functions