pyscf.pbc.gto.pseudo package#
Submodules#
pyscf.pbc.gto.pseudo.pp module#
PP with numeric integration. See also pyscf/pbc/gto/pesudo/pp_int.py
- For GTH/HGH PPs, see:
Goedecker, Teter, Hutter, PRB 54, 1703 (1996) Hartwigsen, Goedecker, and Hutter, PRB 58, 3641 (1998)
- pyscf.pbc.gto.pseudo.pp.get_alphas(cell)[source]#
alpha parameters from the non-divergent Hartree+Vloc G=0 term.
See ewald.pdf
- Returns:
alphas : (natm,) ndarray
- pyscf.pbc.gto.pseudo.pp.get_alphas_gth(cell)[source]#
alpha parameters for the local GTH pseudopotential.
- pyscf.pbc.gto.pseudo.pp.get_gth_projG(cell, Gvs)[source]#
G space projectors from the FT of the real-space projectors.
int e^{iGr} p_j^l(r) Y_{lm}^*(theta,phi) = i^l p_j^l(G) Y_{lm}^*(thetaG, phiG)
See MH Eq.(4.80)
- pyscf.pbc.gto.pseudo.pp.get_gth_vlocG(cell, Gv)[source]#
Local part of the GTH pseudopotential.
See MH (4.79).
- Args:
Gv : (ngrids,3) ndarray
- Returns:
(natm, ngrids) ndarray
- pyscf.pbc.gto.pseudo.pp.get_jvloc_G0(cell, kpt=array([0., 0., 0.]))[source]#
Get the (separately divergent) Hartree + Vloc G=0 contribution.
- pyscf.pbc.gto.pseudo.pp.get_pp(cell, kpt=array([0., 0., 0.]))[source]#
Get the periodic pseudopotential nuc-el AO matrix
- pyscf.pbc.gto.pseudo.pp.get_projG(cell, kpt=array([0., 0., 0.]))[source]#
PP weight and projector for the nonlocal PP in G space.
- Returns:
- hslist( list( np.array( , ) ) )
hs[atm][l][i,j]
- projslist( list( list( list( np.array(ngrids) ) ) ) )
projs[atm][l][m][i][ngrids]
pyscf.pbc.gto.pseudo.pp_int module#
Analytic PP integrals. See also pyscf/pbc/gto/pesudo/pp.py
- For GTH/HGH PPs, see:
Goedecker, Teter, Hutter, PRB 54, 1703 (1996) Hartwigsen, Goedecker, and Hutter, PRB 58, 3641 (1998)
- pyscf.pbc.gto.pseudo.pp_int.fake_cell_vloc(cell, cn=0, atm_id=None)[source]#
Generate fake cell for V_{loc}.
Each term of V_{loc} (erf, C_1, C_2, C_3, C_4) is a gaussian type function. The integral over V_{loc} can be transfered to the 3-center integrals, in which the auxiliary basis is given by the fake cell.
The kwarg cn indiciates which term to generate for the fake cell. If cn = 0, the erf term is generated. C_1,..,C_4 are generated with cn = 1..4
- pyscf.pbc.gto.pseudo.pp_int.fake_cell_vnl(cell)[source]#
Generate fake cell for V_{nl}.
gaussian function p_i^l Y_{lm}
- pyscf.pbc.gto.pseudo.pp_int.get_pp_loc_part1(cell, kpts=None)[source]#
PRB, 58, 3641 Eq (1), integrals associated to erf
- pyscf.pbc.gto.pseudo.pp_int.get_pp_loc_part2(cell, kpts=None)[source]#
PRB, 58, 3641 Eq (1), integrals associated to C1, C2, C3, C4
pyscf.pbc.gto.pseudo.ppnl_velgauge module#
Analytic PP integrals for GTH/HGH PPs in velocity gauge.
- For GTH/HGH PPs, see:
Goedecker, Teter, Hutter, PRB 54, 1703 (1996) Hartwigsen, Goedecker, and Hutter, PRB 58, 3641 (1998)
For the velocity gauge transformation, see: [1] Comparison of Length, Velocity, and Symmetric Gauges for the Calculation
of Absorption and Electric Circular Dichroism Spectra with Real-Time Time-Dependent Density Functional Theory, Johann Mattiat and Sandra Luber Journal of Chemical Theory and Computation 2022 18 (9), 5513-5526, DOI: 10.1021/acs.jctc.2c00644
- class pyscf.pbc.gto.pseudo.ppnl_velgauge.VelGaugePPNLHelper(cell, kpts=None, intors=None, hl_max=3, origin=(0.0, 0.0, 0.0))[source]#
Bases:
objectHelper class for evaluating velocity gauge pseudopotential non-local integrals. Useful to avoid recomputing data that only depends on the cell and k-points.
- pyscf.pbc.gto.pseudo.ppnl_velgauge.ft_aopair_kpts_kern(cell, aosym='s1', kptjs=array([[0., 0., 0.]]), intor='GTO_ft_ovlp', comp=1, bvk_kmesh=None, origin=(0.0, 0.0, 0.0))[source]#
Fourier transform AO pair for a group of k-points sum_T exp(-i k_j * T) int exp(-i(G+q)r) i(r) j(r-T) dr^3
Modified version of pyscf.pbc.df.ft_ao.ft_aopair_kpts that returns the generated ft kernel.
- pyscf.pbc.gto.pseudo.ppnl_velgauge.get_gth_pp_nl_velgauge(cell, q, kpts=None, vgppnl_helper=None)[source]#
Nonlocal part of GTH pseudopotential in velocity gauge.
Parameters#
- cellpyscf.pbc.gto.Cell
System cell
- A_over_cnp.ndarray
Scaled magnetic vector potential. Shape is (3,)
- kptsnp.ndarray, optional
k-point list.
- vgppnl_helperVelGaugePPNLHelper, optional
Helper object for velocity gauge PP integrals. By default None which causes a new helper object to be created and built.
Returns#
- tuple
(ppnl, vgppnl_helper) where ppnl is an ndarray of shape (nkpts, nao, nao) and vgppnl_helper is the VelGaugePPNLHelper object used to compute the integrals.
- pyscf.pbc.gto.pseudo.ppnl_velgauge.get_gth_pp_nl_velgauge_commutator(cell, q, kpts=None, vgppnl_helper=None)[source]#
- pyscf.pbc.gto.pseudo.ppnl_velgauge.prepare_ppnl_ft_data(cell, fakecell, hl_idx, hl_blocks, kpts, intor, origin=(0.0, 0.0, 0.0), comp=1)[source]#
Prepare ft_kernel methods for fast evaluation of velocity gauge ppnl integrals
Parameters#
- cellpyscf.pbc.gto.Cell
System cell
- fakecellpyscf.pbc.gto.Cell
Fake cell containing GTH projectors
- hl_idxint
GTH projector angular momentum index
- hl_blockslist
GTH hl blocks
- kptsnp.ndarray
k-points
- intorstr
GTO ft-ao integral name
- compint, optional
Size of each integral (eg scalar=1, vector=3), by default 1
Returns#
- tuple
shls_slice, ft_kern, cell_conc_fakecell