SCALE-RM
scale_file_netcdf.c
Go to the documentation of this file.
1 #include "scale_file.h"
2 #ifndef MPI_INCLUDED
3 #define MPI_INCLUDED
4 #endif
5 #include "netcdf.h"
6 
7 #define TEPS 1e-6
8 #define NTMAX 102400
9 
10 #define MIN(a,b) ((a)<(b) ? (a) : (b))
11 
12 static int32_t ERROR_SUPPRESS = 0;
13 
14 #define CHECK_ERROR(func) \
15  { \
16  int status_ = (func); \
17  if (status_ != NC_NOERR) { \
18  if ( ! ERROR_SUPPRESS ) { \
19  fprintf(stderr, "Error: at line %d in %s\n", __LINE__, __FILE__); \
20  fprintf(stderr, " %s\n", nc_strerror(status_)); \
21  } \
22  return ERROR_CODE; \
23  } \
24  }
25 
26 #ifdef PNETCDF
27 #include "pnetcdf.h"
28 #define CHECK_PNC_ERROR(func) \
29  { \
30  int status_ = (func); \
31  if (status_ != NC_NOERR) { \
32  if ( ! ERROR_SUPPRESS ) { \
33  fprintf(stderr, "Error: at line %d in %s\n", __LINE__, __FILE__); \
34  fprintf(stderr, " %s\n", ncmpi_strerror(status_)); \
35  } \
36  return ERROR_CODE; \
37  } \
38  }
39 #else
40 #define CHECK_PNC_ERROR(func) \
41  { \
42  fprintf(stderr, "pnetCDF is necessary for shared_mode\n"); \
43  fprintf(stderr, "Please re-compile with pnetCDF\n"); \
44  return ERROR_CODE; \
45  }
46 #define ncmpi_inq_attid(a,b,c,d) NC2_ERR
47 #define ncmpi_inq_varid(a,b,c) NC2_ERR
48 #define ncmpi_inq_dimid(a,b,c) NC2_ERR
49 #endif
50 
51 #define NCTYPE2TYPE(nctype, type) \
52  { \
53  switch ( nctype ) { \
54  case NC_FLOAT: \
55  type = File_REAL4; \
56  break; \
57  case NC_DOUBLE: \
58  type = File_REAL8; \
59  break; \
60  case NC_SHORT: \
61  type = File_INTEGER2; \
62  break; \
63  case NC_INT: \
64  type = File_INTEGER4; \
65  break; \
66  case NC_CHAR: \
67  type = File_TEXT; \
68  break; \
69  default: \
70  fprintf(stderr, "unsupported data type: %d\n", xtype); \
71  return ERROR_CODE; \
72  } \
73  }
74 
75 #define TYPE2NCTYPE(type, nctype) \
76  { \
77  switch ( type ) { \
78  case File_REAL4: \
79  nctype = NC_FLOAT; \
80  break; \
81  case File_REAL8: \
82  nctype = NC_DOUBLE; \
83  break; \
84  default: \
85  fprintf(stderr, "unsupported data type: %d\n", xtype); \
86  return ERROR_CODE; \
87  } \
88  }
89 
90 
91 #define DEFAULT_DEFLATE_LEVEL 2
92 
93 typedef struct {
94  int ncid;
95  char time_units[File_HMID+1];
96  char calendar[File_HSHORT+1];
98 #if defined(NETCDF3) || defined(PNETCDF)
99  int defmode;
100 #endif
102  char fname[256]; // used for debugging, to be deleted
103 } fileinfo_t;
104 
105 typedef struct {
106  int ncid;
107  int dimid;
108  int varid;
109  int bndsid;
110  int count;
111  real64_t t;
112  real64_t tint;
113  real64_t *tval;
114  char name[File_HSHORT+1];
115 } tdim_t;
116 
117 typedef struct {
118  int ncid;
119  int varid;
121  size_t *start;
122  size_t *count;
123  size_t ndims;
124  size_t ndims_t;
125 } varinfo_t;
126 
127 static fileinfo_t *files[FILE_MAX];
128 static int nfile = 0;
129 static varinfo_t *vars[VAR_MAX];
130 static int nvar = 0;
131 static tdim_t *tdims[VAR_MAX];
132 static int nt = 0;
133 
134 
135 static inline int32_t file_enddef( const int32_t fid, const int ncid ) // (in)
136 {
137 #if defined(NETCDF3) || defined(PNETCDF)
138  if (files[fid]->defmode == 1) {
139  if ( files[fid]->shared_mode )
140  CHECK_PNC_ERROR( ncmpi_enddef(ncid) )
141 #ifdef NETCDF3
142  else
143  CHECK_ERROR( nc_enddef(ncid) )
144 #endif
145  files[fid]->defmode = 0;
146  }
147 #endif
148 
149  return SUCCESS_CODE;
150 }
151 
152 static inline int32_t file_redef( const int32_t fid, const int ncid ) // (in)
153 {
154 #if defined(NETCDF3) || defined(PNETCDF)
155  if (files[fid]->defmode == 0) {
156  if (files[fid]->shared_mode)
157  CHECK_PNC_ERROR( ncmpi_redef(ncid) )
158 #ifdef NETCDF3
159  else
160  CHECK_ERROR( nc_redef(ncid) )
161 #endif
162  files[fid]->defmode = 1;
163  }
164 #endif
165 
166  return SUCCESS_CODE;
167 }
168 
169 
170 int32_t file_open_c( int32_t *fid, // (out)
171  const char *fname, // (in)
172  const int32_t mode, // (in)
173  const MPI_Comm comm ) // (in)
174 {
175  int ncid;
176  int len;
177  int shared_mode;
178  char _fname[File_HLONG+4];
179  int add_suffix;
180 
181  if ( nfile >= FILE_MAX ) {
182  fprintf(stderr, "exceed max number of file limit\n");
183  return ERROR_CODE;
184  }
185 
186  len = strlen(fname);
187  strcpy(_fname, fname);
188 
189  if ( mode==File_FREAD || mode==File_FAPPEND ) {
190  FILE *fp = fopen(_fname, "r");
191  if ( fp==NULL ) {
192  add_suffix = 1;
193  } else {
194  fclose(fp);
195  add_suffix = 0;
196  }
197  } else
198  add_suffix = 1;
199 
200  if ( add_suffix )
201  if (fname[len-3] != '.' || fname[len-2] != 'n' || fname[len-1] != 'c' )
202  strcat(_fname, ".nc");
203 
204  if ( comm == MPI_COMM_NULL || comm == MPI_COMM_SELF )
205  shared_mode = 0;
206  else
207  shared_mode = 1;
208 
209  switch ( mode ) {
210  case File_FREAD:
211  if ( shared_mode )
212  CHECK_PNC_ERROR( ncmpi_open(comm, _fname, NC_NOWRITE, MPI_INFO_NULL, &ncid) )
213  else
214  CHECK_ERROR( nc_open(_fname, NC_NOWRITE, &ncid) )
215  break;
216  case File_FWRITE:
217  if ( shared_mode )
218  CHECK_PNC_ERROR( ncmpi_create(comm, _fname, NC_CLOBBER|NC_64BIT_OFFSET, MPI_INFO_NULL, &ncid) )
219  else
220 #ifdef NETCDF3
221  CHECK_ERROR( nc_create(_fname, NC_CLOBBER|NC_64BIT_OFFSET, &ncid) )
222 #else
223  CHECK_ERROR( nc_create(_fname, NC_CLOBBER|NC_NETCDF4, &ncid) )
224 #endif
225  break;
226  case File_FAPPEND:
227  if ( shared_mode )
228  CHECK_PNC_ERROR( ncmpi_open(comm, _fname, NC_WRITE, MPI_INFO_NULL, &ncid) )
229  else
230  CHECK_ERROR( nc_open(_fname, NC_WRITE, &ncid) )
231  break;
232  default:
233  fprintf(stderr, "invalid mode type\n");
234  return ERROR_CODE;
235  }
236 
237  files[nfile] = (fileinfo_t*) malloc(sizeof(fileinfo_t));
238  files[nfile]->ncid = ncid;
239  files[nfile]->deflate_level = DEFAULT_DEFLATE_LEVEL;
240 #if defined(NETCDF3) || defined(PNETCDF)
241  if ( mode == File_FWRITE )
242  files[nfile]->defmode = 1;
243  else
244  files[nfile]->defmode = 0;
245 #endif
246 
247  files[nfile]->shared_mode = shared_mode; /* shared-file I/O mode */
248  strcpy(files[nfile]->fname, fname);
249  *fid = nfile;
250 
251  nfile++;
252 
253  return SUCCESS_CODE;
254 }
255 
256 int32_t file_get_dim_length_c( const int32_t fid, // (in)
257  const char* dimname, // (in)
258  int32_t *len ) // (out)
259 {
260  int ncid, dimid;
261 
262  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
263  ncid = files[fid]->ncid;
264 
265  if ( files[fid]->shared_mode ) {
266  MPI_Offset l;
267  CHECK_PNC_ERROR( ncmpi_inq_dimid(ncid, dimname, &dimid) )
268  CHECK_PNC_ERROR( ncmpi_inq_dimlen(ncid, dimid, &l) )
269  *len = l;
270  } else {
271  size_t l;
272  CHECK_ERROR( nc_inq_dimid(ncid, dimname, &dimid) )
273  CHECK_ERROR( nc_inq_dimlen(ncid, dimid, &l) )
274  *len = l;
275  }
276 
277  return SUCCESS_CODE;
278 }
279 
280 int32_t file_set_option_c( const int32_t fid, // (in)
281  const char* filetype, // (in)
282  const char* key, // (in)
283  const char* val) // (in)
284 {
285  if ( strcmp(filetype, "netcdf") != 0 ) return SUCCESS_CODE;
286 
287  if ( strcmp(key, "deflate_level") == 0 ) {
288  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
289  files[fid]->deflate_level = atoi(val);
290  return SUCCESS_CODE;
291  } else {
292  return ERROR_CODE;
293  }
294 }
295 
296 int32_t file_get_nvars_c( const int32_t fid, // (in)
297  int32_t *nvars )// (out)
298 {
299  int ncid;
300  int ndims, ngatts, unlimdim;
301 
302  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
303  ncid = files[fid]->ncid;
304 
305  if ( files[fid]->shared_mode )
306  CHECK_PNC_ERROR( ncmpi_inq(ncid, &ndims, nvars, &ngatts, &unlimdim) )
307  else
308  CHECK_ERROR( nc_inq(ncid, &ndims, nvars, &ngatts, &unlimdim) )
309 
310  return SUCCESS_CODE;
311 }
312 
313 int32_t file_get_varname_c( const int32_t fid, // (in)
314  const int32_t vid, // (in)
315  char *name, // (out)
316  const int32_t len ) // (in)
317 {
318  int ncid, varid;
319  char buf[MAX_NC_NAME+1];
320  int i;
321 
322  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
323  ncid = files[fid]->ncid;
324  varid = vid-1; // index starts from 1 in fortran space
325 
326  if ( files[fid]->shared_mode )
327  CHECK_PNC_ERROR( ncmpi_inq_varname(ncid, varid, buf) )
328  else
329  CHECK_ERROR( nc_inq_varname(ncid, varid, buf) )
330 
331  for (i=0; i<MIN(len-1,strlen(buf)); i++)
332  name[i] = buf[i];
333  name[i] = '\0';
334 
335  return SUCCESS_CODE;
336 }
337 
338 int32_t file_get_datainfo_c( datainfo_t *dinfo, // (out)
339  const int32_t fid, // (in)
340  const char* varname, // (in)
341  const int32_t step, // (in)
342  const int32_t suppress)// (in)
343 {
344  int ncid, varid;
345  nc_type xtype;
346  int rank;
347  int dimids[RANK_MAX], tdim, uldims[NC_MAX_DIMS];
348  char name[NC_MAX_NAME+1];
349  char *buf;
350  size_t size;
351  int natts;
352  int i, n;
353 
354  ERROR_SUPPRESS = suppress;
355 
356  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
357  ncid = files[fid]->ncid;
358  if ( files[fid]->shared_mode )
359  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, varname, &varid) )
360  else
361  CHECK_ERROR( nc_inq_varid(ncid, varname, &varid) )
362 
363  // fid
364  dinfo->fid = fid;
365  // varname
366  strcpy(dinfo->varname, varname);
367 
368  if ( files[fid]->shared_mode ) {
369  MPI_Offset l;
370  // datatype
371  CHECK_PNC_ERROR( ncmpi_inq_vartype(ncid, varid, &xtype) )
372  NCTYPE2TYPE(xtype, dinfo->datatype);
373  // rank
374  CHECK_PNC_ERROR( ncmpi_inq_varndims(ncid, varid, &rank) )
375  if ( rank > 0 ) {
376 #ifdef PNETCDF
377  // description
378  if ( ncmpi_inq_attlen(ncid, varid, "long_name", &l) == NC_NOERR ) {
379  buf = (char*) malloc(l+1);
380  CHECK_PNC_ERROR( ncmpi_get_att_text(ncid, varid, "long_name", buf) )
381  for (i=0; i<MIN(File_HMID-1,l); i++)
382  dinfo->description[i] = buf[i];
383  dinfo->description[i] = '\0';
384  free(buf);
385  } else
386  dinfo->description[0] = '\0';
387  // units
388  if ( ncmpi_inq_attlen(ncid, varid, "units", &l) == NC_NOERR ) {
389  buf = (char*) malloc(l+1);
390  CHECK_PNC_ERROR( ncmpi_get_att_text(ncid, varid, "units", buf) )
391  for (i=0; i<MIN(File_HSHORT-1,l); i++)
392  dinfo->units[i] = buf[i];
393  dinfo->units[i] = '\0';
394  free(buf);
395  } else
396  dinfo->units[0] = '\0';
397  // standard_name
398  if ( ncmpi_inq_attlen(ncid, varid, "standard_name", &l) == NC_NOERR ) {
399  buf = (char*) malloc(l+1);
400  CHECK_PNC_ERROR( ncmpi_get_att_text(ncid, varid, "standard_name", buf) )
401  for (i=0; i<MIN(File_HMID-1,l); i++)
402  dinfo->standard_name[i] = buf[i];
403  dinfo->standard_name[i] = '\0';
404  free(buf);
405  } else
406  dinfo->standard_name[0] = '\0';
407 #endif
408  // dimensions
409  CHECK_PNC_ERROR( ncmpi_inq_vardimid(ncid, varid, dimids) )
410 #if 1
411  CHECK_PNC_ERROR( ncmpi_inq_unlimdim(ncid, uldims) )
412  n = 1;
413 #else
414  CHECK_PNC_ERROR( ncmpi_inq_unlimdims(ncid, &n, uldims) )
415 #endif
416  }
417  // attributes
418  dinfo->natts = 0;
419  CHECK_PNC_ERROR( ncmpi_inq_varnatts(ncid, varid, &natts) )
420  for (i=0; i<natts; i++) {
421  CHECK_PNC_ERROR( ncmpi_inq_attname(ncid, varid, i, name) )
422  if ( strcmp(name, "long_name") &&
423  strcmp(name, "description") &&
424  strcmp(name, "units") &&
425  strcmp(name, "standard_name") ) {
426  strcpy(dinfo->att_name+File_HSHORT*dinfo->natts, name);
427  CHECK_PNC_ERROR( ncmpi_inq_atttype(ncid, varid, name, &xtype) )
428  NCTYPE2TYPE(xtype, dinfo->att_type[dinfo->natts]);
429  CHECK_PNC_ERROR( ncmpi_inq_attlen(ncid, varid, name, &l) )
430  dinfo->att_len[dinfo->natts] = l;
431  dinfo->natts++;
432  }
433  }
434  } else {
435  size_t l;
436  // datatype
437  CHECK_ERROR( nc_inq_vartype(ncid, varid, &xtype) )
438  NCTYPE2TYPE(xtype, dinfo->datatype);
439  // rank
440  CHECK_ERROR( nc_inq_varndims(ncid, varid, &rank) )
441  if ( rank > 0 ) {
442  // description
443  if ( nc_inq_attlen(ncid, varid, "long_name", &l) == NC_NOERR ) {
444  buf = (char*) malloc(l+1);
445  CHECK_ERROR( nc_get_att_text(ncid, varid, "long_name", buf) )
446  } else if ( nc_inq_attlen(ncid, varid, "description", &l) == NC_NOERR ) { // for WRF file
447  buf = (char*) malloc(l+1);
448  CHECK_ERROR( nc_get_att_text(ncid, varid, "description", buf) )
449  } else
450  l = 0;
451  for (i=0; i<MIN(File_HMID-1,l); i++)
452  dinfo->description[i] = buf[i];
453  dinfo->description[i] = '\0';
454  if ( l>0 ) free(buf);
455  // units
456  if ( nc_inq_attlen(ncid, varid, "units", &l) == NC_NOERR ) {
457  buf = (char*) malloc(l+1);
458  CHECK_ERROR( nc_get_att_text(ncid, varid, "units", buf) )
459  for (i=0; i<MIN(File_HSHORT-1,l); i++)
460  dinfo->units[i] = buf[i];
461  dinfo->units[i] = '\0';
462  free(buf);
463  } else
464  dinfo->units[0] = '\0';
465  // standard_name
466  if ( nc_inq_attlen(ncid, varid, "standard_name", &l) == NC_NOERR ) {
467  buf = (char*) malloc(l+1);
468  CHECK_ERROR( nc_get_att_text(ncid, varid, "standard_name", buf) )
469  for (i=0; i<MIN(File_HMID-1,l); i++)
470  dinfo->standard_name[i] = buf[i];
471  dinfo->standard_name[i] = '\0';
472  free(buf);
473  } else
474  dinfo->standard_name[0] = '\0';
475  // dimension
476  CHECK_ERROR( nc_inq_vardimid(ncid, varid, dimids) )
477 #ifdef NETCDF3
478  CHECK_ERROR( nc_inq_unlimdim(ncid, uldims) )
479  n = 1;
480 #else
481  CHECK_ERROR( nc_inq_unlimdims(ncid, &n, uldims) )
482 #endif
483  }
484  // attributes
485  dinfo->natts = 0;
486  CHECK_ERROR( nc_inq_varnatts(ncid, varid, &natts) )
487  for (i=0; i<natts; i++) {
488  CHECK_ERROR( nc_inq_attname(ncid, varid, i, name) )
489  if ( strcmp(name, "long_name") &&
490  strcmp(name, "description") &&
491  strcmp(name, "units") &&
492  strcmp(name, "standard_name") ) {
493  strcpy(dinfo->att_name+File_HSHORT*dinfo->natts, name);
494  CHECK_ERROR( nc_inq_atttype(ncid, varid, name, &xtype) )
495  NCTYPE2TYPE(xtype, dinfo->att_type[dinfo->natts]);
496  CHECK_ERROR( nc_inq_attlen(ncid, varid, name, &l) )
497  dinfo->att_len[dinfo->natts] = l;
498  dinfo->natts++;
499  }
500  }
501  }
502 
503  if ( rank > 0 ) {
504  tdim = -1;
505  for ( i=0; i<n; i++ ) {
506  if ( uldims[i] == dimids[0] ) {
507  tdim = uldims[i];
508  break;
509  }
510  }
511  if (rank > RANK_MAX) {
512  fprintf(stderr, "rank exceeds limit: %d\n", rank);
513  return ERROR_CODE;
514  }
515  dinfo->rank = tdim >= 0 ? rank -1 : rank; // do not count time dimension
516  // dim_name and dim_size
517  for (i=0; i<dinfo->rank; i++) {
518  // note: C and Fortran orders are opposite
519  if ( files[fid]->shared_mode ) {
520  MPI_Offset size_;
521  CHECK_PNC_ERROR( ncmpi_inq_dim(ncid, dimids[rank-i-1], name, &size_) )
522  size = (size_t)size_;
523  }
524  else
525  CHECK_ERROR( nc_inq_dim(ncid, dimids[rank-i-1], name, &size) )
526  if ( strlen(name) > File_HSHORT-1 ) {
527  fprintf(stderr, "Length of the dimension name (%s) is too long (should be < %d).\n", name, File_HSHORT);
528  return ERROR_CODE;
529  }
530  strncpy(dinfo->dim_name+i*File_HSHORT, name, File_HSHORT);
531  dinfo->dim_size[i] = size;
532  }
533 
534  dinfo->step = step;
535  if ( tdim >= 0 ) {
536  if ( files[fid]->shared_mode ) {
537  MPI_Offset idx[2];
538  MPI_Offset l;
539  // time_end
540  CHECK_PNC_ERROR( ncmpi_inq_dimname(ncid, tdim, name) )
541  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, name, &varid) )
542  idx[0] = step - 1;
543  CHECK_PNC_ERROR( ncmpi_get_var1_double_all(ncid, varid, idx, &(dinfo->time_end)) )
544  // units
545  CHECK_PNC_ERROR( ncmpi_inq_attlen (ncid, varid, "units", &l) )
546  buf = (char*) malloc(l+1);
547  CHECK_PNC_ERROR( ncmpi_get_att_text(ncid, varid, "units", buf) )
548  for (i=0; i<MIN(File_HMID-1,l); i++)
549  dinfo->time_units[i] = buf[i];
550  dinfo->time_units[i] = '\0';
551  free(buf);
552  // calendar
553 #ifdef PNETCDF
554  if ( ncmpi_inq_attlen(ncid, varid, "calendar", &l) == NC_NOERR ) {
555  buf = (char*) malloc(l+1);
556  CHECK_PNC_ERROR( ncmpi_get_att_text(ncid, varid, "calendar", buf) )
557  for (i=0; i<MIN(File_HSHORT-1,l); i++)
558  dinfo->calendar[i] = buf[i];
559  dinfo->calendar[i] = '\0';
560  free(buf);
561  } else
562  dinfo->calendar[0] = '\0';
563 #endif
564  // time_start
565  strcat(name, "_bnds");
566  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, name, &varid) )
567  idx[1] = 0;
568  CHECK_PNC_ERROR( ncmpi_get_var1_double_all(ncid, varid, idx, &(dinfo->time_start)) )
569  } else {
570  size_t idx[2];
571  size_t l;
572  // time_end
573  CHECK_ERROR( nc_inq_dimname(ncid, tdim, name) )
574  if ( nc_inq_varid(ncid, name, &varid) == NC_NOERR ) {
575  idx[0] = step - 1;
576  CHECK_ERROR( nc_get_var1_double(ncid, varid, idx, &(dinfo->time_end)) )
577  // units
578  CHECK_ERROR( nc_inq_attlen (ncid, varid, "units", &l) )
579  buf = (char*) malloc(l+1);
580  CHECK_ERROR( nc_get_att_text(ncid, varid, "units", buf) )
581  for (i=0; i<MIN(File_HMID-1,l); i++)
582  dinfo->time_units[i] = buf[i];
583  dinfo->time_units[i] = '\0';
584  free(buf);
585  // calendar
586  if ( nc_inq_attlen(ncid, varid, "calendar", &l) == NC_NOERR ) {
587  buf = (char*) malloc(l+1);
588  CHECK_ERROR( nc_get_att_text(ncid, varid, "calendar", buf) )
589  for (i=0; i<MIN(File_HSHORT-1,l); i++)
590  dinfo->calendar[i] = buf[i];
591  dinfo->calendar[i] = '\0';
592  free(buf);
593  } else
594  dinfo->calendar[0] = '\0';
595  // time_start
596  strcat(name, "_bnds");
597  CHECK_ERROR( nc_inq_varid(ncid, name, &varid) )
598  idx[1] = 0;
599  CHECK_ERROR( nc_get_var1_double(ncid, varid, idx, &(dinfo->time_start)) )
600  } else {
601  dinfo->time_start = 0.0;
602  dinfo->time_end = 0.0;
603  dinfo->time_units[0] = '\0';
604  dinfo->calendar[0] = '\0';
605  }
606  }
607  } else {
608  if ( step > 1 ) { // if variable does not have time dimention, step > 1 should not exist
609  fprintf(stderr, "requested step is larger than tdim: step=%d tdim=%d\n", step, tdim);
610  return ERROR_CODE;
611  }
612  dinfo->time_start = 0.0;
613  dinfo->time_end = 0.0;
614  dinfo->time_units[0] = '\0';
615  dinfo->calendar[0] = '\0';
616  }
617  } else {
618  dinfo->rank = 0;
619  dinfo->description[0] = '\0';
620  dinfo->units[0] = '\0';
621  dinfo->standard_name[0] = '\0';
622  dinfo->time_start = 0.0;
623  dinfo->time_end = 0.0;
624  dinfo->time_units[0] = '\0';
625  dinfo->calendar[0] = '\0';
626  }
627 
628  ERROR_SUPPRESS = 0;
629 
630  return SUCCESS_CODE;
631 }
632 
633 int32_t file_get_step_size_c( const int32_t fid, // (in)
634  const char* varname, // (in)
635  int32_t *len ) // (out)
636 {
637  int ncid, varid;
638 
639  int dimids[RANK_MAX], uldims[NC_MAX_DIMS], tdim;
640  int n, i;
641 
642  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
643  ncid = files[fid]->ncid;
644  if ( files[fid]->shared_mode )
645  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, varname, &varid) )
646  else
647  CHECK_ERROR( nc_inq_varid(ncid, varname, &varid) )
648 
649  if ( files[fid]->shared_mode ) {
650  CHECK_PNC_ERROR( ncmpi_inq_vardimid(ncid, varid, dimids) )
651  CHECK_PNC_ERROR( ncmpi_inq_unlimdim(ncid, uldims) )
652  n = uldims[0] < 0 ? 0 : 1;
653  } else {
654  CHECK_ERROR( nc_inq_vardimid(ncid, varid, dimids) )
655 #ifdef NETCDF3
656  CHECK_ERROR( nc_inq_unlimdim(ncid, uldims) )
657  n = uldims[0] < 0 ? 0 : 1;
658 #else
659  CHECK_ERROR( nc_inq_unlimdims(ncid, &n, uldims) )
660 #endif
661  }
662 
663  tdim = -1;
664  for ( i=0; i<n; i++ ) {
665  if ( uldims[i] == dimids[0] ) {
666  tdim = uldims[i];
667  break;
668  }
669  }
670 
671  if ( tdim > 0 ) {
672  if ( files[fid]->shared_mode ) {
673  MPI_Offset l;
674  CHECK_PNC_ERROR( ncmpi_inq_dimlen(ncid, tdim, &l) )
675  *len = l;
676  } else {
677  size_t l;
678  CHECK_ERROR( nc_inq_dimlen(ncid, tdim, &l) )
679  *len = l;
680  }
681  } else
682  *len = 0;
683 
684  return SUCCESS_CODE;
685 }
686 
687 int32_t file_read_data_c( void *var, // (out)
688  const datainfo_t *dinfo, // (in)
689  const int32_t precision, // (in)
690  const MPI_Offset ntypes, // (in)
691  const MPI_Datatype dtype, // (in)
692  const int32_t *start, // (in)
693  const int32_t *count ) // (in)
694 {
695  int ncid, varid;
696  int rank;
697  int i;
698  int fid;
699  size_t *str, *cnt;
700  MPI_Offset *strp, *cntp;
701  size_t size;
702  int l_rescale;
703 
704  fid = dinfo->fid;
705  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
706  ncid = files[fid]->ncid;
707  if ( files[fid]->shared_mode ) {
708  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, dinfo->varname, &varid) )
709  CHECK_PNC_ERROR( ncmpi_inq_varndims(ncid, varid, &rank) )
710  strp = (MPI_Offset*) malloc(sizeof(MPI_Offset)*rank);
711  cntp = (MPI_Offset*) malloc(sizeof(MPI_Offset)*rank);
712  } else {
713  CHECK_ERROR( nc_inq_varid(ncid, dinfo->varname, &varid) )
714  CHECK_ERROR( nc_inq_varndims(ncid, varid, &rank) )
715  }
716  str = (size_t*) malloc(sizeof(size_t)*rank);
717  cnt = (size_t*) malloc(sizeof(size_t)*rank);
718 
719  if ( start == NULL || count == NULL ) {
720  for (i=0; i<dinfo->rank; i++) {
721  // note: C and Fortran orders are opposite
722  str[rank -i-1] = 0;
723  cnt[rank -i-1] = dinfo->dim_size[i];
724  }
725  } else {
726  for (i=0; i<dinfo->rank; i++) {
727  // note: C and Fortran orders are opposite
728  str[rank -i-1] = start[i] - 1;
729  cnt[rank -i-1] = count[i];
730  }
731  }
732  if (rank > dinfo->rank) { // have time dimension
733  str[0] = dinfo->step - 1;
734  cnt[0] = 1;
735  }
736 
737  size = 1;
738  for (i=0; i<rank; i++) size *= cnt[i];
739 
740  if ( files[fid]->shared_mode ) {
741 #ifdef PNETCDF
742  for (i=0; i<rank; i++) {
743  strp[i] = (MPI_Offset) str[i];
744  cntp[i] = (MPI_Offset) cnt[i];
745  }
746  free(str);
747  free(cnt);
748  CHECK_PNC_ERROR( ncmpi_iget_vara(ncid, varid, strp, cntp, var, ntypes, dtype, NULL) )
749  free(strp);
750  free(cntp);
751  if ( dtype == MPI_FLOAT ) {
752  float factor, offset;
753  l_rescale = 0;
754  if ( ncmpi_get_att_float(ncid, varid, "scale_factor", &factor) != NC_NOERR )
755  factor = 1.0f;
756  else
757  l_rescale = 1;
758  if ( ncmpi_get_att_float(ncid, varid, "add_offset", &offset) != NC_NOERR )
759  offset = 0.0f;
760  else
761  l_rescale = 1;
762  if ( l_rescale ) for (i=0; i<size; i++) ((float*)var)[i] = ((float*)var)[i] * factor + offset;
763  } else if ( dtype == MPI_DOUBLE ) {
764  double factor, offset;
765  l_rescale = 0;
766  if ( ncmpi_get_att_double(ncid, varid, "scale_factor", &factor) != NC_NOERR )
767  factor = 1.0;
768  else
769  l_rescale = 1;
770  if ( ncmpi_get_att_double(ncid, varid, "add_offset", &offset) != NC_NOERR )
771  offset = 0.0;
772  else
773  l_rescale = 1;
774  if ( l_rescale ) for (i=0; i<size; i++) ((double*)var)[i] = ((double*)var)[i] * factor + offset;
775  } else {
776  float factor, offset;
777  if ( ( ncmpi_get_att_float(ncid, varid, "scale_factor", &factor) == NC_NOERR )
778  || ( ncmpi_get_att_float(ncid, varid, "add_offset", &offset) == NC_NOERR ) ) {
779  fprintf(stderr, "scale_factor and add_offset is not supported with a MPI derived type\n");
780  return ERROR_CODE;
781  }
782  }
783 #else
784  CHECK_PNC_ERROR( dummy )
785 #endif
786  } else {
787  switch ( precision ) {
788  case 8:
789  CHECK_ERROR( nc_get_vara_double(ncid, varid, str, cnt, (double*)var) )
790  {
791  double factor, offset;
792  l_rescale = 0;
793  if ( nc_get_att_double(ncid, varid, "scale_factor", &factor) != NC_NOERR )
794  factor = 1.0;
795  else
796  l_rescale = 1;
797  if ( nc_get_att_double(ncid, varid, "add_offset", &offset) != NC_NOERR )
798  offset = 0.0;
799  else
800  l_rescale = 1;
801  if ( l_rescale ) for (i=0; i<size; i++) ((double*)var)[i] = ((double*)var)[i] * factor + offset;
802  }
803  break;
804  case 4:
805  CHECK_ERROR( nc_get_vara_float(ncid, varid, str, cnt, (float*)var) )
806  {
807  float factor, offset;
808  l_rescale = 0;
809  if ( nc_get_att_float(ncid, varid, "scale_factor", &factor) != NC_NOERR )
810  factor = 1.0f;
811  else
812  l_rescale = 1;
813  if ( nc_get_att_float(ncid, varid, "add_offset", &offset) != NC_NOERR )
814  offset = 0.0f;
815  else
816  l_rescale = 1;
817  if ( l_rescale ) for (i=0; i<size; i++) ((float*)var)[i] = ((float*)var)[i] * factor + offset;
818  }
819  break;
820  default:
821  free(str);
822  free(cnt);
823  fprintf(stderr, "unsupported data precision: %d\n", precision );
824  return ERROR_CODE;
825  }
826  free(str);
827  free(cnt);
828  }
829 
830 
831  return SUCCESS_CODE;
832 }
833 
834 int32_t file_get_attribute_text_c( const int32_t fid, // (in)
835  const char *vname, // (in)
836  const char *key, // (in)
837  const int32_t suppress, // (in)
838  char *value, // (out)
839  const int32_t len) // (in)
840 {
841  int ncid;
842  int varid;
843 
844  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
845  ncid = files[fid]->ncid;
846 
847  ERROR_SUPPRESS = suppress;
848 
849  if ( files[fid]->shared_mode ) {
850  MPI_Offset l;
851  if ( strcmp(vname, "global") == 0 ) {
852  varid = NC_GLOBAL;
853  } else
854  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
855 
856  CHECK_PNC_ERROR( ncmpi_inq_attlen(ncid, varid, key, &l) )
857  if ( len < l ) return ERROR_CODE;
858 
859  CHECK_PNC_ERROR( ncmpi_get_att_text(ncid, varid, key, value) )
860  value[l] = '\0';
861  }
862  else {
863  size_t l;
864  if ( strcmp(vname, "global") == 0 ) {
865  varid = NC_GLOBAL;
866  } else
867  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
868 
869  CHECK_ERROR( nc_inq_attlen(ncid, varid, key, &l) )
870  if ( len < l ) return ERROR_CODE;
871 
872  CHECK_ERROR( nc_get_att_text(ncid, varid, key, value) )
873  value[l] = '\0';
874  }
875 
876  ERROR_SUPPRESS = 0;
877 
878  return SUCCESS_CODE;
879 }
880 
881 int32_t file_get_attribute_int_c( const int32_t fid, // (in)
882  const char *vname, // (in)
883  const char *key, // (in)
884  const int32_t suppress, // (in)
885  int *value, // (out)
886  const size_t len) // (in)
887 {
888  int ncid;
889  int varid;
890 
891  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
892  ncid = files[fid]->ncid;
893 
894  ERROR_SUPPRESS = suppress;
895 
896  if ( files[fid]->shared_mode ) {
897  MPI_Offset l;
898  if ( strcmp(vname, "global") == 0 ) {
899  varid = NC_GLOBAL;
900  } else
901  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
902 
903  CHECK_PNC_ERROR( ncmpi_inq_attlen(ncid, varid, key, &l) )
904  if ( len < l ) return ERROR_CODE;
905  CHECK_PNC_ERROR( ncmpi_get_att_int(ncid, varid, key, value) )
906  }
907  else {
908  size_t l;
909  if ( strcmp(vname, "global") == 0 ) {
910  varid = NC_GLOBAL;
911  } else
912  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
913 
914  CHECK_ERROR( nc_inq_attlen(ncid, varid, key, &l) )
915  if ( len < l ) return ERROR_CODE;
916  CHECK_ERROR( nc_get_att_int(ncid, varid, key, value) )
917  }
918 
919  ERROR_SUPPRESS = 0;
920 
921  return SUCCESS_CODE;
922 }
923 
924 int32_t file_get_attribute_float_c( const int32_t fid, // (in)
925  const char *vname, // (in)
926  const char *key, // (in)
927  const int32_t suppress, // (in)
928  float *value, // (out)
929  const size_t len) // (in)
930 {
931  int ncid;
932  int varid;
933 
934  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
935  ncid = files[fid]->ncid;
936 
937  ERROR_SUPPRESS = suppress;
938 
939  if ( files[fid]->shared_mode ) {
940  MPI_Offset l;
941  if ( strcmp(vname, "global") == 0 ) {
942  varid = NC_GLOBAL;
943  } else
944  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
945 
946  CHECK_PNC_ERROR( ncmpi_inq_attlen(ncid, varid, key, &l) )
947  if ( len < l ) return ERROR_CODE;
948  CHECK_PNC_ERROR( ncmpi_get_att_float(ncid, varid, key, value) )
949  }
950  else {
951  size_t l;
952  if ( strcmp(vname, "global") == 0 ) {
953  varid = NC_GLOBAL;
954  } else
955  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
956 
957  CHECK_ERROR( nc_inq_attlen(ncid, varid, key, &l) )
958  if ( len < l ) return ERROR_CODE;
959  CHECK_ERROR( nc_get_att_float(ncid, varid, key, value) )
960  }
961 
962  ERROR_SUPPRESS = 0;
963 
964  return SUCCESS_CODE;
965 }
966 
967 int32_t file_get_attribute_double_c( const int32_t fid, // (in)
968  const char *vname, // (in)
969  const char *key, // (in)
970  const int32_t suppress, // (in)
971  double *value, // (out)
972  const size_t len) // (in)
973 {
974  int ncid;
975  int varid;
976 
977  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
978  ncid = files[fid]->ncid;
979 
980  ERROR_SUPPRESS = suppress;
981 
982  if ( files[fid]->shared_mode ) {
983  MPI_Offset l;
984  if ( strcmp(vname, "global") == 0 ) {
985  varid = NC_GLOBAL;
986  } else
987  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
988 
989  CHECK_PNC_ERROR( ncmpi_inq_attlen(ncid, varid, key, &l) )
990  if ( len < l ) return ERROR_CODE;
991  CHECK_PNC_ERROR( ncmpi_get_att_double(ncid, varid, key, value) )
992  }
993  else {
994  size_t l;
995  if ( strcmp(vname, "global") == 0 ) {
996  varid = NC_GLOBAL;
997  } else
998  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
999 
1000  CHECK_ERROR( nc_inq_attlen(ncid, varid, key, &l) )
1001  if ( len < l ) return ERROR_CODE;
1002  CHECK_ERROR( nc_get_att_double(ncid, varid, key, value) )
1003  }
1004 
1005  ERROR_SUPPRESS = 0;
1006 
1007  return SUCCESS_CODE;
1008 }
1009 
1010 int32_t file_set_attribute_text_c( const int32_t fid, // (in)
1011  const char *vname, // (in)
1012  const char *key, // (in)
1013  const char *val) // (in)
1014 {
1015  int ncid;
1016  int varid;
1017  int attid;
1018  int32_t ret;
1019 
1020  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1021  ncid = files[fid]->ncid;
1022 
1023  if ( files[fid]->shared_mode ) {
1024  if ( strcmp(vname, "global") == 0 ) {
1025  varid = NC_GLOBAL;
1026  } else
1027  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
1028 
1029  if ( ncmpi_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1030  }
1031  else {
1032  if ( strcmp(vname, "global") == 0 ) {
1033  varid = NC_GLOBAL;
1034  } else
1035  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
1036 
1037  if ( nc_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1038  }
1039 
1040  ret = file_redef(fid, ncid);
1041  if ( ret != SUCCESS_CODE ) return ret;
1042 
1043  if ( files[fid]->shared_mode )
1044  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, key, strlen(val), val) )
1045  else
1046  CHECK_ERROR( nc_put_att_text(ncid, varid, key, strlen(val), val) )
1047 
1048  return SUCCESS_CODE;
1049 }
1050 
1051 int32_t file_set_attribute_int_c( const int32_t fid, // (in)
1052  const char *vname, // (in)
1053  const char *key, // (in)
1054  const int32_t *value, // (in)
1055  const size_t len ) // (in)
1056 {
1057  int ncid;
1058  int varid;
1059  int attid;
1060  int32_t ret;
1061 
1062  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1063  ncid = files[fid]->ncid;
1064 
1065  if ( files[fid]->shared_mode ) {
1066  if ( strcmp(vname, "global") == 0 ) {
1067  varid = NC_GLOBAL;
1068  } else
1069  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
1070 
1071  if ( ncmpi_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1072  }
1073  else {
1074  if ( strcmp(vname, "global") == 0 ) {
1075  varid = NC_GLOBAL;
1076  } else
1077  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
1078 
1079  if ( nc_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1080  }
1081 
1082  ret = file_redef(fid, ncid);
1083  if ( ret != SUCCESS_CODE ) return ret;
1084 
1085  if ( files[fid]->shared_mode )
1086  CHECK_PNC_ERROR( ncmpi_put_att_int(ncid, varid, key, NC_INT, len, value) )
1087  else
1088  CHECK_ERROR( nc_put_att_int(ncid, varid, key, NC_INT, len, value) )
1089 
1090 
1091  return SUCCESS_CODE;
1092 }
1093 
1094 int32_t file_set_attribute_float_c( const int32_t fid, // (in)
1095  const char *vname, // (in)
1096  const char *key, // (in)
1097  const float *value, // (in)
1098  const size_t len ) // (in)
1099 {
1100  int ncid;
1101  int varid;
1102  int attid;
1103  int32_t ret;
1104 
1105  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1106  ncid = files[fid]->ncid;
1107 
1108  if ( files[fid]->shared_mode ) {
1109  if ( strcmp(vname, "global") == 0 ) {
1110  varid = NC_GLOBAL;
1111  } else
1112  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
1113 
1114  if ( ncmpi_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1115  }
1116  else {
1117  if ( strcmp(vname, "global") == 0 ) {
1118  varid = NC_GLOBAL;
1119  } else
1120  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
1121 
1122  if ( nc_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1123  }
1124 
1125  ret = file_redef(fid, ncid);
1126  if ( ret != SUCCESS_CODE ) return ret;
1127 
1128  if ( files[fid]->shared_mode )
1129  CHECK_PNC_ERROR( ncmpi_put_att_float(ncid, varid, key, NC_FLOAT, len, value) )
1130  else
1131  CHECK_ERROR( nc_put_att_float(ncid, varid, key, NC_FLOAT, len, value) )
1132 
1133  return SUCCESS_CODE;
1134 }
1135 
1136 int32_t file_set_attribute_double_c( const int32_t fid, // (in)
1137  const char *vname, // (in)
1138  const char *key, // (in)
1139  const double *value, // (in)
1140  const size_t len ) // (in)
1141 {
1142  int ncid;
1143  int varid;
1144  int attid;
1145 
1146  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1147  ncid = files[fid]->ncid;
1148 
1149  if ( files[fid]->shared_mode ) {
1150  if ( strcmp(vname, "global") == 0 ) {
1151  varid = NC_GLOBAL;
1152  } else
1153  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
1154 
1155  if ( ncmpi_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1156  }
1157  else {
1158  if ( strcmp(vname, "global") == 0 ) {
1159  varid = NC_GLOBAL;
1160  } else
1161  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
1162 
1163  if ( nc_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1164  }
1165 
1166  if ( files[fid]->shared_mode )
1167  CHECK_PNC_ERROR( ncmpi_put_att_double(ncid, varid, key, NC_DOUBLE, len, value) )
1168  else
1169  CHECK_ERROR( nc_put_att_double(ncid, varid, key, NC_DOUBLE, len, value) )
1170 
1171  return SUCCESS_CODE;
1172 }
1173 
1174 int32_t file_add_associatedvariable_c( const int32_t fid, // (in)
1175  const char *vname) // (in)
1176 {
1177  int ncid, varid;
1178  int32_t ret;
1179 
1180  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1181  ncid = files[fid]->ncid;
1182 
1183  if ( nc_inq_varid(ncid, vname, &varid) == NC_NOERR ) // check if existed
1184  return ALREADY_EXISTED_CODE;
1185 
1186  ret = file_redef(fid, ncid);
1187  if ( ret != SUCCESS_CODE ) return ret;
1188 
1189  if ( files[fid]->shared_mode )
1190  CHECK_PNC_ERROR( ncmpi_def_var(ncid, vname, NC_INT, 0, 0, &varid) )
1191  else
1192  CHECK_ERROR( nc_def_var(ncid, vname, NC_INT, 0, 0, &varid) )
1193 
1194  return SUCCESS_CODE;
1195 }
1196 
1197 int32_t file_set_tunits_c( const int32_t fid, // (in)
1198  const char *time_units, // (in)
1199  const char *calendar) // (in)
1200 {
1201  strcpy(files[fid]->time_units, time_units);
1202  strcpy(files[fid]->calendar, calendar);
1203 
1204  return SUCCESS_CODE;
1205 }
1206 
1207 int32_t file_put_axis_c( const int32_t fid, // (in)
1208  const char *name, // (in)
1209  const char *desc, // (in)
1210  const char *units, // (in)
1211  const char *dim_name, // (in)
1212  const int32_t dtype, // (in)
1213  const void* val, // (in)
1214  const int32_t size, // (in)
1215  const int32_t precision) // (in)
1216 {
1217  int ncid, dimid, varid;
1218  nc_type xtype = -1;
1219  int32_t ret;
1220 
1221  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1222  ncid = files[fid]->ncid;
1223 
1224  if ( nc_inq_varid(ncid, name, &varid) == NC_NOERR ) // check if existed
1225  return ALREADY_EXISTED_CODE;
1226 
1227  ret = file_redef(fid, ncid);
1228  if ( ret != SUCCESS_CODE ) return ret;
1229 
1230  if ( nc_inq_dimid(ncid, dim_name, &dimid) != NC_NOERR ) // check if existed
1231  CHECK_ERROR( nc_def_dim(ncid, dim_name, size, &dimid) )
1232 
1233  TYPE2NCTYPE(dtype, xtype);
1234  CHECK_ERROR( nc_def_var(ncid, name, xtype, 1, &dimid, &varid) )
1235  if ( strlen(desc)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1236  if ( strlen(units)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "units", strlen(units), units) )
1237 
1238  ret = file_enddef(fid, ncid);
1239  if ( ret != SUCCESS_CODE ) return ret;
1240 
1241  switch ( precision ) {
1242  case 8:
1243  CHECK_ERROR( nc_put_var_double(ncid, varid, (double*)val) )
1244  break;
1245  case 4:
1246  CHECK_ERROR( nc_put_var_float(ncid, varid, (float*)val) )
1247  break;
1248  default:
1249  fprintf(stderr, "unsupported data precision: %d\n", precision);
1250  return ERROR_CODE;
1251  }
1252 
1253  return SUCCESS_CODE;
1254 }
1255 
1256 int32_t file_def_axis_c( const int32_t fid, // (in)
1257  const char *name, // (in)
1258  const char *desc, // (in)
1259  const char *units, // (in)
1260  const char *dim_name, // (in)
1261  const int32_t dtype, // (in)
1262  const int32_t dim_size, // (in)
1263  const int32_t bounds) // (in)
1264 {
1265  int ncid, dimid, varid;
1266  nc_type xtype = -1;
1267  int dimids[2];
1268  char buf[File_HSHORT+6];
1269  int32_t ret;
1270 
1271  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1272  ncid = files[fid]->ncid;
1273 
1274  if ( files[fid]->shared_mode ) {
1275  if ( ncmpi_inq_varid(ncid, name, &varid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1276  } else {
1277  if ( nc_inq_varid(ncid, name, &varid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1278  }
1279 
1280  ret = file_redef(fid, ncid);
1281  if ( ret != SUCCESS_CODE ) return ret;
1282 
1283  TYPE2NCTYPE(dtype, xtype);
1284  if ( files[fid]->shared_mode ) {
1285  if ( ncmpi_inq_dimid(ncid, dim_name, &dimid) != NC_NOERR ) // check if existed
1286  CHECK_PNC_ERROR( ncmpi_def_dim(ncid, dim_name, dim_size, &dimid) )
1287 
1288  CHECK_PNC_ERROR( ncmpi_def_var(ncid, name, xtype, 1, &dimid, &varid) )
1289  if ( strlen(desc)>0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1290  if ( strlen(units)>0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "units", strlen(units), units) )
1291 
1292  if ( bounds ) {
1293  dimids[0] = dimid;
1294  if ( ncmpi_inq_dimid(ncid, "nv", &(dimids[1])) != NC_NOERR ) // first called
1295  CHECK_PNC_ERROR( ncmpi_def_dim(ncid, "nv", 2, &(dimids[1])) )
1296  sprintf(buf, "%s_bnds", dim_name);
1297  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "bounds", strlen(buf), buf) )
1298  CHECK_PNC_ERROR( ncmpi_def_var(ncid, buf, NC_DOUBLE, 2, dimids, &varid) )
1299  }
1300  } else {
1301  if ( nc_inq_dimid(ncid, dim_name, &dimid) != NC_NOERR ) // check if existed
1302  CHECK_ERROR( nc_def_dim(ncid, dim_name, dim_size, &dimid) )
1303 
1304  CHECK_ERROR( nc_def_var(ncid, name, xtype, 1, &dimid, &varid) )
1305  if ( strlen(desc)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1306  if ( strlen(units)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "units", strlen(units), units) )
1307 
1308  if ( bounds ) {
1309  dimids[0] = dimid;
1310  if ( nc_inq_dimid(ncid, "nv", &(dimids[1])) != NC_NOERR ) // first called
1311  CHECK_ERROR( nc_def_dim(ncid, "nv", 2, &(dimids[1])) )
1312  sprintf(buf, "%s_bnds", dim_name);
1313  CHECK_ERROR( nc_put_att_text(ncid, varid, "bounds", strlen(buf), buf) )
1314  CHECK_ERROR( nc_def_var(ncid, buf, NC_DOUBLE, 2, dimids, &varid) )
1315  }
1316  }
1317 
1318  return SUCCESS_CODE;
1319 }
1320 
1321 int32_t file_write_axis_c( const int32_t fid, // (in)
1322  const char *name, // (in)
1323  const void *val, // (in)
1324  const int32_t precision, // (in)
1325  const MPI_Offset *start, // (in)
1326  const MPI_Offset *count) // (in)
1327 {
1328  int ncid, varid;
1329  int32_t ret;
1330 
1331  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1332  ncid = files[fid]->ncid;
1333 
1334  if ( files[fid]->shared_mode )
1335  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, name, &varid) )
1336  else
1337  CHECK_ERROR( nc_inq_varid(ncid, name, &varid) )
1338 
1339  ret = file_enddef(fid, ncid);
1340  if ( ret != SUCCESS_CODE ) return ret;
1341 
1342  switch ( precision ) {
1343  case 8:
1344  if ( files[fid]->shared_mode )
1345  CHECK_PNC_ERROR( ncmpi_iput_vara_double(ncid, varid, start, count, val, NULL) )
1346  else
1347  CHECK_ERROR( nc_put_var_double(ncid, varid, (double*)val) )
1348  break;
1349  case 4:
1350  if ( files[fid]->shared_mode )
1351  CHECK_PNC_ERROR( ncmpi_iput_vara_float(ncid, varid, start, count, val, NULL) )
1352  else
1353  CHECK_ERROR( nc_put_var_float(ncid, varid, (float*)val) )
1354  break;
1355  default:
1356  fprintf(stderr, "unsupported data precision: %d\n", precision);
1357  return ERROR_CODE;
1358  }
1359 
1360  return SUCCESS_CODE;
1361 }
1362 
1363 int32_t file_put_associatedcoordinate_c( const int32_t fid, // (in)
1364  const char *name, // (in)
1365  const char *desc, // (in)
1366  const char *units, // (in)
1367  const char **dim_names, // (in)
1368  const int32_t ndims, // (in)
1369  const int32_t dtype, // (in)
1370  const void* val, // (in)
1371  const int32_t precision) // (in)
1372 {
1373  int ncid, *dimids, varid;
1374  nc_type xtype = -1;
1375  int i;
1376  int32_t ret;
1377 
1378  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1379  ncid = files[fid]->ncid;
1380 
1381  if ( nc_inq_varid(ncid, name, &varid) == NC_NOERR ) // check if existed
1382  return ALREADY_EXISTED_CODE;
1383 
1384  ret = file_redef(fid, ncid);
1385  if ( ret != SUCCESS_CODE ) return ret;
1386 
1387  dimids = malloc(sizeof(int)*ndims);
1388  for (i=0; i<ndims; i++)
1389  CHECK_ERROR( nc_inq_dimid(ncid, dim_names[i], dimids+ndims-i-1) )
1390 
1391  TYPE2NCTYPE(dtype, xtype);
1392 
1393  CHECK_ERROR( nc_def_var(ncid, name, xtype, ndims, dimids, &varid) )
1394  if ( strlen(desc)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1395  if ( strlen(units)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "units", strlen(units), units) )
1396  free(dimids);
1397 
1398  ret = file_enddef(fid, ncid);
1399  if ( ret != SUCCESS_CODE ) return ret;
1400 
1401  switch ( precision ) {
1402  case 8:
1403  CHECK_ERROR( nc_put_var_double(ncid, varid, (double*)val) )
1404  break;
1405  case 4:
1406  CHECK_ERROR( nc_put_var_float(ncid, varid, (float*)val) )
1407  break;
1408  default:
1409  fprintf(stderr, "unsupported data precision: %d\n", precision);
1410  return ERROR_CODE;
1411  }
1412 
1413  return SUCCESS_CODE;
1414 }
1415 
1416 int32_t file_def_associatedcoordinate_c( const int32_t fid, // (in)
1417  const char *name, // (in)
1418  const char *desc, // (in)
1419  const char *units, // (in)
1420  const char **dim_names, // (in)
1421  const int32_t ndims, // (in)
1422  const int32_t dtype) // (in)
1423 {
1424  int ncid, varid;
1425  nc_type xtype = -1;
1426  int i;
1427  int32_t ret;
1428  int dimids[ndims];
1429 
1430  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1431  ncid = files[fid]->ncid;
1432 
1433  if ( files[fid]->shared_mode ) {
1434  if ( ncmpi_inq_varid(ncid, name, &varid) == NC_NOERR ) // check if existed
1435  return ALREADY_EXISTED_CODE;
1436  } else {
1437  if ( nc_inq_varid(ncid, name, &varid) == NC_NOERR ) // check if existed
1438  return ALREADY_EXISTED_CODE;
1439  }
1440 
1441  ret = file_redef(fid, ncid);
1442  if ( ret != SUCCESS_CODE ) return ret;
1443 
1444  TYPE2NCTYPE(dtype, xtype);
1445 
1446  if ( files[fid]->shared_mode ) {
1447  for (i=0; i<ndims; i++)
1448  CHECK_PNC_ERROR( ncmpi_inq_dimid(ncid, dim_names[i], dimids+ndims-i-1) )
1449 
1450  CHECK_PNC_ERROR( ncmpi_def_var(ncid, name, xtype, ndims, dimids, &varid) )
1451  if ( strlen(desc) >0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1452  if ( strlen(units)>0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "units", strlen(units), units) )
1453  }
1454  else {
1455  for (i=0; i<ndims; i++)
1456  CHECK_ERROR( nc_inq_dimid(ncid, dim_names[i], dimids+ndims-i-1) )
1457 
1458  CHECK_ERROR( nc_def_var(ncid, name, xtype, ndims, dimids, &varid) )
1459  if ( strlen(desc) >0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1460  if ( strlen(units)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "units", strlen(units), units) )
1461  }
1462 
1463  return SUCCESS_CODE;
1464 }
1465 
1466 int32_t file_write_associatedcoordinate_c( const int32_t fid, // (in)
1467  const char *name, // (in)
1468  const void* val, // (in)
1469  const int32_t precision, // (in)
1470  const MPI_Offset *start, // (in)
1471  const MPI_Offset *count) // (in)
1472 {
1473  int ncid, varid;
1474  int32_t ret;
1475 
1476  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1477  ncid = files[fid]->ncid;
1478 
1479  if ( files[fid]->shared_mode )
1480  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, name, &varid) )
1481  else
1482  CHECK_ERROR( nc_inq_varid(ncid, name, &varid) )
1483 
1484  ret = file_enddef(fid, ncid);
1485  if ( ret != SUCCESS_CODE ) return ret;
1486 
1487  switch ( precision ) {
1488  case 8:
1489  if ( files[fid]->shared_mode )
1490  CHECK_PNC_ERROR( ncmpi_iput_vara_double(ncid, varid, start, count, (double*)val, NULL) )
1491  else
1492  CHECK_ERROR( nc_put_var_double(ncid, varid, (double*)val) )
1493  break;
1494  case 4:
1495  if ( files[fid]->shared_mode )
1496  CHECK_PNC_ERROR( ncmpi_iput_vara_float(ncid, varid, start, count, (float*)val, NULL) )
1497  else
1498  CHECK_ERROR( nc_put_var_float(ncid, varid, (float*)val) )
1499  break;
1500  default:
1501  fprintf(stderr, "unsupported data precision: %d\n", precision);
1502  return ERROR_CODE;
1503  }
1504 
1505  return SUCCESS_CODE;
1506 }
1507 
1508 int32_t file_add_variable_c( const int32_t fid, // (in)
1509  const char *varname, // (in)
1510  const char *desc, // (in)
1511  const char *units, // (in)
1512  const char *stdname, // (in)
1513  const char **dims, // (in)
1514  const int32_t ndims, // (in)
1515  const int32_t dtype, // (in)
1516  const real64_t tint, // (in)
1517  const int32_t tavg, // (in)
1518  int32_t *vid) // (out)
1519 {
1520  int ncid, varid, acid, *acdimids;
1521  int dimids[NC_MAX_DIMS], dimid;
1522  char tname[File_HSHORT+1];
1523  int tdimid, tvarid;
1524  nc_type xtype = -1;
1525  char buf[File_HMID+1];
1526  int i, j, k, m, err;
1527  int ndims_t, nndims;
1528  size_t size;
1529  double rmiss = RMISS;
1530  char coord[File_HMID+1];
1531  int has_assoc;
1532  int new;
1533  int32_t ret;
1534 
1535  if ( nvar >= VAR_MAX ) {
1536  fprintf(stderr, "exceed max number of variable limit\n");
1537  return ERROR_CODE;
1538  }
1539 
1540  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1541  ncid = files[fid]->ncid;
1542 
1543  vars[nvar] = (varinfo_t*) malloc(sizeof(varinfo_t));
1544  vars[nvar]->ncid = ncid;
1545  vars[nvar]->t = NULL;
1546  vars[nvar]->start = NULL;
1547  vars[nvar]->count = NULL;
1548  vars[nvar]->ndims = ndims;
1549 
1550  ret = file_redef(fid, ncid);
1551  if ( ret != SUCCESS_CODE ) return ret;
1552 
1553  // get time variable
1554  if ( tint > 0.0 ) {
1555  for ( i=0; i<nt; i++ ) {
1556  if ( tdims[i] != NULL && // still opened
1557  tdims[i]->ncid == ncid && // same file
1558  tdims[i]->tint == tint ) { // same time interval
1559  vars[nvar]->t = tdims[i];
1560  break;
1561  }
1562  }
1563  if ( vars[nvar]->t == NULL ) {
1564  tdims[nt] = (tdim_t*) malloc(sizeof(tdim_t));
1565  tdims[nt]->ncid = ncid;
1566  tdims[nt]->count = -1;
1567  tdims[nt]->tint = tint;
1568  tdims[nt]->tval = (double*) malloc(sizeof(double)*NTMAX);
1569  // generate name
1570  m=0;
1571  for (i=0; i<nt; i++) {
1572  if ( tdims[i] != NULL && tdims[i]->ncid == ncid ) m++;
1573  }
1574  if ( m == 0 ) {
1575  strcpy(tname, "time");
1576  } else {
1577  sprintf(tname, "time%d", m);
1578  }
1579  strcpy(tdims[nt]->name, tname);
1580  // define time dimension and variable
1581  if ( files[fid]->shared_mode ) {
1582  CHECK_PNC_ERROR( ncmpi_def_dim(ncid, tname, 0, &tdimid) )
1583  tdims[nt]->dimid = tdimid;
1584  CHECK_PNC_ERROR( ncmpi_def_var(ncid, tname, NC_DOUBLE, 1, &tdimid, &tvarid) )
1585  tdims[nt]->varid = tvarid;
1586  strcpy(buf, "time");
1587  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, tvarid, "long_name", strlen(buf), buf) )
1588  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, tvarid, "units", strlen(files[fid]->time_units), files[fid]->time_units) )
1589  if ( strlen(files[fid]->calendar) > 0 )
1590  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, tvarid, "calendar", strlen(files[fid]->calendar), files[fid]->calendar) )
1591  // define boundary variable
1592  if ( ncmpi_inq_dimid(ncid, "nv", &(dimids[1])) != NC_NOERR ) // first called
1593  CHECK_PNC_ERROR( ncmpi_def_dim(ncid, "nv", 2, &(dimids[1])) )
1594  sprintf(buf, "%s_bnds", tname);
1595  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, tvarid, "bounds", strlen(buf), buf) )
1596  dimids[0] = tdimid;
1597  CHECK_PNC_ERROR( ncmpi_def_var(ncid, buf, NC_DOUBLE, 2, dimids, &tvarid) )
1598  tdims[nt]->bndsid = tvarid;
1599  //CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, tvarid, "units", strlen(files[fid]->time_units), files[fid]->time_units) )
1600  }
1601  else {
1602  CHECK_ERROR( nc_def_dim(ncid, tname, 0, &tdimid) )
1603  tdims[nt]->dimid = tdimid;
1604  CHECK_ERROR( nc_def_var(ncid, tname, NC_DOUBLE, 1, &tdimid, &tvarid) )
1605  tdims[nt]->varid = tvarid;
1606  strcpy(buf, "time");
1607  CHECK_ERROR( nc_put_att_text(ncid, tvarid, "long_name", strlen(buf), buf) )
1608  CHECK_ERROR( nc_put_att_text(ncid, tvarid, "units", strlen(files[fid]->time_units), files[fid]->time_units) )
1609  if ( strlen(files[fid]->calendar) > 0 )
1610  CHECK_ERROR( nc_put_att_text(ncid, tvarid, "calendar", strlen(files[fid]->calendar), files[fid]->calendar) )
1611  // define boundary variable
1612  if ( nc_inq_dimid(ncid, "nv", &(dimids[1])) != NC_NOERR ) // first called
1613  CHECK_ERROR( nc_def_dim(ncid, "nv", 2, &(dimids[1])) )
1614  sprintf(buf, "%s_bnds", tname);
1615  CHECK_ERROR( nc_put_att_text(ncid, tvarid, "bounds", strlen(buf), buf) )
1616  dimids[0] = tdimid;
1617  CHECK_ERROR( nc_def_var(ncid, buf, NC_DOUBLE, 2, dimids, &tvarid) )
1618  tdims[nt]->bndsid = tvarid;
1619  //CHECK_ERROR( nc_put_att_text(ncid, tvarid, "units", strlen(files[fid]->time_units), files[fid]->time_units) )
1620  }
1621 
1622  vars[nvar]->t = tdims[nt];
1623  nt++;
1624  }
1625  }
1626 
1627  // get dimension IDs
1628  // note: C and Fortran order are opposite
1629  ndims_t = ndims;
1630  if ( tint > 0.0 ) { // add time dimension
1631  dimids[0] = vars[nvar]->t->dimid;
1632  ndims_t++;
1633  }
1634  for (i=ndims_t-ndims; i<ndims_t; i++) dimids[i] = -1;
1635 
1636  has_assoc = 0;
1637  nndims = 0;
1638  for (i=0; i<ndims; i++) {
1639  //printf("%d %s\n", i, dims[i]);
1640  if ( files[fid]->shared_mode )
1641  err = ncmpi_inq_dimid(ncid, dims[i], &dimid);
1642  else
1643  err = nc_inq_dimid(ncid, dims[i], &dimid);
1644  if ( err == NC_NOERR ) {
1645  //printf("not assoc\n");
1646  new = 1;
1647  for (k=0; k<nndims; k++) {
1648  if (dimid == dimids[k]) {
1649  new = 0;
1650  break;
1651  }
1652  }
1653  if (new) {
1654  dimids[ndims_t-(++nndims)] = dimid;
1655  }
1656  } else {
1657  //printf("assoc\n");
1658  if ( files[fid]->shared_mode ) {
1659  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, dims[i], &acid) )
1660  CHECK_PNC_ERROR( ncmpi_inq_varndims(ncid, acid, &m) )
1661  acdimids = (int*) malloc((sizeof(int)*m));
1662  CHECK_PNC_ERROR( ncmpi_inq_vardimid(ncid, acid, acdimids) )
1663  }
1664  else {
1665  CHECK_ERROR( nc_inq_varid(ncid, dims[i], &acid) )
1666  CHECK_ERROR( nc_inq_varndims(ncid, acid, &m) )
1667  acdimids = (int*) malloc((sizeof(int)*m));
1668  CHECK_ERROR( nc_inq_vardimid(ncid, acid, acdimids) )
1669  }
1670  for (j=m-1; j>=0; j--) {
1671  new = 1;
1672  for (k=0; k<ndims_t; k++) {
1673  if (acdimids[j] == dimids[k]) {
1674  new = 0;
1675  break;
1676  }
1677  }
1678  if (new) {
1679  if ( nndims >= ndims_t ) {
1680  fprintf(stderr, "Error: invalid associated coordinates\n");
1681  return ERROR_CODE;
1682  }
1683  dimids[ndims_t-(++nndims)] = acdimids[j];
1684  //nc_inq_dimname(ncid, acdimids[j], tname);
1685  //printf("add %s\n", tname);
1686  }
1687  }
1688  free(acdimids);
1689  has_assoc = 1;
1690  }
1691  }
1692  if (nndims != ndims) {
1693  fprintf(stderr, "Error: invalid associated coordinates: %d %d\n", ndims_t, nndims);
1694  return ERROR_CODE;
1695  }
1696 
1697  TYPE2NCTYPE(dtype, xtype);
1698  if ( files[fid]->shared_mode ) {
1699  CHECK_PNC_ERROR( ncmpi_def_var(ncid, varname, xtype, ndims_t, dimids, &varid) )
1700  // put variable attribute
1701  if ( strlen(desc) >0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1702  if ( strlen(units) >0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "units", strlen(units), units) )
1703  if ( strlen(stdname)>0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "standard_name", strlen(stdname), stdname) )
1704 
1705 // CHECK_PNC_ERROR( ncmpi_put_att_double(ncid, varid, _FillValue, xtype, 1, &rmiss) )
1706  CHECK_PNC_ERROR( ncmpi_put_att_double(ncid, varid, "missing_value", xtype, 1, &rmiss) )
1707  }
1708  else {
1709  CHECK_ERROR( nc_def_var(ncid, varname, xtype, ndims_t, dimids, &varid) )
1710  // put variable attribute
1711  if ( strlen(desc) >0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1712  if ( strlen(units) >0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "units", strlen(units), units) )
1713  if ( strlen(stdname)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "standard_name", strlen(stdname), stdname) )
1714  CHECK_ERROR( nc_put_att_double(ncid, varid, _FillValue, xtype, 1, &rmiss) )
1715  CHECK_ERROR( nc_put_att_double(ncid, varid, "missing_value", xtype, 1, &rmiss) )
1716  }
1717  if ( has_assoc ) {
1718  strcpy(coord, dims[0]);
1719  for(i=1; i<ndims; i++) {
1720  if (strlen(coord)+strlen(dims[i])+1 < File_HMID) {
1721  strcat(coord, " ");
1722  strcat(coord, dims[i]);
1723  }
1724  }
1725  if ( ndims_t > ndims && strlen(coord)+6 < File_HMID) {
1726  strcat(coord, " ");
1727  strcat(coord, vars[nvar]->t->name);
1728  }
1729  if ( files[fid]->shared_mode )
1730  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "coordinates", strlen(coord), coord) )
1731  else
1732  CHECK_ERROR( nc_put_att_text(ncid, varid, "coordinates", strlen(coord), coord) )
1733  }
1734 
1735 
1736  if ( tavg ) {
1737  sprintf(buf, "%s: mean", vars[nvar]->t->name);
1738  if ( files[fid]->shared_mode )
1739  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "cell_methods", strlen(buf), buf) )
1740  else
1741  CHECK_ERROR( nc_put_att_text(ncid, varid, "cell_methods", strlen(buf), buf) )
1742  }
1743 
1744  // set start and count
1745  vars[nvar]->ndims_t = ndims_t;
1746  vars[nvar]->start = (size_t*) malloc(sizeof(size_t)*ndims_t);
1747  vars[nvar]->count = (size_t*) malloc(sizeof(size_t)*ndims_t);
1748  for ( i=0; i<ndims_t; i++ ) {
1749  if ( files[fid]->shared_mode ) {
1750  MPI_Offset dimlen;
1751  CHECK_PNC_ERROR( ncmpi_inq_dimlen(ncid, dimids[i], &dimlen) )
1752  size = (size_t) dimlen;
1753  }
1754  else
1755  CHECK_ERROR( nc_inq_dimlen(ncid, dimids[i], &size) )
1756  vars[nvar]->count[i] = size;
1757  vars[nvar]->start[i] = 0;
1758  }
1759  if ( tint > 0.0 ) vars[nvar]->count[0] = 1;
1760 
1761 #ifndef NETCDF3
1762  // set chunk size and deflate level (NetCDF-4 only)
1763  if ( ! files[fid]->shared_mode && files[fid]->deflate_level > 0 ) {
1764  CHECK_ERROR( nc_def_var_chunking(ncid, varid, NC_CHUNKED, vars[nvar]->count) )
1765  CHECK_ERROR( nc_def_var_deflate(ncid, varid, 0, 1, files[fid]->deflate_level) )
1766  }
1767 #endif
1768 
1769  vars[nvar]->varid = varid;
1770  *vid = nvar;
1771  nvar++;
1772 
1773  return SUCCESS_CODE;
1774 }
1775 
1776 int32_t file_enddef_c( const int32_t fid ) // (in)
1777 {
1778  int ncid;
1779 
1780  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1781 
1782  ncid = files[fid]->ncid;
1783 
1784  return file_enddef(fid, ncid);
1785 }
1786 
1787 int32_t file_attach_buffer_c( const int32_t fid, // (in)
1788  const int64_t buf_amount ) // (in)
1789 {
1790  int ncid;
1791 
1792  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1793  ncid = files[fid]->ncid;
1794 
1795  if ( files[fid]->shared_mode )
1796  CHECK_PNC_ERROR( ncmpi_buffer_attach(ncid, (MPI_Offset)buf_amount) )
1797 
1798  return SUCCESS_CODE;
1799 }
1800 
1801 int32_t file_detach_buffer_c( const int32_t fid ) // (in)
1802 {
1803  int ncid;
1804 
1805  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1806  ncid = files[fid]->ncid;
1807 
1808  if ( files[fid]->shared_mode )
1809  CHECK_PNC_ERROR( ncmpi_buffer_detach(ncid) )
1810 
1811  return SUCCESS_CODE;
1812 }
1813 
1814 int32_t file_flush_c( const int32_t fid ) // (in)
1815 {
1816  int ncid;
1817 
1818  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1819  ncid = files[fid]->ncid;
1820 
1821  if ( files[fid]->shared_mode )
1822  CHECK_PNC_ERROR( ncmpi_wait_all(ncid, NC_REQ_ALL, NULL, NULL) )
1823  else
1824  CHECK_ERROR( nc_sync(ncid) )
1825 
1826  return SUCCESS_CODE;
1827 }
1828 
1829 int32_t file_write_data_c( const int32_t fid, // (in)
1830  const int32_t vid, // (in)
1831  const void *var, // (in)
1832  const real64_t t_start, // (in)
1833  const real64_t t_end, // (in)
1834  const int32_t precision, // (in)
1835  const int32_t ndims, // (in)
1836  const int32_t *start, // (in)
1837  const int32_t *count) // (in)
1838 {
1839  int ncid, varid;
1840  MPI_Offset *str, *cnt;
1841  int32_t ret;
1842 
1843  if ( vars[vid] == NULL ) return ALREADY_CLOSED_CODE;
1844  ncid = vars[vid]->ncid;
1845 
1846  if ( ndims != vars[vid]->ndims ) {
1847  fprintf(stderr, "Error: at line %d in %s\n", __LINE__, __FILE__);
1848  fprintf(stderr, " dimension size %d is not consistent that was added by file_add_variable %d\n", ndims, (int)vars[vid]->ndims );
1849  return ERROR_CODE;
1850  }
1851 
1852  ret = file_enddef(fid, ncid);
1853  if ( ret != SUCCESS_CODE ) return ret;
1854 
1855  varid = vars[vid]->varid;
1856  if ( vars[vid]->t != NULL ) { // have time dimension
1857  if ( vars[vid]->t->count < 0 || // first time
1858  t_end > vars[vid]->t->t + TEPS ) { // time goes next step
1859  vars[vid]->t->count += 1;
1860  vars[vid]->t->t = t_end;
1861  if ( vars[vid]->t->count > NTMAX-1 ) {
1862  fprintf(stderr, "time count exceeds the max limit (%d)\n", NTMAX);
1863  return ERROR_CODE;
1864  }
1865  vars[vid]->t->tval[vars[vid]->t->count] = t_end;
1866  if ( files[fid]->shared_mode ) { // write a new value to variable time
1867  MPI_Offset index[2];
1868  index[0] = (MPI_Offset) vars[vid]->t->count;
1869  CHECK_PNC_ERROR( ncmpi_put_var1_double_all(ncid, vars[vid]->t->varid, index, &t_end) )
1870  index[1] = 0;
1871  CHECK_PNC_ERROR( ncmpi_put_var1_double_all(ncid, vars[vid]->t->bndsid, index, &t_start ) )
1872  index[1] = 1;
1873  CHECK_PNC_ERROR( ncmpi_put_var1_double_all(ncid, vars[vid]->t->bndsid, index, &t_end ) )
1874  } else {
1875  size_t index[2];
1876  index[0] = vars[vid]->t->count;
1877  CHECK_ERROR( nc_put_var1_double(ncid, vars[vid]->t->varid, index, &t_end) )
1878  index[1] = 0;
1879  CHECK_ERROR( nc_put_var1_double(ncid, vars[vid]->t->bndsid, index, &t_start) )
1880  index[1] = 1;
1881  CHECK_ERROR( nc_put_var1_double(ncid, vars[vid]->t->bndsid, index, &t_end) )
1882  }
1883  vars[vid]->start[0] = vars[vid]->t->count;
1884  } else {
1885  size_t nt = vars[vid]->t->count + 1;
1886  int flag, n;
1887  flag = 1;
1888  for(n=nt-1;n>=0;n--) {
1889  if ( fabs(vars[vid]->t->tval[n]-t_end) < TEPS ) {
1890  vars[vid]->start[0] = n;
1891  flag = 0;
1892  break;
1893  }
1894  }
1895  if ( flag ) {
1896  fprintf(stderr, "cannot find time: %f\n", t_end);
1897  fprintf(stderr, " time count is : %d, last time is: %f, diff is: %e\n", vars[vid]->t->count < 0, vars[vid]->t->t, vars[vid]->t->t-t_end);
1898  fprintf(stderr, " time is: ");
1899  for (n=0;n<nt;n++) fprintf(stderr, "%f, ", vars[vid]->t->tval[n]);
1900  fprintf(stderr, "\n");
1901  return ERROR_CODE;
1902  }
1903  }
1904  }
1905 
1906  if ( files[fid]->shared_mode ) {
1907  int i;
1908  int ndims_t = vars[vid]->ndims_t;
1909  str = (MPI_Offset*) malloc(sizeof(MPI_Offset)*(ndims_t));
1910  cnt = (MPI_Offset*) malloc(sizeof(MPI_Offset)*(ndims_t));
1911  if ( vars[vid]->t != NULL ) { // have time dimension
1912  // add time dimension to start[0] and count[0]
1913  str[0] = vars[vid]->start[0]; // start along the time dimension
1914  cnt[0] = vars[vid]->count[0];
1915  for (i=0; i<ndims; i++) {
1916  str[ndims_t-i-1] = start[i] - 1;
1917  cnt[ndims_t-i-1] = count[i];
1918  }
1919  } else {
1920  for (i=0; i<ndims; i++) {
1921  str[ndims-i-1] = start[i] - 1;
1922  cnt[ndims-i-1] = count[i];
1923  }
1924  }
1925  }
1926 
1927  switch (precision) {
1928  case 8:
1929  if ( files[fid]->shared_mode )
1930  CHECK_PNC_ERROR( ncmpi_bput_vara_double(ncid, varid, str, cnt, (double*)var, NULL) )
1931  else
1932  CHECK_ERROR( nc_put_vara_double(ncid, varid, vars[vid]->start, vars[vid]->count, (double*)var) )
1933  break;
1934  case 4:
1935  if ( files[fid]->shared_mode )
1936  CHECK_PNC_ERROR( ncmpi_bput_vara_float(ncid, varid, str, cnt, (float*)var, NULL) )
1937  else
1938  CHECK_ERROR( nc_put_vara_float(ncid, varid, vars[vid]->start, vars[vid]->count, (float*)var) )
1939  break;
1940  default:
1941  fprintf(stderr, "unsupported data precision: %d\n", precision);
1942  return ERROR_CODE;
1943  }
1944 
1945  if ( files[fid]->shared_mode) {
1946  free(str);
1947  free(cnt);
1948  }
1949 
1950  return SUCCESS_CODE;
1951 }
1952 
1953 int32_t file_close_c( const int32_t fid ) // (in)
1954 {
1955  int ncid;
1956  int i;
1957 
1958  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1959  ncid = files[fid]->ncid;
1960 
1961  for (i=0; i<nvar; i++) {
1962  if ( vars[i] != NULL && vars[i]->ncid == ncid ) {
1963  free( vars[i]->start );
1964  free( vars[i]->count );
1965  free( vars[i] );
1966  vars[i] = NULL;
1967  }
1968  }
1969 
1970  for (i=0; i<nt; i++) {
1971  if ( tdims[i] != NULL && tdims[i]->ncid == ncid ) {
1972  free( tdims[i]->tval );
1973  free( tdims[i] );
1974  tdims[i] = NULL;
1975  }
1976  }
1977 
1978  if ( files[fid]->shared_mode )
1979  CHECK_PNC_ERROR( ncmpi_close(ncid) )
1980  else
1981  CHECK_ERROR( nc_close(ncid) )
1982 
1983  free( files[fid] );
1984  files[fid] = NULL;
1985 
1986  return SUCCESS_CODE;
1987 }
char att_name[File_HSHORT *ATT_MAX]
Definition: scale_file.h:34
char standard_name[File_HMID]
Definition: scale_file.h:23
#define RMISS
int32_t file_add_associatedvariable_c(const int32_t fid, const char *vname)
int32_t step
Definition: scale_file.h:28
int32_t file_get_dim_length_c(const int32_t fid, const char *dimname, int32_t *len)
int32_t file_get_step_size_c(const int32_t fid, const char *varname, int32_t *len)
int32_t file_get_varname_c(const int32_t fid, const int32_t vid, char *name, const int32_t len)
int32_t file_get_attribute_double_c(const int32_t fid, const char *vname, const char *key, const int32_t suppress, double *value, const size_t len)
int32_t file_get_attribute_float_c(const int32_t fid, const char *vname, const char *key, const int32_t suppress, float *value, const size_t len)
int32_t file_open_c(int32_t *fid, const char *fname, const int32_t mode, const MPI_Comm comm)
int32_t file_get_attribute_int_c(const int32_t fid, const char *vname, const char *key, const int32_t suppress, int *value, const size_t len)
#define NTMAX
#define SUCCESS_CODE
#define TYPE2NCTYPE(type, nctype)
#define File_FAPPEND
char calendar[File_HSHORT]
Definition: scale_file.h:32
char dim_name[File_HSHORT *RANK_MAX]
Definition: scale_file.h:26
int32_t file_set_attribute_float_c(const int32_t fid, const char *vname, const char *key, const float *value, const size_t len)
int32_t file_detach_buffer_c(const int32_t fid)
#define File_FREAD
int32_t file_put_associatedcoordinate_c(const int32_t fid, const char *name, const char *desc, const char *units, const char **dim_names, const int32_t ndims, const int32_t dtype, const void *val, const int32_t precision)
char description[File_HMID]
Definition: scale_file.h:21
#define File_FWRITE
real64_t time_start
Definition: scale_file.h:29
int32_t file_set_tunits_c(const int32_t fid, const char *time_units, const char *calendar)
#define RANK_MAX
#define File_HLONG
real(rp), allocatable, target, public k
int32_t natts
Definition: scale_file.h:33
real64_t tint
#define MIN(a, b)
#define NCTYPE2TYPE(nctype, type)
int32_t file_set_attribute_int_c(const int32_t fid, const char *vname, const char *key, const int32_t *value, const size_t len)
#define ncmpi_inq_varid(a, b, c)
int32_t file_set_option_c(const int32_t fid, const char *filetype, const char *key, const char *val)
#define ncmpi_inq_dimid(a, b, c)
logical, public i
Definition: scale_file.F90:171
#define VAR_MAX
int32_t dim_size[RANK_MAX]
Definition: scale_file.h:27
int32_t att_len[ATT_MAX]
Definition: scale_file.h:36
#define CHECK_ERROR(func)
#define ALREADY_EXISTED_CODE
#define FILE_MAX
#define File_HSHORT
int32_t fid
Definition: scale_file.h:37
int32_t file_write_axis_c(const int32_t fid, const char *name, const void *val, const int32_t precision, const MPI_Offset *start, const MPI_Offset *count)
int32_t file_enddef_c(const int32_t fid)
int32_t att_type[ATT_MAX]
Definition: scale_file.h:35
int32_t file_set_attribute_double_c(const int32_t fid, const char *vname, const char *key, const double *value, const size_t len)
#define DEFAULT_DEFLATE_LEVEL
#define ALREADY_CLOSED_CODE
int32_t file_get_attribute_text_c(const int32_t fid, const char *vname, const char *key, const int32_t suppress, char *value, const int32_t len)
#define File_HMID
char time_units[File_HMID]
Definition: scale_file.h:31
real64_t * tval
real64_t time_end
Definition: scale_file.h:30
#define ERROR_CODE
int32_t file_close_c(const int32_t fid)
size_t * start
int32_t file_put_axis_c(const int32_t fid, const char *name, const char *desc, const char *units, const char *dim_name, const int32_t dtype, const void *val, const int32_t size, const int32_t precision)
int32_t file_read_data_c(void *var, const datainfo_t *dinfo, const int32_t precision, const MPI_Offset ntypes, const MPI_Datatype dtype, const int32_t *start, const int32_t *count)
char units[File_HSHORT]
Definition: scale_file.h:22
int32_t file_get_datainfo_c(datainfo_t *dinfo, const int32_t fid, const char *varname, const int32_t step, const int32_t suppress)
int32_t file_set_attribute_text_c(const int32_t fid, const char *vname, const char *key, const char *val)
#define TEPS
int32_t file_flush_c(const int32_t fid)
real64_t t
#define CHECK_PNC_ERROR(func)
char varname[File_HSHORT]
Definition: scale_file.h:20
int32_t file_def_associatedcoordinate_c(const int32_t fid, const char *name, const char *desc, const char *units, const char **dim_names, const int32_t ndims, const int32_t dtype)
int32_t file_add_variable_c(const int32_t fid, const char *varname, const char *desc, const char *units, const char *stdname, const char **dims, const int32_t ndims, const int32_t dtype, const real64_t tint, const int32_t tavg, int32_t *vid)
real(rp), allocatable, target, public j
#define ncmpi_inq_attid(a, b, c, d)
int32_t file_def_axis_c(const int32_t fid, const char *name, const char *desc, const char *units, const char *dim_name, const int32_t dtype, const int32_t dim_size, const int32_t bounds)
size_t * count
int32_t file_write_associatedcoordinate_c(const int32_t fid, const char *name, const void *val, const int32_t precision, const MPI_Offset *start, const MPI_Offset *count)
int32_t file_get_nvars_c(const int32_t fid, int32_t *nvars)
int32_t datatype
Definition: scale_file.h:24
int32_t file_write_data_c(const int32_t fid, const int32_t vid, const void *var, const real64_t t_start, const real64_t t_end, const int32_t precision, const int32_t ndims, const int32_t *start, const int32_t *count)
int32_t file_attach_buffer_c(const int32_t fid, const int64_t buf_amount)
int32_t rank
Definition: scale_file.h:25