preCICE v3.1.1
Loading...
Searching...
No Matches
Acceleration.cpp
Go to the documentation of this file.
3#include "utils/Helpers.hpp"
4
6
7void Acceleration::checkDataIDs(const DataMap &cplData) const
8{
9#ifndef NDEBUG
10 for (int id : getDataIDs()) {
11 bool valid = utils::contained(id, cplData);
12 PRECICE_ASSERT(valid, "Data with ID {} unknown.", id);
13 }
14#endif
15}
16
17void Acceleration::applyRelaxation(double omega, DataMap &cplData)
18{
19 for (auto &pair : cplData) {
20 auto &couplingData = *(pair.second);
21
22 if (couplingData.timeStepsStorage().empty()) {
23 continue;
24 }
25 // @todo: This will apply the scaling to the sample at t=0 and t=1 when
26 // calling performAcceleration the first time. Computations at the
27 // previousValuesAtTime/oldValues don't change anything and are unneeded
28 for (auto &stample : couplingData.timeStepsStorage().stamples()) {
29 auto &values = stample.sample.values;
30 auto oldValues = couplingData.getPreviousValuesAtTime(stample.timestamp); // IMPORTANT DETAIL: The interpolation that we use for resampling does not necessarily have to be the same interpolation as the interpolation the user accesses via read-data. (But probably it is easier to just use the same)
31 values = values * omega + oldValues * (1.0 - omega);
32
33 if (couplingData.hasGradient()) {
34 auto &gradients = stample.sample.gradients;
35 auto oldGradients = couplingData.getPreviousGradientsAtTime(stample.timestamp); // IMPORTANT DETAIL: The interpolation that we use for resampling does not necessarily have to be the same interpolation as the interpolation the user accesses via read-data. (But probably it is easier to just use the same)
36 gradients = gradients * omega + oldGradients * (1.0 - omega);
37 }
38 }
39
40 // @todo remove
41 // update the "sample"
42 couplingData.sample() = couplingData.timeStepsStorage().last().sample;
43 }
44}
45
47 const DataMap &cplData, const std::vector<DataID> &dataIDs, Eigen::VectorXd &targetValues, Eigen::VectorXd &targetOldValues) const
48{
49 Eigen::Index offset = 0;
50 for (auto id : dataIDs) {
51 Eigen::Index size = cplData.at(id)->values().size();
52 auto & values = cplData.at(id)->values();
53 const auto & oldValues = cplData.at(id)->previousIteration();
54 PRECICE_ASSERT(targetValues.size() >= offset + size, "Target vector was not initialized.", targetValues.size(), offset + size);
55 PRECICE_ASSERT(targetOldValues.size() >= offset + size, "Target vector was not initialized.");
56 for (Eigen::Index i = 0; i < size; i++) {
57 targetValues(i + offset) = values(i);
58 targetOldValues(i + offset) = oldValues(i);
59 }
60 offset += size;
61 }
62}
63} // namespace precice::acceleration
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
T at(T... args)
void checkDataIDs(const DataMap &cplData) const
Checks if all dataIDs are contained in cplData.
static void applyRelaxation(double omega, DataMap &cplData)
performs a relaxation given a relaxation factor omega
void concatenateCouplingData(const DataMap &cplData, const std::vector< DataID > &dataIDs, Eigen::VectorXd &targetValues, Eigen::VectorXd &targetOldValues) const
Concatenates all coupling data involved into a single vector.
virtual std::vector< int > getDataIDs() const =0
contains implementations of acceleration schemes.
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:39