votca 2024.1-dev
Loading...
Searching...
No Matches
rpa.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2021 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_RPA_H
22#define VOTCA_XTP_RPA_H
23
24// Standard includes
25#include <vector>
26
27// Local VOTCA includes
28#include "eigen.h"
29#include "logger.h"
30
31namespace votca {
32namespace xtp {
33class TCMatrix_gwbse;
34
35class RPA {
36 public:
37 RPA(Logger& log, const TCMatrix_gwbse& Mmn) : log_(log), Mmn_(Mmn) {};
38
39 void configure(Index homo, Index rpamin, Index rpamax) {
40 homo_ = homo;
41 rpamin_ = rpamin;
42 rpamax_ = rpamax;
43 }
44
45 double getEta() const { return eta_; }
46
47 Eigen::MatrixXd calculate_epsilon_i(double frequency) const {
48 return calculate_epsilon<true>(frequency);
49 }
50
51 Eigen::MatrixXd calculate_epsilon_r(double frequency) const {
52 return calculate_epsilon<false>(frequency);
53 }
54
55 Eigen::MatrixXd calculate_epsilon_r(std::complex<double> frequency) const;
56
57 const Eigen::VectorXd& getRPAInputEnergies() const { return energies_; }
58
59 void setRPAInputEnergies(const Eigen::VectorXd& rpaenergies) {
60 energies_ = rpaenergies;
61 }
62
63 // calculates full RPA vector of energies from gwa and dftenergies and qpmin
64 // RPA energies have three parts, lower than qpmin: dftenergies,between qpmin
65 // and qpmax:gwa_energies,above:dftenergies+homo-lumo shift
66 void UpdateRPAInputEnergies(const Eigen::VectorXd& dftenergies,
67 const Eigen::VectorXd& gwaenergies, Index qpmin);
68
70 Eigen::VectorXd omega; // Eigenvalues
71 Eigen::MatrixXd XpY; // Eigenvector components (X + Y)
72 double ERPA_correlation; // total correlation energy
73 };
74
76
77 private:
78 Index homo_; // HOMO index with respect to dft energies
81 const double eta_ = 0.0001;
82
83 Eigen::VectorXd energies_;
84
87
88 template <bool imag>
89 Eigen::MatrixXd calculate_epsilon(double frequency) const;
90
91 Eigen::VectorXd Calculate_H2p_AmB() const;
92 Eigen::MatrixXd Calculate_H2p_ApB() const;
93 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> Diagonalize_H2p_C(
94 const Eigen::MatrixXd& C) const;
95
96 void ShiftUncorrectedEnergies(const Eigen::VectorXd& dftenergies, Index qpmin,
97 Index gwsize);
98
99 double getMaxCorrection(const Eigen::VectorXd& dftenergies, Index min,
100 Index max) const;
101};
102
103} // namespace xtp
104} // namespace votca
105
106#endif // VOTCA_XTP_RPA_H
Logger is used for thread-safe output of messages.
Definition logger.h:164
RPA(Logger &log, const TCMatrix_gwbse &Mmn)
Definition rpa.h:37
Eigen::MatrixXd Calculate_H2p_ApB() const
Definition rpa.cc:220
Eigen::MatrixXd calculate_epsilon_i(double frequency) const
Definition rpa.h:47
void ShiftUncorrectedEnergies(const Eigen::VectorXd &dftenergies, Index qpmin, Index gwsize)
Definition rpa.cc:48
Index rpamax_
Definition rpa.h:80
void setRPAInputEnergies(const Eigen::VectorXd &rpaenergies)
Definition rpa.h:59
Eigen::VectorXd energies_
Definition rpa.h:83
Eigen::MatrixXd calculate_epsilon(double frequency) const
Definition rpa.cc:74
Eigen::MatrixXd calculate_epsilon_r(double frequency) const
Definition rpa.h:51
const TCMatrix_gwbse & Mmn_
Definition rpa.h:86
const double eta_
Definition rpa.h:81
Index rpamin_
Definition rpa.h:79
Eigen::VectorXd Calculate_H2p_AmB() const
Definition rpa.cc:205
void configure(Index homo, Index rpamin, Index rpamax)
Definition rpa.h:39
double getEta() const
Definition rpa.h:45
const Eigen::VectorXd & getRPAInputEnergies() const
Definition rpa.h:57
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > Diagonalize_H2p_C(const Eigen::MatrixXd &C) const
Definition rpa.cc:243
Logger & log_
Definition rpa.h:85
rpa_eigensolution Diagonalize_H2p() const
Definition rpa.cc:157
void UpdateRPAInputEnergies(const Eigen::VectorXd &dftenergies, const Eigen::VectorXd &gwaenergies, Index qpmin)
Definition rpa.cc:30
double getMaxCorrection(const Eigen::VectorXd &dftenergies, Index min, Index max) const
Definition rpa.cc:63
Index homo_
Definition rpa.h:78
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26