votca 2026-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;
37
38 // Normalize legacy and new grid-search settings once at configuration time.
39 // From this point on, the solver operates only on the decoupled canonical
40 // controls, independent of how the user specified them in XML.
42 qptotal_ = opt_.qpmax - opt_.qpmin + 1;
43 rpa_.configure(opt_.homo, opt_.rpamin, opt_.rpamax);
44 sigma_ = SigmaFactory().Create(opt_.sigma_integration, Mmn_, rpa_);
45 Sigma_base::options sigma_opt;
46 sigma_opt.homo = opt_.homo;
47 sigma_opt.qpmax = opt_.qpmax;
48 sigma_opt.qpmin = opt_.qpmin;
49 sigma_opt.rpamin = opt_.rpamin;
50 sigma_opt.rpamax = opt_.rpamax;
51 sigma_opt.eta = opt_.eta;
52 sigma_opt.alpha = opt_.alpha;
53 sigma_opt.quadrature_scheme = opt_.quadrature_scheme;
54 sigma_opt.order = opt_.order;
55 sigma_->configure(sigma_opt);
56 Sigma_x_ = Eigen::MatrixXd::Zero(qptotal_, qptotal_);
57 Sigma_c_ = Eigen::MatrixXd::Zero(qptotal_, qptotal_);
58}
59
60double GW::CalcHomoLumoShift(Eigen::VectorXd frequencies) const {
61 double DFTgap = dft_energies_(opt_.homo + 1) - dft_energies_(opt_.homo);
62 double QPgap = frequencies(opt_.homo + 1 - opt_.qpmin) -
63 frequencies(opt_.homo - opt_.qpmin);
64 return QPgap - DFTgap;
65}
66
67Eigen::MatrixXd GW::getHQP() const {
68 return Sigma_x_ + Sigma_c_ - vxc_ +
69 Eigen::MatrixXd(
70 dft_energies_.segment(opt_.qpmin, qptotal_).asDiagonal());
71}
72
73Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> GW::DiagonalizeQPHamiltonian()
74 const {
75 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> qpdiag(getHQP());
76 PrintQP_Energies(qpdiag.eigenvalues());
77 return qpdiag;
78}
79
81 Eigen::VectorXd gwa_energies = getGWAResults();
82 double shift = CalcHomoLumoShift(gwa_energies);
83
85 << (boost::format(
86 " ====== Perturbative quasiparticle energies (Hartree) ====== "))
87 .str()
88 << std::flush;
90 << (boost::format(" DeltaHLGap = %1$+1.6f Hartree") % shift).str()
91 << std::flush;
92
93 for (Index i = 0; i < qptotal_; i++) {
94 std::string level = " Level";
95 if ((i + opt_.qpmin) == opt_.homo) {
96 level = " HOMO ";
97 } else if ((i + opt_.qpmin) == opt_.homo + 1) {
98 level = " LUMO ";
99 }
100
102 << level
103 << (boost::format(" = %1$4d DFT = %2$+1.4f VXC = %3$+1.4f S-X = "
104 "%4$+1.4f S-C = %5$+1.4f GWA = %6$+1.4f") %
105 (i + opt_.qpmin) % dft_energies_(i + opt_.qpmin) % vxc_(i, i) %
106 Sigma_x_(i, i) % Sigma_c_(i, i) % gwa_energies(i))
107 .str()
108 << std::flush;
109 }
110 return;
111}
112
113void GW::PrintQP_Energies(const Eigen::VectorXd& qp_diag_energies) const {
114 Eigen::VectorXd gwa_energies = getGWAResults();
116 << TimeStamp() << " Full quasiparticle Hamiltonian " << std::flush;
118 << (boost::format(
119 " ====== Diagonalized quasiparticle energies (Hartree) "
120 "====== "))
121 .str()
122 << std::flush;
123 for (Index i = 0; i < qptotal_; i++) {
124 std::string level = " Level";
125 if ((i + opt_.qpmin) == opt_.homo) {
126 level = " HOMO ";
127 } else if ((i + opt_.qpmin) == opt_.homo + 1) {
128 level = " LUMO ";
129 }
131 << level
132 << (boost::format(" = %1$4d PQP = %2$+1.6f DQP = %3$+1.6f ") %
133 (i + opt_.qpmin) % gwa_energies(i) % qp_diag_energies(i))
134 .str()
135 << std::flush;
136 }
137 return;
138}
139
141 const Eigen::VectorXd& dft_energies) const {
142 Eigen::VectorXd shifted_energies = dft_energies;
143 shifted_energies.segment(opt_.homo + 1, dft_energies.size() - opt_.homo - 1)
144 .array() += opt_.shift;
145 return shifted_energies;
146}
147
149
150 Sigma_x_ = (1 - opt_.ScaHFX) * sigma_->CalcExchangeMatrix();
152 << TimeStamp() << " Calculated Hartree exchange contribution"
153 << std::flush;
154 // dftenergies has size aobasissize
155 // rpaenergies/Mmn have size rpatotal
156 // gwaenergies/frequencies have size qptotal
157 // homo index is relative to dft_energies
159 << TimeStamp() << " Scissor shifting DFT energies by: " << opt_.shift
160 << " Hrt" << std::flush;
161 Eigen::VectorXd dft_shifted_energies = ScissorShift_DFTlevel(dft_energies_);
162 rpa_.setRPAInputEnergies(
163 dft_shifted_energies.segment(opt_.rpamin, opt_.rpamax - opt_.rpamin + 1));
164 Eigen::VectorXd frequencies =
165 dft_shifted_energies.segment(opt_.qpmin, qptotal_);
166
167 Anderson mixing_;
168 mixing_.Configure(opt_.gw_mixing_order, opt_.gw_mixing_alpha);
169
170 for (Index i_gw = 0; i_gw < opt_.gw_sc_max_iterations; ++i_gw) {
171 gw_sc_iteration_ = i_gw;
172 if (i_gw % opt_.reset_3c == 0 && i_gw != 0) {
173 Mmn_.Rebuild();
175 << TimeStamp() << " Rebuilding 3c integrals" << std::flush;
176 }
177 sigma_->PrepareScreening();
179 << TimeStamp() << " Calculated screening via RPA" << std::flush;
181 << TimeStamp() << " Solving QP equations " << std::flush;
182 if (opt_.gw_mixing_order > 0 && i_gw > 0) {
183 mixing_.UpdateInput(frequencies);
184 }
185
186 frequencies = SolveQP(frequencies);
187
188 if (opt_.gw_sc_max_iterations > 1) {
189 Eigen::VectorXd rpa_energies_old = rpa_.getRPAInputEnergies();
190 if (opt_.gw_mixing_order > 0 && i_gw > 0) {
191 if (opt_.gw_mixing_order == 1) {
193 << "GWSC using linear mixing with alpha: " << opt_.gw_mixing_alpha
194 << std::flush;
195 } else {
197 << "GWSC using Anderson mixing with history "
198 << opt_.gw_mixing_order << ", alpha: " << opt_.gw_mixing_alpha
199 << std::flush;
200 }
201 mixing_.UpdateOutput(frequencies);
202 Eigen::VectorXd mixed_frequencies = mixing_.MixHistory();
203 rpa_.UpdateRPAInputEnergies(dft_energies_, mixed_frequencies,
204 opt_.qpmin);
205 frequencies = mixed_frequencies;
206
207 } else {
208 XTP_LOG(Log::debug, log_) << "GWSC using plain update " << std::flush;
209 rpa_.UpdateRPAInputEnergies(dft_energies_, frequencies, opt_.qpmin);
210 }
211
212 for (int i = 0; i < frequencies.size(); i++) {
214 << "... GWSC iter " << i_gw << " state " << i << " "
215 << std::setprecision(9) << frequencies(i) << std::flush;
216 }
217
219 << TimeStamp() << " GW_Iteration:" << i_gw
220 << " Shift[Hrt]:" << CalcHomoLumoShift(frequencies) << std::flush;
221 if (Converged(rpa_.getRPAInputEnergies(), rpa_energies_old,
222 opt_.gw_sc_limit)) {
223 XTP_LOG(Log::info, log_) << TimeStamp() << " Converged after "
224 << i_gw + 1 << " GW iterations." << std::flush;
225 break;
226 } else if (i_gw == opt_.gw_sc_max_iterations - 1) {
228 << TimeStamp()
229 << " WARNING! GW-self-consistency cycle not converged after "
230 << opt_.gw_sc_max_iterations << " iterations." << std::flush;
232 << TimeStamp() << " Run continues. Inspect results carefully!"
233 << std::flush;
234 break;
235 }
236 }
237 }
238 Sigma_c_.diagonal() = sigma_->CalcCorrelationDiag(frequencies);
240}
241
242Eigen::VectorXd GW::getGWAResults() const {
243 return Sigma_x_.diagonal() + Sigma_c_.diagonal() - vxc_.diagonal() +
244 dft_energies_.segment(opt_.qpmin, qptotal_);
245}
246
247Eigen::VectorXd GW::SolveQP(const Eigen::VectorXd& frequencies) const {
248 sigma_->ResetDiagEvalCounter();
249 Eigen::VectorXd env = Eigen::VectorXd::Zero(qptotal_);
250
251 const Eigen::VectorXd intercepts =
252 dft_energies_.segment(opt_.qpmin, qptotal_) + Sigma_x_.diagonal() -
253 vxc_.diagonal();
254
255 Eigen::VectorXd frequencies_new = frequencies;
256 Eigen::Array<bool, Eigen::Dynamic, 1> converged =
257 Eigen::Array<bool, Eigen::Dynamic, 1>::Zero(qptotal_);
258
259 QPStats total_stats;
260
261#ifdef _OPENMP
262 Index use_threads =
264#else
265 Index use_threads = 1;
266#endif
267
268#pragma omp parallel for schedule(dynamic) num_threads(use_threads)
269 for (Index gw_level = 0; gw_level < qptotal_; ++gw_level) {
270
271 double initial_f = frequencies[gw_level];
272 double intercept = intercepts[gw_level];
273 boost::optional<double> newf;
274 QPStats local_stats;
275
276 if (opt_.qp_solver == "fixedpoint") {
277 newf = SolveQP_FixedPoint(intercept, initial_f, gw_level, &local_stats);
278 }
279 if (newf) {
280 frequencies_new[gw_level] = newf.value();
281 converged[gw_level] = true;
282 } else {
283 newf = SolveQP_Grid(intercept, initial_f, gw_level, &local_stats);
284 if (newf) {
285 frequencies_new[gw_level] = newf.value();
286 converged[gw_level] = true;
287 } else {
288 newf =
289 SolveQP_Linearisation(intercept, initial_f, gw_level, &local_stats);
290 if (newf) {
291 frequencies_new[gw_level] = newf.value();
292 }
293 }
294 }
295
296#pragma omp critical
297 {
298 total_stats.Add(local_stats);
299 }
300 }
301
302 if (!converged.all()) {
303 std::vector<Index> states;
304 for (Index s = 0; s < converged.size(); s++) {
305 if (!converged[s]) {
306 states.push_back(s);
307 }
308 }
309 IndexParser rp;
310 XTP_LOG(Log::error, log_) << TimeStamp() << " Not converged PQP states are:"
311 << rp.CreateIndexString(states) << std::flush;
313 << TimeStamp() << " Increase the grid search interval" << std::flush;
314 }
315
317 << " Sigma diagonal evaluations in SolveQP: "
318 << sigma_->GetDiagEvalCounter() << std::flush;
319
320 XTP_LOG(Log::info, log_) << TimeStamp() << " QP diagnostics: "
321 << "scan=" << total_stats.sigma_scan_calls
322 << " refine=" << total_stats.sigma_refine_calls
323 << " deriv_sigma="
324 << total_stats.sigma_derivative_calls
325 << " other=" << total_stats.sigma_other_calls
326 << " total_sigma=" << total_stats.TotalSigmaCalls()
327 << " unique_omega="
328 << total_stats.sigma_unique_frequencies
329 << " repeat_sigma=" << total_stats.sigma_repeat_calls
330 << " deriv_calls=" << total_stats.deriv_calls
331 << std::flush;
332
333 return frequencies_new;
334}
335
336boost::optional<double> GW::SolveQP_Linearisation(double intercept0,
337 double frequency0,
338 Index gw_level,
339 QPStats* stats) const {
340 boost::optional<double> newf = boost::none;
341
342 QPFunc fqp(gw_level, *sigma_.get(), intercept0);
343
344 double sigma = fqp.sigma(frequency0, EvalStage::Other);
345 double dsigma_domega = fqp.deriv(frequency0);
346 double Z = 1.0 - dsigma_domega;
347 if (std::abs(Z) > 1e-9) {
348 newf = frequency0 + (intercept0 - frequency0 + sigma) / Z;
349 }
350
351 if (stats != nullptr) {
352 *stats = fqp.GetStats();
353 }
354 return newf;
355}
356
358 double intercept0, double frequency0, Index gw_level, double left_limit,
359 double right_limit, bool allow_rejected_return, QPStats* stats) const {
360
361 QPFunc fqp(gw_level, *sigma_.get(), intercept0);
362
363 qp_solver::SolverOptions solver_opt;
364
365 // Pass only canonical search controls into the shared grid-search utility.
366 // RKS and UKS intentionally use the same search semantics; only the sigma
367 // evaluator differs between the two implementations.
368 solver_opt.g_sc_limit = opt_.g_sc_limit;
369 solver_opt.qp_bisection_max_iter = opt_.g_sc_max_iterations;
370 solver_opt.qp_full_window_half_width = opt_.qp_full_window_half_width;
371 solver_opt.qp_dense_spacing = opt_.qp_dense_spacing;
372 solver_opt.qp_adaptive_shell_width = opt_.qp_adaptive_shell_width;
373 solver_opt.qp_adaptive_shell_count = opt_.qp_adaptive_shell_count;
374
376 std::vector<QPRootCandidate> accepted_roots;
377 std::vector<QPRootCandidate> rejected_roots;
378
379 bool use_brent = false;
380 if (opt_.qp_root_finder == "brent") {
381 use_brent = true;
382 }
383
385 fqp, frequency0, left_limit, right_limit, gw_sc_iteration_, solver_opt,
386 &wdiag, &accepted_roots, &rejected_roots, use_brent);
387
389#pragma omp critical
390 {
391 if (!accepted_roots.empty() || !rejected_roots.empty()) {
393 << " Adaptive scan qplevel:" << gw_level << " center="
394 << std::max(left_limit, std::min(right_limit, frequency0))
395 << " shells=" << wdiag.shells_explored
396 << " first_interval_shell=" << wdiag.first_interval_shell
397 << " first_accepted_shell=" << wdiag.first_accepted_shell
398 << " intervals=" << wdiag.intervals_found << std::flush;
399 }
400 }
401 }
402
403 if (stats != nullptr) {
404 *stats = fqp.GetStats();
405 }
406
407 if (!accepted_roots.empty()) {
408 return result;
409 }
410
411 if (!rejected_roots.empty() && !allow_rejected_return) {
413#pragma omp critical
414 {
416 << " Adaptive scan qplevel:" << gw_level
417 << " produced only rejected roots in [" << left_limit << ", "
418 << right_limit << "], forcing wider retry" << std::flush;
419 }
420 }
421 return boost::none;
422 }
423
424 return result;
425}
426
427boost::optional<double> GW::SolveQP_Grid_Windowed_Dense(
428 double intercept0, double frequency0, Index gw_level, double left_limit,
429 double right_limit, bool allow_rejected_return, QPStats* stats) const {
430
431 QPFunc fqp(gw_level, *sigma_.get(), intercept0);
432
433 qp_solver::SolverOptions solver_opt;
434
435 // The dense path performs a uniform sign-change scan over the requested
436 // interval. Its spacing is qp_dense_spacing, not the adaptive shell width.
437 solver_opt.g_sc_limit = opt_.g_sc_limit;
438 solver_opt.qp_bisection_max_iter = opt_.g_sc_max_iterations;
439 solver_opt.qp_full_window_half_width = opt_.qp_full_window_half_width;
440 solver_opt.qp_dense_spacing = opt_.qp_dense_spacing;
441 solver_opt.qp_adaptive_shell_width = opt_.qp_adaptive_shell_width;
442 solver_opt.qp_adaptive_shell_count = opt_.qp_adaptive_shell_count;
443
444 const bool use_brent = (opt_.qp_root_finder == "brent");
445
446 std::vector<QPRootCandidate> accepted_roots;
447 std::vector<QPRootCandidate> rejected_roots;
448 std::vector<std::pair<double, double>> roots; // omega : Z
449
450 if (left_limit < right_limit) {
451 double freq_prev = left_limit;
452 double targ_prev = fqp.value(freq_prev, EvalStage::Scan);
453
454 // Dense scanning is the robustness path: it uniformly samples the entire
455 // interval to avoid missing sign changes that a coarser adaptive shell
456 // search might step over.
457 const Index n_steps = std::max<Index>(
458 2, static_cast<Index>(
459 std::ceil((right_limit - left_limit) / opt_.qp_dense_spacing)) +
460 1);
461
462 for (Index i_node = 1; i_node < n_steps; ++i_node) {
463 const double freq =
464 (i_node == n_steps - 1)
465 ? right_limit
466 : std::min(right_limit, left_limit + static_cast<double>(i_node) *
467 opt_.qp_dense_spacing);
468
469 const double targ = fqp.value(freq, EvalStage::Scan);
470
471 if (targ_prev * targ < 0.0) {
472 auto cand =
473 qp_solver::RefineQPInterval(freq_prev, targ_prev, freq, targ, fqp,
474 frequency0, solver_opt, use_brent);
475 if (cand) {
476 roots.emplace_back(cand->omega, cand->Z);
477 if (cand->accepted) {
478 accepted_roots.push_back(*cand);
479 } else {
480 rejected_roots.push_back(*cand);
481 }
482 }
483 }
484
485 freq_prev = freq;
486 targ_prev = targ;
487 }
488 }
489
491#pragma omp critical
492 {
493 if (accepted_roots.empty() && rejected_roots.empty()) {
495 << " Dense scan qplevel:" << gw_level << " found no roots in ["
496 << left_limit << ", " << right_limit << "]" << std::flush;
497 } else {
499 << " Dense scan qplevel:" << gw_level << " roots (omega:Z)\n\t\t";
500 for (const auto& root : roots) {
501 XTP_LOG(Log::info, log_) << std::setprecision(5) << root.first << ":"
502 << root.second << " ";
503 }
504 XTP_LOG(Log::info, log_) << std::flush;
505 }
506 }
507 }
508
509 if (stats != nullptr) {
510 *stats = fqp.GetStats();
511 }
512
513 if (!accepted_roots.empty()) {
514 auto best = std::max_element(
515 accepted_roots.begin(), accepted_roots.end(),
516 [](const QPRootCandidate& a, const QPRootCandidate& b) {
517 return qp_solver::ScoreRoot(a) < qp_solver::ScoreRoot(b);
518 });
519 return best->omega;
520 }
521
522 if (!rejected_roots.empty()) {
523 if (!allow_rejected_return) {
525#pragma omp critical
526 {
528 << " Dense scan qplevel:" << gw_level
529 << " produced only rejected roots in [" << left_limit << ", "
530 << right_limit << "], forcing wider retry" << std::flush;
531 }
532 }
533 return boost::none;
534 }
535
536 auto least_bad = std::max_element(
537 rejected_roots.begin(), rejected_roots.end(),
538 [](const QPRootCandidate& a, const QPRootCandidate& b) {
539 return qp_solver::ScoreRoot(a) < qp_solver::ScoreRoot(b);
540 });
541 return least_bad->omega;
542 }
543
544 return boost::none;
545}
546
547boost::optional<double> GW::SolveQP_Grid_Windowed(
548 double intercept0, double frequency0, Index gw_level, double left_limit,
549 double right_limit, bool allow_rejected_return, QPStats* stats) const {
550
551 if (opt_.qp_grid_search_mode == "adaptive") {
552 return SolveQP_Grid_Windowed_Adaptive(intercept0, frequency0, gw_level,
553 left_limit, right_limit,
554 allow_rejected_return, stats);
555 }
556
557 if (opt_.qp_grid_search_mode == "dense") {
558 return SolveQP_Grid_Windowed_Dense(intercept0, frequency0, gw_level,
559 left_limit, right_limit,
560 allow_rejected_return, stats);
561 }
562
563 if (opt_.qp_grid_search_mode == "adaptive_with_dense_fallback") {
564 QPStats total_stats;
565 auto adaptive = SolveQP_Grid_Windowed_Adaptive(
566 intercept0, frequency0, gw_level, left_limit, right_limit,
567 allow_rejected_return, &total_stats);
568 if (adaptive) {
569 if (stats != nullptr) {
570 *stats = total_stats;
571 }
572 return adaptive;
573 }
574
575 QPStats dense_stats;
576 auto dense = SolveQP_Grid_Windowed_Dense(
577 intercept0, frequency0, gw_level, left_limit, right_limit,
578 allow_rejected_return, &dense_stats);
579 total_stats.Add(dense_stats);
580
582#pragma omp critical
583 {
585 << " Adaptive QP scan failed for qplevel:" << gw_level
586 << ", retrying dense grid scan in [" << left_limit << ", "
587 << right_limit << "]" << std::flush;
588 }
589 }
590
591 if (stats != nullptr) {
592 *stats = total_stats;
593 }
594 return dense;
595 }
596
597 throw std::runtime_error("Unknown gw.qp_grid_search_mode '" +
598 opt_.qp_grid_search_mode + "'");
599}
600
601boost::optional<double> GW::SolveQP_Grid(double intercept0, double frequency0,
602 Index gw_level, QPStats* stats) const {
603 // The full QP search window is now controlled explicitly and no longer
604 // inferred from the dense-grid spacing.
605 const double range = opt_.qp_full_window_half_width;
606
607 const double full_left_limit = frequency0 - range;
608 const double full_right_limit = frequency0 + range;
609
610 double restricted_left_limit = full_left_limit;
611 double restricted_right_limit = full_right_limit;
612
613 bool use_restricted_window = false;
614
615 if (opt_.qp_restrict_search) {
616 const Index mo_level = gw_level + opt_.qpmin;
617 const bool is_occupied = (mo_level <= opt_.homo);
618
619 if (is_occupied) {
620 restricted_right_limit = std::min(full_right_limit, -opt_.qp_zero_margin);
621 } else {
622 restricted_left_limit =
623 std::max(full_left_limit, opt_.qp_virtual_min_energy);
624 }
625
626 const double tol = 1e-12;
627 use_restricted_window =
628 (std::abs(restricted_left_limit - full_left_limit) > tol) ||
629 (std::abs(restricted_right_limit - full_right_limit) > tol);
630 }
631
632 if (use_restricted_window && restricted_left_limit < restricted_right_limit) {
633 auto restricted =
634 SolveQP_Grid_Windowed(intercept0, frequency0, gw_level,
635 restricted_left_limit, restricted_right_limit,
636 false, // accepted roots only
637 stats);
638 if (restricted) {
639 return restricted;
640 }
641
643#pragma omp critical
644 {
646 << " Restricted QP search failed for qplevel:" << gw_level
647 << " in window [" << restricted_left_limit << ", "
648 << restricted_right_limit << "], forcing full dense window ["
649 << full_left_limit << ", " << full_right_limit << "]" << std::flush;
650 }
651 }
652
653 // Important: do NOT go back to adaptive here.
654 // The restricted window already failed to produce an accepted root.
655 // Force the old robust full-window dense search.
656 return SolveQP_Grid_Windowed_Dense(intercept0, frequency0, gw_level,
657 full_left_limit, full_right_limit, true,
658 stats);
659 }
660
661 return SolveQP_Grid_Windowed(intercept0, frequency0, gw_level,
662 full_left_limit, full_right_limit, true, stats);
663}
664
665boost::optional<double> GW::SolveQP_FixedPoint(double intercept0,
666 double frequency0,
667 Index gw_level,
668 QPStats* stats) const {
669 boost::optional<double> newf = boost::none;
670 QPFunc f(gw_level, *sigma_.get(), intercept0);
672 opt_.g_sc_max_iterations, opt_.g_sc_limit, opt_.qp_solver_alpha);
673 double freq_new = newton.FindRoot(f, frequency0);
674 if (newton.getInfo() == NewtonRapson<QPFunc>::success) {
675 newf = freq_new;
676 }
677 if (stats != nullptr) {
678 *stats = f.GetStats();
679 }
680 return newf;
681}
682
683bool GW::Converged(const Eigen::VectorXd& e1, const Eigen::VectorXd& e2,
684 double epsilon) const {
685 Index state = 0;
686 bool energies_converged = true;
687 double diff_max = (e1 - e2).cwiseAbs().maxCoeff(&state);
688 if (diff_max > epsilon) {
689 energies_converged = false;
690 }
691 XTP_LOG(Log::info, log_) << TimeStamp() << " E_diff max=" << diff_max
692 << " StateNo:" << state << std::flush;
693 return energies_converged;
694}
695
697 Eigen::VectorXd diag_backup = Sigma_c_.diagonal();
698 Sigma_c_ = sigma_->CalcCorrelationOffDiag(getGWAResults());
699 Sigma_c_.diagonal() = diag_backup;
700}
701
702void GW::PlotSigma(std::string filename, Index steps, double spacing,
703 std::string states) const {
704
705 Eigen::VectorXd frequencies =
706 rpa_.getRPAInputEnergies().segment(opt_.qpmin - opt_.rpamin, qptotal_);
707
708 std::vector<Index> state_inds;
709 IndexParser rp;
710 std::vector<Index> parsed_states = rp.CreateIndexVector(states);
711 for (Index gw_level : parsed_states) {
712 if (gw_level >= opt_.qpmin && gw_level <= opt_.qpmax) {
713 state_inds.push_back(gw_level);
714 }
715 }
717 << TimeStamp() << " PQP(omega) written to '" << filename
718 << "' for states " << rp.CreateIndexString(state_inds) << std::flush;
719
720 const Index num_states = state_inds.size();
721
722 const Eigen::VectorXd intercept =
723 dft_energies_.segment(opt_.qpmin, qptotal_) + Sigma_x_.diagonal() -
724 vxc_.diagonal();
725 Eigen::MatrixXd mat = Eigen::MatrixXd::Zero(steps, 2 * num_states);
726#pragma omp parallel for schedule(dynamic)
727 for (Index grid_point = 0; grid_point < steps; grid_point++) {
728 const double offset =
729 ((double)grid_point - ((double)(steps - 1) / 2.0)) * spacing;
730 for (Index i = 0; i < num_states; i++) {
731 const Index gw_level = state_inds[i];
732 const double omega = frequencies(gw_level) + offset;
733 double sigma = sigma_->CalcCorrelationDiagElement(gw_level, omega);
734 mat(grid_point, 2 * i) = omega;
735 mat(grid_point, 2 * i + 1) = sigma + intercept[gw_level];
736 }
737 }
738
739 std::ofstream out;
740 out.open(filename);
741 for (Index i = 0; i < num_states; i++) {
742 const Index gw_level = state_inds[i];
743 out << boost::format("#%1$somega_%2$d\tE_QP(omega)_%2$d") %
744 (i == 0 ? "" : "\t") % gw_level;
745 }
746 out << std::endl;
747 boost::format numFormat("%+1.6f");
748 Eigen::IOFormat matFormat(Eigen::StreamPrecision, 0, "\t", "\n");
749 out << numFormat % mat.format(matFormat) << std::endl;
750 out.close();
751}
752
753} // namespace xtp
754} // 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 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 QPStats & GetStats() const
Definition gw.h:181
double deriv(double frequency) const
Definition gw.h:174
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
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::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
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
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
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)
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
boost::optional< double > SolveQP_Grid_Windowed(QPFunc &fqp, double frequency0, double left_limit, double right_limit, Index gw_sc_iteration, const SolverOptions &opt, WindowDiagnostics *wdiag=nullptr, std::vector< RootCandidate > *accepted_roots_out=nullptr, std::vector< RootCandidate > *rejected_roots_out=nullptr, bool use_brent=false)
void NormalizeGridSearchOptions(Opt &opt)
boost::optional< RootCandidate > RefineQPInterval(double lowerbound, double f_lowerbound, double upperbound, double f_upperbound, const QPFunc &f, double reference, const SolverOptions &opt, bool use_brent)
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
static Level current_level
Definition globals.h:30
void Add(const Stats &other)
std::size_t TotalSigmaCalls() const