grid.angular module#

Angular grid module for constructing integration grids on the unit sphere.

The Lebedev grid points used were obtained from the theochem/HORTON package. Their calculations were based on F77 translation by Dr. Christoph van Wuellen and these are the comments from that translation.

This subroutine is part of a set of subroutines that generate Lebedev grids [1] [2] [3] [4] [5] [6] for integration on a sphere. The original C-code [1] was kindly provided by Dr. Dmitri N. Laikov and translated into Fortran by Dr. Christoph van Wuellen. This subroutine was translated from C to fortran77 by hand. Users of this code are asked to include reference [1] in their publications, and in the user- and programmers-manuals describing their codes. This code was distributed through CCL (http://www.ccl.net/).

The symmetric spherical t-design were obtained from reference [7].

References#

The following references are for the Lebedev grid points.

The following references are for the symmetric spherical t-design points:

class grid.angular.AngularGrid(degree: int | None = 50, *, size: int | None = None, cache: bool = True, method: str = 'lebedev')#

Bases: Grid

Angular grid for integrating functions on the unit sphere.

This class numerically evaluates the surface integral of a function \(f: S^2 \rightarrow \mathbb{R}\) over the unit-sphere:

\[\int_{S^2} f d\Omega = \int_0^{2\pi} \int_0^\pi f(\theta, \phi) \sin(\phi) d\theta d\phi \approx \sum_{i}^{N} w_i f(\phi_i, \theta_i),\]

where \(S^2\) is the unit-sphere parameterized by the quadrature points \(\theta_i \in [0, 2\pi]\) and \(\phi_i \in [0, \pi)\), and \(w_{i}\) are the weights of the \(N\) quadrature points. The \(4\pi\) normalization factor present in the original quadrature scheme is included in the integration weights.

Two types of angular grids are supported: Lebedev-Laikov grid and symmetric spherical t-design. Specifically, for spherical t-design, the weights are constant value of \(4 \pi / N\), where \(N\) is the number of points in the grid. The weights are chosen so that the spherical harmonics are normalized.

static convert_angular_sizes_to_degrees(sizes: ndarray, method: str)#

Convert given Lebedev/Spherical design grid sizes to degrees.

Parameters#

sizesndarray[int]

Sequence of angular grid sizes (e.g., number of points for each atomic shell).

method: str

Method for constructing the angular grid. Options are “lebedev” (for Lebedev-Laikov) and “spherical” (for symmetric spherical t-design).

Returns#

ndarray[int]

Sequence of the corresponding angular degree of the angular grid corresponding to its size.

property degree#

int: The degree of spherical harmonics that this angular grid can integrate exactly.

property method#

str: Method used for constructing the angular grid.