preCICE v3.2.0
Loading...
Searching...
No Matches
Storage.cpp
Go to the documentation of this file.
1#include <boost/range.hpp>
2
4#include "math/Bspline.hpp"
6#include "time/Storage.hpp"
7#include "time/Time.hpp"
8#include "utils/assertion.hpp"
9
10namespace precice::time {
11
16
18{
19 this->clear();
20 this->_degree = other.getInterpolationDegree();
21 for (const auto &stample : other.stamples()) {
22 this->setSampleAtTime(stample.timestamp, stample.sample);
23 }
24 return *this;
25}
26
28{
29 // The spline has to be recomputed, since the underlying data has changed
30 _bspline.reset();
31
32 if (_stampleStorage.empty()) {
33 _stampleStorage.emplace_back(Stample{time, sample});
34 return;
35 }
36
37 const double currentWindowStart = _stampleStorage.front().timestamp;
38
39 PRECICE_ASSERT(not sample.values.hasNaN());
40 PRECICE_ASSERT(math::smallerEquals(currentWindowStart, time), "Setting sample outside of valid range!", currentWindowStart, time);
41 // check if key "time" exists.
42 auto existingSample = std::find_if(_stampleStorage.begin(), _stampleStorage.end(), [&time](const auto &s) { return math::equals(s.timestamp, time); });
43 if (existingSample == _stampleStorage.end()) { // key does not exist yet
44 PRECICE_ASSERT(math::smaller(maxStoredTime(), time), maxStoredTime(), time, "Trying to write sample with a time that is too small. Please use clear(), if you want to write new samples to the storage.");
45 _stampleStorage.emplace_back(Stample{time, sample});
46 } else {
47 // Overriding sample
48 existingSample->sample = sample;
49 }
50}
51
53{
54 for (auto &stample : _stampleStorage) {
55 stample.sample = sample;
56 }
57}
58
59void Storage::setInterpolationDegree(int interpolationDegree)
60{
61 PRECICE_ASSERT(interpolationDegree >= Time::MIN_WAVEFORM_DEGREE);
62 _degree = interpolationDegree;
63
64 // The spline has to be recomputed, since the underlying data has changed
65 _bspline.reset();
66}
67
69{
70 return _degree;
71}
72
74{
75 if (_stampleStorage.size() == 0) {
76 return -1; // invalid return
77 } else {
78 return _stampleStorage.back().timestamp;
79 }
80}
81
82int Storage::nTimes() const
83{
84 return _stampleStorage.size();
85}
86
87int Storage::nDofs() const
88{
90 return _stampleStorage[0].sample.values.size();
91}
92
94{
95 PRECICE_ASSERT(nTimes() >= 2, "Calling Storage::move() is only allowed, if there is a sample at the beginning and at the end. This ensures that this function is only called at the end of the window.", getTimes());
96 PRECICE_ASSERT(!_stampleStorage.empty(), "Storage does not contain any data!");
97 const double nextWindowStart = _stampleStorage.back().timestamp;
98 _stampleStorage.erase(_stampleStorage.begin(), --_stampleStorage.end());
99 PRECICE_ASSERT(nextWindowStart == _stampleStorage.front().timestamp);
100
101 // The spline has to be recomputed, since the underlying data has changed
102 _bspline.reset();
103}
104
106{
107 PRECICE_ASSERT(!_stampleStorage.empty(), "Storage does not contain any data!");
108 const double thisWindowStart = _stampleStorage.front().timestamp;
109 _stampleStorage.erase(++_stampleStorage.begin(), _stampleStorage.end());
110 PRECICE_ASSERT(_stampleStorage.size() == 1);
111 PRECICE_ASSERT(thisWindowStart == _stampleStorage.front().timestamp);
112
113 // The spline has to be recomputed, since the underlying data has changed
114 _bspline.reset();
115}
116
118{
119 _stampleStorage.clear();
120 PRECICE_ASSERT(_stampleStorage.size() == 0);
121
122 // The spline has to be recomputed, since the underlying data has changed
123 _bspline.reset();
124}
125
127{
128 if (_stampleStorage.empty()) {
129 return;
130 }
131 _stampleStorage.erase(_stampleStorage.begin(), --_stampleStorage.end());
132
133 // The spline has to be recomputed, since the underlying data has changed
134 _bspline.reset();
135}
136
138{
139 auto beforeTime = [time](const auto &s) { return math::smaller(s.timestamp, time); };
140 _stampleStorage.erase(std::remove_if(_stampleStorage.begin(), _stampleStorage.end(), beforeTime), _stampleStorage.end());
141
142 // The spline has to be recomputed, since the underlying data has changed
143 _bspline.reset();
144}
145
147{
148 auto afterTime = [time](const auto &s) { return math::greater(s.timestamp, time); };
149 _stampleStorage.erase(std::remove_if(_stampleStorage.begin(), _stampleStorage.end(), afterTime), _stampleStorage.end());
150
151 // The spline has to be recomputed, since the underlying data has changed
152 _bspline.reset();
153}
154
155const Sample &Storage::getSampleAtOrAfter(double before) const
156{
157 PRECICE_TRACE(before);
158 if (nTimes() == 1) {
159 return _stampleStorage.front().sample; // @todo in this case the name getSampleAtOrAfter does not fit, because _stampleStorage.front().sample is returned for any time before.
160 } else {
161 auto stample = std::find_if(_stampleStorage.begin(), _stampleStorage.end(), [&before](const auto &s) { return math::greaterEquals(s.timestamp, before); });
162 PRECICE_ASSERT(stample != _stampleStorage.end(), "no values found!");
163 return stample->sample;
164 }
165}
166
167Eigen::VectorXd Storage::getTimes() const
168{
169 auto times = Eigen::VectorXd(nTimes());
170 for (int i = 0; i < times.size(); i++) {
171 times[i] = _stampleStorage[i].timestamp;
172 }
173 return times;
174}
175
176bool Storage::empty() const
177{
178 return _stampleStorage.empty();
179}
180
182{
184 return _stampleStorage[_stampleStorage.size() - 1];
185}
186
188{
189 auto times = Eigen::VectorXd(nTimes());
190 auto values = Eigen::MatrixXd(nDofs(), nTimes());
191 for (int i = 0; i < times.size(); i++) {
192 times[i] = _stampleStorage[i].timestamp;
193 values.col(i) = _stampleStorage[i].sample.values;
194 }
195 return std::make_pair(times, values);
196}
197
199{
200 PRECICE_ASSERT(this->nTimes() != 0, "There are no samples available");
201 const int usedDegree = computeUsedDegree(_degree, nTimes());
202
203 if (usedDegree == 0) {
204 return this->getSampleAtOrAfter(time).values;
205 }
206
207 PRECICE_ASSERT(usedDegree >= 1);
208
209 // Find existing samples
210 for (const auto &stample : _stampleStorage) {
211 if (math::equals(stample.timestamp, time)) {
212 return stample.sample.values;
213 }
214 if (math::greater(stample.timestamp, time)) {
215 break;
216 }
217 }
218
219 // Create a new bspline if _bspline does not already contain a spline
220 if (!_bspline.has_value()) {
221 auto [times, values] = getTimesAndValues();
222 _bspline.emplace(times, values, usedDegree);
223 }
224
225 return _bspline.value().interpolateAt(time);
226}
227
228Eigen::MatrixXd Storage::sampleGradients(double time) const
229{
230 const int usedDegree = computeUsedDegree(_degree, nTimes());
231
232 if (usedDegree == 0) {
233 return this->getSampleAtOrAfter(time).gradients;
234 }
235
236 PRECICE_WARN("You specified interpolation degree of {}, but only degree 0 is supported for gradient interpolation", usedDegree); // @todo implement this like for sampleAt
237 return this->getSampleAtOrAfter(time).gradients;
238}
239
240int Storage::computeUsedDegree(int requestedDegree, int numberOfAvailableSamples) const
241{
242 return std::min(requestedDegree, numberOfAvailableSamples - 1);
243}
244
246{
248 return _stampleStorage.back().sample;
249}
250
251} // namespace precice::time
#define PRECICE_WARN(...)
Definition LogMacros.hpp:12
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
void trimAfter(double time)
Definition Storage.cpp:146
const Sample & getSampleAtEnd() const
Returns the last Sample contained in this Storage.
Definition Storage.cpp:245
SampleResult sample(double time) const
Need to use interpolation for the case with changing time grids.
Definition Storage.cpp:198
auto stamples() const
Get the stamples.
Definition Storage.hpp:93
void trimBefore(double time)
Definition Storage.cpp:137
void setSampleAtTime(double time, const Sample &sample)
Store Sample at a specific time.
Definition Storage.cpp:27
Storage()
Stores data samples in time and provides corresponding convenience functions.
Definition Storage.cpp:12
void clear()
Clears this Storage by deleting all values.
Definition Storage.cpp:117
int nDofs() const
Number of Dofs for each values.
Definition Storage.cpp:87
Eigen::MatrixXd sampleGradients(double time) const
Definition Storage.cpp:228
void trim()
Trims this Storage by deleting all values except values associated with the window start.
Definition Storage.cpp:105
std::pair< Eigen::VectorXd, Eigen::MatrixXd > getTimesAndValues() const
Get all normalized dts and values in ascending order (with respect to normalized dts)
Definition Storage.cpp:187
double maxStoredTime() const
Get maximum time that is stored in this Storage.
Definition Storage.cpp:73
std::optional< math::Bspline > _bspline
Definition Storage.hpp:167
int computeUsedDegree(int requestedDegree, int numberOfAvailableSamples) const
Computes which degree may be used for interpolation.
Definition Storage.cpp:240
const Sample & getSampleAtOrAfter(double before) const
Returns the Sample at time following "before" contained in this Storage.
Definition Storage.cpp:155
Storage & operator=(const Storage &other)
Copy assignment operator to assign Storage to this Storage.
Definition Storage.cpp:17
void setInterpolationDegree(int interpolationDegree)
Definition Storage.cpp:59
std::vector< Stample > _stampleStorage
Stores Stamples on the current window.
Definition Storage.hpp:161
void setAllSamples(const Sample &sample)
Overrides all existing samples.
Definition Storage.cpp:52
void move()
Move this Storage by deleting all stamples except the one at the end of the window.
Definition Storage.cpp:93
int nTimes() const
Number of stored times.
Definition Storage.cpp:82
Eigen::VectorXd getTimes() const
Get all normalized dts stored in this Storage sorted ascending.
Definition Storage.cpp:167
const time::Stample & last() const
Definition Storage.cpp:181
int getInterpolationDegree() const
Definition Storage.cpp:68
static const int MIN_WAVEFORM_DEGREE
The minimum required interpolation degree.
Definition Time.hpp:11
T find_if(T... args)
T make_pair(T... args)
T min(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.
std::enable_if< std::is_arithmetic< Scalar >::value, bool >::type smallerEquals(Scalar A, Scalar B, Scalar tolerance=NUMERICAL_ZERO_DIFFERENCE)
std::enable_if< std::is_arithmetic< Scalar >::value, bool >::type smaller(Scalar A, Scalar B, Scalar tolerance=NUMERICAL_ZERO_DIFFERENCE)
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 remove_if(T... args)
Eigen::MatrixXd gradients
The gradients of the data. Use gradients.col(d*i+k) to get the gradient of vertex i,...
Definition Sample.hpp:67
Eigen::VectorXd values
Definition Sample.hpp:64
Stample containing timestampled Sample.
Definition Stample.hpp:7