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_CELL_TET10_H 21 : #define LIBMESH_CELL_TET10_H 22 : 23 : // Local includes 24 : #include "libmesh/cell_tet.h" 25 : 26 : namespace libMesh 27 : { 28 : 29 : /** 30 : * The \p Tet10 is an element in 3D composed of 10 nodes. It is 31 : * numbered like this: 32 : * 33 : * \verbatim 34 : * 3 35 : * TET10: o 36 : * /|\ 37 : * / | \ 38 : * 7 / | \9 39 : * o | o zeta 40 : * / |8 \ ^ 41 : * / o \ | 42 : * / 6 | \ | 43 : * 0 o.....o.|.......o 2 o---> eta 44 : * \ | / \ 45 : * \ | / \ 46 : * \ | / xi (out of page) 47 : * 4 o | o 5 48 : * \ | / 49 : * \ | / 50 : * \|/ 51 : * o 52 : * 1 53 : * \endverbatim 54 : * 55 : * (xi, eta, zeta): { 0 <= xi <= 1 56 : * { 0 <= eta <= 1 57 : * { 0 <= zeta <= 1 58 : * { xi + eta + zeta <= 1 59 : * 60 : * \author Benjamin S. Kirk 61 : * \date 2002 62 : * \brief A 3D tetrahedral element with 10 nodes. 63 : */ 64 : class Tet10 final : public Tet 65 : { 66 : public: 67 : 68 : /** 69 : * Constructor. By default this element has no parent. 70 : */ 71 : explicit 72 1441305 : Tet10 (Elem * p=nullptr) : 73 1441305 : Tet(num_nodes, p, _nodelinks_data) 74 1441305 : {} 75 : 76 : Tet10 (Tet10 &&) = delete; 77 : Tet10 (const Tet10 &) = delete; 78 : Tet10 & operator= (const Tet10 &) = delete; 79 : Tet10 & operator= (Tet10 &&) = delete; 80 1500840 : virtual ~Tet10() = default; 81 : 82 : /** 83 : * \returns \p TET10. 84 : */ 85 1003946696 : virtual ElemType type () const override { return TET10; } 86 : 87 : /** 88 : * \returns 10. 89 : */ 90 115536263 : virtual unsigned int n_nodes() const override { return num_nodes; } 91 : 92 : /** 93 : * \returns 8. 94 : */ 95 0 : virtual unsigned int n_sub_elem() const override { return 8; } 96 : 97 : /** 98 : * \returns \p true if the specified (local) node number is a vertex. 99 : */ 100 : virtual bool is_vertex(const unsigned int i) const override; 101 : 102 : /** 103 : * \returns \p true if the specified (local) node number is an edge. 104 : */ 105 : virtual bool is_edge(const unsigned int i) const override; 106 : 107 : /** 108 : * \returns \p true if the specified (local) node number is a face. 109 : */ 110 : virtual bool is_face(const unsigned int i) const override; 111 : 112 : /** 113 : * \returns \p true if the specified (local) node number is on the 114 : * specified side. 115 : */ 116 : virtual bool is_node_on_side(const unsigned int n, 117 : const unsigned int s) const override; 118 : 119 : virtual std::vector<unsigned int> nodes_on_side(const unsigned int s) const override; 120 : 121 : virtual std::vector<unsigned int> nodes_on_edge(const unsigned int e) const override; 122 : 123 : /** 124 : * \returns \p true if the specified (local) node number is on the 125 : * specified edge. 126 : */ 127 : virtual bool is_node_on_edge(const unsigned int n, 128 : const unsigned int e) const override; 129 : 130 : /** 131 : * \returns \p true if the specified child is on the 132 : * specified side. 133 : */ 134 : virtual bool is_child_on_side(const unsigned int c, 135 : const unsigned int s) const override; 136 : 137 : /** 138 : * \returns \p true if the element map is definitely affine within 139 : * numerical tolerances. 140 : */ 141 : virtual bool has_affine_map () const override; 142 : 143 : /** 144 : * \returns SECOND. 145 : */ 146 : virtual Order default_order() const override; 147 : 148 : /** 149 : * \returns \p Tet10::side_nodes_map[side][side_node] after doing some range checking. 150 : */ 151 : virtual unsigned int local_side_node(unsigned int side, 152 : unsigned int side_node) const override; 153 : 154 : /** 155 : * \returns \p Tet10::edge_nodes_map[edge][edge_node] after doing some range checking. 156 : */ 157 : virtual unsigned int local_edge_node(unsigned int edge, 158 : unsigned int edge_node) const override; 159 : 160 : /** 161 : * Builds a \p TRI6 built coincident with face i. 162 : * The \p std::unique_ptr<Elem> handles the memory aspect. 163 : */ 164 : virtual std::unique_ptr<Elem> build_side_ptr (const unsigned int i) override; 165 : 166 : /** 167 : * Rebuilds a TRI6 built coincident with face i. 168 : */ 169 : virtual void build_side_ptr (std::unique_ptr<Elem> & elem, 170 : const unsigned int i) override; 171 : 172 : // Avoid hiding deprecated version with different signature 173 : using Elem::build_side_ptr; 174 : 175 : /** 176 : * Builds a \p EDGE3 built coincident with edge i. 177 : * The \p std::unique_ptr<Elem> handles the memory aspect. 178 : */ 179 : virtual std::unique_ptr<Elem> build_edge_ptr (const unsigned int i) override; 180 : 181 : /** 182 : * Rebuilds a \p EDGE3 coincident with edge i. 183 : */ 184 : virtual void build_edge_ptr (std::unique_ptr<Elem> & edge, const unsigned int i) override; 185 : 186 : virtual void connectivity(const unsigned int sc, 187 : const IOPackage iop, 188 : std::vector<dof_id_type> & conn) const override; 189 : 190 : /** 191 : * \returns 2 for all \p n. 192 : */ 193 3978288 : virtual unsigned int n_second_order_adjacent_vertices (const unsigned int) const override 194 3978288 : { return 2; } 195 : 196 : /** 197 : * \returns The element-local number of the \f$ v^{th} \f$ vertex 198 : * that defines the \f$ n^{th} \f$ second-order node. 199 : * 200 : * \note \p n is counted as depicted above, \f$ 4 \le n < 10 \f$. 201 : */ 202 : virtual unsigned short int second_order_adjacent_vertex (const unsigned int n, 203 : const unsigned int v) const override; 204 : 205 : /** 206 : * \returns The child number \p c and element-local index \p v of the 207 : * \f$ n^{th} \f$ second-order node on the parent element. See 208 : * elem.h for further details. 209 : */ 210 : virtual std::pair<unsigned short int, unsigned short int> 211 : second_order_child_vertex (const unsigned int n) const override; 212 : 213 : /** 214 : * Geometric constants for Tet10. 215 : */ 216 : static const int num_nodes = 10; 217 : static const int nodes_per_side = 6; 218 : static const int nodes_per_edge = 3; 219 : 220 : /** 221 : * This maps the \f$ j^{th} \f$ node of the \f$ i^{th} \f$ side to 222 : * element node numbers. 223 : */ 224 : static const unsigned int side_nodes_map[num_sides][nodes_per_side]; 225 : 226 : /** 227 : * This maps the \f$ j^{th} \f$ node of the \f$ i^{th} \f$ edge to 228 : * element node numbers. 229 : */ 230 : static const unsigned int edge_nodes_map[num_edges][nodes_per_edge]; 231 : 232 : /** 233 : * A specialization for computing the volume of a Tet10. 234 : */ 235 : virtual Real volume () const override; 236 : 237 : virtual void permute(unsigned int perm_num) override final; 238 : 239 : virtual void flip(BoundaryInfo *) override final; 240 : 241 : ElemType side_type (const unsigned int s) const override final; 242 : 243 : protected: 244 : 245 : /** 246 : * Data for links to nodes. 247 : */ 248 : Node * _nodelinks_data[num_nodes]; 249 : 250 : 251 : 252 : #ifdef LIBMESH_ENABLE_AMR 253 : 254 : /** 255 : * Matrix used to create the elements children. 256 : */ 257 : virtual Real embedding_matrix (const unsigned int i, 258 : const unsigned int j, 259 : const unsigned int k) const override; 260 : 261 : /** 262 : * Matrix that computes new nodal locations/solution values 263 : * from current nodes/solution. 264 : */ 265 : static const Real _embedding_matrix[num_children][num_nodes][num_nodes]; 266 : 267 15164726 : LIBMESH_ENABLE_TOPOLOGY_CACHES; 268 : 269 : #endif // LIBMESH_ENABLE_AMR 270 : 271 : private: 272 : 273 : /** 274 : * Matrix that tells which vertices define the location 275 : * of mid-side (or second-order) nodes 276 : */ 277 : static const unsigned short int _second_order_adjacent_vertices[6][2]; 278 : 279 : /** 280 : * Vector that names a child sharing each second order node. 281 : */ 282 : static const unsigned short int _second_order_vertex_child_number[10]; 283 : 284 : /** 285 : * Vector that names the child vertex index for each second order node. 286 : */ 287 : static const unsigned short int _second_order_vertex_child_index[10]; 288 : }; 289 : 290 : } // namespace libMesh 291 : 292 : 293 : #endif // LIBMESH_CELL_TET10_H