]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
add test for new DataOut assert
authorheister <heister@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 19 Feb 2014 15:01:51 +0000 (15:01 +0000)
committerheister <heister@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 19 Feb 2014 15:01:51 +0000 (15:01 +0000)
git-svn-id: https://svn.dealii.org/trunk@32515 0785d39b-7218-0410-832d-ea1e28bc413d

tests/bits/data_out_09.cc [new file with mode: 0644]
tests/bits/data_out_09.output [new file with mode: 0644]

diff --git a/tests/bits/data_out_09.cc b/tests/bits/data_out_09.cc
new file mode 100644 (file)
index 0000000..25712fa
--- /dev/null
@@ -0,0 +1,131 @@
+// ---------------------------------------------------------------------
+// $Id$
+//
+// Copyright (C) 2003 - 2013 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// check that add_data_vector and type_automatic work correctly and confirm
+// that the Assert is triggered that checks, when we have dof_data and
+// cell_data with the same length.
+
+
+#include "../tests.h"
+#include <deal.II/lac/sparsity_pattern.h>
+#include <deal.II/numerics/data_out.h>
+
+#include <deal.II/base/logstream.h>
+#include <deal.II/lac/vector.h>
+#include <deal.II/lac/block_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_tools.h>
+#include <deal.II/fe/fe_dgq.h>
+
+
+
+
+void test()
+{
+  Triangulation<2> tria;
+  GridGenerator::hyper_cube (tria);
+  tria.refine_global(1);
+  
+  FE_DGQ<2> fe(0);
+  DoFHandler<2> dof_handler(tria);
+  dof_handler.distribute_dofs(fe);
+
+  std::vector<types::global_dof_index> renumbering(dof_handler.n_dofs());
+  renumbering[0]=3;
+  renumbering[1]=2;
+  renumbering[2]=1;
+  renumbering[3]=0;
+  
+  dof_handler.renumber_dofs(renumbering);
+  
+
+  Vector<double> v_node(dof_handler.n_dofs());
+  Vector<double> v_cell(dof_handler.n_dofs());
+  for (unsigned int i=0;i<dof_handler.n_dofs();++i)
+    {
+      v_node[i]=i;
+      v_cell[i]=10+i;
+    }
+
+  //correct one is done if both are possible but specified
+  deallog << "*** Check cell and node data:" << std::endl;
+  {    
+    DataOut<2> data_out;
+    data_out.attach_dof_handler (dof_handler);
+    data_out.add_data_vector (v_node, "node_data", DataOut<2>::type_dof_data);
+    data_out.add_data_vector (v_cell, "cell_data", DataOut<2>::type_cell_data);
+    data_out.build_patches ();
+
+    data_out.write_vtk (deallog.get_file_stream());
+  }
+
+  //only tria, output correctly
+  deallog << "*** Check cell data only:" << std::endl;
+  {  
+    DataOut<2> data_out;
+    data_out.attach_triangulation (tria);
+    data_out.add_data_vector (v_cell, "cell_data");
+    data_out.build_patches ();
+
+    data_out.write_vtk (deallog.get_file_stream());
+  }
+  
+
+  //error if both
+  deallog << "*** should fail:" << std::endl;
+  try
+  {  
+    DataOut<2> data_out;
+    data_out.attach_dof_handler (dof_handler);
+    data_out.add_data_vector (v_cell, "cell_data");
+    data_out.build_patches ();
+
+    data_out.write_vtk (deallog.get_file_stream());
+  }
+  catch (...)
+    {
+      deallog << "exception" << std::endl;
+    }
+
+  // works if we use the other add_data_vector:
+  deallog << "*** Check node data only:" << std::endl;
+  {  
+    DataOut<2> data_out;
+    data_out.attach_dof_handler (dof_handler);
+    data_out.add_data_vector (dof_handler, v_node, "node_data");
+    data_out.build_patches ();
+
+    data_out.write_vtk (deallog.get_file_stream());
+  }
+  
+  
+  
+}
+
+
+int main()
+{
+  deal_II_exceptions::disable_abort_on_exception();
+  initlog();
+  test();
+  
+}
+
diff --git a/tests/bits/data_out_09.output b/tests/bits/data_out_09.output
new file mode 100644 (file)
index 0000000..37af927
--- /dev/null
@@ -0,0 +1,114 @@
+
+DEAL::*** Check cell and node data:
+# vtk DataFile Version 3.0
+#This file was generated 
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 16 double
+0.00000 0.00000 0
+0.500000 0.00000 0
+0.00000 0.500000 0
+0.500000 0.500000 0
+0.500000 0.00000 0
+1.00000 0.00000 0
+0.500000 0.500000 0
+1.00000 0.500000 0
+0.00000 0.500000 0
+0.500000 0.500000 0
+0.00000 1.00000 0
+0.500000 1.00000 0
+0.500000 0.500000 0
+1.00000 0.500000 0
+0.500000 1.00000 0
+1.00000 1.00000 0
+
+CELLS 4 20
+4      0       1       3       2
+4      4       5       7       6
+4      8       9       11      10
+4      12      13      15      14
+
+CELL_TYPES 4
+ 9 9 9 9
+POINT_DATA 16
+SCALARS node_data double 1
+LOOKUP_TABLE default
+3.00000 3.00000 3.00000 3.00000 2.00000 2.00000 2.00000 2.00000 1.00000 1.00000 1.00000 1.00000 0.00000 0.00000 0.00000 0.00000 
+SCALARS cell_data double 1
+LOOKUP_TABLE default
+10.0000 10.0000 10.0000 10.0000 11.0000 11.0000 11.0000 11.0000 12.0000 12.0000 12.0000 12.0000 13.0000 13.0000 13.0000 13.0000 
+DEAL::*** Check cell data only:
+# vtk DataFile Version 3.0
+#This file was generated 
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 16 double
+0.00000 0.00000 0
+0.500000 0.00000 0
+0.00000 0.500000 0
+0.500000 0.500000 0
+0.500000 0.00000 0
+1.00000 0.00000 0
+0.500000 0.500000 0
+1.00000 0.500000 0
+0.00000 0.500000 0
+0.500000 0.500000 0
+0.00000 1.00000 0
+0.500000 1.00000 0
+0.500000 0.500000 0
+1.00000 0.500000 0
+0.500000 1.00000 0
+1.00000 1.00000 0
+
+CELLS 4 20
+4      0       1       3       2
+4      4       5       7       6
+4      8       9       11      10
+4      12      13      15      14
+
+CELL_TYPES 4
+ 9 9 9 9
+POINT_DATA 16
+SCALARS cell_data double 1
+LOOKUP_TABLE default
+10.0000 10.0000 10.0000 10.0000 11.0000 11.0000 11.0000 11.0000 12.0000 12.0000 12.0000 12.0000 13.0000 13.0000 13.0000 13.0000 
+DEAL::*** should fail:
+DEAL::exception
+DEAL::*** Check node data only:
+# vtk DataFile Version 3.0
+#This file was generated 
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 16 double
+0.00000 0.00000 0
+0.500000 0.00000 0
+0.00000 0.500000 0
+0.500000 0.500000 0
+0.500000 0.00000 0
+1.00000 0.00000 0
+0.500000 0.500000 0
+1.00000 0.500000 0
+0.00000 0.500000 0
+0.500000 0.500000 0
+0.00000 1.00000 0
+0.500000 1.00000 0
+0.500000 0.500000 0
+1.00000 0.500000 0
+0.500000 1.00000 0
+1.00000 1.00000 0
+
+CELLS 4 20
+4      0       1       3       2
+4      4       5       7       6
+4      8       9       11      10
+4      12      13      15      14
+
+CELL_TYPES 4
+ 9 9 9 9
+POINT_DATA 16
+SCALARS node_data double 1
+LOOKUP_TABLE default
+3.00000 3.00000 3.00000 3.00000 2.00000 2.00000 2.00000 2.00000 1.00000 1.00000 1.00000 1.00000 0.00000 0.00000 0.00000 0.00000 

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.