votca 2024.2-dev
Loading...
Searching...
No Matches
bse.cc
Go to the documentation of this file.
1/*
2 * Copyright 2009-2020 The VOTCA Development Team
3 * (http://www.votca.org)
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License")
6 *
7 * You may not use this file except in compliance with the License.
8 * You may obtain a copy of the License at
9 *
10 * http://www.apache.org/licenses/LICENSE-2.0
11 *
12 * Unless required by applicable law or agreed to in writing, software
13 * distributed under the License is distributed on an "AS IS" BASIS,
14 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
15 * See the License for the specific language governing permissions and
16 * limitations under the License.
17 *
18 */
19
20// Standard includes
21#include <chrono>
22#include <iostream>
23
24// VOTCA includes
25#include <votca/tools/linalg.h>
26
27// Local VOTCA includes
28#include "votca/xtp/bse.h"
34#include "votca/xtp/rpa.h"
35#include "votca/xtp/vc2index.h"
36
37using boost::format;
38using std::flush;
39
40namespace votca {
41namespace xtp {
42
43void BSE::configure(const options& opt, const Eigen::VectorXd& RPAInputEnergies,
44 const Eigen::MatrixXd& Hqp_in) {
45 opt_ = opt;
47 bse_cmin_ = opt_.homo + 1;
54 Hqp_ = AdjustHqpSize(Hqp_in, RPAInputEnergies);
55 } else {
56 Hqp_ = AdjustHqpSize(Hqp_in, RPAInputEnergies).diagonal().asDiagonal();
57 }
58 SetupDirectInteractionOperator(RPAInputEnergies, 0.0);
59}
60
61Eigen::MatrixXd BSE::AdjustHqpSize(const Eigen::MatrixXd& Hqp,
62 const Eigen::VectorXd& RPAInputEnergies) {
63
64 Index hqp_size = bse_vtotal_ + bse_ctotal_;
65 Index gwsize = opt_.qpmax - opt_.qpmin + 1;
66 Index RPAoffset = opt_.vmin - opt_.rpamin;
67 Eigen::MatrixXd Hqp_BSE = Eigen::MatrixXd::Zero(hqp_size, hqp_size);
68
69 if (opt_.vmin >= opt_.qpmin) {
70 Index start = opt_.vmin - opt_.qpmin;
71 if (opt_.cmax <= opt_.qpmax) {
72 Hqp_BSE = Hqp.block(start, start, hqp_size, hqp_size);
73 } else {
74 Index virtoffset = gwsize - start;
75 Hqp_BSE.topLeftCorner(virtoffset, virtoffset) =
76 Hqp.block(start, start, virtoffset, virtoffset);
77
78 Index virt_extra = opt_.cmax - opt_.qpmax;
79 Hqp_BSE.diagonal().tail(virt_extra) =
80 RPAInputEnergies.segment(RPAoffset + virtoffset, virt_extra);
81 }
82 }
83
84 if (opt_.vmin < opt_.qpmin) {
85 Index occ_extra = opt_.qpmin - opt_.vmin;
86 Hqp_BSE.diagonal().head(occ_extra) =
87 RPAInputEnergies.segment(RPAoffset, occ_extra);
88
89 Hqp_BSE.block(occ_extra, occ_extra, gwsize, gwsize) = Hqp;
90
91 if (opt_.cmax > opt_.qpmax) {
92 Index virtoffset = occ_extra + gwsize;
93 Index virt_extra = opt_.cmax - opt_.qpmax;
94 Hqp_BSE.diagonal().tail(virt_extra) =
95 RPAInputEnergies.segment(RPAoffset + virtoffset, virt_extra);
96 }
97 }
98
99 return Hqp_BSE;
100}
101
103 const Eigen::VectorXd& RPAInputEnergies, double energy) {
104 RPA rpa = RPA(log_, Mmn_);
106 rpa.setRPAInputEnergies(RPAInputEnergies);
107
108 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es(
109 rpa.calculate_epsilon_r(energy));
110 Mmn_.MultiplyRightWithAuxMatrix(es.eigenvectors());
111
112 epsilon_0_inv_ = Eigen::VectorXd::Zero(es.eigenvalues().size());
113 for (Index i = 0; i < es.eigenvalues().size(); ++i) {
114 if (es.eigenvalues()(i) > 1e-8) {
115 epsilon_0_inv_(i) = 1 / es.eigenvalues()(i);
116 }
117 }
118}
119
120template <typename BSE_OPERATOR>
123 opt.cmax = opt_.cmax;
124 opt.homo = opt_.homo;
125 opt.qpmin = opt_.qpmin;
126 opt.rpamin = opt_.rpamin;
127 opt.vmin = opt_.vmin;
128 H.configure(opt);
129}
130
137
140 if (opt_.useTDA) {
142 } else {
144 }
146}
147
150 if (opt_.useTDA) {
152 } else {
154 }
155}
156
158
162 << TimeStamp() << " Setup TDA singlet hamiltonian " << flush;
163 return solve_hermitian(Hs);
164}
165
172
179
180template <typename BSE_OPERATOR>
182
183 std::chrono::time_point<std::chrono::system_clock> start =
184 std::chrono::system_clock::now();
185
186 tools::EigenSystem result;
187
189
195 DS.solve(h, opt_.nmax);
196 result.eigenvalues() = DS.eigenvalues();
197 result.eigenvectors() = DS.eigenvectors();
198
199 std::chrono::time_point<std::chrono::system_clock> end =
200 std::chrono::system_clock::now();
201 std::chrono::duration<double> elapsed_time = end - start;
202
203 XTP_LOG(Log::info, log_) << TimeStamp() << " Diagonalization done in "
204 << elapsed_time.count() << " secs" << flush;
205
206 return result;
207}
208
218
225 << TimeStamp() << " Setup Full triplet hamiltonian " << flush;
226 return Solve_nonhermitian_Davidson(A, B);
227}
228
229template <typename BSE_OPERATOR_A, typename BSE_OPERATOR_B>
231 BSE_OPERATOR_B& Bop) const {
232 std::chrono::time_point<std::chrono::system_clock> start =
233 std::chrono::system_clock::now();
234
235 // operator
237
238 // Davidson solver
245 DS.set_matrix_type("HAM");
246 DS.solve(Hop, opt_.nmax);
247
248 // results
249 tools::EigenSystem result;
250 result.eigenvalues() = DS.eigenvalues();
251 Eigen::MatrixXd tmpX = DS.eigenvectors().topRows(Aop.rows());
252 Eigen::MatrixXd tmpY = DS.eigenvectors().bottomRows(Bop.rows());
253
254 // // normalization so that eigenvector^2 - eigenvector2^2 = 1
255 Eigen::VectorXd normX = tmpX.colwise().squaredNorm();
256 Eigen::VectorXd normY = tmpY.colwise().squaredNorm();
257
258 Eigen::ArrayXd sqinvnorm = (normX - normY).array().inverse().cwiseSqrt();
259
260 result.eigenvectors() = tmpX * sqinvnorm.matrix().asDiagonal();
261 result.eigenvectors2() = tmpY * sqinvnorm.matrix().asDiagonal();
262
263 std::chrono::time_point<std::chrono::system_clock> end =
264 std::chrono::system_clock::now();
265 std::chrono::duration<double> elapsed_time = end - start;
266
267 XTP_LOG(Log::info, log_) << TimeStamp() << " Diagonalization done in "
268 << elapsed_time.count() << " secs" << flush;
269
270 return result;
271}
272
273void BSE::printFragInfo(const std::vector<QMFragment<BSE_Population> >& frags,
274 Index state) const {
275 for (const QMFragment<BSE_Population>& frag : frags) {
276 double dq = frag.value().H[state] + frag.value().E[state];
277 double qeff = dq + frag.value().Gs;
279 << format(
280 " Fragment %1$4d -- hole: %2$5.1f%% electron: "
281 "%3$5.1f%% dQ: %4$+5.2f Qeff: %5$+5.2f") %
282 int(frag.getId()) % (100.0 * frag.value().H[state]) %
283 (-100.0 * frag.value().E[state]) % dq % qeff
284 << flush;
285 }
286 return;
287}
288
289void BSE::PrintWeights(const Eigen::VectorXd& weights) const {
291 for (Index i_bse = 0; i_bse < bse_size_; ++i_bse) {
292 double weight = weights(i_bse);
293 if (weight > opt_.min_print_weight) {
295 << format(" HOMO-%1$-3d -> LUMO+%2$-3d : %3$3.1f%%") %
296 (opt_.homo - vc.v(i_bse)) % (vc.c(i_bse) - opt_.homo - 1) %
297 (100.0 * weight)
298 << flush;
299 }
300 }
301 return;
302}
303
305 const Orbitals& orb) const {
306
308
309 Eigen::VectorXd oscs = orb.Oscillatorstrengths();
310
311 Interaction act = Analyze_eh_interaction(singlet, orb);
312
313 if (fragments.size() > 0) {
314 Lowdin low;
315 low.CalcChargeperFragment(fragments, orb, singlet);
316 }
317
318 const Eigen::VectorXd& energies = orb.BSESinglets().eigenvalues();
319
320 double hrt2ev = tools::conv::hrt2ev;
322 << " ====== singlet energies (eV) ====== " << flush;
323 for (Index i = 0; i < opt_.nmax; ++i) {
324 Eigen::VectorXd weights =
325 orb.BSESinglets().eigenvectors().col(i).cwiseAbs2();
326 if (!orb.getTDAApprox()) {
327 weights -= orb.BSESinglets().eigenvectors2().col(i).cwiseAbs2();
328 }
329
330 double osc = oscs[i];
332 << format(
333 " S = %1$4d Omega = %2$+1.12f eV lamdba = %3$+3.2f nm <FT> "
334 "= %4$+1.4f <K_x> = %5$+1.4f <K_d> = %6$+1.4f") %
335 (i + 1) % (hrt2ev * energies(i)) %
336 (1240.0 / (hrt2ev * energies(i))) %
337 (hrt2ev * act.qp_contrib(i)) %
338 (hrt2ev * act.exchange_contrib(i)) %
339 (hrt2ev * act.direct_contrib(i))
340 << flush;
341
342 const Eigen::Vector3d& trdip = orb.TransitionDipoles()[i];
344 << format(
345 " TrDipole length gauge[e*bohr] dx = %1$+1.4f dy = "
346 "%2$+1.4f dz = %3$+1.4f |d|^2 = %4$+1.4f f = %5$+1.4f") %
347 trdip[0] % trdip[1] % trdip[2] % (trdip.squaredNorm()) % osc
348 << flush;
349
350 PrintWeights(weights);
351 if (fragments.size() > 0) {
352 printFragInfo(fragments, i);
353 }
354
355 XTP_LOG(Log::error, log_) << flush;
356 }
357 return;
358}
359
361 const Orbitals& orb) const {
362
364 Interaction act = Analyze_eh_interaction(triplet, orb);
365
366 if (fragments.size() > 0) {
367 Lowdin low;
368 low.CalcChargeperFragment(fragments, orb, triplet);
369 }
370
371 const Eigen::VectorXd& energies = orb.BSETriplets().eigenvalues();
373 << " ====== triplet energies (eV) ====== " << flush;
374 for (Index i = 0; i < opt_.nmax; ++i) {
375 Eigen::VectorXd weights =
376 orb.BSETriplets().eigenvectors().col(i).cwiseAbs2();
377 if (!orb.getTDAApprox()) {
378 weights -= orb.BSETriplets().eigenvectors2().col(i).cwiseAbs2();
379 }
380
382 << format(
383 " T = %1$4d Omega = %2$+1.12f eV lamdba = %3$+3.2f nm <FT> "
384 "= %4$+1.4f <K_d> = %5$+1.4f") %
385 (i + 1) % (tools::conv::hrt2ev * energies(i)) %
386 (1240.0 / (tools::conv::hrt2ev * energies(i))) %
389 << flush;
390
391 PrintWeights(weights);
392 if (fragments.size() > 0) {
393 printFragInfo(fragments, i);
394 }
395 XTP_LOG(Log::error, log_) << format(" ") << flush;
396 }
397
398 return;
399}
400template <class OP>
401Eigen::VectorXd ExpValue(const Eigen::MatrixXd& state1, OP OPxstate2) {
402 return state1.cwiseProduct(OPxstate2.eval()).colwise().sum().transpose();
403}
404
405Eigen::VectorXd ExpValue(const Eigen::MatrixXd& state1,
406 const Eigen::MatrixXd& OPxstate2) {
407 return state1.cwiseProduct(OPxstate2).colwise().sum().transpose();
408}
409
410template <typename BSE_OPERATOR>
412 const QMStateType& type, const Orbitals& orb, const BSE_OPERATOR& H) const {
413
414 const tools::EigenSystem& BSECoefs =
415 (type == QMStateType::Singlet) ? orb.BSESinglets() : orb.BSETriplets();
416
417 ExpectationValues expectation_values;
418
419 const Eigen::MatrixXd temp = H * BSECoefs.eigenvectors();
420
421 expectation_values.direct_term = ExpValue(BSECoefs.eigenvectors(), temp);
422 if (!orb.getTDAApprox()) {
423 expectation_values.direct_term +=
424 ExpValue(BSECoefs.eigenvectors2(), H * BSECoefs.eigenvectors2());
425 expectation_values.cross_term =
426 2 * ExpValue(BSECoefs.eigenvectors2(), temp);
427 } else {
428 expectation_values.cross_term = Eigen::VectorXd::Zero(0);
429 }
430 return expectation_values;
431}
432
433template <typename BSE_OPERATOR>
435 const QMState& state, const Orbitals& orb, const BSE_OPERATOR& H) const {
436
437 const tools::EigenSystem& BSECoefs = (state.Type() == QMStateType::Singlet)
438 ? orb.BSESinglets()
439 : orb.BSETriplets();
440
441 ExpectationValues expectation_values;
442
443 const Eigen::MatrixXd BSECoefs_state =
444 BSECoefs.eigenvectors().col(state.StateIdx());
445
446 const Eigen::MatrixXd temp = H * BSECoefs_state;
447
448 expectation_values.direct_term = ExpValue(BSECoefs_state, temp);
449
450 if (!orb.getTDAApprox()) {
451 const Eigen::MatrixXd BSECoefs2_state =
452 BSECoefs.eigenvectors2().col(state.StateIdx());
453
454 expectation_values.direct_term +=
455 ExpValue(BSECoefs2_state, H * BSECoefs2_state);
456 expectation_values.cross_term = 2 * ExpValue(BSECoefs2_state, temp);
457 } else {
458 expectation_values.cross_term = Eigen::VectorXd::Zero(0);
459 }
460
461 return expectation_values;
462}
463
464// Composition of the excitation energy in terms of QP, direct (screened),
465// and exchance contributions in the BSE
466// Full BSE:
467//
468// | A* | | H K | | A |
469// | -B* | | -K -H | | B | = A*.H.A + B*.H.B + 2A*.K.B
470//
471// with: H = H_qp + H_d + eta.H_x
472// K = H_d2 + eta.H_x
473//
474// reports composition for FULL BSE as
475// <FT> = A*.H_qp.A + B*.H_qp.B
476// <Kx> = eta.(A*.H_x.A + B*.H_x.B + 2A*.H_x.B)
477// <Kd> = A*.H_d.A + B*.H_d.B + 2A*.H_d2.B
479 const Orbitals& orb) const {
480 Interaction analysis;
481 {
484 ExpectationValues expectation_values =
485 ExpectationValue_Operator(type, orb, hqp);
486 analysis.qp_contrib = expectation_values.direct_term;
487 }
488 {
491 ExpectationValues expectation_values =
492 ExpectationValue_Operator(type, orb, hd);
493 analysis.direct_contrib = expectation_values.direct_term;
494 }
495 if (!orb.getTDAApprox()) {
498 ExpectationValues expectation_values =
499 ExpectationValue_Operator(type, orb, hd2);
500 analysis.direct_contrib += expectation_values.cross_term;
501 }
502
503 if (type == QMStateType::Singlet) {
506 ExpectationValues expectation_values =
507 ExpectationValue_Operator(type, orb, hx);
508 analysis.exchange_contrib = 2.0 * expectation_values.direct_term;
509 if (!orb.getTDAApprox()) {
510 analysis.exchange_contrib += 2.0 * expectation_values.cross_term;
511 }
512 } else {
513 analysis.exchange_contrib = Eigen::VectorXd::Zero(0);
514 }
515
516 return analysis;
517}
518
519// Dynamical Screening in BSE as perturbation to static excitation energies
520// as in Phys. Rev. B 80, 241405 (2009) for the TDA case
522 Orbitals& orb) {
523
524 const tools::EigenSystem& BSECoefs =
525 (type == QMStateType::Singlet) ? orb.BSESinglets() : orb.BSETriplets();
526
527 const Eigen::VectorXd& RPAInputEnergies = orb.RPAInputEnergies();
528
529 // static case as reference
530 SetupDirectInteractionOperator(RPAInputEnergies, 0.0);
531 HdOperator Hd_static(epsilon_0_inv_, Mmn_, Hqp_);
532 configureBSEOperator(Hd_static);
533 ExpectationValues expectation_values =
534 ExpectationValue_Operator(type, orb, Hd_static);
535 Eigen::VectorXd Hd_static_contribution = expectation_values.direct_term;
536 if (!orb.getTDAApprox()) {
537 Hd2Operator Hd2_static(epsilon_0_inv_, Mmn_, Hqp_);
538 configureBSEOperator(Hd2_static);
539 expectation_values = ExpectationValue_Operator(type, orb, Hd2_static);
540 Hd_static_contribution += expectation_values.cross_term;
541 }
542
543 const Eigen::VectorXd& BSEenergies = BSECoefs.eigenvalues();
544
545 // initial copy of static BSE energies to dynamic
546 Eigen::VectorXd BSEenergies_dynamic = BSEenergies;
547
548 // recalculate Hd at the various energies
549 for (Index i_exc = 0; i_exc < BSEenergies.size(); i_exc++) {
550 XTP_LOG(Log::info, log_) << "Dynamical Screening BSE, Excitation " << i_exc
551 << " static " << BSEenergies(i_exc) << flush;
552
553 for (Index iter = 0; iter < max_dyn_iter_; iter++) {
554
555 // setup the direct operator with the last energy as screening frequency
556 double old_energy = BSEenergies_dynamic(i_exc);
557 SetupDirectInteractionOperator(RPAInputEnergies, old_energy);
559 configureBSEOperator(Hd_dyn);
560
561 // get the contribution of Hd for the dynamic case
562 QMState state(type, i_exc, false);
563 expectation_values = ExpectationValue_Operator_State(state, orb, Hd_dyn);
564 Eigen::VectorXd Hd_dynamic_contribution = expectation_values.direct_term;
565 if (!orb.getTDAApprox()) {
567 configureBSEOperator(Hd2_dyn);
568 expectation_values =
569 ExpectationValue_Operator_State(state, orb, Hd2_dyn);
570 Hd_dynamic_contribution += expectation_values.cross_term;
571 }
572
573 // new energy perturbatively
574 BSEenergies_dynamic(i_exc) = BSEenergies(i_exc) +
575 Hd_static_contribution(i_exc) -
576 Hd_dynamic_contribution(0);
577
579 << "Dynamical Screening BSE, excitation " << i_exc << " iteration "
580 << iter << " dynamic " << BSEenergies_dynamic(i_exc) << flush;
581
582 // check tolerance
583 if (std::abs(BSEenergies_dynamic(i_exc) - old_energy) < dyn_tolerance_) {
584 break;
585 }
586 }
587 }
588
589 double hrt2ev = tools::conv::hrt2ev;
590
591 if (type == QMStateType::Singlet) {
592 orb.BSESinglets_dynamic() = BSEenergies_dynamic;
593 XTP_LOG(Log::error, log_) << " ====== singlet energies with perturbative "
594 "dynamical screening (eV) ====== "
595 << flush;
596 Eigen::VectorXd oscs = orb.Oscillatorstrengths();
597 for (Index i = 0; i < opt_.nmax; ++i) {
598 double osc = oscs[i];
600 << format(
601 " S(dynamic) = %1$4d Omega = %2$+1.12f eV lamdba = %3$+3.2f "
602 "nm f "
603 "= %4$+1.4f") %
604 (i + 1) % (hrt2ev * BSEenergies_dynamic(i)) %
605 (1240.0 / (hrt2ev * BSEenergies_dynamic(i))) %
606 (osc * BSEenergies_dynamic(i) / BSEenergies(i))
607 << flush;
608 }
609
610 } else {
611 orb.BSETriplets_dynamic() = BSEenergies_dynamic;
612 XTP_LOG(Log::error, log_) << " ====== triplet energies with perturbative "
613 "dynamical screening (eV) ====== "
614 << flush;
615 for (Index i = 0; i < opt_.nmax; ++i) {
617 << format(
618 " T(dynamic) = %1$4d Omega = %2$+1.12f eV lamdba = %3$+3.2f "
619 "nm ") %
620 (i + 1) % (hrt2ev * BSEenergies_dynamic(i)) %
621 (1240.0 / (hrt2ev * BSEenergies_dynamic(i)))
622 << flush;
623 }
624 }
625}
626
627} // namespace xtp
628} // 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
void Solve_singlets(Orbitals &orb) const
Definition bse.cc:138
void printFragInfo(const std::vector< QMFragment< BSE_Population > > &frags, Index state) const
Definition bse.cc:273
void PrintWeights(const Eigen::VectorXd &weights) const
Definition bse.cc:289
double dyn_tolerance_
Definition bse.h:107
void Analyze_singlets(std::vector< QMFragment< BSE_Population > > fragments, const Orbitals &orb) const
Definition bse.cc:304
tools::EigenSystem Solve_triplets_BTDA() const
Definition bse.cc:219
TripletOperator_TDA getTripletOperator_TDA() const
Definition bse.cc:173
SingletOperator_TDA getSingletOperator_TDA() const
Definition bse.cc:166
ExpectationValues ExpectationValue_Operator_State(const QMState &state, const Orbitals &orb, const BSE_OPERATOR &H) const
Definition bse.cc:434
options opt_
Definition bse.h:86
tools::EigenSystem Solve_singlets_BTDA() const
Definition bse.cc:209
Index bse_vmax_
Definition bse.h:100
Index bse_cmin_
Definition bse.h:101
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 SetupDirectInteractionOperator(const Eigen::VectorXd &DFTenergies, double energy)
Definition bse.cc:102
void configure(const options &opt, const Eigen::VectorXd &RPAEnergies, const Eigen::MatrixXd &Hqp_in)
Definition bse.cc:43
Eigen::MatrixXd Hqp_
Definition bse.h:112
Eigen::VectorXd epsilon_0_inv_
Definition bse.h:109
ExpectationValues ExpectationValue_Operator(const QMStateType &type, const Orbitals &orb, const BSE_OPERATOR &H) const
Definition bse.cc:411
Interaction Analyze_eh_interaction(const QMStateType &type, const Orbitals &orb) const
Definition bse.cc:478
Index bse_ctotal_
Definition bse.h:104
TCMatrix_gwbse & Mmn_
Definition bse.h:111
Index bse_size_
Definition bse.h:102
tools::EigenSystem Solve_singlets_TDA() const
Definition bse.cc:157
void configureBSEOperator(BSE_OPERATOR &H) const
Definition bse.cc:121
Index max_dyn_iter_
Definition bse.h:106
tools::EigenSystem solve_hermitian(BSE_OPERATOR &h) const
Definition bse.cc:181
void Solve_triplets(Orbitals &orb) const
Definition bse.cc:148
Logger & log_
Definition bse.h:99
tools::EigenSystem Solve_nonhermitian_Davidson(BSE_OPERATOR_A &Aop, BSE_OPERATOR_B &Bop) const
Definition bse.cc:230
Eigen::MatrixXd AdjustHqpSize(const Eigen::MatrixXd &Hqp_in, const Eigen::VectorXd &RPAInputEnergies)
Definition bse.cc:61
Index bse_vtotal_
Definition bse.h:103
tools::EigenSystem Solve_triplets_TDA() const
Definition bse.cc:131
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
Definition orbitals.h:46
const tools::EigenSystem & BSETriplets() const
Definition orbitals.h:324
void CalcCoupledTransition_Dipoles()
Definition orbitals.cc:526
const tools::EigenSystem & BSESinglets() const
Definition orbitals.h:334
bool getTDAApprox() const
Definition orbitals.h:269
Eigen::VectorXd Oscillatorstrengths() const
Definition orbitals.cc:428
const std::vector< Eigen::Vector3d > & TransitionDipoles() const
Definition orbitals.h:369
const Eigen::VectorXd & RPAInputEnergies() const
Definition orbitals.h:299
const Eigen::VectorXd & BSETriplets_dynamic() const
Definition orbitals.h:355
void setTDAApprox(bool usedTDA)
Definition orbitals.h:268
const Eigen::VectorXd & BSESinglets_dynamic() const
Definition orbitals.h:343
void CalcChargeperFragment(std::vector< QMFragment< BSE_Population > > &frags, const Orbitals &orbitals, QMStateType type) const
Identifier for QMstates. Strings like S1 are converted into enum +zero indexed int.
Definition qmstate.h:132
Index StateIdx() const
Definition qmstate.h:154
const QMStateType & Type() const
Definition qmstate.h:151
void setRPAInputEnergies(const Eigen::VectorXd &rpaenergies)
Definition rpa.h:59
Eigen::MatrixXd calculate_epsilon_r(double frequency) const
Definition rpa.h:51
void configure(Index homo, Index rpamin, Index rpamax)
Definition rpa.h:39
void MultiplyRightWithAuxMatrix(const Eigen::MatrixXd &matrix)
Timestamp returns the current time as a string Example: cout << TimeStamp()
Definition logger.h:224
This class transforms a pair of indices v,c into a compound index I, via I=ctotal*v+c the fast dimens...
Definition vc2index.h:36
Index c(Index index) const
Definition vc2index.h:48
Index v(Index index) const
Definition vc2index.h:46
#define XTP_LOG(level, log)
Definition logger.h:40
const double hrt2ev
Definition constants.h:53
Eigen::VectorXd ExpValue(const Eigen::MatrixXd &state1, OP OPxstate2)
Definition bse.cc:401
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
Eigen::VectorXd direct_term
Definition bse.h:95
Eigen::VectorXd cross_term
Definition bse.h:96
Eigen::VectorXd qp_contrib
Definition bse.h:91
Eigen::VectorXd exchange_contrib
Definition bse.h:89
Eigen::VectorXd direct_contrib
Definition bse.h:90
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