]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Refactored AdaptationStrategies namespace.
authorMarc Fehling <marc.fehling@gmx.net>
Fri, 24 Apr 2020 21:23:08 +0000 (23:23 +0200)
committerMatthias Maier <tamiko@43-1.org>
Mon, 11 May 2020 20:14:27 +0000 (15:14 -0500)
include/deal.II/distributed/cell_data_transfer.h
include/deal.II/distributed/cell_data_transfer.templates.h
include/deal.II/numerics/adaptation_strategies.h [new file with mode: 0644]
include/deal.II/numerics/cell_data_transfer.h
include/deal.II/numerics/cell_data_transfer.templates.h
include/deal.II/numerics/coarsening_strategies.h [deleted file]
source/distributed/cell_weights.cc
source/distributed/solution_transfer.cc
source/hp/dof_handler.cc
source/numerics/solution_transfer.cc
tests/mpi/cell_data_transfer_03.cc

index aa404813c28f6515043a7f2dd15a80539403d6b5..2e3ec97b346c50085e2c4857a310e317c3229ce2 100644 (file)
@@ -22,7 +22,7 @@
 
 #include <deal.II/distributed/tria.h>
 
-#include <deal.II/numerics/coarsening_strategies.h>
+#include <deal.II/numerics/adaptation_strategies.h>
 
 #include <boost/range/iterator_range.hpp>
 
@@ -122,7 +122,7 @@ namespace parallel
      * matching code snippets for both transfer and serialization.
      *
      * @ingroup distributed
-     * @author Marc Fehling, 2018
+     * @author Marc Fehling, 2018 - 2020
      */
     template <int dim, int spacedim = dim, typename VectorType = Vector<double>>
     class CellDataTransfer
@@ -143,14 +143,14 @@ namespace parallel
        * problem. Such strategies can be passed to the CellDataTransfer and
        * parallel::distributed::CellDataTransfer constructors.
        *
-       * @deprecated Use the namespace dealii::CoarseningStrategies instead.
+       * @deprecated Use the namespace dealii::AdaptationStrategies::Coarsening instead.
        */
       struct DEAL_II_DEPRECATED CoarseningStrategies
       {
         /**
-         * @copydoc dealii::CoarseningStrategies::check_equality()
+         * @copydoc dealii::AdaptationStrategies::Coarsening::check_equality()
          *
-         * @deprecated Use dealii::CoarseningStrategies::check_equality() instead.
+         * @deprecated Use dealii::AdaptationStrategies::Coarsening::check_equality() instead.
          */
         DEAL_II_DEPRECATED static typename VectorType::value_type
         check_equality(const typename parallel::distributed::
@@ -172,9 +172,9 @@ namespace parallel
         }
 
         /**
-         * @copydoc dealii::CoarseningStrategies::sum()
+         * @copydoc dealii::AdaptationStrategies::Coarsening::sum()
          *
-         * @deprecated Use dealii::CoarseningStrategies::sum() instead.
+         * @deprecated Use dealii::AdaptationStrategies::Coarsening::sum() instead.
          */
         DEAL_II_DEPRECATED static typename VectorType::value_type
         sum(const typename parallel::distributed::Triangulation<dim, spacedim>::
@@ -192,9 +192,9 @@ namespace parallel
         }
 
         /**
-         * @copydoc dealii::CoarseningStrategies::mean()
+         * @copydoc dealii::AdaptationStrategies::Coarsening::mean()
          *
-         * @deprecated Use dealii::CoarseningStrategies::mean() instead.
+         * @deprecated Use dealii::AdaptationStrategies::Coarsening::mean() instead.
          */
         DEAL_II_DEPRECATED static typename VectorType::value_type
         mean(const typename parallel::distributed::
@@ -215,19 +215,40 @@ namespace parallel
        *   container stores values that differ in size. A varying amount of data
        *   may be packed per cell, if for example the underlying ValueType of
        *   the VectorType container is a container itself.
+       * @param[in] refinement_strategy Function deciding how data will be
+       * stored on refined cells from its parent cell.
        * @param[in] coarsening_strategy Function deciding which data to store on
        *   a cell whose children will get coarsened into.
        */
       CellDataTransfer(
         const parallel::distributed::Triangulation<dim, spacedim>
-          &        triangulation,
-        const bool transfer_variable_size_data = false,
+          &                               triangulation,
+        const bool                        transfer_variable_size_data = false,
+        const std::function<std::vector<value_type>(
+          const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+            &              parent,
+          const value_type parent_value)> refinement_strategy =
+          &dealii::AdaptationStrategies::Refinement::
+            preserve<dim, spacedim, value_type>,
         const std::function<value_type(
+          const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+            &                            parent,
           const std::vector<value_type> &children_values)> coarsening_strategy =
-          &dealii::CoarseningStrategies::check_equality<value_type>);
+          &dealii::AdaptationStrategies::Coarsening::
+            check_equality<dim, spacedim, value_type>);
 
       /**
-       * @copydoc CellDataTransfer::CellDataTransfer
+       * Constructor.
+       *
+       * @param[in] triangulation The triangulation on which all operations will
+       *   happen. At the time when this constructor is called, the refinement
+       *   in question has not happened yet.
+       * @param[in] transfer_variable_size_data Specify whether your VectorType
+       *   container stores values that differ in size. A varying amount of data
+       *   may be packed per cell, if for example the underlying ValueType of
+       *   the VectorType container is a container itself.
+       * @param[in] coarsening_strategy Function deciding which data to store on
+       *   a cell whose children will get coarsened into.
        *
        * @deprecated Use the above constructor instead.
        */
@@ -320,11 +341,21 @@ namespace parallel
        */
       const bool transfer_variable_size_data;
 
+      /**
+       * Function deciding how data will be stored on refined cells from its
+       * parent cell.
+       */
+      const std::function<std::vector<value_type>(
+        const typename Triangulation<dim, spacedim>::cell_iterator &parent,
+        const value_type parent_value)>
+        refinement_strategy;
+
       /**
        * Function deciding on how to process data from children to be stored on
        * the parent cell.
        */
       const std::function<value_type(
+        const typename Triangulation<dim, spacedim>::cell_iterator &parent,
         const std::vector<value_type> &children_values)>
         coarsening_strategy;
 
index 8c60012c1bf1db3db059bef6e6904c13096d256b..260d97a8e54fa8c9e4bb6bdf9088e01adeaff9e4 100644 (file)
@@ -73,11 +73,18 @@ namespace parallel
     template <int dim, int spacedim, typename VectorType>
     CellDataTransfer<dim, spacedim, VectorType>::CellDataTransfer(
       const parallel::distributed::Triangulation<dim, spacedim> &triangulation,
-      const bool transfer_variable_size_data,
+      const bool                        transfer_variable_size_data,
+      const std::function<std::vector<value_type>(
+        const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+          &              parent,
+        const value_type parent_value)> refinement_strategy,
       const std::function<value_type(
+        const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+          &                            parent,
         const std::vector<value_type> &children_values)> coarsening_strategy)
       : triangulation(&triangulation, typeid(*this).name())
       , transfer_variable_size_data(transfer_variable_size_data)
+      , refinement_strategy(refinement_strategy)
       , coarsening_strategy(coarsening_strategy)
       , handle(numbers::invalid_unsigned_int)
     {}
@@ -93,6 +100,8 @@ namespace parallel
                    const VectorType &input_vector)> coarsening_strategy)
       : triangulation(&triangulation, typeid(*this).name())
       , transfer_variable_size_data(transfer_variable_size_data)
+      , refinement_strategy(&dealii::AdaptationStrategies::Refinement::
+                              preserve<dim, spacedim, value_type>)
       , handle(numbers::invalid_unsigned_int)
     {
       value_type (*const *old_strategy)(
@@ -105,20 +114,23 @@ namespace parallel
                          const VectorType &)>();
 
       if (*old_strategy == CoarseningStrategies::check_equality)
-        const_cast<
-          std::function<value_type(const std::vector<value_type> &)> &>(
-          this->coarsening_strategy) =
-          &dealii::CoarseningStrategies::check_equality<value_type>;
+        const_cast<std::function<value_type(
+          const typename dealii::Triangulation<dim, spacedim>::cell_iterator &,
+          const std::vector<value_type> &)> &>(this->coarsening_strategy) =
+          &dealii::AdaptationStrategies::Coarsening::
+            check_equality<dim, spacedim, value_type>;
       else if (*old_strategy == CoarseningStrategies::sum)
-        const_cast<
-          std::function<value_type(const std::vector<value_type> &)> &>(
-          this->coarsening_strategy) =
-          &dealii::CoarseningStrategies::sum<value_type>;
+        const_cast<std::function<value_type(
+          const typename dealii::Triangulation<dim, spacedim>::cell_iterator &,
+          const std::vector<value_type> &)> &>(this->coarsening_strategy) =
+          &dealii::AdaptationStrategies::Coarsening::
+            sum<dim, spacedim, value_type>;
       else if (*old_strategy == CoarseningStrategies::mean)
-        const_cast<
-          std::function<value_type(const std::vector<value_type> &)> &>(
-          this->coarsening_strategy) =
-          &dealii::CoarseningStrategies::mean<value_type>;
+        const_cast<std::function<value_type(
+          const typename dealii::Triangulation<dim, spacedim>::cell_iterator &,
+          const std::vector<value_type> &)> &>(this->coarsening_strategy) =
+          &dealii::AdaptationStrategies::Coarsening::
+            mean<dim, spacedim, value_type>;
       else
         Assert(
           false,
@@ -313,7 +325,7 @@ namespace parallel
                        child_index < cell->n_children();
                        ++child_index)
                     {
-                      const auto child = cell->child(child_index);
+                      const auto &child = cell->child(child_index);
                       Assert(child->is_active() && child->coarsen_flag_set(),
                              typename dealii::Triangulation<
                                dim>::ExcInconsistentCoarseningFlags());
@@ -322,7 +334,7 @@ namespace parallel
                         (**it_input)[child->active_cell_index()];
                     }
 
-                  *it_output = coarsening_strategy(children_values);
+                  *it_output = coarsening_strategy(cell, children_values);
                   break;
                 }
 
@@ -384,19 +396,23 @@ namespace parallel
                                                       spacedim>::CELL_COARSEN:
               // Cell either persists, or has been coarsened.
               // Thus, cell has no (longer) children.
-              (**it_output)[cell->active_cell_index()] = std::move(*it_input);
+              (**it_output)[cell->active_cell_index()] = *it_input;
               break;
 
             case parallel::distributed::Triangulation<dim,
                                                       spacedim>::CELL_REFINE:
-              // Cell has been refined, and is now parent of its children.
-              // Thus, distribute parent's data on its children.
-              for (unsigned int child_index = 0;
-                   child_index < cell->n_children();
-                   ++child_index)
-                (**it_output)[cell->child(child_index)->active_cell_index()] =
-                  std::move(*it_input);
-              break;
+              {
+                // Cell has been refined, and is now parent of its children.
+                // Thus, distribute parent's data on its children.
+                const std::vector<value_type> children_values =
+                  refinement_strategy(cell, *it_input);
+                for (unsigned int child_index = 0;
+                     child_index < cell->n_children();
+                     ++child_index)
+                  (**it_output)[cell->child(child_index)->active_cell_index()] =
+                    children_values[child_index];
+                break;
+              }
 
             default:
               Assert(false, ExcInternalError());
diff --git a/include/deal.II/numerics/adaptation_strategies.h b/include/deal.II/numerics/adaptation_strategies.h
new file mode 100644 (file)
index 0000000..4517891
--- /dev/null
@@ -0,0 +1,355 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 - 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.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii_adaptation_strategies_h
+#define dealii_adaptation_strategies_h
+
+
+#include <deal.II/base/config.h>
+
+#include <deal.II/base/exceptions.h>
+
+#include <deal.II/grid/tria_accessor.h>
+
+#include <algorithm>
+#include <numeric>
+#include <typeinfo>
+#include <vector>
+
+DEAL_II_NAMESPACE_OPEN
+
+/**
+ * When data is transferred during adaptation, it is not trivial to decide
+ * how to process data from former cells on the old mesh that have been changed
+ * into current cells on the new mesh. Or in other words, how data should be
+ * stored in the cells on the adapted mesh.
+ *
+ * In this namespace, we offer a few strategies that cope with this
+ * problem. Such strategies can be passed to the CellDataTransfer and
+ * parallel::distributed::CellDataTransfer constructors.
+ *
+ * @author Marc Fehling, 2019 - 2020
+ */
+namespace AdaptationStrategies
+{
+  /**
+   * For refinement, all strategies take the parent cell and its associated
+   * data. They return a vector containing data for each individual child that
+   * the parent cell will be refined to.
+   *
+   * The ordering of values in the vector for children data corresponds to the
+   * index when calling TriaAccessor::child_index.
+   */
+  namespace Refinement
+  {
+    /**
+     * Return a vector containing copies of data of the parent cell for each
+     * child.
+     *
+     * @f[
+     *   d_{K_c} = d_{K_p}
+     *   \qquad
+     *   \forall K_c \text{ children of } K_p
+     * @f]
+     */
+    template <int dim, int spacedim, typename value_type>
+    std::vector<value_type>
+    preserve(const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+               &              parent,
+             const value_type parent_value);
+
+    /**
+     * Return a vector which contains data of the parent cell being equally
+     * divided among all children.
+     *
+     * @f[
+     *   d_{K_c} = d_{K_p} / n_\text{children}
+     *   \qquad
+     *   \forall K_c \text{ children of } K_p
+     * @f]
+     *
+     * This strategy preserves the $l_1$-norm of the corresponding global data
+     * Vector before and after adaptation.
+     */
+    template <int dim, int spacedim, typename value_type>
+    std::vector<value_type>
+    split(const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+            &              parent,
+          const value_type parent_value);
+
+    /**
+     * Return a vector which contains squared data of the parent cell being
+     * equally divided among the squares of all children.
+     *
+     * @f[
+     *   d_{K_c}^2 = d_{K_p}^2 / n_\text{children}
+     *   \qquad
+     *   \forall K_c \text{ children of } K_p
+     * @f]
+     *
+     * This strategy preserves the $l_2$-norm of the corresponding global data
+     * Vector before and after adaptation.
+     */
+    template <int dim, int spacedim, typename value_type>
+    std::vector<value_type>
+    l2_norm(const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+              &              parent,
+            const value_type parent_value);
+  } // namespace Refinement
+
+  /**
+   * For coarsening, all strategies take the parent cell and a vector of data
+   * that belonged to its former children. They return the value that will be
+   * assigned to the parent cell.
+   *
+   * The ordering of values in the vector for children data corresponds to the
+   * index when calling TriaAccessor::child_index.
+   */
+  namespace Coarsening
+  {
+    /**
+     * Check if data on all children match, and return value of the first child.
+     *
+     * @f[
+     *   d_{K_p} = d_{K_c}
+     *   \qquad
+     *   \forall K_c \text{ children of } K_p
+     * @f]
+     */
+    template <int dim, int spacedim, typename value_type>
+    value_type
+    check_equality(
+      const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+        &                            parent,
+      const std::vector<value_type> &children_values);
+
+    /**
+     * Return sum of data on all children.
+     *
+     * @f[
+     *   d_{K_p} = \sum d_{K_c}
+     *   \qquad
+     *   \forall K_c \text{ children of } K_p
+     * @f]
+     *
+     * This strategy preserves the $l_1$-norm of the corresponding global data
+     * Vector before and after adaptation.
+     */
+    template <int dim, int spacedim, typename value_type>
+    value_type
+    sum(const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+          &                            parent,
+        const std::vector<value_type> &children_values);
+
+    /**
+     * Return $ l_2 $-norm of data on all children.
+     *
+     * @f[
+     *   d_{K_p}^2 = \sum d_{K_c}^2
+     *   \qquad
+     *   \forall K_c \text{ children of } K_p
+     * @f]
+     *
+     * This strategy preserves the $l_2$-norm of the corresponding global data
+     * Vector before and after adaptation.
+     */
+    template <int dim, int spacedim, typename value_type>
+    value_type
+    l2_norm(const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+              &                            parent,
+            const std::vector<value_type> &children_values);
+
+    /**
+     * Return mean value of data on all children.
+     *
+     * @f[
+     *   d_{K_p} = \sum d_{K_c} / n_\text{children}
+     *   \qquad
+     *   \forall K_c \text{ children of } K_p
+     * @f]
+     */
+    template <int dim, int spacedim, typename value_type>
+    value_type
+    mean(const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+           &                            parent,
+         const std::vector<value_type> &children_values);
+
+    /**
+     * Return maximum value of data on all children.
+     *
+     * @f[
+     *   d_{K_p} = \max \left( d_{K_c} \right)
+     *   \qquad
+     *   \forall K_c \text{ children of } K_p
+     * @f]
+     */
+    template <int dim, int spacedim, typename value_type>
+    value_type
+    max(const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+          &                            parent,
+        const std::vector<value_type> &children_values);
+  } // namespace Coarsening
+} // namespace AdaptationStrategies
+
+
+
+/* ---------------- template functions ---------------- */
+
+#ifndef DOXYGEN
+
+namespace AdaptationStrategies
+{
+  namespace Refinement
+  {
+    template <int dim, int spacedim, typename value_type>
+    std::vector<value_type>
+    preserve(const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+               &              parent,
+             const value_type parent_value)
+    {
+      Assert(parent->n_children() > 0, ExcInternalError());
+      return std::vector<value_type>(parent->n_children(), parent_value);
+    }
+
+
+
+    template <int dim, int spacedim, typename value_type>
+    std::vector<value_type>
+    split(const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+            &              parent,
+          const value_type parent_value)
+    {
+      static_assert(std::is_arithmetic<value_type>::value &&
+                      !std::is_same<value_type, bool>::value,
+                    "The provided value_type may not meet the requirements "
+                    "of this function.");
+
+      Assert(parent->n_children() > 0, ExcInternalError());
+      return std::vector<value_type>(parent->n_children(),
+                                     parent_value / parent->n_children());
+    }
+
+
+
+    template <int dim, int spacedim, typename value_type>
+    std::vector<value_type>
+    l2_norm(const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+              &              parent,
+            const value_type parent_value)
+    {
+      static_assert(std::is_arithmetic<value_type>::value &&
+                      !std::is_same<value_type, bool>::value,
+                    "The provided value_type may not meet the requirements "
+                    "of this function.");
+
+      Assert(parent->n_children() > 0, ExcInternalError());
+      return std::vector<value_type>(parent->n_children(),
+                                     parent_value /
+                                       std::sqrt(parent->n_children()));
+    }
+  } // namespace Refinement
+
+
+
+  namespace Coarsening
+  {
+    template <int dim, int spacedim, typename value_type>
+    value_type
+    check_equality(
+      const typename dealii::Triangulation<dim, spacedim>::cell_iterator &,
+      const std::vector<value_type> &children_values)
+    {
+      Assert(!children_values.empty(), ExcInternalError());
+
+      const auto first_child = children_values.cbegin();
+      for (auto other_child = first_child + 1;
+           other_child != children_values.cend();
+           ++other_child)
+        Assert(*first_child == *other_child,
+               ExcMessage(
+                 "Values on cells that will be coarsened are not equal!"));
+
+      return *first_child;
+    }
+
+
+
+    template <int dim, int spacedim, typename value_type>
+    value_type
+    sum(const typename dealii::Triangulation<dim, spacedim>::cell_iterator &,
+        const std::vector<value_type> &children_values)
+    {
+      static_assert(std::is_arithmetic<value_type>::value &&
+                      !std::is_same<value_type, bool>::value,
+                    "The provided value_type may not meet the requirements "
+                    "of this function.");
+
+      Assert(!children_values.empty(), ExcInternalError());
+      return std::accumulate(children_values.cbegin(),
+                             children_values.cend(),
+                             static_cast<value_type>(0));
+    }
+
+
+
+    template <int dim, int spacedim, typename value_type>
+    value_type
+    l2_norm(
+      const typename dealii::Triangulation<dim, spacedim>::cell_iterator &,
+      const std::vector<value_type> &children_values)
+    {
+      static_assert(std::is_arithmetic<value_type>::value &&
+                      !std::is_same<value_type, bool>::value,
+                    "The provided value_type may not meet the requirements "
+                    "of this function.");
+
+      Assert(!children_values.empty(), ExcInternalError());
+      return std::sqrt(std::inner_product(children_values.cbegin(),
+                                          children_values.cend(),
+                                          children_values.cbegin(),
+                                          static_cast<value_type>(0)));
+    }
+
+
+
+    template <int dim, int spacedim, typename value_type>
+    value_type
+    mean(const typename dealii::Triangulation<dim, spacedim>::cell_iterator
+           &                            parent,
+         const std::vector<value_type> &children_values)
+    {
+      return sum<dim, spacedim, value_type>(parent, children_values) /
+             children_values.size();
+    }
+
+
+
+    template <int dim, int spacedim, typename value_type>
+    value_type
+    max(const typename dealii::Triangulation<dim, spacedim>::cell_iterator &,
+        const std::vector<value_type> &children_values)
+    {
+      Assert(!children_values.empty(), ExcInternalError());
+      return *std::max_element(children_values.cbegin(),
+                               children_values.cend());
+    }
+  } // namespace Coarsening
+} // namespace AdaptationStrategies
+
+#endif
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif /* dealii_adaptation_strategies_h */
index cea9761dd9a3156aeca3180f16e50dc704c5b88f..76d215cd72cde62f706c86cf12ce8885ba518918 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2019 by the deal.II authors
+// Copyright (C) 2019 - 2020 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -22,7 +22,7 @@
 
 #include <deal.II/grid/tria.h>
 
-#include <deal.II/numerics/coarsening_strategies.h>
+#include <deal.II/numerics/adaptation_strategies.h>
 
 #include <algorithm>
 #include <functional>
@@ -102,7 +102,7 @@ DEAL_II_NAMESPACE_OPEN
  *   for transfer.
  *
  * @ingroup numerics
- * @author Marc Fehling, 2019
+ * @author Marc Fehling, 2019 - 2020
  */
 template <int dim, int spacedim = dim, typename VectorType = Vector<double>>
 class CellDataTransfer
@@ -120,14 +120,22 @@ public:
    * @param[in] triangulation The triangulation on which all operations will
    *   happen. At the time when this constructor is called, the refinement
    *   in question has not happened yet.
+   * @param[in] refinement_strategy Function deciding how data will be stored on
+   *   refined cells from its parent cell.
    * @param[in] coarsening_strategy Function deciding which data to store on
    *   a cell whose children will get coarsened into.
    */
   CellDataTransfer(
-    const Triangulation<dim, spacedim> &                triangulation,
+    const Triangulation<dim, spacedim> &triangulation,
+    const std::function<std::vector<value_type>(
+      const typename Triangulation<dim, spacedim>::cell_iterator &parent,
+      const value_type parent_value)>   refinement_strategy =
+      &AdaptationStrategies::Refinement::preserve<dim, spacedim, value_type>,
     const std::function<value_type(
-      const std::vector<value_type> &children_indices)> coarsening_strategy =
-      &CoarseningStrategies::check_equality<value_type>);
+      const typename Triangulation<dim, spacedim>::cell_iterator &parent,
+      const std::vector<value_type> &children_values)> coarsening_strategy =
+      &AdaptationStrategies::Coarsening::
+        check_equality<dim, spacedim, value_type>);
 
   /**
    * Prepare the current object for coarsening and refinement.
@@ -157,11 +165,21 @@ private:
                CellDataTransfer<dim, spacedim, VectorType>>
     triangulation;
 
+  /**
+   * Function deciding how data will be stored on refined cells from its parent
+   * cell.
+   */
+  const std::function<std::vector<value_type>(
+    const typename Triangulation<dim, spacedim>::cell_iterator &parent,
+    const value_type                                            parent_value)>
+    refinement_strategy;
+
   /**
    * Function deciding on how to process data from children to be stored on the
    * parent cell.
    */
   const std::function<value_type(
+    const typename Triangulation<dim, spacedim>::cell_iterator &parent,
     const std::vector<value_type> &children_indices)>
     coarsening_strategy;
 
index 211fe82cca396689c69d54aa41fd0b77ab2bd587..91deef81734f85e2a35bfc12050487e2894ee7b0 100644 (file)
@@ -51,10 +51,15 @@ namespace internal
 
 template <int dim, int spacedim, typename VectorType>
 CellDataTransfer<dim, spacedim, VectorType>::CellDataTransfer(
-  const Triangulation<dim, spacedim> &                triangulation,
+  const Triangulation<dim, spacedim> &               triangulation,
+  const std::function<std::vector<value_type>(
+    const typename Triangulation<dim, spacedim>::cell_iterator &parent,
+    const value_type parent_value)>                  refinement_strategy,
   const std::function<value_type(
-    const std::vector<value_type> &children_indices)> coarsening_strategy)
+    const typename Triangulation<dim, spacedim>::cell_iterator &parent,
+    const std::vector<value_type> &children_values)> coarsening_strategy)
   : triangulation(&triangulation, typeid(*this).name())
+  , refinement_strategy(refinement_strategy)
   , coarsening_strategy(coarsening_strategy)
   , n_active_cells_pre(numbers::invalid_unsigned_int)
 {
@@ -104,7 +109,7 @@ CellDataTransfer<dim, spacedim, VectorType>::
                    child_index < parent->n_children();
                    ++child_index)
                 {
-                  const auto sibling = parent->child(child_index);
+                  const auto &sibling = parent->child(child_index);
                   Assert(sibling->is_active() && sibling->coarsen_flag_set(),
                          typename dealii::Triangulation<
                            dim>::ExcInconsistentCoarseningFlags());
@@ -156,19 +161,27 @@ CellDataTransfer<dim, spacedim, VectorType>::unpack(const VectorType &in,
 
   // Transfer data of the parent cell to all of its children that it has been
   // refined to.
+  std::vector<value_type> children_values;
   for (const auto &refined : refined_cells_active_index)
-    for (unsigned int child_index = 0;
-         child_index < refined.first->n_children();
-         ++child_index)
-      {
-        const auto child = refined.first->child(child_index);
-        Assert(child->is_active(), ExcInternalError());
-        out[child->active_cell_index()] = in[refined.second];
-      }
+    {
+      // Decide how to handle the previous data.
+      children_values = refinement_strategy(refined.first, in[refined.second]);
+      Assert(refined.first->n_children() == children_values.size(),
+             ExcInternalError());
+
+      // Set values for all children cells.
+      for (unsigned int child_index = 0;
+           child_index < refined.first->n_children();
+           ++child_index)
+        {
+          const auto &child = refined.first->child(child_index);
+          Assert(child->is_active(), ExcInternalError());
+          out[child->active_cell_index()] = children_values[child_index];
+        }
+    }
 
   // Transfer data from the former children to the cell that they have been
   // coarsened to.
-  std::vector<value_type> children_values;
   for (const auto &coarsened : coarsened_cells_active_index)
     {
       // Get previous values of former children.
@@ -181,7 +194,8 @@ CellDataTransfer<dim, spacedim, VectorType>::unpack(const VectorType &in,
       Assert(values_it == children_values.end(), ExcInternalError());
 
       // Decide how to handle the previous data.
-      const value_type parent_value = coarsening_strategy(children_values);
+      const value_type parent_value =
+        coarsening_strategy(coarsened.first, children_values);
 
       // Set value for the parent cell.
       Assert(coarsened.first->is_active(), ExcInternalError());
diff --git a/include/deal.II/numerics/coarsening_strategies.h b/include/deal.II/numerics/coarsening_strategies.h
deleted file mode 100644 (file)
index 0a29483..0000000
+++ /dev/null
@@ -1,139 +0,0 @@
-// ---------------------------------------------------------------------
-//
-// Copyright (C) 2019 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.
-//
-// ---------------------------------------------------------------------
-
-#ifndef dealii_coarsening_strategies_h
-#define dealii_coarsening_strategies_h
-
-
-#include <deal.II/base/config.h>
-
-#include <deal.II/base/exceptions.h>
-
-#include <algorithm>
-#include <numeric>
-#include <typeinfo>
-#include <vector>
-
-DEAL_II_NAMESPACE_OPEN
-
-/**
- * When data is transferred during coarsening, it is not trivial to decide
- * how to handle data of child cells which will be coarsened. Or in other
- * words, which data should be stored in the corresponding parent cell.
- *
- * In this namespace, we offer a few strategies that cope with this
- * problem. Such strategies can be passed to the CellDataTransfer and
- * parallel::distributed::CellDataTransfer constructors.
- *
- * @author Marc Fehling, 2019
- */
-namespace CoarseningStrategies
-{
-  /**
-   * Check if data on all children match, and return that value.
-   */
-  template <typename value_type>
-  value_type
-  check_equality(const std::vector<value_type> &children_values);
-
-  /**
-   * Return sum of data on all children.
-   */
-  template <typename value_type>
-  value_type
-  sum(const std::vector<value_type> &children_values);
-
-  /**
-   * Return mean value of data on all children.
-   */
-  template <typename value_type>
-  value_type
-  mean(const std::vector<value_type> &children_values);
-
-  /**
-   * Return maximum value of data on all children.
-   */
-  template <typename value_type>
-  value_type
-  max(const std::vector<value_type> &children_values);
-} // namespace CoarseningStrategies
-
-
-
-/* ---------------- template functions ---------------- */
-
-#ifndef DOXYGEN
-
-namespace CoarseningStrategies
-{
-  template <typename value_type>
-  value_type
-  check_equality(const std::vector<value_type> &children_values)
-  {
-    Assert(!children_values.empty(), ExcInternalError());
-
-    const auto first_child = children_values.cbegin();
-    for (auto other_child = first_child + 1;
-         other_child != children_values.cend();
-         ++other_child)
-      Assert(*first_child == *other_child,
-             ExcMessage(
-               "Values on cells that will be coarsened are not equal!"));
-
-    return *first_child;
-  }
-
-
-
-  template <typename value_type>
-  value_type
-  sum(const std::vector<value_type> &children_values)
-  {
-    static_assert(std::is_arithmetic<value_type>::value &&
-                    !std::is_same<value_type, bool>::value,
-                  "The provided value_type may not meet the requirements "
-                  "of this function.");
-
-    Assert(!children_values.empty(), ExcInternalError());
-    return std::accumulate(children_values.cbegin(),
-                           children_values.cend(),
-                           static_cast<value_type>(0));
-  }
-
-
-
-  template <typename value_type>
-  value_type
-  mean(const std::vector<value_type> &children_values)
-  {
-    return sum<value_type>(children_values) / children_values.size();
-  }
-
-
-
-  template <typename value_type>
-  value_type
-  max(const std::vector<value_type> &children_values)
-  {
-    Assert(!children_values.empty(), ExcInternalError());
-    return *std::max_element(children_values.cbegin(), children_values.cend());
-  }
-} // namespace CoarseningStrategies
-
-#endif
-
-DEAL_II_NAMESPACE_CLOSE
-
-#endif /* dealii_coarsening_strategies_h */
index 978c346d872f196d08410f225e001aafe7ecae2f..2ea3381ff010d969eba12b648999fb57ac097680 100644 (file)
@@ -195,7 +195,7 @@ namespace parallel
             for (unsigned int child_index = 0; child_index < cell->n_children();
                  ++child_index)
               {
-                const auto child = cell->child(child_index);
+                const auto &child = cell->child(child_index);
                 Assert(child->is_active() && child->coarsen_flag_set(),
                        typename dealii::Triangulation<
                          dim>::ExcInconsistentCoarseningFlags());
index 0342befbb0b91cd1956e265b0ab27a74723ba2d3..43ed7e213500a5e2e35364328c8748a72a4aaa46 100644 (file)
@@ -336,7 +336,7 @@ namespace parallel
                        child_index < cell->n_children();
                        ++child_index)
                     {
-                      const auto child = cell->child(child_index);
+                      const auto &child = cell->child(child_index);
                       Assert(child->is_active() && child->coarsen_flag_set(),
                              typename dealii::Triangulation<
                                dim>::ExcInconsistentCoarseningFlags());
index be4e010abb9fe11c71092dbdf31e585b98554c82..16046cb174edbe741ce451eb59b0445aa775208a 100644 (file)
@@ -1126,7 +1126,7 @@ namespace internal
                              child_index < parent->n_children();
                              ++child_index)
                           {
-                            const auto sibling = parent->child(child_index);
+                            const auto &sibling = parent->child(child_index);
                             Assert(sibling->is_active() &&
                                      sibling->coarsen_flag_set(),
                                    typename dealii::Triangulation<
@@ -1229,6 +1229,7 @@ namespace internal
         template <int dim, int spacedim>
         static unsigned int
         determine_fe_from_children(
+          const typename Triangulation<dim, spacedim>::cell_iterator &,
           const std::vector<unsigned int> &        children_fe_indices,
           dealii::hp::FECollection<dim, spacedim> &fe_collection)
         {
@@ -2072,10 +2073,17 @@ namespace hp
             CellDataTransfer<dim, spacedim, std::vector<unsigned int>>>(
           *distributed_tria,
           /*transfer_variable_size_data=*/false,
-          [this](const std::vector<unsigned int> &children_fe_indices) {
+          /*refinement_strategy=*/
+          &dealii::AdaptationStrategies::Refinement::
+            preserve<dim, spacedim, unsigned int>,
+          /*coarsening_strategy=*/
+          [this](
+            const typename Triangulation<dim, spacedim>::cell_iterator &parent,
+            const std::vector<unsigned int> &children_fe_indices)
+            -> unsigned int {
             return dealii::internal::hp::DoFHandlerImplementation::
               Implementation::determine_fe_from_children<dim, spacedim>(
-                children_fe_indices, fe_collection);
+                parent, children_fe_indices, fe_collection);
           });
 
         active_fe_index_transfer->cell_data_transfer
@@ -2193,10 +2201,17 @@ namespace hp
             CellDataTransfer<dim, spacedim, std::vector<unsigned int>>>(
           *distributed_tria,
           /*transfer_variable_size_data=*/false,
-          [this](const std::vector<unsigned int> &children_fe_indices) {
+          /*refinement_strategy=*/
+          &dealii::AdaptationStrategies::Refinement::
+            preserve<dim, spacedim, unsigned int>,
+          /*coarsening_strategy=*/
+          [this](
+            const typename Triangulation<dim, spacedim>::cell_iterator &parent,
+            const std::vector<unsigned int> &children_fe_indices)
+            -> unsigned int {
             return dealii::internal::hp::DoFHandlerImplementation::
               Implementation::determine_fe_from_children<dim, spacedim>(
-                children_fe_indices, fe_collection);
+                parent, children_fe_indices, fe_collection);
           });
 
         // If we work on a p::d::Triangulation, we have to transfer all
@@ -2279,10 +2294,17 @@ namespace hp
             CellDataTransfer<dim, spacedim, std::vector<unsigned int>>>(
           *distributed_tria,
           /*transfer_variable_size_data=*/false,
-          [this](const std::vector<unsigned int> &children_fe_indices) {
+          /*refinement_strategy=*/
+          &dealii::AdaptationStrategies::Refinement::
+            preserve<dim, spacedim, unsigned int>,
+          /*coarsening_strategy=*/
+          [this](
+            const typename Triangulation<dim, spacedim>::cell_iterator &parent,
+            const std::vector<unsigned int> &children_fe_indices)
+            -> unsigned int {
             return dealii::internal::hp::DoFHandlerImplementation::
               Implementation::determine_fe_from_children<dim, spacedim>(
-                children_fe_indices, fe_collection);
+                parent, children_fe_indices, fe_collection);
           });
 
         // Unpack active_fe_indices.
index f07c789dc6c267e14cd4bc14155503d45b29f111..b4d839f53df990091cb8e854651221a16d9f5357 100644 (file)
@@ -363,7 +363,7 @@ SolutionTransfer<dim, VectorType, DoFHandlerType>::
           for (unsigned int child_index = 0; child_index < cell->n_children();
                ++child_index)
             {
-              const auto child = cell->child(child_index);
+              const auto &child = cell->child(child_index);
               Assert(child->is_active() && child->coarsen_flag_set(),
                      typename dealii::Triangulation<
                        dim>::ExcInconsistentCoarseningFlags());
index ef7dc0cc36c337d9d949850ae225fb316ea9248c..75d9d74943dbdb8ce7d1a50fdccf0d5873b67da0 100644 (file)
 #include "../tests.h"
 
 
+template <int dim, int spacedim>
 std::vector<int>
-get_data_of_first_child(const std::vector<std::vector<int>> &children_values)
+get_data_of_first_child(
+  const typename dealii::Triangulation<dim, spacedim>::cell_iterator &,
+  const std::vector<std::vector<int>> &children_values)
 {
   return children_values[0];
 }
@@ -86,9 +89,13 @@ test()
   // ----- transfer -----
   parallel::distributed::
     CellDataTransfer<dim, spacedim, std::vector<std::vector<int>>>
-      cell_data_transfer(tria,
-                         /*transfer_variable_size_data=*/true,
-                         /*coarsening_strategy=*/&get_data_of_first_child);
+    cell_data_transfer(
+      tria,
+      /*transfer_variable_size_data=*/true,
+      /*refinement_strategy=*/
+      &dealii::AdaptationStrategies::Refinement::
+        preserve<dim, spacedim, std::vector<int>>,
+      /*coarsening_strategy=*/&get_data_of_first_child<dim, spacedim>);
 
   cell_data_transfer.prepare_for_coarsening_and_refinement(cell_data);
   tria.execute_coarsening_and_refinement();

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.