preCICE v3.1.2
Loading...
Searching...
No Matches
RadialGeoMultiscaleMappingTest.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <algorithm>
3#include <memory>
5#include "mapping/Mapping.hpp"
7#include "math/constants.hpp"
8#include "mesh/Data.hpp"
9#include "mesh/Mesh.hpp"
11#include "mesh/Utils.hpp"
12#include "mesh/Vertex.hpp"
14#include "testing/Testing.hpp"
15
16using namespace precice;
17using namespace precice::mesh;
18
19BOOST_AUTO_TEST_SUITE(MappingTests)
20BOOST_AUTO_TEST_SUITE(RadialGeoMultiscaleMapping)
21
22BOOST_AUTO_TEST_CASE(ConsistentSpreadX)
23{
24 /* The following test works by creating two dimensionally heterogeneous meshes, namely 1D and 3D, that intersect along the x-axis.
25 Then, the data is mapped from the vertices of the 1D mesh to defined vertices of the 3D mesh (hence, SPREAD).
26 The defined vertices are close to a 1D vertex - in the sense that their projection onto the x-axis is a nearest neighbor of a 1D vertex -
27 and therefore we can predict which value from the 1D vertex should be assigned.
28 Finally, this expected behavior is tested.
29 */
30
31 PRECICE_TEST(1_rank);
32 constexpr int dimensions = 3;
33 using testing::equals;
34
35 // Create mesh to map from
36 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
37 inMesh->createVertex(Eigen::Vector3d::Constant(0.0)); // Point a (1D)
38 inMesh->createVertex(Eigen::Vector3d(3.0, 0.0, 0.0)); // Point b (1D)
39 inMesh->createVertex(Eigen::Vector3d(6.0, 0.0, 0.0)); // Point c (1D)
40 inMesh->allocateDataValues();
41
42 // Create mesh to map to
43 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
44 outMesh->createVertex(Eigen::Vector3d(1.0, 1.0, 0.0)); // Point A (3D): closest to Point a
45 outMesh->createVertex(Eigen::Vector3d(0.5, 0.5, 0.0)); // Point B (3D): closest to Point a
46 outMesh->createVertex(Eigen::Vector3d(4.0, 2.0, 0.0)); // Point C (3D): closest to Point b
47 outMesh->createVertex(Eigen::Vector3d(3.0, 1.0, 0.0)); // Point D (3D): closest to Point b
48 outMesh->createVertex(Eigen::Vector3d(7.0, 3.0, 0.0)); // Point E (3D): closest to Point c
49 outMesh->createVertex(Eigen::Vector3d(6.0, 2.0, 0.0)); // Point F (3D): closest to Point c
50 outMesh->allocateDataValues();
51
52 // Setup mapping with mapping coordinates and geometry used
53 precice::mapping::RadialGeoMultiscaleMapping mapping(mapping::Mapping::CONSISTENT, dimensions, mapping::RadialGeoMultiscaleMapping::MultiscaleType::SPREAD, mapping::RadialGeoMultiscaleMapping::MultiscaleAxis::X);
54 mapping.setMeshes(inMesh, outMesh);
55 BOOST_TEST(mapping.hasComputedMapping() == false);
56
57 // Create data to map
58 Eigen::VectorXd inValues(9);
59 inValues << 1.0, 0.0, 0.0,
60 2.0, 0.0, 0.0,
61 3.0, 0.0, 0.0;
62 const time::Sample inSample{3, inValues};
63 Eigen::VectorXd outValues(18);
64 outValues = Eigen::VectorXd::Zero(18);
65
66 // Map data
67 mapping.computeMapping();
68 mapping.computeMapping(); // Check (only in this case) if calling computeMapping() additional times works.
69 mapping.map(inSample, outValues);
70
71 // Check if data is mapped to closest vertex
72 BOOST_TEST(mapping.hasComputedMapping() == true);
73 mapping.computeMapping(); // Check (only in this case) if calling computeMapping() additional times works.
74
75 // Point A (1D) == Point a (3D)
76 BOOST_TEST(outValues(0) == inSample.values(0));
77 BOOST_TEST(outValues(1) == inSample.values(1));
78 BOOST_TEST(outValues(2) == inSample.values(2));
79
80 // Point B (1D) == Point a (3D)
81 BOOST_TEST(outValues(3) == inSample.values(0));
82 BOOST_TEST(outValues(4) == inSample.values(1));
83 BOOST_TEST(outValues(5) == inSample.values(2));
84
85 // Point C (1D) == Point b (3D)
86 BOOST_TEST(outValues(6) == inSample.values(3));
87 BOOST_TEST(outValues(7) == inSample.values(4));
88 BOOST_TEST(outValues(8) == inSample.values(5));
89
90 // Point D (1D) = Point b (3D)
91 BOOST_TEST(outValues(9) == inSample.values(3));
92 BOOST_TEST(outValues(10) == inSample.values(4));
93 BOOST_TEST(outValues(11) == inSample.values(5));
94
95 // Point E (1D) == Point c (3D)
96 BOOST_TEST(outValues(12) == inSample.values(6));
97 BOOST_TEST(outValues(13) == inSample.values(7));
98 BOOST_TEST(outValues(14) == inSample.values(8));
99
100 // Point F (1D) == Point c (3D)
101 BOOST_TEST(outValues(15) == inSample.values(6));
102 BOOST_TEST(outValues(16) == inSample.values(7));
103 BOOST_TEST(outValues(17) == inSample.values(8));
104}
105
106BOOST_AUTO_TEST_CASE(ConsistentSpreadZ)
107{
108 /* The following test works by creating two dimensionally heterogeneous meshes, namely 1D and 3D, that intersect along the z-axis.
109 Then, the data is mapped from the vertices of the 1D mesh to defined vertices of the 3D mesh (hence, "spread").
110 The defined vertices are close to a 1D vertex - in the sense that their projection onto the z-axis is a nearest neighbor of a 1D vertex -
111 and therefore we can predict which value from the 1D vertex should be assigned.
112 Finally, this expected behavior is tested.
113 */
114
115 PRECICE_TEST(1_rank);
116 constexpr int dimensions = 3;
117 using testing::equals;
118
119 // Create mesh to map from
120 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
121 inMesh->createVertex(Eigen::Vector3d::Constant(0.0)); // Point a (1D)
122 inMesh->createVertex(Eigen::Vector3d(0.0, 0.0, 3.0)); // Point b (1D)
123 inMesh->createVertex(Eigen::Vector3d(0.0, 0.0, 6.0)); // Point c (1D)
124 inMesh->allocateDataValues();
125
126 // Create mesh to map to
127 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
128 outMesh->createVertex(Eigen::Vector3d(0.0, 1.0, 1.0)); // Point A (3D): closest to Point a
129 outMesh->createVertex(Eigen::Vector3d(0.0, 0.5, 0.5)); // Point B (3D): closest to Point a
130 outMesh->createVertex(Eigen::Vector3d(0.0, 2.0, 4.0)); // Point C (3D): closest to Point b
131 outMesh->createVertex(Eigen::Vector3d(0.0, 1.0, 3.0)); // Point D (3D): closest to Point b
132 outMesh->createVertex(Eigen::Vector3d(0.0, 3.0, 7.0)); // Point E (3D): closest to Point c
133 outMesh->createVertex(Eigen::Vector3d(0.0, 2.0, 6.0)); // Point F (3D): closest to Point c
134 outMesh->allocateDataValues();
135
136 // Setup mapping with mapping coordinates and geometry used
137 precice::mapping::RadialGeoMultiscaleMapping mapping(mapping::Mapping::CONSISTENT, dimensions, mapping::RadialGeoMultiscaleMapping::MultiscaleType::SPREAD, mapping::RadialGeoMultiscaleMapping::MultiscaleAxis::Z);
138 mapping.setMeshes(inMesh, outMesh);
139 BOOST_TEST(mapping.hasComputedMapping() == false);
140
141 // Create data to map
142 Eigen::VectorXd inValues(9);
143 inValues << 1.0, 0.0, 0.0,
144 2.0, 0.0, 0.0,
145 3.0, 0.0, 0.0;
146 const time::Sample inSample{3, inValues};
147 Eigen::VectorXd outValues(18);
148 outValues = Eigen::VectorXd::Zero(18);
149
150 // Map data
151 mapping.computeMapping();
152 mapping.map(inSample, outValues);
153
154 // Check if data is mapped to closest vertex
155 BOOST_TEST(mapping.hasComputedMapping() == true);
156
157 // Point A (1D) == Point a (3D)
158 BOOST_TEST(outValues(0) == inSample.values(0));
159 BOOST_TEST(outValues(1) == inSample.values(1));
160 BOOST_TEST(outValues(2) == inSample.values(2));
161
162 // Point B (1D) == Point a (3D)
163 BOOST_TEST(outValues(3) == inSample.values(0));
164 BOOST_TEST(outValues(4) == inSample.values(1));
165 BOOST_TEST(outValues(5) == inSample.values(2));
166
167 // Point C (1D) == Point b (3D)
168 BOOST_TEST(outValues(6) == inSample.values(3));
169 BOOST_TEST(outValues(7) == inSample.values(4));
170 BOOST_TEST(outValues(8) == inSample.values(5));
171
172 // Point D (1D) = Point b (3D)
173 BOOST_TEST(outValues(9) == inSample.values(3));
174 BOOST_TEST(outValues(10) == inSample.values(4));
175 BOOST_TEST(outValues(11) == inSample.values(5));
176
177 // Point E (1D) == Point c (3D)
178 BOOST_TEST(outValues(12) == inSample.values(6));
179 BOOST_TEST(outValues(13) == inSample.values(7));
180 BOOST_TEST(outValues(14) == inSample.values(8));
181
182 // Point F (1D) == Point c (3D)
183 BOOST_TEST(outValues(15) == inSample.values(6));
184 BOOST_TEST(outValues(16) == inSample.values(7));
185 BOOST_TEST(outValues(16) == inSample.values(8));
186}
187
188BOOST_AUTO_TEST_CASE(ConsistentCollectX)
189{
190 /* The following test works by creating two dimensionally heterogeneous meshes, namely 1D and 3D, that intersect along the x-axis.
191 Then, the data is mapped from the vertices of the 3D mesh to defined vertices of the 1D mesh (hence, "collect").
192 The defined vertices are close to a 1D vertex - in the sense that their projection onto the x-axis is a nearest neighbor of a 1D vertex -
193 and therefore we can predict that the 1D vertex value should be the mean of all close 3D vertices.
194 Finally, this expected behavior is tested.
195 */
196
197 PRECICE_TEST(1_rank);
198 constexpr int dimensions = 3;
199 using testing::equals;
200
201 // Create mesh to map from
202 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
203 inMesh->createVertex(Eigen::Vector3d(1.0, 1.0, 0.0)); // Point a (3D)
204 inMesh->createVertex(Eigen::Vector3d(1.0, 2.0, 0.0)); // Point b (3D)
205 inMesh->createVertex(Eigen::Vector3d(2.0, 1.0, 0.0)); // Point c (3D)
206 inMesh->createVertex(Eigen::Vector3d(4.0, 1.0, 0.0)); // Point d (3D)
207 inMesh->createVertex(Eigen::Vector3d(5.0, 2.0, 0.0)); // Point e (3D)
208 inMesh->createVertex(Eigen::Vector3d(7.0, 1.0, 0.0)); // Point f (3D)
209 inMesh->allocateDataValues();
210
211 // Create mesh to map to
212 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
213 outMesh->createVertex(Eigen::Vector3d::Constant(0.0)); // Point A (1D): Averaging from points a,b
214 outMesh->createVertex(Eigen::Vector3d(3.0, 0.0, 0.0)); // Point B (1D): Averaging from points c,d
215 outMesh->createVertex(Eigen::Vector3d(6.0, 0.0, 0.0)); // Point C (1D): Averaging from points e,f
216 outMesh->allocateDataValues();
217
218 // Setup mapping with mapping coordinates and geometry used
219 precice::mapping::RadialGeoMultiscaleMapping mapping(mapping::Mapping::CONSISTENT, dimensions, mapping::RadialGeoMultiscaleMapping::MultiscaleType::COLLECT, mapping::RadialGeoMultiscaleMapping::MultiscaleAxis::X);
220 mapping.setMeshes(inMesh, outMesh);
221 BOOST_TEST(mapping.hasComputedMapping() == false);
222
223 // Create data to map
224 Eigen::VectorXd inValues(18);
225 inValues << 1.0, 0.0, 0.0,
226 3.0, 0.0, 0.0,
227 5.0, 0.0, 0.0,
228 7.0, 0.0, 0.0,
229 2.0, 0.0, 0.0,
230 4.0, 0.0, 0.0;
231 const time::Sample inSample{3, inValues};
232 Eigen::VectorXd outValues(9);
233 outValues = Eigen::VectorXd::Zero(9);
234
235 // Map data
236 mapping.computeMapping();
237 mapping.map(inSample, outValues);
238
239 // Check if data is mapped to closest vertex
240 BOOST_TEST(mapping.hasComputedMapping() == true);
241
242 // Point A (1D) == average of Points a,b (3D)
243 BOOST_TEST(outValues(0) == (inSample.values(0) + inSample.values(3)) / 2);
244 BOOST_TEST(outValues(1) == 0.0);
245 BOOST_TEST(outValues(2) == 0.0);
246
247 // Point B (1D) == average of Points c,d (3D)
248 BOOST_TEST(outValues(3) == (inSample.values(6) + inSample.values(9)) / 2);
249 BOOST_TEST(outValues(4) == 0.0);
250 BOOST_TEST(outValues(5) == 0.0);
251
252 // Point C (1D) == average of Points e,f (3D)
253 BOOST_TEST(outValues(6) == (inSample.values(12) + inSample.values(15)) / 2);
254 BOOST_TEST(outValues(7) == 0.0);
255 BOOST_TEST(outValues(8) == 0.0);
256}
257
258BOOST_AUTO_TEST_CASE(ConsistentCollectZ)
259{
260 /* The following test works by creating two dimensionally heterogeneous meshes, namely 1D and 3D, that intersect along the x-axis.
261 Then, the data is mapped from the vertices of the 3D mesh to defined vertices of the 1D mesh (hence, COLLECT).
262 The defined vertices are close to a 1D vertex and therefore we can predict that the 1D vertex value should be the mean of all close 3D vertices.
263 Finally, this expected behavior is tested.
264 */
265
266 PRECICE_TEST(1_rank);
267 constexpr int dimensions = 3;
268 using testing::equals;
269
270 // Create mesh to map from
271 PtrMesh inMesh(new Mesh("InMesh", dimensions, testing::nextMeshID()));
272 inMesh->createVertex(Eigen::Vector3d(0.0, 1.0, 1.0)); // Point a (3D)
273 inMesh->createVertex(Eigen::Vector3d(0.0, 2.0, 1.0)); // Point b (3D)
274 inMesh->createVertex(Eigen::Vector3d(0.0, 1.0, 2.0)); // Point c (3D)
275 inMesh->createVertex(Eigen::Vector3d(0.0, 1.0, 4.0)); // Point d (3D)
276 inMesh->createVertex(Eigen::Vector3d(0.0, 2.0, 5.0)); // Point e (3D)
277 inMesh->createVertex(Eigen::Vector3d(0.0, 1.0, 7.0)); // Point f (3D)
278 inMesh->allocateDataValues();
279
280 // Create mesh to map to
281 PtrMesh outMesh(new Mesh("OutMesh", dimensions, testing::nextMeshID()));
282 outMesh->createVertex(Eigen::Vector3d::Constant(0.0)); // Point A (1D): Averaging from points a,b
283 outMesh->createVertex(Eigen::Vector3d(0.0, 0.0, 3.0)); // Point B (1D): Averaging from points c,d
284 outMesh->createVertex(Eigen::Vector3d(0.0, 0.0, 6.0)); // Point C (1D): Averaging from points e,f
285 outMesh->allocateDataValues();
286
287 // Setup mapping with mapping coordinates and geometry used
288 precice::mapping::RadialGeoMultiscaleMapping mapping(mapping::Mapping::CONSISTENT, dimensions, mapping::RadialGeoMultiscaleMapping::MultiscaleType::COLLECT, mapping::RadialGeoMultiscaleMapping::MultiscaleAxis::Z);
289 mapping.setMeshes(inMesh, outMesh);
290 BOOST_TEST(mapping.hasComputedMapping() == false);
291
292 // Create data to map
293 Eigen::VectorXd inValues(18);
294 inValues << 1.0, 0.0, 0.0,
295 3.0, 0.0, 0.0,
296 5.0, 0.0, 0.0,
297 7.0, 0.0, 0.0,
298 2.0, 0.0, 0.0,
299 4.0, 0.0, 0.0;
300 const time::Sample inSample{3, inValues};
301 Eigen::VectorXd outValues(9);
302 outValues = Eigen::VectorXd::Zero(9);
303
304 // Map data
305 mapping.computeMapping();
306 mapping.map(inSample, outValues);
307
308 // Check if data is mapped to closest vertex
309 BOOST_TEST(mapping.hasComputedMapping() == true);
310
311 // Point A (1D) == average of Points a,b (3D)
312 BOOST_TEST(outValues(0) == (inSample.values(0) + inSample.values(3)) / 2);
313 BOOST_TEST(outValues(1) == 0.0);
314 BOOST_TEST(outValues(2) == 0.0);
315
316 // Point B (1D) == average of Points c,d (3D)
317 BOOST_TEST(outValues(3) == (inSample.values(6) + inSample.values(9)) / 2);
318 BOOST_TEST(outValues(4) == 0.0);
319 BOOST_TEST(outValues(5) == 0.0);
320
321 // Point C (1D) == average of Points e,f (3D)
322 BOOST_TEST(outValues(6) == (inSample.values(12) + inSample.values(15)) / 2);
323 BOOST_TEST(outValues(7) == 0.0);
324 BOOST_TEST(outValues(8) == 0.0);
325}
326
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
BOOST_AUTO_TEST_CASE(ConsistentSpreadX)
#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
Geometric multiscale mapping in radial direction.
void computeMapping() override
Takes care of compute-heavy operations needed only once to set up the mapping.
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
void allocateDataValues()
Allocates memory for the vertex data values and corresponding gradient values.
Definition Mesh.cpp:233
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.