1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
3 * (C) 2001 by Argonne National Laboratory.
4 * See COPYRIGHT in top-level directory.
7 /* One-Sided MPI 2-D Strided Put Test
9 * Author: James Dinan <dinan@mcs.anl.gov>
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.
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
35 int main(int argc, char **argv) {
36 int i, j, rank, nranks, peer, bufsize, errors;
37 double *win_buf, *src_buf, *dst_buf;
40 MTest_Init(&argc, &argv);
42 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
43 MPI_Comm_size(MPI_COMM_WORLD, &nranks);
45 bufsize = XDIM * YDIM * sizeof(double);
46 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &win_buf);
47 /* Alloc_mem is not required for the origin buffers for RMA operations -
48 just for the Win_create memory */
49 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &src_buf);
50 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &dst_buf);
52 for (i = 0; i < XDIM*YDIM; i++) {
53 *(win_buf + i) = 1.0 + rank;
54 *(src_buf + i) = 1.0 + rank;
57 MPI_Win_create(win_buf, bufsize, 1, MPI_INFO_NULL, MPI_COMM_WORLD, &buf_win);
59 peer = (rank+1) % nranks;
61 /* Perform ITERATIONS strided put operations */
63 for (i = 0; i < ITERATIONS; i++) {
64 MPI_Aint idx_loc[SUB_YDIM];
65 int idx_rem[SUB_YDIM];
66 int blk_len[SUB_YDIM];
67 MPI_Datatype src_type, dst_type;
69 void *base_ptr = dst_buf;
72 MPI_Get_address(base_ptr, &base_int);
74 for (j = 0; j < SUB_YDIM; j++) {
75 MPI_Get_address(&src_buf[j*XDIM], &idx_loc[j]);
76 idx_loc[j] = idx_loc[j] - base_int;
77 idx_rem[j] = j*XDIM*sizeof(double);
78 blk_len[j] = SUB_XDIM*sizeof(double);
81 MPI_Type_create_hindexed(SUB_YDIM, blk_len, idx_loc, MPI_BYTE, &src_type);
82 MPI_Type_create_indexed_block(SUB_YDIM, SUB_XDIM*sizeof(double), idx_rem, MPI_BYTE, &dst_type);
84 MPI_Type_commit(&src_type);
85 MPI_Type_commit(&dst_type);
87 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, peer, 0, buf_win);
88 MPI_Put(base_ptr, 1, src_type, peer, 0, 1, dst_type, buf_win);
89 MPI_Win_unlock(peer, buf_win);
91 MPI_Type_free(&src_type);
92 MPI_Type_free(&dst_type);
95 MPI_Barrier(MPI_COMM_WORLD);
97 /* Verify that the results are correct */
99 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, rank, 0, buf_win);
101 for (i = 0; i < SUB_XDIM; i++) {
102 for (j = 0; j < SUB_YDIM; j++) {
103 const double actual = *(win_buf + i + j*XDIM);
104 const double expected = (1.0 + ((rank+nranks-1)%nranks));
105 if (actual - expected > 1e-10) {
106 SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
107 rank, j, i, expected, actual); );
113 for (i = SUB_XDIM; i < XDIM; i++) {
114 for (j = 0; j < SUB_YDIM; j++) {
115 const double actual = *(win_buf + i + j*XDIM);
116 const double expected = 1.0 + rank;
117 if (actual - expected > 1e-10) {
118 SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
119 rank, j, i, expected, actual); );
125 for (i = 0; i < XDIM; i++) {
126 for (j = SUB_YDIM; j < YDIM; j++) {
127 const double actual = *(win_buf + i + j*XDIM);
128 const double expected = 1.0 + rank;
129 if (actual - expected > 1e-10) {
130 SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
131 rank, j, i, expected, actual); );
137 MPI_Win_unlock(rank, buf_win);
139 MPI_Win_free(&buf_win);
140 MPI_Free_mem(win_buf);
141 MPI_Free_mem(src_buf);
142 MPI_Free_mem(dst_buf);
144 MTest_Finalize( errors );