votca 2024.1-dev
Loading...
Searching...
No Matches
sigma_base.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_BASE_H
22#define VOTCA_XTP_SIGMA_BASE_H
23
24// Local VOTCA includes
25#include "eigen.h"
26
27namespace votca {
28namespace xtp {
29
30class TCMatrix_gwbse;
31class RPA;
32
34 public:
35 Sigma_base(TCMatrix_gwbse& Mmn, const RPA& rpa) : Mmn_(Mmn), rpa_(rpa) {};
36
37 virtual ~Sigma_base() = default;
38
39 struct options {
45 double eta;
46 std::string quadrature_scheme; // Gaussian-quadrature scheme to use in CDA
47 Index order; // used in numerical integration of CDA Sigma
48 double alpha;
49 };
50
51 void configure(options opt) {
52 opt_ = opt;
53 qptotal_ = opt.qpmax - opt.qpmin + 1;
54 rpatotal_ = opt.rpamax - opt.rpamin + 1;
55 }
56
57 // Calculates full exchange matrix
58 Eigen::MatrixXd CalcExchangeMatrix() const;
59 // Calculates correlation diagonal
60 Eigen::VectorXd CalcCorrelationDiag(const Eigen::VectorXd& frequencies) const;
61 // Calculates correlation off-diagonal
62 Eigen::MatrixXd CalcCorrelationOffDiag(
63 const Eigen::VectorXd& frequencies) const;
64
65 // Sets up the screening parametrisation
66 virtual void PrepareScreening() = 0;
67 // Calculates Sigma_c diagonal elements
69 Index gw_level, double frequency) const = 0;
70 virtual double CalcCorrelationDiagElement(Index gw_level,
71 double frequency) const = 0;
72 // Calculates Sigma_c off-diagonal elements
73 virtual double CalcCorrelationOffDiagElement(Index gw_level1, Index gw_level2,
74 double frequency1,
75 double frequency2) const = 0;
76
77 protected:
80 const RPA& rpa_;
81
84};
85} // namespace xtp
86} // namespace votca
87
88#endif // VOTCA_XTP_SIGMA_BASE_H
Sigma_base(TCMatrix_gwbse &Mmn, const RPA &rpa)
Definition sigma_base.h:35
Eigen::MatrixXd CalcCorrelationOffDiag(const Eigen::VectorXd &frequencies) const
Definition sigma_base.cc:65
virtual double CalcCorrelationDiagElementDerivative(Index gw_level, double frequency) const =0
virtual double CalcCorrelationOffDiagElement(Index gw_level1, Index gw_level2, double frequency1, double frequency2) const =0
void configure(options opt)
Definition sigma_base.h:51
virtual double CalcCorrelationDiagElement(Index gw_level, double frequency) const =0
virtual void PrepareScreening()=0
virtual ~Sigma_base()=default
Eigen::VectorXd CalcCorrelationDiag(const Eigen::VectorXd &frequencies) const
Definition sigma_base.cc:54
TCMatrix_gwbse & Mmn_
Definition sigma_base.h:79
Eigen::MatrixXd CalcExchangeMatrix() const
Definition sigma_base.cc:36
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26