preCICE v3.1.2
Loading...
Searching...
No Matches
BarycentricBaseMapping.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <algorithm>
3#include <memory>
4#include <ostream>
5#include <unordered_set>
6#include <utility>
7
10#include "mapping/Mapping.hpp"
11#include "mapping/Polation.hpp"
12#include "math/differences.hpp"
13#include "mesh/Data.hpp"
14#include "mesh/Mesh.hpp"
16#include "mesh/Vertex.hpp"
17#include "profiling/Event.hpp"
18#include "query/Index.hpp"
19#include "utils/IntraComm.hpp"
20#include "utils/Statistics.hpp"
21#include "utils/assertion.hpp"
22
23namespace precice::mapping {
24
26 : Mapping(constraint, dimensions, false, Mapping::InitialGuessRequirement::None)
27{
28}
29
36
37void BarycentricBaseMapping::mapConservative(const time::Sample &inData, Eigen::VectorXd &outData)
38{
40 precice::profiling::Event e("map.bbm.mapData.From" + input()->getName() + "To" + output()->getName(), profiling::Synchronize);
42 PRECICE_DEBUG("Map conservative using {}", getName());
43 PRECICE_ASSERT(_interpolations.size() == input()->nVertices(),
44 _interpolations.size(), input()->nVertices());
45 const int dimensions = inData.dataDims;
46 const Eigen::VectorXd &inValues = inData.values;
47 Eigen::VectorXd & outValues = outData;
48
49 // For each input vertex, distribute the conserved data among the relevant output vertices
50 // Do it for all dimensions (i.e. components if data is a vector)
51 for (size_t i = 0; i < input()->nVertices(); i++) {
52 const size_t inOffset = i * dimensions;
53 const auto & elems = _interpolations[i].getWeightedElements();
54 for (const auto &elem : elems) {
55 size_t outOffset = static_cast<size_t>(elem.vertexID) * dimensions;
56 for (int dim = 0; dim < dimensions; dim++) {
57 PRECICE_ASSERT(outOffset + dim < (size_t) outValues.size());
58 PRECICE_ASSERT(inOffset + dim < (size_t) inValues.size());
59 outValues(outOffset + dim) += elem.weight * inValues(inOffset + dim);
60 }
61 }
62 }
63}
64
65void BarycentricBaseMapping::mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData)
66{
68 precice::profiling::Event e("map.bbm.mapData.From" + input()->getName() + "To" + output()->getName(), profiling::Synchronize);
69 PRECICE_DEBUG("Map {} using {}", (hasConstraint(CONSISTENT) ? "consistent" : "scaled-consistent"), getName());
70 PRECICE_ASSERT(_interpolations.size() == output()->nVertices(),
71 _interpolations.size(), output()->nVertices());
72
73 const int dimensions = inData.dataDims;
74 const Eigen::VectorXd &inValues = inData.values;
75 Eigen::VectorXd & outValues = outData;
76
77 // For each output vertex, compute the linear combination of input vertices
78 // Do it for all dimensions (i.e. components if data is a vector)
79 for (size_t i = 0; i < output()->nVertices(); i++) {
80 const auto &elems = _interpolations[i].getWeightedElements();
81 size_t outOffset = i * dimensions;
82 for (const auto &elem : elems) {
83 const size_t inOffset = static_cast<size_t>(elem.vertexID) * dimensions;
84 for (int dim = 0; dim < dimensions; dim++) {
85 PRECICE_ASSERT(outOffset + dim < (size_t) outValues.size());
86 PRECICE_ASSERT(inOffset + dim < (size_t) inValues.size());
87 outValues(outOffset + dim) += elem.weight * inValues(inOffset + dim);
88 }
89 }
90 }
91}
92
94{
96 precice::profiling::Event e("map.bbm.tagMeshFirstRound.From" + input()->getName() + "To" + output()->getName(), profiling::Synchronize);
97 PRECICE_DEBUG("Compute Mapping for Tagging");
98
100 PRECICE_DEBUG("Tagging First Round");
101
102 // Determine the Mesh to Tag
103 mesh::PtrMesh origins;
105 origins = output();
106 } else {
107 origins = input();
108 }
109
110 // Gather all vertices to be tagged in a first phase.
111 // max_count is used to shortcut if all vertices have been tagged.
113 const std::size_t max_count = origins->nVertices();
114
115 for (const Polation &interpolation : _interpolations) {
116 for (const auto &elem : interpolation.getWeightedElements()) {
117 if (!math::equals(elem.weight, 0.0)) {
118 tagged.insert(elem.vertexID);
119 }
120 }
121 // Shortcut if all vertices are tagged
122 if (tagged.size() == max_count) {
123 break;
124 }
125 }
126
127 // Now tag all vertices to be tagged in the second phase.
128 for (auto &v : origins->vertices()) {
129 if (tagged.count(v.getID()) == 1) {
130 v.tag();
131 }
132 }
133 PRECICE_DEBUG("First Round Tagged {}/{} Vertices", tagged.size(), max_count);
134
135 clear();
136}
137
139{
141 // for NP mapping no operation needed here
142}
143
144} // namespace precice::mapping
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:64
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:95
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
void tagMeshFirstRound() final override
Method used by partition. Tags vertices that could be owned by this rank.
void clear() final override
Removes a computed mapping.
void mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData) override
Maps data using a consistent constraint.
void tagMeshSecondRound() final override
Method used by partition. Tags vertices that can be filtered out.
BarycentricBaseMapping(Constraint constraint, int dimensions)
void mapConservative(const time::Sample &inData, Eigen::VectorXd &outData) override
Maps data using a conservative constraint.
Abstract base class for mapping of data from one mesh to another.
Definition Mapping.hpp:15
mesh::PtrMesh output() const
Returns pointer to output mesh.
Definition Mapping.cpp:91
Constraint
Specifies additional constraints for a mapping.
Definition Mapping.hpp:29
mesh::PtrMesh input() const
Returns pointer to input mesh.
Definition Mapping.cpp:86
bool _hasComputedMapping
Flag to indicate whether computeMapping() has been called.
Definition Mapping.hpp:208
Constraint getConstraint() const
Returns the constraint (consistent/conservative) of the mapping.
Definition Mapping.cpp:45
virtual std::string getName() const =0
Returns the name of the mapping method for logging purpose.
virtual void computeMapping()=0
Computes the mapping coefficients from the in- and output mesh.
virtual bool hasConstraint(const Constraint &constraint) const
Checks whether the mapping has the given constraint or not.
Definition Mapping.cpp:247
InitialGuessRequirement
Specifies whether the mapping requires an initial guess.
Definition Mapping.hpp:63
Calculates the barycentric coordinates of a coordinate on the given vertex/edge/triangle and stores t...
Definition Polation.hpp:24
T count(T... args)
T insert(T... args)
contains data mapping from points to meshes.
constexpr bool equals(const Eigen::MatrixBase< DerivedA > &A, const Eigen::MatrixBase< DerivedB > &B, double tolerance=NUMERICAL_ZERO_DIFFERENCE)
Compares two Eigen::MatrixBase for equality up to tolerance.
static constexpr SynchronizeTag Synchronize
Convenience instance of the SynchronizeTag.
Definition Event.hpp:21
T size(T... args)
int dataDims
The dimensionality of the data.
Definition Sample.hpp:45
Eigen::VectorXd values
Definition Sample.hpp:49