Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
number of elements can be 0 in type constructors
[simgrid.git] / src / smpi / smpi_pmpi.c
index caf671e..0db0452 100644 (file)
@@ -271,6 +271,7 @@ 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;
@@ -906,26 +907,41 @@ int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
   int retval;
 
   smpi_bench_end();
-#ifdef HAVE_TRACING
-  int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
-  int src_traced = smpi_group_index(smpi_comm_group(comm), src);
-  TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
-#endif
+
   if (request == NULL) {
     retval = MPI_ERR_ARG;
   } else if (comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
-  }else if (src == MPI_PROC_NULL) {
+  } else if (src == MPI_PROC_NULL) {
     *request = MPI_REQUEST_NULL;
     retval = MPI_SUCCESS;
+  } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
+    retval = MPI_ERR_COMM;
+  } else if (count < 0) {
+    retval = MPI_ERR_COUNT;
+  } else if (buf==NULL && count > 0) {
+    retval = MPI_ERR_COUNT;
+  } else if (datatype == MPI_DATATYPE_NULL){
+    retval = MPI_ERR_TYPE;
+  } else if(tag<0 && tag !=  MPI_ANY_TAG){
+    retval = MPI_ERR_TAG;
   } else {
+
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
+  int src_traced = smpi_group_index(smpi_comm_group(comm), src);
+  TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
+#endif
+
     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
     retval = MPI_SUCCESS;
-  }
+
 #ifdef HAVE_TRACING
   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
   (*request)->recv = 1;
 #endif
+  }
+
   smpi_bench_begin();
   return retval;
 }
@@ -944,6 +960,16 @@ int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
   } else if (dst == MPI_PROC_NULL) {
     *request = MPI_REQUEST_NULL;
     retval = MPI_SUCCESS;
+  } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
+    retval = MPI_ERR_COMM;
+  } else if (count < 0) {
+    retval = MPI_ERR_COUNT;
+  } else if (buf==NULL && count > 0) {
+    retval = MPI_ERR_COUNT;
+  } else if (datatype == MPI_DATATYPE_NULL){
+    retval = MPI_ERR_TYPE;
+  } else if(tag<0 && tag !=  MPI_ANY_TAG){
+    retval = MPI_ERR_TAG;
   } else {
 
 #ifdef HAVE_TRACING
@@ -983,6 +1009,16 @@ int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
     smpi_empty_status(status);
     status->MPI_SOURCE = MPI_PROC_NULL;
     retval = MPI_SUCCESS;
+  } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
+    retval = MPI_ERR_COMM;
+  } else if (count < 0) {
+    retval = MPI_ERR_COUNT;
+  } else if (buf==NULL && count > 0) {
+    retval = MPI_ERR_COUNT;
+  } else if (datatype == MPI_DATATYPE_NULL){
+    retval = MPI_ERR_TYPE;
+  } else if(tag<0 && tag !=  MPI_ANY_TAG){
+    retval = MPI_ERR_TAG;
   } else {
 #ifdef HAVE_TRACING
   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
@@ -1014,6 +1050,23 @@ int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
   int retval;
 
   smpi_bench_end();
+
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (dst == MPI_PROC_NULL) {
+    retval = MPI_SUCCESS;
+  } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
+    retval = MPI_ERR_COMM;
+  } else if (count < 0) {
+    retval = MPI_ERR_COUNT;
+  } else if (buf==NULL && count > 0) {
+    retval = MPI_ERR_COUNT;
+  } else if (datatype == MPI_DATATYPE_NULL){
+    retval = MPI_ERR_TYPE;
+  } else if(tag<0 && tag !=  MPI_ANY_TAG){
+    retval = MPI_ERR_TAG;
+  } else {
+
 #ifdef HAVE_TRACING
   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
   TRACE_smpi_computing_out(rank);
@@ -1021,18 +1074,16 @@ int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
   TRACE_smpi_send(rank, rank, dst_traced);
 #endif
-  if (comm == MPI_COMM_NULL) {
-    retval = MPI_ERR_COMM;
-  } else if (dst == MPI_PROC_NULL) {
-    retval = MPI_SUCCESS;
-  } else {
+
     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
     retval = MPI_SUCCESS;
-  }
+
 #ifdef HAVE_TRACING
   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
   TRACE_smpi_computing_in(rank);
 #endif
+  }
+
   smpi_bench_begin();
   return retval;
 }
@@ -1045,15 +1096,7 @@ int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
   int retval;
 
   smpi_bench_end();
-#ifdef HAVE_TRACING
-  int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
-  TRACE_smpi_computing_out(rank);
-  int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
-  int src_traced = smpi_group_index(smpi_comm_group(comm), src);
-  TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
-  TRACE_smpi_send(rank, rank, dst_traced);
-  TRACE_smpi_send(rank, src_traced, rank);
-#endif
+
   if (comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
   } else if (sendtype == MPI_DATATYPE_NULL
@@ -1063,18 +1106,41 @@ int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
       smpi_empty_status(status);
       status->MPI_SOURCE = MPI_PROC_NULL;
       retval = MPI_SUCCESS;
+  }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
+      (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
+    retval = MPI_ERR_COMM;
+  } else if (sendcount < 0 || recvcount<0) {
+      retval = MPI_ERR_COUNT;
+  } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
+    retval = MPI_ERR_COUNT;
+  } else if((sendtag<0 && sendtag !=  MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
+    retval = MPI_ERR_TAG;
   } else {
+
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
+  TRACE_smpi_computing_out(rank);
+  int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
+  int src_traced = smpi_group_index(smpi_comm_group(comm), src);
+  TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
+  TRACE_smpi_send(rank, rank, dst_traced);
+  TRACE_smpi_send(rank, src_traced, rank);
+#endif
+
+
     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
                       recvcount, recvtype, src, recvtag, comm, status);
     retval = MPI_SUCCESS;
-  }
+
 #ifdef HAVE_TRACING
   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
   TRACE_smpi_recv(rank, rank, dst_traced);
   TRACE_smpi_recv(rank, src_traced, rank);
   TRACE_smpi_computing_in(rank);
-
 #endif
+
+  }
+
   smpi_bench_begin();
   return retval;
 }
@@ -1278,20 +1344,23 @@ int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * sta
     retval = MPI_SUCCESS;
   }
 #ifdef HAVE_TRACING
-  int src_traced, dst_traced, is_wait_for_receive;
-  xbt_dynar_get_cpy(srcs, *index, &src_traced);
-  xbt_dynar_get_cpy(dsts, *index, &dst_traced);
-  xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
-  if (is_wait_for_receive) {
-    TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
+  if(*index!=MPI_UNDEFINED){
+    int src_traced, dst_traced, is_wait_for_receive;
+    xbt_dynar_get_cpy(srcs, *index, &src_traced);
+    xbt_dynar_get_cpy(dsts, *index, &dst_traced);
+    xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
+    if (is_wait_for_receive) {
+      TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
+    }
+    TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
+    //clean-up of dynars
+    xbt_dynar_free(&srcs);
+    xbt_dynar_free(&dsts);
+    xbt_dynar_free(&recvs);
   }
-  TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
-  //clean-up of dynars
-  xbt_dynar_free(&srcs);
-  xbt_dynar_free(&dsts);
-  xbt_dynar_free(&recvs);
   TRACE_smpi_computing_in(rank_traced);
 
+
 #endif
   smpi_bench_begin();
   return retval;
@@ -1884,7 +1953,7 @@ int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_typ
   smpi_bench_end();
   if (old_type == MPI_DATATYPE_NULL) {
     retval = MPI_ERR_TYPE;
-  } else if (count<=0){
+  } else if (count<0){
     retval = MPI_ERR_COUNT;
   } else {
     retval = smpi_datatype_contiguous(count, old_type, new_type);
@@ -1914,7 +1983,7 @@ int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type,
   smpi_bench_end();
   if (old_type == MPI_DATATYPE_NULL) {
     retval = MPI_ERR_TYPE;
-  } else if (count<=0 || blocklen<=0){
+  } else if (count<0 || blocklen<0){
     retval = MPI_ERR_COUNT;
   } else {
     retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
@@ -1929,7 +1998,7 @@ int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old
   smpi_bench_end();
   if (old_type == MPI_DATATYPE_NULL) {
     retval = MPI_ERR_TYPE;
-  } else if (count<=0 || blocklen<=0){
+  } else if (count<0 || blocklen<0){
     retval = MPI_ERR_COUNT;
   } else {
     retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
@@ -1945,7 +2014,7 @@ int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_
   smpi_bench_end();
   if (old_type == MPI_DATATYPE_NULL) {
     retval = MPI_ERR_TYPE;
-  } else if (count<=0){
+  } else if (count<0){
     retval = MPI_ERR_COUNT;
   } else {
     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
@@ -1960,7 +2029,7 @@ int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatyp
   smpi_bench_end();
   if (old_type == MPI_DATATYPE_NULL) {
     retval = MPI_ERR_TYPE;
-  } else if (count<=0){
+  } else if (count<0){
     retval = MPI_ERR_COUNT;
   } else {
     retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
@@ -1974,7 +2043,7 @@ int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype*
   int retval;
 
   smpi_bench_end();
-  if (count<=0){
+  if (count<0){
     retval = MPI_ERR_COUNT;
   } else {
     retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
@@ -2101,6 +2170,11 @@ int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int
    return not_yet_implemented();
 }
 
+int PMPI_Pcontrol(const int level )
+{
+   return not_yet_implemented();
+}
+
 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
    return not_yet_implemented();
 }