HFMethod
The HFMethod
object specifies the settings for using Hartree-Fock (HF) theory as the method for a calculation.
There is one field required for creating a HFMethod
object:
ao
for specifying the atomic-orbital basis set
All other calculation details such as SCF convergence criteria can be specified in the details
dictionary.
Example
The following examples demonstrate how to create and use HFMethod
.
import sierra
from sierra.inputs import *
# Create a water molecule from pubchem (internet connection required)
water = Molecule(pubchem="water")
# Perform a single-point energy calculation at HF/Def2-SVP level of theory
# using default settings for the other calculation details
basic_input = EnergyInput(
molecule=water,
method=HFMethod(ao="def2-svp"),
)
basic_result = sierra.run(basic_input)
print(f"HF/Def2-SVP energy: {basic_result.energy:.6f} Hartree")
#> HF/Def2-SVP energy: -75.960068 Hartree
# Use more advanced settings for the HF/Def2-SVP method
advanced_input = EnergyInput(
molecule=water,
method=HFMethod(
ao="def2-svp", details={"max_iter": 200, "energy_threshold": 1e-8}
),
)
advanced_result = sierra.run(advanced_input)
print(f"HF/Def2-SVP energy: {advanced_result.energy:.6f} Hartree")
#> HF/Def2-SVP energy: -75.960068 Hartree
Fields
ao
-
Atomic-orbital basis set name.
- Type: str
engine
-
The Entos Engine in which to evaluate the method with.
- Type: EngineIdentifier
- Default: EngineIdentifier(name='qcore', image=None)
kind
-
- Type: Literal['hf', 'HFMethod']
- Default: 'hf'
solvation
-
Include effects of solvation through a continuum model.
- Type: Optional[Solvation]
details
- Details relating to the given method
Available options for details
max_iter
-
Maximum number of SCF iterations.
- The type is int
- The default is 128
- The value must be nonnegative
energy_threshold
-
SCF convergence threshold using the absolute value of the energy difference between iterations.
- The type is float
- The default is 1e-6
- The value must be nonnegative
orbital_grad_threshold
-
SCF convergence threshold using the infinity norm of the orbital gradient.
- The type is float
- The default is 1e-5
- The value must be nonnegative
coulomb_method
-
Method for computing the Coulomb contribution to Fock matrix and energy.
- The type is str
- The default is incore_df
- The value must be one of:
incore_df
-Use conventional incore density-fitting to compute the Coulomb contribution to Fock matrix and energy.
This algorithm is fast but has a memory requirement of about \(8 * N_\text{ao}^2 * N_\text{df}\) Bytes, where \(N_\text{ao}\) and \(N_\text{df}\) are the number of AO and density-fitting basis functions, respectively. Therefore, it may not work for calculations on very large systems or with very large basis; in these cases, use
coulomb_method = direct_df
instead.direct_df
-Use integral-direct density-fitting to compute the Coulomb contribution to Fock matrix and energy.
This algorithm has minimal memory requirement but is generally slower than the
incore_df
method. Use it whencoulomb_method = incore_df
does not work.This method also indicates use of integral screening (see schwarz_threshold) and incremental Fock build (see incremental_fock).
direct_4idx
-Use integral-direct approach to compute the 4-index electron repulsion integrals (ERIs).
Currently, this algorithm is very slow, and is recommended to be used only for testing purposes (on small systems with small basis).
exchange_method
-
Method for computing the exact exchange contribution to Fock matrix and energy.
- The type is str
- The default is pre_transformed_df
- The value must be one of:
incore_df
-Use conventional incore density-fitting to compute the exchange contribution to Fock matrix and energy.
This algorithm is fast but has a memory requirement of about \(8 * N_\text{ao}^2 * N_\text{df}\) Bytes, where \(N_\text{ao}\) and \(N_\text{df}\) are the number of AO and density-fitting basis functions, respectively. Therefore, it may not work for calculations on very large systems or with very large basis; in these cases, use
coulomb_method = direct_df
instead.pre_transformed_df
-Use the pre-transformed density-fitting integrals to compute the exact exchange contribution to Fock matrix and energy.
This algorithm is faster than
incore_df
when SCF converges slowly (ca. with more than 10 SCF steps), and the speepup increases with the number of SCF steps; but is slower when SCF converges within a few steps.This algorithm has the same memory requirement as
incore_df
.occupied_df
-Use the occ-RI-K algorithm (incore version) to compute the exact exchange contribution to Fock matrix and energy.
This algorithm is ~ 2 times faster than
incore_df
for normal calculations (e.g. with double-zeta AO basis), and the speedup increases when larger basis is used.See this paper for a detailed description of the algorithm.
This algorithm has the same memory requirement as
incore_df
.
level_shift
-
Turn on level shifting to improve SCF convergence. Raises the energy of virtual orbitals by level shifting value. The level shift is removed at the end of the calculation.
- The type is float
- The default is 0
- The value must be nonnegative
temperature
-
Specify the temperature for the electrons and perform the SCF calculation in the canonical-ensemble. Convergence is often poor with the default convergence settings, and using the original Pulay DIIS method is recommended (see diis option).
- The type is a temperature quantity, i.e. a string consisting of a number and a temperature unit.
- The default is 0 kelvin
diis
-
Method of DIIS.
- The type is str
- The default is adiis+cdiis
- The value must be one of:
none
- DIIS is not used.cdiis
- This is the original DIIS algorithm by P. Pulay. P. Pulay, Chem. Phys. Lett. 73, 393 (1980).ediis
- Energy-based DIIS which minimizes the Hartree-Fock energy functional. Kudin, Scuseria and Cances, J. Chem. Phys. 116, 8255 (2002).adiis
- DIIS based on the augmented Roothaan--Hall energy. X. Hu and W. Yang, J. Chem. Phys. 132, 054109 (2010).ediis+cdiis
- Combination of these two algorithms; see diis_switch_threshold_gradadiis+cdiis
- Combination of these two algorithms; see diis_switch_threshold_grad
List of avaliable basis sets
- 3-21++G
- 3-21G
- 3-21GSP
- 4-31G
- 6-31++G
- 6-31++G*
- 6-31++G**
- 6-31+G
- 6-31+G*
- 6-31+G**
- 6-311++G
- 6-311++G(2d,2p)
- 6-311++G(3df,3pd)
- 6-311++G*
- 6-311++G**
- 6-311+G
- 6-311+G(2d,p)
- 6-311+G*
- 6-311+G**
- 6-311G
- 6-311G(2df,2pd)
- 6-311G*
- 6-311G**
- 6-31G
- 6-31G(2df,p)
- 6-31G(3df,3pd)
- 6-31G*
- 6-31G**
- DZP
- DZP-DKH
- Def2-QZVP
- Def2-QZVPD
- Def2-QZVPP
- Def2-QZVPPD
- Def2-SVP
- Def2-TZVP
- Def2-TZVPD
- Def2-TZVPP
- Def2-TZVPPD
- STO-2G
- STO-3G
- STO-6G
- aug-cc-pVDZ
- aug-cc-pVTZ
- aug-pc-0
- aug-pc-1
- aug-pc-2
- aug-pc-3
- aug-pcseg-0
- aug-pcseg-1
- aug-pcseg-2
- aug-pcseg-3
- cc-pVDZ
- cc-pVQZ
- cc-pVTZ
- heavy-aug-cc-pVDZ
- heavy-aug-cc-pVTZ
- mTZVP
- pc-0
- pc-1
- pc-2
- pc-3
- pcseg-0
- pcseg-1
- pcseg-2
- pcseg-3