preCICE v3.1.1
Loading...
Searching...
No Matches
Data.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Core>
4#include <stddef.h>
5#include <string>
6
7#include "SharedPointer.hpp"
8#include "logging/Logger.hpp"
10#include "time/Sample.hpp"
11#include "time/Storage.hpp"
12#include "time/Time.hpp"
13#include "time/Waveform.hpp"
14
15namespace precice {
16namespace mesh {
17class Mesh;
18}
19} // namespace precice
20
21// ----------------------------------------------------------- CLASS DEFINITION
22
23namespace precice {
24namespace mesh {
25
29class Data {
30public:
31 // @brief Data dimensions type (scalar/vector)
36
37 // @brief Possible types of data values.
38 // enum DataType {
39 // TYPE_UNDEFINED,
40 // TYPE_DOUBLE,
41 // TYPE_VECTOR
42 // };
43
44 // @brief Name of an undefined data type.
45 //static const std::string TYPE_NAME_UNDEFINED;
46 // @brief Name of a double data type.
47 //static const std::string TYPE_NAME_DOUBLE;
48 // @brief Name of a vector data type.
49 //static const std::string TYPE_NAME_VECTOR;
50
54 Data(
56 DataID id,
57 int dimension,
58 int spatialDimensions = -1,
59 int waveformDegree = time::Time::DEFAULT_WAVEFORM_DEGREE);
60
62 Eigen::VectorXd &values();
63
65 const Eigen::VectorXd &values() const;
66
68 Eigen::MatrixXd &gradients();
69
71 const Eigen::MatrixXd &gradients() const;
72
75
77 const time::Sample &sample() const;
78
85 Eigen::VectorXd sampleAtTime(double time) const;
86
92 int getWaveformDegree() const;
93
96
97 void moveToNextWindow();
98
100 auto stamples() const
101 {
102 return _waveform.stamples();
103 }
104
106 void setSampleAtTime(double time, const time::Sample &sample);
107
109 const std::string &getName() const;
110
112 DataID getID() const;
113
115 void toZero();
116
118 bool hasGradient() const;
119
121 void requireDataGradient();
122
124 int getSpatialDimensions() const;
125
127 int getDimensions() const;
128
134 void allocateValues(int expectedCount);
135
136private:
137 logging::Logger _log{"mesh::Data"};
138
141
144
147
150
153
155 bool _hasGradient = false;
156
158};
159
160} // namespace mesh
161} // namespace precice
std::string name
This class provides a lightweight logger.
Definition Logger.hpp:16
Describes a set of data values belonging to the vertices of a mesh.
Definition Data.hpp:29
time::Sample & sample()
Returns a reference to the _sample.
Definition Data.cpp:48
const std::string & getName() const
Returns the name of the data set, as set in the config file.
Definition Data.cpp:88
std::string _name
Name of the data set.
Definition Data.hpp:143
time::Storage & timeStepsStorage()
Returns a reference to the _timeStepsStorage of _waveform.
Definition Data.cpp:68
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:116
void requireDataGradient()
Set the additional requirement of gradient data.
Definition Data.cpp:111
bool hasGradient() const
Returns if the data contains gradient data.
Definition Data.cpp:106
logging::Logger _log
Definition Data.hpp:137
int _spatialDimensions
Spatial Dimension of one element -> number of rows (only 2, 3 allowed for 2D, 3D).
Definition Data.hpp:152
auto stamples() const
Returns a the stamples from _timeStepsStorage.
Definition Data.hpp:100
void allocateValues(int expectedCount)
Allocates memory for the data values and corresponding gradient values.
Definition Data.cpp:121
int getSpatialDimensions() const
Returns the mesh dimension (i.e., number of rows) of one gradient data value .
Definition Data.cpp:159
int _dimensions
Dimensionality of one data value.
Definition Data.hpp:149
int getWaveformDegree() const
get degree of _waveform.
Definition Data.cpp:63
Eigen::VectorXd sampleAtTime(double time) const
Samples _waveform at given time.
Definition Data.cpp:58
void toZero()
Sets all values to zero.
Definition Data.cpp:98
Eigen::MatrixXd & gradients()
Returns a reference to the gradient data values.
Definition Data.cpp:38
DataID getID() const
Returns the ID of the data set (supposed to be unique).
Definition Data.cpp:93
DataID _id
ID of the data set (supposed to be unique).
Definition Data.hpp:146
void moveToNextWindow()
Definition Data.cpp:73
bool _hasGradient
Whether gradient data is available or not.
Definition Data.hpp:155
Eigen::VectorXd & values()
Returns a reference to the data values.
Definition Data.cpp:28
time::Waveform _waveform
Waveform wrapping this Data.
Definition Data.hpp:140
time::Sample _sample
Definition Data.hpp:157
void setSampleAtTime(double time, const time::Sample &sample)
Add sample at given time to _timeStepsStorage.
Definition Data.cpp:82
static const int DEFAULT_WAVEFORM_DEGREE
To be used, when the interpolation degree is not defined.
Definition Time.hpp:9
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
Main namespace of the precice library.
int DataID
Definition Types.hpp:25