31#if defined PIO_LIB && defined DISTRIBUTE 
   53#if defined PIO_LIB && defined DISTRIBUTE 
   67#if defined PIO_LIB && defined DISTRIBUTE 
 
   76     &                         ifield, ncid, Lmulti, nfiles, S)
 
   81      logical, 
intent(in) :: Lmulti
 
   83      integer, 
intent(in) :: ng, model, job, ifield, nfiles
 
   84      integer, 
intent(in) :: Iout, Irec, Mrec
 
   85      integer, 
intent(inout) :: ncid
 
   87      TYPE (T_IO), 
intent(inout) :: S(nfiles)
 
   91      logical :: CloseFile, Liocycle, Lgridded, Lonerec
 
   92      logical :: foundit, got_var, got_time, special
 
   94      integer :: Fcount, Nrec, Tid, Tindex, Trec, Vid, Vtype
 
   95      integer :: i, ifile, lstr, my_ncid, nvatt, nvdim
 
   98      real(dp) :: Clength, Tend, Tmax, Tmin, Tmono, Tscale, Tstr, Tval
 
  101      character (len=1  ), 
parameter :: blank = 
' ' 
  102      character (len=3  ) :: label
 
  103      character (len=256) :: Fname
 
  105      character (len=*), 
parameter :: MyFile =                          &
 
  106     &  __FILE__//
", inquiry_nf90" 
  124      label=s(1)%label(1:3)
 
  125      vtype=
iinfo(1,ifield,ng)
 
  138          IF (trim(
cinfo(ifield,ng)).eq.trim(s(ifile)%name)) 
THEN 
  140              fcount=s(ifile)%Fcount+1
 
  142              fcount=s(ifile)%Fcount-1
 
  144            IF ((1.gt.fcount).and.(fcount.gt.s(ifile)%Nfiles)) 
THEN 
  147     &                            fcount, s(ifile)%Nfiles
 
  151     &                       __line__, myfile)) 
RETURN 
  153            s(ifile)%Fcount=fcount
 
  154            s(ifile)%name=trim(s(ifile)%files(fcount))
 
  157            IF (label(1:3).eq.
'FRC') 
THEN 
  161            IF (label(1:3).eq.
'BRY') 
THEN 
  165            IF (label(1:3).eq.
'CLM') 
THEN 
  171            fcount=s(ifile)%Fcount   
 
  177      IF (fcount.eq.0) 
THEN 
  179          WRITE (
stdout,20) fcount, label, trim(
vname(1,ifield))
 
  189      query: 
DO ifile=1,nfiles
 
  194        IF (s(ifile)%ncid.eq.-1) 
THEN 
  202          my_ncid=s(ifile)%ncid
 
  216     &                       myvarname = trim(
vname(1,ifield)),         &
 
  217     &                       searchvar = foundit,                       &
 
  229          IF ((nvdim.gt.1).and.(abs(job).gt.1)) 
THEN 
  237          IF (.not.lmulti) 
THEN 
  248            IF (trim(
var_aname(i)).eq.
'scale_factor') 
THEN 
  250              finfo(10,ifield,ng)=scale
 
  251            ELSE IF (trim(
var_aname(i)).eq.
'water_points') 
THEN 
  253              vtype=
iinfo(1,ifield,ng)
 
  261        IF (foundit.and.(abs(job).eq.2)) 
THEN 
  264            IF (index(trim(
var_dname(i)),
'period').ne.0) 
THEN 
  268          linfo(4,ifield,ng)=special
 
  275          IF (len_trim(
vname(5,ifield)).gt.0) 
THEN 
  299              lstr=len_trim(
tname(ifield))
 
  304     &                   
tname(ifield)(1:lstr-1)) 
THEN   
  315          IF (got_time.and.(nrec.eq.0)) 
THEN 
  322          IF (got_time.and.(nrec.eq.0)) 
THEN 
  324     &                                    trim(
vname(1,ifield)),        &
 
  330        IF (abs(job).eq.1) 
THEN 
  331          IF ((iout.eq.1).and.(nrec.gt.mrec)) 
THEN 
  342          CALL get_cycle (ng, model, ifield, job, lmulti,               &
 
  344     &                    
tdays(ng), tid, liocycle, clength, trec,      &
 
  345     &                    tstr, tend, tmin, tmax,  tscale)
 
  355          linfo( 1,ifield,ng)=lgridded
 
  356          linfo( 2,ifield,ng)=liocycle
 
  357          iinfo( 2,ifield,ng)=vid
 
  358          iinfo( 3,ifield,ng)=tid
 
  359          iinfo( 4,ifield,ng)=nrec
 
  360          iinfo( 5,ifield,ng)=vsize(1)
 
  361          iinfo( 6,ifield,ng)=vsize(2)
 
  362          iinfo( 7,ifield,ng)=vsize(3)
 
  363          iinfo(10,ifield,ng)=s(ifile)%Nfiles
 
  364          iinfo(11,ifield,ng)=nvdim
 
  365          finfo(1,ifield,ng)=tmin
 
  366          finfo(2,ifield,ng)=tmax
 
  367          finfo(3,ifield,ng)=tstr
 
  368          finfo(4,ifield,ng)=tend
 
  369          finfo(5,ifield,ng)=clength
 
  370          finfo(6,ifield,ng)=tscale
 
  371          s(ifile)%Nrec(fcount)=nrec
 
  372          s(ifile)%time_min(fcount)=tmin
 
  373          s(ifile)%time_max(fcount)=tmax
 
  388      IF (.not.got_var) 
THEN 
  389        IF ((nfiles.gt.1).and.(label(1:3).eq.
'FRC')) 
THEN 
  393              WRITE (
stdout,
'(15x,a)') trim(s(i)%name)
 
  396        ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.
'BRY')) 
THEN 
  400              WRITE (
stdout,
'(15x,a)') trim(s(i)%name)
 
  403        ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.
'CLM')) 
THEN 
  407              WRITE (
stdout,
'(15x,a)') trim(s(i)%name)
 
  417              WRITE (
stdout,
'(15x,a,a)') 
'file name is blank, ',        &
 
  418     &                                   
'cannot be determined.' 
  425      IF (.not.got_time) 
THEN 
  426        IF ((nfiles.gt.1).and.(label(1:3).eq.
'FRC')) 
THEN 
  430              WRITE (
stdout,
'(15x,a)') trim(s(i)%name)
 
  433        ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.
'BRY')) 
THEN 
  437              WRITE (
stdout,
'(15x,a)') trim(s(i)%name)
 
  440        ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.
'CLM')) 
THEN 
  444              WRITE (
stdout,
'(15x,a)') trim(s(i)%name)
 
  454              WRITE (
stdout,
'(15x,a,a)') 
'file name is blank, ',        &
 
  455     &                                   
'cannot be determined.' 
  467      IF (label(1:3).eq.
'FRC') 
THEN 
  469          IF ((trim(
ncfile).eq.trim(s(ifile)%name)).and.                &
 
  470     &        (
frcids(ifile,ng).ne.-1).and.                             &
 
  471     &        (ncid.ne.
frcids(ifile,ng))) 
THEN 
  477      IF (label(1:3).eq.
'BRY') 
THEN 
  479          IF ((trim(
ncfile).eq.trim(s(ifile)%name)).and.                &
 
  480     &        (
bryids(ifile,ng).ne.-1).and.                             &
 
  481     &        (ncid.ne.
bryids(ifile,ng))) 
THEN 
  487      IF (label(1:3).eq.
'CLM') 
THEN 
  489          IF ((trim(
ncfile).eq.trim(s(ifile)%name)).and.                &
 
  490     &        (
clmids(ifile,ng).ne.-1).and.                             &
 
  491     &        (ncid.ne.
clmids(ifile,ng))) 
THEN 
  512      IF (label(1:3).eq.
'FRC') 
THEN 
  514          IF (trim(
ncfile).eq.trim(s(ifile)%name)) 
THEN 
  520      ELSE IF (label(1:3).eq.
'BRY') 
THEN 
  522          IF (trim(
ncfile).eq.trim(s(ifile)%name)) 
THEN 
  528      ELSE IF (label(1:3).eq.
'CLM') 
THEN 
  530          IF (trim(
ncfile).eq.trim(s(ifile)%name)) 
THEN 
  548      IF (.not.lmulti) 
THEN 
  553            tindex=
iinfo(8,ifield,ng)
 
  556            IF (trec.eq.nrec) 
THEN 
  557              IF (
tdays(ng).lt.tmax) 
THEN 
  563              IF (
tdays(ng).gt.clength) 
THEN 
  574     &                            
rclock%DateNumber, tval,              &
 
  576     &                            start = (/trec-1/),                   &
 
  583            IF (tval.lt.tend) 
THEN 
  594            IF (trec.eq.nrec) 
THEN 
  595              IF (
tdays(ng).lt.tmax) 
THEN 
  598                tmono=
tdays(ng)+(tstr-clength)
 
  599                IF (tstr.eq.tmax) 
THEN 
  600                  tmono=tmono+(tmin-mod(
tdays(ng)+tmin,clength))
 
  602                  tmono=tmono+(tstr-mod(
tdays(ng)+tstr,clength))
 
  606              IF (
tdays(ng).gt.clength) 
THEN 
  618        iinfo(8,ifield,ng)=tindex
 
  619        iinfo(9,ifield,ng)=trec
 
  620        finfo(7,ifield,ng)=tmono
 
  622        iinfo(9,ifield,ng)=trec
 
  628      IF (nrec.eq.1) lonerec=.true.
 
  629      linfo(3,ifield,ng)=lonerec
 
  630      tindex=
iinfo(8,ifield,ng)
 
  637  10  
FORMAT (/,
' INQUIRY_NF90 - out of range multi-files counter for', &
 
  638     &        
' variable: ',a,/,16x,
'Fcount = ',i0,                     &
 
  639     &        
',   Expected range: 1 - ',i0)
 
  640  20  
FORMAT (/,
' INQUIRY_NF90 - unable to assign file counter, ',      &
 
  641     &        
'Fcount = ',i0,/,16x,
'while processing structure: ',a,    &
 
  642     &        /,16x,
'and variable; ',a)
 
  643  30  
FORMAT (/,
' INQUIRY_NF90 - unable to find dimension ',a,          &
 
  644     &        /,16x,
'for variable: ',a,/,16x,
'in file: ',a,             &
 
  645     &        /,16x,
'file is not CF compliant...')
 
  646  40  
FORMAT (/,
' INQUIRY_NF90 - too small dimension for variable ',a,  &
 
  648  50  
FORMAT (/,
' INQUIRY_NF90 - unable to find requested variable: ',  &
 
  650  60  
FORMAT (/,
' INQUIRY_NF90 - unable to open input NetCDF file: ',a)
 
  651  70  
FORMAT (/,
' INQUIRY_NF90 - error while reading variable: ',a,2x,  &
 
  652     &        
' at TIME index = ',i0)
 
 
  657# if defined PIO_LIB && defined DISTRIBUTE 
  661     &                        ifield, pioFile, Lmulti, nfiles, S)
 
  666      logical, 
intent(in) :: Lmulti
 
  668      integer, 
intent(in) :: ng, model, job, ifield, nfiles
 
  669      integer, 
intent(in) :: Iout, Irec, Mrec
 
  671      TYPE (File_desc_t), 
intent(inout) :: pioFile
 
  672      TYPE (T_IO),        
intent(inout) :: S(nfiles)
 
  676      logical :: CloseFile, Liocycle, Lgridded, Lonerec
 
  677      logical :: foundit, got_var, got_time, special
 
  679      integer :: Fcount, Nrec, Tindex, Trec, Vtype
 
  680      integer :: i, ifile, lstr, nvatt, nvdim
 
  683      real(dp) :: Clength, Tend, Tmax, Tmin, Tmono, Tscale, Tstr, Tval
 
  686      character (len=1  ), 
parameter :: blank = 
' ' 
  687      character (len=3  ) :: label
 
  688      character (len=256) :: Fname
 
  690      character (len=*), 
parameter :: MyFile =                          &
 
  691     &  __FILE__//
", inquiry_pio" 
  693      TYPE (File_desc_t) :: my_pioFile
 
  694      TYPE (Var_desc_t)  :: pioVar, TpioVar
 
  712      label=s(1)%label(1:3)
 
  713      vtype=
iinfo(1,ifield,ng)
 
  726          IF (trim(
cinfo(ifield,ng)).eq.trim(s(ifile)%name)) 
THEN 
  728              fcount=s(ifile)%Fcount+1
 
  730              fcount=s(ifile)%Fcount-1
 
  732            IF ((1.gt.fcount).and.(fcount.gt.s(ifile)%Nfiles)) 
THEN 
  735     &                            fcount, s(ifile)%Nfiles
 
  739     &                       __line__, myfile)) 
RETURN 
  741            s(ifile)%Fcount=fcount
 
  742            s(ifile)%name=trim(s(ifile)%files(fcount))
 
  745            IF (label(1:3).eq.
'FRC') 
THEN 
  747              s(ifile)%pioFile%fh=-1
 
  749            IF (label(1:3).eq.
'BRY') 
THEN 
  751              s(ifile)%pioFile%fh=-1
 
  753            IF (label(1:3).eq.
'CLM') 
THEN 
  755              s(ifile)%pioFile%fh=-1
 
  759            fcount=s(ifile)%Fcount   
 
  765      IF (fcount.eq.0) 
THEN 
  767          WRITE (
stdout,20) fcount, label, trim(
vname(1,ifield))
 
  777      query: 
DO ifile=1,nfiles
 
  782        IF (s(ifile)%pioFile%fh.eq.-1) 
THEN 
  790          my_piofile=s(ifile)%pioFile
 
  797     &                             piofile = my_piofile)
 
  803     &                           piofile = my_piofile,                  &
 
  804     &                           myvarname = trim(
vname(1,ifield)),     &
 
  805     &                           searchvar = foundit,                   &
 
  816          IF ((nvdim.gt.1).and.(abs(job).gt.1)) 
THEN 
  824          IF (.not.lmulti) 
THEN 
  835            IF (trim(
var_aname(i)).eq.
'scale_factor') 
THEN 
  837              finfo(10,ifield,ng)=scale
 
  838            ELSE IF (trim(
var_aname(i)).eq.
'water_points') 
THEN 
  840              vtype=
iinfo(1,ifield,ng)
 
  848        IF (foundit.and.(abs(job).eq.2)) 
THEN 
  851            IF (index(trim(
var_dname(i)),
'period').ne.0) 
THEN 
  855          linfo(4,ifield,ng)=special
 
  862          IF (len_trim(
vname(5,ifield)).gt.0) 
THEN 
  886              lstr=len_trim(
tname(ifield))
 
  891     &                   
tname(ifield)(1:lstr-1)) 
THEN   
  902          IF (got_time.and.(nrec.eq.0)) 
THEN 
  909          IF (got_time.and.(nrec.eq.0)) 
THEN 
  911     &                                    trim(
vname(1,ifield)),        &
 
  917        IF (abs(job).eq.1) 
THEN 
  918          IF ((iout.eq.1).and.(nrec.gt.mrec)) 
THEN 
  929          CALL get_cycle (ng, model, ifield, job, lmulti,               &
 
  931     &                    
tdays(ng), tpiovar, liocycle, clength, trec,  &
 
  932     &                    tstr, tend, tmin, tmax,  tscale)
 
  942          dinfo( 1,ifield,ng)%vd=piovar
 
  944          dinfo( 1,ifield,ng)%gtype=vtype
 
  945          dinfo( 2,ifield,ng)%vd=tpiovar
 
  947          dinfo( 2,ifield,ng)%gtype=0
 
  948          linfo( 1,ifield,ng)=lgridded
 
  949          linfo( 2,ifield,ng)=liocycle
 
  950          iinfo( 4,ifield,ng)=nrec
 
  951          iinfo( 5,ifield,ng)=vsize(1)
 
  952          iinfo( 6,ifield,ng)=vsize(2)
 
  953          iinfo( 7,ifield,ng)=vsize(3)
 
  954          iinfo(10,ifield,ng)=s(ifile)%Nfiles
 
  955          iinfo(11,ifield,ng)=nvdim
 
  956          finfo(1,ifield,ng)=tmin
 
  957          finfo(2,ifield,ng)=tmax
 
  958          finfo(3,ifield,ng)=tstr
 
  959          finfo(4,ifield,ng)=tend
 
  960          finfo(5,ifield,ng)=clength
 
  961          finfo(6,ifield,ng)=tscale
 
  962          s(ifile)%Nrec(fcount)=nrec
 
  963          s(ifile)%time_min(fcount)=tmin
 
  964          s(ifile)%time_max(fcount)=tmax
 
  979      IF (.not.got_var) 
THEN 
  980        IF ((nfiles.gt.1).and.(label(1:3).eq.
'FRC')) 
THEN 
  984              WRITE (
stdout,
'(15x,a)') trim(s(i)%name)
 
  987        ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.
'BRY')) 
THEN 
  991              WRITE (
stdout,
'(15x,a)') trim(s(i)%name)
 
  994        ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.
'CLM')) 
THEN 
  998              WRITE (
stdout,
'(15x,a)') trim(s(i)%name)
 
 1008              WRITE (
stdout,
'(15x,a,a)') 
'file name is blank, ',        &
 
 1009     &                                   
'cannot be determined.' 
 1016      IF (.not.got_time) 
THEN 
 1017        IF ((nfiles.gt.1).and.(label(1:3).eq.
'FRC')) 
THEN 
 1021              WRITE (
stdout,
'(15x,a)') trim(s(i)%name)
 
 1024        ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.
'BRY')) 
THEN 
 1028              WRITE (
stdout,
'(15x,a)') trim(s(i)%name)
 
 1031        ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.
'CLM')) 
THEN 
 1035              WRITE (
stdout,
'(15x,a)') trim(s(i)%name)
 
 1045              WRITE (
stdout,
'(15x,a,a)') 
'file name is blank, ',        &
 
 1046     &                                   
'cannot be determined.' 
 1058      IF (label(1:3).eq.
'FRC') 
THEN 
 1060          IF ((trim(
ncfile).eq.trim(s(ifile)%name)).and.                &
 
 1061     &        (
frcdesc(ifile,ng)%fh.ne.-1).and.                         &
 
 1062     &        (piofile%fh.ne.
frcdesc(ifile,ng)%fh)) 
THEN 
 1068      IF (label(1:3).eq.
'BRY') 
THEN 
 1070          IF ((trim(
ncfile).eq.trim(s(ifile)%name)).and.                &
 
 1071     &        (
brydesc(ifile,ng)%fh.ne.-1).and.                         &
 
 1072     &        (piofile%fh.ne.
brydesc(ifile,ng)%fh)) 
THEN 
 1078      IF (label(1:3).eq.
'CLM') 
THEN 
 1080          IF ((trim(
ncfile).eq.trim(s(ifile)%name)).and.                &
 
 1081     &        (
clmdesc(ifile,ng)%fh.ne.-1).and.                         &
 
 1082     &        (piofile%fh.ne.
clmdesc(ifile,ng)%fh)) 
THEN 
 1093      IF (piofile%fh.eq.-1) 
THEN 
 1103      IF (label(1:3).eq.
'FRC') 
THEN 
 1105          IF (trim(
ncfile).eq.trim(s(ifile)%name)) 
THEN 
 1107            s(ifile)%pioFile=piofile
 
 1111      ELSE IF (label(1:3).eq.
'BRY') 
THEN 
 1113          IF (trim(
ncfile).eq.trim(s(ifile)%name)) 
THEN 
 1115            s(ifile)%pioFile=piofile
 
 1119      ELSE IF (label(1:3).eq.
'CLM') 
THEN 
 1121          IF (trim(
ncfile).eq.trim(s(ifile)%name)) 
THEN 
 1123            s(ifile)%pioFile=piofile
 
 1128        s(1)%pioFile=piofile
 
 1139      IF (.not.lmulti) 
THEN 
 1144            tindex=
iinfo(8,ifield,ng)
 
 1147            IF (trec.eq.nrec) 
THEN 
 1148              IF (
tdays(ng).lt.tmax) 
THEN 
 1151                tmono=tend+(
tdays(ng)-mod(
tdays(ng),clength))
 
 1154              IF (
tdays(ng).gt.clength) 
THEN 
 1155                tmono=tend+(
tdays(ng)-mod(
tdays(ng),clength))
 
 1165     &                                
rclock%DateNumber, tval,          &
 
 1166     &                                piofile = piofile,                &
 
 1167     &                                start = (/trec-1/),               &
 
 1174            IF (tval.lt.tend) 
THEN 
 1185            IF (trec.eq.nrec) 
THEN 
 1186              IF (
tdays(ng).lt.tmax) 
THEN 
 1189                tmono=
tdays(ng)+(tstr-clength)
 
 1190                IF (tstr.eq.tmax) 
THEN 
 1191                  tmono=tmono+(tmin-mod(
tdays(ng)+tmin,clength))
 
 1193                  tmono=tmono+(tstr-mod(
tdays(ng)+tstr,clength))
 
 1197              IF (
tdays(ng).gt.clength) 
THEN 
 1209        iinfo(8,ifield,ng)=tindex
 
 1210        iinfo(9,ifield,ng)=trec
 
 1211        finfo(7,ifield,ng)=tmono
 
 1213        iinfo(9,ifield,ng)=trec
 
 1219      IF (nrec.eq.1) lonerec=.true.
 
 1220      linfo(3,ifield,ng)=lonerec
 
 1221      tindex=
iinfo(8,ifield,ng)
 
 1228  10  
FORMAT (/,
' INQUIRY_PIO  - out of range multi-files counter for', &
 
 1229     &        
' variable: ',a,/,16x,
'Fcount = ',i0,                     &
 
 1230     &        
',   Expected range: 1 - ',i0)
 
 1231  20  
FORMAT (/,
' INQUIRY_PIO  - unable to assign file counter, ',      &
 
 1232     &        
'Fcount = ',i0,/,16x,
'while processing structure: ',a,    &
 
 1233     &        /,16x,
'and variable; ',a)
 
 1234  30  
FORMAT (/,
' INQUIRY_PIO  - unable to find dimension ',a,          &
 
 1235     &        /,16x,
'for variable: ',a,/,16x,
'in file: ',a,             &
 
 1236     &        /,16x,
'file is not CF compliant...')
 
 1237  40  
FORMAT (/,
' INQUIRY_PIO  - too small dimension for variable ',a,  &
 
 1239  50  
FORMAT (/,
' INQUIRY_PIO  - unable to find requested variable: ',  &
 
 1241  60  
FORMAT (/,
' INQUIRY_PIO  - unable to open input NetCDF file: ',a)
 
 1242  70  
FORMAT (/,
' INQUIRY_PIO  - error while reading variable: ',a,2x,  &
 
 1243     &        
' at TIME index = ',i0)
 
 
subroutine inquiry_nf90(ng, model, job, iout, irec, mrec, ifield, ncid, lmulti, nfiles, s)
subroutine inquiry_pio(ng, model, job, iout, irec, mrec, ifield, piofile, lmulti, nfiles, s)
integer, dimension(:,:), allocatable clmids
type(file_desc_t), dimension(:,:), allocatable brydesc
integer, dimension(:,:), allocatable bryids
character(len=256) ncfile
type(file_desc_t), dimension(:,:), allocatable clmdesc
integer, dimension(:,:), allocatable frcids
type(file_desc_t), dimension(:,:), allocatable frcdesc
character(len=256) sourcefile
character(len=256), dimension(:,:), allocatable cinfo
logical, dimension(:,:,:), allocatable linfo
type(my_vardesc), dimension(:,:,:), pointer dinfo
real(dp), dimension(:,:,:), allocatable vtime
character(len=46), dimension(0:nv) tname
real(dp), dimension(:,:,:), allocatable finfo
character(len=maxlen), dimension(6, 0:nv) vname
integer, dimension(:,:,:), allocatable iinfo
subroutine, public netcdf_close(ng, model, ncid, ncname, lupdate)
character(len=1024), dimension(nvara) var_achar
subroutine, public netcdf_check_dim(ng, model, ncname, ncid)
character(len=100), dimension(mdims) dim_name
subroutine, public netcdf_open(ng, model, ncname, omode, ncid)
integer, dimension(mdims) dim_size
real(r8), dimension(nvara) var_afloat
integer, dimension(nvard) var_dsize
character(len=100), dimension(nvard) var_dname
character(len=100), dimension(nvara) var_aname
subroutine, public netcdf_inq_var(ng, model, ncname, ncid, myvarname, searchvar, varid, nvardim, nvaratt)
integer, parameter pio_type
subroutine, public pio_netcdf_inq_var(ng, model, ncname, piofile, myvarname, searchvar, piovar, nvardim, nvaratt)
subroutine, public pio_netcdf_close(ng, model, piofile, ncname, lupdate)
subroutine, public pio_netcdf_open(ng, model, ncname, omode, piofile)
subroutine, public pio_netcdf_check_dim(ng, model, ncname, piofile)
integer, parameter pio_tout
real(dp), parameter day2sec
real(dp), dimension(:), allocatable tdays
character(len(sinp)) function, public lowercase(sinp)
logical function, public founderror(flag, noerr, line, routine)