]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add test case mpi/tria_signals_05: tria_signals_02 -> p4est balance
authorLei Qiao <qiaol618@gmail.com>
Thu, 15 Oct 2015 17:08:03 +0000 (12:08 -0500)
committerLei Qiao <qiaol618@gmail.com>
Thu, 22 Oct 2015 14:56:37 +0000 (09:56 -0500)
tests/mpi/tria_signals_05.cc [new file with mode: 0644]
tests/mpi/tria_signals_05.mpirun=1.output [new file with mode: 0644]
tests/mpi/tria_signals_05.mpirun=11.output [new file with mode: 0644]
tests/mpi/tria_signals_05.mpirun=4.output [new file with mode: 0644]

diff --git a/tests/mpi/tria_signals_05.cc b/tests/mpi/tria_signals_05.cc
new file mode 100644 (file)
index 0000000..59d1daf
--- /dev/null
@@ -0,0 +1,169 @@
+// ---------------------------------------------------------------------
+//
+// 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.
+//
+// ---------------------------------------------------------------------
+
+
+#include "../tests.h"
+#include <deal.II/grid/tria.h>
+#include <deal.II/distributed/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/base/std_cxx11/bind.h>
+
+#include <fstream>
+#include <ostream>
+
+// Test on whether signals post_refinement_on_cell and pre_coarsening_on_cell
+// could catch all cell changes.
+// The test is designed to count cell number increase and decrease in signal
+// calls and then compare the result against n_active_cells reported by Tria
+// object. Absolute value change in n_active_cells is not concerned in this test.
+
+// This test is based on tria_signals_02. The difference is in this case we know
+// that p4est is doing mesh smoothing beyond class Triangulation. The case setup
+// is borrowed from tests/distributed_grids/2d_refinement_10.
+
+template<int dim, int spacedim>
+class SignalListener
+{
+public:
+  SignalListener(Triangulation<dim, spacedim> &tria_in)
+    :
+    n_active_cells(tria_in.n_active_cells()),
+    tria(tria_in)
+  {
+    tria_in.signals.post_refinement_on_cell.connect
+    (std_cxx11::bind (&SignalListener<dim, spacedim>::count_on_refine,
+                      this,
+                      std_cxx11::placeholders::_1));
+
+    tria_in.signals.pre_coarsening_on_cell.connect
+    (std_cxx11::bind (&SignalListener<dim, spacedim>::count_on_coarsen,
+                      this,
+                      std_cxx11::placeholders::_1));
+  }
+
+  int n_active_cell_gap()
+  {
+    return (n_active_cells -
+            static_cast<int> (tria.n_active_cells()));
+  }
+
+private:
+  void count_on_refine(const typename Triangulation<dim, spacedim>::cell_iterator &cell)
+  {
+    n_active_cells += cell->n_children();
+    --n_active_cells;
+
+    return;
+  }
+
+  void count_on_coarsen(const typename Triangulation<dim, spacedim>::cell_iterator &cell)
+  {
+    ++n_active_cells;
+    n_active_cells -= cell->n_children();
+
+    return;
+  }
+
+  int n_active_cells;
+  const Triangulation<dim, spacedim> &tria;
+};
+
+
+template<int dim, int spacedim>
+void test()
+{
+  typedef parallel::distributed::Triangulation<dim, spacedim> TriaType;
+
+  {
+    const std::string prefix = Utilities::int_to_string (dim, 1) +
+                               "d-" +
+                               Utilities::int_to_string (spacedim, 1)
+                               + "d";
+    deallog.push(prefix.c_str());
+  }
+
+  TriaType tria(MPI_COMM_WORLD);
+
+  {
+    std::vector<unsigned int> repetitions;
+    Point<dim> p1;
+    Point<dim> p2;
+
+    for (unsigned int d=0; d<dim; ++d)
+      {
+        repetitions.push_back (2);
+        p1[d] = 0.0;
+        p2[d] = 1.0;
+      }
+    GridGenerator::subdivided_hyper_rectangle (tria,
+                                               repetitions,
+                                               p1,
+                                               p2);
+  }
+
+  SignalListener<dim, spacedim> count_cell_via_signal(tria);
+
+  for (unsigned int n_loop = 1; n_loop < 5; ++n_loop)
+    {
+      {
+        Point<dim> p;
+        for (unsigned int d=0; d<dim; ++d)
+          {
+            p[d] = 0.5 - std::pow(0.5, 1.0 + n_loop);
+          }
+        typename TriaType::active_cell_iterator
+        cell = tria.begin_active();
+        for (; cell!= tria.end(); ++cell)
+          if (cell->is_locally_owned() && ((cell->center()).distance(p) < 1e-4))
+            {
+              cell->set_refine_flag();
+            }
+      }
+
+      tria.execute_coarsening_and_refinement ();
+
+      deallog << "n_loop: " << n_loop
+              << ", n_cell_gap: "
+              << count_cell_via_signal.n_active_cell_gap() << std::endl;
+    }
+
+  deallog.pop();
+  return;
+}
+
+int main(int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, /* int max_num_threads */ 1);
+  MPILogInitAll log;
+
+  // parallel::distributed::Triangulation<1, spacedim> is not valid.
+  {
+    const int dim = 2;
+    const int spacedim = 2;
+    test<dim,spacedim> ();
+  }
+
+  // GridGenerator::subdivided_hyper_rectangle do not accept
+  // parallel::distributed::Triangulation<2, 3>.
+
+  {
+    const int dim = 3;
+    const int spacedim = 3;
+    test<dim,spacedim> ();
+  }
+
+  return (0);
+}
+
diff --git a/tests/mpi/tria_signals_05.mpirun=1.output b/tests/mpi/tria_signals_05.mpirun=1.output
new file mode 100644 (file)
index 0000000..3f164b5
--- /dev/null
@@ -0,0 +1,9 @@
+
+DEAL:0:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:0:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:0:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:0:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:0:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:0:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:0:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:0:3d-3d::n_loop: 4, n_cell_gap: 0
diff --git a/tests/mpi/tria_signals_05.mpirun=11.output b/tests/mpi/tria_signals_05.mpirun=11.output
new file mode 100644 (file)
index 0000000..69f2731
--- /dev/null
@@ -0,0 +1,109 @@
+
+DEAL:0:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:0:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:0:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:0:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:0:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:0:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:0:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:0:3d-3d::n_loop: 4, n_cell_gap: 0
+
+DEAL:1:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:1:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:1:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:1:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:1:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:1:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:1:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:1:3d-3d::n_loop: 4, n_cell_gap: 0
+
+
+DEAL:2:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:2:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:2:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:2:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:2:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:2:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:2:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:2:3d-3d::n_loop: 4, n_cell_gap: 0
+
+
+DEAL:3:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:3:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:3:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:3:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:3:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:3:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:3:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:3:3d-3d::n_loop: 4, n_cell_gap: 0
+
+
+DEAL:4:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:4:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:4:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:4:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:4:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:4:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:4:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:4:3d-3d::n_loop: 4, n_cell_gap: 0
+
+
+DEAL:5:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:5:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:5:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:5:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:5:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:5:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:5:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:5:3d-3d::n_loop: 4, n_cell_gap: 0
+
+
+DEAL:6:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:6:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:6:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:6:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:6:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:6:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:6:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:6:3d-3d::n_loop: 4, n_cell_gap: 0
+
+
+DEAL:7:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:7:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:7:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:7:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:7:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:7:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:7:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:7:3d-3d::n_loop: 4, n_cell_gap: 0
+
+
+DEAL:8:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:8:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:8:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:8:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:8:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:8:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:8:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:8:3d-3d::n_loop: 4, n_cell_gap: 0
+
+
+DEAL:9:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:9:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:9:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:9:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:9:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:9:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:9:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:9:3d-3d::n_loop: 4, n_cell_gap: 0
+
+
+DEAL:10:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:10:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:10:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:10:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:10:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:10:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:10:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:10:3d-3d::n_loop: 4, n_cell_gap: 0
+
diff --git a/tests/mpi/tria_signals_05.mpirun=4.output b/tests/mpi/tria_signals_05.mpirun=4.output
new file mode 100644 (file)
index 0000000..0df6a22
--- /dev/null
@@ -0,0 +1,39 @@
+
+DEAL:0:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:0:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:0:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:0:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:0:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:0:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:0:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:0:3d-3d::n_loop: 4, n_cell_gap: 0
+
+DEAL:1:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:1:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:1:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:1:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:1:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:1:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:1:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:1:3d-3d::n_loop: 4, n_cell_gap: 0
+
+
+DEAL:2:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:2:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:2:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:2:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:2:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:2:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:2:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:2:3d-3d::n_loop: 4, n_cell_gap: 0
+
+
+DEAL:3:2d-2d::n_loop: 1, n_cell_gap: 0
+DEAL:3:2d-2d::n_loop: 2, n_cell_gap: 0
+DEAL:3:2d-2d::n_loop: 3, n_cell_gap: 0
+DEAL:3:2d-2d::n_loop: 4, n_cell_gap: 0
+DEAL:3:3d-3d::n_loop: 1, n_cell_gap: 0
+DEAL:3:3d-3d::n_loop: 2, n_cell_gap: 0
+DEAL:3:3d-3d::n_loop: 3, n_cell_gap: 0
+DEAL:3:3d-3d::n_loop: 4, n_cell_gap: 0
+

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.