From: Marc Fehling Date: Thu, 29 Feb 2024 02:07:55 +0000 (-0700) Subject: Fix limit_p_level_difference for p-coarsening. X-Git-Tag: v9.6.0-rc1~511^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F16705%2Fhead;p=dealii.git Fix limit_p_level_difference for p-coarsening. --- diff --git a/doc/news/changes/minor/20240229Fehling b/doc/news/changes/minor/20240229Fehling new file mode 100644 index 0000000000..eba14354b6 --- /dev/null +++ b/doc/news/changes/minor/20240229Fehling @@ -0,0 +1,4 @@ +Fixed: The function hp::Refinement::limit_p_level_difference() now +correctly sets future FE indices in case of p-coarsening. +
+(Marc Fehling, 2024/02/29) diff --git a/source/hp/refinement.cc b/source/hp/refinement.cc index f0eff39172..c9a4b218a4 100644 --- a/source/hp/refinement.cc +++ b/source/hp/refinement.cc @@ -1063,9 +1063,10 @@ namespace hp const unsigned int fe_index = fe_index_for_hierarchy_level[cell_level]; - // only update if necessary if (fe_index != cell->active_fe_index()) cell->set_future_fe_index(fe_index); + else + cell->clear_future_fe_index(); } } diff --git a/tests/hp/limit_p_level_difference_06.cc b/tests/hp/limit_p_level_difference_06.cc new file mode 100644 index 0000000000..583005709c --- /dev/null +++ b/tests/hp/limit_p_level_difference_06.cc @@ -0,0 +1,89 @@ +// ------------------------------------------------------------------------ +// +// SPDX-License-Identifier: LGPL-2.1-or-later +// Copyright (C) 2021 - 2024 by the deal.II authors +// +// This file is part of the deal.II library. +// +// Part of the source code is dual licensed under Apache-2.0 WITH +// LLVM-exception OR LGPL-2.1-or-later. Detailed license information +// governing the source code and code contributions can be found in +// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II. +// +// ------------------------------------------------------------------------ + + +// hp::Refinement::limit_p_level_difference() used to ignore updating future +// FE indices when they coincide with the active FE index of a cell. +// This posed a problem with p-coarsening. + + +#include + +#include + +#include +#include + +#include +#include + +#include "../tests.h" + +#include "../test_grids.h" + + +template +void +test() +{ + // setup FE collection + hp::FECollection fes; + while (fes.size() < 3) + fes.push_back(FE_Q(1)); + + // setup two cell triangulation + Triangulation tria; + TestGrids::hyper_line(tria, 2); + + DoFHandler dofh(tria); + dofh.distribute_dofs(fes); + + // set both cells to the same active FE index, + // and assign a future FE index one higher and one lower + for (const auto &cell : dofh.active_cell_iterators()) + { + cell->set_active_fe_index(1); + + if (cell->id().to_string() == "0_0:") + cell->set_future_fe_index(2); + else if (cell->id().to_string() == "1_0:") + cell->set_future_fe_index(0); + else + Assert(false, ExcInternalError()); + } + + const bool future_fe_indices_changed = + hp::Refinement::limit_p_level_difference(dofh); + + (void)future_fe_indices_changed; + Assert(future_fe_indices_changed, ExcInternalError()); + + // display FE indices for all cells + for (const auto &cell : dofh.active_cell_iterators()) + deallog << cell->id().to_string() << ": active:" << cell->active_fe_index() + << " future:" << cell->future_fe_index() + << " flag:" << cell->future_fe_index_set() << std::endl; + + deallog << "OK" << std::endl; +} + + +int +main() +{ + initlog(); + + test<2>(); + test<3>(); +} diff --git a/tests/hp/limit_p_level_difference_06.output b/tests/hp/limit_p_level_difference_06.output new file mode 100644 index 0000000000..997950631b --- /dev/null +++ b/tests/hp/limit_p_level_difference_06.output @@ -0,0 +1,7 @@ + +DEAL::0_0:: active:1 future:2 flag:1 +DEAL::1_0:: active:1 future:1 flag:0 +DEAL::OK +DEAL::0_0:: active:1 future:2 flag:1 +DEAL::1_0:: active:1 future:1 flag:0 +DEAL::OK