Logo AND Algorithmique Numérique Distribuée

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