libMesh
Classes | Public Types | Public Member Functions | Public Attributes | Private Attributes | List of all members
libMesh::TriangleInterface Class Reference

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

#include <mesh_triangle_interface.h>

Classes

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  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

 TriangleInterface (UnstructuredMesh &mesh)
 The constructor. More...
 
 ~TriangleInterface ()
 Empty destructor. More...
 
void triangulate ()
 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...
 
Realminimum_angle ()
 Sets and/or gets the minimum angle. More...
 
std::string & extra_flags ()
 Sets and/or gets additional flags to be passed to triangle. 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...
 
bool & smooth_after_generating ()
 Sets/gets flag which tells whether to do Delaunay mesh smoothing after generating the grid. More...
 
bool & quiet ()
 Whether or not to make Triangle quiet. More...
 
void attach_hole_list (const std::vector< Hole * > *holes)
 Attaches a vector of Hole* pointers which will be meshed around. More...
 
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...
 

Public Attributes

std::vector< std::pair< unsigned int, unsigned int > > segments
 When constructing a PSLG, it is often not possible to do so implicitly through the ordering of the points. More...
 

Private 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...
 
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...
 
std::string _extra_flags
 Additional flags to be passed to triangle. 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...
 
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 Triangle outputs. More...
 
MeshSerializer _serializer
 Triangle only operates on serial meshes. More...
 

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.

79  {
86 
97  PSLG = 1,
98 
103  };

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),
47  _markers(nullptr),
48  _regions(nullptr),
50  _desired_area(0.1),
51  _minimum_angle(20.0),
52  _extra_flags(""),
54  _insert_extra_points(false),
56  _quiet(true),
58 {}

◆ ~TriangleInterface()

libMesh::TriangleInterface::~TriangleInterface ( )
inline

Empty destructor.

Definition at line 72 of file mesh_triangle_interface.h.

72 {}

Member Function Documentation

◆ attach_boundary_marker()

void libMesh::TriangleInterface::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 194 of file mesh_triangle_interface.h.

194 { _markers = markers; }

References _markers.

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

174 {_holes = holes;}

References _holes.

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

◆ attach_region_list()

void libMesh::TriangleInterface::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 200 of file mesh_triangle_interface.h.

200 { _regions = regions; }

References _regions.

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

136 {return _desired_area;}

References _desired_area.

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

◆ elem_type()

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

Sets and/or gets the desired element type.

Definition at line 130 of file mesh_triangle_interface.h.

130 {return _elem_type;}

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

147 {return _extra_flags;}

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

157 {return _insert_extra_points;}

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

142 {return _minimum_angle;}

References _minimum_angle.

◆ quiet()

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

Whether or not to make Triangle quiet.

Definition at line 168 of file mesh_triangle_interface.h.

168 {return _quiet;}

References _quiet.

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

163 {return _smooth_after_generating;}

References _smooth_after_generating.

Referenced by triangulate_domain().

◆ triangulate()

void libMesh::TriangleInterface::triangulate ( )

This is the main public interface for this function.

Internally, it calls Triangle's triangulate routine.

Definition at line 63 of file mesh_triangle_interface.C.

64 {
65  // Will the triangulation have holes?
66  const bool have_holes = ((_holes != nullptr) && (!_holes->empty()));
67 
68  // If the initial PSLG is really simple, e.g. an L-shaped domain or
69  // a square/rectangle, the resulting triangulation may be very
70  // "structured" looking. Sometimes this is a problem if your
71  // intention is to work with an "unstructured" looking grid. We can
72  // attempt to work around this limitation by inserting midpoints
73  // into the original PSLG. Inserting additional points into a
74  // set of points meant to be a convex hull usually makes less sense.
75 
76  // May or may not need to insert new points ...
78  {
79  // Make a copy of the original points from the Mesh
80  std::vector<Point> original_points;
81  original_points.reserve (_mesh.n_nodes());
82  for (auto & node : _mesh.node_ptr_range())
83  original_points.push_back(*node);
84 
85  // Clear out the mesh
86  _mesh.clear();
87 
88  // Make sure the new Mesh will be 2D
90 
91  // Insert a new point on each PSLG at some random location
92  // np=index into new points vector
93  // n =index into original points vector
94  for (std::size_t np=0, n=0, tops=2*original_points.size(); np<tops; ++np)
95  {
96  // the even entries are the original points
97  if (np%2==0)
98  _mesh.add_point(original_points[n++]);
99 
100  else // the odd entries are the midpoints of the original PSLG segments
101  _mesh.add_point ((original_points[n] + original_points[n-1])/2);
102  }
103  }
104 
105  // Regardless of whether we added additional points, the set of points to
106  // triangulate is now sitting in the mesh.
107 
108  // If the holes vector is non-nullptr (and non-empty) we need to determine
109  // the number of additional points which the holes will add to the
110  // triangulation.
111  // Note that the number of points is always equal to the number of segments
112  // that form the holes.
113  unsigned int n_hole_points = 0;
114 
115  if (have_holes)
116  for (const auto & hole : *_holes)
117  {
118  n_hole_points += hole->n_points();
119  // A hole at least has one enclosure.
120  // Points on enclosures are ordered so that we can add segments implicitly.
121  // Elements in segment_indices() indicates the starting points of all enclosures.
122  // The last element in segment_indices() is the number of total points.
123  libmesh_assert_greater(hole->segment_indices().size(), 1);
124  libmesh_assert_equal_to(hole->segment_indices().back(), hole->n_points());
125  }
126 
127  // Triangle data structure for the mesh
128  TriangleWrapper::triangulateio initial;
129  TriangleWrapper::triangulateio final;
130  TriangleWrapper::triangulateio voronoi;
131 
132  // Pseudo-Constructor for the triangle io structs
133  TriangleWrapper::init(initial);
134  TriangleWrapper::init(final);
135  TriangleWrapper::init(voronoi);
136 
137  initial.numberofpoints = _mesh.n_nodes() + n_hole_points;
138  initial.pointlist = static_cast<REAL*>(std::malloc(initial.numberofpoints * 2 * sizeof(REAL)));
139 
141  {
142  // Implicit segment ordering: One segment per point, including hole points
143  if (this->segments.empty())
144  initial.numberofsegments = initial.numberofpoints;
145 
146  // User-defined segment ordering: One segment per entry in the segments vector
147  else
148  initial.numberofsegments = this->segments.size() + n_hole_points;
149  }
150 
152  initial.numberofsegments = n_hole_points; // One segment for each hole point
153 
154  // Allocate space for the segments (2 int per segment)
155  if (initial.numberofsegments > 0)
156  {
157  initial.segmentlist = static_cast<int *> (std::malloc(initial.numberofsegments * 2 * sizeof(int)));
158  if (_markers)
159  initial.segmentmarkerlist = static_cast<int *> (std::malloc(initial.numberofsegments * sizeof(int)));
160  }
161 
162 
163  // Copy all the holes' points and segments into the triangle struct.
164 
165  // The hole_offset is a constant offset into the points vector which points
166  // past the end of the last hole point added.
167  unsigned int hole_offset=0;
168 
169  if (have_holes)
170  for (const auto & hole : *_holes)
171  {
172  for (unsigned int ctr=0, h=0, i=0, hsism=hole->segment_indices().size()-1; i<hsism; ++i)
173  {
174  unsigned int begp = hole_offset + hole->segment_indices()[i];
175  unsigned int endp = hole->segment_indices()[i+1];
176 
177  for (; h<endp; ctr+=2, ++h)
178  {
179  Point p = hole->point(h);
180 
181  const unsigned int index0 = 2*hole_offset+ctr;
182  const unsigned int index1 = 2*hole_offset+ctr+1;
183 
184  // Save the x,y locations in the triangle struct.
185  initial.pointlist[index0] = p(0);
186  initial.pointlist[index1] = p(1);
187 
188  // Set the points which define the segments
189  initial.segmentlist[index0] = hole_offset+h;
190  initial.segmentlist[index1] = (h == endp - 1) ? begp : hole_offset + h + 1; // wrap around
191  if (_markers)
192  // 1 is reserved for boundaries of holes
193  initial.segmentmarkerlist[hole_offset+h] = 1;
194  }
195  }
196 
197  // Update the hole_offset for the next hole
198  hole_offset += hole->n_points();
199  }
200 
201 
202  // Copy all the non-hole points and segments into the triangle struct.
203  {
204  dof_id_type ctr=0;
205  for (auto & node : _mesh.node_ptr_range())
206  {
207  dof_id_type index = 2*hole_offset + ctr;
208 
209  // Set x,y values in pointlist
210  initial.pointlist[index] = (*node)(0);
211  initial.pointlist[index+1] = (*node)(1);
212 
213  // If the user requested a PSLG, the non-hole points are also segments
215  {
216  // Use implicit ordering to define segments
217  if (this->segments.empty())
218  {
219  dof_id_type n = ctr/2; // ctr is always even
220  initial.segmentlist[index] = hole_offset+n;
221  initial.segmentlist[index+1] = (n==_mesh.n_nodes()-1) ? hole_offset : hole_offset+n+1; // wrap around
222  if (_markers)
223  initial.segmentmarkerlist[hole_offset + n] = (*_markers)[n];
224  }
225  }
226 
227  ctr +=2;
228  }
229  }
230 
231 
232  // If the user provided it, use his ordering to define the segments
233  for (std::size_t ctr=0, s=0, ss=this->segments.size(); s<ss; ctr+=2, ++s)
234  {
235  const unsigned int index0 = 2*hole_offset+ctr;
236  const unsigned int index1 = 2*hole_offset+ctr+1;
237 
238  initial.segmentlist[index0] = hole_offset + this->segments[s].first;
239  initial.segmentlist[index1] = hole_offset + this->segments[s].second;
240  if (_markers)
241  initial.segmentmarkerlist[hole_offset + s] = (*_markers)[s];
242  }
243 
244 
245 
246  // Tell the input struct about the holes
247  if (have_holes)
248  {
249  initial.numberofholes = _holes->size();
250  initial.holelist = static_cast<REAL*>(std::malloc(initial.numberofholes * 2 * sizeof(REAL)));
251  for (std::size_t i=0, ctr=0, hs=_holes->size(); i<hs; ++i, ctr+=2)
252  {
253  Point inside_point = (*_holes)[i]->inside();
254  initial.holelist[ctr] = inside_point(0);
255  initial.holelist[ctr+1] = inside_point(1);
256  }
257  }
258 
259  if (_regions)
260  {
261  initial.numberofregions = _regions->size();
262  initial.regionlist = static_cast<REAL*>(std::malloc(initial.numberofregions * 4 * sizeof(REAL)));
263  for (std::size_t i=0, ctr=0, rs=_regions->size(); i<rs; ++i, ctr+=4)
264  {
265  Point inside_point = (*_regions)[i]->inside();
266  initial.regionlist[ctr] = inside_point(0);
267  initial.regionlist[ctr+1] = inside_point(1);
268  initial.regionlist[ctr+2] = (*_regions)[i]->attribute();
269  initial.regionlist[ctr+3] = (*_regions)[i]->max_area();
270  }
271  }
272 
273  // Set the triangulation flags.
274  // c ~ enclose convex hull with segments
275  // z ~ use zero indexing
276  // B ~ Suppresses boundary markers in the output
277  // Q ~ run in "quiet" mode
278  // p ~ Triangulates a Planar Straight Line Graph
279  // If the `p' switch is used, `segmentlist' must point to a list of
280  // segments, `numberofsegments' must be properly set, and
281  // `segmentmarkerlist' must either be set to nullptr (in which case all
282  // markers default to zero), or must point to a list of markers.
283  // D ~ Conforming Delaunay: use this switch if you want all triangles
284  // in the mesh to be Delaunay, and not just constrained Delaunay
285  // q ~ Quality mesh generation with no angles smaller than 20 degrees.
286  // An alternate minimum angle may be specified after the q
287  // a ~ Imposes a maximum triangle area constraint.
288  // -P Suppresses the output .poly file. Saves disk space, but you lose the ability to maintain
289  // constraining segments on later refinements of the mesh.
290  // -e Outputs (to an .edge file) a list of edges of the triangulation.
291  // -v Outputs the Voronoi diagram associated with the triangulation.
292  // Create the flag strings, depends on element type
293  std::ostringstream flags;
294 
295  // Default flags always used
296  flags << "z";
297 
298  if (_quiet)
299  flags << "QP";
300  else
301  flags << "V";
302 
303  if (_markers)
304  flags << "ev";
305 
306  // Flags which are specific to the type of triangulation
307  switch (_triangulation_type)
308  {
310  {
311  flags << "c";
312  break;
313  }
314 
315  case PSLG:
316  {
317  flags << "p";
318  break;
319  }
320 
322  libmesh_error_msg("ERROR: INVALID_TRIANGULATION_TYPE selected!");
323 
324  default:
325  libmesh_error_msg("Unrecognized _triangulation_type");
326  }
327 
328 
329  // Flags specific to the type of element
330  switch (_elem_type)
331  {
332  case TRI3:
333  {
334  // do nothing.
335  break;
336  }
337 
338  case TRI6:
339  {
340  flags << "o2";
341  break;
342  }
343 
344  default:
345  libmesh_error_msg("ERROR: Unrecognized triangular element type.");
346  }
347 
348 
349  // If we do have holes and the user asked to GENERATE_CONVEX_HULL,
350  // need to add the p flag so the triangulation respects those segments.
351  if ((_triangulation_type==GENERATE_CONVEX_HULL) && (have_holes))
352  flags << "p";
353 
354  // Finally, add the area constraint
355  if (_desired_area > TOLERANCE)
356  flags << "a" << std::fixed << _desired_area;
357 
358  // add minimum angle constraint
360  flags << "q" << std::fixed << _minimum_angle;
361 
362  if (_regions)
363  flags << "Aa";
364 
365  // add user provided extra flags
366  if (_extra_flags.size() > 0)
367  flags << _extra_flags;
368 
369  // Refine the initial output to conform to the area constraint
370  if (_markers)
371  {
372  // need Voronoi to generate boundary information
373  TriangleWrapper::triangulate(const_cast<char *>(flags.str().c_str()),
374  &initial,
375  &final,
376  &voronoi);
377 
378  // Send the information computed by Triangle to the Mesh.
380  _mesh,
381  _elem_type,
382  &voronoi);
383  }
384  else
385  {
386  TriangleWrapper::triangulate(const_cast<char *>(flags.str().c_str()),
387  &initial,
388  &final,
389  nullptr);
390  // Send the information computed by Triangle to the Mesh.
392  _mesh,
393  _elem_type);
394  }
395 
396 
397 
398  // To the naked eye, a few smoothing iterations usually looks better,
399  // so we do this by default unless the user says not to.
400  if (this->_smooth_after_generating)
401  LaplaceMeshSmoother(_mesh).smooth(2);
402 
403 
404  // Clean up.
408 
409  // Prepare the mesh for use before returning. This ensures (among
410  // other things) that it is partitioned and therefore users can
411  // iterate over local elements, etc.
413 }

References _desired_area, _elem_type, _extra_flags, _holes, _insert_extra_points, _markers, _mesh, _minimum_angle, _quiet, _regions, _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(), and triangulate_domain().

◆ triangulation_type()

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

Sets and/or gets the desired triangulation type.

Definition at line 152 of file mesh_triangle_interface.h.

152 {return _triangulation_type;}

References _triangulation_type.

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

Member Data Documentation

◆ _desired_area

Real libMesh::TriangleInterface::_desired_area
private

The desired area for the elements in the resulting mesh.

Definition at line 234 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 229 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 244 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 212 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 258 of file mesh_triangle_interface.h.

Referenced by insert_extra_points(), and triangulate().

◆ _markers

const std::vector<int>* libMesh::TriangleInterface::_markers
private

Boundary markers.

Definition at line 217 of file mesh_triangle_interface.h.

Referenced by attach_boundary_marker(), and triangulate().

◆ _mesh

UnstructuredMesh& libMesh::TriangleInterface::_mesh
private

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

Definition at line 206 of file mesh_triangle_interface.h.

Referenced by triangulate().

◆ _minimum_angle

Real libMesh::TriangleInterface::_minimum_angle
private

Minimum angle in triangles.

Definition at line 239 of file mesh_triangle_interface.h.

Referenced by minimum_angle(), and triangulate().

◆ _quiet

bool libMesh::TriangleInterface::_quiet
private

Flag which tells if we want to suppress Triangle outputs.

Definition at line 269 of file mesh_triangle_interface.h.

Referenced by quiet(), and triangulate().

◆ _regions

const std::vector<Region*>* libMesh::TriangleInterface::_regions
private

A pointer to a vector of Regions*s.

If this is nullptr, there are no regions!

Definition at line 223 of file mesh_triangle_interface.h.

Referenced by attach_region_list(), and triangulate().

◆ _serializer

MeshSerializer libMesh::TriangleInterface::_serializer
private

Triangle only operates on serial meshes.

Definition at line 274 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 264 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 251 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 187 of file mesh_triangle_interface.h.

Referenced by triangulate().


The documentation for this class was generated from the following files:
libMesh::TriangleInterface::_minimum_angle
Real _minimum_angle
Minimum angle in triangles.
Definition: mesh_triangle_interface.h:239
libMesh::TriangleInterface::_holes
const std::vector< Hole * > * _holes
A pointer to a vector of Hole*s.
Definition: mesh_triangle_interface.h:212
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::TriangleInterface::segments
std::vector< std::pair< unsigned int, unsigned int > > segments
When constructing a PSLG, it is often not possible to do so implicitly through the ordering of the po...
Definition: mesh_triangle_interface.h:187
libMesh::TriangleInterface::_serializer
MeshSerializer _serializer
Triangle only operates on serial meshes.
Definition: mesh_triangle_interface.h:274
libMesh::TriangleInterface::_markers
const std::vector< int > * _markers
Boundary markers.
Definition: mesh_triangle_interface.h:217
libMesh::TriangleInterface::_smooth_after_generating
bool _smooth_after_generating
Flag which tells whether we should smooth the mesh after it is generated.
Definition: mesh_triangle_interface.h:264
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::TriangleWrapper::copy_tri_to_mesh
void copy_tri_to_mesh(const triangulateio &triangle_data_input, UnstructuredMesh &mesh_output, const ElemType type, const triangulateio *voronoi=nullptr)
Copies triangulation data computed by triangle from a triangulateio object to a LibMesh mesh.
Definition: mesh_triangle_wrapper.C:103
libMesh::TriangleInterface::_elem_type
ElemType _elem_type
The type of elements to generate.
Definition: mesh_triangle_interface.h:229
libMesh::TriangleWrapper::destroy
void destroy(triangulateio &t, IO_Type)
Frees any memory which has been dynamically allocated by Triangle.
libMesh::MeshBase::node_ptr_range
virtual SimpleRange< node_iterator > node_ptr_range()=0
libMesh::REAL
Real REAL
Definition: mesh_triangle_wrapper.h:44
libMesh::TRI3
Definition: enum_elem_type.h:39
libMesh::TriangleInterface::_desired_area
Real _desired_area
The desired area for the elements in the resulting mesh.
Definition: mesh_triangle_interface.h:234
libMesh::TriangleInterface::GENERATE_CONVEX_HULL
Uses the triangle library to first generate a convex hull from the set of points passed in,...
Definition: mesh_triangle_interface.h:85
libMesh::TriangleInterface::_quiet
bool _quiet
Flag which tells if we want to suppress Triangle outputs.
Definition: mesh_triangle_interface.h:269
libMesh::TriangleInterface::INVALID_TRIANGULATION_TYPE
Does nothing, used as a default value.
Definition: mesh_triangle_interface.h:102
libMesh::MeshBase::n_nodes
virtual dof_id_type n_nodes() const =0
libMesh::TRI6
Definition: enum_elem_type.h:40
libMesh::TriangleWrapper::INPUT
Definition: mesh_triangle_wrapper.h:66
libMesh::TriangleInterface::_mesh
UnstructuredMesh & _mesh
Reference to the mesh which is to be created by triangle.
Definition: mesh_triangle_interface.h:206
libMesh::TriangleInterface::_extra_flags
std::string _extra_flags
Additional flags to be passed to triangle.
Definition: mesh_triangle_interface.h:244
libMesh::TriangleInterface::PSLG
Use the triangle library to triangulate a Planar Straight Line Graph which is defined implicitly by t...
Definition: mesh_triangle_interface.h:97
libMesh::MeshBase::add_point
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::MeshBase::prepare_for_use
void prepare_for_use(const bool skip_renumber_nodes_and_elements=false, const bool skip_find_neighbors=false)
Prepare a newly ecreated (or read) mesh for use.
Definition: mesh_base.C:318
libMesh::MeshBase::clear
virtual void clear()
Deletes all the element and node data that is currently stored.
Definition: mesh_base.C:429
libMesh::TriangleInterface::_regions
const std::vector< Region * > * _regions
A pointer to a vector of Regions*s.
Definition: mesh_triangle_interface.h:223
libMesh::MeshBase::set_mesh_dimension
void set_mesh_dimension(unsigned char d)
Resets the logical dimension of the mesh.
Definition: mesh_base.h:218
libMesh::TriangleInterface::_insert_extra_points
bool _insert_extra_points
Flag which tells whether or not to insert additional nodes before triangulation.
Definition: mesh_triangle_interface.h:258
libMesh::TriangleInterface::_triangulation_type
TriangulationType _triangulation_type
The type of triangulation to perform: choices are: convex hull PSLG.
Definition: mesh_triangle_interface.h:251
libMesh::TriangleWrapper::OUTPUT
Definition: mesh_triangle_wrapper.h:67