preCICE
Loading...
Searching...
No Matches
Waveform.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Core>
4#include <boost/range.hpp>
5#include <optional>
6#include "logging/Logger.hpp"
7#include "math/Bspline.hpp"
9#include "time/Stample.hpp"
10
11namespace precice::time {
12
13class Waveform {
14public:
23 Waveform();
24
25 Waveform(int interpolationDegree);
26
27 Waveform &operator=(const Waveform &other);
28
37 void setSampleAtTime(double time, const Sample &sample);
38
47 void setAllSamples(const Sample &sample);
48
49 void setInterpolationDegree(int interpolationDegree);
50
51 int getInterpolationDegree() const;
52
58 double maxStoredTime() const;
59
68 const Sample &getSampleAtOrAfter(double before) const;
69
75 const Sample &getSampleAtEnd() const;
76
82 Eigen::VectorXd getTimes() const;
83
89 auto stamples() const
90 {
91 return boost::make_iterator_range(_stampleStorage);
92 }
93
94 auto stamples()
95 {
96 return boost::make_iterator_range(_stampleStorage);
97 }
98
99 bool empty() const;
100
101 const time::Stample &last() const;
102
108 std::pair<Eigen::VectorXd, Eigen::MatrixXd> getTimesAndValues() const;
109
115 int nTimes() const;
116
122 int nDofs() const;
123
127 void move();
128
132 void trim();
133
137 void clear();
138
139 void clearExceptLast();
140
141 void trimBefore(double time);
142
143 void trimAfter(double time);
144
151 SampleResult sample(double time) const;
152
153 Eigen::MatrixXd sampleGradients(double time) const;
154
155private:
157 std::vector<Stample> _stampleStorage;
158
159 mutable logging::Logger _log{"time::Storage"};
160
162
163 mutable std::optional<math::Bspline> _bspline;
164
175 int computeUsedDegree(int requestedDegree, int numberOfAvailableSamples) const;
176};
177
178} // namespace precice::time
This class provides a lightweight logger.
Definition Logger.hpp:17
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
logging::Logger _log
Definition Waveform.hpp:159
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
contains the time interpolation logic.
Definition Sample.hpp:8
Stample containing timestampled Sample.
Definition Stample.hpp:7