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_POLYGON_H
21 : #define LIBMESH_FACE_POLYGON_H
22 :
23 :
24 : // Local includes
25 : #include "libmesh/libmesh_common.h"
26 : #include "libmesh/face.h"
27 :
28 : namespace libMesh
29 : {
30 :
31 : /**
32 : * The \p Polygon is an element in 2D with an arbitrary (but fixed)
33 : * number of sides.
34 : *
35 : * \author Roy H. Stogner
36 : * \date 2025
37 : * \brief A 2D element with an arbitrary number of sides.
38 : */
39 : class Polygon : public Face
40 : {
41 : public:
42 :
43 : /**
44 : * Arbitrary polygonal element, takes number of nodes and
45 : * sides, and a (probably unused) parent link. Derived classes
46 : * implement 'true' elements.
47 : */
48 : Polygon (const unsigned int nn,
49 : const unsigned int ns,
50 : Elem * p);
51 :
52 : Polygon (Polygon &&) = delete;
53 : Polygon (const Polygon &) = delete;
54 : Polygon & operator= (const Polygon &) = delete;
55 : Polygon & operator= (Polygon &&) = delete;
56 864 : virtual ~Polygon() = default;
57 :
58 : /**
59 : * \returns true - polygon subclasses can have numbers of sides and
60 : * nodes which vary at runtime.
61 : */
62 12 : virtual bool runtime_topology() const override { return true; }
63 :
64 : /**
65 : * \returns The \p Point associated with local \p Node \p i,
66 : * in master element rather than physical coordinates.
67 : *
68 : * This implementation returns the master vertices; subclasses with
69 : * more points will need to further override.
70 : */
71 : virtual Point master_point (const unsigned int i) const override;
72 :
73 : static const int num_children = 0; // Refinement not yet supported
74 :
75 : /**
76 : * \returns the number of nodes in the polygon.
77 : */
78 3827146 : virtual unsigned int n_nodes() const override { return this->_nodelinks_data.size(); }
79 :
80 : /**
81 : * \returns the number of sides, which is the number of links we had
82 : * to save minus one parent link and minus one interior_parent link.
83 : */
84 14623613 : virtual unsigned int n_sides() const override final { return _elemlinks_data.size()-2; }
85 :
86 : /**
87 : * \returns the number of vertices. All polygons have 1 vertex per side.
88 : */
89 457728 : virtual unsigned int n_vertices() const override final { return _elemlinks_data.size()-2; }
90 :
91 : /**
92 : * \returns the number of edges. All polygons have 1 edge per side.
93 : */
94 28524 : virtual unsigned int n_edges() const override final { return _elemlinks_data.size()-2; }
95 :
96 : /**
97 : * \returns the number of nodes on each side. All polygons have the
98 : * same number of nodes on each side.
99 : */
100 : virtual unsigned int n_nodes_per_side() const = 0;
101 :
102 : /**
103 : * \returns num_children.
104 : */
105 0 : virtual unsigned int n_children() const override final { return num_children; }
106 :
107 : /**
108 : * \returns \p true if the specified child is on the
109 : * specified side.
110 : *
111 : * Not implemented ... indefinitely? I don't think we'll be doing
112 : * hierarchic refinement for general polygons.
113 : */
114 : virtual bool is_child_on_side(const unsigned int c,
115 : const unsigned int s) const override;
116 :
117 : /**
118 : * \returns The side number opposite to \p s (for n_sides() even),
119 : * or throws an error otherwise.
120 : */
121 : virtual unsigned int opposite_side(const unsigned int s) const override final;
122 :
123 : /**
124 : * Don't hide Elem::key() defined in the base class.
125 : */
126 : using Elem::key;
127 :
128 : /**
129 : * \returns An id associated with the global node ids of this
130 : * element. The id is not necessarily unique, but should be
131 : * close.
132 : */
133 : virtual dof_id_type key () const override;
134 :
135 : /**
136 : * \returns An id associated with the \p s side of this element.
137 : * The id is not necessarily unique, but should be close.
138 : */
139 : virtual dof_id_type key (const unsigned int s) const override;
140 :
141 : /**
142 : * \returns An id associated with the \p s side of this element, as
143 : * defined solely by element vertices. The id is not necessarily
144 : * unique, but should be close. This is particularly useful in the
145 : * \p MeshBase::find_neighbors() routine.
146 : */
147 : virtual dof_id_type low_order_key (const unsigned int s) const override;
148 :
149 : /**
150 : * \returns The local node id for node \p side_node on side \p side of
151 : * this Elem.
152 : *
153 : * This implementation assumes a particular node numbering "style"
154 : * and may be overridden in subclasses.
155 : */
156 : virtual unsigned int local_side_node(unsigned int side,
157 : unsigned int side_node) const override;
158 :
159 : /**
160 : * Calls local_side_node(edge, edge_node). For 2D elements, there is an implied
161 : * equivalence between edges and sides, e.g. n_edges() == n_sides(), so we treat
162 : * these two functions the same.
163 : */
164 : virtual unsigned int local_edge_node(unsigned int edge,
165 : unsigned int edge_node) const override;
166 :
167 : /**
168 : * \returns A primitive (2-noded) edge for edge i.
169 : */
170 : virtual std::unique_ptr<Elem> side_ptr (const unsigned int i) override final;
171 :
172 : /**
173 : * Rebuilds an EDGE2 coincident with face i.
174 : */
175 : virtual void side_ptr (std::unique_ptr<Elem> & elem,
176 : const unsigned int i) override final;
177 :
178 : /**
179 : * \returns The suggested quality bounds for
180 : * the hex based on quality measure q.
181 : */
182 : virtual std::pair<Real, Real> qual_bounds (const ElemQuality q) const override;
183 :
184 : /**
185 : * n_sides sides, one orientation each. Not marked final because
186 : * subclasses may have interior nodes with less symmetry
187 : */
188 625 : virtual unsigned int n_permutations() const override { return this->n_sides(); }
189 :
190 : virtual bool is_flipped() const override final;
191 :
192 : virtual std::vector<unsigned int> edges_adjacent_to_node(const unsigned int n) const override;
193 :
194 : /**
195 : * Create a triangulation from the current node locations.
196 : *
197 : * If this is not called by the user, it will be called
198 : * automatically when needed.
199 : *
200 : * If the user moves the polygon nodes, the triangulation may need
201 : * to be regenerated manually, or the changed mapping may be lower
202 : * quality or even inverted.
203 : *
204 : * Pure virtual because this strongly depends on what mid-edge or
205 : * mid-face nodes the subclass might have.
206 : */
207 : virtual void retriangulate() = 0;
208 :
209 : /**
210 : * \returns the size of the triangulation of the polygon
211 : */
212 512927 : unsigned int n_subtriangles() const { return cast_int<unsigned int>(this->_triangulation.size()); }
213 :
214 : /**
215 : * \returns the local indices of points on a subtriangle of the polygon
216 : */
217 2049188 : virtual std::array<int, 3> subtriangle (unsigned int i) const
218 : {
219 174855 : libmesh_assert_less(i, this->_triangulation.size());
220 2214943 : return this->_triangulation[i];
221 : }
222 :
223 : /**
224 : * \returns the master-space points of a subtriangle of the polygon
225 : */
226 : virtual std::array<Point, 3> master_subtriangle (unsigned int i) const;
227 :
228 : /**
229 : * \returns the index of a subtriangle containing the master-space
230 : * point \p p, along with barycentric coordinates for \p, or return
231 : * invalid_uint if no subtriangle contains p to within tolerance \p
232 : * tol.
233 : */
234 : std::tuple<unsigned int, Real, Real>
235 : subtriangle_coordinates (const Point & p,
236 : Real tol = TOLERANCE*TOLERANCE) const;
237 :
238 : virtual bool on_reference_element(const Point & p,
239 : const Real eps = TOLERANCE) const override final;
240 :
241 : protected:
242 :
243 : /**
244 : * Data for links to parent/neighbor/interior_parent elements.
245 : *
246 : * There should be num_sides neighbors, preceded by any parent link
247 : * and followed by any interior_parent link.
248 : *
249 : * We're stuck with a heap allocation here since the number of sides
250 : * in a subclass is chosen at runtime.
251 : */
252 : std::vector<Elem *> _elemlinks_data;
253 :
254 : /**
255 : * Data for links to nodes. We're stuck with a heap allocation here
256 : * since the number of nodes in a subclass is chosen at runtime.
257 : */
258 : std::vector<Node *> _nodelinks_data;
259 :
260 : /**
261 : * Data for a triangulation of the polygon. Our mapping from master
262 : * space to physical space will be based on these subelements.
263 : *
264 : * Positive int values should correspond to local node numbers.
265 : *
266 : * Negative int values may correspond to "special" points, e.g. a
267 : * centroid or skeleton point.
268 : */
269 : std::vector<std::array<int, 3>> _triangulation;
270 : };
271 :
272 : } // namespace libMesh
273 :
274 : #endif // LIBMESH_FACE_POLYGON_H
|