Wildmeshing Toolkit
Loading...
Searching...
No Matches
wildmeshing_old.cpp
Go to the documentation of this file.
1#include "wildmeshing.hpp"
2
4
5#include <wmtk/Mesh.hpp>
6#include <wmtk/Scheduler.hpp>
7#include <wmtk/TetMesh.hpp>
8#include <wmtk/TriMesh.hpp>
10
13#include <wmtk/utils/Logger.hpp>
14
15
19
33
34
38
61
63
65#include <wmtk/utils/orient.hpp>
66
69
70#include <queue>
72#include <wmtk/simplex/link.hpp>
73
74#include <fstream>
75namespace wmtk::components {
76
77using namespace simplex;
78using namespace operations;
79using namespace operations::tri_mesh;
80using namespace operations::tet_mesh;
81using namespace operations::composite;
82using namespace function;
83using namespace invariants;
84
85namespace {
86void write(
87 const std::shared_ptr<Mesh>& mesh,
88 const std::string& out_dir,
89 const std::string& name,
90 const std::string& vname,
91 const int64_t index,
92 const bool intermediate_output)
93{
94 if (intermediate_output) {
95 if (mesh->top_simplex_type() == PrimitiveType::Triangle) {
96 // write trimesh
97 const std::filesystem::path data_dir = "";
99 data_dir / (name + "_" + std::to_string(index)),
100 vname,
101 *mesh,
102 true,
103 true,
104 true,
105 false);
106 mesh->serialize(writer);
107 } else if (mesh->top_simplex_type() == PrimitiveType::Tetrahedron) {
108 // write tetmesh
109 const std::filesystem::path data_dir = "";
111 data_dir / (name + "_" + std::to_string(index)),
112 vname,
113 *mesh,
114 true,
115 true,
116 true,
117 true);
118 mesh->serialize(writer);
119 } else if (mesh->top_simplex_type() == PrimitiveType::Edge) {
120 // write edgemesh
121 const std::filesystem::path data_dir = "";
123 data_dir / (name + "_" + std::to_string(index)),
124 vname,
125 *mesh,
126 true,
127 true,
128 false,
129 false);
130 mesh->serialize(writer);
131 }
132 }
133}
134
135} // namespace
136
138 Mesh& m,
139 const TypedAttributeHandle<Rational>& coordinate_handle,
140 const TypedAttributeHandle<double>& edge_length_handle,
141 const TypedAttributeHandle<double>& sizing_field_scalar_handle,
142 const TypedAttributeHandle<double>& energy_handle,
143 const TypedAttributeHandle<double>& target_edge_length_handle,
144 const TypedAttributeHandle<char>& visited_handle,
145 const double stop_energy,
146 const double current_max_energy,
147 const double initial_target_edge_length,
148 const double min_target_edge_length)
149{
151
152 std::cout << "in here" << std::endl;
153
154 const auto coordinate_accessor = m.create_const_accessor<Rational>(coordinate_handle);
155 const auto edge_length_accessor = m.create_const_accessor<double>(edge_length_handle);
156 const auto energy_accessor = m.create_const_accessor<double>(energy_handle);
157
158 auto sizing_field_scalar_accessor = m.create_accessor<double>(sizing_field_scalar_handle);
159 auto target_edge_length_accessor = m.create_accessor<double>(target_edge_length_handle);
160 auto visited_accessor = m.create_accessor<char>(visited_handle);
161
162 const double stop_filter_energy = stop_energy * 0.8;
163 double filter_energy = std::max(current_max_energy / 100, stop_filter_energy);
164 filter_energy = std::min(filter_energy, 100.0);
165
166 const double recover_scalar = 1.5;
167 const double refine_scalar = 0.5;
168 const double min_refine_scalar = min_target_edge_length / initial_target_edge_length;
169
170 auto vertices_all = m.get_all(PrimitiveType::Vertex);
171 int64_t v_cnt = vertices_all.size();
172
173 std::vector<Vector3d> centroids;
174 std::queue<Tuple> v_queue;
175 // std::vector<double> scale_multipliers(v_cnt, recover_scalar);
176
177 // get centroids and initial v_queue
178 for (const auto& t : m.get_all(PrimitiveType::Tetrahedron)) {
179 // if (std::cbrt(energy_accessor.const_scalar_attribute(t)) < filter_energy) {
180 if (energy_accessor.const_scalar_attribute(t) < filter_energy) {
181 // skip good tets
182 continue;
183 }
184 auto vertices = m.orient_vertices(t);
185 Vector3d c(0, 0, 0);
186 for (int i = 0; i < 4; ++i) {
187 c += coordinate_accessor.const_vector_attribute(vertices[i]).cast<double>();
188 v_queue.emplace(vertices[i]);
189 }
190 centroids.emplace_back(c / 4.0);
191 }
192
193 wmtk::logger().info(
194 "filter energy: {}, low quality tets num: {}",
195 filter_energy,
196 centroids.size());
197
198 const double R = initial_target_edge_length * 1.8; // update field radius
199
200 for (const auto& v : vertices_all) { // reset visited flag
201 visited_accessor.scalar_attribute(v) = char(0);
202 }
203
204 // TODO: use efficient data structure
205 auto get_nearest_dist = [&](const Tuple& v) -> double {
206 Tuple nearest_tuple;
207 double min_dist = std::numeric_limits<double>::max();
208 const Vector3d v_pos = coordinate_accessor.const_vector_attribute(v).cast<double>();
209 for (const auto& c_pos : centroids) {
210 double dist = (c_pos - v_pos).norm();
211 min_dist = std::min(min_dist, dist);
212 }
213 return min_dist;
214 };
215
216 while (!v_queue.empty()) {
217 auto v = v_queue.front();
218 v_queue.pop();
219
220 if (visited_accessor.scalar_attribute(v) == char(1)) continue;
221 visited_accessor.scalar_attribute(v) = char(1);
222
223 double dist = std::max(0., get_nearest_dist(v));
224
225 if (dist > R) {
226 visited_accessor.scalar_attribute(v) = char(0);
227 continue;
228 }
229
230 double scale_multiplier = std::min(
231 recover_scalar,
232 dist / R * (1 - refine_scalar) + refine_scalar); // linear interpolate
233
234 auto new_scale = sizing_field_scalar_accessor.scalar_attribute(v) * scale_multiplier;
235 if (new_scale > 1) {
236 sizing_field_scalar_accessor.scalar_attribute(v) = 1;
237 } else if (new_scale < min_refine_scalar) {
238 sizing_field_scalar_accessor.scalar_attribute(v) = min_refine_scalar;
239 } else {
240 sizing_field_scalar_accessor.scalar_attribute(v) = new_scale;
241 }
242
243 // push one ring vertices into the queue
244 for (const auto& v_one_ring : simplex::link(m, simplex::Simplex::vertex(m, v))
246 if (visited_accessor.scalar_attribute(v_one_ring) == char(1)) continue;
247 v_queue.push(v_one_ring.tuple());
248 }
249 }
250
251 // update the rest
252 for (const auto& v : vertices_all) {
253 if (visited_accessor.scalar_attribute(v) == char(1)) continue;
254 auto new_scale = sizing_field_scalar_accessor.scalar_attribute(v) * 1.5;
255 if (new_scale > 1) {
256 sizing_field_scalar_accessor.scalar_attribute(v) = 1;
257 } else if (new_scale < min_refine_scalar) {
258 sizing_field_scalar_accessor.scalar_attribute(v) = min_refine_scalar;
259 } else {
260 sizing_field_scalar_accessor.scalar_attribute(v) = new_scale;
261 }
262 }
263
264 // update target edge length
265 for (const auto& e : m.get_all(PrimitiveType::Edge)) {
266 target_edge_length_accessor.scalar_attribute(e) =
267 initial_target_edge_length *
268 (sizing_field_scalar_accessor.scalar_attribute(e) +
269 sizing_field_scalar_accessor.scalar_attribute(
271 2;
272 }
273}
274
276 Mesh& m,
277 const TypedAttributeHandle<Rational>& coordinate_handle,
278 const TypedAttributeHandle<double>& energy_handle,
279 const TypedAttributeHandle<char>& energy_filter_handle,
280 const TypedAttributeHandle<char>& visited_handle,
281 const double stop_energy,
282 const double current_max_energy,
283 const double initial_target_edge_length)
284{
285 // two ring version
287
288 const auto coordinate_accessor = m.create_const_accessor<Rational>(coordinate_handle);
289 const auto energy_accessor = m.create_const_accessor<double>(energy_handle);
290
291 auto energy_filter_accessor = m.create_accessor<char>(energy_filter_handle);
292 auto visited_accessor = m.create_accessor<char>(visited_handle);
293
294 const double stop_filter_energy = stop_energy * 0.8;
295 double filter_energy = std::max(current_max_energy / 100, stop_filter_energy);
296 filter_energy = std::min(filter_energy, 100.0);
297
298 auto vertices_all = m.get_all(PrimitiveType::Vertex);
299
300 for (const auto& v : vertices_all) { // reset visited flag
301 visited_accessor.scalar_attribute(v) = char(0);
302 energy_filter_accessor.scalar_attribute(v) = char(0);
303 }
304
305 // get centroids and initial v_queue
306 for (const auto& t : m.get_all(PrimitiveType::Tetrahedron)) {
307 if (energy_accessor.const_scalar_attribute(t) < filter_energy) {
308 // skip good tets
309 continue;
310 }
311 auto vertices = m.orient_vertices(t);
312 for (const auto& v : vertices) {
313 energy_filter_accessor.scalar_attribute(v) = char(1);
314 for (const auto& vv : simplex::k_ring(m, simplex::Simplex::vertex(m, v), 1)
316 energy_filter_accessor.scalar_attribute(vv) = char(1);
317 }
318 }
319 }
320}
321
323 Mesh& m,
324 const TypedAttributeHandle<Rational>& coordinate_handle,
325 const TypedAttributeHandle<double>& energy_handle,
326 const TypedAttributeHandle<char>& energy_filter_handle,
327 const TypedAttributeHandle<char>& visited_handle,
328 const double stop_energy,
329 const double current_max_energy,
330 const double initial_target_edge_length)
331{
333
334 const auto coordinate_accessor = m.create_const_accessor<Rational>(coordinate_handle);
335 const auto energy_accessor = m.create_const_accessor<double>(energy_handle);
336
337 auto energy_filter_accessor = m.create_accessor<char>(energy_filter_handle);
338 auto visited_accessor = m.create_accessor<char>(visited_handle);
339
340 const double stop_filter_energy = stop_energy * 0.8;
341 double filter_energy = std::max(current_max_energy / 100, stop_filter_energy);
342 filter_energy = std::min(filter_energy, 100.0);
343
344 auto vertices_all = m.get_all(PrimitiveType::Vertex);
345 std::vector<Vector3d> centroids;
346 std::queue<Tuple> v_queue;
347
348 // get centroids and initial v_queue
349 for (const auto& t : m.get_all(PrimitiveType::Tetrahedron)) {
350 if (energy_accessor.const_scalar_attribute(t) < filter_energy) {
351 // skip good tets
352 continue;
353 }
354 auto vertices = m.orient_vertices(t);
355 Vector3d c(0, 0, 0);
356 for (int i = 0; i < 4; ++i) {
357 c += coordinate_accessor.const_vector_attribute(vertices[i]).cast<double>();
358 v_queue.emplace(vertices[i]);
359 }
360 centroids.emplace_back(c / 4.0);
361 }
362
363 const double R = initial_target_edge_length * 1.8;
364
365 for (const auto& v : vertices_all) { // reset visited flag
366 visited_accessor.scalar_attribute(v) = char(0);
367 energy_filter_accessor.scalar_attribute(v) = char(0);
368 }
369
370 // TODO: use efficient data structure
371 auto get_nearest_dist = [&](const Tuple& v) -> double {
372 Tuple nearest_tuple;
373 double min_dist = std::numeric_limits<double>::max();
374 const Vector3d v_pos = coordinate_accessor.const_vector_attribute(v).cast<double>();
375 for (const auto& c_pos : centroids) {
376 double dist = (c_pos - v_pos).norm();
377 min_dist = std::min(min_dist, dist);
378 }
379 return min_dist;
380 };
381
382 while (!v_queue.empty()) {
383 auto v = v_queue.front();
384 v_queue.pop();
385
386 if (visited_accessor.scalar_attribute(v) == char(1)) continue;
387 visited_accessor.scalar_attribute(v) = char(1);
388
389 double dist = std::max(0., get_nearest_dist(v));
390
391 if (dist > R) {
392 visited_accessor.scalar_attribute(v) = char(0);
393 continue;
394 }
395
396 energy_filter_accessor.scalar_attribute(v) = char(1);
397
398 // push one ring vertices into the queue
399 for (const auto& v_one_ring : simplex::link(m, simplex::Simplex::vertex(m, v))
401 if (visited_accessor.scalar_attribute(v_one_ring) == char(1)) continue;
402 v_queue.push(v_one_ring.tuple());
403 }
404 }
405}
406
407void wildmeshing(const utils::Paths& paths, const nlohmann::json& j, io::Cache& cache)
408{
409 // opt_logger().set_level(spdlog::level::info);
410
412 // Load mesh from settings
413 WildmeshingOptions options = j.get<WildmeshingOptions>();
414 const std::filesystem::path& file = options.input;
415 auto mesh = cache.read_mesh(options.input);
416
417 if (!mesh->is_connectivity_valid()) {
418 throw std::runtime_error("input mesh for wildmeshing connectivity invalid");
419 }
420
421 wmtk::logger().trace("Getting rational point handle");
422
424 // Retriving vertices
425 //
427 wmtk::logger().trace("Found double attribugte");
428 auto pt_double_attribute =
429 mesh->get_attribute_handle<double>(options.attributes.position, PrimitiveType::Vertex);
430
431 if (!mesh->has_attribute<Rational>(options.attributes.position, PrimitiveType::Vertex)) {
432 wmtk::utils::cast_attribute<wmtk::Rational>(
433 pt_double_attribute,
434 *mesh,
435 options.attributes.position);
436
437
438 } else {
439 auto pt_attribute = mesh->get_attribute_handle<Rational>(
440 options.attributes.position,
442 wmtk::utils::cast_attribute<wmtk::Rational>(pt_double_attribute, pt_attribute);
443 }
444 mesh->delete_attribute(pt_double_attribute);
445 }
446 auto pt_attribute =
447 mesh->get_attribute_handle<Rational>(options.attributes.position, PrimitiveType::Vertex);
448 wmtk::logger().trace("Getting rational point accessor");
449 auto pt_accessor = mesh->create_accessor(pt_attribute.as<Rational>());
450
451 wmtk::logger().trace("Computing bounding box diagonal");
453 // computing bbox diagonal
454 Eigen::VectorXd bmin(mesh->top_cell_dimension());
455 bmin.setConstant(std::numeric_limits<double>::max());
456 Eigen::VectorXd bmax(mesh->top_cell_dimension());
457 bmax.setConstant(std::numeric_limits<double>::lowest());
458
459 const auto vertices = mesh->get_all(PrimitiveType::Vertex);
460 for (const auto& v : vertices) {
461 const auto p = pt_accessor.vector_attribute(v).cast<double>();
462 for (int64_t d = 0; d < bmax.size(); ++d) {
463 bmin[d] = std::min(bmin[d], p[d]);
464 bmax[d] = std::max(bmax[d], p[d]);
465 }
466 }
467
468
469 const double bbdiag = (bmax - bmin).norm();
470
471 wmtk::logger().info("bbox max {}, bbox min {}, diag {}", bmax, bmin, bbdiag);
472
473 const double target_edge_length = options.target_edge_length * bbdiag;
474
475 // const double target_edge_length =
476 // options.target_edge_length * bbdiag /
477 // std::min(
478 // 2.0,
479 // 1 + options.target_edge_length * 2); // min to prevent bad option.target_edge_length
480
481 wmtk::logger().info("target edge length: {}", target_edge_length);
482
484 // store amips
485 auto amips_attribute =
486 mesh->register_attribute<double>("wildmeshing_amips", mesh->top_simplex_type(), 1);
487 auto amips_accessor = mesh->create_accessor(amips_attribute.as<double>());
488 // amips update
489 auto compute_amips = [](const Eigen::MatrixX<Rational>& P) -> Eigen::VectorXd {
490 assert(P.rows() == 2 || P.rows() == 3); // rows --> attribute dimension
491 assert(P.cols() == P.rows() + 1);
492 if (P.cols() == 3) {
493 // triangle
494 assert(P.rows() == 2);
495 std::array<double, 6> pts;
496 for (size_t i = 0; i < 3; ++i) {
497 for (size_t j = 0; j < 2; ++j) {
498 pts[2 * i + j] = P(j, i).to_double();
499 }
500 }
501 const double a = Tri_AMIPS_energy(pts);
502 return Eigen::VectorXd::Constant(1, a);
503 } else {
504 // tet
505 assert(P.rows() == 3);
506 std::array<double, 12> pts;
507 for (size_t i = 0; i < 4; ++i) {
508 for (size_t j = 0; j < 3; ++j) {
509 pts[3 * i + j] = P(j, i).to_double();
510 }
511 }
512 const double a = Tet_AMIPS_energy(pts);
513 return Eigen::VectorXd::Constant(1, a);
514 }
515 };
516 auto amips_update =
517 std::make_shared<wmtk::operations::SingleAttributeTransferStrategy<double, Rational>>(
518 amips_attribute,
519 pt_attribute,
520 compute_amips);
521 amips_update->run_on_all();
522
523 double max_amips = std::numeric_limits<double>::lowest();
524 double min_amips = std::numeric_limits<double>::max();
525
526 for (const auto& t : mesh->get_all(mesh->top_simplex_type())) {
527 // double e = amips->get_value(simplex::Simplex(mesh->top_simplex_type(), t));
528 double e = amips_accessor.scalar_attribute(t);
529 max_amips = std::max(max_amips, e);
530 min_amips = std::min(min_amips, e);
531 }
532
533 logger().info("Initial Max AMIPS Energy: {}, Min AMIPS Energy: {}", max_amips, min_amips);
534
536 // sizing field scalar
538
539 auto sizing_field_scalar_attribute = mesh->register_attribute<double>(
540 "sizing_field_scalar",
542 1,
543 false,
544 1); // defaults to 1
545
546
548 // Storing target edge length
549 auto target_edge_length_attribute = mesh->register_attribute<double>(
550 "wildmeshing_target_edge_length",
552 1,
553 false,
554 target_edge_length); // defaults to target edge length
555
556 // Target edge length update
557 const double min_edge_length = [&]() -> double {
558 if (options.envelopes.empty()) {
559 return 1e-6; // some default value if no envelope exists
560 } else {
561 // use envelope thickness if available
562 double r = 0;
563 for (const auto& e : options.envelopes) {
564 r = std::max(r, e.thickness);
565 }
566 assert(r > 0);
567 return r;
568 }
569 }();
570 const double target_max_amips = options.target_max_amips;
571
572 // Target Edge length update
573 // auto compute_target_edge_length = [&](const Eigen::Vector2d& P) -> double {
574 // return (P[0] + P[1]) / 2 * target_edge_length;
575 // };
576 // auto target_edge_length_update =
577 // std::make_shared<wmtk::operations::SingleAttributeTransferStrategy<double, double>>(
578 // target_edge_length_attribute,
579 // sizing_field_scalar_attribute,
580 // compute_target_edge_length);
581
582 // auto compute_target_edge_length =
583 // [target_edge_length,
584 // target_max_amips,
585 // min_edge_length,
586 // target_edge_length_attribute,
587 // &mesh](const Eigen::MatrixXd& P, const std::vector<Tuple>& neighs) ->
588 // Eigen::VectorXd {
589 // auto target_edge_length_accessor =
590 // mesh->create_accessor(target_edge_length_attribute.as<double>());
591
592 // assert(P.rows() == 1); // rows --> attribute dimension
593 // assert(!neighs.empty());
594 // assert(P.cols() == neighs.size());
595 // const double current_target_edge_length =
596 // target_edge_length_accessor.const_scalar_attribute(neighs[0]);
597 // const double max_amips = P.maxCoeff();
598
599 // double new_target_edge_length = current_target_edge_length;
600 // if (max_amips > target_max_amips) {
601 // new_target_edge_length *= 0.5;
602 // } else {
603 // new_target_edge_length *= 1.5;
604 // }
605 // new_target_edge_length =
606 // std::min(new_target_edge_length, target_edge_length); // upper bound
607 // new_target_edge_length = std::max(new_target_edge_length, min_edge_length); // lower bound
608
609 // return Eigen::VectorXd::Constant(1, new_target_edge_length);
610 // };
611 // auto target_edge_length_update =
612 // std::make_shared<wmtk::operations::SingleAttributeTransferStrategy<double, double>>(
613 // target_edge_length_attribute,
614 // amips_attribute,
615 // compute_target_edge_length);
616
617
619 // auto compute_target_edge_length = [](const Eigen::MatrixX<Rational>& P) ->
620 // Eigen::VectorXd {
621 // assert(P.cols() == 2); // cols --> number of neighbors
622 // assert(P.rows() == 2 || P.rows() == 3); // rows --> attribute dimension
623 // const double x_avg = 0.5 * (P(0, 0) + P(0, 1)).to_double();
624 // const double target_length = x_avg * x_avg;
625 // return Eigen::VectorXd::Constant(1, target_length);
626 //};
627 // auto target_edge_length_update =
628 // std::make_shared<wmtk::operations::SingleAttributeTransferStrategy<double, Rational>>(
629 // target_edge_length_attribute,
630 // rpt_attribute,
631 // compute_target_edge_length);
632 // target_edge_length_update->run_on_all();
633
634
636 // Storing edge lengths
637 auto edge_length_attribute =
638 mesh->register_attribute<double>("edge_length", PrimitiveType::Edge, 1);
639 auto edge_length_accessor = mesh->create_accessor(edge_length_attribute.as<double>());
640 // Edge length update
641 auto compute_edge_length = [](const Eigen::MatrixX<Rational>& P) -> Eigen::VectorXd {
642 assert(P.cols() == 2);
643 assert(P.rows() == 2 || P.rows() == 3);
644 return Eigen::VectorXd::Constant(1, sqrt((P.col(0) - P.col(1)).squaredNorm().to_double()));
645 };
646 auto edge_length_update =
647 std::make_shared<wmtk::operations::SingleAttributeTransferStrategy<double, Rational>>(
648 edge_length_attribute,
649 pt_attribute,
650 compute_edge_length);
651 edge_length_update->run_on_all();
652
654 // default transfer
655 auto pass_through_attributes = utils::get_attributes(cache, *mesh, options.pass_through);
656 pass_through_attributes.push_back(edge_length_attribute);
657 pass_through_attributes.push_back(amips_attribute);
658 // pass_through_attributes.push_back(target_edge_length_attribute);
659
660
662 // Lambdas for priority
664 auto long_edges_first = [&](const simplex::Simplex& s) {
665 assert(s.primitive_type() == PrimitiveType::Edge);
666 return -edge_length_accessor.scalar_attribute(s.tuple());
667 };
668 auto short_edges_first = [&](const simplex::Simplex& s) {
669 assert(s.primitive_type() == PrimitiveType::Edge);
670 return edge_length_accessor.scalar_attribute(s.tuple());
671 };
672
673
675 // envelopes
677
678 using MeshConstrainPair = ProjectOperation::MeshConstrainPair;
679
680 auto envelope_invariant = std::make_shared<InvariantCollection>(*mesh);
681 std::vector<std::shared_ptr<AttributeTransferStrategyBase>> update_child_position,
682 update_parent_position;
683 std::vector<std::shared_ptr<Mesh>> envelopes;
684 std::vector<MeshConstrainPair> mesh_constraint_pairs;
685
686 std::vector<std::shared_ptr<Mesh>> multimesh_meshes;
687
688 // assert(options.envelopes.size() == 4);
689 // four kind of envelopes in tetwild [surface_mesh,
690 // open_boudnary, nonmanifold_edges, is_boundary(bbox)]
691 // TODO: add nonmanifold vertex point mesh
692
693 for (const auto& v : options.envelopes) {
694 auto envelope = cache.read_mesh(v.geometry.mesh);
695 if (envelope == nullptr) {
696 wmtk::logger().info("TetWild: no {} mesh for this mesh", v.geometry.mesh);
697 continue;
698 } else {
699 wmtk::logger().info("TetWild: registered {} mesh", v.geometry.mesh);
700 }
701 envelopes.emplace_back(envelope);
702
703 wmtk::logger().trace(
704 "TetWild: getting constrained position {} from {}",
705 v.constrained_position,
706 v.geometry.mesh);
707 auto constrained = utils::get_attributes(cache, *mesh, v.constrained_position);
708 multimesh_meshes.push_back(constrained.front().mesh().shared_from_this());
709 assert(constrained.size() == 1);
710 pass_through_attributes.emplace_back(constrained.front());
711
712 wmtk::logger().trace(
713 "TetWild: Constructing envelope position handle {} == input = {}",
714 v.geometry.mesh,
715 v.geometry.mesh == "input");
716
717 const bool has_double_pos =
718 envelope->has_attribute<double>(v.geometry.position, PrimitiveType::Vertex);
719 const bool has_rational_pos =
720 envelope->has_attribute<Rational>(v.geometry.position, PrimitiveType::Vertex);
721 assert(has_double_pos || has_rational_pos);
722 assert(!(v.geometry.mesh == "input") || has_double_pos);
723
724 auto envelope_position_handle =
725 has_double_pos
726 ? envelope->get_attribute_handle<double>(v.geometry.position, PrimitiveType::Vertex)
727 : envelope->get_attribute_handle<Rational>(
728 v.geometry.position,
730 // envelope->get_attribute_handle<Rational>(v.geometry.position, PrimitiveType::Vertex);
731
732
733 mesh_constraint_pairs.emplace_back(envelope_position_handle, constrained.front());
734
735 envelope_invariant->add(std::make_shared<EnvelopeInvariant>(
736 envelope_position_handle,
737 v.thickness * bbdiag,
738 constrained.front()));
739
740 update_parent_position.emplace_back(attribute_update::make_cast_attribute_transfer_strategy(
741 /*source=*/constrained.front(),
742 /*target=*/pt_attribute));
743
744 update_child_position.emplace_back(attribute_update::make_cast_attribute_transfer_strategy(
745 /*source=*/pt_attribute,
746 /*target=*/constrained.front()));
747 }
748
750
751
753 // Invariants
755
756 wmtk::logger().trace("Going through invariants");
757 auto inversion_invariant =
758 std::make_shared<SimplexInversionInvariant<Rational>>(*mesh, pt_attribute.as<Rational>());
759
760 std::shared_ptr<function::PerSimplexFunction> amips =
761 std::make_shared<AMIPS>(*mesh, pt_attribute);
762 // auto function_invariant = std::make_shared<FunctionInvariant>(mesh->top_simplex_type(),
763 // amips);
764 auto function_invariant =
765 std::make_shared<MaxFunctionInvariant>(mesh->top_simplex_type(), amips);
766
767 auto link_condition = std::make_shared<MultiMeshLinkConditionInvariant>(*mesh);
768
769 auto todo_larger = std::make_shared<TodoLargerInvariant>(
770 *mesh,
771 edge_length_attribute.as<double>(),
772 target_edge_length_attribute.as<double>(),
773 4.0 / 3.0);
774
775 auto todo_smaller = std::make_shared<TodoSmallerInvariant>(
776 *mesh,
777 edge_length_attribute.as<double>(),
778 target_edge_length_attribute.as<double>(),
779 4.0 / 5.0);
780
781
782 auto interior_edge = std::make_shared<InteriorEdgeInvariant>(*mesh);
783 auto interior_face = std::make_shared<InteriorSimplexInvariant>(*mesh, PrimitiveType::Triangle);
784
785 for (const auto& em : multimesh_meshes) {
786 interior_edge->add_boundary(*em);
787 interior_face->add_boundary(*em);
788 }
789
790 auto valence_3 = std::make_shared<EdgeValenceInvariant>(*mesh, 3);
791 auto valence_4 = std::make_shared<EdgeValenceInvariant>(*mesh, 4);
792 auto swap44_energy_before =
793 std::make_shared<Swap44EnergyBeforeInvariant>(*mesh, pt_attribute.as<Rational>(), 0);
794 auto swap44_2_energy_before =
795 std::make_shared<Swap44EnergyBeforeInvariant>(*mesh, pt_attribute.as<Rational>(), 1);
796 auto swap32_energy_before =
797 std::make_shared<Swap32EnergyBeforeInvariant>(*mesh, pt_attribute.as<Rational>());
798 auto swap23_energy_before =
799 std::make_shared<Swap23EnergyBeforeInvariant>(*mesh, pt_attribute.as<Rational>());
800
801 auto invariant_separate_substructures =
802 std::make_shared<invariants::SeparateSubstructuresInvariant>(*mesh);
803
805 // renew flags
807 auto visited_edge_flag =
808 mesh->register_attribute<char>("visited_edge", PrimitiveType::Edge, 1, false, char(1));
809
810 auto update_flag_func = [](const Eigen::MatrixX<Rational>& P) -> Eigen::VectorX<char> {
811 assert(P.cols() == 2);
812 assert(P.rows() == 2 || P.rows() == 3);
813 return Eigen::VectorX<char>::Constant(1, char(1));
814 };
815 auto tag_update =
816 std::make_shared<wmtk::operations::SingleAttributeTransferStrategy<char, Rational>>(
817 visited_edge_flag,
818 pt_attribute,
819 update_flag_func);
820
822 // sizing field update flags
824 auto visited_vertex_flag =
825 mesh->register_attribute<char>("visited_vertex", PrimitiveType::Vertex, 1, false, char(1));
826 pass_through_attributes.push_back(visited_vertex_flag);
827
829 // energy filter flag
831 auto energy_filter_attribute =
832 mesh->register_attribute<char>("energy_filter", PrimitiveType::Vertex, 1, false, char(1));
833
834 auto energy_filter_accessor = mesh->create_accessor<char>(energy_filter_attribute);
835
836 auto update_energy_filter_func = [](const Eigen::MatrixX<Rational>& P) -> Eigen::VectorX<char> {
837 // assert(P.cols() == 2);
838 // assert(P.rows() == 2 || P.rows() == 3);
839 return Eigen::VectorX<char>::Constant(1, char(1));
840 };
841 auto energy_filter_update =
842 std::make_shared<wmtk::operations::SingleAttributeTransferStrategy<char, Rational>>(
843 energy_filter_attribute,
844 pt_attribute,
845 update_energy_filter_func);
846
847 pass_through_attributes.push_back(energy_filter_attribute);
848
850 // coloring flag for smoothing
852
853 // auto coloring_attribute =
854 // mesh->register_attribute<int64_t>("coloring", PrimitiveType::Vertex, 1, false, -1);
855 // pass_through_attributes.push_back(coloring_attribute);
856
858 // Creation of the 4 ops
860 std::vector<std::shared_ptr<Operation>> ops;
861 std::vector<std::string> ops_name;
862
864 // 0) Rounding
866 auto rounding_pt_attribute = mesh->get_attribute_handle_typed<Rational>(
867 options.attributes.position,
869 auto rounding = std::make_shared<Rounding>(*mesh, rounding_pt_attribute);
870 rounding->add_invariant(
871 std::make_shared<RoundedInvariant>(*mesh, pt_attribute.as<Rational>(), true));
872 rounding->add_invariant(inversion_invariant);
873
874
876 // 1) EdgeSplit
878 auto split = std::make_shared<EdgeSplit>(*mesh);
879 split->set_priority(long_edges_first);
880
881 split->add_invariant(todo_larger);
882 split->add_invariant(inversion_invariant);
883
884 split->set_new_attribute_strategy(pt_attribute);
885 split->set_new_attribute_strategy(sizing_field_scalar_attribute);
886 split->set_new_attribute_strategy(
887 visited_edge_flag,
890 split->set_new_attribute_strategy(
891 target_edge_length_attribute,
894
895 for (const auto& attr : pass_through_attributes) {
896 split->set_new_attribute_strategy(
897 attr,
900 }
901
902 split->add_transfer_strategy(amips_update);
903 split->add_transfer_strategy(edge_length_update);
904 split->add_transfer_strategy(tag_update); // for renew the queue
905 split->add_transfer_strategy(energy_filter_update);
906
907 // split->add_transfer_strategy(target_edge_length_update);
908
909 auto split_then_round = std::make_shared<AndOperationSequence>(*mesh);
910 split_then_round->add_operation(split);
911 split_then_round->add_operation(rounding);
912
913 for (auto& s : update_child_position) {
914 split_then_round->add_transfer_strategy(s);
915 }
916
917 // split unrounded
918 auto split_unrounded = std::make_shared<EdgeSplit>(*mesh);
919 split_unrounded->set_priority(long_edges_first);
920
921 split_unrounded->add_invariant(todo_larger);
922
923 auto split_unrounded_transfer_strategy =
924 std::make_shared<SplitNewAttributeStrategy<Rational>>(pt_attribute);
925 split_unrounded_transfer_strategy->set_strategy(
926 [](const Eigen::VectorX<Rational>& a, const std::bitset<2>&) {
927 return std::array<Eigen::VectorX<Rational>, 2>{{a, a}};
928 });
929 split_unrounded_transfer_strategy->set_rib_strategy(
930 [](const Eigen::VectorX<Rational>& p0_d,
931 const Eigen::VectorX<Rational>& p1_d,
932 const std::bitset<2>& bs) -> Eigen::VectorX<Rational> {
933 Eigen::VectorX<Rational> p0(p0_d.size());
934 Eigen::VectorX<Rational> p1(p1_d.size());
935 for (int i = 0; i < p0_d.size(); ++i) {
936 p0[i] = Rational(p0_d[i], false);
937 p1[i] = Rational(p1_d[i], false);
938 }
939 if (bs[0] == bs[1]) {
940 return (p0 + p1) / Rational(2, false);
941 } else if (bs[0]) {
942 return p0;
943
944 } else {
945 return p1;
946 }
947 });
948
949 split_unrounded->set_new_attribute_strategy(pt_attribute, split_unrounded_transfer_strategy);
950 split_unrounded->set_new_attribute_strategy(sizing_field_scalar_attribute);
951 split_unrounded->set_new_attribute_strategy(
952 visited_edge_flag,
955 for (const auto& attr : pass_through_attributes) {
956 split_unrounded->set_new_attribute_strategy(
957 attr,
960 }
961 split_unrounded->set_new_attribute_strategy(
962 target_edge_length_attribute,
965
966 split_unrounded->add_transfer_strategy(amips_update);
967 split_unrounded->add_transfer_strategy(edge_length_update);
968 split_unrounded->add_transfer_strategy(tag_update); // for renew the queue
969 split_unrounded->add_transfer_strategy(energy_filter_update);
970
971 // split->add_transfer_strategy(target_edge_length_update);
972
973 auto split_sequence = std::make_shared<OrOperationSequence>(*mesh);
974 split_sequence->add_operation(split_then_round);
975 split_sequence->add_operation(split_unrounded);
976 split_sequence->add_invariant(
977 std::make_shared<EnergyFilterInvariant>(*mesh, energy_filter_attribute.as<char>()));
978
979 split_sequence->set_priority(long_edges_first);
980
981 ops.emplace_back(split_sequence);
982 ops_name.emplace_back("SPLIT");
983
984 ops.emplace_back(rounding);
985 ops_name.emplace_back("rounding");
986
987
989 // collapse transfer
991 auto clps_strat1 = std::make_shared<CollapseNewAttributeStrategy<Rational>>(pt_attribute);
992 // clps_strat1->set_simplex_predicate(BasicSimplexPredicate::IsInterior);
993 // clps_strat1->set_strategy(CollapseBasicStrategy::Default);
994 clps_strat1->set_strategy(CollapseBasicStrategy::CopyOther);
995
996 auto clps_strat2 = std::make_shared<CollapseNewAttributeStrategy<Rational>>(pt_attribute);
997 // clps_strat2->set_simplex_predicate(BasicSimplexPredicate::IsInterior);
998 // clps_strat2->set_strategy(CollapseBasicStrategy::Default);
999 clps_strat2->set_strategy(CollapseBasicStrategy::CopyTuple);
1000
1002 // 2) EdgeCollapse
1004
1005 auto setup_collapse = [&](std::shared_ptr<EdgeCollapse>& collapse) {
1006 collapse->add_invariant(invariant_separate_substructures);
1007 collapse->add_invariant(link_condition);
1008 collapse->add_invariant(inversion_invariant);
1009 // collapse->add_invariant(function_invariant);
1010 collapse->add_invariant(envelope_invariant);
1011
1012 collapse->set_new_attribute_strategy(
1013 visited_edge_flag,
1015
1016 collapse->add_transfer_strategy(tag_update);
1017 collapse->add_transfer_strategy(energy_filter_update);
1018 for (const auto& attr : pass_through_attributes) {
1019 collapse->set_new_attribute_strategy(
1020 attr,
1022 }
1023 collapse->set_new_attribute_strategy(
1024 target_edge_length_attribute,
1026 // THis triggers a segfault in release
1027 // solved somehow
1028 // collapse->set_priority(short_edges_first);
1029
1030 collapse->add_transfer_strategy(amips_update);
1031 collapse->add_transfer_strategy(edge_length_update);
1032 // collapse->add_transfer_strategy(target_edge_length_update);
1033
1034 for (auto& s : update_child_position) {
1035 collapse->add_transfer_strategy(s);
1036 }
1037 };
1038
1039 auto collapse1 = std::make_shared<EdgeCollapse>(*mesh);
1040 collapse1->add_invariant(std::make_shared<CollapseEnergyBeforeInvariant>(
1041 *mesh,
1042 pt_attribute.as<Rational>(),
1043 amips_attribute.as<double>(),
1044 1));
1045
1046 collapse1->set_new_attribute_strategy(pt_attribute, clps_strat1);
1047 collapse1->set_new_attribute_strategy(sizing_field_scalar_attribute, clps_strat1);
1048 setup_collapse(collapse1);
1049
1050 auto collapse2 = std::make_shared<EdgeCollapse>(*mesh);
1051 collapse2->add_invariant(std::make_shared<CollapseEnergyBeforeInvariant>(
1052 *mesh,
1053 pt_attribute.as<Rational>(),
1054 amips_attribute.as<double>(),
1055 0));
1056
1057 collapse2->set_new_attribute_strategy(pt_attribute, clps_strat2);
1058 collapse2->set_new_attribute_strategy(sizing_field_scalar_attribute, clps_strat2);
1059 setup_collapse(collapse2);
1060
1061 auto collapse = std::make_shared<OrOperationSequence>(*mesh);
1062 collapse->add_operation(collapse1);
1063 collapse->add_operation(collapse2);
1064 collapse->add_invariant(todo_smaller);
1065
1066 auto collapse_then_round = std::make_shared<AndOperationSequence>(*mesh);
1067 collapse_then_round->add_operation(collapse);
1068 collapse_then_round->add_operation(rounding);
1069
1070 collapse_then_round->set_priority(short_edges_first);
1071 collapse_then_round->add_invariant(
1072 std::make_shared<EnergyFilterInvariant>(*mesh, energy_filter_attribute.as<char>()));
1073
1074
1075 for (auto& s : update_child_position) {
1076 collapse_then_round->add_transfer_strategy(s);
1077 }
1078
1079 ops.emplace_back(collapse_then_round);
1080 ops_name.emplace_back("COLLAPSE");
1081
1082 ops.emplace_back(rounding);
1083 ops_name.emplace_back("rounding");
1084
1086 // 3) Swap
1088
1089 auto setup_swap = [&](Operation& op,
1090 EdgeCollapse& collapse,
1091 EdgeSplit& split,
1092 std::shared_ptr<Invariant> simplex_invariant,
1093 bool is_edge = true) {
1094 if (is_edge) op.set_priority(long_edges_first);
1095
1096 op.add_invariant(simplex_invariant);
1097 op.add_invariant(inversion_invariant);
1098 // op.add_invariant(function_invariant);
1099
1100 op.add_transfer_strategy(amips_update);
1101 op.add_transfer_strategy(edge_length_update);
1102 // op.add_transfer_strategy(target_edge_length_update);
1103 // for (auto& s : update_child_position) {
1104 // op.add_transfer_strategy(s);
1105 // }
1106
1107 collapse.add_invariant(invariant_separate_substructures);
1108 collapse.add_invariant(link_condition);
1109
1110
1111 collapse.set_new_attribute_strategy(pt_attribute, CollapseBasicStrategy::CopyOther);
1112 split.set_new_attribute_strategy(pt_attribute);
1113
1114 split.set_new_attribute_strategy(
1115 target_edge_length_attribute,
1118 collapse.set_new_attribute_strategy(
1119 target_edge_length_attribute,
1121
1122 // this might not be necessary
1123 // for (auto& s : update_child_position) {
1124 // collapse.add_transfer_strategy(s);
1125 // split.add_transfer_strategy(s);
1126 // }
1127
1128
1129 for (const auto& attr : pass_through_attributes) {
1130 split.set_new_attribute_strategy(
1131 attr,
1134 collapse.set_new_attribute_strategy(
1135 attr,
1137 }
1138 };
1139
1140
1141 if (mesh->top_simplex_type() == PrimitiveType::Triangle) {
1142 auto swap = std::make_shared<TriEdgeSwap>(*mesh);
1143 setup_swap(*swap, swap->collapse(), swap->split(), interior_edge);
1144
1145 ops.push_back(swap);
1146 ops_name.push_back("swap");
1147
1148 ops.emplace_back(rounding);
1149 ops_name.emplace_back("rounding");
1150
1151 } else if (mesh->top_simplex_type() == PrimitiveType::Tetrahedron) {
1152 auto swap56 = std::make_shared<MinOperationSequence>(*mesh);
1153 for (int i = 0; i < 5; ++i) {
1154 auto swap = std::make_shared<TetEdgeSwap>(*mesh, i);
1155 swap->collapse().add_invariant(invariant_separate_substructures);
1156 swap->collapse().add_invariant(link_condition);
1157 swap->collapse().set_new_attribute_strategy(
1158 pt_attribute,
1159 CollapseBasicStrategy::CopyOther);
1160 swap->collapse().set_new_attribute_strategy(
1161 sizing_field_scalar_attribute,
1162 CollapseBasicStrategy::CopyOther);
1163 swap->split().set_new_attribute_strategy(pt_attribute);
1164 swap->split().set_new_attribute_strategy(sizing_field_scalar_attribute);
1165 // swap->split().add_transfer_strategy(amips_update);
1166 // swap->collapse().add_transfer_strategy(amips_update);
1167 swap->split().set_new_attribute_strategy(
1168 visited_edge_flag,
1171 swap->collapse().set_new_attribute_strategy(
1172 visited_edge_flag,
1174
1175 swap->split().set_new_attribute_strategy(
1176 target_edge_length_attribute,
1179 swap->collapse().set_new_attribute_strategy(
1180 target_edge_length_attribute,
1182
1183 swap->add_invariant(std::make_shared<Swap56EnergyBeforeInvariant>(
1184 *mesh,
1185 pt_attribute.as<Rational>(),
1186 i));
1187
1188 swap->add_transfer_strategy(amips_update);
1189
1190 // swap->add_invariant(inversion_invariant);
1191 swap->collapse().add_invariant(inversion_invariant);
1192
1193 // swap->collapse().add_invariant(envelope_invariant);
1194
1195 for (const auto& attr : pass_through_attributes) {
1196 swap->split().set_new_attribute_strategy(
1197 attr,
1200 swap->collapse().set_new_attribute_strategy(
1201 attr,
1203 }
1204
1205 swap56->add_operation(swap);
1206 }
1207
1208 auto swap56_energy_check = [&](int64_t idx, const simplex::Simplex& t) -> double {
1209 constexpr static PrimitiveType PV = PrimitiveType::Vertex;
1210 constexpr static PrimitiveType PE = PrimitiveType::Edge;
1211 constexpr static PrimitiveType PF = PrimitiveType::Triangle;
1212 constexpr static PrimitiveType PT = PrimitiveType::Tetrahedron;
1213
1214 auto accessor = mesh->create_const_accessor(pt_attribute.as<Rational>());
1215
1216 const Tuple e0 = t.tuple();
1217 const Tuple e1 = mesh->switch_tuple(e0, PV);
1218
1219 std::array<Tuple, 5> v;
1220 auto iter_tuple = e0;
1221 for (int64_t i = 0; i < 5; ++i) {
1222 v[i] = mesh->switch_tuples(iter_tuple, {PE, PV});
1223 iter_tuple = mesh->switch_tuples(iter_tuple, {PF, PT});
1224 }
1225 if (iter_tuple != e0) return 0;
1226 assert(iter_tuple == e0);
1227
1228 // five iterable vertices remap to 0-4 by m_collapse_index, 0: m_collapse_index, 5:
1229 // e0, 6: e1
1230 std::array<Eigen::Vector3<Rational>, 7> positions = {
1231 {accessor.const_vector_attribute(v[(idx + 0) % 5]),
1232 accessor.const_vector_attribute(v[(idx + 1) % 5]),
1233 accessor.const_vector_attribute(v[(idx + 2) % 5]),
1234 accessor.const_vector_attribute(v[(idx + 3) % 5]),
1235 accessor.const_vector_attribute(v[(idx + 4) % 5]),
1236 accessor.const_vector_attribute(e0),
1237 accessor.const_vector_attribute(e1)}};
1238
1239 std::array<Eigen::Vector3d, 7> positions_double = {
1240 {positions[0].cast<double>(),
1241 positions[1].cast<double>(),
1242 positions[2].cast<double>(),
1243 positions[3].cast<double>(),
1244 positions[4].cast<double>(),
1245 positions[5].cast<double>(),
1246 positions[6].cast<double>()}};
1247
1248 std::array<std::array<int, 4>, 6> new_tets = {
1249 {{{0, 1, 2, 5}},
1250 {{0, 2, 3, 5}},
1251 {{0, 3, 4, 5}},
1252 {{0, 1, 2, 6}},
1253 {{0, 2, 3, 6}},
1254 {{0, 3, 4, 6}}}};
1255
1256 double new_energy_max = std::numeric_limits<double>::lowest();
1257
1258 for (int i = 0; i < 6; ++i) {
1260 positions[new_tets[i][0]],
1261 positions[new_tets[i][1]],
1262 positions[new_tets[i][2]],
1263 positions[new_tets[i][3]]) > 0) {
1265 positions_double[new_tets[i][0]][0],
1266 positions_double[new_tets[i][0]][1],
1267 positions_double[new_tets[i][0]][2],
1268 positions_double[new_tets[i][1]][0],
1269 positions_double[new_tets[i][1]][1],
1270 positions_double[new_tets[i][1]][2],
1271 positions_double[new_tets[i][2]][0],
1272 positions_double[new_tets[i][2]][1],
1273 positions_double[new_tets[i][2]][2],
1274 positions_double[new_tets[i][3]][0],
1275 positions_double[new_tets[i][3]][1],
1276 positions_double[new_tets[i][3]][2],
1277 }});
1278
1279
1280 if (energy > new_energy_max) new_energy_max = energy;
1281 } else {
1283 positions_double[new_tets[i][1]][0],
1284 positions_double[new_tets[i][1]][1],
1285 positions_double[new_tets[i][1]][2],
1286 positions_double[new_tets[i][0]][0],
1287 positions_double[new_tets[i][0]][1],
1288 positions_double[new_tets[i][0]][2],
1289 positions_double[new_tets[i][2]][0],
1290 positions_double[new_tets[i][2]][1],
1291 positions_double[new_tets[i][2]][2],
1292 positions_double[new_tets[i][3]][0],
1293 positions_double[new_tets[i][3]][1],
1294 positions_double[new_tets[i][3]][2],
1295 }});
1296
1297
1298 if (energy > new_energy_max) new_energy_max = energy;
1299 }
1300 }
1301
1302 return new_energy_max;
1303 };
1304
1305 swap56->set_value_function(swap56_energy_check);
1306 swap56->add_invariant(std::make_shared<EdgeValenceInvariant>(*mesh, 5));
1307
1308 // swap44
1309
1310 auto swap44 = std::make_shared<MinOperationSequence>(*mesh);
1311 for (int i = 0; i < 2; ++i) {
1312 auto swap = std::make_shared<TetEdgeSwap>(*mesh, i);
1313 swap->collapse().add_invariant(invariant_separate_substructures);
1314 swap->collapse().add_invariant(link_condition);
1315 swap->collapse().set_new_attribute_strategy(
1316 pt_attribute,
1317 CollapseBasicStrategy::CopyOther);
1318 swap->collapse().set_new_attribute_strategy(
1319 sizing_field_scalar_attribute,
1320 CollapseBasicStrategy::CopyOther);
1321 swap->split().set_new_attribute_strategy(pt_attribute);
1322 swap->split().set_new_attribute_strategy(sizing_field_scalar_attribute);
1323 swap->split().set_new_attribute_strategy(
1324 visited_edge_flag,
1327 swap->collapse().set_new_attribute_strategy(
1328 visited_edge_flag,
1330
1331 // swap->split().add_transfer_strategy(amips_update);
1332 // swap->collapse().add_transfer_strategy(amips_update);
1333
1334 swap->split().set_new_attribute_strategy(
1335 target_edge_length_attribute,
1338 swap->collapse().set_new_attribute_strategy(
1339 target_edge_length_attribute,
1341
1342 swap->add_transfer_strategy(amips_update);
1343
1344 swap->add_invariant(std::make_shared<Swap44EnergyBeforeInvariant>(
1345 *mesh,
1346 pt_attribute.as<Rational>(),
1347 i));
1348
1349 // swap->add_invariant(inversion_invariant);
1350 swap->collapse().add_invariant(inversion_invariant);
1351
1352 // swap->collapse().add_invariant(envelope_invariant);
1353
1354 for (const auto& attr : pass_through_attributes) {
1355 swap->split().set_new_attribute_strategy(
1356 attr,
1359 swap->collapse().set_new_attribute_strategy(
1360 attr,
1362 }
1363
1364 swap44->add_operation(swap);
1365 }
1366
1367 auto swap44_energy_check = [&](int64_t idx, const simplex::Simplex& t) -> double {
1368 constexpr static PrimitiveType PV = PrimitiveType::Vertex;
1369 constexpr static PrimitiveType PE = PrimitiveType::Edge;
1370 constexpr static PrimitiveType PF = PrimitiveType::Triangle;
1371 constexpr static PrimitiveType PT = PrimitiveType::Tetrahedron;
1372
1373 auto accessor = mesh->create_const_accessor(pt_attribute.as<Rational>());
1374
1375 // get the coords of the vertices
1376 // input edge end points
1377 const Tuple e0 = t.tuple();
1378 const Tuple e1 = mesh->switch_tuple(e0, PV);
1379 // other four vertices
1380 std::array<Tuple, 4> v;
1381 auto iter_tuple = e0;
1382 for (int64_t i = 0; i < 4; ++i) {
1383 v[i] = mesh->switch_tuples(iter_tuple, {PE, PV});
1384 iter_tuple = mesh->switch_tuples(iter_tuple, {PF, PT});
1385 }
1386
1387 if (iter_tuple != e0) return 0;
1388 assert(iter_tuple == e0);
1389
1390 std::array<Eigen::Vector3<Rational>, 6> positions = {
1391 {accessor.const_vector_attribute(v[(idx + 0) % 4]),
1392 accessor.const_vector_attribute(v[(idx + 1) % 4]),
1393 accessor.const_vector_attribute(v[(idx + 2) % 4]),
1394 accessor.const_vector_attribute(v[(idx + 3) % 4]),
1395 accessor.const_vector_attribute(e0),
1396 accessor.const_vector_attribute(e1)}};
1397 std::array<Eigen::Vector3d, 6> positions_double = {
1398 {positions[0].cast<double>(),
1399 positions[1].cast<double>(),
1400 positions[2].cast<double>(),
1401 positions[3].cast<double>(),
1402 positions[4].cast<double>(),
1403 positions[5].cast<double>()}};
1404
1405 std::array<std::array<int, 4>, 4> new_tets = {
1406 {{{0, 1, 2, 4}}, {{0, 2, 3, 4}}, {{0, 1, 2, 5}}, {{0, 2, 3, 5}}}};
1407
1408 double new_energy_max = std::numeric_limits<double>::lowest();
1409
1410 for (int i = 0; i < 4; ++i) {
1412 positions[new_tets[i][0]],
1413 positions[new_tets[i][1]],
1414 positions[new_tets[i][2]],
1415 positions[new_tets[i][3]]) > 0) {
1417 positions_double[new_tets[i][0]][0],
1418 positions_double[new_tets[i][0]][1],
1419 positions_double[new_tets[i][0]][2],
1420 positions_double[new_tets[i][1]][0],
1421 positions_double[new_tets[i][1]][1],
1422 positions_double[new_tets[i][1]][2],
1423 positions_double[new_tets[i][2]][0],
1424 positions_double[new_tets[i][2]][1],
1425 positions_double[new_tets[i][2]][2],
1426 positions_double[new_tets[i][3]][0],
1427 positions_double[new_tets[i][3]][1],
1428 positions_double[new_tets[i][3]][2],
1429 }});
1430
1431 if (energy > new_energy_max) new_energy_max = energy;
1432 } else {
1434 positions_double[new_tets[i][1]][0],
1435 positions_double[new_tets[i][1]][1],
1436 positions_double[new_tets[i][1]][2],
1437 positions_double[new_tets[i][0]][0],
1438 positions_double[new_tets[i][0]][1],
1439 positions_double[new_tets[i][0]][2],
1440 positions_double[new_tets[i][2]][0],
1441 positions_double[new_tets[i][2]][1],
1442 positions_double[new_tets[i][2]][2],
1443 positions_double[new_tets[i][3]][0],
1444 positions_double[new_tets[i][3]][1],
1445 positions_double[new_tets[i][3]][2],
1446 }});
1447
1448 if (energy > new_energy_max) new_energy_max = energy;
1449 }
1450 }
1451
1452 return new_energy_max;
1453 };
1454
1455 swap44->set_value_function(swap44_energy_check);
1456 swap44->add_invariant(std::make_shared<EdgeValenceInvariant>(*mesh, 4));
1457
1458 // swap 32
1459 auto swap32 = std::make_shared<TetEdgeSwap>(*mesh, 0);
1460 swap32->add_invariant(std::make_shared<EdgeValenceInvariant>(*mesh, 3));
1461 swap32->add_invariant(
1462 std::make_shared<Swap32EnergyBeforeInvariant>(*mesh, pt_attribute.as<Rational>()));
1463
1464 swap32->collapse().add_invariant(invariant_separate_substructures);
1465 swap32->collapse().add_invariant(link_condition);
1466 swap32->collapse().set_new_attribute_strategy(
1467 pt_attribute,
1468 CollapseBasicStrategy::CopyOther);
1469 swap32->collapse().set_new_attribute_strategy(
1470 sizing_field_scalar_attribute,
1471 CollapseBasicStrategy::CopyOther);
1472 // swap32->add_invariant(inversion_invariant);
1473 swap32->split().set_new_attribute_strategy(pt_attribute);
1474 swap32->split().set_new_attribute_strategy(sizing_field_scalar_attribute);
1475 swap32->split().set_new_attribute_strategy(
1476 visited_edge_flag,
1479 swap32->collapse().set_new_attribute_strategy(
1480 visited_edge_flag,
1482
1483 // swap32->split().add_transfer_strategy(amips_update);
1484 // swap32->collapse().add_transfer_strategy(amips_update);
1485
1486 swap32->split().set_new_attribute_strategy(
1487 target_edge_length_attribute,
1490 swap32->collapse().set_new_attribute_strategy(
1491 target_edge_length_attribute,
1493
1494 swap32->add_transfer_strategy(amips_update);
1495
1496 // hack
1497 swap32->collapse().add_invariant(inversion_invariant);
1498 // swap32->collapse().add_invariant(envelope_invariant);
1499
1500 for (const auto& attr : pass_through_attributes) {
1501 swap32->split().set_new_attribute_strategy(
1502 attr,
1505 swap32->collapse().set_new_attribute_strategy(
1506 attr,
1508 }
1509
1510 // all swaps
1511
1512 auto swap_all = std::make_shared<OrOperationSequence>(*mesh);
1513 swap_all->add_operation(swap32);
1514 swap_all->add_operation(swap44);
1515 swap_all->add_operation(swap56);
1516 swap_all->add_transfer_strategy(tag_update);
1517 swap_all->add_transfer_strategy(energy_filter_update);
1518
1519 auto swap_then_round = std::make_shared<AndOperationSequence>(*mesh);
1520 swap_then_round->add_operation(swap_all);
1521 swap_then_round->add_operation(rounding);
1522 swap_then_round->add_invariant(
1523 std::make_shared<EnergyFilterInvariant>(*mesh, energy_filter_attribute.as<char>()));
1524 swap_then_round->add_invariant(std::make_shared<InteriorEdgeInvariant>(*mesh));
1525 swap_then_round->add_invariant(std::make_shared<NoChildMeshAttachingInvariant>(*mesh));
1526
1527 swap_then_round->set_priority(long_edges_first);
1528
1529
1530 // swap_then_round->add_invariant(inversion_invariant);
1531 for (auto& s : update_child_position) {
1532 swap_then_round->add_transfer_strategy(s);
1533 }
1534
1535 ops.push_back(swap_then_round);
1536 ops_name.push_back("EDGE SWAP");
1537 ops.emplace_back(rounding);
1538 ops_name.emplace_back("rounding");
1539 }
1540
1541 // 4) Smoothing
1542 // //////////////////////////////////////
1543 // // special smoothing on surface
1544 // //////////////////////////////////////
1545
1546 // if (mesh->top_simplex_type() == PrimitiveType::Tetrahedron) {
1547 // for (auto& pair : mesh_constraint_pairs) {
1548 // if (pair.second.mesh().top_simplex_type() != PrimitiveType::Triangle) continue;
1549
1550 // auto& child_mesh = pair.second.mesh();
1551 // auto& child_position_handle = pair.second;
1552
1553 // auto lap_smoothing = std::make_shared<TetWildTangentialLaplacianSmoothing>(
1554 // child_mesh,
1555 // child_position_handle.as<Rational>());
1556 // lap_smoothing->add_invariant(std::make_shared<RoundedInvariant>(
1557 // child_mesh,
1558 // child_position_handle.as<Rational>()));
1559 // lap_smoothing->add_invariant(inversion_invariant);
1560
1561 // auto proj_lap_smoothing =
1562 // std::make_shared<ProjectOperation>(lap_smoothing, mesh_constraint_pairs);
1563 // proj_lap_smoothing->use_random_priority() = true;
1564
1565 // proj_lap_smoothing->add_invariant(envelope_invariant);
1566 // proj_lap_smoothing->add_invariant(inversion_invariant);
1567 // proj_lap_smoothing->add_invariant(
1568 // std::make_shared<EnergyFilterInvariant>(*mesh,
1569 // energy_filter_attribute.as<char>()));
1570
1571 // proj_lap_smoothing->add_transfer_strategy(amips_update);
1572 // proj_lap_smoothing->add_transfer_strategy(edge_length_update);
1573 // for (auto& s : update_parent_position) { // TODO::this should from only one child
1574 // proj_lap_smoothing->add_transfer_strategy(s);
1575 // }
1576 // for (auto& s : update_child_position) {
1577 // proj_lap_smoothing->add_transfer_strategy(s);
1578 // }
1579
1580 // ops.push_back(proj_lap_smoothing);
1581 // ops_name.push_back("LAPLACIAN SMOOTHING");
1582 // }
1583 // }
1584
1585 // auto energy =
1586 // std::make_shared<function::LocalNeighborsSumFunction>(*mesh, pt_attribute, *amips);
1587 // auto smoothing = std::make_shared<OptimizationSmoothing>(energy);
1588 auto smoothing = std::make_shared<AMIPSOptimizationSmoothing>(*mesh, pt_attribute);
1589 smoothing->add_invariant(
1590 std::make_shared<RoundedInvariant>(*mesh, pt_attribute.as<Rational>()));
1591 smoothing->add_invariant(inversion_invariant);
1592 for (auto& s : update_child_position) {
1593 smoothing->add_transfer_strategy(s);
1594 }
1595
1596 // test code
1597 // smoothing->add_invariant(envelope_invariant);
1598 // smoothing->add_transfer_strategy(amips_update);
1599 // smoothing->add_transfer_strategy(edge_length_update);
1600 // ops.push_back(smoothing);
1601 // ops_name.push_back("SMOOTHING");
1602 // ops.emplace_back(rounding);
1603 // ops_name.emplace_back("rounding");
1604 //--------
1605
1606 auto proj_smoothing = std::make_shared<ProjectOperation>(smoothing, mesh_constraint_pairs);
1607 proj_smoothing->use_random_priority() = true;
1608
1609 proj_smoothing->add_invariant(envelope_invariant);
1610 proj_smoothing->add_invariant(inversion_invariant);
1611 proj_smoothing->add_invariant(
1612 std::make_shared<EnergyFilterInvariant>(*mesh, energy_filter_attribute.as<char>()));
1613
1614 proj_smoothing->add_transfer_strategy(amips_update);
1615 proj_smoothing->add_transfer_strategy(edge_length_update);
1616 for (auto& s : update_parent_position) {
1617 proj_smoothing->add_transfer_strategy(s);
1618 }
1619
1620 for (auto& s : update_child_position) {
1621 proj_smoothing->add_transfer_strategy(s);
1622 }
1623 // proj_smoothing->add_transfer_strategy(target_edge_length_update);
1624
1625 for (int i = 0; i < 1; ++i) {
1626 ops.push_back(proj_smoothing);
1627 ops_name.push_back("SMOOTHING");
1628 }
1629
1630 ops.emplace_back(rounding);
1631 ops_name.emplace_back("rounding");
1632
1633 write(
1634 mesh,
1635 paths.output_dir,
1636 options.output,
1637 options.attributes.position,
1638 0,
1639 options.intermediate_output);
1640
1641 auto surface_mesh = mesh->get_child_meshes().front();
1642
1643 write(
1644 surface_mesh,
1645 paths.output_dir,
1646 options.output + "_intermediate_surface",
1647 options.attributes.position,
1648 0,
1649 options.intermediate_output);
1650
1652 // Running all ops in order n times
1653 Scheduler scheduler;
1654 {
1655 const size_t freq = options.scheduler_update_frequency;
1656 scheduler.set_update_frequency(freq == 0 ? std::optional<size_t>{} : freq);
1657 }
1658 int64_t success = 10;
1659
1661 // preprocessing
1663
1664 SchedulerStats pre_stats;
1665
1666 logger().info("----------------------- Preprocess Collapse -----------------------");
1667 for (const auto& t : mesh->get_all(mesh->top_simplex_type())) {
1668 const auto vertices = mesh->orient_vertices(t);
1669 std::vector<Vector3r> pos;
1670 for (int i = 0; i < vertices.size(); ++i) {
1671 pos.push_back(pt_accessor.const_vector_attribute(vertices[i]));
1672 }
1673 if (wmtk::utils::wmtk_orient3d(pos[0], pos[1], pos[2], pos[3]) <= 0) {
1674 wmtk::logger().error("Flipped tet!");
1675 }
1676 }
1677 logger().info("Executing collapse ...");
1678
1679
1680 wmtk::attribute::TypedAttributeHandle<char> visited_edge_flag_t = visited_edge_flag.as<char>();
1681
1682 pre_stats = scheduler.run_operation_on_all(*collapse_then_round, visited_edge_flag_t);
1683 logger().info(
1684 "Executed {}, {} ops (S/F) {}/{}. Time: collecting: {}, sorting: {}, "
1685 "executing: {}",
1686 "preprocessing collapse",
1689 pre_stats.number_of_failed_operations(),
1690 pre_stats.collecting_time,
1691 pre_stats.sorting_time,
1692 pre_stats.executing_time);
1693
1694 success = pre_stats.number_of_successful_operations();
1695
1696 // verbose logger, can be removed
1697 int64_t unrounded = 0;
1698 for (const auto& v : mesh->get_all(PrimitiveType::Vertex)) {
1699 const auto p = pt_accessor.vector_attribute(v);
1700 for (int64_t d = 0; d < 3; ++d) {
1701 if (!p[d].is_rounded()) {
1702 ++unrounded;
1703 break;
1704 }
1705 }
1706 }
1707
1708 logger().info("Mesh has {} unrounded vertices", unrounded);
1709
1710 // compute max energy
1711 double max_energy = std::numeric_limits<double>::lowest();
1712 double min_energy = std::numeric_limits<double>::max();
1713 double avg_energy = 0;
1714 for (const auto& t : mesh->get_all(mesh->top_simplex_type())) {
1715 // double e = amips->get_value(simplex::Simplex(mesh->top_simplex_type(), t));
1716 double e = amips_accessor.scalar_attribute(t);
1717 max_energy = std::max(max_energy, e);
1718 min_energy = std::min(min_energy, e);
1719 avg_energy += e;
1720 }
1721
1722 avg_energy = avg_energy / mesh->get_all(mesh->top_simplex_type()).size();
1723
1724 logger().info(
1725 "Max AMIPS Energy: {}, Min AMIPS Energy: {}, Avg AMIPS Energy: {}",
1726 max_energy,
1727 min_energy,
1728 avg_energy);
1729
1730 // std::ofstream file0("quality_plot_pre.csv");
1731 // file0 << "tid, quality" << std::endl;
1732 // int64_t t_cnt = 0;
1733 // for (const auto& t : mesh->get_all(PrimitiveType::Tetrahedron)) {
1734 // t_cnt++;
1735 // file0 << t_cnt << ", " << amips_accessor.scalar_attribute(t) << std::endl;
1736 // }
1737
1738
1739 double old_max_energy = max_energy;
1740 double old_avg_energy = avg_energy;
1741 int iii = 0;
1742 bool is_double = false;
1743 for (int64_t i = 0; i < options.passes; ++i) {
1744 logger().info("--------------------------- Pass {} ---------------------------", i);
1745
1746 SchedulerStats pass_stats;
1747 int jj = 0;
1748 for (auto& op : ops) {
1749 for (const auto& t : mesh->get_all(mesh->top_simplex_type())) {
1750 const auto vertices = mesh->orient_vertices(t);
1751 std::vector<Vector3r> pos;
1752 for (int i = 0; i < vertices.size(); ++i) {
1753 pos.push_back(pt_accessor.const_vector_attribute(vertices[i]));
1754 }
1755 if (wmtk::utils::wmtk_orient3d(pos[0], pos[1], pos[2], pos[3]) <= 0) {
1756 wmtk::logger().error("Flipped tet!");
1757 }
1758 }
1759 logger().info("Executing {} ...", ops_name[jj]);
1760 SchedulerStats stats;
1761 if (op->primitive_type() == PrimitiveType::Edge) {
1762 stats = scheduler.run_operation_on_all(*op, visited_edge_flag_t);
1763 // } else if (ops_name[jj] == "SMOOTHING") {
1764 // // stats = scheduler.run_operation_on_all(*op);
1765 // stats =
1766 // scheduler.run_operation_on_all_coloring(*op,
1767 // coloring_attribute.as<int64_t>());
1768 } else {
1769 stats = scheduler.run_operation_on_all(*op);
1770 }
1771 pass_stats += stats;
1772 logger().info(
1773 "Executed {}, {} ops (S/F) {}/{}. Time: collecting: {}, sorting: {}, "
1774 "executing: {}",
1775 ops_name[jj],
1779 stats.collecting_time,
1780 stats.sorting_time,
1781 stats.executing_time);
1782
1783 success = stats.number_of_successful_operations();
1784
1785 // verbose logger, can be removed
1786 int64_t unrounded = 0;
1787 for (const auto& v : mesh->get_all(PrimitiveType::Vertex)) {
1788 const auto p = pt_accessor.vector_attribute(v);
1789 for (int64_t d = 0; d < 3; ++d) {
1790 if (!p[d].is_rounded()) {
1791 ++unrounded;
1792 break;
1793 }
1794 }
1795 }
1796
1797 logger().info("Mesh has {} unrounded vertices", unrounded);
1798
1799 avg_energy = 0;
1800
1801 // compute max energy
1802 max_energy = std::numeric_limits<double>::lowest();
1803 min_energy = std::numeric_limits<double>::max();
1804 for (const auto& t : mesh->get_all(mesh->top_simplex_type())) {
1805 // double e = amips->get_value(simplex::Simplex(mesh->top_simplex_type(), t));
1806 double e = amips_accessor.scalar_attribute(t);
1807 max_energy = std::max(max_energy, e);
1808 min_energy = std::min(min_energy, e);
1809 avg_energy += e;
1810 }
1811
1812 avg_energy = avg_energy / mesh->get_all(mesh->top_simplex_type()).size();
1813
1814 logger().info(
1815 "Max AMIPS Energy: {}, Min AMIPS Energy: {}, Avg AMIPS Energy: {}",
1816 max_energy,
1817 min_energy,
1818 avg_energy);
1819
1820
1821 ++jj;
1822 }
1823
1824 logger().info(
1825 "Executed {} ops (S/F) {}/{}. Time: collecting: {}, sorting: {}, executing: {}",
1826 pass_stats.number_of_performed_operations(),
1828 pass_stats.number_of_failed_operations(),
1829 pass_stats.collecting_time,
1830 pass_stats.sorting_time,
1831 pass_stats.executing_time);
1832
1834
1835 write(
1836 mesh,
1837 paths.output_dir,
1838 options.output,
1839 options.attributes.position,
1840 i + 1,
1841 options.intermediate_output);
1842
1843 write(
1844 surface_mesh,
1845 paths.output_dir,
1846 options.output + "_intermediate_surface",
1847 options.attributes.position,
1848 i + 1,
1849 options.intermediate_output);
1850
1851 assert(mesh->is_connectivity_valid());
1852
1853 // compute max energy
1854 max_energy = std::numeric_limits<double>::lowest();
1855 min_energy = std::numeric_limits<double>::max();
1856 avg_energy = 0;
1857 for (const auto& t : mesh->get_all(mesh->top_simplex_type())) {
1858 // double e = amips->get_value(simplex::Simplex(mesh->top_simplex_type(), t));
1859 double e = amips_accessor.scalar_attribute(t);
1860 max_energy = std::max(max_energy, e);
1861 min_energy = std::min(min_energy, e);
1862 avg_energy += e;
1863 }
1864
1865 avg_energy = avg_energy / mesh->get_all(mesh->top_simplex_type()).size();
1866
1867 int64_t unrounded = 0;
1868 if (!is_double) {
1869 bool rational = false;
1870 for (const auto& v : mesh->get_all(PrimitiveType::Vertex)) {
1871 const auto p = pt_accessor.vector_attribute(v);
1872 for (int64_t d = 0; d < bmax.size(); ++d) {
1873 if (!p[d].is_rounded()) {
1874 rational = true;
1875 ++unrounded;
1876 break;
1877 }
1878 }
1879 }
1880
1881 is_double = !rational;
1882 }
1883
1884 logger().info("Mesh has {} unrounded vertices", unrounded);
1885 logger().info(
1886 "Max AMIPS Energy: {}, Min AMIPS Energy: {}, Avg AMIPS Energy: {}",
1887 max_energy,
1888 min_energy,
1889 avg_energy);
1890
1891
1892 // adjust sizing field
1893 if (i > 0 && old_max_energy - max_energy < 5e-1 &&
1894 (old_avg_energy - avg_energy) / avg_energy < 0.1) {
1895 wmtk::logger().info("adjusting sizing field ...");
1896
1898 *mesh,
1899 pt_attribute.as<Rational>(),
1900 edge_length_attribute.as<double>(),
1901 sizing_field_scalar_attribute.as<double>(),
1902 amips_attribute.as<double>(),
1903 target_edge_length_attribute.as<double>(),
1904 visited_vertex_flag.as<char>(),
1905 target_max_amips,
1906 max_energy,
1907 target_edge_length,
1908 min_edge_length);
1909
1910 wmtk::logger().info("adjusting sizing field finished");
1911
1912 // wmtk::logger().info("setting energy filter ...");
1913 // set_operation_energy_filter_after_sizing_field(
1914 // *mesh,
1915 // pt_attribute.as<Rational>(),
1916 // amips_attribute.as<double>(),
1917 // energy_filter_attribute.as<char>(),
1918 // visited_vertex_flag.as<char>(),
1919 // target_max_amips,
1920 // max_energy,
1921 // target_edge_length);
1922 // wmtk::logger().info("setting energy filter finished");
1923
1924 // int64_t e_cnt = 0;
1925 // for (const auto& e : mesh->get_all(PrimitiveType::Edge)) {
1926 // if (energy_filter_accessor.scalar_attribute(e) == char(1) ||
1927 // energy_filter_accessor.scalar_attribute(
1928 // mesh->switch_tuple(e, PrimitiveType::Vertex)) == char(1)) {
1929 // e_cnt++;
1930 // }
1931 // }
1932 // wmtk::logger().info(
1933 // "{} edges are going to be executed out of {}",
1934 // e_cnt,
1935 // mesh->get_all(PrimitiveType::Edge).size());
1936
1937 for (const auto& v : mesh->get_all(PrimitiveType::Vertex)) {
1938 energy_filter_accessor.scalar_attribute(v) = char(1);
1939 }
1940 wmtk::logger().info("reset energy filter");
1941 } else {
1942 wmtk::logger().info("setting energy filter ...");
1944 *mesh,
1945 pt_attribute.as<Rational>(),
1946 amips_attribute.as<double>(),
1947 energy_filter_attribute.as<char>(),
1948 visited_vertex_flag.as<char>(),
1949 target_max_amips,
1950 max_energy,
1951 target_edge_length);
1952 wmtk::logger().info("setting energy filter finished");
1953
1954 int64_t e_cnt = 0;
1955 for (const auto& e : mesh->get_all(PrimitiveType::Edge)) {
1956 if (energy_filter_accessor.scalar_attribute(e) == char(1) ||
1957 energy_filter_accessor.scalar_attribute(
1958 mesh->switch_tuple(e, PrimitiveType::Vertex)) == char(1)) {
1959 e_cnt++;
1960 }
1961 }
1962 wmtk::logger().info(
1963 "{} edges are going to be executed out of {}",
1964 e_cnt,
1965 mesh->get_all(PrimitiveType::Edge).size());
1966 }
1967
1968 old_max_energy = max_energy;
1969 old_avg_energy = avg_energy;
1970
1971 // stop at good quality
1972 if (max_energy <= target_max_amips && is_double) break;
1973 }
1974
1975 // output
1976 cache.write_mesh(*mesh, options.output);
1977}
1978} // namespace wmtk::components
const attribute::Accessor< T, Mesh, D > create_const_accessor(const attribute::MeshAttributeHandle &handle) const
attribute::Accessor< T, Mesh, D > create_accessor(const attribute::MeshAttributeHandle &handle)
std::vector< Tuple > get_all(PrimitiveType type) const
Generate a vector of Tuples from global vertex/edge/triangle/tetrahedron index.
Definition Mesh.cpp:18
virtual Tuple switch_tuple(const Tuple &tuple, PrimitiveType type) const =0
switch the orientation of the Tuple of the given dimension
virtual std::vector< Tuple > orient_vertices(const Tuple &t) const =0
PrimitiveType top_simplex_type() const
Definition Mesh.hpp:982
SchedulerStats run_operation_on_all(operations::Operation &op)
Definition Scheduler.cpp:33
void set_update_frequency(std::optional< size_t > &&freq={})
int64_t number_of_failed_operations() const
Returns the number of failed operations performed by the scheduler.
Definition Scheduler.hpp:23
int64_t number_of_performed_operations() const
Returns the number of performed operations performed by the scheduler.
Definition Scheduler.hpp:30
int64_t number_of_successful_operations() const
Returns the number of successful operations performed by the scheduler.
Definition Scheduler.hpp:16
The Tuple is the basic navigation tool in our mesh data structure.
Definition Tuple.hpp:19
Handle that represents attributes for some mesh.
void write_mesh(const Mesh &m, const std::string &name, const std::map< std::string, std::vector< int64_t > > &multimesh_names={})
Write a mesh to cache.
Definition Cache.cpp:194
std::shared_ptr< Mesh > read_mesh(const std::string &name) const
Load a mesh from cache.
Definition Cache.cpp:171
void add_invariant(std::shared_ptr< Invariant > invariant)
Definition Operation.hpp:48
void set_priority(const std::function< double(const simplex::Simplex &)> &func)
Definition Operation.hpp:50
virtual PrimitiveType primitive_type() const =0
void add_transfer_strategy(const std::shared_ptr< const operations::AttributeTransferStrategyBase > &other)
Definition Operation.cpp:59
std::pair< attribute::MeshAttributeHandle, attribute::MeshAttributeHandle > MeshConstrainPair
const std::vector< Simplex > & simplex_vector() const
Return const reference to the simplex vector.
static Simplex vertex(const Mesh &m, const Tuple &t)
Definition Simplex.hpp:56
constexpr wmtk::PrimitiveType PT
constexpr wmtk::PrimitiveType PF
void set_operation_energy_filter(Mesh &m, const TypedAttributeHandle< double > &coordinate_handle, const TypedAttributeHandle< double > &energy_handle, const TypedAttributeHandle< char > &energy_filter_handle, const TypedAttributeHandle< char > &visited_handle, const double stop_energy, const double current_max_energy, const double initial_target_edge_length)
void adjust_sizing_field(Mesh &m, const TypedAttributeHandle< double > &coordinate_handle, const TypedAttributeHandle< double > &edge_length_handle, const TypedAttributeHandle< double > &sizing_field_scalar_handle, const TypedAttributeHandle< double > &energy_handle, const TypedAttributeHandle< double > &target_edge_length_handle, const TypedAttributeHandle< char > &visited_handle, const double stop_energy, const double current_max_energy, const double initial_target_edge_length, const double min_target_edge_length)
void set_operation_energy_filter_after_sizing_field(Mesh &m, const TypedAttributeHandle< Rational > &coordinate_handle, const TypedAttributeHandle< double > &energy_handle, const TypedAttributeHandle< char > &energy_filter_handle, const TypedAttributeHandle< char > &visited_handle, const double stop_energy, const double current_max_energy, const double initial_target_edge_length)
void write(const Mesh &mesh, const std::string &out_dir, const std::string &name, const std::string &vname, const int64_t index, const bool intermediate_output)
std::vector< attribute::MeshAttributeHandle > get_attributes(const Mesh &m, const nlohmann::json &attribute_names)
std::vector< std::pair< std::shared_ptr< Mesh >, std::string > > wildmeshing(const WildMeshingOptions &option)
double Tet_AMIPS_energy(const std::array< double, 12 > &T)
Definition amips.cpp:380
double Tet_AMIPS_energy(const std::array< double, 12 > &T)
Definition AMIPS.cpp:375
double Tri_AMIPS_energy(const std::array< double, 6 > &T)
Definition AMIPS.cpp:473
void consolidate(Mesh &mesh)
std::shared_ptr< AttributeTransferStrategyBase > make_cast_attribute_transfer_strategy(const wmtk::attribute::MeshAttributeHandle &source, const wmtk::attribute::MeshAttributeHandle &target)
std::vector< Tuple > vertices(const Mesh &m, const Simplex &simplex)
bool link_condition(const EdgeMesh &mesh, const Tuple &edge)
Check if the edge to collapse satisfying the link condition.
SimplexCollection link(const Mesh &mesh, const simplex::Simplex &simplex, const bool sort_and_clean)
Definition link.cpp:84
SimplexCollection k_ring(const Mesh &mesh, const Simplex &simplex, int64_t k)
Definition k_ring.cpp:6
int wmtk_orient3d(const Eigen::Ref< const Eigen::Vector3< Rational > > &p0, const Eigen::Ref< const Eigen::Vector3< Rational > > &p1, const Eigen::Ref< const Eigen::Vector3< Rational > > &p2, const Eigen::Ref< const Eigen::Vector3< Rational > > &p3)
Definition orient.cpp:75
constexpr PrimitiveType PE
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:58
Vector< double, 3 > Vector3d
Definition Types.hpp:39
constexpr PrimitiveType PV
std::vector< nlohmann::json > pass_through
std::vector< WildmeshingOptionsEnvelope > envelopes
const std::filesystem::path data_dir