netcdf-c/nc_test/tst_cdf5format.c

188 lines
4.2 KiB
C
Raw Normal View History

2015-08-16 06:26:35 +08:00
/* This is part of the netCDF package.
2018-12-07 05:56:42 +08:00
Copyright 2018 University Corporation for Atmospheric Research/Unidata
2015-08-16 06:26:35 +08:00
See COPYRIGHT file for conditions of use.
Test fix of bug involving creation of a file with PnetCDF APIs,
2015-08-16 06:26:35 +08:00
then opening and modifying the file with netcdf.
2017-11-06 20:35:34 +08:00
Author: Wei-keng Liao, Ed Hartnett
2015-08-16 06:26:35 +08:00
*/
2017-11-06 20:35:34 +08:00
#include <config.h>
2015-08-16 06:26:35 +08:00
#include <nc_tests.h>
2017-11-06 20:35:34 +08:00
#include <err_macros.h>
2015-08-16 06:26:35 +08:00
#define NVARS 6
2017-11-06 20:35:34 +08:00
#define NX 5
#define NDIM2 2
#define FILENAME "tst_cdf5format.nc"
2015-08-16 06:26:35 +08:00
2017-11-06 21:19:34 +08:00
/* Write a file with 2 dims and 6 vars, including some sample data. */
2015-08-16 06:26:35 +08:00
int
2017-11-06 20:35:34 +08:00
write2(int ncid, int parallel)
2015-08-16 06:26:35 +08:00
{
2017-11-06 20:35:34 +08:00
int dimid[NDIM2];
char str[NC_MAX_NAME + 1];
int varid[NVARS];
2019-02-12 23:10:39 +08:00
int i;
int j;
2018-12-07 05:56:42 +08:00
2017-11-06 20:35:34 +08:00
/* define dimension */
if (nc_def_dim(ncid, "Y", NC_UNLIMITED, &dimid[0])) ERR;
if (nc_def_dim(ncid, "X", NX, &dimid[1])) ERR;
/* Define vars. */
2019-02-12 23:10:39 +08:00
for (i = 0; i < NVARS; i++)
2017-11-06 20:35:34 +08:00
{
if (i % 2)
{
snprintf(str, sizeof(str), "fixed_var_%d",i);
2017-11-06 20:35:34 +08:00
if (nc_def_var(ncid, str, NC_INT, 1, &dimid[1], &varid[i])) ERR;
}
else
{
snprintf(str, sizeof(str), "record_var_%d",i);
2017-11-06 20:35:34 +08:00
if (nc_def_var(ncid, str, NC_INT, 2, dimid, &varid[i])) ERR;
}
}
2018-12-07 05:56:42 +08:00
2017-11-06 20:35:34 +08:00
if (nc_enddef(ncid)) ERR;
2018-12-07 05:56:42 +08:00
2017-11-06 20:35:34 +08:00
/* write all variables */
2019-02-12 23:10:39 +08:00
for (i = 0; i < NVARS; i++)
2017-11-06 20:35:34 +08:00
{
size_t start[NDIM2] = {0, 0};
size_t count[NDIM2];
int buf[NX];
2018-12-07 05:56:42 +08:00
2017-11-06 20:35:34 +08:00
/* Initialize some data. */
2019-02-12 23:10:39 +08:00
for (j = 0; j < NX; j++)
2017-11-06 20:35:34 +08:00
buf[j] = i * 10 + j;
/* Write the data. */
if (i % 2)
{
count[0] = NX;
if (nc_put_vara_int(ncid, varid[i], start, count, buf)) ERR;
}
else
{
count[0] = 1;
count[1] = NX;
if (nc_put_vara_int(ncid, varid[i], start, count, buf)) ERR;
}
}
return 0;
2015-08-16 06:26:35 +08:00
}
2017-11-06 21:19:34 +08:00
/* Add some attributes to the vars of an open file. */
2015-08-16 06:26:35 +08:00
int
extend(int ncid)
{
2017-11-06 20:35:34 +08:00
int i;
char str[32];
if (nc_redef(ncid)) ERR;
2018-12-07 05:56:42 +08:00
2017-11-06 20:35:34 +08:00
/* add attributes to make header grow */
for (i = 0; i < NVARS; i++)
{
snprintf(str, sizeof(str), "annotation_for_var_%d", i);
2017-11-06 20:35:34 +08:00
if (nc_put_att_text(ncid, i, "text_attr", strlen(str), str)) ERR;
}
if (nc_enddef(ncid)) ERR;
return NC_NOERR;
2015-08-16 06:26:35 +08:00
}
2017-11-06 21:19:34 +08:00
/* Read the file and check the data. */
2015-08-16 06:26:35 +08:00
int
2017-11-06 20:35:34 +08:00
read2(int ncid)
2015-08-16 06:26:35 +08:00
{
2019-02-12 23:10:39 +08:00
int i;
int j;
for (i = 0; i < NVARS; i++)
2017-11-06 20:35:34 +08:00
{
int buf[NX];
2018-12-07 05:56:42 +08:00
size_t start[2] = {0, 0}, count[2];
2017-11-06 20:35:34 +08:00
if (i % 2)
{
count[0] = NX;
}
else
{
count[0] = 1;
count[1] = NX;
}
2018-12-07 05:56:42 +08:00
if (nc_get_vara_int(ncid, i, start, count, buf)) ERR;
2019-02-12 23:10:39 +08:00
for (j = 0; j < NX; j++)
2017-11-06 20:35:34 +08:00
{
if (buf[j] != i * 10 + j)
{
printf("unexpected read value var i=%d buf[j=%d]=%d should be %d\n",
i, j, buf[j], i * 10 + j);
ERR;
}
}
}
return 0;
2015-08-16 06:26:35 +08:00
}
int main(int argc, char* argv[])
{
2017-11-06 20:35:34 +08:00
int rank, nprocs, ncid, cmode;
MPI_Comm comm = MPI_COMM_SELF;
MPI_Info info = MPI_INFO_NULL;
MPI_Init(&argc, &argv);
MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
if (rank > 0)
return 2;
printf("\nWrite using PNETCDF; Read using classic netCDF...");
2017-11-06 21:19:34 +08:00
{
/* Create a netCDF classic file with PnetCDF. */
2018-09-23 09:22:34 +08:00
cmode = NC_CLOBBER;
2017-11-06 21:19:34 +08:00
if (nc_create_par(FILENAME, cmode, comm, info, &ncid)) ERR;
if (write2(ncid, 1)) ERR;
if (nc_close(ncid)) ERR;
2018-12-07 05:56:42 +08:00
2017-11-06 21:19:34 +08:00
/* Re-open the file with pnetCDF (parallel) and add var attributes. */
2018-09-23 09:22:34 +08:00
if (nc_open_par(FILENAME, NC_WRITE, comm, info, &ncid)) ERR;
2017-11-06 21:19:34 +08:00
if (extend(ncid)) ERR;
if (nc_close(ncid)) ERR;
2018-12-07 05:56:42 +08:00
2017-11-06 21:19:34 +08:00
/* Open with classic and check. */
if (nc_open(FILENAME, 0, &ncid)) ERR;
if (read2(ncid)) ERR;
if (nc_close(ncid)) ERR;
}
SUMMARIZE_ERR;
2017-11-06 20:35:34 +08:00
2017-11-06 21:19:34 +08:00
printf("\nWrite using CDF-5; Read using PNETCDF...");
{
/* Create a file with CDF5. */
cmode = NC_CDF5 | NC_CLOBBER;
if (nc_create(FILENAME, cmode, &ncid)) ERR;
if (write2(ncid, 0)) ERR;
if (nc_close(ncid)) ERR;
2018-12-07 05:56:42 +08:00
2017-11-06 21:19:34 +08:00
/* Re-open the file with CDF5 and add some atts. */
if (nc_open(FILENAME, NC_WRITE, &ncid)) ERR;
if (extend(ncid)) ERR;
if (nc_close(ncid)) ERR;
/* Re-open with PnetCDF and check. */
2018-09-23 09:22:34 +08:00
cmode = NC_NOCLOBBER;
2017-11-06 21:19:34 +08:00
if (nc_open_par(FILENAME, cmode, comm, info, &ncid)) ERR;
if (read2(ncid)) ERR;
if (nc_close(ncid)) ERR;
}
SUMMARIZE_ERR;
2017-11-06 20:35:34 +08:00
MPI_Finalize();
FINAL_RESULTS;
return 0;
2015-08-16 06:26:35 +08:00
}