preCICE v3.1.2
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"
8#include "time/Stample.hpp"
9
10namespace precice::time {
11
12class Storage {
13public:
22 Storage();
23
30 Storage &operator=(const Storage &other);
31
40 void setSampleAtTime(double time, const Sample &sample);
41
42 void setInterpolationDegree(int interpolationDegree);
43
44 int getInterpolationDegree() const;
45
51 double maxStoredTime() const;
52
61 Sample getSampleAtOrAfter(double before) const;
62
68 Eigen::VectorXd getTimes() const;
69
75 auto stamples() const
76 {
77 return boost::make_iterator_range(_stampleStorage);
78 }
79
80 auto stamples()
81 {
82 return boost::make_iterator_range(_stampleStorage);
83 }
84
85 bool empty() const;
86
87 const time::Stample &last() const;
88
95
101 int nTimes() const;
102
108 int nDofs() const;
109
113 void move();
114
118 void trim();
119
123 void clear();
124
125 void clearExceptLast();
126
127 void trimBefore(double time);
128
129 void trimAfter(double time);
130
137 Eigen::VectorXd sample(double time) const;
138
139 Eigen::MatrixXd sampleGradients(double time) const;
140
141private:
144
145 mutable logging::Logger _log{"time::Storage"};
146
148
150
161 int computeUsedDegree(int requestedDegree, int numberOfAvailableSamples) const;
162
164
166
167 int findTimeId(double time) const;
168};
169
170} // namespace precice::time
This class provides a lightweight logger.
Definition Logger.hpp:16
void trimAfter(double time)
Definition Storage.cpp:139
Eigen::VectorXd sample(double time) const
Need to use interpolation for the case with changing time grids.
Definition Storage.cpp:191
Sample getSampleAtOrAfter(double before) const
Returns the Sample at time following "before" contained in this Storage.
Definition Storage.cpp:148
auto stamples() const
Get the stamples.
Definition Storage.hpp:75
void trimBefore(double time)
Definition Storage.cpp:130
int findTimeId(double time) const
Definition Storage.cpp:244
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:110
int nDofs() const
Number of Dofs for each values.
Definition Storage.cpp:80
Eigen::MatrixXd sampleGradients(double time) const
Definition Storage.cpp:216
void trim()
Trims this Storage by deleting all values except values associated with the window start.
Definition Storage.cpp:98
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:180
double maxStoredTime() const
Get maximum time that is stored in this Storage.
Definition Storage.cpp:66
logging::Logger _log
Definition Storage.hpp:145
std::optional< math::Bspline > _bspline
Definition Storage.hpp:149
int computeUsedDegree(int requestedDegree, int numberOfAvailableSamples) const
Computes which degree may be used for interpolation.
Definition Storage.cpp:228
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:52
std::vector< Stample > _stampleStorage
Stores Stamples on the current window.
Definition Storage.hpp:143
void move()
Move this Storage by deleting all stamples except the one at the end of the window.
Definition Storage.cpp:86
int nTimes() const
Number of stored times.
Definition Storage.cpp:75
Eigen::VectorXd getTimes() const
Get all normalized dts stored in this Storage sorted ascending.
Definition Storage.cpp:160
time::Sample getSampleAtBeginning()
Definition Storage.cpp:234
const time::Stample & last() const
Definition Storage.cpp:174
time::Sample getSampleAtEnd()
Definition Storage.cpp:239
int getInterpolationDegree() const
Definition Storage.cpp:61
contains the time interpolation logic.
Definition Sample.hpp:6
Stample containing timestampled Sample.
Definition Stample.hpp:7