votca 2024.1-dev
Loading...
Searching...
No Matches
sigma_base.cc
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// Standard includes
21#include <cmath>
22
23// Third party includes
24#include <boost/math/constants/constants.hpp>
25
26// VOTCA includes
28
29// Local VOTCA includes
32
33namespace votca {
34namespace xtp {
35
36Eigen::MatrixXd Sigma_base::CalcExchangeMatrix() const {
37 Eigen::MatrixXd result = Eigen::MatrixXd::Zero(qptotal_, qptotal_);
38 Index occlevel = opt_.homo - opt_.rpamin + 1;
39 Index qpmin = opt_.qpmin - opt_.rpamin;
40#pragma omp parallel for schedule(dynamic)
41 for (Index gw_level1 = 0; gw_level1 < qptotal_; gw_level1++) {
42 const Eigen::MatrixXd& Mmn1 = Mmn_[gw_level1 + qpmin];
43 for (Index gw_level2 = gw_level1; gw_level2 < qptotal_; gw_level2++) {
44 const Eigen::MatrixXd& Mmn2 = Mmn_[gw_level2 + qpmin];
45 double sigma_x =
46 -(Mmn1.topRows(occlevel).cwiseProduct(Mmn2.topRows(occlevel))).sum();
47 result(gw_level2, gw_level1) = sigma_x;
48 }
49 }
50 result = result.selfadjointView<Eigen::Lower>();
51 return result;
52}
53
55 const Eigen::VectorXd& frequencies) const {
56 Eigen::VectorXd result = Eigen::VectorXd::Zero(qptotal_);
57#pragma omp parallel for schedule(dynamic)
58 for (Index gw_level = 0; gw_level < qptotal_; gw_level++) {
59 result(gw_level) =
60 CalcCorrelationDiagElement(gw_level, frequencies[gw_level]);
61 }
62 return result;
63}
64
66 const Eigen::VectorXd& frequencies) const {
67 Eigen::MatrixXd result = Eigen::MatrixXd::Zero(qptotal_, qptotal_);
68#pragma omp parallel for schedule(dynamic)
69 for (Index gw_level1 = 0; gw_level1 < qptotal_; gw_level1++) {
70 for (Index gw_level2 = gw_level1 + 1; gw_level2 < qptotal_; gw_level2++) {
71 double sigma_c = CalcCorrelationOffDiagElement(
72 gw_level1, gw_level2, frequencies[gw_level1], frequencies[gw_level2]);
73 result(gw_level2, gw_level1) = sigma_c;
74 }
75 }
76 result = result.selfadjointView<Eigen::Lower>();
77 return result;
78}
79
80} // namespace xtp
81} // namespace votca
Eigen::MatrixXd CalcCorrelationOffDiag(const Eigen::VectorXd &frequencies) const
Definition sigma_base.cc:65
virtual double CalcCorrelationOffDiagElement(Index gw_level1, Index gw_level2, double frequency1, double frequency2) const =0
virtual double CalcCorrelationDiagElement(Index gw_level, double frequency) const =0
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