libMesh
elem.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 
19 // Local includes
20 #include "libmesh/elem.h"
21 
22 #include "libmesh/boundary_info.h"
23 #include "libmesh/fe_type.h"
24 #include "libmesh/fe_interface.h"
25 #include "libmesh/node_elem.h"
26 #include "libmesh/edge_edge2.h"
27 #include "libmesh/edge_edge3.h"
28 #include "libmesh/edge_edge4.h"
29 #include "libmesh/edge_inf_edge2.h"
30 #include "libmesh/face_c0polygon.h"
31 #include "libmesh/face_tri3.h"
32 #include "libmesh/face_tri3_subdivision.h"
33 #include "libmesh/face_tri3_shell.h"
34 #include "libmesh/face_tri6.h"
35 #include "libmesh/face_tri7.h"
36 #include "libmesh/face_quad4.h"
37 #include "libmesh/face_quad4_shell.h"
38 #include "libmesh/face_quad8.h"
39 #include "libmesh/face_quad8_shell.h"
40 #include "libmesh/face_quad9.h"
41 #include "libmesh/face_quad9_shell.h"
42 #include "libmesh/face_inf_quad4.h"
43 #include "libmesh/face_inf_quad6.h"
44 #include "libmesh/cell_tet4.h"
45 #include "libmesh/cell_tet10.h"
46 #include "libmesh/cell_tet14.h"
47 #include "libmesh/cell_hex8.h"
48 #include "libmesh/cell_hex20.h"
49 #include "libmesh/cell_hex27.h"
50 #include "libmesh/cell_inf_hex8.h"
51 #include "libmesh/cell_inf_hex16.h"
52 #include "libmesh/cell_inf_hex18.h"
53 #include "libmesh/cell_prism6.h"
54 #include "libmesh/cell_prism15.h"
55 #include "libmesh/cell_prism18.h"
56 #include "libmesh/cell_prism20.h"
57 #include "libmesh/cell_prism21.h"
58 #include "libmesh/cell_inf_prism6.h"
59 #include "libmesh/cell_inf_prism12.h"
60 #include "libmesh/cell_pyramid5.h"
61 #include "libmesh/cell_pyramid13.h"
62 #include "libmesh/cell_pyramid14.h"
63 #include "libmesh/cell_pyramid18.h"
64 #include "libmesh/fe_base.h"
65 #include "libmesh/mesh_base.h"
66 #include "libmesh/quadrature_nodal.h"
67 #include "libmesh/quadrature_gauss.h"
68 #include "libmesh/remote_elem.h"
69 #include "libmesh/reference_elem.h"
70 #include "libmesh/enum_to_string.h"
71 #include "libmesh/threads.h"
72 #include "libmesh/enum_elem_quality.h"
73 #include "libmesh/enum_io_package.h"
74 #include "libmesh/enum_order.h"
75 #include "libmesh/elem_internal.h"
76 
77 #ifdef LIBMESH_ENABLE_PERIODIC
78 #include "libmesh/mesh.h"
79 #include "libmesh/periodic_boundaries.h"
80 #endif
81 
82 
83 // C++ includes
84 #include <algorithm> // for std::sort
85 #include <array>
86 #include <iterator> // for std::ostream_iterator
87 #include <sstream>
88 #include <limits> // for std::numeric_limits<>
89 #include <cmath> // for std::sqrt()
90 #include <memory>
91 #include <regex> // for exceptions in volume()
92 
93 
94 namespace libMesh
95 {
96 
99 
100 const subdomain_id_type Elem::invalid_subdomain_id = std::numeric_limits<subdomain_id_type>::max();
101 
102 // Initialize static member variables
103 const unsigned int Elem::type_to_dim_map [] =
104  {
105  1, // EDGE2
106  1, // EDGE3
107  1, // EDGE4
108 
109  2, // TRI3
110  2, // TRI6
111 
112  2, // QUAD4
113  2, // QUAD8
114  2, // QUAD9
115 
116  3, // TET4
117  3, // TET10
118 
119  3, // HEX8
120  3, // HEX20
121  3, // HEX27
122 
123  3, // PRISM6
124  3, // PRISM15
125  3, // PRISM18
126 
127  3, // PYRAMID5
128  3, // PYRAMID13
129  3, // PYRAMID14
130 
131  1, // INFEDGE2
132 
133  2, // INFQUAD4
134  2, // INFQUAD6
135 
136  3, // INFHEX8
137  3, // INFHEX16
138  3, // INFHEX18
139 
140  3, // INFPRISM6
141  3, // INFPRISM12
142 
143  0, // NODEELEM
144 
145  0, // REMOTEELEM
146 
147  2, // TRI3SUBDIVISION
148  2, // TRISHELL3
149  2, // QUADSHELL4
150  2, // QUADSHELL8
151 
152  2, // TRI7
153  3, // TET14
154  3, // PRISM20
155  3, // PRISM21
156  3, // PYRAMID18
157 
158  2, // QUADSHELL9
159 
160  2, // C0POLYGON
161  3, // C0POLYHEDRON
162  };
163 
164 const unsigned int Elem::max_n_nodes;
165 
166 const unsigned int Elem::type_to_n_nodes_map [] =
167  {
168  2, // EDGE2
169  3, // EDGE3
170  4, // EDGE4
171 
172  3, // TRI3
173  6, // TRI6
174 
175  4, // QUAD4
176  8, // QUAD8
177  9, // QUAD9
178 
179  4, // TET4
180  10, // TET10
181 
182  8, // HEX8
183  20, // HEX20
184  27, // HEX27
185 
186  6, // PRISM6
187  15, // PRISM15
188  18, // PRISM18
189 
190  5, // PYRAMID5
191  13, // PYRAMID13
192  14, // PYRAMID14
193 
194  2, // INFEDGE2
195 
196  4, // INFQUAD4
197  6, // INFQUAD6
198 
199  8, // INFHEX8
200  16, // INFHEX16
201  18, // INFHEX18
202 
203  6, // INFPRISM6
204  12, // INFPRISM12
205 
206  1, // NODEELEM
207 
208  0, // REMOTEELEM
209 
210  3, // TRI3SUBDIVISION
211  3, // TRISHELL3
212  4, // QUADSHELL4
213  8, // QUADSHELL8
214 
215  7, // TRI7
216  14, // TET14
217  20, // PRISM20
218  21, // PRISM21
219  18, // PYRAMID18
220 
221  9, // QUADSHELL9
222 
223  invalid_uint, // C0POLYGON
224  invalid_uint, // C0POLYHEDRON
225  };
226 
227 const unsigned int Elem::type_to_n_sides_map [] =
228  {
229  2, // EDGE2
230  2, // EDGE3
231  2, // EDGE4
232 
233  3, // TRI3
234  3, // TRI6
235 
236  4, // QUAD4
237  4, // QUAD8
238  4, // QUAD9
239 
240  4, // TET4
241  4, // TET10
242 
243  6, // HEX8
244  6, // HEX20
245  6, // HEX27
246 
247  5, // PRISM6
248  5, // PRISM15
249  5, // PRISM18
250 
251  5, // PYRAMID5
252  5, // PYRAMID13
253  5, // PYRAMID14
254 
255  2, // INFEDGE2
256 
257  3, // INFQUAD4
258  3, // INFQUAD6
259 
260  5, // INFHEX8
261  5, // INFHEX16
262  5, // INFHEX18
263 
264  4, // INFPRISM6
265  4, // INFPRISM12
266 
267  0, // NODEELEM
268 
269  0, // REMOTEELEM
270 
271  3, // TRI3SUBDIVISION
272  3, // TRISHELL3
273  4, // QUADSHELL4
274  4, // QUADSHELL8
275 
276  3, // TRI7
277  4, // TET14
278  5, // PRISM20
279  5, // PRISM21
280  5, // PYRAMID18
281 
282  4, // QUADSHELL9
283 
284  invalid_uint, // C0POLYGON
285  invalid_uint, // C0POLYHEDRON
286  };
287 
288 const unsigned int Elem::type_to_n_edges_map [] =
289  {
290  0, // EDGE2
291  0, // EDGE3
292  0, // EDGE4
293 
294  3, // TRI3
295  3, // TRI6
296 
297  4, // QUAD4
298  4, // QUAD8
299  4, // QUAD9
300 
301  6, // TET4
302  6, // TET10
303 
304  12, // HEX8
305  12, // HEX20
306  12, // HEX27
307 
308  9, // PRISM6
309  9, // PRISM15
310  9, // PRISM18
311 
312  8, // PYRAMID5
313  8, // PYRAMID13
314  8, // PYRAMID14
315 
316  0, // INFEDGE2
317 
318  3, // INFQUAD4
319  3, // INFQUAD6
320 
321  8, // INFHEX8
322  8, // INFHEX16
323  8, // INFHEX18
324 
325  6, // INFPRISM6
326  6, // INFPRISM12
327 
328  0, // NODEELEM
329 
330  0, // REMOTEELEM
331 
332  3, // TRI3SUBDIVISION
333  3, // TRISHELL3
334  4, // QUADSHELL4
335  4, // QUADSHELL8
336 
337  3, // TRI7
338  6, // TET14
339  9, // PRISM20
340  9, // PRISM21
341  8, // PYRAMID18
342 
343  4, // QUADSHELL9
344 
345  invalid_uint, // C0POLYGON
346  invalid_uint, // C0POLYHEDRON
347  };
348 
350  {
351  FIRST, // EDGE2
352  SECOND, // EDGE3
353  THIRD, // EDGE4
354 
355  FIRST, // TRI3
356  SECOND, // TRI6
357 
358  FIRST, // QUAD4
359  SECOND, // QUAD8
360  SECOND, // QUAD9
361 
362  FIRST, // TET4
363  SECOND, // TET10
364 
365  FIRST, // HEX8
366  SECOND, // HEX20
367  SECOND, // HEX27
368 
369  FIRST, // PRISM6
370  SECOND, // PRISM15
371  SECOND, // PRISM18
372 
373  FIRST, // PYRAMID5
374  SECOND, // PYRAMID13
375  SECOND, // PYRAMID14
376 
377  FIRST, // INFEDGE2
378 
379  FIRST, // INFQUAD4
380  SECOND, // INFQUAD6
381 
382  FIRST, // INFHEX8
383  SECOND, // INFHEX16
384  SECOND, // INFHEX18
385 
386  FIRST, // INFPRISM6
387  SECOND, // INFPRISM12
388 
389  CONSTANT, // NODEELEM
390 
391  INVALID_ORDER, // REMOTEELEM
392 
393  FIRST, // TRI3SUBDIVISION
394  FIRST, // TRISHELL3
395  FIRST, // QUADSHELL4
396  SECOND, // QUADSHELL8
397 
398  THIRD, // TRI7
399  THIRD, // TET14
400  THIRD, // PRISM20
401  THIRD, // PRISM21
402  THIRD, // PYRAMID18
403 
404  SECOND, // QUADSHELL9
405 
406  FIRST, // C0POLYGON
407  FIRST, // C0POLYHEDRON
408  };
409 
410 // ------------------------------------------------------------
411 // Elem class member functions
412 std::unique_ptr<Elem> Elem::disconnected_clone() const
413 {
414  std::unique_ptr<Elem> returnval;
415 
416  switch (this->type())
417  {
418  case C0POLYGON:
419  returnval = std::make_unique<C0Polygon>(this->n_sides());
420  break;
421 
422  default:
423  returnval = Elem::build(this->type());
424  }
425 
426  returnval->set_id() = this->id();
427 #ifdef LIBMESH_ENABLE_UNIQUE_ID
428  if (this->valid_unique_id())
429  returnval->set_unique_id(this->unique_id());
430 #endif
431 
432  const auto n_elem_ints = this->n_extra_integers();
433  returnval->add_extra_integers(n_elem_ints);
434  for (unsigned int i = 0; i != n_elem_ints; ++i)
435  returnval->set_extra_integer(i, this->get_extra_integer(i));
436 
437  returnval->inherit_data_from(*this);
438 
439  return returnval;
440 }
441 
442 
443 
444 std::unique_ptr<Elem> Elem::build(const ElemType type,
445  Elem * p)
446 {
447  switch (type)
448  {
449  // 0D elements
450  case NODEELEM:
451  return std::make_unique<NodeElem>(p);
452 
453  // 1D elements
454  case EDGE2:
455  return std::make_unique<Edge2>(p);
456  case EDGE3:
457  return std::make_unique<Edge3>(p);
458  case EDGE4:
459  return std::make_unique<Edge4>(p);
460 
461  // 2D elements
462  case TRI3:
463  return std::make_unique<Tri3>(p);
464  case TRISHELL3:
465  return std::make_unique<TriShell3>(p);
466  case TRI3SUBDIVISION:
467  return std::make_unique<Tri3Subdivision>(p);
468  case TRI6:
469  return std::make_unique<Tri6>(p);
470  case TRI7:
471  return std::make_unique<Tri7>(p);
472  case QUAD4:
473  return std::make_unique<Quad4>(p);
474  case QUADSHELL4:
475  return std::make_unique<QuadShell4>(p);
476  case QUAD8:
477  return std::make_unique<Quad8>(p);
478  case QUADSHELL8:
479  return std::make_unique<QuadShell8>(p);
480  case QUAD9:
481  return std::make_unique<Quad9>(p);
482  case QUADSHELL9:
483  return std::make_unique<QuadShell9>(p);
484 
485  // Well, a hexagon is *a* polygon...
486  case C0POLYGON:
487  return std::make_unique<C0Polygon>(6, p);
488 
489  // Building a polyhedron can't currently be done without creating
490  // its nodes first
491  case C0POLYHEDRON:
492  libmesh_not_implemented_msg
493  ("Polyhedra cannot be built via Elem::build()");
494 
495  // 3D elements
496  case TET4:
497  return std::make_unique<Tet4>(p);
498  case TET10:
499  return std::make_unique<Tet10>(p);
500  case TET14:
501  return std::make_unique<Tet14>(p);
502  case HEX8:
503  return std::make_unique<Hex8>(p);
504  case HEX20:
505  return std::make_unique<Hex20>(p);
506  case HEX27:
507  return std::make_unique<Hex27>(p);
508  case PRISM6:
509  return std::make_unique<Prism6>(p);
510  case PRISM15:
511  return std::make_unique<Prism15>(p);
512  case PRISM18:
513  return std::make_unique<Prism18>(p);
514  case PRISM20:
515  return std::make_unique<Prism20>(p);
516  case PRISM21:
517  return std::make_unique<Prism21>(p);
518  case PYRAMID5:
519  return std::make_unique<Pyramid5>(p);
520  case PYRAMID13:
521  return std::make_unique<Pyramid13>(p);
522  case PYRAMID14:
523  return std::make_unique<Pyramid14>(p);
524  case PYRAMID18:
525  return std::make_unique<Pyramid18>(p);
526 
527 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
528  // 1D infinite elements
529  case INFEDGE2:
530  return std::make_unique<InfEdge2>(p);
531 
532  // 2D infinite elements
533  case INFQUAD4:
534  return std::make_unique<InfQuad4>(p);
535  case INFQUAD6:
536  return std::make_unique<InfQuad6>(p);
537 
538  // 3D infinite elements
539  case INFHEX8:
540  return std::make_unique<InfHex8>(p);
541  case INFHEX16:
542  return std::make_unique<InfHex16>(p);
543  case INFHEX18:
544  return std::make_unique<InfHex18>(p);
545  case INFPRISM6:
546  return std::make_unique<InfPrism6>(p);
547  case INFPRISM12:
548  return std::make_unique<InfPrism12>(p);
549 #endif
550 
551  default:
552  libmesh_error_msg("ERROR: Undefined element type == " << Utility::enum_to_string(type));
553  }
554 }
555 
556 
557 
558 std::unique_ptr<Elem> Elem::build_with_id (const ElemType type,
559  dof_id_type id)
560 {
561  // Call the other build() method with nullptr parent, then set the
562  // required id.
563  auto temp = Elem::build(type, nullptr);
564  temp->set_id(id);
565  return temp;
566 }
567 
568 
569 
570 const Elem * Elem::reference_elem () const
571 {
572  return &(ReferenceElem::get(this->type()));
573 }
574 
575 
576 
577 #ifdef LIBMESH_ENABLE_DEPRECATED
579 {
580  libmesh_do_once(libMesh::err
581  << "Elem::centroid() has been deprecated. Replace with either "
582  << "Elem::vertex_average() to maintain existing behavior, or "
583  << "the more expensive Elem::true_centroid() "
584  << "in cases where the true 'geometric' centroid is required."
585  << std::endl);
586  libmesh_deprecated();
587 
588  return Elem::vertex_average();
589 }
590 #endif // LIBMESH_ENABLE_DEPRECATED
591 
592 
593 
595 {
596  // The base class implementation builds a finite element of the correct
597  // order and computes the centroid, c=(cx, cy, cz), where:
598  //
599  // [cx] [\int x dV]
600  // [cy] := (1/V) * [\int y dV]
601  // [cz] [\int z dV]
602  //
603  // using quadrature. Note that we can expand "x" in the FE space as:
604  //
605  // x = \sum_i x_i \phi_i
606  //
607  // where x_i are the nodal positions of the element and \phi_i are the
608  // associated Lagrange shape functions. This allows us to write the
609  // integrals above as e.g.:
610  //
611  // \int x dV = \sum_i x_i \int \phi_i dV
612  //
613  // Defining:
614  //
615  // V_i := \int \phi_i dV
616  //
617  // we then have:
618  //
619  // [cx] [\sum_i x_i V_i]
620  // [cy] = (1/V) * [\sum_i y_i V_i]
621  // [cz] [\sum_i z_i V_i]
622  //
623  // where:
624  // V = \sum_i V_i
625  //
626  // Derived element types can overload this method to compute
627  // the centroid more efficiently when possible.
628 
629  // If this Elem has an elevated p_level, then we need to generate a
630  // barebones copy of it with zero p_level and call true_centroid()
631  // on that instead. This workaround allows us to avoid issues with
632  // calling FE::reinit() with a default_order() FEType, and then
633  // having that order incorrectly boosted by p_level.
634  if (this->p_level())
635  {
636  auto elem_copy = this->disconnected_clone();
637 #ifdef LIBMESH_ENABLE_AMR
638  elem_copy->set_p_level(0);
639 #endif
640 
641  // Set node pointers
642  for (auto n : this->node_index_range())
643  elem_copy->set_node(n, _nodes[n]);
644 
645  return elem_copy->true_centroid();
646  }
647 
648  const FEFamily mapping_family = FEMap::map_fe_type(*this);
649  const FEType fe_type(this->default_order(), mapping_family);
650 
651  // Build FE and attach quadrature rule. The default quadrature rule
652  // integrates the mass matrix exactly, thus it is overkill to
653  // integrate the basis functions, but this is convenient.
654  std::unique_ptr<FEBase> fe = FEBase::build(this->dim(), fe_type);
655  QGauss qrule (this->dim(), fe_type.default_quadrature_order());
656  fe->attach_quadrature_rule(&qrule);
657 
658  // Pre-request required data
659  const auto & JxW = fe->get_JxW();
660  const auto & phi = fe->get_phi();
661 
662  // Re-compute element-specific values
663  fe->reinit(this);
664 
665  // Number of basis functions
666  auto N = phi.size();
667  libmesh_assert_equal_to(N, this->n_nodes());
668 
669  // Compute V_i
670  std::vector<Real> V(N);
671  for (auto qp : index_range(JxW))
672  for (auto i : make_range(N))
673  V[i] += JxW[qp] * phi[i][qp];
674 
675  // Compute centroid
676  Point cp;
677  Real vol = 0.;
678 
679  for (auto i : make_range(N))
680  {
681  cp += this->point(i) * V[i];
682  vol += V[i];
683  }
684 
685  return cp / vol;
686 }
687 
689 {
690  Point cp;
691 
692  const auto n_vertices = this->n_vertices();
693 
694  for (unsigned int n=0; n<n_vertices; n++)
695  cp.add (this->point(n));
696 
697  return (cp /= static_cast<Real>(n_vertices));
698 }
699 
700 
701 
703 {
704  Real h_min=std::numeric_limits<Real>::max();
705 
706  // Avoid calling a virtual a lot of times
707  const auto n_vertices = this->n_vertices();
708 
709  for (unsigned int n_outer=0; n_outer<n_vertices; n_outer++)
710  for (unsigned int n_inner=n_outer+1; n_inner<n_vertices; n_inner++)
711  {
712  const auto diff = (this->point(n_outer) - this->point(n_inner));
713 
714  h_min = std::min(h_min, diff.norm_sq());
715  }
716 
717  return std::sqrt(h_min);
718 }
719 
720 
721 
723 {
724  Real h_max=0;
725 
726  // Avoid calling a virtual a lot of times
727  const auto n_vertices = this->n_vertices();
728 
729  for (unsigned int n_outer=0; n_outer<n_vertices; n_outer++)
730  for (unsigned int n_inner=n_outer+1; n_inner<n_vertices; n_inner++)
731  {
732  const auto diff = (this->point(n_outer) - this->point(n_inner));
733 
734  h_max = std::max(h_max, diff.norm_sq());
735  }
736 
737  return std::sqrt(h_max);
738 }
739 
740 
741 
742 Real Elem::length(const unsigned int n1,
743  const unsigned int n2) const
744 {
745  libmesh_assert_less ( n1, this->n_vertices() );
746  libmesh_assert_less ( n2, this->n_vertices() );
747 
748  return (this->point(n1) - this->point(n2)).norm();
749 }
750 
751 
752 
754 {
755  const unsigned short n_n = this->n_nodes();
756 
757  std::array<dof_id_type, Elem::max_n_nodes> node_ids;
758 
759  for (unsigned short n=0; n != n_n; ++n)
760  node_ids[n] = this->node_id(n);
761 
762  // Always sort, so that different local node numberings hash to the
763  // same value.
764  std::sort (node_ids.begin(), node_ids.begin()+n_n);
765 
766  return Utility::hashword(node_ids.data(), n_n);
767 }
768 
769 
770 
771 bool Elem::operator == (const Elem & rhs) const
772 {
773  // If the elements aren't the same type, they aren't equal
774  if (this->type() != rhs.type())
775  return false;
776 
777  const unsigned short n_n = this->n_nodes();
778  libmesh_assert_equal_to(n_n, rhs.n_nodes());
779 
780  // Make two sorted arrays of global node ids and compare them for
781  // equality.
782  std::array<dof_id_type, Elem::max_n_nodes> this_ids, rhs_ids;
783 
784  for (unsigned short n = 0; n != n_n; n++)
785  {
786  this_ids[n] = this->node_id(n);
787  rhs_ids[n] = rhs.node_id(n);
788  }
789 
790  // Sort the vectors to rule out different local node numberings.
791  std::sort(this_ids.begin(), this_ids.begin()+n_n);
792  std::sort(rhs_ids.begin(), rhs_ids.begin()+n_n);
793 
794  // If the node ids match, the elements are equal!
795  for (unsigned short n = 0; n != n_n; ++n)
796  if (this_ids[n] != rhs_ids[n])
797  return false;
798  return true;
799 }
800 
801 
802 
803 bool Elem::topologically_equal (const Elem & rhs) const
804 {
805  // If the elements aren't the same type, they aren't equal
806  if (this->type() != rhs.type())
807  return false;
808 
809  libmesh_assert_equal_to(this->n_nodes(), rhs.n_nodes());
810 
811  for (auto n : make_range(this->n_nodes()))
812  if (this->node_id(n) != rhs.node_id(n))
813  return false;
814 
815  for (auto neigh : make_range(this->n_neighbors()))
816  {
817  if (!this->neighbor_ptr(neigh))
818  {
819  if (rhs.neighbor_ptr(neigh))
820  return false;
821  continue;
822  }
823  if (!rhs.neighbor_ptr(neigh) ||
824  this->neighbor_ptr(neigh)->id() !=
825  rhs.neighbor_ptr(neigh)->id())
826  return false;
827  }
828 
829  if (this->parent())
830  {
831  if (!rhs.parent())
832  return false;
833  if (this->parent()->id() != rhs.parent()->id())
834  return false;
835  }
836  else if (rhs.parent())
837  return false;
838 
839  if (this->interior_parent())
840  {
841  if (!rhs.interior_parent())
842  return false;
843  if (this->interior_parent()->id() !=
844  rhs.interior_parent()->id())
845  return false;
846  }
847  else if (rhs.interior_parent())
848  return false;
849 
850  return true;
851 }
852 
853 
854 
855 bool Elem::is_semilocal(const processor_id_type my_pid) const
856 {
857  std::set<const Elem *> point_neighbors;
858 
859  this->find_point_neighbors(point_neighbors);
860 
861  for (const auto & elem : point_neighbors)
862  if (elem->processor_id() == my_pid)
863  return true;
864 
865  return false;
866 }
867 
868 
869 
870 unsigned int Elem::which_side_am_i (const Elem * e) const
871 {
872  libmesh_assert(e);
873 
874  const unsigned int ns = this->n_sides();
875  const unsigned int nn = this->n_nodes();
876 
877  const unsigned int en = e->n_nodes();
878 
879  // e might be on any side until proven otherwise
880  std::vector<bool> might_be_side(ns, true);
881 
882  for (unsigned int i=0; i != en; ++i)
883  {
884  Point side_point = e->point(i);
885  unsigned int local_node_id = libMesh::invalid_uint;
886 
887  // Look for a node of this that's contiguous with node i of
888  // e. Note that the exact floating point comparison of Point
889  // positions is intentional, see the class documentation for
890  // this function.
891  for (unsigned int j=0; j != nn; ++j)
892  if (this->point(j) == side_point)
893  local_node_id = j;
894 
895  // If a node of e isn't contiguous with some node of this, then
896  // e isn't a side of this.
897  if (local_node_id == libMesh::invalid_uint)
898  return libMesh::invalid_uint;
899 
900  // If a node of e isn't contiguous with some node on side s of
901  // this, then e isn't on side s.
902  for (unsigned int s=0; s != ns; ++s)
903  if (!this->is_node_on_side(local_node_id, s))
904  might_be_side[s] = false;
905  }
906 
907  for (unsigned int s=0; s != ns; ++s)
908  if (might_be_side[s])
909  {
910 #ifdef DEBUG
911  for (unsigned int s2=s+1; s2 < ns; ++s2)
912  libmesh_assert (!might_be_side[s2]);
913 #endif
914  return s;
915  }
916 
917  // Didn't find any matching side
918  return libMesh::invalid_uint;
919 }
920 
921 
922 
923 #ifdef LIBMESH_ENABLE_DEPRECATED
924 unsigned int Elem::which_node_am_i(unsigned int side,
925  unsigned int side_node) const
926 {
927  libmesh_deprecated();
928  return local_side_node(side, side_node);
929 }
930 #endif // LIBMESH_ENABLE_DEPRECATED
931 
932 
933 
934 bool Elem::contains_vertex_of(const Elem * e, bool mesh_connection) const
935 {
936  // Our vertices are the first numbered nodes
937  const unsigned int nv = e->n_vertices();
938  const unsigned int my_nv = this->n_vertices();
939 
940  // Check for vertex-to-vertex containment first; contains_point() is
941  // expensive
942  for (auto n : make_range(nv))
943  {
944  const Node * vertex = e->node_ptr(n);
945  for (auto my_n : make_range(my_nv))
946  if (&this->node_ref(my_n) == vertex)
947  return true;
948  }
949 
950  // If e is in our mesh, then we might be done testing
951  if (mesh_connection)
952  {
953  const unsigned int l = this->level();
954  const unsigned int el = e->level();
955 
956  if (l >= el)
957  return false;
958 
959  // We could also return false for l==el-1 iff we knew we had no
960  // triangular faces, but we don't have an API to check that.
961  }
962 
963  // Our vertices are the first numbered nodes
964  for (auto n : make_range(nv))
965  if (this->contains_point(e->point(n)))
966  return true;
967  return false;
968 }
969 
970 
971 
972 bool Elem::contains_edge_of(const Elem * e) const
973 {
974  unsigned int num_contained_edges = 0;
975 
976  // Our vertices are the first numbered nodes
977  for (auto n : make_range(e->n_vertices()))
978  {
979  if (this->contains_point(e->point(n)))
980  {
981  num_contained_edges++;
982  if (num_contained_edges>=2)
983  {
984  return true;
985  }
986  }
987  }
988  return false;
989 }
990 
991 
992 
994  std::set<const Elem *> & neighbor_set) const
995 {
996  libmesh_assert(this->contains_point(p));
997  libmesh_assert(this->active());
998 
999  neighbor_set.clear();
1000  neighbor_set.insert(this);
1001 
1002  std::set<const Elem *> untested_set, next_untested_set;
1003  untested_set.insert(this);
1004 
1005 #ifdef LIBMESH_ENABLE_AMR
1006  std::vector<const Elem *> active_neighbor_children;
1007 #endif // #ifdef LIBMESH_ENABLE_AMR
1008 
1009  while (!untested_set.empty())
1010  {
1011  // Loop over all the elements in the patch that haven't already
1012  // been tested
1013  for (const auto & elem : untested_set)
1014  for (auto current_neighbor : elem->neighbor_ptr_range())
1015  {
1016  if (current_neighbor &&
1017  current_neighbor != remote_elem) // we have a real neighbor on this side
1018  {
1019  if (current_neighbor->active()) // ... if it is active
1020  {
1021  auto it = neighbor_set.lower_bound(current_neighbor);
1022  if ((it == neighbor_set.end() || *it != current_neighbor) &&
1023  current_neighbor->contains_point(p)) // ... don't have and touches p
1024  {
1025  // Add it and test it
1026  next_untested_set.insert(current_neighbor);
1027  neighbor_set.emplace_hint(it, current_neighbor);
1028  }
1029  }
1030 #ifdef LIBMESH_ENABLE_AMR
1031  else // ... the neighbor is *not* active,
1032  { // ... so add *all* neighboring
1033  // active children that touch p
1034  active_neighbor_children.clear();
1035  current_neighbor->active_family_tree_by_neighbor
1036  (active_neighbor_children, elem);
1037 
1038  for (const auto & current_child : active_neighbor_children)
1039  {
1040  auto it = neighbor_set.lower_bound(current_child);
1041  if ((it == neighbor_set.end() || *it != current_child) &&
1042  current_child->contains_point(p))
1043  {
1044  // Add it and test it
1045  next_untested_set.insert(current_child);
1046  neighbor_set.emplace_hint(it, current_child);
1047  }
1048  }
1049  }
1050 #endif // #ifdef LIBMESH_ENABLE_AMR
1051  }
1052  }
1053  untested_set.swap(next_untested_set);
1054  next_untested_set.clear();
1055  }
1056 }
1057 
1058 
1059 
1060 void Elem::find_point_neighbors(std::set<const Elem *> & neighbor_set) const
1061 {
1062  this->find_point_neighbors(neighbor_set, this);
1063 }
1064 
1065 
1066 
1067 void Elem::find_point_neighbors(std::set<const Elem *> & neighbor_set,
1068  const Elem * start_elem) const
1069 {
1070  ElemInternal::find_point_neighbors(this, neighbor_set, start_elem);
1071 }
1072 
1073 
1074 
1075 void Elem::find_point_neighbors(std::set<Elem *> & neighbor_set,
1076  Elem * start_elem)
1077 {
1078  ElemInternal::find_point_neighbors(this, neighbor_set, start_elem);
1079 }
1080 
1081 
1082 
1084  const Point & p2,
1085  std::set<const Elem *> & neighbor_set) const
1086 {
1087  // Simple but perhaps suboptimal code: find elements containing the
1088  // first point, then winnow this set down by removing elements which
1089  // don't also contain the second point
1090 
1091  libmesh_assert(this->contains_point(p2));
1092  this->find_point_neighbors(p1, neighbor_set);
1093 
1094  std::set<const Elem *>::iterator it = neighbor_set.begin();
1095  const std::set<const Elem *>::iterator end = neighbor_set.end();
1096 
1097  while (it != end)
1098  {
1099  // As of C++11, set::erase returns an iterator to the element
1100  // following the erased element, or end.
1101  if (!(*it)->contains_point(p2))
1102  it = neighbor_set.erase(it);
1103  else
1104  ++it;
1105  }
1106 }
1107 
1108 
1109 
1110 void Elem::find_edge_neighbors(std::set<const Elem *> & neighbor_set) const
1111 {
1112  neighbor_set.clear();
1113  neighbor_set.insert(this);
1114 
1115  std::set<const Elem *> untested_set, next_untested_set;
1116  untested_set.insert(this);
1117 
1118  while (!untested_set.empty())
1119  {
1120  // Loop over all the elements in the patch that haven't already
1121  // been tested
1122  for (const auto & elem : untested_set)
1123  {
1124  for (auto current_neighbor : elem->neighbor_ptr_range())
1125  {
1126  if (current_neighbor &&
1127  current_neighbor != remote_elem) // we have a real neighbor on this side
1128  {
1129  if (current_neighbor->active()) // ... if it is active
1130  {
1131  if (this->contains_edge_of(current_neighbor) // ... and touches us
1132  || current_neighbor->contains_edge_of(this))
1133  {
1134  // Make sure we'll test it
1135  if (!neighbor_set.count(current_neighbor))
1136  next_untested_set.insert (current_neighbor);
1137 
1138  // And add it
1139  neighbor_set.insert (current_neighbor);
1140  }
1141  }
1142 #ifdef LIBMESH_ENABLE_AMR
1143  else // ... the neighbor is *not* active,
1144  { // ... so add *all* neighboring
1145  // active children
1146  std::vector<const Elem *> active_neighbor_children;
1147 
1148  current_neighbor->active_family_tree_by_neighbor
1149  (active_neighbor_children, elem);
1150 
1151  for (const auto & current_child : active_neighbor_children)
1152  if (this->contains_edge_of(current_child) || current_child->contains_edge_of(this))
1153  {
1154  // Make sure we'll test it
1155  if (!neighbor_set.count(current_child))
1156  next_untested_set.insert (current_child);
1157 
1158  neighbor_set.insert (current_child);
1159  }
1160  }
1161 #endif // #ifdef LIBMESH_ENABLE_AMR
1162  }
1163  }
1164  }
1165  untested_set.swap(next_untested_set);
1166  next_untested_set.clear();
1167  }
1168 }
1169 
1170 
1171 
1172 void Elem::find_interior_neighbors(std::set<const Elem *> & neighbor_set) const
1173 {
1174  ElemInternal::find_interior_neighbors(this, neighbor_set);
1175 }
1176 
1177 
1178 
1179 void Elem::find_interior_neighbors(std::set<Elem *> & neighbor_set)
1180 {
1181  ElemInternal::find_interior_neighbors(this, neighbor_set);
1182 }
1183 
1184 
1185 
1186 const Elem * Elem::interior_parent () const
1187 {
1188  // interior parents make no sense for full-dimensional elements.
1189  if (this->dim() >= LIBMESH_DIM)
1190  return nullptr;
1191 
1192  // they USED TO BE only good for level-0 elements, but we now
1193  // support keeping interior_parent() valid on refined boundary
1194  // elements.
1195  // if (this->level() != 0)
1196  // return this->parent()->interior_parent();
1197 
1198  // We store the interior_parent pointer after both the parent
1199  // neighbor and neighbor pointers
1200  Elem * interior_p = _elemlinks[1+this->n_sides()];
1201 
1202  // If we have an interior_parent, we USED TO assume it was a
1203  // one-higher-dimensional interior element, but we now allow e.g.
1204  // edge elements to have a 3D interior_parent with no
1205  // intermediate 2D element.
1206  // libmesh_assert (!interior_p ||
1207  // interior_p->dim() == (this->dim()+1));
1208  libmesh_assert (!interior_p ||
1209  (interior_p == remote_elem) ||
1210  (interior_p->dim() > this->dim()));
1211 
1212  // If an element in a multi-dimensional mesh has an interior_parent
1213  // link, it should be at our level or coarser, just like a neighbor
1214  // link. Our collect_families() code relies on this, but it might
1215  // be tempting for users to manually assign something that breaks
1216  // it.
1217  //
1218  // However, we *also* create temporary side elements, and we don't
1219  // bother with creating ancestors for those, so they can be at level
1220  // 0 even when they're sides of non-level-0 elements.
1221  libmesh_assert (!interior_p ||
1222  (interior_p->level() <= this->level()) ||
1223  (this->level() == 0 &&
1224  this->id() == DofObject::invalid_id));
1225 
1226  return interior_p;
1227 }
1228 
1229 
1230 
1232 {
1233  // See the const version for comments
1234  if (this->dim() >= LIBMESH_DIM)
1235  return nullptr;
1236 
1237  Elem * interior_p = _elemlinks[1+this->n_sides()];
1238 
1239  libmesh_assert (!interior_p ||
1240  (interior_p == remote_elem) ||
1241  (interior_p->dim() > this->dim()));
1242 
1243  return interior_p;
1244 }
1245 
1246 
1247 
1249 {
1250  // interior parents make no sense for full-dimensional elements.
1251  libmesh_assert (!p ||
1252  this->dim() < LIBMESH_DIM);
1253 
1254  // If we have an interior_parent, we USED TO assume it was a
1255  // one-higher-dimensional interior element, but we now allow e.g.
1256  // edge elements to have a 3D interior_parent with no
1257  // intermediate 2D element.
1258  // libmesh_assert (!p ||
1259  // p->dim() == (this->dim()+1));
1260  libmesh_assert (!p ||
1261  (p == remote_elem) ||
1262  (p->dim() > this->dim()));
1263 
1264  _elemlinks[1+this->n_sides()] = p;
1265 }
1266 
1267 
1268 
1269 #ifdef LIBMESH_ENABLE_PERIODIC
1270 
1271 Elem * Elem::topological_neighbor (const unsigned int i,
1272  MeshBase & mesh,
1273  const PointLocatorBase & point_locator,
1274  const PeriodicBoundaries * pb)
1275 {
1276  libmesh_assert_less (i, this->n_neighbors());
1277 
1278  Elem * neighbor_i = this->neighbor_ptr(i);
1279  if (neighbor_i != nullptr)
1280  return neighbor_i;
1281 
1282  if (pb)
1283  {
1284  // Since the neighbor is nullptr it must be on a boundary. We need
1285  // see if this is a periodic boundary in which case it will have a
1286  // topological neighbor
1287  std::vector<boundary_id_type> bc_ids;
1288  mesh.get_boundary_info().boundary_ids(this, cast_int<unsigned short>(i), bc_ids);
1289  for (const auto & id : bc_ids)
1290  if (pb->boundary(id))
1291  {
1292  // Since the point locator inside of periodic boundaries
1293  // returns a const pointer we will retrieve the proper
1294  // pointer directly from the mesh object.
1295  const Elem * const cn = pb->neighbor(id, point_locator, this, i);
1296  neighbor_i = const_cast<Elem *>(cn);
1297 
1298  // Since coarse elements do not have more refined
1299  // neighbors we need to make sure that we don't return one
1300  // of these types of neighbors.
1301  if (neighbor_i)
1302  while (level() < neighbor_i->level())
1303  neighbor_i = neighbor_i->parent();
1304  return neighbor_i;
1305  }
1306  }
1307 
1308  return nullptr;
1309 }
1310 
1311 
1312 
1313 const Elem * Elem::topological_neighbor (const unsigned int i,
1314  const MeshBase & mesh,
1315  const PointLocatorBase & point_locator,
1316  const PeriodicBoundaries * pb) const
1317 {
1318  libmesh_assert_less (i, this->n_neighbors());
1319 
1320  const Elem * neighbor_i = this->neighbor_ptr(i);
1321  if (neighbor_i != nullptr)
1322  return neighbor_i;
1323 
1324  if (pb)
1325  {
1326  // Since the neighbor is nullptr it must be on a boundary. We need
1327  // see if this is a periodic boundary in which case it will have a
1328  // topological neighbor
1329  std::vector<boundary_id_type> bc_ids;
1330  mesh.get_boundary_info().boundary_ids(this, cast_int<unsigned short>(i), bc_ids);
1331  for (const auto & id : bc_ids)
1332  if (pb->boundary(id))
1333  {
1334  neighbor_i = pb->neighbor(id, point_locator, this, i);
1335 
1336  // Since coarse elements do not have more refined
1337  // neighbors we need to make sure that we don't return one
1338  // of these types of neighbors.
1339  if (neighbor_i)
1340  while (level() < neighbor_i->level())
1341  neighbor_i = neighbor_i->parent();
1342  return neighbor_i;
1343  }
1344  }
1345 
1346  return nullptr;
1347 }
1348 
1349 
1351  const MeshBase & mesh,
1352  const PointLocatorBase & point_locator,
1353  const PeriodicBoundaries * pb) const
1354 {
1355  // First see if this is a normal "interior" neighbor
1356  if (has_neighbor(elem))
1357  return true;
1358 
1359  for (auto n : this->side_index_range())
1360  if (this->topological_neighbor(n, mesh, point_locator, pb))
1361  return true;
1362 
1363  return false;
1364 }
1365 
1366 
1367 #endif
1368 
1369 #ifndef NDEBUG
1370 
1372 {
1373  libmesh_assert(this->valid_id());
1374  for (auto n : this->node_index_range())
1375  {
1376  libmesh_assert(this->node_ptr(n));
1377  libmesh_assert(this->node_ptr(n)->valid_id());
1378  }
1379 }
1380 
1381 
1382 
1384 {
1385  for (auto n : this->side_index_range())
1386  {
1387  const Elem * neigh = this->neighbor_ptr(n);
1388 
1389  // Any element might have a remote neighbor; checking
1390  // to make sure that's not inaccurate is tough.
1391  if (neigh == remote_elem)
1392  continue;
1393 
1394  if (neigh)
1395  {
1396  // Only subactive elements have subactive neighbors
1397  libmesh_assert (this->subactive() || !neigh->subactive());
1398 
1399  const Elem * elem = this;
1400 
1401  // If we're subactive but our neighbor isn't, its
1402  // return neighbor link will be to our first active
1403  // ancestor OR to our inactive ancestor of the same
1404  // level as neigh,
1405  if (this->subactive() && !neigh->subactive())
1406  {
1407  for (elem = this; !elem->active();
1408  elem = elem->parent())
1409  libmesh_assert(elem);
1410  }
1411  else
1412  {
1413  unsigned int rev = neigh->which_neighbor_am_i(elem);
1414  libmesh_assert_less (rev, neigh->n_neighbors());
1415 
1416  if (this->subactive() && !neigh->subactive())
1417  {
1418  while (neigh->neighbor_ptr(rev) != elem)
1419  {
1420  libmesh_assert(elem->parent());
1421  elem = elem->parent();
1422  }
1423  }
1424  else
1425  {
1426  const Elem * nn = neigh->neighbor_ptr(rev);
1427  libmesh_assert(nn);
1428 
1429  for (; elem != nn; elem = elem->parent())
1430  libmesh_assert(elem);
1431  }
1432  }
1433  }
1434  // If we don't have a neighbor and we're not subactive, our
1435  // ancestors shouldn't have any neighbors in this same
1436  // direction.
1437  else if (!this->subactive())
1438  {
1439  const Elem * my_parent = this->parent();
1440  if (my_parent &&
1441  // A parent with a different dimension isn't really one of
1442  // our ancestors, it means we're on a boundary mesh and this
1443  // is an interior mesh element for which we're on a side.
1444  // Nothing to test for in that case.
1445  (my_parent->dim() == this->dim()))
1446  libmesh_assert (!my_parent->neighbor_ptr(n));
1447  }
1448  }
1449 }
1450 
1451 #endif // !NDEBUG
1452 
1453 
1454 
1455 void Elem::make_links_to_me_local(unsigned int n, unsigned int nn)
1456 {
1457  Elem * neigh = this->neighbor_ptr(n);
1458 
1459  // Don't bother calling this function unless it's necessary
1460  libmesh_assert(neigh);
1461  libmesh_assert(!neigh->is_remote());
1462 
1463  // We never have neighbors more refined than us
1464  libmesh_assert_less_equal (neigh->level(), this->level());
1465 
1466  // We never have subactive neighbors of non subactive elements
1467  libmesh_assert(!neigh->subactive() || this->subactive());
1468 
1469  // If we have a neighbor less refined than us then it must not
1470  // have any more refined descendants we could have pointed to
1471  // instead.
1472  libmesh_assert((neigh->level() == this->level()) ||
1473  (neigh->active() && !this->subactive()) ||
1474  (!neigh->has_children() && this->subactive()));
1475 
1476  // If neigh is at our level, then its family might have
1477  // remote_elem neighbor links which need to point to us
1478  // instead, but if not, then we're done.
1479  if (neigh->level() != this->level())
1480  return;
1481 
1482  // What side of neigh are we on? nn.
1483  //
1484  // We can't use the usual Elem method because we're in the middle of
1485  // restoring topology. We can't compare side_ptr nodes because
1486  // users want to abuse neighbor_ptr to point to
1487  // not-technically-neighbors across mesh slits. We can't compare
1488  // node locations because users want to move those
1489  // not-technically-neighbors until they're
1490  // not-even-geometrically-neighbors.
1491 
1492  // Find any elements that ought to point to elem
1493  std::vector<Elem *> neigh_family;
1494 #ifdef LIBMESH_ENABLE_AMR
1495  if (this->active())
1496  neigh->family_tree_by_side(neigh_family, nn);
1497  else
1498 #endif
1499  neigh_family.push_back(neigh);
1500 
1501  // And point them to elem
1502  for (auto & neigh_family_member : neigh_family)
1503  {
1504  // Only subactive elements point to other subactive elements
1505  if (this->subactive() && !neigh_family_member->subactive())
1506  continue;
1507 
1508  // Ideally, the neighbor link ought to either be correct
1509  // already or ought to be to remote_elem.
1510  //
1511  // However, if we're redistributing a newly created elem,
1512  // after an AMR step but before find_neighbors has fixed up
1513  // neighbor links, we might have an out of date neighbor
1514  // link to elem's parent instead.
1515 #ifdef LIBMESH_ENABLE_AMR
1516  libmesh_assert((neigh_family_member->neighbor_ptr(nn) &&
1517  (neigh_family_member->neighbor_ptr(nn)->active() ||
1518  neigh_family_member->neighbor_ptr(nn)->is_ancestor_of(this))) ||
1519  (neigh_family_member->neighbor_ptr(nn) == remote_elem) ||
1520  ((this->refinement_flag() == JUST_REFINED) &&
1521  (this->parent() != nullptr) &&
1522  (neigh_family_member->neighbor_ptr(nn) == this->parent())));
1523 #else
1524  libmesh_assert((neigh_family_member->neighbor_ptr(nn) == this) ||
1525  (neigh_family_member->neighbor_ptr(nn) == remote_elem));
1526 #endif
1527 
1528  neigh_family_member->set_neighbor(nn, this);
1529  }
1530 }
1531 
1532 
1534 {
1535  libmesh_assert_not_equal_to (this, remote_elem);
1536 
1537  // We need to have handled any children first
1538 #if defined(LIBMESH_ENABLE_AMR) && defined(DEBUG)
1539  if (this->has_children())
1540  for (auto & child : this->child_ref_range())
1541  libmesh_assert_equal_to (&child, remote_elem);
1542 #endif
1543 
1544  // Remotify any neighbor links
1545  for (auto neigh : this->neighbor_ptr_range())
1546  {
1547  if (neigh && neigh != remote_elem)
1548  {
1549  // My neighbor should never be more refined than me; my real
1550  // neighbor would have been its parent in that case.
1551  libmesh_assert_greater_equal (this->level(), neigh->level());
1552 
1553  if (this->level() == neigh->level() &&
1554  neigh->has_neighbor(this))
1555  {
1556 #ifdef LIBMESH_ENABLE_AMR
1557  // My neighbor may have descendants which also consider me a
1558  // neighbor
1559  std::vector<Elem *> family;
1560  neigh->total_family_tree_by_neighbor (family, this);
1561 
1562  // FIXME - There's a lot of ugly const_casts here; we
1563  // may want to make remote_elem non-const
1564  for (auto & n : family)
1565  {
1566  libmesh_assert (n);
1567  if (n == remote_elem)
1568  continue;
1569  unsigned int my_s = n->which_neighbor_am_i(this);
1570  libmesh_assert_less (my_s, n->n_neighbors());
1571  libmesh_assert_equal_to (n->neighbor_ptr(my_s), this);
1572  n->set_neighbor(my_s, const_cast<RemoteElem *>(remote_elem));
1573  }
1574 #else
1575  unsigned int my_s = neigh->which_neighbor_am_i(this);
1576  libmesh_assert_less (my_s, neigh->n_neighbors());
1577  libmesh_assert_equal_to (neigh->neighbor_ptr(my_s), this);
1578  neigh->set_neighbor(my_s, const_cast<RemoteElem *>(remote_elem));
1579 #endif
1580  }
1581 #ifdef LIBMESH_ENABLE_AMR
1582  // Even if my neighbor doesn't link back to me, it might
1583  // have subactive descendants which do
1584  else if (neigh->has_children())
1585  {
1586  // If my neighbor at the same level doesn't have me as a
1587  // neighbor, I must be subactive
1588  libmesh_assert(this->level() > neigh->level() ||
1589  this->subactive());
1590 
1591  // My neighbor must have some ancestor of mine as a
1592  // neighbor
1593  Elem * my_ancestor = this->parent();
1594  libmesh_assert(my_ancestor);
1595  while (!neigh->has_neighbor(my_ancestor))
1596  {
1597  my_ancestor = my_ancestor->parent();
1598  libmesh_assert(my_ancestor);
1599  }
1600 
1601  // My neighbor may have descendants which consider me a
1602  // neighbor
1603  std::vector<Elem *> family;
1604  neigh->total_family_tree_by_subneighbor (family, my_ancestor, this);
1605 
1606  for (auto & n : family)
1607  {
1608  libmesh_assert (n);
1609  if (n->is_remote())
1610  continue;
1611  unsigned int my_s = n->which_neighbor_am_i(this);
1612  libmesh_assert_less (my_s, n->n_neighbors());
1613  libmesh_assert_equal_to (n->neighbor_ptr(my_s), this);
1614  // TODO: we may want to make remote_elem non-const.
1615  n->set_neighbor(my_s, const_cast<RemoteElem *>(remote_elem));
1616  }
1617  }
1618 #endif
1619  }
1620  }
1621 
1622 #ifdef LIBMESH_ENABLE_AMR
1623  // Remotify parent's child link
1624  Elem * my_parent = this->parent();
1625  if (my_parent &&
1626  // As long as it's not already remote
1627  my_parent != remote_elem &&
1628  // And it's a real parent, not an interior parent
1629  this->dim() == my_parent->dim())
1630  {
1631  unsigned int me = my_parent->which_child_am_i(this);
1632  libmesh_assert_equal_to (my_parent->child_ptr(me), this);
1633  my_parent->set_child(me, const_cast<RemoteElem *>(remote_elem));
1634  }
1635 #endif
1636 }
1637 
1638 
1640 {
1641  libmesh_assert_not_equal_to (this, remote_elem);
1642 
1643  // We need to have handled any children first
1644 #ifdef LIBMESH_ENABLE_AMR
1645  libmesh_assert (!this->has_children());
1646 #endif
1647 
1648  // Nullify any neighbor links
1649  for (auto neigh : this->neighbor_ptr_range())
1650  {
1651  if (neigh && neigh != remote_elem)
1652  {
1653  // My neighbor should never be more refined than me; my real
1654  // neighbor would have been its parent in that case.
1655  libmesh_assert_greater_equal (this->level(), neigh->level());
1656 
1657  if (this->level() == neigh->level() &&
1658  neigh->has_neighbor(this))
1659  {
1660 #ifdef LIBMESH_ENABLE_AMR
1661  // My neighbor may have descendants which also consider me a
1662  // neighbor
1663  std::vector<Elem *> family;
1664  neigh->total_family_tree_by_neighbor (family, this);
1665 
1666  for (auto & n : family)
1667  {
1668  libmesh_assert (n);
1669  if (n->is_remote())
1670  continue;
1671  unsigned int my_s = n->which_neighbor_am_i(this);
1672  libmesh_assert_less (my_s, n->n_neighbors());
1673  libmesh_assert_equal_to (n->neighbor_ptr(my_s), this);
1674  n->set_neighbor(my_s, nullptr);
1675  }
1676 #else
1677  unsigned int my_s = neigh->which_neighbor_am_i(this);
1678  libmesh_assert_less (my_s, neigh->n_neighbors());
1679  libmesh_assert_equal_to (neigh->neighbor_ptr(my_s), this);
1680  neigh->set_neighbor(my_s, nullptr);
1681 #endif
1682  }
1683 #ifdef LIBMESH_ENABLE_AMR
1684  // Even if my neighbor doesn't link back to me, it might
1685  // have subactive descendants which do
1686  else if (neigh->has_children())
1687  {
1688  // If my neighbor at the same level doesn't have me as a
1689  // neighbor, I must be subactive
1690  libmesh_assert(this->level() > neigh->level() ||
1691  this->subactive());
1692 
1693  // My neighbor must have some ancestor of mine as a
1694  // neighbor
1695  Elem * my_ancestor = this->parent();
1696  libmesh_assert(my_ancestor);
1697  while (!neigh->has_neighbor(my_ancestor))
1698  {
1699  my_ancestor = my_ancestor->parent();
1700  libmesh_assert(my_ancestor);
1701  }
1702 
1703  // My neighbor may have descendants which consider me a
1704  // neighbor
1705  std::vector<Elem *> family;
1706  neigh->total_family_tree_by_subneighbor (family, my_ancestor, this);
1707 
1708  for (auto & n : family)
1709  {
1710  libmesh_assert (n);
1711  if (n->is_remote())
1712  continue;
1713  unsigned int my_s = n->which_neighbor_am_i(this);
1714  libmesh_assert_less (my_s, n->n_neighbors());
1715  libmesh_assert_equal_to (n->neighbor_ptr(my_s), this);
1716  n->set_neighbor(my_s, nullptr);
1717  }
1718  }
1719 #endif
1720  }
1721  }
1722 
1723 #ifdef LIBMESH_ENABLE_AMR
1724  // We can't currently delete a child with a parent!
1725  libmesh_assert (!this->parent());
1726 #endif
1727 }
1728 
1729 
1730 
1731 void Elem::write_connectivity (std::ostream & out_stream,
1732  const IOPackage iop) const
1733 {
1734  libmesh_assert (out_stream.good());
1736  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
1737 
1738  switch (iop)
1739  {
1740  case TECPLOT:
1741  {
1742  // This connectivity vector will be used repeatedly instead
1743  // of being reconstructed inside the loop.
1744  std::vector<dof_id_type> conn;
1745  for (auto sc : make_range(this->n_sub_elem()))
1746  {
1747  this->connectivity(sc, TECPLOT, conn);
1748 
1749  std::copy(conn.begin(),
1750  conn.end(),
1751  std::ostream_iterator<dof_id_type>(out_stream, " "));
1752 
1753  out_stream << '\n';
1754  }
1755  return;
1756  }
1757 
1758  case UCD:
1759  {
1760  for (auto i : this->node_index_range())
1761  out_stream << this->node_id(i)+1 << "\t";
1762 
1763  out_stream << '\n';
1764  return;
1765  }
1766 
1767  default:
1768  libmesh_error_msg("Unsupported IO package " << iop);
1769  }
1770 }
1771 
1772 
1773 
1775 {
1776  switch (q)
1777  {
1778  // Return the maximum ratio of edge lengths, or zero if that
1779  // maximum would otherwise be infinity.
1780  case EDGE_LENGTH_RATIO:
1781  {
1782  if (this->dim() < 2)
1783  return 1;
1784 
1785  std::vector<Real> edge_lengths(this->n_edges());
1786  for (auto e : index_range(edge_lengths))
1787  edge_lengths[e] = (this->point(this->local_edge_node(e,1)) -
1788  this->point(this->local_edge_node(e,0))).norm();
1789 
1790  auto [min, max] =
1791  std::minmax_element(edge_lengths.begin(), edge_lengths.end());
1792 
1793  if (*min == 0.)
1794  return 0.;
1795  else
1796  return *max / *min;
1797  }
1798 
1799  case MIN_ANGLE:
1800  case MAX_ANGLE:
1801  {
1802  // 1D elements don't have interior angles, so just return some
1803  // dummy value in that case.
1804  if (this->dim() < 2)
1805  return 0.;
1806 
1807  // Initialize return values
1808  Real min_angle = std::numeric_limits<Real>::max();
1809  Real max_angle = -std::numeric_limits<Real>::max();
1810 
1811  for (auto n : this->node_index_range())
1812  {
1813  // Get list of edge ids adjacent to this node.
1814  const auto adjacent_edge_ids = this->edges_adjacent_to_node(n);
1815 
1816  // Skip any nodes with fewer than 2 adjacent edges. You
1817  // need at least two adjacent edges to form an interior
1818  // element angle.
1819  auto N = adjacent_edge_ids.size();
1820  if (N < 2)
1821  continue;
1822 
1823  // Consider all possible pairs of edges adjacent to node n
1824  for (unsigned int first = 0; first < N-1; ++first)
1825  for (unsigned int second = first+1; second < N; ++second)
1826  {
1827  // Get ids of first and second edges
1828  auto first_edge = adjacent_edge_ids[first];
1829  auto second_edge = adjacent_edge_ids[second];
1830 
1831  // Get node ids of first and second edge
1832  auto first_edge_node_0 = this->local_edge_node(first_edge, 0);
1833  auto first_edge_node_1 = this->local_edge_node(first_edge, 1);
1834  auto second_edge_node_0 = this->local_edge_node(second_edge, 0);
1835  auto second_edge_node_1 = this->local_edge_node(second_edge, 1);
1836 
1837  // Orient both edges so that edge node 0 == n
1838  if (first_edge_node_0 != n)
1839  std::swap(first_edge_node_0, first_edge_node_1);
1840  if (second_edge_node_0 != n)
1841  std::swap(second_edge_node_0, second_edge_node_1);
1842 
1843  libmesh_assert_equal_to(first_edge_node_0, n);
1844  libmesh_assert_equal_to(second_edge_node_0, n);
1845 
1846  // Locally oriented edge vectors
1847  Point
1848  first_ev = this->point(first_edge_node_1) - this->point(first_edge_node_0),
1849  second_ev = this->point(second_edge_node_1) - this->point(second_edge_node_0);
1850 
1851  // Angle between them in the range [0, pi]
1852  Real theta = std::acos(first_ev.unit() * second_ev.unit());
1853 
1854  // Track min and max angles seen
1855  min_angle = std::min(theta, min_angle);
1856  max_angle = std::max(theta, max_angle);
1857  }
1858  }
1859 
1860  // Return requested extreme value (in degrees)
1861  return Real(180)/libMesh::pi * ((q == MIN_ANGLE) ? min_angle : max_angle);
1862  }
1863 
1864  case MIN_DIHEDRAL_ANGLE:
1865  case MAX_DIHEDRAL_ANGLE:
1866  {
1867  // For 1D and 2D elements, just call this function with MIN,MAX_ANGLE instead.
1868  if (this->dim() < 3)
1869  return this->quality((q == MIN_DIHEDRAL_ANGLE) ? MIN_ANGLE : MAX_ANGLE);
1870 
1871  // Initialize return values
1872  Real min_angle = std::numeric_limits<Real>::max();
1873  Real max_angle = -std::numeric_limits<Real>::max();
1874 
1875  // Algorithm for computing dihedral angles:
1876  // .) Loop over the edges, use the edge_sides_map to get the
1877  // two sides adjacent to the edge.
1878  // .) For each adjacent side, use the side_nodes_map to get the
1879  // list of three (or more) node ids on that side.
1880  // .) Use the node ids to compute the (approximate) inward
1881  // normal for the side. Note: this is approximate since 3D
1882  // elements with four-sided faces and quadratic tetrahedra
1883  // do not necessarily have planar faces.
1884  // .) Compute the dihedral angle for the current edge, compare
1885  // it to the current min and max dihedral angles.
1886  for (auto e : this->edge_index_range())
1887  {
1888  // Get list of edge ids adjacent to this node.
1889  const auto adjacent_side_ids = this->sides_on_edge(e);
1890 
1891  // All 3D elements should have exactly two sides adjacent to each edge.
1892  libmesh_assert_equal_to(adjacent_side_ids.size(), 2);
1893 
1894  // Get lists of node ids on each side
1895  const auto side_0_node_ids = this->nodes_on_side(adjacent_side_ids[0]);
1896  const auto side_1_node_ids = this->nodes_on_side(adjacent_side_ids[1]);
1897 
1898  // All 3D elements have at least three nodes on each side
1899  libmesh_assert_greater_equal(side_0_node_ids.size(), 3);
1900  libmesh_assert_greater_equal(side_1_node_ids.size(), 3);
1901 
1902  // Construct (approximate) inward normal on each adjacent side
1903  const auto side_0_normal =
1904  (this->point(side_0_node_ids[2]) - this->point(side_0_node_ids[0])).cross
1905  (this->point(side_0_node_ids[1]) - this->point(side_0_node_ids[0])).unit();
1906  const auto side_1_normal =
1907  (this->point(side_1_node_ids[2]) - this->point(side_1_node_ids[0])).cross
1908  (this->point(side_1_node_ids[1]) - this->point(side_1_node_ids[0])).unit();
1909 
1910  // Compute dihedral angle between the planes.
1911  // Using the absolute value ensures that we get the same result
1912  // even if the orientation of one of the planes is flipped, i.e.
1913  // it always gives us a value in the range [0, pi/2].
1914  // https://en.wikipedia.org/wiki/Dihedral_angle
1915  Real theta = std::acos(std::abs(side_0_normal * side_1_normal));
1916 
1917  // Track min and max angles seen
1918  min_angle = std::min(theta, min_angle);
1919  max_angle = std::max(theta, max_angle);
1920  }
1921 
1922  // Return requested extreme value (in degrees)
1923  return Real(180)/libMesh::pi * ((q == MIN_DIHEDRAL_ANGLE) ? min_angle : max_angle);
1924  }
1925 
1926  case JACOBIAN:
1927  case SCALED_JACOBIAN:
1928  {
1929  // 1D elements don't have interior corners, so this metric
1930  // does not really apply to them.
1931  const auto N = this->dim();
1932  if (N < 2)
1933  return 1.;
1934 
1935  // Initialize return value
1936  Real min_node_area = std::numeric_limits<Real>::max();
1937 
1938  for (auto n : this->node_index_range())
1939  {
1940  // Get list of edge ids adjacent to this node.
1941  auto adjacent_edge_ids = this->edges_adjacent_to_node(n);
1942 
1943  // Skip any nodes that don't have dim() adjacent edges. In 2D,
1944  // you need at least two adjacent edges to compute the cross
1945  // product, and in 3D, you need at least three adjacent edges
1946  // to compute the scalar triple product. The only element
1947  // type which has an unusual topology in this regard is the
1948  // Pyramid element in 3D, where 4 edges meet at the apex node.
1949  // For now, we just skip this node when computing the JACOBIAN
1950  // metric for Pyramids.
1951  if (adjacent_edge_ids.size() != N)
1952  continue;
1953 
1954  // Construct oriented edges
1955  std::vector<Point> oriented_edges(N);
1956  for (auto i : make_range(N))
1957  {
1958  auto node_0 = this->local_edge_node(adjacent_edge_ids[i], 0);
1959  auto node_1 = this->local_edge_node(adjacent_edge_ids[i], 1);
1960  if (node_0 != n)
1961  std::swap(node_0, node_1);
1962  oriented_edges[i] = this->point(node_1) - this->point(node_0);
1963  }
1964 
1965  // Compute unscaled area (2D) or volume (3D) using the
1966  // cross product (2D) or scalar triple product (3D) of the
1967  // oriented edges. We take the absolute value so that we
1968  // don't have to worry about the sign of the area.
1969  Real node_area = (N == 2) ?
1970  cross_norm(oriented_edges[0], oriented_edges[1]) :
1971  std::abs(triple_product(oriented_edges[0], oriented_edges[1], oriented_edges[2]));
1972 
1973  // Divide by (non-zero) edge lengths if computing scaled Jacobian.
1974  // If any length is zero, then set node_area to zero. This means that
1975  // e.g. degenerate quadrilaterals will have a zero SCALED_JACOBIAN metric.
1976  if (q == SCALED_JACOBIAN)
1977  for (auto i : make_range(N))
1978  {
1979  Real len_i = oriented_edges[i].norm();
1980  node_area = (len_i == 0.) ? 0. : (node_area / len_i);
1981  }
1982 
1983  // Update minimum
1984  min_node_area = std::min(node_area, min_node_area);
1985  }
1986 
1987  return min_node_area;
1988  }
1989 
1990  // Return 1 if we made it here
1991  default:
1992  {
1993  libmesh_do_once( libmesh_here();
1994 
1995  libMesh::err << "ERROR: quality metric "
1997  << " not implemented on element type "
1998  << Utility::enum_to_string(this->type())
1999  << std::endl
2000  << "Returning 1."
2001  << std::endl; );
2002 
2003  return 1.;
2004  }
2005  }
2006 }
2007 
2008 
2009 
2010 bool Elem::ancestor() const
2011 {
2012 #ifdef LIBMESH_ENABLE_AMR
2013 
2014  // Use a fast, DistributedMesh-safe definition
2015  const bool is_ancestor =
2016  !this->active() && !this->subactive();
2017 
2018  // But check for inconsistencies if we have time
2019 #ifdef DEBUG
2020  if (!is_ancestor && this->has_children())
2021  {
2022  for (auto & c : this->child_ref_range())
2023  {
2024  if (&c != remote_elem)
2025  {
2026  libmesh_assert(!c.active());
2027  libmesh_assert(!c.ancestor());
2028  }
2029  }
2030  }
2031 #endif // DEBUG
2032 
2033  return is_ancestor;
2034 
2035 #else
2036  return false;
2037 #endif
2038 }
2039 
2040 
2041 
2042 #ifdef LIBMESH_ENABLE_AMR
2043 
2044 void Elem::add_child (Elem * elem)
2045 {
2046  const unsigned int nc = this->n_children();
2047 
2048  if (!_children)
2049  {
2050  _children = std::make_unique<Elem *[]>(nc);
2051 
2052  for (unsigned int c = 0; c != nc; c++)
2053  this->set_child(c, nullptr);
2054  }
2055 
2056  for (unsigned int c = 0; c != nc; c++)
2057  {
2058  if (this->_children[c] == nullptr || this->_children[c] == remote_elem)
2059  {
2060  libmesh_assert_equal_to (this, elem->parent());
2061  this->set_child(c, elem);
2062  return;
2063  }
2064  }
2065 
2066  libmesh_error_msg("Error: Tried to add a child to an element with full children array");
2067 }
2068 
2069 
2070 
2071 void Elem::add_child (Elem * elem, unsigned int c)
2072 {
2073  if (!this->has_children())
2074  {
2075  const unsigned int nc = this->n_children();
2076  _children = std::make_unique<Elem *[]>(nc);
2077 
2078  for (unsigned int i = 0; i != nc; i++)
2079  this->set_child(i, nullptr);
2080  }
2081 
2082  libmesh_assert (this->_children[c] == nullptr || this->child_ptr(c) == remote_elem);
2083  libmesh_assert (elem == remote_elem || this == elem->parent());
2084 
2085  this->set_child(c, elem);
2086 }
2087 
2088 
2089 
2090 void Elem::replace_child (Elem * elem, unsigned int c)
2091 {
2092  libmesh_assert(this->has_children());
2093 
2094  libmesh_assert(this->child_ptr(c));
2095 
2096  this->set_child(c, elem);
2097 }
2098 
2099 
2100 
2101 void Elem::family_tree (std::vector<const Elem *> & family,
2102  bool reset) const
2103 {
2104  ElemInternal::family_tree(this, family, reset);
2105 }
2106 
2107 
2108 
2109 void Elem::family_tree (std::vector<Elem *> & family,
2110  bool reset)
2111 {
2112  ElemInternal::family_tree(this, family, reset);
2113 }
2114 
2115 
2116 
2117 void Elem::total_family_tree (std::vector<const Elem *> & family,
2118  bool reset) const
2119 {
2120  ElemInternal::total_family_tree(this, family, reset);
2121 }
2122 
2123 
2124 
2125 void Elem::total_family_tree (std::vector<Elem *> & family,
2126  bool reset)
2127 {
2128  ElemInternal::total_family_tree(this, family, reset);
2129 }
2130 
2131 
2132 
2133 void Elem::active_family_tree (std::vector<const Elem *> & active_family,
2134  bool reset) const
2135 {
2136  ElemInternal::active_family_tree(this, active_family, reset);
2137 }
2138 
2139 
2140 
2141 void Elem::active_family_tree (std::vector<Elem *> & active_family,
2142  bool reset)
2143 {
2144  ElemInternal::active_family_tree(this, active_family, reset);
2145 }
2146 
2147 
2148 
2149 void Elem::family_tree_by_side (std::vector<const Elem *> & family,
2150  unsigned int side,
2151  bool reset) const
2152 {
2153  ElemInternal::family_tree_by_side(this, family, side, reset);
2154 }
2155 
2156 
2157 
2158 void Elem:: family_tree_by_side (std::vector<Elem *> & family,
2159  unsigned int side,
2160  bool reset)
2161 {
2162  ElemInternal::family_tree_by_side(this, family, side, reset);
2163 }
2164 
2165 
2166 
2167 void Elem::active_family_tree_by_side (std::vector<const Elem *> & family,
2168  unsigned int side,
2169  bool reset) const
2170 {
2171  ElemInternal::active_family_tree_by_side(this, family, side, reset);
2172 }
2173 
2174 
2175 
2176 void Elem::active_family_tree_by_side (std::vector<Elem *> & family,
2177  unsigned int side,
2178  bool reset)
2179 {
2180  ElemInternal::active_family_tree_by_side(this, family, side, reset);
2181 }
2182 
2183 
2184 
2185 void Elem::family_tree_by_neighbor (std::vector<const Elem *> & family,
2186  const Elem * neighbor,
2187  bool reset) const
2188 {
2189  ElemInternal::family_tree_by_neighbor(this, family, neighbor, reset);
2190 }
2191 
2192 
2193 
2194 void Elem::family_tree_by_neighbor (std::vector<Elem *> & family,
2195  Elem * neighbor,
2196  bool reset)
2197 {
2198  ElemInternal::family_tree_by_neighbor(this, family, neighbor, reset);
2199 }
2200 
2201 
2202 
2203 void Elem::total_family_tree_by_neighbor (std::vector<const Elem *> & family,
2204  const Elem * neighbor,
2205  bool reset) const
2206 {
2207  ElemInternal::total_family_tree_by_neighbor(this, family, neighbor, reset);
2208 }
2209 
2210 
2211 
2212 void Elem::total_family_tree_by_neighbor (std::vector<Elem *> & family,
2213  Elem * neighbor,
2214  bool reset)
2215 {
2216  ElemInternal::total_family_tree_by_neighbor(this, family, neighbor, reset);
2217 }
2218 
2219 
2220 
2221 void Elem::family_tree_by_subneighbor (std::vector<const Elem *> & family,
2222  const Elem * neighbor,
2223  const Elem * subneighbor,
2224  bool reset) const
2225 {
2226  ElemInternal::family_tree_by_subneighbor(this, family, neighbor, subneighbor, reset);
2227 }
2228 
2229 
2230 
2231 void Elem::family_tree_by_subneighbor (std::vector<Elem *> & family,
2232  Elem * neighbor,
2233  Elem * subneighbor,
2234  bool reset)
2235 {
2236  ElemInternal::family_tree_by_subneighbor(this, family, neighbor, subneighbor, reset);
2237 }
2238 
2239 
2240 
2241 void Elem::total_family_tree_by_subneighbor (std::vector<const Elem *> & family,
2242  const Elem * neighbor,
2243  const Elem * subneighbor,
2244  bool reset) const
2245 {
2246  ElemInternal::total_family_tree_by_subneighbor(this, family, neighbor, subneighbor, reset);
2247 }
2248 
2249 
2250 
2251 void Elem::total_family_tree_by_subneighbor (std::vector<Elem *> & family,
2252  Elem * neighbor,
2253  Elem * subneighbor,
2254  bool reset)
2255 {
2256  ElemInternal::total_family_tree_by_subneighbor(this, family, neighbor, subneighbor, reset);
2257 }
2258 
2259 
2260 
2261 void Elem::active_family_tree_by_neighbor (std::vector<const Elem *> & family,
2262  const Elem * neighbor,
2263  bool reset) const
2264 {
2265  ElemInternal::active_family_tree_by_neighbor(this, family, neighbor, reset);
2266 }
2267 
2268 
2269 
2270 void Elem::active_family_tree_by_neighbor (std::vector<Elem *> & family,
2271  Elem * neighbor,
2272  bool reset)
2273 {
2274  ElemInternal::active_family_tree_by_neighbor(this, family, neighbor, reset);
2275 }
2276 
2277 
2278 
2279 void Elem::active_family_tree_by_topological_neighbor (std::vector<const Elem *> & family,
2280  const Elem * neighbor,
2281  const MeshBase & mesh,
2282  const PointLocatorBase & point_locator,
2283  const PeriodicBoundaries * pb,
2284  bool reset) const
2285 {
2287  mesh, point_locator, pb,
2288  reset);
2289 }
2290 
2291 
2292 
2293 void Elem::active_family_tree_by_topological_neighbor (std::vector<Elem *> & family,
2294  Elem * neighbor,
2295  const MeshBase & mesh,
2296  const PointLocatorBase & point_locator,
2297  const PeriodicBoundaries * pb,
2298  bool reset)
2299 {
2301  mesh, point_locator, pb,
2302  reset);
2303 }
2304 
2305 
2306 bool Elem::is_child_on_edge(const unsigned int c,
2307  const unsigned int e) const
2308 {
2309  libmesh_assert_less (c, this->n_children());
2310  libmesh_assert_less (e, this->n_edges());
2311 
2312  std::unique_ptr<const Elem> my_edge = this->build_edge_ptr(e);
2313  std::unique_ptr<const Elem> child_edge = this->child_ptr(c)->build_edge_ptr(e);
2314 
2315  // We're assuming that an overlapping child edge has the same
2316  // number and orientation as its parent
2317  return (child_edge->node_id(0) == my_edge->node_id(0) ||
2318  child_edge->node_id(1) == my_edge->node_id(1));
2319 }
2320 
2321 
2322 
2323 unsigned int Elem::min_p_level_by_neighbor(const Elem * neighbor_in,
2324  unsigned int current_min) const
2325 {
2326  libmesh_assert(!this->subactive());
2327  libmesh_assert(neighbor_in->active());
2328 
2329  // If we're an active element this is simple
2330  if (this->active())
2331  return std::min(current_min, this->p_level());
2332 
2333  libmesh_assert(has_neighbor(neighbor_in));
2334 
2335  // The p_level() of an ancestor element is already the minimum
2336  // p_level() of its children - so if that's high enough, we don't
2337  // need to examine any children.
2338  if (current_min <= this->p_level())
2339  return current_min;
2340 
2341  unsigned int min_p_level = current_min;
2342 
2343  for (auto & c : this->child_ref_range())
2344  if (&c != remote_elem && c.has_neighbor(neighbor_in))
2345  min_p_level =
2346  c.min_p_level_by_neighbor(neighbor_in, min_p_level);
2347 
2348  return min_p_level;
2349 }
2350 
2351 
2352 unsigned int Elem::min_new_p_level_by_neighbor(const Elem * neighbor_in,
2353  unsigned int current_min) const
2354 {
2355  libmesh_assert(!this->subactive());
2356  libmesh_assert(neighbor_in->active());
2357 
2358  // If we're an active element this is simple
2359  if (this->active())
2360  {
2361  unsigned int new_p_level = this->p_level();
2362  if (this->p_refinement_flag() == Elem::REFINE)
2363  new_p_level += 1;
2364  if (this->p_refinement_flag() == Elem::COARSEN)
2365  {
2366  libmesh_assert_greater (new_p_level, 0);
2367  new_p_level -= 1;
2368  }
2369  return std::min(current_min, new_p_level);
2370  }
2371 
2372  libmesh_assert(has_neighbor(neighbor_in));
2373 
2374  unsigned int min_p_level = current_min;
2375 
2376  for (auto & c : this->child_ref_range())
2377  if (&c != remote_elem && c.has_neighbor(neighbor_in))
2378  min_p_level =
2379  c.min_new_p_level_by_neighbor(neighbor_in, min_p_level);
2380 
2381  return min_p_level;
2382 }
2383 
2384 
2385 
2386 unsigned int Elem::as_parent_node (unsigned int child,
2387  unsigned int child_node) const
2388 {
2389  const unsigned int nc = this->n_children();
2390  libmesh_assert_less(child, nc);
2391 
2392  // Cached return values, indexed first by embedding_matrix version,
2393  // then by child number, then by child node number.
2394  std::vector<std::vector<std::vector<signed char>>> &
2395  cached_parent_indices = this->_get_parent_indices_cache();
2396 
2397  unsigned int em_vers = this->embedding_matrix_version();
2398 
2399  // We may be updating the cache on one thread, and while that
2400  // happens we can't safely access the cache from other threads.
2401  Threads::spin_mutex::scoped_lock lock(parent_indices_mutex);
2402 
2403  if (em_vers >= cached_parent_indices.size())
2404  cached_parent_indices.resize(em_vers+1);
2405 
2406  if (child >= cached_parent_indices[em_vers].size())
2407  {
2408  const signed char nn = cast_int<signed char>(this->n_nodes());
2409 
2410  cached_parent_indices[em_vers].resize(nc);
2411 
2412  for (unsigned int c = 0; c != nc; ++c)
2413  {
2414  const unsigned int ncn = this->n_nodes_in_child(c);
2415  cached_parent_indices[em_vers][c].resize(ncn);
2416  for (unsigned int cn = 0; cn != ncn; ++cn)
2417  {
2418  for (signed char n = 0; n != nn; ++n)
2419  {
2420  const Real em_val = this->embedding_matrix
2421  (c, cn, n);
2422  if (em_val == 1)
2423  {
2424  cached_parent_indices[em_vers][c][cn] = n;
2425  break;
2426  }
2427 
2428  if (em_val != 0)
2429  {
2430  cached_parent_indices[em_vers][c][cn] =
2431  -1;
2432  break;
2433  }
2434 
2435  // We should never see an all-zero embedding matrix
2436  // row
2437  libmesh_assert_not_equal_to (n+1, nn);
2438  }
2439  }
2440  }
2441  }
2442 
2443  const signed char cache_val =
2444  cached_parent_indices[em_vers][child][child_node];
2445  if (cache_val == -1)
2446  return libMesh::invalid_uint;
2447 
2448  return cached_parent_indices[em_vers][child][child_node];
2449 }
2450 
2451 
2452 
2453 const std::vector<std::pair<unsigned char, unsigned char>> &
2455  unsigned int child_node) const
2456 {
2457  // Indexed first by embedding matrix type, then by child id, then by
2458  // child node, then by bracketing pair
2459  std::vector<std::vector<std::vector<std::vector<std::pair<unsigned char, unsigned char>>>>> &
2460  cached_bracketing_nodes = this->_get_bracketing_node_cache();
2461 
2462  const unsigned int em_vers = this->embedding_matrix_version();
2463 
2464  // We may be updating the cache on one thread, and while that
2465  // happens we can't safely access the cache from other threads.
2466  Threads::spin_mutex::scoped_lock lock(parent_bracketing_nodes_mutex);
2467 
2468  if (cached_bracketing_nodes.size() <= em_vers)
2469  cached_bracketing_nodes.resize(em_vers+1);
2470 
2471  const unsigned int nc = this->n_children();
2472 
2473  // If we haven't cached the bracketing nodes corresponding to this
2474  // embedding matrix yet, let's do so now.
2475  if (cached_bracketing_nodes[em_vers].size() < nc)
2476  {
2477  // If we're a second-order element but we're not a full-order
2478  // element, then some of our bracketing nodes may not exist
2479  // except on the equivalent full-order element. Let's build an
2480  // equivalent full-order element and make a copy of its cache to
2481  // use.
2482  if (this->default_order() != FIRST &&
2483  second_order_equivalent_type(this->type(), /*full_ordered=*/ true) != this->type())
2484  {
2485  // Check that we really are the non-full-order type
2486  libmesh_assert_equal_to
2487  (second_order_equivalent_type (this->type(), false),
2488  this->type());
2489 
2490  // Build the full-order type
2491  ElemType full_type =
2492  second_order_equivalent_type(this->type(), /*full_ordered=*/ true);
2493  std::unique_ptr<Elem> full_elem = Elem::build(full_type);
2494 
2495  // This won't work for elements with multiple
2496  // embedding_matrix versions, but every such element is full
2497  // order anyways.
2498  libmesh_assert_equal_to(em_vers, 0);
2499 
2500  // Make sure its cache has been built. We temporarily
2501  // release our mutex lock so that the inner call can
2502  // re-acquire it.
2503  lock.release();
2504  full_elem->parent_bracketing_nodes(0,0);
2505 
2506  // And then we need to lock again, so that if someone *else*
2507  // grabbed our lock before we did we don't risk accessing
2508  // cached_bracketing_nodes while they're working on it.
2509  // Threading is hard.
2510  lock.acquire(parent_bracketing_nodes_mutex);
2511 
2512  // Copy its cache
2513  cached_bracketing_nodes =
2514  full_elem->_get_bracketing_node_cache();
2515 
2516  // Now we don't need to build the cache ourselves.
2517  return cached_bracketing_nodes[em_vers][child][child_node];
2518  }
2519 
2520  cached_bracketing_nodes[em_vers].resize(nc);
2521 
2522  const unsigned int nn = this->n_nodes();
2523 
2524  // We have to examine each child
2525  for (unsigned int c = 0; c != nc; ++c)
2526  {
2527  const unsigned int ncn = this->n_nodes_in_child(c);
2528 
2529  cached_bracketing_nodes[em_vers][c].resize(ncn);
2530 
2531  // We have to examine each node in that child
2532  for (unsigned int n = 0; n != ncn; ++n)
2533  {
2534  // If this child node isn't a vertex or an infinite
2535  // child element's mid-infinite-edge node, then we need
2536  // to find bracketing nodes on the child.
2537  if (!this->is_vertex_on_child(c, n)
2538 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2539  && !this->is_mid_infinite_edge_node(n)
2540 #endif
2541  )
2542  {
2543  // Use the embedding matrix to find the child node
2544  // location in parent master element space
2545  Point bracketed_pt;
2546 
2547  for (unsigned int pn = 0; pn != nn; ++pn)
2548  {
2549  const Real em_val =
2550  this->embedding_matrix(c,n,pn);
2551 
2552  libmesh_assert_not_equal_to (em_val, 1);
2553  if (em_val != 0.)
2554  bracketed_pt.add_scaled(this->master_point(pn), em_val);
2555  }
2556 
2557  // Check each pair of nodes on the child which are
2558  // also both parent nodes
2559  for (unsigned int n1 = 0; n1 != ncn; ++n1)
2560  {
2561  if (n1 == n)
2562  continue;
2563 
2564  unsigned int parent_n1 =
2565  this->as_parent_node(c,n1);
2566 
2567  if (parent_n1 == libMesh::invalid_uint)
2568  continue;
2569 
2570  Point p1 = this->master_point(parent_n1);
2571 
2572  for (unsigned int n2 = n1+1; n2 < nn; ++n2)
2573  {
2574  if (n2 == n)
2575  continue;
2576 
2577  unsigned int parent_n2 =
2578  this->as_parent_node(c,n2);
2579 
2580  if (parent_n2 == libMesh::invalid_uint)
2581  continue;
2582 
2583  Point p2 = this->master_point(parent_n2);
2584 
2585  Point pmid = (p1 + p2)/2;
2586 
2587  if (pmid == bracketed_pt)
2588  {
2589  cached_bracketing_nodes[em_vers][c][n].emplace_back(parent_n1, parent_n2);
2590  break;
2591  }
2592  else
2593  libmesh_assert(!pmid.absolute_fuzzy_equals(bracketed_pt));
2594  }
2595  }
2596  }
2597  // If this child node is a parent node, we need to
2598  // find bracketing nodes on the parent.
2599  else
2600  {
2601  unsigned int parent_node = this->as_parent_node(c,n);
2602 
2603  Point bracketed_pt;
2604 
2605  // If we're not a parent node, use the embedding
2606  // matrix to find the child node location in parent
2607  // master element space
2608  if (parent_node == libMesh::invalid_uint)
2609  {
2610  for (unsigned int pn = 0; pn != nn; ++pn)
2611  {
2612  const Real em_val =
2613  this->embedding_matrix(c,n,pn);
2614 
2615  libmesh_assert_not_equal_to (em_val, 1);
2616  if (em_val != 0.)
2617  bracketed_pt.add_scaled(this->master_point(pn), em_val);
2618  }
2619  }
2620  // If we're a parent node then we need no arithmetic
2621  else
2622  bracketed_pt = this->master_point(parent_node);
2623 
2624  for (unsigned int n1 = 0; n1 != nn; ++n1)
2625  {
2626  if (n1 == parent_node)
2627  continue;
2628 
2629  Point p1 = this->master_point(n1);
2630 
2631  for (unsigned int n2 = n1+1; n2 < nn; ++n2)
2632  {
2633  if (n2 == parent_node)
2634  continue;
2635 
2636  Point pmid = (p1 + this->master_point(n2))/2;
2637 
2638  if (pmid == bracketed_pt)
2639  {
2640  cached_bracketing_nodes[em_vers][c][n].emplace_back(n1, n2);
2641  break;
2642  }
2643  else
2644  libmesh_assert(!pmid.absolute_fuzzy_equals(bracketed_pt));
2645  }
2646  }
2647  }
2648  }
2649  }
2650  }
2651 
2652  return cached_bracketing_nodes[em_vers][child][child_node];
2653 }
2654 
2655 
2656 const std::vector<std::pair<dof_id_type, dof_id_type>>
2657 Elem::bracketing_nodes(unsigned int child,
2658  unsigned int child_node) const
2659 {
2660  std::vector<std::pair<dof_id_type, dof_id_type>> returnval;
2661 
2662  const std::vector<std::pair<unsigned char, unsigned char>> & pbc =
2663  this->parent_bracketing_nodes(child,child_node);
2664 
2665  for (const auto & pb : pbc)
2666  {
2667  const unsigned short n_n = this->n_nodes();
2668  if (pb.first < n_n && pb.second < n_n)
2669  returnval.emplace_back(this->node_id(pb.first), this->node_id(pb.second));
2670  else
2671  {
2672  // We must be on a non-full-order higher order element...
2673  libmesh_assert_not_equal_to(this->default_order(), FIRST);
2674  libmesh_assert_not_equal_to
2675  (second_order_equivalent_type (this->type(), true),
2676  this->type());
2677  libmesh_assert_equal_to
2678  (second_order_equivalent_type (this->type(), false),
2679  this->type());
2680 
2681  // And that's a shame, because this is a nasty search:
2682 
2683  // Build the full-order type
2684  ElemType full_type =
2685  second_order_equivalent_type(this->type(), /*full_ordered=*/ true);
2686  std::unique_ptr<Elem> full_elem = Elem::build(full_type);
2687 
2690 
2691  // Find the bracketing nodes by figuring out what
2692  // already-created children will have them.
2693 
2694  // This only doesn't break horribly because we add children
2695  // and nodes in straightforward + hierarchical orders...
2696  for (unsigned int c=0; c <= child; ++c)
2697  for (auto n : make_range(this->n_nodes_in_child(c)))
2698  {
2699  if (c == child && n == child_node)
2700  break;
2701 
2702  if (pb.first == full_elem->as_parent_node(c,n))
2703  {
2704  // We should be consistent
2705  if (pt1 != DofObject::invalid_id)
2706  libmesh_assert_equal_to(pt1, this->child_ptr(c)->node_id(n));
2707 
2708  pt1 = this->child_ptr(c)->node_id(n);
2709  }
2710 
2711  if (pb.second == full_elem->as_parent_node(c,n))
2712  {
2713  // We should be consistent
2714  if (pt2 != DofObject::invalid_id)
2715  libmesh_assert_equal_to(pt2, this->child_ptr(c)->node_id(n));
2716 
2717  pt2 = this->child_ptr(c)->node_id(n);
2718  }
2719  }
2720 
2721  // We should *usually* find all bracketing nodes by the time
2722  // we query them (again, because of the child & node add
2723  // order)
2724  //
2725  // The exception is if we're a HEX20, in which case we will
2726  // find pairs of vertex nodes and edge nodes bracketing the
2727  // new central node but we *won't* find the pairs of face
2728  // nodes which we would have had on a HEX27. In that case
2729  // we'll still have enough bracketing nodes for a
2730  // topological lookup, but we won't be able to make the
2731  // following assertions.
2732  if (this->type() != HEX20)
2733  {
2734  libmesh_assert_not_equal_to (pt1, DofObject::invalid_id);
2735  libmesh_assert_not_equal_to (pt2, DofObject::invalid_id);
2736  }
2737 
2738  if (pt1 != DofObject::invalid_id &&
2739  pt2 != DofObject::invalid_id)
2740  returnval.emplace_back(pt1, pt2);
2741  }
2742  }
2743 
2744  return returnval;
2745 }
2746 #endif // #ifdef LIBMESH_ENABLE_AMR
2747 
2748 
2749 
2750 
2751 bool Elem::contains_point (const Point & p, Real tol) const
2752 {
2753  // We currently allow the user to enlarge the bounding box by
2754  // providing a tol > TOLERANCE (so this routine is identical to
2755  // Elem::close_to_point()), but print a warning so that the
2756  // user can eventually switch his code over to calling close_to_point()
2757  // instead, which is intended to be used for this purpose.
2758  if (tol > TOLERANCE)
2759  {
2760  libmesh_do_once(libMesh::err
2761  << "WARNING: Resizing bounding box to match user-specified tolerance!\n"
2762  << "In the future, calls to Elem::contains_point() with tol > TOLERANCE\n"
2763  << "will be more optimized, but should not be used\n"
2764  << "to search for points 'close to' elements!\n"
2765  << "Instead, use Elem::close_to_point() for this purpose.\n"
2766  << std::endl;);
2767  return this->point_test(p, tol, tol);
2768  }
2769  else
2770  return this->point_test(p, TOLERANCE, tol);
2771 }
2772 
2773 
2774 
2775 
2776 bool Elem::close_to_point (const Point & p, Real tol) const
2777 {
2778  // This test uses the user's passed-in tolerance for the
2779  // bounding box test as well, thereby allowing the routine to
2780  // find points which are not only "in" the element, but also
2781  // "nearby" to within some tolerance.
2782  return this->point_test(p, tol, tol);
2783 }
2784 
2785 
2786 
2787 
2788 bool Elem::point_test(const Point & p, Real box_tol, Real map_tol) const
2789 {
2790  libmesh_assert_greater (box_tol, 0.);
2791  libmesh_assert_greater (map_tol, 0.);
2792 
2793  // This is a great optimization on first order elements, but it
2794  // could return false negatives on higher orders
2795  if (this->default_order() == FIRST)
2796  {
2797  // Check to make sure the element *could* contain this point, so we
2798  // can avoid an expensive inverse_map call if it doesn't.
2799  bool
2800 #if LIBMESH_DIM > 2
2801  point_above_min_z = false,
2802  point_below_max_z = false,
2803 #endif
2804 #if LIBMESH_DIM > 1
2805  point_above_min_y = false,
2806  point_below_max_y = false,
2807 #endif
2808  point_above_min_x = false,
2809  point_below_max_x = false;
2810 
2811  // For relative bounding box checks in physical space
2812  const Real my_hmax = this->hmax();
2813 
2814  for (auto & n : this->node_ref_range())
2815  {
2816  point_above_min_x = point_above_min_x || (n(0) - my_hmax*box_tol <= p(0));
2817  point_below_max_x = point_below_max_x || (n(0) + my_hmax*box_tol >= p(0));
2818 #if LIBMESH_DIM > 1
2819  point_above_min_y = point_above_min_y || (n(1) - my_hmax*box_tol <= p(1));
2820  point_below_max_y = point_below_max_y || (n(1) + my_hmax*box_tol >= p(1));
2821 #endif
2822 #if LIBMESH_DIM > 2
2823  point_above_min_z = point_above_min_z || (n(2) - my_hmax*box_tol <= p(2));
2824  point_below_max_z = point_below_max_z || (n(2) + my_hmax*box_tol >= p(2));
2825 #endif
2826  }
2827 
2828  if (
2829 #if LIBMESH_DIM > 2
2830  !point_above_min_z ||
2831  !point_below_max_z ||
2832 #endif
2833 #if LIBMESH_DIM > 1
2834  !point_above_min_y ||
2835  !point_below_max_y ||
2836 #endif
2837  !point_above_min_x ||
2838  !point_below_max_x)
2839  return false;
2840  }
2841 
2842  // To be on the safe side, we converge the inverse_map() iteration
2843  // to a slightly tighter tolerance than that requested by the
2844  // user...
2845  const Point mapped_point = FEMap::inverse_map(this->dim(),
2846  this,
2847  p,
2848  0.1*map_tol, // <- this is |dx| tolerance, the Newton residual should be ~ |dx|^2
2849  /*secure=*/ false,
2850  /*extra_checks=*/ false);
2851 
2852  // Check that the refspace point maps back to p! This is only necessary
2853  // for 1D and 2D elements, 3D elements always live in 3D.
2854  //
2855  // TODO: The contains_point() function could most likely be implemented
2856  // more efficiently in the element sub-classes themselves, at least for
2857  // the linear element types.
2858  if (this->dim() < 3)
2859  {
2860  Point xyz = FEMap::map(this->dim(), this, mapped_point);
2861 
2862  // Compute the distance between the original point and the re-mapped point.
2863  // They should be in the same place.
2864  Real dist = (xyz - p).norm();
2865 
2866 
2867  // If dist is larger than some fraction of the tolerance, then return false.
2868  // This can happen when e.g. a 2D element is living in 3D, and
2869  // FEMap::inverse_map() maps p onto the projection of the element,
2870  // effectively "tricking" on_reference_element().
2871  if (dist > this->hmax() * map_tol)
2872  return false;
2873  }
2874 
2875  return this->on_reference_element(mapped_point, map_tol);
2876 }
2877 
2878 
2879 
2880 
2881 bool Elem::has_invertible_map(Real /*tol*/) const
2882 {
2883  QNodal qnodal {this->dim()};
2884  FEMap fe_map;
2885  auto & jac = fe_map.get_jacobian();
2886 
2887  // We have a separate check for is_singular_node() below, so in this
2888  // case its "OK" to do nodal quadrature on pyramids.
2889  qnodal.allow_nodal_pyramid_quadrature = true;
2890  qnodal.init(*this);
2891  auto & qp = qnodal.get_points();
2892  libmesh_assert_equal_to(qp.size(), this->n_nodes());
2893 
2894  std::vector<Point> one_point(1);
2895  std::vector<Real> one_weight(1,1);
2896  for (auto i : index_range(qp))
2897  {
2898  if (this->is_singular_node(i))
2899  continue;
2900 
2901  one_point[0] = qp[i];
2902 
2903  fe_map.init_reference_to_physical_map(this->dim(), one_point, this);
2904  fe_map.compute_map(this->dim(), one_weight, this, false);
2905 
2906  if (jac[0] <= 0)
2907  return false;
2908  }
2909 
2910  return true;
2911 }
2912 
2913 
2914 
2915 void Elem::print_info (std::ostream & os) const
2916 {
2917  os << this->get_info()
2918  << std::endl;
2919 }
2920 
2921 
2922 
2923 std::string Elem::get_info () const
2924 {
2925  std::ostringstream oss;
2926 
2927  oss << " Elem Information" << '\n'
2928  << " id()=";
2929 
2930  if (this->valid_id())
2931  oss << this->id();
2932  else
2933  oss << "invalid";
2934 
2935 #ifdef LIBMESH_ENABLE_UNIQUE_ID
2936  oss << ", unique_id()=";
2937  if (this->valid_unique_id())
2938  oss << this->unique_id();
2939  else
2940  oss << "invalid";
2941 #endif
2942 
2943  oss << ", subdomain_id()=" << this->subdomain_id();
2944  oss << ", processor_id()=" << this->processor_id() << '\n';
2945 
2946  oss << " type()=" << Utility::enum_to_string(this->type()) << '\n'
2947  << " dim()=" << this->dim() << '\n'
2948  << " n_nodes()=" << this->n_nodes() << '\n';
2949 
2950  oss << " mapping=" << Utility::enum_to_string(this->mapping_type()) << '\n';
2951 
2952  for (auto n : this->node_index_range())
2953  {
2954  oss << " " << n << this->node_ref(n);
2955  if (this->mapping_type() == RATIONAL_BERNSTEIN_MAP)
2956  {
2957  const unsigned char datum_index = this->mapping_data();
2958  oss << " weight=" <<
2959  this->node_ref(n).get_extra_datum<Real>(datum_index) << '\n';
2960  }
2961  }
2962 
2963  oss << " n_sides()=" << this->n_sides() << '\n';
2964 
2965  for (auto s : this->side_index_range())
2966  {
2967  oss << " neighbor(" << s << ")=";
2968  if (this->neighbor_ptr(s))
2969  oss << this->neighbor_ptr(s)->id() << '\n';
2970  else
2971  oss << "nullptr\n";
2972  }
2973 
2974  if (!this->infinite())
2975  {
2976  oss << " hmin()=" << this->hmin()
2977  << ", hmax()=" << this->hmax() << '\n'
2978  << " volume()=" << this->volume() << '\n';
2979  }
2980  oss << " active()=" << this->active()
2981  << ", ancestor()=" << this->ancestor()
2982  << ", subactive()=" << this->subactive()
2983  << ", has_children()=" << this->has_children() << '\n'
2984  << " parent()=";
2985  if (this->parent())
2986  oss << this->parent()->id() << '\n';
2987  else
2988  oss << "nullptr\n";
2989  oss << " level()=" << this->level()
2990  << ", p_level()=" << this->p_level() << '\n'
2991 #ifdef LIBMESH_ENABLE_AMR
2992  << " refinement_flag()=" << Utility::enum_to_string(this->refinement_flag()) << '\n'
2993  << " p_refinement_flag()=" << Utility::enum_to_string(this->p_refinement_flag()) << '\n'
2994 #endif
2995 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2996  << " infinite()=" << this->infinite() << '\n';
2997  if (this->infinite())
2998  oss << " origin()=" << this->origin() << '\n'
2999 #endif
3000  ;
3001 
3002  oss << " DoFs=";
3003  for (auto s : make_range(this->n_systems()))
3004  for (auto v : make_range(this->n_vars(s)))
3005  for (auto c : make_range(this->n_comp(s,v)))
3006  oss << '(' << s << '/' << v << '/' << this->dof_number(s,v,c) << ") ";
3007 
3008 
3009  return oss.str();
3010 }
3011 
3012 
3013 
3015 {
3016  // Tell any of my neighbors about my death...
3017  // Looks strange, huh?
3018  for (auto n : this->side_index_range())
3019  {
3020  Elem * current_neighbor = this->neighbor_ptr(n);
3021  if (current_neighbor && current_neighbor != remote_elem)
3022  {
3023  // Note: it is possible that I see the neighbor
3024  // (which is coarser than me)
3025  // but they don't see me, so avoid that case.
3026  if (current_neighbor->level() == this->level())
3027  {
3028  const unsigned int w_n_a_i = current_neighbor->which_neighbor_am_i(this);
3029  libmesh_assert_less (w_n_a_i, current_neighbor->n_neighbors());
3030  current_neighbor->set_neighbor(w_n_a_i, nullptr);
3031  this->set_neighbor(n, nullptr);
3032  }
3033  }
3034  }
3035 }
3036 
3037 
3038 
3039 
3040 unsigned int Elem::n_second_order_adjacent_vertices (const unsigned int) const
3041 {
3042  // for linear elements, always return 0
3043  return 0;
3044 }
3045 
3046 
3047 
3048 unsigned short int Elem::second_order_adjacent_vertex (const unsigned int,
3049  const unsigned int) const
3050 {
3051  // for linear elements, always return 0
3052  return 0;
3053 }
3054 
3055 
3056 
3057 std::pair<unsigned short int, unsigned short int>
3058 Elem::second_order_child_vertex (const unsigned int) const
3059 {
3060  // for linear elements, always return 0
3061  return std::pair<unsigned short int, unsigned short int>(0,0);
3062 }
3063 
3064 
3065 
3067 {
3068  switch (et)
3069  {
3070  case NODEELEM:
3071  return NODEELEM;
3072  case EDGE2:
3073  case EDGE3:
3074  case EDGE4:
3075  return EDGE2;
3076  case TRI3:
3077  case TRI6:
3078  case TRI7:
3079  return TRI3;
3080  case TRISHELL3:
3081  return TRISHELL3;
3082  case QUAD4:
3083  case QUAD8:
3084  case QUAD9:
3085  return QUAD4;
3086  case QUADSHELL4:
3087  case QUADSHELL8:
3088  case QUADSHELL9:
3089  return QUADSHELL4;
3090  case TET4:
3091  case TET10:
3092  case TET14:
3093  return TET4;
3094  case HEX8:
3095  case HEX27:
3096  case HEX20:
3097  return HEX8;
3098  case PRISM6:
3099  case PRISM15:
3100  case PRISM18:
3101  case PRISM20:
3102  case PRISM21:
3103  return PRISM6;
3104  case PYRAMID5:
3105  case PYRAMID13:
3106  case PYRAMID14:
3107  case PYRAMID18:
3108  return PYRAMID5;
3109 
3110 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
3111 
3112  case INFEDGE2:
3113  return INFEDGE2;
3114  case INFQUAD4:
3115  case INFQUAD6:
3116  return INFQUAD4;
3117  case INFHEX8:
3118  case INFHEX16:
3119  case INFHEX18:
3120  return INFHEX8;
3121  case INFPRISM6:
3122  case INFPRISM12:
3123  return INFPRISM6;
3124 
3125 #endif
3126 
3127  default:
3128  // unknown element
3129  return INVALID_ELEM;
3130  }
3131 }
3132 
3133 
3134 
3136  const bool full_ordered)
3137 {
3138  switch (et)
3139  {
3140  case NODEELEM:
3141  return NODEELEM;
3142  case EDGE2:
3143  case EDGE3:
3144  {
3145  // full_ordered not relevant
3146  return EDGE3;
3147  }
3148 
3149  case EDGE4:
3150  {
3151  // full_ordered not relevant
3152  return EDGE4;
3153  }
3154 
3155  case TRI3:
3156  case TRI6:
3157  {
3158  // full_ordered not relevant
3159  return TRI6;
3160  }
3161 
3162  case TRI7:
3163  return TRI7;
3164 
3165  // Currently there is no TRISHELL6, so similarly to other types
3166  // where this is the case, we just return the input.
3167  case TRISHELL3:
3168  return TRISHELL3;
3169 
3170  case QUAD4:
3171  case QUAD8:
3172  {
3173  if (full_ordered)
3174  return QUAD9;
3175  else
3176  return QUAD8;
3177  }
3178 
3179  case QUADSHELL4:
3180  case QUADSHELL8:
3181  {
3182  if (full_ordered)
3183  return QUADSHELL9;
3184  else
3185  return QUADSHELL8;
3186  }
3187 
3188  case QUAD9:
3189  {
3190  // full_ordered not relevant
3191  return QUAD9;
3192  }
3193 
3194  case QUADSHELL9:
3195  {
3196  // full_ordered not relevant
3197  return QUADSHELL9;
3198  }
3199 
3200  case TET4:
3201  case TET10:
3202  {
3203  // full_ordered not relevant
3204  return TET10;
3205  }
3206 
3207  case TET14:
3208  return TET14;
3209 
3210  case HEX8:
3211  case HEX20:
3212  {
3213  // see below how this correlates with INFHEX8
3214  if (full_ordered)
3215  return HEX27;
3216  else
3217  return HEX20;
3218  }
3219 
3220  case HEX27:
3221  {
3222  // full_ordered not relevant
3223  return HEX27;
3224  }
3225 
3226  case PRISM6:
3227  case PRISM15:
3228  {
3229  if (full_ordered)
3230  return PRISM18;
3231  else
3232  return PRISM15;
3233  }
3234 
3235  // full_ordered not relevant, already fully second order
3236  case PRISM18:
3237  return PRISM18;
3238  case PRISM20:
3239  return PRISM20;
3240  case PRISM21:
3241  return PRISM21;
3242  case PYRAMID18:
3243  return PYRAMID18;
3244 
3245  case PYRAMID5:
3246  case PYRAMID13:
3247  {
3248  if (full_ordered)
3249  return PYRAMID14;
3250  else
3251  return PYRAMID13;
3252  }
3253 
3254  case PYRAMID14:
3255  {
3256  // full_ordered not relevant
3257  return PYRAMID14;
3258  }
3259 
3260 
3261 
3262 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
3263 
3264  // infinite elements
3265  case INFEDGE2:
3266  {
3267  return INFEDGE2;
3268  }
3269 
3270  case INFQUAD4:
3271  case INFQUAD6:
3272  {
3273  // full_ordered not relevant
3274  return INFQUAD6;
3275  }
3276 
3277  case INFHEX8:
3278  case INFHEX16:
3279  {
3280  /*
3281  * Note that this matches with \p Hex8:
3282  * For full-ordered, \p InfHex18 and \p Hex27
3283  * belong together, and for not full-ordered,
3284  * \p InfHex16 and \p Hex20 belong together.
3285  */
3286  if (full_ordered)
3287  return INFHEX18;
3288  else
3289  return INFHEX16;
3290  }
3291 
3292  case INFHEX18:
3293  {
3294  // full_ordered not relevant
3295  return INFHEX18;
3296  }
3297 
3298  case INFPRISM6:
3299  case INFPRISM12:
3300  {
3301  // full_ordered not relevant
3302  return INFPRISM12;
3303  }
3304 
3305 #endif
3306 
3307 
3308  default:
3309  {
3310  // what did we miss?
3311  libmesh_error_msg("No second order equivalent element type for et = "
3312  << Utility::enum_to_string(et));
3313  }
3314  }
3315 }
3316 
3317 
3318 
3320 {
3321  switch (et)
3322  {
3323  case NODEELEM:
3324  return NODEELEM;
3325  case EDGE2:
3326  case EDGE3:
3327  return EDGE3;
3328 
3329  case EDGE4:
3330  return EDGE4;
3331 
3332  case TRI3:
3333  case TRI6:
3334  case TRI7:
3335  return TRI7;
3336 
3337  // Currently there is no TRISHELL6, so similarly to other types
3338  // where this is the case, we just return the input.
3339  case TRISHELL3:
3340  return TRISHELL3;
3341 
3342  case QUAD4:
3343  case QUAD8:
3344  case QUAD9:
3345  return QUAD9;
3346 
3347  case QUADSHELL4:
3348  case QUADSHELL8:
3349  case QUADSHELL9:
3350  return QUADSHELL9;
3351 
3352  case TET4:
3353  case TET10:
3354  case TET14:
3355  return TET14;
3356 
3357  case HEX8:
3358  case HEX20:
3359  case HEX27:
3360  return HEX27;
3361 
3362  // We don't strictly need the 21st node for DoFs, but some code
3363  // depends on it for e.g. refinement patterns
3364  case PRISM6:
3365  case PRISM15:
3366  case PRISM18:
3367  case PRISM20:
3368  case PRISM21:
3369  return PRISM21;
3370 
3371  case PYRAMID5:
3372  case PYRAMID13:
3373  case PYRAMID14:
3374  case PYRAMID18:
3375  return PYRAMID18;
3376 
3377 
3378 
3379 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
3380  case INFEDGE2:
3381  return INFEDGE2;
3382 
3383  case INFQUAD4:
3384  case INFQUAD6:
3385  return INFQUAD6;
3386 
3387  case INFHEX8:
3388  case INFHEX16:
3389  case INFHEX18:
3390  return INFHEX18;
3391 
3392  // Probably ought to implement INFPRISM13 and/or 14; until then
3393  // we'll just return the not-complete-order ElemType, in
3394  // accordance which what seems like bad precedent...
3395  case INFPRISM6:
3396  case INFPRISM12:
3397  return INFPRISM12;
3398 #endif // LIBMESH_ENABLE_INFINITE_ELEMENTS
3399 
3400  default:
3401  {
3402  // what did we miss?
3403  libmesh_error_msg("No complete order equivalent element type for et = "
3404  << Utility::enum_to_string(et));
3405  }
3406  }
3407 }
3408 
3409 
3410 
3412 {
3414  return side_iterator(this->_first_side(), this->_last_side(), bsp);
3415 }
3416 
3417 
3418 
3419 
3421 {
3423  return side_iterator(this->_last_side(), this->_last_side(), bsp);
3424 }
3425 
3426 
3427 
3428 
3430 {
3431  // The default implementation builds a finite element of the correct
3432  // order and sums up the JxW contributions. This can be expensive,
3433  // so the various element types can overload this method and compute
3434  // the volume more efficiently.
3435  const FEFamily mapping_family = FEMap::map_fe_type(*this);
3436  const FEType fe_type(this->default_order(), mapping_family);
3437 
3438  std::unique_ptr<FEBase> fe (FEBase::build(this->dim(),
3439  fe_type));
3440 
3441  // Use a map with a negative Jacobian tolerance, in case we're asked
3442  // for the net volume of a tangled element
3443  fe->get_fe_map().set_jacobian_tolerance(std::numeric_limits<Real>::lowest());
3444 
3445  const std::vector<Real> & JxW = fe->get_JxW();
3446 
3447  // The default quadrature rule should integrate the mass matrix,
3448  // thus it should be plenty to compute the area
3449  QGauss qrule (this->dim(), fe_type.default_quadrature_order());
3450 
3451  fe->attach_quadrature_rule(&qrule);
3452 
3453  fe->reinit(this);
3454 
3455  Real vol=0.;
3456  for (auto jxw : JxW)
3457  vol += jxw;
3458 
3459  return vol;
3460 
3461 }
3462 
3463 
3464 
3466 {
3467  Point pmin = this->point(0);
3468  Point pmax = pmin;
3469 
3470  unsigned int n_points = this->n_nodes();
3471  for (unsigned int p=0; p != n_points; ++p)
3472  for (unsigned d=0; d<LIBMESH_DIM; ++d)
3473  {
3474  const Point & pt = this->point(p);
3475  if (pmin(d) > pt(d))
3476  pmin(d) = pt(d);
3477 
3478  if (pmax(d) < pt(d))
3479  pmax(d) = pt(d);
3480  }
3481 
3482  return BoundingBox(pmin, pmax);
3483 }
3484 
3485 
3486 
3487 bool Elem::is_vertex_on_parent(unsigned int c,
3488  unsigned int n) const
3489 {
3490 #ifdef LIBMESH_ENABLE_AMR
3491 
3492  unsigned int my_n_vertices = this->n_vertices();
3493  for (unsigned int n_parent = 0; n_parent != my_n_vertices;
3494  ++n_parent)
3495  if (this->node_ptr(n_parent) == this->child_ptr(c)->node_ptr(n))
3496  return true;
3497  return false;
3498 
3499 #else
3500 
3501  // No AMR?
3502  libmesh_ignore(c,n);
3503  libmesh_error_msg("ERROR: AMR disabled, how did we get here?");
3504  return true;
3505 
3506 #endif
3507 }
3508 
3509 
3510 
3511 unsigned int Elem::opposite_side(const unsigned int /*s*/) const
3512 {
3513  // If the subclass didn't rederive this, using it is an error
3514  libmesh_not_implemented();
3515 }
3516 
3517 
3518 
3519 unsigned int Elem::opposite_node(const unsigned int /*n*/,
3520  const unsigned int /*s*/) const
3521 {
3522  // If the subclass didn't rederive this, using it is an error
3523  libmesh_not_implemented();
3524 }
3525 
3526 
3527 unsigned int Elem::center_node_on_side(const unsigned short libmesh_dbg_var(side)) const
3528 {
3529  libmesh_assert_less (side, this->n_sides());
3530  return invalid_uint;
3531 }
3532 
3533 
3534 void Elem::swap2boundarysides(unsigned short s1,
3535  unsigned short s2,
3536  BoundaryInfo * boundary_info) const
3537 {
3538  std::vector<boundary_id_type> ids1, ids2;
3539  boundary_info->boundary_ids(this, s1, ids1);
3540  boundary_info->boundary_ids(this, s2, ids2);
3541  boundary_info->remove_side(this, s1);
3542  boundary_info->remove_side(this, s2);
3543  if (!ids1.empty())
3544  boundary_info->add_side(this, s2, ids1);
3545  if (!ids2.empty())
3546  boundary_info->add_side(this, s1, ids2);
3547 }
3548 
3549 
3550 void Elem::swap2boundaryedges(unsigned short e1,
3551  unsigned short e2,
3552  BoundaryInfo * boundary_info) const
3553 {
3554  std::vector<boundary_id_type> ids1, ids2;
3555  boundary_info->edge_boundary_ids(this, e1, ids1);
3556  boundary_info->edge_boundary_ids(this, e2, ids2);
3557  boundary_info->remove_edge(this, e1);
3558  boundary_info->remove_edge(this, e2);
3559  if (!ids1.empty())
3560  boundary_info->add_edge(this, e2, ids1);
3561  if (!ids2.empty())
3562  boundary_info->add_edge(this, e1, ids2);
3563 }
3564 
3565 bool
3566 Elem::is_internal(const unsigned int i) const
3567 {
3568  switch (this->dim())
3569  {
3570  case 0:
3571  return false;
3572 
3573  case 1:
3574  return !this->is_vertex(i);
3575 
3576  case 2:
3577  return !this->is_vertex(i) && !this->is_edge(i);
3578 
3579  case 3:
3580  return !this->is_vertex(i) && !this->is_edge(i) && !this->is_face(i);
3581 
3582  default:
3583  libmesh_error_msg("impossible element dimension " << std::to_string(this->dim()));
3584  return 0;
3585  }
3586 }
3587 
3588 bool
3589 Elem::positive_edge_orientation(const unsigned int i) const
3590 {
3591  libmesh_assert_less (i, this->n_edges());
3592 
3593  return this->point(this->local_edge_node(i, 0)) >
3594  this->point(this->local_edge_node(i, 1));
3595 }
3596 
3597 bool
3598 Elem::positive_face_orientation(const unsigned int i) const
3599 {
3600  libmesh_assert_less (i, this->n_faces());
3601 
3602  // Get the number of vertices N of face i. Note that for 3d elements, i.e.
3603  // elements for which this->n_faces() > 0, the number of vertices on any of
3604  // its sides (or faces) is just the number of that face's sides (or edges).
3605  auto side_i = this->side_ptr(i);
3606  const unsigned int N = side_i->n_sides();
3607 
3608  const std::vector<unsigned int> nodes = this->nodes_on_side(i);
3609 
3610  auto cmp = [&](const unsigned int & m, const unsigned int & n) -> bool
3611  { return this->point(m) < this->point(n); };
3612 
3613  const unsigned int V = std::distance(nodes.begin(),
3614  std::min_element(nodes.begin(), nodes.begin() + N, cmp));
3615 
3616  return cmp(nodes[(V + N - 1) % N], nodes[(V + 1) % N]);
3617 }
3618 
3619 } // namespace libMesh
virtual Point true_centroid() const
Definition: elem.C:594
void active_family_tree_by_topological_neighbor(std::vector< const Elem *> &family, const Elem *neighbor, const MeshBase &mesh, const PointLocatorBase &point_locator, const PeriodicBoundaries *pb, bool reset=true) const
Same as the active_family_tree_by_neighbor() member, but the neighbor here may be a topological (e...
Definition: elem.C:2279
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
OStreamProxy err
bool has_neighbor(const Elem *elem) const
Definition: elem.h:2628
unsigned char mapping_data() const
Definition: elem.h:3136
virtual bool is_vertex_on_parent(unsigned int c, unsigned int n) const
Definition: elem.C:3487
static const Order type_to_default_order_map[INVALID_ELEM]
This array maps the integer representation of the ElemType enum to the default approximation order of...
Definition: elem.h:991
RefinementState refinement_flag() const
Definition: elem.h:3210
ElemType
Defines an enum for geometric element types.
void swap2boundaryedges(unsigned short e1, unsigned short e2, BoundaryInfo *boundary_info) const
Swaps two edges in boundary_info, if it is non-null.
Definition: elem.C:3550
void total_family_tree(T elem, std::vector< T > &family, bool reset)
Definition: elem_internal.h:67
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
Definition: dof_object.h:1032
void write_connectivity(std::ostream &out, const IOPackage iop) const
Writes the element connectivity for various IO packages to the passed ostream "out".
Definition: elem.C:1731
const Elem * parent() const
Definition: elem.h:3030
void add_scaled(const TypeVector< T2 > &, const T &)
Add a scaled value to this vector without creating a temporary.
Definition: type_vector.h:629
virtual const std::vector< std::pair< dof_id_type, dof_id_type > > bracketing_nodes(unsigned int c, unsigned int n) const
Definition: elem.C:2657
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the element.
Definition: elem.C:2915
void active_family_tree_by_side(std::vector< const Elem *> &family, unsigned int side, bool reset=true) const
Same as the active_family_tree() member, but only adds elements which are next to side...
Definition: elem.C:2167
void remove_edge(const Elem *elem, const unsigned short int edge)
Removes all boundary conditions associated with edge edge of element elem, if any exist...
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:2245
virtual std::vector< unsigned int > sides_on_edge(const unsigned int) const =0
static ElemType complete_order_equivalent_type(const ElemType et)
Definition: elem.C:3319
virtual Point origin() const
Definition: elem.h:1914
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
unsigned int n_comp(const unsigned int s, const unsigned int var) const
Definition: dof_object.h:1002
virtual void compute_map(const unsigned int dim, const std::vector< Real > &qw, const Elem *elem, bool calculate_d2phi)
Compute the jacobian and some other additional data fields.
Definition: fe_map.C:1439
std::string get_info() const
Prints relevant information about the element to a string.
Definition: elem.C:2923
static const unsigned int type_to_n_sides_map[INVALID_ELEM]
This array maps the integer representation of the ElemType enum to the number of sides on the element...
Definition: elem.h:685
void find_interior_neighbors(T this_elem, std::set< T > &neighbor_set)
const Elem * interior_parent() const
Definition: elem.C:1186
const Elem * topological_neighbor(const unsigned int i, const MeshBase &mesh, const PointLocatorBase &point_locator, const PeriodicBoundaries *pb) const
Definition: elem.C:1313
virtual bool is_face(const unsigned int i) const =0
bool is_semilocal(const processor_id_type my_pid) const
Definition: elem.C:855
void total_family_tree_by_neighbor(T elem, std::vector< T > &family, T neighbor_in, bool reset=true)
This class implements nodal quadrature rules for various element types.
IntRange< unsigned short > side_index_range() const
Definition: elem.h:2710
void libmesh_assert_valid_node_pointers() const
Checks for a valid id and pointers to nodes with valid ids on this element.
Definition: elem.C:1371
void active_family_tree_by_neighbor(T elem, std::vector< T > &family, T neighbor_in, bool reset=true)
virtual dof_id_type key() const
Definition: elem.C:753
static constexpr Real TOLERANCE
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true, const bool extra_checks=true)
Definition: fe_map.C:1628
void active_family_tree(T elem, std::vector< T > &active_family, bool reset=true)
Definition: elem_internal.h:90
We&#39;re using a class instead of a typedef to allow forward declarations and future flexibility...
virtual unsigned int embedding_matrix_version() const
Definition: elem.h:2051
unsigned int which_side_am_i(const Elem *e) const
This function tells you which side the boundary element e is.
Definition: elem.C:870
T cross_norm(const TypeVector< T > &b, const TypeVector< T > &c)
Calls cross_norm_sq() and takes the square root of the result.
Definition: type_vector.h:1096
RefinementState p_refinement_flag() const
Definition: elem.h:3226
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
void find_interior_neighbors(std::set< const Elem *> &neighbor_set) const
This function finds all active elements (not including this one) in the parent manifold of this eleme...
Definition: elem.C:1172
static const unsigned int type_to_dim_map[INVALID_ELEM]
This array maps the integer representation of the ElemType enum to the geometric dimension of the ele...
Definition: elem.h:635
void family_tree_by_subneighbor(std::vector< const Elem *> &family, const Elem *neighbor, const Elem *subneighbor, bool reset=true) const
Same as the family_tree() member, but only adds elements which are next to subneighbor.
Definition: elem.C:2221
Threads::spin_mutex parent_indices_mutex
Definition: elem.C:97
void active_family_tree_by_side(T elem, std::vector< T > &family, unsigned int side, bool reset=true)
void family_tree(std::vector< const Elem *> &family, bool reset=true) const
Fills the vector family with the children of this element, recursively.
Definition: elem.C:2101
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
MeshBase & mesh
PeriodicBoundaryBase * boundary(boundary_id_type id)
void add_child(Elem *elem)
Adds a child pointer to the array of children of this element.
Definition: elem.C:2044
virtual BoundingBox loose_bounding_box() const
Definition: elem.C:3465
virtual bool on_reference_element(const Point &p, const Real eps=TOLERANCE) const =0
unsigned int min_new_p_level_by_neighbor(const Elem *neighbor, unsigned int current_min) const
Definition: elem.C:2352
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const =0
void active_family_tree(std::vector< const Elem *> &active_family, bool reset=true) const
Same as the family_tree() member, but only adds the active children.
Definition: elem.C:2133
unique_id_type unique_id() const
Definition: dof_object.h:844
void swap2boundarysides(unsigned short s1, unsigned short s2, BoundaryInfo *boundary_info) const
Swaps two sides in boundary_info, if it is non-null.
Definition: elem.C:3534
Order default_quadrature_order() const
Definition: fe_type.h:371
virtual unsigned int n_children() const =0
unsigned int p_level() const
Definition: elem.h:3108
void find_edge_neighbors(const Point &p1, const Point &p2, std::set< const Elem *> &neighbor_set) const
This function finds all active elements in the same manifold as this element which touch the current ...
Definition: elem.C:1083
side_iterator boundary_sides_end()
Definition: elem.C:3420
void boundary_ids(const Node *node, std::vector< boundary_id_type > &vec_to_fill) const
Fills a user-provided std::vector with the boundary ids associated with Node node.
void make_links_to_me_remote()
Resets this element&#39;s neighbors&#39; appropriate neighbor pointers and its parent&#39;s and children&#39;s approp...
Definition: elem.C:1533
The libMesh namespace provides an interface to certain functionality in the library.
virtual Real hmax() const
Definition: elem.C:722
Real distance(const Point &p)
unsigned int min_p_level_by_neighbor(const Elem *neighbor, unsigned int current_min) const
Definition: elem.C:2323
void set_interior_parent(Elem *p)
Sets the pointer to the element&#39;s interior_parent.
Definition: elem.C:1248
static const unsigned int type_to_n_nodes_map[INVALID_ELEM]
This array maps the integer representation of the ElemType enum to the number of nodes in the element...
Definition: elem.h:650
bool contains_vertex_of(const Elem *e, bool mesh_connection=false) const
Definition: elem.C:934
uint8_t processor_id_type
Definition: id_types.h:104
This is the MeshBase class.
Definition: mesh_base.h:75
SimpleRange< ChildRefIter > child_ref_range()
Returns a range with all children of a parent element, usable in range-based for loops.
Definition: elem.h:2346
void add(const TypeVector< T2 > &)
Add to this vector without creating a temporary.
Definition: type_vector.h:605
bool ancestor() const
Definition: elem.C:2010
bool has_topological_neighbor(const Elem *elem, const MeshBase &mesh, const PointLocatorBase &point_locator, const PeriodicBoundaries *pb) const
Definition: elem.C:1350
IntRange< unsigned short > edge_index_range() const
Definition: elem.h:2692
side_iterator boundary_sides_begin()
Iterator accessor functions.
Definition: elem.C:3411
virtual Real embedding_matrix(const unsigned int child_num, const unsigned int child_node_num, const unsigned int parent_node_num) const =0
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
T triple_product(const TypeVector< T > &a, const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1029
void family_tree_by_side(std::vector< const Elem *> &family, unsigned int side, bool reset=true) const
Same as the family_tree() member, but only adds elements which are next to side.
Definition: elem.C:2149
void replace_child(Elem *elem, unsigned int c)
Replaces the child pointer at the specified index in the child array.
Definition: elem.C:2090
void total_family_tree_by_subneighbor(std::vector< const Elem *> &family, const Elem *neighbor, const Elem *subneighbor, bool reset=true) const
Same as the family_tree_by_subneighbor() member, but also adds any subactive descendants.
Definition: elem.C:2241
TypeVector< T > unit() const
Definition: type_vector.h:1104
void libmesh_ignore(const Args &...)
static const subdomain_id_type invalid_subdomain_id
A static integral constant representing an invalid subdomain id.
Definition: elem.h:251
void remove_links_to_me()
Resets this element&#39;s neighbors&#39; appropriate neighbor pointers and its parent&#39;s and children&#39;s approp...
Definition: elem.C:1639
virtual std::vector< std::vector< std::vector< signed char > > > & _get_parent_indices_cache() const
Elem subclasses which don&#39;t do their own child-to-parent node calculations will need to supply a stat...
Definition: elem.h:2223
virtual bool contains_point(const Point &p, Real tol=TOLERANCE) const
Definition: elem.C:2751
Real length(const unsigned int n1, const unsigned int n2) const
Definition: elem.C:742
ElemMappingType mapping_type() const
Definition: elem.h:3120
void make_links_to_me_local(unsigned int n, unsigned int neighbor_side)
Resets the neighbor_side pointers of our nth neighbor (and its descendants, if appropriate) to point ...
Definition: elem.C:1455
const Node & node_ref(const unsigned int i) const
Definition: elem.h:2529
void family_tree_by_neighbor(T elem, std::vector< T > &family, T neighbor_in, bool reset=true)
dof_id_type id() const
Definition: dof_object.h:828
virtual Real hmin() const
Definition: elem.C:702
virtual unsigned int n_nodes() const =0
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const =0
unsigned int which_node_am_i(unsigned int side, unsigned int side_node) const
This function is deprecated, call local_side_node(side, side_node) instead.
Definition: elem.C:924
static const unsigned int type_to_n_edges_map[INVALID_ELEM]
This array maps the integer representation of the ElemType enum to the number of edges on the element...
Definition: elem.h:743
unsigned int n_vars(const unsigned int s, const unsigned int vg) const
Definition: dof_object.h:967
bool contains_edge_of(const Elem *e) const
Definition: elem.C:972
unsigned int which_neighbor_am_i(const Elem *e) const
This function tells you which neighbor e is.
Definition: elem.h:2919
virtual unsigned int as_parent_node(unsigned int c, unsigned int n) const
Definition: elem.C:2386
virtual std::vector< unsigned int > edges_adjacent_to_node(const unsigned int) const =0
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
Definition: elem.C:444
unsigned int n_systems() const
Definition: dof_object.h:937
virtual unsigned int opposite_node(const unsigned int n, const unsigned int s) const
Definition: elem.C:3519
virtual std::pair< unsigned short int, unsigned short int > second_order_child_vertex(const unsigned int n) const
Definition: elem.C:3058
Elem ** _elemlinks
Pointers to this element&#39;s parent and neighbors, and for lower-dimensional elements&#39; interior_parent...
Definition: elem.h:2251
void find_point_neighbors(const Point &p, std::set< const Elem *> &neighbor_set) const
This function finds all active elements (including this one) which are in the same manifold as this e...
Definition: elem.C:993
virtual unsigned int local_edge_node(unsigned int edge, unsigned int edge_node) const =0
Similar to Elem::local_side_node(), but instead of a side id, takes an edge id and a node id on that ...
The definition of the struct used for iterating over sides.
Definition: elem.h:3489
bool point_test(const Point &p, Real box_tol, Real map_tol) const
Shared private implementation used by the contains_point() and close_to_point() routines.
Definition: elem.C:2788
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
void family_tree(T elem, std::vector< T > &family, bool reset=true)
Definition: elem_internal.h:41
libmesh_assert(ctx)
static const dof_id_type invalid_id
An invalid id to distinguish an uninitialized DofObject.
Definition: dof_object.h:482
Threads::spin_mutex parent_bracketing_nodes_mutex
Definition: elem.C:98
const Elem * reference_elem() const
Definition: elem.C:570
void libmesh_assert_valid_neighbors() const
Checks for consistent neighbor links on this element.
Definition: elem.C:1383
virtual bool is_remote() const
Definition: elem.h:603
const std::vector< Real > & get_jacobian() const
Definition: fe_map.h:223
bool valid_unique_id() const
Definition: dof_object.h:893
This is the base class for point locators.
virtual unsigned int n_edges() const =0
ElemQuality
Defines an enum for element quality metrics.
virtual std::unique_ptr< Elem > disconnected_clone() const
Definition: elem.C:412
virtual unsigned int n_second_order_adjacent_vertices(const unsigned int n) const
Definition: elem.C:3040
SideIter _last_side()
Definition: elem.h:3478
auto norm(const T &a) -> decltype(std::abs(a))
Definition: tensor_tools.h:74
bool absolute_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
Definition: type_vector.h:972
bool positive_edge_orientation(const unsigned int i) const
Definition: elem.C:3589
void set_neighbor(const unsigned int i, Elem *n)
Assigns n as the neighbor.
Definition: elem.h:2618
static ElemType second_order_equivalent_type(const ElemType et, const bool full_ordered=true)
Definition: elem.C:3135
unsigned int which_child_am_i(const Elem *e) const
Definition: elem.h:3192
std::string enum_to_string(const T e)
void init_reference_to_physical_map(const std::vector< Point > &qp, const Elem *elem)
Definition: fe_map.C:109
void active_family_tree_by_neighbor(std::vector< const Elem *> &family, const Elem *neighbor, bool reset=true) const
Same as the active_family_tree() member, but only adds elements which are next to neighbor...
Definition: elem.C:2261
void total_family_tree_by_neighbor(std::vector< const Elem *> &family, const Elem *neighbor, bool reset=true) const
Same as the family_tree_by_neighbor() member, but also adds any subactive descendants.
Definition: elem.C:2203
SimpleRange< NodeRefIter > node_ref_range()
Returns a range with all nodes of an element, usable in range-based for loops.
Definition: elem.h:2665
Defines a Cartesian bounding box by the two corner extremum.
Definition: bounding_box.h:40
virtual unsigned int n_sides() const =0
bool topologically_equal(const Elem &rhs) const
Definition: elem.C:803
const Elem * neighbor_ptr(unsigned int i) const
Definition: elem.h:2598
void remove_side(const Elem *elem, const unsigned short int side)
Removes all boundary conditions associated with side side of element elem, if any exist...
void family_tree_by_subneighbor(T elem, std::vector< T > &family, T neighbor_in, T subneighbor, bool reset=true)
virtual bool close_to_point(const Point &p, Real tol) const
Definition: elem.C:2776
unsigned int level() const
Definition: elem.h:3074
virtual unsigned int n_vertices() const =0
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual std::unique_ptr< Elem > side_ptr(unsigned int i)=0
subdomain_id_type subdomain_id() const
Definition: elem.h:2582
bool valid_id() const
Definition: dof_object.h:885
virtual bool is_singular_node(unsigned int) const
Definition: elem.h:1838
virtual unsigned short dim() const =0
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:2507
virtual bool is_vertex(const unsigned int i) const =0
static Point map(const unsigned int dim, const Elem *elem, const Point &reference_point)
Definition: fe_map.C:2069
virtual Point master_point(const unsigned int i) const =0
static std::unique_ptr< Elem > build_with_id(const ElemType type, dof_id_type id)
Calls the build() method above with a nullptr parent, and additionally sets the newly-created Elem&#39;s ...
Definition: elem.C:558
unsigned int n_neighbors() const
Definition: elem.h:714
virtual Real quality(const ElemQuality q) const
Definition: elem.C:1774
void add_side(const dof_id_type elem, const unsigned short int side, const boundary_id_type id)
Add side side of element number elem with boundary id id to the boundary information data structure...
virtual unsigned short int second_order_adjacent_vertex(const unsigned int n, const unsigned int v) const
Definition: elem.C:3048
bool subactive() const
Definition: elem.h:2959
virtual Real volume() const
Definition: elem.C:3429
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
const Elem * neighbor(boundary_id_type boundary_id, const PointLocatorBase &point_locator, const Elem *e, unsigned int side, unsigned int *neigh_side=nullptr) const
This class implements specific orders of Gauss quadrature.
void set_child(unsigned int c, Elem *elem)
Sets the pointer to the child for this element.
Definition: elem.h:3182
virtual void connectivity(const unsigned int sc, const IOPackage iop, std::vector< dof_id_type > &conn) const =0
unsigned int size() const
Alias for n_points() to enable use in index_range.
Definition: quadrature.h:142
virtual bool has_invertible_map(Real tol=TOLERANCE *TOLERANCE) const
Definition: elem.C:2881
IntRange< unsigned short > node_index_range() const
Definition: elem.h:2683
bool positive_face_orientation(const unsigned int i) const
Definition: elem.C:3598
unsigned int n_extra_integers() const
Returns how many extra integers are associated to the DofObject.
Definition: dof_object.h:1170
virtual bool is_vertex_on_child(unsigned int, unsigned int n) const
Definition: elem.h:774
void edge_boundary_ids(const Elem *const elem, const unsigned short int edge, std::vector< boundary_id_type > &vec_to_fill) const
void find_point_neighbors(T this_elem, std::set< T > &neighbor_set, T start_elem)
void nullify_neighbors()
Replaces this element with nullptr for all of its neighbors.
Definition: elem.C:3014
SimpleRange< NeighborPtrIter > neighbor_ptr_range()
Returns a range with all neighbors of an element, usable in range-based for loops.
Definition: elem.h:3503
virtual unsigned int n_faces() const =0
T get_extra_datum(const unsigned int index) const
Gets the value on this object of the extra datum associated with index, which should have been obtain...
Definition: dof_object.h:1146
virtual bool infinite() const =0
FEFamily
defines an enum for finite element families.
virtual unsigned int n_sub_elem() const =0
void total_family_tree(std::vector< const Elem *> &family, bool reset=true) const
Same as the family_tree() member, but also adds any subactive descendants.
Definition: elem.C:2117
void family_tree_by_neighbor(std::vector< const Elem *> &family, const Elem *neighbor, bool reset=true) const
Same as the family_tree() member, but only adds elements which are next to neighbor.
Definition: elem.C:2185
virtual std::unique_ptr< Elem > build_edge_ptr(const unsigned int i)=0
virtual Order default_order() const =0
virtual const std::vector< std::pair< unsigned char, unsigned char > > & parent_bracketing_nodes(unsigned int c, unsigned int n) const
Definition: elem.C:2454
virtual std::vector< std::vector< std::vector< std::vector< std::pair< unsigned char, unsigned char > > > > > & _get_bracketing_node_cache() const
Elem subclasses which don&#39;t do their own bracketing node calculations will need to supply a static ca...
Definition: elem.h:2209
virtual unsigned int opposite_side(const unsigned int s) const
Definition: elem.C:3511
bool active() const
Definition: elem.h:2941
void family_tree_by_side(T elem, std::vector< T > &family, unsigned int s, bool reset)
Class contained in FE that encapsulates mapping (i.e.
Definition: fe_map.h:47
static ElemType first_order_equivalent_type(const ElemType et)
Definition: elem.C:3066
processor_id_type processor_id() const
Definition: dof_object.h:905
virtual ElemType type() const =0
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:2475
const Point & point(const unsigned int i) const
Definition: elem.h:2453
virtual Point centroid() const
Calls Elem::vertex_average() for backwards compatibility.
Definition: elem.C:578
bool operator==(const Elem &rhs) const
Definition: elem.C:771
bool has_children() const
Definition: elem.h:2979
virtual unsigned int center_node_on_side(const unsigned short side) const
Definition: elem.C:3527
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
std::unique_ptr< Elem *[]> _children
unique_ptr to array of this element&#39;s children.
Definition: elem.h:2261
virtual bool is_child_on_edge(const unsigned int c, const unsigned int e) const
Definition: elem.C:2306
virtual bool is_edge(const unsigned int i) const =0
const Elem & get(const ElemType type_in)
void total_family_tree_by_subneighbor(T elem, std::vector< T > &family, T neighbor_in, T subneighbor, bool reset=true)
virtual std::vector< unsigned int > nodes_on_side(const unsigned int) const =0
virtual unsigned int n_nodes_in_child(unsigned int) const
Definition: elem.h:674
Point vertex_average() const
Definition: elem.C:688
static const unsigned int max_n_nodes
The maximum number of nodes any element can contain.
Definition: elem.h:661
const Elem * child_ptr(unsigned int i) const
Definition: elem.h:3163
void active_family_tree_by_topological_neighbor(T elem, std::vector< T > &family, T neighbor_in, const MeshBase &mesh, const PointLocatorBase &point_locator, const PeriodicBoundaries *pb, bool reset=true)
uint8_t dof_id_type
Definition: id_types.h:67
void add_edge(const dof_id_type elem, const unsigned short int edge, const boundary_id_type id)
Add edge edge of element number elem with boundary id id to the boundary information data structure...
const Real pi
.
Definition: libmesh.h:299
static FEFamily map_fe_type(const Elem &elem)
Definition: fe_map.C:46
virtual bool is_mid_infinite_edge_node(const unsigned int) const
Definition: elem.h:1905
Used to iterate over the sides of an element which are on the boundary of the Mesh.
bool is_internal(const unsigned int i) const
Definition: elem.C:3566
SideIter _first_side()
Side iterator helper functions.
Definition: elem.h:3470
const RemoteElem * remote_elem
Definition: remote_elem.C:57