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 50 of file mesh_triangle_interface.h.

Member Enumeration Documentation

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 71 of file mesh_triangle_interface.h.

Constructor & Destructor Documentation

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 43 of file mesh_triangle_interface.C.

44  : _mesh(mesh),
47  _desired_area(0.1),
48  _minimum_angle(20.0),
49  _extra_flags(""),
51  _insert_extra_points(false),
54 {}
MeshBase & mesh
const class libmesh_nullptr_t libmesh_nullptr
const std::vector< Hole * > * _holes
libMesh::TriangleInterface::~TriangleInterface ( )
inline

Empty destructor.

Definition at line 65 of file mesh_triangle_interface.h.

65 {}

Member Function Documentation

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 155 of file mesh_triangle_interface.h.

References _holes.

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

155 {_holes = holes;}
const std::vector< Hole * > * _holes
Real& libMesh::TriangleInterface::desired_area ( )
inline

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

Definition at line 122 of file mesh_triangle_interface.h.

References _desired_area.

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

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

Sets and/or gets the desired element type.

Definition at line 116 of file mesh_triangle_interface.h.

References _elem_type.

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

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

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

Definition at line 133 of file mesh_triangle_interface.h.

References _extra_flags.

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

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

Definition at line 143 of file mesh_triangle_interface.h.

References _insert_extra_points.

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

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

Definition at line 128 of file mesh_triangle_interface.h.

References _minimum_angle.

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 149 of file mesh_triangle_interface.h.

References _smooth_after_generating.

void libMesh::TriangleInterface::triangulate ( )

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

Definition at line 59 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_nullptr, 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().

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

Sets and/or gets the desired triangulation type.

Definition at line 138 of file mesh_triangle_interface.h.

References _triangulation_type.

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

138 {return _triangulation_type;}

Member Data Documentation

Real libMesh::TriangleInterface::_desired_area
private

The desired area for the elements in the resulting mesh.

Definition at line 191 of file mesh_triangle_interface.h.

Referenced by desired_area(), and triangulate().

ElemType libMesh::TriangleInterface::_elem_type
private

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

Definition at line 186 of file mesh_triangle_interface.h.

Referenced by elem_type(), and triangulate().

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

Additional flags to be passed to triangle

Definition at line 201 of file mesh_triangle_interface.h.

Referenced by extra_flags(), and triangulate().

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

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

Definition at line 180 of file mesh_triangle_interface.h.

Referenced by attach_hole_list(), and triangulate().

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 215 of file mesh_triangle_interface.h.

Referenced by insert_extra_points(), and triangulate().

UnstructuredMesh& libMesh::TriangleInterface::_mesh
private

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

Definition at line 174 of file mesh_triangle_interface.h.

Referenced by triangulate().

Real libMesh::TriangleInterface::_minimum_angle
private

Minimum angle in triangles

Definition at line 196 of file mesh_triangle_interface.h.

Referenced by minimum_angle(), and triangulate().

MeshSerializer libMesh::TriangleInterface::_serializer
private

Triangle only operates on serial meshes.

Definition at line 226 of file mesh_triangle_interface.h.

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 221 of file mesh_triangle_interface.h.

Referenced by smooth_after_generating(), and triangulate().

TriangulationType libMesh::TriangleInterface::_triangulation_type
private

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

Definition at line 208 of file mesh_triangle_interface.h.

Referenced by triangulate(), and triangulation_type().

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 168 of file mesh_triangle_interface.h.

Referenced by triangulate().


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