preCICE v3.1.2
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
13{
14 PRECICE_ASSERT(data);
15 _providedData = data;
16 PRECICE_ASSERT(mesh);
17 _mesh = mesh;
18}
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 return hasReadMapping() || hasWriteMapping();
88}
89
91{
94
95 int executedMappings{0};
96
97 // Execute the mappings
98 for (auto &context : _mappingContexts) {
99 PRECICE_ASSERT(!context.fromData->stamples().empty(),
100 "There must be samples at this point!");
101
102 // linear lookup should be sufficient here
103 const auto timestampExists = [times = context.toData->timeStepsStorage().getTimes()](double lookup) -> bool {
104 return std::any_of(times.data(), std::next(times.data(), times.size()), [lookup](double time) {
105 return math::equals(time, lookup);
106 });
107 };
108
109 auto &mapping = *context.mapping;
110
111 const auto dataDims = context.fromData->getDimensions();
112
113 for (const auto &stample : context.fromData->stamples()) {
114 // skip stamples before given time
115 if (after && math::smallerEquals(stample.timestamp, *after)) {
116 PRECICE_DEBUG("Skipping stample t={} (not after {})", stample.timestamp, *after);
117 continue;
118 }
119 // skip existing stamples
120 if (timestampExists(stample.timestamp)) {
121 PRECICE_DEBUG("Skipping stample t={} (exists)", stample.timestamp);
122 continue;
123 }
124
125 time::Sample outSample{
126 dataDims,
127 Eigen::VectorXd::Zero(dataDims * mapping.getOutputMesh()->nVertices())};
128
129 // Note that the l2norm is only computed during initialization due to short-circuit evaluation in C++
130 bool skipMapping = skipZero && (utils::IntraComm::l2norm(stample.sample.values) < math::NUMERICAL_ZERO_DIFFERENCE);
131
132 PRECICE_INFO("Mapping \"{}\" for t={} from \"{}\" to \"{}\"{}",
133 getDataName(), stample.timestamp, mapping.getInputMesh()->getName(), mapping.getOutputMesh()->getName(),
134 (skipMapping ? " (skipped zero sample)" : ""));
135 if (!skipMapping) {
136 if (mapping.requiresInitialGuess()) {
137 const FromToDataIDs key{context.fromData->getID(), context.toData->getID()};
138 mapping.map(stample.sample, outSample.values, _initialGuesses[key]);
139 } else {
140 mapping.map(stample.sample, outSample.values);
141 }
142 PRECICE_DEBUG("Mapped values (t={}) = {}", stample.timestamp, utils::previewRange(3, outSample.values));
143 ++executedMappings;
144 }
145
146 // Store data from mapping buffer in storage
147 context.toData->setSampleAtTime(stample.timestamp, std::move(outSample));
148 }
149 }
150 return executedMappings;
151}
152
154{
155 return std::any_of(_mappingContexts.begin(), _mappingContexts.end(), [this](auto &context) { return context.toData == _providedData; });
156}
157
159{
160 return std::any_of(_mappingContexts.begin(), _mappingContexts.end(), [this](auto &context) { return context.fromData == _providedData; });
161}
162
164{
166 return _mesh->isValidVertexID(id);
167}
168
169} // namespace precice::impl
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:64
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:95
#define PRECICE_INFO(...)
Definition LogMacros.hpp:13
T any_of(T... args)
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
std::map< FromToDataIDs, Eigen::VectorXd > _initialGuesses
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...
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)
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.
mesh::PtrData _providedData
Unique data this context is associated with.
void resetInitialGuesses()
Resets initial guesses of transient mappings to zero.
static logging::Logger _log
std::vector< MappingContext > _mappingContexts
Defines all mappings associated to this DataContext. A DataContext may also exist without a mapping.
mesh::PtrMesh _mesh
Unique mesh associated with _providedData.
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
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
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