preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
RTreeAdapterTests.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <algorithm>
3#include <boost/geometry.hpp>
4#include <cmath>
5#include <ostream>
6#include <vector>
7#include "logging/Logger.hpp"
9#include "math/geometry.hpp"
10#include "mesh/Data.hpp"
11#include "mesh/Edge.hpp"
12#include "mesh/Mesh.hpp"
13#include "mesh/Triangle.hpp"
14#include "mesh/Vertex.hpp"
17#include "testing/Testing.hpp"
18
19using namespace precice;
20using namespace precice::mesh;
21using namespace precice::query;
22
23namespace bg = boost::geometry;
24namespace bgi = boost::geometry::index;
25
26BOOST_AUTO_TEST_SUITE(QueryTests)
27BOOST_AUTO_TEST_SUITE(MeshTests)
28BOOST_AUTO_TEST_SUITE(RTreeTests)
29BOOST_AUTO_TEST_SUITE(BGAdapters)
30
33{
35 Eigen::VectorXd vec = Eigen::Vector2d(1, 2);
36 BOOST_TEST(bg::get<0>(vec) == 1);
37 BOOST_TEST(bg::get<1>(vec) == 2);
38 BOOST_TEST(bg::get<2>(vec) == 0);
39}
40
43{
46 auto &v = mesh.createVertex(Eigen::Vector2d(1, 2));
47 BOOST_TEST(bg::get<0>(v) == 1);
48 BOOST_TEST(bg::get<1>(v) == 2);
49 BOOST_TEST(bg::get<2>(v) == 0);
50}
51
53BOOST_AUTO_TEST_CASE(RawCoordinateAdapter)
54{
56 mesh::Vertex::RawCoords v{1, 2, 3};
57 BOOST_TEST(bg::get<0>(v) == 1);
58 BOOST_TEST(bg::get<1>(v) == 2);
59 BOOST_TEST(bg::get<2>(v) == 3);
60 bg::set<1>(v, 5);
61 bg::set<2>(v, 1);
62 BOOST_TEST(bg::get<0>(v) == 1);
63 BOOST_TEST(bg::get<1>(v) == 5);
64 BOOST_TEST(bg::get<2>(v) == 1);
65}
66
69{
72 auto &v1 = mesh.createVertex(Eigen::Vector2d(1, 2));
73 auto &v2 = mesh.createVertex(Eigen::Vector2d(3, 4));
74 auto &e = mesh.createEdge(v1, v2);
75 BOOST_TEST((bg::get<0, 0>(e)) == 1.0);
76 BOOST_TEST((bg::get<0, 1>(e)) == 2.0);
77 BOOST_TEST((bg::get<0, 2>(e)) == 0.0);
78 BOOST_TEST((bg::get<1, 0>(e)) == 3.0);
79 BOOST_TEST((bg::get<1, 1>(e)) == 4.0);
80 BOOST_TEST((bg::get<1, 2>(e)) == 0.0);
81}
82
84BOOST_AUTO_TEST_CASE(TriangleAdapter)
85{
88 auto &v1 = mesh.createVertex(Eigen::Vector3d(0, 2, 0));
89 auto &v2 = mesh.createVertex(Eigen::Vector3d(2, 1, 0));
90 auto &v3 = mesh.createVertex(Eigen::Vector3d(1, 0, 0));
91 auto &e1 = mesh.createEdge(v1, v2);
92 auto &e2 = mesh.createEdge(v2, v3);
93 auto &e3 = mesh.createEdge(v3, v1);
94 auto &t = mesh.createTriangle(e1, e2, e3);
95
96 std::vector<Vertex::RawCoords> vertices(t.begin(), t.end());
97 std::vector<Vertex::RawCoords> refs{v1.rawCoords(), v2.rawCoords(), v3.rawCoords()};
98 BOOST_TEST(vertices.size() == refs.size());
99 BOOST_TEST((std::is_permutation(
100 vertices.begin(), vertices.end(),
101 refs.begin(),
102 [](const auto &lhs, const auto &rhs) {
103 return precice::math::equals(rawToEigen(lhs), rawToEigen(rhs));
104 })));
105}
106
107PRECICE_TEST_SETUP(1_rank)
108BOOST_AUTO_TEST_CASE(DistanceTestFlatSingleTriangle)
109{
110 PRECICE_TEST();
111 precice::mesh::Mesh mesh("MyMesh", 3, testing::nextMeshID());
112 auto &v1 = mesh.createVertex(Eigen::Vector3d(0, 0, 0));
113 auto &v2 = mesh.createVertex(Eigen::Vector3d(0, 1, 0));
114 auto &v3 = mesh.createVertex(Eigen::Vector3d(1, 0, 0));
115 auto &v4 = mesh.createVertex(Eigen::Vector3d(1, 1, 0));
116 auto &v5 = mesh.createVertex(Eigen::Vector3d(0.2, 0.2, 0));
117 auto &e1 = mesh.createEdge(v1, v2);
118 auto &e2 = mesh.createEdge(v2, v3);
119 auto &e3 = mesh.createEdge(v3, v1);
120 auto &t = mesh.createTriangle(e1, e2, e3);
121
122 BOOST_TEST(bg::comparable_distance(v1, v2) > 0.5);
123 BOOST_TEST(bg::comparable_distance(v1, v1) < 0.01);
124 BOOST_TEST(bg::comparable_distance(e1, v2) < 0.01);
125 BOOST_TEST(bg::comparable_distance(e1, v3) > 0.2);
126 BOOST_TEST(bg::comparable_distance(t, v3) < 0.1);
127 BOOST_TEST(bg::comparable_distance(t, v4) > 0.2);
128 BOOST_TEST(bg::comparable_distance(t, v5) < 0.01);
129}
130
131PRECICE_TEST_SETUP(1_rank)
132BOOST_AUTO_TEST_CASE(DistanceTestFlatDoubleTriangle)
133{
134 PRECICE_TEST();
135 precice::mesh::Mesh mesh("MyMesh", 3, testing::nextMeshID());
136 auto &lv1 = mesh.createVertex(Eigen::Vector3d(-1, 1, 0.1));
137 auto &lv2 = mesh.createVertex(Eigen::Vector3d(0, -1, 0));
138 auto &lv3 = mesh.createVertex(Eigen::Vector3d(-2, 0, -0.1));
139 auto &le1 = mesh.createEdge(lv1, lv2);
140 auto &le2 = mesh.createEdge(lv2, lv3);
141 auto &le3 = mesh.createEdge(lv3, lv1);
142 auto &lt = mesh.createTriangle(le1, le2, le3);
143
144 auto &rv1 = mesh.createVertex(Eigen::Vector3d(0, 1, 0.1));
145 auto &rv2 = mesh.createVertex(Eigen::Vector3d(2, 0, -0.1));
146 auto &rv3 = mesh.createVertex(Eigen::Vector3d(1, -1, 0));
147 auto &re1 = mesh.createEdge(rv1, rv2);
148 auto &re2 = mesh.createEdge(rv2, rv3);
149 auto &re3 = mesh.createEdge(rv3, rv1);
150 auto &rt = mesh.createTriangle(re1, re2, re3);
151
152 auto &v1 = mesh.createVertex(Eigen::Vector3d(-2, 1, 0));
153 auto &v2 = mesh.createVertex(Eigen::Vector3d(2, -1, 0));
154 auto &v3 = mesh.createVertex(Eigen::Vector3d(0, 0, 0));
155
156 auto lt_v1 = bg::comparable_distance(lt, v1);
157 auto lt_v2 = bg::comparable_distance(lt, v2);
158
159 auto rt_v1 = bg::comparable_distance(rt, v1);
160 auto rt_v3 = bg::comparable_distance(rt, v3);
161
162 BOOST_TEST(precice::testing::equals(lt_v1, 0.5));
163 BOOST_TEST(lt_v2 > 0);
164 BOOST_TEST(rt_v1 > 1);
165 BOOST_TEST(rt_v3 > 0);
166}
167
168PRECICE_TEST_SETUP(1_rank)
169BOOST_AUTO_TEST_CASE(DistanceTestFlatDoubleTriangleInsideOutside)
170{
171 PRECICE_TEST();
172 precice::mesh::Mesh mesh("MyMesh", 3, testing::nextMeshID());
173 auto &a = mesh.createVertex(Eigen::Vector3d(0, 0, 0));
174 auto &b = mesh.createVertex(Eigen::Vector3d(1, 0, 0));
175 auto &c = mesh.createVertex(Eigen::Vector3d(1, 1, 0));
176 auto &d = mesh.createVertex(Eigen::Vector3d(0, 1, 0));
177
178 auto &ab = mesh.createEdge(a, b);
179 auto &bd = mesh.createEdge(b, d);
180 auto &da = mesh.createEdge(d, a);
181 auto &dc = mesh.createEdge(d, c);
182 auto &cb = mesh.createEdge(c, b);
183
184 auto &lt = mesh.createTriangle(ab, bd, da);
185 auto &rt = mesh.createTriangle(bd, dc, cb);
186 BOOST_TEST_MESSAGE("Left Triangle:" << lt);
187 BOOST_TEST_MESSAGE("Right Triangle:" << rt);
188
189 auto &lv = mesh.createVertex(Eigen::Vector3d(.25, .25, 0));
190 auto &rv = mesh.createVertex(Eigen::Vector3d(.75, .75, 0));
191
192 auto lv_lt = bg::comparable_distance(lv, lt);
193 auto lv_rt = bg::comparable_distance(lv, rt);
194 BOOST_TEST(precice::testing::equals(lv_lt, 0.0), lv_lt << " == 0.0");
195 BOOST_TEST(!precice::testing::equals(lv_rt, 0.0), lv_rt << " != 0.0");
196
197 auto rv_lt = bg::comparable_distance(rv, lt);
198 auto rv_rt = bg::comparable_distance(rv, rt);
199 BOOST_TEST(!precice::testing::equals(rv_lt, 0.0), rv_lt << " != 0.0");
200 BOOST_TEST(precice::testing::equals(rv_rt, 0.0), rv_rt << " == 0.0");
201}
202
203PRECICE_TEST_SETUP(1_rank)
204BOOST_AUTO_TEST_CASE(DistanceTestSlopedTriangle)
205{
206 PRECICE_TEST();
207 precice::mesh::Mesh mesh("MyMesh", 3, testing::nextMeshID());
208 auto &v1 = mesh.createVertex(Eigen::Vector3d(0, 1, 0));
209 auto &v2 = mesh.createVertex(Eigen::Vector3d(1, 1, 1));
210 auto &v3 = mesh.createVertex(Eigen::Vector3d(0, 0, 1));
211 auto &v4 = mesh.createVertex(Eigen::Vector3d(0, 1, 1));
212 auto &v5 = mesh.createVertex(Eigen::Vector3d(1, 0, 0));
213 auto &e1 = mesh.createEdge(v1, v2);
214 auto &e2 = mesh.createEdge(v2, v3);
215 auto &e3 = mesh.createEdge(v3, v1);
216 auto &t = mesh.createTriangle(e1, e2, e3);
217
218 auto t_v4 = bg::comparable_distance(t, v4);
219 auto v4_t = bg::comparable_distance(v4, t);
220 BOOST_TEST(t_v4 > 0.01);
221 BOOST_TEST(v4_t > 0.01);
222 BOOST_TEST(v4_t == t_v4);
223
224 auto t_v5 = bg::comparable_distance(t, v5);
225 auto v5_t = bg::comparable_distance(v5, t);
226 BOOST_TEST(t_v5 > 0.01);
227 BOOST_TEST(v5_t > 0.01);
228 BOOST_TEST(v5_t == t_v5);
229
230 BOOST_TEST(v4_t < t_v5);
231}
232
233PRECICE_TEST_SETUP(1_rank)
234BOOST_AUTO_TEST_CASE(EnvelopeTriangleClockWise)
235{
236 PRECICE_TEST();
238 precice::mesh::Mesh mesh("MyMesh", 3, testing::nextMeshID());
239 auto &v1 = mesh.createVertex(Eigen::Vector3d(0, 1, 0));
240 auto &v2 = mesh.createVertex(Eigen::Vector3d(1, 1, 1));
241 auto &v3 = mesh.createVertex(Eigen::Vector3d(0, 0, 1));
242 auto &e1 = mesh.createEdge(v1, v2);
243 auto &e2 = mesh.createEdge(v2, v3);
244 auto &e3 = mesh.createEdge(v3, v1);
245 auto &t = mesh.createTriangle(e1, e2, e3);
246 auto box = bg::return_envelope<precice::query::RTreeBox>(t);
247 auto min = box.min_corner();
248 BOOST_TEST(min[0] == 0.0);
249 BOOST_TEST(min[1] == 0.0);
250 BOOST_TEST(min[2] == 0.0);
251 auto max = box.max_corner();
252 BOOST_TEST(max[0] == 1.0);
253 BOOST_TEST(max[1] == 1.0);
254 BOOST_TEST(max[2] == 1.0);
255}
256
257PRECICE_TEST_SETUP(1_rank)
258BOOST_AUTO_TEST_CASE(EnvelopeTriangleCounterclockWise)
259{
260 PRECICE_TEST();
262 precice::mesh::Mesh mesh("MyMesh", 3, testing::nextMeshID());
263 auto &v1 = mesh.createVertex(Eigen::Vector3d(0, 1, 0));
264 auto &v2 = mesh.createVertex(Eigen::Vector3d(1, 1, 1));
265 auto &v3 = mesh.createVertex(Eigen::Vector3d(0, 0, 1));
266 auto &e1 = mesh.createEdge(v1, v3);
267 auto &e2 = mesh.createEdge(v3, v2);
268 auto &e3 = mesh.createEdge(v2, v1);
269 auto &t = mesh.createTriangle(e1, e2, e3);
270 auto box = bg::return_envelope<precice::query::RTreeBox>(t);
271 auto min = box.min_corner();
272 BOOST_TEST(min[0] == 0.0);
273 BOOST_TEST(min[1] == 0.0);
274 BOOST_TEST(min[2] == 0.0);
275 auto max = box.max_corner();
276 BOOST_TEST(max[0] == 1.0);
277 BOOST_TEST(max[1] == 1.0);
278 BOOST_TEST(max[2] == 1.0);
279}
280
281BOOST_AUTO_TEST_SUITE_END() // BG Adapters
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
BOOST_AUTO_TEST_CASE(VectorAdapter)
#define PRECICE_TEST()
Definition Testing.hpp:39
#define PRECICE_TEST_SETUP(...)
Creates and attaches a TestSetup to a Boost test case.
Definition Testing.hpp:29
T begin(T... args)
Container and creator for meshes.
Definition Mesh.hpp:38
Triangle & createTriangle(Edge &edgeOne, Edge &edgeTwo, Edge &edgeThree)
Creates and initializes a Triangle object.
Definition Mesh.cpp:121
Edge & createEdge(Vertex &vertexOne, Vertex &vertexTwo)
Creates and initializes an Edge object.
Definition Mesh.cpp:113
Vertex & createVertex(const Eigen::Ref< const Eigen::VectorXd > &coords)
Creates and initializes a Vertex object.
Definition Mesh.cpp:105
T end(T... args)
T is_permutation(T... args)
provides Mesh, Data and primitives.
contains geometrical queries.
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:93
Main namespace of the precice library.
T size(T... args)