preCICE
Loading...
Searching...
No Matches
src
time
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
"
8
#include "
time/SampleResult.hpp
"
9
#include "
time/Stample.hpp
"
10
11
namespace
precice::time
{
12
13
class
Waveform
{
14
public
:
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
155
private
:
157
std::vector<Stample>
_stampleStorage
;
158
159
mutable
logging::Logger
_log
{
"time::Storage"
};
160
161
int
_degree
;
162
163
mutable
std::optional<math::Bspline>
_bspline
;
164
175
int
computeUsedDegree
(
int
requestedDegree,
int
numberOfAvailableSamples)
const
;
176
};
177
178
}
// namespace precice::time
Bspline.hpp
Logger.hpp
SampleResult.hpp
Stample.hpp
precice::logging::Logger
This class provides a lightweight logger.
Definition
Logger.hpp:17
precice::time::SampleResult
Definition
SampleResult.hpp:14
precice::time::Waveform::stamples
auto stamples() const
Get the stamples.
Definition
Waveform.hpp:89
precice::time::Waveform::empty
bool empty() const
Definition
Waveform.cpp:178
precice::time::Waveform::nTimes
int nTimes() const
Number of stored times.
Definition
Waveform.cpp:84
precice::time::Waveform::getTimes
Eigen::VectorXd getTimes() const
Get all normalized dts stored in this Storage sorted ascending.
Definition
Waveform.cpp:169
precice::time::Waveform::stamples
auto stamples()
Definition
Waveform.hpp:94
precice::time::Waveform::_stampleStorage
std::vector< Stample > _stampleStorage
Stores Stamples on the current window.
Definition
Waveform.hpp:157
precice::time::Waveform::clear
void clear()
Clears this Storage by deleting all values.
Definition
Waveform.cpp:119
precice::time::Waveform::setAllSamples
void setAllSamples(const Sample &sample)
Overrides all existing samples.
Definition
Waveform.cpp:54
precice::time::Waveform::trim
void trim()
Trims this Storage by deleting all values except values associated with the window start.
Definition
Waveform.cpp:107
precice::time::Waveform::clearExceptLast
void clearExceptLast()
Definition
Waveform.cpp:128
precice::time::Waveform::Waveform
Waveform()
Stores data samples in time and provides sampling of the Waveform.
Definition
Waveform.cpp:12
precice::time::Waveform::trimBefore
void trimBefore(double time)
Definition
Waveform.cpp:139
precice::time::Waveform::setSampleAtTime
void setSampleAtTime(double time, const Sample &sample)
Store Sample at a specific time.
Definition
Waveform.cpp:29
precice::time::Waveform::getTimesAndValues
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
precice::time::Waveform::setInterpolationDegree
void setInterpolationDegree(int interpolationDegree)
Definition
Waveform.cpp:61
precice::time::Waveform::move
void move()
Move this Storage by deleting all stamples except the one at the end of the window.
Definition
Waveform.cpp:95
precice::time::Waveform::getInterpolationDegree
int getInterpolationDegree() const
Definition
Waveform.cpp:70
precice::time::Waveform::last
const time::Stample & last() const
Definition
Waveform.cpp:183
precice::time::Waveform::sampleGradients
Eigen::MatrixXd sampleGradients(double time) const
Definition
Waveform.cpp:230
precice::time::Waveform::getSampleAtEnd
const Sample & getSampleAtEnd() const
Returns the last Sample contained in this Storage.
Definition
Waveform.cpp:247
precice::time::Waveform::operator=
Waveform & operator=(const Waveform &other)
Definition
Waveform.cpp:19
precice::time::Waveform::computeUsedDegree
int computeUsedDegree(int requestedDegree, int numberOfAvailableSamples) const
Computes which degree may be used for interpolation.
Definition
Waveform.cpp:242
precice::time::Waveform::_degree
int _degree
Definition
Waveform.hpp:161
precice::time::Waveform::_log
logging::Logger _log
Definition
Waveform.hpp:159
precice::time::Waveform::trimAfter
void trimAfter(double time)
Definition
Waveform.cpp:148
precice::time::Waveform::maxStoredTime
double maxStoredTime() const
Get maximum time that is stored in this Storage.
Definition
Waveform.cpp:75
precice::time::Waveform::sample
SampleResult sample(double time) const
Need to use interpolation for the case with changing time grids.
Definition
Waveform.cpp:200
precice::time::Waveform::getSampleAtOrAfter
const Sample & getSampleAtOrAfter(double before) const
Returns the Sample at time following "before" contained in this Storage.
Definition
Waveform.cpp:157
precice::time::Waveform::_bspline
std::optional< math::Bspline > _bspline
Definition
Waveform.hpp:163
precice::time::Waveform::nDofs
int nDofs() const
Number of Dofs for each values.
Definition
Waveform.cpp:89
precice::time
contains the time interpolation logic.
Definition
Sample.hpp:8
precice::time::Sample
Definition
Sample.hpp:15
precice::time::Stample
Stample containing timestampled Sample.
Definition
Stample.hpp:7