Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
dafuk? don't mix seconds and flops!
[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 /* 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_execute_flops(parse_double(action[2]));
30 }
31
32 static void action_send(const char *const *action)
33 {
34   int to = atoi(action[2]);
35   double size=parse_double(action[3]);
36
37   XBT_DEBUG("send %.0f bytes to rank%d (%s)",size, to, action[2]);
38
39   smpi_mpi_send(NULL, size, MPI_BYTE, to , 0, MPI_COMM_WORLD);
40 }
41
42 static void action_Isend(const char *const *action)
43 {
44   int to = atoi(action[2]);
45   double size=parse_double(action[3]);
46
47   MPI_Request request = smpi_mpi_isend(NULL, size, MPI_BYTE, to, 0,
48                                        MPI_COMM_WORLD);
49
50   //TODO do something with request
51   request = NULL;
52 }
53
54 static void action_recv(const char *const *action) {
55   int from = atoi(action[2]);
56   XBT_DEBUG("receive from rank%d (%s)",from, action[2]);
57
58   MPI_Status status;
59   //TODO find a way to get the message size
60   smpi_mpi_recv(NULL, 1e9, MPI_BYTE, from, 0, MPI_COMM_WORLD, &status);
61
62 }
63
64 static void action_Irecv(const char *const *action)
65 {
66   int from = atoi(action[2]);
67   MPI_Request request;
68
69   XBT_DEBUG("Asynchronous receive from rank%d (%s)",from, action[2]);
70
71   //TODO find a way to get the message size
72   request = smpi_mpi_irecv(NULL, 1e9, MPI_BYTE, from, 0, MPI_COMM_WORLD);
73
74   //TODO do something with request
75   request = NULL;
76 }
77
78 static void action_wait(const char *const *action){
79   MPI_Request request;
80   MPI_Status status;
81
82   smpi_mpi_wait(&request, &status);
83 }
84
85 static void action_barrier(const char *const *action){
86   smpi_mpi_barrier(MPI_COMM_WORLD);
87 }
88
89 void smpi_replay_init(int *argc, char***argv){
90   PMPI_Init(argc, argv);
91   _xbt_replay_action_init();
92
93 //  xbt_replay_action_register("init",     action_init);
94 //  xbt_replay_action_register("finalize", action_finalize);
95 //  xbt_replay_action_register("comm_size",action_comm_size);
96   xbt_replay_action_register("send",     action_send);
97   xbt_replay_action_register("Isend",    action_Isend);
98   xbt_replay_action_register("recv",     action_recv);
99   xbt_replay_action_register("Irecv",    action_Irecv);
100   xbt_replay_action_register("wait",     action_wait);
101   xbt_replay_action_register("barrier",  action_barrier);
102 //  xbt_replay_action_register("bcast",    action_bcast);
103 //  xbt_replay_action_register("reduce",   action_reduce);
104 //  xbt_replay_action_register("allReduce",action_allReduce);
105 //  xbt_replay_action_register("sleep",    action_sleep);
106   xbt_replay_action_register("compute",  action_compute);
107
108   xbt_replay_action_runner(*argc, *argv);
109 }
110
111 int smpi_replay_finalize(){
112   return PMPI_Finalize();
113 }