From: Luca Heltai Date: Tue, 9 Apr 2019 21:05:17 +0000 (+0200) Subject: Test ScratchData get_values X-Git-Tag: v9.1.0-rc1~210^2~6 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=c4a8e9cadee0a4e53eeca593bddc52265189d636;p=dealii.git Test ScratchData get_values --- diff --git a/tests/meshworker/scratch_data_04.cc b/tests/meshworker/scratch_data_04.cc new file mode 100644 index 0000000000..4887755155 --- /dev/null +++ b/tests/meshworker/scratch_data_04.cc @@ -0,0 +1,131 @@ +/* --------------------------------------------------------------------- + * + * Copyright (C) 2018 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. + * + * --------------------------------------------------------------------- + */ + +// Compute integral of known function using mesh_loop and ScratchData + +#include +#include +#include + +#include + +#include +#include + +#include +#include +#include +#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; + FE_Q fe(1); + DoFHandler dh(tria); + + FunctionParser integral_function("x"); + GridGenerator::hyper_cube(tria); + tria.refine_global(5); + tria.execute_coarsening_and_refinement(); + dh.distribute_dofs(fe); + + Vector solution(dh.n_dofs()); + + VectorTools::interpolate(dh, integral_function, solution); + + QGauss quad(3); + + UpdateFlags cell_flags = update_values | update_gradients | + update_quadrature_points | update_JxW_values; + + using ScratchData = ScratchData; + struct CopyData + {}; + + + double H1_norm = 0; + FEValuesExtractors::Scalar scalar(0); + + ScratchData scratch(fe, quad, cell_flags); + CopyData copy; + + auto cell = dh.begin_active(); + auto endc = dh.end(); + + typedef decltype(cell) Iterator; + + + auto cell_integrator = [&H1_norm, &solution, &scalar](const Iterator &cell, + ScratchData & s, + CopyData & c) { + const auto &fev = s.reinit(cell); + const auto &JxW = s.get_JxW_values(); + + s.extract_local_dof_values("solution", solution); + + const auto &u = s.get_values("solution", scalar); + const auto &grad_u = s.get_gradients("solution", scalar); + + for (unsigned int q = 0; q < u.size(); ++q) + H1_norm += (u[q] * u[q] + grad_u[q] * grad_u[q]) * JxW[q]; + }; + + const auto empty_copyer = [](const CopyData &) {}; + + mesh_loop(cell, + endc, + cell_integrator, + empty_copyer, + scratch, + copy, + assemble_own_cells); + + deallog << "H1 norm: " << std::sqrt(H1_norm) << std::endl; +} + + +int +main() +{ + initlog(1); + MultithreadInfo::set_thread_limit(1); // to make output deterministic + + test<2, 2>(); +} diff --git a/tests/meshworker/scratch_data_04.with_muparser=on.output b/tests/meshworker/scratch_data_04.with_muparser=on.output new file mode 100644 index 0000000000..b906801d0d --- /dev/null +++ b/tests/meshworker/scratch_data_04.with_muparser=on.output @@ -0,0 +1,2 @@ + +DEAL::H1 norm: 1.15470