Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
simple example of smpi tracing with platform utilization (with three categories)
[simgrid.git] / examples / smpi / mvmul.c
1 /* Copyright (c) 2009, 2010. 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 <stdlib.h>
9 #include <ctype.h>
10 #include <mpi.h>
11
12 #define ITERATIONS         10
13
14 #define USAGE_ERROR        1
15 #define SANITY_ERROR       2
16 #define GETTIMEOFDAY_ERROR 3
17
18 int main(int argc, char *argv[])
19 {
20
21   int size, rank;
22   int N, n, i, j, k, current_iteration, successful_iterations = 0;
23   double *matrix = NULL, *vector = NULL, *vcalc, *vcheck;
24   MPI_Status status;
25   struct timeval *start_time = NULL, *stop_time = NULL;
26   long parallel_usecs, parallel_usecs_total =
27     0, sequential_usecs, sequential_usecs_total = 0;
28
29   MPI_Init(&argc, &argv);
30
31   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
32   MPI_Comm_size(MPI_COMM_WORLD, &size);
33
34   if (0 == rank) {
35
36     // root node parses cmdline args
37     if (2 > argc || !isdigit(*argv[1])) {
38       printf("usage:\n%s <size>\n", argv[0]);
39       MPI_Abort(MPI_COMM_WORLD, USAGE_ERROR);
40       exit(USAGE_ERROR);
41     }
42
43     N = atoi(argv[1]);
44
45     start_time = (struct timeval *) malloc(sizeof(struct timeval));
46     stop_time = (struct timeval *) malloc(sizeof(struct timeval));
47
48   }
49
50   for (current_iteration = 0; current_iteration < ITERATIONS;
51        current_iteration++) {
52
53     if (0 == rank) {
54
55       matrix = (double *) malloc(N * N * sizeof(double));
56       vector = (double *) malloc(N * sizeof(double));
57
58       for (i = 0; i < N * N; i++) {
59         matrix[i] = (double) rand() / ((double) RAND_MAX + 1);
60       }
61
62       for (i = 0; i < N; i++) {
63         vector[i] = (double) rand() / ((double) RAND_MAX + 1);
64       }
65
66       // for the sake of argument, the parallel algorithm begins
67       // when the root node begins to transmit the matrix to the
68       // workers.
69       if (-1 == gettimeofday(start_time, NULL)) {
70         printf("couldn't set start_time on node 0!\n");
71         MPI_Abort(MPI_COMM_WORLD, GETTIMEOFDAY_ERROR);
72         exit(GETTIMEOFDAY_ERROR);
73       }
74
75       for (i = 1; i < size; i++) {
76         MPI_Send(&N, 1, MPI_INT, i, 0, MPI_COMM_WORLD);
77       }
78
79     } else {
80       MPI_Recv(&N, 1, MPI_INT, 0, 0, MPI_COMM_WORLD, &status);
81     }
82
83     // this algorithm uses at most N processors...
84     if (rank < N) {
85
86       if (size > N)
87         size = N;
88       n = N / size + ((rank < (N % size)) ? 1 : 0);
89
90       if (0 == rank) {
91
92         for (i = 1, j = n; i < size && j < N; i++, j += k) {
93           k = N / size + ((i < (N % size)) ? 1 : 0);
94           MPI_Send(matrix + N * j, N * k, MPI_DOUBLE, i, 0, MPI_COMM_WORLD);
95           MPI_Send(vector, N, MPI_DOUBLE, i, 0, MPI_COMM_WORLD);
96         }
97
98         // sanity check
99 #ifdef DEBUG
100         if (i != size || j != N) {
101           printf("index calc error: i = %d, size = %d, j = %d, N = %d\n", i,
102                  size, j, N);
103           MPI_Abort(MPI_COMM_WORLD, SANITY_ERROR);
104           exit(SANITY_ERROR);
105         }
106 #endif
107
108         vcalc = (double *) malloc(N * sizeof(double));
109
110       } else {
111
112         matrix = (double *) malloc(N * n * sizeof(double));
113         vector = (double *) malloc(N * sizeof(double));
114         vcalc = (double *) malloc(n * sizeof(double));
115
116         MPI_Recv(matrix, N * n, MPI_DOUBLE, 0, 0, MPI_COMM_WORLD, &status);
117         MPI_Recv(vector, N, MPI_DOUBLE, 0, 0, MPI_COMM_WORLD, &status);
118
119       }
120
121       for (i = 0; i < n; i++) {
122         vcalc[i] = 0.0;
123         for (j = 0; j < N; j++) {
124           vcalc[i] += matrix[N * i + j] * vector[j];
125         }
126       }
127
128       if (0 != rank) {
129         MPI_Send(vcalc, n, MPI_DOUBLE, 0, 0, MPI_COMM_WORLD);
130       } else {
131
132         for (i = 1, j = n; i < size && j < N; i++, j += k) {
133           k = N / size + ((i < (N % size)) ? 1 : 0);
134           MPI_Recv(vcalc + j, k, MPI_DOUBLE, i, 0, MPI_COMM_WORLD, &status);
135         }
136
137         // sanity check
138 #ifdef DEBUG
139         if (i != size || j != N) {
140           printf("index calc error 2: i = %d, size = %d, j = %d, N = %d\n", i,
141                  size, j, N);
142           MPI_Abort(MPI_COMM_WORLD, SANITY_ERROR);
143           exit(SANITY_ERROR);
144         }
145 #endif
146
147         if (-1 == gettimeofday(stop_time, NULL)) {
148           printf("couldn't set stop_time on node 0!\n");
149           MPI_Abort(MPI_COMM_WORLD, GETTIMEOFDAY_ERROR);
150           exit(GETTIMEOFDAY_ERROR);
151         }
152
153         parallel_usecs =
154           (stop_time->tv_sec * 1000000 + stop_time->tv_usec) -
155           (start_time->tv_sec * 1000000 + start_time->tv_usec);
156
157         if (-1 == gettimeofday(start_time, NULL)) {
158           printf("couldn't set start_time on node 0!\n");
159           MPI_Abort(MPI_COMM_WORLD, GETTIMEOFDAY_ERROR);
160           exit(GETTIMEOFDAY_ERROR);
161         }
162         // calculate serially
163         vcheck = (double *) malloc(N * sizeof(double));
164         for (i = 0; i < N; i++) {
165           vcheck[i] = 0.0;
166           for (j = 0; j < N; j++) {
167             vcheck[i] += matrix[N * i + j] * vector[j];
168           }
169         }
170
171         if (-1 == gettimeofday(stop_time, NULL)) {
172           printf("couldn't set stop_time on node 0!\n");
173           MPI_Abort(MPI_COMM_WORLD, GETTIMEOFDAY_ERROR);
174           exit(GETTIMEOFDAY_ERROR);
175         }
176
177         sequential_usecs =
178           (stop_time->tv_sec * 1000000 + stop_time->tv_usec) -
179           (start_time->tv_sec * 1000000 + start_time->tv_usec);
180
181         // verify correctness
182         for (i = 0; i < N && vcalc[i] == vcheck[i]; i++);
183
184         printf("prog: blocking, i: %d ", current_iteration);
185
186         if (i == N) {
187           printf
188             ("ptime: %ld us, stime: %ld us, speedup: %.3f, nodes: %d, efficiency: %.3f\n",
189              parallel_usecs, sequential_usecs,
190              (double) sequential_usecs / (double) parallel_usecs, size,
191              (double) sequential_usecs / ((double) parallel_usecs *
192                                           (double) size));
193
194           parallel_usecs_total += parallel_usecs;
195           sequential_usecs_total += sequential_usecs;
196           successful_iterations++;
197         } else {
198           printf("parallel calc != serial calc, ");
199         }
200
201         free(vcheck);
202
203       }
204
205       free(matrix);
206       free(vector);
207       free(vcalc);
208     }
209
210   }
211
212   if (0 == rank) {
213     printf("prog: blocking, ");
214     if (0 < successful_iterations) {
215       printf
216         ("iterations: %d, avg. ptime: %.3f us, avg. stime: %.3f us, avg. speedup: %.3f, nodes: %d, avg. efficiency: %.3f\n",
217          successful_iterations,
218          (double) parallel_usecs_total / (double) successful_iterations,
219          (double) sequential_usecs_total / (double) successful_iterations,
220          (double) sequential_usecs_total / (double) parallel_usecs_total,
221          size,
222          (double) sequential_usecs_total / ((double) parallel_usecs_total *
223                                             (double) size));
224     } else {
225       printf("no successful iterations!\n");
226     }
227
228     free(start_time);
229     free(stop_time);
230
231   }
232
233   MPI_Finalize();
234
235   return 0;
236 }