pyscf.nac package#

Submodules#

pyscf.nac.mspdft module#

class pyscf.nac.mspdft.NonAdiabaticCouplings(mc, state=None, mult_ediff=False, use_etfs=False)[source]#

Bases: Gradients

MS-PDFT non-adiabatic couplings (NACs) between states

kwargs/attributes:

statetuple of length 2

The NACs returned are <state[1]|d(state[0])/dR>. In other words, state = (ket, bra).

mult_edifflogical

If True, returns NACs multiplied by the energy difference. Useful near conical intersections to avoid numerical problems.

use_etfslogical

If True, use the ``electron translation factors’’ of Fatehi and Subotnik [JPCL 3, 2039 (2012)], which guarantee conservation of total electron + nuclear momentum when the nuclei are moving (i.e., in non-adiabatic molecular dynamics). This corresponds to omitting the ``model state contribution’’.

get_ham_response(**kwargs)[source]#

write mspdft heff Hellmann-Feynman calculator; sum over diagonal PDFT Hellmann-Feynman terms

get_wfn_response(si_bra=None, si_ket=None, state=None, si=None, verbose=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(*args, **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

nac_model(mo_coeff=None, ci=None, si=None, si_bra=None, si_ket=None, state=None, mf_grad=None, atmlst=None, **kwargs)[source]#
pyscf.nac.mspdft.nac_model(mc_grad, mo_coeff=None, ci=None, si_bra=None, si_ket=None, mf_grad=None, atmlst=None)[source]#

Compute the “model-state contribution” to the MS-PDFT NAC

pyscf.nac.sacasscf module#

class pyscf.nac.sacasscf.NonAdiabaticCouplings(mc, state=None, mult_ediff=False, use_etfs=False)[source]#

Bases: Gradients

SA-CASSCF non-adiabatic couplings (NACs) between states

kwargs/attributes:

statetuple of length 2

The NACs returned are <state[1]|d(state[0])/dR>. In other words, state = (ket, bra).

mult_edifflogical

If True, returns NACs multiplied by the energy difference. Useful near conical intersections to avoid numerical problems.

use_etfslogical

If True, use the ``electron translation factors’’ of Fatehi and Subotnik [JPCL 3, 2039 (2012)], which guarantee conservation of total electron + nuclear momentum when the nuclei are moving (i.e., in non-adiabatic molecular dynamics). This corresponds to omitting the so-called ``CSF contribution’’ [cf. JCTC 12, 3636 (2016)].

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(*args, **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_nacs(state=None, casscf_attr=None, fcisolver_attr=None)[source]#
nac_csf(mo_coeff=None, ci=None, state=None, mf_grad=None, atmlst=None)[source]#
pyscf.nac.sacasscf.gen_g_hop_active(mc, mo, ci0, eris, verbose=None)[source]#

Compute the active-electron part of the orbital rotation gradient by patching out the appropriate block of eris.vhf_c

pyscf.nac.sacasscf.grad_elec_active(mc_grad, mo_coeff=None, ci=None, atmlst=None, eris=None, mf_grad=None, verbose=None)[source]#

Compute the active-electron part of the CASSCF (Hellmann-Feynman) gradient by subtracting the core-electron part.

pyscf.nac.sacasscf.grad_elec_core(mc_grad, mo_coeff=None, atmlst=None, eris=None, mf_grad=None)[source]#

Compute the core-electron part of the CASSCF (Hellmann-Feynman) gradient using a modified RHF grad_elec call.

pyscf.nac.sacasscf.nac_csf(mc_grad, mo_coeff=None, ci=None, state=None, mf_grad=None, atmlst=None)[source]#

Compute the “CSF contribution” to the SA-CASSCF NAC

Module contents#

Analytical Nonadiabatic Coupling Vectors#

Simple usage:

>>> from pyscf import gto, scf, mcscf, nac
>>> mol = gto.M(atom='N 0 0 0; N 0 0 1', basis='ccpvdz')
>>> mf = scf.RHF(mol).run()
>>> mc = mcscf.CASSCF(mf, 2, 2).state_average([0.5, 0.5]).run()
>>> mc_nac = nac.sacasscf.NonAdiabaticCouplings(mc)
>>> mc_nac = mc.nac_method() # Also valid
>>> mc_nac.kernel(state=(0,1), use_etfs=False)