18 #include "libmesh/libmesh_config.h" 19 #ifdef LIBMESH_HAVE_NETGEN 26 #include "libmesh/mesh_netgen_interface.h" 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" 42 _ngmesh = nglib::Ng_NewMesh();
46 nglib::Ng_DeleteMesh(_ngmesh);
50 nglib::Ng_DeleteMesh(_ngmesh);
51 _ngmesh = nglib::Ng_NewMesh();
54 operator nglib::Ng_Mesh* () {
59 nglib::Ng_Mesh * _ngmesh;
79 using namespace nglib;
88 std::vector<MeshSerializer> hole_serializers;
90 for (std::unique_ptr<UnstructuredMesh> & hole : *
_holes)
97 "Found hole with bounding box " << hole_bb <<
98 "\nextending outside of mesh bounding box " << mesh_bb);
100 hole_serializers.emplace_back
112 hole_serializers.clear();
125 Ng_Meshing_Parameters params;
129 params.uselocalh =
false;
131 params.elementsperedge = 1;
132 params.elementspercurve = 1;
133 params.closeedgeenable =
false;
134 params.closeedgefact = 0;
135 params.minedgelenenable =
false;
136 params.minedgelen = 0;
146 params.maxh = std::numeric_limits<double>::max();
152 params.optsteps_3d = 0;
159 std::unordered_map<int, dof_id_type> ng_to_libmesh_id;
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"};
168 std::size_t(result+1) < result_types.size())
170 (result,
"Ng_GenerateVolumeMesh failed: " <<
171 result_types[result+1]);
174 (
"Ng_GenerateVolumeMesh failed with an unknown error code");
177 WrappedNgMesh ngmesh;
186 auto create_surface_component =
187 [
this, &ng_id, &ng_to_libmesh_id, &ngmesh]
195 std::unordered_map<dof_id_type, int> libmesh_to_ng_id;
199 std::unordered_map<dof_id_type, dof_id_type> hole_to_main_mesh_id;
203 std::array<double, 3> point_val;
206 std::array<int, 3> elem_nodes;
208 for (
const auto * elem : srcmesh.element_ptr_range())
212 if (elem->type() ==
TRI6 ||
213 elem->type() ==
TRI7)
214 libmesh_not_implemented_msg
215 (
"Netgen tetrahedralization currently only supports TRI3 boundaries");
218 if (elem->type() !=
TRI3)
226 auto & elem_node = hole_mesh ? elem_nodes[2-ni] : elem_nodes[ni];
228 const Node & n = elem->node_ref(ni);
232 if (
auto it = hole_to_main_mesh_id.find(n_id);
233 it != hole_to_main_mesh_id.end())
241 hole_to_main_mesh_id.emplace(n_id, n_new_id);
246 if (
auto it = libmesh_to_ng_id.find(n_id);
247 it != libmesh_to_ng_id.end())
249 const int existing_ng_id = it->second;
250 elem_node = existing_ng_id;
255 point_val[i] =
double(n(i));
257 Ng_AddPoint(ngmesh, point_val.data());
259 ng_to_libmesh_id[ng_id] = n_id;
260 libmesh_to_ng_id[n_id] = ng_id;
266 Ng_AddSurfaceElement(ngmesh, NG_TRIG, elem_nodes.data());
270 create_surface_component(this->
_mesh,
false);
273 for (
const std::unique_ptr<UnstructuredMesh> & h : *
_holes)
274 create_surface_component(*h,
true);
277 auto result = Ng_GenerateVolumeMesh(ngmesh, ¶ms);
278 handle_ng_result(result);
280 const int n_elem = Ng_GetNE(ngmesh);
282 libmesh_error_msg_if (
n_elem <= 0,
283 "NetGen failed to generate any tetrahedra");
289 if (n_points != old_nodes)
291 std::array<double, 3> point_val;
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;
307 for (
auto i :
make_range(old_nodes, n_points))
310 Ng_GetPoint (ngmesh, i+1, point_val.data());
311 const Point p(point_val[0], point_val[1], point_val[2]);
314 ng_to_libmesh_id[i+1] = n_new_id;
318 for (
auto * elem : this->
_mesh.element_ptr_range())
330 Ng_Volume_Element_Type ngtype =
331 Ng_GetVolumeElement(ngmesh, i+1, ngnodes);
341 libmesh_map_find(ng_to_libmesh_id, ngnodes[n]);
351 hole_serializers.clear();
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.
unsigned check_hull_integrity()
This function checks the integrity of the current set of elements in the Mesh to see if they comprise...
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.
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.
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 delete_elem(Elem *e)=0
Removes element e from the mesh.
This is the MeshCommunication class.
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...
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.
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's ...
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...
bool contains(const BoundingBox &) const
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.
virtual dof_id_type n_nodes() const =0
Class MeshTetInterface provides an abstract interface for tetrahedralization of meshes by subclasses...