Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'pikachuyann/simgrid-stoprofiles'
[simgrid.git] / teshsuite / smpi / mpich3-test / datatype / sendrecvt4.c
1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
2 /*
3  *  (C) 2014 by Argonne National Laboratory.
4  *      See COPYRIGHT in top-level directory.
5  */
6 #include "mpi.h"
7 #include "mpitest.h"
8 #include <stdio.h>
9 #include "dtypes.h"
10
11
12 /*
13    This program is derived from one in the MPICH-1 test suite
14
15    This version sends and receives EVERYTHING from MPI_BOTTOM, by putting
16    the data into a structure.
17  */
18 int main(int argc, char **argv)
19 {
20     MPI_Datatype *types;
21     void **inbufs, **outbufs;
22     int *counts, *bytesize, ntype;
23     MPI_Comm comm;
24     int rank, np, partner, tag, count;
25     int j, k, err, world_rank, errloc;
26     MPI_Status status;
27     char *obuf;
28     MPI_Datatype offsettype;
29     int blen;
30     MPI_Aint displ, extent, natural_extent;
31     char myname[MPI_MAX_OBJECT_NAME];
32     int mynamelen;
33
34     MTest_Init(&argc, &argv);
35
36     MTestDatatype2Allocate(&types, &inbufs, &outbufs, &counts, &bytesize, &ntype);
37     MTestDatatype2Generate(types, inbufs, outbufs, counts, bytesize, &ntype);
38
39     MPI_Comm_rank(MPI_COMM_WORLD, &world_rank);
40
41     /* Test over a wide range of datatypes and communicators */
42     err = 0;
43     tag = 0;
44     while (MTestGetIntracomm(&comm, 2)) {
45         if (comm == MPI_COMM_NULL)
46             continue;
47         MPI_Comm_rank(comm, &rank);
48         MPI_Comm_size(comm, &np);
49         if (np < 2)
50             continue;
51         tag++;
52         for (j = 0; j < ntype; j++) {
53             MPI_Type_get_name(types[j], myname, &mynamelen);
54             if (world_rank == 0)
55                 MTestPrintfMsg(10, "Testing type %s\n", myname);
56             if (rank == 0) {
57                 MPI_Get_address(inbufs[j], &displ);
58                 blen = 1;
59                 MPI_Type_create_struct(1, &blen, &displ, types + j, &offsettype);
60                 MPI_Type_commit(&offsettype);
61                 /* Warning: if the type has an explicit MPI_UB, then using a
62                  * simple shift of the offset won't work.  For now, we skip
63                  * types whose extents are negative; the correct solution is
64                  * to add, where required, an explicit MPI_UB */
65                 MPI_Type_extent(offsettype, &extent);
66                 if (extent < 0) {
67                     if (world_rank == 0)
68                         MTestPrintfMsg(10, "... skipping (appears to have explicit MPI_UB\n");
69                     MPI_Type_free(&offsettype);
70                     continue;
71                 }
72                 MPI_Type_extent(types[j], &natural_extent);
73                 if (natural_extent != extent) {
74                     MPI_Type_free(&offsettype);
75                     continue;
76                 }
77                 partner = np - 1;
78                 MPI_Datatype dup;
79                 MPI_Type_dup(offsettype, &dup);
80                 MPI_Send(MPI_BOTTOM, counts[j], dup, partner, tag, comm);
81                 MPI_Type_free(&offsettype);
82                 MPI_Type_free(&dup);
83             }
84             else if (rank == np - 1) {
85                 partner = 0;
86                 obuf = outbufs[j];
87                 for (k = 0; k < bytesize[j]; k++)
88                     obuf[k] = 0;
89                 MPI_Get_address(outbufs[j], &displ);
90                 blen = 1;
91                 MPI_Type_create_struct(1, &blen, &displ, types + j, &offsettype);
92                 MPI_Type_commit(&offsettype);
93                 /* Warning: if the type has an explicit MPI_UB, then using a
94                  * simple shift of the offset won't work.  For now, we skip
95                  * types whose extents are negative; the correct solution is
96                  * to add, where required, an explicit MPI_UB */
97                 MPI_Type_extent(offsettype, &extent);
98                 if (extent < 0) {
99                     MPI_Type_free(&offsettype);
100                     continue;
101                 }
102                 MPI_Type_extent(types[j], &natural_extent);
103                 if (natural_extent != extent) {
104                     MPI_Type_free(&offsettype);
105                     continue;
106                 }
107                 MPI_Datatype dup;
108                 MPI_Type_dup(offsettype, &dup);
109                 MPI_Recv(MPI_BOTTOM, counts[j], dup, partner, tag, comm, &status);
110                 /* Test for correctness */
111                 MPI_Get_count(&status, types[j], &count);
112                 if (count != counts[j]) {
113                     fprintf(stderr,
114                             "Error in counts (got %d expected %d) with type %s\n",
115                             count, counts[j], myname);
116                     err++;
117                 }
118                 if (status.MPI_SOURCE != partner) {
119                     fprintf(stderr,
120                             "Error in source (got %d expected %d) with type %s\n",
121                             status.MPI_SOURCE, partner, myname);
122                     err++;
123                 }
124                 if ((errloc = MTestDatatype2Check(inbufs[j], outbufs[j], bytesize[j]))) {
125                     fprintf(stderr,
126                             "Error in data with type %s (type %d on %d) at byte %d\n",
127                             myname, j, world_rank, errloc - 1);
128                     if (err < 10) {
129                         /* Give details on only the first 10 errors */
130                         unsigned char *in_p = (unsigned char *) inbufs[j],
131                             *out_p = (unsigned char *) outbufs[j];
132                         int jj;
133                         jj = errloc - 1;
134                         jj &= 0xfffffffc;       /* lop off a few bits */
135                         in_p += jj;
136                         out_p += jj;
137                         fprintf(stderr, "%02x%02x%02x%02x should be %02x%02x%02x%02x\n",
138                                 out_p[0], out_p[1], out_p[2], out_p[3],
139                                 in_p[0], in_p[1], in_p[2], in_p[3]);
140                     }
141                     err++;
142                 }
143                 MPI_Type_free(&offsettype);
144                 MPI_Type_free(&dup);
145             }
146         }
147         MTestFreeComm(&comm);
148     }
149
150     MTestDatatype2Free(types, inbufs, outbufs, counts, bytesize, ntype);
151     MTest_Finalize(err);
152     MPI_Finalize();
153     return MTestReturnValue(err);
154 }