votca 2024-dev
Loading...
Searching...
No Matches
bse.h
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#pragma once
21#ifndef VOTCA_XTP_BSE_H
22#define VOTCA_XTP_BSE_H
23
24// Local VOTCA includes
25#include "logger.h"
26#include "orbitals.h"
27#include "qmstate.h"
28#include "threecenter.h"
29
30namespace votca {
31namespace xtp {
32struct BSE_Population;
33template <Index cqp, Index cx, Index cd, Index cd2>
34class BSE_OPERATOR;
37template <class T>
38class QMFragment;
39
40class BSE {
41
42 public:
43 // BSE(Logger& log, TCMatrix_gwbse& Mmn, const Eigen::MatrixXd& Hqp_in)
44 // : log_(log), Mmn_(Mmn), Hqp_in_(Hqp_in){};
45 BSE(Logger& log, TCMatrix_gwbse& Mmn) : log_(log), Mmn_(Mmn){};
46
47 struct options {
48 bool useTDA;
56 Index nmax; // number of eigenvectors to calculat
58 std::string davidson_tolerance;
59 std::string davidson_update;
61 double min_print_weight; // minimium contribution for state to print it
65 };
66
67 void configure(const options& opt, const Eigen::VectorXd& RPAEnergies,
68 const Eigen::MatrixXd& Hqp_in);
69
70 void Solve_singlets(Orbitals& orb) const;
71 void Solve_triplets(Orbitals& orb) const;
72
73 Eigen::MatrixXd getHqp() const { return Hqp_; };
74
77
78 void Analyze_singlets(std::vector<QMFragment<BSE_Population> > fragments,
79 const Orbitals& orb) const;
80 void Analyze_triplets(std::vector<QMFragment<BSE_Population> > fragments,
81 const Orbitals& orb) const;
82
84
85 private:
87
88 struct Interaction {
89 Eigen::VectorXd exchange_contrib;
90 Eigen::VectorXd direct_contrib;
91 Eigen::VectorXd qp_contrib;
92 };
93
95 Eigen::VectorXd direct_term;
96 Eigen::VectorXd cross_term;
97 };
98
105
108
109 Eigen::VectorXd epsilon_0_inv_;
110
112 Eigen::MatrixXd Hqp_;
113
116
119
120 void PrintWeights(const Eigen::VectorXd& weights) const;
121
122 template <typename BSE_OPERATOR>
124
125 template <typename BSE_OPERATOR>
127
128 template <typename BSE_OPERATOR_ApB, typename BSE_OPERATOR_AmB>
130 BSE_OPERATOR_AmB&) const;
131
132 template <typename BSE_OPERATOR_A, typename BSE_OPERATOR_B>
134 BSE_OPERATOR_B& Bop) const;
135
136 void printFragInfo(const std::vector<QMFragment<BSE_Population> >& frags,
137 Index state) const;
138 void printWeights(Index i_bse, double weight) const;
139 void SetupDirectInteractionOperator(const Eigen::VectorXd& DFTenergies,
140 double energy);
141
142 Eigen::MatrixXd AdjustHqpSize(const Eigen::MatrixXd& Hqp_in,
143 const Eigen::VectorXd& RPAInputEnergies);
144
146 const Orbitals& orb) const;
147 template <typename BSE_OPERATOR>
149 const Orbitals& orb,
150 const BSE_OPERATOR& H) const;
151
152 template <typename BSE_OPERATOR>
154 const QMState& state, const Orbitals& orb, const BSE_OPERATOR& H) const;
155};
156} // namespace xtp
157} // namespace votca
158
159#endif // VOTCA_XTP_BSE_H
void Solve_singlets(Orbitals &orb) const
Definition bse.cc:138
BSE(Logger &log, TCMatrix_gwbse &Mmn)
Definition bse.h:45
Eigen::MatrixXd getHqp() const
Definition bse.h:73
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
void printWeights(Index i_bse, double weight) const
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
tools::EigenSystem Solve_nonhermitian(BSE_OPERATOR_ApB &apb, BSE_OPERATOR_AmB &) const
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
Logger is used for thread-safe output of messages.
Definition logger.h:164
container for molecular orbitals
Definition orbitals.h:46
Identifier for QMstates. Strings like S1 are converted into enum +zero indexed int.
Definition qmstate.h:132
BSE_OPERATOR< 1, 2, 1, 0 > SingletOperator_TDA
Definition bse.h:35
BSE_OPERATOR< 1, 0, 1, 0 > TripletOperator_TDA
Definition bse.h:36
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