Logo AND Algorithmique Numérique Distribuée

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