2#include "Coulomb/YkTable.hpp"
4#include "Physics/PhysConst_constants.hpp"
5#include "Potentials/Parametric_potentials.hpp"
6#include "Potentials/RadPot.hpp"
15class CorrelationPotential;
28 friend bool operator<(
const EpsIts &l,
const EpsIts &r) {
29 return std::abs(l.eps) < std::abs(r.eps);
46std::string parseMethod_short(
const Method &in_method);
54 const std::vector<DiracSpinor> &core,
55 int k_cut = 99,
double lambda_cut = 0.003);
75 std::shared_ptr<const Grid> m_rgrid;
76 std::vector<DiracSpinor> m_core;
77 std::vector<double> m_vnuc;
78 std::optional<QED::RadPot> m_vrad;
79 std::optional<HF::Breit> m_VBr;
83 std::vector<double> m_vdir;
85 int m_max_hf_its = 128;
111 HartreeFock(std::shared_ptr<const Grid> rgrid, std::vector<double>
vnuc,
112 std::vector<DiracSpinor>
core,
113 std::optional<QED::RadPot> vrad = std::nullopt,
116 std::optional<Breit::Params> breit_params = std::nullopt,
118 Parametric::Type potential = Parametric::Type::Green,
119 double H_g = 0.0,
double d_t = 0.0);
130 solve_valence(std::vector<DiracSpinor> *valence,
bool print =
true,
131 const MBPT::CorrelationPotential *
const Sigma =
nullptr)
const;
135 const MBPT::CorrelationPotential *
const Sigma =
nullptr,
136 std::optional<double> eta = std::nullopt,
137 std::optional<int> prev_its = std::nullopt)
const;
141 const MBPT::CorrelationPotential *
const Sigma)
const;
149 return ::HF::vexFa(Fa, m_core, 99);
161 std::shared_ptr<const Grid>
grid_sptr()
const {
return m_rgrid; };
164 const std::vector<double> &
vdir()
const {
return m_vdir; }
165 std::vector<double> &
vdir() {
return m_vdir; }
168 const std::vector<double> &
vnuc()
const {
return m_vnuc; }
169 std::vector<double> &
vnuc() {
return m_vnuc; }
172 std::vector<double>
Hrad_el(
int l = 0)
const;
175 std::vector<double>
Hmag(
int l = 0)
const;
178 std::vector<double>
vlocal(
int l = 0)
const;
188 return m_core.empty() ||
189 !(m_method == Method::HartreeFock || m_method == Method::ApproxHF);
193 const std::vector<DiracSpinor> &
core()
const {
return m_core; }
196 double alpha()
const {
return m_alpha; }
208 double x_Breit()
const {
return m_VBr ? m_VBr->scale_factor() : 0.0; }
216 EpsIts solve_initial_core(
const double eps);
219 EpsIts selfcon_local_core(
const double eps_target_HF);
221 EpsIts hf_approx_core(
const double eps_target_HF);
223 EpsIts hartree_fock_core();
248 void hf_orbital_green(
249 DiracSpinor &Fa,
double en,
const std::vector<double> &vl,
250 const std::vector<double> &H_mag,
const DiracSpinor &VxF,
251 const std::vector<DiracSpinor> &static_core,
252 const std::vector<double> &dv0 = {},
const HF::Breit *
const VBr =
nullptr,
253 const MBPT::CorrelationPotential *
const Sigma =
nullptr)
const;
259 enum class ReScale { yes =
true, no =
false };
261 void update_vdir(ReScale re_scale = ReScale::no);
263 void add_KohnSham_vdir_addition();
267 set_parametric_potential(
bool print =
true,
268 Parametric::Type potential = Parametric::Type::Green,
269 double H_g = 0.0,
double d_t = 0.0);
272 void form_approx_vex_core(std::vector<std::vector<double>> &vex)
const;
273 std::vector<std::vector<double>> form_approx_vex_core()
const;
276 std::vector<double> &vex_a)
const;
277 std::vector<double> form_approx_vex_core_a(
const DiracSpinor &Fa)
const;
280 double enGuessCore(
int n,
int ka)
const;
282 double enGuessVal(
int n,
int ka)
const;
Calculates + stores Hartree Y functions + Angular (w/ look-up), taking advantage of symmetry.
Definition YkTable.hpp:31
Stores radial Dirac spinor: F_nk = (f, g)
Definition DiracSpinor.hpp:42
Non-uniform radial grid with Jacobian, suitable for atomic structure calculations.
Definition Grid.hpp:85
Breit potentials for one- (Hartree-Fock Breit) and two-body Breit integrals.
Definition Breit.hpp:87
Solves relativistic Hartree-Fock equations for core and valence. Optionally includes Breit and QED ef...
Definition HartreeFock.hpp:72
const std::vector< double > & vnuc() const
Returns reference to Vnuc (nuclear potential)
Definition HartreeFock.hpp:168
const Grid & grid() const
Resturns a const reference to the radial grid.
Definition HartreeFock.hpp:158
double x_Breit() const
Breit scale factor (usualy 0 or 1)
Definition HartreeFock.hpp:208
std::vector< double > Hrad_el(int l=0) const
Electric part of radiative potential.
Definition HartreeFock.cpp:1104
double calculateCoreEnergy() const
Calculates the HF core energy (not including Breit?)
Definition HartreeFock.cpp:679
double zion() const
Effective charge at large Z : zion = Z - num_core_electrons.
Definition HartreeFock.cpp:1189
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:200
EpsIts hf_valence_Green(DiracSpinor &Fa, const MBPT::CorrelationPotential *const Sigma) const
Solves HF equation (+ Sigma) for single valence state, alternative method.
Definition HartreeFock.cpp:626
double alpha() const
Value of fine-structure constant used.
Definition HartreeFock.hpp:196
std::vector< double > Hmag(int l=0) const
Magnetic (off-diagonal) part of radiative potential. Doesn't currently depend on l.
Definition HartreeFock.cpp:1107
Method method() const
Which method used to solve HF.
Definition HartreeFock.hpp:181
DiracSpinor vexFa(const DiracSpinor &Fa) const
Calculates exchange term Vex*Fa.
Definition HartreeFock.hpp:147
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:551
const HF::Breit * vBreit() const
pointer to Breit - may be nullptr if no breit
Definition HartreeFock.hpp:206
const std::vector< double > & vdir() const
Returns reference to Vdir (direct HF potential)
Definition HartreeFock.hpp:164
int num_core_electrons() const
Number of electrons in the core.
Definition HartreeFock.cpp:736
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:158
DiracSpinor VBr(const DiracSpinor &Fv) const
Breit interaction V_Br*Fa.
Definition HartreeFock.cpp:728
bool is_localQ() const
Returns true if exchange not included.
Definition HartreeFock.hpp:187
EpsIts solve_core(bool print=true)
Solves HF equations self-consitantly for core orbs. Returns epsilon.
Definition HartreeFock.cpp:91
const std::vector< DiracSpinor > & core() const
vector of core orbitals
Definition HartreeFock.hpp:193
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:161
const QED::RadPot * Vrad() const
Get (const) ptr to Vrad - may be null.
Definition HartreeFock.hpp:202
std::vector< double > vlocal(int l=0) const
vlocal = vnuc + vrad_el + vdir
Definition HartreeFock.cpp:722
Constructs and stores the Flambaum-Ginges QED Radiative Potential.
Definition RadPot.hpp:16
Stores Wavefunction (set of valence orbitals, grid, HF etc.)
Definition Wavefunction.hpp:37
Functions and classes for Hartree-Fock.
Definition CI_Integrals.hpp:13
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:973
Method
Methods available for self-consistant field model.
Definition HartreeFock.hpp:42
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:1069
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:10
constexpr double alpha
Fine-structure constant: alpha = 1/137.035 999 177(21) [CODATA 2022].
Definition PhysConst_constants.hpp:24
Small struct to store: {eps, its, symbol}. eps=convergence; its=iterations; symbol=which state....
Definition HartreeFock.hpp:24