ampsci
c++ program for high-precision atomic structure calculations of single-valence systems
HartreeFock.hpp
1 #pragma once
2 #include "Coulomb/YkTable.hpp"
3 #include "HF/Breit.hpp"
4 #include "Physics/Parametric_potentials.hpp"
5 #include "Physics/PhysConst_constants.hpp"
6 #include "Physics/RadPot.hpp"
7 #include <memory>
8 #include <string>
9 #include <vector>
10 class Wavefunction;
11 class DiracSpinor;
12 class Grid;
13 namespace MBPT {
14 class CorrelationPotential;
15 }
16 
18 namespace HF {
19 
20 //==============================================================================
23 struct EpsIts {
24  double eps{0.0};
25  int its{0};
26  std::string symbol{};
27  friend bool operator<(const EpsIts &l, const EpsIts &r) {
28  return l.eps < r.eps;
29  }
30 };
31 
32 //==============================================================================
34 
41 enum class Method { HartreeFock, ApproxHF, Hartree, KohnSham, Local };
43 Method parseMethod(const std::string &in_method);
44 std::string parseMethod(const Method &in_method);
45 
46 //==============================================================================
51 std::vector<double> vex_approx(const DiracSpinor &Fa,
52  const std::vector<DiracSpinor> &core,
53  int k_cut = 99, double lambda_cut = 0.003);
54 
59 DiracSpinor vexFa(const DiracSpinor &Fa, const std::vector<DiracSpinor> &core,
60  int k_cut = 99);
61 
62 //==============================================================================
63 //==============================================================================
64 //==============================================================================
65 
70 class HartreeFock {
71 
72 private:
73  std::shared_ptr<const Grid> m_rgrid;
74  std::vector<DiracSpinor> m_core;
75  std::vector<double> m_vnuc;
76  std::optional<QED::RadPot> m_vrad;
77  std::optional<HF::Breit> m_VBr;
78  double m_alpha;
79  Method m_method;
80  double m_eps_HF;
81  std::vector<double> m_vdir;
82  Coulomb::YkTable m_Yab;
83  int m_max_hf_its = 128;
84 
85 public:
87 
109  HartreeFock(std::shared_ptr<const Grid> rgrid, std::vector<double> vnuc,
110  std::vector<DiracSpinor> core,
111  std::optional<QED::RadPot> vrad = std::nullopt,
112  double m_alpha = PhysConst::alpha,
113  Method method = Method::HartreeFock, double x_Breit = 0.0,
114  double eps_HF = 0.0,
115  Parametric::Type potential = Parametric::Type::Green,
116  double H_g = 0.0, double d_t = 0.0);
117 
119  EpsIts solve_core(bool print = true);
120 
126  void
127  solve_valence(std::vector<DiracSpinor> *valence, bool print = true,
128  const MBPT::CorrelationPotential *const Sigma = nullptr) const;
129 
132  const MBPT::CorrelationPotential *const Sigma = nullptr,
133  std::optional<double> eta = std::nullopt,
134  std::optional<int> prev_its = std::nullopt) const;
135 
137  double calculateCoreEnergy() const;
138 
140  DiracSpinor vexFa(const DiracSpinor &Fa) const {
141  // calls static version with HF core
142  return ::HF::vexFa(Fa, m_core, 99);
143  }
144 
146  DiracSpinor VBr(const DiracSpinor &Fv) const;
147 
148  //---------------------------
149 
151  const Grid &grid() const { return *m_rgrid; };
154  std::shared_ptr<const Grid> grid_sptr() const { return m_rgrid; };
155 
157  const std::vector<double> &vdir() const { return m_vdir; }
158  std::vector<double> &vdir() { return m_vdir; }
159 
161  const std::vector<double> &vnuc() const { return m_vnuc; }
162  std::vector<double> &vnuc() { return m_vnuc; }
163 
165  std::vector<double> Hrad_el(int l = 0) const;
168  std::vector<double> Hmag(int l = 0) const;
169 
171  std::vector<double> vlocal(int l = 0) const;
172 
174  Method method() const { return m_method; }
175 
177  double zion() const;
178 
180  bool excludeExchangeQ() const {
181  return m_core.empty() ||
182  !(m_method == Method::HartreeFock || m_method == Method::ApproxHF);
183  }
184 
186  const std::vector<DiracSpinor> &core() const { return m_core; }
187 
189  double alpha() const { return m_alpha; }
190 
192  // not core, for testing)
193  void set_Vrad(QED::RadPot in_vrad) { m_vrad = std::move(in_vrad); } // XXX
195  const QED::RadPot *Vrad() const { return m_vrad ? &*m_vrad : nullptr; }
196  QED::RadPot *Vrad() { return m_vrad ? &*m_vrad : nullptr; }
197 
199  const HF::Breit *vBreit() const { return m_VBr ? &*m_VBr : nullptr; }
201  double x_Breit() const { return m_VBr ? m_VBr->scale_factor() : 0.0; }
202 
204  int num_core_electrons() const;
205 
206 private:
207  // Solve Dirac equation for core states (just once) using existing vdir
208  // (usually, vdir set to parametric potential beforehand)
209  EpsIts solve_initial_core(const double eps);
210  // Solve equations self-consistantly for core, using local method (either
211  // core-Hartree or Kohn-Sham)
212  EpsIts selfcon_local_core(const double eps_target_HF);
213  // Solve HF equations self-consistantly for core, using approximate HF method
214  EpsIts hf_approx_core(const double eps_target_HF);
215  // Solve HF equations self-consistantly for core, using Hartree-Fock method
216  EpsIts hartree_fock_core();
217  // Solves HF equation for valence state, assuming local potential (Hartree,
218  // Local, Kohn-Sham or approxHF)
219  EpsIts local_valence(DiracSpinor &Fa) const;
220 
221  /*
222  // same as hf_valence, but uses Green method
223  EpsIts hf_valence_Green(
224  DiracSpinor &Fv,
225  const MBPT::CorrelationPotential *const Sigma = nullptr) const;
226  */
227 
228  // Solves HF equation for given state, using non-local Green's method for
229  // inhomogeneous ODE (used for hartree_fock_core()).
230  // Solve Dirac Equation (Eigenvalue):
231  // (H0 + Vl + Vx)Fa = 0
232  // (H0 + Vl)Fa = -VxFa
233  // Vl is local (e.g., Vnuc + fVdir), Vx is non-local (e.g., (1-f)Vdir + Vex)
234  // where v0 = (1-f)Vdir [f=1 for valence states!, so v0 may be empty]
235  // Vx also includes Breit, and Sigma
236  // Small energy adjustmenets (and wfs), solve:
237  // (Hl - e) dF = de * F -VxFa
238  // e -> e+de, F->F+dF
239  // Core is input so can call in a thread-safe way! (with a 'old_core' copy)
240  // Only used in dE from dF
241  void hf_orbital_green(
242  DiracSpinor &Fa, double en, const std::vector<double> &vl,
243  const std::vector<double> &H_mag, const DiracSpinor &VxF,
244  const std::vector<DiracSpinor> &static_core,
245  const std::vector<double> &dv0 = {}, const HF::Breit *const VBr = nullptr,
246  const MBPT::CorrelationPotential *const Sigma = nullptr) const;
247 
248  // Calc's Vex*Fa, for Fa in the core. Fa must be in the core
249  void vex_Fa_core(const DiracSpinor &Fa, DiracSpinor &vexFa) const;
250 
251  // Option to re-scale diract potential so that V(r)~-zion/r at large r
252  enum class ReScale { yes = true, no = false };
253  // Forms direct potential
254  void update_vdir(ReScale re_scale = ReScale::no);
255  // Adds the additional Kohn-Sham parts to Vdir
256  void add_KohnSham_vdir_addition();
257 
258  // Sets Vdir to be parametric potential. By default, Greens potential
259  void
260  set_parametric_potential(bool print = true,
261  Parametric::Type potential = Parametric::Type::Green,
262  double H_g = 0.0, double d_t = 0.0);
263 
264  // Forms approximate vex for all core states
265  void form_approx_vex_core(std::vector<std::vector<double>> &vex) const;
266  std::vector<std::vector<double>> form_approx_vex_core() const;
267  // Forms approximate vex for given core states
268  void form_approx_vex_core_a(const DiracSpinor &Fa,
269  std::vector<double> &vex_a) const;
270  std::vector<double> form_approx_vex_core_a(const DiracSpinor &Fa) const;
271 
272  // Energy guess for core states
273  double enGuessCore(int n, int ka) const;
274  // Energy guess for valence states
275  double enGuessVal(int n, int ka) const;
276 };
277 
278 } // namespace HF
Calculates + stores Hartree Y functions + Angular (w/ look-up), taking advantage of symmetry.
Definition: YkTable.hpp:26
Stores radial Dirac spinor: F_nk = (f, g)
Definition: DiracSpinor.hpp:41
Holds grid, including type + Jacobian (dr/du)
Definition: Grid.hpp:31
Breit (Hartree-Fock Breit) interaction potential.
Definition: Breit.hpp:52
Solves relativistic Hartree-Fock equations for core and valence. Optionally includes Breit and QED ef...
Definition: HartreeFock.hpp:70
double x_Breit() const
Breit scale factor (usualy 0 or 1)
Definition: HartreeFock.hpp:201
std::vector< double > Hrad_el(int l=0) const
Electric part of radiative potential.
Definition: HartreeFock.cpp:1078
double calculateCoreEnergy() const
Calculates the HF core energy (not including Breit?)
Definition: HartreeFock.cpp:653
const std::vector< double > & vdir() const
Returns reference to Vdir (direct HF potential)
Definition: HartreeFock.hpp:157
std::shared_ptr< const Grid > grid_sptr() const
Resturns copy of shared_ptr to grid [shared resource] - used when we want to construct a new object t...
Definition: HartreeFock.hpp:154
const HF::Breit * vBreit() const
pointer to Breit - may be nullptr if no breit
Definition: HartreeFock.hpp:199
double zion() const
Effective charge at large Z : zion = Z - num_core_electrons.
Definition: HartreeFock.cpp:1162
const std::vector< double > & vnuc() const
Returns reference to Vnuc (nuclear potential)
Definition: HartreeFock.hpp:161
void set_Vrad(QED::RadPot in_vrad)
Update the Vrad used inside HF (only used if we want QED into valence but.
Definition: HartreeFock.hpp:193
double alpha() const
Value of fine-structure constant used.
Definition: HartreeFock.hpp:189
std::vector< double > Hmag(int l=0) const
Magnetic (off-diagonal) part of radiative potential. Doesn't currently depend on l.
Definition: HartreeFock.cpp:1081
HartreeFock(std::shared_ptr< const Grid > rgrid, std::vector< double > vnuc, std::vector< DiracSpinor > core, std::optional< QED::RadPot > vrad=std::nullopt, double m_alpha=PhysConst::alpha, Method method=Method::HartreeFock, double x_Breit=0.0, double eps_HF=0.0, Parametric::Type potential=Parametric::Type::Green, double H_g=0.0, double d_t=0.0)
Method is enum class, eps_HF is convergence goal.
Definition: HartreeFock.cpp:56
Method method() const
Which method used to solve HF.
Definition: HartreeFock.hpp:174
const Grid & grid() const
Resturns a const reference to the radial grid.
Definition: HartreeFock.hpp:151
DiracSpinor vexFa(const DiracSpinor &Fa) const
Calculates exchange term Vex*Fa.
Definition: HartreeFock.hpp:140
const QED::RadPot * Vrad() const
Get (const) ptr to Vrad - may be null.
Definition: HartreeFock.hpp:195
EpsIts hf_valence(DiracSpinor &Fv, const MBPT::CorrelationPotential *const Sigma=nullptr, std::optional< double > eta=std::nullopt, std::optional< int > prev_its=std::nullopt) const
Solves HF equation (+ Sigma) for single valence state.
Definition: HartreeFock.cpp:534
bool excludeExchangeQ() const
Returns true if exchange not included.
Definition: HartreeFock.hpp:180
int num_core_electrons() const
Number of electrons in the core.
Definition: HartreeFock.cpp:710
const std::vector< DiracSpinor > & core() const
vector of core orbitals
Definition: HartreeFock.hpp:186
void solve_valence(std::vector< DiracSpinor > *valence, bool print=true, const MBPT::CorrelationPotential *const Sigma=nullptr) const
Solves HF for given valence list. They need not already be solutions.
Definition: HartreeFock.cpp:144
DiracSpinor VBr(const DiracSpinor &Fv) const
Breit interaction V_Br*Fa.
Definition: HartreeFock.cpp:702
EpsIts solve_core(bool print=true)
Solves HF equations self-consitantly for core orbs. Returns epsilon.
Definition: HartreeFock.cpp:77
std::vector< double > vlocal(int l=0) const
vlocal = vnuc + vrad_el + vdir
Definition: HartreeFock.cpp:696
Class holds Flambaum-Ginges QED Radiative Potential.
Definition: RadPot.hpp:15
Stores Wavefunction (set of valence orbitals, grid, HF etc.)
Definition: Wavefunction.hpp:36
Functions and classes for Hartree-Fock.
Definition: CI_Integrals.hpp:12
std::vector< double > vex_approx(const DiracSpinor &Fa, const std::vector< DiracSpinor > &core, int k_cut, const double lambda_cut)
Forms approx (localised) exchange potential, from scratch.
Definition: HartreeFock.cpp:947
Method
Methods available for self-consistant field model.
Definition: HartreeFock.hpp:41
DiracSpinor vexFa(const DiracSpinor &Fa, const std::vector< DiracSpinor > &core, int k_cut)
Calculates V_exch * Fa, for any orbital Fa (calculates Coulomb integral from scratch).
Definition: HartreeFock.cpp:1043
Method parseMethod(const std::string &in_method)
Convers string (name) of method (HartreeFock, Hartree etc.) to enum.
Definition: HartreeFock.cpp:25
Many-body perturbation theory.
Definition: CI_Integrals.hpp:9
constexpr double alpha
fine-structure constant
Definition: PhysConst_constants.hpp:13
Small struct to store: {eps, its, symbol}. eps=convergence; its=iterations; symbol=which state....
Definition: HartreeFock.hpp:23