Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'master' of github.com:mquinson/simgrid
[simgrid.git] / examples / smpi / smpi_msg_masterslave / masterslave_mailbox_smpi.c
1 /* Copyright (c) 2010-2015. 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 #include <stdio.h>
8 #include "simgrid/msg.h"            /* Yeah! If you want to use msg, you need to include simgrid/msg.h */
9 #include "xbt/sysdep.h"         /* calloc, printf */
10 #include "mpi.h"
11 /* Create a log channel to have nice outputs. */
12 #include "xbt/log.h"
13 #include "xbt/asserts.h"
14 #include "smpi/private.h"
15 XBT_LOG_NEW_DEFAULT_CATEGORY(msg_test,
16                              "Messages specific for this msg example");
17
18 int master(int argc, char *argv[]);
19 int slave(int argc, char *argv[]);
20 int master_mpi(int argc, char *argv[]);
21 int alltoall_mpi(int argc, char *argv[]);
22
23 /** Emitter function  */
24 int master(int argc, char *argv[])
25 {
26   long number_of_tasks = atol(argv[1]);
27   double task_comp_size = atof(argv[2]);
28   double task_comm_size = atof(argv[3]);
29   long slaves_count = atol(argv[4]);
30
31   int i;
32
33   XBT_INFO("Got %ld slaves and %ld tasks to process", slaves_count,
34         number_of_tasks);
35
36   for (i = 0; i < number_of_tasks; i++) {
37     char mailbox[256];
38     char sprintf_buffer[256];
39     msg_task_t task = NULL;
40
41     sprintf(mailbox, "slave-%ld", i % slaves_count);
42     sprintf(sprintf_buffer, "Task_%d", i);
43     task =
44         MSG_task_create(sprintf_buffer, task_comp_size, task_comm_size,
45                         NULL);
46     if (number_of_tasks < 10000 || i % 10000 == 0)
47       XBT_INFO("Sending \"%s\" (of %ld) to mailbox \"%s\"", task->name,
48             number_of_tasks, mailbox);
49
50     MSG_task_send(task, mailbox);
51   }
52
53   XBT_INFO
54       ("All tasks have been dispatched. Let's tell everybody the computation is over.");
55   for (i = 0; i < slaves_count; i++) {
56     char mailbox[80];
57
58     sprintf(mailbox, "slave-%ld", i % slaves_count);
59     msg_task_t finalize = MSG_task_create("finalize", 0, 0, 0);
60     MSG_task_send(finalize, mailbox);
61   }
62
63 //  XBT_INFO("Goodbye now!");
64   return 0;
65 }                               /* end_of_master */
66
67
68 int master_mpi(int argc, char *argv[])
69 {
70   MPI_Init(&argc, &argv);
71
72   int rank;
73   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
74   XBT_INFO("here for rank %d", rank);
75   int test[1000]={rank};
76   if(rank==0)
77         MPI_Send(&test, 1000, 
78                  MPI_INT, 1, 1, MPI_COMM_WORLD); 
79   else
80         MPI_Recv(&test, 1000, 
81                  MPI_INT, 0, 1, MPI_COMM_WORLD, MPI_STATUSES_IGNORE); 
82
83   XBT_INFO("After comm %d", rank);
84   MPI_Finalize();
85
86   XBT_INFO("After finalize %d %d", rank, test[0]);
87   return 0;
88 }                               /* end_of_master */
89
90
91
92 int alltoall_mpi(int argc, char *argv[])
93 {
94   MPI_Init(&argc, &argv);
95
96   int rank, size;
97   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
98   MPI_Comm_size(MPI_COMM_WORLD, &size);
99   XBT_INFO("alltoall for rank %d", rank);
100   int* out=malloc(1000*size*sizeof(int));
101   int* in=malloc(1000*size*sizeof(int));
102   MPI_Alltoall(out, 1000, MPI_INT,in, 1000, MPI_INT, MPI_COMM_WORLD); 
103
104   XBT_INFO("after alltoall %d", rank);
105   free(out);
106   free(in);
107   MPI_Finalize();
108   return 0;
109 }                               /* end_of_master */
110
111
112 /** Receiver function  */
113 int slave(int argc, char *argv[])
114 {
115   msg_task_t task = NULL;
116   XBT_ATTRIB_UNUSED int res;
117   int id = -1;
118   char mailbox[80];
119   XBT_ATTRIB_UNUSED int read;
120
121   read = sscanf(argv[1], "%d", &id);
122   xbt_assert(read, "Invalid argument %s\n", argv[1]);
123
124   sprintf(mailbox, "slave-%d", id);
125
126   while (1) {
127     res = MSG_task_receive(&(task), mailbox);
128     xbt_assert(res == MSG_OK, "MSG_task_get failed");
129
130 //  XBT_INFO("Received \"%s\"", MSG_task_get_name(task));
131     if (!strcmp(MSG_task_get_name(task), "finalize")) {
132       MSG_task_destroy(task);
133       break;
134     }
135 //    XBT_INFO("Processing \"%s\"", MSG_task_get_name(task));
136     MSG_task_execute(task);
137 //    XBT_INFO("\"%s\" done", MSG_task_get_name(task));
138     MSG_task_destroy(task);
139     task = NULL;
140   }
141   XBT_INFO("I'm done. See you!");
142
143   return 0;
144 }                               /* end_of_slave */
145
146 /** Main function */
147 int main(int argc, char *argv[])
148 {
149   msg_error_t res;
150   const char *platform_file;
151   const char *application_file;
152
153   MSG_init(&argc, argv);
154
155   if (argc < 3) {
156     printf("Usage: %s platform_file deployment_file\n", argv[0]);
157     printf("example: %s msg_platform.xml msg_deployment.xml\n", argv[0]);
158     exit(1);
159   }
160   platform_file = argv[1];
161   application_file = argv[2];
162
163   {                             /*  Simulation setting */
164     MSG_create_environment(platform_file);
165   }
166   {                             /*   Application deployment */
167     MSG_function_register("master", master);
168     MSG_function_register("slave", slave);
169     // launch two MPI applications as well, one using master_mpi function as main on 2 nodes
170     SMPI_app_instance_register("master_mpi", master_mpi,2);
171     // the second performing an alltoall on 4 nodes
172     SMPI_app_instance_register("alltoall_mpi", alltoall_mpi,4);
173     MSG_launch_application(application_file);
174     SMPI_init();
175   }
176   res = MSG_main();
177
178   XBT_INFO("Simulation time %g", MSG_get_clock());
179
180
181   SMPI_finalize();
182   if (res == MSG_OK)
183     return 0;
184   else
185     return 1;
186 }                               /* end_of_main */