preCICE
Loading...
Searching...
No Matches
precice::time::Waveform Class Reference

#include <Waveform.hpp>

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

Public Member Functions

 Waveform ()
 Stores data samples in time and provides sampling of the Waveform.
 Waveform (int interpolationDegree)
Waveformoperator= (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 SamplegetSampleAtOrAfter (double before) const
 Returns the Sample at time following "before" contained in this Storage.
const SamplegetSampleAtEnd () 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::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)
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

Detailed Description

Definition at line 13 of file Waveform.hpp.

Constructor & Destructor Documentation

◆ Waveform() [1/2]

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.

◆ Waveform() [2/2]

precice::time::Waveform::Waveform ( int interpolationDegree)

Definition at line 17 of file Waveform.cpp.

Member Function Documentation

◆ clear()

void precice::time::Waveform::clear ( )

Clears this Storage by deleting all values.

Definition at line 119 of file Waveform.cpp.

◆ clearExceptLast()

void precice::time::Waveform::clearExceptLast ( )

Definition at line 128 of file Waveform.cpp.

◆ computeUsedDegree()

int precice::time::Waveform::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 242 of file Waveform.cpp.

◆ empty()

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

Definition at line 178 of file Waveform.cpp.

◆ getInterpolationDegree()

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

Definition at line 70 of file Waveform.cpp.

◆ getSampleAtEnd()

const Sample & precice::time::Waveform::getSampleAtEnd ( ) const

Returns the last Sample contained in this Storage.

Returns
Last Sample in this Storage

Definition at line 247 of file Waveform.cpp.

◆ getSampleAtOrAfter()

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"

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 157 of file Waveform.cpp.

Here is the call graph for this function:

◆ getTimes()

Eigen::VectorXd precice::time::Waveform::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 169 of file Waveform.cpp.

Here is the call graph for this function:

◆ getTimesAndValues()

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)

Returns
std::pair<Eigen::VectorXd, Eigen::MatrixXd> containing all stored times and values in ascending order (with respect to normalized dts).

Definition at line 189 of file Waveform.cpp.

Here is the call graph for this function:

◆ last()

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

Definition at line 183 of file Waveform.cpp.

◆ maxStoredTime()

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

Get maximum time that is stored in this Storage.

Returns
the maximum time from this Storage

Definition at line 75 of file Waveform.cpp.

◆ move()

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.

Here is the call graph for this function:

◆ nDofs()

int precice::time::Waveform::nDofs ( ) const

Number of Dofs for each values.

Returns
int number of dofs

Definition at line 89 of file Waveform.cpp.

◆ nTimes()

int precice::time::Waveform::nTimes ( ) const

Number of stored times.

Returns
int number of stored times

Definition at line 84 of file Waveform.cpp.

◆ operator=()

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

Definition at line 19 of file Waveform.cpp.

Here is the call graph for this function:

◆ sample()

SampleResult precice::time::Waveform::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 200 of file Waveform.cpp.

Here is the call graph for this function:

◆ sampleGradients()

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

Definition at line 230 of file Waveform.cpp.

Here is the call graph for this function:

◆ setAllSamples()

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

Parameters
samplenew sample value

Definition at line 54 of file Waveform.cpp.

Here is the call graph for this function:

◆ setInterpolationDegree()

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

Definition at line 61 of file Waveform.cpp.

◆ setSampleAtTime()

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.

Parameters
timethe time associated with the sample
samplestored sample

Definition at line 29 of file Waveform.cpp.

Here is the call graph for this function:

◆ stamples() [1/2]

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

Definition at line 94 of file Waveform.hpp.

◆ stamples() [2/2]

auto precice::time::Waveform::stamples ( ) const
inline

Get the stamples.

Returns
boost range of stamples

Definition at line 89 of file Waveform.hpp.

◆ trim()

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.

◆ trimAfter()

void precice::time::Waveform::trimAfter ( double time)

Definition at line 148 of file Waveform.cpp.

Here is the call graph for this function:

◆ trimBefore()

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

Definition at line 139 of file Waveform.cpp.

Here is the call graph for this function:

Member Data Documentation

◆ _bspline

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

Definition at line 163 of file Waveform.hpp.

◆ _degree

int precice::time::Waveform::_degree
private

Definition at line 161 of file Waveform.hpp.

◆ _log

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

Definition at line 159 of file Waveform.hpp.

◆ _stampleStorage

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

Stores Stamples on the current window.

Definition at line 157 of file Waveform.hpp.


The documentation for this class was generated from the following files: