]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
improve error checks and comments for parallel::distributed::Triangulation::load...
authorheister <heister@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 24 Jan 2014 15:47:40 +0000 (15:47 +0000)
committerheister <heister@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 24 Jan 2014 15:47:40 +0000 (15:47 +0000)
git-svn-id: https://svn.dealii.org/trunk@32298 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/distributed/tria.h
deal.II/source/distributed/tria.cc

index ce0a862fb224432925eedc2ef5ee205b98c7e90e..b00a01b224483e5e2cbe0825ded255754a79cc50 100644 (file)
@@ -514,17 +514,20 @@ namespace parallel
        * Save the refinement information from the coarse mesh into the given
        * file. This file needs to be reachable from all nodes in the computation
        * on a shared network file system. See the SolutionTransfer class
-       * on how to store solution vectors into this file.
+       * on how to store solution vectors into this file. Additional cell-based data can be saved
+       * using register_data_attach().
        */
       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(). You do not
+       * 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.
+       * 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);
 
index 4149b6e9a1e734b1534880703bdeeadb4b0eebd5..45db88fa607ed019885a64655e1db5048fe2f8f0 100644 (file)
@@ -2046,13 +2046,19 @@ namespace parallel
       if (attached_data_size>0)
         real_data_size = attached_data_size+sizeof(CellStatus);
 
+      Assert(this->n_cells()>0, ExcMessage("Can not save() an empty Triangulation."));
+
       if (my_subdomain==0)
         {
           std::string fname=std::string(filename)+".info";
           std::ofstream f(fname.c_str());
-          f << Utilities::MPI::n_mpi_processes (mpi_communicator) << " "
+          f << "version nproc attached_bytes n_attached_objs n_coarse_cells" << std::endl
+            << 2 << " "
+            << Utilities::MPI::n_mpi_processes (mpi_communicator) << " "
             << real_data_size << " "
-            << attached_data_pack_callbacks.size() << std::endl;
+            << attached_data_pack_callbacks.size() << " "
+            << this->n_cells(0)
+            << std::endl;
         }
 
       if (attached_data_size>0)
@@ -2082,6 +2088,10 @@ namespace parallel
     Triangulation<dim,spacedim>::
     load(const char *filename)
     {
+      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()."));
+
+
       if (parallel_ghost != 0)
         {
           dealii::internal::p4est::functions<dim>::ghost_destroy (parallel_ghost);
@@ -2090,15 +2100,20 @@ namespace parallel
       dealii::internal::p4est::functions<dim>::destroy (parallel_forest);
       parallel_forest = 0;
       dealii::internal::p4est::functions<dim>::connectivity_destroy (connectivity);
-      connectivity=0;
+      connectivity = 0;
 
-      unsigned int numcpus, attached_size, attached_count;
+      unsigned int version, numcpus, attached_size, attached_count, n_coarse_cells;
       {
         std::string fname=std::string(filename)+".info";
         std::ifstream f(fname.c_str());
-        f >> numcpus >> attached_size >> attached_count;
+        std::string firstline;
+        getline(f, firstline); //skip first line
+        f >> version >> numcpus >> attached_size >> attached_count >> n_coarse_cells;
       }
 
+      Assert(version == 2, ExcMessage("Incompatible version found in .info file."));
+      Assert(this->n_cells(0) == n_coarse_cells, ExcMessage("Number of coarse cells differ!"));
+
       attached_data_size = 0;
       n_attached_datas = 0;
       n_attached_deserialize = attached_count;

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.