2#include "HF/HartreeFock.hpp"
3#include "MBPT/RadialMatrix.hpp"
4#include "MBPT/SpinorMatrix.hpp"
5#include "Maths/Grid.hpp"
6#include "Wavefunction/DiracSpinor.hpp"
21struct FeynmanOptions {
24 int max_l_internal{6};
41 std::shared_ptr<const Grid> m_grid;
43 std::size_t m_i0, m_stride, m_subgrid_points;
63 bool m_include_higher_order_hp;
65 bool m_screen_Coulomb;
69 std::vector<ComplexRMatrix> m_qk{};
72 std::vector<ComplexGMatrix> m_pa{};
74 std::vector<GMatrix> m_Vx_kappa{};
81 bool m_Complex_green_method =
false;
93 std::size_t size,
const FeynmanOptions &options,
int n_min_core,
94 bool include_G,
bool verbose =
true,
const std::string &ident =
"");
96 bool screening()
const {
return m_screen_Coulomb; }
97 bool hole_particle()
const {
return m_hole_particle; }
100 std::size_t
stride()
const {
return m_stride; }
102 int n_min()
const {
return m_min_core_n; }
103 double omre()
const {
return m_omre; }
104 double w0()
const {
return m_wgrid.
r0(); }
105 double wratio()
const {
return m_wgrid.
r(1) / m_wgrid.
r(0); }
109 ComplexGMatrix
green(
int kappa, std::complex<double> en,
114 bool hole_particle)
const;
118 std::optional<int> k = {})
const;
136 Grid form_w_grid(
double w0,
double wratio)
const;
140 bool readwrite_qpiq(IO::FRW::RoW rw,
const std::string &fname);
149 const std::complex<double> f)
const;
165 double om_imag)
const;
168 ComplexGMatrix green_core(
int kappa, std::complex<double> en)
const;
181 const double w)
const;
189 const std::complex<double> w)
const;
Stores radial Dirac spinor: F_nk = (f, g)
Definition DiracSpinor.hpp:41
Holds grid, including type + Jacobian (dr/du)
Definition Grid.hpp:31
const std::vector< double > & r() const
Grid points, r.
Definition Grid.hpp:75
auto r0() const
Minium (first) grid point.
Definition Grid.hpp:58
Solves relativistic Hartree-Fock equations for core and valence. Optionally includes Breit and QED ef...
Definition HartreeFock.hpp:70
Matrix class; row-major.
Definition Matrix.hpp:35
Class to construct Feynman diagrams, Green's functions and polarisation op.
Definition Feynman.hpp:36
std::size_t stride() const
Returns stride used for sub-grid.
Definition Feynman.hpp:100
ComplexGMatrix green(int kappa, std::complex< double > en, GreenStates states=GreenStates::both) const
Calculates Green's function for kappa, and complex energy.
Definition Feynman.cpp:227
const ComplexRMatrix & get_qk(int k) const
Returns (reference to) q^k (radial) matrix. Note: includes drj? No?
Definition Feynman.hpp:121
ComplexRMatrix polarisation_k(int k, std::complex< double > omega, bool hole_particle) const
Polarisation operator pi^k(w), for multipolarity k.
Definition Feynman.cpp:539
GMatrix Sigma_direct(int kappa_v, double en_v, std::optional< int > k={}) const
Calculate Direct part of correlation potential.
Definition Feynman.cpp:774
const GMatrix & get_Vx_kappa(int kappa) const
Returns (ref to) radial exchange matrix Vx_kappa. Nb: includes dri*drj.
Definition Feynman.hpp:124
Definition RadialMatrix.hpp:27
constexpr int lFromIndex(int i)
returns l, given kappa_index
Definition Wigner369j.hpp:90
constexpr int indexFromKappa(int ka)
returns kappa_index, given kappa
Definition Wigner369j.hpp:80
Many-body perturbation theory.
Definition CI_Integrals.hpp:9
GreenStates
Which states to include in Green's function:
Definition Feynman.hpp:32
Screening
Options for including Screening.
Definition Feynman.hpp:19
HoleParticle
Options for including hole-particle interaction. include mean all k; include_k0 means k=0 term only.
Definition Feynman.hpp:16