Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Fix HAVE_FOOBAR flags handling
[simgrid.git] / teshsuite / smpi / mpich3-test / datatype / get-struct.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 <string.h>
10
11 /* Communicating a datatype built out of structs
12  * This test was motivated by the failure of an example program for
13  * RMA involving simple operations on a struct that included a struct
14  *
15  * The observed failure was a SEGV in the MPI_Get
16  *
17  * 
18  */
19 #define MAX_KEY_SIZE 64
20 #define MAX_VALUE_SIZE 256
21 typedef struct {
22     MPI_Aint disp;
23     int      rank;
24     void     *lptr;
25 } Rptr;
26 typedef struct {
27     Rptr     next;
28     char     key[MAX_KEY_SIZE], value[MAX_VALUE_SIZE];
29 } ListElm;
30 Rptr nullDptr = {0,-1,0};
31
32 int testCases = -1;
33 #define BYTE_ONLY 0x1
34 #define TWO_STRUCT 0x2
35 int isOneLevel = 0;
36
37 int main(int argc, char **argv)
38 {
39     int  errors=0;
40     Rptr         headDptr;
41     ListElm      *headLptr=0;
42     int          i, wrank;
43     MPI_Datatype dptrType, listelmType;
44     MPI_Win      listwin;
45
46     MTest_Init(&argc, &argv);
47     MPI_Comm_rank(MPI_COMM_WORLD,&wrank);
48
49     for (i=1; i<argc; i++) {
50         if (strcmp(argv[i], "-byteonly") == 0) {
51             testCases = BYTE_ONLY;
52         }
53         else if (strcmp(argv[i], "-twostruct") == 0) {
54             testCases = TWO_STRUCT;
55         }
56         else if (strcmp(argv[i], "-onelevel") == 0) {
57             isOneLevel = 1;
58         }
59         else {
60             printf("Unrecognized argument %s\n", argv[i] );
61         }
62     }
63
64     /* Create the datatypes that we will use to move the data */
65     {
66         int      blens[3];
67         MPI_Aint displ[3];
68         MPI_Datatype dtypes[3];
69         ListElm  sampleElm;
70
71         blens[0] = 1; blens[1] = 1;
72         MPI_Get_address( &nullDptr.disp, &displ[0] );
73         MPI_Get_address( &nullDptr.rank, &displ[1] );
74         displ[1] = displ[1] - displ[0];
75         displ[0] = 0;
76         dtypes[0] = MPI_AINT;
77         dtypes[1] = MPI_INT;
78         MPI_Type_create_struct(2, blens, displ, dtypes, &dptrType);
79         MPI_Type_commit(&dptrType);
80
81         if (isOneLevel) {
82             blens[0]  = sizeof(nullDptr); dtypes[0] = MPI_BYTE;
83         }
84         else {
85             blens[0] = 1;                 dtypes[0] = dptrType;
86         }
87         blens[1] = MAX_KEY_SIZE;   dtypes[1] = MPI_CHAR;
88         blens[2] = MAX_VALUE_SIZE; dtypes[2] = MPI_CHAR;
89         MPI_Get_address(&sampleElm.next,&displ[0]);
90         MPI_Get_address(&sampleElm.key[0],&displ[1]);
91         MPI_Get_address(&sampleElm.value[0],&displ[2]);
92         displ[2] -= displ[0];
93         displ[1] -= displ[0];
94         displ[0] = 0;
95         for (i=0; i<3; i++) {
96             MTestPrintfMsg(0,"%d:count=%d,disp=%ld\n",i, blens[i], displ[i]);
97         }
98         MPI_Type_create_struct(3, blens, displ, dtypes, &listelmType);
99         MPI_Type_commit(&listelmType);
100     }
101
102     MPI_Win_create_dynamic(MPI_INFO_NULL, MPI_COMM_WORLD, &listwin);
103
104     headDptr.rank = 0;
105     if (wrank == 0) {
106         /* Create 1 list element (the head) and initialize it */
107         MPI_Alloc_mem(sizeof(ListElm), MPI_INFO_NULL, &headLptr);
108         MPI_Get_address(headLptr, &headDptr.disp);
109         headLptr->next.rank = -1;
110         headLptr->next.disp = (MPI_Aint)MPI_BOTTOM;
111         headLptr->next.lptr = 0;
112         strncpy(headLptr->key,"key1",MAX_KEY_SIZE);
113         strncpy(headLptr->value,"value1",MAX_VALUE_SIZE);
114         MPI_Win_attach(listwin, headLptr, sizeof(ListElm));
115     }
116     MPI_Bcast(&headDptr.disp, 1, MPI_AINT, 0, MPI_COMM_WORLD);
117
118     MPI_Barrier(MPI_COMM_WORLD);
119
120     if (wrank == 1) {
121         ListElm headcopy;
122
123         MPI_Win_lock_all(0, listwin);
124         /* Get head element with simple get of BYTES */
125         if (testCases & BYTE_ONLY) {
126             headcopy.next.rank=100;
127             headcopy.next.disp=0xefefefef;
128             MPI_Get(&headcopy, sizeof(ListElm), MPI_BYTE,
129                     headDptr.rank, headDptr.disp,
130                     sizeof(ListElm), MPI_BYTE, listwin);
131             MPI_Win_flush(headDptr.rank, listwin);
132             if (headcopy.next.rank != -1 &&
133                 headcopy.next.disp != (MPI_Aint)MPI_BOTTOM) {
134                 errors++;
135                 printf("MPI_BYTE: headcopy contains incorrect next:<%d,%ld>\n",
136                        headcopy.next.rank, (long)headcopy.next.disp);
137             }
138         }
139
140         if (testCases & TWO_STRUCT) {
141             headcopy.next.rank=100;
142             headcopy.next.disp=0xefefefef;
143             /* Get head element using struct of struct type.  This is
144                not an identical get to the simple BYTE one above but should
145                work */
146             MPI_Get(&headcopy, 1, listelmType, headDptr.rank, headDptr.disp,
147                     1, listelmType, listwin);
148             MPI_Win_flush(headDptr.rank, listwin);
149             if (headcopy.next.rank != -1 &&
150                 headcopy.next.disp != (MPI_Aint)MPI_BOTTOM) {
151                 errors++;
152                 printf("ListelmType: headcopy contains incorrect next:<%d,%ld>\n",
153                        headcopy.next.rank, (long)headcopy.next.disp);
154             }
155         }
156
157         MPI_Win_unlock_all(listwin);
158     }
159     MPI_Barrier(MPI_COMM_WORLD);
160     if (wrank == 0) {
161         MPI_Win_detach(listwin,headLptr);
162         MPI_Free_mem(headLptr);
163     }
164     MPI_Win_free(&listwin);
165     MPI_Type_free(&dptrType);
166     MPI_Type_free(&listelmType);
167
168     MTest_Finalize( errors );
169     MPI_Finalize();
170     return MTestReturnValue( errors );
171 }