#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>
* 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
* 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::
}
/**
- * @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>::
}
/**
- * @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::
* 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.
*/
*/
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;
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)
{}
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)(
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,
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());
(**it_input)[child->active_cell_index()];
}
- *it_output = coarsening_strategy(children_values);
+ *it_output = coarsening_strategy(cell, children_values);
break;
}
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());
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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 */
// ---------------------------------------------------------------------
//
-// 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.
//
#include <deal.II/grid/tria.h>
-#include <deal.II/numerics/coarsening_strategies.h>
+#include <deal.II/numerics/adaptation_strategies.h>
#include <algorithm>
#include <functional>
* 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
* @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.
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;
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)
{
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());
// 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.
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());
+++ /dev/null
-// ---------------------------------------------------------------------
-//
-// 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 */
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());
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());
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<
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)
{
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
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
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.
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());
#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];
}
// ----- 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();