preCICE v3.2.0
Loading...
Searching...
No Matches
Storage.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"
9#include "time/Stample.hpp"
10
11namespace precice::time {
12
13class Storage {
14public:
23 Storage();
24
31 Storage &operator=(const Storage &other);
32
41 void setSampleAtTime(double time, const Sample &sample);
42
51 void setAllSamples(const Sample &sample);
52
53 void setInterpolationDegree(int interpolationDegree);
54
55 int getInterpolationDegree() const;
56
62 double maxStoredTime() const;
63
72 const Sample &getSampleAtOrAfter(double before) const;
73
79 const Sample &getSampleAtEnd() const;
80
86 Eigen::VectorXd getTimes() const;
87
93 auto stamples() const
94 {
95 return boost::make_iterator_range(_stampleStorage);
96 }
97
98 auto stamples()
99 {
100 return boost::make_iterator_range(_stampleStorage);
101 }
102
103 bool empty() const;
104
105 const time::Stample &last() const;
106
113
119 int nTimes() const;
120
126 int nDofs() const;
127
131 void move();
132
136 void trim();
137
141 void clear();
142
143 void clearExceptLast();
144
145 void trimBefore(double time);
146
147 void trimAfter(double time);
148
155 SampleResult sample(double time) const;
156
157 Eigen::MatrixXd sampleGradients(double time) const;
158
159private:
162
163 mutable logging::Logger _log{"time::Storage"};
164
166
168
179 int computeUsedDegree(int requestedDegree, int numberOfAvailableSamples) const;
180};
181
182} // namespace precice::time
This class provides a lightweight logger.
Definition Logger.hpp:17
void trimAfter(double time)
Definition Storage.cpp:146
const Sample & getSampleAtEnd() const
Returns the last Sample contained in this Storage.
Definition Storage.cpp:245
SampleResult sample(double time) const
Need to use interpolation for the case with changing time grids.
Definition Storage.cpp:198
auto stamples() const
Get the stamples.
Definition Storage.hpp:93
void trimBefore(double time)
Definition Storage.cpp:137
void setSampleAtTime(double time, const Sample &sample)
Store Sample at a specific time.
Definition Storage.cpp:27
Storage()
Stores data samples in time and provides corresponding convenience functions.
Definition Storage.cpp:12
void clear()
Clears this Storage by deleting all values.
Definition Storage.cpp:117
int nDofs() const
Number of Dofs for each values.
Definition Storage.cpp:87
Eigen::MatrixXd sampleGradients(double time) const
Definition Storage.cpp:228
void trim()
Trims this Storage by deleting all values except values associated with the window start.
Definition Storage.cpp:105
std::pair< Eigen::VectorXd, Eigen::MatrixXd > getTimesAndValues() const
Get all normalized dts and values in ascending order (with respect to normalized dts)
Definition Storage.cpp:187
double maxStoredTime() const
Get maximum time that is stored in this Storage.
Definition Storage.cpp:73
logging::Logger _log
Definition Storage.hpp:163
std::optional< math::Bspline > _bspline
Definition Storage.hpp:167
int computeUsedDegree(int requestedDegree, int numberOfAvailableSamples) const
Computes which degree may be used for interpolation.
Definition Storage.cpp:240
const Sample & getSampleAtOrAfter(double before) const
Returns the Sample at time following "before" contained in this Storage.
Definition Storage.cpp:155
Storage & operator=(const Storage &other)
Copy assignment operator to assign Storage to this Storage.
Definition Storage.cpp:17
void setInterpolationDegree(int interpolationDegree)
Definition Storage.cpp:59
std::vector< Stample > _stampleStorage
Stores Stamples on the current window.
Definition Storage.hpp:161
void setAllSamples(const Sample &sample)
Overrides all existing samples.
Definition Storage.cpp:52
void move()
Move this Storage by deleting all stamples except the one at the end of the window.
Definition Storage.cpp:93
int nTimes() const
Number of stored times.
Definition Storage.cpp:82
Eigen::VectorXd getTimes() const
Get all normalized dts stored in this Storage sorted ascending.
Definition Storage.cpp:167
const time::Stample & last() const
Definition Storage.cpp:181
int getInterpolationDegree() const
Definition Storage.cpp:68
contains the time interpolation logic.
Definition Sample.hpp:8
Stample containing timestampled Sample.
Definition Stample.hpp:7