preCICE
Loading...
Searching...
No Matches
AitkenAcceleration.cpp
Go to the documentation of this file.
2#include <Eigen/Core>
3#include <boost/range/adaptor/map.hpp>
4#include <cmath>
5#include <cstddef>
6#include <limits>
7#include <map>
8#include <memory>
9#include <numeric>
10#include <ostream>
11#include <utility>
12
14#include "logging/LogMacros.hpp"
15#include "math/math.hpp"
17#include "utils/Helpers.hpp"
18#include "utils/IntraComm.hpp"
19#include "utils/assertion.hpp"
20
21namespace precice::acceleration {
22
24 std::vector<int> dataIDs,
25 impl::PtrPreconditioner preconditioner)
26 : _initialRelaxation(initialRelaxation),
27 _primaryDataIDs(std::move(dataIDs)),
28 _aitkenFactor(initialRelaxation),
29 _preconditioner(std::move(preconditioner))
30{
32 "Initial relaxation factor for Aitken acceleration has to "
33 "be larger than zero and smaller or equal to one. "
34 "Current initial relaxation is: {}",
36}
37
39{
40 checkDataIDs(cplData);
41
42 // Accumulate number of entries
43 // Size for each subvector needed for preconditioner
44 std::vector<std::size_t> subVectorSizes;
45 // Gather sizes
46 std::transform(_primaryDataIDs.cbegin(), _primaryDataIDs.cend(), std::back_inserter(subVectorSizes), [&cplData](const auto &d) { return cplData.at(d)->getSize(); });
47 // The accumulated sum
48 Eigen::Index entries = std::accumulate(subVectorSizes.cbegin(), subVectorSizes.cend(), static_cast<Eigen::Index>(0));
49
50 // Allocate memory
51 _oldResiduals = Eigen::VectorXd::Zero(entries);
52 _values = Eigen::VectorXd::Zero(entries);
53 _oldValues = Eigen::VectorXd::Zero(entries);
54 if (_primaryDataIDs.size() > 1) {
55 std::vector<std::string> subVectorNames;
56 std::transform(_primaryDataIDs.cbegin(), _primaryDataIDs.cend(), std::back_inserter(subVectorNames), [&cplData](const auto &d) { return cplData.at(d)->getDataName(); });
57 _preconditioner->initialize(std::move(subVectorSizes), std::move(subVectorNames));
58 }
59}
60
62 DataMap &cplData,
63 double windowStart,
64 double windowEnd)
65{
67
68 // Compute aitken relaxation factor
70
71 concatenateCouplingData(cplData, _primaryDataIDs, _values, _oldValues, windowStart, windowEnd);
72
73 // Compute current residual = values - oldValues
74 Eigen::VectorXd residuals = _values - _oldValues;
75
76 // Compute residual deltas (= residuals - oldResiduals) and store it in _oldResiduals
77 Eigen::VectorXd residualDeltas = residuals - _oldResiduals;
78
79 // We need to update the preconditioner in every iteration
80 if (_primaryDataIDs.size() > 1) {
81 _preconditioner->update(false, _values, residuals);
82 }
83
84 // Select/compute aitken factor depending on current iteration count
85 if (_iterationCounter == 0) {
86 // preconditioner not necessary
88 } else {
89 // If we have more than one data set, we scale the data to get a better approximation
90 // of the Aitken factor
91 if (_primaryDataIDs.size() > 1) {
92 _preconditioner->apply(residualDeltas);
94 }
95 // compute fraction of aitken factor with residuals and residual deltas
96 double nominator = utils::IntraComm::dot(_oldResiduals, residualDeltas);
97 double denominator = utils::IntraComm::dot(residualDeltas, residualDeltas);
98 _aitkenFactor = -_aitkenFactor * (nominator / denominator);
99 }
100
101 PRECICE_DEBUG("AitkenFactor: {}", _aitkenFactor);
102
103 // Perform relaxation with aitken factor
104 applyRelaxation(_aitkenFactor, cplData, windowStart);
105
106 // Store residuals for next iteration
107 _oldResiduals = std::move(residuals);
108
110}
111
113 const DataMap &cplData, double windowStart)
114{
116 if (_primaryDataIDs.size() > 1) {
117 _preconditioner->update(true, _values, _oldResiduals);
118 }
119 _oldResiduals = Eigen::VectorXd::Constant(_oldResiduals.size(), std::numeric_limits<double>::max());
120}
121
123 const DataMap &cplData, const std::vector<DataID> &dataIDs, Eigen::VectorXd &targetValues, Eigen::VectorXd &targetOldValues, double windowStart, double windowEnd) const
124{
125 Eigen::Index offset = 0;
126
127 for (auto id : dataIDs) {
128 Eigen::Index size = cplData.at(id)->getSize();
129
130 auto valuesSample = cplData.at(id)->waveform().sample(windowEnd);
131 auto oldValuesSample = cplData.at(id)->getPreviousValuesAtTime(windowEnd);
132
133 PRECICE_ASSERT(valuesSample.values().size() == size, valuesSample.values().size(), size);
134 PRECICE_ASSERT(oldValuesSample.values().size() == size, oldValuesSample.values().size(), size);
135 PRECICE_ASSERT(targetValues.size() >= offset + size, "Target vector was not initialized.", targetValues.size(), offset + size);
136 PRECICE_ASSERT(targetOldValues.size() >= offset + size, "Target vector was not initialized.");
137 for (Eigen::Index i = 0; i < size; i++) {
138 targetValues(i + offset) = valuesSample.values()(i);
139 targetOldValues(i + offset) = oldValuesSample.values()(i);
140 }
141 offset += size;
142 }
143}
144} // namespace precice::acceleration
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:61
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:32
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
static void applyRelaxation(double omega, DataMap &cplData, double windowStart)
performs a relaxation given a relaxation factor omega
void checkDataIDs(const DataMap &cplData) const
Checks if all dataIDs are contained in cplData.
std::map< int, cplscheme::PtrCouplingData > DataMap
Map from data ID to data values.
void performAcceleration(DataMap &cpldata, double windowStart, double windowEnd) final override
void initialize(const DataMap &cpldata) final override
AitkenAcceleration(double initialRelaxationFactor, std::vector< int > dataIDs, impl::PtrPreconditioner preconditioner)
void concatenateCouplingData(const DataMap &cplData, const std::vector< DataID > &dataIDs, Eigen::VectorXd &targetValues, Eigen::VectorXd &targetOldValues, double windowStart, double windowEnd) const
Concatenates the data and old data in cplData into two long vectors.
impl::PtrPreconditioner _preconditioner
Preconditioner for data vector if multiple data sets are used.
void iterationsConverged(const DataMap &cpldata, double windowStart) final override
static double dot(const Eigen::VectorXd &vec1, const Eigen::VectorXd &vec2)
std::shared_ptr< Preconditioner > PtrPreconditioner
contains implementations of acceleration schemes.
int sign(double number)
Return the sign, one of {-1, 0, 1}.
Definition math.hpp:10
bool contained(const ELEMENT_T &element, const std::vector< ELEMENT_T > &vec)
Returns true, if given element is in vector, otherwise false.
Definition Helpers.hpp:38
STL namespace.