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

pyscf.grad.casci.Grad#

alias of Gradients

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'))
dump_flags(verbose=None)[source]#
grad_elec(mo_coeff=None, ci=None, atmlst=None, verbose=None)#
grad_nuc(mol=None, atmlst=None)[source]#
hcore_generator(mol=None)[source]#
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.casci.grad_elec(mc_grad, mo_coeff=None, ci=None, atmlst=None, verbose=None)[source]#

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

pyscf.grad.casscf.Grad#

alias of Gradients

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.casscf.grad_elec(mc_grad, mo_coeff=None, ci=None, atmlst=None, verbose=None)[source]#

pyscf.grad.ccsd module#

CCSD analytical nuclear gradients

class pyscf.grad.ccsd.CCSD_GradScanner(g)[source]#

Bases: GradScanner

property converged#
pyscf.grad.ccsd.Grad#

alias of Gradients

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)#
grad_nuc(mol=None, atmlst=None)[source]#
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.grad_elec(cc_grad, t1=None, t2=None, l1=None, l2=None, eris=None, atmlst=None, d1=None, d2=None, verbose=4)[source]#

pyscf.grad.ccsd_slow module#

RCCSD

Ref: JCP 90, 1752 (1989); DOI:10.1063/1.456069

pyscf.grad.ccsd_slow.index_frozen_active(cc)[source]#
pyscf.grad.ccsd_slow.kernel(cc, t1, t2, l1, l2, eris=None)[source]#

pyscf.grad.ccsd_t module#

class pyscf.grad.ccsd_t.Gradients(method)[source]#

Bases: Gradients

grad_elec(t1=None, t2=None, l1=None, l2=None, eris=None, atmlst=None, verbose=4)#
pyscf.grad.ccsd_t.grad_elec(cc_grad, t1=None, t2=None, l1=None, l2=None, eris=None, atmlst=None, verbose=4)[source]#

pyscf.grad.cisd module#

CISD analytical nuclear gradients

class pyscf.grad.cisd.CISD_GradScanner(g, state)[source]#

Bases: GradScanner

property converged#
pyscf.grad.cisd.Grad#

alias of Gradients

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'))
dump_flags(verbose=None)[source]#
grad_elec(civec=None, eris=None, atmlst=None, verbose=4)#
grad_nuc(mol=None, atmlst=None)[source]#
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.cisd.grad_elec(cigrad, civec=None, eris=None, atmlst=None, verbose=4)[source]#

pyscf.grad.dhf module#

Relativistic Dirac-Hartree-Fock

pyscf.grad.dhf.Grad#

alias of Gradients

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.

get_veff(mol, dm)[source]#
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.).

make_rdm1e(mo_energy=None, mo_coeff=None, mo_occ=None)[source]#
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

get_hcore(mol=None)[source]#
get_ovlp(mol=None)[source]#
hcore_generator(mol)[source]#
pyscf.grad.dhf.get_coulomb_hf(mol, dm, level='SSSS')[source]#

Dirac-Hartree-Fock Coulomb repulsion

pyscf.grad.dhf.get_hcore(mol)[source]#
pyscf.grad.dhf.get_ovlp(mol)[source]#
pyscf.grad.dhf.get_veff(mol, dm, level='SSSS')#

Dirac-Hartree-Fock Coulomb repulsion

pyscf.grad.dhf.grad_elec(mf_grad, mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)[source]#

pyscf.grad.dispersion module#

gradient of dispersion correction for HF and DFT

pyscf.grad.dispersion.get_dispersion(mf_grad, disp=None, with_3body=None, verbose=None)[source]#

gradient of DFTD3/DFTD4 dispersion correction

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#
debug_lagrange(Lvec, bvec, Aop, Adiag, **kwargs)[source]#
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_init_guess(bvec, Adiag, Aop, precond)[source]#
get_lagrange_callback(Lvec_last, itvec, geff_op)[source]#
get_lagrange_precond(Adiag, level_shift=None, **kwargs)[source]#
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.).

solve_lagrange(Lvec_guess=None, level_shift=None, **kwargs)[source]#
class pyscf.grad.lagrange.LagPrec(Adiag=None, level_shift=None, **kwargs)[source]#

Bases: object

A callable preconditioner for solving the Lagrange equations. Default is 1/(Adiagd+level_shift)

pyscf.grad.mp2 module#

MP2 analytical nuclear gradients

pyscf.grad.mp2.Grad#

alias of 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)#
grad_nuc(mol=None, atmlst=None)[source]#
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.mp2.grad_elec(mp_grad, t2, atmlst=None, verbose=4)[source]#

pyscf.grad.rhf module#

Non-relativistic Hartree-Fock analytical nuclear gradients

pyscf.grad.rhf.Grad#

alias of Gradients

class pyscf.grad.rhf.Gradients(method)[source]#

Bases: GradientsBase

Non-relativistic restricted Hartree-Fock gradients

get_veff(mol=None, dm=None)[source]#
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

make_rdm1e(mo_energy=None, mo_coeff=None, mo_occ=None)[source]#
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'))
dump_flags(verbose=None)[source]#
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_hcore(mol=None)[source]#
get_j(mol=None, dm=None, hermi=0, omega=None)[source]#
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

get_k(mol=None, dm=None, hermi=0, omega=None)[source]#
get_ovlp(mol=None)[source]#
get_veff(mol=None, dm=None)[source]#
grad(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)#
grad_elec()[source]#
grad_nuc(mol=None, atmlst=None)[source]#
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.).

make_rdm1e(mo_energy=None, mo_coeff=None, mo_occ=None)[source]#
optimizer(solver='geometric')[source]#

Geometry optimization solver

Kwargs:

solver (string) : geometry optimization solver, can be “geomeTRIC” (default) or “berny”.

reset(mol=None)[source]#
symmetrize(de, atmlst=None)[source]#

Symmetrize the gradients wrt the point group symmetry of the molecule.

to_gpu()[source]#
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.get_hcore(mol)[source]#

Part of the nuclear gradients of core Hamiltonian

pyscf.grad.rhf.get_jk(mol, dm)[source]#

J = ((-nabla i) j| kl) D_lk K = ((-nabla i) j| kl) D_jk

pyscf.grad.rhf.get_ovlp(mol)[source]#
pyscf.grad.rhf.get_veff(mf_grad, mol, dm)[source]#

NR Hartree-Fock Coulomb repulsion

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.rhf.grad_nuc(mol, atmlst=None)[source]#

Derivatives of nuclear repulsion energy wrt nuclear coordinates

pyscf.grad.rhf.hcore_generator(mf, mol=None)[source]#
pyscf.grad.rhf.make_rdm1e(mo_energy, mo_coeff, mo_occ)[source]#

Energy weighted density matrix

pyscf.grad.rhf.symmetrize(mol, de, atmlst=None)[source]#

Symmetrize the gradients wrt the point group symmetry of the molecule.

pyscf.grad.rks module#

Non-relativistic RKS analytical nuclear gradients

pyscf.grad.rks.Grad#

alias of Gradients

class pyscf.grad.rks.Gradients(mf)[source]#

Bases: Gradients

dump_flags(verbose=None)[source]#
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.rks.get_vxc_full_response(ni, mol, grids, xc_code, dms, relativity=0, hermi=1, max_memory=2000, verbose=None)[source]#

Full response including the response of the grids

pyscf.grad.rks.grids_noresponse_cc(grids)[source]#
pyscf.grad.rks.grids_response_cc(grids)[source]#

pyscf.grad.rohf module#

Non-relativistic ROHF analytical nuclear gradients

pyscf.grad.rohf.Grad#

alias of 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.rohf.make_rdm1e(mf_grad, mo_energy, mo_coeff, mo_occ)[source]#

Energy weighted density matrix

pyscf.grad.roks module#

Non-relativistic ROKS analytical nuclear gradients

pyscf.grad.roks.Grad#

alias of 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'))
debug_lagrange(Lvec, bvec, Aop, Adiag, state=None, mo=None, ci=None, **kwargs)[source]#
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_lagrange_callback(Lvec_last, itvec, geff_op)[source]#
get_lagrange_precond(Adiag, level_shift=None, ci=None, **kwargs)[source]#
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.

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.

pack_uniq_var(xorb, xci)[source]#
project_Aop(Aop, ci, state)[source]#

Wrap the Aop function to project out redundant degrees of freedom for the CI part. What’s redundant changes between SA-CASSCF and MC-PDFT so modify this part in child classes.

unpack_uniq_var(x)[source]#
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!

ci_prec(xci_spins)[source]#
orb_prec(xorb)[source]#
pack_uniq_var(xorb, xci)[source]#
unpack_uniq_var(x)[source]#
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#

pyscf.grad.tdrhf.Grad#

alias of Gradients

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#
dump_flags(verbose=None)[source]#
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.

grad_nuc(mol=None, atmlst=None)[source]#
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#

alias of Gradients

class pyscf.grad.tdrks.Gradients(td)[source]#

Bases: Gradients

grad_elec(xy, singlet, atmlst=None)[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.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#

alias of Gradients

class pyscf.grad.tduhf.Gradients(td)[source]#

Bases: Gradients

grad_elec(xy, singlet=None, atmlst=None)[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.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#

alias of Gradients

class pyscf.grad.tduks.Gradients(td)[source]#

Bases: Gradients

grad_elec(xy, singlet=None, atmlst=None)[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.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.Grad#

alias of Gradients

class pyscf.grad.uccsd.Gradients(method)[source]#

Bases: Gradients

grad_elec(t1=None, t2=None, l1=None, l2=None, eris=None, atmlst=None, d1=None, d2=None, verbose=4)#
pyscf.grad.uccsd.grad_elec(cc_grad, t1=None, t2=None, l1=None, l2=None, eris=None, atmlst=None, d1=None, d2=None, verbose=4)[source]#

pyscf.grad.uccsd_t module#

class pyscf.grad.uccsd_t.Gradients(method)[source]#

Bases: Gradients

grad_elec(t1=None, t2=None, l1=None, l2=None, eris=None, atmlst=None, verbose=4)#
pyscf.grad.uccsd_t.grad_elec(cc_grad, t1=None, t2=None, l1=None, l2=None, eris=None, atmlst=None, verbose=4)[source]#

pyscf.grad.ucisd module#

UCISD analytical nuclear gradients

pyscf.grad.ucisd.Grad#

alias of Gradients

class pyscf.grad.ucisd.Gradients(myci)[source]#

Bases: Gradients

grad_elec(civec=None, eris=None, atmlst=None, verbose=4)#
pyscf.grad.ucisd.grad_elec(cigrad, civec=None, eris=None, atmlst=None, verbose=4)[source]#

pyscf.grad.uhf module#

Non-relativistic unrestricted Hartree-Fock analytical nuclear gradients

pyscf.grad.uhf.Grad#

alias of Gradients

class pyscf.grad.uhf.Gradients(method)[source]#

Bases: GradientsBase

Non-relativistic unrestricted Hartree-Fock gradients

get_veff(mol=None, dm=None)[source]#
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=None, mo_coeff=None, mo_occ=None)[source]#
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.uhf.grad_elec(mf_grad, mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)[source]#

Electronic part of UHF/UKS gradients

Args:

mf_grad : grad.uhf.Gradients or grad.uks.Gradients object

pyscf.grad.uhf.make_rdm1e(mo_energy, mo_coeff, mo_occ)[source]#

Energy weighted density matrix

pyscf.grad.uks module#

Non-relativistic UKS analytical nuclear gradients

pyscf.grad.uks.Grad#

alias of Gradients

class pyscf.grad.uks.Gradients(mf)[source]#

Bases: Gradients

dump_flags(verbose=None)[source]#
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.uks.get_vxc(ni, mol, grids, xc_code, dms, relativity=0, hermi=1, max_memory=2000, verbose=None)[source]#
pyscf.grad.uks.get_vxc_full_response(ni, mol, grids, xc_code, dms, relativity=0, hermi=1, max_memory=2000, verbose=None)[source]#

Full response including the response of the grids

pyscf.grad.ump2 module#

UMP2 analytical nuclear gradients

pyscf.grad.ump2.Grad#

alias of Gradients

class pyscf.grad.ump2.Gradients(method)[source]#

Bases: Gradients

grad_elec(t2, atmlst=None, verbose=4)#
pyscf.grad.ump2.grad_elec(mp_grad, t2, atmlst=None, verbose=4)[source]#

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()