pyscf.grad package#
Submodules#
pyscf.grad.casci module#
CASCI analytical nuclear gradients
Ref. J. Comput. Chem., 5, 589
- class pyscf.grad.casci.CASCI_GradScanner(g, state)[source]#
Bases:
GradScanner
- class pyscf.grad.casci.Gradients(mc)[source]#
Bases:
GradientsBase
Non-relativistic restricted Hartree-Fock gradients
- as_scanner(state=None)#
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.CASCI(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, state=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.grad.casci.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.CASCI(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.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.grad.casscf.CASSCF_GradScanner(g)[source]#
Bases:
GradScanner
- class pyscf.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.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.grad.ccsd module#
CCSD analytical nuclear gradients
- class pyscf.grad.ccsd.CCSD_GradScanner(g)[source]#
Bases:
GradScanner
- property converged#
- class pyscf.grad.ccsd.Gradients(method)[source]#
Bases:
GradientsBase
- 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 total CCSD energy.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the CCSD and the underlying SCF objects (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, cc >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> cc_scanner = cc.CCSD(scf.RHF(mol)).nuc_grad_method().as_scanner() >>> e_tot, grad = cc_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot, grad = cc_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
- grad_elec(t1=None, t2=None, l1=None, l2=None, eris=None, atmlst=None, d1=None, d2=None, verbose=4)#
- kernel(t1=None, t2=None, l1=None, l2=None, eris=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.grad.ccsd.as_scanner(grad_cc)[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 total CCSD energy.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the CCSD and the underlying SCF objects (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, cc >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> cc_scanner = cc.CCSD(scf.RHF(mol)).nuc_grad_method().as_scanner() >>> e_tot, grad = cc_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot, grad = cc_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
pyscf.grad.ccsd_slow module#
RCCSD
Ref: JCP 90, 1752 (1989); DOI:10.1063/1.456069
pyscf.grad.ccsd_t module#
pyscf.grad.ccsd_t_slow module#
pyscf.grad.cisd module#
CISD analytical nuclear gradients
- class pyscf.grad.cisd.CISD_GradScanner(g, state)[source]#
Bases:
GradScanner
- property converged#
- class pyscf.grad.cisd.Gradients(myci)[source]#
Bases:
GradientsBase
- as_scanner(state=0)#
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 total CISD energy.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the CISD and the underlying SCF objects (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, ci >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> ci_scanner = ci.CISD(scf.RHF(mol)).nuc_grad_method().as_scanner() >>> e_tot, grad = ci_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot, grad = ci_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
- grad_elec(civec=None, eris=None, atmlst=None, verbose=4)#
- kernel(civec=None, eris=None, atmlst=None, state=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.grad.cisd.as_scanner(grad_ci, state=0)[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 total CISD energy.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the CISD and the underlying SCF objects (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, ci >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> ci_scanner = ci.CISD(scf.RHF(mol)).nuc_grad_method().as_scanner() >>> e_tot, grad = ci_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot, grad = ci_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
pyscf.grad.cmspdft module#
- pyscf.grad.cmspdft.diab_grad(mc_grad, Lis, atmlst=None, mo=None, ci=None, eris=None, mf_grad=None, **kwargs)[source]#
Computes the partial first derivatives of
Q_a-a = 1/2 sum_I g_pqrs <I|p’q|I> <I|r’s|I>
with respect to geometry perturbation.
- Args:
mc_grad : object of class Gradients (CASSCF or CASCI) Lis : ndarray of shape (nroots*(nroots-1)/2,)
Contains step vector for intermediate state rotations
- Kwargs:
- atmlstlist
List of atoms whose geometries are perturbed. Defaults to all atoms in mc_grad.mol.
- mondarray of shape (nao,nmo)
Contains MO coefficients
- cindarray or list of length (nroots)
Contains intermediate-state CI vectors
- erisobject of class ERIS (CASSCF or CASCI)
Contains (true) ERIs in the MO basis
- mf_grad: object of class Gradients (RHF)
Defaults to mc_grad.base.get_rhf_base ().nuc_grad_method ()
- Returns:
- dendarray of shape (len (atmlst), 3)
Contains gradient vector
- pyscf.grad.cmspdft.diab_grad_o0(mc_grad, Lis, atmlst=None, mo=None, ci=None, eris=None, mf_grad=None, **kwargs)[source]#
Monkeypatch version of diab_grad
- pyscf.grad.cmspdft.diab_response(mc_grad, Lis, mo=None, ci=None, eris=None, **kwargs)[source]#
Computes the Hessian-vector product of
Q_a-a = 1/2 sum_I g_pqrs <I|p’q|I> <I|r’s|I>
where the vector is a vector of intermediate-state rotations and the external derivatives are with respect to orbital rotations and CI transfers.
- Args:
mc_grad : object of class Gradients (CASSCF or CASCI) Lis : ndarray of shape (nroots*(nroots-1)/2,)
Contains step vector for intermediate state rotations
- Kwargs:
- mondarray of shape (nao,nmo)
Contains MO coefficients
- cindarray or list of length (nroots)
Contains intermediate-state CI vectors
- erisobject of class ERIS (CASSCF or CASCI)
Contains (true) ERIs in the MO basis
- Returns:
- Rndarray of shape (mc_grad.ngorb+mc_grad.nci)
Contains Hessian-vector product
pyscf.grad.dhf module#
Relativistic Dirac-Hartree-Fock
- class pyscf.grad.dhf.Gradients(scf_method)[source]#
Bases:
GradientsBase
Unrestricted Dirac-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, grad >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> hf_scanner = scf.RHF(mol).apply(grad.RHF).as_scanner() >>> e_tot, grad = hf_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot, grad = hf_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
- 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.
- grad_elec(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)#
- kernel(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=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.
- class pyscf.grad.dhf.GradientsBase(method)[source]#
Bases:
GradientsBase
Basic nuclear gradient functions for 4C relativistic methods
- pyscf.grad.dhf.get_veff(mol, dm, level='SSSS')#
Dirac-Hartree-Fock Coulomb repulsion
pyscf.grad.dispersion module#
gradient of dispersion correction for HF and DFT
pyscf.grad.lagrange module#
- class pyscf.grad.lagrange.Gradients(method, nlag)[source]#
Bases:
GradientsBase
Dummy parent class for calculating analytical nuclear gradients using the technique of Lagrange multipliers: L = E + sum_i z_i L_i dE/dx = partial L/partial x iff all L_i = 0 for the given wave function I.E., the Lagrange multipliers L_i cancel the direct dependence of the wave function on the nuclear coordinates and allow the Hellmann-Feynman theorem to be used for some non-variational methods.
- property converged#
- get_Aop_Adiag(**kwargs)[source]#
Return a function calculating Lvec . J_wfn, where J_wfn is the Jacobian of the Lagrange cofactors (e.g., in state-averaged CASSCF, the Hessian of the state-averaged energy wrt wfn parameters) along with the diagonal of the Jacobian.
- 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.
- 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_wfn_response(**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.
- kernel(level_shift=None, **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.grad.lpdft module#
- class pyscf.grad.lpdft.Gradients(mc, state=None)[source]#
Bases:
Gradients
- 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(state=None, atmlst=None, verbose=None, mo=None, ci=None, mf_grad=None, feff1=None, feff2=None, **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.grad.mcpdft module#
- class pyscf.grad.mcpdft.Gradients(pdft, state=None)[source]#
Bases:
Gradients
- get_ham_response(state=None, atmlst=None, verbose=None, mo=None, ci=None, eris=None, mf_grad=None, veff1=None, veff2=None, **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.grad.mcpdft.gfock_sym(mc, mo_coeff, casdm1, casdm2, h1e, eris)[source]#
Assume that h2e v_j = v_k
- pyscf.grad.mcpdft.mcpdft_HellmanFeynman_grad(mc, ot, veff1, veff2, mo_coeff=None, ci=None, atmlst=None, mf_grad=None, verbose=None, max_memory=None, auxbasis_response=False)[source]#
Modification of pyscf.grad.casscf.kernel to compute instead the Hellman-Feynman gradient terms of MC-PDFT. From the differentiated Hamiltonian matrix elements, only the core and Coulomb energy parts remain. For the renormalization terms, the effective Fock matrix is as in CASSCF, but with the same Hamiltonian substutition that is used for the energy response terms.
pyscf.grad.mp2 module#
MP2 analytical nuclear gradients
- class pyscf.grad.mp2.Gradients(method)[source]#
Bases:
GradientsBase
- 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 total MP2 energy.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the MP2 and the underlying SCF objects (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, mp >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> mp2_scanner = mp.MP2(scf.RHF(mol)).nuc_grad_method().as_scanner() >>> e_tot, grad = mp2_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot, grad = mp2_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
- grad_elec(t2, atmlst=None, verbose=4)#
- kernel(t2=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.
- class pyscf.grad.mp2.MP2_GradScanner(g)[source]#
Bases:
GradScanner
- property converged#
- pyscf.grad.mp2.as_scanner(grad_mp)[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 total MP2 energy.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the MP2 and the underlying SCF objects (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, mp >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> mp2_scanner = mp.MP2(scf.RHF(mol)).nuc_grad_method().as_scanner() >>> e_tot, grad = mp2_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot, grad = mp2_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
pyscf.grad.mspdft module#
- class pyscf.grad.mspdft.Gradients(mc)[source]#
Bases:
Gradients
- debug_lagrange(Lvec, bvec, Aop, Adiag, state=None, mo=None, ci=None, d2f=None, verbose=None, eris=None, **kwargs)[source]#
- get_Aop_Adiag(verbose=None, mo=None, ci=None, eris=None, level_shift=None, d2f=None, **kwargs)[source]#
Return a function calculating Lvec . J_wfn, where J_wfn is the Jacobian of the Lagrange cofactors (e.g., in state-averaged CASSCF, the Hessian of the state-averaged energy wrt wfn parameters) along with the diagonal of the Jacobian.
- get_LdotJnuc(Lvec, atmlst=None, verbose=None, mo=None, ci=None, eris=None, mf_grad=None, d2f=None, **kwargs)[source]#
Add the IS component
- get_ham_response(si_bra=None, si_ket=None, state=None, mo=None, ci=None, si=None, eris=None, veff1=None, veff2=None, mf_grad=None, atmlst=None, verbose=None, **kwargs)[source]#
write mspdft heff Hellmann-Feynman calculator; sum over diagonal PDFT Hellmann-Feynman terms
- get_init_guess(bvec, Adiag, Aop, precond)[source]#
Initial guess should solve the problem for SA-SA rotations
- get_wfn_response(si_bra=None, si_ket=None, state=None, mo=None, ci=None, si=None, eris=None, veff1=None, veff2=None, _freeze_is=False, d2f=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.
- kernel(state=None, mo=None, ci=None, si=None, _freeze_is=False, **kwargs)[source]#
Cache the Hamiltonian and effective Hamiltonian terms, and pass around the IS hessian
eris, veff1, veff2, and d2f should be available to all top-level functions: get_wfn_response, get_Aop_Adiag, get_ham_response, and get_LdotJnuc
freeze_is == True sets the is component of the response to zero for debugging purposes
- property nis#
- class pyscf.grad.mspdft.MSPDFTLagPrec(Adiag=None, level_shift=None, ci=None, grad_method=None, d2f=None, **kwargs)[source]#
Bases:
SACASLagPrec
Solve IS part exactly, then do everything else the same
- pyscf.grad.mspdft.get_diabfns(obj)[source]#
Interpret the name of the MS-PDFT method as a pair of functions which compute the derivatives of a particular objective function with respect to wave function parameters and geometry perturbations, excluding first and second derivatives wrt intermediate state rotations, which is handled by the energy-class version of this function.
- Args:
- objstring
Specify particular MS-PDFT method. Currently, only “CMS” is supported. Not case-sensitive.
- Returns:
- diab_responsecallable
Computes the orbital-rotation and CI-transfer sectors of the Hessian-vector product of the MS objective function for a vector of intermediate-state rotations
- diab_gradcallable
Computes the gradient of the MS objective function wrt geometry perturbation
- pyscf.grad.mspdft.make_rdm12_heff_offdiag(mc, ci, si_bra, si_ket)[source]#
Compute <bra|O|ket> - sum_i <i|O|i>, where O is the 1- and 2-RDM operator product, and |bra> and |ket> are both states spanning the vector space of |i>, which are multi-determinantal many-electron states in an active space.
- Args:
- mcobject of class CASCI or CASSCF
Only “ncas” and “nelecas” are used, to determine Hilbert of ci
- cindarray or list of length (nroots)
Contains CI vectors spanning a model space
- si_brandarray of shape (nroots)
Coefficients of ci elements for state |bra>
- si_ketndarray of shape (nroots)
Coefficients of ci elements for state |ket>
- Returns:
- casdm1ndarray of shape [ncas,]*2
Contains O = p’q case
- casdm2ndarray of shape [ncas,]*4
Contains O = p’q’sr case
pyscf.grad.rhf module#
Non-relativistic Hartree-Fock analytical nuclear gradients
- class pyscf.grad.rhf.Gradients(method)[source]#
Bases:
GradientsBase
Non-relativistic restricted Hartree-Fock gradients
- grad_elec(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)#
Electronic part of RHF/RKS gradients
- Args:
mf_grad : grad.rhf.Gradients or grad.rks.Gradients object
- class pyscf.grad.rhf.GradientsBase(method)[source]#
Bases:
StreamObject
Basic nuclear gradient functions for non-relativistic methods
- 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, grad >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> hf_scanner = scf.RHF(mol).apply(grad.RHF).as_scanner() >>> e_tot, grad = hf_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot, grad = hf_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
- 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_dispersion(disp=None, with_3body=None, verbose=None)#
gradient of DFTD3/DFTD4 dispersion correction
- get_jk(mol=None, dm=None, hermi=0, omega=None)[source]#
J = ((-nabla i) j| kl) D_lk K = ((-nabla i) j| kl) D_jk
- grad(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)#
- hcore_generator(mol=None)#
- kernel(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=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.).
- optimizer(solver='geometric')[source]#
Geometry optimization solver
- Kwargs:
solver (string) : geometry optimization solver, can be “geomeTRIC” (default) or “berny”.
- pyscf.grad.rhf.GradientsMixin#
alias of
GradientsBase
- class pyscf.grad.rhf.SCF_GradScanner(g)[source]#
Bases:
GradScanner
- pyscf.grad.rhf.as_scanner(mf_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, grad >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> hf_scanner = scf.RHF(mol).apply(grad.RHF).as_scanner() >>> e_tot, grad = hf_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot, grad = hf_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
- pyscf.grad.rhf.grad_elec(mf_grad, mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)[source]#
Electronic part of RHF/RKS gradients
- Args:
mf_grad : grad.rhf.Gradients or grad.rks.Gradients object
pyscf.grad.rks module#
Non-relativistic RKS analytical nuclear gradients
- class pyscf.grad.rks.Gradients(mf)[source]#
Bases:
Gradients
- 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_veff(mol=None, dm=None)#
First order derivative of DFT effective potential matrix (wrt electron coordinates)
- Args:
ks_grad : grad.uhf.Gradients or grad.uks.Gradients object
- grid_response = False#
- pyscf.grad.rks.get_nlc_vxc(ni, mol, grids, xc_code, dm, relativity=0, hermi=1, max_memory=2000, verbose=None)[source]#
- pyscf.grad.rks.get_nlc_vxc_full_response(ni, mol, grids, xc_code, dms, relativity=0, hermi=1, max_memory=2000, verbose=None)[source]#
Full NLC functional response including the response of the grids
- pyscf.grad.rks.get_veff(ks_grad, mol=None, dm=None)[source]#
First order derivative of DFT effective potential matrix (wrt electron coordinates)
- Args:
ks_grad : grad.uhf.Gradients or grad.uks.Gradients object
- pyscf.grad.rks.get_vxc(ni, mol, grids, xc_code, dms, relativity=0, hermi=1, max_memory=2000, verbose=None)[source]#
pyscf.grad.rohf module#
Non-relativistic ROHF analytical nuclear gradients
- class pyscf.grad.rohf.Gradients(method)[source]#
Bases:
Gradients
Non-relativistic restricted open-shell Hartree-Fock gradients
- get_veff(mol, dm)[source]#
First order derivative of HF potential matrix (wrt electron coordinates)
- Args:
mf_grad : grad.uhf.Gradients or grad.uks.Gradients object
- grad_elec(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)#
Electronic part of UHF/UKS gradients
- Args:
mf_grad : grad.uhf.Gradients or grad.uks.Gradients object
- make_rdm1e(mo_energy, mo_coeff, mo_occ)#
Energy weighted density matrix
pyscf.grad.roks module#
Non-relativistic ROKS analytical nuclear gradients
- class pyscf.grad.roks.Gradients(mf)[source]#
Bases:
Gradients
Non-relativistic ROHF gradients
- get_veff(mol=None, dm=None)#
First order derivative of DFT effective potential matrix (wrt electron coordinates)
- Args:
ks_grad : grad.uhf.Gradients or grad.uks.Gradients object
- grad_elec(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)#
Electronic part of UHF/UKS gradients
- Args:
mf_grad : grad.uhf.Gradients or grad.uks.Gradients object
- make_rdm1e(mo_energy, mo_coeff, mo_occ)#
Energy weighted density matrix
pyscf.grad.sacasscf module#
- class pyscf.grad.sacasscf.CASSCF_GradScanner(g, state)[source]#
Bases:
GradScanner
- property converged#
- class pyscf.grad.sacasscf.Gradients(mc, state=None)[source]#
Bases:
Gradients
- as_scanner(state=None)#
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'))
- get_Aop_Adiag(atmlst=None, state=None, verbose=None, mo=None, ci=None, eris=None, level_shift=None, **kwargs)[source]#
Return a function calculating Lvec . J_wfn, where J_wfn is the Jacobian of the Lagrange cofactors (e.g., in state-averaged CASSCF, the Hessian of the state-averaged energy wrt wfn parameters) along with the diagonal of the Jacobian.
- get_LdotJnuc(Lvec, state=None, atmlst=None, verbose=None, mo=None, ci=None, eris=None, mf_grad=None, **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.
- get_ham_response(state=None, atmlst=None, verbose=None, mo=None, ci=None, eris=None, mf_grad=None, **kwargs)[source]#
Return expectation values <dH/dx> where x is nuclear displacement. I.E., the gradient if the method were variational.
- get_wfn_response(atmlst=None, state=None, verbose=None, mo=None, ci=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.
- kernel(state=None, atmlst=None, verbose=None, mo=None, ci=None, eris=None, mf_grad=None, e_states=None, level_shift=None, **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.).
- 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
- make_fcasscf_sa(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 2). Note that penalty methods etc. must be removed, and that child methods such as MC-PDFT which do not reoptimize the orbitals also do not alter this function.
- Kwargs:
- casscf_attrdictionary
Extra attributes to apply to fcasscf. Just in case.
- fcisolver_attrdictionary
Extra attributes to apply to fcasscf.fcisolver. Just in case.
- Returns:
- fcasscfobject of
StateAverageMCSCFSolver
Set up to evaluate second derivatives of SA-CASSCF average energy in the absence of (i.e., spin) penalties.
- fcasscfobject of
- pyscf.grad.sacasscf.Lci_dot_dgci_dx(Lci, weights, mc, mo_coeff=None, ci=None, atmlst=None, mf_grad=None, eris=None, verbose=None)[source]#
Modification of single-state CASSCF electronic energy nuclear gradient to compute instead the CI Lagrange term nuclear gradient:
sum_IJ Lci_IJ d2_Ecas/d_lambda d_PIJ
This involves the effective density matrices ~D_pq = sum_I w_I<L_I|p’q|I> + c.c. ~d_pqrs = sum_I w_I<L_I|p’r’sq|I> + c.c. (NB: All-core terms ~D_ii, ~d_iijj = 0
However, active-core terms ~d_xyii, ~d_xiiy != 0)
- pyscf.grad.sacasscf.Lorb_dot_dgorb_dx(Lorb, mc, mo_coeff=None, ci=None, atmlst=None, mf_grad=None, eris=None, verbose=None)[source]#
Modification of single-state CASSCF electronic energy nuclear gradient to compute instead the orbital Lagrange term nuclear gradient:
sum_pq Lorb_pq d2_Ecas/d_lambda d_kpq
This involves the effective density matrices ~D_pq = L_pr*D_rq + L_qr*D_pr ~d_pqrs = L_pt*d_tqrs + L_rt*d_pqts + L_qt*d_ptrs + L_st*d_pqrt (NB: L_pq = -L_qp)
- class pyscf.grad.sacasscf.SACASLagPrec(Adiag=None, level_shift=None, ci=None, grad_method=None)[source]#
Bases:
LagPrec
A callable preconditioner for solving the Lagrange equations. Based on Mol. Phys. 99, 103 (2001). Attributes:
- nrootsinteger
Number of roots in the SA space
- nlaginteger
Number of Lagrange degrees of freedom
- ngorbinteger
Number of Lagrange degrees of freedom which are orbital rotations
- level_shiftfloat
numerical shift applied to CI rotation Hessian
- cindarray of shape (nroots, ndet or ncscf)
Ci vectors of the SA space
- Rorbndarray of shape (ngorb)
Diagonal inverse Hessian matrix for orbital rotations
- Rcindarray of shape (nroots, ndet or ncsf)
Diagonal inverse Hessian matrix for CI rotations including a level shift
- Rci_sandarray of shape (nroots (I), ndet or ncsf, nroots (K))
First two factors of the inverse diagonal CI Hessian projected into SA space: Rci(I)|J> <J|Rci(I)|K>^{-1} <K|Rci(I) note: right-hand bra and R_I factor not included due to storage considerations Make the operand’s matrix element with <K|Rci(I) before taking the dot product!
- pyscf.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.grad.tdrhf module#
- class pyscf.grad.tdrhf.Gradients(td)[source]#
Bases:
GradientsBase
- as_scanner(state=1)#
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, tdscf, grad >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> td_grad_scanner = scf.RHF(mol).apply(tdscf.TDA).nuc_grad_method().as_scanner() >>> e_tot, grad = td_grad_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot, grad = td_grad_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
- cphf_conv_tol = 1e-08#
- cphf_max_cycle = 20#
- grad_elec(xy, singlet, atmlst=None)[source]#
Electronic part of TDA, TDHF nuclear gradients
- Args:
td_grad : grad.tdrhf.Gradients or grad.tdrks.Gradients object.
- x_ya two-element list of numpy arrays
TDDFT X and Y amplitudes. If Y is set to 0, this function computes TDA energy gradients.
- kernel(xy=None, state=None, singlet=None, atmlst=None)[source]#
- Args:
- stateint
Excited state ID. state = 1 means the first excited state.
- 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.
- class pyscf.grad.tdrhf.TDSCF_GradScanner(g, state)[source]#
Bases:
GradScanner
- property converged#
- pyscf.grad.tdrhf.as_scanner(td_grad, state=1)[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, tdscf, grad >>> mol = gto.M(atom='H 0 0 0; F 0 0 1') >>> td_grad_scanner = scf.RHF(mol).apply(tdscf.TDA).nuc_grad_method().as_scanner() >>> e_tot, grad = td_grad_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1')) >>> e_tot, grad = td_grad_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
- pyscf.grad.tdrhf.grad_elec(td_grad, x_y, singlet=True, atmlst=None, max_memory=2000, verbose=4)[source]#
Electronic part of TDA, TDHF nuclear gradients
- Args:
td_grad : grad.tdrhf.Gradients or grad.tdrks.Gradients object.
- x_ya two-element list of numpy arrays
TDDFT X and Y amplitudes. If Y is set to 0, this function computes TDA energy gradients.
pyscf.grad.tdrks module#
- pyscf.grad.tdrks.grad_elec(td_grad, x_y, singlet=True, atmlst=None, max_memory=2000, verbose=4)[source]#
Electronic part of TDA, TDDFT nuclear gradients
- Args:
td_grad : grad.tdrhf.Gradients or grad.tdrks.Gradients object.
- x_ya two-element list of numpy arrays
TDDFT X and Y amplitudes. If Y is set to 0, this function computes TDA energy gradients.
pyscf.grad.tduhf module#
- pyscf.grad.tduhf.grad_elec(td_grad, x_y, atmlst=None, max_memory=2000, verbose=4)[source]#
Electronic part of TDA, TDHF nuclear gradients
- Args:
td_grad : grad.tduhf.Gradients or grad.tduks.Gradients object.
- x_ya two-element list of numpy arrays
TDDFT X and Y amplitudes. If Y is set to 0, this function computes TDA energy gradients.
pyscf.grad.tduks module#
- pyscf.grad.tduks.grad_elec(td_grad, x_y, atmlst=None, max_memory=2000, verbose=4)[source]#
Electronic part of TDA, TDDFT nuclear gradients
- Args:
td_grad : grad.tdrhf.Gradients or grad.tdrks.Gradients object.
- x_ya two-element list of numpy arrays
TDDFT X and Y amplitudes. If Y is set to 0, this function computes TDA energy gradients.
pyscf.grad.uccsd module#
UCCSD analytical nuclear gradients
pyscf.grad.uccsd_t module#
pyscf.grad.ucisd module#
UCISD analytical nuclear gradients
pyscf.grad.uhf module#
Non-relativistic unrestricted Hartree-Fock analytical nuclear gradients
- class pyscf.grad.uhf.Gradients(method)[source]#
Bases:
GradientsBase
Non-relativistic unrestricted Hartree-Fock gradients
- grad_elec(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)#
Electronic part of UHF/UKS gradients
- Args:
mf_grad : grad.uhf.Gradients or grad.uks.Gradients object
- pyscf.grad.uhf.get_veff(mf_grad, mol, dm)[source]#
First order derivative of HF potential matrix (wrt electron coordinates)
- Args:
mf_grad : grad.uhf.Gradients or grad.uks.Gradients object
pyscf.grad.uks module#
Non-relativistic UKS analytical nuclear gradients
- class pyscf.grad.uks.Gradients(mf)[source]#
Bases:
Gradients
- 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_veff(mol=None, dm=None)#
First order derivative of DFT effective potential matrix (wrt electron coordinates)
- Args:
ks_grad : grad.uhf.Gradients or grad.uks.Gradients object
- grid_response = False#
- pyscf.grad.uks.get_veff(ks_grad, mol=None, dm=None)[source]#
First order derivative of DFT effective potential matrix (wrt electron coordinates)
- Args:
ks_grad : grad.uhf.Gradients or grad.uks.Gradients object
pyscf.grad.ump2 module#
UMP2 analytical nuclear gradients
Module contents#
Analytical nuclear gradients#
Simple usage:
>>> from pyscf import gto, scf, grad
>>> mol = gto.M(atom='N 0 0 0; N 0 0 1', basis='ccpvdz')
>>> mf = scf.RHF(mol).run()
>>> grad.RHF(mf).kernel()