preCICE
Loading...
Searching...
No Matches
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
13 std::string name,
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),
24{
25 PRECICE_ASSERT(dimensions > 0, dimensions);
26 _lowerBound = std::vector<std::optional<double>>(dimensions);
27 _upperBound = std::vector<std::optional<double>>(dimensions);
28}
29
31 std::string name,
32 DataID id,
33 int dimensions,
34 int spatialDimensions,
35 int waveformDegree,
36 std::vector<std::optional<double>> lowerBound,
37 std::vector<std::optional<double>> upperBound)
38 : _waveform(waveformDegree),
39 _lowerBound(lowerBound),
40 _upperBound(upperBound),
41 _name(std::move(name)),
42 _id(id),
43 _dimensions(dimensions),
44 _spatialDimensions(spatialDimensions),
46{
47 PRECICE_ASSERT(dimensions > 0, dimensions);
48 PRECICE_ASSERT(static_cast<int>(lowerBound.size()) == dimensions, lowerBound.size(), dimensions);
49 PRECICE_ASSERT(static_cast<int>(upperBound.size()) == dimensions, upperBound.size(), dimensions);
50}
51
52Eigen::VectorXd &Data::values()
53{
54 return _sample.values;
55}
56
57const Eigen::VectorXd &Data::values() const
58{
59 return _sample.values;
60}
61
62const Eigen::MatrixXd &Data::gradients() const
63{
64 return _sample.gradients;
65}
66
68{
69 return _sample;
70}
71
73{
74 return _waveform.sample(time);
75}
76
78{
79 return _waveform.getInterpolationDegree();
80}
81
86
87std::vector<std::optional<double>> Data::getLowerBound() const
88{
89 return _lowerBound;
90}
91
92std::vector<std::optional<double>> Data::getUpperBound() const
93{
94 return _upperBound;
95}
96
98{
99 if (stamples().size() > 1) { // Needed to avoid CompositionalCouplingScheme callong moveToNextWindow on same Data multiple times. Could be simplified by replacing Waveform::move() with clearBefore(double time). See https://github.com/precice/precice/issues/1821.
100 waveform().move();
101 PRECICE_ASSERT(stamples().size() == 1);
102 setGlobalSample(stamples().back().sample); // @todo at some point we should not need this anymore, when mapping, acceleration ... directly work on _timeStepsWaveform, see https://github.com/precice/precice/issues/1645
103 }
104}
105
107{
108 PRECICE_ASSERT(sample.dataDims == getDimensions(), "Sample has incorrect data dimension");
109 setGlobalSample(sample); // @todo at some point we should not need this anymore, when mapping, acceleration ... directly work on _timeStepsWaveform, see https://github.com/precice/precice/issues/1645
110 _waveform.setSampleAtTime(time, sample);
111}
112
114{
115 PRECICE_ASSERT(not sample.values.hasNaN());
116 _sample = sample;
117}
118
123
124void Data::emplaceSampleAtTime(double time, std::initializer_list<double> values)
125{
127 Eigen::Map<const Eigen::VectorXd>(values.begin(), values.size())});
128}
129
130void Data::emplaceSampleAtTime(double time, std::initializer_list<double> values, std::initializer_list<double> gradients)
131{
132 PRECICE_ASSERT(gradients.size() == values.size() * getSpatialDimensions(), "Gradient isn't correctly sized", values.size(), gradients.size());
133 auto nVertices = values.size() / getDimensions();
135 Eigen::Map<const Eigen::VectorXd>(values.begin(), values.size()),
136 Eigen::Map<const Eigen::MatrixXd>(gradients.begin(), getSpatialDimensions(), nVertices * getDimensions())});
137}
138
139const std::string &Data::getName() const
140{
141 return _name;
142}
143
145{
146 return _id;
147}
148
150{
151 return _hasGradient;
152}
153
155{
156 return !_waveform.empty();
157}
158
160{
161 _hasGradient = true;
162};
163
165{
166 return _dimensions;
167}
168
169void Data::allocateValues(int expectedCount)
170{
171 using SizeType = std::remove_cv<decltype(expectedCount)>::type;
172 // Allocate data values
173 const SizeType expectedSize = expectedCount * _dimensions;
174 const auto actualSize = static_cast<SizeType>(_sample.values.size());
175 // Shrink Buffer
176 if (expectedSize < actualSize) {
177 _sample.values.resize(expectedSize);
178 }
179 // Enlarge Buffer
180 if (expectedSize > actualSize) {
181 const auto leftToAllocate = expectedSize - actualSize;
182 utils::append(_sample.values, Eigen::VectorXd(Eigen::VectorXd::Zero(leftToAllocate)));
183 }
184 PRECICE_DEBUG("Data {} now has {} values", _name, _sample.values.size());
185
186 // Allocate gradient data values
187 if (_hasGradient) {
188 const SizeType spaceDimensions = _spatialDimensions;
189
190 const SizeType expectedColumnSize = expectedCount * _dimensions;
191 const auto actualColumnSize = static_cast<SizeType>(_sample.gradients.cols());
192
193 // Shrink Buffer
194 if (expectedColumnSize < actualColumnSize) {
195 _sample.gradients.resize(spaceDimensions, expectedColumnSize);
196 }
197
198 // Enlarge Buffer
199 if (expectedColumnSize > actualColumnSize) {
200 const auto columnLeftToAllocate = expectedColumnSize - actualColumnSize;
201 utils::append(_sample.gradients, Eigen::MatrixXd(Eigen::MatrixXd::Zero(spaceDimensions, columnLeftToAllocate)));
202 }
203 PRECICE_DEBUG("Gradient Data {} now has {} x {} values", _name, _sample.gradients.rows(), _sample.gradients.cols());
204 }
205}
206
208{
209 return _spatialDimensions;
210}
211
212} // namespace precice::mesh
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:61
#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: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
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
void move()
Move this Storage by deleting all stamples except the one at the end of the window.
Definition Waveform.cpp:95
provides Mesh, Data and primitives.
contains the time interpolation logic.
Definition Sample.hpp:8
void append(Eigen::VectorXd &v, double value)
int DataID
Definition Types.hpp:25
STL namespace.