preCICE
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/Waveform.hpp"
8#include "utils/assertion.hpp"
9
10namespace precice::cplscheme {
11
13public:
14 enum struct Direction : bool { Send,
16
18 mesh::PtrData data,
22 Direction direction);
23
24 int getDimensions() const;
25
26 int getSize() const;
27
28 int nVertices() const;
29
31 const Eigen::VectorXd &values() const;
32
34 const Eigen::MatrixXd &gradients() const;
35
37 int gradientsRows() const;
38
40 int gradientsCols() const;
41
43 const time::Sample &sample() const;
44
47
50
51 Eigen::MatrixXd getPreviousGradientsAtTime(double relativeDt);
52
54 const time::Waveform &waveform() const;
55
57 auto stamples() const
58 {
59 return waveform().stamples();
60 }
61
64
66 void setGlobalSample(const time::Sample &sample); // @todo try to remove this function
67
69 void initializeWithZeroAtTime(double time);
70
72 void emplaceSampleAtTime(double time);
73
75 void emplaceSampleAtTime(double time, std::initializer_list<double> values);
76
78 void emplaceSampleAtTime(double time, std::initializer_list<double> values, std::initializer_list<double> gradients);
79
81 bool hasGradient() const;
82
84 int meshDimensions() const;
85
87 void reinitialize();
88
90 void storeIteration();
91
93 const Eigen::VectorXd &previousIteration() const;
94
96 const Eigen::MatrixXd &previousIterationGradients() const;
97
99 int getPreviousIterationSize() const;
100
102 int getMeshID();
103
105 int getDataID();
106
108 std::string getDataName() const;
109
111 std::string getMeshName() const;
112
114 std::vector<std::optional<double>> getLowerBound() const;
115
117 std::vector<std::optional<double>> getUpperBound() const;
118
120 std::vector<int> getVertexOffsets();
121
123 Direction getDirection() const;
124
127
129 void moveToNextWindow();
130
131 bool exchangeSubsteps() const;
132
133private:
134 logging::Logger _log{"cplscheme::CouplingData"};
135
138
141
144
147
149};
150
151} // namespace precice::cplscheme
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().
This class provides a lightweight logger.
Definition Logger.hpp:17
auto stamples() const
Get the stamples.
Definition Waveform.hpp:89
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