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: GradientsBase

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

get_veff(dm=None, kpts=None)[source]#
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.).

make_rdm1e(mo_energy=None, mo_coeff=None, mo_occ=None)[source]#
class pyscf.pbc.grad.krhf.GradientsBase(method)[source]#

Bases: GradientsBase

Basic nuclear gradient functions for non-relativistic methods

get_hcore(cell=None, kpts=None)[source]#
get_j(dm=None, kpts=None)[source]#
get_jk(dm=None, kpts=None)[source]#

J = ((-nabla i) j| kl) D_lk K = ((-nabla i) j| kl) D_jk

get_k(dm=None, kpts=None, kpts_band=None)[source]#
get_ovlp(cell=None, kpts=None)[source]#
grad_nuc(cell=None, atmlst=None)[source]#
hcore_generator(cell=None, kpts=None)#
property kpts#
optimizer()[source]#

Geometry optimization solver

Kwargs:

solver (string) : geometry optimization solver, can be “geomeTRIC” (default) or “berny”.

reset(cell=None)[source]#

Clean up intermediates

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_j(mf_grad, dm, kpts)[source]#
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.get_k(mf_grad, dm, kpts)[source]#
pyscf.pbc.grad.krhf.get_ovlp(cell, kpts)[source]#
pyscf.pbc.grad.krhf.get_veff(mf_grad, dm, kpts)[source]#

NR Hartree-Fock Coulomb repulsion

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.krhf.grad_nuc(cell, atmlst)[source]#

Derivatives of nuclear repulsion energy wrt nuclear coordinates

Notes:

An optimized version of this function is available in pbc.gto.ewald_methods.ewald_nuc_grad

pyscf.pbc.grad.krhf.hcore_generator(mf_grad, cell=None, kpts=None)[source]#
pyscf.pbc.grad.krhf.make_rdm1e(mo_energy, mo_coeff, mo_occ)[source]#

Energy weighted density matrix

pyscf.pbc.grad.krks module#

Non-relativistic analytical nuclear gradients for restricted Kohn Sham with kpoints sampling

class pyscf.pbc.grad.krks.Gradients(mf)[source]#

Bases: Gradients

dump_flags(verbose=None)[source]#
get_stress()[source]#
get_veff(dm=None, kpts=None)#
reset(cell=None)[source]#

Clean up intermediates

pyscf.pbc.grad.krks.get_veff(ks_grad, dm=None, kpts=None)[source]#
pyscf.pbc.grad.krks.get_vxc(ni, cell, grids, xc_code, dms, kpts, kpts_band=None, relativity=0, hermi=1, max_memory=2000, verbose=None)[source]#

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_kin(cell, kpts)[source]#

Strain derivatives for kinetic matrix

pyscf.pbc.grad.krks_stress.get_ovlp(cell, kpts)[source]#

Strain derivatives for overlap matrix

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

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.

get_veff(dm=None, kpts=None)[source]#
pyscf.pbc.grad.krkspu.generate_first_order_local_orbitals(cell, minao_ref='MINAO', kpts=None)[source]#

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: GradientsBase

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

get_veff(dm=None, kpts=None)[source]#
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.).

make_rdm1e(mo_energy=None, mo_coeff=None, mo_occ=None)[source]#
pyscf.pbc.grad.kuhf.get_veff(mf_grad, dm, kpts)[source]#

NR Hartree-Fock Coulomb repulsion

pyscf.pbc.grad.kuhf.grad_elec(mf_grad, mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None)[source]#
pyscf.pbc.grad.kuhf.make_rdm1e(mo_energy, mo_coeff, mo_occ)[source]#

Energy weighted density matrix

pyscf.pbc.grad.kuks module#

Non-relativistic analytical nuclear gradients for unrestricted Kohn Sham with kpoints sampling

class pyscf.pbc.grad.kuks.Gradients(mf)[source]#

Bases: Gradients

Non-relativistic restricted Hartree-Fock gradients

get_stress()[source]#
get_veff(dm=None, kpts=None)#
pyscf.pbc.grad.kuks.get_veff(ks_grad, dm=None, kpts=None)[source]#
pyscf.pbc.grad.kuks.get_vxc(ni, cell, grids, xc_code, dms, kpts, kpts_band=None, relativity=0, hermi=1, max_memory=2000, verbose=None)[source]#

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

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.

get_veff(dm=None, kpts=None)[source]#

pyscf.pbc.grad.rhf module#

class pyscf.pbc.grad.rhf.Gradients(method)[source]#

Bases: GradientsBase

Non-relativistic Gamma-point restricted Hartree-Fock gradients

get_veff(mol=None, dm=None, kpt=array([0., 0., 0.]))[source]#
grad_elec(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None, kpt=array([0., 0., 0.]))#
make_rdm1e(mo_energy=None, mo_coeff=None, mo_occ=None)[source]#
class pyscf.pbc.grad.rhf.GradientsBase(method)[source]#

Bases: GradientsBase

Base class for Gamma-point nuclear gradient

get_ovlp(mol=None, kpt=array([0., 0., 0.]))[source]#
grad_nuc(mol=None, atmlst=None)[source]#
optimizer()[source]#

Geometry optimization solver

Kwargs:

solver (string) : geometry optimization solver, can be “geomeTRIC” (default) or “berny”.

pyscf.pbc.grad.rhf.get_ovlp(cell, kpt=array([0., 0., 0.]))[source]#
pyscf.pbc.grad.rhf.get_veff(mf_grad, mol, dm, kpt=array([0., 0., 0.]))[source]#
pyscf.pbc.grad.rhf.grad_elec(mf_grad, mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None, kpt=array([0., 0., 0.]))[source]#
pyscf.pbc.grad.rhf.grad_nuc(cell, atmlst=None, ew_eta=None, ew_cut=None)[source]#

pyscf.pbc.grad.rks module#

class pyscf.pbc.grad.rks.Gradients(method)[source]#

Bases: Gradients

Non-relativistic Gamma-point restricted Kohn-Sham DFT gradients

get_stress()[source]#
grids = None#

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.ewald(cell)[source]#
pyscf.pbc.grad.rks_stress.get_kin(cell)[source]#

Strain derivatives for kinetic matrix

pyscf.pbc.grad.rks_stress.get_ovlp(cell)[source]#

Strain derivatives for overlap matrix

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.rks_stress.strain_tensor_dispalcement(x, y, disp)[source]#

pyscf.pbc.grad.uhf module#

class pyscf.pbc.grad.uhf.Gradients(method)[source]#

Bases: GradientsBase

Non-relativistic Gamma-point restricted Hartree-Fock gradients

get_veff(mol=None, dm=None, kpt=array([0., 0., 0.]))[source]#
grad_elec(mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None, kpt=array([0., 0., 0.]))#
make_rdm1e(mo_energy=None, mo_coeff=None, mo_occ=None)[source]#
pyscf.pbc.grad.uhf.get_veff(mf_grad, mol, dm, kpt=array([0., 0., 0.]))[source]#
pyscf.pbc.grad.uhf.grad_elec(mf_grad, mo_energy=None, mo_coeff=None, mo_occ=None, atmlst=None, kpt=array([0., 0., 0.]))[source]#

pyscf.pbc.grad.uks module#

class pyscf.pbc.grad.uks.Gradients(method)[source]#

Bases: Gradients

Non-relativistic Gamma-point unrestricted Kohn-Sham DFT gradients

get_stress()[source]#
grids = None#

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