votca 2024.1-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
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 {
42 const double eta2 = opt_.eta * opt_.eta;
43 const Index lumo = opt_.homo + 1;
44 const Index n_occ = lumo - opt_.rpamin;
45 const Index n_unocc = opt_.rpamax - opt_.homo;
46 double sigma = 0.0;
47 for (Index s = 0; s < rpa_omegas_.size(); s++) {
48 const double eigenvalue = rpa_omegas_(s);
49 const Eigen::ArrayXd res_12 = residues_[gw_level].col(s).cwiseAbs2();
50 Eigen::ArrayXd temp = -rpa_.getRPAInputEnergies().array() + frequency;
51 temp.segment(0, n_occ) += eigenvalue;
52 temp.segment(n_occ, n_unocc) -= eigenvalue;
53 const Eigen::ArrayXd denom = temp.abs2() + eta2;
54 sigma += (res_12 * temp / denom).sum();
55 }
56 return 2 * sigma;
57}
58
60 Index gw_level, double frequency) const {
61 const double eta2 = opt_.eta * opt_.eta;
62 const Index lumo = opt_.homo + 1;
63 const Index n_occ = lumo - opt_.rpamin;
64 const Index n_unocc = opt_.rpamax - opt_.homo;
65 double dsigma_domega = 0.0;
66 for (Index s = 0; s < rpa_omegas_.size(); s++) {
67 const double eigenvalue = rpa_omegas_(s);
68 const Eigen::ArrayXd res_12 = residues_[gw_level].col(s).cwiseAbs2();
69 Eigen::ArrayXd temp = -rpa_.getRPAInputEnergies().array() + frequency;
70 temp.segment(0, n_occ) += eigenvalue;
71 temp.segment(n_occ, n_unocc) -= eigenvalue;
72 const Eigen::ArrayXd denom = temp.abs2() + eta2;
73 dsigma_domega += ((eta2 - temp.abs2()) * res_12 / denom.abs2()).sum();
74 }
75 return 2 * dsigma_domega;
76}
77
79 Index gw_level2,
80 double frequency1,
81 double frequency2) const {
82 const double eta2 = opt_.eta * opt_.eta;
83 const Index lumo = opt_.homo + 1;
84 const Index n_occ = lumo - opt_.rpamin;
85 const Index n_unocc = opt_.rpamax - opt_.homo;
86 const Index rpasize = rpa_omegas_.size();
87 double sigma_c = 0.0;
88 for (Index s = 0; s < rpasize; s++) {
89 const double eigenvalue = rpa_omegas_(s);
90 const Eigen::VectorXd& res1 = residues_[gw_level1].col(s);
91 const Eigen::VectorXd& res2 = residues_[gw_level2].col(s);
92 const Eigen::VectorXd res_12 = res1.cwiseProduct(res2);
93 Eigen::ArrayXd temp1 = -rpa_.getRPAInputEnergies().array();
94 temp1.segment(0, n_occ) += eigenvalue;
95 temp1.segment(n_occ, n_unocc) -= eigenvalue;
96 const Eigen::ArrayXd temp2 = temp1 + frequency2;
97 temp1 += frequency1;
98 const Eigen::ArrayXd numer1 = res_12.array() * temp1;
99 const Eigen::ArrayXd numer2 = res_12.array() * temp2;
100 const Eigen::ArrayXd denom1 = temp1.abs2() + eta2;
101 const Eigen::ArrayXd denom2 = temp2.abs2() + eta2;
102 sigma_c += 0.5 * ((numer1 / denom1) + (numer2 / denom2)).sum();
103 }
104 // Multiply with factor 2.0 to sum over both (identical) spin states
105 return 2.0 * sigma_c;
106}
107
108Eigen::MatrixXd Sigma_Exact::CalcResidues(Index gw_level,
109 const Eigen::MatrixXd& XpY) const {
110 const Index lumo = opt_.homo + 1;
111 const Index n_occ = lumo - opt_.rpamin;
112 const Index n_unocc = opt_.rpamax - opt_.homo;
113 const Index rpasize = n_occ * n_unocc;
114 const Index qpoffset = opt_.qpmin - opt_.rpamin;
115 vc2index vc = vc2index(0, 0, n_unocc);
116 const Eigen::MatrixXd& Mmn_i = Mmn_[gw_level + qpoffset];
117 Eigen::MatrixXd res = Eigen::MatrixXd::Zero(rpatotal_, rpasize);
118 for (Index v = 0; v < n_occ; v++) { // Sum over v
119 auto Mmn_v = Mmn_[v].middleRows(n_occ, n_unocc);
120 auto fc = Mmn_v * Mmn_i.transpose(); // Sum over chi
121 auto XpY_v = XpY.middleRows(vc.I(v, 0), n_unocc);
122 res += fc.transpose() * XpY_v; // Sum over c
123 }
124 return res;
125}
126
127} // namespace xtp
128} // namespace votca
const Eigen::VectorXd & getRPAInputEnergies() const
Definition rpa.h:57
rpa_eigensolution Diagonalize_H2p() const
Definition rpa.cc:157
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:79
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
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26