ampsci
c++ program for high-precision atomic structure calculations of single-valence systems
Loading...
Searching...
No Matches
CorrelationPotential.hpp
1#pragma once
2#include "Angular/CkTable.hpp"
3#include "Angular/SixJTable.hpp"
4#include "Coulomb/YkTable.hpp"
5#include "IO/FRW_fileReadWrite.hpp"
6#include "MBPT/Feynman.hpp"
7#include "MBPT/Goldstone.hpp"
8#include "MBPT/RDMatrix.hpp"
9#include "Wavefunction/DiracSpinor.hpp"
10#include <optional>
11#include <vector>
12
13namespace MBPT {
14
15struct SigmaData {
16 int kappa;
17 double en;
18 RDMatrix<double> Sigma;
19 int n{0};
20 double lambda{1.0};
21};
22
23enum class SigmaMethod { Goldstone, Feynman };
24
25struct rgrid_params {
26 double r0{1.0e-4};
27 double rmax{30.0};
28 std::size_t stride{4};
29};
30
31// method:
32/*
33
34 - Goldstone/Feynman
35 - Screening
36 - hp
37 - Include g?
38
39 - Effective screening: input, or calculate?
40 - If given, use, else: calculate
41
42
43*/
44
45//==============================================================================
46
47class CorrelationPotential {
48 const HF::HartreeFock *m_HF;
49 std::vector<DiracSpinor> m_basis; // so we can delay goldstone construction
50 std::vector<SigmaData> m_Sigmas{};
51 double m_r0, m_rmax;
52 std::size_t m_stride;
53 std::size_t m_i0, m_size; // need?
54
55 SigmaMethod m_method;
56 int m_n_min_core;
57 int m_n_min_core_F;
58 bool m_includeG;
59 bool m_includeBreit;
60 int m_n_max_breit;
61
62 std::optional<Goldstone> m_Gold{};
63
64 FeynmanOptions m_Foptions;
65 bool m_calculate_fk; // if not, need fk and etak
66 std::vector<double> m_fk;
67 std::vector<double> m_etak;
68
69 std::optional<Feynman> m_Fy{};
70
71 // These are only for calculating fk and eta
72 std::optional<Feynman> m_Fy0{};
73 std::optional<Feynman> m_FyX{};
74 std::optional<Feynman> m_FyH{};
75
76public:
77 CorrelationPotential(const std::string &fname, const HF::HartreeFock *vHF,
78 const std::vector<DiracSpinor> &basis, double r0,
79 double rmax, std::size_t stride, int n_min_core,
80 SigmaMethod method, bool include_g = false,
81 bool include_Breit = false, int n_max_breit = 0,
82 const FeynmanOptions &Foptions = {},
83 bool calculate_fk = true,
84 const std::vector<double> &fk = {},
85 const std::vector<double> &etak = {},
86 std::optional<int> n_min_core_F = {});
87
88 // // not thread safe!
89 // void formSigma(int kappa, double en, int n = 0) {}
90 // not thread safe!
91 void formSigma(int kappa, double ev, int n, const DiracSpinor *Fv = nullptr);
92
93 bool empty() const { return m_Sigmas.empty(); }
94
95 const GMatrix *getSigma(int kappa, int n = 0) const;
96
97 double getLambda(int kappa, int n = 0) const;
98
99 void clear() { m_Sigmas.clear(); }
100
103 DiracSpinor SigmaFv(const DiracSpinor &Fv) const;
104 DiracSpinor operator()(const DiracSpinor &Fv) const { return SigmaFv(Fv); }
105
107 void scale_Sigma(const std::vector<double> &lambdas);
108
109 // if n=0, scales _all_
110 void scale_Sigma(double lambda, int kappa, int n = 0);
111
113 void print_scaling() const;
114
116 void print_info() const;
117
119 void print_subGrid() const;
120
121 void write(const std::string &fname) { read_write(fname, IO::FRW::write); }
122
123private:
124 bool read_write(const std::string &fname, IO::FRW::RoW rw);
125 void setup_Feynman();
126 std::vector<double> calculate_fk(double ev, const DiracSpinor &v) const;
127 std::vector<double> calculate_etak(double ev, const DiracSpinor &v) const;
128 const SigmaData *get(int kappa, int n = 0) const;
129
130 GMatrix formSigma_F(int kappa, double ev, const DiracSpinor *Fv = nullptr);
131 GMatrix formSigma_G(int kappa, double ev, const DiracSpinor *Fv = nullptr);
132
133public:
134 CorrelationPotential &operator=(const CorrelationPotential &) = default;
135 CorrelationPotential(const CorrelationPotential &) = default;
136 ~CorrelationPotential() = default;
137
138 //
139};
140
141} // namespace MBPT
Stores radial Dirac spinor: F_nk = (f, g)
Definition DiracSpinor.hpp:41
Solves relativistic Hartree-Fock equations for core and valence. Optionally includes Breit and QED ef...
Definition HartreeFock.hpp:70
Many-body perturbation theory.
Definition CI_Integrals.hpp:9