ampsci
c++ program for high-precision atomic structure calculations of single-valence systems
Loading...
Searching...
No Matches
Feynman.hpp
1#pragma once
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"
7#include <memory>
8#include <vector>
9
10namespace MBPT {
11
12// FeynmanOptions class
13
16enum class HoleParticle { exclude, include, include_k0 };
17
19enum class Screening { exclude, include };
20
21struct FeynmanOptions {
22 Screening screening{Screening::exclude};
23 HoleParticle hole_particle{HoleParticle::exclude};
24 int max_l_internal{6};
25 double omre{-0.3};
26 double w0{0.05};
27 double w_ratio{1.5};
28 // int n_min_core{1};
29};
30
32enum class GreenStates { both, core, excited };
33
34//------------------------------------------------------------------------------
36class Feynman {
37
38 // Pointer to HF potential, and HF core
39 const HF::HartreeFock *m_HF;
40 // Pointer to shared radial grid (full grid)
41 std::shared_ptr<const Grid> m_grid;
42 // Parameters of the sub-grid: initial/final points, stride
43 std::size_t m_i0, m_stride, m_subgrid_points;
44
45 // maximum kappa_index appearing in the core
46 int m_max_ki_core;
47 // maximum kappa index to include for internal lines
48 int m_max_ki;
49 // maximum multipolarity, k
50 int m_max_k;
51 // Lowest n to polarise in polarisation operator
52 int m_min_core_n;
53 // Include relativistic components of the correlation potential?
54 bool m_include_G;
55 // real part of frequency for integration
56 double m_omre;
57 // Stores imaginary frequency grid (setup post-construction)
58 Grid m_wgrid;
59
60 // Option: include hole-particle interaction into polarisation operator
61 bool m_hole_particle;
62 // Option to include higher-k in hp interaction
63 bool m_include_higher_order_hp;
64 // Include screening correction to Coulomb line in Q*Pi*Q => Q*Pi*X*Q
65 bool m_screen_Coulomb;
66
67 // Coulomb operator; include right integration measure
68 // q_ij := (r>^k/r<^{k+1}) * dr_j
69 std::vector<ComplexRMatrix> m_qk{};
70
71 // Core projection operators (one for each core state)
72 std::vector<ComplexGMatrix> m_pa{};
73 // Hartree-Fock Exchange matrix (one for each kappa)
74 std::vector<GMatrix> m_Vx_kappa{};
75
76 // Effective spinless Q*Pi*Q operator: for each imaginary omega, and each k
78
79 // For now, just for testing: switch between complex green methods
80 // This is broken, not sure what changed
81 bool m_Complex_green_method = false;
82
83public:
92 Feynman(const HF::HartreeFock *vHF, std::size_t i0, std::size_t stride,
93 std::size_t size, const FeynmanOptions &options, int n_min_core,
94 bool include_G, bool verbose = true, const std::string &ident = "");
95
96 bool screening() const { return m_screen_Coulomb; }
97 bool hole_particle() const { return m_hole_particle; }
98
100 std::size_t stride() const { return m_stride; }
101
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); }
106 int lmax() const { return Angular::lFromIndex(m_max_ki); }
107
109 ComplexGMatrix green(int kappa, std::complex<double> en,
110 GreenStates states = GreenStates::both) const;
111
113 ComplexRMatrix polarisation_k(int k, std::complex<double> omega,
114 bool hole_particle) const;
115
117 GMatrix Sigma_direct(int kappa_v, double en_v,
118 std::optional<int> k = {}) const;
119
121 const ComplexRMatrix &get_qk(int k) const { return m_qk.at(std::size_t(k)); }
122
124 const GMatrix &get_Vx_kappa(int kappa) const {
125 return m_Vx_kappa.at(std::size_t(Angular::indexFromKappa(kappa)));
126 }
127
128private:
129 // forms Qk matrices, as well as dri, drj
130 void form_qk();
131 // Forms core projection operators
132 void form_pa();
133 // Forms HF exchange potential matrix
134 void form_vx();
135 // Sets up imaginary frequency grid
136 Grid form_w_grid(double w0, double wratio) const;
137 // Constructs the Q*Pi*Q Matrix along w grid, for each k
138 void form_qpiq();
139
140 bool readwrite_qpiq(IO::FRW::RoW rw, const std::string &fname);
141
142 // Screening factor X = [1 + i qk*pik]^-1
143 ComplexRMatrix X_screen(const ComplexRMatrix &pik,
144 const ComplexRMatrix &qdri) const;
145
146 // Forms single "Green's function" contribution f*|ket><bra|
147 // (f is usually 1/(e-en))
148 ComplexGMatrix green_single(const DiracSpinor &ket, const DiracSpinor &bra,
149 const std::complex<double> f) const;
150
151 // Forms full Hartree-Fock green's function. If Fc_hp is given,
152 // includes hole-particle contribution (Fc is polarised core state).
153 // Uses Dyson method to extend to complex energies.
154 ComplexGMatrix green_hf(int kappa, std::complex<double> en,
155 const DiracSpinor *Fc_hp = nullptr) const;
156
157 // Forms full Hartree-Fock green's function. If Fc_hp is given,
158 // includes hole-particle contribution (Fc is polarised core state).
159 // Uses "Complex Dirac" method for complex energies.
160 ComplexGMatrix green_hf_v2(int kappa, std::complex<double> en,
161 const DiracSpinor *Fc_hp = nullptr) const;
162
163 // Given Gr(wr), and wi (wr is real), returns Gr(wr + i*wi)
164 ComplexGMatrix green_to_complex(const ComplexGMatrix &Gr,
165 double om_imag) const;
166
167 // Green's function, only due to core states
168 ComplexGMatrix green_core(int kappa, std::complex<double> en) const;
169
170 // Greens function, only due to excited states (by orthog to core)
171 ComplexGMatrix green_excited(int kappa, std::complex<double> en,
172 const DiracSpinor *Fc_hp = nullptr) const;
173
174 // Forces Green's fn to be orthog. to core
175 ComplexGMatrix orthogonalise_wrt_core(const ComplexGMatrix &g_in,
176 int kappa) const;
177
178 // Given Dirac solutions regular at 0 (x0) and infinity (xI), forms "local"
179 // Green's function, with all four spinor components
180 GMatrix construct_green_g0(const DiracSpinor &x0, const DiracSpinor &xI,
181 const double w) const;
182
183 // Given Dirac solutions regular at 0 (x0 + i*Ix0) and infinity (xI + i*IxI),
184 // forms "local" Green's function (complex).
185 ComplexGMatrix construct_green_g0(const DiracSpinor &x0,
186 const DiracSpinor &Ix0,
187 const DiracSpinor &xI,
188 const DiracSpinor &IxI,
189 const std::complex<double> w) const;
190
191 // Calculates hole-particle potential (1-P)Vhp(1-P)
192 GMatrix calculate_Vhp(int kappa, const DiracSpinor &Fc) const;
193
194 // [[nodiscard]] const ComplexGMatrix &get_dri() const { return m_dri; }
195 // [[nodiscard]] const ComplexGMatrix &get_drj() const { return m_drj; }
196
197public:
198 Feynman &operator=(const Feynman &) = default;
199 Feynman(const Feynman &) = default;
200 ~Feynman() = default;
201};
202
203} // namespace MBPT
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