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 Get Test
9 * Author: James Dinan <dinan@mcs.anl.gov>
10 * Date : December, 2010
12 * This code performs N strided get operations from 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 indexed type.
30 int main(int argc, char **argv) {
31 int i, j, rank, nranks, peer, bufsize, errors;
32 double *win_buf, *loc_buf;
35 int idx_rem[SUB_YDIM];
36 int blk_len[SUB_YDIM];
37 MPI_Datatype loc_type, rem_type;
39 MTest_Init(&argc, &argv);
41 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
42 MPI_Comm_size(MPI_COMM_WORLD, &nranks);
44 bufsize = XDIM * YDIM * sizeof(double);
45 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &win_buf);
46 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &loc_buf);
48 for (i = 0; i < XDIM*YDIM; i++) {
49 *(win_buf + i) = 1.0 + rank;
50 *(loc_buf + i) = -1.0;
53 MPI_Win_create(win_buf, bufsize, 1, MPI_INFO_NULL, MPI_COMM_WORLD, &buf_win);
55 peer = (rank+1) % nranks;
57 /* Build the datatype */
59 for (i = 0; i < SUB_YDIM; i++) {
61 blk_len[i] = SUB_XDIM;
64 MPI_Type_indexed(SUB_YDIM, blk_len, idx_rem, MPI_DOUBLE, &loc_type);
65 MPI_Type_indexed(SUB_YDIM, blk_len, idx_rem, MPI_DOUBLE, &rem_type);
67 MPI_Type_commit(&loc_type);
68 MPI_Type_commit(&rem_type);
70 /* Perform get operation */
72 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, peer, 0, buf_win);
74 MPI_Get(loc_buf, 1, loc_type, peer, 0, 1, rem_type, buf_win);
76 /* Use the datatype only on the remote side (must have SUB_XDIM == XDIM) */
77 /* MPI_Get(loc_buf, SUB_XDIM*SUB_YDIM, MPI_DOUBLE, peer, 0, 1, rem_type, buf_win); */
79 MPI_Win_unlock(peer, buf_win);
81 MPI_Type_free(&loc_type);
82 MPI_Type_free(&rem_type);
84 MPI_Barrier(MPI_COMM_WORLD);
86 /* Verify that the results are correct */
89 for (i = 0; i < SUB_XDIM; i++) {
90 for (j = 0; j < SUB_YDIM; j++) {
91 const double actual = *(loc_buf + i + j*XDIM);
92 const double expected = (1.0 + peer);
93 if (fabs(actual - expected) > 1.0e-10) {
94 SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
95 rank, j, i, expected, actual); );
101 for (i = SUB_XDIM; i < XDIM; i++) {
102 for (j = 0; j < SUB_YDIM; j++) {
103 const double actual = *(loc_buf + i + j*XDIM);
104 const double expected = -1.0;
105 if (fabs(actual - expected) > 1.0e-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 = 0; i < XDIM; i++) {
114 for (j = SUB_YDIM; j < YDIM; j++) {
115 const double actual = *(loc_buf + i + j*XDIM);
116 const double expected = -1.0;
117 if (fabs(actual - expected) > 1.0e-10) {
118 SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
119 rank, j, i, expected, actual); );
126 MPI_Win_free(&buf_win);
127 MPI_Free_mem(win_buf);
128 MPI_Free_mem(loc_buf);
130 MTest_Finalize( errors );
132 return MTestReturnValue( errors );