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_EDGE_EDGE2_H 21 : #define LIBMESH_EDGE_EDGE2_H 22 : 23 : // Local includes 24 : #include "libmesh/libmesh_common.h" 25 : #include "libmesh/edge.h" 26 : 27 : namespace libMesh 28 : { 29 : 30 : /** 31 : * The \p Edge2 is an element in 1D composed of 2 nodes. It is numbered 32 : * like this: 33 : * 34 : * \verbatim 35 : * EDGE2: o--------o o---> xi 36 : * 0 1 37 : * \endverbatim 38 : * 39 : * xi in [-1,1] is the reference element coordinate associated with 40 : * the given numbering. 41 : * 42 : * \author Benjamin S. Kirk 43 : * \date 2002 44 : * \brief A 1D geometric element with 2 nodes. 45 : */ 46 : class Edge2 : public Edge 47 : { 48 : public: 49 : 50 : /** 51 : * Constructor. By default this element has no parent. 52 : */ 53 : explicit 54 38960989 : Edge2 (Elem * p=nullptr) : 55 38960989 : Edge(num_nodes, p, _nodelinks_data) {} 56 : 57 : Edge2 (Edge2 &&) = delete; 58 : Edge2 (const Edge2 &) = delete; 59 : Edge2 & operator= (const Edge2 &) = delete; 60 : Edge2 & operator= (Edge2 &&) = delete; 61 61803993 : virtual ~Edge2() = default; 62 : 63 : /** 64 : * \returns The \p Point associated with local \p Node \p i, 65 : * in master element rather than physical coordinates. 66 : */ 67 1372 : virtual Point master_point (const unsigned int i) const override 68 : { 69 56 : libmesh_assert_less(i, this->n_nodes()); 70 1428 : return Point(2.0f*Real(i)-1.0f,0,0); 71 : } 72 : 73 : /** 74 : * \returns 1. 75 : */ 76 0 : virtual unsigned int n_sub_elem() const override { return 1; } 77 : 78 : /** 79 : * \returns \p true if the specified (local) node number is a vertex. 80 : */ 81 : virtual bool is_vertex(const unsigned int i) const override; 82 : 83 : /** 84 : * \returns \p true if the specified (local) node number is an edge. 85 : */ 86 : virtual bool is_edge(const unsigned int i) const override; 87 : 88 : /** 89 : * \returns \p true if the specified (local) node number is a face. 90 : */ 91 : virtual bool is_face(const unsigned int i) const override; 92 : 93 : /** 94 : * \returns \p true if the specified (local) node number is on the 95 : * specified side. 96 : */ 97 : virtual bool is_node_on_side(const unsigned int n, 98 : const unsigned int s) const override; 99 : 100 : /** 101 : * \returns \p true if the specified (local) node number is on the 102 : * specified edge (always true in 1D). 103 : */ 104 : virtual bool is_node_on_edge(const unsigned int n, 105 : const unsigned int e) const override; 106 : 107 : /** 108 : * \returns \p true if the element map is definitely affine within 109 : * numerical tolerances. 110 : */ 111 44749659 : virtual bool has_affine_map () const override { return true; } 112 : 113 : /** 114 : * \returns \p true if the element has non-zero volume(), false otherwise. 115 : */ 116 : virtual bool has_invertible_map(Real tol) const override; 117 : 118 : /** 119 : * \returns \p true if the Lagrange shape functions on this element 120 : * are linear. 121 : */ 122 2560788 : virtual bool is_linear () const override { return true; } 123 : 124 : /** 125 : * \returns \p EDGE2. 126 : */ 127 897238740 : virtual ElemType type() const override { return EDGE2; } 128 : 129 : /** 130 : * \returns FIRST. 131 : */ 132 : virtual Order default_order() const override; 133 : 134 : virtual void connectivity(const unsigned int sc, 135 : const IOPackage iop, 136 : std::vector<dof_id_type> & conn) const override; 137 : 138 : /** 139 : * An optimized method for computing the centroid of a 140 : * 2-node edge. 141 : */ 142 : virtual Point true_centroid () const override; 143 : 144 : /** 145 : * An optimized method for computing the length of a 2-node edge. 146 : */ 147 : virtual Real volume () const override; 148 : 149 : #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 150 : 151 : /** 152 : * \returns \p false. This is a finite element. 153 : */ 154 13552682 : virtual bool infinite () const override { return false; } 155 : 156 : #endif 157 : 158 : /** 159 : * Don't hide Edge::key(side) defined in the base class. 160 : */ 161 : using Edge::key; 162 : 163 : /** 164 : * \returns An id associated with the global node ids of this 165 : * element. The id is not necessarily unique, but should be 166 : * close. 167 : */ 168 : virtual dof_id_type key () const override; 169 : 170 : /** 171 : * Geometric constants for Edge2. 172 : */ 173 : static const int num_nodes = 2; 174 : 175 : virtual void flip(BoundaryInfo *) override final; 176 : 177 : protected: 178 : 179 : /** 180 : * Data for links to nodes. 181 : */ 182 : Node * _nodelinks_data[num_nodes]; 183 : 184 : 185 : 186 : #ifdef LIBMESH_ENABLE_AMR 187 : 188 : /** 189 : * Matrix used to create the elements children. 190 : */ 191 3814 : virtual Real embedding_matrix (const unsigned int i, 192 : const unsigned int j, 193 : const unsigned int k) const override 194 3814 : { return _embedding_matrix[i][j][k]; } 195 : 196 : /** 197 : * Matrix that computes new nodal locations/solution values 198 : * from current nodes/solution. 199 : */ 200 : static const Real _embedding_matrix[num_children][num_nodes][num_nodes]; 201 : 202 13652 : LIBMESH_ENABLE_TOPOLOGY_CACHES; 203 : 204 : #endif // LIBMESH_ENABLE_AMR 205 : 206 : }; 207 : 208 : } // namespace libMesh 209 : 210 : 211 : #endif // LIBMESH_EDGE_EDGE2_H