preCICE v3.1.2
Loading...
Searching...
No Matches
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
32{
33 PRECICE_TEST(1_rank);
34 Eigen::VectorXd vec = Eigen::Vector2d(1, 2);
35 BOOST_TEST(bg::get<0>(vec) == 1);
36 BOOST_TEST(bg::get<1>(vec) == 2);
37 BOOST_TEST(bg::get<2>(vec) == 0);
38}
39
41{
42 PRECICE_TEST(1_rank);
44 auto & v = mesh.createVertex(Eigen::Vector2d(1, 2));
45 BOOST_TEST(bg::get<0>(v) == 1);
46 BOOST_TEST(bg::get<1>(v) == 2);
47 BOOST_TEST(bg::get<2>(v) == 0);
48}
49
50BOOST_AUTO_TEST_CASE(RawCoordinateAdapter)
51{
52 PRECICE_TEST(1_rank);
53 mesh::Vertex::RawCoords v{1, 2, 3};
54 BOOST_TEST(bg::get<0>(v) == 1);
55 BOOST_TEST(bg::get<1>(v) == 2);
56 BOOST_TEST(bg::get<2>(v) == 3);
57 bg::set<1>(v, 5);
58 bg::set<2>(v, 1);
59 BOOST_TEST(bg::get<0>(v) == 1);
60 BOOST_TEST(bg::get<1>(v) == 5);
61 BOOST_TEST(bg::get<2>(v) == 1);
62}
63
65{
66 PRECICE_TEST(1_rank);
68 auto & v1 = mesh.createVertex(Eigen::Vector2d(1, 2));
69 auto & v2 = mesh.createVertex(Eigen::Vector2d(3, 4));
70 auto & e = mesh.createEdge(v1, v2);
71 BOOST_TEST((bg::get<0, 0>(e)) == 1.0);
72 BOOST_TEST((bg::get<0, 1>(e)) == 2.0);
73 BOOST_TEST((bg::get<0, 2>(e)) == 0.0);
74 BOOST_TEST((bg::get<1, 0>(e)) == 3.0);
75 BOOST_TEST((bg::get<1, 1>(e)) == 4.0);
76 BOOST_TEST((bg::get<1, 2>(e)) == 0.0);
77}
78
79BOOST_AUTO_TEST_CASE(TriangleAdapter)
80{
81 PRECICE_TEST(1_rank);
83 auto & v1 = mesh.createVertex(Eigen::Vector3d(0, 2, 0));
84 auto & v2 = mesh.createVertex(Eigen::Vector3d(2, 1, 0));
85 auto & v3 = mesh.createVertex(Eigen::Vector3d(1, 0, 0));
86 auto & e1 = mesh.createEdge(v1, v2);
87 auto & e2 = mesh.createEdge(v2, v3);
88 auto & e3 = mesh.createEdge(v3, v1);
89 auto & t = mesh.createTriangle(e1, e2, e3);
90
91 std::vector<Vertex::RawCoords> vertices(t.begin(), t.end());
92 std::vector<Vertex::RawCoords> refs{v1.rawCoords(), v2.rawCoords(), v3.rawCoords()};
93 BOOST_TEST(vertices.size() == refs.size());
94 BOOST_TEST((std::is_permutation(
95 vertices.begin(), vertices.end(),
96 refs.begin(),
97 [](const auto &lhs, const auto &rhs) {
98 return precice::math::equals(rawToEigen(lhs), rawToEigen(rhs));
99 })));
100}
101
102BOOST_AUTO_TEST_CASE(DistanceTestFlatSingleTriangle)
103{
104 PRECICE_TEST(1_rank);
105 precice::mesh::Mesh mesh("MyMesh", 3, testing::nextMeshID());
106 auto & v1 = mesh.createVertex(Eigen::Vector3d(0, 0, 0));
107 auto & v2 = mesh.createVertex(Eigen::Vector3d(0, 1, 0));
108 auto & v3 = mesh.createVertex(Eigen::Vector3d(1, 0, 0));
109 auto & v4 = mesh.createVertex(Eigen::Vector3d(1, 1, 0));
110 auto & v5 = mesh.createVertex(Eigen::Vector3d(0.2, 0.2, 0));
111 auto & e1 = mesh.createEdge(v1, v2);
112 auto & e2 = mesh.createEdge(v2, v3);
113 auto & e3 = mesh.createEdge(v3, v1);
114 auto & t = mesh.createTriangle(e1, e2, e3);
115
116 BOOST_TEST(bg::comparable_distance(v1, v2) > 0.5);
117 BOOST_TEST(bg::comparable_distance(v1, v1) < 0.01);
118 BOOST_TEST(bg::comparable_distance(e1, v2) < 0.01);
119 BOOST_TEST(bg::comparable_distance(e1, v3) > 0.2);
120 BOOST_TEST(bg::comparable_distance(t, v3) < 0.1);
121 BOOST_TEST(bg::comparable_distance(t, v4) > 0.2);
122 BOOST_TEST(bg::comparable_distance(t, v5) < 0.01);
123}
124
125BOOST_AUTO_TEST_CASE(DistanceTestFlatDoubleTriangle)
126{
127 PRECICE_TEST(1_rank);
128 precice::mesh::Mesh mesh("MyMesh", 3, testing::nextMeshID());
129 auto & lv1 = mesh.createVertex(Eigen::Vector3d(-1, 1, 0.1));
130 auto & lv2 = mesh.createVertex(Eigen::Vector3d(0, -1, 0));
131 auto & lv3 = mesh.createVertex(Eigen::Vector3d(-2, 0, -0.1));
132 auto & le1 = mesh.createEdge(lv1, lv2);
133 auto & le2 = mesh.createEdge(lv2, lv3);
134 auto & le3 = mesh.createEdge(lv3, lv1);
135 auto & lt = mesh.createTriangle(le1, le2, le3);
136
137 auto &rv1 = mesh.createVertex(Eigen::Vector3d(0, 1, 0.1));
138 auto &rv2 = mesh.createVertex(Eigen::Vector3d(2, 0, -0.1));
139 auto &rv3 = mesh.createVertex(Eigen::Vector3d(1, -1, 0));
140 auto &re1 = mesh.createEdge(rv1, rv2);
141 auto &re2 = mesh.createEdge(rv2, rv3);
142 auto &re3 = mesh.createEdge(rv3, rv1);
143 auto &rt = mesh.createTriangle(re1, re2, re3);
144
145 auto &v1 = mesh.createVertex(Eigen::Vector3d(-2, 1, 0));
146 auto &v2 = mesh.createVertex(Eigen::Vector3d(2, -1, 0));
147 auto &v3 = mesh.createVertex(Eigen::Vector3d(0, 0, 0));
148
149 auto lt_v1 = bg::comparable_distance(lt, v1);
150 auto lt_v2 = bg::comparable_distance(lt, v2);
151
152 auto rt_v1 = bg::comparable_distance(rt, v1);
153 auto rt_v3 = bg::comparable_distance(rt, v3);
154
155 BOOST_TEST(precice::testing::equals(lt_v1, 0.5));
156 BOOST_TEST(lt_v2 > 0);
157 BOOST_TEST(rt_v1 > 1);
158 BOOST_TEST(rt_v3 > 0);
159}
160
161BOOST_AUTO_TEST_CASE(DistanceTestFlatDoubleTriangleInsideOutside)
162{
163 PRECICE_TEST(1_rank);
164 precice::mesh::Mesh mesh("MyMesh", 3, testing::nextMeshID());
165 auto & a = mesh.createVertex(Eigen::Vector3d(0, 0, 0));
166 auto & b = mesh.createVertex(Eigen::Vector3d(1, 0, 0));
167 auto & c = mesh.createVertex(Eigen::Vector3d(1, 1, 0));
168 auto & d = mesh.createVertex(Eigen::Vector3d(0, 1, 0));
169
170 auto &ab = mesh.createEdge(a, b);
171 auto &bd = mesh.createEdge(b, d);
172 auto &da = mesh.createEdge(d, a);
173 auto &dc = mesh.createEdge(d, c);
174 auto &cb = mesh.createEdge(c, b);
175
176 auto &lt = mesh.createTriangle(ab, bd, da);
177 auto &rt = mesh.createTriangle(bd, dc, cb);
178 BOOST_TEST_MESSAGE("Left Triangle:" << lt);
179 BOOST_TEST_MESSAGE("Right Triangle:" << rt);
180
181 auto &lv = mesh.createVertex(Eigen::Vector3d(.25, .25, 0));
182 auto &rv = mesh.createVertex(Eigen::Vector3d(.75, .75, 0));
183
184 auto lv_lt = bg::comparable_distance(lv, lt);
185 auto lv_rt = bg::comparable_distance(lv, rt);
186 BOOST_TEST(precice::testing::equals(lv_lt, 0.0), lv_lt << " == 0.0");
187 BOOST_TEST(!precice::testing::equals(lv_rt, 0.0), lv_rt << " != 0.0");
188
189 auto rv_lt = bg::comparable_distance(rv, lt);
190 auto rv_rt = bg::comparable_distance(rv, rt);
191 BOOST_TEST(!precice::testing::equals(rv_lt, 0.0), rv_lt << " != 0.0");
192 BOOST_TEST(precice::testing::equals(rv_rt, 0.0), rv_rt << " == 0.0");
193}
194
195BOOST_AUTO_TEST_CASE(DistanceTestSlopedTriangle)
196{
197 PRECICE_TEST(1_rank);
198 precice::mesh::Mesh mesh("MyMesh", 3, testing::nextMeshID());
199 auto & v1 = mesh.createVertex(Eigen::Vector3d(0, 1, 0));
200 auto & v2 = mesh.createVertex(Eigen::Vector3d(1, 1, 1));
201 auto & v3 = mesh.createVertex(Eigen::Vector3d(0, 0, 1));
202 auto & v4 = mesh.createVertex(Eigen::Vector3d(0, 1, 1));
203 auto & v5 = mesh.createVertex(Eigen::Vector3d(1, 0, 0));
204 auto & e1 = mesh.createEdge(v1, v2);
205 auto & e2 = mesh.createEdge(v2, v3);
206 auto & e3 = mesh.createEdge(v3, v1);
207 auto & t = mesh.createTriangle(e1, e2, e3);
208
209 auto t_v4 = bg::comparable_distance(t, v4);
210 auto v4_t = bg::comparable_distance(v4, t);
211 BOOST_TEST(t_v4 > 0.01);
212 BOOST_TEST(v4_t > 0.01);
213 BOOST_TEST(v4_t == t_v4);
214
215 auto t_v5 = bg::comparable_distance(t, v5);
216 auto v5_t = bg::comparable_distance(v5, t);
217 BOOST_TEST(t_v5 > 0.01);
218 BOOST_TEST(v5_t > 0.01);
219 BOOST_TEST(v5_t == t_v5);
220
221 BOOST_TEST(v4_t < t_v5);
222}
223
224BOOST_AUTO_TEST_CASE(EnvelopeTriangleClockWise)
225{
226 PRECICE_TEST(1_rank);
228 precice::mesh::Mesh mesh("MyMesh", 3, testing::nextMeshID());
229 auto & v1 = mesh.createVertex(Eigen::Vector3d(0, 1, 0));
230 auto & v2 = mesh.createVertex(Eigen::Vector3d(1, 1, 1));
231 auto & v3 = mesh.createVertex(Eigen::Vector3d(0, 0, 1));
232 auto & e1 = mesh.createEdge(v1, v2);
233 auto & e2 = mesh.createEdge(v2, v3);
234 auto & e3 = mesh.createEdge(v3, v1);
235 auto & t = mesh.createTriangle(e1, e2, e3);
236 auto box = bg::return_envelope<precice::query::RTreeBox>(t);
237 auto min = box.min_corner();
238 BOOST_TEST(min[0] == 0.0);
239 BOOST_TEST(min[1] == 0.0);
240 BOOST_TEST(min[2] == 0.0);
241 auto max = box.max_corner();
242 BOOST_TEST(max[0] == 1.0);
243 BOOST_TEST(max[1] == 1.0);
244 BOOST_TEST(max[2] == 1.0);
245}
246
247BOOST_AUTO_TEST_CASE(EnvelopeTriangleCounterclockWise)
248{
249 PRECICE_TEST(1_rank);
251 precice::mesh::Mesh mesh("MyMesh", 3, testing::nextMeshID());
252 auto & v1 = mesh.createVertex(Eigen::Vector3d(0, 1, 0));
253 auto & v2 = mesh.createVertex(Eigen::Vector3d(1, 1, 1));
254 auto & v3 = mesh.createVertex(Eigen::Vector3d(0, 0, 1));
255 auto & e1 = mesh.createEdge(v1, v3);
256 auto & e2 = mesh.createEdge(v3, v2);
257 auto & e3 = mesh.createEdge(v2, v1);
258 auto & t = mesh.createTriangle(e1, e2, e3);
259 auto box = bg::return_envelope<precice::query::RTreeBox>(t);
260 auto min = box.min_corner();
261 BOOST_TEST(min[0] == 0.0);
262 BOOST_TEST(min[1] == 0.0);
263 BOOST_TEST(min[2] == 0.0);
264 auto max = box.max_corner();
265 BOOST_TEST(max[0] == 1.0);
266 BOOST_TEST(max[1] == 1.0);
267 BOOST_TEST(max[2] == 1.0);
268}
269
270BOOST_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:27
T begin(T... args)
Container and creator for meshes.
Definition Mesh.hpp:39
Triangle & createTriangle(Edge &edgeOne, Edge &edgeTwo, Edge &edgeThree)
Creates and initializes a Triangle object.
Definition Mesh.cpp:119
Edge & createEdge(Vertex &vertexOne, Vertex &vertexTwo)
Creates and initializes an Edge object.
Definition Mesh.cpp:111
Vertex & createVertex(const Eigen::VectorXd &coords)
Creates and initializes a Vertex object.
Definition Mesh.cpp:103
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:65
Main namespace of the precice library.
T size(T... args)