pyscf.pbc.df package#
Submodules#
pyscf.pbc.df.aft module#
Density expansion on plane waves
- class pyscf.pbc.df.aft.AFTDF(cell, kpts=array([[0., 0., 0.]]))[source]#
Bases:
StreamObject
,AFTDFMixin
Density expansion on plane waves
- ao2mo(mo_coeffs, kpts=None, compact=True)#
- ao2mo_7d(mo_coeff_kpts, kpts=None, factor=1, out=None)#
- 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.
- get_ao_eri(kpts=None, compact=True)#
- get_ao_pairs(kpts=array([[0., 0., 0.], [0., 0., 0.]]), q=None, shls_slice=None, compact=False)#
Calculate forward Fourier transform (G|ij) of all AO pairs.
- Returns:
- ao_pairs_G2D complex array
For gamma point, the shape is (ngrids, nao*(nao+1)/2); otherwise the shape is (ngrids, nao*nao)
- get_ao_pairs_G(kpts=array([[0., 0., 0.], [0., 0., 0.]]), q=None, shls_slice=None, compact=False)#
Calculate forward Fourier transform (G|ij) of all AO pairs.
- Returns:
- ao_pairs_G2D complex array
For gamma point, the shape is (ngrids, nao*(nao+1)/2); otherwise the shape is (ngrids, nao*nao)
- get_eri(kpts=None, compact=True)#
- get_jk(dm, hermi=1, kpts=None, kpts_band=None, with_j=True, with_k=True, omega=None, exxdiv=None)[source]#
- get_mo_eri(mo_coeffs, kpts=None, compact=True)#
- get_mo_pairs(mo_coeffs, kpts=array([[0., 0., 0.], [0., 0., 0.]]), q=None, compact=False)#
Calculate forward fourier transform (G|ij) of all MO pairs.
- Args:
- mo_coeff: length-2 list of (nao,nmo) ndarrays
The two sets of MO coefficients to use in calculating the product |ij).
- Returns:
- mo_pairs_G(ngrids, nmoi*nmoj) ndarray
The FFT of the real-space MO pairs.
- get_mo_pairs_G(mo_coeffs, kpts=array([[0., 0., 0.], [0., 0., 0.]]), q=None, compact=False)#
Calculate forward fourier transform (G|ij) of all MO pairs.
- Args:
- mo_coeff: length-2 list of (nao,nmo) ndarrays
The two sets of MO coefficients to use in calculating the product |ij).
- Returns:
- mo_pairs_G(ngrids, nmoi*nmoj) ndarray
The FFT of the real-space MO pairs.
- get_nuc(kpts=None)#
Get the periodic nuc-el AO matrix, with G=0 removed.
- Kwargs:
function _guess_eta from module pbc.df.gdf_builder.
- get_pp(kpts=None)#
Get the periodic pseudopotential nuc-el AO matrix, with G=0 removed.
- Kwargs:
mesh: custom mesh grids. By default mesh is determined by the function _guess_eta from module pbc.df.gdf_builder.
- class pyscf.pbc.df.aft.AFTDFMixin[source]#
Bases:
object
- ft_loop(mesh=None, q=array([0., 0., 0.]), kpts=None, shls_slice=None, max_memory=4000, aosym='s1', intor='GTO_ft_ovlp', comp=1, bvk_kmesh=None, return_complex=True)[source]#
- Fourier transform iterator for all kpti which satisfy
2pi*N = (kpts - kpti - q)*a, N = -1, 0, 1
- pw_loop(mesh=None, kpti_kptj=None, q=None, shls_slice=None, max_memory=2000, aosym='s1', blksize=None, intor='GTO_ft_ovlp', comp=1, bvk_kmesh=None)[source]#
Fourier transform iterator for AO pair
- 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.
- weighted_coulG(kpt=array([0., 0., 0.]), exx=False, mesh=None, omega=None)#
Weighted regular Coulomb kernel, applying cell.omega by default
- pyscf.pbc.df.aft.estimate_eta(cell, cutoff=None)#
Given rcut the boundary of repeated images of the cell, estimates the minimal exponent of the smooth compensated gaussian model charge, requiring that at boundary, density ~ 4pi rmax^2 exp(-eta/2*rmax^2) < cutoff
- pyscf.pbc.df.aft.estimate_eta_for_ke_cutoff(cell, ke_cutoff, precision=None)[source]#
Given ke_cutoff, the upper bound of eta to produce the required precision in AFTDF Coulomb integrals.
- pyscf.pbc.df.aft.estimate_eta_min(cell, cutoff=None)[source]#
Given rcut the boundary of repeated images of the cell, estimates the minimal exponent of the smooth compensated gaussian model charge, requiring that at boundary, density ~ 4pi rmax^2 exp(-eta/2*rmax^2) < cutoff
- pyscf.pbc.df.aft.estimate_ke_cutoff(cell, precision=None)[source]#
Energy cutoff estimation for 4-center Coulomb repulsion integrals
- pyscf.pbc.df.aft.estimate_ke_cutoff_for_eta(cell, eta, precision=None)[source]#
Given eta, the lower bound of ke_cutoff to produce the required precision in AFTDF Coulomb integrals.
- pyscf.pbc.df.aft.estimate_ke_cutoff_for_omega(cell, omega, precision=None)[source]#
Energy cutoff for AFTDF to converge attenuated Coulomb in moment space
- pyscf.pbc.df.aft.estimate_omega(cell, cutoff=None)#
Given cell.rcut the boundary of repeated images of the cell, estimates the minimal omega for the attenuated Coulomb interactions, requiring that at boundary the Coulomb potential of a point charge < cutoff
- pyscf.pbc.df.aft.estimate_omega_for_ke_cutoff(cell, ke_cutoff, precision=None)[source]#
The minimal omega in attenuated Coulombl given energy cutoff
- pyscf.pbc.df.aft.estimate_omega_min(cell, cutoff=None)[source]#
Given cell.rcut the boundary of repeated images of the cell, estimates the minimal omega for the attenuated Coulomb interactions, requiring that at boundary the Coulomb potential of a point charge < cutoff
- pyscf.pbc.df.aft.get_nuc(mydf, kpts=None)[source]#
Get the periodic nuc-el AO matrix, with G=0 removed.
- Kwargs:
function _guess_eta from module pbc.df.gdf_builder.
pyscf.pbc.df.aft_ao2mo module#
Integral transformation with analytic Fourier transformation
- pyscf.pbc.df.aft_ao2mo.get_ao_pairs_G(mydf, kpts=array([[0., 0., 0.], [0., 0., 0.]]), q=None, shls_slice=None, compact=False)[source]#
Calculate forward Fourier transform (G|ij) of all AO pairs.
- Returns:
- ao_pairs_G2D complex array
For gamma point, the shape is (ngrids, nao*(nao+1)/2); otherwise the shape is (ngrids, nao*nao)
- pyscf.pbc.df.aft_ao2mo.get_mo_pairs_G(mydf, mo_coeffs, kpts=array([[0., 0., 0.], [0., 0., 0.]]), q=None, compact=False)[source]#
Calculate forward fourier transform (G|ij) of all MO pairs.
- Args:
- mo_coeff: length-2 list of (nao,nmo) ndarrays
The two sets of MO coefficients to use in calculating the product |ij).
- Returns:
- mo_pairs_G(ngrids, nmoi*nmoj) ndarray
The FFT of the real-space MO pairs.
pyscf.pbc.df.aft_jk module#
JK with analytic Fourier transformation
- pyscf.pbc.df.aft_jk.get_j_for_bands(mydf, dm_kpts, hermi=1, kpts=array([[0., 0., 0.]]), kpts_band=None)[source]#
- pyscf.pbc.df.aft_jk.get_j_kpts(mydf, dm_kpts, hermi=1, kpts=array([[0., 0., 0.]]), kpts_band=None)[source]#
- pyscf.pbc.df.aft_jk.get_jk(mydf, dm, hermi=1, kpt=array([0., 0., 0.]), kpts_band=None, with_j=True, with_k=True, exxdiv=None)[source]#
JK for given k-point
pyscf.pbc.df.df module#
Density fitting
Divide the 3-center Coulomb integrals to two parts. Compute the local part in real space, long range part in reciprocal space.
Note when diffuse functions are used in fitting basis, it is easy to cause linear dependence (non-positive definite) issue under PBC.
Ref: J. Chem. Phys. 147, 164119 (2017)
- class pyscf.pbc.df.df.CDERIArray(data_group, label='j3c')[source]#
Bases:
object
Provide numpy APIs to access cderi tensor. This object can be viewed as an 5-dimension array [kpt-i, kpt-j, aux-index, ao-i, ao-j]
- class pyscf.pbc.df.df.GDF(cell, kpts=array([[0., 0., 0.]]))[source]#
Bases:
StreamObject
,AFTDFMixin
Gaussian density fitting
- ao2mo(mo_coeffs, kpts=None, compact=True)#
- ao2mo_7d(mo_coeff_kpts, kpts=None, factor=1, out=None)#
- property auxbasis#
- blockdim = 240#
- cderi_array(label=None)[source]#
Returns CDERIArray object which provides numpy APIs to access cderi tensor.
- force_dm_kbuild = False#
- get_ao_eri(kpts=None, compact=True)#
- get_eri(kpts=None, compact=True)#
- get_jk(dm, hermi=1, kpts=None, kpts_band=None, with_j=True, with_k=True, omega=None, exxdiv=None)[source]#
- get_mo_eri(mo_coeffs, kpts=None, compact=True)#
- property gs#
- prange(start, stop, step)[source]#
This is a hook for MPI parallelization. DO NOT use it out of the scope of AFTDF/GDF/MDF.
- sr_loop(kpti_kptj=array([[0., 0., 0.], [0., 0., 0.]]), max_memory=2000, compact=True, blksize=None, aux_slice=None)[source]#
Short range part
- to_gpu(out=None)#
Convert a method to its corresponding GPU variant, and recursively converts all attributes of a method to cupy objects or gpu4pyscf objects.
- pyscf.pbc.df.df.make_auxcell(cell, auxbasis=None, drop_eta=None)#
Generate a cell object using the density fitting auxbasis as the basis set. The normalization coefficients of the auxiliary cell are different to the regular (square-norm) convention. To simplify the compensated charge algorithm, they are normalized against int (r^l e^{-ar^2} r^2 dr
- pyscf.pbc.df.df.make_modrho_basis(cell, auxbasis=None, drop_eta=None)[source]#
Generate a cell object using the density fitting auxbasis as the basis set. The normalization coefficients of the auxiliary cell are different to the regular (square-norm) convention. To simplify the compensated charge algorithm, they are normalized against int (r^l e^{-ar^2} r^2 dr
pyscf.pbc.df.df_ao2mo module#
pyscf.pbc.df.df_jk module#
Density fitting with Gaussian basis Ref: J. Chem. Phys. 147, 164119 (2017)
- pyscf.pbc.df.df_jk.density_fit(mf, auxbasis=None, mesh=None, with_df=None)[source]#
Generte density-fitting SCF object
- Args:
- auxbasisstr or basis dict
Same format to the input attribute mol.basis. If auxbasis is None, auxiliary basis based on AO basis (if possible) or even-tempered Gaussian basis will be used.
- meshtuple
number of grids in each direction
with_df : DF object
- pyscf.pbc.df.df_jk.get_j_kpts(mydf, dm_kpts, hermi=1, kpts=array([[0., 0., 0.]]), kpts_band=None)[source]#
- pyscf.pbc.df.df_jk.get_j_kpts_kshift(mydf, dm_kpts, kshift, hermi=0, kpts=array([[0., 0., 0.]]), kpts_band=None)[source]#
- Math:
- J^{k1 k1’}_{pq}
= (1/Nk) sum_{k2} sum_{rs} (p k1 q k1’ |r k2’ s k2) D_{sr}^{k2 k2’}
- where k1’ and k2’ satisfies
(k1 - k1’ - kpts[kshift]) dot a = 2n pi (k2 - k2’ - kpts[kshift]) dot a = 2n pi
For kshift = 0,
get_j_kpts()
is called.
- pyscf.pbc.df.df_jk.get_jk(mydf, dm, hermi=1, kpt=array([0., 0., 0.]), kpts_band=None, with_j=True, with_k=True, exxdiv=None)[source]#
JK for given k-point
- pyscf.pbc.df.df_jk.get_k_kpts(mydf, dm_kpts, hermi=1, kpts=array([[0., 0., 0.]]), kpts_band=None, exxdiv=None)[source]#
- pyscf.pbc.df.df_jk.get_k_kpts_kshift(mydf, dm_kpts, kshift, hermi=0, kpts=array([[0., 0., 0.]]), kpts_band=None, exxdiv=None)[source]#
- Math:
- K^{k1 k1’}_{pq}
= (1/Nk) sum_{k2} sum_{rs} (p k1 s k2 | r k2’ q k1’) D_{sr}^{k2 k2’}
- where k1’ and k2’ satisfies
(k1 - k1’ - kpts[kshift]) dot a = 2n pi (k2 - k2’ - kpts[kshift]) dot a = 2n pi
For kshift = 0,
get_k_kpts()
is called.
pyscf.pbc.df.fft module#
Density expansion on plane waves
- class pyscf.pbc.df.fft.FFTDF(cell, kpts=array([[0., 0., 0.]]))[source]#
Bases:
StreamObject
Density expansion on plane waves
- ao2mo(mo_coeffs, kpts=None, compact=True)#
General MO integral transformation
- ao2mo_7d(mo_coeff_kpts, kpts=None, factor=1, out=None)#
- 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.
- get_ao_eri(kpts=None, compact=True)#
- get_ao_pairs(kpts=array([[0., 0., 0.], [0., 0., 0.]]), q=None, shls_slice=None, compact=False)#
Calculate forward (G|ij) FFT of all AO pairs.
- Returns:
- ao_pairs_G2D complex array
For gamma point, the shape is (ngrids, nao*(nao+1)/2); otherwise the shape is (ngrids, nao*nao)
- get_ao_pairs_G(kpts=array([[0., 0., 0.], [0., 0., 0.]]), q=None, shls_slice=None, compact=False)#
Calculate forward (G|ij) FFT of all AO pairs.
- Returns:
- ao_pairs_G2D complex array
For gamma point, the shape is (ngrids, nao*(nao+1)/2); otherwise the shape is (ngrids, nao*nao)
- get_eri(kpts=None, compact=True)#
- get_jk(dm, hermi=1, kpts=None, kpts_band=None, with_j=True, with_k=True, omega=None, exxdiv=None)[source]#
- get_mo_eri(mo_coeffs, kpts=None, compact=True)#
General MO integral transformation
- get_mo_pairs(mo_coeffs, kpts=array([[0., 0., 0.], [0., 0., 0.]]), q=None, compact=False)#
Calculate forward (G|ij) FFT of all MO pairs.
- Args:
- mo_coeff: length-2 list of (nao,nmo) ndarrays
The two sets of MO coefficients to use in calculating the product |ij).
- Returns:
- mo_pairs_G(ngrids, nmoi*nmoj) ndarray
The FFT of the real-space MO pairs.
- get_mo_pairs_G(mo_coeffs, kpts=array([[0., 0., 0.], [0., 0., 0.]]), q=None, compact=False)#
Calculate forward (G|ij) FFT of all MO pairs.
- Args:
- mo_coeff: length-2 list of (nao,nmo) ndarrays
The two sets of MO coefficients to use in calculating the product |ij).
- Returns:
- mo_pairs_G(ngrids, nmoi*nmoj) ndarray
The FFT of the real-space MO pairs.
- get_nuc(kpts=None)#
- get_pp(kpts=None)#
Get the periodic pseudopotential nuc-el AO matrix, with G=0 removed.
- property mesh#
- range_coulomb(omega)#
Creates a temporary density fitting object for RSH-DF integrals. In this context, only LR or SR integrals for mol and auxmol are computed.
- to_gpu(out=None)#
Convert a method to its corresponding GPU variant, and recursively converts all attributes of a method to cupy objects or gpu4pyscf objects.
pyscf.pbc.df.fft_ao2mo module#
Integral transformation with FFT
- (ij|kl) = int dr1 dr2 i*(r1) j(r1) v(r12) k*(r2) l(r2)
= (ij|G) v(G) (G|kl)
- i*(r) j(r) = 1/N sum_G e^{iGr} (G|ij)
= 1/N sum_G e^{-iGr} (ij|G)
- “forward” FFT:
(G|ij) = sum_r e^{-iGr} i*(r) j(r) = fft[ i*(r) j(r) ]
- “inverse” FFT:
- (ij|G) = sum_r e^{iGr} i*(r) j(r) = N * ifft[ i*(r) j(r) ]
= conj[ sum_r e^{-iGr} j*(r) i(r) ]
- pyscf.pbc.df.fft_ao2mo.general(mydf, mo_coeffs, kpts=None, compact=True)[source]#
General MO integral transformation
- pyscf.pbc.df.fft_ao2mo.get_ao_pairs_G(mydf, kpts=array([[0., 0., 0.], [0., 0., 0.]]), q=None, shls_slice=None, compact=False)[source]#
Calculate forward (G|ij) FFT of all AO pairs.
- Returns:
- ao_pairs_G2D complex array
For gamma point, the shape is (ngrids, nao*(nao+1)/2); otherwise the shape is (ngrids, nao*nao)
- pyscf.pbc.df.fft_ao2mo.get_mo_pairs_G(mydf, mo_coeffs, kpts=array([[0., 0., 0.], [0., 0., 0.]]), q=None, compact=False)[source]#
Calculate forward (G|ij) FFT of all MO pairs.
- Args:
- mo_coeff: length-2 list of (nao,nmo) ndarrays
The two sets of MO coefficients to use in calculating the product |ij).
- Returns:
- mo_pairs_G(ngrids, nmoi*nmoj) ndarray
The FFT of the real-space MO pairs.
pyscf.pbc.df.fft_jk module#
JK with discrete Fourier transformation
- pyscf.pbc.df.fft_jk.get_j(mydf, dm, hermi=1, 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
- kpt(3,) ndarray
The “inner” dummy k-point at which the DM was evaluated (or sampled).
- kpts_band(3,) ndarray or (*,3) ndarray
The “outer” primary k-point at which J and K are evaluated.
- Returns:
The function returns one J matrix, corresponding to the input density matrix (both order and shape).
- pyscf.pbc.df.fft_jk.get_j_e1_kpts(mydf, dm_kpts, kpts=array([[0., 0., 0.]]), kpts_band=None)[source]#
Derivatives of Coulomb (J) AO matrix at sampled k-points.
- pyscf.pbc.df.fft_jk.get_j_kpts(mydf, dm_kpts, hermi=1, kpts=array([[0., 0., 0.]]), 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.
kpts : (nkpts, 3) ndarray
- Kwargs:
- kpts_band(3,) ndarray or (*,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.df.fft_jk.get_jk(mydf, dm, hermi=1, kpt=array([0., 0., 0.]), kpts_band=None, with_j=True, with_k=True, exxdiv=None)[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
- kpt(3,) ndarray
The “inner” dummy k-point at which the DM was evaluated (or sampled).
- kpts_band(3,) ndarray or (*,3) ndarray
The “outer” primary 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.df.fft_jk.get_k(mydf, dm, hermi=1, kpt=array([0., 0., 0.]), kpts_band=None, exxdiv=None)[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
- kpt(3,) ndarray
The “inner” dummy k-point at which the DM was evaluated (or sampled).
- kpts_band(3,) ndarray or (*,3) ndarray
The “outer” primary 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.df.fft_jk.get_k_e1_kpts(mydf, dm_kpts, kpts=array([[0., 0., 0.]]), kpts_band=None, exxdiv=None)[source]#
Derivatives of exchange (K) AO matrix at sampled k-points.
- pyscf.pbc.df.fft_jk.get_k_kpts(mydf, dm_kpts, hermi=1, kpts=array([[0., 0., 0.]]), kpts_band=None, exxdiv=None)[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
kpts : (nkpts, 3) ndarray
- Kwargs:
- hermiint
Whether K matrix is hermitian
0 : not hermitian and not symmetric1 : hermitian- kpts_band(3,) ndarray or (*,3) ndarray
A list of arbitrary “band” k-points 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.df.ft_ao module#
Analytical Fourier transformation AO-pair product for PBC
- class pyscf.pbc.df.ft_ao.ExtendedMole[source]#
Bases:
Mole
An extended Mole object to mimic periodicity
- property bas_map#
- static bas_mask_to_segment(rs_cell, bas_mask, verbose=None)[source]#
bas_mask shape [bvk_ncells, nbas, nimgs]
- gen_ft_kernel(aosym='s1', intor='GTO_ft_ovlp', comp=1, return_complex=False, kpts=None, verbose=None)#
Generate the analytical fourier transform kernel for AO products
sum_T exp(-i k_j * T) int exp(-i(G+q)r) i(r) j(r-T) dr^3
- get_ovlp_mask(cutoff=None)[source]#
integral screening mask for basis product between cell and supmol
- property sh_loc#
- pyscf.pbc.df.ft_ao.estimate_rcut(cell, precision=None)[source]#
Estimate rcut for each basis based on Schwarz inequality Q_ij ~ S_ij * (sqrt(2aij/pi) * aij**(lij*2) * (4*lij-1)!!)**.5
- pyscf.pbc.df.ft_ao.ft_ao(mol, Gv, shls_slice=None, b=None, gxyz=None, Gvbase=None, kpt=array([0., 0., 0.]), verbose=None)[source]#
Analytical FT transform AO int mu(r) exp(-ikr) dr^3
The output tensor has the shape [nGv, nao]
- pyscf.pbc.df.ft_ao.ft_aopair(cell, Gv, shls_slice=None, aosym='s1', b=None, gxyz=None, Gvbase=None, kpti_kptj=array([[0., 0., 0.], [0., 0., 0.]]), q=None, intor='GTO_ft_ovlp', comp=1, verbose=None)[source]#
Fourier transform AO pair for a pair of k-points sum_T exp(-i k_j * T) int exp(-i(G+q)r) i(r) j(r-T) dr^3
- pyscf.pbc.df.ft_ao.ft_aopair_kpts(cell, Gv, shls_slice=None, aosym='s1', b=None, gxyz=None, Gvbase=None, q=array([0., 0., 0.]), kptjs=array([[0., 0., 0.]]), intor='GTO_ft_ovlp', comp=1, bvk_kmesh=None, out=None)[source]#
Fourier transform AO pair for a group of k-points sum_T exp(-i k_j * T) int exp(-i(G+q)r) i(r) j(r-T) dr^3
The return array holds the AO pair corresponding to the kpoints given by kptjs
pyscf.pbc.df.gdf_builder module#
Build GDF tensor with compensated charges
This algorithm can handle the LR-, SR- and regular density fitting integrals with the same framework. The RSGDF algorithms (rsdf.py rsdf_builder.py) are good for regular density fitting and SR-integral density fitting only.
- pyscf.pbc.df.gdf_builder.auxbar(fused_cell)[source]#
Potential average = sum_L V_L*Lpq
The coulomb energy is computed with chargeless density int (rho-C) V, C = (int rho) / vol = Tr(gamma,S)/vol It is equivalent to removing the averaged potential from the short range V vs = vs - (int V)/vol * S
- pyscf.pbc.df.gdf_builder.estimate_eta_for_ke_cutoff(cell, ke_cutoff, precision=None)[source]#
Given ke_cutoff, the upper bound of eta to produce the required precision in AFTDF Coulomb integrals.
- pyscf.pbc.df.gdf_builder.estimate_eta_min(cell, precision=None)[source]#
Given rcut the boundary of repeated images of the cell, estimates the minimal exponent of the smooth compensated gaussian model charge, requiring that at boundary, density ~ 4pi rmax^2 exp(-eta/2*rmax^2) < precision
- pyscf.pbc.df.gdf_builder.estimate_ke_cutoff_for_eta(cell, eta, precision=None)[source]#
Given eta, the lower bound of ke_cutoff to produce the required precision in AFTDF Coulomb integrals.
pyscf.pbc.df.incore module#
- class pyscf.pbc.df.incore.Int3cBuilder(cell, auxcell, kpts=None)[source]#
Bases:
StreamObject
helper functions to compute 3-center integral tensor with double-lattice sum
- gen_int3c_kernel(intor='int3c2e', aosym='s2', comp=None, j_only=False, reindex_k=None, rs_auxcell=None, auxcell=None, supmol=None, return_complex=False)[source]#
Generate function to compute int3c2e with double lattice-sum
rs_auxcell: range-separated auxcell for gdf/rsdf module
reindex_k: an index array to sort the order of k-points in output
- pyscf.pbc.df.incore.aux_e2(cell, auxcell_or_auxbasis, intor='int3c2e', aosym='s1', comp=None, kptij_lst=array([[[0., 0., 0.], [0., 0., 0.]]]), shls_slice=None, **kwargs)[source]#
3-center AO integrals (ij|L) with double lattice sum: sum_{lm} (i[l]j[m]|L[0]), where L is the auxiliary basis.
- Returns:
(nao_pair, naux) array
- pyscf.pbc.df.incore.aux_e2_sum_auxbas(cell, auxcell_or_auxbasis, intor='int3c2e', aosym='s1', comp=None, kptij_lst=array([[[0., 0., 0.], [0., 0., 0.]]]), shls_slice=None, **kwargs)[source]#
Compute \(\sum_{L} (ij|L)\) on the fly.
- Returns:
out : (nao_pair,) array
- pyscf.pbc.df.incore.estimate_rcut(cell, auxcell, precision=None)[source]#
Estimate rcut for 3c2e integrals
- pyscf.pbc.df.incore.fill_2c2e(cell, auxcell_or_auxbasis, intor='int2c2e', hermi=0, kpt=array([0., 0., 0.]))[source]#
2-center 2-electron AO integrals (L|ij), where L is the auxiliary basis.
- pyscf.pbc.df.incore.format_aux_basis(cell, auxbasis='weigend+etb')[source]#
For backward compatibility
- pyscf.pbc.df.incore.int3c1e_nuc_grad(cell, auxcell, dm, intor='int3c1e', aosym='s1', comp=3, kptij_lst=array([[[0., 0., 0.], [0., 0., 0.]]]), shls_slice=None, **kwargs)[source]#
Compute the nuclear gradient contribution to the 2nd local part of PP on the fly. See pbc.gto.pseudo.pp_int.vpploc_part2_nuc_grad.
- Returns:
out : (natm,comp) array
- pyscf.pbc.df.incore.make_auxmol(cell, auxbasis=None)#
See pyscf.df.addons.make_auxmol
- pyscf.pbc.df.incore.wrap_int3c(cell, auxcell, intor='int3c2e', aosym='s1', comp=1, kptij_lst=array([[[0., 0., 0.], [0., 0., 0.]]]), cintopt=None, pbcopt=None)[source]#
Generate a 3-center integral kernel which can be called repeatedly in the incore or outcore driver. The kernel function has a simple function signature f(shls_slice)
pyscf.pbc.df.mdf module#
Gaussian and planewaves mixed density fitting Ref: J. Chem. Phys. 147, 164119 (2017)
- class pyscf.pbc.df.mdf.MDF(cell, kpts=array([[0., 0., 0.]]))[source]#
Bases:
GDF
Gaussian and planewaves mixed density fitting
- ao2mo(mo_coeffs, kpts=None, compact=True)#
- ao2mo_7d(mo_coeff_kpts, kpts=None, factor=1, out=None)#
- get_ao_eri(kpts=None, compact=True)#
- get_eri(kpts=None, compact=True)#
- get_jk(dm, hermi=1, kpts=None, kpts_band=None, with_j=True, with_k=True, omega=None, exxdiv=None)[source]#
- get_mo_eri(mo_coeffs, kpts=None, compact=True)#
- get_nuc(kpts=None)#
Get the periodic nuc-el AO matrix, with G=0 removed.
- get_pp(kpts=None)#
Get the periodic pseudopotential nuc-el AO matrix, with G=0 removed.
pyscf.pbc.df.mdf_ao2mo module#
pyscf.pbc.df.mdf_jk module#
Exact density fitting with Gaussian and planewaves Ref: J. Chem. Phys. 147, 164119 (2017)
- pyscf.pbc.df.mdf_jk.density_fit(mf, auxbasis=None, mesh=None, with_df=None)[source]#
Generte density-fitting SCF object
- Args:
- auxbasisstr or basis dict
Same format to the input attribute mol.basis. If auxbasis is None, auxiliary basis based on AO basis (if possible) or even-tempered Gaussian basis will be used.
- meshtuple
number of grids in each direction
with_df : MDF object
- pyscf.pbc.df.mdf_jk.get_j_kpts(mydf, dm_kpts, hermi=1, kpts=array([[0., 0., 0.]]), kpts_band=None)[source]#
pyscf.pbc.df.outcore module#
- pyscf.pbc.df.outcore.aux_e1(cell, auxcell_or_auxbasis, erifile, intor='int3c2e', aosym='s2ij', comp=None, kptij_lst=None, dataname='eri_mo', shls_slice=None, max_memory=2000, verbose=0)[source]#
3-center AO integrals (L|ij) with double lattice sum: sum_{lm} (L[0]|i[l]j[m]), where L is the auxiliary basis. Three-index integral tensor (kptij_idx, naux, nao_pair) or four-index integral tensor (kptij_idx, comp, naux, nao_pair) are stored on disk.
- Args:
- kptij_lst(*,2,3) array
A list of (kpti, kptj)
pyscf.pbc.df.rsdf module#
Range-separated Gaussian Density Fitting (RSGDF)
rsdf_builder.py uses a different algorithm to build RS-GDF tensors. Note both modules CANNOT handle the long-range operator erf(omega*r12)/r12. It has to be computed with the gdf_builder module.
ref.: [1] For the RSGDF method:
Hong-Zhou Ye and Timothy C. Berkelbach, J. Chem. Phys. 154, 131104 (2021).
- [2] For the SR lattice sum integral screening:
Hong-Zhou Ye and Timothy C. Berkelbach, arXiv:2107.09704.
- In RSGDF, the two-center and three-center Coulomb integrals are calculated in two pars:
j2c = j2c_SR(omega) + j2c_LR(omega) j3c = j3c_SR(omega) + j3c_LR(omega)
- where the SR and LR integrals correspond to using the following potentials
g_SR(r_12;omega) = erfc(omega * r_12) / r_12 g_LR(r_12;omega) = erf(omega * r_12) / r_12
The SR integrals are evaluated in real space using a lattice summation, while the LR integrals are evaluated in reciprocal space with a plane wave basis.
pyscf.pbc.df.rsdf_builder module#
Build GDF tensor with range-separation technique.
rsdf.py is another version of RS-GDF using a different algorithm to compute 3c2e integrals. Note both modules CANNOT handle the long-range operator erf(omega*r12)/r12. It has to be computed with the gdf_builder module.
- Ref:
Sun, arXiv:2012.07929
- pyscf.pbc.df.rsdf_builder.estimate_ft_rcut(rs_cell, precision=None, exclude_dd_block=False)[source]#
Remove less important basis based on Schwarz inequality Q_ij ~ S_ij * (sqrt(2aij/pi) * aij**(lij*2) * (4*lij-1)!!)**.5
- pyscf.pbc.df.rsdf_builder.estimate_ke_cutoff_for_omega(cell, omega, precision=None)[source]#
Energy cutoff for AFTDF to converge attenuated Coulomb in moment space
- pyscf.pbc.df.rsdf_builder.estimate_omega_for_ke_cutoff(cell, ke_cutoff, precision=None)[source]#
The minimal omega in attenuated Coulomb given energy cutoff
- pyscf.pbc.df.rsdf_builder.estimate_omega_min(cell, precision=None)[source]#
Given cell.rcut the boundary of repeated images of the cell, estimates the minimal omega for the attenuated Coulomb interactions, requiring that at boundary the Coulomb potential of a point charge < cutoff
- pyscf.pbc.df.rsdf_builder.estimate_rcut(rs_cell, rs_auxcell, omega, precision=None, exclude_dd_block=False, exclude_d_aux=False)[source]#
Estimate rcut for 3c2e SR-integrals
pyscf.pbc.df.rsdf_helper module#
- class pyscf.pbc.df.rsdf_helper.MoleNoBasSort(**kwargs)[source]#
Bases:
Mole
- build(**kwargs)[source]#
Setup molecule and initialize some control parameters. Whenever you change the value of the attributes of
Mole
, you need call this function to refresh the internal data of Mole.- Kwargs:
- dump_inputbool
whether to dump the contents of input file in the output file
- parse_argbool
whether to read the sys.argv and overwrite the relevant parameters
- verboseint
Print level. If given, overwrite
Mole.verbose
- outputstr or None
Output file. If given, overwrite
Mole.output
- max_memoryint, float
Allowd memory in MB. If given, overwrite
Mole.max_memory
- atomlist or str
To define molecular structure.
- basisdict or str
To define basis set.
- nucmoddict or str
Nuclear model. If given, overwrite
Mole.nucmod
- chargeint
Charge of molecule. It affects the electron numbers If given, overwrite
Mole.charge
- spinint
2S, num. alpha electrons - num. beta electrons to control multiplicity. If setting spin = None , multiplicity will be guessed based on the neutral molecule. If given, overwrite
Mole.spin
- symmetrybool or str
Whether to use symmetry. If given a string of point group name, the given point group symmetry will be used.
- magmomlist
Collinear spin of each atom. Default is [0.0,]*natm
- pyscf.pbc.df.rsdf_helper.estimate_mesh_for_omega(cell, omega, precision=None, kmax=0, round2odd=True)[source]#
Find the appropriate PW mesh size s.t. the LR integrals computed via AFT achieves a given precision.
- Args:
See Args for :func:’estimate_omega_for_npw’.
- pyscf.pbc.df.rsdf_helper.estimate_omega_for_npw(cell, npw_max, precision=None, kmax=0, round2odd=True)[source]#
Find the biggest omega where the appropriate PW mesh for achieving a given precision in the LR integrals computed using AFT does not exceed npw_max in size.
- Args:
cell (PBC Cell object) npw_max (int):
- The returned mesh will stasify
mesh[0]*mesh[1]*mesh[2] <= npw_max
- precision (float; default: None):
Target precision for the integrals. If not provided, cell.precision is used.
- kmax (float; default: 0):
‘Gmin’ in eqn (28) of https://arxiv.org/pdf/2012.07929.pdf See :func:’RSDF._rs_build’ for how to determine it for given a kpt mesh; note the use of a minimum kmax there to improve accuracy for Gamma point.
- round2odd (bool; default: True):
If True, the mesh is required to be odd in each dimension.
- Returns:
omega (float), ke_cutoff (float), mesh (array of size 3).
- pyscf.pbc.df.rsdf_helper.intor_j2c(cell, omega, precision=None, kpts=None, hermi=1, shls_slice=None, no_screening=False)[source]#
Calculate the SR 2c integrals (i| erfc(omega*r12)/r12 |j) via a real- space lattice sum.
- Args:
cell (PBC Cell object) omega (float):
The omega in the SR Coulomb potential (only its absolute value will be used).
- precision (float; default: None):
Target precision for the integrals. If not provided, cell.precision is used.
- kpts (array; default: None):
The kpoint mesh used to sample the Brillouin zone. If not provided, Gamma point is used.
- hermi (int; default: 1):
Hermitian symmetry: 0 for s1, and 1 for s2.
- shls_slice (array of size 4; default: None):
Range of shl indices organized as (ish0, ish1, jsh0, jsh1). If not provided, i/jsh0 = 0 and i/jsh1 = cell.nbas.
- no_screening (bool; default: False):
If set to True, the integral bound is only used in determining a global rcut for the lattice sum and no further (i.e., shlpair-wise) screening is performed. This arg is for debug purpose.
- Returns:
If a single kpt, return the matrix of integrals in the GTO basis. Otherwise, return a list of matrices corresponding to different kpts.
- pyscf.pbc.df.rsdf_helper.remove_exp_basis(basis, amin=None, amax=None)[source]#
Removing primitives with exponents not in (amin,amax).
- Args:
- basis (list or dict):
- PySCF GTO basis format (as in mol/cell._basis), e.g.,:
basis = [[0, (1.5, 1.)], [1, (0.7, 1.)]] basis = {“C”: [[0, (1.5, 1.)], [1, (0.7, 1.)]],
“H”: [[0, (0.3, 1.)]]}
- amin/amax (float; default: None):
Primitives with exponents <= amin or >= amax will be removed. If either is not provided (default), no primitives will be removed by that criterion.
pyscf.pbc.df.rsdf_jk module#
- pyscf.pbc.df.rsdf_jk.density_fit(mf, auxbasis=None, mesh=None, with_df=None)[source]#
Generate density-fitting SCF object
- Args:
- auxbasisstr or basis dict
Same format to the input attribute mol.basis. If auxbasis is None, auxiliary basis based on AO basis (if possible) or even-tempered Gaussian basis will be used.
- meshtuple
number of grids in each direction
with_df : DF object