preCICE
Loading...
Searching...
No Matches
CouplingData.cpp
Go to the documentation of this file.
2
3#include <utility>
4
5#include "mesh/Data.hpp"
6#include "mesh/Mesh.hpp"
8
9namespace precice::cplscheme {
10
30
32{
33 PRECICE_ASSERT(_data != nullptr);
34 return _data->getDimensions();
35}
36
38{
39 return _mesh->nVertices() * getDimensions();
40}
41
42std::vector<std::optional<double>> CouplingData::getLowerBound() const
43{
44 PRECICE_ASSERT(_data != nullptr);
45 return _data->getLowerBound();
46}
47
48std::vector<std::optional<double>> CouplingData::getUpperBound() const
49{
50 PRECICE_ASSERT(_data != nullptr);
51 return _data->getUpperBound();
52}
53
55{
56 return _mesh->nVertices();
57}
58
59const Eigen::VectorXd &CouplingData::values() const
60{
61 return sample().values;
62}
63
64const Eigen::MatrixXd &CouplingData::gradients() const
65{
66 return sample().gradients;
67}
68
70{
71 const int rows = meshDimensions();
72 PRECICE_ASSERT(sample().gradients.rows() == rows, sample().gradients.rows(), rows);
73 return rows;
74}
75
77{
78 const int cols = getSize();
79 PRECICE_ASSERT(sample().gradients.cols() == cols, sample().gradients.cols(), cols);
80 return cols;
81}
82
84{
85 PRECICE_ASSERT(_data != nullptr);
86 return _data->sample();
87}
88
90{
91 PRECICE_ASSERT(_data != nullptr);
92 return _data->waveform();
93}
94
96{
97 PRECICE_ASSERT(_data != nullptr);
98 return _data->waveform();
99}
100
102{
103 return _previousTimeStepsStorage.sample(relativeDt);
104}
105
106Eigen::MatrixXd CouplingData::getPreviousGradientsAtTime(double relativeDt)
107{
108 return _previousTimeStepsStorage.sampleGradients(relativeDt);
109}
110
112{
113 PRECICE_ASSERT(not sample.values.hasNaN());
114 _data->setSampleAtTime(time, sample);
115}
116
118{
119 PRECICE_ASSERT(not sample.values.hasNaN());
120 _data->setGlobalSample(sample);
121}
122
124{
125 if (!hasGradient()) {
126 auto zero = time::Sample(getDimensions(), nVertices());
127 zero.setZero();
128 _data->setSampleAtTime(time, zero);
129 return;
130 }
132 zero.setZero();
133 _data->setSampleAtTime(time, zero);
134}
135
137{
138 _data->emplaceSampleAtTime(time);
139}
140
141void CouplingData::emplaceSampleAtTime(double time, std::initializer_list<double> values)
142{
143 _data->emplaceSampleAtTime(time, values);
144}
145
146void CouplingData::emplaceSampleAtTime(double time, std::initializer_list<double> values, std::initializer_list<double> gradients)
147{
148 _data->emplaceSampleAtTime(time, values, gradients);
149}
150
152{
153 PRECICE_ASSERT(_data != nullptr);
154 return _data->hasGradient();
155}
156
158{
159 return _mesh->getDimensions();
160}
161
163{
164 // TODO port this to subcyling
165
166 // The mesh was reinitialized and new written data will be added later in advance().
167 // Meaning all samples are based on a different mesh.
168 // Without remapping, the best we can do is setting them to zero samples.
169 // We keep the timestamps not to break convergence measures, accelerations, and actions
170 auto zero = time::Sample(getDimensions(), nVertices());
171 zero.setZero();
172
173 _data->waveform().setAllSamples(zero);
174 _previousTimeStepsStorage.setAllSamples(zero);
175}
176
178{
179 const auto &stamples = this->stamples();
180 PRECICE_ASSERT(stamples.size() > 0);
181 _previousTimeStepsStorage = _data->waveform();
182}
183
184const Eigen::VectorXd &CouplingData::previousIteration() const
185{
186 PRECICE_ASSERT(!_previousTimeStepsStorage.stamples().empty());
187 return _previousTimeStepsStorage.stamples().back().sample.values;
188}
189
190const Eigen::MatrixXd &CouplingData::previousIterationGradients() const
191{
192 PRECICE_ASSERT(!_previousTimeStepsStorage.stamples().empty());
193 return _previousTimeStepsStorage.stamples().back().sample.gradients;
194}
195
197{
198 PRECICE_ASSERT(!_previousTimeStepsStorage.stamples().empty());
199 return _previousTimeStepsStorage.stamples().back().sample.values.size();
200}
201
203{
204 return _mesh->getID();
205}
206
208{
209 return _data->getID();
210}
211
212std::string CouplingData::getDataName() const
213{
214 return _data->getName();
215}
216
217std::string CouplingData::getMeshName() const
218{
219 return _mesh->getName();
220}
221
223{
224 return _mesh->getVertexOffsets();
225}
226
231
233{
235 //_data->moveToNextWindow();
236 // _previousTimeStepsStorage = _data->waveform();
237 }
238 _data->moveToNextWindow();
239 _previousTimeStepsStorage = _data->waveform();
240}
241
243{
244 return _exchangeSubsteps;
245}
246} // namespace precice::cplscheme
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
std::vector< int > getVertexOffsets()
get vertex offsets of this CouplingData's mesh. See Mesh::getVertexOffsets().
time::Waveform & waveform()
Returns a reference to the time step storage of the data.
int getMeshID()
get ID of this CouplingData's mesh. See Mesh::getID().
void initializeWithZeroAtTime(double time)
Add sample with zero values at given time to the Waveform.
Eigen::MatrixXd getPreviousGradientsAtTime(double relativeDt)
const Eigen::VectorXd & values() const
Returns a const reference to the data values.
bool _exchangeSubsteps
If true, all substeps will be sent / received for this coupling data.
const Eigen::MatrixXd & previousIterationGradients() const
returns gradient data from previous iteration
void emplaceSampleAtTime(double time)
Creates an empty sample at given time.
mesh::PtrMesh _mesh
Mesh associated with this CouplingData.
auto stamples() const
Returns the stamples in the Waveform.
const Eigen::VectorXd & previousIteration() const
returns data value from previous iteration
std::vector< std::optional< double > > getLowerBound() const
get configured lower bound(s) of the data
int meshDimensions() const
Returns the dimensions of the current mesh (2D or 3D)
void setGlobalSample(const time::Sample &sample)
Set _data::_sample.
Direction getDirection() const
get direction of this coupling data
int getDataID()
get ID of this CouplingData's data. See Data::getID().
time::Waveform _previousTimeStepsStorage
Sample values of previous iteration (end of time window).
int gradientsCols() const
Returns number of columns of the stored gradients.
void setSampleAtTime(double time, time::Sample sample)
Add sample at given time to the Waveform.
const time::Sample & sample() const
Returns a const reference to the data Sample.
const Eigen::MatrixXd & gradients() const
Returns a const reference to the gradient data values.
const bool requiresInitialization
True, if the data values of this CouplingData require to be initialized by this participant.
void storeIteration()
store _data->values() in read-only variable _previousIteration for convergence checks etc.
std::vector< std::optional< double > > getUpperBound() const
get configured upper bound(s) of the data
int getPreviousIterationSize() const
returns size of previous iteration
void moveToNextWindow()
move to next window and initialize data via extrapolation
bool hasGradient() const
Returns if the data contains gradient data.
time::SampleResult getPreviousValuesAtTime(double relativeDt)
returns previous data interpolated to the relativeDt time
void reinitialize()
Reshape the past iterations and initial sample during remeshing.
int gradientsRows() const
Returns number of rows of the stored gradients.
CouplingData(mesh::PtrData data, mesh::PtrMesh mesh, bool requiresInitialization, bool exchangeSubsteps, Direction direction)
std::string getMeshName() const
get name of this CouplingData's mesh. See Mesh::getName().
mesh::PtrData _data
Data associated with this CouplingData.
std::string getDataName() const
get name of this CouplingData's data. See Data::getName().
contains implementations of coupling schemes for coupled simulations.
provides Mesh, Data and primitives.
std::shared_ptr< Data > PtrData
std::shared_ptr< Mesh > PtrMesh
contains the time interpolation logic.
Definition Sample.hpp:8
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
Eigen::VectorXd values
Definition Sample.hpp:64