From: Marc Fehling Date: Fri, 26 Mar 2021 04:09:31 +0000 (-0600) Subject: Moved `dominated_future_fe_on_children` from class member to internal namespace. X-Git-Tag: v9.3.0-rc1~271^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=2dd3e086eb5c79c76d39a46b7ecc48dc6f4a1e2f;p=dealii.git Moved `dominated_future_fe_on_children` from class member to internal namespace. --- diff --git a/doc/news/changes/minor/20201005Fehling b/doc/news/changes/minor/20201005Fehling deleted file mode 100644 index a9d4679274..0000000000 --- a/doc/news/changes/minor/20201005Fehling +++ /dev/null @@ -1,4 +0,0 @@ -New: Helper function DoFCellAccessor::dominated_future_fe_on_children() -to clean up code for hp-coarsening. -
-(Marc Fehling, 2020/10/05) diff --git a/include/deal.II/distributed/cell_weights.h b/include/deal.II/distributed/cell_weights.h index 1693088627..e8c1649178 100644 --- a/include/deal.II/distributed/cell_weights.h +++ b/include/deal.II/distributed/cell_weights.h @@ -20,7 +20,7 @@ #include -#include +#include DEAL_II_NAMESPACE_OPEN diff --git a/include/deal.II/dofs/dof_accessor.h b/include/deal.II/dofs/dof_accessor.h index c3b39e4212..2adb5e7a43 100644 --- a/include/deal.II/dofs/dof_accessor.h +++ b/include/deal.II/dofs/dof_accessor.h @@ -2060,48 +2060,6 @@ public: */ void clear_future_fe_index() const; - - /** - * Return the index of the finite element from the entire hp::FECollection - * that is dominated by those assigned as future finite elements to the - * children of this cell. - * - * We find the corresponding finite element among the future finite elements - * on the children of this cell. If none of them qualify, we extend our search - * on the whole hp::FECollection, which is the element that describes the - * smallest finite element space that includes all future finite elements - * assigned to the children. If the function is not able to find a finite - * element at all, an assertion will be triggered. - * - * In this way, we determine the finite element of the parent cell in case of - * h-coarsening in the hp-context. - * - * @note This function can only be called on direct parent cells, i.e., - * non-active cells whose children are all active. - */ - unsigned int - dominated_future_fe_on_children() const; - /** - * @} - */ - - /** - * @name Exceptions - */ - /** - * @{ - */ - - /** - * Exception - * - * @ingroup Exceptions - */ - DeclExceptionMsg( - ExcNoDominatedFiniteElementOnChildren, - "No FiniteElement has been found in your FECollection that is " - "dominated by all children of a cell you are trying to coarsen!"); - /** * @} */ diff --git a/include/deal.II/dofs/dof_accessor.templates.h b/include/deal.II/dofs/dof_accessor.templates.h index 265f83249e..3ed210610f 100644 --- a/include/deal.II/dofs/dof_accessor.templates.h +++ b/include/deal.II/dofs/dof_accessor.templates.h @@ -2934,45 +2934,6 @@ DoFCellAccessor:: -template -inline unsigned int -DoFCellAccessor:: - dominated_future_fe_on_children() const -{ - Assert(!this->is_active(), - ExcMessage( - "You ask for information on children of this cell which is only " - "available for active cells. This cell has no children.")); - - std::set future_fe_indices_children; - for (const auto &child : this->child_iterators()) - { - Assert( - child->is_active(), - ExcMessage( - "You ask for information on children of this cell which is only " - "available for active cells. One of its children is not active.")); - Assert(child->is_locally_owned(), - ExcMessage( - "You ask for information on children of this cell which is only " - "available for locally owned cells. One of its children is not " - "locally owned.")); - future_fe_indices_children.insert(child->future_fe_index()); - } - Assert(!future_fe_indices_children.empty(), ExcInternalError()); - - const unsigned int future_fe_index = - this->dof_handler->fe_collection.find_dominated_fe_extended( - future_fe_indices_children, /*codim=*/0); - - Assert(future_fe_index != numbers::invalid_unsigned_int, - ExcNoDominatedFiniteElementOnChildren()); - - return future_fe_index; -} - - - template template inline void diff --git a/include/deal.II/dofs/dof_handler.h b/include/deal.II/dofs/dof_handler.h index c4c15ab2e8..6353304739 100644 --- a/include/deal.II/dofs/dof_handler.h +++ b/include/deal.II/dofs/dof_handler.h @@ -1817,6 +1817,46 @@ private: #endif }; +namespace internal +{ + namespace hp + { + namespace DoFHandlerImplementation + { + /** + * Return the index of the finite element from the entire hp::FECollection + * that is dominated by those assigned as future finite elements to the + * children of @p parent. + * + * We find the corresponding finite element among the future finite + * elements on the children of this cell. If none of them qualify, we + * extend our search on the whole hp::FECollection, which is the element + * that describes the smallest finite element space that includes all + * future finite elements assigned to the children. If the function is not + * able to find a finite element at all, an assertion will be triggered. + * + * In this way, we determine the finite element of the parent cell in case + * of h-coarsening in the hp-context. + * + * @note This function can only be called on direct parent cells, i.e., + * non-active cells whose children are all active. + */ + template + unsigned int + dominated_future_fe_on_children( + const typename DoFHandler::cell_iterator &parent); + + /** + * Exception + */ + DeclExceptionMsg( + ExcNoDominatedFiniteElementOnChildren, + "No FiniteElement has been found in your FECollection that is " + "dominated by all children of a cell you are trying to coarsen!"); + } // namespace DoFHandlerImplementation + } // namespace hp +} // namespace internal + #ifndef DOXYGEN /* ----------------------- Inline functions ---------------------------------- diff --git a/source/distributed/cell_weights.cc b/source/distributed/cell_weights.cc index 0688d1db98..d493728065 100644 --- a/source/distributed/cell_weights.cc +++ b/source/distributed/cell_weights.cc @@ -200,7 +200,8 @@ namespace parallel dim>::ExcInconsistentCoarseningFlags()); #endif - fe_index = cell->dominated_future_fe_on_children(); + fe_index = dealii::internal::hp::DoFHandlerImplementation:: + dominated_future_fe_on_children(cell); break; default: diff --git a/source/distributed/error_predictor.cc b/source/distributed/error_predictor.cc index 1fb13c38f9..df2f4eb39c 100644 --- a/source/distributed/error_predictor.cc +++ b/source/distributed/error_predictor.cc @@ -240,7 +240,8 @@ namespace parallel # endif const unsigned int future_fe_index = - cell->dominated_future_fe_on_children(); + dealii::internal::hp::DoFHandlerImplementation:: + dominated_future_fe_on_children(cell); const unsigned int future_fe_degree = dof_handler->get_fe_collection()[future_fe_index].degree; diff --git a/source/distributed/solution_transfer.cc b/source/distributed/solution_transfer.cc index 50104891cd..037d0d70e2 100644 --- a/source/distributed/solution_transfer.cc +++ b/source/distributed/solution_transfer.cc @@ -324,7 +324,10 @@ namespace parallel dim>::ExcInconsistentCoarseningFlags()); # endif - fe_index = cell->dominated_future_fe_on_children(); + fe_index = dealii::internal::hp::DoFHandlerImplementation:: + dominated_future_fe_on_children< + dim, + DoFHandlerType::space_dimension>(cell); break; } diff --git a/source/dofs/dof_handler.cc b/source/dofs/dof_handler.cc index 0fb78f4f63..9cffc4ccc9 100644 --- a/source/dofs/dof_handler.cc +++ b/source/dofs/dof_handler.cc @@ -1000,6 +1000,8 @@ namespace internal }; } // namespace DoFHandlerImplementation + + namespace hp { namespace DoFHandlerImplementation @@ -1872,7 +1874,9 @@ namespace internal #endif const unsigned int fe_index = - parent->dominated_future_fe_on_children(); + dealii::internal::hp::DoFHandlerImplementation:: + dominated_future_fe_on_children( + parent); fe_transfer->coarsened_cells_fe_index.insert( {parent, fe_index}); @@ -1969,12 +1973,74 @@ namespace internal /*codim=*/0); Assert(dominated_fe_index != numbers::invalid_unsigned_int, - (typename dealii::DoFCellAccessor:: - ExcNoDominatedFiniteElementOnChildren())); + ExcNoDominatedFiniteElementOnChildren()); return dominated_fe_index; } + + + /** + * Return the index of the finite element from the entire + * hp::FECollection that is dominated by those assigned as future finite + * elements to the + * children of @p parent. + * + * See documentation in the header file for more information. + */ + template + static unsigned int + dominated_future_fe_on_children( + const typename DoFHandler::cell_iterator &parent) + { + Assert( + !parent->is_active(), + ExcMessage( + "You ask for information on children of this cell which is only " + "available for active cells. This cell has no children.")); + + std::set future_fe_indices_children; + for (const auto &child : parent->child_iterators()) + { + Assert( + child->is_active(), + ExcMessage( + "You ask for information on children of this cell which is only " + "available for active cells. One of its children is not active.")); + Assert( + child->is_locally_owned(), + ExcMessage( + "You ask for information on children of this cell which is only " + "available for locally owned cells. One of its children is not " + "locally owned.")); + future_fe_indices_children.insert(child->future_fe_index()); + } + Assert(!future_fe_indices_children.empty(), ExcInternalError()); + + const unsigned int future_fe_index = + parent->get_dof_handler().fe_collection.find_dominated_fe_extended( + future_fe_indices_children, + /*codim=*/0); + + Assert(future_fe_index != numbers::invalid_unsigned_int, + ExcNoDominatedFiniteElementOnChildren()); + + return future_fe_index; + } }; + + + + /** + * Public wrapper for the above function. + */ + template + unsigned int + dominated_future_fe_on_children( + const typename DoFHandler::cell_iterator &parent) + { + return Implementation::dominated_future_fe_on_children( + parent); + } } // namespace DoFHandlerImplementation } // namespace hp } // namespace internal diff --git a/source/dofs/dof_handler.inst.in b/source/dofs/dof_handler.inst.in index 9d57271784..b6fb6cb840 100644 --- a/source/dofs/dof_handler.inst.in +++ b/source/dofs/dof_handler.inst.in @@ -1,6 +1,6 @@ // --------------------------------------------------------------------- // -// Copyright (C) 1998 - 2020 by the deal.II authors +// Copyright (C) 1998 - 2021 by the deal.II authors // // This file is part of the deal.II library. // @@ -22,6 +22,12 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : DIMENSIONS) template std::string internal::policy_to_string( const dealii::internal::DoFHandlerImplementation::Policy:: PolicyBase &policy); + + template unsigned int internal::hp::DoFHandlerImplementation:: + dominated_future_fe_on_children( + const typename DoFHandler::cell_iterator &); #endif } diff --git a/source/hp/refinement.cc b/source/hp/refinement.cc index def4bb4a8c..5eebda3355 100644 --- a/source/hp/refinement.cc +++ b/source/hp/refinement.cc @@ -573,7 +573,8 @@ namespace hp #endif parent_future_fe_index = - parent->dominated_future_fe_on_children(); + dealii::internal::hp::DoFHandlerImplementation:: + dominated_future_fe_on_children(parent); future_fe_indices_on_coarsened_cells.insert( {parent, parent_future_fe_index}); diff --git a/source/numerics/solution_transfer.cc b/source/numerics/solution_transfer.cc index c67ae4dc41..49cd4a72b6 100644 --- a/source/numerics/solution_transfer.cc +++ b/source/numerics/solution_transfer.cc @@ -397,10 +397,8 @@ SolutionTransfer:: fe_indices_children, /*codim=*/0); Assert(target_fe_index != numbers::invalid_unsigned_int, - (typename dealii::DoFCellAccessor< - dim, - DoFHandlerType::space_dimension, - false>::ExcNoDominatedFiniteElementOnChildren())); + internal::hp::DoFHandlerImplementation:: + ExcNoDominatedFiniteElementOnChildren()); const unsigned int dofs_per_cell = dof_handler->get_fe(target_fe_index).n_dofs_per_cell(); diff --git a/tests/grid/accessor_04.cc b/tests/hp/dominated_future_fe_on_children_01.cc similarity index 90% rename from tests/grid/accessor_04.cc rename to tests/hp/dominated_future_fe_on_children_01.cc index 6593fd965d..0040a04a62 100644 --- a/tests/grid/accessor_04.cc +++ b/tests/hp/dominated_future_fe_on_children_01.cc @@ -16,7 +16,7 @@ // Verify functionality of -// DoFCellAccessor:dominated_future_fe_on_children() +// internal::hp::DoFHandlerImplementation::dominated_future_fe_on_children() #include @@ -52,7 +52,8 @@ test() const auto & parent = dofh.begin(/*level=*/0); const unsigned int parent_future_fe = - parent->dominated_future_fe_on_children(); + internal::hp::DoFHandlerImplementation::dominated_future_fe_on_children< + dim>(parent); deallog << parent_future_fe << std::endl; } diff --git a/tests/grid/accessor_04.output b/tests/hp/dominated_future_fe_on_children_01.output similarity index 100% rename from tests/grid/accessor_04.output rename to tests/hp/dominated_future_fe_on_children_01.output