libMesh
cell_polyhedron.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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  // Do the manual initialization that Elem::Elem and Cell::Cell
97  // couldn't, now that we've resized both our vectors. No need to
98  // manually set nullptr, though, since std::vector does that.
99  this->_elemlinks = _elemlinks_data.data();
100  this->_nodes = _nodelinks_data.data();
101  this->_elemlinks[0] = p;
102 
103  libmesh_assert_equal_to(nn, this->n_nodes());
104 
105  // Figure out the orientation of our sides, now that we've got our
106  // nodes organized enough to find our center. The algorithm below
107  // only works for convex polyhedra, but that's all we're
108  // supporting for now.
109  Point center;
110  for (auto n : make_range(nn))
111  center.add (this->point(n));
112  center /= static_cast<Real>(nn);
113 
114  for (unsigned int s : index_range(sides))
115  {
116  const Polygon & side = *sides[s];
117  const Point x_i = side.point(0);
118  const Point n_i =
119  (side.point(1) - side.point(0)).cross
120  (side.point(0) - side.point(side.n_sides()-1)).unit();
121 
122  bool & inward_normal = std::get<1>(_sidelinks_data[s]);
123  inward_normal = (n_i * (center - x_i) > TOLERANCE);
124  }
125 
126  // We're betting a lot on "our polyhedra are all convex", so let's
127  // check that if we have time.
128 #ifdef DEBUG
129  for (unsigned int s : index_range(sides))
130  {
131  const Polygon & side = *sides[s];
132  const Point x_i = side.point(0);
133  const bool inward_normal = std::get<1>(this->_sidelinks_data[s]);
134 
135  const Point n_i =
136  (side.point(1) - side.point(0)).cross
137  (side.point(0) - side.point(side.n_sides()-1)).unit() *
138  (inward_normal ? -1 : 1);
139 
140  for (const Point & node : this->node_ref_range())
141  {
142  const Point d_n = node - x_i;
143  if (d_n * n_i > TOLERANCE * d_n.norm())
144  libmesh_not_implemented_msg
145  ("Cannot create a non-convex polyhedron");
146  }
147  }
148 #endif
149 
150  // Is this likely to ever be used? We may do refinement with
151  // polyhedra but it's probably not going to have a hierarchy...
152  if (p)
153  {
154  this->subdomain_id() = p->subdomain_id();
155  this->processor_id() = p->processor_id();
156  _map_type = p->mapping_type();
157  _map_data = p->mapping_data();
158 
159 #ifdef LIBMESH_ENABLE_AMR
160  this->set_p_level(p->p_level());
161 #endif
162  }
163 
164  // Make sure the interior parent isn't undefined
165  this->set_interior_parent(nullptr);
166  }
167 
168 
169 
170 Point Polyhedron::master_point (const unsigned int i) const
171 {
172  return this->point(i);
173 }
174 
175 
176 
178 {
179  for (unsigned int s : make_range(this->n_sides()))
180  {
181  const Polygon & side = *std::get<0>(this->_sidelinks_data[s]);
182  const Point x_i = side.point(0);
183  const bool inward_normal = std::get<1>(this->_sidelinks_data[s]);
184 
185  const Point n_i =
186  (side.point(1) - side.point(0)).cross
187  (side.point(0) - side.point(side.n_sides()-1)).unit() *
188  (inward_normal ? -1 : 1);
189 
190  for (const Point & node : this->node_ref_range())
191  {
192  const Point d_n = node - x_i;
193  if (d_n * n_i > TOLERANCE * d_n.norm())
194  return false;
195  }
196  }
197  return true;
198 }
199 
200 
201 
203  const Real eps) const
204 {
205  const unsigned int ns = this->n_sides();
206 
207  // Check that the point is on the same side of all the faces by
208  // testing whether:
209  //
210  // n_i.(p - x_i) <= 0
211  //
212  // for each i, where:
213  // n_i is the outward normal of face i,
214  // x_i is a point on face i.
215 
216  for (auto i : make_range(ns))
217  {
218  const Polygon & face = *std::get<0>(this->_sidelinks_data[i]);
219  const bool inward_normal = std::get<1>(this->_sidelinks_data[i]);
220 
221  const Point x_i = face.point(0);
222 
223  const Point n_i =
224  (face.point(1) - face.point(0)).cross
225  (face.point(0) - face.point(face.n_sides()-1)).unit() *
226  (inward_normal ? -1 : 1);
227 
228  // This only works for polyhedra with flat sides.
229 #ifdef DEBUG
230  for (auto j : make_range(face.n_sides()-1))
231  {
232  const Point x_j = face.point(j+1);
233  const Point d_j = x_j - x_i;
234  if (std::abs(d_j * n_i) > eps * d_j.norm())
235  libmesh_not_implemented_msg
236  ("Polyhedra with non-flat sides are not fully supported.");
237  }
238 #endif
239 
240  if (n_i * (p - x_i) > eps)
241  return false;
242  }
243 
244  return true;
245 }
246 
247 
248 
249 dof_id_type Polyhedron::key (const unsigned int s) const
250 {
251  libmesh_assert_less (s, this->n_sides());
252 
253  const Polygon & face = *std::get<0>(this->_sidelinks_data[s]);
254 
255  return face.key();
256 }
257 
258 
259 
260 dof_id_type Polyhedron::low_order_key (const unsigned int s) const
261 {
262  libmesh_assert_less (s, this->n_sides());
263 
264  const Polygon & face = *std::get<0>(this->_sidelinks_data[s]);
265 
266  const unsigned int nv = face.n_vertices();
267  std::vector<dof_id_type> vertex_ids(nv);
268  for (unsigned int v : make_range(nv))
269  vertex_ids[v] = face.node_id(v);
270 
271  return Utility::hashword(vertex_ids);
272 }
273 
274 
275 
276 unsigned int Polyhedron::local_side_node(unsigned int side,
277  unsigned int side_node) const
278 {
279  libmesh_assert_less (side, this->n_sides());
280 
281  const std::vector<unsigned int> & node_map =
282  std::get<2>(this->_sidelinks_data[side]);
283  libmesh_assert_less (side_node, node_map.size());
284 
285  return node_map[side_node];
286 }
287 
288 
289 
290 unsigned int Polyhedron::local_edge_node(unsigned int edge,
291  unsigned int edge_node) const
292 {
293  libmesh_assert_less (edge, this->n_edges());
294  libmesh_assert_less (edge, _edge_lookup.size());
295 
296  auto [side, edge_of_side] = _edge_lookup[edge];
297 
298  const Polygon & face = *std::get<0>(this->_sidelinks_data[side]);
299 
300  const std::vector<unsigned int> & node_map =
301  std::get<2>(this->_sidelinks_data[side]);
302 
303  return node_map[face.local_edge_node(edge_of_side, edge_node)];
304 }
305 
306 
307 
309 {
310  std::vector<dof_id_type> node_ids;
311  for (const auto & n : this->node_ref_range())
312  node_ids.push_back(n.id());
313 
314  return Utility::hashword(node_ids);
315 }
316 
317 
318 
319 std::unique_ptr<Elem> Polyhedron::side_ptr (const unsigned int i)
320 {
321  return const_cast<const Polyhedron *>(this)->side_ptr(i);
322 }
323 
324 
325 
326 std::unique_ptr<Elem> Polyhedron::side_ptr (const unsigned int i) const
327 {
328  libmesh_assert_less (i, this->n_sides());
329 
330  Polygon & face = *std::get<0>(this->_sidelinks_data[i]);
331  std::unique_ptr<Elem> face_copy = face.disconnected_clone();
332  for (auto n : face.node_index_range())
333  face_copy->set_node(n, face.node_ptr(n));
334 
335  return face_copy;
336 }
337 
338 
339 
340 void Polyhedron::side_ptr (std::unique_ptr<Elem> & side,
341  const unsigned int i)
342 {
343  libmesh_assert_less (i, this->n_sides());
344 
345  // Polyhedra are irregular enough that we're not even going to try
346  // and bother optimizing heap access here.
347  side = this->side_ptr(i);
348 }
349 
350 
351 
352 std::unique_ptr<Elem> Polyhedron::build_side_ptr (const unsigned int i)
353 {
354  auto returnval = this->side_ptr(i);
355  returnval->set_interior_parent(this);
356  returnval->inherit_data_from(*this);
357  return returnval;
358 }
359 
360 
361 
362 void Polyhedron::build_side_ptr (std::unique_ptr<Elem> & side,
363  const unsigned int i)
364 {
365  this->side_ptr(side, i);
366  side->set_interior_parent(this);
367  side->inherit_data_from(*this);
368 }
369 
370 
371 
372 std::unique_ptr<Elem> Polyhedron::build_edge_ptr (const unsigned int i)
373 {
374  auto [s, se] = _edge_lookup[i];
375  Polygon & face = *std::get<0>(_sidelinks_data[s]);
376  return face.build_edge_ptr(se);
377 }
378 
379 
380 
381 void Polyhedron::build_edge_ptr (std::unique_ptr<Elem> & elem,
382  const unsigned int i)
383 {
384  auto [s, se] = _edge_lookup[i];
385  Polygon & face = *std::get<0>(_sidelinks_data[s]);
386  face.build_edge_ptr(elem, se);
387 }
388 
389 
390 
391 bool Polyhedron::is_child_on_side(const unsigned int /*c*/,
392  const unsigned int /*s*/) const
393 {
394  libmesh_not_implemented();
395  return false;
396 }
397 
398 
399 
400 unsigned int Polyhedron::opposite_side(const unsigned int /* side_in */) const
401 {
402  // This is too ambiguous in general.
403  libmesh_not_implemented();
404  return libMesh::invalid_uint;
405 }
406 
407 
408 
409 unsigned int Polyhedron::opposite_node(const unsigned int /* n */,
410  const unsigned int /* s */) const
411 {
412  // This is too ambiguous in general.
413  libmesh_not_implemented();
414  return libMesh::invalid_uint;
415 }
416 
417 
418 
420 {
421  if (this->_triangulation.empty())
422  return false;
423 
424  auto & tet = this->_triangulation[0];
425 
426  const Point v01 = this->point(tet[1]) - this->point(tet[0]);
427  const Point v02 = this->point(tet[2]) - this->point(tet[0]);
428  const Point v03 = this->point(tet[3]) - this->point(tet[0]);
429 
430  return (triple_product(v01, v02, v03) < 0);
431 }
432 
433 
434 std::vector<unsigned int>
435 Polyhedron::edges_adjacent_to_node(const unsigned int n) const
436 {
437  libmesh_assert_less(n, this->n_nodes());
438 
439  // For mid-edge or mid-face nodes, the subclass had better have
440  // overridden this.
441  libmesh_assert_less(n, this->n_vertices());
442 
443  const unsigned int ne = this->n_edges();
444  libmesh_assert_equal_to(ne, _edge_lookup.size());
445 
446  std::vector<unsigned int> adjacent_edges;
447 
448  unsigned int next_edge = 0;
449 
450  // Look for any adjacent edges on each side. Make use of the fact
451  // that we number our edges in order as we encounter them from each
452  // side.
453  for (auto t : index_range(_sidelinks_data))
454  {
455  const Polygon & face = *std::get<0>(_sidelinks_data[t]);
456  const std::vector<unsigned int> & node_map =
457  std::get<2>(_sidelinks_data[t]);
458 
459  while (_edge_lookup[next_edge].first < t)
460  {
461  ++next_edge;
462  if (next_edge == ne)
463  return adjacent_edges;
464  }
465 
466  // If we haven't seen the next edge on this or an earlier side
467  // then we might as well go on to the next.
468  if (_edge_lookup[next_edge].first > t)
469  continue;
470 
471  const unsigned int fnv = face.n_vertices();
472  libmesh_assert_equal_to(fnv, face.n_edges());
473  libmesh_assert_equal_to(fnv, face.n_sides());
474  libmesh_assert_less_equal(fnv, node_map.size());
475 
476  // Polygon faces have one edge per vertex
477  for (auto v : make_range(fnv))
478  {
479  libmesh_assert_equal_to (_edge_lookup[next_edge].first, t);
480 
481  if (_edge_lookup[next_edge].second > v)
482  continue;
483 
484  while (_edge_lookup[next_edge].first == t &&
485  _edge_lookup[next_edge].second < v)
486  {
487  ++next_edge;
488  if (next_edge == ne)
489  return adjacent_edges;
490  }
491 
492  if (_edge_lookup[next_edge].first > t)
493  break;
494 
495  const unsigned int vn = node_map[v];
496  const unsigned int vnp = node_map[(v+1)%fnv];
497  libmesh_assert_less(vn, this->n_vertices());
498  libmesh_assert_less(vnp, this->n_vertices());
499  if (vn == n || vnp == n)
500  adjacent_edges.push_back(next_edge);
501  }
502  }
503 
504  return adjacent_edges;
505 }
506 
507 
508 std::pair<Real, Real> Polyhedron::qual_bounds (const ElemQuality q) const
509 {
510  std::pair<Real, Real> bounds;
511 
512  switch (q)
513  {
514  case EDGE_LENGTH_RATIO:
515  bounds.first = 1.;
516  bounds.second = 4.;
517  break;
518 
519  case MIN_ANGLE:
520  bounds.first = 30.;
521  bounds.second = 180.;
522  break;
523 
524  case MAX_ANGLE:
525  bounds.first = 60.;
526  bounds.second = 180.;
527  break;
528 
529  case JACOBIAN:
530  case SCALED_JACOBIAN:
531  bounds.first = 0.5;
532  bounds.second = 1.;
533  break;
534 
535  default:
536  libMesh::out << "Warning: Invalid quality measure chosen." << std::endl;
537  bounds.first = -1;
538  bounds.second = -1;
539  }
540 
541  return bounds;
542 }
543 
544 
545 
546 std::vector<std::shared_ptr<Polygon>>
548 {
549  const auto ns = this->n_sides();
550 
551  libmesh_assert_equal_to(ns, _sidelinks_data.size());
552 
553  std::vector<std::shared_ptr<Polygon>> cloned_sides(ns);
554 
555  for (auto i : make_range(ns))
556  {
557  const Polygon & face = *std::get<0>(this->_sidelinks_data[i]);
558 
559  Elem * clone = face.disconnected_clone().release();
560  Polygon * polygon_clone = cast_ptr<Polygon *>(clone);
561  cloned_sides[i] = std::shared_ptr<Polygon>(polygon_clone);
562 
563  // We can't actually use a *disconnected* clone to reconstruct
564  // links between sides, so we'll temporarily give the clone our
565  // own nodes; user code that typically replaces the usual
566  // nullptr with permanent nodes will then instead place our
567  // nodes with permanent nodes.
568  for (auto n : make_range(face.n_nodes()))
569  cloned_sides[i]->set_node
570  (n, const_cast<Node *>(face.node_ptr(n)));
571  }
572 
573  return cloned_sides;
574 }
575 
576 
577 
579  unsigned int min_node,
580  unsigned int max_node) const
581 {
582  const Polygon & face = *std::get<0>(_sidelinks_data[s]);
583  const std::vector<unsigned int> & node_map =
584  std::get<2>(this->_sidelinks_data[s]);
585 
586  for (unsigned int e : make_range(face.n_sides()))
587  {
588  std::vector<unsigned int> nodes_on_edge =
589  face.nodes_on_side(e);
590  libmesh_assert_equal_to(nodes_on_edge.size(), 2);
591  nodes_on_edge[0] = node_map[nodes_on_edge[0]];
592  nodes_on_edge[1] = node_map[nodes_on_edge[1]];
593  if ((nodes_on_edge[0] == min_node) &&
594  (nodes_on_edge[1] == max_node))
595  return true;
596  if ((nodes_on_edge[1] == min_node) &&
597  (nodes_on_edge[0] == max_node))
598  return true;
599  }
600 
601  return false;
602 }
603 
604 
605 
606 std::vector<unsigned int>
607 Polyhedron::sides_on_edge(const unsigned int e) const
608 {
609  std::vector<unsigned int> returnval(2);
610  auto [s1, s1e] = _edge_lookup[e];
611  returnval[0] = s1;
612 
613  const Polygon & face1 = *std::get<0>(_sidelinks_data[s1]);
614  const std::vector<unsigned int> & node_map =
615  std::get<2>(this->_sidelinks_data[s1]);
616 
617  std::vector<unsigned int> nodes_on_edge =
618  face1.nodes_on_side(s1e);
619  libmesh_assert_equal_to(nodes_on_edge.size(), 2);
620  nodes_on_edge[0] = node_map[nodes_on_edge[0]];
621  nodes_on_edge[1] = node_map[nodes_on_edge[1]];
622 
623  if (nodes_on_edge[0] > nodes_on_edge[1])
624  std::swap(nodes_on_edge[0], nodes_on_edge[1]);
625 
626  for (unsigned int s2 : make_range(this->n_sides()))
627  {
628  if (s2 == s1)
629  continue;
630 
631  if (this->side_has_edge_nodes(s2, nodes_on_edge[0],
632  nodes_on_edge[1]))
633  {
634  returnval[1] = s2;
635  return returnval;
636  }
637  }
638 
639  libmesh_error();
640 
641  return returnval;
642 }
643 
644 
645 
646 bool Polyhedron::is_edge_on_side(const unsigned int e,
647  const unsigned int s) const
648 {
649  auto [s1, s1e] = _edge_lookup[e];
650 
651  // Did we get lucky with our cache?
652  if (s1 == s)
653  return true;
654 
655  const Polygon & face1 = *std::get<0>(_sidelinks_data[s1]);
656  const std::vector<unsigned int> & node_map =
657  std::get<2>(this->_sidelinks_data[s1]);
658  std::vector<unsigned int> nodes_on_edge1 =
659  face1.nodes_on_side(s1e);
660  libmesh_assert_equal_to(nodes_on_edge1.size(), 2);
661 
662  nodes_on_edge1[0] = node_map[nodes_on_edge1[0]];
663  nodes_on_edge1[1] = node_map[nodes_on_edge1[1]];
664  if (nodes_on_edge1[0] > nodes_on_edge1[1])
665  std::swap(nodes_on_edge1[0], nodes_on_edge1[1]);
666 
667  return this->side_has_edge_nodes(s,
668  nodes_on_edge1[0],
669  nodes_on_edge1[1]);
670 }
671 
672 
673 
674 std::array<Point, 4> Polyhedron::master_subelement (unsigned int i) const
675 {
676  libmesh_assert_less(i, this->_triangulation.size());
677 
678  const auto & tet = this->_triangulation[i];
679 
680  return { this->master_point(tet[0]),
681  this->master_point(tet[1]),
682  this->master_point(tet[2]),
683  this->master_point(tet[3]) };
684 }
685 
686 
687 std::tuple<unsigned int, Real, Real, Real>
689 {
690  std::tuple<unsigned int, Real, Real, Real> returnval =
691  {libMesh::invalid_uint, -1, -1, -1};
692 
693  Real best_bad_coord = -1;
694 
695  for (auto s : make_range(this->n_subelements()))
696  {
697  const std::array<Point, 4> subtet =
698  this->master_subelement(s);
699 
700  // Find barycentric coordinates in subelem
701  const Point v0 = p - subtet[0];
702  // const Point v1 = p - subtet[1];
703 
704  const Point v01 = subtet[1] - subtet[0];
705  const Point v02 = subtet[2] - subtet[0];
706  const Point v03 = subtet[3] - subtet[0];
707 
708  // const Point v12 = subtet[2] - subtet[1];
709  // const Point v13 = subtet[3] - subtet[1];
710 
711  // const Real tp0 = triple_product(v1, v13, v12);
712  const Real tp1 = triple_product(v0, v02, v03);
713  const Real tp2 = triple_product(v0, v03, v01);
714  const Real tp3 = triple_product(v0, v01, v02);
715 
716  const Real six_vol = triple_product(v01, v02, v03);
717 
718  const Real xi = tp1 / six_vol;
719  const Real eta = tp2 / six_vol;
720  const Real zeta = tp3 / six_vol;
721 
722  if (xi>=0 && eta>=0 && zeta>=0 && xi+eta+zeta<=1)
723  return { s, xi, eta, zeta };
724 
725  const Real my_best_bad_coord =
726  std::min(std::min(std::min(xi, eta), zeta), 1-xi-eta-zeta);
727 
728  if (my_best_bad_coord > best_bad_coord)
729  {
730  best_bad_coord = my_best_bad_coord;
731  returnval = { s, xi, eta, zeta };
732  }
733  }
734 
735  if (best_bad_coord > -tol)
736  return returnval;
737 
738  return {libMesh::invalid_uint, -1, -1, -1};
739 }
740 
741 
742 } // 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:3137
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:2298
unsigned char _map_data
Mapping function data; currently used when needed to store the RATIONAL_BERNSTEIN nodal weight data i...
Definition: elem.h:2304
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:2559
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:2246
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 -> decltype(std::norm(T()))
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:310
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:3109
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:1246
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:1029
TypeVector< T > unit() const
Definition: type_vector.h:1104
std::tuple< unsigned int, Real, Real, Real > subelement_coordinates(const Point &p, Real tol=TOLERANCE *TOLERANCE) const
ElemMappingType mapping_type() const
Definition: elem.h:3121
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:2252
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:2666
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:2583
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:2508
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:140
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:2684
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:905
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:2476
const Point & point(const unsigned int i) const
Definition: elem.h:2454
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:117
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