|
preCICE
|
#include <Waveform.hpp>
Public Member Functions | |
| Waveform () | |
| Stores data samples in time and provides sampling of the Waveform. | |
| Waveform (int interpolationDegree) | |
| Waveform & | operator= (const Waveform &other) |
| void | setSampleAtTime (double time, const Sample &sample) |
| Store Sample at a specific time. | |
| void | setAllSamples (const Sample &sample) |
| Overrides all existing samples. | |
| void | setInterpolationDegree (int interpolationDegree) |
| int | getInterpolationDegree () const |
| double | maxStoredTime () const |
| Get maximum time that is stored in this Storage. | |
| const Sample & | getSampleAtOrAfter (double before) const |
| Returns the Sample at time following "before" contained in this Storage. | |
| const Sample & | getSampleAtEnd () const |
| Returns the last Sample contained in this Storage. | |
| Eigen::VectorXd | getTimes () const |
| Get all normalized dts stored in this Storage sorted ascending. | |
| auto | stamples () const |
| Get the stamples. | |
| auto | stamples () |
| bool | empty () const |
| const time::Stample & | last () const |
| std::pair< Eigen::VectorXd, Eigen::MatrixXd > | getTimesAndValues () const |
| Get all normalized dts and values in ascending order (with respect to normalized dts) | |
| int | nTimes () const |
| Number of stored times. | |
| int | nDofs () const |
| Number of Dofs for each values. | |
| void | move () |
| Move this Storage by deleting all stamples except the one at the end of the window. | |
| void | trim () |
| Trims this Storage by deleting all values except values associated with the window start. | |
| void | clear () |
| Clears this Storage by deleting all values. | |
| void | clearExceptLast () |
| void | trimBefore (double time) |
| void | trimAfter (double time) |
| SampleResult | sample (double time) const |
| Need to use interpolation for the case with changing time grids. | |
| Eigen::MatrixXd | sampleGradients (double time) const |
Private Member Functions | |
| int | computeUsedDegree (int requestedDegree, int numberOfAvailableSamples) const |
| Computes which degree may be used for interpolation. | |
Private Attributes | |
| std::vector< Stample > | _stampleStorage |
| Stores Stamples on the current window. | |
| logging::Logger | _log {"time::Storage"} |
| int | _degree |
| std::optional< math::Bspline > | _bspline |
Definition at line 13 of file Waveform.hpp.
| precice::time::Waveform::Waveform | ( | ) |
Stores data samples in time and provides sampling of the Waveform.
The Storage must be initialized before it can be used. Then values can be stored in the Storage. It is only allowed to store samples with increasing times. Overwriting existing samples or writing samples with a time smaller then the maximum stored time is forbidden. The Storage is considered complete, when a sample for the end of the current window is provided. Then one can only sample from the storage. To add further samples one needs to trim the storage or move to the next time window first.
This Storage is used in the context of Waveform relaxation where samples in time are provided.
Definition at line 12 of file Waveform.cpp.
| precice::time::Waveform::Waveform | ( | int | interpolationDegree | ) |
Definition at line 17 of file Waveform.cpp.
| void precice::time::Waveform::clear | ( | ) |
Clears this Storage by deleting all values.
Definition at line 119 of file Waveform.cpp.
| void precice::time::Waveform::clearExceptLast | ( | ) |
Definition at line 128 of file Waveform.cpp.
|
private |
Computes which degree may be used for interpolation.
Actual degree of interpolating B-spline is determined by number of stored samples and maximum degree defined by the user. Example: If only two samples are available, the maximum degree we may use is 1, even if the user demands degree 2.
| requestedDegree | B-spline degree requested by the user. |
| numberOfAvailableSamples | Samples available for interpolation. |
Definition at line 242 of file Waveform.cpp.
| bool precice::time::Waveform::empty | ( | ) | const |
Definition at line 178 of file Waveform.cpp.
| int precice::time::Waveform::getInterpolationDegree | ( | ) | const |
Definition at line 70 of file Waveform.cpp.
| const Sample & precice::time::Waveform::getSampleAtEnd | ( | ) | const |
Returns the last Sample contained in this Storage.
Definition at line 247 of file Waveform.cpp.
| const Sample & precice::time::Waveform::getSampleAtOrAfter | ( | double | before | ) | const |
Returns the Sample at time following "before" contained in this Storage.
The stored normalized dt is larger or equal than "before". If "before" is a normalized dt stored in this Storage, this function returns the Sample at "before"
| before | a double, where we want to find a normalized dt that comes directly after this one |
Definition at line 157 of file Waveform.cpp.
| Eigen::VectorXd precice::time::Waveform::getTimes | ( | ) | const |
Get all normalized dts stored in this Storage sorted ascending.
Definition at line 169 of file Waveform.cpp.
| std::pair< Eigen::VectorXd, Eigen::MatrixXd > precice::time::Waveform::getTimesAndValues | ( | ) | const |
Get all normalized dts and values in ascending order (with respect to normalized dts)
Definition at line 189 of file Waveform.cpp.
| const time::Stample & precice::time::Waveform::last | ( | ) | const |
Definition at line 183 of file Waveform.cpp.
| double precice::time::Waveform::maxStoredTime | ( | ) | const |
Get maximum time that is stored in this Storage.
Definition at line 75 of file Waveform.cpp.
| void precice::time::Waveform::move | ( | ) |
Move this Storage by deleting all stamples except the one at the end of the window.
Definition at line 95 of file Waveform.cpp.
| int precice::time::Waveform::nDofs | ( | ) | const |
Number of Dofs for each values.
Definition at line 89 of file Waveform.cpp.
| int precice::time::Waveform::nTimes | ( | ) | const |
Number of stored times.
Definition at line 84 of file Waveform.cpp.
| SampleResult precice::time::Waveform::sample | ( | double | time | ) | const |
Need to use interpolation for the case with changing time grids.
| time | a double, where we want to sample the waveform |
Definition at line 200 of file Waveform.cpp.
| Eigen::MatrixXd precice::time::Waveform::sampleGradients | ( | double | time | ) | const |
| void precice::time::Waveform::setAllSamples | ( | const Sample & | sample | ) |
Overrides all existing samples.
This keeps time stamps, but overwrites all samples with a given one. Required in reinitialization
| sample | new sample value |
Definition at line 54 of file Waveform.cpp.
| void precice::time::Waveform::setInterpolationDegree | ( | int | interpolationDegree | ) |
Definition at line 61 of file Waveform.cpp.
| void precice::time::Waveform::setSampleAtTime | ( | double | time, |
| const Sample & | sample ) |
Store Sample at a specific time.
It is only allowed to store a Sample in time that comes after a Sample that was already stored. Therefore, time has to be larger than maxStoredTime. The function trim() should be used before providing new samples.
| time | the time associated with the sample |
| sample | stored sample |
Definition at line 29 of file Waveform.cpp.
|
inline |
Definition at line 94 of file Waveform.hpp.
|
inline |
| void precice::time::Waveform::trim | ( | ) |
Trims this Storage by deleting all values except values associated with the window start.
Definition at line 107 of file Waveform.cpp.
| void precice::time::Waveform::trimAfter | ( | double | time | ) |
| void precice::time::Waveform::trimBefore | ( | double | time | ) |
|
mutableprivate |
Definition at line 163 of file Waveform.hpp.
|
private |
Definition at line 161 of file Waveform.hpp.
|
mutableprivate |
Definition at line 159 of file Waveform.hpp.
|
private |
Stores Stamples on the current window.
Definition at line 157 of file Waveform.hpp.