votca 2026-dev
Loading...
Searching...
No Matches
gwbse.cc
Go to the documentation of this file.
1
2
3/*
4 * Copyright 2009-2023 The VOTCA Development Team
5 * (http://www.votca.org)
6 *
7 * Licensed under the Apache License, Version 2.0 (the "License")
8 *
9 * You may not use this file except in compliance with the License.
10 * You may obtain a copy of the License at
11 *
12 * http://www.apache.org/licenses/LICENSE-2.0
13 *
14 * Unless required by applicable law or agreed to in writing, software
15 * distributed under the License is distributed on an "AS IS" BASIS,
16 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
17 * See the License for the specific language governing permissions and
18 * limitations under the License.
19 *
20 */
21
22// Third party includes
23#include <boost/algorithm/string.hpp>
24#include <boost/format.hpp>
25
26// VOTCA includes
27#include <stdexcept>
29
30// Local VOTCA includes
31#include "votca/xtp/basisset.h"
32#include "votca/xtp/bse.h"
34#include "votca/xtp/gwbse.h"
35#include "votca/xtp/logger.h"
37#include "votca/xtp/orbitals.h"
38#include "votca/xtp/rpa_uks.h"
39#include "votca/xtp/vxc_grid.h"
41
42using std::flush;
43namespace votca {
44namespace xtp {
45
47 Index ignored_corelevels = 0;
48 if (!orbitals_.hasECPName()) {
49 ECPBasisSet basis;
50 basis.Load("corelevels");
51 Index coreElectrons = 0;
52 for (const auto& atom : orbitals_.QMAtoms()) {
53 coreElectrons += basis.getElement(atom.getElement()).getNcore();
54 }
55 ignored_corelevels = coreElectrons / 2;
56 }
57 return ignored_corelevels;
58}
59
61
62 // getting level ranges
63 Index rpamax = 0;
64 Index rpamin = 0; // never changes
65 Index qpmin = 0;
66 Index qpmax = 0;
67 Index bse_vmin = 0;
68 Index bse_cmax = 0;
69
70 const bool is_uks = orbitals_.hasUnrestrictedOrbitals();
71 Index homo = orbitals_.getHomo(); // indexed from 0
72 if (is_uks) {
73 homo = std::max(orbitals_.getHomoAlpha(), orbitals_.getHomoBeta());
74 }
75 Index num_of_levels = orbitals_.getBasisSetSize();
76 Index num_of_occlevels = orbitals_.getNumberOfAlphaElectrons();
77 if (is_uks) {
78 num_of_occlevels = std::max(orbitals_.getNumberOfAlphaElectrons(),
79 orbitals_.getNumberOfBetaElectrons());
80 }
81
82 std::string ranges = options.get(".ranges").as<std::string>();
83
84 // now check validity, and get rpa, qp, and bse level ranges accordingly
85
86 if (ranges == "factor") {
87
88 double rpamaxfactor = options.get(".rpamax").as<double>();
89 rpamax = Index(rpamaxfactor * double(num_of_levels)) -
90 1; // total number of levels
91
92 double qpminfactor = options.get(".qpmin").as<double>();
93 qpmin =
94 num_of_occlevels - Index(qpminfactor * double(num_of_occlevels)) - 1;
95
96 double qpmaxfactor = options.get(".qpmax").as<double>();
97 qpmax =
98 num_of_occlevels + Index(qpmaxfactor * double(num_of_occlevels)) - 1;
99
100 double bseminfactor = options.get(".bsemin").as<double>();
101 bse_vmin =
102 num_of_occlevels - Index(bseminfactor * double(num_of_occlevels)) - 1;
103
104 double bsemaxfactor = options.get(".bsemax").as<double>();
105 bse_cmax =
106 num_of_occlevels + Index(bsemaxfactor * double(num_of_occlevels)) - 1;
107
108 } else if (ranges == "explicit") {
109 // get explicit numbers
110 rpamax = options.get(".rpamax").as<Index>();
111 qpmin = options.get(".qpmin").as<Index>();
112 qpmax = options.get(".qpmax").as<Index>();
113 bse_vmin = options.get(".bsemin").as<Index>();
114 bse_cmax = options.get(".bsemax").as<Index>();
115 } else if (ranges == "default") {
116 rpamax = num_of_levels - 1;
117 qpmin = 0;
118 qpmax = 3 * homo + 1;
119 bse_vmin = 0;
120 bse_cmax = 3 * homo + 1;
121 } else if (ranges == "full") {
122 rpamax = num_of_levels - 1;
123 qpmin = 0;
124 qpmax = num_of_levels - 1;
125 bse_vmin = 0;
126 bse_cmax = num_of_levels - 1;
127 }
128 std::string ignore_corelevels =
129 options.get(".ignore_corelevels").as<std::string>();
130
131 if (ignore_corelevels == "RPA" || ignore_corelevels == "GW" ||
132 ignore_corelevels == "BSE") {
133 Index ignored_corelevels = CountCoreLevels();
134 if (ignore_corelevels == "RPA") {
135 rpamin = ignored_corelevels;
136 }
137 if (ignore_corelevels == "GW" || ignore_corelevels == "RPA") {
138 if (qpmin < ignored_corelevels) {
139 qpmin = ignored_corelevels;
140 }
141 }
142 if (ignore_corelevels == "GW" || ignore_corelevels == "RPA" ||
143 ignore_corelevels == "BSE") {
144 if (bse_vmin < ignored_corelevels) {
145 bse_vmin = ignored_corelevels;
146 }
147 }
148
150 << TimeStamp() << " Ignoring " << ignored_corelevels
151 << " core levels for " << ignore_corelevels << " and beyond." << flush;
152 }
153
154 // check maximum and minimum sizes
155 if (rpamax >= num_of_levels) {
156 rpamax = num_of_levels - 1;
157 }
158 if (qpmax >= num_of_levels) {
159 qpmax = num_of_levels - 1;
160 }
161 if (bse_cmax >= num_of_levels) {
162 bse_cmax = num_of_levels - 1;
163 }
164 if (bse_vmin < 0) {
165 bse_vmin = 0;
166 }
167 if (qpmin < 0) {
168 qpmin = 0;
169 }
170
171 Index bse_vmax = homo;
172 Index bse_cmin = homo + 1;
173
174 if (rpamin < 0 || rpamax < 0 || rpamin > rpamax) {
175 throw std::runtime_error("Invalid RPA level range after setup/clamping.");
176 }
177 if (qpmin < 0 || qpmax < 0 || qpmin > qpmax) {
178 throw std::runtime_error("Invalid GW level range after setup/clamping.");
179 }
180 if (bse_vmin < 0 || bse_vmax < 0 || bse_vmin > bse_vmax) {
181 throw std::runtime_error(
182 "Invalid BSE occupied level range after setup/clamping.");
183 }
184 if (bse_cmin < 0 || bse_cmax < 0 || bse_cmin > bse_cmax) {
185 throw std::runtime_error(
186 "Invalid BSE virtual level range after setup/clamping.");
187 }
188 if (bse_vmax >= num_of_levels || bse_cmax >= num_of_levels) {
189 throw std::runtime_error(
190 "BSE level range exceeds available orbital indices.");
191 }
192 if (bse_vmax >= bse_cmin) {
193 throw std::runtime_error("BSE occupied and virtual level ranges overlap.");
194 }
195
196 gwopt_.homo = homo;
197 gwopt_.qpmin = qpmin;
198 gwopt_.qpmax = qpmax;
199 gwopt_.rpamin = rpamin;
200 gwopt_.rpamax = rpamax;
201
202 bseopt_.vmin = bse_vmin;
203 bseopt_.cmax = bse_cmax;
204 bseopt_.homo = homo;
205 bseopt_.qpmin = qpmin;
206 bseopt_.qpmax = qpmax;
207 bseopt_.rpamin = rpamin;
208 bseopt_.rpamax = rpamax;
209
210 orbitals_.setRPAindices(rpamin, rpamax);
211 orbitals_.setGWindices(qpmin, qpmax);
212 orbitals_.setBSEindices(bse_vmin, bse_cmax);
213
214 // Index bse_vmax = homo;
215 // Index bse_cmin = homo + 1;
216 Index bse_vtotal = bse_vmax - bse_vmin + 1;
217 Index bse_ctotal = bse_cmax - bse_cmin + 1;
218 Index bse_size = bse_vtotal * bse_ctotal;
219
220 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " RPA level range [" << rpamin
221 << ":" << rpamax << "]" << flush;
222 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " GW level range [" << qpmin
223 << ":" << qpmax << "]" << flush;
225 << TimeStamp() << " BSE level range occ[" << bse_vmin << ":" << bse_vmax
226 << "] virt[" << bse_cmin << ":" << bse_cmax << "]" << flush;
227
228 gwopt_.reset_3c = options.get(".gw.rebuild_3c_freq").as<Index>();
229
230 bseopt_.nmax = options.get(".bse.exctotal").as<Index>();
231 if (bseopt_.nmax > bse_size || bseopt_.nmax < 0) {
232 bseopt_.nmax = bse_size;
233 }
234
235 bseopt_.davidson_correction =
236 options.get("bse.davidson.correction").as<std::string>();
237
238 bseopt_.davidson_tolerance =
239 options.get("bse.davidson.tolerance").as<std::string>();
240
241 bseopt_.davidson_update =
242 options.get("bse.davidson.update").as<std::string>();
243
244 bseopt_.davidson_maxiter = options.get("bse.davidson.maxiter").as<Index>();
245
246 bseopt_.useTDA = options.get("bse.useTDA").as<bool>();
247 orbitals_.setTDAApprox(bseopt_.useTDA);
248 if (!bseopt_.useTDA) {
249 XTP_LOG(Log::error, *pLog_) << " BSE type: full" << flush;
250 } else {
251 XTP_LOG(Log::error, *pLog_) << " BSE type: TDA" << flush;
252 }
253
254 Index full_bse_size = (bseopt_.useTDA) ? bse_size : 2 * bse_size;
255 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " BSE Hamiltonian has size "
256 << full_bse_size << "x" << full_bse_size << flush;
257
258 bseopt_.use_Hqp_offdiag = options.get("bse.use_Hqp_offdiag").as<bool>();
259 orbitals_.SetFlagUseHqpOffdiag(bseopt_.use_Hqp_offdiag);
260
261 if (!bseopt_.use_Hqp_offdiag) {
263 << " BSE without Hqp offdiagonal elements" << flush;
264 } else {
266 << " BSE with Hqp offdiagonal elements" << flush;
267 }
268
269 bseopt_.max_dyn_iter = options.get("bse.dyn_screen_max_iter").as<Index>();
270 bseopt_.dyn_tolerance = options.get("bse.dyn_screen_tol").as<double>();
271 if (bseopt_.max_dyn_iter > 0) {
273 }
274
275 functional_ = orbitals_.getXCFunctionalName();
276 grid_ = orbitals_.getXCGrid();
277
278 dftbasis_name_ = orbitals_.getDFTbasisName();
279 if (orbitals_.hasAuxbasisName()) {
280 auxbasis_name_ = orbitals_.getAuxbasisName();
281 } else if (options.exists(".auxbasisset")) {
282 auxbasis_name_ = options.get(".auxbasisset").as<std::string>();
283 } else {
285 try {
286 BasisSet b;
288 } catch (std::runtime_error&) {
289 std::runtime_error(
290 "There is no auxbasis from the dftcalculation nor did you specify an "
291 "auxbasisset for the gwbse calculation. Also no auxiliary basisset "
292 "for basisset " +
293 dftbasis_name_ + " could be found!");
294 }
296 << " Could not find an auxbasisset using " << auxbasis_name_ << flush;
297 }
298
299 std::string mode = options.get("gw.mode").as<std::string>();
300 if (mode == "G0W0") {
301 gwopt_.gw_sc_max_iterations = 1;
302 } else if (mode == "evGW") {
303 gwopt_.g_sc_limit = 0.1 * gwopt_.gw_sc_limit;
304 gwopt_.eta = 0.1;
305 }
306
307 XTP_LOG(Log::error, *pLog_) << " Running GW as: " << mode << flush;
308 gwopt_.ScaHFX = orbitals_.getScaHFX();
309
310 gwopt_.shift = options.get("gw.scissor_shift").as<double>();
311 gwopt_.g_sc_limit =
312 options.get("gw.qp_sc_limit").as<double>(); // convergence criteria
313 // for qp iteration
314 // [Hartree]]
315 gwopt_.g_sc_max_iterations =
316 options.get("gw.qp_sc_max_iter").as<Index>(); // convergence
317 // criteria for qp
318 // iteration
319 // [Hartree]]
320
321 if (mode == "evGW") {
322 gwopt_.gw_sc_max_iterations = options.get("gw.sc_max_iter").as<Index>();
323 }
324
325 gwopt_.gw_sc_limit =
326 options.get("gw.sc_limit").as<double>(); // convergence criteria
327 // for shift it
329 << " qp_sc_limit [Hartree]: " << gwopt_.g_sc_limit << flush;
330 if (gwopt_.gw_sc_max_iterations > 1) {
332 << " gw_sc_limit [Hartree]: " << gwopt_.gw_sc_limit << flush;
333 }
334 bseopt_.min_print_weight = options.get("bse.print_weight").as<double>();
335 // print exciton WF composition weight larger than minimum
336
337 // possible tasks
338 std::string tasks_string = options.get(".tasks").as<std::string>();
339 boost::algorithm::to_lower(tasks_string);
340
341 if (tasks_string.find("gw") != std::string::npos) {
342 do_gw_ = true;
343 }
344
345 if (is_uks) {
346 if (tasks_string.find("excitons") != std::string::npos ||
347 tasks_string.find("exciton_uks") != std::string::npos) {
348 do_gw_ = true;
349 do_bse_exciton_uks_ = true;
350 }
351
352 if (tasks_string.find("singlets") != std::string::npos ||
353 tasks_string.find("triplets") != std::string::npos) {
354 throw std::runtime_error(
355 "Invalid gwbse task for UKS reference: 'singlets' and 'triplets' "
356 "are not defined for open-shell systems.\n"
357 "Use 'exciton_uks' (or 'excitons') instead.");
358 }
359 } else {
360 if (tasks_string.find("all") != std::string::npos) {
361 if (is_uks) {
362 do_gw_ = true;
363 do_bse_exciton_uks_ = true;
364 } else {
365 do_gw_ = true;
366 do_bse_singlets_ = true;
367 do_bse_triplets_ = true;
368 }
369 }
370 if (tasks_string.find("singlets") != std::string::npos) {
371 do_bse_singlets_ = true;
372 }
373 if (tasks_string.find("triplets") != std::string::npos) {
374 do_bse_triplets_ = true;
375 }
376 }
377
378 XTP_LOG(Log::error, *pLog_) << " Tasks: " << flush;
379 if (do_gw_) {
380 XTP_LOG(Log::error, *pLog_) << " GW " << flush;
381 }
382 if (do_bse_singlets_) {
383 XTP_LOG(Log::error, *pLog_) << " singlets " << flush;
384 }
385 if (do_bse_triplets_) {
386 XTP_LOG(Log::error, *pLog_) << " triplets " << flush;
387 }
389 XTP_LOG(Log::error, *pLog_) << " exciton_uks " << flush;
390 }
391
392 if (options.exists("bse.fragments")) {
393 std::vector<tools::Property*> prop_region =
394 options.Select("bse.fragments.fragment");
395 Index index = 0;
396 for (tools::Property* prop : prop_region) {
397 std::string indices = prop->get("indices").as<std::string>();
398 fragments_.push_back(QMFragment<BSE_Population>(index, indices));
399 index++;
400 }
401 }
402
403 gwopt_.sigma_integration =
404 options.get("gw.sigma_integrator").as<std::string>();
406 << " Sigma integration: " << gwopt_.sigma_integration << flush;
407 gwopt_.eta = options.get("gw.eta").as<double>();
408 XTP_LOG(Log::error, *pLog_) << " eta: " << gwopt_.eta << flush;
409 if (gwopt_.sigma_integration == "exact") {
411 << " RPA Hamiltonian size: " << (homo + 1 - rpamin) * (rpamax - homo)
412 << flush;
413 }
414 if (gwopt_.sigma_integration == "cda") {
415 gwopt_.order = options.get("gw.quadrature_order").as<Index>();
417 << " Quadrature integration order : " << gwopt_.order << flush;
418 gwopt_.quadrature_scheme =
419 options.get("gw.quadrature_scheme").as<std::string>();
421 << " Quadrature integration scheme : " << gwopt_.quadrature_scheme
422 << flush;
423 gwopt_.alpha = options.get("gw.alpha").as<double>();
425 << " Alpha smoothing parameter : " << gwopt_.alpha << flush;
426 }
427 gwopt_.qp_solver = options.get("gw.qp_solver").as<std::string>();
428
429 XTP_LOG(Log::error, *pLog_) << " QP solver: " << gwopt_.qp_solver << flush;
430 if (gwopt_.qp_solver == "grid") {
431 // Grid-based QP solving now exposes the full search extent, dense scan,
432 // and adaptive shell scan as separate controls. This avoids the old
433 // overloading of qp_grid_steps / qp_grid_spacing, where one pair of values
434 // simultaneously determined the window width and both scan resolutions.
435 //
436 // New-style options are read first. Deprecated legacy aliases are then
437 // translated only for values that are still unset.
438 // New decoupled QP search options
439 if (options.exists("gw.qp_full_window_half_width")) {
440 gwopt_.qp_full_window_half_width =
441 options.get("gw.qp_full_window_half_width").as<double>();
442 }
443 if (options.exists("gw.qp_dense_spacing")) {
444 gwopt_.qp_dense_spacing = options.get("gw.qp_dense_spacing").as<double>();
445 }
446 if (options.exists("gw.qp_adaptive_shell_width")) {
447 gwopt_.qp_adaptive_shell_width =
448 options.get("gw.qp_adaptive_shell_width").as<double>();
449 }
450 if (options.exists("gw.qp_adaptive_shell_count")) {
451 gwopt_.qp_adaptive_shell_count =
452 options.get("gw.qp_adaptive_shell_count").as<Index>();
453 }
454
455 // Deprecated legacy aliases. They are accepted for backward
456 // compatibility and mapped onto the canonical controls below.
457 const bool has_legacy_steps = options.exists("gw.qp_grid_steps");
458 const bool has_legacy_spacing = options.exists("gw.qp_grid_spacing");
459
460 if (has_legacy_steps != has_legacy_spacing) {
461 throw std::runtime_error(
462 "Deprecated gw.qp_grid_steps and gw.qp_grid_spacing must be given "
463 "together if either is used.");
464 }
465
466 if (has_legacy_steps && has_legacy_spacing) {
467 gwopt_.qp_grid_steps = options.get("gw.qp_grid_steps").as<Index>();
468 gwopt_.qp_grid_spacing = options.get("gw.qp_grid_spacing").as<double>();
469
470 const double legacy_half_width =
472 const double legacy_shell_width =
474
475 if (gwopt_.qp_full_window_half_width <= 0.0) {
476 gwopt_.qp_full_window_half_width = legacy_half_width;
477 }
478 if (gwopt_.qp_dense_spacing <= 0.0) {
479 gwopt_.qp_dense_spacing = gwopt_.qp_grid_spacing;
480 }
481 if (gwopt_.qp_adaptive_shell_count <= 0 &&
482 gwopt_.qp_adaptive_shell_width <= 0.0) {
483 gwopt_.qp_adaptive_shell_width = legacy_shell_width;
484 }
485
487 << " Deprecated GW QP options detected: "
488 << "qp_grid_steps=" << gwopt_.qp_grid_steps
489 << " qp_grid_spacing=" << gwopt_.qp_grid_spacing
490 << " -> mapped to qp_full_window_half_width="
491 << gwopt_.qp_full_window_half_width
492 << " qp_dense_spacing=" << gwopt_.qp_dense_spacing
493 << " qp_adaptive_shell_width=" << gwopt_.qp_adaptive_shell_width
494 << flush;
495 }
496
497 // Final normalization produces one canonical option set used by both
498 // GW and GW_UKS implementations.
500
502 << " QP full window half-width: " << gwopt_.qp_full_window_half_width
503 << flush;
505 << " QP dense spacing: " << gwopt_.qp_dense_spacing << flush;
507 << " QP adaptive shell width: " << gwopt_.qp_adaptive_shell_width
508 << flush;
510 << " QP adaptive shell count: " << gwopt_.qp_adaptive_shell_count
511 << flush;
512 }
513 gwopt_.qp_root_finder = options.get("gw.qp_root_finder").as<std::string>();
515 << " QP root finder: "
516 << gwopt_.qp_root_finder << flush;
517
518
519 gwopt_.gw_mixing_order =
520 options.get("gw.mixing_order").as<Index>(); // max history in
521 // mixing (0: plain,
522 // 1: linear, >1
523 // Anderson)
524
525 gwopt_.gw_mixing_alpha = options.get("gw.mixing_alpha").as<double>();
526 gwopt_.qp_grid_search_mode =
527 options.get("gw.qp_grid_search_mode").as<std::string>();
528 gwopt_.qp_restrict_search = options.get("gw.qp_restrict_search").as<bool>();
529 gwopt_.qp_zero_margin = options.get("gw.qp_zero_margin").as<double>();
530 gwopt_.qp_virtual_min_energy =
531 options.get("gw.qp_virtual_min_energy").as<double>();
532 if (mode == "evGW") {
533 if (gwopt_.gw_mixing_order == 0) {
534 XTP_LOG(Log::error, *pLog_) << " evGW with plain update " << std::flush;
535 } else if (gwopt_.gw_mixing_order == 1) {
536 XTP_LOG(Log::error, *pLog_) << " evGW with linear update using alpha "
537 << gwopt_.gw_mixing_alpha << std::flush;
538 } else {
539 XTP_LOG(Log::error, *pLog_) << " evGW with Anderson update with history "
540 << gwopt_.gw_mixing_order << " using alpha "
541 << gwopt_.gw_mixing_alpha << std::flush;
542 }
543 }
544
545 if (options.exists("gw.sigma_plot")) {
546 sigma_plot_states_ = options.get("gw.sigma_plot.states").as<std::string>();
547 sigma_plot_steps_ = options.get("gw.sigma_plot.steps").as<Index>();
548 sigma_plot_spacing_ = options.get("gw.sigma_plot.spacing").as<double>();
550 options.get("gw.sigma_plot.filename").as<std::string>();
551
553 << " Sigma plot states: " << sigma_plot_states_ << flush;
555 << " Sigma plot steps: " << sigma_plot_steps_ << flush;
557 << " Sigma plot spacing: " << sigma_plot_spacing_ << flush;
559 << " Sigma plot filename: " << sigma_plot_filename_ << flush;
560 }
561}
562
564
565 const double hrt2ev = tools::conv::hrt2ev;
566 tools::Property& gwbse_summary = summary.add("GWBSE", "");
567 if (do_gw_) {
568 gwbse_summary.setAttribute("units", "eV");
569 gwbse_summary.setAttribute(
570 "DFTEnergy",
571 (boost::format("%1$+1.6f ") % (orbitals_.getDFTTotalEnergy() * hrt2ev))
572 .str());
573
574 if (orbitals_.hasUnrestrictedOrbitals()) {
575 gwbse_summary.add("dft_alpha", "");
576 gwbse_summary.add("dft_beta", "");
577
578 tools::Property& alpha_summary = gwbse_summary.get("dft_alpha");
579 tools::Property& beta_summary = gwbse_summary.get("dft_beta");
580
581 alpha_summary.setAttribute("HOMO", orbitals_.getHomoAlpha());
582 alpha_summary.setAttribute("LUMO", orbitals_.getLumoAlpha());
583 beta_summary.setAttribute("HOMO", orbitals_.getHomoBeta());
584 beta_summary.setAttribute("LUMO", orbitals_.getLumoBeta());
585 for (Index state = 0; state < gwopt_.qpmax + 1 - gwopt_.qpmin; state++) {
586 tools::Property& level_alpha = alpha_summary.add("level", "");
587 level_alpha.setAttribute("number", state + gwopt_.qpmin);
588 level_alpha.add(
589 "dft_energy",
590 (boost::format("%1$+1.6f ") %
591 (orbitals_.MOs().eigenvalues()(state + gwopt_.qpmin) * hrt2ev))
592 .str());
593 level_alpha.add("gw_energy",
594 (boost::format("%1$+1.6f ") %
595 (orbitals_.QPpertEnergiesAlpha()(state) * hrt2ev))
596 .str());
597 level_alpha.add(
598 "qp_energy",
599 (boost::format("%1$+1.6f ") %
600 (orbitals_.QPdiagAlpha().eigenvalues()(state) * hrt2ev))
601 .str());
602
603 tools::Property& level_beta = beta_summary.add("level", "");
604 level_beta.setAttribute("number", state + gwopt_.qpmin);
605 level_beta.add(
606 "dft_energy",
607 (boost::format("%1$+1.6f ") %
608 (orbitals_.MOs_beta().eigenvalues()(state + gwopt_.qpmin) *
609 hrt2ev))
610 .str());
611 level_beta.add("gw_energy",
612 (boost::format("%1$+1.6f ") %
613 (orbitals_.QPpertEnergiesBeta()(state) * hrt2ev))
614 .str());
615 level_beta.add("qp_energy",
616 (boost::format("%1$+1.6f ") %
617 (orbitals_.QPdiagBeta().eigenvalues()(state) * hrt2ev))
618 .str());
619 }
620 } else {
621 tools::Property& dft_summary = gwbse_summary.add("dft", "");
622 dft_summary.setAttribute("HOMO", gwopt_.homo);
623 dft_summary.setAttribute("LUMO", gwopt_.homo + 1);
624
625 for (Index state = 0; state < gwopt_.qpmax + 1 - gwopt_.qpmin; state++) {
626 tools::Property& level_summary = dft_summary.add("level", "");
627 level_summary.setAttribute("number", state + gwopt_.qpmin);
628 level_summary.add(
629 "dft_energy",
630 (boost::format("%1$+1.6f ") %
631 (orbitals_.MOs().eigenvalues()(state + gwopt_.qpmin) * hrt2ev))
632 .str());
633 level_summary.add("gw_energy",
634 (boost::format("%1$+1.6f ") %
635 (orbitals_.QPpertEnergies()(state) * hrt2ev))
636 .str());
637 level_summary.add("qp_energy",
638 (boost::format("%1$+1.6f ") %
639 (orbitals_.QPdiag().eigenvalues()(state) * hrt2ev))
640 .str());
641 }
642 }
643 }
644 if (do_bse_singlets_) {
645 tools::Property& singlet_summary = gwbse_summary.add("singlets", "");
646 for (Index state = 0; state < bseopt_.nmax; ++state) {
647 tools::Property& level_summary = singlet_summary.add("level", "");
648 level_summary.setAttribute("number", state + 1);
649 level_summary.add(
650 "omega", (boost::format("%1$+1.6f ") %
651 (orbitals_.BSESinglets().eigenvalues()(state) * hrt2ev))
652 .str());
653 if (orbitals_.hasTransitionDipoles()) {
654
655 const Eigen::Vector3d& dipoles = (orbitals_.TransitionDipoles())[state];
656 double f = 2 * dipoles.squaredNorm() *
657 orbitals_.BSESinglets().eigenvalues()(state) / 3.0;
658
659 level_summary.add("f", (boost::format("%1$+1.6f ") % f).str());
660 tools::Property& dipol_summary = level_summary.add(
661 "Trdipole", (boost::format("%1$+1.4f %2$+1.4f %3$+1.4f") %
662 dipoles.x() % dipoles.y() % dipoles.z())
663 .str());
664 dipol_summary.setAttribute("unit", "e*bohr");
665 dipol_summary.setAttribute("gauge", "length");
666 }
667 }
668 }
669 if (do_bse_triplets_) {
670 tools::Property& triplet_summary = gwbse_summary.add("triplets", "");
671 for (Index state = 0; state < bseopt_.nmax; ++state) {
672 tools::Property& level_summary = triplet_summary.add("level", "");
673 level_summary.setAttribute("number", state + 1);
674 level_summary.add(
675 "omega", (boost::format("%1$+1.6f ") %
676 (orbitals_.BSETriplets().eigenvalues()(state) * hrt2ev))
677 .str());
678 }
679 }
681 tools::Property& exciton_uks_summary = gwbse_summary.add("exciton_uks", "");
682
683 const bool has_dipoles = orbitals_.hasTransitionDipoles();
684 Eigen::VectorXd oscs = Eigen::VectorXd::Zero(0);
685 if (has_dipoles) {
686 oscs =
687 orbitals_.Oscillatorstrengths(QMStateType(QMStateType::ExcitonUKS));
688 }
689
690 for (Index state = 0;
691 state <
692 std::min<Index>(bseopt_.nmax, orbitals_.BSEUKS().eigenvalues().size());
693 ++state) {
694 tools::Property& level_summary = exciton_uks_summary.add("level", "");
695 level_summary.setAttribute("number", state + 1);
696 level_summary.add("omega",
697 (boost::format("%1$+1.6f ") %
698 (orbitals_.BSEUKS().eigenvalues()(state) * hrt2ev))
699 .str());
700 const Index ndip =
701 static_cast<Index>(orbitals_.TransitionDipoles().size());
702 const Index nosc = static_cast<Index>(oscs.size());
703
704 if (has_dipoles && state < ndip && state < nosc) {
705 const Eigen::Vector3d& dipoles = orbitals_.TransitionDipoles()[state];
706
707 level_summary.add("f",
708 (boost::format("%1$+1.6f ") % oscs(state)).str());
709
710 tools::Property& dipol_summary = level_summary.add(
711 "Trdipole", (boost::format("%1$+1.4f %2$+1.4f %3$+1.4f") %
712 dipoles.x() % dipoles.y() % dipoles.z())
713 .str());
714 dipol_summary.setAttribute("unit", "e*bohr");
715 dipol_summary.setAttribute("gauge", "length");
716 }
717 }
718 }
719
720 return;
721}
722
723/*
724 * Many-body Green's fuctions theory implementation
725 *
726 * data required from orbitals file
727 * - atomic coordinates
728 * - DFT molecular orbitals (energies and coeffcients)
729 * - DFT exchange-correlation potential matrix in atomic orbitals
730 * - number of electrons, number of levels
731 */
732
733Eigen::MatrixXd GWBSE::CalculateVXC(const AOBasis& dftbasis) {
734 if (orbitals_.getXCFunctionalName().empty()) {
735 orbitals_.setXCFunctionalName(functional_);
736 } else {
737 if (!(functional_ == orbitals_.getXCFunctionalName())) {
738 throw std::runtime_error("Functionals from DFT " +
739 orbitals_.getXCFunctionalName() + " GWBSE " +
740 functional_ + " differ!");
741 }
742 }
743
744 Vxc_Grid grid;
745 grid.GridSetup(grid_, orbitals_.QMAtoms(), dftbasis);
747 << TimeStamp() << " Setup grid for integration with gridsize: " << grid_
748 << " with " << grid.getGridSize() << " points, divided into "
749 << grid.getBoxesSize() << " boxes" << flush;
750 Vxc_Potential<Vxc_Grid> vxcpotential(grid);
751 vxcpotential.setXCfunctional(functional_);
753 << TimeStamp() << " Integrating Vxc with functional " << functional_
754 << flush;
755 Eigen::MatrixXd DMAT = orbitals_.DensityMatrixGroundState();
756 Mat_p_Energy e_vxc_ao = vxcpotential.IntegrateVXC(DMAT);
757 XTP_LOG(Log::info, *pLog_) << TimeStamp() << " Calculated Vxc" << flush;
759 << TimeStamp() << " Set hybrid exchange factor: " << orbitals_.getScaHFX()
760 << flush;
761 Index qptotal = gwopt_.qpmax - gwopt_.qpmin + 1;
762 Eigen::MatrixXd mos =
763 orbitals_.MOs().eigenvectors().middleCols(gwopt_.qpmin, qptotal);
764
765 Eigen::MatrixXd vxc = mos.transpose() * e_vxc_ao.matrix() * mos;
767 << TimeStamp() << " Calculated exchange-correlation expectation values "
768 << flush;
769
770 return vxc;
771}
772
773std::pair<Eigen::MatrixXd, Eigen::MatrixXd> GWBSE::CalculateVXCSpinResolved(
774 const AOBasis& dftbasis) {
775 if (orbitals_.getXCFunctionalName().empty()) {
776 orbitals_.setXCFunctionalName(functional_);
777 } else if (!(functional_ == orbitals_.getXCFunctionalName())) {
778 throw std::runtime_error("Functionals from DFT " +
779 orbitals_.getXCFunctionalName() + " GWBSE " +
780 functional_ + " differ!");
781 }
782
783 Vxc_Grid grid;
784 grid.GridSetup(grid_, orbitals_.QMAtoms(), dftbasis);
785 Vxc_Potential<Vxc_Grid> vxcpotential(grid);
786 vxcpotential.setXCfunctional(functional_);
787
788 auto dmats = orbitals_.DensityMatrixGroundStateSpinResolved();
789 auto e_vxc_ao = vxcpotential.IntegrateVXCSpin(dmats[0], dmats[1]);
790
791 Index qptotal = gwopt_.qpmax - gwopt_.qpmin + 1;
792 Eigen::MatrixXd mos_alpha =
793 orbitals_.MOs().eigenvectors().middleCols(gwopt_.qpmin, qptotal);
794 Eigen::MatrixXd mos_beta =
795 orbitals_.MOs_beta().eigenvectors().middleCols(gwopt_.qpmin, qptotal);
796
797 Eigen::MatrixXd vxc_alpha =
798 mos_alpha.transpose() * e_vxc_ao.vxc_alpha * mos_alpha;
799 Eigen::MatrixXd vxc_beta =
800 mos_beta.transpose() * e_vxc_ao.vxc_beta * mos_beta;
801
802 return std::make_pair(vxc_alpha, vxc_beta);
803}
804
806
807 // set the parallelization
808 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " Using "
809 << OPENMP::getMaxThreads() << " threads" << flush;
810
811 if (XTP_HAS_MKL_OVERLOAD()) {
813 << TimeStamp() << " Using MKL overload for Eigen " << flush;
814 } else {
816 << TimeStamp()
817 << " Using native Eigen implementation, no BLAS overload " << flush;
818 }
819 Index nogpus = OpenMP_CUDA::UsingGPUs();
820 if (nogpus > 0) {
822 << TimeStamp() << " Using CUDA support for tensor multiplication with "
823 << nogpus << " GPUs." << flush;
824 }
825
827 << TimeStamp() << " Molecule Coordinates [A] " << flush;
828 for (QMAtom& atom : orbitals_.QMAtoms()) {
829 std::string output = (boost::format("%5d"
830 "%5s"
831 " %1.4f %1.4f %1.4f") %
832 atom.getId() % atom.getElement() %
833 (atom.getPos().x() * tools::conv::bohr2ang) %
834 (atom.getPos().y() * tools::conv::bohr2ang) %
835 (atom.getPos().z() * tools::conv::bohr2ang))
836 .str();
837
838 XTP_LOG(Log::error, *pLog_) << output << flush;
839 }
840
841 std::string dft_package = orbitals_.getQMpackage();
843 << TimeStamp() << " DFT data was created by " << dft_package << flush;
844
845 BasisSet dftbs;
846 dftbs.Load(dftbasis_name_);
847
849 << TimeStamp() << " Loaded DFT Basis Set " << dftbasis_name_ << flush;
850
851 AOBasis dftbasis = orbitals_.getDftBasis();
852 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " Filled DFT Basis of size "
853 << dftbasis.AOBasisSize() << flush;
855 << TimeStamp() << " Loaded Auxbasis Set " << auxbasis_name_ << flush;
856
857 orbitals_.SetupAuxBasis(auxbasis_name_);
858 AOBasis auxbasis = orbitals_.getAuxBasis();
859 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " Filled Auxbasis of size "
860 << auxbasis.AOBasisSize() << flush;
861
863 fragments_.size() > 0) {
864 for (const auto& frag : fragments_) {
865 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " Fragment " << frag.getId()
866 << " size:" << frag.size() << flush;
867 }
868 }
869
870 if (!do_gw_) {
871 if ((!orbitals_.hasUnrestrictedOrbitals() && !orbitals_.hasQPdiag()) ||
872 (orbitals_.hasUnrestrictedOrbitals() &&
873 (!orbitals_.hasQPdiagAlpha() || !orbitals_.hasQPdiagBeta()))) {
874 throw std::runtime_error(
875 "You want no GW calculation but the orb file has no matching QP "
876 "coefficients.");
877 }
878 }
879 const bool is_uks = orbitals_.hasUnrestrictedOrbitals();
880 TCMatrix_gwbse Mmn;
881 TCMatrix_gwbse_spin Mmn_spin;
882 // rpamin here, because RPA needs till rpamin
883 Index max_3c = std::max(bseopt_.cmax, gwopt_.qpmax);
884 if (is_uks) {
885 Mmn_spin.alpha.Initialize(auxbasis.AOBasisSize(), gwopt_.rpamin, max_3c,
886 gwopt_.rpamin, gwopt_.rpamax);
887 Mmn_spin.beta.Initialize(auxbasis.AOBasisSize(), gwopt_.rpamin, max_3c,
888 gwopt_.rpamin, gwopt_.rpamax);
890 << TimeStamp() << " Calculating alpha spin Mmn " << flush;
891 Mmn_spin.alpha.Fill(auxbasis, dftbasis, orbitals_.MOs().eigenvectors());
893 << TimeStamp() << " Calculating beta spin Mmn " << flush;
894 Mmn_spin.beta.Fill(auxbasis, dftbasis, orbitals_.MOs_beta().eigenvectors());
895 } else {
896 Mmn.Initialize(auxbasis.AOBasisSize(), gwopt_.rpamin, max_3c, gwopt_.rpamin,
897 gwopt_.rpamax);
899 << TimeStamp() << " Calculating Mmn (3-center-repulsion x orbitals) "
900 << flush;
901 Mmn.Fill(auxbasis, dftbasis, orbitals_.MOs().eigenvectors());
903 << TimeStamp() << " Removed " << Mmn.Removedfunctions()
904 << " functions from Aux Coulomb matrix to avoid near linear "
905 "dependencies"
906 << flush;
908 << TimeStamp() << " Calculated Mmn (3-center-repulsion x orbitals) "
909 << flush;
910 }
911
912 Eigen::MatrixXd Hqp;
913 Eigen::MatrixXd Hqp_alpha;
914 Eigen::MatrixXd Hqp_beta;
915 if (do_gw_) {
916
917 std::chrono::time_point<std::chrono::system_clock> start =
918 std::chrono::system_clock::now();
919 if (is_uks) {
920 auto vxc = CalculateVXCSpinResolved(dftbasis);
921 GW_UKS::options gwopt_uks;
922 gwopt_uks.homo_alpha = orbitals_.getHomoAlpha();
923 gwopt_uks.homo_beta = orbitals_.getHomoBeta();
924 gwopt_uks.qpmin = gwopt_.qpmin;
925 gwopt_uks.qpmax = gwopt_.qpmax;
926 gwopt_uks.rpamin = gwopt_.rpamin;
927 gwopt_uks.rpamax = gwopt_.rpamax;
928 gwopt_uks.eta = gwopt_.eta;
929 gwopt_uks.g_sc_limit = gwopt_.g_sc_limit;
930 gwopt_uks.g_sc_max_iterations = gwopt_.g_sc_max_iterations;
931 gwopt_uks.gw_sc_limit = gwopt_.gw_sc_limit;
932 gwopt_uks.gw_sc_max_iterations = gwopt_.gw_sc_max_iterations;
933 gwopt_uks.shift = gwopt_.shift;
934 gwopt_uks.ScaHFX = gwopt_.ScaHFX;
935 gwopt_uks.sigma_integration = gwopt_.sigma_integration;
936 gwopt_uks.reset_3c = gwopt_.reset_3c;
937 gwopt_uks.qp_solver = gwopt_.qp_solver;
938 gwopt_uks.qp_solver_alpha = gwopt_.qp_solver_alpha;
939 gwopt_uks.qp_grid_steps = gwopt_.qp_grid_steps;
940 gwopt_uks.qp_grid_spacing = gwopt_.qp_grid_spacing;
941
942 // Forward the already-normalized canonical grid-search controls so that
943 // RKS and UKS use the same search design and only differ in the
944 // spin-dependent sigma evaluation.
945 gwopt_uks.qp_full_window_half_width = gwopt_.qp_full_window_half_width;
946 gwopt_uks.qp_dense_spacing = gwopt_.qp_dense_spacing;
947 gwopt_uks.qp_adaptive_shell_width = gwopt_.qp_adaptive_shell_width;
948 gwopt_uks.qp_adaptive_shell_count = gwopt_.qp_adaptive_shell_count;
949 gwopt_uks.gw_mixing_order = gwopt_.gw_mixing_order;
950 gwopt_uks.gw_mixing_alpha = gwopt_.gw_mixing_alpha;
951 gwopt_uks.quadrature_scheme = gwopt_.quadrature_scheme;
952 gwopt_uks.order = gwopt_.order;
953 gwopt_uks.alpha = gwopt_.alpha;
954 gwopt_uks.qp_restrict_search = gwopt_.qp_restrict_search;
955 gwopt_uks.qp_zero_margin = gwopt_.qp_zero_margin;
956 gwopt_uks.qp_virtual_min_energy = gwopt_.qp_virtual_min_energy;
957 gwopt_uks.qp_root_finder = gwopt_.qp_root_finder;
958 gwopt_uks.qp_grid_search_mode = gwopt_.qp_grid_search_mode;
959
960 GW_UKS gw = GW_UKS(*pLog_, Mmn_spin, vxc.first, vxc.second,
961 orbitals_.MOs().eigenvalues(),
962 orbitals_.MOs_beta().eigenvalues());
963 gw.configure(gwopt_uks);
965 orbitals_.QPpertEnergiesAlpha() = gw.getGWAResultsAlpha();
966 orbitals_.QPpertEnergiesBeta() = gw.getGWAResultsBeta();
967 orbitals_.RPAInputEnergiesAlpha() = gw.RPAInputEnergiesAlpha();
968 orbitals_.RPAInputEnergiesBeta() = gw.RPAInputEnergiesBeta();
969 gw.CalculateHQP();
970 Hqp_alpha = gw.getHQPAlpha();
971 Hqp_beta = gw.getHQPBeta();
972 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es_alpha =
974 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es_beta =
976 orbitals_.QPdiagAlpha().eigenvectors() = es_alpha.eigenvectors();
977 orbitals_.QPdiagAlpha().eigenvalues() = es_alpha.eigenvalues();
978 orbitals_.QPdiagBeta().eigenvectors() = es_beta.eigenvectors();
979 orbitals_.QPdiagBeta().eigenvalues() = es_beta.eigenvalues();
980 std::chrono::duration<double> elapsed_time =
981 std::chrono::system_clock::now() - start;
983 << TimeStamp() << " UKS GW calculation took " << elapsed_time.count()
984 << " seconds." << flush;
985 } else {
986 Eigen::MatrixXd vxc = CalculateVXC(dftbasis);
987 GW gw = GW(*pLog_, Mmn, vxc, orbitals_.MOs().eigenvalues());
988 gw.configure(gwopt_);
989
991
992 if (!sigma_plot_states_.empty()) {
995 }
996
997 orbitals_.QPpertEnergies() = gw.getGWAResults();
998 orbitals_.RPAInputEnergies() = gw.RPAInputEnergies();
999
1001 << TimeStamp() << " Calculating offdiagonal part of Sigma " << flush;
1002 gw.CalculateHQP();
1004 << TimeStamp() << " Calculated offdiagonal part of Sigma " << flush;
1005
1006 Hqp = gw.getHQP();
1007
1008 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es =
1010 if (es.info() == Eigen::ComputationInfo::Success) {
1012 << TimeStamp() << " Diagonalized QP Hamiltonian " << flush;
1013 }
1014
1015 orbitals_.QPdiag().eigenvectors() = es.eigenvectors();
1016 orbitals_.QPdiag().eigenvalues() = es.eigenvalues();
1017 std::chrono::duration<double> elapsed_time =
1018 std::chrono::system_clock::now() - start;
1020 << TimeStamp() << " GW calculation took " << elapsed_time.count()
1021 << " seconds." << flush;
1022 }
1023
1024 } else {
1025 if (orbitals_.getGWAmax() != gwopt_.qpmax ||
1026 orbitals_.getGWAmin() != gwopt_.qpmin ||
1027 orbitals_.getRPAmax() != gwopt_.rpamax ||
1028 orbitals_.getRPAmin() != gwopt_.rpamin) {
1029 throw std::runtime_error(
1030 "The ranges for GW and RPA do not agree with the ranges from the "
1031 ".orb file, rerun your GW calculation");
1032 }
1033
1034 if (is_uks) {
1035 const Eigen::MatrixXd& qpcoeff_alpha =
1036 orbitals_.QPdiagAlpha().eigenvectors();
1037 const Eigen::MatrixXd& qpcoeff_beta =
1038 orbitals_.QPdiagBeta().eigenvectors();
1039
1040 Hqp_alpha = qpcoeff_alpha *
1041 orbitals_.QPdiagAlpha().eigenvalues().asDiagonal() *
1042 qpcoeff_alpha.transpose();
1043 Hqp_beta = qpcoeff_beta *
1044 orbitals_.QPdiagBeta().eigenvalues().asDiagonal() *
1045 qpcoeff_beta.transpose();
1046 } else {
1047 const Eigen::MatrixXd& qpcoeff = orbitals_.QPdiag().eigenvectors();
1048
1049 Hqp = qpcoeff * orbitals_.QPdiag().eigenvalues().asDiagonal() *
1050 qpcoeff.transpose();
1051 }
1052 }
1053
1054 // proceed only if BSE requested
1056
1057 std::chrono::time_point<std::chrono::system_clock> start =
1058 std::chrono::system_clock::now();
1059
1060 if (is_uks) {
1061
1062 // Build a single spin-summed dielectric screening from the unrestricted
1063 // RPA response and reuse it for both exciton spin channels.
1064 RPA_UKS rpa_uks_bse(*pLog_, Mmn_spin);
1065 rpa_uks_bse.configure(orbitals_.getHomoAlpha(), orbitals_.getHomoBeta(),
1066 bseopt_.rpamin, bseopt_.rpamax);
1067 rpa_uks_bse.setRPAInputEnergies(orbitals_.RPAInputEnergiesAlpha(),
1068 orbitals_.RPAInputEnergiesBeta());
1069
1070 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es_bse(
1071 rpa_uks_bse.calculate_epsilon_r(0.0));
1072
1073 Eigen::VectorXd epsilon_0_inv_bse =
1074 Eigen::VectorXd::Zero(es_bse.eigenvalues().size());
1075 for (Index i = 0; i < es_bse.eigenvalues().size(); ++i) {
1076 if (es_bse.eigenvalues()(i) > 1e-8) {
1077 epsilon_0_inv_bse(i) = 1.0 / es_bse.eigenvalues()(i);
1078 }
1079 }
1080
1081 if (do_bse_exciton_uks_) {
1082 BSE_UKS::options bseopt_uks;
1083 bseopt_uks.useTDA = bseopt_.useTDA;
1084 bseopt_uks.rpamin = bseopt_.rpamin;
1085 bseopt_uks.rpamax = bseopt_.rpamax;
1086 bseopt_uks.qpmin = bseopt_.qpmin;
1087 bseopt_uks.qpmax = bseopt_.qpmax;
1088 bseopt_uks.vmin = bseopt_.vmin;
1089 bseopt_uks.cmax = bseopt_.cmax;
1090 bseopt_uks.nmax = bseopt_.nmax;
1091 bseopt_uks.davidson_correction = bseopt_.davidson_correction;
1092 bseopt_uks.davidson_tolerance = bseopt_.davidson_tolerance;
1093 bseopt_uks.davidson_update = bseopt_.davidson_update;
1094 bseopt_uks.davidson_maxiter = bseopt_.davidson_maxiter;
1095 bseopt_uks.min_print_weight = bseopt_.min_print_weight;
1096 bseopt_uks.use_Hqp_offdiag = bseopt_.use_Hqp_offdiag;
1097 bseopt_uks.max_dyn_iter = bseopt_.max_dyn_iter;
1098 bseopt_uks.dyn_tolerance = bseopt_.dyn_tolerance;
1099
1100 BSE_UKS bse_uks(*pLog_, Mmn_spin);
1102 bseopt_uks, orbitals_.getHomoAlpha(), orbitals_.getHomoBeta(),
1103 orbitals_.RPAInputEnergiesAlpha(), orbitals_.RPAInputEnergiesBeta(),
1104 Hqp_alpha, Hqp_beta, epsilon_0_inv_bse, es_bse.eigenvectors());
1105
1108 << TimeStamp() << " Solved combined UKS BSE exciton problem "
1109 << flush;
1111
1114 }
1115 }
1116 } else {
1117 BSE bse = BSE(*pLog_, Mmn);
1118 bse.configure(bseopt_, orbitals_.RPAInputEnergies(), Hqp);
1119
1120 if (do_bse_triplets_) {
1123 << TimeStamp() << " Solved BSE for triplets " << flush;
1125 }
1126
1127 if (do_bse_singlets_) {
1130 << TimeStamp() << " Solved BSE for singlets " << flush;
1132 }
1133
1135 if (do_bse_triplets_) {
1137 orbitals_);
1138 }
1139
1140 if (do_bse_singlets_) {
1142 orbitals_);
1143 }
1144 }
1145 }
1146
1147 std::chrono::duration<double> elapsed_time =
1148 std::chrono::system_clock::now() - start;
1149 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " BSE calculation took "
1150 << elapsed_time.count() << " seconds." << flush;
1151 }
1153 << TimeStamp() << " GWBSE calculation finished " << flush;
1154 return true;
1155}
1156
1157} // namespace xtp
1158} // namespace votca
class to manage program options with xml serialization functionality
Definition property.h:55
Property & add(const std::string &key, const std::string &value)
add a new property to structure
Definition property.cc:108
Property & get(const std::string &key)
get existing property
Definition property.cc:79
bool exists(const std::string &key) const
check whether property exists
Definition property.cc:122
T as() const
return value as type
Definition property.h:283
std::vector< Property * > Select(const std::string &filter)
select property based on a filter
Definition property.cc:185
void setAttribute(const std::string &attribute, const T &value)
set an attribute
Definition property.h:314
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
Index AOBasisSize() const
Definition aobasis.h:46
void configure_with_precomputed_screening(const options &opt, Index homo_alpha, Index homo_beta, const Eigen::VectorXd &RPAInputEnergiesAlpha, const Eigen::VectorXd &RPAInputEnergiesBeta, const Eigen::MatrixXd &Hqp_alpha_in, const Eigen::MatrixXd &Hqp_beta_in, const Eigen::VectorXd &epsilon_0_inv, const Eigen::MatrixXd &epsilon_eigenvectors)
Definition bse_uks.cc:51
void Perturbative_DynamicalScreening(Orbitals &orb)
Definition bse_uks.cc:640
void Solve_excitons_uks(Orbitals &orb) const
Definition bse_uks.cc:442
void Analyze_excitons_uks(std::vector< QMFragment< BSE_Population > > fragments, const Orbitals &orb) const
Definition bse_uks.cc:597
void Solve_singlets(Orbitals &orb) const
Definition bse.cc:224
void Analyze_singlets(std::vector< QMFragment< BSE_Population > > fragments, const Orbitals &orb) const
Definition bse.cc:394
void Perturbative_DynamicalScreening(const QMStateType &type, Orbitals &orb)
Definition bse.cc:610
void Analyze_triplets(std::vector< QMFragment< BSE_Population > > fragments, const Orbitals &orb) const
Definition bse.cc:448
void configure(const options &opt, const Eigen::VectorXd &RPAEnergies, const Eigen::MatrixXd &Hqp_in)
Definition bse.cc:43
void Solve_triplets(Orbitals &orb) const
Definition bse.cc:234
void Load(const std::string &name)
Definition basisset.cc:149
void Load(const std::string &name)
const ECPElement & getElement(std::string element_type) const
Index getNcore() const
Definition ecpbasisset.h:90
bool Evaluate()
Definition gwbse.cc:805
bool do_bse_exciton_uks_
Definition gwbse.h:84
std::pair< Eigen::MatrixXd, Eigen::MatrixXd > CalculateVXCSpinResolved(const AOBasis &dftbasis)
Definition gwbse.cc:773
std::vector< QMFragment< BSE_Population > > fragments_
Definition gwbse.h:102
std::string sigma_plot_states_
Definition gwbse.h:93
BSE::options bseopt_
Definition gwbse.h:91
void addoutput(tools::Property &summary)
Definition gwbse.cc:563
bool do_bse_triplets_
Definition gwbse.h:82
std::string grid_
Definition gwbse.h:88
std::string auxbasis_name_
Definition gwbse.h:99
bool do_bse_singlets_
Definition gwbse.h:81
Orbitals & orbitals_
Definition gwbse.h:77
Eigen::MatrixXd CalculateVXC(const AOBasis &dftbasis)
Definition gwbse.cc:733
Index sigma_plot_steps_
Definition gwbse.h:94
std::string functional_
Definition gwbse.h:87
std::string dftbasis_name_
Definition gwbse.h:100
void Initialize(tools::Property &options)
Definition gwbse.cc:60
Index CountCoreLevels()
Definition gwbse.cc:46
GW::options gwopt_
Definition gwbse.h:90
double sigma_plot_spacing_
Definition gwbse.h:95
Logger * pLog_
Definition gwbse.h:76
std::string sigma_plot_filename_
Definition gwbse.h:96
bool do_dynamical_screening_bse_
Definition gwbse.h:83
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > DiagonalizeQPHamiltonianAlpha() const
Definition gw_uks.cc:779
Eigen::VectorXd getGWAResultsAlpha() const
Definition gw_uks.cc:297
void CalculateHQP()
Definition gw_uks.cc:758
const Eigen::VectorXd & RPAInputEnergiesAlpha() const
Definition gw_uks.cc:307
const Eigen::VectorXd & RPAInputEnergiesBeta() const
Definition gw_uks.cc:310
Eigen::MatrixXd getHQPBeta() const
Definition gw_uks.cc:772
Eigen::VectorXd getGWAResultsBeta() const
Definition gw_uks.cc:302
Eigen::MatrixXd getHQPAlpha() const
Definition gw_uks.cc:767
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > DiagonalizeQPHamiltonianBeta() const
Definition gw_uks.cc:786
void configure(const options &opt)
Definition gw_uks.cc:38
void CalculateGWPerturbation()
Definition gw_uks.cc:188
void PlotSigma(std::string filename, Index steps, double spacing, std::string states) const
Definition gw.cc:702
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > DiagonalizeQPHamiltonian() const
Definition gw.cc:73
void configure(const options &opt)
Definition gw.cc:35
Eigen::MatrixXd getHQP() const
Definition gw.cc:67
void CalculateHQP()
Definition gw.cc:696
Eigen::VectorXd getGWAResults() const
Definition gw.cc:242
Eigen::VectorXd RPAInputEnergies() const
Definition gw.h:120
void CalculateGWPerturbation()
Definition gw.cc:148
Eigen::MatrixXd & matrix()
Definition eigen.h:80
static Index UsingGPUs()
container for QM atoms
Definition qmatom.h:37
Unrestricted RPA helper for spin-resolved GW screening.
Definition rpa_uks.h:68
void setRPAInputEnergies(const Eigen::VectorXd &rpaenergies_alpha, const Eigen::VectorXd &rpaenergies_beta)
Set alpha and beta RPA input energies directly.
Definition rpa_uks.h:158
void configure(Index homo_alpha, Index homo_beta, Index rpamin, Index rpamax)
Configure orbital window and spin-resolved HOMO indices.
Definition rpa_uks.h:91
Eigen::MatrixXd calculate_epsilon_r(double frequency) const
Dielectric matrix on the real frequency axis.
Definition rpa_uks.h:121
void Fill(const AOBasis &auxbasis, const AOBasis &dftbasis, const Eigen::MatrixXd &dft_orbitals)
void Initialize(Index basissize, Index mmin, Index mmax, Index nmin, Index nmax)
Index Removedfunctions() const
Definition threecenter.h:45
Timestamp returns the current time as a string Example: cout << TimeStamp()
Definition logger.h:224
Index getGridSize() const
Definition vxc_grid.h:41
void GridSetup(const std::string &type, const QMMolecule &atoms, const AOBasis &basis)
Definition vxc_grid.cc:229
Index getBoxesSize() const
Definition vxc_grid.h:42
Mat_p_Energy IntegrateVXC(const Eigen::MatrixXd &density_matrix) const
void setXCfunctional(const std::string &functional)
SpinResult IntegrateVXCSpin(const Eigen::MatrixXd &dmat_alpha, const Eigen::MatrixXd &dmat_beta) const
#define XTP_LOG(level, log)
Definition logger.h:40
const double hrt2ev
Definition constants.h:53
const double bohr2ang
Definition constants.h:49
Index getMaxThreads()
Definition eigen.h:128
double LegacyAdaptiveShellWidth(const Opt &opt)
void NormalizeGridSearchOptions(Opt &opt)
double LegacyFullWindowHalfWidth(const Opt &opt)
Charge transport classes.
Definition ERIs.h:28
bool XTP_HAS_MKL_OVERLOAD()
Definition eigen.h:41
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
std::string davidson_update
Definition bse_uks.h:47
std::string davidson_tolerance
Definition bse_uks.h:46
std::string davidson_correction
Definition bse_uks.h:45
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