1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2020 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 
20 #ifndef LIBMESH_ELEM_H
21 #define LIBMESH_ELEM_H
22 
23 // Local includes
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/bounding_box.h"
26 #include "libmesh/dof_object.h"
27 #include "libmesh/id_types.h"
28 #include "libmesh/reference_counted_object.h"
29 #include "libmesh/node.h"
30 #include "libmesh/enum_elem_type.h" // INVALID_ELEM
31 #include "libmesh/multi_predicates.h"
32 #include "libmesh/pointer_to_pointer_iter.h"
33 #include "libmesh/int_range.h"
34 #include "libmesh/simple_range.h"
35 #include "libmesh/variant_filter_iterator.h"
36 #include "libmesh/hashword.h" // Used in compute_key() functions
37 
38 #ifdef LIBMESH_FORWARD_DECLARE_ENUMS
39 namespace libMesh
40 {
41 enum ElemQuality : int;
42 enum IOPackage : int;
43 enum Order : int;
44 }
45 #else
46 #include "libmesh/enum_elem_quality.h"
47 #include "libmesh/enum_io_package.h"
48 #include "libmesh/enum_order.h"
49 #endif
50 
51 // C++ includes
52 #include <algorithm>
53 #include <cstddef>
54 #include <iostream>
55 #include <limits.h> // CHAR_BIT
56 #include <set>
57 #include <vector>
58 #include <memory>
59 #include <array>
60 
61 namespace libMesh
62 {
63 
64 // Forward declarations
65 class MeshBase;
66 class MeshRefinement;
67 class Elem;
68 #ifdef LIBMESH_ENABLE_PERIODIC
69 class PeriodicBoundaries;
70 class PointLocatorBase;
71 #endif
72 template <class SideType, class ParentType>
73 class Side;
74 
75 
76 /**
77  * This is the base class from which all geometric element types are
78  * derived.  The \p Elem class provides standard information such as
79  * the number of nodes, edges, faces, vertices, children, and
80  * neighbors it has, as well as access to (or the ability to
81  * construct) these entities.
82  *
83  * An \p Elem has pointers to its \p Node objects.  Some of these
84  * nodes live at the vertices of the element, while others may live on
85  * edges (and faces in 3D), or interior to the element.  The number of
86  * nodes in a given element, \p n_nodes(), is encoded into the name of
87  * the class.  For example, a \p Tri3 has three nodes which correspond
88  * to the vertices, while a \p Tri6 has six nodes, three of which are
89  * located at vertices, and three which are located at the midpoint of
90  * each edge.  Nodes on edges, faces, and element interiors are called
91  * second-order nodes.
92  *
93  * A 1D Elem is an \p Edge, a 2D Elem is a \p Face, and a 3D Elem is a
94  * \p Cell.  An \p Elem is composed of a number of sides, which can
95  * be accessed as dim-1 dimensional \p Elem types.  For example, a \p
96  * Hex8 is a 3D hexahedral element. A \p Hex8 has 6 sides, which are
97  * \p Faces of type Quad4.
98  *
99  * \author Benjamin S. Kirk
100  * \date 2002-2007
101  * \brief The base class for all geometric element types.
102  */
103 class Elem : public ReferenceCountedObject<Elem>,
104              public DofObject
105 {
106 protected:
107 
108   /**
109    * Constructor.  Creates an element with \p n_nodes nodes,
110    * \p n_sides sides, \p n_children possible children, and
111    * parent \p p.  The constructor allocates the memory necessary
112    * to support this data.
113    */
114   Elem (const unsigned int n_nodes,
115         const unsigned int n_sides,
116         Elem * parent,
117         Elem ** elemlinkdata,
118         Node ** nodelinkdata);
119 
120 public:
121 
122   /**
123    * Elems are responsible for allocating and deleting their children
124    * during refinement, so they cannot be (default) copied or
125    * assigned. We therefore explicitly delete these operations.  Note
126    * that because _children is a C-style array, an Elem cannot even be
127    * safely default move-constructed (we would have to maintain a
128    * custom move constructor that explicitly sets _children to nullptr
129    * to do this safely).
130    */
131   Elem (Elem &&) = delete;
132   Elem (const Elem &) = delete;
133   Elem & operator= (const Elem &) = delete;
134   Elem & operator= (Elem &&) = delete;
135 
136   /**
137    * Destructor.  Frees all the memory associated with the element.
138    */
139   virtual ~Elem();
140 
141   /**
142    * \returns The \p Point associated with local \p Node \p i.
143    */
144   const Point & point (const unsigned int i) const;
145 
146   /**
147    * \returns The \p Point associated with local \p Node \p i
148    * as a writable reference.
149    */
150   Point & point (const unsigned int i);
151 
152   /**
153    * \returns The \p Point associated with local \p Node \p i,
154    * in master element rather than physical coordinates.
155    */
156   virtual Point master_point (const unsigned int i) const = 0;
157 
158   /**
159    * \returns The global id number of local \p Node \p i.
160    */
161   dof_id_type node_id (const unsigned int i) const;
162 
163   /**
164    * \returns The local id number of global \p Node id \p i,
165    * or \p invalid_uint if Node id \p i is not local.
166    */
167   unsigned int local_node (const dof_id_type i) const;
168 
169   /**
170    * \returns The local index for the \p Node pointer \p node_ptr,
171    * or \p invalid_id if \p node_ptr is not a local node.
172    */
173   unsigned int get_node_index (const Node * node_ptr) const;
174 
175   /**
176    * \returns A pointer to an array of local node pointers.
177    */
178   const Node * const * get_nodes () const;
179 
180   /**
181    * \returns A const pointer to local \p Node \p i.
182    */
183   const Node * node_ptr (const unsigned int i) const;
184 
185   /**
186    * \returns A non-const pointer to local \p Node \p i.
187    */
188   Node * node_ptr (const unsigned int i);
189 
190   /**
191    * \returns A const reference to local \p Node \p i.
192    */
193   const Node & node_ref (const unsigned int i) const;
194 
195   /**
196    * \returns A writable reference to local \p Node \p i.
197    */
198   Node & node_ref (const unsigned int i);
199 
200   /**
201    * \returns The pointer to local \p Node \p i as a writable reference.
202    */
203   virtual Node * & set_node (const unsigned int i);
204 
205   /**
206    * Nested classes for use iterating over all nodes of an element.
207    */
208   class NodeRefIter;
209   class ConstNodeRefIter;
210 
211   /**
212    * Returns a range with all nodes of an element, usable in
213    * range-based for loops.  The exact type of the return value here
214    * may be subject to change in future libMesh releases, but the
215    * iterators will always dereference to produce a reference to a
216    * Node.
217    */
218   SimpleRange<NodeRefIter> node_ref_range();
219 
220   SimpleRange<ConstNodeRefIter> node_ref_range() const;
221 
222   /**
223    * \returns The subdomain that this element belongs to.
224    */
225   subdomain_id_type subdomain_id () const;
226 
227   /**
228    * \returns The subdomain that this element belongs to as a
229    * writable reference.
230    */
231   subdomain_id_type & subdomain_id ();
232 
233   /**
234    * A static integral constant representing an invalid subdomain id.
235    * See also DofObject::{invalid_id, invalid_unique_id, invalid_processor_id}.
236    *
237    * \note We don't use the static_cast(-1) trick here since
238    * \p subdomain_id_type is sometimes a *signed* integer for
239    * compatibility reasons (see libmesh/id_types.h).
240    *
241    * \note Normally you can declare static const integral types
242    * directly in the header file (C++ standard, 9.4.2/4) but
243    * std::numeric_limits<T>::max() is not considered a "constant
244    * expression".  This one is therefore defined in elem.C.
245    * http://stackoverflow.com/questions/2738435/using-numeric-limitsmax-in-constant-expressions
246    */
247   static const subdomain_id_type invalid_subdomain_id;
248 
249   /**
250    * \returns A pointer to the "reference element" associated
251    * with this element.  The reference element is the image of this
252    * element in reference parametric space. Importantly, it is *not*
253    * an actual element in the mesh, but rather a Singleton-type
254    * object, so for example all \p Quad4 elements share the same
255    * \p reference_elem().
256    */
257   const Elem * reference_elem () const;
258 
259   /**
260    * \returns An id associated with the \p s side of this element.
261    * The id is not necessarily unique, but should be close.  This is
262    * particularly useful in the \p MeshBase::find_neighbors() routine.
263    */
264   virtual dof_id_type key (const unsigned int s) const = 0;
265 
266   /**
267    * \returns An id associated with the global node ids of this
268    * element.  The id is not necessarily unique, but should be
269    * close. Uses the same hash as the key(s) function, so for example
270    * if "tri3" is side 0 of "tet4", then tri3->key()==tet4->key(0).
271    */
272   virtual dof_id_type key () const;
273 
274   /**
275    * \returns \p true if two elements are identical, false otherwise.
276    * This is true if the elements are connected to identical global
277    * nodes, regardless of how those nodes might be numbered local
278    * to the elements.
279    */
280   bool operator == (const Elem & rhs) const;
281 
282   /**
283    * \returns A const pointer to the \f$ i^{th} \f$ neighbor of this
284    * element, or \p nullptr if \p MeshBase::find_neighbors() has not been
285    * called.
286    *
287    * \note If \p MeshBase::find_neighbors() has been called and this
288    * function still returns \p nullptr, then the side is on a boundary of
289    * the domain.
290    */
291   const Elem * neighbor_ptr (unsigned int i) const;
292 
293   /**
294    * \returns A non-const pointer to the \f$ i^{th} \f$ neighbor of this element.
295    */
296   Elem * neighbor_ptr (unsigned int i);
297 
298   /**
299    * Nested "classes" for use iterating over all neighbors of an element.
300    */
301   typedef Elem * const *       NeighborPtrIter;
302   typedef const Elem * const * ConstNeighborPtrIter;
303 
304   /**
305    * Returns a range with all neighbors of an element, usable in
306    * range-based for loops.  The exact type of the return value here
307    * may be subject to change in future libMesh releases, but the
308    * iterators will always dereference to produce a pointer to a
309    * neighbor element (or a null pointer, for sides which have no
310    * neighbors).
311    */
312   SimpleRange<NeighborPtrIter> neighbor_ptr_range();
313 
314   SimpleRange<ConstNeighborPtrIter> neighbor_ptr_range() const;
315 
316 #ifdef LIBMESH_ENABLE_PERIODIC
317   /**
318    * \returns A pointer to the \f$ i^{th} \f$ neighbor of this element
319    * for interior elements.  If an element is on a periodic
320    * boundary, it will return a corresponding element on the opposite
321    * side.
322    */
323   const Elem * topological_neighbor (const unsigned int i,
324                                      const MeshBase & mesh,
325                                      const PointLocatorBase & point_locator,
326                                      const PeriodicBoundaries * pb) const;
327 
328   /**
329    * \returns A writable pointer to the \f$ i^{th} \f$ neighbor of
330    * this element for interior elements.  If an element is on a
331    * periodic boundary, it will return a corresponding element on the
332    * opposite side.
333    */
334   Elem * topological_neighbor (const unsigned int i,
335                                MeshBase & mesh,
336                                const PointLocatorBase & point_locator,
337                                const PeriodicBoundaries * pb);
338 
339   /**
340    * \returns \p true if the element \p elem in question is a neighbor or
341    * topological neighbor of this element, \p false otherwise.
342    */
343   bool has_topological_neighbor (const Elem * elem,
344                                  const MeshBase & mesh,
345                                  const PointLocatorBase & point_locator,
346                                  const PeriodicBoundaries * pb) const;
347 #endif
348 
349   /**
350    * Assigns \p n as the \f$ i^{th} \f$ neighbor.
351    */
352   void set_neighbor (const unsigned int i, Elem * n);
353 
354   /**
355    * \returns \p true if the element \p elem in question is a neighbor
356    * of this element, \p false otherwise.
357    */
358   bool has_neighbor (const Elem * elem) const;
359 
360   /**
361    * \returns If \p elem is a neighbor of a child of this element, a
362    * pointer to that child, otherwise \p nullptr.
363    */
364   Elem * child_neighbor (Elem * elem);
365 
366   /**
367    * \returns If \p elem is a neighbor of a child of this element, a
368    * pointer to that child, otherwise \p nullptr.
369    */
370   const Elem * child_neighbor (const Elem * elem) const;
371 
372   /**
373    * \returns \p true if this element has a side coincident
374    * with a boundary (indicated by a \p nullptr neighbor), \p false
375    * otherwise.
376    */
377   bool on_boundary () const;
378 
379   /**
380    * \returns \p true if this element is semilocal to the calling
381    * processor, which must specify its rank.
382    */
383   bool is_semilocal (const processor_id_type my_pid) const;
384 
385   /**
386    * This function tells you which neighbor \p e is.
387    * I.e. if s = a->which_neighbor_am_i(e); then
388    * a->neighbor(s) will be an ancestor of e.
389    */
390   unsigned int which_neighbor_am_i(const Elem * e) const;
391 
392   /**
393    * This function tells you which side the boundary element \p e is.
394    * I.e. if e = a->build_side_ptr(s) or e = a->side_ptr(s); then
395    * a->which_side_am_i(e) will be s.
396    *
397    * \note An \e exact floating point comparison of the nodal
398    * positions of \p e is made with the nodal positions of \p this in
399    * order to perform this test. The idea is that the test will return
400    * a valid side id if \p e either directly shares Node pointers with
401    * \p this, or was created by exactly copying some of the nodes of
402    * \p this (e.g. through BoundaryMesh::sync()). In these
403    * circumstances, non-fuzzy floating point equality is expected.
404    *
405    * \returns The side of \p this the element which \p e is, otherwise
406    * \p invalid_uint.
407    */
408   unsigned int which_side_am_i(const Elem * e) const;
409 
410   /**
411    * \returns The local node id for node \p side_node on side \p side of
412    * this Elem. Simply relies on the \p side_nodes_map for each of the
413    * derived types. For example,
414    * Tri3::local_side_node(0, 0) -> 0
415    * Tri3::local_side_node(0, 1) -> 1
416    * Tri3::local_side_node(1, 0) -> 1
417    * Tri3::local_side_node(1, 1) -> 2
418    * etc...
419    */
420   virtual unsigned int local_side_node(unsigned int side,
421                                        unsigned int side_node) const = 0;
422 
423   /**
424    * Similar to Elem::local_side_node(), but instead of a side id, takes
425    * an edge id and a node id on that edge and returns a local node number
426    * for the Elem. The implementation relies on the "edge_nodes_map" tables
427    * for 3D elements. For 2D elements, calls local_side_node(). Throws an
428    * error if called on 1D elements.
429    */
430   virtual unsigned int local_edge_node(unsigned int edge,
431                                        unsigned int edge_node) const = 0;
432 
433   /**
434    * This function is deprecated, call local_side_node(side, side_node) instead.
435    */
436   unsigned int which_node_am_i(unsigned int side,
437                                unsigned int side_node) const;
438 
439   /**
440    * \returns \p true if a vertex of \p e is contained
441    * in this element.
442    */
443   bool contains_vertex_of(const Elem * e) const;
444 
445   /**
446    * \returns \p true if an edge of \p e is contained in
447    * this element.  (Internally, this is done by checking whether at
448    * least two vertices of \p e are contained in this element).
449    */
450   bool contains_edge_of(const Elem * e) const;
451 
452   /**
453    * This function finds all active elements (including this one)
454    * which are in the same manifold as this element and which touch
455    * the current active element at the specified point, which should
456    * be a point in the current element.
457    *
458    * Elements which are not "in the same manifold" (e.g. the
459    * interior_parent of a boundary element) will not be found with
460    * this method.
461    *
462    * Elements which overlap the specified point but which are only
463    * connected to the current element via elements which do not
464    * overlap that point (e.g. in a folded or tangled mesh) are not
465    * considered to "touch" the current element and will not be found
466    * with this method.
467    */
468   void find_point_neighbors(const Point & p,
469                             std::set<const Elem *> & neighbor_set) const;
470 
471   /**
472    * This function finds all active elements (including this one) in
473    * the same manifold as this element which touch this active element
474    * at any point.
475    */
476   void find_point_neighbors(std::set<const Elem *> & neighbor_set) const;
477 
478   /**
479    * This function finds all active elements (including this one) in
480    * the same manifold as start_elem (which must be active and must
481    * touch this element) which touch this element at any point.
482    */
483   void find_point_neighbors(std::set<const Elem *> & neighbor_set,
484                             const Elem * start_elem) const;
485 
486   /**
487    * Non-const version of function above. Fills a set of non-const Elem pointers.
488    */
489   void find_point_neighbors(std::set<Elem *> & neighbor_set,
490                             Elem * start_elem);
491 
492   /**
493    * This function finds all active elements in the same manifold as
494    * this element which touch the current active element along the
495    * whole edge defined by the two points \p p1 and \p p2.
496    */
497   void find_edge_neighbors(const Point & p1,
498                            const Point & p2,
499                            std::set<const Elem *> & neighbor_set) const;
500 
501   /**
502    * This function finds all active elements in the same manifold as
503    * this element which touch the current active element along any
504    * edge (more precisely, at at least two points).
505    *
506    * In this case, elements are included even if they do not touch a
507    * *whole* edge of this element.
508    */
509   void find_edge_neighbors(std::set<const Elem *> & neighbor_set) const;
510 
511   /**
512    * This function finds all active elements (*not* including this
513    * one) in the parent manifold of this element whose intersection
514    * with this element has non-zero measure.
515    */
516   void find_interior_neighbors(std::set<const Elem *> & neighbor_set) const;
517 
518   /**
519    * Non-const version of function above that fills up a vector of
520    * non-const Elem pointers instead.
521    */
522   void find_interior_neighbors(std::set<Elem *> & neighbor_set);
523 
524   /**
525    * Resets this element's neighbors' appropriate neighbor pointers
526    * and its parent's and children's appropriate pointers
527    * to point to null instead of to this.
528    *
529    * Used by the library before an element is deleted from a mesh.
530    */
531   void remove_links_to_me ();
532 
533   /**
534    * Resets this element's neighbors' appropriate neighbor pointers
535    * and its parent's and children's appropriate pointers
536    * to point to the global remote_elem instead of this.
537    * Used by the library before an element becomes remote on the
538    * local processor.
539    */
540   void make_links_to_me_remote ();
541 
542   /**
543    * Resets the \p neighbor_side pointers of our nth neighbor (and
544    * its descendants, if appropriate) to point to this Elem instead of
545    * to the global remote_elem.  Used by the library when a formerly
546    * remote element is being added to the local processor.
547    */
548   void make_links_to_me_local (unsigned int n, unsigned int neighbor_side);
549 
550   /**
551    * \returns \p true if this element is remote, false otherwise.
552    *
553    * A remote element (see \p RemoteElem) is a syntactic convenience --
554    * it is a placeholder for an element which exists on some other
555    * processor.  Local elements are required to have valid neighbors,
556    * and these ghost elements may have remote neighbors for data
557    * structure consistency.  The use of remote elements helps ensure
558    * that any element we may access has a \p nullptr neighbor only if it
559    * lies on the physical boundary of the domain.
560    */
is_remote()561   virtual bool is_remote () const
562   { return false; }
563 
564   /**
565    * \returns The connectivity for this element in a specific
566    * format, which is specified by the IOPackage tag.
567    */
568   virtual void connectivity(const unsigned int sc,
569                             const IOPackage iop,
570                             std::vector<dof_id_type> & conn) const = 0;
571 
572   /**
573    * Writes the element connectivity for various IO packages
574    * to the passed ostream "out".  Not virtual, since it is
575    * implemented in the base class.
576    */
577   void write_connectivity (std::ostream & out,
578                            const IOPackage iop) const;
579 
580   /**
581    * \returns The type of element that has been derived from this
582    * base class.
583    */
584   virtual ElemType type () const = 0;
585 
586   /**
587    * \returns The dimensionality of the object.
588    */
589   virtual unsigned short dim () const = 0;
590 
591   /**
592    * This array maps the integer representation of the \p ElemType enum
593    * to the number of nodes in the element.
594    */
595   static const unsigned int type_to_n_nodes_map[INVALID_ELEM];
596 
597   /**
598    * \returns The number of nodes this element contains.
599    */
600   virtual unsigned int n_nodes () const = 0;
601 
602   /**
603    * The maximum number of nodes *any* element can contain.
604    * This is useful for replacing heap vectors with stack arrays.
605    */
606   static const unsigned int max_n_nodes = 27;
607 
608   /**
609    * \returns An integer range from 0 up to (but not including)
610    * the number of nodes this element contains.
611    */
612   IntRange<unsigned short> node_index_range () const;
613 
614   /**
615    * \returns The number of nodes the given child of this element
616    * contains.  Except in odd cases like pyramid refinement this will
617    * be the same as the number of nodes in the parent element.
618    */
n_nodes_in_child(unsigned int)619   virtual unsigned int n_nodes_in_child (unsigned int /*c*/) const
620   { return this->n_nodes(); }
621 
622   /**
623    * This array maps the integer representation of the \p ElemType enum
624    * to the number of sides on the element.
625    */
626   static const unsigned int type_to_n_sides_map[INVALID_ELEM];
627 
628   /**
629    * \returns The number of sides the element that has been derived
630    * from this class has. In 2D the number of sides is the number
631    * of edges, in 3D the number of sides is the number of faces.
632    */
633   virtual unsigned int n_sides () const = 0;
634 
635   /**
636    * \returns An integer range from 0 up to (but not including)
637    * the number of sides this element has.
638    */
639   IntRange<unsigned short> side_index_range () const;
640 
641   /**
642    * \returns The number of neighbors the element that has been derived
643    * from this class has.
644    *
645    * Only face (or edge in 2D) neighbors are stored, so this method
646    * returns n_sides().  At one point we intended to allow derived
647    * classes to override this, but too much current libMesh code
648    * assumes n_neighbors==n_sides.
649    */
n_neighbors()650   unsigned int n_neighbors () const
651   { return this->n_sides(); }
652 
653   /**
654    * \returns The number of vertices the element that has been derived
655    * from this class has.
656    */
657   virtual unsigned int n_vertices () const = 0;
658 
659   /**
660    * \returns The number of edges the element that has been derived
661    * from this class has.
662    */
663   virtual unsigned int n_edges () const = 0;
664 
665   /**
666    * \returns An integer range from 0 up to (but not including)
667    * the number of edges this element has.
668    */
669   IntRange<unsigned short> edge_index_range () const;
670 
671   /**
672    * This array maps the integer representation of the \p ElemType enum
673    * to the number of edges on the element.
674    */
675   static const unsigned int type_to_n_edges_map[INVALID_ELEM];
676 
677   /**
678    * \returns The number of faces the element that has been derived
679    * from this class has.
680    */
681   virtual unsigned int n_faces () const = 0;
682 
683   /**
684    * \returns The number of children the element that has been derived
685    * from this class may have.
686    */
687   virtual unsigned int n_children () const = 0;
688 
689   /**
690    * \returns \p true if the specified (local) node number is a vertex.
691    */
692   virtual bool is_vertex(const unsigned int i) const = 0;
693 
694   /**
695    * \returns \p true if the specified child has a vertex at the
696    * specified (child-local) node number.
697    * Except in odd cases like pyramid refinement the child will have
698    * the same local structure as the parent element.
699    */
is_vertex_on_child(unsigned int,unsigned int n)700   virtual unsigned int is_vertex_on_child (unsigned int /*c*/,
701                                            unsigned int n) const
702   { return this->is_vertex(n); }
703 
704   /**
705    * \returns \p true if this element has a vertex at the specified
706    * (child-local) node number \p n of the specified child \p c.
707    */
708   virtual bool is_vertex_on_parent(unsigned int c,
709                                    unsigned int n) const;
710 
711   /**
712    * \returns \p true if the specified (local) node number is an edge.
713    */
714   virtual bool is_edge(const unsigned int i) const = 0;
715 
716   /**
717    * \returns \p true if the specified (local) node number is a face.
718    */
719   virtual bool is_face(const unsigned int i) const = 0;
720 
721   /**
722    * \returns \p true if the specified (local) node number is on the
723    * specified side.
724    */
725   virtual bool is_node_on_side(const unsigned int n,
726                                const unsigned int s) const = 0;
727 
728   /**
729    * \returns the (local) node numbers on the specified side
730    */
731   virtual std::vector<unsigned int> nodes_on_side(const unsigned int /*s*/) const = 0;
732 
733   /**
734    * \returns the (local) node numbers on the specified edge
735    */
736   virtual std::vector<unsigned int> nodes_on_edge(const unsigned int /*e*/) const = 0;
737 
738   /**
739    * \returns the (local) side numbers that touch the specified edge
740    */
741   virtual std::vector<unsigned int> sides_on_edge(const unsigned int /*e*/) const = 0;
742 
743   /**
744    * \returns \p true if the specified (local) node number is on the
745    * specified edge.
746    */
747   virtual bool is_node_on_edge(const unsigned int n,
748                                const unsigned int e) const = 0;
749 
750   /**
751    * \returns \p true if the specified edge is on the specified side.
752    */
753   virtual bool is_edge_on_side(const unsigned int e,
754                                const unsigned int s) const = 0;
755 
756   /**
757    * \returns The side number opposite to \p s (for a tensor product
758    * element), or throws an error otherwise.
759    */
760   virtual unsigned int opposite_side(const unsigned int s) const;
761 
762   /**
763    * \returns The local node number for the node opposite to node n
764    * on side \p opposite_side(s) (for a tensor product element), or
765    * throws an error otherwise.
766    */
767   virtual unsigned int opposite_node(const unsigned int n,
768                                      const unsigned int s) const;
769 
770   /**
771    * \returns The number of sub-elements this element may be broken
772    * down into for visualization purposes.  For example, 1 for a
773    * linear triangle, 4 for a quadratic (6-noded) triangle, etc...
774    */
775   virtual unsigned int n_sub_elem () const = 0;
776 
777   /**
778    * \returns A proxy element coincident with side \p i.
779    *
780    * \note This method returns the _minimum_ element necessary to
781    * uniquely identify the side.  For example, the side of a
782    * hexahedron is always returned as a 4-noded quadrilateral,
783    * regardless of what type of hex you are dealing with.  If you want
784    * the full-ordered face (i.e. a 9-noded quad face for a 27-noded
785    * hexahedron) use the build_side method.
786    *
787    * \note The const version of this function is non-virtual; it
788    * simply calls the virtual non-const version and const_casts the
789    * return type.
790    */
791   virtual std::unique_ptr<Elem> side_ptr (unsigned int i) = 0;
792   std::unique_ptr<const Elem> side_ptr (unsigned int i) const;
793 
794   /**
795    * Resets the loose element \p side, which may currently point to a
796    * different side than \p i or even a different element than \p
797    * this, to point to side \p i on \p this.  If \p side is currently
798    * an element of the wrong type, it will be freed and a new element
799    * allocated; otherwise no memory allocation will occur.
800    *
801    * This should not be called with proxy Side elements.  This will
802    * cause \p side to be a minimum-ordered element, even if it is
803    * handed a higher-ordered element that must be replaced.
804    *
805    * The const version of this function is non-virtual; it simply
806    * calls the virtual non-const version and const_casts the return
807    * type.
808    */
809   virtual void side_ptr (std::unique_ptr<Elem> & side, const unsigned int i) = 0;
810   void side_ptr (std::unique_ptr<const Elem> & side, const unsigned int i) const;
811 
812   /**
813    * \returns An element coincident with side \p i wrapped in a smart pointer.
814    *
815    * The element returned is full-ordered, in contrast to the side
816    * method.  For example, calling build_side_ptr(0) on a 20-noded hex
817    * will build a 8-noded quadrilateral coincident with face 0 and
818    * pass back the pointer.  A \p std::unique_ptr<Elem> is returned to
819    * prevent a memory leak.  This way the user need not remember to
820    * delete the object.
821    *
822    * The second argument, which is true by default, specifies that a
823    * "proxy" element (of type Side) will be returned.  This type of
824    * value is useful because it does not allocate additional
825    * memory, and is usually sufficient for FE calculation purposes.
826    * If you really need a full-ordered, non-proxy side object, call
827    * this function with proxy=false.
828    *
829    * The const version of this function is non-virtual; it simply
830    * calls the virtual non-const version and const_casts the return
831    * type.
832    */
833   virtual std::unique_ptr<Elem> build_side_ptr (const unsigned int i, bool proxy=true) = 0;
834   std::unique_ptr<const Elem> build_side_ptr (const unsigned int i, bool proxy=true) const;
835 
836   /**
837    * Resets the loose element \p side, which may currently point to a
838    * different side than \p i or even a different element than \p
839    * this, to point to side \p i on \p this.  If \p side is currently
840    * an element of the wrong type, it will be freed and a new element
841    * allocated; otherwise no memory allocation will occur.
842    *
843    * This should not be called with proxy Side elements.  This will
844    * cause \p side to be a full-ordered element, even if it is handed
845    * a lower-ordered element that must be replaced.
846    *
847    * The const version of this function is non-virtual; it simply
848    * calls the virtual non-const version and const_casts the return
849    * type.
850    */
851   virtual void build_side_ptr (std::unique_ptr<Elem> & side, const unsigned int i) = 0;
852   void build_side_ptr (std::unique_ptr<const Elem> & side, const unsigned int i) const;
853 
854   /**
855    * \returns An element coincident with edge \p i wrapped in a smart pointer.
856    *
857    * The element returned is full-ordered.  For example, calling
858    * build_edge_ptr(0) on a 20-noded hex will build a 3-noded edge
859    * coincident with edge 0 and pass back the pointer.  A \p
860    * std::unique_ptr<Elem> is returned to prevent a memory leak.  This way
861    * the user need not remember to delete the object.
862    *
863    * The const version of this function is non-virtual; it simply
864    * calls the virtual non-const version and const_casts the return
865    * type.
866    */
867   virtual std::unique_ptr<Elem> build_edge_ptr (const unsigned int i) = 0;
868   std::unique_ptr<const Elem> build_edge_ptr (const unsigned int i) const;
869 
870   /**
871    * \returns The default approximation order for this element type.
872    * This is the order that will be used to compute the map to the
873    * reference element.
874    */
875   virtual Order default_order () const = 0;
876 
877   /**
878    * \returns The centroid of the element. The centroid is
879    * computed as the average of all the element vertices.
880    *
881    * This method is virtual since some derived elements
882    * might want to use shortcuts to compute their centroid.
883    */
884   virtual Point centroid () const;
885 
886   /**
887    * \returns The minimum vertex separation for the element.
888    */
889   virtual Real hmin () const;
890 
891   /**
892    * \returns The maximum vertex separation for the element.
893    */
894   virtual Real hmax () const;
895 
896   /**
897    * \returns The (length/area/volume) of the geometric element.
898    */
899   virtual Real volume () const;
900 
901   /**
902    * \returns A bounding box (not necessarily the minimal bounding box)
903    * containing the geometric element.
904    *
905    * The base class implementation determines a bounding box for the
906    * element *nodes*, which should be sufficient for first order
907    * finite elements.  Higher order geometric elements will need to
908    * override with an implementation which takes curved elements into
909    * account.
910    */
911   virtual BoundingBox loose_bounding_box () const;
912 
913   /**
914    * \returns A quantitative assessment of element quality based on
915    * the quality metric \p q specified by the user.
916    */
917   virtual Real quality (const ElemQuality q) const;
918 
919   /**
920    * \returns The suggested quality bounds for the Elem based on
921    * quality measure \p q.
922    *
923    * These are the values suggested by the CUBIT User's Manual.  Since
924    * this function can have no possible meaning for an abstract Elem,
925    * it is an error in the base class.
926    */
qual_bounds(const ElemQuality)927   virtual std::pair<Real,Real> qual_bounds (const ElemQuality) const
928   { libmesh_not_implemented(); return std::make_pair(0.,0.); }
929 
930   /**
931    * \returns \p true if the point p is contained in this element,
932    * false otherwise.
933    *
934    * For linear elements, performs an initial tight bounding box check
935    * (as an optimization step) and (if that passes) then uses the
936    * user-defined tolerance "tol" in a call to inverse_map() to actually
937    * test if the point is in the element.  For quadratic elements, the
938    * bounding box optimization is skipped, and only the inverse_map()
939    * steps are performed.
940    *
941    * \note This routine should not be used to determine if a point
942    * is merely "nearby" an element to within some tolerance. For that,
943    * use Elem::close_to_point() instead.
944    */
945   virtual bool contains_point (const Point & p, Real tol=TOLERANCE) const;
946 
947   /**
948    * \returns \p true if this element is "close" to the point p, where
949    * "close" is determined by the tolerance tol.
950    */
951   virtual bool close_to_point(const Point & p, Real tol) const;
952 
953 private:
954   /**
955    * Shared private implementation used by the contains_point()
956    * and close_to_point() routines.  The box_tol tolerance is
957    * used in the bounding box optimization, the map_tol tolerance is used
958    * in the calls to inverse_map() and on_reference_element().
959    */
960   bool point_test(const Point & p, Real box_tol, Real map_tol) const;
961 
962 public:
963   /**
964    * \returns \p true if the element map is definitely affine (i.e. the same at
965    * every quadrature point) within numerical tolerances.
966    */
has_affine_map()967   virtual bool has_affine_map () const { return false; }
968 
969   /**
970    * \returns \p true if the Lagrange shape functions on this element
971    * are linear.
972    */
is_linear()973   virtual bool is_linear () const { return false; }
974 
975   /**
976    * Prints relevant information about the element.
977    */
978   void print_info (std::ostream & os=libMesh::out) const;
979 
980   /**
981    * Prints relevant information about the element to a string.
982    */
983   std::string get_info () const;
984 
985   /**
986    * \returns \p true if the element is active (i.e. has no active
987    * descendants) or AMR is disabled, \p false otherwise.
988    *
989    * \note It suffices to check the first child only.
990    */
991   bool active () const;
992 
993   /**
994    * \returns \p true if the element is an ancestor (i.e. has an
995    * active child or ancestor child), \p false otherwise or when AMR
996    * is disabled.
997    */
998   bool ancestor () const;
999 
1000   /**
1001    * \returns \p true if the element is subactive (i.e. has no active
1002    * descendants), \p false otherwise or if AMR is disabled.
1003    */
1004   bool subactive () const;
1005 
1006   /**
1007    * \returns \p true if the element has any children (active or not),
1008    * \p false otherwise, or if AMR is disabled.
1009    */
1010   bool has_children () const;
1011 
1012   /**
1013    * \returns \p true if the element has any descendants other than
1014    * its immediate children, \p false otherwise, or if AMR is disabled.
1015    */
1016   bool has_ancestor_children () const;
1017 
1018   /**
1019    * \returns \p true if \p descendant is a child of \p this, or a
1020    * child of a child of \p this, etc., \p false otherwise or if AMR
1021    * is disabled.
1022    */
1023   bool is_ancestor_of(const Elem * descendant) const;
1024 
1025   /**
1026    * \returns A const pointer to the element's parent, or \p nullptr if
1027    * the element was not created via refinement.
1028    */
1029   const Elem * parent () const;
1030 
1031   /**
1032    * \returns A pointer to the element's parent, or \p nullptr if
1033    * the element was not created via refinement.
1034    */
1035   Elem * parent ();
1036 
1037   /**
1038    * Sets the pointer to the element's parent.
1039    * Dangerous! Only use this if you know what you are doing!
1040    */
1041   void set_parent (Elem * p);
1042 
1043   /**
1044    * \returns A pointer to the element's top-most (i.e. level-0) parent.
1045    *
1046    * That is, \p this if this is a level-0 element, this element's parent
1047    * if this is a level-1 element, this element's grandparent if this is
1048    * a level-2 element, etc...
1049    */
1050   const Elem * top_parent () const;
1051 
1052   /**
1053    * \returns The higher-dimensional Elem for which this Elem is a face.
1054    *
1055    * In some cases it is desirable to extract the boundary (or a subset thereof)
1056    * of a D-dimensional mesh as a (D-1)-dimensional manifold.  In this case
1057    * we may want to know the 'parent' element from which the manifold elements
1058    * were extracted.  We can easily do that for the level-0 manifold elements
1059    * by storing the D-dimensional parent.  This method provides access to that
1060    * element.
1061    *
1062    * This method returns nullptr if this->dim() == LIBMESH_DIM; in
1063    * such cases no data storage for an interior parent pointer has
1064    * been allocated.
1065    */
1066   const Elem * interior_parent () const;
1067 
1068   Elem * interior_parent ();
1069 
1070   /**
1071    * Sets the pointer to the element's interior_parent.
1072    * Dangerous! Only use this if you know what you are doing!
1073    */
1074   void set_interior_parent (Elem * p);
1075 
1076   /**
1077    * \returns The distance between nodes n1 and n2.
1078    *
1079    * Useful for computing the lengths of the sides of elements.
1080    */
1081   Real length (const unsigned int n1,
1082                const unsigned int n2) const;
1083 
1084   /**
1085    * \returns The number of adjacent vertices that uniquely define the
1086    * location of the \f$ n^{th} \f$ second-order node, or 0 for linear
1087    * elements.
1088    *
1089    * This method is useful when converting linear elements to quadratic
1090    * elements.
1091    *
1092    * \note \p n has to be greater than or equal to \p this->n_vertices().
1093    */
1094   virtual unsigned int n_second_order_adjacent_vertices (const unsigned int n) const;
1095 
1096   /**
1097    * \returns The element-local number of the \f$ v^{th} \f$ vertex
1098    * that defines the \f$ n^{th} \f$ second-order node, or 0 for
1099    * linear elements.
1100    *
1101    * \note The value is always less than \p this->n_vertices(), while
1102    * \p n has to be greater than or equal to \p this->n_vertices().
1103    */
1104   virtual unsigned short int second_order_adjacent_vertex (const unsigned int n,
1105                                                            const unsigned int v) const;
1106 
1107   /**
1108    * \returns A pair (c,v), where
1109    * c == child index, and
1110    * v == element-local index of the \p \f$ n^{th} \f$
1111    *      second-order node on the parent element.
1112    * For linear elements, (0,0) is returned.
1113    *
1114    * \note The return values are always less than \p this->n_children()
1115    * and \p this->child_ptr(c)->n_vertices().
1116    *
1117    * \note \p n has to be greater than or equal to \p this->n_vertices().
1118    *
1119    * \note On refined second-order elements, the return value will
1120    * satisfy \p this->node_ptr(n) == this->child_ptr(c)->node_ptr(v).
1121    */
1122   virtual std::pair<unsigned short int, unsigned short int>
1123   second_order_child_vertex (const unsigned int n) const;
1124 
1125   /**
1126    * \returns The element type of the associated second-order element,
1127    * or INVALID_ELEM for second-order or other elements that cannot be
1128    * converted into higher order equivalents.
1129    *
1130    * For example, when \p this is a \p TET4, then \p TET10 is returned.
1131    *
1132    * For some elements, there exist two second-order equivalents, e.g.
1133    * for \p Quad4 there is \p Quad8 and \p Quad9.  When the optional
1134    * \p full_ordered is \p true, then \p QUAD9 is returned.  When
1135    * \p full_ordered is \p false, then \p QUAD8 is returned.
1136    */
1137   static ElemType second_order_equivalent_type (const ElemType et,
1138                                                 const bool full_ordered=true);
1139 
1140   /**
1141    * \returns The element type of the associated first-order element,
1142    * or \p INVALID_ELEM for first-order or other elements that cannot be
1143    * converted into lower order equivalents.
1144    *
1145    * For example, when \p this is a \p TET10, then \p TET4 is returned.
1146    */
1147   static ElemType first_order_equivalent_type (const ElemType et);
1148 
1149 
1150   /**
1151    * \returns The refinement level of the current element.
1152    *
1153    * If the element's parent is \p nullptr then by convention it is at
1154    * level 0, otherwise it is simply at one level greater than its
1155    * parent.
1156    */
1157   unsigned int level () const;
1158 
1159   /**
1160    * \returns The value of the p refinement level of an active
1161    * element, or the minimum value of the p refinement levels
1162    * of an ancestor element's descendants.
1163    */
1164   unsigned int p_level () const;
1165 
1166   /**
1167    * \returns \p true if the specified child is on the specified side.
1168    */
1169   virtual bool is_child_on_side(const unsigned int c,
1170                                 const unsigned int s) const = 0;
1171 
1172   /**
1173    * \returns The value of the mapping type for the element.
1174    */
1175   ElemMappingType mapping_type () const;
1176 
1177   /**
1178    * Sets the value of the mapping type for the element.
1179    */
1180   void set_mapping_type (const ElemMappingType type);
1181 
1182   /**
1183    * \returns The value of the mapping data for the element.
1184    */
1185   unsigned char mapping_data () const;
1186 
1187   /**
1188    * Sets the value of the mapping data for the element.
1189    */
1190   void set_mapping_data (const unsigned char data);
1191 
1192 
1193 #ifdef LIBMESH_ENABLE_AMR
1194 
1195   /**
1196    * Enumeration of possible element refinement states.
1197    */
1198   enum RefinementState { COARSEN = 0,
1199                          DO_NOTHING,
1200                          REFINE,
1201                          JUST_REFINED,
1202                          JUST_COARSENED,
1203                          INACTIVE,
1204                          COARSEN_INACTIVE,
1205                          INVALID_REFINEMENTSTATE };
1206 
1207   /**
1208    * \returns A constant pointer to the \f$ i^{th} \f$ child for this element.
1209    * For internal use only - skips assertions about null pointers.
1210    */
1211   const Elem * raw_child_ptr (unsigned int i) const;
1212 
1213   /**
1214    * \returns A constant pointer to the \f$ i^{th} \f$ child for this element.
1215    * Do not call if this element has no children, i.e. is active.
1216    */
1217   const Elem * child_ptr (unsigned int i) const;
1218 
1219   /**
1220    * \returns A non-constant pointer to the \f$ i^{th} \f$ child for this element.
1221    * Do not call if this element has no children, i.e. is active.
1222    */
1223   Elem * child_ptr (unsigned int i);
1224 
1225   /**
1226    * Nested classes for use iterating over all children of a parent
1227    * element.
1228    */
1229   class ChildRefIter;
1230   class ConstChildRefIter;
1231 
1232   /**
1233    * Returns a range with all children of a parent element, usable in
1234    * range-based for loops.  The exact type of the return value here
1235    * may be subject to change in future libMesh releases, but the
1236    * iterators will always dereference to produce a reference to a
1237    * child element.
1238    */
1239   SimpleRange<ChildRefIter> child_ref_range();
1240 
1241   SimpleRange<ConstChildRefIter> child_ref_range() const;
1242 
1243 private:
1244   /**
1245    * Sets the pointer to the \f$ i^{th} \f$ child for this element.
1246    * Do not call if this element has no children, i.e. is active.
1247    */
1248   void set_child (unsigned int c, Elem * elem);
1249 
1250 public:
1251   /**
1252    * \returns The child index which \p e corresponds to.
1253    *
1254    * I.e. if c = a->which_child_am_i(e); then a->child_ptr(c) will be
1255    * e.
1256    */
1257   unsigned int which_child_am_i(const Elem * e) const;
1258 
1259   /**
1260    * \returns \p true if the specified child is on the specified edge.
1261    */
1262   virtual bool is_child_on_edge(const unsigned int c,
1263                                 const unsigned int e) const;
1264 
1265   /**
1266    * Adds a child pointer to the array of children of this element.
1267    * If this is the first child to be added, this method allocates
1268    * memory in the parent's _children array, otherwise, it just sets
1269    * the pointer.
1270    */
1271   void add_child (Elem * elem);
1272 
1273   /**
1274    * Adds a new child pointer to the specified index in the array of
1275    * children of this element.  If this is the first child to be added,
1276    * this method allocates memory in the parent's _children array,
1277    * otherwise, it just sets the pointer.
1278    */
1279   void add_child (Elem * elem, unsigned int c);
1280 
1281   /**
1282    * Replaces the child pointer at the specified index in the child array.
1283    */
1284   void replace_child (Elem * elem, unsigned int c);
1285 
1286   /**
1287    * Fills the vector \p family with the children of this element,
1288    * recursively.  Calling this method on a twice-refined element
1289    * will give you the element itself, its direct children, and their
1290    * children, etc...  When the optional parameter \p reset is
1291    * true, the vector will be cleared before the element and its
1292    * descendants are added.
1293    *
1294    * The family tree only includes ancestor and active elements. To
1295    * include subactive elements as well, use total_family_tree().
1296    */
1297   void family_tree (std::vector<const Elem *> & family,
1298                     bool reset = true) const;
1299 
1300   /**
1301    * Non-const version of function above; fills a vector of non-const pointers.
1302    */
1303   void family_tree (std::vector<Elem *> & family,
1304                     bool reset = true);
1305 
1306   /**
1307    * Same as the \p family_tree() member, but also adds any subactive
1308    * descendants.
1309    */
1310   void total_family_tree (std::vector<const Elem *> & family,
1311                           bool reset = true) const;
1312 
1313   /**
1314    * Non-const version of function above; fills a vector of non-const pointers.
1315    */
1316   void total_family_tree (std::vector<Elem *> & family,
1317                           bool reset = true);
1318 
1319   /**
1320    * Same as the \p family_tree() member, but only adds the active
1321    * children.  Can be thought of as removing all the inactive
1322    * elements from the vector created by \p family_tree, but is
1323    * implemented more efficiently.
1324    */
1325   void active_family_tree (std::vector<const Elem *> & active_family,
1326                            bool reset = true) const;
1327 
1328   /**
1329    * Non-const version of function above; fills a vector of non-const pointers.
1330    */
1331   void active_family_tree (std::vector<Elem *> & active_family,
1332                            bool reset = true);
1333 
1334   /**
1335    * Same as the \p family_tree() member, but only adds elements
1336    * which are next to \p side.
1337    */
1338   void family_tree_by_side (std::vector<const Elem *> & family,
1339                             unsigned int side,
1340                             bool reset = true) const;
1341 
1342   /**
1343    * Non-const version of function above; fills a vector of non-const pointers.
1344    */
1345   void family_tree_by_side (std::vector<Elem *> & family,
1346                             unsigned int side,
1347                             bool reset = true);
1348 
1349   /**
1350    * Same as the \p active_family_tree() member, but only adds elements
1351    * which are next to \p side.
1352    */
1353   void active_family_tree_by_side (std::vector<const Elem *> & family,
1354                                    unsigned int side,
1355                                    bool reset = true) const;
1356 
1357   /**
1358    * Non-const version of function above; fills a vector of non-const pointers.
1359    */
1360   void active_family_tree_by_side (std::vector<Elem *> & family,
1361                                    unsigned int side,
1362                                    bool reset = true);
1363 
1364   /**
1365    * Same as the \p family_tree() member, but only adds elements
1366    * which are next to \p neighbor.
1367    */
1368   void family_tree_by_neighbor (std::vector<const Elem *> & family,
1369                                 const Elem * neighbor,
1370                                 bool reset = true) const;
1371 
1372   /**
1373    * Non-const version of function above; fills a vector of non-const pointers.
1374    */
1375   void family_tree_by_neighbor (std::vector<Elem *> & family,
1376                                 Elem * neighbor,
1377                                 bool reset = true);
1378 
1379   /**
1380    * Same as the \p family_tree_by_neighbor() member, but also adds
1381    * any subactive descendants.
1382    */
1383   void total_family_tree_by_neighbor (std::vector<const Elem *> & family,
1384                                       const Elem * neighbor,
1385                                       bool reset = true) const;
1386 
1387   /**
1388    * Non-const version of function above; fills a vector of non-const pointers.
1389    */
1390   void total_family_tree_by_neighbor (std::vector<Elem *> & family,
1391                                       Elem * neighbor,
1392                                       bool reset = true);
1393 
1394   /**
1395    * Same as the \p family_tree() member, but only adds elements
1396    * which are next to \p subneighbor.  Only applicable when
1397    * \p this->has_neighbor(neighbor) and
1398    * \p neighbor->is_ancestor(subneighbor)
1399    */
1400   void family_tree_by_subneighbor (std::vector<const Elem *> & family,
1401                                    const Elem * neighbor,
1402                                    const Elem * subneighbor,
1403                                    bool reset = true) const;
1404 
1405   /**
1406    * Non-const version of function above; fills a vector of non-const pointers.
1407    */
1408   void family_tree_by_subneighbor (std::vector<Elem *> & family,
1409                                    Elem * neighbor,
1410                                    Elem * subneighbor,
1411                                    bool reset = true);
1412 
1413   /**
1414    * Same as the \p family_tree_by_subneighbor() member, but also adds
1415    * any subactive descendants.
1416    */
1417   void total_family_tree_by_subneighbor (std::vector<const Elem *> & family,
1418                                          const Elem * neighbor,
1419                                          const Elem * subneighbor,
1420                                          bool reset = true) const;
1421 
1422   /**
1423    * Non-const version of function above; fills a vector of non-const pointers.
1424    */
1425   void total_family_tree_by_subneighbor (std::vector<Elem *> & family,
1426                                          Elem * neighbor,
1427                                          Elem * subneighbor,
1428                                          bool reset = true);
1429 
1430   /**
1431    * Same as the \p active_family_tree() member, but only adds elements
1432    * which are next to \p neighbor.
1433    */
1434   void active_family_tree_by_neighbor (std::vector<const Elem *> & family,
1435                                        const Elem * neighbor,
1436                                        bool reset = true) const;
1437 
1438   /**
1439    * Non-const version of function above; fills a vector of non-const pointers.
1440    */
1441   void active_family_tree_by_neighbor (std::vector<Elem *> & family,
1442                                        Elem * neighbor,
1443                                        bool reset = true);
1444 
1445   /**
1446    * Same as the \p active_family_tree_by_neighbor() member, but the
1447    * \p neighbor here may be a topological (e.g. periodic boundary
1448    * condition) neighbor, not just a local neighbor.
1449    */
1450   void active_family_tree_by_topological_neighbor (std::vector<const Elem *> & family,
1451                                                    const Elem * neighbor,
1452                                                    const MeshBase & mesh,
1453                                                    const PointLocatorBase & point_locator,
1454                                                    const PeriodicBoundaries * pb,
1455                                                    bool reset = true) const;
1456 
1457   /**
1458    * Non-const version of function above; fills a vector of non-const pointers.
1459    */
1460   void active_family_tree_by_topological_neighbor (std::vector<Elem *> & family,
1461                                                    Elem * neighbor,
1462                                                    const MeshBase & mesh,
1463                                                    const PointLocatorBase & point_locator,
1464                                                    const PeriodicBoundaries * pb,
1465                                                    bool reset = true);
1466 
1467   /**
1468    * \returns The value of the refinement flag for the element.
1469    */
1470   RefinementState refinement_flag () const;
1471 
1472   /**
1473    * Sets the value of the refinement flag for the element.
1474    */
1475   void set_refinement_flag (const RefinementState rflag);
1476 
1477   /**
1478    * \returns The value of the p-refinement flag for the element.
1479    */
1480   RefinementState p_refinement_flag () const;
1481 
1482   /**
1483    * Sets the value of the p-refinement flag for the element.
1484    */
1485   void set_p_refinement_flag (const RefinementState pflag);
1486 
1487   /**
1488    * \returns The maximum value of the p-refinement levels of
1489    * an ancestor element's descendants.
1490    */
1491   unsigned int max_descendant_p_level () const;
1492 
1493   /**
1494    * \returns The minimum p-refinement level of elements which are
1495    * descended from this element, and which share a side with the
1496    * active \p neighbor.
1497    */
1498   unsigned int min_p_level_by_neighbor (const Elem * neighbor,
1499                                         unsigned int current_min) const;
1500 
1501   /**
1502    * \returns The minimum new p-refinement level (i.e. after refinement
1503    * and coarsening is done) of elements which are descended from this
1504    * element and which share a side with the active \p neighbor.
1505    */
1506   unsigned int min_new_p_level_by_neighbor (const Elem * neighbor,
1507                                             unsigned int current_min) const;
1508 
1509   /**
1510    * Sets the value of the p-refinement level for the element.
1511    *
1512    * \note The maximum p-refinement level is currently 255.
1513    */
1514   void set_p_level (const unsigned int p);
1515 
1516   /**
1517    * Sets the value of the p-refinement level for the element
1518    * without altering the p-level of its ancestors
1519    */
1520   void hack_p_level (const unsigned int p);
1521 
1522   /**
1523    * Refine the element.
1524    */
1525   virtual void refine (MeshRefinement & mesh_refinement);
1526 
1527   /**
1528    * Coarsen the element.  This function is non-virtual since it is the same
1529    * for all element types.
1530    */
1531   void coarsen ();
1532 
1533   /**
1534    * Contract an active element, i.e. remove pointers to any
1535    * subactive children.  This should only be called via
1536    * MeshRefinement::contract, which will also remove subactive
1537    * children from the mesh.
1538    */
1539   void contract ();
1540 
1541 #endif
1542 
1543 #ifdef DEBUG
1544   /**
1545    * Checks for consistent neighbor links on this element.
1546    */
1547   void libmesh_assert_valid_neighbors() const;
1548 
1549   /**
1550    * Checks for a valid id and pointers to nodes with valid ids on
1551    * this element.
1552    */
1553   void libmesh_assert_valid_node_pointers() const;
1554 #endif // DEBUG
1555 
1556 protected:
1557 
1558   /**
1559    * The protected nested SideIter class is used to iterate over the
1560    * sides of this Elem.  It is a specially-designed class since
1561    * no sides are actually stored by the element.  This iterator-like
1562    * class has to provide the following three operations
1563    * 1) operator*
1564    * 2) operator++
1565    * 3) operator==
1566    * The definition can be found at the end of this header file.
1567    */
1568   class SideIter;
1569 
1570 public:
1571   /**
1572    * Useful iterator typedefs
1573    */
1574   typedef Predicates::multi_predicate Predicate;
1575 
1576   /**
1577    * Data structure for iterating over sides.  Defined at the end of
1578    * this header file.
1579    */
1580   struct side_iterator;
1581 
1582   /**
1583    * Iterator accessor functions
1584    */
1585   side_iterator boundary_sides_begin();
1586   side_iterator boundary_sides_end();
1587 
1588 private:
1589   /**
1590    * Side iterator helper functions.  Used to replace the begin()
1591    * and end() functions of the STL containers.
1592    */
1593   SideIter _first_side();
1594   SideIter _last_side();
1595 
1596 public:
1597 
1598 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1599 
1600   /**
1601    * \returns \p true if the element is an infinite element,
1602    * \p false otherwise.
1603    */
1604   virtual bool infinite () const = 0;
1605 
1606   /**
1607    * \returns \p true if the specified (local) node number is a
1608    * "mid-edge" node on an infinite element edge.
1609    *
1610    * This is false for all nodes on non-infinite elements, so we won't
1611    * make it pure virtual, to simplify their code.
1612    */
is_mid_infinite_edge_node(const unsigned int)1613   virtual bool is_mid_infinite_edge_node(const unsigned int /* n */) const
1614   { libmesh_assert (!this->infinite()); return false; }
1615 
1616   /**
1617    * \returns The origin for an infinite element.
1618    *
1619    * Currently, all infinite elements used in a mesh share the same
1620    * origin.  Override this in infinite element classes.
1621    */
origin()1622   virtual Point origin () const { libmesh_not_implemented(); return Point(); }
1623 
1624 #endif
1625 
1626 
1627 
1628 
1629   /**
1630    * \returns An Elem of type \p type wrapped in a smart pointer.
1631    */
1632   static std::unique_ptr<Elem> build (const ElemType type,
1633                                       Elem * p=nullptr);
1634   /**
1635    * Calls the build() method above with a nullptr parent, and
1636    * additionally sets the newly-created Elem's id. This can be useful
1637    * when adding pre-numbered Elems to a Mesh via add_elem() calls.
1638    */
1639   static std::unique_ptr<Elem> build_with_id (const ElemType type,
1640                                               dof_id_type id);
1641 
1642 #ifdef LIBMESH_ENABLE_AMR
1643 
1644   /**
1645    * \returns The local node id on the parent which corresponds to node
1646    * \p n of child \p c, or \p invalid_uint if no such parent
1647    * node exists.
1648    */
1649   virtual unsigned int as_parent_node (unsigned int c,
1650                                        unsigned int n) const;
1651 
1652   /**
1653    * \returns All the pairs of nodes (indexed by local node id) which
1654    * should bracket node \p n of child \p c.
1655    */
1656   virtual
1657   const std::vector<std::pair<unsigned char, unsigned char>> &
1658   parent_bracketing_nodes(unsigned int c,
1659                           unsigned int n) const;
1660 
1661   /**
1662    * \returns All the pairs of nodes (indexed by global node id) which
1663    * should bracket node \p n of child \p c.
1664    */
1665   virtual
1666   const std::vector<std::pair<dof_id_type, dof_id_type>>
1667   bracketing_nodes(unsigned int c,
1668                    unsigned int n) const;
1669 
1670 
1671   /**
1672    * \returns The embedding matrix entry for the requested child.
1673    */
1674   virtual float embedding_matrix (const unsigned int child_num,
1675                                   const unsigned int child_node_num,
1676                                   const unsigned int parent_node_num) const = 0;
1677 
1678   /**
1679    * \returns A "version number" that identifies which embedding
1680    * matrix is in use.
1681    *
1682    * Some element types may use a different embedding matrix depending
1683    * on their geometric characteristics.
1684    */
embedding_matrix_version()1685   virtual unsigned int embedding_matrix_version () const { return 0; }
1686 
1687 #endif // LIBMESH_ENABLE_AMR
1688 
1689 
1690 protected:
1691 
1692   /**
1693    * \returns A hash key computed from a single node id.
1694    */
1695   static dof_id_type compute_key (dof_id_type n0);
1696 
1697   /**
1698    * \returns A hash key computed from two node ids.
1699    */
1700   static dof_id_type compute_key (dof_id_type n0,
1701                                   dof_id_type n1);
1702 
1703   /**
1704    * \returns A hash key computed from three node ids.
1705    */
1706   static dof_id_type compute_key (dof_id_type n0,
1707                                   dof_id_type n1,
1708                                   dof_id_type n2);
1709 
1710   /**
1711    * \returns A hash key computed from four node ids.
1712    */
1713   static dof_id_type compute_key (dof_id_type n0,
1714                                   dof_id_type n1,
1715                                   dof_id_type n2,
1716                                   dof_id_type n3);
1717 
1718   /**
1719    * An implementation for simple (all sides equal) elements
1720    */
1721   template <typename Sideclass, typename Subclass>
1722   std::unique_ptr<Elem>
1723   simple_build_side_ptr(const unsigned int i,
1724                         bool proxy);
1725 
1726   /**
1727    * An implementation for simple (all sides equal) elements
1728    */
1729   template <typename Subclass>
1730   void simple_build_side_ptr(std::unique_ptr<Elem> & side,
1731                              const unsigned int i,
1732                              ElemType sidetype);
1733 
1734   /**
1735    * An implementation for simple (all sides equal) elements
1736    */
1737   template <typename Subclass, typename Mapclass>
1738   void simple_side_ptr(std::unique_ptr<Elem> & side,
1739                        const unsigned int i,
1740                        ElemType sidetype);
1741 
1742 #ifdef LIBMESH_ENABLE_AMR
1743 
1744   /**
1745    * Elem subclasses which don't do their own bracketing node
1746    * calculations will need to supply a static cache, since the
1747    * default calculation is slow.
1748    */
1749   virtual
1750   std::vector<std::vector<std::vector<std::vector<std::pair<unsigned char, unsigned char>>>>> &
_get_bracketing_node_cache()1751   _get_bracketing_node_cache() const
1752   {
1753     static std::vector<std::vector<std::vector<std::vector<std::pair<unsigned char, unsigned char>>>>> c;
1754     libmesh_error();
1755     return c;
1756   }
1757 
1758   /**
1759    * Elem subclasses which don't do their own child-to-parent node
1760    * calculations will need to supply a static cache, since the
1761    * default calculation is slow.
1762    */
1763   virtual
1764   std::vector<std::vector<std::vector<signed char>>> &
_get_parent_indices_cache()1765   _get_parent_indices_cache() const
1766   {
1767     static std::vector<std::vector<std::vector<signed char>>> c;
1768     libmesh_error();
1769     return c;
1770   }
1771 
1772 #endif // LIBMESH_ENABLE_AMR
1773 
1774 public:
1775 
1776   /**
1777    * Replaces this element with \p nullptr for all of its neighbors.
1778    * This is useful when deleting an element.
1779    */
1780   void nullify_neighbors ();
1781 
1782 protected:
1783 
1784   /**
1785    * Pointers to the nodes we are connected to.
1786    */
1787   Node ** _nodes;
1788 
1789   /**
1790    * Pointers to this element's parent and neighbors, and for
1791    * lower-dimensional elements' interior_parent.
1792    */
1793   Elem ** _elemlinks;
1794 
1795 #ifdef LIBMESH_ENABLE_AMR
1796   /**
1797    * Pointers to this element's children.
1798    */
1799   Elem ** _children;
1800 #endif
1801 
1802   /**
1803    * The subdomain to which this element belongs.
1804    */
1805   subdomain_id_type _sbd_id;
1806 
1807 #ifdef LIBMESH_ENABLE_AMR
1808   /**
1809    * h refinement flag. This is stored as an unsigned char
1810    * to save space.
1811    */
1812   unsigned char _rflag;
1813 
1814   /**
1815    * p refinement flag. This is stored as an unsigned char
1816    * to save space.
1817    */
1818   unsigned char _pflag;
1819 
1820   /**
1821    * p refinement level - the difference between the
1822    * polynomial degree on this element and the minimum
1823    * polynomial degree on the mesh.
1824    * This is stored as an unsigned char to save space.
1825    * In theory, these last four bytes might have
1826    * been padding anyway.
1827    */
1828   unsigned char _p_level;
1829 #endif
1830 
1831   /**
1832    * Mapping function type; currently either 0 (LAGRANGE) or 1
1833    * (RATIONAL_BERNSTEIN).
1834    */
1835   unsigned char _map_type;
1836 
1837   /**
1838    * Mapping function data; currently used when needed to store the
1839    * RATIONAL_BERNSTEIN nodal weight data index.
1840    */
1841   unsigned char _map_data;
1842 };
1843 
1844 
1845 
1846 // ------------------------------------------------------------
1847 // Elem helper classes
1848 //
1849 class
1850 Elem::NodeRefIter : public PointerToPointerIter<Node>
1851 {
1852 public:
NodeRefIter(Node * const * nodepp)1853   NodeRefIter (Node * const * nodepp) : PointerToPointerIter<Node>(nodepp) {}
1854 };
1855 
1856 
1857 class
1858 Elem::ConstNodeRefIter : public PointerToPointerIter<const Node>
1859 {
1860 public:
ConstNodeRefIter(const Node * const * nodepp)1861   ConstNodeRefIter (const Node * const * nodepp) : PointerToPointerIter<const Node>(nodepp) {}
1862 };
1863 
1864 
1865 #ifdef LIBMESH_ENABLE_AMR
1866 class
1867 Elem::ChildRefIter : public PointerToPointerIter<Elem>
1868 {
1869 public:
ChildRefIter(Elem * const * childpp)1870   ChildRefIter (Elem * const * childpp) : PointerToPointerIter<Elem>(childpp) {}
1871 };
1872 
1873 
1874 class
1875 Elem::ConstChildRefIter : public PointerToPointerIter<const Elem>
1876 {
1877 public:
ConstChildRefIter(const Elem * const * childpp)1878   ConstChildRefIter (const Elem * const * childpp) : PointerToPointerIter<const Elem>(childpp) {}
1879 };
1880 
1881 
1882 inline
child_ref_range()1883 SimpleRange<Elem::ChildRefIter> Elem::child_ref_range()
1884 {
1885   libmesh_assert(_children);
1886   return {_children, _children + this->n_children()};
1887 }
1888 
1889 
1890 inline
child_ref_range()1891 SimpleRange<Elem::ConstChildRefIter> Elem::child_ref_range() const
1892 {
1893   libmesh_assert(_children);
1894   return {_children, _children + this->n_children()};
1895 }
1896 #endif // LIBMESH_ENABLE_AMR
1897 
1898 
1899 
1900 
1901 // ------------------------------------------------------------
1902 // global Elem functions
1903 
1904 inline
1905 std::ostream & operator << (std::ostream & os, const Elem & e)
1906 {
1907   e.print_info(os);
1908   return os;
1909 }
1910 
1911 
1912 // ------------------------------------------------------------
1913 // Elem class member functions
1914 inline
Elem(const unsigned int nn,const unsigned int ns,Elem * p,Elem ** elemlinkdata,Node ** nodelinkdata)1915 Elem::Elem(const unsigned int nn,
1916            const unsigned int ns,
1917            Elem * p,
1918            Elem ** elemlinkdata,
1919            Node ** nodelinkdata) :
1920   _nodes(nodelinkdata),
1921   _elemlinks(elemlinkdata),
1922 #ifdef LIBMESH_ENABLE_AMR
1923   _children(nullptr),
1924 #endif
1925   _sbd_id(0),
1926 #ifdef LIBMESH_ENABLE_AMR
1927   _rflag(Elem::DO_NOTHING),
1928   _pflag(Elem::DO_NOTHING),
1929   _p_level(0),
1930 #endif
1931   _map_type(p ? p->mapping_type() : 0),
1932   _map_data(p ? p->mapping_data() : 0)
1933 {
1934   this->processor_id() = DofObject::invalid_processor_id;
1935 
1936   // If this ever legitimately fails we need to increase max_n_nodes
1937   libmesh_assert_less_equal(nn, max_n_nodes);
1938 
1939   // Initialize the nodes data structure
1940   if (_nodes)
1941     {
1942       for (unsigned int n=0; n<nn; n++)
1943         _nodes[n] = nullptr;
1944     }
1945 
1946   // Initialize the neighbors/parent data structure
1947   // _elemlinks = new Elem *[ns+1];
1948 
1949   if (_elemlinks)
1950     {
1951       _elemlinks[0] = p;
1952 
1953       for (unsigned int n=1; n<ns+1; n++)
1954         _elemlinks[n] = nullptr;
1955     }
1956 
1957   // Optionally initialize data from the parent
1958   if (this->parent() != nullptr)
1959     {
1960       this->subdomain_id() = this->parent()->subdomain_id();
1961       this->processor_id() = this->parent()->processor_id();
1962       _map_type = this->parent()->_map_type;
1963       _map_data = this->parent()->_map_data;
1964     }
1965 
1966 #ifdef LIBMESH_ENABLE_AMR
1967   if (this->parent())
1968     this->set_p_level(this->parent()->p_level());
1969 #endif
1970 }
1971 
1972 
1973 
1974 inline
~Elem()1975 Elem::~Elem()
1976 {
1977   // Deleting my parent/neighbor/nodes storage isn't necessary since it's
1978   // handled by the subclass
1979 
1980   // if (_nodes != nullptr)
1981   //   delete [] _nodes;
1982   // _nodes = nullptr;
1983 
1984   // delete [] _elemlinks;
1985 
1986 #ifdef LIBMESH_ENABLE_AMR
1987 
1988   // Delete my children's storage
1989   if (_children != nullptr)
1990     delete [] _children;
1991   _children = nullptr;
1992 
1993 #endif
1994 }
1995 
1996 
1997 
1998 inline
point(const unsigned int i)1999 const Point & Elem::point (const unsigned int i) const
2000 {
2001   libmesh_assert_less (i, this->n_nodes());
2002   libmesh_assert(_nodes[i]);
2003   libmesh_assert_not_equal_to (_nodes[i]->id(), Node::invalid_id);
2004 
2005   return *_nodes[i];
2006 }
2007 
2008 
2009 
2010 inline
point(const unsigned int i)2011 Point & Elem::point (const unsigned int i)
2012 {
2013   libmesh_assert_less (i, this->n_nodes());
2014 
2015   return *_nodes[i];
2016 }
2017 
2018 
2019 
2020 inline
node_id(const unsigned int i)2021 dof_id_type Elem::node_id (const unsigned int i) const
2022 {
2023   libmesh_assert_less (i, this->n_nodes());
2024   libmesh_assert(_nodes[i]);
2025   libmesh_assert_not_equal_to (_nodes[i]->id(), Node::invalid_id);
2026 
2027   return _nodes[i]->id();
2028 }
2029 
2030 
2031 
2032 inline
local_node(const dof_id_type i)2033 unsigned int Elem::local_node (const dof_id_type i) const
2034 {
2035   for (auto n : make_range(this->n_nodes()))
2036     if (this->node_id(n) == i)
2037       return n;
2038 
2039   return libMesh::invalid_uint;
2040 }
2041 
2042 
2043 
2044 inline
get_nodes()2045 const Node * const * Elem::get_nodes () const
2046 {
2047   return _nodes;
2048 }
2049 
2050 
2051 
2052 inline
node_ptr(const unsigned int i)2053 const Node * Elem::node_ptr (const unsigned int i) const
2054 {
2055   libmesh_assert_less (i, this->n_nodes());
2056   libmesh_assert(_nodes[i]);
2057 
2058   return _nodes[i];
2059 }
2060 
2061 
2062 
2063 inline
node_ptr(const unsigned int i)2064 Node * Elem::node_ptr (const unsigned int i)
2065 {
2066   libmesh_assert_less (i, this->n_nodes());
2067   libmesh_assert(_nodes[i]);
2068 
2069   return _nodes[i];
2070 }
2071 
2072 
2073 
2074 inline
node_ref(const unsigned int i)2075 const Node & Elem::node_ref (const unsigned int i) const
2076 {
2077   return *this->node_ptr(i);
2078 }
2079 
2080 
2081 
2082 inline
node_ref(const unsigned int i)2083 Node & Elem::node_ref (const unsigned int i)
2084 {
2085   return *this->node_ptr(i);
2086 }
2087 
2088 
2089 
2090 inline
get_node_index(const Node * node_ptr)2091 unsigned int Elem::get_node_index (const Node * node_ptr) const
2092 {
2093   for (auto n : make_range(this->n_nodes()))
2094     if (this->_nodes[n] == node_ptr)
2095       return n;
2096 
2097   return libMesh::invalid_uint;
2098 }
2099 
2100 
2101 
2102 inline
set_node(const unsigned int i)2103 Node * & Elem::set_node (const unsigned int i)
2104 {
2105   libmesh_assert_less (i, this->n_nodes());
2106 
2107   return _nodes[i];
2108 }
2109 
2110 
2111 
2112 inline
subdomain_id()2113 subdomain_id_type Elem::subdomain_id () const
2114 {
2115   return _sbd_id;
2116 }
2117 
2118 
2119 
2120 inline
subdomain_id()2121 subdomain_id_type & Elem::subdomain_id ()
2122 {
2123   return _sbd_id;
2124 }
2125 
2126 
2127 
2128 inline
neighbor_ptr(unsigned int i)2129 const Elem * Elem::neighbor_ptr (unsigned int i) const
2130 {
2131   libmesh_assert_less (i, this->n_neighbors());
2132 
2133   return _elemlinks[i+1];
2134 }
2135 
2136 
2137 
2138 inline
neighbor_ptr(unsigned int i)2139 Elem * Elem::neighbor_ptr (unsigned int i)
2140 {
2141   libmesh_assert_less (i, this->n_neighbors());
2142 
2143   return _elemlinks[i+1];
2144 }
2145 
2146 
2147 
2148 inline
set_neighbor(const unsigned int i,Elem * n)2149 void Elem::set_neighbor (const unsigned int i, Elem * n)
2150 {
2151   libmesh_assert_less (i, this->n_neighbors());
2152 
2153   _elemlinks[i+1] = n;
2154 }
2155 
2156 
2157 
2158 inline
has_neighbor(const Elem * elem)2159 bool Elem::has_neighbor (const Elem * elem) const
2160 {
2161   for (auto n : this->neighbor_ptr_range())
2162     if (n == elem)
2163       return true;
2164 
2165   return false;
2166 }
2167 
2168 
2169 
2170 inline
child_neighbor(Elem * elem)2171 Elem * Elem::child_neighbor (Elem * elem)
2172 {
2173   for (auto n : elem->neighbor_ptr_range())
2174     if (n && n->parent() == this)
2175       return n;
2176 
2177   return nullptr;
2178 }
2179 
2180 
2181 
2182 inline
child_neighbor(const Elem * elem)2183 const Elem * Elem::child_neighbor (const Elem * elem) const
2184 {
2185   for (auto n : elem->neighbor_ptr_range())
2186     if (n && n->parent() == this)
2187       return n;
2188 
2189   return nullptr;
2190 }
2191 
2192 
2193 
2194 inline
2195 SimpleRange<Elem::NodeRefIter>
node_ref_range()2196 Elem::node_ref_range()
2197 {
2198   return {_nodes, _nodes+this->n_nodes()};
2199 }
2200 
2201 
2202 
2203 inline
2204 SimpleRange<Elem::ConstNodeRefIter>
node_ref_range()2205 Elem::node_ref_range() const
2206 {
2207   return {_nodes, _nodes+this->n_nodes()};
2208 }
2209 
2210 
2211 
2212 inline
2213 IntRange<unsigned short>
node_index_range()2214 Elem::node_index_range() const
2215 {
2216   return {0, cast_int<unsigned short>(this->n_nodes())};
2217 }
2218 
2219 
2220 
2221 inline
2222 IntRange<unsigned short>
edge_index_range()2223 Elem::edge_index_range() const
2224 {
2225   return {0, cast_int<unsigned short>(this->n_edges())};
2226 }
2227 
2228 
2229 
2230 inline
2231 IntRange<unsigned short>
side_index_range()2232 Elem::side_index_range() const
2233 {
2234   return {0, cast_int<unsigned short>(this->n_sides())};
2235 }
2236 
2237 
2238 
2239 
2240 inline
side_ptr(unsigned int i)2241 std::unique_ptr<const Elem> Elem::side_ptr (unsigned int i) const
2242 {
2243   // Call the non-const version of this function, return the result as
2244   // a std::unique_ptr<const Elem>.
2245   Elem * me = const_cast<Elem *>(this);
2246   const Elem * s = const_cast<const Elem *>(me->side_ptr(i).release());
2247   return std::unique_ptr<const Elem>(s);
2248 }
2249 
2250 
2251 
2252 inline
2253 void
side_ptr(std::unique_ptr<const Elem> & elem,const unsigned int i)2254 Elem::side_ptr (std::unique_ptr<const Elem> & elem,
2255                 const unsigned int i) const
2256 {
2257   // Hand off to the non-const version of this function
2258   Elem * me = const_cast<Elem *>(this);
2259   std::unique_ptr<Elem> e {const_cast<Elem *>(elem.release())};
2260   me->side_ptr(e, i);
2261   elem.reset(e.release());
2262 }
2263 
2264 
2265 
2266 inline
2267 std::unique_ptr<const Elem>
build_side_ptr(const unsigned int i,bool proxy)2268 Elem::build_side_ptr (const unsigned int i, bool proxy) const
2269 {
2270   // Call the non-const version of this function, return the result as
2271   // a std::unique_ptr<const Elem>.
2272   Elem * me = const_cast<Elem *>(this);
2273   const Elem * s = const_cast<const Elem *>(me->build_side_ptr(i, proxy).release());
2274   return std::unique_ptr<const Elem>(s);
2275 }
2276 
2277 
2278 
2279 inline
2280 void
build_side_ptr(std::unique_ptr<const Elem> & elem,const unsigned int i)2281 Elem::build_side_ptr (std::unique_ptr<const Elem> & elem,
2282                       const unsigned int i) const
2283 {
2284   // Hand off to the non-const version of this function
2285   Elem * me = const_cast<Elem *>(this);
2286   std::unique_ptr<Elem> e {const_cast<Elem *>(elem.release())};
2287   me->build_side_ptr(e, i);
2288   elem.reset(e.release());
2289 }
2290 
2291 
2292 
2293 template <typename Sideclass, typename Subclass>
2294 inline
2295 std::unique_ptr<Elem>
simple_build_side_ptr(const unsigned int i,bool proxy)2296 Elem::simple_build_side_ptr (const unsigned int i,
2297                              bool proxy)
2298 {
2299   libmesh_assert_less (i, this->n_sides());
2300 
2301   std::unique_ptr<Elem> face;
2302   if (proxy)
2303     face = libmesh_make_unique<Side<Sideclass,Subclass>>(this,i);
2304   else
2305     {
2306       face = libmesh_make_unique<Sideclass>(this);
2307       for (auto n : face->node_index_range())
2308         face->set_node(n) = this->node_ptr(Subclass::side_nodes_map[i][n]);
2309     }
2310 
2311 #ifdef LIBMESH_ENABLE_DEPRECATED
2312   if (!proxy) // proxy sides used to leave parent() set
2313 #endif
2314     face->set_parent(nullptr);
2315   face->set_interior_parent(this);
2316 
2317   return face;
2318 }
2319 
2320 
2321 
2322 template <typename Subclass>
2323 inline
2324 void
simple_build_side_ptr(std::unique_ptr<Elem> & side,const unsigned int i,ElemType sidetype)2325 Elem::simple_build_side_ptr (std::unique_ptr<Elem> & side,
2326                              const unsigned int i,
2327                              ElemType sidetype)
2328 {
2329   libmesh_assert_less (i, this->n_sides());
2330 
2331   if (!side.get() || side->type() != sidetype)
2332     {
2333       Subclass & real_me = cast_ref<Subclass&>(*this);
2334       side = real_me.Subclass::build_side_ptr(i, false);
2335     }
2336   else
2337     {
2338       side->subdomain_id() = this->subdomain_id();
2339 #ifdef LIBMESH_ENABLE_AMR
2340       side->set_p_level(this->p_level());
2341 #endif
2342       for (auto n : side->node_index_range())
2343         side->set_node(n) = this->node_ptr(Subclass::side_nodes_map[i][n]);
2344     }
2345 }
2346 
2347 
2348 
2349 template <typename Subclass, typename Mapclass>
2350 inline
2351 void
simple_side_ptr(std::unique_ptr<Elem> & side,const unsigned int i,ElemType sidetype)2352 Elem::simple_side_ptr (std::unique_ptr<Elem> & side,
2353                        const unsigned int i,
2354                        ElemType sidetype)
2355 {
2356   libmesh_assert_less (i, this->n_sides());
2357 
2358   if (!side.get() || side->type() != sidetype)
2359     {
2360       Subclass & real_me = cast_ref<Subclass&>(*this);
2361       side = real_me.Subclass::side_ptr(i);
2362     }
2363   else
2364     {
2365       side->subdomain_id() = this->subdomain_id();
2366 
2367       for (auto n : side->node_index_range())
2368         side->set_node(n) = this->node_ptr(Mapclass::side_nodes_map[i][n]);
2369     }
2370 }
2371 
2372 
2373 
2374 inline
2375 std::unique_ptr<const Elem>
build_edge_ptr(const unsigned int i)2376 Elem::build_edge_ptr (const unsigned int i) const
2377 {
2378   // Call the non-const version of this function, return the result as
2379   // a std::unique_ptr<const Elem>.
2380   Elem * me = const_cast<Elem *>(this);
2381   const Elem * e = const_cast<const Elem *>(me->build_edge_ptr(i).release());
2382   return std::unique_ptr<const Elem>(e);
2383 }
2384 
2385 
2386 
2387 inline
on_boundary()2388 bool Elem::on_boundary () const
2389 {
2390   // By convention, the element is on the boundary
2391   // if it has a nullptr neighbor.
2392   return this->has_neighbor(nullptr);
2393 }
2394 
2395 
2396 
2397 inline
which_neighbor_am_i(const Elem * e)2398 unsigned int Elem::which_neighbor_am_i (const Elem * e) const
2399 {
2400   libmesh_assert(e);
2401 
2402   const Elem * eparent = e;
2403 
2404   while (eparent->level() > this->level())
2405     {
2406       eparent = eparent->parent();
2407       libmesh_assert(eparent);
2408     }
2409 
2410   for (auto s : make_range(this->n_sides()))
2411     if (this->neighbor_ptr(s) == eparent)
2412       return s;
2413 
2414   return libMesh::invalid_uint;
2415 }
2416 
2417 
2418 
2419 inline
active()2420 bool Elem::active() const
2421 {
2422 #ifdef LIBMESH_ENABLE_AMR
2423   if ((this->refinement_flag() == INACTIVE) ||
2424       (this->refinement_flag() == COARSEN_INACTIVE))
2425     return false;
2426   else
2427     return true;
2428 #else
2429   return true;
2430 #endif
2431 }
2432 
2433 
2434 
2435 
2436 
2437 inline
subactive()2438 bool Elem::subactive() const
2439 {
2440 #ifdef LIBMESH_ENABLE_AMR
2441   if (this->active())
2442     return false;
2443   if (!this->has_children())
2444     return true;
2445   for (const Elem * my_ancestor = this->parent();
2446        my_ancestor != nullptr;
2447        my_ancestor = my_ancestor->parent())
2448     if (my_ancestor->active())
2449       return true;
2450 #endif
2451 
2452   return false;
2453 }
2454 
2455 
2456 
2457 inline
has_children()2458 bool Elem::has_children() const
2459 {
2460 #ifdef LIBMESH_ENABLE_AMR
2461   if (_children == nullptr)
2462     return false;
2463   else
2464     return true;
2465 #else
2466   return false;
2467 #endif
2468 }
2469 
2470 
2471 inline
has_ancestor_children()2472 bool Elem::has_ancestor_children() const
2473 {
2474 #ifdef LIBMESH_ENABLE_AMR
2475   if (_children == nullptr)
2476     return false;
2477   else
2478     for (auto & c : child_ref_range())
2479       if (c.has_children())
2480         return true;
2481 #endif
2482   return false;
2483 }
2484 
2485 
2486 
2487 inline
is_ancestor_of(const Elem * descendant)2488 bool Elem::is_ancestor_of(const Elem *
2489 #ifdef LIBMESH_ENABLE_AMR
2490                           descendant
2491 #endif
2492                           ) const
2493 {
2494 #ifdef LIBMESH_ENABLE_AMR
2495   const Elem * e = descendant;
2496   while (e)
2497     {
2498       if (this == e)
2499         return true;
2500       e = e->parent();
2501     }
2502 #endif
2503   return false;
2504 }
2505 
2506 
2507 
2508 inline
parent()2509 const Elem * Elem::parent () const
2510 {
2511   return _elemlinks[0];
2512 }
2513 
2514 
2515 
2516 inline
parent()2517 Elem * Elem::parent ()
2518 {
2519   return _elemlinks[0];
2520 }
2521 
2522 
2523 
2524 inline
set_parent(Elem * p)2525 void Elem::set_parent (Elem * p)
2526 {
2527   // We no longer support using parent() as interior_parent()
2528   libmesh_assert_equal_to(this->dim(), p ? p->dim() : this->dim());
2529   _elemlinks[0] = p;
2530 }
2531 
2532 
2533 
2534 inline
top_parent()2535 const Elem * Elem::top_parent () const
2536 {
2537   const Elem * tp = this;
2538 
2539   // Keep getting the element's parent
2540   // until that parent is at level-0
2541   while (tp->parent() != nullptr)
2542     tp = tp->parent();
2543 
2544   libmesh_assert(tp);
2545   libmesh_assert_equal_to (tp->level(), 0);
2546 
2547   return tp;
2548 }
2549 
2550 
2551 
2552 inline
level()2553 unsigned int Elem::level() const
2554 {
2555 #ifdef LIBMESH_ENABLE_AMR
2556 
2557   // if I don't have a parent I was
2558   // created directly from file
2559   // or by the user, so I am a
2560   // level-0 element
2561   if (this->parent() == nullptr)
2562     return 0;
2563 
2564   // if the parent and this element are of different
2565   // dimensionality we are at the same level as
2566   // the parent (e.g. we are the 2D side of a
2567   // 3D element)
2568   if (this->dim() != this->parent()->dim())
2569     return this->parent()->level();
2570 
2571   // otherwise we are at a level one
2572   // higher than our parent
2573   return (this->parent()->level() + 1);
2574 
2575 #else
2576 
2577   // Without AMR all elements are
2578   // at level 0.
2579   return 0;
2580 
2581 #endif
2582 }
2583 
2584 
2585 
2586 inline
p_level()2587 unsigned int Elem::p_level() const
2588 {
2589 #ifdef LIBMESH_ENABLE_AMR
2590   return _p_level;
2591 #else
2592   return 0;
2593 #endif
2594 }
2595 
2596 
2597 
2598 inline
mapping_type()2599 ElemMappingType Elem::mapping_type () const
2600 {
2601   return static_cast<ElemMappingType>(_map_type);
2602 }
2603 
2604 
2605 
2606 inline
set_mapping_type(const ElemMappingType type)2607 void Elem::set_mapping_type(const ElemMappingType type)
2608 {
2609   _map_type = cast_int<unsigned char>(type);
2610 }
2611 
2612 
2613 
2614 inline
mapping_data()2615 unsigned char Elem::mapping_data () const
2616 {
2617   return _map_data;
2618 }
2619 
2620 
2621 
2622 inline
set_mapping_data(const unsigned char data)2623 void Elem::set_mapping_data(const unsigned char data)
2624 {
2625   _map_data = data;
2626 }
2627 
2628 
2629 
2630 #ifdef LIBMESH_ENABLE_AMR
2631 
2632 inline
raw_child_ptr(unsigned int i)2633 const Elem * Elem::raw_child_ptr (unsigned int i) const
2634 {
2635   if (!_children)
2636     return nullptr;
2637 
2638   return _children[i];
2639 }
2640 
2641 inline
child_ptr(unsigned int i)2642 const Elem * Elem::child_ptr (unsigned int i) const
2643 {
2644   libmesh_assert(_children);
2645   libmesh_assert(_children[i]);
2646 
2647   return _children[i];
2648 }
2649 
2650 inline
child_ptr(unsigned int i)2651 Elem * Elem::child_ptr (unsigned int i)
2652 {
2653   libmesh_assert(_children);
2654   libmesh_assert(_children[i]);
2655 
2656   return _children[i];
2657 }
2658 
2659 
2660 inline
set_child(unsigned int c,Elem * elem)2661 void Elem::set_child (unsigned int c, Elem * elem)
2662 {
2663   libmesh_assert (this->has_children());
2664 
2665   _children[c] = elem;
2666 }
2667 
2668 
2669 
2670 inline
which_child_am_i(const Elem * e)2671 unsigned int Elem::which_child_am_i (const Elem * e) const
2672 {
2673   libmesh_assert(e);
2674   libmesh_assert (this->has_children());
2675 
2676   unsigned int nc = this->n_children();
2677   for (unsigned int c=0; c != nc; c++)
2678     if (this->child_ptr(c) == e)
2679       return c;
2680 
2681   libmesh_error_msg("ERROR:  which_child_am_i() was called with a non-child!");
2682 
2683   return libMesh::invalid_uint;
2684 }
2685 
2686 
2687 
2688 inline
refinement_flag()2689 Elem::RefinementState Elem::refinement_flag () const
2690 {
2691   return static_cast<RefinementState>(_rflag);
2692 }
2693 
2694 
2695 
2696 inline
set_refinement_flag(RefinementState rflag)2697 void Elem::set_refinement_flag(RefinementState rflag)
2698 {
2699   _rflag = cast_int<unsigned char>(rflag);
2700 }
2701 
2702 
2703 
2704 inline
p_refinement_flag()2705 Elem::RefinementState Elem::p_refinement_flag () const
2706 {
2707   return static_cast<RefinementState>(_pflag);
2708 }
2709 
2710 
2711 
2712 inline
set_p_refinement_flag(RefinementState pflag)2713 void Elem::set_p_refinement_flag(RefinementState pflag)
2714 {
2715   if (this->p_level() == 0)
2716     libmesh_assert_not_equal_to
2717       (pflag, Elem::JUST_REFINED);
2718 
2719   _pflag = cast_int<unsigned char>(pflag);
2720 }
2721 
2722 
2723 
2724 inline
max_descendant_p_level()2725 unsigned int Elem::max_descendant_p_level () const
2726 {
2727   // This is undefined for subactive elements,
2728   // which have no active descendants
2729   libmesh_assert (!this->subactive());
2730   if (this->active())
2731     return this->p_level();
2732 
2733   unsigned int max_p_level = _p_level;
2734   for (auto & c : child_ref_range())
2735     max_p_level = std::max(max_p_level,
2736                            c.max_descendant_p_level());
2737   return max_p_level;
2738 }
2739 
2740 
2741 
2742 inline
hack_p_level(unsigned int p)2743 void Elem::hack_p_level(unsigned int p)
2744 {
2745   if (p == 0)
2746     libmesh_assert_not_equal_to
2747       (this->p_refinement_flag(), Elem::JUST_REFINED);
2748 
2749   _p_level = cast_int<unsigned char>(p);
2750 }
2751 
2752 
2753 
2754 #endif // ifdef LIBMESH_ENABLE_AMR
2755 
2756 
2757 inline
compute_key(dof_id_type n0)2758 dof_id_type Elem::compute_key (dof_id_type n0)
2759 {
2760   return n0;
2761 }
2762 
2763 
2764 
2765 inline
compute_key(dof_id_type n0,dof_id_type n1)2766 dof_id_type Elem::compute_key (dof_id_type n0,
2767                                dof_id_type n1)
2768 {
2769   // Order the two so that n0 < n1
2770   if (n0 > n1) std::swap (n0, n1);
2771 
2772   return Utility::hashword2(n0, n1);
2773 }
2774 
2775 
2776 
2777 inline
compute_key(dof_id_type n0,dof_id_type n1,dof_id_type n2)2778 dof_id_type Elem::compute_key (dof_id_type n0,
2779                                dof_id_type n1,
2780                                dof_id_type n2)
2781 {
2782   std::array<dof_id_type, 3> array = {{n0, n1, n2}};
2783   std::sort(array.begin(), array.end());
2784   return Utility::hashword(array);
2785 }
2786 
2787 
2788 
2789 inline
compute_key(dof_id_type n0,dof_id_type n1,dof_id_type n2,dof_id_type n3)2790 dof_id_type Elem::compute_key (dof_id_type n0,
2791                                dof_id_type n1,
2792                                dof_id_type n2,
2793                                dof_id_type n3)
2794 {
2795   std::array<dof_id_type, 4> array = {{n0, n1, n2, n3}};
2796   std::sort(array.begin(), array.end());
2797   return Utility::hashword(array);
2798 }
2799 
2800 
2801 
2802 /**
2803  * The definition of the protected nested SideIter class.
2804  */
2805 class Elem::SideIter
2806 {
2807 public:
2808   // Constructor with arguments.
SideIter(const unsigned int side_number,Elem * parent)2809   SideIter(const unsigned int side_number,
2810            Elem * parent)
2811     : _side(),
2812       _side_ptr(nullptr),
2813       _parent(parent),
2814       _side_number(side_number)
2815   {}
2816 
2817 
2818   // Empty constructor.
SideIter()2819   SideIter()
2820     : _side(),
2821       _side_ptr(nullptr),
2822       _parent(nullptr),
2823       _side_number(libMesh::invalid_uint)
2824   {}
2825 
2826 
2827   // Copy constructor
SideIter(const SideIter & other)2828   SideIter(const SideIter & other)
2829     : _side(),
2830       _side_ptr(nullptr),
2831       _parent(other._parent),
2832       _side_number(other._side_number)
2833   {}
2834 
2835 
2836   // op=
2837   SideIter & operator=(const SideIter & other)
2838   {
2839     this->_parent      = other._parent;
2840     this->_side_number = other._side_number;
2841     return *this;
2842   }
2843 
2844   // unary op*
2845   Elem *& operator*() const
2846   {
2847     // Set the std::unique_ptr
2848     this->_update_side_ptr();
2849 
2850     // Return a reference to _side_ptr
2851     return this->_side_ptr;
2852   }
2853 
2854   // op++
2855   SideIter & operator++()
2856   {
2857     ++_side_number;
2858     return *this;
2859   }
2860 
2861   // op==  Two side iterators are equal if they have
2862   // the same side number and the same parent element.
2863   bool operator == (const SideIter & other) const
2864   {
2865     return (this->_side_number == other._side_number &&
2866             this->_parent      == other._parent);
2867   }
2868 
2869 
2870   // Consults the parent Elem to determine if the side
2871   // is a boundary side.  Note: currently side N is a
2872   // boundary side if neighbor N is nullptr.  Be careful,
2873   // this could possibly change in the future?
side_on_boundary()2874   bool side_on_boundary() const
2875   {
2876     return this->_parent->neighbor_ptr(_side_number) == nullptr;
2877   }
2878 
2879 private:
2880   // Update the _side pointer by building the correct side.
2881   // This has to be called before dereferencing.
_update_side_ptr()2882   void _update_side_ptr() const
2883   {
2884     // Construct new side, store in std::unique_ptr
2885     this->_side = this->_parent->build_side_ptr(this->_side_number);
2886 
2887     // Also set our internal naked pointer.  Memory is still owned
2888     // by the std::unique_ptr.
2889     this->_side_ptr = _side.get();
2890   }
2891 
2892   // std::unique_ptr to the actual side, handles memory management for
2893   // the sides which are created during the course of iteration.
2894   mutable std::unique_ptr<Elem> _side;
2895 
2896   // Raw pointer needed to facilitate passing back to the user a
2897   // reference to a non-temporary raw pointer in order to conform to
2898   // the variant_filter_iterator interface.  It points to the same
2899   // thing the std::unique_ptr "_side" above holds.  What happens if the user
2900   // calls delete on the pointer passed back?  Well, this is an issue
2901   // which is not addressed by the iterators in libMesh.  Basically it
2902   // is a bad idea to ever call delete on an iterator from the library.
2903   mutable Elem * _side_ptr;
2904 
2905   // Pointer to the parent Elem class which generated this iterator
2906   Elem * _parent;
2907 
2908   // A counter variable which keeps track of the side number
2909   unsigned int _side_number;
2910 };
2911 
2912 
2913 
2914 
2915 
2916 
2917 // Private implementation functions in the Elem class for the side iterators.
2918 // They have to come after the definition of the SideIter class.
2919 inline
_first_side()2920 Elem::SideIter Elem::_first_side()
2921 {
2922   return SideIter(0, this);
2923 }
2924 
2925 
2926 
2927 inline
_last_side()2928 Elem::SideIter Elem::_last_side()
2929 {
2930   return SideIter(this->n_neighbors(), this);
2931 }
2932 
2933 
2934 
2935 
2936 /**
2937  * The definition of the struct used for iterating over sides.
2938  */
2939 struct
2940 Elem::side_iterator : variant_filter_iterator<Elem::Predicate, Elem *>
2941 {
2942   // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
2943   template <typename PredType, typename IterType>
side_iteratorside_iterator2944   side_iterator (const IterType & d,
2945                  const IterType & e,
2946                  const PredType & p ) :
2947     variant_filter_iterator<Elem::Predicate, Elem *>(d,e,p) {}
2948 };
2949 
2950 
2951 
2952 inline
neighbor_ptr_range()2953 SimpleRange<Elem::NeighborPtrIter> Elem::neighbor_ptr_range()
2954 {
2955   return {_elemlinks+1, _elemlinks + 1 + this->n_neighbors()};
2956 }
2957 
2958 
2959 inline
neighbor_ptr_range()2960 SimpleRange<Elem::ConstNeighborPtrIter> Elem::neighbor_ptr_range() const
2961 {
2962   return {_elemlinks+1, _elemlinks + 1 + this->n_neighbors()};
2963 }
2964 
2965 } // namespace libMesh
2966 
2967 
2968 // Helper function for default caches in Elem subclasses
2969 
2970 #define LIBMESH_ENABLE_TOPOLOGY_CACHES                                  \
2971   virtual                                                               \
2972   std::vector<std::vector<std::vector<std::vector<std::pair<unsigned char, unsigned char>>>>> & \
2973   _get_bracketing_node_cache() const override                   \
2974   {                                                                     \
2975     static std::vector<std::vector<std::vector<std::vector<std::pair<unsigned char, unsigned char>>>>> c; \
2976     return c;                                                           \
2977   }                                                                     \
2978                                                                         \
2979   virtual                                                               \
2980   std::vector<std::vector<std::vector<signed char>>> &                  \
2981   _get_parent_indices_cache() const override                    \
2982   {                                                                     \
2983     static std::vector<std::vector<std::vector<signed char>>> c;        \
2984     return c;                                                           \
2985   }
2986 
2987 
2988 
2989 
2990 
2991 
2992 #endif // LIBMESH_ELEM_H
2993