libMesh
face_polygon.h
Go to the documentation of this file.
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 
39 class Polygon : public Face
40 {
41 public:
42 
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  virtual ~Polygon() = default;
57 
62  virtual bool runtime_topology() const override { return true; }
63 
71  virtual Point master_point (const unsigned int i) const override;
72 
73  static const int num_children = 0; // Refinement not yet supported
74 
78  virtual unsigned int n_nodes() const override { return this->_nodelinks_data.size(); }
79 
84  virtual unsigned int n_sides() const override final { return _elemlinks_data.size()-2; }
85 
89  virtual unsigned int n_vertices() const override final { return _elemlinks_data.size()-2; }
90 
94  virtual unsigned int n_edges() const override final { return _elemlinks_data.size()-2; }
95 
100  virtual unsigned int n_nodes_per_side() const = 0;
101 
105  virtual unsigned int n_children() const override final { return num_children; }
106 
114  virtual bool is_child_on_side(const unsigned int c,
115  const unsigned int s) const override;
116 
121  virtual unsigned int opposite_side(const unsigned int s) const override final;
122 
126  using Elem::key;
127 
133  virtual dof_id_type key () const override;
134 
139  virtual dof_id_type key (const unsigned int s) const override;
140 
147  virtual dof_id_type low_order_key (const unsigned int s) const override;
148 
156  virtual unsigned int local_side_node(unsigned int side,
157  unsigned int side_node) const override;
158 
164  virtual unsigned int local_edge_node(unsigned int edge,
165  unsigned int edge_node) const override;
166 
170  virtual std::unique_ptr<Elem> side_ptr (const unsigned int i) override final;
171 
175  virtual void side_ptr (std::unique_ptr<Elem> & elem,
176  const unsigned int i) override final;
177 
182  virtual std::pair<Real, Real> qual_bounds (const ElemQuality q) const override;
183 
188  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 
207  virtual void retriangulate() = 0;
208 
212  unsigned int n_subtriangles() const { return cast_int<unsigned int>(this->_triangulation.size()); }
213 
217  virtual std::array<int, 3> subtriangle (unsigned int i) const
218  {
219  libmesh_assert_less(i, this->_triangulation.size());
220  return this->_triangulation[i];
221  }
222 
226  virtual std::array<Point, 3> master_subtriangle (unsigned int i) const;
227 
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 
252  std::vector<Elem *> _elemlinks_data;
253 
258  std::vector<Node *> _nodelinks_data;
259 
269  std::vector<std::array<int, 3>> _triangulation;
270 };
271 
272 } // namespace libMesh
273 
274 #endif // LIBMESH_FACE_POLYGON_H
virtual unsigned int n_nodes() const override
Definition: face_polygon.h:78
static const int num_children
Definition: face_polygon.h:73
std::vector< Node * > _nodelinks_data
Data for links to nodes.
Definition: face_polygon.h:258
virtual void retriangulate()=0
Create a triangulation from the current node locations.
The Face is an abstract element type that lives in two dimensions.
Definition: face.h:37
unsigned int n_subtriangles() const
Definition: face_polygon.h:212
virtual dof_id_type key() const override
Definition: face_polygon.C:200
virtual dof_id_type key() const
Definition: elem.C:753
std::tuple< unsigned int, Real, Real > subtriangle_coordinates(const Point &p, Real tol=TOLERANCE *TOLERANCE) const
Definition: face_polygon.C:356
virtual bool is_flipped() const override final
Definition: face_polygon.C:263
static constexpr Real TOLERANCE
std::vector< std::array< int, 3 > > _triangulation
Data for a triangulation of the polygon.
Definition: face_polygon.h:269
virtual unsigned int n_nodes_per_side() const =0
virtual Point master_point(const unsigned int i) const override
Definition: face_polygon.C:80
virtual bool on_reference_element(const Point &p, const Real eps=TOLERANCE) const override final
Definition: face_polygon.C:109
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
virtual std::unique_ptr< Elem > side_ptr(const unsigned int i) override final
Definition: face_polygon.C:211
virtual unsigned int n_children() const override final
Definition: face_polygon.h:105
virtual unsigned int opposite_side(const unsigned int s) const override final
Definition: face_polygon.C:252
The libMesh namespace provides an interface to certain functionality in the library.
virtual unsigned int local_edge_node(unsigned int edge, unsigned int edge_node) const override
Calls local_side_node(edge, edge_node).
Definition: face_polygon.C:192
Polygon & operator=(const Polygon &)=delete
virtual std::array< int, 3 > subtriangle(unsigned int i) const
Definition: face_polygon.h:217
virtual ~Polygon()=default
virtual unsigned int n_edges() const override final
Definition: face_polygon.h:94
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const override
Definition: face_polygon.C:243
virtual unsigned int n_permutations() const override
n_sides sides, one orientation each.
Definition: face_polygon.h:188
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const override
Definition: face_polygon.C:174
The Polygon is an element in 2D with an arbitrary (but fixed) number of sides.
Definition: face_polygon.h:39
virtual std::array< Point, 3 > master_subtriangle(unsigned int i) const
Definition: face_polygon.C:343
Polygon(const unsigned int nn, const unsigned int ns, Elem *p)
Arbitrary polygonal element, takes number of nodes and sides, and a (probably unused) parent link...
Definition: face_polygon.C:46
virtual bool runtime_topology() const override
Definition: face_polygon.h:62
virtual std::pair< Real, Real > qual_bounds(const ElemQuality q) const override
Definition: face_polygon.C:306
virtual std::vector< unsigned int > edges_adjacent_to_node(const unsigned int n) const override
Definition: face_polygon.C:288
ElemQuality
Defines an enum for element quality metrics.
std::vector< Elem * > _elemlinks_data
Data for links to parent/neighbor/interior_parent elements.
Definition: face_polygon.h:252
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual unsigned int n_vertices() const override final
Definition: face_polygon.h:89
virtual unsigned int n_sides() const override final
Definition: face_polygon.h:84
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
virtual dof_id_type low_order_key(const unsigned int s) const override
Definition: face_polygon.C:163
uint8_t dof_id_type
Definition: id_types.h:67