15.7. Jacobian
15.7.1. Overview
This module provides the implementation of the Jacobian matrix used for computing singlet excitation energies based on the pair Coupled Cluster Doubles (pCCD) reference wavefunction.
Unlike full EOM solvers, these Jacobian classes construct and diagonalize the true LR Jacobian matrix by removing ground-state coupling terms (\(\langle 0|H|X \rangle\) and \(\langle X|H|0 \rangle\)). This gives excitation energies directly from the similarity-transformed Hamiltonian.
These modules support both iterative (Davidson) and direct (full) diagonalization methods. If you use this module, please cite [ahmadkhani2024].
15.7.2. Implemented Classes
JacobianpCCD
: Jacobian based on pCCD referenceJacobianpCCDS
: Jacobian based on pCCD with singles
Each class is derived from its corresponding EOM base class and reuses most of the functionality, modifying only the construction of the Hamiltonian to match the actual LR structure.
15.7.3. Key Features
Excitation energies from a pCCD reference
Davidson and exact diagonalization supported
Efficient handling of single and pair excitations
Compatible with PyBEST’s modular tensor contraction engine
15.7.4. Methodology
These solvers are based on the standard LR equation:
Here, \(\mathbf{J}\) is the Jacobian matrix formed from the similarity-transformed Hamiltonian, and \(\omega_k\) are excitation energies. The ground state is excluded from the eigenproblem by zeroing the first row and column of \(\mathbf{J}\).
15.7.5. Quick Guide: Jacobian for pCCD and pCCD+S
This is a short introduction to the Jacobian module. More information on the input and output structure can be found in the following sections. Similar to the previous modules, we will assume the following names for all PyBEST objects
- lf
a
LinalgFactory
instance (see Preliminaries).- occ_model
an Aufbau occupation model of the
AufbauOccModel
class- kin
the kinetic energy integrals
- ne
the nucleus-electron attraction integrals
- eri
the two-electron repulsion integrals
The current version of PyBEST supports two different Jacobian flavors, both based on a
pCCD reference function (The pCCD module).
The JacobianpCCD
class allows you to determine
the Jacobian for pair-excitations only. A code snippet can be found below.
The output of a previous pCCD calculation is stored in the pccd_output
container.
from pybest.ee_jacobian import JacobianpCCD
# For pCCD Jacobian (only electron pairs)
jac = JacobianpCCD(lf, occ_model)
result = jac(kin, ne, eri, pccd_output)
# For full diagonalization instead of Davidson:
result = jac(kin, ne, eri, pccd_output, davidson=False)
To include single excitations in the Jacobian, you can use the corresponding
JacobianpCCDS
class.
from pybest.ee_jacobian import JacobianpCCDS
# For pCCD+S Jacobian (singles and electron pairs)
jac = JacobianpCCDS(lf, occ_model)
result = jac(kin, ne, eri, pccd_output)
# For full diagonalization instead of Davidson:
result = jac(kin, ne, eri, pccd_output, davidson=False)
from pybest.ee_jacobian import JacobianpCCD
from pybest.ee_jacobian import JacobianpCCDS
# For pCCD Jacobian (no singles)
jac = JacobianpCCD(lf, occ_model)
result = jac(kin, ne, eri, pccd_output)
# For pCCDS Jacobian (includes singles)
jac = JacobianpCCDS(lf, occ_model)
result = jac(kin, ne, eri, pccds_output)
# For full diagonalization instead of Davidson:
result = jac(kin, ne, eri, pccd_output, davidson=False)
The result includes excitation energies and right eigenvectors of the Jacobian matrix.
The results are returned as an IOData
container,
while all results are saved to the pybest-results/checkpoint_LRpCCD.h5
file.
Specifically, the IOData
container contains
the following attributes
- e_ref
Total reference energy (pCCD)
- e_ee
Excitation energies (first element is 0.0 by convention)
- civ_ee
Eigenvectors of the Jacobian matrix
- orb_a
Molecular orbitals
15.7.6. Keyword Arguments
These classes support the following keyword arguments:
- nroot
(int) Number of excited states (default: Hamiltonian dimension - 1 (ground state))
- davidson
(bool) Use Davidson diagonalization (default: True)
- tolerance
(float) Convergence threshold (default: 1e-6)
- maxiter
(int) Maximum Davidson iterations (default: 200)
Note
For Jacobian-based LR methods, keep the default value for nroot
(equal to dimension - 1
).
This is required to construct the full Jacobian matrix from all eigenvectors and eigenvalues.
Changing nroot
will lead to a dimension error.
15.7.7. Limitations
We do not yet support open-shell electronic structures. The accuracy of the implemented LR models strongly depends on the quality of the pCCD reference function and the +S excited state extension.