votca 2026-dev
Loading...
Searching...
No Matches
sigma_cda_uks.cc
Go to the documentation of this file.
1/*
2 * Copyright 2009-2026 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#include "sigma_cda_uks.h"
21#include "votca/xtp/gw.h"
22
23#include <algorithm>
24#include <iostream>
25#include <limits>
26
27namespace votca {
28namespace xtp {
29
30namespace {
31
32constexpr double kResidueFactorTol = 1e-10;
33constexpr double kCdaWarnMinDelta = 1e-4; // Ha
34constexpr double kCdaWarnKappa = 1e3;
35constexpr double kTiny = 1e-14;
36
37const char* SpinName(TCMatrix::SpinChannel spin) {
38 return (spin == TCMatrix::SpinChannel::Alpha) ? "alpha" : "beta";
39}
40
41} // namespace
42
45 opt.homo = opt_.homo;
46 opt.order = opt_.order;
47 opt.qptotal = qptotal_;
48 opt.qpmin = opt_.qpmin;
49 opt.rpamax = opt_.rpamax;
50 opt.rpamin = opt_.rpamin;
51 opt.alpha = opt_.alpha;
52 opt.quadrature_scheme = opt_.quadrature_scheme;
53
55 rpa_.calculate_epsilon_r(std::complex<double>(0.0, 0.0)).inverse();
56 kDielMxInv_zero_.diagonal().array() -= 1.0;
57 gq_.configure(opt, rpa_, kDielMxInv_zero_);
58}
59
61 const Eigen::MatrixXd::ConstRowXpr& Imx_row, double delta,
62 double eta) const {
63 std::complex<double> delta_eta(delta, eta);
64
65 Eigen::MatrixXd DielMxInv = rpa_.calculate_epsilon_r(delta_eta);
66 Eigen::VectorXd x =
67 DielMxInv.partialPivLu().solve(Imx_row.transpose()) - Imx_row.transpose();
68
69 return x.dot(Imx_row.transpose());
70}
71
72double Sigma_CDA_UKS::CalcResiduePrefactor(double e_f, double e_m,
73 double frequency) const {
74 double factor = 0.0;
75 double tolerance = 1e-10;
76 if (e_f < e_m && e_m < frequency) {
77 factor = 1.0;
78 } else if (e_f > e_m && e_m > frequency) {
79 factor = -1.0;
80 } else if (std::abs(e_m - frequency) < tolerance && e_f > e_m) {
81 factor = -0.5;
82 } else if (std::abs(e_m - frequency) < tolerance && e_f < e_m) {
83 factor = 0.5;
84 }
85 return factor;
86}
87
89 Index gw_level) const {
90 const Eigen::VectorXd& rpa_energies = getSpinRPAInputEnergies();
91 Index rpatotal = rpa_energies.size();
92 Index gw_level_offset = gw_level + opt_.qpmin - opt_.rpamin;
93
94 double sigma_c = 0.0;
95 double sigma_c_tail = 0.0;
96 Index homo = opt_.homo - opt_.rpamin;
97 Index lumo = homo + 1;
98 double fermi_rpa = (rpa_energies(lumo) + rpa_energies(homo)) / 2.0;
99
100 double min_abs_delta = std::numeric_limits<double>::infinity();
101 double sum_abs_contributions = 0.0;
102 double sum_signed_contributions = 0.0;
103
104 const Eigen::MatrixXd& Imx = Mmn_[gw_level_offset];
105
106 for (Index i = 0; i < rpatotal; ++i) {
107 double delta = rpa_energies(i) - frequency;
108 double abs_delta = std::abs(delta);
109 min_abs_delta = std::min(min_abs_delta, abs_delta);
110
111 double factor = CalcResiduePrefactor(fermi_rpa, rpa_energies(i), frequency);
112
113 if (std::abs(factor) > kResidueFactorTol) {
114 double diag = CalcDiagContribution(Imx.row(i), abs_delta, rpa_.getEta());
115 double contribution = factor * diag;
116 sigma_c += contribution;
117
118 sum_abs_contributions += std::abs(contribution);
119 sum_signed_contributions += contribution;
120 }
121
122 if (abs_delta > kResidueFactorTol) {
123 double tail =
124 CalcDiagContributionValue_tail(Imx.row(i), delta, opt_.alpha);
125 sigma_c_tail += tail;
126 }
127 }
128
129 double kappa =
130 sum_abs_contributions / (std::abs(sum_signed_contributions) + kTiny);
131
132 if (min_abs_delta < kCdaWarnMinDelta || kappa > kCdaWarnKappa) {
133 std::cout << "\nWarning: CDA may be unreliable for " << SpinName(spin_)
134 << " GW level " << (gw_level + opt_.qpmin)
135 << " at omega = " << frequency
136 << " Ha: min |e_m - omega| = " << min_abs_delta
137 << " Ha, cancellation metric = " << kappa << std::flush;
138 }
139
140 return sigma_c + sigma_c_tail;
141}
142
144 double frequency) const {
145 double sigma_c_residue = CalcResidueContribution(frequency, gw_level);
146 double sigma_c_integral = gq_.SigmaGQDiag(frequency, gw_level, rpa_.getEta());
147 return sigma_c_residue + sigma_c_integral;
148}
149
151 const Eigen::MatrixXd::ConstRowXpr& Imx_row, double delta,
152 double alpha) const {
153 double erfc_factor = 0.5 * std::copysign(1.0, delta) *
154 std::exp(std::pow(alpha * delta, 2)) *
155 std::erfc(std::abs(alpha * delta));
156
157 double value = (Imx_row * kDielMxInv_zero_).dot(Imx_row);
158 return value * erfc_factor;
159}
160
161} // namespace xtp
162} // namespace votca
double CalcCorrelationDiagElement(Index gw_level, double frequency) const final
double CalcResiduePrefactor(double e_f, double e_m, double frequency) const
double CalcResidueContribution(double frequency, Index gw_level) const
double CalcDiagContributionValue_tail(const Eigen::MatrixXd::ConstRowXpr &Imx_row, double delta, double alpha) const
ImaginaryAxisIntegration gq_
Eigen::MatrixXd kDielMxInv_zero_
double CalcDiagContribution(const Eigen::MatrixXd::ConstRowXpr &Imx_row, double delta, double eta) const
const Eigen::VectorXd & getSpinRPAInputEnergies() const
TCMatrix::SpinChannel spin_
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