Wildmeshing Toolkit
TetMesh.cpp
Go to the documentation of this file.
1 #include "TetMesh.hpp"
2 
3 
5 #include <numeric>
12 #include <wmtk/utils/Logger.hpp>
13 
14 namespace wmtk {
15 
16 using namespace autogen;
17 TetMesh::~TetMesh() = default;
18 
20  : MeshCRTP<TetMesh>(3)
21  , m_vt_handle(register_attribute_typed<int64_t>("m_vt", PrimitiveType::Vertex, 1, false, -1))
22  , m_et_handle(register_attribute_typed<int64_t>("m_et", PrimitiveType::Edge, 1, false, -1))
23  , m_ft_handle(register_attribute_typed<int64_t>("m_ft", PrimitiveType::Triangle, 1, false, -1))
24  , m_tv_handle(
25  register_attribute_typed<int64_t>("m_tv", PrimitiveType::Tetrahedron, 4, false, -1))
26  , m_te_handle(
27  register_attribute_typed<int64_t>("m_te", PrimitiveType::Tetrahedron, 6, false, -1))
28  , m_tf_handle(
29  register_attribute_typed<int64_t>("m_tf", PrimitiveType::Tetrahedron, 4, false, -1))
30  , m_tt_handle(
31  register_attribute_typed<int64_t>("m_tt", PrimitiveType::Tetrahedron, 4, false, -1))
32 {
34 }
35 
36 
38  : MeshCRTP<TetMesh>(std::move(o))
39 {
40  m_vt_handle = o.m_vt_handle;
41  m_et_handle = o.m_et_handle;
42  m_ft_handle = o.m_ft_handle;
43  m_tv_handle = o.m_tv_handle;
44  m_te_handle = o.m_te_handle;
45  m_tf_handle = o.m_tf_handle;
46  m_tt_handle = o.m_tt_handle;
47 
49 }
51 {
52  Mesh::operator=(std::move(o));
53  m_vt_handle = o.m_vt_handle;
54  m_et_handle = o.m_et_handle;
55  m_ft_handle = o.m_ft_handle;
56  m_tv_handle = o.m_tv_handle;
57  m_te_handle = o.m_te_handle;
58  m_tf_handle = o.m_tf_handle;
59  m_tt_handle = o.m_tt_handle;
60 
62  return *this;
63 }
64 
65 
67 {
68  m_vt_accessor = std::make_unique<attribute::Accessor<int64_t, TetMesh>>(*this, m_vt_handle);
69  m_et_accessor = std::make_unique<attribute::Accessor<int64_t, TetMesh>>(*this, m_et_handle);
70  m_ft_accessor = std::make_unique<attribute::Accessor<int64_t, TetMesh>>(*this, m_ft_handle);
71 
72  m_tv_accessor = std::make_unique<attribute::Accessor<int64_t, TetMesh>>(*this, m_tv_handle);
73  m_te_accessor = std::make_unique<attribute::Accessor<int64_t, TetMesh>>(*this, m_te_handle);
74  m_tf_accessor = std::make_unique<attribute::Accessor<int64_t, TetMesh>>(*this, m_tf_handle);
75  m_tt_accessor = std::make_unique<attribute::Accessor<int64_t, TetMesh>>(*this, m_tt_handle);
76 }
77 
78 
80  Eigen::Ref<const RowVectors4l> TV,
81  Eigen::Ref<const RowVectors6l> TE,
82  Eigen::Ref<const RowVectors4l> TF,
83  Eigen::Ref<const RowVectors4l> TT,
84  Eigen::Ref<const VectorXl> VT,
85  Eigen::Ref<const VectorXl> ET,
86  Eigen::Ref<const VectorXl> FT)
87 
88 {
89  // reserve memory for attributes
90 
91  std::vector<int64_t> cap{
92  static_cast<int64_t>(VT.rows()),
93  static_cast<int64_t>(ET.rows()),
94  static_cast<int64_t>(FT.rows()),
95  static_cast<int64_t>(TT.rows())};
96  set_capacities(cap);
97 
98  // get Accessors for topology
99  auto vt_accessor = create_accessor<int64_t>(m_vt_handle);
100  auto et_accessor = create_accessor<int64_t>(m_et_handle);
101  auto ft_accessor = create_accessor<int64_t>(m_ft_handle);
102  auto tv_accessor = create_accessor<int64_t>(m_tv_handle);
103  auto te_accessor = create_accessor<int64_t>(m_te_handle);
104  auto tf_accessor = create_accessor<int64_t>(m_tf_handle);
105  auto tt_accessor = create_accessor<int64_t>(m_tt_handle);
109  attribute::FlagAccessor<TetMesh> t_flag_accessor =
111 
112  // iterate over the matrices and fill attributes
113  for (int64_t i = 0; i < capacity(PrimitiveType::Tetrahedron); ++i) {
114  tv_accessor.index_access().vector_attribute<4>(i) = TV.row(i).transpose();
115  te_accessor.index_access().vector_attribute<6>(i) = TE.row(i).transpose();
116  tf_accessor.index_access().vector_attribute<4>(i) = TF.row(i).transpose();
117  tt_accessor.index_access().vector_attribute<4>(i) = TT.row(i).transpose();
118  t_flag_accessor.index_access().activate(i);
119  e_flag_accessor.index_access().activate(i);
120  }
121  // m_vt
122  for (int64_t i = 0; i < capacity(PrimitiveType::Vertex); ++i) {
123  vt_accessor.index_access().scalar_attribute(i) = VT(i);
124  v_flag_accessor.index_access().activate(i);
125  }
126  // m_et
127  for (int64_t i = 0; i < capacity(PrimitiveType::Edge); ++i) {
128  et_accessor.index_access().scalar_attribute(i) = ET(i);
129  e_flag_accessor.index_access().activate(i);
130  }
131  // m_ft
132  for (int64_t i = 0; i < capacity(PrimitiveType::Triangle); ++i) {
133  ft_accessor.index_access().scalar_attribute(i) = FT(i);
134  f_flag_accessor.index_access().activate(i);
135  }
136 }
137 
138 
139 void TetMesh::initialize(Eigen::Ref<const RowVectors4l> T, bool is_free)
140 {
141  this->m_is_free = is_free;
142  auto [TE, TF, TT, VT, ET, FT] = tetmesh_topology_initialization(T);
143  if (is_free) {
144  TT.setConstant(-1);
145  }
146  initialize(T, TE, TF, TT, VT, ET, FT);
147 }
148 void TetMesh::initialize_free(int64_t count)
149 {
150  RowVectors4l S(count, 4);
151  std::iota(S.data(), S.data() + S.size(), int64_t(0));
152  initialize(S, true);
153 }
154 
156 {
157  int64_t t = m_vt_accessor->index_access().const_scalar_attribute(id);
158  auto tv = m_tv_accessor->index_access().const_vector_attribute<4>(t);
159  int64_t lvid = -1;
160 
161  for (int64_t i = 0; i < 4; ++i) {
162  if (tv(i) == id) {
163  lvid = i;
164  break;
165  }
166  }
167 
169  assert(is_valid(v_tuple));
170  return v_tuple;
171 }
172 
174 {
175  int64_t t = m_et_accessor->index_access().const_scalar_attribute(id);
176  auto te = m_te_accessor->index_access().const_vector_attribute<6>(t);
177 
178  int64_t leid = -1;
179 
180  for (int64_t i = 0; i < 6; ++i) {
181  if (te(i) == id) {
182  leid = i;
183  break;
184  }
185  }
187  assert(is_valid(e_tuple));
188  return e_tuple;
189 }
190 
192 {
193  int64_t t = m_ft_accessor->index_access().const_scalar_attribute(id);
194  auto tf = m_tf_accessor->index_access().const_vector_attribute<4>(t);
195 
196  int64_t lfid = -1;
197 
198  for (int64_t i = 0; i < 4; ++i) {
199  if (tf(i) == id) {
200  lfid = i;
201  break;
202  }
203  }
204 
205 
206  assert(lfid >= 0);
208  assert(is_valid(f_tuple));
209  return f_tuple;
210 }
211 
213 {
214  const int64_t lvid = 0;
215  const auto [nlvid, leid, lfid] = autogen::tet_mesh::auto_3d_table_complete_vertex[lvid];
216  assert(lvid == nlvid);
217 
218  Tuple t_tuple = Tuple(lvid, leid, lfid, id);
219  assert(is_ccw(t_tuple));
220  assert(is_valid(t_tuple));
221  return t_tuple;
222 }
223 
224 Tuple TetMesh::tuple_from_id(const PrimitiveType type, const int64_t gid) const
225 {
226  switch (type) {
227  case PrimitiveType::Vertex: {
228  return vertex_tuple_from_id(gid);
229  break;
230  }
231  case PrimitiveType::Edge: {
232  return edge_tuple_from_id(gid);
233  break;
234  }
236  return face_tuple_from_id(gid);
237  break;
238  }
240  return tet_tuple_from_id(gid);
241  break;
242  }
243  default: assert(false); // "Invalid primitive type"
244  }
245 
246  return Tuple();
247 }
248 
249 
251 {
252  assert(is_valid(tuple));
253  switch (type) {
254  // bool ccw = is_ccw(tuple);
256  assert(!is_boundary_face(tuple));
257  // need test
258  const int64_t gvid = id(tuple, PrimitiveType::Vertex);
259  const int64_t geid = id(tuple, PrimitiveType::Edge);
260  const int64_t gfid = id(tuple, PrimitiveType::Triangle);
261 
262  auto tt = m_tt_accessor->const_vector_attribute<4>(tuple);
263 
264  int64_t gcid_new = tt(tuple.m_local_fid);
265 
266  /*handle exception here*/
267  assert(gcid_new != -1);
268  // check if is_boundary allows removing this exception in 3d cases
269  // if (gcid_new == -1) {
270  // return Tuple(-1, -1, -1, -1, -1);
271  // }
272  /*handle exception end*/
273 
274  int64_t lvid_new = -1, leid_new = -1, lfid_new = -1;
275 
276  auto tv = m_tv_accessor->index_access().const_vector_attribute<4>(gcid_new);
277 
278  auto te = m_te_accessor->index_access().const_vector_attribute<6>(gcid_new);
279 
280  auto tf = m_tf_accessor->index_access().const_vector_attribute<4>(gcid_new);
281 
282  for (int64_t i = 0; i < 4; ++i) {
283  if (tv(i) == gvid) {
284  lvid_new = i;
285  }
286  if (tf(i) == gfid) {
287  lfid_new = i;
288  }
289  }
290 
291  for (int64_t i = 0; i < 6; ++i) {
292  if (te(i) == geid) {
293  leid_new = i;
294  break; // check if the break is correct
295  }
296  }
297 
298 
299  assert(lvid_new != -1);
300  assert(leid_new != -1);
301  assert(lfid_new != -1);
302 
303  const Tuple res(lvid_new, leid_new, lfid_new, gcid_new);
304  assert(is_valid(res));
305  return res;
306  }
308  case PrimitiveType::Edge:
310  default: return autogen::tet_mesh::local_switch_tuple(tuple, type);
311  }
312 }
313 
314 bool TetMesh::is_ccw(const Tuple& tuple) const
315 {
316  assert(is_valid(tuple));
317  return autogen::tet_mesh::is_ccw(tuple);
318 }
319 
320 bool TetMesh::is_valid(const Tuple& tuple) const
321 {
322  if (!Mesh::is_valid(tuple)) {
323  return false;
324  }
325  const bool is_connectivity_valid = tuple.m_local_vid >= 0 && tuple.m_local_eid >= 0 &&
326  tuple.m_local_fid >= 0 && tuple.m_global_cid >= 0 &&
328 
329  if (!is_connectivity_valid) {
330 #if !defined(NDEBUG)
331  logger().trace(
332  "tuple.m_local_vid={} >= 0 && tuple.m_local_eid={} >= 0 &&"
333  "tuple.m_local_fid={} >= 0 &&"
334  " tuple.m_global_cid={} >= 0 &&"
335  " autogen::tet_mesh::tuple_is_valid_for_ccw(tuple)={}",
336  tuple.m_local_vid,
337  tuple.m_local_eid,
338  tuple.m_local_fid,
339  tuple.m_global_cid,
341  assert(tuple.m_local_vid >= 0);
342  assert(tuple.m_local_eid >= 0);
343  assert(tuple.m_local_fid >= 0);
344  assert(tuple.m_global_cid >= 0);
346 #endif
347  return false;
348  }
349 
350  return true;
351 }
352 
353 bool TetMesh::is_boundary(PrimitiveType pt, const Tuple& tuple) const
354 {
355  switch (pt) {
356  case PrimitiveType::Vertex: return is_boundary_vertex(tuple);
357  case PrimitiveType::Edge: return is_boundary_edge(tuple);
358  case PrimitiveType::Triangle: return is_boundary_face(tuple);
360  default: break;
361  }
362  assert(
363  false); // "tried to compute the boundary of an tet mesh for an invalid simplex dimension"
364  return false;
365 }
366 
367 
368 bool TetMesh::is_boundary_face(const Tuple& tuple) const
369 {
370  const attribute::Accessor<int64_t> tt_accessor = create_const_accessor<int64_t>(m_tt_handle);
371  return tt_accessor.const_vector_attribute<4>(tuple)(tuple.m_local_fid) < 0;
372 }
373 
374 bool TetMesh::is_boundary_edge(const Tuple& edge) const
375 {
377  *this,
378  simplex::Simplex::edge(*this, edge),
380  if (is_boundary_face(f)) {
381  return true;
382  }
383  }
384  return false;
385 }
386 bool TetMesh::is_boundary_vertex(const Tuple& vertex) const
387 {
388  // go through all faces and check if they are boundary
389  const simplex::SimplexCollection neigh =
390  wmtk::simplex::open_star(*this, simplex::Simplex::vertex(*this, vertex));
392  if (is_boundary(s)) {
393  return true;
394  }
395  }
396 
397  return false;
398 }
399 
401 {
402  // get Accessors for topology
403  const attribute::Accessor<int64_t> tv_accessor = create_const_accessor<int64_t>(m_tv_handle);
404  const attribute::Accessor<int64_t> te_accessor = create_const_accessor<int64_t>(m_te_handle);
405  const attribute::Accessor<int64_t> tf_accessor = create_const_accessor<int64_t>(m_tf_handle);
406  const attribute::Accessor<int64_t> tt_accessor = create_const_accessor<int64_t>(m_tt_handle);
407  const attribute::Accessor<int64_t> vt_accessor = create_const_accessor<int64_t>(m_vt_handle);
408  const attribute::Accessor<int64_t> et_accessor = create_const_accessor<int64_t>(m_et_handle);
409  const attribute::Accessor<int64_t> ft_accessor = create_const_accessor<int64_t>(m_ft_handle);
410  const attribute::FlagAccessor<TetMesh> v_flag_accessor =
413  const attribute::FlagAccessor<TetMesh> f_flag_accessor =
415  const attribute::FlagAccessor<TetMesh> t_flag_accessor =
417 
418 
419  for (int64_t i = 0; i < capacity(PrimitiveType::Tetrahedron); ++i) {
420  if (!t_flag_accessor.index_access().is_active(i)) {
421  continue;
422  }
423  auto tf = tf_accessor.index_access().const_vector_attribute<4>(i);
424  auto te = te_accessor.index_access().const_vector_attribute<6>(i);
425  auto tv = tv_accessor.index_access().const_vector_attribute<4>(i);
426 
427  bool bad_face = false;
428  for (int64_t j = 0; j < 6; ++j) {
429  int64_t ei = te(j);
430  if (!e_flag_accessor.index_access().is_active(ei)) {
431  wmtk::logger().error(
432  "Tet {} refers to edge {} at local index {} which was deleted",
433  i,
434  ei,
435  j);
436  bad_face = true;
437  }
438  }
439 
440  for (int64_t j = 0; j < 4; ++j) {
441  int64_t vi = tv(j);
442  int64_t fi = tf(j);
443  if (!v_flag_accessor.index_access().is_active(vi)) {
444  wmtk::logger().error(
445  "Tet {} refers to vertex{} at local index {} which was deleted",
446  i,
447  vi,
448  j);
449  bad_face = true;
450  }
451  if (!f_flag_accessor.index_access().is_active(fi)) {
452  wmtk::logger()
453  .error("Tet {} refers to face{} at local index {} which was deleted", i, fi, j);
454  bad_face = true;
455  }
456  }
457  if (bad_face) {
458  return false;
459  }
460  }
461  // VT and TV
462  for (int64_t i = 0; i < capacity(PrimitiveType::Vertex); ++i) {
463  if (!v_flag_accessor.index_access().is_active(i)) {
464  continue;
465  }
466  int cnt = 0;
467  for (int j = 0; j < 4; ++j) {
468  if (tv_accessor.index_access().const_vector_attribute<4>(
469  vt_accessor.index_access().const_scalar_attribute(i))[j] == i) {
470  cnt++;
471  }
472  }
473  if (cnt != 1) {
474  wmtk::logger().info("fail VT and TV");
475  return false;
476  }
477  }
478 
479  // ET and TE
480  for (int64_t i = 0; i < capacity(PrimitiveType::Edge); ++i) {
481  if (!e_flag_accessor.index_access().is_active(i)) {
482  continue;
483  }
484  int cnt = 0;
485  for (int j = 0; j < 6; ++j) {
486  if (te_accessor.index_access().const_vector_attribute<6>(
487  et_accessor.index_access().const_scalar_attribute(i))[j] == i) {
488  cnt++;
489  }
490  }
491  if (cnt != 1) {
492  wmtk::logger().info("fail ET and TE");
493  return false;
494  }
495  }
496 
497  // FT and TF
498  for (int64_t i = 0; i < capacity(PrimitiveType::Triangle); ++i) {
499  if (!f_flag_accessor.index_access().is_active(i)) {
500  continue;
501  }
502  int cnt = 0;
503  for (int j = 0; j < 4; ++j) {
504  if (tf_accessor.index_access().const_vector_attribute<4>(
505  ft_accessor.index_access().const_scalar_attribute(i))[j] == i) {
506  cnt++;
507  }
508  }
509  if (cnt != 1) {
510  wmtk::logger().info("fail FT and TF");
511  return false;
512  }
513  }
514 
515  // TF and TT
516  for (int64_t i = 0; i < capacity(PrimitiveType::Tetrahedron); ++i) {
517  if (!t_flag_accessor.index_access().is_active(i)) {
518  continue;
519  }
520 
521  for (int j = 0; j < 4; ++j) {
522  int64_t nb = tt_accessor.index_access().const_vector_attribute<4>(i)(j);
523  if (nb == -1) {
524  if (ft_accessor.index_access().const_scalar_attribute(
525  tf_accessor.index_access().const_vector_attribute<4>(i)(j)) != i) {
526  wmtk::logger().error("FT[TF[{},{}]] != {}", i, j, i);
527  return false;
528  }
529  continue;
530  }
531 
532  int cnt = 0;
533  int id_in_nb;
534  for (int k = 0; k < 4; ++k) {
535  if (tt_accessor.index_access().const_vector_attribute<4>(nb)(k) == i) {
536  cnt++;
537  id_in_nb = k;
538  }
539  }
540  if (cnt != 1) {
541  wmtk::logger().error("Tet {} was adjacent to tet {} {} <= 1 times", nb, i, cnt);
542  return false;
543  }
544 
545  if (tf_accessor.index_access().const_vector_attribute<4>(i)(j) !=
546  tf_accessor.index_access().const_vector_attribute<4>(nb)(id_in_nb)) {
547  wmtk::logger().error(
548  "TF[{},{}] = {} != {} = TF[{},{}] even though TT[{},{}] == {}",
549  i,
550  j,
551  tf_accessor.index_access().const_vector_attribute<4>(i)(j),
552  tf_accessor.index_access().const_vector_attribute<4>(nb)(id_in_nb),
553  nb,
554  id_in_nb,
555  i,
556  j,
557  nb);
558 
559  return false;
560  }
561  }
562  }
563 
564  return true;
565 }
566 
567 std::vector<std::vector<TypedAttributeHandle<int64_t>>> TetMesh::connectivity_attributes() const
568 {
569  std::vector<std::vector<TypedAttributeHandle<int64_t>>> handles(4);
570 
571  handles[0].push_back(m_tv_handle);
572  handles[1].push_back(m_te_handle);
573  handles[2].push_back(m_tf_handle);
574 
575  handles[3].push_back(m_tt_handle);
576  handles[3].push_back(m_vt_handle);
577  handles[3].push_back(m_et_handle);
578  handles[3].push_back(m_ft_handle);
579 
580  return handles;
581 }
582 
583 Tuple TetMesh::tuple_from_global_ids(int64_t tid, int64_t fid, int64_t eid, int64_t vid) const
584 {
585  auto tv = m_tv_accessor->index_access().const_vector_attribute<4>(tid);
586  auto te = m_te_accessor->index_access().const_vector_attribute<6>(tid);
587  auto tf = m_tf_accessor->index_access().const_vector_attribute<4>(tid);
588 
589  int64_t lvid = -1, leid = -1, lfid = -1;
590 
591  for (int j = 0; j < 4; ++j) {
592  if (tv(j) == vid) {
593  lvid = j;
594  }
595  if (tf(j) == fid) {
596  lfid = j;
597  }
598  }
599 
600  for (int j = 0; j < 6; ++j) {
601  if (te(j) == eid) {
602  leid = j;
603  break;
604  }
605  }
606 
607  assert(lvid != -1);
608  assert(leid != -1);
609  assert(lfid != -1);
610 
611  return Tuple(lvid, leid, lfid, tid);
612 }
613 
614 std::vector<Tuple> TetMesh::orient_vertices(const Tuple& tuple) const
615 {
616  int64_t cid = tuple.m_global_cid;
617  return {Tuple(0, 0, 2, cid), Tuple(1, 0, 3, cid), Tuple(2, 1, 1, cid), Tuple(3, 2, 2, cid)};
618 }
619 
620 
621 } // namespace wmtk
A Curiously Recurring Template Pattern shim to enable generic specialization of functions.
Definition: MeshCRTP.hpp:24
int64_t id(const Tuple &tuple, PrimitiveType type) const
return the global id of the Tuple of the given dimension
Definition: Mesh.hpp:1021
bool m_is_free
Definition: Mesh.hpp:847
void set_capacities(std::vector< int64_t > capacities)
int64_t capacity(PrimitiveType type) const
read in the m_capacities return the upper bound for the number of entities of the given dimension
bool is_free() const
Definition: Mesh.hpp:988
virtual bool is_valid(const Tuple &tuple) const
check validity of tuple including its hash
Definition: Mesh.cpp:113
Mesh & operator=(const Mesh &other)=delete
const attribute::FlagAccessor< Mesh > get_flag_accessor(PrimitiveType type) const
Definition: Mesh.cpp:159
TypedAttributeHandle< int64_t > m_tt_handle
Definition: TetMesh.hpp:95
std::unique_ptr< attribute::Accessor< int64_t, TetMesh > > m_tv_accessor
Definition: TetMesh.hpp:101
TetMesh & operator=(const TetMesh &o)=delete
TypedAttributeHandle< int64_t > m_vt_handle
Definition: TetMesh.hpp:87
Tuple tuple_from_id(const PrimitiveType type, const int64_t gid) const final override
internal function that returns the tuple of requested type, and has the global index cid
Definition: TetMesh.cpp:224
bool is_boundary_vertex(const Tuple &tuple) const
Definition: TetMesh.cpp:386
std::unique_ptr< attribute::Accessor< int64_t, TetMesh > > m_et_accessor
Definition: TetMesh.hpp:98
void make_cached_accessors()
Definition: TetMesh.cpp:66
std::vector< std::vector< TypedAttributeHandle< int64_t > > > connectivity_attributes() const final override
Returns a vector of vectors of attribute handles.
Definition: TetMesh.cpp:567
std::unique_ptr< attribute::Accessor< int64_t, TetMesh > > m_tt_accessor
Definition: TetMesh.hpp:104
Tuple switch_tuple(const Tuple &tuple, PrimitiveType type) const final override
switch the orientation of the Tuple of the given dimension
Definition: TetMesh.cpp:250
Tuple tuple_from_global_ids(int64_t tid, int64_t fid, int64_t eid, int64_t vid) const
Definition: TetMesh.cpp:583
bool is_boundary(const simplex::Simplex &tuple) const
check if a simplex lies on a boundary or not
Definition: Mesh.cpp:107
TypedAttributeHandle< int64_t > m_tv_handle
Definition: TetMesh.hpp:92
bool is_boundary_edge(const Tuple &tuple) const
Definition: TetMesh.cpp:374
Tuple edge_tuple_from_id(int64_t id) const
Definition: TetMesh.cpp:173
Tuple face_tuple_from_id(int64_t id) const
Definition: TetMesh.cpp:191
void initialize_free(int64_t count)
Definition: TetMesh.cpp:148
Tuple tet_tuple_from_id(int64_t id) const
Definition: TetMesh.cpp:212
std::unique_ptr< attribute::Accessor< int64_t, TetMesh > > m_tf_accessor
Definition: TetMesh.hpp:103
TypedAttributeHandle< int64_t > m_tf_handle
Definition: TetMesh.hpp:94
bool is_ccw(const Tuple &tuple) const final override
TODO this needs dimension?
Definition: TetMesh.cpp:314
Tuple vertex_tuple_from_id(int64_t id) const
Definition: TetMesh.cpp:155
void initialize(Eigen::Ref< const RowVectors4l > TV, Eigen::Ref< const RowVectors6l > TE, Eigen::Ref< const RowVectors4l > TF, Eigen::Ref< const RowVectors4l > TT, Eigen::Ref< const VectorXl > VT, Eigen::Ref< const VectorXl > ET, Eigen::Ref< const VectorXl > FT)
Definition: TetMesh.cpp:79
TypedAttributeHandle< int64_t > m_te_handle
Definition: TetMesh.hpp:93
bool is_connectivity_valid() const final override
Definition: TetMesh.cpp:400
TypedAttributeHandle< int64_t > m_et_handle
Definition: TetMesh.hpp:89
std::unique_ptr< attribute::Accessor< int64_t, TetMesh > > m_ft_accessor
Definition: TetMesh.hpp:99
TypedAttributeHandle< int64_t > m_ft_handle
Definition: TetMesh.hpp:90
~TetMesh() override
std::unique_ptr< attribute::Accessor< int64_t, TetMesh > > m_vt_accessor
Definition: TetMesh.hpp:97
std::unique_ptr< attribute::Accessor< int64_t, TetMesh > > m_te_accessor
Definition: TetMesh.hpp:102
bool is_valid(const Tuple &tuple) const final override
check validity of tuple including its hash
Definition: TetMesh.cpp:320
bool is_boundary_face(const Tuple &tuple) const
Definition: TetMesh.cpp:368
std::vector< Tuple > orient_vertices(const Tuple &t) const override
Definition: TetMesh.cpp:614
int8_t m_local_vid
Definition: Tuple.hpp:47
int8_t m_local_eid
Definition: Tuple.hpp:48
int64_t m_global_cid
Definition: Tuple.hpp:46
int8_t m_local_fid
Definition: Tuple.hpp:49
CachingBaseType & index_access()
Definition: Accessor.hpp:95
ConstMapResult< D > const_vector_attribute(const ArgType &t) const
T const_scalar_attribute(const int64_t index) const
ConstMapResult< D > const_vector_attribute(const int64_t index) const
const IndexBaseType & index_access() const
bool is_active(int64_t t) const
const std::vector< Simplex > & simplex_vector() const
Return const reference to the simplex vector.
static Simplex edge(const Mesh &m, const Tuple &t)
Definition: Simplex.hpp:61
static Simplex vertex(const Mesh &m, const Tuple &t)
Definition: Simplex.hpp:56
Definition: autodiff.h:995
bool tuple_is_valid_for_ccw(const Tuple &t)
Definition: is_ccw.hxx:17
Tuple local_switch_tuple(const Tuple &t, PrimitiveType pt)
bool is_ccw(const Tuple &t)
Definition: is_ccw.hxx:10
Tuple get_tuple_from_simplex_local_edge_id(int8_t local_id, int64_t global_id)
Tuple get_tuple_from_simplex_local_face_id(int8_t local_id, int64_t global_id)
const int64_t auto_3d_table_complete_vertex[4][3]
Tuple get_tuple_from_simplex_local_vertex_id(int8_t local_id, int64_t global_id)
SimplexCollection open_star(const Mesh &mesh, const Simplex &simplex, const bool sort_and_clean)
Definition: open_star.cpp:12
CofacesSingleDimensionIterable cofaces_single_dimension_iterable(const Mesh &mesh, const Simplex &simplex, const PrimitiveType cofaces_type)
Definition: Accessor.hpp:6
RowVectors< int64_t, 4 > RowVectors4l
Definition: Types.hpp:48
spdlog::logger & logger()
Retrieves the current logger.
Definition: Logger.cpp:58
std::tuple< RowVectors6l, RowVectors4l, RowVectors4l, VectorXl, VectorXl, VectorXl > tetmesh_topology_initialization(Eigen::Ref< const RowVectors4l > T)
Given the mesh connectivity in matrix format, finds unique edges and faces and their relations.