votca 2026-dev
Loading...
Searching...
No Matches
sigma_base_uks.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
21
22namespace votca {
23namespace xtp {
24
25Eigen::MatrixXd Sigma_base_UKS::CalcExchangeMatrix() const {
26 Eigen::MatrixXd result = Eigen::MatrixXd::Zero(qptotal_, qptotal_);
27 Index occlevel = opt_.homo - opt_.rpamin + 1;
28 Index qpmin = opt_.qpmin - opt_.rpamin;
29
30#pragma omp parallel for schedule(dynamic)
31 for (Index gw_level1 = 0; gw_level1 < qptotal_; gw_level1++) {
32 const Eigen::MatrixXd& Mmn1 = Mmn_[gw_level1 + qpmin];
33 for (Index gw_level2 = gw_level1; gw_level2 < qptotal_; gw_level2++) {
34 const Eigen::MatrixXd& Mmn2 = Mmn_[gw_level2 + qpmin];
35 double sigma_x =
36 -(Mmn1.topRows(occlevel).cwiseProduct(Mmn2.topRows(occlevel))).sum();
37 result(gw_level2, gw_level1) = sigma_x;
38 }
39 }
40
41 result = result.selfadjointView<Eigen::Lower>();
42 return result;
43}
44
46 const Eigen::VectorXd& frequencies) const {
47 Eigen::VectorXd result = Eigen::VectorXd::Zero(qptotal_);
48#pragma omp parallel for schedule(dynamic)
49 for (Index gw_level = 0; gw_level < qptotal_; gw_level++) {
50 result(gw_level) =
51 CalcCorrelationDiagElement(gw_level, frequencies[gw_level]);
52 }
53 return result;
54}
55
57 const Eigen::VectorXd& frequencies) const {
58 Eigen::MatrixXd result = Eigen::MatrixXd::Zero(qptotal_, qptotal_);
59#pragma omp parallel for schedule(dynamic)
60 for (Index gw_level1 = 0; gw_level1 < qptotal_; gw_level1++) {
61 for (Index gw_level2 = gw_level1 + 1; gw_level2 < qptotal_; gw_level2++) {
62 double sigma_c = CalcCorrelationOffDiagElement(
63 gw_level1, gw_level2, frequencies[gw_level1], frequencies[gw_level2]);
64 result(gw_level2, gw_level1) = sigma_c;
65 }
66 }
67 result = result.selfadjointView<Eigen::Lower>();
68 return result;
69}
70
71} // namespace xtp
72} // namespace votca
virtual double CalcCorrelationDiagElement(Index gw_level, double frequency) const =0
Eigen::MatrixXd CalcCorrelationOffDiag(const Eigen::VectorXd &frequencies) const
virtual double CalcCorrelationOffDiagElement(Index gw_level1, Index gw_level2, double frequency1, double frequency2) const =0
Eigen::MatrixXd CalcExchangeMatrix() const
Eigen::VectorXd CalcCorrelationDiag(const Eigen::VectorXd &frequencies) const
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