preCICE v3.2.0
Loading...
Searching...
No Matches
Preconditioner.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Core>
4#include <numeric>
5#include <vector>
6
9#include "logging/Logger.hpp"
10#include "utils/assertion.hpp"
11
13
24public:
25 Preconditioner(int maxNonConstTimeWindows)
26 : _maxNonConstTimeWindows(maxNonConstTimeWindows)
27 {
28 }
29
31 virtual ~Preconditioner() = default;
32
38 {
40
41 _subVectorSizes = svs;
42
43 size_t N = std::accumulate(_subVectorSizes.begin(), _subVectorSizes.end(), static_cast<std::size_t>(0));
44
45 // cannot do this already in the constructor as the size is unknown at that point
46 _weights.resize(N, 1.0);
47 _invWeights.resize(N, 1.0);
48 }
49
54 void apply(Eigen::MatrixXd &M, bool transpose)
55 {
57 if (transpose) {
58 PRECICE_DEBUG_IF((int) _weights.size() != M.cols(), "The number of columns of the matrix {} and weights size {} mismatched.", M.rows(), _weights.size());
59
60 int validCols = std::min(static_cast<int>(M.cols()), (int) _weights.size());
61 for (int i = 0; i < validCols; i++) {
62 for (int j = 0; j < M.rows(); j++) {
63 M(j, i) *= _weights[i];
64 }
65 }
66 } else {
67 PRECICE_DEBUG_IF((int) _weights.size() != M.rows(), "The number of rows of the matrix {} and weights size {} mismatched.", M.rows(), _weights.size());
68
69 int validRows = std::min(static_cast<int>(M.rows()), (int) _weights.size());
70 for (int i = 0; i < M.cols(); i++) {
71 for (int j = 0; j < validRows; j++) {
72 M(j, i) *= _weights[j];
73 }
74 }
75 }
76 }
77
82 void revert(Eigen::MatrixXd &M, bool transpose)
83 {
85 // PRECICE_ASSERT(_needsGlobalWeights);
86 if (transpose) {
87 PRECICE_DEBUG_IF((int) _weights.size() != M.cols(), "The number of columns of the matrix {} and weights size {} mismatched.", M.cols(), _weights.size());
88
89 int validCols = std::min(static_cast<int>(M.cols()), (int) _weights.size());
90 for (int i = 0; i < validCols; i++) {
91 for (int j = 0; j < M.rows(); j++) {
92 M(j, i) *= _invWeights[i];
93 }
94 }
95 } else {
96 PRECICE_DEBUG_IF((int) _weights.size() != M.rows(), "The number of rows of the matrix {} and weights size {} mismatched.", M.rows(), _weights.size());
97
98 int validRows = std::min(static_cast<int>(M.rows()), (int) _weights.size());
99 for (int i = 0; i < M.cols(); i++) {
100 for (int j = 0; j < validRows; j++) {
101 M(j, i) *= _invWeights[j];
102 }
103 }
104 }
105 }
106
108 void apply(Eigen::MatrixXd &M)
109 {
111 PRECICE_DEBUG_IF((int) _weights.size() != M.rows(), "The number of rows of the matrix {} and weights size {} mismatched.", M.rows(), _weights.size());
112
113 // scale matrix M
114 int validRows = std::min(static_cast<int>(M.rows()), (int) _weights.size());
115 for (int i = 0; i < M.cols(); i++) {
116 for (int j = 0; j < validRows; j++) {
117 M(j, i) *= _weights[j];
118 }
119 }
120 }
121
123 void apply(Eigen::VectorXd &v)
124 {
126 PRECICE_DEBUG_IF((int) _weights.size() != v.size(), "The vector size {} and weights size {} mismatched.", v.size(), _weights.size());
127
128 // scale vector
129 int validSize = std::min(static_cast<int>(v.size()), (int) _weights.size());
130 for (int j = 0; j < validSize; j++) {
131 v[j] *= _weights[j];
132 }
133 }
134
136 void revert(Eigen::MatrixXd &M)
137 {
139 PRECICE_DEBUG_IF((int) _weights.size() != M.rows(), "The number of rows of the matrix {} and weights size {} mismatched.", M.rows(), _weights.size());
140
141 // scale matrix M
142 int validRows = std::min(static_cast<int>(M.rows()), (int) _weights.size());
143 for (int i = 0; i < M.cols(); i++) {
144 for (int j = 0; j < validRows; j++) {
145 M(j, i) *= _invWeights[j];
146 }
147 }
148 }
149
151 void revert(Eigen::VectorXd &v)
152 {
154 PRECICE_DEBUG_IF((int) _weights.size() != v.size(), "The vector size {} and weights size {} mismatched.", v.size(), _weights.size());
155
156 // revert vector scaling
157 int validSize = std::min(static_cast<int>(v.size()), (int) _weights.size());
158 for (int j = 0; j < validSize; j++) {
159 v[j] *= _invWeights[j];
160 }
161 }
162
168 void update(bool timeWindowComplete, const Eigen::VectorXd &oldValues, const Eigen::VectorXd &res)
169 {
171
172 // if number of allowed non-const time windows is exceeded, do not update weights
173 if (_frozen)
174 return;
175
176 // increment number of time windows that has been scaled with changing preconditioning weights
177 if (timeWindowComplete) {
180 _frozen = true;
181 }
182
183 // type specific update functionality
184 _update_(timeWindowComplete, oldValues, res);
185 }
186
189 {
191 return _requireNewQR;
192 }
193
196 {
197 _requireNewQR = false;
198 }
199
201 {
202 return _weights;
203 }
204
205 bool isConst()
206 {
207 return _frozen;
208 }
209
210protected:
213
216
219
224
227
229 bool _requireNewQR = false;
230
232 bool _frozen = false;
233
239 virtual void _update_(bool timeWindowComplete, const Eigen::VectorXd &oldValues, const Eigen::VectorXd &res) = 0;
240
241private:
242 logging::Logger _log{"acceleration::Preconditioner"};
243};
244
245} // namespace precice::acceleration::impl
#define PRECICE_DEBUG_IF(condition,...)
Definition LogMacros.hpp:63
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
T accumulate(T... args)
void apply(Eigen::MatrixXd &M, bool transpose)
Apply preconditioner to matrix.
void revert(Eigen::MatrixXd &M, bool transpose)
Apply inverse preconditioner to matrix.
int _nbNonConstTimeWindows
Counts the number of completed time windows with a non-const weighting.
std::vector< size_t > _subVectorSizes
Sizes of each sub-vector, i.e. each coupling data.
void newQRfulfilled()
to tell the preconditioner that QR-decomposition has been recomputed
void update(bool timeWindowComplete, const Eigen::VectorXd &oldValues, const Eigen::VectorXd &res)
Update the scaling after every FSI iteration and require a new QR decomposition (if necessary)
virtual void _update_(bool timeWindowComplete, const Eigen::VectorXd &oldValues, const Eigen::VectorXd &res)=0
Update the scaling after every FSI iteration and require a new QR decomposition (if necessary)
virtual void initialize(std::vector< size_t > &svs)
initialize the preconditioner
virtual ~Preconditioner()=default
Destructor, empty.
void revert(Eigen::MatrixXd &M)
To transform balanced values back to physical values. Matrix version.
std::vector< double > _invWeights
Inverse weights (for efficiency reasons)
bool requireNewQR()
returns true if a QR decomposition from scratch is necessary
int _maxNonConstTimeWindows
maximum number of non-const time windows, i.e., after this number of time windows,...
std::vector< double > _weights
Weights used to scale the matrix V and the residual.
void revert(Eigen::VectorXd &v)
To transform balanced values back to physical values. Vector version.
bool _requireNewQR
True if a QR decomposition from scratch is necessary.
void apply(Eigen::VectorXd &v)
To transform physical values to balanced values. Vector version.
void apply(Eigen::MatrixXd &M)
To transform physical values to balanced values. Matrix version.
bool _frozen
True if _nbNonConstTimeWindows >= _maxNonConstTimeWindows, i.e., preconditioner is not updated any mo...
This class provides a lightweight logger.
Definition Logger.hpp:17
T min(T... args)