]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Clean up hp-coarsening code. 11011/head
authorMarc Fehling <mafehling.git@gmail.com>
Tue, 6 Oct 2020 01:34:25 +0000 (19:34 -0600)
committerMarc Fehling <mafehling.git@gmail.com>
Wed, 7 Oct 2020 21:43:40 +0000 (15:43 -0600)
12 files changed:
doc/news/changes/minor/20201005Fehling [new file with mode: 0644]
include/deal.II/dofs/dof_accessor.h
include/deal.II/dofs/dof_accessor.templates.h
include/deal.II/hp/fe_collection.h
source/distributed/cell_weights.cc
source/distributed/error_predictor.cc
source/distributed/solution_transfer.cc
source/dofs/dof_handler.cc
source/hp/refinement.cc
source/numerics/solution_transfer.cc
tests/grid/accessor_04.cc [new file with mode: 0644]
tests/grid/accessor_04.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20201005Fehling b/doc/news/changes/minor/20201005Fehling
new file mode 100644 (file)
index 0000000..a9d4679
--- /dev/null
@@ -0,0 +1,4 @@
+New: Helper function DoFCellAccessor::dominated_future_fe_on_children()
+to clean up code for hp-coarsening.
+<br>
+(Marc Fehling, 2020/10/05)
index 14fe641575df0d3052e654ae47ec57d089bc3319..75b130f0f91c752e33788d261c9659c3032cdc2c 100644 (file)
@@ -2039,6 +2039,48 @@ public:
    */
   void
   clear_future_fe_index() const;
+
+  /**
+   * Return the index of the finite element from the entire hp::FECollection
+   * that is dominated by those assigned as future finite elements to the
+   * children of this cell.
+   *
+   * We find the corresponding finite element among the future finite elements
+   * on the children of this cell. If none of them qualify, we extend our search
+   * on the whole hp::FECollection, which is the element that describes the
+   * smallest finite element space that includes all future finite elements
+   * assigned to the children. If the function is not able to find a finite
+   * element at all, an assertion will be triggered.
+   *
+   * In this way, we determine the finite element of the parent cell in case of
+   * h-coarsening in the hp-context.
+   *
+   * @note This function can only be called on direct parent cells, i.e.,
+   * non-active cells whose children are all active.
+   */
+  unsigned int
+  dominated_future_fe_on_children() const;
+  /**
+   * @}
+   */
+
+  /**
+   * @name Exceptions
+   */
+  /**
+   * @{
+   */
+
+  /**
+   * Exception
+   *
+   * @ingroup Exceptions
+   */
+  DeclExceptionMsg(
+    ExcNoDominatedFiniteElementOnChildren,
+    "No FiniteElement has been found in your FECollection that is "
+    "dominated by all children of a cell you are trying to coarsen!");
+
   /**
    * @}
    */
index 3284b305fbe5318101188d7dd987312dd6ab632f..40f1038ad0a3d3b09178e52efb38ebb06566d609 100644 (file)
@@ -2920,6 +2920,40 @@ DoFCellAccessor<dimension_, space_dimension_, level_dof_access>::
 
 
 
+template <int dimension_, int space_dimension_, bool level_dof_access>
+inline unsigned int
+DoFCellAccessor<dimension_, space_dimension_, level_dof_access>::
+  dominated_future_fe_on_children() const
+{
+  Assert(!this->is_active(),
+         ExcMessage(
+           "You ask for information on children of this cell which is only "
+           "available for active cells. This cell has no children."));
+
+  std::set<unsigned int> future_fe_indices_children;
+  for (const auto &child : this->child_iterators())
+    {
+      Assert(
+        child->is_active(),
+        ExcMessage(
+          "You ask for information on children of this cell which is only "
+          "available for active cells. One of its children is not active."));
+      future_fe_indices_children.insert(child->future_fe_index());
+    }
+  Assert(!future_fe_indices_children.empty(), ExcInternalError());
+
+  const unsigned int future_fe_index =
+    this->dof_handler->get_fe_collection().find_dominated_fe_extended(
+      future_fe_indices_children, /*codim=*/0);
+
+  Assert(future_fe_index != numbers::invalid_unsigned_int,
+         ExcNoDominatedFiniteElementOnChildren());
+
+  return future_fe_index;
+}
+
+
+
 template <int dimension_, int space_dimension_, bool level_dof_access>
 template <typename number, typename OutputVector>
 inline void
index c6e6f87956611e6faf5ccf883fe9e28e8cd2bddf..5c319f378be9cb2faffff320bf4e293daebfeb62 100644 (file)
@@ -735,16 +735,6 @@ namespace hp
      */
     DeclException0(ExcNoFiniteElements);
 
-    /**
-     * Exception
-     *
-     * @ingroup Exceptions
-     */
-    DeclExceptionMsg(
-      ExcNoDominatedFiniteElementAmongstChildren,
-      "No FiniteElement has been found in your FECollection that is "
-      "dominated by all children of a cell you are trying to coarsen!");
-
     /**
      * @}
      */
index 72f0c0850f6ed79712beed1f4fcac0ee180b685d..8aaefbcd4ee6c841ec5344c8c1295100a186acb7 100644 (file)
@@ -192,28 +192,14 @@ namespace parallel
           break;
 
         case Triangulation<dim, spacedim>::CELL_COARSEN:
-          {
-            std::set<unsigned int> fe_indices_children;
-            for (unsigned int child_index = 0; child_index < cell->n_children();
-                 ++child_index)
-              {
-                const auto &child = cell->child(child_index);
-                Assert(child->is_active() && child->coarsen_flag_set(),
-                       typename dealii::Triangulation<
-                         dim>::ExcInconsistentCoarseningFlags());
-
-                fe_indices_children.insert(child->future_fe_index());
-              }
-            Assert(!fe_indices_children.empty(), ExcInternalError());
-
-            fe_index =
-              dof_handler.get_fe_collection().find_dominated_fe_extended(
-                fe_indices_children, /*codim=*/0);
-
-            Assert(fe_index != numbers::invalid_unsigned_int,
-                   typename dealii::hp::FECollection<
-                     dim>::ExcNoDominatedFiniteElementAmongstChildren());
-          }
+#ifdef DEBUG
+          for (const auto &child : cell->child_iterators())
+            Assert(child->is_active() && child->coarsen_flag_set(),
+                   typename dealii::Triangulation<
+                     dim>::ExcInconsistentCoarseningFlags());
+#endif
+
+          fe_index = cell->dominated_future_fe_on_children();
           break;
 
         default:
index ef12125b818aeb50d30da79ced9ac8f895d5788f..65d8bbba29183ec75690b014ee6f9e8c9c5b50b9 100644 (file)
@@ -233,22 +233,15 @@ namespace parallel
                 // parent cell after h-adaptation analogously to
                 // dealii::internal::hp::DoFHandlerImplementation::Implementation::
                 //   collect_fe_indices_on_cells_to_be_refined()
-                std::set<unsigned int> fe_indices_children;
-                for (unsigned int child_index = 0;
-                     child_index < cell->n_children();
-                     ++child_index)
-                  {
-                    const auto child = cell->child(child_index);
-                    Assert(child->is_active() && child->coarsen_flag_set(),
-                           typename dealii::Triangulation<
-                             dim>::ExcInconsistentCoarseningFlags());
-
-                    fe_indices_children.insert(child->future_fe_index());
-                  }
+#  ifdef DEBUG
+                for (const auto &child : cell->child_iterators())
+                  Assert(child->is_active() && child->coarsen_flag_set(),
+                         typename dealii::Triangulation<
+                           dim>::ExcInconsistentCoarseningFlags());
+#  endif
 
                 const unsigned int future_fe_index =
-                  dof_handler->get_fe_collection().find_dominated_fe_extended(
-                    fe_indices_children, /*codim=*/0);
+                  cell->dominated_future_fe_on_children();
 
                 const unsigned int future_fe_degree =
                   dof_handler->get_fe_collection()[future_fe_index].degree;
index 481990c56d8550e58a1ec93f59ca8e6a28e6ebe1..c1862d0126a2567161640e1b1ea29b2d5d5e31ef 100644 (file)
@@ -311,28 +311,14 @@ namespace parallel
                   // In case of coarsening, we need to find a suitable fe index
                   // for the parent cell. We choose the 'least dominant fe'
                   // on all children from the associated FECollection.
-                  std::set<unsigned int> fe_indices_children;
-                  for (unsigned int child_index = 0;
-                       child_index < cell->n_children();
-                       ++child_index)
-                    {
-                      const auto &child = cell->child(child_index);
-                      Assert(child->is_active() && child->coarsen_flag_set(),
-                             typename dealii::Triangulation<
-                               dim>::ExcInconsistentCoarseningFlags());
-
-                      fe_indices_children.insert(child->future_fe_index());
-                    }
-                  Assert(!fe_indices_children.empty(), ExcInternalError());
-
-                  fe_index =
-                    dof_handler->get_fe_collection().find_dominated_fe_extended(
-                      fe_indices_children, /*codim=*/0);
-
-                  Assert(fe_index != numbers::invalid_unsigned_int,
-                         typename dealii::hp::FECollection<
-                           dim>::ExcNoDominatedFiniteElementAmongstChildren());
-
+#  ifdef DEBUG
+                  for (const auto &child : cell->child_iterators())
+                    Assert(child->is_active() && child->coarsen_flag_set(),
+                           typename dealii::Triangulation<
+                             dim>::ExcInconsistentCoarseningFlags());
+#  endif
+
+                  fe_index = cell->dominated_future_fe_on_children();
                   break;
                 }
 
index 22a0c26f5acb848c048948fdffd9209e5bf71c5b..9f75ec85c4f01b8e62be65cdf69601ab17bff346 100644 (file)
@@ -1879,30 +1879,16 @@ namespace internal
                         // based on the 'least dominant finite element' of its
                         // children. Consider the childrens' hypothetical future
                         // index when they have been flagged for p-refinement.
-                        std::set<unsigned int> fe_indices_children;
-                        for (unsigned int child_index = 0;
-                             child_index < parent->n_children();
-                             ++child_index)
-                          {
-                            const auto sibling = parent->child(child_index);
-                            Assert(sibling->is_active() &&
-                                     sibling->coarsen_flag_set(),
-                                   typename dealii::Triangulation<
-                                     dim>::ExcInconsistentCoarseningFlags());
-
-                            fe_indices_children.insert(
-                              sibling->future_fe_index());
-                          }
-                        Assert(!fe_indices_children.empty(),
-                               ExcInternalError());
+#ifdef DEBUG
+                        for (const auto &child : parent->child_iterators())
+                          Assert(child->is_active() &&
+                                   child->coarsen_flag_set(),
+                                 typename dealii::Triangulation<
+                                   dim>::ExcInconsistentCoarseningFlags());
+#endif
 
                         const unsigned int fe_index =
-                          dof_handler.fe_collection.find_dominated_fe_extended(
-                            fe_indices_children, /*codim=*/0);
-
-                        Assert(fe_index != numbers::invalid_unsigned_int,
-                               typename dealii::hp::FECollection<dim>::
-                                 ExcNoDominatedFiniteElementAmongstChildren());
+                          parent->dominated_future_fe_on_children();
 
                         fe_transfer->coarsened_cells_fe_index.insert(
                           {parent, fe_index});
@@ -1999,8 +1985,8 @@ namespace internal
                                                      /*codim=*/0);
 
           Assert(dominated_fe_index != numbers::invalid_unsigned_int,
-                 typename dealii::hp::FECollection<
-                   dim>::ExcNoDominatedFiniteElementAmongstChildren());
+                 (typename dealii::DoFCellAccessor<dim, spacedim, false>::
+                    ExcNoDominatedFiniteElementOnChildren()));
 
           return dominated_fe_index;
         }
index 3e80a37954317f0454e67ae08abb2374d5275d5d..7c9c02a4949df228a18626f1f5ebac0d8ff371f8 100644 (file)
@@ -542,29 +542,15 @@ namespace hp
                 if (future_fe_indices_on_coarsened_cells.find(parent) ==
                     future_fe_indices_on_coarsened_cells.end())
                   {
-                    std::set<unsigned int> fe_indices_children;
-                    for (unsigned int child_index = 0;
-                         child_index < parent->n_children();
-                         ++child_index)
-                      {
-                        const auto &child = parent->child(child_index);
-                        Assert(child->is_active() && child->coarsen_flag_set(),
-                               typename dealii::Triangulation<
-                                 dim>::ExcInconsistentCoarseningFlags());
-
-                        fe_indices_children.insert(child->future_fe_index());
-                      }
-                    Assert(!fe_indices_children.empty(), ExcInternalError());
+#ifdef DEBUG
+                    for (const auto &child : parent->child_iterators())
+                      Assert(child->is_active() && child->coarsen_flag_set(),
+                             typename dealii::Triangulation<
+                               dim>::ExcInconsistentCoarseningFlags());
+#endif
 
                     parent_future_fe_index =
-                      dof_handler.get_fe_collection()
-                        .find_dominated_fe_extended(fe_indices_children,
-                                                    /*codim=*/0);
-
-                    Assert(
-                      parent_future_fe_index != numbers::invalid_unsigned_int,
-                      typename dealii::hp::FECollection<
-                        dim>::ExcNoDominatedFiniteElementAmongstChildren());
+                      parent->dominated_future_fe_on_children();
 
                     future_fe_indices_on_coarsened_cells.insert(
                       {parent, parent_future_fe_index});
index 67c25059ae9550be4be1bd0f701a828723a33d40..4ea1fa8aa6820900308a89b8733ea623be71375c 100644 (file)
@@ -372,8 +372,10 @@ SolutionTransfer<dim, VectorType, DoFHandlerType>::
               fe_indices_children, /*codim=*/0);
 
           Assert(target_fe_index != numbers::invalid_unsigned_int,
-                 typename dealii::hp::FECollection<
-                   dim>::ExcNoDominatedFiniteElementAmongstChildren());
+                 (typename dealii::DoFCellAccessor<
+                   dim,
+                   DoFHandlerType::space_dimension,
+                   false>::ExcNoDominatedFiniteElementOnChildren()));
 
           const unsigned int dofs_per_cell =
             dof_handler->get_fe(target_fe_index).n_dofs_per_cell();
diff --git a/tests/grid/accessor_04.cc b/tests/grid/accessor_04.cc
new file mode 100644 (file)
index 0000000..578c2d9
--- /dev/null
@@ -0,0 +1,75 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Verify functionality of
+// DoFCellAccessor:dominated_future_fe_on_children()
+
+#include <deal.II/base/logstream.h>
+
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/hp/fe_collection.h>
+
+#include "../tests.h"
+
+
+
+template <int dim>
+void
+test()
+{
+  Triangulation<dim> tria;
+  GridGenerator::hyper_cube(tria);
+  tria.refine_global(1);
+
+  hp::FECollection<dim> fes;
+  fes.push_back(FE_Q<dim>(1));
+  fes.push_back(FE_Q<dim>(2));
+
+  DoFHandler<dim> dofh(tria, /*hp_capability_enabled=*/true);
+  dofh.set_fe(fes);
+  dofh.begin_active()->set_active_fe_index(1);
+
+  const auto &       parent = dofh.begin(/*level=*/0);
+  const unsigned int parent_future_fe =
+    parent->dominated_future_fe_on_children();
+
+  deallog << parent_future_fe << std::endl;
+}
+
+
+int
+main()
+{
+  initlog();
+
+  deallog.push("1d");
+  test<1>();
+  deallog.pop();
+  deallog.push("2d");
+  test<2>();
+  deallog.pop();
+  deallog.push("3d");
+  test<3>();
+  deallog.pop();
+}
diff --git a/tests/grid/accessor_04.output b/tests/grid/accessor_04.output
new file mode 100644 (file)
index 0000000..8b5ba12
--- /dev/null
@@ -0,0 +1,4 @@
+
+DEAL:1d::1
+DEAL:2d::1
+DEAL:3d::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.