preCICE v3.1.2
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
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 sample().values.size();
40}
41
42Eigen::VectorXd &CouplingData::values()
43{
44 return sample().values;
45}
46
47const Eigen::VectorXd &CouplingData::values() const
48{
49 return sample().values;
50}
51
52Eigen::MatrixXd &CouplingData::gradients()
53{
54 return sample().gradients;
55}
56
57const Eigen::MatrixXd &CouplingData::gradients() const
58{
59 return sample().gradients;
60}
61
63{
64 PRECICE_ASSERT(_data != nullptr);
65 return _data->timeStepsStorage();
66}
67
69{
70 PRECICE_ASSERT(_data != nullptr);
71 return _data->timeStepsStorage();
72}
73
74Eigen::VectorXd CouplingData::getPreviousValuesAtTime(double relativeDt)
75{
76 return _previousTimeStepsStorage.sample(relativeDt);
77}
78
79Eigen::MatrixXd CouplingData::getPreviousGradientsAtTime(double relativeDt)
80{
82}
83
85{
86 PRECICE_ASSERT(not sample.values.hasNaN());
87 this->sample() = sample; // @todo at some point we should not need this anymore, when mapping, acceleration ... directly work on _timeStepsStorage
88 _data->setSampleAtTime(time, sample);
89}
90
92{
93 PRECICE_ASSERT(_data != nullptr);
94 return _data->hasGradient();
95}
96
98{
99 return _mesh->getDimensions();
100}
101
103{
104 const auto &stamples = this->stamples();
105 PRECICE_ASSERT(stamples.size() > 0);
106 this->sample() = stamples.back().sample;
107 _previousTimeStepsStorage = _data->timeStepsStorage();
108}
109
110const Eigen::VectorXd &CouplingData::previousIteration() const
111{
113 return _previousTimeStepsStorage.stamples().back().sample.values;
114}
115
116const Eigen::MatrixXd &CouplingData::previousIterationGradients() const
117{
119 return _previousTimeStepsStorage.stamples().back().sample.gradients;
120}
121
123{
125 return _previousTimeStepsStorage.stamples().back().sample.values.size();
126}
127
129{
130 return _mesh->getID();
131}
132
134{
135 return _data->getID();
136}
137
139{
140 return _data->getName();
141}
142
144{
145 return _mesh->getVertexOffsets();
146}
147
152
154{
156 //_data->moveToNextWindow();
157 // _previousTimeStepsStorage = _data->timeStepsStorage();
158 }
159 _data->moveToNextWindow();
160 _previousTimeStepsStorage = _data->timeStepsStorage();
161}
162
164{
165 PRECICE_ASSERT(_data != nullptr);
166 return _data->sample();
167}
168
170{
171 PRECICE_ASSERT(_data != nullptr);
172 return _data->sample();
173}
174
176{
177 return _exchangeSubsteps;
178}
179} // namespace precice::cplscheme
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
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().
Eigen::MatrixXd getPreviousGradientsAtTime(double relativeDt)
bool _exchangeSubsteps
If true, all substeps will be sent / received for this coupling data.
Eigen::MatrixXd & gradients()
Returns a reference to the gradient data values.
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
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)
Eigen::VectorXd & values()
Returns a reference to the data values.
Direction getDirection() const
get direction of this coupling data
int getDataID()
get ID of this CouplingData's data. See Data::getID().
void setSampleAtTime(double time, time::Sample sample)
Add sample at given time to _timeStepsStorage.
std::string getDataName()
get name of this CouplingData's data. See Data::getName().
Eigen::VectorXd getPreviousValuesAtTime(double relativeDt)
returns previous data interpolated to the relativeDt time
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).
time::Sample & sample()
Returns a reference to the gradient data Sample.
bool hasGradient() const
Returns if the data contains gradient data.
CouplingData(mesh::PtrData data, mesh::PtrMesh mesh, bool requiresInitialization, bool exchangeSubsteps, Direction direction)
mesh::PtrData _data
Data associated with this CouplingData.
Eigen::VectorXd sample(double time) const
Need to use interpolation for the case with changing time grids.
Definition Storage.cpp:191
auto stamples() const
Get the stamples.
Definition Storage.hpp:75
Eigen::MatrixXd sampleGradients(double time) const
Definition Storage.cpp:216
contains implementations of coupling schemes for coupled simulations.
STL namespace.
Eigen::MatrixXd gradients
The gradients of the data. Use gradients.col(i) to get the gradient at vertex i.
Definition Sample.hpp:52
Eigen::VectorXd values
Definition Sample.hpp:49
T use_count(T... args)