preCICE v3.2.0
Loading...
Searching...
No Matches
WriteDataContext.cpp
Go to the documentation of this file.
2
3namespace precice::impl {
4
5logging::Logger WriteDataContext::_log{"impl::WriteDataContext"};
6
13
15{
17 _writeDataBuffer.values.setZero();
18 _writeDataBuffer.gradients.setZero();
19}
20
22{
23 _providedData->timeStepsStorage().trimAfter(time);
24
25 // reset all toData
26 PRECICE_ASSERT(!hasReadMapping(), "Read mapping is not allowed for WriteDataContext.");
27 if (hasWriteMapping()) {
28 std::for_each(_mappingContexts.begin(), _mappingContexts.end(), [time](auto &context) { context.toData->timeStepsStorage().trimAfter(time); });
29 }
30}
31
33{
36 // finalize mapping the data stored in the cache and transfer it to the _writeDataBuffer
37 Eigen::Map<Eigen::MatrixXd> map(_writeDataBuffer.values.data(), _providedData->getDimensions(), _writeDataBuffer.values.size() / _providedData->getDimensions());
38 justInTimeMapping->completeJustInTimeMapping(*mappingCache, map);
39 }
40}
41
43{
47 PRECICE_ASSERT((coordinates.size() / getSpatialDimensions()) * getDataDimensions() == values.size());
48 PRECICE_ASSERT(_writeDataBuffer.values.data());
49
50 // We forward both, the _writeDataBuffer and the cache to the justInTimeMapping
51 Eigen::Map<const Eigen::MatrixXd> coords(coordinates.data(), getSpatialDimensions(), coordinates.size() / getSpatialDimensions());
52 Eigen::Map<const Eigen::MatrixXd> inputData(values.data(), getDataDimensions(), coordinates.size() / getDataDimensions());
53 Eigen::Map<Eigen::MatrixXd> localData(_writeDataBuffer.values.data(), getDataDimensions(), getMeshVertexCount());
54
55 // Function to fill the localData
56 justInTimeMapping->mapConservativeAt(coords, inputData, *mappingCache, localData);
57}
58
60{
61 PRECICE_ASSERT(vertices.size() * getDataDimensions() == values.size());
62 PRECICE_ASSERT(_writeDataBuffer.values.data());
63
64 Eigen::Map<const Eigen::MatrixXd> inputData(values.data(), getDataDimensions(), vertices.size());
65 Eigen::Map<Eigen::MatrixXd> localData(_writeDataBuffer.values.data(), getDataDimensions(), getMeshVertexCount());
66
67 for (int i = 0; i < static_cast<int>(vertices.size()); ++i) {
68 PRECICE_ASSERT(vertices[i] < localData.cols());
69 localData.col(vertices[i]) = inputData.col(i);
70 }
71}
72
74{
75 const auto gradientComponents = getSpatialDimensions() * getDataDimensions();
76
77 PRECICE_ASSERT(gradientComponents * vertices.size() == gradients.size());
78 PRECICE_ASSERT(_writeDataBuffer.gradients.data());
79
80 Eigen::Map<const Eigen::MatrixXd> inputGradients(gradients.data(), gradientComponents, vertices.size());
81 Eigen::Map<Eigen::MatrixXd> localGradients(_writeDataBuffer.gradients.data(), gradientComponents, getMeshVertexCount());
82
83 for (int i = 0; i < static_cast<int>(vertices.size()); ++i) {
84 PRECICE_ASSERT(vertices[i] < localGradients.cols());
85 localGradients.col(vertices[i]) = inputGradients.col(i);
86 }
87}
88
90{
91 using SizeType = std::remove_cv<decltype(nVertices)>::type;
92
93 // Allocate data values
94 const SizeType expectedSize = nVertices * getDataDimensions();
95 const auto actualSize = static_cast<SizeType>(_writeDataBuffer.values.size());
96 const auto change = expectedSize - actualSize;
97
98 _writeDataBuffer.values.conservativeResize(expectedSize);
99 // Enlarge Buffer
100 if (change > 0) {
101 _writeDataBuffer.values.tail(change).setZero();
102 }
103 PRECICE_DEBUG("Data {} now has {} values", getDataName(), _writeDataBuffer.values.size());
104
105 if (!_providedData->hasGradient()) {
106 return;
107 }
108
109 // Allocate gradient data values
110 _writeDataBuffer.gradients.conservativeResize(getSpatialDimensions(), expectedSize);
111
112 // Enlarge Buffer
113 if (change > 0) {
114 _writeDataBuffer.gradients.rightCols(change).setZero();
115 }
116 PRECICE_DEBUG("Gradient Data {} now has {} x {} values", getDataName(), _writeDataBuffer.gradients.rows(), _writeDataBuffer.gradients.cols());
117}
118
120{
121 _providedData->setSampleAtTime(currentTime, _writeDataBuffer);
122}
123
125{
126 PRECICE_ASSERT(meshContext.mesh->hasDataName(getDataName()));
127 mesh::PtrData data = meshContext.mesh->data(getDataName());
128 PRECICE_ASSERT(data != _providedData, "Data the write mapping is mapping to needs to be different from _providedData");
129 mappingContext.fromData = _providedData;
130 mappingContext.toData = data;
131 appendMapping(mappingContext);
133}
134
135} // namespace precice::impl
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:61
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
int getMeshVertexCount() const
Get the number of vertices of mesh.
bool hasWriteMapping() const
Informs the user whether this DataContext has any write mapping.
DataContext(mesh::PtrData data, mesh::PtrMesh mesh)
Construct a new DataContext without a mapping. Protected, because only ReadDataContext and WriteDataC...
mapping::PtrMapping justInTimeMapping
The just-in-time mapping for this data context.
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.
std::unique_ptr< mapping::impl::MappingDataCache > mappingCache
Cache for just-in-time mapping.
int getSpatialDimensions() const
Get the spatial dimensions of _providedData.
bool hasReadMapping() const
Informs the user whether this DataContext has any read mapping.
void invalidateMappingCacheAndResetData()
Resets the time stamp of the MappingDataCache and resets the data it holds.
std::vector< MappingContext > _mappingContexts
Defines all mappings associated to this DataContext. A DataContext may also exist without a mapping.
void completeJustInTimeMapping()
Evaluates the MappingDataCache and stores the result in the _writeDataBuffer.
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 resetBufferedData()
Resets the writeDataBuffer and the mapping data cache.
void writeAndMapValues(::precice::span< const double > coordinates, ::precice::span< const double > values)
Forwards the just-in-time mapping API call for writing data to the data context.
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)
provides Mesh, Data and primitives.
std::shared_ptr< Data > PtrData
std::shared_ptr< Mesh > PtrMesh
contains the time interpolation logic.
Definition Sample.hpp:8
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.