preCICE v3.2.0
Loading...
Searching...
No Matches
DataContext.cpp
Go to the documentation of this file.
1#include <iterator>
2#include <memory>
3#include <utility>
4
7#include "utils/IntraComm.hpp"
8namespace precice::impl {
9
10logging::Logger DataContext::_log{"impl::DataContext"};
11
19
25
27{
28 for (auto &kv : _initialGuesses) {
29 kv.second.setZero();
30 }
31}
32
34{
36 return _providedData->getDimensions();
37}
38
40{
42 return _providedData->getSpatialDimensions();
43}
44
46{
48 return _mesh->getName();
49}
50
52{
54 return _mesh->nVertices();
55}
56
58{
60 return _mesh->getID();
61}
62
64{
66 return _providedData->hasGradient();
67}
68
70{
71 PRECICE_ASSERT(mappingContext.fromData);
72 PRECICE_ASSERT(mappingContext.toData);
73 // Make sure we don't append a mapping twice
74#ifndef NDEBUG
75 for (auto &context : _mappingContexts) {
76 PRECICE_ASSERT(!((context.mapping == mappingContext.mapping) && (context.fromData == mappingContext.fromData) && (context.fromData == mappingContext.toData)), "The appended mapping already exists.");
77 }
78#endif
79 _mappingContexts.emplace_back(mappingContext);
80 PRECICE_ASSERT(mappingContext.fromData == _providedData || mappingContext.toData == _providedData, "Either fromData or toData has to equal _providedData.");
81 PRECICE_ASSERT(mappingContext.fromData->getName() == getDataName());
82 PRECICE_ASSERT(mappingContext.toData->getName() == getDataName());
83}
84
86{
87 PRECICE_ASSERT(meshContext.mesh->hasDataName(getDataName()));
89 PRECICE_ASSERT(mappingCache == nullptr);
90 mesh::PtrData data = meshContext.mesh->data(getDataName());
91 // the mapping itself has even for just-in-time mapping no notion about the data
92 // maybe remove the data pointer here or set them to nullptr
93 // the data access happens through the API functions
94 mappingContext.toData = data;
95 mappingContext.fromData = data;
96
97 PRECICE_ASSERT(mappingContext.fromData);
98 PRECICE_ASSERT(mappingContext.toData);
99
100 PRECICE_ASSERT(mappingContext.fromData == _providedData || mappingContext.toData == _providedData, "Either fromData or toData has to equal _providedData.");
102 justInTimeMapping = mappingContext.mapping;
103}
104
106{
107 return hasReadMapping() || hasWriteMapping();
108}
109
111{
114
115 int executedMappings{0};
116
117 // Execute the mappings
118 for (auto &context : _mappingContexts) {
119 PRECICE_CHECK(!context.fromData->stamples().empty(),
120 "Data {0} on mesh {1} didn't contain any data samples while attempting to map to mesh {2}. "
121 "Check your exchange tags to ensure your coupling scheme exchanges the data or the participant produces it using an action. "
122 "The expected exchange tag should look like this: <exchange data=\"{0}\" mesh=\"{1}\" from=... to=... />.",
123 context.fromData->getName(), context.mapping->getInputMesh()->getName(), context.mapping->getOutputMesh()->getName());
124
125 // linear lookup should be sufficient here
126 const auto timestampExists = [times = context.toData->timeStepsStorage().getTimes()](double lookup) -> bool {
127 return std::any_of(times.data(), std::next(times.data(), times.size()), [lookup](double time) {
128 return math::equals(time, lookup);
129 });
130 };
131
132 auto &mapping = *context.mapping;
133
134 const auto dataDims = context.fromData->getDimensions();
135
136 for (const auto &stample : context.fromData->stamples()) {
137 // skip stamples before given time
138 if (after && math::smallerEquals(stample.timestamp, *after)) {
139 PRECICE_DEBUG("Skipping stample t={} (not after {})", stample.timestamp, *after);
140 continue;
141 }
142 // skip existing stamples
143 if (timestampExists(stample.timestamp)) {
144 PRECICE_DEBUG("Skipping stample t={} (exists)", stample.timestamp);
145 continue;
146 }
147
148 time::Sample outSample{
149 dataDims,
150 Eigen::VectorXd::Zero(dataDims * mapping.getOutputMesh()->nVertices())};
151
152 // Note that the l2norm is only computed during initialization due to short-circuit evaluation in C++
153 bool skipMapping = skipZero && (utils::IntraComm::l2norm(stample.sample.values) < math::NUMERICAL_ZERO_DIFFERENCE);
154
155 PRECICE_INFO("Mapping \"{}\" for t={} from \"{}\" to \"{}\"{}",
156 getDataName(), stample.timestamp, mapping.getInputMesh()->getName(), mapping.getOutputMesh()->getName(),
157 (skipMapping ? " (skipped zero sample)" : ""));
158 if (!skipMapping) {
159 if (mapping.requiresInitialGuess()) {
160 const FromToDataIDs key{context.fromData->getID(), context.toData->getID()};
161 mapping.map(stample.sample, outSample.values, _initialGuesses[key]);
162 } else {
163 mapping.map(stample.sample, outSample.values);
164 }
165 PRECICE_DEBUG("Mapped values (t={}) = {}", stample.timestamp, utils::previewRange(3, outSample.values));
166 ++executedMappings;
167 }
168
169 // Store data from mapping buffer in storage
170 context.toData->setSampleAtTime(stample.timestamp, std::move(outSample));
171 }
172 }
173 return executedMappings;
174}
175
177{
178 return std::any_of(_mappingContexts.begin(), _mappingContexts.end(), [this](auto &context) { return context.toData == _providedData; });
179}
180
182{
183 return std::any_of(_mappingContexts.begin(), _mappingContexts.end(), [this](auto &context) { return context.fromData == _providedData; });
184}
185
187{
189 return _mesh->isValidVertexID(id);
190}
191
192} // namespace precice::impl
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:61
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_INFO(...)
Definition LogMacros.hpp:14
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:32
T any_of(T... args)
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
std::map< FromToDataIDs, Eigen::VectorXd > _initialGuesses
int getMeshVertexCount() const
Get the number of vertices of mesh.
std::pair< int, int > FromToDataIDs
bool hasWriteMapping() const
Informs the user whether this DataContext has any write mapping.
void addJustInTimeMapping(MappingContext &mappingContext, MeshContext &meshContext)
Attach a just-in-time mapping to this data context and setup a corresponding MappingDataCache.
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.
MeshID getMeshID() const
Get the ID 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.
bool hasGradient() const
Returns whether _providedData has gradient.
int getDataDimensions() const
Get the dimensions of _providedData.
int mapData(std::optional< double > after=std::nullopt, bool skipZero=false)
Perform the mapping for mapping contexts and the corresponding data context (from and to data)
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.
bool isValidVertexID(const VertexID id) const
Returns true if the given vertexID is valid.
void resetInitialGuesses()
Resets initial guesses of transient mappings to zero.
static logging::Logger _log
Unique mesh associated with _providedData.
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.
static double l2norm(const Eigen::VectorXd &vec)
The l2 norm of a vector is calculated on distributed data.
Definition IntraComm.cpp:67
T make_unique(T... args)
contains data mapping from points to meshes.
std::enable_if< std::is_arithmetic< Scalar >::value, bool >::type smallerEquals(Scalar A, Scalar B, Scalar tolerance=NUMERICAL_ZERO_DIFFERENCE)
constexpr double NUMERICAL_ZERO_DIFFERENCE
provides Mesh, Data and primitives.
std::shared_ptr< Data > PtrData
std::shared_ptr< Mesh > PtrMesh
contains the time interpolation logic.
Definition Sample.hpp:8
const RangePreview< Iter > previewRange(Size n, const Range &range)
int MeshID
Definition Types.hpp:30
int VertexID
Definition Types.hpp:13
T next(T... args)
Holds a data mapping and related information.
mesh::PtrData fromData
data which is mapped from mesh
mapping::PtrMapping mapping
Data mapping.
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::VectorXd values
Definition Sample.hpp:64