From bbf712c9716ff8054816d09495f2e156749995d9 Mon Sep 17 00:00:00 2001 From: bangerth Date: Thu, 10 Aug 2006 19:24:22 +0000 Subject: [PATCH] Same test, but for DG again git-svn-id: https://svn.dealii.org/trunk@13638 0785d39b-7218-0410-832d-ea1e28bc413d --- tests/deal.II/interpolate_dgq_02.cc | 127 +++++++++++++++++++ tests/deal.II/interpolate_dgq_02/cmp/generic | 47 +++++++ tests/deal.II/interpolate_q_02.cc | 2 +- 3 files changed, 175 insertions(+), 1 deletion(-) create mode 100644 tests/deal.II/interpolate_dgq_02.cc create mode 100644 tests/deal.II/interpolate_dgq_02/cmp/generic diff --git a/tests/deal.II/interpolate_dgq_02.cc b/tests/deal.II/interpolate_dgq_02.cc new file mode 100644 index 0000000000..164bf3b716 --- /dev/null +++ b/tests/deal.II/interpolate_dgq_02.cc @@ -0,0 +1,127 @@ +//---------------------------- interpolate_dgq_02.cc --------------------------- +// $Id: interpolate_dgq_02.cc 12732 2006-03-28 23:15:45Z wolf $ +// Version: $Name$ +// +// Copyright (C) 2006 by the deal.II authors +// +// This file is subject to QPL and may not be distributed +// without copyright and license information. Please refer +// to the file deal.II/doc/license.html for the text and +// further information on this license. +// +//---------------------------- interpolate_dgq_02.cc --------------------------- + + +// check that VectorTools::interpolate works for FE_DGQ(p) elements correctly on +// an adaptively refined mesh for functions of degree q + +#include "../tests.h" +#include +#include +#include +#include + +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include + +#include +#include + + +template +class F : public Function +{ + public: + F (const unsigned int q) : q(q) {} + + virtual double value (const Point &p, + const unsigned int) const + { + double v=0; + for (unsigned int d=0; d +void test () +{ + Triangulation triangulation; + GridGenerator::hyper_cube (triangulation); + triangulation.refine_global (1); + triangulation.begin_active()->set_refine_flag (); + triangulation.execute_coarsening_and_refinement (); + triangulation.refine_global (1); + + for (unsigned int p=1; p<6-dim; ++p) + { + FE_DGQ fe(p); + DoFHandler dof_handler(triangulation); + dof_handler.distribute_dofs (fe); + + ConstraintMatrix constraints; + DoFTools::make_hanging_node_constraints (dof_handler, constraints); + constraints.close (); + + Vector interpolant (dof_handler.n_dofs()); + Vector error (triangulation.n_active_cells()); + for (unsigned int q=0; q<=p+2; ++q) + { + // interpolate the function + VectorTools::interpolate (dof_handler, + F (q), + interpolant); + constraints.distribute (interpolant); + + // then compute the interpolation error + VectorTools::integrate_difference (dof_handler, + interpolant, + F (q), + error, + QGauss(q+2), + VectorTools::L2_norm); + if (q<=p) + Assert (error.l2_norm() < 1e-12*interpolant.l2_norm(), + ExcInternalError()); + + deallog << fe.get_name() << ", P_" << q + << ", rel. error=" << error.l2_norm() / interpolant.l2_norm() + << std::endl; + } + } +} + + + +int main () +{ + std::ofstream logfile("interpolate_dgq_02/output"); + logfile.precision (3); + + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + test<1>(); + test<2>(); + test<3>(); +} + diff --git a/tests/deal.II/interpolate_dgq_02/cmp/generic b/tests/deal.II/interpolate_dgq_02/cmp/generic new file mode 100644 index 0000000000..ad68d0fd2b --- /dev/null +++ b/tests/deal.II/interpolate_dgq_02/cmp/generic @@ -0,0 +1,47 @@ + +DEAL::FE_DGQ<1>(1), P_0, rel. error=0 +DEAL::FE_DGQ<1>(1), P_1, rel. error=0 +DEAL::FE_DGQ<1>(1), P_2, rel. error=0.00243 +DEAL::FE_DGQ<1>(1), P_3, rel. error=0.00676 +DEAL::FE_DGQ<1>(2), P_0, rel. error=0 +DEAL::FE_DGQ<1>(2), P_1, rel. error=0 +DEAL::FE_DGQ<1>(2), P_2, rel. error=0 +DEAL::FE_DGQ<1>(2), P_3, rel. error=8.78e-05 +DEAL::FE_DGQ<1>(2), P_4, rel. error=0.000315 +DEAL::FE_DGQ<1>(3), P_0, rel. error=0 +DEAL::FE_DGQ<1>(3), P_1, rel. error=0 +DEAL::FE_DGQ<1>(3), P_2, rel. error=0 +DEAL::FE_DGQ<1>(3), P_3, rel. error=0 +DEAL::FE_DGQ<1>(3), P_4, rel. error=3.99e-06 +DEAL::FE_DGQ<1>(3), P_5, rel. error=1.74e-05 +DEAL::FE_DGQ<1>(4), P_0, rel. error=0 +DEAL::FE_DGQ<1>(4), P_1, rel. error=0 +DEAL::FE_DGQ<1>(4), P_2, rel. error=0 +DEAL::FE_DGQ<1>(4), P_3, rel. error=0 +DEAL::FE_DGQ<1>(4), P_4, rel. error=0 +DEAL::FE_DGQ<1>(4), P_5, rel. error=2.07e-07 +DEAL::FE_DGQ<1>(4), P_6, rel. error=1.06e-06 +DEAL::FE_DGQ<2>(1), P_0, rel. error=0 +DEAL::FE_DGQ<2>(1), P_1, rel. error=0 +DEAL::FE_DGQ<2>(1), P_2, rel. error=0.00100 +DEAL::FE_DGQ<2>(1), P_3, rel. error=0.00249 +DEAL::FE_DGQ<2>(2), P_0, rel. error=0 +DEAL::FE_DGQ<2>(2), P_1, rel. error=0 +DEAL::FE_DGQ<2>(2), P_2, rel. error=0 +DEAL::FE_DGQ<2>(2), P_3, rel. error=2.41e-05 +DEAL::FE_DGQ<2>(2), P_4, rel. error=7.74e-05 +DEAL::FE_DGQ<2>(3), P_0, rel. error=0 +DEAL::FE_DGQ<2>(3), P_1, rel. error=0 +DEAL::FE_DGQ<2>(3), P_2, rel. error=0 +DEAL::FE_DGQ<2>(3), P_3, rel. error=0 +DEAL::FE_DGQ<2>(3), P_4, rel. error=1.07e-06 +DEAL::FE_DGQ<2>(3), P_5, rel. error=4.06e-06 +DEAL::FE_DGQ<3>(1), P_0, rel. error=0 +DEAL::FE_DGQ<3>(1), P_1, rel. error=0 +DEAL::FE_DGQ<3>(1), P_2, rel. error=0.000375 +DEAL::FE_DGQ<3>(1), P_3, rel. error=0.000898 +DEAL::FE_DGQ<3>(2), P_0, rel. error=0 +DEAL::FE_DGQ<3>(2), P_1, rel. error=0 +DEAL::FE_DGQ<3>(2), P_2, rel. error=0 +DEAL::FE_DGQ<3>(2), P_3, rel. error=6.40e-06 +DEAL::FE_DGQ<3>(2), P_4, rel. error=1.98e-05 diff --git a/tests/deal.II/interpolate_q_02.cc b/tests/deal.II/interpolate_q_02.cc index 469126cfd9..c60cda04f6 100644 --- a/tests/deal.II/interpolate_q_02.cc +++ b/tests/deal.II/interpolate_q_02.cc @@ -12,7 +12,7 @@ //---------------------------- interpolate_q_02.cc --------------------------- -// check that VectorTools::interpolate works for FE_DGQ(p) elements correctly on +// check that VectorTools::interpolate works for FE_Q(p) elements correctly on // an adaptively refined mesh for functions of degree q #include "../tests.h" -- 2.39.5