pyscf.df.grad package#
Submodules#
pyscf.df.grad.casdm2_util module#
- pyscf.df.grad.casdm2_util.energy_elec_dferi(mc, mo_cas=None, ci=None, dfcasdm2=None, casdm2=None)[source]#
Evaluate E2 = (P|ij) d_Pij / 2, where d_Pij is the DF-2rdm obtained by solve_df_rdm2. For testing purposes. Note that the only index permutation this function understands is (P|ij) = (P|ji) if i and j span the same range of MOs. The caller has to handle everything else, including, for instance, multiplication by 2 if a nonsymmetric slice of the 2RDM is used.
- Args:
mc: MC-SCF energy method object
- Kwargs:
- mo_cas: ndarray, list, or tuple containing active-space MO coefficients
If a tuple of length 2, the same pair of MO sets are assumed to apply to the internally-contracted and externally-contracted indices of the DF-2rdm: (P|Q)d_Qij = (P|kl)d_ijkl -> (P|Q)d_Qij = (P|ij)d_ijij If a tuple of length 4, the 4 MO sets are applied to ijkl above in that order (first two external, last two internal).
- ci: ndarray, tuple, or list containing CI coefficients in mo_cas basis.
Not used if dfcasdm2 is provided.
- dfcasdm2: ndarray, tuple, or list containing DF-2rdm in mo_cas basis.
Computed by solve_df_rdm2 if omitted.
- casdm2: ndarray, tuple, or list containing rdm2 in mo_cas basis.
Computed by mc_or_mc_grad.fcisolver.make_rdm12 (ci,…) if omitted.
- Returns:
- energy: list
List of energies corresponding to the dfcasdm2s, E = (P|ij) d_Pij / 2 = (P|ij) (P|Q)^-1 (Q|kl) d_ijkl / 2
- pyscf.df.grad.casdm2_util.get_int3c_mo(mol, auxmol, mo_coeff, compact=True, max_memory=None)[source]#
Evaluate (P|uv) c_ui c_vj -> (P|ij)
- Args:
mol: gto.Mole auxmol: gto.Mole, contains auxbasis mo_coeff: ndarray, list, or tuple containing MO coefficients
if two ndarrays mo_coeff = (mo0, mo1) are provided, mo0 and mo1 are used for the two AO dimensions
- Kwargs:
- compact: bool
If true, will return only unique ERIs along the two MO dimensions. Does nothing if mo_coeff contains two different sets of orbitals.
- max_memory: int
Maximum memory consumption in MB
- Returns:
int3c: ndarray of shape (naux, nmo0, nmo1) or (naux, nmo*(nmo+1)//2)
- pyscf.df.grad.casdm2_util.gfock_dferi(mc, mo_cas=None, ci=None, dfcasdm2=None, casdm2=None, max_memory=None, ao_basis=True)[source]#
Evaluate F_ij = (P|ik) d_Pjk - this was a giant reinvention of the wheel that didn’t need to happen because with_df._cderi is plenty good enough to calculate gfock. Oh well.
- Args:
mc_grad: MC-SCF gradients method object
- Kwargs:
- mc_cas: ndarray, list, or tuple containing active-space MO coefficients
If a tuple of length 2, the same pair of MO sets are assumed to apply to the internally-contracted and externally-contracted indices of the DF-2rdm: (P|Q)d_Qij = (P|kl)d_ijkl -> (P|Q)d_Qij = (P|ij)d_ijij If a tuple of length 4, the 4 MO sets are applied to ijkl above in that order (first two external, last two internal).
- ci: ndarray, tuple, or list containing CI coefficients in mo_cas basis.
Not used if dfcasdm2 is provided.
- dfcasdm2: ndarray, tuple, or list containing DF-2rdm in mo_cas basis.
Computed by solve_df_rdm2 if omitted.
- casdm2: ndarray, tuple, or list containing rdm2 in mo_cas basis.
Computed by mc_grad.fcisolver.make_rdm12 (ci,…) if omitted.
- max_memory: int
Maximum memory usage in MB
- ao_basis: bool
If true, return gfock in AO basis
- Returns:
gfock: ndarray of shape (nset, nmo[0], nmo[1]) or (nset, nao, nao)
- pyscf.df.grad.casdm2_util.grad_elec_auxresponse_dferi(mc_grad, mo_cas=None, ci=None, dfcasdm2=None, casdm2=None, atmlst=None, max_memory=None, dferi=None, incl_2c=True)[source]#
Evaluate the [(P’|ij) + (P’|Q) g_Qij] d_Pij contribution to the electronic gradient, where d_Pij is the DF-2RDM obtained by solve_df_rdm2 and g_Qij solves (P|Q) g_Qij = (P|ij). The caller must symmetrize if necessary (i.e., (P|Q) d_Qij = (P|kl) d_ijkl <-> (P|Q) d_Qkl = (P|ij) d_ijkl in order to get at Q’). Args:
mc_grad: MC-SCF gradients method object
- Kwargs:
- mc_cas: ndarray, list, or tuple containing active-space MO coefficients
If a tuple of length 2, the same pair of MO sets are assumed to apply to the internally-contracted and externally-contracted indices of the DF-2rdm: (P|Q)d_Qij = (P|kl)d_ijkl -> (P|Q)d_Qij = (P|ij)d_ijij If a tuple of length 4, the 4 MO sets are applied to ijkl above in that order (first two external, last two internal).
- ci: ndarray, tuple, or list containing CI coefficients in mo_cas basis.
Not used if dfcasdm2 is provided.
- dfcasdm2: ndarray, tuple, or list containing DF-2rdm in mo_cas basis.
Computed by solve_df_rdm2 if omitted.
- casdm2: ndarray, tuple, or list containing rdm2 in mo_cas basis.
Computed by mc_grad.fcisolver.make_rdm12 (ci,…) if omitted.
- atmlst: list of integers
List of nonfrozen atoms, as in grad_elec functions. Defaults to list (range (mol.natm))
- max_memory: int
Maximum memory usage in MB
dferi: ndarray containing g_Pij for optional precalculation incl_2c: bool
If False, omit the terms depending on (P’|Q)
- Returns:
dE: list of ndarray of shape (len (atmlst), 3)
- pyscf.df.grad.casdm2_util.grad_elec_dferi(mc_grad, mo_cas=None, ci=None, dfcasdm2=None, casdm2=None, atmlst=None, max_memory=None)[source]#
Evaluate the (P|i’j) d_Pij contribution to the electronic gradient, where d_Pij is the DF-2RDM obtained by solve_df_rdm2. The caller must symmetrize (i.e., [(P|i’j) + (P|ij’)] d_Pij / 2) if necessary.
- Args:
mc_grad: MC-SCF gradients method object
- Kwargs:
- mc_cas: ndarray, list, or tuple containing active-space MO coefficients
If a tuple of length 2, the same pair of MO sets are assumed to apply to the internally-contracted and externally-contracted indices of the DF-2rdm: (P|Q)d_Qij = (P|kl)d_ijkl -> (P|Q)d_Qij = (P|ij)d_ijij If a tuple of length 4, the 4 MO sets are applied to ijkl above in that order (first two external, last two internal).
- ci: ndarray, tuple, or list containing CI coefficients in mo_cas basis.
Not used if dfcasdm2 is provided.
- dfcasdm2: ndarray, tuple, or list containing DF-2rdm in mo_cas basis.
Computed by solve_df_rdm2 if omitted.
- casdm2: ndarray, tuple, or list containing rdm2 in mo_cas basis.
Computed by mc_grad.fcisolver.make_rdm12 (ci,…) if omitted.
- atmlst: list of integers
List of nonfrozen atoms, as in grad_elec functions. Defaults to list (range (mol.natm))
- max_memory: int
Maximum memory usage in MB
- Returns:
dE: ndarray of shape (len (dfcasdm2), len (atmlst), 3)
- pyscf.df.grad.casdm2_util.solve_df_eri(mc_or_mc_grad, mo_cas=None, compact=True)[source]#
Solve (P|Q) g_Qij = (P|ij) for g_Qij using MOs i,j. I mean this should be a basic function but whatever.
- pyscf.df.grad.casdm2_util.solve_df_rdm2(mc_or_mc_grad, mo_cas=None, ci=None, casdm2=None)[source]#
Solve (P|Q)d_Qij = (P|kl)d_ijkl for d_Qij in the MO basis.
- Args:
mc_or_mc_grad: DF-MCSCF energy or gradients method object.
- Kwargs:
- mo_cas: ndarray, tuple, or list containing active mo coefficients.
if two ndarrays mo_cas = (mo0, mo1) are provided, mo0 and mo1 are assumed to correspond to casdm2’s LAST two dimensions in that order, regardless of len (ci) or len (casdm2). (This will facilitate SA-CASSCF gradients at some point. Note the difference from grad_elec_dferi!)
- ci: ndarray, tuple, or list containing CI coefficients in mo_cas basis.
Not used if casdm2 is provided.
- casdm2: ndarray, tuple, or list containing rdm2 in mo_cas basis.
Computed by mc_or_mc_grad.fcisolver.make_rdm12 (ci,…) if omitted.
- compact: bool
If true, tries to return d_Pqr in lower-triangular form if possible
- Returns:
- dfcasdm2: ndarray or list containing 3-center 2RDM, d_Pqr, where P is
auxbasis index and q, r are mo_cas basis indices.
pyscf.df.grad.casscf module#
CASSCF analytical nuclear gradients
Ref. J. Comput. Chem., 5, 589
MRH: copied from pyscf.grad.casscf.py on 12/07/2019 Contains my modifications for SA-CASSCF gradients 1. Generalized Fock has nonzero i->a and u->a 2. Memory footprint for differentiated eris bugfix
- class pyscf.df.grad.casscf.CASSCF_GradScanner(g)[source]#
Bases:
GradScanner
- class pyscf.df.grad.casscf.Gradients(mc)[source]#
Bases:
Gradients
Non-relativistic restricted Hartree-Fock gradients
- as_scanner()#
Generating a nuclear gradients scanner/solver (for geometry optimizer).
The returned solver is a function. This function requires one argument “mol” as input and returns energy and first order nuclear derivatives.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the nuc-grad object and SCF object (DIIS, conv_tol, max_memory etc) are automatically applied in the solver.
Note scanner has side effects. It may change many underlying objects (_scf, with_df, with_x2c, …) during calculation.
Examples:
>>> from pyscf import gto, scf, mcscf >>> mol = gto.M(atom='N 0 0 0; N 0 0 1.1', verbose=0) >>> mc_grad_scanner = mcscf.CASSCF(scf.RHF(mol), 4, 4).nuc_grad_method().as_scanner() >>> etot, grad = mc_grad_scanner(gto.M(atom='N 0 0 0; N 0 0 1.1')) >>> etot, grad = mc_grad_scanner(gto.M(atom='N 0 0 0; N 0 0 1.5'))
- grad_elec(mo_coeff=None, ci=None, atmlst=None, verbose=None)#
- kernel(mo_coeff=None, ci=None, atmlst=None, verbose=None)[source]#
Kernel function is the main driver of a method. Every method should define the kernel function as the entry of the calculation. Note the return value of kernel function is not strictly defined. It can be anything related to the method (such as the energy, the wave-function, the DFT mesh grids etc.).
- to_gpu(out=None)#
Convert a method to its corresponding GPU variant, and recursively converts all attributes of a method to cupy objects or gpu4pyscf objects.
- pyscf.df.grad.casscf.as_scanner(mcscf_grad)[source]#
Generating a nuclear gradients scanner/solver (for geometry optimizer).
The returned solver is a function. This function requires one argument “mol” as input and returns energy and first order nuclear derivatives.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the nuc-grad object and SCF object (DIIS, conv_tol, max_memory etc) are automatically applied in the solver.
Note scanner has side effects. It may change many underlying objects (_scf, with_df, with_x2c, …) during calculation.
Examples:
>>> from pyscf import gto, scf, mcscf >>> mol = gto.M(atom='N 0 0 0; N 0 0 1.1', verbose=0) >>> mc_grad_scanner = mcscf.CASSCF(scf.RHF(mol), 4, 4).nuc_grad_method().as_scanner() >>> etot, grad = mc_grad_scanner(gto.M(atom='N 0 0 0; N 0 0 1.1')) >>> etot, grad = mc_grad_scanner(gto.M(atom='N 0 0 0; N 0 0 1.5'))
pyscf.df.grad.lpdft module#
- class pyscf.df.grad.lpdft.Gradients(pdft, state=None)[source]#
-
- get_Aop_Adiag(verbose=None, mo=None, ci=None, eris=None, state=None, **kwargs)[source]#
This function accounts for the fact that the CI vectors are no longer eigenstates of the CAS Hamiltonian. It adds back in the necessary values to the Hessian.
- get_ham_response(**kwargs)[source]#
Return expectation values <dH/dx> where x is nuclear displacement. I.E., the gradient if the method were variational.
- get_otp_gradient_response(mo=None, ci=None, state=0)[source]#
Generate the 1- and 2-body on-top gradient response terms which have been partially contracted with the Delta density generated from state.
- Args:
- mondarray of shape (nao,nmo)
A full set of molecular orbital coefficients. Taken from self if not provided.
- cilist of ndarrays of length nroots
CI vectors should be from a converged L-PDFT calculation
- stateint
State to generate the Delta density with
- Returns:
- feff1ndarray of shape (nao, nao) 1-particle On-top gradient response which as been contracted with the
Delta density generated from state. Should include the Coulomb term as well.
- feff2pyscf.mcscf.mc_ao2mo._ERIS instance Relevant 2-body on-top gradient response terms in the MO
basis. Also, partially contracted with the Delta density.
- get_wfn_response(state=None, verbose=None, mo=None, ci=None, feff1=None, feff2=None, **kwargs)[source]#
Returns the derivative of the L-PDFT energy for the given state with respect to MO parameters and CI parameters. Care is take to account for the implicit and explicit response terms, and make sure the CI vectors are properly projected out.
- Args:
- stateint
Which state energy to get response of.
- mondarray of shape (nao, nmo)
A full set of molecular orbital coefficients. Taken from self if not provided.
- cilist of ndarrays of length nroots
CI vectors should be from a converged L-PDFT calculation.
- feff1ndarray of shape (nao, nao) 1-particle On-top gradient response which as been contracted with the
Delta density generated from state. Should include the Coulomb term as well.
- feff2pyscf.mcscf.mc_ao2mo._ERIS instance Relevant 2-body on-top gradient response terms in the MO
basis. Also, partially contracted with the Delta density.
Returns: g_all : ndarray of shape nlag First sector [:self.ngorb] contains the derivatives with respect to MO parameters. Second sector [self.ngorb:] contains the derivatives with respect to CI parameters.
- kernel(**kwargs)[source]#
Kernel function is the main driver of a method. Every method should define the kernel function as the entry of the calculation. Note the return value of kernel function is not strictly defined. It can be anything related to the method (such as the energy, the wave-function, the DFT mesh grids etc.).
pyscf.df.grad.mcpdft module#
- class pyscf.df.grad.mcpdft.Gradients(pdft, state=None)[source]#
-
- get_ham_response(**kwargs)[source]#
Return expectation values <dH/dx> where x is nuclear displacement. I.E., the gradient if the method were variational.
- get_init_guess(bvec, Adiag, Aop, precond)[source]#
Initial guess should solve the problem for SA-SA rotations
- get_wfn_response(state=None, verbose=None, mo=None, ci=None, veff1=None, veff2=None, nlag=None, **kwargs)[source]#
Return first derivative of the energy wrt wave function parameters conjugate to the Lagrange multipliers. Used to calculate the value of the Lagrange multipliers.
pyscf.df.grad.mspdft module#
- class pyscf.df.grad.mspdft.Gradients(pdft)[source]#
Bases:
Gradients
- get_ham_response(**kwargs)[source]#
write mspdft heff Hellmann-Feynman calculator; sum over diagonal PDFT Hellmann-Feynman terms
- make_fcasscf(state=None, casscf_attr={}, fcisolver_attr={})[source]#
SA-CASSCF nuclear gradients require 1) first derivatives wrt wave function variables and nuclear shifts of the target state’s energy, AND 2) first and second derivatives of the objective function used to determine the MO coefficients and CI vectors. This function addresses 1).
- Kwargs:
- stateinteger
The specific state whose energy is being differentiated. This kwarg is necessary in the context of state_average_mix, where the number of electrons and the make_rdm* functions differ from state to state.
- casscf_attrdictionary
Extra attributes to apply to fcasscf. Relevant to child methods (i.e., MC-PDFT; NACs)
- fcisolver_attrdictionary
Extra attributes to apply to fcasscf.fcisolver. Relevant to child methods (i.e., MC-PDFT; NACs)
- Returns:
- fcasscfobject of
mc1step.CASSCF
Set up to evaluate first derivatives of state “state”. Only functions, classes, and the nelecas variable are set up; the caller should assign MO coefficients and CI vectors explicitly post facto.
- fcasscfobject of
pyscf.df.grad.rhf module#
- class pyscf.df.grad.rhf.Gradients(mf)[source]#
Bases:
Gradients
Restricted density-fitting Hartree-Fock gradients
- auxbasis_response = True#
- check_sanity()[source]#
Check input of class/object attributes, check whether a class method is overwritten. It does not check the attributes which are prefixed with “_”. The return value of method set is the object itself. This allows a series of functions/methods to be executed in pipe.
- extra_force(atom_id, envs)[source]#
Hook for extra contributions in analytical gradients.
Contributions like the response of auxiliary basis in density fitting method, the grid response in DFT numerical integration can be put in this function.
- pyscf.df.grad.rhf.get_jk(mf_grad, mol=None, dm=None, hermi=0, with_j=True, with_k=True, decompose_j2c='CD', lindep=1e-07)[source]#
Computes J and K derivatives with density fitting
- Args:
mf_grad : instance of
Gradients
mol : instance of
gto.Mole
- dm: numpy.ndarray
Zeroth order density matrix
- hermi: int
Is the density matrix hermitian or not
- with_j: boolean
Whether to compute J matrix
- with_k: boolean
Whether to compute K matrix
- decompose_j2c: string
The method to decompose the metric defined by int2c. It can be set to CD (cholesky decomposition) or ED (eigenvalue decomposition).
- lindepfloat
The threshold to discard linearly dependent basis when decompose_j2c is set to ED.
pyscf.df.grad.rks module#
- class pyscf.df.grad.rks.Gradients(mf)[source]#
Bases:
Gradients
- auxbasis_response = True#
- extra_force(atom_id, envs)[source]#
Hook for extra contributions in analytical gradients.
Contributions like the response of auxiliary basis in density fitting method, the grid response in DFT numerical integration can be put in this function.
- get_jk(mol=None, dm=None, hermi=0, with_j=True, with_k=True, omega=None)[source]#
J = ((-nabla i) j| kl) D_lk K = ((-nabla i) j| kl) D_jk
- get_veff(mol=None, dm=None)#
Coulomb + XC functional
pyscf.df.grad.rohf module#
pyscf.df.grad.roks module#
pyscf.df.grad.sacasscf module#
- class pyscf.df.grad.sacasscf.CASSCF_GradScanner(g, state)[source]#
Bases:
GradScanner
- class pyscf.df.grad.sacasscf.Gradients(mc, state=None)[source]#
Bases:
Gradients
- get_LdotJnuc(Lvec, **kwargs)[source]#
Return Lvec . J_nuc, where J_nuc is the Jacobian of the Lagrange cofactors wrt nuclear displacement. This is the second term of the final gradient expectation value.
- kernel(**kwargs)[source]#
Kernel function is the main driver of a method. Every method should define the kernel function as the entry of the calculation. Note the return value of kernel function is not strictly defined. It can be anything related to the method (such as the energy, the wave-function, the DFT mesh grids etc.).
- to_gpu(out=None)#
Convert a method to its corresponding GPU variant, and recursively converts all attributes of a method to cupy objects or gpu4pyscf objects.
- pyscf.df.grad.sacasscf.Lci_dot_dgci_dx(Lci, weights, mc, mo_coeff=None, ci=None, atmlst=None, mf_grad=None, eris=None, verbose=None, auxbasis_response=True)[source]#
Modification of pyscf.grad.casscf.kernel to compute instead the CI Lagrange term nuclear gradient (sum_IJ Lci_IJ d2_Ecas/d_lambda d_PIJ) This involves removing all core-core and nuclear-nuclear terms and making the substitution sum_I w_I<L_I|p’q|I> + c.c. -> <0|p’q|0> sum_I w_I<L_I|p’r’sq|I> + c.c. -> <0|p’r’sq|0> The active-core terms (sum_I w_I<L_I|x’iyi|I>, sum_I w_I <L_I|x’iiy|I>, c.c.) must be retained.
- pyscf.df.grad.sacasscf.Lorb_dot_dgorb_dx(Lorb, mc, mo_coeff=None, ci=None, atmlst=None, mf_grad=None, eris=None, verbose=None, auxbasis_response=True)[source]#
Modification of pyscf.grad.casscf.kernel to compute instead the orbital Lagrange term nuclear gradient (sum_pq Lorb_pq d2_Ecas/d_lambda d_kpq) This involves removing nuclear-nuclear terms and making the substitution (D_[p]q + D_p[q]) -> D_pq (d_[p]qrs + d_pq[r]s + d_p[q]rs + d_pqr[s]) -> d_pqrs Where [] around an index implies contraction with Lorb from the left, so that the external index (regardless of whether the index on the rdm is bra or ket) is always the first index of Lorb.
- pyscf.df.grad.sacasscf.as_scanner(mcscf_grad, state=None)[source]#
Generating a nuclear gradients scanner/solver (for geometry optimizer).
The returned solver is a function. This function requires one argument “mol” as input and returns energy and first order nuclear derivatives.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the nuc-grad object and SCF object (DIIS, conv_tol, max_memory etc) are automatically applied in the solver.
Note scanner has side effects. It may change many underlying objects (_scf, with_df, with_x2c, …) during calculation.
Examples:
>>> from pyscf import gto, scf, mcscf >>> mol = gto.M(atom='N 0 0 0; N 0 0 1.1', verbose=0) >>> mc_grad_scanner = mcscf.CASSCF(scf.RHF(mol), 4, 4).nuc_grad_method().as_scanner() >>> etot, grad = mc_grad_scanner(gto.M(atom='N 0 0 0; N 0 0 1.1')) >>> etot, grad = mc_grad_scanner(gto.M(atom='N 0 0 0; N 0 0 1.5'))
pyscf.df.grad.uhf module#
- class pyscf.df.grad.uhf.Gradients(mf)[source]#
Bases:
Gradients
Unrestricted density-fitting Hartree-Fock gradients
- auxbasis_response = True#
- extra_force(atom_id, envs)[source]#
Hook for extra contributions in analytical gradients.
Contributions like the response of auxiliary basis in density fitting method, the grid response in DFT numerical integration can be put in this function.
pyscf.df.grad.uks module#
- class pyscf.df.grad.uks.Gradients(mf)[source]#
Bases:
Gradients
- auxbasis_response = True#
- extra_force(atom_id, envs)[source]#
Hook for extra contributions in analytical gradients.
Contributions like the response of auxiliary basis in density fitting method, the grid response in DFT numerical integration can be put in this function.
- get_jk(mol=None, dm=None, hermi=0, with_j=True, with_k=True, omega=None)[source]#
J = ((-nabla i) j| kl) D_lk K = ((-nabla i) j| kl) D_jk
- get_veff(mol=None, dm=None)#
Coulomb + XC functional