SCALE-RM
Classes | Macros | Functions
scale_file_netcdf.c File Reference
#include "scale_file.h"
#include "netcdf.h"
Include dependency graph for scale_file_netcdf.c:

Go to the source code of this file.

Classes

struct  fileinfo_t
 
struct  tdim_t
 
struct  varinfo_t
 

Macros

#define MPI_INCLUDED
 
#define TEPS   1e-6
 
#define NTMAX   102400
 
#define MIN(a, b)   ((a)<(b) ? (a) : (b))
 
#define CHECK_ERROR(func)
 
#define CHECK_PNC_ERROR(func)
 
#define ncmpi_inq_attid(a, b, c, d)   NC2_ERR
 
#define ncmpi_inq_varid(a, b, c)   NC2_ERR
 
#define ncmpi_inq_dimid(a, b, c)   NC2_ERR
 
#define NCTYPE2TYPE(nctype, type)
 
#define TYPE2NCTYPE(type, nctype)
 
#define DEFAULT_DEFLATE_LEVEL   2
 

Functions

int file_open_c (int *fid, const char *fname, const int mode, const int icomm)
 
int file_get_dim_length_c (int *len, const int fid, const char *dimname, const bool suppress)
 
int file_set_option_c (const int fid, const char *filetype, const char *key, const char *val)
 
int file_get_nvars_c (int *nvars, const int fid)
 
int file_get_varname_c (char *name, const int fid, const int vid, const int len)
 
int file_get_datainfo_c (datainfo_t *dinfo, const int fid, const char *varname, const int step, const bool suppress)
 
int file_get_step_size_c (int *len, const int fid, const char *varname)
 
int file_read_data_c (void *var, const datainfo_t *dinfo, const int precision, const int ntypes, const int dtype, const int *start, const int *count)
 
int file_get_attribute_text_c (char *value, const int fid, const char *vname, const char *key, const bool suppress, const int len)
 
int file_get_attribute_int_c (int *value, const int fid, const char *vname, const char *key, const bool suppress, const int len)
 
int file_get_attribute_float_c (float *value, const int fid, const char *vname, const char *key, const bool suppress, const int len)
 
int file_get_attribute_double_c (double *value, const int fid, const char *vname, const char *key, const bool suppress, const int len)
 
int file_set_attribute_text_c (const int fid, const char *vname, const char *key, const char *val)
 
int file_set_attribute_int_c (const int fid, const char *vname, const char *key, const int *value, const int len)
 
int file_set_attribute_float_c (const int fid, const char *vname, const char *key, const float *value, const int len)
 
int file_set_attribute_double_c (const int fid, const char *vname, const char *key, const double *value, const int len)
 
int file_add_associatedvariable_c (const int fid, const char *vname)
 
int file_set_tunits_c (const int fid, const char *time_units, const char *calendar)
 
int file_put_axis_c (const int fid, const char *name, const char *desc, const char *units, const char *dim_name, const int dtype, const void *val, const int size, const int precision)
 
int file_def_axis_c (const int fid, const char *name, const char *desc, const char *units, const char *dim_name, const int dtype, const int dim_size, const int bounds)
 
int file_write_axis_c (const int fid, const char *name, const void *val, const int precision, const int *start, const int *count)
 
int file_put_associatedcoordinate_c (const int fid, const char *name, const char *desc, const char *units, const char **dim_names, const int ndims, const int dtype, const void *val, const int precision)
 
int file_def_associatedcoordinate_c (const int fid, const char *name, const char *desc, const char *units, const char **dim_names, const int ndims, const int dtype)
 
int file_write_associatedcoordinate_c (const int fid, const char *name, const void *val, const int ndims, const int precision, const int *start, const int *count)
 
int file_add_variable_c (int *vid, const int fid, const char *varname, const char *desc, const char *units, const char *stdname, const char **dims, const int ndims, const int dtype, const double tint, const char *tstats)
 
int file_enddef_c (const int fid)
 
int file_redef_c (const int fid)
 
int file_attach_buffer_c (const int fid, const int64_t buf_amount)
 
int file_detach_buffer_c (const int fid)
 
int file_flush_c (const int fid)
 
int file_write_data_c (const int fid, const int vid, const void *var, const double t_start, const double t_end, const int ndims, const int precision, const int *start, const int *count)
 
int file_close_c (const int fid, const bool abort)
 

Macro Definition Documentation

◆ MPI_INCLUDED

#define MPI_INCLUDED

Definition at line 3 of file scale_file_netcdf.c.

◆ TEPS

#define TEPS   1e-6

Definition at line 7 of file scale_file_netcdf.c.

◆ NTMAX

#define NTMAX   102400

Definition at line 8 of file scale_file_netcdf.c.

◆ MIN

#define MIN (   a,
 
)    ((a)<(b) ? (a) : (b))

Definition at line 10 of file scale_file_netcdf.c.

◆ CHECK_ERROR

#define CHECK_ERROR (   func)
Value:
{ \
int status_ = (func); \
if (status_ != NC_NOERR) { \
if ( ! ERROR_SUPPRESS ) { \
fprintf(stderr, "Error: at line %d in %s\n", __LINE__, __FILE__); \
fprintf(stderr, " %s\n", nc_strerror(status_)); \
} \
return ERROR_CODE; \
} \
}

Definition at line 14 of file scale_file_netcdf.c.

◆ CHECK_PNC_ERROR

#define CHECK_PNC_ERROR (   func)
Value:
{ \
fprintf(stderr, "pnetCDF is necessary for shared_mode\n"); \
fprintf(stderr, "Please re-compile with pnetCDF\n"); \
return ERROR_CODE; \
}

Definition at line 40 of file scale_file_netcdf.c.

◆ ncmpi_inq_attid

#define ncmpi_inq_attid (   a,
  b,
  c,
 
)    NC2_ERR

Definition at line 46 of file scale_file_netcdf.c.

◆ ncmpi_inq_varid

#define ncmpi_inq_varid (   a,
  b,
 
)    NC2_ERR

Definition at line 47 of file scale_file_netcdf.c.

◆ ncmpi_inq_dimid

#define ncmpi_inq_dimid (   a,
  b,
 
)    NC2_ERR

Definition at line 48 of file scale_file_netcdf.c.

◆ NCTYPE2TYPE

#define NCTYPE2TYPE (   nctype,
  type 
)
Value:
{ \
switch ( nctype ) { \
case NC_FLOAT: \
type = File_REAL4; \
break; \
case NC_DOUBLE: \
type = File_REAL8; \
break; \
case NC_SHORT: \
type = File_INTEGER2; \
break; \
case NC_INT: \
type = File_INTEGER4; \
break; \
case NC_CHAR: \
type = File_TEXT; \
break; \
default: \
fprintf(stderr, "unsupported data type: %d\n", xtype); \
return ERROR_CODE; \
} \
}

Definition at line 51 of file scale_file_netcdf.c.

◆ TYPE2NCTYPE

#define TYPE2NCTYPE (   type,
  nctype 
)
Value:
{ \
switch ( type ) { \
case File_REAL4: \
nctype = NC_FLOAT; \
break; \
case File_REAL8: \
nctype = NC_DOUBLE; \
break; \
default: \
fprintf(stderr, "unsupported data type: %d\n", xtype); \
return ERROR_CODE; \
} \
}

Definition at line 75 of file scale_file_netcdf.c.

◆ DEFAULT_DEFLATE_LEVEL

#define DEFAULT_DEFLATE_LEVEL   2

Definition at line 91 of file scale_file_netcdf.c.

Function Documentation

◆ file_open_c()

int file_open_c ( int *  fid,
const char *  fname,
const int  mode,
const int  icomm 
)

Definition at line 170 of file scale_file_netcdf.c.

174 {
175  int ncid;
176  int len;
177  int shared_mode;
178  char _fname[File_HLONG+4];
179  int add_suffix;
180  MPI_Comm comm;
181 
182  if ( nfile >= FILE_MAX ) {
183  fprintf(stderr, "exceed max number of file limit\n");
184  return ERROR_CODE;
185  }
186 
187  len = strlen(fname);
188  strcpy(_fname, fname);
189 
190  if ( mode==File_FREAD || mode==File_FAPPEND ) {
191  FILE *fp = fopen(_fname, "r");
192  if ( fp==NULL ) {
193  add_suffix = 1;
194  } else {
195  fclose(fp);
196  add_suffix = 0;
197  }
198  } else
199  add_suffix = 1;
200 
201  if ( add_suffix )
202  if (fname[len-3] != '.' || fname[len-2] != 'n' || fname[len-1] != 'c' )
203  strcat(_fname, ".nc");
204 
205  comm = MPI_Comm_f2c(icomm);
206  if ( comm == MPI_COMM_NULL || comm == MPI_COMM_SELF )
207  shared_mode = 0;
208  else
209  shared_mode = 1;
210 
211  switch ( mode ) {
212  case File_FREAD:
213  if ( shared_mode )
214  CHECK_PNC_ERROR( ncmpi_open(comm, _fname, NC_NOWRITE, MPI_INFO_NULL, &ncid) )
215  else
216  CHECK_ERROR( nc_open(_fname, NC_NOWRITE, &ncid) )
217  break;
218  case File_FWRITE:
219  if ( shared_mode )
220  CHECK_PNC_ERROR( ncmpi_create(comm, _fname, NC_CLOBBER|NC_64BIT_OFFSET, MPI_INFO_NULL, &ncid) )
221  else
222 #ifdef NETCDF3
223  CHECK_ERROR( nc_create(_fname, NC_CLOBBER|NC_64BIT_OFFSET, &ncid) )
224 #else
225  CHECK_ERROR( nc_create(_fname, NC_CLOBBER|NC_NETCDF4, &ncid) )
226 #endif
227  break;
228  case File_FAPPEND:
229  if ( shared_mode )
230  CHECK_PNC_ERROR( ncmpi_open(comm, _fname, NC_WRITE, MPI_INFO_NULL, &ncid) )
231  else
232  CHECK_ERROR( nc_open(_fname, NC_WRITE, &ncid) )
233  break;
234  default:
235  fprintf(stderr, "invalid mode type\n");
236  return ERROR_CODE;
237  }
238 
239  files[nfile] = (fileinfo_t*) malloc(sizeof(fileinfo_t));
240  files[nfile]->ncid = ncid;
241  files[nfile]->deflate_level = DEFAULT_DEFLATE_LEVEL;
242 #if defined(NETCDF3) || defined(PNETCDF)
243  if ( mode == File_FWRITE )
244  files[nfile]->defmode = 1;
245  else
246  files[nfile]->defmode = 0;
247 #endif
248 
249  files[nfile]->shared_mode = shared_mode; /* shared-file I/O mode */
250  strcpy(files[nfile]->fname, fname);
251  *fid = nfile;
252 
253  nfile++;
254 
255  return SUCCESS_CODE;
256 }

References File_HLONG.

Referenced by scale_file::file_get_aggregate().

Here is the caller graph for this function:

◆ file_get_dim_length_c()

int file_get_dim_length_c ( int *  len,
const int  fid,
const char *  dimname,
const bool  suppress 
)

Definition at line 258 of file scale_file_netcdf.c.

262 {
263  int ncid, dimid;
264 
265  ERROR_SUPPRESS = suppress;
266 
267  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
268  ncid = files[fid]->ncid;
269 
270  if ( files[fid]->shared_mode ) {
271  MPI_Offset l;
272  CHECK_PNC_ERROR( ncmpi_inq_dimid(ncid, dimname, &dimid) )
273  CHECK_PNC_ERROR( ncmpi_inq_dimlen(ncid, dimid, &l) )
274  *len = l;
275  } else {
276  size_t l;
277  CHECK_ERROR( nc_inq_dimid(ncid, dimname, &dimid) )
278  CHECK_ERROR( nc_inq_dimlen(ncid, dimid, &l) )
279  *len = l;
280  }
281 
282  ERROR_SUPPRESS = 0;
283 
284  return SUCCESS_CODE;
285 }

Referenced by scale_file::file_get_dimlength().

Here is the caller graph for this function:

◆ file_set_option_c()

int file_set_option_c ( const int  fid,
const char *  filetype,
const char *  key,
const char *  val 
)

Definition at line 287 of file scale_file_netcdf.c.

291 {
292  if ( strcmp(filetype, "netcdf") != 0 ) return SUCCESS_CODE;
293 
294  if ( strcmp(key, "deflate_level") == 0 ) {
295  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
296  files[fid]->deflate_level = atoi(val);
297  return SUCCESS_CODE;
298  } else {
299  return ERROR_CODE;
300  }
301 }

References SUCCESS_CODE.

Referenced by scale_file::file_set_option().

Here is the caller graph for this function:

◆ file_get_nvars_c()

int file_get_nvars_c ( int *  nvars,
const int  fid 
)

Definition at line 303 of file scale_file_netcdf.c.

305 {
306  int ncid;
307  int ndims, ngatts, unlimdim;
308 
309  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
310  ncid = files[fid]->ncid;
311 
312  if ( files[fid]->shared_mode )
313  CHECK_PNC_ERROR( ncmpi_inq(ncid, &ndims, nvars, &ngatts, &unlimdim) )
314  else
315  CHECK_ERROR( nc_inq(ncid, &ndims, nvars, &ngatts, &unlimdim) )
316 
317  return SUCCESS_CODE;
318 }

Referenced by scale_file::file_create().

Here is the caller graph for this function:

◆ file_get_varname_c()

int file_get_varname_c ( char *  name,
const int  fid,
const int  vid,
const int  len 
)

Definition at line 320 of file scale_file_netcdf.c.

324 {
325  int ncid, varid;
326  char buf[MAX_NC_NAME+1];
327  int i;
328 
329  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
330  ncid = files[fid]->ncid;
331  varid = vid-1; // index starts from 1 in fortran space
332 
333  if ( files[fid]->shared_mode )
334  CHECK_PNC_ERROR( ncmpi_inq_varname(ncid, varid, buf) )
335  else
336  CHECK_ERROR( nc_inq_varname(ncid, varid, buf) )
337 
338  for (i=0; i<MIN(len-1,strlen(buf)); i++)
339  name[i] = buf[i];
340  name[i] = '\0';
341 
342  return SUCCESS_CODE;
343 }

References scale_file::i.

Referenced by scale_file::file_create().

Here is the caller graph for this function:

◆ file_get_datainfo_c()

int file_get_datainfo_c ( datainfo_t dinfo,
const int  fid,
const char *  varname,
const int  step,
const bool  suppress 
)

Definition at line 345 of file scale_file_netcdf.c.

350 {
351  int ncid, varid;
352  nc_type xtype;
353  int rank;
354  int dimids[RANK_MAX], tdim, uldims[NC_MAX_DIMS];
355  char name[NC_MAX_NAME+1];
356  char *buf;
357  size_t size;
358  int natts;
359  int i, n;
360 
361  ERROR_SUPPRESS = suppress;
362 
363  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
364  ncid = files[fid]->ncid;
365  if ( files[fid]->shared_mode )
366  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, varname, &varid) )
367  else
368  CHECK_ERROR( nc_inq_varid(ncid, varname, &varid) )
369 
370  // fid
371  dinfo->fid = fid;
372  // varname
373  strcpy(dinfo->varname, varname);
374 
375  if ( files[fid]->shared_mode ) {
376  MPI_Offset l;
377  // datatype
378  CHECK_PNC_ERROR( ncmpi_inq_vartype(ncid, varid, &xtype) )
379  NCTYPE2TYPE(xtype, dinfo->datatype);
380  // rank
381  CHECK_PNC_ERROR( ncmpi_inq_varndims(ncid, varid, &rank) )
382  if ( rank > 0 ) {
383 #ifdef PNETCDF
384  // description
385  if ( ncmpi_inq_attlen(ncid, varid, "long_name", &l) == NC_NOERR ) {
386  buf = (char*) malloc(l+1);
387  CHECK_PNC_ERROR( ncmpi_get_att_text(ncid, varid, "long_name", buf) )
388  for (i=0; i<MIN(File_HMID-1,l); i++)
389  dinfo->description[i] = buf[i];
390  dinfo->description[i] = '\0';
391  free(buf);
392  } else
393  dinfo->description[0] = '\0';
394  // units
395  if ( ncmpi_inq_attlen(ncid, varid, "units", &l) == NC_NOERR ) {
396  buf = (char*) malloc(l+1);
397  CHECK_PNC_ERROR( ncmpi_get_att_text(ncid, varid, "units", buf) )
398  for (i=0; i<MIN(File_HSHORT-1,l); i++)
399  dinfo->units[i] = buf[i];
400  dinfo->units[i] = '\0';
401  free(buf);
402  } else
403  dinfo->units[0] = '\0';
404  // standard_name
405  if ( ncmpi_inq_attlen(ncid, varid, "standard_name", &l) == NC_NOERR ) {
406  buf = (char*) malloc(l+1);
407  CHECK_PNC_ERROR( ncmpi_get_att_text(ncid, varid, "standard_name", buf) )
408  for (i=0; i<MIN(File_HMID-1,l); i++)
409  dinfo->standard_name[i] = buf[i];
410  dinfo->standard_name[i] = '\0';
411  free(buf);
412  } else
413  dinfo->standard_name[0] = '\0';
414 #endif
415  // dimensions
416  CHECK_PNC_ERROR( ncmpi_inq_vardimid(ncid, varid, dimids) )
417 #if 1
418  CHECK_PNC_ERROR( ncmpi_inq_unlimdim(ncid, uldims) )
419  n = 1;
420 #else
421  CHECK_PNC_ERROR( ncmpi_inq_unlimdims(ncid, &n, uldims) )
422 #endif
423  }
424  // attributes
425  dinfo->natts = 0;
426  CHECK_PNC_ERROR( ncmpi_inq_varnatts(ncid, varid, &natts) )
427  for (i=0; i<natts; i++) {
428  CHECK_PNC_ERROR( ncmpi_inq_attname(ncid, varid, i, name) )
429  if ( strcmp(name, "long_name") &&
430  strcmp(name, "description") &&
431  strcmp(name, "units") &&
432  strcmp(name, "standard_name") ) {
433  strcpy(dinfo->att_name+File_HSHORT*dinfo->natts, name);
434  CHECK_PNC_ERROR( ncmpi_inq_atttype(ncid, varid, name, &xtype) )
435  NCTYPE2TYPE(xtype, dinfo->att_type[dinfo->natts]);
436  CHECK_PNC_ERROR( ncmpi_inq_attlen(ncid, varid, name, &l) )
437  dinfo->att_len[dinfo->natts] = l;
438  dinfo->natts++;
439  }
440  }
441  } else {
442  size_t l;
443  // datatype
444  CHECK_ERROR( nc_inq_vartype(ncid, varid, &xtype) )
445  NCTYPE2TYPE(xtype, dinfo->datatype);
446  // rank
447  CHECK_ERROR( nc_inq_varndims(ncid, varid, &rank) )
448  if ( rank > 0 ) {
449  // description
450  if ( nc_inq_attlen(ncid, varid, "long_name", &l) == NC_NOERR ) {
451  buf = (char*) malloc(l+1);
452  CHECK_ERROR( nc_get_att_text(ncid, varid, "long_name", buf) )
453  } else if ( nc_inq_attlen(ncid, varid, "description", &l) == NC_NOERR ) { // for WRF file
454  buf = (char*) malloc(l+1);
455  CHECK_ERROR( nc_get_att_text(ncid, varid, "description", buf) )
456  } else
457  l = 0;
458  for (i=0; i<MIN(File_HMID-1,l); i++)
459  dinfo->description[i] = buf[i];
460  dinfo->description[i] = '\0';
461  if ( l>0 ) free(buf);
462  // units
463  if ( nc_inq_attlen(ncid, varid, "units", &l) == NC_NOERR ) {
464  buf = (char*) malloc(l+1);
465  CHECK_ERROR( nc_get_att_text(ncid, varid, "units", buf) )
466  for (i=0; i<MIN(File_HSHORT-1,l); i++)
467  dinfo->units[i] = buf[i];
468  dinfo->units[i] = '\0';
469  free(buf);
470  } else
471  dinfo->units[0] = '\0';
472  // standard_name
473  if ( nc_inq_attlen(ncid, varid, "standard_name", &l) == NC_NOERR ) {
474  buf = (char*) malloc(l+1);
475  CHECK_ERROR( nc_get_att_text(ncid, varid, "standard_name", buf) )
476  for (i=0; i<MIN(File_HMID-1,l); i++)
477  dinfo->standard_name[i] = buf[i];
478  dinfo->standard_name[i] = '\0';
479  free(buf);
480  } else
481  dinfo->standard_name[0] = '\0';
482  // dimension
483  CHECK_ERROR( nc_inq_vardimid(ncid, varid, dimids) )
484 #ifdef NETCDF3
485  CHECK_ERROR( nc_inq_unlimdim(ncid, uldims) )
486  n = 1;
487 #else
488  CHECK_ERROR( nc_inq_unlimdims(ncid, &n, uldims) )
489 #endif
490  }
491  // attributes
492  dinfo->natts = 0;
493  CHECK_ERROR( nc_inq_varnatts(ncid, varid, &natts) )
494  for (i=0; i<natts; i++) {
495  CHECK_ERROR( nc_inq_attname(ncid, varid, i, name) )
496  if ( strcmp(name, "long_name") &&
497  strcmp(name, "description") &&
498  strcmp(name, "units") &&
499  strcmp(name, "standard_name") ) {
500  strcpy(dinfo->att_name+File_HSHORT*dinfo->natts, name);
501  CHECK_ERROR( nc_inq_atttype(ncid, varid, name, &xtype) )
502  NCTYPE2TYPE(xtype, dinfo->att_type[dinfo->natts]);
503  CHECK_ERROR( nc_inq_attlen(ncid, varid, name, &l) )
504  dinfo->att_len[dinfo->natts] = l;
505  dinfo->natts++;
506  }
507  }
508  }
509 
510  if ( rank > 0 ) {
511  tdim = -1;
512  for ( i=0; i<n; i++ ) {
513  if ( uldims[i] == dimids[0] ) {
514  tdim = uldims[i];
515  break;
516  }
517  }
518  if (rank > RANK_MAX) {
519  fprintf(stderr, "rank exceeds limit: %d\n", rank);
520  return ERROR_CODE;
521  }
522  dinfo->rank = tdim >= 0 ? rank -1 : rank; // do not count time dimension
523  // dim_name and dim_size
524  for (i=0; i<dinfo->rank; i++) {
525  // note: C and Fortran orders are opposite
526  if ( files[fid]->shared_mode ) {
527  MPI_Offset size_;
528  CHECK_PNC_ERROR( ncmpi_inq_dim(ncid, dimids[rank-i-1], name, &size_) )
529  size = (size_t)size_;
530  }
531  else
532  CHECK_ERROR( nc_inq_dim(ncid, dimids[rank-i-1], name, &size) )
533  if ( strlen(name) > File_HSHORT-1 ) {
534  fprintf(stderr, "Length of the dimension name (%s) is too long (should be < %d).\n", name, File_HSHORT);
535  return ERROR_CODE;
536  }
537  strncpy(dinfo->dim_name+i*File_HSHORT, name, File_HSHORT);
538  dinfo->dim_size[i] = size;
539  }
540 
541  dinfo->step = step;
542  if ( tdim >= 0 ) {
543  if ( files[fid]->shared_mode ) {
544  MPI_Offset idx[2];
545  MPI_Offset l;
546  // time_end
547  CHECK_PNC_ERROR( ncmpi_inq_dimname(ncid, tdim, name) )
548  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, name, &varid) )
549  idx[0] = step - 1;
550  CHECK_PNC_ERROR( ncmpi_get_var1_double_all(ncid, varid, idx, &(dinfo->time_end)) )
551  // units
552  CHECK_PNC_ERROR( ncmpi_inq_attlen (ncid, varid, "units", &l) )
553  buf = (char*) malloc(l+1);
554  CHECK_PNC_ERROR( ncmpi_get_att_text(ncid, varid, "units", buf) )
555  for (i=0; i<MIN(File_HMID-1,l); i++)
556  dinfo->time_units[i] = buf[i];
557  dinfo->time_units[i] = '\0';
558  free(buf);
559  // calendar
560 #ifdef PNETCDF
561  if ( ncmpi_inq_attlen(ncid, varid, "calendar", &l) == NC_NOERR ) {
562  buf = (char*) malloc(l+1);
563  CHECK_PNC_ERROR( ncmpi_get_att_text(ncid, varid, "calendar", buf) )
564  for (i=0; i<MIN(File_HSHORT-1,l); i++)
565  dinfo->calendar[i] = buf[i];
566  dinfo->calendar[i] = '\0';
567  free(buf);
568  } else
569  dinfo->calendar[0] = '\0';
570 #endif
571  // time_start
572  strcat(name, "_bnds");
573  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, name, &varid) )
574  idx[1] = 0;
575  CHECK_PNC_ERROR( ncmpi_get_var1_double_all(ncid, varid, idx, &(dinfo->time_start)) )
576  } else {
577  size_t idx[2];
578  size_t l;
579  // time_end
580  CHECK_ERROR( nc_inq_dimname(ncid, tdim, name) )
581  if ( nc_inq_varid(ncid, name, &varid) == NC_NOERR ) {
582  idx[0] = step - 1;
583  CHECK_ERROR( nc_get_var1_double(ncid, varid, idx, &(dinfo->time_end)) )
584  // units
585  CHECK_ERROR( nc_inq_attlen (ncid, varid, "units", &l) )
586  buf = (char*) malloc(l+1);
587  CHECK_ERROR( nc_get_att_text(ncid, varid, "units", buf) )
588  for (i=0; i<MIN(File_HMID-1,l); i++)
589  dinfo->time_units[i] = buf[i];
590  dinfo->time_units[i] = '\0';
591  free(buf);
592  // calendar
593  if ( nc_inq_attlen(ncid, varid, "calendar", &l) == NC_NOERR ) {
594  buf = (char*) malloc(l+1);
595  CHECK_ERROR( nc_get_att_text(ncid, varid, "calendar", buf) )
596  for (i=0; i<MIN(File_HSHORT-1,l); i++)
597  dinfo->calendar[i] = buf[i];
598  dinfo->calendar[i] = '\0';
599  free(buf);
600  } else
601  dinfo->calendar[0] = '\0';
602  // time_start
603  strcat(name, "_bnds");
604  CHECK_ERROR( nc_inq_varid(ncid, name, &varid) )
605  idx[1] = 0;
606  CHECK_ERROR( nc_get_var1_double(ncid, varid, idx, &(dinfo->time_start)) )
607  } else {
608  dinfo->time_start = 0.0;
609  dinfo->time_end = 0.0;
610  dinfo->time_units[0] = '\0';
611  dinfo->calendar[0] = '\0';
612  }
613  }
614  dinfo->has_tdim = 1;
615  } else {
616  if ( step > 1 ) { // if variable does not have time dimention, step > 1 should not exist
617  fprintf(stderr, "requested step is larger than tdim: step=%d tdim=%d\n", step, tdim);
618  return ERROR_CODE;
619  }
620  dinfo->has_tdim = 0;
621  dinfo->time_start = 0.0;
622  dinfo->time_end = 0.0;
623  dinfo->time_units[0] = '\0';
624  dinfo->calendar[0] = '\0';
625  }
626  } else {
627  dinfo->rank = 0;
628  dinfo->description[0] = '\0';
629  dinfo->units[0] = '\0';
630  dinfo->standard_name[0] = '\0';
631  dinfo->has_tdim = 0;
632  dinfo->time_start = 0.0;
633  dinfo->time_end = 0.0;
634  dinfo->time_units[0] = '\0';
635  dinfo->calendar[0] = '\0';
636  }
637 
638  ERROR_SUPPRESS = 0;
639 
640  return SUCCESS_CODE;
641 }

Referenced by scale_file::file_def_variable(), and scale_file::file_get_stepsize().

Here is the caller graph for this function:

◆ file_get_step_size_c()

int file_get_step_size_c ( int *  len,
const int  fid,
const char *  varname 
)

Definition at line 643 of file scale_file_netcdf.c.

646 {
647  int ncid, varid;
648 
649  int dimids[RANK_MAX], uldims[NC_MAX_DIMS], tdim;
650  int n, i;
651 
652  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
653  ncid = files[fid]->ncid;
654  if ( files[fid]->shared_mode )
655  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, varname, &varid) )
656  else
657  CHECK_ERROR( nc_inq_varid(ncid, varname, &varid) )
658 
659  if ( files[fid]->shared_mode ) {
660  CHECK_PNC_ERROR( ncmpi_inq_vardimid(ncid, varid, dimids) )
661  CHECK_PNC_ERROR( ncmpi_inq_unlimdim(ncid, uldims) )
662  n = uldims[0] < 0 ? 0 : 1;
663  } else {
664  CHECK_ERROR( nc_inq_vardimid(ncid, varid, dimids) )
665 #ifdef NETCDF3
666  CHECK_ERROR( nc_inq_unlimdim(ncid, uldims) )
667  n = uldims[0] < 0 ? 0 : 1;
668 #else
669  CHECK_ERROR( nc_inq_unlimdims(ncid, &n, uldims) )
670 #endif
671  }
672 
673  tdim = -1;
674  for ( i=0; i<n; i++ ) {
675  if ( uldims[i] == dimids[0] ) {
676  tdim = uldims[i];
677  break;
678  }
679  }
680 
681  if ( tdim > 0 ) {
682  if ( files[fid]->shared_mode ) {
683  MPI_Offset l;
684  CHECK_PNC_ERROR( ncmpi_inq_dimlen(ncid, tdim, &l) )
685  *len = l;
686  } else {
687  size_t l;
688  CHECK_ERROR( nc_inq_dimlen(ncid, tdim, &l) )
689  *len = l;
690  }
691  } else
692  *len = 0;
693 
694  return SUCCESS_CODE;
695 }

References scale_file::i, and RANK_MAX.

Referenced by scale_file::file_get_stepsize().

Here is the caller graph for this function:

◆ file_read_data_c()

int file_read_data_c ( void *  var,
const datainfo_t dinfo,
const int  precision,
const int  ntypes,
const int  dtype,
const int *  start,
const int *  count 
)

Definition at line 697 of file scale_file_netcdf.c.

704 {
705  int ncid, varid;
706  int rank;
707  int i;
708  int fid;
709  size_t *str, *cnt;
710  MPI_Datatype dtype_;
711  MPI_Offset *strp, *cntp;
712  size_t size;
713 
714  fid = dinfo->fid;
715  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
716  ncid = files[fid]->ncid;
717  if ( files[fid]->shared_mode ) {
718  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, dinfo->varname, &varid) )
719  CHECK_PNC_ERROR( ncmpi_inq_varndims(ncid, varid, &rank) )
720  strp = (MPI_Offset*) malloc(sizeof(MPI_Offset)*rank);
721  cntp = (MPI_Offset*) malloc(sizeof(MPI_Offset)*rank);
722  } else {
723  CHECK_ERROR( nc_inq_varid(ncid, dinfo->varname, &varid) )
724  CHECK_ERROR( nc_inq_varndims(ncid, varid, &rank) )
725  }
726  str = (size_t*) malloc(sizeof(size_t)*rank);
727  cnt = (size_t*) malloc(sizeof(size_t)*rank);
728 
729  if ( start==NULL || start[0]<=0 || count==NULL || count[0]<=0 ) {
730  for (i=0; i<dinfo->rank; i++) {
731  // note: C and Fortran orders are opposite
732  str[rank -i-1] = 0;
733  cnt[rank -i-1] = dinfo->dim_size[i];
734  }
735  } else {
736  for (i=0; i<dinfo->rank; i++) {
737  // note: C and Fortran orders are opposite
738  str[rank -i-1] = start[i] - 1;
739  cnt[rank -i-1] = count[i];
740  }
741  }
742  if (rank > dinfo->rank) { // have time dimension
743  str[0] = dinfo->step - 1;
744  cnt[0] = 1;
745  }
746 
747  size = 1;
748  for (i=0; i<rank; i++) size *= cnt[i];
749 
750  if ( files[fid]->shared_mode ) {
751 #ifdef PNETCDF
752  for (i=0; i<rank; i++) {
753  strp[i] = (MPI_Offset) str[i];
754  cntp[i] = (MPI_Offset) cnt[i];
755  }
756  free(str);
757  free(cnt);
758  dtype_ = dtype == 0 ? MPI_DATATYPE_NULL : MPI_Type_f2c(dtype);
759  CHECK_PNC_ERROR( ncmpi_iget_vara(ncid, varid, strp, cntp, var, (MPI_Offset) ntypes, dtype_, NULL) )
760  free(strp);
761  free(cntp);
762  if ( dtype_ == MPI_FLOAT ) {
763  float factor, offset, misval;
764  int l_rescale = 0;
765  if ( ncmpi_get_att_float(ncid, varid, "missing_value", &misval) == NC_NOERR )
766  if ( misval != RMISS )
767  for (i=0; i<size; i++) if ( ((float*)var)[i] == misval ) ((float*)var)[i] = RMISS;
768  if ( ncmpi_get_att_float(ncid, varid, "scale_factor", &factor) != NC_NOERR )
769  factor = 1.0f;
770  else
771  l_rescale = 1;
772  if ( ncmpi_get_att_float(ncid, varid, "add_offset", &offset) != NC_NOERR )
773  offset = 0.0f;
774  else
775  l_rescale = 1;
776  if ( l_rescale ) for (i=0; i<size; i++) ((float*)var)[i] = ((float*)var)[i] * factor + offset;
777  } else if ( dtype_ == MPI_DOUBLE ) {
778  double factor, offset, misval;
779  int l_rescale = 0;
780  if ( ncmpi_get_att_double(ncid, varid, "missing_value", &misval) == NC_NOERR )
781  if ( (float)misval != (float)RMISS )
782  for (i=0; i<size; i++) if ( ((double*)var)[i] == misval ) ((double*)var)[i] = RMISS;
783  if ( ncmpi_get_att_double(ncid, varid, "scale_factor", &factor) != NC_NOERR )
784  factor = 1.0;
785  else
786  l_rescale = 1;
787  if ( ncmpi_get_att_double(ncid, varid, "add_offset", &offset) != NC_NOERR )
788  offset = 0.0;
789  else
790  l_rescale = 1;
791  if ( l_rescale ) for (i=0; i<size; i++) ((double*)var)[i] = ((double*)var)[i] * factor + offset;
792  } else {
793  float factor, offset, misval;
794  if ( ( ncmpi_get_att_float(ncid, varid, "missing_value", &misval) == NC_NOERR ) )
795  if ( misval != (float)RMISS ) {
796  fprintf(stderr, "missing_value (!=UNDEF) is not supported with a MPI derived type\n");
797  return ERROR_CODE;
798  }
799  if ( ( ncmpi_get_att_float(ncid, varid, "scale_factor", &factor) == NC_NOERR )
800  || ( ncmpi_get_att_float(ncid, varid, "add_offset", &offset) == NC_NOERR ) ) {
801  fprintf(stderr, "scale_factor and add_offset is not supported with a MPI derived type\n");
802  return ERROR_CODE;
803  }
804  }
805 #else
806  CHECK_PNC_ERROR( dummy )
807 #endif
808  } else {
809  switch ( precision ) {
810  case 8:
811  CHECK_ERROR( nc_get_vara_double(ncid, varid, str, cnt, (double*)var) )
812  {
813  double factor, offset, misval;
814  float a;
815  int l_rescale = 0;
816  nc_get_att_float(ncid, varid, "missing_value", &a);
817  if ( nc_get_att_double(ncid, varid, "missing_value", &misval) == NC_NOERR )
818  if ( (float)misval != (float)RMISS )
819  for (i=0; i<size; i++) if ( ((double*)var)[i] == misval ) ((double*)var)[i] = RMISS;
820  if ( nc_get_att_double(ncid, varid, "scale_factor", &factor) != NC_NOERR )
821  factor = 1.0;
822  else
823  l_rescale = 1;
824  if ( nc_get_att_double(ncid, varid, "add_offset", &offset) != NC_NOERR )
825  offset = 0.0;
826  else
827  l_rescale = 1;
828  if ( l_rescale ) for (i=0; i<size; i++) ((double*)var)[i] = ((double*)var)[i] * factor + offset;
829  }
830  break;
831  case 4:
832  CHECK_ERROR( nc_get_vara_float(ncid, varid, str, cnt, (float*)var) )
833  {
834  float factor, offset, misval;
835  int l_rescale = 0;
836  if ( nc_get_att_float(ncid, varid, "missing_value", &misval) == NC_NOERR )
837  if ( misval != (float)RMISS )
838  for (i=0; i<size; i++) if ( ((float*)var)[i] == misval ) ((float*)var)[i] = RMISS;
839  if ( nc_get_att_float(ncid, varid, "scale_factor", &factor) != NC_NOERR )
840  factor = 1.0f;
841  else
842  l_rescale = 1;
843  if ( nc_get_att_float(ncid, varid, "add_offset", &offset) != NC_NOERR )
844  offset = 0.0f;
845  else
846  l_rescale = 1;
847  if ( l_rescale ) for (i=0; i<size; i++) ((float*)var)[i] = ((float*)var)[i] * factor + offset;
848  }
849  break;
850  default:
851  free(str);
852  free(cnt);
853  fprintf(stderr, "unsupported data precision: %d\n", precision );
854  return ERROR_CODE;
855  }
856  free(str);
857  free(cnt);
858  }
859 
860 
861  return SUCCESS_CODE;
862 }

Referenced by scale_file::file_get_stepsize().

Here is the caller graph for this function:

◆ file_get_attribute_text_c()

int file_get_attribute_text_c ( char *  value,
const int  fid,
const char *  vname,
const char *  key,
const bool  suppress,
const int  len 
)

Definition at line 864 of file scale_file_netcdf.c.

870 {
871  int ncid;
872  int varid;
873 
874  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
875  ncid = files[fid]->ncid;
876 
877  ERROR_SUPPRESS = suppress;
878 
879  if ( files[fid]->shared_mode ) {
880  MPI_Offset l;
881  if ( strcmp(vname, "global") == 0 ) {
882  varid = NC_GLOBAL;
883  } else
884  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
885 
886  CHECK_PNC_ERROR( ncmpi_inq_attlen(ncid, varid, key, &l) )
887  if ( len-1 < l ) return ERROR_CODE;
888 
889  CHECK_PNC_ERROR( ncmpi_get_att_text(ncid, varid, key, value) )
890  value[l] = '\0';
891  }
892  else {
893  size_t l;
894  if ( strcmp(vname, "global") == 0 ) {
895  varid = NC_GLOBAL;
896  } else
897  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
898 
899  CHECK_ERROR( nc_inq_attlen(ncid, varid, key, &l) )
900  if ( len-1 < l ) return ERROR_CODE;
901 
902  CHECK_ERROR( nc_get_att_text(ncid, varid, key, value) )
903  value[l] = '\0';
904  }
905 
906  ERROR_SUPPRESS = 0;
907 
908  return SUCCESS_CODE;
909 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_get_attribute_int_c()

int file_get_attribute_int_c ( int *  value,
const int  fid,
const char *  vname,
const char *  key,
const bool  suppress,
const int  len 
)

Definition at line 911 of file scale_file_netcdf.c.

917 {
918  int ncid;
919  int varid;
920 
921  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
922  ncid = files[fid]->ncid;
923 
924  ERROR_SUPPRESS = suppress;
925 
926  if ( files[fid]->shared_mode ) {
927  MPI_Offset l;
928  if ( strcmp(vname, "global") == 0 ) {
929  varid = NC_GLOBAL;
930  } else
931  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
932 
933  CHECK_PNC_ERROR( ncmpi_inq_attlen(ncid, varid, key, &l) )
934  if ( len != l ) return ERROR_CODE;
935  CHECK_PNC_ERROR( ncmpi_get_att_int(ncid, varid, key, value) )
936  }
937  else {
938  size_t l;
939  if ( strcmp(vname, "global") == 0 ) {
940  varid = NC_GLOBAL;
941  } else
942  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
943 
944  CHECK_ERROR( nc_inq_attlen(ncid, varid, key, &l) )
945  if ( len != l ) return ERROR_CODE;
946  CHECK_ERROR( nc_get_att_int(ncid, varid, key, value) )
947  }
948 
949  ERROR_SUPPRESS = 0;
950 
951  return SUCCESS_CODE;
952 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_get_attribute_float_c()

int file_get_attribute_float_c ( float value,
const int  fid,
const char *  vname,
const char *  key,
const bool  suppress,
const int  len 
)

Definition at line 954 of file scale_file_netcdf.c.

960 {
961  int ncid;
962  int varid;
963 
964  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
965  ncid = files[fid]->ncid;
966 
967  ERROR_SUPPRESS = suppress;
968 
969  if ( files[fid]->shared_mode ) {
970  MPI_Offset l;
971  if ( strcmp(vname, "global") == 0 ) {
972  varid = NC_GLOBAL;
973  } else
974  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
975 
976  CHECK_PNC_ERROR( ncmpi_inq_attlen(ncid, varid, key, &l) )
977  if ( len != l ) return ERROR_CODE;
978  CHECK_PNC_ERROR( ncmpi_get_att_float(ncid, varid, key, value) )
979  }
980  else {
981  size_t l;
982  if ( strcmp(vname, "global") == 0 ) {
983  varid = NC_GLOBAL;
984  } else
985  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
986 
987  CHECK_ERROR( nc_inq_attlen(ncid, varid, key, &l) )
988  if ( len != l ) return ERROR_CODE;
989  CHECK_ERROR( nc_get_att_float(ncid, varid, key, value) )
990  }
991 
992  ERROR_SUPPRESS = 0;
993 
994  return SUCCESS_CODE;
995 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_get_attribute_double_c()

int file_get_attribute_double_c ( double value,
const int  fid,
const char *  vname,
const char *  key,
const bool  suppress,
const int  len 
)

Definition at line 997 of file scale_file_netcdf.c.

1003 {
1004  int ncid;
1005  int varid;
1006 
1007  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1008  ncid = files[fid]->ncid;
1009 
1010  ERROR_SUPPRESS = suppress;
1011 
1012  if ( files[fid]->shared_mode ) {
1013  MPI_Offset l;
1014  if ( strcmp(vname, "global") == 0 ) {
1015  varid = NC_GLOBAL;
1016  } else
1017  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
1018 
1019  CHECK_PNC_ERROR( ncmpi_inq_attlen(ncid, varid, key, &l) )
1020  if ( len != l ) return ERROR_CODE;
1021  CHECK_PNC_ERROR( ncmpi_get_att_double(ncid, varid, key, value) )
1022  }
1023  else {
1024  size_t l;
1025  if ( strcmp(vname, "global") == 0 ) {
1026  varid = NC_GLOBAL;
1027  } else
1028  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
1029 
1030  CHECK_ERROR( nc_inq_attlen(ncid, varid, key, &l) )
1031  if ( len != l ) return ERROR_CODE;
1032  CHECK_ERROR( nc_get_att_double(ncid, varid, key, value) )
1033  }
1034 
1035  ERROR_SUPPRESS = 0;
1036 
1037  return SUCCESS_CODE;
1038 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_set_attribute_text_c()

int file_set_attribute_text_c ( const int  fid,
const char *  vname,
const char *  key,
const char *  val 
)

Definition at line 1040 of file scale_file_netcdf.c.

1044 {
1045  int ncid;
1046  int varid;
1047  int attid;
1048  int ret;
1049 
1050  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1051  ncid = files[fid]->ncid;
1052 
1053  if ( files[fid]->shared_mode ) {
1054  if ( strcmp(vname, "global") == 0 ) {
1055  varid = NC_GLOBAL;
1056  } else
1057  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
1058 
1059  if ( ncmpi_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1060  }
1061  else {
1062  if ( strcmp(vname, "global") == 0 ) {
1063  varid = NC_GLOBAL;
1064  } else
1065  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
1066 
1067  if ( nc_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1068  }
1069 
1070  ret = file_redef(fid, ncid);
1071  if ( ret != SUCCESS_CODE ) return ret;
1072 
1073  if ( files[fid]->shared_mode )
1074  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, key, strlen(val), val) )
1075  else
1076  CHECK_ERROR( nc_put_att_text(ncid, varid, key, strlen(val), val) )
1077 
1078  return SUCCESS_CODE;
1079 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_set_attribute_int_c()

int file_set_attribute_int_c ( const int  fid,
const char *  vname,
const char *  key,
const int *  value,
const int  len 
)

Definition at line 1081 of file scale_file_netcdf.c.

1086 {
1087  int ncid;
1088  int varid;
1089  int attid;
1090  int ret;
1091 
1092  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1093  ncid = files[fid]->ncid;
1094 
1095  if ( files[fid]->shared_mode ) {
1096  if ( strcmp(vname, "global") == 0 ) {
1097  varid = NC_GLOBAL;
1098  } else
1099  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
1100 
1101  if ( ncmpi_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1102  }
1103  else {
1104  if ( strcmp(vname, "global") == 0 ) {
1105  varid = NC_GLOBAL;
1106  } else
1107  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
1108 
1109  if ( nc_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1110  }
1111 
1112  ret = file_redef(fid, ncid);
1113  if ( ret != SUCCESS_CODE ) return ret;
1114 
1115  if ( files[fid]->shared_mode )
1116  CHECK_PNC_ERROR( ncmpi_put_att_int(ncid, varid, key, NC_INT, len, value) )
1117  else
1118  CHECK_ERROR( nc_put_att_int(ncid, varid, key, NC_INT, len, value) )
1119 
1120 
1121  return SUCCESS_CODE;
1122 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_set_attribute_float_c()

int file_set_attribute_float_c ( const int  fid,
const char *  vname,
const char *  key,
const float value,
const int  len 
)

Definition at line 1124 of file scale_file_netcdf.c.

1129 {
1130  int ncid;
1131  int varid;
1132  int attid;
1133  int ret;
1134 
1135  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1136  ncid = files[fid]->ncid;
1137 
1138  if ( files[fid]->shared_mode ) {
1139  if ( strcmp(vname, "global") == 0 ) {
1140  varid = NC_GLOBAL;
1141  } else
1142  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
1143 
1144  if ( ncmpi_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1145  }
1146  else {
1147  if ( strcmp(vname, "global") == 0 ) {
1148  varid = NC_GLOBAL;
1149  } else
1150  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
1151 
1152  if ( nc_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1153  }
1154 
1155  ret = file_redef(fid, ncid);
1156  if ( ret != SUCCESS_CODE ) return ret;
1157 
1158  if ( files[fid]->shared_mode )
1159  CHECK_PNC_ERROR( ncmpi_put_att_float(ncid, varid, key, NC_FLOAT, len, value) )
1160  else
1161  CHECK_ERROR( nc_put_att_float(ncid, varid, key, NC_FLOAT, len, value) )
1162 
1163  return SUCCESS_CODE;
1164 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_set_attribute_double_c()

int file_set_attribute_double_c ( const int  fid,
const char *  vname,
const char *  key,
const double value,
const int  len 
)

Definition at line 1166 of file scale_file_netcdf.c.

1171 {
1172  int ncid;
1173  int varid;
1174  int attid;
1175 
1176  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1177  ncid = files[fid]->ncid;
1178 
1179  if ( files[fid]->shared_mode ) {
1180  if ( strcmp(vname, "global") == 0 ) {
1181  varid = NC_GLOBAL;
1182  } else
1183  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
1184 
1185  if ( ncmpi_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1186  }
1187  else {
1188  if ( strcmp(vname, "global") == 0 ) {
1189  varid = NC_GLOBAL;
1190  } else
1191  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
1192 
1193  if ( nc_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1194  }
1195 
1196  if ( files[fid]->shared_mode )
1197  CHECK_PNC_ERROR( ncmpi_put_att_double(ncid, varid, key, NC_DOUBLE, len, value) )
1198  else
1199  CHECK_ERROR( nc_put_att_double(ncid, varid, key, NC_DOUBLE, len, value) )
1200 
1201  return SUCCESS_CODE;
1202 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_add_associatedvariable_c()

int file_add_associatedvariable_c ( const int  fid,
const char *  vname 
)

Definition at line 1204 of file scale_file_netcdf.c.

1206 {
1207  int ncid, varid;
1208  int ret;
1209 
1210  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1211  ncid = files[fid]->ncid;
1212 
1213  if ( nc_inq_varid(ncid, vname, &varid) == NC_NOERR ) // check if existed
1214  return ALREADY_EXISTED_CODE;
1215 
1216  ret = file_redef(fid, ncid);
1217  if ( ret != SUCCESS_CODE ) return ret;
1218 
1219  if ( files[fid]->shared_mode )
1220  CHECK_PNC_ERROR( ncmpi_def_var(ncid, vname, NC_INT, 0, 0, &varid) )
1221  else
1222  CHECK_ERROR( nc_def_var(ncid, vname, NC_INT, 0, 0, &varid) )
1223 
1224  return SUCCESS_CODE;
1225 }

Referenced by scale_file::file_add_associatedvariable().

Here is the caller graph for this function:

◆ file_set_tunits_c()

int file_set_tunits_c ( const int  fid,
const char *  time_units,
const char *  calendar 
)

Definition at line 1227 of file scale_file_netcdf.c.

1230 {
1231  strcpy(files[fid]->time_units, time_units);
1232  strcpy(files[fid]->calendar, calendar);
1233 
1234  return SUCCESS_CODE;
1235 }

Referenced by scale_file::file_create().

Here is the caller graph for this function:

◆ file_put_axis_c()

int file_put_axis_c ( const int  fid,
const char *  name,
const char *  desc,
const char *  units,
const char *  dim_name,
const int  dtype,
const void *  val,
const int  size,
const int  precision 
)

Definition at line 1237 of file scale_file_netcdf.c.

1246 {
1247  int ncid, dimid, varid;
1248  nc_type xtype = -1;
1249  int ret;
1250 
1251  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1252  ncid = files[fid]->ncid;
1253 
1254  if ( nc_inq_varid(ncid, name, &varid) == NC_NOERR ) // check if existed
1255  return ALREADY_EXISTED_CODE;
1256 
1257  ret = file_redef(fid, ncid);
1258  if ( ret != SUCCESS_CODE ) return ret;
1259 
1260  if ( nc_inq_dimid(ncid, dim_name, &dimid) != NC_NOERR ) // check if existed
1261  CHECK_ERROR( nc_def_dim(ncid, dim_name, size, &dimid) )
1262 
1263  TYPE2NCTYPE(dtype, xtype);
1264  CHECK_ERROR( nc_def_var(ncid, name, xtype, 1, &dimid, &varid) )
1265  if ( strlen(desc)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1266  if ( strlen(units)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "units", strlen(units), units) )
1267 
1268  ret = file_enddef(fid, ncid);
1269  if ( ret != SUCCESS_CODE ) return ret;
1270 
1271  switch ( precision ) {
1272  case 8:
1273  CHECK_ERROR( nc_put_var_double(ncid, varid, (double*)val) )
1274  break;
1275  case 4:
1276  CHECK_ERROR( nc_put_var_float(ncid, varid, (float*)val) )
1277  break;
1278  default:
1279  fprintf(stderr, "unsupported data precision: %d\n", precision);
1280  return ERROR_CODE;
1281  }
1282 
1283  return SUCCESS_CODE;
1284 }

Referenced by scale_file::file_get_dimlength().

Here is the caller graph for this function:

◆ file_def_axis_c()

int file_def_axis_c ( const int  fid,
const char *  name,
const char *  desc,
const char *  units,
const char *  dim_name,
const int  dtype,
const int  dim_size,
const int  bounds 
)

Definition at line 1286 of file scale_file_netcdf.c.

1294 {
1295  int ncid, dimid, varid;
1296  nc_type xtype = -1;
1297  int dimids[2];
1298  char buf[File_HSHORT+6];
1299  int ret;
1300 
1301  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1302  ncid = files[fid]->ncid;
1303 
1304  if ( files[fid]->shared_mode ) {
1305  if ( ncmpi_inq_varid(ncid, name, &varid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1306  } else {
1307  if ( nc_inq_varid(ncid, name, &varid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1308  }
1309 
1310  ret = file_redef(fid, ncid);
1311  if ( ret != SUCCESS_CODE ) return ret;
1312 
1313  TYPE2NCTYPE(dtype, xtype);
1314  if ( files[fid]->shared_mode ) {
1315  if ( ncmpi_inq_dimid(ncid, dim_name, &dimid) != NC_NOERR ) // check if existed
1316  CHECK_PNC_ERROR( ncmpi_def_dim(ncid, dim_name, dim_size, &dimid) )
1317 
1318  CHECK_PNC_ERROR( ncmpi_def_var(ncid, name, xtype, 1, &dimid, &varid) )
1319  if ( strlen(desc)>0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1320  if ( strlen(units)>0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "units", strlen(units), units) )
1321 
1322  if ( bounds ) {
1323  dimids[0] = dimid;
1324  if ( ncmpi_inq_dimid(ncid, "nv", &(dimids[1])) != NC_NOERR ) // first called
1325  CHECK_PNC_ERROR( ncmpi_def_dim(ncid, "nv", 2, &(dimids[1])) )
1326  sprintf(buf, "%s_bnds", dim_name);
1327  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "bounds", strlen(buf), buf) )
1328  CHECK_PNC_ERROR( ncmpi_def_var(ncid, buf, NC_DOUBLE, 2, dimids, &varid) )
1329  }
1330  } else {
1331  if ( nc_inq_dimid(ncid, dim_name, &dimid) != NC_NOERR ) // check if existed
1332  CHECK_ERROR( nc_def_dim(ncid, dim_name, dim_size, &dimid) )
1333 
1334  CHECK_ERROR( nc_def_var(ncid, name, xtype, 1, &dimid, &varid) )
1335  if ( strlen(desc)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1336  if ( strlen(units)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "units", strlen(units), units) )
1337 
1338  if ( bounds ) {
1339  dimids[0] = dimid;
1340  if ( nc_inq_dimid(ncid, "nv", &(dimids[1])) != NC_NOERR ) // first called
1341  CHECK_ERROR( nc_def_dim(ncid, "nv", 2, &(dimids[1])) )
1342  sprintf(buf, "%s_bnds", dim_name);
1343  CHECK_ERROR( nc_put_att_text(ncid, varid, "bounds", strlen(buf), buf) )
1344  CHECK_ERROR( nc_def_var(ncid, buf, NC_DOUBLE, 2, dimids, &varid) )
1345  }
1346  }
1347 
1348  return SUCCESS_CODE;
1349 }

References File_HSHORT.

Referenced by scale_file::file_def_axis().

Here is the caller graph for this function:

◆ file_write_axis_c()

int file_write_axis_c ( const int  fid,
const char *  name,
const void *  val,
const int  precision,
const int *  start,
const int *  count 
)

Definition at line 1351 of file scale_file_netcdf.c.

1357 {
1358  int ncid, varid;
1359  int ret;
1360  MPI_Offset start_[1], count_[1];
1361 
1362  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1363  ncid = files[fid]->ncid;
1364 
1365  if ( files[fid]->shared_mode )
1366  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, name, &varid) )
1367  else
1368  CHECK_ERROR( nc_inq_varid(ncid, name, &varid) )
1369 
1370  ret = file_enddef(fid, ncid);
1371  if ( ret != SUCCESS_CODE ) return ret;
1372 
1373  start_[0] = (MPI_Offset) start[0];
1374  count_[0] = (MPI_Offset) count[0];
1375 
1376  switch ( precision ) {
1377  case 8:
1378  if ( files[fid]->shared_mode )
1379  CHECK_PNC_ERROR( ncmpi_iput_vara_double(ncid, varid, start_, count_, val, NULL) )
1380  else
1381  CHECK_ERROR( nc_put_var_double(ncid, varid, (double*)val) )
1382  break;
1383  case 4:
1384  if ( files[fid]->shared_mode )
1385  CHECK_PNC_ERROR( ncmpi_iput_vara_float(ncid, varid, start_, count_, val, NULL) )
1386  else
1387  CHECK_ERROR( nc_put_var_float(ncid, varid, (float*)val) )
1388  break;
1389  default:
1390  fprintf(stderr, "unsupported data precision: %d\n", precision);
1391  return ERROR_CODE;
1392  }
1393 
1394  return SUCCESS_CODE;
1395 }

Referenced by scale_file::file_def_axis().

Here is the caller graph for this function:

◆ file_put_associatedcoordinate_c()

int file_put_associatedcoordinate_c ( const int  fid,
const char *  name,
const char *  desc,
const char *  units,
const char **  dim_names,
const int  ndims,
const int  dtype,
const void *  val,
const int  precision 
)

Definition at line 1397 of file scale_file_netcdf.c.

1406 {
1407  int ncid, *dimids, varid;
1408  nc_type xtype = -1;
1409  int i;
1410  int ret;
1411 
1412  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1413  ncid = files[fid]->ncid;
1414 
1415  if ( nc_inq_varid(ncid, name, &varid) == NC_NOERR ) // check if existed
1416  return ALREADY_EXISTED_CODE;
1417 
1418  ret = file_redef(fid, ncid);
1419  if ( ret != SUCCESS_CODE ) return ret;
1420 
1421  dimids = malloc(sizeof(int)*ndims);
1422  for (i=0; i<ndims; i++)
1423  CHECK_ERROR( nc_inq_dimid(ncid, dim_names[i], dimids+ndims-i-1) )
1424 
1425  TYPE2NCTYPE(dtype, xtype);
1426 
1427  CHECK_ERROR( nc_def_var(ncid, name, xtype, ndims, dimids, &varid) )
1428  if ( strlen(desc)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1429  if ( strlen(units)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "units", strlen(units), units) )
1430  free(dimids);
1431 
1432  ret = file_enddef(fid, ncid);
1433  if ( ret != SUCCESS_CODE ) return ret;
1434 
1435  switch ( precision ) {
1436  case 8:
1437  CHECK_ERROR( nc_put_var_double(ncid, varid, (double*)val) )
1438  break;
1439  case 4:
1440  CHECK_ERROR( nc_put_var_float(ncid, varid, (float*)val) )
1441  break;
1442  default:
1443  fprintf(stderr, "unsupported data precision: %d\n", precision);
1444  return ERROR_CODE;
1445  }
1446 
1447  return SUCCESS_CODE;
1448 }

References scale_file::i.

Referenced by scale_file::file_def_axis().

Here is the caller graph for this function:

◆ file_def_associatedcoordinate_c()

int file_def_associatedcoordinate_c ( const int  fid,
const char *  name,
const char *  desc,
const char *  units,
const char **  dim_names,
const int  ndims,
const int  dtype 
)

Definition at line 1450 of file scale_file_netcdf.c.

1457 {
1458  int ncid, varid;
1459  nc_type xtype = -1;
1460  int i;
1461  int ret;
1462  int dimids[ndims];
1463 
1464  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1465  ncid = files[fid]->ncid;
1466 
1467  if ( files[fid]->shared_mode ) {
1468  if ( ncmpi_inq_varid(ncid, name, &varid) == NC_NOERR ) // check if existed
1469  return ALREADY_EXISTED_CODE;
1470  } else {
1471  if ( nc_inq_varid(ncid, name, &varid) == NC_NOERR ) // check if existed
1472  return ALREADY_EXISTED_CODE;
1473  }
1474 
1475  ret = file_redef(fid, ncid);
1476  if ( ret != SUCCESS_CODE ) return ret;
1477 
1478  TYPE2NCTYPE(dtype, xtype);
1479 
1480  if ( files[fid]->shared_mode ) {
1481  for (i=0; i<ndims; i++)
1482  CHECK_PNC_ERROR( ncmpi_inq_dimid(ncid, dim_names[i], dimids+ndims-i-1) )
1483 
1484  CHECK_PNC_ERROR( ncmpi_def_var(ncid, name, xtype, ndims, dimids, &varid) )
1485  if ( strlen(desc) >0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1486  if ( strlen(units)>0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "units", strlen(units), units) )
1487  }
1488  else {
1489  for (i=0; i<ndims; i++)
1490  CHECK_ERROR( nc_inq_dimid(ncid, dim_names[i], dimids+ndims-i-1) )
1491 
1492  CHECK_ERROR( nc_def_var(ncid, name, xtype, ndims, dimids, &varid) )
1493  if ( strlen(desc) >0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1494  if ( strlen(units)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "units", strlen(units), units) )
1495  }
1496 
1497  return SUCCESS_CODE;
1498 }

References scale_file::i.

Referenced by scale_file::file_def_associatedcoordinate().

Here is the caller graph for this function:

◆ file_write_associatedcoordinate_c()

int file_write_associatedcoordinate_c ( const int  fid,
const char *  name,
const void *  val,
const int  ndims,
const int  precision,
const int *  start,
const int *  count 
)

Definition at line 1500 of file scale_file_netcdf.c.

1507 {
1508  int ncid, varid;
1509  int ret;
1510  MPI_Offset start_[ndims], count_[ndims];
1511  int i;
1512 
1513  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1514  ncid = files[fid]->ncid;
1515 
1516  if ( files[fid]->shared_mode )
1517  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, name, &varid) )
1518  else
1519  CHECK_ERROR( nc_inq_varid(ncid, name, &varid) )
1520 
1521  ret = file_enddef(fid, ncid);
1522  if ( ret != SUCCESS_CODE ) return ret;
1523 
1524  for (i=0; i<ndims; i++) {
1525  start_[i] = (MPI_Offset) start[i];
1526  count_[i] = (MPI_Offset) count[i];
1527  }
1528 
1529  switch ( precision ) {
1530  case 8:
1531  if ( files[fid]->shared_mode )
1532  CHECK_PNC_ERROR( ncmpi_iput_vara_double(ncid, varid, start_, count_, (double*)val, NULL) )
1533  else
1534  CHECK_ERROR( nc_put_var_double(ncid, varid, (double*)val) )
1535  break;
1536  case 4:
1537  if ( files[fid]->shared_mode )
1538  CHECK_PNC_ERROR( ncmpi_iput_vara_float(ncid, varid, start_, count_, (float*)val, NULL) )
1539  else
1540  CHECK_ERROR( nc_put_var_float(ncid, varid, (float*)val) )
1541  break;
1542  default:
1543  fprintf(stderr, "unsupported data precision: %d\n", precision);
1544  return ERROR_CODE;
1545  }
1546 
1547  return SUCCESS_CODE;
1548 }

Referenced by scale_file::file_def_associatedcoordinate().

Here is the caller graph for this function:

◆ file_add_variable_c()

int file_add_variable_c ( int *  vid,
const int  fid,
const char *  varname,
const char *  desc,
const char *  units,
const char *  stdname,
const char **  dims,
const int  ndims,
const int  dtype,
const double  tint,
const char *  tstats 
)

Definition at line 1550 of file scale_file_netcdf.c.

1561 {
1562  int ncid, varid, acid, *acdimids;
1563  int dimids[NC_MAX_DIMS], dimid;
1564  char tname[File_HSHORT+1];
1565  int tdimid, tvarid;
1566  nc_type xtype = -1;
1567  char buf[File_HMID+1];
1568  int i, j, k, m, err;
1569  int ndims_t, nndims;
1570  size_t size;
1571  double rmiss = RMISS;
1572  char coord[File_HMID+1];
1573  int has_assoc;
1574  int new;
1575  int ret;
1576 
1577  if ( nvar >= VAR_MAX ) {
1578  fprintf(stderr, "exceed max number of variable limit\n");
1579  return ERROR_CODE;
1580  }
1581 
1582  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1583  ncid = files[fid]->ncid;
1584 
1585  vars[nvar] = (varinfo_t*) malloc(sizeof(varinfo_t));
1586  vars[nvar]->ncid = ncid;
1587  vars[nvar]->t = NULL;
1588  vars[nvar]->start = NULL;
1589  vars[nvar]->count = NULL;
1590  vars[nvar]->ndims = ndims;
1591 
1592  ret = file_redef(fid, ncid);
1593  if ( ret != SUCCESS_CODE ) return ret;
1594 
1595  // get time variable
1596  if ( tint > 0.0 ) {
1597  for ( i=0; i<nt; i++ ) {
1598  if ( tdims[i] != NULL && // still opened
1599  tdims[i]->ncid == ncid && // same file
1600  tdims[i]->tint == tint ) { // same time interval
1601  vars[nvar]->t = tdims[i];
1602  break;
1603  }
1604  }
1605  if ( vars[nvar]->t == NULL ) {
1606  tdims[nt] = (tdim_t*) malloc(sizeof(tdim_t));
1607  tdims[nt]->ncid = ncid;
1608  tdims[nt]->count = -1;
1609  tdims[nt]->tint = tint;
1610  tdims[nt]->tval = (double*) malloc(sizeof(double)*NTMAX);
1611  // generate name
1612  m=0;
1613  for (i=0; i<nt; i++) {
1614  if ( tdims[i] != NULL && tdims[i]->ncid == ncid ) m++;
1615  }
1616  if ( m == 0 ) {
1617  strcpy(tname, "time");
1618  } else {
1619  sprintf(tname, "time%d", m);
1620  }
1621  strcpy(tdims[nt]->name, tname);
1622  // define time dimension and variable
1623  if ( files[fid]->shared_mode ) {
1624  CHECK_PNC_ERROR( ncmpi_def_dim(ncid, tname, 0, &tdimid) )
1625  tdims[nt]->dimid = tdimid;
1626  CHECK_PNC_ERROR( ncmpi_def_var(ncid, tname, NC_DOUBLE, 1, &tdimid, &tvarid) )
1627  tdims[nt]->varid = tvarid;
1628  strcpy(buf, "time");
1629  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, tvarid, "long_name", strlen(buf), buf) )
1630  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, tvarid, "units", strlen(files[fid]->time_units), files[fid]->time_units) )
1631  if ( strlen(files[fid]->calendar) > 0 )
1632  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, tvarid, "calendar", strlen(files[fid]->calendar), files[fid]->calendar) )
1633  // define boundary variable
1634  if ( ncmpi_inq_dimid(ncid, "nv", &(dimids[1])) != NC_NOERR ) // first called
1635  CHECK_PNC_ERROR( ncmpi_def_dim(ncid, "nv", 2, &(dimids[1])) )
1636  sprintf(buf, "%s_bnds", tname);
1637  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, tvarid, "bounds", strlen(buf), buf) )
1638  dimids[0] = tdimid;
1639  CHECK_PNC_ERROR( ncmpi_def_var(ncid, buf, NC_DOUBLE, 2, dimids, &tvarid) )
1640  tdims[nt]->bndsid = tvarid;
1641  //CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, tvarid, "units", strlen(files[fid]->time_units), files[fid]->time_units) )
1642  }
1643  else {
1644  CHECK_ERROR( nc_def_dim(ncid, tname, 0, &tdimid) )
1645  tdims[nt]->dimid = tdimid;
1646  CHECK_ERROR( nc_def_var(ncid, tname, NC_DOUBLE, 1, &tdimid, &tvarid) )
1647  tdims[nt]->varid = tvarid;
1648  strcpy(buf, "time");
1649  CHECK_ERROR( nc_put_att_text(ncid, tvarid, "long_name", strlen(buf), buf) )
1650  CHECK_ERROR( nc_put_att_text(ncid, tvarid, "units", strlen(files[fid]->time_units), files[fid]->time_units) )
1651  if ( strlen(files[fid]->calendar) > 0 )
1652  CHECK_ERROR( nc_put_att_text(ncid, tvarid, "calendar", strlen(files[fid]->calendar), files[fid]->calendar) )
1653  // define boundary variable
1654  if ( nc_inq_dimid(ncid, "nv", &(dimids[1])) != NC_NOERR ) // first called
1655  CHECK_ERROR( nc_def_dim(ncid, "nv", 2, &(dimids[1])) )
1656  sprintf(buf, "%s_bnds", tname);
1657  CHECK_ERROR( nc_put_att_text(ncid, tvarid, "bounds", strlen(buf), buf) )
1658  dimids[0] = tdimid;
1659  CHECK_ERROR( nc_def_var(ncid, buf, NC_DOUBLE, 2, dimids, &tvarid) )
1660  tdims[nt]->bndsid = tvarid;
1661  //CHECK_ERROR( nc_put_att_text(ncid, tvarid, "units", strlen(files[fid]->time_units), files[fid]->time_units) )
1662  }
1663 
1664  vars[nvar]->t = tdims[nt];
1665  nt++;
1666  }
1667  }
1668 
1669  // get dimension IDs
1670  // note: C and Fortran order are opposite
1671  ndims_t = ndims;
1672  if ( tint > 0.0 ) { // add time dimension
1673  dimids[0] = vars[nvar]->t->dimid;
1674  ndims_t++;
1675  }
1676  for (i=ndims_t-ndims; i<ndims_t; i++) dimids[i] = -1;
1677 
1678  has_assoc = 0;
1679  nndims = 0;
1680  for (i=0; i<ndims; i++) {
1681  //printf("%d %s\n", i, dims[i]);
1682  if ( files[fid]->shared_mode )
1683  err = ncmpi_inq_dimid(ncid, dims[i], &dimid);
1684  else
1685  err = nc_inq_dimid(ncid, dims[i], &dimid);
1686  if ( err == NC_NOERR ) {
1687  //printf("not assoc\n");
1688  new = 1;
1689  for (k=0; k<nndims; k++) {
1690  if (dimid == dimids[k]) {
1691  new = 0;
1692  break;
1693  }
1694  }
1695  if (new) {
1696  dimids[ndims_t-(++nndims)] = dimid;
1697  }
1698  } else {
1699  //printf("assoc\n");
1700  if ( files[fid]->shared_mode ) {
1701  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, dims[i], &acid) )
1702  CHECK_PNC_ERROR( ncmpi_inq_varndims(ncid, acid, &m) )
1703  acdimids = (int*) malloc((sizeof(int)*m));
1704  CHECK_PNC_ERROR( ncmpi_inq_vardimid(ncid, acid, acdimids) )
1705  }
1706  else {
1707  CHECK_ERROR( nc_inq_varid(ncid, dims[i], &acid) )
1708  CHECK_ERROR( nc_inq_varndims(ncid, acid, &m) )
1709  acdimids = (int*) malloc((sizeof(int)*m));
1710  CHECK_ERROR( nc_inq_vardimid(ncid, acid, acdimids) )
1711  }
1712  for (j=m-1; j>=0; j--) {
1713  new = 1;
1714  for (k=0; k<ndims_t; k++) {
1715  if (acdimids[j] == dimids[k]) {
1716  new = 0;
1717  break;
1718  }
1719  }
1720  if (new) {
1721  if ( nndims >= ndims_t ) {
1722  fprintf(stderr, "Error: invalid associated coordinates\n");
1723  return ERROR_CODE;
1724  }
1725  dimids[ndims_t-(++nndims)] = acdimids[j];
1726  //nc_inq_dimname(ncid, acdimids[j], tname);
1727  //printf("add %s\n", tname);
1728  }
1729  }
1730  free(acdimids);
1731  has_assoc = 1;
1732  }
1733  }
1734  if (nndims != ndims) {
1735  fprintf(stderr, "Error: invalid associated coordinates: %d %d %d\n", ndims_t, nndims, ndims);
1736  return ERROR_CODE;
1737  }
1738 
1739  TYPE2NCTYPE(dtype, xtype);
1740  if ( files[fid]->shared_mode ) {
1741  CHECK_PNC_ERROR( ncmpi_def_var(ncid, varname, xtype, ndims_t, dimids, &varid) )
1742  // put variable attribute
1743  if ( strlen(desc) >0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1744  if ( strlen(units) >0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "units", strlen(units), units) )
1745  if ( strlen(stdname)>0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "standard_name", strlen(stdname), stdname) )
1746 
1747 // CHECK_PNC_ERROR( ncmpi_put_att_double(ncid, varid, _FillValue, xtype, 1, &rmiss) )
1748  CHECK_PNC_ERROR( ncmpi_put_att_double(ncid, varid, "missing_value", xtype, 1, &rmiss) )
1749  }
1750  else {
1751  CHECK_ERROR( nc_def_var(ncid, varname, xtype, ndims_t, dimids, &varid) )
1752  // put variable attribute
1753  if ( strlen(desc) >0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1754  if ( strlen(units) >0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "units", strlen(units), units) )
1755  if ( strlen(stdname)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "standard_name", strlen(stdname), stdname) )
1756  CHECK_ERROR( nc_put_att_double(ncid, varid, _FillValue, xtype, 1, &rmiss) )
1757  CHECK_ERROR( nc_put_att_double(ncid, varid, "missing_value", xtype, 1, &rmiss) )
1758  }
1759  if ( has_assoc ) {
1760  strcpy(coord, dims[0]);
1761  for(i=1; i<ndims; i++) {
1762  if (strlen(coord)+strlen(dims[i])+1 < File_HMID) {
1763  strcat(coord, " ");
1764  strcat(coord, dims[i]);
1765  }
1766  }
1767  if ( ndims_t > ndims && strlen(coord)+6 < File_HMID) {
1768  strcat(coord, " ");
1769  strcat(coord, vars[nvar]->t->name);
1770  }
1771  if ( files[fid]->shared_mode )
1772  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "coordinates", strlen(coord), coord) )
1773  else
1774  CHECK_ERROR( nc_put_att_text(ncid, varid, "coordinates", strlen(coord), coord) )
1775  }
1776 
1777 
1778  if ( strcmp(tstats,"none") != 0 ) {
1779  sprintf(buf, "%s: %s", vars[nvar]->t->name, tstats);
1780  if ( files[fid]->shared_mode )
1781  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "cell_methods", strlen(buf), buf) )
1782  else
1783  CHECK_ERROR( nc_put_att_text(ncid, varid, "cell_methods", strlen(buf), buf) )
1784  }
1785 
1786  // set start and count
1787  vars[nvar]->ndims_t = ndims_t;
1788  vars[nvar]->start = (size_t*) malloc(sizeof(size_t)*ndims_t);
1789  vars[nvar]->count = (size_t*) malloc(sizeof(size_t)*ndims_t);
1790  for ( i=0; i<ndims_t; i++ ) {
1791  if ( files[fid]->shared_mode ) {
1792  MPI_Offset dimlen;
1793  CHECK_PNC_ERROR( ncmpi_inq_dimlen(ncid, dimids[i], &dimlen) )
1794  size = (size_t) dimlen;
1795  }
1796  else
1797  CHECK_ERROR( nc_inq_dimlen(ncid, dimids[i], &size) )
1798  vars[nvar]->count[i] = size;
1799  vars[nvar]->start[i] = 0;
1800  }
1801  if ( tint > 0.0 ) vars[nvar]->count[0] = 1;
1802 
1803 #ifndef NETCDF3
1804  // set chunk size and deflate level (NetCDF-4 only)
1805  if ( ! files[fid]->shared_mode && files[fid]->deflate_level > 0 ) {
1806  CHECK_ERROR( nc_def_var_chunking(ncid, varid, NC_CHUNKED, vars[nvar]->count) )
1807  CHECK_ERROR( nc_def_var_deflate(ncid, varid, 0, 1, files[fid]->deflate_level) )
1808  }
1809 #endif
1810 
1811  vars[nvar]->varid = varid;
1812  *vid = nvar;
1813  nvar++;
1814 
1815  return SUCCESS_CODE;
1816 }

Referenced by scale_file::file_def_associatedcoordinate(), and scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_enddef_c()

int file_enddef_c ( const int  fid)

Definition at line 1818 of file scale_file_netcdf.c.

1819 {
1820  int ncid;
1821 
1822  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1823 
1824  ncid = files[fid]->ncid;
1825 
1826  return file_enddef(fid, ncid);
1827 }

Referenced by scale_file::file_enddef().

Here is the caller graph for this function:

◆ file_redef_c()

int file_redef_c ( const int  fid)

Definition at line 1829 of file scale_file_netcdf.c.

1830 {
1831  int ncid;
1832 
1833  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1834 
1835  ncid = files[fid]->ncid;
1836 
1837  return file_redef(fid, ncid);
1838 }

Referenced by scale_file::file_redef().

Here is the caller graph for this function:

◆ file_attach_buffer_c()

int file_attach_buffer_c ( const int  fid,
const int64_t  buf_amount 
)

Definition at line 1840 of file scale_file_netcdf.c.

1842 {
1843  int ncid;
1844 
1845  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1846  ncid = files[fid]->ncid;
1847 
1848  if ( files[fid]->shared_mode )
1849  CHECK_PNC_ERROR( ncmpi_buffer_attach(ncid, (MPI_Offset)buf_amount) )
1850 
1851  return SUCCESS_CODE;
1852 }

Referenced by scale_file::file_attach_buffer().

Here is the caller graph for this function:

◆ file_detach_buffer_c()

int file_detach_buffer_c ( const int  fid)

Definition at line 1854 of file scale_file_netcdf.c.

1855 {
1856  int ncid;
1857 
1858  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1859  ncid = files[fid]->ncid;
1860 
1861  if ( files[fid]->shared_mode )
1862  CHECK_PNC_ERROR( ncmpi_buffer_detach(ncid) )
1863 
1864  return SUCCESS_CODE;
1865 }

Referenced by scale_file::file_detach_buffer().

Here is the caller graph for this function:

◆ file_flush_c()

int file_flush_c ( const int  fid)

Definition at line 1867 of file scale_file_netcdf.c.

1868 {
1869  int ncid;
1870 
1871  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1872  ncid = files[fid]->ncid;
1873 
1874  if ( files[fid]->shared_mode )
1875  CHECK_PNC_ERROR( ncmpi_wait_all(ncid, NC_REQ_ALL, NULL, NULL) )
1876  else
1877  CHECK_ERROR( nc_sync(ncid) )
1878 
1879  return SUCCESS_CODE;
1880 }

Referenced by scale_file::file_flush().

Here is the caller graph for this function:

◆ file_write_data_c()

int file_write_data_c ( const int  fid,
const int  vid,
const void *  var,
const double  t_start,
const double  t_end,
const int  ndims,
const int  precision,
const int *  start,
const int *  count 
)

Definition at line 1882 of file scale_file_netcdf.c.

1891 {
1892  int ncid, varid;
1893  MPI_Offset *str, *cnt;
1894  int ret;
1895 
1896  if ( vars[vid] == NULL ) return ALREADY_CLOSED_CODE;
1897  ncid = vars[vid]->ncid;
1898 
1899  if ( ndims != vars[vid]->ndims ) {
1900  fprintf(stderr, "Error: at line %d in %s\n", __LINE__, __FILE__);
1901  fprintf(stderr, " dimension size %d is not consistent that was added by file_add_variable %d\n", ndims, (int)vars[vid]->ndims );
1902  return ERROR_CODE;
1903  }
1904 
1905  ret = file_enddef(fid, ncid);
1906  if ( ret != SUCCESS_CODE ) return ret;
1907 
1908  varid = vars[vid]->varid;
1909  if ( vars[vid]->t != NULL ) { // have time dimension
1910  if ( vars[vid]->t->count < 0 || // first time
1911  t_end > vars[vid]->t->t + TEPS ) { // time goes next step
1912  vars[vid]->t->count += 1;
1913  vars[vid]->t->t = t_end;
1914  if ( vars[vid]->t->count > NTMAX-1 ) {
1915  fprintf(stderr, "time count exceeds the max limit (%d)\n", NTMAX);
1916  return ERROR_CODE;
1917  }
1918  vars[vid]->t->tval[vars[vid]->t->count] = t_end;
1919  if ( files[fid]->shared_mode ) { // write a new value to variable time
1920  MPI_Offset index[2];
1921  index[0] = (MPI_Offset) vars[vid]->t->count;
1922  CHECK_PNC_ERROR( ncmpi_put_var1_double_all(ncid, vars[vid]->t->varid, index, &t_end) )
1923  index[1] = 0;
1924  CHECK_PNC_ERROR( ncmpi_put_var1_double_all(ncid, vars[vid]->t->bndsid, index, &t_start ) )
1925  index[1] = 1;
1926  CHECK_PNC_ERROR( ncmpi_put_var1_double_all(ncid, vars[vid]->t->bndsid, index, &t_end ) )
1927  } else {
1928  size_t index[2];
1929  index[0] = vars[vid]->t->count;
1930  CHECK_ERROR( nc_put_var1_double(ncid, vars[vid]->t->varid, index, &t_end) )
1931  index[1] = 0;
1932  CHECK_ERROR( nc_put_var1_double(ncid, vars[vid]->t->bndsid, index, &t_start) )
1933  index[1] = 1;
1934  CHECK_ERROR( nc_put_var1_double(ncid, vars[vid]->t->bndsid, index, &t_end) )
1935  }
1936  vars[vid]->start[0] = vars[vid]->t->count;
1937  } else {
1938  size_t nt = vars[vid]->t->count + 1;
1939  int flag, n;
1940  flag = 1;
1941  for(n=nt-1;n>=0;n--) {
1942  if ( fabs(vars[vid]->t->tval[n]-t_end) < TEPS ) {
1943  vars[vid]->start[0] = n;
1944  flag = 0;
1945  break;
1946  }
1947  }
1948  if ( flag ) {
1949  fprintf(stderr, "cannot find time: %f\n", t_end);
1950  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);
1951  fprintf(stderr, " time is: ");
1952  for (n=0;n<nt;n++) fprintf(stderr, "%f, ", vars[vid]->t->tval[n]);
1953  fprintf(stderr, "\n");
1954  return ERROR_CODE;
1955  }
1956  }
1957  }
1958 
1959  if ( files[fid]->shared_mode ) {
1960  int i;
1961  int ndims_t = vars[vid]->ndims_t;
1962  str = (MPI_Offset*) malloc(sizeof(MPI_Offset)*(ndims_t));
1963  cnt = (MPI_Offset*) malloc(sizeof(MPI_Offset)*(ndims_t));
1964  if ( vars[vid]->t != NULL ) { // have time dimension
1965  // add time dimension to start[0] and count[0]
1966  str[0] = vars[vid]->start[0]; // start along the time dimension
1967  cnt[0] = vars[vid]->count[0];
1968  for (i=0; i<ndims; i++) {
1969  str[ndims_t-i-1] = start[i] - 1;
1970  cnt[ndims_t-i-1] = count[i];
1971  }
1972  } else {
1973  for (i=0; i<ndims; i++) {
1974  str[ndims-i-1] = start[i] - 1;
1975  cnt[ndims-i-1] = count[i];
1976  }
1977  }
1978  }
1979 
1980  switch (precision) {
1981  case 8:
1982  if ( files[fid]->shared_mode )
1983  CHECK_PNC_ERROR( ncmpi_bput_vara_double(ncid, varid, str, cnt, (double*)var, NULL) )
1984  else
1985  CHECK_ERROR( nc_put_vara_double(ncid, varid, vars[vid]->start, vars[vid]->count, (double*)var) )
1986  break;
1987  case 4:
1988  if ( files[fid]->shared_mode )
1989  CHECK_PNC_ERROR( ncmpi_bput_vara_float(ncid, varid, str, cnt, (float*)var, NULL) )
1990  else
1991  CHECK_ERROR( nc_put_vara_float(ncid, varid, vars[vid]->start, vars[vid]->count, (float*)var) )
1992  break;
1993  default:
1994  fprintf(stderr, "unsupported data precision: %d\n", precision);
1995  return ERROR_CODE;
1996  }
1997 
1998  if ( files[fid]->shared_mode) {
1999  free(str);
2000  free(cnt);
2001  }
2002 
2003  return SUCCESS_CODE;
2004 }

Referenced by scale_file::file_get_stepsize().

Here is the caller graph for this function:

◆ file_close_c()

int file_close_c ( const int  fid,
const bool  abort 
)

Definition at line 2006 of file scale_file_netcdf.c.

2008 {
2009  int ncid;
2010  int i;
2011 
2012  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
2013  ncid = files[fid]->ncid;
2014 
2015  for (i=0; i<nvar; i++) {
2016  if ( vars[i] != NULL && vars[i]->ncid == ncid ) {
2017  free( vars[i]->start );
2018  free( vars[i]->count );
2019  free( vars[i] );
2020  vars[i] = NULL;
2021  }
2022  }
2023 
2024  for (i=0; i<nt; i++) {
2025  if ( tdims[i] != NULL && tdims[i]->ncid == ncid ) {
2026  free( tdims[i]->tval );
2027  free( tdims[i] );
2028  tdims[i] = NULL;
2029  }
2030  }
2031 
2032  if ( files[fid]->shared_mode ) {
2033  if ( ! abort ) CHECK_PNC_ERROR( ncmpi_close(ncid) )
2034  } else
2035  CHECK_ERROR( nc_close(ncid) )
2036 
2037  free( files[fid] );
2038  files[fid] = NULL;
2039 
2040  return SUCCESS_CODE;
2041 }

Referenced by scale_file::file_close().

Here is the caller graph for this function:
CHECK_PNC_ERROR
#define CHECK_PNC_ERROR(func)
Definition: scale_file_netcdf.c:40
varinfo_t
Definition: scale_file_netcdf.c:117
tdim_t::tint
real64_t tint
Definition: scale_file_netcdf.c:112
File_INTEGER2
#define File_INTEGER2
Definition: scale_file_const.h:9
datainfo_t::time_end
real64_t time_end
Definition: scale_file.h:32
varinfo_t::ndims
size_t ndims
Definition: scale_file_netcdf.c:123
tdim_t
Definition: scale_file_netcdf.c:105
TEPS
#define TEPS
Definition: scale_file_netcdf.c:7
tdim_t::count
int count
Definition: scale_file_netcdf.c:110
scale_file::file_enddef
subroutine, public file_enddef(fid)
Definition: scale_file.F90:6061
mod_atmos_bnd_driver::index
integer, allocatable, public index
Definition: mod_atmos_bnd_driver.F90:43
File_FAPPEND
#define File_FAPPEND
Definition: scale_file_const.h:17
MIN
#define MIN(a, b)
Definition: scale_file_netcdf.c:10
File_TEXT
#define File_TEXT
Definition: scale_file_const.h:12
ncmpi_inq_varid
#define ncmpi_inq_varid(a, b, c)
Definition: scale_file_netcdf.c:47
SUCCESS_CODE
#define SUCCESS_CODE
Definition: scale_file_const.h:21
datainfo_t::att_type
int32_t att_type[ATT_MAX]
Definition: scale_file.h:37
tdim_t::t
real64_t t
Definition: scale_file_netcdf.c:111
varinfo_t::varid
int varid
Definition: scale_file_netcdf.c:119
FILE_MAX
#define FILE_MAX
Definition: scale_file_const.h:26
datainfo_t::att_len
int32_t att_len[ATT_MAX]
Definition: scale_file.h:38
NCTYPE2TYPE
#define NCTYPE2TYPE(nctype, type)
Definition: scale_file_netcdf.c:51
File_FREAD
#define File_FREAD
Definition: scale_file_const.h:15
tdim_t::bndsid
int bndsid
Definition: scale_file_netcdf.c:109
varinfo_t::ncid
int ncid
Definition: scale_file_netcdf.c:118
File_REAL4
#define File_REAL4
Definition: scale_file_const.h:7
ncmpi_inq_dimid
#define ncmpi_inq_dimid(a, b, c)
Definition: scale_file_netcdf.c:48
float
typedef float(real32_t)
datainfo_t::calendar
char calendar[File_HSHORT]
Definition: scale_file.h:34
tdim_t::varid
int varid
Definition: scale_file_netcdf.c:108
datainfo_t::has_tdim
_Bool has_tdim
Definition: scale_file.h:30
NTMAX
#define NTMAX
Definition: scale_file_netcdf.c:8
mod_atmos_vars::j
real(rp), allocatable, target, public j
Definition: mod_atmos_vars.F90:141
datainfo_t::dim_size
int32_t dim_size[RANK_MAX]
Definition: scale_file.h:28
varinfo_t::t
tdim_t * t
Definition: scale_file_netcdf.c:120
RANK_MAX
#define RANK_MAX
Definition: scale_file_const.h:30
datainfo_t::att_name
char att_name[File_HSHORT *ATT_MAX]
Definition: scale_file.h:36
datainfo_t::fid
int32_t fid
Definition: scale_file.h:39
datainfo_t::units
char units[File_HSHORT]
Definition: scale_file.h:23
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:45
VAR_MAX
#define VAR_MAX
Definition: scale_file_const.h:28
File_HLONG
#define File_HLONG
Definition: scale_file_const.h:4
scale_file::i
logical, public i
Definition: scale_file.F90:196
datainfo_t::standard_name
char standard_name[File_HMID]
Definition: scale_file.h:24
datainfo_t::datatype
int32_t datatype
Definition: scale_file.h:25
tdim_t::ncid
int ncid
Definition: scale_file_netcdf.c:106
datainfo_t::description
char description[File_HMID]
Definition: scale_file.h:22
datainfo_t::time_units
char time_units[File_HMID]
Definition: scale_file.h:33
scale_coriolis::type
character(len=h_short), public type
Definition: scale_coriolis.F90:37
File_FWRITE
#define File_FWRITE
Definition: scale_file_const.h:16
File_HSHORT
#define File_HSHORT
Definition: scale_file_const.h:2
datainfo_t::rank
int32_t rank
Definition: scale_file.h:26
tdim_t::tval
real64_t * tval
Definition: scale_file_netcdf.c:113
ALREADY_CLOSED_CODE
#define ALREADY_CLOSED_CODE
Definition: scale_file_const.h:22
tdim_t::dimid
int dimid
Definition: scale_file_netcdf.c:107
scale_tracer::offset
real(rp), public offset
Definition: scale_tracer.F90:45
datainfo_t::natts
int32_t natts
Definition: scale_file.h:35
ERROR_CODE
#define ERROR_CODE
Definition: scale_file_const.h:20
File_HMID
#define File_HMID
Definition: scale_file_const.h:3
datainfo_t::dim_name
char dim_name[File_HSHORT *RANK_MAX]
Definition: scale_file.h:27
datainfo_t::time_start
real64_t time_start
Definition: scale_file.h:31
varinfo_t::start
size_t * start
Definition: scale_file_netcdf.c:121
scale_file::file_redef
subroutine, public file_redef(fid)
Definition: scale_file.F90:6093
ALREADY_EXISTED_CODE
#define ALREADY_EXISTED_CODE
Definition: scale_file_const.h:23
CHECK_ERROR
#define CHECK_ERROR(func)
Definition: scale_file_netcdf.c:14
ncmpi_inq_attid
#define ncmpi_inq_attid(a, b, c, d)
Definition: scale_file_netcdf.c:46
fileinfo_t
Definition: scale_file_netcdf.c:93
datainfo_t::step
int32_t step
Definition: scale_file.h:29
datainfo_t::varname
char varname[File_HSHORT]
Definition: scale_file.h:21
File_INTEGER4
#define File_INTEGER4
Definition: scale_file_const.h:10
DEFAULT_DEFLATE_LEVEL
#define DEFAULT_DEFLATE_LEVEL
Definition: scale_file_netcdf.c:91
scale_tracer::name
character(len=h_short), public name
Definition: scale_tracer.F90:39
varinfo_t::ndims_t
size_t ndims_t
Definition: scale_file_netcdf.c:124
TYPE2NCTYPE
#define TYPE2NCTYPE(type, nctype)
Definition: scale_file_netcdf.c:75
varinfo_t::count
size_t * count
Definition: scale_file_netcdf.c:122
RMISS
#define RMISS
Definition: scale_file_const.h:35
File_REAL8
#define File_REAL8
Definition: scale_file_const.h:8