preCICE v3.1.2
Loading...
Searching...
No Matches
NearestNeighborBaseMapping.cpp
Go to the documentation of this file.
2
3#include <boost/container/flat_set.hpp>
4#include <functional>
5#include <iostream>
7#include "mapping/Mapping.hpp"
9#include "mesh/Vertex.hpp"
10#include "profiling/Event.hpp"
11#include "utils/IntraComm.hpp"
12#include "utils/Parallel.hpp"
13#include "utils/Statistics.hpp"
14#include "utils/assertion.hpp"
15
16namespace precice::mapping {
17
19 Constraint constraint,
20 int dimensions,
21 bool requiresGradientData,
22 std::string mappingName,
23 std::string mappingNameShort)
24 : Mapping(constraint, dimensions, requiresGradientData, Mapping::InitialGuessRequirement::None),
25 mappingName(std::move(mappingName)),
26 mappingNameShort(std::move(mappingNameShort))
27{
28}
29
31{
32 PRECICE_TRACE(input()->nVertices());
33
34 PRECICE_ASSERT(input().get() != nullptr);
35 PRECICE_ASSERT(output().get() != nullptr);
36
37 const std::string baseEvent = "map." + mappingNameShort + ".computeMapping.From" + input()->getName() + "To" + output()->getName();
39
40 // Setup Direction of Mapping
41 mesh::PtrMesh origins, searchSpace;
43 PRECICE_DEBUG("Compute conservative mapping");
44 origins = input();
45 searchSpace = output();
46 } else {
47 PRECICE_DEBUG((hasConstraint(CONSISTENT) ? "Compute consistent mapping" : "Compute scaled-consistent mapping"));
48 origins = output();
49 searchSpace = input();
50 }
51
52 // Set up of output arrays
53 const size_t verticesSize = origins->nVertices();
54 const auto & sourceVertices = origins->vertices();
55 _vertexIndices.resize(verticesSize);
56
57 // Needed for error calculations
59
60 auto &index = searchSpace->index();
61 for (size_t i = 0; i < verticesSize; ++i) {
62 const auto &sourceCoords = sourceVertices[i].getCoords();
63 const auto &matchedVertex = index.getClosestVertex(sourceCoords);
64 _vertexIndices[i] = matchedVertex.index;
65
66 // Compute distance between input and output vertiex for the stats
67 const auto &matchCoords = searchSpace->vertex(matchedVertex.index).getCoords();
68 auto distance = (sourceCoords - matchCoords).norm();
69 distanceStatistics(distance);
70 }
71
72 // For gradient mapping, the calculation of offsets between source and matched vertex necessary
73 onMappingComputed(origins, searchSpace);
74
75 // This is the distance object between the coordinates of the vertices and its match in the mesh.
76 // This prints min, max, average and count of the distances.
77 if (distanceStatistics.empty()) {
78 PRECICE_INFO("Mapping distance not available due to empty partition.");
79 } else {
80 PRECICE_INFO("Mapping distance {}", distanceStatistics);
81 }
82
84}
85
87{
90 _hasComputedMapping = false;
91
94
95 if (getConstraint() == CONSISTENT) {
96 input()->index().clear();
97 } else {
98 output()->index().clear();
99 }
100}
101
103{
104 // Does nothing by default
105}
106
108{
110 precice::profiling::Event e("map." + mappingNameShort + ".tagMeshFirstRound.From" + input()->getName() + "To" + output()->getName(), profiling::Synchronize);
111
113
114 // Lookup table of all indices used in the mapping
115 const boost::container::flat_set<int> indexSet(_vertexIndices.begin(), _vertexIndices.end());
116
117 // Get the source mesh depending on the constraint
118 const mesh::PtrMesh &source = hasConstraint(CONSERVATIVE) ? output() : input();
119
120 // Tag all vertices used in the mapping
121 for (mesh::Vertex &v : source->vertices()) {
122 if (indexSet.count(v.getID()) != 0) {
123 v.tag();
124 }
125 }
126
127 clear();
128}
129
131{
133 // for NN mapping no operation needed here
134}
135
136} // namespace precice::mapping
std::size_t distance
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:64
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:95
#define PRECICE_INFO(...)
Definition LogMacros.hpp:13
unsigned int index
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
T begin(T... args)
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
bool requiresGradientData() const
Returns whether the mapping requires gradient data.
Definition Mapping.cpp:113
virtual std::string getName() const =0
Returns the name of the mapping method for logging purpose.
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
void clear() final override
Removes a computed mapping.
virtual void onMappingComputed(mesh::PtrMesh origins, mesh::PtrMesh searchSpace)
std::vector< int > _vertexIndices
Computed output vertex indices to map data from input vertices to.
void tagMeshSecondRound() final override
Method used by partition. Tags vertices that can be filtered out.
void computeMapping() final override
Computes the mapping coefficients from the in- and output mesh.
NearestNeighborBaseMapping(Constraint constraint, int dimensions, bool hasGradient, std::string mappingName, std::string mappingNameShort)
Constructor.
void tagMeshFirstRound() final override
Method used by partition. Tags vertices that could be owned by this rank.
Vertex of a mesh.
Definition Vertex.hpp:16
bool empty() const
Returns count == 0.
T clear(T... args)
T end(T... args)
contains data mapping from points to meshes.
static constexpr SynchronizeTag Synchronize
Convenience instance of the SynchronizeTag.
Definition Event.hpp:21
constexpr auto get(span< E, S > s) -> decltype(s[N])
Definition span.hpp:602
STL namespace.
T resize(T... args)