]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add test distributed_grids/2d_refinement_10
authorLei Qiao <qiaol618@gmail.com>
Thu, 24 Sep 2015 23:19:19 +0000 (18:19 -0500)
committerLei Qiao <qiaol618@gmail.com>
Fri, 25 Sep 2015 16:29:14 +0000 (11:29 -0500)
tests/distributed_grids/2d_refinement_10.cc [new file with mode: 0644]
tests/distributed_grids/2d_refinement_10.output [new file with mode: 0644]

diff --git a/tests/distributed_grids/2d_refinement_10.cc b/tests/distributed_grids/2d_refinement_10.cc
new file mode 100644 (file)
index 0000000..1f96d4d
--- /dev/null
@@ -0,0 +1,335 @@
+// ---------------------------------------------------------------------
+//
+// 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/grid_generator.h>
+#include <deal.II/distributed/grid_refinement.h>
+#include <deal.II/numerics/data_out.h>
+#include <deal.II/distributed/tria.h>
+#include <deal.II/lac/vector.h>
+#include <fstream>
+
+using namespace dealii;
+
+// Purpose of this test:
+// Demonstrate that mesh smoothing option limit_level_difference_at_vertices
+// is always effective across different initial cells (different trees in p4est context)
+// not matter it is set or not.
+//
+// Design of this test:
+// The initial mesh is set to 4 cells in unite square. Then the center most cell
+// in each level is refined consecutively for 4 times. Finally, total active
+// cell number in each refine level and cell center coordinates of the final
+// active cells are compared for cases with and without flag mesh smoothing
+// limit_level_difference_at_vertices.
+//
+// Expected result:
+// In principle, the result of the two refinement runs should be different.
+// However, under the current implementation, their results are the same.
+// So, if some day this test failed, it might be a good news.
+
+// Now the output of mesh and refine indicators in vtu format is disabled to
+// prevent unnecessary outputs during regression test. It you want to see the
+// mesh, uncomment the following macro definition.
+// #define __WRITE_VTU__
+
+template<int dim> class Location;
+
+template<int dim>
+class TriaTest
+{
+private:
+  typedef typename parallel::distributed::Triangulation<dim> TypeTria;
+public:
+  TriaTest(const typename dealii::Triangulation<dim>::MeshSmoothing smoothing_option = dealii::Triangulation<dim>::none);
+  void run(std::vector<unsigned int> &n_cell,
+           std::set<Location<dim> > &position_list);
+private:
+  void write_vtu(const unsigned int counter) const;
+
+  const MPI_Comm mpi_communicator;
+  TypeTria triangulation;
+  const unsigned int myid;
+  const bool I_am_host;
+  const std::string case_name;
+};
+
+// Equip class Point<dim> with operator < in order to sort the cell center coordinates.
+template<int dim>
+class Location: public Point<dim>
+{
+public:
+  Location(const Point<dim> &in)
+    :
+    Point<dim>(in)
+  {}
+
+  inline
+  bool operator<(const Location<dim> &op) const
+  {
+    for (unsigned int d=0; d<dim; ++d)
+      {
+        if ((*this)[d] == op[d])
+          {
+            continue;
+          }
+        return ((*this)[d] < op[d]);
+      }
+    return (false);
+  }
+};
+
+int main (int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, /* int max_num_threads */ 1);
+  const bool I_am_host = (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0);
+  // Although most part of this test is designed to run in parallel, there is still
+  // one place that doesn't work perfectly in parallel. Now the sorting of final cell
+  // centers is just a local operation and won't produce any reasonable result in parallel.
+  AssertThrow(I_am_host, ExcMessage("Current code works properly only with one process."));
+
+  const unsigned int dim = 2;
+  std::ofstream logfile;
+  if (I_am_host)
+    {
+      logfile.open("output");
+      deallog.attach(logfile);
+      deallog.depth_console(0);
+      deallog.threshold_double(1.e-10);
+      deallog.push("2d");
+    }
+
+  std::vector<unsigned int> n_cell_smooth;
+  std::vector<unsigned int> n_cell_no_smooth;
+
+  std::set<Location<dim> > final_cell_center_loactions_smooth;
+  std::set<Location<dim> > final_cell_center_loactions_no_smooth;
+  if (I_am_host)
+    {
+      deallog << "Flag limit_level_difference_at_vertices set:" << std::endl;
+    }
+  {
+    TriaTest<dim> tria_test(dealii::Triangulation<dim>::limit_level_difference_at_vertices);
+    tria_test.run(n_cell_smooth, final_cell_center_loactions_smooth);
+  }
+  if (I_am_host)
+    {
+      deallog << "Flag limit_level_difference_at_vertices unset:" << std::endl;
+    }
+  {
+    TriaTest<dim> tria_test(dealii::Triangulation<dim>::none);
+    tria_test.run(n_cell_no_smooth, final_cell_center_loactions_no_smooth);
+  }
+  if (I_am_host)
+    {
+      deallog << "Compare result:" << std::endl;
+    }
+  {
+    bool n_cells_are_same = (n_cell_smooth.size() == n_cell_no_smooth.size());
+
+    for (unsigned int i = 0; n_cells_are_same && (i<n_cell_smooth.size()); ++i)
+      {
+        n_cells_are_same = n_cells_are_same &&
+                           (n_cell_smooth[i] == n_cell_no_smooth[i]);
+      }
+    if (I_am_host)
+      {
+        deallog << "n_cells_are_same = " << n_cells_are_same << std::endl;
+      }
+  }
+  {
+    bool cell_center_loactions_are_same =
+      (final_cell_center_loactions_smooth.size() == final_cell_center_loactions_no_smooth.size());
+
+    std::set<Location<dim> >::const_iterator it1=final_cell_center_loactions_smooth.begin();
+    std::set<Location<dim> >::const_iterator it2=final_cell_center_loactions_no_smooth.begin();
+
+    for (; cell_center_loactions_are_same && (it1!=final_cell_center_loactions_smooth.end()); ++it1, ++it2)
+      {
+        cell_center_loactions_are_same = cell_center_loactions_are_same &&
+                                         (*it1 == *it2);
+      }
+    if (I_am_host)
+      {
+        deallog << "cell_center_loactions_are_same = " << cell_center_loactions_are_same << std::endl;
+      }
+  }
+
+  deallog.pop();
+
+  return (0);
+}
+
+template<int dim>
+TriaTest<dim>::TriaTest(const typename dealii::Triangulation<dim>::MeshSmoothing smoothing_option)
+  :
+  mpi_communicator(MPI_COMM_WORLD),
+  triangulation(mpi_communicator, smoothing_option),
+  myid (Utilities::MPI::this_mpi_process (mpi_communicator)),
+  I_am_host(myid == 0),
+  case_name(smoothing_option?"smooth-":"no_smooth-")
+{
+  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 (triangulation,
+                                             repetitions,
+                                             p1,
+                                             p2);
+
+  // GridGenerator::hyper_cube(triangulation);
+  // triangulation.refine_global(1);
+}
+
+template<int dim>
+void TriaTest<dim>::run(std::vector<unsigned int> &n_cell,
+                        std::set<Location<dim> > &position_list)
+{
+  n_cell.clear();
+  position_list.clear();
+
+  unsigned int counter = 0;
+  if (I_am_host)
+    {
+      deallog << "n_loop  n_cell" << std::endl;
+      deallog << counter << "       "
+              << triangulation.n_global_active_cells() << std::endl;
+    }
+  ++counter;
+
+  for (; counter < 5; ++counter)
+    {
+      {
+        Point<dim> p;
+        for (unsigned int d=0; d<dim; ++d)
+          {
+            p[d] = 0.5 - std::pow(0.5, 1.0 + counter);
+          }
+        typename TypeTria::active_cell_iterator
+        cell = triangulation.begin_active();
+        for (; cell!= triangulation.end(); ++cell)
+          if (cell->is_locally_owned() && ((cell->center()).distance(p) < 1e-4))
+            {
+              cell->set_refine_flag();
+            }
+      }
+
+      triangulation.prepare_coarsening_and_refinement();
+      write_vtu(counter);
+      triangulation.execute_coarsening_and_refinement();
+
+      if (I_am_host)
+        {
+          deallog << counter << "       "
+                  << triangulation.n_global_active_cells() << std::endl;
+        }
+      n_cell.push_back(triangulation.n_global_active_cells());
+    }
+
+  // Output cell center coordinates in a sorted manner to prevent change in
+  // cell ordering causing failure of this test.
+  // This is just a local operation, you have to re-implement this part if you
+  // want to run this test in distributed parallel.
+  {
+    deallog << " position of cell centers:" << std::endl;
+
+    for (typename TypeTria::active_cell_iterator cell = triangulation.begin_active();
+         cell!= triangulation.end();
+         ++cell)
+      if (cell->is_locally_owned())
+        {
+          Location<dim> loc(cell->center());
+          position_list.insert(loc);
+        }
+
+    for (typename std::set<Location<dim> >::const_iterator it = position_list.begin();
+         it != position_list.end();
+         ++it)
+      {
+        deallog << *it << std::endl;
+      }
+  }
+
+  return;
+}
+
+template<int dim>
+void TriaTest<dim>::write_vtu(const unsigned int counter) const
+{
+#ifdef __WRITE_VTU__
+  // save refine flag
+  Vector<float> refine_mark;
+  {
+    refine_mark.reinit (triangulation.n_active_cells());
+
+    typename TypeTria::active_cell_iterator
+    cell = triangulation.begin_active();
+    const typename TypeTria::active_cell_iterator
+    endc = triangulation.end();
+    for (; cell != endc; ++cell)
+      if (cell->is_locally_owned())
+        {
+          refine_mark[cell->active_cell_index()] =
+            cell->refine_flag_set();
+        }
+  }
+
+  DataOut<dim> data_out;
+  data_out.attach_triangulation (triangulation);
+  {
+    const std::string data_name ("refine_flag");
+    data_out.add_data_vector (refine_mark,
+                              data_name,
+                              DataOut<dim>::type_cell_data);
+  }
+
+  data_out.build_patches();
+
+  const std::string output_tag = case_name +
+                                 Utilities::int_to_string (counter, 4);
+  const std::string slot_itag = ".slot-" + Utilities::int_to_string (myid, 4);
+
+  std::ofstream output ((output_tag + slot_itag + ".vtu").c_str());
+  data_out.write_vtu (output);
+
+  if (I_am_host)
+    {
+      std::vector<std::string> filenames;
+      for (unsigned int i=0;
+           i<Utilities::MPI::n_mpi_processes (mpi_communicator);
+           ++i)
+        {
+          filenames.push_back (output_tag +
+                               ".slot-" +
+                               Utilities::int_to_string (i, 4) +
+                               ".vtu");
+        }
+      std::ofstream master_output ((output_tag + ".pvtu").c_str());
+      data_out.write_pvtu_record (master_output, filenames);
+    }
+#else
+  (void)counter;
+#endif
+  return;
+}
+
diff --git a/tests/distributed_grids/2d_refinement_10.output b/tests/distributed_grids/2d_refinement_10.output
new file mode 100644 (file)
index 0000000..fa349ed
--- /dev/null
@@ -0,0 +1,106 @@
+
+DEAL:2d::Flag limit_level_difference_at_vertices set:
+DEAL:2d::n_loop  n_cell
+DEAL:2d::0       4
+DEAL:2d::1       7
+DEAL:2d::2       19
+DEAL:2d::3       31
+DEAL:2d::4       43
+DEAL:2d:: position of cell centers:
+DEAL:2d::0.125000 0.125000
+DEAL:2d::0.125000 0.375000
+DEAL:2d::0.125000 0.625000
+DEAL:2d::0.125000 0.875000
+DEAL:2d::0.312500 0.312500
+DEAL:2d::0.312500 0.437500
+DEAL:2d::0.312500 0.562500
+DEAL:2d::0.312500 0.687500
+DEAL:2d::0.375000 0.125000
+DEAL:2d::0.375000 0.875000
+DEAL:2d::0.406250 0.406250
+DEAL:2d::0.406250 0.468750
+DEAL:2d::0.406250 0.531250
+DEAL:2d::0.406250 0.593750
+DEAL:2d::0.437500 0.312500
+DEAL:2d::0.437500 0.687500
+DEAL:2d::0.453125 0.453125
+DEAL:2d::0.453125 0.484375
+DEAL:2d::0.468750 0.406250
+DEAL:2d::0.468750 0.531250
+DEAL:2d::0.468750 0.593750
+DEAL:2d::0.484375 0.453125
+DEAL:2d::0.484375 0.484375
+DEAL:2d::0.531250 0.406250
+DEAL:2d::0.531250 0.468750
+DEAL:2d::0.531250 0.531250
+DEAL:2d::0.531250 0.593750
+DEAL:2d::0.562500 0.312500
+DEAL:2d::0.562500 0.687500
+DEAL:2d::0.593750 0.406250
+DEAL:2d::0.593750 0.468750
+DEAL:2d::0.593750 0.531250
+DEAL:2d::0.593750 0.593750
+DEAL:2d::0.625000 0.125000
+DEAL:2d::0.625000 0.875000
+DEAL:2d::0.687500 0.312500
+DEAL:2d::0.687500 0.437500
+DEAL:2d::0.687500 0.562500
+DEAL:2d::0.687500 0.687500
+DEAL:2d::0.875000 0.125000
+DEAL:2d::0.875000 0.375000
+DEAL:2d::0.875000 0.625000
+DEAL:2d::0.875000 0.875000
+DEAL:2d::Flag limit_level_difference_at_vertices unset:
+DEAL:2d::n_loop  n_cell
+DEAL:2d::0       4
+DEAL:2d::1       7
+DEAL:2d::2       19
+DEAL:2d::3       31
+DEAL:2d::4       43
+DEAL:2d:: position of cell centers:
+DEAL:2d::0.125000 0.125000
+DEAL:2d::0.125000 0.375000
+DEAL:2d::0.125000 0.625000
+DEAL:2d::0.125000 0.875000
+DEAL:2d::0.312500 0.312500
+DEAL:2d::0.312500 0.437500
+DEAL:2d::0.312500 0.562500
+DEAL:2d::0.312500 0.687500
+DEAL:2d::0.375000 0.125000
+DEAL:2d::0.375000 0.875000
+DEAL:2d::0.406250 0.406250
+DEAL:2d::0.406250 0.468750
+DEAL:2d::0.406250 0.531250
+DEAL:2d::0.406250 0.593750
+DEAL:2d::0.437500 0.312500
+DEAL:2d::0.437500 0.687500
+DEAL:2d::0.453125 0.453125
+DEAL:2d::0.453125 0.484375
+DEAL:2d::0.468750 0.406250
+DEAL:2d::0.468750 0.531250
+DEAL:2d::0.468750 0.593750
+DEAL:2d::0.484375 0.453125
+DEAL:2d::0.484375 0.484375
+DEAL:2d::0.531250 0.406250
+DEAL:2d::0.531250 0.468750
+DEAL:2d::0.531250 0.531250
+DEAL:2d::0.531250 0.593750
+DEAL:2d::0.562500 0.312500
+DEAL:2d::0.562500 0.687500
+DEAL:2d::0.593750 0.406250
+DEAL:2d::0.593750 0.468750
+DEAL:2d::0.593750 0.531250
+DEAL:2d::0.593750 0.593750
+DEAL:2d::0.625000 0.125000
+DEAL:2d::0.625000 0.875000
+DEAL:2d::0.687500 0.312500
+DEAL:2d::0.687500 0.437500
+DEAL:2d::0.687500 0.562500
+DEAL:2d::0.687500 0.687500
+DEAL:2d::0.875000 0.125000
+DEAL:2d::0.875000 0.375000
+DEAL:2d::0.875000 0.625000
+DEAL:2d::0.875000 0.875000
+DEAL:2d::Compare result:
+DEAL:2d::n_cells_are_same = 1
+DEAL:2d::cell_center_loactions_are_same = 1

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.