From: Jean-Paul Pelteret Date: Fri, 7 Jan 2022 20:58:01 +0000 (+0100) Subject: Add a test X-Git-Tag: v9.4.0-rc1~398^2~6 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=f5000c61d91dad2b07b6e85df530095cd8b96a59;p=dealii.git Add a test --- diff --git a/tests/meshworker/scratch_data_01_hp.cc b/tests/meshworker/scratch_data_01_hp.cc new file mode 100644 index 0000000000..84a361a255 --- /dev/null +++ b/tests/meshworker/scratch_data_01_hp.cc @@ -0,0 +1,166 @@ +/* --------------------------------------------------------------------- + * + * Copyright (C) 2019 - 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. + * + * --------------------------------------------------------------------- + */ + +// Test ScratchData on simple volume computation +// - hp variant + +#include + +#include +#include +#include +#include +#include + +#include +#include + +#include +#include +#include + +#include +#include +#include +#include + +#include "../tests.h" + +using namespace MeshWorker; + +template +void +test() +{ + Triangulation tria; + hp::FECollection fe(FE_Q(1), + FE_Q(2)); + DoFHandler dh(tria); + GridGenerator::hyper_cube(tria); + tria.refine_global(1); + tria.begin_active()->set_refine_flag(); + tria.execute_coarsening_and_refinement(); + + const auto bounding_box = GridTools::compute_bounding_box(tria); + for (auto &cell : dh.active_cell_iterators()) + { + bool use_alt_index = true; + for (unsigned int d = 0; d < spacedim; ++d) + if (cell->center()[d] < bounding_box.center()[d]) + { + use_alt_index = false; + break; + } + + if (use_alt_index) + cell->set_active_fe_index(1); + else + cell->set_active_fe_index(0); + } + + dh.distribute_dofs(fe); + + hp::QCollection cell_quad(QGauss(3), QGauss(4)); + hp::QCollection face_quad(QGauss(3), QGauss(4)); + + UpdateFlags cell_flags = update_JxW_values; + UpdateFlags face_flags = update_JxW_values; + + using ScratchData = MeshWorker::ScratchData; + using CopyData = MeshWorker::CopyData<0, 1, 0>; + + ScratchData scratch(fe, cell_quad, cell_flags, face_quad, face_flags); + CopyData copy(3); // store measures + + auto cell = dh.begin_active(); + auto endc = dh.end(); + + using Iterator = decltype(cell); + + auto measures = std::make_tuple(0.0, 0.0, 0.0); + + + auto cell_worker = [](const Iterator &cell, ScratchData &s, CopyData &c) { + const auto &fev = s.reinit(cell); + const auto &JxW = s.get_JxW_values(); + for (auto w : JxW) + c.vectors[0][0] += w; + }; + + auto boundary_worker = [](const Iterator & cell, + const unsigned int &f, + ScratchData & s, + CopyData & c) { + const auto &fev = s.reinit(cell, f); + const auto &JxW = s.get_JxW_values(); + for (auto w : JxW) + c.vectors[0][1] += w; + }; + + auto face_worker = [](const Iterator & cell, + const unsigned int &f, + const unsigned int &sf, + const Iterator & ncell, + const unsigned int &nf, + const unsigned int &nsf, + ScratchData & s, + CopyData & c) { + const auto &fev = s.reinit(cell, f, sf); + const auto &nfev = s.reinit_neighbor(ncell, nf, nsf); + + const auto &JxW = s.get_JxW_values(); + const auto &nJxW = s.get_neighbor_JxW_values(); + for (auto w : JxW) + c.vectors[0][2] += w; + }; + + auto copier = [&measures](const CopyData &c) { + std::get<0>(measures) += c.vectors[0][0]; + std::get<1>(measures) += c.vectors[0][1]; + std::get<2>(measures) += c.vectors[0][2]; + }; + + mesh_loop(cell, + endc, + cell_worker, + copier, + scratch, + copy, + assemble_own_cells | assemble_boundary_faces | + assemble_own_interior_faces_once, + boundary_worker, + face_worker); + + deallog << "Testing <" << dim << "," << spacedim << ">" << std::endl; + + deallog << "Volume: " << std::get<0>(measures) << std::endl; + + deallog << "Boundary surface: " << std::get<1>(measures) << std::endl; + + deallog << "Interior surface: " << std::get<2>(measures) << std::endl; +} + + + +int +main() +{ + initlog(); + + test<2, 2>(); + test<2, 3>(); + test<3, 3>(); +} diff --git a/tests/meshworker/scratch_data_01_hp.output b/tests/meshworker/scratch_data_01_hp.output new file mode 100644 index 0000000000..b9c1ecd3d3 --- /dev/null +++ b/tests/meshworker/scratch_data_01_hp.output @@ -0,0 +1,13 @@ + +DEAL::Testing <2,2> +DEAL::Volume: 1.00000 +DEAL::Boundary surface: 4.00000 +DEAL::Interior surface: 3.00000 +DEAL::Testing <2,3> +DEAL::Volume: 1.00000 +DEAL::Boundary surface: 4.00000 +DEAL::Interior surface: 3.00000 +DEAL::Testing <3,3> +DEAL::Volume: 1.00000 +DEAL::Boundary surface: 6.00000 +DEAL::Interior surface: 3.75000