preCICE v3.1.2
Loading...
Searching...
No Matches
MeshTest.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <algorithm>
3#include <deque>
4#include <iosfwd>
5#include <memory>
6#include <string>
7#include <vector>
8
9#include "logging/Logger.hpp"
10#include "mesh/BoundingBox.hpp"
11#include "mesh/Data.hpp"
12#include "mesh/Edge.hpp"
13#include "mesh/Mesh.hpp"
15#include "mesh/Triangle.hpp"
16#include "mesh/Utils.hpp"
17#include "mesh/Vertex.hpp"
19#include "testing/Testing.hpp"
20#include "utils/algorithm.hpp"
21
22using namespace precice;
23using namespace precice::mesh;
24using Eigen::Vector2d;
25using Eigen::Vector3d;
27
28BOOST_AUTO_TEST_SUITE(MeshTests)
29
30BOOST_AUTO_TEST_SUITE(MeshTests)
31
32BOOST_AUTO_TEST_CASE(BoundingBoxCOG_2D)
33{
34 PRECICE_TEST(1_rank);
35 Eigen::Vector2d coords0(2, 0);
36 Eigen::Vector2d coords1(-1, 4);
37 Eigen::Vector2d coords2(0, 1);
38
39 mesh::Mesh mesh("2D Testmesh", 2, testing::nextMeshID());
40 mesh.createVertex(coords0);
41 mesh.createVertex(coords1);
42 mesh.createVertex(coords2);
43
44 mesh.computeBoundingBox();
45
47 auto cog = bBox.center();
48
49 mesh::BoundingBox referenceBox({-1.0, 2.0,
50 0.0, 4.0});
51
52 std::vector<double> referenceCOG = {0.5, 2.0};
53
54 BOOST_TEST(bBox.getDimension() == 2);
55 BOOST_TEST(cog.size() == 2);
56 BOOST_TEST(referenceBox == bBox);
57
58 for (decltype(cog.size()) d = 0; d < cog.size(); d++) {
59 BOOST_TEST(referenceCOG.at(d) == cog(d));
60 }
61}
62
63BOOST_AUTO_TEST_CASE(BoundingBoxCOG_3D)
64{
65 PRECICE_TEST(1_rank);
66 Eigen::Vector3d coords0(2, 0, -3);
67 Eigen::Vector3d coords1(-1, 4, 8);
68 Eigen::Vector3d coords2(0, 1, -2);
69 Eigen::Vector3d coords3(3.5, 2, -2);
70
71 mesh::Mesh mesh("3D Testmesh", 3, testing::nextMeshID());
72 mesh.createVertex(coords0);
73 mesh.createVertex(coords1);
74 mesh.createVertex(coords2);
75 mesh.createVertex(coords3);
76
77 mesh.computeBoundingBox();
78
80 auto cog = bBox.center();
81
82 mesh::BoundingBox referenceBox({-1.0, 3.5,
83 0.0, 4.0,
84 -3.0, 8.0});
85
86 std::vector<double> referenceCOG = {1.25, 2.0, 2.5};
87
88 BOOST_TEST(bBox.getDimension() == 3);
89 BOOST_TEST(cog.size() == 3);
90 BOOST_TEST(referenceBox == bBox);
91
92 for (decltype(cog.size()) d = 0; d < cog.size(); d++) {
93 BOOST_TEST(referenceCOG.at(d) == cog(d));
94 }
95}
96
98{
99 PRECICE_TEST(1_rank);
100 for (int dim = 2; dim <= 3; dim++) {
101 // Create mesh object
102 std::string meshName("MyMesh");
103 precice::mesh::Mesh mesh(meshName, dim, testing::nextMeshID());
104
105 // Validate mesh object state
106 BOOST_TEST(mesh.getName() == meshName);
107
108 // Create mesh vertices
109 Eigen::VectorXd coords0(dim);
110 Eigen::VectorXd coords1(dim);
111 Eigen::VectorXd coords2(dim);
112 if (dim == 2) {
113 coords0 << 0.0, 0.0;
114 coords1 << 1.0, 0.0;
115 coords2 << 0.0, 1.0;
116 } else {
117 coords0 << 0.0, 0.0, 0.0;
118 coords1 << 1.0, 0.0, 0.0;
119 coords2 << 0.0, 0.0, 1.0;
120 }
121 Vertex &v0 = mesh.createVertex(coords0);
122 Vertex &v1 = mesh.createVertex(coords1);
123 Vertex &v2 = mesh.createVertex(coords2);
124
125 // Validate mesh vertices state
126 // This is the preferred way to iterate over elements in a mesh, it hides
127 // the details of the vertex container in class Mesh.
128 size_t index = 0;
129 for (Vertex &vertex : mesh.vertices()) {
130 if (index == 0) {
131 BOOST_TEST(vertex.getID() == v0.getID());
132 } else if (index == 1) {
133 BOOST_TEST(vertex.getID() == v1.getID());
134 } else if (index == 2) {
135 BOOST_TEST(vertex.getID() == v2.getID());
136 } else {
137 BOOST_TEST(false);
138 }
139 index++;
140 }
141
142 BOOST_TEST(!mesh.hasEdges());
143 BOOST_TEST(!mesh.hasTriangles());
144
145 // Create mesh edges
146 Edge &e0 = mesh.createEdge(v0, v1);
147 Edge &e1 = mesh.createEdge(v1, v2);
148 Edge &e2 = mesh.createEdge(v2, v0);
149
150 BOOST_TEST(mesh.hasEdges());
151 BOOST_TEST(mesh.edges().size() == 3);
152 BOOST_TEST(!mesh.hasTriangles());
153
154 Triangle *t = nullptr;
155 if (dim == 3) {
156 // Create triangle
157 t = &mesh.createTriangle(e0, e1, e2);
158
159 BOOST_TEST(mesh.hasTriangles());
160 } else {
161 BOOST_TEST(!mesh.hasTriangles());
162 }
163
164 BOOST_TEST(mesh.hasEdges());
165
166 // Create vertex data
167 std::string dataName("MyData");
168 int dataDimensions = dim;
169 // Add a data set to the mesh. Every data value is associated to a vertex in
170 // the mesh via the vertex ID.
171 PtrData data = mesh.createData(dataName, dataDimensions, 0_dataID);
172
173 // Validate data state
174 BOOST_TEST(data->getName() == dataName);
175 BOOST_TEST(data->getDimensions() == dataDimensions);
176
177 // Validate state of mesh with data
178 BOOST_TEST(mesh.data().size() == 1);
179 BOOST_TEST(mesh.data(0)->getName() == dataName);
180
181 // Allocate memory for the data values of set data. Before data value access
182 // leads to assertions.
183 mesh.allocateDataValues();
184
185 // Access data values
186 Eigen::VectorXd &dataValues = data->values();
187 BOOST_TEST(dataValues.size() == 3 * dim);
188 BOOST_TEST(v0.getID() == 0);
189 Eigen::VectorXd value = Eigen::VectorXd::Zero(dim);
190 for (int i = 0; i < dim; i++) {
191 value(i) = dataValues(v0.getID() * dim + i);
192 }
193 }
194}
195
197{
198 PRECICE_TEST(1_rank);
199 int dim = 3;
200 Mesh mesh1("Mesh1", dim, testing::nextMeshID());
201 Mesh mesh2("Mesh2", dim, testing::nextMeshID());
202 std::array<Mesh *, 2> meshes = {&mesh1, &mesh2};
203 for (auto ptr : meshes) {
204 auto & mesh = *ptr;
205 Eigen::VectorXd coords0(dim);
206 Eigen::VectorXd coords1(dim);
207 Eigen::VectorXd coords2(dim);
208 Eigen::VectorXd coords3(dim);
209 coords0 << 0.0, 0.0, 0.0;
210 coords1 << 1.0, 0.0, 0.0;
211 coords2 << 0.0, 0.0, 1.0;
212 coords3 << 1.0, 0.0, 1.0;
213 Vertex &v0 = mesh.createVertex(coords0);
214 Vertex &v1 = mesh.createVertex(coords1);
215 Vertex &v2 = mesh.createVertex(coords2);
216 Vertex &v3 = mesh.createVertex(coords3);
217 Edge & e0 = mesh.createEdge(v0, v1); // LINESTRING (0 0 0, 1 0 0)
218 Edge & e1 = mesh.createEdge(v1, v2); // LINESTRING (1 0 0, 0 0 1)
219 Edge & e2 = mesh.createEdge(v2, v0); // LINESTRING (0 0 1, 0 0 0)
220 mesh.createEdge(v1, v3); // LINESTRING (1 0 0, 1 0 1)
221 mesh.createEdge(v3, v2); // LINESTRING (1 0 1, 0 0 1)
222 mesh.createTriangle(e0, e1, e2);
223 }
224 BOOST_TEST(mesh1 == mesh2);
225}
226
228{
229 PRECICE_TEST(1_rank);
230 Mesh mesh("WKTMesh", 3, testing::nextMeshID());
231 Vertex &v0 = mesh.createVertex(Eigen::Vector3d(0., 0., 0.));
232 Vertex &v1 = mesh.createVertex(Eigen::Vector3d(1., 0., 0.));
233 Vertex &v2 = mesh.createVertex(Eigen::Vector3d(0., 0., 1.));
234 Vertex &v3 = mesh.createVertex(Eigen::Vector3d(1., 0., 1.));
235 Edge & e0 = mesh.createEdge(v0, v1); // LINESTRING (0 0 0, 1 0 0)
236 Edge & e1 = mesh.createEdge(v1, v2); // LINESTRING (1 0 0, 0 0 1)
237 Edge & e2 = mesh.createEdge(v2, v0); // LINESTRING (0 0 1, 0 0 0)
238 mesh.createEdge(v1, v3); // LINESTRING (1 0 0, 1 0 1)
239 mesh.createEdge(v3, v2); // LINESTRING (1 0 1, 0 0 1)
240 mesh.createTriangle(e0, e1, e2);
241 std::stringstream sstream;
242 sstream << mesh;
243 std::string reference(
244 "Mesh \"WKTMesh\", dimensionality = 3:\n"
245 "GEOMETRYCOLLECTION(\n"
246 "POINT (0 0 0), POINT (1 0 0), POINT (0 0 1), POINT (1 0 1),\n"
247 "LINESTRING (0 0 0, 1 0 0), LINESTRING (1 0 0, 0 0 1), LINESTRING (0 0 0, 0 0 1), LINESTRING (1 0 0, 1 0 1), LINESTRING (0 0 1, 1 0 1),\n"
248 "POLYGON ((0 0 0, 1 0 0, 0 0 1, 0 0 0))\n"
249 ")");
250 BOOST_TEST(reference == sstream.str());
251}
252
253BOOST_AUTO_TEST_CASE(ResizeDataGrow)
254{
255 PRECICE_TEST(1_rank);
256 precice::mesh::Mesh mesh("MyMesh", 3, testing::nextMeshID());
257 const auto & values = mesh.createData("Data", 1, 0_dataID)->values();
258
259 // Create mesh
260 mesh.createVertex(Vector3d(0.0, 0.0, 0.0));
261 mesh.createVertex(Vector3d(1.0, 0.0, 1.0));
262
263 BOOST_TEST(mesh.nVertices() == 2);
264 mesh.allocateDataValues();
265 BOOST_TEST(values.size() == 2);
266
267 mesh.createVertex(Vector3d(1.0, 1.0, 1.0));
268 mesh.createVertex(Vector3d(2.0, 0.0, 2.0));
269 mesh.createVertex(Vector3d(2.0, 0.0, 2.1));
270
271 BOOST_TEST(mesh.nVertices() == 5);
272 mesh.allocateDataValues();
273 BOOST_TEST(values.size() == 5);
274}
275
276BOOST_AUTO_TEST_CASE(ResizeDataShrink)
277{
278 PRECICE_TEST(1_rank);
279 precice::mesh::Mesh mesh("MyMesh", 3, testing::nextMeshID());
280 const auto & values = mesh.createData("Data", 1, 0_dataID)->values();
281
282 // Create mesh
283 mesh.createVertex(Vector3d(0.0, 0.0, 0.0));
284 mesh.createVertex(Vector3d(1.0, 0.0, 1.0));
285 mesh.createVertex(Vector3d(1.0, 1.0, 1.0));
286 mesh.createVertex(Vector3d(2.0, 2.0, 2.0));
287
288 BOOST_TEST(mesh.nVertices() == 4);
289 mesh.allocateDataValues();
290 BOOST_TEST(values.size() == 4);
291
292 mesh.clear();
293 mesh.createVertex(Vector3d(0.0, 0.0, 0.0));
294 mesh.createVertex(Vector3d(1.0, 0.0, 1.0));
295
296 BOOST_TEST(mesh.nVertices() == 2);
297 mesh.allocateDataValues();
298 BOOST_TEST(values.size() == 2);
299}
300
302
304{
305 PRECICE_TEST(1_rank);
306 Mesh mesh("Mesh1", 3, testing::nextMeshID());
307 mesh.createData("Data", 1, 0_dataID);
308
309 Eigen::Vector3d coords0;
310 Eigen::Vector3d coords1;
311 Eigen::Vector3d coords2;
312 Eigen::Vector3d coords3;
313 coords0 << 1.0, 0.0, 0.0;
314 coords1 << 0.0, 0.0, 0.0;
315 coords2 << 2.0, 2.0, 0.0;
316 coords3 << 0.0, 1.0, 0.0;
317 Vertex &v0 = mesh.createVertex(coords0);
318 Vertex &v1 = mesh.createVertex(coords1);
319 Vertex &v2 = mesh.createVertex(coords2);
320 Vertex &v3 = mesh.createVertex(coords3);
321
322 BOOST_TEST(mesh.nVertices() == 4);
323
324 Edge &e0 = mesh.createEdge(v0, v1);
325 Edge &e1 = mesh.createEdge(v3, v2);
326 Edge &e2 = mesh.createEdge(v3, v1);
327 Edge &e3 = mesh.createEdge(v0, v2);
328
329 BOOST_TEST(mesh.edges().size() == 4);
330
331 auto chain = asChain(utils::make_array(&e0, &e1, &e2, &e3));
332
333 BOOST_REQUIRE(chain.connected);
334
335 BOOST_TEST(chain.edges.at(0) == &e0);
336 BOOST_TEST(chain.edges.at(1) == &e2);
337 BOOST_TEST(chain.edges.at(2) == &e1);
338 BOOST_TEST(chain.edges.at(3) == &e3);
339
340 BOOST_TEST(chain.vertices.at(0) == &v1);
341 BOOST_TEST(chain.vertices.at(1) == &v3);
342 BOOST_TEST(chain.vertices.at(2) == &v2);
343 BOOST_TEST(chain.vertices.at(3) == &v0);
344}
345
347{
348 PRECICE_TEST(1_rank);
349 Mesh mesh("Mesh1", 3, testing::nextMeshID());
350 mesh.createData("Data", 1, 0_dataID);
351
352 Eigen::Vector3d coords0;
353 Eigen::Vector3d coords1;
354 Eigen::Vector3d coords2;
355 Eigen::Vector3d coords3;
356 coords0 << 1.0, 0.0, 0.0;
357 coords1 << 0.0, 0.0, 0.0;
358 coords2 << 2.0, 2.0, 0.0;
359 coords3 << 0.0, 1.0, 0.0;
360 Vertex *v0 = &mesh.createVertex(coords0);
361 Vertex *v1 = &mesh.createVertex(coords1);
362 Vertex *v2 = &mesh.createVertex(coords2);
363 Vertex *v3 = &mesh.createVertex(coords3);
364 BOOST_REQUIRE(mesh.nVertices() == 4);
365
366 Edge &e0 = mesh.createEdge(*v0, *v1);
367 Edge &e1 = mesh.createEdge(*v1, *v2);
368 Edge &e2 = mesh.createEdge(*v2, *v3);
369 Edge &e3 = mesh.createEdge(*v3, *v0);
370 BOOST_REQUIRE(mesh.edges().size() == 4);
371
372 BOOST_TEST(sharedVertex(e0, e1) == v1);
373 BOOST_TEST(sharedVertex(e1, e2) == v2);
374 BOOST_TEST(sharedVertex(e2, e3) == v3);
375 BOOST_TEST(sharedVertex(e3, e0) == v0);
376
377 // not connected
378 BOOST_TEST((sharedVertex(e0, e2) == nullptr));
379 BOOST_TEST((sharedVertex(e1, e3) == nullptr));
380
381 // check symmetry
382 BOOST_TEST(sharedVertex(e0, e1) == sharedVertex(e1, e0));
383 BOOST_TEST(sharedVertex(e0, e2) == sharedVertex(e2, e0));
384}
385
387{
388 PRECICE_TEST(1_rank);
389
390 Eigen::Vector3d coords0;
391 Eigen::Vector3d coords1;
392 coords0 << 1.0, 0.0, 0.0;
393 coords1 << 0.0, 1.0, 0.0;
394 Vertex v0{coords0, 0};
395 Vertex v1{coords1, 1};
396 Edge e(v0, v1);
397 BOOST_TEST(edgeLength(e) == std::sqrt(2));
398}
399
401{
402 PRECICE_TEST(1_rank);
403 Mesh mesh("Mesh1", 3, testing::nextMeshID());
404 mesh.createData("Data", 1, 0_dataID);
405
406 Eigen::Vector3d coords0;
407 Eigen::Vector3d coords1;
408 Eigen::Vector3d coords2;
409 Eigen::Vector3d coords3;
410 coords0 << 1.0, 0.0, 0.0;
411 coords1 << 0.0, 0.0, 0.0;
412 coords2 << 2.0, 2.0, 0.0;
413 coords3 << 0.0, 1.0, 0.0;
414 Vertex &v0 = mesh.createVertex(coords0);
415 Vertex &v1 = mesh.createVertex(coords1);
416 Vertex &v2 = mesh.createVertex(coords2);
417 mesh.createVertex(coords3);
418 BOOST_TEST(mesh.nVertices() == 4);
419
420 std::array<int, 3> ids{v0.getID(), v2.getID(), v1.getID()};
421 std::array<Vertex *, 3> expected{&v0, &v2, &v1};
422
423 auto result = vertexPtrsFor(mesh, ids);
424 BOOST_REQUIRE(result.size() == 3);
425 BOOST_TEST(result == expected);
426}
427
429{
430 PRECICE_TEST(1_rank);
431 Mesh mesh("Mesh1", 3, testing::nextMeshID());
432 mesh.createData("Data", 1, 0_dataID);
433
434 Eigen::Vector3d coords0;
435 Eigen::Vector3d coords1;
436 Eigen::Vector3d coords2;
437 Eigen::Vector3d coords3;
438 coords0 << 1.0, 0.0, 0.0;
439 coords1 << 0.0, 0.0, 0.0;
440 coords2 << 2.0, 2.0, 0.0;
441 coords3 << 0.0, 1.0, 0.0;
442 Vertex &v0 = mesh.createVertex(coords0);
443 Vertex &v1 = mesh.createVertex(coords1);
444 Vertex &v2 = mesh.createVertex(coords2);
445 mesh.createVertex(coords3);
446 BOOST_TEST(mesh.nVertices() == 4);
447
448 std::array<int, 3> ids{v0.getID(), v2.getID(), v1.getID()};
449 std::array<Eigen::VectorXd, 3> expected{coords0, coords2, coords1};
450
451 auto result = coordsFor(mesh, ids);
452 BOOST_REQUIRE(result.size() == 3);
453 BOOST_TEST(result == expected);
454}
455
457{
458 PRECICE_TEST(1_rank);
459 Mesh mesh("Mesh1", 3, testing::nextMeshID());
460 mesh.createData("Data", 1, 0_dataID);
461
462 Eigen::Vector3d coords0;
463 Eigen::Vector3d coords1;
464 Eigen::Vector3d coords2;
465 Eigen::Vector3d coords3;
466 coords0 << 1.0, 0.0, 0.0;
467 coords1 << 0.0, 0.0, 0.0;
468 coords2 << 2.0, 2.0, 0.0;
469 coords3 << 0.0, 1.0, 0.0;
470 Vertex &v0 = mesh.createVertex(coords0);
471 Vertex &v1 = mesh.createVertex(coords1);
472 Vertex &v2 = mesh.createVertex(coords2);
473 mesh.createVertex(coords3);
474 BOOST_TEST(mesh.nVertices() == 4);
475
476 std::array<Vertex *, 3> ptrs{&v0, &v2, &v1};
477 std::array<Eigen::VectorXd, 3> expected{coords0, coords2, coords1};
478
479 auto result = coordsFor(ptrs);
480 BOOST_REQUIRE(result.size() == 3);
481 BOOST_TEST(result == expected);
482}
483
484BOOST_AUTO_TEST_CASE(Integrate2DScalarData)
485{
486 PRECICE_TEST(1_rank);
487 PtrMesh mesh = std::make_shared<Mesh>("Mesh1", 2, testing::nextMeshID());
488 mesh->createData("Data", 1, 0_dataID);
489
490 auto &v1 = mesh->createVertex(Eigen::Vector2d(0.0, 0.0));
491 auto &v2 = mesh->createVertex(Eigen::Vector2d(0.0, 0.5));
492 auto &v3 = mesh->createVertex(Eigen::Vector2d(0.0, 1.5));
493 auto &v4 = mesh->createVertex(Eigen::Vector2d(0.0, 3.5));
494 mesh->allocateDataValues();
495
496 mesh->createEdge(v1, v2); // Length = 0.5
497 mesh->createEdge(v2, v3); // Length = 1.0
498 mesh->createEdge(v3, v4); // Length = 2.0
499
500 mesh->data(0)->values()(0) = 1.0;
501 mesh->data(0)->values()(1) = 3.0;
502 mesh->data(0)->values()(2) = 5.0;
503 mesh->data(0)->values()(3) = 7.0;
504
505 auto result = mesh::integrateSurface(mesh, mesh->data(0)->values());
506 double expected = 17.0;
507 BOOST_REQUIRE(result.size() == 1);
508 BOOST_TEST(result(0) == expected);
509}
510
511BOOST_AUTO_TEST_CASE(Integrate2DVectorData)
512{
513 PRECICE_TEST(1_rank);
514 PtrMesh mesh = std::make_shared<Mesh>("Mesh1", 2, testing::nextMeshID());
515 mesh->createData("Data", 2, 0_dataID);
516
517 auto &v1 = mesh->createVertex(Eigen::Vector2d(0.0, 0.0));
518 auto &v2 = mesh->createVertex(Eigen::Vector2d(0.0, 0.5));
519 auto &v3 = mesh->createVertex(Eigen::Vector2d(0.0, 1.5));
520 auto &v4 = mesh->createVertex(Eigen::Vector2d(0.0, 3.5));
521 mesh->allocateDataValues();
522
523 mesh->createEdge(v1, v2); // Length = 0.5
524 mesh->createEdge(v2, v3); // Length = 1.0
525 mesh->createEdge(v3, v4); // Length = 2.0
526
527 mesh->data(0)->values()(0) = 1.0;
528 mesh->data(0)->values()(1) = 2.0;
529 mesh->data(0)->values()(2) = 3.0;
530 mesh->data(0)->values()(3) = 4.0;
531 mesh->data(0)->values()(4) = 5.0;
532 mesh->data(0)->values()(5) = 6.0;
533 mesh->data(0)->values()(6) = 7.0;
534 mesh->data(0)->values()(7) = 8.0;
535
536 auto result = mesh::integrateSurface(mesh, mesh->data(0)->values());
537 Eigen::Vector2d expected(17.0, 20.5);
538 BOOST_REQUIRE(result.size() == 2);
539 BOOST_TEST(result(0) == expected(0));
540 BOOST_TEST(result(1) == expected(1));
541}
542
543BOOST_AUTO_TEST_CASE(Integrate3DScalarData)
544{
545 PRECICE_TEST(1_rank);
546 PtrMesh mesh = std::make_shared<Mesh>("Mesh1", 3, testing::nextMeshID());
547 mesh->createData("Data", 1, 0_dataID);
548
549 auto &v1 = mesh->createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
550 auto &v2 = mesh->createVertex(Eigen::Vector3d(3.0, 0.0, 0.0));
551 auto &v3 = mesh->createVertex(Eigen::Vector3d(3.0, 4.0, 0.0));
552 auto &v4 = mesh->createVertex(Eigen::Vector3d(0.0, 8.0, 0.0));
553 mesh->allocateDataValues();
554
555 auto &e1 = mesh->createEdge(v1, v2);
556 auto &e2 = mesh->createEdge(v2, v3);
557 auto &e3 = mesh->createEdge(v3, v4);
558 auto &e4 = mesh->createEdge(v1, v4);
559 auto &e5 = mesh->createEdge(v1, v3);
560
561 mesh->createTriangle(e1, e2, e5); // Area = 6.0
562 mesh->createTriangle(e3, e4, e5); // Area = 12.0
563
564 mesh->data(0)->values()(0) = 1.0;
565 mesh->data(0)->values()(1) = 3.0;
566 mesh->data(0)->values()(2) = 5.0;
567 mesh->data(0)->values()(3) = 7.0;
568
569 auto result = mesh::integrateSurface(mesh, mesh->data(0)->values());
570 double expected = 70.0;
571 BOOST_REQUIRE(result.size() == 1);
572 BOOST_TEST(result(0) == expected);
573}
574
575BOOST_AUTO_TEST_CASE(Integrate3DVectorData)
576{
577 PRECICE_TEST(1_rank);
578 PtrMesh mesh = std::make_shared<Mesh>("Mesh1", 3, testing::nextMeshID());
579 mesh->createData("Data", 2, 0_dataID);
580
581 auto &v1 = mesh->createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
582 auto &v2 = mesh->createVertex(Eigen::Vector3d(3.0, 0.0, 0.0));
583 auto &v3 = mesh->createVertex(Eigen::Vector3d(3.0, 4.0, 0.0));
584 auto &v4 = mesh->createVertex(Eigen::Vector3d(0.0, 8.0, 0.0));
585 mesh->allocateDataValues();
586
587 auto &e1 = mesh->createEdge(v1, v2);
588 auto &e2 = mesh->createEdge(v2, v3);
589 auto &e3 = mesh->createEdge(v3, v4);
590 auto &e4 = mesh->createEdge(v1, v4);
591 auto &e5 = mesh->createEdge(v1, v3);
592
593 mesh->createTriangle(e1, e2, e5); // Area = 6.0
594 mesh->createTriangle(e3, e4, e5); // Area = 12.0
595
596 mesh->data(0)->values()(0) = 1.0;
597 mesh->data(0)->values()(1) = 2.0;
598 mesh->data(0)->values()(2) = 3.0;
599 mesh->data(0)->values()(3) = 4.0;
600 mesh->data(0)->values()(4) = 5.0;
601 mesh->data(0)->values()(5) = 6.0;
602 mesh->data(0)->values()(6) = 7.0;
603 mesh->data(0)->values()(7) = 8.0;
604
605 auto result = mesh::integrateSurface(mesh, mesh->data(0)->values());
606 Eigen::Vector2d expected(70.0, 88.0);
607 BOOST_REQUIRE(result.size() == 2);
608 BOOST_TEST(result(0) == expected(0));
609 BOOST_TEST(result(1) == expected(1));
610}
611
613
614BOOST_AUTO_TEST_SUITE(VolumeIntegrals)
615
617 Eigen::Vector2d x0{0.0, 0.0};
618 Eigen::Vector2d x1{1.0, 0.0};
619 Eigen::Vector2d x2{1.0, 1.0};
620 Eigen::Vector2d x3{0.0, 1.0};
621};
622
623BOOST_FIXTURE_TEST_CASE(Integrate2DScalarDataVolume, UnitSquareFixture)
624{
625 PRECICE_TEST(1_rank);
626 PtrMesh mesh = std::make_shared<Mesh>("Mesh1", 2, testing::nextMeshID());
627 mesh->createData("Data", 1, 0_dataID);
628
629 auto &v1 = mesh->createVertex(x0);
630 auto &v2 = mesh->createVertex(x1);
631 auto &v3 = mesh->createVertex(x2);
632 auto &v4 = mesh->createVertex(x3);
633 mesh->allocateDataValues();
634
635 auto &e1 = mesh->createEdge(v1, v2);
636 auto &e2 = mesh->createEdge(v2, v3);
637 auto &e3 = mesh->createEdge(v3, v4);
638 auto &e4 = mesh->createEdge(v1, v4);
639 auto &e5 = mesh->createEdge(v1, v3);
640
641 // 2 triangles as halves of a unit square
642 mesh->createTriangle(e1, e2, e5);
643 mesh->createTriangle(e3, e4, e5);
644
645 BOOST_REQUIRE(mesh->triangles()[0].getArea() == 0.5);
646 BOOST_REQUIRE(mesh->triangles()[1].getArea() == 0.5);
647 BOOST_REQUIRE(mesh->triangles().size() == 2);
648
649 // Integrand is 1 + 4x + 2y integrated over unit square: 4
650 mesh->data(0)->values()(0) = 1.0;
651 mesh->data(0)->values()(1) = 3.0;
652 mesh->data(0)->values()(2) = 7.0;
653 mesh->data(0)->values()(3) = 5.0;
654
655 auto result = mesh::integrateVolume(mesh, mesh->data(0)->values());
656 double expected = 4.0;
657 BOOST_REQUIRE(result.size() == 1);
658 BOOST_TEST(result(0) == expected);
659}
660
661BOOST_FIXTURE_TEST_CASE(Integrate2DVectorDataVolume, UnitSquareFixture)
662{
663 PRECICE_TEST(1_rank);
664 PtrMesh mesh = std::make_shared<Mesh>("Mesh1", 2, testing::nextMeshID());
665 mesh->createData("Data", 2, 0_dataID);
666
667 auto &v1 = mesh->createVertex(x0);
668 auto &v2 = mesh->createVertex(x1);
669 auto &v3 = mesh->createVertex(x2);
670 auto &v4 = mesh->createVertex(x3);
671 mesh->allocateDataValues();
672
673 auto &e1 = mesh->createEdge(v1, v2);
674 auto &e2 = mesh->createEdge(v2, v3);
675 auto &e3 = mesh->createEdge(v3, v4);
676 auto &e4 = mesh->createEdge(v1, v4);
677 auto &e5 = mesh->createEdge(v1, v3);
678
679 // 2 triangles as halves of a unit square
680 mesh->createTriangle(e1, e2, e5);
681 mesh->createTriangle(e3, e4, e5);
682
683 BOOST_REQUIRE(mesh->triangles()[0].getArea() == 0.5);
684 BOOST_REQUIRE(mesh->triangles()[1].getArea() == 0.5);
685 BOOST_REQUIRE(mesh->triangles().size() == 2);
686
687 // Integrand of 1st component is 1 + 4x + 2y integrated over unit square: 4
688 mesh->data(0)->values()(0) = 1.0;
689 mesh->data(0)->values()(2) = 3.0;
690 mesh->data(0)->values()(4) = 7.0;
691 mesh->data(0)->values()(6) = 5.0;
692 // Integrand of 1st component is 1 + 4x + 2y integrated over unit square: 4
693 mesh->data(0)->values()(1) = 2.0;
694 mesh->data(0)->values()(3) = 4.0;
695 mesh->data(0)->values()(5) = 8.0;
696 mesh->data(0)->values()(7) = 6.0;
697
698 auto result = mesh::integrateVolume(mesh, mesh->data(0)->values());
699 Eigen::Vector2d expected(4.0, 5.0);
700 BOOST_REQUIRE(result.size() == 2);
701 BOOST_TEST(result(0) == expected(0));
702 BOOST_TEST(result(1) == expected(1));
703}
704
706 Eigen::Vector3d x1{0.0, 0.0, 0.0};
707 Eigen::Vector3d x2{1.0, 0.0, 0.0};
708 Eigen::Vector3d x3{0.0, 1.0, 0.0};
709 Eigen::Vector3d x4{0.0, 0.0, 1.0};
710};
711
712BOOST_FIXTURE_TEST_CASE(Integrate3DScalarDataVolume, OneTetraFixture)
713{
714 PRECICE_TEST(1_rank);
715 PtrMesh mesh = std::make_shared<Mesh>("Mesh1", 3, testing::nextMeshID());
716 mesh->createData("Data", 1, 0_dataID);
717
718 auto &v1 = mesh->createVertex(x1);
719 auto &v2 = mesh->createVertex(x2);
720 auto &v3 = mesh->createVertex(x3);
721 auto &v4 = mesh->createVertex(x4);
722
723 mesh->allocateDataValues();
724
725 mesh->createTetrahedron(v1, v2, v3, v4);
726
727 BOOST_REQUIRE(mesh->tetrahedra()[0].getVolume() == 1. / 6);
728 BOOST_REQUIRE(mesh->tetrahedra().size() == 1);
729
730 // Integrand is 1 + 2x + 4y + 6z integrated over one tetra
731 mesh->data(0)->values()(0) = 1.0;
732 mesh->data(0)->values()(1) = 3.0;
733 mesh->data(0)->values()(2) = 5.0;
734 mesh->data(0)->values()(3) = 7.0;
735
736 auto result = mesh::integrateVolume(mesh, mesh->data(0)->values());
737 double expected = 4.0 / 6;
738 BOOST_REQUIRE(result.size() == 1);
739 BOOST_TEST(result(0) == expected);
740}
741
742BOOST_FIXTURE_TEST_CASE(Integrate3DVectorDataVolume, OneTetraFixture)
743{
744 PRECICE_TEST(1_rank);
745 PtrMesh mesh = std::make_shared<Mesh>("Mesh1", 3, testing::nextMeshID());
746 mesh->createData("Data", 3, 0_dataID);
747
748 auto &v1 = mesh->createVertex(x1);
749 auto &v2 = mesh->createVertex(x2);
750 auto &v3 = mesh->createVertex(x3);
751 auto &v4 = mesh->createVertex(x4);
752
753 mesh->allocateDataValues();
754
755 mesh->createTetrahedron(v1, v2, v3, v4);
756
757 BOOST_REQUIRE(mesh->tetrahedra()[0].getVolume() == 1. / 6);
758 BOOST_REQUIRE(mesh->tetrahedra().size() == 1);
759
760 // Integrand is (1 + 2x + 4y + 6z, 1, 1-x-y-z) integrated over one tetra
761 mesh->data(0)->values()(0 * 3) = 1.0;
762 mesh->data(0)->values()(1 * 3) = 3.0;
763 mesh->data(0)->values()(2 * 3) = 5.0;
764 mesh->data(0)->values()(3 * 3) = 7.0;
765
766 mesh->data(0)->values()(0 * 3 + 1) = 1.0;
767 mesh->data(0)->values()(1 * 3 + 1) = 1.0;
768 mesh->data(0)->values()(2 * 3 + 1) = 1.0;
769 mesh->data(0)->values()(3 * 3 + 1) = 1.0;
770
771 mesh->data(0)->values()(0 * 3 + 2) = 1.0;
772 mesh->data(0)->values()(1 * 3 + 2) = 0.0;
773 mesh->data(0)->values()(2 * 3 + 2) = 0.0;
774 mesh->data(0)->values()(3 * 3 + 2) = 0.0;
775
776 auto result = mesh::integrateVolume(mesh, mesh->data(0)->values());
777 Eigen::Vector3d expected(4.0 / 6, 1.0 / 6, 1.0 / 24);
778 BOOST_REQUIRE(result.size() == 3);
779 BOOST_TEST(result(0) == expected(0));
780 BOOST_TEST(result(1) == expected(1));
781 BOOST_TEST(result(2) == expected(2));
782}
783
784BOOST_AUTO_TEST_SUITE_END() // VolumeIntegrals
785
787{
788 PRECICE_TEST(1_rank);
789 PtrMesh globalMesh = std::make_shared<Mesh>("Mesh1", 3, testing::nextMeshID());
790 PtrMesh subMesh = std::make_shared<Mesh>("Mesh2", 3, testing::nextMeshID());
791
792 // Fill globalMesh, then subMesh, then add and check the result is the sul
793 auto &v01 = globalMesh->createVertex(Eigen::Vector3d{0.0, 0.0, 0.0});
794 auto &v02 = globalMesh->createVertex(Eigen::Vector3d{1.0, 0.0, 0.0});
795 auto &v03 = globalMesh->createVertex(Eigen::Vector3d{0.0, 1.0, 0.0});
796 auto &v04 = globalMesh->createVertex(Eigen::Vector3d{0.0, 0.0, 1.0});
797
798 // Add some elements of all kind to global
799 globalMesh->createEdge(v01, v02);
800 globalMesh->createTriangle(v01, v02, v03);
801 globalMesh->createTetrahedron(v01, v02, v03, v04);
802
803 auto &v11 = subMesh->createVertex(Eigen::Vector3d{0.0, 0.0, 0.0});
804 auto &v12 = subMesh->createVertex(Eigen::Vector3d{2.0, 0.0, 0.0});
805 auto &v13 = subMesh->createVertex(Eigen::Vector3d{0.0, 4.0, 0.0});
806 auto &v14 = subMesh->createVertex(Eigen::Vector3d{0.0, 0.0, 3.0});
807 auto &v15 = subMesh->createVertex(Eigen::Vector3d{5.0, 0.0, 1.0});
808
809 // Add some elements of all kind to sub mesh
810 subMesh->createEdge(v11, v12);
811 subMesh->createEdge(v14, v15);
812 subMesh->createTriangle(v11, v13, v15);
813 subMesh->createTriangle(v11, v13, v14);
814 subMesh->createTetrahedron(v11, v12, v13, v14);
815 subMesh->createTetrahedron(v15, v12, v13, v14);
816
817 globalMesh->addMesh(*subMesh);
818 BOOST_TEST(globalMesh->nVertices() == 9);
819 BOOST_TEST(globalMesh->edges().size() == 3);
820 BOOST_TEST(globalMesh->triangles().size() == 3);
821 BOOST_TEST(globalMesh->tetrahedra().size() == 3);
822}
823
825
826BOOST_AUTO_TEST_CASE(DuplicateEdges)
827{
828 PRECICE_TEST(1_rank);
829 Mesh mesh{"Mesh1", 3, 0};
830
831 auto &v1 = mesh.createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
832 auto &v2 = mesh.createVertex(Eigen::Vector3d(3.0, 0.0, 0.0));
833 auto &v3 = mesh.createVertex(Eigen::Vector3d(3.0, 4.0, 0.0));
834 auto &v4 = mesh.createVertex(Eigen::Vector3d(0.0, 8.0, 0.0));
835
836 mesh.createEdge(v1, v2);
837 mesh.createEdge(v2, v3);
838 mesh.createEdge(v3, v4);
839
840 // add single duplicate
841 mesh.createEdge(v3, v4);
842
843 // add 1000 duplicates
844 for (int i = 0; i < 1000; ++i) {
845 mesh.createEdge(v2, v3);
846 };
847 BOOST_TEST(mesh.edges().size() == 1004);
848
849 mesh.preprocess();
850
851 BOOST_TEST(mesh.edges().size() == 3);
852 BOOST_TEST(mesh.triangles().empty());
853 BOOST_TEST(mesh.tetrahedra().empty());
854
855 std::vector<Edge> expectedEdges{{v1, v2}, {v2, v3}, {v3, v4}};
856 for (auto &e : expectedEdges) {
857 auto cnt = std::count(mesh.edges().begin(), mesh.edges().end(), e);
858 BOOST_TEST(cnt == 1);
859 }
860}
861
862BOOST_AUTO_TEST_CASE(SingleTriangle)
863{
864 PRECICE_TEST(1_rank);
865 Mesh mesh{"Mesh1", 3, 0};
866
867 auto &v1 = mesh.createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
868 auto &v2 = mesh.createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
869 auto &v3 = mesh.createVertex(Eigen::Vector3d(1.0, 1.0, 0.0));
870
871 mesh.createTriangle(v1, v2, v3);
872
873 mesh.preprocess();
874
875 BOOST_TEST(mesh.edges().size() == 3);
876 BOOST_TEST(mesh.triangles().size() == 1);
877 BOOST_TEST(mesh.tetrahedra().size() == 0);
878}
879
880BOOST_AUTO_TEST_CASE(SingleTriangulatedQuad)
881{
882 PRECICE_TEST(1_rank);
883 Mesh mesh{"Mesh1", 3, 0};
884
885 auto &v1 = mesh.createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
886 auto &v2 = mesh.createVertex(Eigen::Vector3d(1.0, 0.0, 0.0));
887 auto &v3 = mesh.createVertex(Eigen::Vector3d(1.0, 1.0, 0.0));
888 auto &v4 = mesh.createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
889
890 mesh.createTriangle(v1, v2, v3);
891 mesh.createTriangle(v2, v3, v4);
892
893 mesh.preprocess();
894
895 BOOST_TEST(mesh.edges().size() == 5);
896 BOOST_TEST(mesh.triangles().size() == 2);
897 BOOST_TEST(mesh.tetrahedra().size() == 0);
898}
899
900BOOST_AUTO_TEST_CASE(SingleTetrahedron)
901{
902 PRECICE_TEST(1_rank);
903 Mesh mesh{"Mesh1", 3, 0};
904
905 auto &v1 = mesh.createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
906 auto &v2 = mesh.createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
907 auto &v3 = mesh.createVertex(Eigen::Vector3d(1.0, 1.0, 0.0));
908 auto &v4 = mesh.createVertex(Eigen::Vector3d(0.3, 0.3, 1.0));
909
910 mesh.createTetrahedron(v1, v2, v3, v4);
911
912 mesh.preprocess();
913
914 BOOST_TEST(mesh.edges().size() == 6);
915 BOOST_TEST(mesh.triangles().size() == 4);
916 BOOST_TEST(mesh.tetrahedra().size() == 1);
917}
918
919BOOST_AUTO_TEST_CASE(TouchingTetrahedra)
920{
921 PRECICE_TEST(1_rank);
922 Mesh mesh{"Mesh1", 3, 0};
923
924 auto &v1 = mesh.createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
925 auto &v2 = mesh.createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
926 auto &v3 = mesh.createVertex(Eigen::Vector3d(1.0, 1.0, 0.0));
927
928 auto &v4 = mesh.createVertex(Eigen::Vector3d(0.3, 0.3, 1.0));
929 auto &v5 = mesh.createVertex(Eigen::Vector3d(0.3, 0.3, -1.0));
930
931 mesh.createTetrahedron(v1, v2, v3, v4);
932 mesh.createTetrahedron(v1, v2, v3, v5);
933
934 mesh.preprocess();
935
936 BOOST_TEST(mesh.edges().size() == 9);
937 BOOST_TEST(mesh.triangles().size() == 7);
938 BOOST_TEST(mesh.tetrahedra().size() == 2);
939}
940
941BOOST_AUTO_TEST_CASE(SingleTriangle2ImplicitEdges)
942{
943 PRECICE_TEST(1_rank);
944 Mesh mesh{"Mesh1", 3, 0};
945
946 auto &v1 = mesh.createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
947 auto &v2 = mesh.createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
948 auto &v3 = mesh.createVertex(Eigen::Vector3d(1.0, 1.0, 0.0));
949
950 mesh.createTriangle(v1, v2, v3);
951 mesh.createEdge(v2, v3);
952
953 mesh.preprocess();
954
955 BOOST_TEST(mesh.edges().size() == 3);
956 BOOST_TEST(mesh.triangles().size() == 1);
957 BOOST_TEST(mesh.tetrahedra().size() == 0);
958}
959
960BOOST_AUTO_TEST_CASE(SingleTriangle1ImplicitEdges)
961{
962 PRECICE_TEST(1_rank);
963 Mesh mesh{"Mesh1", 3, 0};
964
965 auto &v1 = mesh.createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
966 auto &v2 = mesh.createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
967 auto &v3 = mesh.createVertex(Eigen::Vector3d(1.0, 1.0, 0.0));
968
969 mesh.createTriangle(v1, v2, v3);
970 mesh.createEdge(v1, v2);
971 mesh.createEdge(v2, v3);
972
973 mesh.preprocess();
974
975 BOOST_TEST(mesh.edges().size() == 3);
976 BOOST_TEST(mesh.triangles().size() == 1);
977 BOOST_TEST(mesh.tetrahedra().size() == 0);
978}
979
980BOOST_AUTO_TEST_CASE(SingleTetrahedron3ImplicitTriangles)
981{
982 PRECICE_TEST(1_rank);
983 Mesh mesh{"Mesh1", 3, 0};
984
985 auto &v1 = mesh.createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
986 auto &v2 = mesh.createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
987 auto &v3 = mesh.createVertex(Eigen::Vector3d(1.0, 1.0, 0.0));
988 auto &v4 = mesh.createVertex(Eigen::Vector3d(0.3, 0.3, 1.0));
989
990 mesh.createTetrahedron(v1, v2, v3, v4);
991 mesh.createTriangle(v1, v2, v3);
992
993 mesh.preprocess();
994
995 BOOST_TEST(mesh.edges().size() == 6);
996 BOOST_TEST(mesh.triangles().size() == 4);
997 BOOST_TEST(mesh.tetrahedra().size() == 1);
998}
999
1000BOOST_AUTO_TEST_CASE(SingleTetrahedronNoImplicitTriangles)
1001{
1002 PRECICE_TEST(1_rank);
1003 Mesh mesh{"Mesh1", 3, 0};
1004
1005 auto &v1 = mesh.createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
1006 auto &v2 = mesh.createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
1007 auto &v3 = mesh.createVertex(Eigen::Vector3d(1.0, 1.0, 0.0));
1008 auto &v4 = mesh.createVertex(Eigen::Vector3d(0.3, 0.3, 1.0));
1009
1010 mesh.createTetrahedron(v1, v2, v3, v4);
1011 mesh.createTriangle(v1, v2, v3);
1012 mesh.createTriangle(v1, v2, v4);
1013 mesh.createTriangle(v3, v4, v1);
1014 mesh.createTriangle(v3, v4, v2);
1015
1016 mesh.preprocess();
1017
1018 BOOST_TEST(mesh.edges().size() == 6);
1019 BOOST_TEST(mesh.triangles().size() == 4);
1020 BOOST_TEST(mesh.tetrahedra().size() == 1);
1021}
1022
1024{
1025 PRECICE_TEST(1_rank);
1026 Mesh mesh{"Mesh1", 3, 0};
1027
1028 auto &v1 = mesh.createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
1029 auto &v2 = mesh.createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
1030 auto &v3 = mesh.createVertex(Eigen::Vector3d(1.0, 1.0, 0.0));
1031 auto &v4 = mesh.createVertex(Eigen::Vector3d(0.3, 0.3, 1.0));
1032
1033 mesh.createTetrahedron(v1, v2, v3, v4);
1034 mesh.createTriangle(v1, v2, v3);
1035 mesh.createEdge(v3, v4);
1036
1037 mesh.preprocess();
1038
1039 BOOST_TEST(mesh.edges().size() == 6);
1040 BOOST_TEST(mesh.triangles().size() == 4);
1041 BOOST_TEST(mesh.tetrahedra().size() == 1);
1042}
1043
1045{
1046 PRECICE_TEST(1_rank);
1047 Mesh mesh{"Mesh1", 3, 0};
1048
1049 auto &v1 = mesh.createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
1050 auto &v2 = mesh.createVertex(Eigen::Vector3d(3.0, 0.0, 0.0));
1051 auto &v3 = mesh.createVertex(Eigen::Vector3d(3.0, 4.0, 0.0));
1052 auto &v4 = mesh.createVertex(Eigen::Vector3d(0.0, 8.0, 0.0));
1053
1054 mesh.createEdge(v1, v2);
1055 mesh.createEdge(v2, v3);
1056 mesh.createEdge(v3, v4);
1057
1058 // add single duplicate
1059 mesh.createEdge(v3, v4);
1060
1061 // add 1000 duplicates
1062 for (int i = 0; i < 1000; ++i) {
1063 mesh.createEdge(v2, v3);
1064 };
1065 BOOST_TEST(mesh.edges().size() == 1004);
1066
1067 mesh.createTriangle(v1, v2, v3);
1068 mesh.createTriangle(v2, v3, v4);
1069 for (int i = 0; i < 1000; ++i) {
1070 mesh.createTriangle(v2, v3, v4);
1071 };
1072 BOOST_TEST(mesh.triangles().size() == 1002);
1073
1074 // creates edges 1-3 and 2-4
1075 mesh.preprocess();
1076
1077 BOOST_TEST(mesh.edges().size() == 5);
1078 BOOST_TEST(mesh.triangles().size() == 2);
1079
1080 std::vector<Edge> expectedEdges{{v1, v2}, {v1, v3}, {v2, v3}, {v3, v4}, {v2, v4}};
1081 for (auto &e : expectedEdges) {
1082 auto cnt = std::count(mesh.edges().begin(), mesh.edges().end(), e);
1083 BOOST_TEST(cnt == 1);
1084 }
1085
1086 std::vector<Triangle> expectedTriangles{{v1, v2, v3}, {v2, v3, v4}};
1087 for (auto &t : expectedTriangles) {
1088 auto cnt = std::count(mesh.triangles().begin(), mesh.triangles().end(), t);
1089 BOOST_TEST(cnt == 1);
1090 }
1091}
1092
double const * ptr
Definition ExportCSV.cpp:23
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_CASE(BoundingBoxCOG_2D)
Definition MeshTest.cpp:32
BOOST_FIXTURE_TEST_CASE(Integrate2DScalarDataVolume, UnitSquareFixture)
Definition MeshTest.cpp:623
BOOST_AUTO_TEST_SUITE_END()
unsigned int index
#define PRECICE_TEST(...)
Definition Testing.hpp:27
An axis-aligned bounding box around a (partition of a) mesh.
int getDimension() const
Getter dimension of the bounding box.
Eigen::VectorXd center() const
Returns the Center Of Gravity of the mesh.
Eigen::VectorXd & values()
Returns a reference to the data values.
Definition Data.cpp:28
Linear edge of a mesh, defined by two Vertex objects.
Definition Edge.hpp:16
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
VertexContainer & vertices()
Returns modifieable container holding all vertices.
Definition Mesh.cpp:53
void clear()
Removes all mesh elements and data values (does not remove data or the bounding boxes).
Definition Mesh.cpp:256
bool hasEdges() const
Definition Mesh.hpp:99
void addMesh(Mesh &deltaMesh)
Definition Mesh.cpp:309
const std::string & getName() const
Returns the name of the mesh, as set in the config file.
Definition Mesh.cpp:218
std::size_t nVertices() const
Returns the number of vertices.
Definition Mesh.cpp:63
TetraContainer & tetrahedra()
Returns modifiable container holding all tetrahedra.
Definition Mesh.cpp:93
PtrData & createData(const std::string &name, int dimension, DataID id, int waveformDegree=time::Time::DEFAULT_WAVEFORM_DEGREE)
Create only data for vertex.
Definition Mesh.cpp:151
void computeBoundingBox()
Computes the boundingBox for the vertices.
Definition Mesh.cpp:242
Tetrahedron & createTetrahedron(Vertex &vertexOne, Vertex &vertexTwo, Vertex &vertexThree, Vertex &vertexFour)
Creates and initializes a Tetrahedron object.
Definition Mesh.cpp:141
const DataContainer & data() const
Allows access to all data.
Definition Mesh.cpp:170
TriangleContainer & triangles()
Returns modifiable container holding all triangles.
Definition Mesh.cpp:78
const BoundingBox & getBoundingBox() const
Returns the bounding box of the mesh.
Definition Mesh.cpp:365
Edge & createEdge(Vertex &vertexOne, Vertex &vertexTwo)
Creates and initializes an Edge object.
Definition Mesh.cpp:111
EdgeContainer & edges()
Returns modifiable container holding all edges.
Definition Mesh.cpp:68
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
bool hasTriangles() const
Definition Mesh.hpp:110
Triangle of a mesh, defined by three vertices.
Definition Triangle.hpp:27
Vertex of a mesh.
Definition Vertex.hpp:16
VertexID getID() const
Returns the unique (among vertices of one mesh on one processor) ID of the vertex.
Definition Vertex.hpp:111
T count(T... args)
provides Mesh, Data and primitives.
Eigen::VectorXd integrateSurface(const PtrMesh &mesh, const Eigen::VectorXd &input)
Given the data and the mesh, this function returns the surface integral. Assumes no overlap exists fo...
Definition Utils.cpp:11
Eigen::VectorXd integrateVolume(const PtrMesh &mesh, const Eigen::VectorXd &input)
Given the data and the mesh, this function returns the volume integral. Assumes no overlap exists for...
Definition Utils.cpp:41
std::array< Eigen::VectorXd, n > coordsFor(const Mesh &mesh, const std::array< int, n > &vertexIDs)
Given a mesh and an array of vertexIDS, this function returns an array of coordinates of the vertices...
Definition Utils.hpp:123
Vertex * sharedVertex(Edge &a, Edge &b)
Definition Utils.hpp:26
std::array< Vertex *, n > vertexPtrsFor(Mesh &mesh, const std::array< int, n > &vertexIDs)
Given a mesh and an array of vertexIDS, this function returns an array of pointers to vertices.
Definition Utils.hpp:112
double edgeLength(const Edge &e)
Definition Utils.hpp:47
Chain< n > asChain(std::array< mesh::Edge *, n > edges)
Definition Utils.hpp:73
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
auto make_array(Elements &&... elements) -> std::array< typename std::common_type< Elements... >::type, sizeof...(Elements)>
Function that generates an array from given elements.
Definition algorithm.hpp:51
Main namespace of the precice library.
T size(T... args)
T sqrt(T... args)
T str(T... args)
Eigen::Vector3d x3
Definition MeshTest.cpp:708
Eigen::Vector3d x2
Definition MeshTest.cpp:707
Eigen::Vector3d x4
Definition MeshTest.cpp:709
Eigen::Vector3d x1
Definition MeshTest.cpp:706