Line data Source code
1 : // The libMesh Finite Element Library.
2 : // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 :
4 : // This library is free software; you can redistribute it and/or
5 : // modify it under the terms of the GNU Lesser General Public
6 : // License as published by the Free Software Foundation; either
7 : // version 2.1 of the License, or (at your option) any later version.
8 :
9 : // This library is distributed in the hope that it will be useful,
10 : // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 : // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 : // Lesser General Public License for more details.
13 :
14 : // You should have received a copy of the GNU Lesser General Public
15 : // License along with this library; if not, write to the Free Software
16 : // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 :
18 :
19 :
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 : // C++ includes
39 : #include <algorithm>
40 : #include <cstddef>
41 : #include <iostream>
42 : #include <limits.h> // CHAR_BIT
43 : #include <set>
44 : #include <vector>
45 : #include <memory>
46 : #include <array>
47 :
48 : namespace libMesh
49 : {
50 :
51 : // Forward declarations
52 : class BoundaryInfo;
53 : class Elem;
54 : class MeshBase;
55 : class MeshRefinement;
56 : #ifdef LIBMESH_ENABLE_PERIODIC
57 : class PeriodicBoundaries;
58 : class PointLocatorBase;
59 : #endif
60 : template <class SideType, class ParentType>
61 : class Side;
62 : enum ElemQuality : int;
63 : enum IOPackage : int;
64 : enum Order : int;
65 :
66 :
67 : /**
68 : * This is the base class from which all geometric element types are
69 : * derived. The \p Elem class provides standard information such as
70 : * the number of nodes, edges, faces, vertices, children, and
71 : * neighbors it has, as well as access to (or the ability to
72 : * construct) these entities.
73 : *
74 : * An \p Elem has pointers to its \p Node objects. Some of these
75 : * nodes live at the vertices of the element, while others may live on
76 : * edges (and faces in 3D), or interior to the element. The number of
77 : * nodes in a given element, \p n_nodes(), is encoded into the name of
78 : * the class. For example, a \p Tri3 has three nodes which correspond
79 : * to the vertices, while a \p Tri6 has six nodes, three of which are
80 : * located at vertices, and three which are located at the midpoint of
81 : * each edge. Nodes on edges, faces, and element interiors are called
82 : * second-order nodes.
83 : *
84 : * A 1D Elem is an \p Edge, a 2D Elem is a \p Face, and a 3D Elem is a
85 : * \p Cell. An \p Elem is composed of a number of sides, which can
86 : * be accessed as dim-1 dimensional \p Elem types. For example, a \p
87 : * Hex8 is a 3D hexahedral element. A \p Hex8 has 6 sides, which are
88 : * \p Faces of type Quad4.
89 : *
90 : * \author Benjamin S. Kirk
91 : * \date 2002-2007
92 : * \brief The base class for all geometric element types.
93 : */
94 : class Elem : public ReferenceCountedObject<Elem>,
95 : public DofObject
96 : {
97 : protected:
98 :
99 : /**
100 : * Constructor. Creates an element with \p n_nodes nodes,
101 : * \p n_sides sides, \p n_children possible children, and
102 : * parent \p p. The constructor allocates the memory necessary
103 : * to support this data.
104 : */
105 : Elem (const unsigned int n_nodes,
106 : const unsigned int n_sides,
107 : Elem * parent,
108 : Elem ** elemlinkdata,
109 : Node ** nodelinkdata);
110 :
111 : public:
112 :
113 : /**
114 : * Elems are responsible for allocating and deleting space for
115 : * storing pointers to their children during refinement, so they
116 : * cannot currently be (default) copy-constructed or copy-
117 : * assigned. We therefore explicitly delete these operations. In
118 : * addition, the DofObject base class currently has private copy
119 : * construction and assignment operators, so that prevents us from
120 : * copying Elems as well.
121 : */
122 : Elem (Elem &&) = delete;
123 : Elem (const Elem &) = delete;
124 : Elem & operator= (const Elem &) = delete;
125 : Elem & operator= (Elem &&) = delete;
126 :
127 : /**
128 : * Destructor.
129 : */
130 4293778012 : virtual ~Elem() = default;
131 :
132 : /**
133 : * \returns The \p Point associated with local \p Node \p i.
134 : */
135 : const Point & point (const unsigned int i) const;
136 :
137 : /**
138 : * \returns The \p Point associated with local \p Node \p i
139 : * as a writable reference.
140 : */
141 : Point & point (const unsigned int i);
142 :
143 : /**
144 : * \returns The \p Point associated with local \p Node \p i,
145 : * in master element rather than physical coordinates.
146 : */
147 : virtual Point master_point (const unsigned int i) const = 0;
148 :
149 : /**
150 : * \returns The global id number of local \p Node \p i.
151 : */
152 : dof_id_type node_id (const unsigned int i) const;
153 :
154 : /**
155 : * \returns The local id number of global \p Node id \p i,
156 : * or \p invalid_uint if Node id \p i is not local.
157 : */
158 : unsigned int local_node (const dof_id_type i) const;
159 :
160 : /**
161 : * \returns The local index for the \p Node pointer \p node_ptr,
162 : * or \p invalid_uint if \p node_ptr is not a local node.
163 : */
164 : unsigned int get_node_index (const Node * node_ptr) const;
165 :
166 : /**
167 : * \returns A pointer to an array of local node pointers.
168 : */
169 : const Node * const * get_nodes () const;
170 :
171 : /**
172 : * \returns A const pointer to local \p Node \p i.
173 : */
174 : const Node * node_ptr (const unsigned int i) const;
175 :
176 : /**
177 : * \returns A non-const pointer to local \p Node \p i.
178 : */
179 : Node * node_ptr (const unsigned int i);
180 :
181 : /**
182 : * \returns A const reference to local \p Node \p i.
183 : */
184 : const Node & node_ref (const unsigned int i) const;
185 :
186 : /**
187 : * \returns A writable reference to local \p Node \p i.
188 : */
189 : Node & node_ref (const unsigned int i);
190 :
191 : #ifdef LIBMESH_ENABLE_DEPRECATED
192 : /**
193 : * \returns The pointer to the \p Node with local number \p i as a
194 : * writable reference.
195 : *
196 : * \deprecated This setter cannot update the multiple node pointers
197 : * used in a general polyhedron; use the \p set_node overload that
198 : * takes an argument.
199 : */
200 : virtual Node * & set_node (const unsigned int i);
201 : #endif // LIBMESH_ENABLE_DEPRECATED
202 :
203 : /**
204 : * Sets local \p Node \p i to refer to \p node.
205 : */
206 : virtual void set_node (const unsigned int i,
207 : Node * node);
208 :
209 : /**
210 : * Nested classes for use iterating over all nodes of an element.
211 : */
212 : class NodeRefIter;
213 : class ConstNodeRefIter;
214 :
215 : /**
216 : * Returns a range with all nodes of an element, usable in
217 : * range-based for loops. The exact type of the return value here
218 : * may be subject to change in future libMesh releases, but the
219 : * iterators will always dereference to produce a reference to a
220 : * Node.
221 : */
222 : SimpleRange<NodeRefIter> node_ref_range();
223 :
224 : SimpleRange<ConstNodeRefIter> node_ref_range() const;
225 :
226 : /**
227 : * \returns The subdomain that this element belongs to.
228 : */
229 : subdomain_id_type subdomain_id () const;
230 :
231 : /**
232 : * \returns The subdomain that this element belongs to as a
233 : * writable reference.
234 : */
235 : subdomain_id_type & subdomain_id ();
236 :
237 : /**
238 : * A static integral constant representing an invalid subdomain id.
239 : * See also DofObject::{invalid_id, invalid_unique_id, invalid_processor_id}.
240 : *
241 : * \note We don't use the static_cast(-1) trick here since
242 : * \p subdomain_id_type is sometimes a *signed* integer for
243 : * compatibility reasons (see libmesh/id_types.h).
244 : */
245 : static constexpr subdomain_id_type invalid_subdomain_id
246 : = std::numeric_limits<subdomain_id_type>::max();
247 :
248 : /**
249 : * \returns true iff this element type can vary in topology (e.g.
250 : * have different numbers of sides and/or nodes) at runtime. For
251 : * such general polygons or polyhedra, APIs which assume a fixed
252 : * topology are not safe to use.
253 : */
254 290152832 : virtual bool runtime_topology() const { return false; }
255 :
256 : /**
257 : * \returns A pointer to the "reference element" associated
258 : * with this element. The reference element is the image of this
259 : * element in reference parametric space. Importantly, it is *not*
260 : * an actual element in the mesh, but rather a Singleton-type
261 : * object, so for example all \p Quad4 elements share the same
262 : * \p reference_elem().
263 : *
264 : * If the element is of a type that can admit multiple topologies,
265 : * such as a Polygon subtype, then there is no reference element;
266 : * for such types this method should not be used.
267 : */
268 : const Elem * reference_elem () const;
269 :
270 : /**
271 : * \returns An id associated with the \p s side of this element.
272 : * The id is not necessarily unique, but should be close.
273 : */
274 : virtual dof_id_type key (const unsigned int s) const = 0;
275 :
276 : /**
277 : * \returns An id associated with the \p s side of this element, as
278 : * defined solely by element vertices. The id is not necessarily
279 : * unique, but should be close. This is particularly useful in the
280 : * \p MeshBase::find_neighbors() routine.
281 : */
282 : virtual dof_id_type low_order_key (const unsigned int s) const = 0;
283 :
284 : /**
285 : * \returns An id associated with the global node ids of this
286 : * element. The id is not necessarily unique, but should be
287 : * close. Uses the same hash as the key(s) function, so for example
288 : * if "tri3" is side 0 of "tet4", then tri3->key()==tet4->key(0).
289 : */
290 : virtual dof_id_type key () const;
291 :
292 : /**
293 : * \returns \p true if two elements are equivalent, false otherwise.
294 : * This is true if the elements are connected to identical global
295 : * nodes, regardless of how those nodes might be numbered local
296 : * to the elements.
297 : */
298 : bool operator == (const Elem & rhs) const;
299 :
300 : /**
301 : * \returns \p true if two elements have equal topologies, false
302 : * otherwise.
303 : * This is true if the elements connect to nodes of the same id in
304 : * the same order, and neighbors of the same id on each side, the
305 : * same id on any parent and/or interior_parent link, etc.
306 : */
307 : bool topologically_equal (const Elem & rhs) const;
308 :
309 : /**
310 : * \returns A const pointer to the \f$ i^{th} \f$ neighbor of this
311 : * element, or \p nullptr if \p MeshBase::find_neighbors() has not been
312 : * called.
313 : *
314 : * \note If \p MeshBase::find_neighbors() has been called and this
315 : * function still returns \p nullptr, then the side is on a boundary of
316 : * the domain.
317 : */
318 : const Elem * neighbor_ptr (unsigned int i) const;
319 :
320 : /**
321 : * \returns A non-const pointer to the \f$ i^{th} \f$ neighbor of this element.
322 : */
323 : Elem * neighbor_ptr (unsigned int i);
324 :
325 : /**
326 : * Nested "classes" for use iterating over all neighbors of an element.
327 : */
328 : typedef Elem * const * NeighborPtrIter;
329 : typedef const Elem * const * ConstNeighborPtrIter;
330 :
331 : /**
332 : * Returns a range with all neighbors of an element, usable in
333 : * range-based for loops. The exact type of the return value here
334 : * may be subject to change in future libMesh releases, but the
335 : * iterators will always dereference to produce a pointer to a
336 : * neighbor element (or a null pointer, for sides which have no
337 : * neighbors).
338 : */
339 : SimpleRange<NeighborPtrIter> neighbor_ptr_range();
340 :
341 : SimpleRange<ConstNeighborPtrIter> neighbor_ptr_range() const;
342 :
343 : #ifdef LIBMESH_ENABLE_PERIODIC
344 : /**
345 : * \returns A pointer to the \f$ i^{th} \f$ neighbor of this element
346 : * for interior elements. If an element is on a periodic
347 : * boundary, it will return a corresponding element on the opposite
348 : * side.
349 : */
350 : const Elem * topological_neighbor (const unsigned int i,
351 : const MeshBase & mesh,
352 : const PointLocatorBase & point_locator,
353 : const PeriodicBoundaries * pb) const;
354 :
355 : /**
356 : * \returns A writable pointer to the \f$ i^{th} \f$ neighbor of
357 : * this element for interior elements. If an element is on a
358 : * periodic boundary, it will return a corresponding element on the
359 : * opposite side.
360 : */
361 : Elem * topological_neighbor (const unsigned int i,
362 : MeshBase & mesh,
363 : const PointLocatorBase & point_locator,
364 : const PeriodicBoundaries * pb);
365 :
366 : /**
367 : * \returns \p true if the element \p elem in question is a neighbor or
368 : * topological neighbor of this element, \p false otherwise.
369 : */
370 : bool has_topological_neighbor (const Elem * elem,
371 : const MeshBase & mesh,
372 : const PointLocatorBase & point_locator,
373 : const PeriodicBoundaries * pb) const;
374 : #endif
375 :
376 : /**
377 : * Assigns \p n as the \f$ i^{th} \f$ neighbor.
378 : */
379 : void set_neighbor (const unsigned int i, Elem * n);
380 :
381 : /**
382 : * \returns \p true if the element \p elem in question is a neighbor
383 : * of this element, \p false otherwise.
384 : */
385 : bool has_neighbor (const Elem * elem) const;
386 :
387 : /**
388 : * \returns If \p elem is a neighbor of a child of this element, a
389 : * pointer to that child, otherwise \p nullptr.
390 : */
391 : Elem * child_neighbor (Elem * elem);
392 :
393 : /**
394 : * \returns If \p elem is a neighbor of a child of this element, a
395 : * pointer to that child, otherwise \p nullptr.
396 : */
397 : const Elem * child_neighbor (const Elem * elem) const;
398 :
399 : /**
400 : * \returns \p true if this element has a side coincident
401 : * with a boundary (indicated by a \p nullptr neighbor), \p false
402 : * otherwise.
403 : */
404 : bool on_boundary () const;
405 :
406 : /**
407 : * \returns \p true if this element is "semilocal" to the calling
408 : * processor, which must specify its rank.
409 : *
410 : * This method is discouraged, as it uses the *old* definition of
411 : * semilocal (elements which are not local but which are point
412 : * neighbors of something local) rather than any of the new
413 : * definitions discussed in ghosting_functor.h
414 : */
415 : bool is_semilocal (const processor_id_type my_pid) const;
416 :
417 : /**
418 : * This function tells you which neighbor \p e is.
419 : * I.e. if s = a->which_neighbor_am_i(e); then
420 : * a->neighbor(s) will be an ancestor of e.
421 : */
422 : unsigned int which_neighbor_am_i(const Elem * e) const;
423 :
424 : /**
425 : * This function tells you which side the boundary element \p e is.
426 : * I.e. if e = a->build_side_ptr(s) or e = a->side_ptr(s); then
427 : * a->which_side_am_i(e) will be s.
428 : *
429 : * \note An \e exact floating point comparison of the nodal
430 : * positions of \p e is made with the nodal positions of \p this in
431 : * order to perform this test. The idea is that the test will return
432 : * a valid side id if \p e either directly shares Node pointers with
433 : * \p this, or was created by exactly copying some of the nodes of
434 : * \p this (e.g. through BoundaryMesh::sync()). In these
435 : * circumstances, non-fuzzy floating point equality is expected.
436 : *
437 : * \returns The side of \p this the element which \p e is, otherwise
438 : * \p invalid_uint.
439 : */
440 : unsigned int which_side_am_i(const Elem * e) const;
441 :
442 : /**
443 : * \returns The local node id for node \p side_node on side \p side of
444 : * this Elem. Simply relies on the \p side_nodes_map for each of the
445 : * derived types. For example,
446 : * Tri3::local_side_node(0, 0) -> 0
447 : * Tri3::local_side_node(0, 1) -> 1
448 : * Tri3::local_side_node(1, 0) -> 1
449 : * Tri3::local_side_node(1, 1) -> 2
450 : * etc...
451 : */
452 : virtual unsigned int local_side_node(unsigned int side,
453 : unsigned int side_node) const = 0;
454 :
455 : /**
456 : * Similar to Elem::local_side_node(), but instead of a side id, takes
457 : * an edge id and a node id on that edge and returns a local node number
458 : * for the Elem. The implementation relies on the "edge_nodes_map" tables
459 : * for 3D elements. For 2D elements, calls local_side_node(). Throws an
460 : * error if called on 1D elements.
461 : */
462 : virtual unsigned int local_edge_node(unsigned int edge,
463 : unsigned int edge_node) const = 0;
464 :
465 : /**
466 : * \returns \p true if a vertex of \p e is contained
467 : * in this element. If \p mesh_connection is true, looks
468 : * specifically for containment possibilities of an element \p e
469 : * that is connected to \p this via membership in the same manifold
470 : * of the same mesh.
471 : */
472 : bool contains_vertex_of(const Elem * e, bool mesh_connection=false) const;
473 :
474 : /**
475 : * \returns \p true if an edge of \p e is contained in
476 : * this element. (Internally, this is done by checking whether at
477 : * least two vertices of \p e are contained in this element).
478 : */
479 : bool contains_edge_of(const Elem * e) const;
480 :
481 : /**
482 : * This function finds all active elements (including this one)
483 : * which are in the same manifold as this element and which touch
484 : * the current active element at the specified point, which should
485 : * be a point in the current element.
486 : *
487 : * Elements which are not "in the same manifold" (e.g. the
488 : * interior_parent of a boundary element) will not be found with
489 : * this method.
490 : *
491 : * Elements which overlap the specified point but which are only
492 : * connected to the current element via elements which do not
493 : * overlap that point (e.g. in a folded or tangled mesh) are not
494 : * considered to "touch" the current element and will not be found
495 : * with this method.
496 : */
497 : void find_point_neighbors(const Point & p,
498 : std::set<const Elem *> & neighbor_set) const;
499 :
500 : /**
501 : * This function finds all active elements (including this one) in
502 : * the same manifold as this element which touch this active element
503 : * at any point.
504 : */
505 : void find_point_neighbors(std::set<const Elem *> & neighbor_set) const;
506 :
507 : /**
508 : * This function finds all active elements (including this one) in
509 : * the same manifold as start_elem (which must be active and must
510 : * touch this element) which touch this element at any point.
511 : */
512 : void find_point_neighbors(std::set<const Elem *> & neighbor_set,
513 : const Elem * start_elem) const;
514 :
515 : /**
516 : * Non-const version of function above. Fills a set of non-const Elem pointers.
517 : */
518 : void find_point_neighbors(std::set<Elem *> & neighbor_set,
519 : Elem * start_elem);
520 :
521 : /**
522 : * This function finds all active elements in the same manifold as
523 : * this element which touch the current active element along the
524 : * whole edge defined by the two points \p p1 and \p p2.
525 : */
526 : void find_edge_neighbors(const Point & p1,
527 : const Point & p2,
528 : std::set<const Elem *> & neighbor_set) const;
529 :
530 : /**
531 : * This function finds all active elements in the same manifold as
532 : * this element which touch the current active element along any
533 : * edge (more precisely, at at least two points).
534 : *
535 : * In this case, elements are included even if they do not touch a
536 : * *whole* edge of this element.
537 : */
538 : void find_edge_neighbors(std::set<const Elem *> & neighbor_set) const;
539 :
540 : /**
541 : * This function finds all active elements (*not* including this
542 : * one) in the parent manifold of this element whose intersection
543 : * with this element has non-zero measure.
544 : */
545 : void find_interior_neighbors(std::set<const Elem *> & neighbor_set) const;
546 :
547 : /**
548 : * Non-const version of function above that fills up a vector of
549 : * non-const Elem pointers instead.
550 : */
551 : void find_interior_neighbors(std::set<Elem *> & neighbor_set);
552 :
553 : /**
554 : * Resets this element's neighbors' appropriate neighbor pointers
555 : * and its parent's and children's appropriate pointers
556 : * to point to null instead of to this.
557 : *
558 : * To be used before an element is deleted from a mesh.
559 : */
560 : void remove_links_to_me ();
561 :
562 : /**
563 : * Resets this element's neighbors' appropriate neighbor pointers
564 : * and its parent's and children's appropriate pointers
565 : * to point to the global remote_elem instead of this.
566 : * Used by the library before an element becomes remote on the
567 : * local processor.
568 : */
569 : void make_links_to_me_remote ();
570 :
571 : /**
572 : * Resets the \p neighbor_side pointers of our nth neighbor (and
573 : * its descendants, if appropriate) to point to this Elem instead of
574 : * to the global remote_elem. Used by the library when a formerly
575 : * remote element is being added to the local processor.
576 : */
577 : void make_links_to_me_local (unsigned int n, unsigned int neighbor_side);
578 :
579 : /**
580 : * \returns \p true if this element is remote, false otherwise.
581 : *
582 : * A remote element (see \p RemoteElem) is a syntactic convenience --
583 : * it is a placeholder for an element which exists on some other
584 : * processor. Local elements are required to have valid neighbors,
585 : * and these ghost elements may have remote neighbors for data
586 : * structure consistency. The use of remote elements helps ensure
587 : * that any element we may access has a \p nullptr neighbor only if it
588 : * lies on the physical boundary of the domain.
589 : */
590 7280126 : virtual bool is_remote () const
591 7280126 : { return false; }
592 :
593 : /**
594 : * \returns The connectivity for this element in a specific
595 : * format, which is specified by the IOPackage tag.
596 : */
597 : virtual void connectivity(const unsigned int sc,
598 : const IOPackage iop,
599 : std::vector<dof_id_type> & conn) const = 0;
600 :
601 : /**
602 : * Writes the element connectivity for various IO packages
603 : * to the passed ostream "out". Not virtual, since it is
604 : * implemented in the base class.
605 : */
606 : void write_connectivity (std::ostream & out,
607 : const IOPackage iop) const;
608 :
609 : /**
610 : * \returns The type of element that has been derived from this
611 : * base class.
612 : */
613 : virtual ElemType type () const = 0;
614 :
615 : /**
616 : * This array maps the integer representation of the \p ElemType enum
617 : * to the geometric dimension of the element.
618 : *
619 : * This is currently usable even for complicated subclasses with
620 : * runtime-varying topology.
621 : */
622 : static const unsigned int type_to_dim_map[INVALID_ELEM];
623 :
624 : /**
625 : * \returns The dimensionality of the object.
626 : */
627 : virtual unsigned short dim () const = 0;
628 :
629 : /**
630 : * This array maps the integer representation of the \p ElemType enum
631 : * to the number of nodes in the element.
632 : *
633 : * This is only usable for simple types for which the node number
634 : * is fixed; for more general types like Polygon subclasses an actual
635 : * instantiated Elem must be queried.
636 : */
637 : static const unsigned int type_to_n_nodes_map[INVALID_ELEM];
638 :
639 : /**
640 : * \returns The number of nodes this element contains.
641 : */
642 : virtual unsigned int n_nodes () const = 0;
643 :
644 : /**
645 : * The maximum number of nodes *any* element can contain.
646 : * This is useful for replacing heap vectors with stack arrays.
647 : */
648 : static const unsigned int max_n_nodes = 27;
649 :
650 : /**
651 : * \returns An integer range from 0 up to (but not including)
652 : * the number of nodes this element contains.
653 : */
654 : IntRange<unsigned short> node_index_range () const;
655 :
656 : /**
657 : * \returns The number of nodes the given child of this element
658 : * contains. Except in odd cases like pyramid refinement this will
659 : * be the same as the number of nodes in the parent element.
660 : */
661 18635211 : virtual unsigned int n_nodes_in_child (unsigned int /*c*/) const
662 18635211 : { return this->n_nodes(); }
663 :
664 : /**
665 : * This array maps the integer representation of the \p ElemType enum
666 : * to the number of sides on the element.
667 : *
668 : * This is only usable for simple types for which the node number
669 : * is fixed; for more general types like Polygon subclasses an actual
670 : * instantiated Elem must be queried.
671 : */
672 : static const unsigned int type_to_n_sides_map[INVALID_ELEM];
673 :
674 : /**
675 : * \returns The number of sides the element that has been derived
676 : * from this class has. In 2D the number of sides is the number
677 : * of edges, in 3D the number of sides is the number of faces.
678 : */
679 : virtual unsigned int n_sides () const = 0;
680 :
681 : /**
682 : * \returns The type of element for side \p s.
683 : */
684 : virtual ElemType side_type (const unsigned int s) const = 0;
685 :
686 : /**
687 : * \returns the normal (outwards-facing) of the side of the element at the vertex-average of the side
688 : * @param s the side of interest
689 : */
690 : virtual Point side_vertex_average_normal(const unsigned int s) const;
691 :
692 : /**
693 : * \returns An integer range from 0 up to (but not including)
694 : * the number of sides this element has.
695 : */
696 : IntRange<unsigned short> side_index_range () const;
697 :
698 : /**
699 : * \returns The number of neighbors the element that has been derived
700 : * from this class has.
701 : *
702 : * Only face (or edge in 2D) neighbors are stored, so this method
703 : * returns n_sides(). At one point we intended to allow derived
704 : * classes to override this, but too much current libMesh code
705 : * assumes n_neighbors==n_sides.
706 : */
707 85299598 : unsigned int n_neighbors () const
708 152254805 : { return this->n_sides(); }
709 :
710 : /**
711 : * \returns The number of vertices the element that has been derived
712 : * from this class has.
713 : */
714 : virtual unsigned int n_vertices () const = 0;
715 :
716 : /**
717 : * \returns The number of edges the element that has been derived
718 : * from this class has.
719 : */
720 : virtual unsigned int n_edges () const = 0;
721 :
722 : /**
723 : * \returns An integer range from 0 up to (but not including)
724 : * the number of edges this element has.
725 : */
726 : IntRange<unsigned short> edge_index_range () const;
727 :
728 : /**
729 : * This array maps the integer representation of the \p ElemType enum
730 : * to the number of edges on the element.
731 : *
732 : * This is only usable for simple types for which the node number
733 : * is fixed; for more general types like Polygon subclasses an actual
734 : * instantiated Elem must be queried.
735 : */
736 : static const unsigned int type_to_n_edges_map[INVALID_ELEM];
737 :
738 : /**
739 : * \returns The number of faces the element that has been derived
740 : * from this class has.
741 : */
742 : virtual unsigned int n_faces () const = 0;
743 :
744 : /**
745 : * \returns An integer range from 0 up to (but not including)
746 : * the number of faces this element has.
747 : */
748 : IntRange<unsigned short> face_index_range () const;
749 :
750 : /**
751 : * \returns The number of children the element that has been derived
752 : * from this class may have.
753 : */
754 : virtual unsigned int n_children () const = 0;
755 :
756 : /**
757 : * \returns \p true if the specified (local) node number is a vertex node.
758 : */
759 : virtual bool is_vertex(const unsigned int i) const = 0;
760 :
761 : /**
762 : * \returns \p true if the specified child has a vertex at the
763 : * specified (child-local) node number.
764 : * Except in odd cases like pyramid refinement the child will have
765 : * the same local structure as the parent element.
766 : */
767 169263 : virtual bool is_vertex_on_child (unsigned int /*c*/,
768 : unsigned int n) const
769 169263 : { return this->is_vertex(n); }
770 :
771 : /**
772 : * \returns \p true if this element has a vertex at the specified
773 : * (child-local) node number \p n of the specified child \p c.
774 : */
775 : virtual bool is_vertex_on_parent(unsigned int c,
776 : unsigned int n) const;
777 :
778 : /**
779 : * \returns \p true if the specified (local) node number is an edge node.
780 : * For 1D elements, is_edge() is equivalent to is_internal().
781 : */
782 : virtual bool is_edge(const unsigned int i) const = 0;
783 :
784 : /**
785 : * \returns \p true if the specified (local) node number is a face node.
786 : * For 2D elements, is_face() is equivalent to is_internal().
787 : * For 1D elements, is_face() == false.
788 : */
789 : virtual bool is_face(const unsigned int i) const = 0;
790 :
791 : /**
792 : * \returns \p true if the specified (local) node number is an internal node.
793 : */
794 : bool is_internal(const unsigned int i) const;
795 :
796 : /**
797 : * \returns \p true if the specified (local) node number is on the
798 : * specified side.
799 : */
800 : virtual bool is_node_on_side(const unsigned int n,
801 : const unsigned int s) const = 0;
802 :
803 : /**
804 : * \returns the (local) node numbers on the specified side
805 : */
806 : virtual std::vector<unsigned int> nodes_on_side(const unsigned int /*s*/) const = 0;
807 :
808 : /**
809 : * \returns the (local) node numbers on the specified edge
810 : */
811 : virtual std::vector<unsigned int> nodes_on_edge(const unsigned int /*e*/) const = 0;
812 :
813 : /**
814 : * \returns the (local) side numbers that touch the specified edge
815 : */
816 : virtual std::vector<unsigned int> sides_on_edge(const unsigned int /*e*/) const = 0;
817 :
818 : /**
819 : * \returns the (local) edge numbers that touch the specified node
820 : */
821 : virtual std::vector<unsigned int> edges_adjacent_to_node(const unsigned int /*n*/) const = 0;
822 :
823 : /**
824 : * \returns \p true if the specified (local) node number is on the
825 : * specified edge.
826 : */
827 : virtual bool is_node_on_edge(const unsigned int n,
828 : const unsigned int e) const = 0;
829 :
830 : /**
831 : * \returns \p true if the specified edge is on the specified side.
832 : */
833 : virtual bool is_edge_on_side(const unsigned int e,
834 : const unsigned int s) const = 0;
835 :
836 : /**
837 : * \returns The side number opposite to \p s (for a tensor product
838 : * element), or throws an error otherwise.
839 : */
840 : virtual unsigned int opposite_side(const unsigned int s) const;
841 :
842 : /**
843 : * \returns The local node number for the node opposite to node n
844 : * on side \p opposite_side(s) (for a tensor product element), or
845 : * throws an error otherwise.
846 : */
847 : virtual unsigned int opposite_node(const unsigned int n,
848 : const unsigned int s) const;
849 :
850 : /**
851 : * \returns The number of sub-elements this element may be broken
852 : * down into for visualization purposes. For example, 1 for a
853 : * linear triangle, 4 for a quadratic (6-noded) triangle, etc...
854 : */
855 : virtual unsigned int n_sub_elem () const = 0;
856 :
857 : /**
858 : * \returns A temporary element coincident with side \p i.
859 : *
860 : * This method returns the _minimum_ element necessary to uniquely
861 : * identify the side. For example, the side of a hexahedron is
862 : * always returned as a 4-noded quadrilateral, regardless of what
863 : * type of hex you are dealing with. Important data like subdomain
864 : * id, p level, or mapping type may be omitted from the temporary
865 : * element. If you want a first-class full-ordered face (i.e. a
866 : * 9-noded quad face for a 27-noded hexahedron), use the
867 : * build_side_ptr method.
868 : *
869 : * \note The const version of this function is non-virtual; it
870 : * simply calls the virtual non-const version and const_casts the
871 : * return type.
872 : */
873 : virtual std::unique_ptr<Elem> side_ptr (unsigned int i) = 0;
874 : std::unique_ptr<const Elem> side_ptr (unsigned int i) const;
875 :
876 : /**
877 : * Resets the loose element \p side, which may currently point to a
878 : * different side than \p i or even a different element than \p
879 : * this, to point to side \p i on \p this. If \p side is currently
880 : * an element of the wrong type, it will be freed and a new element
881 : * allocated; otherwise no memory allocation will occur.
882 : *
883 : * This will cause \p side to be a minimum-ordered element, even if
884 : * it is handed a higher-ordered element that must be replaced.
885 : *
886 : * The const version of this function is non-virtual; it simply
887 : * calls the virtual non-const version and const_casts the return
888 : * type.
889 : */
890 : virtual void side_ptr (std::unique_ptr<Elem> & side, const unsigned int i) = 0;
891 : void side_ptr (std::unique_ptr<const Elem> & side, const unsigned int i) const;
892 :
893 : /**
894 : * \returns An temporary element coincident with side \p i wrapped
895 : * in a smart pointer.
896 : *
897 : * The element returned is full-ordered and full-featured, in
898 : * contrast to the side method. For example, calling
899 : * build_side_ptr(0) on a 20-noded hex in subdomain 5 will build a
900 : * 8-noded quadrilateral coincident with face 0, assign it subdomain
901 : * id 5, and pass back the pointer.
902 : *
903 : * The side element's id() is undefined; it is a temporary element
904 : * not added to any mesh.
905 : *
906 : * A \p std::unique_ptr<Elem> is returned to prevent a memory leak.
907 : * This way the user need not remember to delete the object.
908 : *
909 : * The const version of this function is non-virtual; it simply
910 : * calls the virtual non-const version and const_casts the return
911 : * type.
912 : */
913 : virtual std::unique_ptr<Elem> build_side_ptr (const unsigned int i) = 0;
914 : std::unique_ptr<const Elem> build_side_ptr (const unsigned int i) const;
915 :
916 : #ifdef LIBMESH_ENABLE_DEPRECATED
917 : /*
918 : * Older versions of libMesh supported a "proxy" option here.
919 : */
920 0 : virtual std::unique_ptr<Elem> build_side_ptr (const unsigned int i, bool proxy)
921 0 : { if (proxy) libmesh_error(); libmesh_deprecated(); return this->build_side_ptr(i); }
922 :
923 : std::unique_ptr<const Elem> build_side_ptr (const unsigned int i, bool proxy) const
924 : { if (proxy) libmesh_error(); libmesh_deprecated(); return this->build_side_ptr(i); }
925 : #endif
926 :
927 : /**
928 : * Resets the loose element \p side, which may currently point to a
929 : * different side than \p i or even a different element than \p
930 : * this, to point to side \p i on \p this. If \p side is currently
931 : * an element of the wrong type, it will be freed and a new element
932 : * allocated; otherwise no memory allocation will occur.
933 : *
934 : * This will cause \p side to be a full-ordered element, even if it
935 : * is handed a lower-ordered element that must be replaced.
936 : *
937 : * The const version of this function is non-virtual; it simply
938 : * calls the virtual non-const version and const_casts the return
939 : * type.
940 : */
941 : virtual void build_side_ptr (std::unique_ptr<Elem> & side, const unsigned int i) = 0;
942 : void build_side_ptr (std::unique_ptr<const Elem> & side, const unsigned int i) const;
943 :
944 : /**
945 : * \returns An element coincident with edge \p i wrapped in a smart pointer.
946 : *
947 : * The element returned is full-ordered. For example, calling
948 : * build_edge_ptr(0) on a 20-noded hex will build a 3-noded edge
949 : * coincident with edge 0 and pass back the pointer. A \p
950 : * std::unique_ptr<Elem> is returned to prevent a memory leak. This way
951 : * the user need not remember to delete the object.
952 : *
953 : * The const version of this function is non-virtual; it simply
954 : * calls the virtual non-const version and const_casts the return
955 : * type.
956 : */
957 : virtual std::unique_ptr<Elem> build_edge_ptr (const unsigned int i) = 0;
958 : std::unique_ptr<const Elem> build_edge_ptr (const unsigned int i) const;
959 :
960 : /**
961 : * Resets the loose element \p edge, which may currently point to a
962 : * different edge than \p i or even a different element than \p
963 : * this, to point to edge \p i on \p this. If \p edge is currently
964 : * an element of the wrong type, it will be freed and a new element
965 : * allocated; otherwise no memory allocation will occur.
966 : *
967 : * This will cause \p edge to be a full-ordered element, even if it
968 : * is handed a lower-ordered element that must be replaced.
969 : *
970 : * The const version of this function is non-virtual; it simply
971 : * calls the virtual non-const version and const_casts the return
972 : * type.
973 : */
974 : virtual void build_edge_ptr (std::unique_ptr<Elem> & edge, const unsigned int i) = 0;
975 : void build_edge_ptr (std::unique_ptr<const Elem> & edge, const unsigned int i) const;
976 :
977 : /**
978 : * This array maps the integer representation of the \p ElemType enum
979 : * to the default approximation order of elements of that type.
980 : *
981 : * This is currently usable even for complicated subclasses with
982 : * runtime-varying topology.
983 : */
984 : static const Order type_to_default_order_map[INVALID_ELEM];
985 :
986 : /**
987 : * \returns The default approximation order for this element. This
988 : * is the order that will be used to compute the map to the
989 : * reference element.
990 : */
991 : virtual Order default_order () const = 0;
992 :
993 : /**
994 : * \returns The maximum supported approximation order for nodal
995 : * (Lagrange or Rational Bezier-Bernstein) variables on this element
996 : * type. This is usually the same as the default order.
997 : */
998 38858696 : virtual Order supported_nodal_order() const { return default_order(); }
999 :
1000 : /**
1001 : * \returns The default approximation order for side elements of
1002 : * this element type. This may be lower for elements with 'bubble
1003 : * functions' in the Lagrange basis.
1004 : */
1005 481660 : virtual Order default_side_order () const { return default_order(); }
1006 :
1007 : /**
1008 : * \returns The "true" geometric centroid of the element, c=(cx, cy,
1009 : * cz), where:
1010 : *
1011 : * [cx] [\int x dV]
1012 : * [cy] := (1/V) * [\int y dV]
1013 : * [cz] [\int z dV]
1014 : *
1015 : * This method is virtual since some derived elements might want to
1016 : * use shortcuts to compute their centroid. For most element types,
1017 : * this method is more expensive than calling vertex_average(), so
1018 : * if you only need a point which is located "somewhere" in the
1019 : * interior of the element, consider calling vertex_average() instead.
1020 : */
1021 : virtual Point true_centroid () const;
1022 :
1023 : /**
1024 : * \returns A Point at the average of the elment's vertices.
1025 : *
1026 : * \note This used to be the base class centroid() implementation, but
1027 : * the centroid is only equal to the vertex average in some special cases.
1028 : * The centroid() implementation now returns the "true" centroid of the
1029 : * element (up to quadrature error).
1030 : */
1031 : Point vertex_average () const;
1032 :
1033 : /**
1034 : * \returns The "circumcenter of mass" (area-weighted average of
1035 : * triangulation circumcenters) of the element.
1036 : *
1037 : * Not implemented for infinite elements, not currently implemented
1038 : * for 3D elements, currently ignores curvature of element edges.
1039 : */
1040 0 : virtual Point quasicircumcenter () const
1041 0 : { libmesh_not_implemented(); }
1042 :
1043 : /**
1044 : * \returns The minimum vertex separation for the element.
1045 : */
1046 : virtual Real hmin () const;
1047 :
1048 : /**
1049 : * \returns The maximum vertex separation for the element.
1050 : */
1051 : virtual Real hmax () const;
1052 :
1053 : /**
1054 : * \returns The (length/area/volume) of the geometric element.
1055 : *
1056 : * If the element is twisted or inverted such that the mapping
1057 : * Jacobian is singular at any point, implementations of this method
1058 : * may return a "net" volume or may simply return NaN.
1059 : */
1060 : virtual Real volume () const;
1061 :
1062 : /**
1063 : * \returns A bounding box (not necessarily the minimal bounding box)
1064 : * containing the geometric element.
1065 : *
1066 : * The base class implementation determines a bounding box for the
1067 : * element *nodes*, which should be sufficient for first order
1068 : * finite elements. Higher order geometric elements will need to
1069 : * override with an implementation which takes curved elements into
1070 : * account.
1071 : */
1072 : virtual BoundingBox loose_bounding_box () const;
1073 :
1074 : /**
1075 : * \returns A quantitative assessment of element quality based on
1076 : * the quality metric \p q specified by the user. Not all ElemQuality
1077 : * metrics are supported for all Elem types; consult the Elem::quality()
1078 : * overrides for specific Elem types to determine which quality metrics
1079 : * are supported. The ElemQuality metrics with generic support for all
1080 : * Elems with dimension > 1 are:
1081 : * .) EDGE_LENGTH_RATIO - ratio of maximum to minimum edge (in 2D,
1082 : * side) length, where the min/max is taken over all Elem edges.
1083 : * .) MIN,MAX_ANGLE - The minimum (respectively maximum) angle
1084 : * between all pairs of adjacent Elem edges, in degrees. Note
1085 : * that, in 3D, these are *not* the dihedral angles (angle
1086 : * between adjacent planar faces of the element), which we plan to
1087 : * add support for in the future. In 2D, we compute the angle
1088 : * between adjacent sides for this metric.
1089 : */
1090 : virtual Real quality (const ElemQuality q) const;
1091 :
1092 : /**
1093 : * \returns The suggested quality bounds for the Elem based on
1094 : * quality measure \p q.
1095 : *
1096 : * These are the values suggested by the CUBIT User's Manual. Since
1097 : * this function can have no possible meaning for an abstract Elem,
1098 : * it is an error in the base class.
1099 : */
1100 0 : virtual std::pair<Real,Real> qual_bounds (const ElemQuality) const
1101 0 : { libmesh_not_implemented(); return std::make_pair(0.,0.); }
1102 :
1103 : /**
1104 : * \returns \p true if the physical point p is contained in this
1105 : * element, false otherwise.
1106 : *
1107 : * For linear elements, performs an initial tight bounding box check
1108 : * (as an optimization step) and (if that passes) then uses the
1109 : * user-defined tolerance "tol" in a call to inverse_map() to actually
1110 : * test if the point is in the element. For quadratic elements, the
1111 : * bounding box optimization is skipped, and only the inverse_map()
1112 : * steps are performed.
1113 : *
1114 : * \note This routine should not be used to determine if a point
1115 : * is merely "nearby" an element to within some tolerance. For that,
1116 : * use Elem::close_to_point() instead.
1117 : */
1118 : virtual bool contains_point (const Point & p, Real tol=TOLERANCE) const;
1119 :
1120 : /**
1121 : * \returns \p true if the master-space point p is contained in the
1122 : * reference element corresponding to this element, false otherwise.
1123 : *
1124 : * Since we are doing floating point comparisons here the parameter
1125 : * \p eps can be specified to indicate a tolerance. For example,
1126 : * \f$ x \le 1 \f$ becomes \f$ x \le 1 + \epsilon \f$.
1127 : */
1128 : virtual bool on_reference_element(const Point & p,
1129 : const Real eps = TOLERANCE) const = 0;
1130 :
1131 : /**
1132 : * \returns \p true if this element is "close" to the point p, where
1133 : * "close" is determined by the tolerance tol.
1134 : */
1135 : virtual bool close_to_point(const Point & p, Real tol) const;
1136 :
1137 : /**
1138 : * \returns \p true if edge \p i is positively oriented. An edge is
1139 : * positively oriented iff its first vertex (i.e. zeroth node) is
1140 : * lexicographically greater than its second vertex (i.e. first node).
1141 : */
1142 : bool positive_edge_orientation(const unsigned int i) const;
1143 :
1144 : /**
1145 : * \returns \p true if face \p i is positively oriented. A face is
1146 : * positively oriented iff the triangle defined by the lexicographically
1147 : * least vertex and its two adjacent vertices on the same face is
1148 : * positively oriented. Said triangle is positively oriented iff its
1149 : * vertices are an odd permutation of their lexicographic ordering.
1150 : */
1151 : bool positive_face_orientation(const unsigned int i) const;
1152 :
1153 :
1154 : /**
1155 : * A helper function for copying generic element data (mapping,
1156 : * subdomain, processor) from an element to a derived (child, side,
1157 : * edge) element. Useful for forwards compatibility when new data
1158 : * is added.
1159 : */
1160 : void inherit_data_from(const Elem & src);
1161 :
1162 :
1163 : private:
1164 : /**
1165 : * Shared private implementation used by the contains_point()
1166 : * and close_to_point() routines. The box_tol tolerance is
1167 : * used in the bounding box optimization, the map_tol tolerance is used
1168 : * in the calls to inverse_map() and on_reference_element().
1169 : */
1170 : bool point_test(const Point & p, Real box_tol, Real map_tol) const;
1171 :
1172 : public:
1173 : /**
1174 : * \returns \p true if the element map is definitely affine (i.e. the same at
1175 : * every quadrature point) within numerical tolerances.
1176 : */
1177 0 : virtual bool has_affine_map () const { return false; }
1178 :
1179 : /**
1180 : * \returns \p true if the element map is invertible everywhere on
1181 : * the element, to within a user-specified tolerance. The tolerance
1182 : * is generally used in comparisons against zero, so it should be an
1183 : * absolute rather than a relative tolerance. Throws a
1184 : * libmesh_not_implemented() error unless specialized by derived
1185 : * classes.
1186 : */
1187 : virtual bool has_invertible_map(Real tol = TOLERANCE*TOLERANCE) const;
1188 :
1189 : /**
1190 : * \returns \p true if the Lagrange shape functions on this element
1191 : * are linear.
1192 : */
1193 0 : virtual bool is_linear () const { return false; }
1194 :
1195 : /**
1196 : * Prints relevant information about the element.
1197 : */
1198 : void print_info (std::ostream & os=libMesh::out) const;
1199 :
1200 : /**
1201 : * Prints relevant information about the element to a string.
1202 : */
1203 : std::string get_info () const;
1204 :
1205 : /**
1206 : * \returns \p true if the element is active (i.e. has no active
1207 : * descendants) or AMR is disabled, \p false otherwise.
1208 : *
1209 : * \note It suffices to check the first child only.
1210 : */
1211 : bool active () const;
1212 :
1213 : /**
1214 : * \returns \p true if the element is an ancestor (i.e. has an
1215 : * active child or ancestor child), \p false otherwise or when AMR
1216 : * is disabled.
1217 : */
1218 : bool ancestor () const;
1219 :
1220 : /**
1221 : * \returns \p true if the element is subactive (i.e. has no active
1222 : * descendants), \p false otherwise or if AMR is disabled.
1223 : */
1224 : bool subactive () const;
1225 :
1226 : /**
1227 : * \returns \p true if the element has any children (active or not),
1228 : * \p false otherwise, or if AMR is disabled.
1229 : */
1230 : bool has_children () const;
1231 :
1232 : /**
1233 : * \returns \p true if the element has any descendants other than
1234 : * its immediate children, \p false otherwise, or if AMR is disabled.
1235 : */
1236 : bool has_ancestor_children () const;
1237 :
1238 : /**
1239 : * \returns \p true if \p descendant is a child of \p this, or a
1240 : * child of a child of \p this, etc., \p false otherwise or if AMR
1241 : * is disabled.
1242 : */
1243 : bool is_ancestor_of(const Elem * descendant) const;
1244 :
1245 : /**
1246 : * \returns A const pointer to the element's parent, or \p nullptr if
1247 : * the element was not created via refinement.
1248 : */
1249 : const Elem * parent () const;
1250 :
1251 : /**
1252 : * \returns A pointer to the element's parent, or \p nullptr if
1253 : * the element was not created via refinement.
1254 : */
1255 : Elem * parent ();
1256 :
1257 : /**
1258 : * Sets the pointer to the element's parent.
1259 : * Dangerous! Only use this if you know what you are doing!
1260 : */
1261 : void set_parent (Elem * p);
1262 :
1263 : /**
1264 : * \returns A pointer to the element's top-most (i.e. level-0) parent.
1265 : *
1266 : * That is, \p this if this is a level-0 element, this element's parent
1267 : * if this is a level-1 element, this element's grandparent if this is
1268 : * a level-2 element, etc...
1269 : */
1270 : const Elem * top_parent () const;
1271 :
1272 : /**
1273 : * \returns The higher-dimensional Elem for which this Elem is a face.
1274 : *
1275 : * In some cases it is desirable to extract the boundary (or a subset thereof)
1276 : * of a D-dimensional mesh as a (D-1)-dimensional manifold. In this case
1277 : * we may want to know the 'parent' element from which the manifold elements
1278 : * were extracted. We can easily do that for the level-0 manifold elements
1279 : * by storing the D-dimensional parent. This method provides access to that
1280 : * element.
1281 : *
1282 : * This method returns nullptr if this->dim() == LIBMESH_DIM; in
1283 : * such cases no data storage for an interior parent pointer has
1284 : * been allocated.
1285 : */
1286 : const Elem * interior_parent () const;
1287 :
1288 : Elem * interior_parent ();
1289 :
1290 : /**
1291 : * Sets the pointer to the element's interior_parent.
1292 : * Dangerous! Only use this if you know what you are doing!
1293 : */
1294 : void set_interior_parent (Elem * p);
1295 :
1296 : /**
1297 : * \returns The distance between nodes n1 and n2.
1298 : *
1299 : * Useful for computing the lengths of the sides of elements.
1300 : */
1301 : Real length (const unsigned int n1,
1302 : const unsigned int n2) const;
1303 :
1304 : /**
1305 : * \returns The number of adjacent vertices that uniquely define the
1306 : * location of the \f$ n^{th} \f$ second-order node, or 0 for linear
1307 : * elements.
1308 : *
1309 : * This method is useful when converting linear elements to quadratic
1310 : * elements.
1311 : *
1312 : * \note \p n has to be greater than or equal to \p this->n_vertices().
1313 : */
1314 : virtual unsigned int n_second_order_adjacent_vertices (const unsigned int n) const;
1315 :
1316 : /**
1317 : * \returns The element-local number of the \f$ v^{th} \f$ vertex
1318 : * that defines the \f$ n^{th} \f$ second-order node, or 0 for
1319 : * linear elements.
1320 : *
1321 : * \note The value is always less than \p this->n_vertices(), while
1322 : * \p n has to be greater than or equal to \p this->n_vertices().
1323 : */
1324 : virtual unsigned short int second_order_adjacent_vertex (const unsigned int n,
1325 : const unsigned int v) const;
1326 :
1327 : /**
1328 : * \returns A pair (c,v), where
1329 : * c == child index, and
1330 : * v == element-local index of the \p \f$ n^{th} \f$
1331 : * second-order node on the parent element.
1332 : * For linear elements, (0,0) is returned.
1333 : *
1334 : * \note The return values are always less than \p this->n_children()
1335 : * and \p this->child_ptr(c)->n_vertices().
1336 : *
1337 : * \note \p n has to be greater than or equal to \p this->n_vertices().
1338 : *
1339 : * \note On refined second-order elements, the return value will
1340 : * satisfy \p this->node_ptr(n) == this->child_ptr(c)->node_ptr(v).
1341 : */
1342 : virtual std::pair<unsigned short int, unsigned short int>
1343 : second_order_child_vertex (const unsigned int n) const;
1344 :
1345 : /**
1346 : * \returns The ElemType of the associated second-order element
1347 : * (which will be the same as the input if the input is already a
1348 : * second-order ElemType) or INVALID_ELEM for elements that cannot be
1349 : * converted into higher order equivalents.
1350 : *
1351 : * For example, when \p this is a \p TET4, then \p TET10 is returned.
1352 : *
1353 : * For some elements, there exist two second-order equivalents, e.g.
1354 : * for \p Quad4 there is \p Quad8 and \p Quad9. When the optional
1355 : * \p full_ordered is \p true, then \p QUAD9 is returned. When
1356 : * \p full_ordered is \p false, then \p QUAD8 is returned.
1357 : */
1358 : static ElemType second_order_equivalent_type (const ElemType et,
1359 : const bool full_ordered=true);
1360 :
1361 : /**
1362 : * \returns The element type of the associated first-order element,
1363 : * or \p INVALID_ELEM for first-order or other elements that cannot be
1364 : * converted into lower order equivalents.
1365 : *
1366 : * For example, when \p this is a \p TET10, then \p TET4 is returned.
1367 : */
1368 : static ElemType first_order_equivalent_type (const ElemType et);
1369 :
1370 : /**
1371 : * \returns The ElemType of the associated "complete" order element
1372 : * (which will be the same as the input if the input is already a
1373 : * complete-order ElemType), or INVALID_ELEM for elements that cannot be
1374 : * converted into complete-order equivalents.
1375 : *
1376 : * The "complete" version of an element is an element which can
1377 : * represent the same geometry but which has nodes available to
1378 : * restore degrees of freedom on any vertex, edge, or face.
1379 : *
1380 : * For example, when \p this is a \p TET4, then \p TET14 is returned.
1381 : */
1382 : static ElemType complete_order_equivalent_type (const ElemType et);
1383 :
1384 : /**
1385 : * \returns The refinement level of the current element.
1386 : *
1387 : * If the element's parent is \p nullptr then by convention it is at
1388 : * level 0, otherwise it is simply at one level greater than its
1389 : * parent.
1390 : */
1391 : unsigned int level () const;
1392 :
1393 : /**
1394 : * \returns The value of the p refinement level of an active
1395 : * element, or the minimum value of the p refinement levels
1396 : * of an ancestor element's descendants.
1397 : */
1398 : unsigned int p_level () const;
1399 :
1400 : /**
1401 : * \returns \p true if the specified child is on the specified side.
1402 : */
1403 : virtual bool is_child_on_side(const unsigned int c,
1404 : const unsigned int s) const = 0;
1405 :
1406 : /**
1407 : * \returns The value of the mapping type for the element.
1408 : */
1409 : ElemMappingType mapping_type () const;
1410 :
1411 : /**
1412 : * Sets the value of the mapping type for the element.
1413 : */
1414 : void set_mapping_type (const ElemMappingType type);
1415 :
1416 : /**
1417 : * \returns The value of the mapping data for the element.
1418 : */
1419 : unsigned char mapping_data () const;
1420 :
1421 : /**
1422 : * Sets the value of the mapping data for the element.
1423 : */
1424 : void set_mapping_data (const unsigned char data);
1425 :
1426 :
1427 : #ifdef LIBMESH_ENABLE_AMR
1428 :
1429 : /**
1430 : * Enumeration of possible element refinement states.
1431 : */
1432 : enum RefinementState { COARSEN = 0,
1433 : DO_NOTHING,
1434 : REFINE,
1435 : JUST_REFINED,
1436 : JUST_COARSENED,
1437 : INACTIVE,
1438 : COARSEN_INACTIVE,
1439 : INVALID_REFINEMENTSTATE };
1440 :
1441 : /**
1442 : * \returns A constant pointer to the \f$ i^{th} \f$ child for this element.
1443 : * For internal use only - skips assertions about null pointers.
1444 : */
1445 : const Elem * raw_child_ptr (unsigned int i) const;
1446 :
1447 : /**
1448 : * \returns A constant pointer to the \f$ i^{th} \f$ child for this element.
1449 : * Do not call if this element has no children, i.e. is active.
1450 : */
1451 : const Elem * child_ptr (unsigned int i) const;
1452 :
1453 : /**
1454 : * \returns A non-constant pointer to the \f$ i^{th} \f$ child for this element.
1455 : * Do not call if this element has no children, i.e. is active.
1456 : */
1457 : Elem * child_ptr (unsigned int i);
1458 :
1459 : /**
1460 : * Nested classes for use iterating over all children of a parent
1461 : * element.
1462 : */
1463 : class ChildRefIter;
1464 : class ConstChildRefIter;
1465 :
1466 : /**
1467 : * Returns a range with all children of a parent element, usable in
1468 : * range-based for loops. The exact type of the return value here
1469 : * may be subject to change in future libMesh releases, but the
1470 : * iterators will always dereference to produce a reference to a
1471 : * child element.
1472 : */
1473 : SimpleRange<ChildRefIter> child_ref_range();
1474 :
1475 : SimpleRange<ConstChildRefIter> child_ref_range() const;
1476 :
1477 : private:
1478 : /**
1479 : * Sets the pointer to the \f$ i^{th} \f$ child for this element.
1480 : * Do not call if this element has no children, i.e. is active.
1481 : */
1482 : void set_child (unsigned int c, Elem * elem);
1483 :
1484 : public:
1485 : /**
1486 : * \returns The child index which \p e corresponds to.
1487 : *
1488 : * I.e. if c = a->which_child_am_i(e); then a->child_ptr(c) will be
1489 : * e.
1490 : */
1491 : unsigned int which_child_am_i(const Elem * e) const;
1492 :
1493 : /**
1494 : * \returns \p true if the specified child is on the specified edge.
1495 : */
1496 : virtual bool is_child_on_edge(const unsigned int c,
1497 : const unsigned int e) const;
1498 :
1499 : /**
1500 : * Adds a child pointer to the array of children of this element.
1501 : * If this is the first child to be added, this method allocates
1502 : * memory in the parent's _children array, otherwise, it just sets
1503 : * the pointer.
1504 : */
1505 : void add_child (Elem * elem);
1506 :
1507 : /**
1508 : * Adds a new child pointer to the specified index in the array of
1509 : * children of this element. If this is the first child to be added,
1510 : * this method allocates memory in the parent's _children array,
1511 : * otherwise, it just sets the pointer.
1512 : */
1513 : void add_child (Elem * elem, unsigned int c);
1514 :
1515 : /**
1516 : * Replaces the child pointer at the specified index in the child array.
1517 : */
1518 : void replace_child (Elem * elem, unsigned int c);
1519 :
1520 : /**
1521 : * Fills the vector \p family with the children of this element,
1522 : * recursively. Calling this method on a twice-refined element
1523 : * will give you the element itself, its direct children, and their
1524 : * children, etc... When the optional parameter \p reset is
1525 : * true, the vector will be cleared before the element and its
1526 : * descendants are added.
1527 : *
1528 : * The family tree only includes ancestor and active elements. To
1529 : * include subactive elements as well, use total_family_tree().
1530 : */
1531 : void family_tree (std::vector<const Elem *> & family,
1532 : bool reset = true) const;
1533 :
1534 : /**
1535 : * Non-const version of function above; fills a vector of non-const pointers.
1536 : */
1537 : void family_tree (std::vector<Elem *> & family,
1538 : bool reset = true);
1539 :
1540 : /**
1541 : * Same as the \p family_tree() member, but also adds any subactive
1542 : * descendants.
1543 : */
1544 : void total_family_tree (std::vector<const Elem *> & family,
1545 : bool reset = true) const;
1546 :
1547 : /**
1548 : * Non-const version of function above; fills a vector of non-const pointers.
1549 : */
1550 : void total_family_tree (std::vector<Elem *> & family,
1551 : bool reset = true);
1552 :
1553 : /**
1554 : * Same as the \p family_tree() member, but only adds the active
1555 : * children. Can be thought of as removing all the inactive
1556 : * elements from the vector created by \p family_tree, but is
1557 : * implemented more efficiently.
1558 : */
1559 : void active_family_tree (std::vector<const Elem *> & active_family,
1560 : bool reset = true) const;
1561 :
1562 : /**
1563 : * Non-const version of function above; fills a vector of non-const pointers.
1564 : */
1565 : void active_family_tree (std::vector<Elem *> & active_family,
1566 : bool reset = true);
1567 :
1568 : /**
1569 : * Same as the \p family_tree() member, but only adds elements
1570 : * which are next to \p side.
1571 : */
1572 : void family_tree_by_side (std::vector<const Elem *> & family,
1573 : unsigned int side,
1574 : bool reset = true) const;
1575 :
1576 : /**
1577 : * Non-const version of function above; fills a vector of non-const pointers.
1578 : */
1579 : void family_tree_by_side (std::vector<Elem *> & family,
1580 : unsigned int side,
1581 : bool reset = true);
1582 :
1583 : /**
1584 : * Same as the \p active_family_tree() member, but only adds elements
1585 : * which are next to \p side.
1586 : */
1587 : void active_family_tree_by_side (std::vector<const Elem *> & family,
1588 : unsigned int side,
1589 : bool reset = true) const;
1590 :
1591 : /**
1592 : * Non-const version of function above; fills a vector of non-const pointers.
1593 : */
1594 : void active_family_tree_by_side (std::vector<Elem *> & family,
1595 : unsigned int side,
1596 : bool reset = true);
1597 :
1598 : /**
1599 : * Same as the \p family_tree() member, but only adds elements
1600 : * which are next to \p neighbor.
1601 : */
1602 : void family_tree_by_neighbor (std::vector<const Elem *> & family,
1603 : const Elem * neighbor,
1604 : bool reset = true) const;
1605 :
1606 : /**
1607 : * Non-const version of function above; fills a vector of non-const pointers.
1608 : */
1609 : void family_tree_by_neighbor (std::vector<Elem *> & family,
1610 : Elem * neighbor,
1611 : bool reset = true);
1612 :
1613 : /**
1614 : * Same as the \p family_tree_by_neighbor() member, but also adds
1615 : * any subactive descendants.
1616 : */
1617 : void total_family_tree_by_neighbor (std::vector<const Elem *> & family,
1618 : const Elem * neighbor,
1619 : bool reset = true) const;
1620 :
1621 : /**
1622 : * Non-const version of function above; fills a vector of non-const pointers.
1623 : */
1624 : void total_family_tree_by_neighbor (std::vector<Elem *> & family,
1625 : Elem * neighbor,
1626 : bool reset = true);
1627 :
1628 : /**
1629 : * Same as the \p family_tree() member, but only adds elements
1630 : * which are next to \p subneighbor. Only applicable when
1631 : * \p this->has_neighbor(neighbor) and
1632 : * \p neighbor->is_ancestor(subneighbor)
1633 : */
1634 : void family_tree_by_subneighbor (std::vector<const Elem *> & family,
1635 : const Elem * neighbor,
1636 : const Elem * subneighbor,
1637 : bool reset = true) const;
1638 :
1639 : /**
1640 : * Non-const version of function above; fills a vector of non-const pointers.
1641 : */
1642 : void family_tree_by_subneighbor (std::vector<Elem *> & family,
1643 : Elem * neighbor,
1644 : Elem * subneighbor,
1645 : bool reset = true);
1646 :
1647 : /**
1648 : * Same as the \p family_tree_by_subneighbor() member, but also adds
1649 : * any subactive descendants.
1650 : */
1651 : void total_family_tree_by_subneighbor (std::vector<const Elem *> & family,
1652 : const Elem * neighbor,
1653 : const Elem * subneighbor,
1654 : bool reset = true) const;
1655 :
1656 : /**
1657 : * Non-const version of function above; fills a vector of non-const pointers.
1658 : */
1659 : void total_family_tree_by_subneighbor (std::vector<Elem *> & family,
1660 : Elem * neighbor,
1661 : Elem * subneighbor,
1662 : bool reset = true);
1663 :
1664 : /**
1665 : * Same as the \p active_family_tree() member, but only adds elements
1666 : * which are next to \p neighbor.
1667 : */
1668 : void active_family_tree_by_neighbor (std::vector<const Elem *> & family,
1669 : const Elem * neighbor,
1670 : bool reset = true) const;
1671 :
1672 : /**
1673 : * Non-const version of function above; fills a vector of non-const pointers.
1674 : */
1675 : void active_family_tree_by_neighbor (std::vector<Elem *> & family,
1676 : Elem * neighbor,
1677 : bool reset = true);
1678 :
1679 : /**
1680 : * Same as the \p active_family_tree_by_neighbor() member, but the
1681 : * \p neighbor here may be a topological (e.g. periodic boundary
1682 : * condition) neighbor, not just a local neighbor.
1683 : */
1684 : void active_family_tree_by_topological_neighbor (std::vector<const Elem *> & family,
1685 : const Elem * neighbor,
1686 : const MeshBase & mesh,
1687 : const PointLocatorBase & point_locator,
1688 : const PeriodicBoundaries * pb,
1689 : bool reset = true) const;
1690 :
1691 : /**
1692 : * Non-const version of function above; fills a vector of non-const pointers.
1693 : */
1694 : void active_family_tree_by_topological_neighbor (std::vector<Elem *> & family,
1695 : Elem * neighbor,
1696 : const MeshBase & mesh,
1697 : const PointLocatorBase & point_locator,
1698 : const PeriodicBoundaries * pb,
1699 : bool reset = true);
1700 :
1701 : /**
1702 : * \returns The value of the refinement flag for the element.
1703 : */
1704 : RefinementState refinement_flag () const;
1705 :
1706 : /**
1707 : * Sets the value of the refinement flag for the element.
1708 : */
1709 : void set_refinement_flag (const RefinementState rflag);
1710 :
1711 : /**
1712 : * \returns The value of the p-refinement flag for the element.
1713 : */
1714 : RefinementState p_refinement_flag () const;
1715 :
1716 : /**
1717 : * Sets the value of the p-refinement flag for the element.
1718 : */
1719 : void set_p_refinement_flag (const RefinementState pflag);
1720 :
1721 : /**
1722 : * \returns The maximum value of the p-refinement levels of
1723 : * an ancestor element's descendants.
1724 : */
1725 : unsigned int max_descendant_p_level () const;
1726 :
1727 : /**
1728 : * \returns The minimum p-refinement level of elements which are
1729 : * descended from this element, and which share a side with the
1730 : * active \p neighbor.
1731 : */
1732 : unsigned int min_p_level_by_neighbor (const Elem * neighbor,
1733 : unsigned int current_min) const;
1734 :
1735 : /**
1736 : * \returns The minimum new p-refinement level (i.e. after refinement
1737 : * and coarsening is done) of elements which are descended from this
1738 : * element and which share a side with the active \p neighbor.
1739 : */
1740 : unsigned int min_new_p_level_by_neighbor (const Elem * neighbor,
1741 : unsigned int current_min) const;
1742 :
1743 : /**
1744 : * Sets the value of the p-refinement level for the element.
1745 : *
1746 : * \note The maximum p-refinement level is currently 255.
1747 : */
1748 : void set_p_level (const unsigned int p);
1749 :
1750 : /**
1751 : * Sets the value of the p-refinement level for the element
1752 : * without altering the p-level of its ancestors
1753 : */
1754 : void hack_p_level (const unsigned int p);
1755 :
1756 : /**
1757 : * Sets the value of the p-refinement level for the element
1758 : * without altering the p-level of its ancestors; also sets the
1759 : * p_refinement_flag, simultaneously so that they can be safely
1760 : * checked for mutual consistency
1761 : */
1762 : void hack_p_level_and_refinement_flag (const unsigned int p,
1763 : RefinementState pflag);
1764 :
1765 : /**
1766 : * Refine the element.
1767 : */
1768 : virtual void refine (MeshRefinement & mesh_refinement);
1769 :
1770 : /**
1771 : * Coarsen the element. This function is non-virtual since it is the same
1772 : * for all element types.
1773 : */
1774 : void coarsen ();
1775 :
1776 : /**
1777 : * Contract an active element, i.e. remove pointers to any
1778 : * subactive children. This should only be called via
1779 : * MeshRefinement::contract, which will also remove subactive
1780 : * children from the mesh.
1781 : */
1782 : void contract ();
1783 :
1784 : #endif
1785 :
1786 : #ifndef NDEBUG
1787 : /**
1788 : * Checks for consistent neighbor links on this element.
1789 : */
1790 : void libmesh_assert_valid_neighbors() const;
1791 :
1792 : /**
1793 : * Checks for a valid id and pointers to nodes with valid ids on
1794 : * this element.
1795 : */
1796 : void libmesh_assert_valid_node_pointers() const;
1797 : #endif // !NDEBUG
1798 :
1799 : /**
1800 : * \returns The local node index of the given point IF said node
1801 : * has a singular Jacobian for this element. If the given point
1802 : * is not a node or is a node and does not have a singular Jacobian,
1803 : * this will return invalid_uint.
1804 : *
1805 : * The intention is for this to be overridden in derived element
1806 : * classes that do have nodes that have singular Jacobians. When
1807 : * mapping failures are caught, we can check this to see if the
1808 : * failed physical point is actually a singular point and
1809 : * return the correct master point.
1810 : */
1811 0 : virtual unsigned int local_singular_node(const Point & /* p */, const Real /* tol */ = TOLERANCE*TOLERANCE) const
1812 0 : { return invalid_uint; }
1813 :
1814 : /**
1815 : * \returns true iff the node at the given index has a singular
1816 : * mapping; i.e. is the degree-4 node on a Pyramid.
1817 : */
1818 0 : virtual bool is_singular_node(unsigned int /* node_i */) const { return false; }
1819 :
1820 : /**
1821 : * \returns The local index of the center node on the side \p side.
1822 : *
1823 : * A center node is a node that is located at the centroid of the given side.
1824 : * If the given side does not have a center node, this will return invalid_uint.
1825 : */
1826 : virtual unsigned int center_node_on_side(const unsigned short side) const;
1827 :
1828 : protected:
1829 :
1830 : /**
1831 : * The protected nested SideIter class is used to iterate over the
1832 : * sides of this Elem. It is a specially-designed class since
1833 : * no sides are actually stored by the element. This iterator-like
1834 : * class has to provide the following three operations
1835 : * 1) operator*
1836 : * 2) operator++
1837 : * 3) operator==
1838 : * The definition can be found at the end of this header file.
1839 : */
1840 : class SideIter;
1841 :
1842 : public:
1843 : /**
1844 : * Useful iterator typedefs
1845 : */
1846 : typedef Predicates::multi_predicate Predicate;
1847 :
1848 : /**
1849 : * Data structure for iterating over sides. Defined at the end of
1850 : * this header file.
1851 : */
1852 : struct side_iterator;
1853 :
1854 : /**
1855 : * Iterator accessor functions
1856 : */
1857 : side_iterator boundary_sides_begin();
1858 : side_iterator boundary_sides_end();
1859 :
1860 : private:
1861 : /**
1862 : * Side iterator helper functions. Used to replace the begin()
1863 : * and end() functions of the STL containers.
1864 : */
1865 : SideIter _first_side();
1866 : SideIter _last_side();
1867 :
1868 : public:
1869 :
1870 : #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1871 :
1872 : /**
1873 : * \returns \p true if the element is an infinite element,
1874 : * \p false otherwise.
1875 : */
1876 : virtual bool infinite () const = 0;
1877 :
1878 : /**
1879 : * \returns \p true if the specified (local) node number is a
1880 : * "mid-edge" node on an infinite element edge.
1881 : *
1882 : * This is false for all nodes on non-infinite elements, so we won't
1883 : * make it pure virtual, to simplify their code.
1884 : */
1885 0 : virtual bool is_mid_infinite_edge_node(const unsigned int /* n */) const
1886 0 : { libmesh_assert (!this->infinite()); return false; }
1887 :
1888 : /**
1889 : * \returns The origin for an infinite element.
1890 : *
1891 : * Currently, all infinite elements used in a mesh share the same
1892 : * origin. Override this in infinite element classes.
1893 : */
1894 0 : virtual Point origin () const { libmesh_not_implemented(); return Point(); }
1895 :
1896 : #else
1897 :
1898 : static constexpr bool infinite () { return false; }
1899 :
1900 : #endif
1901 :
1902 : /**
1903 : * \returns An Elem of type \p type wrapped in a smart pointer.
1904 : */
1905 : static std::unique_ptr<Elem> build (const ElemType type,
1906 : Elem * p=nullptr);
1907 :
1908 : /**
1909 : * Calls the build() method above with a nullptr parent, and
1910 : * additionally sets the newly-created Elem's id. This can be useful
1911 : * when adding pre-numbered Elems to a Mesh via add_elem() calls.
1912 : */
1913 : static std::unique_ptr<Elem> build_with_id (const ElemType type,
1914 : dof_id_type id);
1915 :
1916 : /**
1917 : * \returns An Elem of the same type as \p this, wrapped in a smart
1918 : * pointer.
1919 : *
1920 : * This is not a complete clone() method (since e.g. it does not set
1921 : * node pointers; the standard use case reassigns node pointers from
1922 : * a different mesh), but it is necessary to use this instead of
1923 : * build() for runtime-polymorphic elements like Polygon subtypes
1924 : * whose "type" depends on more than their type(), and it is useful
1925 : * to use this for elements whose id, unique_id, extra integers,
1926 : * etc. should be preserved in the near-clone.
1927 : */
1928 : virtual std::unique_ptr<Elem> disconnected_clone () const;
1929 :
1930 : /**
1931 : * Returns the number of independent permutations of element nodes -
1932 : * e.g. a cube can be reoriented to put side 0 where side N is (for
1933 : * 0 <= N < 6) and then rotated in one of four ways, giving 24
1934 : * possible permutations.
1935 : *
1936 : * Permutations which change the mapping Jacobian of an element
1937 : * (i.e. flipping the element) are not allowed in this definition.
1938 : */
1939 : virtual unsigned int n_permutations() const = 0;
1940 :
1941 : /**
1942 : * Permutes the element (by swapping node and neighbor pointers)
1943 : * according to the specified index.
1944 : *
1945 : * This is useful for regression testing, by making it easy to make
1946 : * a structured mesh behave more like an arbitrarily unstructured
1947 : * mesh.
1948 : *
1949 : * This is so far *only* used for regression testing, so we do
1950 : * not currently provide a way to permute any boundary side/edge ids
1951 : * along with the element permutation.
1952 : */
1953 : virtual void permute(unsigned int perm_num) = 0;
1954 :
1955 : /**
1956 : * Flips the element (by swapping node and neighbor pointers) to
1957 : * have a mapping Jacobian of opposite sign.
1958 : *
1959 : * This is useful for automatically fixing up elements that have
1960 : * been newly created (e.g. from extrusions) with a negative
1961 : * Jacobian.
1962 : *
1963 : * If \p boundary_info is not null, swap boundary side/edge ids
1964 : * consistently.
1965 : */
1966 : virtual void flip(BoundaryInfo * boundary_info) = 0;
1967 :
1968 : /**
1969 : * \returns Whether the element is flipped compared to standard
1970 : * libMesh (e.g. clockwise for 2D elements) node orientations.
1971 : *
1972 : * Always returns \p false if a 2D element is not in the XY plane or
1973 : * a 1D element is not on the X axis; user code designed to work for
1974 : * embedded manifolds should handle any consistent orientation, and
1975 : * determining whether an orientation is consistent is not a local
1976 : * operation.
1977 : */
1978 : virtual bool is_flipped() const = 0;
1979 :
1980 : /**
1981 : * Flips the element (by swapping node and neighbor pointers) to
1982 : * have a mapping Jacobian of opposite sign, iff we find a negative
1983 : * orientation. This only fixes flipped elements; for tangled
1984 : * elements the only fixes possible are non-local.
1985 : */
1986 : void orient(BoundaryInfo * boundary_info);
1987 :
1988 : #ifdef LIBMESH_ENABLE_AMR
1989 :
1990 : /**
1991 : * \returns The local node id on the parent which corresponds to node
1992 : * \p n of child \p c, or \p invalid_uint if no such parent
1993 : * node exists.
1994 : */
1995 : virtual unsigned int as_parent_node (unsigned int c,
1996 : unsigned int n) const;
1997 :
1998 : /**
1999 : * \returns All the pairs of nodes (indexed by local node id) which
2000 : * should bracket node \p n of child \p c.
2001 : */
2002 : virtual
2003 : const std::vector<std::pair<unsigned char, unsigned char>> &
2004 : parent_bracketing_nodes(unsigned int c,
2005 : unsigned int n) const;
2006 :
2007 : /**
2008 : * \returns All the pairs of nodes (indexed by global node id) which
2009 : * should bracket node \p n of child \p c.
2010 : */
2011 : virtual
2012 : const std::vector<std::pair<dof_id_type, dof_id_type>>
2013 : bracketing_nodes(unsigned int c,
2014 : unsigned int n) const;
2015 :
2016 :
2017 : /**
2018 : * \returns The embedding matrix entry for the requested child.
2019 : */
2020 : virtual Real embedding_matrix (const unsigned int child_num,
2021 : const unsigned int child_node_num,
2022 : const unsigned int parent_node_num) const = 0;
2023 :
2024 : /**
2025 : * \returns A "version number" that identifies which embedding
2026 : * matrix is in use.
2027 : *
2028 : * Some element types may use a different embedding matrix depending
2029 : * on their geometric characteristics.
2030 : */
2031 0 : virtual unsigned int embedding_matrix_version () const { return 0; }
2032 :
2033 : #endif // LIBMESH_ENABLE_AMR
2034 :
2035 :
2036 : protected:
2037 :
2038 : /**
2039 : * Default tolerance to use in has_affine_map().
2040 : */
2041 : static constexpr Real affine_tol = TOLERANCE*TOLERANCE;
2042 :
2043 : /**
2044 : * \returns A hash key computed from a single node id.
2045 : */
2046 : static dof_id_type compute_key (dof_id_type n0);
2047 :
2048 : /**
2049 : * \returns A hash key computed from two node ids.
2050 : */
2051 : static dof_id_type compute_key (dof_id_type n0,
2052 : dof_id_type n1);
2053 :
2054 : /**
2055 : * \returns A hash key computed from three node ids.
2056 : */
2057 : static dof_id_type compute_key (dof_id_type n0,
2058 : dof_id_type n1,
2059 : dof_id_type n2);
2060 :
2061 : /**
2062 : * \returns A hash key computed from four node ids.
2063 : */
2064 : static dof_id_type compute_key (dof_id_type n0,
2065 : dof_id_type n1,
2066 : dof_id_type n2,
2067 : dof_id_type n3);
2068 :
2069 : /**
2070 : * Swaps two node_ptrs
2071 : */
2072 28794629 : void swap2nodes(unsigned int n1, unsigned int n2)
2073 : {
2074 4334448 : Node * temp = this->node_ptr(n1);
2075 30961853 : this->set_node(n1, this->node_ptr(n2));
2076 28794629 : this->set_node(n2, temp);
2077 28794629 : }
2078 :
2079 : /**
2080 : * Swaps two neighbor_ptrs
2081 : */
2082 6786404 : void swap2neighbors(unsigned int n1, unsigned int n2)
2083 : {
2084 793040 : Elem * temp = this->neighbor_ptr(n1);
2085 548184 : this->set_neighbor(n1, this->neighbor_ptr(n2));
2086 548184 : this->set_neighbor(n2, temp);
2087 6861068 : }
2088 :
2089 : /**
2090 : * Swaps two sides in \p boundary_info, if it is non-null.
2091 : */
2092 : void swap2boundarysides(unsigned short s1, unsigned short s2,
2093 : BoundaryInfo * boundary_info) const;
2094 :
2095 : /**
2096 : * Swaps two edges in \p boundary_info, if it is non-null.
2097 : */
2098 : void swap2boundaryedges(unsigned short e1, unsigned short e2,
2099 : BoundaryInfo * boundary_info) const;
2100 :
2101 : /**
2102 : * Swaps three node_ptrs, "rotating" them.
2103 : */
2104 11312324 : void swap3nodes(unsigned int n1, unsigned int n2, unsigned int n3)
2105 : {
2106 12189036 : swap2nodes(n1, n2);
2107 12189036 : swap2nodes(n2, n3);
2108 11312324 : }
2109 :
2110 : /**
2111 : * Swaps three neighbor_ptrs, "rotating" them.
2112 : */
2113 3071806 : void swap3neighbors(unsigned int n1, unsigned int n2,
2114 : unsigned int n3)
2115 : {
2116 3071806 : swap2neighbors(n1, n2);
2117 3071806 : swap2neighbors(n2, n3);
2118 3071806 : }
2119 :
2120 : /**
2121 : * Swaps four node_ptrs, "rotating" them.
2122 : */
2123 3278918 : void swap4nodes(unsigned int n1, unsigned int n2, unsigned int n3,
2124 : unsigned int n4)
2125 : {
2126 2952118 : swap3nodes(n1, n2, n3);
2127 3278918 : swap2nodes(n3, n4);
2128 3278918 : }
2129 :
2130 : /**
2131 : * Swaps four neighbor_ptrs, "rotating" them.
2132 : */
2133 618012 : void swap4neighbors(unsigned int n1, unsigned int n2,
2134 : unsigned int n3, unsigned int n4)
2135 : {
2136 618012 : swap3neighbors(n1, n2, n3);
2137 618012 : swap2neighbors(n3, n4);
2138 618012 : }
2139 :
2140 :
2141 : /**
2142 : * An implementation for simple (all sides equal) elements
2143 : */
2144 : template <typename Sideclass, typename Subclass>
2145 : std::unique_ptr<Elem>
2146 : simple_build_side_ptr(const unsigned int i);
2147 :
2148 : /**
2149 : * An implementation for simple (all sides equal) elements
2150 : */
2151 : template <typename Subclass>
2152 : void simple_build_side_ptr(std::unique_ptr<Elem> & side,
2153 : const unsigned int i,
2154 : ElemType sidetype);
2155 :
2156 : /**
2157 : * An implementation for simple (all sides equal) elements
2158 : */
2159 : template <typename Subclass, typename Mapclass>
2160 : void simple_side_ptr(std::unique_ptr<Elem> & side,
2161 : const unsigned int i,
2162 : ElemType sidetype);
2163 :
2164 : /**
2165 : * An implementation for simple (all edges equal) elements
2166 : */
2167 : template <typename Edgeclass, typename Subclass>
2168 : std::unique_ptr<Elem>
2169 : simple_build_edge_ptr(const unsigned int i);
2170 :
2171 : /**
2172 : * An implementation for simple (all edges equal) elements
2173 : */
2174 : template <typename Subclass>
2175 : void simple_build_edge_ptr(std::unique_ptr<Elem> & edge,
2176 : const unsigned int i,
2177 : ElemType edgetype);
2178 :
2179 :
2180 : #ifdef LIBMESH_ENABLE_AMR
2181 :
2182 : /**
2183 : * Elem subclasses which don't do their own bracketing node
2184 : * calculations will need to supply a static cache, since the
2185 : * default calculation is slow.
2186 : */
2187 : virtual
2188 : std::vector<std::vector<std::vector<std::vector<std::pair<unsigned char, unsigned char>>>>> &
2189 0 : _get_bracketing_node_cache() const
2190 : {
2191 0 : static std::vector<std::vector<std::vector<std::vector<std::pair<unsigned char, unsigned char>>>>> c;
2192 0 : libmesh_error();
2193 : return c;
2194 : }
2195 :
2196 : /**
2197 : * Elem subclasses which don't do their own child-to-parent node
2198 : * calculations will need to supply a static cache, since the
2199 : * default calculation is slow.
2200 : */
2201 : virtual
2202 : std::vector<std::vector<std::vector<signed char>>> &
2203 0 : _get_parent_indices_cache() const
2204 : {
2205 0 : static std::vector<std::vector<std::vector<signed char>>> c;
2206 0 : libmesh_error();
2207 : return c;
2208 : }
2209 :
2210 : #endif // LIBMESH_ENABLE_AMR
2211 :
2212 : public:
2213 :
2214 : /**
2215 : * Replaces this element with \p nullptr for all of its neighbors.
2216 : * This is useful when deleting an element.
2217 : */
2218 : void nullify_neighbors ();
2219 :
2220 : protected:
2221 :
2222 : /**
2223 : * Pointers to the nodes we are connected to.
2224 : */
2225 : Node ** _nodes;
2226 :
2227 : /**
2228 : * Pointers to this element's parent and neighbors, and for
2229 : * lower-dimensional elements' interior_parent.
2230 : */
2231 : Elem ** _elemlinks;
2232 :
2233 : #ifdef LIBMESH_ENABLE_AMR
2234 : /**
2235 : * unique_ptr to array of this element's children.
2236 : *
2237 : * A Mesh ultimately owns the child Elems so we are not responsible
2238 : * for deleting them, but we are responsible for cleaning up the
2239 : * array allocated to hold those Elems, hence the unique_ptr.
2240 : */
2241 : std::unique_ptr<Elem *[]> _children;
2242 : #endif
2243 :
2244 : /**
2245 : * The subdomain to which this element belongs.
2246 : */
2247 : subdomain_id_type _sbd_id;
2248 :
2249 : #ifdef LIBMESH_ENABLE_AMR
2250 : /**
2251 : * h refinement flag. This is stored as an unsigned char
2252 : * to save space.
2253 : */
2254 : unsigned char _rflag;
2255 :
2256 : /**
2257 : * p refinement flag. This is stored as an unsigned char
2258 : * to save space.
2259 : */
2260 : unsigned char _pflag;
2261 :
2262 : /**
2263 : * p refinement level - the difference between the
2264 : * polynomial degree on this element and the minimum
2265 : * polynomial degree on the mesh.
2266 : * This is stored as an unsigned char to save space.
2267 : * In theory, these last four bytes might have
2268 : * been padding anyway.
2269 : */
2270 : unsigned char _p_level;
2271 : #endif
2272 :
2273 : /**
2274 : * Mapping function type; currently either 0 (LAGRANGE) or 1
2275 : * (RATIONAL_BERNSTEIN).
2276 : */
2277 : unsigned char _map_type;
2278 :
2279 : /**
2280 : * Mapping function data; currently used when needed to store the
2281 : * RATIONAL_BERNSTEIN nodal weight data index.
2282 : */
2283 : unsigned char _map_data;
2284 : };
2285 :
2286 :
2287 :
2288 : // ------------------------------------------------------------
2289 : // Elem helper classes
2290 : //
2291 : class
2292 : Elem::NodeRefIter : public PointerToPointerIter<Node>
2293 : {
2294 : public:
2295 22573684 : NodeRefIter (Node * const * nodepp) : PointerToPointerIter<Node>(nodepp) {}
2296 : };
2297 :
2298 :
2299 : class
2300 : Elem::ConstNodeRefIter : public PointerToPointerIter<const Node>
2301 : {
2302 : public:
2303 17714810 : ConstNodeRefIter (const Node * const * nodepp) : PointerToPointerIter<const Node>(nodepp) {}
2304 : };
2305 :
2306 :
2307 : #ifdef LIBMESH_ENABLE_AMR
2308 : class
2309 : Elem::ChildRefIter : public PointerToPointerIter<Elem>
2310 : {
2311 : public:
2312 6437588 : ChildRefIter (Elem * const * childpp) : PointerToPointerIter<Elem>(childpp) {}
2313 : };
2314 :
2315 :
2316 : class
2317 : Elem::ConstChildRefIter : public PointerToPointerIter<const Elem>
2318 : {
2319 : public:
2320 2483114 : ConstChildRefIter (const Elem * const * childpp) : PointerToPointerIter<const Elem>(childpp) {}
2321 : };
2322 :
2323 :
2324 :
2325 : inline
2326 57316574 : SimpleRange<Elem::ChildRefIter> Elem::child_ref_range()
2327 : {
2328 3218794 : libmesh_assert(_children);
2329 59997391 : return {_children.get(), _children.get() + this->n_children()};
2330 : }
2331 :
2332 :
2333 : inline
2334 6674250 : SimpleRange<Elem::ConstChildRefIter> Elem::child_ref_range() const
2335 : {
2336 1241557 : libmesh_assert(_children);
2337 6901684 : return {_children.get(), _children.get() + this->n_children()};
2338 : }
2339 : #endif // LIBMESH_ENABLE_AMR
2340 :
2341 :
2342 :
2343 :
2344 : // ------------------------------------------------------------
2345 : // global Elem functions
2346 :
2347 : inline
2348 0 : std::ostream & operator << (std::ostream & os, const Elem & e)
2349 : {
2350 0 : e.print_info(os);
2351 0 : return os;
2352 : }
2353 :
2354 :
2355 : // ------------------------------------------------------------
2356 : // Elem class member functions
2357 : inline
2358 1131203854 : Elem::Elem(const unsigned int nn,
2359 : const unsigned int ns,
2360 : Elem * p,
2361 : Elem ** elemlinkdata,
2362 1131203854 : Node ** nodelinkdata) :
2363 658408385 : _nodes(nodelinkdata),
2364 658408385 : _elemlinks(elemlinkdata),
2365 658408385 : _sbd_id(0),
2366 : #ifdef LIBMESH_ENABLE_AMR
2367 658408385 : _rflag(Elem::DO_NOTHING),
2368 658408385 : _pflag(Elem::DO_NOTHING),
2369 658408385 : _p_level(0),
2370 : #endif
2371 658990269 : _map_type(p ? p->mapping_type() : 0),
2372 1160006321 : _map_data(p ? p->mapping_data() : 0)
2373 : {
2374 1131203854 : this->processor_id() = DofObject::invalid_processor_id;
2375 :
2376 : // If this ever legitimately fails we need to increase max_n_nodes
2377 57049312 : libmesh_assert_less_equal(nn, max_n_nodes);
2378 :
2379 : // We currently only support refinement of elements into child
2380 : // elements of the same type. We can't test elem->type() here,
2381 : // because that's virtual and we're still in the base class
2382 : // constructor, but we can at least usually verify constency with
2383 : // the arguments we were handed.
2384 : #ifndef NDEBUG
2385 57049312 : if (p && !p->runtime_topology())
2386 : {
2387 287180 : libmesh_assert_equal_to(nn, p->n_nodes());
2388 287180 : libmesh_assert_equal_to(ns, p->n_sides());
2389 : }
2390 : #endif
2391 :
2392 : // Initialize the nodes data structure if we're given a pointer to
2393 : // memory for it.
2394 1131203854 : if (_nodes)
2395 : {
2396 5364372364 : for (unsigned int n=0; n<nn; n++)
2397 4233219058 : _nodes[n] = nullptr;
2398 : }
2399 :
2400 : // Initialize the neighbors/parent data structure
2401 : // _elemlinks = new Elem *[ns+1];
2402 :
2403 : // Initialize the elements data structure if we're given a pointer
2404 : // to memory for it. If we *weren't* given memory for it, e.g.
2405 : // because a subclass like an arbitrary Polygon needs to
2406 : // heap-allocate this memory, then that subclass will have to handle
2407 : // this initialization too.
2408 1131203854 : if (_elemlinks)
2409 : {
2410 1131169685 : _elemlinks[0] = p;
2411 :
2412 5171741450 : for (unsigned int n=1; n<ns+1; n++)
2413 4040571765 : _elemlinks[n] = nullptr;
2414 :
2415 : // Optionally initialize data from the parent
2416 1131169685 : if (this->parent())
2417 : {
2418 28802467 : this->subdomain_id() = this->parent()->subdomain_id();
2419 28802467 : this->processor_id() = this->parent()->processor_id();
2420 28802467 : _map_type = this->parent()->_map_type;
2421 28802467 : _map_data = this->parent()->_map_data;
2422 :
2423 : #ifdef LIBMESH_ENABLE_AMR
2424 29097171 : this->set_p_level(this->parent()->p_level());
2425 : #endif
2426 : }
2427 : }
2428 1131203854 : }
2429 :
2430 :
2431 :
2432 : inline
2433 3417328926 : const Point & Elem::point (const unsigned int i) const
2434 : {
2435 3417328926 : libmesh_assert_less (i, this->n_nodes());
2436 3417328926 : libmesh_assert(_nodes[i]);
2437 3417328926 : libmesh_assert_not_equal_to (_nodes[i]->id(), Node::invalid_id);
2438 :
2439 47431942560 : return *_nodes[i];
2440 : }
2441 :
2442 :
2443 :
2444 : inline
2445 2660766 : Point & Elem::point (const unsigned int i)
2446 : {
2447 2660766 : libmesh_assert_less (i, this->n_nodes());
2448 :
2449 167012950 : return *_nodes[i];
2450 : }
2451 :
2452 :
2453 :
2454 : inline
2455 87772341 : dof_id_type Elem::node_id (const unsigned int i) const
2456 : {
2457 87772341 : libmesh_assert_less (i, this->n_nodes());
2458 87772341 : libmesh_assert(_nodes[i]);
2459 87772341 : libmesh_assert_not_equal_to (_nodes[i]->id(), Node::invalid_id);
2460 :
2461 979106051 : return _nodes[i]->id();
2462 : }
2463 :
2464 :
2465 :
2466 : inline
2467 1374823 : unsigned int Elem::local_node (const dof_id_type i) const
2468 : {
2469 6447702 : for (auto n : make_range(this->n_nodes()))
2470 6447702 : if (this->node_id(n) == i)
2471 2328 : return n;
2472 :
2473 0 : return libMesh::invalid_uint;
2474 : }
2475 :
2476 :
2477 :
2478 : inline
2479 25757310 : const Node * const * Elem::get_nodes () const
2480 : {
2481 140029722 : return _nodes;
2482 : }
2483 :
2484 :
2485 :
2486 : inline
2487 113359507 : const Node * Elem::node_ptr (const unsigned int i) const
2488 : {
2489 113359507 : libmesh_assert_less (i, this->n_nodes());
2490 113359507 : libmesh_assert(_nodes[i]);
2491 :
2492 24244585811 : return _nodes[i];
2493 : }
2494 :
2495 :
2496 :
2497 : inline
2498 206792455 : Node * Elem::node_ptr (const unsigned int i)
2499 : {
2500 206792455 : libmesh_assert_less (i, this->n_nodes());
2501 206792455 : libmesh_assert(_nodes[i]);
2502 :
2503 4716922410 : return _nodes[i];
2504 : }
2505 :
2506 :
2507 :
2508 : inline
2509 19108478 : const Node & Elem::node_ref (const unsigned int i) const
2510 : {
2511 175166503 : return *this->node_ptr(i);
2512 : }
2513 :
2514 :
2515 :
2516 : inline
2517 26066851 : Node & Elem::node_ref (const unsigned int i)
2518 : {
2519 54581947 : return *this->node_ptr(i);
2520 : }
2521 :
2522 :
2523 :
2524 : inline
2525 1679854 : unsigned int Elem::get_node_index (const Node * node_ptr) const
2526 : {
2527 9014632 : for (auto n : make_range(this->n_nodes()))
2528 9002456 : if (this->_nodes[n] == node_ptr)
2529 160640 : return n;
2530 :
2531 0 : return libMesh::invalid_uint;
2532 : }
2533 :
2534 :
2535 :
2536 : #ifdef LIBMESH_ENABLE_DEPRECATED
2537 : inline
2538 0 : Node * & Elem::set_node (const unsigned int i)
2539 : {
2540 0 : libmesh_assert_less (i, this->n_nodes());
2541 :
2542 : libmesh_deprecated();
2543 :
2544 0 : return _nodes[i];
2545 : }
2546 : #endif // LIBMESH_ENABLE_DEPRECATED
2547 :
2548 :
2549 :
2550 : inline
2551 6250271598 : void Elem::set_node (const unsigned int i,
2552 : Node * node)
2553 : {
2554 177702 : libmesh_assert_less (i, this->n_nodes());
2555 :
2556 6538461113 : _nodes[i] = node;
2557 5886353898 : }
2558 :
2559 :
2560 :
2561 : inline
2562 57629362 : subdomain_id_type Elem::subdomain_id () const
2563 : {
2564 544443252 : return _sbd_id;
2565 : }
2566 :
2567 :
2568 :
2569 : inline
2570 42933511 : subdomain_id_type & Elem::subdomain_id ()
2571 : {
2572 124773900 : return _sbd_id;
2573 : }
2574 :
2575 :
2576 :
2577 : inline
2578 98443825 : const Elem * Elem::neighbor_ptr (unsigned int i) const
2579 : {
2580 98443825 : libmesh_assert_less (i, this->n_neighbors());
2581 :
2582 467161131 : return _elemlinks[i+1];
2583 : }
2584 :
2585 :
2586 :
2587 : inline
2588 18202435 : Elem * Elem::neighbor_ptr (unsigned int i)
2589 : {
2590 18202435 : libmesh_assert_less (i, this->n_neighbors());
2591 :
2592 1227440702 : return _elemlinks[i+1];
2593 : }
2594 :
2595 :
2596 :
2597 : inline
2598 9747194 : void Elem::set_neighbor (const unsigned int i, Elem * n)
2599 : {
2600 9747194 : libmesh_assert_less (i, this->n_neighbors());
2601 :
2602 406478512 : _elemlinks[i+1] = n;
2603 509642149 : }
2604 :
2605 :
2606 :
2607 : inline
2608 91040518 : bool Elem::has_neighbor (const Elem * elem) const
2609 : {
2610 272334224 : for (auto n : this->neighbor_ptr_range())
2611 261615379 : if (n == elem)
2612 80624522 : return true;
2613 :
2614 10415996 : return false;
2615 : }
2616 :
2617 :
2618 :
2619 : inline
2620 : Elem * Elem::child_neighbor (Elem * elem)
2621 : {
2622 : for (auto n : elem->neighbor_ptr_range())
2623 : if (n && n->parent() == this)
2624 : return n;
2625 :
2626 : return nullptr;
2627 : }
2628 :
2629 :
2630 :
2631 : inline
2632 : const Elem * Elem::child_neighbor (const Elem * elem) const
2633 : {
2634 : for (auto n : elem->neighbor_ptr_range())
2635 : if (n && n->parent() == this)
2636 : return n;
2637 :
2638 : return nullptr;
2639 : }
2640 :
2641 :
2642 :
2643 : inline
2644 : SimpleRange<Elem::NodeRefIter>
2645 193464098 : Elem::node_ref_range()
2646 : {
2647 202994074 : return {_nodes, _nodes+this->n_nodes()};
2648 : }
2649 :
2650 :
2651 :
2652 : inline
2653 : SimpleRange<Elem::ConstNodeRefIter>
2654 261960791 : Elem::node_ref_range() const
2655 : {
2656 279186192 : return {_nodes, _nodes+this->n_nodes()};
2657 : }
2658 :
2659 :
2660 :
2661 : inline
2662 : IntRange<unsigned short>
2663 60956284 : Elem::node_index_range() const
2664 : {
2665 1424868200 : return {0, cast_int<unsigned short>(this->n_nodes())};
2666 : }
2667 :
2668 :
2669 :
2670 : inline
2671 : IntRange<unsigned short>
2672 410664 : Elem::edge_index_range() const
2673 : {
2674 41967376 : return {0, cast_int<unsigned short>(this->n_edges())};
2675 : }
2676 :
2677 :
2678 :
2679 : inline
2680 : IntRange<unsigned short>
2681 308695 : Elem::face_index_range() const
2682 : {
2683 3655747 : return {0, cast_int<unsigned short>(this->n_faces())};
2684 : }
2685 :
2686 :
2687 :
2688 : inline
2689 : IntRange<unsigned short>
2690 8891180 : Elem::side_index_range() const
2691 : {
2692 357962690 : return {0, cast_int<unsigned short>(this->n_sides())};
2693 : }
2694 :
2695 :
2696 :
2697 :
2698 : inline
2699 826716870 : std::unique_ptr<const Elem> Elem::side_ptr (unsigned int i) const
2700 : {
2701 : // Call the non-const version of this function, return the result as
2702 : // a std::unique_ptr<const Elem>.
2703 75970818 : Elem * me = const_cast<Elem *>(this);
2704 902687688 : return me->side_ptr(i);
2705 : }
2706 :
2707 :
2708 :
2709 : inline
2710 : void
2711 15336 : Elem::side_ptr (std::unique_ptr<const Elem> & elem,
2712 : const unsigned int i) const
2713 : {
2714 : // Hand off to the non-const version of this function
2715 432 : Elem * me = const_cast<Elem *>(this);
2716 864 : std::unique_ptr<Elem> e {const_cast<Elem *>(elem.release())};
2717 15336 : me->side_ptr(e, i);
2718 14904 : elem = std::move(e);
2719 15336 : }
2720 :
2721 :
2722 :
2723 : inline
2724 : std::unique_ptr<const Elem>
2725 105961679 : Elem::build_side_ptr (const unsigned int i) const
2726 : {
2727 : // Call the non-const version of this function, return the result as
2728 : // a std::unique_ptr<const Elem>.
2729 2828115 : Elem * me = const_cast<Elem *>(this);
2730 116109785 : return me->build_side_ptr(i);
2731 : }
2732 :
2733 :
2734 :
2735 : inline
2736 : void
2737 38825034 : Elem::build_side_ptr (std::unique_ptr<const Elem> & elem,
2738 : const unsigned int i) const
2739 : {
2740 : // Hand off to the non-const version of this function
2741 827596 : Elem * me = const_cast<Elem *>(this);
2742 1655192 : std::unique_ptr<Elem> e {const_cast<Elem *>(elem.release())};
2743 38825034 : me->build_side_ptr(e, i);
2744 37156276 : elem = std::move(e);
2745 38825034 : }
2746 :
2747 :
2748 :
2749 : template <typename Sideclass, typename Subclass>
2750 : inline
2751 : std::unique_ptr<Elem>
2752 155887313 : Elem::simple_build_side_ptr (const unsigned int i)
2753 : {
2754 33129175 : libmesh_assert_less (i, this->n_sides());
2755 :
2756 155887313 : std::unique_ptr<Elem> face = std::make_unique<Sideclass>();
2757 1201217579 : for (auto n : face->node_index_range())
2758 1045330266 : face->set_node(n, this->node_ptr(Subclass::side_nodes_map[i][n]));
2759 :
2760 155887313 : face->set_interior_parent(this);
2761 144901583 : face->inherit_data_from(*this);
2762 :
2763 155887313 : return face;
2764 0 : }
2765 :
2766 :
2767 :
2768 : template <typename Subclass>
2769 : inline
2770 : void
2771 38530810 : Elem::simple_build_side_ptr (std::unique_ptr<Elem> & side,
2772 : const unsigned int i,
2773 : ElemType sidetype)
2774 : {
2775 1679657 : libmesh_assert_less (i, this->n_sides());
2776 :
2777 38530810 : if (!side.get() || side->type() != sidetype)
2778 : {
2779 147145 : Subclass & real_me = cast_ref<Subclass&>(*this);
2780 1450515 : side = real_me.Subclass::build_side_ptr(i);
2781 : }
2782 : else
2783 : {
2784 37731980 : side->set_interior_parent(this);
2785 36199437 : side->inherit_data_from(*this);
2786 267457308 : for (auto n : side->node_index_range())
2787 229725328 : side->set_node(n, this->node_ptr(Subclass::side_nodes_map[i][n]));
2788 : }
2789 38530810 : }
2790 :
2791 :
2792 :
2793 : template <typename Subclass, typename Mapclass>
2794 : inline
2795 : void
2796 237870206 : Elem::simple_side_ptr (std::unique_ptr<Elem> & side,
2797 : const unsigned int i,
2798 : ElemType sidetype)
2799 : {
2800 4784348 : libmesh_assert_less (i, this->n_sides());
2801 :
2802 237870206 : if (!side.get() || side->type() != sidetype)
2803 : {
2804 16242 : Subclass & real_me = cast_ref<Subclass&>(*this);
2805 1174252 : side = real_me.Subclass::side_ptr(i);
2806 : }
2807 : else
2808 : {
2809 242000609 : side->subdomain_id() = this->subdomain_id();
2810 :
2811 808368527 : for (auto n : side->node_index_range())
2812 571093568 : side->set_node(n, this->node_ptr(Mapclass::side_nodes_map[i][n]));
2813 : }
2814 237870206 : }
2815 :
2816 :
2817 :
2818 : inline
2819 : std::unique_ptr<const Elem>
2820 46700338 : Elem::build_edge_ptr (const unsigned int i) const
2821 : {
2822 : // Call the non-const version of this function, return the result as
2823 : // a std::unique_ptr<const Elem>.
2824 27214118 : Elem * me = const_cast<Elem *>(this);
2825 48809379 : return me->build_edge_ptr(i);
2826 : }
2827 :
2828 :
2829 :
2830 : inline
2831 : void
2832 99164 : Elem::build_edge_ptr (std::unique_ptr<const Elem> & elem,
2833 : const unsigned int i) const
2834 : {
2835 : // Hand off to the non-const version of this function
2836 3020 : Elem * me = const_cast<Elem *>(this);
2837 6040 : std::unique_ptr<Elem> e {const_cast<Elem *>(elem.release())};
2838 99164 : me->build_edge_ptr(e, i);
2839 96144 : elem = std::move(e);
2840 99164 : }
2841 :
2842 :
2843 : template <typename Edgeclass, typename Subclass>
2844 : inline
2845 : std::unique_ptr<Elem>
2846 27878431 : Elem::simple_build_edge_ptr (const unsigned int i)
2847 : {
2848 9812036 : libmesh_assert_less (i, this->n_edges());
2849 :
2850 27878431 : std::unique_ptr<Elem> edge = std::make_unique<Edgeclass>();
2851 :
2852 97685766 : for (auto n : edge->node_index_range())
2853 69807335 : edge->set_node(n, this->node_ptr(Subclass::edge_nodes_map[i][n]));
2854 :
2855 27878431 : edge->set_interior_parent(this);
2856 26092187 : edge->inherit_data_from(*this);
2857 :
2858 27878431 : return edge;
2859 0 : }
2860 :
2861 :
2862 :
2863 :
2864 : template <typename Subclass>
2865 : inline
2866 : void
2867 137744 : Elem::simple_build_edge_ptr (std::unique_ptr<Elem> & edge,
2868 : const unsigned int i,
2869 : ElemType edgetype)
2870 : {
2871 4196 : libmesh_assert_less (i, this->n_edges());
2872 :
2873 137744 : if (!edge.get() || edge->type() != edgetype)
2874 : {
2875 25 : Subclass & real_me = cast_ref<Subclass&>(*this);
2876 655 : edge = real_me.Subclass::build_edge_ptr(i);
2877 : }
2878 : else
2879 : {
2880 133233 : edge->inherit_data_from(*this);
2881 414504 : for (auto n : edge->node_index_range())
2882 277100 : edge->set_node(n, this->node_ptr(Subclass::edge_nodes_map[i][n]));
2883 : }
2884 137744 : }
2885 :
2886 :
2887 :
2888 : inline
2889 39 : bool Elem::on_boundary () const
2890 : {
2891 : // By convention, the element is on the boundary
2892 : // if it has a nullptr neighbor.
2893 351 : return this->has_neighbor(nullptr);
2894 : }
2895 :
2896 :
2897 :
2898 : inline
2899 155935029 : unsigned int Elem::which_neighbor_am_i (const Elem * e) const
2900 : {
2901 8215679 : libmesh_assert(e);
2902 :
2903 8215679 : const Elem * eparent = e;
2904 :
2905 156726376 : while (eparent->level() > this->level())
2906 : {
2907 262972 : eparent = eparent->parent();
2908 260896 : libmesh_assert(eparent);
2909 : }
2910 :
2911 408456411 : for (auto s : make_range(this->n_sides()))
2912 408448964 : if (this->neighbor_ptr(s) == eparent)
2913 8215679 : return s;
2914 :
2915 0 : return libMesh::invalid_uint;
2916 : }
2917 :
2918 :
2919 :
2920 : inline
2921 180556814 : bool Elem::active() const
2922 : {
2923 : #ifdef LIBMESH_ENABLE_AMR
2924 3780373335 : if ((this->refinement_flag() == INACTIVE) ||
2925 126867312 : (this->refinement_flag() == COARSEN_INACTIVE))
2926 53798188 : return false;
2927 : else
2928 126758626 : return true;
2929 : #else
2930 : return true;
2931 : #endif
2932 : }
2933 :
2934 :
2935 :
2936 :
2937 :
2938 : inline
2939 418195435 : bool Elem::subactive() const
2940 : {
2941 : #ifdef LIBMESH_ENABLE_AMR
2942 46093207 : if (this->active())
2943 34007824 : return false;
2944 12085383 : if (!this->has_children())
2945 3662996 : return true;
2946 85393393 : for (const Elem * my_ancestor = this->parent();
2947 297528150 : my_ancestor != nullptr;
2948 28658029 : my_ancestor = my_ancestor->parent())
2949 24846628 : if (my_ancestor->active())
2950 10516 : return true;
2951 : #endif
2952 :
2953 8411871 : return false;
2954 : }
2955 :
2956 :
2957 :
2958 : inline
2959 46254826 : bool Elem::has_children() const
2960 : {
2961 : #ifdef LIBMESH_ENABLE_AMR
2962 358749024 : if (!_children)
2963 9120437 : return false;
2964 : else
2965 37134389 : return true;
2966 : #else
2967 : return false;
2968 : #endif
2969 : }
2970 :
2971 :
2972 : inline
2973 : bool Elem::has_ancestor_children() const
2974 : {
2975 : #ifdef LIBMESH_ENABLE_AMR
2976 : if (!_children)
2977 : return false;
2978 : else
2979 : for (auto & c : child_ref_range())
2980 : if (c.has_children())
2981 : return true;
2982 : #endif
2983 : return false;
2984 : }
2985 :
2986 :
2987 :
2988 : inline
2989 0 : bool Elem::is_ancestor_of(const Elem *
2990 : #ifdef LIBMESH_ENABLE_AMR
2991 : descendant
2992 : #endif
2993 : ) const
2994 : {
2995 : #ifdef LIBMESH_ENABLE_AMR
2996 0 : const Elem * e = descendant;
2997 10626 : while (e)
2998 : {
2999 8638 : if (this == e)
3000 0 : return true;
3001 0 : e = e->parent();
3002 : }
3003 : #endif
3004 0 : return false;
3005 : }
3006 :
3007 :
3008 :
3009 : inline
3010 1175985739 : const Elem * Elem::parent () const
3011 : {
3012 2541422960 : return _elemlinks[0];
3013 : }
3014 :
3015 :
3016 :
3017 : inline
3018 98478222 : Elem * Elem::parent ()
3019 : {
3020 1585373382 : return _elemlinks[0];
3021 : }
3022 :
3023 :
3024 :
3025 : inline
3026 344184 : void Elem::set_parent (Elem * p)
3027 : {
3028 : // We no longer support using parent() as interior_parent()
3029 344184 : libmesh_assert_equal_to(this->dim(), p ? p->dim() : this->dim());
3030 10101067 : _elemlinks[0] = p;
3031 527167 : }
3032 :
3033 :
3034 :
3035 : inline
3036 824242 : const Elem * Elem::top_parent () const
3037 : {
3038 824242 : const Elem * tp = this;
3039 :
3040 : // Keep getting the element's parent
3041 : // until that parent is at level-0
3042 4369837 : while (tp->parent() != nullptr)
3043 2044533 : tp = tp->parent();
3044 :
3045 824242 : libmesh_assert(tp);
3046 824242 : libmesh_assert_equal_to (tp->level(), 0);
3047 :
3048 824242 : return tp;
3049 : }
3050 :
3051 :
3052 :
3053 : inline
3054 7214879204 : unsigned int Elem::level() const
3055 : {
3056 : #ifdef LIBMESH_ENABLE_AMR
3057 :
3058 : // if I don't have a parent I was
3059 : // created directly from file
3060 : // or by the user, so I am a
3061 : // level-0 element
3062 19956978999 : if (this->parent() == nullptr)
3063 112102499 : return 0;
3064 :
3065 : // if the parent and this element are of different
3066 : // dimensionality we are at the same level as
3067 : // the parent (e.g. we are the 2D side of a
3068 : // 3D element)
3069 13041120849 : if (this->dim() != this->parent()->dim())
3070 0 : return this->parent()->level();
3071 :
3072 : // otherwise we are at a level one
3073 : // higher than our parent
3074 13164307529 : return (this->parent()->level() + 1);
3075 :
3076 : #else
3077 :
3078 : // Without AMR all elements are
3079 : // at level 0.
3080 : return 0;
3081 :
3082 : #endif
3083 : }
3084 :
3085 :
3086 :
3087 : inline
3088 2172863273 : unsigned int Elem::p_level() const
3089 : {
3090 : #ifdef LIBMESH_ENABLE_AMR
3091 53856786316 : return _p_level;
3092 : #else
3093 : return 0;
3094 : #endif
3095 : }
3096 :
3097 :
3098 :
3099 : inline
3100 100868534 : ElemMappingType Elem::mapping_type () const
3101 : {
3102 1263113540 : return static_cast<ElemMappingType>(_map_type);
3103 : }
3104 :
3105 :
3106 :
3107 : inline
3108 46788722 : void Elem::set_mapping_type(const ElemMappingType type)
3109 : {
3110 302180741 : _map_type = cast_int<unsigned char>(type);
3111 46788722 : }
3112 :
3113 :
3114 :
3115 : inline
3116 33329523 : unsigned char Elem::mapping_data () const
3117 : {
3118 138145373 : return _map_data;
3119 : }
3120 :
3121 :
3122 :
3123 : inline
3124 46788715 : void Elem::set_mapping_data(const unsigned char data)
3125 : {
3126 301917967 : _map_data = data;
3127 47051822 : }
3128 :
3129 :
3130 :
3131 : #ifdef LIBMESH_ENABLE_AMR
3132 :
3133 : inline
3134 0 : const Elem * Elem::raw_child_ptr (unsigned int i) const
3135 : {
3136 4240 : if (!_children)
3137 0 : return nullptr;
3138 :
3139 2880 : return _children[i];
3140 : }
3141 :
3142 : inline
3143 90036455 : const Elem * Elem::child_ptr (unsigned int i) const
3144 : {
3145 90036455 : libmesh_assert(_children);
3146 90036455 : libmesh_assert(_children[i]);
3147 :
3148 341161385 : return _children[i];
3149 : }
3150 :
3151 : inline
3152 2095804 : Elem * Elem::child_ptr (unsigned int i)
3153 : {
3154 2095804 : libmesh_assert(_children);
3155 2095804 : libmesh_assert(_children[i]);
3156 :
3157 23606159 : return _children[i];
3158 : }
3159 :
3160 :
3161 : inline
3162 133312 : void Elem::set_child (unsigned int c, Elem * elem)
3163 : {
3164 133312 : libmesh_assert (this->has_children());
3165 :
3166 52216408 : _children[c] = elem;
3167 27078398 : }
3168 :
3169 :
3170 :
3171 : inline
3172 109569346 : unsigned int Elem::which_child_am_i (const Elem * e) const
3173 : {
3174 27330342 : libmesh_assert(e);
3175 27330342 : libmesh_assert (this->has_children());
3176 :
3177 109569346 : unsigned int nc = this->n_children();
3178 292065503 : for (unsigned int c=0; c != nc; c++)
3179 292065503 : if (this->child_ptr(c) == e)
3180 27330342 : return c;
3181 :
3182 0 : libmesh_error_msg("ERROR: which_child_am_i() was called with a non-child!");
3183 :
3184 : return libMesh::invalid_uint;
3185 : }
3186 :
3187 :
3188 :
3189 : inline
3190 352184148 : Elem::RefinementState Elem::refinement_flag () const
3191 : {
3192 4622297333 : return static_cast<RefinementState>(_rflag);
3193 : }
3194 :
3195 :
3196 :
3197 : inline
3198 3563040 : void Elem::set_refinement_flag(RefinementState rflag)
3199 : {
3200 60281177 : _rflag = cast_int<unsigned char>(rflag);
3201 11435438 : }
3202 :
3203 :
3204 :
3205 : inline
3206 86988644 : Elem::RefinementState Elem::p_refinement_flag () const
3207 : {
3208 234388605 : return static_cast<RefinementState>(_pflag);
3209 : }
3210 :
3211 :
3212 :
3213 : inline
3214 2523117 : void Elem::set_p_refinement_flag(RefinementState pflag)
3215 : {
3216 2522768 : if (this->p_level() == 0)
3217 2517656 : libmesh_assert_not_equal_to
3218 : (pflag, Elem::JUST_REFINED);
3219 :
3220 40630022 : _pflag = cast_int<unsigned char>(pflag);
3221 25205360 : }
3222 :
3223 :
3224 :
3225 : inline
3226 0 : unsigned int Elem::max_descendant_p_level () const
3227 : {
3228 : // This is undefined for subactive elements,
3229 : // which have no active descendants
3230 0 : libmesh_assert (!this->subactive());
3231 0 : if (this->active())
3232 0 : return this->p_level();
3233 :
3234 0 : unsigned int max_p_level = _p_level;
3235 0 : for (auto & c : child_ref_range())
3236 0 : max_p_level = std::max(max_p_level,
3237 0 : c.max_descendant_p_level());
3238 0 : return max_p_level;
3239 : }
3240 :
3241 :
3242 :
3243 : inline
3244 46748627 : void Elem::hack_p_level(unsigned int p)
3245 : {
3246 46748627 : if (p == 0)
3247 46664767 : libmesh_assert_not_equal_to
3248 : (this->p_refinement_flag(), Elem::JUST_REFINED);
3249 :
3250 356952146 : _p_level = cast_int<unsigned char>(p);
3251 296188329 : }
3252 :
3253 :
3254 : inline
3255 3602 : void Elem::hack_p_level_and_refinement_flag (unsigned int p,
3256 : RefinementState pflag)
3257 : {
3258 34518747 : _pflag = cast_int<unsigned char>(pflag);
3259 3602 : this->hack_p_level(p);
3260 34499920 : }
3261 :
3262 : #endif // ifdef LIBMESH_ENABLE_AMR
3263 :
3264 :
3265 : inline
3266 290889 : void Elem::orient(BoundaryInfo * boundary_info)
3267 : {
3268 304149 : if (this->is_flipped())
3269 138978 : this->flip(boundary_info);
3270 290889 : }
3271 :
3272 :
3273 : inline
3274 43566 : dof_id_type Elem::compute_key (dof_id_type n0)
3275 : {
3276 43566 : return n0;
3277 : }
3278 :
3279 :
3280 :
3281 : inline
3282 6484395 : dof_id_type Elem::compute_key (dof_id_type n0,
3283 : dof_id_type n1)
3284 : {
3285 : // Order the two so that n0 < n1
3286 177667983 : if (n0 > n1) std::swap (n0, n1);
3287 :
3288 177667983 : return Utility::hashword2(n0, n1);
3289 : }
3290 :
3291 :
3292 :
3293 : inline
3294 49970180 : dof_id_type Elem::compute_key (dof_id_type n0,
3295 : dof_id_type n1,
3296 : dof_id_type n2)
3297 : {
3298 49970180 : std::array<dof_id_type, 3> array = {{n0, n1, n2}};
3299 1792288 : std::sort(array.begin(), array.end());
3300 51762468 : return Utility::hashword(array);
3301 : }
3302 :
3303 :
3304 :
3305 : inline
3306 35914773 : dof_id_type Elem::compute_key (dof_id_type n0,
3307 : dof_id_type n1,
3308 : dof_id_type n2,
3309 : dof_id_type n3)
3310 : {
3311 35914773 : std::array<dof_id_type, 4> array = {{n0, n1, n2, n3}};
3312 1077366 : std::sort(array.begin(), array.end());
3313 36992139 : return Utility::hashword(array);
3314 : }
3315 :
3316 :
3317 :
3318 : inline
3319 228006763 : void Elem::inherit_data_from (const Elem & src)
3320 : {
3321 60834302 : this->set_mapping_type(src.mapping_type());
3322 60834302 : this->set_mapping_data(src.mapping_data());
3323 243331595 : this->subdomain_id() = src.subdomain_id();
3324 60834302 : this->processor_id(src.processor_id());
3325 : #ifdef LIBMESH_ENABLE_AMR
3326 243331595 : this->set_p_level(src.p_level());
3327 : #endif
3328 228006763 : }
3329 :
3330 :
3331 :
3332 : /**
3333 : * The definition of the protected nested SideIter class.
3334 : */
3335 0 : class Elem::SideIter
3336 : {
3337 : public:
3338 : // Constructor with arguments.
3339 0 : SideIter(const unsigned int side_number,
3340 : Elem * parent)
3341 0 : : _side(),
3342 0 : _side_ptr(nullptr),
3343 0 : _parent(parent),
3344 0 : _side_number(side_number)
3345 0 : {}
3346 :
3347 :
3348 : // Empty constructor.
3349 : SideIter()
3350 : : _side(),
3351 : _side_ptr(nullptr),
3352 : _parent(nullptr),
3353 : _side_number(libMesh::invalid_uint)
3354 : {}
3355 :
3356 :
3357 : // Copy constructor
3358 0 : SideIter(const SideIter & other)
3359 0 : : _side(),
3360 0 : _side_ptr(nullptr),
3361 0 : _parent(other._parent),
3362 0 : _side_number(other._side_number)
3363 0 : {}
3364 :
3365 :
3366 : // op=
3367 : SideIter & operator=(const SideIter & other)
3368 : {
3369 : this->_parent = other._parent;
3370 : this->_side_number = other._side_number;
3371 : return *this;
3372 : }
3373 :
3374 : // unary op*
3375 0 : Elem *& operator*() const
3376 : {
3377 : // Set the std::unique_ptr
3378 0 : this->_update_side_ptr();
3379 :
3380 : // Return a reference to _side_ptr
3381 0 : return this->_side_ptr;
3382 : }
3383 :
3384 : // op++
3385 0 : SideIter & operator++()
3386 : {
3387 0 : ++_side_number;
3388 0 : return *this;
3389 : }
3390 :
3391 : // op== Two side iterators are equal if they have
3392 : // the same side number and the same parent element.
3393 0 : bool operator == (const SideIter & other) const
3394 : {
3395 0 : return (this->_side_number == other._side_number &&
3396 0 : this->_parent == other._parent);
3397 : }
3398 :
3399 :
3400 : // Consults the parent Elem to determine if the side
3401 : // is a boundary side. Note: currently side N is a
3402 : // boundary side if neighbor N is nullptr. Be careful,
3403 : // this could possibly change in the future?
3404 0 : bool side_on_boundary() const
3405 : {
3406 0 : return this->_parent->neighbor_ptr(_side_number) == nullptr;
3407 : }
3408 :
3409 : private:
3410 : // Update the _side pointer by building the correct side.
3411 : // This has to be called before dereferencing.
3412 0 : void _update_side_ptr() const
3413 : {
3414 : // Construct new side, store in std::unique_ptr
3415 0 : this->_side = this->_parent->build_side_ptr(this->_side_number);
3416 :
3417 : // Also set our internal naked pointer. Memory is still owned
3418 : // by the std::unique_ptr.
3419 0 : this->_side_ptr = _side.get();
3420 0 : }
3421 :
3422 : // std::unique_ptr to the actual side, handles memory management for
3423 : // the sides which are created during the course of iteration.
3424 : mutable std::unique_ptr<Elem> _side;
3425 :
3426 : // Raw pointer needed to facilitate passing back to the user a
3427 : // reference to a non-temporary raw pointer in order to conform to
3428 : // the variant_filter_iterator interface. It points to the same
3429 : // thing the std::unique_ptr "_side" above holds. What happens if the user
3430 : // calls delete on the pointer passed back? Well, this is an issue
3431 : // which is not addressed by the iterators in libMesh. Basically it
3432 : // is a bad idea to ever call delete on an iterator from the library.
3433 : mutable Elem * _side_ptr;
3434 :
3435 : // Pointer to the parent Elem class which generated this iterator
3436 : Elem * _parent;
3437 :
3438 : // A counter variable which keeps track of the side number
3439 : unsigned int _side_number;
3440 : };
3441 :
3442 :
3443 :
3444 :
3445 :
3446 :
3447 : // Private implementation functions in the Elem class for the side iterators.
3448 : // They have to come after the definition of the SideIter class.
3449 : inline
3450 0 : Elem::SideIter Elem::_first_side()
3451 : {
3452 0 : return SideIter(0, this);
3453 : }
3454 :
3455 :
3456 :
3457 : inline
3458 0 : Elem::SideIter Elem::_last_side()
3459 : {
3460 0 : return SideIter(this->n_neighbors(), this);
3461 : }
3462 :
3463 :
3464 :
3465 :
3466 : /**
3467 : * The definition of the struct used for iterating over sides.
3468 : */
3469 : struct
3470 : Elem::side_iterator : variant_filter_iterator<Elem::Predicate, Elem *>
3471 : {
3472 : // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
3473 : template <typename PredType, typename IterType>
3474 0 : side_iterator (const IterType & d,
3475 : const IterType & e,
3476 : const PredType & p ) :
3477 0 : variant_filter_iterator<Elem::Predicate, Elem *>(d,e,p) {}
3478 : };
3479 :
3480 :
3481 :
3482 : inline
3483 77875859 : SimpleRange<Elem::NeighborPtrIter> Elem::neighbor_ptr_range()
3484 : {
3485 80955372 : return {_elemlinks+1, _elemlinks + 1 + this->n_neighbors()};
3486 : }
3487 :
3488 :
3489 : inline
3490 386729663 : SimpleRange<Elem::ConstNeighborPtrIter> Elem::neighbor_ptr_range() const
3491 : {
3492 388150744 : return {_elemlinks+1, _elemlinks + 1 + this->n_neighbors()};
3493 : }
3494 :
3495 : } // namespace libMesh
3496 :
3497 :
3498 : // Helper function for default caches in Elem subclasses
3499 :
3500 : #define LIBMESH_ENABLE_TOPOLOGY_CACHES \
3501 : virtual \
3502 : std::vector<std::vector<std::vector<std::vector<std::pair<unsigned char, unsigned char>>>>> & \
3503 : _get_bracketing_node_cache() const override \
3504 : { \
3505 : static std::vector<std::vector<std::vector<std::vector<std::pair<unsigned char, unsigned char>>>>> c; \
3506 : return c; \
3507 : } \
3508 : \
3509 : virtual \
3510 : std::vector<std::vector<std::vector<signed char>>> & \
3511 : _get_parent_indices_cache() const override \
3512 : { \
3513 : static std::vector<std::vector<std::vector<signed char>>> c; \
3514 : return c; \
3515 : }
3516 :
3517 :
3518 :
3519 :
3520 :
3521 :
3522 : #endif // LIBMESH_ELEM_H
|