preCICE
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 <string>
6#include <vector>
7
10#include "logging/Logger.hpp"
11#include "utils/assertion.hpp"
12
14
25public:
26 Preconditioner(int maxNonConstTimeWindows)
27 : _maxNonConstTimeWindows(maxNonConstTimeWindows)
28 {
29 }
30
32 virtual ~Preconditioner() = default;
33
39 virtual void initialize(std::vector<size_t> svs, std::vector<std::string> svNames)
40 {
42 PRECICE_ASSERT(svs.size() == svNames.size());
43
44 _subVectorSizes = std::move(svs);
45 _subVectorNames = std::move(svNames);
46
47 // Compute offsets of each subvector
48 _subVectorOffsets.resize(_subVectorSizes.size(), 0);
49 std::partial_sum(_subVectorSizes.begin(), --_subVectorSizes.end(), ++_subVectorOffsets.begin());
50
51 size_t N = std::accumulate(_subVectorSizes.begin(), _subVectorSizes.end(), static_cast<std::size_t>(0));
52
53 // cannot do this already in the constructor as the size is unknown at that point
54 _weights.resize(N, 1.0);
55 _invWeights.resize(N, 1.0);
56 }
57
62 void apply(Eigen::MatrixXd &M, bool transpose)
63 {
65 if (transpose) {
66 PRECICE_DEBUG_IF((int) _weights.size() != M.cols(), "The number of columns of the matrix {} and weights size {} mismatched.", M.rows(), _weights.size());
67
68 int validCols = std::min(static_cast<int>(M.cols()), (int) _weights.size());
69 for (int i = 0; i < validCols; i++) {
70 for (int j = 0; j < M.rows(); j++) {
71 M(j, i) *= _weights[i];
72 }
73 }
74 } else {
75 PRECICE_DEBUG_IF((int) _weights.size() != M.rows(), "The number of rows of the matrix {} and weights size {} mismatched.", M.rows(), _weights.size());
76
77 int validRows = std::min(static_cast<int>(M.rows()), (int) _weights.size());
78 for (int i = 0; i < M.cols(); i++) {
79 for (int j = 0; j < validRows; j++) {
80 M(j, i) *= _weights[j];
81 }
82 }
83 }
84 }
85
90 void revert(Eigen::MatrixXd &M, bool transpose)
91 {
93 // PRECICE_ASSERT(_needsGlobalWeights);
94 if (transpose) {
95 PRECICE_DEBUG_IF((int) _weights.size() != M.cols(), "The number of columns of the matrix {} and weights size {} mismatched.", M.cols(), _weights.size());
96
97 int validCols = std::min(static_cast<int>(M.cols()), (int) _weights.size());
98 for (int i = 0; i < validCols; i++) {
99 for (int j = 0; j < M.rows(); j++) {
100 M(j, i) *= _invWeights[i];
101 }
102 }
103 } else {
104 PRECICE_DEBUG_IF((int) _weights.size() != M.rows(), "The number of rows of the matrix {} and weights size {} mismatched.", M.rows(), _weights.size());
105
106 int validRows = std::min(static_cast<int>(M.rows()), (int) _weights.size());
107 for (int i = 0; i < M.cols(); i++) {
108 for (int j = 0; j < validRows; j++) {
109 M(j, i) *= _invWeights[j];
110 }
111 }
112 }
113 }
114
116 void apply(Eigen::MatrixXd &M)
117 {
119 PRECICE_DEBUG_IF((int) _weights.size() != M.rows(), "The number of rows of the matrix {} and weights size {} mismatched.", M.rows(), _weights.size());
120
121 // scale matrix M
122 int validRows = std::min(static_cast<int>(M.rows()), (int) _weights.size());
123 for (int i = 0; i < M.cols(); i++) {
124 for (int j = 0; j < validRows; j++) {
125 M(j, i) *= _weights[j];
126 }
127 }
128 }
129
131 void apply(Eigen::VectorXd &v)
132 {
134 PRECICE_DEBUG_IF((int) _weights.size() != v.size(), "The vector size {} and weights size {} mismatched.", v.size(), _weights.size());
135
136 // scale vector
137 int validSize = std::min(static_cast<int>(v.size()), (int) _weights.size());
138 for (int j = 0; j < validSize; j++) {
139 v[j] *= _weights[j];
140 }
141 }
142
144 void revert(Eigen::MatrixXd &M)
145 {
147 PRECICE_DEBUG_IF((int) _weights.size() != M.rows(), "The number of rows of the matrix {} and weights size {} mismatched.", M.rows(), _weights.size());
148
149 // scale matrix M
150 int validRows = std::min(static_cast<int>(M.rows()), (int) _weights.size());
151 for (int i = 0; i < M.cols(); i++) {
152 for (int j = 0; j < validRows; j++) {
153 M(j, i) *= _invWeights[j];
154 }
155 }
156 }
157
159 void revert(Eigen::VectorXd &v)
160 {
162 PRECICE_DEBUG_IF((int) _weights.size() != v.size(), "The vector size {} and weights size {} mismatched.", v.size(), _weights.size());
163
164 // revert vector scaling
165 int validSize = std::min(static_cast<int>(v.size()), (int) _weights.size());
166 for (int j = 0; j < validSize; j++) {
167 v[j] *= _invWeights[j];
168 }
169 }
170
176 void update(bool timeWindowComplete, const Eigen::VectorXd &oldValues, const Eigen::VectorXd &res)
177 {
179
180 // if number of allowed non-const time windows is exceeded, do not update weights
181 if (_frozen)
182 return;
183
184 // increment number of time windows that has been scaled with changing preconditioning weights
185 if (timeWindowComplete) {
188 _frozen = true;
189 }
190
191 // type specific update functionality
192 _update_(timeWindowComplete, oldValues, res);
193 }
194
197 {
199 return _requireNewQR;
200 }
201
204 {
205 _requireNewQR = false;
206 }
207
208 std::vector<double> &getWeights()
209 {
210 return _weights;
211 }
212
213 bool isConst()
214 {
215 return _frozen;
216 }
217
218protected:
220 std::vector<double> _weights;
221
223 std::vector<double> _invWeights;
224
226 std::vector<size_t> _subVectorSizes;
227
229 std::vector<std::string> _subVectorNames;
230
232 std::vector<size_t> _subVectorOffsets;
233
238
241
243 bool _requireNewQR = false;
244
246 bool _frozen = false;
247
253 virtual void _update_(bool timeWindowComplete, const Eigen::VectorXd &oldValues, const Eigen::VectorXd &res) = 0;
254
255private:
256 logging::Logger _log{"acceleration::Preconditioner"};
257};
258
259} // namespace precice::acceleration::impl
#define PRECICE_DEBUG_IF(condition,...)
Definition LogMacros.hpp:63
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
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, std::vector< std::string > svNames)
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.
std::vector< std::string > _subVectorNames
Names of the coupling data corresponding to each sub-vector.
void revert(Eigen::VectorXd &v)
To transform balanced values back to physical values. Vector version.
std::vector< size_t > _subVectorOffsets
Offsets of each sub-vector in concatenated data, i.e. each coupling data.
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