API Reference

Representations

class anisoap.representations.ellipsoidal_density_projection.EllipsoidalDensityProjection(max_angular, radial_basis_name, cutoff_radius, compute_gradients=False, subtract_center_contribution=False, radial_gaussian_width=None, max_radial=None, rotation_key='quaternion', rotation_type='quaternion', basis_rcond=0, basis_tol=1e-08)

Bases: object

Class for computing spherical projection coefficients.

Compute the spherical projection coefficients for a system of ellipsoids assuming a multivariate Gaussian density.

Initialize the calculator using the hyperparameters.

Parameters:
  • max_angular (int) – Number of angular functions

  • radial_basis (_RadialBasis) – The radial basis. Currently implemented are ‘gto’, ‘monomial’.

  • compute_gradients (bool) – Compute gradients

  • subtract_center_contribution (bool) – Subtract contribution from the central atom.

  • rotation_key (string) – Key under which rotations are stored in ase frames arrays

  • rotation_type (string) – Type of rotation object being passed. Currently implemented are ‘quaternion’ and ‘matrix’

  • max_radial (None, int, list of int) – Number of radial bases to use. Can either correspond to number of bases per spherical harmonic or a value to use with every harmonic. If None, then for every l, (max_angular - l) // 2 + 1 will be used.

features
Type:

numpy.ndarray

feature_gradients
Type:

numpy.ndarray

power_spectrum(frames, mean_over_samples=True, show_progress=False, rust_moments=True)

Helper function to compute the power spectrum of AniSOAP

Computes the power spectrum of AniSOAP with the inputs of AniSOAP hyperparameters and ellipsoidal frames using ellipsoidal density projection. It checks if each ellipsoidal frame contains all required attributes and processes Clebsch-Gordan coefficients for the angular component of the AniSOAP descriptors.

Parameters:
  • frames (list) – A list of ellipsoidal frames, where each frame contains attributes: ‘c_diameter[1]’, ‘c_diameter[2]’, ‘c_diameter[3]’, ‘c_q’, ‘positions’, and ‘numbers’. It only accepts c_q for the angular attribute of each frame.

  • mean_over_samples (bool) – A function that returns the sum of coefficients of the frames in the sample.

Returns:

  • x_asoap_raw when kwarg mean_over_samples=True or mvg_nu2 when mean_over_samples=False

  • x_asoap_raw (A 2-dimensional np.array with shape (n_samples, n_features). This AniSOAP power spectrum aggregates (averages) over each sample.)

  • mvg_nu2 (a TensorMap of unaggregated power spectrum features.)

transform(frames, show_progress=False, normalize=True, rust_moments=True, return_pef=False)

Computes features and gradients for frames

Computes the features and (if compute_gradients == True) gradients for all the provided frames. The features and gradients are stored in features and feature_gradients attribute.

Parameters:
  • frames (ase.Atoms) – List containing all ase.Atoms types

  • show_progress (bool) – Show progress bar for frame analysis and feature generation

  • normalize (bool) – Whether to perform Lowdin Symmetric Orthonormalization or not.

  • rust_moments (bool) – Use the ported rust code, which should result in increased speed. Default = True. In the future, once we ensure integrity checks with the original python code, this kwarg will be deprecated, and the rust version will always be used.

Returns:

  • None, but stores the projection coefficients and (if desired)

  • gradients as arrays as features and features_gradients.

anisoap.representations.ellipsoidal_density_projection.contract_pairwise_feat(pair_ellip_feat, types, show_progress=False)

Function to sum over the pairwise expansion

\[\sum_{j \in a} \langle anlm|\rho_{ij} \rangle = \langle anlm|\rho_i \rangle\]
Parameters:
  • pair_ellip_feat (metatensor.TensorMap) – TensorMap returned from “pairwise_ellip_expansion()” with keys (types_1, types_2,l) enumerating the possible species pairs and the angular channels.

  • types (list of ints) – List of atomic numbers present across the data frames

  • show_progress (bool) – Show progress bar for frame analysis and feature generation

Returns:

A metatensor TensorMap with keys (types, l) “types” takes the value of the atomic numbers present in the dataset and “l” is the angular channel. Each block of this tensormap has as samples (“type”, “center”), yielding the indices of the frames and atoms that correspond to “species” are present. block.value is a 3D array of the form (num_samples, num_components, properties) where num_components take on the same values as in the pair_ellip_feat_feat.block. block.properties now has an additional index for neighbor_species that corresponds to “a” in \(\langle anlm|rho_i \rangle\)

Return type:

TensorMap

anisoap.representations.ellipsoidal_density_projection.pairwise_ellip_expansion(lmax, neighbor_list, types, frame_to_global_atom_idx, rotation_matrices, ellipsoid_lengths, sph_to_cart, radial_basis, show_progress=False, rust_moments=True, normalize=True)

Computes pairwise expansion

Function to compute the pairwise expansion \(\langle anlm|\rho_{ij} \rangle\) by combining the moments and the spherical to Cartesian transformation.

Parameters:
  • lmax (int) – Maximum angular

  • neighbor_list (metatensor.TensorMap) – Full neighborlist with keys (types_1, types_2) enumerating the possible species pairs. Each block contains as samples, the atom indices of (first_atom, second_atom) that correspond to the key, and block.value is a 3D array of the form (num_samples, num_components,properties), with num_components=3 corresponding to the (x,y,z) components of the vector from first_atom to second_atom. Depending on the cutoff some species pairs may not appear. Self pairs are not present but in PBC, pairs between copies of the same atom are accounted for.

  • types (list of ints) – List of atomic numbers present across the data frames

  • frame_to_global_atom_idx (list of ints) – The length of the list equals the number of frames, each entry enumerating the number of atoms in corresponding frame

  • rotation_matrices (np.array of dimension ((num_atoms,3,3))) –

  • ellipsoid_lengths (np.array of dimension ((num_atoms,3))) –

  • radial_basis (Instance of the RadialBasis Class) – anisoap.representations.radial_basis.RadialBasis that has been instantiated appropriately with the cutoff radius, radial basis type.

  • show_progress (bool) – Show progress bar for frame analysis and feature generation

  • rust_moments (bool) – Use the ported rust code, which should result in increased speed. Default = True. In the future, once we ensure integrity checks with the original python code, this kwarg will be deprecated, and the rust version will always be used.

Returns:

A metatensor TensorMap with keys (types_1, types_2, l) where (“types_1”, “types_2”) is key in the neighbor_list and “l” is the angular channel. Each block of this tensormap has the same samples as the corresponding block of the neighbor_list. block.value is a 3D array of the form (num_samples, num_components, properties) where num_components form the 2*l+1 values for the corresponding angular channel and the properties dimension corresponds to the radial channel.

Return type:

TensorMap

class anisoap.representations.radial_basis.GTORadialBasis(max_angular, cutoff_radius, *, radial_gaussian_width, max_radial=None, rcond=1e-08, tol=0.001)

Bases: _RadialBasis

A subclass of _RadialBasis that contains attributes and methods required for the GTO basis.

The GTO basis of order n is defined to be \(R(r) = r^n e^{(-r^2/(2\sigma^2))}\).

For GTO basis with defined nmax and lmax, our radial basis set consists of GTO of order n=0 to n=lmax + 2nmax.

For GTO basis with coupled nmax and lmax, our radial basis set consists of GTO of order n=0 to n=max(lmax + 2nmax)

calc_overlap_matrix()

Computes the overlap matrix for GTOs.

The overlap matrix is a Gram matrix whose entries are the overlap:

\[S_{ij} = \int_0^\infty dr r^2 \phi_i \phi_j\]

The overlap has an analytic solution (see above functions).

The overlap matrix is the first step to generating an orthonormal basis set of functions (Lodwin Symmetric Orthonormalization). The actual orthonormalization matrix cannot be fully precomputed because each tensor block uses a different set of GTOs. Hence, we precompute the full overlap matrix of dim l_max, and while orthonormalizing each tensor block, we generate the respective orthonormal matrices from slices of the full overlap matrix.

Returns:

The overlap matrix

Return type:

2D array

get_basis(rs)

Evaluate orthonormalized GTO basis functions on mesh rs.

If lmax and nmax defined, then the number of functions outputted is lmax*(nmax+1)

If lmax and nmax coupled, then the number of functions outputted is \(\sum_{l=0}^{\text{lmax}} (\text{number_of_radial_functions_per_l})\)

Parameters:

rs (np.array) – a mesh to evaluate the basis functions

Returns:

a matrix containing orthonormalized GTO basis functions evaluated on rs

Return type:

np.array

get_num_radial_functions()

Output the number of radial basis functions considered per value of l. If max_angular and max_radial are both specified, then the list will contain repeated values of max_radial Otherwise, the outputted list will specify the number of radial basis functions per l, which may be automatically calculated if max_radial=None. If a custom list of max_radial is specified when initializing, then it will return the same inputted list.

Returns:

The attribute self.num_radial_functions, which is a list containing the number of radial basis functions considered per l.

Return type:

array_like

orthonormalize_basis(features: TensorMap)

Applies in-place orthonormalization on the features.

Apply an in-place orthonormalization on the features, using Lodwin Symmetric Orthonormalization. Each block in the features TensorMap uses a GTO set of l + 2n, so we must take the appropriate slices of the overlap matrix to compute the orthonormalization matrix. An instructive example of Lodwin Symmetric Orthonormalization of a 2-element basis set is found here: https://booksite.elsevier.com/9780444594365/downloads/16755_10030.pdf

Parameters:
  • features (TensorMap) – contains blocks whose values we wish to orthonormalize. Note that features is modified in place, so a copy of features must be made before the function if you wish to retain the unnormalized values.

  • radial_basis (RadialBasis) –

Returns:

features containing values multiplied by normalization factors.

Return type:

TensorMap

plot_basis(n_r=100)

Plot the normalized basis functions used in calculating the expansion coefficients

Parameters:

n_r (int) – number of mesh points. Default: 100

class anisoap.representations.radial_basis.MonomialBasis(max_angular, cutoff_radius, max_radial=None, rcond=1e-08, tol=0.001)

Bases: _RadialBasis

A subclass of _RadialBasis that contains attributes and methods required for the Monomial basis.

The monomial basis of order n is defined to be \(R(r) = r^n\).

For monomial basis with defined \(n_{\text{max}}\) and \(l_{\text{max}}\), our radial basis set consists of monomials of order \(n=0\) to \(n=l_{\text{max}} + 2n_{\text{max}}\).

For monomial basis with coupled \(n_{\text{max}}\) and \(l_{\text{max}}\), our radial basis set consists of monomials of order \(n=0\) to \(n=\max{l_{\text{max}} + 2n_{\text{max}}}\)

Monomials are not square-integrable from \([0, \infty]\), so we orthonormalize by integrating up to the cutoff radius

calc_overlap_matrix()

Computes the overlap matrix for Monomnials over a fixed interval.

The overlap matrix is a Gram matrix whose entries are the overlap:

\[S_{ij} = \int_{0}^{r_{cut}} dr r^2 \phi_i \phi_j\]

The overlap has an analytic solution (see above functions).

The overlap matrix is the first step to generating an orthonormal basis set of functions (Lodwin Symmetric Orthonormalization). The actual orthonormalization matrix cannot be fully precomputed because each tensor block use a different set of bases. Hence, we precompute the full overlap matrix of dim l_max, and while orthonormalizing each tensor block, we generate the respective orthonormal matrices from slices of the full overlap matrix.

Returns:

The overlap matrix

Return type:

2D array

get_basis(rs)

Evaluate orthonormalized monomial basis functions on mesh rs.

If lmax and nmax defined, then the number of functions outputted is lmax*(nmax+1)

If lmax and nmax coupled, then the number of functions outputted is \(\sum_{l=0}^{lmax} (number\_of\_radial\_functions\_per\_l)\)

Parameters:

rs (np.array) – a mesh to evaluate the basis functions

Returns:

a matrix containing orthonormalized monomial basis functions evaluated on rs

Return type:

np.array

get_num_radial_functions()

Output the number of radial basis functions considered per value of l. If max_angular and max_radial are both specified, then the list will contain repeated values of max_radial Otherwise, the outputted list will specify the number of radial basis functions per l, which may be automatically calculated if max_radial=None. If a custom list of max_radial is specified when initializing, then it will return the same inputted list.

Returns:

The attribute self.num_radial_functions, which is a list containing the number of radial basis functions considered per l.

Return type:

array_like

orthonormalize_basis(features: TensorMap)

Apply an in-place orthonormalization on the features, using Lodwin Symmetric Orthonormalization. Each block in the features TensorMap uses a basis set of l + 2n, so we must take the appropriate slices of the overlap matrix to compute the orthonormalization matrix. An instructive example of Lodwin Symmetric Orthonormalization of a 2-element basis set is found here: https://booksite.elsevier.com/9780444594365/downloads/16755_10030.pdf

Parameters:
  • features (TensorMap) – A TensorMap whose blocks’ values we wish to orthonormalize. Note that features is modified in place, so a copy of features must be made before the function if you wish to retain the unnormalized values.

  • radial_basis (_RadialBasis) –

Returns:

features containing values multiplied by proper normalization factors.

Return type:

TensorMap

plot_basis(n_r=100)

Plot the normalized basis functions used in calculating the expansion coefficients

Parameters:

n_r (int) – number of mesh points. Default: 100

anisoap.representations.radial_basis.gto_overlap(n, m, sigma_n, sigma_m)

Compute overlap of two normalized GTOs

Note that the overlap of two GTOs can be modeled as the square norm of one GTO, with an effective \(n\) and \(\sigma\). All we need to do is to calculate those effective parameters, then compute the normalization prefactor.

\[\begin{split}\langle \phi_n, \phi_m \rangle &= \int_0^\infty dr \, r^2 r^n e^{-r^2/(2\sigma_n^2)} \, r^m e^{-r^2/(2\sigma_m^2)} \\ &= \int_0^\infty dr \, r^2 \lvert r^{(n+m)/2} e^{-r^2/4 (1/\sigma_n^2 + 1/\sigma_m^2)}\rvert^2 \\ &= \int_0^\infty dr \, r^2 r^n_\text{eff} e^{-r^2/(2\sigma_\text{eff}^2)}\end{split}\]
Parameters:
  • n – order of the first GTO

  • m – order of the second GTO

  • sigma_n – sigma parameter of the first GTO

  • sigma_m – sigma parameter of the second GTO

Returns:

overlap of the two normalized GTOs

Return type:

float

anisoap.representations.radial_basis.gto_prefactor(n, sigma)

Computes the normalization prefactor of an unnormalized GTO.

This prefactor is simply \(\frac{1}{\sqrt{\text{square_norm_area}}}\). Scaling a GTO by this prefactor will ensure that the GTO has square norm equal to 1.

Parameters:
  • n – order of GTO

  • sigma – width of GTO

Returns:

The normalization constant

Return type:

float

anisoap.representations.radial_basis.gto_square_norm(n, sigma)

Compute the square norm of GTOs (inner product of itself over \(R^3\)).

An unnormalized GTO of order n is \(\phi_n = r^n e^{-r^2/(2*\sigma^2)}\) The square norm of the unnormalized GTO has an analytic solution:

\[\begin{split}\braket{\phi_n | \phi_n} &= \int_0^\infty dr \, r^2 \lvert\phi_n\rvert^2 \\ &= \frac{1}{2} \sigma^{2n+3} \Gamma(n+\frac{3}{2})\end{split}\]

This function uses the above expression.

Parameters:
  • n – order of the GTO

  • sigma – width of the GTO

Returns:

The square norm of the unnormalized GTO

Return type:

float

anisoap.representations.radial_basis.inverse_matrix_sqrt(matrix: array, rcond=1e-08, tol=0.001)

Returns the inverse matrix square root.

The inverse square root of the overlap matrix (or slices of the overlap matrix) yields the orthonormalization matrix

Parameters:
  • matrix (np.array) – Symmetric square matrix to find the inverse square root of

  • rcond (float) – Lower bound for eigenvalues for inverse square root

  • tol (float) – Tolerance for differences between original matrix and reconstruction via inverse square root

Returns:

\(S^{-1/2}\)

Return type:

inverse_sqrt_matrix

anisoap.representations.radial_basis.monomial_overlap(n, m, r_cut)

Compute overlap of two normalized monomials

Note that the overlap of two monomials can be modeled as the square norm of one monomial, with an effective n. All we need to do is to calculate those effective parameters, then compute the normalization:

\[\begin{split}\langle \phi_n, \phi_m \rangle &= \int_0^\infty dr r^2 r^n r^m \\ &= \int_0^{r_{cut}} dr r^2 |r^{(n+m)/2}|^2 \\ &= \int_0^{r_{cut}} dr r^2 |r^{n_{eff}}|^2\end{split}\]
Parameters:
  • n – order of the first monomial

  • m – order of the second monomial

  • r_cut (float) – cutoff radius

Returns:

overlap of the two normalized monomial

Return type:

float

anisoap.representations.radial_basis.monomial_prefactor(n, r_cut)

Computes the normalization prefactor of an unnormalized monomial basis. This prefactor is simply \(1/\sqrt{square\_norm\_area}\). Scaling a basis by this prefactor will ensure that the basis has square norm equal to 1.

Parameters:

n – order of the basis

Returns:

normalization constant

Return type:

float

anisoap.representations.radial_basis.monomial_square_norm(n, r_cut)

Compute the square norm of monomials (inner product of itself over R^3).

Parameters:

n – order of the basis

Returns:

The square norm of the unnormalized basis

Return type:

float

Utilities

class anisoap.utils.cyclic_list.CGRCacheList(size: int)

Bases: object

This is a simple class that only exists to be used as a “private” cache for computed ClebschGordanReal matrices (self._cg)

It only stores last n CG matrices computed for distinct l_max values. Since it is specialized to store (l_max, cg matrix) pair, it is NOT supposed to be used to cache any other things. While type of “value” (in (key, value) pair) does not matter as much, the type of key MUST BE A 32-BIT INTEGER with most significant bit unused, as it performs bitwise operations on the most significant bit to store and check the replacement flags.

It will replace the entries using something called “clock algorithm” for page replacement, which mimics the “least recently used” algorithm but with less computation requirement.

[Link to clock algorithm](https://en.wikipedia.org/wiki/Page_replacement_algorithm#Clock)

clear_cache() None

Resets the cache (makes list empty)

get_val(key: int) dict

Obtains the value for given key. Raises IndexError if the given key is not found in the list.

insert(key: int, value: dict) None

This insert algorithm mimics the “clock algorithm” for page replacement technique. The algorithm works like this: 1. Get the entry that the “clock hand” (self._ins_index) points to. 2. Keep advancing the “clock hand” until key with replacement flag = 0 is found. Replacement flag, in this case, is the most significant bit 3. Insert the new [key, value] pair into the list, with replacement flag = 1. Advance “clock hand” by one position. Note that clock hand wraps around the list, should it reach the end of the list.

keys() List[int]

Return list of keys, with keys

anisoap.utils.metatensor_utils.cg_combine(x_a, x_b, feature_names=None, clebsch_gordan=None, lcut=None, other_keys_match=None)

Performs a CG product of two sets of equivariants.

The only requirement is that sparse indices are labeled as (“inversion_sigma”, “spherical_harmonics_l”, “order_nu”).

Parameters:
  • x_a – First set of equivariants

  • x_b – Second set of equivariants

  • feature_names (list, optional) – Overrides automatically-generated names of output features. By default, all other key labels are combined via their outer product, i.e. if there is a key-side neighbor-species in both x_a and x_b, the returned keys will have two neighbor_species labels, corresponding to the parent features.

  • clebsch_gordan (ClebschGordanReal, optional) –

  • lcut (np.ndarray()) – Cutoff in new features

  • other_keys_match (list, optional) – List of keys that should match. These will not need to have their outer product taken, but will instead be merged into a new key. For instance, passing [“types center”] will combine the keys with the same type center, yielding a single key with the same types_center in the results.

Returns:

The Clebsch-Gordan product of x_a and x_b

Return type:

TensorMap

anisoap.utils.metatensor_utils.standardize_keys(descriptor)

Standardize the naming scheme of density expansion coefficient blocks (nu=1)

anisoap.utils.moment_generator.compute_moments_diagonal_inefficient_implementation(principal_components, a, maxdeg)

Computes all moments for a diagonal dilation matrix.

The implementation focuses on conceptual simplicity, while sacrificing memory efficiency. To be specific, the moments array allows access to the value of the moment \(\langle x^{n_0} y^{n_1} z^{n_2} \rangle\) simply as moments[n0, n1, n2]. This leads to more intuitive code, at the cost of wasting around a third of the memory in the array to store zeros.

Parameters:
  • principal_components (np.ndarray of shape (3,)) – Array containing the three principal components

  • a (np.ndarray of shape (3,)) – Vectorial center of the trivariate Gaussian

  • maxdeg (int) – Maximum degree for which the moments need to be computed.

Returns:

Moments calculated. moments[n0,n1,n2] is the \(\left(n_0,n_1,n_2\right)^\text{th}\) moment of the Gaussian defined as

\[\langle x^{n_0} y^{n_1} z^{n_2} \rangle = \int(x^{n_0} y^{n_1} z^{n_2}) e^{-\frac{1}{2}(r-a)^T \Sigma (r-a)} dx\,dy\,dz\, \sum_{i=1}^{\infty} x_{i}\]

where \(\Sigma\) is the covariance matrix.

Return type:

np.ndarray of shape (3, maxdeg + 1)

Note

The term “moments” in probability theory are defined for normalized Gaussian distributions. Here, we take the Gaussian without prefactor, meaning that all moments are scaled by a global factor.

anisoap.utils.moment_generator.compute_moments_inefficient_implementation(A, a, maxdeg)

Computes all moments for a general dilation matrix.

Parameters:
  • A (np.ndarray of shape (3,3)) – Dilation matrix of the Gaussian that determines its shape. It can be written as \(\mathbf{A} = \mathbf{R}\mathbf{D}\mathbf{R}^T\), where R is a rotation matrix that specifies the orientation of the three principal axes, while \(\mathbf{D}\) is a diagonal matrix whose three diagonal elements are the lengths of the principal axes.

  • a (np.ndarray of shape (3,)) – Vectorial center of the trivariate Gaussian.

  • maxdeg (int) – Maximum degree for which the moments need to be computed.

Returns:

The list of moments defined as

\[\langle x^{n_0} y^{n_1} z^{n_2} \rangle = \int(x^{n_0} y^{n_1} z^{n_2}) \exp(-0.5(r-a)^T \Sigma (r-a)) dx\,dy\,dz\, \sum_{i=1}^{\infty} x_{i}\]

where \(\Sigma\) is the covariance matrix.

Return type:

np.ndarray

Note

Note that the term “moments” in probability theory are defined for normalized Gaussian distributions. These recursive relations find the moments of a normalized Gaussian, but we actually want to find the moments of the unnormalized underlying gaussian (as seen in the equation above) because calculating expansion coefficients of any density should not involve normalization. Therefore, we scale the end-result by global_factor, which is the reciprocal of the normalizing constant in front of any gaussian.

anisoap.utils.moment_generator.compute_moments_single_variable(A, a, maxdeg)

Computes all moments for a single variable Gaussian.

Parameters:
  • A (float) – inverse of variance (see below for exact mathematical form)

  • a (float) – Center of Gaussian function

  • maxdeg (int) – Maximum degree for which the moments need to be computed.

Returns:

A numpy array of size (maxdeg + 1,) containing the moments defined as

\[\langle x^n \rangle = \int x^n e^{-A(x-a)^2/2} dx\]

Return type:

np.array

Note

The Gaussian is not normalized, meaning that we need to multiply all results by a global factor if we wish to interpret these as moments of a probability density.

class anisoap.utils.monomial_iterator.TrivariateMonomialIndices(deg)

Bases: object

Class for generating an iterator object over trivariate monomials.

Generates an iterator object over trivariate monomials of the form \(f(x,y,z) = x^{n_0} y^{n_1} z^{n_2}\) sorted in lexicographical order.

Without this class, iterating over all monomials at some fixed degree requires the use of a double loop of the form:

idx = 0
for n0 in range(deg, -1, -1):
        for n1 in range(deg-n0, -1, -1):
        n2 = deg - n0 - n1

        ... # do something with exponents (n0,n1,n2)

        idx += 1

Instead, with this class, these lines can be reduced to:

myiter = iter(TrivariateMonomialIndices(deg=2))
for idx, n0, n1, n2 in myiter:
    ... # do something with exponents (n0, n1, n2)
find_idx(exponents)

Finds the index of a given tuple of exponents.

Parameters:

exponents (3-tuple of the form (n0, n1, n2)) – The exponents of the monomial \(x^{n_0} y^{n_1} z^{n_2}\)

Return type:

The index of the tuple in lexicographical order

anisoap.utils.spherical_to_cartesian.prefact_minus1(l)

Computes the prefactor that multiplies the \(T_{l-1}^\text{th}\) term in the iteration.

For \(m \in \left[-l, -l+2, ..., l \right]\), compute the factor as

\[\left( \frac{\sqrt{(l+1-m)!}}{(l+1+m)!} \right) \left( \frac{\sqrt{(l+m)!}}{(l-m)!} \right) \left( \frac{2l+1}{l+1-m} \right)\]
Parameters:

l (int) – Term immediately proceeding the term for which the prefactor is computed.

Returns:

corresponds to the prefactor that multiplies the \(T_{l-1}^\text{th}\) term in the iteration

Return type:

list of size (2l + 1)

anisoap.utils.spherical_to_cartesian.prefact_minus2(l)

Computes the prefactor that multiplies the \(T_{l-2}^\text{th}\) term in the iteration.

For \(m \in \left[-l+1, -l+2, ..., l-1\right]\), compute the factor as

\[\left( \frac{\sqrt{(l+1-m)!}}{(l+1+m)!} \right) \left(\frac{\sqrt{(l-1+m)!}}{(l-1-m)!} \right) \left( \frac{l+m}{l-m+1} \right)\]
Parameters:

l (int) – Term two places after the term for which the prefactor is computed

Returns:

Corresponds to the prefactor that multiplies the term in question

Return type:

list of size (2l - 1)

anisoap.utils.spherical_to_cartesian.spherical_to_cartesian(lmax, num_ns)

Finds the coefficients for the cartesian polynomial form of solid harmonics \(R_{lm} = \sqrt{(4\pi)/(2l+1)}*r^l*Y_{lm}\). Note that our AniSOAP expansion does not contain the \(\sqrt{(4\pi)/(2l+1)}\) prefactor, so in calculating expansion coefficients, we need to divide by that coefficient.

Parameters:
  • lmax (int) –

  • num_ns (int) –