--- /dev/null
+New: Helper function DoFCellAccessor::dominated_future_fe_on_children()
+to clean up code for hp-coarsening.
+<br>
+(Marc Fehling, 2020/10/05)
*/
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!");
+
/**
* @}
*/
+template <int dimension_, int space_dimension_, bool level_dof_access>
+inline unsigned int
+DoFCellAccessor<dimension_, space_dimension_, level_dof_access>::
+ dominated_future_fe_on_children() const
+{
+ Assert(!this->is_active(),
+ ExcMessage(
+ "You ask for information on children of this cell which is only "
+ "available for active cells. This cell has no children."));
+
+ std::set<unsigned int> future_fe_indices_children;
+ for (const auto &child : this->child_iterators())
+ {
+ Assert(
+ child->is_active(),
+ ExcMessage(
+ "You ask for information on children of this cell which is only "
+ "available for active cells. One of its children is not active."));
+ future_fe_indices_children.insert(child->future_fe_index());
+ }
+ Assert(!future_fe_indices_children.empty(), ExcInternalError());
+
+ const unsigned int future_fe_index =
+ this->dof_handler->get_fe_collection().find_dominated_fe_extended(
+ future_fe_indices_children, /*codim=*/0);
+
+ Assert(future_fe_index != numbers::invalid_unsigned_int,
+ ExcNoDominatedFiniteElementOnChildren());
+
+ return future_fe_index;
+}
+
+
+
template <int dimension_, int space_dimension_, bool level_dof_access>
template <typename number, typename OutputVector>
inline void
*/
DeclException0(ExcNoFiniteElements);
- /**
- * Exception
- *
- * @ingroup Exceptions
- */
- DeclExceptionMsg(
- ExcNoDominatedFiniteElementAmongstChildren,
- "No FiniteElement has been found in your FECollection that is "
- "dominated by all children of a cell you are trying to coarsen!");
-
/**
* @}
*/
break;
case Triangulation<dim, spacedim>::CELL_COARSEN:
- {
- std::set<unsigned int> fe_indices_children;
- for (unsigned int child_index = 0; child_index < cell->n_children();
- ++child_index)
- {
- const auto &child = cell->child(child_index);
- Assert(child->is_active() && child->coarsen_flag_set(),
- typename dealii::Triangulation<
- dim>::ExcInconsistentCoarseningFlags());
-
- fe_indices_children.insert(child->future_fe_index());
- }
- Assert(!fe_indices_children.empty(), ExcInternalError());
-
- fe_index =
- dof_handler.get_fe_collection().find_dominated_fe_extended(
- fe_indices_children, /*codim=*/0);
-
- Assert(fe_index != numbers::invalid_unsigned_int,
- typename dealii::hp::FECollection<
- dim>::ExcNoDominatedFiniteElementAmongstChildren());
- }
+#ifdef DEBUG
+ for (const auto &child : cell->child_iterators())
+ Assert(child->is_active() && child->coarsen_flag_set(),
+ typename dealii::Triangulation<
+ dim>::ExcInconsistentCoarseningFlags());
+#endif
+
+ fe_index = cell->dominated_future_fe_on_children();
break;
default:
// parent cell after h-adaptation analogously to
// dealii::internal::hp::DoFHandlerImplementation::Implementation::
// collect_fe_indices_on_cells_to_be_refined()
- std::set<unsigned int> fe_indices_children;
- for (unsigned int child_index = 0;
- child_index < cell->n_children();
- ++child_index)
- {
- const auto child = cell->child(child_index);
- Assert(child->is_active() && child->coarsen_flag_set(),
- typename dealii::Triangulation<
- dim>::ExcInconsistentCoarseningFlags());
-
- fe_indices_children.insert(child->future_fe_index());
- }
+# ifdef DEBUG
+ for (const auto &child : cell->child_iterators())
+ Assert(child->is_active() && child->coarsen_flag_set(),
+ typename dealii::Triangulation<
+ dim>::ExcInconsistentCoarseningFlags());
+# endif
const unsigned int future_fe_index =
- dof_handler->get_fe_collection().find_dominated_fe_extended(
- fe_indices_children, /*codim=*/0);
+ cell->dominated_future_fe_on_children();
const unsigned int future_fe_degree =
dof_handler->get_fe_collection()[future_fe_index].degree;
// In case of coarsening, we need to find a suitable fe index
// for the parent cell. We choose the 'least dominant fe'
// on all children from the associated FECollection.
- std::set<unsigned int> fe_indices_children;
- for (unsigned int child_index = 0;
- child_index < cell->n_children();
- ++child_index)
- {
- const auto &child = cell->child(child_index);
- Assert(child->is_active() && child->coarsen_flag_set(),
- typename dealii::Triangulation<
- dim>::ExcInconsistentCoarseningFlags());
-
- fe_indices_children.insert(child->future_fe_index());
- }
- Assert(!fe_indices_children.empty(), ExcInternalError());
-
- fe_index =
- dof_handler->get_fe_collection().find_dominated_fe_extended(
- fe_indices_children, /*codim=*/0);
-
- Assert(fe_index != numbers::invalid_unsigned_int,
- typename dealii::hp::FECollection<
- dim>::ExcNoDominatedFiniteElementAmongstChildren());
-
+# ifdef DEBUG
+ for (const auto &child : cell->child_iterators())
+ Assert(child->is_active() && child->coarsen_flag_set(),
+ typename dealii::Triangulation<
+ dim>::ExcInconsistentCoarseningFlags());
+# endif
+
+ fe_index = cell->dominated_future_fe_on_children();
break;
}
// based on the 'least dominant finite element' of its
// children. Consider the childrens' hypothetical future
// index when they have been flagged for p-refinement.
- std::set<unsigned int> fe_indices_children;
- for (unsigned int child_index = 0;
- child_index < parent->n_children();
- ++child_index)
- {
- const auto sibling = parent->child(child_index);
- Assert(sibling->is_active() &&
- sibling->coarsen_flag_set(),
- typename dealii::Triangulation<
- dim>::ExcInconsistentCoarseningFlags());
-
- fe_indices_children.insert(
- sibling->future_fe_index());
- }
- Assert(!fe_indices_children.empty(),
- ExcInternalError());
+#ifdef DEBUG
+ for (const auto &child : parent->child_iterators())
+ Assert(child->is_active() &&
+ child->coarsen_flag_set(),
+ typename dealii::Triangulation<
+ dim>::ExcInconsistentCoarseningFlags());
+#endif
const unsigned int fe_index =
- dof_handler.fe_collection.find_dominated_fe_extended(
- fe_indices_children, /*codim=*/0);
-
- Assert(fe_index != numbers::invalid_unsigned_int,
- typename dealii::hp::FECollection<dim>::
- ExcNoDominatedFiniteElementAmongstChildren());
+ parent->dominated_future_fe_on_children();
fe_transfer->coarsened_cells_fe_index.insert(
{parent, fe_index});
/*codim=*/0);
Assert(dominated_fe_index != numbers::invalid_unsigned_int,
- typename dealii::hp::FECollection<
- dim>::ExcNoDominatedFiniteElementAmongstChildren());
+ (typename dealii::DoFCellAccessor<dim, spacedim, false>::
+ ExcNoDominatedFiniteElementOnChildren()));
return dominated_fe_index;
}
if (future_fe_indices_on_coarsened_cells.find(parent) ==
future_fe_indices_on_coarsened_cells.end())
{
- std::set<unsigned int> fe_indices_children;
- for (unsigned int child_index = 0;
- child_index < parent->n_children();
- ++child_index)
- {
- const auto &child = parent->child(child_index);
- Assert(child->is_active() && child->coarsen_flag_set(),
- typename dealii::Triangulation<
- dim>::ExcInconsistentCoarseningFlags());
-
- fe_indices_children.insert(child->future_fe_index());
- }
- Assert(!fe_indices_children.empty(), ExcInternalError());
+#ifdef DEBUG
+ for (const auto &child : parent->child_iterators())
+ Assert(child->is_active() && child->coarsen_flag_set(),
+ typename dealii::Triangulation<
+ dim>::ExcInconsistentCoarseningFlags());
+#endif
parent_future_fe_index =
- dof_handler.get_fe_collection()
- .find_dominated_fe_extended(fe_indices_children,
- /*codim=*/0);
-
- Assert(
- parent_future_fe_index != numbers::invalid_unsigned_int,
- typename dealii::hp::FECollection<
- dim>::ExcNoDominatedFiniteElementAmongstChildren());
+ parent->dominated_future_fe_on_children();
future_fe_indices_on_coarsened_cells.insert(
{parent, parent_future_fe_index});
fe_indices_children, /*codim=*/0);
Assert(target_fe_index != numbers::invalid_unsigned_int,
- typename dealii::hp::FECollection<
- dim>::ExcNoDominatedFiniteElementAmongstChildren());
+ (typename dealii::DoFCellAccessor<
+ dim,
+ DoFHandlerType::space_dimension,
+ false>::ExcNoDominatedFiniteElementOnChildren()));
const unsigned int dofs_per_cell =
dof_handler->get_fe(target_fe_index).n_dofs_per_cell();
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Verify functionality of
+// DoFCellAccessor:dominated_future_fe_on_children()
+
+#include <deal.II/base/logstream.h>
+
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/hp/fe_collection.h>
+
+#include "../tests.h"
+
+
+
+template <int dim>
+void
+test()
+{
+ Triangulation<dim> tria;
+ GridGenerator::hyper_cube(tria);
+ tria.refine_global(1);
+
+ hp::FECollection<dim> fes;
+ fes.push_back(FE_Q<dim>(1));
+ fes.push_back(FE_Q<dim>(2));
+
+ DoFHandler<dim> dofh(tria, /*hp_capability_enabled=*/true);
+ dofh.set_fe(fes);
+ dofh.begin_active()->set_active_fe_index(1);
+
+ const auto & parent = dofh.begin(/*level=*/0);
+ const unsigned int parent_future_fe =
+ parent->dominated_future_fe_on_children();
+
+ deallog << parent_future_fe << std::endl;
+}
+
+
+int
+main()
+{
+ initlog();
+
+ deallog.push("1d");
+ test<1>();
+ deallog.pop();
+ deallog.push("2d");
+ test<2>();
+ deallog.pop();
+ deallog.push("3d");
+ test<3>();
+ deallog.pop();
+}
--- /dev/null
+
+DEAL:1d::1
+DEAL:2d::1
+DEAL:3d::1