]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Add a boolean autopartition to
authormaier <maier@0785d39b-7218-0410-832d-ea1e28bc413d>
Sat, 26 Apr 2014 00:38:26 +0000 (00:38 +0000)
committermaier <maier@0785d39b-7218-0410-832d-ea1e28bc413d>
Sat, 26 Apr 2014 00:38:26 +0000 (00:38 +0000)
parallel::distributed::Triangulation<...>::load to control automatic
reparitioning behaviour of p4est

git-svn-id: https://svn.dealii.org/trunk@32839 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/include/deal.II/distributed/tria.h
deal.II/source/distributed/tria.cc
tests/mpi/p4est_save_01.cc
tests/mpi/p4est_save_02.cc
tests/mpi/p4est_save_03.cc

index 899dd3ce74e21d2a483e582ed4813f688e356b1d..4af2ee48b45cb9eeaeaa1cb65b79b6f6ed2c1744 100644 (file)
@@ -148,6 +148,15 @@ inconvenience this causes.
 <h3>Specific improvements</h3>
 
 <ol>
+  <li> Fixed: parallel::distributed::Triangulation<...>::load now has an
+  additional parameter autopartition to control p4est's behaviour of
+  rebalancing triangulations between mpi nodes upon reading. It is
+  particularly useful to disable this behaviour when data is stored
+  seperately (see tests mpi/p4est_save_0?).
+  <br>
+  (Alexander Grayver, Matthias Maier, 2014/04/26)
+  </li>
+
   <li> Fixed: GridTools::find_active_cell_around_point() could get into an infinite
   loop if the point we are looking for is in fact not within the domain. This is now
   fixed.
index 17dbdee03466ee08d283976bca6ece4a8cd20dae..673b0dc13202bbb51579bbd34921e895955650e5 100644 (file)
@@ -553,16 +553,24 @@ namespace parallel
       void save(const char *filename) const;
 
       /**
-       * Load the refinement information saved with save() back in. The mesh
-       * must contain the same coarse mesh that was used in save() before calling
-       * this function. You do not
-       * need to load with the same number of MPI processes that you saved
-       * with. Rather, if a mesh is loaded with a different number of MPI
-       * processes than used at the time of saving, the mesh is repartitioned
-       * appropriately. Cell-based data that was saved with register_data_attach()
-       * can be read in with notify_ready_to_unpack() after calling load().
-       */
-      void load(const char *filename);
+       * Load the refinement information saved with save() back in. The
+       * mesh must contain the same coarse mesh that was used in save()
+       * before calling this function.
+       *
+       * You do not need to load with the same number of MPI processes that
+       * you saved with. Rather, if a mesh is loaded with a different
+       * number of MPI processes than used at the time of saving, the mesh
+       * is repartitioned appropriately. Cell-based data that was saved
+       * with register_data_attach() can be read in with
+       * notify_ready_to_unpack() after calling load().
+       *
+       * If you use p4est version > 0.3.4.2 the @p autopartition flag tells
+       * p4est to ignore the partitioning that the triangulation had when
+       * it was saved and make it uniform upon loading. If @p autopartition
+       * set to false, the triangulation is only repartitioned if needed
+       * (i.e. if a different number of MPI processes is encountered).
+       */
+      void load(const char *filename, bool autopartition = true);
 
       /**
        * Used to inform in the callbacks of register_data_attach() and
index 6a22f877f596ce437f6191281bd9a6e42a3ef01e..48044be5a94d0304ab29cd7bfd37f94899050fc2 100644 (file)
@@ -2214,7 +2214,7 @@ namespace parallel
     template <int dim, int spacedim>
     void
     Triangulation<dim,spacedim>::
-    load(const char *filename)
+    load(const char *filename, bool autopartition)
     {
       Assert(this->n_cells()>0, ExcMessage("load() only works if Triangulation already contains the same coarse mesh!"));
       Assert(this->n_levels()==1, ExcMessage("Triangulation may only contain coarse cells when calling load()."));
@@ -2257,7 +2257,7 @@ namespace parallel
       parallel_forest = dealii::internal::p4est::functions<dim>::load_ext (
                           filename, mpi_communicator,
                           attached_size, attached_size>0,
-                          1, 0,
+                          autopartition, 0,
                           this,
                           &connectivity);
 #else
index aa26490503847670b31e9eca783fbab49ab6d9f0..6b03577a5b773e47b1141641e03dfc7aabec496c 100644 (file)
@@ -87,7 +87,7 @@ void test()
     parallel::distributed::Triangulation<dim> tr(MPI_COMM_WORLD);
 
     GridGenerator::hyper_cube(tr);
-    tr.load(filename.c_str());
+    tr.load(filename.c_str(), false);
 
     if (myid == 0)
       {
index 20f71ba0707d192307cde9e9dc0c4ec42125885a..177242bcf1b93e99ba947fd811753c68d79a39cb 100644 (file)
@@ -115,7 +115,7 @@ void test()
     parallel::distributed::Triangulation<dim> tr(MPI_COMM_WORLD);
 
     GridGenerator::hyper_cube(tr);
-    tr.load(filename.c_str());
+    tr.load(filename.c_str(), false);
     FE_Q<dim> fe(1);
     DoFHandler<dim> dh(tr);
 
index b427832fc8805b2bbbaed1faf2deeb5ef8a805aa..9a224f7d14ee9e48b5db3c84db62111298e698db 100644 (file)
@@ -124,7 +124,7 @@ void test()
     parallel::distributed::Triangulation<dim> tr (MPI_COMM_WORLD);
 
     GridGenerator::hyper_cube (tr);
-    tr.load (filename.c_str());
+    tr.load (filename.c_str(), false);
     FE_Q<dim> fe (1);
     DoFHandler<dim> dh (tr);
 

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.