libMesh::TriangleInterface Class Reference

#include <mesh_triangle_interface.h>

List of all members.

Classes

class  ArbitraryHole
class  Hole
class  PolygonHole

Public Types

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

Public Member Functions

 TriangleInterface (UnstructuredMesh &mesh)
 ~TriangleInterface ()
void triangulate ()
ElemType & elem_type ()
Realdesired_area ()
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
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

Definition at line 49 of file mesh_triangle_interface.h.


Member Enumeration Documentation

The TriangulationType is used with the general triangulate function defind 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 70 of file mesh_triangle_interface.h.

00071       {
00077         GENERATE_CONVEX_HULL = 0,
00078 
00089         PSLG = 1,
00090 
00094         INVALID_TRIANGULATION_TYPE
00095       };


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.

00044     : _mesh(mesh),
00045       _holes(NULL),
00046       _elem_type(TRI3),
00047       _desired_area(0.1),
00048       _triangulation_type(GENERATE_CONVEX_HULL),
00049       _insert_extra_points(false),
00050       _smooth_after_generating(true),
00051       _serializer(_mesh)
00052   {}

libMesh::TriangleInterface::~TriangleInterface (  )  [inline]

Empty destructor.

Definition at line 64 of file mesh_triangle_interface.h.

00064 {}


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

References _holes.

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

00142 {_holes = holes;}

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

Sets and/or gets the desired triangle area.

Definition at line 120 of file mesh_triangle_interface.h.

References _desired_area.

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

00120 {return _desired_area;}

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

Sets and/or gets the desired element type.

Definition at line 115 of file mesh_triangle_interface.h.

References _elem_type.

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

00115 {return _elem_type;}

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

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

Definition at line 130 of file mesh_triangle_interface.h.

References _insert_extra_points.

00130 {return _insert_extra_points;}

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

References _smooth_after_generating.

00136 {return _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 57 of file mesh_triangle_interface.C.

References _desired_area, _elem_type, _holes, _insert_extra_points, _mesh, _smooth_after_generating, _triangulation_type, libMesh::MeshBase::add_point(), libMesh::MeshBase::clear(), libMesh::TriangleWrapper::copy_tri_to_mesh(), libMesh::TriangleWrapper::destroy(), end, libMesh::err, GENERATE_CONVEX_HULL, libMesh::init(), libMesh::TriangleWrapper::INPUT, INVALID_TRIANGULATION_TYPE, libMesh::MeshBase::n_nodes(), libMesh::MeshBase::nodes_begin(), libMesh::MeshBase::nodes_end(), libMesh::TriangleWrapper::OUTPUT, PSLG, segments, libMesh::MeshBase::set_mesh_dimension(), libMeshEnums::TRI3, and libMeshEnums::TRI6.

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

00058   {
00059     // Will the triangulation have holes?
00060     const bool have_holes = ((_holes != NULL) && (!_holes->empty()));
00061 
00062     // If the initial PSLG is really simple, e.g. an L-shaped domain or
00063     // a square/rectangle, the resulting triangulation may be very
00064     // "structured" looking.  Sometimes this is a problem if your
00065     // intention is to work with an "unstructured" looking grid.  We can
00066     // attempt to work around this limitation by inserting midpoints
00067     // into the original PSLG.  Inserting additional points into a
00068     // set of points meant to be a convex hull usually makes less sense.
00069 
00070     // May or may not need to insert new points ...
00071     if ((_triangulation_type==PSLG) && (_insert_extra_points))
00072       {
00073         // Make a copy of the original points from the Mesh
00074         std::vector<Point> original_points (_mesh.n_nodes());
00075 
00076         MeshBase::node_iterator       node_it  = _mesh.nodes_begin();
00077         const MeshBase::node_iterator node_end = _mesh.nodes_end();
00078 
00079         for (unsigned int ctr=0; node_it != node_end; ++node_it)
00080           original_points[ctr++] = **node_it;
00081 
00082         // Clear out the mesh
00083         _mesh.clear();
00084 
00085         // Make sure the new Mesh will be 2D
00086         _mesh.set_mesh_dimension(2);
00087 
00088         // Insert a new point on each PSLG at some random location
00089         // np=index into new points vector
00090         // n =index into original points vector
00091         for (unsigned int np=0, n=0; np<2*original_points.size(); ++np)
00092           {
00093             // the even entries are the original points
00094             if (np%2==0)
00095               _mesh.add_point(original_points[n++]);
00096 
00097             else // the odd entries are the midpoints of the original PSLG segments
00098               _mesh.add_point ((original_points[n] + original_points[n-1])/2);
00099           }
00100       }
00101 
00102     // Regardless of whether we added additional points, the set of points to
00103     // triangulate is now sitting in the mesh.
00104 
00105     // If the holes vector is non-NULL (and non-empty) we need to determine
00106     // the number of additional points which the holes will add to the
00107     // triangulation.
00108     unsigned int n_hole_points = 0;
00109 
00110     if (have_holes)
00111       {
00112         for (unsigned int i=0; i<_holes->size(); ++i)
00113           n_hole_points += (*_holes)[i]->n_points();
00114       }
00115 
00116     // Triangle data structure for the mesh
00117     TriangleWrapper::triangulateio initial;
00118     TriangleWrapper::triangulateio final;
00119 
00120     // Pseudo-Constructor for the triangle io structs
00121     TriangleWrapper::init(initial);
00122     TriangleWrapper::init(final);
00123 
00124     initial.numberofpoints = _mesh.n_nodes() + n_hole_points;
00125     initial.pointlist      = static_cast<REAL*>(std::malloc(initial.numberofpoints * 2 * sizeof(REAL)));
00126 
00127     if (_triangulation_type==PSLG)
00128       {
00129         // Implicit segment ordering: One segment per point, including hole points
00130         if (this->segments.empty())
00131           initial.numberofsegments = initial.numberofpoints;
00132 
00133         // User-defined segment ordering: One segment per entry in the segments vector
00134         else
00135           initial.numberofsegments = this->segments.size();
00136       }
00137 
00138     else if (_triangulation_type==GENERATE_CONVEX_HULL)
00139       initial.numberofsegments = n_hole_points; // One segment for each hole point
00140 
00141     // Debugging
00142     // libMesh::out << "Number of segments set to: " << initial.numberofsegments << std::endl;
00143 
00144     // Allocate space for the segments (2 int per segment)
00145     if (initial.numberofsegments > 0)
00146       {
00147         initial.segmentlist = static_cast<int*> (std::malloc(initial.numberofsegments * 2 * sizeof(int)));
00148       }
00149 
00150 
00151     // Copy all the holes' points and segments into the triangle struct.
00152 
00153     // The hole_offset is a constant offset into the points vector which points
00154     // past the end of the last hole point added.
00155     unsigned int hole_offset=0;
00156 
00157     if (have_holes)
00158       for (unsigned int i=0; i<_holes->size(); ++i)
00159         {
00160           for (unsigned int ctr=0, h=0; h<(*_holes)[i]->n_points(); ctr+=2, ++h)
00161             {
00162               Point p = (*_holes)[i]->point(h);
00163 
00164               const unsigned int index0 = 2*hole_offset+ctr;
00165               const unsigned int index1 = 2*hole_offset+ctr+1;
00166 
00167               // Save the x,y locations in the triangle struct.
00168               initial.pointlist[index0] = p(0);
00169               initial.pointlist[index1] = p(1);
00170 
00171               // Set the points which define the segments
00172               initial.segmentlist[index0] = hole_offset+h;
00173               initial.segmentlist[index1] = (h==(*_holes)[i]->n_points()-1) ? hole_offset : hole_offset+h+1; // wrap around
00174             }
00175 
00176           // Update the hole_offset for the next hole
00177           hole_offset += (*_holes)[i]->n_points();
00178         }
00179 
00180 
00181     // Copy all the non-hole points and segments into the triangle struct.
00182     {
00183       MeshBase::node_iterator it = _mesh.nodes_begin();
00184       const MeshBase::node_iterator end = _mesh.nodes_end();
00185 
00186       for (dof_id_type ctr=0; it != end; ctr+=2, ++it)
00187         {
00188           dof_id_type index = 2*hole_offset + ctr;
00189 
00190           // Get pointer to the current node
00191           Node* node = *it;
00192 
00193           // Set x,y values in pointlist
00194           initial.pointlist[index] = (*node)(0);
00195           initial.pointlist[index+1] = (*node)(1);
00196 
00197           // If the user requested a PSLG, the non-hole points are also segments
00198           if (_triangulation_type==PSLG)
00199             {
00200               // Use implicit ordering to define segments
00201               if (this->segments.empty())
00202                 {
00203                   dof_id_type n = ctr/2; // ctr is always even
00204                   initial.segmentlist[index] = hole_offset+n;
00205                   initial.segmentlist[index+1] = (n==_mesh.n_nodes()-1) ? hole_offset : hole_offset+n+1; // wrap around
00206                 }
00207             }
00208         }
00209     }
00210 
00211 
00212     // If the user provided it, use his ordering to define the segments
00213     for (unsigned int ctr=0, s=0; s<this->segments.size(); ctr+=2, ++s)
00214       {
00215         const unsigned int index0 = 2*hole_offset+ctr;
00216         const unsigned int index1 = 2*hole_offset+ctr+1;
00217 
00218         initial.segmentlist[index0] = hole_offset + this->segments[s].first;
00219         initial.segmentlist[index1] = hole_offset + this->segments[s].second;
00220       }
00221 
00222 
00223 
00224     // Tell the input struct about the holes
00225     if (have_holes)
00226       {
00227         initial.numberofholes = _holes->size();
00228         initial.holelist      = static_cast<REAL*>(std::malloc(initial.numberofholes * 2 * sizeof(REAL)));
00229         for (unsigned int i=0, ctr=0; i<_holes->size(); ++i, ctr+=2)
00230           {
00231             Point inside_point = (*_holes)[i]->inside();
00232             initial.holelist[ctr]   = inside_point(0);
00233             initial.holelist[ctr+1] = inside_point(1);
00234           }
00235       }
00236 
00237     // Set the triangulation flags.
00238     // c ~ enclose convex hull with segments
00239     // z ~ use zero indexing
00240     // B ~ Suppresses boundary markers in the output
00241     // Q ~ run in "quiet" mode
00242     // p ~ Triangulates a Planar Straight Line Graph
00243     //     If the `p' switch is used, `segmentlist' must point to a list of
00244     //     segments, `numberofsegments' must be properly set, and
00245     //     `segmentmarkerlist' must either be set to NULL (in which case all
00246     //     markers default to zero), or must point to a list of markers.
00247     // D ~ Conforming Delaunay: use this switch if you want all triangles
00248     //     in the mesh to be Delaunay, and not just constrained Delaunay
00249     // q ~  Quality mesh generation with no angles smaller than 20 degrees.
00250     //      An alternate minimum angle may be specified after the q
00251     // a ~ Imposes a maximum triangle area constraint.
00252     // -P  Suppresses the output .poly file. Saves disk space, but you lose the ability to maintain
00253     //     constraining segments on later refinements of the mesh.
00254     // Create the flag strings, depends on element type
00255     std::ostringstream flags;
00256 
00257     // Default flags always used
00258     flags << "zBPQq";
00259 
00260     // Flags which are specific to the type of triangulation
00261     switch (_triangulation_type)
00262       {
00263       case GENERATE_CONVEX_HULL:
00264         {
00265           flags << "c";
00266           break;
00267         }
00268 
00269       case PSLG:
00270         {
00271           flags << "p";
00272           break;
00273         }
00274 
00275       case INVALID_TRIANGULATION_TYPE:
00276         {
00277           libmesh_error();
00278           break;
00279         }
00280 
00281       default:
00282         {
00283           libmesh_error();
00284         }
00285       }
00286 
00287 
00288     // Flags specific to the type of element
00289     switch (_elem_type)
00290       {
00291       case TRI3:
00292         {
00293           // do nothing.
00294           break;
00295         }
00296 
00297       case TRI6:
00298         {
00299           flags << "o2";
00300           break;
00301         }
00302 
00303       default:
00304         {
00305           libMesh::err << "ERROR: Unrecognized triangular element type." << std::endl;
00306           libmesh_error();
00307         }
00308       }
00309 
00310 
00311     // If we do have holes and the user asked to GENERATE_CONVEX_HULL,
00312     // need to add the p flag so the triangulation respects those segments.
00313     if ((_triangulation_type==GENERATE_CONVEX_HULL) && (have_holes))
00314       flags << "p";
00315 
00316     // Finally, add the area constraint
00317     flags << "a" << std::fixed << _desired_area;
00318 
00319     // Refine the initial output to conform to the area constraint
00320     TriangleWrapper::triangulate(const_cast<char*>(flags.str().c_str()),
00321                                  &initial,
00322                                  &final,
00323                                  NULL); // voronoi ouput -- not used
00324 
00325 
00326     // Send the information computed by Triangle to the Mesh.
00327     TriangleWrapper::copy_tri_to_mesh(final,
00328                                       _mesh,
00329                                       _elem_type);
00330 
00331     // To the naked eye, a few smoothing iterations usually looks better,
00332     // so we do this by default unless the user says not to.
00333     if (this->_smooth_after_generating)
00334       LaplaceMeshSmoother(_mesh).smooth(2);
00335 
00336 
00337     // Clean up.
00338     TriangleWrapper::destroy(initial,      TriangleWrapper::INPUT);
00339     TriangleWrapper::destroy(final,        TriangleWrapper::OUTPUT);
00340 
00341   }

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

Sets and/or gets the desired triangulation type.

Definition at line 125 of file mesh_triangle_interface.h.

References _triangulation_type.

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

00125 {return _triangulation_type;}


Member Data Documentation

The desired area for the elements in the resulting mesh.

Definition at line 178 of file mesh_triangle_interface.h.

Referenced by desired_area(), and triangulate().

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

Definition at line 173 of file mesh_triangle_interface.h.

Referenced by elem_type(), 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 167 of file mesh_triangle_interface.h.

Referenced by attach_hole_list(), and triangulate().

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

Referenced by insert_extra_points(), and triangulate().

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

Definition at line 161 of file mesh_triangle_interface.h.

Referenced by triangulate().

Triangle only operates on serial meshes.

Definition at line 203 of file mesh_triangle_interface.h.

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

Definition at line 198 of file mesh_triangle_interface.h.

Referenced by smooth_after_generating(), and triangulate().

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

Definition at line 185 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) Note: for this case, you could use the implicit ordering!

Definition at line 155 of file mesh_triangle_interface.h.

Referenced by triangulate().


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

Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:29 UTC

Hosted By:
SourceForge.net Logo