From: heltai Date: Wed, 29 Aug 2007 20:15:32 +0000 (+0000) Subject: Added ParsedFunction and FEFieldFunction tests. X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=910e2f6ed2d8876769a519265d1cbf97ce8c40ce;p=dealii-svn.git Added ParsedFunction and FEFieldFunction tests. git-svn-id: https://svn.dealii.org/trunk@15090 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/tests/bits/Makefile b/tests/bits/Makefile index 164ba8b138..ea9420c0f2 100644 --- a/tests/bits/Makefile +++ b/tests/bits/Makefile @@ -34,12 +34,14 @@ tests_x = grid_generator_?? \ fe_tools_[0-9]* \ fe_tools_cpfqpm* \ fe_face_to_cell_index \ + fe_field_function \ anna_? \ gerold_1 \ roy_1 \ denis_1 \ unit_support_points \ parameter_handler_* \ + parsed_function \ sparse_lu_decomposition_1 \ block_sparse_matrix_* \ block_vector_* \ diff --git a/tests/bits/fe_field_function.cc b/tests/bits/fe_field_function.cc new file mode 100644 index 0000000000..b98636b327 --- /dev/null +++ b/tests/bits/fe_field_function.cc @@ -0,0 +1,94 @@ +//---------------------------- template.cc --------------------------- +// $Id: testsuite.html 13373 2006-07-13 13:12:08Z kanschat $ +// Version: $Name$ +// +// Copyright (C) 2005 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 --------------------------- + + +// a short (a few lines) description of what the program does + +#include "../tests.h" +#include + +// all include files you need here +#include +#include +#include +#include +#include +#include +#include +#include +#include + +double abs_zero(double a) { + if( std::abs(a) < 1e-10) + return 0; + else + return a; +} + +template +void test() { + Triangulation tria; + GridGenerator::hyper_cube(tria); + tria.refine_global(8/dim); + + FE_Q fe(1); + DoFHandler dh(tria); + + dh.distribute_dofs(fe); + deallog << "Ndofs :" << dh.n_dofs() << std::endl; + + Functions::CosineFunction ff; + + Vector v1(dh.n_dofs()), v2(dh.n_dofs()); + + VectorTools::interpolate(dh, ff, v1); + deallog << "V norm: " << v1.l2_norm() << std::endl; + + Functions::FEFieldFunction, Vector > + fef(dh, v1); + + VectorTools::interpolate(dh, fef, v2); + + v2.add(-1, v1); + deallog << "Interpolation error: " << abs_zero(v2.l2_norm()) + << std::endl; + + Vector error(tria.n_active_cells()); + QGauss quad(2); + VectorTools::integrate_difference(dh, v1, ff, + error, quad, + VectorTools::H1_norm); + deallog << "H1 Interpolation error: " + << abs_zero(error.l2_norm()) << std::endl; + error = 0; + VectorTools::integrate_difference(dh, v1, fef, + error, quad, + VectorTools::H1_norm); + deallog << "H1 Interpolation error with fef: " + << abs_zero(error.l2_norm()) << std::endl; + +} + +int main () +{ + std::ofstream logfile("fe_field_function/output"); + deallog.attach(logfile); + deallog.depth_console(0); + + test<1>(); + test<2>(); + test<3>(); + + return 0; +} + diff --git a/tests/bits/fe_field_function/cmp/generic b/tests/bits/fe_field_function/cmp/generic new file mode 100644 index 0000000000..bb356f97aa --- /dev/null +++ b/tests/bits/fe_field_function/cmp/generic @@ -0,0 +1,16 @@ + +DEAL::Ndofs :257 +DEAL::V norm: 11.3358 +DEAL::Interpolation error: 0 +DEAL::H1 Interpolation error: 0.00196741 +DEAL::H1 Interpolation error with fef: 0 +DEAL::Ndofs :289 +DEAL::V norm: 8.50000 +DEAL::Interpolation error: 0 +DEAL::H1 Interpolation error: 0.0314971 +DEAL::H1 Interpolation error with fef: 0 +DEAL::Ndofs :125 +DEAL::V norm: 3.95285 +DEAL::Interpolation error: 0 +DEAL::H1 Interpolation error: 0.112288 +DEAL::H1 Interpolation error with fef: 0 diff --git a/tests/bits/parsed_function.cc b/tests/bits/parsed_function.cc new file mode 100644 index 0000000000..fc4bda239f --- /dev/null +++ b/tests/bits/parsed_function.cc @@ -0,0 +1,102 @@ +//---------------------------- function_parser.cc --------------------------- +// $Id: function_parser.cc 11749 2005-11-09 19:11:20Z wolf $ +// Version: $Name$ +// +// Copyright (C) 2005 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. +// +//---------------------------- function_parser.cc --------------------------- + + +// This program tests the functionality of the function parser +// wrapper. + +#include "../tests.h" +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include + +template +void Test() { + // A parameter handler + ParameterHandler prm; + + // A triangulation + Triangulation tria; + GridGenerator::hyper_cube(tria,0,1); + tria.refine_global(3); + + // Vertices + const std::vector > &vertices = tria.get_vertices(); + + // Test vector declaration + for(unsigned int i=0; i::declare_parameters(prm, i+1); + prm.set("Function constants", "f=" + Utilities::int_to_string(i+1)); + + // It is cos(pi f x_i) *t^i, numbering from zero + std::string expr = "cos(pi*f*x)*t"; + if(i>0) expr+= "; cos(pi*f*y)*t"; + if(i>1) expr+= "; cos(pi*f*z)*t"; + + prm.set("Function expression", expr); + + Functions::ParsedFunction function(i+1); + function.parse_parameters(prm); + + prm.leave_subsection(); + + // Now test the difference from t=0 to t=1 + for(double t=0.; t<1; t+= .1) { + function.set_time(t); + Point p; + std::vector > values(vertices.size(), + Vector(i+1)); + function.vector_value_list(vertices, values); + for(unsigned int j=0; j 1e-10) + deallog << "p(" << di << "): " << vertices[j] << ", delta: " + << delta << std::endl; + } + } + } + } + deallog << "Tested on: " << std::endl; + prm.log_parameters(deallog); +} + +int main () +{ + std::ofstream logfile("parsed_function/output"); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + Test<1>(); + Test<2>(); + Test<3>(); +} + + + + diff --git a/tests/bits/parsed_function/cmp/generic b/tests/bits/parsed_function/cmp/generic new file mode 100644 index 0000000000..30938838a3 --- /dev/null +++ b/tests/bits/parsed_function/cmp/generic @@ -0,0 +1,22 @@ + +DEAL::Tested on: +DEAL:parameters:Function 1 - 0::Function constants: f=1 +DEAL:parameters:Function 1 - 0::Function expression: cos(pi*f*x)*t +DEAL:parameters:Function 1 - 0::Variable names: x,t +DEAL::Tested on: +DEAL:parameters:Function 2 - 0::Function constants: f=1 +DEAL:parameters:Function 2 - 0::Function expression: cos(pi*f*x)*t +DEAL:parameters:Function 2 - 0::Variable names: x,y,t +DEAL:parameters:Function 2 - 1::Function constants: f=2 +DEAL:parameters:Function 2 - 1::Function expression: cos(pi*f*x)*t; cos(pi*f*y)*t +DEAL:parameters:Function 2 - 1::Variable names: x,y,t +DEAL::Tested on: +DEAL:parameters:Function 3 - 0::Function constants: f=1 +DEAL:parameters:Function 3 - 0::Function expression: cos(pi*f*x)*t +DEAL:parameters:Function 3 - 0::Variable names: x,y,z,t +DEAL:parameters:Function 3 - 1::Function constants: f=2 +DEAL:parameters:Function 3 - 1::Function expression: cos(pi*f*x)*t; cos(pi*f*y)*t +DEAL:parameters:Function 3 - 1::Variable names: x,y,z,t +DEAL:parameters:Function 3 - 2::Function constants: f=3 +DEAL:parameters:Function 3 - 2::Function expression: cos(pi*f*x)*t; cos(pi*f*y)*t; cos(pi*f*z)*t +DEAL:parameters:Function 3 - 2::Variable names: x,y,z,t