libMesh
Functions
miscellaneous_ex6.C File Reference

Go to the source code of this file.

Functions

void triangulate_domain (const Parallel::Communicator &comm)
 
void tetrahedralize_domain (const Parallel::Communicator &comm)
 
void add_cube_convex_hull_to_mesh (MeshBase &mesh, Point lower_limit, Point upper_limit)
 
int main (int argc, char **argv)
 

Function Documentation

◆ add_cube_convex_hull_to_mesh()

void add_cube_convex_hull_to_mesh ( MeshBase mesh,
Point  lower_limit,
Point  upper_limit 
)

Definition at line 250 of file miscellaneous_ex6.C.

References libMesh::MeshBase::add_elem(), libMesh::MeshBase::add_point(), libMesh::Elem::build(), libMesh::MeshTools::Generation::build_cube(), libMesh::ParallelObject::comm(), libMesh::HEX8, libMesh::DofObject::id(), libMesh::libmesh_ignore(), mesh, libMesh::MeshTools::n_elem(), libMesh::MeshBase::node_ptr(), libMesh::TetGenMeshInterface::pointset_convexhull(), libMesh::Elem::set_node(), and libMesh::TRI3.

Referenced by tetrahedralize_domain().

253 {
254 #ifdef LIBMESH_HAVE_TETGEN
255  ReplicatedMesh cube_mesh(mesh.comm(), 3);
256 
257  unsigned n_elem = 1;
258 
260  n_elem, n_elem, n_elem, // n. elements in each direction
261  lower_limit(0), upper_limit(0),
262  lower_limit(1), upper_limit(1),
263  lower_limit(2), upper_limit(2),
264  HEX8);
265 
266  // The pointset_convexhull() algorithm will ignore the Hex8s
267  // in the Mesh, and just construct the triangulation
268  // of the convex hull.
269  TetGenMeshInterface t(cube_mesh);
270  t.pointset_convexhull();
271 
272  // Now add all nodes from the boundary of the cube_mesh to the input mesh.
273 
274  // Map from "node id in cube_mesh" -> "node id in mesh". Initially inserted
275  // with a dummy value, later to be assigned a value by the input mesh.
276  std::map<unsigned, unsigned> node_id_map;
277  typedef std::map<unsigned, unsigned>::iterator iterator;
278 
279  for (auto & elem : cube_mesh.element_ptr_range())
280  for (auto s : elem->side_index_range())
281  if (elem->neighbor_ptr(s) == nullptr)
282  {
283  // Add the node IDs of this side to the set
284  std::unique_ptr<Elem> side = elem->side_ptr(s);
285 
286  for (auto n : side->node_index_range())
287  node_id_map.emplace(side->node_id(n), /*dummy_value=*/0);
288  }
289 
290  // For each node in the map, insert it into the input mesh and keep
291  // track of the ID assigned.
292  for (iterator it=node_id_map.begin(); it != node_id_map.end(); ++it)
293  {
294  // Id of the node in the cube mesh
295  unsigned id = (*it).first;
296 
297  // Pointer to node in the cube mesh
298  Node & old_node = cube_mesh.node_ref(id);
299 
300  // Add geometric point to input mesh
301  Node * new_node = mesh.add_point (old_node);
302 
303  // Track ID value of new_node in map
304  (*it).second = new_node->id();
305  }
306 
307  // With the points added and the map data structure in place, we are
308  // ready to add each TRI3 element of the cube_mesh to the input Mesh
309  // with proper node assignments
310  for (auto & old_elem : cube_mesh.element_ptr_range())
311  if (old_elem->type() == TRI3)
312  {
313  Elem * new_elem = mesh.add_elem(Elem::build(TRI3));
314 
315  // Assign nodes in new elements. Since this is an example,
316  // we'll do it in several steps.
317  for (auto i : old_elem->node_index_range())
318  {
319  // Mapping to node ID in input mesh
320  unsigned int new_node_id = libmesh_map_find(node_id_map, old_elem->node_id(i));
321 
322  // Node pointer assigned from input mesh
323  new_elem->set_node(i) = mesh.node_ptr(new_node_id);
324  }
325  }
326 #else
327  // Avoid compiler warnings
328  libmesh_ignore(mesh, lower_limit, upper_limit);
329 #endif // LIBMESH_HAVE_TETGEN
330 }
Class TetGenMeshInterface provides an interface for tetrahedralization of meshes using the TetGen lib...
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
virtual Node *& set_node(const unsigned int i)
Definition: elem.h:2381
A Node is like a Point, but with more information.
Definition: node.h:52
dof_id_type n_elem(const MeshBase::const_element_iterator &begin, const MeshBase::const_element_iterator &end)
Count up the number of elements of a specific type (as defined by an iterator range).
Definition: mesh_tools.C:850
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
MeshBase & mesh
const Parallel::Communicator & comm() const
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.
void libmesh_ignore(const Args &...)
dof_id_type id() const
Definition: dof_object.h:823
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
virtual const Node * node_ptr(const dof_id_type i) const =0
void build_cube(UnstructuredMesh &mesh, const unsigned int nx=0, const unsigned int ny=0, const unsigned int nz=0, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const Real zmin=0., const Real zmax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
Builds a (elements) cube.

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 55 of file miscellaneous_ex6.C.

References libMesh::TriangleWrapper::init(), libMesh::out, tetrahedralize_domain(), and triangulate_domain().

56 {
57  // Initialize libMesh and any dependent libraries, like in example 2.
58  LibMeshInit init (argc, argv);
59 
60  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
61 
62  libMesh::out << "Triangulating an L-shaped domain with holes" << std::endl;
63 
64  // 1.) 2D triangulation of L-shaped domain with three holes of different shape
65  triangulate_domain(init.comm());
66 
67  libmesh_example_requires(3 <= LIBMESH_DIM, "3D support");
68 
69  libMesh::out << "Tetrahedralizing a prismatic domain with a hole" << std::endl;
70 
71  // 2.) 3D tetrahedralization of rectangular domain with hole.
73 
74  return 0;
75 }
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
OStreamProxy out
void tetrahedralize_domain(const Parallel::Communicator &comm)
void triangulate_domain(const Parallel::Communicator &comm)

◆ tetrahedralize_domain()

void tetrahedralize_domain ( const Parallel::Communicator comm)

Definition at line 169 of file miscellaneous_ex6.C.

References add_cube_convex_hull_to_mesh(), libMesh::MeshBase::find_neighbors(), libMesh::libmesh_ignore(), mesh, libMesh::MeshBase::prepare_for_use(), libMesh::TetGenMeshInterface::set_switches(), libMesh::TetGenMeshInterface::triangulate_conformingDelaunayMesh_carvehole(), and libMesh::MeshBase::write().

Referenced by main().

170 {
171 #ifdef LIBMESH_HAVE_TETGEN
172  // The algorithm is broken up into several steps:
173  // 1.) A convex hull is constructed for a rectangular hole.
174  // 2.) A convex hull is constructed for the domain exterior.
175  // 3.) Neighbor information is updated so TetGen knows there is a convex hull
176  // 4.) A vector of hole points is created.
177  // 5.) The domain is tetrahedralized, the mesh is written out, etc.
178 
179  // The mesh we will eventually generate
180  ReplicatedMesh mesh(comm, 3);
181 
182  // Lower and Upper bounding box limits for a rectangular hole within the unit cube.
183  Point hole_lower_limit(0.2, 0.2, 0.4);
184  Point hole_upper_limit(0.8, 0.8, 0.6);
185 
186  // 1.) Construct a convex hull for the hole
187  add_cube_convex_hull_to_mesh(mesh, hole_lower_limit, hole_upper_limit);
188 
189  // 2.) Generate elements comprising the outer boundary of the domain.
191  Point(0., 0., 0.),
192  Point(1., 1., 1.));
193 
194  // 3.) Update neighbor information so that TetGen can verify there is a convex hull.
196 
197  // 4.) Set up vector of hole points
198  std::vector<Point> hole(1);
199  hole[0] = Point(0.5*(hole_lower_limit + hole_upper_limit));
200 
201  // 5.) Set parameters and tetrahedralize the domain
202 
203  // 0 means "use TetGen default value"
204  double quality_constraint = 2.0;
205 
206  // The volume constraint determines the max-allowed tetrahedral
207  // volume in the Mesh. TetGen will split cells which are larger than
208  // this size
209  double volume_constraint = 0.001;
210 
211  // Construct the Delaunay tetrahedralization
213 
214  // In debug mode, set the tetgen switches
215  // -V (verbose) and
216  // -CC (check consistency and constrained Delaunay)
217  // In optimized mode, only switch Q (quiet) is set.
218  // For more options, see tetgen website:
219  // (http://wias-berlin.de/software/tetgen/1.5/doc/manual/manual005.html#cmd-m)
220 #ifdef DEBUG
221  t.set_switches("VCC");
222 #endif
223  t.triangulate_conformingDelaunayMesh_carvehole(hole,
224  quality_constraint,
225  volume_constraint);
226 
227  // Find neighbors, etc in preparation for writing out the Mesh
229 
230  // Finally, write out the result
231  mesh.write("hole_3D.e");
232 #else
233  // Avoid compiler warnings
234  libmesh_ignore(comm);
235 #endif // LIBMESH_HAVE_TETGEN
236 }
Class TetGenMeshInterface provides an interface for tetrahedralization of meshes using the TetGen lib...
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
void add_cube_convex_hull_to_mesh(MeshBase &mesh, Point lower_limit, Point upper_limit)
void prepare_for_use(const bool skip_renumber_nodes_and_elements, const bool skip_find_neighbors)
Prepare a newly ecreated (or read) mesh for use.
Definition: mesh_base.C:710
MeshBase & mesh
void libmesh_ignore(const Args &...)
virtual void find_neighbors(const bool reset_remote_elements=false, const bool reset_current_list=true)=0
Locate element face (edge in 2D) neighbors.
virtual void write(const std::string &name)=0
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39

◆ triangulate_domain()

void triangulate_domain ( const Parallel::Communicator comm)

Definition at line 80 of file miscellaneous_ex6.C.

References libMesh::MeshBase::add_point(), libMesh::TriangulatorInterface::attach_hole_list(), libMesh::TriangulatorInterface::desired_area(), libMesh::libmesh_ignore(), mesh, libMesh::pi, libMesh::TriangulatorInterface::PSLG, libMesh::Real, libMesh::TriangulatorInterface::smooth_after_generating(), libMesh::TriangleInterface::triangulate(), libMesh::TriangulatorInterface::triangulation_type(), and libMesh::MeshBase::write().

Referenced by main().

81 {
82 #ifdef LIBMESH_HAVE_TRIANGLE
83  // Use typedefs for slightly less typing.
84  typedef TriangleInterface::Hole Hole;
85  typedef TriangleInterface::PolygonHole PolygonHole;
86  typedef TriangleInterface::ArbitraryHole ArbitraryHole;
87 
88  // Libmesh mesh that will eventually be created.
89  Mesh mesh(comm, 2);
90 
91  // The points which make up the L-shape:
92  mesh.add_point(Point(0., 0.));
93  mesh.add_point(Point(0., -1.));
94  mesh.add_point(Point(-1., -1.));
95  mesh.add_point(Point(-1., 1.));
96  mesh.add_point(Point(1., 1.));
97  mesh.add_point(Point(1., 0.));
98 
99  // Declare the TriangleInterface object. This is where
100  // we can set parameters of the triangulation and where the
101  // actual triangulate function lives.
103 
104  // Customize the variables for the triangulation
105  t.desired_area() = .01;
106 
107  // A Planar Straight Line Graph (PSLG) is essentially a list
108  // of segments which have to exist in the final triangulation.
109  // For an L-shaped domain, Triangle will compute the convex
110  // hull of boundary points if we do not specify the PSLG.
111  // The PSLG algorithm is also required for triangulating domains
112  // containing holes
113  t.triangulation_type() = TriangleInterface::PSLG;
114 
115  // Turn on/off Laplacian mesh smoothing after generation.
116  // By default this is on.
117  t.smooth_after_generating() = true;
118 
119  // Define holes...
120 
121  // hole_1 is a circle (discretized by 50 points)
122  PolygonHole hole_1(Point(-0.5, 0.5), // center
123  0.25, // radius
124  50); // n. points
125 
126  // hole_2 is itself a triangle
127  PolygonHole hole_2(Point(0.5, 0.5), // center
128  0.1, // radius
129  3); // n. points
130 
131  // hole_3 is an ellipse of 100 points which we define here
132  Point ellipse_center(-0.5, -0.5);
133  const unsigned int n_ellipse_points=100;
134  std::vector<Point> ellipse_points(n_ellipse_points);
135  const Real
136  dtheta = 2*libMesh::pi / static_cast<Real>(n_ellipse_points),
137  a = .1,
138  b = .2;
139 
140  for (unsigned int i=0; i<n_ellipse_points; ++i)
141  ellipse_points[i]= Point(ellipse_center(0)+a*cos(i*dtheta),
142  ellipse_center(1)+b*sin(i*dtheta));
143 
144  ArbitraryHole hole_3(ellipse_center, ellipse_points);
145 
146  // Create the vector of Hole*'s ...
147  std::vector<Hole *> holes;
148  holes.push_back(&hole_1);
149  holes.push_back(&hole_2);
150  holes.push_back(&hole_3);
151 
152  // ... and attach it to the triangulator object
153  t.attach_hole_list(&holes);
154 
155  // Triangulate!
156  t.triangulate();
157 
158  // Write the result to file
159  mesh.write("delaunay_l_shaped_hole.e");
160 
161 #else
162  // Avoid compiler warnings
163  libmesh_ignore(comm);
164 #endif // LIBMESH_HAVE_TRIANGLE
165 }
MeshBase & mesh
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.
void libmesh_ignore(const Args &...)
virtual void write(const std::string &name)=0
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
A C++ interface between LibMesh and the Triangle library written by J.R.
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
const Real pi
.
Definition: libmesh.h:274