$extrastylesheet
#include <mesh_triangle_interface.h>
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 () |
| Real & | desired_area () |
| Real & | minimum_angle () |
| std::string & | extra_flags () |
| TriangulationType & | triangulation_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 |
A C++ interface between LibMesh and the Triangle library written by J.R. Shewchuk.
Definition at line 49 of file mesh_triangle_interface.h.
The TriangulationType is used with the general triangulate function defind below.
| 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.
{
GENERATE_CONVEX_HULL = 0,
PSLG = 1,
INVALID_TRIANGULATION_TYPE
};
| 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.
: _mesh(mesh), _holes(NULL), _elem_type(TRI3), _desired_area(0.1), _minimum_angle(20.0), _extra_flags(""), _triangulation_type(GENERATE_CONVEX_HULL), _insert_extra_points(false), _smooth_after_generating(true), _serializer(_mesh) {}
| libMesh::TriangleInterface::~TriangleInterface | ( | ) | [inline] |
| 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 154 of file mesh_triangle_interface.h.
References _holes.
Referenced by libMesh::MeshTools::Generation::build_delaunay_square().
{_holes = 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 121 of file mesh_triangle_interface.h.
References _desired_area.
Referenced by libMesh::MeshTools::Generation::build_delaunay_square().
{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().
{return _elem_type;}
| std::string& libMesh::TriangleInterface::extra_flags | ( | ) | [inline] |
Sets and/or gets additional flags to be passed to triangle
Definition at line 132 of file mesh_triangle_interface.h.
References _extra_flags.
{return _extra_flags;}
| bool& libMesh::TriangleInterface::insert_extra_points | ( | ) | [inline] |
Sets and/or gets the flag for inserting add'l points.
Definition at line 142 of file mesh_triangle_interface.h.
References _insert_extra_points.
{return _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 127 of file mesh_triangle_interface.h.
References _minimum_angle.
{return _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 148 of file mesh_triangle_interface.h.
References _smooth_after_generating.
{return _smooth_after_generating;}
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(), end, GENERATE_CONVEX_HULL, libMesh::TriangleWrapper::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(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::TOLERANCE, libMesh::TRI3, and libMesh::TRI6.
Referenced by libMesh::MeshTools::Generation::build_delaunay_square().
{
// Will the triangulation have holes?
const bool have_holes = ((_holes != NULL) && (!_holes->empty()));
// If the initial PSLG is really simple, e.g. an L-shaped domain or
// a square/rectangle, the resulting triangulation may be very
// "structured" looking. Sometimes this is a problem if your
// intention is to work with an "unstructured" looking grid. We can
// attempt to work around this limitation by inserting midpoints
// into the original PSLG. Inserting additional points into a
// set of points meant to be a convex hull usually makes less sense.
// May or may not need to insert new points ...
if ((_triangulation_type==PSLG) && (_insert_extra_points))
{
// Make a copy of the original points from the Mesh
std::vector<Point> original_points (_mesh.n_nodes());
MeshBase::node_iterator node_it = _mesh.nodes_begin();
const MeshBase::node_iterator node_end = _mesh.nodes_end();
for (unsigned int ctr=0; node_it != node_end; ++node_it)
original_points[ctr++] = **node_it;
// Clear out the mesh
_mesh.clear();
// Make sure the new Mesh will be 2D
_mesh.set_mesh_dimension(2);
// Insert a new point on each PSLG at some random location
// np=index into new points vector
// n =index into original points vector
for (unsigned int np=0, n=0; np<2*original_points.size(); ++np)
{
// the even entries are the original points
if (np%2==0)
_mesh.add_point(original_points[n++]);
else // the odd entries are the midpoints of the original PSLG segments
_mesh.add_point ((original_points[n] + original_points[n-1])/2);
}
}
// Regardless of whether we added additional points, the set of points to
// triangulate is now sitting in the mesh.
// If the holes vector is non-NULL (and non-empty) we need to determine
// the number of additional points which the holes will add to the
// triangulation.
unsigned int n_hole_points = 0;
if (have_holes)
{
for (unsigned int i=0; i<_holes->size(); ++i)
n_hole_points += (*_holes)[i]->n_points();
}
// Triangle data structure for the mesh
TriangleWrapper::triangulateio initial;
TriangleWrapper::triangulateio final;
// Pseudo-Constructor for the triangle io structs
TriangleWrapper::init(initial);
TriangleWrapper::init(final);
initial.numberofpoints = _mesh.n_nodes() + n_hole_points;
initial.pointlist = static_cast<REAL*>(std::malloc(initial.numberofpoints * 2 * sizeof(REAL)));
if (_triangulation_type==PSLG)
{
// Implicit segment ordering: One segment per point, including hole points
if (this->segments.empty())
initial.numberofsegments = initial.numberofpoints;
// User-defined segment ordering: One segment per entry in the segments vector
else
initial.numberofsegments = this->segments.size();
}
else if (_triangulation_type==GENERATE_CONVEX_HULL)
initial.numberofsegments = n_hole_points; // One segment for each hole point
// Debugging
// libMesh::out << "Number of segments set to: " << initial.numberofsegments << std::endl;
// Allocate space for the segments (2 int per segment)
if (initial.numberofsegments > 0)
{
initial.segmentlist = static_cast<int*> (std::malloc(initial.numberofsegments * 2 * sizeof(int)));
}
// Copy all the holes' points and segments into the triangle struct.
// The hole_offset is a constant offset into the points vector which points
// past the end of the last hole point added.
unsigned int hole_offset=0;
if (have_holes)
for (unsigned int i=0; i<_holes->size(); ++i)
{
for (unsigned int ctr=0, h=0; h<(*_holes)[i]->n_points(); ctr+=2, ++h)
{
Point p = (*_holes)[i]->point(h);
const unsigned int index0 = 2*hole_offset+ctr;
const unsigned int index1 = 2*hole_offset+ctr+1;
// Save the x,y locations in the triangle struct.
initial.pointlist[index0] = p(0);
initial.pointlist[index1] = p(1);
// Set the points which define the segments
initial.segmentlist[index0] = hole_offset+h;
initial.segmentlist[index1] = (h==(*_holes)[i]->n_points()-1) ? hole_offset : hole_offset+h+1; // wrap around
}
// Update the hole_offset for the next hole
hole_offset += (*_holes)[i]->n_points();
}
// Copy all the non-hole points and segments into the triangle struct.
{
MeshBase::node_iterator it = _mesh.nodes_begin();
const MeshBase::node_iterator end = _mesh.nodes_end();
for (dof_id_type ctr=0; it != end; ctr+=2, ++it)
{
dof_id_type index = 2*hole_offset + ctr;
// Get pointer to the current node
Node* node = *it;
// Set x,y values in pointlist
initial.pointlist[index] = (*node)(0);
initial.pointlist[index+1] = (*node)(1);
// If the user requested a PSLG, the non-hole points are also segments
if (_triangulation_type==PSLG)
{
// Use implicit ordering to define segments
if (this->segments.empty())
{
dof_id_type n = ctr/2; // ctr is always even
initial.segmentlist[index] = hole_offset+n;
initial.segmentlist[index+1] = (n==_mesh.n_nodes()-1) ? hole_offset : hole_offset+n+1; // wrap around
}
}
}
}
// If the user provided it, use his ordering to define the segments
for (unsigned int ctr=0, s=0; s<this->segments.size(); ctr+=2, ++s)
{
const unsigned int index0 = 2*hole_offset+ctr;
const unsigned int index1 = 2*hole_offset+ctr+1;
initial.segmentlist[index0] = hole_offset + this->segments[s].first;
initial.segmentlist[index1] = hole_offset + this->segments[s].second;
}
// Tell the input struct about the holes
if (have_holes)
{
initial.numberofholes = _holes->size();
initial.holelist = static_cast<REAL*>(std::malloc(initial.numberofholes * 2 * sizeof(REAL)));
for (unsigned int i=0, ctr=0; i<_holes->size(); ++i, ctr+=2)
{
Point inside_point = (*_holes)[i]->inside();
initial.holelist[ctr] = inside_point(0);
initial.holelist[ctr+1] = inside_point(1);
}
}
// Set the triangulation flags.
// c ~ enclose convex hull with segments
// z ~ use zero indexing
// B ~ Suppresses boundary markers in the output
// Q ~ run in "quiet" mode
// p ~ Triangulates a Planar Straight Line Graph
// If the `p' switch is used, `segmentlist' must point to a list of
// segments, `numberofsegments' must be properly set, and
// `segmentmarkerlist' must either be set to NULL (in which case all
// markers default to zero), or must point to a list of markers.
// D ~ Conforming Delaunay: use this switch if you want all triangles
// in the mesh to be Delaunay, and not just constrained Delaunay
// q ~ Quality mesh generation with no angles smaller than 20 degrees.
// An alternate minimum angle may be specified after the q
// a ~ Imposes a maximum triangle area constraint.
// -P Suppresses the output .poly file. Saves disk space, but you lose the ability to maintain
// constraining segments on later refinements of the mesh.
// Create the flag strings, depends on element type
std::ostringstream flags;
// Default flags always used
flags << "zBPQ";
// Flags which are specific to the type of triangulation
switch (_triangulation_type)
{
case GENERATE_CONVEX_HULL:
{
flags << "c";
break;
}
case PSLG:
{
flags << "p";
break;
}
case INVALID_TRIANGULATION_TYPE:
libmesh_error_msg("ERROR: INVALID_TRIANGULATION_TYPE selected!");
default:
libmesh_error_msg("Unrecognized _triangulation_type");
}
// Flags specific to the type of element
switch (_elem_type)
{
case TRI3:
{
// do nothing.
break;
}
case TRI6:
{
flags << "o2";
break;
}
default:
libmesh_error_msg("ERROR: Unrecognized triangular element type.");
}
// If we do have holes and the user asked to GENERATE_CONVEX_HULL,
// need to add the p flag so the triangulation respects those segments.
if ((_triangulation_type==GENERATE_CONVEX_HULL) && (have_holes))
flags << "p";
// Finally, add the area constraint
if (_desired_area > TOLERANCE)
flags << "a" << std::fixed << _desired_area;
// add minimum angle constraint
if (_minimum_angle > TOLERANCE)
flags << "q" << std::fixed << _minimum_angle;
// add user provided extra flags
if (_extra_flags.size() > 0)
flags << _extra_flags;
// Refine the initial output to conform to the area constraint
TriangleWrapper::triangulate(const_cast<char*>(flags.str().c_str()),
&initial,
&final,
NULL); // voronoi ouput -- not used
// Send the information computed by Triangle to the Mesh.
TriangleWrapper::copy_tri_to_mesh(final,
_mesh,
_elem_type);
// To the naked eye, a few smoothing iterations usually looks better,
// so we do this by default unless the user says not to.
if (this->_smooth_after_generating)
LaplaceMeshSmoother(_mesh).smooth(2);
// Clean up.
TriangleWrapper::destroy(initial, TriangleWrapper::INPUT);
TriangleWrapper::destroy(final, TriangleWrapper::OUTPUT);
}
Sets and/or gets the desired triangulation type.
Definition at line 137 of file mesh_triangle_interface.h.
References _triangulation_type.
Referenced by libMesh::MeshTools::Generation::build_delaunay_square().
{return _triangulation_type;}
The desired area for the elements in the resulting mesh.
Definition at line 190 of file mesh_triangle_interface.h.
Referenced by desired_area(), and triangulate().
The type of elements to generate. (Defaults to TRI3).
Definition at line 185 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 200 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 179 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 214 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 173 of file mesh_triangle_interface.h.
Referenced by triangulate().
Minimum angle in triangles
Definition at line 195 of file mesh_triangle_interface.h.
Referenced by minimum_angle(), and triangulate().
Triangle only operates on serial meshes.
Definition at line 225 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 220 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 207 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 167 of file mesh_triangle_interface.h.
Referenced by triangulate().