PolyFEM
Loading...
Searching...
No Matches
WildRemesher.cpp
Go to the documentation of this file.
1#include "WildRemesher.hpp"
2
6
7#include <wmtk/utils/TupleUtils.hpp>
8
9#include <unordered_map>
10
11#define VERTEX_ATTRIBUTE_GETTER(name, attribute) \
12 template <class WMTKMesh> \
13 Eigen::MatrixXd WildRemesher<WMTKMesh>::name() const \
14 { \
15 Eigen::MatrixXd attributes = Eigen::MatrixXd::Constant(WMTKMesh::vert_capacity(), DIM, NaN); \
16 for (const Tuple &t : WMTKMesh::get_vertices()) \
17 attributes.row(t.vid(*this)) = vertex_attrs[t.vid(*this)].attribute; \
18 return attributes; \
19 }
20
21#define VERTEX_ATTRIBUTE_SETTER(name, attribute) \
22 template <class WMTKMesh> \
23 void WildRemesher<WMTKMesh>::name(const Eigen::MatrixXd &attributes) \
24 { \
25 for (const Tuple &t : WMTKMesh::get_vertices()) \
26 vertex_attrs[t.vid(*this)].attribute = attributes.row(t.vid(*this)); \
27 }
28
29namespace polyfem::mesh
30{
31 static constexpr double NaN = std::numeric_limits<double>::quiet_NaN();
32
33 template <class WMTKMesh>
35 const State &state,
36 const Eigen::MatrixXd &obstacle_displacements,
37 const Eigen::MatrixXd &obstacle_vals,
38 const double current_time,
39 const double starting_energy)
40 : Remesher(state, obstacle_displacements, obstacle_vals, current_time, starting_energy),
41 WMTKMesh()
42 {
43 }
44
45 template <typename WMTKMesh>
47 rank_edges(const Remesher::EdgeMap<double> &edge_energy, const json &args)
48 {
49 if (edge_energy.empty())
51
52 double min_energy = std::numeric_limits<double>::infinity();
53 double max_energy = -std::numeric_limits<double>::infinity();
54 std::vector<double> sorted_energies;
55 sorted_energies.reserve(edge_energy.size());
56 for (const auto &[edge, energy] : edge_energy)
57 {
58 min_energy = std::min(min_energy, energy);
59 max_energy = std::max(max_energy, energy);
60 sorted_energies.push_back(energy);
61 }
62 std::sort(sorted_energies.begin(), sorted_energies.end());
63
64 const double split_threshold = args["split"]["culling_threshold"];
65 const double split_tolerance = args["split"]["acceptance_tolerance"];
66 const double top_energy_threshold = (max_energy - min_energy) * split_threshold + min_energy;
67 const double top_element_threshold = sorted_energies[int(sorted_energies.size() * split_threshold)];
68 const double top_threshold = std::max(std::min(top_energy_threshold, top_element_threshold), std::min(1e-12, split_tolerance));
69
70 const double collapse_threshold = args["collapse"]["culling_threshold"];
71 const double bottom_energy_threshold = (max_energy - min_energy) * collapse_threshold + min_energy;
72 // const double bottom_element_threshold = sorted_energies[int(sorted_energies.size() * collapse_threshold)];
73 const double bottom_threshold = bottom_energy_threshold;
75 logger().info("min energy: {}, max energy: {}, thresholds: {}, {}", min_energy, max_energy, bottom_threshold, top_threshold);
76
78 for (const auto &[edge, energy] : edge_energy)
79 {
80 if (energy >= top_threshold)
82 else if (energy <= bottom_threshold)
84 else
86 }
87
88 return edge_ranks;
89 }
91 template <class WMTKMesh>
93 const Eigen::MatrixXd &rest_positions,
94 const Eigen::MatrixXd &positions,
95 const Eigen::MatrixXi &elements,
96 const Eigen::MatrixXd &projection_quantities,
97 const BoundaryMap<int> &boundary_to_id,
98 const std::vector<int> &body_ids,
99 const EdgeMap<double> &elastic_energy,
100 const EdgeMap<double> &contact_energy)
101 {
103 rest_positions, positions, elements, projection_quantities,
104 boundary_to_id, body_ids, elastic_energy, contact_energy);
105
106 total_volume = 0;
107 for (const Tuple &t : get_elements())
108 total_volume += element_volume(t);
109 assert(total_volume > 0);
110
111#ifndef NDEBUG
112 assert(get_elements().size() == elements.rows());
113 for (const Tuple &t : get_elements())
114 assert(!is_inverted(t));
115#endif
116
117 const auto edge_elastic_ranks = rank_edges<WMTKMesh>(elastic_energy, args);
118 const auto edge_contact_ranks = rank_edges<WMTKMesh>(contact_energy, args);
119
121 for (const auto &[edge, elastic_rank] : edge_elastic_ranks)
122 {
123 const auto contact_rank = edge_contact_ranks.empty() ? elastic_rank : edge_contact_ranks.at(edge);
124 if (elastic_rank == EdgeAttributes::EnergyRank::TOP || contact_rank == EdgeAttributes::EnergyRank::TOP)
125 edge_ranks[edge] = EdgeAttributes::EnergyRank::TOP;
126 else if (elastic_rank == EdgeAttributes::EnergyRank::BOTTOM && contact_rank == EdgeAttributes::EnergyRank::BOTTOM)
127 edge_ranks[edge] = EdgeAttributes::EnergyRank::BOTTOM;
128 else
129 edge_ranks[edge] = EdgeAttributes::EnergyRank::MIDDLE;
130 }
131
132 for (const Tuple &edge : WMTKMesh::get_edges())
133 {
134 const size_t e0 = edge.vid(*this);
135 const size_t e1 = edge.switch_vertex(*this).vid(*this);
136 edge_attr(edge.eid(*this)).energy_rank = edge_ranks.at({{e0, e1}});
137 }
138
139 // write_edge_ranks_mesh(edge_elastic_ranks, edge_contact_ranks);
140 }
142 // -------------------------------------------------------------------------
143 // Getters
144
145 VERTEX_ATTRIBUTE_GETTER(rest_positions, rest_position)
146 VERTEX_ATTRIBUTE_GETTER(positions, position)
147 VERTEX_ATTRIBUTE_GETTER(displacements, displacement())
148
149 template <class WMTKMesh>
150 Eigen::MatrixXi WildRemesher<WMTKMesh>::edges() const
151 {
152 const std::vector<Tuple> edges = WMTKMesh::get_edges();
153 Eigen::MatrixXi E = Eigen::MatrixXi::Constant(edges.size(), 2, -1);
154 for (int i = 0; i < edges.size(); ++i)
155 {
156 const Tuple &e = edges[i];
157 E(i, 0) = e.vid(*this);
158 E(i, 1) = e.switch_vertex(*this).vid(*this);
159 }
160 return E;
162
163 template <class WMTKMesh>
164 Eigen::MatrixXi WildRemesher<WMTKMesh>::elements() const
165 {
166 const std::vector<Tuple> elements = get_elements();
167 Eigen::MatrixXi F = Eigen::MatrixXi::Constant(elements.size(), dim() + 1, -1);
168 for (size_t i = 0; i < elements.size(); i++)
170 const auto vids = element_vids(elements[i]);
171
172 for (int j = 0; j < vids.size(); j++)
174 F(i, j) = vids[j];
176 }
177 return F;
178 }
179
180 template <class WMTKMesh>
182 {
183 Eigen::MatrixXd projection_quantities =
184 Eigen::MatrixXd::Constant(dim() * WMTKMesh::vert_capacity(), n_quantities(), NaN);
185
186 for (const Tuple &t : WMTKMesh::get_vertices())
187 {
188 const int vi = t.vid(*this);
189 projection_quantities.middleRows(dim() * vi, dim()) = vertex_attrs[vi].projection_quantities;
190 }
191
192 return projection_quantities;
194
195 template <class WMTKMesh>
196 std::vector<int> WildRemesher<WMTKMesh>::body_ids() const
198 const std::vector<Tuple> elements = get_elements();
199 std::vector<int> body_ids(elements.size(), -1);
200 for (size_t i = 0; i < elements.size(); i++)
201 {
202 body_ids[i] = element_attrs[element_id(elements[i])].body_id;
203 }
204 return body_ids;
207 template <class WMTKMesh>
209 const Eigen::VectorXi &vertex_to_basis) const
211 std::vector<int> boundary_nodes;
213 const std::unordered_map<int, std::array<bool, 3>> bcs =
214 state.boundary_conditions_ids("dirichlet_boundary");
215 std::vector<int> boundary_ids;
216 const std::vector<Tuple> boundary_facets = this->boundary_facets(&boundary_ids);
217 for (int i = 0; i < boundary_facets.size(); ++i)
218 {
219 const Tuple &t = boundary_facets[i];
220 const auto bc = bcs.find(boundary_ids[i]);
221
222 if (bc == bcs.end())
223 continue;
224
225 for (int d = 0; d < this->dim(); ++d)
226 {
227 if (!bc->second[d])
228 continue;
229
230 for (const size_t vid : facet_vids(t))
231 boundary_nodes.push_back(dim() * vertex_to_basis[vid] + d);
232 }
233 }
234
235 // Sort and remove the duplicate boundary_nodes.
236 std::sort(boundary_nodes.begin(), boundary_nodes.end());
237 auto new_end = std::unique(boundary_nodes.begin(), boundary_nodes.end());
238 boundary_nodes.erase(new_end, boundary_nodes.end());
239
240 return boundary_nodes;
241 }
242
243 template <class WMTKMesh>
244 std::vector<typename WMTKMesh::Tuple>
245 WildRemesher<WMTKMesh>::boundary_facets(std::vector<int> *boundary_ids) const
246 {
247 POLYFEM_REMESHER_SCOPED_TIMER("boundary_facets");
248
249 size_t element_capacity;
250 if constexpr (std::is_same_v<wmtk::TriMesh, WMTKMesh>)
251 element_capacity = WMTKMesh::tri_capacity();
252 else
253 element_capacity = WMTKMesh::tet_capacity();
254
255 const size_t max_tid = std::min(boundary_attrs.size() / FACETS_PER_ELEMENT, element_capacity);
256
257 std::vector<Tuple> boundary_facets;
258 for (int tid = 0; tid < max_tid; ++tid)
259 {
260 const Tuple t = tuple_from_element(tid);
261 if (!t.is_valid(*this))
262 continue;
263
264 for (int local_fid = 0; local_fid < FACETS_PER_ELEMENT; ++local_fid)
265 {
266 const Tuple facet_tuple = tuple_from_facet(tid, local_fid);
267 const int boundary_id = boundary_attrs[facet_id(facet_tuple)].boundary_id;
268 assert((boundary_id >= 0) == is_boundary_facet(facet_tuple));
269 if (boundary_id >= 0)
270 {
271 boundary_facets.push_back(facet_tuple);
272 if (boundary_ids)
273 boundary_ids->push_back(boundary_id);
274 }
275 }
277
278#ifndef NDEBUG
279 {
280 std::vector<size_t> boundary_facet_ids;
281 for (const Tuple &f : boundary_facets)
282 boundary_facet_ids.push_back(facet_id(f));
283 std::sort(boundary_facet_ids.begin(), boundary_facet_ids.end());
284
285 std::vector<size_t> gt_boundary_facet_ids;
286 for (const Tuple &f : get_facets())
287 if (is_boundary_facet(f))
288 gt_boundary_facet_ids.push_back(facet_id(f));
289 std::sort(gt_boundary_facet_ids.begin(), gt_boundary_facet_ids.end());
290
291 assert(boundary_facets.size() == gt_boundary_facet_ids.size());
292 assert(boundary_facet_ids == gt_boundary_facet_ids);
293 }
294#endif
295
296 return boundary_facets;
297 }
298
299 // -------------------------------------------------------------------------
300 // Setters
301
302 VERTEX_ATTRIBUTE_SETTER(set_rest_positions, rest_position)
303 VERTEX_ATTRIBUTE_SETTER(set_positions, position)
304
305 template <class WMTKMesh>
307 const Eigen::MatrixXd &projection_quantities)
308 {
309 assert(projection_quantities.rows() == dim() * WMTKMesh::vert_capacity());
310 m_n_quantities = projection_quantities.cols();
311
312 for (const Tuple &t : WMTKMesh::get_vertices())
313 {
314 const int vi = t.vid(*this);
315 vertex_attrs[vi].projection_quantities =
316 projection_quantities.middleRows(dim() * vi, dim());
317 }
318 }
319
320 template <class WMTKMesh>
322 const std::vector<bool> &fixed)
323 {
324 assert(fixed.size() == WMTKMesh::vert_capacity());
325 for (const Tuple &t : WMTKMesh::get_vertices())
326 vertex_attrs[t.vid(*this)].fixed = fixed[t.vid(*this)];
327 }
328
329 template <class WMTKMesh>
331 const std::vector<int> &body_ids)
332 {
333 const std::vector<Tuple> elements = get_elements();
334 for (int i = 0; i < elements.size(); ++i)
335 {
336 element_attrs[element_id(elements[i])].body_id = body_ids.at(i);
337 }
338 }
339
340 template <class WMTKMesh>
342 {
343 for (const Tuple &face : get_facets())
344 {
345 if constexpr (std::is_same_v<wmtk::TriMesh, WMTKMesh>)
346 {
347 assert(std::holds_alternative<EdgeMap<int>>(boundary_to_id));
348 const EdgeMap<int> &edge_to_boundary_id = std::get<EdgeMap<int>>(boundary_to_id);
349 boundary_attrs[facet_id(face)].boundary_id = edge_to_boundary_id.at(facet_vids(face));
350 }
351 else
352 {
353 assert(std::holds_alternative<FaceMap<int>>(boundary_to_id));
354 const FaceMap<int> &face_to_boundary_id = std::get<FaceMap<int>>(boundary_to_id);
355 boundary_attrs[facet_id(face)].boundary_id = face_to_boundary_id.at(facet_vids(face));
356 }
357
358#ifndef NDEBUG
359 if (is_boundary_facet(face))
360 assert(boundary_attrs[facet_id(face)].boundary_id >= 0);
361 else
362 assert(boundary_attrs[facet_id(face)].boundary_id == -1);
363#endif
364 }
365 }
366
367 // -------------------------------------------------------------------------
368
369 template <class WMTKMesh>
370 bool WildRemesher<WMTKMesh>::invariants(const std::vector<Tuple> &new_tris)
371 {
372 POLYFEM_REMESHER_SCOPED_TIMER("WildRemesher::invariants");
373 // for (auto &t : new_tris)
374 for (auto &t : get_elements())
375 {
376 if (is_inverted(t))
377 {
378 static int inversion_cnt = 0;
379 write_mesh(state.resolve_output_path(fmt::format("inversion_{:04d}.vtu", inversion_cnt++)));
380 log_and_throw_error("Inverted element found, invariants violated!");
381 return false;
382 }
383 }
384 return true;
385 }
386
387 // -------------------------------------------------------------------------
388 // Utils
389
390 template <class WMTKMesh>
392 {
393 const auto &e0 = vertex_attrs[e.vid(*this)].rest_position;
394 const auto &e1 = vertex_attrs[e.switch_vertex(*this).vid(*this)].rest_position;
395 return (e1 - e0).norm();
396 }
397
398 template <class WMTKMesh>
400 {
401 const auto &e0 = vertex_attrs[e.vid(*this)].position;
402 const auto &e1 = vertex_attrs[e.switch_vertex(*this).vid(*this)].position;
403 return (e1 - e0).norm();
404 }
405
406 template <class WMTKMesh>
409 {
410 const VectorNd &e0 = vertex_attrs[e.vid(*this)].rest_position;
411 const VectorNd &e1 = vertex_attrs[e.switch_vertex(*this).vid(*this)].rest_position;
412 return (e1 + e0) / 2.0;
413 }
414
415 template <class WMTKMesh>
418 {
419 const VectorNd &e0 = vertex_attrs[e.vid(*this)].position;
420 const VectorNd &e1 = vertex_attrs[e.switch_vertex(*this).vid(*this)].position;
421 return (e1 + e0) / 2.0;
422 }
423
424 template <class WMTKMesh>
425 std::vector<typename WMTKMesh::Tuple> WildRemesher<WMTKMesh>::get_edges_for_elements(
426 const std::vector<Tuple> &elements) const
427 {
428 std::vector<Tuple> edges;
429 for (auto t : elements)
430 for (int j = 0; j < EDGES_PER_ELEMENT; ++j)
431 edges.push_back(WMTKMesh::tuple_from_edge(element_id(t), j));
432 wmtk::unique_edge_tuples(*this, edges);
433 return edges;
434 }
435
436 template <class WMTKMesh>
438 {
439 double vol_tol;
440 std::vector<Tuple> adjacent_elements;
441 if constexpr (std::is_same_v<wmtk::TriMesh, WMTKMesh>)
442 {
443 adjacent_elements = {{e}};
444 const std::optional<Tuple> f = e.switch_face(*this);
445 if (f.has_value())
446 adjacent_elements.push_back(f.value());
447 }
448 else
449 {
450 adjacent_elements = get_incident_elements_for_edge(e);
451 }
452
453 Eigen::VectorXd adjacent_element_volumes(adjacent_elements.size());
454 for (int i = 0; i < adjacent_elements.size(); ++i)
455 adjacent_element_volumes[i] = element_volume(adjacent_elements[i]);
456
457 return adjacent_element_volumes;
458 }
459
460 template <class WMTKMesh>
461 void WildRemesher<WMTKMesh>::extend_local_patch(std::vector<Tuple> &patch) const
462 {
463 const size_t starting_size = patch.size();
464
465 std::unordered_set<size_t> element_ids;
466 for (const Tuple &t : patch)
467 element_ids.insert(element_id(t));
468
469 for (size_t i = 0; i < starting_size; ++i)
470 {
471 const size_t id = element_id(patch[i]);
472 for (int j = 0; j < EDGES_PER_ELEMENT; ++j)
473 {
474 const Tuple e = WMTKMesh::tuple_from_edge(id, j);
475
476 for (const Tuple &t : get_incident_elements_for_edge(e))
477 {
478 const size_t t_id = element_id(t);
479 if (element_ids.find(t_id) == element_ids.end())
480 {
481 patch.push_back(t);
482 element_ids.insert(t_id);
483 }
484 }
485 }
486 }
487 }
488
489 template <class WMTKMesh>
491 {
492 min.setConstant(dim(), std::numeric_limits<double>::infinity());
493 max.setConstant(dim(), -std::numeric_limits<double>::infinity());
494
495 for (const size_t vid : element_vids(t))
496 {
497 min = min.cwiseMin(vertex_attrs[vid].rest_position);
498 max = max.cwiseMax(vertex_attrs[vid].rest_position);
499 }
500 }
501
502 template <class WMTKMesh>
505 const EdgeMap<typename EdgeAttributes::EnergyRank> &contact_ranks) const
506 {
507 const std::vector<Tuple> edges = WMTKMesh::get_edges();
508
509 // Create two vertices per edge to get per edge values.
510 const int n_vertices = 2 * edges.size();
511
512 std::vector<std::vector<int>> elements(edges.size(), std::vector<int>(2));
513 Eigen::MatrixXd rest_positions(n_vertices, dim());
514 Eigen::MatrixXd displacements(n_vertices, dim());
515 Eigen::VectorXd energy_ranks(n_vertices);
516 Eigen::VectorXd elastic_energy_ranks(n_vertices);
517 Eigen::VectorXd contact_energy_ranks(n_vertices);
518
519 for (int ei = 0; ei < edges.size(); ei++)
520 {
521 const std::array<size_t, 2> vids = {{
522 edges[ei].vid(*this),
523 edges[ei].switch_vertex(*this).vid(*this),
524 }};
525
526 for (int vi = 0; vi < vids.size(); ++vi)
527 {
528 elements[ei][vi] = 2 * ei + vi;
529 rest_positions.row(elements[ei][vi]) = vertex_attrs[vids[vi]].rest_position;
530 displacements.row(elements[ei][vi]) = vertex_attrs[vids[vi]].displacement();
531 energy_ranks(elements[ei][vi]) = int(edge_attr(edges[ei].eid(*this)).energy_rank);
532 elastic_energy_ranks(elements[ei][vi]) = int(elastic_ranks.at(vids));
533 if (!contact_ranks.empty())
534 contact_energy_ranks(elements[ei][vi]) = int(contact_ranks.at(vids));
535 }
536 }
537
538 const double t0 = state.args["time"]["t0"];
539 const double dt = state.args["time"]["dt"];
540 const double save_dt = dt / 3;
541 // current_time = t0 + t * dt for t = 1, 2, 3, ...
542 // t = (current_time - t0) / dt
543 const int time_steps = int(round((current_time - t0) / save_dt));
544
545 // 0 -> 0 * dt + save_dt
546 // 1 -> 1 * dt + save_dt
547 // 2 -> 2 * dt + save_dt
548
549 const auto vtu_name = [&](int i) -> std::string {
550 return state.resolve_output_path(fmt::format("edge_ranks_{:d}.vtu", i));
551 };
552
553 paraviewo::VTUWriter writer;
554 writer.add_field("displacement", displacements);
555 writer.add_field("elastic_energy_rank", elastic_energy_ranks);
556 if (!contact_ranks.empty())
557 writer.add_field("contact_energy_rank", contact_energy_ranks);
558 writer.add_field("energy_rank", energy_ranks);
559 writer.write_mesh(vtu_name(time_steps - 3), rest_positions, elements, /*is_simplicial=*/true, /*has_poly=*/false);
560
561 writer.add_field("displacement", displacements);
562 writer.add_field("elastic_energy_rank", elastic_energy_ranks);
563 if (!contact_ranks.empty())
564 writer.add_field("contact_energy_rank", contact_energy_ranks);
565 writer.add_field("energy_rank", energy_ranks);
566 writer.write_mesh(vtu_name(time_steps - 2), rest_positions, elements, /*is_simplicial=*/true, /*has_poly=*/false);
567
568 writer.add_field("displacement", displacements);
569 writer.add_field("elastic_energy_rank", elastic_energy_ranks);
570 if (!contact_ranks.empty())
571 writer.add_field("contact_energy_rank", contact_energy_ranks);
572 writer.add_field("energy_rank", energy_ranks);
573 writer.write_mesh(vtu_name(time_steps - 1), rest_positions, elements, /*is_simplicial=*/true, /*has_poly=*/false);
574
575 state.out_geom.save_pvd(
576 state.resolve_output_path("edge_ranks.pvd"), vtu_name,
577 time_steps - 1, /*t0=*/save_dt, save_dt);
578 }
579
580 // -------------------------------------------------------------------------
581 // Template specializations
582
583 template class WildRemesher<wmtk::TriMesh>;
584 template class WildRemesher<wmtk::TetMesh>;
585
586} // namespace polyfem::mesh
#define POLYFEM_REMESHER_SCOPED_TIMER(name)
Definition Remesher.hpp:15
#define VERTEX_ATTRIBUTE_SETTER(name, attribute)
#define VERTEX_ATTRIBUTE_GETTER(name, attribute)
main class that contains the polyfem solver and all its state
Definition State.hpp:79
std::unordered_map< std::array< size_t, 2 >, T, polyfem::utils::HashUnorderedArray< size_t, 2 >, polyfem::utils::EqualUnorderedArray< size_t, 2 > > EdgeMap
Map from a (sorted) edge to an integer (ID)
Definition Remesher.hpp:29
std::variant< EdgeMap< T >, FaceMap< T > > BoundaryMap
Definition Remesher.hpp:42
std::unordered_map< std::array< size_t, 3 >, T, polyfem::utils::HashUnorderedArray< size_t, 3 >, polyfem::utils::EqualUnorderedArray< size_t, 3 > > FaceMap
Map from a (sorted) edge to an integer (ID)
Definition Remesher.hpp:35
virtual void init(const Eigen::MatrixXd &rest_positions, const Eigen::MatrixXd &positions, const Eigen::MatrixXi &elements, const Eigen::MatrixXd &projection_quantities, const BoundaryMap< int > &boundary_to_id, const std::vector< int > &body_ids, const EdgeMap< double > &elastic_energy, const EdgeMap< double > &contact_energy)
Initialize the mesh.
Definition Remesher.cpp:35
Eigen::MatrixXi elements() const override
Exports elements of the stored mesh.
Eigen::MatrixXd projection_quantities() const override
Exports projected quantities of the stored mesh.
double deformed_edge_length(const Tuple &e) const
void set_projection_quantities(const Eigen::MatrixXd &projection_quantities) override
Set projected quantities of the stored mesh.
void write_edge_ranks_mesh(const EdgeMap< typename EdgeAttributes::EnergyRank > &elastic_ranks, const EdgeMap< typename EdgeAttributes::EnergyRank > &contact_ranks) const
std::vector< Tuple > boundary_facets(std::vector< int > *boundary_ids=nullptr) const
Get the boundary facets of the mesh.
void set_body_ids(const std::vector< int > &body_ids) override
Set the body IDs of all elements.
virtual void init(const Eigen::MatrixXd &rest_positions, const Eigen::MatrixXd &positions, const Eigen::MatrixXi &elements, const Eigen::MatrixXd &projection_quantities, const BoundaryMap< int > &boundary_to_id, const std::vector< int > &body_ids, const EdgeMap< double > &elastic_energy, const EdgeMap< double > &contact_energy) override
Initialize the mesh.
typename WMTKMesh::Tuple Tuple
WildRemesher(const State &state, const Eigen::MatrixXd &obstacle_displacements, const Eigen::MatrixXd &obstacle_vals, const double current_time, const double starting_energy)
Construct a new WildRemesher object.
void set_boundary_ids(const BoundaryMap< int > &boundary_to_id) override
Set the boundary IDs of all edges.
std::vector< Tuple > get_edges_for_elements(const std::vector< Tuple > &elements) const
double rest_edge_length(const Tuple &e) const
Compute the length of an edge.
VectorNd rest_edge_center(const Tuple &e) const
Compute the center of the edge.
std::vector< int > body_ids() const override
Exports body ids of the stored mesh.
void element_aabb(const Tuple &t, polyfem::VectorNd &el_min, polyfem::VectorNd &el_max) const
Get a AABB for an element.
std::vector< int > boundary_nodes(const Eigen::VectorXi &vertex_to_basis) const override
Get the boundary nodes of the stored mesh.
Eigen::VectorXd edge_adjacent_element_volumes(const Tuple &e) const
void extend_local_patch(std::vector< Tuple > &patch) const
Extend the local patch by including neighboring elements.
Eigen::Matrix< double, DIM, 1 > VectorNd
void set_fixed(const std::vector< bool > &fixed) override
Set if a vertex is fixed.
bool invariants(const std::vector< Tuple > &new_tris) override
Check if invariants are satisfied.
VectorNd deformed_edge_center(const Tuple &e) const
WildTetRemesher::Tuple Tuple
Remesher::EdgeMap< typename WildRemesher< WMTKMesh >::EdgeAttributes::EnergyRank > rank_edges(const Remesher::EdgeMap< double > &edge_energy, const json &args)
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > VectorNd
Definition Types.hpp:9
nlohmann::json json
Definition Common.hpp:9
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71