libMesh::TetGenMeshInterface Class Reference

#include <mesh_tetgen_interface.h>

List of all members.

Public Member Functions

 TetGenMeshInterface (UnstructuredMesh &mesh)
 ~TetGenMeshInterface ()
void triangulate_pointset ()
void pointset_convexhull ()
void triangulate_conformingDelaunayMesh (double quality_constraint=0., double volume_constraint=0.)
void triangulate_conformingDelaunayMesh_carvehole (const std::vector< Point > &holes, double quality_constraint=0., double volume_constraint=0.)

Protected Member Functions

void fill_pointlist (TetGenWrapper &wrapper)
void assign_nodes_to_elem (unsigned *node_labels, Elem *elem)
unsigned check_hull_integrity ()
void process_hull_integrity_result (unsigned result)
void delete_2D_hull_elements ()

Protected Attributes

UnstructuredMesh_mesh
std::vector< unsigned > _sequential_to_libmesh_node_map
MeshSerializer _serializer

Detailed Description

Class TetGenMeshInterface provides an interface for tetrahedrization of meshes using the TetGen library. For information about TetGen cf. TetGen home page.

Author:
, Steffen Petersen, 2004 Refactoring, John W. Peterson, 2011

Definition at line 53 of file mesh_tetgen_interface.h.


Constructor & Destructor Documentation

libMesh::TetGenMeshInterface::TetGenMeshInterface ( UnstructuredMesh mesh  )  [explicit]

Constructor. Takes a reference to the mesh.

Definition at line 38 of file mesh_tetgen_interface.C.

00038                                                                   :
00039     _mesh(mesh),
00040     _serializer(_mesh)
00041   {
00042   }

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

Empty destructor.

Definition at line 66 of file mesh_tetgen_interface.h.

00066 {};


Member Function Documentation

void libMesh::TetGenMeshInterface::assign_nodes_to_elem ( unsigned *  node_labels,
Elem elem 
) [protected]

Assigns the node IDs contained in the 'node_labels' array to 'elem'.

Definition at line 368 of file mesh_tetgen_interface.C.

References _mesh, _sequential_to_libmesh_node_map, libMesh::Elem::n_nodes(), libMesh::MeshBase::node_ptr(), and libMesh::Elem::set_node().

Referenced by pointset_convexhull(), triangulate_conformingDelaunayMesh_carvehole(), and triangulate_pointset().

00369   {
00370     for (unsigned int j=0; j<elem->n_nodes(); ++j)
00371       {
00372         // Get the mapped node index to ask the Mesh for
00373         unsigned mapped_node_id = _sequential_to_libmesh_node_map[ node_labels[j] ];
00374 
00375         // Parallel mesh can return NULL pointers, this is bad...
00376         Node* current_node = this->_mesh.node_ptr( mapped_node_id );
00377 
00378         if (current_node == NULL)
00379           {
00380             std::cerr << "Error! Mesh returned NULL node pointer!" << std::endl;
00381             libmesh_error();
00382           }
00383 
00384         elem->set_node(j) = current_node;
00385       }
00386   }

unsigned libMesh::TetGenMeshInterface::check_hull_integrity (  )  [protected]

This function checks the integrity of the current set of elements in the Mesh to see if they comprise a convex hull, that is: .) If they are all TRI3 elements .) They all have non-NULL neighbors

Returns: .) 0 if the mesh forms a valid convex hull .) 1 if a non-TRI3 element is found .) 2 if an element with a NULL-neighbor is found

Definition at line 392 of file mesh_tetgen_interface.C.

References _mesh, libMesh::MeshBase::elements_begin(), libMesh::MeshBase::elements_end(), end, libMesh::MeshBase::n_elem(), libMesh::Elem::n_neighbors(), libMesh::Elem::neighbor(), libMeshEnums::TRI3, and libMesh::Elem::type().

Referenced by triangulate_conformingDelaunayMesh_carvehole().

00393   {
00394     // Check for easy return: if the Mesh is empty (i.e. if
00395     // somebody called triangulate_conformingDelaunayMesh on
00396     // a Mesh with no elements, then hull integrity check must
00397     // fail...
00398     if (_mesh.n_elem() == 0)
00399       return 3;
00400 
00401     MeshBase::element_iterator it        = this->_mesh.elements_begin();
00402     const MeshBase::element_iterator end = this->_mesh.elements_end();
00403 
00404     for (; it != end ; ++it)
00405       {
00406         Elem* elem = *it;
00407 
00408         // Check for proper element type
00409         if (elem->type() != TRI3)
00410           {
00411             //libMesh::err << "ERROR: Some of the elements in the original mesh were not TRI3!" << std::endl;
00412             //libmesh_error();
00413             return 1;
00414           }
00415 
00416         for (unsigned int i=0; i<elem->n_neighbors(); ++i)
00417           {
00418             if (elem->neighbor(i) == NULL)
00419               {
00420                 // libMesh::err << "ERROR: Non-convex hull, cannot be tetrahedralized." << std::endl;
00421                 // libmesh_error();
00422                 return 2;
00423               }
00424           }
00425       }
00426 
00427     // If we made it here, return success!
00428     return 0;
00429   }

void libMesh::TetGenMeshInterface::delete_2D_hull_elements (  )  [protected]

Delete original convex hull elements from the Mesh after performing a Delaunay tetrahedralization.

Definition at line 457 of file mesh_tetgen_interface.C.

References _mesh, libMesh::MeshBase::delete_elem(), libMesh::MeshBase::elements_begin(), libMesh::MeshBase::elements_end(), end, libMeshEnums::TRI3, and libMesh::Elem::type().

Referenced by triangulate_conformingDelaunayMesh_carvehole().

00458   {
00459     MeshBase::element_iterator it        = this->_mesh.elements_begin();
00460     const MeshBase::element_iterator end = this->_mesh.elements_end();
00461 
00462     for (; it != end ; ++it)
00463       {
00464         Elem* elem = *it;
00465 
00466         // Check for proper element type. Yes, we legally delete elements while
00467         // iterating over them because no entries from the underlying container
00468         // are actually erased.
00469         if (elem->type() == TRI3)
00470           _mesh.delete_elem(elem);
00471       }
00472   }

void libMesh::TetGenMeshInterface::fill_pointlist ( TetGenWrapper wrapper  )  [protected]

This function copies nodes from the _mesh into TetGen's pointlist. Takes some pains to ensure that non-sequential node numberings (which can happen with e.g. ParallelMesh) are handled.

Definition at line 341 of file mesh_tetgen_interface.C.

References _mesh, _sequential_to_libmesh_node_map, libMesh::TetGenWrapper::allocate_pointlist(), end, libMesh::MeshBase::n_nodes(), libMesh::MeshBase::nodes_begin(), libMesh::MeshBase::nodes_end(), and libMesh::TetGenWrapper::set_node().

Referenced by pointset_convexhull(), triangulate_conformingDelaunayMesh_carvehole(), and triangulate_pointset().

00342   {
00343     // fill input structure with point set data:
00344     wrapper.allocate_pointlist( this->_mesh.n_nodes() );
00345 
00346     // Make enough space to store a mapping between the implied sequential
00347     // node numbering used in tetgen and libmesh's (possibly) non-sequential
00348     // numbering scheme.
00349     _sequential_to_libmesh_node_map.clear();
00350     _sequential_to_libmesh_node_map.resize( this->_mesh.n_nodes() );
00351 
00352     {
00353       unsigned index = 0;
00354       MeshBase::node_iterator it  = this->_mesh.nodes_begin();
00355       const MeshBase::node_iterator end = this->_mesh.nodes_end();
00356       for ( ; it != end; ++it)
00357         {
00358           _sequential_to_libmesh_node_map[index] = (*it)->id();
00359           wrapper.set_node(index++, (**it)(0), (**it)(1), (**it)(2));
00360         }
00361     }
00362   }

void libMesh::TetGenMeshInterface::pointset_convexhull (  ) 

Method invokes TetGen library to compute a Delaunay tetrahedrization from the nodes point set. Stores only 2D hull surface elements.

Definition at line 100 of file mesh_tetgen_interface.C.

References _mesh, libMesh::MeshBase::add_elem(), assign_nodes_to_elem(), libMesh::MeshBase::delete_elem(), libMesh::MeshBase::elements_begin(), libMesh::MeshBase::elements_end(), end, fill_pointlist(), libMesh::TetGenWrapper::get_numberoftrifaces(), libMesh::TetGenWrapper::get_triface_node(), libMesh::Elem::n_nodes(), libMesh::TetGenWrapper::run_tetgen(), and libMesh::TetGenWrapper::set_switches().

00101   {
00102     // class tetgen_wrapper allows library access on a basic level
00103     TetGenWrapper tetgen_wrapper;
00104 
00105     // Copy Mesh's node points into TetGen data structure
00106     this->fill_pointlist(tetgen_wrapper);
00107 
00108     // Run TetGen triangulation method:
00109     // Q = quiet, no terminal output
00110     // Note: if no switch is used, the input must be a list of 3D points
00111     // (.node file) and the Delaunay tetrahedralization of this point set
00112     // will be generated.  In this particular function, we are throwing
00113     // away the tetrahedra generated by TetGen, and keeping only the
00114     // convex hull...
00115     tetgen_wrapper.set_switches("Q");
00116     tetgen_wrapper.run_tetgen();
00117     unsigned int num_elements   = tetgen_wrapper.get_numberoftrifaces();
00118 
00119     // Delete *all* old elements.  Yes, we legally delete elements while
00120     // iterating over them because no entries from the underlying container
00121     // are actually erased.
00122     {
00123       MeshBase::element_iterator       it  = this->_mesh.elements_begin();
00124       const MeshBase::element_iterator end = this->_mesh.elements_end();
00125       for ( ; it != end; ++it)
00126         this->_mesh.delete_elem (*it);
00127     }
00128 
00129 
00130     // Add the 2D elements which comprise the convex hull back to the mesh.
00131     // Vector that temporarily holds the node labels defining element.
00132     unsigned int node_labels[3];
00133 
00134     for (unsigned int i=0; i<num_elements; ++i)
00135       {
00136         Elem* elem = new Tri3;
00137 
00138         // Get node labels associated with this element
00139         for (unsigned int j=0; j<elem->n_nodes(); ++j)
00140           node_labels[j] = tetgen_wrapper.get_triface_node(i,j);
00141 
00142         this->assign_nodes_to_elem(node_labels, elem);
00143 
00144         // Finally, add this element to the mesh.
00145         this->_mesh.add_elem(elem);
00146       }
00147   }

void libMesh::TetGenMeshInterface::process_hull_integrity_result ( unsigned  result  )  [protected]

This function prints an informative message and crashes based on the output of the check_hull_integrity() function. It is a separate function so that you can check hull integrity without crashing if you desire.

Definition at line 435 of file mesh_tetgen_interface.C.

References libMesh::err.

Referenced by triangulate_conformingDelaunayMesh_carvehole().

00436   {
00437     if (result != 0)
00438       {
00439         libMesh::err << "Error! Conforming Delaunay mesh tetrahedralization requires a convex hull." << std::endl;
00440 
00441         if (result==1)
00442           {
00443             libMesh::err << "Non-TRI3 elements were found in the input Mesh.  ";
00444             libMesh::err << "A constrained Delaunay triangulation requires a convex hull of TRI3 elements." << std::endl;
00445           }
00446 
00447         libMesh::err << "Consider calling TetGenMeshInterface::pointset_convexhull() followed " << std::endl;
00448         libMesh::err << "by Mesh::find_neighbors() first." << std::endl;
00449 
00450         libmesh_error();
00451       }
00452   }

void libMesh::TetGenMeshInterface::triangulate_conformingDelaunayMesh ( double  quality_constraint = 0.,
double  volume_constraint = 0. 
)

Method invokes TetGen library to compute a Delaunay tetrahedrization from the nodes point set. Boundary constraints are taken from elements array.

Definition at line 153 of file mesh_tetgen_interface.C.

References triangulate_conformingDelaunayMesh_carvehole().

00155   {
00156     // start triangulation method with empty holes list:
00157     std::vector<Point> noholes;
00158     triangulate_conformingDelaunayMesh_carvehole(noholes, quality_constraint, volume_constraint);
00159   }

void libMesh::TetGenMeshInterface::triangulate_conformingDelaunayMesh_carvehole ( const std::vector< Point > &  holes,
double  quality_constraint = 0.,
double  volume_constraint = 0. 
)

Method invokes TetGen library to compute a Delaunay tetrahedrization from the nodes point set. Boundary constraints are taken from elements array. Include carve-out functionality.

Definition at line 163 of file mesh_tetgen_interface.C.

References _mesh, _sequential_to_libmesh_node_map, libMesh::MeshBase::add_elem(), libMesh::MeshBase::add_point(), libMesh::TetGenWrapper::allocate_facet_polygonlist(), libMesh::TetGenWrapper::allocate_facetlist(), libMesh::TetGenWrapper::allocate_polygon_vertexlist(), assign_nodes_to_elem(), libMesh::Utility::binary_find(), check_hull_integrity(), delete_2D_hull_elements(), libMesh::MeshBase::elements_begin(), libMesh::MeshBase::elements_end(), end, libMesh::err, fill_pointlist(), libMesh::TetGenWrapper::get_element_node(), libMesh::TetGenWrapper::get_numberofpoints(), libMesh::TetGenWrapper::get_numberoftetrahedra(), libMesh::TetGenWrapper::get_output_node(), libMesh::DofObject::id(), libMesh::MeshBase::n_elem(), libMesh::MeshBase::n_nodes(), libMesh::Elem::n_nodes(), libMesh::Elem::node(), process_hull_integrity_result(), libMesh::TetGenWrapper::run_tetgen(), libMesh::TetGenWrapper::set_hole(), libMesh::TetGenWrapper::set_switches(), and libMesh::TetGenWrapper::set_vertex().

Referenced by triangulate_conformingDelaunayMesh().

00166   {
00167     // Before calling this function, the Mesh must contain a convex hull
00168     // of TRI3 elements which define the boundary.
00169     unsigned hull_integrity_check = check_hull_integrity();
00170 
00171     // Possibly die if hull integrity check failed
00172     this->process_hull_integrity_result(hull_integrity_check);
00173 
00174     // class tetgen_wrapper allows library access on a basic level
00175     TetGenWrapper tetgen_wrapper;
00176 
00177     // Copy Mesh's node points into TetGen data structure
00178     this->fill_pointlist(tetgen_wrapper);
00179 
00180     // >>> fill input structure "tetgenio" with facet data:
00181     int facet_num = this->_mesh.n_elem();
00182 
00183     // allocate memory in "tetgenio" structure:
00184     tetgen_wrapper.allocate_facetlist(facet_num, holes.size());
00185 
00186 
00187     // Set up tetgen data structures with existing facet information
00188     // from the convex hull.
00189     {
00190       int insertnum = 0;
00191       MeshBase::element_iterator it        = this->_mesh.elements_begin();
00192       const MeshBase::element_iterator end = this->_mesh.elements_end();
00193       for (; it != end ; ++it)
00194         {
00195           tetgen_wrapper.allocate_facet_polygonlist(insertnum, 1);
00196           tetgen_wrapper.allocate_polygon_vertexlist(insertnum, 0, 3);
00197 
00198           Elem* elem = *it;
00199 
00200           for (unsigned int j=0; j<elem->n_nodes(); ++j)
00201             {
00202               // We need to get the sequential index of elem->get_node(j), but
00203               // it should already be stored in _sequential_to_libmesh_node_map...
00204               unsigned libmesh_node_id = elem->node(j);
00205 
00206               // The libmesh node IDs may not be sequential, but can we assume
00207               // they are at least in order???  We will do so here.
00208               std::vector<unsigned>::iterator node_iter =
00209                 Utility::binary_find(_sequential_to_libmesh_node_map.begin(),
00210                                      _sequential_to_libmesh_node_map.end(),
00211                                      libmesh_node_id);
00212 
00213               // Check to see if not found: this could also indicate the sequential
00214               // node map is not sorted...
00215               if (node_iter == _sequential_to_libmesh_node_map.end())
00216                 {
00217                   libMesh::err << "Global node " << libmesh_node_id << " not found in sequential node map!"  << std::endl;
00218                   libmesh_error();
00219                 }
00220 
00221               std::vector<unsigned>::difference_type
00222                 sequential_index = std::distance(_sequential_to_libmesh_node_map.begin(), node_iter);
00223 
00224               // Debugging:
00225               //            std::cout << "libmesh_node_id=" << libmesh_node_id
00226               //                      << ", sequential_index=" << sequential_index
00227               //                      << std::endl;
00228 
00229               tetgen_wrapper.set_vertex(insertnum, // facet number
00230                                         0,         // polygon (always 0)
00231                                         j,         // local vertex index in tetgen input
00232                                         sequential_index);
00233             }
00234 
00235           // Go to next facet in polygonlist
00236           insertnum++;
00237         }
00238     }
00239 
00240 
00241 
00242     // fill hole list (if there are holes):
00243     if (holes.size() > 0)
00244       {
00245         std::vector<Point>::const_iterator ihole;
00246         unsigned hole_index = 0;
00247         for (ihole=holes.begin(); ihole!=holes.end(); ++ihole)
00248           tetgen_wrapper.set_hole(hole_index++, (*ihole)(0), (*ihole)(1), (*ihole)(2));
00249       }
00250 
00251 
00252     // Run TetGen triangulation method:
00253     // p = tetrahedralizes a piecewise linear complex (see definition in user manual)
00254     // Q = quiet, no terminal output
00255     // q = specify a minimum radius/edge ratio
00256     // a = tetrahedron volume constraint
00257 
00258     // assemble switches:
00259     std::ostringstream oss; // string holding switches
00260     oss << "pQ";
00261 
00262     if (quality_constraint != 0)
00263       oss << "q" << std::fixed << quality_constraint;
00264 
00265     if (volume_constraint != 0)
00266       oss << "a" << std::fixed << volume_constraint;
00267 
00268     std::string params = oss.str();
00269 
00270     tetgen_wrapper.set_switches(params); // TetGen switches: Piecewise linear complex, Quiet mode
00271     tetgen_wrapper.run_tetgen();
00272 
00273     // => nodes:
00274     unsigned int old_nodesnum = this->_mesh.n_nodes();
00275     REAL x=0., y=0., z=0.;
00276     const unsigned int num_nodes = tetgen_wrapper.get_numberofpoints();
00277 
00278     // Debugging:
00279     // std::cout << "Original mesh had " << old_nodesnum << " nodes." << std::endl;
00280     // std::cout << "Reserving space for " << num_nodes << " total nodes." << std::endl;
00281 
00282     // Reserve space for additional nodes in the node map
00283     _sequential_to_libmesh_node_map.reserve(num_nodes);
00284 
00285     // Add additional nodes to the Mesh.
00286     // Original code had i<=num_nodes here (Note: the indexing is:
00287     // foo[3*i], [3*i+1], [3*i+2]) But according to the TetGen docs, "In
00288     // all cases, the first item in any array is stored starting at
00289     // index [0]."
00290     for (unsigned int i=old_nodesnum; i<num_nodes; i++)
00291       {
00292         // Fill in x, y, z values
00293         tetgen_wrapper.get_output_node(i, x,y,z);
00294 
00295         // Catch the node returned by add_point()... this will tell us the ID
00296         // assigned by the Mesh.
00297         Node* new_node = this->_mesh.add_point ( Point(x,y,z) );
00298 
00299         // Store this new ID in our sequential-to-libmesh node mapping array
00300         _sequential_to_libmesh_node_map.push_back( new_node->id() );
00301       }
00302 
00303     // Debugging:
00304     //  std::copy(_sequential_to_libmesh_node_map.begin(),
00305     //      _sequential_to_libmesh_node_map.end(),
00306     //      std::ostream_iterator<unsigned>(std::cout, " "));
00307     //  std::cout << std::endl;
00308 
00309 
00310     // => tetrahedra:
00311     const unsigned int num_elements = tetgen_wrapper.get_numberoftetrahedra();
00312 
00313     // Vector that temporarily holds the node labels defining element connectivity.
00314     unsigned int node_labels[4];
00315 
00316     for (unsigned int i=0; i<num_elements; i++)
00317       {
00318         // TetGen only supports Tet4 elements.
00319         Elem* elem = new Tet4;
00320 
00321         // Fill up the the node_labels vector
00322         for (unsigned int j=0; j<elem->n_nodes(); j++)
00323           node_labels[j] = tetgen_wrapper.get_element_node(i,j);
00324 
00325         // Associate nodes with this element
00326         this->assign_nodes_to_elem(node_labels, elem);
00327 
00328         // Finally, add this element to the mesh
00329         this->_mesh.add_elem(elem);
00330       }
00331 
00332     // Delete original convex hull elements.  Is there ever a case where
00333     // we should not do this?
00334     this->delete_2D_hull_elements();
00335   }

void libMesh::TetGenMeshInterface::triangulate_pointset (  ) 

Method invokes TetGen library to compute a Delaunay tetrahedrization from the nodes point set.

Definition at line 46 of file mesh_tetgen_interface.C.

References _mesh, libMesh::MeshBase::add_elem(), assign_nodes_to_elem(), fill_pointlist(), libMesh::TetGenWrapper::get_element_node(), libMesh::TetGenWrapper::get_numberoftetrahedra(), libMesh::Elem::n_nodes(), libMesh::TetGenWrapper::run_tetgen(), and libMesh::TetGenWrapper::set_switches().

00047   {
00048     // class tetgen_wrapper allows library access on a basic level:
00049     TetGenWrapper tetgen_wrapper;
00050 
00051     // fill input structure with point set data:
00052     this->fill_pointlist(tetgen_wrapper);
00053 
00054     // Run TetGen triangulation method:
00055     // Q = quiet, no terminal output
00056     // V = verbose, more terminal output
00057     // Note: if no switch is used, the input must be a list of 3D points
00058     // (.node file) and the Delaunay tetrahedralization of this point set
00059     // will be generated.
00060 
00061     // Can we apply quality and volume constraints in
00062     // triangulate_pointset()?.  On at least one test problem,
00063     // specifying any quality or volume constraints here causes tetgen
00064     // to segfault down in the insphere method: a NULL pointer is passed
00065     // to the routine.
00066     std::ostringstream oss;
00067     oss << "Q"; // quiet operation
00068     // oss << "V"; // verbose operation
00069     //oss  << "q" << std::fixed << 2.0;  // quality constraint
00070     //oss  << "a" << std::fixed << 100.; // volume constraint
00071     tetgen_wrapper.set_switches(oss.str());
00072 
00073     // Run tetgen
00074     tetgen_wrapper.run_tetgen();
00075 
00076     // save elements to mesh structure, nodes will not be changed:
00077     const unsigned int num_elements   = tetgen_wrapper.get_numberoftetrahedra();
00078 
00079     // Vector that temporarily holds the node labels defining element.
00080     unsigned int node_labels[4];
00081 
00082     for (unsigned int i=0; i<num_elements; ++i)
00083       {
00084         Elem* elem = new Tet4;
00085 
00086         // Get the nodes associated with this element
00087         for (unsigned int j=0; j<elem->n_nodes(); ++j)
00088           node_labels[j] = tetgen_wrapper.get_element_node(i,j);
00089 
00090         // Associate the nodes with this element
00091         this->assign_nodes_to_elem(node_labels, elem);
00092 
00093         // Finally, add this element to the mesh.
00094         this->_mesh.add_elem(elem);
00095       }
00096   }


Member Data Documentation

We should not assume libmesh nodes are numbered sequentially... This is not the default behavior of ParallelMesh, for example, unless you specify node IDs explicitly. So this array allows us to keep a mapping between the sequential numbering in tetgen_data.pointlist.

Definition at line 154 of file mesh_tetgen_interface.h.

Referenced by assign_nodes_to_elem(), fill_pointlist(), and triangulate_conformingDelaunayMesh_carvehole().

Tetgen only operates on serial meshes.

Definition at line 159 of file mesh_tetgen_interface.h.


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