libMesh
face_c0polygon.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 // Local includes
19 #include "libmesh/face_c0polygon.h"
20 
21 #include "libmesh/edge_edge2.h"
22 #include "libmesh/enum_order.h"
23 #include "libmesh/tensor_value.h"
24 
25 // C++ headers
26 #include <numeric> // std::iota
27 
28 namespace libMesh
29 {
30 
31 
32 C0Polygon::C0Polygon (const unsigned int num_sides, Elem * p) :
33  Polygon(num_sides, num_sides, p)
34 {
35  // A default triangulation is better than nothing
36  for (int i : make_range(num_sides-2))
37  this->_triangulation.push_back({0, i+1, i+2});
38 }
39 
40 
41 unsigned int C0Polygon::opposite_node(const unsigned int node_in,
42  const unsigned int side_in) const
43 {
44  const auto ns = this->n_sides();
45  if (ns % 2)
46  libmesh_error();
47 
48  libmesh_assert_less (node_in, ns);
49  libmesh_assert_less (node_in, this->n_nodes());
50  libmesh_assert_less (side_in, this->n_sides());
51  libmesh_assert(this->is_node_on_side(node_in, side_in));
52 
53  if (node_in == side_in)
54  return (node_in + ns/2 + 1) % ns;
55 
56  libmesh_assert(node_in == side_in + 1);
57  return (node_in + ns/2 - 1) % ns;
58 }
59 
60 
61 
62 // ------------------------------------------------------------
63 // C0Polygon class member functions
64 
65 bool C0Polygon::is_vertex(const unsigned int libmesh_dbg_var(i)) const
66 {
67  libmesh_assert (i < this->n_nodes());
68 
69  return true;
70 }
71 
72 bool C0Polygon::is_edge(const unsigned int libmesh_dbg_var(i)) const
73 {
74  libmesh_assert (i < this->n_nodes());
75 
76  return false;
77 }
78 
79 bool C0Polygon::is_face(const unsigned int libmesh_dbg_var(i)) const
80 {
81  libmesh_assert (i < this->n_nodes());
82 
83  return false;
84 }
85 
86 bool C0Polygon::is_node_on_side(const unsigned int n,
87  const unsigned int s) const
88 {
89  const auto ns = this->n_sides();
90  libmesh_assert_less (s, ns);
91  libmesh_assert_less (n, this->n_nodes());
92 
93  return ((n % ns) == s) ||
94  ((n < ns) && ((s+1)%ns) == n);
95 }
96 
97 std::vector<unsigned int>
98 C0Polygon::nodes_on_side(const unsigned int s) const
99 {
100  const auto ns = this->n_sides();
101  libmesh_assert(!(this->n_nodes() % ns));
102  libmesh_assert_less(s, ns);
103 
104  std::vector<unsigned int> returnval(2);
105  returnval[0] = s;
106  returnval[1] = (s+1)%ns;
107 
108  return returnval;
109 }
110 
111 std::vector<unsigned int>
112 C0Polygon::nodes_on_edge(const unsigned int e) const
113 {
114  return this->nodes_on_side(e);
115 }
116 
118 {
119  const unsigned int ns = this->n_sides();
120 
121  // Get a good tolerance to use
122  Real perimeter_l1 = 0;
123  for (auto s : make_range(ns))
124  perimeter_l1 += (this->point((s+1)%ns) - this->point(s)).l1_norm();
125  const Real tol = perimeter_l1 * affine_tol;
126 
127  // Find the map we expect to have
128  const Point veci = this->point(1) - this->point(0);
129  const Point vec12 = this->point(2) - this->point(1);
130 
131  // Exterior angle
132  const Real theta = 2 * libMesh::pi / ns;
133  const Real costheta = cos(theta);
134  const Real sintheta = sin(theta);
135 
136  RealTensor map;
137  map(0, 0) = veci(0);
138  map(1, 0) = veci(1);
139  map(2, 0) = veci(2);
140  map(0, 1) = (vec12(0) - veci(0)*costheta)/sintheta;
141  map(1, 1) = (vec12(1) - veci(1)*costheta)/sintheta;
142  map(2, 1) = (vec12(2) - veci(2)*costheta)/sintheta;
143 
144  libmesh_assert_less((map * this->master_point(1) -
145  (this->point(1) - this->point(0))).l1_norm(),
146  tol);
147  libmesh_assert_less((map * this->master_point(2) -
148  (this->point(2) - this->point(0))).l1_norm(),
149  tol);
150 
151  for (auto i : make_range(3u, ns))
152  if ((map * this->master_point(i) -
153  (this->point(i) - this->point(0))).l1_norm() >
154  tol)
155  return false;
156 
157  return true;
158 }
159 
160 
161 
163 {
164  return FIRST;
165 }
166 
167 
168 
169 std::unique_ptr<Elem> C0Polygon::build_side_ptr (const unsigned int i)
170 {
171  const auto ns = this->n_sides();
172  libmesh_assert_less (i, ns);
173 
174  std::unique_ptr<Elem> sidep = std::make_unique<Edge2>();
175  sidep->set_node(0, this->node_ptr(i));
176  sidep->set_node(1, this->node_ptr((i+1)%ns));
177 
178  sidep->set_interior_parent(this);
179  sidep->inherit_data_from(*this);
180 
181  return sidep;
182 }
183 
184 
185 
186 void C0Polygon::build_side_ptr (std::unique_ptr<Elem> & side,
187  const unsigned int i)
188 {
189  const auto ns = this->n_sides();
190  libmesh_assert_less (i, ns);
191 
192  if (!side.get() || side->type() != EDGE2)
193  {
194  side = this->build_side_ptr(i);
195  }
196  else
197  {
198  side->inherit_data_from(*this);
199 
200  side->set_node(0, this->node_ptr(i));
201  side->set_node(1, this->node_ptr((i+1)%ns));
202  }
203 }
204 
205 
206 
207 void C0Polygon::connectivity(const unsigned int /*sf*/,
208  const IOPackage /*iop*/,
209  std::vector<dof_id_type> & /*conn*/) const
210 {
211  libmesh_not_implemented();
212 }
213 
214 
215 
217 {
218  // This specialization is good for Lagrange mappings only
219  if (this->mapping_type() != LAGRANGE_MAP)
220  return this->Elem::volume();
221 
222  // We use a triangulation to calculate here; assert that it's as
223  // consistent as possible.
224  libmesh_assert_equal_to (this->_triangulation.size(),
225  this->n_nodes() - 2);
226 
227  Real double_area = 0;
228  for (const auto & triangle : this->_triangulation)
229  {
230  Point v01 = this->point(triangle[1]) -
231  this->point(triangle[0]);
232  Point v02 = this->point(triangle[2]) -
233  this->point(triangle[0]);
234  double_area += std::sqrt(cross_norm_sq(v01, v02));
235  }
236 
237  return double_area/2;
238 }
239 
240 
241 
243 {
244  // This specialization is good for Lagrange mappings only
245  if (this->mapping_type() != LAGRANGE_MAP)
246  return this->Elem::true_centroid();
247 
248  // We use a triangulation to calculate here; assert that it's as
249  // consistent as possible.
250  libmesh_assert_equal_to (this->_triangulation.size(),
251  this->n_nodes() - 2);
252 
253  Real double_area = 0;
254  Point double_area_weighted_centroid;
255  for (const auto & triangle : this->_triangulation)
256  {
257  Point v01 = this->point(triangle[1]) -
258  this->point(triangle[0]);
259  Point v02 = this->point(triangle[2]) -
260  this->point(triangle[0]);
261 
262  const Real double_tri_area = std::sqrt(cross_norm_sq(v01, v02));
263  const Point tri_centroid = (this->point(triangle[0]) +
264  this->point(triangle[1]) +
265  this->point(triangle[2]))/3;
266 
267  double_area += double_tri_area;
268 
269  double_area_weighted_centroid += double_tri_area * tri_centroid;
270  }
271 
272  return double_area_weighted_centroid / double_area;
273 }
274 
275 
276 
277 std::pair<unsigned short int, unsigned short int>
278 C0Polygon::second_order_child_vertex (const unsigned int /*n*/) const
279 {
280  libmesh_not_implemented();
281  return std::pair<unsigned short int, unsigned short int> (0, 0);
282 }
283 
284 
285 void C0Polygon::permute(unsigned int perm_num)
286 {
287  const auto ns = this->n_sides();
288  libmesh_assert_less (perm_num, ns);
289 
290  // This is mostly for unit testing so I'll just make it O(N).
291  for (unsigned int p : make_range(perm_num))
292  {
293  libmesh_ignore(p);
294 
295  Node * tempnode = this->node_ptr(0);
296  Elem * tempneigh = this->neighbor_ptr(0);
297  for (unsigned int s : make_range(ns-1))
298  {
299  this->set_node(s, this->node_ptr((s+1)%ns));
300  this->set_neighbor(s, this->neighbor_ptr(s+1));
301  }
302  this->set_node(ns-1, tempnode);
303  this->set_neighbor(ns-1, tempneigh);
304 
305  // Keep the same triangulation, just with permuted node order,
306  // so we can really expect this to act like the *same* polygon
307  for (auto & triangle : this->_triangulation)
308  for (auto i : make_range(3))
309  triangle[i] = (triangle[i]+1)%ns;
310  }
311 }
312 
313 
314 void C0Polygon::flip(BoundaryInfo * boundary_info)
315 {
316  libmesh_assert(boundary_info);
317 
318  const auto ns = this->n_sides();
319 
320  // Reorder nodes in such a way as to keep side ns-1 the same
321  for (auto s : make_range(0u, ns/2))
322  {
323  swap2nodes(s, ns-1-s);
324  swap2neighbors(s, ns-2-s);
325  swap2boundarysides(s, ns-2-s, boundary_info);
326  swap2boundaryedges(s, ns-2-s, boundary_info);
327  }
328 }
329 
330 
331 
332 ElemType C0Polygon::side_type (const unsigned int libmesh_dbg_var(s)) const
333 {
334  libmesh_assert_less (s, this->n_sides());
335  return EDGE2;
336 }
337 
338 
339 
341 {
342  this->_triangulation.clear();
343 
344  // Start with the full polygon
345  std::vector<int> remaining_nodes(this->n_nodes());
346  std::iota(remaining_nodes.begin(), remaining_nodes.end(), 0);
347 
348  const auto ns = this->n_sides();
349 
350  const Point vavg = this->vertex_average();
351 
352  // Find out what plane we're in, if we're in 3D
353 #if LIBMESH_DIM > 2
354  Point plane_normal;
355  for (auto i : make_range(ns))
356  {
357  const Point vi = this->point(i) - vavg;
358  const Point viplus = this->point((i+1)%ns) - vavg;
359  plane_normal += vi.cross(viplus);
360  }
361  plane_normal = plane_normal.unit();
362 #endif
363 
364  // Greedy heuristic: find the node with the most acutely convex
365  // angle and strip it out as a triangle.
366  while (remaining_nodes.size() > 2)
367  {
368  Real min_cos_angle = 1;
369  int best_vertex = -1;
370  for (auto n : index_range(remaining_nodes))
371  {
372  const Point & pn = this->point(remaining_nodes[n]);
373  const Point & pnext = this->point(remaining_nodes[(n+1)%ns]);
374  const Point & pprev = this->point(remaining_nodes[(n+ns-1)%ns]);
375  const Point vprev = (pn - pprev).unit();
376  const Point vnext = (pnext - pn).unit();
377 
378  const Real sign_check = (vprev.cross(vnext)) * plane_normal;
379  if (sign_check <= 0)
380  continue;
381 
382  const Real cos_angle = vprev * vnext;
383  if (cos_angle < min_cos_angle)
384  {
385  min_cos_angle = cos_angle;
386  best_vertex = n;
387  }
388  }
389 
390  libmesh_assert(best_vertex >= 0);
391 
392  this->_triangulation.push_back({remaining_nodes[(best_vertex+ns-1)%ns],
393  remaining_nodes[best_vertex],
394  remaining_nodes[(best_vertex+1)%ns]});
395  remaining_nodes.erase(remaining_nodes.begin()+best_vertex);
396  }
397 }
398 
399 
400 
401 } // namespace libMesh
virtual Point true_centroid() const
Definition: elem.C:594
virtual Real volume() const override
An optimized method for calculating the area of a C0Polygon.
virtual unsigned int n_nodes() const override
Definition: face_polygon.h:78
ElemType
Defines an enum for geometric element types.
void swap2boundaryedges(unsigned short e1, unsigned short e2, BoundaryInfo *boundary_info) const
Swaps two edges in boundary_info, if it is non-null.
Definition: elem.C:3550
virtual std::pair< unsigned short int, unsigned short int > second_order_child_vertex(const unsigned int n) const override
Element refinement is not implemented for polygons.
virtual void permute(unsigned int perm_num) override final
Permutes the element (by swapping node and neighbor pointers) according to the specified index...
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
virtual Node *& set_node(const unsigned int i)
Definition: elem.h:2558
A Node is like a Point, but with more information.
Definition: node.h:52
virtual Order default_order() const override
virtual bool is_edge(const unsigned int i) const override
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
std::vector< std::array< int, 3 > > _triangulation
Data for a triangulation of the polygon.
Definition: face_polygon.h:269
virtual Point master_point(const unsigned int i) const override
Definition: face_polygon.C:80
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int e) const override
void swap2boundarysides(unsigned short s1, unsigned short s2, BoundaryInfo *boundary_info) const
Swaps two sides in boundary_info, if it is non-null.
Definition: elem.C:3534
The libMesh namespace provides an interface to certain functionality in the library.
ElemType side_type(const unsigned int s) const override final
void iota(ForwardIter first, ForwardIter last, T value)
Utility::iota was created back when std::iota was just an SGI STL extension.
Definition: utility.h:229
T cross_norm_sq(const TypeVector< T > &b, const TypeVector< T > &c)
Compute |b x c|^2 without creating the extra temporary produced by calling b.cross(c).norm_sq().
Definition: type_vector.h:1075
virtual bool is_vertex(const unsigned int i) const override
TypeVector< T > unit() const
Definition: type_vector.h:1104
void libmesh_ignore(const Args &...)
virtual unsigned int opposite_node(const unsigned int n, const unsigned int s) const override final
virtual bool has_affine_map() const override
ElemMappingType mapping_type() const
Definition: elem.h:3120
void swap2nodes(unsigned int n1, unsigned int n2)
Swaps two node_ptrs.
Definition: elem.h:2092
auto l1_norm(const NumericVector< T > &vec)
The Polygon is an element in 2D with an arbitrary (but fixed) number of sides.
Definition: face_polygon.h:39
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
libmesh_assert(ctx)
virtual bool is_face(const unsigned int i) const override
virtual void connectivity(const unsigned int sf, const IOPackage iop, std::vector< dof_id_type > &conn) const override
TypeVector< typename CompareTypes< T, T2 >::supertype > cross(const TypeVector< T2 > &v) const
Definition: type_vector.h:884
static constexpr Real affine_tol
Default tolerance to use in has_affine_map().
Definition: elem.h:2061
void set_neighbor(const unsigned int i, Elem *n)
Assigns n as the neighbor.
Definition: elem.h:2618
void swap2neighbors(unsigned int n1, unsigned int n2)
Swaps two neighbor_ptrs.
Definition: elem.h:2102
const Elem * neighbor_ptr(unsigned int i) const
Definition: elem.h:2598
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i) override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:2507
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
virtual Real volume() const
Definition: elem.C:3429
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
C0Polygon(const unsigned int num_sides, Elem *p=nullptr)
Constructor.
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
const Point & point(const unsigned int i) const
Definition: elem.h:2453
virtual void retriangulate() override final
Create a triangulation from the current node locations.
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117
virtual void flip(BoundaryInfo *) override final
Flips the element (by swapping node and neighbor pointers) to have a mapping Jacobian of opposite sig...
Point vertex_average() const
Definition: elem.C:688
virtual Point true_centroid() const override
An optimized method for calculating the centroid of a C0Polygon.
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
const Real pi
.
Definition: libmesh.h:299