votca 2024.2-dev
Loading...
Searching...
No Matches
gw.cc
Go to the documentation of this file.
1/*
2 * Copyright 2009-2023 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// Standard includes
21#include <fstream>
22#include <iostream>
23
24// Local VOTCA includes
27#include "votca/xtp/gw.h"
29#include "votca/xtp/rpa.h"
31
32namespace votca {
33namespace xtp {
34
35void GW::configure(const options& opt) {
36 opt_ = opt;
40 Sigma_base::options sigma_opt;
41 sigma_opt.homo = opt_.homo;
42 sigma_opt.qpmax = opt_.qpmax;
43 sigma_opt.qpmin = opt_.qpmin;
44 sigma_opt.rpamin = opt_.rpamin;
45 sigma_opt.rpamax = opt_.rpamax;
46 sigma_opt.eta = opt_.eta;
47 sigma_opt.alpha = opt_.alpha;
49 sigma_opt.order = opt_.order;
50 sigma_->configure(sigma_opt);
51 Sigma_x_ = Eigen::MatrixXd::Zero(qptotal_, qptotal_);
52 Sigma_c_ = Eigen::MatrixXd::Zero(qptotal_, qptotal_);
53}
54
55double GW::CalcHomoLumoShift(Eigen::VectorXd frequencies) const {
56 double DFTgap = dft_energies_(opt_.homo + 1) - dft_energies_(opt_.homo);
57 double QPgap = frequencies(opt_.homo + 1 - opt_.qpmin) -
58 frequencies(opt_.homo - opt_.qpmin);
59 return QPgap - DFTgap;
60}
61
62Eigen::MatrixXd GW::getHQP() const {
63 return Sigma_x_ + Sigma_c_ - vxc_ +
64 Eigen::MatrixXd(
65 dft_energies_.segment(opt_.qpmin, qptotal_).asDiagonal());
66}
67
68Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> GW::DiagonalizeQPHamiltonian()
69 const {
70 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> qpdiag(getHQP());
71 PrintQP_Energies(qpdiag.eigenvalues());
72 return qpdiag;
73}
74
76 Eigen::VectorXd gwa_energies = getGWAResults();
77 double shift = CalcHomoLumoShift(gwa_energies);
78
80 << (boost::format(
81 " ====== Perturbative quasiparticle energies (Hartree) ====== "))
82 .str()
83 << std::flush;
85 << (boost::format(" DeltaHLGap = %1$+1.6f Hartree") % shift).str()
86 << std::flush;
87
88 for (Index i = 0; i < qptotal_; i++) {
89 std::string level = " Level";
90 if ((i + opt_.qpmin) == opt_.homo) {
91 level = " HOMO ";
92 } else if ((i + opt_.qpmin) == opt_.homo + 1) {
93 level = " LUMO ";
94 }
95
97 << level
98 << (boost::format(" = %1$4d DFT = %2$+1.4f VXC = %3$+1.4f S-X = "
99 "%4$+1.4f S-C = %5$+1.4f GWA = %6$+1.4f") %
100 (i + opt_.qpmin) % dft_energies_(i + opt_.qpmin) % vxc_(i, i) %
101 Sigma_x_(i, i) % Sigma_c_(i, i) % gwa_energies(i))
102 .str()
103 << std::flush;
104 }
105 return;
106}
107
108void GW::PrintQP_Energies(const Eigen::VectorXd& qp_diag_energies) const {
109 Eigen::VectorXd gwa_energies = getGWAResults();
111 << TimeStamp() << " Full quasiparticle Hamiltonian " << std::flush;
113 << (boost::format(
114 " ====== Diagonalized quasiparticle energies (Hartree) "
115 "====== "))
116 .str()
117 << std::flush;
118 for (Index i = 0; i < qptotal_; i++) {
119 std::string level = " Level";
120 if ((i + opt_.qpmin) == opt_.homo) {
121 level = " HOMO ";
122 } else if ((i + opt_.qpmin) == opt_.homo + 1) {
123 level = " LUMO ";
124 }
126 << level
127 << (boost::format(" = %1$4d PQP = %2$+1.6f DQP = %3$+1.6f ") %
128 (i + opt_.qpmin) % gwa_energies(i) % qp_diag_energies(i))
129 .str()
130 << std::flush;
131 }
132 return;
133}
134
136 const Eigen::VectorXd& dft_energies) const {
137 Eigen::VectorXd shifted_energies = dft_energies;
138 shifted_energies.segment(opt_.homo + 1, dft_energies.size() - opt_.homo - 1)
139 .array() += opt_.shift;
140 return shifted_energies;
141}
142
144
145 Sigma_x_ = (1 - opt_.ScaHFX) * sigma_->CalcExchangeMatrix();
147 << TimeStamp() << " Calculated Hartree exchange contribution"
148 << std::flush;
149 // dftenergies has size aobasissize
150 // rpaenergies/Mmn have size rpatotal
151 // gwaenergies/frequencies have size qptotal
152 // homo index is relative to dft_energies
154 << TimeStamp() << " Scissor shifting DFT energies by: " << opt_.shift
155 << " Hrt" << std::flush;
156 Eigen::VectorXd dft_shifted_energies = ScissorShift_DFTlevel(dft_energies_);
158 dft_shifted_energies.segment(opt_.rpamin, opt_.rpamax - opt_.rpamin + 1));
159 Eigen::VectorXd frequencies =
160 dft_shifted_energies.segment(opt_.qpmin, qptotal_);
161
162 Anderson mixing_;
164
165 for (Index i_gw = 0; i_gw < opt_.gw_sc_max_iterations; ++i_gw) {
166
167 if (i_gw % opt_.reset_3c == 0 && i_gw != 0) {
168 Mmn_.Rebuild();
170 << TimeStamp() << " Rebuilding 3c integrals" << std::flush;
171 }
172 sigma_->PrepareScreening();
174 << TimeStamp() << " Calculated screening via RPA" << std::flush;
176 << TimeStamp() << " Solving QP equations " << std::flush;
177 if (opt_.gw_mixing_order > 0 && i_gw > 0) {
178 mixing_.UpdateInput(frequencies);
179 }
180
181 frequencies = SolveQP(frequencies);
182
183 if (opt_.gw_sc_max_iterations > 1) {
184 Eigen::VectorXd rpa_energies_old = rpa_.getRPAInputEnergies();
185 if (opt_.gw_mixing_order > 0 && i_gw > 0) {
186 if (opt_.gw_mixing_order == 1) {
188 << "GWSC using linear mixing with alpha: " << opt_.gw_mixing_alpha
189 << std::flush;
190 } else {
192 << "GWSC using Anderson mixing with history "
193 << opt_.gw_mixing_order << ", alpha: " << opt_.gw_mixing_alpha
194 << std::flush;
195 }
196 mixing_.UpdateOutput(frequencies);
197 Eigen::VectorXd mixed_frequencies = mixing_.MixHistory();
198 rpa_.UpdateRPAInputEnergies(dft_energies_, mixed_frequencies,
199 opt_.qpmin);
200 frequencies = mixed_frequencies;
201
202 } else {
203 XTP_LOG(Log::debug, log_) << "GWSC using plain update " << std::flush;
205 }
206
207 for (int i = 0; i < frequencies.size(); i++) {
209 << "... GWSC iter " << i_gw << " state " << i << " "
210 << std::setprecision(9) << frequencies(i) << std::flush;
211 }
212
214 << TimeStamp() << " GW_Iteration:" << i_gw
215 << " Shift[Hrt]:" << CalcHomoLumoShift(frequencies) << std::flush;
216 if (Converged(rpa_.getRPAInputEnergies(), rpa_energies_old,
217 opt_.gw_sc_limit)) {
218 XTP_LOG(Log::info, log_) << TimeStamp() << " Converged after "
219 << i_gw + 1 << " GW iterations." << std::flush;
220 break;
221 } else if (i_gw == opt_.gw_sc_max_iterations - 1) {
223 << TimeStamp()
224 << " WARNING! GW-self-consistency cycle not converged after "
225 << opt_.gw_sc_max_iterations << " iterations." << std::flush;
227 << TimeStamp() << " Run continues. Inspect results carefully!"
228 << std::flush;
229 break;
230 }
231 }
232 }
233 Sigma_c_.diagonal() = sigma_->CalcCorrelationDiag(frequencies);
235}
236
237Eigen::VectorXd GW::getGWAResults() const {
238 return Sigma_x_.diagonal() + Sigma_c_.diagonal() - vxc_.diagonal() +
240}
241
242Eigen::VectorXd GW::SolveQP(const Eigen::VectorXd& frequencies) const {
243
244 Eigen::VectorXd env = Eigen::VectorXd::Zero(qptotal_);
245
246 const Eigen::VectorXd intercepts =
247 dft_energies_.segment(opt_.qpmin, qptotal_) + Sigma_x_.diagonal() -
248 vxc_.diagonal();
249
250 Eigen::VectorXd frequencies_new = frequencies;
251 Eigen::Array<bool, Eigen::Dynamic, 1> converged =
252 Eigen::Array<bool, Eigen::Dynamic, 1>::Zero(qptotal_);
253#ifdef _OPENMP
254 Index use_threads =
256#endif
257#pragma omp parallel for schedule(dynamic) num_threads(use_threads)
258 for (Index gw_level = 0; gw_level < qptotal_; ++gw_level) {
259
260 double initial_f = frequencies[gw_level];
261 double intercept = intercepts[gw_level];
262 boost::optional<double> newf;
263 if (opt_.qp_solver == "fixedpoint") {
264 newf = SolveQP_FixedPoint(intercept, initial_f, gw_level);
265 }
266 if (newf) {
267 frequencies_new[gw_level] = newf.value();
268 converged[gw_level] = true;
269 } else {
270 newf = SolveQP_Grid(intercept, initial_f, gw_level);
271 if (newf) {
272 frequencies_new[gw_level] = newf.value();
273 converged[gw_level] = true;
274 } else {
275 newf = SolveQP_Linearisation(intercept, initial_f, gw_level);
276 if (newf) {
277 frequencies_new[gw_level] = newf.value();
278 }
279 }
280 }
281 }
282
283 if (!converged.all()) {
284 std::vector<Index> states;
285 for (Index s = 0; s < converged.size(); s++) {
286 if (!converged[s]) {
287 states.push_back(s);
288 }
289 }
290 IndexParser rp;
291 XTP_LOG(Log::error, log_) << TimeStamp() << " Not converged PQP states are:"
292 << rp.CreateIndexString(states) << std::flush;
294 << TimeStamp() << " Increase the grid search interval" << std::flush;
295 }
296 return frequencies_new;
297}
298
299boost::optional<double> GW::SolveQP_Linearisation(double intercept0,
300 double frequency0,
301 Index gw_level) const {
302 boost::optional<double> newf = boost::none;
303
304 double sigma = sigma_->CalcCorrelationDiagElement(gw_level, frequency0);
305 double dsigma_domega =
306 sigma_->CalcCorrelationDiagElementDerivative(gw_level, frequency0);
307 double Z = 1.0 - dsigma_domega;
308 if (std::abs(Z) > 1e-9) {
309 newf = frequency0 + (intercept0 - frequency0 + sigma) / Z;
310 }
311 return newf;
312}
313
314boost::optional<double> GW::SolveQP_Grid(double intercept0, double frequency0,
315 Index gw_level) const {
316 std::vector<std::pair<double, double>> roots;
317 const double range =
318 opt_.qp_grid_spacing * double(opt_.qp_grid_steps - 1) / 2.0;
319 boost::optional<double> newf = boost::none;
320 double freq_prev = frequency0 - range;
321 QPFunc fqp(gw_level, *sigma_.get(), intercept0);
322 double targ_prev = fqp.value(freq_prev);
323 double qp_energy = 0.0;
324 double gradient_max = std::numeric_limits<double>::max();
325 bool pole_found = false;
326 for (Index i_node = 1; i_node < opt_.qp_grid_steps; ++i_node) {
327 double freq = freq_prev + opt_.qp_grid_spacing;
328 double targ = fqp.value(freq);
329 if (targ_prev * targ < 0.0) { // Sign change
330 double f = SolveQP_Bisection(freq_prev, targ_prev, freq, targ, fqp);
331 double gradient = fqp.deriv(f);
332 double qp_weight = -1.0 / gradient;
333 roots.push_back(std::make_pair(f, qp_weight));
334 if (std::abs(gradient) < gradient_max) {
335 gradient_max = std::abs(gradient);
336 qp_energy = f;
337 pole_found = true;
338 }
339 }
340 freq_prev = freq;
341 targ_prev = targ;
342 }
344#pragma omp critical
345 {
346 if (!pole_found) {
348 << " No roots found for qplevel:" << gw_level << std::flush;
349 } else {
350 XTP_LOG(Log::info, log_) << " Roots found for qplevel:" << gw_level
351 << " (qpenergy:qpweight)\n\t\t";
352 for (auto& root : roots) {
353 XTP_LOG(Log::info, log_) << std::setprecision(5) << root.first << ":"
354 << root.second << " ";
355 }
356 XTP_LOG(Log::info, log_) << "Root chosen " << qp_energy << std::flush;
357 }
358 }
359 }
360
361 if (pole_found) {
362 newf = qp_energy;
363 }
364 return newf;
365}
366
367boost::optional<double> GW::SolveQP_FixedPoint(double intercept0,
368 double frequency0,
369 Index gw_level) const {
370 boost::optional<double> newf = boost::none;
371 QPFunc f(gw_level, *sigma_.get(), intercept0);
374 double freq_new = newton.FindRoot(f, frequency0);
375 if (newton.getInfo() == NewtonRapson<QPFunc>::success) {
376 newf = freq_new;
377 }
378 return newf;
379}
380
381// https://en.wikipedia.org/wiki/Bisection_method
382double GW::SolveQP_Bisection(double lowerbound, double f_lowerbound,
383 double upperbound, double f_upperbound,
384 const QPFunc& f) const {
385
386 if (f_lowerbound * f_upperbound > 0) {
387 throw std::runtime_error(
388 "Bisection needs a postive and negative function value");
389 }
390 double zero = 0.0;
391 while (true) {
392 double c = 0.5 * (lowerbound + upperbound);
393 if (std::abs(upperbound - lowerbound) < opt_.g_sc_limit) {
394 zero = c;
395 break;
396 }
397 double y_c = f.value(c);
398 if (std::abs(y_c) < opt_.g_sc_limit) {
399 zero = c;
400 break;
401 }
402 if (y_c * f_lowerbound > 0) {
403 lowerbound = c;
404 f_lowerbound = y_c;
405 } else {
406 upperbound = c;
407 f_upperbound = y_c;
408 }
409 }
410 return zero;
411}
412
413bool GW::Converged(const Eigen::VectorXd& e1, const Eigen::VectorXd& e2,
414 double epsilon) const {
415 Index state = 0;
416 bool energies_converged = true;
417 double diff_max = (e1 - e2).cwiseAbs().maxCoeff(&state);
418 if (diff_max > epsilon) {
419 energies_converged = false;
420 }
421 XTP_LOG(Log::info, log_) << TimeStamp() << " E_diff max=" << diff_max
422 << " StateNo:" << state << std::flush;
423 return energies_converged;
424}
425
427 Eigen::VectorXd diag_backup = Sigma_c_.diagonal();
428 Sigma_c_ = sigma_->CalcCorrelationOffDiag(getGWAResults());
429 Sigma_c_.diagonal() = diag_backup;
430}
431
432void GW::PlotSigma(std::string filename, Index steps, double spacing,
433 std::string states) const {
434
435 Eigen::VectorXd frequencies =
437
438 std::vector<Index> state_inds;
439 IndexParser rp;
440 std::vector<Index> parsed_states = rp.CreateIndexVector(states);
441 for (Index gw_level : parsed_states) {
442 if (gw_level >= opt_.qpmin && gw_level <= opt_.qpmax) {
443 state_inds.push_back(gw_level);
444 }
445 }
447 << TimeStamp() << " PQP(omega) written to '" << filename
448 << "' for states " << rp.CreateIndexString(state_inds) << std::flush;
449
450 const Index num_states = state_inds.size();
451
452 const Eigen::VectorXd intercept =
453 dft_energies_.segment(opt_.qpmin, qptotal_) + Sigma_x_.diagonal() -
454 vxc_.diagonal();
455 Eigen::MatrixXd mat = Eigen::MatrixXd::Zero(steps, 2 * num_states);
456#pragma omp parallel for schedule(dynamic)
457 for (Index grid_point = 0; grid_point < steps; grid_point++) {
458 const double offset =
459 ((double)grid_point - ((double)(steps - 1) / 2.0)) * spacing;
460 for (Index i = 0; i < num_states; i++) {
461 const Index gw_level = state_inds[i];
462 const double omega = frequencies(gw_level) + offset;
463 double sigma = sigma_->CalcCorrelationDiagElement(gw_level, omega);
464 mat(grid_point, 2 * i) = omega;
465 mat(grid_point, 2 * i + 1) = sigma + intercept[gw_level];
466 }
467 }
468
469 std::ofstream out;
470 out.open(filename);
471 for (Index i = 0; i < num_states; i++) {
472 const Index gw_level = state_inds[i];
473 out << boost::format("#%1$somega_%2$d\tE_QP(omega)_%2$d") %
474 (i == 0 ? "" : "\t") % gw_level;
475 }
476 out << std::endl;
477 boost::format numFormat("%+1.6f");
478 Eigen::IOFormat matFormat(Eigen::StreamPrecision, 0, "\t", "\n");
479 out << numFormat % mat.format(matFormat) << std::endl;
480 out.close();
481}
482
483} // namespace xtp
484} // namespace votca
virtual std::unique_ptr< T > Create(const key_t &key, args_t &&...arguments)
Anderson mixing as convergence acceleration in SCF/fixed point problems.
void UpdateInput(const Eigen::VectorXd &newInput)
const Eigen::VectorXd MixHistory()
void UpdateOutput(const Eigen::VectorXd &newOutput)
void Configure(const Index order, const double alpha)
double deriv(double frequency) const
Definition gw.h:126
double value(double frequency) const
Definition gw.h:122
void PlotSigma(std::string filename, Index steps, double spacing, std::string states) const
Definition gw.cc:432
Index qptotal_
Definition gw.h:94
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > DiagonalizeQPHamiltonian() const
Definition gw.cc:68
void PrintGWA_Energies() const
Definition gw.cc:75
RPA rpa_
Definition gw.h:107
double SolveQP_Bisection(double lowerbound, double f_lowerbound, double upperbound, double f_upperbound, const QPFunc &f) const
Definition gw.cc:382
boost::optional< double > SolveQP_FixedPoint(double intercept0, double frequency0, Index gw_level) const
Definition gw.cc:367
bool Converged(const Eigen::VectorXd &e1, const Eigen::VectorXd &e2, double epsilon) const
Definition gw.cc:413
void configure(const options &opt)
Definition gw.cc:35
options opt_
Definition gw.h:99
Eigen::MatrixXd getHQP() const
Definition gw.cc:62
Eigen::VectorXd SolveQP(const Eigen::VectorXd &frequencies) const
Definition gw.cc:242
TCMatrix_gwbse & Mmn_
Definition gw.h:103
Eigen::MatrixXd Sigma_c_
Definition gw.h:97
boost::optional< double > SolveQP_Grid(double intercept0, double frequency0, Index gw_level) const
Definition gw.cc:314
const Eigen::VectorXd & dft_energies_
Definition gw.h:105
std::unique_ptr< Sigma_base > sigma_
Definition gw.h:101
void CalculateHQP()
Definition gw.cc:426
Logger & log_
Definition gw.h:102
Eigen::VectorXd getGWAResults() const
Definition gw.cc:237
double CalcHomoLumoShift(Eigen::VectorXd frequencies) const
Definition gw.cc:55
const Eigen::MatrixXd & vxc_
Definition gw.h:104
void PrintQP_Energies(const Eigen::VectorXd &qp_diag_energies) const
Definition gw.cc:108
boost::optional< double > SolveQP_Linearisation(double intercept0, double frequency0, Index gw_level) const
Definition gw.cc:299
Eigen::VectorXd ScissorShift_DFTlevel(const Eigen::VectorXd &dft_energies) const
Definition gw.cc:135
Eigen::MatrixXd Sigma_x_
Definition gw.h:96
void CalculateGWPerturbation()
Definition gw.cc:143
std::vector< Index > CreateIndexVector(const std::string &Ids) const
std::string CreateIndexString(const std::vector< Index > &indeces) const
Newton Rapson rootfinder for 1d functions.
double FindRoot(const Func &f, double x0)
void setRPAInputEnergies(const Eigen::VectorXd &rpaenergies)
Definition rpa.h:59
void configure(Index homo, Index rpamin, Index rpamax)
Definition rpa.h:39
const Eigen::VectorXd & getRPAInputEnergies() const
Definition rpa.h:57
void UpdateRPAInputEnergies(const Eigen::VectorXd &dftenergies, const Eigen::VectorXd &gwaenergies, Index qpmin)
Definition rpa.cc:30
Timestamp returns the current time as a string Example: cout << TimeStamp()
Definition logger.h:224
#define XTP_LOG(level, log)
Definition logger.h:40
Index getMaxThreads()
Definition eigen.h:128
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
static Level current_level
Definition globals.h:30
double qp_solver_alpha
Definition gw.h:61
double gw_sc_limit
Definition gw.h:53
Index gw_mixing_order
Definition gw.h:64
Index gw_sc_max_iterations
Definition gw.h:54
Index qp_grid_steps
Definition gw.h:62
std::string quadrature_scheme
Definition gw.h:66
std::string qp_solver
Definition gw.h:60
Index g_sc_max_iterations
Definition gw.h:52
double qp_grid_spacing
Definition gw.h:63
double g_sc_limit
Definition gw.h:51
double gw_mixing_alpha
Definition gw.h:65
std::string sigma_integration
Definition gw.h:57