#include <Storage.hpp>
Definition at line 12 of file Storage.hpp.
◆ Storage()
precice::time::Storage::Storage |
( |
| ) |
|
Stores data samples in time and provides corresponding convenience functions.
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 Storage.cpp.
◆ clear()
void precice::time::Storage::clear |
( |
| ) |
|
◆ clearExceptLast()
void precice::time::Storage::clearExceptLast |
( |
| ) |
|
◆ computeUsedDegree()
int precice::time::Storage::computeUsedDegree |
( |
int | requestedDegree, |
|
|
int | numberOfAvailableSamples ) const |
|
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.
- Parameters
-
requestedDegree | B-spline degree requested by the user. |
numberOfAvailableSamples | Samples available for interpolation. |
- Returns
- B-spline degree that may be used.
Definition at line 228 of file Storage.cpp.
◆ empty()
bool precice::time::Storage::empty |
( |
| ) |
const |
◆ findTimeId()
int precice::time::Storage::findTimeId |
( |
double | time | ) |
const |
|
private |
◆ getInterpolationDegree()
int precice::time::Storage::getInterpolationDegree |
( |
| ) |
const |
◆ getSampleAtBeginning()
time::Sample precice::time::Storage::getSampleAtBeginning |
( |
| ) |
|
|
private |
◆ getSampleAtEnd()
◆ getSampleAtOrAfter()
Sample precice::time::Storage::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"
- Parameters
-
before | a double, where we want to find a normalized dt that comes directly after this one |
- Returns
- Sample in this Storage at or directly after "before"
Definition at line 148 of file Storage.cpp.
◆ getTimes()
Eigen::VectorXd precice::time::Storage::getTimes |
( |
| ) |
const |
Get all normalized dts stored in this Storage sorted ascending.
- Returns
- Eigen::VectorXd containing all stored normalized dts in ascending order.
Definition at line 160 of file Storage.cpp.
◆ getTimesAndValues()
std::pair< Eigen::VectorXd, Eigen::MatrixXd > precice::time::Storage::getTimesAndValues |
( |
| ) |
const |
Get all normalized dts and values in ascending order (with respect to normalized dts)
- Returns
- std::pair<Eigen::VectorXd, Eigen::MatrixXd> containing all stored times and values in ascending order (with respect to normalized dts).
Definition at line 180 of file Storage.cpp.
◆ last()
◆ maxStoredTime()
double precice::time::Storage::maxStoredTime |
( |
| ) |
const |
◆ move()
void precice::time::Storage::move |
( |
| ) |
|
Move this Storage by deleting all stamples except the one at the end of the window.
Definition at line 86 of file Storage.cpp.
◆ nDofs()
int precice::time::Storage::nDofs |
( |
| ) |
const |
Number of Dofs for each values.
- Returns
- int number of dofs
Definition at line 80 of file Storage.cpp.
◆ nTimes()
int precice::time::Storage::nTimes |
( |
| ) |
const |
Number of stored times.
- Returns
- int number of stored times
Definition at line 75 of file Storage.cpp.
◆ operator=()
◆ sample()
Eigen::VectorXd precice::time::Storage::sample |
( |
double | time | ) |
const |
Need to use interpolation for the case with changing time grids.
- Parameters
-
time | a double, where we want to sample the waveform |
- Returns
- Eigen::VectorXd values in this Storage at or directly after "before"
Definition at line 191 of file Storage.cpp.
◆ sampleGradients()
Eigen::MatrixXd precice::time::Storage::sampleGradients |
( |
double | time | ) |
const |
◆ setInterpolationDegree()
void precice::time::Storage::setInterpolationDegree |
( |
int | interpolationDegree | ) |
|
◆ setSampleAtTime()
void precice::time::Storage::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.
- Parameters
-
time | the time associated with the sample |
sample | stored sample |
Definition at line 27 of file Storage.cpp.
◆ stamples() [1/2]
auto precice::time::Storage::stamples |
( |
| ) |
|
|
inline |
◆ stamples() [2/2]
auto precice::time::Storage::stamples |
( |
| ) |
const |
|
inline |
Get the stamples.
- Returns
- boost range of stamples
Definition at line 75 of file Storage.hpp.
◆ trim()
void precice::time::Storage::trim |
( |
| ) |
|
Trims this Storage by deleting all values except values associated with the window start.
Definition at line 98 of file Storage.cpp.
◆ trimAfter()
void precice::time::Storage::trimAfter |
( |
double | time | ) |
|
◆ trimBefore()
void precice::time::Storage::trimBefore |
( |
double | time | ) |
|
◆ _bspline
◆ _degree
int precice::time::Storage::_degree |
|
private |
◆ _log
◆ _stampleStorage
Stores Stamples on the current window.
Definition at line 143 of file Storage.hpp.
The documentation for this class was generated from the following files: