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