]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix unintended terminate in gmsh-API version of read_msh 16714/head
authorSimon Sticko <simon@sticko.se>
Mon, 4 Mar 2024 18:10:11 +0000 (19:10 +0100)
committerSimon Sticko <simon@sticko.se>
Mon, 11 Mar 2024 17:02:23 +0000 (18:02 +0100)
The gmsh-api version of GridIn::read_msh uses a throw statement without
specifing an exception. The intention is cleary to throw any exception
but the standard specifies that doing this outside of a catch block
should lead to std::terminate being called. Set the boundary_id directly,
without using throw, and add a test that covers this part of the function.

source/grid/grid_in.cc
tests/gmsh/gmsh_api_06.cc [new file with mode: 0644]
tests/gmsh/gmsh_api_06.with_gmsh_with_api=on.output [new file with mode: 0644]

index 4521749cc468a2266cc0b5f3ed7eee35f53ae139..71687f1984bc44470426b5cbc3f3f8d9d6edad57 100644 (file)
@@ -2788,51 +2788,54 @@ GridIn<dim, spacedim>::read_msh(const std::string &fname)
             std::string name;
             gmsh::model::getPhysicalName(entity_dim, physical_tag, name);
             if (!name.empty())
-              try
-                {
-                  std::map<std::string, int> id_names;
-                  Patterns::Tools::to_value(name, id_names);
-                  bool throw_anyway      = false;
-                  bool found_boundary_id = false;
-                  // If the above did not throw, we keep going, and retrieve
-                  // all the information that we know how to translate.
-                  for (const auto &it : id_names)
-                    {
-                      const auto &name = it.first;
-                      const auto &id   = it.second;
-                      if (entity_dim == dim && name == "MaterialID")
-                        {
-                          boundary_id = static_cast<types::boundary_id>(id);
-                          found_boundary_id = true;
-                        }
-                      else if (entity_dim < dim && name == "BoundaryID")
-                        {
-                          boundary_id = static_cast<types::boundary_id>(id);
-                          found_boundary_id = true;
-                        }
-                      else if (name == "ManifoldID")
-                        manifold_id = static_cast<types::manifold_id>(id);
-                      else
-                        // We did not recognize one of the keys. We'll fall
-                        // back to setting the boundary id to the physical tag
-                        // after reading all strings.
-                        throw_anyway = true;
-                    }
-                  // If we didn't find a BoundaryID:XX or MaterialID:XX, and
-                  // something was found but not recognized, then we set the
-                  // material id or boundary id in the catch block below,
-                  // using directly the physical tag
-                  if (throw_anyway && !found_boundary_id)
-                    throw;
-                }
-              catch (...)
-                {
-                  // When the above didn't work, we revert to the old
-                  // behaviour: the physical tag itself is interpreted either
-                  // as a material_id or a boundary_id, and no manifold id is
-                  // known
-                  boundary_id = physical_tag;
-                }
+              {
+                // Patterns::Tools::to_value throws an exception, if it can not
+                // convert name to a map from string to int.
+                try
+                  {
+                    std::map<std::string, int> id_names;
+                    Patterns::Tools::to_value(name, id_names);
+                    bool found_unrecognized_tag = false;
+                    bool found_boundary_id      = false;
+                    // If the above did not throw, we keep going, and retrieve
+                    // all the information that we know how to translate.
+                    for (const auto &it : id_names)
+                      {
+                        const auto &name = it.first;
+                        const auto &id   = it.second;
+                        if (entity_dim == dim && name == "MaterialID")
+                          {
+                            boundary_id = static_cast<types::boundary_id>(id);
+                            found_boundary_id = true;
+                          }
+                        else if (entity_dim < dim && name == "BoundaryID")
+                          {
+                            boundary_id = static_cast<types::boundary_id>(id);
+                            found_boundary_id = true;
+                          }
+                        else if (name == "ManifoldID")
+                          manifold_id = static_cast<types::manifold_id>(id);
+                        else
+                          // We did not recognize one of the keys. We'll fall
+                          // back to setting the boundary id to the physical tag
+                          // after reading all strings.
+                          found_unrecognized_tag = true;
+                      }
+                    // If we didn't find a BoundaryID:XX or MaterialID:XX, and
+                    // something was found but not recognized, then we set the
+                    // material id or boundary using the physical tag directly.
+                    if (found_unrecognized_tag && !found_boundary_id)
+                      boundary_id = physical_tag;
+                  }
+                catch (...)
+                  {
+                    // When the above didn't work, we revert to the old
+                    // behaviour: the physical tag itself is interpreted either
+                    // as a material_id or a boundary_id, and no manifold id is
+                    // known
+                    boundary_id = physical_tag;
+                  }
+              }
           }
 
       // Get the mesh elements for the entity (dim, tag):
diff --git a/tests/gmsh/gmsh_api_06.cc b/tests/gmsh/gmsh_api_06.cc
new file mode 100644 (file)
index 0000000..b43e4fd
--- /dev/null
@@ -0,0 +1,95 @@
+// ------------------------------------------------------------------------
+//
+// SPDX-License-Identifier: LGPL-2.1-or-later
+// Copyright (C) 2021 - 2023 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// Part of the source code is dual licensed under Apache-2.0 WITH
+// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
+// governing the source code and code contributions can be found in
+// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
+//
+// ------------------------------------------------------------------------
+
+// Check that gmsh api correctly reads and writes a mesh with manifold
+// information, in all coordinate dimensions
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_in.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/grid/tria.h>
+
+#include <boost/algorithm/string.hpp>
+
+#include <fstream>
+
+#include "../tests.h"
+
+/*
+ * Test that we can use the gmsh-API version of GridIn::read_msh, when the
+ * description of the "Physical name" does not contain MaterialID, BoundaryID,
+ * or ManifoldID.
+ *
+ * Create a hypercube triangulation, write it to file so that we get a valid
+ * msh-file. Read the filecontent back in and replace MaterialId with a string
+ * that will not be recognized. Write it back to the file and try to read the
+ * triangulation from the msh-file.
+ */
+template <int dim>
+void
+test()
+{
+  // Create a mesh, write it to file.
+  Triangulation<dim> original_mesh;
+  const double       left     = 0.;
+  const double       right    = 1.;
+  const bool         colorize = true;
+  GridGenerator::hyper_cube(original_mesh, left, right, colorize);
+  // Set materialID to something nonzero so that we
+  // get a physical name "MaterialID: 1" in the msh-file.
+  auto cell_on_original = original_mesh.begin_active();
+  cell_on_original->set_material_id(17);
+
+  const std::string mesh_filename = "output.msh";
+  GridOut           grid_out;
+  grid_out.write_msh(original_mesh, mesh_filename);
+
+  // Read the file back and change the physical group description to be
+  // unrecognizable.
+  std::ifstream     input(mesh_filename, std::ifstream::in);
+  std::stringstream buffert;
+  buffert << input.rdbuf();
+  std::string filecontent = buffert.str();
+  input.close();
+
+  boost::replace_all(filecontent, "MaterialID", "Unrecognizable");
+  boost::replace_all(filecontent, "BoundaryID", "Unrecognizable");
+  boost::replace_all(filecontent, "ManifoldID", "Unrecognizable");
+
+  // Write the mesh to file and read the triangulation back in.
+  std::ofstream output(mesh_filename);
+  output << filecontent;
+  output.close();
+
+  Triangulation<dim> triangulation_from_file;
+  GridIn<dim>        grid_in(triangulation_from_file);
+  grid_in.read_msh(mesh_filename);
+
+  // The function should set material_id and boundary_id based on the physical
+  // tag. Write these to log to check them.
+  const auto cell = triangulation_from_file.begin_active();
+  deallog << "MaterialID = " << cell->material_id() << std::endl;
+  for (unsigned int f : cell->face_indices())
+    deallog << "BoundaryID = " << cell->face(f)->boundary_id() << std::endl;
+}
+
+int
+main(int argc, char **argv)
+{
+  // gmsh might be build with mpi support enabled.
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv);
+  initlog();
+
+  test<2>();
+}
diff --git a/tests/gmsh/gmsh_api_06.with_gmsh_with_api=on.output b/tests/gmsh/gmsh_api_06.with_gmsh_with_api=on.output
new file mode 100644 (file)
index 0000000..acc1037
--- /dev/null
@@ -0,0 +1,6 @@
+
+DEAL::MaterialID = 4
+DEAL::BoundaryID = 0
+DEAL::BoundaryID = 1
+DEAL::BoundaryID = 2
+DEAL::BoundaryID = 3

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.