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