]> https://gitweb.dealii.org/ - dealii.git/commitdiff
- Fixed a bug in grid_tools compute_mesh_predicate_bounding_box that 9053/head
authorblaisb <blais.bruno@gmail.com>
Fri, 15 Nov 2019 11:53:11 +0000 (06:53 -0500)
committerblaisb <blais.bruno@gmail.com>
Fri, 15 Nov 2019 12:59:45 +0000 (07:59 -0500)
lead to a crash when the triangulation was <dim, spacedim> with dim <
spacedim (e.g. <2,3>).
- Adapted the tests to have the right <dim,spacedim> template parameters
for all classes
NOTE: The commented parts in the test  will still crash if you have dim<spacedim (e.g. <2,3>)
AND if the maxbbox < number of cell in the refinement level. This is
because the bounding box merge does not work when the cells are part of
a Triangulation<dim,spacedim> with dim!=spacedim. This part should be
further investigated.

source/grid/grid_tools.cc
tests/grid/grid_tools_compute_mesh_predicate_bounding_box_1.cc

index a328fcfe6460028a197bba820e27dc2bf69150e3..d75bda3261e428bd124929b7d751c235d16f608f 100644 (file)
@@ -2021,7 +2021,7 @@ namespace GridTools
         for (; i < active_cells.size(); ++i)
           if (predicate(active_cells[i]))
             for (unsigned int v = 0;
-                 v < GeometryInfo<spacedim>::vertices_per_cell;
+                 v < GeometryInfo<MeshType::dimension>::vertices_per_cell;
                  ++v)
               for (unsigned int d = 0; d < spacedim; ++d)
                 {
index 1376f748eecb49c2e44f0a956d37bb963ec5573c..fb85be1dcf5c6495ba2c24f6aadaf632460ab16b 100644 (file)
 #include "../tests.h"
 
 // predicate: test if the cell is locally owned
-template <int spacedim>
+template <int dim, int spacedim>
 bool
-pred_locally_owned(
-  const typename parallel::distributed::Triangulation<spacedim>::cell_iterator
-    &cell)
+pred_locally_owned(const typename parallel::distributed::
+                     Triangulation<dim, spacedim>::cell_iterator &cell)
 {
   return cell->is_locally_owned();
 }
@@ -53,19 +52,19 @@ test_hypercube(unsigned int ref, unsigned int max_bbox)
           << " refinement: " << ref << " max number of boxes: " << max_bbox
           << std::endl;
 
-  parallel::distributed::Triangulation<spacedim> tria(mpi_communicator);
+  parallel::distributed::Triangulation<dim, spacedim> tria(mpi_communicator);
   GridGenerator::hyper_cube(tria);
   tria.refine_global(4);
 
-  typedef typename parallel::distributed::Triangulation<
-    spacedim>::active_cell_iterator cell_iterator;
+  typedef typename parallel::distributed::Triangulation<dim, spacedim>::
+    active_cell_iterator cell_iterator;
 
   std::function<bool(const cell_iterator &)> predicate =
-    pred_locally_owned<spacedim>;
+    pred_locally_owned<dim, spacedim>;
 
   std::vector<BoundingBox<spacedim>> local_bbox =
     GridTools::compute_mesh_predicate_bounding_box<
-      parallel::distributed::Triangulation<spacedim>>(
+      parallel::distributed::Triangulation<dim, spacedim>>(
       tria, predicate, ref, true, max_bbox);
 
   if (local_bbox.size() > 0)
@@ -94,8 +93,7 @@ test_hypercube(unsigned int ref, unsigned int max_bbox)
   // Looking if every point is at least inside a bounding box
   for (; cell < endc; ++cell)
     if (cell->is_locally_owned())
-      for (unsigned int v = 0; v < GeometryInfo<spacedim>::vertices_per_cell;
-           ++v)
+      for (unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell; ++v)
         {
           bool inside_a_box = false;
           for (unsigned int i = 0; i < local_bbox.size(); ++i)
@@ -126,6 +124,10 @@ main(int argc, char *argv[])
   test_hypercube<2>(1, 4);
   test_hypercube<2>(2, 4);
   test_hypercube<2>(2, 2);
+  // test_hypercube<2, 3>(0, 4);
+  // test_hypercube<2, 3>(1, 4);
+  // test_hypercube<2, 3>(2, 4);
+  // test_hypercube<2, 3>(2, 2);
   test_hypercube<3>(0, 4);
   test_hypercube<3>(1, 4);
   test_hypercube<3>(2, 4);

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.