File: put_varn_float.c

package info (click to toggle)
pnetcdf 1.14.1-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 13,812 kB
  • sloc: ansic: 85,298; f90: 10,707; fortran: 9,283; cpp: 8,864; makefile: 3,084; perl: 2,833; sh: 2,538; yacc: 1,227; lex: 216
file content (252 lines) | stat: -rw-r--r-- 8,713 bytes parent folder | download | duplicates (2)
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
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
/*********************************************************************
 *
 *  Copyright (C) 2013, Northwestern University and Argonne National Laboratory
 *  See COPYRIGHT notice in top-level directory.
 *
 *********************************************************************/
/* $Id$ */

/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
 * This example shows how to use a single call of ncmpi_put_varn_float_all()
 * to write a sequence of one-element requests with arbitrary array indices.
 *
 * The compile and run commands are given below, together with an ncmpidump of
 * the output file.
 *
 *    % mpicc -O2 -o put_varn_float put_varn_float.c -lpnetcdf
 *    % mpiexec -n 4 ./put_varn_float /pvfs2/wkliao/testfile.nc
 *    % ncmpidump /pvfs2/wkliao/testfile.nc
 *    netcdf testfile {
 *    // file format: CDF-5 (big variables)
 *    dimensions:
 *             Y = 4 ;
 *             X = 10 ;
 *    variables:
 *             int var(Y, X) ;
 *    data:
 *
 *     var =
 *       3, 3, 3, 1, 1, 0, 0, 2, 1, 1,
 *       0, 2, 2, 2, 3, 1, 1, 2, 2, 2,
 *       1, 1, 2, 3, 3, 3, 0, 0, 1, 1,
 *       0, 0, 0, 2, 1, 1, 1, 3, 3, 3 ;
 *    }
 * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */

#include <stdio.h>
#include <stdlib.h>
#include <string.h> /* strcpy(), strncpy() */
#include <unistd.h> /* getopt() */
#include <mpi.h>
#include <pnetcdf.h>

#define NY 4
#define NX 10
#define NDIMS 2

static int verbose;

#define ERR {if(err!=NC_NOERR){printf("Error at %s:%d : %s\n", __FILE__,__LINE__, ncmpi_strerror(err));nerrs++;}}

static void
usage(char *argv0)
{
    char *help =
    "Usage: %s [-h] | [-q] [file_name]\n"
    "       [-h] Print help\n"
    "       [-q] Quiet mode (reports when fail)\n"
    "       [filename] output netCDF file name\n";
    fprintf(stderr, help, argv0);
}

/*----< pnetcdf_check_mem_usage() >------------------------------------------*/
/* check PnetCDF library internal memory usage */
static int
pnetcdf_check_mem_usage(MPI_Comm comm)
{
    int err, nerrs=0, rank;
    MPI_Offset malloc_size, sum_size;

    MPI_Comm_rank(comm, &rank);

    /* print info about PnetCDF internal malloc usage */
    err = ncmpi_inq_malloc_max_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 && verbose)
            printf("maximum heap memory allocated by PnetCDF internally is %lld bytes\n",
                   sum_size);

        /* check if there is any PnetCDF internal malloc residue */
        err = ncmpi_inq_malloc_size(&malloc_size);
        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);
    }
    else if (err != NC_ENOTENABLED) {
        printf("Error at %s:%d: %s\n", __FILE__,__LINE__,ncmpi_strerror(err));
        nerrs++;
    }
    return nerrs;
}

int main(int argc, char** argv)
{
    extern int optind;
    char filename[256];
    int i, rank, nprocs, err, nerrs=0;
    int ncid, cmode, varid, dimid[2], num_reqs;
    float *buffer;
    MPI_Offset **starts=NULL;

    MPI_Init(&argc, &argv);
    MPI_Comm_rank(MPI_COMM_WORLD, &rank);
    MPI_Comm_size(MPI_COMM_WORLD, &nprocs);

    verbose = 1;

    /* 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]);

    if (verbose && nprocs != 4 && rank == 0)
        printf("Warning: this program is intended to run on 4 processes\n");

    /* create a new file for writing ----------------------------------------*/
    cmode = NC_CLOBBER;
    err = ncmpi_create(MPI_COMM_WORLD, filename, cmode, MPI_INFO_NULL, &ncid);
    ERR

    /* create a global array of size NY * NX */
    err = ncmpi_def_dim(ncid, "Y", NY, &dimid[0]);
    ERR
    err = ncmpi_def_dim(ncid, "X", NX, &dimid[1]);
    ERR
    err = ncmpi_def_var(ncid, "var", NC_FLOAT, NDIMS, dimid, &varid);
    ERR
    if (nprocs < 4) { /* need 4 processes to fill the variables */
        err = ncmpi_set_fill(ncid, NC_FILL, NULL);
        ERR
    }
    err = ncmpi_enddef(ncid);
    ERR

    /* pick arbitrary numbers of requests for 4 processes */
    num_reqs = 0;
    if (rank == 0)      num_reqs = 8;
    else if (rank == 1) num_reqs = 13;
    else if (rank == 2) num_reqs = 9;
    else if (rank == 3) num_reqs = 10;

    if (num_reqs > 0) {
        starts    = (MPI_Offset**) malloc(sizeof(MPI_Offset*) * num_reqs);
        starts[0] = (MPI_Offset*)  calloc(num_reqs*NDIMS, sizeof(MPI_Offset));
        for (i=1; i<num_reqs; i++)
            starts[i] = starts[i-1] + NDIMS;
    }

    /* assign arbitrary starts */
    const int y=0, x=1;
    if (rank == 0) {
        starts[0][y] = 0; starts[0][x] = 5;
        starts[1][y] = 1; starts[1][x] = 0;
        starts[2][y] = 2; starts[2][x] = 6;
        starts[3][y] = 3; starts[3][x] = 0;
        starts[4][y] = 0; starts[4][x] = 6;
        starts[5][y] = 2; starts[5][x] = 7;
        starts[6][y] = 3; starts[6][x] = 1;
        starts[7][y] = 3; starts[7][x] = 2;
        /* rank 0 is writing the following locations: ("-" means skip)
                  -  -  -  -  -  0  0  -  -  -
                  0  -  -  -  -  -  -  -  -  -
                  -  -  -  -  -  -  0  0  -  -
                  0  0  0  -  -  -  -  -  -  -
         */
    } else if (rank ==1) {
        starts[ 0][y] = 0; starts[ 0][x] = 3;
        starts[ 1][y] = 0; starts[ 1][x] = 8;
        starts[ 2][y] = 1; starts[ 2][x] = 5;
        starts[ 3][y] = 2; starts[ 3][x] = 0;
        starts[ 4][y] = 2; starts[ 4][x] = 8;
        starts[ 5][y] = 3; starts[ 5][x] = 4;
        starts[ 6][y] = 0; starts[ 6][x] = 4;
        starts[ 7][y] = 0; starts[ 7][x] = 9;
        starts[ 8][y] = 1; starts[ 8][x] = 6;
        starts[ 9][y] = 2; starts[ 9][x] = 1;
        starts[10][y] = 2; starts[10][x] = 9;
        starts[11][y] = 3; starts[11][x] = 5;
        starts[12][y] = 3; starts[12][x] = 6;
        /* rank 1 is writing the following locations: ("-" means skip)
                  -  -  -  1  1  -  -  -  1  1
                  -  -  -  -  -  1  1  -  -  -
                  1  1  -  -  -  -  -  -  1  1
                  -  -  -  -  1  1  1  -  -  -
         */
    } else if (rank ==2) {
        starts[0][y] = 0; starts[0][x] = 7;
        starts[1][y] = 1; starts[1][x] = 1;
        starts[2][y] = 1; starts[2][x] = 7;
        starts[3][y] = 2; starts[3][x] = 2;
        starts[4][y] = 3; starts[4][x] = 3;
        starts[5][y] = 1; starts[5][x] = 2;
        starts[6][y] = 1; starts[6][x] = 8;
        starts[7][y] = 1; starts[7][x] = 3;
        starts[8][y] = 1; starts[8][x] = 9;
        /* rank 2 is writing the following locations: ("-" means skip)
                  -  -  -  -  -  -  -  2  -  -
                  -  2  2  2  -  -  -  2  2  2
                  -  -  2  -  -  -  -  -  -  -
                  -  -  -  2  -  -  -  -  -  -
         */
    } else if (rank ==3) {
        starts[0][y] = 0; starts[0][x] = 0;
        starts[1][y] = 1; starts[1][x] = 4;
        starts[2][y] = 2; starts[2][x] = 3;
        starts[3][y] = 3; starts[3][x] = 7;
        starts[4][y] = 0; starts[4][x] = 1;
        starts[5][y] = 2; starts[5][x] = 4;
        starts[6][y] = 3; starts[6][x] = 8;
        starts[7][y] = 0; starts[7][x] = 2;
        starts[8][y] = 2; starts[8][x] = 5;
        starts[9][y] = 3; starts[9][x] = 9;
        /* rank 3 is writing the following locations: ("-" means skip)
                  3  3  3  -  -  -  -  -  -  -
                  -  -  -  -  3  -  -  -  -  -
                  -  -  -  3  3  3  -  -  -  -
                  -  -  -  -  -  -  -  3  3  3
         */
    }

    /* allocate I/O buffer and initialize its contents */
    buffer = (float*) malloc(sizeof(float) * num_reqs);
    for (i=0; i<num_reqs; i++) buffer[i] = (float)rank;

    /* set the buffer pointers to different offsets to the I/O buffer */
    err = ncmpi_put_varn_float_all(ncid, varid, num_reqs, starts, NULL, buffer);
    ERR

    err = ncmpi_close(ncid);
    ERR

    free(buffer);
    if (num_reqs > 0) {
        free(starts[0]);
        free(starts);
    }

    nerrs += pnetcdf_check_mem_usage(MPI_COMM_WORLD);

    MPI_Finalize();
    return (nerrs > 0);
}