Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Update copyright lines with new year.
[simgrid.git] / src / smpi / colls / bcast / bcast-scatter-LR-allgather.cpp
index c9d2b20..7148af8 100644 (file)
@@ -1,10 +1,10 @@
-/* Copyright (c) 2013-2014. The SimGrid Team.
- * All rights reserved.                                                     */
+/* Copyright (c) 2013-2020. The SimGrid Team. All rights reserved.          */
 
 /* This program is free software; you can redistribute it and/or modify it
  * under the terms of the license (GNU LGPL) which comes with this package. */
 
-#include "../colls_private.h"
+#include "../colls_private.hpp"
+#include "smpi_status.hpp"
 
 /*****************************************************************************
 
@@ -64,19 +64,20 @@ SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
 
  * Descrp: broadcasts using a scatter followed by LR allgather.
 
- * Auther: MPIH / modified by Ahmad Faraj
+ * Author: MPIH / modified by Ahmad Faraj
 
  ****************************************************************************/
-int
-smpi_coll_tuned_bcast_scatter_LR_allgather(void *buff, int count,
-                                           MPI_Datatype data_type, int root,
-                                           MPI_Comm comm)
+namespace simgrid {
+namespace smpi {
+int bcast__scatter_LR_allgather(void *buff, int count,
+                                MPI_Datatype data_type, int root,
+                                MPI_Comm comm)
 {
   MPI_Aint extent;
   MPI_Status status;
   int i, src, dst, rank, num_procs;
   int mask, relative_rank, curr_size, recv_size, send_size, nbytes;
-  int scatter_size, left, right, next_src, *recv_counts, *disps;
+  int scatter_size, left, right, next_src;
   int tag = COLL_TAG_BCAST;
 
   rank = comm->rank();
@@ -85,7 +86,7 @@ smpi_coll_tuned_bcast_scatter_LR_allgather(void *buff, int count,
 
 
   nbytes = extent * count;
-  scatter_size = (nbytes + num_procs - 1) / num_procs;  // ceiling division 
+  scatter_size = (nbytes + num_procs - 1) / num_procs;  // ceiling division
   curr_size = (rank == root) ? nbytes : 0;      // root starts with all the data
   relative_rank = (rank >= root) ? rank - root : rank - root + num_procs;
 
@@ -101,11 +102,11 @@ smpi_coll_tuned_bcast_scatter_LR_allgather(void *buff, int count,
       //  allows you to post a larger recv.
       if (recv_size <= 0)
         curr_size = 0;          // this process doesn't receive any data
-      // because of uneven division 
+      // because of uneven division
       else {
         Request::recv((char *) buff + relative_rank * scatter_size, recv_size,
                  MPI_BYTE, src, tag, comm, &status);
-        curr_size = smpi_mpi_get_count(&status, MPI_BYTE);
+        curr_size = Status::get_count(&status, MPI_BYTE);
       }
       break;
     }
@@ -121,7 +122,7 @@ smpi_coll_tuned_bcast_scatter_LR_allgather(void *buff, int count,
   while (mask > 0) {
     if (relative_rank + mask < num_procs) {
       send_size = curr_size - scatter_size * mask;
-      // mask is also the size of this process's subtree 
+      // mask is also the size of this process's subtree
 
       if (send_size > 0) {
         dst = rank + mask;
@@ -137,8 +138,8 @@ smpi_coll_tuned_bcast_scatter_LR_allgather(void *buff, int count,
   }
 
   // done scatter now do allgather
-  recv_counts = (int *) xbt_malloc(sizeof(int) * num_procs);
-  disps = (int *) xbt_malloc(sizeof(int) * num_procs);
+  int* recv_counts = new int[num_procs];
+  int* disps       = new int[num_procs];
 
   for (i = 0; i < num_procs; i++) {
     recv_counts[i] = nbytes - i * scatter_size;
@@ -170,9 +171,11 @@ smpi_coll_tuned_bcast_scatter_LR_allgather(void *buff, int count,
     next_src = (num_procs + next_src - 1) % num_procs;
   }
 
-
-  free(recv_counts);
-  free(disps);
+  delete[] recv_counts;
+  delete[] disps;
 
   return MPI_SUCCESS;
 }
+
+}
+}