votca 2026-dev
Loading...
Searching...
No Matches
sigma_cda_uks.h
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#pragma once
21#ifndef VOTCA_XTP_SIGMA_CDA_UKS_H
22#define VOTCA_XTP_SIGMA_CDA_UKS_H
23
25#include "votca/xtp/logger.h"
26#include "votca/xtp/rpa_uks.h"
28#include <complex>
29
30namespace votca {
31namespace xtp {
32
34 public:
37 : Sigma_base_UKS(Mmn, rpa, spin),
38 gq_(getSpinRPAInputEnergies(), Mmn[spin]) {}
39
40 ~Sigma_CDA_UKS() override = default;
41
42 void PrepareScreening() final;
43
44 double CalcCorrelationDiagElement(Index gw_level,
45 double frequency) const final;
46
48 double frequency) const final {
49 double h = 1e-3;
50 double plus = CalcCorrelationDiagElement(gw_level, frequency + h);
51 double minus = CalcCorrelationDiagElement(gw_level, frequency - h);
52 return (plus - minus) / (2 * h);
53 }
54
56 double) const final {
57 return 0.0;
58 }
59
60 private:
61 double CalcResiduePrefactor(double e_f, double e_m, double frequency) const;
62 double CalcResidueContribution(double frequency, Index gw_level) const;
63 double CalcDiagContribution(const Eigen::MatrixXd::ConstRowXpr& Imx_row,
64 double delta, double eta) const;
66 const Eigen::MatrixXd::ConstRowXpr& Imx_row, double delta,
67 double alpha) const;
68
70 Eigen::MatrixXd kDielMxInv_zero_;
71};
72
73} // namespace xtp
74} // namespace votca
75
76#endif
Unrestricted RPA helper for spin-resolved GW screening.
Definition rpa_uks.h:68
double CalcCorrelationDiagElement(Index gw_level, double frequency) const final
double CalcResiduePrefactor(double e_f, double e_m, double frequency) const
~Sigma_CDA_UKS() override=default
double CalcResidueContribution(double frequency, Index gw_level) const
double CalcCorrelationDiagElementDerivative(Index gw_level, double frequency) const final
double CalcDiagContributionValue_tail(const Eigen::MatrixXd::ConstRowXpr &Imx_row, double delta, double alpha) const
ImaginaryAxisIntegration gq_
Sigma_CDA_UKS(TCMatrix_gwbse_spin &Mmn, RPA_UKS &rpa, TCMatrix::SpinChannel spin)
double CalcCorrelationOffDiagElement(Index, Index, double, double) const final
Eigen::MatrixXd kDielMxInv_zero_
double CalcDiagContribution(const Eigen::MatrixXd::ConstRowXpr &Imx_row, double delta, double eta) const
const Eigen::VectorXd & getSpinRPAInputEnergies() const
Sigma_base_UKS(TCMatrix_gwbse_spin &Mmn, const RPA_UKS &rpa, TCMatrix::SpinChannel spin)
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26