Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Code smells spotted by static analyzers.
[simgrid.git] / teshsuite / smpi / gh-139 / gh-139.c
1 /* Copyright (c) 2019. 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 <mpi.h>
8 #include <simgrid/msg.h>
9 #include <simgrid/simix.h>
10 #include <stdio.h>
11
12 XBT_LOG_NEW_DEFAULT_CATEGORY(smpi_test, "Messages specific for this SMPI example");
13
14 struct param {
15   MPI_Request* req;
16   int rank;
17 };
18
19 struct threadwrap {
20   void* father_data;
21   void* (*f)(void*);
22   void* param;
23 };
24
25 static int global_rank;
26 void* req_wait(void* bar);
27
28 // Thread creation helper
29
30 static int thread_create_wrapper(XBT_ATTRIB_UNUSED int argc, XBT_ATTRIB_UNUSED char* argv[])
31 {
32   int the_global_rank  = global_rank;
33   struct threadwrap* t = (struct threadwrap*)sg_actor_self_data();
34   XBT_INFO("new thread has parameter rank %d and global variable rank %d", ((struct param*)(t->param))->rank,
35            the_global_rank);
36   sg_actor_self_data_set(t->father_data);
37   t->f(t->param);
38   sg_actor_self_data_set(NULL);
39   free(t);
40   return 0;
41 }
42
43 static void mpi_thread_create(const char* name, void* (*f)(void*), void* param)
44 {
45   struct threadwrap* threadwrap = (struct threadwrap*)malloc(sizeof(*threadwrap));
46   threadwrap->father_data       = sg_actor_self_data();
47   threadwrap->f                 = f;
48   threadwrap->param             = param;
49   sg_actor_t actor              = sg_actor_init(name, sg_host_self());
50   sg_actor_data_set(actor, threadwrap);
51   sg_actor_start(actor, thread_create_wrapper, 0, NULL);
52 }
53
54 static void mythread_create(const char* name, MPI_Request* req, int rank)
55 {
56   struct param* param = (struct param*)malloc(sizeof(*param));
57   param->req          = req;
58   param->rank         = rank;
59   mpi_thread_create(name, req_wait, param);
60 }
61
62 // Actual application
63
64 void* req_wait(void* bar)
65 {
66   struct param* param = (struct param*)bar;
67   int rank;
68   MPI_Status status;
69
70   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
71
72   XBT_INFO("%d has MPI rank %d and global variable rank %d", param->rank, rank, global_rank);
73   XBT_INFO("%d waiting request", rank);
74   MPI_Wait(param->req, &status);
75   XBT_INFO("%d request done", rank);
76   XBT_INFO("%d still has MPI rank %d and global variable %d", param->rank, rank, global_rank);
77   free(param);
78   return NULL;
79 }
80
81 int main(int argc, char* argv[])
82 {
83   int rank;
84   int size;
85   char c = 0;
86   MPI_Request req;
87   MPI_Init(&argc, &argv);
88   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
89   MPI_Comm_size(MPI_COMM_WORLD, &size);
90   XBT_INFO("I'm %d/%d", rank, size);
91   global_rank = rank;
92
93   if (rank == 0) {
94     c = 42;
95     MPI_Isend(&c, 1, MPI_CHAR, 1, 0, MPI_COMM_WORLD, &req);
96     mythread_create("wait send", &req, rank);
97   } else if (rank == 1) {
98     MPI_Irecv(&c, 1, MPI_CHAR, 0, 0, MPI_COMM_WORLD, &req);
99     mythread_create("wait recv", &req, rank);
100   }
101
102   sg_actor_sleep_for(1 + rank);
103   MPI_Finalize();
104   XBT_INFO("finally %d", c);
105   return 0;
106 }