Skip to content

Commit

Permalink
Enable get_face_interpolation_matrix() in 1D on FE_Q and FE_Bernstein
Browse files Browse the repository at this point in the history
Currently, get_face_interpolation_matrix() on FE_Q and FE_Bernstein
throws an ExcImpossibleInDim(1) exception in 1D. This is strange
because the interpolation matrix should just be a 1-by-1 unit matrix.
Remove the dim > 1 assert in FE_Q_base and FE_Bernstein to fix this
and modify two tests to check that it works.
  • Loading branch information
simonsticko committed Jun 8, 2021
1 parent 0cb66ef commit ecd685f
Show file tree
Hide file tree
Showing 6 changed files with 634 additions and 31 deletions.
5 changes: 2 additions & 3 deletions source/fe/fe_bernstein.cc
Original file line number Diff line number Diff line change
Expand Up @@ -98,7 +98,6 @@ FE_Bernstein<dim, spacedim>::get_face_interpolation_matrix(
FullMatrix<double> & interpolation_matrix,
const unsigned int face_no) const
{
Assert(dim > 1, ExcImpossibleInDim(1));
get_subface_interpolation_matrix(source_fe,
numbers::invalid_unsigned_int,
interpolation_matrix,
Expand Down Expand Up @@ -145,8 +144,8 @@ FE_Bernstein<dim, spacedim>::get_subface_interpolation_matrix(
// Rule of thumb for FP accuracy, that can be expected for a given
// polynomial degree. This value is used to cut off values close to
// zero.
double eps =
2e-13 * std::max(this->degree, source_fe->degree) * (dim - 1);
const double eps = 2e-13 * std::max(this->degree, source_fe->degree) *
std::max(dim - 1, 1);

// compute the interpolation matrix by simply taking the value at the
// support points.
Expand Down
3 changes: 1 addition & 2 deletions source/fe/fe_q_base.cc
Original file line number Diff line number Diff line change
Expand Up @@ -611,7 +611,6 @@ FE_Q_Base<dim, spacedim>::get_face_interpolation_matrix(
FullMatrix<double> & interpolation_matrix,
const unsigned int face_no) const
{
Assert(dim > 1, ExcImpossibleInDim(1));
get_subface_interpolation_matrix(source_fe,
numbers::invalid_unsigned_int,
interpolation_matrix,
Expand Down Expand Up @@ -660,7 +659,7 @@ FE_Q_Base<dim, spacedim>::get_subface_interpolation_matrix(
// Rule of thumb for FP accuracy, that can be expected for a given
// polynomial degree. This value is used to cut off values close to
// zero.
double eps = 2e-13 * this->q_degree * (dim - 1);
const double eps = 2e-13 * this->q_degree * std::max(dim - 1, 1);

// compute the interpolation matrix by simply taking the value at the
// support points.
Expand Down
3 changes: 0 additions & 3 deletions tests/bits/face_interpolation.cc
Original file line number Diff line number Diff line change
Expand Up @@ -27,9 +27,6 @@ template <int dim>
void
check_this(const FiniteElement<dim> &fe1, const FiniteElement<dim> &fe2)
{
if (dim == 1)
return;

// check all combinations of fe1 and fe2
FullMatrix<double> face_constraints;
try
Expand Down
Loading

0 comments on commit ecd685f

Please sign in to comment.