]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix hp::Refinement::choose_p_over_h() for hp-coarsening on p:s:Triangulation. 11980/head
authorMarc Fehling <mafehling.git@gmail.com>
Mon, 29 Mar 2021 20:12:34 +0000 (14:12 -0600)
committerMarc Fehling <mafehling.git@gmail.com>
Mon, 29 Mar 2021 20:32:56 +0000 (14:32 -0600)
source/hp/refinement.cc
tests/sharedtria/hp_choose_p_over_h.cc

index 5eebda33550b251fc6115545874f462ec1b8a956..fb91081755ca71c7f5ea3e8f218268603c026d7b 100644 (file)
@@ -22,6 +22,7 @@
 #include <deal.II/distributed/shared_tria.h>
 #include <deal.II/distributed/tria_base.h>
 
+#include <deal.II/dofs/dof_accessor.templates.h>
 #include <deal.II/dofs/dof_handler.h>
 
 #include <deal.II/grid/grid_refinement.h>
@@ -655,33 +656,14 @@ namespace hp
         dof_handler.has_hp_capabilities(),
         (typename dealii::DoFHandler<dim, spacedim>::ExcOnlyAvailableWithHP()));
 
-      // Siblings of cells to be coarsened may not be owned by the same
-      // processor. We will exchange coarsening flags on ghost cells and
-      // temporarily store them.
-      std::map<CellId, std::pair<bool, bool>> ghost_buffer;
-
-      if (dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
-            &dof_handler.get_triangulation()))
-        {
-          auto pack =
-            [](const typename dealii::DoFHandler<dim,
-                                                 spacedim>::active_cell_iterator
-                 &cell) -> std::pair<bool, bool> {
-            return {cell->coarsen_flag_set(), cell->future_fe_index_set()};
-          };
-
-          auto unpack =
-            [&ghost_buffer](const typename dealii::DoFHandler<dim, spacedim>::
-                              active_cell_iterator &    cell,
-                            const std::pair<bool, bool> pair) -> void {
-            ghost_buffer.emplace(cell->id(), pair);
-          };
-
-          GridTools::exchange_cell_data_to_ghosts<
-            std::pair<bool, bool>,
-            dealii::DoFHandler<dim, spacedim>>(dof_handler, pack, unpack);
-        }
-
+      // Ghost siblings might occur on parallel::shared::Triangulation objects.
+      // We need information about future FE indices on all locally relevant
+      // cells here, and thus communicate them.
+      if (dynamic_cast<const parallel::shared::Triangulation<dim, spacedim> *>(
+            &dof_handler.get_triangulation()) != nullptr)
+        dealii::internal::hp::DoFHandlerImplementation::
+          communicate_future_fe_indices(
+            const_cast<dealii::DoFHandler<dim, spacedim> &>(dof_handler));
 
       for (const auto &cell : dof_handler.active_cell_iterators())
         if (cell->is_locally_owned() && cell->future_fe_index_set())
@@ -711,12 +693,27 @@ namespace hp
                           }
                         else if (child->is_ghost())
                           {
-                            const std::pair<bool, bool> &flags =
-                              ghost_buffer[child->id()];
+                            // The case of siblings being owned by different
+                            // processors can only occur for
+                            // parallel::shared::Triangulation objects.
+                            Assert(
+                              (dynamic_cast<const parallel::shared::
+                                              Triangulation<dim, spacedim> *>(
+                                 &dof_handler.get_triangulation()) != nullptr),
+                              ExcInternalError());
 
-                            if (flags.first)
+                            if (child->coarsen_flag_set())
                               ++h_flagged_children;
-                            if (flags.second)
+                            // The public interface does not allow to access
+                            // future FE indices on ghost cells. However, we
+                            // need this information here and thus call the
+                            // internal function that does not check for cell
+                            // ownership.
+                            if (dealii::internal::
+                                  DoFCellAccessorImplementation::
+                                    Implementation::
+                                      future_fe_index_set<dim, spacedim, false>(
+                                        *child))
                               ++p_flagged_children;
                           }
                         else
index 5ad88c3d4bd50f4cd72b0ffef9b9d7390cbd184e..dce9a701e05f015fae7a90aee59f991cac8ace4e 100644 (file)
@@ -69,11 +69,11 @@ test()
           for (unsigned int i = 0; i < cell->n_children(); ++i)
             {
               const auto &child = cell->child(i);
+
+              child->set_coarsen_flag();
+
               if (child->is_locally_owned())
-                {
-                  child->set_future_fe_index(1);
-                  child->set_coarsen_flag();
-                }
+                child->set_future_fe_index(1);
             }
         }
       else if (cell->id().to_string() == "1_0:")
@@ -84,12 +84,12 @@ test()
           for (unsigned int i = 0; i < cell->n_children(); ++i)
             {
               const auto &child = cell->child(i);
+
+              child->set_coarsen_flag();
+
               if (child->is_locally_owned())
-                {
-                  if (i == 0)
-                    child->set_future_fe_index(1);
-                  child->set_coarsen_flag();
-                }
+                if (i == 0)
+                  child->set_future_fe_index(1);
             }
         }
     }

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.