grid.basegrid module
Construct basic grid data structure.
- class Grid(points, weights)[source]
Bases:
objectBasic Grid class for grid information storage.
- __init__(points, weights)[source]
Construct Grid instance.
- Parameters:
points (np.ndarray(N,) or np.ndarray(N, M)) – An array with positions of the grid points.
weights (np.ndarray(N,)) – An array of weights associated with each point on the grid.
- get_localgrid(center, radius)[source]
Create a grid containing 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:
- integrate(*value_arrays)[source]
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)[source]
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 function, 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\) is evaluated at 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 at the mth 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)[source]
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 weights
the weights of each grid point.
- Type:
np.ndarray(N,)
- class LocalGrid(points, weights, center, indices=None)[source]
Bases:
GridLocal portion of a grid, containing all points within a sphere.
- __init__(points, weights, center, indices=None)[source]
Initialize a local grid.
- Parameters:
points (np.ndarray(N,) or np.ndarray(N,M)) – Cartesian coordinates of \(N\) grid points in 1D or M-D space.
weights (np.ndarray(N)) – Integration weight of \(N\) grid points
center (float or np.ndarray(M,)) – Cartesian coordinates of the center of the local grid in 3D space.
indices (np.ndarray(N,), optional) – Indices of \(N\) grid points and weights in the parent grid.
- property center
Cartesian coordinates of the center of the local grid.
- Type:
np.ndarray(3,)
- get_localgrid(center, radius)[source]
Create a grid containing 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:
- property indices
Indices of grid points and weights in the parent grid.
- Type:
np.ndarray(N,)
- integrate(*value_arrays)[source]
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)[source]
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 function, 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\) is evaluated at 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 at the mth 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)[source]
Save the points, indices 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 weights
the weights of each grid point.
- Type:
np.ndarray(N,)
- class OneDGrid(points, weights, domain=None)[source]
Bases:
GridOne-Dimensional Grid.
- __init__(points, weights, domain=None)[source]
Construct grid.
- Parameters:
points (np.ndarray(N,)) – A 1-D array of coordinates of \(N\) points in one-dimension.
weights (np.ndarray(N,)) – A 1-D array of integration weights of \(N\) points in one-dimension.
domain (tuple(float, float), optional) – Lower and upper bounds for which the grid can carry out numerical integration. This does not always coincide with the positions of the first and last grid point. For example, in case of the Gauss-Chebyshev quadrature the domain is [-1,1] but all grid points lie in (-1, 1).
- property domain
the range of grid points.
- Type:
(float, float)
- get_localgrid(center, radius)[source]
Create a grid containing 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:
- integrate(*value_arrays)[source]
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)[source]
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 function, 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\) is evaluated at 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 at the mth 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)[source]
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 weights
the weights of each grid point.
- Type:
np.ndarray(N,)