CUDA C++ API#
Defines the CUDA API for sphericart
.
-
namespace cuda#
The
sphericart::cuda
namespace contains the CUDA API forsphericart
.-
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.
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. Iftrue
, 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 T *xyz, size_t nsamples, bool compute_with_gradients, bool compute_with_hessian, T *sph, T *dsph = nullptr, T *ddsph = nullptr, void *cuda_stream = nullptr)#
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 isn_samples
long, accounting for different samples, while the inner dimension has size 3 and it represents the x, y, and z coordinates respectively.nsamples – 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
dsph –
On entry, nullptr or a preallocated device-side array of size . 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 . 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.
-
SphericalHarmonics(size_t l_max, bool normalized = false)#
-
template<typename T>