preCICE v3.1.2
Loading...
Searching...
No Matches
CouplingData.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Core>
4#include <vector>
7#include "time/Storage.hpp"
8#include "utils/assertion.hpp"
9
10namespace precice {
11namespace cplscheme {
12
14public:
15 enum struct Direction : bool { Send,
16 Receive };
17
19 mesh::PtrData data,
20 mesh::PtrMesh mesh,
23 Direction direction);
24
25 int getDimensions() const;
26
27 int getSize() const;
28
30 Eigen::VectorXd &values();
31
33 const Eigen::VectorXd &values() const;
34
36 Eigen::MatrixXd &gradients();
37
39 const Eigen::MatrixXd &gradients() const;
40
43
45 const time::Sample &sample() const;
46
49
51 Eigen::VectorXd getPreviousValuesAtTime(double relativeDt);
52
53 Eigen::MatrixXd getPreviousGradientsAtTime(double relativeDt);
54
56 const time::Storage &timeStepsStorage() const;
57
59 auto stamples() const
60 {
61 return timeStepsStorage().stamples();
62 }
63
65 void setSampleAtTime(double time, time::Sample sample);
66
68 bool hasGradient() const;
69
71 int meshDimensions() const;
72
74 void storeIteration();
75
77 const Eigen::VectorXd &previousIteration() const;
78
80 const Eigen::MatrixXd &previousIterationGradients() const;
81
83 int getPreviousIterationSize() const;
84
86 int getMeshID();
87
89 int getDataID();
90
93
96
98 Direction getDirection() const;
99
102
104 void moveToNextWindow();
105
106 bool exchangeSubsteps() const;
107
108private:
109 logging::Logger _log{"cplscheme::CouplingData"};
110
113
116
119
122
124};
125
126} // namespace cplscheme
127} // namespace precice
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().
const bool requiresInitialization
True, if the data values of this CouplingData require to be initialized by this participant.
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.
This class provides a lightweight logger.
Definition Logger.hpp:16
auto stamples() const
Get the stamples.
Definition Storage.hpp:75
Main namespace of the precice library.