Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
don't need to synchronize processes anymore
[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   /* start a simulated timer */
43   smpi_process_simulated_start();
44 }
45
46 static void action_finalize(const char *const *action)
47 {
48   double sim_time= 1.;
49   smpi_replay_globals_t globals =
50       (smpi_replay_globals_t) smpi_process_get_user_data();
51
52   if (globals){
53     XBT_DEBUG("There are %lu isends and %lu irecvs in the dynars",
54          xbt_dynar_length(globals->isends),xbt_dynar_length(globals->irecvs));
55     xbt_dynar_free_container(&(globals->isends));
56     xbt_dynar_free_container(&(globals->irecvs));
57   }
58   free(globals);
59   /* end the simulated timer */
60   sim_time = smpi_process_simulated_elapsed();
61   if (!smpi_process_index())
62     XBT_INFO("Simulation time %g", sim_time);
63
64 }
65
66 static void action_comm_size(const char *const *action)
67 {
68   double clock = smpi_process_simulated_elapsed();
69
70   communicator_size = parse_double(action[2]);
71
72   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
73            smpi_process_simulated_elapsed()-clock);
74 }
75
76
77 static void action_compute(const char *const *action)
78 {
79   double clock = smpi_process_simulated_elapsed();
80   smpi_execute_flops(parse_double(action[2]));
81
82   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
83            smpi_process_simulated_elapsed()-clock);
84 }
85
86 static void action_send(const char *const *action)
87 {
88   int to = atoi(action[2]);
89   double size=parse_double(action[3]);
90   double clock = smpi_process_simulated_elapsed();
91
92   smpi_mpi_send(NULL, size, MPI_BYTE, to , 0, MPI_COMM_WORLD);
93   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
94            smpi_process_simulated_elapsed()-clock);
95 }
96
97 static void action_Isend(const char *const *action)
98 {
99   int to = atoi(action[2]);
100   double size=parse_double(action[3]);
101   double clock = smpi_process_simulated_elapsed();
102   smpi_replay_globals_t globals =
103      (smpi_replay_globals_t) smpi_process_get_user_data();
104
105   MPI_Request request = smpi_mpi_isend(NULL, size, MPI_BYTE, to, 0,
106                                        MPI_COMM_WORLD);
107
108   xbt_dynar_push(globals->isends,&request);
109
110   //TODO do the asynchronous cleanup
111   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
112            smpi_process_simulated_elapsed()-clock);
113 }
114
115 static void action_recv(const char *const *action) {
116   int from = atoi(action[2]);
117   double size=parse_double(action[3]);
118   double clock = smpi_process_simulated_elapsed();
119   MPI_Status status;
120
121   smpi_mpi_recv(NULL, size, MPI_BYTE, from, 0, MPI_COMM_WORLD, &status);
122   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
123            smpi_process_simulated_elapsed()-clock);
124 }
125
126 static void action_Irecv(const char *const *action)
127 {
128   int from = atoi(action[2]);
129   double size=parse_double(action[3]);
130   double clock = smpi_process_simulated_elapsed();
131   MPI_Request request;
132   smpi_replay_globals_t globals =
133      (smpi_replay_globals_t) smpi_process_get_user_data();
134
135   XBT_DEBUG("Asynchronous receive of %.0f from rank%d (%s)",size, from,
136             action[2]);
137
138   request = smpi_mpi_irecv(NULL, size, MPI_BYTE, from, 0, MPI_COMM_WORLD);
139   xbt_dynar_push(globals->irecvs,&request);
140
141   //TODO do the asynchronous cleanup
142   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
143            smpi_process_simulated_elapsed()-clock);
144 }
145
146 static void action_wait(const char *const *action){
147   double clock = smpi_process_simulated_elapsed();
148   MPI_Request request;
149   MPI_Status status;
150   smpi_replay_globals_t globals =
151       (smpi_replay_globals_t) smpi_process_get_user_data();
152
153   xbt_assert(xbt_dynar_length(globals->irecvs),
154       "action wait not preceded by any irecv: %s",
155       xbt_str_join_array(action," "));
156   request = xbt_dynar_pop_as(globals->irecvs,MPI_Request);
157   smpi_mpi_wait(&request, &status);
158   XBT_DEBUG("MPI_Wait Status is : (source=%d tag=%d error=%d count=%d)",
159             status.MPI_SOURCE, status.MPI_TAG, status.MPI_ERROR, status.count);
160   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
161            smpi_process_simulated_elapsed()-clock);
162 }
163
164 static void action_barrier(const char *const *action){
165   double clock = smpi_process_simulated_elapsed();
166   smpi_mpi_barrier(MPI_COMM_WORLD);
167   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
168            smpi_process_simulated_elapsed()-clock);
169 }
170
171 static void action_bcast(const char *const *action)
172 {
173   double size = parse_double(action[2]);
174   double clock = smpi_process_simulated_elapsed();
175
176   smpi_mpi_bcast(NULL, size, MPI_BYTE, 0, MPI_COMM_WORLD);
177
178   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
179            smpi_process_simulated_elapsed()-clock);
180 }
181
182 static void action_reduce(const char *const *action)
183 {
184   double size = parse_double(action[2]);
185   double clock = smpi_process_simulated_elapsed();
186   smpi_mpi_reduce(NULL, NULL, size, MPI_BYTE, MPI_OP_NULL, 0, MPI_COMM_WORLD);
187
188   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
189            smpi_process_simulated_elapsed()-clock);
190 }
191
192 static void action_allReduce(const char *const *action) {
193   double comm_size = parse_double(action[2]);
194   double comp_size = parse_double(action[3]);
195   double clock = smpi_process_simulated_elapsed();
196   smpi_mpi_reduce(NULL, NULL, comm_size, MPI_BYTE, MPI_OP_NULL, 0, MPI_COMM_WORLD);
197   smpi_execute_flops(comp_size);
198   smpi_mpi_bcast(NULL, comm_size, MPI_BYTE, 0, MPI_COMM_WORLD);
199   XBT_VERB("%s %f", xbt_str_join_array(action, " "),
200            smpi_process_simulated_elapsed()-clock);
201 }
202
203 void smpi_replay_init(int *argc, char***argv){
204   PMPI_Init(argc, argv);
205   _xbt_replay_action_init();
206
207   xbt_replay_action_register("init",     action_init);
208   xbt_replay_action_register("finalize", action_finalize);
209   xbt_replay_action_register("comm_size",action_comm_size);
210   xbt_replay_action_register("send",     action_send);
211   xbt_replay_action_register("Isend",    action_Isend);
212   xbt_replay_action_register("recv",     action_recv);
213   xbt_replay_action_register("Irecv",    action_Irecv);
214   xbt_replay_action_register("wait",     action_wait);
215   xbt_replay_action_register("barrier",  action_barrier);
216   xbt_replay_action_register("bcast",    action_bcast);
217   xbt_replay_action_register("reduce",   action_reduce);
218   xbt_replay_action_register("allReduce",action_allReduce);
219   xbt_replay_action_register("compute",  action_compute);
220
221   xbt_replay_action_runner(*argc, *argv);
222 }
223
224 int smpi_replay_finalize(){
225   return PMPI_Finalize();
226 }