preCICE v3.1.2
Loading...
Searching...
No Matches
WriteDataContext.cpp
Go to the documentation of this file.
2
4
5namespace precice::impl {
6
7logging::Logger WriteDataContext::_log{"impl::WriteDataContext"};
8
10 mesh::PtrMesh mesh)
11 : DataContext(data, mesh),
12 _writeDataBuffer(data->getDimensions())
13{
14}
15
21
23{
24 _providedData->timeStepsStorage().trimAfter(time);
25
26 // reset all toData
27 PRECICE_ASSERT(!hasReadMapping(), "Read mapping is not allowed for WriteDataContext.");
28 if (hasWriteMapping()) {
29 std::for_each(_mappingContexts.begin(), _mappingContexts.end(), [time](auto &context) { context.toData->timeStepsStorage().trimAfter(time); });
30 }
31}
32
34{
35 PRECICE_ASSERT(vertices.size() * getDataDimensions() == values.size());
37
38 Eigen::Map<const Eigen::MatrixXd> inputData(values.data(), getDataDimensions(), vertices.size());
39 Eigen::Map<Eigen::MatrixXd> localData(_writeDataBuffer.values.data(), getDataDimensions(), getMeshVertexCount());
40
41 for (int i = 0; i < static_cast<int>(vertices.size()); ++i) {
42 PRECICE_ASSERT(vertices[i] < localData.cols());
43 localData.col(vertices[i]) = inputData.col(i);
44 }
45}
46
48{
49 const auto gradientComponents = getSpatialDimensions() * getDataDimensions();
50
51 PRECICE_ASSERT(gradientComponents * vertices.size() == gradients.size());
53
54 Eigen::Map<const Eigen::MatrixXd> inputGradients(gradients.data(), gradientComponents, vertices.size());
55 Eigen::Map<Eigen::MatrixXd> localGradients(_writeDataBuffer.gradients.data(), gradientComponents, getMeshVertexCount());
56
57 for (int i = 0; i < static_cast<int>(vertices.size()); ++i) {
58 PRECICE_ASSERT(vertices[i] < localGradients.cols());
59 localGradients.col(vertices[i]) = inputGradients.col(i);
60 }
61}
62
64{
65 using SizeType = std::remove_cv<decltype(nVertices)>::type;
66
67 // Allocate data values
68 const SizeType expectedSize = nVertices * getDataDimensions();
69 const auto actualSize = static_cast<SizeType>(_writeDataBuffer.values.size());
70 // Shrink Buffer
71 if (expectedSize < actualSize) {
72 _writeDataBuffer.values.resize(expectedSize);
73 }
74 // Enlarge Buffer
75 if (expectedSize > actualSize) {
76 const auto leftToAllocate = expectedSize - actualSize;
77 utils::append(_writeDataBuffer.values, Eigen::VectorXd(Eigen::VectorXd::Zero(leftToAllocate)));
78 }
79 PRECICE_DEBUG("Data {} now has {} values", getDataName(), _writeDataBuffer.values.size());
80
81 // Allocate gradient data values
82 if (_providedData->hasGradient()) {
83 const SizeType spaceDimensions = getSpatialDimensions();
84
85 const SizeType expectedColumnSize = expectedSize * getDataDimensions();
86 const auto actualColumnSize = static_cast<SizeType>(_writeDataBuffer.gradients.cols());
87
88 // Shrink Buffer
89 if (expectedColumnSize < actualColumnSize) {
90 _writeDataBuffer.gradients.resize(spaceDimensions, expectedColumnSize);
91 }
92
93 // Enlarge Buffer
94 if (expectedColumnSize > actualColumnSize) {
95 const auto columnLeftToAllocate = expectedColumnSize - actualColumnSize;
96 utils::append(_writeDataBuffer.gradients, Eigen::MatrixXd(Eigen::MatrixXd::Zero(spaceDimensions, columnLeftToAllocate)));
97 }
98 PRECICE_DEBUG("Gradient Data {} now has {} x {} values", getDataName(), _writeDataBuffer.gradients.rows(), _writeDataBuffer.gradients.cols());
99 }
100}
101
103{
104 _providedData->setSampleAtTime(currentTime, _writeDataBuffer);
105}
106
108{
109 PRECICE_ASSERT(meshContext.mesh->hasDataName(getDataName()));
110 mesh::PtrData data = meshContext.mesh->data(getDataName());
111 PRECICE_ASSERT(data != _providedData, "Data the write mapping is mapping to needs to be different from _providedData");
112 mappingContext.fromData = _providedData;
113 mappingContext.toData = data;
114 appendMapping(mappingContext);
116}
117
118} // namespace precice::impl
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:64
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
Stores one Data object with related mesh.
int getMeshVertexCount() const
Get the number of vertices of mesh.
bool hasWriteMapping() const
Informs the user whether this DataContext has any write mapping.
void appendMapping(MappingContext mappingContext)
Helper to append a mappingContext, fromData and toData to the corresponding data containers.
std::string getDataName() const
Get the Name of _providedData.
int getDataDimensions() const
Get the dimensions of _providedData.
int getSpatialDimensions() const
Get the spatial dimensions of _providedData.
bool hasReadMapping() const
Informs the user whether this DataContext has any read mapping.
mesh::PtrData _providedData
Unique data this context is associated with.
std::vector< MappingContext > _mappingContexts
Defines all mappings associated to this DataContext. A DataContext may also exist without a mapping.
void trimAfter(double time)
Removes stample before time and (if mapping exists) fromData or toData.
time::Sample _writeDataBuffer
Buffer to store written data until it is copied to _providedData->timeStepsStorage()
void writeGradientsIntoDataBuffer(::precice::span< const VertexID > vertices, ::precice::span< const double > gradients)
Store gradients in _writeDataBuffer.
void storeBufferedData(double currentTime)
Store data from _writeDataBuffer in persistent storage.
void resetBuffer()
Resets writeDataBuffer.
WriteDataContext(mesh::PtrData data, mesh::PtrMesh mesh)
Construct a new WriteDataContext object without a mapping.
void writeValuesIntoDataBuffer(::precice::span< const VertexID > vertices, ::precice::span< const double > values)
Store values in _writeDataBuffer.
void appendMappingConfiguration(MappingContext &mappingContext, const MeshContext &meshContext) override
Adds a MappingContext and the MeshContext required by the write mapping to the corresponding WriteDat...
A C++ 11 implementation of the non-owning C++20 std::span type.
Definition span.hpp:284
constexpr pointer data() const noexcept
Definition span.hpp:500
constexpr size_type size() const noexcept
Definition span.hpp:469
T for_each(T... args)
void append(Eigen::VectorXd &v, double value)
Holds a data mapping and related information.
mesh::PtrData fromData
data which is mapped from mesh
mesh::PtrData toData
data which is mapped to mesh
Stores a mesh and related objects and data.
mesh::PtrMesh mesh
Mesh holding the geometry data structure.
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