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_TET4_H 21 : #define LIBMESH_CELL_TET4_H 22 : 23 : // Local includes 24 : #include "libmesh/cell_tet.h" 25 : 26 : namespace libMesh 27 : { 28 : 29 : /** 30 : * The \p Tet4 is an element in 3D composed of 4 nodes. It is 31 : * numbered like this: 32 : * 33 : * \verbatim 34 : * TET4: 35 : * 3 36 : * o zeta 37 : * /|\ ^ 38 : * / | \ | 39 : * / | \ | 40 : * 0 o...|...o 2 o---> eta 41 : * \ | / \ 42 : * \ | / \ 43 : * \|/ xi (out of page) 44 : * o 45 : * 1 46 : * \endverbatim 47 : * 48 : * (xi, eta, zeta): { 0 <= xi <= 1 49 : * { 0 <= eta <= 1 50 : * { 0 <= zeta <= 1 51 : * { xi + eta + zeta <= 1 52 : * are the reference element coordinates associated with the given 53 : * numbering. 54 : * 55 : * \author Benjamin S. Kirk 56 : * \date 2002 57 : * \brief A 3D tetrahedral element with 4 nodes. 58 : */ 59 : class Tet4 final : public Tet 60 : { 61 : public: 62 : 63 : /** 64 : * Constructor. By default this element has no parent. 65 : */ 66 : explicit 67 14819362 : Tet4 (Elem * p=nullptr) : 68 14819362 : Tet(num_nodes, p, _nodelinks_data) 69 14819362 : {} 70 : 71 : Tet4 (Tet4 &&) = delete; 72 : Tet4 (const Tet4 &) = delete; 73 : Tet4 & operator= (const Tet4 &) = delete; 74 : Tet4 & operator= (Tet4 &&) = delete; 75 15669437 : virtual ~Tet4() = default; 76 : 77 : /** 78 : * \returns \p TET4. 79 : */ 80 189028297 : virtual ElemType type () const override { return TET4; } 81 : 82 : /** 83 : * \returns 4. 84 : */ 85 131548740 : virtual unsigned int n_nodes() const override { return num_nodes; } 86 : 87 : /** 88 : * \returns 1. 89 : */ 90 8685 : virtual unsigned int n_sub_elem() const override { return 1; } 91 : 92 : /** 93 : * \returns \p true if the specified (local) node number is a vertex. 94 : */ 95 : virtual bool is_vertex(const unsigned int i) const override; 96 : 97 : /** 98 : * \returns \p true if the specified (local) node number is an edge. 99 : */ 100 : virtual bool is_edge(const unsigned int i) const override; 101 : 102 : /** 103 : * \returns \p true if the specified (local) node number is a face. 104 : */ 105 : virtual bool is_face(const unsigned int i) const override; 106 : 107 : /** 108 : * \returns \p true if the specified (local) node number is on the 109 : * specified side. 110 : */ 111 : virtual bool is_node_on_side(const unsigned int n, 112 : const unsigned int s) const override; 113 : 114 : virtual std::vector<unsigned int> nodes_on_side(const unsigned int s) const override; 115 : 116 : virtual std::vector<unsigned int> nodes_on_edge(const unsigned int e) const override; 117 : 118 : /** 119 : * \returns \p true if the specified (local) node number is on the 120 : * specified edge. 121 : */ 122 : virtual bool is_node_on_edge(const unsigned int n, 123 : const unsigned int e) const override; 124 : 125 : /** 126 : * \returns \p true if the specified child is on the 127 : * specified side. 128 : */ 129 : virtual bool is_child_on_side(const unsigned int c, 130 : const unsigned int s) const override; 131 : 132 : /** 133 : * \returns \p true if the element map is definitely affine within 134 : * numerical tolerances. 135 : */ 136 2106316 : virtual bool has_affine_map () const override { return true; } 137 : 138 : /** 139 : * \returns \p true if the element has non-zero volume(), false otherwise. 140 : */ 141 : virtual bool has_invertible_map(Real tol) const override; 142 : 143 : /** 144 : * \returns \p true if the Lagrange shape functions on this element 145 : * are linear. 146 : */ 147 156914 : virtual bool is_linear () const override { return true; } 148 : 149 : /** 150 : * \returns FIRST. 151 : */ 152 : virtual Order default_order() const override; 153 : 154 : /** 155 : * Builds a \p TRI3 built coincident with face i. 156 : * The \p std::unique_ptr<Elem> handles the memory aspect. 157 : */ 158 : virtual std::unique_ptr<Elem> build_side_ptr (const unsigned int i) override; 159 : 160 : /** 161 : * Rebuilds a TRI3 built coincident with face i. 162 : */ 163 : virtual void build_side_ptr (std::unique_ptr<Elem> & elem, 164 : const unsigned int i) override; 165 : 166 : // Avoid hiding deprecated version with different signature 167 : using Elem::build_side_ptr; 168 : 169 : /** 170 : * Builds a \p EDGE2 built coincident with face i. 171 : * The \p std::unique_ptr<Elem> handles the memory aspect. 172 : */ 173 : virtual std::unique_ptr<Elem> build_edge_ptr (const unsigned int i) override; 174 : 175 : /** 176 : * Rebuilds a \p EDGE2 coincident with edge i. 177 : */ 178 : virtual void build_edge_ptr (std::unique_ptr<Elem> & edge, const unsigned int i) override; 179 : 180 : virtual void connectivity(const unsigned int sc, 181 : const IOPackage iop, 182 : std::vector<dof_id_type> & conn) const override; 183 : 184 : /** 185 : * Geometric constants for Tet4. 186 : */ 187 : static const int num_nodes = 4; 188 : static const int nodes_per_side = 3; 189 : static const int nodes_per_edge = 2; 190 : 191 : /** 192 : * This maps the \f$ j^{th} \f$ node of the \f$ i^{th} \f$ side to 193 : * element node numbers. 194 : */ 195 : static const unsigned int side_nodes_map[num_sides][nodes_per_side]; 196 : 197 : /** 198 : * This maps the \f$ j^{th} \f$ node of the \f$ i^{th} \f$ edge to 199 : * element node numbers. 200 : */ 201 : static const unsigned int edge_nodes_map[num_edges][nodes_per_edge]; 202 : 203 : /** 204 : * The centroid of a 4-node tetrahedron is simply given by the 205 : * average of its vertex positions. 206 : */ 207 : virtual Point true_centroid () const override; 208 : 209 : /** 210 : * An optimized method for computing the area of a 211 : * 4-node tetrahedron. 212 : */ 213 : virtual Real volume () const override; 214 : 215 : /** 216 : * \returns The min and max *dihedral* angles for the tetrahedron. 217 : * 218 : * \note There are 6 dihedral angles (angles between the planar 219 : * faces) for the Tet4. Dihedral angles near 180 deg. are generally 220 : * bad for interpolation. Small dihedral angles are not necessarily 221 : * bad for interpolation, but they can affect the stiffness matrix 222 : * condition number. 223 : */ 224 : std::pair<Real, Real> min_and_max_angle() const; 225 : 226 : /** 227 : * Don't hide Tet::key(side) defined in the base class. 228 : */ 229 : using Tet::key; 230 : 231 : /** 232 : * \returns An id associated with the global node ids of this 233 : * element. The id is not necessarily unique, but should be 234 : * close. 235 : */ 236 : virtual dof_id_type key () const override; 237 : 238 : /** 239 : * Uses simple geometric tests to determine if the point p is inside 240 : * the tetrahedron. 241 : */ 242 : virtual bool contains_point (const Point & p, Real tol) const override; 243 : 244 : virtual void permute(unsigned int perm_num) override final; 245 : 246 : virtual void flip(BoundaryInfo *) override final; 247 : 248 : ElemType side_type (const unsigned int s) const override final; 249 : 250 : protected: 251 : 252 : /** 253 : * Data for links to nodes. 254 : */ 255 : Node * _nodelinks_data[num_nodes]; 256 : 257 : 258 : 259 : #ifdef LIBMESH_ENABLE_AMR 260 : 261 : /** 262 : * Matrix used to create the elements children. 263 : */ 264 : virtual Real embedding_matrix (const unsigned int i, 265 : const unsigned int j, 266 : const unsigned int k) const override; 267 : 268 : /** 269 : * Matrix that computes new nodal locations/solution values 270 : * from current nodes/solution. 271 : */ 272 : static const Real _embedding_matrix[num_children][num_nodes][num_nodes]; 273 : 274 7014604 : LIBMESH_ENABLE_TOPOLOGY_CACHES; 275 : 276 : #endif 277 : 278 : }; 279 : 280 : } // namespace libMesh 281 : 282 : 283 : #endif // LIBMESH_CELL_TET4_H