From: Martin Kronbichler Date: Tue, 16 Mar 2021 15:41:42 +0000 (+0100) Subject: New test case X-Git-Tag: v9.3.0-rc1~320^2~1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=6401762ddbedeef2afc810308e985522cce7f13d;p=dealii.git New test case --- diff --git a/tests/matrix_free/boundary_id.cc b/tests/matrix_free/boundary_id.cc new file mode 100644 index 0000000000..b57e5e93d2 --- /dev/null +++ b/tests/matrix_free/boundary_id.cc @@ -0,0 +1,79 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2021 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. +// +// --------------------------------------------------------------------- + + + +// Check MatrixFree::boundary_id() for some large boundary id numbers +// (exceeding a very old case of 8 bit integers) + + +#include + +#include +#include + +#include +#include + +#include + +#include "../tests.h" + + +template +void +test() +{ + Triangulation tria; + GridGenerator::hyper_cube(tria, 0, 1); + unsigned int count = 0; + for (unsigned int face = 0; face < 2 * dim; ++face) + { + tria.begin()->face(face)->set_boundary_id(count); + count += 100000; + } + + MappingQGeneric mapping(1); + FE_Q fe(1); + DoFHandler dof(tria); + dof.distribute_dofs(fe); + AffineConstraints constraints; + constraints.close(); + + MatrixFree mf_data; + typename MatrixFree::AdditionalData data; + data.tasks_parallel_scheme = MatrixFree::AdditionalData::none; + data.mapping_update_flags_inner_faces = + (update_gradients | update_JxW_values); + data.mapping_update_flags_boundary_faces = + (update_gradients | update_JxW_values); + + mf_data.reinit(mapping, dof, constraints, QGauss<1>(2), data); + for (unsigned int i = 0; i < mf_data.n_boundary_face_batches(); ++i) + deallog << "Face " + << static_cast( + mf_data.get_face_info(i).interior_face_no) + << " boundary id " << mf_data.get_boundary_id(i) << std::endl; +} + + +int +main() +{ + initlog(); + + test<2>(); + test<3>(); +} diff --git a/tests/matrix_free/boundary_id.output b/tests/matrix_free/boundary_id.output new file mode 100644 index 0000000000..2ff9f1cb8b --- /dev/null +++ b/tests/matrix_free/boundary_id.output @@ -0,0 +1,11 @@ + +DEAL::Face 0 boundary id 0 +DEAL::Face 1 boundary id 100000 +DEAL::Face 2 boundary id 200000 +DEAL::Face 3 boundary id 300000 +DEAL::Face 0 boundary id 0 +DEAL::Face 1 boundary id 100000 +DEAL::Face 2 boundary id 200000 +DEAL::Face 3 boundary id 300000 +DEAL::Face 4 boundary id 400000 +DEAL::Face 5 boundary id 500000