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_FACE_TRI7_H 21 : #define LIBMESH_FACE_TRI7_H 22 : 23 : // Local includes 24 : #include "libmesh/libmesh_common.h" 25 : #include "libmesh/face_tri.h" 26 : 27 : namespace libMesh 28 : { 29 : 30 : /** 31 : * The \p Tri7 is an element in 2D composed of 7 nodes. It is 32 : * numbered like this: 33 : * 34 : * \verbatim 35 : * TRI7: 36 : * 2 37 : * o 38 : * |\ 39 : * | \ 40 : * | \ eta 41 : * | \ ^ 42 : * | \ | 43 : * 5 o o 4 | 44 : * | o \ o---> xi 45 : * | 6 \ 46 : * | \ 47 : * o----o----o 48 : * 0 3 1 49 : * \endverbatim 50 : * 51 : * (xi, eta): { 0 <= xi <= 1 52 : * { 0 <= eta <= 1 53 : * { xi + eta <= 1 54 : * are the reference element coordinates associated with the given 55 : * numbering. 56 : * 57 : * \author Roy H Stogner 58 : * \date 2021 59 : * \brief A 2D triangular element with 7 nodes. 60 : */ 61 : class Tri7 : public Tri 62 : { 63 : public: 64 : 65 : /** 66 : * Constructor. By default this element has no parent. 67 : */ 68 : explicit 69 19961433 : Tri7 (Elem * p=nullptr) : 70 19961433 : Tri(num_nodes, p, _nodelinks_data) {} 71 : 72 : Tri7 (Tri7 &&) = delete; 73 : Tri7 (const Tri7 &) = delete; 74 : Tri7 & operator= (const Tri7 &) = delete; 75 : Tri7 & operator= (Tri7 &&) = delete; 76 21147512 : virtual ~Tri7() = default; 77 : 78 : /** 79 : * \returns \p TRI7. 80 : */ 81 1750959351 : virtual ElemType type () const override { return TRI7; } 82 : 83 : /** 84 : * \returns 7. 85 : */ 86 151909013 : virtual unsigned int n_nodes() const override { return num_nodes; } 87 : 88 : /** 89 : * \returns 4. 90 : */ 91 0 : virtual unsigned int n_sub_elem() const override { return 4; } 92 : 93 : /** 94 : * \returns \p true if the specified (local) node number is a vertex. 95 : */ 96 : virtual bool is_vertex(const unsigned int i) const override; 97 : 98 : /** 99 : * \returns \p true if the specified (local) node number is an edge. 100 : */ 101 : virtual bool is_edge(const unsigned int i) const override; 102 : 103 : /** 104 : * \returns \p true if the specified (local) node number is a face. 105 : */ 106 : virtual bool is_face(const unsigned int i) const override; 107 : 108 : /** 109 : * \returns \p true if the specified (local) node number is on the 110 : * specified side. 111 : */ 112 : virtual bool is_node_on_side(const unsigned int n, 113 : const unsigned int s) const override; 114 : 115 : virtual std::vector<unsigned int> nodes_on_side(const unsigned int s) const override; 116 : 117 : virtual std::vector<unsigned int> nodes_on_edge(const unsigned int e) const override; 118 : 119 : /** 120 : * \returns \p true if the specified (local) node number is on the 121 : * specified edge (== is_node_on_side in 2D). 122 : */ 123 792 : virtual bool is_node_on_edge(const unsigned int n, 124 : const unsigned int e) const override 125 792 : { return this->is_node_on_side(n,e); } 126 : 127 : /** 128 : * \returns \p true if the element map is definitely affine within 129 : * numerical tolerances. 130 : */ 131 : virtual bool has_affine_map () const override; 132 : 133 : /** 134 : * \returns THIRD. 135 : */ 136 : virtual Order default_order() const override; 137 : 138 : /** 139 : * \returns SECOND - the sides are EDGE3 140 : */ 141 : virtual Order default_side_order() const override; 142 : 143 : /** 144 : * Don't hide Elem::key() defined in the base class. 145 : */ 146 : using Elem::key; 147 : 148 : /** 149 : * \returns An id associated with the \p s side of this element. 150 : * The id is not necessarily unique, but should be close. 151 : * 152 : * We reimplement this method here for the \p Tri7 since we can 153 : * use the center node of each edge to provide a perfect (unique) 154 : * key. 155 : */ 156 : virtual dof_id_type key (const unsigned int s) const override; 157 : 158 : /** 159 : * \returns \p Tri7::side_nodes_map[side][side_node] after doing some range checking. 160 : */ 161 : virtual unsigned int local_side_node(unsigned int side, 162 : unsigned int side_node) const override; 163 : 164 : virtual std::unique_ptr<Elem> build_side_ptr (const unsigned int i) override; 165 : 166 : /** 167 : * Rebuilds an EDGE2 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 : virtual void connectivity(const unsigned int sf, 176 : const IOPackage iop, 177 : std::vector<dof_id_type> & conn) const override; 178 : 179 : /** 180 : * \returns 2 for edge nodes and 3 for the interior node. 181 : */ 182 : virtual unsigned int n_second_order_adjacent_vertices (const unsigned int) const override; 183 : 184 : /** 185 : * \returns The element-local number of the \f$ v^{th} \f$ vertex 186 : * that defines the \f$ n^{th} \f$ second-order node. 187 : * 188 : * \note \p n is counted as depicted above, \f$ 3 \le n < 7 \f$. 189 : */ 190 : virtual unsigned short int second_order_adjacent_vertex (const unsigned int n, 191 : const unsigned int v) const override; 192 : 193 : /** 194 : * \returns The child number \p c and element-local index \p v of the 195 : * \f$ n^{th} \f$ second-order node on the parent element. See 196 : * elem.h for further details. 197 : */ 198 : virtual std::pair<unsigned short int, unsigned short int> 199 : second_order_child_vertex (const unsigned int n) const override; 200 : 201 : /** 202 : * Geometric constants for Tri7. 203 : */ 204 : static const int num_nodes = 7; 205 : static const int nodes_per_side = 3; 206 : 207 : /** 208 : * This maps the \f$ j^{th} \f$ node of the \f$ i^{th} \f$ side to 209 : * element node numbers. 210 : */ 211 : static const unsigned int side_nodes_map[num_sides][nodes_per_side]; 212 : 213 : /** 214 : * \returns A bounding box (not necessarily the minimal bounding box) 215 : * containing the geometric element. 216 : */ 217 : virtual BoundingBox loose_bounding_box () const override; 218 : 219 : virtual void permute(unsigned int perm_num) override final; 220 : 221 : virtual void flip(BoundaryInfo *) override final; 222 : 223 : #ifdef LIBMESH_ENABLE_AMR 224 : virtual 225 : const std::vector<std::pair<unsigned char, unsigned char>> & 226 72805 : parent_bracketing_nodes(unsigned int c, 227 : unsigned int n) const override 228 72805 : { return _parent_bracketing_nodes[c][n]; } 229 : #endif 230 : 231 : unsigned int center_node_on_side(const unsigned short side) const override final; 232 : 233 : ElemType side_type (const unsigned int s) const override final; 234 : 235 : protected: 236 : 237 : /** 238 : * Data for links to nodes. 239 : */ 240 : Node * _nodelinks_data[num_nodes]; 241 : 242 : 243 : 244 : #ifdef LIBMESH_ENABLE_AMR 245 : 246 : /** 247 : * Matrix used to create the elements children. 248 : */ 249 217122 : virtual Real embedding_matrix (const unsigned int i, 250 : const unsigned int j, 251 : const unsigned int k) const override 252 217122 : { return _embedding_matrix[i][j][k]; } 253 : 254 : /** 255 : * Matrix that computes new nodal locations/solution values 256 : * from current nodes/solution. 257 : */ 258 : static const Real _embedding_matrix[num_children][num_nodes][num_nodes]; 259 : 260 : /** 261 : * Pairs of nodes that bracket child nodes when doing mesh 262 : * refinement. 263 : */ 264 : static const std::vector<std::pair<unsigned char, unsigned char>> 265 : _parent_bracketing_nodes[num_children][num_nodes]; 266 : 267 112336 : 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[num_sides][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[num_nodes]; 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[num_nodes]; 288 : }; 289 : 290 : 291 : } // namespace libMesh 292 : 293 : #endif // LIBMESH_FACE_TRI6_H