preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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/Storage.hpp"
13#include "time/Time.hpp"
14#include "time/Waveform.hpp"
15
16namespace precice::mesh {
17class Mesh;
18}
19
20// ----------------------------------------------------------- CLASS DEFINITION
21
22namespace precice::mesh {
23
27class Data {
28public:
29 // @brief Data dimensions type (scalar/vector)
34
35 // @brief Possible types of data values.
36 // enum DataType {
37 // TYPE_UNDEFINED,
38 // TYPE_DOUBLE,
39 // TYPE_VECTOR
40 // };
41
42 // @brief Name of an undefined data type.
43 // static const std::string TYPE_NAME_UNDEFINED;
44 // @brief Name of a double data type.
45 // static const std::string TYPE_NAME_DOUBLE;
46 // @brief Name of a vector data type.
47 // static const std::string TYPE_NAME_VECTOR;
48
52 Data(
54 DataID id,
55 int dimension,
56 int spatialDimensions = -1,
57 int waveformDegree = time::Time::DEFAULT_WAVEFORM_DEGREE);
58
60 Eigen::VectorXd &values();
61
63 const Eigen::VectorXd &values() const;
64
66 const Eigen::MatrixXd &gradients() const;
67
69 const time::Sample &sample() const;
70
77 time::SampleResult sampleAtTime(double time) const;
78
84 int getWaveformDegree() const;
85
88
89 void moveToNextWindow();
90
92 auto stamples() const
93 {
94 return _waveform.stamples();
95 }
96
98 void setSampleAtTime(double time, const time::Sample &sample);
99
101 void setGlobalSample(const time::Sample &sample); // @todo try to remove this function
102
104 void emplaceSampleAtTime(double time);
105
108
111
113 const std::string &getName() const;
114
116 DataID getID() const;
117
119 bool hasGradient() const;
120
122 bool hasSamples() const;
123
125 void requireDataGradient();
126
128 int getSpatialDimensions() const;
129
131 int getDimensions() const;
132
138 void allocateValues(int expectedCount);
139
140private:
141 logging::Logger _log{"mesh::Data"};
142
145
148
151
154
157
159 bool _hasGradient = false;
160
162};
163
164} // namespace precice::mesh
std::string name
This class provides a lightweight logger.
Definition Logger.hpp:17
Describes a set of data values belonging to the vertices of a mesh.
Definition Data.hpp:27
const std::string & getName() const
Returns the name of the data set, as set in the config file.
Definition Data.cpp:105
std::string _name
Name of the data set.
Definition Data.hpp:147
time::Storage & timeStepsStorage()
Returns a reference to the _timeStepsStorage of _waveform.
Definition Data.cpp:58
bool hasSamples() const
Returns if there are sample of this data.
Definition Data.cpp:120
Data(std::string name, DataID id, int dimension, int spatialDimensions=-1, int waveformDegree=time::Time::DEFAULT_WAVEFORM_DEGREE)
Constructor.
Definition Data.cpp:12
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:130
void requireDataGradient()
Set the additional requirement of gradient data.
Definition Data.cpp:125
bool hasGradient() const
Returns if the data contains gradient data.
Definition Data.cpp:115
logging::Logger _log
Definition Data.hpp:141
int _spatialDimensions
Spatial Dimension of one element -> number of rows (only 2, 3 allowed for 2D, 3D).
Definition Data.hpp:156
auto stamples() const
Returns a the stamples from _timeStepsStorage.
Definition Data.hpp:92
void allocateValues(int expectedCount)
Allocates memory for the data values and corresponding gradient values.
Definition Data.cpp:135
int getSpatialDimensions() const
Returns the mesh dimension (i.e., number of rows) of one gradient data value .
Definition Data.cpp:173
int _dimensions
Dimensionality of one data value.
Definition Data.hpp:153
int getWaveformDegree() const
get degree of _waveform.
Definition Data.cpp:53
void setGlobalSample(const time::Sample &sample)
Set _sample.
Definition Data.cpp:79
time::SampleResult sampleAtTime(double time) const
Samples _waveform at given time.
Definition Data.cpp:48
DataID getID() const
Returns the ID of the data set (supposed to be unique).
Definition Data.cpp:110
DataID _id
ID of the data set (supposed to be unique).
Definition Data.hpp:150
void moveToNextWindow()
Definition Data.cpp:63
bool _hasGradient
Whether gradient data is available or not.
Definition Data.hpp:159
Eigen::VectorXd & values()
Returns a reference to the data values.
Definition Data.cpp:28
const time::Sample & sample() const
Returns a const reference to the _sample.
Definition Data.cpp:43
time::Waveform _waveform
Waveform wrapping this Data.
Definition Data.hpp:144
time::Sample _sample
Definition Data.hpp:161
const Eigen::MatrixXd & gradients() const
Returns a const reference to the gradient data values.
Definition Data.cpp:38
void setSampleAtTime(double time, const time::Sample &sample)
Add sample at given time to _timeStepsStorage.
Definition Data.cpp:72
void emplaceSampleAtTime(double time)
Creates an empty sample at given time.
Definition Data.cpp:85
static const int DEFAULT_WAVEFORM_DEGREE
To be used, when the interpolation degree is not defined.
Definition Time.hpp:8
Allows to perform interpolation on samples in storage of given data.
Definition Waveform.hpp:25
auto stamples() const
Returns a the stamples from _timeStepsStorage.
Definition Waveform.hpp:42
provides Mesh, Data and primitives.
int DataID
Definition Types.hpp:25