1 /* -*- Mode: C; c-basic-offset:4 ; -*- */
2 /*
3  *  (C) 2001 by Argonne National Laboratory.
4  *      See COPYRIGHT in top-level directory.
5  */
6 
7 /* One-Sided MPI 2-D Strided Put Test
8  *
9  * Author: James Dinan <dinan@mcs.anl.gov>
10  * Date  : March, 2011
11  *
12  * This code performs N strided put operations into a 2d patch of a shared
13  * array.  The array has dimensions [X, Y] and the subarray has dimensions
14  * [SUB_X, SUB_Y] and begins at index [0, 0].  The input and output buffers are
15  * specified using an MPI datatype.
16  *
17  * This test generates a datatype that is relative to an arbitrary base address
18  * in memory and tests the RMA implementation's ability to perform the correct
19  * transfer.
20  */
21 
22 #include <stdio.h>
23 #include <stdlib.h>
24 #include <stdint.h>
25 #include <mpi.h>
26 #include "mpitest.h"
27 
28 #define XDIM 1024
29 #define YDIM 1024
30 #define SUB_XDIM 1024
31 #define SUB_YDIM 1024
32 #define ITERATIONS 10
33 
34 static const int SQ_LIMIT = 10;
35 static       int SQ_COUNT = 0;
36 
37 #define SQUELCH(X)                      \
38   do {                                  \
39     if (SQ_COUNT < SQ_LIMIT) {          \
40       SQ_COUNT++;                       \
41       X                                 \
42     }                                   \
43   } while (0)
44 
45 static int verbose = 0;
46 
main(int argc,char ** argv)47 int main(int argc, char **argv) {
48     int i, j, rank, nranks, peer, bufsize, errors;
49     double  *win_buf, *src_buf, *dst_buf;
50     MPI_Win buf_win;
51 
52     MTest_Init(&argc, &argv);
53 
54     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
55     MPI_Comm_size(MPI_COMM_WORLD, &nranks);
56 
57     bufsize = XDIM * YDIM * sizeof(double);
58     MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &win_buf);
59     /* Alloc_mem is not required for the origin buffers for RMA operations -
60        just for the Win_create memory */
61     MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &src_buf);
62     MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &dst_buf);
63 
64     if (rank == 0)
65         if (verbose) printf("MPI RMA Strided Put Test:\n");
66 
67     for (i = 0; i < XDIM*YDIM; i++) {
68         *(win_buf  + i) = 1.0 + rank;
69         *(src_buf + i) = 1.0 + rank;
70     }
71 
72     MPI_Win_create(win_buf, bufsize, 1, MPI_INFO_NULL, MPI_COMM_WORLD, &buf_win);
73 
74     peer = (rank+1) % nranks;
75 
76     /* Perform ITERATIONS strided put operations */
77 
78     for (i = 0; i < ITERATIONS; i++) {
79       MPI_Aint idx_loc[SUB_YDIM];
80       int idx_rem[SUB_YDIM];
81       int blk_len[SUB_YDIM];
82       MPI_Datatype src_type, dst_type;
83 
84       void *base_ptr = dst_buf;
85       MPI_Aint base_int;
86 
87       MPI_Get_address(base_ptr, &base_int);
88 
89       if (rank == 0)
90         if (verbose) printf(" + iteration %d\n", i);
91 
92       for (j = 0; j < SUB_YDIM; j++) {
93         MPI_Get_address(&src_buf[j*XDIM], &idx_loc[j]);
94         idx_loc[j] = idx_loc[j] - base_int;
95         idx_rem[j] = j*XDIM*sizeof(double);
96         blk_len[j] = SUB_XDIM*sizeof(double);
97       }
98 
99       MPI_Type_create_hindexed(SUB_YDIM, blk_len, idx_loc, MPI_BYTE, &src_type);
100       MPI_Type_create_indexed_block(SUB_YDIM, SUB_XDIM*sizeof(double), idx_rem, MPI_BYTE, &dst_type);
101 
102       MPI_Type_commit(&src_type);
103       MPI_Type_commit(&dst_type);
104 
105       MPI_Win_lock(MPI_LOCK_EXCLUSIVE, peer, 0, buf_win);
106       MPI_Put(base_ptr, 1, src_type, peer, 0, 1, dst_type, buf_win);
107       MPI_Win_unlock(peer, buf_win);
108 
109       MPI_Type_free(&src_type);
110       MPI_Type_free(&dst_type);
111     }
112 
113     MPI_Barrier(MPI_COMM_WORLD);
114 
115     /* Verify that the results are correct */
116 
117     MPI_Win_lock(MPI_LOCK_EXCLUSIVE, rank, 0, buf_win);
118     errors = 0;
119     for (i = 0; i < SUB_XDIM; i++) {
120       for (j = 0; j < SUB_YDIM; j++) {
121         const double actual   = *(win_buf + i + j*XDIM);
122         const double expected = (1.0 + ((rank+nranks-1)%nranks));
123         if (actual - expected > 1e-10) {
124           SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
125               rank, j, i, expected, actual); );
126           errors++;
127           fflush(stdout);
128         }
129       }
130     }
131     for (i = SUB_XDIM; i < XDIM; i++) {
132       for (j = 0; j < SUB_YDIM; j++) {
133         const double actual   = *(win_buf + i + j*XDIM);
134         const double expected = 1.0 + rank;
135         if (actual - expected > 1e-10) {
136           SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
137               rank, j, i, expected, actual); );
138           errors++;
139           fflush(stdout);
140         }
141       }
142     }
143     for (i = 0; i < XDIM; i++) {
144       for (j = SUB_YDIM; j < YDIM; j++) {
145         const double actual   = *(win_buf + i + j*XDIM);
146         const double expected = 1.0 + rank;
147         if (actual - expected > 1e-10) {
148           SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
149               rank, j, i, expected, actual); );
150           errors++;
151           fflush(stdout);
152         }
153       }
154     }
155     MPI_Win_unlock(rank, buf_win);
156 
157     MPI_Win_free(&buf_win);
158     MPI_Free_mem(win_buf);
159     MPI_Free_mem(src_buf);
160     MPI_Free_mem(dst_buf);
161 
162     MTest_Finalize( errors );
163     MPI_Finalize();
164 
165     return 0;
166 }
167