Wildmeshing Toolkit
Loading...
Searching...
No Matches
test_component_isotropic_remeshing.cpp
Go to the documentation of this file.
1#include <catch2/catch_test_macros.hpp>
2#include <nlohmann/json.hpp>
3#include <tools/DEBUG_EdgeMesh.hpp>
4#include <tools/DEBUG_TriMesh.hpp>
5#include <tools/EdgeMesh_examples.hpp>
6#include <tools/TriMesh_examples.hpp>
7#include <wmtk/Scheduler.hpp>
8#include <wmtk/Types.hpp>
10#include <wmtk/components/isotropic_remeshing/internal/IsotropicRemeshing.hpp>
11#include <wmtk/components/isotropic_remeshing/internal/IsotropicRemeshingOptions.hpp>
14#include <wmtk/components/utils/Paths.hpp>
20#include <wmtk/io/Cache.hpp>
33#include <wmtk/simplex/link.hpp>
35
36using namespace wmtk::components::utils;
37
38using json = nlohmann::json;
39using namespace wmtk;
40using namespace wmtk::tests;
41using namespace wmtk::simplex;
42
43const std::filesystem::path data_dir = WMTK_DATA_DIR;
44
45void print_tuple_map_iso(const DEBUG_TriMesh& parent, const DEBUG_MultiMeshManager& p_mul_manager)
46{
47 int64_t child_id = 0;
48 for (auto& child_data : p_mul_manager.children()) {
49 PrimitiveType map_ptype = child_data.mesh->top_simplex_type();
50 auto parent_to_child_accessor =
51 parent.create_const_accessor<int64_t>(child_data.map_handle);
52 for (int64_t parent_gid = 0; parent_gid < parent.capacity(map_ptype); ++parent_gid) {
53 auto parent_to_child_data = parent_to_child_accessor.const_vector_attribute(
54 parent.tuple_from_id(map_ptype, parent_gid));
55 auto [parent_tuple, child_tuple] =
57 }
58 }
59}
60
61TEST_CASE("smoothing_mesh", "[components][isotropic_remeshing][2D]")
62{
63 using namespace operations;
64
65 wmtk::io::Cache cache("wmtk_cache", ".");
66
67 // input
68 json input_component_json = {
69 {"type", "input"},
70 {"name", "input_mesh"},
71 {"cell_dimension", 2},
72 {"file", (data_dir / "bumpyDice.msh").string()},
73 {"ignore_z", false},
74 {"tetrahedron_attributes", json::array()}};
75 wmtk::components::input(Paths(), input_component_json, cache);
76
77
78 auto mesh_in = cache.read_mesh(input_component_json["name"]);
79
80 TriMesh& m = static_cast<TriMesh&>(*mesh_in);
81
83 m.get_attribute_handle<double>("vertices", PrimitiveType::Vertex);
84
86 op.set_function(VertexLaplacianSmooth(pos_attribute));
87
88 op.add_invariant(
89 std::make_shared<invariants::InteriorSimplexInvariant>(m, PrimitiveType::Vertex));
90
91 Scheduler scheduler;
92
93 for (int i = 0; i < 3; ++i) {
94 scheduler.run_operation_on_all(op);
95 }
96
97 // output
98 {
99 ParaviewWriter writer(
100 cache.get_cache_path() / "mesh_smooth",
101 "vertices",
102 m,
103 false,
104 false,
105 true,
106 false);
107 m.serialize(writer);
108 }
109}
110
111TEST_CASE("smoothing_simple_examples", "[components][isotropic_remeshing][2D]")
112{
113 using namespace operations;
114 using namespace tri_mesh;
115
116 SECTION("hex_plus_two")
117 {
118 DEBUG_TriMesh mesh = wmtk::tests::hex_plus_two_with_position();
119
121 mesh.get_attribute_handle<double>("vertices", PrimitiveType::Vertex);
122
123 // offset interior vertex
124 auto pos = mesh.create_accessor<double>(pos_attribute);
125 Tuple v4 = mesh.tuple_from_id(PrimitiveType::Vertex, 4);
126 pos.vector_attribute(v4) = Eigen::Vector3d{0.6, 0.9, 0};
127
129 op.set_function(VertexLaplacianSmooth(pos_attribute));
130 op.add_invariant(
131 std::make_shared<invariants::InteriorSimplexInvariant>(mesh, PrimitiveType::Vertex));
132
133 Scheduler scheduler;
134 scheduler.run_operation_on_all(op);
135
136 v4 = mesh.tuple_from_id(PrimitiveType::Vertex, 4);
137 Eigen::Vector3d after_smooth = pos.vector_attribute(v4);
138 CHECK((after_smooth - Eigen::Vector3d{1, 0, 0}).squaredNorm() == 0);
139 }
140
141 SECTION("edge_region")
142 {
143 DEBUG_TriMesh mesh = wmtk::tests::edge_region_with_position();
144
146 mesh.get_attribute_handle<double>("vertices", PrimitiveType::Vertex);
147
148 // offset interior vertex
149 auto pos = mesh.create_accessor<double>(pos_attribute);
150 Tuple v4 = mesh.tuple_from_id(PrimitiveType::Vertex, 4);
151 Tuple v5 = mesh.tuple_from_id(PrimitiveType::Vertex, 5);
152 pos.vector_attribute(v4) = Eigen::Vector3d{0.6, 0.9, 0};
153 pos.vector_attribute(v5) = Eigen::Vector3d{1.4, -0.9, 0};
154
156 op.set_function(VertexLaplacianSmooth(pos_attribute));
157 op.add_invariant(
158 std::make_shared<invariants::InteriorSimplexInvariant>(mesh, PrimitiveType::Vertex));
159
160 Scheduler scheduler;
161
162 for (size_t i = 0; i < 10; ++i) {
163 scheduler.run_operation_on_all(op);
164 // v4 = mesh.tuple_from_id(PrimitiveType::Vertex, 4);
165 // /*Eigen::Vector3d p4_after_smooth =*/pos.vector_attribute(v4);
166 }
167
168 v4 = mesh.tuple_from_id(PrimitiveType::Vertex, 4);
169 Eigen::Vector3d p4_after_smooth = pos.vector_attribute(v4);
170 CHECK((p4_after_smooth - Eigen::Vector3d{1, 0, 0}).squaredNorm() < 1e-10);
171
172 v5 = mesh.tuple_from_id(PrimitiveType::Vertex, 5);
173 Eigen::Vector3d p5_after_smooth = pos.vector_attribute(v5);
174 CHECK((p5_after_smooth - Eigen::Vector3d{2, 0, 0}).squaredNorm() < 1e-10);
175 }
176}
177
178TEST_CASE("tangential_smoothing", "[components][isotropic_remeshing][2D]")
179{
180 using namespace operations;
181
182 DEBUG_TriMesh mesh = wmtk::tests::hex_plus_two_with_position();
183
185 mesh.get_attribute_handle<double>("vertices", PrimitiveType::Vertex);
186
187 // offset interior vertex
188 auto pos = mesh.create_accessor<double>(pos_attribute);
189 Tuple v4 = mesh.tuple_from_id(PrimitiveType::Vertex, 4);
190
191 Eigen::Vector3d p_init;
192 SECTION("1_0_1")
193 {
194 p_init = Eigen::Vector3d{1, 0, 1};
195 }
196 SECTION("0.5_0.5_1")
197 {
198 p_init = Eigen::Vector3d{0.5, 0.5, 1};
199 }
200 SECTION("0_0_7")
201 {
202 p_init = Eigen::Vector3d{0, 0, 7};
203 }
204
205 pos.vector_attribute(v4) = p_init;
206
209 op.add_invariant(
210 std::make_shared<invariants::InteriorSimplexInvariant>(mesh, PrimitiveType::Vertex));
211
212 Scheduler scheduler;
213 scheduler.run_operation_on_all(op);
214
215 v4 = mesh.tuple_from_id(PrimitiveType::Vertex, 4);
216 Eigen::Vector3d after_smooth = pos.vector_attribute(v4);
217 Eigen::Vector3d target = Eigen::Vector3d{1, 0, p_init[2]};
218 CHECK((after_smooth - target).squaredNorm() == 0);
219}
220
221TEST_CASE("tangential_smoothing_boundary", "[components][isotropic_remeshing][2D]")
222{
223 using namespace operations;
224 using namespace tri_mesh;
225
226 DEBUG_TriMesh mesh = wmtk::tests::hex_plus_two_with_position();
227
229 mesh.get_attribute_handle<double>("vertices", PrimitiveType::Vertex);
230
231 // offset interior vertex
232 auto pos = mesh.create_accessor<double>(pos_attribute);
233 Tuple v1 = mesh.tuple_from_id(PrimitiveType::Vertex, 1);
234
235 Eigen::Vector3d p_init;
236 SECTION("1.7_1.1_0")
237 {
238 p_init = Eigen::Vector3d{1.7, 1.1, 0};
239 }
240 SECTION("2.2_2_0")
241 {
242 p_init = Eigen::Vector3d{2.2, 2, 0};
243 }
244 SECTION("2.2_2_5")
245 {
246 p_init = Eigen::Vector3d{2.2, 2, 5};
247 }
248
249 pos.vector_attribute(v1) = p_init;
250
253
254 const bool success = !op(Simplex::vertex(mesh, v1)).empty();
255 REQUIRE(success);
256
257 v1 = mesh.tuple_from_id(PrimitiveType::Vertex, 1);
258 Eigen::Vector3d after_smooth = pos.vector_attribute(v1);
259 CHECK((after_smooth - Eigen::Vector3d{1.5, p_init[1], p_init[2]}).squaredNorm() == 0);
260}
261
262TEST_CASE("split_long_edges", "[components][isotropic_remeshing][split][2D]")
263{
264 using namespace operations;
265
266 // This test does not fully work yet
267
268 DEBUG_TriMesh mesh = wmtk::tests::edge_region_with_position();
269
271 mesh.get_attribute_handle<double>("vertices", PrimitiveType::Vertex);
272
273 EdgeSplit op(mesh);
274 op.set_new_attribute_strategy(
275 pos_attribute,
276 SplitBasicStrategy::None,
277 SplitRibBasicStrategy::Mean);
278
279 {
280 auto pos = mesh.create_accessor<double>(pos_attribute);
281 const Tuple v4 = mesh.tuple_from_id(PrimitiveType::Vertex, 4);
282 const Tuple v5 = mesh.tuple_from_id(PrimitiveType::Vertex, 5);
283 // reposition interior vertices
284 pos.vector_attribute(v4) = Eigen::Vector3d{0.6, 0.9, 0};
285 pos.vector_attribute(v5) = Eigen::Vector3d{2.4, -0.9, 0};
286 // << std::endl;
287 }
288
289 double min_split_length_squared = -1;
290
291 SECTION("6.4")
292 {
293 min_split_length_squared = 6.4;
294 op.add_invariant(std::make_shared<MinEdgeLengthInvariant>(
295 mesh,
296 pos_attribute.as<double>(),
297 min_split_length_squared));
298
299 Scheduler scheduler;
300
301 size_t n_vertices = mesh.get_all(PrimitiveType::Vertex).size();
302 size_t n_iterations = 0;
303 for (; n_iterations < 10; ++n_iterations) {
304 scheduler.run_operation_on_all(op);
305
306 const size_t n_vertices_new = mesh.get_all(PrimitiveType::Vertex).size();
307 if (n_vertices_new == n_vertices) {
308 break;
309 } else {
310 n_vertices = n_vertices_new;
311 }
312 }
313
314 CHECK(n_iterations == 1);
315 REQUIRE(n_vertices == 11);
316
317 // check position of new vertex
318 auto pos = mesh.create_accessor<double>(pos_attribute);
319 const Tuple v10 = mesh.tuple_from_id(PrimitiveType::Vertex, 10);
320 CHECK((pos.vector_attribute(v10) - Eigen::Vector3d{1.5, 0, 0}).squaredNorm() == 0);
321 }
322 SECTION("3.5")
323 {
324 min_split_length_squared = 3.5;
325 op.add_invariant(std::make_shared<MinEdgeLengthInvariant>(
326 mesh,
327 pos_attribute.as<double>(),
328 min_split_length_squared));
329
330 Scheduler scheduler;
331
332 size_t n_vertices = mesh.get_all(PrimitiveType::Vertex).size();
333 size_t n_iterations = 0;
334 for (; n_iterations < 10; ++n_iterations) {
335 scheduler.run_operation_on_all(op);
336
337 const size_t n_vertices_new = mesh.get_all(PrimitiveType::Vertex).size();
338 if (n_vertices_new == n_vertices) {
339 break;
340 } else {
341 n_vertices = n_vertices_new;
342 }
343 }
344
345 CHECK(n_iterations < 5);
346 CHECK(n_vertices == 15);
347 }
348
349 // check edge lengths
350 auto pos = mesh.create_accessor<double>(pos_attribute);
351 for (const Tuple& e : mesh.get_all(PrimitiveType::Edge)) {
352 const Eigen::Vector3d p0 = pos.vector_attribute(e);
353 const Eigen::Vector3d p1 = pos.vector_attribute(mesh.switch_vertex(e));
354 const double l_squared = (p1 - p0).squaredNorm();
355 CHECK(l_squared < min_split_length_squared);
356 }
357}
358
359TEST_CASE("collapse_short_edges", "[components][isotropic_remeshing][collapse][2D]")
360{
361 using namespace operations;
362
363 DEBUG_TriMesh mesh = wmtk::tests::edge_region_with_position();
364
365 auto pos_attribute = mesh.get_attribute_handle<double>("vertices", PrimitiveType::Vertex);
366
367
368 EdgeCollapse op(mesh);
369 op.add_invariant(std::make_shared<MultiMeshLinkConditionInvariant>(mesh));
370 op.set_new_attribute_strategy(pos_attribute, CollapseBasicStrategy::Mean);
371
372 SECTION("interior")
373 {
374 {
375 auto pos = mesh.create_accessor<double>(pos_attribute);
376 const Tuple v4 = mesh.tuple_from_id(PrimitiveType::Vertex, 4);
377 const Tuple v5 = mesh.tuple_from_id(PrimitiveType::Vertex, 5);
378 // reposition interior vertices
379 pos.vector_attribute(v4) = Eigen::Vector3d{1.4, 0, 0};
380 pos.vector_attribute(v5) = Eigen::Vector3d{1.6, 0, 0};
381 }
382
383 op.add_invariant(
384 std::make_shared<MaxEdgeLengthInvariant>(mesh, pos_attribute.as<double>(), 0.1));
385
386 Scheduler scheduler;
387
388 size_t n_vertices = mesh.get_all(PrimitiveType::Vertex).size();
389 size_t n_iterations = 0;
390 for (; n_iterations < 10; ++n_iterations) {
391 scheduler.run_operation_on_all(op);
392
393 const size_t n_vertices_new = mesh.get_all(PrimitiveType::Vertex).size();
394 if (n_vertices_new == n_vertices) {
395 break;
396 } else {
397 n_vertices = n_vertices_new;
398 }
399 }
400
401 REQUIRE(n_iterations == 1);
402 REQUIRE(n_vertices == 9);
403
404 // CHECK_THROWS(mesh.tuple_from_id(PrimitiveType::Vertex, 4));
405 const Tuple v5 = mesh.tuple_from_id(PrimitiveType::Vertex, 5);
406#if defined(WMTK_ENABLE_HASH_UPDATE)
407 REQUIRE(mesh.is_valid_with_hash(v5));
408#else
409 REQUIRE(mesh.is_valid(v5));
410#endif
411
412 auto pos = mesh.create_accessor<double>(pos_attribute);
413 Eigen::Vector3d p5 = pos.vector_attribute(v5);
414 CHECK((p5 - Eigen::Vector3d{1.5, 0, 0}).squaredNorm() == 0);
415 }
416 SECTION("towards_boundary_true")
417 {
418 {
419 auto pos = mesh.create_accessor<double>(pos_attribute);
420 const Tuple v4 = mesh.tuple_from_id(PrimitiveType::Vertex, 4);
421 // reposition vertex
422 pos.vector_attribute(v4) = Eigen::Vector3d{0.6, 0.9, 0};
423 }
424
425 // set collapse towards boundary
426 {
427 auto tmp = std::make_shared<CollapseNewAttributeStrategy<double>>(pos_attribute);
428 tmp->set_strategy(CollapseBasicStrategy::Mean);
429 tmp->set_simplex_predicate(BasicSimplexPredicate::IsInterior);
430 op.set_new_attribute_strategy(pos_attribute, tmp);
431 }
432
433 op.add_invariant(
434 std::make_shared<MaxEdgeLengthInvariant>(mesh, pos_attribute.as<double>(), 0.1));
435
436 Scheduler scheduler;
437
438 size_t n_vertices = mesh.get_all(PrimitiveType::Vertex).size();
439 size_t n_iterations = 0;
440 for (; n_iterations < 10; ++n_iterations) {
441 scheduler.run_operation_on_all(op);
442
443 const size_t n_vertices_new = mesh.get_all(PrimitiveType::Vertex).size();
444 if (n_vertices_new == n_vertices) {
445 break;
446 } else {
447 n_vertices = n_vertices_new;
448 }
449 }
450
451 REQUIRE(n_iterations == 1);
452 REQUIRE(n_vertices == 9);
453
454 const Tuple v0 = mesh.tuple_from_id(PrimitiveType::Vertex, 0);
455#if defined(WMTK_ENABLE_HASH_UPDATE)
456 REQUIRE(mesh.is_valid_with_hash(v0));
457#else
458 REQUIRE(mesh.is_valid(v0));
459#endif
460
461 auto pos = mesh.create_accessor<double>(pos_attribute);
462 Eigen::Vector3d p0 = pos.vector_attribute(v0);
463 CHECK((p0 - Eigen::Vector3d{0.5, 1, 0}).squaredNorm() == 0);
464 }
465 SECTION("towards_boundary_false")
466 {
467 {
468 auto pos = mesh.create_accessor<double>(pos_attribute);
469 const Tuple v4 = mesh.tuple_from_id(PrimitiveType::Vertex, 4);
470 // reposition vertex
471 pos.vector_attribute(v4) = Eigen::Vector3d{0.6, 0.9, 0};
472 }
473
474 op.add_invariant(
475 std::make_shared<MaxEdgeLengthInvariant>(mesh, pos_attribute.as<double>(), 0.1));
476 // op_settings.collapse_towards_boundary = false; <-- invariant missing
477
478
479 Scheduler scheduler;
480
481 size_t n_vertices = mesh.get_all(PrimitiveType::Vertex).size();
482 size_t n_iterations = 0;
483 for (; n_iterations < 10; ++n_iterations) {
484 scheduler.run_operation_on_all(op);
485
486 const size_t n_vertices_new = mesh.get_all(PrimitiveType::Vertex).size();
487 if (n_vertices_new == n_vertices) {
488 break;
489 } else {
490 n_vertices = n_vertices_new;
491 }
492 }
493
494 REQUIRE(n_iterations == 1);
495 REQUIRE(n_vertices == 9);
496
497 const Tuple v0 = mesh.tuple_from_id(PrimitiveType::Vertex, 0);
498#if defined(WMTK_ENABLE_HASH_UPDATE)
499 REQUIRE(mesh.is_valid_with_hash(v0));
500#else
501 REQUIRE(mesh.is_valid(v0));
502#endif
503
504 auto pos = mesh.create_accessor<double>(pos_attribute);
505 Eigen::Vector3d p0 = pos.vector_attribute(v0);
506 CHECK((p0 - Eigen::Vector3d{0.55, 0.95, 0}).squaredNorm() == 0);
507 }
508 SECTION("collapse_boundary_true")
509 {
510 {
511 auto pos = mesh.create_accessor<double>(pos_attribute);
512 const Tuple v1 = mesh.tuple_from_id(PrimitiveType::Vertex, 1);
513 // reposition vertex
514 pos.vector_attribute(v1) = Eigen::Vector3d{0.6, 1, 0};
515 }
516
517 op.add_invariant(
518 std::make_shared<MaxEdgeLengthInvariant>(mesh, pos_attribute.as<double>(), 0.1));
519
520 Scheduler scheduler;
521
522 scheduler.run_operation_on_all(op);
523
524 REQUIRE(mesh.get_all(PrimitiveType::Vertex).size() == 9);
525
526 const Tuple v0 = mesh.tuple_from_id(PrimitiveType::Vertex, 0);
527#if defined(WMTK_ENABLE_HASH_UPDATE)
528 REQUIRE(mesh.is_valid_with_hash(v0));
529#else
530 REQUIRE(mesh.is_valid(v0));
531#endif
532
533 auto pos = mesh.create_accessor<double>(pos_attribute);
534 Eigen::Vector3d p0 = pos.vector_attribute(v0);
535 CHECK((p0 - Eigen::Vector3d{0.55, 1, 0}).squaredNorm() == 0);
536 }
537 SECTION("collapse_boundary_false")
538 {
539 {
540 auto pos = mesh.create_accessor<double>(pos_attribute);
541 const Tuple v1 = mesh.tuple_from_id(PrimitiveType::Vertex, 1);
542 // reposition vertex
543 pos.vector_attribute(v1) = Eigen::Vector3d{0.6, 1, 0};
544 }
545
546 op.add_invariant(
547 std::make_shared<MaxEdgeLengthInvariant>(mesh, pos_attribute.as<double>(), 0.1));
548 op.add_invariant(
549 std::make_shared<invariants::InteriorSimplexInvariant>(mesh, PrimitiveType::Edge));
550
551 Scheduler scheduler;
552
553 scheduler.run_operation_on_all(op);
554
555 REQUIRE(mesh.get_all(PrimitiveType::Vertex).size() == 10);
556 }
557}
558
559TEST_CASE("swap_edge_for_valence", "[components][isotropic_remeshing][swap][2D]")
560{
561 using namespace operations;
562
563 DEBUG_TriMesh mesh = wmtk::tests::embedded_diamond();
564
565 composite::TriEdgeSwap op(mesh);
566 op.add_invariant(
567 std::make_shared<invariants::InteriorSimplexInvariant>(mesh, PrimitiveType::Edge));
568 op.collapse().add_invariant(std::make_shared<MultiMeshLinkConditionInvariant>(mesh));
569
570 Tuple swap_edge = mesh.edge_tuple_with_vs_and_t(6, 7, 5);
571
572 auto vertex_one_ring = [](TriMesh& m, const Tuple& t) {
574 .simplex_vector(PrimitiveType::Vertex);
575 };
576
577 SECTION("single_op_fail")
578 {
579 op.add_invariant(std::make_shared<invariants::ValenceImprovementInvariant>(mesh));
580 CHECK(op(Simplex::edge(mesh, swap_edge)).empty());
581 }
582 SECTION("swap_success")
583 {
584 // swap edge to create inbalence in valence
585 {
586 const auto ret = op(Simplex::edge(mesh, swap_edge));
587 REQUIRE(!ret.empty());
588 const Tuple& return_tuple = ret[0].tuple();
589 swap_edge = return_tuple;
590 long id0 = mesh.id_vertex(swap_edge);
591 long id1 = mesh.id_vertex(mesh.switch_vertex(swap_edge));
592
593 // check valence
594 const Tuple v3 = mesh.tuple_from_id(PrimitiveType::Vertex, 3);
595 const Tuple v6 = mesh.tuple_from_id(PrimitiveType::Vertex, 6);
596 const Tuple v7 = mesh.tuple_from_id(PrimitiveType::Vertex, 7);
597 const Tuple v10 = mesh.tuple_from_id(PrimitiveType::Vertex, 10);
598 CHECK(vertex_one_ring(mesh, v3).size() == 7);
599 CHECK(vertex_one_ring(mesh, v10).size() == 7);
600 CHECK(vertex_one_ring(mesh, v6).size() == 5);
601 CHECK(vertex_one_ring(mesh, v7).size() == 5);
602 }
603
604 op.add_invariant(std::make_shared<invariants::ValenceImprovementInvariant>(mesh));
605
606
607 SECTION("single_op")
608 {
609 const auto ret = op(Simplex::edge(mesh, swap_edge));
610 REQUIRE(!ret.empty());
611 const Tuple& return_tuple = ret[0].tuple();
612 swap_edge = return_tuple;
613 CHECK(mesh.id(Simplex::vertex(mesh, swap_edge)) == 7);
614 CHECK(mesh.id(Simplex::vertex(mesh, mesh.switch_vertex(swap_edge))) == 6);
615 }
616 SECTION("with_scheduler")
617 {
618 Scheduler scheduler;
619 scheduler.run_operation_on_all(op);
620 }
621
622
623 // check valence
624 {
625 const Tuple v3 = mesh.tuple_from_id(PrimitiveType::Vertex, 3);
626 const Tuple v6 = mesh.tuple_from_id(PrimitiveType::Vertex, 6);
627 const Tuple v7 = mesh.tuple_from_id(PrimitiveType::Vertex, 7);
628 const Tuple v10 = mesh.tuple_from_id(PrimitiveType::Vertex, 10);
629 CHECK(vertex_one_ring(mesh, v3).size() == 6);
630 CHECK(vertex_one_ring(mesh, v10).size() == 6);
631 CHECK(vertex_one_ring(mesh, v6).size() == 6);
632 CHECK(vertex_one_ring(mesh, v7).size() == 6);
633 }
634 }
635 SECTION("swap_fail")
636 {
637 const Tuple e = mesh.edge_tuple_with_vs_and_t(6, 7, 5);
638 op.add_invariant(std::make_shared<invariants::ValenceImprovementInvariant>(mesh));
639 const bool success = !op(Simplex::edge(mesh, e)).empty();
640 CHECK(!success);
641 }
642}
643
644TEST_CASE("component_isotropic_remeshing", "[components][isotropic_remeshing][2D]")
645{
646 io::Cache cache("wmtk_cache", ".");
647
648 if constexpr (true) {
649 const std::filesystem::path input_file = data_dir / "small.msh";
650 json input_component_json = {
651 {"name", "input_mesh"},
652 {"file", input_file.string()},
653 {"ignore_z", false},
654 {"tetrahedron_attributes", json::array()}};
655 REQUIRE_NOTHROW(wmtk::components::input(Paths(), input_component_json, cache));
656 }
657
658 if constexpr (true) {
659 json mesh_isotropic_remeshing_json = R"({
660 "input": "input_mesh",
661 "output": "output_mesh",
662 "attributes": {"position": "vertices", "inversion_position": [], "other_positions": [], "update_other_positions": false},
663 "pass_through": [],
664 "iterations": 1,
665 "length_abs": 0.003,
666 "length_rel": -1,
667 "lock_boundary": true,
668 "use_for_periodic": false,
669 "dont_disable_split": false
670 })"_json;
671 REQUIRE_NOTHROW(
672 wmtk::components::isotropic_remeshing(Paths(), mesh_isotropic_remeshing_json, cache));
673 }
674 if constexpr (true) {
675 json component_json = R"({
676 "input": "output_mesh",
677 "attributes": {"position": "vertices"},
678 "file": "bunny_isotropic_remeshing"
679 })"_json;
680
681 CHECK_NOTHROW(wmtk::components::output(Paths(), component_json, cache));
682
683 // auto mesh_in = cache.read_mesh("output_mesh");
684 // TriMesh& m = static_cast<TriMesh&>(*mesh_in);
685 // ParaviewWriter writer(
686 // cache.get_cache_path() / "isotropic_remeshing_output",
687 // "vertices",
688 // m,
689 // false,
690 // false,
691 // true,
692 // false);
693 // m.serialize(writer);
694 }
695}
696
697TEST_CASE("remeshing_tetrahedron", "[components][isotropic_remeshing][2D]")
698{
699 using namespace wmtk::components::internal;
700
701 io::Cache cache("wmtk_cache", ".");
702
703 // input
704 DEBUG_TriMesh mesh = tetrahedron_with_position();
705
706 auto pos_handle = mesh.get_attribute_handle<double>("vertices", PrimitiveType::Vertex);
707 std::vector<attribute::MeshAttributeHandle> pass_through_attributes;
708 pass_through_attributes.emplace_back(
709 mesh.register_attribute<int64_t>("dummy_tag", PrimitiveType::Triangle, 1));
710
711 {
712 auto acc = mesh.create_accessor<int64_t>(pass_through_attributes[0]);
713 acc.scalar_attribute(mesh.face_tuple_from_vids(0, 1, 2)) = 1;
714 acc.scalar_attribute(mesh.face_tuple_from_vids(0, 1, 3)) = 1;
715 }
716
717
718 CHECK_NOTHROW(wmtk::components::internal::isotropic_remeshing(
719 pos_handle,
720 pass_through_attributes,
721 0.1,
722 true,
723 false,
724 false,
725 false,
726 10));
727
728 {
729 ParaviewWriter writer(
730 cache.get_cache_path() / "tet_remeshing",
731 "vertices",
732 mesh,
733 false,
734 false,
735 true,
736 false);
737 mesh.serialize(writer);
738 }
739}
740
741TEST_CASE("remeshing_with_boundary", "[components][isotropic_remeshing][2D]")
742{
743 using namespace wmtk::components::internal;
744
745 io::Cache cache("wmtk_cache", ".");
746
747 // input
748 TriMesh mesh = edge_region_with_position();
749 auto pos_handle = mesh.get_attribute_handle<double>("vertices", PrimitiveType::Vertex);
750 std::vector<attribute::MeshAttributeHandle> pass_through_attributes;
751
752 SECTION("lock_boundary_false")
753 {
754 wmtk::components::internal::isotropic_remeshing(
755 pos_handle,
756 pass_through_attributes,
757 0.5,
758 false,
759 false,
760 false,
761 false,
762 5);
763
764 size_t n_boundary_edges = 0;
765 for (const Tuple& e : mesh.get_all(PrimitiveType::Edge)) {
766 if (mesh.is_boundary_edge(e)) {
767 ++n_boundary_edges;
768 }
769 }
770 CHECK(n_boundary_edges > 8);
771
772 {
773 ParaviewWriter writer(
774 cache.get_cache_path() / "w_bd_remeshing_lock_false",
775 "vertices",
776 mesh,
777 false,
778 false,
779 true,
780 false);
781 mesh.serialize(writer);
782 }
783 }
784
785 SECTION("lock_boundary_true")
786 {
787 wmtk::components::internal::isotropic_remeshing(
788 pos_handle,
789 pass_through_attributes,
790 0.5,
791 true,
792 false,
793 false,
794 false,
795 5);
796
797
798 size_t n_boundary_edges = 0;
799 for (const Tuple& e : mesh.get_all(PrimitiveType::Edge)) {
800 if (mesh.is_boundary_edge(e)) {
801 ++n_boundary_edges;
802 }
803 }
804 CHECK(n_boundary_edges == 8);
805
806 {
807 ParaviewWriter writer(
808 cache.get_cache_path() / "w_bd_remeshing_lock_true",
809 "vertices",
810 mesh,
811 false,
812 false,
813 true,
814 false);
815 mesh.serialize(writer);
816 }
817 }
818}
819
820TEST_CASE("remeshing_preserve_topology", "[components][isotropic_remeshing][2D][.]")
821{
822 using namespace wmtk::components::internal;
823
824 // input
825 DEBUG_TriMesh mesh = edge_region_with_position();
826 // DEBUG_TriMesh mesh = hex_plus_two_with_position();
827
828 auto pos_handle = mesh.get_attribute_handle<double>("vertices", PrimitiveType::Vertex);
829 std::vector<attribute::MeshAttributeHandle> pass_through_attributes;
830
831 auto tag_handle = mesh.register_attribute<int64_t>("is_boundary", wmtk::PrimitiveType::Edge, 1);
832 auto tag_accessor = mesh.create_accessor<int64_t>(tag_handle);
833 for (const Tuple& e : mesh.get_all(PrimitiveType::Edge)) {
834 if (mesh.is_boundary_edge(e)) {
835 tag_accessor.scalar_attribute(e) = 1;
836 } else {
837 tag_accessor.scalar_attribute(e) = 0;
838 }
839 }
840 std::shared_ptr<Mesh> child_ptr =
842 mesh,
843 "is_boundary",
844 1,
845 PrimitiveType::Edge);
846
847 REQUIRE(mesh.get_child_meshes().size() == 1);
848 mesh.multi_mesh_manager().check_map_valid(mesh);
849 const auto& child_mesh = *child_ptr;
850 CHECK(child_mesh.get_all(PrimitiveType::Edge).size() == 8);
851 CHECK(child_mesh.get_all(PrimitiveType::Vertex).size() == 8);
852
853
854 wmtk::components::internal::isotropic_remeshing(
855 pos_handle,
856 pass_through_attributes,
857 0.5,
858 /*lock_boundary*/ false,
859 true,
860 true,
861 false,
862 5);
863
864 REQUIRE(mesh.is_connectivity_valid());
865 mesh.multi_mesh_manager().check_map_valid(mesh);
866
867
868 size_t n_boundary_edges = 0;
869 for (const Tuple& e : mesh.get_all(PrimitiveType::Edge)) {
870 if (mesh.is_boundary_edge(e)) {
871 ++n_boundary_edges;
872 }
873 }
874 // CHECK(n_boundary_edges > 8);
875
876 // output
877 {
878 ParaviewWriter writer("remeshing_test", "vertices", mesh, true, true, true, false);
879 mesh.serialize(writer);
880 }
881}
882
883TEST_CASE("remeshing_preserve_topology_realmesh", "[components][isotropic_remeshing][2D][.]")
884{
885 using namespace wmtk::components::internal;
886 using namespace operations;
887
888 wmtk::io::Cache cache("wmtk_cache", ".");
889
890 // input
891 // TODO: What is the default attribute for "vertices". From the reader it seems to be
892 // "vertices". need change "vertices" to "vertices" isotropic_remeshing.hpp
893 json input_component_json = {
894 {"type", "input"},
895 {"name", "input_mesh"},
896 {"cell_dimension", 2},
897 {"file", (data_dir / "circle.msh").string()},
898 {"ignore_z", false}};
899 wmtk::components::input(Paths(), input_component_json, cache);
900
901 auto m = cache.read_mesh(input_component_json["name"]);
902 tests::DEBUG_TriMesh& mesh = static_cast<tests::DEBUG_TriMesh&>(*m);
903
904 auto pos_handle = mesh.get_attribute_handle<double>("vertices", PrimitiveType::Vertex);
905 std::vector<attribute::MeshAttributeHandle> pass_through_attributes;
906
907 auto tag_handle = mesh.register_attribute<int64_t>("is_boundary", wmtk::PrimitiveType::Edge, 1);
908 auto tag_accessor = mesh.create_accessor<int64_t>(tag_handle);
909 for (const Tuple& e : mesh.get_all(PrimitiveType::Edge)) {
910 if (mesh.is_boundary_edge(e)) {
911 tag_accessor.scalar_attribute(e) = 1;
912 } else {
913 tag_accessor.scalar_attribute(e) = 0;
914 }
915 }
916 std::shared_ptr<Mesh> child_ptr =
918 mesh,
919 "is_boundary",
920 1,
921 PrimitiveType::Edge);
922
923 REQUIRE(mesh.get_child_meshes().size() == 1);
924 // mesh.multi_mesh_manager().check_map_valid(mesh);
925 // const auto& child_mesh = *child_ptr;
926
927
928 // IsotropicRemeshing isotropicRemeshing(mesh, 0.5, false, false, false);
929
930 for (int i = 0; i < 25; i++) {
931 wmtk::components::internal::isotropic_remeshing(
932 pos_handle,
933 pass_through_attributes,
934 0.05,
935 false,
936 true,
937 false,
938 false,
939 1);
940 REQUIRE(mesh.is_connectivity_valid());
941 mesh.multi_mesh_manager().check_map_valid(mesh);
942 }
943
944
945 auto child_vertex_handle =
946 child_ptr->register_attribute<double>("vertices", wmtk::PrimitiveType::Vertex, 3);
947 auto child_vertex_accessor = child_ptr->create_accessor<int64_t>(child_vertex_handle);
948
949 auto parent_vertex_handle =
950 mesh.get_attribute_handle<double>("vertices", PrimitiveType::Vertex);
951 auto parent_vertex_accessor = mesh.create_accessor<int64_t>(parent_vertex_handle);
952
953 for (const auto& v : child_ptr->get_all(PrimitiveType::Vertex)) {
954 auto parent_v = child_ptr->map_to_root_tuple(Simplex(*child_ptr, PrimitiveType::Vertex, v));
955 child_vertex_accessor.vector_attribute(v) =
956 parent_vertex_accessor.vector_attribute(parent_v);
957 }
958
959 // output
960 {
961 ParaviewWriter writer(
962 cache.get_cache_path() / "remeshing_test_circle_final",
963 "vertices",
964 mesh,
965 true,
966 true,
967 true,
968 false);
969 mesh.serialize(writer);
970
971 ParaviewWriter writer2(
972 cache.get_cache_path() / "remeshing_test_circle_child_mesh_final",
973 "vertices",
974 *child_ptr,
975 true,
976 true,
977 false,
978 false);
979 child_ptr->serialize(writer2);
980 }
981}
982
983TEST_CASE("remeshing_realmesh", "[components][isotropic_remeshing][2D][.]")
984{
985 using namespace wmtk::components::internal;
986 using namespace operations;
987
988 wmtk::io::Cache cache("wmtk_cache", ".");
989
990 // input
991 json input_component_json = {
992 {"type", "input"},
993 {"name", "input_mesh"},
994 {"cell_dimension", 2},
995 {"file", (data_dir / "circle.msh").string()},
996 {"ignore_z", false}};
997 wmtk::components::input(Paths(), input_component_json, cache);
998
999 auto m = cache.read_mesh(input_component_json["name"]);
1000 TriMesh& mesh = static_cast<TriMesh&>(*m);
1001
1002 auto pos_handle = mesh.get_attribute_handle<double>("vertices", PrimitiveType::Vertex);
1003 std::vector<attribute::MeshAttributeHandle> pass_through_attributes;
1004
1005 // auto tag_handle = mesh.register_attribute<int64_t>("is_boundary", wmtk::PrimitiveType::Edge,
1006 // 1); auto tag_accessor = mesh.create_accessor<int64_t>(tag_handle); for (const Tuple& e :
1007 // mesh.get_all(PrimitiveType::Edge)) {
1008 // if (mesh.is_boundary_edge(e)) {
1009 // tag_accessor.scalar_attribute(e) = 1;
1010 // } else {
1011 // tag_accessor.scalar_attribute(e) = 0;
1012 // }
1013 // }
1014 // std::shared_ptr<Mesh> child_ptr =
1015 // wmtk::multimesh::utils::extract_and_register_child_mesh_from_tag(
1016 // mesh,
1017 // "is_boundary",
1018 // 1,
1019 // PrimitiveType::Edge);
1020
1021 // REQUIRE(mesh.get_child_meshes().size() == 1);
1022 // mesh.multi_mesh_manager().check_map_valid(mesh);
1023 // const auto& child_mesh = *child_ptr;
1024
1025 wmtk::components::internal::isotropic_remeshing(
1026 pos_handle,
1027 pass_through_attributes,
1028 0.5,
1029 false,
1030 true,
1031 false,
1032 false,
1033 25);
1034
1035 REQUIRE(mesh.is_connectivity_valid());
1036 // mesh.multi_mesh_manager().check_map_valid(mesh);
1037
1038 // output
1039 {
1040 ParaviewWriter writer(
1041 cache.get_cache_path() / "remeshing_test_circle_no_nultimesh",
1042 "vertices",
1043 mesh,
1044 true,
1045 true,
1046 true,
1047 false);
1048 mesh.serialize(writer);
1049 }
1050}
void serialize(MeshWriter &writer, const Mesh *local_root=nullptr) const
Definition Mesh.cpp:93
attribute::MeshAttributeHandle get_attribute_handle(const std::string &name, const PrimitiveType ptype) const
Definition Mesh.hpp:903
std::vector< Tuple > get_all(PrimitiveType type) const
Generate a vector of Tuples from global vertex/edge/triangle/tetrahedron index.
Definition Mesh.cpp:18
SchedulerStats run_operation_on_all(operations::Operation &op)
Definition Scheduler.cpp:33
bool is_connectivity_valid() const final override
Definition TriMesh.cpp:359
bool is_boundary_edge(const Tuple &tuple) const
Definition TriMesh.cpp:72
The Tuple is the basic navigation tool in our mesh data structure.
Definition Tuple.hpp:19
auto as() const -> const held_handle_type< held_type_from_primitive< T >()> &
std::filesystem::path get_cache_path() const
Get the path of the cache folder.
Definition Cache.cpp:146
std::shared_ptr< Mesh > read_mesh(const std::string &name) const
Load a mesh from cache.
Definition Cache.cpp:171
void set_function(const UpdateFunction &func)
void add_invariant(std::shared_ptr< Invariant > invariant)
Definition Operation.hpp:48
Performs an edge swap, implemented as a combination of swap and collapse.
const std::vector< Simplex > & simplex_vector() const
Return const reference to the simplex vector.
static Simplex edge(const Mesh &m, const Tuple &t)
Definition Simplex.hpp:61
static Simplex vertex(const Mesh &m, const Tuple &t)
Definition Simplex.hpp:56
std::shared_ptr< Mesh > extract_and_register_child_mesh_from_tag(Mesh &m, const std::string &tag, const int64_t tag_value, const PrimitiveType pt, bool child_is_free)
extract a child mesh based on the given tag and tag value, and register it to the input (parent) mesh
std::tuple< Tuple, Tuple > vectors_to_tuples(const Eigen::Ref< const TwoTupleVector > &v)
SimplexCollection link(const Mesh &mesh, const simplex::Simplex &simplex, const bool sort_and_clean)
Definition link.cpp:84
void print_tuple_map_iso(const DEBUG_TriMesh &parent, const DEBUG_MultiMeshManager &p_mul_manager)
const std::filesystem::path data_dir
TEST_CASE("smoothing_mesh", "[components][isotropic_remeshing][2D]")
nlohmann::json json
Definition input.cpp:9