pyscf.sgx package#
Subpackages#
Submodules#
pyscf.sgx.sgx module#
Pseudo-spectral methods (COSX, PS, SN-K)
- class pyscf.sgx.sgx.SGX(mol, auxbasis=None)[source]#
Bases:
StreamObjectObject to perform seminumerical evaluation of J/K matrix.
- Basic attributes:
- grids_thrdfloat
Remove grids with small weights using this threshold. Default is 1e-10.
- grids_level_iint
Initial grids level (coarse grid for early SCF iterations). Default is 2.
- grids_level_fint
Final grids level (dense grid for final SCF iterations). Default is 2.
- use_opt_gridsbool
Use optimized grids for SGX based on ORCA. Default is True.
- fit_ovlpbool
Numerically fit overlap matrix to improve numerical precision. Default is True.
- grids_switch_thrdfloat
Threshold of density matrix convergence to switch from the coarse initial grid (grids_level_i) to the denser final grid (grids_level_f). Ignored if grids_level_i == grids_level_f. Default is 0.03.
- dfjbool
Compute J matrix using DF and K matrix using SGX, like in the RIJCOSX method in ORCA. Default is True.
- Optimized SGX-K Attributes:
- optkbool
Turn on optimization for evaluating the K-matrix with SGX. Only has an effect if dfj is True. Default is True.
- sgx_tol_energyfloat or str
DM screening error tolerance for the energy. “auto” means set automatically using direct_scf_tol. Default=”auto”.
- sgx_tol_potentialfloat or str
DM screening error tolerance for the potential (K-matrix). “auto” means set to sqrt(sgx_tol_energy), or sqrt(direct_scf_tol) if sgx_tol_energy is None. Default=”auto”.
- bound_algostr
Bound algo determines how the three-center integral upper bounds are estimated. Default is “sample_pos”.
- If (sgx_tol_energy, sgx_tol_potential) is
(float/”auto”, float/”auto”): Bound energy and potential error; (float/”auto”, None): Bound energy error only; (None, float/”auto”): Bound potential error only; (None, None): Turn off DM screening (no energy or potential error).
It is recommended to bound both energy and potential error for numerical stability.
- Attribute bound_algo can be
- “ovlp”:
Screen integrals based on overlap of orbital pairs. Overlap serves as a rough approximation of the maximum ESP integral.
- “sample”:
Provide an approximate but accurate upper bound for the ESP integrals by sampling _nquad points for each shell pair.
- “sample_pos”:
Same as sample, but the ESP bounds are position-dependent, which gives a slight speed increase for large systems and a significant speed increase for short-range hybrids.
Default is “sample_pos” and is recommended for most cases.
- Debugging Attributes:
- debugbool
Whether to use debug mode. Default is False.
- blockdimint
Max block size for grids when debug=True. Default 1200.
- direct_jbool
Run the calculation with direct integral J-matrix. Default False.
- _symm_ovlp_fit: bool
Default True. Perform a “symmetric” overlap fit when optk is True. When fit_ovlp=True, _symm_ovlp_fit=True is required for exact analytical gradients. Note that symmetric overlap fitting is not “exact” overlap fitting. It fits the overlap correctly to first order in the difference between the exact and numerical overlap matrix. This difference is quite small for any reasonable grid, so this approximation typically works well. Set to False (not typically recommended) to obtain “exact” overlap fitting and abandon exact analytical gradients.
- property auxbasis#
- get_jk(dm, hermi=1, vhfopt=None, with_j=True, with_k=True, direct_scf_tol=1e-13, omega=None)[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.).
- 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.sgx.sgx.sgx_fit(mf, auxbasis=None, with_df=None, pjs=False)[source]#
For the given SCF object, update the J, K matrix constructor with corresponding SGX or 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.
- with_dfSGX
Existing SGX object for the system.
- pjsbool
If True, the SGX object is set up to screen negligible integrals using the density matrix (i.e. P-junction screening), and density fitting is used for the J-matrix. If False, no P-junction screening is performed, and SGX is used for the J-matrix. Screening settings can also be adjusted after initialization. See dfj, optk, sgx_tol_energy, and sgx_tol_potential for details.
- Returns:
An SCF object with a modified J, K matrix constructor which uses density fitting and seminumerical integrals to compute J and K
Examples:
>>> mol = gto.M(atom='H 0 0 0; F 0 0 1', basis='ccpvdz', verbose=0) >>> mf = sgx_fit(scf.RHF(mol)) >>> mf.scf() -100.00978770917165
>>> mol.symmetry = 1 >>> mol.build(0, 0) >>> mf = sgx_fit(scf.UHF(mol)) >>> mf.scf() -100.00978770951018
pyscf.sgx.sgx_jk module#
semi-grid Coulomb and eXchange without differential density matrix
To lower the scaling of coulomb and exchange matrix construction for large system, one coordinate is analytical and the other is grid. The traditional two electron integrals turn to analytical one electron integrals and numerical integration based on grid.(see Friesner, R. A. Chem. Phys. Lett. 1985, 116, 39)
Minimizing numerical errors using overlap fitting correction.(see Lzsak, R. et. al. J. Chem. Phys. 2011, 135, 144105) Grid screening for weighted AO value and DktXkg. Two SCF steps: coarse grid then fine grid.
See the docs of pyscf.sgx.sgx for details on adjustable parameters.
- class pyscf.sgx.sgx_jk.SGXData(mol, grids, fit_ovlp=True, sym_ovlp=False, max_memory=2000, hermi=0, bound_algo='sample_pos', sgxopt=None, etol=None, vtol=None, direct_scf_tol=0.0)[source]#
Bases:
object- get_loop_data(nset=1, with_pair_mask=True, grad=False)[source]#
Get all the info we need to loop over grids, i.e. nbins, screen_index (mask), pair_mask, ao_loc, and the block size.
- property grids#
- property mol#
- reduce_ao_to_shl_mat(ao, ao_loc, wt)[source]#
For F-ordered orbitals ao[g, u] = X_ug, compute the sum of X_ug^2 * wt_g within each grid block and AO shell.
- reduce_mask_ni2sgx(xni_bi)[source]#
Compute xsgx_bi (block-reduced array) in blocks of SGX_BLKSIZE, using sni_bi in blocks of BLKSIZE
- set_block(b0, weights, full_f_bi=None)[source]#
Se the grid block information in the C sgxopt object. This must be called right before SGXnr_direct_k_drv. TODO might be better to hav a safer way of managing this data.
- Args:
b0: Initial block idx weights: grid weights for this set of blocks full_f_bi: block-reduced estimate (upper bound)
for the full density matrix f_gu. In direct_scf, f_gu passed directly to SGXnr_direct_k_drv is a DM difference, and we also need the “full” DM for energy-based screening. If None, we assume DM difference and full DM are the same.
- unset_block()[source]#
Reset the grid block data a the C level. Should be called after a call to SGXnr_direct_k_drv.
- property use_dm_screening#
- pyscf.sgx.sgx_jk.get_jk(sgx, dm, hermi=1, with_j=True, with_k=True, direct_scf_tol=1e-13)#
- pyscf.sgx.sgx_jk.get_jk_favorj(sgx, dm, hermi=1, with_j=True, with_k=True, direct_scf_tol=1e-13)[source]#