]> https://gitweb.dealii.org/ - dealii.git/commitdiff
fix variable transfer as well
authorTimo Heister <timo.heister@gmail.com>
Fri, 15 Oct 2021 14:56:39 +0000 (10:56 -0400)
committerTimo Heister <timo.heister@gmail.com>
Sat, 23 Oct 2021 13:42:30 +0000 (09:42 -0400)
source/distributed/tria_base.cc

index 0cf076328d886a74b421194a329b6039de531c41..e34cd0557a4a3ff4bbfeb9045afbde4f3d5a13e9 100644 (file)
@@ -1393,6 +1393,8 @@ namespace parallel
 
     const int myrank = Utilities::MPI::this_mpi_process(mpi_communicator);
 
+    const unsigned int bytes_per_cell = sizes_fixed_cumulative.back();
+
     //
     // ---------- Fixed size data ----------
     //
@@ -1444,8 +1446,8 @@ namespace parallel
       // Make sure we do the following computation in 64bit integers to be able
       // to handle 4GB+ files:
       const MPI_Offset my_global_file_position =
-        size_header + static_cast<MPI_Offset>(global_first_cell) *
-                        sizes_fixed_cumulative.back();
+        size_header +
+        static_cast<MPI_Offset>(global_first_cell) * bytes_per_cell;
 
       ierr = MPI_File_write_at(fh,
                                my_global_file_position,
@@ -1490,41 +1492,42 @@ namespace parallel
 
         // Write sizes of each cell into file simultaneously.
         {
+          const MPI_Offset my_global_file_position =
+            static_cast<MPI_Offset>(global_first_cell) * sizeof(unsigned int);
+
           const int *data = src_sizes_variable.data();
-          ierr =
-            MPI_File_write_at(fh,
-                              global_first_cell *
-                                sizeof(unsigned int), // global position in file
-                              DEAL_II_MPI_CONST_CAST(data),
-                              src_sizes_variable.size(), // local buffer
-                              MPI_INT,
-                              MPI_STATUS_IGNORE);
+
+          ierr = MPI_File_write_at(fh,
+                                   my_global_file_position,
+                                   DEAL_II_MPI_CONST_CAST(data),
+                                   src_sizes_variable.size(), // local buffer
+                                   MPI_INT,
+                                   MPI_STATUS_IGNORE);
           AssertThrowMPI(ierr);
         }
 
-
-        const unsigned int offset_variable =
-          global_num_cells * sizeof(unsigned int);
-
-        // Gather size of data in bytes we want to store from this processor.
-        const unsigned int size_on_proc = src_data_variable.size();
-
-        // Compute prefix sum
-        unsigned int prefix_sum = 0;
+        // Gather size of data in bytes we want to store from this
+        // processor and compute the prefix sum. We do this in 64 bit
+        // to avoid overflow for files larger than 4GB:
+        const unsigned long size_on_proc = src_data_variable.size();
+        unsigned long       prefix_sum   = 0;
         ierr = MPI_Exscan(DEAL_II_MPI_CONST_CAST(&size_on_proc),
                           &prefix_sum,
                           1,
-                          MPI_UNSIGNED,
+                          MPI_UNSIGNED_LONG,
                           MPI_SUM,
                           mpi_communicator);
         AssertThrowMPI(ierr);
 
+        const MPI_Offset my_global_file_position =
+          static_cast<MPI_Offset>(global_num_cells) * sizeof(unsigned int) +
+          prefix_sum;
+
         const char *data = src_data_variable.data();
 
         // Write data consecutively into file.
         ierr = MPI_File_write_at(fh,
-                                 offset_variable +
-                                   prefix_sum, // global position in file
+                                 my_global_file_position,
                                  DEAL_II_MPI_CONST_CAST(data),
                                  src_data_variable.size(), // local buffer
                                  MPI_CHAR,
@@ -1601,7 +1604,8 @@ namespace parallel
       AssertThrowMPI(ierr);
 
       // Allocate sufficient memory.
-      dest_data_fixed.resize(local_num_cells * sizes_fixed_cumulative.back());
+      const unsigned int bytes_per_cell = sizes_fixed_cumulative.back();
+      dest_data_fixed.resize(local_num_cells * bytes_per_cell);
 
       // Read packed data from file simultaneously.
       const MPI_Offset size_header =
@@ -1610,8 +1614,8 @@ namespace parallel
       // Make sure we do the following computation in 64bit integers to be able
       // to handle 4GB+ files:
       const MPI_Offset my_global_file_position =
-        size_header + static_cast<MPI_Offset>(global_first_cell) *
-                        sizes_fixed_cumulative.back();
+        size_header +
+        static_cast<MPI_Offset>(global_first_cell) * bytes_per_cell;
 
       ierr = MPI_File_read_at(fh,
                               my_global_file_position,
@@ -1650,34 +1654,42 @@ namespace parallel
 
         // Read sizes of all locally owned cells.
         dest_sizes_variable.resize(local_num_cells);
+
+        const MPI_Offset my_global_file_position_sizes =
+          static_cast<MPI_Offset>(global_first_cell) * sizeof(unsigned int);
+
         ierr = MPI_File_read_at(fh,
-                                global_first_cell * sizeof(unsigned int),
+                                my_global_file_position_sizes,
                                 dest_sizes_variable.data(),
                                 dest_sizes_variable.size(),
                                 MPI_INT,
                                 MPI_STATUS_IGNORE);
         AssertThrowMPI(ierr);
 
-        const unsigned int offset = global_num_cells * sizeof(unsigned int);
 
-        const unsigned int size_on_proc =
+        // Compute my data size in bytes and compute prefix sum. We do this in
+        // 64 bit to avoid overflow for files larger than 4 GB:
+        const unsigned long size_on_proc =
           std::accumulate(dest_sizes_variable.begin(),
                           dest_sizes_variable.end(),
                           0);
 
-        // share information among all processors by prefix sum
-        unsigned int prefix_sum = 0;
+        unsigned long prefix_sum = 0;
         ierr = MPI_Exscan(DEAL_II_MPI_CONST_CAST(&size_on_proc),
                           &prefix_sum,
                           1,
-                          MPI_UNSIGNED,
+                          MPI_UNSIGNED_LONG,
                           MPI_SUM,
                           mpi_communicator);
         AssertThrowMPI(ierr);
 
+        const MPI_Offset my_global_file_position =
+          static_cast<MPI_Offset>(global_num_cells) * sizeof(unsigned int) +
+          prefix_sum;
+
         dest_data_variable.resize(size_on_proc);
         ierr = MPI_File_read_at(fh,
-                                offset + prefix_sum,
+                                my_global_file_position,
                                 dest_data_variable.data(),
                                 dest_data_variable.size(),
                                 MPI_CHAR,

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.