Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Increase required width for times in tesh file.
[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 <simix/simix.h>
12 #include <xbt/replay.h>
13
14 XBT_LOG_NEW_DEFAULT_CATEGORY(smpi_replay,"Trace Replay with SMPI");
15
16 /* Helper function */
17 static double parse_double(const char *string)
18 {
19   double value;
20   char *endptr;
21   value = strtod(string, &endptr);
22   if (*endptr != '\0')
23     THROWF(unknown_error, 0, "%s is not a double", string);
24   return value;
25 }
26
27 static void action_compute(const char *const *action)
28 {
29   XBT_DEBUG("Start to compute %.0f flops", parse_double(action[2]));
30   smpi_sample_flops(parse_double(action[2]));
31 }
32
33 static void action_send(const char *const *action)
34 {
35   int to = atoi(action[2]);
36   const char *size_str = action[3];
37   double size=parse_double(size_str);
38   char *buffer = (char*) malloc (size);
39
40   XBT_DEBUG("send %.0f bytes to rank%d (%s)",size, to, action[2]);
41
42   MPI_Send(buffer, size, MPI_BYTE, to, 0, MPI_COMM_WORLD);
43
44   free(buffer);
45 }
46
47 static void action_Isend(const char *const *action)
48 {
49   int to = atoi(action[2]);
50   const char *size_str = action[3];
51   double size=parse_double(size_str);
52   char *buffer = (char*) malloc (size);
53
54   MPI_Request request;
55
56   MPI_Isend(buffer, size, MPI_BYTE, to, 0, MPI_COMM_WORLD, &request);
57
58   free(buffer);
59 }
60
61 static void action_recv(const char *const *action)
62 {
63   int from = atoi(action[2]);
64   XBT_DEBUG("receive from rank%d (%s)",from, action[2]);
65   char *buffer = (char*) malloc (1e10);
66
67   MPI_Status status;
68
69   MPI_Recv(buffer, 1e9, MPI_BYTE, from, 0, MPI_COMM_WORLD, &status);
70
71   free(buffer);
72 }
73
74 int main(int argc, char *argv[])
75 {
76   int i;
77
78   MPI_Init(&argc, &argv);
79   MPI_Comm_rank(MPI_COMM_WORLD, &i);
80   _xbt_replay_action_init();
81
82 //  xbt_replay_action_register("init",     action_init);
83 //  xbt_replay_action_register("finalize", action_finalize);
84 //  xbt_replay_action_register("comm_size",action_comm_size);
85   xbt_replay_action_register("send",     action_send);
86   xbt_replay_action_register("Isend",    action_Isend);
87   xbt_replay_action_register("recv",     action_recv);
88 //  xbt_replay_action_register("Irecv",    action_Irecv);
89 //  xbt_replay_action_register("wait",     action_wait);
90 //  xbt_replay_action_register("barrier",  action_barrier);
91 //  xbt_replay_action_register("bcast",    action_bcast);
92 //  xbt_replay_action_register("reduce",   action_reduce);
93 //  xbt_replay_action_register("allReduce",action_allReduce);
94 //  xbt_replay_action_register("sleep",    action_sleep);
95   xbt_replay_action_register("compute",  action_compute);
96
97   xbt_replay_action_runner(argc, argv);
98   /* Actually do the simulation using MSG_action_trace_run */
99   smpi_action_trace_run(argv[1]);  // it's ok to pass a NULL argument here
100
101   XBT_DEBUG("rank%d: I'm done!",i);
102   MPI_Finalize();
103   return 0;
104 }