]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Added ParsedFunction and FEFieldFunction tests.
authorheltai <heltai@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 29 Aug 2007 20:15:32 +0000 (20:15 +0000)
committerheltai <heltai@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 29 Aug 2007 20:15:32 +0000 (20:15 +0000)
git-svn-id: https://svn.dealii.org/trunk@15090 0785d39b-7218-0410-832d-ea1e28bc413d

tests/bits/Makefile
tests/bits/fe_field_function.cc [new file with mode: 0644]
tests/bits/fe_field_function/cmp/generic [new file with mode: 0644]
tests/bits/parsed_function.cc [new file with mode: 0644]
tests/bits/parsed_function/cmp/generic [new file with mode: 0644]

index 164ba8b138eebc2e054e9bc42bf046e601b46ce4..ea9420c0f29c69b484609c55f9691e6c79f6cb8a 100644 (file)
@@ -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 (file)
index 0000000..b98636b
--- /dev/null
@@ -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 <fstream>
+
+// all include files you need here
+#include <numerics/fe_field_function.h>
+#include <grid/tria.h>
+#include <grid/grid_generator.h>
+#include <dofs/dof_handler.h>
+#include <dofs/dof_tools.h>
+#include <fe/fe_q.h>
+#include <base/function_lib.h>
+#include <base/quadrature_lib.h>
+#include <numerics/vectors.h>
+
+double abs_zero(double a) {
+    if( std::abs(a) < 1e-10)
+       return 0;
+    else
+       return a;
+}
+
+template <int dim> 
+void test() {
+  Triangulation<dim> tria;
+  GridGenerator::hyper_cube(tria);
+  tria.refine_global(8/dim);
+
+  FE_Q<dim> fe(1);
+  DoFHandler<dim> dh(tria);
+  
+  dh.distribute_dofs(fe);
+  deallog << "Ndofs :" << dh.n_dofs() << std::endl;
+
+  Functions::CosineFunction<dim> ff;
+  
+  Vector<double> v1(dh.n_dofs()), v2(dh.n_dofs());
+  
+  VectorTools::interpolate(dh, ff, v1);
+  deallog << "V norm: " << v1.l2_norm() << std::endl;
+  
+  Functions::FEFieldFunction<dim, DoFHandler<dim>, Vector<double> > 
+    fef(dh, v1);
+  
+  VectorTools::interpolate(dh, fef, v2);
+  
+  v2.add(-1, v1);
+  deallog << "Interpolation error: " << abs_zero(v2.l2_norm())
+         << std::endl;
+  
+  Vector<double> error(tria.n_active_cells());
+  QGauss<dim> 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 (file)
index 0000000..bb356f9
--- /dev/null
@@ -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 (file)
index 0000000..fc4bda2
--- /dev/null
@@ -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 <fstream>
+#include <iostream>
+#include <map>
+#include <base/logstream.h>
+#include <base/point.h>
+#include <lac/vector.h>
+#include <base/parsed_function.h>
+#include <base/parameter_handler.h>
+#include <base/utilities.h>
+#include <grid/tria.h>
+#include <grid/grid_generator.h>
+
+template <int dim> 
+void Test() {
+  // A parameter handler
+  ParameterHandler prm;
+  
+  // A triangulation
+  Triangulation<dim> tria;
+  GridGenerator::hyper_cube(tria,0,1);
+  tria.refine_global(3);
+  
+  // Vertices
+  const std::vector<Point<dim> > &vertices = tria.get_vertices();
+  
+  // Test vector declaration
+  for(unsigned int i=0; i<dim;++i) {
+    std::string id = 
+      "Function " + Utilities::int_to_string(dim) +
+      " - " + Utilities::int_to_string(i);
+    prm.enter_subsection(id);
+    
+    Functions::ParsedFunction<dim>::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<dim> 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<dim> p;
+      std::vector<Vector<double> > values(vertices.size(),
+                                         Vector<double>(i+1));
+      function.vector_value_list(vertices, values);
+      for(unsigned int j=0; j<vertices.size(); ++j) {
+       double delta=0.;
+       for(unsigned int di=0; di<i; ++di) {
+         delta = values[j](di) - std::cos(M_PI*(i+1)*vertices[j][di])*t;
+         if(std::abs(delta) > 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 (file)
index 0000000..3093883
--- /dev/null
@@ -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

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.