openquad.S2#
- class openquad.S2(method_specs, polar_sampling='cos')#
Two-dimensional unit sphere.
This class provides an interface to construct quadrature methods for the integral
\[\int_{0}^{\pi} \int_{0}^{2\pi} f(\theta, \phi) \, \sin(\theta) \, d\phi \, d\theta\]The innermost integral is evaluated first.
The default coordinates of this geometry are the spherical polar angles [3] in the order \((\theta, \phi)\).
- Parameters:
- method_specslist of tuples
Specification of the quadrature methods. Each entry of the list contains a tuple corresponding to a quadrature method and its options.
The first entry of each tuple needs to be a string denoting the quadrature method. See this list for an overview of available methods. The other entries of the tuple need to be keyword arguments for the individual methods. At least,
size
needs to be given, ordegree
if applicable. Other options are tabulated here.Individual quadrature methods are applied to the default coordinates in the order mentioned above. See below, for examples.
- polar_sampling{‘angle’, ‘cos’}, optional
Controls if the second integral is sampled along \(\theta\) or \(\cos\theta\).
References
- [1]: Wikipedia page:
- [2]: Wikipedia page:
- [3]: Wikipedia page:
Examples
A quadrature composed of one-dimensional quadratures. Gauss-Legendre quadrature is applied to the first coordinate, i.e., \(\theta\). The composite trapezoidal rule is applied to the second coordinate, \(\phi\):
>>> openquad.S2([ ... ('GaussLegendre', dict(degree=6)), ... ('trapezoid', dict(size=7)), ... ])
A quadrature composed of a two-dimensional \(\mathrm{S}^2\) Gauss quadrature:
>>> openquad.S2([ ... ('S2-Gauss-LebedevLaikov', dict(degree=5)), ... ])
- Attributes:
- dimint
Dimension of the geometry.
- sizeint
Total number of sampling points.
- shapetuple of int
The shape of points,
(dim, size)
.- pointsndarray
Sampling points in the default coordinates. The array has shape
(2, size)
.- weightsndarray
Quadrature weights. The array has shape
(size,)
.- anglesndarray
Equivalent to points
- xyzndarray
Sampling points in cartesian coordinates, \((x, y, z)\). This array has shape
(3, size)
.
- integrate(f, f_kwargs={}, axis=-1)#
Integrate the given data or callable.
- Parameters:
- fcallable or array_like
Representation of the integrand. Either an array-like or an callable.
- f_kwargsdict, optional
Additional keyword arguments to pass to f, if f is a callable.
- axissequence of ints, optional
Axes along which to perform the integration. The number of axes provided must match the number of integration methods hold by the quadrature object. By default, start with the last axis of f.
Caution
Experimental feature
- Returns:
- resultscalar or ndarray
The result of the integration. If the dimension of f equals the number of integration methods, the result is a scalar.
- savetxt(filename, weights=True)#
Save sampling points and weights to a text file.
This function wraps
numpy.savetxt()
.- Parameters:
- filenamefilename, file handle or pathlib.Path
Name of the text file.
- weightslogical, optional
If
False
, only save sampling points and not the quadrature weights.
- angles#
- property dim#
int: Dimension of this geometry.
- points#
- weights#
- xyz#