C API#

Defines the C API for sphericart.

Typedefs

typedef struct sphericart_calculator_t sphericart_calculator_t#

A type referring to the sphericart_calculator_t struct.

typedef struct sphericart_calculator_f_t sphericart_calculator_f_t#

A type referring to the sphericart_calculator_f_t struct.

Functions

sphericart_calculator_t *sphericart_new(size_t l_max, bool normalized)#

Initializes a spherical harmonics calculator and returns a pointer that can then be used by functions that evaluate spherical harmonics over arrays or individual samples.

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

  • normalized – If false, computes the scaled spherical harmonics, which are 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.

Returns:

A pointer to a sphericart_calculator_t object

sphericart_calculator_f_t *sphericart_new_f(size_t l_max, bool normalized)#

Similar to sphericart_new, but it returns a sphericart_calculator_f_t, which performs calculations on the float type.

void sphericart_delete(sphericart_calculator_t *calculator)#

Deletes a previously allocated sphericart_calculator_t calculator.

void sphericart_delete_f(sphericart_calculator_f_t *calculator)#

Deletes a previously allocated sphericart_calculator_f_t calculator.

void sphericart_compute_array(sphericart_calculator_t *calculator, const double *xyz, size_t xyz_length, double *sph, size_t sph_length)#

This function calculates the spherical harmonics and, optionally, their derivatives for an array of 3D points.

Parameters:
  • calculator – A pointer to a sphericart_calculator_t struct that holds prefactors and options to compute the spherical harmonics.

  • 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 – size of the xyz allocation, i.e, 3 * n_samples

  • sph – pointer to the first element of an array containing n_samples

  • sph_lengthsize of the sph allocation, should be

void sphericart_compute_array_with_gradients(sphericart_calculator_t *calculator, const double *xyz, size_t xyz_length, double *sph, size_t sph_length, double *dsph, size_t dsph_length)#

This function calculates the spherical harmonics and their derivatives for an array of 3D points.

Parameters:
  • calculator – A pointer to a sphericart_calculator_t struct that holds prefactors and options to compute the spherical harmonics.

  • 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 – size of the xyz allocation, i.e, 3 * n_samples

  • sph – pointer to the first element of an array containing n_samples

  • sph_lengthsize of the sph allocation, should be

  • dsphpointer to the first element of an array containing n_samples * 3 * (l_max + 1) * (l_max + 1) elements. On exit, this array will contain the spherical harmonics derivatives organized along three dimensions. As for the parameter, the leading dimension represents the different samples, while the inner-most dimension size is , 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.

  • dsph_lengthsize of the dsph allocation, which should be

void sphericart_compute_array_with_hessians(sphericart_calculator_t *calculator, const double *xyz, size_t xyz_length, double *sph, size_t sph_length, double *dsph, size_t dsph_length, double *ddsph, size_t ddsph_length)#

This function calculates the spherical harmonics, their derivatives and second derivatives for an array of 3D points.

Parameters:
  • calculator – A pointer to a sphericart_calculator_t struct that holds prefactors and options to compute the spherical harmonics.

  • 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 – size of the xyz allocation, i.e, 3 * n_samples

  • sph – pointer to the first element of an array containing n_samples

  • sph_lengthsize of the sph allocation, should be

  • dsphpointer to the first element of an array containing elements. On exit, this array will contain the spherical harmonics' derivatives organized along three dimensions. As for the parameter, the leading dimension represents the different samples, while the inner-most dimension size is , 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.

  • dsph_lengthsize of the dsph allocation, which should be

  • ddsph

  • ddsph_length

void sphericart_compute_sample(sphericart_calculator_t *calculator, const double *xyz, size_t xyz_length, double *sph, size_t sph_length)#

Similar to :func:sphericart_compute_array, but it computes the spherical harmonics for a single 3D point in space.

void sphericart_compute_sample_with_gradients(sphericart_calculator_t *calculator, const double *xyz, size_t xyz_length, double *sph, size_t sph_length, double *dsph, size_t dsph_length)#

Similar to :func:sphericart_compute_array_with_gradients, but it computes the spherical harmonics for a single 3D point in space.

void sphericart_compute_sample_with_hessians(sphericart_calculator_t *calculator, const double *xyz, size_t xyz_length, double *sph, size_t sph_length, double *dsph, size_t dsph_length, double *ddsph, size_t ddsph_length)#

Similar to :func:sphericart_compute_array_with_hessians, but it computes the spherical harmonics for a single 3D point in space.

void sphericart_compute_array_f(sphericart_calculator_f_t *calculator, const float *xyz, size_t xyz_length, float *sph, size_t sph_length)#

Similar to :func:sphericart_compute_array, but using the float data type.

void sphericart_compute_array_with_gradients_f(sphericart_calculator_f_t *calculator, const float *xyz, size_t xyz_length, float *sph, size_t sph_length, float *dsph, size_t dsph_length)#

Similar to :func:sphericart_compute_array_with_gradients, but using the float data type.

void sphericart_compute_array_with_hessians_f(sphericart_calculator_f_t *calculator, const float *xyz, size_t xyz_length, float *sph, size_t sph_length, float *dsph, size_t dsph_length, float *ddsph, size_t ddsph_length)#

Similar to :func:sphericart_compute_array_with_hessians, but using the float data type.

void sphericart_compute_sample_f(sphericart_calculator_f_t *calculator, const float *xyz, size_t xyz_length, float *sph, size_t sph_length)#

Similar to :func:sphericart_compute_sample, but using the float data type.

void sphericart_compute_sample_with_gradients_f(sphericart_calculator_f_t *calculator, const float *xyz, size_t xyz_length, float *sph, size_t sph_length, float *dsph, size_t dsph_length)#

Similar to :func:sphericart_compute_sample_with_gradients, but using the float data type.

void sphericart_compute_sample_with_hessians_f(sphericart_calculator_f_t *calculator, const float *xyz, size_t xyz_length, float *sph, size_t sph_length, float *dsph, size_t dsph_length, float *ddsph, size_t ddsph_length)#

Similar to :func:sphericart_compute_sample_with_hessians, but using the float data type.

int sphericart_omp_num_threads(sphericart_calculator_t *calculator)#

Get the number of OpenMP threads used by a calculator. If sphericart is computed without OpenMP support returns 1.

int sphericart_omp_num_threads_f(sphericart_calculator_f_t *calculator)#