pyscf.pbc.dft package#
Subpackages#
- pyscf.pbc.dft.multigrid package
- Submodules
- pyscf.pbc.dft.multigrid.multigrid module
- pyscf.pbc.dft.multigrid.multigrid_pair module
- pyscf.pbc.dft.multigrid.pp module
- pyscf.pbc.dft.multigrid.utils module
- Module contents
Submodules#
pyscf.pbc.dft.cdft module#
- pyscf.pbc.dft.cdft.cdft(mf, cell, offset, orbital, basis=None)[source]#
- 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.gen_grid module#
- pyscf.pbc.dft.gen_grid.AtomicGrids#
alias of
BeckeGrids
- class pyscf.pbc.dft.gen_grid.BeckeGrids(cell)[source]#
Bases:
Grids
Atomic grids for all-electron calculation.
- make_mask(cell=None, coords=None, relativity=0, shls_slice=None, verbose=None)[source]#
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.
- 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.
- class pyscf.pbc.dft.gen_grid.UniformGrids(cell)[source]#
Bases:
StreamObject
Uniform Grid class.
- property coords#
- cutoff = 1e-15#
- kernel(cell=None, with_non0tab=False)[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.).
- make_mask(cell=None, coords=None, relativity=0, shls_slice=None, verbose=None)[source]#
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.
- 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.
- 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>)[source]#
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, cutoff=None, verbose=None)[source]#
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.gks module#
Generalized collinear Kohn-Sham in the spin-orbital basis for periodic systems at a single k-point
- See Also:
- pyscf.pbc.dft.kgks.pyGeneral spin-orbital Kohn-Sham for periodic
systems with k-point sampling
- class pyscf.pbc.dft.gks.GKS(cell, kpt=array([0., 0., 0.]), xc='LDA,VWN', exxdiv='ewald')[source]#
Bases:
KohnShamDFT
,GHF
GKS class adapted for PBCs at a single k-point.
This is a literal duplication of the molecular GKS class with some mol variables replaced by cell.
- property collinear#
- 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-particle density matrix
- h1e2D ndarray
Core hamiltonian
- Returns:
RKS electronic energy and the 2-electron contribution
- get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)#
Coulomb + XC functional for GKS.
- property spin_samples#
- 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.
- x2c()#
Adds spin-orbit coupling effects to H0 through the x2c1e approximation
pyscf.pbc.dft.kgks module#
Generalized collinear Kohn-Sham in the spin-orbital basis for periodic systems with k-point sampling
- See Also:
- pyscf.pbc.dft.rks.pyGeneral spin-orbital Kohn-Sham for periodic
systems at a single k-point
- class pyscf.pbc.dft.kgks.KGKS(cell, kpts=array([[0., 0., 0.]]), xc='LDA,VWN', exxdiv='ewald')[source]#
Bases:
KohnShamDFT
,KGHF
GKS class adapted for PBCs with k-point sampling.
- property collinear#
- 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 KGKS
- Args:
- ksan instance of
GKS
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
- ksan instance of
- Returns:
Veff : (nkpts, 2*nao, 2*nao) or (*, nkpts, 2*nao, 2*nao) ndarray Veff = J + Vxc.
- property spin_samples#
- 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.
- x2c()#
Adds spin-orbit coupling effects to H0 through the x2c1e approximation
- pyscf.pbc.dft.kgks.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)[source]#
Coulomb + XC functional for KGKS
- Args:
- ksan instance of
GKS
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
- ksan instance of
- Returns:
Veff : (nkpts, 2*nao, 2*nao) or (*, nkpts, 2*nao, 2*nao) ndarray Veff = J + Vxc.
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')[source]#
Bases:
KohnShamDFT
,KRHF
RKS class adapted for PBCs with k-point sampling.
- 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
- ksan instance of
- Returns:
Veff : (nkpts, nao, nao) or (*, nkpts, nao, nao) ndarray Veff = J + Vxc.
- 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.pbc.dft.krks.get_rho(mf, dm=None, grids=None, kpts=None)[source]#
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)[source]#
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
- ksan instance of
- Returns:
Veff : (nkpts, nao, nao) or (*, nkpts, nao, nao) ndarray Veff = J + Vxc.
pyscf.pbc.dft.krks_ksymm module#
- pyscf.pbc.dft.krks_ksymm.KRKS#
alias of
KsymAdaptedKRKS
- class pyscf.pbc.dft.krks_ksymm.KsymAdaptedKRKS(cell, kpts=<pyscf.pbc.lib.kpts.KPoints object>, xc='LDA, VWN', exxdiv='ewald', **kwargs)[source]#
Bases:
KRKS
,KsymAdaptedKRHF
- dump_chk(envs)#
Serialize the SCF object and save it to the specified chkfile.
- Args:
- envs_or_file:
If this argument is a file path, the serialized SCF object is saved to the file specified by this argument. If this attribute is a dict (created by locals()), the necessary variables are saved to the file specified by the attribute mf.chkfile.
- eig(h_kpts, s_kpts)#
Solver for generalized eigenvalue problem
\[HC = SCE\]
- get_hcore(cell=None, kpts=None)#
Get the core Hamiltonian AO matrices at sampled k-points.
- Args:
kpts : (nkpts, 3) ndarray
- Returns:
hcore : (nkpts, nao, nao) ndarray
- get_init_guess(cell=None, key='minao', s1e=None)#
- get_jk(cell=None, dm_kpts=None, hermi=1, kpts=None, kpts_band=None, with_j=True, with_k=True, omega=None, **kwargs)#
Get the Coulomb (J) and exchange (K) AO matrices at sampled k-points.
- Args:
- dm_kpts(nkpts, nao, nao) ndarray
Density matrix at each k-point. It needs to be Hermitian.
- Kwargs:
- kpts_band(3,) ndarray
A list of arbitrary “band” k-point at which to evalute the matrix.
- Returns:
vj : (nkpts, nao, nao) ndarray vk : (nkpts, nao, nao) ndarray or list of vj and vk if the input dm_kpts is a list of DMs
- get_occ(mo_energy_kpts=None, mo_coeff_kpts=None)#
Label the occupancies for each orbital for sampled k-points.
This is a k-point version of scf.hf.SCF.get_occ
- get_orbsym(mo_coeff=None, s=None)#
- get_ovlp(cell=None, kpts=None)#
Get the overlap AO matrices at sampled k-points.
- Args:
kpts : (nkpts, 3) ndarray
- Returns:
ovlp_kpts : (nkpts, nao, nao) ndarray
- 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
- ksan instance of
- Returns:
Veff : (nkpts, nao, nao) or (*, nkpts, nao, nao) ndarray Veff = J + Vxc.
- init_guess_by_chkfile(chk=None, project=None, kpts=None)#
Read the HF results from checkpoint file, then project it to the basis defined by
mol
- Kwargs:
- projectNone or bool
Whether to project chkfile’s orbitals to the new basis. Note when the geometry of the chkfile and the given molecule are very different, this projection can produce very poor initial guess. In PES scanning, it is recommended to switch off project.
If project is set to None, the projection is only applied when the basis sets of the chkfile’s molecule are different to the basis sets of the given molecule (regardless whether the geometry of the two molecules are different). Note the basis sets are considered to be different if the two molecules are derived from the same molecule with different ordering of atoms.
- Returns:
Density matrix, 2D ndarray
- property kpts#
- property orbsym#
- pyscf.pbc.dft.krks_ksymm.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)[source]#
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
- ksan instance of
- Returns:
Veff : (nkpts, nao, nao) or (*, nkpts, nao, nao) ndarray Veff = J + Vxc.
pyscf.pbc.dft.krkspu module#
Restricted DFT+U with kpoint sampling. Based on KRHF routine.
Refs: PRB, 1998, 57, 1505.
- class pyscf.pbc.dft.krkspu.KRKSpU(cell, kpts=array([[0., 0., 0.]]), xc='LDA,VWN', exxdiv='ewald', U_idx=[], U_val=[], C_ao_lo='minao', minao_ref='MINAO', **kwargs)[source]#
Bases:
KRKS
RKSpU class adapted for PBCs with k-point sampling.
- energy_elec(dm_kpts=None, h1e_kpts=None, vhf=None)#
Electronic energy for KRKSpU.
- get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)#
Coulomb + XC functional + Hubbard U terms.
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
- ksan instance of
- Returns:
Veff : (nkpts, nao, nao) or (*, nkpts, nao, nao) ndarray Veff = J + Vxc + V_U.
- to_hf(*args, **kwargs)#
Convert to KRHF object.
- pyscf.pbc.dft.krkspu.energy_elec(ks, dm_kpts=None, h1e_kpts=None, vhf=None)[source]#
Electronic energy for KRKSpU.
- pyscf.pbc.dft.krkspu.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)[source]#
Coulomb + XC functional + Hubbard U terms.
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
- ksan instance of
- Returns:
Veff : (nkpts, nao, nao) or (*, nkpts, nao, nao) ndarray Veff = J + Vxc + V_U.
- pyscf.pbc.dft.krkspu.mdot(*args)[source]#
Compute the dot product of a list of arrays in a single function call.
pyscf.pbc.dft.krkspu_ksymm module#
- pyscf.pbc.dft.krkspu_ksymm.KRKSpU#
alias of
KsymAdaptedKRKSpU
- class pyscf.pbc.dft.krkspu_ksymm.KsymAdaptedKRKSpU(cell, kpts=<pyscf.pbc.lib.kpts.KPoints object>, xc='LDA, VWN', exxdiv='ewald', U_idx=[], U_val=[], C_ao_lo='minao', minao_ref='MINAO', **kwargs)[source]#
Bases:
KsymAdaptedKRKS
RKSpU class adapted for PBCs with k-point sampling.
- energy_elec(dm_kpts=None, h1e_kpts=None, vhf=None)#
Electronic energy for KRKSpU.
- get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)#
Coulomb + XC functional + Hubbard U terms.
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
- ksan instance of
- Returns:
Veff : (nkpts, nao, nao) or (*, nkpts, nao, nao) ndarray Veff = J + Vxc + V_U.
- to_hf(*args, **kwargs)#
Convert to KRHF object.
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')[source]#
Bases:
KohnShamDFT
,KROHF
RKS class adapted for PBCs with k-point sampling.
- 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.
- 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.pbc.dft.kroks.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)[source]#
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')[source]#
Bases:
KohnShamDFT
,KUHF
UKS class adapted for PBCs with k-point sampling.
- 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.
- 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.pbc.dft.kuks.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)[source]#
Coulomb + XC functional for UKS. See pyscf/pbc/dft/uks.py
get_veff()
fore more details.
pyscf.pbc.dft.kuks_ksymm module#
- pyscf.pbc.dft.kuks_ksymm.KUKS#
alias of
KsymAdaptedKUKS
- class pyscf.pbc.dft.kuks_ksymm.KsymAdaptedKUKS(cell, kpts=<pyscf.pbc.lib.kpts.KPoints object>, xc='LDA, VWN', exxdiv='ewald', **kwargs)[source]#
Bases:
KUKS
,KsymAdaptedKUHF
- dump_chk(envs)#
Serialize the SCF object and save it to the specified chkfile.
- Args:
- envs_or_file:
If this argument is a file path, the serialized SCF object is saved to the file specified by this argument. If this attribute is a dict (created by locals()), the necessary variables are saved to the file specified by the attribute mf.chkfile.
- eig(h_kpts, s_kpts)#
Solver for generalized eigenvalue problem
\[HC = SCE\]
- get_hcore(cell=None, kpts=None)#
Get the core Hamiltonian AO matrices at sampled k-points.
- Args:
kpts : (nkpts, 3) ndarray
- Returns:
hcore : (nkpts, nao, nao) ndarray
- get_init_guess(cell=None, key='minao', s1e=None)#
- get_jk(cell=None, dm_kpts=None, hermi=1, kpts=None, kpts_band=None, with_j=True, with_k=True, omega=None, **kwargs)#
Get the Coulomb (J) and exchange (K) AO matrices at sampled k-points.
- Args:
- dm_kpts(nkpts, nao, nao) ndarray
Density matrix at each k-point. It needs to be Hermitian.
- Kwargs:
- kpts_band(3,) ndarray
A list of arbitrary “band” k-point at which to evalute the matrix.
- Returns:
vj : (nkpts, nao, nao) ndarray vk : (nkpts, nao, nao) ndarray or list of vj and vk if the input dm_kpts is a list of DMs
- get_occ(mo_energy_kpts=None, mo_coeff_kpts=None)#
Label the occupancies for each orbital for sampled k-points.
This is a k-point version of scf.hf.SCF.get_occ
- get_orbsym(mo_coeff=None, s=None)#
- get_ovlp(cell=None, kpts=None)#
Get the overlap AO matrices at sampled k-points.
- Args:
kpts : (nkpts, 3) ndarray
- Returns:
ovlp_kpts : (nkpts, nao, nao) ndarray
- 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.
- init_guess_by_chkfile(chk=None, project=None, kpts=None)#
Read the HF results from checkpoint file, then project it to the basis defined by
mol
- Kwargs:
- projectNone or bool
Whether to project chkfile’s orbitals to the new basis. Note when the geometry of the chkfile and the given molecule are very different, this projection can produce very poor initial guess. In PES scanning, it is recommended to switch off project.
If project is set to None, the projection is only applied when the basis sets of the chkfile’s molecule are different to the basis sets of the given molecule (regardless whether the geometry of the two molecules are different). Note the basis sets are considered to be different if the two molecules are derived from the same molecule with different ordering of atoms.
- Returns:
Density matrix, 2D ndarray
- property kpts#
- property nelec#
- property orbsym#
- pyscf.pbc.dft.kuks_ksymm.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)[source]#
Coulomb + XC functional for UKS. See pyscf/pbc/dft/uks.py
get_veff()
fore more details.
pyscf.pbc.dft.kukspu module#
Unrestricted DFT+U with kpoint sampling. Based on KUHF routine.
Refs: PRB, 1998, 57, 1505.
- class pyscf.pbc.dft.kukspu.KUKSpU(cell, kpts=array([[0., 0., 0.]]), xc='LDA,VWN', exxdiv='ewald', U_idx=[], U_val=[], C_ao_lo='minao', minao_ref='MINAO', **kwargs)[source]#
Bases:
KUKS
UKSpU class adapted for PBCs with k-point sampling.
- energy_elec(dm_kpts=None, h1e_kpts=None, vhf=None)#
Electronic energy for KUKSpU.
- get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)#
Coulomb + XC functional + (Hubbard - double counting) for KUKSpU.
- to_hf(*args, **kwargs)#
Convert to KUHF object.
pyscf.pbc.dft.kukspu_ksymm module#
- pyscf.pbc.dft.kukspu_ksymm.KUKSpU#
alias of
KsymAdaptedKUKSpU
- class pyscf.pbc.dft.kukspu_ksymm.KsymAdaptedKUKSpU(cell, kpts=<pyscf.pbc.lib.kpts.KPoints object>, xc='LDA, VWN', exxdiv='ewald', U_idx=[], U_val=[], C_ao_lo='minao', minao_ref='MINAO', **kwargs)[source]#
Bases:
KsymAdaptedKUKS
UKSpU class adapted for PBCs with k-point sampling.
- energy_elec(dm_kpts=None, h1e_kpts=None, vhf=None)#
Electronic energy for KUKSpU.
- get_veff(cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpts=None, kpts_band=None)#
Coulomb + XC functional + (Hubbard - double counting) for KUKSpU.
- to_hf(*args, **kwargs)#
Convert to KRHF object.
pyscf.pbc.dft.numint module#
- class pyscf.pbc.dft.numint.KNumInt(kpts=array([[0., 0., 0.]]))[source]#
Bases:
StreamObject
,LibXCMixin
Generalization of pyscf’s NumInt class for k-point sampling and periodic images.
- block_loop(cell, grids, nao=None, deriv=0, kpts=None, kpts_band=None, max_memory=2000, non0tab=None, blksize=None)[source]#
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.
- cache_xc_kernel1(cell, grids, xc_code, dm, 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.
- static eval_ao(cell, coords, kpts=None, deriv=0, relativity=0, shls_slice=None, non0tab=None, cutoff=None, out=None, verbose=None, **kwargs)#
- Returns:
- ao_kpts: (nkpts, [comp], ngrids, nao) ndarray
AO values at each k-point
- eval_rho(cell, ao_kpts, dm_kpts, non0tab=None, xctype='LDA', hermi=0, with_lapl=True, verbose=None)[source]#
Collocate the 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_rho1(cell, ao_kpts, dm_kpts, non0tab=None, xctype='LDA', hermi=0, with_lapl=True, cutoff=1e-15, grids=None, verbose=None)[source]#
- eval_rho2(cell, ao_kpts, mo_coeff_kpts, mo_occ_kpts, non0tab=None, xctype='LDA', with_lapl=True, verbose=None)[source]#
- get_rho(cell, dm, grids, kpts=array([[0., 0., 0.]]), max_memory=2000)#
Density in real space
- get_vxc(cell, grids, xc_code, dms, spin=0, relativity=0, hermi=1, kpts=None, kpts_band=None, max_memory=2000, verbose=None)#
Evaluate RKS/UKS XC functional and potential matrix. See
nr_rks()
andnr_uks()
for more details.
- make_mask(cell, coords, relativity=0, shls_slice=None, verbose=None)[source]#
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, relativity=0, hermi=1, kpts=None, kpts_band=None, max_memory=2000, verbose=None, **kwargs)[source]#
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.
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- spinint
spin polarized if spin = 1
- relativityint
No effects.
- hermiint
Whether the input density matrix is hermitian
- 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
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- hermiint
Whether the input density matrix is hermitian
- 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_rks_fxc_st(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
- nr_uks(cell, grids, xc_code, dms, relativity=0, hermi=1, kpts=None, kpts_band=None, max_memory=2000, verbose=None, **kwargs)[source]#
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.
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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
- gridsan instance of
- Kwargs:
- spinint
spin polarized if spin = 1
- relativityint
No effects.
- hermiint
Whether the input density matrix is hermitian
- 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
orKNumInt
cell : instance of
Mole
orCell
- 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
- gridsan instance of
- 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=1, kpts=None, kpts_band=None, max_memory=2000, verbose=None)[source]#
Evaluate RKS/UKS XC functional and potential matrix. See
nr_rks()
andnr_uks()
for more details.
- 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.
- class pyscf.pbc.dft.numint.NumInt[source]#
Bases:
StreamObject
,LibXCMixin
Generalization of pyscf’s NumInt class for a single k-point shift and periodic images.
- block_loop(cell, grids, nao=None, deriv=0, kpt=None, kpts_band=None, max_memory=2000, non0tab=None, blksize=None)[source]#
Define this macro to loop over grids by blocks.
- cutoff = 1e-13#
- static eval_ao(cell, coords, kpt=array([0., 0., 0.]), deriv=0, relativity=0, shls_slice=None, non0tab=None, cutoff=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
- static eval_rho(cell, ao, dm, non0tab=None, xctype='LDA', hermi=0, with_lapl=True, verbose=None)#
Collocate the density (opt. gradients) on the real-space grid.
- Args:
cell : instance of
Mole
orCell
- 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_rho1(cell, ao, dm, screen_index=None, xctype='LDA', hermi=0, with_lapl=True, cutoff=None, ao_cutoff=None, verbose=None)[source]#
- static eval_rho2(cell, ao, mo_coeff, mo_occ, non0tab=None, xctype='LDA', with_lapl=True, verbose=None)#
Refer to pyscf.dft.numint.eval_rho2 for full documentation.
- get_fxc(cell, grids, xc_code, dm0, dms, spin=0, relativity=0, hermi=0, rho0=None, vxc=None, fxc=None, kpt=None, max_memory=2000, verbose=None)#
Contract RKS XC kernel matrix with given density matrices
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- hermiint
Whether the input density matrix is hermitian
- 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:
- get_rho(cell, dm, grids, kpts=array([[0., 0., 0.]]), max_memory=2000)#
Density in real space
- get_vxc(cell, grids, xc_code, dms, spin=0, relativity=0, hermi=1, kpt=None, kpts_band=None, max_memory=2000, verbose=None)#
Evaluate RKS/UKS XC functional and potential matrix. See
nr_rks()
andnr_uks()
for more details.
- static make_mask(cell, coords, relativity=0, shls_slice=None, cutoff=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_fxc(cell, grids, xc_code, dm0, dms, spin=0, relativity=0, hermi=0, rho0=None, vxc=None, fxc=None, kpt=None, max_memory=2000, verbose=None)[source]#
Contract RKS XC kernel matrix with given density matrices
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- hermiint
Whether the input density matrix is hermitian
- 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_nlc_vxc(cell, grids, xc_code, dms, relativity=0, hermi=1, kpts=None, kpts_band=None, max_memory=2000, verbose=None)#
Calculate NLC functional and potential matrix on given grids
- Args:
ni : an instance of
NumInt
cell : an instance of
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.- dm2D array
Density matrix or multiple density matrices
- gridsan instance of
- Kwargs:
- hermiint
Input density matrices symmetric or not. It also indicates whether the potential matrices in return are symmetric or not.
- max_memoryint or float
The maximum size of cache to use (in MB).
- 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(cell, grids, xc_code, dms, relativity=0, hermi=1, kpt=array([0., 0., 0.]), kpts_band=None, max_memory=2000, verbose=None)[source]#
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.
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- spinint
spin polarized if spin = 1
- relativityint
No effects.
- hermiint
Whether the input density matrix is hermitian
- 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
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- hermiint
Whether the input density matrix is hermitian
- 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_rks_fxc_st(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
- nr_uks(cell, grids, xc_code, dms, relativity=0, hermi=1, kpt=array([0., 0., 0.]), kpts_band=None, max_memory=2000, verbose=None)[source]#
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.
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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
- gridsan instance of
- Kwargs:
- spinint
spin polarized if spin = 1
- relativityint
No effects.
- hermiint
Whether the input density matrix is hermitian
- 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
orKNumInt
cell : instance of
Mole
orCell
- 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
- gridsan instance of
- 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=1, kpt=None, kpts_band=None, max_memory=2000, verbose=None)[source]#
Evaluate RKS/UKS XC functional and potential matrix. See
nr_rks()
andnr_uks()
for more details.
- 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.pbc.dft.numint.cache_xc_kernel(ni, cell, grids, xc_code, mo_coeff, mo_occ, spin=0, kpts=None, max_memory=2000)[source]#
Compute the 0th order density, Vxc and fxc. They can be used in TDDFT, DFT hessian module etc.
- pyscf.pbc.dft.numint.cache_xc_kernel1(ni, cell, grids, xc_code, dm, spin=0, kpts=None, max_memory=2000)[source]#
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, cutoff=None, out=None, verbose=None)[source]#
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, cutoff=None, out=None, verbose=None, **kwargs)[source]#
- 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, with_lapl=True, verbose=None)[source]#
Collocate the density (opt. gradients) on the real-space grid.
- Args:
cell : instance of
Mole
orCell
- 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', with_lapl=True, verbose=None)[source]#
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)[source]#
Density in real space
- pyscf.pbc.dft.numint.nr_nlc_vxc(ni, cell, grids, xc_code, dms, relativity=0, hermi=1, kpts=None, kpts_band=None, max_memory=2000, verbose=None)[source]#
Calculate NLC functional and potential matrix on given grids
- Args:
ni : an instance of
NumInt
cell : an instance of
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.- dm2D array
Density matrix or multiple density matrices
- gridsan instance of
- Kwargs:
- hermiint
Input density matrices symmetric or not. It also indicates whether the potential matrices in return are symmetric or not.
- max_memoryint or float
The maximum size of cache to use (in MB).
- 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(ni, cell, grids, xc_code, dms, spin=0, relativity=0, hermi=1, kpts=None, kpts_band=None, max_memory=2000, verbose=None)[source]#
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.
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- spinint
spin polarized if spin = 1
- relativityint
No effects.
- hermiint
Whether the input density matrix is hermitian
- 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)[source]#
Contract RKS XC kernel matrix with given density matrices
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- hermiint
Whether the input density matrix is hermitian
- 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)[source]#
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=1, 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.
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- spinint
spin polarized if spin = 1
- relativityint
No effects.
- hermiint
Whether the input density matrix is hermitian
- 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=1, kpts=None, kpts_band=None, max_memory=2000, verbose=None)[source]#
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.
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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
- gridsan instance of
- Kwargs:
- spinint
spin polarized if spin = 1
- relativityint
No effects.
- hermiint
Whether the input density matrix is hermitian
- 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)[source]#
Contract UKS XC kernel matrix with given density matrices
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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
- gridsan instance of
- 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=1, 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.
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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
- gridsan instance of
- Kwargs:
- spinint
spin polarized if spin = 1
- relativityint
No effects.
- hermiint
Whether the input density matrix is hermitian
- 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.numint2c module#
Numerical integration functions for (2-component) GKS and KGKS
- Ref:
Phys. Rev. Research 5, 013036
- class pyscf.pbc.dft.numint2c.KNumInt2C(kpts=array([[0., 0., 0.]]))[source]#
Bases:
StreamObject
,LibXCMixin
- block_loop(cell, grids, nao=None, deriv=0, kpts=None, 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_kpts, mo_occ_kpts, spin=0, kpts=None, max_memory=2000)[source]#
Compute the 0th order density, Vxc and fxc. They can be used in TDDFT, DFT hessian module etc.
- cache_xc_kernel1(cell, grids, xc_code, dm_kpts, spin=0, kpts=None, max_memory=2000)[source]#
Compute the 0th order density, Vxc and fxc. They can be used in TDDFT, DFT hessian module etc.
- collinear = 'col'#
- collinear_samples = 200#
- collinear_thrd = 0.99#
- static eval_ao(cell, coords, kpts=None, deriv=0, relativity=0, shls_slice=None, non0tab=None, cutoff=None, out=None, verbose=None, **kwargs)#
- Returns:
- ao_kpts: (nkpts, [comp], ngrids, nao) ndarray
AO values at each k-point
- eval_rho(cell, ao_kpts, dm_kpts, non0tab=None, xctype='LDA', hermi=0, with_lapl=True, verbose=None)[source]#
Collocate the 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_rho1(cell, ao_kpts, dm_kpts, screen_index=None, xctype='LDA', hermi=0, with_lapl=True, cutoff=None, ao_cutoff=None, pair_mask=None, verbose=None)[source]#
- eval_rho2(cell, ao_kpts, mo_coeff_kpts, mo_occ_kpts, non0tab=None, xctype='LDA', with_lapl=True, verbose=None)[source]#
- eval_xc_eff(xc_code, rho, deriv=1, omega=None, xctype=None, verbose=None)#
Returns the derivative tensor against the density parameters
- [[density_a, (nabla_x)_a, (nabla_y)_a, (nabla_z)_a, tau_a],
[density_b, (nabla_x)_b, (nabla_y)_b, (nabla_z)_b, tau_b]].
It differs from the eval_xc method in the derivatives of non-local part. The eval_xc method returns the XC functional derivatives to sigma (|nabla rho|^2)
- Args:
- rho: 2-dimensional or 3-dimensional array
Total density and spin density (and their derivatives if GGA or MGGA functionals) on grids
- Kwargs:
- deriv: int
derivative orders
- omega: float
define the exponent in the attenuated Coulomb for RSH functional
- get_fxc(cell, grids, xc_code, dm0, dms, spin=0, relativity=0, hermi=0, rho0=None, vxc=None, fxc=None, kpt=None, max_memory=2000, verbose=None)#
Contract RKS XC kernel matrix with given density matrices
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- hermiint
Whether the input density matrix is hermitian
- 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:
- get_rho(cell, dm, grids, kpts=array([[0., 0., 0.]]), max_memory=2000)[source]#
Density in real space
- get_vxc(cell, grids, xc_code, dms, spin=0, relativity=0, hermi=1, kpts=None, kpts_band=None, max_memory=2000, verbose=None)#
Evaluate RKS/UKS XC functional and potential matrix. See
nr_rks()
andnr_uks()
for more details.
- make_mask(*args, **kwargs)#
- mcfun_eval_xc_adapter(xc_code)#
Wrapper to generate the eval_xc function required by mcfun
- nr_fxc(cell, grids, xc_code, dm0, dms, spin=0, relativity=0, hermi=0, rho0=None, vxc=None, fxc=None, kpt=None, max_memory=2000, verbose=None)#
Contract RKS XC kernel matrix with given density matrices
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- hermiint
Whether the input density matrix is hermitian
- 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_gks_fxc(cell, grids, xc_code, dm0, dms, spin=0, relativity=0, hermi=0, rho0=None, vxc=None, fxc=None, kpt=None, max_memory=2000, verbose=None)#
Contract RKS XC kernel matrix with given density matrices
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- hermiint
Whether the input density matrix is hermitian
- 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_gks_vxc(cell, grids, xc_code, dms, spin=0, relativity=0, hermi=1, kpts=None, kpts_band=None, max_memory=2000, verbose=None)#
Evaluate RKS/UKS XC functional and potential matrix. See
nr_rks()
andnr_uks()
for more details.
- nr_vxc(cell, grids, xc_code, dms, spin=0, relativity=0, hermi=1, kpts=None, kpts_band=None, max_memory=2000, verbose=None)[source]#
Evaluate RKS/UKS XC functional and potential matrix. See
nr_rks()
andnr_uks()
for more details.
- spin_samples = 770#
- class pyscf.pbc.dft.numint2c.NumInt2C[source]#
Bases:
StreamObject
,LibXCMixin
Numerical integration methods for 2-component basis (used by GKS)
- block_loop(cell, grids, nao=None, deriv=0, kpt=None, 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, kpt=None, max_memory=2000)[source]#
Compute the 0th order density, Vxc and fxc. They can be used in TDDFT, DFT hessian module etc.
- cache_xc_kernel1(cell, grids, xc_code, dm, spin=0, kpt=None, max_memory=2000)[source]#
Compute the 0th order density, Vxc and fxc. They can be used in TDDFT, DFT hessian module etc.
- collinear = 'col'#
- collinear_samples = 200#
- collinear_thrd = 0.99#
- static eval_ao(cell, coords, kpt=array([0., 0., 0.]), deriv=0, relativity=0, shls_slice=None, non0tab=None, cutoff=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
- static eval_rho(mol, ao, dm, non0tab=None, xctype='LDA', hermi=0, with_lapl=True, verbose=None)#
Calculate the electron density and magnetization spin density in the framework of 2-component real basis. Calculate the electron density for LDA functional, and the density derivatives for GGA and MGGA functionals.
- Args:
mol : an instance of
Mole
- ao2D array of shape (N,nao) for LDA, 3D array of shape (4,N,nao) for GGA
or meta-GGA. N is the number of grids, nao is the number of AO functions. If xctype is GGA (MGGA), ao[0] is AO value and ao[1:3] are the AO gradients. ao[4:10] are second derivatives of ao values if applicable.
- 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.
- hermibool
dm is hermitian or not
- 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; For meta-GGA, returns can be a (6,N) (with_lapl=True) array where last two rows are nabla^2 rho and tau = 1/2(nabla f)^2 or (5,N) (with_lapl=False) where the last row is tau = 1/2(nabla f)^2
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, deriv=0) >>> dm = numpy.random.random((mol.nao_nr(),mol.nao_nr())) >>> dm = dm + dm.T >>> rho, dx_rho, dy_rho, dz_rho = eval_rho(mol, ao, dm, xctype='LDA')
- eval_rho1(cell, ao, dm, screen_index=None, xctype='LDA', hermi=0, with_lapl=True, cutoff=None, ao_cutoff=None, pair_mask=None, verbose=None)[source]#
- eval_rho2(mol, ao, mo_coeff, mo_occ, non0tab=None, xctype='LDA', with_lapl=True, verbose=None)[source]#
Calculate the electron density for LDA functional and the density derivatives for GGA functional in the framework of 2-component basis.
- eval_xc_eff(xc_code, rho, deriv=1, omega=None, xctype=None, verbose=None)#
Returns the derivative tensor against the density parameters
- [[density_a, (nabla_x)_a, (nabla_y)_a, (nabla_z)_a, tau_a],
[density_b, (nabla_x)_b, (nabla_y)_b, (nabla_z)_b, tau_b]].
It differs from the eval_xc method in the derivatives of non-local part. The eval_xc method returns the XC functional derivatives to sigma (|nabla rho|^2)
- Args:
- rho: 2-dimensional or 3-dimensional array
Total density and spin density (and their derivatives if GGA or MGGA functionals) on grids
- Kwargs:
- deriv: int
derivative orders
- omega: float
define the exponent in the attenuated Coulomb for RSH functional
- get_fxc(cell, grids, xc_code, dm0, dms, spin=0, relativity=0, hermi=0, rho0=None, vxc=None, fxc=None, kpt=None, max_memory=2000, verbose=None)#
Contract RKS XC kernel matrix with given density matrices
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- hermiint
Whether the input density matrix is hermitian
- 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:
- get_vxc(cell, grids, xc_code, dms, spin=0, relativity=0, hermi=1, kpt=None, kpts_band=None, max_memory=2000, verbose=None)#
Evaluate RKS/UKS XC functional and potential matrix. See
nr_rks()
andnr_uks()
for more details.
- make_mask(*args, **kwargs)#
- mcfun_eval_xc_adapter(xc_code)#
Wrapper to generate the eval_xc function required by mcfun
- nr_fxc(cell, grids, xc_code, dm0, dms, spin=0, relativity=0, hermi=0, rho0=None, vxc=None, fxc=None, kpt=None, max_memory=2000, verbose=None)[source]#
Contract RKS XC kernel matrix with given density matrices
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- hermiint
Whether the input density matrix is hermitian
- 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_gks_fxc(cell, grids, xc_code, dm0, dms, spin=0, relativity=0, hermi=0, rho0=None, vxc=None, fxc=None, kpt=None, max_memory=2000, verbose=None)#
Contract RKS XC kernel matrix with given density matrices
- Args:
ni : an instance of
NumInt
orKNumInt
cell : instance of
Mole
orCell
- 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)
- gridsan instance of
- Kwargs:
- hermiint
Whether the input density matrix is hermitian
- 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_gks_vxc(cell, grids, xc_code, dms, spin=0, relativity=0, hermi=1, kpt=None, kpts_band=None, max_memory=2000, verbose=None)#
Evaluate RKS/UKS XC functional and potential matrix. See
nr_rks()
andnr_uks()
for more details.
- nr_vxc(cell, grids, xc_code, dms, spin=0, relativity=0, hermi=1, kpt=None, kpts_band=None, max_memory=2000, verbose=None)[source]#
Evaluate RKS/UKS XC functional and potential matrix. See
nr_rks()
andnr_uks()
for more details.
- spin_samples = 770#
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')[source]#
Bases:
KohnShamDFT
PBC-KS
- density_fit(auxbasis=None, with_df=None, *args, **kwargs)#
- get_rho(dm=None, grids=None, kpt=None)#
Compute density in real space
- initialize_grids(cell, dm, kpts, ground_state=True)[source]#
Initialize self.grids the first time call get_veff
- jk_method(J='FFTDF', K=None)[source]#
Set up the schemes to evaluate Coulomb and exchange matrix
FFTDF: planewave density fitting using Fast Fourier Transform AFTDF: planewave density fitting using analytic Fourier Transform GDF: Gaussian density fitting MDF: Gaussian and planewave mix density fitting RS: range-separation JK builder RSDF: range-separation density fitting
- mix_density_fit(auxbasis=None, with_df=None, *args, **kwargs)#
- rs_density_fit(auxbasis=None, with_df=None, *args, **kwargs)#
- to_gks(xc=None)[source]#
Convert the input mean-field object to a GKS object.
Note this conversion only changes the class of the mean-field object. The total energy and wave-function are the same as them in the input mean-field object.
- to_hf()[source]#
Convert the input KS object to the associated HF object.
Note this conversion only changes the class of the mean-field object. The total energy and wave-function are the same as them in the input mean-field object.
- class pyscf.pbc.dft.rks.RKS(cell, kpt=array([0., 0., 0.]), xc='LDA,VWN', exxdiv='ewald')[source]#
Bases:
KohnShamDFT
,RHF
RKS class adapted for PBCs.
This is a literal duplication of the molecular RKS class with some mol variables replaced by cell.
- 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-particle density matrix
- h1e2D ndarray
Core hamiltonian
- Returns:
RKS electronic energy and the 2-electron contribution
- 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
- ksan instance of
- Returns:
matrix Veff = J + Vxc. Veff can be a list matrices, if the input dm is a list of density matrices.
- get_vsap(mol=None)#
Superposition of atomic potentials
S. Lehtola, Assessment of initial guesses for self-consistent field calculations. Superposition of Atomic Potentials: simple yet efficient, J. Chem. Theory Comput. 15, 1593 (2019). DOI: 10.1021/acs.jctc.8b01089. arXiv:1810.11659.
This function evaluates the effective charge of a neutral atom, given by exchange-only LDA on top of spherically symmetric unrestricted Hartree-Fock calculations as described in
S. Lehtola, L. Visscher, E. Engel, Efficient implementation of the superposition of atomic potentials initial guess for electronic structure calculations in Gaussian basis sets, J. Chem. Phys., in press (2020).
The potentials have been calculated for the ground-states of spherically symmetric atoms at the non-relativistic level of theory as described in
S. Lehtola, “Fully numerical calculations on atoms with fractional occupations and range-separated exchange functionals”, Phys. Rev. A 101, 012516 (2020). DOI: 10.1103/PhysRevA.101.012516
using accurate finite-element calculations as described in
S. Lehtola, “Fully numerical Hartree-Fock and density functional calculations. I. Atoms”, Int. J. Quantum Chem. e25945 (2019). DOI: 10.1002/qua.25945
Note
This function will modify the input ks object.
- Args:
- ksan instance of
RKS
XC functional are controlled by ks.xc attribute. Attribute ks.grids might be initialized.
- ksan instance of
- Returns:
matrix Vsap = Vnuc + J + Vxc.
- init_guess_by_vsap(mol=None)#
Form SAP guess
- 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.pbc.dft.rks.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)[source]#
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
- ksan instance of
- Returns:
matrix Veff = J + Vxc. Veff can be a list matrices, if the input dm is a list of density matrices.
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')[source]#
Bases:
KohnShamDFT
,ROHF
UKS class adapted for PBCs.
This is a literal duplication of the molecular UKS class with some mol variables replaced by cell.
- 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-particle 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_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.
- get_vsap(mol=None)#
Superposition of atomic potentials
S. Lehtola, Assessment of initial guesses for self-consistent field calculations. Superposition of Atomic Potentials: simple yet efficient, J. Chem. Theory Comput. 15, 1593 (2019). DOI: 10.1021/acs.jctc.8b01089. arXiv:1810.11659.
This function evaluates the effective charge of a neutral atom, given by exchange-only LDA on top of spherically symmetric unrestricted Hartree-Fock calculations as described in
S. Lehtola, L. Visscher, E. Engel, Efficient implementation of the superposition of atomic potentials initial guess for electronic structure calculations in Gaussian basis sets, J. Chem. Phys., in press (2020).
The potentials have been calculated for the ground-states of spherically symmetric atoms at the non-relativistic level of theory as described in
S. Lehtola, “Fully numerical calculations on atoms with fractional occupations and range-separated exchange functionals”, Phys. Rev. A 101, 012516 (2020). DOI: 10.1103/PhysRevA.101.012516
using accurate finite-element calculations as described in
S. Lehtola, “Fully numerical Hartree-Fock and density functional calculations. I. Atoms”, Int. J. Quantum Chem. e25945 (2019). DOI: 10.1002/qua.25945
Note
This function will modify the input ks object.
- Args:
- ksan instance of
RKS
XC functional are controlled by ks.xc attribute. Attribute ks.grids might be initialized.
- ksan instance of
- Returns:
matrix Vsap = Vnuc + J + Vxc.
- init_guess_by_vsap(mol=None)#
Form SAP guess
- 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.pbc.dft.roks.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)[source]#
Coulomb + XC functional for UKS. See pyscf/pbc/dft/uks.py
get_veff()
fore more details.
pyscf.pbc.dft.uks module#
Non-relativistic unrestricted 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')[source]#
Bases:
KohnShamDFT
,UHF
UKS class adapted for PBCs.
This is a literal duplication of the molecular UKS class with some mol variables replaced by cell.
- 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.
- get_vsap(mol=None)#
Superposition of atomic potentials
S. Lehtola, Assessment of initial guesses for self-consistent field calculations. Superposition of Atomic Potentials: simple yet efficient, J. Chem. Theory Comput. 15, 1593 (2019). DOI: 10.1021/acs.jctc.8b01089. arXiv:1810.11659.
This function evaluates the effective charge of a neutral atom, given by exchange-only LDA on top of spherically symmetric unrestricted Hartree-Fock calculations as described in
S. Lehtola, L. Visscher, E. Engel, Efficient implementation of the superposition of atomic potentials initial guess for electronic structure calculations in Gaussian basis sets, J. Chem. Phys., in press (2020).
The potentials have been calculated for the ground-states of spherically symmetric atoms at the non-relativistic level of theory as described in
S. Lehtola, “Fully numerical calculations on atoms with fractional occupations and range-separated exchange functionals”, Phys. Rev. A 101, 012516 (2020). DOI: 10.1103/PhysRevA.101.012516
using accurate finite-element calculations as described in
S. Lehtola, “Fully numerical Hartree-Fock and density functional calculations. I. Atoms”, Int. J. Quantum Chem. e25945 (2019). DOI: 10.1002/qua.25945
Note
This function will modify the input ks object.
- Args:
- ksan instance of
RKS
XC functional are controlled by ks.xc attribute. Attribute ks.grids might be initialized.
- ksan instance of
- Returns:
matrix Vsap = Vnuc + J + Vxc.
- init_guess_by_vsap(mol=None)#
Form SAP guess
- 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.pbc.dft.uks.get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None)[source]#
Coulomb + XC functional for UKS. See pyscf/pbc/dft/uks.py
get_veff()
fore more details.
Module contents#
Kohn-Sham DFT for periodic systems
- pyscf.pbc.dft.KKS(cell, *args, **kwargs)[source]#
A wrap function to create DFT object with k-point sampling (KRKS or KUKS).
RKS class adapted for PBCs with k-point sampling.