votca 2026-dev
Loading...
Searching...
No Matches
gw_uks.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2026 The VOTCA Development Team
3 * (http://www.votca.org)
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License")
6 */
7#pragma once
8#ifndef VOTCA_XTP_GW_UKS_H
9#define VOTCA_XTP_GW_UKS_H
10
11#include <memory>
12#include <unordered_set>
13
14#include "logger.h"
15#include "orbitals.h"
16#include "qp_solver_utils.h"
17#include "rpa_uks.h"
18#include "sigma_base_uks.h"
19#include "threecenter.h"
20#include "votca/xtp/ppm.h"
21
22namespace votca {
23namespace xtp {
24
25class GW_UKS {
26
31
32 public:
33 struct options {
40 double eta;
41 double g_sc_limit;
45 double shift = 0;
46 double ScaHFX = 0.0;
47 std::string sigma_integration;
49 std::string qp_solver;
50 double qp_solver_alpha = 0.75;
51 // Legacy aliases from the original grid-search implementation.
52 // They used to define both the total QP window and the scan resolution at
53 // the same time. The new implementation keeps them only so that old XML
54 // files and existing tests still work; configuration code maps them onto
55 // the decoupled controls below.
57 double qp_grid_spacing = 0.0;
58
59 // Decoupled QP search controls:
60 // - qp_full_window_half_width defines the outer full search interval
61 // - qp_dense_spacing defines the dense scan used for robust sign changes
62 // - qp_adaptive_shell_width / count define the outward shell search
63 // All values are normalized in configure(); negative values mean "unset"
64 // and trigger either legacy mapping or robust defaults.
65 double qp_full_window_half_width = -1.0; // Ha; <=0 means "unset"
66 double qp_dense_spacing = -1.0; // Ha; <=0 means "unset"
67 double qp_adaptive_shell_width = -1.0; // Ha; <=0 means "unset"
68 Index qp_adaptive_shell_count = 0; // 0 means "use shell width"
71 std::string quadrature_scheme;
73 double alpha;
74 bool qp_restrict_search = true;
75 double qp_zero_margin = 1e-6;
76 double qp_virtual_min_energy = -0.1;
77 std::string qp_root_finder = "bisection";
78 std::string qp_grid_search_mode = "adaptive_with_dense_fallback";
79 };
80
82 const Eigen::MatrixXd& vxc_alpha, const Eigen::MatrixXd& vxc_beta,
83 const Eigen::VectorXd& dft_energies_alpha,
84 const Eigen::VectorXd& dft_energies_beta);
85
86 void configure(const options& opt);
88 void CalculateHQP();
89
90 Eigen::VectorXd getGWAResultsAlpha() const;
91 Eigen::VectorXd getGWAResultsBeta() const;
92 const Eigen::VectorXd& RPAInputEnergiesAlpha() const;
93 const Eigen::VectorXd& RPAInputEnergiesBeta() const;
94 Eigen::MatrixXd getHQPAlpha() const;
95 Eigen::MatrixXd getHQPBeta() const;
96 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> DiagonalizeQPHamiltonianAlpha()
97 const;
98 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> DiagonalizeQPHamiltonianBeta()
99 const;
100
101 private:
102 enum class Spin { Alpha, Beta };
103
104 class QPFunc {
105 public:
106 QPFunc(Index gw_level, const Sigma_base_UKS& sigma, double offset)
107 : gw_level_(gw_level), offset_(offset), sigma_c_func_(sigma) {}
108
109 std::pair<double, double> operator()(double frequency) const {
110 std::pair<double, double> result;
111 result.first = value(frequency, EvalStage::Other);
112 result.second = deriv(frequency);
113 return result;
114 }
115
116 double sigma(double frequency, EvalStage stage = EvalStage::Other) const {
117 const std::uint64_t key = FrequencyKey(frequency);
118
119 auto insert_result = seen_frequencies_.insert(key);
120 if (!insert_result.second) {
121 ++stats_.sigma_repeat_calls;
122 } else {
123 ++stats_.sigma_unique_frequencies;
124 }
125
126 CountSigmaStage(stage);
127 return sigma_c_func_.CalcCorrelationDiagElement(gw_level_, frequency);
128 }
129
130 double value(double frequency, EvalStage stage = EvalStage::Other) const {
131 return sigma(frequency, stage) + offset_ - frequency;
132 }
133
134 double deriv(double frequency) const {
135 ++stats_.deriv_calls;
136 return sigma_c_func_.CalcCorrelationDiagElementDerivative(gw_level_,
137 frequency) -
138 1.0;
139 }
140
141 const QPStats& GetStats() const { return stats_; }
142
143 private:
144 static std::uint64_t FrequencyKey(double x) {
145 std::uint64_t key = 0;
146 static_assert(sizeof(double) == sizeof(std::uint64_t),
147 "Unexpected double size");
148 std::memcpy(&key, &x, sizeof(double));
149 return key;
150 }
151
152 void CountSigmaStage(EvalStage stage) const {
153 switch (stage) {
154 case EvalStage::Scan:
155 ++stats_.sigma_scan_calls;
156 break;
157 case EvalStage::Refine:
158 ++stats_.sigma_refine_calls;
159 break;
160 case EvalStage::Derivative:
161 ++stats_.sigma_derivative_calls;
162 break;
163 case EvalStage::Other:
164 default:
165 ++stats_.sigma_other_calls;
166 break;
167 }
168 }
169
171 double offset_;
173
174 mutable std::unordered_set<std::uint64_t> seen_frequencies_;
176 };
177
179 const Eigen::VectorXd& DftEnergies(Spin spin) const;
180 const Eigen::MatrixXd& Vxc(Spin spin) const;
181 Eigen::MatrixXd& SigmaX(Spin spin);
182 Eigen::MatrixXd& SigmaC(Spin spin);
183 const Eigen::MatrixXd& SigmaX(Spin spin) const;
184 const Eigen::MatrixXd& SigmaC(Spin spin) const;
186 const Sigma_base_UKS& SigmaEvaluator(Spin spin) const;
187 Index Homo(Spin spin) const;
188 const char* SpinName(Spin spin) const;
189
190 Eigen::VectorXd ScissorShift_DFTlevel(const Eigen::VectorXd& dft_energies,
191 Index homo) const;
192 double CalcSpinHomoLumoShift(const Eigen::VectorXd& frequencies,
193 Spin spin) const;
194 void PrintGWA_Energies(Spin spin) const;
195 void PrintQP_Energies(Spin spin,
196 const Eigen::VectorXd& qp_diag_energies) const;
197 Eigen::VectorXd SolveQP(Spin spin, const Eigen::VectorXd& frequencies) const;
198
199 boost::optional<double> SolveQP_Grid(Spin spin, double intercept0,
200 double frequency0, Index gw_level,
201 QPStats* stats = nullptr) const;
202
203 boost::optional<double> SolveQP_Grid_Windowed(
204 Spin spin, double intercept0, double frequency0, Index gw_level,
205 double left_limit, double right_limit, bool allow_rejected_return = true,
206 QPStats* stats = nullptr) const;
207
208 boost::optional<double> SolveQP_Grid_Windowed_Adaptive(
209 Spin spin, double intercept0, double frequency0, Index gw_level,
210 double left_limit, double right_limit, bool allow_rejected_return = true,
211 QPStats* stats = nullptr) const;
212
213 boost::optional<double> SolveQP_Grid_Windowed_Dense(
214 Spin spin, double intercept0, double frequency0, Index gw_level,
215 double left_limit, double right_limit, bool allow_rejected_return = true,
216 QPStats* stats = nullptr) const;
217
218 boost::optional<double> SolveQP_FixedPoint(Spin spin, double intercept0,
219 double frequency0, Index gw_level,
220 QPStats* stats = nullptr) const;
221
222 boost::optional<double> SolveQP_Linearisation(Spin spin, double intercept0,
223 double frequency0,
224 Index gw_level,
225 QPStats* stats = nullptr) const;
226
227 bool Converged(const Eigen::VectorXd& e1, const Eigen::VectorXd& e2,
228 double epsilon) const;
229
230 boost::optional<QPRootCandidate> RefineQPInterval(
231 double lowerbound, double f_lowerbound, double upperbound,
232 double f_upperbound, const QPFunc& f, double reference) const;
233
234 std::string LevelLabel(Spin spin, Index level) const;
235 const char* OccupationTag(Spin spin, Index level) const;
236
242 const Eigen::MatrixXd& vxc_alpha_;
243 const Eigen::MatrixXd& vxc_beta_;
244 const Eigen::VectorXd& dft_energies_alpha_;
245 const Eigen::VectorXd& dft_energies_beta_;
247 std::unique_ptr<Sigma_base_UKS> sigma_alpha_;
248 std::unique_ptr<Sigma_base_UKS> sigma_beta_;
249 Eigen::MatrixXd Sigma_x_alpha_;
250 Eigen::MatrixXd Sigma_x_beta_;
251 Eigen::MatrixXd Sigma_c_alpha_;
252 Eigen::MatrixXd Sigma_c_beta_;
253};
254
255} // namespace xtp
256} // namespace votca
257
258#endif
void CountSigmaStage(EvalStage stage) const
Definition gw_uks.h:152
std::pair< double, double > operator()(double frequency) const
Definition gw_uks.h:109
std::unordered_set< std::uint64_t > seen_frequencies_
Definition gw_uks.h:174
double sigma(double frequency, EvalStage stage=EvalStage::Other) const
Definition gw_uks.h:116
const QPStats & GetStats() const
Definition gw_uks.h:141
double deriv(double frequency) const
Definition gw_uks.h:134
double value(double frequency, EvalStage stage=EvalStage::Other) const
Definition gw_uks.h:130
static std::uint64_t FrequencyKey(double x)
Definition gw_uks.h:144
QPFunc(Index gw_level, const Sigma_base_UKS &sigma, double offset)
Definition gw_uks.h:106
const Sigma_base_UKS & sigma_c_func_
Definition gw_uks.h:172
boost::optional< double > SolveQP_Grid_Windowed(Spin spin, double intercept0, double frequency0, Index gw_level, double left_limit, double right_limit, bool allow_rejected_return=true, QPStats *stats=nullptr) const
Definition gw_uks.cc:671
TCMatrix_gwbse_spin & Mmn_
Definition gw_uks.h:241
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > DiagonalizeQPHamiltonianAlpha() const
Definition gw_uks.cc:779
bool Converged(const Eigen::VectorXd &e1, const Eigen::VectorXd &e2, double epsilon) const
Definition gw_uks.cc:749
void PrintQP_Energies(Spin spin, const Eigen::VectorXd &qp_diag_energies) const
Definition gw_uks.cc:832
Eigen::MatrixXd & SigmaC(Spin spin)
Definition gw_uks.cc:134
const Eigen::VectorXd & DftEnergies(Spin spin) const
Definition gw_uks.cc:125
std::unique_ptr< Sigma_base_UKS > sigma_alpha_
Definition gw_uks.h:247
boost::optional< double > SolveQP_Linearisation(Spin spin, double intercept0, double frequency0, Index gw_level, QPStats *stats=nullptr) const
Definition gw_uks.cc:385
std::string LevelLabel(Spin spin, Index level) const
Definition gw_uks.cc:156
Eigen::MatrixXd Sigma_x_alpha_
Definition gw_uks.h:249
std::unique_ptr< Sigma_base_UKS > sigma_beta_
Definition gw_uks.h:248
boost::optional< double > SolveQP_Grid_Windowed_Adaptive(Spin spin, double intercept0, double frequency0, Index gw_level, double left_limit, double right_limit, bool allow_rejected_return=true, QPStats *stats=nullptr) const
Definition gw_uks.cc:473
boost::optional< double > SolveQP_Grid(Spin spin, double intercept0, double frequency0, Index gw_level, QPStats *stats=nullptr) const
Definition gw_uks.cc:409
Eigen::MatrixXd Sigma_c_beta_
Definition gw_uks.h:252
Eigen::MatrixXd & SigmaX(Spin spin)
Definition gw_uks.cc:131
Eigen::VectorXd getGWAResultsAlpha() const
Definition gw_uks.cc:297
boost::optional< QPRootCandidate > RefineQPInterval(double lowerbound, double f_lowerbound, double upperbound, double f_upperbound, const QPFunc &f, double reference) const
const Eigen::MatrixXd & vxc_alpha_
Definition gw_uks.h:242
Eigen::VectorXd ScissorShift_DFTlevel(const Eigen::VectorXd &dft_energies, Index homo) const
Definition gw_uks.cc:170
void CalculateHQP()
Definition gw_uks.cc:758
const Eigen::VectorXd & RPAInputEnergiesAlpha() const
Definition gw_uks.cc:307
qp_solver::WindowDiagnostics QPWindowDiagnostics
Definition gw_uks.h:30
const char * OccupationTag(Spin spin, Index level) const
Definition gw_uks.cc:166
qp_solver::Stats QPStats
Definition gw_uks.h:28
Sigma_base_UKS & SigmaEvaluator(Spin spin)
Definition gw_uks.cc:143
Eigen::VectorXd SolveQP(Spin spin, const Eigen::VectorXd &frequencies) const
Definition gw_uks.cc:314
const Eigen::VectorXd & RPAInputEnergiesBeta() const
Definition gw_uks.cc:310
Eigen::MatrixXd getHQPBeta() const
Definition gw_uks.cc:772
void PrintGWA_Energies(Spin spin) const
Definition gw_uks.cc:792
Eigen::VectorXd getGWAResultsBeta() const
Definition gw_uks.cc:302
Eigen::MatrixXd getHQPAlpha() const
Definition gw_uks.cc:767
Eigen::MatrixXd Sigma_x_beta_
Definition gw_uks.h:250
Index gw_sc_iteration_
Definition gw_uks.h:240
const char * SpinName(Spin spin) const
Definition gw_uks.cc:152
GW_UKS(Logger &log, TCMatrix_gwbse_spin &Mmn, const Eigen::MatrixXd &vxc_alpha, const Eigen::MatrixXd &vxc_beta, const Eigen::VectorXd &dft_energies_alpha, const Eigen::VectorXd &dft_energies_beta)
Definition gw_uks.cc:25
const Eigen::VectorXd & dft_energies_alpha_
Definition gw_uks.h:244
double CalcSpinHomoLumoShift(const Eigen::VectorXd &frequencies, Spin spin) const
Definition gw_uks.cc:178
boost::optional< double > SolveQP_Grid_Windowed_Dense(Spin spin, double intercept0, double frequency0, Index gw_level, double left_limit, double right_limit, bool allow_rejected_return=true, QPStats *stats=nullptr) const
Definition gw_uks.cc:548
const Eigen::MatrixXd & vxc_beta_
Definition gw_uks.h:243
const Eigen::VectorXd & dft_energies_beta_
Definition gw_uks.h:245
qp_solver::EvalStage EvalStage
Definition gw_uks.h:27
Logger & log_
Definition gw_uks.h:239
qp_solver::RootCandidate QPRootCandidate
Definition gw_uks.h:29
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > DiagonalizeQPHamiltonianBeta() const
Definition gw_uks.cc:786
Eigen::MatrixXd Sigma_c_alpha_
Definition gw_uks.h:251
const Eigen::MatrixXd & Vxc(Spin spin) const
Definition gw_uks.cc:128
void configure(const options &opt)
Definition gw_uks.cc:38
void CalculateGWPerturbation()
Definition gw_uks.cc:188
boost::optional< double > SolveQP_FixedPoint(Spin spin, double intercept0, double frequency0, Index gw_level, QPStats *stats=nullptr) const
Definition gw_uks.cc:727
Index Homo(Spin spin) const
Definition gw_uks.cc:149
Logger is used for thread-safe output of messages.
Definition logger.h:164
Unrestricted RPA helper for spin-resolved GW screening.
Definition rpa_uks.h:68
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
std::string quadrature_scheme
Definition gw_uks.h:71
std::string qp_grid_search_mode
Definition gw_uks.h:78
std::string sigma_integration
Definition gw_uks.h:47
std::string qp_root_finder
Definition gw_uks.h:77