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)

Bases: pyscf.pbc.grad.krhf.GradientsBasics

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)

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

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)
class pyscf.pbc.grad.krhf.GradientsBasics(method)

Bases: pyscf.grad.rhf.GradientsBasics

Basic nuclear gradient functions for non-relativistic methods

get_hcore(cell=None, kpts=None)
get_j(dm=None, kpts=None)
get_jk(dm=None, kpts=None)

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

get_k(dm=None, kpts=None, kpts_band=None)
get_ovlp(cell=None, kpts=None)
grad_nuc(cell=None, atmlst=None)
hcore_generator(cell=None, kpts=None)
pyscf.pbc.grad.krhf.as_scanner(mf_grad)

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)

Part of the nuclear gradients of core Hamiltonian

pyscf.pbc.grad.krhf.get_j(mf_grad, dm, kpts)
pyscf.pbc.grad.krhf.get_jk(mf_grad, dm, kpts)

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

pyscf.pbc.grad.krhf.get_k(mf_grad, dm, kpts)
pyscf.pbc.grad.krhf.get_ovlp(cell, kpts)
pyscf.pbc.grad.krhf.get_veff(mf_grad, dm, kpts)

NR Hartree-Fock Coulomb repulsion

pyscf.pbc.grad.krhf.grad_elec(mf_grad, 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

pyscf.pbc.grad.krhf.grad_nuc(cell, atmlst)

Derivatives of nuclear repulsion energy wrt nuclear coordinates

pyscf.pbc.grad.krhf.hcore_generator(mf, cell=None, kpts=None)
pyscf.pbc.grad.krhf.make_rdm1e(mo_energy, mo_coeff, mo_occ)

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)

Bases: pyscf.pbc.grad.krhf.Gradients

dump_flags(verbose=None)
get_veff(dm=None, kpts=None)
pyscf.pbc.grad.krks.get_veff(ks_grad, dm=None, kpts=None)
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)

pyscf.pbc.grad.kuhf module

Non-relativistic analytical nuclear gradients for unrestricted Hartree Fock with kpoints sampling

class pyscf.pbc.grad.kuhf.Gradients(method)

Bases: pyscf.pbc.grad.krhf.Gradients

Non-relativistic restricted Hartree-Fock gradients

get_veff(dm=None, kpts=None)
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

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

NR Hartree-Fock Coulomb repulsion

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

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)

Bases: pyscf.pbc.grad.kuhf.Gradients

Non-relativistic restricted Hartree-Fock gradients

get_veff(dm=None, kpts=None)
pyscf.pbc.grad.kuks.get_veff(ks_grad, dm=None, kpts=None)
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)

Module contents

Analytical nuclear gradients for PBC