1 /*********************************************************************
2  *
3  *  Copyright (C) 2014, Northwestern University and Argonne National Laboratory
4  *  See COPYRIGHT notice in top-level directory.
5  *
6  *********************************************************************/
7 /* $Id: bput_varn_int64.c 2717 2016-12-18 01:20:47Z wkliao $ */
8 
9 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
10  * This example tests nonblocking buffered write varn APIs, including
11  * ncmpi_bput_varn_longlong() and ncmpi_bput_varn(),
12  * It first writes a sequence of requests with arbitrary array indices and
13  * lengths to four variables of type NC_INT64, and reads back.
14  *
15  * The compile and run commands are given below, together with an ncmpidump of
16  * the output file.
17  *
18  *    % mpicc -O2 -o i_varn_int64 i_varn_int64.c -lpnetcdf
19  *    % mpiexec -n 4 ./i_varn_int64 /pvfs2/wkliao/testfile.nc
20  *    % ncmpidump /pvfs2/wkliao/testfile.nc
21  *    netcdf testfile {
22  *    // file format: CDF-5 (big variables)
23  *    dimensions:
24  *             Y = 10 ;
25  *             X = 4 ;
26  *    variables:
27  *            int64 var0(Y, X) ;
28  *            int64 var1(Y, X) ;
29  *            int64 var2(Y, X) ;
30  *            int64 var3(Y, X) ;
31  *    data:
32  *
33  *     var0 =
34  *      1, 0, 3, 0,
35  *      1, 2, 3, 0,
36  *      1, 2, 2, 0,
37  *      3, 2, 1, 2,
38  *      3, 1, 1, 3,
39  *      0, 3, 1, 3,
40  *      0, 3, 0, 3,
41  *      2, 2, 0, 1,
42  *      3, 2, 3, 1,
43  *      3, 2, 3, 1 ;
44  *
45  *     var1 =
46  *      2, 1, 0, 1,
47  *      2, 3, 0, 1,
48  *      2, 3, 3, 1,
49  *      0, 3, 2, 3,
50  *      0, 2, 2, 0,
51  *      1, 0, 2, 0,
52  *      1, 0, 1, 0,
53  *      3, 3, 1, 2,
54  *      0, 3, 0, 2,
55  *      0, 3, 0, 2 ;
56  *
57  *     var2 =
58  *      3, 2, 1, 2,
59  *      3, 0, 1, 2,
60  *      3, 0, 0, 2,
61  *      1, 0, 3, 0,
62  *      1, 3, 3, 1,
63  *      2, 1, 3, 1,
64  *      2, 1, 2, 1,
65  *      0, 0, 2, 3,
66  *      1, 0, 1, 3,
67  *      1, 0, 1, 3 ;
68  *
69  *     var3 =
70  *      0, 3, 2, 3,
71  *      0, 1, 2, 3,
72  *      0, 1, 1, 3,
73  *      2, 1, 0, 1,
74  *      2, 0, 0, 2,
75  *      3, 2, 0, 2,
76  *      3, 2, 3, 2,
77  *      1, 1, 3, 0,
78  *      2, 1, 2, 0,
79  *      2, 1, 2, 0 ;
80  *    }
81  * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
82 
83 #include <stdio.h>
84 #include <stdlib.h>
85 #include <string.h> /* strcpy(), strncpy() */
86 #include <unistd.h> /* getopt() */
87 #include <mpi.h>
88 #include <pnetcdf.h>
89 
90 #ifndef MPI_OFFSET
91 #define MPI_OFFSET MPI_LONG_LONG_INT
92 #endif
93 
94 #define NY 10
95 #define NX 4
96 #define NDIMS 2
97 
98 #define ERR {if(err!=NC_NOERR){printf("Error at line=%d: %s\n", __LINE__, ncmpi_strerror(err));nerrs++;}}
99 
100 #define ERRS(n,a) { \
101     int _i; \
102     for (_i=0; _i<(n); _i++) { \
103         if ((a)[_i] != NC_NOERR) { \
104             printf("Error at line=%d: err[%d] %s\n", __LINE__, _i, \
105                    ncmpi_strerror((a)[_i])); \
106                    nerrs++; \
107         } \
108     } \
109 }
110 
calloc_3D(size_t z,size_t y,size_t x)111 static MPI_Offset *** calloc_3D(size_t z, size_t y, size_t x)
112 {
113     if (z*y*x == 0) return NULL;
114     int _j, _k;
115     MPI_Offset ***buf  = (MPI_Offset***) malloc(z     * sizeof(MPI_Offset**));
116     MPI_Offset  **bufy = (MPI_Offset**)  malloc(z*y   * sizeof(MPI_Offset*));
117     MPI_Offset   *bufx = (MPI_Offset*)   calloc(z*y*x,  sizeof(MPI_Offset));
118     for (_k=0; _k<z; _k++, bufy+=y) {
119         buf[_k] = bufy;
120         for (_j=0; _j<y; _j++, bufx+=x)
121             buf[_k][_j] = bufx;
122     }
123     return buf;
124 }
125 
free_3D(MPI_Offset *** buf)126 static void free_3D(MPI_Offset ***buf)
127 {
128     free(buf[0][0]);
129     free(buf[0]);
130     free(buf);
131 }
132 
133 static void
usage(char * argv0)134 usage(char *argv0)
135 {
136     char *help =
137     "Usage: %s [-h] | [-q] [file_name]\n"
138     "       [-h] Print help\n"
139     "       [-q] Quiet mode (reports when fail)\n"
140     "       [filename] output netCDF file name\n";
141     fprintf(stderr, help, argv0);
142 }
143 
check_contents(int ncid,int * varid)144 static int check_contents(int ncid, int *varid)
145 {
146     int i, j, err, nerrs=0;
147     long long expected[4][NY*NX] = {{1, 0, 3, 0, 1, 2, 3, 0, 1, 2,
148                                      2, 0, 3, 2, 1, 2, 3, 1, 1, 3,
149                                      0, 3, 1, 3, 0, 3, 0, 3, 2, 2,
150                                      0, 1, 3, 2, 3, 1, 3, 2, 3, 1},
151                                     {2, 1, 0, 1, 2, 3, 0, 1, 2, 3,
152                                      3, 1, 0, 3, 2, 3, 0, 2, 2, 0,
153                                      1, 0, 2, 0, 1, 0, 1, 0, 3, 3,
154                                      1, 2, 0, 3, 0, 2, 0, 3, 0, 2},
155                                     {3, 2, 1, 2, 3, 0, 1, 2, 3, 0,
156                                      0, 2, 1, 0, 3, 0, 1, 3, 3, 1,
157                                      2, 1, 3, 1, 2, 1, 2, 1, 0, 0,
158                                      2, 3, 1, 0, 1, 3, 1, 0, 1, 3},
159                                     {0, 3, 2, 3, 0, 1, 2, 3, 0, 1,
160                                      1, 3, 2, 1, 0, 1, 2, 0, 0, 2,
161                                      3, 2, 0, 2, 3, 2, 3, 2, 1, 1,
162                                      3, 0, 2, 1, 2, 0, 2, 1, 2, 0}};
163 
164     long long *r_buffer = (long long*) malloc(NY*NX*sizeof(long long));
165 
166     for (i=0; i<4; i++) {
167         for (j=0; j<NY*NX; j++) r_buffer[j] = 99999;
168         err = ncmpi_get_var_longlong_all(ncid, varid[i], r_buffer);
169         ERR
170 
171         /* check if the contents of buf are expected */
172         for (j=0; j<NY*NX; j++)
173             if (r_buffer[j] != expected[i][j]) {
174                 printf("Expected file contents [%d][%d]=%lld, but got %lld\n",
175                        i,j,expected[i][j],r_buffer[j]);
176                 nerrs++;
177             }
178     }
179     free(r_buffer);
180     return nerrs;
181 }
182 
main(int argc,char ** argv)183 int main(int argc, char** argv)
184 {
185     extern int optind;
186     char filename[256], *exec;
187     int i, j, k, n, rank, nprocs, verbose=1, err, nerrs=0;
188     int ncid, cmode, varid[4], dimid[2], nreqs, reqs[4], sts[4];
189     long long *buffer[4];
190     int num_segs[4], req_lens[4];
191     MPI_Offset ***starts, ***counts;
192     MPI_Info info;
193 
194     MPI_Init(&argc, &argv);
195     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
196     MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
197     exec = argv[0];
198 
199     /* get command-line arguments */
200     while ((i = getopt(argc, argv, "hq")) != EOF)
201         switch(i) {
202             case 'q': verbose = 0;
203                       break;
204             case 'h':
205             default:  if (rank==0) usage(argv[0]);
206                       MPI_Finalize();
207                       return 0;
208         }
209     argc -= optind;
210     argv += optind;
211     if (argc == 1) snprintf(filename, 256, "%s", argv[0]);
212     else           strcpy (filename, "testfile.nc");
213 
214     if (nprocs != 4 && rank == 0 && verbose)
215         printf("Warning: %s is intended to run on 4 processes\n",exec);
216 
217     /* set an MPI-IO hint to disable file offset alignment for fixed-size
218      * variables */
219     MPI_Info_create(&info);
220     MPI_Info_set(info, "nc_var_align_size", "1");
221 
222     /* create a new file for writing ----------------------------------------*/
223     cmode = NC_CLOBBER | NC_64BIT_DATA;
224     err = ncmpi_create(MPI_COMM_WORLD, filename, cmode, info, &ncid);
225     ERR
226 
227     MPI_Info_free(&info);
228 
229     /* create a global array of size NY * NX */
230     err = ncmpi_def_dim(ncid, "Y", NY, &dimid[0]); ERR
231     err = ncmpi_def_dim(ncid, "X", NX, &dimid[1]); ERR
232     err = ncmpi_def_var(ncid, "var0", NC_INT64, NDIMS, dimid, &varid[0]); ERR
233     err = ncmpi_def_var(ncid, "var1", NC_INT64, NDIMS, dimid, &varid[1]); ERR
234     err = ncmpi_def_var(ncid, "var2", NC_INT64, NDIMS, dimid, &varid[2]); ERR
235     err = ncmpi_def_var(ncid, "var3", NC_INT64, NDIMS, dimid, &varid[3]); ERR
236     err = ncmpi_enddef(ncid); ERR
237 
238     /* allocate space for starts and counts */
239     starts = calloc_3D(4, 13, NDIMS);
240     counts = calloc_3D(4, 13, NDIMS);
241 
242     n = rank % 4;
243     num_segs[n] = 8; /* number of segments for this request */
244     starts[n][0][0]=0; starts[n][0][1]=1; counts[n][0][0]=1; counts[n][0][1]=1;
245     starts[n][1][0]=0; starts[n][1][1]=3; counts[n][1][0]=1; counts[n][1][1]=1;
246     starts[n][2][0]=1; starts[n][2][1]=3; counts[n][2][0]=1; counts[n][2][1]=1;
247     starts[n][3][0]=2; starts[n][3][1]=3; counts[n][3][0]=1; counts[n][3][1]=1;
248     starts[n][4][0]=5; starts[n][4][1]=0; counts[n][4][0]=1; counts[n][4][1]=1;
249     starts[n][5][0]=6; starts[n][5][1]=0; counts[n][5][0]=1; counts[n][5][1]=1;
250     starts[n][6][0]=6; starts[n][6][1]=2; counts[n][6][0]=1; counts[n][6][1]=1;
251     starts[n][7][0]=7; starts[n][7][1]=2; counts[n][7][0]=1; counts[n][7][1]=1;
252     /* starts[n][][] n_counts[n][][] indicate the following: ("-" means skip)
253               _  X  _  X
254               _  _  _  X
255               _  _  _  X
256               _  _  _  _
257               _  _  _  _
258               X  _  _  _
259               X  _  X  _
260               _  _  X  _
261               _  _  _  _
262               _  _  _  _
263      */
264     n = (rank+1) % 4;
265     num_segs[n] = 13; /* number of segments for this request */
266     starts[n][ 0][0]=0;starts[n][ 0][1]=2;counts[n][ 0][0]=1;counts[n][ 0][1]=1;
267     starts[n][ 1][0]=1;starts[n][ 1][1]=2;counts[n][ 1][0]=1;counts[n][ 1][1]=1;
268     starts[n][ 2][0]=3;starts[n][ 2][1]=0;counts[n][ 2][0]=1;counts[n][ 2][1]=1;
269     starts[n][ 3][0]=4;starts[n][ 3][1]=0;counts[n][ 3][0]=1;counts[n][ 3][1]=1;
270     starts[n][ 4][0]=4;starts[n][ 4][1]=3;counts[n][ 4][0]=1;counts[n][ 4][1]=1;
271     starts[n][ 5][0]=5;starts[n][ 5][1]=1;counts[n][ 5][0]=1;counts[n][ 5][1]=1;
272     starts[n][ 6][0]=5;starts[n][ 6][1]=3;counts[n][ 6][0]=1;counts[n][ 6][1]=1;
273     starts[n][ 7][0]=6;starts[n][ 7][1]=1;counts[n][ 7][0]=1;counts[n][ 7][1]=1;
274     starts[n][ 8][0]=6;starts[n][ 8][1]=3;counts[n][ 8][0]=1;counts[n][ 8][1]=1;
275     starts[n][ 9][0]=8;starts[n][ 9][1]=0;counts[n][ 9][0]=1;counts[n][ 9][1]=1;
276     starts[n][10][0]=8;starts[n][10][1]=2;counts[n][10][0]=1;counts[n][10][1]=1;
277     starts[n][11][0]=9;starts[n][11][1]=0;counts[n][11][0]=1;counts[n][11][1]=1;
278     starts[n][12][0]=9;starts[n][12][1]=2;counts[n][12][0]=1;counts[n][12][1]=1;
279     /* starts[n][][] counts[n][][] indicate the following pattern.
280               _  _  X  _
281               _  _  X  _
282               _  _  _  _
283               X  _  _  _
284               X  _  _  X
285               _  X  _  X
286               _  X  _  X
287               _  _  _  _
288               X  _  X  _
289               X  _  X  _
290      */
291     n = (rank+2) % 4;
292     num_segs[n] = 7; /* number of segments for this request */
293     starts[n][0][0]=1; starts[n][0][1]=1; counts[n][0][0]=1; counts[n][0][1]=1;
294     starts[n][1][0]=2; starts[n][1][1]=1; counts[n][1][0]=1; counts[n][1][1]=2;
295     starts[n][2][0]=3; starts[n][2][1]=1; counts[n][2][0]=1; counts[n][2][1]=1;
296     starts[n][3][0]=3; starts[n][3][1]=3; counts[n][3][0]=1; counts[n][3][1]=1;
297     starts[n][4][0]=7; starts[n][4][1]=0; counts[n][4][0]=1; counts[n][4][1]=2;
298     starts[n][5][0]=8; starts[n][5][1]=1; counts[n][5][0]=1; counts[n][5][1]=1;
299     starts[n][6][0]=9; starts[n][6][1]=1; counts[n][6][0]=1; counts[n][6][1]=1;
300     /* starts[n][][] counts[n][][] indicate the following pattern.
301               _  _  _  _
302               _  X  _  _
303               _  X  X  _
304               _  X  _  X
305               _  _  _  _
306               _  _  _  _
307               _  _  _  _
308               X  X  _  _
309               _  X  _  _
310               _  X  _  _
311      */
312     n = (rank+3) % 4;
313     num_segs[n] = 10; /* number of segments for this request */
314     starts[n][0][0]=0; starts[n][0][1]=0; counts[n][0][0]=1; counts[n][0][1]=1;
315     starts[n][1][0]=1; starts[n][1][1]=0; counts[n][1][0]=1; counts[n][1][1]=1;
316     starts[n][2][0]=2; starts[n][2][1]=0; counts[n][2][0]=1; counts[n][2][1]=1;
317     starts[n][3][0]=3; starts[n][3][1]=2; counts[n][3][0]=1; counts[n][3][1]=1;
318     starts[n][4][0]=4; starts[n][4][1]=1; counts[n][4][0]=1; counts[n][4][1]=2;
319     starts[n][5][0]=5; starts[n][5][1]=2; counts[n][5][0]=1; counts[n][5][1]=1;
320     starts[n][6][0]=7; starts[n][6][1]=3; counts[n][6][0]=1; counts[n][6][1]=1;
321     starts[n][7][0]=8; starts[n][7][1]=3; counts[n][7][0]=1; counts[n][7][1]=1;
322     starts[n][8][0]=9; starts[n][8][1]=3; counts[n][8][0]=1; counts[n][8][1]=1;
323      /*starts[n][][] counts[n][][] indicate the following pattern.
324               X  _  _  _
325               X  _  _  _
326               X  _  _  _
327               _  _  X  _
328               _  X  X  _
329               _  _  X  _
330               _  _  _  _
331               _  _  _  X
332               _  _  _  X
333               _  _  _  X
334      */
335 
336     /* only rank 0, 1, 2, and 3 do I/O:
337      * each of ranks 0 to 3 write 4 nonblocking requests */
338     nreqs = 4;
339     if (rank >= 4) nreqs = 0;
340 
341     /* bufsize must be max of data type converted before and after */
342     MPI_Offset bufsize = 0;
343 
344     /* calculate length of each varn request, number of segments in each
345      * varn request, and allocate write buffer */
346     for (i=0; i<nreqs; i++) {
347         req_lens[i] = 0; /* total length this request */
348         for (j=0; j<num_segs[i]; j++) {
349             MPI_Offset req_len=1;
350             for (k=0; k<NDIMS; k++)
351                 req_len *= counts[i][j][k];
352             req_lens[i] += req_len;
353         }
354         if (verbose) printf("req_lens[%d]=%d\n",i,req_lens[i]);
355 
356         /* allocate I/O buffer and initialize its contents */
357         buffer[i] = (long long*) malloc(req_lens[i] * sizeof(long long));
358         for (j=0; j<req_lens[i]; j++) buffer[i][j] = rank;
359 
360         bufsize += req_lens[i];
361     }
362     bufsize *= sizeof(long long);
363     if (verbose) printf("%d: Attach buffer size %lld\n", rank, bufsize);
364 
365     /* give PnetCDF a space to buffer the nonblocking requests */
366     if (bufsize > 0) {
367         err = ncmpi_buffer_attach(ncid, bufsize); ERR
368     }
369 
370     /* write using varn API */
371     for (i=0; i<nreqs; i++) {
372         err = ncmpi_bput_varn_longlong(ncid, varid[i], num_segs[i], starts[i],
373                                        counts[i], buffer[i], &reqs[i]);
374         ERR
375     }
376     err = ncmpi_wait_all(ncid, nreqs, reqs, sts);
377     ERRS(nreqs, sts)
378 
379     /* check file contents */
380     if (nprocs >= 4) nerrs += check_contents(ncid, varid);
381 
382     for (i=0; i<nreqs; i++) free(buffer[i]);
383 
384     /* test flexible put API, using a noncontiguous buftype */
385     for (i=0; i<nreqs; i++) {
386         MPI_Datatype buftype;
387         MPI_Type_vector(req_lens[i], 1, 2, MPI_UNSIGNED, &buftype);
388         MPI_Type_commit(&buftype);
389         buffer[i] = (long long*) malloc(req_lens[i] * 2 * sizeof(long long));
390         for (j=0; j<req_lens[i]*2; j++) buffer[i][j] = rank;
391 
392         err = ncmpi_bput_varn(ncid, varid[i], num_segs[i], starts[i],
393                               counts[i], buffer[i], 1, buftype, &reqs[i]);
394         ERR
395         MPI_Type_free(&buftype);
396     }
397     err = ncmpi_wait_all(ncid, nreqs, reqs, sts);
398     ERRS(nreqs, sts)
399 
400     /* check file contents */
401     if (nprocs >= 4) nerrs += check_contents(ncid, varid);
402 
403     /* free the buffer space for bput */
404     if (bufsize > 0) {
405         err = ncmpi_buffer_detach(ncid); ERR
406     }
407 
408     err = ncmpi_close(ncid);
409     ERR
410 
411     for (i=0; i<nreqs; i++) free(buffer[i]);
412     free_3D(starts);
413     free_3D(counts);
414 
415     /* check if there is any PnetCDF internal malloc residue */
416     MPI_Offset malloc_size, sum_size;
417     err = ncmpi_inq_malloc_size(&malloc_size);
418     if (err == NC_NOERR) {
419         MPI_Reduce(&malloc_size, &sum_size, 1, MPI_OFFSET, MPI_SUM, 0, MPI_COMM_WORLD);
420         if (rank == 0 && sum_size > 0)
421             printf("heap memory allocated by PnetCDF internally has %lld bytes yet to be freed\n",
422                    sum_size);
423     }
424 
425     MPI_Finalize();
426     return nerrs;
427 }
428 
429