From: Marc Fehling Date: Fri, 24 Apr 2020 21:23:08 +0000 (+0200) Subject: Refactored AdaptationStrategies namespace. X-Git-Tag: v9.2.0-rc1~15^2~3 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=ca86e5f7f6d122ff52de6312c2194d180c691abc;p=dealii.git Refactored AdaptationStrategies namespace. --- diff --git a/include/deal.II/distributed/cell_data_transfer.h b/include/deal.II/distributed/cell_data_transfer.h index aa404813c2..2e3ec97b34 100644 --- a/include/deal.II/distributed/cell_data_transfer.h +++ b/include/deal.II/distributed/cell_data_transfer.h @@ -22,7 +22,7 @@ #include -#include +#include #include @@ -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 > 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:: @@ -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 - & triangulation, - const bool transfer_variable_size_data = false, + & triangulation, + const bool transfer_variable_size_data = false, + const std::function( + const typename dealii::Triangulation::cell_iterator + & parent, + const value_type parent_value)> refinement_strategy = + &dealii::AdaptationStrategies::Refinement:: + preserve, const std::function::cell_iterator + & parent, const std::vector &children_values)> coarsening_strategy = - &dealii::CoarseningStrategies::check_equality); + &dealii::AdaptationStrategies::Coarsening:: + check_equality); /** - * @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( + const typename Triangulation::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::cell_iterator &parent, const std::vector &children_values)> coarsening_strategy; diff --git a/include/deal.II/distributed/cell_data_transfer.templates.h b/include/deal.II/distributed/cell_data_transfer.templates.h index 8c60012c1b..260d97a8e5 100644 --- a/include/deal.II/distributed/cell_data_transfer.templates.h +++ b/include/deal.II/distributed/cell_data_transfer.templates.h @@ -73,11 +73,18 @@ namespace parallel template CellDataTransfer::CellDataTransfer( const parallel::distributed::Triangulation &triangulation, - const bool transfer_variable_size_data, + const bool transfer_variable_size_data, + const std::function( + const typename dealii::Triangulation::cell_iterator + & parent, + const value_type parent_value)> refinement_strategy, const std::function::cell_iterator + & parent, const std::vector &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) , 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 &)> &>( - this->coarsening_strategy) = - &dealii::CoarseningStrategies::check_equality; + const_cast::cell_iterator &, + const std::vector &)> &>(this->coarsening_strategy) = + &dealii::AdaptationStrategies::Coarsening:: + check_equality; else if (*old_strategy == CoarseningStrategies::sum) - const_cast< - std::function &)> &>( - this->coarsening_strategy) = - &dealii::CoarseningStrategies::sum; + const_cast::cell_iterator &, + const std::vector &)> &>(this->coarsening_strategy) = + &dealii::AdaptationStrategies::Coarsening:: + sum; else if (*old_strategy == CoarseningStrategies::mean) - const_cast< - std::function &)> &>( - this->coarsening_strategy) = - &dealii::CoarseningStrategies::mean; + const_cast::cell_iterator &, + const std::vector &)> &>(this->coarsening_strategy) = + &dealii::AdaptationStrategies::Coarsening:: + mean; 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::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 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 index 0000000000..451789150c --- /dev/null +++ b/include/deal.II/numerics/adaptation_strategies.h @@ -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 + +#include + +#include + +#include +#include +#include +#include + +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 + std::vector + preserve(const typename dealii::Triangulation::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 + std::vector + split(const typename dealii::Triangulation::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 + std::vector + l2_norm(const typename dealii::Triangulation::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 + value_type + check_equality( + const typename dealii::Triangulation::cell_iterator + & parent, + const std::vector &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 + value_type + sum(const typename dealii::Triangulation::cell_iterator + & parent, + const std::vector &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 + value_type + l2_norm(const typename dealii::Triangulation::cell_iterator + & parent, + const std::vector &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 + value_type + mean(const typename dealii::Triangulation::cell_iterator + & parent, + const std::vector &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 + value_type + max(const typename dealii::Triangulation::cell_iterator + & parent, + const std::vector &children_values); + } // namespace Coarsening +} // namespace AdaptationStrategies + + + +/* ---------------- template functions ---------------- */ + +#ifndef DOXYGEN + +namespace AdaptationStrategies +{ + namespace Refinement + { + template + std::vector + preserve(const typename dealii::Triangulation::cell_iterator + & parent, + const value_type parent_value) + { + Assert(parent->n_children() > 0, ExcInternalError()); + return std::vector(parent->n_children(), parent_value); + } + + + + template + std::vector + split(const typename dealii::Triangulation::cell_iterator + & parent, + const value_type parent_value) + { + static_assert(std::is_arithmetic::value && + !std::is_same::value, + "The provided value_type may not meet the requirements " + "of this function."); + + Assert(parent->n_children() > 0, ExcInternalError()); + return std::vector(parent->n_children(), + parent_value / parent->n_children()); + } + + + + template + std::vector + l2_norm(const typename dealii::Triangulation::cell_iterator + & parent, + const value_type parent_value) + { + static_assert(std::is_arithmetic::value && + !std::is_same::value, + "The provided value_type may not meet the requirements " + "of this function."); + + Assert(parent->n_children() > 0, ExcInternalError()); + return std::vector(parent->n_children(), + parent_value / + std::sqrt(parent->n_children())); + } + } // namespace Refinement + + + + namespace Coarsening + { + template + value_type + check_equality( + const typename dealii::Triangulation::cell_iterator &, + const std::vector &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 + value_type + sum(const typename dealii::Triangulation::cell_iterator &, + const std::vector &children_values) + { + static_assert(std::is_arithmetic::value && + !std::is_same::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(0)); + } + + + + template + value_type + l2_norm( + const typename dealii::Triangulation::cell_iterator &, + const std::vector &children_values) + { + static_assert(std::is_arithmetic::value && + !std::is_same::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(0))); + } + + + + template + value_type + mean(const typename dealii::Triangulation::cell_iterator + & parent, + const std::vector &children_values) + { + return sum(parent, children_values) / + children_values.size(); + } + + + + template + value_type + max(const typename dealii::Triangulation::cell_iterator &, + const std::vector &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 */ diff --git a/include/deal.II/numerics/cell_data_transfer.h b/include/deal.II/numerics/cell_data_transfer.h index cea9761dd9..76d215cd72 100644 --- a/include/deal.II/numerics/cell_data_transfer.h +++ b/include/deal.II/numerics/cell_data_transfer.h @@ -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 -#include +#include #include #include @@ -102,7 +102,7 @@ DEAL_II_NAMESPACE_OPEN * for transfer. * * @ingroup numerics - * @author Marc Fehling, 2019 + * @author Marc Fehling, 2019 - 2020 */ template > 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 & triangulation, + const Triangulation &triangulation, + const std::function( + const typename Triangulation::cell_iterator &parent, + const value_type parent_value)> refinement_strategy = + &AdaptationStrategies::Refinement::preserve, const std::function &children_indices)> coarsening_strategy = - &CoarseningStrategies::check_equality); + const typename Triangulation::cell_iterator &parent, + const std::vector &children_values)> coarsening_strategy = + &AdaptationStrategies::Coarsening:: + check_equality); /** * Prepare the current object for coarsening and refinement. @@ -157,11 +165,21 @@ private: CellDataTransfer> triangulation; + /** + * Function deciding how data will be stored on refined cells from its parent + * cell. + */ + const std::function( + const typename Triangulation::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::cell_iterator &parent, const std::vector &children_indices)> coarsening_strategy; diff --git a/include/deal.II/numerics/cell_data_transfer.templates.h b/include/deal.II/numerics/cell_data_transfer.templates.h index 211fe82cca..91deef8173 100644 --- a/include/deal.II/numerics/cell_data_transfer.templates.h +++ b/include/deal.II/numerics/cell_data_transfer.templates.h @@ -51,10 +51,15 @@ namespace internal template CellDataTransfer::CellDataTransfer( - const Triangulation & triangulation, + const Triangulation & triangulation, + const std::function( + const typename Triangulation::cell_iterator &parent, + const value_type parent_value)> refinement_strategy, const std::function &children_indices)> coarsening_strategy) + const typename Triangulation::cell_iterator &parent, + const std::vector &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:: 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::unpack(const VectorType &in, // Transfer data of the parent cell to all of its children that it has been // refined to. + std::vector 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 children_values; for (const auto &coarsened : coarsened_cells_active_index) { // Get previous values of former children. @@ -181,7 +194,8 @@ CellDataTransfer::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 index 0a29483962..0000000000 --- a/include/deal.II/numerics/coarsening_strategies.h +++ /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 - -#include - -#include -#include -#include -#include - -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 - value_type - check_equality(const std::vector &children_values); - - /** - * Return sum of data on all children. - */ - template - value_type - sum(const std::vector &children_values); - - /** - * Return mean value of data on all children. - */ - template - value_type - mean(const std::vector &children_values); - - /** - * Return maximum value of data on all children. - */ - template - value_type - max(const std::vector &children_values); -} // namespace CoarseningStrategies - - - -/* ---------------- template functions ---------------- */ - -#ifndef DOXYGEN - -namespace CoarseningStrategies -{ - template - value_type - check_equality(const std::vector &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 - value_type - sum(const std::vector &children_values) - { - static_assert(std::is_arithmetic::value && - !std::is_same::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(0)); - } - - - - template - value_type - mean(const std::vector &children_values) - { - return sum(children_values) / children_values.size(); - } - - - - template - value_type - max(const std::vector &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 */ diff --git a/source/distributed/cell_weights.cc b/source/distributed/cell_weights.cc index 978c346d87..2ea3381ff0 100644 --- a/source/distributed/cell_weights.cc +++ b/source/distributed/cell_weights.cc @@ -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()); diff --git a/source/distributed/solution_transfer.cc b/source/distributed/solution_transfer.cc index 0342befbb0..43ed7e2135 100644 --- a/source/distributed/solution_transfer.cc +++ b/source/distributed/solution_transfer.cc @@ -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()); diff --git a/source/hp/dof_handler.cc b/source/hp/dof_handler.cc index be4e010abb..16046cb174 100644 --- a/source/hp/dof_handler.cc +++ b/source/hp/dof_handler.cc @@ -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 static unsigned int determine_fe_from_children( + const typename Triangulation::cell_iterator &, const std::vector & children_fe_indices, dealii::hp::FECollection &fe_collection) { @@ -2072,10 +2073,17 @@ namespace hp CellDataTransfer>>( *distributed_tria, /*transfer_variable_size_data=*/false, - [this](const std::vector &children_fe_indices) { + /*refinement_strategy=*/ + &dealii::AdaptationStrategies::Refinement:: + preserve, + /*coarsening_strategy=*/ + [this]( + const typename Triangulation::cell_iterator &parent, + const std::vector &children_fe_indices) + -> unsigned int { return dealii::internal::hp::DoFHandlerImplementation:: Implementation::determine_fe_from_children( - 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>>( *distributed_tria, /*transfer_variable_size_data=*/false, - [this](const std::vector &children_fe_indices) { + /*refinement_strategy=*/ + &dealii::AdaptationStrategies::Refinement:: + preserve, + /*coarsening_strategy=*/ + [this]( + const typename Triangulation::cell_iterator &parent, + const std::vector &children_fe_indices) + -> unsigned int { return dealii::internal::hp::DoFHandlerImplementation:: Implementation::determine_fe_from_children( - 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>>( *distributed_tria, /*transfer_variable_size_data=*/false, - [this](const std::vector &children_fe_indices) { + /*refinement_strategy=*/ + &dealii::AdaptationStrategies::Refinement:: + preserve, + /*coarsening_strategy=*/ + [this]( + const typename Triangulation::cell_iterator &parent, + const std::vector &children_fe_indices) + -> unsigned int { return dealii::internal::hp::DoFHandlerImplementation:: Implementation::determine_fe_from_children( - children_fe_indices, fe_collection); + parent, children_fe_indices, fe_collection); }); // Unpack active_fe_indices. diff --git a/source/numerics/solution_transfer.cc b/source/numerics/solution_transfer.cc index f07c789dc6..b4d839f53d 100644 --- a/source/numerics/solution_transfer.cc +++ b/source/numerics/solution_transfer.cc @@ -363,7 +363,7 @@ SolutionTransfer:: 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()); diff --git a/tests/mpi/cell_data_transfer_03.cc b/tests/mpi/cell_data_transfer_03.cc index ef7dc0cc36..75d9d74943 100644 --- a/tests/mpi/cell_data_transfer_03.cc +++ b/tests/mpi/cell_data_transfer_03.cc @@ -29,8 +29,11 @@ #include "../tests.h" +template std::vector -get_data_of_first_child(const std::vector> &children_values) +get_data_of_first_child( + const typename dealii::Triangulation::cell_iterator &, + const std::vector> &children_values) { return children_values[0]; } @@ -86,9 +89,13 @@ test() // ----- transfer ----- parallel::distributed:: CellDataTransfer>> - 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>, + /*coarsening_strategy=*/&get_data_of_first_child); cell_data_transfer.prepare_for_coarsening_and_refinement(cell_data); tria.execute_coarsening_and_refinement();