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