preCICE v3.1.2
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 PRECICE_ASSERT(meshContext.mesh->hasDataName(getDataName()));
18 mesh::PtrData data = meshContext.mesh->data(getDataName());
19 PRECICE_ASSERT(data != _providedData, "Data the read mapping is mapping from needs to be different from _providedData");
20 mappingContext.fromData = data;
21 mappingContext.toData = _providedData;
22 appendMapping(mappingContext);
24}
25
27{
28 Eigen::Map<Eigen::MatrixXd> outputData(values.data(), getDataDimensions(), values.size());
29 const Eigen::MatrixXd sample{_providedData->sampleAtTime(readTime)};
30 Eigen::Map<const Eigen::MatrixXd> localData(sample.data(), getDataDimensions(), getMeshVertexCount());
31 for (int i = 0; i < static_cast<int>(vertices.size()); ++i) {
32 outputData.col(i) = localData.col(vertices[i]);
33 }
34}
35
37{
38 return _providedData->getWaveformDegree();
39}
40
42{
45 for (auto &context : _mappingContexts) {
46 auto id = context.fromData->getID();
47 if (from.contains(id)) {
48 if (from.toKeep(id)) {
49 context.toData->timeStepsStorage().clearExceptLast();
50 } else {
51 context.toData->timeStepsStorage().clear();
52 }
53 }
54 }
55}
56
58{
61 for (auto &context : _mappingContexts) {
62 if (from.contains(context.fromData->getID())) {
63 context.toData->timeStepsStorage().trimAfter(t);
64 }
65 }
66}
67
68} // namespace precice::impl
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:95
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
Stores one Data object with related mesh.
int getMeshVertexCount() const
Get the number of vertices of mesh.
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.
bool hasMapping() const
Informs the user whether this DataContext has any _mappingContext.
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.
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.
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 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
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.