preCICE v3.3.1
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>
9#include "profiling/Event.hpp"
11#include "utils/IntraComm.hpp"
12#include "utils/assertion.hpp"
13
14namespace precice::mapping {
15
17 Constraint constraint,
18 int dimensions)
19 : NearestNeighborBaseMapping(constraint, dimensions, true, "NearestNeighborGradientMapping", "nng")
20{
22
24 "The scaled-consistent mapping hasn't been specifically tested with nearest-neighbor-gradient. Please avoid using it or choose another mapping method. ");
25
26 if (isScaledConsistent()) {
29 } else {
32 }
33}
34
36{
37
38 // Initialize the offsets list
39 _offsetsMatched.resize(_vertexIndices.size());
40
41 // Calculate offsets
42 for (size_t i = 0; i < _vertexIndices.size(); ++i) {
43
44 const auto &matchedVertexCoords = searchSpace->vertex(_vertexIndices[i]).getCoords();
45 const auto &sourceVertexCoords = origins->vertex(i).getCoords();
46
47 // We calculate the distances uniformly for consistent mapping constraint as the difference (output - input)
48 // For consistent mapping: the source is the output vertex and the matched vertex is the input since we iterate over all outputs
49 // and assign each exactly one vertex form the search space, which are our origins vertices.
50 _offsetsMatched[i] = sourceVertexCoords - matchedVertexCoords;
51 }
52};
53
54void NearestNeighborGradientMapping::mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData)
55{
57 precice::profiling::Event e("map." + mappingNameShort + ".mapData.From" + input()->getName() + "To" + output()->getName(), profiling::Synchronize);
58
59 PRECICE_ASSERT(inData.values.size() == 0 || inData.gradients.size() != 0,
60 "Mesh \"{}\" does not contain gradient data. Using Nearest Neighbor Gradient mapping requires gradient data.",
61 input()->getName());
62
64 PRECICE_WARN_IF(input()->empty(), "The mesh doesn't contain any vertices.");
65
66 const int valueDimensions = inData.dataDims;
67 const Eigen::VectorXd &inputValues = inData.values;
68 Eigen::VectorXd &outputValues = outData;
69 const Eigen::MatrixXd &gradients = inData.gradients;
70
71 // Consistent mapping
72 PRECICE_DEBUG("Map {} using {}", (hasConstraint(CONSISTENT) ? "consistent" : "scaled-consistent"), getName());
73 const size_t outSize = output()->nVertices();
74
75 for (size_t i = 0; i < outSize; i++) {
76 int inputIndex = _vertexIndices[i] * valueDimensions;
77
78 for (int dim = 0; dim < valueDimensions; dim++) {
79
80 const int mapOutputIndex = (i * valueDimensions) + dim;
81 const int mapInputIndex = inputIndex + dim;
82
83 outputValues(mapOutputIndex) = inputValues(mapInputIndex) + _offsetsMatched[i].transpose() * gradients.col(mapInputIndex);
84 }
85 }
86
87 PRECICE_DEBUG("Mapped values (with gradient) = {}", utils::previewRange(3, outputValues));
88}
89
90void NearestNeighborGradientMapping::mapConservative(const time::Sample & /* inData */, Eigen::VectorXd & /* outData */)
91{
92 PRECICE_ASSERT(false, "Not implemented.");
93}
94
96{
97 return "nearest-neighbor-gradient";
98}
99
100} // namespace precice::mapping
#define PRECICE_WARN_IF(condition,...)
Definition LogMacros.hpp:18
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:61
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
std::size_t nVertices() const
Returns the number of vertices.
Definition Mesh.cpp:65
Vertex & vertex(VertexID id)
Mutable access to a vertex by VertexID.
Definition Mesh.cpp:43
mesh::PtrMesh output() const
Returns pointer to output mesh.
Definition Mapping.cpp:92
Constraint
Specifies additional constraints for a mapping.
Definition Mapping.hpp:30
mesh::PtrMesh input() const
Returns pointer to input mesh.
Definition Mapping.cpp:87
bool isScaledConsistent() const
Returns true if mapping is a form of scaled consistent mapping.
Definition Mapping.cpp:258
void setInputRequirement(MeshRequirement requirement)
Sets the mesh requirement for the input mesh.
Definition Mapping.cpp:97
void setOutputRequirement(MeshRequirement requirement)
Sets the mesh requirement for the output mesh.
Definition Mapping.cpp:103
virtual bool hasConstraint(const Constraint &constraint) const
Checks whether the mapping has the given constraint or not.
Definition Mapping.cpp:248
std::vector< int > _vertexIndices
Computed output vertex indices to map data from input vertices to.
NearestNeighborBaseMapping(Constraint constraint, int dimensions, bool hasGradient, std::string mappingName, std::string mappingNameShort)
Constructor.
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
Eigen::VectorXd getCoords() const
Returns the coordinates of the vertex.
Definition Vertex.hpp:114
contains data mapping from points to meshes.
std::shared_ptr< Mesh > PtrMesh
static constexpr SynchronizeTag Synchronize
Convenience instance of the SynchronizeTag.
Definition Event.hpp:28
const RangePreview< Iter > previewRange(Size n, const Range &range)
Eigen::MatrixXd gradients
The gradients of the data. Use gradients.col(d*i+k) to get the gradient of vertex i,...
Definition Sample.hpp:67
int dataDims
The dimensionality of the data.
Definition Sample.hpp:60
Eigen::VectorXd values
Definition Sample.hpp:64