Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
MPI_Topo -> c++
authordegomme <augustin.degomme@unibas.ch>
Mon, 6 Mar 2017 17:51:16 +0000 (18:51 +0100)
committerdegomme <augustin.degomme@unibas.ch>
Tue, 7 Mar 2017 10:10:15 +0000 (11:10 +0100)
include/smpi/forward.hpp
include/smpi/smpi.h
src/smpi/private.h
src/smpi/smpi_comm.cpp
src/smpi/smpi_pmpi.cpp
src/smpi/smpi_topo.cpp
src/smpi/smpi_topo.hpp [new file with mode: 0644]
tools/cmake/DefinePackages.cmake

index e24b1f0..81292da 100644 (file)
@@ -15,17 +15,29 @@ namespace SMPI {
 
 class Group;
 class Comm;
 
 class Group;
 class Comm;
+class Topo;
+class Cart;
+class Graph;
+class Dist_Graph;
 
 }
 }
 
 typedef simgrid::SMPI::Group SMPI_Group;
 typedef simgrid::SMPI::Comm SMPI_Comm;
 
 }
 }
 
 typedef simgrid::SMPI::Group SMPI_Group;
 typedef simgrid::SMPI::Comm SMPI_Comm;
+typedef simgrid::SMPI::Topo SMPI_Topology;
+typedef simgrid::SMPI::Graph SMPI_Graph_topology;
+typedef simgrid::SMPI::Cart SMPI_Cart_topology;
+typedef simgrid::SMPI::Dist_Graph SMPI_Dist_Graph_topology;
 
 #else
 
 typedef struct SMPI_Group SMPI_Group;
 typedef struct SMPI_Comm SMPI_Comm;
 
 #else
 
 typedef struct SMPI_Group SMPI_Group;
 typedef struct SMPI_Comm SMPI_Comm;
+typedef struct SMPI_Topology SMPI_Topology;
+typedef struct SMPI_Graph_topology SMPI_Graph_topology;
+typedef struct SMPI_Cart_topology SMPI_Cart_topology;
+typedef struct SMPI_Dist_Graph_topology SMPI_Dist_Graph_topology;
 
 #endif
 
 
 #endif
 
index 6214488..c058e3a 100644 (file)
@@ -363,8 +363,7 @@ XBT_PUBLIC_DATA( MPI_Op ) MPI_BXOR;
 //For accumulate
 XBT_PUBLIC_DATA( MPI_Op ) MPI_REPLACE;
 
 //For accumulate
 XBT_PUBLIC_DATA( MPI_Op ) MPI_REPLACE;
 
-struct s_smpi_mpi_topology;
-typedef struct s_smpi_mpi_topology *MPI_Topology;
+typedef SMPI_Topology *MPI_Topology;
 
 typedef SMPI_Group* MPI_Group;
 
 
 typedef SMPI_Group* MPI_Group;
 
index d600eb7..e2e68fc 100644 (file)
@@ -11,6 +11,7 @@
 #include "smpi/smpi.h"
 #include "src/smpi/smpi_group.hpp"
 #include "src/smpi/smpi_comm.hpp"
 #include "smpi/smpi.h"
 #include "src/smpi/smpi_group.hpp"
 #include "src/smpi/smpi_comm.hpp"
+#include "src/smpi/smpi_topo.hpp"
 #include "src/include/smpi/smpi_interface.h"
 #include "src/instr/instr_private.h"
 #include "src/internal_config.h"
 #include "src/include/smpi/smpi_interface.h"
 #include "src/instr/instr_private.h"
 #include "src/internal_config.h"
@@ -135,32 +136,11 @@ XBT_PRIVATE int smpi_process_finalized();
 XBT_PRIVATE int smpi_process_initialized();
 XBT_PRIVATE void smpi_process_mark_as_initialized();
 
 XBT_PRIVATE int smpi_process_initialized();
 XBT_PRIVATE void smpi_process_mark_as_initialized();
 
-struct s_smpi_mpi_cart_topology;
-typedef struct s_smpi_mpi_cart_topology *MPIR_Cart_Topology;
-
-struct s_smpi_mpi_graph_topology;
-typedef struct s_smpi_mpi_graph_topology *MPIR_Graph_Topology;
-
-struct s_smpi_dist_graph_topology;
-typedef struct s_smpi_dist_graph_topology *MPIR_Dist_Graph_Topology;
-
-// MPI_Topology defined in smpi.h, as it is public
-
-XBT_PRIVATE void smpi_topo_destroy(MPI_Topology topo);
-XBT_PRIVATE MPI_Topology smpi_topo_create(MPIR_Topo_type kind);
-XBT_PRIVATE void smpi_cart_topo_destroy(MPIR_Cart_Topology cart);
-XBT_PRIVATE MPI_Topology smpi_cart_topo_create(int ndims);
-XBT_PRIVATE int smpi_mpi_cart_create(MPI_Comm comm_old, int ndims, int dims[], int periods[], int reorder,
-                                     MPI_Comm *comm_cart);
-XBT_PRIVATE int smpi_mpi_cart_sub(MPI_Comm comm, const int remain_dims[], MPI_Comm *newcomm);
-XBT_PRIVATE int smpi_mpi_cart_coords(MPI_Comm comm, int rank, int maxdims, int coords[]);
-XBT_PRIVATE int smpi_mpi_cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords);
-XBT_PRIVATE int smpi_mpi_cart_rank(MPI_Comm comm, int* coords, int* rank);
-XBT_PRIVATE int smpi_mpi_cart_shift(MPI_Comm comm, int direction, int disp, int *rank_source, int *rank_dest);
-XBT_PRIVATE int smpi_mpi_cartdim_get(MPI_Comm comm, int *ndims);
-XBT_PRIVATE int smpi_mpi_dims_create(int nnodes, int ndims, int dims[]);
-XBT_PRIVATE void smpi_graph_topo_destroy(MPIR_Graph_Topology cart);
-XBT_PRIVATE void smpi_dist_graph_topo_destroy(MPIR_Dist_Graph_Topology cart);
+typedef SMPI_Cart_topology *MPIR_Cart_Topology;
+
+typedef SMPI_Graph_topology *MPIR_Graph_Topology;
+
+typedef SMPI_Dist_Graph_topology *MPIR_Dist_Graph_Topology;
 
 XBT_PRIVATE smpi_process_data_t smpi_process_data();
 XBT_PRIVATE smpi_process_data_t smpi_process_remote_data(int index);
 
 XBT_PRIVATE smpi_process_data_t smpi_process_data();
 XBT_PRIVATE smpi_process_data_t smpi_process_remote_data(int index);
index ec4d78b..19e7a8e 100644 (file)
@@ -68,7 +68,7 @@ void Comm::destroy()
 {
   if (this == MPI_COMM_UNINITIALIZED)
     return smpi_process_comm_world()->destroy();
 {
   if (this == MPI_COMM_UNINITIALIZED)
     return smpi_process_comm_world()->destroy();
-  smpi_topo_destroy(m_topo); // there's no use count on topos
+  delete m_topo; // there's no use count on topos
   this->unuse();
 }
 
   this->unuse();
 }
 
index 57075df..6002014 100644 (file)
@@ -2306,7 +2306,8 @@ int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periodic, int
   } else if (ndims < 0 || (ndims > 0 && (dims == nullptr || periodic == nullptr)) || comm_cart == nullptr) {
     return MPI_ERR_ARG;
   } else{
   } else if (ndims < 0 || (ndims > 0 && (dims == nullptr || periodic == nullptr)) || comm_cart == nullptr) {
     return MPI_ERR_ARG;
   } else{
-    return smpi_mpi_cart_create(comm_old, ndims, dims, periodic, reorder, comm_cart);
+    new simgrid::SMPI::Cart(comm_old, ndims, dims, periodic, reorder, comm_cart);
+    return MPI_SUCCESS;
   }
 }
 
   }
 }
 
@@ -2317,7 +2318,11 @@ int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
   if (coords == nullptr) {
     return MPI_ERR_ARG;
   }
   if (coords == nullptr) {
     return MPI_ERR_ARG;
   }
-  return smpi_mpi_cart_rank(comm, coords, rank);
+  simgrid::SMPI::Cart* topo = static_cast<simgrid::SMPI::Cart*>(comm->topo());
+  if (topo==nullptr) {
+    return MPI_ERR_ARG;
+  }
+  return topo->rank(coords, rank);
 }
 
 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
 }
 
 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
@@ -2327,7 +2332,11 @@ int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* d
   if (source == nullptr || dest == nullptr || direction < 0 ) {
     return MPI_ERR_ARG;
   }
   if (source == nullptr || dest == nullptr || direction < 0 ) {
     return MPI_ERR_ARG;
   }
-  return smpi_mpi_cart_shift(comm, direction, displ, source, dest);
+  simgrid::SMPI::Cart* topo = static_cast<simgrid::SMPI::Cart*>(comm->topo());
+  if (topo==nullptr) {
+    return MPI_ERR_ARG;
+  }
+  return topo->shift(direction, displ, source, dest);
 }
 
 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
 }
 
 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
@@ -2343,7 +2352,11 @@ int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
   if(coords == nullptr) {
     return MPI_ERR_ARG;
   }
   if(coords == nullptr) {
     return MPI_ERR_ARG;
   }
-  return smpi_mpi_cart_coords(comm, rank, maxdims, coords);
+  simgrid::SMPI::Cart* topo = static_cast<simgrid::SMPI::Cart*>(comm->topo());
+  if (topo==nullptr) {
+    return MPI_ERR_ARG;
+  }
+  return topo->coords(rank, maxdims, coords);
 }
 
 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
 }
 
 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
@@ -2353,7 +2366,11 @@ int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coor
   if(maxdims <= 0 || dims == nullptr || periods == nullptr || coords == nullptr) {
     return MPI_ERR_ARG;
   }
   if(maxdims <= 0 || dims == nullptr || periods == nullptr || coords == nullptr) {
     return MPI_ERR_ARG;
   }
-  return smpi_mpi_cart_get(comm, maxdims, dims, periods, coords);
+  simgrid::SMPI::Cart* topo = static_cast<simgrid::SMPI::Cart*>(comm->topo());
+  if (topo==nullptr) {
+    return MPI_ERR_ARG;
+  }
+  return topo->get(maxdims, dims, periods, coords);
 }
 
 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
 }
 
 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
@@ -2363,7 +2380,11 @@ int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
   if (ndims == nullptr) {
     return MPI_ERR_ARG;
   }
   if (ndims == nullptr) {
     return MPI_ERR_ARG;
   }
-  return smpi_mpi_cartdim_get(comm, ndims);
+  simgrid::SMPI::Cart* topo = static_cast<simgrid::SMPI::Cart*>(comm->topo());
+  if (topo==nullptr) {
+    return MPI_ERR_ARG;
+  }
+  return topo->dim_get(ndims);
 }
 
 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
 }
 
 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
@@ -2373,8 +2394,7 @@ int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
   if (ndims < 1 || nnodes < 1) {
     return MPI_ERR_DIMS;
   }
   if (ndims < 1 || nnodes < 1) {
     return MPI_ERR_DIMS;
   }
-
-  return smpi_mpi_dims_create(nnodes, ndims, dims);
+  return simgrid::SMPI::Dims_create(nnodes, ndims, dims);
 }
 
 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
 }
 
 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
@@ -2384,7 +2404,14 @@ int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
   if (comm_new == nullptr) {
     return MPI_ERR_ARG;
   }
   if (comm_new == nullptr) {
     return MPI_ERR_ARG;
   }
-  return smpi_mpi_cart_sub(comm, remain_dims, comm_new);
+  simgrid::SMPI::Cart* topo = static_cast<simgrid::SMPI::Cart*>(comm->topo());
+  if (topo==nullptr) {
+    return MPI_ERR_ARG;
+  }
+  simgrid::SMPI::Cart* cart = topo->sub(remain_dims, comm_new);
+  if(cart==nullptr)
+    return  MPI_ERR_ARG;
+  return MPI_SUCCESS;
 }
 
 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
 }
 
 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
index d7912c7..39a1a62 100644 (file)
 #include <vector>
 #include <math.h>
 
 #include <vector>
 #include <math.h>
 
-typedef struct s_smpi_mpi_cart_topology {
-  int nnodes;
-  int ndims;
-  int *dims;
-  int *periodic;
-  int *position;
-} s_smpi_mpi_cart_topology_t;
-
-typedef struct s_smpi_mpi_graph_topology {
-  int nnodes;
-  int nedges;
-  int *index;
-  int *edges;
-} s_smpi_mpi_graph_topology_t;
-
-typedef struct s_smpi_dist_graph_topology {
-  int indegree;
-  int *in;
-  int *in_weights;
-  int outdegree;
-  int *out;
-  int *out_weights;
-  int is_weighted;
-} s_smpi_mpi_dist_graph_topology_t;
-
-typedef struct s_smpi_mpi_topology { 
-  MPIR_Topo_type kind;
-  union topo {
-    MPIR_Graph_Topology graph;
-    MPIR_Cart_Topology  cart;
-    MPIR_Dist_Graph_Topology dist_graph;
-  } topo;
-} s_smpi_mpi_topology_t;
-
-void smpi_topo_destroy(MPI_Topology topo) {
-  if(topo == nullptr) {
-    return;
-  }
-  switch (topo->kind) {
-  case MPI_CART:
-    smpi_cart_topo_destroy(topo->topo.cart);
-    break;
-  case MPI_GRAPH:
-    // This topology is not supported by SMPI yet
-    smpi_graph_topo_destroy(topo->topo.graph);
-    break;
-  case MPI_DIST_GRAPH:
-    // This topology is not supported by SMPI yet
-     smpi_dist_graph_topo_destroy(topo->topo.dist_graph);
-    break;
-  default:
-    return;
-  }
-  xbt_free(topo);
-}
+/* static functions */
+static int assignnodes(int ndim, int nfactor, int *pfacts,int **pdims);
+static int getfactors(int num, int *nfators, int **factors);
 
 
-MPI_Topology smpi_topo_create(MPIR_Topo_type kind) {
-  MPI_Topology newTopo = static_cast<MPI_Topology>(xbt_malloc(sizeof(*newTopo)));
-  newTopo->kind = kind;
-  // Allocate and initialize the right topo should be done by the caller
-  return newTopo;
-}
 
 
-void smpi_graph_topo_destroy(MPIR_Graph_Topology graph) {
-  if (graph) {
-    delete[] graph->index;
-    delete[] graph->edges;
-    xbt_free(graph);
-  }
+namespace simgrid{
+namespace SMPI{
+
+
+Graph::~Graph() 
+{
+  delete[] m_index;
+  delete[] m_edges;
 }
 
 }
 
-void smpi_dist_graph_topo_destroy(MPIR_Dist_Graph_Topology dist_graph) {
-  if (dist_graph) {
-    delete[] dist_graph->in;
-    delete[] dist_graph->in_weights;
-    delete[] dist_graph->out;
-    delete[] dist_graph->out_weights;
-    xbt_free(dist_graph);
-  }
+Dist_Graph::~Dist_Graph() 
+{
+  delete[] m_in;
+  delete[] m_in_weights;
+  delete[] m_out;
+  delete[] m_out_weights;
 }
 
 /*******************************************************************************
  * Cartesian topologies
  ******************************************************************************/
 }
 
 /*******************************************************************************
  * Cartesian topologies
  ******************************************************************************/
-void smpi_cart_topo_destroy(MPIR_Cart_Topology cart) {
-  if (cart) {
-    delete[] cart->dims;
-    delete[] cart->periodic;
-    delete[] cart->position;
-    xbt_free(cart);
-  }
+Cart::~Cart() 
+{
+  delete[] m_dims;
+  delete[] m_periodic;
+  delete[] m_position;
 }
 
 }
 
-MPI_Topology smpi_cart_topo_create(int ndims) {
-    MPI_Topology newTopo = smpi_topo_create(MPI_CART);
-    MPIR_Cart_Topology newCart = static_cast<MPIR_Cart_Topology>(xbt_malloc(sizeof(*newCart)));
-    newCart->nnodes = 0;
-    newCart->ndims = ndims;
-    newCart->dims = new int[ndims];
-    newCart->periodic = new int[ndims];
-    newCart->position = new int[ndims];
-    newTopo->topo.cart = newCart;
-    return newTopo;
+Cart::Cart(int ndims)
+{
+  m_nnodes = 0;
+  m_ndims = ndims;
+  m_dims = new int[ndims];
+  m_periodic = new int[ndims];
+  m_position = new int[ndims];
 }
 
 /* reorder is ignored, don't know what would be the consequences of a dumb reordering but neither do I see the point of
  * reordering*/
 }
 
 /* reorder is ignored, don't know what would be the consequences of a dumb reordering but neither do I see the point of
  * reordering*/
-int smpi_mpi_cart_create(MPI_Comm comm_old, int ndims, int dims[], int periods[], int reorder, MPI_Comm *comm_cart) {
-  int retval = MPI_SUCCESS;
-  MPI_Topology newCart;
+Cart::Cart(MPI_Comm comm_old, int ndims, int dims[], int periods[], int reorder, MPI_Comm *comm_cart) : Cart(ndims) {
   MPI_Group newGroup;
   MPI_Group oldGroup;
   int nranks;
   MPI_Group newGroup;
   MPI_Group oldGroup;
   int nranks;
@@ -133,49 +68,46 @@ int smpi_mpi_cart_create(MPI_Comm comm_old, int ndims, int dims[], int periods[]
     }
     if(rank >= newSize) {
       *comm_cart = MPI_COMM_NULL;
     }
     if(rank >= newSize) {
       *comm_cart = MPI_COMM_NULL;
-      return retval;
+      return;
     }
     }
-    newCart = smpi_cart_topo_create(ndims);
     oldGroup = comm_old->group();
     newGroup = new simgrid::SMPI::Group(newSize);
     for (int i = 0 ; i < newSize ; i++) {
       newGroup->set_mapping(oldGroup->index(i), i);
     }
 
     oldGroup = comm_old->group();
     newGroup = new simgrid::SMPI::Group(newSize);
     for (int i = 0 ; i < newSize ; i++) {
       newGroup->set_mapping(oldGroup->index(i), i);
     }
 
-    newCart->topo.cart->nnodes = newSize;
+    m_nnodes = newSize;
 
 
-    //  FIXME : code duplication... See smpi_mpi_cart_coords
+    //  FIXME : code duplication... See coords
     nranks = newSize;
     for (int i=0; i<ndims; i++) {
     nranks = newSize;
     for (int i=0; i<ndims; i++) {
-      newCart->topo.cart->dims[i] = dims[i];
-      newCart->topo.cart->periodic[i] = periods[i];
+      m_dims[i] = dims[i];
+     m_periodic[i] = periods[i];
       nranks = nranks / dims[i];
       /* FIXME: nranks could be zero (?) */
       nranks = nranks / dims[i];
       /* FIXME: nranks could be zero (?) */
-      newCart->topo.cart->position[i] = rank / nranks;
+      m_position[i] = rank / nranks;
       rank = rank % nranks;
     }
 
       rank = rank % nranks;
     }
 
-    *comm_cart = new simgrid::SMPI::Comm(newGroup, newCart);
+    *comm_cart = new simgrid::SMPI::Comm(newGroup, this);
   } else {
     if (rank == 0) {
   } else {
     if (rank == 0) {
-      newCart = smpi_cart_topo_create(ndims);
-      *comm_cart = new simgrid::SMPI::Comm(new simgrid::SMPI::Group(MPI_COMM_SELF->group()), newCart);
+      *comm_cart = new simgrid::SMPI::Comm(new simgrid::SMPI::Group(MPI_COMM_SELF->group()), this);
     } else {
       *comm_cart = MPI_COMM_NULL;
     }
   }
     } else {
       *comm_cart = MPI_COMM_NULL;
     }
   }
-  return retval;
+  m_comm=*comm_cart;
 }
 
 }
 
-int smpi_mpi_cart_sub(MPI_Comm comm, const int remain_dims[], MPI_Comm *newcomm) {
-  MPI_Topology oldTopo = comm->topo();
-  int oldNDims = oldTopo->topo.cart->ndims;
+Cart* Cart::sub(const int remain_dims[], MPI_Comm *newcomm) {
+  int oldNDims = m_ndims;
   int j = 0;
   int *newDims = nullptr;
   int *newPeriodic = nullptr;
 
   if (remain_dims == nullptr && oldNDims != 0) {
   int j = 0;
   int *newDims = nullptr;
   int *newPeriodic = nullptr;
 
   if (remain_dims == nullptr && oldNDims != 0) {
-    return MPI_ERR_ARG;
+    return nullptr;
   }
   int newNDims = 0;
   for (int i = 0 ; i < oldNDims ; i++) {
   }
   int newNDims = 0;
   for (int i = 0 ; i < oldNDims ; i++) {
@@ -190,40 +122,37 @@ int smpi_mpi_cart_sub(MPI_Comm comm, const int remain_dims[], MPI_Comm *newcomm)
     // that should not segfault
     for (int i = 0 ; j < newNDims ; i++) {
       if(remain_dims[i]) {
     // that should not segfault
     for (int i = 0 ; j < newNDims ; i++) {
       if(remain_dims[i]) {
-        newDims[j] = oldTopo->topo.cart->dims[i];
-        newPeriodic[j] = oldTopo->topo.cart->periodic[i];
+        newDims[j] =m_dims[i];
+        newPeriodic[j] =m_periodic[i];
         j++;
       }
     }
   }
         j++;
       }
     }
   }
-  return smpi_mpi_cart_create(comm, newNDims, newDims, newPeriodic, 0, newcomm);
+  return new Cart(m_comm, newNDims, newDims, newPeriodic, 0, newcomm);
 }
 
 }
 
-int smpi_mpi_cart_coords(MPI_Comm comm, int rank, int maxdims, int coords[]) {
-  MPI_Topology topo = comm->topo();
-  int nnodes = topo->topo.cart->nnodes;
-  for (int i = 0; i< topo->topo.cart->ndims; i++ ) {
-    nnodes    = nnodes / topo->topo.cart->dims[i];
+int Cart::coords(int rank, int maxdims, int coords[]) {
+  int nnodes = m_nnodes;
+  for (int i = 0; i< m_ndims; i++ ) {
+    nnodes    = nnodes /m_dims[i];
     coords[i] = rank / nnodes;
     rank      = rank % nnodes;
   }
   return MPI_SUCCESS;
 }
 
     coords[i] = rank / nnodes;
     rank      = rank % nnodes;
   }
   return MPI_SUCCESS;
 }
 
-int smpi_mpi_cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
-  MPI_Topology topo = comm->topo();
-  int ndims=topo->topo.cart->ndims < maxdims ? topo->topo.cart->ndims : maxdims;
+int Cart::get(int maxdims, int* dims, int* periods, int* coords) {
+  int ndims=m_ndims < maxdims ?m_ndims : maxdims;
   for(int i = 0 ; i < ndims ; i++) {
   for(int i = 0 ; i < ndims ; i++) {
-    dims[i] = topo->topo.cart->dims[i];
-    periods[i] = topo->topo.cart->periodic[i];
-    coords[i] = topo->topo.cart->position[i];
+    dims[i] =m_dims[i];
+    periods[i] =m_periodic[i];
+    coords[i] =m_position[i];
   }
   return MPI_SUCCESS;
 }
 
   }
   return MPI_SUCCESS;
 }
 
-int smpi_mpi_cart_rank(MPI_Comm comm, int* coords, int* rank) {
-  MPI_Topology topo = comm->topo();
-  int ndims = topo->topo.cart->ndims;
+int Cart::rank(int* coords, int* rank) {
+  int ndims =m_ndims;
   int coord;
   *rank = 0;
   int multiplier = 1;
   int coord;
   *rank = 0;
   int multiplier = 1;
@@ -234,19 +163,19 @@ int smpi_mpi_cart_rank(MPI_Comm comm, int* coords, int* rank) {
     /* The user can give us whatever coordinates he wants. If one of them is out of range, either this dimension is
      * periodic, and we consider the equivalent coordinate inside the bounds, or it's not and then it's an error
      */
     /* The user can give us whatever coordinates he wants. If one of them is out of range, either this dimension is
      * periodic, and we consider the equivalent coordinate inside the bounds, or it's not and then it's an error
      */
-    if (coord >= topo->topo.cart->dims[i]) {
-      if ( topo->topo.cart->periodic[i] ) {
-        coord = coord % topo->topo.cart->dims[i];
+    if (coord >=m_dims[i]) {
+      if (m_periodic[i] ) {
+        coord = coord %m_dims[i];
       } else {
         // Should I do that ?
         *rank = -1;
         return MPI_ERR_ARG;
       }
     } else if (coord <  0) {
       } else {
         // Should I do that ?
         *rank = -1;
         return MPI_ERR_ARG;
       }
     } else if (coord <  0) {
-      if(topo->topo.cart->periodic[i]) {
-        coord = coord % topo->topo.cart->dims[i];
+      if(m_periodic[i]) {
+        coord = coord %m_dims[i];
         if (coord)
         if (coord)
-          coord = topo->topo.cart->dims[i] + coord;
+          coord =m_dims[i] + coord;
       } else {
         *rank = -1;
         return MPI_ERR_ARG;
       } else {
         *rank = -1;
         return MPI_ERR_ARG;
@@ -254,56 +183,54 @@ int smpi_mpi_cart_rank(MPI_Comm comm, int* coords, int* rank) {
     }
 
     *rank += multiplier * coord;
     }
 
     *rank += multiplier * coord;
-    multiplier *= topo->topo.cart->dims[i];
+    multiplier *=m_dims[i];
   }
   return MPI_SUCCESS;
 }
 
   }
   return MPI_SUCCESS;
 }
 
-int smpi_mpi_cart_shift(MPI_Comm comm, int direction, int disp, int *rank_source, int *rank_dest) {
-  MPI_Topology topo = comm->topo();
-  int position[topo->topo.cart->ndims];
+int Cart::shift(int direction, int disp, int *rank_source, int *rank_dest) {
 
 
-  if(topo->topo.cart->ndims == 0) {
+  int position[m_ndims];
+
+  if(m_ndims == 0) {
     return MPI_ERR_ARG;
   }
     return MPI_ERR_ARG;
   }
-  if (topo->topo.cart->ndims < direction) {
+  if (m_ndims < direction) {
     return MPI_ERR_DIMS;
   }
 
     return MPI_ERR_DIMS;
   }
 
-  smpi_mpi_cart_coords(comm, comm->rank(), topo->topo.cart->ndims, position);
+  this->coords(m_comm->rank(),m_ndims, position);
   position[direction] += disp;
 
   if(position[direction] < 0 ||
   position[direction] += disp;
 
   if(position[direction] < 0 ||
-      position[direction] >= topo->topo.cart->dims[direction]) {
-    if(topo->topo.cart->periodic[direction]) {
-      position[direction] %= topo->topo.cart->dims[direction];
-      smpi_mpi_cart_rank(comm, position, rank_dest);
+      position[direction] >=m_dims[direction]) {
+    if(m_periodic[direction]) {
+      position[direction] %=m_dims[direction];
+      this->rank(position, rank_dest);
     } else {
       *rank_dest = MPI_PROC_NULL;
     }
   } else {
     } else {
       *rank_dest = MPI_PROC_NULL;
     }
   } else {
-    smpi_mpi_cart_rank(comm, position, rank_dest);
+    this->rank(position, rank_dest);
   }
 
   }
 
-  position[direction] =  topo->topo.cart->position[direction] - disp;
-  if(position[direction] < 0 || position[direction] >= topo->topo.cart->dims[direction]) {
-    if(topo->topo.cart->periodic[direction]) {
-      position[direction] %= topo->topo.cart->dims[direction];
-      smpi_mpi_cart_rank(comm, position, rank_source);
+  position[direction] = m_position[direction] - disp;
+  if(position[direction] < 0 || position[direction] >=m_dims[direction]) {
+    if(m_periodic[direction]) {
+      position[direction] %=m_dims[direction];
+      this->rank(position, rank_source);
     } else {
       *rank_source = MPI_PROC_NULL;
     }
   } else {
     } else {
       *rank_source = MPI_PROC_NULL;
     }
   } else {
-    smpi_mpi_cart_rank(comm, position, rank_source);
+    this->rank(position, rank_source);
   }
 
   return MPI_SUCCESS;
 }
 
   }
 
   return MPI_SUCCESS;
 }
 
-int smpi_mpi_cartdim_get(MPI_Comm comm, int *ndims) {
-  MPI_Topology topo = comm->topo();
-
-  *ndims = topo->topo.cart->ndims;
+int Cart::dim_get(int *ndims) {
+  *ndims =m_ndims;
   return MPI_SUCCESS;
 }
 
   return MPI_SUCCESS;
 }
 
@@ -330,14 +257,11 @@ int smpi_mpi_cartdim_get(MPI_Comm comm, int *ndims) {
  * $HEADER$
  */
 
  * $HEADER$
  */
 
-/* static functions */
-static int assignnodes(int ndim, int nfactor, int *pfacts,int **pdims);
-static int getfactors(int num, int *nfators, int **factors);
 
 /*
  * This is a utility function, no need to have anything in the lower layer for this at all
  */
 
 /*
  * This is a utility function, no need to have anything in the lower layer for this at all
  */
-int smpi_mpi_dims_create(int nnodes, int ndims, int dims[])
+int Dims_create(int nnodes, int ndims, int dims[])
 {
   /* Get # of free-to-be-assigned processes and # of free dimensions */
   int freeprocs = nnodes;
 {
   /* Get # of free-to-be-assigned processes and # of free dimensions */
   int freeprocs = nnodes;
@@ -401,6 +325,9 @@ int smpi_mpi_dims_create(int nnodes, int ndims, int dims[])
   return MPI_SUCCESS;
 }
 
   return MPI_SUCCESS;
 }
 
+}
+}
+
 /*
  *  assignnodes
  *
 /*
  *  assignnodes
  *
@@ -506,3 +433,4 @@ static int getfactors(int num, int *nfactors, int **factors) {
   (*nfactors) = i;
   return MPI_SUCCESS;
 }
   (*nfactors) = i;
   return MPI_SUCCESS;
 }
+
diff --git a/src/smpi/smpi_topo.hpp b/src/smpi/smpi_topo.hpp
new file mode 100644 (file)
index 0000000..6fadde1
--- /dev/null
@@ -0,0 +1,76 @@
+/* Copyright (c) 2010-2015. The SimGrid Team.
+ * All rights reserved.                                                     */
+
+/* This program is free software; you can redistribute it and/or modify it
+ * under the terms of the license (GNU LGPL) which comes with this package. */
+
+#ifndef SMPI_TOPO_HPP_INCLUDED
+#define SMPI_TOPO_HPP_INCLUDED
+
+#include "private.h"
+
+
+namespace simgrid{
+namespace SMPI{
+
+class Topo {
+  protected:
+  MPI_Comm m_comm;
+};
+
+
+class Cart: public Topo {
+  private:
+    int m_nnodes;
+    int m_ndims;
+    int *m_dims;
+    int *m_periodic;
+    int *m_position;
+  public:
+    Cart(int ndims);
+    ~Cart();
+    Cart(MPI_Comm comm_old, int ndims, int dims[], int periods[], int reorder, MPI_Comm *comm_cart);
+    Cart* sub(const int remain_dims[], MPI_Comm *newcomm) ;
+    int coords(int rank, int maxdims, int coords[]) ;
+    int get(int maxdims, int* dims, int* periods, int* coords);
+    int rank(int* coords, int* rank);
+    int shift(int direction, int disp, int *rank_source, int *rank_dest) ;
+    int dim_get(int *ndims);
+};
+
+
+class Graph: public Topo {
+  private:
+    int m_nnodes;
+    int m_nedges;
+    int *m_index;
+    int *m_edges;
+  public:
+    Graph();
+    ~Graph();
+};
+
+class Dist_Graph: public Topo {
+  private:
+    int m_indegree;
+    int *m_in;
+    int *m_in_weights;
+    int m_outdegree;
+    int *m_out;
+    int *m_out_weights;
+    int m_is_weighted;
+  public:
+    Dist_Graph();
+    ~Dist_Graph();
+};
+
+/*
+ * This is a utility function, no need to have anything in the lower layer for this at all
+ */
+extern int Dims_create(int nnodes, int ndims, int dims[]);
+
+}
+}
+
+
+#endif
index 2a28a03..02e408d 100644 (file)
@@ -224,6 +224,7 @@ set(SMPI_SRC
   src/smpi/smpi_replay.cpp
   src/smpi/smpi_rma.cpp
   src/smpi/smpi_topo.cpp
   src/smpi/smpi_replay.cpp
   src/smpi/smpi_rma.cpp
   src/smpi/smpi_topo.cpp
+  src/smpi/smpi_topo.hpp
   src/smpi/smpi_utils.cpp
   src/smpi/smpi_f77.cpp
   
   src/smpi/smpi_utils.cpp
   src/smpi/smpi_f77.cpp