1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2020 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
17 
18 
19 
20 #ifndef LIBMESH_CELL_TET_H
21 #define LIBMESH_CELL_TET_H
22 
23 // Local includes
24 #include "libmesh/cell.h"
25 
26 namespace libMesh
27 {
28 
29 /**
30  * The \p Tet is an element in 3D composed of 4 sides.
31  *
32  * \author Benjamin S. Kirk
33  * \date 2002
34  * \brief The base class for all tetrahedral element types.
35  */
36 class Tet : public Cell
37 {
38 public:
39 
40   /**
41    * Default tetrahedral element, takes number of nodes and
42    * parent. Derived classes implement 'true' elements.
43    */
Tet(const unsigned int nn,Elem * p,Node ** nodelinkdata)44   Tet (const unsigned int nn, Elem * p, Node ** nodelinkdata) :
45     Cell(nn, Tet::n_sides(), p, _elemlinks_data, nodelinkdata),
46     _diagonal_selection(INVALID_DIAG)
47   {
48     // Make sure the interior parent isn't undefined
49     if (LIBMESH_DIM > 3)
50       this->set_interior_parent(nullptr);
51   }
52 
53   Tet (Tet &&) = delete;
54   Tet (const Tet &) = delete;
55   Tet & operator= (const Tet &) = delete;
56   Tet & operator= (Tet &&) = delete;
57   virtual ~Tet() = default;
58 
59   /**
60    * \returns The \p Point associated with local \p Node \p i,
61    * in master element rather than physical coordinates.
62    */
master_point(const unsigned int i)63   virtual Point master_point (const unsigned int i) const override final
64   {
65     libmesh_assert_less(i, this->n_nodes());
66     return Point(_master_points[i][0],
67                  _master_points[i][1],
68                  _master_points[i][2]);
69   }
70 
71   /**
72    * \returns 4.
73    */
n_sides()74   virtual unsigned int n_sides() const override final { return 4; }
75 
76   /**
77    * \returns 4.  All tetrahedra have 4 vertices.
78    */
n_vertices()79   virtual unsigned int n_vertices() const override final { return 4; }
80 
81   /**
82    * \returns 6.  All tetrahedra have 6 edges.
83    */
n_edges()84   virtual unsigned int n_edges() const override final { return 6; }
85 
86   /**
87    * \returns 4.  All tetrahedra have 4 faces.
88    */
n_faces()89   virtual unsigned int n_faces() const override final { return 4; }
90 
91   /**
92    * \returns 8.
93    */
n_children()94   virtual unsigned int n_children() const override final { return 8; }
95 
96   /**
97    * \returns \p true if the specified edge is on the specified side.
98    */
99   virtual bool is_edge_on_side(const unsigned int e,
100                                const unsigned int s) const override final;
101 
102   /**
103    * Don't hide Elem::key() defined in the base class.
104    */
105   using Elem::key;
106 
107   /**
108    * \returns An id associated with the \p s side of this element.
109    * The id is not necessarily unique, but should be close.  This is
110    * particularly useful in the \p MeshBase::find_neighbors() routine.
111    */
112   virtual dof_id_type key (const unsigned int s) const override;
113 
114   /**
115    * \returns \p Tet4::side_nodes_map[side][side_node] after doing some range checking.
116    */
117   virtual unsigned int local_side_node(unsigned int side,
118                                        unsigned int side_node) const override;
119 
120   /**
121    * \returns \p Tet4::edge_nodes_map[edge][edge_node] after doing some range checking.
122    */
123   virtual unsigned int local_edge_node(unsigned int edge,
124                                        unsigned int edge_node) const override;
125 
126   /**
127    * \returns A primitive (3-noded) triangle for face i.
128    */
129   virtual std::unique_ptr<Elem> side_ptr (const unsigned int i) override final;
130 
131   /**
132    * Rebuilds a primitive (3-noded) triangle for face i.
133    */
134   virtual void side_ptr (std::unique_ptr<Elem> & side, const unsigned int i) override final;
135 
136   /**
137    * \returns A quantitative assessment of element quality based on
138    * the quality metric \p q specified by the user.
139    */
140   virtual Real quality (const ElemQuality q) const override;
141 
142   /**
143    * \returns The suggested quality bounds for the hex based on quality
144    * measure \p q.  These are the values suggested by the CUBIT User's
145    * Manual.
146    */
147   virtual std::pair<Real, Real> qual_bounds (const ElemQuality q) const override;
148 
149   /**
150    * This enumeration keeps track of which diagonal is selected during
151    * refinement.  In general there are three possible diagonals to
152    * choose when splitting the octahedron, and by choosing the shortest
153    * one we obtain the best element shape.
154    */
155   enum Diagonal
156     {
157       DIAG_02_13=0,    // diagonal between edges (0,2) and (1,3)
158       DIAG_03_12=1,    // diagonal between edges (0,3) and (1,2)
159       DIAG_01_23=2,    // diagonal between edges (0,1) and (2,3)
160       INVALID_DIAG=99  // diagonal not yet selected
161     };
162 
163   /**
164    * \returns The diagonal that has been selected during refinement.
165    */
diagonal_selection()166   Diagonal diagonal_selection () const { return _diagonal_selection; }
167 
168   /**
169    * Allows the user to select the diagonal for the refinement.  This
170    * function may only be called before the element is ever refined.
171    */
172   void select_diagonal (const Diagonal diag) const;
173 
174   virtual std::vector<unsigned int> sides_on_edge(const unsigned int e) const override final;
175 
176 
177 #ifdef LIBMESH_ENABLE_AMR
178 
179 
180   /**
181    * Tetrahedral elements permute the embedding matrix depending on which
182    * interior diagonal is used to subdivide into child elements.
183    * But we want to cache topology data based on that matrix.  So we return a
184    * "version number" based on the diagonal selection.
185    */
embedding_matrix_version()186   virtual unsigned int embedding_matrix_version () const override final
187   {
188     this->choose_diagonal();
189     return this->diagonal_selection();
190   }
191 
192 #endif // LIBMESH_ENABLE_AMR
193 
194 
195 
196 protected:
197 
198   /**
199    * Data for links to parent/neighbor/interior_parent elements.
200    */
201   Elem * _elemlinks_data[5+(LIBMESH_DIM>3)];
202 
203   /**
204    * Master element node locations
205    */
206   static const Real _master_points[10][3];
207 
208   /**
209    * Called by descendant classes with appropriate data to determine
210    * if child c is on side s.  Only works if LIBMESH_ENABLE_AMR.
211    */
212   bool is_child_on_side_helper(const unsigned int c,
213                                const unsigned int s,
214                                const unsigned int checked_nodes[][3] ) const;
215 
216   /**
217    * The currently-selected diagonal used during refinement.
218    * Initialized to INVALID_DIAG.
219    */
220   mutable Diagonal _diagonal_selection;
221 
222   /**
223    * Derived classes use this function to select an initial
224    * diagonal during refinement. The optimal choice is the shortest
225    * of the three.
226    */
227   void choose_diagonal() const;
228 };
229 
230 } // namespace libMesh
231 
232 #endif // LIBMESH_CELL_TET_H
233