CAPHamiltonian
This section briefly describes how to use the CAPHamiltonian object to generate eigenvalue trajectories.
Initialization
CAPHamiltonian
objects can be initialized in one of two ways.
The first is to pass H0 and W as numpy arrays:
from pyopencap.analysis.CAPTrajectory import CAPHamiltonian
CAPH = CAPHamiltonian(H0=h0,W=mat)
The other is to read them in from an OpenCAP output file, or from a Q-Chem output file generated by a Projected CAP-EOM-CC or Projected CAP-ADC calculation.
CAPH = CAPHamiltonian(output="path/to/output.out")
If one wishes to exclude some of the states from the analysis, this can be accomplished through the by placing their indices in a list (starting from 0) and passing it into the exclude_states keyword argument:
exclude_states = [2,5,7]
CAPH = CAPHamiltonian(H0=h0,W=mat,exclude_states=exclude_states)
Similarly, the include_states argument includes only the desired states. Note that these two keywords are incompatible.
include_states = [0,1,2,3,4]
CAPH = CAPHamiltonian(H0=h0,W=mat,include_states=include_states)
Importantly, in all cases, the W matrix is assumed to be pre-multiplied by a factor of -1.0.
Diagonalization
The run_trajectory()
function diagonalizes the
CAP Hamiltonian over a range of eta values (and at a specified value of the cap lambda parameter
if using a CR-CAP).
eta_list = np.linspace(0,2000,101)
eta_list = eta_list * 1E-5
CAPH.run_trajectory(eta_list,cap_lambda=0.0)
Since the W matrix is assumed to be multiplied by a factor of -1.0 upon instantiation, the following matrix is actually diagonalized at each step:
\(H^{CAP}=H_0+(i\eta- \lambda) W\)
and each eigenpair is stored in a Root
object.
After all of the diagonalizations are finished, individual states can be tracked
using the track_state()
function:
traj = CAPH.track_state(4,tracking="overlap")
The traj variable is a EigenvalueTrajectory
object, which
contains helpful functions for analysis. Indices for states start from 0, and there are
two options for tracking states: “overlap” (the default), and “energy”. See EigenvalueTrajectory
for more details.
Visualization
The energies of all states computed are stored in the total_energies instance variable of the
CAPHamiltonian
object. This can
very useful for graphical searches e.g.
import matplotlib.pyplot as plt
import numpy as np
plt.plot(np.real(CAPH.total_energies),np.imag(CAPH.total_energies),'ro')
plt.show()
There is also a function energies_ev()
which returns
the excitation energies in eV with respect to specified reference energy.
E_ev = CAPH.energies_ev(ref_energy)
plt.plot(np.real(E_ev),np.imag(CAPH.energies_ev(E_ev),'ro')
plt.show()
- class pyopencap.analysis.CAPHamiltonian(pc=None, H0=None, W=None, output=None, irrep='', onset='')[source]
Projected CAP Hamiltonian handler for generating eigenvalue trajectories.
The instance variables H0,W etc. are only set after run_trajectory is executed. The original matrices passed/parsed when the object is constructed are stored in _H0, _W, etc. This makes it easy to run multiple trajectories with different states included in the projection scheme without having to construct a new object.
- H0
Zeroth order Hamiltonian in state basis
- Type
np.ndarray of float: default=None
- W
CAP matrix in state basis. -1 prefactor is assumed.
- Type
np.ndarray of float: default=None
- nstates
Number of states
- Type
int
- total_energies
Energies of all states found by repeated diagonalization of CAP Hamiltonian
- Type
list of float
- etas
List of CAP strengths in trajectory.
- Type
list of float
- cap_lambda
Real CAP strength used for continuum remover CAP. Set to 0.0 by default.
- Type
float
- __init__(pc=None, H0=None, W=None, output=None, irrep='', onset='')[source]
Initializes CAPHamiltonian object from H0 and W matrix in state basis.
Object can be initialized in one of two ways. The user can pass H0 and W directly as numpy matrices, or they can specify a path to a properly formatted output file (either OpenCAP or Q-Chem) which contains these two matrices.
W matrix is assumed to already have a -1 prefactor, as that is how OpenCAP output is formatted.
- Parameters
pc (
CAP
: default=None) – PyOpenCAP CAP objectH0 (np.ndarray of float: default=None) – Zeroth order Hamiltonian in state basis
W (np.ndarray of float: default=None) – CAP matrix in state basis. -1 prefactor is assumed.
output (str: default=None) – Path to Q-Chem or OpenCAP output file.
irrep (str: default=None) – Title of irreducible representation of state of interest. Only compatible with Q-Chem projected CAP-EOM-CC outputs. Set to ‘all’ to include all symmetries in CAP projection.
onset (str: default=None) – Title of CAP onset. Only compatible with Q-Chem projected CAP-ADC outputs.
- energies_ev(ref_energy)[source]
Returns excitation energies of all calculated states in eV with respect to specified reference energy.
- Parameters
ref_energy (float) – Reference energy
- Returns
E_eV – Excitation energies in eV with respect to specified reference energy.
- Return type
list of floats
- export(finame)[source]
Exports Zeroth order Hamiltonian and CAP matrix to an OpenCAP formatted output file for further analysis. Useful for saving the results of an expensive electronic structure calculation performed in a python environment.
- Parameters
finame (str) – File handle to export data.
- run_trajectory(eta_list, cap_lambda=0.0, exclude_states=None, include_states=None, biorthogonalize=False)[source]
Diagonalizes CAP Hamiltonian over range of eta values.
CAP Hamiltonian is defined as H^CAP = H0 + i*eta*W - cap_lambda * W. W matrix is assumed to already have a -1 prefactor, as that is how OpenCAP output is formatted. Recommended range for eta_list is between 1E-5 and 1E-2, though this can vary widely based on system and CAP shape.
- Parameters
eta_list (iterable object) – List of eta values to use
cap_lambda (float, default=0.0) – Real CAP strength to use for continuum remover CAP
exclude_states (list of int, default=None) – List of states to exclude from subspace projection. Not compatible with include_states parameter.
include_states (list of int, default=None) – List of states to include in subspace projection. Not compatible with exclude_states parameter.
biorthogonalize (bool, default=False) – Biorthogonalize left and right eigenvectors. If false, left eigenvectors are assumed to equal right eigenvectors for density matrix correction
- track_state(state_idx, tracking='overlap', correction='derivative')[source]
Tracks eigenvalue trajectory over range of eta values.
- Parameters
state_idx (int) – Index of state to track
tracking (str, default="overlap") – Method to use to track the state. Options are “overlap”, which tracks based on eigenvector overlap, and “energy” which tracks based on energy.
correction (str, default="derivative") – Choice of correction scheme. Either “density” or “derivative”.
- Returns
traj – Eigenvalue trajectory for further analysis.
- Return type