votca 2024.2-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/vxc_grid.h"
40
41using boost::format;
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 Index homo = orbitals_.getHomo(); // indexed from 0
71 Index num_of_levels = orbitals_.getBasisSetSize();
72 Index num_of_occlevels = orbitals_.getNumberOfAlphaElectrons();
73
74 std::string ranges = options.get(".ranges").as<std::string>();
75
76 // now check validity, and get rpa, qp, and bse level ranges accordingly
77
78 if (ranges == "factor") {
79
80 double rpamaxfactor = options.get(".rpamax").as<double>();
81 rpamax = Index(rpamaxfactor * double(num_of_levels)) -
82 1; // total number of levels
83
84 double qpminfactor = options.get(".qpmin").as<double>();
85 qpmin =
86 num_of_occlevels - Index(qpminfactor * double(num_of_occlevels)) - 1;
87
88 double qpmaxfactor = options.get(".qpmax").as<double>();
89 qpmax =
90 num_of_occlevels + Index(qpmaxfactor * double(num_of_occlevels)) - 1;
91
92 double bseminfactor = options.get(".bsemin").as<double>();
93 bse_vmin =
94 num_of_occlevels - Index(bseminfactor * double(num_of_occlevels)) - 1;
95
96 double bsemaxfactor = options.get(".bsemax").as<double>();
97 bse_cmax =
98 num_of_occlevels + Index(bsemaxfactor * double(num_of_occlevels)) - 1;
99
100 } else if (ranges == "explicit") {
101 // get explicit numbers
102 rpamax = options.get(".rpamax").as<Index>();
103 qpmin = options.get(".qpmin").as<Index>();
104 qpmax = options.get(".qpmax").as<Index>();
105 bse_vmin = options.get(".bsemin").as<Index>();
106 bse_cmax = options.get(".bsemax").as<Index>();
107 } else if (ranges == "default") {
108 rpamax = num_of_levels - 1;
109 qpmin = 0;
110 qpmax = 3 * homo + 1;
111 bse_vmin = 0;
112 bse_cmax = 3 * homo + 1;
113 } else if (ranges == "full") {
114 rpamax = num_of_levels - 1;
115 qpmin = 0;
116 qpmax = num_of_levels - 1;
117 bse_vmin = 0;
118 bse_cmax = num_of_levels - 1;
119 }
120 std::string ignore_corelevels =
121 options.get(".ignore_corelevels").as<std::string>();
122
123 if (ignore_corelevels == "RPA" || ignore_corelevels == "GW" ||
124 ignore_corelevels == "BSE") {
125 Index ignored_corelevels = CountCoreLevels();
126 if (ignore_corelevels == "RPA") {
127 rpamin = ignored_corelevels;
128 }
129 if (ignore_corelevels == "GW" || ignore_corelevels == "RPA") {
130 if (qpmin < ignored_corelevels) {
131 qpmin = ignored_corelevels;
132 }
133 }
134 if (ignore_corelevels == "GW" || ignore_corelevels == "RPA" ||
135 ignore_corelevels == "BSE") {
136 if (bse_vmin < ignored_corelevels) {
137 bse_vmin = ignored_corelevels;
138 }
139 }
140
142 << TimeStamp() << " Ignoring " << ignored_corelevels
143 << " core levels for " << ignore_corelevels << " and beyond." << flush;
144 }
145
146 // check maximum and minimum sizes
147 if (rpamax > num_of_levels) {
148 rpamax = num_of_levels - 1;
149 }
150 if (qpmax > num_of_levels) {
151 qpmax = num_of_levels - 1;
152 }
153 if (bse_cmax > num_of_levels) {
154 bse_cmax = num_of_levels - 1;
155 }
156 if (bse_vmin < 0) {
157 bse_vmin = 0;
158 }
159 if (qpmin < 0) {
160 qpmin = 0;
161 }
162
163 gwopt_.homo = homo;
164 gwopt_.qpmin = qpmin;
165 gwopt_.qpmax = qpmax;
166 gwopt_.rpamin = rpamin;
167 gwopt_.rpamax = rpamax;
168
169 bseopt_.vmin = bse_vmin;
170 bseopt_.cmax = bse_cmax;
171 bseopt_.homo = homo;
172 bseopt_.qpmin = qpmin;
173 bseopt_.qpmax = qpmax;
174 bseopt_.rpamin = rpamin;
175 bseopt_.rpamax = rpamax;
176
177 orbitals_.setRPAindices(rpamin, rpamax);
178 orbitals_.setGWindices(qpmin, qpmax);
179 orbitals_.setBSEindices(bse_vmin, bse_cmax);
181
182 Index bse_vmax = homo;
183 Index bse_cmin = homo + 1;
184 Index bse_vtotal = bse_vmax - bse_vmin + 1;
185 Index bse_ctotal = bse_cmax - bse_cmin + 1;
186 Index bse_size = bse_vtotal * bse_ctotal;
187
188 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " RPA level range [" << rpamin
189 << ":" << rpamax << "]" << flush;
190 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " GW level range [" << qpmin
191 << ":" << qpmax << "]" << flush;
193 << TimeStamp() << " BSE level range occ[" << bse_vmin << ":" << bse_vmax
194 << "] virt[" << bse_cmin << ":" << bse_cmax << "]" << flush;
195
196 gwopt_.reset_3c = options.get(".gw.rebuild_3c_freq").as<Index>();
197
198 bseopt_.nmax = options.get(".bse.exctotal").as<Index>();
199 if (bseopt_.nmax > bse_size || bseopt_.nmax < 0) {
200 bseopt_.nmax = bse_size;
201 }
202
204 options.get("bse.davidson.correction").as<std::string>();
205
207 options.get("bse.davidson.tolerance").as<std::string>();
208
210 options.get("bse.davidson.update").as<std::string>();
211
212 bseopt_.davidson_maxiter = options.get("bse.davidson.maxiter").as<Index>();
213
214 bseopt_.useTDA = options.get("bse.useTDA").as<bool>();
216 if (!bseopt_.useTDA) {
217 XTP_LOG(Log::error, *pLog_) << " BSE type: full" << flush;
218 } else {
219 XTP_LOG(Log::error, *pLog_) << " BSE type: TDA" << flush;
220 }
221
222 Index full_bse_size = (bseopt_.useTDA) ? bse_size : 2 * bse_size;
223 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " BSE Hamiltonian has size "
224 << full_bse_size << "x" << full_bse_size << flush;
225
226 bseopt_.use_Hqp_offdiag = options.get("bse.use_Hqp_offdiag").as<bool>();
227
230 << " BSE without Hqp offdiagonal elements" << flush;
231 } else {
233 << " BSE with Hqp offdiagonal elements" << flush;
234 }
235
236 bseopt_.max_dyn_iter = options.get("bse.dyn_screen_max_iter").as<Index>();
237 bseopt_.dyn_tolerance = options.get("bse.dyn_screen_tol").as<double>();
238 if (bseopt_.max_dyn_iter > 0) {
240 }
241
244
248 } else if (options.exists(".auxbasisset")) {
249 auxbasis_name_ = options.get(".auxbasisset").as<std::string>();
250 } else {
252 try {
253 BasisSet b;
255 } catch (std::runtime_error&) {
256 std::runtime_error(
257 "There is no auxbasis from the dftcalculation nor did you specify an "
258 "auxbasisset for the gwbse calculation. Also no auxiliary basisset "
259 "for basisset " +
260 dftbasis_name_ + " could be found!");
261 }
263 << " Could not find an auxbasisset using " << auxbasis_name_ << flush;
264 }
265
266 std::string mode = options.get("gw.mode").as<std::string>();
267 if (mode == "G0W0") {
269 } else if (mode == "evGW") {
271 gwopt_.eta = 0.1;
272 }
273
274 XTP_LOG(Log::error, *pLog_) << " Running GW as: " << mode << flush;
276
277 gwopt_.shift = options.get("gw.scissor_shift").as<double>();
279 options.get("gw.qp_sc_limit").as<double>(); // convergence criteria
280 // for qp iteration
281 // [Hartree]]
283 options.get("gw.qp_sc_max_iter").as<Index>(); // convergence
284 // criteria for qp
285 // iteration
286 // [Hartree]]
287
288 if (mode == "evGW") {
289 gwopt_.gw_sc_max_iterations = options.get("gw.sc_max_iter").as<Index>();
290 }
291
293 options.get("gw.sc_limit").as<double>(); // convergence criteria
294 // for shift it
296 << " qp_sc_limit [Hartree]: " << gwopt_.g_sc_limit << flush;
299 << " gw_sc_limit [Hartree]: " << gwopt_.gw_sc_limit << flush;
300 }
301 bseopt_.min_print_weight = options.get("bse.print_weight").as<double>();
302 // print exciton WF composition weight larger than minimum
303
304 // possible tasks
305 std::string tasks_string = options.get(".tasks").as<std::string>();
306 boost::algorithm::to_lower(tasks_string);
307 if (tasks_string.find("all") != std::string::npos) {
308 do_gw_ = true;
309 do_bse_singlets_ = true;
310 do_bse_triplets_ = true;
311 }
312 if (tasks_string.find("gw") != std::string::npos) {
313 do_gw_ = true;
314 }
315 if (tasks_string.find("singlets") != std::string::npos) {
316 do_bse_singlets_ = true;
317 }
318 if (tasks_string.find("triplets") != std::string::npos) {
319 do_bse_triplets_ = true;
320 }
321
322 XTP_LOG(Log::error, *pLog_) << " Tasks: " << flush;
323 if (do_gw_) {
324 XTP_LOG(Log::error, *pLog_) << " GW " << flush;
325 }
326 if (do_bse_singlets_) {
327 XTP_LOG(Log::error, *pLog_) << " singlets " << flush;
328 }
329 if (do_bse_triplets_) {
330 XTP_LOG(Log::error, *pLog_) << " triplets " << flush;
331 }
332 XTP_LOG(Log::error, *pLog_) << " Store: " << flush;
333 if (do_gw_) {
334 XTP_LOG(Log::error, *pLog_) << " GW " << flush;
335 }
336
337 if (options.exists("bse.fragments")) {
338 std::vector<tools::Property*> prop_region =
339 options.Select("bse.fragments.fragment");
340 Index index = 0;
341 for (tools::Property* prop : prop_region) {
342 std::string indices = prop->get("indices").as<std::string>();
343 fragments_.push_back(QMFragment<BSE_Population>(index, indices));
344 index++;
345 }
346 }
347
349 options.get("gw.sigma_integrator").as<std::string>();
351 << " Sigma integration: " << gwopt_.sigma_integration << flush;
352 gwopt_.eta = options.get("gw.eta").as<double>();
353 XTP_LOG(Log::error, *pLog_) << " eta: " << gwopt_.eta << flush;
354 if (gwopt_.sigma_integration == "exact") {
356 << " RPA Hamiltonian size: " << (homo + 1 - rpamin) * (rpamax - homo)
357 << flush;
358 }
359 if (gwopt_.sigma_integration == "cda") {
360 gwopt_.order = options.get("gw.quadrature_order").as<Index>();
362 << " Quadrature integration order : " << gwopt_.order << flush;
364 options.get("gw.quadrature_scheme").as<std::string>();
366 << " Quadrature integration scheme : " << gwopt_.quadrature_scheme
367 << flush;
368 gwopt_.alpha = options.get("gw.alpha").as<double>();
370 << " Alpha smoothing parameter : " << gwopt_.alpha << flush;
371 }
372 gwopt_.qp_solver = options.get("gw.qp_solver").as<std::string>();
373
374 XTP_LOG(Log::error, *pLog_) << " QP solver: " << gwopt_.qp_solver << flush;
375 if (gwopt_.qp_solver == "grid") {
376 gwopt_.qp_grid_steps = options.get("gw.qp_grid_steps").as<Index>();
377 gwopt_.qp_grid_spacing = options.get("gw.qp_grid_spacing").as<double>();
379 << " QP grid steps: " << gwopt_.qp_grid_steps << flush;
381 << " QP grid spacing: " << gwopt_.qp_grid_spacing << flush;
382 }
384 options.get("gw.mixing_order").as<Index>(); // max history in
385 // mixing (0: plain,
386 // 1: linear, >1
387 // Anderson)
388
389 gwopt_.gw_mixing_alpha = options.get("gw.mixing_alpha").as<double>();
390
391 if (mode == "evGW") {
392 if (gwopt_.gw_mixing_order == 0) {
393 XTP_LOG(Log::error, *pLog_) << " evGW with plain update " << std::flush;
394 } else if (gwopt_.gw_mixing_order == 1) {
395 XTP_LOG(Log::error, *pLog_) << " evGW with linear update using alpha "
396 << gwopt_.gw_mixing_alpha << std::flush;
397 } else {
398 XTP_LOG(Log::error, *pLog_) << " evGW with Anderson update with history "
399 << gwopt_.gw_mixing_order << " using alpha "
400 << gwopt_.gw_mixing_alpha << std::flush;
401 }
402 }
403
404 if (options.exists("gw.sigma_plot")) {
405 sigma_plot_states_ = options.get("gw.sigma_plot.states").as<std::string>();
406 sigma_plot_steps_ = options.get("gw.sigma_plot.steps").as<Index>();
407 sigma_plot_spacing_ = options.get("gw.sigma_plot.spacing").as<double>();
409 options.get("gw.sigma_plot.filename").as<std::string>();
410
412 << " Sigma plot states: " << sigma_plot_states_ << flush;
414 << " Sigma plot steps: " << sigma_plot_steps_ << flush;
416 << " Sigma plot spacing: " << sigma_plot_spacing_ << flush;
418 << " Sigma plot filename: " << sigma_plot_filename_ << flush;
419 }
420}
421
423
424 const double hrt2ev = tools::conv::hrt2ev;
425 tools::Property& gwbse_summary = summary.add("GWBSE", "");
426 if (do_gw_) {
427 gwbse_summary.setAttribute("units", "eV");
428 gwbse_summary.setAttribute(
429 "DFTEnergy",
430 (format("%1$+1.6f ") % (orbitals_.getDFTTotalEnergy() * hrt2ev)).str());
431
432 tools::Property& dft_summary = gwbse_summary.add("dft", "");
433 dft_summary.setAttribute("HOMO", gwopt_.homo);
434 dft_summary.setAttribute("LUMO", gwopt_.homo + 1);
435
436 for (Index state = 0; state < gwopt_.qpmax + 1 - gwopt_.qpmin; state++) {
437 tools::Property& level_summary = dft_summary.add("level", "");
438 level_summary.setAttribute("number", state + gwopt_.qpmin);
439 level_summary.add(
440 "dft_energy",
441 (format("%1$+1.6f ") %
442 (orbitals_.MOs().eigenvalues()(state + gwopt_.qpmin) * hrt2ev))
443 .str());
444 level_summary.add(
445 "gw_energy",
446 (format("%1$+1.6f ") % (orbitals_.QPpertEnergies()(state) * hrt2ev))
447 .str());
448
449 level_summary.add("qp_energy",
450 (format("%1$+1.6f ") %
451 (orbitals_.QPdiag().eigenvalues()(state) * hrt2ev))
452 .str());
453 }
454 }
455 if (do_bse_singlets_) {
456 tools::Property& singlet_summary = gwbse_summary.add("singlets", "");
457 for (Index state = 0; state < bseopt_.nmax; ++state) {
458 tools::Property& level_summary = singlet_summary.add("level", "");
459 level_summary.setAttribute("number", state + 1);
460 level_summary.add(
461 "omega", (format("%1$+1.6f ") %
462 (orbitals_.BSESinglets().eigenvalues()(state) * hrt2ev))
463 .str());
465
466 const Eigen::Vector3d& dipoles = (orbitals_.TransitionDipoles())[state];
467 double f = 2 * dipoles.squaredNorm() *
468 orbitals_.BSESinglets().eigenvalues()(state) / 3.0;
469
470 level_summary.add("f", (format("%1$+1.6f ") % f).str());
471 tools::Property& dipol_summary = level_summary.add(
472 "Trdipole", (format("%1$+1.4f %2$+1.4f %3$+1.4f") % dipoles.x() %
473 dipoles.y() % dipoles.z())
474 .str());
475 dipol_summary.setAttribute("unit", "e*bohr");
476 dipol_summary.setAttribute("gauge", "length");
477 }
478 }
479 }
480 if (do_bse_triplets_) {
481 tools::Property& triplet_summary = gwbse_summary.add("triplets", "");
482 for (Index state = 0; state < bseopt_.nmax; ++state) {
483 tools::Property& level_summary = triplet_summary.add("level", "");
484 level_summary.setAttribute("number", state + 1);
485 level_summary.add(
486 "omega", (format("%1$+1.6f ") %
487 (orbitals_.BSETriplets().eigenvalues()(state) * hrt2ev))
488 .str());
489 }
490 }
491 return;
492}
493
494/*
495 * Many-body Green's fuctions theory implementation
496 *
497 * data required from orbitals file
498 * - atomic coordinates
499 * - DFT molecular orbitals (energies and coeffcients)
500 * - DFT exchange-correlation potential matrix in atomic orbitals
501 * - number of electrons, number of levels
502 */
503
504Eigen::MatrixXd GWBSE::CalculateVXC(const AOBasis& dftbasis) {
505 if (orbitals_.getXCFunctionalName().empty()) {
507 } else {
509 throw std::runtime_error("Functionals from DFT " +
510 orbitals_.getXCFunctionalName() + " GWBSE " +
511 functional_ + " differ!");
512 }
513 }
514
515 Vxc_Grid grid;
516 grid.GridSetup(grid_, orbitals_.QMAtoms(), dftbasis);
518 << TimeStamp() << " Setup grid for integration with gridsize: " << grid_
519 << " with " << grid.getGridSize() << " points, divided into "
520 << grid.getBoxesSize() << " boxes" << flush;
521 Vxc_Potential<Vxc_Grid> vxcpotential(grid);
522 vxcpotential.setXCfunctional(functional_);
524 << TimeStamp() << " Integrating Vxc with functional " << functional_
525 << flush;
526 Eigen::MatrixXd DMAT = orbitals_.DensityMatrixGroundState();
527 Mat_p_Energy e_vxc_ao = vxcpotential.IntegrateVXC(DMAT);
528 XTP_LOG(Log::info, *pLog_) << TimeStamp() << " Calculated Vxc" << flush;
530 << TimeStamp() << " Set hybrid exchange factor: " << orbitals_.getScaHFX()
531 << flush;
532 Index qptotal = gwopt_.qpmax - gwopt_.qpmin + 1;
533 Eigen::MatrixXd mos =
534 orbitals_.MOs().eigenvectors().middleCols(gwopt_.qpmin, qptotal);
535
536 Eigen::MatrixXd vxc = mos.transpose() * e_vxc_ao.matrix() * mos;
538 << TimeStamp() << " Calculated exchange-correlation expectation values "
539 << flush;
540
541 return vxc;
542}
543
545
546 // set the parallelization
547 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " Using "
548 << OPENMP::getMaxThreads() << " threads" << flush;
549
550 if (XTP_HAS_MKL_OVERLOAD()) {
552 << TimeStamp() << " Using MKL overload for Eigen " << flush;
553 } else {
555 << TimeStamp()
556 << " Using native Eigen implementation, no BLAS overload " << flush;
557 }
558 Index nogpus = OpenMP_CUDA::UsingGPUs();
559 if (nogpus > 0) {
561 << TimeStamp() << " Using CUDA support for tensor multiplication with "
562 << nogpus << " GPUs." << flush;
563 }
564
566 << TimeStamp() << " Molecule Coordinates [A] " << flush;
567 for (QMAtom& atom : orbitals_.QMAtoms()) {
568 std::string output = (boost::format("%5d"
569 "%5s"
570 " %1.4f %1.4f %1.4f") %
571 atom.getId() % atom.getElement() %
572 (atom.getPos().x() * tools::conv::bohr2ang) %
573 (atom.getPos().y() * tools::conv::bohr2ang) %
574 (atom.getPos().z() * tools::conv::bohr2ang))
575 .str();
576
577 XTP_LOG(Log::error, *pLog_) << output << flush;
578 }
579
580 std::string dft_package = orbitals_.getQMpackage();
582 << TimeStamp() << " DFT data was created by " << dft_package << flush;
583
584 BasisSet dftbs;
585 dftbs.Load(dftbasis_name_);
586
588 << TimeStamp() << " Loaded DFT Basis Set " << dftbasis_name_ << flush;
589
590 AOBasis dftbasis = orbitals_.getDftBasis();
591 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " Filled DFT Basis of size "
592 << dftbasis.AOBasisSize() << flush;
594 << TimeStamp() << " Loaded Auxbasis Set " << auxbasis_name_ << flush;
595
597 AOBasis auxbasis = orbitals_.getAuxBasis();
598 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " Filled Auxbasis of size "
599 << auxbasis.AOBasisSize() << flush;
600
601 if ((do_bse_singlets_ || do_bse_triplets_) && fragments_.size() > 0) {
602 for (const auto& frag : fragments_) {
603 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " Fragment " << frag.getId()
604 << " size:" << frag.size() << flush;
605 }
606 }
607
608 if (!do_gw_ && !orbitals_.hasQPdiag()) {
609 throw std::runtime_error(
610 "You want no GW calculation but the orb file has no QPcoefficients for "
611 "BSE");
612 }
613 TCMatrix_gwbse Mmn;
614 // rpamin here, because RPA needs till rpamin
615 Index max_3c = std::max(bseopt_.cmax, gwopt_.qpmax);
616 Mmn.Initialize(auxbasis.AOBasisSize(), gwopt_.rpamin, max_3c, gwopt_.rpamin,
617 gwopt_.rpamax);
619 << TimeStamp()
620 << " Calculating Mmn_beta (3-center-repulsion x orbitals) " << flush;
621 Mmn.Fill(auxbasis, dftbasis, orbitals_.MOs().eigenvectors());
623 << TimeStamp() << " Removed " << Mmn.Removedfunctions()
624 << " functions from Aux Coulomb matrix to avoid near linear dependencies"
625 << flush;
627 << TimeStamp() << " Calculated Mmn_beta (3-center-repulsion x orbitals) "
628 << flush;
629
630 Eigen::MatrixXd Hqp;
631 if (do_gw_) {
632
633 std::chrono::time_point<std::chrono::system_clock> start =
634 std::chrono::system_clock::now();
635 Eigen::MatrixXd vxc = CalculateVXC(dftbasis);
636 GW gw = GW(*pLog_, Mmn, vxc, orbitals_.MOs().eigenvalues());
637 gw.configure(gwopt_);
638
640
641 if (!sigma_plot_states_.empty()) {
644 }
645
646 // store perturbative QP energy data in orbitals object (DFT, S_x,S_c, V_xc,
647 // E_qp)
650
652 << TimeStamp() << " Calculating offdiagonal part of Sigma " << flush;
653 gw.CalculateHQP();
655 << TimeStamp() << " Calculated offdiagonal part of Sigma " << flush;
656
657 Hqp = gw.getHQP();
658
659 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es =
661 if (es.info() == Eigen::ComputationInfo::Success) {
663 << TimeStamp() << " Diagonalized QP Hamiltonian " << flush;
664 }
665
666 orbitals_.QPdiag().eigenvectors() = es.eigenvectors();
667 orbitals_.QPdiag().eigenvalues() = es.eigenvalues();
668 std::chrono::duration<double> elapsed_time =
669 std::chrono::system_clock::now() - start;
670 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " GW calculation took "
671 << elapsed_time.count() << " seconds." << flush;
672
673 } else {
674 if (orbitals_.getGWAmax() != gwopt_.qpmax ||
678 throw std::runtime_error(
679 "The ranges for GW and RPA do not agree with the ranges from the "
680 ".orb file, rerun your GW calculation");
681 }
682 const Eigen::MatrixXd& qpcoeff = orbitals_.QPdiag().eigenvectors();
683
684 Hqp = qpcoeff * orbitals_.QPdiag().eigenvalues().asDiagonal() *
685 qpcoeff.transpose();
686 }
687
688 // proceed only if BSE requested
690
691 std::chrono::time_point<std::chrono::system_clock> start =
692 std::chrono::system_clock::now();
693
694 BSE bse = BSE(*pLog_, Mmn);
696
697 // store the direct contribution to the static BSE results
698 Eigen::VectorXd Hd_static_contrib_triplet;
699 Eigen::VectorXd Hd_static_contrib_singlet;
700
701 if (do_bse_triplets_) {
704 << TimeStamp() << " Solved BSE for triplets " << flush;
706 }
707
708 if (do_bse_singlets_) {
711 << TimeStamp() << " Solved BSE for singlets " << flush;
713 }
714
715 // do perturbative dynamical screening in BSE
717
718 if (do_bse_triplets_) {
720 orbitals_);
721 }
722
723 if (do_bse_singlets_) {
725 orbitals_);
726 }
727 }
728
729 std::chrono::duration<double> elapsed_time =
730 std::chrono::system_clock::now() - start;
731 XTP_LOG(Log::error, *pLog_) << TimeStamp() << " BSE calculation took "
732 << elapsed_time.count() << " seconds." << flush;
733 }
735 << TimeStamp() << " GWBSE calculation finished " << flush;
736 return true;
737}
738
739} // namespace xtp
740} // namespace votca
const Eigen::VectorXd & eigenvalues() const
Definition eigensystem.h:30
const Eigen::MatrixXd & eigenvectors() const
Definition eigensystem.h:33
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 Solve_singlets(Orbitals &orb) const
Definition bse.cc:138
void Analyze_singlets(std::vector< QMFragment< BSE_Population > > fragments, const Orbitals &orb) const
Definition bse.cc:304
void Perturbative_DynamicalScreening(const QMStateType &type, Orbitals &orb)
Definition bse.cc:521
void Analyze_triplets(std::vector< QMFragment< BSE_Population > > fragments, const Orbitals &orb) const
Definition bse.cc:360
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:148
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:544
std::vector< QMFragment< BSE_Population > > fragments_
Definition gwbse.h:97
std::string sigma_plot_states_
Definition gwbse.h:88
BSE::options bseopt_
Definition gwbse.h:86
void addoutput(tools::Property &summary)
Definition gwbse.cc:422
bool do_bse_triplets_
Definition gwbse.h:78
std::string grid_
Definition gwbse.h:83
std::string auxbasis_name_
Definition gwbse.h:94
bool do_bse_singlets_
Definition gwbse.h:77
Orbitals & orbitals_
Definition gwbse.h:73
Eigen::MatrixXd CalculateVXC(const AOBasis &dftbasis)
Definition gwbse.cc:504
Index sigma_plot_steps_
Definition gwbse.h:89
std::string functional_
Definition gwbse.h:82
std::string dftbasis_name_
Definition gwbse.h:95
void Initialize(tools::Property &options)
Definition gwbse.cc:60
Index CountCoreLevels()
Definition gwbse.cc:46
GW::options gwopt_
Definition gwbse.h:85
double sigma_plot_spacing_
Definition gwbse.h:90
Logger * pLog_
Definition gwbse.h:72
std::string sigma_plot_filename_
Definition gwbse.h:91
bool do_dynamical_screening_bse_
Definition gwbse.h:79
void PlotSigma(std::string filename, Index steps, double spacing, std::string states) const
Definition gw.cc:432
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > DiagonalizeQPHamiltonian() const
Definition gw.cc:68
void configure(const options &opt)
Definition gw.cc:35
Eigen::MatrixXd getHQP() const
Definition gw.cc:62
void CalculateHQP()
Definition gw.cc:426
Eigen::VectorXd getGWAResults() const
Definition gw.cc:237
Eigen::VectorXd RPAInputEnergies() const
Definition gw.h:89
void CalculateGWPerturbation()
Definition gw.cc:143
Eigen::MatrixXd & matrix()
Definition eigen.h:80
static Index UsingGPUs()
const tools::EigenSystem & BSETriplets() const
Definition orbitals.h:324
bool hasQPdiag() const
Definition orbitals.h:314
const tools::EigenSystem & QPdiag() const
Definition orbitals.h:317
Index getHomo() const
Definition orbitals.h:68
double getDFTTotalEnergy() const
Definition orbitals.h:193
bool hasAuxbasisName() const
Definition orbitals.h:234
const tools::EigenSystem & BSESinglets() const
Definition orbitals.h:334
const Eigen::VectorXd & QPpertEnergies() const
Definition orbitals.h:308
const AOBasis & getAuxBasis() const
Definition orbitals.h:218
Index getNumberOfAlphaElectrons() const
Definition orbitals.h:93
const std::string & getQMpackage() const
Definition orbitals.h:112
void setRPAindices(Index rpamin, Index rpamax)
Definition orbitals.h:257
bool hasTransitionDipoles() const
Definition orbitals.h:365
void SetupAuxBasis(std::string aux_basis_name)
Definition orbitals.cc:99
Index getRPAmax() const
Definition orbitals.h:264
void setBSEindices(Index vmin, Index cmax)
Definition orbitals.h:273
void setGWindices(Index qpmin, Index qpmax)
Definition orbitals.h:244
const std::vector< Eigen::Vector3d > & TransitionDipoles() const
Definition orbitals.h:369
const std::string getAuxbasisName() const
Definition orbitals.h:238
double getScaHFX() const
Definition orbitals.h:292
Index getBasisSetSize() const
Definition orbitals.h:64
const Eigen::VectorXd & RPAInputEnergies() const
Definition orbitals.h:299
const std::string & getXCGrid() const
Definition orbitals.h:188
Index getGWAmin() const
Definition orbitals.h:249
bool hasECPName() const
Definition orbitals.h:102
const tools::EigenSystem & MOs() const
Definition orbitals.h:122
Eigen::MatrixXd DensityMatrixGroundState() const
Definition orbitals.cc:183
const QMMolecule & QMAtoms() const
Definition orbitals.h:172
const std::string & getXCFunctionalName() const
Definition orbitals.h:185
Index getGWAmax() const
Definition orbitals.h:251
void SetFlagUseHqpOffdiag(bool flag)
Definition orbitals.h:409
const std::string & getDFTbasisName() const
Definition orbitals.h:202
void setTDAApprox(bool usedTDA)
Definition orbitals.h:268
Index getRPAmin() const
Definition orbitals.h:262
const AOBasis & getDftBasis() const
Definition orbitals.h:208
void setXCFunctionalName(std::string functionalname)
Definition orbitals.h:182
container for QM atoms
Definition qmatom.h:37
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)
#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
bool XTP_HAS_MKL_OVERLOAD()
Definition eigen.h:41
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
std::string davidson_tolerance
Definition bse.h:58
std::string davidson_update
Definition bse.h:59
double min_print_weight
Definition bse.h:61
std::string davidson_correction
Definition bse.h:57
double gw_sc_limit
Definition gw.h:53
Index gw_mixing_order
Definition gw.h:64
Index gw_sc_max_iterations
Definition gw.h:54
Index qp_grid_steps
Definition gw.h:62
std::string quadrature_scheme
Definition gw.h:66
std::string qp_solver
Definition gw.h:60
Index g_sc_max_iterations
Definition gw.h:52
double qp_grid_spacing
Definition gw.h:63
double g_sc_limit
Definition gw.h:51
double gw_mixing_alpha
Definition gw.h:65
std::string sigma_integration
Definition gw.h:57