libMesh::TriangleInterface Class Reference

#include <mesh_triangle_interface.h>

Classes

class  ArbitraryHole
 
class  Hole
 Class for parameterizing 2D holes to be meshed with Triangle. More...
 
class  PolygonHole
 

Public Types

enum  TriangulationType { GENERATE_CONVEX_HULL = 0, PSLG = 1, INVALID_TRIANGULATION_TYPE }
 

Public Member Functions

 TriangleInterface (UnstructuredMesh &mesh)
 
 ~TriangleInterface ()
 
void triangulate ()
 
ElemTypeelem_type ()
 
Realdesired_area ()
 
Realminimum_angle ()
 
std::string & extra_flags ()
 
TriangulationTypetriangulation_type ()
 
bool & insert_extra_points ()
 
bool & smooth_after_generating ()
 
void attach_hole_list (const std::vector< Hole *> *holes)
 

Public Attributes

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

Private Attributes

UnstructuredMesh_mesh
 
const std::vector< Hole * > * _holes
 
ElemType _elem_type
 
Real _desired_area
 
Real _minimum_angle
 
std::string _extra_flags
 
TriangulationType _triangulation_type
 
bool _insert_extra_points
 
bool _smooth_after_generating
 
MeshSerializer _serializer
 

Detailed Description

A C++ interface between LibMesh and the Triangle library written by J.R. Shewchuk.

Author
John W. Peterson
Date
2011

Definition at line 57 of file mesh_triangle_interface.h.

Member Enumeration Documentation

◆ TriangulationType

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

Enumerator
GENERATE_CONVEX_HULL 

Uses the triangle library to 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 

Use the triangle library to triangulate 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. In case your triangulation is a little too "structured" looking (which can happen when the initial PSLG is really simple) you can try to make the resulting triangulation a little more "unstructured" looking by setting insert_points to true in the triangulate() function.

INVALID_TRIANGULATION_TYPE 

Does nothing, used as a default value.

Definition at line 78 of file mesh_triangle_interface.h.

Constructor & Destructor Documentation

◆ TriangleInterface()

libMesh::TriangleInterface::TriangleInterface ( 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 44 of file mesh_triangle_interface.C.

45  : _mesh(mesh),
46  _holes(nullptr),
48  _desired_area(0.1),
49  _minimum_angle(20.0),
50  _extra_flags(""),
52  _insert_extra_points(false),
55 {}
MeshBase & mesh
const std::vector< Hole * > * _holes

◆ ~TriangleInterface()

libMesh::TriangleInterface::~TriangleInterface ( )
inline

Empty destructor.

Definition at line 72 of file mesh_triangle_interface.h.

72 {}

Member Function Documentation

◆ attach_hole_list()

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

Attaches a vector of Hole* pointers which will be meshed around.

Definition at line 162 of file mesh_triangle_interface.h.

References _holes.

Referenced by libMesh::MeshTools::Generation::build_delaunay_square().

162 {_holes = holes;}
const std::vector< Hole * > * _holes

◆ desired_area()

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

Sets and/or gets the desired triangle area. Set to zero to disable area constraint.

Definition at line 129 of file mesh_triangle_interface.h.

References _desired_area.

Referenced by libMesh::MeshTools::Generation::build_delaunay_square().

◆ elem_type()

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

Sets and/or gets the desired element type.

Definition at line 123 of file mesh_triangle_interface.h.

References _elem_type.

Referenced by libMesh::MeshTools::Generation::build_delaunay_square().

◆ extra_flags()

std::string& libMesh::TriangleInterface::extra_flags ( )
inline

Sets and/or gets additional flags to be passed to triangle

Definition at line 140 of file mesh_triangle_interface.h.

References _extra_flags.

◆ insert_extra_points()

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

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

Definition at line 150 of file mesh_triangle_interface.h.

References _insert_extra_points.

◆ minimum_angle()

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

Sets and/or gets the minimum angle. Set to zero to disable area constraint.

Definition at line 135 of file mesh_triangle_interface.h.

References _minimum_angle.

◆ smooth_after_generating()

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

Sets/gets flag which tells whether to do Delaunay mesh smoothing after generating the grid.

Definition at line 156 of file mesh_triangle_interface.h.

References _smooth_after_generating.

◆ triangulate()

void libMesh::TriangleInterface::triangulate ( )

This is the main public interface for this function. Internally, it calls Triangle's triangulate routine.

Definition at line 60 of file mesh_triangle_interface.C.

References _desired_area, _elem_type, _extra_flags, _holes, _insert_extra_points, _mesh, _minimum_angle, _smooth_after_generating, _triangulation_type, libMesh::MeshBase::add_point(), libMesh::MeshBase::clear(), libMesh::TriangleWrapper::copy_tri_to_mesh(), libMesh::TriangleWrapper::destroy(), GENERATE_CONVEX_HULL, libMesh::TriangleWrapper::init(), libMesh::TriangleWrapper::INPUT, INVALID_TRIANGULATION_TYPE, libMesh::MeshBase::n_nodes(), libMesh::MeshBase::node_ptr_range(), libMesh::TriangleWrapper::OUTPUT, libMesh::MeshBase::prepare_for_use(), PSLG, segments, libMesh::MeshBase::set_mesh_dimension(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::TOLERANCE, libMesh::TRI3, and libMesh::TRI6.

Referenced by libMesh::MeshTools::Generation::build_delaunay_square().

61 {
62  // Will the triangulation have holes?
63  const bool have_holes = ((_holes != nullptr) && (!_holes->empty()));
64 
65  // If the initial PSLG is really simple, e.g. an L-shaped domain or
66  // a square/rectangle, the resulting triangulation may be very
67  // "structured" looking. Sometimes this is a problem if your
68  // intention is to work with an "unstructured" looking grid. We can
69  // attempt to work around this limitation by inserting midpoints
70  // into the original PSLG. Inserting additional points into a
71  // set of points meant to be a convex hull usually makes less sense.
72 
73  // May or may not need to insert new points ...
75  {
76  // Make a copy of the original points from the Mesh
77  std::vector<Point> original_points;
78  original_points.reserve (_mesh.n_nodes());
79  for (auto & node : _mesh.node_ptr_range())
80  original_points.push_back(*node);
81 
82  // Clear out the mesh
83  _mesh.clear();
84 
85  // Make sure the new Mesh will be 2D
87 
88  // Insert a new point on each PSLG at some random location
89  // np=index into new points vector
90  // n =index into original points vector
91  for (std::size_t np=0, n=0; np<2*original_points.size(); ++np)
92  {
93  // the even entries are the original points
94  if (np%2==0)
95  _mesh.add_point(original_points[n++]);
96 
97  else // the odd entries are the midpoints of the original PSLG segments
98  _mesh.add_point ((original_points[n] + original_points[n-1])/2);
99  }
100  }
101 
102  // Regardless of whether we added additional points, the set of points to
103  // triangulate is now sitting in the mesh.
104 
105  // If the holes vector is non-nullptr (and non-empty) we need to determine
106  // the number of additional points which the holes will add to the
107  // triangulation.
108  unsigned int n_hole_points = 0;
109 
110  if (have_holes)
111  {
112  for (std::size_t i=0; i<_holes->size(); ++i)
113  n_hole_points += (*_holes)[i]->n_points();
114  }
115 
116  // Triangle data structure for the mesh
117  TriangleWrapper::triangulateio initial;
118  TriangleWrapper::triangulateio final;
119 
120  // Pseudo-Constructor for the triangle io structs
121  TriangleWrapper::init(initial);
122  TriangleWrapper::init(final);
123 
124  initial.numberofpoints = _mesh.n_nodes() + n_hole_points;
125  initial.pointlist = static_cast<REAL*>(std::malloc(initial.numberofpoints * 2 * sizeof(REAL)));
126 
128  {
129  // Implicit segment ordering: One segment per point, including hole points
130  if (this->segments.empty())
131  initial.numberofsegments = initial.numberofpoints;
132 
133  // User-defined segment ordering: One segment per entry in the segments vector
134  else
135  initial.numberofsegments = this->segments.size();
136  }
137 
139  initial.numberofsegments = n_hole_points; // One segment for each hole point
140 
141  // Debugging
142  // libMesh::out << "Number of segments set to: " << initial.numberofsegments << std::endl;
143 
144  // Allocate space for the segments (2 int per segment)
145  if (initial.numberofsegments > 0)
146  {
147  initial.segmentlist = static_cast<int *> (std::malloc(initial.numberofsegments * 2 * sizeof(int)));
148  }
149 
150 
151  // Copy all the holes' points and segments into the triangle struct.
152 
153  // The hole_offset is a constant offset into the points vector which points
154  // past the end of the last hole point added.
155  unsigned int hole_offset=0;
156 
157  if (have_holes)
158  for (std::size_t i=0; i<_holes->size(); ++i)
159  {
160  for (unsigned int ctr=0, h=0; h<(*_holes)[i]->n_points(); ctr+=2, ++h)
161  {
162  Point p = (*_holes)[i]->point(h);
163 
164  const unsigned int index0 = 2*hole_offset+ctr;
165  const unsigned int index1 = 2*hole_offset+ctr+1;
166 
167  // Save the x,y locations in the triangle struct.
168  initial.pointlist[index0] = p(0);
169  initial.pointlist[index1] = p(1);
170 
171  // Set the points which define the segments
172  initial.segmentlist[index0] = hole_offset+h;
173  initial.segmentlist[index1] = (h==(*_holes)[i]->n_points()-1) ? hole_offset : hole_offset+h+1; // wrap around
174  }
175 
176  // Update the hole_offset for the next hole
177  hole_offset += (*_holes)[i]->n_points();
178  }
179 
180 
181  // Copy all the non-hole points and segments into the triangle struct.
182  {
183  dof_id_type ctr=0;
184  for (auto & node : _mesh.node_ptr_range())
185  {
186  dof_id_type index = 2*hole_offset + ctr;
187 
188  // Set x,y values in pointlist
189  initial.pointlist[index] = (*node)(0);
190  initial.pointlist[index+1] = (*node)(1);
191 
192  // If the user requested a PSLG, the non-hole points are also segments
194  {
195  // Use implicit ordering to define segments
196  if (this->segments.empty())
197  {
198  dof_id_type n = ctr/2; // ctr is always even
199  initial.segmentlist[index] = hole_offset+n;
200  initial.segmentlist[index+1] = (n==_mesh.n_nodes()-1) ? hole_offset : hole_offset+n+1; // wrap around
201  }
202  }
203 
204  ctr +=2;
205  }
206  }
207 
208 
209  // If the user provided it, use his ordering to define the segments
210  for (std::size_t ctr=0, s=0; s<this->segments.size(); ctr+=2, ++s)
211  {
212  const unsigned int index0 = 2*hole_offset+ctr;
213  const unsigned int index1 = 2*hole_offset+ctr+1;
214 
215  initial.segmentlist[index0] = hole_offset + this->segments[s].first;
216  initial.segmentlist[index1] = hole_offset + this->segments[s].second;
217  }
218 
219 
220 
221  // Tell the input struct about the holes
222  if (have_holes)
223  {
224  initial.numberofholes = _holes->size();
225  initial.holelist = static_cast<REAL*>(std::malloc(initial.numberofholes * 2 * sizeof(REAL)));
226  for (std::size_t i=0, ctr=0; i<_holes->size(); ++i, ctr+=2)
227  {
228  Point inside_point = (*_holes)[i]->inside();
229  initial.holelist[ctr] = inside_point(0);
230  initial.holelist[ctr+1] = inside_point(1);
231  }
232  }
233 
234  // Set the triangulation flags.
235  // c ~ enclose convex hull with segments
236  // z ~ use zero indexing
237  // B ~ Suppresses boundary markers in the output
238  // Q ~ run in "quiet" mode
239  // p ~ Triangulates a Planar Straight Line Graph
240  // If the `p' switch is used, `segmentlist' must point to a list of
241  // segments, `numberofsegments' must be properly set, and
242  // `segmentmarkerlist' must either be set to nullptr (in which case all
243  // markers default to zero), or must point to a list of markers.
244  // D ~ Conforming Delaunay: use this switch if you want all triangles
245  // in the mesh to be Delaunay, and not just constrained Delaunay
246  // q ~ Quality mesh generation with no angles smaller than 20 degrees.
247  // An alternate minimum angle may be specified after the q
248  // a ~ Imposes a maximum triangle area constraint.
249  // -P Suppresses the output .poly file. Saves disk space, but you lose the ability to maintain
250  // constraining segments on later refinements of the mesh.
251  // Create the flag strings, depends on element type
252  std::ostringstream flags;
253 
254  // Default flags always used
255  flags << "zBPQ";
256 
257  // Flags which are specific to the type of triangulation
258  switch (_triangulation_type)
259  {
261  {
262  flags << "c";
263  break;
264  }
265 
266  case PSLG:
267  {
268  flags << "p";
269  break;
270  }
271 
273  libmesh_error_msg("ERROR: INVALID_TRIANGULATION_TYPE selected!");
274 
275  default:
276  libmesh_error_msg("Unrecognized _triangulation_type");
277  }
278 
279 
280  // Flags specific to the type of element
281  switch (_elem_type)
282  {
283  case TRI3:
284  {
285  // do nothing.
286  break;
287  }
288 
289  case TRI6:
290  {
291  flags << "o2";
292  break;
293  }
294 
295  default:
296  libmesh_error_msg("ERROR: Unrecognized triangular element type.");
297  }
298 
299 
300  // If we do have holes and the user asked to GENERATE_CONVEX_HULL,
301  // need to add the p flag so the triangulation respects those segments.
302  if ((_triangulation_type==GENERATE_CONVEX_HULL) && (have_holes))
303  flags << "p";
304 
305  // Finally, add the area constraint
306  if (_desired_area > TOLERANCE)
307  flags << "a" << std::fixed << _desired_area;
308 
309  // add minimum angle constraint
311  flags << "q" << std::fixed << _minimum_angle;
312 
313  // add user provided extra flags
314  if (_extra_flags.size() > 0)
315  flags << _extra_flags;
316 
317  // Refine the initial output to conform to the area constraint
318  TriangleWrapper::triangulate(const_cast<char *>(flags.str().c_str()),
319  &initial,
320  &final,
321  nullptr); // voronoi ouput -- not used
322 
323 
324  // Send the information computed by Triangle to the Mesh.
326  _mesh,
327  _elem_type);
328 
329  // To the naked eye, a few smoothing iterations usually looks better,
330  // so we do this by default unless the user says not to.
331  if (this->_smooth_after_generating)
332  LaplaceMeshSmoother(_mesh).smooth(2);
333 
334 
335  // Clean up.
338 
339  // Prepare the mesh for use before returning. This ensures (among
340  // other things) that it is partitioned and therefore users can
341  // iterate over local elements, etc.
343 }
static const Real TOLERANCE
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
const std::vector< Hole * > * _holes
void copy_tri_to_mesh(const triangulateio &triangle_data_input, UnstructuredMesh &mesh_output, const ElemType type)
virtual SimpleRange< node_iterator > node_ptr_range()=0
void init(triangulateio &t)
void prepare_for_use(const bool skip_renumber_nodes_and_elements=false, const bool skip_find_neighbors=false)
Definition: mesh_base.C:152
void set_mesh_dimension(unsigned char d)
Definition: mesh_base.h:213
std::vector< std::pair< unsigned int, unsigned int > > segments
virtual void clear()
Definition: mesh_base.C:260
void destroy(triangulateio &t, IO_Type)
virtual dof_id_type n_nodes() const =0
uint8_t dof_id_type
Definition: id_types.h:64

◆ triangulation_type()

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

Sets and/or gets the desired triangulation type.

Definition at line 145 of file mesh_triangle_interface.h.

References _triangulation_type.

Referenced by libMesh::MeshTools::Generation::build_delaunay_square().

145 {return _triangulation_type;}

Member Data Documentation

◆ _desired_area

Real libMesh::TriangleInterface::_desired_area
private

The desired area for the elements in the resulting mesh.

Definition at line 198 of file mesh_triangle_interface.h.

Referenced by desired_area(), and triangulate().

◆ _elem_type

ElemType libMesh::TriangleInterface::_elem_type
private

The type of elements to generate. (Defaults to TRI3).

Definition at line 193 of file mesh_triangle_interface.h.

Referenced by elem_type(), and triangulate().

◆ _extra_flags

std::string libMesh::TriangleInterface::_extra_flags
private

Additional flags to be passed to triangle

Definition at line 208 of file mesh_triangle_interface.h.

Referenced by extra_flags(), and triangulate().

◆ _holes

const std::vector<Hole*>* libMesh::TriangleInterface::_holes
private

A pointer to a vector of Hole*s. If this is nullptr, there are no holes!

Definition at line 187 of file mesh_triangle_interface.h.

Referenced by attach_hole_list(), and triangulate().

◆ _insert_extra_points

bool libMesh::TriangleInterface::_insert_extra_points
private

Flag which tells whether or not to insert additional nodes before triangulation. This can sometimes be used to "de-regularize" the resulting triangulation.

Definition at line 222 of file mesh_triangle_interface.h.

Referenced by insert_extra_points(), and triangulate().

◆ _mesh

UnstructuredMesh& libMesh::TriangleInterface::_mesh
private

Reference to the mesh which is to be created by triangle.

Definition at line 181 of file mesh_triangle_interface.h.

Referenced by triangulate().

◆ _minimum_angle

Real libMesh::TriangleInterface::_minimum_angle
private

Minimum angle in triangles

Definition at line 203 of file mesh_triangle_interface.h.

Referenced by minimum_angle(), and triangulate().

◆ _serializer

MeshSerializer libMesh::TriangleInterface::_serializer
private

Triangle only operates on serial meshes.

Definition at line 233 of file mesh_triangle_interface.h.

◆ _smooth_after_generating

bool libMesh::TriangleInterface::_smooth_after_generating
private

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

Definition at line 228 of file mesh_triangle_interface.h.

Referenced by smooth_after_generating(), and triangulate().

◆ _triangulation_type

TriangulationType libMesh::TriangleInterface::_triangulation_type
private

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

Definition at line 215 of file mesh_triangle_interface.h.

Referenced by triangulate(), and triangulation_type().

◆ segments

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

When constructing a PSLG, it is often not possible to do so 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 this case you could actually use the implicit ordering!

Definition at line 175 of file mesh_triangle_interface.h.

Referenced by triangulate().


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