CheMPS2
|
#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... | |
FCI class.
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.
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.
Ham | The Hamiltonian matrix elements |
Nel_up | The number of up (alpha) electrons |
Nel_down | The number of down (beta) electrons |
TargetIrrep | The targeted point group irrep |
maxMemWorkMB | Maximum workspace size in MB to be used for matrix vector product (this does not include the FCI vectors as stored for example in GSDavidson!!) |
FCIverbose | The FCI verbose level: 0 print nothing, 1 print start and solution, 2 print everything |
Definition at line 37 of file FCI.cpp.
|
protected |
Set resultVector to the number operator of a specific site acting on sourceVector.
orbIndex | Orbital index of the number operator |
resultVector | Vector with length getVecLength(0) where the result of the operation should be stored |
sourceVector | Vector with length getVecLength(0) on which the number operator acts |
Definition at line 1957 of file FCI.cpp.
|
protected |
Set thisVector to a creator/annihilator acting on otherVector.
whichOperator | With which operator should be acted on the other FCI state: C means creator and A means annihilator |
isUp | Boolean which denotes if the operator corresponds to an up (alpha) or down (beta) electron |
orbIndex | Orbital index on which the operator acts |
thisVector | Vector with length getVecLength(0) where the result of the operation should be stored |
otherFCI | FCI instance which corresponds to the FCI vector otherVector on which is acted |
otherVector | Vector with length otherFCI->getVecLength(0) which contains the FCI vector on which is acted |
Definition at line 1975 of file FCI.cpp.
|
protected |
Apply E_{crea,anni} to |orig_vector> and store the result in |result_vector>
orig_vector | The original vector |
result_vector | The result vector |
crea | The orbital index of the creator |
anni | The orbital index of the annihilator |
orig_target_irrep | The irrep of the orig_vector |
Definition at line 772 of file FCI.cpp.
|
staticprotected |
Convertor between two representations of a same spin-projection Slater determinant.
Lvalue | The number of orbitals |
bits | Array with the Lvalue bits which should be combined to a single integer |
Definition at line 397 of file FCI.cpp.
double CheMPS2::FCI::CalcSpinSquared | ( | double * | vector | ) | const |
|
protected |
Calculate out = (alpha + beta * Hamiltonian) * in (Econstant is taken into account!!)
alpha | The parameter alpha of the operator |
beta | The parameter beta of the operator |
in | On entry the vector of size getVecLength(0) on which the operator should be applied; unchanged on exit |
out | Array of size getVecLength(0), which contains on exit (alpha + beta * Hamiltonian) * in |
Definition at line 2142 of file FCI.cpp.
|
protected |
Calculate the solution of the system Operator |Sol> = |RESID> with Operator = precon * [ ( alpha + beta * H )^2 + eta^2 ] * precon.
alpha | The parameter alpha of the operator |
beta | The parameter beta of the operator |
eta | The parameter eta of the operator |
precon | The diagonal preconditioner |
Sol | On entry this array of size getVecLength(0) contains the initial guess; on exit the solution is stored here |
RESID | On entry the RHS of the equation is stored in this array of size getVecLength(0), on exit it is overwritten |
PVEC | Workspace of size getVecLength(0) |
OxPVEC | Workspace of size getVecLength(0) |
temp | Workspace of size getVecLength(0) |
temp2 | Workspace of size getVecLength(0) |
|
protected |
Calculate (without approximation) diagonal = diag[ (alpha + beta * Hamiltonian)^2 + eta^2 ] (Econstant is taken into account!!)
alpha | The parameter alpha of the operator |
beta | The parameter beta of the operator |
eta | The parameter eta of the operator |
diagonal | Array of size getVecLength(0), which contains on exit diag[ (alpha + beta * Hamiltonian)^2 + eta^2 ] |
workspace | Workspace of size getVecLength(0) |
Definition at line 2162 of file FCI.cpp.
|
protected |
Calculate out = [(alpha + beta * Hamiltonian)^2 + eta^2] * in (Econstant is taken into account!!)
alpha | The parameter alpha of the operator |
beta | The parameter beta of the operator |
eta | The parameter eta of the operator |
in | On entry the vector of size getVecLength(0) on which the operator should be applied; unchanged on exit |
out | Array of size getVecLength(0), which contains on exit [(alpha + beta * Hamiltonian)^2 + eta^2] * in |
temp | Workspace of size getVecLength(0) |
Definition at line 2153 of file FCI.cpp.
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.
alpha | The real part of the scalar in the operator |
beta | The real-valued prefactor of the Hamiltonian in the operator |
eta | The imaginary part of the scalar in the operator |
RHS | The real-valued right-hand side of the equation with length getVecLength(0) |
RealSol | On exit this array of length getVecLength(0) contains the real part of the solution |
ImagSol | On exit this array of length getVecLength(0) contains the imaginary part of the solution |
checkError | If true, the RMS error without preconditioner will be calculated and printed after convergence |
Definition at line 2066 of file FCI.cpp.
|
static |
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)
omega | The frequency value |
eta | The regularization parameter (... + I*eta in the denominator) |
orb_alpha | The first orbital index |
orb_beta | The second orbital index |
GSenergy | The ground state energy returned by GSDavidson |
GSvector | The ground state vector as calculated by GSDavidson |
RePartGF | On exit RePartGF[0] contains the real part of the density response Green's function |
ImPartGF | On exit ImPartGF[0] contains the imaginary part of the density response Green's function |
Definition at line 2398 of file FCI.cpp.
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>
omega | The frequency value |
eta | The regularization parameter |
orb_alpha | The first orbital index |
orb_beta | The second orbital index |
GSenergy | The ground state energy returned by GSDavidson |
GSvector | The ground state vector as calculated by GSDavidson |
RePartGF | On exit RePartGF[0] contains the real part of the backward propagating part of the density response Green's function |
ImPartGF | On exit ImPartGF[0] contains the imaginary part of the backward propagating part of the density response Green's function |
TwoRDMreal | If not NULL, on exit the 2-RDM of Re{ [ omega + Ham - GSenergy + I*eta ]^{-1} ( n_alpha - <GSvector| n_alpha |GSvector> ) |GSvector> } |
TwoRDMimag | If not NULL, on exit the 2-RDM of Im{ [ omega + Ham - GSenergy + I*eta ]^{-1} ( n_alpha - <GSvector| n_alpha |GSvector> ) |GSvector> } |
TwoRDMdens | If not NULL, on exit the 2-RDM of ( n_alpha - <GSvector| n_alpha |GSvector> ) |GSvector> |
Definition at line 2458 of file FCI.cpp.
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>
omega | The frequency value |
eta | The regularization parameter |
orb_alpha | The first orbital index |
orb_beta | The second orbital index |
GSenergy | The ground state energy returned by GSDavidson |
GSvector | The ground state vector as calculated by GSDavidson |
RePartGF | On exit RePartGF[0] contains the real part of the forward propagating part of the density response Green's function |
ImPartGF | On exit ImPartGF[0] contains the imaginary part of the forward propagating part of the density response Green's function |
TwoRDMreal | If not NULL, on exit the 2-RDM of Re{ [ omega - Ham + GSenergy + I*eta ]^{-1} ( n_beta - <GSvector| n_beta |GSvector> ) |GSvector> } |
TwoRDMimag | If not NULL, on exit the 2-RDM of Im{ [ omega - Ham + GSenergy + I*eta ]^{-1} ( n_beta - <GSvector| n_beta |GSvector> ) |GSvector> } |
TwoRDMdens | If not NULL, on exit the 2-RDM of ( n_beta - <GSvector| n_beta |GSvector> ) |GSvector> |
Definition at line 2422 of file FCI.cpp.
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)
vector | The FCI vector of length getVecLength(0) |
three_rdm | The spin-summed 3-RDM as calculated by Fill3RDM |
orbz | The orbital z which is fixed in Gamma^4(i,j,k,z,p,q,r,z) |
output | To 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.
|
protected |
Function which returns the diagonal elements of the FCI Hamiltonian (without Econstant!!) = Slater determinant energies.
diag | Vector with getVecLength(0) variables which contains on exit the diagonal elements of the FCI Hamiltonian |
Definition at line 1467 of file FCI.cpp.
|
protected |
Function which returns the diagonal elements of the FCI Hamiltonian squared (without Econstant!!)
output | Vector with getVecLength(0) variables which contains on exit the diagonal elements of the FCI Hamiltonian squared |
Definition at line 1504 of file FCI.cpp.
|
staticprotected |
|
staticprotected |
|
staticprotected |
|
staticprotected |
|
staticprotected |
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 ) ) ].
vector | The FCI vector of length getVecLength(0) |
TwoRDM | To store the 2-RDM; needs to be of size getL()^4; point group symmetry shows in 2-RDM elements being zero |
Definition at line 805 of file FCI.cpp.
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 ) ) ) ) ].
vector | The FCI vector of length getVecLength(0) |
ThreeRDM | To 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.
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 ) ) ) ) ) ) ].
vector | The FCI vector of length getVecLength(0) |
FourRDM | To 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.
|
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.
vecLength | The length of the vector; when used for GSDavidson it should be getVecLength(0) |
vec | The vector to fill with random numbers |
Definition at line 1914 of file FCI.cpp.
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)
vector | The FCI vector of length getVecLength(0) |
ThreeRDM | The spin-summed 3-RDM as calculated by Fill3RDM |
Fock | The symmetric Fock operator Fock(i,j) = Fock[ i + L * j ] = Fock[ j + L * i ] |
output | To 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.
|
protected |
Find the bit representation of a global counter corresponding to " E_ij | FCI vector > " ; where irrep_center = I_i x I_j.
irrep_center | The single electron excitation irrep I_i x I_j |
counter | The given global counter corresponding to " E_ij | FCI vector > " |
bits_up | Array of length L to store the bit representation of the up (alpha) electrons in |
bits_down | Array of length L to store the bit representation of the down (beta) electrons in |
Definition at line 417 of file FCI.cpp.
|
inline |
|
inline |
Function which returns the electron repulsion integral (in chemists notation: integral dr1 dr2 orb1(r1) * orb2(r1) * orb3(r2) * orb4(r2) / |r1-r2|)
orb1 | First orbital index (electron at position r1) |
orb2 | Second orbital index (electron at position r1) |
orb3 | Third orbital index (electron at position r2) |
orb4 | Fourth orbital index (electron at position r2) |
Definition at line 82 of file FCI.h.
double CheMPS2::FCI::getFCIcoeff | ( | int * | bits_up, |
int * | bits_down, | ||
double * | vector | ||
) | const |
Function which returns a FCI coefficient.
bits_up | The bit string representation of the up or alpha electron Slater determinant |
bits_down | The bit string representation of the down or beta electron Slater determinant |
vector | The FCI vector with getVecLength(0) variables from which a coefficient is desired |
Definition at line 435 of file FCI.cpp.
|
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) )
orb1 | First orbital index |
orb2 | Second orbital index |
Definition at line 88 of file FCI.h.
|
inline |
|
protected |
Sandwich the Hamiltonian between two Slater determinants (return a specific element) (without Econstant!!)
bits_bra_up | Bit representation of the <bra| Slater determinant of the up (alpha) electrons (length L) |
bits_bra_down | Bit representation of the <bra| Slater determinant of the down (beta) electrons (length L) |
bits_ket_up | Bit representation of the |ket> Slater determinant of the up (alpha) electrons (length L) |
bits_ket_down | Bit representation of the |ket> Slater determinant of the down (beta) electrons (length L) |
work | Work array of length 8 |
Definition at line 1722 of file FCI.cpp.
|
inline |
|
inline |
|
inline |
|
inline |
|
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.
irrep_center | The single electron excitation irrep I_i x I_j |
counter | The given global counter corresponding to " E_ij | FCI vector > " |
Definition at line 409 of file FCI.cpp.
|
inline |
Getter for the number of variables in the vector " E_ij | FCI vector > " ; where irrep_center = I_i x I_j.
irrep_center | The single electron excitation irrep I_i x I_j |
Definition at line 65 of file FCI.h.
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>
alpha | Constant parameter in the resolvent |
beta | Prefector of the Hamiltonian in the resolvent |
eta | The regularization parameter |
orbsLeft | The left orbital indices |
numLeft | The number of left orbital indices |
orbsRight | The right orbital indices |
numRight | The number of right orbital indices |
isUp | If true, the spin projection value of the second quantized operators is up, otherwise it will be down |
GSvector | The ground state vector as calculated by GSDavidson |
Ham | The Hamiltonian, which contains the matrix elements |
RePartsGF | On exit RePartsGF[i+numLeft*j] contains the real part of the addition Green's function |
ImPartsGF | On exit ImPartsGF[i+numLeft*j] contains the imaginary part of the addition Green's function |
TwoRDMreal | If not NULL, TwoRDMreal[j] contains on exit the 2-RDM of Re{ [ alpha + beta * Ham + I*eta ]^{-1} a^+_{orbsRight[j], spin} | GSvector > } |
TwoRDMimag | If not NULL, TwoRDMimag[j] contains on exit the 2-RDM of Im{ [ alpha + beta * Ham + I*eta ]^{-1} a^+_{orbsRight[j], spin} | GSvector > } |
TwoRDMadd | If not NULL, TwoRDMadd[j] contains on exit the 2-RDM of a^+_{orbsRight[j], spin} | GSvector > |
Definition at line 2215 of file FCI.cpp.
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>
alpha | Constant parameter in the resolvent |
beta | Prefector of the Hamiltonian in the resolvent |
eta | The regularization parameter |
orbsLeft | The left orbital indices |
numLeft | The number of left orbital indices |
orbsRight | The right orbital indices |
numRight | The number of right orbital indices |
isUp | If true, the spin projection value of the second quantized operators is up, otherwise it will be down |
GSvector | The ground state vector as calculated by GSDavidson |
Ham | The Hamiltonian, which contains the matrix elements |
RePartsGF | On exit RePartsGF[i+numLeft*j] contains the real part of the removal Green's function GF[i+numLeft*j] |
ImPartsGF | On exit ImPartsGF[i+numLeft*j] contains the imaginary part of the removal Green's function GF[i+numLeft*j] |
TwoRDMreal | If not NULL, TwoRDMreal[j] contains on exit the 2-RDM of Re{ [ alpha + beta * Ham + I*eta ]^{-1} a_{orbsRight[j], spin} | GSvector > } |
TwoRDMimag | If not NULL, TwoRDMimag[j] contains on exit the 2-RDM of Im{ [ alpha + beta * Ham + I*eta ]^{-1} a_{orbsRight[j], spin} | GSvector > } |
TwoRDMrem | If not NULL, TwoRDMrem[j] contains on exit the 2-RDM of a_{orbsRight[j], spin} | GSvector > |
Definition at line 2306 of file FCI.cpp.
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.
inoutput | If 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_VEC | The maximum number of vectors to use in Davidson's algorithm; adjustable in case memory becomes an issue |
Definition at line 1920 of file FCI.cpp.
unsigned int CheMPS2::FCI::LowestEnergyDeterminant | ( | ) | 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.
input | The vector of length getVecLength(0) on which the Hamiltonian should act |
output | Vector of length getVecLength(0) which contains on exit the Hamiltonian times input |
Definition at line 618 of file FCI.cpp.
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)
omega | The frequency value |
eta | The regularization parameter (... + I*eta in the denominator) |
orb_alpha | The first orbital index |
orb_beta | The second orbital index |
isUp | If true, the spin projection value of the second quantized operators is up, otherwise it will be down |
GSenergy | The ground state energy returned by GSDavidson |
GSvector | The ground state vector as calculated by GSDavidson |
Ham | The Hamiltonian, which contains the matrix elements |
RePartGF | On exit RePartGF[0] contains the real part of the retarded Green's function |
ImPartGF | On exit ImPartGF[0] contains the imaginary part of the retarded Green's function |
Definition at line 2191 of file FCI.cpp.
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>
omega | The frequency value |
eta | The regularization parameter |
orb_alpha | The first orbital index |
orb_beta | The second orbital index |
isUp | If true, the spin projection value of the second quantized operators is up, otherwise it will be down |
GSenergy | The ground state energy returned by GSDavidson |
GSvector | The ground state vector as calculated by GSDavidson |
Ham | The Hamiltonian, which contains the matrix elements |
RePartGF | On exit RePartGF[0] contains the real part of the addition part of the retarded Green's function |
ImPartGF | On exit ImPartGF[0] contains the imaginary part of the addition part of the retarded Green's function |
TwoRDMreal | If not NULL, on exit the 2-RDM of Re{ [ omega - Ham + GSenergy + I*eta ]^{-1} a^+_{beta, spin(isUp)} | GSvector > } |
TwoRDMimag | If not NULL, on exit the 2-RDM of Im{ [ omega - Ham + GSenergy + I*eta ]^{-1} a^+_{beta, spin(isUp)} | GSvector > } |
TwoRDMadd | If not NULL, on exit the 2-RDM of a^+_{beta, spin(isUp)} | GSvector > |
Definition at line 2287 of file FCI.cpp.
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>
omega | The frequency value |
eta | The regularization parameter |
orb_alpha | The first orbital index |
orb_beta | The second orbital index |
isUp | If true, the spin projection value of the second quantized operators is up, otherwise it will be down |
GSenergy | The ground state energy returned by GSDavidson |
GSvector | The ground state vector as calculated by GSDavidson |
Ham | The Hamiltonian, which contains the matrix elements |
RePartGF | On exit RePartGF[0] contains the real part of the removal part of the retarded Green's function |
ImPartGF | On exit ImPartGF[0] contains the imaginary part of the removal part of the retarded Green's function |
TwoRDMreal | If not NULL, on exit the 2-RDM of Re{ [ omega + Ham - GSenergy + I*eta ]^{-1} a_{alpha, spin(isUp)} | GSvector > } |
TwoRDMimag | If not NULL, on exit the 2-RDM of Im{ [ omega + Ham - GSenergy + I*eta ]^{-1} a_{alpha, spin(isUp)} | GSvector > } |
TwoRDMrem | If not NULL, on exit the 2-RDM of a_{alpha, spin(isUp)} | GSvector > |
Definition at line 2378 of file FCI.cpp.
|
staticprotected |
Convertor between two representations of a same spin-projection Slater determinant.
Lvalue | The number of orbitals |
bitstring | The input integer, whos bits are the occupation numbers of the orbitals |
bits | Contains on exit the Lvalue bits of bitstring |
Definition at line 391 of file FCI.cpp.