Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
1fe801ae11f4116b87241d4137f5ab802cdf6b88
[simgrid.git] / examples / smpi / smpi_msg_masterslave / masterslave_mailbox_smpi.c
1 /* Copyright (c) 2010-2019. The SimGrid Team. All rights reserved.          */
2
3 /* This program is free software; you can redistribute it and/or modify it
4  * under the terms of the license (GNU LGPL) which comes with this package. */
5
6 #include "simgrid/msg.h"
7 #include "mpi.h"
8
9 #include <stdio.h> /* snprintf */
10
11 XBT_LOG_NEW_DEFAULT_CATEGORY(msg_test, "Messages specific for this msg example");
12
13 static int master(int argc, char *argv[])
14 {
15   long number_of_tasks = xbt_str_parse_int(argv[1], "Invalid amount of tasks: %s");
16   double task_comp_size = xbt_str_parse_double(argv[2], "Invalid computational size: %s");
17   double task_comm_size = xbt_str_parse_double(argv[3], "Invalid communication size: %s");
18   long slaves_count = xbt_str_parse_int(argv[4], "Invalid amount of slaves: %s");
19
20   XBT_INFO("Got %ld slaves and %ld tasks to process", slaves_count, number_of_tasks);
21
22   for (int i = 0; i < number_of_tasks; i++) {
23     char mailbox[256];
24     char sprintf_buffer[256];
25
26     snprintf(mailbox,256, "slave-%ld", i % slaves_count);
27     snprintf(sprintf_buffer,256, "Task_%d", i);
28     msg_task_t task = MSG_task_create(sprintf_buffer, task_comp_size, task_comm_size, NULL);
29     if (number_of_tasks < 10000 || i % 10000 == 0)
30       XBT_INFO("Sending \"%s\" (of %ld) to mailbox \"%s\"", MSG_task_get_name(task), number_of_tasks, mailbox);
31
32     MSG_task_send(task, mailbox);
33   }
34
35   XBT_INFO("All tasks have been dispatched. Let's tell everybody the computation is over.");
36   for (int i = 0; i < slaves_count; i++) {
37     char mailbox[80];
38
39     snprintf(mailbox,80, "slave-%ld", i % slaves_count);
40     msg_task_t finalize = MSG_task_create("finalize", 0, 0, 0);
41     MSG_task_send(finalize, mailbox);
42   }
43
44   return 0;
45 }
46
47 static int master_mpi(int argc, char *argv[])
48 {
49   MPI_Init(&argc, &argv);
50
51   int rank;
52   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
53   XBT_INFO("here for rank %d", rank);
54   int test[1000]={rank};
55   if(rank==0)
56     MPI_Send(&test, 1000, MPI_INT, 1, 1, MPI_COMM_WORLD);
57   else
58     MPI_Recv(&test, 1000, MPI_INT, 0, 1, MPI_COMM_WORLD, MPI_STATUSES_IGNORE);
59
60   XBT_INFO("After comm %d", rank);
61   MPI_Finalize();
62
63   XBT_INFO("After finalize %d %d", rank, test[0]);
64   return 0;
65 }
66
67 static int alltoall_mpi(int argc, char *argv[])
68 {
69   MPI_Init(&argc, &argv);
70
71   int rank;
72   int size;
73   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
74   MPI_Comm_size(MPI_COMM_WORLD, &size);
75   XBT_INFO("alltoall for rank %d", rank);
76   int* out=malloc(1000*size*sizeof(int));
77   int* in=malloc(1000*size*sizeof(int));
78   MPI_Alltoall(out, 1000, MPI_INT,in, 1000, MPI_INT, MPI_COMM_WORLD);
79
80   XBT_INFO("after alltoall %d", rank);
81   free(out);
82   free(in);
83   MPI_Finalize();
84   return 0;
85 }
86
87 static int slave(int argc, char *argv[])
88 {
89   msg_task_t task = NULL;
90   XBT_ATTRIB_UNUSED int res;
91   int id = -1;
92   char mailbox[80];
93   XBT_ATTRIB_UNUSED int read;
94
95   read = sscanf(argv[1], "%d", &id);
96   xbt_assert(read, "Invalid argument %s\n", argv[1]);
97
98   snprintf(mailbox,80, "slave-%d", id);
99
100   while (1) {
101     res = MSG_task_receive(&(task), mailbox);
102     xbt_assert(res == MSG_OK, "MSG_task_get failed");
103
104     if (strcmp(MSG_task_get_name(task), "finalize")==0) {
105       MSG_task_destroy(task);
106       break;
107     }
108     MSG_task_execute(task);
109     MSG_task_destroy(task);
110     task = NULL;
111   }
112   XBT_INFO("I'm done. See you!");
113
114   return 0;
115 }
116
117 int main(int argc, char *argv[])
118 {
119   msg_error_t res;
120
121   MSG_init(&argc, argv);
122   SMPI_init();
123
124   xbt_assert(argc > 2,"Usage: %s platform_file deployment_file\n"
125              "\nexample: %s msg_platform.xml msg_deployment.xml\n", argv[0], argv[0]);
126
127   MSG_create_environment(argv[1]);
128
129   MSG_function_register("master", master);
130   MSG_function_register("slave", slave);
131   // launch two MPI applications as well, one using master_mpi function as main on 2 nodes
132   SMPI_app_instance_register("master_mpi", master_mpi,2);
133   // the second performing an alltoall on 4 nodes
134   SMPI_app_instance_register("alltoall_mpi", alltoall_mpi,4);
135   MSG_launch_application(argv[2]);
136
137   res = MSG_main();
138
139   XBT_INFO("Simulation time %g", MSG_get_clock());
140
141   SMPI_finalize();
142   return res != MSG_OK;
143 }