preCICE
Loading...
Searching...
No Matches
Data.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Core>
4#include <initializer_list>
5#include <stddef.h>
6#include <string>
7
8#include "SharedPointer.hpp"
9#include "logging/Logger.hpp"
11#include "time/Sample.hpp"
12#include "time/Time.hpp"
13#include "time/Waveform.hpp"
14
15namespace precice::mesh {
16class Mesh;
17}
18
19// ----------------------------------------------------------- CLASS DEFINITION
20
21namespace precice::mesh {
22
26class Data {
27public:
28 // @brief Data dimensions type (scalar/vector)
33
34 // @brief Possible types of data values.
35 // enum DataType {
36 // TYPE_UNDEFINED,
37 // TYPE_DOUBLE,
38 // TYPE_VECTOR
39 // };
40
41 // @brief Name of an undefined data type.
42 // static const std::string TYPE_NAME_UNDEFINED;
43 // @brief Name of a double data type.
44 // static const std::string TYPE_NAME_DOUBLE;
45 // @brief Name of a vector data type.
46 // static const std::string TYPE_NAME_VECTOR;
47
51 Data(
52 std::string name,
53 DataID id,
54 int dimension,
55 int spatialDimensions = -1,
56 int waveformDegree = time::Time::DEFAULT_WAVEFORM_DEGREE);
57
61 Data(
62 std::string name,
63 DataID id,
64 int dimension,
65 int spatialDimensions,
66 int waveformDegree,
67 std::vector<std::optional<double>>
68 lowerBound,
69 std::vector<std::optional<double>>
70 upperBound);
71
73 Eigen::VectorXd &values();
74
76 const Eigen::VectorXd &values() const;
77
79 const Eigen::MatrixXd &gradients() const;
80
82 const time::Sample &sample() const;
83
91
97 int getWaveformDegree() const;
98
101
102 std::vector<std::optional<double>> getLowerBound() const;
103
104 std::vector<std::optional<double>> getUpperBound() const;
105
106 void moveToNextWindow();
107
109 auto stamples() const
110 {
111 return _waveform.stamples();
112 }
113
115 void setSampleAtTime(double time, const time::Sample &sample);
116
118 void setGlobalSample(const time::Sample &sample); // @todo try to remove this function
119
121 void emplaceSampleAtTime(double time);
122
124 void emplaceSampleAtTime(double time, std::initializer_list<double> values);
125
127 void emplaceSampleAtTime(double time, std::initializer_list<double> values, std::initializer_list<double> gradients);
128
130 const std::string &getName() const;
131
133 DataID getID() const;
134
136 bool hasGradient() const;
137
139 bool hasSamples() const;
140
142 void requireDataGradient();
143
145 int getSpatialDimensions() const;
146
148 int getDimensions() const;
149
155 void allocateValues(int expectedCount);
156
157private:
158 logging::Logger _log{"mesh::Data"};
159
162
164 std::vector<std::optional<double>> _lowerBound;
165
167 std::vector<std::optional<double>> _upperBound;
168
170 std::string _name;
171
174
177
180
182 bool _hasGradient = false;
183
185};
186
187} // namespace precice::mesh
This class provides a lightweight logger.
Definition Logger.hpp:17
const std::string & getName() const
Returns the name of the data set, as set in the config file.
Definition Data.cpp:139
std::string _name
Name of the data set.
Definition Data.hpp:170
bool hasSamples() const
Returns if there are sample of this data.
Definition Data.cpp:154
Data(std::string name, DataID id, int dimension, int spatialDimensions=-1, int waveformDegree=time::Time::DEFAULT_WAVEFORM_DEGREE)
Constructor.
Definition Data.cpp:12
std::vector< std::optional< double > > getUpperBound() const
Definition Data.cpp:92
int getDimensions() const
Returns the dimension (i.e., number of components) of one data value (i.e number of columns of one gr...
Definition Data.cpp:164
void requireDataGradient()
Set the additional requirement of gradient data.
Definition Data.cpp:159
bool hasGradient() const
Returns if the data contains gradient data.
Definition Data.cpp:149
logging::Logger _log
Definition Data.hpp:158
int _spatialDimensions
Spatial Dimension of one element -> number of rows (only 2, 3 allowed for 2D, 3D).
Definition Data.hpp:179
auto stamples() const
Returns a the stamples from the waveform.
Definition Data.hpp:109
void allocateValues(int expectedCount)
Allocates memory for the data values and corresponding gradient values.
Definition Data.cpp:169
int getSpatialDimensions() const
Returns the mesh dimension (i.e., number of rows) of one gradient data value .
Definition Data.cpp:207
int _dimensions
Dimensionality of one data value.
Definition Data.hpp:176
int getWaveformDegree() const
get degree of _waveform.
Definition Data.cpp:77
void setGlobalSample(const time::Sample &sample)
Set _sample.
Definition Data.cpp:113
time::SampleResult sampleAtTime(double time) const
Samples _waveform at given time.
Definition Data.cpp:72
std::vector< std::optional< double > > getLowerBound() const
Definition Data.cpp:87
DataID getID() const
Returns the ID of the data set (supposed to be unique).
Definition Data.cpp:144
DataID _id
ID of the data set (supposed to be unique).
Definition Data.hpp:173
void moveToNextWindow()
Definition Data.cpp:97
bool _hasGradient
Whether gradient data is available or not.
Definition Data.hpp:182
Eigen::VectorXd & values()
Returns a reference to the data values.
Definition Data.cpp:52
std::vector< std::optional< double > > _lowerBound
Lower bound for data values. This vector with optional elements has size of the data dimension....
Definition Data.hpp:164
const time::Sample & sample() const
Returns a const reference to the _sample.
Definition Data.cpp:67
time::Waveform _waveform
Sample storage of this Data.
Definition Data.hpp:161
std::vector< std::optional< double > > _upperBound
Upper bound for data values. This vector with optional elements has size of the data dimension....
Definition Data.hpp:167
time::Sample _sample
Definition Data.hpp:184
time::Waveform & waveform()
Returns a reference to the waveform.
Definition Data.cpp:82
const Eigen::MatrixXd & gradients() const
Returns a const reference to the gradient data values.
Definition Data.cpp:62
void setSampleAtTime(double time, const time::Sample &sample)
Add sample at given time to the waveform.
Definition Data.cpp:106
void emplaceSampleAtTime(double time)
Creates an empty sample at given time.
Definition Data.cpp:119
Container and creator for meshes.
Definition Mesh.hpp:38
static const int DEFAULT_WAVEFORM_DEGREE
To be used, when the interpolation degree is not defined.
Definition Time.hpp:8
provides Mesh, Data and primitives.
contains the time interpolation logic.
Definition Sample.hpp:8
int DataID
Definition Types.hpp:25