]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Tests for DoFRenumbering::hierarchical().
authorheister <heister@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 24 Aug 2011 18:24:25 +0000 (18:24 +0000)
committerheister <heister@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 24 Aug 2011 18:24:25 +0000 (18:24 +0000)
git-svn-id: https://svn.dealii.org/trunk@24180 0785d39b-7218-0410-832d-ea1e28bc413d

tests/deal.II/dof_renumbering_zorder_01.cc [new file with mode: 0644]
tests/deal.II/dof_renumbering_zorder_01/cmp/generic [new file with mode: 0644]
tests/deal.II/dof_renumbering_zorder_02.cc [new file with mode: 0644]
tests/deal.II/dof_renumbering_zorder_02/cmp/generic [new file with mode: 0644]

diff --git a/tests/deal.II/dof_renumbering_zorder_01.cc b/tests/deal.II/dof_renumbering_zorder_01.cc
new file mode 100644 (file)
index 0000000..fac2b15
--- /dev/null
@@ -0,0 +1,113 @@
+//----------------------------------------------------------------------
+//    $Id: dof_renumbering_02.cc 23710 2011-05-17 04:50:10Z bangerth $
+//    Version: $Name$
+//
+//    Copyright (C) 2000, 2001, 2003, 2004, 2007, 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.
+//
+//----------------------------------------------------------------------
+
+// Check DoFRenumbering::hierachical changes nothing for a regular refined mesh
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/function_lib.h>
+#include <deal.II/lac/vector.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/dofs/dof_renumbering.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_dgq.h>
+#include <deal.II/fe/fe_dgp.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/fe_values.h>
+
+#include <fstream>
+#include <sstream>
+
+
+template <int dim, class stream>
+void
+print_dofs (const DoFHandler<dim> &dof, stream & out)
+{
+  const FiniteElement<dim>& fe = dof.get_fe();
+  std::vector<unsigned int> v (fe.dofs_per_cell);
+  std_cxx1x::shared_ptr<FEValues<dim> > fevalues;
+
+  if (fe.has_support_points())
+    {
+      Quadrature<dim> quad(fe.get_unit_support_points());
+      fevalues = std_cxx1x::shared_ptr<FEValues<dim> >(new FEValues<dim>(fe, quad, update_q_points));
+    }
+
+  for (typename DoFHandler<dim>::active_cell_iterator cell=dof.begin_active();
+       cell != dof.end(); ++cell)
+    {
+      Point<dim> p = cell->center();
+      if (fevalues.get() != 0)
+       fevalues->reinit(cell);
+
+      cell->get_dof_indices (v);
+      for (unsigned int i=0; i<v.size(); ++i)
+       if (fevalues.get() != 0)
+         out << fevalues->quadrature_point(i) << '\t' << v[i] << std::endl;
+       else
+         out << p << '\t' << v[i] << std::endl;
+      out << std::endl;
+    }
+}
+
+
+
+
+template <int dim>
+void
+check ()
+{
+  Triangulation<dim> tr;
+  GridGenerator::hyper_cube(tr, -1., 1.);
+  tr.refine_global (1);
+
+  FE_Q<dim> fe(1);
+  DoFHandler<dim> dof(tr);
+  dof.distribute_dofs(fe);
+
+  std::ostringstream o1, o2;
+  
+  print_dofs(dof, o1);
+  deallog << "**" << endl;
+
+  DoFRenumbering::hierarchical(dof);
+
+  print_dofs(dof, o2);
+
+  if (o1.str()==o2.str())
+       deallog << "OK" << endl;
+  
+}
+
+int main ()
+{
+  std::ofstream logfile ("dof_renumbering_zorder_01/output");
+  deallog << std::setprecision (2);
+  deallog << std::fixed;
+  deallog.attach(logfile);
+  deallog.depth_console (0);
+
+  deallog.push ("1d");
+  check<1> ();
+  deallog.pop ();
+  deallog.push ("2d");
+  check<2> ();
+  deallog.pop ();
+  deallog.push ("3d");
+  check<3> ();
+  deallog.pop ();
+}
diff --git a/tests/deal.II/dof_renumbering_zorder_01/cmp/generic b/tests/deal.II/dof_renumbering_zorder_01/cmp/generic
new file mode 100644 (file)
index 0000000..16b2558
--- /dev/null
@@ -0,0 +1,7 @@
+
+DEAL:1d::**
+DEAL:1d::OK
+DEAL:2d::**
+DEAL:2d::OK
+DEAL:3d::**
+DEAL:3d::OK
diff --git a/tests/deal.II/dof_renumbering_zorder_02.cc b/tests/deal.II/dof_renumbering_zorder_02.cc
new file mode 100644 (file)
index 0000000..bb6af0d
--- /dev/null
@@ -0,0 +1,131 @@
+//----------------------------------------------------------------------
+//    $Id: dof_renumbering_02.cc 23710 2011-05-17 04:50:10Z bangerth $
+//    Version: $Name$
+//
+//    Copyright (C) 2000, 2001, 2003, 2004, 2007, 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.
+//
+//----------------------------------------------------------------------
+
+// Check that DoFRenumbering::hierarchical works in the most simple case.
+// checked by hand.
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/function_lib.h>
+#include <deal.II/lac/vector.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/dofs/dof_renumbering.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_dgq.h>
+#include <deal.II/fe/fe_dgp.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/fe_values.h>
+
+#include <fstream>
+#include <sstream>
+
+
+template <int dim, class stream>
+void
+print_dofs (const DoFHandler<dim> &dof, stream & out)
+{
+  out << std::setprecision (2);
+  out << std::fixed;
+  const FiniteElement<dim>& fe = dof.get_fe();
+  std::vector<unsigned int> v (fe.dofs_per_cell);
+  std_cxx1x::shared_ptr<FEValues<dim> > fevalues;
+
+  if (fe.has_support_points())
+    {
+      Quadrature<dim> quad(fe.get_unit_support_points());
+      fevalues = std_cxx1x::shared_ptr<FEValues<dim> >(new FEValues<dim>(fe, quad, update_q_points));
+    }
+
+  for (typename DoFHandler<dim>::active_cell_iterator cell=dof.begin_active();
+       cell != dof.end(); ++cell)
+    {
+      Point<dim> p = cell->center();
+      if (fevalues.get() != 0)
+       fevalues->reinit(cell);
+
+      cell->get_dof_indices (v);
+      for (unsigned int i=0; i<v.size(); ++i)
+       if (fevalues.get() != 0)
+         out << fevalues->quadrature_point(i) << '\t' << v[i] << std::endl;
+       else
+         out << p << '\t' << v[i] << std::endl;
+    }
+}
+
+
+
+
+template <int dim>
+void
+check ()
+{
+  Triangulation<dim> tr;
+  GridGenerator::hyper_cube(tr, -1., 1.);
+  tr.refine_global (1);
+  tr.begin_active()->set_refine_flag ();
+  tr.execute_coarsening_and_refinement ();
+
+  FE_Q<dim> fe(1);
+  DoFHandler<dim> dof(tr);
+  dof.distribute_dofs(fe);
+
+  std::ostringstream o1, o2, o3;
+  
+  print_dofs(dof, o1);
+  deallog << o1.str();
+  deallog << "**" << endl;
+
+  DoFRenumbering::hierarchical(dof);
+
+  print_dofs(dof, o2);
+  deallog << o2.str();
+
+  if (o1.str()!=o2.str())
+       deallog << "OK" << endl;
+  else
+       Assert(false, ExcInternalError());
+
+  DoFRenumbering::hierarchical(dof);
+  print_dofs(dof, o3);
+  
+
+  // doing renumbering twice does not change the result?!
+  if (o2.str()==o3.str())
+       deallog << "OK" << endl;
+  else
+       Assert(false, ExcInternalError());
+
+}
+
+int main ()
+{
+  std::ofstream logfile ("dof_renumbering_zorder_02/output");
+  deallog << std::setprecision (2);
+  deallog << std::fixed;
+  deallog.attach(logfile);
+  deallog.depth_console (0);
+
+  deallog.push ("1d");
+  check<1> ();
+  deallog.pop ();
+  deallog.push ("2d");
+  check<2> ();
+  deallog.pop ();
+  deallog.push ("3d");
+  check<3> ();
+  deallog.pop ();
+}
diff --git a/tests/deal.II/dof_renumbering_zorder_02/cmp/generic b/tests/deal.II/dof_renumbering_zorder_02/cmp/generic
new file mode 100644 (file)
index 0000000..011277e
--- /dev/null
@@ -0,0 +1,318 @@
+
+DEAL:1d::0.00  0
+1.00   1
+-1.00  2
+-0.50  3
+-0.50  3
+0.00   0
+**
+DEAL:1d::0.00  2
+1.00   3
+-1.00  0
+-0.50  1
+-0.50  1
+0.00   2
+OK
+DEAL:1d::OK
+DEAL:2d::0.00 -1.00    0
+1.00 -1.00     1
+0.00 0.00      2
+1.00 0.00      3
+-1.00 0.00     4
+0.00 0.00      2
+-1.00 1.00     5
+0.00 1.00      6
+0.00 0.00      2
+1.00 0.00      3
+0.00 1.00      6
+1.00 1.00      7
+-1.00 -1.00    8
+-0.50 -1.00    9
+-1.00 -0.50    10
+-0.50 -0.50    11
+-0.50 -1.00    9
+0.00 -1.00     0
+-0.50 -0.50    11
+0.00 -0.50     12
+-1.00 -0.50    10
+-0.50 -0.50    11
+-1.00 0.00     4
+-0.50 0.00     13
+-0.50 -0.50    11
+0.00 -0.50     12
+-0.50 0.00     13
+0.00 0.00      2
+**
+DEAL:2d::0.00 -1.00    4
+1.00 -1.00     9
+0.00 0.00      8
+1.00 0.00      10
+-1.00 0.00     6
+0.00 0.00      8
+-1.00 1.00     11
+0.00 1.00      12
+0.00 0.00      8
+1.00 0.00      10
+0.00 1.00      12
+1.00 1.00      13
+-1.00 -1.00    0
+-0.50 -1.00    1
+-1.00 -0.50    2
+-0.50 -0.50    3
+-0.50 -1.00    1
+0.00 -1.00     4
+-0.50 -0.50    3
+0.00 -0.50     5
+-1.00 -0.50    2
+-0.50 -0.50    3
+-1.00 0.00     6
+-0.50 0.00     7
+-0.50 -0.50    3
+0.00 -0.50     5
+-0.50 0.00     7
+0.00 0.00      8
+OK
+DEAL:2d::OK
+DEAL:3d::0.00 -1.00 -1.00      0
+1.00 -1.00 -1.00       1
+0.00 0.00 -1.00        2
+1.00 0.00 -1.00        3
+0.00 -1.00 0.00        4
+1.00 -1.00 0.00        5
+0.00 0.00 0.00 6
+1.00 0.00 0.00 7
+-1.00 0.00 -1.00       8
+0.00 0.00 -1.00        2
+-1.00 1.00 -1.00       9
+0.00 1.00 -1.00        10
+-1.00 0.00 0.00        11
+0.00 0.00 0.00 6
+-1.00 1.00 0.00        12
+0.00 1.00 0.00 13
+0.00 0.00 -1.00        2
+1.00 0.00 -1.00        3
+0.00 1.00 -1.00        10
+1.00 1.00 -1.00        14
+0.00 0.00 0.00 6
+1.00 0.00 0.00 7
+0.00 1.00 0.00 13
+1.00 1.00 0.00 15
+-1.00 -1.00 0.00       16
+0.00 -1.00 0.00        4
+-1.00 0.00 0.00        11
+0.00 0.00 0.00 6
+-1.00 -1.00 1.00       17
+0.00 -1.00 1.00        18
+-1.00 0.00 1.00        19
+0.00 0.00 1.00 20
+0.00 -1.00 0.00        4
+1.00 -1.00 0.00        5
+0.00 0.00 0.00 6
+1.00 0.00 0.00 7
+0.00 -1.00 1.00        18
+1.00 -1.00 1.00        21
+0.00 0.00 1.00 20
+1.00 0.00 1.00 22
+-1.00 0.00 0.00        11
+0.00 0.00 0.00 6
+-1.00 1.00 0.00        12
+0.00 1.00 0.00 13
+-1.00 0.00 1.00        19
+0.00 0.00 1.00 20
+-1.00 1.00 1.00        23
+0.00 1.00 1.00 24
+0.00 0.00 0.00 6
+1.00 0.00 0.00 7
+0.00 1.00 0.00 13
+1.00 1.00 0.00 15
+0.00 0.00 1.00 20
+1.00 0.00 1.00 22
+0.00 1.00 1.00 24
+1.00 1.00 1.00 25
+-1.00 -1.00 -1.00      26
+-0.50 -1.00 -1.00      27
+-1.00 -0.50 -1.00      28
+-0.50 -0.50 -1.00      29
+-1.00 -1.00 -0.50      30
+-0.50 -1.00 -0.50      31
+-1.00 -0.50 -0.50      32
+-0.50 -0.50 -0.50      33
+-0.50 -1.00 -1.00      27
+0.00 -1.00 -1.00       0
+-0.50 -0.50 -1.00      29
+0.00 -0.50 -1.00       34
+-0.50 -1.00 -0.50      31
+0.00 -1.00 -0.50       35
+-0.50 -0.50 -0.50      33
+0.00 -0.50 -0.50       36
+-1.00 -0.50 -1.00      28
+-0.50 -0.50 -1.00      29
+-1.00 0.00 -1.00       8
+-0.50 0.00 -1.00       37
+-1.00 -0.50 -0.50      32
+-0.50 -0.50 -0.50      33
+-1.00 0.00 -0.50       38
+-0.50 0.00 -0.50       39
+-0.50 -0.50 -1.00      29
+0.00 -0.50 -1.00       34
+-0.50 0.00 -1.00       37
+0.00 0.00 -1.00        2
+-0.50 -0.50 -0.50      33
+0.00 -0.50 -0.50       36
+-0.50 0.00 -0.50       39
+0.00 0.00 -0.50        40
+-1.00 -1.00 -0.50      30
+-0.50 -1.00 -0.50      31
+-1.00 -0.50 -0.50      32
+-0.50 -0.50 -0.50      33
+-1.00 -1.00 0.00       16
+-0.50 -1.00 0.00       41
+-1.00 -0.50 0.00       42
+-0.50 -0.50 0.00       43
+-0.50 -1.00 -0.50      31
+0.00 -1.00 -0.50       35
+-0.50 -0.50 -0.50      33
+0.00 -0.50 -0.50       36
+-0.50 -1.00 0.00       41
+0.00 -1.00 0.00        4
+-0.50 -0.50 0.00       43
+0.00 -0.50 0.00        44
+-1.00 -0.50 -0.50      32
+-0.50 -0.50 -0.50      33
+-1.00 0.00 -0.50       38
+-0.50 0.00 -0.50       39
+-1.00 -0.50 0.00       42
+-0.50 -0.50 0.00       43
+-1.00 0.00 0.00        11
+-0.50 0.00 0.00        45
+-0.50 -0.50 -0.50      33
+0.00 -0.50 -0.50       36
+-0.50 0.00 -0.50       39
+0.00 0.00 -0.50        40
+-0.50 -0.50 0.00       43
+0.00 -0.50 0.00        44
+-0.50 0.00 0.00        45
+0.00 0.00 0.00 6
+**
+DEAL:3d::0.00 -1.00 -1.00      8
+1.00 -1.00 -1.00       27
+0.00 0.00 -1.00        16
+1.00 0.00 -1.00        28
+0.00 -1.00 0.00        22
+1.00 -1.00 0.00        29
+0.00 0.00 0.00 26
+1.00 0.00 0.00 30
+-1.00 0.00 -1.00       12
+0.00 0.00 -1.00        16
+-1.00 1.00 -1.00       31
+0.00 1.00 -1.00        32
+-1.00 0.00 0.00        24
+0.00 0.00 0.00 26
+-1.00 1.00 0.00        33
+0.00 1.00 0.00 34
+0.00 0.00 -1.00        16
+1.00 0.00 -1.00        28
+0.00 1.00 -1.00        32
+1.00 1.00 -1.00        35
+0.00 0.00 0.00 26
+1.00 0.00 0.00 30
+0.00 1.00 0.00 34
+1.00 1.00 0.00 36
+-1.00 -1.00 0.00       18
+0.00 -1.00 0.00        22
+-1.00 0.00 0.00        24
+0.00 0.00 0.00 26
+-1.00 -1.00 1.00       37
+0.00 -1.00 1.00        38
+-1.00 0.00 1.00        39
+0.00 0.00 1.00 40
+0.00 -1.00 0.00        22
+1.00 -1.00 0.00        29
+0.00 0.00 0.00 26
+1.00 0.00 0.00 30
+0.00 -1.00 1.00        38
+1.00 -1.00 1.00        41
+0.00 0.00 1.00 40
+1.00 0.00 1.00 42
+-1.00 0.00 0.00        24
+0.00 0.00 0.00 26
+-1.00 1.00 0.00        33
+0.00 1.00 0.00 34
+-1.00 0.00 1.00        39
+0.00 0.00 1.00 40
+-1.00 1.00 1.00        43
+0.00 1.00 1.00 44
+0.00 0.00 0.00 26
+1.00 0.00 0.00 30
+0.00 1.00 0.00 34
+1.00 1.00 0.00 36
+0.00 0.00 1.00 40
+1.00 0.00 1.00 42
+0.00 1.00 1.00 44
+1.00 1.00 1.00 45
+-1.00 -1.00 -1.00      0
+-0.50 -1.00 -1.00      1
+-1.00 -0.50 -1.00      2
+-0.50 -0.50 -1.00      3
+-1.00 -1.00 -0.50      4
+-0.50 -1.00 -0.50      5
+-1.00 -0.50 -0.50      6
+-0.50 -0.50 -0.50      7
+-0.50 -1.00 -1.00      1
+0.00 -1.00 -1.00       8
+-0.50 -0.50 -1.00      3
+0.00 -0.50 -1.00       9
+-0.50 -1.00 -0.50      5
+0.00 -1.00 -0.50       10
+-0.50 -0.50 -0.50      7
+0.00 -0.50 -0.50       11
+-1.00 -0.50 -1.00      2
+-0.50 -0.50 -1.00      3
+-1.00 0.00 -1.00       12
+-0.50 0.00 -1.00       13
+-1.00 -0.50 -0.50      6
+-0.50 -0.50 -0.50      7
+-1.00 0.00 -0.50       14
+-0.50 0.00 -0.50       15
+-0.50 -0.50 -1.00      3
+0.00 -0.50 -1.00       9
+-0.50 0.00 -1.00       13
+0.00 0.00 -1.00        16
+-0.50 -0.50 -0.50      7
+0.00 -0.50 -0.50       11
+-0.50 0.00 -0.50       15
+0.00 0.00 -0.50        17
+-1.00 -1.00 -0.50      4
+-0.50 -1.00 -0.50      5
+-1.00 -0.50 -0.50      6
+-0.50 -0.50 -0.50      7
+-1.00 -1.00 0.00       18
+-0.50 -1.00 0.00       19
+-1.00 -0.50 0.00       20
+-0.50 -0.50 0.00       21
+-0.50 -1.00 -0.50      5
+0.00 -1.00 -0.50       10
+-0.50 -0.50 -0.50      7
+0.00 -0.50 -0.50       11
+-0.50 -1.00 0.00       19
+0.00 -1.00 0.00        22
+-0.50 -0.50 0.00       21
+0.00 -0.50 0.00        23
+-1.00 -0.50 -0.50      6
+-0.50 -0.50 -0.50      7
+-1.00 0.00 -0.50       14
+-0.50 0.00 -0.50       15
+-1.00 -0.50 0.00       20
+-0.50 -0.50 0.00       21
+-1.00 0.00 0.00        24
+-0.50 0.00 0.00        25
+-0.50 -0.50 -0.50      7
+0.00 -0.50 -0.50       11
+-0.50 0.00 -0.50       15
+0.00 0.00 -0.50        17
+-0.50 -0.50 0.00       21
+0.00 -0.50 0.00        23
+-0.50 0.00 0.00        25
+0.00 0.00 0.00 26
+OK
+DEAL:3d::OK

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.