votca 2026-dev
Loading...
Searching...
No Matches
rpa_uks.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_UKS_H
22#define VOTCA_XTP_RPA_UKS_H
23
24#include <complex>
25#include <vector>
26
27#include "eigen.h"
28#include "logger.h"
30
31namespace votca {
32namespace xtp {
33
68class RPA_UKS {
69 public:
77 RPA_UKS(Logger& log, const TCMatrix_gwbse_spin& Mmn) : log_(log), Mmn_(Mmn) {}
78
91 void configure(Index homo_alpha, Index homo_beta, Index rpamin,
92 Index rpamax) {
93 homo_alpha_ = homo_alpha;
94 homo_beta_ = homo_beta;
95 rpamin_ = rpamin;
96 rpamax_ = rpamax;
98 }
99
103 double getEta() const { return eta_; }
104
111 Eigen::MatrixXd calculate_epsilon_i(double frequency) const {
112 return calculate_epsilon<true>(frequency);
113 }
114
121 Eigen::MatrixXd calculate_epsilon_r(double frequency) const {
122 return calculate_epsilon<false>(frequency);
123 }
124
131 Eigen::MatrixXd calculate_epsilon_r(std::complex<double> frequency) const;
132
140 const Eigen::VectorXd& getRPAInputEnergiesAlpha() const {
141 return energies_alpha_;
142 }
143
147 const Eigen::VectorXd& getRPAInputEnergiesBeta() const {
148 return energies_beta_;
149 }
150
158 void setRPAInputEnergies(const Eigen::VectorXd& rpaenergies_alpha,
159 const Eigen::VectorXd& rpaenergies_beta) {
160 energies_alpha_ = rpaenergies_alpha;
161 energies_beta_ = rpaenergies_beta;
163 }
164
173 void UpdateRPAInputEnergies(const Eigen::VectorXd& dftenergies_alpha,
174 const Eigen::VectorXd& dftenergies_beta,
175 const Eigen::VectorXd& gwaenergies_alpha,
176 const Eigen::VectorXd& gwaenergies_beta,
177 Index qpmin);
178
230
243 Eigen::VectorXd omega;
244 Eigen::MatrixXd XpY;
246 };
247
256 const rpa_eigensolution& Diagonalize_H2p() const;
257
258 Index getHomoAlpha() const { return homo_alpha_; }
259 Index getHomoBeta() const { return homo_beta_; }
260
262 const Eigen::VectorXd*& omegas,
263 const std::vector<Eigen::VectorXd>*& modes) const;
264
265 private:
275 template <bool imag>
276 Eigen::MatrixXd calculate_epsilon(double frequency) const;
277
278 void InvalidateH2pCache();
279
291 Eigen::VectorXd Calculate_H2p_AmB() const;
292
305 Eigen::MatrixXd Calculate_H2p_ApB() const;
306
313 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> Diagonalize_H2p_C(
314 const Eigen::MatrixXd& C) const;
315
324 void ShiftUncorrectedEnergies(Eigen::VectorXd& energies,
325 const Eigen::VectorXd& dftenergies, Index homo,
326 Index qpmin, Index gwsize);
327
331 double getMaxCorrection(const Eigen::VectorXd& energies,
332 const Eigen::VectorXd& dftenergies, Index min,
333 Index max) const;
334
335 // Spin-resolved HOMO indices in the full orbital numbering
338
339 // Shared orbital window used to build the RPA quantities
342
343 // Small broadening parameter used in real-frequency expressions
344 const double eta_ = 0.0001;
345
346 // Current alpha/beta energies used in the response over [rpamin_, rpamax_]
347 Eigen::VectorXd energies_alpha_;
348 Eigen::VectorXd energies_beta_;
349
351
352 // Spin-resolved three-center integrals:
353 // Mmn_.alpha for alpha MOs, Mmn_.beta for beta MOs
355
356 mutable bool h2p_cached_ = false;
358
359 // Cached screening modes (shared across alpha/beta exact GW)
360 mutable bool screening_cached_ = false;
361 mutable Eigen::VectorXd screening_omegas_;
362 mutable std::vector<Eigen::VectorXd> screening_modes_;
363 void BuildCachedScreeningModes() const;
364};
365
366} // namespace xtp
367} // namespace votca
368
369#endif
Logger is used for thread-safe output of messages.
Definition logger.h:164
void BuildCachedScreeningModes() const
Definition rpa_uks.cc:92
Eigen::VectorXd screening_omegas_
Definition rpa_uks.h:361
Eigen::VectorXd energies_beta_
Definition rpa_uks.h:348
const Eigen::VectorXd & getRPAInputEnergiesBeta() const
Access the current beta RPA input energies used in the response.
Definition rpa_uks.h:147
Eigen::MatrixXd calculate_epsilon_i(double frequency) const
Dielectric matrix on the imaginary frequency axis.
Definition rpa_uks.h:111
const TCMatrix_gwbse_spin & Mmn_
Definition rpa_uks.h:354
Eigen::MatrixXd calculate_epsilon(double frequency) const
Internal implementation of epsilon(w).
Definition rpa_uks.cc:204
Index getHomoBeta() const
Definition rpa_uks.h:259
Index getHomoAlpha() const
Definition rpa_uks.h:258
const double eta_
Definition rpa_uks.h:344
rpa_eigensolution h2p_solution_cache_
Definition rpa_uks.h:357
void InvalidateH2pCache()
Definition rpa_uks.cc:73
double getMaxCorrection(const Eigen::VectorXd &energies, const Eigen::VectorXd &dftenergies, Index min, Index max) const
Return the largest absolute GW correction in a given index range.
Definition rpa_uks.cc:186
std::vector< Eigen::VectorXd > screening_modes_
Definition rpa_uks.h:362
void GetCachedScreeningModes(const Eigen::VectorXd *&omegas, const std::vector< Eigen::VectorXd > *&modes) const
Definition rpa_uks.cc:82
Eigen::VectorXd Calculate_H2p_AmB() const
Construct the diagonal AmB block of the two-particle Hamiltonian.
Definition rpa_uks.cc:425
Eigen::VectorXd energies_alpha_
Definition rpa_uks.h:347
double getEta() const
Small positive broadening used in the real-frequency response.
Definition rpa_uks.h:103
const Eigen::VectorXd & getRPAInputEnergiesAlpha() const
Access the current alpha RPA input energies used in the response.
Definition rpa_uks.h:140
void ShiftUncorrectedEnergies(Eigen::VectorXd &energies, const Eigen::VectorXd &dftenergies, Index homo, Index qpmin, Index gwsize)
Shift uncorrected states outside the GW window.
Definition rpa_uks.cc:163
void setRPAInputEnergies(const Eigen::VectorXd &rpaenergies_alpha, const Eigen::VectorXd &rpaenergies_beta)
Set alpha and beta RPA input energies directly.
Definition rpa_uks.h:158
Eigen::MatrixXd Calculate_H2p_ApB() const
Construct the ApB block of the unrestricted two-particle Hamiltonian.
Definition rpa_uks.cc:464
Eigen::SelfAdjointEigenSolver< Eigen::MatrixXd > Diagonalize_H2p_C(const Eigen::MatrixXd &C) const
Diagonalize the symmetrized two-particle matrix C.
Definition rpa_uks.cc:560
RPA_UKS(Logger &log, const TCMatrix_gwbse_spin &Mmn)
Definition rpa_uks.h:77
const rpa_eigensolution & Diagonalize_H2p() const
Diagonalize the unrestricted two-particle Hamiltonian.
Definition rpa_uks.cc:369
void configure(Index homo_alpha, Index homo_beta, Index rpamin, Index rpamax)
Configure orbital window and spin-resolved HOMO indices.
Definition rpa_uks.h:91
Eigen::MatrixXd calculate_epsilon_r(double frequency) const
Dielectric matrix on the real frequency axis.
Definition rpa_uks.h:121
void UpdateRPAInputEnergies(const Eigen::VectorXd &dftenergies_alpha, const Eigen::VectorXd &dftenergies_beta, const Eigen::VectorXd &gwaenergies_alpha, const Eigen::VectorXd &gwaenergies_beta, Index qpmin)
Update alpha and beta RPA input energies from DFT and GW energies.
Definition rpa_uks.cc:41
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
Solution of the two-particle Hamiltonian in the unrestricted basis.
Definition rpa_uks.h:242