3#include "Coulomb/QkTable.hpp"
4#include "Coulomb/meTable.hpp"
5#include "LinAlg/Matrix.hpp"
10class CorrelationPotential;
22 DiracSpinor::Index w, DiracSpinor::Index x,
23 DiracSpinor::Index y,
int twoJ);
28 DiracSpinor::Index w, DiracSpinor::Index x,
29 DiracSpinor::Index y,
int twoJ);
34 DiracSpinor::Index w, DiracSpinor::Index x,
35 DiracSpinor::Index y,
int twoJ);
54 const std::vector<DiracSpinor> &s1_basis_core,
55 const std::vector<DiracSpinor> &s1_basis_excited,
60 const MBPT::CorrelationPotential &Sigma,
65 const std::vector<DiracSpinor> &cis2_basis,
66 const std::vector<DiracSpinor> &s2_basis_core,
67 const std::vector<DiracSpinor> &s2_basis_excited,
69 bool exclude_wrong_parity_box,
70 bool no_new_integrals =
false);
75 const std::vector<DiracSpinor> &ci_basis,
80std::vector<DiracSpinor>
82 const std::string &subset_string,
83 const std::string &frozen_core_string =
"");
90 const std::vector<CI::CSF2> &CSFAs,
int twoJA,
92 const std::vector<CI::CSF2> &CSFBs,
int twoJB,
102 int K_rank,
int Parity) {
104 Bs.
twoJ(), h, K_rank, Parity);
112std::pair<int, int>
Term_S_L(
int l1,
int l2,
int twoJ,
double gJ_target);
115std::string
Term_Symbol(
int two_J,
int L,
int two_S,
int parity);
117std::string
Term_Symbol(
int L,
int two_S,
int parity);
Very basic two-electron CSF. Only two-electron is implemented.
Definition CSF.hpp:12
Stores the CI Solutions for given J and parity (only two-electron).
Definition CSF.hpp:62
LinAlg::View< const double > coefs(std::size_t i) const
List of CI expansion coefs for the ith CI solution.
Definition CSF.cpp:163
const std::vector< CSF2 > & CSFs() const
Full list of CSFs.
Definition CSF.cpp:151
int twoJ() const
2J for the CI solutions
Definition CSF.cpp:178
Base (pure virtual) class to store Coulomb integrals, and similar. 3 derived classes,...
Definition QkTable.hpp:57
Breit (Hartree-Fock Breit) interaction potential.
Definition Breit.hpp:52
Matrix class; row-major.
Definition Matrix.hpp:35
Proved a "view" onto an array.
Definition Matrix.ipp:7
Functions and classes for Configuration Interaction calculations.
Definition CI_Integrals.cpp:11
double CSF2_Breit(const Coulomb::WkTable &Bk, DiracSpinor::Index v, DiracSpinor::Index w, DiracSpinor::Index x, DiracSpinor::Index y, int twoJ)
Calculates the anti-symmetrised Breit integral for 2-particle states: C1*C2*(b_abcd-b_abdc),...
Definition CI_Integrals.cpp:131
Coulomb::meTable< double > calculate_h1_table(const std::vector< DiracSpinor > &ci_basis, const std::vector< DiracSpinor > &s1_basis_core, const std::vector< DiracSpinor > &s1_basis_excited, const Coulomb::QkTable &qk, bool include_Sigma1)
Calculates table of single-particle matrix elements of one-body Hamiltonian. Note: assumes basis are ...
Definition CI_Integrals.cpp:231
double Sigma2_AB(const CI::CSF2 &A, const CI::CSF2 &B, int twoJ, const Coulomb::LkTable &Sk)
Calculates the Sigma(2) correction to Hab()
Definition CI_Integrals.cpp:179
LinAlg::Matrix< double > construct_Hci(const PsiJPi &psi, const Coulomb::meTable< double > &h1, const Coulomb::QkTable &qk, const Coulomb::WkTable *Bk, const Coulomb::LkTable *Sk)
Constructs the CI matrix, optionally including Sigma corrections.
Definition CI_Integrals.cpp:573
double Breit_AB(const CI::CSF2 &A, const CI::CSF2 &B, int twoJ, const Coulomb::WkTable &Bk)
Calculates the Breit correction to Hab()
Definition CI_Integrals.cpp:187
double RME_CSF2(const CI::CSF2 &X, int twoJX, const CI::CSF2 &V, int twoJV, const Coulomb::meTable< double > &h, int K_rank)
Calculate reduce ME between two 2-particle CSFs - XXX not quite right??
Definition CI_Integrals.cpp:470
std::string Term_Symbol(int two_J, int L, int two_S, int parity)
Returns Term_Symbol as string.
Definition CI_Integrals.cpp:559
Coulomb::WkTable calculate_Bk(const std::string &bk_filename, const HF::Breit *const pBr, const std::vector< DiracSpinor > &ci_basis, int max_k)
Calculates table of single-particle matrix elements of two-body Breit operator.
Definition CI_Integrals.cpp:351
std::pair< int, int > Term_S_L(int l1, int l2, int twoJ, double gJ_target)
Determines best-fit for S and L for two-electron state by matching g-factor.
Definition CI_Integrals.cpp:518
double Hab(const CI::CSF2 &X, const CI::CSF2 &V, int twoJ, const Coulomb::meTable< double > &h1, const Coulomb::QkTable &qk)
Determines CI Hamiltonian matrix element for two 2-particle CSFs, a and b. Does NOT include Sigma(2) ...
Definition CI_Integrals.cpp:196
double ReducedME(const LinAlg::View< const double > &cA, const std::vector< CI::CSF2 > &CSFAs, int twoJA, const LinAlg::View< const double > &cB, const std::vector< CI::CSF2 > &CSFBs, int twoJB, const Coulomb::meTable< double > &h, int K_rank, int Parity)
Calculate reduced matrix elements between two CI states.
Definition CI_Integrals.cpp:431
Coulomb::LkTable calculate_Sk(const std::string &filename, const std::vector< DiracSpinor > &cis2_basis, const std::vector< DiracSpinor > &s2_basis_core, const std::vector< DiracSpinor > &s2_basis_excited, const Coulomb::QkTable &qk, int max_k, bool exclude_wrong_parity_box, bool no_new_integrals)
Calculates table of single-particle matrix elements of two-body Sigma_2 operator.
Definition CI_Integrals.cpp:302
std::vector< DiracSpinor > basis_subset(const std::vector< DiracSpinor > &basis, const std::string &subset_string, const std::string &frozen_core_string)
Takes a subset of input basis according to subset_string. Only states not included in frozen_core_str...
Definition CI_Integrals.cpp:391
double CSF2_Sigma2(const Coulomb::LkTable &Sk, DiracSpinor::Index v, DiracSpinor::Index w, DiracSpinor::Index x, DiracSpinor::Index y, int twoJ)
Calculates the correlation [Sigma(2)] correction to CSF2_Coulomb(). Sk is the table of two-body singl...
Definition CI_Integrals.cpp:75
double CSF2_Coulomb(const Coulomb::QkTable &qk, DiracSpinor::Index v, DiracSpinor::Index w, DiracSpinor::Index x, DiracSpinor::Index y, int twoJ)
Calculates the anti-symmetrised Coulomb integral for 2-particle states: C1*C2*(g_abcd-g_abdc),...
Definition CI_Integrals.cpp:16
Functions and classes for Hartree-Fock.
Definition CI_Integrals.hpp:12
Many-body perturbation theory.
Definition CI_Integrals.hpp:9