CheMPS2
CheMPS2::FCI Class Reference

#include <FCI.h>

Public Member Functions

 FCI (CheMPS2::Hamiltonian *Ham, const unsigned int Nel_up, const unsigned int Nel_down, const int TargetIrrep, const double maxMemWorkMB=100.0, const int FCIverbose=2)
 Constructor. More...
 
virtual ~FCI ()
 Destructor.
 
unsigned int getL () const
 Getter for the number of orbitals. More...
 
unsigned int getNel_up () const
 Getter for the number of up or alpha electrons. More...
 
unsigned int getNel_down () const
 Getter for the number of down or beta electrons. More...
 
unsigned int getVecLength (const int irrep_center) const
 Getter for the number of variables in the vector " E_ij | FCI vector > " ; where irrep_center = I_i x I_j. More...
 
int getTargetIrrep () const
 Get the target irrep. More...
 
int getOrb2Irrep (const int orb) const
 Get an orbital irrep. More...
 
double getERI (const int orb1, const int orb2, const int orb3, const int orb4) const
 Function which returns the electron repulsion integral (in chemists notation: integral dr1 dr2 orb1(r1) * orb2(r1) * orb3(r2) * orb4(r2) / |r1-r2|) More...
 
double getGmat (const int orb1, const int orb2) const
 Function which returns the AUGMENTED one-body matrix elements ( G_{ij} = T_{ij} - 0.5 * sum_k ERI_{ikkj} ; see Chem. Phys. Lett. 111 (4-5), 315-321 (1984) ) More...
 
double getEconst () const
 Function which returns the nuclear repulsion energy. More...
 
double GSDavidson (double *inoutput=NULL, const int DVDSN_NUM_VEC=CheMPS2::DAVIDSON_NUM_VEC) const
 Calculates the FCI ground state with Davidson's algorithm. More...
 
unsigned int LowestEnergyDeterminant () const
 Return the global counter of the Slater determinant with the lowest energy. More...
 
double Fill2RDM (double *vector, double *TwoRDM) const
 Construct the (spin-summed) 2-RDM of a FCI vector: Gamma^2(i,j,k,l) = sum_sigma,tau < a^+_i,sigma a^+_j,tau a_l,tau a_k,sigma > = TwoRDM[ i + L * ( j + L * ( k + L * l ) ) ]. More...
 
void Fill3RDM (double *vector, double *ThreeRDM) const
 Construct the (spin-summed) 3-RDM of a FCI vector: Gamma^3(i,j,k,l,m,n) = sum_sigma,tau,s < a^+_{i,sigma} a^+_{j,tau} a^+_{k,s} a_{n,s} a_{m,tau} a_{l,sigma} > = ThreeRDM[ i + L * ( j + L * ( k + L * ( l + L * ( m + L * n ) ) ) ) ]. More...
 
void Fill4RDM (double *vector, double *FourRDM) const
 Construct the (spin-summed) 4-RDM of a FCI vector: Gamma^4(i,j,k,l,p,q,r,t) = sum_sigma,tau,s,z < a^+_{i,sigma} a^+_{j,tau} a^+_{k,s} a^+_{l,z} a_{t,z} a_{r,s} a_{q,tau} a_{p,sigma} > = FourRDM[ i + L * ( j + L * ( k + L * ( l + L * ( p + L * ( q + L * ( r + L * t ) ) ) ) ) ) ]. More...
 
void Fock4RDM (double *vector, double *ThreeRDM, double *Fock, double *output) const
 Construct the (spin-summed) contraction of the 4-RDM with the Fock operator: output(i,j,k,p,q,r) = sum_{l,t} Fock(l,t) * Gamma^4(i,j,k,l,p,q,r,t) More...
 
void Diag4RDM (double *vector, double *three_rdm, const unsigned int orbz, double *output) const
 Construct part of the 4-RDM: output(i,j,k,p,q,r) = Gamma^4(i,j,k,z,p,q,r,z) More...
 
double CalcSpinSquared (double *vector) const
 Measure S(S+1) (spin squared) More...
 
void RetardedGF (const double omega, const double eta, const unsigned int orb_alpha, const unsigned int orb_beta, const bool isUp, const double GSenergy, double *GSvector, CheMPS2::Hamiltonian *Ham, double *RePartGF, double *ImPartGF) const
 Calculate the retarded Green's function (= addition + removal amplitude) More...
 
void RetardedGF_addition (const double omega, const double eta, const unsigned int orb_alpha, const unsigned int orb_beta, const bool isUp, const double GSenergy, double *GSvector, CheMPS2::Hamiltonian *Ham, double *RePartGF, double *ImPartGF, double *TwoRDMreal=NULL, double *TwoRDMimag=NULL, double *TwoRDMadd=NULL) const
 Calculate the addition part of the retarded Green's function: <GSvector| a_{alpha, spin(isUp)} [ omega - Ham + GSenergy + I*eta ]^{-1} a^+_{beta, spin(isUp)} |GSvector> More...
 
void GFmatrix_addition (const double alpha, const double beta, const double eta, int *orbsLeft, const unsigned int numLeft, int *orbsRight, const unsigned int numRight, const bool isUp, double *GSvector, CheMPS2::Hamiltonian *Ham, double *RePartsGF, double *ImPartsGF, double **TwoRDMreal=NULL, double **TwoRDMimag=NULL, double **TwoRDMadd=NULL) const
 Calculate the addition Green's function: GF[i+numLeft*j] = <GSvector| a_{orbsLeft[i], spin} [ alpha + beta * Ham + I*eta ]^{-1} a^+_{orbsRight[j], spin} |GSvector> More...
 
void RetardedGF_removal (const double omega, const double eta, const unsigned int orb_alpha, const unsigned int orb_beta, const bool isUp, const double GSenergy, double *GSvector, CheMPS2::Hamiltonian *Ham, double *RePartGF, double *ImPartGF, double *TwoRDMreal=NULL, double *TwoRDMimag=NULL, double *TwoRDMrem=NULL) const
 Calculate the removal part of the retarded Green's function: <GSvector| a^+_{beta, spin(isUp)} [ omega + Ham - GSenergy + I*eta ]^{-1} a_{alpha, spin(isUp)} |GSvector> More...
 
void GFmatrix_removal (const double alpha, const double beta, const double eta, int *orbsLeft, const unsigned int numLeft, int *orbsRight, const unsigned int numRight, const bool isUp, double *GSvector, CheMPS2::Hamiltonian *Ham, double *RePartsGF, double *ImPartsGF, double **TwoRDMreal=NULL, double **TwoRDMimag=NULL, double **TwoRDMrem=NULL) const
 Calculate the removal Green's function: GF[i+numLeft*j] = <GSvector| a^+_{orbsLeft[i], spin} [ alpha + beta * Ham + I*eta ]^{-1} a_{orbsRight[j], spin} |GSvector> More...
 
void DensityResponseGF (const double omega, const double eta, const unsigned int orb_alpha, const unsigned int orb_beta, const double GSenergy, double *GSvector, double *RePartGF, double *ImPartGF) const
 Calculate the density response Green's function (= forward - backward propagating part) More...
 
void DensityResponseGF_forward (const double omega, const double eta, const unsigned int orb_alpha, const unsigned int orb_beta, const double GSenergy, double *GSvector, double *RePartGF, double *ImPartGF, double *TwoRDMreal=NULL, double *TwoRDMimag=NULL, double *TwoRDMdens=NULL) const
 Calculate the forward propagating part of the density response Green's function: <GSvector| ( n_alpha - <GSvector| n_alpha |GSvector> ) [ omega - Ham + GSenergy + I*eta ]^{-1} ( n_beta - <GSvector| n_beta |GSvector> ) |GSvector> More...
 
void DensityResponseGF_backward (const double omega, const double eta, const unsigned int orb_alpha, const unsigned int orb_beta, const double GSenergy, double *GSvector, double *RePartGF, double *ImPartGF, double *TwoRDMreal=NULL, double *TwoRDMimag=NULL, double *TwoRDMdens=NULL) const
 Calculate the backward propagating part of the density response Green's function: <GSvector| ( n_beta - <GSvector| n_beta |GSvector> ) [ omega + Ham - GSenergy + I*eta ]^{-1} ( n_alpha - <GSvector| n_alpha |GSvector> ) |GSvector> More...
 
void CGSolveSystem (const double alpha, const double beta, const double eta, double *RHS, double *RealSol, double *ImagSol, const bool checkError=true) const
 Calculate the solution of the equation ( alpha + beta * Hamiltonian + I * eta ) Solution = RHS with conjugate gradient. More...
 
double getFCIcoeff (int *bits_up, int *bits_down, double *vector) const
 Function which returns a FCI coefficient. More...
 

Static Public Member Functions

static void FillRandom (const unsigned int vecLength, double *vec)
 Fill a vector with random numbers in the interval [-1,1[; used when output for GSDavidson is desired but no specific input can be given. More...
 
static void ClearVector (const unsigned int vecLength, double *vec)
 Set the entries of a vector to zero. More...
 

Protected Member Functions

void DiagHam (double *diag) const
 Function which returns the diagonal elements of the FCI Hamiltonian (without Econstant!!) = Slater determinant energies. More...
 
void DiagHamSquared (double *output) const
 Function which returns the diagonal elements of the FCI Hamiltonian squared (without Econstant!!) More...
 
void matvec (double *input, double *output) const
 Function which performs the Hamiltonian times Vector product (without Econstant!!) making use of the (ij|kl) = (ji|kl) = (ij|lk) = (ji|lk) symmetry of the electron repulsion integrals. More...
 
double GetMatrixElement (int *bits_bra_up, int *bits_bra_down, int *bits_ket_up, int *bits_ket_down, int *work) const
 Sandwich the Hamiltonian between two Slater determinants (return a specific element) (without Econstant!!) More...
 
void getBitsOfCounter (const int irrep_center, const unsigned int counter, int *bits_up, int *bits_down) const
 Find the bit representation of a global counter corresponding to " E_ij | FCI vector > " ; where irrep_center = I_i x I_j. More...
 
int getUpIrrepOfCounter (const int irrep_center, const unsigned int counter) const
 Find the irrep of the up Slater determinant of a global counter corresponding to " E_ij | FCI vector > " ; where irrep_center = I_i x I_j. More...
 
void ActWithSecondQuantizedOperator (const char whichOperator, const bool isUp, const unsigned int orbIndex, double *thisVector, const FCI *otherFCI, double *otherVector) const
 Set thisVector to a creator/annihilator acting on otherVector. More...
 
void ActWithNumberOperator (const unsigned int orbIndex, double *resultVector, double *sourceVector) const
 Set resultVector to the number operator of a specific site acting on sourceVector. More...
 
void CGCoreSolver (const double alpha, const double beta, const double eta, double *precon, double *Sol, double *RESID, double *PVEC, double *OxPVEC, double *temp, double *temp2) const
 Calculate the solution of the system Operator |Sol> = |RESID> with Operator = precon * [ ( alpha + beta * H )^2 + eta^2 ] * precon. More...
 
void CGAlphaPlusBetaHAM (const double alpha, const double beta, double *in, double *out) const
 Calculate out = (alpha + beta * Hamiltonian) * in (Econstant is taken into account!!) More...
 
void CGoperator (const double alpha, const double beta, const double eta, double *in, double *temp, double *out) const
 Calculate out = [(alpha + beta * Hamiltonian)^2 + eta^2] * in (Econstant is taken into account!!) More...
 
void CGdiagonal (const double alpha, const double beta, const double eta, double *diagonal, double *workspace) const
 Calculate (without approximation) diagonal = diag[ (alpha + beta * Hamiltonian)^2 + eta^2 ] (Econstant is taken into account!!) More...
 
void apply_excitation (double *orig_vector, double *result_vector, const int crea, const int anni, const int orig_target_irrep) const
 Apply E_{crea,anni} to |orig_vector> and store the result in |result_vector> More...
 

Static Protected Member Functions

static void str2bits (const unsigned int Lvalue, const unsigned int bitstring, int *bits)
 Convertor between two representations of a same spin-projection Slater determinant. More...
 
static unsigned int bits2str (const unsigned int Lvalue, int *bits)
 Convertor between two representations of a same spin-projection Slater determinant. More...
 
static double FCIddot (const unsigned int vecLength, double *vec1, double *vec2)
 Take the inproduct of two vectors. More...
 
static void FCIdcopy (const unsigned int vecLength, double *origin, double *target)
 Copy a vector. More...
 
static double FCIfrobeniusnorm (const unsigned int vecLength, double *vec)
 Calculate the 2-norm of a vector. More...
 
static void FCIdaxpy (const unsigned int vecLength, const double alpha, double *vec_x, double *vec_y)
 Do lapack's daxpy vec_y += alpha * vec_x. More...
 
static void FCIdscal (const unsigned int vecLength, const double alpha, double *vec)
 Do lapack's dscal vec *= alpha. More...
 

Detailed Description

FCI class.

Author
Sebastian Wouters sebas.nosp@m.tian.nosp@m.woute.nosp@m.rs@g.nosp@m.mail..nosp@m.com
Date
November 6, 2014

The FCI class performs the full configuration interaction ground state calculation in a given particle number, irrep, and spin projection sector of a given Hamiltonian. It also contains the functionality to calculate Green's functions.

Definition at line 32 of file FCI.h.

Constructor & Destructor Documentation

CheMPS2::FCI::FCI ( CheMPS2::Hamiltonian Ham,
const unsigned int  Nel_up,
const unsigned int  Nel_down,
const int  TargetIrrep,
const double  maxMemWorkMB = 100.0,
const int  FCIverbose = 2 
)

Constructor.

Parameters
HamThe Hamiltonian matrix elements
Nel_upThe number of up (alpha) electrons
Nel_downThe number of down (beta) electrons
TargetIrrepThe targeted point group irrep
maxMemWorkMBMaximum workspace size in MB to be used for matrix vector product (this does not include the FCI vectors as stored for example in GSDavidson!!)
FCIverboseThe FCI verbose level: 0 print nothing, 1 print start and solution, 2 print everything

Definition at line 37 of file FCI.cpp.

+ Here is the call graph for this function:

Member Function Documentation

void CheMPS2::FCI::ActWithNumberOperator ( const unsigned int  orbIndex,
double *  resultVector,
double *  sourceVector 
) const
protected

Set resultVector to the number operator of a specific site acting on sourceVector.

Parameters
orbIndexOrbital index of the number operator
resultVectorVector with length getVecLength(0) where the result of the operation should be stored
sourceVectorVector with length getVecLength(0) on which the number operator acts

Definition at line 1957 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::ActWithSecondQuantizedOperator ( const char  whichOperator,
const bool  isUp,
const unsigned int  orbIndex,
double *  thisVector,
const FCI otherFCI,
double *  otherVector 
) const
protected

Set thisVector to a creator/annihilator acting on otherVector.

Parameters
whichOperatorWith which operator should be acted on the other FCI state: C means creator and A means annihilator
isUpBoolean which denotes if the operator corresponds to an up (alpha) or down (beta) electron
orbIndexOrbital index on which the operator acts
thisVectorVector with length getVecLength(0) where the result of the operation should be stored
otherFCIFCI instance which corresponds to the FCI vector otherVector on which is acted
otherVectorVector with length otherFCI->getVecLength(0) which contains the FCI vector on which is acted

Definition at line 1975 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::apply_excitation ( double *  orig_vector,
double *  result_vector,
const int  crea,
const int  anni,
const int  orig_target_irrep 
) const
protected

Apply E_{crea,anni} to |orig_vector> and store the result in |result_vector>

Parameters
orig_vectorThe original vector
result_vectorThe result vector
creaThe orbital index of the creator
anniThe orbital index of the annihilator
orig_target_irrepThe irrep of the orig_vector

Definition at line 772 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

unsigned int CheMPS2::FCI::bits2str ( const unsigned int  Lvalue,
int *  bits 
)
staticprotected

Convertor between two representations of a same spin-projection Slater determinant.

Parameters
LvalueThe number of orbitals
bitsArray with the Lvalue bits which should be combined to a single integer
Returns
The single integer which represents the bits in bits

Definition at line 397 of file FCI.cpp.

+ Here is the caller graph for this function:

double CheMPS2::FCI::CalcSpinSquared ( double *  vector) const

Measure S(S+1) (spin squared)

Parameters
vectorThe FCI vector of length getVecLength(0)
Returns
Measured value of S(S+1)

Definition at line 1405 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::CGAlphaPlusBetaHAM ( const double  alpha,
const double  beta,
double *  in,
double *  out 
) const
protected

Calculate out = (alpha + beta * Hamiltonian) * in (Econstant is taken into account!!)

Parameters
alphaThe parameter alpha of the operator
betaThe parameter beta of the operator
inOn entry the vector of size getVecLength(0) on which the operator should be applied; unchanged on exit
outArray of size getVecLength(0), which contains on exit (alpha + beta * Hamiltonian) * in

Definition at line 2142 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::CGCoreSolver ( const double  alpha,
const double  beta,
const double  eta,
double *  precon,
double *  Sol,
double *  RESID,
double *  PVEC,
double *  OxPVEC,
double *  temp,
double *  temp2 
) const
protected

Calculate the solution of the system Operator |Sol> = |RESID> with Operator = precon * [ ( alpha + beta * H )^2 + eta^2 ] * precon.

Parameters
alphaThe parameter alpha of the operator
betaThe parameter beta of the operator
etaThe parameter eta of the operator
preconThe diagonal preconditioner
SolOn entry this array of size getVecLength(0) contains the initial guess; on exit the solution is stored here
RESIDOn entry the RHS of the equation is stored in this array of size getVecLength(0), on exit it is overwritten
PVECWorkspace of size getVecLength(0)
OxPVECWorkspace of size getVecLength(0)
tempWorkspace of size getVecLength(0)
temp2Workspace of size getVecLength(0)

+ Here is the caller graph for this function:

void CheMPS2::FCI::CGdiagonal ( const double  alpha,
const double  beta,
const double  eta,
double *  diagonal,
double *  workspace 
) const
protected

Calculate (without approximation) diagonal = diag[ (alpha + beta * Hamiltonian)^2 + eta^2 ] (Econstant is taken into account!!)

Parameters
alphaThe parameter alpha of the operator
betaThe parameter beta of the operator
etaThe parameter eta of the operator
diagonalArray of size getVecLength(0), which contains on exit diag[ (alpha + beta * Hamiltonian)^2 + eta^2 ]
workspaceWorkspace of size getVecLength(0)

Definition at line 2162 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::CGoperator ( const double  alpha,
const double  beta,
const double  eta,
double *  in,
double *  temp,
double *  out 
) const
protected

Calculate out = [(alpha + beta * Hamiltonian)^2 + eta^2] * in (Econstant is taken into account!!)

Parameters
alphaThe parameter alpha of the operator
betaThe parameter beta of the operator
etaThe parameter eta of the operator
inOn entry the vector of size getVecLength(0) on which the operator should be applied; unchanged on exit
outArray of size getVecLength(0), which contains on exit [(alpha + beta * Hamiltonian)^2 + eta^2] * in
tempWorkspace of size getVecLength(0)

Definition at line 2153 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::CGSolveSystem ( const double  alpha,
const double  beta,
const double  eta,
double *  RHS,
double *  RealSol,
double *  ImagSol,
const bool  checkError = true 
) const

Calculate the solution of the equation ( alpha + beta * Hamiltonian + I * eta ) Solution = RHS with conjugate gradient.

Parameters
alphaThe real part of the scalar in the operator
betaThe real-valued prefactor of the Hamiltonian in the operator
etaThe imaginary part of the scalar in the operator
RHSThe real-valued right-hand side of the equation with length getVecLength(0)
RealSolOn exit this array of length getVecLength(0) contains the real part of the solution
ImagSolOn exit this array of length getVecLength(0) contains the imaginary part of the solution
checkErrorIf true, the RMS error without preconditioner will be calculated and printed after convergence

Definition at line 2066 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::ClearVector ( const unsigned int  vecLength,
double *  vec 
)
static

Set the entries of a vector to zero.

Parameters
vecLengthThe vector length
vecThe vector which has to be set to zero

Definition at line 1908 of file FCI.cpp.

+ Here is the caller graph for this function:

void CheMPS2::FCI::DensityResponseGF ( const double  omega,
const double  eta,
const unsigned int  orb_alpha,
const unsigned int  orb_beta,
const double  GSenergy,
double *  GSvector,
double *  RePartGF,
double *  ImPartGF 
) const

Calculate the density response Green's function (= forward - backward propagating part)

Parameters
omegaThe frequency value
etaThe regularization parameter (... + I*eta in the denominator)
orb_alphaThe first orbital index
orb_betaThe second orbital index
GSenergyThe ground state energy returned by GSDavidson
GSvectorThe ground state vector as calculated by GSDavidson
RePartGFOn exit RePartGF[0] contains the real part of the density response Green's function
ImPartGFOn exit ImPartGF[0] contains the imaginary part of the density response Green's function

Definition at line 2398 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::DensityResponseGF_backward ( const double  omega,
const double  eta,
const unsigned int  orb_alpha,
const unsigned int  orb_beta,
const double  GSenergy,
double *  GSvector,
double *  RePartGF,
double *  ImPartGF,
double *  TwoRDMreal = NULL,
double *  TwoRDMimag = NULL,
double *  TwoRDMdens = NULL 
) const

Calculate the backward propagating part of the density response Green's function: <GSvector| ( n_beta - <GSvector| n_beta |GSvector> ) [ omega + Ham - GSenergy + I*eta ]^{-1} ( n_alpha - <GSvector| n_alpha |GSvector> ) |GSvector>

Parameters
omegaThe frequency value
etaThe regularization parameter
orb_alphaThe first orbital index
orb_betaThe second orbital index
GSenergyThe ground state energy returned by GSDavidson
GSvectorThe ground state vector as calculated by GSDavidson
RePartGFOn exit RePartGF[0] contains the real part of the backward propagating part of the density response Green's function
ImPartGFOn exit ImPartGF[0] contains the imaginary part of the backward propagating part of the density response Green's function
TwoRDMrealIf not NULL, on exit the 2-RDM of Re{ [ omega + Ham - GSenergy + I*eta ]^{-1} ( n_alpha - <GSvector| n_alpha |GSvector> ) |GSvector> }
TwoRDMimagIf not NULL, on exit the 2-RDM of Im{ [ omega + Ham - GSenergy + I*eta ]^{-1} ( n_alpha - <GSvector| n_alpha |GSvector> ) |GSvector> }
TwoRDMdensIf not NULL, on exit the 2-RDM of ( n_alpha - <GSvector| n_alpha |GSvector> ) |GSvector>

Definition at line 2458 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::DensityResponseGF_forward ( const double  omega,
const double  eta,
const unsigned int  orb_alpha,
const unsigned int  orb_beta,
const double  GSenergy,
double *  GSvector,
double *  RePartGF,
double *  ImPartGF,
double *  TwoRDMreal = NULL,
double *  TwoRDMimag = NULL,
double *  TwoRDMdens = NULL 
) const

Calculate the forward propagating part of the density response Green's function: <GSvector| ( n_alpha - <GSvector| n_alpha |GSvector> ) [ omega - Ham + GSenergy + I*eta ]^{-1} ( n_beta - <GSvector| n_beta |GSvector> ) |GSvector>

Parameters
omegaThe frequency value
etaThe regularization parameter
orb_alphaThe first orbital index
orb_betaThe second orbital index
GSenergyThe ground state energy returned by GSDavidson
GSvectorThe ground state vector as calculated by GSDavidson
RePartGFOn exit RePartGF[0] contains the real part of the forward propagating part of the density response Green's function
ImPartGFOn exit ImPartGF[0] contains the imaginary part of the forward propagating part of the density response Green's function
TwoRDMrealIf not NULL, on exit the 2-RDM of Re{ [ omega - Ham + GSenergy + I*eta ]^{-1} ( n_beta - <GSvector| n_beta |GSvector> ) |GSvector> }
TwoRDMimagIf not NULL, on exit the 2-RDM of Im{ [ omega - Ham + GSenergy + I*eta ]^{-1} ( n_beta - <GSvector| n_beta |GSvector> ) |GSvector> }
TwoRDMdensIf not NULL, on exit the 2-RDM of ( n_beta - <GSvector| n_beta |GSvector> ) |GSvector>

Definition at line 2422 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::Diag4RDM ( double *  vector,
double *  three_rdm,
const unsigned int  orbz,
double *  output 
) const

Construct part of the 4-RDM: output(i,j,k,p,q,r) = Gamma^4(i,j,k,z,p,q,r,z)

Parameters
vectorThe FCI vector of length getVecLength(0)
three_rdmThe spin-summed 3-RDM as calculated by Fill3RDM
orbzThe orbital z which is fixed in Gamma^4(i,j,k,z,p,q,r,z)
outputTo store part of the 4-RDM output(i,j,k,p,q,r) = output[ i + L * ( j + L * ( k + L * ( p + L * ( q + L * r ) ) ) ) ]; needs to be of size getL()^6; point group symmetry shows in elements being zero; has 12-fold permutation symmetry just like 3-RDM

Definition at line 1176 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::DiagHam ( double *  diag) const
protected

Function which returns the diagonal elements of the FCI Hamiltonian (without Econstant!!) = Slater determinant energies.

Parameters
diagVector with getVecLength(0) variables which contains on exit the diagonal elements of the FCI Hamiltonian

Definition at line 1467 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::DiagHamSquared ( double *  output) const
protected

Function which returns the diagonal elements of the FCI Hamiltonian squared (without Econstant!!)

Parameters
outputVector with getVecLength(0) variables which contains on exit the diagonal elements of the FCI Hamiltonian squared

Definition at line 1504 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::FCIdaxpy ( const unsigned int  vecLength,
const double  alpha,
double *  vec_x,
double *  vec_y 
)
staticprotected

Do lapack's daxpy vec_y += alpha * vec_x.

Parameters
vecLengthThe vector length
alphaThe scalar factor
vec_xThe vector which has to be added to vec_y in rescaled form
vec_yThe target vector

Definition at line 1890 of file FCI.cpp.

+ Here is the caller graph for this function:

void CheMPS2::FCI::FCIdcopy ( const unsigned int  vecLength,
double *  origin,
double *  target 
)
staticprotected

Copy a vector.

Parameters
vecLengthThe vector length
originVector to be copied
targetWhere to copy the vector to

Definition at line 1868 of file FCI.cpp.

+ Here is the caller graph for this function:

double CheMPS2::FCI::FCIddot ( const unsigned int  vecLength,
double *  vec1,
double *  vec2 
)
staticprotected

Take the inproduct of two vectors.

Parameters
vecLengthThe vector length
vec1The first vector
vec2The second vector
Returns
The inproduct < vec1 | vec2 >

Definition at line 1876 of file FCI.cpp.

+ Here is the caller graph for this function:

void CheMPS2::FCI::FCIdscal ( const unsigned int  vecLength,
const double  alpha,
double *  vec 
)
staticprotected

Do lapack's dscal vec *= alpha.

Parameters
vecLengthThe vector length
alphaThe scalar factor
vecThe vector which has to be rescaled

Definition at line 1899 of file FCI.cpp.

+ Here is the caller graph for this function:

double CheMPS2::FCI::FCIfrobeniusnorm ( const unsigned int  vecLength,
double *  vec 
)
staticprotected

Calculate the 2-norm of a vector.

Parameters
vecLengthThe vector length
vecThe vector
Returns
The 2-norm of vec

Definition at line 1884 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

double CheMPS2::FCI::Fill2RDM ( double *  vector,
double *  TwoRDM 
) const

Construct the (spin-summed) 2-RDM of a FCI vector: Gamma^2(i,j,k,l) = sum_sigma,tau < a^+_i,sigma a^+_j,tau a_l,tau a_k,sigma > = TwoRDM[ i + L * ( j + L * ( k + L * l ) ) ].

Parameters
vectorThe FCI vector of length getVecLength(0)
TwoRDMTo store the 2-RDM; needs to be of size getL()^4; point group symmetry shows in 2-RDM elements being zero
Returns
The energy of the given FCI vector, calculated by contraction of the 2-RDM with Gmat and ERI

Definition at line 805 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::Fill3RDM ( double *  vector,
double *  ThreeRDM 
) const

Construct the (spin-summed) 3-RDM of a FCI vector: Gamma^3(i,j,k,l,m,n) = sum_sigma,tau,s < a^+_{i,sigma} a^+_{j,tau} a^+_{k,s} a_{n,s} a_{m,tau} a_{l,sigma} > = ThreeRDM[ i + L * ( j + L * ( k + L * ( l + L * ( m + L * n ) ) ) ) ].

Parameters
vectorThe FCI vector of length getVecLength(0)
ThreeRDMTo store the 3-RDM; needs to be of size getL()^6; point group symmetry shows in 3-RDM elements being zero

Definition at line 1168 of file FCI.cpp.

+ Here is the caller graph for this function:

void CheMPS2::FCI::Fill4RDM ( double *  vector,
double *  FourRDM 
) const

Construct the (spin-summed) 4-RDM of a FCI vector: Gamma^4(i,j,k,l,p,q,r,t) = sum_sigma,tau,s,z < a^+_{i,sigma} a^+_{j,tau} a^+_{k,s} a^+_{l,z} a_{t,z} a_{r,s} a_{q,tau} a_{p,sigma} > = FourRDM[ i + L * ( j + L * ( k + L * ( l + L * ( p + L * ( q + L * ( r + L * t ) ) ) ) ) ) ].

Parameters
vectorThe FCI vector of length getVecLength(0)
FourRDMTo store the 4-RDM; needs to be of size getL()^8; point group symmetry shows in 4-RDM elements being zero

Definition at line 896 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::FillRandom ( const unsigned int  vecLength,
double *  vec 
)
static

Fill a vector with random numbers in the interval [-1,1[; used when output for GSDavidson is desired but no specific input can be given.

Parameters
vecLengthThe length of the vector; when used for GSDavidson it should be getVecLength(0)
vecThe vector to fill with random numbers

Definition at line 1914 of file FCI.cpp.

+ Here is the caller graph for this function:

void CheMPS2::FCI::Fock4RDM ( double *  vector,
double *  ThreeRDM,
double *  Fock,
double *  output 
) const

Construct the (spin-summed) contraction of the 4-RDM with the Fock operator: output(i,j,k,p,q,r) = sum_{l,t} Fock(l,t) * Gamma^4(i,j,k,l,p,q,r,t)

Parameters
vectorThe FCI vector of length getVecLength(0)
ThreeRDMThe spin-summed 3-RDM as calculated by Fill3RDM
FockThe symmetric Fock operator Fock(i,j) = Fock[ i + L * j ] = Fock[ j + L * i ]
outputTo store the contraction output(i,j,k,p,q,r) = output[ i + L * ( j + L * ( k + L * ( p + L * ( q + L * r ) ) ) ) ]; needs to be of size getL()^6; point group symmetry shows in elements being zero; has 12-fold permutation symmetry just like 3-RDM

Definition at line 1160 of file FCI.cpp.

+ Here is the caller graph for this function:

void CheMPS2::FCI::getBitsOfCounter ( const int  irrep_center,
const unsigned int  counter,
int *  bits_up,
int *  bits_down 
) const
protected

Find the bit representation of a global counter corresponding to " E_ij | FCI vector > " ; where irrep_center = I_i x I_j.

Parameters
irrep_centerThe single electron excitation irrep I_i x I_j
counterThe given global counter corresponding to " E_ij | FCI vector > "
bits_upArray of length L to store the bit representation of the up (alpha) electrons in
bits_downArray of length L to store the bit representation of the down (beta) electrons in

Definition at line 417 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

double CheMPS2::FCI::getEconst ( ) const
inline

Function which returns the nuclear repulsion energy.

Returns
The nuclear repulsion energy

Definition at line 92 of file FCI.h.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

double CheMPS2::FCI::getERI ( const int  orb1,
const int  orb2,
const int  orb3,
const int  orb4 
) const
inline

Function which returns the electron repulsion integral (in chemists notation: integral dr1 dr2 orb1(r1) * orb2(r1) * orb3(r2) * orb4(r2) / |r1-r2|)

Parameters
orb1First orbital index (electron at position r1)
orb2Second orbital index (electron at position r1)
orb3Third orbital index (electron at position r2)
orb4Fourth orbital index (electron at position r2)
Returns
The desired electron repulsion integral

Definition at line 82 of file FCI.h.

+ Here is the caller graph for this function:

double CheMPS2::FCI::getFCIcoeff ( int *  bits_up,
int *  bits_down,
double *  vector 
) const

Function which returns a FCI coefficient.

Parameters
bits_upThe bit string representation of the up or alpha electron Slater determinant
bits_downThe bit string representation of the down or beta electron Slater determinant
vectorThe FCI vector with getVecLength(0) variables from which a coefficient is desired
Returns
The corresponding FCI coefficient; 0.0 if bits_up and bits_down do not form a valid FCI determinant

Definition at line 435 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

double CheMPS2::FCI::getGmat ( const int  orb1,
const int  orb2 
) const
inline

Function which returns the AUGMENTED one-body matrix elements ( G_{ij} = T_{ij} - 0.5 * sum_k ERI_{ikkj} ; see Chem. Phys. Lett. 111 (4-5), 315-321 (1984) )

Parameters
orb1First orbital index
orb2Second orbital index
Returns
The desired AUGMENTED one-body matrix element

Definition at line 88 of file FCI.h.

+ Here is the caller graph for this function:

unsigned int CheMPS2::FCI::getL ( ) const
inline

Getter for the number of orbitals.

Returns
The number of orbitals

Definition at line 52 of file FCI.h.

+ Here is the caller graph for this function:

double CheMPS2::FCI::GetMatrixElement ( int *  bits_bra_up,
int *  bits_bra_down,
int *  bits_ket_up,
int *  bits_ket_down,
int *  work 
) const
protected

Sandwich the Hamiltonian between two Slater determinants (return a specific element) (without Econstant!!)

Parameters
bits_bra_upBit representation of the <bra| Slater determinant of the up (alpha) electrons (length L)
bits_bra_downBit representation of the <bra| Slater determinant of the down (beta) electrons (length L)
bits_ket_upBit representation of the |ket> Slater determinant of the up (alpha) electrons (length L)
bits_ket_downBit representation of the |ket> Slater determinant of the down (beta) electrons (length L)
workWork array of length 8
Returns
The FCI Hamiltonian element which connects the given two Slater determinants

Definition at line 1722 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

unsigned int CheMPS2::FCI::getNel_down ( ) const
inline

Getter for the number of down or beta electrons.

Returns
The number of down or beta electrons

Definition at line 60 of file FCI.h.

+ Here is the caller graph for this function:

unsigned int CheMPS2::FCI::getNel_up ( ) const
inline

Getter for the number of up or alpha electrons.

Returns
The number of up or alpha electrons

Definition at line 56 of file FCI.h.

+ Here is the caller graph for this function:

int CheMPS2::FCI::getOrb2Irrep ( const int  orb) const
inline

Get an orbital irrep.

Parameters
orbThe orbital index
Returns
The irrep of orbital orb

Definition at line 74 of file FCI.h.

+ Here is the caller graph for this function:

int CheMPS2::FCI::getTargetIrrep ( ) const
inline

Get the target irrep.

Returns
The target irrep

Definition at line 69 of file FCI.h.

+ Here is the caller graph for this function:

int CheMPS2::FCI::getUpIrrepOfCounter ( const int  irrep_center,
const unsigned int  counter 
) const
protected

Find the irrep of the up Slater determinant of a global counter corresponding to " E_ij | FCI vector > " ; where irrep_center = I_i x I_j.

Parameters
irrep_centerThe single electron excitation irrep I_i x I_j
counterThe given global counter corresponding to " E_ij | FCI vector > "
Returns
The corresponding irrep of the up Slater determinant

Definition at line 409 of file FCI.cpp.

+ Here is the caller graph for this function:

unsigned int CheMPS2::FCI::getVecLength ( const int  irrep_center) const
inline

Getter for the number of variables in the vector " E_ij | FCI vector > " ; where irrep_center = I_i x I_j.

Parameters
irrep_centerThe single electron excitation irrep I_i x I_j
Returns
The number of variables in the corresponding vector

Definition at line 65 of file FCI.h.

+ Here is the caller graph for this function:

void CheMPS2::FCI::GFmatrix_addition ( const double  alpha,
const double  beta,
const double  eta,
int *  orbsLeft,
const unsigned int  numLeft,
int *  orbsRight,
const unsigned int  numRight,
const bool  isUp,
double *  GSvector,
CheMPS2::Hamiltonian Ham,
double *  RePartsGF,
double *  ImPartsGF,
double **  TwoRDMreal = NULL,
double **  TwoRDMimag = NULL,
double **  TwoRDMadd = NULL 
) const

Calculate the addition Green's function: GF[i+numLeft*j] = <GSvector| a_{orbsLeft[i], spin} [ alpha + beta * Ham + I*eta ]^{-1} a^+_{orbsRight[j], spin} |GSvector>

Parameters
alphaConstant parameter in the resolvent
betaPrefector of the Hamiltonian in the resolvent
etaThe regularization parameter
orbsLeftThe left orbital indices
numLeftThe number of left orbital indices
orbsRightThe right orbital indices
numRightThe number of right orbital indices
isUpIf true, the spin projection value of the second quantized operators is up, otherwise it will be down
GSvectorThe ground state vector as calculated by GSDavidson
HamThe Hamiltonian, which contains the matrix elements
RePartsGFOn exit RePartsGF[i+numLeft*j] contains the real part of the addition Green's function
ImPartsGFOn exit ImPartsGF[i+numLeft*j] contains the imaginary part of the addition Green's function
TwoRDMrealIf not NULL, TwoRDMreal[j] contains on exit the 2-RDM of Re{ [ alpha + beta * Ham + I*eta ]^{-1} a^+_{orbsRight[j], spin} | GSvector > }
TwoRDMimagIf not NULL, TwoRDMimag[j] contains on exit the 2-RDM of Im{ [ alpha + beta * Ham + I*eta ]^{-1} a^+_{orbsRight[j], spin} | GSvector > }
TwoRDMaddIf not NULL, TwoRDMadd[j] contains on exit the 2-RDM of a^+_{orbsRight[j], spin} | GSvector >

Definition at line 2215 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::GFmatrix_removal ( const double  alpha,
const double  beta,
const double  eta,
int *  orbsLeft,
const unsigned int  numLeft,
int *  orbsRight,
const unsigned int  numRight,
const bool  isUp,
double *  GSvector,
CheMPS2::Hamiltonian Ham,
double *  RePartsGF,
double *  ImPartsGF,
double **  TwoRDMreal = NULL,
double **  TwoRDMimag = NULL,
double **  TwoRDMrem = NULL 
) const

Calculate the removal Green's function: GF[i+numLeft*j] = <GSvector| a^+_{orbsLeft[i], spin} [ alpha + beta * Ham + I*eta ]^{-1} a_{orbsRight[j], spin} |GSvector>

Parameters
alphaConstant parameter in the resolvent
betaPrefector of the Hamiltonian in the resolvent
etaThe regularization parameter
orbsLeftThe left orbital indices
numLeftThe number of left orbital indices
orbsRightThe right orbital indices
numRightThe number of right orbital indices
isUpIf true, the spin projection value of the second quantized operators is up, otherwise it will be down
GSvectorThe ground state vector as calculated by GSDavidson
HamThe Hamiltonian, which contains the matrix elements
RePartsGFOn exit RePartsGF[i+numLeft*j] contains the real part of the removal Green's function GF[i+numLeft*j]
ImPartsGFOn exit ImPartsGF[i+numLeft*j] contains the imaginary part of the removal Green's function GF[i+numLeft*j]
TwoRDMrealIf not NULL, TwoRDMreal[j] contains on exit the 2-RDM of Re{ [ alpha + beta * Ham + I*eta ]^{-1} a_{orbsRight[j], spin} | GSvector > }
TwoRDMimagIf not NULL, TwoRDMimag[j] contains on exit the 2-RDM of Im{ [ alpha + beta * Ham + I*eta ]^{-1} a_{orbsRight[j], spin} | GSvector > }
TwoRDMremIf not NULL, TwoRDMrem[j] contains on exit the 2-RDM of a_{orbsRight[j], spin} | GSvector >

Definition at line 2306 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

double CheMPS2::FCI::GSDavidson ( double *  inoutput = NULL,
const int  DVDSN_NUM_VEC = CheMPS2::DAVIDSON_NUM_VEC 
) const

Calculates the FCI ground state with Davidson's algorithm.

Parameters
inoutputIf inoutput!=NULL, vector with getVecLength(0) variables which contains the initial guess at the start, and on exit the solution of the FCI calculation
DVDSN_NUM_VECThe maximum number of vectors to use in Davidson's algorithm; adjustable in case memory becomes an issue
Returns
The ground state energy

Definition at line 1920 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

unsigned int CheMPS2::FCI::LowestEnergyDeterminant ( ) const

Return the global counter of the Slater determinant with the lowest energy.

Returns
The global counter of the Slater determinant with the lowest energy

Definition at line 1700 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::matvec ( double *  input,
double *  output 
) const
protected

Function which performs the Hamiltonian times Vector product (without Econstant!!) making use of the (ij|kl) = (ji|kl) = (ij|lk) = (ji|lk) symmetry of the electron repulsion integrals.

Parameters
inputThe vector of length getVecLength(0) on which the Hamiltonian should act
outputVector of length getVecLength(0) which contains on exit the Hamiltonian times input

Definition at line 618 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::RetardedGF ( const double  omega,
const double  eta,
const unsigned int  orb_alpha,
const unsigned int  orb_beta,
const bool  isUp,
const double  GSenergy,
double *  GSvector,
CheMPS2::Hamiltonian Ham,
double *  RePartGF,
double *  ImPartGF 
) const

Calculate the retarded Green's function (= addition + removal amplitude)

Parameters
omegaThe frequency value
etaThe regularization parameter (... + I*eta in the denominator)
orb_alphaThe first orbital index
orb_betaThe second orbital index
isUpIf true, the spin projection value of the second quantized operators is up, otherwise it will be down
GSenergyThe ground state energy returned by GSDavidson
GSvectorThe ground state vector as calculated by GSDavidson
HamThe Hamiltonian, which contains the matrix elements
RePartGFOn exit RePartGF[0] contains the real part of the retarded Green's function
ImPartGFOn exit ImPartGF[0] contains the imaginary part of the retarded Green's function

Definition at line 2191 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::RetardedGF_addition ( const double  omega,
const double  eta,
const unsigned int  orb_alpha,
const unsigned int  orb_beta,
const bool  isUp,
const double  GSenergy,
double *  GSvector,
CheMPS2::Hamiltonian Ham,
double *  RePartGF,
double *  ImPartGF,
double *  TwoRDMreal = NULL,
double *  TwoRDMimag = NULL,
double *  TwoRDMadd = NULL 
) const

Calculate the addition part of the retarded Green's function: <GSvector| a_{alpha, spin(isUp)} [ omega - Ham + GSenergy + I*eta ]^{-1} a^+_{beta, spin(isUp)} |GSvector>

Parameters
omegaThe frequency value
etaThe regularization parameter
orb_alphaThe first orbital index
orb_betaThe second orbital index
isUpIf true, the spin projection value of the second quantized operators is up, otherwise it will be down
GSenergyThe ground state energy returned by GSDavidson
GSvectorThe ground state vector as calculated by GSDavidson
HamThe Hamiltonian, which contains the matrix elements
RePartGFOn exit RePartGF[0] contains the real part of the addition part of the retarded Green's function
ImPartGFOn exit ImPartGF[0] contains the imaginary part of the addition part of the retarded Green's function
TwoRDMrealIf not NULL, on exit the 2-RDM of Re{ [ omega - Ham + GSenergy + I*eta ]^{-1} a^+_{beta, spin(isUp)} | GSvector > }
TwoRDMimagIf not NULL, on exit the 2-RDM of Im{ [ omega - Ham + GSenergy + I*eta ]^{-1} a^+_{beta, spin(isUp)} | GSvector > }
TwoRDMaddIf not NULL, on exit the 2-RDM of a^+_{beta, spin(isUp)} | GSvector >

Definition at line 2287 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::RetardedGF_removal ( const double  omega,
const double  eta,
const unsigned int  orb_alpha,
const unsigned int  orb_beta,
const bool  isUp,
const double  GSenergy,
double *  GSvector,
CheMPS2::Hamiltonian Ham,
double *  RePartGF,
double *  ImPartGF,
double *  TwoRDMreal = NULL,
double *  TwoRDMimag = NULL,
double *  TwoRDMrem = NULL 
) const

Calculate the removal part of the retarded Green's function: <GSvector| a^+_{beta, spin(isUp)} [ omega + Ham - GSenergy + I*eta ]^{-1} a_{alpha, spin(isUp)} |GSvector>

Parameters
omegaThe frequency value
etaThe regularization parameter
orb_alphaThe first orbital index
orb_betaThe second orbital index
isUpIf true, the spin projection value of the second quantized operators is up, otherwise it will be down
GSenergyThe ground state energy returned by GSDavidson
GSvectorThe ground state vector as calculated by GSDavidson
HamThe Hamiltonian, which contains the matrix elements
RePartGFOn exit RePartGF[0] contains the real part of the removal part of the retarded Green's function
ImPartGFOn exit ImPartGF[0] contains the imaginary part of the removal part of the retarded Green's function
TwoRDMrealIf not NULL, on exit the 2-RDM of Re{ [ omega + Ham - GSenergy + I*eta ]^{-1} a_{alpha, spin(isUp)} | GSvector > }
TwoRDMimagIf not NULL, on exit the 2-RDM of Im{ [ omega + Ham - GSenergy + I*eta ]^{-1} a_{alpha, spin(isUp)} | GSvector > }
TwoRDMremIf not NULL, on exit the 2-RDM of a_{alpha, spin(isUp)} | GSvector >

Definition at line 2378 of file FCI.cpp.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void CheMPS2::FCI::str2bits ( const unsigned int  Lvalue,
const unsigned int  bitstring,
int *  bits 
)
staticprotected

Convertor between two representations of a same spin-projection Slater determinant.

Parameters
LvalueThe number of orbitals
bitstringThe input integer, whos bits are the occupation numbers of the orbitals
bitsContains on exit the Lvalue bits of bitstring

Definition at line 391 of file FCI.cpp.

+ Here is the caller graph for this function:


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