Inverter

class n2v.Inverter(wfn, pbs='same', debug=False)[source]

Bases: n2v.methods.direct.Direct, n2v.methods.wuyang.WuYang, n2v.methods.zmp.ZMP, n2v.methods.mrks.MRKS, n2v.methods.oucarter.OC, n2v.methods.pdeco.PDECO, n2v.grid.grider.Grider

wfnpsi4.core.{RHF, UHF, RKS, UKS, Wavefunction, CCWavefuncion…}

Psi4 wavefunction object

molpsi4.core.Molecule

Psi4 molecule object

basispsi4.core.BasisSet

Psi4 basis set object

basis_strstr

Basis set

nbfint

Number of basis functions for main calculation

nalphaint

Number of alpha electrons

nbetaint

Number of beta electrons

ref{1,2}

Reference calculation 1 -> Restricted 2 -> Unrestricted

DtList

List of psi4.core.Matrix for target density matrices (on AO).

ctList

List of psi4.core.Matrix for input occupied orbitals. This might not be correct for post-HartreeFock methods.

pbs_str: string

name of Potential basis set

pbspsi4.core.BasisSet

Potential basis set.

npbsint

the length of pbs

v_pbsnp.ndarray shape (npbs, ) for ref==1 and (2*npbs, ) for ref==2.

potential vector on the Potential Baiss Set. If the potential is not represented on the basis set, this should remain 0. It will be initialized to a 0 array. One can set this value for initial guesses before Wu-Yang method (WY) or PDE-Constrained Optimization method (PDE-CO). For example, if PDE-CO is ran after a WY calculation, the initial for PDE-CO will be the result of WY if v_pbs is not zeroed.

S2np.ndarray

The ao overlap matrix (i.e. S matrix)

S3np.ndarray

The three ao overlap matrix (ao, ao, pbs)

S4np.ndarray

The four ao overlap matrix, the size should be (ao, ao, ao, ao)

jkpsi4.core.JK

Psi4 jk object. Built if wfn has no jk, otherwise use wfn.jk

Tnp.ndarray

kinetic matrix on ao

Vnp.ndarray

external potential matrix on ao

T_pbs: np.ndarray

kinetic matrix on pbs. Useful for regularization.

guide_potential_components: list of string

guide potential components name

va, vb: np.ndarray of shape (nbasis, nbasis)

guide potential Fock matrix.

Methods Summary

diagonalize(matrix, ndocc)

Diagonalizes Fock Matrix

form_jk(Cocc_a, Cocc_b)

Generates Coulomb and Exchange matrices from occupied orbitals

generate_components(guide_potential_components)

Generates exact potential components to be added to the Hamiltonian to aide in the inversion procedure.

generate_jk([gen_K])

Creates jk object for generation of Coulomb and Exchange matrices 1.0e9 B -> 1.0 GB

generate_mints_matrices()

Generates matrices that are methods of a mints object

invert(method[, guide_potential_components, …])

Handler to all available inversion methods

Methods Documentation

diagonalize(matrix, ndocc)[source]

Diagonalizes Fock Matrix

Parameters
  • marrix (np.ndarray) – Matrix to be diagonalized

  • ndocc (int) – Number of occupied orbitals

Returns

  • C (np.ndarray) – Orbital Matrix

  • Cocc (np.ndarray) – Occupied Orbital Matrix

  • D (np.ndarray) – Density Matrix

  • eigves (np.ndarray) – Eigenvalues

form_jk(Cocc_a, Cocc_b)[source]

Generates Coulomb and Exchange matrices from occupied orbitals

generate_components(guide_potential_components)[source]

Generates exact potential components to be added to the Hamiltonian to aide in the inversion procedure.

guide_potential_components: list

Components added as to guide inversion. Can be chosen from {“fermi_amandi”, “svwn”}

generate_jk(gen_K=False)[source]

Creates jk object for generation of Coulomb and Exchange matrices 1.0e9 B -> 1.0 GB

generate_mints_matrices()[source]

Generates matrices that are methods of a mints object

invert(method, guide_potential_components=['fermi_amaldi'], opt_max_iter=50, **keywords)[source]

Handler to all available inversion methods

method: str

Method used to invert density. Can be chosen from {wuyang, zmp, mrks, oc}. See documentation below for each method.

guide_potential_components: list, opt

Components added as to guide inversion. Can be chosen from {“fermi_amandi”, “svwn”} Default: [“fermi_amaldi”]

opt_max_iter: int, opt

Maximum number of iterations inside the chosen inversion. Default: 50

Direct inversion of a set of Kohn-Sham equations. $$v_{xc}(r) =

rac{1}{n(r)} sum_i^N [phi_i^{*} (r) abla^2 phi_i(r) + arepsilon_i | phi_i(r)|^2] $$

grid: np.ndarray, opt

Grid where result will be expressed in. If not provided, dft grid will be used instead.

the Wu-Yang method: The Journal of chemical physics 118.6 (2003): 2498-2509.

opt_max_iter: int

maximum iteration

opt_method: string, opt

Method for scipy optimizer Currently only used by wuyang and pdeco method. Defaul: ‘trust-krylov’ https://docs.scipy.org/doc/scipy/reference/generated/scipy.optimize.minimize.html

regfloat, opt

Regularization constant for Wuyant Inversion. Default: None -> No regularization is added. Becomes attribute of inverter -> inverter.lambda_reg

tol: float

tol for scipy.optimize.minimize

gtol: float

gtol for scipy.optimize.minimize: the gradient norm for convergence

opt: dict
options for scipy.optimize.minimize

Notice that opt has lower priorities than opt_max_iter and gtol.

return:

the result are stored in self.v_pbs

The Zhao-Morrison-Parr Method: Phys. Rev. A 50, 2138

lambda_list: list

List of Lamda parameters used as a coefficient for Hartree difference in SCF cycle.

zmp_mixing: float, optional

mixing in [0,1]. How much of the new potential is added in. For example, zmp_mixing = 0 means the traditional ZMP, i.e. all the potentials from previous smaller lambda are ignored. Zmp_mixing = 1 means that all the potentials of previous lambdas are accumulated, the larger lambda potential are meant to fix the wrong/inaccurate region of the potential of the sum of the previous potentials instead of providing an entire new potentials. default: 1

opt_max_iter: float

Maximum number of iterations for scf cycle

opt_tol: float

Convergence criteria set for Density Difference and DIIS error.

return:

The result will be stored in self.proto_density_a and self.proto_density_b For zmp_mixing==1, restricted (ref==1):

self.proto_density_a = sum_i lambda_i * (Da_i - Dt[0]) - 1/N * (Dt[0]) self.proto_density_b = sum_i lambda_i * (Db_i - Dt[1]) - 1/N * (Dt[1]);

unrestricted (ref==1):

self.proto_density_a = sum_i lambda_i * (Da_i - Dt[0]) - 1/N * (Dt[0] + Dt[1]) self.proto_density_b = sum_i lambda_i * (Db_i - Dt[1]) - 1/N * (Dt[0] + Dt[1]);

For restricted (ref==1):

vxc = int dr’

rac{self.proto_density_a + self.proto_density_b}{|r-r'|}

= 2 * int dr’

rac{self.proto_density_a}{|r-r'|};
for unrestricted (ref==2):

vxc_up = int dr’

rac{self.proto_density_a}{|r-r'|}

vxc_down = int dr’

rac{self.proto_density_b}{|r-r'|}.
To get potential on grid, one needs to do

vxc = self.on_grid_esp(Da=self.proto_density_a, Db=self.proto_density_b, grid=grid) for restricted; vxc_up = self.on_grid_esp(Da=self.proto_density_a, Db=np.zeros_like(self.proto_density_a),

grid=grid) for unrestricted;

the modified Ryabinkin-Kohut-Staroverov method: Phys. Rev. Lett. 115, 083001 J. Chem. Phys. 146, 084103p

maxiter: int

same as opt_max_iter

vxc_grid: np.ndarray of shape (3, num_grid_points), opt

When this is given, the final result will be represented

v_tol: float, opt

convergence criteria for vxc Fock matrices. default: 1e-4

D_tol: float, opt

convergence criteria for density matrices. default: 1e-7

eig_tol: float, opt

convergence criteria for occupied eigenvalue spectrum. default: 1e-4

frac_old: float, opt

Linear mixing parameter for current vxc and old vxc. If 0, no old vxc is mixed in. Should be in [0,1) default: 0.5.

init: string or psi4.core.Wavefunction, opt

Initial guess method. default: “SCAN” 1) If None, input wfn info will be used as initial guess. 2) If “continue” is given, then it will not initialize but use the densities and orbitals stored. Meaningly, one can run a quick WY calculation as the initial guess. This can also be used to user speficified initial guess by setting Da, Coca, eigvec_a. 3) If it’s not continue, it would be expecting a method name string that works for psi4. A separate psi4 calculation would be performed.

sing: tuple of float of length 4, opt.

Singularity parameter for _vxc_hole_quadrature() default: (1e-5, 1e-4, 1e-5, 1e-4) [0]: atol, [1]: atol1 for dft_spherical grid calculation. [2]: atol, [3]: atol1 for vxc_grid calculation.

return:

The result will be stored in self.grid.vxc

Ou-Carter method J. Chem. Theory Comput. 2018, 14, 5680−5689

maxiter: int

same as opt_max_iter

vxc_grid: np.ndarray of shape (3, num_grid_points)

The final result will be represented on this grid default: 1e-4

D_tol: float, opt

convergence criteria for density matrices. default: 1e-7

eig_tol: float, opt

convergence criteria for occupied eigenvalue spectrum. default: 1e-4

frac_old: float, opt

Linear mixing parameter for current vxc and old vxc. If 0, no old vxc is mixed in. Should be in [0,1) default: 0.5.

init: string, opt

Initial guess method. default: “SCAN” 1) If None, input wfn info will be used as initial guess. 2) If “continue” is given, then it will not initialize but use the densities and orbitals stored. Meaningly, one can run a quick WY calculation as the initial guess. This can also be used to user speficified initial guess by setting Da, Coca, eigvec_a. 3) If it’s not continue, it would be expecting a method name string that works for psi4. A separate psi4 calculation would be performed. wuyang

the PDE-Constrained Optimization method: Int J Quantum Chem. 2018;118:e25425; Nat Commun 10, 4497 (2019).

opt_max_iter: int

maximum iteration

opt_method: string, opt

Method for scipy optimizer Currently only used by wuyang and pdeco method. Defaul: ‘L-BFGS-B’ Options: [‘L-BFGS-B’, ‘BFGS’] https://docs.scipy.org/doc/scipy/reference/generated/scipy.optimize.minimize.html

regfloat, opt

Regularization constant for Wuyant Inversion. Default: None -> No regularization is added. Becomes attribute of inverter -> inverter.lambda_reg

gtol: float

gtol for scipy.optimize.minimize: the gradient norm for convergence

opt: dict

options for scipy.optimize.minimize Notice that opt has lower priorities than opt_max_iter and gtol.

return:

the result are stored in self.v_pbs