libMesh
Classes | Public Types | Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
libMesh::TriangulatorInterface Class Referenceabstract

#include <triangulator_interface.h>

Inheritance diagram for libMesh::TriangulatorInterface:
[legend]

Classes

class  AffineHole
 A way to translate and/or rotate an existing hole; perhaps to tile it in many places or to put it at an angle that the underlying hole doesn't support. More...
 
class  ArbitraryHole
 Another concrete instantiation of the hole, this one should be sufficiently general for most non-polygonal purposes. More...
 
class  Hole
 An abstract class for defining a 2-dimensional hole. More...
 
class  MeshedHole
 Another concrete instantiation of the hole, as general as ArbitraryHole, but based on an existing 1D or 2D mesh. More...
 
class  PolygonHole
 A concrete instantiation of the Hole class that describes polygonal (triangular, square, pentagonal, ...) holes. More...
 
class  Region
 A class for defining a 2-dimensional region for Triangle. More...
 

Public Types

enum  TriangulationType { GENERATE_CONVEX_HULL = 0, PSLG = 1, INVALID_TRIANGULATION_TYPE }
 The TriangulationType is used with the general triangulate function defined below. More...
 

Public Member Functions

 TriangulatorInterface (UnstructuredMesh &mesh)
 The constructor. More...
 
virtual ~TriangulatorInterface ()=default
 Empty destructor. More...
 
virtual void triangulate ()=0
 This is the main public interface for this function. More...
 
ElemTypeelem_type ()
 Sets and/or gets the desired element type. More...
 
Realdesired_area ()
 Sets and/or gets the desired triangle area. More...
 
virtual void set_desired_area_function (FunctionBase< Real > *)
 Set a function giving desired triangle area as a function of position. More...
 
virtual FunctionBase< Real > * get_desired_area_function ()
 Get the function giving desired triangle area as a function of position, or nullptr if no such function has been set. More...
 
Realminimum_angle ()
 Sets and/or gets the minimum desired angle. More...
 
TriangulationTypetriangulation_type ()
 Sets and/or gets the desired triangulation type. More...
 
bool & insert_extra_points ()
 Sets and/or gets the flag for inserting add'l points. More...
 
void set_interpolate_boundary_points (int n_points)
 Complicated setter, for compatibility with insert_extra_points() More...
 
int get_interpolate_boundary_points () const
 Complicated getter, for compatibility with insert_extra_points() More...
 
virtual void set_refine_boundary_allowed (bool)
 Set whether or not the triangulation is allowed to refine the mesh boundary when refining the interior. More...
 
virtual bool refine_boundary_allowed () const
 Get whether or not the triangulation is allowed to refine the mesh boundary when refining the interior. More...
 
bool & smooth_after_generating ()
 Sets/gets flag which tells whether to do two steps of Laplace mesh smoothing after generating the grid. More...
 
bool & quiet ()
 Whether not to silence internal messages to stdout. More...
 
void attach_hole_list (const std::vector< Hole *> *holes)
 Attaches a vector of Hole* pointers which will be meshed around. More...
 
void set_verify_hole_boundaries (bool v)
 Verifying that hole boundaries don't cross the outer boundary or each other is something like O(N_bdys^2*N_points_per_bdy^2), so we only do it if requested. More...
 
bool get_verify_hole_boundaries () const
 
void attach_boundary_marker (const std::vector< int > *markers)
 Attaches boundary markers. More...
 
void attach_region_list (const std::vector< Region *> *regions)
 Attaches regions for using attribute to set subdomain IDs and better controlling the triangle sizes within the regions. More...
 
void set_auto_area_function (const Parallel::Communicator &comm, const unsigned int num_nearest_pts, const unsigned int power, const Real background_value, const Real background_eff_dist)
 Generate an auto area function based on spacing of boundary points. More...
 
bool has_auto_area_function ()
 Whether or not an auto area function has been set. More...
 
FunctionBase< Real > * get_auto_area_function ()
 Get the auto area function. More...
 
void calculate_auto_desired_area_samples (std::vector< Point > &function_points, std::vector< Real > &function_sizes, const Real &area_factor=1.5)
 The external boundary and all hole boundaries are collected. More...
 
void set_outer_boundary_ids (std::set< std::size_t > bdy_ids)
 A set of ids to allow on the outer boundary loop: interpreted as boundary ids of 2D elements and/or subdomain ids of 1D edges. More...
 
const std::set< std::size_t > & get_outer_boundary_ids () const
 

Public Attributes

std::vector< std::pair< unsigned int, unsigned int > > segments
 When constructing a PSLG, if the node numbers do not define the desired boundary segments implicitly through the ordering of the points, you can use the segments vector to specify the segments explicitly, Ex: unit square numbered counter-clockwise starting from origin segments[0] = (0,1) segments[1] = (1,2) segments[2] = (2,3) segments[3] = (3,0) (For the above case you could actually use the implicit ordering!) More...
 
std::vector< Pointsegment_midpoints
 When constructing a second-order triangulation from a second-order boundary, we may do the triangulation using first-order elements, in which case we need to save midpoint location data in order to reconstruct curvature along boundaries. More...
 

Protected Member Functions

void elems_to_segments ()
 Helper function to create PSLG segments from our other boundary-defining options (1D mesh edges, 2D mesh boundary sides), if no segments already exist. More...
 
void nodes_to_segments (dof_id_type max_node_id)
 Helper function to create PSLG segments from our node ordering, up to the maximum node id, if no segments already exist. More...
 
void insert_any_extra_boundary_points ()
 Helper function to add extra points (midpoints of initial segments) to a PSLG triangulation. More...
 
void increase_triangle_order ()
 Helper function to upconvert Tri3 to any higher order triangle type if requested via _elem_type. More...
 
void verify_holes (const Hole &outer_bdy)
 Helper function to check holes for intersections if requested. More...
 
unsigned int total_hole_points ()
 Helper function to count points in and verify holes. More...
 

Protected Attributes

UnstructuredMesh_mesh
 Reference to the mesh which is to be created by triangle. More...
 
const std::vector< Hole * > * _holes
 A pointer to a vector of Hole*s. More...
 
const std::vector< int > * _markers
 Boundary markers. More...
 
const std::vector< Region * > * _regions
 A pointer to a vector of Regions*s. More...
 
std::set< std::size_t > _bdy_ids
 A set of ids to allow on the outer boundary loop. More...
 
ElemType _elem_type
 The type of elements to generate. More...
 
Real _desired_area
 The desired area for the elements in the resulting mesh. More...
 
Real _minimum_angle
 Minimum angle in triangles. More...
 
TriangulationType _triangulation_type
 The type of triangulation to perform: choices are: convex hull PSLG. More...
 
bool _insert_extra_points
 Flag which tells whether or not to insert additional nodes before triangulation. More...
 
int _interpolate_boundary_points
 Flag which tells how many additional nodes should be inserted between each pair of original mesh points. More...
 
bool _smooth_after_generating
 Flag which tells whether we should smooth the mesh after it is generated. More...
 
bool _quiet
 Flag which tells if we want to suppress stdout outputs. More...
 
bool _verify_hole_boundaries
 Flag which tells if we want to check hole geometry. More...
 
std::unique_ptr< AutoAreaFunction_auto_area_function
 The auto area function based on the spacing of boundary points. More...
 

Detailed Description

Definition at line 90 of file triangulator_interface.h.

Member Enumeration Documentation

◆ TriangulationType

The TriangulationType is used with the general triangulate function defined below.

Enumerator
GENERATE_CONVEX_HULL 

First generate a convex hull from the set of points passed in, and then triangulate this set of points.

This is probably the most common type of usage.

PSLG 

Triangulate the interior of a Planar Straight Line Graph, which is defined implicitly by the order of the "points" vector: a straight line is assumed to lie between each successive pair of points, with an additional line joining the final and first points.

Explicitly telling the triangulator to add additional points may be important for this option.

INVALID_TRIANGULATION_TYPE 

Does nothing, used as a "null" value.

Definition at line 111 of file triangulator_interface.h.

112  {
119 
130  PSLG = 1,
131 
136  };
First generate a convex hull from the set of points passed in, and then triangulate this set of point...
Triangulate the interior of a Planar Straight Line Graph, which is defined implicitly by the order of...

Constructor & Destructor Documentation

◆ TriangulatorInterface()

libMesh::TriangulatorInterface::TriangulatorInterface ( UnstructuredMesh mesh)
explicit

The constructor.

A reference to the mesh containing the points which are to be triangulated must be provided. Unless otherwise specified, a convex hull will be computed for the set of input points and the convex hull will be meshed.

Definition at line 107 of file triangulator_interface.C.

108  : _mesh(mesh),
109  _holes(nullptr),
110  _markers(nullptr),
111  _regions(nullptr),
112  _elem_type(TRI3),
113  _desired_area(0.1),
114  _minimum_angle(20.0),
116  _insert_extra_points(false),
118  _quiet(true),
119  _auto_area_function(nullptr)
120 {}
bool _insert_extra_points
Flag which tells whether or not to insert additional nodes before triangulation.
const std::vector< int > * _markers
Boundary markers.
Real _minimum_angle
Minimum angle in triangles.
const std::vector< Region * > * _regions
A pointer to a vector of Regions*s.
MeshBase & mesh
UnstructuredMesh & _mesh
Reference to the mesh which is to be created by triangle.
Real _desired_area
The desired area for the elements in the resulting mesh.
ElemType _elem_type
The type of elements to generate.
bool _quiet
Flag which tells if we want to suppress stdout outputs.
const std::vector< Hole * > * _holes
A pointer to a vector of Hole*s.
First generate a convex hull from the set of points passed in, and then triangulate this set of point...
bool _smooth_after_generating
Flag which tells whether we should smooth the mesh after it is generated.
std::unique_ptr< AutoAreaFunction > _auto_area_function
The auto area function based on the spacing of boundary points.
TriangulationType _triangulation_type
The type of triangulation to perform: choices are: convex hull PSLG.

◆ ~TriangulatorInterface()

virtual libMesh::TriangulatorInterface::~TriangulatorInterface ( )
virtualdefault

Empty destructor.

Member Function Documentation

◆ attach_boundary_marker()

void libMesh::TriangulatorInterface::attach_boundary_marker ( const std::vector< int > *  markers)
inline

Attaches boundary markers.

If segments is set, the number of markers must be equal to the size of segments, otherwise, it is equal to the number of points.

Definition at line 294 of file triangulator_interface.h.

References _markers.

294 { _markers = markers; }
const std::vector< int > * _markers
Boundary markers.

◆ attach_hole_list()

void libMesh::TriangulatorInterface::attach_hole_list ( const std::vector< Hole *> *  holes)
inline

◆ attach_region_list()

void libMesh::TriangulatorInterface::attach_region_list ( const std::vector< Region *> *  regions)
inline

Attaches regions for using attribute to set subdomain IDs and better controlling the triangle sizes within the regions.

Definition at line 300 of file triangulator_interface.h.

References _regions.

300 { _regions = regions; }
const std::vector< Region * > * _regions
A pointer to a vector of Regions*s.

◆ calculate_auto_desired_area_samples()

void libMesh::TriangulatorInterface::calculate_auto_desired_area_samples ( std::vector< Point > &  function_points,
std::vector< Real > &  function_sizes,
const Real area_factor = 1.5 
)

The external boundary and all hole boundaries are collected.

The centroid of each EDGE element is used as the point position and the desired area is calculated as the area of the equilateral triangle with the edge length as the length of the EDGE element times an area_factor (default is 1.5).

Definition at line 500 of file triangulator_interface.C.

References _bdy_ids, _holes, _mesh, and libMesh::Real.

Referenced by get_auto_area_function().

503 {
504  // Get the hole mesh of the outer boundary
505  // Holes should already be attached if applicable when this function is called
506  const TriangulatorInterface::MeshedHole bdry_mh { _mesh, this->_bdy_ids };
507  // Collect all the centroid points of the outer boundary segments
508  // and the corresponding element sizes
509  for (unsigned int i = 0; i < bdry_mh.n_points(); i++)
510  {
511  function_points.push_back((bdry_mh.point(i) + bdry_mh.point((i + 1) % bdry_mh.n_points())) /
512  Real(2.0));
513  function_sizes.push_back(
514  (bdry_mh.point(i) - bdry_mh.point((i + 1) % bdry_mh.n_points())).norm());
515  }
516  // If holes are present, do the same for the hole boundaries
517  if(_holes)
518  for (const Hole * hole : *_holes)
519  {
520  for (unsigned int i = 0; i < hole->n_points(); i++)
521  {
522  function_points.push_back(
523  (hole->point(i) + hole->point((i + 1) % hole->n_points())) / Real(2.0));
524  function_sizes.push_back(
525  (hole->point(i) - hole->point((i + 1) % hole->n_points())).norm());
526  }
527  }
528 
529  std::for_each(
530  function_sizes.begin(), function_sizes.end(), [&area_factor](Real & a) { a = a * a * area_factor * std::sqrt(3.0) / 4.0; });
531 
532 }
std::set< std::size_t > _bdy_ids
A set of ids to allow on the outer boundary loop.
UnstructuredMesh & _mesh
Reference to the mesh which is to be created by triangle.
const std::vector< Hole * > * _holes
A pointer to a vector of Hole*s.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ desired_area()

Real& libMesh::TriangulatorInterface::desired_area ( )
inline

Sets and/or gets the desired triangle area.

Set to zero to disable area constraint.

If a desired_area_function is set, then desired_area() should be used to set a minimum desired area; this will reduce "false negatives" by suggesting how finely to sample desired_area_function inside large triangles, where ideally the desired_area_function will be satisfied in the triangle interior and not just at the triangle vertices.

Definition at line 177 of file triangulator_interface.h.

References _desired_area.

Referenced by libMesh::MeshTools::Generation::build_delaunay_square(), MeshTriangulationTest::commonSettings(), libMesh::Poly2TriTriangulator::insert_refinement_points(), main(), libMesh::Poly2TriTriangulator::should_refine_elem(), MeshTriangulationTest::testPoly2TriRefinementBase(), MeshTriangulationTest::testTriangulatorInterp(), and triangulate_domain().

177 {return _desired_area;}
Real _desired_area
The desired area for the elements in the resulting mesh.

◆ elem_type()

ElemType& libMesh::TriangulatorInterface::elem_type ( )
inline

Sets and/or gets the desired element type.

Definition at line 164 of file triangulator_interface.h.

References _elem_type.

Referenced by libMesh::MeshTools::Generation::build_delaunay_square(), and MeshTriangulationTest::testTriangulatorRoundHole().

164 {return _elem_type;}
ElemType _elem_type
The type of elements to generate.

◆ elems_to_segments()

void libMesh::TriangulatorInterface::elems_to_segments ( )
protected

Helper function to create PSLG segments from our other boundary-defining options (1D mesh edges, 2D mesh boundary sides), if no segments already exist.

Definition at line 149 of file triangulator_interface.C.

References _bdy_ids, _holes, _mesh, _verify_hole_boundaries, libMesh::MeshBase::clear_elems(), libMesh::MeshBase::delete_node(), libMesh::make_range(), libMesh::MeshBase::n_elem(), libMesh::MeshBase::node_ptr(), segment_midpoints, segments, and verify_holes().

Referenced by libMesh::TriangleInterface::triangulate(), and libMesh::Poly2TriTriangulator::triangulate().

150 {
151  // Don't try to override manually specified segments
152  if (!this->segments.empty())
153  return;
154 
155  // If we have edges, they should form the polyline with the ordering
156  // we want. Let's turn them into segments for later use, because
157  // we're going to delete the original elements to replace with our
158  // triangulation.
159  if (_mesh.n_elem())
160  {
161  // Mapping from points to node ids, to back those out from
162  // MeshedHole results later
163  std::map<Point, dof_id_type> point_id_map;
164 
165  for (Node * node : _mesh.node_ptr_range())
166  {
167  // We're not going to support overlapping nodes on the boundary
168  libmesh_error_msg_if
169  (point_id_map.count(*node),
170  "TriangulatorInterface does not support overlapping nodes found at "
171  << static_cast<Point&>(*node));
172 
173  point_id_map.emplace(*node, node->id());
174  }
175 
176  // We don't support directly generating Tri6, so for
177  // compatibility with future stitching we need to be working
178  // with first-order elements. Let's get rid of any non-vertex
179  // nodes we just added.
180  for (Elem * elem : _mesh.element_ptr_range())
181  for (auto n : make_range(elem->n_vertices(), elem->n_nodes()))
182  point_id_map.erase(elem->point(n));
183 
184  // We'll steal the ordering calculation from
185  // the MeshedHole code
186  const TriangulatorInterface::MeshedHole mh { _mesh, this->_bdy_ids };
187 
188  // If we've specified only a subset of the mesh as our outer
189  // boundary, then we may have nodes that don't actually fall
190  // inside that boundary. Triangulator code doesn't like Steiner
191  // points that aren't inside the triangulation domain, so we
192  // need to get rid of them.
193  //
194  // Also, if we're using Edge3 elements to define our outer
195  // boundary, we're only dealing with their 2 end nodes and we'll
196  // need to get rid of their central nodes.
197  std::unordered_set<Node *> nodes_to_delete;
198 
199  for (Elem * elem : _mesh.element_ptr_range())
200  for (auto n : make_range(elem->n_vertices(), elem->n_nodes()))
201  nodes_to_delete.insert(elem->node_ptr(n));
202 
203  if (!this->_bdy_ids.empty())
204  {
205  for (auto & node : _mesh.node_ptr_range())
206  if (!mh.contains(*node))
207  nodes_to_delete.insert(node);
208  }
209 
210  // And now we're done with elements. Delete them lest they have
211  // dangling pointers to nodes we'll be deleting.
212  _mesh.clear_elems();
213 
214  // Make segments from boundary nodes; also make sure we don't
215  // delete them.
216  const std::size_t np = mh.n_points();
217  for (auto i : make_range(np))
218  {
219  const Point pt = mh.point(i);
220  const dof_id_type id0 = libmesh_map_find(point_id_map, pt);
221  nodes_to_delete.erase(_mesh.node_ptr(id0));
222  const Point next_pt = mh.point((np+i+1)%np);
223  const dof_id_type id1 = libmesh_map_find(point_id_map, next_pt);
224  this->segments.emplace_back(id0, id1);
225  for (auto m : make_range(mh.n_midpoints()))
226  this->segment_midpoints.emplace_back(mh.midpoint(m, i));
227  }
228 
229  for (Node * node : nodes_to_delete)
230  _mesh.delete_node(node);
231 
232  if (this->_verify_hole_boundaries && _holes)
233  this->verify_holes(mh);
234  }
235 }
bool _verify_hole_boundaries
Flag which tells if we want to check hole geometry.
std::set< std::size_t > _bdy_ids
A set of ids to allow on the outer boundary loop.
UnstructuredMesh & _mesh
Reference to the mesh which is to be created by triangle.
const std::vector< Hole * > * _holes
A pointer to a vector of Hole*s.
virtual void delete_node(Node *n)=0
Removes the Node n from the mesh.
std::vector< Point > segment_midpoints
When constructing a second-order triangulation from a second-order boundary, we may do the triangulat...
virtual void clear_elems()=0
Deletes all the element data that is currently stored.
std::vector< std::pair< unsigned int, unsigned int > > segments
When constructing a PSLG, if the node numbers do not define the desired boundary segments implicitly ...
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
void verify_holes(const Hole &outer_bdy)
Helper function to check holes for intersections if requested.
virtual dof_id_type n_elem() const =0
virtual const Node * node_ptr(const dof_id_type i) const =0
uint8_t dof_id_type
Definition: id_types.h:67

◆ get_auto_area_function()

FunctionBase< Real > * libMesh::TriangulatorInterface::get_auto_area_function ( )

Get the auto area function.

Definition at line 487 of file triangulator_interface.C.

References _auto_area_function, and calculate_auto_desired_area_samples().

Referenced by libMesh::Poly2TriTriangulator::should_refine_elem().

488 {
489  if (!_auto_area_function->initialized())
490  {
491  // Points and target element sizes for the interpolation
492  std::vector<Point> function_points;
493  std::vector<Real> function_sizes;
494  calculate_auto_desired_area_samples(function_points, function_sizes);
495  _auto_area_function->init_mfi(function_points, function_sizes);
496  }
497  return _auto_area_function.get();
498 }
void calculate_auto_desired_area_samples(std::vector< Point > &function_points, std::vector< Real > &function_sizes, const Real &area_factor=1.5)
The external boundary and all hole boundaries are collected.
std::unique_ptr< AutoAreaFunction > _auto_area_function
The auto area function based on the spacing of boundary points.

◆ get_desired_area_function()

virtual FunctionBase<Real>* libMesh::TriangulatorInterface::get_desired_area_function ( )
inlinevirtual

Get the function giving desired triangle area as a function of position, or nullptr if no such function has been set.

Reimplemented in libMesh::Poly2TriTriangulator.

Definition at line 193 of file triangulator_interface.h.

194  { return nullptr; }

◆ get_interpolate_boundary_points()

int libMesh::TriangulatorInterface::get_interpolate_boundary_points ( ) const

Complicated getter, for compatibility with insert_extra_points()

Definition at line 137 of file triangulator_interface.C.

References _insert_extra_points, and _interpolate_boundary_points.

Referenced by insert_any_extra_boundary_points().

138 {
139  // backwards compatibility - someone might have turned this off via
140  // the old API
142  return 0;
143 
145 }
bool _insert_extra_points
Flag which tells whether or not to insert additional nodes before triangulation.
int _interpolate_boundary_points
Flag which tells how many additional nodes should be inserted between each pair of original mesh poin...

◆ get_outer_boundary_ids()

const std::set<std::size_t>& libMesh::TriangulatorInterface::get_outer_boundary_ids ( ) const
inline

Definition at line 352 of file triangulator_interface.h.

References _bdy_ids.

352 { return _bdy_ids; }
std::set< std::size_t > _bdy_ids
A set of ids to allow on the outer boundary loop.

◆ get_verify_hole_boundaries()

bool libMesh::TriangulatorInterface::get_verify_hole_boundaries ( ) const
inline

Definition at line 264 of file triangulator_interface.h.

References _verify_hole_boundaries.

264 {return _verify_hole_boundaries;}
bool _verify_hole_boundaries
Flag which tells if we want to check hole geometry.

◆ has_auto_area_function()

bool libMesh::TriangulatorInterface::has_auto_area_function ( )
inline

Whether or not an auto area function has been set.

Definition at line 329 of file triangulator_interface.h.

References _auto_area_function.

Referenced by libMesh::Poly2TriTriangulator::insert_refinement_points(), and libMesh::Poly2TriTriangulator::should_refine_elem().

329 {return _auto_area_function != nullptr;}
std::unique_ptr< AutoAreaFunction > _auto_area_function
The auto area function based on the spacing of boundary points.

◆ increase_triangle_order()

void libMesh::TriangulatorInterface::increase_triangle_order ( )
protected

Helper function to upconvert Tri3 to any higher order triangle type if requested via _elem_type.

Should be called at the end of triangulate()

Definition at line 332 of file triangulator_interface.C.

References _elem_type, _holes, _mesh, libMesh::MeshBase::all_complete_order(), libMesh::MeshBase::all_second_order(), libMesh::FIRST, libMesh::libmesh_assert(), libMesh::make_range(), libMesh::MeshBase::point(), segment_midpoints, segments, libMesh::TRI3, libMesh::TRI6, and libMesh::TRI7.

Referenced by libMesh::TriangleInterface::triangulate(), and libMesh::Poly2TriTriangulator::triangulate().

333 {
334  switch (_elem_type)
335  {
336  case TRI3:
337  // Nothing to do if we're not requested to increase order
338  return;
339  case TRI6:
341  break;
342  case TRI7:
344  break;
345  default:
346  libmesh_not_implemented();
347  }
348 
349  // If we have any midpoint location data, we'll want to look it up
350  // by point. all_midpoints[{p, m}] will be the mth midpoint
351  // location following after point p (when traversing a triangle
352  // counter-clockwise)
353  std::map<std::pair<Point, unsigned int>, Point> all_midpoints;
354  unsigned int n_midpoints =
355  this->segment_midpoints.size() / this->segments.size();
356  libmesh_assert_equal_to(this->segments.size() * n_midpoints,
357  this->segment_midpoints.size());
358  for (auto m : make_range(n_midpoints))
359  for (auto i : make_range(this->segments.size()))
360  {
361  const Point & p = _mesh.point(this->segments[i].first);
362  all_midpoints[{p,m}] =
363  this->segment_midpoints[i*n_midpoints+m];
364  }
365 
366  if (_holes)
367  for (const Hole * hole : *_holes)
368  {
369  if (!hole->n_midpoints())
370  continue;
371  if (!n_midpoints)
372  n_midpoints = hole->n_midpoints();
373  else if (hole->n_midpoints() != n_midpoints)
374  libmesh_not_implemented_msg
375  ("Differing boundary midpoint counts " <<
376  hole->n_midpoints() << " and " << n_midpoints);
377 
378  // Our inner holes are expected to have points in
379  // counter-clockwise order, which is backwards from how we
380  // want to traverse them when iterating in counter-clockwise
381  // order over a triangle, so we'll need to reverse our maps
382  // carefully here.
383  const auto n_hole_points = hole->n_points();
384  libmesh_assert(n_hole_points);
385  for (auto m : make_range(n_midpoints))
386  {
387  for (auto i : make_range(n_hole_points-1))
388  {
389  const Point & p = hole->point(i+1);
390  all_midpoints[{p,m}] = hole->midpoint(n_midpoints-m-1, i);
391  }
392  const Point & p = hole->point(0);
393  all_midpoints[{p,m}] = hole->midpoint(n_midpoints-m-1, n_hole_points-1);
394  }
395  }
396 
397  // The n_midpoints > 1 case is for future proofing, but in the
398  // present we have EDGE4 and no TRI10 yet.
399  if (n_midpoints > 1)
400  libmesh_not_implemented_msg
401  ("Cannot construct triangles with more than 1 midpoint per edge");
402 
403  if (!n_midpoints)
404  return;
405 
406  for (Elem * elem : _mesh.element_ptr_range())
407  {
408  // This should only be called right after we've finished
409  // converting a triangulation to higher order
410  libmesh_assert_equal_to(elem->n_vertices(), 3);
411  libmesh_assert_not_equal_to(elem->default_order(), FIRST);
412 
413  for (auto n : make_range(3))
414  {
415  // Only hole/outer boundary segments need adjusted midpoints
416  if (elem->neighbor_ptr(n))
417  continue;
418 
419  const Point & p = elem->point(n);
420 
421  if (const auto it = all_midpoints.find({p,0});
422  it != all_midpoints.end())
423  elem->point(n+3) = it->second;
424  }
425  }
426 }
UnstructuredMesh & _mesh
Reference to the mesh which is to be created by triangle.
ElemType _elem_type
The type of elements to generate.
const std::vector< Hole * > * _holes
A pointer to a vector of Hole*s.
libmesh_assert(ctx)
std::vector< Point > segment_midpoints
When constructing a second-order triangulation from a second-order boundary, we may do the triangulat...
std::vector< std::pair< unsigned int, unsigned int > > segments
When constructing a PSLG, if the node numbers do not define the desired boundary segments implicitly ...
virtual void all_complete_order()
Calls the range-based version of this function with a range consisting of all elements in the mesh...
Definition: mesh_base.C:1613
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
void all_second_order(const bool full_ordered=true)
Calls the range-based version of this function with a range consisting of all elements in the mesh...
Definition: mesh_base.C:1608
virtual const Point & point(const dof_id_type i) const =0

◆ insert_any_extra_boundary_points()

void libMesh::TriangulatorInterface::insert_any_extra_boundary_points ( )
protected

Helper function to add extra points (midpoints of initial segments) to a PSLG triangulation.

Definition at line 278 of file triangulator_interface.C.

References _mesh, _triangulation_type, libMesh::MeshBase::add_point(), get_interpolate_boundary_points(), libMesh::DofObject::id(), libMesh::libmesh_assert(), libMesh::make_range(), libMesh::MeshBase::max_node_id(), libMesh::MeshBase::node_ptr(), PSLG, and segments.

Referenced by libMesh::TriangleInterface::triangulate(), and libMesh::Poly2TriTriangulator::triangulate().

279 {
280  // If the initial PSLG is really simple, e.g. an L-shaped domain or
281  // a square/rectangle, the resulting triangulation may be very
282  // "structured" looking. Sometimes this is a problem if your
283  // intention is to work with an "unstructured" looking grid. We can
284  // attempt to work around this limitation by inserting midpoints
285  // into the original PSLG. Inserting additional points into a
286  // set of points meant to be a convex hull usually makes less sense.
287 
288  const int n_interpolated = this->get_interpolate_boundary_points();
289  if ((_triangulation_type==PSLG) && n_interpolated)
290  {
291  // If we were lucky enough to start with contiguous node ids,
292  // let's keep them that way.
294 
295  std::vector<std::pair<unsigned int, unsigned int>> old_segments =
296  std::move(this->segments);
297 
298  // We expect to have converted any elems and/or nodes into
299  // segments by now.
300  libmesh_assert(!old_segments.empty());
301 
302  this->segments.clear();
303 
304  // Insert a new point on each segment at evenly spaced locations
305  // between existing boundary points.
306  // np=index into new points vector
307  // n =index into original points vector
308  for (auto old_segment : old_segments)
309  {
310  Node * begin_node = _mesh.node_ptr(old_segment.first);
311  Node * end_node = _mesh.node_ptr(old_segment.second);
312  dof_id_type current_id = begin_node->id();
313  for (auto i : make_range(n_interpolated))
314  {
315  // new points are equispaced along the original segments
316  const Point new_point =
317  ((n_interpolated-i) * *(Point *)(begin_node) +
318  (i+1) * *(Point *)(end_node)) /
319  (n_interpolated + 1);
320  Node * next_node = _mesh.add_point(new_point, nn++);
321  this->segments.emplace_back(current_id,
322  next_node->id());
323  current_id = next_node->id();
324  }
325  this->segments.emplace_back(current_id,
326  end_node->id());
327  }
328  }
329 }
UnstructuredMesh & _mesh
Reference to the mesh which is to be created by triangle.
virtual Node * add_point(const Point &p, const dof_id_type id=DofObject::invalid_id, const processor_id_type proc_id=DofObject::invalid_processor_id)=0
Add a new Node at Point p to the end of the vertex array, with processor_id procid.
libmesh_assert(ctx)
Triangulate the interior of a Planar Straight Line Graph, which is defined implicitly by the order of...
std::vector< std::pair< unsigned int, unsigned int > > segments
When constructing a PSLG, if the node numbers do not define the desired boundary segments implicitly ...
int get_interpolate_boundary_points() const
Complicated getter, for compatibility with insert_extra_points()
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 dof_id_type max_node_id() const =0
virtual const Node * node_ptr(const dof_id_type i) const =0
TriangulationType _triangulation_type
The type of triangulation to perform: choices are: convex hull PSLG.
uint8_t dof_id_type
Definition: id_types.h:67

◆ insert_extra_points()

bool& libMesh::TriangulatorInterface::insert_extra_points ( )
inline

Sets and/or gets the flag for inserting add'l points.

Definition at line 210 of file triangulator_interface.h.

References _insert_extra_points.

210 {return _insert_extra_points;}
bool _insert_extra_points
Flag which tells whether or not to insert additional nodes before triangulation.

◆ minimum_angle()

Real& libMesh::TriangulatorInterface::minimum_angle ( )
inline

Sets and/or gets the minimum desired angle.

Set to zero to disable angle constraint.

Definition at line 200 of file triangulator_interface.h.

References _minimum_angle.

Referenced by MeshTriangulationTest::commonSettings(), libMesh::Poly2TriTriangulator::insert_refinement_points(), and main().

200 {return _minimum_angle;}
Real _minimum_angle
Minimum angle in triangles.

◆ nodes_to_segments()

void libMesh::TriangulatorInterface::nodes_to_segments ( dof_id_type  max_node_id)
protected

Helper function to create PSLG segments from our node ordering, up to the maximum node id, if no segments already exist.

Definition at line 239 of file triangulator_interface.C.

References _holes, _mesh, _triangulation_type, _verify_hole_boundaries, libMesh::DofObject::id(), libMesh::MeshBase::point(), PSLG, segments, and verify_holes().

Referenced by libMesh::TriangleInterface::triangulate(), and libMesh::Poly2TriTriangulator::triangulate().

240 {
241  // Don't try to override manually specified segments, or try to add
242  // segments if we're doing a convex hull
243  if (!this->segments.empty() || _triangulation_type != PSLG)
244  return;
245 
246  for (auto node_it = _mesh.nodes_begin(),
247  node_end = _mesh.nodes_end();
248  node_it != node_end;)
249  {
250  Node * node = *node_it;
251 
252  // If we're out of boundary nodes, the rest are going to be
253  // Steiner points or hole points
254  if (node->id() >= max_node_id)
255  break;
256 
257  ++node_it;
258 
259  Node * next_node = (node_it == node_end) ?
260  *_mesh.nodes_begin() : *node_it;
261 
262  this->segments.emplace_back(node->id(), next_node->id());
263  }
264 
265  if (this->_verify_hole_boundaries && _holes)
266  {
267  std::vector<Point> outer_pts;
268  for (auto segment : this->segments)
269  outer_pts.push_back(_mesh.point(segment.first));
270 
271  ArbitraryHole ah(outer_pts);
272  this->verify_holes(ah);
273  }
274 }
bool _verify_hole_boundaries
Flag which tells if we want to check hole geometry.
UnstructuredMesh & _mesh
Reference to the mesh which is to be created by triangle.
const std::vector< Hole * > * _holes
A pointer to a vector of Hole*s.
Triangulate the interior of a Planar Straight Line Graph, which is defined implicitly by the order of...
std::vector< std::pair< unsigned int, unsigned int > > segments
When constructing a PSLG, if the node numbers do not define the desired boundary segments implicitly ...
virtual const Point & point(const dof_id_type i) const =0
void verify_holes(const Hole &outer_bdy)
Helper function to check holes for intersections if requested.
TriangulationType _triangulation_type
The type of triangulation to perform: choices are: convex hull PSLG.

◆ quiet()

bool& libMesh::TriangulatorInterface::quiet ( )
inline

Whether not to silence internal messages to stdout.

Definition at line 249 of file triangulator_interface.h.

References _quiet.

249 {return _quiet;}
bool _quiet
Flag which tells if we want to suppress stdout outputs.

◆ refine_boundary_allowed()

virtual bool libMesh::TriangulatorInterface::refine_boundary_allowed ( ) const
inlinevirtual

Get whether or not the triangulation is allowed to refine the mesh boundary when refining the interior.

True by default.

Reimplemented in libMesh::Poly2TriTriangulator.

Definition at line 237 of file triangulator_interface.h.

238  { return true; }

◆ set_auto_area_function()

void libMesh::TriangulatorInterface::set_auto_area_function ( const Parallel::Communicator comm,
const unsigned int  num_nearest_pts,
const unsigned int  power,
const Real  background_value,
const Real  background_eff_dist 
)

Generate an auto area function based on spacing of boundary points.

The external boundary as well as the hole boundaries are taken into consideration to generate the auto area function based on inverse distance interpolation. For each EDGE element on these boundaries, its centroid (midpoint) is used as the point position and the desired area is calculated as 1.5 times of the area of the equilateral triangle with the edge length as the length of the EDGE element. For a given position, the inverse distance interpolation only considers a number of nearest points (set by num_nearest_pts) to calculate the desired area. The weight of the value at each point is calculated as 1/distance^power.

In addition to these conventional inverse distance interpolation features, a concept of "background value" and "background effective distance" is introduced. The background value belongs to a virtual point located at a constant distance (background effective distance) from the given position. The weight of the value at this virtual point is calculated as 1/background_effective_distance^power. Effectively, the background value is the value when the given position is far away from the boundary points.

Definition at line 478 of file triangulator_interface.C.

References _auto_area_function.

483 {
484  _auto_area_function = std::make_unique<AutoAreaFunction>(comm, num_nearest_pts, power, background_value, background_eff_dist);
485 }
std::unique_ptr< AutoAreaFunction > _auto_area_function
The auto area function based on the spacing of boundary points.

◆ set_desired_area_function()

virtual void libMesh::TriangulatorInterface::set_desired_area_function ( FunctionBase< Real > *  )
inlinevirtual

Set a function giving desired triangle area as a function of position.

Set this to nullptr to disable position-dependent area constraint (falling back on desired_area()).

This may not be implemented in all subclasses.

Reimplemented in libMesh::Poly2TriTriangulator.

Definition at line 186 of file triangulator_interface.h.

187  { libmesh_not_implemented(); }

◆ set_interpolate_boundary_points()

void libMesh::TriangulatorInterface::set_interpolate_boundary_points ( int  n_points)

Complicated setter, for compatibility with insert_extra_points()

Definition at line 123 of file triangulator_interface.C.

References _insert_extra_points, _interpolate_boundary_points, and libMesh::libmesh_assert().

Referenced by MeshTriangulationTest::testTriangulatorInterp().

124 {
125  // Maybe we'll reserve a meaning for negatives later?
126  libmesh_assert(n_points >= 0);
127 
128  _interpolate_boundary_points = n_points;
129 
130  // backwards compatibility - someone (including us) might want to
131  // query this via the old API.
132  _insert_extra_points = n_points;
133 }
bool _insert_extra_points
Flag which tells whether or not to insert additional nodes before triangulation.
int _interpolate_boundary_points
Flag which tells how many additional nodes should be inserted between each pair of original mesh poin...
libmesh_assert(ctx)

◆ set_outer_boundary_ids()

void libMesh::TriangulatorInterface::set_outer_boundary_ids ( std::set< std::size_t >  bdy_ids)
inline

A set of ids to allow on the outer boundary loop: interpreted as boundary ids of 2D elements and/or subdomain ids of 1D edges.

If this is empty, then the outer boundary may be constructed from boundary edges of any id!

Definition at line 351 of file triangulator_interface.h.

References _bdy_ids.

Referenced by MeshTriangulationTest::testHalfDomain().

351 { _bdy_ids = std::move(bdy_ids); }
std::set< std::size_t > _bdy_ids
A set of ids to allow on the outer boundary loop.

◆ set_refine_boundary_allowed()

virtual void libMesh::TriangulatorInterface::set_refine_boundary_allowed ( bool  )
inlinevirtual

Set whether or not the triangulation is allowed to refine the mesh boundary when refining the interior.

This is true by default, but may be set to false to make the mesh boundary more predictable (and so easier to stitch to other meshes) later.

This may not be implemented in all subclasses.

Reimplemented in libMesh::Poly2TriTriangulator.

Definition at line 230 of file triangulator_interface.h.

231  { libmesh_not_implemented(); }

◆ set_verify_hole_boundaries()

void libMesh::TriangulatorInterface::set_verify_hole_boundaries ( bool  v)
inline

Verifying that hole boundaries don't cross the outer boundary or each other is something like O(N_bdys^2*N_points_per_bdy^2), so we only do it if requested.

Definition at line 262 of file triangulator_interface.h.

References _verify_hole_boundaries.

Referenced by MeshTriangulationTest::commonSettings().

bool _verify_hole_boundaries
Flag which tells if we want to check hole geometry.

◆ smooth_after_generating()

bool& libMesh::TriangulatorInterface::smooth_after_generating ( )
inline

Sets/gets flag which tells whether to do two steps of Laplace mesh smoothing after generating the grid.

Definition at line 244 of file triangulator_interface.h.

References _smooth_after_generating.

Referenced by MeshTriangulationTest::commonSettings(), and triangulate_domain().

244 {return _smooth_after_generating;}
bool _smooth_after_generating
Flag which tells whether we should smooth the mesh after it is generated.

◆ total_hole_points()

unsigned int libMesh::TriangulatorInterface::total_hole_points ( )
protected

Helper function to count points in and verify holes.

Definition at line 454 of file triangulator_interface.C.

References _holes.

Referenced by libMesh::TriangleInterface::triangulate().

455 {
456  // If the holes vector is non-nullptr (and non-empty) we need to determine
457  // the number of additional points which the holes will add to the
458  // triangulation.
459  // Note that the number of points is always equal to the number of segments
460  // that form the holes.
461  unsigned int n_hole_points = 0;
462 
463  if (_holes)
464  for (const auto & hole : *_holes)
465  {
466  n_hole_points += hole->n_points();
467  // A hole at least has one enclosure.
468  // Points on enclosures are ordered so that we can add segments implicitly.
469  // Elements in segment_indices() indicates the starting points of all enclosures.
470  // The last element in segment_indices() is the number of total points.
471  libmesh_assert_greater(hole->segment_indices().size(), 1);
472  libmesh_assert_equal_to(hole->segment_indices().back(), hole->n_points());
473  }
474 
475  return n_hole_points;
476 }
const std::vector< Hole * > * _holes
A pointer to a vector of Hole*s.

◆ triangulate()

virtual void libMesh::TriangulatorInterface::triangulate ( )
pure virtual

◆ triangulation_type()

TriangulationType& libMesh::TriangulatorInterface::triangulation_type ( )
inline

Sets and/or gets the desired triangulation type.

Definition at line 205 of file triangulator_interface.h.

References _triangulation_type.

Referenced by libMesh::MeshTools::Generation::build_delaunay_square(), MeshTriangulationTest::commonSettings(), main(), and triangulate_domain().

205 {return _triangulation_type;}
TriangulationType _triangulation_type
The type of triangulation to perform: choices are: convex hull PSLG.

◆ verify_holes()

void libMesh::TriangulatorInterface::verify_holes ( const Hole outer_bdy)
protected

Helper function to check holes for intersections if requested.

Definition at line 429 of file triangulator_interface.C.

References _holes, libMesh::TriangulatorInterface::Hole::contains(), and libMesh::make_range().

Referenced by elems_to_segments(), and nodes_to_segments().

430 {
431  for (const Hole * hole : *_holes)
432  {
433  for (const Hole * hole2 : *_holes)
434  {
435  if (hole == hole2)
436  continue;
437 
438  for (auto i : make_range(hole2->n_points()))
439  if (hole->contains(hole2->point(i)))
440  libmesh_error_msg
441  ("Found point " << hole2->point(i) <<
442  " on one hole boundary and another's interior");
443  }
444 
445  for (auto i : make_range(hole->n_points()))
446  if (!outer_bdy.contains(hole->point(i)))
447  libmesh_error_msg
448  ("Found point " << hole->point(i) <<
449  " on hole boundary but outside outer boundary");
450  }
451 }
const std::vector< Hole * > * _holes
A pointer to a vector of Hole*s.
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

Member Data Documentation

◆ _auto_area_function

std::unique_ptr<AutoAreaFunction> libMesh::TriangulatorInterface::_auto_area_function
protected

The auto area function based on the spacing of boundary points.

Definition at line 476 of file triangulator_interface.h.

Referenced by get_auto_area_function(), has_auto_area_function(), and set_auto_area_function().

◆ _bdy_ids

std::set<std::size_t> libMesh::TriangulatorInterface::_bdy_ids
protected

A set of ids to allow on the outer boundary loop.

Definition at line 416 of file triangulator_interface.h.

Referenced by calculate_auto_desired_area_samples(), elems_to_segments(), get_outer_boundary_ids(), and set_outer_boundary_ids().

◆ _desired_area

Real libMesh::TriangulatorInterface::_desired_area
protected

The desired area for the elements in the resulting mesh.

Definition at line 427 of file triangulator_interface.h.

Referenced by desired_area(), and libMesh::TriangleInterface::triangulate().

◆ _elem_type

ElemType libMesh::TriangulatorInterface::_elem_type
protected

The type of elements to generate.

(Defaults to TRI3).

Definition at line 422 of file triangulator_interface.h.

Referenced by elem_type(), increase_triangle_order(), libMesh::TriangleInterface::triangulate(), and libMesh::Poly2TriTriangulator::triangulate().

◆ _holes

const std::vector<Hole*>* libMesh::TriangulatorInterface::_holes
protected

◆ _insert_extra_points

bool libMesh::TriangulatorInterface::_insert_extra_points
protected

Flag which tells whether or not to insert additional nodes before triangulation.

This can sometimes be used to "de-regularize" the resulting triangulation.

This flag is supported for backwards compatibility; setting _interpolate_boundary_points = 1 is equivalent.

Definition at line 449 of file triangulator_interface.h.

Referenced by get_interpolate_boundary_points(), insert_extra_points(), and set_interpolate_boundary_points().

◆ _interpolate_boundary_points

int libMesh::TriangulatorInterface::_interpolate_boundary_points
protected

Flag which tells how many additional nodes should be inserted between each pair of original mesh points.

Definition at line 455 of file triangulator_interface.h.

Referenced by get_interpolate_boundary_points(), and set_interpolate_boundary_points().

◆ _markers

const std::vector<int>* libMesh::TriangulatorInterface::_markers
protected

◆ _mesh

UnstructuredMesh& libMesh::TriangulatorInterface::_mesh
protected

◆ _minimum_angle

Real libMesh::TriangulatorInterface::_minimum_angle
protected

Minimum angle in triangles.

Definition at line 432 of file triangulator_interface.h.

Referenced by minimum_angle(), and libMesh::TriangleInterface::triangulate().

◆ _quiet

bool libMesh::TriangulatorInterface::_quiet
protected

Flag which tells if we want to suppress stdout outputs.

Definition at line 466 of file triangulator_interface.h.

Referenced by quiet(), and libMesh::TriangleInterface::triangulate().

◆ _regions

const std::vector<Region*>* libMesh::TriangulatorInterface::_regions
protected

A pointer to a vector of Regions*s.

If this is nullptr, there are no regions!

Definition at line 411 of file triangulator_interface.h.

Referenced by attach_region_list(), libMesh::TriangleInterface::triangulate(), and libMesh::Poly2TriTriangulator::triangulate().

◆ _smooth_after_generating

bool libMesh::TriangulatorInterface::_smooth_after_generating
protected

Flag which tells whether we should smooth the mesh after it is generated.

True by default.

Definition at line 461 of file triangulator_interface.h.

Referenced by smooth_after_generating(), libMesh::TriangleInterface::triangulate(), and libMesh::Poly2TriTriangulator::triangulate().

◆ _triangulation_type

TriangulationType libMesh::TriangulatorInterface::_triangulation_type
protected

The type of triangulation to perform: choices are: convex hull PSLG.

Definition at line 439 of file triangulator_interface.h.

Referenced by insert_any_extra_boundary_points(), nodes_to_segments(), libMesh::TriangleInterface::triangulate(), libMesh::Poly2TriTriangulator::triangulate(), and triangulation_type().

◆ _verify_hole_boundaries

bool libMesh::TriangulatorInterface::_verify_hole_boundaries
protected

Flag which tells if we want to check hole geometry.

Definition at line 471 of file triangulator_interface.h.

Referenced by elems_to_segments(), get_verify_hole_boundaries(), nodes_to_segments(), and set_verify_hole_boundaries().

◆ segment_midpoints

std::vector<Point> libMesh::TriangulatorInterface::segment_midpoints

When constructing a second-order triangulation from a second-order boundary, we may do the triangulation using first-order elements, in which case we need to save midpoint location data in order to reconstruct curvature along boundaries.

Definition at line 287 of file triangulator_interface.h.

Referenced by elems_to_segments(), and increase_triangle_order().

◆ segments

std::vector<std::pair<unsigned int, unsigned int> > libMesh::TriangulatorInterface::segments

When constructing a PSLG, if the node numbers do not define the desired boundary segments implicitly through the ordering of the points, you can use the segments vector to specify the segments explicitly, Ex: unit square numbered counter-clockwise starting from origin segments[0] = (0,1) segments[1] = (1,2) segments[2] = (2,3) segments[3] = (3,0) (For the above case you could actually use the implicit ordering!)

Definition at line 279 of file triangulator_interface.h.

Referenced by elems_to_segments(), increase_triangle_order(), insert_any_extra_boundary_points(), libMesh::Poly2TriTriangulator::insert_refinement_points(), nodes_to_segments(), MeshTriangulationTest::testTriangulatorSegments(), libMesh::TriangleInterface::triangulate(), and libMesh::Poly2TriTriangulator::triangulate_current_points().


The documentation for this class was generated from the following files: