Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add support for various MPI_Type functions, to handle creation of new MPI types
authorAugustin Degomme <degomme@idpann.imag.fr>
Fri, 21 Sep 2012 16:07:08 +0000 (18:07 +0200)
committerAugustin Degomme <degomme@idpann.imag.fr>
Fri, 21 Sep 2012 16:07:44 +0000 (18:07 +0200)
include/smpi/smpi.h
src/smpi/private.h
src/smpi/smpi_mpi_dt.c
src/smpi/smpi_pmpi.c

index 3071425..253b117 100644 (file)
@@ -173,7 +173,28 @@ MPI_CALL(XBT_PUBLIC(int), MPI_Type_get_extent,
 MPI_CALL(XBT_PUBLIC(int), MPI_Type_extent, (MPI_Datatype datatype, MPI_Aint * extent));
 MPI_CALL(XBT_PUBLIC(int), MPI_Type_lb, (MPI_Datatype datatype, MPI_Aint * disp));
 MPI_CALL(XBT_PUBLIC(int), MPI_Type_ub, (MPI_Datatype datatype, MPI_Aint * disp));
-
+MPI_CALL(XBT_PUBLIC(int), MPI_Type_commit, (MPI_Datatype* datatype));
+MPI_CALL(XBT_PUBLIC(int), MPI_Type_hindexed,
+                            (int count, int* blocklens, MPI_Aint* indices,
+                            MPI_Datatype old_type, MPI_Datatype* newtype));
+MPI_CALL(XBT_PUBLIC(int), MPI_Type_hvector,
+                            (int count, int blocklen, MPI_Aint stride,
+                             MPI_Datatype old_type, MPI_Datatype* newtype));
+MPI_CALL(XBT_PUBLIC(int), MPI_Type_indexed,
+                            (int count, int* blocklens, int* indices,
+                             MPI_Datatype old_type, MPI_Datatype* newtype));
+MPI_CALL(XBT_PUBLIC(int), MPI_Type_struct,
+                            (int count, int* blocklens, MPI_Aint* indices,
+                             MPI_Datatype* old_types, MPI_Datatype* newtype));
+MPI_CALL(XBT_PUBLIC(int), MPI_Type_vector,
+                            (int count, int blocklen, int stride,
+                             MPI_Datatype old_type, MPI_Datatype* newtype));
+MPI_CALL(XBT_PUBLIC(int), MPI_Type_contiguous,
+                            (int count, MPI_Datatype old_type,
+                             MPI_Datatype* newtype));
+MPI_CALL(XBT_PUBLIC(int), MPI_Testall,
+                            (int count, MPI_Request* requests, int* flag,
+                             MPI_Status* statuses));
 MPI_CALL(XBT_PUBLIC(int), MPI_Op_create,
                             (MPI_User_function * function, int commute,
                              MPI_Op * op));
@@ -352,6 +373,7 @@ MPI_CALL(XBT_PUBLIC(int), MPI_Probe,
                             (int source, int tag, MPI_Comm comm,
                              MPI_Status* status));
 
+
 //FIXME: these are not yet implemented
 typedef void MPI_Handler_function(MPI_Comm*, int*, ...);
 typedef void* MPI_Errhandler;
@@ -381,19 +403,12 @@ MPI_CALL(XBT_PUBLIC(int), MPI_Errhandler_free, (MPI_Errhandler* errhandler));
 MPI_CALL(XBT_PUBLIC(int), MPI_Errhandler_get, (MPI_Comm comm, MPI_Errhandler* errhandler));
 MPI_CALL(XBT_PUBLIC(int), MPI_Error_string, (int errorcode, char* string, int* resultlen));
 MPI_CALL(XBT_PUBLIC(int), MPI_Errhandler_set, (MPI_Comm comm, MPI_Errhandler errhandler));
-MPI_CALL(XBT_PUBLIC(int), MPI_Type_contiguous, (int count, MPI_Datatype old_type, MPI_Datatype* newtype));
 MPI_CALL(XBT_PUBLIC(int), MPI_Cancel, (MPI_Request* request));
 MPI_CALL(XBT_PUBLIC(int), MPI_Buffer_attach, (void* buffer, int size));
 MPI_CALL(XBT_PUBLIC(int), MPI_Buffer_detach, (void* buffer, int* size));
 MPI_CALL(XBT_PUBLIC(int), MPI_Testsome, (int incount, MPI_Request* requests, int* outcount, int* indices, MPI_Status* statuses));
 MPI_CALL(XBT_PUBLIC(int), MPI_Comm_test_inter, (MPI_Comm comm, int* flag));
 MPI_CALL(XBT_PUBLIC(int), MPI_Unpack, (void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm));
-MPI_CALL(XBT_PUBLIC(int), MPI_Type_commit, (MPI_Datatype* datatype));
-MPI_CALL(XBT_PUBLIC(int), MPI_Type_hindexed, (int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* newtype));
-MPI_CALL(XBT_PUBLIC(int), MPI_Type_hvector, (int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* newtype));
-MPI_CALL(XBT_PUBLIC(int), MPI_Type_indexed, (int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* newtype));
-MPI_CALL(XBT_PUBLIC(int), MPI_Type_struct, (int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* newtype));
-MPI_CALL(XBT_PUBLIC(int), MPI_Type_vector, (int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* newtype));
 MPI_CALL(XBT_PUBLIC(int), MPI_Ssend, (void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm));
 MPI_CALL(XBT_PUBLIC(int), MPI_Ssend_init, (void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request));
 MPI_CALL(XBT_PUBLIC(int), MPI_Intercomm_create, (MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out));
@@ -414,7 +429,6 @@ MPI_CALL(XBT_PUBLIC(int), MPI_Keyval_create, (MPI_Copy_function* copy_fn, MPI_De
 MPI_CALL(XBT_PUBLIC(int), MPI_Keyval_free, (int* keyval));
 MPI_CALL(XBT_PUBLIC(int), MPI_Test_cancelled, (MPI_Status* status, int* flag));
 MPI_CALL(XBT_PUBLIC(int), MPI_Pack, (void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm));
-MPI_CALL(XBT_PUBLIC(int), MPI_Testall, (int count, MPI_Request* requests, int* flag, MPI_Status* statuses));
 MPI_CALL(XBT_PUBLIC(int), MPI_Get_elements, (MPI_Status* status, MPI_Datatype datatype, int* elements));
 MPI_CALL(XBT_PUBLIC(int), MPI_Dims_create, (int nnodes, int ndims, int* dims));
 MPI_CALL(XBT_PUBLIC(int), MPI_Initialized, (int* flag));
index b0a8f5f..c657011 100644 (file)
@@ -68,6 +68,24 @@ int smpi_datatype_extent(MPI_Datatype datatype, MPI_Aint * lb,
 int smpi_datatype_copy(void *sendbuf, int sendcount, MPI_Datatype sendtype,
                        void *recvbuf, int recvcount,
                        MPI_Datatype recvtype);
+int smpi_datatype_contiguous(int count, MPI_Datatype old_type,
+                       MPI_Datatype* new_type);
+int smpi_datatype_vector(int count, int blocklen, int stride,
+                      MPI_Datatype old_type, MPI_Datatype* new_type);
+
+int smpi_datatype_hvector(int count, int blocklen, MPI_Aint stride,
+                      MPI_Datatype old_type, MPI_Datatype* new_type);
+int smpi_datatype_indexed(int count, int* blocklens, int* indices,
+                     MPI_Datatype old_type, MPI_Datatype* new_type);
+int smpi_datatype_hindexed(int count, int* blocklens, MPI_Aint* indices,
+                     MPI_Datatype old_type, MPI_Datatype* new_type);
+int smpi_datatype_struct(int count, int* blocklens, MPI_Aint* indices,
+                    MPI_Datatype* old_types, MPI_Datatype* new_type);
+void smpi_datatype_create(MPI_Datatype* new_type, int size, int flags);
+void smpi_datatype_free(MPI_Datatype* type);
+void smpi_datatype_commit(MPI_Datatype* datatype);
+
+
 
 MPI_Op smpi_op_new(MPI_User_function * function, int commute);
 void smpi_op_destroy(MPI_Op op);
index 700870d..ca60a35 100644 (file)
@@ -151,6 +151,110 @@ int smpi_datatype_copy(void *sendbuf, int sendcount, MPI_Datatype sendtype,
   return retval;
 }
 
+void smpi_datatype_create(MPI_Datatype* new_type, int size, int flags){
+  MPI_Datatype new_t= xbt_new(s_smpi_mpi_datatype_t,1);
+  new_t->size=size;
+  new_t->lb=0;
+  new_t->ub=size;
+  new_t->flags=flags;
+  *new_type = new_t;
+}
+
+void smpi_datatype_free(MPI_Datatype* type){
+  xbt_free(*type);
+}
+
+int smpi_datatype_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type)
+{
+  int retval;
+  if ((old_type->flags & DT_FLAG_COMMITED) != DT_FLAG_COMMITED) {
+     retval = MPI_ERR_TYPE;
+  } else {
+    smpi_datatype_create(new_type, count * smpi_datatype_size(old_type), DT_FLAG_CONTIGUOUS);
+    retval=MPI_SUCCESS;
+  }
+  return retval;
+}
+
+int smpi_datatype_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type)
+{
+  int retval;
+  if (blocklen<=0)return MPI_ERR_ARG;
+  if ((old_type->flags & DT_FLAG_COMMITED) != DT_FLAG_COMMITED) {
+     retval = MPI_ERR_TYPE;
+  } else {
+    smpi_datatype_create(new_type, count * (blocklen+stride) * smpi_datatype_size(old_type), DT_FLAG_VECTOR);
+    retval=MPI_SUCCESS;
+  }
+  return retval;
+}
+
+int smpi_datatype_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type)
+{
+  int retval;
+  if (blocklen<=0)return MPI_ERR_ARG;
+  if ((old_type->flags & DT_FLAG_COMMITED) != DT_FLAG_COMMITED) {
+     retval = MPI_ERR_TYPE;
+  } else {
+    smpi_datatype_create(new_type, count * ((blocklen * smpi_datatype_size(old_type))+stride), DT_FLAG_VECTOR);
+    retval=MPI_SUCCESS;
+  }
+  return retval;
+}
+
+
+int smpi_datatype_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type)
+{
+  int i;
+  int retval;
+  for(i=0; i< count; i++){
+    if   (blocklens[i]<=0)
+      return MPI_ERR_ARG;
+  }
+  if ((old_type->flags & DT_FLAG_COMMITED) != DT_FLAG_COMMITED) {
+     retval = MPI_ERR_TYPE;
+  } else {
+    smpi_datatype_create(new_type,  (blocklens[count-1] + indices[count-1]) * smpi_datatype_size(old_type), DT_FLAG_DATA);
+    retval=MPI_SUCCESS;
+  }
+  return retval;
+}
+
+int smpi_datatype_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type)
+{
+  int i;
+  int retval;
+  for(i=0; i< count; i++){
+    if   (blocklens[i]<=0)
+      return MPI_ERR_ARG;
+  }
+  if ((old_type->flags & DT_FLAG_COMMITED) != DT_FLAG_COMMITED) {
+     retval = MPI_ERR_TYPE;
+  } else {
+    smpi_datatype_create(new_type,indices[count-1] + (blocklens[count-1]  * smpi_datatype_size(old_type)), DT_FLAG_DATA);
+    retval=MPI_SUCCESS;
+  }
+  return retval;
+}
+
+int smpi_datatype_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type)
+{
+  int i;
+  for(i=0; i< count; i++){
+    if (blocklens[i]<=0)
+      return MPI_ERR_ARG;
+    if ((old_types[i]->flags & DT_FLAG_COMMITED) != DT_FLAG_COMMITED)
+      return MPI_ERR_TYPE;
+  }
+  smpi_datatype_create(new_type,indices[count-1] + (blocklens[count-1]  * smpi_datatype_size(old_types[count-1])), DT_FLAG_DATA);
+  return MPI_SUCCESS;
+}
+
+void smpi_datatype_commit(MPI_Datatype* datatype)
+{
+  (*datatype)->flags= ( (*datatype)->flags | DT_FLAG_COMMITED);
+}
+
 typedef struct s_smpi_mpi_op {
   MPI_User_function *func;
 } s_smpi_mpi_op_t;
index 5f7bd69..16fccc0 100644 (file)
@@ -138,8 +138,8 @@ int PMPI_Type_free(MPI_Datatype * datatype)
   if (!datatype) {
     retval = MPI_ERR_ARG;
   } else {
-    // FIXME: always fail for now
-    retval = MPI_ERR_TYPE;
+    smpi_datatype_free(datatype);
+    retval = MPI_SUCCESS;
   }
   smpi_bench_begin();
   return retval;
@@ -922,6 +922,7 @@ int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
   return retval;
 }
 
+
 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
               int tag, MPI_Comm comm, MPI_Request * request)
 {
@@ -952,6 +953,8 @@ int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
   return retval;
 }
 
+
+
 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
              MPI_Comm comm, MPI_Status * status)
 {
@@ -1806,12 +1809,117 @@ int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
   return retval;
 }
 
+int PMPI_Type_contiguous(int count, 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_contiguous(count, old_type, new_type);
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+int PMPI_Type_commit(MPI_Datatype* datatype) {
+  int retval;
+
+  smpi_bench_end();
+  if (datatype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else {
+    smpi_datatype_commit(datatype);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+
+int PMPI_Type_vector(int count, int blocklen, int stride, 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 || blocklen<=0){
+    retval = MPI_ERR_COUNT;
+  } else {
+    retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, 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 || blocklen<=0){
+    retval = MPI_ERR_COUNT;
+  } else {
+    retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+
+int PMPI_Type_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_hindexed(int count, int* blocklens, MPI_Aint* 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_hindexed(count, blocklens, indices, old_type, new_type);
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+
+int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
+  int retval;
+
+  smpi_bench_end();
+  if (count<=0){
+    retval = MPI_ERR_COUNT;
+  } else {
+    retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
+  }
+  smpi_bench_begin();
+  return retval;}
+
+
 /* The following calls are not yet implemented and will fail at runtime. */
 /* Once implemented, please move them above this notice. */
 
 static int not_yet_implemented(void) {
-   xbt_die("Not yet implemented");
-   return MPI_ERR_OTHER;
+         XBT_WARN("Not yet implemented");
+   return MPI_SUCCESS;
 }
 
 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
@@ -1902,9 +2010,6 @@ int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
    return not_yet_implemented();
 }
 
-int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* newtype) {
-   return not_yet_implemented();
-}
 
 int PMPI_Cancel(MPI_Request* request) {
    return not_yet_implemented();
@@ -1930,30 +2035,6 @@ int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcou
    return not_yet_implemented();
 }
 
-int PMPI_Type_commit(MPI_Datatype* datatype) {
-   return not_yet_implemented();
-}
-
-int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* newtype) {
-   return not_yet_implemented();
-}
-
-int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* newtype) {
-   return not_yet_implemented();
-}
-
-int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* newtype) {
-   return not_yet_implemented();
-}
-
-int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* newtype) {
-   return not_yet_implemented();
-}
-
-int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* newtype) {
-   return not_yet_implemented();
-}
-
 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
    return not_yet_implemented();
 }
@@ -1995,7 +2076,6 @@ int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dest, int tag,
 }
 
 
-
 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
    return not_yet_implemented();
 }