From: Marc Fehling Date: Fri, 14 Apr 2023 04:35:10 +0000 (-0600) Subject: Add a test for FEInterfaceValues::reinit with no dominating FE. X-Git-Tag: v9.5.0-rc1~326^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F15103%2Fhead;p=dealii.git Add a test for FEInterfaceValues::reinit with no dominating FE. --- diff --git a/tests/feinterface/fe_interface_values_13.cc b/tests/feinterface/fe_interface_values_13.cc new file mode 100644 index 0000000000..57ef0acdbe --- /dev/null +++ b/tests/feinterface/fe_interface_values_13.cc @@ -0,0 +1,95 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2023 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. +// +// --------------------------------------------------------------------- + + +// Set up a two cells with different elements, such that neither element +// dominates. Create an FEInterfaceValues object and call reinit on it. + +#include + +#include +#include +#include + +#include + +#include "../tests.h" + +#include "../test_grids.h" + + +template +void +test() +{ + Triangulation tria; + TestGrids::hyper_line(tria, 2); + + const FESystem fe_left(FE_Q(1), FE_Q(2)); + const QGauss q_left(3); + + const FESystem fe_right(FE_Q(3), FE_Q(1)); + const QGauss q_right(4); + + const hp::FECollection fe_collection(fe_left, fe_right); + const hp::QCollection q_collection(q_left, q_right); + + + const unsigned int face_index = 1; + + + DoFHandler dofh(tria); + + auto cell = dofh.begin(); + auto neighbor = cell->neighbor(face_index); + + cell->set_active_fe_index(0); + neighbor->set_active_fe_index(1); + + dofh.distribute_dofs(fe_collection); + + + const UpdateFlags update_flags = update_values | update_gradients | + update_quadrature_points | update_JxW_values; + FEInterfaceValues fiv(fe_collection, q_collection, update_flags); + + try + { + fiv.reinit(cell, + face_index, + numbers::invalid_unsigned_int, + neighbor, + cell->neighbor_of_neighbor(face_index), + numbers::invalid_unsigned_int); + } + catch (ExceptionBase &exc) + { + deallog << exc.what() << std::endl; + } + + deallog << "OK" << std::endl; +} + + + +int +main() +{ + deal_II_exceptions::disable_abort_on_exception(); + + initlog(); + + test<2>(); +} diff --git a/tests/feinterface/fe_interface_values_13.debug.output b/tests/feinterface/fe_interface_values_13.debug.output new file mode 100644 index 0000000000..06ec85d5ad --- /dev/null +++ b/tests/feinterface/fe_interface_values_13.debug.output @@ -0,0 +1,15 @@ + +DEAL:: +-------------------------------------------------------- +An error occurred in file in function + void dealii::FEInterfaceValues::reinit(const CellIteratorType&, unsigned int, unsigned int, const CellNeighborIteratorType&, unsigned int, unsigned int, unsigned int, unsigned int, unsigned int) [with CellIteratorType = dealii::TriaIterator >; CellNeighborIteratorType = dealii::TriaIterator >; int dim = 2; int spacedim = 2] +The violated condition was: + dominated_fe_index != numbers::invalid_fe_index +Additional information: + You called this function with 'q_index' left at its default value, but + this can only work if one of the two finite elements adjacent to this + face dominates the other. See the documentation of this function for + more information of how to deal with this situation. +-------------------------------------------------------- + +DEAL::OK