votca 2025.1-dev
Loading...
Searching...
No Matches
symmetric_matrix.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 <cassert>
22#include <iostream>
23
24// Local VOTCA includes
26
27namespace votca {
28namespace xtp {
29
30Symmetric_Matrix::Symmetric_Matrix(const Eigen::MatrixXd& full)
31 : dimension(full.rows()) {
32 assert(full.rows() == full.cols() && "Input matrix not quadratic");
33 data.resize((dimension + 1) * dimension / 2);
34 for (Index i = 0; i < full.rows(); ++i) {
35 for (Index j = 0; j <= i; ++j) {
36 this->operator()(i, j) = full(i, j);
37 }
38 }
39}
40
41std::ostream& operator<<(std::ostream& out, const Symmetric_Matrix& a) {
42
43 out << "[" << a.dimension << "," << a.dimension << "]\n";
44 for (Index i = 0; i < a.dimension; ++i) {
45 for (Index j = 0; j <= i; ++j) {
46 out << a(i, j);
47 if (i == j) {
48 out << "\n";
49 } else {
50 out << " ";
51 }
52 }
53 }
54 return out;
55}
56
58 assert(data.size() == a.data.size() && "Matrices do not have same size");
59 double result = 0.0;
60
61 for (Index i = 0; i < dimension; ++i) {
62 const Index index = (i * (i + 1)) / 2 + i;
63 result += +data[index] * a.data[index];
64 }
65 for (Index i = 0; i < dimension; ++i) {
66 const Index start = (i * (i + 1)) / 2;
67 for (Index j = 0; j < i; ++j) {
68 result += 2 * data[start + j] * a.data[start + j];
69 }
70 }
71 return result;
72}
73
74void Symmetric_Matrix::AddtoEigenMatrix(Eigen::MatrixXd& full,
75 double factor) const {
76 for (Index j = 0; j < full.cols(); ++j) {
77 const Index start = (j * (j + 1)) / 2;
78 for (Index i = 0; i <= j; ++i) {
79 full(i, j) += factor * data[start + i];
80 }
81
82 for (Index i = j + 1; i < full.rows(); ++i) {
83 const Index index = (i * (i + 1)) / 2 + j;
84 full(i, j) += factor * data[index];
85 }
86 }
87 return;
88}
89
91 Eigen::SelfAdjointView<Eigen::MatrixXd, Eigen::Upper>& upper,
92 double factor) const {
93 for (Index j = 0; j < upper.cols(); ++j) {
94 const Index start = (j * (j + 1)) / 2;
95 for (Index i = 0; i <= j; ++i) {
96 upper(i, j) += factor * data[start + i];
97 }
98 }
99 return;
100}
101
102Eigen::MatrixXd Symmetric_Matrix::FullMatrix() const {
103 Eigen::MatrixXd result = Eigen::MatrixXd(dimension, dimension);
104 for (Index j = 0; j < result.cols(); ++j) {
105 const Index start = (j * (j + 1)) / 2;
106 for (Index i = 0; i <= j; ++i) {
107 result(i, j) = data[start + i];
108 }
109
110 for (Index i = j + 1; i < result.rows(); ++i) {
111 const Index index = (i * (i + 1)) / 2 + j;
112 result(i, j) = data[index];
113 }
114 }
115 return result;
116}
117
118Eigen::MatrixXd Symmetric_Matrix::UpperMatrix() const {
119 Eigen::MatrixXd result = Eigen::MatrixXd::Zero(dimension, dimension);
120 for (Index j = 0; j < result.cols(); ++j) {
121 const Index start = (j * (j + 1)) / 2;
122 for (Index i = 0; i <= j; ++i) {
123 result(i, j) = data[start + i];
124 }
125 }
126 return result;
127}
128
130 Index index;
131 if (i >= j) {
132 index = (i * (i + 1)) / 2 + j;
133 } else {
134 index = (j * (j + 1)) / 2 + i;
135 }
136 return index;
137}
138
139} // namespace xtp
140} // namespace votca
Index index(Index i, Index j) const
void AddtoEigenUpperMatrix(Eigen::SelfAdjointView< Eigen::MatrixXd, Eigen::Upper > &upper, double factor=1.0) const
double TraceofProd(const Symmetric_Matrix &a) const
Eigen::MatrixXd FullMatrix() const
double & operator()(Index i, Index j)
std::vector< double > data
void AddtoEigenMatrix(Eigen::MatrixXd &full, double factor=1.0) const
Eigen::MatrixXd UpperMatrix() const
std::ostream & operator<<(std::ostream &out, const Correlate &c)
Definition correlate.h:53
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