Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add a Q&D implementation of trace replay with SMPI. Just to get the
authorsuter <frederic.suter@cc.in2p3.fr>
Mon, 9 Jul 2012 12:33:19 +0000 (14:33 +0200)
committersuter <frederic.suter@cc.in2p3.fr>
Mon, 9 Jul 2012 12:33:44 +0000 (14:33 +0200)
flavor of it, some issues have to be solved and code has to be polished
(we don't know, marion might read it ...)

examples/smpi/CMakeLists.txt
examples/smpi/actions.txt [new file with mode: 0644]
examples/smpi/replay.c [new file with mode: 0644]

index 1c950f3..b98c20f 100644 (file)
@@ -28,6 +28,7 @@ if(enable_smpi)
   add_executable(ttest01 ttest01.c)
   add_executable(mc_bugged1 mc_bugged1.c)
   add_executable(mc_bugged2 mc_bugged2.c)
+  add_executable(replay replay.c)
 
   target_link_libraries(alltoall2 m simgrid smpi )
   target_link_libraries(alltoall_basic m simgrid smpi )
@@ -50,6 +51,7 @@ if(enable_smpi)
   target_link_libraries(ttest01 m simgrid smpi )
   target_link_libraries(mc_bugged1 m simgrid smpi )
   target_link_libraries(mc_bugged2 m simgrid smpi )
+  target_link_libraries(replay m simgrid smpi )
 
   set_target_properties(smpi_sendrecv PROPERTIES RENAME sendrecv)
 endif(enable_smpi)
@@ -76,6 +78,7 @@ set(examples_src
   ${CMAKE_CURRENT_SOURCE_DIR}/alltoall_basic.c
   ${CMAKE_CURRENT_SOURCE_DIR}/sendrecv.c
   ${CMAKE_CURRENT_SOURCE_DIR}/mc_bugged1.c
+  ${CMAKE_CURRENT_SOURCE_DIR}/replay.c
   ${CMAKE_CURRENT_SOURCE_DIR}/reduce.c
   ${CMAKE_CURRENT_SOURCE_DIR}/mvmul.c
   ${CMAKE_CURRENT_SOURCE_DIR}/compute2.c
diff --git a/examples/smpi/actions.txt b/examples/smpi/actions.txt
new file mode 100644 (file)
index 0000000..5456907
--- /dev/null
@@ -0,0 +1,5 @@
+0 send 1 1e6
+1 recv 0
+1 compute 1e9
+1 Isend 0 1e6
+0 recv 1
diff --git a/examples/smpi/replay.c b/examples/smpi/replay.c
new file mode 100644 (file)
index 0000000..d6b7b43
--- /dev/null
@@ -0,0 +1,104 @@
+/* Copyright (c) 2009, 2010, 2011, 2012. 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. */
+
+
+#include <stdio.h>
+#include <mpi.h>
+#include <xbt.h>
+#include <simix/simix.h>
+#include <xbt/replay.h>
+
+XBT_LOG_NEW_DEFAULT_CATEGORY(smpi_replay,"Trace Replay with SMPI");
+
+/* Helper function */
+static double parse_double(const char *string)
+{
+  double value;
+  char *endptr;
+  value = strtod(string, &endptr);
+  if (*endptr != '\0')
+    THROWF(unknown_error, 0, "%s is not a double", string);
+  return value;
+}
+
+static void action_compute(const char *const *action)
+{
+  XBT_DEBUG("Start to compute %.0f flops", parse_double(action[2]));
+  smpi_sample_flops(parse_double(action[2]));
+}
+
+static void action_send(const char *const *action)
+{
+  int to = atoi(action[2]);
+  const char *size_str = action[3];
+  double size=parse_double(size_str);
+  char *buffer = (char*) malloc (size);
+
+  XBT_DEBUG("send %.0f bytes to rank%d (%s)",size, to, action[2]);
+
+  MPI_Send(buffer, size, MPI_BYTE, to, 0, MPI_COMM_WORLD);
+
+  free(buffer);
+}
+
+static void action_Isend(const char *const *action)
+{
+  int to = atoi(action[2]);
+  const char *size_str = action[3];
+  double size=parse_double(size_str);
+  char *buffer = (char*) malloc (size);
+
+  MPI_Request request;
+
+  MPI_Isend(buffer, size, MPI_BYTE, to, 0, MPI_COMM_WORLD, &request);
+
+  free(buffer);
+}
+
+static void action_recv(const char *const *action)
+{
+  int from = atoi(action[2]);
+  XBT_DEBUG("receive from rank%d (%s)",from, action[2]);
+  char *buffer = (char*) malloc (1e10);
+
+  MPI_Status status;
+
+  MPI_Recv(buffer, 1e9, MPI_BYTE, from, 0, MPI_COMM_WORLD, &status);
+
+  free(buffer);
+}
+
+int main(int argc, char *argv[])
+{
+  int i;
+
+  MPI_Init(&argc, &argv);
+  MPI_Comm_rank(MPI_COMM_WORLD, &i);
+  _xbt_replay_action_init();
+
+//  xbt_replay_action_register("init",     action_init);
+//  xbt_replay_action_register("finalize", action_finalize);
+//  xbt_replay_action_register("comm_size",action_comm_size);
+  xbt_replay_action_register("send",     action_send);
+  xbt_replay_action_register("Isend",    action_Isend);
+  xbt_replay_action_register("recv",     action_recv);
+//  xbt_replay_action_register("Irecv",    action_Irecv);
+//  xbt_replay_action_register("wait",     action_wait);
+//  xbt_replay_action_register("barrier",  action_barrier);
+//  xbt_replay_action_register("bcast",    action_bcast);
+//  xbt_replay_action_register("reduce",   action_reduce);
+//  xbt_replay_action_register("allReduce",action_allReduce);
+//  xbt_replay_action_register("sleep",    action_sleep);
+  xbt_replay_action_register("compute",  action_compute);
+
+  xbt_replay_action_runner(argc, argv);
+  /* Actually do the simulation using MSG_action_trace_run */
+  smpi_action_trace_run(argv[1]);  // it's ok to pass a NULL argument here
+
+  XBT_DEBUG("rank%d: I'm done!",i);
+  MPI_Finalize();
+  return 0;
+}