CheMPS2
CheMPS2::CASSCF Class Reference

#include <CASSCF.h>

Public Member Functions

 CASSCF (Hamiltonian *ham_in, int *docc, int *socc, int *nocc, int *ndmrg, int *nvirt, const string tmp_folder=CheMPS2::defaultTMPpath)
 Constructor. More...
 
virtual ~CASSCF ()
 Destructor.
 
int get_num_irreps ()
 Get the number of irreps. More...
 
double solve (const int Nelectrons, const int TwoS, const int Irrep, ConvergenceScheme *OptScheme, const int rootNum, DMRGSCFoptions *scf_options)
 Do the CASSCF cycles with the augmented Hessian Newton-Raphson method. More...
 
double caspt2 (const int Nelectrons, const int TwoS, const int Irrep, ConvergenceScheme *OptScheme, const int rootNum, DMRGSCFoptions *scf_options, const double IPEA, const double IMAG, const bool PSEUDOCANONICAL, const bool CHECKPOINT=false, const bool CUMULANT=false)
 Calculate the caspt2 correction energy for a converged casscf wavefunction. More...
 

Static Public Member Functions

static void deleteStoredUnitary (const string filename=CheMPS2::DMRGSCF_unitary_storage_name)
 CASSCF unitary rotation remove call.
 
static void deleteStoredDIIS (const string filename=CheMPS2::DMRGSCF_diis_storage_name)
 CASSCF DIIS vectors remove call.
 
static void buildFmat (DMRGSCFmatrix *localFmat, const DMRGSCFmatrix *localTmat, const DMRGSCFmatrix *localJKocc, const DMRGSCFmatrix *localJKact, const DMRGSCFindices *localIdx, const DMRGSCFintegrals *theInts, double *local2DM, double *local1DM)
 Build the F-matrix (Eq. (11) in the Siegbahn paper [CAS3]) More...
 
static void buildWtilde (DMRGSCFwtilde *localwtilde, const DMRGSCFmatrix *localTmat, const DMRGSCFmatrix *localJKocc, const DMRGSCFmatrix *localJKact, const DMRGSCFindices *localIdx, const DMRGSCFintegrals *theInts, double *local2DM, double *local1DM)
 Build the Wtilde-matrix (Eq. (20b) in the Siegbahn paper [CAS3]) More...
 
static void augmentedHessianNR (DMRGSCFmatrix *localFmat, DMRGSCFwtilde *localwtilde, const DMRGSCFindices *localIdx, const DMRGSCFunitary *localUmat, double *theupdate, double *updateNorm, double *gradNorm)
 Calculate the augmented Hessian Newton-Raphson update for the orthogonal orbital rotation matrix. More...
 
static void copy2DMover (TwoDM *theDMRG2DM, const int LAS, double *two_dm)
 Copy over the DMRG 2-RDM. More...
 
static void setDMRG1DM (const int num_elec, const int LAS, double *one_dm, double *two_dm)
 Construct the 1-RDM from the 2-RDM. More...
 
static void copy_active (double *origin, DMRGSCFmatrix *result, const DMRGSCFindices *idx, const bool one_rdm)
 Copy a one-orbital quantity from array format to DMRGSCFmatrix format. More...
 
static void copy_active (const DMRGSCFmatrix *origin, double *result, const DMRGSCFindices *idx)
 Copy a one-orbital quantity from DMRGSCFmatrix format to array format. More...
 
static void fillLocalizedOrbitalRotations (DMRGSCFunitary *umat, DMRGSCFindices *idx, double *eigenvecs)
 From an Edmiston-Ruedenberg active space rotation, fetch the eigenvectors and store them in eigenvecs. More...
 
static void block_diagonalize (const char space, const DMRGSCFmatrix *Mat, DMRGSCFunitary *Umat, double *work1, double *work2, const DMRGSCFindices *idx, const bool invert, double *two_dm, double *three_dm, double *contract)
 Block-diagonalize Mat. More...
 
static void construct_fock (DMRGSCFmatrix *Fock, const DMRGSCFmatrix *Tmat, const DMRGSCFmatrix *Qocc, const DMRGSCFmatrix *Qact, const DMRGSCFindices *idx)
 Construct the Fock matrix. More...
 
static double deviation_from_blockdiag (DMRGSCFmatrix *matrix, const DMRGSCFindices *idx)
 Return the RMS deviation from block-diagonal. More...
 
static void write_f4rdm_checkpoint (const string f4rdm_file, int *hamorb1, int *hamorb2, const int tot_dmrg_power6, double *contract)
 Write the checkpoint file for the contraction of the generalized Fock operator with the 4-RDM to disk. More...
 
static bool read_f4rdm_checkpoint (const string f4rdm_file, int *hamorb1, int *hamorb2, const int tot_dmrg_power6, double *contract)
 Read the checkpoint file for the contraction of the generalized Fock operator with the 4-RDM from disk. More...
 
static void fock_dot_4rdm (double *fockmx, CheMPS2::DMRG *dmrgsolver, CheMPS2::Hamiltonian *ham, int next_orb1, int next_orb2, double *work, double *result, const bool CHECKPOINT, const bool PSEUDOCANONICAL)
 Build the contraction of the fock matrix with the 4-RDM. More...
 

Detailed Description

CASSCF class.

Author
Sebastian Wouters sebas.nosp@m.tian.nosp@m.woute.nosp@m.rs@g.nosp@m.mail..nosp@m.com
Date
June 18, 2013

Introduction

In methods which use a FCI solver, this solver can be replaced by DMRG. DMRG allows for an efficient extraction of the 2-RDM [CAS1, CAS2]. The 2-RDM of the active space is required in the complete active space self-consistent field (CASSCF) method to compute the gradient and the Hessian with respect to orbital rotations [CAS3]. It is therefore natural to introduce a CASSCF variant with DMRG as active space solver, called DMRG-SCF [CAS2, CAS4, CAS5], which allows to treat static correlation in large active spaces. In CheMPS2, I have implemented the augmented Hessian Newton-Raphson DMRG-SCF method, with exact Hessian [CAS5, CAS6]. It can be called with the function CheMPS2::CASSCF::doCASSCFnewtonraphson.

Orbital gradient and Hessian

The calculation of the orbital gradient and Hessian for DMRG-SCF is based on [CAS3]. The basic idea is to express the energy with the unitary group generators:

\begin{eqnarray*} \hat{E}_{pq} & = & \sum\limits_{\sigma} \hat{a}^{\dagger}_{p \sigma} \hat{a}_{q \sigma} \\ \left[ \hat{E}_{pq} , \hat{E}_{rs} \right] & = & \delta_{qr} \hat{E}_{ps} - \delta_{ps} \hat{E}_{rq} \\ \hat{E}^{-}_{pq} & = & \hat{E}_{pq} - \hat{E}_{qp} \\ \hat{T} & = & \sum\limits_{p<q} x_{pq} \hat{E}^{-}_{pq} \\ E(\vec{x}) & = & \braket{0 | e^{\hat{T}} \hat{H} e^{-\hat{T}} | 0 } \\ \left. \frac{\partial E(\vec{x})}{\partial x_{ij}} \right|_{0} & = & \braket{ 0 | \left[ \hat{E}_{ij}^{-}, \hat{H} \right] | 0 } \\ \left. \frac{\partial^2 E(\vec{x})}{\partial x_{ij} \partial x_{kl}} \right|_{0} & = & \frac{1}{2} \braket{ 0 | \left[ \hat{E}_{ij}^{-}, \left[ \hat{E}_{kl}^{-}, \hat{H} \right] \right] | 0 } + \frac{1}{2} \braket{ 0 | \left[ \hat{E}_{kl}^{-}, \left[ \hat{E}_{ij}^{-}, \hat{H} \right] \right] | 0 } \end{eqnarray*}

The variables $x_{pq}$ only connect orbitals with the same irrep ( $I_p=I_q$). Assuming that DMRG is exact, $x_{pq}$ in addition only connects orbitals when they belong to different occupation blocks: occupied, active, virtual. With some algebra, the derivatives can be rewritten. Real-valued symmetric one-electron integrals $h_{ij}$ and real-valued eightfold permutation symmetric two-electron integrals $(ij | kl)$ are assumed (chemical notation for the latter).

\begin{eqnarray*} \Gamma^{2A}_{ijkl} & = & \sum\limits_{\sigma \tau} \braket{ 0 | \hat{a}^{\dagger}_{i \sigma} \hat{a}_{j \tau}^{\dagger} \hat{a}_{l \tau} \hat{a}_{k \sigma} | 0} \\ \Gamma^1_{ij} & = & \sum\limits_{\sigma} \braket{ 0 | \hat{a}^{\dagger}_{i \sigma} \hat{a}_{j \sigma} | 0} \\ \left. \frac{\partial E(\vec{x})}{\partial x_{ij}} \right|_{0} & = & 2 \left( F_{ij} - F_{ji} \right) \\ F_{pq} & = & \sum\limits_{r} \Gamma_{pr}^{1} h_{qr} + \sum\limits_{rst} \Gamma^{2A}_{psrt} (qr | st) \\ \left. \frac{\partial^2 E(\vec{x})}{\partial x_{ij} \partial x_{kl}} \right|_{0} & = & w_{ijkl} - w_{jikl} - w_{ijlk} + w_{jilk} \\ w_{pqrs} & = & \delta_{qr} \left( F_{ps} + F_{sp} \right) + \tilde{w}_{pqrs} \\ \tilde{w}_{pqrs} & = & 2 \Gamma^1_{pr} h_{qs} + 2 \sum\limits_{\alpha \beta} \left( \Gamma^{2A}_{r \alpha p \beta} (qs | \alpha \beta ) + \left( \Gamma^{2A}_{r \alpha \beta p} + \Gamma^{2A}_{r p \beta \alpha} \right) (q \alpha | s \beta ) \right) \end{eqnarray*}

In the calculation of $F_{pq}$, the indices $prst$ can only be occupied or active due to their appearance in the density matrices, and the only index which can be virtual is hence $q$. Moreover, due to the irrep symmetry of the integrals and density matrices, $F_{pq}$ is diagonal in the irreps: $I_p = I_q$. Alternatively, this can be understood by the fact that $x_{pq}$ only connects orbitals with the same irrep.

In the calculation of $\tilde{w}_{pqrs}$, the indices $pr\alpha\beta$ can only be occupied or active due to their appearance in the density matrices, and the only indices which can be virtual are hence $qs$. Together with the remark for $F_{pq}$, this can save time for the two-electron integral rotation. Moreover, as $x_{pq}$ only connects orbitals with the same irrep, $I_p = I_q$ and $I_r = I_s$ in $\tilde{w}_{pqrs}$.

By rewriting the density matrices, the calculation of $F_{pq}$ and $\tilde{w}_{pqrs}$ can be simplified. In the following, occ and act denote the doubly occupied and active orbital spaces, respectively.

\begin{eqnarray*} \Gamma^1_{ij} & = & 2 \delta_{ij}^{occ} + \Gamma^{1,act}_{ij} \\ \Gamma^{2A}_{ijkl} & = & \Gamma^{2A,act}_{ijkl} + 2 \delta_{ik}^{occ} \Gamma^{1,act}_{jl} - \delta_{il}^{occ} \Gamma^{1,act}_{jk} \\ & + & 2 \delta_{jl}^{occ} \Gamma^{1,act}_{ik} - \delta_{jk}^{occ} \Gamma^{1,act}_{il} + 4 \delta_{ik}^{occ} \delta_{jl}^{occ} - 2 \delta_{il}^{occ} \delta_{jk}^{occ} \end{eqnarray*}

Define the following symmetric charge (Coulomb + exchange) matrices:

\begin{eqnarray*} Q^{occ}_{ij} & = & \sum\limits_{s \in occ} \left[ 2 (ij | ss) - (is | js ) \right] \\ Q^{act}_{ij} & = & \sum\limits_{st \in act} \frac{1}{2} \Gamma^{1,act}_{st} \left[ 2 (ij | st) - (is | jt ) \right] \end{eqnarray*}

They can be calculated efficiently by (1) rotating the occupied and active density matrices from the current basis to the original basis, (2) contracting the rotated density matrices with the two-electron integrals in the original basis, and (3) rotating these contractions to the current basis. The constant part and the one-electron integrals of the active space Hamiltonian are:

\begin{eqnarray*} \tilde{E}_{const} & = & E_{const} + \sum\limits_{s \in occ} \left( 2 h_{ss} + Q_{ss}^{occ} \right) \\ \tilde{h}_{ij} & = & h_{ij} + Q_{ij}^{occ} \end{eqnarray*}

The calculation of $F_{pq}$ boils down to:

\begin{eqnarray*} p \in occ & : & F_{pq} = 2 \left( h_{qp} + Q^{occ}_{qp} + Q^{act}_{qp} \right) \\ p \in act & : & F_{pq} = \sum\limits_{r \in act} \Gamma^{1,act}_{pr} \left[ h_{qr} + Q^{occ}_{qr} \right] + \sum\limits_{rst \in act} \Gamma^{2A,act}_{psrt} (qr | st) \end{eqnarray*}

And the calculation of $\tilde{w}_{pqrs}$ (remember that $I_p = I_q$ and $I_r = I_s$):

\begin{eqnarray*} (p,r) & \in & (occ,occ) : \tilde{w}_{pqrs} \\ & = & 4 \delta_{pr}^{occ} \left[ h_{qs} + Q^{occ}_{qs} + Q^{act}_{qs} \right] \\ & + & 4 \left[ 4 (qp | sr) - ( qs | pr ) - ( qr | sp ) \right] \\ (p,r) & \in & (act,act) : \tilde{w}_{pqrs} = 2 \Gamma^{1,act}_{rp} \left[ h_{qs} + Q^{occ}_{qs} \right] \\ & + & 2 \sum\limits_{\alpha\beta \in act} \left[ \Gamma^{2A,act}_{r \alpha p \beta} (qs | \alpha \beta ) + \left( \Gamma^{2A,act}_{r \alpha \beta p} + \Gamma^{2A,act}_{r p \beta \alpha} \right) (q \alpha | s \beta ) \right] \\ (p,r) & \in & (act,occ) : \tilde{w}_{pqrs} \\ & = & 2 \sum\limits_{\alpha \in act} \Gamma^{1,act}_{\alpha p} \left[ 4 (q \alpha | s r) - (qs | \alpha r) - (qr | s \alpha) \right] \\ (p,r) & \in & (occ,act) : \tilde{w}_{pqrs} \\ & = & 2 \sum\limits_{\beta \in act} \Gamma^{1,act}_{r \beta} \left[ 4 (q p | s \beta) - (qs | p \beta) - (q \beta | sp) \right] \end{eqnarray*}

Augmented Hessian Newton-Raphson DMRG-SCF

The CASSCF energy is a function of $\vec{x}$. Up to second order, the energy is given by

\[ E(\vec{x}) = \braket{0 | e^{\hat{T}(\vec{x})} \hat{H} e^{-\hat{T}(\vec{x})} | 0 } \approx E(0) + \vec{x}^T \vec{g} + \frac{1}{2} \vec{x}^T \mathbf{H} \vec{x}. \]

The vector $\vec{g}$ is the gradient and the matrix $\mathbf{H}$ the Hessian for orbital rotations [CAS3]. They have been described in the previous section. The minimum of $E(\vec{x})$ is found at $\vec{x} = - \mathbf{H}^{-1} \vec{g}$. The variables $\vec{x}$ parametrize an additional orbital rotation $\mathbf{U}_{add} = \exp(\mathbf{X}(\vec{x}))$, with $\mathbf{X}(\vec{x}) = -\mathbf{X}^T(\vec{x})$ a real-valued skew-symmetric matrix. The additional orbital rotation $\mathbf{U}_{add}$ transforms the current orbitals $\mathbf{U}(n)$ to the new orbitals

\[ \mathbf{U}(n+1) = \mathbf{U}_{add} \mathbf{U}(n) = \exp(\mathbf{X}(\vec{x}(n))) \mathbf{U}(n). \]

This updating scheme is called the Newton-Raphson method [CAS3]. If the Hessian is positive definite, these updates are stable. For saddle points in the energy landscape, the Hessian has negative eigenvalues, and these updates can be unstable. It is therefore better to use the augmented Hessian Newton-Raphson method [CAS7]:

\[ \left[ \begin{array}{cc} \mathbf{H} & \vec{g} \\ \vec{g}^T & 0 \end{array} \right] \left[ \begin{array}{c} \vec{x} \\ 1 \end{array} \right] = \alpha \left[ \begin{array}{c} \vec{x} \\ 1 \end{array} \right]. \]

The eigenvector with smallest algebraic eigenvalue determines a stable update $\vec{x}$, as is well explained in Ref. [CAS7].

As a final remark in this section, I would like to say that orbitals have gauge freedom. One can always multiply them with a phase factor. It is therefore possible to choose the orbital gauges so that all $\mathbf{U}$ are always special orthogonal: $\det(\mathbf{U})=+1$.

Direct inversion of the iterative subspace (DIIS)

When the update norm $\|\vec{x}\|_2$ is small enough, the convergence can be accelerated by the direct inversion of the iterative subspace (DIIS) [CAS5, CAS8, CAS9, CAS10]. For a given set of orbitals $\mathbf{U}(n)$, the update $\vec{x}(n)$ is calculated with the augmented Hessian Newton-Raphson method. This update defines the next set of orbitals:

\[ \mathbf{U}(n+1) = \mathbf{U}_{add} \mathbf{U}(n) = \exp(\mathbf{X}(\vec{x}(n))) \mathbf{U}(n). \]

In DIIS, the error vector $\vec{x}(n)$ and the state vector $\mathbf{Y}(n+1) = \log(\mathbf{U}(n+1))$ are added to a list. The error

\[ e = \left\| \sum\limits_{i = 1}^{\kappa} c_i \vec{x}(n - \kappa + i) \right\|_2 \]

is minimized under the constraint $\sum_{i} c_i = 1$. $\kappa$ is the size of the list memory, i.e. the number of retained vectors. The minimization of the error $e$ can be performed with Lagrangian calculus:

\[ \left[ \begin{array}{cc} \mathbf{B} & \vec{1} \\ \vec{1}^T & 0 \end{array} \right] \left[ \begin{array}{c} \vec{c} \\ \lambda \end{array} \right] = \left[ \begin{array}{c} \vec{0} \\ 1 \end{array} \right], \]

where $2\lambda$ is the Lagrangian multiplier and

\[ \left[\mathbf{B}\right]_{ij} = \vec{x}^T(n - \kappa + i) \vec{x}(n - \kappa + j) = \left[\mathbf{B}\right]_{ji}. \]

The new state vector is then defined as

\[ \mathbf{Y}_{new} = \sum\limits_{i = 1}^{\kappa} c_i \mathbf{Y}(n+1 - \kappa + i). \]

The new state vector $\mathbf{Y}_{new}$ is calculated by the function CheMPS2::DIIS::calculateParam. The current orbitals are then set to $\mathbf{U}(n+1) = \exp(\mathbf{Y}_{new})$.

References

[CAS1] D. Zgid and M. Nooijen, Journal of Chemical Physics 128, 144115 (2008). http://dx.doi.org/10.1063/1.2883980
[CAS2] D. Ghosh, J. Hachmann, T. Yanai and G.K.-L. Chan, Journal of Chemical Physics 128, 144117 (2008). http://dx.doi.org/10.1063/1.2883976
[CAS3] P.E.M. Siegbahn, J. Almlof, A. Heiberg and B.O. Roos, Journal of Chemical Physics 74, 2384-2396 (1981). http://dx.doi.org/10.1063/1.441359
[CAS4] D. Zgid and M. Nooijen, Journal of Chemical Physics 128, 144116 (2008). http://dx.doi.org/10.1063/1.2883981
[CAS5] T. Yanai, Y. Kurashige, D. Ghosh and G.K.-L. Chan, International Journal of Quantum Chemistry 109, 2178-2190 (2009). http://dx.doi.org/10.1002/qua.22099
[CAS6] S. Wouters, W. Poelmans, P.W. Ayers and D. Van Neck, Computer Physics Communications 185, 1501-1514 (2014). http://dx.doi.org/10.1016/j.cpc.2014.01.019
[CAS7] A. Banerjee, N. Adams, J. Simons and R. Shepard, Journal of Physical Chemistry 89, 52-57 (1985). http://dx.doi.org/10.1021/j100247a015
[CAS8] P. Pulay, Chemical Physics Letters 73, 393-398 (1980). http://dx.doi.org/10.1016/0009-2614(80)80396-4
[CAS9] C.D. Sherrill, Programming DIIS, http://vergil.chemistry.gatech.edu/notes/diis/node3.html (2000).
[CAS10] T. Rohwedder and R. Schneider, Journal of Mathematical Chemistry 49, 1889-1914 (2011). http://dx.doi.org/10.1007/s10910-011-9863-y

Definition at line 164 of file CASSCF.h.

Constructor & Destructor Documentation

CheMPS2::CASSCF::CASSCF ( Hamiltonian ham_in,
int *  docc,
int *  socc,
int *  nocc,
int *  ndmrg,
int *  nvirt,
const string  tmp_folder = CheMPS2::defaultTMPpath 
)

Constructor.

Parameters
ham_inHamiltonian containing the matrix elements of the Hamiltonian for which a CASSCF calculation is desired
doccArray containing the number of doubly occupied HF orbitals per irrep
soccArray containing the number of singly occupied HF orbitals per irrep
noccArray containing the number of doubly occupied (inactive) orbitals per irrep
ndmrgArray containing the number of active orbitals per irrep
nvirtArray containing the number of virtual (secondary) orbitals per irrep
tmp_folderTemporary work folder for the DMRG renormalized operators and the ERI rotations

Definition at line 39 of file CASSCF.cpp.

+ Here is the call graph for this function:

Member Function Documentation

void CheMPS2::CASSCF::augmentedHessianNR ( DMRGSCFmatrix localFmat,
DMRGSCFwtilde localwtilde,
const DMRGSCFindices localIdx,
const DMRGSCFunitary localUmat,
double *  theupdate,
double *  updateNorm,
double *  gradNorm 
)
static

Calculate the augmented Hessian Newton-Raphson update for the orthogonal orbital rotation matrix.

Parameters
localFmatMatrix which contains the Fock operator (Eq. (11) in the Siegbahn paper [CAS3])
localwtildeObject which contains the second order derivative of the energy with respect to the unitary (Eq. (20b) in the Siegbahn paper [CAS3])
localIdxOrbital index bookkeeper for the CASSCF calculations
localUmatThe unitary matrix for CASSCF calculations (in this function it is used to fetch the orbital ordering convention of the skew-symmetric parametrization)
theupdateWhere the augmented Hessian Newton-Raphson update will be stored
updateNormPointer to one double to store the update norm
gradNormPointer to one double to store the gradient norm

Definition at line 316 of file CASSCFnewtonraphson.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::CASSCF::block_diagonalize ( const char  space,
const DMRGSCFmatrix Mat,
DMRGSCFunitary Umat,
double *  work1,
double *  work2,
const DMRGSCFindices idx,
const bool  invert,
double *  two_dm,
double *  three_dm,
double *  contract 
)
static

Block-diagonalize Mat.

Parameters
spaceCan be 'O', 'A', or 'V' and denotes which block of Mat should be considered
MatMatrix to block-diagonalize
UmatThe unitary rotation will be updated so that Mat is block-diagonal in the orbitals 'space'
work1Workspace
work2Workspace
idxObject which handles the index conventions for CASSCF
invertIf true, the eigenvectors are sorted from large to small instead of the other way around
two_dmIf not NULL, this 4-index array will be rotated to the new eigenvecs if space == 'A'
three_dmIf not NULL, this 6-index array will be rotated to the new eigenvecs if space == 'A'
contractIf not NULL, this 6-index array will be rotated to the new eigenvecs if space == 'A'

Definition at line 433 of file CASSCF.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::CASSCF::buildFmat ( DMRGSCFmatrix localFmat,
const DMRGSCFmatrix localTmat,
const DMRGSCFmatrix localJKocc,
const DMRGSCFmatrix localJKact,
const DMRGSCFindices localIdx,
const DMRGSCFintegrals theInts,
double *  local2DM,
double *  local1DM 
)
static

Build the F-matrix (Eq. (11) in the Siegbahn paper [CAS3])

Parameters
localFmatMatrix where the result should be stored
localTmatMatrix which contains the one-electron integrals
localJKoccMatrix which contains the Coulomb and exchange interaction due to the frozen core orbitals
localJKactMatrix which contains the Coulomb and exchange interaction due to the active space
localIdxOrbital index bookkeeper for the CASSCF calculations
theIntsThe rotated two-electron integrals (at most 2 virtual indices)
local2DMThe DMRG 2-RDM
local1DMThe DMRG 1-RDM

Definition at line 823 of file CASSCFnewtonraphson.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::CASSCF::buildWtilde ( DMRGSCFwtilde localwtilde,
const DMRGSCFmatrix localTmat,
const DMRGSCFmatrix localJKocc,
const DMRGSCFmatrix localJKact,
const DMRGSCFindices localIdx,
const DMRGSCFintegrals theInts,
double *  local2DM,
double *  local1DM 
)
static

Build the Wtilde-matrix (Eq. (20b) in the Siegbahn paper [CAS3])

Parameters
localwtildeWhere the result should be stored
localTmatMatrix which contains the one-electron integrals
localJKoccMatrix which contains the Coulomb and exchange interaction due to the frozen core orbitals
localJKactMatrix which contains the Coulomb and exchange interaction due to the active space
localIdxOrbital index bookkeeper for the CASSCF calculations
theIntsThe rotated two-electron integrals (at most 2 virtual indices)
local2DMThe DMRG 2-RDM
local1DMThe DMRG 1-RDM

Definition at line 611 of file CASSCFnewtonraphson.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

double CheMPS2::CASSCF::caspt2 ( const int  Nelectrons,
const int  TwoS,
const int  Irrep,
ConvergenceScheme OptScheme,
const int  rootNum,
DMRGSCFoptions scf_options,
const double  IPEA,
const double  IMAG,
const bool  PSEUDOCANONICAL,
const bool  CHECKPOINT = false,
const bool  CUMULANT = false 
)

Calculate the caspt2 correction energy for a converged casscf wavefunction.

Parameters
NelectronsTotal number of electrons in the system: occupied HF orbitals + active space
TwoSTwice the targeted spin
IrrepDesired wave-function irrep
OptSchemeThe optimization scheme to run the inner DMRG loop. If NULL: use FCI instead of DMRG.
rootNumDenotes the targeted state in state-specific CASSCF; 1 means ground state, 2 first excited state etc.
scf_optionsContains the DMRGSCF options
IPEAThe CASPT2 IPEA shift from Ghigo, Roos and Malmqvist, Chemical Physics Letters 396, 142-149 (2004)
IMAGThe CASPT2 imaginary shift from Forsberg and Malmqvist, Chemical Physics Letters 274, 196-204 (1997)
PSEUDOCANONICALIf true, use the exact DMRG 4-RDM in the pseudocanonical basis. If false, use the cumulant approximated DMRG 4-RDM in the unrotated basis.
CHECKPOINTIf true, write checkpoints to disk and read them back in again in order to perform the contraction of the generalized Fock operator with the 4-RDM in multiple runs.
CUMULANTIf true, a cumulant approximation is used for the 4-RDM and CHECKPOINT is overwritten to false. If false, the full 4-RDM is used.
Returns
The CASPT2 variational correction energy

Definition at line 189 of file CASSCFpt2.cpp.

+ Here is the call graph for this function:

void CheMPS2::CASSCF::construct_fock ( DMRGSCFmatrix Fock,
const DMRGSCFmatrix Tmat,
const DMRGSCFmatrix Qocc,
const DMRGSCFmatrix Qact,
const DMRGSCFindices idx 
)
static

Construct the Fock matrix.

Parameters
FockMatrix to store the Fock operator in
TmatMatrix with the one-electron integrals
QoccMatrix with the Coulomb and exchange contributions of the occupied (inactive) orbitals
QactMatrix with the Coulomb and exchange contributions of the active space orbitals
idxObject which handles the index conventions for CASSCF

Definition at line 395 of file CASSCFpt2.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::CASSCF::copy2DMover ( TwoDM theDMRG2DM,
const int  LAS,
double *  two_dm 
)
static

Copy over the DMRG 2-RDM.

Parameters
theDMRG2DMThe 2-RDM from the DMRG object
LASThe total number of DMRG orbitals
two_dmThe CASSCF 2-RDM

Definition at line 119 of file CASSCF.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::CASSCF::copy_active ( double *  origin,
DMRGSCFmatrix result,
const DMRGSCFindices idx,
const bool  one_rdm 
)
static

Copy a one-orbital quantity from array format to DMRGSCFmatrix format.

Parameters
originArray to copy
resultDMRGSCFmatrix to store the copy
idxObject which handles the index conventions for CASSCF
one_rdmIf true, the occupied orbitals get occupation 2

Definition at line 375 of file CASSCF.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::CASSCF::copy_active ( const DMRGSCFmatrix origin,
double *  result,
const DMRGSCFindices idx 
)
static

Copy a one-orbital quantity from DMRGSCFmatrix format to array format.

Parameters
originDMRGSCFmatrix to copy
resultArray to store the copy
idxObject which handles the index conventions for CASSCF

Definition at line 407 of file CASSCF.cpp.

+ Here is the call graph for this function:

double CheMPS2::CASSCF::deviation_from_blockdiag ( DMRGSCFmatrix matrix,
const DMRGSCFindices idx 
)
static

Return the RMS deviation from block-diagonal.

Parameters
matrixMatrix to be assessed
idxObject which handles the index conventions for CASSCF
Returns
RMS deviation from block-diagonal

Definition at line 411 of file CASSCFpt2.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::CASSCF::fillLocalizedOrbitalRotations ( DMRGSCFunitary umat,
DMRGSCFindices idx,
double *  eigenvecs 
)
static

From an Edmiston-Ruedenberg active space rotation, fetch the eigenvectors and store them in eigenvecs.

Parameters
umatThe Edmiston-Ruedenberg active space rotation
idxObject which handles the index conventions for CASSCF
eigenvecsWhere the eigenvectors are stored

Definition at line 160 of file CASSCF.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::CASSCF::fock_dot_4rdm ( double *  fockmx,
CheMPS2::DMRG dmrgsolver,
CheMPS2::Hamiltonian ham,
int  next_orb1,
int  next_orb2,
double *  work,
double *  result,
const bool  CHECKPOINT,
const bool  PSEUDOCANONICAL 
)
static

Build the contraction of the fock matrix with the 4-RDM.

Parameters
fockmxArray of size ham->getL() x ham->getL() containing the Fock matrix elements
dmrgsolverDMRG object which is solved, and for which the 2-RDM and 3-RDM have been calculated as well
hamActive space Hamiltonian, which is needed for the size of the active space and the orbital irreps
next_orb1The next first orbital for which the contraction should be continued
next_orb2The next second orbital for which the contraction should be continued
workWork array of size ham->getL()**6
resultOn entry, contains the partial contraction corresponding to (next_orb1, next_orb2). On exit, contains the full contraction.
CHECKPOINTWhether or not the standard CheMPS2 F.4-RDM checkpoint should be created/updated to continue the contraction at later times
PSEUDOCANONICALWhether or not pseudocanonical orbitals are used in the active space

Definition at line 140 of file CASSCFpt2.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

int CheMPS2::CASSCF::get_num_irreps ( )

Get the number of irreps.

Returns
The number of irreps

Definition at line 117 of file CASSCF.cpp.

bool CheMPS2::CASSCF::read_f4rdm_checkpoint ( const string  f4rdm_file,
int *  hamorb1,
int *  hamorb2,
const int  tot_dmrg_power6,
double *  contract 
)
static

Read the checkpoint file for the contraction of the generalized Fock operator with the 4-RDM from disk.

Parameters
f4rdm_fileThe filename
hamorb1The next hamiltonian orbital 1
hamorb2The next hamiltonian orbital 2
tot_dmrg_power6The size of the array contract
contractThe current partial contraction
Returns
Whether the file was found and read

Definition at line 87 of file CASSCFpt2.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::CASSCF::setDMRG1DM ( const int  num_elec,
const int  LAS,
double *  one_dm,
double *  two_dm 
)
static

Construct the 1-RDM from the 2-RDM.

Parameters
num_elecThe number of DMRG active space electrons
LASThe total number of DMRG orbitals
one_dmThe CASSCF 1-RDM
two_dmThe CASSCF 2-RDM

Definition at line 134 of file CASSCF.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

double CheMPS2::CASSCF::solve ( const int  Nelectrons,
const int  TwoS,
const int  Irrep,
ConvergenceScheme OptScheme,
const int  rootNum,
DMRGSCFoptions scf_options 
)

Do the CASSCF cycles with the augmented Hessian Newton-Raphson method.

Parameters
NelectronsTotal number of electrons in the system: occupied HF orbitals + active space
TwoSTwice the targeted spin
IrrepDesired wave-function irrep
OptSchemeThe optimization scheme to run the inner DMRG loop. If NULL: use FCI instead of DMRG.
rootNumDenotes the targeted state in state-specific CASSCF; 1 means ground state, 2 first excited state etc.
scf_optionsContains the DMRGSCF options
Returns
The converged DMRGSCF energy

Definition at line 56 of file CASSCFnewtonraphson.cpp.

+ Here is the call graph for this function:

void CheMPS2::CASSCF::write_f4rdm_checkpoint ( const string  f4rdm_file,
int *  hamorb1,
int *  hamorb2,
const int  tot_dmrg_power6,
double *  contract 
)
static

Write the checkpoint file for the contraction of the generalized Fock operator with the 4-RDM to disk.

Parameters
f4rdm_fileThe filename
hamorb1The next hamiltonian orbital 1
hamorb2The next hamiltonian orbital 2
tot_dmrg_power6The size of the array contract
contractThe current partial contraction

Definition at line 44 of file CASSCFpt2.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:


The documentation for this class was generated from the following files: