18 #include "libmesh/libmesh_config.h" 19 #ifdef LIBMESH_HAVE_TETGEN 26 #include "libmesh/mesh_tetgen_interface.h" 28 #include "libmesh/boundary_info.h" 29 #include "libmesh/cell_tet4.h" 30 #include "libmesh/face_tri3.h" 31 #include "libmesh/mesh_smoother_laplace.h" 32 #include "libmesh/unstructured_mesh.h" 33 #include "libmesh/utility.h" 34 #include "libmesh/mesh_tetgen_wrapper.h" 88 std::ostringstream oss;
107 unsigned int node_labels[4];
109 for (
unsigned int i=0; i<num_elements; ++i)
114 for (
auto j : elem->node_index_range())
154 for (
auto & elem : this->
_mesh.element_ptr_range())
163 unsigned int node_labels[3];
165 for (
unsigned int i=0; i<num_elements; ++i)
170 for (
auto j : elem->node_index_range())
190 double volume_constraint)
193 std::vector<Point> noholes;
200 double quality_constraint,
201 double volume_constraint)
221 (facet_num, cast_int<int>(holes.size()));
228 for (
auto & elem : this->
_mesh.element_ptr_range())
233 for (
auto j : elem->node_index_range())
237 unsigned libmesh_node_id = elem->node_id(j);
241 std::vector<unsigned>::iterator node_iter =
249 "Global node " << libmesh_node_id <<
" not found in sequential node map!");
251 int sequential_index = cast_int<int>
274 if (holes.size() > 0)
276 unsigned hole_index = 0;
277 for (
Point ihole : holes)
278 tetgen_wrapper.
set_hole(hole_index++,
292 std::ostringstream oss;
296 if (quality_constraint != 0)
297 oss <<
"q" << std::fixed << quality_constraint;
299 if (volume_constraint != 0)
300 oss <<
"a" << std::fixed << volume_constraint;
302 std::string params = oss.str();
309 REAL x=0., y=0., z=0.;
324 for (
unsigned int i=old_nodesnum; i<num_nodes; i++)
348 unsigned int node_labels[4];
350 for (
unsigned int i=0; i<num_elements; i++)
356 for (
auto j : elem->node_index_range())
393 for (
auto & node : this->
_mesh.node_ptr_range())
424 #endif // #ifdef LIBMESH_HAVE_TETGEN virtual Node *& set_node(const unsigned int i)
A Node is like a Point, but with more information.
virtual void smooth() override
Redefinition of the smooth function from the base class.
void set_switches(std::string new_switches)
Method to set switches to tetgen, allowing for different behaviours.
int get_element_node(unsigned i, unsigned j)
unsigned check_hull_integrity()
This function checks the integrity of the current set of elements in the Mesh to see if they comprise...
void get_output_node(unsigned i, REAL &x, REAL &y, REAL &z)
int get_triface_node(unsigned i, unsigned j)
std::string _switches
Parameter controlling the behaviour of tetgen.
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...
This is the base class from which all geometric element types are derived.
std::vector< unsigned > _sequential_to_libmesh_node_map
We should not assume libmesh nodes are numbered sequentially...
void set_vertex(unsigned i, unsigned j, unsigned k, int nodeindex)
Sets index of ith facet, jth polygon, kth vertex in the TetGen input.
TetGenMeshInterface(UnstructuredMesh &mesh)
Constructor.
The libMesh namespace provides an interface to certain functionality in the library.
void allocate_polygon_vertexlist(unsigned i, unsigned j, int numofvertices)
Allocates memory, sets number of vertices for polygon j, facet i in the TetGen input.
void allocate_pointlist(int numofpoints)
Allocates memory, sets number of nodes in the TetGen input.
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
int get_numberoftetrahedra()
Real distance(const Point &p)
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.
This class defines the data structures necessary for Laplace smoothing.
void triangulate_pointset()
Method invokes TetGen library to compute a Delaunay tetrahedralization from the nodes point set...
The TetGenWrapper provides an interface for basic access to TetGen data structures and methods...
virtual void delete_elem(Elem *e)=0
Removes element e from the mesh.
void set_node(unsigned i, REAL x, REAL y, REAL z)
Sets coordinates of point i in the TetGen input.
The UnstructuredMesh class is derived from the MeshBase class.
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...
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
void assign_nodes_to_elem(unsigned *node_labels, Elem *elem)
Assigns the node IDs contained in the 'node_labels' array to 'elem'.
void run_tetgen()
Starts the triangulation.
void regenerate_id_sets()
Clears and regenerates the cached sets of ids.
void allocate_facet_polygonlist(unsigned i, int numofpolygons)
Allocates memory, sets number of polygons for facet i in the TetGen input.
void set_switches(std::string_view s)
Method to set TetGen commandline switches -p Tetrahedralizes a piecewise linear complex (...
void set_hole(unsigned i, REAL x, REAL y, REAL z)
Sets coordinates of hole i in the TetGen input.
Real _desired_volume
The desired volume for the elements in the resulting mesh.
void process_hull_integrity_result(unsigned result)
This function prints an informative message and crashes based on the output of the check_hull_integri...
int get_numberoftrifaces()
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...
IntRange< unsigned short > node_index_range() const
void allocate_facetlist(int numoffacets, int numofholes)
Allocates memory, sets number of facets, holes in the TetGen input.
void fill_pointlist(TetGenWrapper &wrapper)
This function copies nodes from the _mesh into TetGen's pointlist.
void pointset_convexhull()
Method invokes TetGen library to compute a Delaunay tetrahedralization from the nodes point set...
virtual dof_id_type n_elem() const =0
virtual const Node * node_ptr(const dof_id_type i) const =0
virtual void triangulate() override
Method invokes TetGen library to compute a Delaunay tetrahedralization.
UnstructuredMesh & _mesh
Local reference to the mesh we are working with.
A Point defines a location in LIBMESH_DIM dimensional Real space.
void delete_2D_hull_elements()
Delete original convex hull elements from the Mesh after performing a Delaunay tetrahedralization.
virtual dof_id_type n_nodes() const =0
bool _smooth_after_generating
Flag which tells whether we should smooth the mesh after it is generated.
Class MeshTetInterface provides an abstract interface for tetrahedralization of meshes by subclasses...