/
NearestNeighborMapping.cpp
88 lines (69 loc) · 2.99 KB
/
NearestNeighborMapping.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
#include "NearestNeighborMapping.hpp"
#include <Eigen/Core>
#include <boost/container/flat_set.hpp>
#include <functional>
#include "logging/LogMacros.hpp"
#include "profiling/Event.hpp"
#include "profiling/EventUtils.hpp"
#include "utils/EigenHelperFunctions.hpp"
#include "utils/IntraComm.hpp"
#include "utils/assertion.hpp"
namespace precice::mapping {
NearestNeighborMapping::NearestNeighborMapping(
Constraint constraint,
int dimensions)
: NearestNeighborBaseMapping(constraint, dimensions, false, "NearestNeighborMapping", "nn")
{
if (isScaledConsistent()) {
setInputRequirement(Mapping::MeshRequirement::FULL);
setOutputRequirement(Mapping::MeshRequirement::FULL);
} else {
setInputRequirement(Mapping::MeshRequirement::VERTEX);
setOutputRequirement(Mapping::MeshRequirement::VERTEX);
}
}
void NearestNeighborMapping::mapConservative(const time::Sample &inData, Eigen::VectorXd &outData)
{
PRECICE_TRACE();
precice::profiling::Event e("map." + mappingNameShort + ".mapData.From" + input()->getName() + "To" + output()->getName(), profiling::Synchronize);
PRECICE_DEBUG("Map conservative using {}", getName());
const Eigen::VectorXd &inputValues = inData.values;
Eigen::VectorXd & outputValues = outData;
// Data dimensions (for scalar = 1, for vectors > 1)
const size_t inSize = input()->nVertices();
const int valueDimensions = inData.dataDims;
for (size_t i = 0; i < inSize; i++) {
int const outputIndex = _vertexIndices[i] * valueDimensions;
for (int dim = 0; dim < valueDimensions; dim++) {
const int mapOutputIndex = outputIndex + dim;
const int mapInputIndex = (i * valueDimensions) + dim;
outputValues(mapOutputIndex) += inputValues(mapInputIndex);
}
}
PRECICE_DEBUG("Mapped values = {}", utils::previewRange(3, outputValues));
}
void NearestNeighborMapping::mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData)
{
PRECICE_TRACE();
precice::profiling::Event e("map." + mappingNameShort + ".mapData.From" + input()->getName() + "To" + output()->getName(), profiling::Synchronize);
PRECICE_DEBUG("Map {} using {}", (hasConstraint(CONSISTENT) ? "consistent" : "scaled-consistent"), getName());
const Eigen::VectorXd &inputValues = inData.values;
Eigen::VectorXd & outputValues = outData;
// Data dimensions (for scalar = 1, for vectors > 1)
const size_t outSize = output()->nVertices();
const int valueDimensions = inData.dataDims;
for (size_t i = 0; i < outSize; i++) {
int inputIndex = _vertexIndices[i] * valueDimensions;
for (int dim = 0; dim < valueDimensions; dim++) {
const int mapOutputIndex = (i * valueDimensions) + dim;
const int mapInputIndex = inputIndex + dim;
outputValues(mapOutputIndex) = inputValues(mapInputIndex);
}
}
PRECICE_DEBUG("Mapped values = {}", utils::previewRange(3, outputValues));
}
std::string NearestNeighborMapping::getName() const
{
return "nearest-neighbor";
}
} // namespace precice::mapping