grid.atomgrid module#
Atomic Grid Module.
- class grid.atomgrid.AtomGrid(rgrid: OneDGrid, degrees: ndarray | list | None = [50], *, sizes: ndarray | list = None, center: ndarray = None, rotate: int = 0, method: str = 'lebedev')#
Bases:
GridAtomic grid construction class for integrating three-dimensional functions.
Atomic grid is composed of a radial grid \(\{(r_i, w_i)\}_{i=1}^{N}\) meant to integrate the radius coordinate in spherical coordinates. Further, each radial point is associated with an Angular (Lebedev or Symmetric spherical t-design) grid \(\{(\theta^i_j, \phi^i_j, w_j^i)\}_{j=1}^{M_i}\) that integrates over a sphere (angles in spherical coordinates). The atomic grid points can also be centered at a given location.
- property basis#
ndarray(N, 3): Generate spherical harmonic basis evaluated on atomic grid points.
- property center#
ndarray(3,): Cartesian coordinates of the grid center.
- convert_cartesian_to_spherical(points: ndarray = None, center: ndarray = None)#
Convert a set of points from Cartesian to spherical coordinates.
\[\begin{split}\begin{align} r &= \sqrt{x^2 + y^2 + z^2}\\ \theta &= arc\tan (\frac{y}{x})\\ \phi &= arc\cos(\frac{z}{r}) \end{align}\end{split}\]with the canonical choice \(r=0\), then \(\theta,\phi = 0\). If the points attribute is not specified, then atomic grid points are used and the canonical choice when \(r=0\), is the points \((r=0, \theta_j, \phi_j)\) where \((\theta_j, \phi_j)\) come from the Angular grid with the degree at \(r=0\).
Parameters#
- pointsndarray(N, 3), optional
Three-dimensions Cartesian coordinates of points in atomic units. If None, the atomic grid points will be used.
- centerndarray(3,), optional
Three-dimensional Cartesian coordinates of the center of coordinate system in atomic units. If None, the atomic grid center will be used.
Returns#
- ndarray(N, 3)
Spherical coordinates of points respect to the center denoted as (radius \(r\), azimuthal \(\theta\), polar \(\phi\)).
- property degrees#
ndarray(L,): Return the degree of each angular grid at each radial point.
- classmethod from_preset(atnum: int, preset: str, rgrid: OneDGrid = None, center: ndarray = None, rotate: int = 0, method: str = 'lebedev')#
Construct an atomic grid with pre-defined angular grids for a given atomic number.
Parameters#
- atnumint
The atomic number.
- presetstr, optional
The name of pre-defined grid specifying the radial sectors and their corresponding number of angular grid points. Supported options include our custom presets: ‘coarse’, ‘medium’, ‘fine’, ‘veryfine’, ‘ultrafine’, and ‘insane’, as well as, other standard pre-defined grids including: ‘sg_0’, ‘sg_1’, ‘sg_2’, and ‘sg_3’, and the Ochsenfeld grids: ‘g1’, ‘g2’, ‘g3’, ‘g4’, ‘g5’, ‘g6’, and ‘g7’, with higher number indicating greater accuracy but denser grid. See Notes for more information.
- rgridOneDGrid, optional
The 1D radial grid representing the radius of spherical grids. If None, a default radial grid (PowerRTransform of UniformInteger grid) for the give atnum is constructed.
- centerndarray(3,), optional
Cartesian coordinates of the grid center. If None, the origin is used.
- rotateint, optional
Integer used as a seed for generating random rotation matrices to rotate the angular spherical grids at each radial grid point. If 0, then no rotate is made.
- method: str, optional
Method for constructing the angular grid. Options are “lebedev” (for Lebedev-Laikov) and “spherical” (for symmetric spherical t-design).
Notes#
The standard and Ochsenfeld presets were not designed with symmetric spherical t-design in mind.
The “standard grids” [1] “SG-0” and “SG-1” are designed for large molecules with LDA (GGA) functionals, whereas “SG-2” and “SG-3” are designed for Meta-GGA functionals and B95/Minnesota functionals, respectively.
The Ochsenfeld pruned grids [2] are obtained based on the paper.
References#
- classmethod from_pruned(rgrid: OneDGrid, radius: float, r_sectors: list | ndarray, d_sectors: list | ndarray | None = None, *, s_sectors: list | ndarray = None, center: ndarray = None, rotate: int = 0, method: str = 'lebedev')#
Initialize AtomGrid class that splits radial sections into sectors which specified degrees.
Given a sequence of radial sectors \(\{a_i\}_{i=1}^Q\), a radius number \(R\) and angular degree sectors \(\{L_i \}_{i=1}^{Q+1}\). This assigned the degrees to the following radial points:
\[\begin{split}\begin{align*} &L_1 \text{ when } r < R a_1 \\ &L_2 \text{ when } R a_1 \leq r < R a_2 \vdots \\ &L_{Q+1} \text{ when } R a_{Q} < r. \end{align*}\end{split}\]Parameters#
- rgridOneDGrid
The (one-dimensional) radial grid representing the radius of spherical grids.
- radiusfloat
The atomic radius to be multiplied with r_sectors in atomic units (to make the radial sectors atom specific).
- r_sectorslist or ndarray(S,)
Sequence of boundary radius (in atomic units) specifying sectors of the pruned radial grid. The first sector is
(0, radius*r_sectors[0]), then(radius*r_sectors[0], radius*r_sectors[1]), and so on.- d_sectorslist or ndarray(S + 1, dtype=int) or None
Sequence of angular degrees for each radial sector of r_sectors in the pruned grid. If None, then s_sectors should be given.
- s_sectorslist or ndarray(S + 1, dtype=int) or None, optional, keyword-only
Sequence of angular sizes for each radial sector of r_sectors in the pruned grid. If both d_sectors and s_sectors are given, s_sectors is used.
- centerndarray(3,), optional
Three-dimensional Cartesian coordinates of the grid center in atomic units. If None, the origin is used.
- rotateint, optional
Integer used as a seed for generating random rotation matrices to rotate the angular spherical grids at each radial grid point. If the integer is zero, then no rotate is used.
- method: str, optional
Method for constructing the angular grid. Options are “lebedev” (for Lebedev-Laikov) and “spherical” (for symmetric spherical t-design).
Returns#
- AtomGrid
Generated AtomGrid instance for this special init method.
- get_shell_grid(index: int, r_sq: bool = True)#
Get the spherical integral grid at radial point from specified index.
The spherical integration grid has points scaled with the ith radial point and weights multiplied by the ith weight of the radial grid.
Note that when \(r=0\) then the Cartesian points are all zeros.
Parameters#
- indexint
Index of radial points.
- r_sqbool, default True
If true, multiplies the angular grid weights with r**2.
Returns#
- AngularGrid
Instance of AngularGrid for the given radial index position and weights.
- property indices#
ndarray(M+1,): Indices saved for each spherical shell.
- integrate_angular_coordinates(func_vals: ndarray)#
Integrate the angular coordinates of a sequence of functions.
Given a series of functions \(f_k \in L^2(\mathbb{R}^3)\), this returns the values
\[f_k(r_i) = \int \int f(r_i, \theta, \phi) sin(\phi) d\theta d\phi\]on each radial point \(r_i\) in the atomic grid.
Parameters#
- func_valsndarray(…, N)
The function values evaluated on all \(N\) points on the atomic grid for many types of functions. This can also be one-dimensional.
Returns#
- ndarray(…, M) :
The function \(f_{...}(r_i)\) on each \(M\) radial points.
- interpolate(func_vals: ndarray)#
Return function that interpolates (and its derivatives) from function values.
Any real-valued function \(f(r, \theta, \phi)\) can be decomposed as
\[f(r, \theta, \phi) = \sum_l \sum_{m=-l}^l \sum_i \rho_{ilm}(r) Y_{lm}(\theta, \phi)\]A cubic spline is used to interpolate the radial functions \(\sum_i \rho_{ilm}(r)\). This is then multiplied by the corresponding spherical harmonics at all \((\theta_j, \phi_j)\) angles and summed to obtain the equation above.
Parameters#
- func_valsndarray(N,)
The function values evaluated on all \(N\) points on the atomic grid.
Returns#
- Callable[[ndarray(M, 3), int] -> ndarray(M)]:
Callable function that interpolates the function and its derivative provided. The function takes the following attributes:
- pointsndarray(N, 3)
Cartesian coordinates of \(N\) points to evaluate the splines on.
- derivint, optional
If deriv is zero, then only returns function values. If it is one, then returns the first derivative of the interpolated function with respect to either Cartesian or spherical coordinates. Only higher-order derivatives (deriv=2,3) are supported for the derivatives wrt to radial components.
- deriv_sphericalbool
If True, then returns the derivatives with respect to spherical coordinates \((r, \theta, \phi)\). Default False.
- only_radial_derivbool
If true, then the derivative wrt to radius \(r\) is returned.
This function returns the following.
- ndarray(M,…):
The interpolated function values or its derivatives with respect to Cartesian \((x,y,z)\) or if deriv_spherical then \((r, \theta, \phi)\) or if only_radial_derivs then derivative wrt to \(r\) is only returned.
- property l_max#
int: Largest angular degree L value in angular grids.
- property method#
str: Method used for constructing an angular grid.
- property n_shells#
int: Number of shells in radial points.
- property points#
ndarray(N, 3): 3D Cartesian coordinates of the grid points in atomic units.
- radial_component_splines(func_vals: ndarray)#
Return spline to interpolate radial components wrt to expansion in real spherical harmonics.
For each pt \(r_i\) of the atomic grid with associated angular degree \(l_i\), the function \(f(r_i, \theta, \phi)\) is projected onto the spherical harmonic expansion:
\[f(r_i, \theta, \phi) \approx \sum_{l=0}^{l_i} \sum_{m=-l}^l \rho^{lm}(r_i) Y^m_l(\theta, \phi)\]where \(Y^m_l\) is the real Spherical harmonic of degree \(l\) and order \(m\). The radial components \(\rho^{lm}(r_i)\) are calculated via integration on the \(i\)th Lebedev/angular grid of the atomic grid:
\[\rho^{lm}(r_i) = \int \int f(r_i, \theta, \phi) Y^m_l(\theta, \phi) \sin(\phi) d\theta d\phi,\]and then interpolated using a cubic spline over all radial points of the atomic grid.
Parameters#
- func_valsndarray(N,)
The function values evaluated on all \(N\) points on the atomic grid.
Returns#
- list[scipy.PPoly]
A list of size \((l_{max}/2 + 1)^2\) of CubicSpline object for interpolating the coefficients \(\rho^{lm}(r)\). The input of spline is array of \(N\) points on \([0, \infty)\) and the output is \(\{\rho^{lm}(r_i)\}\). The list starts with degrees \(l=0,\cdots l_{max}\), and the for each degree, the zeroth order spline is first, followed by positive orders then negative.
- property rgrid#
OneDGrid: The radial grid representing the radius of spherical grids.
- property rotate#
int: Integer representing the seed for rotating the angular grid.
- save(filename)#
Save atomic grid attributes as a npz file.
Parameters#
- filename: str
The path/name of the .npz file.
- spherical_average(func_vals: ndarray)#
Return spline of the spherical average of a function.
This function takes a function \(f\) evaluated on the atomic grid points and returns the spherical average of it defined as:
\[f_{avg}(r) := \frac{\int \int f(r, \theta, \phi) \sin(\phi) d\theta d\phi}{4 \pi}.\]The definition is chosen such that \(\int f_{avg}(r) 4\pi r^2 dr\) matches the full integral \(\int \int \int f(x,y,z)dxdydz\).
Parameters#
- func_valsndarray(N,)
The function values evaluated on all \(N\) points on the atomic grid.
Returns#
- CubicSpline:
Cubic spline with input r in the positive real axis and output \(f_{avg}(r)\).