Molecular grid class for integration of three-dimensional functions.
Molecular grid is defined here to be a weighted average of \(M\) atomic grids
(see AtomGrid). This is defined by a atom in molecule weights (or nuclear weight functions)
\(w_n(r)\) for each center n such that \(\sum^M_n w_n(r) = 1\) for all points
\(r\in\mathbb{R}^3\). [1]
atnums (np.ndarray(M, 3)) – Atomic number of \(M\) atoms in molecule.
atgrids (list[AtomGrid]) – List of atomic grids of size \(M\) for each atom in molecule.
aim_weights (Callable or np.ndarray(sum^M_n N_n,)) – Atoms in molecule weights \({ {w_n(r_k)}_k^{N_i}}_n^{M}\), where
\(N_i\) is the number of points in the ith atomic grid.
store (bool) – If true, then the atomic grids atgrids are stored as attribute atgrids.
atnums (np.ndarray(M,)) – Array of atomic numbers.
atcoords (np.ndarray(M, 3)) – Atomic coordinates of atoms.
preset ((str, list[str], dict[int: str])) – Preset grid accuracy scheme. If string is provided, then preset is used
for all atoms, either it is specified by a list, or a dictionary whose keys
are from atnums. These predefined grid specify the radial sectors and
their corresponding number of Lebedev grid points. Supported preset options include:
‘coarse’, ‘medium’, ‘fine’, ‘veryfine’, ‘ultrafine’, and ‘insane’.
Other options include the “standard grids”:
‘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.
rgrid ((OneDGrid, list[OneDGrid], dict[int: OneDGrid], None), optional) – One dimensional radial grid. If of type OneDGrid then this radial grid is used for
all atoms. If a list is provided, then ith grid correspond to the ith atom. If
dictionary is provided, then the keys correspond to the atnums[i] attribute.
If None, a default radial grid (PowerRTransform of UniformInteger grid) is constructed
based on the given atomic numbers.
aim_weights (Callable or np.ndarray(K,), optional) – Atoms in molecule weights. If None, then aim_weights is Becke weights with order=3.
rotate (bool or int, optional) – Flag to set auto rotation for atomic grid, if given int, the number
will be used as a seed to generate rantom matrix.
store (bool, optional) – Store atomic grid as a class attribute.
Notes
The standard and Ochsenfeld presets were not designed with symmetric spherical t-design
in mind.
The “standard grids” [2] “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 [3] are obtained based on the paper.
Initialize a MolGrid instance with pruned method from AtomGrid.
Parameters:
atnums (ndarray(M,)) – List of atomic numbers for each atom.
atcoords (np.ndarray(M, 3)) – Three-dimensional Cartesian coordinates for each atoms in atomic units.
radius (float, List[float]) – The atomic radius to be multiplied with r_sectors in atomic units (to make the
radial sectors atom specific). If float, then the same atomic radius is used for all
atoms, otherwise a list with \(M\) elements is used, where \(M\) is the number
of atoms in the molecule. If list, then the ith element is used for the ith atom.
r_sectors (list of List[float]) – List of sequences of the boundary radius (in atomic units) specifying sectors of
the pruned radial grid of \(M\) atoms. For the first atom, the first
sector is (0,radius*r_sectors[0][0]), then (radius*r_sectors[0][0],radius*r_sectors[0][1]), and so on. See AtomGrid.from_pruned for more information.
d_sectors (int or list of List[int], optional) – List of sequences of the angular degrees for radial sectors of \(M\) atoms.
If a number is given, then the same number of degrees is used for all sectors of
all atoms.
s_sectors (int or list of List[int] or None, optional, keyword-only) – List of sequences of angular sizes for each radial sector of of \(M\) atoms.
If both d_sectors and s_sectors are given, s_sectors is used.
rgrid (OneDGrid or List[OneDGrid] or Dict[int: OneDGrid], optional) – One dimensional grid for the radial component. If a list is provided,then ith
grid correspond to the ith atom. If dictionary is provided, then the keys are
correspond to the atnums[i] attribute. If None, then using atomic numbers it will
generate a default radial grid (PowerRTransform of UniformInteger grid).
aim_weights (Callable or np.ndarray(sum^M_n N_n,), optional) – Atoms in molecule/nuclear weights \({ {w_n(r_k)}_k^{N_i}}_n^{M}\), where
\(N_i\) is the number of points in the ith atomic grid. If None, then aim_weights
is Becke weights with order=3.
rotate (bool or int , optional) – Flag to set auto rotation for atomic grid, if given int, the number
will be used as a seed to generate random matrix.
store (bool, optional) – Flag to store each original atomic grid information.
Initialize a MolGrid instance with Horton Style input.
Parameters:
atnums (np.ndarray(M, 3)) – Atomic number of \(M\) atoms in molecule.
atcoords (np.ndarray(N, 3)) – Cartesian coordinates for each atoms
size (int) – Number of points on each shell of angular grid.
rgrid (OneDGrid or None, optional) – One-dimensional grid to construct the atomic grid. If none, then
default radial grid is generated based on atomic numbers.
aim_weights (Callable or np.ndarray(K,), optional) – Atoms in molecule weights. If None, then aim_weights is Becke weights with order=3.
rotate (bool or int , optional) – Flag to set auto rotation for atomic grid, if given int, the number
will be used as a seed to generate rantom matrix.
store (bool, optional) – Flag to store each original atomic grid information
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.
Return function that interpolates (and its derivatives) from function values.
Consider a real-valued function \(f(r, \theta, \phi) = \sum_A f_A(r, \theta, \phi)\)
written as a sum of atomic centers \(f_A(r, \theta, \phi)\). Each of these
functions can be further decomposed based on the atom grid centered at \(A\):
A cubic spline is used to interpolate the radial functions \(\sum_i \rho^A_{ilm}(r)\)
based on the atomic grid centered at \(A\).
This is then multipled by the corresponding spherical harmonics at all
\((\theta_j, \phi_j)\) angles and summed to obtain approximation to \(f_A\). This
is then further summed over all centers to get \(f\).
Parameters:
func_vals (ndarray(sum_i N_i,)) – The function values evaluated on all \(N_i\) points on the \(i\)th atomic grid.
Returns:
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.
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) ,…].