pyscf.pbc.dft.multigrid package#
Submodules#
pyscf.pbc.dft.multigrid.multigrid module#
Multigrid to compute DFT integrals
- class pyscf.pbc.dft.multigrid.multigrid.MultiGridFFTDF(cell, kpts=array([[0., 0., 0.]]))[source]#
Bases:
FFTDF
- get_jk(dm, hermi=1, kpts=None, kpts_band=None, with_j=True, with_k=True, exxdiv='ewald', **kwargs)[source]#
- get_nuc(kpts=None)#
- get_pp(kpts=None, max_memory=4000)#
Get the periodic pseudopotential nuc-el AO matrix, with G=0 removed.
- get_rho(dm, kpts=array([[0., 0., 0.]]))#
Density in real space
- pyscf.pbc.dft.multigrid.multigrid.cache_xc_kernel(mydf, xc_code, mo_coeff, mo_occ, spin=0, kpts=None)[source]#
- pyscf.pbc.dft.multigrid.multigrid.cache_xc_kernel1(mydf, xc_code, dm, spin=0, kpts=None)[source]#
Compute the 0th order density, Vxc and fxc. They can be used in TDDFT, DFT hessian module etc.
- pyscf.pbc.dft.multigrid.multigrid.eval_mat(cell, weights, shls_slice=None, comp=1, hermi=0, xctype='LDA', kpts=None, mesh=None, offset=None, submesh=None)[source]#
- pyscf.pbc.dft.multigrid.multigrid.eval_rho(cell, dm, shls_slice=None, hermi=0, xctype='LDA', kpts=None, mesh=None, offset=None, submesh=None, ignore_imag=False, out=None)[source]#
Collocate the real density (opt. gradients) on the real-space grid.
- Kwargs:
- ignore_image :
The output density is assumed to be real if ignore_imag=True.
- pyscf.pbc.dft.multigrid.multigrid.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.dft.multigrid.multigrid.get_pp(mydf, kpts=None, max_memory=4000)[source]#
Get the periodic pseudopotential nuc-el AO matrix, with G=0 removed.
- pyscf.pbc.dft.multigrid.multigrid.get_rho(mydf, dm, kpts=array([[0., 0., 0.]]))[source]#
Density in real space
- pyscf.pbc.dft.multigrid.multigrid.multi_grids_tasks_for_ke_cut(cell, fft_mesh=None, verbose=None)[source]#
- pyscf.pbc.dft.multigrid.multigrid.multi_grids_tasks_for_rcut(cell, fft_mesh=None, verbose=None)[source]#
- pyscf.pbc.dft.multigrid.multigrid.multigrid(mf)#
Use MultiGridFFTDF to replace the default FFTDF integration method in the DFT object.
- pyscf.pbc.dft.multigrid.multigrid.multigrid_fftdf(mf)[source]#
Use MultiGridFFTDF to replace the default FFTDF integration method in the DFT object.
- pyscf.pbc.dft.multigrid.multigrid.nr_rks(mydf, xc_code, dm_kpts, hermi=1, kpts=None, kpts_band=None, with_j=False, return_j=False, verbose=None)[source]#
Compute the XC energy and RKS XC matrix at sampled k-points. multigrid version of function pbc.dft.numint.nr_rks.
- Args:
- dm_kpts(nkpts, nao, nao) ndarray or a list of (nkpts,nao,nao) ndarray
Density matrix at each k-point.
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:
exc : XC energy nelec : number of electrons obtained from the numerical integration veff : (nkpts, nao, nao) ndarray
or list of veff if the input dm_kpts is a list of DMs
- vj(nkpts, nao, nao) ndarray
or list of vj if the input dm_kpts is a list of DMs
- pyscf.pbc.dft.multigrid.multigrid.nr_rks_fxc(mydf, xc_code, dm0, dms, hermi=0, with_j=False, rho0=None, vxc=None, fxc=None, kpts=None, verbose=None)[source]#
multigrid version of function pbc.dft.numint.nr_rks_fxc
- pyscf.pbc.dft.multigrid.multigrid.nr_rks_fxc_st(mydf, xc_code, dm0, dms_alpha, singlet=True, rho0=None, vxc=None, fxc=None, kpts=None, verbose=None)[source]#
multigrid version of function pbc.dft.numint.nr_rks_fxc_st
- pyscf.pbc.dft.multigrid.multigrid.nr_uks(mydf, xc_code, dm_kpts, hermi=1, kpts=None, kpts_band=None, with_j=False, return_j=False, verbose=None)[source]#
Compute the XC energy and UKS XC matrix at sampled k-points. multigrid version of function pbc.dft.numint.nr_uks
- Args:
- dm_kpts(nkpts, nao, nao) ndarray or a list of (nkpts,nao,nao) ndarray
Density matrix at each k-point.
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:
exc : XC energy nelec : number of electrons obtained from the numerical integration veff : (2, nkpts, nao, nao) ndarray
or list of veff if the input dm_kpts is a list of DMs
- vj(nkpts, nao, nao) ndarray
or list of vj if the input dm_kpts is a list of DMs
pyscf.pbc.dft.multigrid.multigrid_pair module#
- class pyscf.pbc.dft.multigrid.multigrid_pair.GridLevel_Info[source]#
Bases:
Structure
Info about the grid levels.
- cutoff#
Structure/Union member
- mesh#
Structure/Union member
- nlevels#
Structure/Union member
- rel_cutoff#
Structure/Union member
- class pyscf.pbc.dft.multigrid.multigrid_pair.MultiGridFFTDF2(cell, kpts=array([[0., 0., 0.]]))[source]#
Bases:
MultiGridFFTDF
Base class for multigrid DFT (version 2).
- Attributes:
- task_listTaskList instance
Task list recording which primitive basis function pairs need to be considered.
- vpplocG_part1arrary
Short-range part of the local pseudopotential represented in the reciprocal space. It is cached to reduce cost.
- rhoGarray
Electronic density represented in the reciprocal space. It is cached in nuclear gradient calculations to reduce cost.
- get_pp(kpts=None)[source]#
Compute the GTH pseudopotential matrix, which includes the second part of the local potential and the non-local potential. The first part of the local potential is cached as vpplocG_part1, which is the reciprocal space representation, to be added to the electron density for computing the Coulomb matrix. In order to get the full PP matrix, the potential due to vpplocG_part1 needs to be added.
- ke_ratio = 3.0#
- ngrids = 4#
- rel_cutoff = 20.0#
- vpploc_part1_nuc_grad(dm, kpts=array([[0., 0., 0.]]), atm_id=None, precision=None)#
- class pyscf.pbc.dft.multigrid.multigrid_pair.PGFPair[source]#
Bases:
Structure
A primitive Gaussian function pair.
- iL#
Structure/Union member
- ipgf#
Structure/Union member
- ish#
Structure/Union member
- jpgf#
Structure/Union member
- jsh#
Structure/Union member
- radius#
Structure/Union member
- class pyscf.pbc.dft.multigrid.multigrid_pair.RS_Grid[source]#
Bases:
Structure
Values on real space multigrid.
- comp#
Structure/Union member
- data#
Structure/Union member
- gridlevel_info#
Structure/Union member
- nlevels#
Structure/Union member
- class pyscf.pbc.dft.multigrid.multigrid_pair.Task[source]#
Bases:
Structure
A single task.
- buf_size#
Structure/Union member
- ntasks#
Structure/Union member
- pgfpairs#
Structure/Union member
- radius#
Structure/Union member
- class pyscf.pbc.dft.multigrid.multigrid_pair.TaskList[source]#
Bases:
Structure
A task list.
- gridlevel_info#
Structure/Union member
- hermi#
Structure/Union member
- nlevels#
Structure/Union member
- tasks#
Structure/Union member
- pyscf.pbc.dft.multigrid.multigrid_pair.build_task_list(cell, gridlevel_info, cell1=None, Ls=None, hermi=0, precision=None)[source]#
Build the task list for multigrid DFT calculations.
- Arguments:
- cell
pbc.gto.cell.Cell
The
Cell
instance for the bra basis functions.- gridlevel_info
ctypes.POINTER
The C pointer of the
GridLevel_Info
structure.- cell1
pbc.gto.cell.Cell
, optional The
Cell
instance for the ket basis functions. If not given, both bra and ket basis functions come from cell.- Ls(*,3) array, optional
The cartesian coordinates of the periodic images. Default is calculated by
cell.get_lattice_Ls()
.- hermiint, optional
If \(hermi=1\), the task list is built only for the upper triangle of the matrix. Default is 0.
- precisionfloat, optional
The integral precision. Default is
cell.precision
.
- cell
- Returns:
ctypes.POINTER
The C pointer of the
TaskList
structure.
- pyscf.pbc.dft.multigrid.multigrid_pair.eval_mat(cell, weights, task_list, shls_slice=None, comp=1, hermi=0, deriv=0, xctype='LDA', kpts=None, grid_level=None, dimension=None, mesh=None, cell1=None, shls_slice1=None, Ls=None, a=None)[source]#
- pyscf.pbc.dft.multigrid.multigrid_pair.eval_rho(cell, dm, task_list, shls_slice=None, hermi=0, xctype='LDA', kpts=None, dimension=None, cell1=None, shls_slice1=None, Ls=None, a=None, ignore_imag=False)[source]#
Collocate density (opt. gradients) on the real-space grid. The two sets of Gaussian functions can be different.
- Returns:
- rho: RS_Grid object
Densities on real space multigrids.
- pyscf.pbc.dft.multigrid.multigrid_pair.free_task_list(task_list)[source]#
- Note:
This will also free task_list.contents.gridlevel_info.
- pyscf.pbc.dft.multigrid.multigrid_pair.get_veff_ip1(mydf, dm_kpts, xc_code=None, kpts=array([[0., 0., 0.]]), kpts_band=None, spin=0)[source]#
- pyscf.pbc.dft.multigrid.multigrid_pair.gradient_gs(f_gs, Gv)[source]#
Compute the G-space components of \(\nabla f(r)\) given \(f(G)\) and \(G\), which is equivalent to einsum(‘np,px->nxp’, f_gs, 1j*Gv)
- pyscf.pbc.dft.multigrid.multigrid_pair.init_rs_grid(gridlevel_info, comp)[source]#
Initialize values on real space multigrid
- pyscf.pbc.dft.multigrid.multigrid_pair.multi_grids_tasks(cell, ke_cutoff=None, hermi=0, ngrids=4, ke_ratio=3.0, rel_cutoff=20.0)[source]#