1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
3 * (C) 2014 by Argonne National Laboratory.
4 * See COPYRIGHT in top-level directory.
15 #define CHECK(fn) {int errcode; errcode = (fn); if (errcode != MPI_SUCCESS) handle_error(errcode, #fn);}
17 static void handle_error(int errcode, const char *str)
19 char msg[MPI_MAX_ERROR_STRING];
21 MPI_Error_string(errcode, msg, &resultlen);
22 fprintf(stderr, "%s: %s\n", str, msg);
23 MPI_Abort(MPI_COMM_WORLD, 1);
27 int main(int argc, char *argv[])
29 int i, j, nerrors = 0, total_errors = 0;
39 /* Define array distribution
40 * A 2x2 block size works with ROMIO, a 3x3 block size breaks it. */
41 int distrib[2] = { MPI_DISTRIBUTE_CYCLIC, MPI_DISTRIBUTE_CYCLIC };
42 int bsize[2] = { NBLOCK, NBLOCK };
43 int gsize[2] = { NSIDE, NSIDE };
44 int psize[2] = { NPROC, NPROC };
46 double data[NSIDE * NSIDE];
47 double *ldata, *pdata;
54 filename = (argc > 1) ? argv[1] : "testfile";
56 MPI_Init(&argc, &argv);
58 MPI_Comm_size(MPI_COMM_WORLD, &size);
59 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
62 CHECK(MPI_Type_create_darray(size, rank, 2, gsize, distrib,
63 bsize, psize, MPI_ORDER_FORTRAN, MPI_DOUBLE, &darray));
64 CHECK(MPI_Type_commit(&darray));
65 CHECK(MPI_Type_size(darray, &tsize));
66 nelem = tsize / sizeof(double);
68 for (i = 0; i < (NSIDE * NSIDE); i++)
72 CHECK(MPI_File_open(MPI_COMM_SELF, filename,
73 MPI_MODE_CREATE | MPI_MODE_WRONLY, MPI_INFO_NULL, &dfile));
74 CHECK(MPI_File_write(dfile, data, NSIDE * NSIDE, MPI_DOUBLE, &status));
75 CHECK(MPI_File_close(&dfile));
77 MPI_Barrier(MPI_COMM_WORLD);
80 ldata = (double *) malloc(tsize);
81 pdata = (double *) malloc(tsize);
83 /* Use Pack to pull out array */
85 CHECK(MPI_Pack(data, 1, darray, pdata, tsize, &bpos, MPI_COMM_WORLD));
87 MPI_Barrier(MPI_COMM_WORLD);
89 /* Read in array from file. */
90 CHECK(MPI_File_open(MPI_COMM_WORLD, filename, MPI_MODE_RDONLY, MPI_INFO_NULL, &mpi_fh));
91 CHECK(MPI_File_set_view(mpi_fh, 0, MPI_DOUBLE, darray, "native", MPI_INFO_NULL));
92 CHECK(MPI_File_iread_all(mpi_fh, ldata, nelem, MPI_DOUBLE, &request));
93 CHECK(MPI_Wait(&request, &status));
94 CHECK(MPI_File_close(&mpi_fh));
96 for (i = 0; i < size; i++) {
98 MPI_Barrier(MPI_COMM_WORLD);
100 printf("=== Rank %i === (%i elements) \nPacked: ", rank, nelem);
101 for (j = 0; j < nelem; j++) {
102 printf("%4.1f ", pdata[j]);
106 for (j = 0; j < nelem; j++) {
107 printf("%4.1f ", ldata[j]);
115 for (j = 0; j < nelem; j++) {
116 if (pdata[j] != ldata[j]) {
117 fprintf(stderr, "rank %d at index %d: packbuf %4.1f filebuf %4.1f\n",
118 rank, j, pdata[j], ldata[j]);
124 MPI_Allreduce(&nerrors, &total_errors, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
125 if (rank == 0 && total_errors == 0)
126 printf(" No Errors\n");
130 MPI_Type_free(&darray);