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 07f5101..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;
@@ -914,10 +915,16 @@ int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
   } 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
@@ -953,10 +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
@@ -996,10 +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;
@@ -1036,10 +1055,16 @@ int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
     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
@@ -1081,8 +1106,15 @@ 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
@@ -1921,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);
@@ -1951,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);
@@ -1966,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);
@@ -1982,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);
@@ -1997,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);
@@ -2011,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);