libMesh
cell_polyhedron.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2026 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 #include "libmesh/cell_polyhedron.h"
19 
20 // Local includes
21 #include "libmesh/face_polygon.h"
22 #include "libmesh/enum_elem_quality.h"
23 #include "libmesh/hashword.h"
24 
25 // C++ includes
26 #include <array>
27 #include <unordered_map>
28 #include <unordered_set>
29 
30 
31 namespace libMesh
32 {
33 
34 // ------------------------------------------------------------
35 // Polyhedron class static member initializations
36 const int Polyhedron::num_children;
37 
38 // ------------------------------------------------------------
39 // Polyhedron class member functions
40 
41 
42 Polyhedron::Polyhedron (const std::vector<std::shared_ptr<Polygon>> & sides,
43  Elem * p) :
44  Cell(/* unused here */ 0, sides.size(), p, nullptr, nullptr),
45  _elemlinks_data(sides.size()+2), // neighbors + parent + interior_parent
46  _nodelinks_data(0), // We'll have to resize *this* later too!
47  _sidelinks_data(sides.size())
48  {
49  // Set our sides, and while we're at it figure out our node maps
50  // and side normal directions and edge lookup table, and count our
51  // sides' nodes. If we have internal nodes too then the subclass
52  // will append those afterward.
53  unsigned int nn = 0;
54  std::unordered_map<Node *, unsigned int> local_node_number;
55  std::unordered_set<std::pair<const Node *, const Node *>,
56  libMesh::hash> edges_seen;
57  std::unique_ptr<const Elem> edge;
58  for (unsigned int s : index_range(sides))
59  {
60  libmesh_assert(sides[s].get());
61  auto & side_tuple = _sidelinks_data[s];
62  std::get<0>(side_tuple) = sides[s];
63 
64  Polygon & side = *sides[s]; // not const, for writeable nodes
65  for (auto n : make_range(side.n_nodes()))
66  {
67  Node * node = side.node_ptr(n);
68  if (auto it = local_node_number.find(node);
69  it != local_node_number.end())
70  {
71  std::get<2>(side_tuple).push_back(it->second);
72  }
73  else
74  {
75  std::get<2>(side_tuple).push_back(nn);
76  local_node_number[node] = nn++;
77  _nodelinks_data.push_back(node);
78  }
79  }
80 
81  for (unsigned int e : make_range(side.n_edges()))
82  {
83  side.build_edge_ptr(edge, e);
84  auto edge_vertices = std::make_pair(edge->node_ptr(0), edge->node_ptr(1));
85  if (edge_vertices.first > edge_vertices.second)
86  std::swap(edge_vertices.first, edge_vertices.second);
87 
88  if (!edges_seen.count(edge_vertices))
89  {
90  edges_seen.insert(edge_vertices);
91  _edge_lookup.emplace_back(s, e);
92  }
93  }
94  }
95 
96  // Plan room for an extra node so we don't invalidate this->_nodes when we add to
97  // nodelinks_data in the derived class
98  _nodelinks_data.reserve(_nodelinks_data.size() + 1);
99 
100  // Do the manual initialization that Elem::Elem and Cell::Cell
101  // couldn't, now that we've resized both our vectors. No need to
102  // manually set nullptr, though, since std::vector does that.
103  this->_elemlinks = _elemlinks_data.data();
104  this->_nodes = _nodelinks_data.data();
105  this->_elemlinks[0] = p;
106 
107  libmesh_assert_equal_to(nn, this->n_nodes());
108 
109  // Figure out the orientation of our sides, now that we've got our
110  // nodes organized enough to find our center. The algorithm below
111  // only works for convex polyhedra, but that's all we're
112  // supporting for now.
113  Point center;
114  for (auto n : make_range(nn))
115  center.add (this->point(n));
116  center /= static_cast<Real>(nn);
117 
118  for (unsigned int s : index_range(sides))
119  {
120  const Polygon & side = *sides[s];
121  const Point x_i = side.point(0);
122  const Point n_i =
123  (side.point(1) - side.point(0)).cross
124  (side.point(0) - side.point(side.n_sides()-1)).unit();
125 
126  bool & inward_normal = std::get<1>(_sidelinks_data[s]);
127  inward_normal = (n_i * (center - x_i) > TOLERANCE);
128  }
129 
130  // We're betting a lot on "our polyhedra are all convex", so let's
131  // check that if we have time.
132 #ifdef DEBUG
133  for (unsigned int s : index_range(sides))
134  {
135  const Polygon & side = *sides[s];
136  const Point x_i = side.point(0);
137  const bool inward_normal = std::get<1>(this->_sidelinks_data[s]);
138 
139  const Point n_i =
140  (side.point(1) - side.point(0)).cross
141  (side.point(0) - side.point(side.n_sides()-1)).unit() *
142  (inward_normal ? -1 : 1);
143 
144  for (const Point & node : this->node_ref_range())
145  {
146  const Point d_n = node - x_i;
147  if (d_n * n_i > TOLERANCE * d_n.norm())
148  libmesh_not_implemented_msg
149  ("Cannot create a non-convex polyhedron");
150  }
151  }
152 #endif
153 
154  // Is this likely to ever be used? We may do refinement with
155  // polyhedra but it's probably not going to have a hierarchy...
156  if (p)
157  {
158  this->subdomain_id() = p->subdomain_id();
159  this->processor_id() = p->processor_id();
160  _map_type = p->mapping_type();
161  _map_data = p->mapping_data();
162 
163 #ifdef LIBMESH_ENABLE_AMR
164  this->set_p_level(p->p_level());
165 #endif
166  }
167 
168  // Make sure the interior parent isn't undefined
169  this->set_interior_parent(nullptr);
170  }
171 
172 
173 
174 Point Polyhedron::master_point (const unsigned int i) const
175 {
176  return this->point(i);
177 }
178 
179 
180 
182 {
183  for (unsigned int s : make_range(this->n_sides()))
184  {
185  const Polygon & side = *std::get<0>(this->_sidelinks_data[s]);
186  const Point x_i = side.point(0);
187  const bool inward_normal = std::get<1>(this->_sidelinks_data[s]);
188 
189  const Point n_i =
190  (side.point(1) - side.point(0)).cross
191  (side.point(0) - side.point(side.n_sides()-1)).unit() *
192  (inward_normal ? -1 : 1);
193 
194  for (const Point & node : this->node_ref_range())
195  {
196  const Point d_n = node - x_i;
197  if (d_n * n_i > TOLERANCE * d_n.norm())
198  return false;
199  }
200  }
201  return true;
202 }
203 
204 
205 
207  const Real eps) const
208 {
209  const unsigned int ns = this->n_sides();
210 
211  // Check that the point is on the same side of all the faces by
212  // testing whether:
213  //
214  // n_i.(p - x_i) <= 0
215  //
216  // for each i, where:
217  // n_i is the outward normal of face i,
218  // x_i is a point on face i.
219 
220  for (auto i : make_range(ns))
221  {
222  const Polygon & face = *std::get<0>(this->_sidelinks_data[i]);
223  const bool inward_normal = std::get<1>(this->_sidelinks_data[i]);
224 
225  const Point x_i = face.point(0);
226 
227  const Point n_i =
228  (face.point(1) - face.point(0)).cross
229  (face.point(0) - face.point(face.n_sides()-1)).unit() *
230  (inward_normal ? -1 : 1);
231 
232  // This only works for polyhedra with flat sides.
233 #ifdef DEBUG
234  for (auto j : make_range(face.n_sides()-1))
235  {
236  const Point x_j = face.point(j+1);
237  const Point d_j = x_j - x_i;
238  if (std::abs(d_j * n_i) > eps * d_j.norm())
239  libmesh_not_implemented_msg
240  ("Polyhedra with non-flat sides are not fully supported.");
241  }
242 #endif
243 
244  if (n_i * (p - x_i) > eps)
245  return false;
246  }
247 
248  return true;
249 }
250 
251 
252 
253 dof_id_type Polyhedron::key (const unsigned int s) const
254 {
255  libmesh_assert_less (s, this->n_sides());
256 
257  const Polygon & face = *std::get<0>(this->_sidelinks_data[s]);
258 
259  return face.key();
260 }
261 
262 
263 
264 dof_id_type Polyhedron::low_order_key (const unsigned int s) const
265 {
266  libmesh_assert_less (s, this->n_sides());
267 
268  const Polygon & face = *std::get<0>(this->_sidelinks_data[s]);
269 
270  const unsigned int nv = face.n_vertices();
271  std::vector<dof_id_type> vertex_ids(nv);
272  for (unsigned int v : make_range(nv))
273  vertex_ids[v] = face.node_id(v);
274 
275  return Utility::hashword(vertex_ids);
276 }
277 
278 
279 
280 unsigned int Polyhedron::local_side_node(unsigned int side,
281  unsigned int side_node) const
282 {
283  libmesh_assert_less (side, this->n_sides());
284 
285  const std::vector<unsigned int> & node_map =
286  std::get<2>(this->_sidelinks_data[side]);
287  libmesh_assert_less (side_node, node_map.size());
288 
289  return node_map[side_node];
290 }
291 
292 
293 
294 unsigned int Polyhedron::local_edge_node(unsigned int edge,
295  unsigned int edge_node) const
296 {
297  libmesh_assert_less (edge, this->n_edges());
298  libmesh_assert_less (edge, _edge_lookup.size());
299 
300  auto [side, edge_of_side] = _edge_lookup[edge];
301 
302  const Polygon & face = *std::get<0>(this->_sidelinks_data[side]);
303 
304  const std::vector<unsigned int> & node_map =
305  std::get<2>(this->_sidelinks_data[side]);
306 
307  return node_map[face.local_edge_node(edge_of_side, edge_node)];
308 }
309 
310 
311 
313 {
314  std::vector<dof_id_type> node_ids;
315  for (const auto & n : this->node_ref_range())
316  node_ids.push_back(n.id());
317 
318  return Utility::hashword(node_ids);
319 }
320 
321 
322 
323 std::unique_ptr<Elem> Polyhedron::side_ptr (const unsigned int i)
324 {
325  return const_cast<const Polyhedron *>(this)->side_ptr(i);
326 }
327 
328 
329 
330 std::unique_ptr<Elem> Polyhedron::side_ptr (const unsigned int i) const
331 {
332  libmesh_assert_less (i, this->n_sides());
333 
334  Polygon & face = *std::get<0>(this->_sidelinks_data[i]);
335  std::unique_ptr<Elem> face_copy = face.disconnected_clone();
336  for (auto n : face.node_index_range())
337  face_copy->set_node(n, face.node_ptr(n));
338 
339  return face_copy;
340 }
341 
342 
343 
344 void Polyhedron::side_ptr (std::unique_ptr<Elem> & side,
345  const unsigned int i)
346 {
347  libmesh_assert_less (i, this->n_sides());
348 
349  // Polyhedra are irregular enough that we're not even going to try
350  // and bother optimizing heap access here.
351  side = this->side_ptr(i);
352 }
353 
354 
355 
356 std::unique_ptr<Elem> Polyhedron::build_side_ptr (const unsigned int i)
357 {
358  auto returnval = this->side_ptr(i);
359  returnval->set_interior_parent(this);
360  returnval->inherit_data_from(*this);
361  return returnval;
362 }
363 
364 
365 
366 void Polyhedron::build_side_ptr (std::unique_ptr<Elem> & side,
367  const unsigned int i)
368 {
369  this->side_ptr(side, i);
370  side->set_interior_parent(this);
371  side->inherit_data_from(*this);
372 }
373 
374 
375 
376 std::unique_ptr<Elem> Polyhedron::build_edge_ptr (const unsigned int i)
377 {
378  auto [s, se] = _edge_lookup[i];
379  Polygon & face = *std::get<0>(_sidelinks_data[s]);
380  return face.build_edge_ptr(se);
381 }
382 
383 
384 
385 void Polyhedron::build_edge_ptr (std::unique_ptr<Elem> & elem,
386  const unsigned int i)
387 {
388  auto [s, se] = _edge_lookup[i];
389  Polygon & face = *std::get<0>(_sidelinks_data[s]);
390  face.build_edge_ptr(elem, se);
391 }
392 
393 
394 
395 bool Polyhedron::is_child_on_side(const unsigned int /*c*/,
396  const unsigned int /*s*/) const
397 {
398  libmesh_not_implemented();
399  return false;
400 }
401 
402 
403 
404 unsigned int Polyhedron::opposite_side(const unsigned int /* side_in */) const
405 {
406  // This is too ambiguous in general.
407  libmesh_not_implemented();
408  return libMesh::invalid_uint;
409 }
410 
411 
412 
413 unsigned int Polyhedron::opposite_node(const unsigned int /* n */,
414  const unsigned int /* s */) const
415 {
416  // This is too ambiguous in general.
417  libmesh_not_implemented();
418  return libMesh::invalid_uint;
419 }
420 
421 
422 
424 {
425  if (this->_triangulation.empty())
426  return false;
427 
428  auto & tet = this->_triangulation[0];
429 
430  const Point v01 = this->point(tet[1]) - this->point(tet[0]);
431  const Point v02 = this->point(tet[2]) - this->point(tet[0]);
432  const Point v03 = this->point(tet[3]) - this->point(tet[0]);
433 
434  return (triple_product(v01, v02, v03) < 0);
435 }
436 
437 
438 std::vector<unsigned int>
439 Polyhedron::edges_adjacent_to_node(const unsigned int n) const
440 {
441  libmesh_assert_less(n, this->n_nodes());
442 
443  // For mid-edge or mid-face nodes, the subclass had better have
444  // overridden this.
445  libmesh_assert_less(n, this->n_vertices());
446 
447  const unsigned int ne = this->n_edges();
448  libmesh_assert_equal_to(ne, _edge_lookup.size());
449 
450  std::vector<unsigned int> adjacent_edges;
451 
452  unsigned int next_edge = 0;
453 
454  // Look for any adjacent edges on each side. Make use of the fact
455  // that we number our edges in order as we encounter them from each
456  // side.
457  for (auto t : index_range(_sidelinks_data))
458  {
459  const Polygon & face = *std::get<0>(_sidelinks_data[t]);
460  const std::vector<unsigned int> & node_map =
461  std::get<2>(_sidelinks_data[t]);
462 
463  while (_edge_lookup[next_edge].first < t)
464  {
465  ++next_edge;
466  if (next_edge == ne)
467  return adjacent_edges;
468  }
469 
470  // If we haven't seen the next edge on this or an earlier side
471  // then we might as well go on to the next.
472  if (_edge_lookup[next_edge].first > t)
473  continue;
474 
475  const unsigned int fnv = face.n_vertices();
476  libmesh_assert_equal_to(fnv, face.n_edges());
477  libmesh_assert_equal_to(fnv, face.n_sides());
478  libmesh_assert_less_equal(fnv, node_map.size());
479 
480  // Polygon faces have one edge per vertex
481  for (auto v : make_range(fnv))
482  {
483  libmesh_assert_equal_to (_edge_lookup[next_edge].first, t);
484 
485  if (_edge_lookup[next_edge].second > v)
486  continue;
487 
488  while (_edge_lookup[next_edge].first == t &&
489  _edge_lookup[next_edge].second < v)
490  {
491  ++next_edge;
492  if (next_edge == ne)
493  return adjacent_edges;
494  }
495 
496  if (_edge_lookup[next_edge].first > t)
497  break;
498 
499  const unsigned int vn = node_map[v];
500  const unsigned int vnp = node_map[(v+1)%fnv];
501  libmesh_assert_less(vn, this->n_vertices());
502  libmesh_assert_less(vnp, this->n_vertices());
503  if (vn == n || vnp == n)
504  adjacent_edges.push_back(next_edge);
505  }
506  }
507 
508  return adjacent_edges;
509 }
510 
511 
512 std::pair<Real, Real> Polyhedron::qual_bounds (const ElemQuality q) const
513 {
514  std::pair<Real, Real> bounds;
515 
516  switch (q)
517  {
518  case EDGE_LENGTH_RATIO:
519  bounds.first = 1.;
520  bounds.second = 4.;
521  break;
522 
523  case MIN_ANGLE:
524  bounds.first = 30.;
525  bounds.second = 180.;
526  break;
527 
528  case MAX_ANGLE:
529  bounds.first = 60.;
530  bounds.second = 180.;
531  break;
532 
533  case JACOBIAN:
534  case SCALED_JACOBIAN:
535  bounds.first = 0.5;
536  bounds.second = 1.;
537  break;
538 
539  default:
540  libMesh::out << "Warning: Invalid quality measure chosen." << std::endl;
541  bounds.first = -1;
542  bounds.second = -1;
543  }
544 
545  return bounds;
546 }
547 
548 
549 
550 std::vector<std::shared_ptr<Polygon>>
552 {
553  const auto ns = this->n_sides();
554 
555  libmesh_assert_equal_to(ns, _sidelinks_data.size());
556 
557  std::vector<std::shared_ptr<Polygon>> cloned_sides(ns);
558 
559  for (auto i : make_range(ns))
560  {
561  const Polygon & face = *std::get<0>(this->_sidelinks_data[i]);
562 
563  Elem * clone = face.disconnected_clone().release();
564  Polygon * polygon_clone = cast_ptr<Polygon *>(clone);
565  cloned_sides[i] = std::shared_ptr<Polygon>(polygon_clone);
566 
567  // We can't actually use a *disconnected* clone to reconstruct
568  // links between sides, so we'll temporarily give the clone our
569  // own nodes; user code that typically replaces the usual
570  // nullptr with permanent nodes will then instead place our
571  // nodes with permanent nodes.
572  for (auto n : make_range(face.n_nodes()))
573  cloned_sides[i]->set_node
574  (n, const_cast<Node *>(face.node_ptr(n)));
575  }
576 
577  return cloned_sides;
578 }
579 
580 
581 
583  unsigned int min_node,
584  unsigned int max_node) const
585 {
586  const Polygon & face = *std::get<0>(_sidelinks_data[s]);
587  const std::vector<unsigned int> & node_map =
588  std::get<2>(this->_sidelinks_data[s]);
589 
590  for (unsigned int e : make_range(face.n_sides()))
591  {
592  std::vector<unsigned int> nodes_on_edge =
593  face.nodes_on_side(e);
594  libmesh_assert_equal_to(nodes_on_edge.size(), 2);
595  nodes_on_edge[0] = node_map[nodes_on_edge[0]];
596  nodes_on_edge[1] = node_map[nodes_on_edge[1]];
597  if ((nodes_on_edge[0] == min_node) &&
598  (nodes_on_edge[1] == max_node))
599  return true;
600  if ((nodes_on_edge[1] == min_node) &&
601  (nodes_on_edge[0] == max_node))
602  return true;
603  }
604 
605  return false;
606 }
607 
608 
609 
610 std::vector<unsigned int>
611 Polyhedron::sides_on_edge(const unsigned int e) const
612 {
613  std::vector<unsigned int> returnval(2);
614  auto [s1, s1e] = _edge_lookup[e];
615  returnval[0] = s1;
616 
617  const Polygon & face1 = *std::get<0>(_sidelinks_data[s1]);
618  const std::vector<unsigned int> & node_map =
619  std::get<2>(this->_sidelinks_data[s1]);
620 
621  std::vector<unsigned int> nodes_on_edge =
622  face1.nodes_on_side(s1e);
623  libmesh_assert_equal_to(nodes_on_edge.size(), 2);
624  nodes_on_edge[0] = node_map[nodes_on_edge[0]];
625  nodes_on_edge[1] = node_map[nodes_on_edge[1]];
626 
627  if (nodes_on_edge[0] > nodes_on_edge[1])
628  std::swap(nodes_on_edge[0], nodes_on_edge[1]);
629 
630  for (unsigned int s2 : make_range(this->n_sides()))
631  {
632  if (s2 == s1)
633  continue;
634 
635  if (this->side_has_edge_nodes(s2, nodes_on_edge[0],
636  nodes_on_edge[1]))
637  {
638  returnval[1] = s2;
639  return returnval;
640  }
641  }
642 
643  libmesh_error();
644 
645  return returnval;
646 }
647 
648 
649 
650 bool Polyhedron::is_edge_on_side(const unsigned int e,
651  const unsigned int s) const
652 {
653  auto [s1, s1e] = _edge_lookup[e];
654 
655  // Did we get lucky with our cache?
656  if (s1 == s)
657  return true;
658 
659  const Polygon & face1 = *std::get<0>(_sidelinks_data[s1]);
660  const std::vector<unsigned int> & node_map =
661  std::get<2>(this->_sidelinks_data[s1]);
662  std::vector<unsigned int> nodes_on_edge1 =
663  face1.nodes_on_side(s1e);
664  libmesh_assert_equal_to(nodes_on_edge1.size(), 2);
665 
666  nodes_on_edge1[0] = node_map[nodes_on_edge1[0]];
667  nodes_on_edge1[1] = node_map[nodes_on_edge1[1]];
668  if (nodes_on_edge1[0] > nodes_on_edge1[1])
669  std::swap(nodes_on_edge1[0], nodes_on_edge1[1]);
670 
671  return this->side_has_edge_nodes(s,
672  nodes_on_edge1[0],
673  nodes_on_edge1[1]);
674 }
675 
676 
677 
678 std::array<Point, 4> Polyhedron::master_subelement (unsigned int i) const
679 {
680  libmesh_assert_less(i, this->_triangulation.size());
681 
682  const auto & tet = this->_triangulation[i];
683 
684  return { this->master_point(tet[0]),
685  this->master_point(tet[1]),
686  this->master_point(tet[2]),
687  this->master_point(tet[3]) };
688 }
689 
690 
691 
692 std::tuple<unsigned int, Real, Real, Real>
694 {
695  std::tuple<unsigned int, Real, Real, Real> returnval =
696  {libMesh::invalid_uint, -1, -1, -1};
697 
698  Real best_bad_coord = -1;
699 
700  for (auto s : make_range(this->n_subelements()))
701  {
702  const std::array<Point, 4> subtet =
703  this->master_subelement(s);
704 
705  // Find barycentric coordinates in subelem
706  const Point v0 = p - subtet[0];
707  // const Point v1 = p - subtet[1];
708 
709  const Point v01 = subtet[1] - subtet[0];
710  const Point v02 = subtet[2] - subtet[0];
711  const Point v03 = subtet[3] - subtet[0];
712 
713  // const Point v12 = subtet[2] - subtet[1];
714  // const Point v13 = subtet[3] - subtet[1];
715 
716  // const Real tp0 = triple_product(v1, v13, v12);
717  const Real tp1 = triple_product(v0, v02, v03);
718  const Real tp2 = triple_product(v0, v03, v01);
719  const Real tp3 = triple_product(v0, v01, v02);
720 
721  const Real six_vol = triple_product(v01, v02, v03);
722 
723  const Real xi = tp1 / six_vol;
724  const Real eta = tp2 / six_vol;
725  const Real zeta = tp3 / six_vol;
726 
727  if (xi>=0 && eta>=0 && zeta>=0 && xi+eta+zeta<=1)
728  return { s, xi, eta, zeta };
729 
730  const Real my_best_bad_coord =
731  std::min(std::min(std::min(xi, eta), zeta), 1-xi-eta-zeta);
732 
733  if (my_best_bad_coord > best_bad_coord)
734  {
735  best_bad_coord = my_best_bad_coord;
736  returnval = { s, xi, eta, zeta };
737  }
738  }
739 
740  if (best_bad_coord > -tol)
741  return returnval;
742 
743  return {libMesh::invalid_uint, -1, -1, -1};
744 }
745 
746 
747 } // namespace libMesh
void set_p_level(const unsigned int p)
Sets the value of the p-refinement level for the element.
unsigned char mapping_data() const
Definition: elem.h:3150
virtual unsigned int n_nodes() const override
virtual unsigned int n_nodes() const override
Definition: face_polygon.h:78
unsigned char _map_type
Mapping function type; currently either 0 (LAGRANGE) or 1 (RATIONAL_BERNSTEIN).
Definition: elem.h:2303
unsigned char _map_data
Mapping function data; currently used when needed to store the RATIONAL_BERNSTEIN nodal weight data i...
Definition: elem.h:2309
std::vector< Elem * > _elemlinks_data
Data for links to parent/neighbor/interior_parent elements.
virtual Node *& set_node(const unsigned int i)
Definition: elem.h:2564
A Node is like a Point, but with more information.
Definition: node.h:52
Node ** _nodes
Pointers to the nodes we are connected to.
Definition: elem.h:2251
The Polyhedron is an element in 3D with an arbitrary number of polygonal faces.
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const override
auto norm() const
Definition: type_vector.h:907
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
Definition: libmesh.h:303
std::vector< Node * > _nodelinks_data
Data for links to nodes.
virtual unsigned int opposite_side(const unsigned int s) const override final
Throws an error.
virtual dof_id_type key() const override
Definition: face_polygon.C:200
virtual unsigned int local_edge_node(unsigned int edge, unsigned int edge_node) const override
Similar to Elem::local_side_node(), but instead of a side id, takes an edge id and a node id on that ...
static constexpr Real TOLERANCE
virtual std::array< Point, 4 > master_subelement(unsigned int i) const
std::vector< std::tuple< std::shared_ptr< Polygon >, bool, std::vector< unsigned int > > > _sidelinks_data
Data for links to sides.
virtual std::pair< Real, Real > qual_bounds(const ElemQuality q) const override
virtual bool on_reference_element(const Point &p, const Real eps=TOLERANCE) const override final
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
static const int num_children
unsigned int p_level() const
Definition: elem.h:3122
Polyhedron(const std::vector< std::shared_ptr< Polygon >> &sides, Elem *p)
Arbitrary polyhedral element, takes a vector of shared pointers to sides (which should already be con...
The libMesh namespace provides an interface to certain functionality in the library.
virtual unsigned int local_edge_node(unsigned int edge, unsigned int edge_node) const override
Calls local_side_node(edge, edge_node).
Definition: face_polygon.C:192
virtual unsigned int n_sides() const override final
virtual std::unique_ptr< Elem > build_edge_ptr(const unsigned int i) override final
void set_interior_parent(Elem *p)
Sets the pointer to the element&#39;s interior_parent.
Definition: elem.C:1218
void add(const TypeVector< T2 > &)
Add to this vector without creating a temporary.
Definition: type_vector.h:605
std::vector< std::pair< unsigned int, unsigned int > > _edge_lookup
One entry for each polyhedron edge, a pair indicating the side number and the edge-of-side number whi...
virtual unsigned int n_edges() const override final
Definition: face_polygon.h:94
uint32_t hashword(const uint32_t *k, size_t length, uint32_t initval=0)
The hashword function takes an array of uint32_t&#39;s of length &#39;length&#39; and computes a single key from ...
Definition: hashword.h:158
bool side_has_edge_nodes(unsigned int side, unsigned int min_node, unsigned int max_node) const
Helper method for finding the non-cached side that shares an edge, by examining the local node ids th...
T triple_product(const TypeVector< T > &a, const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1031
TypeVector< T > unit() const
Definition: type_vector.h:1109
std::tuple< unsigned int, Real, Real, Real > subelement_coordinates(const Point &p, Real tol=TOLERANCE *TOLERANCE) const
ElemMappingType mapping_type() const
Definition: elem.h:3134
virtual bool is_flipped() const override final
virtual bool is_edge_on_side(const unsigned int e, const unsigned int s) const override final
unsigned int n_subelements() const
virtual dof_id_type low_order_key(const unsigned int s) const override
The Polygon is an element in 2D with an arbitrary (but fixed) number of sides.
Definition: face_polygon.h:39
Elem ** _elemlinks
Pointers to this element&#39;s parent and neighbors, and for lower-dimensional elements&#39; interior_parent...
Definition: elem.h:2257
The Cell is an abstract element type that lives in three dimensions.
Definition: cell.h:38
libmesh_assert(ctx)
virtual unsigned int opposite_node(const unsigned int n, const unsigned int s) const override final
Throws an error - opposite_side(s) is too hard to define in general on polyhedra. ...
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i) override
Copies the Polygon side coincident with side i.
ElemQuality
Defines an enum for element quality metrics.
virtual std::unique_ptr< Elem > disconnected_clone() const
Definition: elem.C:410
virtual dof_id_type key() const override
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int) const =0
SimpleRange< NodeRefIter > node_ref_range()
Returns a range with all nodes of an element, usable in range-based for loops.
Definition: elem.h:2679
std::vector< std::shared_ptr< Polygon > > side_clones() const
virtual Point master_point(const unsigned int i) const override
virtual unsigned int n_edges() const override final
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const override
virtual unsigned int n_vertices() const =0
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual unsigned int n_vertices() const override final
Definition: face_polygon.h:89
subdomain_id_type subdomain_id() const
Definition: elem.h:2588
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:2513
OStreamProxy out
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:173
virtual std::vector< unsigned int > edges_adjacent_to_node(const unsigned int n) const override
IntRange< unsigned short > node_index_range() const
Definition: elem.h:2697
virtual unsigned int n_sides() const override final
Definition: face_polygon.h:84
virtual std::unique_ptr< Elem > side_ptr(const unsigned int i) override final
std::vector< std::array< int, 4 > > _triangulation
Data for a triangulation (tetrahedralization) of the polyhedron.
processor_id_type processor_id() const
Definition: dof_object.h:881
virtual std::vector< unsigned int > sides_on_edge(const unsigned int e) const override final
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
dof_id_type node_id(const unsigned int i) const
Definition: elem.h:2481
const Point & point(const unsigned int i) const
Definition: elem.h:2459
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:150
virtual std::vector< unsigned int > nodes_on_side(const unsigned int) const =0
uint8_t dof_id_type
Definition: id_types.h:67
virtual std::unique_ptr< Elem > build_edge_ptr(const unsigned int i) override final
build_side and build_edge are identical for faces.
Definition: face.h:72