votca 2026-dev
Loading...
Searching...
No Matches
gw.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2020 The VOTCA Development Team
3 * (http://www.votca.org)
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License")
6 *
7 * You may not use this file except in compliance with the License.
8 * You may obtain a copy of the License at
9 *
10 * http://www.apache.org/licenses/LICENSE-2.0
11 *
12 * Unless required by applicable law or agreed to in writing, software
13 * distributed under the License is distributed on an "AS IS" BASIS,
14 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
15 * See the License for the specific language governing permissions and
16 * limitations under the License.
17 *
18 */
19
20#pragma once
21#ifndef VOTCA_XTP_GW_H
22#define VOTCA_XTP_GW_H
23
24#include <unordered_set>
25
26// Local VOTCA includes
27#include "logger.h"
28#include "orbitals.h"
29#include "qp_solver_utils.h"
30#include "rpa.h"
31#include "sigma_base.h"
32#include "threecenter.h"
33
34namespace votca {
35namespace xtp {
36
37class GW {
38
43
44 public:
45 GW(Logger& log, TCMatrix_gwbse& Mmn, const Eigen::MatrixXd& vxc,
46 const Eigen::VectorXd& dft_energies)
47 : log_(log),
48 Mmn_(Mmn),
49 vxc_(vxc),
50 dft_energies_(dft_energies),
51 rpa_(log, Mmn) {};
52
53 struct options {
59 double eta;
60 double g_sc_limit;
64 double shift = 0;
65 double ScaHFX = 0.0;
66 std::string sigma_integration;
67 Index reset_3c; // how often the 3c integrals in iterate should be
68 // rebuilt
69 std::string qp_solver;
70 double qp_solver_alpha = 0.75;
71 // Legacy aliases from the original grid-search implementation.
72 // They used to define both the total QP window and the scan resolution at
73 // the same time. The new implementation keeps them only so that old XML
74 // files and existing tests still work; configuration code maps them onto
75 // the decoupled controls below.
77 double qp_grid_spacing = 0.0;
78
79 // Decoupled QP search controls:
80 // - qp_full_window_half_width defines the outer full search interval
81 // - qp_dense_spacing defines the dense scan used for robust sign changes
82 // - qp_adaptive_shell_width / count define the outward shell search
83 // All values are normalized in configure(); negative values mean "unset"
84 // and trigger either legacy mapping or robust defaults.
85 double qp_full_window_half_width = -1.0; // Ha; <=0 means "unset"
86 double qp_dense_spacing = -1.0; // Ha; <=0 means "unset"
87 double qp_adaptive_shell_width = -1.0; // Ha; <=0 means "unset"
88 Index qp_adaptive_shell_count = 0; // 0 means "use shell width"
89
90 Index gw_mixing_order; // mixing order
91 double gw_mixing_alpha; // mixing alpha, also linear mixing
92 std::string quadrature_scheme; // Kind of Gaussian-quadrature scheme to use
93 Index order; // only needed for complex integration sigma CDA
94 double alpha; // smooth tail in complex integration sigma CDA
95 bool qp_restrict_search = true;
96 double qp_zero_margin = 1e-6;
97 double qp_virtual_min_energy = -0.1;
98 std::string qp_root_finder = "bisection";
99 std::string qp_grid_search_mode = "adaptive_with_dense_fallback";
100 };
101
102 void configure(const options& opt);
103
104 Eigen::VectorXd getGWAResults() const;
105 // Calculates the diagonal elements up to self consistency
107
108 // Calculated offdiagonal elements as well
109 void CalculateHQP();
110
111 Eigen::MatrixXd getHQP() const;
112
113 // Diagonalize QP particle Hamiltonian
114 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> DiagonalizeQPHamiltonian()
115 const;
116
117 void PlotSigma(std::string filename, Index steps, double spacing,
118 std::string states) const;
119
120 Eigen::VectorXd RPAInputEnergies() const {
121 return rpa_.getRPAInputEnergies();
122 }
123
124 private:
126
127 Eigen::MatrixXd Sigma_x_;
128 Eigen::MatrixXd Sigma_c_;
129
131
132 std::unique_ptr<Sigma_base> sigma_ = nullptr;
135 const Eigen::MatrixXd& vxc_;
136 const Eigen::VectorXd& dft_energies_;
137
139
141 // small class which calculates f(w) with and df/dw(w)
142 // f=Sigma_c(w)+offset-w
143 // offset= e_dft+Sigma_x-Vxc
144 class QPFunc {
145 public:
146 QPFunc(Index gw_level, const Sigma_base& sigma, double offset)
147 : gw_level_(gw_level), offset_(offset), sigma_c_func_(sigma) {}
148
149 std::pair<double, double> operator()(double frequency) const {
150 std::pair<double, double> result;
151 result.first = value(frequency, EvalStage::Other);
152 result.second = deriv(frequency);
153 return result;
154 }
155
156 double sigma(double frequency, EvalStage stage = EvalStage::Other) const {
157 const std::uint64_t key = FrequencyKey(frequency);
158
159 auto insert_result = seen_frequencies_.insert(key);
160 if (!insert_result.second) {
161 ++stats_.sigma_repeat_calls;
162 } else {
163 ++stats_.sigma_unique_frequencies;
164 }
165
166 CountSigmaStage(stage);
167 return sigma_c_func_.CalcCorrelationDiagElement(gw_level_, frequency);
168 }
169
170 double value(double frequency, EvalStage stage = EvalStage::Other) const {
171 return sigma(frequency, stage) + offset_ - frequency;
172 }
173
174 double deriv(double frequency) const {
175 ++stats_.deriv_calls;
176 return sigma_c_func_.CalcCorrelationDiagElementDerivative(gw_level_,
177 frequency) -
178 1.0;
179 }
180
181 const QPStats& GetStats() const { return stats_; }
182
183 private:
184 static std::uint64_t FrequencyKey(double x) {
185 std::uint64_t key = 0;
186 static_assert(sizeof(double) == sizeof(std::uint64_t),
187 "Unexpected double size");
188 std::memcpy(&key, &x, sizeof(double));
189 return key;
190 }
191
192 void CountSigmaStage(EvalStage stage) const {
193 switch (stage) {
194 case EvalStage::Scan:
195 ++stats_.sigma_scan_calls;
196 break;
197 case EvalStage::Refine:
198 ++stats_.sigma_refine_calls;
199 break;
200 case EvalStage::Derivative:
201 ++stats_.sigma_derivative_calls;
202 break;
203 case EvalStage::Other:
204 default:
205 ++stats_.sigma_other_calls;
206 break;
207 }
208 }
209
211 double offset_;
213
214 mutable std::unordered_set<std::uint64_t> seen_frequencies_;
216 };
217
218 double CalcHomoLumoShift(Eigen::VectorXd frequencies) const;
219 Eigen::VectorXd ScissorShift_DFTlevel(
220 const Eigen::VectorXd& dft_energies) const;
221 void PrintQP_Energies(const Eigen::VectorXd& qp_diag_energies) const;
222 void PrintGWA_Energies() const;
223
224 Eigen::VectorXd SolveQP(const Eigen::VectorXd& frequencies) const;
225 boost::optional<double> SolveQP_Grid(double intercept0, double frequency0,
226 Index gw_level,
227 QPStats* stats = nullptr) const;
228
229 boost::optional<double> SolveQP_Grid_Windowed(
230 double intercept0, double frequency0, Index gw_level, double left_limit,
231 double right_limit, bool allow_rejected_return = true,
232 QPStats* stats = nullptr) const;
233
234 boost::optional<double> SolveQP_Grid_Windowed_Adaptive(
235 double intercept0, double frequency0, Index gw_level, double left_limit,
236 double right_limit, bool allow_rejected_return = true,
237 QPStats* stats = nullptr) const;
238
239 boost::optional<double> SolveQP_Grid_Windowed_Dense(
240 double intercept0, double frequency0, Index gw_level, double left_limit,
241 double right_limit, bool allow_rejected_return = true,
242 QPStats* stats = nullptr) const;
243
244 boost::optional<double> SolveQP_FixedPoint(double intercept0,
245 double frequency0, Index gw_level,
246 QPStats* stats = nullptr) const;
247 boost::optional<double> SolveQP_Linearisation(double intercept0,
248 double frequency0,
249 Index gw_level,
250 QPStats* stats = nullptr) const;
251 bool Converged(const Eigen::VectorXd& e1, const Eigen::VectorXd& e2,
252 double epsilon) const;
253
254 boost::optional<QPRootCandidate> RefineQPInterval(
255 double lowerbound, double f_lowerbound, double upperbound,
256 double f_upperbound, const QPFunc& f, double reference) const;
257};
258} // namespace xtp
259} // namespace votca
260
261#endif // VOTCA_XTP_GW_H
void CountSigmaStage(EvalStage stage) const
Definition gw.h:192
double value(double frequency, EvalStage stage=EvalStage::Other) const
Definition gw.h:170
double sigma(double frequency, EvalStage stage=EvalStage::Other) const
Definition gw.h:156
const Sigma_base & sigma_c_func_
Definition gw.h:212
const QPStats & GetStats() const
Definition gw.h:181
std::unordered_set< std::uint64_t > seen_frequencies_
Definition gw.h:214
static std::uint64_t FrequencyKey(double x)
Definition gw.h:184
QPFunc(Index gw_level, const Sigma_base &sigma, double offset)
Definition gw.h:146
double deriv(double frequency) const
Definition gw.h:174
std::pair< double, double > operator()(double frequency) const
Definition gw.h:149
boost::optional< double > SolveQP_Grid_Windowed(double intercept0, double frequency0, Index gw_level, double left_limit, double right_limit, bool allow_rejected_return=true, QPStats *stats=nullptr) const
Definition gw.cc:547
void PlotSigma(std::string filename, Index steps, double spacing, std::string states) const
Definition gw.cc:702
Index qptotal_
Definition gw.h:125
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > DiagonalizeQPHamiltonian() const
Definition gw.cc:73
void PrintGWA_Energies() const
Definition gw.cc:80
RPA rpa_
Definition gw.h:140
boost::optional< QPRootCandidate > RefineQPInterval(double lowerbound, double f_lowerbound, double upperbound, double f_upperbound, const QPFunc &f, double reference) const
bool Converged(const Eigen::VectorXd &e1, const Eigen::VectorXd &e2, double epsilon) const
Definition gw.cc:683
void configure(const options &opt)
Definition gw.cc:35
boost::optional< double > SolveQP_Grid_Windowed_Adaptive(double intercept0, double frequency0, Index gw_level, double left_limit, double right_limit, bool allow_rejected_return=true, QPStats *stats=nullptr) const
Definition gw.cc:357
options opt_
Definition gw.h:130
Eigen::MatrixXd getHQP() const
Definition gw.cc:67
Eigen::VectorXd SolveQP(const Eigen::VectorXd &frequencies) const
Definition gw.cc:247
TCMatrix_gwbse & Mmn_
Definition gw.h:134
Index gw_sc_iteration_
Definition gw.h:138
qp_solver::EvalStage EvalStage
Definition gw.h:39
qp_solver::Stats QPStats
Definition gw.h:40
boost::optional< double > SolveQP_Grid_Windowed_Dense(double intercept0, double frequency0, Index gw_level, double left_limit, double right_limit, bool allow_rejected_return=true, QPStats *stats=nullptr) const
Definition gw.cc:427
Eigen::MatrixXd Sigma_c_
Definition gw.h:128
qp_solver::RootCandidate QPRootCandidate
Definition gw.h:41
const Eigen::VectorXd & dft_energies_
Definition gw.h:136
std::unique_ptr< Sigma_base > sigma_
Definition gw.h:132
void CalculateHQP()
Definition gw.cc:696
Logger & log_
Definition gw.h:133
Eigen::VectorXd getGWAResults() const
Definition gw.cc:242
Eigen::VectorXd RPAInputEnergies() const
Definition gw.h:120
qp_solver::WindowDiagnostics QPWindowDiagnostics
Definition gw.h:42
double CalcHomoLumoShift(Eigen::VectorXd frequencies) const
Definition gw.cc:60
boost::optional< double > SolveQP_Grid(double intercept0, double frequency0, Index gw_level, QPStats *stats=nullptr) const
Definition gw.cc:601
const Eigen::MatrixXd & vxc_
Definition gw.h:135
void PrintQP_Energies(const Eigen::VectorXd &qp_diag_energies) const
Definition gw.cc:113
boost::optional< double > SolveQP_FixedPoint(double intercept0, double frequency0, Index gw_level, QPStats *stats=nullptr) const
Definition gw.cc:665
Eigen::VectorXd ScissorShift_DFTlevel(const Eigen::VectorXd &dft_energies) const
Definition gw.cc:140
GW(Logger &log, TCMatrix_gwbse &Mmn, const Eigen::MatrixXd &vxc, const Eigen::VectorXd &dft_energies)
Definition gw.h:45
Eigen::MatrixXd Sigma_x_
Definition gw.h:127
boost::optional< double > SolveQP_Linearisation(double intercept0, double frequency0, Index gw_level, QPStats *stats=nullptr) const
Definition gw.cc:336
void CalculateGWPerturbation()
Definition gw.cc:148
Logger is used for thread-safe output of messages.
Definition logger.h:164
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
double qp_full_window_half_width
Definition gw.h:85
Index qp_adaptive_shell_count
Definition gw.h:88
double qp_adaptive_shell_width
Definition gw.h:87
double qp_solver_alpha
Definition gw.h:70
double qp_zero_margin
Definition gw.h:96
double gw_sc_limit
Definition gw.h:62
std::string qp_grid_search_mode
Definition gw.h:99
double qp_virtual_min_energy
Definition gw.h:97
Index gw_mixing_order
Definition gw.h:90
Index gw_sc_max_iterations
Definition gw.h:63
Index qp_grid_steps
Definition gw.h:76
std::string quadrature_scheme
Definition gw.h:92
std::string qp_solver
Definition gw.h:69
Index g_sc_max_iterations
Definition gw.h:61
double qp_grid_spacing
Definition gw.h:77
double g_sc_limit
Definition gw.h:60
std::string qp_root_finder
Definition gw.h:98
double gw_mixing_alpha
Definition gw.h:91
std::string sigma_integration
Definition gw.h:66
double qp_dense_spacing
Definition gw.h:86
bool qp_restrict_search
Definition gw.h:95