Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
add better support for MPI datatypes extent values (to correct behavior of gather...
[simgrid.git] / src / smpi / smpi_mpi_dt.c
index d75b7fa..d277579 100644 (file)
@@ -28,6 +28,16 @@ XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_mpi_dt, smpi,
   };                                          \
 MPI_Datatype name = &mpi_##name;
 
+#define CREATE_MPI_DATATYPE_NULL(name)       \
+  static s_smpi_mpi_datatype_t mpi_##name = { \
+    0,  /* size */                 \
+    0,             /*was 1 has_subtype*/             \
+    0,             /* lb */                   \
+    0,  /* ub = lb + size */       \
+    DT_FLAG_BASIC,  /* flags */              \
+    NULL           /* pointer on extended struct*/ \
+  };                                          \
+MPI_Datatype name = &mpi_##name;
 
 //The following are datatypes for the MPI functions MPI_MAXLOC and MPI_MINLOC.
 typedef struct {
@@ -93,6 +103,8 @@ CREATE_MPI_DATATYPE(MPI_SHORT_INT, short_int);
 CREATE_MPI_DATATYPE(MPI_2INT, int_int);
 CREATE_MPI_DATATYPE(MPI_LONG_DOUBLE_INT, long_double_int);
 
+CREATE_MPI_DATATYPE_NULL(MPI_UB);
+CREATE_MPI_DATATYPE_NULL(MPI_LB);
 // Internal use only
 CREATE_MPI_DATATYPE(MPI_PTR, void*);
 
@@ -190,13 +202,18 @@ void serialize_vector( const void *noncontiguous_vector,
 {
   s_smpi_mpi_vector_t* type_c = (s_smpi_mpi_vector_t*)type;
   int i;
+  char* contiguous_vector_char = (char*)contiguous_vector;
+  char* noncontiguous_vector_char = (char*)noncontiguous_vector;
+
   for (i = 0; i < type_c->block_count * count; i++) {
-    memcpy(contiguous_vector,
-           noncontiguous_vector, type_c->block_length * type_c->size_oldtype);
-    contiguous_vector += type_c->block_length*type_c->size_oldtype;
-    noncontiguous_vector += type_c->block_stride*type_c->size_oldtype;
+    memcpy(contiguous_vector_char,
+           noncontiguous_vector_char, type_c->block_length * type_c->size_oldtype);
+
+    contiguous_vector_char += type_c->block_length*type_c->size_oldtype;
+    noncontiguous_vector_char += type_c->block_stride*type_c->size_oldtype;
   }
 }
+
 /*
  *  Copies contiguous data into noncontiguous memory.
  *  @param noncontiguous_vector - output vector
@@ -213,17 +230,22 @@ void unserialize_vector( const void *contiguous_vector,
 {
   s_smpi_mpi_vector_t* type_c = (s_smpi_mpi_vector_t*)type;
   int i;
+
+  char* contiguous_vector_char = (char*)contiguous_vector;
+  char* noncontiguous_vector_char = (char*)noncontiguous_vector;
+
   for (i = 0; i < type_c->block_count * count; i++) {
-    memcpy(noncontiguous_vector,
-           contiguous_vector, type_c->block_length * type_c->size_oldtype);
-    contiguous_vector += type_c->block_length*type_c->size_oldtype;
-    noncontiguous_vector += type_c->block_stride*type_c->size_oldtype;
+    memcpy(noncontiguous_vector_char,
+           contiguous_vector_char, type_c->block_length * type_c->size_oldtype);
+
+    contiguous_vector_char += type_c->block_length*type_c->size_oldtype;
+    noncontiguous_vector_char += type_c->block_stride*type_c->size_oldtype;
   }
 }
 
 /*
  * Create a Sub type vector to be able to serialize and unserialize it
- * the structre s_smpi_mpi_vector_t is derived from s_smpi_subtype which
+ * the structure s_smpi_mpi_vector_t is derived from s_smpi_subtype which
  * required the functions unserialize and serialize
  *
  */
@@ -235,6 +257,7 @@ s_smpi_mpi_vector_t* smpi_datatype_vector_create( int block_stride,
   s_smpi_mpi_vector_t *new_t= xbt_new(s_smpi_mpi_vector_t,1);
   new_t->base.serialize = &serialize_vector;
   new_t->base.unserialize = &unserialize_vector;
+  new_t->base.subtype_free = &free_vector;
   new_t->block_stride = block_stride;
   new_t->block_length = block_length;
   new_t->block_count = block_count;
@@ -243,19 +266,22 @@ s_smpi_mpi_vector_t* smpi_datatype_vector_create( int block_stride,
   return new_t;
 }
 
-void smpi_datatype_create(MPI_Datatype* new_type, int size, int has_subtype,
+void smpi_datatype_create(MPI_Datatype* new_type, int size,int extent, int has_subtype,
                           void *struct_type, int flags){
   MPI_Datatype new_t= xbt_new(s_smpi_mpi_datatype_t,1);
-  new_t->size=size;
-  new_t->has_subtype=has_subtype;
-  new_t->lb=0;
-  new_t->ub=size;
-  new_t->flags=flags;
-  new_t->substruct=struct_type;
+  new_t->size = size;
+  new_t->has_subtype = has_subtype;
+  new_t->lb = 0;
+  new_t->ub = extent;
+  new_t->flags = flags;
+  new_t->substruct = struct_type;
   *new_type = new_t;
 }
 
 void smpi_datatype_free(MPI_Datatype* type){
+  if ((*type)->has_subtype == 1){
+    ((s_smpi_subtype_t *)(*type)->substruct)->subtype_free(type);  
+  }
   xbt_free(*type);
 }
 
@@ -266,7 +292,8 @@ int smpi_datatype_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new
     retval = MPI_ERR_TYPE;
   } else {
     smpi_datatype_create(new_type, count *
-                         smpi_datatype_size(old_type),1,NULL, DT_FLAG_CONTIGUOUS);
+                         smpi_datatype_size(old_type),count *
+                         smpi_datatype_size(old_type),0,NULL, DT_FLAG_CONTIGUOUS);
     retval=MPI_SUCCESS;
   }
   return retval;
@@ -280,6 +307,8 @@ int smpi_datatype_vector(int count, int blocklen, int stride, MPI_Datatype old_t
     retval = MPI_ERR_TYPE;
   } else {
     if(stride != blocklen){
+if (old_type->has_subtype == 1)
+      XBT_WARN("vector contains a complex type - not yet handled");
       s_smpi_mpi_vector_t* subtype = smpi_datatype_vector_create( stride,
                                                                   blocklen,
                                                                   count,
@@ -288,6 +317,7 @@ int smpi_datatype_vector(int count, int blocklen, int stride, MPI_Datatype old_t
 
       smpi_datatype_create(new_type, count * (blocklen) *
                            smpi_datatype_size(old_type),
+                          ((count -1) * stride + blocklen) * smpi_datatype_size(old_type),
                            1,
                            subtype,
                            DT_FLAG_VECTOR);
@@ -295,17 +325,109 @@ int smpi_datatype_vector(int count, int blocklen, int stride, MPI_Datatype old_t
     }else{
       /* in this situation the data are contignous thus it's not
        * required to serialize and unserialize it*/
-      smpi_datatype_create(new_type, count * (blocklen) *
+      smpi_datatype_create(new_type, count * blocklen *
+                           smpi_datatype_size(old_type), ((count -1) * stride + blocklen)*
                            smpi_datatype_size(old_type),
                            0,
                            NULL,
-                           DT_FLAG_VECTOR);
+                           DT_FLAG_VECTOR|DT_FLAG_CONTIGUOUS);
       retval=MPI_SUCCESS;
     }
   }
   return retval;
 }
 
+void free_vector(MPI_Datatype* d){
+}
+
+/*
+Hvector Implementation - Vector with stride in bytes
+*/
+
+
+/*
+ *  Copies noncontiguous data into contiguous memory.
+ *  @param contiguous_hvector - output hvector
+ *  @param noncontiguous_hvector - input hvector
+ *  @param type - pointer contening :
+ *      - stride - stride of between noncontiguous data, in bytes
+ *      - block_length - the width or height of blocked matrix
+ *      - count - the number of rows of matrix
+ */
+void serialize_hvector( const void *noncontiguous_hvector,
+                       void *contiguous_hvector,
+                       size_t count,
+                       void *type)
+{
+  s_smpi_mpi_hvector_t* type_c = (s_smpi_mpi_hvector_t*)type;
+  int i;
+  char* contiguous_vector_char = (char*)contiguous_hvector;
+  char* noncontiguous_vector_char = (char*)noncontiguous_hvector;
+
+  for (i = 0; i < type_c->block_count * count; i++) {
+    memcpy(contiguous_vector_char,
+           noncontiguous_vector_char, type_c->block_length * type_c->size_oldtype);
+
+    contiguous_vector_char += type_c->block_length*type_c->size_oldtype;
+    noncontiguous_vector_char += type_c->block_stride;
+  }
+}
+/*
+ *  Copies contiguous data into noncontiguous memory.
+ *  @param noncontiguous_vector - output hvector
+ *  @param contiguous_vector - input hvector
+ *  @param type - pointer contening :
+ *      - stride - stride of between noncontiguous data, in bytes
+ *      - block_length - the width or height of blocked matrix
+ *      - count - the number of rows of matrix
+ */
+void unserialize_hvector( const void *contiguous_vector,
+                         void *noncontiguous_vector,
+                         size_t count,
+                         void *type)
+{
+  s_smpi_mpi_hvector_t* type_c = (s_smpi_mpi_hvector_t*)type;
+  int i;
+
+  char* contiguous_vector_char = (char*)contiguous_vector;
+  char* noncontiguous_vector_char = (char*)noncontiguous_vector;
+
+  for (i = 0; i < type_c->block_count * count; i++) {
+    memcpy(noncontiguous_vector_char,
+           contiguous_vector_char, type_c->block_length * type_c->size_oldtype);
+
+    contiguous_vector_char += type_c->block_length*type_c->size_oldtype;
+    noncontiguous_vector_char += type_c->block_stride;
+  }
+}
+
+/*
+ * Create a Sub type vector to be able to serialize and unserialize it
+ * the structure s_smpi_mpi_vector_t is derived from s_smpi_subtype which
+ * required the functions unserialize and serialize
+ *
+ */
+s_smpi_mpi_hvector_t* smpi_datatype_hvector_create( MPI_Aint block_stride,
+                                                  int block_length,
+                                                  int block_count,
+                                                  MPI_Datatype old_type,
+                                                  int size_oldtype){
+  s_smpi_mpi_hvector_t *new_t= xbt_new(s_smpi_mpi_hvector_t,1);
+  new_t->base.serialize = &serialize_hvector;
+  new_t->base.unserialize = &unserialize_hvector;
+  new_t->base.subtype_free = &free_hvector;
+  new_t->block_stride = block_stride;
+  new_t->block_length = block_length;
+  new_t->block_count = block_count;
+  new_t->old_type = old_type;
+  new_t->size_oldtype = size_oldtype;
+  return new_t;
+}
+
+//do nothing for vector types
+void free_hvector(MPI_Datatype* d){
+}
+
 int smpi_datatype_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type)
 {
   int retval;
@@ -313,83 +435,447 @@ int smpi_datatype_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype
   if ((old_type->flags & DT_FLAG_COMMITED) != DT_FLAG_COMMITED) {
     retval = MPI_ERR_TYPE;
   } else {
-    /*FIXME: as for the vector the data should be serialized and
-     * unserialized moreover a structure derived from s_smpi_subtype should
-     * be created*/
-    smpi_datatype_create(new_type, count * ((blocklen *
-                                             smpi_datatype_size(old_type))+stride),
-                         0,
-                         NULL,
-                         DT_FLAG_VECTOR);
-    retval=MPI_SUCCESS;
+if (old_type->has_subtype == 1)
+      XBT_WARN("hvector contains a complex type - not yet handled");
+    if(stride != blocklen*smpi_datatype_size(old_type)){
+      s_smpi_mpi_hvector_t* subtype = smpi_datatype_hvector_create( stride,
+                                                                    blocklen,
+                                                                    count,
+                                                                    old_type,
+                                                                    smpi_datatype_size(old_type));
+
+      smpi_datatype_create(new_type, count * blocklen *
+                           smpi_datatype_size(old_type), (count-1) * stride + blocklen *
+                           smpi_datatype_size(old_type),
+                           1,
+                           subtype,
+                           DT_FLAG_VECTOR);
+      retval=MPI_SUCCESS;
+    }else{
+      smpi_datatype_create(new_type, count * blocklen *
+                                               smpi_datatype_size(old_type),count * blocklen *
+                                               smpi_datatype_size(old_type),
+                                              0,
+                                              NULL,
+                                              DT_FLAG_VECTOR|DT_FLAG_CONTIGUOUS);
+      retval=MPI_SUCCESS;
+    }
   }
   return retval;
 }
 
 
+/*
+Indexed Implementation
+*/
+
+/*
+ *  Copies noncontiguous data into contiguous memory.
+ *  @param contiguous_indexed - output indexed
+ *  @param noncontiguous_indexed - input indexed
+ *  @param type - pointer contening :
+ *      - block_lengths - the width or height of blocked matrix
+ *      - block_indices - indices of each data, in element
+ *      - count - the number of rows of matrix
+ */
+void serialize_indexed( const void *noncontiguous_indexed,
+                       void *contiguous_indexed,
+                       size_t count,
+                       void *type)
+{
+  s_smpi_mpi_indexed_t* type_c = (s_smpi_mpi_indexed_t*)type;
+  int i,j;
+  char* contiguous_indexed_char = (char*)contiguous_indexed;
+  char* noncontiguous_indexed_char = (char*)noncontiguous_indexed;
+  for(j=0; j<count;j++){
+    for (i = 0; i < type_c->block_count; i++) {
+      memcpy(contiguous_indexed_char,
+             noncontiguous_indexed_char, type_c->block_lengths[i] * type_c->size_oldtype);
+
+      contiguous_indexed_char += type_c->block_lengths[i]*type_c->size_oldtype;
+      if (i<type_c->block_count-1)noncontiguous_indexed_char = (char*)noncontiguous_indexed + type_c->block_indices[i+1]*type_c->size_oldtype;
+      else noncontiguous_indexed_char += type_c->block_lengths[i]*type_c->size_oldtype;
+    }
+    noncontiguous_indexed=(void*)noncontiguous_indexed_char;
+  }
+}
+/*
+ *  Copies contiguous data into noncontiguous memory.
+ *  @param noncontiguous_indexed - output indexed
+ *  @param contiguous_indexed - input indexed
+ *  @param type - pointer contening :
+ *      - block_lengths - the width or height of blocked matrix
+ *      - block_indices - indices of each data, in element
+ *      - count - the number of rows of matrix
+ */
+void unserialize_indexed( const void *contiguous_indexed,
+                         void *noncontiguous_indexed,
+                         size_t count,
+                         void *type)
+{
+  s_smpi_mpi_indexed_t* type_c = (s_smpi_mpi_indexed_t*)type;
+  int i,j;
+
+  char* contiguous_indexed_char = (char*)contiguous_indexed;
+  char* noncontiguous_indexed_char = (char*)noncontiguous_indexed;
+  for(j=0; j<count;j++){
+    for (i = 0; i < type_c->block_count; i++) {
+      memcpy(noncontiguous_indexed_char,
+             contiguous_indexed_char, type_c->block_lengths[i] * type_c->size_oldtype);
+
+      contiguous_indexed_char += type_c->block_lengths[i]*type_c->size_oldtype;
+      if (i<type_c->block_count-1)noncontiguous_indexed_char = (char*)noncontiguous_indexed + type_c->block_indices[i+1]*type_c->size_oldtype;
+      else noncontiguous_indexed_char += type_c->block_lengths[i]*type_c->size_oldtype;
+    }
+    noncontiguous_indexed=(void*)noncontiguous_indexed_char;
+  }
+}
+
+void free_indexed(MPI_Datatype* type){
+  xbt_free(((s_smpi_mpi_indexed_t *)(*type)->substruct)->block_lengths);
+  xbt_free(((s_smpi_mpi_indexed_t *)(*type)->substruct)->block_indices);
+}
+
+/*
+ * Create a Sub type indexed to be able to serialize and unserialize it
+ * the structure s_smpi_mpi_indexed_t is derived from s_smpi_subtype which
+ * required the functions unserialize and serialize
+ */
+s_smpi_mpi_indexed_t* smpi_datatype_indexed_create( int* block_lengths,
+                                                  int* block_indices,
+                                                  int block_count,
+                                                  MPI_Datatype old_type,
+                                                  int size_oldtype){
+  s_smpi_mpi_indexed_t *new_t= xbt_new(s_smpi_mpi_indexed_t,1);
+  new_t->base.serialize = &serialize_indexed;
+  new_t->base.unserialize = &unserialize_indexed;
+  new_t->base.subtype_free = &free_indexed;
+ //TODO : add a custom function for each time to clean these 
+  new_t->block_lengths= xbt_new(int, block_count);
+  new_t->block_indices= xbt_new(int, block_count);
+  int i;
+  for(i=0;i<block_count;i++){
+    new_t->block_lengths[i]=block_lengths[i];
+    new_t->block_indices[i]=block_indices[i];
+  }
+  new_t->block_count = block_count;
+  new_t->old_type = old_type;
+  new_t->size_oldtype = size_oldtype;
+  return new_t;
+}
+
+
 int smpi_datatype_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type)
 {
   int i;
   int retval;
   int size = 0;
+  int contiguous=1;
   for(i=0; i< count; i++){
     if   (blocklens[i]<=0)
       return MPI_ERR_ARG;
     size += blocklens[i];
+
+    if ( (i< count -1) && (indices[i]+blocklens[i] != indices[i+1]) )contiguous=0;
   }
   if ((old_type->flags & DT_FLAG_COMMITED) != DT_FLAG_COMMITED) {
     retval = MPI_ERR_TYPE;
   } else {
-    /*FIXME: as for the vector the data should be serialized and
-     * unserialized moreover a structure derived from s_smpi_subtype should
-     * be created*/
-    smpi_datatype_create(new_type,  (size) *
-                         smpi_datatype_size(old_type),0, NULL, DT_FLAG_DATA);
+
+    if (old_type->has_subtype == 1)
+      XBT_WARN("indexed contains a complex type - not yet handled");
+
+    if(!contiguous){
+      s_smpi_mpi_indexed_t* subtype = smpi_datatype_indexed_create( blocklens,
+                                                                    indices,
+                                                                    count,
+                                                                    old_type,
+                                                                    smpi_datatype_size(old_type));
+
+      smpi_datatype_create(new_type,  size *
+                           smpi_datatype_size(old_type),(indices[count-1]+blocklens[count-1])*smpi_datatype_size(old_type),1, subtype, DT_FLAG_DATA);
+}else{
+      smpi_datatype_create(new_type,  size *
+                           smpi_datatype_size(old_type),size *
+                           smpi_datatype_size(old_type),0, NULL, DT_FLAG_DATA|DT_FLAG_CONTIGUOUS);
+}
     retval=MPI_SUCCESS;
   }
   return retval;
 }
 
+
+/*
+Hindexed Implementation - Indexed with indices in bytes 
+*/
+
+/*
+ *  Copies noncontiguous data into contiguous memory.
+ *  @param contiguous_hindexed - output hindexed
+ *  @param noncontiguous_hindexed - input hindexed
+ *  @param type - pointer contening :
+ *      - block_lengths - the width or height of blocked matrix
+ *      - block_indices - indices of each data, in bytes
+ *      - count - the number of rows of matrix
+ */
+void serialize_hindexed( const void *noncontiguous_hindexed,
+                       void *contiguous_hindexed,
+                       size_t count,
+                       void *type)
+{
+  s_smpi_mpi_hindexed_t* type_c = (s_smpi_mpi_hindexed_t*)type;
+  int i,j;
+  char* contiguous_hindexed_char = (char*)contiguous_hindexed;
+  char* noncontiguous_hindexed_char = (char*)noncontiguous_hindexed;
+  for(j=0; j<count;j++){
+    for (i = 0; i < type_c->block_count; i++) {
+      memcpy(contiguous_hindexed_char,
+             noncontiguous_hindexed_char, type_c->block_lengths[i] * type_c->size_oldtype);
+
+      contiguous_hindexed_char += type_c->block_lengths[i]*type_c->size_oldtype;
+      if (i<type_c->block_count-1)noncontiguous_hindexed_char = (char*)noncontiguous_hindexed + type_c->block_indices[i+1];
+      else noncontiguous_hindexed_char += type_c->block_lengths[i]*type_c->size_oldtype;
+    }
+    noncontiguous_hindexed=(void*)noncontiguous_hindexed_char;
+  }
+}
+/*
+ *  Copies contiguous data into noncontiguous memory.
+ *  @param noncontiguous_hindexed - output hindexed
+ *  @param contiguous_hindexed - input hindexed
+ *  @param type - pointer contening :
+ *      - block_lengths - the width or height of blocked matrix
+ *      - block_indices - indices of each data, in bytes
+ *      - count - the number of rows of matrix
+ */
+void unserialize_hindexed( const void *contiguous_hindexed,
+                         void *noncontiguous_hindexed,
+                         size_t count,
+                         void *type)
+{
+  s_smpi_mpi_hindexed_t* type_c = (s_smpi_mpi_hindexed_t*)type;
+  int i,j;
+
+  char* contiguous_hindexed_char = (char*)contiguous_hindexed;
+  char* noncontiguous_hindexed_char = (char*)noncontiguous_hindexed;
+  for(j=0; j<count;j++){
+    for (i = 0; i < type_c->block_count; i++) {
+      memcpy(noncontiguous_hindexed_char,
+             contiguous_hindexed_char, type_c->block_lengths[i] * type_c->size_oldtype);
+
+      contiguous_hindexed_char += type_c->block_lengths[i]*type_c->size_oldtype;
+      if (i<type_c->block_count-1)noncontiguous_hindexed_char = (char*)noncontiguous_hindexed + type_c->block_indices[i+1];
+      else noncontiguous_hindexed_char += type_c->block_lengths[i]*type_c->size_oldtype;
+    }
+    noncontiguous_hindexed=(void*)noncontiguous_hindexed_char;
+  }
+}
+
+void free_hindexed(MPI_Datatype* type){
+  xbt_free(((s_smpi_mpi_hindexed_t *)(*type)->substruct)->block_lengths);
+  xbt_free(((s_smpi_mpi_hindexed_t *)(*type)->substruct)->block_indices);
+}
+
+/*
+ * Create a Sub type hindexed to be able to serialize and unserialize it
+ * the structure s_smpi_mpi_hindexed_t is derived from s_smpi_subtype which
+ * required the functions unserialize and serialize
+ */
+s_smpi_mpi_hindexed_t* smpi_datatype_hindexed_create( int* block_lengths,
+                                                  MPI_Aint* block_indices,
+                                                  int block_count,
+                                                  MPI_Datatype old_type,
+                                                  int size_oldtype){
+  s_smpi_mpi_hindexed_t *new_t= xbt_new(s_smpi_mpi_hindexed_t,1);
+  new_t->base.serialize = &serialize_hindexed;
+  new_t->base.unserialize = &unserialize_hindexed;
+  new_t->base.subtype_free = &free_hindexed;
+ //TODO : add a custom function for each time to clean these 
+  new_t->block_lengths= xbt_new(int, block_count);
+  new_t->block_indices= xbt_new(MPI_Aint, block_count);
+  int i;
+  for(i=0;i<block_count;i++){
+    new_t->block_lengths[i]=block_lengths[i];
+    new_t->block_indices[i]=block_indices[i];
+  }
+  new_t->block_count = block_count;
+  new_t->old_type = old_type;
+  new_t->size_oldtype = size_oldtype;
+  return new_t;
+}
+
+
 int smpi_datatype_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type)
 {
   int i;
   int retval;
   int size = 0;
+  int contiguous=1;
   for(i=0; i< count; i++){
     if   (blocklens[i]<=0)
       return MPI_ERR_ARG;
     size += blocklens[i];
+
+
+    if ( (i< count -1) && (indices[i]+blocklens[i]*smpi_datatype_size(old_type) != indices[i+1]) )contiguous=0;
   }
   if ((old_type->flags & DT_FLAG_COMMITED) != DT_FLAG_COMMITED) {
     retval = MPI_ERR_TYPE;
   } else {
-    /*FIXME: as for the vector the data should be serialized and
-     * unserialized moreover a structure derived from s_smpi_subtype should
-     * be created*/
-    smpi_datatype_create(new_type,(size * smpi_datatype_size(old_type)), 0,NULL, DT_FLAG_DATA);
+    if (old_type->has_subtype == 1)
+      XBT_WARN("hindexed contains a complex type - not yet handled");
+
+    if(!contiguous){
+      s_smpi_mpi_hindexed_t* subtype = smpi_datatype_hindexed_create( blocklens,
+                                                                    indices,
+                                                                    count,
+                                                                    old_type,
+                                                                    smpi_datatype_size(old_type));
+
+      smpi_datatype_create(new_type,  size *
+                           smpi_datatype_size(old_type),indices[count-1]+blocklens[count-1]*smpi_datatype_size(old_type)
+                           ,1, subtype, DT_FLAG_DATA);
+    }else{
+      smpi_datatype_create(new_type,  size *
+                           smpi_datatype_size(old_type),size *
+                           smpi_datatype_size(old_type),0, NULL, DT_FLAG_DATA|DT_FLAG_CONTIGUOUS);
+    }
     retval=MPI_SUCCESS;
   }
   return retval;
 }
 
-int smpi_datatype_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type)
+
+/*
+struct Implementation - Indexed with indices in bytes 
+*/
+
+/*
+ *  Copies noncontiguous data into contiguous memory.
+ *  @param contiguous_struct - output struct
+ *  @param noncontiguous_struct - input struct
+ *  @param type - pointer contening :
+ *      - stride - stride of between noncontiguous data
+ *      - block_length - the width or height of blocked matrix
+ *      - count - the number of rows of matrix
+ */
+void serialize_struct( const void *noncontiguous_struct,
+                       void *contiguous_struct,
+                       size_t count,
+                       void *type)
+{
+  s_smpi_mpi_struct_t* type_c = (s_smpi_mpi_struct_t*)type;
+  int i,j;
+  char* contiguous_struct_char = (char*)contiguous_struct;
+  char* noncontiguous_struct_char = (char*)noncontiguous_struct;
+  for(j=0; j<count;j++){
+    for (i = 0; i < type_c->block_count; i++) {
+      memcpy(contiguous_struct_char,
+             noncontiguous_struct_char, type_c->block_lengths[i] * smpi_datatype_size(type_c->old_types[i]));
+      contiguous_struct_char += type_c->block_lengths[i]*smpi_datatype_size(type_c->old_types[i]);
+      if (i<type_c->block_count-1)noncontiguous_struct_char = (char*)noncontiguous_struct + type_c->block_indices[i+1];
+      else noncontiguous_struct_char += type_c->block_lengths[i]*smpi_datatype_size(type_c->old_types[i]);//let's hope this is MPI_UB ?
+    }
+    noncontiguous_struct=(void*)noncontiguous_struct_char;
+  }
+}
+/*
+ *  Copies contiguous data into noncontiguous memory.
+ *  @param noncontiguous_struct - output struct
+ *  @param contiguous_struct - input struct
+ *  @param type - pointer contening :
+ *      - stride - stride of between noncontiguous data
+ *      - block_length - the width or height of blocked matrix
+ *      - count - the number of rows of matrix
+ */
+void unserialize_struct( const void *contiguous_struct,
+                         void *noncontiguous_struct,
+                         size_t count,
+                         void *type)
 {
+  s_smpi_mpi_struct_t* type_c = (s_smpi_mpi_struct_t*)type;
+  int i,j;
+
+  char* contiguous_struct_char = (char*)contiguous_struct;
+  char* noncontiguous_struct_char = (char*)noncontiguous_struct;
+  for(j=0; j<count;j++){
+    for (i = 0; i < type_c->block_count; i++) {
+      memcpy(noncontiguous_struct_char,
+             contiguous_struct_char, type_c->block_lengths[i] * smpi_datatype_size(type_c->old_types[i]));
+      contiguous_struct_char += type_c->block_lengths[i]*smpi_datatype_size(type_c->old_types[i]);
+      if (i<type_c->block_count-1)noncontiguous_struct_char =  (char*)noncontiguous_struct + type_c->block_indices[i+1];
+      else noncontiguous_struct_char += type_c->block_lengths[i]*smpi_datatype_size(type_c->old_types[i]);
+    }
+    noncontiguous_struct=(void*)noncontiguous_struct_char;
+    
+  }
+}
+
+void free_struct(MPI_Datatype* type){
+  xbt_free(((s_smpi_mpi_struct_t *)(*type)->substruct)->block_lengths);
+  xbt_free(((s_smpi_mpi_struct_t *)(*type)->substruct)->block_indices);
+  xbt_free(((s_smpi_mpi_struct_t *)(*type)->substruct)->old_types);
+}
+
+/*
+ * Create a Sub type struct to be able to serialize and unserialize it
+ * the structure s_smpi_mpi_struct_t is derived from s_smpi_subtype which
+ * required the functions unserialize and serialize
+ */
+s_smpi_mpi_struct_t* smpi_datatype_struct_create( int* block_lengths,
+                                                  MPI_Aint* block_indices,
+                                                  int block_count,
+                                                  MPI_Datatype* old_types){
+  s_smpi_mpi_struct_t *new_t= xbt_new(s_smpi_mpi_struct_t,1);
+  new_t->base.serialize = &serialize_struct;
+  new_t->base.unserialize = &unserialize_struct;
+  new_t->base.subtype_free = &free_struct;
+ //TODO : add a custom function for each time to clean these 
+  new_t->block_lengths= xbt_new(int, block_count);
+  new_t->block_indices= xbt_new(MPI_Aint, block_count);
+  new_t->old_types=  xbt_new(MPI_Datatype, block_count);
   int i;
-  size_t size; //Khalid added this
+  for(i=0;i<block_count;i++){
+    new_t->block_lengths[i]=block_lengths[i];
+    new_t->block_indices[i]=block_indices[i];
+    new_t->old_types[i]=old_types[i];
+  }
+  //new_t->block_lengths = block_lengths;
+  //new_t->block_indices = block_indices;
+  new_t->block_count = block_count;
+  //new_t->old_types = old_types;
+  return new_t;
+}
+
 
+int smpi_datatype_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type)
+{
+  int i;
+  size_t size = 0;
+  int contiguous=1;
+  size = 0;
   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;
+    if (old_types[i]->has_subtype == 1)
+      XBT_WARN("Struct contains a complex type - not yet handled");
     size += blocklens[i]*smpi_datatype_size(old_types[i]);
+
+    if ( (i< count -1) && (indices[i]+blocklens[i]*smpi_datatype_size(old_types[i]) != indices[i+1]) )contiguous=0;
+  }
+
+  if(!contiguous){
+    s_smpi_mpi_struct_t* subtype = smpi_datatype_struct_create( blocklens,
+                                                              indices,
+                                                              count,
+                                                              old_types);
+
+    smpi_datatype_create(new_type,  size, indices[count-1] + blocklens[count-1]*smpi_datatype_size(old_types[count-1]),1, subtype, DT_FLAG_DATA);
+  }else{
+    smpi_datatype_create(new_type,  size, indices[count-1] + blocklens[count-1]*smpi_datatype_size(old_types[count-1]),0, NULL, DT_FLAG_DATA|DT_FLAG_CONTIGUOUS);
   }
-  /*FIXME: as for the vector the data should be serialized and
-   * unserialized moreover a structure derived from s_smpi_subtype should
-   * be created*/
-  smpi_datatype_create(new_type, size,
-                       0, NULL,
-                       DT_FLAG_DATA);
   return MPI_SUCCESS;
 }