CUDA C++ API

Defines the CUDA C++ API for sphericart. Two classes are available: SphericalHarmonics and SolidHarmonics. The former calculates the real spherical harmonics \( Y^m_l \) as defined on Wikipedia, which are homogeneous polynomials of (x/r, y/r, z/r). The latter calculates the same polynomials but as a function of the Cartesian coordinates (x, y, z), or, equivalently, \( r^l\,Y^m_l \).

namespace cuda

The sphericart::cuda namespace contains the CUDA API for sphericart.

template<typename T>
class SphericalHarmonics
#include <sphericart_cuda.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.

Subclassed by cuda::SolidHarmonics< T >

Public Functions

SphericalHarmonics(size_t l_max)

Initialize the SphericalHarmonics class setting maximum degree and normalization

Parameters:

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

void compute(const T *xyz, const size_t n_samples, T *sph, void *cuda_stream = nullptr)

Default constructor Required so sphericart_torch can conditionally instantiate this class depending on if cuda is available. Computes the spherical harmonics for one or more 3D points, using pre-allocated device-side pointers

Parameters:
  • xyz – A pre-allocated device-side 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.

  • n_samples – Number of samples contained within xyz.

  • sph – On entry, a preallocated device-side array of size n_samples x (l_max + 1)^2. 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)^2 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.

  • cuda_stream – Pointer to a cudaStream_t or nullptr. If this is nullptr, the kernel launch will be performed on the default stream.

void compute_with_gradients(const T *xyz, const size_t n_samples, T *sph, T *dsph, void *cuda_stream = nullptr)

Computes the spherical harmonics and their first derivatives for one or more 3D points, using pre-allocated device-side pointers

Parameters:
  • xyz – A pre-allocated device-side 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.

  • n_samples – Number of samples contained within xyz.

  • sph – On entry, a preallocated device-side array of size n_samples x (l_max + 1)^2. 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)^2 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.

  • dsph – On entry, nullptr or a preallocated device-side array of size n_samples x 3 x (l_max + 1)^2. If the pointer is not nullptr, then compute_with_gradients must also be true in order for gradients to be computed.

  • cuda_stream – Pointer to a cudaStream_t or nullptr. If this is nullptr, the kernel launch will be performed on the default stream.

void compute_with_hessians(const T *xyz, const size_t n_samples, T *sph, T *dsph, T *ddsph, void *cuda_stream = nullptr)

Computes the spherical harmonics and their first and second derivatives for one or more 3D points, using pre-allocated device-side pointers

Parameters:
  • xyz – A pre-allocated device-side 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.

  • n_samples – Number of samples contained within xyz.

  • compute_with_gradients – Whether we should compute dsph. If true, the pointer dsph must also be allocated on device.

  • compute_with_hessians – Whether we should compute ddsph. If true, the pointer ddsph must also be allocated on device.

  • sph – On entry, a preallocated device-side array of size n_samples x (l_max + 1)^2. 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)^2 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.

  • dsph – On entry, nullptr or a preallocated device-side array of size n_samples x 3 x (l_max + 1)^2. If the pointer is not nullptr, then compute_with_gradients must also be true in order for gradients to be computed.

  • ddsph – On entry, nullptr or a preallocated device-side array of size n_samples x 3 x 3 x (l_max + 1)^2. If the pointer is not nullptr, then compute_with_hessians must also be true in order for gradients to be computed.

  • cuda_stream – Pointer to a cudaStream_t or nullptr. If this is nullptr, the kernel launch will be performed on the default stream.

Friends

friend class SolidHarmonics
template<typename T>
class SolidHarmonics : public cuda::SphericalHarmonics<T>
#include <sphericart_cuda.hpp>

A solid harmonics calculator.

Its interface is the same as that of the SphericalHarmonics class, but it calculates the solid harmonics \( r^l\,Y^m_l \) instead of the real spherical harmonics \( Y^m_l \), allowing for faster computations.

Public Functions

SolidHarmonics(size_t l_max)

Initialize the SolidHarmonics class setting its maximum degree

Parameters:

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