NCEPLIBS-ncio  1.0.0
All Namespaces Files Functions
module_ncio.f90
Go to the documentation of this file.
1 
4 
9 module module_ncio
10 
11  use netcdf
12  use mpi
13 
14  implicit none
15  private
16 
17  type variable
18  integer varid
19  integer ndims
20  integer dtype
21  integer natts
22  integer deflate_level
23  logical shuffle
24  logical hasunlim
25  character(len=nf90_max_name) :: name
26  integer, allocatable, dimension(:) :: dimids
27  integer, allocatable, dimension(:) :: dimindxs
28  character(len=nf90_max_name), allocatable, dimension(:) :: dimnames
29  integer, allocatable, dimension(:) :: dimlens
30  integer, allocatable, dimension(:) :: chunksizes
31  end type variable
32 
33  type dimension
34  integer dimid
35  integer len
36  logical isunlimited
37  character(len=nf90_max_name) :: name
38  end type dimension
39 
40  type dataset
41  integer :: ncid
42  integer :: nvars
43  integer :: ndims
44  integer :: natts
45  integer :: nunlimdim
46  logical :: ishdf5
47  logical :: isparallel
48  character(len=500) filename
49  type(Variable), allocatable, dimension(:) :: variables
50  type(Dimension), allocatable, dimension(:) :: dimensions
51  end type dataset
52 
72  interface read_vardata
73  module procedure read_vardata_1d_r4, read_vardata_2d_r4, read_vardata_3d_r4,&
74  read_vardata_4d_r4, read_vardata_5d_r4, &
75  read_vardata_1d_r8, read_vardata_2d_r8, read_vardata_3d_r8,&
76  read_vardata_4d_r8, read_vardata_5d_r8, &
77  read_vardata_1d_int, read_vardata_2d_int, &
78  read_vardata_3d_int, read_vardata_4d_int, read_vardata_5d_int, &
79  read_vardata_1d_short, read_vardata_2d_short, &
80  read_vardata_3d_short, read_vardata_4d_short, read_vardata_5d_short , &
81  read_vardata_1d_byte, read_vardata_2d_byte, &
82  read_vardata_3d_byte, read_vardata_4d_byte, read_vardata_5d_byte, &
83  read_vardata_1d_char, read_vardata_2d_char, &
84  read_vardata_3d_char, read_vardata_4d_char, read_vardata_5d_char
85  end interface read_vardata
86 
87 
107  interface write_vardata
108  module procedure write_vardata_1d_r4, write_vardata_2d_r4, write_vardata_3d_r4,&
109  write_vardata_4d_r4, write_vardata_1d_r8, write_vardata_2d_r8, write_vardata_3d_r8,&
110  write_vardata_4d_r8, write_vardata_1d_int, write_vardata_2d_int, &
111  write_vardata_3d_int, write_vardata_4d_int, &
112  write_vardata_5d_int, write_vardata_5d_r4, write_vardata_5d_r8, &
113  write_vardata_1d_short, write_vardata_2d_short, write_vardata_3d_short, &
114  write_vardata_4d_short, write_vardata_5d_short, &
115  write_vardata_1d_byte, write_vardata_2d_byte, write_vardata_3d_byte, &
116  write_vardata_4d_byte, write_vardata_5d_byte, &
117  write_vardata_1d_char, write_vardata_2d_char, write_vardata_3d_char, &
118  write_vardata_4d_char, write_vardata_5d_char
119  end interface write_vardata
120 
126  interface read_attribute
127  module procedure read_attribute_r4_scalar, read_attribute_int_scalar,&
128  read_attribute_r8_scalar, read_attribute_r4_1d,&
129  read_attribute_int_1d, read_attribute_r8_1d, read_attribute_char, &
130  read_attribute_short_scalar, read_attribute_short_1d, &
131  read_attribute_byte_scalar, read_attribute_byte_1d
132  end interface read_attribute
133 
134  !! Write attribute 'attname' with data in 'values'. If optional
135  !! argument 'varname' is given, a variable attribute is written.
136  !! values can be a real(4), real(8), integer, string or 1d array.
137  !! @author jeff whitaker
138  interface write_attribute
139  module procedure write_attribute_r4_scalar, write_attribute_int_scalar,&
140  write_attribute_r8_scalar, write_attribute_r4_1d,&
141  write_attribute_int_1d, write_attribute_r8_1d, write_attribute_char, &
142  write_attribute_short_scalar, write_attribute_short_1d, &
143  write_attribute_byte_scalar, write_attribute_byte_1d
144  end interface write_attribute
145 
146  !! Quantize data.
147  !!
148  !! @author Jeff Whitaker, Cory Martin
149  interface quantize_data
150  module procedure quantize_data_2d, quantize_data_3d, &
151  quantize_data_4d, quantize_data_5d
152  end interface quantize_data
153  public :: open_dataset, create_dataset, close_dataset, dataset, variable, dimension, &
154  read_vardata, read_attribute, write_vardata, write_attribute, get_ndim, &
157 
158 contains
165  subroutine nccheck(status,halt,fname)
166  implicit none
167  integer, intent (in) :: status
168  logical, intent(in), optional :: halt
169  character(len=500), intent(in), optional :: fname
170  logical stopit
171  if (present(halt)) then
172  stopit = halt
173  else
174  stopit = .true.
175  endif
176  if (status /= nf90_noerr) then
177  write(0,*) status, trim(nf90_strerror(status))
178  if (present(fname)) then
179  write(0,*) trim(fname)
180  end if
181  if (stopit) stop 99
182  end if
183  end subroutine nccheck
189  function get_dim(dset, dimname) result(dim)
190  type(dataset) :: dset
191  type(dimension) :: dim
192  character(len=*), intent(in) :: dimname
193  integer ndim
194  ndim = get_ndim(dset, dimname)
195  dim = dset%dimensions(ndim)
196  end function get_dim
205  integer function get_ndim(dset, dimname)
206  type(dataset), intent(in) :: dset
207  character(len=*), intent(in) :: dimname
208  integer ndim
209  get_ndim = -1
210  do ndim=1,dset%ndims
211  if (trim(dset%dimensions(ndim)%name) == trim(dimname)) then
212  get_ndim = ndim
213  exit
214  endif
215  enddo
216  end function get_ndim
224  function get_var(dset, varname) result (var)
225  type(dataset) :: dset
226  type(variable) :: var
227  character(len=*) :: varname
228  integer nvar
229  nvar = get_nvar(dset, varname)
230  var = dset%variables(nvar)
231  end function get_var
238  logical function has_var(dset, varname)
239  type(dataset) :: dset
240  character(len=*) :: varname
241  integer nvar
242  nvar = get_nvar(dset, varname)
243  if (nvar > 0) then
244  has_var=.true.
245  else
246  has_var=.false.
247  endif
248  end function has_var
258  logical function has_attr(dset, attname, varname)
259  type(dataset) :: dset
260  character(len=*) :: attname
261  character(len=*), optional :: varname
262  integer nvar, varid, ncerr
263  nvar = get_nvar(dset, varname)
264  if(present(varname))then
265  nvar = get_nvar(dset,varname)
266  if (nvar < 0) then
267  has_attr = .false.
268  return
269  endif
270  varid = dset%variables(nvar)%varid
271  else
272  varid = nf90_global
273  endif
274  ncerr = nf90_inquire_attribute(dset%ncid, varid, attname)
275  if (ncerr /= 0) then
276  has_attr=.false.
277  else
278  has_attr=.true.
279  endif
280  end function has_attr
288  integer function get_nvar(dset,varname)
289  type(dataset), intent(in) :: dset
290  character(len=*), intent(in) :: varname
291  integer nvar
292  get_nvar = -1
293  do nvar=1,dset%nvars
294  if (trim(dset%variables(nvar)%name) == trim(varname)) then
295  get_nvar = nvar
296  exit
297  endif
298  enddo
299  end function get_nvar
307  subroutine set_varunlimdimlens_(dset,errcode)
308  type(dataset), intent(inout) :: dset
309  integer, intent(out), optional :: errcode
310  integer ndim,n,nvar,ncerr
311  logical return_errcode
312  if(present(errcode)) then
313  return_errcode=.true.
314  errcode = 0
315  else
316  return_errcode=.false.
317  endif
318  ! loop over all vars
319  do nvar=1,dset%nvars
320  ! does var have unlim dimension?
321  if (dset%variables(nvar)%hasunlim) then
322  ! loop over all var dimensions
323  do ndim=1,dset%variables(nvar)%ndims
324  n = dset%variables(nvar)%dimindxs(ndim)
325  ! n is the dimension index for this variable dimension
326  ! if this dim is unlimited, update dimlens entry
327  if (dset%dimensions(n)%isunlimited) then
328  ncerr = nf90_inquire_dimension(dset%ncid,&
329  dset%dimensions(n)%dimid, &
330  len=dset%variables(nvar)%dimlens(ndim))
331  if (return_errcode) then
332  call nccheck(ncerr,halt=.false.)
333  errcode=ncerr
334  return
335  else
336  call nccheck(ncerr)
337  endif
338  ! also update len attribute of Dimension object
339  dset%dimensions(n)%len = dset%variables(nvar)%dimlens(ndim)
340  endif
341  enddo
342  endif
343  enddo
344  end subroutine set_varunlimdimlens_
358  function open_dataset(filename,errcode,paropen, mpicomm) result(dset)
359  implicit none
360  character(len=*), intent(in) :: filename
361  type(dataset) :: dset
362  integer, intent(out), optional :: errcode
363  logical, intent(in), optional :: paropen
364  integer, intent(in), optional :: mpicomm
365  integer ncerr,nunlimdim,ndim,nvar,n,formatnum
366  logical return_errcode
367  if(present(errcode)) then
368  return_errcode=.true.
369  errcode = 0
370  else
371  return_errcode=.false.
372  endif
373  if (present(paropen)) then
374  if (paropen) then
375  dset%isparallel = .true.
376  else
377  dset%isparallel = .false.
378  end if
379  else
380  dset%isparallel = .false.
381  end if
382  ! open netcdf file, get info, populate Dataset object.
383  if (dset%isparallel) then
384  if (present(mpicomm)) then
385  ncerr = nf90_open(trim(filename), ior(nf90_nowrite, nf90_mpiio), &
386  comm=mpicomm, info = mpi_info_null, ncid=dset%ncid)
387  else
388  ncerr = nf90_open(trim(filename), ior(nf90_nowrite, nf90_mpiio), &
389  comm=mpi_comm_world, info = mpi_info_null, ncid=dset%ncid)
390  end if
391  else
392  ncerr = nf90_open(trim(filename), nf90_nowrite, ncid=dset%ncid)
393  end if
394  if (return_errcode) then
395  call nccheck(ncerr,halt=.false.,fname=filename)
396  errcode=ncerr
397  if (ncerr /= 0) return
398  else
399  call nccheck(ncerr,fname=filename)
400  endif
401  ncerr = nf90_inquire(dset%ncid, dset%ndims, dset%nvars, dset%natts, nunlimdim, formatnum=formatnum)
402  if (return_errcode) then
403  errcode=ncerr
404  call nccheck(ncerr,halt=.false.,fname=filename)
405  if (ncerr /= 0) return
406  else
407  call nccheck(ncerr,fname=filename)
408  endif
409  if (formatnum == nf90_format_netcdf4 .or. formatnum == nf90_format_netcdf4_classic) then
410  dset%ishdf5 = .true.
411  else
412  dset%ishdf5 = .false.
413  endif
414  dset%filename = trim(filename)
415  allocate(dset%variables(dset%nvars))
416  allocate(dset%dimensions(dset%ndims))
417  do ndim=1,dset%ndims
418  dset%dimensions(ndim)%dimid = ndim
419  ncerr = nf90_inquire_dimension(dset%ncid, ndim, name=dset%dimensions(ndim)%name, &
420  len=dset%dimensions(ndim)%len)
421  if (return_errcode) then
422  errcode=ncerr
423  call nccheck(ncerr,halt=.false.,fname=filename)
424  if (ncerr /= 0) return
425  else
426  call nccheck(ncerr,fname=filename)
427  endif
428  if (ndim == nunlimdim) then
429  dset%dimensions(ndim)%isunlimited = .true.
430  else
431  dset%dimensions(ndim)%isunlimited = .false.
432  endif
433  enddo
434  do nvar=1,dset%nvars
435  dset%variables(nvar)%hasunlim = .false.
436  dset%variables(nvar)%varid = nvar
437  ncerr = nf90_inquire_variable(dset%ncid, nvar,&
438  name=dset%variables(nvar)%name,&
439  natts=dset%variables(nvar)%natts,&
440  xtype=dset%variables(nvar)%dtype,&
441  ndims=dset%variables(nvar)%ndims)
442  if (return_errcode) then
443  errcode=ncerr
444  call nccheck(ncerr,halt=.false.,fname=filename)
445  if (ncerr /= 0) return
446  else
447  call nccheck(ncerr,fname=filename)
448  endif
449  allocate(dset%variables(nvar)%dimids(dset%variables(nvar)%ndims))
450  allocate(dset%variables(nvar)%dimindxs(dset%variables(nvar)%ndims))
451  allocate(dset%variables(nvar)%dimlens(dset%variables(nvar)%ndims))
452  allocate(dset%variables(nvar)%chunksizes(dset%variables(nvar)%ndims))
453  allocate(dset%variables(nvar)%dimnames(dset%variables(nvar)%ndims))
454  if (dset%ishdf5) then
455  ncerr = nf90_inquire_variable(dset%ncid, nvar,&
456  dimids=dset%variables(nvar)%dimids,&
457  deflate_level=dset%variables(nvar)%deflate_level,&
458  chunksizes=dset%variables(nvar)%chunksizes,&
459  shuffle=dset%variables(nvar)%shuffle)
460  else
461  ncerr = nf90_inquire_variable(dset%ncid, nvar,&
462  dimids=dset%variables(nvar)%dimids)
463  endif
464  if (return_errcode) then
465  errcode=ncerr
466  call nccheck(ncerr,halt=.false.,fname=filename)
467  if (ncerr /= 0) return
468  else
469  call nccheck(ncerr,fname=filename)
470  endif
471  do ndim=1,dset%variables(nvar)%ndims
472  do n=1,dset%ndims
473  if (dset%variables(nvar)%dimids(ndim) == dset%dimensions(n)%dimid) then
474  exit
475  endif
476  enddo
477  dset%variables(nvar)%dimindxs(ndim) = n
478  dset%variables(nvar)%dimlens(ndim) = dset%dimensions(n)%len
479  dset%variables(nvar)%dimnames(ndim) = dset%dimensions(n)%name
480  if (dset%dimensions(n)%isunlimited) then
481  dset%variables(nvar)%hasunlim = .true.
482  endif
483  enddo
484  enddo
485  end function open_dataset
486 !
506  function create_dataset(filename, dsetin, copy_vardata, paropen, nocompress, mpicomm, errcode) result(dset)
507  implicit none
508  character(len=*), intent(in) :: filename
509  character(len=nf90_max_name) :: attname, varname
510  logical, intent(in), optional :: copy_vardata
511  type(dataset) :: dset
512  type(dataset), intent(in) :: dsetin
513  logical, intent(in), optional :: paropen
514  integer, intent(in), optional :: mpicomm
515  logical, intent(in), optional :: nocompress
516  integer, intent(out), optional :: errcode
517  integer ncerr,ndim,nvar,n,ishuffle,natt
518  logical copyd, coordvar, compress
519  real(8), allocatable, dimension(:) :: values_1d
520  real(8), allocatable, dimension(:,:) :: values_2d
521  real(8), allocatable, dimension(:,:,:) :: values_3d
522  real(8), allocatable, dimension(:,:,:,:) :: values_4d
523  real(8), allocatable, dimension(:,:,:,:,:) :: values_5d
524  integer, allocatable, dimension(:) :: ivalues_1d
525  integer, allocatable, dimension(:,:) :: ivalues_2d
526  integer, allocatable, dimension(:,:,:) :: ivalues_3d
527  integer, allocatable, dimension(:,:,:,:) :: ivalues_4d
528  integer, allocatable, dimension(:,:,:,:,:) :: ivalues_5d
529  character, allocatable, dimension(:) :: cvalues_1d
530  character, allocatable, dimension(:,:) :: cvalues_2d
531  character, allocatable, dimension(:,:,:) :: cvalues_3d
532  character, allocatable, dimension(:,:,:,:) :: cvalues_4d
533  character, allocatable, dimension(:,:,:,:,:) :: cvalues_5d
534  logical return_errcode
535  if(present(errcode)) then
536  return_errcode=.true.
537  errcode = 0
538  else
539  return_errcode=.false.
540  endif
541  if (present(copy_vardata)) then
542  ! copy all variable data
543  copyd = .true.
544  else
545  ! only copy coordinate variable data
546  copyd = .false.
547  endif
548  if (present(paropen)) then
549  if (paropen) then
550  dset%isparallel = .true.
551  else
552  dset%isparallel = .false.
553  end if
554  else
555  dset%isparallel = .false.
556  end if
557  compress = .true.
558  if (present(nocompress)) then
559  if (nocompress) then
560  compress = .false.
561  end if
562  end if
563  ! create netcdf file
564  if (dsetin%ishdf5) then
565  if (dset%isparallel) then
566  if (present(mpicomm)) then
567  ncerr = nf90_create(trim(filename), &
568  cmode=ior(nf90_clobber,nf90_netcdf4), ncid=dset%ncid, &
569  comm = mpicomm, info = mpi_info_null)
570  else
571  ncerr = nf90_create(trim(filename), &
572  cmode=ior(nf90_clobber,nf90_netcdf4), ncid=dset%ncid, &
573  comm = mpi_comm_world, info = mpi_info_null)
574  end if
575  else
576  ncerr = nf90_create(trim(filename), &
577  cmode=ior(nf90_clobber,nf90_netcdf4), &
578  ncid=dset%ncid)
579  end if
580  dset%ishdf5 = .true.
581  else
582  ncerr = nf90_create(trim(filename), &
583  cmode=ior(nf90_clobber,nf90_cdf5), &
584  ncid=dset%ncid)
585  dset%ishdf5 = .false.
586  endif
587  if (return_errcode) then
588  errcode=ncerr
589  call nccheck(ncerr,halt=.false.,fname=filename)
590  if (ncerr /= 0) return
591  else
592  call nccheck(ncerr,fname=filename)
593  endif
594  ! copy global attributes
595  do natt=1,dsetin%natts
596  ncerr = nf90_inq_attname(dsetin%ncid, nf90_global, natt, attname)
597  if (return_errcode) then
598  errcode=ncerr
599  call nccheck(ncerr,halt=.false.)
600  if (ncerr /= 0) return
601  else
602  call nccheck(ncerr)
603  endif
604  ncerr = nf90_copy_att(dsetin%ncid, nf90_global, attname, dset%ncid, nf90_global)
605  if (return_errcode) then
606  errcode=ncerr
607  call nccheck(ncerr,halt=.false.)
608  if (ncerr /= 0) return
609  else
610  call nccheck(ncerr)
611  endif
612  enddo
613  dset%natts = dsetin%natts
614  dset%filename = trim(filename)
615  dset%ndims = dsetin%ndims
616  dset%nvars = dsetin%nvars
617  allocate(dset%variables(dsetin%nvars))
618  allocate(dset%dimensions(dsetin%ndims))
619  ! create dimensions
620  do ndim=1,dsetin%ndims
621  if (dsetin%dimensions(ndim)%isunlimited) then
622  ncerr = nf90_def_dim(dset%ncid, trim(dsetin%dimensions(ndim)%name), &
623  nf90_unlimited, &
624  dset%dimensions(ndim)%dimid)
625  if (return_errcode) then
626  errcode=ncerr
627  call nccheck(ncerr,halt=.false.)
628  if (ncerr /= 0) return
629  else
630  call nccheck(ncerr)
631  endif
632  dset%dimensions(ndim)%isunlimited = .true.
633  dset%nunlimdim = ndim
634  dset%dimensions(ndim)%len = 0
635  dset%dimensions(ndim)%name = trim(dsetin%dimensions(ndim)%name)
636  else
637  ncerr = nf90_def_dim(dset%ncid, trim(dsetin%dimensions(ndim)%name),&
638  dsetin%dimensions(ndim)%len, &
639  dset%dimensions(ndim)%dimid)
640  if (return_errcode) then
641  errcode=ncerr
642  call nccheck(ncerr,halt=.false.)
643  if (ncerr /= 0) return
644  else
645  call nccheck(ncerr)
646  endif
647  dset%dimensions(ndim)%len = dsetin%dimensions(ndim)%len
648  dset%dimensions(ndim)%isunlimited = .false.
649  dset%dimensions(ndim)%name = trim(dsetin%dimensions(ndim)%name)
650  endif
651  enddo
652  ! create variables
653  do nvar=1,dsetin%nvars
654  dset%variables(nvar)%hasunlim = .false.
655  dset%variables(nvar)%ndims = dsetin%variables(nvar)%ndims
656  allocate(dset%variables(nvar)%dimids(dset%variables(nvar)%ndims))
657  allocate(dset%variables(nvar)%dimindxs(dset%variables(nvar)%ndims))
658  allocate(dset%variables(nvar)%dimnames(dset%variables(nvar)%ndims))
659  allocate(dset%variables(nvar)%dimlens(dset%variables(nvar)%ndims))
660  allocate(dset%variables(nvar)%chunksizes(dset%variables(nvar)%ndims))
661  dset%variables(nvar)%chunksizes = dsetin%variables(nvar)%chunksizes
662  do ndim=1,dset%variables(nvar)%ndims
663  do n=1,dset%ndims
664  if (trim(dsetin%variables(nvar)%dimnames(ndim)) == &
665  trim(dset%dimensions(n)%name)) then
666  exit
667  endif
668  enddo
669  dset%variables(nvar)%dimindxs(ndim) = n
670  dset%variables(nvar)%dimids(ndim) = dset%dimensions(n)%dimid
671  dset%variables(nvar)%dimlens(ndim) = dset%dimensions(n)%len
672  dset%variables(nvar)%dimnames(ndim) = dset%dimensions(n)%name
673  if (dset%dimensions(n)%isunlimited) then
674  dset%variables(nvar)%hasunlim = .true.
675  endif
676  enddo
677  dset%variables(nvar)%name = dsetin%variables(nvar)%name
678  dset%variables(nvar)%dtype = dsetin%variables(nvar)%dtype
679  if (maxval(dset%variables(nvar)%chunksizes) > 0 .and. dset%ishdf5) then
680  ! workaround for older versions of netcdf-fortran that don't
681  ! like zero chunksize to be specified.
682  ncerr = nf90_def_var(dset%ncid, &
683  trim(dset%variables(nvar)%name),&
684  dset%variables(nvar)%dtype, &
685  dset%variables(nvar)%dimids, &
686  dset%variables(nvar)%varid, &
687  chunksizes=dset%variables(nvar)%chunksizes)
688  else
689  ncerr = nf90_def_var(dset%ncid, &
690  trim(dset%variables(nvar)%name),&
691  dset%variables(nvar)%dtype, &
692  dset%variables(nvar)%dimids, &
693  dset%variables(nvar)%varid)
694  endif
695  if (return_errcode) then
696  errcode=ncerr
697  call nccheck(ncerr,halt=.false.)
698  if (ncerr /= 0) return
699  else
700  call nccheck(ncerr)
701  endif
702  if (dsetin%variables(nvar)%deflate_level > 0 .and. dset%ishdf5 .and. compress) then
703  if (dsetin%variables(nvar)%shuffle) then
704  ishuffle=1
705  else
706  ishuffle=0
707  endif
708  ncerr = nf90_def_var_deflate(dset%ncid, dset%variables(nvar)%varid,&
709  ishuffle,1,dsetin%variables(nvar)%deflate_level)
710  if (return_errcode) then
711  errcode=ncerr
712  call nccheck(ncerr,halt=.false.)
713  if (ncerr /= 0) return
714  else
715  call nccheck(ncerr)
716  endif
717  dset%variables(nvar)%shuffle = dsetin%variables(nvar)%shuffle
718  dset%variables(nvar)%deflate_level = &
719  dsetin%variables(nvar)%deflate_level
720  endif
721  ! copy variable attributes
722  do natt=1,dsetin%variables(nvar)%natts
723  ncerr = nf90_inq_attname(dsetin%ncid, dsetin%variables(nvar)%varid, natt, attname)
724  if (return_errcode) then
725  errcode=ncerr
726  call nccheck(ncerr,halt=.false.)
727  if (ncerr /= 0) return
728  else
729  call nccheck(ncerr)
730  endif
731  if (.not. compress) then
732  if (trim(attname) == 'max_abs_compression_error' &
733  .or. trim(attname) == 'nbits') then
734  cycle
735  end if
736  end if
737  ncerr = nf90_copy_att(dsetin%ncid, dsetin%variables(nvar)%varid, attname, dset%ncid, dset%variables(nvar)%varid)
738  if (return_errcode) then
739  errcode=ncerr
740  call nccheck(ncerr,halt=.false.)
741  if (ncerr /= 0) return
742  else
743  call nccheck(ncerr)
744  endif
745  enddo
746  enddo
747  ncerr = nf90_enddef(dset%ncid)
748  if (return_errcode) then
749  errcode=ncerr
750  call nccheck(ncerr,halt=.false.)
751  if (ncerr /= 0) return
752  else
753  call nccheck(ncerr)
754  endif
755  ! copy variable data
756  ! assumes data is real (32 or 64 bit), or integer (16 or 32 bit) and 1-4d.
757  do nvar=1,dsetin%nvars
758  varname = trim(dsetin%variables(nvar)%name)
759  ! is this variable a coordinate variable?
760  coordvar = .false.
761  if (trim(varname) == 'lats' .or. trim(varname) == 'lons' .or. &
762  trim(varname) == 'lat' .or. trim(varname) == 'lon') then
763  coordvar = .true.
764  else
765  do ndim=1,dset%ndims
766  if (trim(varname) == trim(dset%dimensions(ndim)%name)) then
767  coordvar = .true.
768  endif
769  enddo
770  endif
771  ! if copy_data flag not given, and not a coordinate var,
772  ! skip to next var.
773  if (.not. coordvar .and. .not. copyd) cycle
774  ! real variable
775  if (dsetin%variables(nvar)%dtype == nf90_float .or.&
776  dsetin%variables(nvar)%dtype == nf90_double) then
777  if (dsetin%variables(nvar)%ndims == 1) then
778  call read_vardata(dsetin, varname, values_1d)
779  call write_vardata(dset, varname, values_1d)
780  else if (dsetin%variables(nvar)%ndims == 2) then
781  call read_vardata(dsetin, varname, values_2d)
782  call write_vardata(dset, varname, values_2d)
783  else if (dsetin%variables(nvar)%ndims == 3) then
784  call read_vardata(dsetin, varname, values_3d)
785  call write_vardata(dset, varname, values_3d)
786  else if (dsetin%variables(nvar)%ndims == 4) then
787  call read_vardata(dsetin, varname, values_4d)
788  call write_vardata(dset, varname, values_4d)
789  else if (dsetin%variables(nvar)%ndims == 5) then
790  call read_vardata(dsetin, varname, values_5d)
791  call write_vardata(dset, varname, values_5d)
792  endif
793  ! integer var
794  elseif (dsetin%variables(nvar)%dtype == nf90_int .or.&
795  dsetin%variables(nvar)%dtype == nf90_byte .or.&
796  dsetin%variables(nvar)%dtype == nf90_short) then
797  ! TODO: support NF90_UBYTE, USHORT, UINT, INT64, UINT64
798  if (dsetin%variables(nvar)%ndims == 1) then
799  call read_vardata(dsetin, varname, ivalues_1d)
800  call write_vardata(dset, varname, ivalues_1d)
801  else if (dsetin%variables(nvar)%ndims == 2) then
802  call read_vardata(dsetin, varname, ivalues_2d)
803  call write_vardata(dset, varname, ivalues_2d)
804  else if (dsetin%variables(nvar)%ndims == 3) then
805  call read_vardata(dsetin, varname, ivalues_3d)
806  call write_vardata(dset, varname, ivalues_3d)
807  else if (dsetin%variables(nvar)%ndims == 4) then
808  call read_vardata(dsetin, varname, ivalues_4d)
809  call write_vardata(dset, varname, ivalues_4d)
810  else if (dsetin%variables(nvar)%ndims == 5) then
811  call read_vardata(dsetin, varname, ivalues_5d)
812  call write_vardata(dset, varname, ivalues_5d)
813  endif
814  elseif (dsetin%variables(nvar)%dtype == nf90_char) then
815  if (dsetin%variables(nvar)%ndims == 1) then
816  call read_vardata(dsetin, varname, cvalues_1d)
817  call write_vardata(dset, varname, cvalues_1d)
818  else if (dsetin%variables(nvar)%ndims == 2) then
819  call read_vardata(dsetin, varname, cvalues_2d)
820  call write_vardata(dset, varname, cvalues_2d)
821  else if (dsetin%variables(nvar)%ndims == 3) then
822  call read_vardata(dsetin, varname, cvalues_3d)
823  call write_vardata(dset, varname, cvalues_3d)
824  else if (dsetin%variables(nvar)%ndims == 4) then
825  call read_vardata(dsetin, varname, cvalues_4d)
826  call write_vardata(dset, varname, cvalues_4d)
827  else if (dsetin%variables(nvar)%ndims == 5) then
828  call read_vardata(dsetin, varname, cvalues_5d)
829  call write_vardata(dset, varname, cvalues_5d)
830  endif
831  else
832  print *,'not copying variable ',trim(adjustl(varname)),&
833  ' (unsupported data type or rank)'
834  endif
835  enddo
836  end function create_dataset
845  subroutine close_dataset(dset,errcode)
846  type(dataset), intent(inout) :: dset
847  integer, intent(out), optional :: errcode
848  integer ncerr, nvar
849  logical return_errcode
850  if(present(errcode)) then
851  return_errcode=.true.
852  errcode = 0
853  else
854  return_errcode=.false.
855  endif
856  ncerr = nf90_close(ncid=dset%ncid)
857  if (return_errcode) then
858  errcode=ncerr
859  call nccheck(ncerr,halt=.false.)
860  if (ncerr /= 0) return
861  else
862  call nccheck(ncerr)
863  endif
864  do nvar=1,dset%nvars
865  deallocate(dset%variables(nvar)%dimids)
866  deallocate(dset%variables(nvar)%dimindxs)
867  deallocate(dset%variables(nvar)%dimlens)
868  deallocate(dset%variables(nvar)%chunksizes)
869  deallocate(dset%variables(nvar)%dimnames)
870  enddo
871  deallocate(dset%variables,dset%dimensions)
872  end subroutine close_dataset
873 
874  subroutine read_vardata_1d_r4(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
875  real(4), allocatable, dimension(:), intent(inout) :: values
876  include "read_vardata_code_1d.f90"
877  end subroutine read_vardata_1d_r4
878 
879  subroutine read_vardata_2d_r4(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
880  real(4), allocatable, dimension(:,:), intent(inout) :: values
881  include "read_vardata_code_2d.f90"
882  end subroutine read_vardata_2d_r4
883 
884  subroutine read_vardata_3d_r4(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
885  real(4), allocatable, dimension(:,:,:), intent(inout) :: values
886  include "read_vardata_code_3d.f90"
887  end subroutine read_vardata_3d_r4
888 
889  subroutine read_vardata_4d_r4(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
890  real(4), allocatable, dimension(:,:,:,:), intent(inout) :: values
891  include "read_vardata_code_4d.f90"
892  end subroutine read_vardata_4d_r4
893 
894  subroutine read_vardata_5d_r4(dset, varname, values, errcode)
895  real(4), allocatable, dimension(:,:,:,:,:), intent(inout) :: values
896  include "read_vardata_code_5d.f90"
897  end subroutine read_vardata_5d_r4
898 
899  subroutine read_vardata_1d_r8(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
900  real(8), allocatable, dimension(:), intent(inout) :: values
901  include "read_vardata_code_1d.f90"
902  end subroutine read_vardata_1d_r8
903 
904  subroutine read_vardata_2d_r8(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
905  real(8), allocatable, dimension(:,:), intent(inout) :: values
906  include "read_vardata_code_2d.f90"
907  end subroutine read_vardata_2d_r8
908 
909  subroutine read_vardata_3d_r8(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
910  real(8), allocatable, dimension(:,:,:), intent(inout) :: values
911  include "read_vardata_code_3d.f90"
912  end subroutine read_vardata_3d_r8
913 
914  subroutine read_vardata_4d_r8(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
915  real(8), allocatable, dimension(:,:,:,:), intent(inout) :: values
916  include "read_vardata_code_4d.f90"
917  end subroutine read_vardata_4d_r8
918 
919  subroutine read_vardata_5d_r8(dset, varname, values, errcode)
920  real(8), allocatable, dimension(:,:,:,:,:), intent(inout) :: values
921  include "read_vardata_code_5d.f90"
922  end subroutine read_vardata_5d_r8
923 
924  subroutine read_vardata_1d_int(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
925  integer, allocatable, dimension(:), intent(inout) :: values
926  include "read_vardata_code_1d.f90"
927  end subroutine read_vardata_1d_int
928 
929  subroutine read_vardata_2d_int(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
930  integer, allocatable, dimension(:,:), intent(inout) :: values
931  include "read_vardata_code_2d.f90"
932  end subroutine read_vardata_2d_int
933 
934  subroutine read_vardata_3d_int(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
935  integer, allocatable, dimension(:,:,:), intent(inout) :: values
936  include "read_vardata_code_3d.f90"
937  end subroutine read_vardata_3d_int
938 
939  subroutine read_vardata_4d_int(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
940  integer, allocatable, dimension(:,:,:,:), intent(inout) :: values
941  include "read_vardata_code_4d.f90"
942  end subroutine read_vardata_4d_int
943 
944  subroutine read_vardata_5d_int(dset, varname, values, errcode)
945  integer, allocatable, dimension(:,:,:,:,:), intent(inout) :: values
946  include "read_vardata_code_5d.f90"
947  end subroutine read_vardata_5d_int
948 
949  subroutine read_vardata_1d_short(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
950  integer(2), allocatable, dimension(:), intent(inout) :: values
951  include "read_vardata_code_1d.f90"
952  end subroutine read_vardata_1d_short
953 
954  subroutine read_vardata_2d_short(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
955  integer(2), allocatable, dimension(:,:), intent(inout) :: values
956  include "read_vardata_code_2d.f90"
957  end subroutine read_vardata_2d_short
958 
959  subroutine read_vardata_3d_short(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
960  integer(2), allocatable, dimension(:,:,:), intent(inout) :: values
961  include "read_vardata_code_3d.f90"
962  end subroutine read_vardata_3d_short
963 
964  subroutine read_vardata_4d_short(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
965  integer(2), allocatable, dimension(:,:,:,:), intent(inout) :: values
966  include "read_vardata_code_4d.f90"
967  end subroutine read_vardata_4d_short
968 
969  subroutine read_vardata_5d_short(dset, varname, values, errcode)
970  integer(2), allocatable, dimension(:,:,:,:,:), intent(inout) :: values
971  include "read_vardata_code_5d.f90"
972  end subroutine read_vardata_5d_short
973 
974  subroutine read_vardata_1d_byte(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
975  integer(1), allocatable, dimension(:), intent(inout) :: values
976  include "read_vardata_code_1d.f90"
977  end subroutine read_vardata_1d_byte
978 
979  subroutine read_vardata_2d_byte(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
980  integer(1), allocatable, dimension(:,:), intent(inout) :: values
981  include "read_vardata_code_2d.f90"
982  end subroutine read_vardata_2d_byte
983 
984  subroutine read_vardata_3d_byte(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
985  integer(1), allocatable, dimension(:,:,:), intent(inout) :: values
986  include "read_vardata_code_3d.f90"
987  end subroutine read_vardata_3d_byte
988 
989  subroutine read_vardata_4d_byte(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
990  integer(1), allocatable, dimension(:,:,:,:), intent(inout) :: values
991  include "read_vardata_code_4d.f90"
992  end subroutine read_vardata_4d_byte
993 
994  subroutine read_vardata_5d_byte(dset, varname, values, errcode)
995  integer(1), allocatable, dimension(:,:,:,:,:), intent(inout) :: values
996  include "read_vardata_code_5d.f90"
997  end subroutine read_vardata_5d_byte
998 
999  subroutine read_vardata_1d_char(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1000  character, allocatable, dimension(:), intent(inout) :: values
1001  include "read_vardata_code_1d.f90"
1002  end subroutine read_vardata_1d_char
1003 
1004  subroutine read_vardata_2d_char(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1005  character, allocatable, dimension(:,:), intent(inout) :: values
1006  include "read_vardata_code_2d.f90"
1007  end subroutine read_vardata_2d_char
1008 
1009  subroutine read_vardata_3d_char(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1010  character, allocatable, dimension(:,:,:), intent(inout) :: values
1011  include "read_vardata_code_3d.f90"
1012  end subroutine read_vardata_3d_char
1013 
1014  subroutine read_vardata_4d_char(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1015  character, allocatable, dimension(:,:,:,:), intent(inout) :: values
1016  include "read_vardata_code_4d.f90"
1017  end subroutine read_vardata_4d_char
1018 
1019  subroutine read_vardata_5d_char(dset, varname, values, errcode)
1020  character, allocatable, dimension(:,:,:,:,:), intent(inout) :: values
1021  include "read_vardata_code_5d.f90"
1022  end subroutine read_vardata_5d_char
1023 
1024  subroutine write_vardata_1d_r4(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1025  real(4), dimension(:), intent(in) :: values
1026  integer, intent(in), optional :: ncstart(1)
1027  integer, intent(in), optional :: nccount(1)
1028  include "write_vardata_code.f90"
1029  end subroutine write_vardata_1d_r4
1030 
1031  subroutine write_vardata_2d_r4(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1032  real(4), dimension(:,:), intent(in) :: values
1033  integer, intent(in), optional :: ncstart(2)
1034  integer, intent(in), optional :: nccount(2)
1035  include "write_vardata_code.f90"
1036  end subroutine write_vardata_2d_r4
1037 
1038  subroutine write_vardata_3d_r4(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1039  real(4), dimension(:,:,:), intent(in) :: values
1040  integer, intent(in), optional :: ncstart(3)
1041  integer, intent(in), optional :: nccount(3)
1042  include "write_vardata_code.f90"
1043  end subroutine write_vardata_3d_r4
1044 
1045  subroutine write_vardata_4d_r4(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1046  real(4), dimension(:,:,:,:), intent(in) :: values
1047  integer, intent(in), optional :: ncstart(4)
1048  integer, intent(in), optional :: nccount(4)
1049  include "write_vardata_code.f90"
1050  end subroutine write_vardata_4d_r4
1051 
1052  subroutine write_vardata_5d_r4(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1053  real(4), dimension(:,:,:,:,:), intent(in) :: values
1054  integer, intent(in), optional :: ncstart(5)
1055  integer, intent(in), optional :: nccount(5)
1056  include "write_vardata_code.f90"
1057  end subroutine write_vardata_5d_r4
1058 
1059  subroutine write_vardata_1d_r8(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1060  real(8), dimension(:), intent(in) :: values
1061  integer, intent(in), optional :: ncstart(1)
1062  integer, intent(in), optional :: nccount(1)
1063  include "write_vardata_code.f90"
1064  end subroutine write_vardata_1d_r8
1065 
1066  subroutine write_vardata_2d_r8(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1067  real(8), dimension(:,:), intent(in) :: values
1068  integer, intent(in), optional :: ncstart(2)
1069  integer, intent(in), optional :: nccount(2)
1070  include "write_vardata_code.f90"
1071  end subroutine write_vardata_2d_r8
1072 
1073  subroutine write_vardata_3d_r8(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1074  real(8), dimension(:,:,:), intent(in) :: values
1075  integer, intent(in), optional :: ncstart(3)
1076  integer, intent(in), optional :: nccount(3)
1077  include "write_vardata_code.f90"
1078  end subroutine write_vardata_3d_r8
1079 
1080  subroutine write_vardata_4d_r8(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1081  real(8), dimension(:,:,:,:), intent(in) :: values
1082  integer, intent(in), optional :: ncstart(4)
1083  integer, intent(in), optional :: nccount(4)
1084  include "write_vardata_code.f90"
1085  end subroutine write_vardata_4d_r8
1086 
1087  subroutine write_vardata_5d_r8(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1088  real(8), dimension(:,:,:,:,:), intent(in) :: values
1089  integer, intent(in), optional :: ncstart(5)
1090  integer, intent(in), optional :: nccount(5)
1091  include "write_vardata_code.f90"
1092  end subroutine write_vardata_5d_r8
1093 
1094  subroutine write_vardata_1d_int(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1095  integer, dimension(:), intent(in) :: values
1096  integer, intent(in), optional :: ncstart(1)
1097  integer, intent(in), optional :: nccount(1)
1098  include "write_vardata_code.f90"
1099  end subroutine write_vardata_1d_int
1100 
1101  subroutine write_vardata_2d_int(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1102  integer, dimension(:,:), intent(in) :: values
1103  integer, intent(in), optional :: ncstart(2)
1104  integer, intent(in), optional :: nccount(2)
1105  include "write_vardata_code.f90"
1106  end subroutine write_vardata_2d_int
1107 
1108  subroutine write_vardata_3d_int(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1109  integer, dimension(:,:,:), intent(in) :: values
1110  integer, intent(in), optional :: ncstart(3)
1111  integer, intent(in), optional :: nccount(3)
1112  include "write_vardata_code.f90"
1113  end subroutine write_vardata_3d_int
1114 
1115  subroutine write_vardata_4d_int(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1116  integer, dimension(:,:,:,:), intent(in) :: values
1117  integer, intent(in), optional :: ncstart(4)
1118  integer, intent(in), optional :: nccount(4)
1119  include "write_vardata_code.f90"
1120  end subroutine write_vardata_4d_int
1121 
1122  subroutine write_vardata_5d_int(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1123  integer, dimension(:,:,:,:,:), intent(in) :: values
1124  integer, intent(in), optional :: ncstart(5)
1125  integer, intent(in), optional :: nccount(5)
1126  include "write_vardata_code.f90"
1127  end subroutine write_vardata_5d_int
1128 
1129  subroutine write_vardata_1d_short(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1130  integer(2), dimension(:), intent(in) :: values
1131  integer, intent(in), optional :: ncstart(1)
1132  integer, intent(in), optional :: nccount(1)
1133  include "write_vardata_code.f90"
1134  end subroutine write_vardata_1d_short
1135 
1136  subroutine write_vardata_2d_short(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1137  integer(2), dimension(:,:), intent(in) :: values
1138  integer, intent(in), optional :: ncstart(2)
1139  integer, intent(in), optional :: nccount(2)
1140  include "write_vardata_code.f90"
1141  end subroutine write_vardata_2d_short
1142 
1143  subroutine write_vardata_3d_short(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1144  integer(2), dimension(:,:,:), intent(in) :: values
1145  integer, intent(in), optional :: ncstart(3)
1146  integer, intent(in), optional :: nccount(3)
1147  include "write_vardata_code.f90"
1148  end subroutine write_vardata_3d_short
1149 
1150  subroutine write_vardata_4d_short(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1151  integer(2), dimension(:,:,:,:), intent(in) :: values
1152  integer, intent(in), optional :: ncstart(4)
1153  integer, intent(in), optional :: nccount(4)
1154  include "write_vardata_code.f90"
1155  end subroutine write_vardata_4d_short
1156 
1157  subroutine write_vardata_5d_short(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1158  integer(2), dimension(:,:,:,:,:), intent(in) :: values
1159  integer, intent(in), optional :: ncstart(5)
1160  integer, intent(in), optional :: nccount(5)
1161  include "write_vardata_code.f90"
1162  end subroutine write_vardata_5d_short
1163 
1164  subroutine write_vardata_1d_byte(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1165  integer(1), dimension(:), intent(in) :: values
1166  integer, intent(in), optional :: ncstart(1)
1167  integer, intent(in), optional :: nccount(1)
1168  include "write_vardata_code.f90"
1169  end subroutine write_vardata_1d_byte
1170 
1171  subroutine write_vardata_2d_byte(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1172  integer(1), dimension(:,:), intent(in) :: values
1173  integer, intent(in), optional :: ncstart(2)
1174  integer, intent(in), optional :: nccount(2)
1175  include "write_vardata_code.f90"
1176  end subroutine write_vardata_2d_byte
1177 
1178  subroutine write_vardata_3d_byte(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1179  integer(1), dimension(:,:,:), intent(in) :: values
1180  integer, intent(in), optional :: ncstart(3)
1181  integer, intent(in), optional :: nccount(3)
1182  include "write_vardata_code.f90"
1183  end subroutine write_vardata_3d_byte
1184 
1185  subroutine write_vardata_4d_byte(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1186  integer(1), dimension(:,:,:,:), intent(in) :: values
1187  integer, intent(in), optional :: ncstart(4)
1188  integer, intent(in), optional :: nccount(4)
1189  include "write_vardata_code.f90"
1190  end subroutine write_vardata_4d_byte
1191 
1192  subroutine write_vardata_5d_byte(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1193  integer(1), dimension(:,:,:,:,:), intent(in) :: values
1194  integer, intent(in), optional :: ncstart(5)
1195  integer, intent(in), optional :: nccount(5)
1196  include "write_vardata_code.f90"
1197  end subroutine write_vardata_5d_byte
1198 
1199  subroutine write_vardata_1d_char(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1200  character, dimension(:), intent(in) :: values
1201  integer, intent(in), optional :: ncstart(1)
1202  integer, intent(in), optional :: nccount(1)
1203  include "write_vardata_code.f90"
1204  end subroutine write_vardata_1d_char
1205 
1206  subroutine write_vardata_2d_char(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1207  character, dimension(:,:), intent(in) :: values
1208  integer, intent(in), optional :: ncstart(2)
1209  integer, intent(in), optional :: nccount(2)
1210  include "write_vardata_code.f90"
1211  end subroutine write_vardata_2d_char
1212 
1213  subroutine write_vardata_3d_char(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1214  character, dimension(:,:,:), intent(in) :: values
1215  integer, intent(in), optional :: ncstart(3)
1216  integer, intent(in), optional :: nccount(3)
1217  include "write_vardata_code.f90"
1218  end subroutine write_vardata_3d_char
1219 
1220  subroutine write_vardata_4d_char(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1221  character, dimension(:,:,:,:), intent(in) :: values
1222  integer, intent(in), optional :: ncstart(4)
1223  integer, intent(in), optional :: nccount(4)
1224  include "write_vardata_code.f90"
1225  end subroutine write_vardata_4d_char
1226 
1227  subroutine write_vardata_5d_char(dset, varname, values, nslice, slicedim, ncstart, nccount, errcode)
1228  character, dimension(:,:,:,:,:), intent(in) :: values
1229  integer, intent(in), optional :: ncstart(5)
1230  integer, intent(in), optional :: nccount(5)
1231  include "write_vardata_code.f90"
1232  end subroutine write_vardata_5d_char
1233 
1234  subroutine read_attribute_int_scalar(dset, attname, values, varname, errcode)
1235  integer, intent(inout) :: values
1236  include "read_scalar_attribute_code.f90"
1237  end subroutine read_attribute_int_scalar
1238 
1239  subroutine read_attribute_short_scalar(dset, attname, values, varname, errcode)
1240  integer(2), intent(inout) :: values
1241  include "read_scalar_attribute_code.f90"
1242  end subroutine read_attribute_short_scalar
1243 
1244  subroutine read_attribute_byte_scalar(dset, attname, values, varname, errcode)
1245  integer(1), intent(inout) :: values
1246  include "read_scalar_attribute_code.f90"
1247  end subroutine read_attribute_byte_scalar
1248 
1249  subroutine read_attribute_r4_scalar(dset, attname, values, varname, errcode)
1250  real(4), intent(inout) :: values
1251  include "read_scalar_attribute_code.f90"
1252  end subroutine read_attribute_r4_scalar
1253 
1254  subroutine read_attribute_r8_scalar(dset, attname, values, varname, errcode)
1255  real(8), intent(inout) :: values
1256  include "read_scalar_attribute_code.f90"
1257  end subroutine read_attribute_r8_scalar
1258 
1259  subroutine read_attribute_r4_1d(dset, attname, values, varname, errcode)
1260  real(4), intent(inout), allocatable, dimension(:) :: values
1261  include "read_attribute_code.f90"
1262  end subroutine read_attribute_r4_1d
1263 
1264  subroutine read_attribute_r8_1d(dset, attname, values, varname, errcode)
1265  real(8), intent(inout), allocatable, dimension(:) :: values
1266  include "read_attribute_code.f90"
1267  end subroutine read_attribute_r8_1d
1268 
1269  subroutine read_attribute_int_1d(dset, attname, values, varname, errcode)
1270  integer, intent(inout), allocatable, dimension(:) :: values
1271  include "read_attribute_code.f90"
1272  end subroutine read_attribute_int_1d
1273 
1274  subroutine read_attribute_short_1d(dset, attname, values, varname, errcode)
1275  integer(2), intent(inout), allocatable, dimension(:) :: values
1276  include "read_attribute_code.f90"
1277  end subroutine read_attribute_short_1d
1278 
1279  subroutine read_attribute_byte_1d(dset, attname, values, varname, errcode)
1280  integer(1), intent(inout), allocatable, dimension(:) :: values
1281  include "read_attribute_code.f90"
1282  end subroutine read_attribute_byte_1d
1283 
1284  subroutine read_attribute_char(dset, attname, values, varname, errcode)
1285  character(len=*), intent(inout) :: values
1286  include "read_scalar_attribute_code.f90"
1287  end subroutine read_attribute_char
1288 
1289  subroutine write_attribute_int_scalar(dset, attname, values, varname, errcode)
1290  integer, intent(in) :: values
1291  include "write_attribute_code.f90"
1292  end subroutine write_attribute_int_scalar
1293 
1294  subroutine write_attribute_short_scalar(dset, attname, values, varname, errcode)
1295  integer(2), intent(in) :: values
1296  include "write_attribute_code.f90"
1297  end subroutine write_attribute_short_scalar
1298 
1299  subroutine write_attribute_byte_scalar(dset, attname, values, varname, errcode)
1300  integer(1), intent(in) :: values
1301  include "write_attribute_code.f90"
1302  end subroutine write_attribute_byte_scalar
1303 
1304  subroutine write_attribute_r4_scalar(dset, attname, values, varname, errcode)
1305  real(4), intent(in) :: values
1306  include "write_attribute_code.f90"
1307  end subroutine write_attribute_r4_scalar
1308 
1309  subroutine write_attribute_r8_scalar(dset, attname, values, varname, errcode)
1310  real(8), intent(in) :: values
1311  include "write_attribute_code.f90"
1312  end subroutine write_attribute_r8_scalar
1313 
1314  subroutine write_attribute_r4_1d(dset, attname, values, varname, errcode)
1315  real(4), intent(in), allocatable, dimension(:) :: values
1316  include "write_attribute_code.f90"
1317  end subroutine write_attribute_r4_1d
1318 
1319  subroutine write_attribute_r8_1d(dset, attname, values, varname, errcode)
1320  real(8), intent(in), allocatable, dimension(:) :: values
1321  include "write_attribute_code.f90"
1322  end subroutine write_attribute_r8_1d
1323 
1324  subroutine write_attribute_int_1d(dset, attname, values, varname, errcode)
1325  integer, intent(in), allocatable, dimension(:) :: values
1326  include "write_attribute_code.f90"
1327  end subroutine write_attribute_int_1d
1328 
1329  subroutine write_attribute_short_1d(dset, attname, values, varname, errcode)
1330  integer(2), intent(in), allocatable, dimension(:) :: values
1331  include "write_attribute_code.f90"
1332  end subroutine write_attribute_short_1d
1333 
1334  subroutine write_attribute_byte_1d(dset, attname, values, varname, errcode)
1335  integer(1), intent(in), allocatable, dimension(:) :: values
1336  include "write_attribute_code.f90"
1337  end subroutine write_attribute_byte_1d
1338 
1339  subroutine write_attribute_char(dset, attname, values, varname, errcode)
1340  character(len=*), intent(in) :: values
1341  include "write_attribute_code.f90"
1342  end subroutine write_attribute_char
1343 
1346  function get_idate_from_time_units(dset) result(idate)
1347  type(dataset), intent(in) :: dset
1348  integer idate(6)
1349  character(len=nf90_max_name) :: time_units
1350  integer ipos1,ipos2
1351  call read_attribute(dset, 'units', time_units, 'time')
1352  ipos1 = scan(time_units,"since",back=.true.)+1
1353  ipos2 = scan(time_units,"-",back=.false.)-1
1354  read(time_units(ipos1:ipos2),*) idate(1)
1355  ipos1 = ipos2+2; ipos2=ipos1+1
1356  read(time_units(ipos1:ipos2),*) idate(2)
1357  ipos1 = ipos2+2; ipos2=ipos1+1
1358  read(time_units(ipos1:ipos2),*) idate(3)
1359  ipos1 = scan(time_units,":")-2
1360  ipos2 = ipos1+1
1361  read(time_units(ipos1:ipos2),*) idate(4)
1362  ipos1 = ipos2+2
1363  ipos2 = ipos1+1
1364  read(time_units(ipos1:ipos2),*) idate(5)
1365  ipos1 = ipos2+2
1366  ipos2 = ipos1+1
1367  read(time_units(ipos1:ipos2),*) idate(6)
1368  end function get_idate_from_time_units
1373  function get_time_units_from_idate(idate, time_measure) result(time_units)
1374  character(len=*), intent(in), optional :: time_measure
1375  integer, intent(in) :: idate(6)
1376  character(len=12) :: timechar
1377  character(len=nf90_max_name) :: time_units
1378  if (present(time_measure)) then
1379  timechar = trim(time_measure)
1380  else
1381  timechar = 'hours'
1382  endif
1383  write(time_units,101) idate
1384 101 format(' since ',i4.4,'-',i2.2,'-',i2.2,' ',&
1385  i2.2,':',i2.2,':',i2.2)
1386  time_units = trim(adjustl(timechar))//time_units
1387  end function get_time_units_from_idate
1388 
1389  subroutine quantize_data_2d(dataIn, dataOut, nbits, compress_err)
1390  real(4), intent(in) :: datain(:,:)
1391  real(4), intent(out) :: dataout(:,:)
1392  include "quantize_data_code.f90"
1393  end subroutine quantize_data_2d
1394 
1395  subroutine quantize_data_3d(dataIn, dataOut, nbits, compress_err)
1396  real(4), intent(in) :: datain(:,:,:)
1397  real(4), intent(out) :: dataout(:,:,:)
1398  include "quantize_data_code.f90"
1399  end subroutine quantize_data_3d
1400 
1401  subroutine quantize_data_4d(dataIn, dataOut, nbits, compress_err)
1402  real(4), intent(in) :: datain(:,:,:,:)
1403  real(4), intent(out) :: dataout(:,:,:,:)
1404  include "quantize_data_code.f90"
1405  end subroutine quantize_data_4d
1406 
1407  subroutine quantize_data_5d(dataIn, dataOut, nbits, compress_err)
1408  real(4), intent(in) :: datain(:,:,:,:,:)
1409  real(4), intent(out) :: dataout(:,:,:,:,:)
1410  include "quantize_data_code.f90"
1411  end subroutine quantize_data_5d
1412 
1413 end module module_ncio
module_ncio::has_attr
logical function, public has_attr(dset, attname, varname)
Definition: module_ncio.f90:259
module_ncio::has_var
logical function, public has_var(dset, varname)
Definition: module_ncio.f90:239
module_ncio::get_var
type(variable) function, public get_var(dset, varname)
Get Variable object given name.
Definition: module_ncio.f90:225
module_ncio::get_time_units_from_idate
character(len=nf90_max_name) function, public get_time_units_from_idate(idate, time_measure)
create time units attribute of form 'hours since YYYY-MM-DD HH:MM:SS' from integer array with year,...
Definition: module_ncio.f90:1374
module_ncio::close_dataset
subroutine, public close_dataset(dset, errcode)
Close a netcdf file, deallocate members of dataset object.
Definition: module_ncio.f90:846
module_ncio::open_dataset
type(dataset) function, public open_dataset(filename, errcode, paropen, mpicomm)
Open existing dataset, create dataset object for reading netcdf file.
Definition: module_ncio.f90:359
module_ncio::get_dim
type(dimension) function, public get_dim(dset, dimname)
Get Dimension object given name.
Definition: module_ncio.f90:190
module_ncio::get_idate_from_time_units
integer function, dimension(6), public get_idate_from_time_units(dset)
return integer array with year,month,day,hour,minute,second parsed from time units attribute.
Definition: module_ncio.f90:1347
module_ncio::create_dataset
type(dataset) function, public create_dataset(filename, dsetin, copy_vardata, paropen, nocompress, mpicomm, errcode)
Create new dataset, using an existing dataset object to define.
Definition: module_ncio.f90:507
module_ncio::get_nvar
integer function, public get_nvar(dset, varname)
Get Variable index given name.
Definition: module_ncio.f90:289
module_ncio
writing requires a template file.
Definition: module_ncio.f90:9
module_ncio::get_ndim
integer function, public get_ndim(dset, dimname)
Get Dimension index given name.
Definition: module_ncio.f90:206