1 /*
2   Copyright 2007, UCAR/Unidata
3   See COPYRIGHT file for copying and redistribution conditions.
4 
5   This is part of netCDF.
6 
7   This program tests for a bug discovered with nofill mode that failed
8   only on file systems with block size in a particular range.  It fails
9   when invoked with the blksize argument between 2091953 and 2150032,
10   inclusive, and succeeds for other blksizes.
11 
12   $Id: tst_nofill.c 2744 2016-12-28 16:25:22Z wkliao $
13 */
14 
15 #include <stdio.h>
16 #include <string.h>
17 #include <libgen.h> /* basename() */
18 #include <limits.h>
19 #include <stdlib.h>
20 #include <errno.h>
21 #include <pnetcdf.h>
22 
23 #include <testutils.h>
24 
25 #define ERR {if (err != NC_NOERR) {printf("Error at %s line %d: %s\n",__func__,__LINE__,ncmpi_strerror(err)); return 1;}}
26 
27 static void
check_err(const int stat,const int line,const char * file)28 check_err(const int stat, const int line, const char *file) {
29     if (stat != NC_NOERR) {
30 	(void) fprintf(stderr, "line %d of %s: %s\n", line, file, ncmpi_strerror(stat));
31         fflush(stderr);
32 	exit(1);
33     }
34 }
35 
36 #define LON_LEN 240
37 #define LAT_LEN 121
38 #define LVL_LEN 31
39 #define TIME_LEN 1
40 
41 static int
create_file(char * file_name,int fill_mode)42 create_file(char *file_name, int fill_mode)
43 {
44    int i;
45    int  stat;			/* return status */
46    int  ncid;			/* netCDF id */
47 
48    /* dimension ids */
49    int lon_dim;
50    int lat_dim;
51    int lvl_dim;
52    int time_dim;
53 
54    /* dimension lengths */
55    MPI_Offset lon_len = LON_LEN;
56    MPI_Offset lat_len = LAT_LEN;
57    MPI_Offset lvl_len = LVL_LEN;
58    MPI_Offset time_len = TIME_LEN;
59 
60    /* variable ids */
61    int time_id;
62    int lat_id;
63    int lon_id;
64    int lvl_id;
65    int sfc_pres_id;
66    int temp_scrn_id;
67    int qsair_scrn_id;
68    int topog_id;
69    int mslp_id;
70    int sfc_temp_id;
71    int zonal_wnd_id;
72 
73    /* rank (number of dimensions) for each variable */
74 #  define RANK_time 1
75 #  define RANK_lat 1
76 #  define RANK_lon 1
77 #  define RANK_lvl 1
78 #  define RANK_sfc_pres 3
79 #  define RANK_temp_scrn 3
80 #  define RANK_qsair_scrn 3
81 #  define RANK_topog 3
82 #  define RANK_mslp 3
83 #  define RANK_sfc_temp 3
84 #  define RANK_zonal_wnd 4
85 
86    /* variable shapes */
87    int time_dims[RANK_time];
88    int lat_dims[RANK_lat];
89    int lon_dims[RANK_lon];
90    int lvl_dims[RANK_lvl];
91    int sfc_pres_dims[RANK_sfc_pres];
92    int temp_scrn_dims[RANK_temp_scrn];
93    int qsair_scrn_dims[RANK_qsair_scrn];
94    int topog_dims[RANK_topog];
95    int mslp_dims[RANK_mslp];
96    int sfc_temp_dims[RANK_sfc_temp];
97    int zonal_wnd_dims[RANK_zonal_wnd];
98 
99    MPI_Offset zonal_wnd_start[RANK_zonal_wnd];
100    MPI_Offset zonal_wnd_count[RANK_zonal_wnd];
101    float zonal_wnd[LON_LEN*LAT_LEN*TIME_LEN];
102    int ii;
103 
104    int old_fill_mode;
105 
106    stat = ncmpi_create(MPI_COMM_WORLD, file_name,NC_CLOBBER, MPI_INFO_NULL, &ncid);
107 
108    check_err(stat,__LINE__,__FILE__);
109    stat = ncmpi_set_fill(ncid, fill_mode, &old_fill_mode);
110    check_err(stat,__LINE__,__FILE__);
111 
112    /* define dimensions */
113    stat = ncmpi_def_dim(ncid, "lon", lon_len, &lon_dim);
114    check_err(stat,__LINE__,__FILE__);
115    stat = ncmpi_def_dim(ncid, "lat", lat_len, &lat_dim);
116    check_err(stat,__LINE__,__FILE__);
117    stat = ncmpi_def_dim(ncid, "lvl", lvl_len, &lvl_dim);
118    check_err(stat,__LINE__,__FILE__);
119    stat = ncmpi_def_dim(ncid, "time", time_len, &time_dim);
120    check_err(stat,__LINE__,__FILE__);
121 
122    /* define variables */
123    time_dims[0] = time_dim;
124    stat = ncmpi_def_var(ncid, "time", NC_DOUBLE, RANK_time, time_dims, &time_id);
125    check_err(stat,__LINE__,__FILE__);
126 
127    lat_dims[0] = lat_dim;
128    stat = ncmpi_def_var(ncid, "lat", NC_FLOAT, RANK_lat, lat_dims, &lat_id);
129    check_err(stat,__LINE__,__FILE__);
130 
131    lon_dims[0] = lon_dim;
132    stat = ncmpi_def_var(ncid, "lon", NC_FLOAT, RANK_lon, lon_dims, &lon_id);
133    check_err(stat,__LINE__,__FILE__);
134 
135    lvl_dims[0] = lvl_dim;
136    stat = ncmpi_def_var(ncid, "lvl", NC_FLOAT, RANK_lvl, lvl_dims, &lvl_id);
137    check_err(stat,__LINE__,__FILE__);
138 
139    sfc_pres_dims[0] = time_dim;
140    sfc_pres_dims[1] = lat_dim;
141    sfc_pres_dims[2] = lon_dim;
142    stat = ncmpi_def_var(ncid, "sfc_pres", NC_FLOAT, RANK_sfc_pres, sfc_pres_dims, &sfc_pres_id);
143    check_err(stat,__LINE__,__FILE__);
144 
145    temp_scrn_dims[0] = time_dim;
146    temp_scrn_dims[1] = lat_dim;
147    temp_scrn_dims[2] = lon_dim;
148    stat = ncmpi_def_var(ncid, "temp_scrn", NC_FLOAT, RANK_temp_scrn, temp_scrn_dims, &temp_scrn_id);
149    check_err(stat,__LINE__,__FILE__);
150 
151    qsair_scrn_dims[0] = time_dim;
152    qsair_scrn_dims[1] = lat_dim;
153    qsair_scrn_dims[2] = lon_dim;
154    stat = ncmpi_def_var(ncid, "qsair_scrn", NC_FLOAT, RANK_qsair_scrn, qsair_scrn_dims, &qsair_scrn_id);
155    check_err(stat,__LINE__,__FILE__);
156 
157    topog_dims[0] = time_dim;
158    topog_dims[1] = lat_dim;
159    topog_dims[2] = lon_dim;
160    stat = ncmpi_def_var(ncid, "topog", NC_FLOAT, RANK_topog, topog_dims, &topog_id);
161    check_err(stat,__LINE__,__FILE__);
162 
163    mslp_dims[0] = time_dim;
164    mslp_dims[1] = lat_dim;
165    mslp_dims[2] = lon_dim;
166    stat = ncmpi_def_var(ncid, "mslp", NC_FLOAT, RANK_mslp, mslp_dims, &mslp_id);
167    check_err(stat,__LINE__,__FILE__);
168 
169    sfc_temp_dims[0] = time_dim;
170    sfc_temp_dims[1] = lat_dim;
171    sfc_temp_dims[2] = lon_dim;
172    stat = ncmpi_def_var(ncid, "sfc_temp", NC_FLOAT, RANK_sfc_temp, sfc_temp_dims, &sfc_temp_id);
173    check_err(stat,__LINE__,__FILE__);
174 
175    zonal_wnd_dims[0] = time_dim;
176    zonal_wnd_dims[1] = lvl_dim;
177    zonal_wnd_dims[2] = lat_dim;
178    zonal_wnd_dims[3] = lon_dim;
179    stat = ncmpi_def_var(ncid, "zonal_wnd", NC_FLOAT, RANK_zonal_wnd, zonal_wnd_dims, &zonal_wnd_id);
180    check_err(stat,__LINE__,__FILE__);
181 
182    /* leave define mode */
183    stat = ncmpi_enddef (ncid);
184    check_err(stat,__LINE__,__FILE__);
185 
186    {				/* store time */
187        MPI_Offset time_start[RANK_time];
188        MPI_Offset time_count[RANK_time];
189        double time[TIME_LEN] = {1.};
190        time_len = 1;
191        time_start[0] = 0;
192        time_count[0] = time_len;
193        stat = ncmpi_put_vara_double_all(ncid, time_id, time_start, time_count, time);
194        check_err(stat,__LINE__,__FILE__);
195    }
196 
197    {				/* store lat */
198        float lat[] = {90, 88.5, 87, 85.5, 84, 82.5, 81, 79.5, 78, 76.5, 75, 73.5, 72, 70.5, 69, 67.5, 66, 64.5, 63, 61.5, 60, 58.5, 57, 55.5, 54, 52.5, 51, 49.5, 48, 46.5, 45, 43.5, 42, 40.5, 39, 37.5, 36, 34.5, 33, 31.5, 30, 28.5, 27, 25.5, 24, 22.5, 21, 19.5, 18, 16.5, 15, 13.5, 12, 10.5, 9, 7.5, 6, 4.5, 3, 1.5, 0, -1.5, -3, -4.5, -6, -7.5, -9, -10.5, -12, -13.5, -15, -16.5, -18, -19.5, -21, -22.5, -24, -25.5, -27, -28.5, -30, -31.5, -33, -34.5, -36, -37.5, -39, -40.5, -42, -43.5, -45, -46.5, -48, -49.5, -51, -52.5, -54, -55.5, -57, -58.5, -60, -61.5, -63, -64.5, -66, -67.5, -69, -70.5, -72, -73.5, -75, -76.5, -78, -79.5, -81, -82.5, -84, -85.5, -87, -88.5, -90};
199        stat = ncmpi_put_var_float_all(ncid, lat_id, lat);
200        check_err(stat,__LINE__,__FILE__);
201    }
202 
203    {				/* store lon */
204        float lon[] = {0, 1.5, 3, 4.5, 6, 7.5, 9, 10.5, 12, 13.5, 15, 16.5, 18, 19.5, 21, 22.5, 24, 25.5, 27, 28.5, 30, 31.5, 33, 34.5, 36, 37.5, 39, 40.5, 42, 43.5, 45, 46.5, 48, 49.5, 51, 52.5, 54, 55.5, 57, 58.5, 60, 61.5, 63, 64.5, 66, 67.5, 69, 70.5, 72, 73.5, 75, 76.5, 78, 79.5, 81, 82.5, 84, 85.5, 87, 88.5, 90, 91.5, 93, 94.5, 96, 97.5, 99, 100.5, 102, 103.5, 105, 106.5, 108, 109.5, 111, 112.5, 114, 115.5, 117, 118.5, 120, 121.5, 123, 124.5, 126, 127.5, 129, 130.5, 132, 133.5, 135, 136.5, 138, 139.5, 141, 142.5, 144, 145.5, 147, 148.5, 150, 151.5, 153, 154.5, 156, 157.5, 159, 160.5, 162, 163.5, 165, 166.5, 168, 169.5, 171, 172.5, 174, 175.5, 177, 178.5, 180, 181.5, 183, 184.5, 186, 187.5, 189, 190.5, 192, 193.5, 195, 196.5, 198, 199.5, 201, 202.5, 204, 205.5, 207, 208.5, 210, 211.5, 213, 214.5, 216, 217.5, 219, 220.5, 222, 223.5, 225, 226.5, 228, 229.5, 231, 232.5, 234, 235.5, 237, 238.5, 240, 241.5, 243, 244.5, 246, 247.5, 249, 250.5, 252, 253.5, 255, 256.5, 258, 259.5, 261, 262.5, 264, 265.5, 267, 268.5, 270, 271.5, 273, 274.5, 276, 277.5, 279, 280.5, 282, 283.5, 285, 286.5, 288, 289.5, 291, 292.5, 294, 295.5, 297, 298.5, 300, 301.5, 303, 304.5, 306, 307.5, 309, 310.5, 312, 313.5, 315, 316.5, 318, 319.5, 321, 322.5, 324, 325.5, 327, 328.5, 330, 331.5, 333, 334.5, 336, 337.5, 339, 340.5, 342, 343.5, 345, 346.5, 348, 349.5, 351, 352.5, 354, 355.5, 357, 358.5};
205        stat = ncmpi_put_var_float_all(ncid, lon_id, lon);
206        check_err(stat,__LINE__,__FILE__);
207    }
208 
209    {				/* store lvl */
210        float lvl[] = {1000, 995, 990, 985, 975, 950, 925, 900, 875, 850, 800, 750, 700, 600, 500, 450, 400, 350, 300, 275, 250, 225, 200, 175, 150, 100, 70, 50, 30, 20, 10};
211        stat = ncmpi_put_var_float_all(ncid, lvl_id, lvl);
212        check_err(stat,__LINE__,__FILE__);
213    }
214 
215    {				/* store sfc_pres */
216        MPI_Offset sfc_pres_start[RANK_sfc_pres];
217        MPI_Offset sfc_pres_count[RANK_sfc_pres];
218        float sfc_pres[LON_LEN*LAT_LEN];
219 
220        for(ii = 0; ii < LAT_LEN * LON_LEN; ii++) {
221 	   sfc_pres[ii] = 6;
222        }
223        sfc_pres_start[0] = 0;
224        sfc_pres_start[1] = 0;
225        sfc_pres_start[2] = 0;
226        sfc_pres_count[0] = time_len;
227        sfc_pres_count[1] = lat_len;
228        sfc_pres_count[2] = lon_len;
229        stat = ncmpi_put_vara_float_all(ncid, sfc_pres_id, sfc_pres_start, sfc_pres_count, sfc_pres);
230        check_err(stat,__LINE__,__FILE__);
231    }
232 
233    {				/* store temp_scrn */
234        MPI_Offset temp_scrn_start[RANK_temp_scrn];
235        MPI_Offset temp_scrn_count[RANK_temp_scrn];
236        float temp_scrn[LON_LEN*LAT_LEN];
237 
238        for(ii = 0; ii < LAT_LEN * LON_LEN; ii++) {
239 	   temp_scrn[ii] = 11;
240        }
241        temp_scrn_start[0] = 0;
242        temp_scrn_start[1] = 0;
243        temp_scrn_start[2] = 0;
244        temp_scrn_count[0] = time_len;
245        temp_scrn_count[1] = lat_len;
246        temp_scrn_count[2] = lon_len;
247        stat = ncmpi_put_vara_float_all(ncid, temp_scrn_id, temp_scrn_start, temp_scrn_count, temp_scrn);
248        check_err(stat,__LINE__,__FILE__);
249    }
250 
251    {				/* store qsair_scrn */
252        MPI_Offset qsair_scrn_start[RANK_qsair_scrn];
253        MPI_Offset qsair_scrn_count[RANK_qsair_scrn];
254        float qsair_scrn[LON_LEN*LAT_LEN];
255 
256        for(ii = 0; ii < LAT_LEN * LON_LEN; ii++) {
257 	   qsair_scrn[ii] = 22;
258        }
259        qsair_scrn_start[0] = 0;
260        qsair_scrn_start[1] = 0;
261        qsair_scrn_start[2] = 0;
262        qsair_scrn_count[0] = time_len;
263        qsair_scrn_count[1] = lat_len;
264        qsair_scrn_count[2] = lon_len;
265        stat = ncmpi_put_vara_float_all(ncid, qsair_scrn_id, qsair_scrn_start, qsair_scrn_count, qsair_scrn);
266        check_err(stat,__LINE__,__FILE__);
267    }
268 
269    {				/* store topog */
270        MPI_Offset topog_start[RANK_topog];
271        MPI_Offset topog_count[RANK_topog];
272        float topog[LON_LEN*LAT_LEN];
273 
274        for(ii = 0; ii < LAT_LEN * LON_LEN; ii++) {
275 	   topog[ii] = 33;
276        }
277        topog_start[0] = 0;
278        topog_start[1] = 0;
279        topog_start[2] = 0;
280        topog_count[0] = time_len;
281        topog_count[1] = lat_len;
282        topog_count[2] = lon_len;
283        stat = ncmpi_put_vara_float_all(ncid, topog_id, topog_start, topog_count, topog);
284        check_err(stat,__LINE__,__FILE__);
285    }
286 
287    {				/* store mslp */
288        MPI_Offset mslp_start[RANK_mslp];
289        MPI_Offset mslp_count[RANK_mslp];
290        float mslp[LON_LEN*LAT_LEN];
291 
292        for(ii = 0; ii < LAT_LEN * LON_LEN; ii++) {
293 	   mslp[ii] = 44;
294        }
295        mslp_start[0] = 0;
296        mslp_start[1] = 0;
297        mslp_start[2] = 0;
298        mslp_count[0] = time_len;
299        mslp_count[1] = lat_len;
300        mslp_count[2] = lon_len;
301        stat = ncmpi_put_vara_float_all(ncid, mslp_id, mslp_start, mslp_count, mslp);
302        check_err(stat,__LINE__,__FILE__);
303    }
304 
305    {				/* store sfc_temp */
306        MPI_Offset sfc_temp_start[RANK_sfc_temp];
307        MPI_Offset sfc_temp_count[RANK_sfc_temp];
308        float sfc_temp[LON_LEN*LAT_LEN];
309 
310        for(ii = 0; ii < LAT_LEN * LON_LEN; ii++) {
311 	   sfc_temp[ii] = 55;
312        }
313        sfc_temp_start[0] = 0;
314        sfc_temp_start[1] = 0;
315        sfc_temp_start[2] = 0;
316        sfc_temp_count[0] = time_len;
317        sfc_temp_count[1] = lat_len;
318        sfc_temp_count[2] = lon_len;
319        stat = ncmpi_put_vara_float_all(ncid, sfc_temp_id, sfc_temp_start, sfc_temp_count, sfc_temp);
320        check_err(stat,__LINE__,__FILE__);
321    }
322 
323    {		      /* store zonal_wnd */
324        /* Bug exposed when written in reverse order. */
325        for(i = LVL_LEN - 1; i>=0; i--)
326        /* for(i = 0; i < LVL_LEN; i++) */
327        {
328 	   int izw;
329 	   for(izw = 0; izw < TIME_LEN * LAT_LEN * LON_LEN; izw++) {
330 	       zonal_wnd[izw] = 100 + i;
331 	   }
332 	   zonal_wnd_start[0] = 0;
333 	   zonal_wnd_start[1] = i;
334 	   zonal_wnd_start[2] = 0;
335 	   zonal_wnd_start[3] = 0;
336 	   zonal_wnd_count[0] = time_len;
337 	   zonal_wnd_count[1] = 1;
338 	   zonal_wnd_count[2] = lat_len;
339 	   zonal_wnd_count[3] = lon_len;
340 	   stat = ncmpi_put_vara_float_all(ncid, zonal_wnd_id, zonal_wnd_start, zonal_wnd_count, zonal_wnd);
341 	   check_err(stat,__LINE__,__FILE__);
342        }
343    }
344    stat = ncmpi_close(ncid);
345    check_err(stat,__LINE__,__FILE__);
346    return 0;
347 }
348 
349 int
main(int argc,char ** argv)350 main(int argc, char **argv)
351 {
352     char filename[256], fill_filename[256], nofill_filename[256];
353     int rank, nprocs, err, nerrs=0;
354 
355     MPI_Init(&argc, &argv);
356     MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
357     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
358 
359     if (argc > 2) {
360         if (!rank) printf("Usage: %s [filename]\n",argv[0]);
361         MPI_Finalize();
362         return 0;
363     }
364     if (argc == 2) snprintf(filename, 256, "%s", argv[1]);
365     else           strcpy(filename, "testfile.nc");
366     MPI_Bcast(filename, 256, MPI_CHAR, 0, MPI_COMM_WORLD);
367 
368     char *cmd_str = (char*)malloc(strlen(argv[0]) + 256);
369     sprintf(cmd_str, "*** TESTING C   %s for fill/nofill modes ", basename(argv[0]));
370     if (rank == 0) printf("%-66s ------ ", cmd_str);
371     free(cmd_str);
372 
373     sprintf(fill_filename, "%s.fill", filename);
374     sprintf(nofill_filename, "%s.nofill", filename);
375     nerrs += create_file(nofill_filename, NC_NOFILL);
376     nerrs += create_file(fill_filename, NC_FILL);
377     {
378        int ncid1, ncid2;
379        int nvars1, nvars2;
380        int varid;
381        int badvars;
382 
383        /* compare data in two files created with nofill mode and fill
384 	* mode, which should be identical if all the data were written */
385        err = ncmpi_open(MPI_COMM_WORLD, nofill_filename, NC_NOWRITE, MPI_INFO_NULL, &ncid1); ERR
386        err = ncmpi_open(MPI_COMM_WORLD, fill_filename, NC_NOWRITE, MPI_INFO_NULL, &ncid2); ERR
387 
388        err = ncmpi_inq_nvars(ncid1, &nvars1); ERR
389        err = ncmpi_inq_nvars(ncid2, &nvars2); ERR
390        if (nvars1 != nvars2) {
391            printf("Number of variables disagree %d != %d\n",nvars1,nvars2);
392            nerrs++;
393        }
394        badvars = 0;
395        for(varid = 0; varid < nvars1; varid++) {
396 	   MPI_Offset nvals, nn;
397 	   int ndims, *dimids, dim;
398 	   nc_type vtype;
399 	   char varname1[NC_MAX_NAME];
400 	   char varname2[NC_MAX_NAME];
401 	   /* How many values in this variable to compare? */
402 	   err = ncmpi_inq_varndims(ncid1, varid, &ndims); ERR
403 	   dimids = malloc((ndims + 1) * sizeof(int));
404 	   if (!dimids) printf("Error in file %s line %d\n",__FILE__,__LINE__);
405 	   err = ncmpi_inq_vardimid (ncid1, varid, dimids); ERR
406 	   nvals = 1;
407 	   for(dim = 0; dim < ndims; dim++) {
408 	       MPI_Offset len;
409 	       err = ncmpi_inq_dimlen(ncid1, dimids[dim], &len); ERR
410 	       nvals *= len;
411 	   }
412 	   err = ncmpi_inq_vartype(ncid1, varid, &vtype); ERR
413 	   err = ncmpi_inq_varname(ncid1, varid, varname1); ERR
414 	   err = ncmpi_inq_varname(ncid1, varid, varname2); ERR
415 
416 	   if (vtype != NC_CHAR) {  /* numeric data, just read in as doubles */
417 	       double *data1, *data2;
418 	       /* Allocate space to hold values in both files */
419 	       data1 = malloc((nvals + 1) * sizeof(double));
420 	       if (!data1) printf("Error in file %s line %d\n",__FILE__,__LINE__);
421 	       data2 = malloc((nvals + 1) * sizeof(double));
422 	       if (!data2) printf("Error in file %s line %d\n",__FILE__,__LINE__);
423 	       /* Read in values */
424 	       err = ncmpi_get_var_double_all(ncid1, varid, data1); ERR
425 	       err = ncmpi_get_var_double_all(ncid2, varid, data2); ERR
426 	       /* Compare values */
427 	       for(nn = 0; nn < nvals; nn++) {
428 		   if (data1[nn] != data2[nn]) {
429 		       badvars++;
430 		       fprintf(stderr,
431 			       "\tFrom nofill file, %s[%lu] = %.15g\tFrom fill file, %s[%lu] = %.15g\n",
432 			       varname1, (unsigned long)nn, data1[nn], varname2, (unsigned long)nn, data2[nn]);
433 		       break;
434 		   };
435 	       }
436 	       free(data1);
437 	       free(data2);
438 	   } else {		/* character data */
439 	       char *data1, *data2;
440 	       /* Allocate space to hold values in both files */
441 	       data1 = malloc((nvals + 1) * sizeof(char));
442 	       if (!data1) printf("Error in file %s line %d\n",__FILE__,__LINE__);
443 	       data2 = malloc((nvals + 1) * sizeof(char));
444 	       if (!data2) printf("Error in file %s line %d\n",__FILE__,__LINE__);
445 	       /* Read in values */
446 	       err = ncmpi_get_var_text_all(ncid1, varid, data1); ERR
447 	       err = ncmpi_get_var_text_all(ncid2, varid, data2); ERR
448 	       /* Compare values */
449 	       for(nn = 0; nn < nvals; nn++) {
450 		   if (data1[nn] != data2[nn]) {
451 		       badvars++;
452 		       fprintf(stderr,
453 			       "\tFrom nofill file, %s[%lu] = %d\tFrom fill file, %s[%lu] = %d\n",
454 			       varname1, (unsigned long)nn, data1[nn], varname2, (unsigned long)nn, data2[nn]);
455 		       break;
456 		   };
457 	       }
458 	       free(data1);
459 	       free(data2);
460 	   }
461 	   free(dimids);
462        }
463        if(badvars > 0) printf("Error in file %s line %d\n",__FILE__,__LINE__);
464        err = ncmpi_close(ncid1); ERR
465        err = ncmpi_close(ncid2); ERR
466    }
467 
468     /* check if PnetCDF freed all internal malloc */
469     MPI_Offset malloc_size, sum_size;
470     err = ncmpi_inq_malloc_size(&malloc_size);
471     if (err == NC_NOERR) {
472         MPI_Reduce(&malloc_size, &sum_size, 1, MPI_OFFSET, MPI_SUM, 0, MPI_COMM_WORLD);
473         if (rank == 0 && sum_size > 0)
474             printf("heap memory allocated by PnetCDF internally has %lld bytes yet to be freed\n",
475                    sum_size);
476     }
477 
478     MPI_Allreduce(MPI_IN_PLACE, &nerrs, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
479     if (rank == 0) {
480         if (nerrs) printf(FAIL_STR,nerrs);
481         else       printf(PASS_STR);
482     }
483 
484     MPI_Finalize();
485     return 0;
486 }
487