]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
More tests, this time for systems of Q elements
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 10 Aug 2006 19:32:53 +0000 (19:32 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 10 Aug 2006 19:32:53 +0000 (19:32 +0000)
git-svn-id: https://svn.dealii.org/trunk@13639 0785d39b-7218-0410-832d-ea1e28bc413d

tests/deal.II/interpolate_dgq_01.cc
tests/deal.II/interpolate_q_01.cc
tests/deal.II/interpolate_q_system_01.cc [new file with mode: 0644]
tests/deal.II/interpolate_q_system_01/cmp/generic [new file with mode: 0644]
tests/deal.II/interpolate_q_system_02.cc [new file with mode: 0644]
tests/deal.II/interpolate_q_system_02/cmp/generic [new file with mode: 0644]

index 1d7c6b6cfba0d7953549558ab0054458388fceb8..18e88840c16148376ebf0ac33189581585e103e3 100644 (file)
@@ -12,7 +12,7 @@
 //----------------------------  interpolate_dgq_01.cc  ---------------------------
 
 
-// check that VectorTools::interpolate works for FE_Q(p) elements correctly on
+// check that VectorTools::interpolate works for FE_DGQ(p) elements correctly on
 // a uniformly refined mesh for functions of degree q
 
 #include "../tests.h"
index 7a2b5de7d03209dd90bf0ebddb12ab9f974a9294..2689c4f983cabf7e322def59ca7c6f8906298cc0 100644 (file)
@@ -12,7 +12,7 @@
 //----------------------------  interpolate_q_01.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
 // a uniformly refined mesh for functions of degree q
 
 #include "../tests.h"
diff --git a/tests/deal.II/interpolate_q_system_01.cc b/tests/deal.II/interpolate_q_system_01.cc
new file mode 100644 (file)
index 0000000..35d5a53
--- /dev/null
@@ -0,0 +1,121 @@
+//----------------------------  interpolate_q_system_01.cc  ---------------------------
+//    $Id: interpolate_q_system_01.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_q_system_01.cc  ---------------------------
+
+
+// check that VectorTools::interpolate works for FE_System(FE_Q(p)) elements correctly on
+// a uniformly refined mesh for functions of degree q
+
+#include "../tests.h"
+#include <base/function.h>
+#include <base/logstream.h>
+#include <base/quadrature_lib.h>
+#include <lac/vector.h>
+
+#include <grid/tria.h>
+#include <dofs/dof_handler.h>
+#include <grid/grid_generator.h>
+#include <grid/grid_refinement.h>
+#include <grid/tria_accessor.h>
+#include <grid/tria_iterator.h>
+#include <grid/tria_boundary_lib.h>
+#include <dofs/dof_accessor.h>
+#include <dofs/dof_tools.h>
+#include <numerics/vectors.h>
+#include <fe/fe_q.h>
+#include <fe/fe_system.h>
+
+#include <fstream>
+#include <vector>
+
+
+template <int dim>
+class F :  public Function<dim>
+{
+  public:
+    F (const unsigned int q) : Function<dim>(3), q(q) {}
+    
+    virtual void vector_value (const Point<dim> &p,
+                              Vector<double>   &v) const
+      {
+       for (unsigned int c=0; c<v.size(); ++c)
+         {
+           v(c) = 0;
+           for (unsigned int d=0; d<dim; ++d)
+             for (unsigned int i=0; i<=q; ++i)
+               v(c) += (d+1)*(i+1)*std::pow (p[d], 1.*i)+c;
+         }
+      }
+
+  private:
+    const unsigned int q;
+};
+
+
+
+template <int dim>
+void test ()
+{
+  Triangulation<dim>     triangulation;
+  GridGenerator::hyper_cube (triangulation);
+  triangulation.refine_global (3);
+
+  for (unsigned int p=1; p<6-dim; ++p)
+    {
+      FESystem<dim>              fe(FE_Q<dim>(p), 2,
+                                   FE_Q<dim>(p+1), 1);
+      DoFHandler<dim>        dof_handler(triangulation);
+      dof_handler.distribute_dofs (fe);
+
+      Vector<double> interpolant (dof_handler.n_dofs());
+      Vector<float>  error (triangulation.n_active_cells());
+      for (unsigned int q=0; q<=p+2; ++q)
+       {
+                                          // interpolate the function
+         VectorTools::interpolate (dof_handler,
+                                   F<dim> (q),
+                                   interpolant);
+      
+                                          // then compute the interpolation error
+         VectorTools::integrate_difference (dof_handler,
+                                            interpolant,
+                                            F<dim> (q),
+                                            error,
+                                            QGauss<dim>(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_q_system_01/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_q_system_01/cmp/generic b/tests/deal.II/interpolate_q_system_01/cmp/generic
new file mode 100644 (file)
index 0000000..35cf8fe
--- /dev/null
@@ -0,0 +1,47 @@
+
+DEAL::FESystem<1>[FE_Q<1>(1)^2-FE_Q<1>(2)], P_0, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(1)^2-FE_Q<1>(2)], P_1, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(1)^2-FE_Q<1>(2)], P_2, rel. error=0.000278
+DEAL::FESystem<1>[FE_Q<1>(1)^2-FE_Q<1>(2)], P_3, rel. error=0.000654
+DEAL::FESystem<1>[FE_Q<1>(2)^2-FE_Q<1>(3)], P_0, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(2)^2-FE_Q<1>(3)], P_1, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(2)^2-FE_Q<1>(3)], P_2, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(2)^2-FE_Q<1>(3)], P_3, rel. error=5.15e-06
+DEAL::FESystem<1>[FE_Q<1>(2)^2-FE_Q<1>(3)], P_4, rel. error=1.52e-05
+DEAL::FESystem<1>[FE_Q<1>(3)^2-FE_Q<1>(4)], P_0, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(3)^2-FE_Q<1>(4)], P_1, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(3)^2-FE_Q<1>(4)], P_2, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(3)^2-FE_Q<1>(4)], P_3, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(3)^2-FE_Q<1>(4)], P_4, rel. error=1.20e-07
+DEAL::FESystem<1>[FE_Q<1>(3)^2-FE_Q<1>(4)], P_5, rel. error=4.27e-07
+DEAL::FESystem<1>[FE_Q<1>(4)^2-FE_Q<1>(5)], P_0, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(4)^2-FE_Q<1>(5)], P_1, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(4)^2-FE_Q<1>(5)], P_2, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(4)^2-FE_Q<1>(5)], P_3, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(4)^2-FE_Q<1>(5)], P_4, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(4)^2-FE_Q<1>(5)], P_5, rel. error=3.21e-09
+DEAL::FESystem<1>[FE_Q<1>(4)^2-FE_Q<1>(5)], P_6, rel. error=1.33e-08
+DEAL::FESystem<2>[FE_Q<2>(1)^2-FE_Q<2>(2)], P_0, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(1)^2-FE_Q<2>(2)], P_1, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(1)^2-FE_Q<2>(2)], P_2, rel. error=8.74e-05
+DEAL::FESystem<2>[FE_Q<2>(1)^2-FE_Q<2>(2)], P_3, rel. error=0.000201
+DEAL::FESystem<2>[FE_Q<2>(2)^2-FE_Q<2>(3)], P_0, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(2)^2-FE_Q<2>(3)], P_1, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(2)^2-FE_Q<2>(3)], P_2, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(2)^2-FE_Q<2>(3)], P_3, rel. error=1.02e-06
+DEAL::FESystem<2>[FE_Q<2>(2)^2-FE_Q<2>(3)], P_4, rel. error=3.01e-06
+DEAL::FESystem<2>[FE_Q<2>(3)^2-FE_Q<2>(4)], P_0, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(3)^2-FE_Q<2>(4)], P_1, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(3)^2-FE_Q<2>(4)], P_2, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(3)^2-FE_Q<2>(4)], P_3, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(3)^2-FE_Q<2>(4)], P_4, rel. error=2.22e-08
+DEAL::FESystem<2>[FE_Q<2>(3)^2-FE_Q<2>(4)], P_5, rel. error=7.82e-08
+DEAL::FESystem<3>[FE_Q<3>(1)^2-FE_Q<3>(2)], P_0, rel. error=0
+DEAL::FESystem<3>[FE_Q<3>(1)^2-FE_Q<3>(2)], P_1, rel. error=0
+DEAL::FESystem<3>[FE_Q<3>(1)^2-FE_Q<3>(2)], P_2, rel. error=2.52e-05
+DEAL::FESystem<3>[FE_Q<3>(1)^2-FE_Q<3>(2)], P_3, rel. error=5.74e-05
+DEAL::FESystem<3>[FE_Q<3>(2)^2-FE_Q<3>(3)], P_0, rel. error=0
+DEAL::FESystem<3>[FE_Q<3>(2)^2-FE_Q<3>(3)], P_1, rel. error=0
+DEAL::FESystem<3>[FE_Q<3>(2)^2-FE_Q<3>(3)], P_2, rel. error=0
+DEAL::FESystem<3>[FE_Q<3>(2)^2-FE_Q<3>(3)], P_3, rel. error=2.04e-07
+DEAL::FESystem<3>[FE_Q<3>(2)^2-FE_Q<3>(3)], P_4, rel. error=6.07e-07
diff --git a/tests/deal.II/interpolate_q_system_02.cc b/tests/deal.II/interpolate_q_system_02.cc
new file mode 100644 (file)
index 0000000..d89dcf6
--- /dev/null
@@ -0,0 +1,131 @@
+//----------------------------  interpolate_q_system_01.cc  ---------------------------
+//    $Id: interpolate_q_system_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_q_system_02.cc  ---------------------------
+
+
+// check that VectorTools::interpolate works for FE_System(FE_Q(p)) elements correctly on
+// an adaptively refined mesh for functions of degree q
+
+#include "../tests.h"
+#include <base/function.h>
+#include <base/logstream.h>
+#include <base/quadrature_lib.h>
+#include <lac/vector.h>
+
+#include <grid/tria.h>
+#include <dofs/dof_handler.h>
+#include <dofs/dof_tools.h>
+#include <dofs/dof_constraints.h>
+#include <grid/grid_generator.h>
+#include <grid/grid_refinement.h>
+#include <grid/tria_accessor.h>
+#include <grid/tria_iterator.h>
+#include <grid/tria_boundary_lib.h>
+#include <dofs/dof_accessor.h>
+#include <dofs/dof_tools.h>
+#include <numerics/vectors.h>
+#include <fe/fe_q.h>
+#include <fe/fe_system.h>
+
+#include <fstream>
+#include <vector>
+
+
+template <int dim>
+class F :  public Function<dim>
+{
+  public:
+    F (const unsigned int q) : Function<dim>(3), q(q) {}
+    
+    virtual void vector_value (const Point<dim> &p,
+                              Vector<double>   &v) const
+      {
+       for (unsigned int c=0; c<v.size(); ++c)
+         {
+           v(c) = 0;
+           for (unsigned int d=0; d<dim; ++d)
+             for (unsigned int i=0; i<=q; ++i)
+               v(c) += (d+1)*(i+1)*std::pow (p[d], 1.*i)+c;
+         }
+      }
+
+  private:
+    const unsigned int q;
+};
+
+
+
+template <int dim>
+void test ()
+{
+  Triangulation<dim>     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)
+    {
+      FESystem<dim>              fe(FE_Q<dim>(p), 2,
+                                   FE_Q<dim>(p+1), 1);
+      DoFHandler<dim>        dof_handler(triangulation);
+      dof_handler.distribute_dofs (fe);
+
+      ConstraintMatrix constraints;
+      DoFTools::make_hanging_node_constraints (dof_handler, constraints);
+      constraints.close ();
+
+      Vector<double> interpolant (dof_handler.n_dofs());
+      Vector<float>  error (triangulation.n_active_cells());
+      for (unsigned int q=0; q<=p+2; ++q)
+       {
+                                          // interpolate the function
+         VectorTools::interpolate (dof_handler,
+                                   F<dim> (q),
+                                   interpolant);
+         constraints.distribute (interpolant);
+      
+                                          // then compute the interpolation error
+         VectorTools::integrate_difference (dof_handler,
+                                            interpolant,
+                                            F<dim> (q),
+                                            error,
+                                            QGauss<dim>(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_q_system_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_q_system_02/cmp/generic b/tests/deal.II/interpolate_q_system_02/cmp/generic
new file mode 100644 (file)
index 0000000..1cc5567
--- /dev/null
@@ -0,0 +1,47 @@
+
+DEAL::FESystem<1>[FE_Q<1>(1)^2-FE_Q<1>(2)], P_0, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(1)^2-FE_Q<1>(2)], P_1, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(1)^2-FE_Q<1>(2)], P_2, rel. error=0.000970
+DEAL::FESystem<1>[FE_Q<1>(1)^2-FE_Q<1>(2)], P_3, rel. error=0.00282
+DEAL::FESystem<1>[FE_Q<1>(2)^2-FE_Q<1>(3)], P_0, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(2)^2-FE_Q<1>(3)], P_1, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(2)^2-FE_Q<1>(3)], P_2, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(2)^2-FE_Q<1>(3)], P_3, rel. error=3.59e-05
+DEAL::FESystem<1>[FE_Q<1>(2)^2-FE_Q<1>(3)], P_4, rel. error=0.000134
+DEAL::FESystem<1>[FE_Q<1>(3)^2-FE_Q<1>(4)], P_0, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(3)^2-FE_Q<1>(4)], P_1, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(3)^2-FE_Q<1>(4)], P_2, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(3)^2-FE_Q<1>(4)], P_3, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(3)^2-FE_Q<1>(4)], P_4, rel. error=1.69e-06
+DEAL::FESystem<1>[FE_Q<1>(3)^2-FE_Q<1>(4)], P_5, rel. error=7.67e-06
+DEAL::FESystem<1>[FE_Q<1>(4)^2-FE_Q<1>(5)], P_0, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(4)^2-FE_Q<1>(5)], P_1, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(4)^2-FE_Q<1>(5)], P_2, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(4)^2-FE_Q<1>(5)], P_3, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(4)^2-FE_Q<1>(5)], P_4, rel. error=0
+DEAL::FESystem<1>[FE_Q<1>(4)^2-FE_Q<1>(5)], P_5, rel. error=9.09e-08
+DEAL::FESystem<1>[FE_Q<1>(4)^2-FE_Q<1>(5)], P_6, rel. error=4.83e-07
+DEAL::FESystem<2>[FE_Q<2>(1)^2-FE_Q<2>(2)], P_0, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(1)^2-FE_Q<2>(2)], P_1, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(1)^2-FE_Q<2>(2)], P_2, rel. error=0.000478
+DEAL::FESystem<2>[FE_Q<2>(1)^2-FE_Q<2>(2)], P_3, rel. error=0.00120
+DEAL::FESystem<2>[FE_Q<2>(2)^2-FE_Q<2>(3)], P_0, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(2)^2-FE_Q<2>(3)], P_1, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(2)^2-FE_Q<2>(3)], P_2, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(2)^2-FE_Q<2>(3)], P_3, rel. error=1.14e-05
+DEAL::FESystem<2>[FE_Q<2>(2)^2-FE_Q<2>(3)], P_4, rel. error=3.74e-05
+DEAL::FESystem<2>[FE_Q<2>(3)^2-FE_Q<2>(4)], P_0, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(3)^2-FE_Q<2>(4)], P_1, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(3)^2-FE_Q<2>(4)], P_2, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(3)^2-FE_Q<2>(4)], P_3, rel. error=0
+DEAL::FESystem<2>[FE_Q<2>(3)^2-FE_Q<2>(4)], P_4, rel. error=5.13e-07
+DEAL::FESystem<2>[FE_Q<2>(3)^2-FE_Q<2>(4)], P_5, rel. error=1.99e-06
+DEAL::FESystem<3>[FE_Q<3>(1)^2-FE_Q<3>(2)], P_0, rel. error=0
+DEAL::FESystem<3>[FE_Q<3>(1)^2-FE_Q<3>(2)], P_1, rel. error=0
+DEAL::FESystem<3>[FE_Q<3>(1)^2-FE_Q<3>(2)], P_2, rel. error=0.000199
+DEAL::FESystem<3>[FE_Q<3>(1)^2-FE_Q<3>(2)], P_3, rel. error=0.000477
+DEAL::FESystem<3>[FE_Q<3>(2)^2-FE_Q<3>(3)], P_0, rel. error=0
+DEAL::FESystem<3>[FE_Q<3>(2)^2-FE_Q<3>(3)], P_1, rel. error=0
+DEAL::FESystem<3>[FE_Q<3>(2)^2-FE_Q<3>(3)], P_2, rel. error=0
+DEAL::FESystem<3>[FE_Q<3>(2)^2-FE_Q<3>(3)], P_3, rel. error=3.39e-06
+DEAL::FESystem<3>[FE_Q<3>(2)^2-FE_Q<3>(3)], P_4, rel. error=1.06e-05

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.