18 #include "libmesh/face_polygon.h" 21 #include "libmesh/edge_edge2.h" 22 #include "libmesh/enum_elem_quality.h" 23 #include "libmesh/hashword.h" 47 const unsigned int ns,
49 Face(nn, ns, p, nullptr, nullptr),
50 _elemlinks_data(ns+2),
68 #ifdef LIBMESH_ENABLE_AMR 82 const unsigned int ns = this->
n_sides();
87 libmesh_not_implemented();
93 const Real min_r = 0.5/tan(pi_over_ns);
96 const Real max_r = std::sqrt(min_r*min_r + 0.25);
98 const Point center(0.5, min_r);
100 libmesh_assert_less(i, this->
n_nodes());
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));
110 const Real eps)
const 112 const unsigned int ns = this->
n_sides();
118 const Real min_r = 0.5/tan(pi_over_ns);
121 const Real max_r = std::sqrt(min_r*min_r + 0.25);
123 const Point center(0.5, min_r);
137 center +
Point(max_r*sin((
int(i)*2-1)*pi_over_ns),
138 -max_r*cos((
int(i)*2-1)*pi_over_ns));
140 Point(sin(
int(i)*2*pi_over_ns),
141 -cos(
int(i)*2*pi_over_ns));
143 if (n_i * (p - x_i) > eps)
154 const int ns = this->
n_sides();
155 libmesh_assert_less (s, ns);
165 const int ns = this->
n_sides();
166 libmesh_assert_less (s, ns);
175 unsigned int side_node)
const 177 const auto ns = this->
n_sides();
178 libmesh_assert_less (side, ns);
185 return (side + 1) % ns;
187 return (side_node - 1) * ns + side;
193 unsigned int edge_node)
const 202 std::vector<dof_id_type> node_ids;
204 node_ids.push_back(n.id());
213 const auto ns = this->
n_sides();
214 libmesh_assert_less (i, ns);
216 std::unique_ptr<Elem> edge = std::make_unique<Edge2>();
218 edge->set_node(0, this->
node_ptr(i));
219 edge->set_node(1, this->
node_ptr((i+1)%ns));
227 const unsigned int i)
229 const auto ns = this->
n_sides();
230 libmesh_assert_less (i, ns);
232 if (!side.get() || side->type() !=
EDGE2)
233 side = std::make_unique<Edge2>();
237 side->set_node(0, this->
node_ptr(i));
238 side->set_node(1, this->
node_ptr((i+1)%ns));
244 const unsigned int )
const 246 libmesh_not_implemented();
254 const auto ns = this->
n_sides();
258 return (ns / 2 + side_in) % ns;
266 const auto ns = this->
n_sides();
271 if (this->
point(n)(2))
275 Real twice_signed_area = 0;
280 twice_signed_area += p1(0)*p2(1)-p1(1)*p2(0);
283 return (twice_signed_area < 0);
287 std::vector<unsigned int>
290 libmesh_assert_less(n, this->
n_nodes());
295 const unsigned int ns = this->
n_sides();
300 return {n, (n+ns-1)%ns};
308 std::pair<Real, Real> bounds;
318 bounds.first = 90. - (180./this->
n_sides());
319 bounds.second = 180. - (360./this->
n_sides());
323 bounds.first = 180. - (360./this->
n_sides());
324 bounds.second = 180. - (180./this->
n_sides());
334 libMesh::out <<
"Warning: Invalid quality measure chosen." << std::endl;
355 std::tuple<unsigned int, Real, Real>
360 Real best_bad_coord = -1;
364 const std::array<Point, 3> subtri =
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));
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;
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;
387 if (xi>=0 && eta>=0 && xi+eta<=1)
388 return { s, xi, eta };
390 Real my_best_bad_coord = std::min(std::min(xi, eta), 1-xi-eta);
392 if (my_best_bad_coord > best_bad_coord)
394 best_bad_coord = my_best_bad_coord;
395 returnval = { s, xi, eta };
399 if (best_bad_coord > -tol)
void set_p_level(const unsigned int p)
Sets the value of the p-refinement level for the element.
unsigned char mapping_data() const
virtual unsigned int n_nodes() const override
unsigned char _map_type
Mapping function type; currently either 0 (LAGRANGE) or 1 (RATIONAL_BERNSTEIN).
unsigned char _map_data
Mapping function data; currently used when needed to store the RATIONAL_BERNSTEIN nodal weight data i...
Node ** _nodes
Pointers to the nodes we are connected to.
static const int num_children
std::vector< Node * > _nodelinks_data
Data for links to nodes.
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
The Face is an abstract element type that lives in two dimensions.
unsigned int n_subtriangles() const
virtual dof_id_type key() const override
std::tuple< unsigned int, Real, Real > subtriangle_coordinates(const Point &p, Real tol=TOLERANCE *TOLERANCE) const
virtual bool is_flipped() const override final
std::vector< std::array< int, 3 > > _triangulation
Data for a triangulation of the polygon.
virtual unsigned int n_nodes_per_side() const =0
virtual Point master_point(const unsigned int i) const override
virtual bool on_reference_element(const Point &p, const Real eps=TOLERANCE) const override final
This is the base class from which all geometric element types are derived.
unsigned int p_level() const
virtual std::unique_ptr< Elem > side_ptr(const unsigned int i) override final
virtual unsigned int opposite_side(const unsigned int s) const override final
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).
void set_interior_parent(Elem *p)
Sets the pointer to the element's interior_parent.
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const override
uint32_t hashword(const uint32_t *k, size_t length, uint32_t initval=0)
The hashword function takes an array of uint32_t's of length 'length' and computes a single key from ...
ElemMappingType mapping_type() const
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const override
Elem ** _elemlinks
Pointers to this element's parent and neighbors, and for lower-dimensional elements' interior_parent...
virtual std::array< Point, 3 > master_subtriangle(unsigned int i) const
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...
virtual std::pair< Real, Real > qual_bounds(const ElemQuality q) const override
virtual std::vector< unsigned int > edges_adjacent_to_node(const unsigned int n) const override
ElemQuality
Defines an enum for element quality metrics.
std::vector< Elem * > _elemlinks_data
Data for links to parent/neighbor/interior_parent elements.
SimpleRange< NodeRefIter > node_ref_range()
Returns a range with all nodes of an element, usable in range-based for loops.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
subdomain_id_type subdomain_id() const
const Node * node_ptr(const unsigned int i) const
virtual bool is_vertex(const unsigned int i) const =0
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...
virtual unsigned int n_sides() const override final
static dof_id_type compute_key(dof_id_type n0)
processor_id_type processor_id() const
A Point defines a location in LIBMESH_DIM dimensional Real space.
dof_id_type node_id(const unsigned int i) const
const Point & point(const unsigned int i) const
virtual dof_id_type low_order_key(const unsigned int s) const override