libMesh
mesh_netgen_interface.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2023 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 #include "libmesh/libmesh_config.h"
19 #ifdef LIBMESH_HAVE_NETGEN
20 
21 
22 // C++ includes
23 #include <sstream>
24 
25 // Local includes
26 #include "libmesh/mesh_netgen_interface.h"
27 
28 #include "libmesh/boundary_info.h"
29 #include "libmesh/cell_tet4.h"
30 #include "libmesh/face_tri3.h"
31 #include "libmesh/mesh_communication.h"
32 #include "libmesh/unstructured_mesh.h"
33 #include "libmesh/utility.h" // libmesh_map_find
34 
35 namespace {
36 
37 // RAII for exception safety
38 class WrappedNgMesh
39 {
40 public:
41  WrappedNgMesh() {
42  _ngmesh = nglib::Ng_NewMesh();
43  }
44 
45  ~WrappedNgMesh() {
46  nglib::Ng_DeleteMesh(_ngmesh);
47  }
48 
49  void clear() {
50  nglib::Ng_DeleteMesh(_ngmesh);
51  _ngmesh = nglib::Ng_NewMesh();
52  }
53 
54  operator nglib::Ng_Mesh* () {
55  return _ngmesh;
56  }
57 
58 private:
59  nglib::Ng_Mesh * _ngmesh;
60 };
61 
62 }
63 
64 namespace libMesh
65 {
66 
67 //----------------------------------------------------------------------
68 // NetGenMeshInterface class members
71  _serializer(mesh)
72 {
73 }
74 
75 
76 
78 {
79  using namespace nglib;
80 
81  // We're hoping to do volume_to_surface_mesh in parallel at least,
82  // but then we'll need to serialize any hole meshes to rank 0 so it
83  // can use them in serial.
84 
85  const BoundingBox mesh_bb =
87 
88  std::vector<MeshSerializer> hole_serializers;
89  if (_holes)
90  for (std::unique_ptr<UnstructuredMesh> & hole : *_holes)
91  {
92  const BoundingBox hole_bb =
94 
95  libmesh_error_msg_if
96  (!mesh_bb.contains(hole_bb),
97  "Found hole with bounding box " << hole_bb <<
98  "\nextending outside of mesh bounding box " << mesh_bb);
99 
100  hole_serializers.emplace_back
101  (*hole, /* need_serial */ true,
102  /* serial_only_needed_on_proc_0 */ true);
103  }
104 
105  // If we're not rank 0, we're just going to wait for rank 0 to call
106  // Netgen, then receive its data afterward, we're not going to hope
107  // that Netgen does the exact same thing on every processor.
108  if (this->_mesh.processor_id() != 0)
109  {
110  // We don't need our holes anymore. Delete their serializers
111  // first to avoid dereferencing dangling pointers.
112  hole_serializers.clear();
113  if (_holes)
114  _holes->clear();
115 
116  // Receive the mesh data rank 0 will send later, then fix it up
117  // together
119  this->_mesh.prepare_for_use();
120  return;
121  }
122 
123  this->check_hull_integrity();
124 
125  Ng_Meshing_Parameters params;
126 
127  // Override any default parameters we might need to, to avoid
128  // inserting nodes we don't want.
129  params.uselocalh = false;
130  params.minh = 0;
131  params.elementsperedge = 1;
132  params.elementspercurve = 1;
133  params.closeedgeenable = false;
134  params.closeedgefact = 0;
135  params.minedgelenenable = false;
136  params.minedgelen = 0;
137 
138  // Try to get a no-extra-nodes mesh if we're asked to, or try to
139  // translate our desired volume into NetGen terms otherwise.
140  //
141  // Spoiler alert: all we can do is try; NetGen uses a marching front
142  // algorithm that can insert extra nodes despite all my best
143  // efforts.
144  if (_desired_volume == 0) // shorthand for "no refinement"
145  {
146  params.maxh = std::numeric_limits<double>::max();
147  params.fineness = 0; // "coarse" in the docs
148  params.grading = 1; // "aggressive local grading" to avoid smoothing??
149 
150  // Turning off optimization steps avoids another opportunity for
151  // Netgen to try to add more nodes.
152  params.optsteps_3d = 0;
153  }
154  else
155  params.maxh = double(std::pow(_desired_volume, 1./3.));
156 
157  // Keep track of how NetGen copies of nodes map back to our original
158  // nodes, so we can connect new elements to nodes correctly.
159  std::unordered_map<int, dof_id_type> ng_to_libmesh_id;
160 
161  auto handle_ng_result = [](Ng_Result result) {
162  static const std::vector<std::string> result_types =
163  {"Netgen error", "Netgen success", "Netgen surface input error",
164  "Netgen volume failure", "Netgen STL input error",
165  "Netgen surface failure", "Netgen file not found"};
166 
167  if (result+1 >= 0 &&
168  std::size_t(result+1) < result_types.size())
169  libmesh_error_msg_if
170  (result, "Ng_GenerateVolumeMesh failed: " <<
171  result_types[result+1]);
172  else
173  libmesh_error_msg
174  ("Ng_GenerateVolumeMesh failed with an unknown error code");
175  };
176 
177  // Keep track of what boundary ids we want to assign to each new
178  // triangle. We'll give the outer boundary BC 0, and give holes ids
179  // starting from 1.
180  // We key on sorted tuples of node ids to identify a side.
181  std::unordered_map<std::array<dof_id_type,3>,
182  boundary_id_type, libMesh::hash> side_boundary_id;
183 
184  auto insert_id = []
185  (std::array<dof_id_type,3> & array,
186  dof_id_type n_id)
187  {
188  libmesh_assert_less(n_id, DofObject::invalid_id);
189  unsigned int i=0;
190  while (array[i] < n_id)
191  ++i;
192  while (i < 3)
193  std::swap(array[i++], n_id);
194  };
195 
196  WrappedNgMesh ngmesh;
197 
198  // Create surface mesh in the WrappedNgMesh
199  {
200  // NetGen appears to use ONE-BASED numbering for its nodes, and
201  // since it doesn't return an id when adding nodes we'll have to
202  // track the numbering ourselves.
203  int ng_id = 1;
204 
205  auto create_surface_component =
206  [this, &ng_id, &ng_to_libmesh_id, &ngmesh, &side_boundary_id, &insert_id]
207  (UnstructuredMesh & srcmesh, bool hole_mesh,
208  boundary_id_type bcid)
209  {
210  // Keep track of what nodes we've already added to the Netgen
211  // mesh vs what nodes we need to add. We'll keep track by id,
212  // not by point location. I don't know if Netgen can handle
213  // multiple nodes with the same point location, but if they can
214  // it's not going to be *us* who breaks that feature.
215  std::unordered_map<dof_id_type, int> libmesh_to_ng_id;
216 
217  // Keep track of what nodes we've already added to the main
218  // mesh from a hole mesh.
219  std::unordered_map<dof_id_type, dof_id_type> hole_to_main_mesh_id;
220 
221  // Use a separate array for passing points to NetGen, just in case
222  // we're not using double-precision ourselves.
223  std::array<double, 3> point_val;
224 
225  // And an array for element vertices
226  std::array<int, 3> elem_nodes;
227 
228  for (const auto * elem : srcmesh.element_ptr_range())
229  {
230  // If someone has triangles we can't triangulate, we have a
231  // problem
232  if (elem->type() == TRI6 ||
233  elem->type() == TRI7)
234  libmesh_not_implemented_msg
235  ("Netgen tetrahedralization currently only supports TRI3 boundaries");
236 
237  // If someone has non-triangles, let's just ignore them.
238  if (elem->type() != TRI3)
239  continue;
240 
241  std::array<dof_id_type,3> sorted_ids =
244 
245  for (int ni : make_range(3))
246  {
247  // Just using the "invert_trigs" option in NetGen params
248  // doesn't work for me, so we'll have to have properly
249  // oriented the tris earlier.
250  auto & elem_node = hole_mesh ? elem_nodes[2-ni] : elem_nodes[ni];
251 
252  const Node & n = elem->node_ref(ni);
253  auto n_id = n.id();
254  if (hole_mesh)
255  {
256  if (auto it = hole_to_main_mesh_id.find(n_id);
257  it != hole_to_main_mesh_id.end())
258  {
259  n_id = it->second;
260  }
261  else
262  {
263  Node * n_new = this->_mesh.add_point(n);
264  const dof_id_type n_new_id = n_new->id();
265  hole_to_main_mesh_id.emplace(n_id, n_new_id);
266  n_id = n_new_id;
267  }
268  }
269 
270  if (auto it = libmesh_to_ng_id.find(n_id);
271  it != libmesh_to_ng_id.end())
272  {
273  const int existing_ng_id = it->second;
274  elem_node = existing_ng_id;
275  }
276  else
277  {
278  for (auto i : make_range(3))
279  point_val[i] = double(n(i));
280 
281  Ng_AddPoint(ngmesh, point_val.data());
282 
283  ng_to_libmesh_id[ng_id] = n_id;
284  libmesh_to_ng_id[n_id] = ng_id;
285  elem_node = ng_id;
286  ++ng_id;
287  }
288 
289  insert_id(sorted_ids, n_id);
290  }
291 
292  side_boundary_id[sorted_ids] = bcid;
293 
294  Ng_AddSurfaceElement(ngmesh, NG_TRIG, elem_nodes.data());
295  }
296  };
297 
298  // Number the outer boundary 0, and the holes starting from 1
299  boundary_id_type bcid = 0;
300 
301  create_surface_component(this->_mesh, false, bcid);
302 
303  if (_holes)
304  for (const std::unique_ptr<UnstructuredMesh> & h : *_holes)
305  create_surface_component(*h, true, ++bcid);
306  }
307 
308  auto result = Ng_GenerateVolumeMesh(ngmesh, &params);
309  handle_ng_result(result);
310 
311  const int n_elem = Ng_GetNE(ngmesh);
312 
313  libmesh_error_msg_if (n_elem <= 0,
314  "NetGen failed to generate any tetrahedra");
315 
316  const dof_id_type n_points = Ng_GetNP(ngmesh);
317  const dof_id_type old_nodes = this->_mesh.n_nodes();
318 
319  // Netgen may have generated new interior nodes
320  if (n_points != old_nodes)
321  {
322  std::array<double, 3> point_val;
323 
324  // We should only be getting new nodes if we asked for them
325  if (!_desired_volume)
326  {
327  std::cout <<
328  "NetGen output " << n_points <<
329  " points when we gave it " <<
330  old_nodes << " and disabled refinement\n" <<
331  "If new interior points are acceptable in your mesh, please set\n" <<
332  "a non-zero desired_volume to indicate that. If new interior\n" <<
333  "points are not acceptable in your mesh, you may need a different\n" <<
334  "(non-advancing-front?) mesh generator." << std::endl;
335  libmesh_error();
336  }
337  else
338  for (auto i : make_range(old_nodes, n_points))
339  {
340  // i+1 since ng uses ONE-BASED numbering
341  Ng_GetPoint (ngmesh, i+1, point_val.data());
342  const Point p(point_val[0], point_val[1], point_val[2]);
343  Node * n_new = this->_mesh.add_point(p);
344  const dof_id_type n_new_id = n_new->id();
345  ng_to_libmesh_id[i+1] = n_new_id;
346  }
347  }
348 
349  for (auto * elem : this->_mesh.element_ptr_range())
350  this->_mesh.delete_elem(elem);
351 
352  BoundaryInfo * bi = & this->_mesh.get_boundary_info();
353 
354  for (auto i : make_range(n_elem))
355  {
356  // Enough data to return even a Tet10 without a segfault if nglib
357  // went nuts
358  int ngnodes[11];
359 
360  // i+1 since we must be 1-based with these ids too...
361  Ng_Volume_Element_Type ngtype =
362  Ng_GetVolumeElement(ngmesh, i+1, ngnodes);
363 
364  // But really nglib shouldn't go nuts
365  libmesh_assert(ngtype == NG_TET);
366  libmesh_ignore(ngtype);
367 
368  auto elem = this->_mesh.add_elem(Elem::build_with_id(TET4, i));
369  for (auto n : make_range(4))
370  {
371  const dof_id_type node_id =
372  libmesh_map_find(ng_to_libmesh_id, ngnodes[n]);
373  elem->set_node(n, this->_mesh.node_ptr(node_id));
374  }
375 
376  // NetGen and we disagree about node numbering orientation
377  elem->orient(bi);
378 
379  for (auto s : make_range(4))
380  {
381  std::array<dof_id_type,3> sorted_ids =
384 
385  std::vector<unsigned int> nos = elem->nodes_on_side(s);
386  for (auto n : nos)
387  insert_id(sorted_ids, elem->node_id(n));
388 
389  if (auto it = side_boundary_id.find(sorted_ids);
390  it != side_boundary_id.end())
391  bi->add_side(elem, s, it->second);
392  }
393  }
394 
395  // We don't need our holes anymore. Delete their serializers
396  // first to avoid dereferencing dangling pointers.
397  hole_serializers.clear();
398  if (_holes)
399  _holes->clear();
400 
401  // Send our data to other ranks, then fix it up together
403  this->_mesh.prepare_for_use();
404 }
405 
406 
407 
408 } // namespace libMesh
409 
410 
411 #endif // #ifdef LIBMESH_HAVE_NETGEN
NetGenMeshInterface(UnstructuredMesh &mesh)
Constructor.
std::unique_ptr< std::vector< std::unique_ptr< UnstructuredMesh > > > _holes
A pointer to a vector of meshes each defining a hole.
A Node is like a Point, but with more information.
Definition: node.h:52
unsigned check_hull_integrity()
This function checks the integrity of the current set of elements in the Mesh to see if they comprise...
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:969
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:759
MeshBase & mesh
virtual void triangulate() override
Method invokes NetGen library to compute a tetrahedralization.
The libMesh namespace provides an interface to certain functionality in the library.
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:165
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 &...)
T pow(const T &x)
Definition: utility.h:328
int8_t boundary_id_type
Definition: id_types.h:51
virtual void delete_elem(Elem *e)=0
Removes element e from the mesh.
This is the MeshCommunication class.
dof_id_type id() const
Definition: dof_object.h:828
The UnstructuredMesh class is derived from the MeshBase class.
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
libmesh_assert(ctx)
static const dof_id_type invalid_id
An invalid id to distinguish an uninitialized DofObject.
Definition: dof_object.h:482
static BoundingBox volume_to_surface_mesh(UnstructuredMesh &mesh)
Remove volume elements from the given mesh, after converting their outer boundary faces to surface el...
Defines a Cartesian bounding box by the two corner extremum.
Definition: bounding_box.h:40
Real _desired_volume
The desired volume for the elements in the resulting mesh.
static std::unique_ptr< Elem > build_with_id(const ElemType type, dof_id_type id)
Calls the build() method above with a nullptr parent, and additionally sets the newly-created Elem&#39;s ...
Definition: elem.C:558
void broadcast(MeshBase &) const
Finds all the processors that may contain elements that neighbor my elements.
void add_side(const dof_id_type elem, const unsigned short int side, const boundary_id_type id)
Add side side of element number elem with boundary id id to the boundary information data structure...
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
bool contains(const BoundingBox &) const
Definition: bounding_box.h:209
virtual const Node * node_ptr(const dof_id_type i) const =0
processor_id_type processor_id() const
UnstructuredMesh & _mesh
Local reference to the mesh we are working with.
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
virtual dof_id_type n_nodes() const =0
Class MeshTetInterface provides an abstract interface for tetrahedralization of meshes by subclasses...
uint8_t dof_id_type
Definition: id_types.h:67