preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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
12 mesh::PtrData data,
13 mesh::PtrMesh mesh,
14 bool requiresInitialization,
15 bool exchangeSubsteps,
16 Direction direction)
17 : requiresInitialization(requiresInitialization),
18 _mesh(std::move(mesh)),
19 _data(std::move(data)),
20 _previousTimeStepsStorage(),
21 _exchangeSubsteps(exchangeSubsteps),
22 _direction(direction)
23{
24 PRECICE_ASSERT(_data != nullptr);
25 _previousTimeStepsStorage = _data->timeStepsStorage();
26
27 PRECICE_ASSERT(_mesh != nullptr);
29}
30
32{
33 PRECICE_ASSERT(_data != nullptr);
34 return _data->getDimensions();
35}
36
38{
39 return _mesh->nVertices() * getDimensions();
40}
41
43{
44 return _mesh->nVertices();
45}
46
47const Eigen::VectorXd &CouplingData::values() const
48{
49 return sample().values;
50}
51
52const Eigen::MatrixXd &CouplingData::gradients() const
53{
54 return sample().gradients;
55}
56
58{
59 const int rows = meshDimensions();
60 PRECICE_ASSERT(sample().gradients.rows() == rows, sample().gradients.rows(), rows);
61 return rows;
62}
63
65{
66 const int cols = getSize();
67 PRECICE_ASSERT(sample().gradients.cols() == cols, sample().gradients.cols(), cols);
68 return cols;
69}
70
72{
73 PRECICE_ASSERT(_data != nullptr);
74 return _data->sample();
75}
76
78{
79 PRECICE_ASSERT(_data != nullptr);
80 return _data->timeStepsStorage();
81}
82
84{
85 PRECICE_ASSERT(_data != nullptr);
86 return _data->timeStepsStorage();
87}
88
93
94Eigen::MatrixXd CouplingData::getPreviousGradientsAtTime(double relativeDt)
95{
97}
98
100{
101 PRECICE_ASSERT(not sample.values.hasNaN());
102 _data->setSampleAtTime(time, sample);
103}
104
106{
107 PRECICE_ASSERT(not sample.values.hasNaN());
108 _data->setGlobalSample(sample);
109}
110
112{
113 if (!hasGradient()) {
114 auto zero = time::Sample(getDimensions(), nVertices());
115 zero.setZero();
116 _data->setSampleAtTime(time, zero);
117 return;
118 }
120 zero.setZero();
121 _data->setSampleAtTime(time, zero);
122}
123
125{
126 _data->emplaceSampleAtTime(time);
127}
128
130{
131 _data->emplaceSampleAtTime(time, values);
132}
133
135{
136 _data->emplaceSampleAtTime(time, values, gradients);
137}
138
140{
141 PRECICE_ASSERT(_data != nullptr);
142 return _data->hasGradient();
143}
144
146{
147 return _mesh->getDimensions();
148}
149
151{
152 // TODO port this to subcyling
153
154 // The mesh was reinitialized and new written data will be added later in advance().
155 // Meaning all samples are based on a different mesh.
156 // Without remapping, the best we can do is setting them to zero samples.
157 // We keep the timestamps not to break convergence measures, accelerations, and actions
158 auto zero = time::Sample(getDimensions(), nVertices());
159 zero.setZero();
160
161 _data->timeStepsStorage().setAllSamples(zero);
163}
164
166{
167 const auto &stamples = this->stamples();
168 PRECICE_ASSERT(stamples.size() > 0);
169 _previousTimeStepsStorage = _data->timeStepsStorage();
170}
171
172const Eigen::VectorXd &CouplingData::previousIteration() const
173{
175 return _previousTimeStepsStorage.stamples().back().sample.values;
176}
177
178const Eigen::MatrixXd &CouplingData::previousIterationGradients() const
179{
181 return _previousTimeStepsStorage.stamples().back().sample.gradients;
182}
183
185{
187 return _previousTimeStepsStorage.stamples().back().sample.values.size();
188}
189
191{
192 return _mesh->getID();
193}
194
196{
197 return _data->getID();
198}
199
201{
202 return _data->getName();
203}
204
206{
207 return _mesh->getName();
208}
209
211{
212 return _mesh->getVertexOffsets();
213}
214
219
221{
223 //_data->moveToNextWindow();
224 // _previousTimeStepsStorage = _data->timeStepsStorage();
225 }
226 _data->moveToNextWindow();
227 _previousTimeStepsStorage = _data->timeStepsStorage();
228}
229
231{
232 return _exchangeSubsteps;
233}
234} // 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().
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 _timeStepsStorage.
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.
time::Storage & timeStepsStorage()
Returns a reference to the time step storage of the 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 _timeStepsStorage.
const Eigen::VectorXd & previousIteration() const
returns data value from previous iteration
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().
int gradientsCols() const
Returns number of columns of the stored gradients.
void setSampleAtTime(double time, time::Sample sample)
Add sample at given time to _timeStepsStorage.
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.
void storeIteration()
store _data->values() in read-only variable _previousIteration for convergence checks etc.
int getPreviousIterationSize() const
returns size of previous iteration
void moveToNextWindow()
move to next window and initialize data via extrapolation
time::Storage _previousTimeStepsStorage
Sample values of previous iteration (end of time window).
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().
SampleResult sample(double time) const
Need to use interpolation for the case with changing time grids.
Definition Storage.cpp:198
auto stamples() const
Get the stamples.
Definition Storage.hpp:93
Eigen::MatrixXd sampleGradients(double time) const
Definition Storage.cpp:228
void setAllSamples(const Sample &sample)
Overrides all existing samples.
Definition Storage.cpp:52
contains implementations of coupling schemes for coupled simulations.
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
T use_count(T... args)