votca 2026-dev
Loading...
Searching...
No Matches
sigma_ppm.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// VOTCA includes
22#include <votca/tools/globals.h>
23
24// Local VOTCA includes
25#include "sigma_ppm.h"
26#include "votca/xtp/ppm.h"
28
29namespace votca {
30namespace xtp {
31
33 ppm_.PPM_construct_parameters(rpa_);
34 Mmn_.MultiplyRightWithAuxMatrix(ppm_.getPpm_phi());
35}
36
38 double frequency) const {
40 const Index lumo = opt_.homo + 1;
41 const double eta2 = opt_.eta * opt_.eta;
42 const Index levelsum = Mmn_.nsize(); // total number of bands
43 const Index qpmin_offset = opt_.qpmin - opt_.rpamin;
44 double sigma = 0.0;
45 for (Index i_aux = 0; i_aux < Mmn_.auxsize(); i_aux++) {
46 // the ppm_weights smaller 1.e-5 are set to zero in rpa.cc
47 // PPM_construct_parameters
48 if (ppm_.getPpm_weight()(i_aux) < 1.e-9) {
49 continue;
50 }
51 const double ppm_freq = ppm_.getPpm_freq()(i_aux);
52 const double fac = 0.5 * ppm_.getPpm_weight()(i_aux) * ppm_freq;
53 const Eigen::ArrayXd Mmn2 =
54 Mmn_[gw_level + qpmin_offset].col(i_aux).cwiseAbs2();
55 Eigen::ArrayXd temp = frequency - rpa_.getRPAInputEnergies().array();
56 temp.segment(0, lumo) += ppm_freq;
57 temp.segment(lumo, levelsum - lumo) -= ppm_freq;
58 Eigen::ArrayXd denom = temp.abs2() + eta2;
59 sigma += fac * (Mmn2 * temp / denom).sum();
60 }
61 return sigma;
62}
63
65 double frequency) const {
66 const Index lumo = opt_.homo + 1;
67 const double eta2 = opt_.eta * opt_.eta;
68 const Index levelsum = Mmn_.nsize(); // total number of bands
69 const Index qpmin_offset = opt_.qpmin - opt_.rpamin;
70 double dsigma_domega = 0.0;
71 for (Index i_aux = 0; i_aux < Mmn_.auxsize(); i_aux++) {
72 // the ppm_weights smaller 1.e-5 are set to zero in rpa.cc
73 // PPM_construct_parameters
74 if (ppm_.getPpm_weight()(i_aux) < 1.e-9) {
75 continue;
76 }
77 const double ppm_freq = ppm_.getPpm_freq()(i_aux);
78 const double fac = 0.5 * ppm_.getPpm_weight()(i_aux) * ppm_freq;
79 const Eigen::ArrayXd Mmn2 =
80 Mmn_[gw_level + qpmin_offset].col(i_aux).cwiseAbs2();
81 Eigen::ArrayXd temp = frequency - rpa_.getRPAInputEnergies().array();
82 temp.segment(0, lumo) += ppm_freq;
83 temp.segment(lumo, levelsum - lumo) -= ppm_freq;
84 Eigen::ArrayXd denom = temp.abs2() + eta2;
85 dsigma_domega += fac * ((eta2 - temp.abs2()) * Mmn2 / denom.abs2()).sum();
86 }
87 return dsigma_domega;
88}
89
91 Index gw_level2,
92 double frequency1,
93 double frequency2) const {
94 const Index lumo = opt_.homo + 1;
95 const double eta2 = opt_.eta * opt_.eta;
96 const Index levelsum = Mmn_.nsize(); // total number of bands
97 const Index auxsize = Mmn_.auxsize(); // size of the GW basis
98 const Eigen::VectorXd ppm_weight = ppm_.getPpm_weight();
99 const Eigen::VectorXd ppm_freqs = ppm_.getPpm_freq();
100 const Index qpmin_offset = opt_.qpmin - opt_.rpamin;
101 const Eigen::VectorXd RPAEnergies = rpa_.getRPAInputEnergies();
102 double sigma_c = 0;
103 for (Index i_aux = 0; i_aux < auxsize; i_aux++) {
104 // the ppm_weights smaller 1.e-5 are set to zero in rpa.cc
105 // PPM_construct_parameters
106 if (ppm_weight(i_aux) < 1.e-9) {
107 continue;
108 }
109 const double ppm_freq = ppm_freqs(i_aux);
110 const double fac = 0.25 * ppm_weight(i_aux) * ppm_freq;
111 const Eigen::MatrixXd& Mmn1 = Mmn_[gw_level1 + qpmin_offset];
112 const Eigen::MatrixXd& Mmn2 = Mmn_[gw_level2 + qpmin_offset];
113 const Eigen::ArrayXd Mmn1xMmn2 =
114 Mmn1.col(i_aux).cwiseProduct(Mmn2.col(i_aux));
115 Eigen::ArrayXd temp1 = RPAEnergies;
116 temp1.segment(0, lumo) -= ppm_freq;
117 temp1.segment(lumo, levelsum - lumo) += ppm_freq;
118 Eigen::ArrayXd temp2 = (frequency2 - temp1);
119 temp1 = (frequency1 - temp1);
120 sigma_c +=
121 fac * ((temp1 / (temp1.abs2() + eta2) + temp2 / (temp2.abs2() + eta2)) *
122 Mmn1xMmn2)
123 .sum();
124 }
125 return sigma_c;
126}
127
128} // namespace xtp
129} // namespace votca
double CalcCorrelationOffDiagElement(Index gw_level1, Index gw_level2, double frequency1, double frequency2) const final
Definition sigma_ppm.cc:90
double CalcCorrelationDiagElement(Index gw_level, double frequency) const final
Definition sigma_ppm.cc:37
void PrepareScreening() final
Definition sigma_ppm.cc:32
double CalcCorrelationDiagElementDerivative(Index gw_level, double frequency) const final
Definition sigma_ppm.cc:64
TCMatrix_gwbse & Mmn_
Definition sigma_base.h:82
void CountDiagEval() const
Definition sigma_base.h:88
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