pyscf.pbc.scf package

Submodules

pyscf.pbc.scf.addons module

pyscf.pbc.scf.addons.canonical_occ(mf, nelec=None)

Label the occupancies for each orbital for sampled k-points. This is for KUHF objects. Each k-point has a fixed number of up and down electrons in this, which results in a finite size error for metallic systems but can accelerate convergence.

pyscf.pbc.scf.addons.canonical_occ_(mf, nelec=None)[source]

Label the occupancies for each orbital for sampled k-points. This is for KUHF objects. Each k-point has a fixed number of up and down electrons in this, which results in a finite size error for metallic systems but can accelerate convergence.

pyscf.pbc.scf.addons.convert_to_ghf(mf, out=None)[source]

Convert the given mean-field object to the generalized HF/KS object

Args:

mf : SCF object

Returns:

An generalized SCF object

pyscf.pbc.scf.addons.convert_to_khf(mf, out=None)

Convert gamma point SCF object to k-point SCF object

pyscf.pbc.scf.addons.convert_to_kscf(mf, out=None)[source]

Convert gamma point SCF object to k-point SCF object

pyscf.pbc.scf.addons.convert_to_rhf(mf, out=None)[source]

Convert the given mean-field object to the corresponding restricted HF/KS object

pyscf.pbc.scf.addons.convert_to_uhf(mf, out=None)[source]

Convert the given mean-field object to the corresponding unrestricted HF/KS object

pyscf.pbc.scf.addons.mo_energy_with_exxdiv_none(mf, mo_coeff=None)[source]

compute mo energy from the diagonal of fock matrix with exxdiv=None

pyscf.pbc.scf.addons.project_dm_k2k(cell, dm, kpts1, kpts2)[source]

Project density matrix from k-point mesh 1 to k-point mesh 2

pyscf.pbc.scf.addons.project_mo_nr2nr(cell1, mo1, cell2, kpts=None)[source]

Project orbital coefficients

\[ \begin{align}\begin{aligned}|\psi1> = |AO1> C1\\|\psi2> = P |\psi1> = |AO2>S^{-1}<AO2| AO1> C1 = |AO2> C2\\C2 = S^{-1}<AO2|AO1> C1\end{aligned}\end{align} \]
pyscf.pbc.scf.addons.smearing(mf, sigma=None, method='fermi', mu0=None, fix_spin=False)[source]

Fermi-Dirac or Gaussian smearing

pyscf.pbc.scf.addons.smearing_(mf, *args, **kwargs)[source]

pyscf.pbc.scf.chkfile module

pyscf.pbc.scf.chkfile.load_scf(chkfile)[source]

pyscf.pbc.scf.cphf module

Restricted coupled pertubed Hartree-Fock solver Modified from pyscf.scf.cphf

pyscf.pbc.scf.cphf.kernel(fvind, mo_energy, mo_occ, h1, s1=None, max_cycle=20, tol=1e-09, hermi=False, verbose=2)
Args:
fvindfunction

Given density matrix, compute (ij|kl)D_{lk}*2 - (ij|kl)D_{jk}

Kwargs:
hermiboolean

Whether the matrix defined by fvind is Hermitian or not.

pyscf.pbc.scf.cphf.solve(fvind, mo_energy, mo_occ, h1, s1=None, max_cycle=20, tol=1e-09, hermi=False, verbose=2)[source]
Args:
fvindfunction

Given density matrix, compute (ij|kl)D_{lk}*2 - (ij|kl)D_{jk}

Kwargs:
hermiboolean

Whether the matrix defined by fvind is Hermitian or not.

pyscf.pbc.scf.cphf.solve_nos1(fvind, mo_energy, mo_occ, h1, max_cycle=20, tol=1e-09, hermi=False, verbose=2)[source]

For field independent basis. First order overlap matrix is zero

pyscf.pbc.scf.cphf.solve_withs1(fvind, mo_energy, mo_occ, h1, s1, max_cycle=20, tol=1e-09, hermi=False, verbose=2)[source]

For field dependent basis. First order overlap matrix is non-zero. The first order orbitals are set to C^1_{ij} = -1/2 S1 e1 = h1 - s1*e0 + (e0_j-e0_i)*c1 + vhf[c1]

Kwargs:
hermiboolean

Whether the matrix defined by fvind is Hermitian or not.

Returns:

First order orbital coefficients (in MO basis) and first order orbital energy matrix

pyscf.pbc.scf.ghf module

Generalized Hartree-Fock for periodic systems at a single k-point

class pyscf.pbc.scf.ghf.GHF(cell, kpt=array([0., 0., 0.]), exxdiv='ewald')[source]

Bases: SCF, GHF

GHF class for PBCs.

CCSD(*args, **kwargs)

restricted CCSD

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

conv_tolfloat

converge threshold. Default is 1e-7.

conv_tol_normtfloat

converge threshold for norm(t1,t2). Default is 1e-5.

max_cycleint

max number of iterations. Default is 50.

diis_spaceint

DIIS space size. Default is 6.

diis_start_cycleint

The step to start DIIS. Default is 0.

iterative_dampingfloat

The self consistent damping parameter.

directbool

AO-direct CCSD. Default is False.

async_iobool

Allow for asynchronous function execution. Default is True.

incore_completebool

Avoid all I/O (also for DIIS). Default is False.

level_shiftfloat

A shift on virtual orbital energies to stablize the CCSD iteration

frozenint or list

If integer is given, the inner-most orbitals are frozen from CC amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in CC calculation.

>>> mol = gto.M(atom = 'H 0 0 0; F 0 0 1.1', basis = 'ccpvdz')
>>> mf = scf.RHF(mol).run()
>>> # freeze 2 core orbitals
>>> mycc = cc.CCSD(mf).set(frozen = 2).run()
>>> # auto-generate the number of core orbitals to be frozen (1 in this case)
>>> mycc = cc.CCSD(mf).set_frozen().run()
>>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals
>>> mycc.set(frozen = [0,1,16,17,18]).run()
callbackfunction(envs_dict) => None

callback function takes one dict as the argument which is generated by the builtin function locals(), so that the callback function can access all local variables in the current environment.

Saved results

convergedbool

CCSD converged or not

e_corrfloat

CCSD correlation correction

e_totfloat

Total CCSD energy (HF + correlation)

t1, t2 :

T amplitudes t1[i,a], t2[i,j,a,b] (i,j in occ, a,b in virt)

l1, l2 :

Lambda amplitudes l1[i,a], l2[i,j,a,b] (i,j in occ, a,b in virt)

CISD(*args, **kwargs)
MP2(*args, **kwargs)
analyze(*args, **kwargs)

Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Diople moment.

convert_from_(mf)[source]

Convert given mean-field object to GHF

gen_response(mo_coeff=None, mo_occ=None, with_j=True, hermi=0, max_memory=None)

Generate a function to compute the product of KGHF response function and KGHF density matrices.

get_bands(kpts_band, cell=None, dm=None, kpt=None)[source]

Get energy bands at the given (arbitrary) ‘band’ k-points.

Returns:
mo_energy(nmo,) ndarray or a list of (nmo,) ndarray

Bands energies E_n(k)

mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray

Band orbitals psi_n(k)

get_grad(mo_coeff, mo_occ, fock=None)[source]

RHF orbital gradients

Args:
mo_coeff2D ndarray

Obital coefficients

mo_occ1D ndarray

Orbital occupancy

fock_ao2D ndarray

Fock matrix in AO representation

Returns:

Gradients in MO representation. It’s a num_occ*num_vir vector.

get_hcore(cell=None, kpt=None)[source]
get_jk(cell=None, dm=None, hermi=0, kpt=None, kpts_band=None, with_j=True, with_k=True, **kwargs)

Get Coulomb (J) and exchange (K) following scf.hf.RHF.get_jk_(). for particular k-point (kpt).

When kpts_band is given, the J, K matrices on kpts_band are evaluated.

J_{pq} = sum_{rs} (pq|rs) dm[s,r] K_{pq} = sum_{rs} (pr|sq) dm[r,s]

where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.

get_occ(mo_energy=None, mo_coeff=None)

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])
get_ovlp(cell=None, kpt=None)[source]
get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)[source]

Hartree-Fock potential matrix for the given density matrix. See scf.hf.get_veff() and scf.hf.RHF.get_veff()

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

Returns:

Density matrix, 2D ndarray

init_guess_by_huckel(mol=None)[source]

Generate initial guess density matrix from a Huckel calculation based on occupancy averaged atomic RHF calculations, doi:10.1021/acs.jctc.8b01089

Returns:

Density matrix, 2D ndarray

init_guess_by_minao(mol=None)[source]

Generate initial guess density matrix based on ANO basis, then project the density matrix to the basis set defined by mol

Returns:

Density matrix, 2D ndarray

Examples:

>>> from pyscf import gto, scf
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1')
>>> scf.hf.init_guess_by_minao(mol)
array([[ 0.94758917,  0.09227308],
       [ 0.09227308,  0.94758917]])
mulliken_meta(*args, **kwargs)

Mulliken population analysis, based on meta-Lowdin AOs. In the meta-lowdin, the AOs are grouped in three sets: core, valence and Rydberg, the orthogonalization are carried out within each subsets.

Args:

mol : an instance of Mole

dmndarray or 2-item list of ndarray

Density matrix. ROHF dm is a 2-item list of 2D array

Kwargs:

verbose : int or instance of lib.logger.Logger

pre_orth_methodstr

Pre-orthogonalization, which localized GTOs for each atom. To obtain the occupied and unoccupied atomic shells, there are three methods

‘ano’ : Project GTOs to ANO basis
‘minao’ : Project GTOs to MINAO basis
‘scf’ : Symmetry-averaged fractional occupation atomic RHF
Returns:

A list : pop, charges

popnparray

Mulliken population on each atomic orbitals

chargesnparray

Mulliken charges

mulliken_pop(*args, **kwargs)

Mulliken population analysis

\[M_{ij} = D_{ij} S_{ji}\]

Mulliken charges

\[\delta_i = \sum_j M_{ij}\]
Returns:

A list : pop, charges

popnparray

Mulliken population on each atomic orbitals

chargesnparray

Mulliken charges

sfx2c1e()

X2C with spin-orbit coupling effects in spin-orbital basis

spin_square(mo_coeff=None, s=None)[source]

Spin of the GHF wavefunction

\[S^2 = \frac{1}{2}(S_+ S_- + S_- S_+) + S_z^2\]

where \(S_+ = \sum_i S_{i+}\) is effective for all beta occupied orbitals; \(S_- = \sum_i S_{i-}\) is effective for all alpha occupied orbitals.

  1. There are two possibilities for \(S_+ S_-\)
    1. same electron \(S_+ S_- = \sum_i s_{i+} s_{i-}\),

    \[\sum_i \langle UHF|s_{i+} s_{i-}|UHF\rangle = \sum_{pq}\langle p|s_+s_-|q\rangle \gamma_{qp} = n_\alpha\]

    2) different electrons \(S_+ S_- = \sum s_{i+} s_{j-}, (i\neq j)\). There are in total \(n(n-1)\) terms. As a two-particle operator,

    \[\langle S_+ S_- \rangle =\sum_{ij}(\langle i^\alpha|i^\beta\rangle \langle j^\beta|j^\alpha\rangle - \langle i^\alpha|j^\beta\rangle \langle j^\beta|i^\alpha\rangle)\]
  2. Similarly, for \(S_- S_+\)
    1. same electron

    \[\sum_i \langle s_{i-} s_{i+}\rangle = n_\beta\]
    1. different electrons

    \[\langle S_- S_+ \rangle =\sum_{ij}(\langle i^\beta|i^\alpha\rangle \langle j^\alpha|j^\beta\rangle - \langle i^\beta|j^\alpha\rangle \langle j^\alpha|i^\beta\rangle)\]
  3. For \(S_z^2\)
    1. same electron

    \[\langle s_z^2\rangle = \frac{1}{4}(n_\alpha + n_\beta)\]
    1. different electrons

    \[\begin{split}&\sum_{ij}(\langle ij|s_{z1}s_{z2}|ij\rangle -\langle ij|s_{z1}s_{z2}|ji\rangle) \\ &=\frac{1}{4}\sum_{ij}(\langle i^\alpha|i^\alpha\rangle \langle j^\alpha|j^\alpha\rangle - \langle i^\alpha|i^\alpha\rangle \langle j^\beta|j^\beta\rangle - \langle i^\beta|i^\beta\rangle \langle j^\alpha|j^\alpha\rangle + \langle i^\beta|i^\beta\rangle \langle j^\beta|j^\beta\rangle) \\ &-\frac{1}{4}\sum_{ij}(\langle i^\alpha|j^\alpha\rangle \langle j^\alpha|i^\alpha\rangle - \langle i^\alpha|j^\alpha\rangle \langle j^\beta|i^\beta\rangle - \langle i^\beta|j^\beta\rangle \langle j^\alpha|i^\alpha\rangle + \langle i^\beta|j^\beta\rangle\langle j^\beta|i^\beta\rangle) \\ &=\frac{1}{4}\sum_{ij}|\langle i^\alpha|i^\alpha\rangle - \langle i^\beta|i^\beta\rangle|^2 -\frac{1}{4}\sum_{ij}|\langle i^\alpha|j^\alpha\rangle - \langle i^\beta|j^\beta\rangle|^2 \\ &=\frac{1}{4}(n_\alpha - n_\beta)^2 -\frac{1}{4}\sum_{ij}|\langle i^\alpha|j^\alpha\rangle - \langle i^\beta|j^\beta\rangle|^2\end{split}\]
Args:
moa list of 2 ndarrays

Occupied alpha and occupied beta orbitals

Kwargs:
sndarray

AO overlap

Returns:

A list of two floats. The first is the expectation value of S^2. The second is the corresponding 2S+1

Examples:

>>> mol = gto.M(atom='O 0 0 0; H 0 0 1; H 0 1 0', basis='ccpvdz', charge=1, spin=1, verbose=0)
>>> mf = scf.UHF(mol)
>>> mf.kernel()
-75.623975516256706
>>> mo = (mf.mo_coeff[0][:,mf.mo_occ[0]>0], mf.mo_coeff[1][:,mf.mo_occ[1]>0])
>>> print('S^2 = %.7f, 2S+1 = %.7f' % spin_square(mo, mol.intor('int1e_ovlp_sph')))
S^2 = 0.7570150, 2S+1 = 2.0070027
stability(internal=None, external=None, verbose=None, return_status=False)[source]
to_ks(xc='HF')[source]

Convert to RKS object.

x2c()

X2C with spin-orbit coupling effects in spin-orbital basis

x2c1e()[source]

X2C with spin-orbit coupling effects in spin-orbital basis

pyscf.pbc.scf.ghf.get_jk(mf, cell=None, dm=None, hermi=0, kpt=None, kpts_band=None, with_j=True, with_k=True, **kwargs)[source]

pyscf.pbc.scf.hf module

Hartree-Fock for periodic systems at a single k-point

See Also:

pyscf.pbc.scf.khf.py : Hartree-Fock for periodic systems with k-point sampling

class pyscf.pbc.scf.hf.RHF(cell, kpt=array([0., 0., 0.]), exxdiv='ewald')[source]

Bases: SCF, RHF

CCSD(*args, **kwargs)

restricted CCSD

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

conv_tolfloat

converge threshold. Default is 1e-7.

conv_tol_normtfloat

converge threshold for norm(t1,t2). Default is 1e-5.

max_cycleint

max number of iterations. Default is 50.

diis_spaceint

DIIS space size. Default is 6.

diis_start_cycleint

The step to start DIIS. Default is 0.

iterative_dampingfloat

The self consistent damping parameter.

directbool

AO-direct CCSD. Default is False.

async_iobool

Allow for asynchronous function execution. Default is True.

incore_completebool

Avoid all I/O (also for DIIS). Default is False.

level_shiftfloat

A shift on virtual orbital energies to stablize the CCSD iteration

frozenint or list

If integer is given, the inner-most orbitals are frozen from CC amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in CC calculation.

>>> mol = gto.M(atom = 'H 0 0 0; F 0 0 1.1', basis = 'ccpvdz')
>>> mf = scf.RHF(mol).run()
>>> # freeze 2 core orbitals
>>> mycc = cc.CCSD(mf).set(frozen = 2).run()
>>> # auto-generate the number of core orbitals to be frozen (1 in this case)
>>> mycc = cc.CCSD(mf).set_frozen().run()
>>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals
>>> mycc.set(frozen = [0,1,16,17,18]).run()
callbackfunction(envs_dict) => None

callback function takes one dict as the argument which is generated by the builtin function locals(), so that the callback function can access all local variables in the current environment.

Saved results

convergedbool

CCSD converged or not

e_corrfloat

CCSD correlation correction

e_totfloat

Total CCSD energy (HF + correlation)

t1, t2 :

T amplitudes t1[i,a], t2[i,j,a,b] (i,j in occ, a,b in virt)

l1, l2 :

Lambda amplitudes l1[i,a], l2[i,j,a,b] (i,j in occ, a,b in virt)

CISD(*args, **kwargs)
MP2(*args, **kwargs)

restricted MP2 with canonical HF and non-canonical HF reference

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

conv_tolfloat

For non-canonical MP2, converge threshold for MP2 correlation energy. Default value is 1e-7.

conv_tol_normtfloat

For non-canonical MP2, converge threshold for norm(t2). Default value is 1e-5.

max_cycleint

For non-canonical MP2, max number of MP2 iterations. Default value is 50.

diis_spaceint

For non-canonical MP2, DIIS space size in MP2 iterations. Default is 6.

level_shiftfloat

A shift on virtual orbital energies to stablize the MP2 iterations.

frozenint or list

If integer is given, the inner-most orbitals are excluded from MP2 amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in MP2 calculation.

>>> mol = gto.M(atom = 'H 0 0 0; F 0 0 1.1', basis = 'ccpvdz')
>>> mf = scf.RHF(mol).run()
>>> # freeze 2 core orbitals
>>> pt = mp.MP2(mf).set(frozen = 2).run()
>>> # auto-generate the number of core orbitals to be frozen (1 in this case)
>>> pt = mp.MP2(mf).set_frozen().run()
>>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals
>>> pt.set(frozen = [0,1,16,17,18]).run()

Saved results

e_corrfloat

MP2 correlation correction

e_corr_ss/osfloat

Same-spin and opposite-spin component of the MP2 correlation energy

e_totfloat

Total MP2 energy (HF + correlation)

t2 :

T amplitudes t2[i,j,a,b] (i,j in occ, a,b in virt)

TDA(*args, **kwargs)

Tamm-Dancoff approximation

Attributes:
conv_tolfloat

Diagonalization convergence tolerance. Default is 1e-9.

nstatesint

Number of TD states to be computed. Default is 3.

Saved results:

convergedbool

Diagonalization converged or not

e1D array

excitation energy for each excited state.

xyA list of two 2D arrays

The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.

TDHF(*args, **kwargs)

Time-dependent Hartree-Fock

Attributes:
conv_tolfloat

Diagonalization convergence tolerance. Default is 1e-9.

nstatesint

Number of TD states to be computed. Default is 3.

Saved results:

convergedbool

Diagonalization converged or not

e1D array

excitation energy for each excited state.

xyA list of two 2D arrays

The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.

analyze(verbose=None, with_meta_lowdin=True, **kwargs)

Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Diople moment.

convert_from_(mf)[source]

Convert given mean-field object to RHF/ROHF

gen_response(mo_coeff=None, mo_occ=None, singlet=None, hermi=0, max_memory=None)

Generate a function to compute the product of RHF response function and RHF density matrices.

Kwargs:
singlet (None or boolean)If singlet is None, response function for

orbital hessian or CPHF will be generated. If singlet is boolean, it is used in TDDFT response kernel.

nuc_grad_method()[source]

Hook to create object for analytical nuclear gradients.

spin_square(mo_coeff=None, s=None)[source]

Spin square and multiplicity of RHF determinant

stability(internal=True, external=False, verbose=None, return_status=False)[source]

RHF/RKS stability analysis.

See also pyscf.scf.stability.rhf_stability function.

Kwargs:
internalbool

Internal stability, within the RHF optimization space.

externalbool

External stability. Including the RHF -> UHF and real -> complex stability analysis.

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.

to_ks(xc='HF')[source]

Convert to RKS object.

class pyscf.pbc.scf.hf.SCF(cell, kpt=array([0., 0., 0.]), exxdiv='ewald')[source]

Bases: SCF

SCF base class adapted for PBCs.

Attributes:
kpt(3,) ndarray

The AO k-point in Cartesian coordinates, in units of 1/Bohr.

exxdivstr

Exchange divergence treatment, can be one of

None : ignore G=0 contribution in exchange
‘ewald’ : Ewald probe charge correction [JCP 122, 234102 (2005); DOI:10.1063/1.1926272]
with_dfdensity fitting object

Default is the FFT based DF model. For all-electron calculation, MDF model is favored for better accuracy. See also pyscf.pbc.df.

analyze(verbose=None, with_meta_lowdin=None, **kwargs)[source]

Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Diople moment.

build(cell=None)[source]
check_sanity()[source]

Check input of class/object attributes, check whether a class method is overwritten. It does not check the attributes which are prefixed with “_”. The return value of method set is the object itself. This allows a series of functions/methods to be executed in pipe.

density_fit(auxbasis=None, with_df=None)[source]
dip_moment(cell=None, dm=None, unit='Debye', verbose=3, **kwargs)[source]

Dipole moment in the cell (is it well defined)?

Args:

cell : an instance of Cell

dm (ndarray) : density matrix

Return:

A list: the dipole moment on x, y and z components

dump_chk(envs)[source]
dump_flags(verbose=None)[source]
energy_nuc()[source]
from_chk(chk=None, project=None, kpt=None)[source]

Read the HF results from checkpoint file, then project it to the basis defined by mol

Returns:

Density matrix, 2D ndarray

get_bands(kpts_band, cell=None, dm=None, kpt=None)

Get energy bands at the given (arbitrary) ‘band’ k-points.

Returns:
mo_energy(nmo,) ndarray or a list of (nmo,) ndarray

Bands energies E_n(k)

mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray

Band orbitals psi_n(k)

get_hcore(cell=None, kpt=None)[source]
get_init_guess(cell=None, key='minao')[source]
get_j(cell=None, dm=None, hermi=1, kpt=None, kpts_band=None, omega=None)[source]

Compute J matrix for the given density matrix and k-point (kpt). When kpts_band is given, the J matrices on kpts_band are evaluated.

J_{pq} = sum_{rs} (pq|rs) dm[s,r]

where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.

get_jk(cell=None, dm=None, hermi=1, kpt=None, kpts_band=None, with_j=True, with_k=True, omega=None, **kwargs)[source]

Get Coulomb (J) and exchange (K) following scf.hf.RHF.get_jk_(). for particular k-point (kpt).

When kpts_band is given, the J, K matrices on kpts_band are evaluated.

J_{pq} = sum_{rs} (pq|rs) dm[s,r] K_{pq} = sum_{rs} (pr|sq) dm[r,s]

where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.

get_jk_incore(cell=None, dm=None, hermi=1, kpt=None, omega=None, **kwargs)[source]

Get Coulomb (J) and exchange (K) following scf.hf.RHF.get_jk_().

Incore version of Coulomb and exchange build only. Currently RHF always uses PBC AO integrals (unlike RKS), since exchange is currently computed by building PBC AO integrals.

get_k(cell=None, dm=None, hermi=1, kpt=None, kpts_band=None, omega=None)[source]

Compute K matrix for the given density matrix.

get_ovlp(cell=None, kpt=None)[source]
get_rho(dm=None, grids=None, kpt=None)

Compute density in real space

get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)[source]

Hartree-Fock potential matrix for the given density matrix. See scf.hf.get_veff() and scf.hf.RHF.get_veff()

init_direct_scf(*args, **kwargs)
init_guess_by_1e(cell=None)[source]

Generate initial guess density matrix from core hamiltonian

Returns:

Density matrix, 2D ndarray

init_guess_by_chkfile(chk=None, project=None, kpt=None)[source]

Read the HF results from checkpoint file, then project it to the basis defined by mol

Returns:

Density matrix, 2D ndarray

jk_method(J='FFTDF', K=None)[source]

Set up the schemes to evaluate Coulomb and exchange matrix

FFTDF: planewave density fitting using Fast Fourier Transform AFTDF: planewave density fitting using analytic Fourier Transform GDF: Gaussian density fitting MDF: Gaussian and planewave mix density fitting RS: range-separation JK builder RSDF: range-separation density fitting

property kpt
property kpts
mix_density_fit(auxbasis=None, with_df=None)[source]
property mol
mulliken_pop()[source]

Mulliken population analysis

\[M_{ij} = D_{ij} S_{ji}\]

Mulliken charges

\[\delta_i = \sum_j M_{ij}\]
Returns:

A list : pop, charges

popnparray

Mulliken population on each atomic orbitals

chargesnparray

Mulliken charges

reset(cell=None)[source]

Reset cell and relevant attributes associated to the old cell object

rs_density_fit(auxbasis=None, with_df=None)[source]
sfx2c1e()[source]
to_ghf()[source]

Convert the input mean-field object to a GHF/GKS object

to_kscf()[source]

Convert gamma point SCF object to k-point SCF object

to_rhf()[source]

Convert the input mean-field object to a RHF/ROHF/RKS/ROKS object

to_uhf()[source]

Convert the input mean-field object to a UHF/UKS object

x2c()
x2c1e()
pyscf.pbc.scf.hf.dip_moment(cell, dm, unit='Debye', verbose=3, grids=None, rho=None, kpt=array([0., 0., 0.]), origin=None)[source]

Dipole moment in the cell (is it well defined)?

Args:

cell : an instance of Cell

dm (ndarray) : density matrix

Return:

A list: the dipole moment on x, y and z components

pyscf.pbc.scf.hf.get_bands(mf, kpts_band, cell=None, dm=None, kpt=None)[source]

Get energy bands at the given (arbitrary) ‘band’ k-points.

Returns:
mo_energy(nmo,) ndarray or a list of (nmo,) ndarray

Bands energies E_n(k)

mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray

Band orbitals psi_n(k)

pyscf.pbc.scf.hf.get_hcore(cell, kpt=array([0., 0., 0.]))[source]

Get the core Hamiltonian AO matrix.

pyscf.pbc.scf.hf.get_j(cell, dm, hermi=1, vhfopt=None, kpt=array([0., 0., 0.]), kpts_band=None)[source]

Get the Coulomb (J) AO matrix for the given density matrix.

Args:
dmndarray or list of ndarrays

A density matrix or a list of density matrices

Kwargs:
hermiint

Whether J, K matrix is hermitian | 0 : no hermitian or symmetric | 1 : hermitian | 2 : anti-hermitian

vhfopt :

A class which holds precomputed quantities to optimize the computation of J, K matrices

kpt(3,) ndarray

The “inner” dummy k-point at which the DM was evaluated (or sampled).

kpts_band(3,) ndarray or (*,3) ndarray

An arbitrary “band” k-point at which J is evaluated.

Returns:

The function returns one J matrix, corresponding to the input density matrix (both order and shape).

pyscf.pbc.scf.hf.get_jk(mf, cell, dm, hermi=1, vhfopt=None, kpt=array([0., 0., 0.]), kpts_band=None, with_j=True, with_k=True, omega=None, **kwargs)[source]

Get the Coulomb (J) and exchange (K) AO matrices for the given density matrix.

Args:
dmndarray or list of ndarrays

A density matrix or a list of density matrices

Kwargs:
hermiint

Whether J, K matrix is hermitian | 0 : no hermitian or symmetric | 1 : hermitian | 2 : anti-hermitian

vhfopt :

A class which holds precomputed quantities to optimize the computation of J, K matrices

kpt(3,) ndarray

The “inner” dummy k-point at which the DM was evaluated (or sampled).

kpts_band(3,) ndarray or (*,3) ndarray

An arbitrary “band” k-point at which J and K are evaluated.

Returns:

The function returns one J and one K matrix, corresponding to the input density matrix (both order and shape).

pyscf.pbc.scf.hf.get_nuc(cell, kpt=array([0., 0., 0.]))[source]

Get the bare periodic nuc-el AO matrix, with G=0 removed.

See Martin (12.16)-(12.21).

pyscf.pbc.scf.hf.get_ovlp(cell, kpt=array([0., 0., 0.]))[source]

Get the overlap AO matrix.

pyscf.pbc.scf.hf.get_rho(mf, dm=None, grids=None, kpt=None)[source]

Compute density in real space

pyscf.pbc.scf.hf.get_t(cell, kpt=array([0., 0., 0.]))[source]

Get the kinetic energy AO matrix.

pyscf.pbc.scf.hf.init_guess_by_chkfile(cell, chkfile_name, project=None, kpt=None)[source]

Read the HF results from checkpoint file, then project it to the basis defined by cell

Returns:

Density matrix, (nao,nao) ndarray

pyscf.pbc.scf.hf.makov_payne_correction(mf)[source]

Makov-Payne correction (Phys. Rev. B, 51, 4014)

pyscf.pbc.scf.hf.normalize_dm_(mf, dm)[source]

Scale density matrix to make it produce the correct number of electrons.

pyscf.pbc.scf.kghf module

Generalized Hartree-Fock for periodic systems with k-point sampling

class pyscf.pbc.scf.kghf.KGHF(cell, kpts=array([[0., 0., 0.]]), exxdiv='ewald')[source]

Bases: KSCF, GHF

GHF class for PBCs.

CCSD(*args, **kwargs)

restricted CCSD

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

conv_tolfloat

converge threshold. Default is 1e-7.

conv_tol_normtfloat

converge threshold for norm(t1,t2). Default is 1e-5.

max_cycleint

max number of iterations. Default is 50.

diis_spaceint

DIIS space size. Default is 6.

diis_start_cycleint

The step to start DIIS. Default is 0.

iterative_dampingfloat

The self consistent damping parameter.

directbool

AO-direct CCSD. Default is False.

async_iobool

Allow for asynchronous function execution. Default is True.

incore_completebool

Avoid all I/O (also for DIIS). Default is False.

level_shiftfloat

A shift on virtual orbital energies to stablize the CCSD iteration

frozenint or list

If integer is given, the inner-most orbitals are frozen from CC amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in CC calculation.

>>> mol = gto.M(atom = 'H 0 0 0; F 0 0 1.1', basis = 'ccpvdz')
>>> mf = scf.RHF(mol).run()
>>> # freeze 2 core orbitals
>>> mycc = cc.CCSD(mf).set(frozen = 2).run()
>>> # auto-generate the number of core orbitals to be frozen (1 in this case)
>>> mycc = cc.CCSD(mf).set_frozen().run()
>>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals
>>> mycc.set(frozen = [0,1,16,17,18]).run()
callbackfunction(envs_dict) => None

callback function takes one dict as the argument which is generated by the builtin function locals(), so that the callback function can access all local variables in the current environment.

Saved results

convergedbool

CCSD converged or not

e_corrfloat

CCSD correlation correction

e_totfloat

Total CCSD energy (HF + correlation)

t1, t2 :

T amplitudes t1[i,a], t2[i,j,a,b] (i,j in occ, a,b in virt)

l1, l2 :

Lambda amplitudes l1[i,a], l2[i,j,a,b] (i,j in occ, a,b in virt)

MP2 = None
analyze(verbose=None, with_meta_lowdin=True, **kwargs)

Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Dipole moment

convert_from_(mf)

Convert given mean-field object to GHF

gen_response(mo_coeff=None, mo_occ=None, with_j=True, hermi=0, max_memory=None)

Generate a function to compute the product of KGHF response function and KGHF density matrices.

get_bands(kpts_band, cell=None, dm_kpts=None, kpts=None)[source]

Get energy bands at the given (arbitrary) ‘band’ k-points.

Returns:
mo_energy(nmo,) ndarray or a list of (nmo,) ndarray

Bands energies E_n(k)

mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray

Band orbitals psi_n(k)

get_grad(mo_coeff_kpts, mo_occ_kpts, fock=None)[source]

returns 1D array of gradients, like non K-pt version note that occ and virt indices of different k pts now occur in sequential patches of the 1D array

get_hcore(cell=None, kpts=None)[source]

Get the core Hamiltonian AO matrices at sampled k-points.

Args:

kpts : (nkpts, 3) ndarray

Returns:

hcore : (nkpts, nao, nao) ndarray

get_init_guess(cell=None, key='minao')
get_j(cell=None, dm_kpts=None, hermi=0, kpts=None, kpts_band=None)[source]

Compute J matrix for the given density matrix and k-point (kpt). When kpts_band is given, the J matrices on kpts_band are evaluated.

J_{pq} = sum_{rs} (pq|rs) dm[s,r]

where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.

get_jk(cell=None, dm_kpts=None, hermi=0, kpts=None, kpts_band=None, with_j=True, with_k=True, **kwargs)

Get Coulomb (J) and exchange (K) following scf.hf.RHF.get_jk_(). for particular k-point (kpt).

When kpts_band is given, the J, K matrices on kpts_band are evaluated.

J_{pq} = sum_{rs} (pq|rs) dm[s,r] K_{pq} = sum_{rs} (pr|sq) dm[r,s]

where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.

get_k(cell=None, dm_kpts=None, hermi=0, kpts=None, kpts_band=None)[source]

Compute K matrix for the given density matrix.

get_occ(mo_energy_kpts=None, mo_coeff_kpts=None)

Label the occupancies for each orbital for sampled k-points.

This is a k-point version of scf.hf.SCF.get_occ

get_ovlp(cell=None, kpts=None)[source]

Get the overlap AO matrices at sampled k-points.

Args:

kpts : (nkpts, 3) ndarray

Returns:

ovlp_kpts : (nkpts, nao, nao) ndarray

get_veff(cell=None, dm_kpts=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)[source]

Hartree-Fock potential matrix for the given density matrix. See scf.hf.get_veff() and scf.hf.RHF.get_veff()

init_guess_by_atom(cell=None, kpts=None)

Generates initial guess density matrix and the orbitals of the initial guess DM 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_name, project=None)

Read SCF chkfile and make the density matrix for GHF initial guess.

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.

init_guess_by_minao(cell=None, kpts=None)

Generates initial guess density matrix and the orbitals of the initial guess DM Generate initial guess density matrix based on ANO basis, then project the density matrix to the basis set defined by mol

Returns:

Density matrix, 2D ndarray

Examples:

>>> from pyscf import gto, scf
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1')
>>> scf.hf.init_guess_by_minao(mol)
array([[ 0.94758917,  0.09227308],
       [ 0.09227308,  0.94758917]])
mulliken_meta(cell=None, dm=None, verbose=5, pre_orth_method='ANO', s=None)[source]

A modified Mulliken population analysis, based on meta-Lowdin AOs.

Note this function only computes the Mulliken population for the gamma point density matrix.

mulliken_pop()[source]

Mulliken population analysis

\[M_{ij} = D_{ij} S_{ji}\]

Mulliken charges

\[\delta_i = \sum_j M_{ij}\]
Returns:

A list : pop, charges

popnparray

Mulliken population on each atomic orbitals

chargesnparray

Mulliken charges

to_ks(xc='HF')[source]

Convert to RKS object.

x2c()

X2C with spin-orbit coupling effects in spin-orbital basis

x2c1e()[source]

X2C with spin-orbit coupling effects in spin-orbital basis

pyscf.pbc.scf.kghf.get_jk(mf, cell=None, dm_kpts=None, hermi=0, kpts=None, kpts_band=None, with_j=True, with_k=True, **kwargs)[source]
pyscf.pbc.scf.kghf.get_occ(mf, mo_energy_kpts=None, mo_coeff_kpts=None)[source]

Label the occupancies for each orbital for sampled k-points.

This is a k-point version of scf.hf.SCF.get_occ

pyscf.pbc.scf.kghf.mulliken_meta(cell, dm_ao_kpts, verbose=5, pre_orth_method='ANO', s=None)[source]

A modified Mulliken population analysis, based on meta-Lowdin AOs.

Note this function only computes the Mulliken population for the gamma point density matrix.

pyscf.pbc.scf.kghf_ksymm module

pyscf.pbc.scf.kghf_ksymm.KGHF

alias of KsymAdaptedKGHF

class pyscf.pbc.scf.kghf_ksymm.KsymAdaptedKGHF(cell, kpts=<pyscf.pbc.lib.kpts.KPoints object>, exxdiv='ewald', use_ao_symmetry=True)[source]

Bases: KsymAdaptedKSCF, KGHF

KGHF with k-point symmetry

convert_from_(mf)

Convert given mean-field object to GHF

eig(h_kpts, s_kpts)[source]

Solver for generalized eigenvalue problem

\[HC = SCE\]
energy_elec(dm_kpts=None, h1e_kpts=None, vhf_kpts=None)

Following pyscf.scf.hf.energy_elec()

get_hcore(cell=None, kpts=None)[source]

Get the core Hamiltonian AO matrices at sampled k-points.

Args:

kpts : (nkpts, 3) ndarray

Returns:

hcore : (nkpts, nao, nao) ndarray

get_init_guess(cell=None, key='minao')
get_jk(cell=None, dm_kpts=None, hermi=0, kpts=None, kpts_band=None, with_j=True, with_k=True, **kwargs)

Get the Coulomb (J) and exchange (K) AO matrices at sampled k-points.

Args:
dm_kpts(nkpts, nao, nao) ndarray

Density matrix at each k-point. It needs to be Hermitian.

Kwargs:
kpts_band(3,) ndarray

A list of arbitrary “band” k-point at which to evalute the matrix.

Returns:

vj : (nkpts, nao, nao) ndarray vk : (nkpts, nao, nao) ndarray or list of vj and vk if the input dm_kpts is a list of DMs

get_occ(mo_energy_kpts=None, mo_coeff_kpts=None)

Label the occupancies for each orbital for sampled k-points.

This is a k-point version of scf.hf.SCF.get_occ

get_orbsym(mo_coeff=None, s=None)[source]
get_ovlp(cell=None, kpts=None)[source]

Get the overlap AO matrices at sampled k-points.

Args:

kpts : (nkpts, 3) ndarray

Returns:

ovlp_kpts : (nkpts, nao, nao) ndarray

init_guess_by_atom(cell=None, kpts=None)

Generates initial guess density matrix and the orbitals of the initial guess DM 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_name, project=None)

Read SCF chkfile and make the density matrix for GHF initial guess.

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.

init_guess_by_minao(cell=None, kpts=None)

Generates initial guess density matrix and the orbitals of the initial guess DM Generate initial guess density matrix based on ANO basis, then project the density matrix to the basis set defined by mol

Returns:

Density matrix, 2D ndarray

Examples:

>>> from pyscf import gto, scf
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1')
>>> scf.hf.init_guess_by_minao(mol)
array([[ 0.94758917,  0.09227308],
       [ 0.09227308,  0.94758917]])
property orbsym
to_ks(xc='HF')

Convert to RKS object.

pyscf.pbc.scf.kghf_ksymm.eig(kmf, h_kpts, s_kpts)[source]
pyscf.pbc.scf.kghf_ksymm.get_jk(mf, cell=None, dm_kpts=None, hermi=0, kpts=None, kpts_band=None, with_j=True, with_k=True, **kwargs)[source]
pyscf.pbc.scf.kghf_ksymm.get_occ(mf, mo_energy_kpts=None, mo_coeff_kpts=None)[source]

Label the occupancies for each orbital for sampled k-points.

This is a k-point version of scf.hf.SCF.get_occ

pyscf.pbc.scf.khf module

Hartree-Fock for periodic systems with k-point sampling

See Also:

hf.py : Hartree-Fock for periodic systems at a single k-point

class pyscf.pbc.scf.khf.KRHF(cell, kpts=array([[0., 0., 0.]]), exxdiv='ewald')[source]

Bases: KSCF, RHF

CCSD(*args, **kwargs)

restricted CCSD

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

conv_tolfloat

converge threshold. Default is 1e-7.

conv_tol_normtfloat

converge threshold for norm(t1,t2). Default is 1e-5.

max_cycleint

max number of iterations. Default is 50.

diis_spaceint

DIIS space size. Default is 6.

diis_start_cycleint

The step to start DIIS. Default is 0.

iterative_dampingfloat

The self consistent damping parameter.

directbool

AO-direct CCSD. Default is False.

async_iobool

Allow for asynchronous function execution. Default is True.

incore_completebool

Avoid all I/O (also for DIIS). Default is False.

level_shiftfloat

A shift on virtual orbital energies to stablize the CCSD iteration

frozenint or list

If integer is given, the inner-most orbitals are frozen from CC amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in CC calculation.

>>> mol = gto.M(atom = 'H 0 0 0; F 0 0 1.1', basis = 'ccpvdz')
>>> mf = scf.RHF(mol).run()
>>> # freeze 2 core orbitals
>>> mycc = cc.CCSD(mf).set(frozen = 2).run()
>>> # auto-generate the number of core orbitals to be frozen (1 in this case)
>>> mycc = cc.CCSD(mf).set_frozen().run()
>>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals
>>> mycc.set(frozen = [0,1,16,17,18]).run()
callbackfunction(envs_dict) => None

callback function takes one dict as the argument which is generated by the builtin function locals(), so that the callback function can access all local variables in the current environment.

Saved results

convergedbool

CCSD converged or not

e_corrfloat

CCSD correlation correction

e_totfloat

Total CCSD energy (HF + correlation)

t1, t2 :

T amplitudes t1[i,a], t2[i,j,a,b] (i,j in occ, a,b in virt)

l1, l2 :

Lambda amplitudes l1[i,a], l2[i,j,a,b] (i,j in occ, a,b in virt)

MP2(*args, **kwargs)

restricted MP2 with canonical HF and non-canonical HF reference

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

conv_tolfloat

For non-canonical MP2, converge threshold for MP2 correlation energy. Default value is 1e-7.

conv_tol_normtfloat

For non-canonical MP2, converge threshold for norm(t2). Default value is 1e-5.

max_cycleint

For non-canonical MP2, max number of MP2 iterations. Default value is 50.

diis_spaceint

For non-canonical MP2, DIIS space size in MP2 iterations. Default is 6.

level_shiftfloat

A shift on virtual orbital energies to stablize the MP2 iterations.

frozenint or list

If integer is given, the inner-most orbitals are excluded from MP2 amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in MP2 calculation.

>>> mol = gto.M(atom = 'H 0 0 0; F 0 0 1.1', basis = 'ccpvdz')
>>> mf = scf.RHF(mol).run()
>>> # freeze 2 core orbitals
>>> pt = mp.MP2(mf).set(frozen = 2).run()
>>> # auto-generate the number of core orbitals to be frozen (1 in this case)
>>> pt = mp.MP2(mf).set_frozen().run()
>>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals
>>> pt.set(frozen = [0,1,16,17,18]).run()

Saved results

e_corrfloat

MP2 correlation correction

e_corr_ss/osfloat

Same-spin and opposite-spin component of the MP2 correlation energy

e_totfloat

Total MP2 energy (HF + correlation)

t2 :

T amplitudes t2[i,j,a,b] (i,j in occ, a,b in virt)

TDA(*args, **kwargs)

Tamm-Dancoff approximation

Attributes:
conv_tolfloat

Diagonalization convergence tolerance. Default is 1e-9.

nstatesint

Number of TD states to be computed. Default is 3.

Saved results:

convergedbool

Diagonalization converged or not

e1D array

excitation energy for each excited state.

xyA list of two 2D arrays

The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.

TDHF(*args, **kwargs)

Time-dependent Hartree-Fock

Attributes:
conv_tolfloat

Diagonalization convergence tolerance. Default is 1e-9.

nstatesint

Number of TD states to be computed. Default is 3.

Saved results:

convergedbool

Diagonalization converged or not

e1D array

excitation energy for each excited state.

xyA list of two 2D arrays

The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.

analyze(verbose=None, with_meta_lowdin=True, **kwargs)

Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Dipole moment

check_sanity()[source]

Check input of class/object attributes, check whether a class method is overwritten. It does not check the attributes which are prefixed with “_”. The return value of method set is the object itself. This allows a series of functions/methods to be executed in pipe.

convert_from_(mf)[source]

Convert given mean-field object to KRHF/KROHF

gen_response(mo_coeff=None, mo_occ=None, singlet=None, hermi=0, max_memory=None)

Generate a function to compute the product of RHF response function and RHF density matrices.

Kwargs:
singlet (None or boolean)If singlet is None, response function for

orbital hessian or CPHF will be generated. If singlet is boolean, it is used in TDDFT response kernel.

get_init_guess(cell=None, key='minao')[source]
mulliken_meta(cell=None, dm=None, verbose=5, pre_orth_method='ANO', s=None)[source]

A modified Mulliken population analysis, based on meta-Lowdin AOs.

Note this function only computes the Mulliken population for the gamma point density matrix.

nuc_grad_method()[source]

Hook to create object for analytical nuclear gradients.

spin_square(mo_coeff=None, s=None)

Spin square and multiplicity of RHF determinant

stability(internal=True, external=False, verbose=None)[source]

RHF/RKS stability analysis.

See also pyscf.scf.stability.rhf_stability function.

Kwargs:
internalbool

Internal stability, within the RHF optimization space.

externalbool

External stability. Including the RHF -> UHF and real -> complex stability analysis.

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.

to_ks(xc='HF')[source]

Convert to RKS object.

class pyscf.pbc.scf.khf.KSCF(cell, kpts=array([[0., 0., 0.]]), exxdiv='ewald')[source]

Bases: SCF

SCF base class with k-point sampling.

Compared to molecular SCF, some members such as mo_coeff, mo_occ now have an additional first dimension for the k-points, e.g. mo_coeff is (nkpts, nao, nao) ndarray

Attributes:
kpts(nks,3) ndarray

The sampling k-points in Cartesian coordinates, in units of 1/Bohr.

analyze(verbose=None, with_meta_lowdin=True, **kwargs)[source]

Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Diople moment.

canonicalize(mo_coeff_kpts, mo_occ_kpts, fock=None)

Canonicalization diagonalizes the Fock matrix within occupied, open, virtual subspaces separatedly (without change occupancy).

check_sanity()

Check input of class/object attributes, check whether a class method is overwritten. It does not check the attributes which are prefixed with “_”. The return value of method set is the object itself. This allows a series of functions/methods to be executed in pipe.

conv_tol_grad = None
convert_from_(mf)[source]

Convert the abinput mean-field object to the associated KS object.

density_fit(auxbasis=None, with_df=None)[source]
dip_moment(cell=None, dm=None, unit='Debye', verbose=3, **kwargs)[source]

Dipole moment in the cell (is it well defined)?

Args:

cell : an instance of Cell

dm_kpts (a list of ndarrays) : density matrices of k-points

Return:

A list: the dipole moment on x, y and z components

dump_chk(envs)[source]
dump_flags(verbose=None)[source]
eig(h_kpts, s_kpts)[source]

Solver for generalized eigenvalue problem

\[HC = SCE\]
energy_elec(dm_kpts=None, h1e_kpts=None, vhf_kpts=None)

Following pyscf.scf.hf.energy_elec()

energy_nuc()
from_chk(chk=None, project=None, kpts=None)[source]

Read the HF results from checkpoint file, then project it to the basis defined by mol

Returns:

Density matrix, 2D ndarray

get_bands(kpts_band, cell=None, dm_kpts=None, kpts=None)[source]

Get energy bands at the given (arbitrary) ‘band’ k-points.

Returns:
mo_energy(nmo,) ndarray or a list of (nmo,) ndarray

Bands energies E_n(k)

mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray

Band orbitals psi_n(k)

get_fermi(mo_energy_kpts=None, mo_occ_kpts=None)

Fermi level

get_fock(h1e=None, s1e=None, vhf=None, dm=None, cycle=-1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)

F = h^{core} + V^{HF}

Special treatment (damping, DIIS, or level shift) will be applied to the Fock matrix if diis and cycle is specified (The two parameters are passed to get_fock function during the SCF iteration)

Kwargs:
h1e2D ndarray

Core hamiltonian

s1e2D ndarray

Overlap matrix, for DIIS

vhf2D ndarray

HF potential matrix

dm2D ndarray

Density matrix, for DIIS

cycleint

Then present SCF iteration step, for DIIS

diisan object of SCF.DIIS class

DIIS object to hold intermediate Fock and error vectors

diis_start_cycleint

The step to start DIIS. Default is 0.

level_shift_factorfloat or int

Level shift (in AU) for virtual space. Default is 0.

get_grad(mo_coeff_kpts, mo_occ_kpts, fock=None)[source]

returns 1D array of gradients, like non K-pt version note that occ and virt indices of different k pts now occur in sequential patches of the 1D array

get_hcore(cell=None, kpts=None)

Get the core Hamiltonian AO matrices at sampled k-points.

Args:

kpts : (nkpts, 3) ndarray

Returns:

hcore : (nkpts, nao, nao) ndarray

get_init_guess(cell=None, key='minao')[source]
get_j(cell=None, dm_kpts=None, hermi=1, kpts=None, kpts_band=None, omega=None)[source]

Compute J matrix for the given density matrix and k-point (kpt). When kpts_band is given, the J matrices on kpts_band are evaluated.

J_{pq} = sum_{rs} (pq|rs) dm[s,r]

where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.

get_jk(cell=None, dm_kpts=None, hermi=1, kpts=None, kpts_band=None, with_j=True, with_k=True, omega=None, **kwargs)[source]

Get Coulomb (J) and exchange (K) following scf.hf.RHF.get_jk_(). for particular k-point (kpt).

When kpts_band is given, the J, K matrices on kpts_band are evaluated.

J_{pq} = sum_{rs} (pq|rs) dm[s,r] K_{pq} = sum_{rs} (pr|sq) dm[r,s]

where r,s are orbitals on kpt. p and q are orbitals on kpts_band if kpts_band is given otherwise p and q are orbitals on kpt.

get_jk_incore(*args, **kwargs)

Get Coulomb (J) and exchange (K) following scf.hf.RHF.get_jk_().

Incore version of Coulomb and exchange build only. Currently RHF always uses PBC AO integrals (unlike RKS), since exchange is currently computed by building PBC AO integrals.

get_k(cell=None, dm_kpts=None, hermi=1, kpts=None, kpts_band=None, omega=None)[source]

Compute K matrix for the given density matrix.

get_occ(mo_energy_kpts=None, mo_coeff_kpts=None)

Label the occupancies for each orbital for sampled k-points.

This is a k-point version of scf.hf.SCF.get_occ

get_ovlp(cell=None, kpts=None)

Get the overlap AO matrices at sampled k-points.

Args:

kpts : (nkpts, 3) ndarray

Returns:

ovlp_kpts : (nkpts, nao, nao) ndarray

get_rho(dm=None, grids=None, kpts=None)

Compute density in real space

get_veff(cell=None, dm_kpts=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)[source]

Hartree-Fock potential matrix for the given density matrix. See scf.hf.get_veff() and scf.hf.RHF.get_veff()

init_direct_scf(*args, **kwargs)
init_guess_by_1e(cell=None)[source]

Generate initial guess density matrix from core hamiltonian

Returns:

Density matrix, 2D ndarray

init_guess_by_atom(cell=None, kpts=None)

Generates initial guess density matrix and the orbitals of the initial guess DM 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(chk=None, project=None, kpts=None)[source]

Read the HF results from checkpoint file, then project it to the basis defined by mol

Returns:

Density matrix, 2D ndarray

init_guess_by_minao(cell=None, kpts=None)

Generates initial guess density matrix and the orbitals of the initial guess DM Generate initial guess density matrix based on ANO basis, then project the density matrix to the basis set defined by mol

Returns:

Density matrix, 2D ndarray

Examples:

>>> from pyscf import gto, scf
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1')
>>> scf.hf.init_guess_by_minao(mol)
array([[ 0.94758917,  0.09227308],
       [ 0.09227308,  0.94758917]])
make_rdm1(mo_coeff_kpts=None, mo_occ_kpts=None, **kwargs)[source]

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

make_rdm2(mo_coeff_kpts, mo_occ_kpts, **kwargs)[source]

Two-particle density matrix in AO representation

Args:
mo_coeff2D ndarray

Orbital coefficients. Each column is one orbital.

mo_occ1D ndarray

Occupancy

Returns:

Two-particle density matrix, 4D ndarray

mix_density_fit(auxbasis=None, with_df=None)[source]
property mo_coeff_kpts
property mo_energy_kpts
property mo_occ_kpts
property mol
mulliken_meta(cell=None, dm=None, verbose=5, pre_orth_method='ANO', s=None)[source]

Mulliken population analysis, based on meta-Lowdin AOs. In the meta-lowdin, the AOs are grouped in three sets: core, valence and Rydberg, the orthogonalization are carried out within each subsets.

Args:

mol : an instance of Mole

dmndarray or 2-item list of ndarray

Density matrix. ROHF dm is a 2-item list of 2D array

Kwargs:

verbose : int or instance of lib.logger.Logger

pre_orth_methodstr

Pre-orthogonalization, which localized GTOs for each atom. To obtain the occupied and unoccupied atomic shells, there are three methods

‘ano’ : Project GTOs to ANO basis
‘minao’ : Project GTOs to MINAO basis
‘scf’ : Symmetry-averaged fractional occupation atomic RHF
Returns:

A list : pop, charges

popnparray

Mulliken population on each atomic orbitals

chargesnparray

Mulliken charges

mulliken_pop()[source]

Mulliken population analysis

\[M_{ij} = D_{ij} S_{ji}\]

Mulliken charges

\[\delta_i = \sum_j M_{ij}\]
Returns:

A list : pop, charges

popnparray

Mulliken population on each atomic orbitals

chargesnparray

Mulliken charges

newton()[source]

Create an SOSCF object based on the mean-field object

remove_soscf()[source]

Remove the SOSCF decorator

reset(cell=None)

Reset cell and relevant attributes associated to the old cell object

rs_density_fit(auxbasis=None, with_df=None)[source]
sfx2c1e()[source]
to_ghf()[source]

Convert the input mean-field object to a KGHF/KGKS object

to_khf()[source]

Disable point group symmetry

to_kscf()[source]

Convert to k-point SCF object

to_rhf()[source]

Convert the input mean-field object to a KRHF/KROHF/KRKS/KROKS object

to_uhf()[source]

Convert the input mean-field object to a KUHF/KUKS object

x2c()
x2c1e()
pyscf.pbc.scf.khf.analyze(mf, verbose=None, with_meta_lowdin=True, **kwargs)[source]

Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Dipole moment

pyscf.pbc.scf.khf.canonicalize(mf, mo_coeff_kpts, mo_occ_kpts, fock=None)[source]
pyscf.pbc.scf.khf.dip_moment(cell, dm_kpts, unit='Debye', verbose=3, grids=None, rho=None, kpts=array([[0., 0., 0.]]))[source]

Dipole moment in the cell (is it well defined)?

Args:

cell : an instance of Cell

dm_kpts (a list of ndarrays) : density matrices of k-points

Return:

A list: the dipole moment on x, y and z components

pyscf.pbc.scf.khf.energy_elec(mf, dm_kpts=None, h1e_kpts=None, vhf_kpts=None)[source]

Following pyscf.scf.hf.energy_elec()

pyscf.pbc.scf.khf.get_fermi(mf, mo_energy_kpts=None, mo_occ_kpts=None)[source]

Fermi level

pyscf.pbc.scf.khf.get_fock(mf, h1e=None, s1e=None, vhf=None, dm=None, cycle=-1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)[source]
pyscf.pbc.scf.khf.get_grad(mo_coeff_kpts, mo_occ_kpts, fock)[source]

returns 1D array of gradients, like non K-pt version note that occ and virt indices of different k pts now occur in sequential patches of the 1D array

pyscf.pbc.scf.khf.get_hcore(mf, cell=None, kpts=None)[source]

Get the core Hamiltonian AO matrices at sampled k-points.

Args:

kpts : (nkpts, 3) ndarray

Returns:

hcore : (nkpts, nao, nao) ndarray

pyscf.pbc.scf.khf.get_j(mf, cell, dm_kpts, kpts, kpts_band=None)[source]

Get the Coulomb (J) AO matrix at sampled k-points.

Args:
dm_kpts(nkpts, nao, nao) ndarray or a list of (nkpts,nao,nao) ndarray

Density matrix at each k-point. If a list of k-point DMs, eg, UHF alpha and beta DM, the alpha and beta DMs are contracted separately. It needs to be Hermitian.

Kwargs:
kpts_band(k,3) ndarray

A list of arbitrary “band” k-points at which to evalute the matrix.

Returns:

vj : (nkpts, nao, nao) ndarray or list of vj if the input dm_kpts is a list of DMs

pyscf.pbc.scf.khf.get_jk(mf, cell, dm_kpts, kpts, kpts_band=None, with_j=True, with_k=True, omega=None, **kwargs)[source]

Get the Coulomb (J) and exchange (K) AO matrices at sampled k-points.

Args:
dm_kpts(nkpts, nao, nao) ndarray

Density matrix at each k-point. It needs to be Hermitian.

Kwargs:
kpts_band(3,) ndarray

A list of arbitrary “band” k-point at which to evalute the matrix.

Returns:

vj : (nkpts, nao, nao) ndarray vk : (nkpts, nao, nao) ndarray or list of vj and vk if the input dm_kpts is a list of DMs

pyscf.pbc.scf.khf.get_occ(mf, mo_energy_kpts=None, mo_coeff_kpts=None)[source]

Label the occupancies for each orbital for sampled k-points.

This is a k-point version of scf.hf.SCF.get_occ

pyscf.pbc.scf.khf.get_ovlp(mf, cell=None, kpts=None)[source]

Get the overlap AO matrices at sampled k-points.

Args:

kpts : (nkpts, 3) ndarray

Returns:

ovlp_kpts : (nkpts, nao, nao) ndarray

pyscf.pbc.scf.khf.get_rho(mf, dm=None, grids=None, kpts=None)[source]

Compute density in real space

pyscf.pbc.scf.khf.init_guess_by_atom(cell, kpts=None)[source]

Generates initial guess density matrix and the orbitals of the initial guess DM based on the superposition of atomic HF density matrix.

pyscf.pbc.scf.khf.init_guess_by_chkfile(cell, chkfile_name, project=None, kpts=None)[source]

Read the KHF results from checkpoint file, then project it to the basis defined by cell

Returns:

Density matrix, 3D ndarray

pyscf.pbc.scf.khf.init_guess_by_minao(cell, kpts=None)[source]

Generates initial guess density matrix and the orbitals of the initial guess DM based on ANO basis.

pyscf.pbc.scf.khf.make_rdm1(mo_coeff_kpts, mo_occ_kpts, **kwargs)[source]

One particle density matrices for all k-points.

Returns:

dm_kpts : (nkpts, nao, nao) ndarray

pyscf.pbc.scf.khf.mulliken_meta(cell, dm_ao_kpts, verbose=5, pre_orth_method='ANO', s=None)[source]

A modified Mulliken population analysis, based on meta-Lowdin AOs.

Note this function only computes the Mulliken population for the gamma point density matrix.

pyscf.pbc.scf.khf_ksymm module

pyscf.pbc.scf.khf_ksymm.KRHF

alias of KsymAdaptedKRHF

class pyscf.pbc.scf.khf_ksymm.KsymAdaptedKRHF(cell, kpts=<pyscf.pbc.lib.kpts.KPoints object>, exxdiv='ewald', use_ao_symmetry=True)[source]

Bases: KsymAdaptedKSCF, KRHF

MP2(*args, **kwargs)

restricted MP2 with canonical HF and non-canonical HF reference

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

conv_tolfloat

For non-canonical MP2, converge threshold for MP2 correlation energy. Default value is 1e-7.

conv_tol_normtfloat

For non-canonical MP2, converge threshold for norm(t2). Default value is 1e-5.

max_cycleint

For non-canonical MP2, max number of MP2 iterations. Default value is 50.

diis_spaceint

For non-canonical MP2, DIIS space size in MP2 iterations. Default is 6.

level_shiftfloat

A shift on virtual orbital energies to stablize the MP2 iterations.

frozenint or list

If integer is given, the inner-most orbitals are excluded from MP2 amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in MP2 calculation.

>>> mol = gto.M(atom = 'H 0 0 0; F 0 0 1.1', basis = 'ccpvdz')
>>> mf = scf.RHF(mol).run()
>>> # freeze 2 core orbitals
>>> pt = mp.MP2(mf).set(frozen = 2).run()
>>> # auto-generate the number of core orbitals to be frozen (1 in this case)
>>> pt = mp.MP2(mf).set_frozen().run()
>>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals
>>> pt.set(frozen = [0,1,16,17,18]).run()

Saved results

e_corrfloat

MP2 correlation correction

e_corr_ss/osfloat

Same-spin and opposite-spin component of the MP2 correlation energy

e_totfloat

Total MP2 energy (HF + correlation)

t2 :

T amplitudes t2[i,j,a,b] (i,j in occ, a,b in virt)

convert_from_(mf)

Convert given mean-field object to KRHF/KROHF

get_init_guess(cell=None, key='minao')[source]
to_ks(xc='HF')

Convert to RKS object.

class pyscf.pbc.scf.khf_ksymm.KsymAdaptedKSCF(cell, kpts=<pyscf.pbc.lib.kpts.KPoints object>, exxdiv='ewald', use_ao_symmetry=True)[source]

Bases: KSCF

KRHF with k-point symmetry

dump_chk(envs)[source]
dump_flags(verbose=None)[source]
eig(h_kpts, s_kpts)[source]

Solver for generalized eigenvalue problem

\[HC = SCE\]
energy_elec(dm_kpts=None, h1e_kpts=None, vhf_kpts=None)

Following pyscf.scf.hf.energy_elec()

get_hcore(cell=None, kpts=None)[source]

Get the core Hamiltonian AO matrices at sampled k-points.

Args:

kpts : (nkpts, 3) ndarray

Returns:

hcore : (nkpts, nao, nao) ndarray

get_jk(cell=None, dm_kpts=None, hermi=1, kpts=None, kpts_band=None, with_j=True, with_k=True, omega=None, **kwargs)[source]

Get the Coulomb (J) and exchange (K) AO matrices at sampled k-points.

Args:
dm_kpts(nkpts, nao, nao) ndarray

Density matrix at each k-point. It needs to be Hermitian.

Kwargs:
kpts_band(3,) ndarray

A list of arbitrary “band” k-point at which to evalute the matrix.

Returns:

vj : (nkpts, nao, nao) ndarray vk : (nkpts, nao, nao) ndarray or list of vj and vk if the input dm_kpts is a list of DMs

get_occ(mo_energy_kpts=None, mo_coeff_kpts=None)

Label the occupancies for each orbital for sampled k-points.

This is a k-point version of scf.hf.SCF.get_occ

get_orbsym(mo_coeff=None, s=None)[source]
get_ovlp(cell=None, kpts=None)[source]

Get the overlap AO matrices at sampled k-points.

Args:

kpts : (nkpts, 3) ndarray

Returns:

ovlp_kpts : (nkpts, nao, nao) ndarray

get_rho(dm=None, grids=None, kpts=None)

Compute density in real space

init_guess_by_chkfile(chk=None, project=None, kpts=None)[source]

Read the HF results from checkpoint file, then project it to the basis defined by mol

Returns:

Density matrix, 2D ndarray

property kpts
property orbsym
sfx2c1e()[source]
to_khf()[source]

transform to non-symmetry object

x2c()
x2c1e()
pyscf.pbc.scf.khf_ksymm.eig(kmf, h_kpts, s_kpts)[source]
pyscf.pbc.scf.khf_ksymm.energy_elec(mf, dm_kpts=None, h1e_kpts=None, vhf_kpts=None)[source]

Following pyscf.scf.hf.energy_elec()

pyscf.pbc.scf.khf_ksymm.get_occ(mf, mo_energy_kpts=None, mo_coeff_kpts=None)[source]

Label the occupancies for each orbital for sampled k-points.

This is a k-point version of scf.hf.SCF.get_occ

pyscf.pbc.scf.khf_ksymm.get_rho(mf, dm=None, grids=None, kpts=None)[source]

Compute density in real space

pyscf.pbc.scf.khf_ksymm.ksymm_scf_common_init(kmf, cell, kpts, use_ao_symmetry=True)[source]

pyscf.pbc.scf.krohf module

Restricted open-shell Hartree-Fock for periodic systems with k-point sampling

class pyscf.pbc.scf.krohf.KROHF(cell, kpts=array([[0., 0., 0.]]), exxdiv='ewald')[source]

Bases: KRHF, ROHF

UHF class with k-point sampling.

CCSD = None
MP2 = None
TDA = None
TDHF = None
analyze(verbose=None, with_meta_lowdin=True, **kwargs)

Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Dipole moment

canonicalize(mo_coeff_kpts, mo_occ_kpts, fock=None)

Canonicalization diagonalizes the ROHF Fock matrix within occupied, virtual subspaces separatedly (without change occupancy).

conv_tol_grad = None
dump_flags(verbose=None)[source]
eig(fock, s)[source]

Solver for generalized eigenvalue problem

\[HC = SCE\]
energy_elec(dm_kpts=None, h1e_kpts=None, vhf_kpts=None)

Following pyscf.scf.hf.energy_elec()

gen_response(mo_coeff=None, mo_occ=None, with_j=True, hermi=0, max_memory=None)

Generate a function to compute the product of UHF response function and UHF density matrices.

get_fock(h1e=None, s1e=None, vhf=None, dm=None, cycle=-1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)

Build fock matrix based on Roothaan’s effective fock. See also get_roothaan_fock()

get_grad(mo_coeff_kpts, mo_occ_kpts, fock=None)[source]

returns 1D array of gradients, like non K-pt version note that occ and virt indices of different k pts now occur in sequential patches of the 1D array

get_init_guess(cell=None, key='minao')
get_occ(mo_energy_kpts=None, mo_coeff_kpts=None)

Label the occupancies for each orbital for sampled k-points.

This is a k-point version of scf.hf.SCF.get_occ

get_rho(dm=None, grids=None, kpts=None)

Compute density in real space

get_veff(cell=None, dm_kpts=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)[source]

Hartree-Fock potential matrix for the given density matrix. See scf.hf.get_veff() and scf.hf.RHF.get_veff()

init_guess_by_atom(mol=None)

Generates initial guess density matrix and the orbitals of the initial guess DM 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(chk=None, project=True, kpts=None)[source]

Read the HF results from checkpoint file, then project it to the basis defined by mol

Returns:

Density matrix, 2D ndarray

init_guess_by_huckel(mol=None)

Generate initial guess density matrix from a Huckel calculation based on occupancy averaged atomic RHF calculations, doi:10.1021/acs.jctc.8b01089

Returns:

Density matrix, 2D ndarray

init_guess_by_minao(mol=None)

Generates initial guess density matrix and the orbitals of the initial guess DM Generate initial guess density matrix based on ANO basis, then project the density matrix to the basis set defined by mol

Returns:

Density matrix, 2D ndarray

Examples:

>>> from pyscf import gto, scf
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1')
>>> scf.hf.init_guess_by_minao(mol)
array([[ 0.94758917,  0.09227308],
       [ 0.09227308,  0.94758917]])
init_guess_by_mod_huckel(mol=None)

Generate initial guess density matrix from a Huckel calculation based on occupancy averaged atomic RHF calculations, doi:10.1021/acs.jctc.8b01089

In contrast to init_guess_by_huckel, this routine employs the updated GWH rule from doi:10.1021/ja00480a005 to form the guess.

Returns:

Density matrix, 2D ndarray

make_rdm1(mo_coeff_kpts=None, mo_occ_kpts=None, **kwargs)[source]

One-particle density matrix. mo_occ is a 1D array, with occupancy 1 or 2.

mulliken_meta(cell=None, dm=None, verbose=5, pre_orth_method='ANO', s=None)[source]

A modified Mulliken population analysis, based on meta-Lowdin AOs.

Note this function only computes the Mulliken population for the gamma point density matrix.

property nelec
spin_square(mo_coeff=None, s=None)

Spin square and multiplicity of RHF determinant

stability(internal=True, external=False, verbose=None)[source]

ROHF/ROKS stability analysis.

See also pyscf.scf.stability.rohf_stability function.

Kwargs:
internalbool

Internal stability, within the RHF optimization space.

externalbool

External stability. It is not available in current version.

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.

to_ks(xc='HF')[source]

Convert to RKS object.

pyscf.pbc.scf.krohf.canonicalize(mf, mo_coeff_kpts, mo_occ_kpts, fock=None)[source]

Canonicalization diagonalizes the ROHF Fock matrix within occupied, virtual subspaces separatedly (without change occupancy).

pyscf.pbc.scf.krohf.get_fock(mf, h1e=None, s1e=None, vhf=None, dm=None, cycle=-1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)[source]
pyscf.pbc.scf.krohf.get_occ(mf, mo_energy_kpts=None, mo_coeff_kpts=None)[source]

Label the occupancies for each orbital for sampled k-points.

This is a k-point version of scf.hf.SCF.get_occ

pyscf.pbc.scf.krohf.get_roothaan_fock(focka_fockb, dma_dmb, s)[source]

Roothaan’s effective fock.

space

closed

open

virtual

closed

Fc

Fb

Fc

open

Fb

Fc

Fa

virtual

Fc

Fa

Fc

where Fc = (Fa + Fb) / 2

Returns:

Roothaan effective Fock matrix

pyscf.pbc.scf.krohf.make_rdm1(mo_coeff_kpts, mo_occ_kpts, **kwargs)[source]

Alpha and beta spin one particle density matrices for all k-points.

Returns:

dm_kpts : (2, nkpts, nao, nao) ndarray

pyscf.pbc.scf.krohf.mulliken_meta(cell, dm_ao_kpts, verbose=5, pre_orth_method='ANO', s=None)[source]

A modified Mulliken population analysis, based on meta-Lowdin AOs.

Note this function only computes the Mulliken population for the gamma point density matrix.

pyscf.pbc.scf.kuhf module

Hartree-Fock for periodic systems with k-point sampling

See Also:

hf.py : Hartree-Fock for periodic systems at a single k-point

class pyscf.pbc.scf.kuhf.KUHF(cell, kpts=array([[0., 0., 0.]]), exxdiv='ewald')[source]

Bases: KSCF, UHF

UHF class with k-point sampling.

CCSD(*args, **kwargs)

restricted CCSD

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

conv_tolfloat

converge threshold. Default is 1e-7.

conv_tol_normtfloat

converge threshold for norm(t1,t2). Default is 1e-5.

max_cycleint

max number of iterations. Default is 50.

diis_spaceint

DIIS space size. Default is 6.

diis_start_cycleint

The step to start DIIS. Default is 0.

iterative_dampingfloat

The self consistent damping parameter.

directbool

AO-direct CCSD. Default is False.

async_iobool

Allow for asynchronous function execution. Default is True.

incore_completebool

Avoid all I/O (also for DIIS). Default is False.

level_shiftfloat

A shift on virtual orbital energies to stablize the CCSD iteration

frozenint or list

If integer is given, the inner-most orbitals are frozen from CC amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in CC calculation.

>>> mol = gto.M(atom = 'H 0 0 0; F 0 0 1.1', basis = 'ccpvdz')
>>> mf = scf.RHF(mol).run()
>>> # freeze 2 core orbitals
>>> mycc = cc.CCSD(mf).set(frozen = 2).run()
>>> # auto-generate the number of core orbitals to be frozen (1 in this case)
>>> mycc = cc.CCSD(mf).set_frozen().run()
>>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals
>>> mycc.set(frozen = [0,1,16,17,18]).run()
callbackfunction(envs_dict) => None

callback function takes one dict as the argument which is generated by the builtin function locals(), so that the callback function can access all local variables in the current environment.

Saved results

convergedbool

CCSD converged or not

e_corrfloat

CCSD correlation correction

e_totfloat

Total CCSD energy (HF + correlation)

t1, t2 :

T amplitudes t1[i,a], t2[i,j,a,b] (i,j in occ, a,b in virt)

l1, l2 :

Lambda amplitudes l1[i,a], l2[i,j,a,b] (i,j in occ, a,b in virt)

MP2(*args, **kwargs)
TDA(*args, **kwargs)

Tamm-Dancoff approximation

Attributes:
conv_tolfloat

Diagonalization convergence tolerance. Default is 1e-9.

nstatesint

Number of TD states to be computed. Default is 3.

Saved results:

convergedbool

Diagonalization converged or not

e1D array

excitation energy for each excited state.

xyA list of two 2D arrays

The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.

TDHF(*args, **kwargs)
analyze(verbose=None, with_meta_lowdin=True, **kwargs)

Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Dipole moment

canonicalize(mo_coeff_kpts, mo_occ_kpts, fock=None)

Canonicalization diagonalizes the UHF Fock matrix within occupied, virtual subspaces separatedly (without change occupancy).

conv_tol_grad = None
convert_from_(mf)[source]

Convert given mean-field object to KUHF

dip_moment(cell=None, dm=None, unit='Debye', verbose=3, **kwargs)[source]

Dipole moment in the cell.

Args:

cell : an instance of Cell

dm_kpts (two lists of ndarrays) : KUHF density matrices of k-points

Return:

A list: the dipole moment on x, y and z components

dump_flags(verbose=None)[source]
eig(h_kpts, s_kpts)[source]

Solver for generalized eigenvalue problem

\[HC = SCE\]
energy_elec(dm_kpts=None, h1e_kpts=None, vhf_kpts=None)

Following pyscf.scf.hf.energy_elec()

gen_response(mo_coeff=None, mo_occ=None, with_j=True, hermi=0, max_memory=None)

Generate a function to compute the product of UHF response function and UHF density matrices.

get_bands(kpts_band, cell=None, dm_kpts=None, kpts=None)[source]

Get energy bands at the given (arbitrary) ‘band’ k-points.

Returns:
mo_energy(nmo,) ndarray or a list of (nmo,) ndarray

Bands energies E_n(k)

mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray

Band orbitals psi_n(k)

get_fermi(mo_energy_kpts=None, mo_occ_kpts=None)

A pair of Fermi level for spin-up and spin-down orbitals

get_fock(h1e=None, s1e=None, vhf=None, dm=None, cycle=-1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)

F = h^{core} + V^{HF}

Special treatment (damping, DIIS, or level shift) will be applied to the Fock matrix if diis and cycle is specified (The two parameters are passed to get_fock function during the SCF iteration)

Kwargs:
h1e2D ndarray

Core hamiltonian

s1e2D ndarray

Overlap matrix, for DIIS

vhf2D ndarray

HF potential matrix

dm2D ndarray

Density matrix, for DIIS

cycleint

Then present SCF iteration step, for DIIS

diisan object of SCF.DIIS class

DIIS object to hold intermediate Fock and error vectors

diis_start_cycleint

The step to start DIIS. Default is 0.

level_shift_factorfloat or int

Level shift (in AU) for virtual space. Default is 0.

get_grad(mo_coeff_kpts, mo_occ_kpts, fock=None)[source]

returns 1D array of gradients, like non K-pt version note that occ and virt indices of different k pts now occur in sequential patches of the 1D array

get_init_guess(cell=None, key='minao')[source]
get_occ(mo_energy_kpts=None, mo_coeff_kpts=None)

Label the occupancies for each orbital for sampled k-points.

This is a k-point version of scf.hf.SCF.get_occ

get_rho(dm=None, grids=None, kpts=None)

Compute density in real space

get_veff(cell=None, dm_kpts=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)[source]

Hartree-Fock potential matrix for the given density matrix. See scf.hf.get_veff() and scf.hf.RHF.get_veff()

init_guess_by_1e(cell=None)

Generate initial guess density matrix from core hamiltonian

Returns:

Density matrix, 2D ndarray

init_guess_by_atom(mol=None, breaksym=True)

Generates initial guess density matrix and the orbitals of the initial guess DM 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(chk=None, project=True, kpts=None)[source]

Read the HF results from checkpoint file, then project it to the basis defined by mol

Returns:

Density matrix, 2D ndarray

init_guess_by_huckel(mol=None, breaksym=True)

Generate initial guess density matrix from a Huckel calculation based on occupancy averaged atomic RHF calculations, doi:10.1021/acs.jctc.8b01089

Returns:

Density matrix, 2D ndarray

init_guess_by_minao(mol=None, breaksym=True)

Initial guess in terms of the overlap to minimal basis.

init_guess_by_mod_huckel(mol=None, breaksym=True)

Generate initial guess density matrix from a Huckel calculation based on occupancy averaged atomic RHF calculations, doi:10.1021/acs.jctc.8b01089

In contrast to init_guess_by_huckel, this routine employs the updated GWH rule from doi:10.1021/ja00480a005 to form the guess.

Returns:

Density matrix, 2D ndarray

make_rdm1(mo_coeff_kpts=None, mo_occ_kpts=None, **kwargs)[source]

One-particle density matrix in AO representation

Args:
mo_coefftuple of 2D ndarrays

Orbital coefficients for alpha and beta spins. Each column is one orbital.

mo_occtuple of 1D ndarrays

Occupancies for alpha and beta spins.

Returns:

A list of 2D ndarrays for alpha and beta spins

mulliken_meta(cell=None, dm=None, verbose=5, pre_orth_method='ANO', s=None)[source]

A modified Mulliken population analysis, based on meta-Lowdin AOs.

Note this function only computes the Mulliken population for the gamma point density matrix.

mulliken_pop()[source]

Mulliken population analysis

\[M_{ij} = D_{ij} S_{ji}\]

Mulliken charges

\[\delta_i = \sum_j M_{ij}\]
Returns:

A list : pop, charges

popnparray

Mulliken population on each atomic orbitals

chargesnparray

Mulliken charges

property nelec
nuc_grad_method()[source]

Hook to create object for analytical nuclear gradients.

spin_square(mo_coeff=None, s=None)[source]

Spin square and multiplicity of UHF determinant

\[S^2 = \frac{1}{2}(S_+ S_- + S_- S_+) + S_z^2\]

where \(S_+ = \sum_i S_{i+}\) is effective for all beta occupied orbitals; \(S_- = \sum_i S_{i-}\) is effective for all alpha occupied orbitals.

  1. There are two possibilities for \(S_+ S_-\)
    1. same electron \(S_+ S_- = \sum_i s_{i+} s_{i-}\),

    \[\sum_i \langle UHF|s_{i+} s_{i-}|UHF\rangle = \sum_{pq}\langle p|s_+s_-|q\rangle \gamma_{qp} = n_\alpha\]

    2) different electrons \(S_+ S_- = \sum s_{i+} s_{j-}, (i\neq j)\). There are in total \(n(n-1)\) terms. As a two-particle operator,

    \[\langle S_+ S_- \rangle = \langle ij|s_+ s_-|ij\rangle - \langle ij|s_+ s_-|ji\rangle = -\langle i^\alpha|j^\beta\rangle \langle j^\beta|i^\alpha\rangle\]
  2. Similarly, for \(S_- S_+\)
    1. same electron

    \[\sum_i \langle s_{i-} s_{i+}\rangle = n_\beta\]
    1. different electrons

    \[\langle S_- S_+ \rangle = -\langle i^\beta|j^\alpha\rangle \langle j^\alpha|i^\beta\rangle\]
  3. For \(S_z^2\)
    1. same electron

    \[\langle s_z^2\rangle = \frac{1}{4}(n_\alpha + n_\beta)\]
    1. different electrons

    \[\begin{split}&\frac{1}{2}\sum_{ij}(\langle ij|2s_{z1}s_{z2}|ij\rangle -\langle ij|2s_{z1}s_{z2}|ji\rangle) \\ &=\frac{1}{4}(\langle i^\alpha|i^\alpha\rangle \langle j^\alpha|j^\alpha\rangle - \langle i^\alpha|i^\alpha\rangle \langle j^\beta|j^\beta\rangle - \langle i^\beta|i^\beta\rangle \langle j^\alpha|j^\alpha\rangle + \langle i^\beta|i^\beta\rangle \langle j^\beta|j^\beta\rangle) \\ &-\frac{1}{4}(\langle i^\alpha|j^\alpha\rangle \langle j^\alpha|i^\alpha\rangle + \langle i^\beta|j^\beta\rangle\langle j^\beta|i^\beta\rangle) \\ &=\frac{1}{4}(n_\alpha^2 - n_\alpha n_\beta - n_\beta n_\alpha + n_\beta^2) -\frac{1}{4}(n_\alpha + n_\beta) \\ &=\frac{1}{4}((n_\alpha-n_\beta)^2 - (n_\alpha+n_\beta))\end{split}\]

In total

\[\begin{split}\langle S^2\rangle &= \frac{1}{2} (n_\alpha-\sum_{ij}\langle i^\alpha|j^\beta\rangle \langle j^\beta|i^\alpha\rangle +n_\beta -\sum_{ij}\langle i^\beta|j^\alpha\rangle\langle j^\alpha|i^\beta\rangle) + \frac{1}{4}(n_\alpha-n_\beta)^2 \\\end{split}\]
Args:
moa list of 2 ndarrays

Occupied alpha and occupied beta orbitals

Kwargs:
sndarray

AO overlap

Returns:

A list of two floats. The first is the expectation value of S^2. The second is the corresponding 2S+1

Examples:

>>> mol = gto.M(atom='O 0 0 0; H 0 0 1; H 0 1 0', basis='ccpvdz', charge=1, spin=1, verbose=0)
>>> mf = scf.UHF(mol)
>>> mf.kernel()
-75.623975516256706
>>> mo = (mf.mo_coeff[0][:,mf.mo_occ[0]>0], mf.mo_coeff[1][:,mf.mo_occ[1]>0])
>>> print('S^2 = %.7f, 2S+1 = %.7f' % spin_square(mo, mol.intor('int1e_ovlp_sph')))
S^2 = 0.7570150, 2S+1 = 2.0070027
stability(internal=True, external=False, verbose=None)[source]

Stability analysis for UHF/UKS method.

See also pyscf.scf.stability.uhf_stability function.

Args:

mf : UHF or UKS object

Kwargs:
internalbool

Internal stability, within the UHF space.

externalbool

External stability. Including the UHF -> GHF and real -> complex stability analysis.

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.

to_ks(xc='HF')[source]

Convert to RKS object.

pyscf.pbc.scf.kuhf.canonicalize(mf, mo_coeff_kpts, mo_occ_kpts, fock=None)[source]

Canonicalization diagonalizes the UHF Fock matrix within occupied, virtual subspaces separatedly (without change occupancy).

pyscf.pbc.scf.kuhf.dip_moment(cell, dm_kpts, unit='Debye', verbose=3, grids=None, rho=None, kpts=array([[0., 0., 0.]]))[source]

Dipole moment in the cell.

Args:

cell : an instance of Cell

dm_kpts (two lists of ndarrays) : KUHF density matrices of k-points

Return:

A list: the dipole moment on x, y and z components

pyscf.pbc.scf.kuhf.energy_elec(mf, dm_kpts=None, h1e_kpts=None, vhf_kpts=None)[source]

Following pyscf.scf.hf.energy_elec()

pyscf.pbc.scf.kuhf.get_fermi(mf, mo_energy_kpts=None, mo_occ_kpts=None)[source]

A pair of Fermi level for spin-up and spin-down orbitals

pyscf.pbc.scf.kuhf.get_fock(mf, h1e=None, s1e=None, vhf=None, dm=None, cycle=-1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)[source]
pyscf.pbc.scf.kuhf.get_occ(mf, mo_energy_kpts=None, mo_coeff_kpts=None)[source]

Label the occupancies for each orbital for sampled k-points.

This is a k-point version of scf.hf.SCF.get_occ

pyscf.pbc.scf.kuhf.init_guess_by_chkfile(cell, chkfile_name, project=None, kpts=None)[source]

Read the KHF results from checkpoint file, then project it to the basis defined by cell

Returns:

Density matrix, 3D ndarray

pyscf.pbc.scf.kuhf.make_rdm1(mo_coeff_kpts, mo_occ_kpts, **kwargs)[source]

Alpha and beta spin one particle density matrices for all k-points.

Returns:

dm_kpts : (2, nkpts, nao, nao) ndarray

pyscf.pbc.scf.kuhf.mulliken_meta(cell, dm_ao_kpts, verbose=5, pre_orth_method='ANO', s=None)[source]

A modified Mulliken population analysis, based on meta-Lowdin AOs.

Note this function only computes the Mulliken population for the gamma point density matrix.

pyscf.pbc.scf.kuhf_ksymm module

pyscf.pbc.scf.kuhf_ksymm.KUHF

alias of KsymAdaptedKUHF

class pyscf.pbc.scf.kuhf_ksymm.KsymAdaptedKUHF(cell, kpts=<pyscf.pbc.lib.kpts.KPoints object>, exxdiv='ewald', use_ao_symmetry=True)[source]

Bases: KsymAdaptedKSCF, KUHF

KUHF with k-point symmetry

convert_from_(mf)

Convert given mean-field object to KUHF

dump_flags(verbose=None)[source]
eig(h_kpts, s_kpts)[source]

Solver for generalized eigenvalue problem

\[HC = SCE\]
energy_elec(dm_kpts=None, h1e_kpts=None, vhf_kpts=None)

Following pyscf.scf.hf.energy_elec()

get_init_guess(cell=None, key='minao')[source]
get_occ(mo_energy_kpts=None, mo_coeff_kpts=None)

Label the occupancies for each orbital for sampled k-points.

This is a k-point version of scf.hf.SCF.get_occ

get_orbsym(mo_coeff=None, s=None)[source]
get_rho(dm=None, grids=None, kpts=None)

Compute density in real space

property nelec
property orbsym
to_ks(xc='HF')

Convert to RKS object.

pyscf.pbc.scf.kuhf_ksymm.energy_elec(mf, dm_kpts=None, h1e_kpts=None, vhf_kpts=None)[source]

Following pyscf.scf.hf.energy_elec()

pyscf.pbc.scf.kuhf_ksymm.get_occ(mf, mo_energy_kpts=None, mo_coeff_kpts=None)[source]

Label the occupancies for each orbital for sampled k-points.

This is a k-point version of scf.hf.SCF.get_occ

pyscf.pbc.scf.newton_ah module

Co-iterative augmented hessian second order SCF solver (CIAH-SOSCF)

pyscf.pbc.scf.newton_ah.gen_g_hop_rhf(mf, mo_coeff, mo_occ, fock_ao=None, h1e=None)[source]
pyscf.pbc.scf.newton_ah.gen_g_hop_rohf(mf, mo_coeff, mo_occ, fock_ao=None, h1e=None)[source]
pyscf.pbc.scf.newton_ah.gen_g_hop_uhf(mf, mo_coeff, mo_occ, fock_ao=None, h1e=None)[source]
pyscf.pbc.scf.newton_ah.newton(mf)[source]

pyscf.pbc.scf.rohf module

Restricted open-shell Hartree-Fock for periodic systems at a single k-point

See Also:

pyscf/pbc/scf/khf.py : Hartree-Fock for periodic systems with k-point sampling

class pyscf.pbc.scf.rohf.ROHF(cell, kpt=array([0., 0., 0.]), exxdiv='ewald')[source]

Bases: RHF, ROHF

ROHF class for PBCs.

CCSD = None
CISD = None
MP2 = None
TDA = None
TDHF = None
analyze(verbose=None, with_meta_lowdin=True, **kwargs)[source]

Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis

canonicalize(mo_coeff, mo_occ, fock=None)

Canonicalization diagonalizes the Fock matrix within occupied, open, virtual subspaces separatedly (without change occupancy).

dip_moment(cell=None, dm=None, unit='Debye', verbose=3, **kwargs)

Dipole moment in the cell (is it well defined)?

Args:

cell : an instance of Cell

dm (ndarray) : density matrix

Return:

A list: the dipole moment on x, y and z components

dump_flags(verbose=None)[source]
eig(fock, s)[source]

Solver for generalized eigenvalue problem

\[HC = SCE\]
energy_elec(dm=None, h1e=None, vhf=None)

Electronic part of Hartree-Fock energy, for given core hamiltonian and HF potential

… math:

E = \sum_{ij}h_{ij} \gamma_{ji}
  + \frac{1}{2}\sum_{ijkl} \gamma_{ji}\gamma_{lk} \langle ik||jl\rangle

Note this function has side effects which cause mf.scf_summary updated.

Args:

mf : an instance of SCF class

Kwargs:
dm2D ndarray

one-partical density matrix

h1e2D ndarray

Core hamiltonian

vhf2D ndarray

HF potential

Returns:

Hartree-Fock electronic energy and the Coulomb energy

Examples:

>>> from pyscf import gto, scf
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1')
>>> mf = scf.RHF(mol)
>>> mf.scf()
>>> dm = mf.make_rdm1()
>>> scf.hf.energy_elec(mf, dm)
(-1.5176090667746334, 0.60917167853723675)
>>> mf.energy_elec(dm)
(-1.5176090667746334, 0.60917167853723675)
gen_response(mo_coeff=None, mo_occ=None, with_j=True, hermi=0, max_memory=None)

Generate a function to compute the product of UHF response function and UHF density matrices.

get_bands(kpts_band, cell=None, dm=None, kpt=None)[source]

Get energy bands at the given (arbitrary) ‘band’ k-points.

Returns:
mo_energy(nmo,) ndarray or a list of (nmo,) ndarray

Bands energies E_n(k)

mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray

Band orbitals psi_n(k)

get_fock(h1e=None, s1e=None, vhf=None, dm=None, cycle=-1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)

Build fock matrix based on Roothaan’s effective fock. See also get_roothaan_fock()

get_grad(mo_coeff, mo_occ, fock=None)[source]

ROHF gradients is the off-diagonal block [co + cv + ov], where [ cc co cv ] [ oc oo ov ] [ vc vo vv ]

get_init_guess(cell=None, key='minao')
get_occ(mo_energy=None, mo_coeff=None)

Label the occupancies for each orbital. NOTE the occupancies are not assigned based on the orbital energy ordering. The first N orbitals are assigned to be occupied orbitals.

Examples:

>>> mol = gto.M(atom='H 0 0 0; O 0 0 1.1', spin=1)
>>> mf = scf.hf.SCF(mol)
>>> energy = numpy.array([-10., -1., 1, -2., 0, -3])
>>> mf.get_occ(energy)
array([2, 2, 2, 2, 1, 0])
get_rho(dm=None, grids=None, kpt=None)

Compute density in real space

get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)[source]

Hartree-Fock potential matrix for the given density matrix. See scf.hf.get_veff() and scf.hf.RHF.get_veff()

init_guess_by_1e(cell=None)[source]

Generate initial guess density matrix from core hamiltonian

Returns:

Density matrix, 2D ndarray

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(chk=None, project=True, kpt=None)

Read the HF results from checkpoint file, then project it to the basis defined by mol

Returns:

Density matrix, 2D ndarray

init_guess_by_huckel(mol=None)[source]

Generate initial guess density matrix from a Huckel calculation based on occupancy averaged atomic RHF calculations, doi:10.1021/acs.jctc.8b01089

Returns:

Density matrix, 2D ndarray

init_guess_by_minao(mol=None)[source]

Generate initial guess density matrix based on ANO basis, then project the density matrix to the basis set defined by mol

Returns:

Density matrix, 2D ndarray

Examples:

>>> from pyscf import gto, scf
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1')
>>> scf.hf.init_guess_by_minao(mol)
array([[ 0.94758917,  0.09227308],
       [ 0.09227308,  0.94758917]])
init_guess_by_mod_huckel(mol=None)[source]

Generate initial guess density matrix from a Huckel calculation based on occupancy averaged atomic RHF calculations, doi:10.1021/acs.jctc.8b01089

In contrast to init_guess_by_huckel, this routine employs the updated GWH rule from doi:10.1021/ja00480a005 to form the guess.

Returns:

Density matrix, 2D ndarray

make_rdm1(mo_coeff=None, mo_occ=None, **kwargs)[source]

One-particle density matrix. mo_occ is a 1D array, with occupancy 1 or 2.

property nelec
spin_square(mo_coeff=None, s=None)[source]

Spin square and multiplicity of RHF determinant

stability(internal=True, external=False, verbose=None, return_status=False)[source]

ROHF/ROKS stability analysis.

See also pyscf.scf.stability.rohf_stability function.

Kwargs:
internalbool

Internal stability, within the RHF optimization space.

externalbool

External stability. It is not available in current version.

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.

to_ks(xc='HF')[source]

Convert to RKS object.

pyscf.pbc.scf.rsjk module

Range separation JK builder

Ref:

Q. Sun, arXiv:2012.07929 Q. Sun, arXiv:2302.11307

class pyscf.pbc.scf.rsjk.RangeSeparatedJKBuilder(cell, kpts=array([[0., 0., 0.]]))[source]

Bases: StreamObject

build(omega=None, intor='int2e')[source]
dump_flags(verbose=None)[source]
get_jk(dm_kpts, hermi=1, kpts=None, kpts_band=None, with_j=True, with_k=True, omega=None, exxdiv=None)[source]
has_long_range()[source]

Whether to add the long-range part computed with AFT/FFT integrals

property qindex
reset(cell=None)[source]
weighted_coulG(kpt=array([0., 0., 0.]), exx=False, mesh=None, omega=None)

Weighted regular Coulomb kernel, applying cell.omega by default

weighted_coulG_LR(kpt=array([0., 0., 0.]), exx=False, mesh=None)[source]
weighted_coulG_SR(kpt=array([0., 0., 0.]), exx=False, mesh=None)[source]
pyscf.pbc.scf.rsjk.RangeSeparationJKBuilder

alias of RangeSeparatedJKBuilder

pyscf.pbc.scf.rsjk.estimate_ke_cutoff_for_omega(cell, omega, precision=None)[source]

Energy cutoff for FFTDF to converge attenuated Coulomb in moment space

pyscf.pbc.scf.rsjk.estimate_omega_for_ke_cutoff(cell, ke_cutoff, precision=None)[source]

The minimal omega in attenuated Coulombl given energy cutoff

pyscf.pbc.scf.rsjk.estimate_rcut(rs_cell, omega, precision=None, exclude_dd_block=True)[source]

Estimate rcut for 2e SR-integrals

pyscf.pbc.scf.scfint module

SCF (Hartree-Fock and DFT) tools for periodic systems at a single k-point,

using analytical GTO integrals instead of PWs.

See Also:

kscf.py : SCF tools for periodic systems with k-point sampling.

pyscf.pbc.scf.scfint.get_hcore(cell, kpt=None)[source]

Get the core Hamiltonian AO matrix, following dft.rks.get_veff_().

pyscf.pbc.scf.scfint.get_int1e(intor, cell, kpt=None)[source]

Get the one-electron integral defined by intor using lattice sums.

pyscf.pbc.scf.scfint.get_int1e_cross(intor, cell1, cell2, kpt=None, comp=1)[source]

1-electron integrals from two molecules like

\[\langle \mu | intor | \nu \rangle, \mu \in cell1, \nu \in cell2\]
pyscf.pbc.scf.scfint.get_ovlp(cell, kpt=None)[source]

Get the overlap AO matrix.

pyscf.pbc.scf.scfint.get_t(cell, kpt=None)[source]

Get the kinetic energy AO matrix.

pyscf.pbc.scf.stability module

Wave Function Stability Analysis

Ref. JCP 66, 3045 (1977); DOI:10.1063/1.434318 JCP 104, 9047 (1996); DOI:10.1063/1.471637

See also tddft/rhf.py and scf/newton_ah.py

pyscf.pbc.scf.stability.rhf_external(mf, verbose=None)[source]
pyscf.pbc.scf.stability.rhf_internal(mf, verbose=None)[source]
pyscf.pbc.scf.stability.rhf_stability(mf, internal=True, external=False, verbose=None)[source]
pyscf.pbc.scf.stability.uhf_external(mf, verbose=None)[source]
pyscf.pbc.scf.stability.uhf_internal(mf, verbose=None)[source]
pyscf.pbc.scf.stability.uhf_stability(mf, internal=True, external=False, verbose=None)[source]

pyscf.pbc.scf.uhf module

Unrestricted Hartree-Fock for periodic systems at a single k-point

See Also:

pyscf/pbc/scf/khf.py : Hartree-Fock for periodic systems with k-point sampling

class pyscf.pbc.scf.uhf.UHF(cell, kpt=array([0., 0., 0.]), exxdiv='ewald')[source]

Bases: SCF, UHF

UHF class for PBCs.

CCSD(*args, **kwargs)

restricted CCSD

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

conv_tolfloat

converge threshold. Default is 1e-7.

conv_tol_normtfloat

converge threshold for norm(t1,t2). Default is 1e-5.

max_cycleint

max number of iterations. Default is 50.

diis_spaceint

DIIS space size. Default is 6.

diis_start_cycleint

The step to start DIIS. Default is 0.

iterative_dampingfloat

The self consistent damping parameter.

directbool

AO-direct CCSD. Default is False.

async_iobool

Allow for asynchronous function execution. Default is True.

incore_completebool

Avoid all I/O (also for DIIS). Default is False.

level_shiftfloat

A shift on virtual orbital energies to stablize the CCSD iteration

frozenint or list

If integer is given, the inner-most orbitals are frozen from CC amplitudes. Given the orbital indices (0-based) in a list, both occupied and virtual orbitals can be frozen in CC calculation.

>>> mol = gto.M(atom = 'H 0 0 0; F 0 0 1.1', basis = 'ccpvdz')
>>> mf = scf.RHF(mol).run()
>>> # freeze 2 core orbitals
>>> mycc = cc.CCSD(mf).set(frozen = 2).run()
>>> # auto-generate the number of core orbitals to be frozen (1 in this case)
>>> mycc = cc.CCSD(mf).set_frozen().run()
>>> # freeze 2 core orbitals and 3 high lying unoccupied orbitals
>>> mycc.set(frozen = [0,1,16,17,18]).run()
callbackfunction(envs_dict) => None

callback function takes one dict as the argument which is generated by the builtin function locals(), so that the callback function can access all local variables in the current environment.

Saved results

convergedbool

CCSD converged or not

e_corrfloat

CCSD correlation correction

e_totfloat

Total CCSD energy (HF + correlation)

t1, t2 :

T amplitudes t1[i,a], t2[i,j,a,b] (i,j in occ, a,b in virt)

l1, l2 :

Lambda amplitudes l1[i,a], l2[i,j,a,b] (i,j in occ, a,b in virt)

CISD(*args, **kwargs)
MP2(*args, **kwargs)
TDA(*args, **kwargs)

Tamm-Dancoff approximation

Attributes:
conv_tolfloat

Diagonalization convergence tolerance. Default is 1e-9.

nstatesint

Number of TD states to be computed. Default is 3.

Saved results:

convergedbool

Diagonalization converged or not

e1D array

excitation energy for each excited state.

xyA list of two 2D arrays

The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.

TDHF(*args, **kwargs)
analyze(verbose=None, with_meta_lowdin=True, **kwargs)[source]

Analyze the given SCF object: print orbital energies, occupancies; print orbital coefficients; Mulliken population analysis; Diople moment.

canonicalize(mo_coeff, mo_occ, fock=None)

Canonicalization diagonalizes the UHF Fock matrix within occupied, virtual subspaces separatedly (without change occupancy).

convert_from_(mf)[source]

Convert given mean-field object to UHF

dip_moment(cell=None, dm=None, unit='Debye', verbose=3, **kwargs)[source]

Dipole moment in the cell.

Args:

cell : an instance of Cell

dm_kpts (a list of ndarrays) : density matrices of k-points

Return:

A list: the dipole moment on x, y and z components

dump_flags(verbose=None)[source]
eig(fock, s)[source]

Solver for generalized eigenvalue problem

\[HC = SCE\]
energy_elec(dm=None, h1e=None, vhf=None)

Electronic energy of Unrestricted Hartree-Fock

Note this function has side effects which cause mf.scf_summary updated.

Returns:

Hartree-Fock electronic energy and the 2-electron part contribution

gen_response(mo_coeff=None, mo_occ=None, with_j=True, hermi=0, max_memory=None)

Generate a function to compute the product of UHF response function and UHF density matrices.

get_bands(kpts_band, cell=None, dm=None, kpt=None)[source]

Get energy bands at the given (arbitrary) ‘band’ k-points.

Returns:
mo_energy(nmo,) ndarray or a list of (nmo,) ndarray

Bands energies E_n(k)

mo_coeff(nao, nmo) ndarray or a list of (nao,nmo) ndarray

Band orbitals psi_n(k)

get_fock(h1e=None, s1e=None, vhf=None, dm=None, cycle=-1, diis=None, diis_start_cycle=None, level_shift_factor=None, damp_factor=None)

F = h^{core} + V^{HF}

Special treatment (damping, DIIS, or level shift) will be applied to the Fock matrix if diis and cycle is specified (The two parameters are passed to get_fock function during the SCF iteration)

Kwargs:
h1e2D ndarray

Core hamiltonian

s1e2D ndarray

Overlap matrix, for DIIS

vhf2D ndarray

HF potential matrix

dm2D ndarray

Density matrix, for DIIS

cycleint

Then present SCF iteration step, for DIIS

diisan object of SCF.DIIS class

DIIS object to hold intermediate Fock and error vectors

diis_start_cycleint

The step to start DIIS. Default is 0.

level_shift_factorfloat or int

Level shift (in AU) for virtual space. Default is 0.

get_grad(mo_coeff, mo_occ, fock=None)[source]

RHF orbital gradients

Args:
mo_coeff2D ndarray

Obital coefficients

mo_occ1D ndarray

Orbital occupancy

fock_ao2D ndarray

Fock matrix in AO representation

Returns:

Gradients in MO representation. It’s a num_occ*num_vir vector.

get_init_guess(cell=None, key='minao')[source]
get_occ(mo_energy=None, mo_coeff=None)

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])
get_rho(dm=None, grids=None, kpt=None)

Compute density in real space

get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)[source]

Hartree-Fock potential matrix for the given density matrix. See scf.hf.get_veff() and scf.hf.RHF.get_veff()

init_guess_by_1e(cell=None)[source]

Generate initial guess density matrix from core hamiltonian

Returns:

Density matrix, 2D ndarray

init_guess_by_atom(mol=None, breaksym=True)[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(chk=None, project=True, kpt=None)[source]

Read the HF results from checkpoint file, then project it to the basis defined by mol

Returns:

Density matrix, 2D ndarray

init_guess_by_huckel(mol=None, breaksym=True)[source]

Generate initial guess density matrix from a Huckel calculation based on occupancy averaged atomic RHF calculations, doi:10.1021/acs.jctc.8b01089

Returns:

Density matrix, 2D ndarray

init_guess_by_minao(mol=None, breaksym=True)[source]

Initial guess in terms of the overlap to minimal basis.

init_guess_by_mod_huckel(mol=None, breaksym=True)[source]

Generate initial guess density matrix from a Huckel calculation based on occupancy averaged atomic RHF calculations, doi:10.1021/acs.jctc.8b01089

In contrast to init_guess_by_huckel, this routine employs the updated GWH rule from doi:10.1021/ja00480a005 to form the guess.

Returns:

Density matrix, 2D ndarray

make_rdm1(mo_coeff=None, mo_occ=None, **kwargs)[source]

One-particle density matrix in AO representation

Args:
mo_coefftuple of 2D ndarrays

Orbital coefficients for alpha and beta spins. Each column is one orbital.

mo_occtuple of 1D ndarrays

Occupancies for alpha and beta spins.

Returns:

A list of 2D ndarrays for alpha and beta spins

make_rdm2(mo_coeff=None, mo_occ=None, **kwargs)[source]

Two-particle density matrix in AO representation

Args:
mo_coefftuple of 2D ndarrays

Orbital coefficients for alpha and beta spins. Each column is one orbital.

mo_occtuple of 1D ndarrays

Occupancies for alpha and beta spins.

Returns:

A tuple of three 4D ndarrays for alpha,alpha and alpha,beta and beta,beta spins

mulliken_meta(mol=None, dm=None, verbose=5, pre_orth_method='ANO', s=None)[source]

Mulliken population analysis, based on meta-Lowdin AOs. In the meta-lowdin, the AOs are grouped in three sets: core, valence and Rydberg, the orthogonalization are carried out within each subsets.

Args:

mol : an instance of Mole

dmndarray or 2-item list of ndarray

Density matrix. ROHF dm is a 2-item list of 2D array

Kwargs:

verbose : int or instance of lib.logger.Logger

pre_orth_methodstr

Pre-orthogonalization, which localized GTOs for each atom. To obtain the occupied and unoccupied atomic shells, there are three methods

‘ano’ : Project GTOs to ANO basis
‘minao’ : Project GTOs to MINAO basis
‘scf’ : Symmetry-averaged fractional occupation atomic RHF
Returns:

A list : pop, charges

popnparray

Mulliken population on each atomic orbitals

chargesnparray

Mulliken charges

mulliken_pop(mol=None, dm=None, s=None, verbose=5)[source]

Mulliken population analysis

\[M_{ij} = D_{ij} S_{ji}\]

Mulliken charges

\[\delta_i = \sum_j M_{ij}\]
Returns:

A list : pop, charges

popnparray

Mulliken population on each atomic orbitals

chargesnparray

Mulliken charges

property nelec
spin_square(mo_coeff=None, s=None)[source]

Spin square and multiplicity of UHF determinant

\[S^2 = \frac{1}{2}(S_+ S_- + S_- S_+) + S_z^2\]

where \(S_+ = \sum_i S_{i+}\) is effective for all beta occupied orbitals; \(S_- = \sum_i S_{i-}\) is effective for all alpha occupied orbitals.

  1. There are two possibilities for \(S_+ S_-\)
    1. same electron \(S_+ S_- = \sum_i s_{i+} s_{i-}\),

    \[\sum_i \langle UHF|s_{i+} s_{i-}|UHF\rangle = \sum_{pq}\langle p|s_+s_-|q\rangle \gamma_{qp} = n_\alpha\]

    2) different electrons \(S_+ S_- = \sum s_{i+} s_{j-}, (i\neq j)\). There are in total \(n(n-1)\) terms. As a two-particle operator,

    \[\langle S_+ S_- \rangle = \langle ij|s_+ s_-|ij\rangle - \langle ij|s_+ s_-|ji\rangle = -\langle i^\alpha|j^\beta\rangle \langle j^\beta|i^\alpha\rangle\]
  2. Similarly, for \(S_- S_+\)
    1. same electron

    \[\sum_i \langle s_{i-} s_{i+}\rangle = n_\beta\]
    1. different electrons

    \[\langle S_- S_+ \rangle = -\langle i^\beta|j^\alpha\rangle \langle j^\alpha|i^\beta\rangle\]
  3. For \(S_z^2\)
    1. same electron

    \[\langle s_z^2\rangle = \frac{1}{4}(n_\alpha + n_\beta)\]
    1. different electrons

    \[\begin{split}&\frac{1}{2}\sum_{ij}(\langle ij|2s_{z1}s_{z2}|ij\rangle -\langle ij|2s_{z1}s_{z2}|ji\rangle) \\ &=\frac{1}{4}(\langle i^\alpha|i^\alpha\rangle \langle j^\alpha|j^\alpha\rangle - \langle i^\alpha|i^\alpha\rangle \langle j^\beta|j^\beta\rangle - \langle i^\beta|i^\beta\rangle \langle j^\alpha|j^\alpha\rangle + \langle i^\beta|i^\beta\rangle \langle j^\beta|j^\beta\rangle) \\ &-\frac{1}{4}(\langle i^\alpha|j^\alpha\rangle \langle j^\alpha|i^\alpha\rangle + \langle i^\beta|j^\beta\rangle\langle j^\beta|i^\beta\rangle) \\ &=\frac{1}{4}(n_\alpha^2 - n_\alpha n_\beta - n_\beta n_\alpha + n_\beta^2) -\frac{1}{4}(n_\alpha + n_\beta) \\ &=\frac{1}{4}((n_\alpha-n_\beta)^2 - (n_\alpha+n_\beta))\end{split}\]

In total

\[\begin{split}\langle S^2\rangle &= \frac{1}{2} (n_\alpha-\sum_{ij}\langle i^\alpha|j^\beta\rangle \langle j^\beta|i^\alpha\rangle +n_\beta -\sum_{ij}\langle i^\beta|j^\alpha\rangle\langle j^\alpha|i^\beta\rangle) + \frac{1}{4}(n_\alpha-n_\beta)^2 \\\end{split}\]
Args:
moa list of 2 ndarrays

Occupied alpha and occupied beta orbitals

Kwargs:
sndarray

AO overlap

Returns:

A list of two floats. The first is the expectation value of S^2. The second is the corresponding 2S+1

Examples:

>>> mol = gto.M(atom='O 0 0 0; H 0 0 1; H 0 1 0', basis='ccpvdz', charge=1, spin=1, verbose=0)
>>> mf = scf.UHF(mol)
>>> mf.kernel()
-75.623975516256706
>>> mo = (mf.mo_coeff[0][:,mf.mo_occ[0]>0], mf.mo_coeff[1][:,mf.mo_occ[1]>0])
>>> print('S^2 = %.7f, 2S+1 = %.7f' % spin_square(mo, mol.intor('int1e_ovlp_sph')))
S^2 = 0.7570150, 2S+1 = 2.0070027
stability(internal=True, external=False, verbose=None, return_status=False)[source]

Stability analysis for UHF/UKS method.

See also pyscf.scf.stability.uhf_stability function.

Args:

mf : UHF or UKS object

Kwargs:
internalbool

Internal stability, within the UHF space.

externalbool

External stability. Including the UHF -> GHF and real -> complex stability analysis.

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.

to_ks(xc='HF')[source]

Convert to RKS object.

pyscf.pbc.scf.uhf.dip_moment(cell, dm, unit='Debye', verbose=3, grids=None, rho=None, kpt=array([0., 0., 0.]))[source]

Dipole moment in the cell.

Args:

cell : an instance of Cell

dm_kpts (a list of ndarrays) : density matrices of k-points

Return:

A list: the dipole moment on x, y and z components

pyscf.pbc.scf.uhf.init_guess_by_chkfile(cell, chkfile_name, project=None, kpt=None)[source]

Read the HF results from checkpoint file and make the density matrix for UHF initial guess.

Returns:

Density matrix, (nao,nao) ndarray

Module contents

Hartree-Fock for periodic systems

pyscf.pbc.scf.HF(cell, *args, **kwargs)[source]
pyscf.pbc.scf.KGHF(cell, *args, **kwargs)[source]
pyscf.pbc.scf.KHF(cell, *args, **kwargs)[source]
pyscf.pbc.scf.KKS(cell, *args, **kwargs)[source]
pyscf.pbc.scf.KRHF(cell, *args, **kwargs)[source]
pyscf.pbc.scf.KRKS(cell, *args, **kwargs)[source]
pyscf.pbc.scf.KROKS(cell, *args, **kwargs)[source]
pyscf.pbc.scf.KS(cell, *args, **kwargs)[source]
pyscf.pbc.scf.KUHF(cell, *args, **kwargs)[source]
pyscf.pbc.scf.KUKS(cell, *args, **kwargs)[source]
pyscf.pbc.scf.RHF(cell, *args, **kwargs)[source]
pyscf.pbc.scf.RKS(cell, *args, **kwargs)[source]
pyscf.pbc.scf.ROKS(cell, *args, **kwargs)[source]
pyscf.pbc.scf.UKS(cell, *args, **kwargs)[source]