Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
1e26cc249bd482c83d95f9de5369e51730747a36
[simgrid.git] / src / smpi / 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 "private.h"
9 #include <stdio.h>
10 #include <xbt.h>
11 #include <xbt/replay.h>
12
13 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_replay,smpi,"Trace Replay with SMPI");
14
15 int communicator_size = 0;
16
17 typedef struct {
18   xbt_dynar_t isends; /* of MPI_Request */
19   xbt_dynar_t irecvs; /* of MPI_Request */
20 } s_smpi_replay_globals_t, *smpi_replay_globals_t;
21
22 /* Helper function */
23 static double parse_double(const char *string)
24 {
25   double value;
26   char *endptr;
27   value = strtod(string, &endptr);
28   if (*endptr != '\0')
29     THROWF(unknown_error, 0, "%s is not a double", string);
30   return value;
31 }
32
33 static void action_init(const char *const *action)
34 {
35   XBT_DEBUG("Initialize the counters");
36   smpi_replay_globals_t globals =  xbt_new(s_smpi_replay_globals_t, 1);
37   globals->isends = xbt_dynar_new(sizeof(MPI_Request),NULL);
38   globals->irecvs = xbt_dynar_new(sizeof(MPI_Request),NULL);
39
40   smpi_process_set_user_data((void*) globals);
41
42   /* be sure that everyone has initiated the counters before proceeding */
43   smpi_mpi_barrier(MPI_COMM_WORLD);
44   /* start a simulated timer */
45   smpi_process_simulated_start();
46 }
47
48 static void action_finalize(const char *const *action)
49 {
50   double sim_time= 1.;
51   smpi_replay_globals_t globals =
52       (smpi_replay_globals_t) smpi_process_get_user_data();
53
54   if (globals){
55     XBT_DEBUG("There are %lu isends and %lu irecvs in the dynars",
56          xbt_dynar_length(globals->isends),xbt_dynar_length(globals->irecvs));
57     xbt_dynar_free_container(&(globals->isends));
58     xbt_dynar_free_container(&(globals->irecvs));
59   }
60   free(globals);
61   /* end the simulated timer */
62   sim_time = smpi_process_simulated_elapsed();
63   if (!smpi_process_index())
64     XBT_INFO("Simulation time %g", sim_time);
65
66 }
67
68 static void action_comm_size(const char *const *action)
69 {
70   double clock = smpi_process_simulated_elapsed();
71
72   communicator_size = parse_double(action[2]);
73
74   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
75            smpi_process_simulated_elapsed()-clock);
76 }
77
78
79 static void action_compute(const char *const *action)
80 {
81   double clock = smpi_process_simulated_elapsed();
82   smpi_execute_flops(parse_double(action[2]));
83
84   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
85            smpi_process_simulated_elapsed()-clock);
86 }
87
88 static void action_send(const char *const *action)
89 {
90   int to = atoi(action[2]);
91   double size=parse_double(action[3]);
92   double clock = smpi_process_simulated_elapsed();
93
94   smpi_mpi_send(NULL, size, MPI_BYTE, to , 0, MPI_COMM_WORLD);
95   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
96            smpi_process_simulated_elapsed()-clock);
97 }
98
99 static void action_Isend(const char *const *action)
100 {
101   int to = atoi(action[2]);
102   double size=parse_double(action[3]);
103   double clock = smpi_process_simulated_elapsed();
104   smpi_replay_globals_t globals =
105      (smpi_replay_globals_t) smpi_process_get_user_data();
106
107   MPI_Request request = smpi_mpi_isend(NULL, size, MPI_BYTE, to, 0,
108                                        MPI_COMM_WORLD);
109
110   xbt_dynar_push(globals->isends,&request);
111
112   //TODO do the asynchronous cleanup
113   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
114            smpi_process_simulated_elapsed()-clock);
115 }
116
117 static void action_recv(const char *const *action) {
118   int from = atoi(action[2]);
119   double size=parse_double(action[3]);
120   double clock = smpi_process_simulated_elapsed();
121   MPI_Status status;
122
123   smpi_mpi_recv(NULL, size, MPI_BYTE, from, 0, MPI_COMM_WORLD, &status);
124   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
125            smpi_process_simulated_elapsed()-clock);
126 }
127
128 static void action_Irecv(const char *const *action)
129 {
130   int from = atoi(action[2]);
131   double size=parse_double(action[3]);
132   double clock = smpi_process_simulated_elapsed();
133   MPI_Request request;
134   smpi_replay_globals_t globals =
135      (smpi_replay_globals_t) smpi_process_get_user_data();
136
137   XBT_DEBUG("Asynchronous receive of %.0f from rank%d (%s)",size, from,
138             action[2]);
139
140   request = smpi_mpi_irecv(NULL, size, MPI_BYTE, from, 0, MPI_COMM_WORLD);
141   xbt_dynar_push(globals->irecvs,&request);
142
143   //TODO do the asynchronous cleanup
144   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
145            smpi_process_simulated_elapsed()-clock);
146 }
147
148 static void action_wait(const char *const *action){
149   double clock = smpi_process_simulated_elapsed();
150   MPI_Request request;
151   MPI_Status status;
152   smpi_replay_globals_t globals =
153       (smpi_replay_globals_t) smpi_process_get_user_data();
154
155   xbt_assert(xbt_dynar_length(globals->irecvs),
156       "action wait not preceded by any irecv: %s",
157       xbt_str_join_array(action," "));
158   request = xbt_dynar_pop_as(globals->irecvs,MPI_Request);
159   smpi_mpi_wait(&request, &status);
160   XBT_DEBUG("MPI_Wait Status is : (source=%d tag=%d error=%d count=%d)",
161             status.MPI_SOURCE, status.MPI_TAG, status.MPI_ERROR, status.count);
162   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
163            smpi_process_simulated_elapsed()-clock);
164 }
165
166 static void action_barrier(const char *const *action){
167   double clock = smpi_process_simulated_elapsed();
168   smpi_mpi_barrier(MPI_COMM_WORLD);
169   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
170            smpi_process_simulated_elapsed()-clock);
171 }
172
173 static void action_bcast(const char *const *action)
174 {
175   double size = parse_double(action[2]);
176   double clock = smpi_process_simulated_elapsed();
177
178   smpi_mpi_bcast(NULL, size, MPI_BYTE, 0, MPI_COMM_WORLD);
179
180   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
181            smpi_process_simulated_elapsed()-clock);
182 }
183
184 static void action_reduce(const char *const *action)
185 {
186   double size = parse_double(action[2]);
187   double clock = smpi_process_simulated_elapsed();
188   smpi_mpi_reduce(NULL, NULL, size, MPI_BYTE, MPI_OP_NULL, 0, MPI_COMM_WORLD);
189
190   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
191            smpi_process_simulated_elapsed()-clock);
192 }
193
194 static void action_allReduce(const char *const *action) {
195   double comm_size = parse_double(action[2]);
196   double comp_size = parse_double(action[3]);
197   double clock = smpi_process_simulated_elapsed();
198   smpi_mpi_reduce(NULL, NULL, comm_size, MPI_BYTE, MPI_OP_NULL, 0, MPI_COMM_WORLD);
199   smpi_execute_flops(comp_size);
200   smpi_mpi_bcast(NULL, comm_size, MPI_BYTE, 0, MPI_COMM_WORLD);
201   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
202            smpi_process_simulated_elapsed()-clock);
203 }
204
205 void smpi_replay_init(int *argc, char***argv){
206   PMPI_Init(argc, argv);
207   _xbt_replay_action_init();
208
209   xbt_replay_action_register("init",     action_init);
210   xbt_replay_action_register("finalize", action_finalize);
211   xbt_replay_action_register("comm_size",action_comm_size);
212   xbt_replay_action_register("send",     action_send);
213   xbt_replay_action_register("Isend",    action_Isend);
214   xbt_replay_action_register("recv",     action_recv);
215   xbt_replay_action_register("Irecv",    action_Irecv);
216   xbt_replay_action_register("wait",     action_wait);
217   xbt_replay_action_register("barrier",  action_barrier);
218   xbt_replay_action_register("bcast",    action_bcast);
219   xbt_replay_action_register("reduce",   action_reduce);
220   xbt_replay_action_register("allReduce",action_allReduce);
221   xbt_replay_action_register("compute",  action_compute);
222
223   xbt_replay_action_runner(*argc, *argv);
224 }
225
226 int smpi_replay_finalize(){
227   return PMPI_Finalize();
228 }