preCICE v3.3.0
Loading...
Searching...
No Matches
ResidualPreconditioner.cpp
Go to the documentation of this file.
2#include <cstddef>
3#include <vector>
5#include "utils/IntraComm.hpp"
6#include "utils/assertion.hpp"
7
9
11 : Preconditioner(maxNonConstTimeWindows)
12{
13}
14
15void ResidualPreconditioner::_update_(bool timeWindowComplete,
16 const Eigen::VectorXd &oldValues,
17 const Eigen::VectorXd &res)
18{
19 if (not timeWindowComplete) {
20 std::vector<double> norms(_subVectorSizes.size(), 0.0);
21
22 for (size_t k = 0; k < _subVectorSizes.size(); k++) {
23 Eigen::VectorXd part = res.segment(_subVectorOffsets[k], _subVectorSizes[k]);
24 norms[k] = utils::IntraComm::l2norm(part);
25 }
26
27 for (size_t k = 0; k < _subVectorSizes.size(); k++) {
28 if (norms[k] < math::NUMERICAL_ZERO_DIFFERENCE) {
29 PRECICE_WARN("A sub-vector in the residual preconditioner became numerically zero. "
30 "If this occurred in the second iteration and the initial-relaxation factor is equal to 1.0, "
31 "check if the coupling data values of one solver is zero in the first iteration. "
32 "The preconditioner scaling factors were not applied for this iteration.");
33 } else {
34 for (size_t i = 0; i < _subVectorSizes[k]; i++) {
35 PRECICE_ASSERT(norms[k] > 0.0);
36 auto offset = _subVectorOffsets[k];
37 _weights[i + offset] = 1.0 / norms[k];
38 _invWeights[i + offset] = norms[k];
39 }
40 }
41 }
42
43 _requireNewQR = true;
44 }
45}
46
47} // namespace precice::acceleration::impl
#define PRECICE_WARN(...)
Definition LogMacros.hpp:12
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
std::vector< size_t > _subVectorSizes
Sizes of each sub-vector, i.e. each coupling data.
std::vector< double > _invWeights
Inverse weights (for efficiency reasons)
std::vector< double > _weights
Weights used to scale the matrix V and the residual.
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 _update_(bool timeWindowComplete, const Eigen::VectorXd &oldValues, const Eigen::VectorXd &res) override
Update the scaling after every FSI iteration.
static double l2norm(const Eigen::VectorXd &vec)
The l2 norm of a vector is calculated on distributed data.
Definition IntraComm.cpp:67
constexpr double NUMERICAL_ZERO_DIFFERENCE