horton_part.alisa module
Module for Linear Iterative Stockholder Analysis (L-ISA) partitioning scheme.
- class LinearISAWPart(coordinates, numbers, pseudo_numbers, grid, moldens, spindens=None, lmax=3, logger=None, threshold=1e-06, maxiter=500, inner_threshold=1e-08, radius_cutoff=inf, solver='cvxopt', solver_options=None, basis_func='gauss', basis_type='analytic', grid_type=1, **kwargs)
Bases:
GaussianISAWPart
Implements the Linear Iterative Stockholder Analysis (L-ISA) partitioning scheme.
This class extends
GaussianISAWPart
and specializes in performing electron density partitioning in molecules using various L-ISA schemes. L-ISA is a method for dividing the electron density of a molecule into atomic contributions. This class offers a variety of schemes for this partitioning, both at local [1] and global optimization levels.Optimization Problem Schemes
Convex optimization (`solver`=”cvxopt”)
See also
horton_part.lisa_g
Global Linear approximation of ISA.
References
- __init__(coordinates, numbers, pseudo_numbers, grid, moldens, spindens=None, lmax=3, logger=None, threshold=1e-06, maxiter=500, inner_threshold=1e-08, radius_cutoff=inf, solver='cvxopt', solver_options=None, basis_func='gauss', basis_type='analytic', grid_type=1, **kwargs)
LISA initial function.
Optional arguments: (that are not defined in
GaussianISAWPart
)- Parameters:
basis_func (string, optional) – The type of basis functions, and Gaussian is the default.
- property bs_helper
A basis function helper.
- builtin_solvers = {'cdiis': <function solver_cdiis>, 'cvxopt': <function solver_cvxopt>, 'diis': <function solver_diis>, 'm-newton': <function solver_m_newton>, 'newton': <function solver_newton>, 'quasi-newton': <function solver_quasi_newton>, 'sc': <function solver_sc>, 'sc-1-iter': <function solver_sc_1_iter>, 'sc-plus-convex': <function solver_sc_plus_cvxopt>, 'trust-region': <function solver_trust_region>}
- property cache
Cache.
- calc_radial_distances()
Calculate radial distance w.r.t coordinates.
- clear()
Discard all cached results, e.g. because wfn changed
- compute_change(propars1, propars2)
Compute the difference between an old and a new proatoms
- compute_pseudo_population(index)
Compute pseudo population
- property coordinates
Center/Atomic coordinates.
- Type:
ndarray(M, 3)
- property density_cutoff
Get the density cutoff value.
Density values below this cutoff are considered to be invalid.
- Returns:
The cutoff value for density.
- Return type:
float
- do_all()
Computes all properties and return a list of their keys.
- do_charges()
Compute atomic charges.
- do_density_decomposition()
Compute density decomposition.
- do_moments()
Compute atomic multiple moments.
Calculates various types of multipoles, including Cartesian, Spherical, and Radial moments. The order of the moments is determined by the lmax parameter.
- do_partitioning()
Do partitioning
- do_populations()
Compute atomic populations.
- do_prosplines()
Do pro-atom splines
- do_spin_charges()
Compute atomic spin charges.
- eval_proatom(index, output, grid)
Evaluate function on the molecular grid.
The size of the local grid is specified by the radius of the sphere where the local grid is considered. For example, when the radius is np.inf, the grid corresponds to the whole molecular grid.
- Parameters:
index (int) – The index of an atom in the molecule.
output (1D np.array) – The size of output should be the same as the size of the local grid.
grid (2D np.array) – The molecular grid.
- eval_spline(index, spline, output, grid, label='noname')
Evaluate a given spline at radial distances from a specified atom center and store the results in the output array.
This method calculates the radial distances from the specified atom center to each point in the provided grid. It then evaluates the provided spline function at these distances, storing the results in the given output array.
- Parameters:
index (int) – The index of the atom whose center is used for calculating radial distances.
spline (callable) – The spline function to be evaluated. This should be a function that takes an array of radial distances and returns the corresponding spline values.
output (1D ndarray) – The array where the evaluated spline values will be stored. This array is modified in-place.
grid (Grid) – An object representing the grid points. It should have an attribute points which is an array of grid point coordinates.
label (str, optional) – A label for identification purposes, defaults to “noname”.
Notes
The method computes the Euclidean norm (radial distance) from the atom center, specified by index, to each point in the grid. The spline function is then evaluated at these distances. The results are stored directly in the output array, overwriting any existing data.
- fix_proatom_rho(index, rho, deriv)
Check if the radial density for the proatom is correct and fix as needed.
- Parameters:
index (
int
) – The atom for which this proatom rho is created.rho (1D np.ndarray) – The radial density
deriv (
int
) – the derivative of the radial density or None.
- get_grid(index=None)
Return an integration grid
Optional arguments:
- index
The index of the atom. If not given, a grid for the entire system is returned. If self.local is False, a full system grid is always returned.
- get_moldens(index=None, output=None)
Retrieves the molecular electron density (moldens) on the atomic grid.
This method converts the molecular electron density to the atomic grid specified by the index. If an output array is provided, the result is stored in that array.
- Parameters:
index (int or None, optional) – The index of the atom for which the electron density is required. If None, electron density for all atoms is considered. Default is None.
output (np.ndarray or None, optional) – An optional array to store the resulting electron density. If provided, the result is saved in this array. Default is None.
- Returns:
The molecular electron density on the atomic grid.
- Return type:
np.ndarray
- get_proatom_rho(iatom, propars=None, **kwargs)
Get pro-atom density for atom iatom.
If propars is None, the cache values are used; otherwise, the propars are used.
- Parameters:
iatom (int) – The index of atom iatom.
propars (1D np.array) – The pro-atom parameters.
- get_proatom_spline(index, *args, **kwargs)
Create and return a spline representation of the radial density for a given atomic index.
This method first retrieves the radial density and its derivatives for the specified atomic index. It then ensures the correctness of these values and constructs a spline representation based on the radial grid points.
- Parameters:
index (int) – The index of the atom for which the radial density spline is to be calculated.
*args – Variable length argument list, used for passing non-keyworded arguments.
**kwargs – Arbitrary keyword arguments, used for passing additional data.
- Returns:
A spline representation of the radial density. If derivatives are available, a CubicHermiteSpline is returned. Otherwise, a CubicSpline is used.
- Return type:
CubicSpline or CubicHermiteSpline
Notes
The method internally calls get_proatom_rho to obtain the radial density (rho) and its derivatives (deriv), and fix_proatom_rho to validate and potentially correct these values. It also uses get_rgrid to acquire the radial grid points (rgrid.points). The spline is constructed with these grid points and density values, with the type of spline depending on the availability of derivative information.
- get_rgrid(index)
Load radial grid.
- get_spindens(index=None, output=None)
Retrieves the spin density (spindens) on the atomic grid.
This method converts the spin density to the atomic grid specified by the index. If an output array is provided, the result is stored in that array.
- Parameters:
index (int or None, optional) – The index of the atom for which the spin density is required. If None, spin density for all atoms is considered. Default is None.
output (np.ndarray or None, optional) – An optional array to store the resulting spin density. If provided, the result is saved in this array. Default is None.
- Returns:
The spin density on the atomic grid.
- Return type:
np.ndarray
- get_wcor(index)
Load correction of weights.
- property grid
Molecular grid.
- property grid_type
Get the type of grids used for partitioning density.
- Returns:
The type of grids used in the partitioning process.
- Return type:
str
- property lmax
The maximum angular momentum index for moment calculations.
- property local
Whether local grids are included.
- name = 'lisa'
- property natom
The number of atoms in the molecule.
- property negative_cutoff
Get the negative cutoff value.
Values less than this threshold are treated as negative in computations.
- Returns:
The negative cutoff value.
- Return type:
float
- property nelec
The number of electrons in the molecule.
- property numbers
Atomic numbers
- property on_molgrid
Check whether quantities are computed on molecular grids.
These grids are used for evaluating various properties, including:
AIM (Atoms-in-Molecule) weight functions: \(w_a(\mathbf{r})\).
Pro-atom density: :math:`
- ho_a^0(mathbf{r})`.
Pro-molecule density: :math:`
- ho^0(mathbf{r})`.
Mean-square deviation during computations.
- bool
True if quantities are computed on molecular grids, False otherwise.
- property only_use_molgrid
Check whether intermediate values are computed exclusively using molecular grids.
When set to True, all quantities are computed solely on the molecular grid.
- Returns:
True if intermediate values are computed only using the molecular grid, False otherwise.
- Return type:
bool
- property population_cutoff
Get the population cutoff criterion.
This represents the allowed difference between the sum of proatom parameters and the reference population for determining accuracy of methods.
- Returns:
The cutoff value for population differences.
- Return type:
float
- property pseudo_numbers
Atomic charges.
- property radial_distances
Get the radial distances of points from the atomic coordinates.
The radial distances are calculated as the L2 norm (Euclidean distance) of the points relative to the atomic coordinates.
Notes
Accessing this property triggers the calculation of radial distances via the calc_radial_distances method.
- Returns:
A list containing the radial distances of points for each atom.
- Return type:
list
- setup_grids()
Setup grids used in partitioning.
# 1. atom_grids + mol_grid, use atoms for everything and mol_grid is only used for weights calculation # 2. atom_grids + mol_grid, use mol_grid for everything but atom_grids are used applied contratins. # 3. mol_grid, use mol_grid for everything, like in gLISA+.
- to_atomic_grid(index, data)
Load atomic contribution of molecular properties.
- update_at_weights(force_on_molgrid=False)
See
Part.update_at_weights
.
- update_pro(index, proatdens, promoldens, force_on_molgrid=False)
Update propars.
- Parameters:
index (int) – The index of the atom for which the pro-atom and pro-molecule densities are to be updated.
proatdens (1D np.ndarray) – The array representing the pro-atom density. This array is updated with the new density values for the specified atom.
promoldens (1D np.ndarray) – The array representing the pro-molecule density. This array accumulates the density contributions from each atom, including the one specified by index.
- variables_stored_in_cache()
The properties stored in cache obj.
- setup_bs_helper(part)
Setup basis function helper.
- solver_cdiis(bs_funcs, rho, propars, points, weights, threshold, logger, density_cutoff, negative_cutoff, population_cutoff, check_mono=True, **cdiis_options)
Optimize parameters for proatom density functions using CDIIS algorithm.
- Parameters:
bs_funcs (2D np.ndarray) – Basis functions array with shape (M, N), where ‘M’ is the number of basis functions and ‘N’ is the number of grid points.
rho (1D np.ndarray) – Spherically-averaged atomic density as a function of radial distance, with shape (N,).
propars (1D np.ndarray) – Pro-atom parameters with shape (M). ‘M’ is the number of basis functions.
points (1D np.ndarray) – Radial coordinates of grid points, with shape (N,).
weights (1D np.ndarray) – Weights for integration, including the angular part (4πr²), with shape (N,).
threshold (float, optional) – default value : 1e-08 tolerence parameter for convergence test on residual (commutator)
logger (logging.Logger or None) – The log object.
density_cutoff (float) – Density values below this cutoff are considered invalid.
negative_cutoff (float) – The value less than negative_cutoff is treated as negative.
population_cutoff (float) – The criteria for the difference between the sum of propars and reference population.
cdiis_options (dict) – Other options for the solver.
check_mono (bool, opotionl) – Check if the density is monotonically decaying.
- Returns:
conv – if convergence, True, else, False
- Return type:
boolean
- solver_cvxopt(bs_funcs, rho, propars, points, weights, threshold, logger, density_cutoff, negative_cutoff, population_cutoff, allow_neg_params=False, **cvxopt_options)
Optimize parameters for pro-atom density functions using convex optimization method.
- Parameters:
bs_funcs (2D np.ndarray) – Basis functions array with shape (M, N), where ‘M’ is the number of basis functions and ‘N’ is the number of grid points.
rho (1D np.ndarray) – Spherically-averaged atomic density as a function of radial distance, with shape (N,).
propars (1D np.ndarray) – Pro-atom parameters with shape (M). ‘M’ is the number of basis functions.
points (1D np.ndarray) – Radial coordinates of grid points, with shape (N,).
weights (1D np.ndarray) – Weights for integration, including the angular part (4πr²), with shape (N,).
threshold (float) – Convergence threshold for the iterative process.
logger (logging.Logger or None) – The log object.
density_cutoff (float) – Density values below this cutoff are considered invalid.
negative_cutoff (float) – The value less than negative_cutoff is treated as negative.
population_cutoff (float) – The criteria for the difference between the sum of propars and reference population.
allow_neg_params (bool (optional)) – Whether negative parameters are allowed. The default is False.
cvxopt_options (dict) – Other options for the cvxopt solver.
- Returns:
Optimized proatom parameters.
- Return type:
1D np.ndarray
- Raises:
RuntimeError – If the inner iteration does not converge.
- solver_diis(bs_funcs, rho, propars, points, weights, threshold, logger, density_cutoff, negative_cutoff, population_cutoff, check_mono=True, **diis_options)
Optimize parameters for proatom density functions using direct inversion in an iterative space (DIIS).
- Parameters:
bs_funcs (2D np.ndarray) – Basis functions array with shape (M, N), where ‘M’ is the number of basis functions and ‘N’ is the number of grid points.
rho (1D np.ndarray) – Spherically-averaged atomic density as a function of radial distance, with shape (N,).
propars (1D np.ndarray) – Pro-atom parameters with shape (M). ‘M’ is the number of basis functions.
points (1D np.ndarray) – Radial coordinates of grid points, with shape (N,).
weights (1D np.ndarray) – Weights for integration, including the angular part (4πr²), with shape (N,).
threshold (float) – Convergence threshold for the iterative process.
logger (logging.Logger or None) – The log object.
density_cutoff (float) – Density values below this cutoff are considered invalid.
negative_cutoff (float) – The value less than negative_cutoff is treated as negative.
population_cutoff (float) – The criteria for the difference between the sum of propars and reference population.
check_mono (bool, opotionl) – Check if the density is monotonically decaying.
diis_options (dict) – The other options for the solver.
- Returns:
Optimized proatom parameters.
- Return type:
1D np.ndarray
- Raises:
RuntimeError – If the inner iteration does not converge.
- solver_m_newton(bs_funcs, rho, propars, points, weights, threshold, logger, density_cutoff, negative_cutoff, population_cutoff, **newton_options)
Optimize parameters for pro-atom density functions using Newton method
- Parameters:
bs_funcs (2D np.ndarray) – Basis functions array with shape (M, N), where ‘M’ is the number of basis functions and ‘N’ is the number of grid points.
rho (1D np.ndarray) – Spherically-averaged atomic density as a function of radial distance, with shape (N,).
propars (1D np.ndarray) – Pro-atom parameters with shape (M). ‘M’ is the number of basis functions.
points (1D np.ndarray) – Radial coordinates of grid points, with shape (N,).
weights (1D np.ndarray) – Weights for integration, including the angular part (4πr²), with shape (N,).
threshold (float) – Convergence threshold for the iterative process.
logger (logging.Logger or None) – The log object.
density_cutoff (float) – Density values below this cutoff are considered invalid.
negative_cutoff (float) – The value less than negative_cutoff is treated as negative.
population_cutoff (float) – The criteria for the difference between the sum of propars and reference population.
newton_options (dict) – Other options for the solver.
- Returns:
Optimized proatom parameters.
- Return type:
1D np.ndarray
- Raises:
RuntimeError – If the inner iteration does not converge.
- solver_newton(bs_funcs, rho, propars, points, weights, threshold, logger, density_cutoff, negative_cutoff, population_cutoff, **newton_options)
Optimize parameters for pro-atom density functions using Newton method
- Parameters:
bs_funcs (2D np.ndarray) – Basis functions array with shape (M, N), where ‘M’ is the number of basis functions and ‘N’ is the number of grid points.
rho (1D np.ndarray) – Spherically-averaged atomic density as a function of radial distance, with shape (N,).
propars (1D np.ndarray) – Pro-atom parameters with shape (M). ‘M’ is the number of basis functions.
points (1D np.ndarray) – Radial coordinates of grid points, with shape (N,).
weights (1D np.ndarray) – Weights for integration, including the angular part (4πr²), with shape (N,).
threshold (float) – Convergence threshold for the iterative process.
logger (logging.Logger or None) – The log object.
density_cutoff (float) – Density values below this cutoff are considered invalid.
negative_cutoff (float) – The value less than negative_cutoff is treated as negative.
population_cutoff (float) – The criteria for the difference between the sum of propars and reference population.
newton_options (dict) – Other options for the solver.
- Returns:
Optimized proatom parameters.
- Return type:
1D np.ndarray
- Raises:
RuntimeError – If the inner iteration does not converge.
- solver_quasi_newton(bs_funcs, rho, propars, points, weights, threshold, logger, density_cutoff, negative_cutoff, population_cutoff, **newton_options)
Optimize parameters for pro-atom density functions using Quasi-Newton method
- Parameters:
bs_funcs (2D np.ndarray) – Basis functions array with shape (M, N), where ‘M’ is the number of basis functions and ‘N’ is the number of grid points.
rho (1D np.ndarray) – Spherically-averaged atomic density as a function of radial distance, with shape (N,).
propars (1D np.ndarray) – Pro-atom parameters with shape (M). ‘M’ is the number of basis functions.
points (1D np.ndarray) – Radial coordinates of grid points, with shape (N,).
weights (1D np.ndarray) – Weights for integration, including the angular part (4πr²), with shape (N,).
threshold (float) – Convergence threshold for the iterative process.
logger (logging.Logger or None) – The log object.
density_cutoff (float) – Density values below this cutoff are considered invalid.
negative_cutoff (float) – The value less than negative_cutoff is treated as negative.
population_cutoff (float) – The criteria for the difference between the sum of propars and reference population.
newton_options (dict) – Other options for the solver.
- Returns:
Optimized proatom parameters.
- Return type:
1D np.ndarray
- Raises:
RuntimeError – If the inner iteration does not converge.
- solver_sc(bs_funcs, rho, propars, points, weights, threshold, logger, density_cutoff, negative_cutoff, population_cutoff, max_niter_inner=100000)
Optimize parameters for proatom density functions using a self-consistent (SC) method.
This approach analytically computes the parameters, aiming to yield results comparable to those obtained via L-ISA algorithms, which require non-negative parameters.
\[N_{Ai} = \int \rho_A(r) \frac{\rho_{Ai}^0(r)}{\rho_A^0(r)} dr\]- Parameters:
bs_funcs (2D np.ndarray) – Basis functions array with shape (M, N), where ‘M’ is the number of basis functions and ‘N’ is the number of grid points.
rho (1D np.ndarray) – Spherically-averaged atomic density as a function of radial distance, with shape (N,).
propars (1D np.ndarray) – Pro-atom parameters with shape (M). ‘M’ is the number of basis functions.
points (1D np.ndarray) – Radial coordinates of grid points, with shape (N,).
weights (1D np.ndarray) – Weights for integration, including the angular part (4πr²), with shape (N,).
threshold (float) – Convergence threshold for the iterative process.
logger (logging.Logger or None) – The log object.
density_cutoff (float) – Density values below this cutoff are considered invalid.
negative_cutoff (float) – The value less than negative_cutoff is treated as negative.
population_cutoff (float) – The criteria for the difference between the sum of propars and reference population.
max_niter_inner (int) – The maximum number of inner iterations.
- Returns:
Optimized proatom parameters.
- Return type:
1D np.ndarray
- Raises:
RuntimeError – If the inner iteration does not converge.
Notes
The method iteratively optimizes the proatom density function parameters. In each iteration, the basis functions and current parameters are used to compute updated parameters, assessing convergence against the specified threshold.
- solver_sc_1_iter(bs_funcs, rho, propars, points, weights, threshold, logger, density_cutoff, negative_cutoff, population_cutoff)
Optimize parameters for proatom density functions using a self-consistent (SC) method.
\[N_{Ai} = \int \rho_A(r) \frac{\rho_{Ai}^0(r)}{\rho_A^0(r)} dr\]- Parameters:
bs_funcs (2D np.ndarray) – Basis functions array with shape (M, N), where ‘M’ is the number of basis functions and ‘N’ is the number of grid points.
rho (1D np.ndarray) – Spherically-averaged atomic density as a function of radial distance, with shape (N,).
propars (1D np.ndarray) – Pro-atom parameters with shape (M). ‘M’ is the number of basis functions.
points (1D np.ndarray) – Radial coordinates of grid points, with shape (N,).
weights (1D np.ndarray) – Weights for integration, including the angular part (4πr²), with shape (N,).
threshold (float) – Convergence threshold for the iterative process.
logger (logging.Logger or None) – The log object.
density_cutoff (float) – Density values below this cutoff are considered invalid.
negative_cutoff (float) – The value less than negative_cutoff is treated as negative.
population_cutoff (float) – The criteria for the difference between the sum of propars and reference population.
- Returns:
Optimized proatom parameters.
- Return type:
1D np.ndarray
- Raises:
RuntimeError – If the inner iteration does not converge.
- solver_sc_plus_cvxopt(bs_funcs, rho, propars, points, weights, threshold, logger, density_cutoff, negative_cutoff, population_cutoff, sc_iter_limit=1000, **cvxopt_options)
Optimize parameters for proatom density functions using a mixing of self-consistent (SC) method and convex optimization method.
- Parameters:
bs_funcs (2D np.ndarray) – Basis functions array with shape (M, N), where ‘M’ is the number of basis functions and ‘N’ is the number of grid points.
rho (1D np.ndarray) – Spherically-averaged atomic density as a function of radial distance, with shape (N,).
propars (1D np.ndarray) – Pro-atom parameters with shape (M). ‘M’ is the number of basis functions.
points (1D np.ndarray) – Radial coordinates of grid points, with shape (N,).
weights (1D np.ndarray) – Weights for integration, including the angular part (4πr²), with shape (N,).
threshold (float) – Convergence threshold for the iterative process.
logger (logging.Logger or None) – The log object.
density_cutoff (float) – Density values below this cutoff are considered invalid.
negative_cutoff (float) – The value less than negative_cutoff is treated as negative.
population_cutoff (float) – The criteria for the difference between the sum of propars and reference population.
sc_iter_limit (int) – The number of iteration steps of self-consistent method.
cvxopt_options (dict) – The options for cvxopt solver. See solver_cvxopt.
- Returns:
Optimized proatom parameters.
- Return type:
1D np.ndarray
- Raises:
RuntimeError – If the inner iteration does not converge.
- solver_trust_region(bs_funcs, rho, propars, points, weights, threshold, logger, density_cutoff, negative_cutoff, population_cutoff, explicit_constr=True, trust_region_options=None)
Optimize parameters for pro-atom density functions using trust-region method.
- Parameters:
bs_funcs (2D np.ndarray) – Basis functions array with shape (M, N), where ‘M’ is the number of basis functions and ‘N’ is the number of grid points.
rho (1D np.ndarray) – Spherically-averaged atomic density as a function of radial distance, with shape (N,).
propars (1D np.ndarray) – Pro-atom parameters with shape (M). ‘M’ is the number of basis functions.
points (1D np.ndarray) – Radial coordinates of grid points, with shape (N,).
weights (1D np.ndarray) – Weights for integration, including the angular part (4πr²), with shape (N,).
threshold (float) – Convergence threshold for the iterative process.
logger (logging.Logger or None) – The log object.
density_cutoff (float) – Density values below this cutoff are considered invalid.
negative_cutoff (float) – The value less than negative_cutoff is treated as negative.
population_cutoff (float) – The criteria for the difference between the sum of propars and reference population.
explicit_constr (bool (optional)) – Whether adding an explicit constraint for the atomic population. The default is False.
trust_region_options (dict) – Other options for trust_region solver.
- Returns:
Optimized proatom parameters.
- Return type:
1D np.ndarray
- Raises:
RuntimeError – If the inner iteration does not converge.