pyscf.pbc.dft package

Submodules

pyscf.pbc.dft.cdft module

pyscf.pbc.dft.cdft.cdft(mf, cell, offset, orbital, basis=None)
Input:

mf – a mean field object for DFT or (in principle) HF (doesn’t really matter) shift – float – a semi aribitrary energy which displaces the selected orbitals by the diagonal orbital – int – indicating which orbital are shifted in the selected basis basis – 2D numpy array – the working basis in the basis of AOs from ‘cell’ (Defaults to AO basis)

Returns:

mf – converged mean field object (with AO basis)

pyscf.pbc.dft.cdft.fast_iao_mullikan_pop(mf, cell, a=None)

Input: mf – a preconverged mean fild object Returns: mullikan populaion analysis in the basisIAO a

pyscf.pbc.dft.gen_grid module

pyscf.pbc.dft.gen_grid.AtomicGrids

alias of pyscf.pbc.dft.gen_grid.BeckeGrids

class pyscf.pbc.dft.gen_grid.BeckeGrids(cell)

Bases: pyscf.dft.gen_grid.Grids

Atomic grids for all-electron calculation.

build(cell=None, with_non0tab=False)
make_mask(cell=None, coords=None, relativity=0, shls_slice=None, verbose=None)

Mask to indicate whether a shell is zero on grid. The resultant mask array is an extension to the mask array used in molecular code (see also pyscf.dft.numint.make_mask function). For given shell ID and block ID, the value of the extended mask array means the number of images in Ls that does not vanish.

class pyscf.pbc.dft.gen_grid.UniformGrids(cell)

Bases: pyscf.lib.misc.StreamObject

Uniform Grid class.

build(cell=None, with_non0tab=False)
property coords
dump_flags(verbose=None)
kernel(cell=None, with_non0tab=False)

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.).

make_mask(cell=None, coords=None, relativity=0, shls_slice=None, verbose=None)

Mask to indicate whether a shell is zero on grid. The resultant mask array is an extension to the mask array used in molecular code (see also pyscf.dft.numint.make_mask function). For given shell ID and block ID, the value of the extended mask array means the number of images in Ls that does not vanish.

reset(cell=None)
property weights
pyscf.pbc.dft.gen_grid.gen_becke_grids(cell, atom_grid={}, radi_method=<function gauss_chebyshev>, level=3, prune=<function nwchem_prune>)

real-space grids using Becke scheme

Args:

cell : instance of Cell

Returns:
coords(N, 3) ndarray

The real-space grid point coordinates.

weights : (N) ndarray

pyscf.pbc.dft.gen_grid.get_becke_grids(cell, atom_grid={}, radi_method=<function gauss_chebyshev>, level=3, prune=<function nwchem_prune>)

real-space grids using Becke scheme

Args:

cell : instance of Cell

Returns:
coords(N, 3) ndarray

The real-space grid point coordinates.

weights : (N) ndarray

pyscf.pbc.dft.gen_grid.make_mask(cell, coords, relativity=0, shls_slice=None, verbose=None)

Mask to indicate whether a shell is zero on grid. The resultant mask array is an extension to the mask array used in molecular code (see also pyscf.dft.numint.make_mask function). For given shell ID and block ID, the value of the extended mask array means the number of images in Ls that does not vanish.

pyscf.pbc.dft.krks module

Non-relativistic Restricted Kohn-Sham for periodic systems with k-point sampling

See Also:
pyscf.pbc.dft.rks.pyNon-relativistic Restricted Kohn-Sham for periodic

systems at a single k-point

class pyscf.pbc.dft.krks.KRKS(cell, kpts=array([[0., 0., 0.]]), xc='LDA,VWN', exxdiv='ewald')

Bases: pyscf.pbc.dft.rks.KohnShamDFT, pyscf.pbc.scf.khf.KRHF

RKS class adapted for PBCs with k-point sampling.

density_fit(auxbasis=None, with_df=None, *args, **kwargs)
dump_flags(verbose=None)
energy_elec(dm_kpts=None, h1e_kpts=None, vhf=None)

Following pyscf.scf.hf.energy_elec()

get_rho(dm=None, grids=None, kpts=None)

Compute density in real space

get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)

Coulomb + XC functional

Note

This is a replica of pyscf.dft.rks.get_veff with kpts added. This function will change the ks object.

Args:
ksan instance of RKS

XC functional are controlled by ks.xc attribute. Attribute ks.grids might be initialized.

dmndarray or list of ndarrays

A density matrix or a list of density matrices

Returns:

Veff : (nkpts, nao, nao) or (*, nkpts, nao, nao) ndarray Veff = J + Vxc.

mix_density_fit(auxbasis=None, with_df=None, *args, **kwargs)
nuc_grad_method()

Hook to create object for analytical nuclear gradients.

pyscf.pbc.dft.krks.get_rho(mf, dm=None, grids=None, kpts=None)

Compute density in real space

pyscf.pbc.dft.krks.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)

Coulomb + XC functional

Note

This is a replica of pyscf.dft.rks.get_veff with kpts added. This function will change the ks object.

Args:
ksan instance of RKS

XC functional are controlled by ks.xc attribute. Attribute ks.grids might be initialized.

dmndarray or list of ndarrays

A density matrix or a list of density matrices

Returns:

Veff : (nkpts, nao, nao) or (*, nkpts, nao, nao) ndarray Veff = J + Vxc.

pyscf.pbc.dft.kroks module

Restricted open-shell Kohn-Sham for periodic systems with k-point sampling

class pyscf.pbc.dft.kroks.KROKS(cell, kpts=array([[0., 0., 0.]]), xc='LDA,VWN', exxdiv='ewald')

Bases: pyscf.pbc.dft.rks.KohnShamDFT, pyscf.pbc.scf.krohf.KROHF

RKS class adapted for PBCs with k-point sampling.

density_fit(auxbasis=None, with_df=None, *args, **kwargs)
dump_flags(verbose=None)
energy_elec(dm_kpts=None, h1e_kpts=None, vhf=None)

Following pyscf.scf.hf.energy_elec()

get_rho(dm=None, grids=None, kpts=None)

Compute density in real space

get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)

Coulomb + XC functional for UKS. See pyscf/pbc/dft/uks.py get_veff() fore more details.

mix_density_fit(auxbasis=None, with_df=None, *args, **kwargs)
pyscf.pbc.dft.kroks.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)

Coulomb + XC functional for UKS. See pyscf/pbc/dft/uks.py get_veff() fore more details.

pyscf.pbc.dft.kuks module

Non-relativistic Restricted Kohn-Sham for periodic systems with k-point sampling

See Also:
pyscf.pbc.dft.rks.pyNon-relativistic Restricted Kohn-Sham for periodic

systems at a single k-point

class pyscf.pbc.dft.kuks.KUKS(cell, kpts=array([[0., 0., 0.]]), xc='LDA,VWN', exxdiv='ewald')

Bases: pyscf.pbc.dft.rks.KohnShamDFT, pyscf.pbc.scf.kuhf.KUHF

RKS class adapted for PBCs with k-point sampling.

density_fit(auxbasis=None, with_df=None, *args, **kwargs)
dump_flags(verbose=None)
energy_elec(dm_kpts=None, h1e_kpts=None, vhf=None)

Following pyscf.scf.hf.energy_elec()

get_rho(dm=None, grids=None, kpts=None)

Compute density in real space

get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)

Coulomb + XC functional for UKS. See pyscf/pbc/dft/uks.py get_veff() fore more details.

mix_density_fit(auxbasis=None, with_df=None, *args, **kwargs)
nuc_grad_method()

Hook to create object for analytical nuclear gradients.

pyscf.pbc.dft.kuks.energy_elec(mf, dm_kpts=None, h1e_kpts=None, vhf=None)
pyscf.pbc.dft.kuks.get_rho(mf, dm=None, grids=None, kpts=None)

Compute density in real space

pyscf.pbc.dft.kuks.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)

Coulomb + XC functional for UKS. See pyscf/pbc/dft/uks.py get_veff() fore more details.

pyscf.pbc.dft.multigrid module

Multigrid to compute DFT integrals

class pyscf.pbc.dft.multigrid.MultiGridFFTDF(cell, kpts=array([[0., 0., 0.]]))

Bases: pyscf.pbc.df.fft.FFTDF

build()
get_jk(dm, hermi=1, kpts=None, kpts_band=None, with_j=True, with_k=True, exxdiv='ewald', **kwargs)
get_nuc(kpts=None)
get_pp(kpts=None)

Get the periodic pseudotential nuc-el AO matrix, with G=0 removed.

get_rho(dm, kpts=array([[0., 0., 0.]]))

Density in real space

reset(cell=None)
pyscf.pbc.dft.multigrid.cache_xc_kernel(mydf, xc_code, dm, spin=0, kpts=None)

Compute the 0th order density, Vxc and fxc. They can be used in TDDFT, DFT hessian module etc.

pyscf.pbc.dft.multigrid.eval_mat(cell, weights, shls_slice=None, comp=1, hermi=0, xctype='LDA', kpts=None, mesh=None, offset=None, submesh=None)
pyscf.pbc.dft.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)

Collocate the real density (opt. gradients) on the real-space grid.

pyscf.pbc.dft.multigrid.get_j_kpts(mydf, dm_kpts, hermi=1, kpts=array([[0., 0., 0.]]), kpts_band=None)

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.get_nuc(mydf, kpts=None)
pyscf.pbc.dft.multigrid.get_pp(mydf, kpts=None)

Get the periodic pseudotential nuc-el AO matrix, with G=0 removed.

pyscf.pbc.dft.multigrid.get_rho(mydf, dm, kpts=array([[0., 0., 0.]]))

Density in real space

pyscf.pbc.dft.multigrid.multi_grids_tasks(cell, fft_mesh=None, verbose=None)
pyscf.pbc.dft.multigrid.multi_grids_tasks_for_ke_cut(cell, fft_mesh=None, verbose=None)
pyscf.pbc.dft.multigrid.multi_grids_tasks_for_rcut(cell, fft_mesh=None, verbose=None)
pyscf.pbc.dft.multigrid.multigrid(mf)

Use MultiGridFFTDF to replace the default FFTDF integration method in the DFT object.

pyscf.pbc.dft.multigrid.nr_rks(mydf, xc_code, dm_kpts, hermi=1, kpts=None, kpts_band=None, with_j=False, return_j=False, verbose=None)

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.nr_rks_fxc(mydf, xc_code, dm0, dms, hermi=1, with_j=False, rho0=None, vxc=None, fxc=None, kpts=None, verbose=None)

multigrid version of function pbc.dft.numint.nr_rks_fxc

pyscf.pbc.dft.multigrid.nr_rks_fxc_st(mydf, xc_code, dm0, dms_alpha, hermi=1, singlet=True, with_j=False, rho0=None, vxc=None, fxc=None, kpts=None, verbose=None)

multigrid version of function pbc.dft.numint.nr_rks_fxc_st

pyscf.pbc.dft.multigrid.nr_uks(mydf, xc_code, dm_kpts, hermi=1, kpts=None, kpts_band=None, with_j=False, return_j=False, verbose=None)

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.nr_uks_fxc(mydf, xc_code, dm0, dms, hermi=1, with_j=False, rho0=None, vxc=None, fxc=None, kpts=None, verbose=None)

multigrid version of function pbc.dft.numint.nr_uks_fxc

pyscf.pbc.dft.numint module

class pyscf.pbc.dft.numint.KNumInt(kpts=array([[0., 0., 0.]]))

Bases: pyscf.dft.numint.NumInt

Generalization of pyscf’s NumInt class for k-point sampling and periodic images.

block_loop(cell, grids, nao=None, deriv=0, kpts=array([[0., 0., 0.]]), kpts_band=None, max_memory=2000, non0tab=None, blksize=None)

Define this macro to loop over grids by blocks.

cache_xc_kernel(cell, grids, xc_code, mo_coeff, mo_occ, spin=0, kpts=None, max_memory=2000)

Compute the 0th order density, Vxc and fxc. They can be used in TDDFT, DFT hessian module etc.

eval_ao(cell, coords, kpts=array([[0., 0., 0.]]), deriv=0, relativity=0, shls_slice=None, non0tab=None, out=None, verbose=None, **kwargs)

Evaluate AO function value on the given grids.

Args:

mol : an instance of Mole

coords2D array, shape (N,3)

The coordinates of the grids.

Kwargs:
derivint

AO derivative order. It affects the shape of the return array. If deriv=0, the returned AO values are stored in a (N,nao) array. Otherwise the AO values are stored in an array of shape (M,N,nao). Here N is the number of grids, nao is the number of AO functions, M is the size associated to the derivative deriv.

relativitybool

No effects.

shls_slice2-element list

(shl_start, shl_end). If given, only part of AOs (shl_start <= shell_id < shl_end) are evaluated. By default, all shells defined in mol will be evaluated.

non0tab2D bool array

mask array to indicate whether the AO values are zero. The mask array can be obtained by calling make_mask()

outndarray

If provided, results are written into this array.

verboseint or object of Logger

No effects.

Returns:

2D array of shape (N,nao) for AO values if deriv = 0. Or 3D array of shape (:,N,nao) for AO values and AO derivatives if deriv > 0. In the 3D array, the first (N,nao) elements are the AO values, followed by (3,N,nao) for x,y,z compoents; Then 2nd derivatives (6,N,nao) for xx, xy, xz, yy, yz, zz; Then 3rd derivatives (10,N,nao) for xxx, xxy, xxz, xyy, xyz, xzz, yyy, yyz, yzz, zzz; …

Examples:

>>> mol = gto.M(atom='O 0 0 0; H 0 0 1; H 0 1 0', basis='ccpvdz')
>>> coords = numpy.random.random((100,3))  # 100 random points
>>> ao_value = eval_ao(mol, coords)
>>> print(ao_value.shape)
(100, 24)
>>> ao_value = eval_ao(mol, coords, deriv=1, shls_slice=(1,4))
>>> print(ao_value.shape)
(4, 100, 7)
>>> ao_value = eval_ao(mol, coords, deriv=2, shls_slice=(1,4))
>>> print(ao_value.shape)
(10, 100, 7)
eval_mat(cell, ao_kpts, weight, rho, vxc, non0tab=None, xctype='LDA', spin=0, verbose=None)
eval_rho(cell, ao_kpts, dm_kpts, non0tab=None, xctype='LDA', hermi=0, verbose=None)

Collocate the real density (opt. gradients) on the real-space grid.

Args:

cell : Mole or Cell object ao_kpts : (nkpts, ngrids, nao) ndarray

AO values at each k-point

dm_kpts: (nkpts, nao, nao) ndarray

Density matrix at each k-point

Returns:

rhoR : (ngrids,) ndarray

eval_rho2(cell, ao_kpts, mo_coeff_kpts, mo_occ_kpts, non0tab=None, xctype='LDA', verbose=None)

Calculate the electron density for LDA functional, and the density derivatives for GGA functional. This function has the same functionality as eval_rho() except that the density are evaluated based on orbital coefficients and orbital occupancy. It is more efficient than eval_rho() in most scenario.

Args:

mol : an instance of Mole

ao2D array of shape (N,nao) for LDA, 3D array of shape (4,N,nao) for GGA

or (5,N,nao) for meta-GGA. N is the number of grids, nao is the number of AO functions. If xctype is GGA, ao[0] is AO value and ao[1:3] are the AO gradients. If xctype is meta-GGA, ao[4:10] are second derivatives of ao values.

dm2D array

Density matrix

Kwargs:
non0tab2D bool array

mask array to indicate whether the AO values are zero. The mask array can be obtained by calling make_mask()

xctypestr

LDA/GGA/mGGA. It affects the shape of the return density.

verboseint or object of Logger

No effects.

Returns:

1D array of size N to store electron density if xctype = LDA; 2D array of (4,N) to store density and “density derivatives” for x,y,z components if xctype = GGA; (6,N) array for meta-GGA, where last two rows are nabla^2 rho and tau = 1/2(nabla f)^2

get_rho(cell, dm, grids, kpts=array([[0., 0., 0.]]), max_memory=2000)

Density in real space

make_mask(cell, coords, relativity=0, shls_slice=None, verbose=None)

Mask to indicate whether a shell is zero on grid. The resultant mask array is an extension to the mask array used in molecular code (see also pyscf.dft.numint.make_mask function). For given shell ID and block ID, the value of the extended mask array means the number of images in Ls that does not vanish.

nr_rks(cell, grids, xc_code, dms, hermi=0, kpts=None, kpts_band=None, max_memory=2000, verbose=None, **kwargs)

Calculate RKS XC functional and potential matrix for given meshgrids and density matrix

Note: This is a replica of pyscf.dft.numint.nr_rks_vxc with kpts added. This implemented uses slow function in numint, which only calls eval_rho, eval_mat. Faster function uses eval_rho2 which is not yet implemented.

Args:

ni : an instance of NumInt or KNumInt

cell : instance of Mole or Cell

gridsan instance of Grids

grids.coords and grids.weights are needed for coordinates and weights of meshgrids.

xc_codestr

XC functional description. See parse_xc() of pyscf/dft/libxc.py for more details.

dms2D/3D array or a list of 2D/3D arrays

Density matrices (2D) / density matrices for k-points (3D)

Kwargs:
spinint

spin polarized if spin = 1

relativityint

No effects.

hermiint

No effects

max_memoryint or float

The maximum size of cache to use (in MB).

verboseint or object of Logger

No effects.

kpts(3,) ndarray or (nkpts,3) ndarray

Single or multiple k-points sampled for the DM. Default is gamma point.

kpts_band(3,) ndarray or (*,3) ndarray

A list of arbitrary “band” k-points at which to evaluate the XC matrix.

Returns:

nelec, excsum, vmat. nelec is the number of electrons generated by numerical integration. excsum is the XC functional value. vmat is the XC potential matrix in 2D array of shape (nao,nao) where nao is the number of AO functions.

nr_rks_fxc(cell, grids, xc_code, dm0, dms, relativity=0, hermi=0, rho0=None, vxc=None, fxc=None, kpts=None, max_memory=2000, verbose=None)

Contract RKS XC kernel matrix with given density matrices

Args:

ni : an instance of NumInt or KNumInt

cell : instance of Mole or Cell

gridsan instance of Grids

grids.coords and grids.weights are needed for coordinates and weights of meshgrids.

xc_codestr

XC functional description. See parse_xc() of pyscf/dft/libxc.py for more details.

dms2D/3D array or a list of 2D/3D arrays

Density matrices (2D) / density matrices for k-points (3D)

Kwargs:
hermiint

Input density matrices symmetric or not

max_memoryint or float

The maximum size of cache to use (in MB).

rho0float array

Zero-order density (and density derivative for GGA). Giving kwargs rho0, vxc and fxc to improve better performance.

vxcfloat array

First order XC derivatives

fxcfloat array

Second order XC derivatives

Examples:

nr_uks(cell, grids, xc_code, dms, hermi=0, kpts=None, kpts_band=None, max_memory=2000, verbose=None, **kwargs)

Calculate UKS XC functional and potential matrix for given meshgrids and density matrix

Note: This is a replica of pyscf.dft.numint.nr_rks_vxc with kpts added. This implemented uses slow function in numint, which only calls eval_rho, eval_mat. Faster function uses eval_rho2 which is not yet implemented.

Args:

ni : an instance of NumInt or KNumInt

cell : instance of Mole or Cell

gridsan instance of Grids

grids.coords and grids.weights are needed for coordinates and weights of meshgrids.

xc_codestr

XC functional description. See parse_xc() of pyscf/dft/libxc.py for more details.

dms :

Density matrices

Kwargs:
spinint

spin polarized if spin = 1

relativityint

No effects.

hermiint

Input density matrices symmetric or not

max_memoryint or float

The maximum size of cache to use (in MB).

verboseint or object of Logger

No effects.

kpts(3,) ndarray or (nkpts,3) ndarray

Single or multiple k-points sampled for the DM. Default is gamma point. kpts_band : (3,) ndarray or (*,3) ndarray A list of arbitrary “band” k-points at which to evaluate the XC matrix.

Returns:

nelec, excsum, vmat. nelec is the number of electrons generated by numerical integration. excsum is the XC functional value. vmat is the XC potential matrix in 2D array of shape (nao,nao) where nao is the number of AO functions.

nr_uks_fxc(cell, grids, xc_code, dm0, dms, relativity=0, hermi=0, rho0=None, vxc=None, fxc=None, kpts=None, max_memory=2000, verbose=None)

Contract UKS XC kernel matrix with given density matrices

Args:

ni : an instance of NumInt or KNumInt

cell : instance of Mole or Cell

gridsan instance of Grids

grids.coords and grids.weights are needed for coordinates and weights of meshgrids.

xc_codestr

XC functional description. See parse_xc() of pyscf/dft/libxc.py for more details.

dms2D array a list of 2D arrays

Density matrix or multiple density matrices

Kwargs:
hermiint

Input density matrices symmetric or not

max_memoryint or float

The maximum size of cache to use (in MB).

rho0float array

Zero-order density (and density derivative for GGA). Giving kwargs rho0, vxc and fxc to improve better performance.

vxcfloat array

First order XC derivatives

fxcfloat array

Second order XC derivatives

Returns:

nelec, excsum, vmat. nelec is the number of electrons generated by numerical integration. excsum is the XC functional value. vmat is the XC potential matrix in 2D array of shape (nao,nao) where nao is the number of AO functions.

Examples:

nr_vxc(cell, grids, xc_code, dms, spin=0, relativity=0, hermi=0, kpts=None, kpts_band=None, max_memory=2000, verbose=None)

Evaluate RKS/UKS XC functional and potential matrix. See nr_rks() and nr_uks() for more details.

class pyscf.pbc.dft.numint.NumInt

Bases: pyscf.dft.numint.NumInt

Generalization of pyscf’s NumInt class for a single k-point shift and periodic images.

block_loop(cell, grids, nao=None, deriv=0, kpt=array([0., 0., 0.]), kpts_band=None, max_memory=2000, non0tab=None, blksize=None)

Define this macro to loop over grids by blocks.

cache_xc_kernel(cell, grids, xc_code, mo_coeff, mo_occ, spin=0, kpts=None, max_memory=2000)

Compute the 0th order density, Vxc and fxc. They can be used in TDDFT, DFT hessian module etc.

eval_ao(cell, coords, kpt=array([0., 0., 0.]), deriv=0, relativity=0, shls_slice=None, non0tab=None, out=None, verbose=None)

Evaluate AO function value on the given grids.

Args:

mol : an instance of Mole

coords2D array, shape (N,3)

The coordinates of the grids.

Kwargs:
derivint

AO derivative order. It affects the shape of the return array. If deriv=0, the returned AO values are stored in a (N,nao) array. Otherwise the AO values are stored in an array of shape (M,N,nao). Here N is the number of grids, nao is the number of AO functions, M is the size associated to the derivative deriv.

relativitybool

No effects.

shls_slice2-element list

(shl_start, shl_end). If given, only part of AOs (shl_start <= shell_id < shl_end) are evaluated. By default, all shells defined in mol will be evaluated.

non0tab2D bool array

mask array to indicate whether the AO values are zero. The mask array can be obtained by calling make_mask()

outndarray

If provided, results are written into this array.

verboseint or object of Logger

No effects.

Returns:

2D array of shape (N,nao) for AO values if deriv = 0. Or 3D array of shape (:,N,nao) for AO values and AO derivatives if deriv > 0. In the 3D array, the first (N,nao) elements are the AO values, followed by (3,N,nao) for x,y,z compoents; Then 2nd derivatives (6,N,nao) for xx, xy, xz, yy, yz, zz; Then 3rd derivatives (10,N,nao) for xxx, xxy, xxz, xyy, xyz, xzz, yyy, yyz, yzz, zzz; …

Examples:

>>> mol = gto.M(atom='O 0 0 0; H 0 0 1; H 0 1 0', basis='ccpvdz')
>>> coords = numpy.random.random((100,3))  # 100 random points
>>> ao_value = eval_ao(mol, coords)
>>> print(ao_value.shape)
(100, 24)
>>> ao_value = eval_ao(mol, coords, deriv=1, shls_slice=(1,4))
>>> print(ao_value.shape)
(4, 100, 7)
>>> ao_value = eval_ao(mol, coords, deriv=2, shls_slice=(1,4))
>>> print(ao_value.shape)
(10, 100, 7)
eval_mat(cell, ao, weight, rho, vxc, non0tab=None, xctype='LDA', spin=0, verbose=None)
eval_rho(cell, ao, dm, non0tab=None, xctype='LDA', hermi=0, verbose=None)

Collocate the real density (opt. gradients) on the real-space grid.

Args:

cell : instance of Mole or Cell

ao([4,] nx*ny*nz, nao=cell.nao_nr()) ndarray

The value of the AO crystal orbitals on the real-space grid by default. If xctype=’GGA’, also contains the value of the gradient in the x, y, and z directions.

Returns:
rho([4,] nx*ny*nz) ndarray

The value of the density on the real-space grid. If xctype=’GGA’, also contains the value of the gradient in the x, y, and z directions.

See Also:

pyscf.dft.numint.eval_rho

eval_rho2(cell, ao, mo_coeff, mo_occ, non0tab=None, xctype='LDA', verbose=None)

Calculate the electron density for LDA functional, and the density derivatives for GGA functional. This function has the same functionality as eval_rho() except that the density are evaluated based on orbital coefficients and orbital occupancy. It is more efficient than eval_rho() in most scenario.

Args:

mol : an instance of Mole

ao2D array of shape (N,nao) for LDA, 3D array of shape (4,N,nao) for GGA

or (5,N,nao) for meta-GGA. N is the number of grids, nao is the number of AO functions. If xctype is GGA, ao[0] is AO value and ao[1:3] are the AO gradients. If xctype is meta-GGA, ao[4:10] are second derivatives of ao values.

dm2D array

Density matrix

Kwargs:
non0tab2D bool array

mask array to indicate whether the AO values are zero. The mask array can be obtained by calling make_mask()

xctypestr

LDA/GGA/mGGA. It affects the shape of the return density.

verboseint or object of Logger

No effects.

Returns:

1D array of size N to store electron density if xctype = LDA; 2D array of (4,N) to store density and “density derivatives” for x,y,z components if xctype = GGA; (6,N) array for meta-GGA, where last two rows are nabla^2 rho and tau = 1/2(nabla f)^2

get_rho(cell, dm, grids, kpts=array([[0., 0., 0.]]), max_memory=2000)

Density in real space

make_mask(cell, coords, relativity=0, shls_slice=None, verbose=None)

Mask to indicate whether a shell is zero on grid. The resultant mask array is an extension to the mask array used in molecular code (see also pyscf.dft.numint.make_mask function). For given shell ID and block ID, the value of the extended mask array means the number of images in Ls that does not vanish.

nr_rks(cell, grids, xc_code, dms, hermi=0, kpt=array([0., 0., 0.]), kpts_band=None, max_memory=2000, verbose=None)

Calculate RKS XC functional and potential matrix for given meshgrids and density matrix

Note: This is a replica of pyscf.dft.numint.nr_rks_vxc with kpts added. This implemented uses slow function in numint, which only calls eval_rho, eval_mat. Faster function uses eval_rho2 which is not yet implemented.

Args:

ni : an instance of NumInt or KNumInt

cell : instance of Mole or Cell

gridsan instance of Grids

grids.coords and grids.weights are needed for coordinates and weights of meshgrids.

xc_codestr

XC functional description. See parse_xc() of pyscf/dft/libxc.py for more details.

dms2D/3D array or a list of 2D/3D arrays

Density matrices (2D) / density matrices for k-points (3D)

Kwargs:
spinint

spin polarized if spin = 1

relativityint

No effects.

hermiint

No effects

max_memoryint or float

The maximum size of cache to use (in MB).

verboseint or object of Logger

No effects.

kpts(3,) ndarray or (nkpts,3) ndarray

Single or multiple k-points sampled for the DM. Default is gamma point.

kpts_band(3,) ndarray or (*,3) ndarray

A list of arbitrary “band” k-points at which to evaluate the XC matrix.

Returns:

nelec, excsum, vmat. nelec is the number of electrons generated by numerical integration. excsum is the XC functional value. vmat is the XC potential matrix in 2D array of shape (nao,nao) where nao is the number of AO functions.

nr_rks_fxc(cell, grids, xc_code, dm0, dms, relativity=0, hermi=0, rho0=None, vxc=None, fxc=None, kpts=None, max_memory=2000, verbose=None)

Contract RKS XC kernel matrix with given density matrices

Args:

ni : an instance of NumInt or KNumInt

cell : instance of Mole or Cell

gridsan instance of Grids

grids.coords and grids.weights are needed for coordinates and weights of meshgrids.

xc_codestr

XC functional description. See parse_xc() of pyscf/dft/libxc.py for more details.

dms2D/3D array or a list of 2D/3D arrays

Density matrices (2D) / density matrices for k-points (3D)

Kwargs:
hermiint

Input density matrices symmetric or not

max_memoryint or float

The maximum size of cache to use (in MB).

rho0float array

Zero-order density (and density derivative for GGA). Giving kwargs rho0, vxc and fxc to improve better performance.

vxcfloat array

First order XC derivatives

fxcfloat array

Second order XC derivatives

Examples:

nr_uks(cell, grids, xc_code, dms, hermi=0, kpt=array([0., 0., 0.]), kpts_band=None, max_memory=2000, verbose=None)

Calculate UKS XC functional and potential matrix for given meshgrids and density matrix

Note: This is a replica of pyscf.dft.numint.nr_rks_vxc with kpts added. This implemented uses slow function in numint, which only calls eval_rho, eval_mat. Faster function uses eval_rho2 which is not yet implemented.

Args:

ni : an instance of NumInt or KNumInt

cell : instance of Mole or Cell

gridsan instance of Grids

grids.coords and grids.weights are needed for coordinates and weights of meshgrids.

xc_codestr

XC functional description. See parse_xc() of pyscf/dft/libxc.py for more details.

dms :

Density matrices

Kwargs:
spinint

spin polarized if spin = 1

relativityint

No effects.

hermiint

Input density matrices symmetric or not

max_memoryint or float

The maximum size of cache to use (in MB).

verboseint or object of Logger

No effects.

kpts(3,) ndarray or (nkpts,3) ndarray

Single or multiple k-points sampled for the DM. Default is gamma point. kpts_band : (3,) ndarray or (*,3) ndarray A list of arbitrary “band” k-points at which to evaluate the XC matrix.

Returns:

nelec, excsum, vmat. nelec is the number of electrons generated by numerical integration. excsum is the XC functional value. vmat is the XC potential matrix in 2D array of shape (nao,nao) where nao is the number of AO functions.

nr_uks_fxc(cell, grids, xc_code, dm0, dms, relativity=0, hermi=0, rho0=None, vxc=None, fxc=None, kpts=None, max_memory=2000, verbose=None)

Contract UKS XC kernel matrix with given density matrices

Args:

ni : an instance of NumInt or KNumInt

cell : instance of Mole or Cell

gridsan instance of Grids

grids.coords and grids.weights are needed for coordinates and weights of meshgrids.

xc_codestr

XC functional description. See parse_xc() of pyscf/dft/libxc.py for more details.

dms2D array a list of 2D arrays

Density matrix or multiple density matrices

Kwargs:
hermiint

Input density matrices symmetric or not

max_memoryint or float

The maximum size of cache to use (in MB).

rho0float array

Zero-order density (and density derivative for GGA). Giving kwargs rho0, vxc and fxc to improve better performance.

vxcfloat array

First order XC derivatives

fxcfloat array

Second order XC derivatives

Returns:

nelec, excsum, vmat. nelec is the number of electrons generated by numerical integration. excsum is the XC functional value. vmat is the XC potential matrix in 2D array of shape (nao,nao) where nao is the number of AO functions.

Examples:

nr_vxc(cell, grids, xc_code, dms, spin=0, relativity=0, hermi=0, kpt=None, kpts_band=None, max_memory=2000, verbose=None)

Evaluate RKS/UKS XC functional and potential matrix. See nr_rks() and nr_uks() for more details.

pyscf.pbc.dft.numint.cache_xc_kernel(ni, cell, grids, xc_code, mo_coeff, mo_occ, spin=0, kpts=None, max_memory=2000)

Compute the 0th order density, Vxc and fxc. They can be used in TDDFT, DFT hessian module etc.

pyscf.pbc.dft.numint.eval_ao(cell, coords, kpt=array([0., 0., 0.]), deriv=0, relativity=0, shls_slice=None, non0tab=None, out=None, verbose=None)

Collocate AO crystal orbitals (opt. gradients) on the real-space grid.

Args:

cell : instance of Cell

coords(nx*ny*nz, 3) ndarray

The real-space grid point coordinates.

Kwargs:
kpt(3,) ndarray

The k-point corresponding to the crystal AO.

derivint

AO derivative order. It affects the shape of the return array. If deriv=0, the returned AO values are stored in a (N,nao) array. Otherwise the AO values are stored in an array of shape (M,N,nao). Here N is the number of grids, nao is the number of AO functions, M is the size associated to the derivative deriv.

Returns:
aoR([4,] nx*ny*nz, nao=cell.nao_nr()) ndarray

The value of the AO crystal orbitals on the real-space grid by default. If deriv=1, also contains the value of the orbitals gradient in the x, y, and z directions. It can be either complex or float array, depending on the kpt argument. If kpt is not given (gamma point), aoR is a float array.

See Also:

pyscf.dft.numint.eval_ao

pyscf.pbc.dft.numint.eval_ao_kpts(cell, coords, kpts=None, deriv=0, relativity=0, shls_slice=None, non0tab=None, out=None, verbose=None, **kwargs)
Returns:
ao_kpts: (nkpts, [comp], ngrids, nao) ndarray

AO values at each k-point

pyscf.pbc.dft.numint.eval_rho(cell, ao, dm, non0tab=None, xctype='LDA', hermi=0, verbose=None)

Collocate the real density (opt. gradients) on the real-space grid.

Args:

cell : instance of Mole or Cell

ao([4,] nx*ny*nz, nao=cell.nao_nr()) ndarray

The value of the AO crystal orbitals on the real-space grid by default. If xctype=’GGA’, also contains the value of the gradient in the x, y, and z directions.

Returns:
rho([4,] nx*ny*nz) ndarray

The value of the density on the real-space grid. If xctype=’GGA’, also contains the value of the gradient in the x, y, and z directions.

See Also:

pyscf.dft.numint.eval_rho

pyscf.pbc.dft.numint.eval_rho2(cell, ao, mo_coeff, mo_occ, non0tab=None, xctype='LDA', verbose=None)

Refer to pyscf.dft.numint.eval_rho2 for full documentation.

pyscf.pbc.dft.numint.get_rho(ni, cell, dm, grids, kpts=array([[0., 0., 0.]]), max_memory=2000)

Density in real space

pyscf.pbc.dft.numint.nr_rks(ni, cell, grids, xc_code, dms, spin=0, relativity=0, hermi=0, kpts=None, kpts_band=None, max_memory=2000, verbose=None)

Calculate RKS XC functional and potential matrix for given meshgrids and density matrix

Note: This is a replica of pyscf.dft.numint.nr_rks_vxc with kpts added. This implemented uses slow function in numint, which only calls eval_rho, eval_mat. Faster function uses eval_rho2 which is not yet implemented.

Args:

ni : an instance of NumInt or KNumInt

cell : instance of Mole or Cell

gridsan instance of Grids

grids.coords and grids.weights are needed for coordinates and weights of meshgrids.

xc_codestr

XC functional description. See parse_xc() of pyscf/dft/libxc.py for more details.

dms2D/3D array or a list of 2D/3D arrays

Density matrices (2D) / density matrices for k-points (3D)

Kwargs:
spinint

spin polarized if spin = 1

relativityint

No effects.

hermiint

No effects

max_memoryint or float

The maximum size of cache to use (in MB).

verboseint or object of Logger

No effects.

kpts(3,) ndarray or (nkpts,3) ndarray

Single or multiple k-points sampled for the DM. Default is gamma point.

kpts_band(3,) ndarray or (*,3) ndarray

A list of arbitrary “band” k-points at which to evaluate the XC matrix.

Returns:

nelec, excsum, vmat. nelec is the number of electrons generated by numerical integration. excsum is the XC functional value. vmat is the XC potential matrix in 2D array of shape (nao,nao) where nao is the number of AO functions.

pyscf.pbc.dft.numint.nr_rks_fxc(ni, cell, grids, xc_code, dm0, dms, relativity=0, hermi=0, rho0=None, vxc=None, fxc=None, kpts=None, max_memory=2000, verbose=None)

Contract RKS XC kernel matrix with given density matrices

Args:

ni : an instance of NumInt or KNumInt

cell : instance of Mole or Cell

gridsan instance of Grids

grids.coords and grids.weights are needed for coordinates and weights of meshgrids.

xc_codestr

XC functional description. See parse_xc() of pyscf/dft/libxc.py for more details.

dms2D/3D array or a list of 2D/3D arrays

Density matrices (2D) / density matrices for k-points (3D)

Kwargs:
hermiint

Input density matrices symmetric or not

max_memoryint or float

The maximum size of cache to use (in MB).

rho0float array

Zero-order density (and density derivative for GGA). Giving kwargs rho0, vxc and fxc to improve better performance.

vxcfloat array

First order XC derivatives

fxcfloat array

Second order XC derivatives

Examples:

pyscf.pbc.dft.numint.nr_rks_fxc_st(ni, cell, grids, xc_code, dm0, dms_alpha, relativity=0, singlet=True, rho0=None, vxc=None, fxc=None, kpts=None, max_memory=2000, verbose=None)

Associated to singlet or triplet Hessian Note the difference to nr_rks_fxc, dms_alpha is the response density matrices of alpha spin, alpha+/-beta DM is applied due to singlet/triplet coupling

Ref. CPL, 256, 454

pyscf.pbc.dft.numint.nr_rks_vxc(ni, cell, grids, xc_code, dms, spin=0, relativity=0, hermi=0, kpts=None, kpts_band=None, max_memory=2000, verbose=None)

Calculate RKS XC functional and potential matrix for given meshgrids and density matrix

Note: This is a replica of pyscf.dft.numint.nr_rks_vxc with kpts added. This implemented uses slow function in numint, which only calls eval_rho, eval_mat. Faster function uses eval_rho2 which is not yet implemented.

Args:

ni : an instance of NumInt or KNumInt

cell : instance of Mole or Cell

gridsan instance of Grids

grids.coords and grids.weights are needed for coordinates and weights of meshgrids.

xc_codestr

XC functional description. See parse_xc() of pyscf/dft/libxc.py for more details.

dms2D/3D array or a list of 2D/3D arrays

Density matrices (2D) / density matrices for k-points (3D)

Kwargs:
spinint

spin polarized if spin = 1

relativityint

No effects.

hermiint

No effects

max_memoryint or float

The maximum size of cache to use (in MB).

verboseint or object of Logger

No effects.

kpts(3,) ndarray or (nkpts,3) ndarray

Single or multiple k-points sampled for the DM. Default is gamma point.

kpts_band(3,) ndarray or (*,3) ndarray

A list of arbitrary “band” k-points at which to evaluate the XC matrix.

Returns:

nelec, excsum, vmat. nelec is the number of electrons generated by numerical integration. excsum is the XC functional value. vmat is the XC potential matrix in 2D array of shape (nao,nao) where nao is the number of AO functions.

pyscf.pbc.dft.numint.nr_uks(ni, cell, grids, xc_code, dms, spin=1, relativity=0, hermi=0, kpts=None, kpts_band=None, max_memory=2000, verbose=None)

Calculate UKS XC functional and potential matrix for given meshgrids and density matrix

Note: This is a replica of pyscf.dft.numint.nr_rks_vxc with kpts added. This implemented uses slow function in numint, which only calls eval_rho, eval_mat. Faster function uses eval_rho2 which is not yet implemented.

Args:

ni : an instance of NumInt or KNumInt

cell : instance of Mole or Cell

gridsan instance of Grids

grids.coords and grids.weights are needed for coordinates and weights of meshgrids.

xc_codestr

XC functional description. See parse_xc() of pyscf/dft/libxc.py for more details.

dms :

Density matrices

Kwargs:
spinint

spin polarized if spin = 1

relativityint

No effects.

hermiint

Input density matrices symmetric or not

max_memoryint or float

The maximum size of cache to use (in MB).

verboseint or object of Logger

No effects.

kpts(3,) ndarray or (nkpts,3) ndarray

Single or multiple k-points sampled for the DM. Default is gamma point. kpts_band : (3,) ndarray or (*,3) ndarray A list of arbitrary “band” k-points at which to evaluate the XC matrix.

Returns:

nelec, excsum, vmat. nelec is the number of electrons generated by numerical integration. excsum is the XC functional value. vmat is the XC potential matrix in 2D array of shape (nao,nao) where nao is the number of AO functions.

pyscf.pbc.dft.numint.nr_uks_fxc(ni, cell, grids, xc_code, dm0, dms, relativity=0, hermi=0, rho0=None, vxc=None, fxc=None, kpts=None, max_memory=2000, verbose=None)

Contract UKS XC kernel matrix with given density matrices

Args:

ni : an instance of NumInt or KNumInt

cell : instance of Mole or Cell

gridsan instance of Grids

grids.coords and grids.weights are needed for coordinates and weights of meshgrids.

xc_codestr

XC functional description. See parse_xc() of pyscf/dft/libxc.py for more details.

dms2D array a list of 2D arrays

Density matrix or multiple density matrices

Kwargs:
hermiint

Input density matrices symmetric or not

max_memoryint or float

The maximum size of cache to use (in MB).

rho0float array

Zero-order density (and density derivative for GGA). Giving kwargs rho0, vxc and fxc to improve better performance.

vxcfloat array

First order XC derivatives

fxcfloat array

Second order XC derivatives

Returns:

nelec, excsum, vmat. nelec is the number of electrons generated by numerical integration. excsum is the XC functional value. vmat is the XC potential matrix in 2D array of shape (nao,nao) where nao is the number of AO functions.

Examples:

pyscf.pbc.dft.numint.nr_uks_vxc(ni, cell, grids, xc_code, dms, spin=1, relativity=0, hermi=0, kpts=None, kpts_band=None, max_memory=2000, verbose=None)

Calculate UKS XC functional and potential matrix for given meshgrids and density matrix

Note: This is a replica of pyscf.dft.numint.nr_rks_vxc with kpts added. This implemented uses slow function in numint, which only calls eval_rho, eval_mat. Faster function uses eval_rho2 which is not yet implemented.

Args:

ni : an instance of NumInt or KNumInt

cell : instance of Mole or Cell

gridsan instance of Grids

grids.coords and grids.weights are needed for coordinates and weights of meshgrids.

xc_codestr

XC functional description. See parse_xc() of pyscf/dft/libxc.py for more details.

dms :

Density matrices

Kwargs:
spinint

spin polarized if spin = 1

relativityint

No effects.

hermiint

Input density matrices symmetric or not

max_memoryint or float

The maximum size of cache to use (in MB).

verboseint or object of Logger

No effects.

kpts(3,) ndarray or (nkpts,3) ndarray

Single or multiple k-points sampled for the DM. Default is gamma point. kpts_band : (3,) ndarray or (*,3) ndarray A list of arbitrary “band” k-points at which to evaluate the XC matrix.

Returns:

nelec, excsum, vmat. nelec is the number of electrons generated by numerical integration. excsum is the XC functional value. vmat is the XC potential matrix in 2D array of shape (nao,nao) where nao is the number of AO functions.

pyscf.pbc.dft.rks module

Non-relativistic Restricted Kohn-Sham for periodic systems at a single k-point

See Also:
pyscf.pbc.dft.krks.pyNon-relativistic Restricted Kohn-Sham for periodic

systems with k-point sampling

class pyscf.pbc.dft.rks.KohnShamDFT(xc='LDA,VWN')

Bases: pyscf.dft.rks.KohnShamDFT

PBC-KS

dump_flags(verbose=None)
reset(mol=None)
class pyscf.pbc.dft.rks.RKS(cell, kpt=array([0., 0., 0.]), xc='LDA,VWN', exxdiv='ewald')

Bases: pyscf.pbc.dft.rks.KohnShamDFT, pyscf.pbc.scf.hf.RHF

RKS class adapted for PBCs.

This is a literal duplication of the molecular RKS class with some mol variables replaced by cell.

density_fit(auxbasis=None, with_df=None, *args, **kwargs)
dump_flags(verbose=None)
energy_elec(dm=None, h1e=None, vhf=None)

Electronic part of RKS energy.

Note this function has side effects which cause mf.scf_summary updated.

Args:

ks : an instance of DFT class

dm2D ndarray

one-partical density matrix

h1e2D ndarray

Core hamiltonian

Returns:

RKS electronic energy and the 2-electron contribution

get_rho(dm=None, grids=None, kpt=None)

Compute density in real space

get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)

Coulomb + XC functional

Note

This function will change the ks object.

Args:
ksan instance of RKS

XC functional are controlled by ks.xc attribute. Attribute ks.grids might be initialized.

dmndarray or list of ndarrays

A density matrix or a list of density matrices

Returns:

matrix Veff = J + Vxc. Veff can be a list matrices, if the input dm is a list of density matrices.

mix_density_fit(auxbasis=None, with_df=None, *args, **kwargs)
pyscf.pbc.dft.rks.get_rho(mf, dm=None, grids=None, kpt=None)

Compute density in real space

pyscf.pbc.dft.rks.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)

Coulomb + XC functional

Note

This function will change the ks object.

Args:
ksan instance of RKS

XC functional are controlled by ks.xc attribute. Attribute ks.grids might be initialized.

dmndarray or list of ndarrays

A density matrix or a list of density matrices

Returns:

matrix Veff = J + Vxc. Veff can be a list matrices, if the input dm is a list of density matrices.

pyscf.pbc.dft.rks.prune_small_rho_grids_(ks, cell, dm, grids, kpts)

pyscf.pbc.dft.roks module

Restricted open-shell Kohn-Sham for periodic systems at a single k-point

class pyscf.pbc.dft.roks.ROKS(cell, kpt=array([0., 0., 0.]), xc='LDA,VWN', exxdiv='ewald')

Bases: pyscf.pbc.dft.rks.KohnShamDFT, pyscf.pbc.scf.rohf.ROHF

UKS class adapted for PBCs.

This is a literal duplication of the molecular UKS class with some mol variables replaced by cell.

density_fit(auxbasis=None, with_df=None, *args, **kwargs)
dump_flags(verbose=None)
energy_elec(dm=None, h1e=None, vhf=None)

Electronic part of Hartree-Fock energy, for given core hamiltonian and HF potential

… math:

E = \sum_{ij}h_{ij} \gamma_{ji}
  + \frac{1}{2}\sum_{ijkl} \gamma_{ji}\gamma_{lk} \langle ik||jl\rangle

Note this function has side effects which cause mf.scf_summary updated.

Args:

mf : an instance of SCF class

Kwargs:
dm2D ndarray

one-partical density matrix

h1e2D ndarray

Core hamiltonian

vhf2D ndarray

HF potential

Returns:

Hartree-Fock electronic energy and the Coulomb energy

Examples:

>>> from pyscf import gto, scf
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1')
>>> mf = scf.RHF(mol)
>>> mf.scf()
>>> dm = mf.make_rdm1()
>>> scf.hf.energy_elec(mf, dm)
(-1.5176090667746334, 0.60917167853723675)
>>> mf.energy_elec(dm)
(-1.5176090667746334, 0.60917167853723675)
get_rho(dm=None, grids=None, kpt=None)

Compute density in real space

get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)

Coulomb + XC functional for UKS. See pyscf/pbc/dft/uks.py get_veff() fore more details.

mix_density_fit(auxbasis=None, with_df=None, *args, **kwargs)
pyscf.pbc.dft.roks.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)

Coulomb + XC functional for UKS. See pyscf/pbc/dft/uks.py get_veff() fore more details.

pyscf.pbc.dft.uks module

Non-relativistic Restricted Kohn-Sham for periodic systems at a single k-point

See Also:
pyscf.pbc.dft.krks.pyNon-relativistic Restricted Kohn-Sham for periodic

systems with k-point sampling

class pyscf.pbc.dft.uks.UKS(cell, kpt=array([0., 0., 0.]), xc='LDA,VWN', exxdiv='ewald')

Bases: pyscf.pbc.dft.rks.KohnShamDFT, pyscf.pbc.scf.uhf.UHF

UKS class adapted for PBCs.

This is a literal duplication of the molecular UKS class with some mol variables replaced by cell.

density_fit(auxbasis=None, with_df=None, *args, **kwargs)
dump_flags(verbose=None)
energy_elec(dm=None, h1e=None, vhf=None)

Electronic energy of Unrestricted Hartree-Fock

Note this function has side effects which cause mf.scf_summary updated.

Returns:

Hartree-Fock electronic energy and the 2-electron part contribution

get_rho(dm=None, grids=None, kpt=None)

Compute density in real space

get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)

Coulomb + XC functional for UKS. See pyscf/pbc/dft/uks.py get_veff() fore more details.

mix_density_fit(auxbasis=None, with_df=None, *args, **kwargs)
pyscf.pbc.dft.uks.get_rho(mf, dm=None, grids=None, kpt=None)

Compute density in real space

pyscf.pbc.dft.uks.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)

Coulomb + XC functional for UKS. See pyscf/pbc/dft/uks.py get_veff() fore more details.

Module contents

pyscf.pbc.dft.KKS(cell, *args, **kwargs)

A wrap function to create DFT object with k-point sampling (KRKS or KUKS).

RKS class adapted for PBCs with k-point sampling.

pyscf.pbc.dft.KS(cell, *args, **kwargs)

A wrap function to create DFT object (RKS or UKS) for PBC systems.

RKS class adapted for PBCs.

This is a literal duplication of the molecular RKS class with some mol variables replaced by cell.

pyscf.pbc.dft.RKS(cell, *args, **kwargs)

RKS class adapted for PBCs.

This is a literal duplication of the molecular RKS class with some mol variables replaced by cell.