Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Do the same as before to replay traces with SMPI, but one layer lower.
[simgrid.git] / examples / smpi / replay.c
1 /* Copyright (c) 2009, 2010, 2011, 2012. The SimGrid Team.
2  * All rights reserved.                                                     */
3
4 /* This program is free software; you can redistribute it and/or modify it
5  * under the terms of the license (GNU LGPL) which comes with this package. */
6
7
8 #include <stdio.h>
9 #include <mpi.h>
10 #include <xbt.h>
11 #include <xbt/replay.h>
12
13 XBT_LOG_NEW_DEFAULT_CATEGORY(smpi_replay,"Trace Replay with SMPI");
14
15 /* Helper function */
16 static double parse_double(const char *string)
17 {
18   double value;
19   char *endptr;
20   value = strtod(string, &endptr);
21   if (*endptr != '\0')
22     THROWF(unknown_error, 0, "%s is not a double", string);
23   return value;
24 }
25
26 static void action_compute(const char *const *action)
27 {
28   XBT_DEBUG("Start to compute %.0f flops", parse_double(action[2]));
29   smpi_sample_flops(parse_double(action[2]));
30 }
31
32 static void action_send(const char *const *action)
33 {
34   int to = atoi(action[2]);
35   const char *size_str = action[3];
36   double size=parse_double(size_str);
37   char *buffer = (char*) malloc (size);
38
39   XBT_DEBUG("send %.0f bytes to rank%d (%s)",size, to, action[2]);
40
41   MPI_Send(buffer, size, MPI_BYTE, to, 0, MPI_COMM_WORLD);
42
43   free(buffer);
44 }
45
46 static void action_Isend(const char *const *action)
47 {
48   int to = atoi(action[2]);
49   const char *size_str = action[3];
50   double size=parse_double(size_str);
51   char *buffer = (char*) malloc (size);
52
53   MPI_Request request;
54
55   MPI_Isend(buffer, size, MPI_BYTE, to, 0, MPI_COMM_WORLD, &request);
56
57   free(buffer);
58 }
59
60 static void action_recv(const char *const *action)
61 {
62   int from = atoi(action[2]);
63   XBT_DEBUG("receive from rank%d (%s)",from, action[2]);
64   char *buffer = (char*) malloc (1e9);
65
66   MPI_Status status;
67
68   MPI_Recv(buffer, 1e9, MPI_BYTE, from, 0, MPI_COMM_WORLD, &status);
69
70   free(buffer);
71 }
72
73 static void action_Irecv(const char *const *action)
74 {
75   int from = atoi(action[2]);
76   char *buffer = (char*) malloc (1e9);
77   MPI_Request request;
78
79   XBT_DEBUG("Asynchronous receive from rank%d (%s)",from, action[2]);
80
81   MPI_Irecv(buffer, 1e9, MPI_BYTE, from, 0, MPI_COMM_WORLD, &request);
82
83   free(buffer);
84 }
85
86 static void action_wait(const char *const *action)
87 {
88   MPI_Request request;
89   MPI_Status status;
90
91   MPI_Wait(&request, &status);
92 }
93
94
95
96 int main(int argc, char *argv[])
97 {
98   int i;
99
100   MPI_Init(&argc, &argv);
101   MPI_Comm_rank(MPI_COMM_WORLD, &i);
102   _xbt_replay_action_init();
103
104 //  xbt_replay_action_register("init",     action_init);
105 //  xbt_replay_action_register("finalize", action_finalize);
106 //  xbt_replay_action_register("comm_size",action_comm_size);
107   xbt_replay_action_register("send",     action_send);
108   xbt_replay_action_register("Isend",    action_Isend);
109   xbt_replay_action_register("recv",     action_recv);
110   xbt_replay_action_register("Irecv",    action_Irecv);
111   xbt_replay_action_register("wait",     action_wait);
112 //  xbt_replay_action_register("barrier",  action_barrier);
113 //  xbt_replay_action_register("bcast",    action_bcast);
114 //  xbt_replay_action_register("reduce",   action_reduce);
115 //  xbt_replay_action_register("allReduce",action_allReduce);
116 //  xbt_replay_action_register("sleep",    action_sleep);
117   xbt_replay_action_register("compute",  action_compute);
118
119   xbt_replay_action_runner(argc, argv);
120   /* Actually do the simulation using MSG_action_trace_run */
121   smpi_action_trace_run(argv[1]);  // it's ok to pass a NULL argument here
122
123   XBT_DEBUG("rank%d: I'm done!",i);
124   MPI_Finalize();
125   return 0;
126 }