votca 2026-dev
Loading...
Searching...
No Matches
rpa.cc
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// Local VOTCA includes
21#include "votca/xtp/rpa.h"
22#include "votca/xtp/aomatrix.h"
25#include "votca/xtp/vc2index.h"
26#include <iomanip>
27#include <sstream>
28
29namespace votca {
30namespace xtp {
31
32void RPA::UpdateRPAInputEnergies(const Eigen::VectorXd& dftenergies,
33 const Eigen::VectorXd& gwaenergies,
34 Index qpmin) {
35 Index rpatotal = rpamax_ - rpamin_ + 1;
36 energies_ = dftenergies.segment(rpamin_, rpatotal);
37 Index gwsize = Index(gwaenergies.size());
38
39 energies_.segment(qpmin - rpamin_, gwsize) = gwaenergies;
40
41 ShiftUncorrectedEnergies(dftenergies, qpmin, gwsize);
42}
43
44// Shifts energies of levels that are not QP corrected but
45// used in the RPA:
46// between rpamin and qpmin: by maximum abs of explicit QP corrections
47// from qpmin to HOMO
48// between qpmax and rpamax: by maximum abs of explicit QP corrections
49// from LUMO to qpmax
50void RPA::ShiftUncorrectedEnergies(const Eigen::VectorXd& dftenergies,
51 Index qpmin, Index gwsize) {
52
53 Index lumo = homo_ + 1;
54 Index qpmax = qpmin + gwsize - 1;
55
56 // get max abs QP corrections for occupied/virtual levels
57 double max_correction_occ = getMaxCorrection(dftenergies, qpmin, homo_);
58 double max_correction_virt = getMaxCorrection(dftenergies, lumo, qpmax);
59
60 // shift energies
61 energies_.head(qpmin).array() -= max_correction_occ;
62 energies_.tail(rpamax_ - qpmax).array() += max_correction_virt;
63}
64
65double RPA::getMaxCorrection(const Eigen::VectorXd& dftenergies, Index min,
66 Index max) const {
67
68 Index range = max - min + 1;
69 Eigen::VectorXd corrections =
70 energies_.segment(min, range) - dftenergies.segment(min - rpamin_, range);
71
72 return (corrections.cwiseAbs()).maxCoeff();
73}
74
75template <bool imag>
76Eigen::MatrixXd RPA::calculate_epsilon(double frequency) const {
77 const Index size = Mmn_.auxsize();
78
79 const Index lumo = homo_ + 1;
80 const Index n_occ = lumo - rpamin_;
81 const Index n_unocc = rpamax_ - lumo + 1;
82 const double freq2 = frequency * frequency;
83 const double eta2 = eta_ * eta_;
84
85 OpenMP_CUDA transform;
86 transform.createTemporaries(n_unocc, size);
87
88#pragma omp parallel
89 {
90 Index threadid = OPENMP::getThreadId();
91#pragma omp for schedule(dynamic)
92 for (Index m_level = 0; m_level < n_occ; m_level++) {
93 const double qp_energy_m = energies_(m_level);
94
95 Eigen::MatrixXd Mmn_RPA = Mmn_[m_level].bottomRows(n_unocc);
96 transform.PushMatrix(Mmn_RPA, threadid);
97 const Eigen::ArrayXd deltaE =
98 energies_.tail(n_unocc).array() - qp_energy_m;
99 Eigen::VectorXd denom;
100 if (imag) {
101 denom = 4 * deltaE / (deltaE.square() + freq2);
102 } else {
103 Eigen::ArrayXd deltEf = deltaE - frequency;
104 Eigen::ArrayXd sum = deltEf / (deltEf.square() + eta2);
105 deltEf = deltaE + frequency;
106 sum += deltEf / (deltEf.square() + eta2);
107 denom = 2 * sum;
108 }
109
110 transform.A_TDA(denom, threadid);
111 }
112 }
113 Eigen::MatrixXd result = transform.getReductionVar();
114 result.diagonal().array() += 1.0;
115 return result;
116}
117
118template Eigen::MatrixXd RPA::calculate_epsilon<true>(double frequency) const;
119template Eigen::MatrixXd RPA::calculate_epsilon<false>(double frequency) const;
120
121Eigen::MatrixXd RPA::calculate_epsilon_r(std::complex<double> frequency) const {
122
123 const Index size = Mmn_.auxsize();
124
125 const Index lumo = homo_ + 1;
126 const Index n_occ = lumo - rpamin_;
127 const Index n_unocc = rpamax_ - lumo + 1;
128 OpenMP_CUDA transform;
129 transform.createTemporaries(n_unocc, size);
130#pragma omp parallel
131 {
132 Index threadid = OPENMP::getThreadId();
133#pragma omp for schedule(dynamic)
134 for (Index m_level = 0; m_level < n_occ; m_level++) {
135
136 const double qp_energy_m = energies_(m_level);
137 Eigen::MatrixXd Mmn_RPA = Mmn_[m_level].bottomRows(n_unocc);
138 transform.PushMatrix(Mmn_RPA, threadid);
139 const Eigen::ArrayXd deltaE =
140 energies_.tail(n_unocc).array() - qp_energy_m;
141
142 Eigen::ArrayXd deltaEm = frequency.real() - deltaE;
143 Eigen::ArrayXd deltaEp = frequency.real() + deltaE;
144
145 double sigma_1 = std::pow(frequency.imag() + eta_, 2);
146 double sigma_2 = std::pow(frequency.imag() - eta_, 2);
147
148 Eigen::VectorXd chi =
149 deltaEm * (deltaEm.cwiseAbs2() + sigma_1).cwiseInverse() -
150 deltaEp * (deltaEp.cwiseAbs2() + sigma_2).cwiseInverse();
151 transform.A_TDA(chi, threadid);
152 }
153 }
154 Eigen::MatrixXd result = -2 * transform.getReductionVar();
155 result.diagonal().array() += 1.0;
156 return result;
157}
158
160 const Index lumo = homo_ + 1;
161 const Index n_occ = lumo - rpamin_;
162 const Index n_unocc = rpamax_ - lumo + 1;
163 const Index rpasize = n_occ * n_unocc;
164
165 Eigen::VectorXd AmB = Calculate_H2p_AmB();
166 Eigen::MatrixXd ApB = Calculate_H2p_ApB();
167
169 sol.ERPA_correlation = -0.25 * (ApB.trace() + AmB.sum());
170
171 // C = AmB^1/2 * ApB * AmB^1/2
172 Eigen::MatrixXd& C = ApB;
173 C.applyOnTheLeft(AmB.cwiseSqrt().asDiagonal());
174 C.applyOnTheRight(AmB.cwiseSqrt().asDiagonal());
175
176 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es = Diagonalize_H2p_C(C);
177
178 // Do not remove this line! It has to be there for MKL to not crash
179 sol.omega = Eigen::VectorXd::Zero(es.eigenvalues().size());
180 sol.omega = es.eigenvalues().cwiseSqrt();
181 sol.ERPA_correlation += 0.5 * sol.omega.sum();
182
183 {
184 std::ostringstream oss;
185 oss << TimeStamp() << " RKS H2p poles (first 20, eV): ";
186 const Index nprint = std::min<Index>(20, sol.omega.size());
187 for (Index i = 0; i < nprint; ++i) {
188 oss << std::fixed << std::setprecision(6)
189 << tools::conv::hrt2ev * sol.omega(i);
190 if (i + 1 < nprint) {
191 oss << ", ";
192 }
193 }
194 XTP_LOG(Log::error, log_) << oss.str() << std::flush;
195 }
196
198 << " Lowest neutral excitation energy (eV): "
199 << tools::conv::hrt2ev * sol.omega.minCoeff()
200 << std::flush;
201
202 // RPA correlation energy calculated from Eq.9 of J. Chem. Phys. 132, 234114
203 // (2010)
205 << TimeStamp()
206 << " RPA correlation energy (Hartree): " << sol.ERPA_correlation
207 << std::flush;
208
209 sol.XpY = Eigen::MatrixXd(rpasize, rpasize);
210
211 Eigen::VectorXd AmB_sqrt = AmB.cwiseSqrt();
212 Eigen::VectorXd Omega_sqrt_inv = sol.omega.cwiseSqrt().cwiseInverse();
213 for (int s = 0; s < rpasize; s++) {
214 sol.XpY.col(s) =
215 Omega_sqrt_inv(s) * AmB_sqrt.cwiseProduct(es.eigenvectors().col(s));
216 }
217
218 return sol;
219}
220
221Eigen::VectorXd RPA::Calculate_H2p_AmB() const {
222 const Index lumo = homo_ + 1;
223 const Index n_occ = lumo - rpamin_;
224 const Index n_unocc = rpamax_ - lumo + 1;
225 const Index rpasize = n_occ * n_unocc;
226 vc2index vc = vc2index(0, 0, n_unocc);
227 Eigen::VectorXd AmB = Eigen::VectorXd::Zero(rpasize);
228 for (Index v = 0; v < n_occ; v++) {
229 Index i = vc.I(v, 0);
230 AmB.segment(i, n_unocc) =
231 energies_.segment(n_occ, n_unocc).array() - energies_(v);
232 }
233 return AmB;
234}
235
236Eigen::MatrixXd RPA::Calculate_H2p_ApB() const {
237 const Index lumo = homo_ + 1;
238 const Index n_occ = lumo - rpamin_;
239 const Index n_unocc = rpamax_ - lumo + 1;
240 const Index rpasize = n_occ * n_unocc;
241 vc2index vc = vc2index(0, 0, n_unocc);
242 Eigen::MatrixXd ApB = Eigen::MatrixXd::Zero(rpasize, rpasize);
243#pragma omp parallel for schedule(guided)
244 for (Index v2 = 0; v2 < n_occ; v2++) {
245 Index i2 = vc.I(v2, 0);
246 const Eigen::MatrixXd Mmn_v2T =
247 Mmn_[v2].middleRows(n_occ, n_unocc).transpose();
248 for (Index v1 = v2; v1 < n_occ; v1++) {
249 Index i1 = vc.I(v1, 0);
250 // Multiply with factor 2 to sum over both (identical) spin states
251 ApB.block(i1, i2, n_unocc, n_unocc) =
252 2 * 2 * Mmn_[v1].middleRows(n_occ, n_unocc) * Mmn_v2T;
253 }
254 }
255 ApB.diagonal() += Calculate_H2p_AmB();
256 return ApB;
257}
258
259Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> RPA::Diagonalize_H2p_C(
260 const Eigen::MatrixXd& C) const {
262 << TimeStamp() << " Diagonalizing two-particle Hamiltonian "
263 << std::flush;
264 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es(C); // Uses lower triangle
266 << TimeStamp() << " Diagonalization done " << std::flush;
267 double minCoeff = es.eigenvalues().minCoeff();
268 if (minCoeff <= 0.0) {
270 << TimeStamp() << " Detected non-positive eigenvalue: " << minCoeff
271 << std::flush;
272 throw std::runtime_error("Detected non-positive eigenvalue.");
273 }
274 return es;
275}
276
277} // namespace xtp
278} // namespace votca
void A_TDA(const Eigen::VectorXd &vec, Index OpenmpThread)
Eigen::MatrixXd getReductionVar()
void createTemporaries(Index rows, Index cols)
void PushMatrix(const Eigen::MatrixXd &mat, Index OpenmpThread)
Eigen::MatrixXd Calculate_H2p_ApB() const
Definition rpa.cc:236
void ShiftUncorrectedEnergies(const Eigen::VectorXd &dftenergies, Index qpmin, Index gwsize)
Definition rpa.cc:50
Index rpamax_
Definition rpa.h:80
Eigen::VectorXd energies_
Definition rpa.h:83
Eigen::MatrixXd calculate_epsilon(double frequency) const
Definition rpa.cc:76
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:221
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > Diagonalize_H2p_C(const Eigen::MatrixXd &C) const
Definition rpa.cc:259
Logger & log_
Definition rpa.h:85
rpa_eigensolution Diagonalize_H2p() const
Definition rpa.cc:159
void UpdateRPAInputEnergies(const Eigen::VectorXd &dftenergies, const Eigen::VectorXd &gwaenergies, Index qpmin)
Definition rpa.cc:32
double getMaxCorrection(const Eigen::VectorXd &dftenergies, Index min, Index max) const
Definition rpa.cc:65
Index homo_
Definition rpa.h:78
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 I(Index v, Index c) const
Definition vc2index.h:42
#define XTP_LOG(level, log)
Definition logger.h:40
const double hrt2ev
Definition constants.h:53
Index getThreadId()
Definition eigen.h:143
Charge transport classes.
Definition ERIs.h:28
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26