pyscf.pbc.grad package#
Submodules#
pyscf.pbc.grad.krhf module#
Non-relativistic analytical nuclear gradients for restricted Hartree Fock with kpoints sampling
- class pyscf.pbc.grad.krhf.Gradients(method)[source]#
Bases:
GradientsBaseNon-relativistic restricted Hartree-Fock gradients
- as_scanner()#
Generating a nuclear gradients scanner/solver (for geometry optimizer).
The returned solver is a function. This function requires one argument “cell” as input and returns energy and first order nuclear derivatives.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the nuc-grad object and SCF object (DIIS, conv_tol, max_memory etc) are automatically applied in the solver.
Note scanner has side effects. It may change many underlying objects (_scf, with_df, with_x2c, …) during calculation.
- extra_force(atom_id, envs)[source]#
Hook for extra contributions in analytical gradients.
Contributions like the response of auxiliary basis in density fitting method, the grid response in DFT numerical integration can be put in this function.
- grad_elec(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)#
Electronic part of KRHF/KRKS gradients Args:
mf_grad : pbc.grad.krhf.Gradients or pbc.grad.krks.Gradients object
- kernel(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)[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.).
- class pyscf.pbc.grad.krhf.GradientsBase(method)[source]#
Bases:
GradientsBaseBasic nuclear gradient functions for non-relativistic methods
- hcore_generator(cell=None, kpts=None)#
- property kpts#
- class pyscf.pbc.grad.krhf.SCF_GradScanner(g)[source]#
Bases:
GradScanner
- pyscf.pbc.grad.krhf.as_scanner(mf_grad)[source]#
Generating a nuclear gradients scanner/solver (for geometry optimizer).
The returned solver is a function. This function requires one argument “cell” as input and returns energy and first order nuclear derivatives.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the nuc-grad object and SCF object (DIIS, conv_tol, max_memory etc) are automatically applied in the solver.
Note scanner has side effects. It may change many underlying objects (_scf, with_df, with_x2c, …) during calculation.
- pyscf.pbc.grad.krhf.get_hcore(cell, kpts)[source]#
Part of the nuclear gradients of core Hamiltonian
- pyscf.pbc.grad.krhf.get_jk(mf_grad, dm, kpts)[source]#
J = ((-nabla i) j| kl) D_lk K = ((-nabla i) j| kl) D_jk
- pyscf.pbc.grad.krhf.grad_elec(mf_grad, mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)[source]#
Electronic part of KRHF/KRKS gradients Args:
mf_grad : pbc.grad.krhf.Gradients or pbc.grad.krks.Gradients object
pyscf.pbc.grad.krks module#
Non-relativistic analytical nuclear gradients for restricted Kohn Sham with kpoints sampling
pyscf.pbc.grad.krks_stress module#
The energy derivatives for the strain tensor e_ij is
1 d E
- sigma_ij = — ——
V d e_ij
The strain tesnor e_ij describes the transformation for real space coordinates in the crystal
sum_j [deta_ij + e_ij] R_j [for j = x, y, z]
The strain tensor is generally not a symmetric tensor. Symmetrization
[e1 e6/2 e5/2] [e6/2 e2 e4/2] [e5/2 e4/2 e3 ]
is applied to form 6 independent component.
e1 = e_11 e2 = e_22 e3 = e_33 e6 = e_12 + e_21 e5 = e_13 + e_31 e4 = e_32 + e_23
The 6 component strain is then used to define the symmetric stress tensor.
1 d E
- sigma_i = — —— for i = 1 .. 6
V d e_i
The symmetric stress tensor represented in the 6 Voigt notation can be transformed from the asymmetric stress tensor sigma_ij
sigma1 = sigma_11 sigma2 = sigma_22 sigma3 = sigma_33 sigma6 = (sigma_12 + sigma_21)/2 sigma5 = (sigma_13 + sigma_31)/2 sigma4 = (sigma_23 + sigma_32)/2
See K. Doll, Mol Phys (2010), 108, 223
- pyscf.pbc.grad.krks_stress.get_vxc(ks_grad, cell, dm_kpts, kpts, with_j=False, with_nuc=False)[source]#
Strain derivatives for Coulomb and XC at gamma point
- Kwargs:
with_j : Whether to include the electron-electron Coulomb interactions with_nuc : Whether to include the electron-nuclear Coulomb interactions
- pyscf.pbc.grad.krks_stress.kernel(mf_grad)[source]#
Compute the energy derivatives for strain tensor (e_ij)
1 d E
- sigma_ij = — ——
V d e_ij
sigma is a asymmetric 3x3 matrix. The symmetric stress tensor in the 6 Voigt notation can be transformed from the asymmetric stress tensor
sigma1 = sigma_11 sigma2 = sigma_22 sigma3 = sigma_33 sigma6 = (sigma_12 + sigma_21)/2 sigma5 = (sigma_13 + sigma_31)/2 sigma4 = (sigma_23 + sigma_32)/2
See K. Doll, Mol Phys (2010), 108, 223
pyscf.pbc.grad.krkspu module#
Analytical derivatives for DFT+U with kpoints sampling
- class pyscf.pbc.grad.krkspu.Gradients(mf)[source]#
Bases:
Gradients
pyscf.pbc.grad.kuhf module#
Non-relativistic analytical nuclear gradients for unrestricted Hartree Fock with kpoints sampling
- class pyscf.pbc.grad.kuhf.Gradients(method)[source]#
Bases:
GradientsBaseNon-relativistic restricted Hartree-Fock gradients
- as_scanner()#
Generating a nuclear gradients scanner/solver (for geometry optimizer).
The returned solver is a function. This function requires one argument “cell” as input and returns energy and first order nuclear derivatives.
The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the nuc-grad object and SCF object (DIIS, conv_tol, max_memory etc) are automatically applied in the solver.
Note scanner has side effects. It may change many underlying objects (_scf, with_df, with_x2c, …) during calculation.
- extra_force(atom_id, envs)[source]#
Hook for extra contributions in analytical gradients.
Contributions like the response of auxiliary basis in density fitting method, the grid response in DFT numerical integration can be put in this function.
- grad_elec(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)#
- kernel(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)[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.).
pyscf.pbc.grad.kuks module#
Non-relativistic analytical nuclear gradients for unrestricted Kohn Sham with kpoints sampling
pyscf.pbc.grad.kuks_stress module#
- pyscf.pbc.grad.kuks_stress.get_vxc(ks_grad, cell, dm_kpts, kpts, with_j=False, with_nuc=False)[source]#
Strain derivatives for Coulomb and XC at gamma point
- Kwargs:
with_j : Whether to include the electron-electron Coulomb interactions with_nuc : Whether to include the electron-nuclear Coulomb interactions
- pyscf.pbc.grad.kuks_stress.kernel(mf_grad)[source]#
Compute the energy derivatives for strain tensor (e_ij)
1 d E
- sigma_ij = — ——
V d e_ij
sigma is a asymmetric 3x3 matrix. The symmetric stress tensor in the 6 Voigt notation can be transformed from the asymmetric stress tensor
sigma1 = sigma_11 sigma2 = sigma_22 sigma3 = sigma_33 sigma6 = (sigma_12 + sigma_21)/2 sigma5 = (sigma_13 + sigma_31)/2 sigma4 = (sigma_23 + sigma_32)/2
See K. Doll, Mol Phys (2010), 108, 223
pyscf.pbc.grad.kukspu module#
Analytical derivatives for DFT+U with kpoints sampling
- class pyscf.pbc.grad.kukspu.Gradients(mf)[source]#
Bases:
Gradients
pyscf.pbc.grad.rhf module#
- class pyscf.pbc.grad.rhf.Gradients(method)[source]#
Bases:
GradientsBaseNon-relativistic Gamma-point restricted Hartree-Fock gradients
- grad_elec(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None, kpt=array([0., 0., 0.]))#
- class pyscf.pbc.grad.rhf.GradientsBase(method)[source]#
Bases:
GradientsBaseBase class for Gamma-point nuclear gradient
pyscf.pbc.grad.rks module#
pyscf.pbc.grad.rks_stress module#
The energy derivatives for the strain tensor e_ij is
1 d E
- sigma_ij = — ——
V d e_ij
The strain tesnor e_ij describes the transformation for real space coordinates in the crystal
sum_j [deta_ij + e_ij] R_j [for j = x, y, z]
Due to numerical errors, the strain tensor may slightly break the symmetry within the stress tensor. The 6 independent components of the stress tensor
[e1 e6/2 e5/2] [e6/2 e2 e4/2] [e5/2 e4/2 e3 ]
is constructed by symmetrizing the strain tensor as follows:
e1 = e_11 e2 = e_22 e3 = e_33 e6 = e_12 + e_21 e5 = e_13 + e_31 e4 = e_32 + e_23
See K. Doll, Mol Phys (2010), 108, 223
- pyscf.pbc.grad.rks_stress.get_vxc(ks_grad, cell, dm, with_j=False, with_nuc=False)[source]#
Strain derivatives for Coulomb and XC at gamma point
- Kwargs:
with_j : Whether to include the electron-electron Coulomb interactions with_nuc : Whether to include the electron-nuclear Coulomb interactions
- pyscf.pbc.grad.rks_stress.kernel(mf_grad)[source]#
Compute the energy derivatives for strain tensor (e_ij)
1 d E
- sigma_ij = — ——
V d e_ij
sigma is a asymmetric 3x3 matrix. The symmetric stress tensor in the 6 Voigt notation can be transformed from the asymmetric stress tensor
sigma1 = sigma_11 sigma2 = sigma_22 sigma3 = sigma_33 sigma6 = (sigma_12 + sigma_21)/2 sigma5 = (sigma_13 + sigma_31)/2 sigma4 = (sigma_23 + sigma_32)/2
See K. Doll, Mol Phys (2010), 108, 223
pyscf.pbc.grad.uhf module#
- class pyscf.pbc.grad.uhf.Gradients(method)[source]#
Bases:
GradientsBaseNon-relativistic Gamma-point restricted Hartree-Fock gradients
- grad_elec(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None, kpt=array([0., 0., 0.]))#
pyscf.pbc.grad.uks module#
pyscf.pbc.grad.uks_stress module#
- pyscf.pbc.grad.uks_stress.get_vxc(ks_grad, cell, dm, with_j=False, with_nuc=False)[source]#
Strain derivatives for Coulomb and XC at gamma point
- Kwargs:
with_j : Whether to include the electron-electron Coulomb interactions with_nuc : Whether to include the electron-nuclear Coulomb interactions
- pyscf.pbc.grad.uks_stress.kernel(mf_grad)[source]#
Compute the energy derivatives for strain tensor (e_ij)
1 d E
- sigma_ij = — ——
V d e_ij
sigma is a asymmetric 3x3 matrix. The symmetric stress tensor in the 6 Voigt notation can be transformed from the asymmetric stress tensor
sigma1 = sigma_11 sigma2 = sigma_22 sigma3 = sigma_33 sigma6 = (sigma_12 + sigma_21)/2 sigma5 = (sigma_13 + sigma_31)/2 sigma4 = (sigma_23 + sigma_32)/2
See K. Doll, Mol Phys (2010), 108, 223
Module contents#
Analytical nuclear gradients for PBC