preCICE v3.1.2
Loading...
Searching...
No Matches
NearestNeighborGradientMapping.cpp
Go to the documentation of this file.
2
3#include <Eigen/Core>
4#include <Eigen/src/Core/Matrix.h>
5#include <boost/container/flat_set.hpp>
6#include <functional>
7#include <iostream>
8#include <strings.h>
10#include "profiling/Event.hpp"
12#include "utils/IntraComm.hpp"
13#include "utils/assertion.hpp"
14
15namespace precice::mapping {
16
18 Constraint constraint,
19 int dimensions)
20 : NearestNeighborBaseMapping(constraint, dimensions, true, "NearestNeighborGradientMapping", "nng")
21{
23
25 "The scaled-consistent mapping hasn't been specifically tested with nearest-neighbor-gradient. Please avoid using it or choose another mapping method. ");
26
27 if (isScaledConsistent()) {
30 } else {
33 }
34}
35
37{
38
39 // Initialize the offsets list
41
42 // Calculate offsets
43 for (size_t i = 0; i < _vertexIndices.size(); ++i) {
44
45 const auto &matchedVertexCoords = searchSpace->vertex(_vertexIndices[i]).getCoords();
46 const auto &sourceVertexCoords = origins->vertex(i).getCoords();
47
48 // We calculate the distances uniformly for consistent mapping constraint as the difference (output - input)
49 // For consistent mapping: the source is the output vertex and the matched vertex is the input since we iterate over all outputs
50 // and assign each exactly one vertex form the search space, which are our origins vertices.
51 _offsetsMatched[i] = sourceVertexCoords - matchedVertexCoords;
52 }
53};
54
55void NearestNeighborGradientMapping::mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData)
56{
58 precice::profiling::Event e("map." + mappingNameShort + ".mapData.From" + input()->getName() + "To" + output()->getName(), profiling::Synchronize);
59
60 PRECICE_ASSERT(inData.values.size() == 0 || inData.gradients.size() != 0,
61 "Mesh \"{}\" does not contain gradient data. Using Nearest Neighbor Gradient mapping requires gradient data.",
62 input()->getName());
63
65 PRECICE_WARN_IF(input()->empty(), "The mesh doesn't contain any vertices.");
66
67 const int valueDimensions = inData.dataDims;
68 const Eigen::VectorXd &inputValues = inData.values;
69 Eigen::VectorXd & outputValues = outData;
70 const Eigen::MatrixXd &gradients = inData.gradients;
71
72 // Consistent mapping
73 PRECICE_DEBUG("Map {} using {}", (hasConstraint(CONSISTENT) ? "consistent" : "scaled-consistent"), getName());
74 const size_t outSize = output()->nVertices();
75
76 for (size_t i = 0; i < outSize; i++) {
77 int inputIndex = _vertexIndices[i] * valueDimensions;
78
79 for (int dim = 0; dim < valueDimensions; dim++) {
80
81 const int mapOutputIndex = (i * valueDimensions) + dim;
82 const int mapInputIndex = inputIndex + dim;
83
84 outputValues(mapOutputIndex) = inputValues(mapInputIndex) + _offsetsMatched[i].transpose() * gradients.col(mapInputIndex);
85 }
86 }
87
88 PRECICE_DEBUG("Mapped values (with gradient) = {}", utils::previewRange(3, outputValues));
89}
90
91void NearestNeighborGradientMapping::mapConservative(const time::Sample & /* inData */, Eigen::VectorXd & /* outData */)
92{
93 PRECICE_ASSERT(false, "Not implemented.");
94}
95
97{
98 return "nearest-neighbor-gradient";
99}
100
101} // namespace precice::mapping
#define PRECICE_WARN_IF(condition,...)
Definition LogMacros.hpp:21
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:64
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:95
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
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 isScaledConsistent() const
Returns true if mapping is a form of scaled consistent mapping.
Definition Mapping.cpp:257
void setInputRequirement(MeshRequirement requirement)
Sets the mesh requirement for the input mesh.
Definition Mapping.cpp:96
void setOutputRequirement(MeshRequirement requirement)
Sets the mesh requirement for the output mesh.
Definition Mapping.cpp:102
virtual bool hasConstraint(const Constraint &constraint) const
Checks whether the mapping has the given constraint or not.
Definition Mapping.cpp:247
std::vector< int > _vertexIndices
Computed output vertex indices to map data from input vertices to.
NearestNeighborGradientMapping(Constraint constraint, int dimensions)
Constructor.
void mapConservative(const time::Sample &inData, Eigen::VectorXd &outData) final override
Maps data using a conservative constraint.
void onMappingComputed(mesh::PtrMesh origins, mesh::PtrMesh searchSpace) final override
Calculates the offsets needed for the gradient mappings after calculating the matched vertices.
void mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData) final override
Maps data using a consistent constraint.
std::string getName() const final override
name of the nng mapping
contains data mapping from points to meshes.
static constexpr SynchronizeTag Synchronize
Convenience instance of the SynchronizeTag.
Definition Event.hpp:21
const RangePreview< Iter > previewRange(Size n, const Range &range)
T resize(T... args)
T size(T... args)
Eigen::MatrixXd gradients
The gradients of the data. Use gradients.col(i) to get the gradient at vertex i.
Definition Sample.hpp:52
int dataDims
The dimensionality of the data.
Definition Sample.hpp:45
Eigen::VectorXd values
Definition Sample.hpp:49