Loading [MathJax]/extensions/tex2jax.js
libMesh
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
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  WrappedNgMesh ngmesh;
178 
179  // Create surface mesh in the WrappedNgMesh
180  {
181  // NetGen appears to use ONE-BASED numbering for its nodes, and
182  // since it doesn't return an id when adding nodes we'll have to
183  // track the numbering ourselves.
184  int ng_id = 1;
185 
186  auto create_surface_component =
187  [this, &ng_id, &ng_to_libmesh_id, &ngmesh]
188  (UnstructuredMesh & srcmesh, bool hole_mesh)
189  {
190  // Keep track of what nodes we've already added to the Netgen
191  // mesh vs what nodes we need to add. We'll keep track by id,
192  // not by point location. I don't know if Netgen can handle
193  // multiple nodes with the same point location, but if they can
194  // it's not going to be *us* who breaks that feature.
195  std::unordered_map<dof_id_type, int> libmesh_to_ng_id;
196 
197  // Keep track of what nodes we've already added to the main
198  // mesh from a hole mesh.
199  std::unordered_map<dof_id_type, dof_id_type> hole_to_main_mesh_id;
200 
201  // Use a separate array for passing points to NetGen, just in case
202  // we're not using double-precision ourselves.
203  std::array<double, 3> point_val;
204 
205  // And an array for element vertices
206  std::array<int, 3> elem_nodes;
207 
208  for (const auto * elem : srcmesh.element_ptr_range())
209  {
210  // If someone has triangles we can't triangulate, we have a
211  // problem
212  if (elem->type() == TRI6 ||
213  elem->type() == TRI7)
214  libmesh_not_implemented_msg
215  ("Netgen tetrahedralization currently only supports TRI3 boundaries");
216 
217  // If someone has non-triangles, let's just ignore them.
218  if (elem->type() != TRI3)
219  continue;
220 
221  for (int ni : make_range(3))
222  {
223  // Just using the "invert_trigs" option in NetGen params
224  // doesn't work for me, so we'll have to have properly
225  // oriented the tris earlier.
226  auto & elem_node = hole_mesh ? elem_nodes[2-ni] : elem_nodes[ni];
227 
228  const Node & n = elem->node_ref(ni);
229  auto n_id = n.id();
230  if (hole_mesh)
231  {
232  if (auto it = hole_to_main_mesh_id.find(n_id);
233  it != hole_to_main_mesh_id.end())
234  {
235  n_id = it->second;
236  }
237  else
238  {
239  Node * n_new = this->_mesh.add_point(n);
240  const dof_id_type n_new_id = n_new->id();
241  hole_to_main_mesh_id.emplace(n_id, n_new_id);
242  n_id = n_new_id;
243  }
244  }
245 
246  if (auto it = libmesh_to_ng_id.find(n_id);
247  it != libmesh_to_ng_id.end())
248  {
249  const int existing_ng_id = it->second;
250  elem_node = existing_ng_id;
251  }
252  else
253  {
254  for (auto i : make_range(3))
255  point_val[i] = double(n(i));
256 
257  Ng_AddPoint(ngmesh, point_val.data());
258 
259  ng_to_libmesh_id[ng_id] = n_id;
260  libmesh_to_ng_id[n_id] = ng_id;
261  elem_node = ng_id;
262  ++ng_id;
263  }
264  }
265 
266  Ng_AddSurfaceElement(ngmesh, NG_TRIG, elem_nodes.data());
267  }
268  };
269 
270  create_surface_component(this->_mesh, false);
271 
272  if (_holes)
273  for (const std::unique_ptr<UnstructuredMesh> & h : *_holes)
274  create_surface_component(*h, true);
275  }
276 
277  auto result = Ng_GenerateVolumeMesh(ngmesh, &params);
278  handle_ng_result(result);
279 
280  const int n_elem = Ng_GetNE(ngmesh);
281 
282  libmesh_error_msg_if (n_elem <= 0,
283  "NetGen failed to generate any tetrahedra");
284 
285  const dof_id_type n_points = Ng_GetNP(ngmesh);
286  const dof_id_type old_nodes = this->_mesh.n_nodes();
287 
288  // Netgen may have generated new interior nodes
289  if (n_points != old_nodes)
290  {
291  std::array<double, 3> point_val;
292 
293  // We should only be getting new nodes if we asked for them
294  if (!_desired_volume)
295  {
296  std::cout <<
297  "NetGen output " << n_points <<
298  " points when we gave it " <<
299  old_nodes << " and disabled refinement\n" <<
300  "If new interior points are acceptable in your mesh, please set\n" <<
301  "a non-zero desired_volume to indicate that. If new interior\n" <<
302  "points are not acceptable in your mesh, you may need a different\n" <<
303  "(non-advancing-front?) mesh generator." << std::endl;
304  libmesh_error();
305  }
306  else
307  for (auto i : make_range(old_nodes, n_points))
308  {
309  // i+1 since ng uses ONE-BASED numbering
310  Ng_GetPoint (ngmesh, i+1, point_val.data());
311  const Point p(point_val[0], point_val[1], point_val[2]);
312  Node * n_new = this->_mesh.add_point(p);
313  const dof_id_type n_new_id = n_new->id();
314  ng_to_libmesh_id[i+1] = n_new_id;
315  }
316  }
317 
318  for (auto * elem : this->_mesh.element_ptr_range())
319  this->_mesh.delete_elem(elem);
320 
321  BoundaryInfo * bi = & this->_mesh.get_boundary_info();
322 
323  for (auto i : make_range(n_elem))
324  {
325  // Enough data to return even a Tet10 without a segfault if nglib
326  // went nuts
327  int ngnodes[11];
328 
329  // i+1 since we must be 1-based with these ids too...
330  Ng_Volume_Element_Type ngtype =
331  Ng_GetVolumeElement(ngmesh, i+1, ngnodes);
332 
333  // But really nglib shouldn't go nuts
334  libmesh_assert(ngtype == NG_TET);
335  libmesh_ignore(ngtype);
336 
337  auto elem = this->_mesh.add_elem(Elem::build_with_id(TET4, i));
338  for (auto n : make_range(4))
339  {
340  const dof_id_type node_id =
341  libmesh_map_find(ng_to_libmesh_id, ngnodes[n]);
342  elem->set_node(n) = this->_mesh.node_ptr(node_id);
343  }
344 
345  // NetGen and we disagree about node numbering orientation
346  elem->orient(bi);
347  }
348 
349  // We don't need our holes anymore. Delete their serializers
350  // first to avoid dereferencing dangling pointers.
351  hole_serializers.clear();
352  if (_holes)
353  _holes->clear();
354 
355  // Send our data to other ranks, then fix it up together
357  this->_mesh.prepare_for_use();
358 }
359 
360 
361 
362 } // namespace libMesh
363 
364 
365 #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:963
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:741
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:160
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
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 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:429
void broadcast(MeshBase &) const
Finds all the processors that may contain elements that neighbor my elements.
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