preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
ResidualSumPreconditioner.cpp
Go to the documentation of this file.
2#include <algorithm>
3#include <cmath>
6#include "utils/IntraComm.hpp"
7#include "utils/assertion.hpp"
8
10
12 int maxNonConstTimeWindows,
13 bool preconditionerUpdateOnThreshold)
14 : Preconditioner(maxNonConstTimeWindows),
15 _preconditionerUpdateOnThreshold(preconditionerUpdateOnThreshold)
16{
17}
18
27
28void ResidualSumPreconditioner::_update_(bool timeWindowComplete,
29 const Eigen::VectorXd &oldValues,
30 const Eigen::VectorXd &res)
31{
32 if (timeWindowComplete) {
33 _firstTimeWindow = false;
35 return;
36 }
37
39
40 double sum = 0.0;
41
42 int offset = 0;
43 for (size_t k = 0; k < _subVectorSizes.size(); k++) {
44 Eigen::VectorXd part = Eigen::VectorXd::Zero(_subVectorSizes[k]);
45 for (size_t i = 0; i < _subVectorSizes[k]; i++) {
46 part(i) = res(i + offset);
47 }
48 norms[k] = utils::IntraComm::dot(part, part);
49 sum += norms[k];
50 offset += _subVectorSizes[k];
51 norms[k] = std::sqrt(norms[k]);
52 }
53 sum = std::sqrt(sum);
55 math::equals(sum, 0.0),
56 "All residual sub-vectors in the residual-sum preconditioner are numerically zero ( sum = {}). "
57 "This indicates that the data values exchanged between two successive iterations did not change. "
58 "The simulation may be unstable, e.g. produces NAN values. Please check the data values exchanged "
59 "between the solvers is not identical between iterations. The preconditioner scaling factors were "
60 "not updated in this iteration and the scaling factors determined in the previous iteration were used.",
61 sum);
62
63 for (size_t k = 0; k < _subVectorSizes.size(); k++) {
65 _residualSum[k] += norms[k] / sum;
66
69 "A sub-vector in the residual-sum preconditioner became numerically zero ( sub-vector = {}). "
70 "If this occurred in the second iteration and the initial-relaxation factor is equal to 1.0, "
71 "check if the coupling data values of one solver is zero in the first iteration. "
72 "The preconditioner scaling factors were not updated for this iteration and the scaling factors "
73 "determined in the previous iteration were used.",
74 _residualSum[k]);
75 }
76
77 // Determine if weights needs to be reset
78 // if _preconditionerUpdateOnThreshold is true, the weights are reset only if the ratio of the new scaling weight to the previous residual sum has changed significantly
80 if (!resetWeights) {
81 for (size_t k = 0; k < _subVectorSizes.size(); k++) {
82 double resSum = _residualSum[k];
83 if (math::equals(resSum, 0.0)) {
84 continue; // These will be ignored when resetting the weights
85 }
86 // Check if the ratio of the new scaling weight to the previous residual sum
87 // has changed significantly, either exceeding the threshold of 10.0
88 // or dropping below its inverse.
89 double factor = _previousResidualSum[k] / resSum;
90 if ((factor > 10.0) || (factor < 0.1)) {
91 resetWeights = true;
92 PRECICE_DEBUG("Significant scaling weight change is detected. The pre-scaling weights will be reset.");
93 break;
94 }
95 }
96 }
97
98 if (!resetWeights) {
99 return;
100 }
101
102 // Reset the weights for non-zero residual sums
103 offset = 0;
104 for (size_t k = 0; k < _subVectorSizes.size(); k++) {
105 double resSum = _residualSum[k];
106 if (not math::equals(resSum, 0.0)) {
107 for (size_t i = 0; i < _subVectorSizes[k]; i++) {
108 _weights[i + offset] = 1 / resSum;
109 _invWeights[i + offset] = resSum;
110 }
111 PRECICE_DEBUG("preconditioner scaling factor[{}] = {}", k, 1 / resSum);
112 }
113 _previousResidualSum[k] = resSum;
114 offset += _subVectorSizes[k];
115 }
116 _requireNewQR = true;
117}
118
119} // namespace precice::acceleration::impl
#define PRECICE_WARN_IF(condition,...)
Definition LogMacros.hpp:18
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:61
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
T begin(T... args)
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.
virtual void initialize(std::vector< size_t > &svs)
initialize the preconditioner
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.
ResidualSumPreconditioner(int maxNonConstTimeWindows, bool preconditionerUpdateOnThreshold)
void initialize(std::vector< size_t > &svs) override
initialize the preconditioner
void _update_(bool timeWindowComplete, const Eigen::VectorXd &oldValues, const Eigen::VectorXd &res) override
Update the scaling after every FSI iteration.
static double dot(const Eigen::VectorXd &vec1, const Eigen::VectorXd &vec2)
T end(T... args)
T fill(T... args)
constexpr bool equals(const Eigen::MatrixBase< DerivedA > &A, const Eigen::MatrixBase< DerivedB > &B, double tolerance=NUMERICAL_ZERO_DIFFERENCE)
Compares two Eigen::MatrixBase for equality up to tolerance.
constexpr double NUMERICAL_ZERO_DIFFERENCE
T resize(T... args)
T size(T... args)
T sqrt(T... args)