votca 2026-dev
Loading...
Searching...
No Matches
sigma_exact.cc
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// Local VOTCA includes
21#include "sigma_exact.h"
22#include "votca/xtp/rpa.h"
24#include "votca/xtp/vc2index.h"
25
26namespace votca {
27namespace xtp {
28
30 RPA::rpa_eigensolution rpa_solution = rpa_.Diagonalize_H2p();
31 rpa_omegas_ = rpa_solution.omega;
32 residues_ = std::vector<Eigen::MatrixXd>(qptotal_);
33#pragma omp parallel for schedule(dynamic)
34 for (Index gw_level = 0; gw_level < qptotal_; gw_level++) {
35 residues_[gw_level] = CalcResidues(gw_level, rpa_solution.XpY);
36 }
37 return;
38}
39
41 double frequency) const {
43 const double eta2 = opt_.eta * opt_.eta;
44 const Index lumo = opt_.homo + 1;
45 const Index n_occ = lumo - opt_.rpamin;
46 const Index n_unocc = opt_.rpamax - opt_.homo;
47 double sigma = 0.0;
48 for (Index s = 0; s < rpa_omegas_.size(); s++) {
49 const double eigenvalue = rpa_omegas_(s);
50 const Eigen::ArrayXd res_12 = residues_[gw_level].col(s).cwiseAbs2();
51 Eigen::ArrayXd temp = -rpa_.getRPAInputEnergies().array() + frequency;
52 temp.segment(0, n_occ) += eigenvalue;
53 temp.segment(n_occ, n_unocc) -= eigenvalue;
54 const Eigen::ArrayXd denom = temp.abs2() + eta2;
55 sigma += (res_12 * temp / denom).sum();
56 }
57 return 2 * sigma;
58}
59
61 Index gw_level, double frequency) const {
62 const double eta2 = opt_.eta * opt_.eta;
63 const Index lumo = opt_.homo + 1;
64 const Index n_occ = lumo - opt_.rpamin;
65 const Index n_unocc = opt_.rpamax - opt_.homo;
66 double dsigma_domega = 0.0;
67 for (Index s = 0; s < rpa_omegas_.size(); s++) {
68 const double eigenvalue = rpa_omegas_(s);
69 const Eigen::ArrayXd res_12 = residues_[gw_level].col(s).cwiseAbs2();
70 Eigen::ArrayXd temp = -rpa_.getRPAInputEnergies().array() + frequency;
71 temp.segment(0, n_occ) += eigenvalue;
72 temp.segment(n_occ, n_unocc) -= eigenvalue;
73 const Eigen::ArrayXd denom = temp.abs2() + eta2;
74 dsigma_domega += ((eta2 - temp.abs2()) * res_12 / denom.abs2()).sum();
75 }
76 return 2 * dsigma_domega;
77}
78
80 Index gw_level2,
81 double frequency1,
82 double frequency2) const {
83 const double eta2 = opt_.eta * opt_.eta;
84 const Index lumo = opt_.homo + 1;
85 const Index n_occ = lumo - opt_.rpamin;
86 const Index n_unocc = opt_.rpamax - opt_.homo;
87 const Index rpasize = rpa_omegas_.size();
88 double sigma_c = 0.0;
89 for (Index s = 0; s < rpasize; s++) {
90 const double eigenvalue = rpa_omegas_(s);
91 const Eigen::VectorXd& res1 = residues_[gw_level1].col(s);
92 const Eigen::VectorXd& res2 = residues_[gw_level2].col(s);
93 const Eigen::VectorXd res_12 = res1.cwiseProduct(res2);
94 Eigen::ArrayXd temp1 = -rpa_.getRPAInputEnergies().array();
95 temp1.segment(0, n_occ) += eigenvalue;
96 temp1.segment(n_occ, n_unocc) -= eigenvalue;
97 const Eigen::ArrayXd temp2 = temp1 + frequency2;
98 temp1 += frequency1;
99 const Eigen::ArrayXd numer1 = res_12.array() * temp1;
100 const Eigen::ArrayXd numer2 = res_12.array() * temp2;
101 const Eigen::ArrayXd denom1 = temp1.abs2() + eta2;
102 const Eigen::ArrayXd denom2 = temp2.abs2() + eta2;
103 sigma_c += 0.5 * ((numer1 / denom1) + (numer2 / denom2)).sum();
104 }
105 // Multiply with factor 2.0 to sum over both (identical) spin states
106 return 2.0 * sigma_c;
107}
108
109Eigen::MatrixXd Sigma_Exact::CalcResidues(Index gw_level,
110 const Eigen::MatrixXd& XpY) const {
111 const Index lumo = opt_.homo + 1;
112 const Index n_occ = lumo - opt_.rpamin;
113 const Index n_unocc = opt_.rpamax - opt_.homo;
114 const Index rpasize = n_occ * n_unocc;
115 const Index qpoffset = opt_.qpmin - opt_.rpamin;
116 vc2index vc = vc2index(0, 0, n_unocc);
117 const Eigen::MatrixXd& Mmn_i = Mmn_[gw_level + qpoffset];
118 Eigen::MatrixXd res = Eigen::MatrixXd::Zero(rpatotal_, rpasize);
119 for (Index v = 0; v < n_occ; v++) { // Sum over v
120 auto Mmn_v = Mmn_[v].middleRows(n_occ, n_unocc);
121 auto fc = Mmn_v * Mmn_i.transpose(); // Sum over chi
122 auto XpY_v = XpY.middleRows(vc.I(v, 0), n_unocc);
123 res += fc.transpose() * XpY_v; // Sum over c
124 }
125 return res;
126}
127
128} // namespace xtp
129} // namespace votca
Eigen::MatrixXd CalcResidues(Index gw_level, const Eigen::MatrixXd &XpY) const
std::vector< Eigen::MatrixXd > residues_
Definition sigma_exact.h:53
double CalcCorrelationDiagElementDerivative(Index gw_level, double frequency) const final
void PrepareScreening() final
double CalcCorrelationDiagElement(Index gw_level, double frequency) const final
Eigen::VectorXd rpa_omegas_
Definition sigma_exact.h:52
double CalcCorrelationOffDiagElement(Index gw_level1, Index gw_level2, double frequency1, double frequency2) const final
TCMatrix_gwbse & Mmn_
Definition sigma_base.h:82
void CountDiagEval() const
Definition sigma_base.h:88
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
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