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