preCICE v3.2.0
Loading...
Searching...
No Matches
ReadDataContext.cpp
Go to the documentation of this file.
1#include "ReadDataContext.hpp"
2
3namespace precice::impl {
4
5logging::Logger ReadDataContext::_log{"impl::ReadDataContext"};
6
13
15{
16 PRECICE_ASSERT(!hasReadMapping(), "The read data context must be unique. Otherwise we would have an ambiguous read data operation on the user side.");
17 // The read mapping must be unique, but having read and write in the same context is not possible either
19 PRECICE_ASSERT(meshContext.mesh->hasDataName(getDataName()));
20 mesh::PtrData data = meshContext.mesh->data(getDataName());
21 PRECICE_ASSERT(data != _providedData, "Data the read mapping is mapping from needs to be different from _providedData");
22 mappingContext.fromData = data;
23 mappingContext.toData = _providedData;
24 appendMapping(mappingContext);
26}
27
29{
30 return _providedData->hasSamples();
31}
32
34{
35 Eigen::Map<Eigen::MatrixXd> outputData(values.data(), getDataDimensions(), values.size());
36 auto sampleResult = _providedData->sampleAtTime(readTime);
37 auto localData = sampleResult.values().reshaped(getDataDimensions(), getMeshVertexCount());
38 for (int i = 0; i < static_cast<int>(vertices.size()); ++i) {
39 outputData.col(i) = localData.col(vertices[i]);
40 }
41}
42
44{
45 PRECICE_TRACE(getMeshName(), getDataName(), coordinates.size(), values.size(), readTime);
48
49 // First, check if we have the current readTime already in our MappingDataCache
50 if (!mappingCache->hasDataAtTimeStamp(readTime)) {
51 // if not, sample our waveform and update the cache
52 justInTimeMapping->updateMappingDataCache(*mappingCache.get(), _providedData->sampleAtTime(readTime).values());
53 mappingCache->setTimeStamp(readTime);
54 }
55 // Now we are certain that our cache contains the data at readTime
56 Eigen::Map<const Eigen::MatrixXd> coords(coordinates.data(), getSpatialDimensions(), coordinates.size() / getSpatialDimensions());
57 Eigen::Map<Eigen::MatrixXd> target(values.data(), getDataDimensions(), values.size() / getDataDimensions());
58
59 // ...hence, we forward the coordinates, cache and the target to the just-in-time mapping
60 justInTimeMapping->mapConsistentAt(coords, *mappingCache.get(), target);
61}
62
64{
65 return _providedData->getWaveformDegree();
66}
67
69{
72 for (auto &context : _mappingContexts) {
73 auto id = context.fromData->getID();
74 if (from.contains(id)) {
75 if (from.toKeep(id)) {
76 context.toData->timeStepsStorage().clearExceptLast();
77 } else {
78 context.toData->timeStepsStorage().clear();
79 }
80 }
81 }
82}
83
85{
88 for (auto &context : _mappingContexts) {
89 if (from.contains(context.fromData->getID())) {
90 context.toData->timeStepsStorage().trimAfter(t);
91 }
92 }
93}
94
95} // namespace precice::impl
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
int getMeshVertexCount() const
Get the number of vertices of mesh.
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.
bool hasMapping() const
Informs the user whether this DataContext has any _mappingContext.
int getSpatialDimensions() const
Get the spatial dimensions of _providedData.
bool hasReadMapping() const
Informs the user whether this DataContext has any read mapping.
std::vector< MappingContext > _mappingContexts
Defines all mappings associated to this DataContext. A DataContext may also exist without a mapping.
std::string getMeshName() const
Get the name of _mesh.
void trimToDataAfterFor(const cplscheme::ImplicitData &from, double t)
Trims all toData of associated mappings after the given t.
bool hasSamples() const
Are there samples to read from?
int getWaveformDegree() const
Gets degree of waveform.
ReadDataContext(mesh::PtrData data, mesh::PtrMesh mesh)
Construct a new ReadDataContext object without a mapping.
void readValues(::precice::span< const VertexID > vertices, double time, ::precice::span< double > values) const
Samples data at a given point in time within the current time window for given indices.
void clearToDataFor(const cplscheme::ImplicitData &from)
Removes all toData samples from mappings.
void mapAndReadValues(::precice::span< const double > coordinates, double readTime, ::precice::span< double > values)
Forwards the just-in-time mapping API call for reading data to the data context.
void appendMappingConfiguration(MappingContext &mappingContext, const MeshContext &meshContext) override
Adds a MappingContext and the MeshContext required by the read mapping to the corresponding ReadDataC...
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
provides Mesh, Data and primitives.
std::shared_ptr< Data > PtrData
std::shared_ptr< Mesh > PtrMesh
bool contains(DataID did) const
bool toKeep(DataID did) const
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.