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_PRISM6_H 21 : #define LIBMESH_CELL_PRISM6_H 22 : 23 : // Local includes 24 : #include "libmesh/cell_prism.h" 25 : 26 : namespace libMesh 27 : { 28 : 29 : /** 30 : * The \p Prism6 is an element in 3D composed of 6 nodes. It is 31 : * numbered like this: 32 : * 33 : * \verbatim 34 : * PRISM6: 35 : * 5 36 : * o 37 : * /:\ 38 : * / : \ zeta 39 : * / o \ ^ eta (into page) 40 : * 3 o-------o 4 | / 41 : * | . 2 . | |/ 42 : * |. .| o---> xi 43 : * o-------o 44 : * 0 1 45 : * \endverbatim 46 : * 47 : * (xi, eta, zeta): { 0 <= xi <= 1 48 : * { 0 <= eta <= 1 49 : * { -1 <= zeta <= 1 50 : * { xi + eta <= 1 51 : * are the reference element coordinates associated with the given 52 : * numbering. 53 : * 54 : * \author Benjamin S. Kirk 55 : * \date 2002 56 : * \brief A 3D prismatic element with 6 nodes. 57 : */ 58 : class Prism6 final : public Prism 59 : { 60 : public: 61 : 62 : /** 63 : * Constructor. By default this element has no parent. 64 : */ 65 : explicit 66 271273 : Prism6 (Elem * p=nullptr) : 67 271273 : Prism(num_nodes, p, _nodelinks_data) 68 271273 : {} 69 : 70 : Prism6 (Prism6 &&) = delete; 71 : Prism6 (const Prism6 &) = delete; 72 : Prism6 & operator= (const Prism6 &) = delete; 73 : Prism6 & operator= (Prism6 &&) = delete; 74 319096 : virtual ~Prism6() = default; 75 : 76 : /** 77 : * \returns \p PRISM6. 78 : */ 79 248555164 : virtual ElemType type () const override { return PRISM6; } 80 : 81 : /** 82 : * \returns 1. 83 : */ 84 184320 : virtual unsigned int n_sub_elem() const override { return 1; } 85 : 86 : /** 87 : * \returns \p true if the specified (local) node number is a vertex. 88 : */ 89 : virtual bool is_vertex(const unsigned int i) const override; 90 : 91 : /** 92 : * \returns \p true if the specified (local) node number is an edge. 93 : */ 94 : virtual bool is_edge(const unsigned int i) const override; 95 : 96 : /** 97 : * \returns \p true if the specified (local) node number is a face. 98 : */ 99 : virtual bool is_face(const unsigned int i) const override; 100 : 101 : /** 102 : * \returns \p true if the specified (local) node number is on the 103 : * specified side. 104 : */ 105 : virtual bool is_node_on_side(const unsigned int n, 106 : const unsigned int s) const override; 107 : 108 : virtual std::vector<unsigned int> nodes_on_side(const unsigned int s) const override; 109 : 110 : virtual std::vector<unsigned int> nodes_on_edge(const unsigned int e) const override; 111 : 112 : /** 113 : * \returns \p true if the specified (local) node number is on the 114 : * specified edge. 115 : */ 116 : virtual bool is_node_on_edge(const unsigned int n, 117 : const unsigned int e) const override; 118 : 119 : /** 120 : * \returns \p true if the element map is definitely affine within 121 : * numerical tolerances. 122 : */ 123 : virtual bool has_affine_map () const override; 124 : 125 : /** 126 : * \returns FIRST. 127 : */ 128 : virtual Order default_order() const override; 129 : 130 : /** 131 : * Builds a \p QUAD4 or \p TRI3 built coincident with face i. 132 : * The \p std::unique_ptr<Elem> handles the memory aspect. 133 : */ 134 : virtual std::unique_ptr<Elem> build_side_ptr (const unsigned int i) override; 135 : 136 : /** 137 : * Rebuilds a \p QUAD4 or \p TRI3 built coincident with face i. 138 : */ 139 : virtual void build_side_ptr (std::unique_ptr<Elem> & elem, 140 : const unsigned int i) override; 141 : 142 : // Avoid hiding deprecated version with different signature 143 : using Elem::build_side_ptr; 144 : 145 : /** 146 : * Builds a \p EDGE2 or \p INFEDGE2 built coincident with face i. 147 : * The \p std::unique_ptr<Elem> handles the memory aspect. 148 : */ 149 : virtual std::unique_ptr<Elem> build_edge_ptr (const unsigned int i) override; 150 : 151 : /** 152 : * Rebuilds a \p EDGE2 built coincident with edges 0 to 2, or \p 153 : * INFEDGE2 built coincident with edges 3 to 5. 154 : */ 155 : virtual void build_edge_ptr (std::unique_ptr<Elem> & edge, const unsigned int i) override; 156 : 157 : virtual void connectivity(const unsigned int sc, 158 : const IOPackage iop, 159 : std::vector<dof_id_type> & conn) const override; 160 : 161 : /** 162 : * Geometric constants for Prism6. 163 : */ 164 : static const int num_nodes = 6; 165 : static const int nodes_per_side = 4; 166 : static const int nodes_per_edge = 2; 167 : 168 : /** 169 : * This maps the \f$ j^{th} \f$ node of the \f$ i^{th} \f$ side to 170 : * element node numbers. 171 : */ 172 : static const unsigned int side_nodes_map[num_sides][nodes_per_side]; 173 : 174 : /** 175 : * This maps the child elements with the associated side of the parent element 176 : */ 177 : static const unsigned int side_elems_map[num_sides][nodes_per_side]; 178 : 179 : /** 180 : * This maps the \f$ j^{th} \f$ node of the \f$ i^{th} \f$ edge to 181 : * element node numbers. 182 : */ 183 : static const unsigned int edge_nodes_map[num_edges][nodes_per_edge]; 184 : 185 : /** 186 : * An Optimized numerical quadrature approach for computing the 187 : * centroid of the Prism6. 188 : */ 189 : virtual Point true_centroid () const override; 190 : 191 : /** 192 : * Specialized function for computing the element volume. 193 : */ 194 : virtual Real volume () const override; 195 : 196 : /** 197 : * Builds a bounding box out of the nodal positions 198 : */ 199 : virtual BoundingBox loose_bounding_box () const override; 200 : 201 : virtual void permute(unsigned int perm_num) override final; 202 : 203 : virtual void flip(BoundaryInfo *) override final; 204 : 205 : ElemType side_type (const unsigned int s) const override final; 206 : 207 : protected: 208 : 209 : /** 210 : * Data for links to nodes. 211 : */ 212 : Node * _nodelinks_data[num_nodes]; 213 : 214 : 215 : 216 : #ifdef LIBMESH_ENABLE_AMR 217 : 218 : /** 219 : * Matrix used to create the elements children. 220 : */ 221 491922 : virtual Real embedding_matrix (const unsigned int i, 222 : const unsigned int j, 223 : const unsigned int k) const override 224 491922 : { return _embedding_matrix[i][j][k]; } 225 : 226 : /** 227 : * Matrix that computes new nodal locations/solution values 228 : * from current nodes/solution. 229 : */ 230 : static const Real _embedding_matrix[num_children][num_nodes][num_nodes]; 231 : 232 1655574 : LIBMESH_ENABLE_TOPOLOGY_CACHES; 233 : 234 : #endif // LIBMESH_ENABLE_AMR 235 : 236 : }; 237 : 238 : } // namespace libMesh 239 : 240 : #endif // LIBMESH_CELL_PRISM6_H