preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
NearestNeighborMapping.cpp
Go to the documentation of this file.
2
3#include <Eigen/Core>
4#include <boost/container/flat_set.hpp>
5#include <functional>
7#include "profiling/Event.hpp"
10#include "utils/IntraComm.hpp"
11#include "utils/assertion.hpp"
12
13namespace precice::mapping {
14
28
29void NearestNeighborMapping::mapConsistentAt(const Eigen::Ref<const Eigen::MatrixXd> &coordinates, const impl::MappingDataCache &cache, Eigen::Ref<Eigen::MatrixXd> values)
30{
31 precice::profiling::Event e("map.nn.mapConsistentAt.From" + input()->getName());
32 auto &index = input()->index();
33
34 // Set up of output arrays
35 Eigen::Map<const Eigen::MatrixXd> localData(cache.inData.sample.values.data(), cache.getDataDimensions(), cache.inData.sample.values.size() / cache.getDataDimensions());
36 for (Eigen::Index i = 0; i < coordinates.cols(); ++i) {
37 values.col(i) = localData.col(index.getClosestVertex(coordinates.col(i)).index);
38 }
39}
40
41void NearestNeighborMapping::mapConservativeAt(const Eigen::Ref<const Eigen::MatrixXd> &coordinates, const Eigen::Ref<const Eigen::MatrixXd> &source, impl::MappingDataCache &, Eigen::Ref<Eigen::MatrixXd> target)
42{
43 precice::profiling::Event e("map.nn.mapConservativeAt.From" + input()->getName());
44 auto &index = output()->index();
45 for (Eigen::Index i = 0; i < coordinates.cols(); ++i) {
46 target.col(index.getClosestVertex(coordinates.col(i)).index) += source.col(i);
47 }
48}
49
50void NearestNeighborMapping::mapConservative(const time::Sample &inData, Eigen::VectorXd &outData)
51{
53 precice::profiling::Event e("map." + mappingNameShort + ".mapData.From" + input()->getName() + "To" + output()->getName(), profiling::Synchronize);
54 PRECICE_DEBUG("Map conservative using {}", getName());
55
56 // scalar = 1, vector > 1
57 const int valueDimensions = inData.dataDims;
58 // Create Eigen::Map to interpret data as a matrix (enables col-wise operations)
59 Eigen::Map<const Eigen::MatrixXd> inMap(inData.values.data(), valueDimensions, input()->nVertices());
60 Eigen::Map<Eigen::MatrixXd> outMap(outData.data(), valueDimensions, output()->nVertices());
61
62 // Apply mapping
63 for (Eigen::Index i = 0; i < inMap.cols(); ++i) {
64 outMap.col(_vertexIndices[i]) += inMap.col(i);
65 }
66
67 PRECICE_DEBUG("Mapped values = {}", utils::previewRange(3, outData));
68}
69
70void NearestNeighborMapping::mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData)
71{
73 precice::profiling::Event e("map." + mappingNameShort + ".mapData.From" + input()->getName() + "To" + output()->getName(), profiling::Synchronize);
74 PRECICE_DEBUG("Map {} using {}", (hasConstraint(CONSISTENT) ? "consistent" : "scaled-consistent"), getName());
75
76 const int valueDimensions = inData.dataDims;
77 // Map the raw pointers as (valueDimensions x vertices) in column-major Eigen style
78 Eigen::Map<const Eigen::MatrixXd> inputMap(inData.values.data(), valueDimensions, input()->nVertices());
79 Eigen::Map<Eigen::MatrixXd> outputMap(outData.data(), valueDimensions, output()->nVertices());
80
81 for (Eigen::Index i = 0; i < outputMap.cols(); ++i) {
82 // _vertexIndices[i] is the solution of our mapping
83 outputMap.col(i) = inputMap.col(_vertexIndices[i]);
84 }
85 PRECICE_DEBUG("Mapped values = {}", utils::previewRange(3, outData));
86}
87
89{
90 return "nearest-neighbor";
91}
92
93} // namespace precice::mapping
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:61
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
unsigned int index
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.
void mapConservativeAt(const Eigen::Ref< const Eigen::MatrixXd > &coordinates, const Eigen::Ref< const Eigen::MatrixXd > &source, impl::MappingDataCache &cache, Eigen::Ref< Eigen::MatrixXd > target) final override
For writing data just-in-time (only conservative at the moment)
void mapConsistentAt(const Eigen::Ref< const Eigen::MatrixXd > &coordinates, const impl::MappingDataCache &cache, Eigen::Ref< Eigen::MatrixXd > values) final override
For reading data just-in-time (only consistent at the moment)
std::string getName() const final override
name of the nn mapping
void mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData) final override
Maps data using a consistent constraint.
NearestNeighborMapping(Constraint constraint, int dimensions)
Constructor.
void mapConservative(const time::Sample &inData, Eigen::VectorXd &outData) final override
Maps data using a conservative constraint.
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)
int getDataDimensions() const
Returns the number of data components.
int dataDims
The dimensionality of the data.
Definition Sample.hpp:60
Eigen::VectorXd values
Definition Sample.hpp:64