SCALE-RM
Classes | Functions
scale_file.h File Reference
#include <stdio.h>
#include <stdint.h>
#include <stdlib.h>
#include <string.h>
#include <math.h>
#include <mpi.h>
#include "scale_file_const.h"
Include dependency graph for scale_file.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

struct  datainfo_t
 

Functions

typedef float (real32_t)
 
typedef double (real64_t)
 
int32_t file_open_c (int32_t *fid, const char *fname, const int32_t mode, const MPI_Comm comm)
 
int32_t file_get_dim_length_c (const int32_t fid, const char *dimname, int32_t *len)
 
int32_t file_set_option_c (const int32_t fid, const char *filetype, const char *key, const char *val)
 
int32_t file_get_nvars_c (const int32_t fid, int32_t *nvars)
 
int32_t file_get_varname_c (const int32_t fid, const int32_t vid, char *name, const int32_t len)
 
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_get_step_size_c (const int32_t fid, const char *varname, int32_t *len)
 
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)
 
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)
 
int32_t file_get_attribute_int_c (const int32_t fid, const char *vname, const char *key, const int32_t suppress, int32_t *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_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_set_attribute_text_c (const int32_t fid, const char *vname, const char *key, const char *value)
 
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)
 
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_set_attribute_double_c (const int32_t fid, const char *vname, const char *key, const double *value, const size_t len)
 
int32_t file_add_associatedvariable_c (const int32_t fid, const char *vname)
 
int32_t file_set_tunits_c (const int32_t fid, const char *time_units, const char *calendar)
 
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_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)
 
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_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)
 
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_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_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)
 
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_enddef_c (const int32_t fid)
 
int32_t file_redef_c (const int32_t fid)
 
int32_t file_attach_buffer_c (const int32_t fid, const int64_t buf_amount)
 
int32_t file_detach_buffer_c (const int32_t fid)
 
int32_t file_flush_c (const int32_t fid)
 
int32_t file_close_c (const int32_t fid, const int32_t abort)
 

Function Documentation

◆ float()

typedef float ( real32_t  )

◆ double()

typedef double ( real64_t  )

◆ file_open_c()

int32_t file_open_c ( int32_t *  fid,
const char *  fname,
const int32_t  mode,
const MPI_Comm  comm 
)

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 
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 }

References File_HLONG.

Referenced by scale_file::file_get_aggregate().

Here is the caller graph for this function:

◆ file_get_dim_length_c()

int32_t file_get_dim_length_c ( const int32_t  fid,
const char *  dimname,
int32_t *  len 
)

Definition at line 256 of file scale_file_netcdf.c.

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 }

Referenced by scale_file::file_get_dimlength().

Here is the caller graph for this function:

◆ file_set_option_c()

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

Definition at line 280 of file scale_file_netcdf.c.

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 }

References SUCCESS_CODE.

Referenced by scale_file::file_set_option().

Here is the caller graph for this function:

◆ file_get_nvars_c()

int32_t file_get_nvars_c ( const int32_t  fid,
int32_t *  nvars 
)

Definition at line 296 of file scale_file_netcdf.c.

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 }

Referenced by scale_file::file_create().

Here is the caller graph for this function:

◆ file_get_varname_c()

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

Definition at line 313 of file scale_file_netcdf.c.

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 }

References scale_file::i.

Referenced by scale_file::file_create().

Here is the caller graph for this function:

◆ file_get_datainfo_c()

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

Definition at line 338 of file scale_file_netcdf.c.

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 }

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()

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

Definition at line 633 of file scale_file_netcdf.c.

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 }

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()

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 
)

Definition at line 687 of file scale_file_netcdf.c.

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

Referenced by scale_file::file_get_stepsize().

Here is the caller graph for this function:

◆ file_get_attribute_text_c()

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 
)

Definition at line 852 of file scale_file_netcdf.c.

858 {
859  int ncid;
860  int varid;
861 
862  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
863  ncid = files[fid]->ncid;
864 
865  ERROR_SUPPRESS = suppress;
866 
867  if ( files[fid]->shared_mode ) {
868  MPI_Offset l;
869  if ( strcmp(vname, "global") == 0 ) {
870  varid = NC_GLOBAL;
871  } else
872  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
873 
874  CHECK_PNC_ERROR( ncmpi_inq_attlen(ncid, varid, key, &l) )
875  if ( len < l ) return ERROR_CODE;
876 
877  CHECK_PNC_ERROR( ncmpi_get_att_text(ncid, varid, key, value) )
878  value[l] = '\0';
879  }
880  else {
881  size_t l;
882  if ( strcmp(vname, "global") == 0 ) {
883  varid = NC_GLOBAL;
884  } else
885  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
886 
887  CHECK_ERROR( nc_inq_attlen(ncid, varid, key, &l) )
888  if ( len < l ) return ERROR_CODE;
889 
890  CHECK_ERROR( nc_get_att_text(ncid, varid, key, value) )
891  value[l] = '\0';
892  }
893 
894  ERROR_SUPPRESS = 0;
895 
896  return SUCCESS_CODE;
897 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_get_attribute_int_c()

int32_t file_get_attribute_int_c ( const int32_t  fid,
const char *  vname,
const char *  key,
const int32_t  suppress,
int32_t *  value,
const size_t  len 
)

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_get_attribute_float_c()

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 
)

Definition at line 942 of file scale_file_netcdf.c.

948 {
949  int ncid;
950  int varid;
951 
952  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
953  ncid = files[fid]->ncid;
954 
955  ERROR_SUPPRESS = suppress;
956 
957  if ( files[fid]->shared_mode ) {
958  MPI_Offset l;
959  if ( strcmp(vname, "global") == 0 ) {
960  varid = NC_GLOBAL;
961  } else
962  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
963 
964  CHECK_PNC_ERROR( ncmpi_inq_attlen(ncid, varid, key, &l) )
965  if ( len < l ) return ERROR_CODE;
966  CHECK_PNC_ERROR( ncmpi_get_att_float(ncid, varid, key, value) )
967  }
968  else {
969  size_t l;
970  if ( strcmp(vname, "global") == 0 ) {
971  varid = NC_GLOBAL;
972  } else
973  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
974 
975  CHECK_ERROR( nc_inq_attlen(ncid, varid, key, &l) )
976  if ( len < l ) return ERROR_CODE;
977  CHECK_ERROR( nc_get_att_float(ncid, varid, key, value) )
978  }
979 
980  ERROR_SUPPRESS = 0;
981 
982  return SUCCESS_CODE;
983 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_get_attribute_double_c()

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 
)

Definition at line 985 of file scale_file_netcdf.c.

991 {
992  int ncid;
993  int varid;
994 
995  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
996  ncid = files[fid]->ncid;
997 
998  ERROR_SUPPRESS = suppress;
999 
1000  if ( files[fid]->shared_mode ) {
1001  MPI_Offset l;
1002  if ( strcmp(vname, "global") == 0 ) {
1003  varid = NC_GLOBAL;
1004  } else
1005  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
1006 
1007  CHECK_PNC_ERROR( ncmpi_inq_attlen(ncid, varid, key, &l) )
1008  if ( len < l ) return ERROR_CODE;
1009  CHECK_PNC_ERROR( ncmpi_get_att_double(ncid, varid, key, value) )
1010  }
1011  else {
1012  size_t l;
1013  if ( strcmp(vname, "global") == 0 ) {
1014  varid = NC_GLOBAL;
1015  } else
1016  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
1017 
1018  CHECK_ERROR( nc_inq_attlen(ncid, varid, key, &l) )
1019  if ( len < l ) return ERROR_CODE;
1020  CHECK_ERROR( nc_get_att_double(ncid, varid, key, value) )
1021  }
1022 
1023  ERROR_SUPPRESS = 0;
1024 
1025  return SUCCESS_CODE;
1026 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_set_attribute_text_c()

int32_t file_set_attribute_text_c ( const int32_t  fid,
const char *  vname,
const char *  key,
const char *  value 
)

Definition at line 1028 of file scale_file_netcdf.c.

1032 {
1033  int ncid;
1034  int varid;
1035  int attid;
1036  int32_t ret;
1037 
1038  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1039  ncid = files[fid]->ncid;
1040 
1041  if ( files[fid]->shared_mode ) {
1042  if ( strcmp(vname, "global") == 0 ) {
1043  varid = NC_GLOBAL;
1044  } else
1045  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
1046 
1047  if ( ncmpi_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1048  }
1049  else {
1050  if ( strcmp(vname, "global") == 0 ) {
1051  varid = NC_GLOBAL;
1052  } else
1053  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
1054 
1055  if ( nc_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1056  }
1057 
1058  ret = file_redef(fid, ncid);
1059  if ( ret != SUCCESS_CODE ) return ret;
1060 
1061  if ( files[fid]->shared_mode )
1062  CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, key, strlen(val), val) )
1063  else
1064  CHECK_ERROR( nc_put_att_text(ncid, varid, key, strlen(val), val) )
1065 
1066  return SUCCESS_CODE;
1067 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_set_attribute_int_c()

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 
)

Definition at line 1069 of file scale_file_netcdf.c.

1074 {
1075  int ncid;
1076  int varid;
1077  int attid;
1078  int32_t ret;
1079 
1080  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1081  ncid = files[fid]->ncid;
1082 
1083  if ( files[fid]->shared_mode ) {
1084  if ( strcmp(vname, "global") == 0 ) {
1085  varid = NC_GLOBAL;
1086  } else
1087  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
1088 
1089  if ( ncmpi_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1090  }
1091  else {
1092  if ( strcmp(vname, "global") == 0 ) {
1093  varid = NC_GLOBAL;
1094  } else
1095  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
1096 
1097  if ( nc_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1098  }
1099 
1100  ret = file_redef(fid, ncid);
1101  if ( ret != SUCCESS_CODE ) return ret;
1102 
1103  if ( files[fid]->shared_mode )
1104  CHECK_PNC_ERROR( ncmpi_put_att_int(ncid, varid, key, NC_INT, len, value) )
1105  else
1106  CHECK_ERROR( nc_put_att_int(ncid, varid, key, NC_INT, len, value) )
1107 
1108 
1109  return SUCCESS_CODE;
1110 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_set_attribute_float_c()

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

Definition at line 1112 of file scale_file_netcdf.c.

1117 {
1118  int ncid;
1119  int varid;
1120  int attid;
1121  int32_t ret;
1122 
1123  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1124  ncid = files[fid]->ncid;
1125 
1126  if ( files[fid]->shared_mode ) {
1127  if ( strcmp(vname, "global") == 0 ) {
1128  varid = NC_GLOBAL;
1129  } else
1130  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
1131 
1132  if ( ncmpi_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1133  }
1134  else {
1135  if ( strcmp(vname, "global") == 0 ) {
1136  varid = NC_GLOBAL;
1137  } else
1138  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
1139 
1140  if ( nc_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1141  }
1142 
1143  ret = file_redef(fid, ncid);
1144  if ( ret != SUCCESS_CODE ) return ret;
1145 
1146  if ( files[fid]->shared_mode )
1147  CHECK_PNC_ERROR( ncmpi_put_att_float(ncid, varid, key, NC_FLOAT, len, value) )
1148  else
1149  CHECK_ERROR( nc_put_att_float(ncid, varid, key, NC_FLOAT, len, value) )
1150 
1151  return SUCCESS_CODE;
1152 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_set_attribute_double_c()

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

Definition at line 1154 of file scale_file_netcdf.c.

1159 {
1160  int ncid;
1161  int varid;
1162  int attid;
1163 
1164  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1165  ncid = files[fid]->ncid;
1166 
1167  if ( files[fid]->shared_mode ) {
1168  if ( strcmp(vname, "global") == 0 ) {
1169  varid = NC_GLOBAL;
1170  } else
1171  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, vname, &varid) )
1172 
1173  if ( ncmpi_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1174  }
1175  else {
1176  if ( strcmp(vname, "global") == 0 ) {
1177  varid = NC_GLOBAL;
1178  } else
1179  CHECK_ERROR( nc_inq_varid(ncid, vname, &varid) )
1180 
1181  if ( nc_inq_attid(ncid, varid, key, &attid) == NC_NOERR ) return ALREADY_EXISTED_CODE;
1182  }
1183 
1184  if ( files[fid]->shared_mode )
1185  CHECK_PNC_ERROR( ncmpi_put_att_double(ncid, varid, key, NC_DOUBLE, len, value) )
1186  else
1187  CHECK_ERROR( nc_put_att_double(ncid, varid, key, NC_DOUBLE, len, value) )
1188 
1189  return SUCCESS_CODE;
1190 }

Referenced by scale_file::file_def_variable().

Here is the caller graph for this function:

◆ file_add_associatedvariable_c()

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

Definition at line 1192 of file scale_file_netcdf.c.

1194 {
1195  int ncid, varid;
1196  int32_t ret;
1197 
1198  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1199  ncid = files[fid]->ncid;
1200 
1201  if ( nc_inq_varid(ncid, vname, &varid) == NC_NOERR ) // check if existed
1202  return ALREADY_EXISTED_CODE;
1203 
1204  ret = file_redef(fid, ncid);
1205  if ( ret != SUCCESS_CODE ) return ret;
1206 
1207  if ( files[fid]->shared_mode )
1208  CHECK_PNC_ERROR( ncmpi_def_var(ncid, vname, NC_INT, 0, 0, &varid) )
1209  else
1210  CHECK_ERROR( nc_def_var(ncid, vname, NC_INT, 0, 0, &varid) )
1211 
1212  return SUCCESS_CODE;
1213 }

Referenced by scale_file::file_add_associatedvariable().

Here is the caller graph for this function:

◆ file_set_tunits_c()

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

Definition at line 1215 of file scale_file_netcdf.c.

1218 {
1219  strcpy(files[fid]->time_units, time_units);
1220  strcpy(files[fid]->calendar, calendar);
1221 
1222  return SUCCESS_CODE;
1223 }

Referenced by scale_file::file_create().

Here is the caller graph for this function:

◆ file_put_axis_c()

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 
)

Definition at line 1225 of file scale_file_netcdf.c.

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

Referenced by scale_file::file_get_dimlength().

Here is the caller graph for this function:

◆ file_def_axis_c()

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 
)

Definition at line 1274 of file scale_file_netcdf.c.

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

References File_HSHORT.

Referenced by scale_file::file_def_axis().

Here is the caller graph for this function:

◆ file_write_axis_c()

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 
)

Definition at line 1339 of file scale_file_netcdf.c.

1345 {
1346  int ncid, varid;
1347  int32_t ret;
1348 
1349  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1350  ncid = files[fid]->ncid;
1351 
1352  if ( files[fid]->shared_mode )
1353  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, name, &varid) )
1354  else
1355  CHECK_ERROR( nc_inq_varid(ncid, name, &varid) )
1356 
1357  ret = file_enddef(fid, ncid);
1358  if ( ret != SUCCESS_CODE ) return ret;
1359 
1360  switch ( precision ) {
1361  case 8:
1362  if ( files[fid]->shared_mode )
1363  CHECK_PNC_ERROR( ncmpi_iput_vara_double(ncid, varid, start, count, val, NULL) )
1364  else
1365  CHECK_ERROR( nc_put_var_double(ncid, varid, (double*)val) )
1366  break;
1367  case 4:
1368  if ( files[fid]->shared_mode )
1369  CHECK_PNC_ERROR( ncmpi_iput_vara_float(ncid, varid, start, count, val, NULL) )
1370  else
1371  CHECK_ERROR( nc_put_var_float(ncid, varid, (float*)val) )
1372  break;
1373  default:
1374  fprintf(stderr, "unsupported data precision: %d\n", precision);
1375  return ERROR_CODE;
1376  }
1377 
1378  return SUCCESS_CODE;
1379 }

Referenced by scale_file::file_def_axis().

Here is the caller graph for this function:

◆ file_put_associatedcoordinate_c()

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 
)

Definition at line 1381 of file scale_file_netcdf.c.

1390 {
1391  int ncid, *dimids, varid;
1392  nc_type xtype = -1;
1393  int i;
1394  int32_t ret;
1395 
1396  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1397  ncid = files[fid]->ncid;
1398 
1399  if ( nc_inq_varid(ncid, name, &varid) == NC_NOERR ) // check if existed
1400  return ALREADY_EXISTED_CODE;
1401 
1402  ret = file_redef(fid, ncid);
1403  if ( ret != SUCCESS_CODE ) return ret;
1404 
1405  dimids = malloc(sizeof(int)*ndims);
1406  for (i=0; i<ndims; i++)
1407  CHECK_ERROR( nc_inq_dimid(ncid, dim_names[i], dimids+ndims-i-1) )
1408 
1409  TYPE2NCTYPE(dtype, xtype);
1410 
1411  CHECK_ERROR( nc_def_var(ncid, name, xtype, ndims, dimids, &varid) )
1412  if ( strlen(desc)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1413  if ( strlen(units)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "units", strlen(units), units) )
1414  free(dimids);
1415 
1416  ret = file_enddef(fid, ncid);
1417  if ( ret != SUCCESS_CODE ) return ret;
1418 
1419  switch ( precision ) {
1420  case 8:
1421  CHECK_ERROR( nc_put_var_double(ncid, varid, (double*)val) )
1422  break;
1423  case 4:
1424  CHECK_ERROR( nc_put_var_float(ncid, varid, (float*)val) )
1425  break;
1426  default:
1427  fprintf(stderr, "unsupported data precision: %d\n", precision);
1428  return ERROR_CODE;
1429  }
1430 
1431  return SUCCESS_CODE;
1432 }

References scale_file::i.

Referenced by scale_file::file_def_axis().

Here is the caller graph for this function:

◆ file_def_associatedcoordinate_c()

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 
)

Definition at line 1434 of file scale_file_netcdf.c.

1441 {
1442  int ncid, varid;
1443  nc_type xtype = -1;
1444  int i;
1445  int32_t ret;
1446  int dimids[ndims];
1447 
1448  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1449  ncid = files[fid]->ncid;
1450 
1451  if ( files[fid]->shared_mode ) {
1452  if ( ncmpi_inq_varid(ncid, name, &varid) == NC_NOERR ) // check if existed
1453  return ALREADY_EXISTED_CODE;
1454  } else {
1455  if ( nc_inq_varid(ncid, name, &varid) == NC_NOERR ) // check if existed
1456  return ALREADY_EXISTED_CODE;
1457  }
1458 
1459  ret = file_redef(fid, ncid);
1460  if ( ret != SUCCESS_CODE ) return ret;
1461 
1462  TYPE2NCTYPE(dtype, xtype);
1463 
1464  if ( files[fid]->shared_mode ) {
1465  for (i=0; i<ndims; i++)
1466  CHECK_PNC_ERROR( ncmpi_inq_dimid(ncid, dim_names[i], dimids+ndims-i-1) )
1467 
1468  CHECK_PNC_ERROR( ncmpi_def_var(ncid, name, xtype, ndims, dimids, &varid) )
1469  if ( strlen(desc) >0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1470  if ( strlen(units)>0 ) CHECK_PNC_ERROR( ncmpi_put_att_text(ncid, varid, "units", strlen(units), units) )
1471  }
1472  else {
1473  for (i=0; i<ndims; i++)
1474  CHECK_ERROR( nc_inq_dimid(ncid, dim_names[i], dimids+ndims-i-1) )
1475 
1476  CHECK_ERROR( nc_def_var(ncid, name, xtype, ndims, dimids, &varid) )
1477  if ( strlen(desc) >0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "long_name", strlen(desc), desc) )
1478  if ( strlen(units)>0 ) CHECK_ERROR( nc_put_att_text(ncid, varid, "units", strlen(units), units) )
1479  }
1480 
1481  return SUCCESS_CODE;
1482 }

References scale_file::i.

Referenced by scale_file::file_def_associatedcoordinate().

Here is the caller graph for this function:

◆ file_write_associatedcoordinate_c()

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 
)

Definition at line 1484 of file scale_file_netcdf.c.

1490 {
1491  int ncid, varid;
1492  int32_t ret;
1493 
1494  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1495  ncid = files[fid]->ncid;
1496 
1497  if ( files[fid]->shared_mode )
1498  CHECK_PNC_ERROR( ncmpi_inq_varid(ncid, name, &varid) )
1499  else
1500  CHECK_ERROR( nc_inq_varid(ncid, name, &varid) )
1501 
1502  ret = file_enddef(fid, ncid);
1503  if ( ret != SUCCESS_CODE ) return ret;
1504 
1505  switch ( precision ) {
1506  case 8:
1507  if ( files[fid]->shared_mode )
1508  CHECK_PNC_ERROR( ncmpi_iput_vara_double(ncid, varid, start, count, (double*)val, NULL) )
1509  else
1510  CHECK_ERROR( nc_put_var_double(ncid, varid, (double*)val) )
1511  break;
1512  case 4:
1513  if ( files[fid]->shared_mode )
1514  CHECK_PNC_ERROR( ncmpi_iput_vara_float(ncid, varid, start, count, (float*)val, NULL) )
1515  else
1516  CHECK_ERROR( nc_put_var_float(ncid, varid, (float*)val) )
1517  break;
1518  default:
1519  fprintf(stderr, "unsupported data precision: %d\n", precision);
1520  return ERROR_CODE;
1521  }
1522 
1523  return SUCCESS_CODE;
1524 }

Referenced by scale_file::file_def_associatedcoordinate().

Here is the caller graph for this function:

◆ file_add_variable_c()

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 
)

Definition at line 1526 of file scale_file_netcdf.c.

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

References File_HMID, File_HSHORT, scale_file::i, mod_atmos_vars::j, and scale_tracer::k.

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

Here is the caller graph for this function:

◆ file_write_data_c()

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 
)

Definition at line 1858 of file scale_file_netcdf.c.

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

Referenced by scale_file::file_get_stepsize().

Here is the caller graph for this function:

◆ file_enddef_c()

int32_t file_enddef_c ( const int32_t  fid)

Definition at line 1794 of file scale_file_netcdf.c.

1795 {
1796  int ncid;
1797 
1798  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1799 
1800  ncid = files[fid]->ncid;
1801 
1802  return file_enddef(fid, ncid);
1803 }

Referenced by scale_file::file_enddef().

Here is the caller graph for this function:

◆ file_redef_c()

int32_t file_redef_c ( const int32_t  fid)

Definition at line 1805 of file scale_file_netcdf.c.

1806 {
1807  int ncid;
1808 
1809  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1810 
1811  ncid = files[fid]->ncid;
1812 
1813  return file_redef(fid, ncid);
1814 }

Referenced by scale_file::file_redef().

Here is the caller graph for this function:

◆ file_attach_buffer_c()

int32_t file_attach_buffer_c ( const int32_t  fid,
const int64_t  buf_amount 
)

Definition at line 1816 of file scale_file_netcdf.c.

1818 {
1819  int ncid;
1820 
1821  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1822  ncid = files[fid]->ncid;
1823 
1824  if ( files[fid]->shared_mode )
1825  CHECK_PNC_ERROR( ncmpi_buffer_attach(ncid, (MPI_Offset)buf_amount) )
1826 
1827  return SUCCESS_CODE;
1828 }

Referenced by scale_file::file_attach_buffer().

Here is the caller graph for this function:

◆ file_detach_buffer_c()

int32_t file_detach_buffer_c ( const int32_t  fid)

Definition at line 1830 of file scale_file_netcdf.c.

1831 {
1832  int ncid;
1833 
1834  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1835  ncid = files[fid]->ncid;
1836 
1837  if ( files[fid]->shared_mode )
1838  CHECK_PNC_ERROR( ncmpi_buffer_detach(ncid) )
1839 
1840  return SUCCESS_CODE;
1841 }

Referenced by scale_file::file_detach_buffer().

Here is the caller graph for this function:

◆ file_flush_c()

int32_t file_flush_c ( const int32_t  fid)

Definition at line 1843 of file scale_file_netcdf.c.

1844 {
1845  int ncid;
1846 
1847  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1848  ncid = files[fid]->ncid;
1849 
1850  if ( files[fid]->shared_mode )
1851  CHECK_PNC_ERROR( ncmpi_wait_all(ncid, NC_REQ_ALL, NULL, NULL) )
1852  else
1853  CHECK_ERROR( nc_sync(ncid) )
1854 
1855  return SUCCESS_CODE;
1856 }

Referenced by scale_file::file_flush().

Here is the caller graph for this function:

◆ file_close_c()

int32_t file_close_c ( const int32_t  fid,
const int32_t  abort 
)

Definition at line 1982 of file scale_file_netcdf.c.

1984 {
1985  int ncid;
1986  int i;
1987 
1988  if ( files[fid] == NULL ) return ALREADY_CLOSED_CODE;
1989  ncid = files[fid]->ncid;
1990 
1991  for (i=0; i<nvar; i++) {
1992  if ( vars[i] != NULL && vars[i]->ncid == ncid ) {
1993  free( vars[i]->start );
1994  free( vars[i]->count );
1995  free( vars[i] );
1996  vars[i] = NULL;
1997  }
1998  }
1999 
2000  for (i=0; i<nt; i++) {
2001  if ( tdims[i] != NULL && tdims[i]->ncid == ncid ) {
2002  free( tdims[i]->tval );
2003  free( tdims[i] );
2004  tdims[i] = NULL;
2005  }
2006  }
2007 
2008  if ( files[fid]->shared_mode ) {
2009  if ( ! abort ) CHECK_PNC_ERROR( ncmpi_close(ncid) )
2010  } else
2011  CHECK_ERROR( nc_close(ncid) )
2012 
2013  free( files[fid] );
2014  files[fid] = NULL;
2015 
2016  return SUCCESS_CODE;
2017 }

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
datainfo_t::time_end
real64_t time_end
Definition: scale_file.h:30
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:4585
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
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:35
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:36
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
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:32
tdim_t::varid
int varid
Definition: scale_file_netcdf.c:108
NTMAX
#define NTMAX
Definition: scale_file_netcdf.c:8
mod_atmos_vars::j
real(rp), allocatable, target, public j
Definition: mod_atmos_vars.F90:140
datainfo_t::dim_size
int32_t dim_size[RANK_MAX]
Definition: scale_file.h:27
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:34
datainfo_t::fid
int32_t fid
Definition: scale_file.h:37
datainfo_t::units
char units[File_HSHORT]
Definition: scale_file.h:22
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:44
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:182
datainfo_t::standard_name
char standard_name[File_HMID]
Definition: scale_file.h:23
datainfo_t::datatype
int32_t datatype
Definition: scale_file.h:24
tdim_t::ncid
int ncid
Definition: scale_file_netcdf.c:106
datainfo_t::description
char description[File_HMID]
Definition: scale_file.h:21
datainfo_t::time_units
char time_units[File_HMID]
Definition: scale_file.h:31
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:25
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:44
datainfo_t::natts
int32_t natts
Definition: scale_file.h:33
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:26
datainfo_t::time_start
real64_t time_start
Definition: scale_file.h:29
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:4613
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:28
datainfo_t::varname
char varname[File_HSHORT]
Definition: scale_file.h:20
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:38
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