preCICE v3.1.2
Loading...
Searching...
No Matches
ValuePreconditioner.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 int maxNonConstTimeWindows)
12 : Preconditioner(maxNonConstTimeWindows)
13{
14}
15
16void ValuePreconditioner::_update_(bool timeWindowComplete,
17 const Eigen::VectorXd &oldValues,
18 const Eigen::VectorXd &res)
19{
20 if (timeWindowComplete || _firstTimeWindow) {
21
23
24 int offset = 0;
25 for (size_t k = 0; k < _subVectorSizes.size(); k++) {
26 Eigen::VectorXd part = Eigen::VectorXd::Zero(_subVectorSizes[k]);
27 for (size_t i = 0; i < _subVectorSizes[k]; i++) {
28 part(i) = oldValues(i + offset);
29 }
30 norms[k] = utils::IntraComm::l2norm(part);
31 offset += _subVectorSizes[k];
32 }
33
34 offset = 0;
35 for (size_t k = 0; k < _subVectorSizes.size(); k++) {
36 if (norms[k] < math::NUMERICAL_ZERO_DIFFERENCE) {
37 PRECICE_WARN("A sub-vector in the residual preconditioner became numerically zero. "
38 "If this occurred in the second iteration and the initial-relaxation factor is equal to 1.0, "
39 "check if the coupling data values of one solver is zero in the first iteration. "
40 "The preconditioner scaling factors were not applied for this iteration.");
41 } else {
42 for (size_t i = 0; i < _subVectorSizes[k]; i++) {
43 PRECICE_ASSERT(norms[k] > 0.0);
44 _weights[i + offset] = 1.0 / norms[k];
45 _invWeights[i + offset] = norms[k];
46 }
47 }
48 offset += _subVectorSizes[k];
49 }
50
51 _requireNewQR = true;
52 _firstTimeWindow = false;
53 }
54}
55
56} // namespace precice::acceleration::impl
#define PRECICE_WARN(...)
Definition LogMacros.hpp:11
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
Interface for preconditioner variants that can be applied to quasi-Newton acceleration schemes.
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.
bool _requireNewQR
True if a QR decomposition from scratch is necessary.
virtual void _update_(bool timeWindowComplete, const Eigen::VectorXd &oldValues, const Eigen::VectorXd &res)
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
T size(T... args)