Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Update copyright lines with new year.
[simgrid.git] / src / smpi / colls / reduce / reduce-binomial.cpp
index ac0b789..065c434 100644 (file)
@@ -1,16 +1,18 @@
-/* Copyright (c) 2013-2014. The SimGrid Team.
+/* 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 <algorithm>
 
 //#include <star-reduction.c>
-
-int Coll_reduce_binomial::reduce(void *sendbuf, void *recvbuf, int count,
-                                    MPI_Datatype datatype, MPI_Op op, int root,
-                                    MPI_Comm comm)
+namespace simgrid{
+namespace smpi{
+int reduce__binomial(const void *sendbuf, void *recvbuf, int count,
+                     MPI_Datatype datatype, MPI_Op op, int root,
+                     MPI_Comm comm)
 {
   MPI_Status status;
   int comm_size, rank;
@@ -18,7 +20,6 @@ int Coll_reduce_binomial::reduce(void *sendbuf, void *recvbuf, int count,
   int dst;
   int tag = COLL_TAG_REDUCE;
   MPI_Aint extent;
-  void *tmp_buf;
   MPI_Aint true_lb, true_extent;
   if (count == 0)
     return 0;
@@ -27,12 +28,12 @@ int Coll_reduce_binomial::reduce(void *sendbuf, void *recvbuf, int count,
 
   extent = datatype->get_extent();
 
-  tmp_buf = (void *) smpi_get_tmp_sendbuffer(count * extent);
-  int is_commutative =  (op==MPI_OP_NULL || op->is_commutative());
+  unsigned char* tmp_buf = smpi_get_tmp_sendbuffer(count * extent);
+  bool is_commutative    = (op == MPI_OP_NULL || op->is_commutative());
   mask = 1;
-  
+
   int lroot;
-  if (is_commutative) 
+  if (is_commutative)
         lroot   = root;
   else
         lroot   = 0;
@@ -41,12 +42,12 @@ int Coll_reduce_binomial::reduce(void *sendbuf, void *recvbuf, int count,
   datatype->extent(&true_lb, &true_extent);
 
   /* adjust for potential negative lower bound in datatype */
-  tmp_buf = (void *)((char*)tmp_buf - true_lb);
-    
+  tmp_buf = tmp_buf - true_lb;
+
   /* If I'm not the root, then my recvbuf may not be valid, therefore
      I have to allocate a temporary one */
   if (rank != root) {
-      recvbuf = (void *) smpi_get_tmp_recvbuffer(count*(MAX(extent,true_extent)));
+      recvbuf = (void*)smpi_get_tmp_recvbuffer(count * std::max(extent, true_extent));
       recvbuf = (void *)((char*)recvbuf - true_lb);
   }
    if ((rank != root) || (sendbuf != MPI_IN_PLACE)) {
@@ -60,7 +61,7 @@ int Coll_reduce_binomial::reduce(void *sendbuf, void *recvbuf, int count,
       if (source < comm_size) {
         source = (source + lroot) % comm_size;
         Request::recv(tmp_buf, count, datatype, source, tag, comm, &status);
-        
+
         if (is_commutative) {
           if(op!=MPI_OP_NULL) op->apply( tmp_buf, recvbuf, &count, datatype);
         } else {
@@ -76,7 +77,7 @@ int Coll_reduce_binomial::reduce(void *sendbuf, void *recvbuf, int count,
     mask <<= 1;
   }
 
-  if (!is_commutative && (root != 0)){
+  if (not is_commutative && (root != 0)) {
     if (rank == 0){
       Request::send(recvbuf, count, datatype, root,tag, comm);
     }else if (rank == root){
@@ -85,9 +86,11 @@ int Coll_reduce_binomial::reduce(void *sendbuf, void *recvbuf, int count,
   }
 
   if (rank != root) {
-         smpi_free_tmp_buffer(recvbuf);
+    smpi_free_tmp_buffer(static_cast<unsigned char*>(recvbuf));
   }
   smpi_free_tmp_buffer(tmp_buf);
 
   return 0;
 }
+}
+}