libMesh
face_polygon.C
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 #include "libmesh/face_polygon.h"
19 
20 // Local includes
21 #include "libmesh/edge_edge2.h"
22 #include "libmesh/enum_elem_quality.h"
23 #include "libmesh/hashword.h"
24 
25 // C++ includes
26 #include <array>
27 
28 
29 namespace libMesh
30 {
31 
32 // ------------------------------------------------------------
33 // Polygon class static member initializations
34 
35 const int Polygon::num_children;
36 
37 // ------------------------------------------------------------
38 // Polygon class member functions
39 
40 
46 Polygon::Polygon (const unsigned int nn,
47  const unsigned int ns,
48  Elem * p) :
49  Face(nn, ns, p, nullptr, nullptr),
50  _elemlinks_data(ns+2), // neighbors + parent + interior_parent
51  _nodelinks_data(nn)
52 {
53  // Do the manual initialization that Elem::Elem couldn't. No need
54  // to manually set nullptr, though, since std::vector does that.
55  this->_elemlinks = _elemlinks_data.data();
56  this->_nodes = _nodelinks_data.data();
57  this->_elemlinks[0] = p;
58 
59  // Is this likely to ever be used? We may do refinement with
60  // polygons but it's probably not going to have a hierarchy...
61  if (p)
62  {
63  this->subdomain_id() = p->subdomain_id();
64  this->processor_id() = p->processor_id();
65  _map_type = p->mapping_type();
66  _map_data = p->mapping_data();
67 
68 #ifdef LIBMESH_ENABLE_AMR
69  this->set_p_level(p->p_level());
70 #endif
71  }
72 
73  // Make sure the interior parent isn't undefined
74  if (LIBMESH_DIM > 2)
75  this->set_interior_parent(nullptr);
76 }
77 
78 
79 
80 Point Polygon::master_point (const unsigned int i) const
81 {
82  const unsigned int ns = this->n_sides();
83 
84  // Non-vertex points need to be handled by a subclass overriding
85  // this.
86  if (i > ns)
87  libmesh_not_implemented();
88 
89  // Center-to-midside to center-to-vertex angle
90  const Real pi_over_ns = libMesh::pi / ns;
91 
92  // Center-to-midside distance
93  const Real min_r = 0.5/tan(pi_over_ns);
94 
95  // Center-to-vertex distance
96  const Real max_r = std::sqrt(min_r*min_r + 0.25);
97 
98  const Point center(0.5, min_r);
99 
100  libmesh_assert_less(i, this->n_nodes());
101 
102  return center + Point(max_r*sin((int(i)*2-1)*pi_over_ns),
103  -max_r*cos((int(i)*2-1)*pi_over_ns));
104 
105 }
106 
107 
108 
110  const Real eps) const
111 {
112  const unsigned int ns = this->n_sides();
113 
114  // Center-to-midside to center-to-vertex angle
115  const Real pi_over_ns = libMesh::pi / ns;
116 
117  // Center-to-midside distance
118  const Real min_r = 0.5/tan(pi_over_ns);
119 
120  // Center-to-vertex distance
121  const Real max_r = std::sqrt(min_r*min_r + 0.25);
122 
123  const Point center(0.5, min_r);
124 
125  // Check that the point is on the same side of all the faces by
126  // testing whether:
127  //
128  // n_i.(p - x_i) <= 0
129  //
130  // for each i, where:
131  // n_i is the outward normal of face i,
132  // x_i is a point on face i.
133 
134  for (auto i : make_range(ns))
135  {
136  const Point x_i =
137  center + Point(max_r*sin((int(i)*2-1)*pi_over_ns),
138  -max_r*cos((int(i)*2-1)*pi_over_ns));
139  const Point n_i =
140  Point(sin(int(i)*2*pi_over_ns),
141  -cos(int(i)*2*pi_over_ns));
142 
143  if (n_i * (p - x_i) > eps)
144  return false;
145  }
146 
147  return true;
148 }
149 
150 
151 
152 dof_id_type Polygon::key (const unsigned int s) const
153 {
154  const int ns = this->n_sides();
155  libmesh_assert_less (s, ns);
156 
157  return this->compute_key(this->node_id(s),
158  this->node_id((s+1)%ns));
159 }
160 
161 
162 
163 dof_id_type Polygon::low_order_key (const unsigned int s) const
164 {
165  const int ns = this->n_sides();
166  libmesh_assert_less (s, ns);
167 
168  return this->compute_key(this->node_id(s),
169  this->node_id((s+1)%ns));
170 }
171 
172 
173 
174 unsigned int Polygon::local_side_node(unsigned int side,
175  unsigned int side_node) const
176 {
177  const auto ns = this->n_sides();
178  libmesh_assert_less (side, ns);
179  libmesh_assert_less (side_node, this->n_nodes_per_side());
180 
181  if (!side_node)
182  return side;
183 
184  if (side_node == 1)
185  return (side + 1) % ns;
186 
187  return (side_node - 1) * ns + side;
188 }
189 
190 
191 
192 unsigned int Polygon::local_edge_node(unsigned int edge,
193  unsigned int edge_node) const
194 {
195  return local_side_node(edge, edge_node);
196 }
197 
198 
199 
201 {
202  std::vector<dof_id_type> node_ids;
203  for (const auto & n : this->node_ref_range())
204  node_ids.push_back(n.id());
205 
206  return Utility::hashword(node_ids);
207 }
208 
209 
210 
211 std::unique_ptr<Elem> Polygon::side_ptr (const unsigned int i)
212 {
213  const auto ns = this->n_sides();
214  libmesh_assert_less (i, ns);
215 
216  std::unique_ptr<Elem> edge = std::make_unique<Edge2>();
217 
218  edge->set_node(0, this->node_ptr(i));
219  edge->set_node(1, this->node_ptr((i+1)%ns));
220 
221  return edge;
222 }
223 
224 
225 
226 void Polygon::side_ptr (std::unique_ptr<Elem> & side,
227  const unsigned int i)
228 {
229  const auto ns = this->n_sides();
230  libmesh_assert_less (i, ns);
231 
232  if (!side.get() || side->type() != EDGE2)
233  side = std::make_unique<Edge2>();
234 
235  side->subdomain_id() = this->subdomain_id();
236 
237  side->set_node(0, this->node_ptr(i));
238  side->set_node(1, this->node_ptr((i+1)%ns));
239 }
240 
241 
242 
243 bool Polygon::is_child_on_side(const unsigned int /*c*/,
244  const unsigned int /*s*/) const
245 {
246  libmesh_not_implemented();
247  return false;
248 }
249 
250 
251 
252 unsigned int Polygon::opposite_side(const unsigned int side_in) const
253 {
254  const auto ns = this->n_sides();
255  if (ns % 2)
256  libmesh_error();
257 
258  return (ns / 2 + side_in) % ns;
259 }
260 
261 
262 
264 {
265  // Just check based on vertices, the first ns nodes
266  const auto ns = this->n_sides();
267 
268 #if LIBMESH_DIM > 2
269  // Don't bother outside the XY plane
270  for (auto n : make_range(ns))
271  if (this->point(n)(2))
272  return false;
273 #endif
274 
275  Real twice_signed_area = 0;
276  for (auto n : make_range(ns))
277  {
278  const Point & p1 = this->point(n);
279  const Point & p2 = this->point((n+1)%ns);
280  twice_signed_area += p1(0)*p2(1)-p1(1)*p2(0);
281  }
282 
283  return (twice_signed_area < 0);
284 }
285 
286 
287 std::vector<unsigned int>
288 Polygon::edges_adjacent_to_node(const unsigned int n) const
289 {
290  libmesh_assert_less(n, this->n_nodes());
291 
292  // For mid-face nodes, the subclass had better have overridden this.
293  libmesh_assert_less(n, this->n_sides() * (this->n_nodes_per_side() + 1));
294 
295  const unsigned int ns = this->n_sides();
296 
297  // For vertices, we have two adjacent edges; otherwise each of the
298  // mid-edge nodes is adjacent only to the edge it is on.
299  if (this->is_vertex(n))
300  return {n, (n+ns-1)%ns};
301 
302  return {n % ns};
303 }
304 
305 
306 std::pair<Real, Real> Polygon::qual_bounds (const ElemQuality q) const
307 {
308  std::pair<Real, Real> bounds;
309 
310  switch (q)
311  {
312  case EDGE_LENGTH_RATIO:
313  bounds.first = 1.;
314  bounds.second = 4.;
315  break;
316 
317  case MIN_ANGLE:
318  bounds.first = 90. - (180./this->n_sides());
319  bounds.second = 180. - (360./this->n_sides());
320  break;
321 
322  case MAX_ANGLE:
323  bounds.first = 180. - (360./this->n_sides());
324  bounds.second = 180. - (180./this->n_sides());
325  break;
326 
327  case JACOBIAN:
328  case SCALED_JACOBIAN:
329  bounds.first = 0.5;
330  bounds.second = 1.;
331  break;
332 
333  default:
334  libMesh::out << "Warning: Invalid quality measure chosen." << std::endl;
335  bounds.first = -1;
336  bounds.second = -1;
337  }
338 
339  return bounds;
340 }
341 
342 
343 std::array<Point, 3> Polygon::master_subtriangle (unsigned int i) const
344 {
345  libmesh_assert_less(i, this->_triangulation.size());
346 
347  const auto & tri = this->_triangulation[i];
348 
349  return { this->master_point(tri[0]),
350  this->master_point(tri[1]),
351  this->master_point(tri[2]) };
352 }
353 
354 
355 std::tuple<unsigned int, Real, Real>
357 {
358  std::tuple<unsigned int, Real, Real> returnval = {libMesh::invalid_uint, -1, -1};
359 
360  Real best_bad_coord = -1;
361 
362  for (auto s : make_range(this->n_subtriangles()))
363  {
364  const std::array<Point, 3> subtri =
365  this->master_subtriangle (s);
366 
367  // Find barycentric coordinates in subtri.
368  // Hand coding these since we don't need a full 3D cross product
369  // in 2D
370  const Real twice_area = (- subtri[1](1) * subtri[2](0)
371  - subtri[0](1) * subtri[1](0)
372  + subtri[0](1) * subtri[2](0)
373  + subtri[0](0) * subtri[1](1)
374  - subtri[0](0) * subtri[2](1)
375  + subtri[1](0) * subtri[2](1));
376 
377  const Real xi = ( subtri[0](1)*subtri[2](0)
378  - subtri[0](0)*subtri[2](1)
379  + (subtri[2](1)-subtri[0](1))*p(0)
380  + (subtri[0](0)-subtri[2](0))*p(1)) / twice_area;
381 
382  const Real eta = ( subtri[0](0)*subtri[1](1)
383  - subtri[0](1)*subtri[1](0)
384  + (subtri[0](1)-subtri[1](1))*p(0)
385  + (subtri[1](0)-subtri[0](0))*p(1)) / twice_area;
386 
387  if (xi>=0 && eta>=0 && xi+eta<=1)
388  return { s, xi, eta };
389 
390  Real my_best_bad_coord = std::min(std::min(xi, eta), 1-xi-eta);
391 
392  if (my_best_bad_coord > best_bad_coord)
393  {
394  best_bad_coord = my_best_bad_coord;
395  returnval = { s, xi, eta };
396  }
397  }
398 
399  if (best_bad_coord > -tol)
400  return returnval;
401 
402  return {libMesh::invalid_uint, -1, -1};
403 }
404 
405 
406 } // namespace libMesh
void set_p_level(const unsigned int p)
Sets the value of the p-refinement level for the element.
unsigned char mapping_data() const
Definition: elem.h:3136
virtual unsigned int n_nodes() const override
Definition: face_polygon.h:78
unsigned char _map_type
Mapping function type; currently either 0 (LAGRANGE) or 1 (RATIONAL_BERNSTEIN).
Definition: elem.h:2297
unsigned char _map_data
Mapping function data; currently used when needed to store the RATIONAL_BERNSTEIN nodal weight data i...
Definition: elem.h:2303
Node ** _nodes
Pointers to the nodes we are connected to.
Definition: elem.h:2245
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
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
Definition: libmesh.h:310
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
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
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
unsigned int p_level() const
Definition: elem.h:3108
virtual std::unique_ptr< Elem > side_ptr(const unsigned int i) override final
Definition: face_polygon.C:211
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
void set_interior_parent(Elem *p)
Sets the pointer to the element&#39;s interior_parent.
Definition: elem.C:1248
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const override
Definition: face_polygon.C:243
uint32_t hashword(const uint32_t *k, size_t length, uint32_t initval=0)
The hashword function takes an array of uint32_t&#39;s of length &#39;length&#39; and computes a single key from ...
Definition: hashword.h:158
ElemMappingType mapping_type() const
Definition: elem.h:3120
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const override
Definition: face_polygon.C:174
Elem ** _elemlinks
Pointers to this element&#39;s parent and neighbors, and for lower-dimensional elements&#39; interior_parent...
Definition: elem.h:2251
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 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
SimpleRange< NodeRefIter > node_ref_range()
Returns a range with all nodes of an element, usable in range-based for loops.
Definition: elem.h:2665
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
subdomain_id_type subdomain_id() const
Definition: elem.h:2582
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:2507
virtual bool is_vertex(const unsigned int i) const =0
OStreamProxy out
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
virtual unsigned int n_sides() const override final
Definition: face_polygon.h:84
static dof_id_type compute_key(dof_id_type n0)
Definition: elem.h:3294
processor_id_type processor_id() const
Definition: dof_object.h:905
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
dof_id_type node_id(const unsigned int i) const
Definition: elem.h:2475
const Point & point(const unsigned int i) const
Definition: elem.h:2453
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
const Real pi
.
Definition: libmesh.h:299