ampsci
c++ program for high-precision atomic structure calculations of single-valence systems
|
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 verbose=true) | |
Construct Feynman diagram. More... | |
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) const |
Calculates Green's function for kappa, and complex energy. | |
ComplexGMatrix | 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 ComplexGMatrix & | get_qk (int k) const |
Returns (reference to) q^k (radial) matrix. Note: includes drj. | |
const GMatrix & | get_Vx_kappa (int kappa) const |
Returns (ref to) radial exchange matrix Vx_kappa. Nb: includes dri*drj. | |
Feynman & | operator= (const Feynman &)=default |
Feynman (const Feynman &)=default | |
Class to construct Feynman diagrams, Green's functions and polarisation op.
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 | verbose = true |
||
) |
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?