Logo AND Algorithmique Numérique Distribuée

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