Python API#
- class sphericart.SphericalHarmonics(l_max: int, normalized: bool = False)#
Spherical harmonics calculator, up to degree
l_max
.By default, this class computes a non-normalized form of the real spherical harmonics, i.e. \(r^l Y^l_m\). These scaled spherical harmonics are polynomials in the Cartesian coordinates of the input points.
normalize=True
can be set to compute \(Y^l_m\).In order to minimize the cost of each call, the SphericalHarmonics object computes prefactors and initializes buffers upon creation
>>> import numpy as np >>> import sphericart as sc >>> sh = sc.SphericalHarmonics(l_max=8, normalized=False)
Then, the
compute()
method can be called on an array of 3D Cartesian points to compute the spherical harmonics>>> xyz = np.random.normal(size=(10,3)) >>> sh_values = sh.compute(xyz) >>> sh_values.shape (10, 81)
In order to also compute derivatives, you can use
>>> sh_values, sh_grads = sh.compute_with_gradients(xyz) >>> sh_grads.shape (10, 3, 81)
which returns the gradient as a tensor with size (n_samples, 3, (l_max+1)**2).
- Parameters:
l_max – the maximum degree of the spherical harmonics to be calculated
normalized – whether to normalize the spherical harmonics (default: False)
- Returns:
a calculator, in the form of a SphericalHarmonics object
- compute(xyz: ndarray) ndarray #
Calculates the spherical harmonics for a set of 3D points, whose coordinates are in the
xyz
array.>>> import numpy as np >>> import sphericart as sc >>> sh = sc.SphericalHarmonics(l_max=8, normalized=False) >>> xyz = np.random.normal(size=(10,3)) >>> sh_values = sh.compute(xyz) >>> sh_values.shape (10, 81)
- Parameters:
xyz – The Cartesian coordinates of the 3D points, as an array with shape
(n_samples, 3)
- Returns:
An array of shape
(n_samples, (l_max+1)**2)
containing all the spherical harmonics up to degree l_max in lexicographic order. For example, ifl_max = 2
, The last axis will correspond to spherical harmonics with(l, m) = (0, 0), (1, -1), (1, 0), (1, 1), (2, -2), (2, -1), (2, 0), (2, 1), (2, 2)
, in this order.
- compute_with_gradients(xyz: ndarray) Tuple[ndarray, ndarray] #
Calculates the spherical harmonics for a set of 3D points, whose coordinates are in the
xyz
array, together with their Cartesian derivatives.>>> import numpy as np >>> import sphericart as sc >>> sh = sc.SphericalHarmonics(l_max=8, normalized=False) >>> xyz = np.random.normal(size=(10,3)) >>> sh_values, sh_grads = sh.compute_with_gradients(xyz) >>> sh_grads.shape (10, 3, 81)
- Parameters:
xyz – The Cartesian coordinates of the 3D points, as an array with shape
(n_samples, 3)
.- Returns:
A tuple containing: * an array of shape
(n_samples, (l_max+1)**2)
containing all the spherical harmonics up to degree l_max in lexicographic order. For example, ifl_max = 2
, The last axis will correspond to spherical harmonics with(l, m) = (0, 0), (1, -1), (1, 0), (1, 1), (2, -2), (2, -1), (2, 0), (2, 1), (2, 2)
, in this order. * An array of shape(n_samples, 3, (l_max+1)**2)
containing all the spherical harmonics’ derivatives up to degreel_max
. The last axis is organized in the same way as in the spherical harmonics return array, while the second-to-last axis refers to derivatives in the the x, y, and z directions, respectively.
- compute_with_hessians(xyz: ndarray) Tuple[ndarray, ndarray, ndarray] #
Calculates the spherical harmonics for a set of 3D points, whose coordinates are in the
xyz
array, together with their Cartesian derivatives and second derivatives.>>> import numpy as np >>> import sphericart as sc >>> sh = sc.SphericalHarmonics(l_max=8, normalized=False) >>> xyz = np.random.normal(size=(10,3)) >>> sh_values, sh_grads, sh_hessians = sh.compute_with_hessians(xyz) >>> sh_values.shape (10, 81) >>> sh_grads.shape (10, 3, 81) >>> sh_hessians.shape (10, 3, 3, 81)
- Parameters:
xyz – The Cartesian coordinates of the 3D points, as an array with shape
(n_samples, 3)
.- Returns:
A tuple containing: * an array of shape
(n_samples, (l_max+1)**2)
containing all the spherical harmonics up to degree l_max in lexicographic order. For example, ifl_max = 2
, The last axis will correspond to spherical harmonics with(l, m) = (0, 0), (1, -1), (1, 0), (1, 1), (2, -2), (2, -1), (2, 0), (2, 1), (2, 2)
, in this order. * An array of shape(n_samples, 3, (l_max+1)**2)
containing all the spherical harmonics’ derivatives up to degreel_max
. The last axis is organized in the same way as in the spherical harmonics return array, while the second-to-last axis refers to derivatives in the the x, y, and z directions, respectively. * An array of shape(n_samples, 3, 3, (l_max+1)**2)
containing all the spherical harmonics’ second derivatives up to degreel_max
. The last axis is organized in the same way as in the spherical harmonics return array, while the two intermediate axes represent the Hessian dimensions.