votca 2026-dev
Loading...
Searching...
No Matches
bse_uks.cc
Go to the documentation of this file.
1/*
2 * Copyright 2009-2026 The VOTCA Development Team
3 * (http://www.votca.org)
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License")
6 *
7 * You may not use this file except in compliance with the License.
8 * You may obtain a copy of the License at
9 *
10 * http://www.apache.org/licenses/LICENSE-2.0
11 *
12 * Unless required by applicable law or agreed to in writing, software
13 * distributed under the License is distributed on an "AS IS" BASIS,
14 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
15 * See the License for the specific language governing permissions and
16 * limitations under the License.
17 *
18 */
19
20#include <chrono>
21#include <iostream>
22
23#include <votca/tools/linalg.h>
24
27#include "votca/xtp/bse_uks.h"
30#include "votca/xtp/rpa_uks.h"
31
32using std::flush;
33
34namespace {
35
36template <class OP>
37Eigen::VectorXd ExpValue(const Eigen::MatrixXd& state1, OP OPxstate2) {
38 return state1.cwiseProduct(OPxstate2.eval()).colwise().sum().transpose();
39}
40
41Eigen::VectorXd ExpValue(const Eigen::MatrixXd& state1,
42 const Eigen::MatrixXd& OPxstate2) {
43 return state1.cwiseProduct(OPxstate2).colwise().sum().transpose();
44}
45
46} // namespace
47
48namespace votca {
49namespace xtp {
50
52 const options& opt, Index homo_alpha, Index homo_beta,
53 const Eigen::VectorXd& RPAInputEnergiesAlpha,
54 const Eigen::VectorXd& RPAInputEnergiesBeta,
55 const Eigen::MatrixXd& Hqp_alpha_in, const Eigen::MatrixXd& Hqp_beta_in,
56 const Eigen::VectorXd& epsilon_0_inv,
57 const Eigen::MatrixXd& epsilon_eigenvectors) {
58
59 opt_ = opt;
60 homo_alpha_ = homo_alpha;
61 homo_beta_ = homo_beta;
62
63 alpha_vtotal_ = homo_alpha_ - opt_.vmin + 1;
64 alpha_ctotal_ = opt_.cmax - (homo_alpha_ + 1) + 1;
66
67 beta_vtotal_ = homo_beta_ - opt_.vmin + 1;
68 beta_ctotal_ = opt_.cmax - (homo_beta_ + 1) + 1;
70
71 if (opt_.use_Hqp_offdiag) {
73 AdjustHqpSize(Hqp_alpha_in, RPAInputEnergiesAlpha, homo_alpha_);
74 Hqp_beta_ = AdjustHqpSize(Hqp_beta_in, RPAInputEnergiesBeta, homo_beta_);
75 } else {
76 Hqp_alpha_ = AdjustHqpSize(Hqp_alpha_in, RPAInputEnergiesAlpha, homo_alpha_)
77 .diagonal()
78 .asDiagonal();
79 Hqp_beta_ = AdjustHqpSize(Hqp_beta_in, RPAInputEnergiesBeta, homo_beta_)
80 .diagonal()
81 .asDiagonal();
82 }
83
84 // Store the statically screened interaction in the working copy Mmn_.
85 Mmn_ = Mmn_raw_;
86 Mmn_.alpha.MultiplyRightWithAuxMatrix(epsilon_eigenvectors);
87 Mmn_.beta.MultiplyRightWithAuxMatrix(epsilon_eigenvectors);
88
89 epsilon_0_inv_ = epsilon_0_inv;
90}
91
92Eigen::MatrixXd BSE_UKS::AdjustHqpSize(const Eigen::MatrixXd& Hqp,
93 const Eigen::VectorXd& RPAInputEnergies,
94 Index homo) const {
95 Index bse_vtotal = homo - opt_.vmin + 1;
96 Index bse_ctotal = opt_.cmax - (homo + 1) + 1;
97 Index hqp_size = bse_vtotal + bse_ctotal;
98 Index gwsize = opt_.qpmax - opt_.qpmin + 1;
99 Index RPAoffset = opt_.vmin - opt_.rpamin;
100 Eigen::MatrixXd Hqp_BSE = Eigen::MatrixXd::Zero(hqp_size, hqp_size);
101
102 if (opt_.vmin >= opt_.qpmin) {
103 Index start = opt_.vmin - opt_.qpmin;
104 if (opt_.cmax <= opt_.qpmax) {
105 Hqp_BSE = Hqp.block(start, start, hqp_size, hqp_size);
106 } else {
107 Index virtoffset = gwsize - start;
108 Hqp_BSE.topLeftCorner(virtoffset, virtoffset) =
109 Hqp.block(start, start, virtoffset, virtoffset);
110
111 Index virt_extra = opt_.cmax - opt_.qpmax;
112 Hqp_BSE.diagonal().tail(virt_extra) =
113 RPAInputEnergies.segment(RPAoffset + virtoffset, virt_extra);
114 }
115 }
116
117 if (opt_.vmin < opt_.qpmin) {
118 Index occ_extra = opt_.qpmin - opt_.vmin;
119 Hqp_BSE.diagonal().head(occ_extra) =
120 RPAInputEnergies.segment(RPAoffset, occ_extra);
121
122 Hqp_BSE.block(occ_extra, occ_extra, gwsize, gwsize) = Hqp;
123
124 if (opt_.cmax > opt_.qpmax) {
125 Index virtoffset = occ_extra + gwsize;
126 Index virt_extra = opt_.cmax - opt_.qpmax;
127 Hqp_BSE.diagonal().tail(virt_extra) =
128 RPAInputEnergies.segment(RPAoffset + virtoffset, virt_extra);
129 }
130 }
131
132 return Hqp_BSE;
133}
134
136 const Eigen::VectorXd& RPAInputEnergiesAlpha,
137 const Eigen::VectorXd& RPAInputEnergiesBeta, double energy) {
138
139 RPA_UKS rpa(log_, Mmn_raw_);
140 rpa.configure(homo_alpha_, homo_beta_, opt_.rpamin, opt_.rpamax);
141 rpa.setRPAInputEnergies(RPAInputEnergiesAlpha, RPAInputEnergiesBeta);
142
143 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es(
144 rpa.calculate_epsilon_r(energy));
145
146 Mmn_ = Mmn_raw_;
147 Mmn_.alpha.MultiplyRightWithAuxMatrix(es.eigenvectors());
148 Mmn_.beta.MultiplyRightWithAuxMatrix(es.eigenvectors());
149
150 epsilon_0_inv_ = Eigen::VectorXd::Zero(es.eigenvalues().size());
151 for (Index i = 0; i < es.eigenvalues().size(); ++i) {
152 if (es.eigenvalues()(i) > 1e-8) {
153 epsilon_0_inv_(i) = 1.0 / es.eigenvalues()(i);
154 }
155 }
156}
157
158template <typename BSE_OPERATOR>
162 opt.homo_beta = homo_beta_;
163 opt.rpamin = opt_.rpamin;
164 opt.qpmin = opt_.qpmin;
165 opt.vmin = opt_.vmin;
166 opt.cmax = opt_.cmax;
167 H.configure(opt);
168}
169
170template <typename BSE_OPERATOR>
172 const Orbitals& orb, const BSE_OPERATOR& H) const {
173
174 const tools::EigenSystem& BSECoefs = orb.BSEUKS();
175
176 ExpectationValues expectation_values;
177
178 const Eigen::MatrixXd temp = H * BSECoefs.eigenvectors();
179
180 expectation_values.direct_term = ExpValue(BSECoefs.eigenvectors(), temp);
181
182 if (!orb.getTDAApprox()) {
183 expectation_values.direct_term +=
184 ExpValue(BSECoefs.eigenvectors2(), H * BSECoefs.eigenvectors2());
185 expectation_values.cross_term =
186 2.0 * ExpValue(BSECoefs.eigenvectors2(), temp);
187 } else {
188 expectation_values.cross_term = Eigen::VectorXd::Zero(0);
189 }
190
191 return expectation_values;
192}
193
194template <typename BSE_OPERATOR>
196 Index state, const Orbitals& orb, const BSE_OPERATOR& H) const {
197
198 const tools::EigenSystem& BSECoefs = orb.BSEUKS();
199
200 ExpectationValues expectation_values;
201
202 const Eigen::MatrixXd Xstate = BSECoefs.eigenvectors().col(state);
203 const Eigen::MatrixXd temp = H * Xstate;
204
205 expectation_values.direct_term = ExpValue(Xstate, temp);
206
207 if (!orb.getTDAApprox()) {
208 const Eigen::MatrixXd Ystate = BSECoefs.eigenvectors2().col(state);
209 expectation_values.direct_term += ExpValue(Ystate, H * Ystate);
210 expectation_values.cross_term = 2.0 * ExpValue(Ystate, temp);
211 } else {
212 expectation_values.cross_term = Eigen::VectorXd::Zero(0);
213 }
214
215 return expectation_values;
216}
217
222 opt.homo_beta = homo_beta_;
223 opt.rpamin = opt_.rpamin;
224 opt.qpmin = opt_.qpmin;
225 opt.vmin = opt_.vmin;
226 opt.cmax = opt_.cmax;
227 H.configure(opt);
228
230 << TimeStamp() << " Setup combined UKS TDA Hamiltonian " << flush;
231 return solve_hermitian(H);
232}
233
234/* tools::EigenSystem BSE_UKS::Solve_excitons_uks_BTDA() const {
235 ExcitonUKSOperator_TDA A(epsilon_0_inv_, Mmn_, Hqp_alpha_, Hqp_beta_);
236 ExcitonUKSOperator_BTDA_B B(epsilon_0_inv_, Mmn_, Hqp_alpha_, Hqp_beta_);
237
238 BSEOperatorUKS_Options opt;
239 opt.homo_alpha = homo_alpha_;
240 opt.homo_beta = homo_beta_;
241 opt.rpamin = opt_.rpamin;
242 opt.qpmin = opt_.qpmin;
243 opt.vmin = opt_.vmin;
244 opt.cmax = opt_.cmax;
245
246 A.configure(opt);
247 B.configure(opt);
248
249 XTP_LOG(Log::error, log_)
250 << TimeStamp() << " Setup combined UKS full BSE exciton hamiltonian "
251 << flush;
252
253 // Cheap production initial guess:
254 // rank excitations by local diagonal full-BSE estimate
255 // omega_i = sqrt(max(0, A_ii^2 - B_ii^2))
256 // but keep the actual start vectors X-only for robustness with the
257 // current HAM Davidson implementation.
258 const Eigen::VectorXd adiag = A.diagonal();
259 const Eigen::VectorXd bdiag = B.diagonal();
260 Eigen::MatrixXd initial_guess =
261 BuildFullBSEXRankedInitialGuess(adiag, bdiag, opt_.nmax);
262
263 HamiltonianOperator<ExcitonUKSOperator_TDA, ExcitonUKSOperator_BTDA_B> Hop(A,
264 B);
265
266 DavidsonSolver DS(log_);
267 DS.set_correction(opt_.davidson_correction);
268 DS.set_tolerance(opt_.davidson_tolerance);
269 DS.set_size_update(opt_.davidson_update);
270 DS.set_iter_max(opt_.davidson_maxiter);
271 DS.set_max_search_space(10 * opt_.nmax);
272 DS.set_matrix_type("HAM");
273 DS.solve(Hop, opt_.nmax, initial_guess);
274
275 tools::EigenSystem result;
276 result.eigenvalues() = DS.eigenvalues();
277
278 Eigen::MatrixXd tmpX = DS.eigenvectors().topRows(A.rows());
279 Eigen::MatrixXd tmpY = DS.eigenvectors().bottomRows(B.rows());
280
281 Eigen::VectorXd normX = tmpX.colwise().squaredNorm();
282 Eigen::VectorXd normY = tmpY.colwise().squaredNorm();
283 Eigen::ArrayXd sqinvnorm = (normX - normY).array().inverse().cwiseSqrt();
284
285 result.eigenvectors() = tmpX * sqinvnorm.matrix().asDiagonal();
286 result.eigenvectors2() = tmpY * sqinvnorm.matrix().asDiagonal();
287
288 return result;
289}*/
290
294
297 opt.homo_beta = homo_beta_;
298 opt.rpamin = opt_.rpamin;
299 opt.qpmin = opt_.qpmin;
300 opt.vmin = opt_.vmin;
301 opt.cmax = opt_.cmax;
302
303 A.configure(opt);
304 B.configure(opt);
305
307 << TimeStamp() << " Setup combined UKS full BSE exciton hamiltonian "
308 << flush;
309
310 // Dense fallback for small systems:
311 // this avoids architecture-sensitive Davidson/HAM behavior in tiny
312 // regression tests while keeping the iterative solver for production sizes.
313 constexpr Index dense_threshold = 128;
314
315 if (A.rows() <= dense_threshold) {
317 << TimeStamp()
318 << " Using dense full UKS-BSE solve for small system (dim=" << A.rows()
319 << ")" << flush;
320 std::chrono::time_point<std::chrono::system_clock> start =
321 std::chrono::system_clock::now();
322
323 const Eigen::MatrixXd Ad = A.dense_matrix();
324 const Eigen::MatrixXd Bd = B.dense_matrix();
325 const Index n = Ad.rows();
326
327 Eigen::MatrixXd Hfull = Eigen::MatrixXd::Zero(2 * n, 2 * n);
328 Hfull.topLeftCorner(n, n) = Ad;
329 Hfull.topRightCorner(n, n) = Bd;
330 Hfull.bottomLeftCorner(n, n) = -Bd;
331 Hfull.bottomRightCorner(n, n) = -Ad;
332
333 Eigen::EigenSolver<Eigen::MatrixXd> es(Hfull, true);
334 if (es.info() != Eigen::Success) {
335 throw std::runtime_error("Dense full UKS-BSE diagonalization failed.");
336 }
337
338 using RootInfo =
339 std::tuple<double, Index, double>; // (eval, col, imag_abs)
340 std::vector<RootInfo> positive_roots;
341 positive_roots.reserve(2 * n);
342
343 for (Index i = 0; i < 2 * n; ++i) {
344 const std::complex<double> ev = es.eigenvalues()(i);
345 if (std::abs(ev.imag()) < 1e-8 && ev.real() > 0.0) {
346 positive_roots.emplace_back(ev.real(), i, std::abs(ev.imag()));
347 }
348 }
349
350 std::sort(positive_roots.begin(), positive_roots.end(),
351 [](const RootInfo& a, const RootInfo& b) {
352 return std::get<0>(a) < std::get<0>(b);
353 });
354
355 if (positive_roots.empty()) {
356 throw std::runtime_error(
357 "Dense full UKS-BSE diagonalization produced no positive real "
358 "roots.");
359 }
360
361 const Index nroots =
362 std::min<Index>(opt_.nmax, static_cast<Index>(positive_roots.size()));
363
364 tools::EigenSystem result;
365 result.eigenvalues() = Eigen::VectorXd::Zero(nroots);
366 result.eigenvectors() = Eigen::MatrixXd::Zero(n, nroots);
367 result.eigenvectors2() = Eigen::MatrixXd::Zero(n, nroots);
368
369 for (Index iroot = 0; iroot < nroots; ++iroot) {
370 const Index col = std::get<1>(positive_roots[iroot]);
371 Eigen::VectorXd vec = es.eigenvectors().col(col).real();
372
373 Eigen::VectorXd X = vec.head(n);
374 Eigen::VectorXd Y = vec.tail(n);
375
376 // Deterministic phase convention: make the largest |X_k| positive.
377 Eigen::Index imax = 0;
378 X.cwiseAbs().maxCoeff(&imax);
379 if (X(imax) < 0.0) {
380 X *= -1.0;
381 Y *= -1.0;
382 }
383
384 const double norm = X.squaredNorm() - Y.squaredNorm();
385 if (std::abs(norm) < 1e-12) {
386 throw std::runtime_error(
387 "Dense full UKS-BSE eigenvector has near-zero (X^2-Y^2) norm.");
388 }
389
390 const double invnorm = 1.0 / std::sqrt(std::abs(norm));
391
392 result.eigenvalues()(iroot) = std::get<0>(positive_roots[iroot]);
393 result.eigenvectors().col(iroot) = X * invnorm;
394 result.eigenvectors2().col(iroot) = Y * invnorm;
395 }
396
397 std::chrono::time_point<std::chrono::system_clock> end =
398 std::chrono::system_clock::now();
399 std::chrono::duration<double> elapsed_time = end - start;
400
402 << TimeStamp() << " Dense full UKS-BSE diagonalization done in "
403 << elapsed_time.count() << " secs" << flush;
404
405 return result;
406 }
407
408 // Default Davidson/HAM path for larger systems
409 const Eigen::VectorXd adiag = A.diagonal();
410 const Eigen::VectorXd bdiag = B.diagonal();
411 Eigen::MatrixXd initial_guess =
412 BuildFullBSEXRankedInitialGuess(adiag, bdiag, opt_.nmax);
413
415 B);
416
418 DS.set_correction(opt_.davidson_correction);
419 DS.set_tolerance(opt_.davidson_tolerance);
420 DS.set_size_update(opt_.davidson_update);
421 DS.set_iter_max(opt_.davidson_maxiter);
422 DS.set_max_search_space(10 * opt_.nmax);
423 DS.set_matrix_type("HAM");
424 DS.solve(Hop, opt_.nmax, initial_guess);
425
426 tools::EigenSystem result;
427 result.eigenvalues() = DS.eigenvalues();
428
429 Eigen::MatrixXd tmpX = DS.eigenvectors().topRows(A.rows());
430 Eigen::MatrixXd tmpY = DS.eigenvectors().bottomRows(B.rows());
431
432 Eigen::VectorXd normX = tmpX.colwise().squaredNorm();
433 Eigen::VectorXd normY = tmpY.colwise().squaredNorm();
434 Eigen::ArrayXd sqinvnorm = (normX - normY).array().inverse().cwiseSqrt();
435
436 result.eigenvectors() = tmpX * sqinvnorm.matrix().asDiagonal();
437 result.eigenvectors2() = tmpY * sqinvnorm.matrix().asDiagonal();
438
439 return result;
440}
441
443 orb.setTDAApprox(opt_.useTDA);
444 if (opt_.useTDA) {
446 } else {
448 }
450}
451
452template <typename BSE_OPERATOR>
454 std::chrono::time_point<std::chrono::system_clock> start =
455 std::chrono::system_clock::now();
456
457 tools::EigenSystem result;
458
460 DS.set_correction(opt_.davidson_correction);
461 DS.set_tolerance(opt_.davidson_tolerance);
462 DS.set_size_update(opt_.davidson_update);
463 DS.set_iter_max(opt_.davidson_maxiter);
464 DS.set_max_search_space(10 * opt_.nmax);
465 DS.solve(h, opt_.nmax);
466
467 result.eigenvalues() = DS.eigenvalues();
468 result.eigenvectors() = DS.eigenvectors();
469
470 std::chrono::time_point<std::chrono::system_clock> end =
471 std::chrono::system_clock::now();
472 std::chrono::duration<double> elapsed_time = end - start;
473
474 XTP_LOG(Log::info, log_) << TimeStamp() << " Diagonalization done in "
475 << elapsed_time.count() << " secs" << flush;
476
477 return result;
478}
479
480template <typename BSE_OPERATOR_A, typename BSE_OPERATOR_B>
482 BSE_OPERATOR_A& Aop, BSE_OPERATOR_B& Bop) const {
483 std::chrono::time_point<std::chrono::system_clock> start =
484 std::chrono::system_clock::now();
485
487
489 DS.set_correction(opt_.davidson_correction);
490 DS.set_tolerance(opt_.davidson_tolerance);
491 DS.set_size_update(opt_.davidson_update);
492 DS.set_iter_max(opt_.davidson_maxiter);
493 DS.set_max_search_space(10 * opt_.nmax);
494 DS.set_matrix_type("HAM");
495 DS.solve(Hop, opt_.nmax);
496
497 tools::EigenSystem result;
498 result.eigenvalues() = DS.eigenvalues();
499
500 Eigen::MatrixXd tmpX = DS.eigenvectors().topRows(Aop.rows());
501 Eigen::MatrixXd tmpY = DS.eigenvectors().bottomRows(Bop.rows());
502
503 Eigen::VectorXd normX = tmpX.colwise().squaredNorm();
504 Eigen::VectorXd normY = tmpY.colwise().squaredNorm();
505 Eigen::ArrayXd sqinvnorm = (normX - normY).array().inverse().cwiseSqrt();
506
507 result.eigenvectors() = tmpX * sqinvnorm.matrix().asDiagonal();
508 result.eigenvectors2() = tmpY * sqinvnorm.matrix().asDiagonal();
509
510 std::chrono::time_point<std::chrono::system_clock> end =
511 std::chrono::system_clock::now();
512 std::chrono::duration<double> elapsed_time = end - start;
513
514 XTP_LOG(Log::info, log_) << TimeStamp() << " Diagonalization done in "
515 << elapsed_time.count() << " secs" << flush;
516
517 return result;
518}
519
520void BSE_UKS::PrintWeightsUKS(const Eigen::VectorXd& weights) const {
521 struct Contribution {
522 double weight;
523 bool is_alpha;
524 Index v;
525 Index c;
526 };
527
528 constexpr Index max_print = 8;
529
530 double alpha_weight = weights.head(alpha_size_).sum();
531 double beta_weight = weights.tail(beta_size_).sum();
532
534 << boost::format(
535 " alpha-sector: %1$+6.2f%% beta-sector: %2$+6.2f%%") %
536 (100.0 * alpha_weight) % (100.0 * beta_weight)
537 << flush;
538
539 std::vector<Contribution> contributions;
540 contributions.reserve(alpha_size_ + beta_size_);
541
542 for (Index i = 0; i < alpha_size_; ++i) {
543 if (std::abs(weights(i)) > opt_.min_print_weight) {
544 contributions.push_back(
545 {weights(i), true, i / alpha_ctotal_, i % alpha_ctotal_});
546 }
547 }
548
549 for (Index i = 0; i < beta_size_; ++i) {
550 const double w = weights(alpha_size_ + i);
551 if (std::abs(w) > opt_.min_print_weight) {
552 contributions.push_back({w, false, i / beta_ctotal_, i % beta_ctotal_});
553 }
554 }
555
556 std::sort(contributions.begin(), contributions.end(),
557 [](const Contribution& a, const Contribution& b) {
558 return std::abs(a.weight) > std::abs(b.weight);
559 });
560
561 const Index nprint =
562 std::min<Index>(max_print, static_cast<Index>(contributions.size()));
563
564 for (Index i = 0; i < nprint; ++i) {
565 const Contribution& c = contributions[i];
566
567 if (c.is_alpha) {
569 << boost::format(
570 " [alpha] HOMO-%1$-3d -> LUMO+%2$-3d : "
571 "%3$+6.2f%%") %
572 (homo_alpha_ - (opt_.vmin + c.v)) %
573 ((homo_alpha_ + 1 + c.c) - (homo_alpha_ + 1)) %
574 (100.0 * c.weight)
575 << flush;
576 } else {
578 << boost::format(
579 " [beta ] HOMO-%1$-3d -> LUMO+%2$-3d : "
580 "%3$+6.2f%%") %
581 (homo_beta_ - (opt_.vmin + c.v)) %
582 ((homo_beta_ + 1 + c.c) - (homo_beta_ + 1)) %
583 (100.0 * c.weight)
584 << flush;
585 }
586 }
587
588 if (static_cast<Index>(contributions.size()) > nprint) {
590 << boost::format(
591 " ... %1$d more contributions above threshold") %
592 (static_cast<Index>(contributions.size()) - nprint)
593 << flush;
594 }
595}
596
598 std::vector<QMFragment<BSE_Population>> fragments,
599 const Orbitals& orb) const {
600 (void)fragments;
601
602 const tools::EigenSystem& es = orb.BSEUKS();
603 const Eigen::VectorXd oscs =
605
606 if (orb.getTDAApprox()) {
608 << " ====== combined UKS TDA exciton energies (eV) ====== " << flush;
609 } else {
611 << " ====== combined UKS full-BSE exciton energies (eV) ====== "
612 << flush;
613 }
614
615 for (Index i = 0; i < std::min<Index>(opt_.nmax, es.eigenvalues().size());
616 ++i) {
617 XTP_LOG(Log::error, log_) << boost::format(" XU%-4d %+1.6f") % (i + 1) %
619 << flush;
620
621 const Eigen::Vector3d& trdip = orb.TransitionDipoles()[i];
622 const double osc = (i < oscs.size()) ? oscs(i) : 0.0;
623
625 << boost::format(
626 " TrDipole length gauge[e*bohr] dx = %1$+1.4f dy = "
627 "%2$+1.4f dz = %3$+1.4f |d|^2 = %4$+1.4f f = %5$+1.4f") %
628 trdip[0] % trdip[1] % trdip[2] % (trdip.squaredNorm()) % osc
629 << flush;
630
631 Eigen::VectorXd weights = es.eigenvectors().col(i).cwiseAbs2();
632 if (!orb.getTDAApprox()) {
633 weights -= es.eigenvectors2().col(i).cwiseAbs2();
634 }
635
636 PrintWeightsUKS(weights);
637 }
638}
639
641
642 const tools::EigenSystem& BSECoefs = orb.BSEUKS();
643 const Eigen::VectorXd& BSEenergies = BSECoefs.eigenvalues();
644
645 const Eigen::VectorXd& RPAInputEnergiesAlpha = orb.RPAInputEnergiesAlpha();
646 const Eigen::VectorXd& RPAInputEnergiesBeta = orb.RPAInputEnergiesBeta();
647
648 // Static reference contribution of the direct kernel.
649 SetupDirectInteractionOperator(RPAInputEnergiesAlpha, RPAInputEnergiesBeta,
650 0.0);
651
653 configureBSEOperator(Hd_static);
654
655 ExpectationValues expectation_values =
656 ExpectationValue_Operator(orb, Hd_static);
657 Eigen::VectorXd Hd_static_contribution = expectation_values.direct_term;
658
659 if (!orb.getTDAApprox()) {
661 configureBSEOperator(Hd2_static);
662 expectation_values = ExpectationValue_Operator(orb, Hd2_static);
663 Hd_static_contribution += expectation_values.cross_term;
664 }
665
666 Eigen::VectorXd BSEenergies_dynamic = BSEenergies;
667
668 for (Index i_exc = 0; i_exc < BSEenergies.size(); ++i_exc) {
670 << "Dynamical Screening UKS BSE, Excitation " << i_exc << " static "
671 << BSEenergies(i_exc) << flush;
672
673 for (Index iter = 0; iter < opt_.max_dyn_iter; ++iter) {
674 const double old_energy = BSEenergies_dynamic(i_exc);
675
676 SetupDirectInteractionOperator(RPAInputEnergiesAlpha,
677 RPAInputEnergiesBeta, old_energy);
678
680 configureBSEOperator(Hd_dyn);
681
682 expectation_values = ExpectationValue_Operator_State(i_exc, orb, Hd_dyn);
683 Eigen::VectorXd Hd_dynamic_contribution = expectation_values.direct_term;
684
685 if (!orb.getTDAApprox()) {
687 configureBSEOperator(Hd2_dyn);
688 expectation_values =
689 ExpectationValue_Operator_State(i_exc, orb, Hd2_dyn);
690 Hd_dynamic_contribution += expectation_values.cross_term;
691 }
692
693 BSEenergies_dynamic(i_exc) = BSEenergies(i_exc) +
694 Hd_static_contribution(i_exc) -
695 Hd_dynamic_contribution(0);
696
697 XTP_LOG(Log::info, log_) << "Dynamical Screening UKS BSE, excitation "
698 << i_exc << " iteration " << iter << " dynamic "
699 << BSEenergies_dynamic(i_exc) << flush;
700
701 if (std::abs(BSEenergies_dynamic(i_exc) - old_energy) <
702 opt_.dyn_tolerance) {
703 break;
704 }
705 }
706 }
707
708 orb.BSEUKS_dynamic() = BSEenergies_dynamic;
709
710 const double hrt2ev = tools::conv::hrt2ev;
711 const bool has_dipoles = orb.hasTransitionDipoles();
712
713 Eigen::VectorXd oscs = Eigen::VectorXd::Zero(0);
714 if (has_dipoles) {
716 }
717
718 if (orb.getTDAApprox()) {
720 << " ====== combined UKS TDA exciton energies with perturbative "
721 "dynamical screening (eV) ====== "
722 << flush;
723 } else {
725 << " ====== combined UKS full-BSE exciton energies with perturbative "
726 "dynamical screening (eV) ====== "
727 << flush;
728 }
729
730 const Index nprint = std::min<Index>(opt_.nmax, BSEenergies_dynamic.size());
731
732 for (Index i = 0; i < nprint; ++i) {
733 double osc = 0.0;
734 if (has_dipoles && i < oscs.size()) {
735 osc = oscs(i) * BSEenergies_dynamic(i) / BSEenergies(i);
736 }
737
739 << boost::format(
740 " XU(dynamic) = %1$4d Omega = %2$+1.12f eV lambda = %3$+3.2f "
741 "nm"
742 " f = %4$+1.4f") %
743 (i + 1) % (hrt2ev * BSEenergies_dynamic(i)) %
744 (1240.0 / (hrt2ev * BSEenergies_dynamic(i))) % osc
745 << flush;
746 }
747}
748
751
755
761 HdUKSOperator&) const;
763 Hd2UKSOperator&) const;
764
767 const Orbitals&, const ExcitonUKSOperator_TDA&) const;
770 const Orbitals&, const ExcitonUKSOperator_BTDA_B&) const;
772 HdUKSOperator>(const Orbitals&, const HdUKSOperator&) const;
774 Hd2UKSOperator>(const Orbitals&, const Hd2UKSOperator&) const;
775
778 Index, const Orbitals&, const ExcitonUKSOperator_TDA&) const;
781 Index, const Orbitals&, const ExcitonUKSOperator_BTDA_B&) const;
783 HdUKSOperator>(Index, const Orbitals&, const HdUKSOperator&) const;
785 Hd2UKSOperator>(Index, const Orbitals&, const Hd2UKSOperator&) const;
786
787} // namespace xtp
788} // namespace votca
const Eigen::MatrixXd & eigenvectors2() const
Definition eigensystem.h:36
const Eigen::VectorXd & eigenvalues() const
Definition eigensystem.h:30
const Eigen::MatrixXd & eigenvectors() const
Definition eigensystem.h:33
Eigen::VectorXd diagonal() const override
void configure(BSEOperatorUKS_Options opt)
Eigen::MatrixXd dense_matrix() const
const TCMatrix_gwbse_spin & Mmn_raw_
Definition bse_uks.h:80
ExpectationValues ExpectationValue_Operator_State(Index state, const Orbitals &orb, const BSE_OPERATOR &H) const
Definition bse_uks.cc:195
tools::EigenSystem solve_hermitian(BSE_OPERATOR &h) const
Definition bse_uks.cc:453
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
TCMatrix_gwbse_spin Mmn_
Definition bse_uks.h:81
void configureBSEOperator(BSE_OPERATOR &H) const
Definition bse_uks.cc:159
void PrintWeightsUKS(const Eigen::VectorXd &coeffs) const
Definition bse_uks.cc:520
tools::EigenSystem Solve_nonhermitian_Davidson(BSE_OPERATOR_A &Aop, BSE_OPERATOR_B &Bop) const
Definition bse_uks.cc:481
Eigen::MatrixXd Hqp_beta_
Definition bse_uks.h:96
Eigen::VectorXd epsilon_0_inv_
Definition bse_uks.h:94
void Perturbative_DynamicalScreening(Orbitals &orb)
Definition bse_uks.cc:640
ExpectationValues ExpectationValue_Operator(const Orbitals &orb, const BSE_OPERATOR &H) const
Definition bse_uks.cc:171
Eigen::MatrixXd AdjustHqpSize(const Eigen::MatrixXd &Hqp_in, const Eigen::VectorXd &RPAInputEnergies, Index homo) const
Definition bse_uks.cc:92
tools::EigenSystem Solve_excitons_uks_BTDA() const
Definition bse_uks.cc:291
void SetupDirectInteractionOperator(const Eigen::VectorXd &RPAInputEnergiesAlpha, const Eigen::VectorXd &RPAInputEnergiesBeta, double energy)
Definition bse_uks.cc:135
void Solve_excitons_uks(Orbitals &orb) const
Definition bse_uks.cc:442
Eigen::MatrixXd Hqp_alpha_
Definition bse_uks.h:95
void Analyze_excitons_uks(std::vector< QMFragment< BSE_Population > > fragments, const Orbitals &orb) const
Definition bse_uks.cc:597
tools::EigenSystem Solve_excitons_uks_TDA() const
Definition bse_uks.cc:218
Use Davidson algorithm to solve A*V=E*V.
void set_max_search_space(Index N)
Eigen::VectorXd eigenvalues() const
void solve(const MatrixReplacement &A, Index neigen, Index size_initial_guess=0)
void set_matrix_type(std::string mt)
void set_correction(std::string method)
void set_size_update(std::string update_size)
void set_tolerance(std::string tol)
Eigen::MatrixXd eigenvectors() const
Container for molecular orbitals and derived one-particle data.
Definition orbitals.h:47
const Eigen::VectorXd & BSEUKS_dynamic() const
Definition orbitals.h:696
void CalcCoupledTransition_Dipoles()
Definition orbitals.cc:760
bool getTDAApprox() const
Return whether the Tamm-Dancoff approximation is enabled.
Definition orbitals.h:388
const tools::EigenSystem & BSEUKS() const
Definition orbitals.h:690
Eigen::VectorXd Oscillatorstrengths() const
Definition orbitals.cc:643
bool hasTransitionDipoles() const
Report whether transition dipole moments have been computed.
Definition orbitals.h:516
const std::vector< Eigen::Vector3d > & TransitionDipoles() const
Return the stored transition dipole moments.
Definition orbitals.h:521
const Eigen::VectorXd & RPAInputEnergiesBeta() const
Definition orbitals.h:658
void setTDAApprox(bool usedTDA)
Enable or disable the Tamm-Dancoff approximation flag.
Definition orbitals.h:386
const Eigen::VectorXd & RPAInputEnergiesAlpha() const
Definition orbitals.h:655
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
Timestamp returns the current time as a string Example: cout << TimeStamp()
Definition logger.h:224
#define XTP_LOG(level, log)
Definition logger.h:40
const double hrt2ev
Definition constants.h:53
Charge transport classes.
Definition ERIs.h:28
BSE_OPERATOR_UKS< 0, 1, 0, 1 > ExcitonUKSOperator_BTDA_B
BSE_OPERATOR_UKS< 1, 1, 1, 0 > ExcitonUKSOperator_TDA
Eigen::VectorXd ExpValue(const Eigen::MatrixXd &state1, OP OPxstate2)
Definition bse.cc:491
BSE_OPERATOR_UKS< 0, 0, 0, 1 > Hd2UKSOperator
Eigen::MatrixXd BuildFullBSEXRankedInitialGuess(const Eigen::VectorXd &adiag, const Eigen::VectorXd &bdiag, Index nroots)
Build an initial guess for full-BSE Davidson diagonalization.
BSE_OPERATOR_UKS< 0, 0, 1, 0 > HdUKSOperator
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26