votca 2026-dev
Loading...
Searching...
No Matches
gw_uks.cc
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
8#include <fstream>
9#include <iostream>
10#include <limits>
11
14#include "votca/xtp/gw_uks.h"
17#include "votca/xtp/rpa_uks.h"
19
21
22namespace votca {
23namespace xtp {
24
26 const Eigen::MatrixXd& vxc_alpha,
27 const Eigen::MatrixXd& vxc_beta,
28 const Eigen::VectorXd& dft_energies_alpha,
29 const Eigen::VectorXd& dft_energies_beta)
30 : log_(log),
31 Mmn_(Mmn),
32 vxc_alpha_(vxc_alpha),
33 vxc_beta_(vxc_beta),
34 dft_energies_alpha_(dft_energies_alpha),
35 dft_energies_beta_(dft_energies_beta),
36 rpa_(log, Mmn) {}
37
38void GW_UKS::configure(const options& opt) {
39 opt_ = opt;
40 qptotal_ = opt_.qpmax - opt_.qpmin + 1;
41 rpa_.configure(opt_.homo_alpha, opt_.homo_beta, opt_.rpamin, opt_.rpamax);
42
43 sigma_alpha_ = SigmaFactory_UKS().Create(opt_.sigma_integration, Mmn_, rpa_,
45 sigma_beta_ = SigmaFactory_UKS().Create(opt_.sigma_integration, Mmn_, rpa_,
47
48 if (!sigma_alpha_ || !sigma_beta_) {
49 throw std::runtime_error(
50 "GW_UKS currently supports only implemented unrestricted sigma "
51 "evaluators. At present this is 'ppm', 'cda', and 'exact'.");
52 }
53
55 sigma_opt.homo = opt_.homo_alpha;
56 sigma_opt.qpmax = opt_.qpmax;
57 sigma_opt.qpmin = opt_.qpmin;
58 sigma_opt.rpamin = opt_.rpamin;
59 sigma_opt.rpamax = opt_.rpamax;
60 sigma_opt.eta = opt_.eta;
61 sigma_opt.alpha = opt_.alpha;
62 sigma_opt.quadrature_scheme = opt_.quadrature_scheme;
63 sigma_opt.order = opt_.order;
64 sigma_alpha_->configure(sigma_opt);
65 sigma_opt.homo = opt_.homo_beta;
66 sigma_beta_->configure(sigma_opt);
67
68 Sigma_x_alpha_ = Eigen::MatrixXd::Zero(qptotal_, qptotal_);
69 Sigma_x_beta_ = Eigen::MatrixXd::Zero(qptotal_, qptotal_);
70 Sigma_c_alpha_ = Eigen::MatrixXd::Zero(qptotal_, qptotal_);
71 Sigma_c_beta_ = Eigen::MatrixXd::Zero(qptotal_, qptotal_);
72
73 auto print_spin_ranges = [&](Spin spin, Index homo) {
74 const Index lumo = homo + 1;
75
76 const Index n_occ_rpa =
77 (opt_.rpamin <= homo) ? (homo - opt_.rpamin + 1) : 0;
78 const Index n_virt_rpa =
79 (lumo <= opt_.rpamax) ? (opt_.rpamax - lumo + 1) : 0;
80
81 const Index qp_occ_min = opt_.qpmin;
82 const Index qp_occ_max = std::min(opt_.qpmax, homo);
83 const Index n_occ_qp =
84 (qp_occ_min <= qp_occ_max) ? (qp_occ_max - qp_occ_min + 1) : 0;
85
86 const Index qp_virt_min = std::max(opt_.qpmin, lumo);
87 const Index qp_virt_max = opt_.qpmax;
88 const Index n_virt_qp =
89 (qp_virt_min <= qp_virt_max) ? (qp_virt_max - qp_virt_min + 1) : 0;
90
92 << TimeStamp() << " UKS " << SpinName(spin)
93 << " effective ranges: HOMO=" << homo << " LUMO=" << lumo
94 << " | RPA occ[" << opt_.rpamin << ":" << homo << "]"
95 << " virt[" << lumo << ":" << opt_.rpamax << "]"
96 << " (#occ=" << n_occ_rpa << ", #virt=" << n_virt_rpa << ")"
97 << " | GW occ[" << qp_occ_min << ":" << qp_occ_max << "]"
98 << " virt[" << qp_virt_min << ":" << qp_virt_max << "]"
99 << " (#occ=" << n_occ_qp << ", #virt=" << n_virt_qp << ")"
100 << std::flush;
101 };
102
103 print_spin_ranges(Spin::Alpha, opt_.homo_alpha);
104 print_spin_ranges(Spin::Beta, opt_.homo_beta);
105
106 if (opt_.sigma_integration == "ppm") {
107 auto* ppm_a = dynamic_cast<Sigma_PPM_UKS*>(sigma_alpha_.get());
108 auto* ppm_b = dynamic_cast<Sigma_PPM_UKS*>(sigma_beta_.get());
109
110 if (ppm_a == nullptr || ppm_b == nullptr) {
111 throw std::runtime_error(
112 "GW_UKS: expected Sigma_PPM_UKS evaluators for "
113 "sigma_integration=ppm");
114 }
115
116 ppm_a->SetSharedPPM(ppm_);
117 ppm_b->SetSharedPPM(ppm_);
118 }
119}
120
121const Eigen::VectorXd& GW_UKS::DftEnergies(Spin spin) const {
123}
124const Eigen::MatrixXd& GW_UKS::Vxc(Spin spin) const {
125 return (spin == Spin::Alpha) ? vxc_alpha_ : vxc_beta_;
126}
127Eigen::MatrixXd& GW_UKS::SigmaX(Spin spin) {
128 return (spin == Spin::Alpha) ? Sigma_x_alpha_ : Sigma_x_beta_;
129}
130Eigen::MatrixXd& GW_UKS::SigmaC(Spin spin) {
131 return (spin == Spin::Alpha) ? Sigma_c_alpha_ : Sigma_c_beta_;
132}
133const Eigen::MatrixXd& GW_UKS::SigmaX(Spin spin) const {
134 return (spin == Spin::Alpha) ? Sigma_x_alpha_ : Sigma_x_beta_;
135}
136const Eigen::MatrixXd& GW_UKS::SigmaC(Spin spin) const {
137 return (spin == Spin::Alpha) ? Sigma_c_alpha_ : Sigma_c_beta_;
138}
143 return (spin == Spin::Alpha) ? *sigma_alpha_ : *sigma_beta_;
144}
146 return (spin == Spin::Alpha) ? opt_.homo_alpha : opt_.homo_beta;
147}
148const char* GW_UKS::SpinName(Spin spin) const {
149 return (spin == Spin::Alpha) ? "alpha" : "beta";
150}
151
152std::string GW_UKS::LevelLabel(Spin spin, Index level) const {
153 if (level == Homo(spin)) {
154 return " HOMO ";
155 } else if (level == Homo(spin) + 1) {
156 return " LUMO ";
157 } else {
158 return " Level";
159 }
160}
161
162const char* GW_UKS::OccupationTag(Spin spin, Index level) const {
163 return (level <= Homo(spin)) ? "occ" : "virt";
164}
165
167 const Eigen::VectorXd& dft_energies, Index homo) const {
168 Eigen::VectorXd shifted_energies = dft_energies;
169 shifted_energies.segment(homo + 1, dft_energies.size() - homo - 1).array() +=
170 opt_.shift;
171 return shifted_energies;
172}
173
174double GW_UKS::CalcSpinHomoLumoShift(const Eigen::VectorXd& frequencies,
175 Spin spin) const {
176 const Index homo = Homo(spin);
177 const Eigen::VectorXd& dft = DftEnergies(spin);
178 const double DFTgap = dft(homo + 1) - dft(homo);
179 const double QPgap =
180 frequencies(homo + 1 - opt_.qpmin) - frequencies(homo - opt_.qpmin);
181 return QPgap - DFTgap;
182}
183
185 Sigma_x_alpha_ = (1 - opt_.ScaHFX) * sigma_alpha_->CalcExchangeMatrix();
186 Sigma_x_beta_ = (1 - opt_.ScaHFX) * sigma_beta_->CalcExchangeMatrix();
188 << TimeStamp()
189 << " Calculated spin-resolved Hartree exchange contribution"
190 << std::flush;
191
192 Eigen::VectorXd dft_shifted_alpha =
194 Eigen::VectorXd dft_shifted_beta =
196
198 << TimeStamp()
199 << " Scissor shifting alpha/beta DFT energies by: " << opt_.shift
200 << " Hrt" << std::flush;
201
202 rpa_.setRPAInputEnergies(
203 dft_shifted_alpha.segment(opt_.rpamin, opt_.rpamax - opt_.rpamin + 1),
204 dft_shifted_beta.segment(opt_.rpamin, opt_.rpamax - opt_.rpamin + 1));
205
206 Eigen::VectorXd frequencies_alpha =
207 dft_shifted_alpha.segment(opt_.qpmin, qptotal_);
208 Eigen::VectorXd frequencies_beta =
209 dft_shifted_beta.segment(opt_.qpmin, qptotal_);
210
211 Anderson mixing_alpha;
212 Anderson mixing_beta;
213 mixing_alpha.Configure(opt_.gw_mixing_order, opt_.gw_mixing_alpha);
214 mixing_beta.Configure(opt_.gw_mixing_order, opt_.gw_mixing_alpha);
215
216 for (Index i_gw = 0; i_gw < opt_.gw_sc_max_iterations; ++i_gw) {
217 if (i_gw % opt_.reset_3c == 0 && i_gw != 0) {
218 Mmn_.alpha.Rebuild();
219 Mmn_.beta.Rebuild();
221 << TimeStamp() << " Rebuilding alpha/beta 3c integrals" << std::flush;
222 }
223
224 if (opt_.sigma_integration == "ppm") {
225 ppm_.PPM_construct_parameters(rpa_);
226 sigma_alpha_->PrepareScreening();
227 sigma_beta_->PrepareScreening();
228 } else {
229 sigma_alpha_->PrepareScreening();
230 sigma_beta_->PrepareScreening();
231 }
232
234 << TimeStamp() << " Calculated unrestricted screening via RPA"
235 << std::flush;
236
237 if (opt_.gw_mixing_order > 0 && i_gw > 0) {
238 mixing_alpha.UpdateInput(frequencies_alpha);
239 mixing_beta.UpdateInput(frequencies_beta);
240 }
241
242 frequencies_alpha = SolveQP(Spin::Alpha, frequencies_alpha);
243 frequencies_beta = SolveQP(Spin::Beta, frequencies_beta);
244
245 if (opt_.gw_sc_max_iterations > 1) {
246 Eigen::VectorXd rpa_alpha_old = rpa_.getRPAInputEnergiesAlpha();
247 Eigen::VectorXd rpa_beta_old = rpa_.getRPAInputEnergiesBeta();
248
249 if (opt_.gw_mixing_order > 0 && i_gw > 0) {
250 mixing_alpha.UpdateOutput(frequencies_alpha);
251 mixing_beta.UpdateOutput(frequencies_beta);
252 frequencies_alpha = mixing_alpha.MixHistory();
253 frequencies_beta = mixing_beta.MixHistory();
254 }
255
256 rpa_.UpdateRPAInputEnergies(dft_energies_alpha_, dft_energies_beta_,
257 frequencies_alpha, frequencies_beta,
258 opt_.qpmin);
259
261 << TimeStamp() << " GW_Iteration:" << i_gw << " Shift_alpha[Hrt]:"
262 << CalcSpinHomoLumoShift(frequencies_alpha, Spin::Alpha)
263 << " Shift_beta[Hrt]:"
264 << CalcSpinHomoLumoShift(frequencies_beta, Spin::Beta) << std::flush;
265
266 const bool converged_alpha = Converged(rpa_.getRPAInputEnergiesAlpha(),
267 rpa_alpha_old, opt_.gw_sc_limit);
268 const bool converged_beta = Converged(rpa_.getRPAInputEnergiesBeta(),
269 rpa_beta_old, opt_.gw_sc_limit);
270 if (converged_alpha && converged_beta) {
272 << TimeStamp() << " Converged after " << i_gw + 1
273 << " unrestricted GW iterations." << std::flush;
274 break;
275 } else if (i_gw == opt_.gw_sc_max_iterations - 1) {
277 << TimeStamp()
278 << " WARNING! UKS GW-self-consistency cycle not converged after "
279 << opt_.gw_sc_max_iterations << " iterations." << std::flush;
280 break;
281 }
282 }
283 }
284
285 Sigma_c_alpha_.diagonal() =
286 sigma_alpha_->CalcCorrelationDiag(frequencies_alpha);
287 Sigma_c_beta_.diagonal() = sigma_beta_->CalcCorrelationDiag(frequencies_beta);
290}
291
292Eigen::VectorXd GW_UKS::getGWAResultsAlpha() const {
293 return Sigma_x_alpha_.diagonal() + Sigma_c_alpha_.diagonal() -
294 vxc_alpha_.diagonal() +
295 dft_energies_alpha_.segment(opt_.qpmin, qptotal_);
296}
297Eigen::VectorXd GW_UKS::getGWAResultsBeta() const {
298 return Sigma_x_beta_.diagonal() + Sigma_c_beta_.diagonal() -
299 vxc_beta_.diagonal() +
300 dft_energies_beta_.segment(opt_.qpmin, qptotal_);
301}
302const Eigen::VectorXd& GW_UKS::RPAInputEnergiesAlpha() const {
303 return rpa_.getRPAInputEnergiesAlpha();
304}
305const Eigen::VectorXd& GW_UKS::RPAInputEnergiesBeta() const {
306 return rpa_.getRPAInputEnergiesBeta();
307}
308
309Eigen::VectorXd GW_UKS::SolveQP(Spin spin,
310 const Eigen::VectorXd& frequencies) const {
311 const Eigen::VectorXd intercepts =
312 DftEnergies(spin).segment(opt_.qpmin, qptotal_) +
313 SigmaX(spin).diagonal() - Vxc(spin).diagonal();
314 Eigen::VectorXd frequencies_new = frequencies;
315 Eigen::Array<bool, Eigen::Dynamic, 1> converged =
316 Eigen::Array<bool, Eigen::Dynamic, 1>::Zero(qptotal_);
317 Index use_threads = qptotal_;
318#ifdef _OPENMP
319 use_threads =
321#endif
322#pragma omp parallel for schedule(dynamic) num_threads(use_threads)
323 for (Index gw_level = 0; gw_level < qptotal_; ++gw_level) {
324 const double initial_f = frequencies[gw_level];
325 const double intercept = intercepts[gw_level];
326 boost::optional<double> newf;
327 if (opt_.qp_solver == "fixedpoint") {
328 newf = SolveQP_FixedPoint(spin, intercept, initial_f, gw_level);
329 }
330 if (newf) {
331 frequencies_new[gw_level] = newf.value();
332 converged[gw_level] = true;
333 } else {
334 newf = SolveQP_Grid(spin, intercept, initial_f, gw_level);
335 if (newf) {
336 frequencies_new[gw_level] = newf.value();
337 converged[gw_level] = true;
338 } else {
339 newf = SolveQP_Linearisation(spin, intercept, initial_f, gw_level);
340 if (newf) {
341 frequencies_new[gw_level] = newf.value();
342 }
343 }
344 }
345 }
346 return frequencies_new;
347}
348
349boost::optional<double> GW_UKS::SolveQP_Linearisation(Spin spin,
350 double intercept0,
351 double frequency0,
352 Index gw_level) const {
353 boost::optional<double> newf = boost::none;
354 const auto& sigma = SigmaEvaluator(spin);
355 const double sig = sigma.CalcCorrelationDiagElement(gw_level, frequency0);
356 const double dsigma_domega =
357 sigma.CalcCorrelationDiagElementDerivative(gw_level, frequency0);
358 const double Z = 1.0 - dsigma_domega;
359 if (std::abs(Z) > 1e-9) {
360 newf = frequency0 + (intercept0 - frequency0 + sig) / Z;
361 }
362 return newf;
363}
364
365boost::optional<double> GW_UKS::SolveQP_Grid(Spin spin, double intercept0,
366 double frequency0,
367 Index gw_level) const {
368 const double range =
369 opt_.qp_grid_spacing * double(opt_.qp_grid_steps - 1) / 2.0;
370 boost::optional<double> newf = boost::none;
371 double freq_prev = frequency0 - range;
372 QPFunc fqp(gw_level, SigmaEvaluator(spin), intercept0);
373 double targ_prev = fqp.value(freq_prev);
374 double qp_energy = 0.0;
375 double gradient_max = std::numeric_limits<double>::max();
376 bool pole_found = false;
377 for (Index i_node = 1; i_node < opt_.qp_grid_steps; ++i_node) {
378 double freq = freq_prev + opt_.qp_grid_spacing;
379 double targ = fqp.value(freq);
380 if (targ_prev * targ < 0.0) {
381 double f = SolveQP_Bisection(freq_prev, targ_prev, freq, targ, fqp);
382 double gradient = fqp.deriv(f);
383 if (std::abs(gradient) < gradient_max) {
384 gradient_max = std::abs(gradient);
385 qp_energy = f;
386 pole_found = true;
387 }
388 }
389 freq_prev = freq;
390 targ_prev = targ;
391 }
392 if (pole_found) {
393 newf = qp_energy;
394 }
395 return newf;
396}
397
398boost::optional<double> GW_UKS::SolveQP_FixedPoint(Spin spin, double intercept0,
399 double frequency0,
400 Index gw_level) const {
401 boost::optional<double> newf = boost::none;
402 QPFunc f(gw_level, SigmaEvaluator(spin), intercept0);
404 opt_.g_sc_max_iterations, opt_.g_sc_limit, opt_.qp_solver_alpha);
405 double freq_new = newton.FindRoot(f, frequency0);
406 if (newton.getInfo() == NewtonRapson<QPFunc>::success) {
407 newf = freq_new;
408 }
409 return newf;
410}
411
412double GW_UKS::SolveQP_Bisection(double lowerbound, double f_lowerbound,
413 double upperbound, double f_upperbound,
414 const QPFunc& f) const {
415 if (f_lowerbound * f_upperbound > 0) {
416 throw std::runtime_error(
417 "Bisection needs a positive and negative function value");
418 }
419 while (true) {
420 const double c = 0.5 * (lowerbound + upperbound);
421 if (std::abs(upperbound - lowerbound) < opt_.g_sc_limit) {
422 return c;
423 }
424 const double y_c = f.value(c);
425 if (std::abs(y_c) < opt_.g_sc_limit) {
426 return c;
427 }
428 if (y_c * f_lowerbound > 0) {
429 lowerbound = c;
430 f_lowerbound = y_c;
431 } else {
432 upperbound = c;
433 f_upperbound = y_c;
434 }
435 }
436}
437
438bool GW_UKS::Converged(const Eigen::VectorXd& e1, const Eigen::VectorXd& e2,
439 double epsilon) const {
440 Index state = 0;
441 const double diff_max = (e1 - e2).cwiseAbs().maxCoeff(&state);
442 XTP_LOG(Log::info, log_) << TimeStamp() << " E_diff max=" << diff_max
443 << " StateNo:" << state << std::flush;
444 return diff_max <= epsilon;
445}
446
448 Eigen::VectorXd diag_backup_alpha = Sigma_c_alpha_.diagonal();
449 Eigen::VectorXd diag_backup_beta = Sigma_c_beta_.diagonal();
450 Sigma_c_alpha_ = sigma_alpha_->CalcCorrelationOffDiag(getGWAResultsAlpha());
451 Sigma_c_beta_ = sigma_beta_->CalcCorrelationOffDiag(getGWAResultsBeta());
452 Sigma_c_alpha_.diagonal() = diag_backup_alpha;
453 Sigma_c_beta_.diagonal() = diag_backup_beta;
454}
455
456Eigen::MatrixXd GW_UKS::getHQPAlpha() const {
458 Eigen::MatrixXd(
459 dft_energies_alpha_.segment(opt_.qpmin, qptotal_).asDiagonal());
460}
461Eigen::MatrixXd GW_UKS::getHQPBeta() const {
463 Eigen::MatrixXd(
464 dft_energies_beta_.segment(opt_.qpmin, qptotal_).asDiagonal());
465}
466
467Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd>
469 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> qpdiag(getHQPAlpha());
470 PrintQP_Energies(Spin::Alpha, qpdiag.eigenvalues());
471 return qpdiag;
472}
473
474Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd>
476 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> qpdiag(getHQPBeta());
477 PrintQP_Energies(Spin::Beta, qpdiag.eigenvalues());
478 return qpdiag;
479}
480
482 const Eigen::VectorXd gwa =
484
485 const Index homo = Homo(spin);
486 const Eigen::VectorXd& dft = DftEnergies(spin);
487
489 << (boost::format(" ====== Perturbative %1% quasiparticle energies "
490 "(Hartree) ======") %
491 SpinName(spin))
492 .str()
493 << std::flush;
494
495 if (opt_.qpmin <= homo && homo + 1 <= opt_.qpmax) {
496 const double dft_gap = dft(homo + 1) - dft(homo);
497 const double qp_gap = gwa(homo + 1 - opt_.qpmin) - gwa(homo - opt_.qpmin);
498
500 << (boost::format(" %1% DFT gap = %2$+1.6f Hartree %1% GWA gap = "
501 "%3$+1.6f Hartree Delta = %4$+1.6f Hartree") %
502 SpinName(spin) % dft_gap % qp_gap % (qp_gap - dft_gap))
503 .str()
504 << std::flush;
505 }
506
507 for (Index i = 0; i < qptotal_; i++) {
508 const Index level = i + opt_.qpmin;
510 << LevelLabel(spin, level)
511 << (boost::format(" = %1$4d (%2%) %3% DFT = %4$+1.4f VXC = %5$+1.4f "
512 "S-X = %6$+1.4f S-C = %7$+1.4f GWA = %8$+1.4f") %
513 level % OccupationTag(spin, level) % SpinName(spin) %
514 DftEnergies(spin)(level) % Vxc(spin)(i, i) % SigmaX(spin)(i, i) %
515 SigmaC(spin)(i, i) % gwa(i))
516 .str()
517 << std::flush;
518 }
519}
520
522 const Eigen::VectorXd& qp_diag_energies) const {
523 const Eigen::VectorXd gwa =
525
526 const Index homo = Homo(spin);
527 const Eigen::VectorXd& dft = DftEnergies(spin);
528
529 XTP_LOG(Log::error, log_) << TimeStamp() << " Full " << SpinName(spin)
530 << " quasiparticle Hamiltonian" << std::flush;
531
532 if (opt_.qpmin <= homo && homo + 1 <= opt_.qpmax) {
533 const double dft_gap = dft(homo + 1) - dft(homo);
534 const double pqp_gap = gwa(homo + 1 - opt_.qpmin) - gwa(homo - opt_.qpmin);
535 const double dqp_gap = qp_diag_energies(homo + 1 - opt_.qpmin) -
536 qp_diag_energies(homo - opt_.qpmin);
537
539 << (boost::format(" %1% DFT gap = %2$+1.6f Hartree %1% PQP gap = "
540 "%3$+1.6f Hartree %1% DQP gap = %4$+1.6f Hartree") %
541 SpinName(spin) % dft_gap % pqp_gap % dqp_gap)
542 .str()
543 << std::flush;
544 }
545
546 for (Index i = 0; i < qptotal_; i++) {
547 const Index level = i + opt_.qpmin;
549 << LevelLabel(spin, level)
550 << (boost::format(" = %1$4d (%2%) %3% PQP = %4$+1.6f DQP = %5$+1.6f") %
551 level % OccupationTag(spin, level) % SpinName(spin) % gwa(i) %
552 qp_diag_energies(i))
553 .str()
554 << std::flush;
555 }
556}
557
558} // namespace xtp
559} // 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_uks.h:92
double value(double frequency) const
Definition gw_uks.h:87
boost::optional< double > SolveQP_FixedPoint(Spin spin, double intercept0, double frequency0, Index gw_level) const
Definition gw_uks.cc:398
TCMatrix_gwbse_spin & Mmn_
Definition gw_uks.h:144
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > DiagonalizeQPHamiltonianAlpha() const
Definition gw_uks.cc:468
bool Converged(const Eigen::VectorXd &e1, const Eigen::VectorXd &e2, double epsilon) const
Definition gw_uks.cc:438
void PrintQP_Energies(Spin spin, const Eigen::VectorXd &qp_diag_energies) const
Definition gw_uks.cc:521
Eigen::MatrixXd & SigmaC(Spin spin)
Definition gw_uks.cc:130
const Eigen::VectorXd & DftEnergies(Spin spin) const
Definition gw_uks.cc:121
std::unique_ptr< Sigma_base_UKS > sigma_alpha_
Definition gw_uks.h:150
std::string LevelLabel(Spin spin, Index level) const
Definition gw_uks.cc:152
Eigen::MatrixXd Sigma_x_alpha_
Definition gw_uks.h:152
std::unique_ptr< Sigma_base_UKS > sigma_beta_
Definition gw_uks.h:151
Eigen::MatrixXd Sigma_c_beta_
Definition gw_uks.h:155
Eigen::MatrixXd & SigmaX(Spin spin)
Definition gw_uks.cc:127
Eigen::VectorXd getGWAResultsAlpha() const
Definition gw_uks.cc:292
const Eigen::MatrixXd & vxc_alpha_
Definition gw_uks.h:145
Eigen::VectorXd ScissorShift_DFTlevel(const Eigen::VectorXd &dft_energies, Index homo) const
Definition gw_uks.cc:166
void CalculateHQP()
Definition gw_uks.cc:447
const Eigen::VectorXd & RPAInputEnergiesAlpha() const
Definition gw_uks.cc:302
const char * OccupationTag(Spin spin, Index level) const
Definition gw_uks.cc:162
Sigma_base_UKS & SigmaEvaluator(Spin spin)
Definition gw_uks.cc:139
Eigen::VectorXd SolveQP(Spin spin, const Eigen::VectorXd &frequencies) const
Definition gw_uks.cc:309
boost::optional< double > SolveQP_Linearisation(Spin spin, double intercept0, double frequency0, Index gw_level) const
Definition gw_uks.cc:349
const Eigen::VectorXd & RPAInputEnergiesBeta() const
Definition gw_uks.cc:305
Eigen::MatrixXd getHQPBeta() const
Definition gw_uks.cc:461
boost::optional< double > SolveQP_Grid(Spin spin, double intercept0, double frequency0, Index gw_level) const
Definition gw_uks.cc:365
void PrintGWA_Energies(Spin spin) const
Definition gw_uks.cc:481
Eigen::VectorXd getGWAResultsBeta() const
Definition gw_uks.cc:297
Eigen::MatrixXd getHQPAlpha() const
Definition gw_uks.cc:456
Eigen::MatrixXd Sigma_x_beta_
Definition gw_uks.h:153
const char * SpinName(Spin spin) const
Definition gw_uks.cc:148
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:147
double CalcSpinHomoLumoShift(const Eigen::VectorXd &frequencies, Spin spin) const
Definition gw_uks.cc:174
const Eigen::MatrixXd & vxc_beta_
Definition gw_uks.h:146
const Eigen::VectorXd & dft_energies_beta_
Definition gw_uks.h:148
Logger & log_
Definition gw_uks.h:143
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > DiagonalizeQPHamiltonianBeta() const
Definition gw_uks.cc:475
Eigen::MatrixXd Sigma_c_alpha_
Definition gw_uks.h:154
const Eigen::MatrixXd & Vxc(Spin spin) const
Definition gw_uks.cc:124
void configure(const options &opt)
Definition gw_uks.cc:38
void CalculateGWPerturbation()
Definition gw_uks.cc:184
Index Homo(Spin spin) const
Definition gw_uks.cc:145
double SolveQP_Bisection(double lowerbound, double f_lowerbound, double upperbound, double f_upperbound, const QPFunc &f) const
Definition gw_uks.cc:412
Logger is used for thread-safe output of messages.
Definition logger.h:164
Newton Rapson rootfinder for 1d functions.
double FindRoot(const Func &f, double x0)
void SetSharedPPM(const PPM &ppm)
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
Charge transport classes.
Definition ERIs.h:28
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26