preCICE v3.1.2
Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Private Attributes | List of all members
precice::time::Storage Class Reference

#include <Storage.hpp>

Collaboration diagram for precice::time::Storage:
[legend]

Public Member Functions

 Storage ()
 Stores data samples in time and provides corresponding convenience functions.
 
Storageoperator= (const Storage &other)
 Copy assignment operator to assign Storage to this Storage.
 
void setSampleAtTime (double time, const Sample &sample)
 Store Sample at a specific time.
 
void setInterpolationDegree (int interpolationDegree)
 
int getInterpolationDegree () const
 
double maxStoredTime () const
 Get maximum time that is stored in this Storage.
 
Sample getSampleAtOrAfter (double before) const
 Returns the Sample at time following "before" 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::Stamplelast () 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)
 
Eigen::VectorXd 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.
 
time::Sample getSampleAtBeginning ()
 
time::Sample getSampleAtEnd ()
 
int findTimeId (double time) const
 

Private Attributes

std::vector< Stample_stampleStorage
 Stores Stamples on the current window.
 
logging::Logger _log {"time::Storage"}
 
int _degree
 
std::optional< math::Bspline_bspline
 

Detailed Description

Definition at line 12 of file Storage.hpp.

Constructor & Destructor Documentation

◆ 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.

Member Function Documentation

◆ clear()

void precice::time::Storage::clear ( )

Clears this Storage by deleting all values.

Definition at line 110 of file Storage.cpp.

◆ clearExceptLast()

void precice::time::Storage::clearExceptLast ( )

Definition at line 119 of file Storage.cpp.

◆ 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
requestedDegreeB-spline degree requested by the user.
numberOfAvailableSamplesSamples available for interpolation.
Returns
B-spline degree that may be used.

Definition at line 228 of file Storage.cpp.

Here is the call graph for this function:

◆ empty()

bool precice::time::Storage::empty ( ) const

Definition at line 169 of file Storage.cpp.

◆ findTimeId()

int precice::time::Storage::findTimeId ( double time) const
private

Definition at line 244 of file Storage.cpp.

Here is the call graph for this function:

◆ getInterpolationDegree()

int precice::time::Storage::getInterpolationDegree ( ) const

Definition at line 61 of file Storage.cpp.

◆ getSampleAtBeginning()

time::Sample precice::time::Storage::getSampleAtBeginning ( )
private

Definition at line 234 of file Storage.cpp.

◆ getSampleAtEnd()

time::Sample precice::time::Storage::getSampleAtEnd ( )
private

Definition at line 239 of file Storage.cpp.

◆ 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
beforea 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.

Here is the call graph for this function:

◆ 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.

Here is the call graph for this function:

◆ 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.

Here is the call graph for this function:

◆ last()

const time::Stample & precice::time::Storage::last ( ) const

Definition at line 174 of file Storage.cpp.

◆ maxStoredTime()

double precice::time::Storage::maxStoredTime ( ) const

Get maximum time that is stored in this Storage.

Returns
the maximum time from this Storage

Definition at line 66 of file Storage.cpp.

◆ 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.

Here is the call graph for this function:

◆ 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=()

Storage & precice::time::Storage::operator= ( const Storage & other)

Copy assignment operator to assign Storage to this Storage.

Parameters
otherStorage
Returns
Storage&

Definition at line 17 of file Storage.cpp.

Here is the call graph for this function:

◆ sample()

Eigen::VectorXd precice::time::Storage::sample ( double time) const

Need to use interpolation for the case with changing time grids.

Parameters
timea 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.

Here is the call graph for this function:

◆ sampleGradients()

Eigen::MatrixXd precice::time::Storage::sampleGradients ( double time) const

Definition at line 216 of file Storage.cpp.

Here is the call graph for this function:

◆ setInterpolationDegree()

void precice::time::Storage::setInterpolationDegree ( int interpolationDegree)

Definition at line 52 of file Storage.cpp.

◆ 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
timethe time associated with the sample
samplestored sample

Definition at line 27 of file Storage.cpp.

Here is the call graph for this function:

◆ stamples() [1/2]

auto precice::time::Storage::stamples ( )
inline

Definition at line 80 of file Storage.hpp.

◆ 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)

Definition at line 139 of file Storage.cpp.

Here is the call graph for this function:

◆ trimBefore()

void precice::time::Storage::trimBefore ( double time)

Definition at line 130 of file Storage.cpp.

Here is the call graph for this function:

Member Data Documentation

◆ _bspline

std::optional<math::Bspline> precice::time::Storage::_bspline
mutableprivate

Definition at line 149 of file Storage.hpp.

◆ _degree

int precice::time::Storage::_degree
private

Definition at line 147 of file Storage.hpp.

◆ _log

logging::Logger precice::time::Storage::_log {"time::Storage"}
mutableprivate

Definition at line 145 of file Storage.hpp.

◆ _stampleStorage

std::vector<Stample> precice::time::Storage::_stampleStorage
private

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: