1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200
|
/*********************************************************************
*
* Copyright (C) 2014, Northwestern University and Argonne National Laboratory
* See COPYRIGHT notice in top-level directory.
*
*********************************************************************/
/* $Id$ */
/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
* This example makes a number of nonblocking API calls, each writes a single
* column of a 2D integer array. Each process writes NX columns and any two
* consecutive columns are of nprocs columns distance apart from each other.
* In this case, the fileview of each process interleaves with all other
* processes.
* If simply concatenating fileviews of all the nonblocking calls will result
* in a fileview that violates the MPI-IO requirement on the fileview of which
* flattened file offsets must be monotonically non-decreasing. PnetCDF handles
* this case by breaking down each nonblocking call into a list of offset-length
* pairs, merging the pairs across multiple nonblocking calls, and sorting
* them into an increasing order. The sorted pairs are used to construct a
* fileview that meets the monotonically non-decreasing offset requirement,
* and thus the nonblocking requests can be serviced by a single MPI-IO call.
*
* The compile and run commands are given below, together with an ncmpidump of
* the output file.
*
* % mpicxx -O2 -o column_wise column_wise.cpp -lpnetcdf
* % mpiexec -l -n 4 ./column_wise /pvfs2/wkliao/testfile.nc
* 0: 0: myOff= 0 myNX= 4
* 1: 1: myOff= 4 myNX= 4
* 2: 2: myOff= 8 myNX= 4
* 3: 3: myOff= 12 myNX= 4
* 0: 0: start= 0 0 count= 10 1
* 1: 1: start= 0 1 count= 10 1
* 2: 2: start= 0 2 count= 10 1
* 3: 3: start= 0 3 count= 10 1
*
* % ncmpidump /pvfs2/wkliao/testfile.nc
* netcdf testfile {
* // file format: CDF-5 (big variables)
* dimensions:
* Y = 10 ;
* X = 16 ;
* variables:
* int var(Y, X) ;
* data:
*
* var =
* 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3,
* 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3,
* 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3,
* 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3,
* 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3,
* 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3,
* 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3,
* 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3,
* 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3,
* 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3, 0, 1, 2, 3 ;
* }
* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
#include <stdio.h>
#include <stdlib.h>
#include <iostream>
using namespace std;
#include <string.h> /* strcpy(), strncpy() */
#include <unistd.h> /* getopt() */
#include <pnetcdf>
using namespace PnetCDF;
using namespace PnetCDF::exceptions;
#define NY 10
#define NX 4
static void
usage(char *argv0)
{
cerr <<
"Usage: %s [-h] | [-q] [file_name]\n"
" [-h] Print help\n"
" [-q] Quiet mode (reports when fail)\n"
" [filename] output netCDF file name\n"
<< argv0;
}
int main(int argc, char** argv)
{
extern int optind;
char filename[256];
int i, j, verbose=1, rank, nprocs, myNX, G_NX, myOff, num_reqs;
int *reqs, *sts, **buf;
vector<MPI_Offset> start(2), count(2);
MPI_Init(&argc, &argv);
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
/* get command-line arguments */
while ((i = getopt(argc, argv, "hq")) != EOF)
switch(i) {
case 'q': verbose = 0;
break;
case 'h':
default: if (rank==0) usage(argv[0]);
MPI_Finalize();
return 1;
}
if (argv[optind] == NULL) strcpy(filename, "testfile.nc");
else snprintf(filename, 256, "%s", argv[optind]);
try {
NcmpiFile nc(MPI_COMM_WORLD, filename, NcmpiFile::replace,
NcmpiFile::classic5);
/* the global array is NY * (NX * nprocs) */
G_NX = NX * nprocs;
myOff = NX * rank;
myNX = NX;
if (verbose) printf("%2d: myOff=%3d myNX=%3d\n",rank,myOff,myNX);
/* define dimensions Y and X */
vector<NcmpiDim> dimid(2);
dimid[0] = nc.addDim("Y", NY);
dimid[1] = nc.addDim("X", G_NX);
/* define a 2D variable of integer type */
NcmpiVar var = nc.addVar("var", ncmpiInt, dimid);
/* First, fill the entire array with zeros, using a blocking I/O.
Every process writes a subarray of size NY * myNX */
buf = (int**) malloc(sizeof(int*) * myNX);
buf[0] = (int*) calloc(NY * myNX, sizeof(int));
start[0] = 0; start[1] = myOff;
count[0] = NY; count[1] = myNX;
var.putVar_all(start, count, &buf[0][0]);
free(buf[0]);
nc.flush();
/* initialize the buffer with rank ID. Also make the case interesting,
by allocating buffers separately */
for (i=0; i<myNX; i++) {
buf[i] = (int*) malloc(sizeof(int) * NY);
for (j=0; j<NY; j++) buf[i][j] = rank;
}
reqs = (int*) malloc(sizeof(int) * myNX);
sts = (int*) malloc(sizeof(int) * myNX);
/* each proc writes myNX single columns of the 2D array */
start[0] = 0; start[1] = rank;
count[0] = NY; count[1] = 1;
if (verbose)
printf("%2d: start=%3lld %3lld count=%3lld %3lld\n",
rank, start[0],start[1], count[0],count[1]);
num_reqs = 0;
for (i=0; i<myNX; i++) {
var.iputVar(start, count, &buf[0][0], &reqs[num_reqs++]);
start[1] += nprocs;
}
nc.Wait_all(num_reqs, reqs, sts);
/* check status of all requests */
for (i=0; i<num_reqs; i++)
if (sts[i] != NC_NOERR)
printf("Error: nonblocking write fails on request %d (%s)\n",
i, ncmpi_strerror(sts[i]));
free(sts);
free(reqs);
for (i=0; i<myNX; i++) free(buf[i]);
free(buf);
/* file is close implicitly */
}
catch(NcmpiException& e) {
cout << e.what() << " error code=" << e.errorCode() << " Error!\n";
return 1;
}
/* check if there is any PnetCDF internal malloc residue */
MPI_Offset malloc_size, sum_size;
int err = ncmpi_inq_malloc_size(&malloc_size);
if (err == NC_NOERR) {
MPI_Reduce(&malloc_size, &sum_size, 1, MPI_OFFSET, MPI_SUM, 0, MPI_COMM_WORLD);
if (rank == 0 && sum_size > 0)
printf("heap memory allocated by PnetCDF internally has %lld bytes yet to be freed\n",
sum_size);
}
MPI_Finalize();
return 0;
}
|