preCICE v3.1.2
Loading...
Searching...
No Matches
NearestNeighborGradientMappingTest.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <algorithm>
4#include "mapping/Mapping.hpp"
6#include "math/constants.hpp"
7#include "mesh/Mesh.hpp"
9#include "mesh/Vertex.hpp"
11#include "testing/Testing.hpp"
12
13using namespace precice;
14using namespace precice::mesh;
15
16BOOST_AUTO_TEST_SUITE(MappingTests)
17BOOST_AUTO_TEST_SUITE(NearestNeighborGradientMapping)
18
19BOOST_AUTO_TEST_CASE(ConsistentNonIncremental)
20{
21 PRECICE_TEST(1_rank)
22 int dimensions = 2;
23 using testing::equals;
24
25 // Create mesh to map from
26 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
27 Vertex &inVertex0 = inMesh->createVertex(Eigen::Vector2d::Constant(0.0));
28 Vertex &inVertex1 = inMesh->createVertex(Eigen::Vector2d::Constant(1.0));
29
30 // Create data
31 Eigen::VectorXd inValuesScalar = Eigen::VectorXd::Zero(2);
32 Eigen::VectorXd inValuesVector = Eigen::VectorXd::Zero(4);
33 inValuesScalar << 1.0, 2.0;
34 inValuesVector << 1.0, 2.0, 3.0, 4.0;
35
36 // Create corresponding gradient data (all gradient values = const = 1)
37 Eigen::MatrixXd inGradientsScalar(dimensions, 2);
38 Eigen::MatrixXd inGradientsVector(dimensions, 4);
39 inGradientsScalar.setOnes();
40 inGradientsVector.setOnes();
41
42 // Create mesh to map to
43 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
44 Vertex &outVertex0 = outMesh->createVertex(Eigen::Vector2d::Constant(0.0));
45 Vertex &outVertex1 = outMesh->createVertex(Eigen::Vector2d::Constant(1.0));
46
47 // Setup mapping with mapping coordinates and geometry used
49 mapping.setMeshes(inMesh, outMesh);
50 BOOST_TEST(mapping.hasComputedMapping() == false);
51
52 // Map data with coinciding vertices, has to result in equal values.
53 // Distance between in and out vertices is zero
54 mapping.computeMapping();
55 Eigen::VectorXd outValuesScalar = Eigen::VectorXd::Zero(2);
56 time::Sample inSampleScalar(1, inValuesScalar, inGradientsScalar);
57 mapping.map(inSampleScalar, outValuesScalar);
58 BOOST_TEST(mapping.hasComputedMapping() == true);
59
60 BOOST_TEST(outValuesScalar(0) == inValuesScalar(0));
61 BOOST_TEST(outValuesScalar(1) == inValuesScalar(1));
62
63 Eigen::VectorXd outValuesVector = Eigen::VectorXd::Zero(4);
64 time::Sample inSampleVector(2, inValuesVector, inGradientsVector);
65 mapping.map(inSampleVector, outValuesVector);
66 BOOST_CHECK(equals(inValuesVector, outValuesVector));
67
68 // Map data with almost coinciding vertices, with a null gradient, has to result in equal values
69 inGradientsScalar.setZero();
70 inGradientsVector.setZero();
71
72 outVertex0.setCoords(inVertex0.getCoords() + Eigen::Vector2d::Constant(0.1));
73 outVertex1.setCoords(inVertex1.getCoords() + Eigen::Vector2d::Constant(0.1));
74
75 mapping.computeMapping();
76 inSampleScalar = time::Sample(1, inValuesScalar, inGradientsScalar);
77 mapping.map(inSampleScalar, outValuesScalar);
78 BOOST_TEST(mapping.hasComputedMapping() == true);
79 BOOST_TEST(outValuesScalar(0) == inValuesScalar(0));
80 BOOST_TEST(outValuesScalar(1) == inValuesScalar(1));
81 inSampleVector = time::Sample(2, inValuesVector, inGradientsVector);
82 mapping.map(inSampleVector, outValuesVector);
83 Eigen::Vector4d expected(1.0, 2.0, 3.0, 4.0);
84 BOOST_CHECK(equals(expected, outValuesVector));
85
86 // Map data with almost coinciding vertices, should be a little different with the gradient optimization
87 inGradientsScalar.setOnes();
88 inGradientsVector.setOnes();
89 outVertex0.setCoords(inVertex0.getCoords() + Eigen::Vector2d::Constant(0.1));
90 outVertex1.setCoords(inVertex1.getCoords() + Eigen::Vector2d::Constant(0.1));
91
92 mapping.computeMapping();
93 inSampleScalar = time::Sample(1, inValuesScalar, inGradientsScalar);
94 mapping.map(inSampleScalar, outValuesScalar);
95 BOOST_TEST(mapping.hasComputedMapping() == true);
96 BOOST_TEST(outValuesScalar(0) == inValuesScalar(0) + 0.2);
97 BOOST_TEST(outValuesScalar(1) == inValuesScalar(1) + 0.2);
98 inSampleVector = time::Sample(2, inValuesVector, inGradientsVector);
99 mapping.map(inSampleVector, outValuesVector);
100 expected << 1.2, 2.2, 3.2, 4.2;
101 BOOST_CHECK(equals(expected, outValuesVector));
102}
103
104BOOST_AUTO_TEST_CASE(ConsistentGradientNotConstant)
105{
106
107 PRECICE_TEST(1_rank)
108 int dimensions = 2;
109 using testing::equals;
110
111 // Create mesh to map from
112 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
113 inMesh->createVertex(Eigen::Vector2d::Constant(0.0));
114 inMesh->createVertex(Eigen::Vector2d::Constant(1.0));
115
116 // Create data
117 Eigen::VectorXd inValuesScalar = Eigen::VectorXd::Zero(2);
118 Eigen::VectorXd inValuesVector = Eigen::VectorXd::Zero(4);
119 inValuesScalar << 1.0, 2.0;
120 inValuesVector << 1.0, 2.0, 3.0, 4.0;
121
122 // Create corresponding gradient data (all gradient values = const = 1)
123 Eigen::MatrixXd inGradientsScalar(dimensions, 2);
124 Eigen::MatrixXd inGradientsVector(dimensions, 4);
125
126 inGradientsScalar.col(0) << 2.0, 3.0;
127 inGradientsScalar.col(1) << 2.0, 3.0;
128
129 inGradientsVector.col(0) << 2.0, 3.0;
130 inGradientsVector.col(1) << 4.0, 5.0;
131 inGradientsVector.col(2) << 2.0, 3.0;
132 inGradientsVector.col(3) << 4.0, 5.0;
133
134 // Create mesh to map to
135 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
136 outMesh->createVertex(Eigen::Vector2d::Constant(0.1));
137 outMesh->createVertex(Eigen::Vector2d::Constant(1.1));
138
139 // Setup mapping with mapping coordinates and geometry used
141 mapping.setMeshes(inMesh, outMesh);
142 BOOST_TEST(mapping.hasComputedMapping() == false);
143
144 mapping.computeMapping();
145 Eigen::VectorXd outValuesScalar = Eigen::VectorXd::Zero(2);
146 time::Sample inSampleScalar(1, inValuesScalar, inGradientsScalar);
147 mapping.map(inSampleScalar, outValuesScalar);
148 BOOST_TEST(mapping.hasComputedMapping() == true);
149
150 BOOST_TEST(outValuesScalar(0) == inValuesScalar(0) + 0.5);
151 BOOST_TEST(outValuesScalar(1) == inValuesScalar(1) + 0.5);
152
153 Eigen::VectorXd outValuesVector = Eigen::VectorXd::Zero(4);
154 time::Sample inSampleVector(2, inValuesVector, inGradientsVector);
155 mapping.map(inSampleVector, outValuesVector);
156 Eigen::Vector4d expected(1.5, 2.9, 3.5, 4.9);
157 BOOST_CHECK(equals(expected, outValuesVector));
158}
159
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
BOOST_AUTO_TEST_CASE(ConsistentNonIncremental)
#define PRECICE_TEST(...)
Definition Testing.hpp:27
void setMeshes(const mesh::PtrMesh &input, const mesh::PtrMesh &output)
Sets input and output meshes carrying data to be mapped.
Definition Mapping.cpp:27
bool hasComputedMapping() const
Returns true, if the mapping has been computed.
Definition Mapping.cpp:252
void map(int inputDataID, int outputDataID)
Definition Mapping.cpp:126
void computeMapping() final override
Computes the mapping coefficients from the in- and output mesh.
Mapping using nearest neighboring vertices and their local gradient values.
Container and creator for meshes.
Definition Mesh.hpp:39
Vertex & createVertex(const Eigen::VectorXd &coords)
Creates and initializes a Vertex object.
Definition Mesh.cpp:103
Vertex of a mesh.
Definition Vertex.hpp:16
void setCoords(const VECTOR_T &coordinates)
Sets the coordinates of the vertex.
Definition Vertex.hpp:102
Eigen::VectorXd getCoords() const
Returns the coordinates of the vertex.
Definition Vertex.hpp:116
provides Mesh, Data and primitives.
boost::test_tools::predicate_result equals(const std::vector< float > &VectorA, const std::vector< float > &VectorB, float tolerance)
equals to be used in tests. Compares two std::vectors using a given tolerance. Prints both operands o...
Definition Testing.cpp:65
Main namespace of the precice library.