ROMS
Loading...
Searching...
No Matches
inquiry_mod::inquiry Interface Reference

Public Member Functions

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)
 

Detailed Description

Definition at line 65 of file inquiry.F.

Member Function/Subroutine Documentation

◆ inquiry_nf90()

subroutine inquiry_mod::inquiry::inquiry_nf90 ( integer, intent(in) ng,
integer, intent(in) model,
integer, intent(in) job,
integer, intent(in) iout,
integer, intent(in) irec,
integer, intent(in) mrec,
integer, intent(in) ifield,
integer, intent(inout) ncid,
logical, intent(in) lmulti,
integer, intent(in) nfiles,
type (t_io), dimension(nfiles), intent(inout) s )

Definition at line 75 of file inquiry.F.

77!***********************************************************************
78!
79! Imported variable declarations.
80!
81 logical, intent(in) :: Lmulti
82!
83 integer, intent(in) :: ng, model, job, ifield, nfiles
84 integer, intent(in) :: Iout, Irec, Mrec
85 integer, intent(inout) :: ncid
86!
87 TYPE (T_IO), intent(inout) :: S(nfiles)
88!
89! Local variable declarations.
90!
91 logical :: CloseFile, Liocycle, Lgridded, Lonerec
92 logical :: foundit, got_var, got_time, special
93!
94 integer :: Fcount, Nrec, Tid, Tindex, Trec, Vid, Vtype
95 integer :: i, ifile, lstr, my_ncid, nvatt, nvdim
96 integer :: Vsize(4)
97!
98 real(dp) :: Clength, Tend, Tmax, Tmin, Tmono, Tscale, Tstr, Tval
99 real(dp) :: scale
100!
101 character (len=1 ), parameter :: blank = ' '
102 character (len=3 ) :: label
103 character (len=256) :: Fname
104
105 character (len=*), parameter :: MyFile = &
106 & __FILE__//", inquiry_nf90"
107!
108 sourcefile=myfile
109!
110!-----------------------------------------------------------------------
111! On first call, inquire about the contents of input NetCDF file.
112!-----------------------------------------------------------------------
113!
114! Intialize local variables.
115!
116 vid=-1
117 tid=-1
118 nrec=0
119 liocycle=.false.
120 lgridded=.false.
121 lonerec=.false.
122 got_var=.false.
123 got_time=.false.
124 label=s(1)%label(1:3)
125 vtype=iinfo(1,ifield,ng)
126!
127! Clear ncfile module variable.
128!
129 DO i=1,len(ncfile)
130 ncfile(i:i)=blank
131 END DO
132!
133! If multi-files, increase (decrease if backward logic) file counter
134! and set new file names.
135!
136 IF (lmulti) THEN
137 DO ifile=1,nfiles
138 IF (trim(cinfo(ifield,ng)).eq.trim(s(ifile)%name)) THEN
139 IF (job.gt.0) THEN
140 fcount=s(ifile)%Fcount+1
141 ELSE
142 fcount=s(ifile)%Fcount-1
143 END IF
144 IF ((1.gt.fcount).and.(fcount.gt.s(ifile)%Nfiles)) THEN
145 IF (master) THEN
146 WRITE (stdout,10) trim(vname(1,ifield)), &
147 & fcount, s(ifile)%Nfiles
148 END IF
149 exit_flag=4
150 IF (founderror(exit_flag, noerror, &
151 & __line__, myfile)) RETURN
152 END IF
153 s(ifile)%Fcount=fcount
154 s(ifile)%name=trim(s(ifile)%files(fcount))
155 CALL netcdf_close (ng, model, ncid)
156 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
157 IF (label(1:3).eq.'FRC') THEN
158 frcids(ifile,ng)=-1
159 s(ifile)%ncid=-1
160 END IF
161 IF (label(1:3).eq.'BRY') THEN
162 bryids(ifile,ng)=-1
163 s(ifile)%ncid=-1
164 END IF
165 IF (label(1:3).eq.'CLM') THEN
166 clmids(ifile,ng)=-1
167 s(ifile)%ncid=-1
168 END IF
169 EXIT
170 ELSE ! S(ifile)%name and Fcount already
171 fcount=s(ifile)%Fcount ! updated in first field processed
172 END IF ! for current new file
173 END DO
174 ELSE
175 fcount=s(1)%Fcount
176 END IF
177 IF (fcount.eq.0) THEN
178 IF (master) THEN
179 WRITE (stdout,20) fcount, label, trim(vname(1,ifield))
180 END IF
181 exit_flag=4
182 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
183 END IF
184!
185! If several input NetCDF files (nfiles>1), scan files until the
186! requested variable is found.
187!
188 foundit=.false.
189 query: DO ifile=1,nfiles
190 fname=s(ifile)%name
191!
192! Open NetCDF file for reading.
193!
194 IF (s(ifile)%ncid.eq.-1) THEN
195 CALL netcdf_open (ng, model, fname, 0, my_ncid)
196 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
197 IF (master) WRITE (stdout,60) trim(fname)
198 RETURN
199 END IF
200 closefile=.true.
201 ELSE
202 my_ncid=s(ifile)%ncid
203 closefile=.false.
204 END IF
205!
206! Inquire about the dimensions and check for consistency.
207!
208 CALL netcdf_check_dim (ng, model, fname, &
209 & ncid = my_ncid)
210 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
211!
212! Inquire about requested variable.
213!
214 CALL netcdf_inq_var (ng, model, fname, &
215 & ncid = my_ncid, &
216 & myvarname = trim(vname(1,ifield)), &
217 & searchvar = foundit, &
218 & varid = vid, &
219 & nvardim = nvdim, &
220 & nvaratt = nvatt)
221
222 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
223!
224! Determine if gridded or point data. Set variable dimensions.
225!
226 IF (foundit) THEN
227 got_var=.true.
228 ncfile=fname
229 IF ((nvdim.gt.1).and.(abs(job).gt.1)) THEN
230 lgridded=.true.
231 END IF
232 vsize=0
233 DO i=1,nvdim
234 vsize(i)=var_dsize(i)
235 END DO
236 ELSE
237 IF (.not.lmulti) THEN
238 ncfile=fname ! need for error report
239 END IF
240 END IF
241!
242! If "scale_factor" attribute is present for a variable, the data are
243! to be multiplied by this factor. Check if only water points are
244! available.
245!
246 IF (foundit) THEN
247 DO i=1,nvatt
248 IF (trim(var_aname(i)).eq.'scale_factor') THEN
249 scale=var_afloat(i)
250 finfo(10,ifield,ng)=scale
251 ELSE IF (trim(var_aname(i)).eq.'water_points') THEN
252 iinfo(1,ifield,ng)=-abs(iinfo(1,ifield,ng))
253 vtype=iinfo(1,ifield,ng)
254 END IF
255 END DO
256 END IF
257!
258! Check if processing special 2D fields with additional dimensions
259! (for example, tides).
260!
261 IF (foundit.and.(abs(job).eq.2)) THEN
262 special=.false.
263 DO i=1,nvdim
264 IF (index(trim(var_dname(i)),'period').ne.0) THEN
265 special=.true.
266 END IF
267 END DO
268 linfo(4,ifield,ng)=special
269 END IF
270!
271! Inquire about associated time dimension, if any, and get number of
272! available time records.
273!
274 IF (foundit) THEN
275 IF (len_trim(vname(5,ifield)).gt.0) THEN
276 tname(ifield)=trim(vname(5,ifield))
277 DO i=1,nvdim
278 IF (var_dname(i).eq.trim(tname(ifield))) THEN
279 nrec=var_dsize(i)
280 got_time=.true.
281 END IF
282 END DO
283 END IF
284!
285! If the associated time variable is different to that specified in
286! "varinfo.yaml" (Nrec still zero), reset associated time-variable name
287! to that specified in the "time" attribute. If founded, also check
288! for a C-language null termination character in the time variable
289! string, which is possible in NetCDF files generated from C-programs.
290!
291 IF (nrec.eq.0) THEN
292 DO i=1,nvatt
293 IF (trim(var_aname(i)).eq.'time') THEN
294 tname(ifield)=trim(var_achar(i))
295 got_time=.true.
296 END IF
297 END DO
298 IF (got_time) THEN
299 lstr=len_trim(tname(ifield))
300 DO i=1,n_dim
301 IF (trim(dim_name(i)).eq.tname(ifield)(1:lstr)) THEN
302 nrec=dim_size(i)
303 ELSE IF (trim(dim_name(i)).eq. &
304 & tname(ifield)(1:lstr-1)) THEN ! There is a C
305 nrec=dim_size(i) ! language null
306 tname(ifield)=tname(ifield)(1:lstr-1) ! termination
307 END IF ! character
308 END DO
309 END IF
310 END IF
311!
312! If Nrec=0, input file is not CF compliant, check variable dimension
313! to see if the dimension contains the "time" string.
314!
315 IF (got_time.and.(nrec.eq.0)) THEN
316 DO i=1,n_vdim
317 IF (index(trim(lowercase(var_dname(i))),'time').ne.0) THEN
318 nrec=var_dsize(i)
319 END IF
320 END DO
321 END IF
322 IF (got_time.and.(nrec.eq.0)) THEN
323 IF (master) WRITE (stdout,30) trim(tname(ifield)), &
324 & trim(vname(1,ifield)), &
325 & trim(fname)
326 exit_flag=4
327 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
328 END IF
329 END IF
330 IF (abs(job).eq.1) THEN
331 IF ((iout.eq.1).and.(nrec.gt.mrec)) THEN
332 IF (master) WRITE (stdout,40) trim(vname(1,ifield)), &
333 & mrec, nrec
334 exit_flag=4
335 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
336 END IF
337 END IF
338!
339! Determine initial time record to read and cycling switch.
340!
341 IF (foundit) THEN
342 CALL get_cycle (ng, model, ifield, job, lmulti, &
343 & ncfile, my_ncid, tname(ifield), nrec, &
344 & tdays(ng), tid, liocycle, clength, trec, &
345 & tstr, tend, tmin, tmax, tscale)
346 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
347 sourcefile=myfile
348 END IF
349!
350! Store variable information into global information arrays. Also
351! fill some structure values.
352!
353 IF (foundit) THEN
354 ncid=my_ncid
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
374 EXIT query
375 END IF
376!
377! Close input NetCDF file if opened during the query. Files opened
378! outside the query loop remain open. This is done to avoid opening
379! too many files.
380!
381 IF (closefile) THEN
382 CALL netcdf_close (ng, model, my_ncid, fname, .false.)
383 END IF
384 END DO query
385!
386! Terminate execution requested variables are not found.
387!
388 IF (.not.got_var) THEN
389 IF ((nfiles.gt.1).and.(label(1:3).eq.'FRC')) THEN
390 IF (master) THEN
391 WRITE (stdout,50) trim(vname(1,ifield)), 'files:'
392 DO i=1,nfiles
393 WRITE (stdout,'(15x,a)') trim(s(i)%name)
394 END DO
395 END IF
396 ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.'BRY')) THEN
397 IF (master) THEN
398 WRITE (stdout,50) trim(vname(1,ifield)), 'files:'
399 DO i=1,nfiles
400 WRITE (stdout,'(15x,a)') trim(s(i)%name)
401 END DO
402 END IF
403 ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.'CLM')) THEN
404 IF (master) THEN
405 WRITE (stdout,50) trim(vname(1,ifield)), 'files:'
406 DO i=1,nfiles
407 WRITE (stdout,'(15x,a)') trim(s(i)%name)
408 END DO
409 END IF
410 ELSE
411 lstr=len_trim(ncfile)
412 IF (master) THEN
413 WRITE (stdout,50) trim(vname(1,ifield)), 'file:'
414 IF (lstr.gt.0) THEN
415 WRITE (stdout,'(15x,a)') trim(ncfile)
416 ELSE
417 WRITE (stdout,'(15x,a,a)') 'file name is blank, ', &
418 & 'cannot be determined.'
419 END IF
420 END IF
421 END IF
422 exit_flag=2
423 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
424 END IF
425 IF (.not.got_time) THEN
426 IF ((nfiles.gt.1).and.(label(1:3).eq.'FRC')) THEN
427 IF (master) THEN
428 WRITE (stdout,50) trim(tname(ifield)), 'files:'
429 DO i=1,nfiles
430 WRITE (stdout,'(15x,a)') trim(s(i)%name)
431 END DO
432 END IF
433 ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.'BRY')) THEN
434 IF (master) THEN
435 WRITE (stdout,50) trim(tname(ifield)), 'files:'
436 DO i=1,nfiles
437 WRITE (stdout,'(15x,a)') trim(s(i)%name)
438 END DO
439 END IF
440 ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.'CLM')) THEN
441 IF (master) THEN
442 WRITE (stdout,50) trim(tname(ifield)), 'files:'
443 DO i=1,nfiles
444 WRITE (stdout,'(15x,a)') trim(s(i)%name)
445 END DO
446 END IF
447 ELSE
448 lstr=len_trim(ncfile)
449 IF (master) THEN
450 WRITE (stdout,50) trim(tname(ifield)), 'file:'
451 IF (lstr.gt.0) THEN
452 WRITE (stdout,'(15x,a)') trim(ncfile)
453 ELSE
454 WRITE (stdout,'(15x,a,a)') 'file name is blank, ', &
455 & 'cannot be determined.'
456 END IF
457 END IF
458 END IF
459 exit_flag=2
460 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
461 END IF
462!
463! If processing model forcing (multiple files allowed), check if file
464! for requested field (based on ncid) has been already opened and get
465! its ID from the association table (FRCids).
466!
467 IF (label(1:3).eq.'FRC') THEN
468 DO ifile=1,nfiles
469 IF ((trim(ncfile).eq.trim(s(ifile)%name)).and. &
470 & (frcids(ifile,ng).ne.-1).and. &
471 & (ncid.ne.frcids(ifile,ng))) THEN
472 ncid=frcids(ifile,ng)
473 EXIT
474 END IF
475 END DO
476 END IF
477 IF (label(1:3).eq.'BRY') THEN
478 DO ifile=1,nfiles
479 IF ((trim(ncfile).eq.trim(s(ifile)%name)).and. &
480 & (bryids(ifile,ng).ne.-1).and. &
481 & (ncid.ne.bryids(ifile,ng))) THEN
482 ncid=bryids(ifile,ng)
483 EXIT
484 END IF
485 END DO
486 END IF
487 IF (label(1:3).eq.'CLM') THEN
488 DO ifile=1,nfiles
489 IF ((trim(ncfile).eq.trim(s(ifile)%name)).and. &
490 & (clmids(ifile,ng).ne.-1).and. &
491 & (ncid.ne.clmids(ifile,ng))) THEN
492 ncid=clmids(ifile,ng)
493 EXIT
494 END IF
495 END DO
496 END IF
497!
498! If appropriate, open input NetCDF file for reading.
499! (HGA: it is not needed since file was already open and
500! ncid=my_ncid was assigned above.
501!
502 IF (ncid.eq.-1) THEN
503 CALL netcdf_open (ng, model, ncfile, 0, ncid)
504 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
505 IF (master) WRITE (stdout,60) trim(ncfile)
506 RETURN
507 END IF
508 END IF
509!
510! If processing model forcing, load ID into association table (FRCids).
511!
512 IF (label(1:3).eq.'FRC') THEN
513 DO ifile=1,nfiles
514 IF (trim(ncfile).eq.trim(s(ifile)%name)) THEN
515 frcids(ifile,ng)=ncid
516 s(ifile)%ncid=ncid
517 EXIT
518 END IF
519 END DO
520 ELSE IF (label(1:3).eq.'BRY') THEN
521 DO ifile=1,nfiles
522 IF (trim(ncfile).eq.trim(s(ifile)%name)) THEN
523 bryids(ifile,ng)=ncid
524 s(ifile)%ncid=ncid
525 EXIT
526 END IF
527 END DO
528 ELSE IF (label(1:3).eq.'CLM') THEN
529 DO ifile=1,nfiles
530 IF (trim(ncfile).eq.trim(s(ifile)%name)) THEN
531 clmids(ifile,ng)=ncid
532 s(ifile)%ncid=ncid
533 EXIT
534 END IF
535 END DO
536 ELSE
537 s(1)%ncid=ncid
538 END IF
539 cinfo(ifield,ng)=trim(ncfile)
540!
541! The strategy here is to create a local, monotonically increasing
542! time variable so the interpolation between snapshots is trivial
543! when cycling forcing fields. Subtract one to time record counter
544! "Trec" to avoid doing special case at initialization. If the
545! job task is negative, the logic is backaward in time for the
546! adjoint model.
547!
548 IF (.not.lmulti) THEN
549 IF (job.lt.0) THEN ! backward time logic, adjoint
550 IF (irec.eq.1) THEN
551 tindex=iout
552 ELSE
553 tindex=iinfo(8,ifield,ng)
554 END IF
555 IF (liocycle) THEN
556 IF (trec.eq.nrec) THEN
557 IF (tdays(ng).lt.tmax) THEN
558 tmono=tend
559 ELSE
560 tmono=tend+(tdays(ng)-mod(tdays(ng),clength))
561 END IF
562 ELSE
563 IF (tdays(ng).gt.clength) THEN
564 tmono=tend+(tdays(ng)-mod(tdays(ng),clength))
565 ELSE
566 tmono=tend
567 END IF
568 END IF
569 trec=trec+2
570 ELSE
571 tmono=tend
572 trec=trec+1
573 CALL netcdf_get_time (ng, model, ncfile, tname(ifield), &
574 & rclock%DateNumber, tval, &
575 & ncid = ncid, &
576 & start = (/trec-1/), &
577 & total = (/1/))
578 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
579 IF (master) WRITE (stdout,70) trim(tname(ifield)), trec
580 RETURN
581 END IF
582 tval=tval*tscale
583 IF (tval.lt.tend) THEN
584 trec=trec+1
585 END IF
586 END IF
587 ELSE ! forward time logic
588 IF (irec.eq.1) THEN
589 tindex=iout
590 ELSE
591 tindex=1
592 END IF
593 IF (liocycle) THEN
594 IF (trec.eq.nrec) THEN
595 IF (tdays(ng).lt.tmax) THEN
596 tmono=tstr-clength
597 ELSE
598 tmono=tdays(ng)+(tstr-clength)
599 IF (tstr.eq.tmax) THEN
600 tmono=tmono+(tmin-mod(tdays(ng)+tmin,clength))
601 ELSE
602 tmono=tmono+(tstr-mod(tdays(ng)+tstr,clength))
603 END IF
604 END IF
605 ELSE
606 IF (tdays(ng).gt.clength) THEN
607 tmono=tdays(ng)-mod(tdays(ng)-tstr,clength)
608 ELSE
609 tmono=tstr
610 END IF
611 END IF
612 ELSE
613 tmono=tstr
614 END IF
615 trec=trec-1
616 END IF
617 tmono=tmono*day2sec
618 iinfo(8,ifield,ng)=tindex
619 iinfo(9,ifield,ng)=trec
620 finfo(7,ifield,ng)=tmono
621 ELSE
622 iinfo(9,ifield,ng)=trec
623 END IF
624!
625! Set switch for one time record dataset. In this case, the input
626! data is always the same and time interpolation is not performed.
627!
628 IF (nrec.eq.1) lonerec=.true.
629 linfo(3,ifield,ng)=lonerec
630 tindex=iinfo(8,ifield,ng)
631 IF (job.lt.0) THEN
632 vtime(tindex,ifield,ng)=finfo(4,ifield,ng)
633 ELSE
634 vtime(tindex,ifield,ng)=finfo(3,ifield,ng)
635 END IF
636!
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, &
647 & ': ',i0,2x,i0)
648 50 FORMAT (/,' INQUIRY_NF90 - unable to find requested variable: ', &
649 & a,/,16x,'in ',a)
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)
653!
654 RETURN

References mod_iounits::bryids, mod_ncparam::cinfo, mod_iounits::clmids, mod_scalars::day2sec, mod_netcdf::dim_name, mod_netcdf::dim_size, mod_scalars::exit_flag, mod_ncparam::finfo, strings_mod::founderror(), mod_iounits::frcids, mod_ncparam::iinfo, mod_ncparam::linfo, strings_mod::lowercase(), mod_parallel::master, mod_netcdf::n_dim, mod_netcdf::n_vdim, mod_iounits::ncfile, mod_netcdf::netcdf_check_dim(), mod_netcdf::netcdf_close(), mod_netcdf::netcdf_inq_var(), mod_netcdf::netcdf_open(), mod_scalars::noerror, mod_scalars::rclock, mod_iounits::sourcefile, mod_iounits::stdout, mod_scalars::tdays, mod_ncparam::tname, mod_netcdf::var_achar, mod_netcdf::var_afloat, mod_netcdf::var_aname, mod_netcdf::var_dname, mod_netcdf::var_dsize, mod_ncparam::vname, and mod_ncparam::vtime.

Here is the call graph for this function:

◆ inquiry_pio()

subroutine inquiry_mod::inquiry::inquiry_pio ( integer, intent(in) ng,
integer, intent(in) model,
integer, intent(in) job,
integer, intent(in) iout,
integer, intent(in) irec,
integer, intent(in) mrec,
integer, intent(in) ifield,
type (file_desc_t), intent(inout) piofile,
logical, intent(in) lmulti,
integer, intent(in) nfiles,
type (t_io), dimension(nfiles), intent(inout) s )

Definition at line 660 of file inquiry.F.

662!***********************************************************************
663!
664! Imported variable declarations.
665!
666 logical, intent(in) :: Lmulti
667!
668 integer, intent(in) :: ng, model, job, ifield, nfiles
669 integer, intent(in) :: Iout, Irec, Mrec
670!
671 TYPE (File_desc_t), intent(inout) :: pioFile
672 TYPE (T_IO), intent(inout) :: S(nfiles)
673!
674! Local variable declarations.
675!
676 logical :: CloseFile, Liocycle, Lgridded, Lonerec
677 logical :: foundit, got_var, got_time, special
678!
679 integer :: Fcount, Nrec, Tindex, Trec, Vtype
680 integer :: i, ifile, lstr, nvatt, nvdim
681 integer :: Vsize(4)
682!
683 real(dp) :: Clength, Tend, Tmax, Tmin, Tmono, Tscale, Tstr, Tval
684 real(dp) :: scale
685!
686 character (len=1 ), parameter :: blank = ' '
687 character (len=3 ) :: label
688 character (len=256) :: Fname
689
690 character (len=*), parameter :: MyFile = &
691 & __FILE__//", inquiry_pio"
692!
693 TYPE (File_desc_t) :: my_pioFile
694 TYPE (Var_desc_t) :: pioVar, TpioVar
695!
696 sourcefile=myfile
697!
698!-----------------------------------------------------------------------
699! On first call, inquire about the contents of input NetCDF file.
700!-----------------------------------------------------------------------
701!
702! Intialize local variables.
703!
704 piovar%varID=-1
705 tpiovar%varID=-1
706 nrec=0
707 liocycle=.false.
708 lgridded=.false.
709 lonerec=.false.
710 got_var=.false.
711 got_time=.false.
712 label=s(1)%label(1:3)
713 vtype=iinfo(1,ifield,ng)
714!
715! Clear ncfile module variable.
716!
717 DO i=1,len(ncfile)
718 ncfile(i:i)=blank
719 END DO
720!
721! If multi-files, increase (decrease if backward logic) file counter
722! and set new file names.
723!
724 IF (lmulti) THEN
725 DO ifile=1,nfiles
726 IF (trim(cinfo(ifield,ng)).eq.trim(s(ifile)%name)) THEN
727 IF (job.gt.0) THEN
728 fcount=s(ifile)%Fcount+1
729 ELSE
730 fcount=s(ifile)%Fcount-1
731 END IF
732 IF ((1.gt.fcount).and.(fcount.gt.s(ifile)%Nfiles)) THEN
733 IF (master) THEN
734 WRITE (stdout,10) trim(vname(1,ifield)), &
735 & fcount, s(ifile)%Nfiles
736 END IF
737 exit_flag=4
738 IF (founderror(exit_flag, noerror, &
739 & __line__, myfile)) RETURN
740 END IF
741 s(ifile)%Fcount=fcount
742 s(ifile)%name=trim(s(ifile)%files(fcount))
743 CALL pio_netcdf_close (ng, model, piofile)
744 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
745 IF (label(1:3).eq.'FRC') THEN
746 frcdesc(ifile,ng)%fh=-1
747 s(ifile)%pioFile%fh=-1
748 END IF
749 IF (label(1:3).eq.'BRY') THEN
750 brydesc(ifile,ng)%fh=-1
751 s(ifile)%pioFile%fh=-1
752 END IF
753 IF (label(1:3).eq.'CLM') THEN
754 clmdesc(ifile,ng)%fh=-1
755 s(ifile)%pioFile%fh=-1
756 END IF
757 EXIT
758 ELSE ! S(ifile)%name and Fcount already
759 fcount=s(ifile)%Fcount ! updated in first field processed
760 END IF ! for current new file
761 END DO
762 ELSE
763 fcount=s(1)%Fcount
764 END IF
765 IF (fcount.eq.0) THEN
766 IF (master) THEN
767 WRITE (stdout,20) fcount, label, trim(vname(1,ifield))
768 END IF
769 exit_flag=4
770 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
771 END IF
772!
773! If several input NetCDF files (nfiles>1), scan files until the
774! requested variable is found.
775!
776 foundit=.false.
777 query: DO ifile=1,nfiles
778 fname=s(ifile)%name
779!
780! Open NetCDF file for reading.
781!
782 IF (s(ifile)%pioFile%fh.eq.-1) THEN
783 CALL pio_netcdf_open (ng, model, fname, 0, my_piofile)
784 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
785 IF (master) WRITE (stdout,60) trim(fname)
786 RETURN
787 END IF
788 closefile=.true.
789 ELSE
790 my_piofile=s(ifile)%pioFile
791 closefile=.false.
792 END IF
793!
794! Inquire about the dimensions and check for consistency.
795!
796 CALL pio_netcdf_check_dim (ng, model, fname, &
797 & piofile = my_piofile)
798 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
799!
800! Inquire about requested variable.
801!
802 CALL pio_netcdf_inq_var (ng, model, fname, &
803 & piofile = my_piofile, &
804 & myvarname = trim(vname(1,ifield)), &
805 & searchvar = foundit, &
806 & piovar = piovar, &
807 & nvardim = nvdim, &
808 & nvaratt = nvatt)
809 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
810!
811! Determine if gridded or point data. Set variable dimensions.
812!
813 IF (foundit) THEN
814 got_var=.true.
815 ncfile=fname
816 IF ((nvdim.gt.1).and.(abs(job).gt.1)) THEN
817 lgridded=.true.
818 END IF
819 vsize=0
820 DO i=1,nvdim
821 vsize(i)=var_dsize(i)
822 END DO
823 ELSE
824 IF (.not.lmulti) THEN
825 ncfile=fname ! need for error report
826 END IF
827 END IF
828!
829! If "scale_factor" attribute is present for a variable, the data are
830! to be multiplied by this factor. Check if only water points are
831! available.
832!
833 IF (foundit) THEN
834 DO i=1,nvatt
835 IF (trim(var_aname(i)).eq.'scale_factor') THEN
836 scale=var_afloat(i)
837 finfo(10,ifield,ng)=scale
838 ELSE IF (trim(var_aname(i)).eq.'water_points') THEN
839 iinfo(1,ifield,ng)=-abs(iinfo(1,ifield,ng))
840 vtype=iinfo(1,ifield,ng)
841 END IF
842 END DO
843 END IF
844!
845! Check if processing special 2D fields with additional dimensions
846! (for example, tides).
847!
848 IF (foundit.and.(abs(job).eq.2)) THEN
849 special=.false.
850 DO i=1,nvdim
851 IF (index(trim(var_dname(i)),'period').ne.0) THEN
852 special=.true.
853 END IF
854 END DO
855 linfo(4,ifield,ng)=special
856 END IF
857!
858! Inquire about associated time dimension, if any, and get number of
859! available time records.
860!
861 IF (foundit) THEN
862 IF (len_trim(vname(5,ifield)).gt.0) THEN
863 tname(ifield)=trim(vname(5,ifield))
864 DO i=1,nvdim
865 IF (var_dname(i).eq.trim(tname(ifield))) THEN
866 nrec=var_dsize(i)
867 got_time=.true.
868 END IF
869 END DO
870 END IF
871!
872! If the associated time variable is different to that specified in
873! "varinfo.yaml" (Nrec still zero), reset associated time-variable name
874! to that specified in the "time" attribute. If founded, also check
875! for a C-language null termination character in the time variable
876! string, which is possible in NetCDF files generated from C-programs.
877!
878 IF (nrec.eq.0) THEN
879 DO i=1,nvatt
880 IF (trim(var_aname(i)).eq.'time') THEN
881 tname(ifield)=trim(var_achar(i))
882 got_time=.true.
883 END IF
884 END DO
885 IF (got_time) THEN
886 lstr=len_trim(tname(ifield))
887 DO i=1,n_dim
888 IF (trim(dim_name(i)).eq.tname(ifield)(1:lstr)) THEN
889 nrec=dim_size(i)
890 ELSE IF (trim(dim_name(i)).eq. &
891 & tname(ifield)(1:lstr-1)) THEN ! There is a C
892 nrec=dim_size(i) ! language null
893 tname(ifield)=tname(ifield)(1:lstr-1) ! termination
894 END IF ! character
895 END DO
896 END IF
897 END IF
898!
899! If Nrec=0, input file is not CF compliant, check variable dimension
900! to see if the dimension contains the "time" string.
901!
902 IF (got_time.and.(nrec.eq.0)) THEN
903 DO i=1,n_vdim
904 IF (index(trim(lowercase(var_dname(i))),'time').ne.0) THEN
905 nrec=var_dsize(i)
906 END IF
907 END DO
908 END IF
909 IF (got_time.and.(nrec.eq.0)) THEN
910 IF (master) WRITE (stdout,30) trim(tname(ifield)), &
911 & trim(vname(1,ifield)), &
912 & trim(fname)
913 exit_flag=4
914 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
915 END IF
916 END IF
917 IF (abs(job).eq.1) THEN
918 IF ((iout.eq.1).and.(nrec.gt.mrec)) THEN
919 IF (master) WRITE (stdout,40) trim(vname(1,ifield)), &
920 & mrec, nrec
921 exit_flag=4
922 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
923 END IF
924 END IF
925!
926! Determine initial time record to read and cycling switch.
927!
928 IF (foundit) THEN
929 CALL get_cycle (ng, model, ifield, job, lmulti, &
930 & ncfile, my_piofile, tname(ifield), nrec, &
931 & tdays(ng), tpiovar, liocycle, clength, trec, &
932 & tstr, tend, tmin, tmax, tscale)
933 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
934 sourcefile=myfile
935 END IF
936!
937! Store variable information into global information arrays. Also
938! fill some structure values.
939!
940 IF (foundit) THEN
941 piofile=my_piofile
942 dinfo( 1,ifield,ng)%vd=piovar
943 dinfo( 1,ifield,ng)%dkind=pio_type
944 dinfo( 1,ifield,ng)%gtype=vtype
945 dinfo( 2,ifield,ng)%vd=tpiovar
946 dinfo( 1,ifield,ng)%dkind=pio_tout
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
965 EXIT query
966 END IF
967!
968! Close input NetCDF file if opened during the query. Files opened
969! outside the query loop remain open. This is done to avoid opening
970! too many files.
971!
972 IF (closefile) THEN
973 CALL pio_netcdf_close (ng, model, my_piofile, fname, .false.)
974 END IF
975 END DO query
976!
977! Terminate execution requested variables are not found.
978!
979 IF (.not.got_var) THEN
980 IF ((nfiles.gt.1).and.(label(1:3).eq.'FRC')) THEN
981 IF (master) THEN
982 WRITE (stdout,50) trim(vname(1,ifield)), 'files:'
983 DO i=1,nfiles
984 WRITE (stdout,'(15x,a)') trim(s(i)%name)
985 END DO
986 END IF
987 ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.'BRY')) THEN
988 IF (master) THEN
989 WRITE (stdout,50) trim(vname(1,ifield)), 'files:'
990 DO i=1,nfiles
991 WRITE (stdout,'(15x,a)') trim(s(i)%name)
992 END DO
993 END IF
994 ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.'CLM')) THEN
995 IF (master) THEN
996 WRITE (stdout,50) trim(vname(1,ifield)), 'files:'
997 DO i=1,nfiles
998 WRITE (stdout,'(15x,a)') trim(s(i)%name)
999 END DO
1000 END IF
1001 ELSE
1002 lstr=len_trim(ncfile)
1003 IF (master) THEN
1004 WRITE (stdout,50) trim(vname(1,ifield)), 'file:'
1005 IF (lstr.gt.0) THEN
1006 WRITE (stdout,'(15x,a)') trim(ncfile)
1007 ELSE
1008 WRITE (stdout,'(15x,a,a)') 'file name is blank, ', &
1009 & 'cannot be determined.'
1010 END IF
1011 END IF
1012 END IF
1013 exit_flag=2
1014 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1015 END IF
1016 IF (.not.got_time) THEN
1017 IF ((nfiles.gt.1).and.(label(1:3).eq.'FRC')) THEN
1018 IF (master) THEN
1019 WRITE (stdout,50) trim(tname(ifield)), 'files:'
1020 DO i=1,nfiles
1021 WRITE (stdout,'(15x,a)') trim(s(i)%name)
1022 END DO
1023 END IF
1024 ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.'BRY')) THEN
1025 IF (master) THEN
1026 WRITE (stdout,50) trim(tname(ifield)), 'files:'
1027 DO i=1,nfiles
1028 WRITE (stdout,'(15x,a)') trim(s(i)%name)
1029 END DO
1030 END IF
1031 ELSE IF ((nfiles.gt.1).and.(label(1:3).eq.'CLM')) THEN
1032 IF (master) THEN
1033 WRITE (stdout,50) trim(tname(ifield)), 'files:'
1034 DO i=1,nfiles
1035 WRITE (stdout,'(15x,a)') trim(s(i)%name)
1036 END DO
1037 END IF
1038 ELSE
1039 lstr=len_trim(ncfile)
1040 IF (master) THEN
1041 WRITE (stdout,50) trim(tname(ifield)), 'file:'
1042 IF (lstr.gt.0) THEN
1043 WRITE (stdout,'(15x,a)') trim(ncfile)
1044 ELSE
1045 WRITE (stdout,'(15x,a,a)') 'file name is blank, ', &
1046 & 'cannot be determined.'
1047 END IF
1048 END IF
1049 END IF
1050 exit_flag=2
1051 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1052 END IF
1053!
1054! If processing model forcing (multiple files allowed), check if file
1055! for requested field (based on pioFile) has been already opened and
1056! get its file descriptor from the association table (FRCdesc).
1057!
1058 IF (label(1:3).eq.'FRC') THEN
1059 DO ifile=1,nfiles
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
1063 piofile=frcdesc(ifile,ng)
1064 EXIT
1065 END IF
1066 END DO
1067 END IF
1068 IF (label(1:3).eq.'BRY') THEN
1069 DO ifile=1,nfiles
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
1073 piofile=brydesc(ifile,ng)
1074 EXIT
1075 END IF
1076 END DO
1077 END IF
1078 IF (label(1:3).eq.'CLM') THEN
1079 DO ifile=1,nfiles
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
1083 piofile=clmdesc(ifile,ng)
1084 EXIT
1085 END IF
1086 END DO
1087 END IF
1088!
1089! If appropriate, open input NetCDF file for reading.
1090! (HGA: it is not needed since file was already open and
1091! pioFile=my_pioFile was assigned above.
1092!
1093 IF (piofile%fh.eq.-1) THEN
1094 CALL pio_netcdf_open (ng, model, ncfile, 0, piofile)
1095 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
1096 IF (master) WRITE (stdout,60) trim(ncfile)
1097 RETURN
1098 END IF
1099 END IF
1100!
1101! If processing model forcing, load ID into association table (FRCids).
1102!
1103 IF (label(1:3).eq.'FRC') THEN
1104 DO ifile=1,nfiles
1105 IF (trim(ncfile).eq.trim(s(ifile)%name)) THEN
1106 frcdesc(ifile,ng)=piofile
1107 s(ifile)%pioFile=piofile
1108 EXIT
1109 END IF
1110 END DO
1111 ELSE IF (label(1:3).eq.'BRY') THEN
1112 DO ifile=1,nfiles
1113 IF (trim(ncfile).eq.trim(s(ifile)%name)) THEN
1114 brydesc(ifile,ng)=piofile
1115 s(ifile)%pioFile=piofile
1116 EXIT
1117 END IF
1118 END DO
1119 ELSE IF (label(1:3).eq.'CLM') THEN
1120 DO ifile=1,nfiles
1121 IF (trim(ncfile).eq.trim(s(ifile)%name)) THEN
1122 clmdesc(ifile,ng)=piofile
1123 s(ifile)%pioFile=piofile
1124 EXIT
1125 END IF
1126 END DO
1127 ELSE
1128 s(1)%pioFile=piofile
1129 END IF
1130 cinfo(ifield,ng)=trim(ncfile)
1131!
1132! The strategy here is to create a local, monotonically increasing
1133! time variable so the interpolation between snapshots is trivial
1134! when cycling forcing fields. Subtract one to time record counter
1135! "Trec" to avoid doing special case at initialization. If the
1136! job task is negative, the logic is backaward in time for the
1137! adjoint model.
1138!
1139 IF (.not.lmulti) THEN
1140 IF (job.lt.0) THEN ! backward time logic, adjoint
1141 IF (irec.eq.1) THEN
1142 tindex=iout
1143 ELSE
1144 tindex=iinfo(8,ifield,ng)
1145 END IF
1146 IF (liocycle) THEN
1147 IF (trec.eq.nrec) THEN
1148 IF (tdays(ng).lt.tmax) THEN
1149 tmono=tend
1150 ELSE
1151 tmono=tend+(tdays(ng)-mod(tdays(ng),clength))
1152 END IF
1153 ELSE
1154 IF (tdays(ng).gt.clength) THEN
1155 tmono=tend+(tdays(ng)-mod(tdays(ng),clength))
1156 ELSE
1157 tmono=tend
1158 END IF
1159 END IF
1160 trec=trec+2
1161 ELSE
1162 tmono=tend
1163 trec=trec+1
1164 CALL pio_netcdf_get_time (ng, model, ncfile, tname(ifield), &
1165 & rclock%DateNumber, tval, &
1166 & piofile = piofile, &
1167 & start = (/trec-1/), &
1168 & total = (/1/))
1169 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
1170 IF (master) WRITE (stdout,70) trim(tname(ifield)), trec
1171 RETURN
1172 END IF
1173 tval=tval*tscale
1174 IF (tval.lt.tend) THEN
1175 trec=trec+1
1176 END IF
1177 END IF
1178 ELSE ! forward time logic
1179 IF (irec.eq.1) THEN
1180 tindex=iout
1181 ELSE
1182 tindex=1
1183 END IF
1184 IF (liocycle) THEN
1185 IF (trec.eq.nrec) THEN
1186 IF (tdays(ng).lt.tmax) THEN
1187 tmono=tstr-clength
1188 ELSE
1189 tmono=tdays(ng)+(tstr-clength)
1190 IF (tstr.eq.tmax) THEN
1191 tmono=tmono+(tmin-mod(tdays(ng)+tmin,clength))
1192 ELSE
1193 tmono=tmono+(tstr-mod(tdays(ng)+tstr,clength))
1194 END IF
1195 END IF
1196 ELSE
1197 IF (tdays(ng).gt.clength) THEN
1198 tmono=tdays(ng)-mod(tdays(ng)-tstr,clength)
1199 ELSE
1200 tmono=tstr
1201 END IF
1202 END IF
1203 ELSE
1204 tmono=tstr
1205 END IF
1206 trec=trec-1
1207 END IF
1208 tmono=tmono*day2sec
1209 iinfo(8,ifield,ng)=tindex
1210 iinfo(9,ifield,ng)=trec
1211 finfo(7,ifield,ng)=tmono
1212 ELSE
1213 iinfo(9,ifield,ng)=trec
1214 END IF
1215!
1216! Set switch for one time record dataset. In this case, the input
1217! data is always the same and time interpolation is not performed.
1218!
1219 IF (nrec.eq.1) lonerec=.true.
1220 linfo(3,ifield,ng)=lonerec
1221 tindex=iinfo(8,ifield,ng)
1222 IF (job.lt.0) THEN
1223 vtime(tindex,ifield,ng)=finfo(4,ifield,ng)
1224 ELSE
1225 vtime(tindex,ifield,ng)=finfo(3,ifield,ng)
1226 END IF
1227!
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, &
1238 & ': ',i0,2x,i0)
1239 50 FORMAT (/,' INQUIRY_PIO - unable to find requested variable: ', &
1240 & a,/,16x,'in ',a)
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)
1244!
1245 RETURN

References mod_iounits::brydesc, mod_ncparam::cinfo, mod_iounits::clmdesc, mod_scalars::day2sec, mod_netcdf::dim_name, mod_netcdf::dim_size, mod_ncparam::dinfo, mod_scalars::exit_flag, mod_ncparam::finfo, strings_mod::founderror(), mod_iounits::frcdesc, mod_ncparam::iinfo, mod_ncparam::linfo, strings_mod::lowercase(), mod_parallel::master, mod_netcdf::n_dim, mod_netcdf::n_vdim, mod_iounits::ncfile, mod_scalars::noerror, mod_pio_netcdf::pio_netcdf_check_dim(), mod_pio_netcdf::pio_netcdf_close(), mod_pio_netcdf::pio_netcdf_inq_var(), mod_pio_netcdf::pio_netcdf_open(), mod_pio_netcdf::pio_tout, mod_pio_netcdf::pio_type, mod_scalars::rclock, mod_iounits::sourcefile, mod_iounits::stdout, mod_scalars::tdays, mod_ncparam::tname, mod_netcdf::var_achar, mod_netcdf::var_afloat, mod_netcdf::var_aname, mod_netcdf::var_dname, mod_netcdf::var_dsize, mod_ncparam::vname, and mod_ncparam::vtime.

Here is the call graph for this function:

The documentation for this interface was generated from the following file: