preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
TimeGrids.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <numeric>
3#include <vector>
4
5#include "TimeGrids.hpp"
7#include "utils/assertion.hpp"
8
9namespace precice::time {
10
11TimeGrids::TimeGrids(const DataMap &cplData, std::vector<int> dataIDs, bool reducedTimeGrid)
12{
13 for (int dataID : dataIDs) {
14 Eigen::VectorXd timeGrid = cplData.at(dataID)->timeStepsStorage().getTimes();
15 if (reducedTimeGrid) {
16 _timeGrids.insert(std::pair<int, Eigen::VectorXd>(dataID, timeGrid.tail<1>()));
17 } else {
19 }
20 }
21}
22
23Eigen::VectorXd TimeGrids::getTimeGridAfter(int dataID, double time) const
24{
25 PRECICE_ASSERT(_timeGrids.count(dataID), "there does not exists a stored time grid corresponding to this dataID");
26
27 std::vector<double> reduced;
28 for (double d : _timeGrids.at(dataID)) {
29 if (math::greater(d, time)) {
30 reduced.push_back(d);
31 }
32 }
33
34 return Eigen::Map<Eigen::VectorXd>(reduced.data(), reduced.size());
35}
36
38{
39 for (auto &pair : _timeGrids) {
40 if (pair.second.size() == 1) {
41 _timeGrids.at(pair.first) = cplData.at(pair.first)->timeStepsStorage().getTimes().tail<1>();
42 } else {
43 int dataID = pair.first;
44 // Only way to access the first time stamp is through the whole vector
45 Eigen::VectorXd newtimeGrid = cplData.at(dataID)->timeStepsStorage().getTimes();
46 double newTimesMin = newtimeGrid(0);
47 double newTimesMax = newtimeGrid(newtimeGrid.size() - 1);
48
49 Eigen::VectorXd timeGrid = pair.second;
50 double oldTimesMin = _timeGrids.at(dataID)(0);
51 double oldTimesMax = _timeGrids.at(dataID)(timeGrid.size() - 1);
52
53 // Linearly transforms the time grid from the old time window [t_{N-1}, t_N] to the new time window [t_N, t_{N+1}] by translating and scaling the time grid
54 auto transformNewTime = [oldTimesMin, oldTimesMax, newTimesMin, newTimesMax](double t) -> double { return (t - oldTimesMin) / (oldTimesMax - oldTimesMin) * (newTimesMax - newTimesMin) + newTimesMin; };
55 timeGrid = timeGrid.unaryExpr(transformNewTime);
56 _timeGrids.at(dataID) = timeGrid;
57 }
58 }
59}
60} // namespace precice::time
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
T at(T... args)
Eigen::VectorXd getTimeGridAfter(int dataID, double time) const
Definition TimeGrids.cpp:23
TimeGrids(const DataMap &cplData, std::vector< int > dataIDs, bool reduced)
saves the current time grid of the couplig data
Definition TimeGrids.cpp:11
std::map< int, Eigen::VectorXd > _timeGrids
List of the time grid to which all the data will be interpolated to Stored in a map,...
Definition TimeGrids.hpp:46
void moveTimeGridToNewWindow(const DataMap &cplData)
Definition TimeGrids.cpp:37
T count(T... args)
T data(T... args)
T insert(T... args)
std::enable_if< std::is_arithmetic< Scalar >::value, bool >::type greater(Scalar A, Scalar B, Scalar tolerance=NUMERICAL_ZERO_DIFFERENCE)
contains the time interpolation logic.
Definition Sample.hpp:8
T push_back(T... args)
T size(T... args)