1 /*
2  *  Copyright (C) 2015, Northwestern University and Argonne National Laboratory
3  *  See COPYRIGHT notice in top-level directory.
4  *
5  *  $Id: add_var.c 2744 2016-12-28 16:25:22Z wkliao $
6  */
7 
8 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
9  *
10  * This program adds two new variables to an existing netCDF file.
11  * It is used to test if PnetCDF can correctly calculate the file offsets
12  * for the two new variables, in particular for files that align the
13  * fix-size variables to a boundary larger than 4 bytes, for instance
14  * a file created by PnetCDF with defaut alignment of 512 bytes.
15  *
16  * The compile and run commands are given below.
17  *
18  *    % mpicc -g -o add_var add_var.c -lpnetcdf
19  *
20  *    % mpiexec -l -n 1 add_var testfile.nc
21  *
22  * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
23 #include <stdio.h>
24 #include <stdlib.h>
25 #include <string.h>
26 #include <libgen.h> /* basename() */
27 #include <pnetcdf.h>
28 
29 #include <testutils.h>
30 
31 #define ERR {if(err!=NC_NOERR){printf("Error at line=%d: %s\n", __LINE__, ncmpi_strerror(err));}}
32 
main(int argc,char ** argv)33 int main(int argc, char** argv) {
34     char filename[256], var_name[NC_MAX_NAME];
35     int i, nvars, rank, nprocs, err, nerrs=0;
36     int ncid, varid, dimid[2];
37     MPI_Offset prev_off, off;
38 
39     MPI_Init(&argc, &argv);
40     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
41     MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
42 
43     if (argc > 2) {
44         if (!rank) printf("Usage: %s [filename]\n",argv[0]);
45         MPI_Finalize();
46         return 0;
47     }
48     if (argc == 2) snprintf(filename, 256, "%s", argv[1]);
49     else           strcpy(filename, "testfile.nc");
50     MPI_Bcast(filename, 256, MPI_CHAR, 0, MPI_COMM_WORLD);
51 
52     if (rank == 0) {
53         char *cmd_str = (char*)malloc(strlen(argv[0]) + 256);
54         sprintf(cmd_str, "*** TESTING C   %s for checking offsets of new variables ", basename(argv[0]));
55         printf("%-66s ------ ", cmd_str); fflush(stdout);
56         free(cmd_str);
57     }
58 
59     err = ncmpi_create(MPI_COMM_WORLD, filename, NC_CLOBBER, MPI_INFO_NULL, &ncid); ERR
60 
61     /* define dimensions */
62     err = ncmpi_def_dim(ncid, "dim_1", 5, &dimid[0]); ERR
63     err = ncmpi_def_dim(ncid, "dim_2", 4, &dimid[1]); ERR
64 
65     /* define a bunch of variables */
66     for (i=0; i<10; i++) {
67         sprintf(var_name, "var_%d", i);
68         err = ncmpi_def_var(ncid, var_name, NC_INT, 2, dimid, &varid); ERR
69     }
70     err = ncmpi_enddef(ncid); ERR
71 
72     /* re-enter define mode */
73     err = ncmpi_redef(ncid); ERR
74 
75     /* add 2 new dimensions */
76     err = ncmpi_def_dim(ncid, "new_dim_1", 5, &dimid[0]); ERR
77     err = ncmpi_def_dim(ncid, "new_dim_2", 4, &dimid[1]); ERR
78 
79     /* add 2 new variables */
80     err = ncmpi_def_var(ncid, "new_var1", NC_INT,   2, dimid, &varid); ERR
81     err = ncmpi_def_var(ncid, "new_var2", NC_FLOAT, 2, dimid, &varid); ERR
82     err = ncmpi_enddef(ncid); ERR
83 
84     err = ncmpi_inq_nvars(ncid, &nvars); ERR
85     err = ncmpi_inq_varoffset(ncid, 0, &prev_off); ERR
86     for (i=1; i<nvars; i++) {
87         err = ncmpi_inq_varoffset(ncid, i, &off); ERR
88         if (off < prev_off + 5*4*4) { /* each variable is of size 5*4*4 bytes */
89             err = ncmpi_inq_varname(ncid, i, var_name); ERR
90             printf("Error in %s line %d: variable %s offset is set incorrectly\n",
91                    __FILE__,__LINE__,var_name);
92             nerrs++;
93         }
94         prev_off = off;
95     }
96 
97     err = ncmpi_close(ncid); ERR
98 
99     /* check if PnetCDF freed all internal malloc */
100     MPI_Offset malloc_size, sum_size;
101     err = ncmpi_inq_malloc_size(&malloc_size);
102     if (err == NC_NOERR) {
103         MPI_Reduce(&malloc_size, &sum_size, 1, MPI_OFFSET, MPI_SUM, 0, MPI_COMM_WORLD);
104         if (rank == 0 && sum_size > 0)
105             printf("heap memory allocated by PnetCDF internally has %lld bytes yet to be freed\n",
106                    sum_size);
107     }
108 
109     MPI_Allreduce(MPI_IN_PLACE, &nerrs, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
110     if (rank == 0) {
111         if (nerrs) printf(FAIL_STR,nerrs);
112         else       printf(PASS_STR);
113     }
114 
115     MPI_Finalize();
116     return 0;
117 }
118 
119