Skip to content

Commit

Permalink
fix bug in FESystem degreee with 0 multiplicity
Browse files Browse the repository at this point in the history
Bug was introduced in dealii#17062
  • Loading branch information
tjhei committed May 30, 2024
1 parent 16d3914 commit bb3d3eb
Show file tree
Hide file tree
Showing 3 changed files with 70 additions and 1 deletion.
3 changes: 2 additions & 1 deletion include/deal.II/fe/fe_tools.templates.h
Original file line number Diff line number Diff line change
Expand Up @@ -157,7 +157,8 @@ namespace FETools

multiplied_n_components += fes[i]->n_components() * multiplicities[i];

degree = std::max(degree, fes[i]->tensor_degree());
if (multiplicities[i] > 0)
degree = std::max(degree, fes[i]->tensor_degree());
}

// assume conformity of the first finite element and then take away
Expand Down
55 changes: 55 additions & 0 deletions tests/fe/system_05.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,55 @@
// ------------------------------------------------------------------------
//
// SPDX-License-Identifier: LGPL-2.1-or-later
// Copyright (C) 2012 - 2018 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.
//
// ------------------------------------------------------------------------


// Document wrong behavior with multiplicity 0 FESystem and incorrectly reported
// degree

#include <deal.II/fe/fe_q.h>
#include <deal.II/fe/fe_system.h>
#include <deal.II/fe/fe_tools.h>

#include <string>

#include "../tests.h"

template <int dim>
void
check(const FiniteElement<dim> &fe)
{
deallog << fe.get_name() << std::endl;
deallog << "components: " << fe.n_components() << std::endl;
deallog << "blocks: " << fe.n_blocks() << std::endl;
deallog << "conforms H1: " << fe.conforms(FiniteElementData<dim>::H1)
<< std::endl;
deallog << "n_base_elements: " << fe.n_base_elements() << std::endl;

deallog << "degree: " << fe.degree << std::endl;
}

int
main()
{
initlog();

{
FESystem<2> fe(FE_Q<2>(1) ^ 2, FE_Q<2>(1), FE_Q<2>(1), FE_Q<2>(2) ^ 0);
check<2>(fe);
Assert(fe.degree == 1);
auto copy = fe.clone();
check<2>(*copy.get());
Assert(fe.degree == 1);
}
return 0;
}
13 changes: 13 additions & 0 deletions tests/fe/system_05.output
Original file line number Diff line number Diff line change
@@ -0,0 +1,13 @@

DEAL::FESystem<2>[FE_Q<2>(1)^2-FE_Q<2>(1)-FE_Q<2>(1)]
DEAL::components: 4
DEAL::blocks: 4
DEAL::conforms H1: 1
DEAL::n_base_elements: 3
DEAL::degree: 1
DEAL::FESystem<2>[FE_Q<2>(1)^2-FE_Q<2>(1)-FE_Q<2>(1)]
DEAL::components: 4
DEAL::blocks: 4
DEAL::conforms H1: 1
DEAL::n_base_elements: 3
DEAL::degree: 1

0 comments on commit bb3d3eb

Please sign in to comment.