Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
add two more functions
[simgrid.git] / src / smpi / smpi_pmpi.c
index f073ed7..f3269fd 100644 (file)
@@ -27,6 +27,7 @@ void TRACE_smpi_set_category(const char *category)
 int PMPI_Init(int *argc, char ***argv)
 {
   smpi_process_init(argc, argv);
+  smpi_process_mark_as_initialized();
 #ifdef HAVE_TRACING
   int rank = smpi_process_index();
   TRACE_smpi_init(rank);
@@ -303,7 +304,6 @@ int PMPI_Group_free(MPI_Group * group)
   if (group == NULL) {
     retval = MPI_ERR_ARG;
   } else {
-    if(*group!= smpi_comm_group(MPI_COMM_WORLD))// do not free the group of the comm_world
     smpi_group_destroy(*group);
     *group = MPI_GROUP_NULL;
     retval = MPI_SUCCESS;
@@ -419,7 +419,6 @@ int PMPI_Group_union(MPI_Group group1, MPI_Group group2,
         smpi_group_set_mapping(*newgroup, proc2, i);
       }
     }
-    smpi_group_use(*newgroup);
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin();
@@ -459,7 +458,6 @@ int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
         }
       }
     }
-    smpi_group_use(*newgroup);
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin();
@@ -496,7 +494,6 @@ int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgro
         }
       }
     }
-    smpi_group_use(*newgroup);
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin();
@@ -517,6 +514,11 @@ int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
       *newgroup = MPI_GROUP_EMPTY;
     } else if (n == smpi_group_size(group)) {
       *newgroup = group;
+      if(group!= smpi_comm_group(MPI_COMM_WORLD)
+                && group != MPI_GROUP_NULL
+                && group != smpi_comm_group(MPI_COMM_SELF)
+                && group != MPI_GROUP_EMPTY)
+      smpi_group_use(group);
     } else {
       *newgroup = smpi_group_new(n);
       for (i = 0; i < n; i++) {
@@ -524,7 +526,6 @@ int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
         smpi_group_set_mapping(*newgroup, index, i);
       }
     }
-    smpi_group_use(*newgroup);
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin();
@@ -543,6 +544,11 @@ int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
   } else {
     if (n == 0) {
       *newgroup = group;
+      if(group!= smpi_comm_group(MPI_COMM_WORLD)
+                && group != MPI_GROUP_NULL
+                && group != smpi_comm_group(MPI_COMM_SELF)
+                && group != MPI_GROUP_EMPTY)
+      smpi_group_use(group);
     } else if (n == smpi_group_size(group)) {
       *newgroup = MPI_GROUP_EMPTY;
     } else {
@@ -567,7 +573,6 @@ int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
 
       xbt_free(to_exclude);
     }
-    smpi_group_use(*newgroup);
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin();
@@ -626,7 +631,6 @@ int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3],
         }
       }
     }
-    smpi_group_use(*newgroup);
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin();
@@ -646,6 +650,11 @@ int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
   } else {
     if (n == 0) {
       *newgroup = group;
+      if(group!= smpi_comm_group(MPI_COMM_WORLD)
+                && group != MPI_GROUP_NULL
+                && group != smpi_comm_group(MPI_COMM_SELF)
+                && group != MPI_GROUP_EMPTY)
+      smpi_group_use(group);
     } else {
       size = smpi_group_size(group);
       for (i = 0; i < n; i++) {
@@ -699,7 +708,6 @@ int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
         }
       }
     }
-    smpi_group_use(*newgroup);
 
     retval = MPI_SUCCESS;
   }
@@ -769,6 +777,11 @@ int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
     retval = MPI_ERR_ARG;
   } else {
     *group = smpi_comm_group(comm);
+    if(*group!= smpi_comm_group(MPI_COMM_WORLD)
+              && *group != MPI_GROUP_NULL
+              && *group != smpi_comm_group(MPI_COMM_SELF)
+              && *group != MPI_GROUP_EMPTY)
+    smpi_group_use(*group);
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin();
@@ -1853,9 +1866,9 @@ int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
              || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
     retval = MPI_ERR_TYPE;
   } else {
-
-    if(recvbuf==MPI_IN_PLACE){
-       recvcount=0;
+    if (recvbuf == MPI_IN_PLACE) {
+        recvtype=sendtype;
+        recvcount=sendcount;
     }
     mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
                      recvtype, root, comm);
@@ -1890,11 +1903,10 @@ int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
              || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
     retval = MPI_ERR_TYPE;
   } else {
-
-    if(recvbuf==MPI_IN_PLACE){
-       recvcount=0;
+    if (recvbuf == MPI_IN_PLACE) {
+        recvtype=sendtype;
+        recvcount=sendcounts[smpi_comm_rank(comm)];
     }
-
     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
                       recvcount, recvtype, root, comm);
     retval = MPI_SUCCESS;
@@ -1925,17 +1937,7 @@ int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
     retval = MPI_ERR_ARG;
   } else {
 
-    char* sendtmpbuf = (char*) sendbuf;
-    if( sendbuf == MPI_IN_PLACE ) {
-      sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
-      smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
-    }
-
-    mpi_coll_reduce_fun(sendtmpbuf, recvbuf, count, datatype, op, root, comm);
-
-    if( sendbuf == MPI_IN_PLACE ) {
-      xbt_free(sendtmpbuf);
-    }
+    mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
 
     retval = MPI_SUCCESS;
   }
@@ -2315,6 +2317,21 @@ int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_
   return retval;
 }
 
+int PMPI_Type_create_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
+  int retval;
+
+  smpi_bench_end();
+  if (old_type == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else if (count<0){
+    retval = MPI_ERR_COUNT;
+  } else {
+    retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
   int retval,i;
 
@@ -2349,6 +2366,10 @@ int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatyp
   return retval;
 }
 
+int PMPI_Type_create_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
+  return PMPI_Type_hindexed(count, blocklens,indices,old_type,new_type);
+}
+
 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
   int retval,i;
 
@@ -2394,7 +2415,7 @@ int PMPI_Error_class(int errorcode, int* errorclass) {
 
 
 int PMPI_Initialized(int* flag) {
-   *flag=(smpi_process_data()!=NULL);
+   *flag=smpi_process_initialized();
    return MPI_SUCCESS;
 }