pyscf.x2c package#
Submodules#
pyscf.x2c.dft module#
X2C 2-component DFT methods
- class pyscf.x2c.dft.RKS(mol, xc='LDA,VWN')[source]#
Bases:
KohnShamDFT
,RHF
- class pyscf.x2c.dft.UKS(mol, xc='LDA,VWN')[source]#
Bases:
KohnShamDFT
,UHF
pyscf.x2c.newton_ah module#
Second order X2C-SCF solver
- class pyscf.x2c.newton_ah.SecondOrderX2CRHF(mf)[source]#
Bases:
_CIAH_SOSCF
- gen_g_hop(mo_coeff, mo_occ, fock_ao=None, h1e=None, with_symmetry=False)#
pyscf.x2c.sfx2c1e module#
1-electron Spin-free X2C approximation
- class pyscf.x2c.sfx2c1e.SFX2C1E_SCF(mf)[source]#
Bases:
_X2C_SCF
- Attributes for spin-free X2C:
with_x2c : X2C object
- dip_moment(mol=None, dm=None, unit='Debye', verbose=3, picture_change=True, **kwargs)[source]#
Dipole moment calculation with picture change correction
- Args:
mol: an instance of
Mole
dm : a 2D ndarrays density matrices- Kwarg:
picture_chang (bool) : Whether to compute the dipole moment with picture change correction.
- Return:
A list: the dipole moment on x, y and z component
- pyscf.x2c.sfx2c1e.SpinFreeX2C#
alias of
SpinFreeX2CHelper
- class pyscf.x2c.sfx2c1e.SpinFreeX2CHelper(mol)[source]#
Bases:
X2CHelperBase
1-component X2c (spin-free part only)
- picture_change(even_operator=(None, None), odd_operator=None)[source]#
Picture change for even_operator + odd_operator
even_operator has two terms at diagonal blocks [ v 0 ] [ 0 w ]
odd_operator has the term at off-diagonal blocks [ 0 p ] [ p^T 0 ]
v, w, and p can be strings (integral name) or matrices.
- pyscf.x2c.sfx2c1e.sfx2c(mf)#
Spin-free X2C. For the given SCF object, it updates the hcore constructor. All integrals are computed in the real spherical GTO basis.
- Args:
mf : an SCF object
- Returns:
An SCF object
Examples:
>>> mol = gto.M(atom='H 0 0 0; F 0 0 1', basis='ccpvdz', verbose=0) >>> mf = scf.RHF(mol).sfx2c1e() >>> mf.scf()
>>> import pyscf.x2c.sfx2c1e >>> mol.symmetry = 1 >>> mol.build(0, 0) >>> mf = pyscf.x2c.sfx2c1e.sfx2c1e(scf.UHF(mol)) >>> mf.scf()
- pyscf.x2c.sfx2c1e.sfx2c1e(mf)[source]#
Spin-free X2C. For the given SCF object, it updates the hcore constructor. All integrals are computed in the real spherical GTO basis.
- Args:
mf : an SCF object
- Returns:
An SCF object
Examples:
>>> mol = gto.M(atom='H 0 0 0; F 0 0 1', basis='ccpvdz', verbose=0) >>> mf = scf.RHF(mol).sfx2c1e() >>> mf.scf()
>>> import pyscf.x2c.sfx2c1e >>> mol.symmetry = 1 >>> mol.build(0, 0) >>> mf = pyscf.x2c.sfx2c1e.sfx2c1e(scf.UHF(mol)) >>> mf.scf()
pyscf.x2c.sfx2c1e_grad module#
Analytical nuclear gradients for 1-electron spin-free x2c method
Ref. JCP 135, 084114 (2011); DOI:10.1063/1.3624397
pyscf.x2c.sfx2c1e_hess module#
Analytical nuclear hessian for 1-electron spin-free x2c method
Ref. JCP 135, 244104 (2011); DOI:10.1063/1.3667202 JCTC 8, 2617 (2012); DOI:10.1021/ct300127e
pyscf.x2c.stability module#
Generate X2C-SCF response functions
- pyscf.x2c.stability.x2chf_stability(mf, verbose=None, return_status=False)[source]#
Stability analysis for X2C-HF/X2C-KS method.
- Args:
mf : DHF or DKS object
- Kwargs:
- return_status: bool
Whether to return stable_i and stable_e
- Returns:
If return_status is False (default), the return value includes a new set of orbitals, which are more close to the stable condition.
Else, another one boolean variable (indicating current status: stable or unstable) is returned.
pyscf.x2c.tdscf module#
- class pyscf.x2c.tdscf.TDBase(mf)[source]#
Bases:
TDBase
- analyze(verbose=None)#
- get_ab(mf=None)[source]#
A and B matrices for TDDFT response function.
A[i,a,j,b] = delta_{ab}delta_{ij}(E_a - E_i) + (ia||bj) B[i,a,j,b] = (ia||jb)
- get_nto(state=1, threshold=0.3, verbose=None)#
Natural transition orbital analysis.
The natural transition density matrix between ground state and excited state \(Tia = \langle \Psi_{ex} | i a^\dagger | \Psi_0 \rangle\) can be transformed to diagonal form through SVD \(T = O \sqrt{\lambda} V^\dagger\). O and V are occupied and virtual natural transition orbitals. The diagonal elements \(\lambda\) are the weights of the occupied-virtual orbital pair in the excitation.
Ref: Martin, R. L., JCP, 118, 4775-4777
Note in the TDHF/TDDFT calculations, the excitation part (X) is interpreted as the CIS coefficients and normalized to 1. The de-excitation part (Y) is ignored.
- Args:
tdobj : TDA, or TDHF, or TDDFT object
- stateint
Excited state ID. state = 1 means the first excited state. If state < 0, state ID is counted from the last excited state.
- Kwargs:
- thresholdfloat
Above which the NTO coefficients will be printed in the output.
- Returns:
A list (weights, NTOs). NTOs are natural orbitals represented in AO basis. The first N_occ NTOs are occupied NTOs and the rest are virtual NTOs.
- pyscf.x2c.tdscf.gen_tda_hop(mf, fock_ao=None)#
A x
- pyscf.x2c.tdscf.gen_tdhf_operation(mf, fock_ao=None)[source]#
Generate function to compute
[ A B][X] [-B -A][Y]
pyscf.x2c.x2c module#
X2C 2-component HF methods
- class pyscf.x2c.x2c.RHF(mol)[source]#
Bases:
SCF
- TDA(*args, **kwargs)#
- TDHF(*args, **kwargs)#
- gen_response(mo_coeff=None, mo_occ=None, with_j=True, hermi=0, max_memory=None)#
Generate a function to compute the product of X2C-HF response function and density matrices.
- 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.x2c.x2c.SCF(mol)[source]#
Bases:
SCF
The full X2C problem (scaler + soc terms) in j-adapted spinor basis
- analyze(verbose=None)[source]#
Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Diople moment.
- dip_moment(mol=None, dm=None, unit='Debye', verbose=3, picture_change=True, **kwargs)[source]#
Dipole moment calculation with picture change correction
- Args:
mol: an instance of
Mole
dm : a 2D ndarrays density matrices- Kwarg:
picture_change (bool) : Whether to compute the dipole moment with picture change correction.
- Return:
A list: the dipole moment on x, y and z component
- get_jk(mol=None, dm=None, hermi=1, with_j=True, with_k=True, omega=None)[source]#
Compute J, K matrices for all input density matrices
- Args:
mol : an instance of
Mole
- dmndarray or list of ndarrays
A density matrix or a list of density matrices
- Kwargs:
- hermiint
Whether J, K matrix is hermitian
0 : not hermitian and not symmetric1 : hermitian or symmetric2 : anti-hermitian- vhfopt :
A class which holds precomputed quantities to optimize the computation of J, K matrices
- with_jboolean
Whether to compute J matrices
- with_kboolean
Whether to compute K matrices
- omegafloat
Parameter of range-separated Coulomb operator: erf( omega * r12 ) / r12. If specified, integration are evaluated based on the long-range part of the range-separated Coulomb operator.
- Returns:
Depending on the given dm, the function returns one J and one K matrix, or a list of J matrices and a list of K matrices, corresponding to the input density matrices.
Examples:
>>> from pyscf import gto, scf >>> from pyscf.scf import _vhf >>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1') >>> dms = numpy.random.random((3,mol.nao_nr(),mol.nao_nr())) >>> j, k = scf.hf.get_jk(mol, dms, hermi=0) >>> print(j.shape) (3, 2, 2)
- get_occ(mo_energy=None, mo_coeff=None)[source]#
Label the occupancies for each orbital
- Kwargs:
- mo_energy1D ndarray
Obital energies
- mo_coeff2D ndarray
Obital coefficients
Examples:
>>> from pyscf import gto, scf >>> mol = gto.M(atom='H 0 0 0; F 0 0 1.1') >>> mf = scf.hf.SCF(mol) >>> energy = numpy.array([-10., -1., 1, -2., 0, -3]) >>> mf.get_occ(energy) array([2, 2, 0, 2, 2, 2])
- init_guess_by_atom(mol=None)[source]#
Generate initial guess density matrix from superposition of atomic HF density matrix. The atomic HF is occupancy averaged RHF
- Returns:
Density matrix, 2D ndarray
- init_guess_by_chkfile(chkfile=None, project=None)[source]#
Read the HF results from checkpoint file, then project it to the basis defined by
mol
- Kwargs:
- projectNone or bool
Whether to project chkfile’s orbitals to the new basis. Note when the geometry of the chkfile and the given molecule are very different, this projection can produce very poor initial guess. In PES scanning, it is recommended to switch off project.
If project is set to None, the projection is only applied when the basis sets of the chkfile’s molecule are different to the basis sets of the given molecule (regardless whether the geometry of the two molecules are different). Note the basis sets are considered to be different if the two molecules are derived from the same molecule with different ordering of atoms.
- Returns:
Density matrix, 2D ndarray
- make_rdm1(mo_coeff=None, mo_occ=None, **kwargs)#
One-particle density matrix in AO representation
- Args:
- mo_coeff2D ndarray
Orbital coefficients. Each column is one orbital.
- mo_occ1D ndarray
Occupancy
- Returns:
One-particle density matrix, 2D ndarray
- stability(internal=None, external=None, verbose=None, return_status=False)[source]#
X2C-HF/X2C-KS stability analysis.
See also pyscf.scf.stability.rhf_stability function.
- Kwargs:
- return_status: bool
Whether to return stable_i and stable_e
- Returns:
If return_status is False (default), the return value includes two set of orbitals, which are more close to the stable condition. The first corresponds to the internal stability and the second corresponds to the external stability.
Else, another two boolean variables (indicating current status: stable or unstable) are returned. The first corresponds to the internal stability and the second corresponds to the external stability.
- x2c()#
- class pyscf.x2c.x2c.SpinOrbitalX2CHelper(mol)[source]#
Bases:
X2CHelperBase
2-component X2c (including spin-free and spin-dependent terms) in the Gaussian type spin-orbital basis (as the spin-orbital basis in GHF)
- get_hcore(mol=None)[source]#
2-component X2c Foldy-Wouthuysen (FW) Hamiltonian (including spin-free and spin-dependent terms) in the j-adapted spinor basis.
- picture_change(even_operator=(None, None), odd_operator=None)[source]#
Picture change for even_operator + odd_operator
even_operator has two terms at diagonal blocks [ v 0 ] [ 0 w ]
odd_operator has the term at off-diagonal blocks [ 0 p ] [ p^T 0 ]
v, w, and p can be strings (integral name) or matrices.
- class pyscf.x2c.x2c.SpinorX2CHelper(mol)[source]#
Bases:
X2CHelperBase
2-component X2c (including spin-free and spin-dependent terms) in the j-adapted spinor basis.
- class pyscf.x2c.x2c.UHF(mol)[source]#
Bases:
SCF
- TDA(*args, **kwargs)#
- TDHF(*args, **kwargs)#
- gen_response(mo_coeff=None, mo_occ=None, with_j=True, hermi=0, max_memory=None)#
Generate a function to compute the product of X2C-HF response function and density matrices.
- 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.x2c.x2c.X2C#
alias of
SpinorX2CHelper
- class pyscf.x2c.x2c.X2C1E_GSCF(mf)[source]#
Bases:
_X2C_SCF
- Attributes for spin-orbital X2C:
with_x2c : X2C object
- dip_moment(mol=None, dm=None, unit='Debye', verbose=3, picture_change=True, **kwargs)[source]#
Dipole moment calculation with picture change correction
- Args:
mol: an instance of
Mole
dm : a 2D ndarrays density matrices- Kwarg:
picture_change (bool) : Whether to compute the dipole moment with picture change correction.
- Return:
A list: the dipole moment on x, y and z component
- class pyscf.x2c.x2c.X2CHelperBase(mol)[source]#
Bases:
StreamObject
2-component X2c (including spin-free and spin-dependent terms) in the j-adapted spinor basis.
- approx = '1e'#
- basis = None#
- get_hcore(mol=None)[source]#
2-component X2c Foldy-Wouthuysen (FW) Hamiltonian (including spin-free and spin-dependent terms) in the j-adapted spinor basis.
- picture_change(even_operator=(None, None), odd_operator=None)[source]#
Picture change for even_operator + odd_operator
even_operator has two terms at diagonal blocks [ v 0 ] [ 0 w ]
odd_operator has the term at off-diagonal blocks [ 0 p ] [ p^T 0 ]
v, w, and p can be strings (integral name) or matrices.
- xuncontract = True#
- pyscf.x2c.x2c.get_hcore(mol)[source]#
2-component X2c hcore Hamiltonian (including spin-free and spin-dependent terms) in the j-adapted spinor basis.
- pyscf.x2c.x2c.get_jk(mol, dm, hermi=1, mf_opt=None, with_j=True, with_k=True, omega=None)[source]#
non-relativistic J/K matrices (without SSO,SOO etc) in the j-adapted spinor basis.
- pyscf.x2c.x2c.init_guess_by_minao(mol)[source]#
Initial guess in terms of the overlap to minimal basis.
- pyscf.x2c.x2c.x2c1e_ghf(mf)[source]#
For the given GHF object, generate X2C-GSCF object in GHF spin-orbital basis. Note the orbital basis of X2C_GSCF is different to the X2C_RHF and X2C_UHF objects. X2C_RHF and X2C_UHF use spinor basis.
- Args:
mf : an GHF/GKS object
- Returns:
An GHF/GKS object
Examples:
>>> mol = pyscf.M(atom='H 0 0 0; F 0 0 1', basis='ccpvdz', verbose=0) >>> mf = scf.GHF(mol).x2c1e().run()