Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
[mc] Move MC_ignore_local_variable and MC_ignore_global_variable to initialisation...
[simgrid.git] / teshsuite / smpi / mpich3-test / datatype / unusual-noncontigs.c
1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
2 /*
3  *  (C) 2001 by Argonne National Laboratory.
4  *      See COPYRIGHT in top-level directory.
5  */
6 #include <math.h>
7 #include <assert.h>
8 #include <stdio.h>
9 #include <stdlib.h>
10 #include <string.h>
11 #include "mpi.h"
12
13 /* 
14    The default behavior of the test routines should be to briefly indicate
15    the cause of any errors - in this test, that means that verbose needs
16    to be set. Verbose should turn on output that is independent of error
17    levels.
18 */
19 static int verbose = 1;
20
21 int parse_args(int argc, char **argv);
22 int struct_negdisp_test(void);
23 int vector_negstride_test(void);
24 int indexed_negdisp_test(void);
25 int struct_struct_test(void);
26 int flatten_test(void);
27
28 int build_array_section_type(MPI_Aint aext, MPI_Aint astart, MPI_Aint aend, MPI_Datatype *datatype);
29
30 int main(int argc, char *argv[])
31 {
32     int err, errs = 0;
33
34     /* Initialize MPI */
35     MPI_Init(&argc, &argv);
36     parse_args(argc, argv);
37
38     /* To improve reporting of problems about operations, we
39        change the error handler to errors return */
40     MPI_Comm_set_errhandler( MPI_COMM_WORLD, MPI_ERRORS_RETURN );
41
42     err = struct_negdisp_test();
43     if (verbose && err) fprintf(stderr, "error in struct_negdisp_test\n");
44     errs += err;
45
46     err = vector_negstride_test();
47     if (verbose && err) fprintf(stderr, "error in vector_negstride_test\n");
48     errs += err;
49
50     err = indexed_negdisp_test();
51     if (verbose && err) fprintf(stderr, "error in indexed_negdisp_test\n");
52     errs += err;
53
54     err = struct_struct_test();
55     if (verbose && err) fprintf(stderr, "error in struct_struct_test\n");
56     errs += err;
57
58     err = flatten_test();
59     if (verbose && err) fprintf(stderr, "error in flatten_test\n");
60     errs += err;
61
62     /* print message and exit */
63     if (errs) {
64         fprintf(stderr, "Found %d errors\n", errs);
65     }
66     else {
67         printf(" No Errors\n");
68     }
69     MPI_Finalize();
70     return 0;
71 }
72
73 /* test uses a struct type that describes data that is contiguous,
74  * but processed in a noncontiguous way.
75  */
76 int struct_negdisp_test(void)
77 {
78     int err, errs = 0;
79     int sendbuf[6] = { 1, 2, 3, 4, 5, 6 };
80     int recvbuf[6] = { -1, -2, -3, -4, -5, -6 };
81     MPI_Datatype mystruct;
82     MPI_Request request;
83     MPI_Status status;
84
85     MPI_Aint disps[2]     = { 0,       -1*((int) sizeof(int)) };
86     int blks[2]           = { 1,       1, };
87     MPI_Datatype types[2] = { MPI_INT, MPI_INT };
88
89     err = MPI_Type_struct(2, blks, disps, types, &mystruct);
90     if (err != MPI_SUCCESS) {
91         errs++;
92         if (verbose) {
93             fprintf(stderr, "MPI_Type_struct returned error\n");
94         }
95     }
96
97     MPI_Type_commit(&mystruct);
98
99     err = MPI_Irecv(recvbuf+1, 4, MPI_INT, 0, 0, MPI_COMM_SELF, &request);
100     if (err != MPI_SUCCESS) {
101         errs++;
102         if (verbose) {
103             fprintf(stderr, "MPI_Irecv returned error\n");
104         }
105     }
106
107     err = MPI_Send(sendbuf+2, 2, mystruct, 0, 0, MPI_COMM_SELF);
108     if (err != MPI_SUCCESS) {
109         errs++;
110         if (verbose) {
111             fprintf(stderr, "MPI_Send returned error\n");
112         }
113     }
114
115     err = MPI_Wait(&request, &status);
116     if (err != MPI_SUCCESS) {
117         errs++;
118         if (verbose) {
119             fprintf(stderr, "MPI_Wait returned error\n");
120         }
121     }
122
123     /* verify data */
124     if (recvbuf[0] != -1) {
125         errs++;
126         if (verbose) {
127             fprintf(stderr, "recvbuf[0] = %d; should be %d\n", recvbuf[0], -1);
128         }
129     }
130     if (recvbuf[1] != 3) {
131         errs++;
132         if (verbose) {
133             fprintf(stderr, "recvbuf[1] = %d; should be %d\n", recvbuf[1], 3);
134         }
135     }
136     if (recvbuf[2] != 2) {
137         errs++;
138         if (verbose) {
139             fprintf(stderr, "recvbuf[2] = %d; should be %d\n", recvbuf[2], 2);
140         }
141     }
142     if (recvbuf[3] != 5) {
143         errs++;
144         if (verbose) {
145             fprintf(stderr, "recvbuf[3] = %d; should be %d\n", recvbuf[3], 5);
146         }
147     }
148     if (recvbuf[4] != 4) {
149         errs++;
150         if (verbose) {
151             fprintf(stderr, "recvbuf[4] = %d; should be %d\n", recvbuf[4], 4);
152         }
153     }
154     if (recvbuf[5] != -6) {
155         errs++;
156         if (verbose) {
157             fprintf(stderr, "recvbuf[5] = %d; should be %d\n", recvbuf[5], -6);
158         }
159     }
160
161     MPI_Type_free(&mystruct);
162
163     return errs;
164 }
165
166 /* test uses a vector type that describes data that is contiguous,
167  * but processed in a noncontiguous way.  this is effectively the
168  * same type as in the struct_negdisp_test above.
169  */
170 int vector_negstride_test(void)
171 {
172     int err, errs = 0;
173     int sendbuf[6] = { 1, 2, 3, 4, 5, 6 };
174     int recvbuf[6] = { -1, -2, -3, -4, -5, -6 };
175     MPI_Datatype myvector;
176     MPI_Request request;
177     MPI_Status status;
178
179     err = MPI_Type_vector(2, 1, -1, MPI_INT, &myvector);
180     if (err != MPI_SUCCESS) {
181         errs++;
182         if (verbose) {
183             fprintf(stderr, "MPI_Type_vector returned error\n");
184         }
185     }
186
187     MPI_Type_commit(&myvector);
188
189     err = MPI_Irecv(recvbuf+1, 4, MPI_INT, 0, 0, MPI_COMM_SELF, &request);
190     if (err != MPI_SUCCESS) {
191         errs++;
192         if (verbose) {
193             fprintf(stderr, "MPI_Irecv returned error\n");
194         }
195     }
196
197     err = MPI_Send(sendbuf+2, 2, myvector, 0, 0, MPI_COMM_SELF);
198     if (err != MPI_SUCCESS) {
199         errs++;
200         if (verbose) {
201             fprintf(stderr, "MPI_Send returned error\n");
202         }
203     }
204
205     err = MPI_Wait(&request, &status);
206     if (err != MPI_SUCCESS) {
207         errs++;
208         if (verbose) {
209             fprintf(stderr, "MPI_Wait returned error\n");
210         }
211     }
212
213     /* verify data */
214     if (recvbuf[0] != -1) {
215         errs++;
216         if (verbose) {
217             fprintf(stderr, "recvbuf[0] = %d; should be %d\n", recvbuf[0], -1);
218         }
219     }
220     if (recvbuf[1] != 3) {
221         errs++;
222         if (verbose) {
223             fprintf(stderr, "recvbuf[1] = %d; should be %d\n", recvbuf[1], 3);
224         }
225     }
226     if (recvbuf[2] != 2) {
227         errs++;
228         if (verbose) {
229             fprintf(stderr, "recvbuf[2] = %d; should be %d\n", recvbuf[2], 2);
230         }
231     }
232     if (recvbuf[3] != 5) {
233         errs++;
234         if (verbose) {
235             fprintf(stderr, "recvbuf[3] = %d; should be %d\n", recvbuf[3], 5);
236         }
237     }
238     if (recvbuf[4] != 4) {
239         errs++;
240         if (verbose) {
241             fprintf(stderr, "recvbuf[4] = %d; should be %d\n", recvbuf[4], 4);
242         }
243     }
244     if (recvbuf[5] != -6) {
245         errs++;
246         if (verbose) {
247             fprintf(stderr, "recvbuf[5] = %d; should be %d\n", recvbuf[5], -6);
248         }
249     }
250
251     MPI_Type_free(&myvector);
252
253     return errs;
254 }
255
256 /* test uses a indexed type that describes data that is contiguous,
257  * but processed in a noncontiguous way.  this is effectively the same
258  * type as in the two tests above.
259  */
260 int indexed_negdisp_test(void)
261 {
262     int err, errs = 0;
263     int sendbuf[6] = { 1, 2, 3, 4, 5, 6 };
264     int recvbuf[6] = { -1, -2, -3, -4, -5, -6 };
265     MPI_Datatype myindexed;
266     MPI_Request request;
267     MPI_Status status;
268
269     int disps[2]     = { 0, -1 };
270     int blks[2]           = { 1, 1 };
271
272     err = MPI_Type_indexed(2, blks, disps, MPI_INT, &myindexed);
273     if (err != MPI_SUCCESS) {
274         errs++;
275         if (verbose) {
276             fprintf(stderr, "MPI_Type_indexed returned error\n");
277         }
278     }
279
280     MPI_Type_commit(&myindexed);
281
282     err = MPI_Irecv(recvbuf+1, 4, MPI_INT, 0, 0, MPI_COMM_SELF, &request);
283     if (err != MPI_SUCCESS) {
284         errs++;
285         if (verbose) {
286             fprintf(stderr, "MPI_Irecv returned error\n");
287         }
288     }
289
290     err = MPI_Send(sendbuf+2, 2, myindexed, 0, 0, MPI_COMM_SELF);
291     if (err != MPI_SUCCESS) {
292         errs++;
293         if (verbose) {
294             fprintf(stderr, "MPI_Send returned error\n");
295         }
296     }
297
298     err = MPI_Wait(&request, &status);
299     if (err != MPI_SUCCESS) {
300         errs++;
301         if (verbose) {
302             fprintf(stderr, "MPI_Wait returned error\n");
303         }
304     }
305
306     /* verify data */
307     if (recvbuf[0] != -1) {
308         errs++;
309         if (verbose) {
310             fprintf(stderr, "recvbuf[0] = %d; should be %d\n", recvbuf[0], -1);
311         }
312     }
313     if (recvbuf[1] != 3) {
314         errs++;
315         if (verbose) {
316             fprintf(stderr, "recvbuf[1] = %d; should be %d\n", recvbuf[1], 3);
317         }
318     }
319     if (recvbuf[2] != 2) {
320         errs++;
321         if (verbose) {
322             fprintf(stderr, "recvbuf[2] = %d; should be %d\n", recvbuf[2], 2);
323         }
324     }
325     if (recvbuf[3] != 5) {
326         errs++;
327         if (verbose) {
328             fprintf(stderr, "recvbuf[3] = %d; should be %d\n", recvbuf[3], 5);
329         }
330     }
331     if (recvbuf[4] != 4) {
332         errs++;
333         if (verbose) {
334             fprintf(stderr, "recvbuf[4] = %d; should be %d\n", recvbuf[4], 4);
335         }
336     }
337     if (recvbuf[5] != -6) {
338         errs++;
339         if (verbose) {
340             fprintf(stderr, "recvbuf[5] = %d; should be %d\n", recvbuf[5], -6);
341         }
342     }
343
344     MPI_Type_free(&myindexed);
345
346     return errs;
347 }
348
349 #define check_err(fn_name_)                                                   \
350     do {                                                                      \
351         if (err != MPI_SUCCESS) {                                             \
352             errs++;                                                           \
353             if (verbose) {                                                    \
354                 int len_;                                                     \
355                 char err_str_[MPI_MAX_ERROR_STRING];                          \
356                 MPI_Error_string(err, err_str_, &len_);                       \
357                 fprintf(stderr, #fn_name_ " failed at line %d, err=%d: %s\n", \
358                         __LINE__, err, err_str_);                             \
359             }                                                                 \
360         }                                                                     \
361     } while (0)
362 /* test case from tt#1030 ported to C
363  *
364  * Thanks to Matthias Lieber for reporting the bug and providing a good test
365  * program. */
366 int struct_struct_test(void)
367 {
368     int err, errs = 0;
369     int i, j, dt_size = 0;
370     MPI_Request req[2];
371
372
373 #define COUNT (2)
374     MPI_Aint displ[COUNT];
375     int blens[COUNT];
376     MPI_Datatype types[COUNT];
377     MPI_Datatype datatype;
378
379     /* A slight difference from the F90 test: F90 arrays are column-major, C
380      * arrays are row-major.  So we invert the order of dimensions. */
381 #define N (2)
382 #define M (4)
383     int array[N][M] =    { {-1, -1, -1, -1}, {-1, -1, -1, -1} };
384     int expected[N][M] = { {-1,  1,  2,  5}, {-1,  3,  4,  6} };
385     int seq_array[N*M];
386     MPI_Aint astart, aend;
387     MPI_Aint size_exp = 0;
388
389     /* 1st section selects elements 1 and 2 out of 2nd dimension, complete 1st dim.
390      * should receive the values 1, 2, 3, 4 */
391     astart = 1;
392     aend   = 2;
393     err = build_array_section_type(M, astart, aend, &types[0]);
394     if (err) {
395         errs++;
396         if (verbose) fprintf(stderr, "build_array_section_type failed\n");
397         return errs;
398     }
399     blens[0] = N;
400     displ[0] = 0;
401     size_exp = size_exp + N * (aend-astart+1) * sizeof(int);
402
403     /* 2nd section selects last element of 2nd dimension, complete 1st dim.
404      * should receive the values 5, 6 */
405     astart = 3;
406     aend   = 3;
407     err = build_array_section_type(M, astart, aend, &types[1]);
408     if (err) {
409         errs++;
410         if (verbose) fprintf(stderr, "build_array_section_type failed\n");
411         return errs;
412     }
413     blens[1] = N;
414     displ[1] = 0;
415     size_exp = size_exp + N * (aend-astart+1) * sizeof(int);
416
417     /* create type */
418     err = MPI_Type_create_struct(COUNT, blens, displ, types, &datatype);
419     check_err(MPI_Type_create_struct);
420     err = MPI_Type_commit(&datatype);
421     check_err(MPI_Type_commit);
422
423     err = MPI_Type_size(datatype, &dt_size);
424     check_err(MPI_Type_size);
425     if (dt_size != size_exp) {
426         errs++;
427         if (verbose) fprintf(stderr, "unexpected type size\n");
428     }
429
430
431     /* send the type to ourselves to make sure that the type describes data correctly */
432     for (i = 0; i < (N*M) ; ++i)
433         seq_array[i] = i + 1; /* source values 1..(N*M) */
434     err = MPI_Isend(&seq_array[0], dt_size/sizeof(int), MPI_INT, 0, 42, MPI_COMM_SELF, &req[0]);
435     check_err(MPI_Isend);
436     err = MPI_Irecv(&array[0][0], 1, datatype, 0, 42, MPI_COMM_SELF, &req[1]);
437     check_err(MPI_Irecv);
438     err = MPI_Waitall(2, req, MPI_STATUSES_IGNORE);
439     check_err(MPI_Waitall);
440
441     /* check against expected */
442     for (i = 0; i < N; ++i) {
443         for (j = 0; j < M; ++j) {
444             if (array[i][j] != expected[i][j]) {
445                 errs++;
446                 if (verbose)
447                     fprintf(stderr, "array[%d][%d]=%d, should be %d\n", i, j, array[i][j], expected[i][j]);
448             }
449         }
450     }
451
452     err = MPI_Type_free(&datatype);
453     check_err(MPI_Type_free);
454     err = MPI_Type_free(&types[0]);
455     check_err(MPI_Type_free);
456     err = MPI_Type_free(&types[1]);
457     check_err(MPI_Type_free);
458
459     return errs;
460 #undef M
461 #undef N
462 #undef COUNT
463 }
464
465 /*   create a datatype for a 1D int array subsection
466
467      - a subsection of the first dimension is defined via astart, aend
468      - indexes are assumed to start with 0, that means:
469        - 0 <= astart <= aend < aext
470      - astart and aend are inclusive
471
472      example:
473
474      aext = 8, astart=2, aend=4 would produce:
475
476      index     | 0 | 1 | 2 | 3 | 4 | 5 | 6 | 7 | 8 |
477      1D array   ###############################
478      datatype   LB      ###########             UB
479  */
480 int build_array_section_type(MPI_Aint aext, MPI_Aint astart, MPI_Aint aend, MPI_Datatype *datatype)
481 {
482 #define COUNT (3)
483     int err, errs = 0;
484     MPI_Aint displ[COUNT];
485     int blens[COUNT];
486     MPI_Datatype types[COUNT];
487
488     *datatype = MPI_DATATYPE_NULL;
489
490     /* lower bound marker */
491     types[0] = MPI_LB;
492     displ[0] = 0;
493     blens[0] = 1;
494
495     /* subsection starting at astart */
496     displ[1] = astart * sizeof(int);
497     types[1] = MPI_INT;
498     blens[1] = aend - astart + 1;
499
500     /* upper bound marker */
501     types[2] = MPI_UB;
502     displ[2] = aext * sizeof(int);
503     blens[2] = 1;
504
505     err = MPI_Type_create_struct(COUNT, blens, displ, types, datatype);
506     if (err != MPI_SUCCESS) {
507         errs++;
508         if (verbose) {
509             fprintf(stderr, "MPI_Type_create_struct failed, err=%d\n", err);
510         }
511     }
512
513     return errs;
514 #undef COUNT
515 }
516
517 /* start_idx is the "zero" point for the unpack */
518 static int pack_and_check_expected(MPI_Datatype type, const char *name,
519                                    int start_idx, int size,
520                                    int *array, int *expected)
521 {
522     int i;
523     int err, errs = 0;
524     int pack_size = -1;
525     int *pack_buf = NULL;
526     int pos;
527     int type_size = -1;
528     int sendbuf[8] = {0,1,2,3,4,5,6,7};
529
530     err = MPI_Type_size(type, &type_size);
531     check_err(MPI_Type_size);
532     assert(sizeof(sendbuf) >= type_size);
533
534     err = MPI_Pack_size(type_size/sizeof(int), MPI_INT, MPI_COMM_SELF, &pack_size);
535     check_err(MPI_Pack_size);
536     pack_buf = malloc(pack_size);
537     assert(pack_buf);
538
539     pos = 0;
540     err = MPI_Pack(&sendbuf[0], type_size/sizeof(int), MPI_INT, pack_buf, pack_size, &pos, MPI_COMM_SELF);
541     check_err(MPI_Pack);
542     pos = 0;
543     err = MPI_Unpack(pack_buf, pack_size, &pos, &array[start_idx], 1, type, MPI_COMM_SELF);
544     check_err(MPI_Unpack);
545     free(pack_buf);
546
547     /* check against expected */
548     for (i = 0; i < size; ++i) {
549         if (array[i] != expected[i]) {
550             errs++;
551             if (verbose)
552                 fprintf(stderr, "%s: array[%d]=%d, should be %d\n", name, i, array[i], expected[i]);
553         }
554     }
555
556     return errs;
557 }
558
559 /* regression for tt#1030, checks for bad offset math in the
560  * blockindexed and indexed dataloop flattening code */
561 int flatten_test(void)
562 {
563     int err, errs = 0;
564 #define ARR_SIZE (9)
565     /* real indices              0  1  2  3  4  5  6  7  8
566      * indices w/ &array[3]     -3 -2 -1  0  1  2  3  4  5 */
567     int array[ARR_SIZE]      = {-1,-1,-1,-1,-1,-1,-1,-1,-1};
568     int expected[ARR_SIZE]   = {-1, 0, 1,-1, 2,-1, 3,-1, 4};
569     MPI_Datatype idx_type = MPI_DATATYPE_NULL;
570     MPI_Datatype blkidx_type = MPI_DATATYPE_NULL;
571     MPI_Datatype combo = MPI_DATATYPE_NULL;
572 #define COUNT (2)
573     int displ[COUNT];
574     MPI_Aint adispl[COUNT];
575     int blens[COUNT];
576     MPI_Datatype types[COUNT];
577
578     /* indexed type layout:
579      * XX_X
580      * 2101  <-- pos (left of 0 is neg)
581      *
582      * different blens to prevent optimization into a blockindexed
583      */
584     blens[0] = 2;
585     displ[0] = -2; /* elements, puts byte after block end at 0 */
586     blens[1] = 1;
587     displ[1] = 1; /*elements*/
588
589     err = MPI_Type_indexed(COUNT, blens, displ, MPI_INT, &idx_type);
590     check_err(MPI_Type_indexed);
591     err = MPI_Type_commit(&idx_type);
592     check_err(MPI_Type_commit);
593
594     /* indexed type layout:
595      * _X_X
596      * 2101  <-- pos (left of 0 is neg)
597      */
598     displ[0] = -1;
599     displ[1] = 1;
600     err = MPI_Type_create_indexed_block(COUNT, 1, displ, MPI_INT, &blkidx_type);
601     check_err(MPI_Type_indexed_block);
602     err = MPI_Type_commit(&blkidx_type);
603     check_err(MPI_Type_commit);
604
605     /* struct type layout:
606      * II_I_B_B  (I=idx_type, B=blkidx_type)
607      * 21012345  <-- pos (left of 0 is neg)
608      */
609     blens[0]  = 1;
610     adispl[0] = 0; /*bytes*/
611     types[0]  = idx_type;
612
613     blens[1]  = 1;
614     adispl[1] = 4 * sizeof(int); /* bytes */
615     types[1]  = blkidx_type;
616
617     /* must be a struct in order to trigger flattening code */
618     err = MPI_Type_create_struct(COUNT, blens, adispl, types, &combo);
619     check_err(MPI_Type_indexed);
620     err = MPI_Type_commit(&combo);
621     check_err(MPI_Type_commit);
622
623     /* pack/unpack with &array[3] */
624     errs += pack_and_check_expected(combo, "combo", 3, ARR_SIZE, array, expected);
625
626     MPI_Type_free(&combo);
627     MPI_Type_free(&idx_type);
628     MPI_Type_free(&blkidx_type);
629
630     return errs;
631 #undef COUNT
632 }
633 #undef check_err
634
635 int parse_args(int argc, char **argv)
636 {
637     /*
638     int ret;
639
640     while ((ret = getopt(argc, argv, "v")) >= 0)
641     {
642         switch (ret) {
643             case 'v':
644                 verbose = 1;
645                 break;
646         }
647     }
648     */
649     if (argc > 1 && strcmp(argv[1], "-v") == 0)
650         verbose = 1;
651     return 0;
652 }