Skip to content

Latest commit

 

History

History
147 lines (95 loc) · 6.08 KB

df.rst

File metadata and controls

147 lines (95 loc) · 6.08 KB

Density fitting

Modules: :mod:`df`, :mod:`pbc.df`

.. module:: df
   :synopsis: Density fitting and RI approximation
.. sectionauthor:: Qiming Sun <[email protected]>.

Introduction

Density fitting (DF), sometimes also called the resolution of identity (RI) approximation, is a method to approximate the four-index electron repulsion integrals (ERIs) by two- and three-index tensors. :cite:`Dunlap00PCCP` It achieves so by expanding the atomic orbital (AO) pair density in a second, auxiliary Gaussian basis set. PySCF has extensive supports for using DF in both molecular and crystalline calculations. This section of the user document covers the basic usage of DF in molecular applications. See :ref:`theory_pbc_df` for the use of DF in periodic calculations.

Using DF in calculations

DF is not used by default when initializing a :class:`SCF` or :class:`MCSCF` object, but can be invoked via the :func:`density_fit` method as follows:

from pyscf import gto, scf, mcscf
mol = gto.M(atom='N 0 0 0; N 0 0 1.2', basis='def2-tzvp')
mf = scf.RHF(mol).density_fit().run()
mc = mcscf.CASSCF(mf, 8, 10).density_fit().run()

Alternatively, one can initialize a :class:`DF` object and overwrite the :attr:`with_df` attribute of the :class:`SCF` object:

from pyscf import gto, scf, df
mol = gto.M(atom='N 0 0 0; N 0 0 1.2', basis='def2-tzvp')
mydf = df.DF(mol).build()
mf = scf.RHF(mol).density_fit()
mf.with_df = mydf
mf.kernel()

Setting :attr:`with_df` to None turns DF off.

Once the SCF calculations using DF are done, all post-SCF methods will use DF automatically. The example below performs a TDDFT/TDA calculation using DF:

from pyscf import gto, dft, tddft
mol = gto.M(atom='N 0 0 0; N 0 0 1.2', basis='def2-tzvp')
mf = dft.RKS(mol, xc="b3lyp").density_fit().run()
td = tddft.TDA(mf).run()
print(td.e)

Choice of auxiliary basis

General consideration

The choice of auxiliary basis depends on the AO basis. By default, PySCF uses the pre-defined auxiliary basis set optimized for the AO basis set if the former exists. This includes many commonly used AO basis sets in electronic structure calculations, e.g., the Ahlrichs' def2 family, :cite:`Hellweg07TCA,Weigend98CPL` the Dunning's cc family, :cite:`Weigend02PCCP` etc.

When pre-defined auxiliary basis set is not available, an even-tempered basis (ETB) set is generated by the following rule

\varphi &= r^l \exp(-\zeta_{il} r^2), \quad i = 0, 1, \cdots, n \\
\zeta_{il} &= \alpha_l \times \beta^i

where \alpha_l and N are inferred automatically from the AO basis, and \beta = 2.0 by default.

The user can overwrite the default choice by setting :attr:`auxbasis` upon initialization or at a later stage:

from pyscf import gto, scf, df
mol = gto.M(atom='N 0 0 0; N 0 0 1.2', basis='cc-pvdz')
mf = scf.RHF(mol).density_fit(auxbasis="weigend")
mf.kernel() # -108.910953335055
mf.with_df.auxbasis = "cc-pvdz-jkfit"   # this is the default
mf.kernel() # -108.913710743723
mf.with_df.auxbasis = df.aug_etb(mol, beta=1.7) # ETB with beta = 1.7
mf.kernel() # -108.914059329528

The user can check what auxiliary basis set is used at any time by:

print(mf.with_df.auxmol.basis)

More examples on using customized auxiliary basis can be found in :source:`examples/df/01-auxbasis.py`.

Special consideration for DFT

For DFT calculations with pure exchange-correlation functionals (i.e., LDA and GGA), the default auxiliary basis, which is designed for fitting both the Coulomb and the Hartree-Fock exchange integrals, may be unnecessarily too large. We recommend using the weigend family :cite:`Weigend06PCCP` for a more cost-effective choice as shown in the following example:

from pyscf import gto, dft
mol = gto.M(atom='N 0 0 0; N 0 0 1.2', basis='def2-tzvpp')
mf = dft.RKS(mol, xc="pbe").density_fit().run() # -109.432329411505
print(mf.with_df.auxmol.basis)                  # default: def2-tzvpp-jkfit
print(mf.with_df.auxmol.nao_nr())               # 154 aux basis functions
mf.with_df.auxbasis = "weigend"
mf.kernel()                                     # -109.432334646585
print(mf.with_df.auxmol.nao_nr())               # 98 aux basis functions

Saving and reusing DF integrals

By default, the DF integrals are discarded after the calculations are done. Sometimes it is beneficial to save those integrals on disk and re-use them for later calculations. This can be done by specifying a HDF5 file where the DF integrals will be saved via setting :attr:`_cderi_to_save` either at the SCF stage:

mf = scf.RHF(mol).density_fit()
mf.with_df._cderi_to_save = 'saved_cderi.h5'
mf.kernel()

or initializing a :class:`DF` object separately:

mydf = df.DF(mol)
mydf.auxbasis = df.make_auxbasis(mol)
mydf._cderi_to_save = 'saved_cderi.h5'
mydf.build()

The save integrals can be used later by setting :attr:`_cderi` to the HDF5 file holding the DF integrals:

mf = scf.RHF(mol).density_fit()
mf.with_df._cderi = 'saved_cderi.h5'
mf.kernel()

More examples on saving and using DF integrals can be found in :source:`examples/df/10-access_df_integrals.py`, :source:`examples/df/11-access_df_tensor.py`, and :source:`examples/df/40-precompute_df_integrals.py`.

Advanced examples

More examples on advanced topic of using the :mod:`df` module include

References

.. bibliography:: ref_df.bib
   :style: unsrt