ampsci
High-precision calculations for one- and two-valence atomic systems
Loading...
Searching...
No Matches
MBPT::Feynman

Class to construct Feynman diagrams, Green's functions and polarisation op. More...

#include <Feynman.hpp>

Public Member Functions

 Feynman (const HF::HartreeFock *vHF, std::size_t i0, std::size_t stride, std::size_t size, const FeynmanOptions &options, int n_min_core, bool include_G, bool verbose=true, const std::string &ident="")
 Construct Feynman diagram.
 
bool screening () const
 
bool hole_particle () const
 
std::size_t stride () const
 Returns stride used for sub-grid.
 
int n_min () const
 
double omre () const
 
double w0 () const
 
double wratio () const
 
int lmax () const
 
ComplexGMatrix green (int kappa, std::complex< double > en, GreenStates states=GreenStates::both) const
 Calculates Green's function for kappa, and complex energy.
 
ComplexRMatrix polarisation_k (int k, std::complex< double > omega, bool hole_particle) const
 Polarisation operator pi^k(w), for multipolarity k.
 
GMatrix Sigma_direct (int kappa_v, double en_v, std::optional< int > k={}) const
 Calculate Direct part of correlation potential.
 
const ComplexRMatrixget_qk (int k) const
 Returns (reference to) q^k (radial) matrix. Note: includes drj? No?
 
const GMatrixget_Vx_kappa (int kappa) const
 Returns (ref to) radial exchange matrix Vx_kappa. Nb: includes dri*drj.
 
Feynmanoperator= (const Feynman &)=default
 
 Feynman (const Feynman &)=default
 

Detailed Description

Class to construct Feynman diagrams, Green's functions and polarisation op.

Constructor & Destructor Documentation

◆ Feynman()

MBPT::Feynman::Feynman ( const HF::HartreeFock vHF,
std::size_t  i0,
std::size_t  stride,
std::size_t  size,
const FeynmanOptions &  options,
int  n_min_core,
bool  include_G,
bool  verbose = true,
const std::string &  ident = "" 
)

Construct Feynman diagram.

rgrid_params={r0, rmax, stride}; omre is real part of frequency, w_0 is initial point along imaginary freq. axis; w_ratio is ratio used for logarithic omega grid (integration); scr_option and hp_option are screening and hole-particle interactions; max_l is maximum l to include for internal lines (Green's functions); n_min_core is minimum n to include in polarisation loop ** ** Currently have issue: polarising deep n leads to failure?

Member Function Documentation

◆ stride()

std::size_t MBPT::Feynman::stride ( ) const
inline

Returns stride used for sub-grid.

◆ green()

ComplexGMatrix MBPT::Feynman::green ( int  kappa,
std::complex< double >  en,
GreenStates  states = GreenStates::both 
) const

Calculates Green's function for kappa, and complex energy.

◆ polarisation_k()

ComplexRMatrix MBPT::Feynman::polarisation_k ( int  k,
std::complex< double >  omega,
bool  hole_particle 
) const

Polarisation operator pi^k(w), for multipolarity k.

◆ Sigma_direct()

GMatrix MBPT::Feynman::Sigma_direct ( int  kappa_v,
double  en_v,
std::optional< int >  k = {} 
) const

Calculate Direct part of correlation potential.

◆ get_qk()

const ComplexRMatrix & MBPT::Feynman::get_qk ( int  k) const
inline

Returns (reference to) q^k (radial) matrix. Note: includes drj? No?

◆ get_Vx_kappa()

const GMatrix & MBPT::Feynman::get_Vx_kappa ( int  kappa) const
inline

Returns (ref to) radial exchange matrix Vx_kappa. Nb: includes dri*drj.


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