votca 2026-dev
Loading...
Searching...
No Matches
ImaginaryAxisIntegration.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_IMAGINARYAXISINTEGRATION_H
21#define VOTCA_XTP_IMAGINARYAXISINTEGRATION_H
22
23#include "eigen.h"
24#include "quadrature_factory.h"
25#include "rpa.h"
26#include "rpa_uks.h"
27#include <memory>
28
29// Computes the contribution from the Gauss-Laguerre quadrature to the
30// self-energy expectation matrix for given RPA and frequencies
31namespace votca {
32namespace xtp {
33
35
36 public:
47
48 ImaginaryAxisIntegration(const Eigen::VectorXd& energies,
49 const TCMatrix_gwbse& Mmn);
50
51 void configure(options opt, const RPA& rpa,
52 const Eigen::MatrixXd& kDielMxInv_zero);
53 void configure(options opt, const RPA_UKS& rpa,
54 const Eigen::MatrixXd& kDielMxInv_zero);
55
56 double SigmaGQDiag(double frequency, Index gw_level, double eta) const;
57
58 private:
60
61 std::unique_ptr<GaussianQuadratureBase> gq_ = nullptr;
62
63 // This function calculates and stores inverses of the microscopic dielectric
64 // matrix in a matrix vector
65 void CalcDielInvVector(const RPA& rpa,
66 const Eigen::MatrixXd& kDielMxInv_zero);
67 void CalcDielInvVector(const RPA_UKS& rpa,
68 const Eigen::MatrixXd& kDielMxInv_zero);
69 const Eigen::VectorXd& energies_;
70 std::vector<Eigen::MatrixXd> dielinv_matrices_r_;
72};
73} // namespace xtp
74} // namespace votca
75#endif // VOTCA_XTP_IMAGINARYAXISINTEGRATION_H
double SigmaGQDiag(double frequency, Index gw_level, double eta) const
void CalcDielInvVector(const RPA &rpa, const Eigen::MatrixXd &kDielMxInv_zero)
std::vector< Eigen::MatrixXd > dielinv_matrices_r_
ImaginaryAxisIntegration(const Eigen::VectorXd &energies, const TCMatrix_gwbse &Mmn)
void configure(options opt, const RPA &rpa, const Eigen::MatrixXd &kDielMxInv_zero)
std::unique_ptr< GaussianQuadratureBase > gq_
Unrestricted RPA helper for spin-resolved GW screening.
Definition rpa_uks.h:68
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