votca 2024.1-dev
Loading...
Searching...
No Matches
sigma_cda.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#ifndef VOTCA_XTP_SIGMA_CDA_H
21#define VOTCA_XTP_SIGMA_CDA_H
23#include "votca/xtp/logger.h"
24#include "votca/xtp/rpa.h"
26#include <complex>
27
28// This computes the whole expectation matrix for the correlational part of the
29// self-energy with the Contour Deformation Approach according to Eqns 28 and 29
30// of JCP 152, 114103 (2020). There are two contributions:
31// - from a numerical integration using Gaussian quadratures along the imaginary
32// frequncy axis (Eq. 28)
33// - from the residues included in the contours (Eq.29)
34// Both contributions contain term from a Gaussian tail with parameter alpha.
35namespace votca {
36namespace xtp {
37
38class Sigma_CDA : public Sigma_base {
39
40 public:
42 : Sigma_base(Mmn, rpa), gq_(rpa.getRPAInputEnergies(), Mmn) {};
43
44 ~Sigma_CDA() = default;
45
46 // Prepares the zero and imaginary frequency kappa matrices with
47 // kappa(omega) = epsilon^-1(omega) - 1 needed in numerical
48 // integration and for the Gaussian tail
49 void PrepareScreening() final;
50
51 // calculates the diagonal elements of the self-energy correlation part
52 double CalcCorrelationDiagElement(Index gw_level,
53 double frequency) const final;
54
55 // numerical derivatice of the self-energy
57 double frequency) const final {
58 double h = 1e-3;
59 double plus = CalcCorrelationDiagElement(gw_level, frequency + h);
60 double minus = CalcCorrelationDiagElement(gw_level, frequency - h);
61 return (plus - minus) / (2 * h);
62 }
63 // Calculates Sigma_c off-diagonal elements
65 double) const final {
66 return 0;
67 }
68
69 private:
70 // Theta-function weight of a residue
71 double CalcResiduePrefactor(double e_f, double e_m, double frequency) const;
72
73 // Sigma_c from all possible residues for given gw_level and frequency
74 double CalcResidueContribution(double frequency, Index gw_level) const;
75
76 // Sigma_c part from a single residue for a given gw_level and frequency
77 double CalcDiagContribution(const Eigen::MatrixXd::ConstRowXpr& Imx_row,
78 double delta, double eta) const;
79
80 // Sigma_c part from Gaussian tail correction
82 const Eigen::MatrixXd::ConstRowXpr& Imx_row, double delta,
83 double alpha) const;
84
86 Eigen::MatrixXd kDielMxInv_zero_; // kappa = eps^-1 - 1 matrix
87};
88
89} // namespace xtp
90
91} // namespace votca
92
93#endif // VOTCA_XTP_SIGMA_CDA_H
Sigma_CDA(TCMatrix_gwbse &Mmn, RPA &rpa)
Definition sigma_cda.h:41
double CalcCorrelationDiagElementDerivative(Index gw_level, double frequency) const final
Definition sigma_cda.h:56
double CalcDiagContribution(const Eigen::MatrixXd::ConstRowXpr &Imx_row, double delta, double eta) const
Definition sigma_cda.cc:52
double CalcCorrelationOffDiagElement(Index, Index, double, double) const final
Definition sigma_cda.h:64
double CalcCorrelationDiagElement(Index gw_level, double frequency) const final
Definition sigma_cda.cc:121
ImaginaryAxisIntegration gq_
Definition sigma_cda.h:85
double CalcResidueContribution(double frequency, Index gw_level) const
Definition sigma_cda.cc:82
double CalcDiagContributionValue_tail(const Eigen::MatrixXd::ConstRowXpr &Imx_row, double delta, double alpha) const
Definition sigma_cda.cc:131
double CalcResiduePrefactor(double e_f, double e_m, double frequency) const
Definition sigma_cda.cc:64
void PrepareScreening() final
Definition sigma_cda.cc:30
Eigen::MatrixXd kDielMxInv_zero_
Definition sigma_cda.h:86
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26