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_INF_EDGE2_H 21 : #define LIBMESH_EDGE_INF_EDGE2_H 22 : 23 : #include "libmesh/libmesh_common.h" 24 : 25 : #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 26 : 27 : // Local includes 28 : #include "libmesh/edge.h" 29 : 30 : namespace libMesh 31 : { 32 : 33 : /** 34 : * The \p InfEdge2 is an infinite element in 1D composed of 2 nodes. 35 : * It is numbered like this: 36 : * 37 : * \verbatim 38 : * INFEDGE2: 39 : * o closer to infinity 40 : * | 1 41 : * | 42 : * | 43 : * | 44 : * o base node 45 : * 0 46 : * \endverbatim 47 : * 48 : * \author Daniel Dreyer 49 : * \date 2003 50 : * \brief A 1D infinite element with 2 nodes. 51 : */ 52 : class InfEdge2 : public Edge 53 : { 54 : public: 55 : 56 : /** 57 : * Constructor. By default this element has no parent. 58 : */ 59 : explicit 60 21754 : InfEdge2 (Elem * p=nullptr) : 61 21754 : Edge(num_nodes, p, _nodelinks_data) {} 62 : 63 : InfEdge2 (InfEdge2 &&) = delete; 64 : InfEdge2 (const InfEdge2 &) = delete; 65 : InfEdge2 & operator= (const InfEdge2 &) = delete; 66 : InfEdge2 & operator= (InfEdge2 &&) = delete; 67 40186 : virtual ~InfEdge2() = default; 68 : 69 : /** 70 : * Geometric constants for InfEdge2. 71 : */ 72 : static const int num_nodes = 2; 73 : 74 : /** 75 : * \returns The \p Point associated with local \p Node \p i, 76 : * in master element rather than physical coordinates. 77 : */ 78 0 : virtual Point master_point (const unsigned int i) const override 79 : { 80 0 : libmesh_assert_less(i, this->n_nodes()); 81 0 : return Point(0,i,0); 82 : } 83 : 84 : /** 85 : * \returns 1. 86 : */ 87 0 : virtual unsigned int n_sub_elem() const override { return 1; } 88 : 89 : /** 90 : * \returns \p true if the specified (local) node number is a vertex. 91 : */ 92 : virtual bool is_vertex(const unsigned int i) const override; 93 : 94 : /** 95 : * \returns \p true if the specified (local) node number is an edge. 96 : */ 97 : virtual bool is_edge(const unsigned int i) const override; 98 : 99 : /** 100 : * \returns \p true if the specified (local) node number is a face. 101 : */ 102 : virtual bool is_face(const unsigned int i) const override; 103 : 104 : /** 105 : * \returns \p true if the specified (local) node number is on the 106 : * specified side. 107 : */ 108 : virtual bool is_node_on_side(const unsigned int n, 109 : const unsigned int s) const override; 110 : 111 : virtual std::vector<unsigned int> nodes_on_side(const unsigned int s) const override; 112 : 113 : /** 114 : * \returns \p true if the specified (local) node number is on the 115 : * specified edge (always true in 1D). 116 : */ 117 : virtual bool is_node_on_edge(const unsigned int n, 118 : const unsigned int e) const override; 119 : 120 : /** 121 : * \returns \p INFEDGE2. 122 : */ 123 794 : virtual ElemType type() const override { return INFEDGE2; } 124 : 125 : /** 126 : * \returns FIRST. 127 : */ 128 : virtual Order default_order() const override; 129 : 130 : virtual void connectivity(const unsigned int se, 131 : const IOPackage iop, 132 : std::vector<dof_id_type> & conn) const override; 133 : 134 : /** 135 : * No such thing as a misoriented InfEdge 136 : */ 137 0 : virtual void flip(BoundaryInfo *) override final {}; 138 : 139 : #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 140 : 141 : /** 142 : * \returns \p true. This is an infinite element. 143 : */ 144 582 : virtual bool infinite () const override { return true; } 145 : 146 : /** 147 : * \returns The origin of this infinite element. 148 : */ 149 : virtual Point origin () const override; 150 : 151 : /** 152 : * \returns \p true if the specified (local) node number is a 153 : * "mid-edge" node on an infinite element edge. 154 : */ 155 0 : virtual bool is_mid_infinite_edge_node(const unsigned int i) const 156 0 : override { return (i > 0); } 157 : 158 : #endif 159 : 160 : 161 : protected: 162 : 163 : /** 164 : * Data for links to nodes. 165 : */ 166 : Node * _nodelinks_data[2]; 167 : 168 : 169 : 170 : #ifdef LIBMESH_ENABLE_AMR 171 : 172 : /** 173 : * Matrix used to create the elements children. 174 : */ 175 0 : virtual Real embedding_matrix (const unsigned int, 176 : const unsigned int, 177 : const unsigned int) const override 178 0 : { libmesh_not_implemented(); return 0.; } 179 : 180 0 : LIBMESH_ENABLE_TOPOLOGY_CACHES; 181 : 182 : #endif // LIBMESH_ENABLE_AMR 183 : 184 : }; 185 : 186 : 187 : 188 : 189 : // ------------------------------------------------------------ 190 : // InfEdge2 class member functions 191 : inline 192 0 : Point InfEdge2::origin () const 193 : { 194 0 : return ( this->point(0)*2 - this->point(1) ); 195 : } 196 : 197 : 198 : } // namespace libMesh 199 : 200 : #endif 201 : 202 : #endif // LIBMESH_EDGE_INF_EDGE2_H