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