ampsci
High-precision calculations for one- and two-valence atomic systems
CI Namespace Reference

Detailed Description

Functions and classes for Configuration Interaction calculations.

Main functions are:

Main Classes are:

Classes

struct  ConfigInfo
 Configuration metadata for a single CI level. More...
 
class  CSF2
 Two-electron configuration state function (CSF). More...
 
class  PsiJPi
 Container for CI solutions in a single (J, parity) sector. More...
 

Functions

double CSF2_Coulomb (const Coulomb::QkTable &qk, DiracSpinor::Index v, DiracSpinor::Index w, DiracSpinor::Index x, DiracSpinor::Index y, int twoJ)
 Antisymmetrised two-body Coulomb matrix element in the coupled CSF basis.
 
double CSF2_Sigma2 (const Coulomb::LkTable &Sk, DiracSpinor::Index v, DiracSpinor::Index w, DiracSpinor::Index x, DiracSpinor::Index y, int twoJ)
 Two-body \( \Sigma_2 \) (MBPT) correction to CSF2_Coulomb().
 
double CSF2_Breit (const Coulomb::WkTable &Bk, DiracSpinor::Index v, DiracSpinor::Index w, DiracSpinor::Index x, DiracSpinor::Index y, int twoJ)
 Antisymmetrised two-body Breit matrix element in the coupled CSF basis.
 
double Sigma2_AB (const CI::CSF2 &A, const CI::CSF2 &B, int twoJ, const Coulomb::LkTable &Sk)
 Two-body \( \Sigma_2 \) correction to Hab().
 
double Breit_AB (const CI::CSF2 &A, const CI::CSF2 &B, int twoJ, const Coulomb::WkTable &Bk)
 Breit correction to Hab().
 
double Hab (const CI::CSF2 &A, const CI::CSF2 &B, int twoJ, const Coulomb::meTable< double > &h1, const Coulomb::QkTable &qk)
 CI Hamiltonian matrix element between two two-electron CSFs.
 
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)
 Builds the one-body Hamiltonian matrix element table for the CI basis.
 
Coulomb::meTable< double > calculate_h1_table (const std::vector< DiracSpinor > &ci_basis, const MBPT::CorrelationPotential &Sigma, bool include_Sigma1)
 Builds the one-body Hamiltonian table using a precomputed CorrelationPotential.
 
Coulomb::WkTable calculate_Bk (const std::string &bk_filename, const HF::Breit *const pBr, const std::vector< DiracSpinor > &ci_basis, int max_k, bool no_new_integralsQ=false)
 Builds or loads the two-body Breit integral table.
 
std::vector< DiracSpinorbasis_subset (const std::vector< DiracSpinor > &basis, const std::string &subset_string, const std::string &frozen_core_string="")
 Returns the subset of basis matching subset_string, excluding states in frozen_core_string.
 
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)
 Reduced matrix element between two CI states (low-level overload).
 
double RME_CSF2 (const CI::CSF2 &X, int twoJX, const CI::CSF2 &V, int twoJV, const Coulomb::meTable< double > &h, int K_rank)
 Reduced matrix element between two two-electron CSFs.
 
std::pair< int, int > Term_S_L (int l1, int l2, int twoJ, double gJ_target)
 Determines the best-fit (S, L) term for a two-electron state by matching the g-factor.
 
std::string Term_Symbol (int two_J, int L, int two_S, int parity)
 Returns spectroscopic term symbol string, e.g. "3P_1".
 
std::string Term_Symbol (int L, int two_S, int parity)
 Returns term symbol without the J subscript, e.g. "3P".
 
LinAlg::Matrix< double > construct_Hci (const PsiJPi &psi, const Coulomb::meTable< double > &h1, const Coulomb::QkTable &qk, const Coulomb::WkTable *Bk=nullptr, const Coulomb::LkTable *Sk=nullptr)
 Constructs the full CI Hamiltonian matrix in the CSF basis.
 
double ReducedME (const PsiJPi &As, std::size_t iA, const PsiJPi &Bs, std::size_t iB, const Coulomb::meTable< double > &h, int K_rank, int Parity)
 Reduced matrix element between two CI states (PsiJPi overload).
 
std::vector< PsiJPiconfiguration_interaction (const IO::InputBlock &input, const Wavefunction &wf)
 Runs Configuration Interation: returns CI solutions for all requested J and parity values.
 
PsiJPi run_CI (const std::vector< DiracSpinor > &ci_sp_basis, int twoJ, int parity, int num_solutions, std::optional< double > all_below_cm, const Coulomb::meTable< double > &h1, const Coulomb::QkTable &qk, const Coulomb::WkTable &Bk, const Coulomb::LkTable &Sk, bool include_Sigma2, bool print_details, std::ostream &outstream=std::cout)
 Constructs and solves the CI eigenvalue problem for a single J,pi.
 
bool operator== (const CSF2 &A, const CSF2 &B)
 
bool operator!= (const CSF2 &A, const CSF2 &B)
 
std::vector< CSF2form_CSFs (int twoJ, int parity, const std::vector< DiracSpinor > &cisp_basis)
 Forms all two-electron CSFs with given total J and parity.
 

Class Documentation

◆ CI::ConfigInfo

struct CI::ConfigInfo
Class Members
string config {} Dominant configuration label (typically non-relativistic notation)
double ci2 {0.0} Squared CI coefficient of the dominant configuration (or sum over non-rel degenerates)
double gJ {0.0}
double L {-1.0} Approximate orbital angular momentum L (-1 if not assigned)
double twoS {-1.0} Twice the approximate spin S (-1 if not assigned)

Function Documentation

◆ CSF2_Coulomb()

double CI::CSF2_Coulomb ( const Coulomb::QkTable qk,
DiracSpinor::Index  v,
DiracSpinor::Index  w,
DiracSpinor::Index  x,
DiracSpinor::Index  y,
int  twoJ 
)

Antisymmetrised two-body Coulomb matrix element in the coupled CSF basis.

Evaluates the angular-reduced, antisymmetrised Coulomb interaction between two two-electron CSFs \( |vw; J\rangle \) and \( |xy; J\rangle \):

\[ \langle vw; J \| g \| xy; J \rangle = \eta_{vw}\eta_{xy} \sum_k (-1)^{j_v+j_x+k+J} \begin{Bmatrix} j_v & j_w & J \\ j_y & j_x & k \end{Bmatrix} Q^k_{vwxy} + \text{exchange}, \]

where \( \eta_{ab} = 1/\sqrt{2} \) if \( a = b \) (identical-particle normalisation) and 1 otherwise, and \( Q^k \) are the Coulomb integrals stored in qk.

Parameters
qkTable of Coulomb \( Q^k \) integrals.
v,wIndices of the bra single-particle states.
x,yIndices of the ket single-particle states.
twoJTwice the total angular momentum 2J of the coupled pair.
Returns
Antisymmetrised, angular-reduced two-body Coulomb matrix element.

◆ CSF2_Sigma2()

double CI::CSF2_Sigma2 ( const Coulomb::LkTable Sk,
DiracSpinor::Index  v,
DiracSpinor::Index  w,
DiracSpinor::Index  x,
DiracSpinor::Index  y,
int  twoJ 
)

Two-body \( \Sigma_2 \) (MBPT) correction to CSF2_Coulomb().

Evaluates the same angular reduction as CSF2_Coulomb(), but using the two-body \( \Sigma_2 \) integrals \( S^k \) stored in Sk in place of the Coulomb \( Q^k \) integrals. Adds the second-order MBPT correction to the two-electron interaction.

Parameters
SkTable of two-body \( \Sigma_2 \) ( \( L^k \)) integrals.
v,wIndices of the bra single-particle states.
x,yIndices of the ket single-particle states.
twoJTwice the total angular momentum 2J of the coupled pair.
Returns
Antisymmetrised two-body \( \Sigma_2 \) matrix element.

◆ CSF2_Breit()

double CI::CSF2_Breit ( const Coulomb::WkTable Bk,
DiracSpinor::Index  v,
DiracSpinor::Index  w,
DiracSpinor::Index  x,
DiracSpinor::Index  y,
int  twoJ 
)

Antisymmetrised two-body Breit matrix element in the coupled CSF basis.

Evaluates the same angular reduction as CSF2_Coulomb(), but using the Breit \( B^k \) integrals stored in Bk.

Parameters
BkTable of Breit \( W^k \) integrals.
v,wIndices of the bra single-particle states.
x,yIndices of the ket single-particle states.
twoJTwice the total angular momentum 2J of the coupled pair.
Returns
Antisymmetrised two-body Breit matrix element.

◆ Sigma2_AB()

double CI::Sigma2_AB ( const CI::CSF2 A,
const CI::CSF2 B,
int  twoJ,
const Coulomb::LkTable Sk 
)

Two-body \( \Sigma_2 \) correction to Hab().

Evaluates the MBPT \( \Sigma_2 \) contribution to the CI matrix element using CSF2_Sigma2(). Add to Hab() to form the full CI+MBPT Hamiltonian matrix element.

Parameters
A,BThe two CSFs.
twoJTwice the total angular momentum 2J.
SkTable of \( \Sigma_2 \) ( \( L^k \)) integrals.
Returns
\( \Sigma_2 \) correction to \( H_{AB} \).

◆ Breit_AB()

double CI::Breit_AB ( const CI::CSF2 A,
const CI::CSF2 B,
int  twoJ,
const Coulomb::WkTable Bk 
)

Breit correction to Hab().

Evaluates the two-body Breit contribution to the CI matrix element using CSF2_Breit(). Add to Hab() to include the Breit interaction.

Parameters
A,BThe two CSFs.
twoJTwice the total angular momentum 2J.
BkTable of Breit \( W^k \) integrals.
Returns
Breit correction to \( H_{AB} \).

◆ Hab()

double CI::Hab ( const CI::CSF2 A,
const CI::CSF2 B,
int  twoJ,
const Coulomb::meTable< double > &  h1,
const Coulomb::QkTable qk 
)

CI Hamiltonian matrix element between two two-electron CSFs.

Computes \( H_{AB} = \langle A | \hat{H} | B \rangle \) using the Slater-Condon rules, including one-body terms from h1 (which may already incorporate \( \Sigma_1 \) corrections) and the two-body Coulomb interaction via CSF2_Coulomb().

Does NOT include \( \Sigma_2 \) or Breit corrections; add those via Sigma2_AB() and Breit_AB() respectively.

Parameters
A,BThe two CSFs.
twoJTwice the total angular momentum 2J.
h1Table of one-body matrix elements \( \langle a | h_1 | b \rangle \).
qkTable of Coulomb \( Q^k \) integrals.
Returns
CI Hamiltonian matrix element \( H_{AB} \).

◆ calculate_h1_table() [1/2]

Coulomb::meTable< double > CI::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 
)

Builds the one-body Hamiltonian matrix element table for the CI basis.

Constructs a lookup table of single-particle matrix elements \( \langle a | h_1 | b \rangle \) for all pairs \( a, b \) in ci_basis. The diagonal elements are the HF single-particle energies.

If include_Sigma1 is true, the one-body MBPT \( \Sigma_1 \) correction is computed from the Coulomb integrals in qk using s1_basis_core and s1_basis_excited as the internal lines of the MBPT diagrams and added to the diagonal.

Parameters
ci_basisBasis states for which table entries are needed.
s1_basis_coreCore states used as internal lines for \( \Sigma_1 \).
s1_basis_excitedExcited states used as internal lines for \( \Sigma_1 \).
qkTable of Coulomb \( Q^k \) integrals.
include_Sigma1If true, add one-body MBPT \( \Sigma_1 \) corrections.
Returns
Table of \( \langle a | h_1 | b \rangle \) matrix elements.
Warning
Assumes ci_basis states are Hartree-Fock eigenstates, so off-diagonal HF terms vanish.

◆ calculate_h1_table() [2/2]

Coulomb::meTable< double > CI::calculate_h1_table ( const std::vector< DiracSpinor > &  ci_basis,
const MBPT::CorrelationPotential &  Sigma,
bool  include_Sigma1 
)

Builds the one-body Hamiltonian table using a precomputed CorrelationPotential.

Overload of calculate_h1_table() that uses a CorrelationPotential object (i.e., a precomputed \( \Sigma_1 \) operator) instead of computing MBPT diagrams on the fly. Preferred when a CorrelationPotential is available, as it is generally faster and more complete.

Parameters
ci_basisBasis states for which table entries are needed.
SigmaPrecomputed one-body correlation potential \( \Sigma_1 \).
include_Sigma1If true, include \( \Sigma_1 \) corrections from Sigma.
Returns
Table of \( \langle a | h_1 | b \rangle \) matrix elements.

◆ calculate_Bk()

Coulomb::WkTable CI::calculate_Bk ( const std::string &  bk_filename,
const HF::Breit *const  pBr,
const std::vector< DiracSpinor > &  ci_basis,
int  max_k,
bool  no_new_integralsQ = false 
)

Builds or loads the two-body Breit integral table.

Computes Breit \( W^k \) integrals for all pairs in ci_basis using the Breit operator pBr. Results are cached to/from bk_filename.

If pBr is nullptr or no_new_integralsQ is true, no new integrals are computed; only cached values are loaded.

Parameters
bk_filenameFilename for caching the \( W^k \) table.
pBrPointer to Breit operator; if nullptr, returns empty table.
ci_basisBasis for which Breit integrals are needed.
max_kMaximum multipolarity k to include.
no_new_integralsQIf true, skip computing any new integrals.
Returns
Table of Breit \( W^k \) integrals.

◆ basis_subset()

std::vector< DiracSpinor > CI::basis_subset ( const std::vector< DiracSpinor > &  basis,
const std::string &  subset_string,
const std::string &  frozen_core_string = "" 
)

Returns the subset of basis matching subset_string, excluding states in frozen_core_string.

Filters basis to retain only states described by the ampsci basis-string notation (e.g., "20spdf") that are not part of the frozen core.

Parameters
basisFull single-particle basis to filter.
subset_stringBasis-string specifying which states to keep.
frozen_core_stringBasis-string specifying core states to exclude.
Returns
Filtered basis vector.

◆ ReducedME() [1/2]

double CI::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 
)

Reduced matrix element between two CI states (low-level overload).

Evaluates the reduced matrix element of a rank-K_rank tensor operator between two CI states:

\[ \redmatel{A}{T^K}{B} = \sum_{ij} c_i^A \, c_j^B \, \redmatel{\text{CSF}_i}{T^K}{\text{CSF}_j}, \]

where the single-particle reduced matrix elements are looked up from h.

Parameters
cA,cBCI expansion coefficient vectors for states A and B.
CSFAs,CSFBsCSF bases for states A and B respectively.
twoJA,twoJBTwice the total angular momentum of states A and B.
hLookup table of single-particle reduced matrix elements.
K_rankRank of the tensor operator.
ParityParity of the operator (+1 or -1).
Returns
Reduced matrix element \( \redmatel{A}{T^K}{B} \).

◆ RME_CSF2()

double CI::RME_CSF2 ( const CI::CSF2 X,
int  twoJX,
const CI::CSF2 V,
int  twoJV,
const Coulomb::meTable< double > &  h,
int  K_rank 
)

Reduced matrix element between two two-electron CSFs.

Evaluates \( \redmatel{X; J_X}{T^K}{V; J_V} \) for a rank-K_rank one-body tensor operator using the standard 6j angular reduction, accounting for identical-particle normalisation factors.

Warning
This function may not handle all cases correctly; results should be verified for non-trivial configurations.

◆ Term_S_L()

std::pair< int, int > CI::Term_S_L ( int  l1,
int  l2,
int  twoJ,
double  gJ_target 
)

Determines the best-fit (S, L) term for a two-electron state by matching the g-factor.

Iterates over all allowed (S, L) combinations for given orbital angular momenta l1, l2 and total twoJ /2, and returns the pair whose Lande g-factor is closest to gJ_target.

Parameters
l1,l2Orbital angular momenta of the two electrons.
twoJTwice the total angular momentum 2J.
gJ_targetTarget g-factor to match.
Returns
Best-fit {2S, L} pair.

◆ Term_Symbol() [1/2]

std::string CI::Term_Symbol ( int  two_J,
int  L,
int  two_S,
int  parity 
)

Returns spectroscopic term symbol string, e.g. "3P_1".

◆ Term_Symbol() [2/2]

std::string CI::Term_Symbol ( int  L,
int  two_S,
int  parity 
)

Returns term symbol without the J subscript, e.g. "3P".

◆ construct_Hci()

LinAlg::Matrix< double > CI::construct_Hci ( const PsiJPi psi,
const Coulomb::meTable< double > &  h1,
const Coulomb::QkTable qk,
const Coulomb::WkTable Bk = nullptr,
const Coulomb::LkTable Sk = nullptr 
)

Constructs the full CI Hamiltonian matrix in the CSF basis.

Builds the symmetric matrix \( H_{AB} \) for all CSF pairs in psi, calling Hab() for each element and optionally adding Breit and \( \Sigma_2 \) corrections.

Parameters
psiCI solution container holding the CSF basis and J/parity.
h1One-body matrix element table (may include \( \Sigma_1 \)).
qkCoulomb \( Q^k \) integral table.
BkPointer to Breit \( W^k \) table; ignored if nullptr.
SkPointer to \( \Sigma_2 \) \( L^k \) table; ignored if nullptr.
Returns
Full CI Hamiltonian matrix in the CSF basis.

◆ ReducedME() [2/2]

double CI::ReducedME ( const PsiJPi As,
std::size_t  iA,
const PsiJPi Bs,
std::size_t  iB,
const Coulomb::meTable< double > &  h,
int  K_rank,
int  Parity 
)
inline

Reduced matrix element between two CI states (PsiJPi overload).

Convenience wrapper around the low-level ReducedME() overload. Extracts expansion coefficients and CSF lists from As and Bs for the requested solution indices iA and iB.

Parameters
As,BsCI solution containers for the two states.
iA,iBSolution indices within As and Bs.
hLookup table of single-particle reduced matrix elements.
K_rankRank of the tensor operator.
ParityParity of the operator (+1 or -1).
Returns
Reduced matrix element \( \redmatel{A}{T^K}{B} \).

◆ configuration_interaction()

std::vector< PsiJPi > CI::configuration_interaction ( const IO::InputBlock input,
const Wavefunction wf 
)

Runs Configuration Interation: returns CI solutions for all requested J and parity values.

Reads options from input, the CI Input Block, and constructs the CI basis from wf, computes the required Coulomb (and optionally Breit and two-body MBPT) integrals, then calls run_CI() for each requested (J, parity) pair.

The returned vector contains one PsiJPi per {J, parity} combination, each holding the eigenvalues and CI expansion coefficients for the requested number of solutions.

As of writing, options are:

// Available CI options/blocks
ci_basis;
// Basis used for CI expansion; must be a sub-set of
// full ampsci basis [default: 20spdf]
J;
// List of total angular momentum J for CI solutions
// (comma separated). Must be integers (two-electron
// only). []
J+;
// As above, but for EVEN CSFs only (takes precedence
// over J).
J-;
// As above, but for ODD CSFs (takes precedence over J).
num_solutions;
// Number of CI solutions to find (for each J/pi) [5]
all_below_cm;
// Find all CI solutions for energies below this
// threshold, in inverse cm. Note that this is the total
// energy, not the excitation energy. If set,
// num_solutions is ignored.
sigma1;
// Include one-body MBPT correlations? [false]
sigma2;
// Include two-body MBPT correlations? [false]
Brueckner;
// Use Brueckner (spectrum) states for CI basis? Must
// have Spectrum and sigma1. [false]
cis2_basis;
// The subset of ci_basis for which the two-body MBPT
// corrections are calculated. Must be a subset of
// ci_basis. If existing sk file has more integrals,
// they will be used. [default: Nspdf, where N is
// maximum n for core + 3]
Breit2;
// Include two-body Breit? Default is true if Breit
// included in HF. Ignored if Breit not included in HF.
// [true]
Breit_basis;
// Subset of ci_basis used to include two-body Breit
// corrections into CI matrix. Large basis is slow, uses
// huge memory, and makes small contribution. [default:
// Nspdf, where N is maximum n for core + 6]
s1_basis;
// Usually should be left as default. Basis used for the
// one-body MBPT diagrams (Sigma^1) internal lines.
// These are the most important, so in general the
// default (all basis states) should be used. Must be a
// subset of full ampsci basis. [default: full basis]
// - Note: if CorrelationPotential is available, it
// will be used instead of calculating the Sigma_1
// integrals
s2_basis;
// Usually should be left blank. Basis used for internal
// lines of the two-body MBPT diagrams (Sigma^2)
// internal lines. Must be a subset of s1_basis.
// [default: s1_basis]
n_min_core;
// Minimum n for core to be included in MBPT [1]
max_k;
// Maximum k (multipolarity) to include when calculating
// new Coulomb integrals. Higher k often contribute
// negligably. Note: if qk file already has higher-k
// terms, they will be included. Set negative (or very
// large) to include all k. [8]
denominators;
// 'RS', 'Fermi', 'Fermi0'. Denominators used in Sigma2
// matrix elements. RS uses actual states for external
// legs, Fermi uses the lowest excited state for each
// kappa, Fermi0 uses lowest excited state for all
// kappas (and thus cancels in all except diagram 'd').
// [Fermi0]
qk_file;
// Filename for storing two-body Coulomb integrals. By
// default, is ~ At.qk, where At is atomic symbol +
// 'identity'.
sk_file;
// Filename for storing two-body Sigma_2 integrals. By
// default, is At_n_b_k.sk, where At is atomic symbol, n
// is n_min_core, b is cis2_basis, k is max_k.
bk_file;
// Filename for storing two-body Breit integrals. By
// default, is ~ At.bk, where At is atomic symbol +
// 'identity'.
no_new_integrals;
// Usually false. If set to true, ampsci will not
// calculate any new Coulomb or Sigma_2 integrals, even
// if they are implied by the above settings. This saves
// time when we know all required integrals already
// exist, since the code doesn't need to check. [true]
exclude_wrong_parity_box;
// Excludes the Sigma_2 box corrections that have
// 'wrong' parity when calculating Sigma2 matrix
// elements. Note: If existing sk file already has
// these, they will be included [false]
sort_output;
// Sort output by energy? Default is to sort by J and Pi
// first. [false]
print_details;
// Condition to print details of each CI solution
// (otherwise just prints summary) [true]
parallel_ci;
// Run CI in parallel (solve each J/Pi in parallel).
// Faster, uses slightly more memory [true]
}
Functions and classes for Configuration Interaction calculations.
Definition CI_Integrals.cpp:11

Always check for up-to-date options from command line: $ ampsci -i CI See also run_CI, which this function calls

Parameters
inputInput block containing CI options.
wfFully initialised Wavefunction object supplying the orbital basis and radial grid.
Returns
Vector of PsiJPi, one entry per (J, parity) combination requested.

◆ run_CI()

PsiJPi CI::run_CI ( const std::vector< DiracSpinor > &  ci_sp_basis,
int  twoJ,
int  parity,
int  num_solutions,
std::optional< double >  all_below_cm,
const Coulomb::meTable< double > &  h1,
const Coulomb::QkTable qk,
const Coulomb::WkTable Bk,
const Coulomb::LkTable Sk,
bool  include_Sigma2,
bool  print_details,
std::ostream &  outstream = std::cout 
)

Constructs and solves the CI eigenvalue problem for a single J,pi.

Builds the CI+MBPT Hamiltonian matrix in the basis of two-electron configuration state functions (CSFs) with total angular momentum twoJ /2 and parity parity, then solves the eigenvalue problem to obtain CI energies and expansion coefficients.

The Hamiltonian includes:

  • one-body terms from h1
  • two-body Coulomb interaction via the \( Q^k \) integrals in qk
  • optionally, two-body Breit corrections via \( B^k \) integrals in Bk (used if Bk is non-empty)
  • optionally, two-body MBPT \( \Sigma_2 \) corrections via \( S^k \) integrals in Sk (used if include_Sigma2 is true, and Sk is non-empty)

The number of solutions returned is controlled by num_solutions and all_below: if all_below is set it takes precedence and all eigenstates with total energy below the threshold are found.

Parameters
ci_sp_basisSingle-particle basis states spanning the CI space.
twoJTwice the total angular momentum, 2J (must be a non-negative even integer for two-electron systems).
parityParity of the sector: +1 (even) or -1 (odd).
num_solutionsNumber of lowest eigenstates to find. Ignored if all_below_cm is set. Pass 0 to find all solutions.
all_below_cmIf set, find all eigenstates with total energy below this value (in cm^-1). Overrides num_solutions.
h1Table of one-body Hamiltonian matrix elements between single-particle basis states.
qkTable of two-body Coulomb \( Q^k \) integrals.
BkTable of two-body Breit \( B^k \) integrals. Ignored (treated as absent) if the table is empty.
SkTable of two-body MBPT \( \Sigma_2 \) ( \( S^k \)) integrals. Only used when include_Sigma2 is true.
include_Sigma2If true, add two-body MBPT corrections from Sk to the CI Hamiltonian.
print_detailsIf true, print a breakdown of the leading configurations for each solution. Leads to very large output if num_solutions is large
outstreamOutput stream for progress and results [default: stdout].
Returns
PsiJPi (PsiJPi) containing the CI eigenvalues and expansion coefficients for the requested solutions.

◆ form_CSFs()

std::vector< CSF2 > CI::form_CSFs ( int  twoJ,
int  parity,
const std::vector< DiracSpinor > &  cisp_basis 
)

Forms all two-electron CSFs with given total J and parity.

Iterates over all pairs of single-particle states in cisp_basis and retains those whose angular momenta can be coupled to total \( J = \) twoJ /2 and whose combined parity equals parity. Duplicate pairs are excluded by construction.

Parameters
twoJTwice the total angular momentum 2J.
parityTotal parity: +1 (even) or -1 (odd).
cisp_basisSingle-particle basis from which CSFs are constructed.
Returns
Sorted list of all valid two-electron CSFs for the given J and parity.