preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Data.cpp
Go to the documentation of this file.
1#include "Data.hpp"
2#include <algorithm>
3#include <utility>
4
8#include "utils/assertion.hpp"
9
10namespace precice::mesh {
11
14 DataID id,
15 int dimensions,
16 int spatialDimensions,
17 int waveformDegree)
18 : _waveform(waveformDegree),
19 _name(std::move(name)),
20 _id(id),
21 _dimensions(dimensions),
22 _spatialDimensions(spatialDimensions),
23 _sample(_dimensions)
24{
25 PRECICE_ASSERT(dimensions > 0, dimensions);
26}
27
28Eigen::VectorXd &Data::values()
29{
30 return _sample.values;
31}
32
33const Eigen::VectorXd &Data::values() const
34{
35 return _sample.values;
36}
37
38const Eigen::MatrixXd &Data::gradients() const
39{
40 return _sample.gradients;
41}
42
44{
45 return _sample;
46}
47
49{
50 return _waveform.sample(time);
51}
52
57
62
64{
65 if (stamples().size() > 1) { // Needed to avoid CompositionalCouplingScheme callong moveToNextWindow on same Data multiple times. Could be simplified by replacing Storage::move() with clearBefore(double time). See https://github.com/precice/precice/issues/1821.
67 PRECICE_ASSERT(stamples().size() == 1);
68 setGlobalSample(stamples().back().sample); // @todo at some point we should not need this anymore, when mapping, acceleration ... directly work on _timeStepsStorage, see https://github.com/precice/precice/issues/1645
69 }
70}
71
72void Data::setSampleAtTime(double time, const time::Sample &sample)
73{
74 PRECICE_ASSERT(sample.dataDims == getDimensions(), "Sample has incorrect data dimension");
75 setGlobalSample(sample); // @todo at some point we should not need this anymore, when mapping, acceleration ... directly work on _timeStepsStorage, see https://github.com/precice/precice/issues/1645
77}
78
80{
81 PRECICE_ASSERT(not sample.values.hasNaN());
83}
84
89
91{
93 Eigen::Map<const Eigen::VectorXd>(values.begin(), values.size())});
94}
95
97{
98 PRECICE_ASSERT(gradients.size() == values.size() * getSpatialDimensions(), "Gradient isn't correctly sized", values.size(), gradients.size());
99 auto nVertices = values.size() / getDimensions();
101 Eigen::Map<const Eigen::VectorXd>(values.begin(), values.size()),
102 Eigen::Map<const Eigen::MatrixXd>(gradients.begin(), getSpatialDimensions(), nVertices * getDimensions())});
103}
104
106{
107 return _name;
108}
109
111{
112 return _id;
113}
114
116{
117 return _hasGradient;
118}
119
121{
122 return !_waveform.stamples().empty();
123}
124
126{
127 _hasGradient = true;
128};
129
131{
132 return _dimensions;
133}
134
135void Data::allocateValues(int expectedCount)
136{
137 using SizeType = std::remove_cv<decltype(expectedCount)>::type;
138 // Allocate data values
139 const SizeType expectedSize = expectedCount * _dimensions;
140 const auto actualSize = static_cast<SizeType>(_sample.values.size());
141 // Shrink Buffer
142 if (expectedSize < actualSize) {
143 _sample.values.resize(expectedSize);
144 }
145 // Enlarge Buffer
146 if (expectedSize > actualSize) {
147 const auto leftToAllocate = expectedSize - actualSize;
148 utils::append(_sample.values, Eigen::VectorXd(Eigen::VectorXd::Zero(leftToAllocate)));
149 }
150 PRECICE_DEBUG("Data {} now has {} values", _name, _sample.values.size());
151
152 // Allocate gradient data values
153 if (_hasGradient) {
154 const SizeType spaceDimensions = _spatialDimensions;
155
156 const SizeType expectedColumnSize = expectedCount * _dimensions;
157 const auto actualColumnSize = static_cast<SizeType>(_sample.gradients.cols());
158
159 // Shrink Buffer
160 if (expectedColumnSize < actualColumnSize) {
161 _sample.gradients.resize(spaceDimensions, expectedColumnSize);
162 }
163
164 // Enlarge Buffer
165 if (expectedColumnSize > actualColumnSize) {
166 const auto columnLeftToAllocate = expectedColumnSize - actualColumnSize;
167 utils::append(_sample.gradients, Eigen::MatrixXd(Eigen::MatrixXd::Zero(spaceDimensions, columnLeftToAllocate)));
168 }
169 PRECICE_DEBUG("Gradient Data {} now has {} x {} values", _name, _sample.gradients.rows(), _sample.gradients.cols());
170 }
171}
172
174{
175 return _spatialDimensions;
176}
177
178} // namespace precice::mesh
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:61
std::string name
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
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
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
void setSampleAtTime(double time, const Sample &sample)
Store Sample at a specific time.
Definition Storage.cpp:27
void move()
Move this Storage by deleting all stamples except the one at the end of the window.
Definition Storage.cpp:93
int getInterpolationDegree() const
Definition Storage.cpp:68
auto stamples() const
Returns a the stamples from _timeStepsStorage.
Definition Waveform.hpp:42
time::Storage & timeStepsStorage()
Returns a reference to the _timeStepsStorage.
Definition Waveform.cpp:16
SampleResult sample(const double time) const
Evaluate waveform at specific point in time. Uses interpolation if necessary.
Definition Waveform.cpp:26
provides Mesh, Data and primitives.
void append(Eigen::VectorXd &v, double value)
int DataID
Definition Types.hpp:25
STL namespace.
Eigen::MatrixXd gradients
The gradients of the data. Use gradients.col(d*i+k) to get the gradient of vertex i,...
Definition Sample.hpp:67
int dataDims
The dimensionality of the data.
Definition Sample.hpp:60
Eigen::VectorXd values
Definition Sample.hpp:64