grid.angular module

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

The Lebedev grid points here were obtained and calculated 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-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.

1(1,2)

V.I. Lebedev, and D.N. Laikov “A quadrature formula for the sphere of the 131st algebraic order of accuracy” Doklady Mathematics, Vol. 59, No. 3, 1999, pp. 477-481.

2

V.I. Lebedev “A quadrature formula for the sphere of 59th algebraic order of accuracy” Russian Acad. Sci. Dokl. Math., Vol. 50, 1995, pp. 283-286.

3

V.I. Lebedev, and A.L. Skorokhodov “Quadrature formulas of orders 41, 47, and 53 for the sphere” Russian Acad. Sci. Dokl. Math., Vol. 45, 1992, pp. 587-592.

4

V.I. Lebedev “Spherical quadrature formulas exact to orders 25-29” Siberian Mathematical Journal, Vol. 18, 1977, pp. 99-107.

5

V.I. Lebedev “Quadratures on a sphere” Computational Mathematics and Mathematical Physics, Vol. 16, 1976, pp. 10-24.

6

V.I. Lebedev “Values of the nodes and weights of ninth to seventeenth order Gauss-Markov quadrature formulae invariant under the octahedron group with inversion” Computational Mathematics and Mathematical Physics, Vol. 15, 1975, pp. 44-51.

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

7

R. S. Womersley, Efficient Spherical Designs with Good Geometric Properties. In: Dick J., Kuo F., Wozniakowski H. (eds) Contemporary Computational Mathematics - A Celebration of the 80th Birthday of Ian Sloan. Springer (2018) pp. 1243-1285 https://doi.org/10.1007/978-3-319-72456-0_57

class AngularGrid(points=None, weights=None, *, degree=None, size=None, cache=True, use_spherical=False)

Bases: Grid

Angular grid for integrating functions on the unit sphere.

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

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

where \(S^2\) is the unit-sphere, \(w^{ang}_i\) are the quadrature weights that includes \(4\pi\).

Two types of angular grids are provided: Lebedev-Laikov grid (default) or the 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 integrate to \(4\pi\):

\[\int_0^\pi \int_0^{2\pi} Y^l_m(\theta, \phi) \sin(\theta) d\theta d\phi = 4\pi.\]
__init__(points=None, weights=None, *, degree=None, size=None, cache=True, use_spherical=False)

Generate angular grid for a given degree and/or size.

Three choices to construct this grid:

  1. Add your own points and weights.

  2. Specify maximum degree of spherical harmonics that the angular grid can integrate accurately on a unit sphere.

  3. The number of points (size) wanted in the grid. The size corresponds to a maximum degree.

Parameters
  • points (ndarray(N, 3), optional) – The Cartesian coordinates of integral grid points. The attribute weights should also be provided.

  • weights (ndarray(N,), optional) – The weights of each point on the integral quadrature grid. The attribute points should also be applied.

  • degree (int, optional) – Maximum angular degree \(l\) of spherical harmonics that the Lebedev grid can integrate accurately. If the Lebedev grid corresponding to the given angular degree is not supported, the next largest degree is used.

  • size (int, optional) – Number of Lebedev grid points. If the Lebedev grid corresponding to the given size is not supported, the next largest size is used. If both degree and size are given, degree is used for constructing the grid.

  • cache (bool, optional) – If true, then store the points and weights of the AngularGrid in cache to avoid duplicate grids that have the same degree. Default to True.

  • use_spherical (bool, optional) – If true, then it loads/uses the symmetric spherical t-design, rather than Lebedev-Laikov grid.

Returns

An angular grid with points and weights (including \(4\pi\)) on a unit sphere.

Return type

AngularGrid

Examples

>>> # Initialize with degree
>>> AngularGrid(degree=3)
>>>
>>> # Initialize with size
>>> AngularGrid(size=6)
>>>
>>> # Initialize with points and weights
>>> pts = np.array([[0, 1, 0], [1, 0, 1], [0, 1, 1]], dtype=float)
>>> wts = np.array([0.5, 0.4, 0.3])
>>> AngularGrid(pts, wts)

Notes

  • Sometimes the weights for Lebedev-Laikov grids can be negative. Choosing degrees that have positive weights can mitigate round-off errors. Degrees equal to 13, 25 or 27 have negative weights.

static convert_angular_sizes_to_degrees(sizes, use_spherical=False)

Convert given Lebedev/Spherical design grid sizes to degrees.

Parameters
  • sizes (ndarray[int]) – Sequence of angular grid sizes (e.g., number of points for each atomic shell).

  • use_spherical (bool, optional) – If true, loads the symmetric spherical t-design grid rather than the Lebedev-Laikov grid.

Returns

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

Return type

ndarray[int]

property degree

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

Type

int

get_localgrid(center, radius)

Create a grid contain points within the given radius of center.

Parameters
  • center (float or np.array(M,)) – Cartesian coordinates of the center of the local grid.

  • radius (float) – Radius of sphere around the center. When equal to np.inf, the local grid coincides with the whole grid, which can be useful for debugging.

Returns

Instance of LocalGrid.

Return type

LocalGrid

integrate(*value_arrays)

Integrate over the whole grid for given multiple value arrays.

Product of all value_arrays will be computed element-wise then integrated on the grid with its weights:

\[\int w(x) \prod_i f_i(x) dx.\]
Parameters

*value_arrays (np.ndarray(N, )) – One or multiple value array to integrate.

Returns

The calculated integral over given integrand or function

Return type

float

moments(orders, centers, func_vals, type_mom='cartesian', return_orders=False)

Compute the multipole moment integral of a function over centers.

The Cartesian type moments are:

\[m_{n_x, n_y, n_z} = \int (x - X_c)^{n_x} (y - Y_c)^{n_y} (z - Z_c)^{n_z} f(r) dr,\]

where \(\textbf{R}_c = (X_c, Y_c, Z_c)\) is the center of the moment, \(\f(r)\) is the density, and \((n_x, n_y, n_z)\) are the Cartesian orders.

The spherical/pure moments with \((l, m)\) parameter are:

\[m_{lm} = \int | \textbf{r} - \textbf{R}_c|^l S_l^m(\theta, \phi) f(\textbf{r}) d\textbf{r},\]

where \(S_l^m\) is a regular, real solid harmonic.

The radial moments with \(n\) parameter are:

\[m_n = \int | \textbf{r} - \textbf{R}_c|^{n} f(\textbf{r}) d\textbf{r}\]

The radial combined with spherical/pure moments \((n, l, m)\) are:

\[m_{nlm} = \int | \textbf{r} - \textbf{R}_c|^{n+1} S_l^m(\theta, \phi) f(\textbf{r}) d\textbf{r}\]
Parameters
  • orders (int) – Generates all orders with Horton order depending on the type of the multipole moment type_mom.

  • centers (ndarray(M, 3)) – The centers \(\textbf{R}_c\) of the moments to compute from.

  • func_vals (ndarray(N,)) – The function \(f\) values evaluated on all \(N\) points on the integration grid.

  • type_mom (str) – The type of multipole moments: “cartesian”, “pure”, “radial” and “pure-radial”.

  • return_orders (bool) – If true, it will also return a list of size \(L\) of the orders corresponding to each integral/row of the output.

Returns

Computes the moment integral of the function on the m`th center for all orders. If `return_orders is true, then this also returns a list that describes what each row/order is, e.g. for Cartesian, [(0, 0, 0), (1, 0, 0) ,…].

Return type

ndarray(L, M), or (ndarray(L, M), list)

property points

Positions of the grid points.

Type

np.ndarray(N,) or np.ndarray(N, M)

save(filename)

Save the points and weights as a npz file.

Parameters

filename (str) – The path/name of the .npz file.

property size

the total number of points on the grid.

Type

int

property use_spherical

If true, then symmetric spherical t-design was used else Lebedev-Laikov grid.

Type

bool

property weights

the weights of each grid point.

Type

np.ndarray(N,)