Skip to content

Commit

Permalink
Merge pull request #16705 from marcfehling/fix_limit-p-level
Browse files Browse the repository at this point in the history
Fix `hp::Refinement::limit_p_level_difference` for p-coarsening
  • Loading branch information
bangerth committed Mar 3, 2024
2 parents 60de7a6 + 5eea2e1 commit 3fc2fae
Show file tree
Hide file tree
Showing 4 changed files with 102 additions and 1 deletion.
4 changes: 4 additions & 0 deletions doc/news/changes/minor/20240229Fehling
Original file line number Diff line number Diff line change
@@ -0,0 +1,4 @@
Fixed: The function hp::Refinement::limit_p_level_difference() now
correctly sets future FE indices in case of p-coarsening.
<br>
(Marc Fehling, 2024/02/29)
3 changes: 2 additions & 1 deletion source/hp/refinement.cc
Original file line number Diff line number Diff line change
Expand Up @@ -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();
}
}

Expand Down
89 changes: 89 additions & 0 deletions tests/hp/limit_p_level_difference_06.cc
Original file line number Diff line number Diff line change
@@ -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 <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 <deal.II/hp/refinement.h>

#include "../tests.h"

#include "../test_grids.h"


template <int dim>
void
test()
{
// setup FE collection
hp::FECollection<dim> fes;
while (fes.size() < 3)
fes.push_back(FE_Q<dim>(1));

// setup two cell triangulation
Triangulation<dim> tria;
TestGrids::hyper_line(tria, 2);

DoFHandler<dim> 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>();
}
7 changes: 7 additions & 0 deletions tests/hp/limit_p_level_difference_06.output
Original file line number Diff line number Diff line change
@@ -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

0 comments on commit 3fc2fae

Please sign in to comment.