]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add 3D test for vertex_to_cell_map function. 1694/head
authorBruno Turcksin <bruno.turcksin@gmail.com>
Wed, 30 Sep 2015 23:28:05 +0000 (18:28 -0500)
committerBruno Turcksin <bruno.turcksin@gmail.com>
Wed, 30 Sep 2015 23:28:05 +0000 (18:28 -0500)
tests/grid/grid_tools_vertex_to_cell_map_01.cc [moved from tests/grid/grid_tools_vertex_to_cell_map.cc with 97% similarity]
tests/grid/grid_tools_vertex_to_cell_map_01.mpirun=2.output [moved from tests/grid/grid_tools_vertex_to_cell_map.mpirun=2.output with 100% similarity]
tests/grid/grid_tools_vertex_to_cell_map_02.cc [new file with mode: 0644]
tests/grid/grid_tools_vertex_to_cell_map_02.output [new file with mode: 0644]

similarity index 97%
rename from tests/grid/grid_tools_vertex_to_cell_map.cc
rename to tests/grid/grid_tools_vertex_to_cell_map_01.cc
index 3e66bf990297f9efaf0509fe34beab09bcc8605f..6eb8dfe385014e8141d565f84dfec0e261b133f3 100644 (file)
@@ -12,7 +12,7 @@
 //
 // ---------------------------------------------------------------------
 
-
+// Test vertex_to_cell_map for 2D problem with mpi and hanging nodes.
 
 #include "../tests.h"
 #include <deal.II/base/logstream.h>
diff --git a/tests/grid/grid_tools_vertex_to_cell_map_02.cc b/tests/grid/grid_tools_vertex_to_cell_map_02.cc
new file mode 100644 (file)
index 0000000..fdef460
--- /dev/null
@@ -0,0 +1,81 @@
+
+//
+// Copyright (C) 2015 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.
+//
+// ---------------------------------------------------------------------
+
+
+// Test vertex_to_cell_map for 3D problem with hanging nodes.
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/cell_id.h>
+
+#include <vector>
+
+std::ofstream logfile("output");
+
+void test()
+{
+  Triangulation<3> tria;
+  Point<3> a(0.,0.,0.);
+  Point<3> b(1.,1.,1.);
+  std::vector<unsigned int> repetitions(3);
+  repetitions[0] = 2;
+  repetitions[1] = 2;
+  repetitions[2] = 1;
+  GridGenerator::subdivided_hyper_rectangle(tria,repetitions,a,b);
+  typename Triangulation<3>::active_cell_iterator cell = tria.begin_active();
+  cell->set_refine_flag();
+  tria.execute_coarsening_and_refinement();
+
+  std::vector<std::set<typename Triangulation<3>::active_cell_iterator> >
+  vertex_to_cell = GridTools::vertex_to_cell_map(tria);
+
+  AssertThrow(tria.n_vertices()==vertex_to_cell.size(),
+              ExcMessage("Wrong number of vertices"));
+
+  std::vector<unsigned int> n_cells;
+  for (unsigned int i=0; i<vertex_to_cell.size(); ++i)
+    n_cells.push_back(vertex_to_cell[i].size());
+
+  std::vector<unsigned int> histogram(9,0);
+  for (unsigned int i=0; i<n_cells.size(); ++i)
+    histogram[n_cells[i]] += 1;
+
+  AssertThrow(histogram[0]==0, ExcMessage("Wrong cell distribution"));
+  AssertThrow(histogram[1]==8, ExcMessage("Wrong cell distribution"));
+  AssertThrow(histogram[2]==13, ExcMessage("Wrong cell distribution"));
+  AssertThrow(histogram[3]==6, ExcMessage("Wrong cell distribution"));
+  AssertThrow(histogram[4]==6, ExcMessage("Wrong cell distribution"));
+  AssertThrow(histogram[5]==3, ExcMessage("Wrong cell distribution"));
+  AssertThrow(histogram[6]==0, ExcMessage("Wrong cell distribution"));
+  AssertThrow(histogram[7]==0, ExcMessage("Wrong cell distribution"));
+  AssertThrow(histogram[8]==1, ExcMessage("Wrong cell distribution"));
+}
+
+int main (int argc, char *argv[])
+{
+  deallog << std::setprecision(4);
+  logfile << std::setprecision(4);
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+
+  test();
+
+  deallog << "OK" << std::endl;
+
+  return 0;
+}
diff --git a/tests/grid/grid_tools_vertex_to_cell_map_02.output b/tests/grid/grid_tools_vertex_to_cell_map_02.output
new file mode 100644 (file)
index 0000000..0fd8fc1
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL::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.