libMesh
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
libMesh::TetGenMeshInterface Class Reference

Class TetGenMeshInterface provides an interface for tetrahedralization of meshes using the TetGen library. More...

#include <mesh_tetgen_interface.h>

Public Member Functions

 TetGenMeshInterface (UnstructuredMesh &mesh)
 Constructor. More...
 
 ~TetGenMeshInterface ()
 Empty destructor. More...
 
void set_switches (const std::string &)
 Method to set switches to tetgen, allowing for different behaviours. More...
 
void triangulate_pointset ()
 Method invokes TetGen library to compute a Delaunay tetrahedralization from the nodes point set. More...
 
void pointset_convexhull ()
 Method invokes TetGen library to compute a Delaunay tetrahedralization from the nodes point set. More...
 
void triangulate_conformingDelaunayMesh (double quality_constraint=0., double volume_constraint=0.)
 Method invokes TetGen library to compute a Delaunay tetrahedralization from the nodes point set. More...
 
void triangulate_conformingDelaunayMesh_carvehole (const std::vector< Point > &holes, double quality_constraint=0., double volume_constraint=0.)
 Method invokes TetGen library to compute a Delaunay tetrahedralization from the nodes point set. More...
 

Protected Member Functions

void fill_pointlist (TetGenWrapper &wrapper)
 This function copies nodes from the _mesh into TetGen's pointlist. More...
 
void assign_nodes_to_elem (unsigned *node_labels, Elem *elem)
 Assigns the node IDs contained in the 'node_labels' array to 'elem'. More...
 
unsigned check_hull_integrity ()
 This function checks the integrity of the current set of elements in the Mesh to see if they comprise a convex hull, that is: More...
 
void process_hull_integrity_result (unsigned result)
 This function prints an informative message and crashes based on the output of the check_hull_integrity() function. More...
 
void delete_2D_hull_elements ()
 Delete original convex hull elements from the Mesh after performing a Delaunay tetrahedralization. More...
 

Protected Attributes

UnstructuredMesh_mesh
 Local reference to the mesh we are working with. More...
 
std::vector< unsigned > _sequential_to_libmesh_node_map
 We should not assume libmesh nodes are numbered sequentially... More...
 
MeshSerializer _serializer
 Tetgen only operates on serial meshes. More...
 
std::string _switches
 Parameter controlling the behaviour of tetgen. More...
 

Detailed Description

Class TetGenMeshInterface provides an interface for tetrahedralization of meshes using the TetGen library.

For information about TetGen cf. TetGen home page.

Author
Steffen Petersen
Date
2004
Author
John W. Peterson
Date
2011

Definition at line 54 of file mesh_tetgen_interface.h.

Constructor & Destructor Documentation

◆ TetGenMeshInterface()

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

Constructor.

Takes a reference to the mesh.

Definition at line 38 of file mesh_tetgen_interface.C.

38  :
39  _mesh(mesh),
41  _switches("Q")
42 {
43 }

◆ ~TetGenMeshInterface()

libMesh::TetGenMeshInterface::~TetGenMeshInterface ( )
inline

Empty destructor.

Definition at line 67 of file mesh_tetgen_interface.h.

67 {}

Member Function Documentation

◆ assign_nodes_to_elem()

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 377 of file mesh_tetgen_interface.C.

378 {
379  for (auto j : elem->node_index_range())
380  {
381  // Get the mapped node index to ask the Mesh for
382  unsigned mapped_node_id = _sequential_to_libmesh_node_map[ node_labels[j] ];
383 
384  Node * current_node = this->_mesh.node_ptr( mapped_node_id );
385 
386  elem->set_node(j) = current_node;
387  }
388 }

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

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

◆ check_hull_integrity()

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-nullptr 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 nullptr-neighbor is found

Definition at line 394 of file mesh_tetgen_interface.C.

395 {
396  // Check for easy return: if the Mesh is empty (i.e. if
397  // somebody called triangulate_conformingDelaunayMesh on
398  // a Mesh with no elements, then hull integrity check must
399  // fail...
400  if (_mesh.n_elem() == 0)
401  return 3;
402 
403  for (auto & elem : this->_mesh.element_ptr_range())
404  {
405  // Check for proper element type
406  if (elem->type() != TRI3)
407  {
408  //libmesh_error_msg("ERROR: Some of the elements in the original mesh were not TRI3!");
409  return 1;
410  }
411 
412  for (auto neigh : elem->neighbor_ptr_range())
413  {
414  if (neigh == nullptr)
415  {
416  // libmesh_error_msg("ERROR: Non-convex hull, cannot be tetrahedralized.");
417  return 2;
418  }
419  }
420  }
421 
422  // If we made it here, return success!
423  return 0;
424 }

References _mesh, libMesh::MeshBase::element_ptr_range(), libMesh::MeshBase::n_elem(), and libMesh::TRI3.

Referenced by triangulate_conformingDelaunayMesh_carvehole().

◆ delete_2D_hull_elements()

void libMesh::TetGenMeshInterface::delete_2D_hull_elements ( )
protected

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

Definition at line 449 of file mesh_tetgen_interface.C.

450 {
451  for (auto & elem : this->_mesh.element_ptr_range())
452  {
453  // Check for proper element type. Yes, we legally delete elements while
454  // iterating over them because no entries from the underlying container
455  // are actually erased.
456  if (elem->type() == TRI3)
457  _mesh.delete_elem(elem);
458  }
459 
460  // We just removed any boundary info associated with hull element
461  // edges, so let's update the boundary id caches.
463 }

References _mesh, libMesh::MeshBase::delete_elem(), libMesh::MeshBase::element_ptr_range(), libMesh::MeshBase::get_boundary_info(), libMesh::BoundaryInfo::regenerate_id_sets(), and libMesh::TRI3.

Referenced by triangulate_conformingDelaunayMesh_carvehole().

◆ fill_pointlist()

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. DistributedMesh) are handled.

Definition at line 349 of file mesh_tetgen_interface.C.

350 {
351  // fill input structure with point set data:
352  wrapper.allocate_pointlist( this->_mesh.n_nodes() );
353 
354  // Make enough space to store a mapping between the implied sequential
355  // node numbering used in tetgen and libmesh's (possibly) non-sequential
356  // numbering scheme.
359 
360  {
361  unsigned index = 0;
362  for (auto & node : this->_mesh.node_ptr_range())
363  {
364  _sequential_to_libmesh_node_map[index] = node->id();
365  wrapper.set_node(index++,
366  REAL((*node)(0)),
367  REAL((*node)(1)),
368  REAL((*node)(2)));
369  }
370  }
371 }

References _mesh, _sequential_to_libmesh_node_map, libMesh::TetGenWrapper::allocate_pointlist(), libMesh::MeshBase::n_nodes(), libMesh::MeshBase::node_ptr_range(), and libMesh::TetGenWrapper::set_node().

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

◆ pointset_convexhull()

void libMesh::TetGenMeshInterface::pointset_convexhull ( )

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

Stores only 2D hull surface elements.

Definition at line 113 of file mesh_tetgen_interface.C.

114 {
115  // class tetgen_wrapper allows library access on a basic level
116  TetGenWrapper tetgen_wrapper;
117 
118  // Copy Mesh's node points into TetGen data structure
119  this->fill_pointlist(tetgen_wrapper);
120 
121  // Run TetGen triangulation method:
122  // Q = quiet, no terminal output
123  // Note: if no switch is used, the input must be a list of 3D points
124  // (.node file) and the Delaunay tetrahedralization of this point set
125  // will be generated. In this particular function, we are throwing
126  // away the tetrahedra generated by TetGen, and keeping only the
127  // convex hull...
128  tetgen_wrapper.set_switches(_switches);
129  tetgen_wrapper.run_tetgen();
130  unsigned int num_elements = tetgen_wrapper.get_numberoftrifaces();
131 
132  // Delete *all* old elements. Yes, we legally delete elements while
133  // iterating over them because no entries from the underlying container
134  // are actually erased.
135  for (auto & elem : this->_mesh.element_ptr_range())
136  this->_mesh.delete_elem (elem);
137 
138  // We just removed any boundary info associated with element faces
139  // or edges, so let's update the boundary id caches.
141 
142  // Add the 2D elements which comprise the convex hull back to the mesh.
143  // Vector that temporarily holds the node labels defining element.
144  unsigned int node_labels[3];
145 
146  for (unsigned int i=0; i<num_elements; ++i)
147  {
148  Elem * elem = new Tri3;
149 
150  // Get node labels associated with this element
151  for (auto j : elem->node_index_range())
152  node_labels[j] = tetgen_wrapper.get_triface_node(i,j);
153 
154  this->assign_nodes_to_elem(node_labels, elem);
155 
156  // Finally, add this element to the mesh.
157  this->_mesh.add_elem(elem);
158  }
159 }

References _mesh, _switches, libMesh::MeshBase::add_elem(), assign_nodes_to_elem(), libMesh::MeshBase::delete_elem(), libMesh::MeshBase::element_ptr_range(), fill_pointlist(), libMesh::MeshBase::get_boundary_info(), libMesh::TetGenWrapper::get_numberoftrifaces(), libMesh::TetGenWrapper::get_triface_node(), libMesh::Elem::node_index_range(), libMesh::BoundaryInfo::regenerate_id_sets(), libMesh::TetGenWrapper::run_tetgen(), and libMesh::TetGenWrapper::set_switches().

Referenced by add_cube_convex_hull_to_mesh().

◆ process_hull_integrity_result()

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 430 of file mesh_tetgen_interface.C.

431 {
432  if (result != 0)
433  {
434  libMesh::err << "Error! Conforming Delaunay mesh tetrahedralization requires a convex hull." << std::endl;
435 
436  if (result==1)
437  {
438  libMesh::err << "Non-TRI3 elements were found in the input Mesh. ";
439  libMesh::err << "A constrained Delaunay triangulation requires a convex hull of TRI3 elements." << std::endl;
440  }
441 
442  libmesh_error_msg("Consider calling TetGenMeshInterface::pointset_convexhull() followed by Mesh::find_neighbors() first.");
443  }
444 }

References libMesh::err.

Referenced by triangulate_conformingDelaunayMesh_carvehole().

◆ set_switches()

void libMesh::TetGenMeshInterface::set_switches ( const std::string &  switches)

Method to set switches to tetgen, allowing for different behaviours.

Definition at line 45 of file mesh_tetgen_interface.C.

46 {
47  // set the tetgen switch manually:
48  // p = tetrahedralizes a piecewise linear complex (see definition in user manual)
49  // Q = quiet, no terminal output
50  // q = specify a minimum radius/edge ratio
51  // a = tetrahedron volume constraint
52  // V = verbose output
53  // for full list of options and their meaning: see the tetgen manual
54  // (http://wias-berlin.de/software/tetgen/1.5/doc/manual/manual005.html)
55  _switches = switches;
56 }

References _switches.

Referenced by tetrahedralize_domain().

◆ triangulate_conformingDelaunayMesh()

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

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

Boundary constraints are taken from elements array.

Definition at line 165 of file mesh_tetgen_interface.C.

167 {
168  // start triangulation method with empty holes list:
169  std::vector<Point> noholes;
170  triangulate_conformingDelaunayMesh_carvehole(noholes, quality_constraint, volume_constraint);
171 }

References triangulate_conformingDelaunayMesh_carvehole().

◆ triangulate_conformingDelaunayMesh_carvehole()

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 tetrahedralization from the nodes point set.

Boundary constraints are taken from elements array. Include carve-out functionality.

Definition at line 175 of file mesh_tetgen_interface.C.

178 {
179  // Before calling this function, the Mesh must contain a convex hull
180  // of TRI3 elements which define the boundary.
181  unsigned hull_integrity_check = check_hull_integrity();
182 
183  // Possibly die if hull integrity check failed
184  this->process_hull_integrity_result(hull_integrity_check);
185 
186  // class tetgen_wrapper allows library access on a basic level
187  TetGenWrapper tetgen_wrapper;
188 
189  // Copy Mesh's node points into TetGen data structure
190  this->fill_pointlist(tetgen_wrapper);
191 
192  // >>> fill input structure "tetgenio" with facet data:
193  int facet_num = this->_mesh.n_elem();
194 
195  // allocate memory in "tetgenio" structure:
196  tetgen_wrapper.allocate_facetlist
197  (facet_num, cast_int<int>(holes.size()));
198 
199 
200  // Set up tetgen data structures with existing facet information
201  // from the convex hull.
202  {
203  int insertnum = 0;
204  for (auto & elem : this->_mesh.element_ptr_range())
205  {
206  tetgen_wrapper.allocate_facet_polygonlist(insertnum, 1);
207  tetgen_wrapper.allocate_polygon_vertexlist(insertnum, 0, 3);
208 
209  for (auto j : elem->node_index_range())
210  {
211  // We need to get the sequential index of elem->node_ptr(j), but
212  // it should already be stored in _sequential_to_libmesh_node_map...
213  unsigned libmesh_node_id = elem->node_id(j);
214 
215  // The libmesh node IDs may not be sequential, but can we assume
216  // they are at least in order??? We will do so here.
217  std::vector<unsigned>::iterator node_iter =
220  libmesh_node_id);
221 
222  // Check to see if not found: this could also indicate the sequential
223  // node map is not sorted...
224  if (node_iter == _sequential_to_libmesh_node_map.end())
225  libmesh_error_msg("Global node " << libmesh_node_id << " not found in sequential node map!");
226 
227  int sequential_index = cast_int<int>
229  node_iter));
230 
231  // Debugging:
232  // libMesh::out << "libmesh_node_id=" << libmesh_node_id
233  // << ", sequential_index=" << sequential_index
234  // << std::endl;
235 
236  tetgen_wrapper.set_vertex(insertnum, // facet number
237  0, // polygon (always 0)
238  j, // local vertex index in tetgen input
239  sequential_index);
240  }
241 
242  // Go to next facet in polygonlist
243  insertnum++;
244  }
245  }
246 
247 
248 
249  // fill hole list (if there are holes):
250  if (holes.size() > 0)
251  {
252  unsigned hole_index = 0;
253  for (Point ihole : holes)
254  tetgen_wrapper.set_hole(hole_index++,
255  REAL(ihole(0)),
256  REAL(ihole(1)),
257  REAL(ihole(2)));
258  }
259 
260 
261  // Run TetGen triangulation method
262 
263  // Assemble switches: we append the user's switches (if any) to 'p',
264  // which tetrahedralizes a piecewise linear complex (see definition
265  // in user manual)
266  std::ostringstream oss;
267  oss << "p";
268  oss << _switches;
269 
270  if (quality_constraint != 0)
271  oss << "q" << std::fixed << quality_constraint;
272 
273  if (volume_constraint != 0)
274  oss << "a" << std::fixed << volume_constraint;
275 
276  std::string params = oss.str();
277 
278  tetgen_wrapper.set_switches(params); // TetGen switches: Piecewise linear complex, Quiet mode
279  tetgen_wrapper.run_tetgen();
280 
281  // => nodes:
282  unsigned int old_nodesnum = this->_mesh.n_nodes();
283  REAL x=0., y=0., z=0.;
284  const unsigned int num_nodes = tetgen_wrapper.get_numberofpoints();
285 
286  // Debugging:
287  // libMesh::out << "Original mesh had " << old_nodesnum << " nodes." << std::endl;
288  // libMesh::out << "Reserving space for " << num_nodes << " total nodes." << std::endl;
289 
290  // Reserve space for additional nodes in the node map
291  _sequential_to_libmesh_node_map.reserve(num_nodes);
292 
293  // Add additional nodes to the Mesh.
294  // Original code had i<=num_nodes here (Note: the indexing is:
295  // foo[3*i], [3*i+1], [3*i+2]) But according to the TetGen docs, "In
296  // all cases, the first item in any array is stored starting at
297  // index [0]."
298  for (unsigned int i=old_nodesnum; i<num_nodes; i++)
299  {
300  // Fill in x, y, z values
301  tetgen_wrapper.get_output_node(i, x,y,z);
302 
303  // Catch the node returned by add_point()... this will tell us the ID
304  // assigned by the Mesh.
305  Node * new_node = this->_mesh.add_point ( Point(x,y,z) );
306 
307  // Store this new ID in our sequential-to-libmesh node mapping array
308  _sequential_to_libmesh_node_map.push_back( new_node->id() );
309  }
310 
311  // Debugging:
312  // std::copy(_sequential_to_libmesh_node_map.begin(),
313  // _sequential_to_libmesh_node_map.end(),
314  // std::ostream_iterator<unsigned>(std::cout, " "));
315  // std::cout << std::endl;
316 
317 
318  // => tetrahedra:
319  const unsigned int num_elements = tetgen_wrapper.get_numberoftetrahedra();
320 
321  // Vector that temporarily holds the node labels defining element connectivity.
322  unsigned int node_labels[4];
323 
324  for (unsigned int i=0; i<num_elements; i++)
325  {
326  // TetGen only supports Tet4 elements.
327  Elem * elem = new Tet4;
328 
329  // Fill up the the node_labels vector
330  for (auto j : elem->node_index_range())
331  node_labels[j] = tetgen_wrapper.get_element_node(i,j);
332 
333  // Associate nodes with this element
334  this->assign_nodes_to_elem(node_labels, elem);
335 
336  // Finally, add this element to the mesh
337  this->_mesh.add_elem(elem);
338  }
339 
340  // Delete original convex hull elements. Is there ever a case where
341  // we should not do this?
342  this->delete_2D_hull_elements();
343 }

References _mesh, _sequential_to_libmesh_node_map, _switches, 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(), distance(), libMesh::MeshBase::element_ptr_range(), 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::node_index_range(), process_hull_integrity_result(), libMesh::TetGenWrapper::run_tetgen(), libMesh::TetGenWrapper::set_hole(), libMesh::TetGenWrapper::set_switches(), and libMesh::TetGenWrapper::set_vertex().

Referenced by tetrahedralize_domain(), and triangulate_conformingDelaunayMesh().

◆ triangulate_pointset()

void libMesh::TetGenMeshInterface::triangulate_pointset ( )

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

Definition at line 59 of file mesh_tetgen_interface.C.

60 {
61  // class tetgen_wrapper allows library access on a basic level:
62  TetGenWrapper tetgen_wrapper;
63 
64  // fill input structure with point set data:
65  this->fill_pointlist(tetgen_wrapper);
66 
67  // Run TetGen triangulation method:
68  // Q = quiet, no terminal output
69  // V = verbose, more terminal output
70  // Note: if no switch is used, the input must be a list of 3D points
71  // (.node file) and the Delaunay tetrahedralization of this point set
72  // will be generated.
73 
74  // Can we apply quality and volume constraints in
75  // triangulate_pointset()?. On at least one test problem,
76  // specifying any quality or volume constraints here causes tetgen
77  // to segfault down in the insphere method: a nullptr is passed
78  // to the routine.
79  std::ostringstream oss;
80  oss << _switches;
81  // oss << "V"; // verbose operation
82  //oss << "q" << std::fixed << 2.0; // quality constraint
83  //oss << "a" << std::fixed << 100.; // volume constraint
84  tetgen_wrapper.set_switches(oss.str());
85 
86  // Run tetgen
87  tetgen_wrapper.run_tetgen();
88 
89  // save elements to mesh structure, nodes will not be changed:
90  const unsigned int num_elements = tetgen_wrapper.get_numberoftetrahedra();
91 
92  // Vector that temporarily holds the node labels defining element.
93  unsigned int node_labels[4];
94 
95  for (unsigned int i=0; i<num_elements; ++i)
96  {
97  Elem * elem = new Tet4;
98 
99  // Get the nodes associated with this element
100  for (auto j : elem->node_index_range())
101  node_labels[j] = tetgen_wrapper.get_element_node(i,j);
102 
103  // Associate the nodes with this element
104  this->assign_nodes_to_elem(node_labels, elem);
105 
106  // Finally, add this element to the mesh.
107  this->_mesh.add_elem(elem);
108  }
109 }

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

Member Data Documentation

◆ _mesh

UnstructuredMesh& libMesh::TetGenMeshInterface::_mesh
protected

◆ _sequential_to_libmesh_node_map

std::vector<unsigned> libMesh::TetGenMeshInterface::_sequential_to_libmesh_node_map
protected

We should not assume libmesh nodes are numbered sequentially...

This is not the default behavior of DistributedMesh, 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 160 of file mesh_tetgen_interface.h.

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

◆ _serializer

MeshSerializer libMesh::TetGenMeshInterface::_serializer
protected

Tetgen only operates on serial meshes.

Definition at line 165 of file mesh_tetgen_interface.h.

◆ _switches

std::string libMesh::TetGenMeshInterface::_switches
protected

Parameter controlling the behaviour of tetgen.

By default quiet.

Definition at line 171 of file mesh_tetgen_interface.h.

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


The documentation for this class was generated from the following files:
libMesh::MeshBase::get_boundary_info
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:132
libMesh::MeshBase::delete_elem
virtual void delete_elem(Elem *e)=0
Removes element e from the mesh.
libMesh::MeshBase::n_elem
virtual dof_id_type n_elem() const =0
libMesh::TetGenMeshInterface::fill_pointlist
void fill_pointlist(TetGenWrapper &wrapper)
This function copies nodes from the _mesh into TetGen's pointlist.
Definition: mesh_tetgen_interface.C:349
libMesh::TetGenMeshInterface::check_hull_integrity
unsigned check_hull_integrity()
This function checks the integrity of the current set of elements in the Mesh to see if they comprise...
Definition: mesh_tetgen_interface.C:394
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::node_ptr
virtual const Node * node_ptr(const dof_id_type i) const =0
libMesh::TetGenMeshInterface::triangulate_conformingDelaunayMesh_carvehole
void triangulate_conformingDelaunayMesh_carvehole(const std::vector< Point > &holes, double quality_constraint=0., double volume_constraint=0.)
Method invokes TetGen library to compute a Delaunay tetrahedralization from the nodes point set.
Definition: mesh_tetgen_interface.C:175
libMesh::TetGenMeshInterface::_mesh
UnstructuredMesh & _mesh
Local reference to the mesh we are working with.
Definition: mesh_tetgen_interface.h:151
libMesh::TetGenMeshInterface::_sequential_to_libmesh_node_map
std::vector< unsigned > _sequential_to_libmesh_node_map
We should not assume libmesh nodes are numbered sequentially...
Definition: mesh_tetgen_interface.h:160
libMesh::MeshBase::element_ptr_range
virtual SimpleRange< element_iterator > element_ptr_range()=0
libMesh::TetGenMeshInterface::process_hull_integrity_result
void process_hull_integrity_result(unsigned result)
This function prints an informative message and crashes based on the output of the check_hull_integri...
Definition: mesh_tetgen_interface.C:430
libMesh::MeshBase::node_ptr_range
virtual SimpleRange< node_iterator > node_ptr_range()=0
libMesh::TetGenMeshInterface::_switches
std::string _switches
Parameter controlling the behaviour of tetgen.
Definition: mesh_tetgen_interface.h:171
libMesh::REAL
Real REAL
Definition: mesh_triangle_wrapper.h:44
libMesh::TetGenMeshInterface::_serializer
MeshSerializer _serializer
Tetgen only operates on serial meshes.
Definition: mesh_tetgen_interface.h:165
libMesh::TetGenMeshInterface::delete_2D_hull_elements
void delete_2D_hull_elements()
Delete original convex hull elements from the Mesh after performing a Delaunay tetrahedralization.
Definition: mesh_tetgen_interface.C:449
libMesh::TRI3
Definition: enum_elem_type.h:39
libMesh::Utility::binary_find
ForwardIterator binary_find(ForwardIterator first, ForwardIterator last, const T &value)
The STL provides std::binary_search() which returns true or false depending on whether the searched-f...
Definition: utility.h:183
libMesh::BoundaryInfo::regenerate_id_sets
void regenerate_id_sets()
Clears and regenerates the cached sets of ids.
Definition: boundary_info.C:159
libMesh::TetGenMeshInterface::assign_nodes_to_elem
void assign_nodes_to_elem(unsigned *node_labels, Elem *elem)
Assigns the node IDs contained in the 'node_labels' array to 'elem'.
Definition: mesh_tetgen_interface.C:377
libMesh::MeshBase::n_nodes
virtual dof_id_type n_nodes() const =0
distance
Real distance(const Point &p)
Definition: subdomains_ex3.C:50
libMesh::MeshBase::add_elem
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
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::err
OStreamProxy err