From: pauletti Date: Fri, 17 Dec 2010 20:33:18 +0000 (+0000) Subject: failing interpolate boundary 1d in codim==1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=17562ef1b2a8814304f54a583cd5db639654e888;p=dealii-svn.git failing interpolate boundary 1d in codim==1 git-svn-id: https://svn.dealii.org/trunk@22996 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/tests/codim_one/interpolate_boundary_values_03.cc b/tests/codim_one/interpolate_boundary_values_03.cc new file mode 100644 index 0000000000..8868f57b85 --- /dev/null +++ b/tests/codim_one/interpolate_boundary_values_03.cc @@ -0,0 +1,100 @@ + +//---------------------------- template.cc --------------------------- +// $Id$ +// Version: $Name$ +// +// Copyright (C) 2005, 2008, 2010 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. +// +//---------------------------- template.cc --------------------------- + + +// test VectorTools::interpolate_boundary_values for codim=1. like _01 +// but for 1d triangulations + +#include "../tests.h" +#include +#include +#include +#include +#include +#include +#include +#include +#include + +#include + +std::ofstream logfile("interpolate_boundary_values_03/output"); + +void test() { + const int dim = 1; + const int spacedim = 2; + + Triangulation tria; + Triangulation volume_mesh; + GridGenerator::half_hyper_ball(volume_mesh); + std::set boundary_ids; + boundary_ids.insert(0); + + GridTools::extract_boundary_mesh (volume_mesh, tria,boundary_ids); + + deallog << tria.n_active_cells() << " active cells" << std::endl; + + FE_Q fe(1); + DoFHandler dof_handler (tria); + dof_handler.distribute_dofs (fe); + + deallog << dof_handler.n_dofs() << " degrees of freedom" << std::endl; + + // test left and right boundary + // separately + for (unsigned int boundary_id=0; boundary_id<2; ++boundary_id) + { + std::map bv; + VectorTools::interpolate_boundary_values (dof_handler, + boundary_id, + Functions::SquareFunction(), + bv); + deallog << bv.size() << " boundary degrees of freedom" << std::endl; + + for (std::map::const_iterator i = bv.begin(); + i != bv.end(); ++i) + deallog << i->first << ' ' << i->second << std::endl; + + for (DoFHandler::active_cell_iterator + cell = dof_handler.begin_active(); cell != dof_handler.end(); ++cell) + for (unsigned int f=0; f::faces_per_cell; ++f) + if (cell->at_boundary(f) && + (cell->face(f)->boundary_indicator() == boundary_id)) + for (unsigned int v=0; v::vertices_per_face; ++v) + for (unsigned int i=0; iface(f)->vertex_dof_index(v,i)) + != bv.end(), + ExcInternalError()); + Assert (bv[cell->face(f)->vertex_dof_index(v,i)] + == + Functions::SquareFunction() + .value(cell->face(f)->vertex(v),i), + ExcInternalError()); + } + } +} + + + +int main () +{ + deallog.attach(logfile); + deallog.depth_console(0); + + test(); + + return 0; +} + diff --git a/tests/codim_one/interpolate_boundary_values_03/cmp/generic b/tests/codim_one/interpolate_boundary_values_03/cmp/generic new file mode 100644 index 0000000000..84e32afe29 --- /dev/null +++ b/tests/codim_one/interpolate_boundary_values_03/cmp/generic @@ -0,0 +1,7 @@ + +DEAL::1 active cells +DEAL::3 degrees of freedom +DEAL::1 boundary degrees of freedom +DEAL::0 0 +DEAL::1 boundary degrees of freedom +DEAL::1 1.00000