pyscf.df package

Subpackages

Submodules

pyscf.df.addons module

pyscf.df.addons.aug_etb(mol, beta=2.0)[source]

To generate the even-tempered auxiliary Gaussian basis

pyscf.df.addons.aug_etb_for_dfbasis(mol, dfbasis='weigend', beta=2.0, start_at=36)[source]

augment weigend basis with even-tempered gaussian basis exps = alpha*beta^i for i = 1..N

class pyscf.df.addons.load(eri, dataname='j3c')[source]

Bases: load

load 3c2e integrals from hdf5 file. It can be used in the context manager:

with load(cderifile) as eri:

print(eri.shape)

pyscf.df.addons.make_auxbasis(mol, mp2fit=False)[source]

Depending on the orbital basis, generating even-tempered Gaussians or the optimized auxiliary basis defined in DEFAULT_AUXBASIS

pyscf.df.addons.make_auxmol(mol, auxbasis=None)[source]

Generate a fake Mole object which uses the density fitting auxbasis as the basis sets. If auxbasis is not specified, the optimized auxiliary fitting basis set will be generated according to the rules recorded in pyscf.df.addons.DEFAULT_AUXBASIS. If the optimized auxiliary basis is not available (either not specified in DEFAULT_AUXBASIS or the basis set of the required elements not defined in the optimized auxiliary basis), even-tempered Gaussian basis set will be generated.

See also the paper JCTC, 13, 554 about generating auxiliary fitting basis.

pyscf.df.df module

J-metric density fitting

class pyscf.df.df.DF(mol, auxbasis=None)[source]

Bases: StreamObject

Object to hold 3-index tensor

Attributes:
auxbasisstr or dict

Same input format as Mole.basis

auxmolMole object

Read only Mole object to hold the auxiliary basis. auxmol is generated automatically in the initialization step based on the given auxbasis. It is used in the rest part of the code to determine the problem size, the integral batches etc. This object should NOT be modified.

_cderi_to_savestr

If _cderi_to_save is specified, the DF integral tensor will be saved in this file.

_cderistr or numpy array

If _cderi is specified, the DF integral tensor will be read from this HDF5 file (or numpy array). When the DF integral tensor is provided from the HDF5 file, its dataset name should be consistent with DF._dataname, which is ‘j3c’ by default. The DF integral tensor \(V_{x,ij}\) should be a 2D array in C (row-major) convention, where x corresponds to index of auxiliary basis, and the composite index ij is the orbital pair index. The hermitian symmetry is assumed for the composite ij index, ie the elements of \(V_{x,i,j}\) with \(i\geq j\) are existed in the DF integral tensor. Thus the shape of DF integral tensor is (M,N*(N+1)/2), where M is the number of auxbasis functions and N is the number of basis functions of the orbital basis.

blockdimint

When reading DF integrals from disk the chunk size to load. It is used to improve IO performance.

ao2mo(mo_coeffs, compact=True)[source]
property auxbasis
blockdim = 240
build()[source]
dump_flags(verbose=None)[source]
get_ao_eri()
get_eri()[source]
get_jk(dm, hermi=1, with_j=True, with_k=True, direct_scf_tol=1e-13, omega=None)[source]
get_mo_eri(mo_coeffs, compact=True)
get_naoaux()[source]
kernel(*args, **kwargs)[source]

Kernel function is the main driver of a method. Every method should define the kernel function as the entry of the calculation. Note the return value of kernel function is not strictly defined. It can be anything related to the method (such as the energy, the wave-function, the DFT mesh grids etc.).

loop(blksize=None)[source]
prange(start, end, step)[source]
range_coulomb(omega)[source]

Creates a temporary density fitting object for RSH-DF integrals. In this context, only LR or SR integrals for mol and auxmol are computed.

reset(mol=None)[source]

Reset mol and clean up relevant attributes for scanner mode

class pyscf.df.df.DF4C(mol, auxbasis=None)[source]

Bases: DF

Relativistic 4-component

ao2mo(mo_coeffs)[source]
build()[source]
get_jk(dm, hermi=1, with_j=True, with_k=True, direct_scf_tol=1e-13, omega=None)[source]
loop(blksize=None)[source]
pyscf.df.df.GDF

alias of DF

pyscf.df.df.GDF4C

alias of DF4C

pyscf.df.df_jk module

pyscf.df.df_jk.density_fit(mf, auxbasis=None, with_df=None, only_dfj=False)[source]

For the given SCF object, update the J, K matrix constructor with corresponding density fitting integrals.

Args:

mf : an SCF object

Kwargs:
auxbasisstr or basis dict

Same format to the input attribute mol.basis. If auxbasis is None, optimal auxiliary basis based on AO basis (if possible) or even-tempered Gaussian basis will be used.

only_dfjstr

Compute Coulomb integrals only and no approximation for HF exchange. Same to RIJONX in ORCA

Returns:

An SCF object with a modified J, K matrix constructor which uses density fitting integrals to compute J and K

Examples:

>>> mol = gto.M(atom='H 0 0 0; F 0 0 1', basis='ccpvdz', verbose=0)
>>> mf = scf.density_fit(scf.RHF(mol))
>>> mf.scf()
-100.005306000435510
>>> mol.symmetry = 1
>>> mol.build(0, 0)
>>> mf = scf.density_fit(scf.UHF(mol))
>>> mf.scf()
-100.005306000435510
pyscf.df.df_jk.get_j(dfobj, dm, hermi=1, direct_scf_tol=1e-13)[source]
pyscf.df.df_jk.get_jk(dfobj, dm, hermi=1, with_j=True, with_k=True, direct_scf_tol=1e-13)[source]
pyscf.df.df_jk.r_get_jk(dfobj, dms, hermi=1, with_j=True, with_k=True)[source]

Relativistic density fitting JK

pyscf.df.incore module

pyscf.df.incore.aux_e1(mol, auxmol_or_auxbasis, intor='int3c2e', aosym='s1', comp=None, out=None)[source]

3-center 2-electron AO integrals (L|ij), where L is the auxiliary basis.

Note aux_e1 is basically analogous to aux_e2 function. It can be viewed as the version of transposed aux_e2 tensor: if comp == 1:

aux_e1 = aux_e2().T

else:

aux_e1 = aux_e2().transpose(0,2,1)

The same arguments as function aux_e2 can be input to aux_e1.

pyscf.df.incore.aux_e2(mol, auxmol_or_auxbasis, intor='int3c2e', aosym='s1', comp=None, out=None, cintopt=None)[source]

3-center AO integrals (ij|L), where L is the auxiliary basis.

Kwargs:
cintoptLibcint-3.14 and newer version support to compute int3c2e

without the opt for the 3rd index. It can be precomputed to reduce the overhead of cintopt initialization repeatedly.

cintopt = gto.moleintor.make_cintopt(mol._atm, mol._bas, mol._env, ‘int3c2e’)

pyscf.df.incore.cholesky_eri(mol, auxbasis='weigend+etb', auxmol=None, int3c='int3c2e', aosym='s2ij', int2c='int2c2e', comp=1, max_memory=2000, decompose_j2c='cd', lindep=1e-12, verbose=0, fauxe2=<function aux_e2>)[source]
Returns:

2D array of (naux,nao*(nao+1)/2) in C-contiguous

pyscf.df.incore.cholesky_eri_debug(mol, auxbasis='weigend+etb', auxmol=None, int3c='int3c2e', aosym='s2ij', int2c='int2c2e', comp=1, verbose=0, fauxe2=<function aux_e2>)[source]
Returns:

2D array of (naux,nao*(nao+1)/2) in C-contiguous

pyscf.df.incore.fill_2c2e(mol, auxmol_or_auxbasis, intor='int2c2e', comp=None, hermi=1, out=None)[source]

2-center 2-electron AO integrals for auxiliary basis (auxmol)

pyscf.df.outcore module

pyscf.df.outcore.cholesky_eri(mol, erifile, auxbasis='weigend+etb', dataname='j3c', tmpdir=None, int3c='int3c2e', aosym='s2ij', int2c='int2c2e', comp=1, max_memory=2000, auxmol=None, verbose=3)[source]

3-index density-fitting tensor.

pyscf.df.outcore.cholesky_eri_b(mol, erifile, auxbasis='weigend+etb', dataname='j3c', int3c='int3c2e', aosym='s2ij', int2c='int2c2e', comp=1, max_memory=2000, auxmol=None, decompose_j2c='CD', lindep=1e-12, verbose=3)[source]

3-center 2-electron DF tensor. Similar to cholesky_eri while this function stores DF tensor in blocks.

Args:
dataname: string

Dataset label of the DF tensor in HDF5 file.

decompose_j2c: string

The method to decompose the metric defined by int2c. It can be set to CD (cholesky decomposition) or ED (eigenvalue decomposition).

lindepfloat

The threshold to discard linearly dependent basis when decompose_j2c is set to ED.

pyscf.df.outcore.general(mol, mo_coeffs, erifile, auxbasis='weigend+etb', dataname='eri_mo', tmpdir=None, int3c='int3c2e', aosym='s2ij', int2c='int2c2e', comp=1, max_memory=2000, verbose=0, compact=True)[source]

Transform ij of (ij|L) to MOs.

pyscf.df.r_incore module

pyscf.df.r_incore.aux_e1(mol, auxmol, intor='int3c2e_spinor', aosym='s1', comp=1, hermi=0)[source]
pyscf.df.r_incore.aux_e2(mol, auxmol, intor='int3c2e_spinor', aosym='s1', comp=None, hermi=0)[source]
pyscf.df.r_incore.cholesky_eri(mol, auxbasis='weigend+etb', auxmol=None, int3c='int3c2e_spinor', aosym='s1', int2c='int2c2e_sph', comp=1, verbose=0)[source]

Module contents

Density fitting

This module provides the fundamental functions to handle the 3-index tensors (including the 3-center 2-electron AO and MO integrals, the Cholesky decomposed integrals) required by the density fitting method or the RI (resolution of identity) approximation.

Simple usage:

>>> from pyscf import gto, dft
>>> mol = gto.M(atom='N 0 0 0; N 0 0 1', basis='ccpvdz')
>>> mf = dft.RKS(mol).density_fit().run()
pyscf.df.density_fit(obj, *args, **kwargs)[source]

Given SCF/MCSCF or post-HF object, use density fitting technique to approximate the 2e integrals.