C++ API#

Defines the C++ API for sphericart.

namespace sphericart#
template<typename T>
class SphericalHarmonics#
#include <sphericart.hpp>

A spherical harmonics calculator.

It handles initialization of the prefactors upon initialization and it stores the buffers that are necessary to compute the spherical harmonics efficiently.

Public Functions

SphericalHarmonics(size_t l_max, bool normalized = false)#

Initialize the SphericalHarmonics class setting maximum degree and normalization

Parameters:
  • l_max – The maximum degree of the spherical harmonics to be calculated.

  • normalized – If false (default) computes the scaled spherical harmonics, which are homogeneous polynomials in the Cartesian coordinates of the input points. If true, computes the normalized spherical harmonics that are evaluated on the unit sphere. In practice, this simply computes the scaled harmonics at the normalized coordinates \((x/r, y/r, z/r)\), and adapts the derivatives accordingly.

void compute(const std::vector<T> &xyz, std::vector<T> &sph)#

Computes the spherical harmonics for one or more 3D points, using std::vectors.

Parameters:
  • xyz – A std::vector array of size n_samples x 3. It contains the Cartesian coordinates of the 3D points for which the spherical harmonics are to be computed, organized along two dimensions. The outer dimension is n_samples long, accounting for different samples, while the inner dimension has size 3 and it represents the x, y, and z coordinates respectively. If xyz it contains a single point, the class will call a simpler function that directly evaluates the point, without a loop.

  • sph – On entry, a (possibly uninitialized) std::vector, which, if needed, will be resized to n_samples * (l_max + 1) * (l_max + 1). On exit, this array will contain the spherical harmonics organized along two dimensions. The leading dimension is n_samples long and it represents the different samples, while the inner dimension is (l_max + 1) * (l_max + 1) long and it contains the spherical harmonics. These are laid out in lexicographic order. For example, if l_max=2, it will contain (l, m) = (0, 0), (1, -1), (1, 0), (1, 1), (2, -2), (2, -1), (2, 0), (2, 1), (2, 2), in this order.

void compute_with_gradients(const std::vector<T> &xyz, std::vector<T> &sph, std::vector<T> &dsph)#

Computes the spherical harmonics and their derivatives with respect to the Cartesian coordinates of one or more 3D points, using std::vectors.

Parameters:
  • xyz – A std::vector array of size n_samples x 3. It contains the Cartesian coordinates of the 3D points for which the spherical harmonics are to be computed, organized along two dimensions. The outer dimension is n_samples long, accounting for different samples, while the inner dimension has size 3 and it represents the x, y, and z coordinates respectively. If xyz it contains a single point, the class will call a simpler functions that directly evaluates the point, without a loop.

  • sph – On entry, a (possibly uninitialized) std::vector, which, if needed, will be resized to n_samples * (l_max + 1) * (l_max + 1). On exit, this array will contain the spherical harmonics organized along two dimensions. The leading dimension is n_samples long and it represents the different samples, while the inner dimension size is (l_max + 1) * (l_max + 1) long and it contains the spherical harmonics. These are laid out in lexicographic order. For example, if l_max=2, it will contain (l, m) = (0, 0), (1, -1), (1, 0), (1, 1), (2, -2), (2, -1), (2, 0), (2, 1), (2, 2), in this order.

  • dsphstd::vector for spherical harmonics derivatives. It is a (possibly uninitialized) std::vector, which, if needed, will be resized to n_samples * 3 * (l_max + 1) * (l_max + 1). On exit, this array will contain the derivatives of the spherical harmonics organized along three dimensions. As for the sph parameter, the leading dimension represents the different samples, while the inner-most dimension size is (l_max + 1) * (l_max + 1), and it represents the degree and order of the spherical harmonics (again, organized in lexicographic order). The intermediate dimension corresponds to different spatial derivatives of the spherical harmonics: x, y, and z, respectively.

void compute_with_hessians(const std::vector<T> &xyz, std::vector<T> &sph, std::vector<T> &dsph, std::vector<T> &ddsph)#

Computes the spherical harmonics, their derivatives and second derivatives with respect to the Cartesian coordinates of one or more 3D points, using std::vectors.

Parameters:
  • xyz – A std::vector array of size n_samples x 3. It contains the Cartesian coordinates of the 3D points for which the spherical harmonics are to be computed, organized along two dimensions. The outer dimension is n_samples long, accounting for different samples, while the inner dimension has size 3 and it represents the x, y, and z coordinates respectively. If xyz it contains a single point, the class will call a simpler functions that directly evaluates the point, without a loop.

  • sph – On entry, a (possibly uninitialized) std::vector, which, if needed, will be resized to n_samples * (l_max + 1) * (l_max + 1). On exit, this array will contain the spherical harmonics organized along two dimensions. The leading dimension is n_samples long and it represents the different samples, while the inner dimension size is (l_max + 1) * (l_max + 1) long and it contains the spherical harmonics. These are laid out in lexicographic order. For example, if l_max=2, it will contain (l, m) = (0, 0), (1, -1), (1, 0), (1, 1), (2, -2), (2, -1), (2, 0), (2, 1), (2, 2), in this order.

  • dsphstd::vector for spherical harmonics derivatives. It is a (possibly uninitialized) std::vector, which, if needed, will be resized to n_samples * 3 * (l_max + 1) * (l_max + 1). On exit, this array will contain the derivatives of the spherical harmonics organized along three dimensions. As for the sph parameter, the leading dimension represents the different samples, while the inner-most dimension size is (l_max + 1) * (l_max + 1), and it represents the degree and order of the spherical harmonics (again, organized in lexicographic order). The intermediate dimension corresponds to different spatial derivatives of the spherical harmonics: x, y, and z, respectively.

  • ddsphstd::vector for spherical harmonics second derivatives. It is a (possibly uninitialized) std::vector, which, if needed, will be resized to n_samples * 3 * 3 * (l_max + 1) * (l_max + 1). On exit, this array will contain the second derivatives of the spherical harmonics organized along four dimensions. As for the sph parameter, the leading dimension represents the different samples, while the inner-most dimension size is (l_max + 1) * (l_max + 1), and it represents the degree and order of the spherical harmonics (again, organized in lexicographic order). The intermediate dimensions correspond to the different spatial second derivatives of the spherical harmonics, i.e., to the dimensions of the Hessian matrix.

void compute_array(const T *xyz, size_t xyz_length, T *sph, size_t sph_length)#

Computes the spherical harmonics for a set of 3D points using bare arrays.

Parameters:
  • xyz – An array of size n_samples x 3. It contains the Cartesian coordinates of the 3D points for which the spherical harmonics are to be computed, organized along two dimensions. The outer dimension is n_samples long, accounting for different samples, while the inner dimension has size 3 and it represents the x, y, and z coordinates respectively.

  • xyz_length – Total length of the xyz array: n_samples * 3.

void compute_array_with_gradients(const T *xyz, size_t xyz_length, T *sph, size_t sph_length, T *dsph, size_t dsph_length)#

Computes the spherical harmonics and their derivatives for a set of 3D points using bare arrays.

Parameters:
  • xyz – An array of size n_samples x 3. It contains the Cartesian coordinates of the 3D points for which the spherical harmonics are to be computed, organized along two dimensions. The outer dimension is n_samples long, accounting for different samples, while the inner dimension has size 3 and it represents the x, y, and z coordinates respectively.

  • xyz_length – Total length of the xyz array: n_samples * 3.

  • sph – On entry, an array of size n_samples * (l_max + 1) * (l_max + 1). On exit, this array will contain the spherical harmonics organized along two dimensions. The leading dimension is n_samples long and it represents the different samples, while the inner dimension size is (l_max + 1) * (l_max + 1) long and it contains the spherical harmonics. These are laid out in lexicographic order. For example, if l_max=2, it will contain (l, m) = (0, 0), (1, -1), (1, 0), (1, 1), (2, -2), (2, -1), (2, 0), (2, 1), (2, 2), in this order.

  • sph_length – Total length of the sph array: `n_samples * (l_max + 1)

    • (l_max + 1)`.

  • dsph – Array for spherical harmonics derivatives. It is an array of size n_samples * 3 * (l_max + 1) * (l_max + 1). On exit, this array will contain the derivatives of the spherical harmonics organized along three dimensions. As for the sph parameter, the leading dimension represents the different samples, while the inner-most dimension size is (l_max + 1) * (l_max + 1), and it represents the degree and order of the spherical harmonics (again, organized in lexicographic order). The intermediate dimension corresponds to the different spatial derivatives of the spherical harmonics: x, y, and z, respectively.

  • dsph_length – Total length of the dsph array: n_samples * 3 * (l_max + 1) * (l_max + 1).

void compute_array_with_hessians(const T *xyz, size_t xyz_length, T *sph, size_t sph_length, T *dsph, size_t dsph_length, T *ddsph, size_t ddsph_length)#

Computes the spherical harmonics, their derivatives and second derivatives for a set of 3D points using bare arrays.

Parameters:
  • xyz – An array of size n_samples x 3. It contains the Cartesian coordinates of the 3D points for which the spherical harmonics are to be computed, organized along two dimensions. The outer dimension is n_samples long, accounting for different samples, while the inner dimension has size 3 and it represents the x, y, and z coordinates respectively.

  • xyz_length – Total length of the xyz array: n_samples * 3.

  • sph – On entry, an array of size n_samples * (l_max + 1) * (l_max + 1). On exit, this array will contain the spherical harmonics organized along two dimensions. The leading dimension is n_samples long and it represents the different samples, while the inner dimension size is (l_max + 1) * (l_max + 1) long and it contains the spherical harmonics. These are laid out in lexicographic order. For example, if l_max=2, it will contain (l, m) = (0, 0), (1, -1), (1, 0), (1, 1), (2, -2), (2, -1), (2, 0), (2, 1), (2, 2), in this order.

  • sph_length – Total length of the sph array: `n_samples * (l_max + 1)

    • (l_max + 1)`.

  • dsph – Array for spherical harmonics derivatives. It is an array of size n_samples * 3 * (l_max + 1) * (l_max + 1). On exit, this array will contain the derivatives of the spherical harmonics organized along three dimensions. As for the sph parameter, the leading dimension represents the different samples, while the inner-most dimension size is (l_max + 1) * (l_max + 1), and it represents the degree and order of the spherical harmonics (again, organized in lexicographic order). The intermediate dimension corresponds to the different spatial derivatives of the spherical harmonics: x, y, and z, respectively.

  • dsph_length – Total length of the dsph array: n_samples * 3 * (l_max + 1) * (l_max + 1).

  • ddsph – Array for spherical harmonics second derivatives. It is an array of size n_samples * 3 * 3 * (l_max + 1) * (l_max + 1). On exit, this array will contain the second derivatives of the spherical harmonics organized along four dimensions. As for the sph parameter, the leading dimension represents the different samples, while the inner-most dimension size is (l_max + 1) * (l_max + 1), and it represents the degree and order of the spherical harmonics (again, organized in lexicographic order). The intermediate dimensions correspond to the different spatial second derivatives of the spherical harmonics, i.e., to the dimensions of the Hessian matrix.

  • ddsph_length – Total length of the ddsph array: n_samples * 9 * (l_max + 1) * (l_max + 1).

void compute_sample(const T *xyz, size_t xyz_length, T *sph, size_t sph_length)#

Computes the spherical harmonics for a single 3D point using bare arrays.

Parameters:
  • xyz – An array of size 3. It contains the Cartesian coordinates of the 3D point for which the spherical harmonics are to be computed. x, y, and z coordinates respectively.

  • xyz_length – Length of the xyz array: 3.

  • sph – On entry, an array of size (l_max + 1) * (l_max + 1). On exit, this array will contain the spherical harmonics laid out in lexicographic order. For example, if l_max=2, it will contain the spherical harmonics in the following order: (l, m) = (0, 0), (1, -1), (1, 0), (1, 1), (2, -2), (2, -1), (2, 0), (2, 1), (2, 2).

  • sph_length – Total length of the sph array: (l_max + 1) * (l_max

void compute_sample_with_gradients(const T *xyz, size_t xyz_length, T *sph, size_t sph_length, T *dsph, size_t dsph_length)#

Computes the spherical harmonics and their derivatives for a single 3D point using bare arrays.

Parameters:
  • xyz – An array of size 3. It contains the Cartesian coordinates of the 3D point for which the spherical harmonics are to be computed. x, y, and z coordinates respectively.

  • xyz_length – Length of the xyz array: 3.

  • sph – On entry, an array of size (l_max + 1) * (l_max + 1). On exit, this array will contain the spherical harmonics laid out in lexicographic order. For example, if l_max=2, it will contain the spherical harmonics in the following order: (l, m) = (0, 0), (1, -1), (1, 0), (1, 1), (2, -2), (2, -1), (2, 0), (2, 1), (2, 2).

  • sph_length – Total length of the sph array: (l_max + 1) * (l_max

  • dsphArray for spherical harmonics derivatives. It is an array of size . On exit, this array will contain the spherical harmonics' derivatives organized along two dimensions. The second dimension's size is , and it represents the degree and order of the spherical harmonics (again, organized in lexicographic order). The first dimension corresponds to the different spatial derivatives of the spherical harmonics: x, y, and z, respectively.

  • dsph_lengthTotal length of the array: .

void compute_sample_with_hessians(const T *xyz, size_t xyz_length, T *sph, size_t sph_length, T *dsph, size_t dsph_length, T *ddsph, size_t ddsph_length)#

Computes the spherical harmonics, their derivatives and second derivatives for a single 3D point using bare arrays.

Parameters:
  • xyz – An array of size 3. It contains the Cartesian coordinates of the 3D point for which the spherical harmonics are to be computed. x, y, and z coordinates respectively.

  • xyz_length – Length of the xyz array: 3.

  • sph – On entry, an array of size (l_max + 1) * (l_max + 1). On exit, this array will contain the spherical harmonics laid out in lexicographic order. For example, if l_max=2, it will contain the spherical harmonics in the following order: (l, m) = (0, 0), (1, -1), (1, 0), (1, 1), (2, -2), (2, -1), (2, 0), (2, 1), (2, 2).

  • sph_length – Total length of the sph array: (l_max + 1) * (l_max

  • dsphArray for spherical harmonics derivatives. It is an array of size . On exit, this array will contain the spherical harmonics' derivatives organized along two dimensions. The second dimension's size is , and it represents the degree and order of the spherical harmonics (again, organized in lexicographic order). The first dimension corresponds to the different spatial derivatives of the spherical harmonics: x, y, and z, respectively.

  • dsph_lengthTotal length of the array: .

  • ddsphArray for spherical harmonics second derivatives. It is an array of size . On exit, this array will contain the second derivatives of the spherical harmonics organized along three dimensions. As for the parameter, the inner-most dimension size is , and it represents the degree and order of the spherical harmonics (again, organized in lexicographic order). The first two dimensions correspond to the different spatial second derivatives of the spherical harmonics, i.e., to the dimensions of the Hessian matrix.

  • ddsph_lengthTotal length of the array: .

inline int get_omp_num_threads()#

Returns the number of threads used in the calculation