ROMS
Loading...
Searching...
No Matches
get_ngfldr.F File Reference
#include "cppdefs.h"
Include dependency graph for get_ngfldr.F:

Go to the source code of this file.

Functions/Subroutines

subroutine get_ngfldr (ng, model, ifield, ncid, piofile, nfiles, s, recordless, update, lbi, ubi, ubj, ubk, istr, iend, jrec, fout)
 
subroutine get_ngfldr_nf90 (ng, model, ifield, ncid, nfiles, s, recordless, update, lbi, ubi, ubj, ubk, istr, iend, jrec, fout)
 
subroutine get_ngfldr_pio (ng, model, ifield, piofile, nfiles, s, recordless, update, lbi, ubi, ubj, ubk, istr, iend, jrec, fout)
 

Function/Subroutine Documentation

◆ get_ngfldr()

subroutine get_ngfldr ( integer, intent(in) ng,
integer, intent(in) model,
integer, intent(in) ifield,
integer, intent(inout) ncid,
type (file_desc_t), intent(inout) piofile,
integer, intent(in) nfiles,
type(t_io), dimension(nfiles), intent(inout) s,
logical, intent(in) recordless,
logical, intent(out) update,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) ubj,
integer, intent(in) ubk,
integer, intent(in) istr,
integer, intent(in) iend,
integer, intent(in) jrec,
real(r8), dimension(lbi:ubi,ubj,ubk), intent(inout) fout )

Definition at line 3 of file get_ngfldr.F.

10!
11!git $Id$
12!================================================== Hernan G. Arango ===
13! Copyright (c) 2002-2025 The ROMS Group !
14! Licensed under a MIT/X style license !
15! See License_ROMS.md !
16!=======================================================================
17! !
18! This routine reads in requested non-grided field from specified !
19! NetCDF file. A non-grided field has different dimensions than !
20! model spatial dimensions. Backward time processing. !
21! !
22! On Input: !
23! !
24! ng Nested grid number. !
25! model Calling model identifier. !
26! ifield Field ID. !
27! ncid NetCDF file ID. !
28! nfiles Number of input NetCDF files. !
29! S I/O derived type structure, TYPE(T_IO). !
30! recordless Switch for time invariant field (logical). !
31! fname NetCDF file name(s). !
32! LBi "Fout" 1st dimension lower-bound value. !
33! UBi "Fout" 1st dimension upper-bound value. !
34! UBj "Fout" 2nd dimension upper-bound value, if any. !
35! Otherwise, a value of one is expected. !
36! UBk "Fout" time dimension upper-bound value, if any. !
37! Otherwise, a value of one is expected. !
38! Istr Starting location of read data in the 1st dimension. !
39! Iend Ending location of read data in the 1st dimension; !
40! Number of records read is: Iend-Istr+1. !
41! Jrec Number of records read in the 2st dimension. !
42! !
43! On Output: !
44! !
45! Fout Read field. !
46! update Switch indicating reading of the requested field !
47! the current time step. !
48! !
49!=======================================================================
50!
51 USE mod_param
52 USE mod_parallel
53 USE mod_iounits
54 USE mod_ncparam
55 USE mod_scalars
56!
57 USE strings_mod, ONLY : founderror
58!
59 implicit none
60!
61! Imported variable declarations.
62!
63 logical, intent(in) :: recordless
64 logical, intent(out) :: update
65!
66 integer, intent(in) :: ng, model, ifield, nfiles
67 integer, intent(in) :: LBi, UBi, UBj, UBk, Istr, Iend, Jrec
68 integer, intent(inout) :: ncid
69!
70# if defined PIO_LIB && defined DISTRIBUTE
71 TYPE (File_desc_t), intent(inout) :: pioFile
72# endif
73 TYPE(T_IO), intent(inout) :: S(nfiles)
74!
75 real(r8), intent(inout) :: Fout(LBi:UBi,UBj,UBk)
76!
77! Local variable declarations.
78!
79 character (len=*), parameter :: MyFile = &
80 & __FILE__
81!
82!-----------------------------------------------------------------------
83! Reads in requested non-grided field according to IO type.
84!-----------------------------------------------------------------------
85!
86 SELECT CASE (s(1)%IOtype)
87 CASE (io_nf90)
88 CALL get_ngfldr_nf90 (ng, model, ifield, ncid, &
89 & nfiles, s, recordless, update, &
90 & lbi, ubi, ubj, ubk, istr, iend, jrec, &
91 & fout)
92
93# if defined PIO_LIB && defined DISTRIBUTE
94 CASE (io_pio)
95 CALL get_ngfldr_pio (ng, model, ifield, piofile, &
96 & nfiles, s, recordless, update, &
97 & lbi, ubi, ubj, ubk, istr, iend, jrec, &
98 & fout)
99# endif
100 CASE DEFAULT
101 IF (master) WRITE (stdout,10) s(1)%IOtype
102 exit_flag=2
103 END SELECT
104 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
105!
106 10 FORMAT (' GET_NGFLDR - Illegal input file type, io_type = ',i0, &
107 & /,12x,'Check KeyWord ''INP_LIB'' in ''roms.in''.')
108!
109 RETURN
subroutine get_ngfldr_pio(ng, model, ifield, piofile, nfiles, s, recordless, update, lbi, ubi, ubj, ubk, istr, iend, jrec, fout)
Definition get_ngfldr.F:415
subroutine get_ngfldr_nf90(ng, model, ifield, ncid, nfiles, s, recordless, update, lbi, ubi, ubj, ubk, istr, iend, jrec, fout)
Definition get_ngfldr.F:117
integer stdout
integer, parameter io_nf90
Definition mod_ncparam.F:95
integer, parameter io_pio
Definition mod_ncparam.F:96
logical master
integer exit_flag
integer noerror
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52

References mod_scalars::exit_flag, strings_mod::founderror(), get_ngfldr_nf90(), get_ngfldr_pio(), mod_ncparam::io_nf90, mod_ncparam::io_pio, mod_parallel::master, mod_scalars::noerror, and mod_iounits::stdout.

Referenced by ad_get_data(), and ad_get_idata().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_ngfldr_nf90()

subroutine get_ngfldr_nf90 ( integer, intent(in) ng,
integer, intent(in) model,
integer, intent(in) ifield,
integer, intent(inout) ncid,
integer, intent(in) nfiles,
type(t_io), dimension(nfiles), intent(inout) s,
logical, intent(in) recordless,
logical, intent(out) update,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) ubj,
integer, intent(in) ubk,
integer, intent(in) istr,
integer, intent(in) iend,
integer, intent(in) jrec,
real(r8), dimension(lbi:ubi,ubj,ubk), intent(inout) fout )

Definition at line 113 of file get_ngfldr.F.

117!***********************************************************************
118!
119 USE mod_param
120 USE mod_parallel
121 USE mod_iounits
122 USE mod_ncparam
123 USE mod_netcdf
124 USE mod_scalars
125!
126 USE dateclock_mod, ONLY : time_string
127# ifdef CHECKSUM
128 USE get_hash_mod, ONLY : get_hash
129# endif
130 USE inquiry_mod, ONLY : inquiry
131 USE strings_mod, ONLY : founderror
132!
133 implicit none
134!
135! Imported variable declarations.
136!
137 logical, intent(in) :: recordless
138 logical, intent(out) :: update
139!
140 integer, intent(in) :: ng, model, ifield, nfiles
141 integer, intent(in) :: LBi, UBi, UBj, UBk, Istr, Iend, Jrec
142 integer, intent(inout) :: ncid
143!
144 TYPE(T_IO), intent(inout) :: S(nfiles)
145!
146 real(r8), intent(inout) :: Fout(LBi:UBi,UBj,UBk)
147!
148! Local variable declarations.
149!
150 logical :: Linquire, Liocycle, Lmulti, Lonerec
151!
152 integer :: Nrec, Tid, Tindex, Trec, Vid, Vtype
153 integer :: i, ic, j, job, lend, lstr, npts, nvdim, status
154# ifdef CHECKSUM
155 integer(i8b) :: Fhash
156# endif
157!
158 real(r8) :: Aval, Fmax, Fmin
159
160 real(dp) :: Clength, Tdelta, Tend
161 real(dp) :: Tmax, Tmin, Tmono, Tscale, Tstr
162 real(dp) :: Tsec, Tval
163
164 real(r8), dimension((UBi-LBi+1)*UBj) :: Awrk
165!
166 character (len=22) :: t_code
167
168 character (len=*), parameter :: MyFile = &
169 & __FILE__//", get_ngfldr_nf90"
170!
171 sourcefile=myfile
172!
173!-----------------------------------------------------------------------
174! Initialize.
175!-----------------------------------------------------------------------
176!
177 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
178!
179! Determine if inquiring about the field to process in input NetCDF
180! file(s). This usually happens on first call or when the field
181! time records are split (saved) in several multi-files.
182!
183 linquire=.false.
184 lmulti=.false.
185 IF (iic(ng).eq.0) linquire=.true.
186 IF (.not.linquire.and. &
187 & ((iinfo(10,ifield,ng).gt.1).and. &
188 & (linfo( 5,ifield,ng).or. &
189 & (finfo( 1,ifield,ng)*day2sec.gt.time(ng))))) THEN
190 linquire=.true.
191 lmulti=.true.
192 END IF
193!
194! If appropriate, inquire about the contents of input NetCDF file and
195! fill information arrays.
196!
197! Also, if appropriate, deactivate the Linfo(5,ifield,ng) switch after
198! the query for the LOWER snapshot interpolant from previous multifile
199! in the list. The switch was activated previously to indicate the
200! processing of the FIRST record of the file for the UPPER snapshot
201! interpolant.
202!
203 IF (linquire) THEN
204 job=-1
205 CALL inquiry (ng, model, job, ubk, iend, ubi, ifield, ncid, &
206 & lmulti, nfiles, s)
207 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
208 IF (linfo(5,ifield,ng)) THEN
209 linfo(5,ifield,ng)=.false.
210 END IF
211 IF (linfo(6,ifield,ng)) THEN
212 linfo(6,ifield,ng)=.false.
213 END IF
214 END IF
215!
216!-----------------------------------------------------------------------
217! If appropriate, read in new data.
218!-----------------------------------------------------------------------
219!
220 update=.false.
221 tmono=finfo(7,ifield,ng)
222!
223 IF ((tmono.gt.time(ng)).or.(iic(ng).eq.0).or. &
224 & (iic(ng).eq.ntstart(ng))) THEN
225!
226! Load properties for requested field from information arrays.
227!
228 liocycle=linfo( 2,ifield,ng)
229 lonerec =linfo( 3,ifield,ng)
230 vtype =iinfo( 1,ifield,ng)
231 vid =iinfo( 2,ifield,ng)
232 tid =iinfo( 3,ifield,ng)
233 nrec =iinfo( 4,ifield,ng)
234 tindex =iinfo( 8,ifield,ng)
235 trec =iinfo( 9,ifield,ng)
236 nvdim =iinfo(11,ifield,ng)
237 tmin =finfo( 1,ifield,ng)
238 tmax =finfo( 2,ifield,ng)
239 clength =finfo( 5,ifield,ng)
240 tscale =finfo( 6,ifield,ng)
241 ncfile =cinfo(ifield,ng)
242!
243 IF (liocycle) THEN
244 trec=mod(trec,nrec)-1
245 IF (trec.le.0) trec=nrec+trec
246 ELSE
247 trec=trec-1
248 END IF
249 iinfo(9,ifield,ng)=trec
250!
251 IF ((1.le.trec).and.(trec.le.nrec)) THEN
252!
253! Set rolling index for two-time record storage of input data. If
254! "UBk" is unity, input data is stored in recordless array by the
255! calling program.
256!
257 IF (ubk.eq.1) THEN
258 tindex=1
259 ELSE
260 tindex=3-tindex
261 END IF
262 iinfo(8,ifield,ng)=tindex
263!
264! Read in time coordinate.
265!
266 IF (.not.recordless.and.(tid.ge.0).and.(tid.ne.vid)) THEN
267 CALL netcdf_get_time (ng, model, ncfile, tname(ifield), &
268 & rclock%DateNumber, tval, &
269 & ncid = ncid, &
270 & start = (/trec/), &
271 & total = (/1/))
272 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
273 IF (master) WRITE (stdout,50) trim(tname(ifield)), trec
274 RETURN
275 END IF
276 tval=tval*tscale
277 vtime(tindex,ifield,ng)=tval
278!
279! Activate switch Linfo(5,ifield,ng) if processing the FIRST record of
280! the file for the UPPER time snapshot. We need to get the LOWER time
281! snapshot from previous multifile.
282!
283 IF ((trec.eq.1).and.(tval*day2sec.ge.time(ng))) THEN
284 linfo(5,ifield,ng)=.true.
285 END IF
286 END IF
287!
288! Read in non-grided data. The conditional statement on Jrec is to
289! differentiate between reading a 3D and 2D array.
290!
291 IF (vid.ge.0) THEN
292 fmin=0.0_r8
293 fmax=0.0_r8
294 IF (nvdim.eq.1) THEN
295 npts=iend-istr+1
296 CALL netcdf_get_fvar (ng, model, ncfile, &
297 & vname(1,ifield), awrk, &
298 & ncid = ncid, &
299 & start = (/1/), &
300 & total = (/iend-istr+1/))
301 ELSE IF (nvdim.eq.2) THEN
302 IF (recordless) THEN
303 npts=(iend-istr+1)*jrec
304 CALL netcdf_get_fvar (ng, model, ncfile, &
305 & vname(1,ifield), awrk, &
306 & ncid = ncid, &
307 & start = (/1,1/), &
308 & total = (/iend-istr+1,jrec/))
309 ELSE
310 npts=iend-istr+1
311 CALL netcdf_get_fvar (ng, model, ncfile, &
312 & vname(1,ifield), awrk, &
313 & ncid = ncid, &
314 & start = (/1,trec/), &
315 & total = (/iend-istr+1,1/))
316 END IF
317 ELSE IF (nvdim.eq.3) THEN
318 npts=(iend-istr+1)*jrec
319 CALL netcdf_get_fvar (ng, model, ncfile, &
320 & vname(1,ifield), awrk, &
321 & ncid = ncid, &
322 & start = (/1,1,trec/), &
323 & total = (/iend-istr+1,jrec,1/))
324 END IF
325 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
326 IF (master) WRITE (stdout,50) trim(vname(1,ifield)), trec
327 RETURN
328 END IF
329 fmin=awrk(1)*fscale(ifield,ng)
330 fmax=awrk(1)*fscale(ifield,ng)
331 ic=0
332 DO j=1,jrec
333 DO i=istr,iend
334 ic=ic+1
335 aval=awrk(ic)*fscale(ifield,ng)
336 fmin=min(fmin,aval)
337 fmax=max(fmax,aval)
338 fout(i,j,tindex)=aval
339 END DO
340 END DO
341 finfo(8,ifield,ng)=fmin
342 finfo(9,ifield,ng)=fmax
343#ifdef CHECKSUM
344 CALL get_hash (awrk, npts, fhash)
345#endif
346 IF (master) THEN
347 IF (ubk.eq.1) THEN
348 WRITE (stdout,60) trim(vname(2,ifield)), ng, fmin, fmax
349 ELSE
350 lstr=scan(ncfile,'/',back=.true.)+1
351 lend=len_trim(ncfile)
352 tsec=tval*day2sec
353 CALL time_string (tsec, t_code)
354 WRITE (stdout,70) trim(vname(2,ifield)), t_code, &
355 & ng, trec, tindex, ncfile(lstr:lend), &
356 & tmin, tmax, tval, fmin, fmax
357 END IF
358# ifdef CHECKSUM
359 WRITE (stdout,80) fhash
360# endif
361 END IF
362 update=.true.
363 END IF
364 END IF
365!
366! Decrease the local time variable "Tmono" by the interval between
367! snapshots. If the interval is negative, indicating cycling, add in
368! a cycle length. Load time value (sec) into "Tintrp" which used
369! during interpolation between snapshots.
370!
371 IF (.not.lonerec.and.(.not.recordless)) THEN
372 tdelta=vtime(3-tindex,ifield,ng)-vtime(tindex,ifield,ng)
373 IF (liocycle.and.(tdelta.lt.0.0_r8)) THEN
374 tdelta=tdelta+clength
375 END IF
376 tmono=tmono-tdelta*day2sec
377 finfo(7,ifield,ng)=tmono
378 tintrp(tindex,ifield,ng)=tmono
379 END IF
380 END IF
381!
382 10 FORMAT (/,' GET_NGFLDR_NF90 - unable to find dimension ',a, &
383 & /,19x,'for variable: ',a,/,19x,'in file: ',a, &
384 & /,19x,'file is not CF compliant...')
385 20 FORMAT (/,' GET_NGFLDR_NF90 - too small dimension for variable', &
386 & 1x,a,': ',i0,2x,i0)
387 30 FORMAT (/,' GET_NGFLDR_NF90 - unable to find requested', &
388 & ' variable: ',a,/,19x,'in file: ',a)
389 40 FORMAT (/,' GET_NGFLDR_NF90 - unable to open input NetCDF', &
390 & ' file: ',a)
391 50 FORMAT (/,' GET_NGFLDR_NF90 - error while reading variable: ',a, &
392 & 2x,' at TIME index = ',i0)
393 60 FORMAT (2x,'GET_NGFLDR_NF90 - ',a,/,22x,'(Grid = ',i2.2, &
394 & ', Min = ',1pe15.8,' Max = ',1pe15.8,')')
395 70 FORMAT (2x,'GET_NGFLDR_NF90 - ',a,',',t75,a,/,22x, &
396 & '(Grid=',i2.2,', Rec=',i0,', Index=',i1, &
397 & ', File: ',a,')',/,22x, &
398 & '(Tmin= ', f15.4, ' Tmax= ', f15.4,')', &
399 & t71, 't = ', f15.4 ,/,22x, &
400 & '(Min = ', 1p,e15.8,0p,' Max = ',1p,e15.8,0p,')')
401# ifdef CHECKSUM
402 80 FORMAT (22x,'(CheckSum = ',i0,')')
403# endif
404!
405 RETURN
subroutine, public time_string(mytime, date_string)
Definition dateclock.F:1272
subroutine, public get_hash(a, asize, hash, lreduce)
Definition get_hash.F:72
character(len=256) ncfile
character(len=256) sourcefile
character(len=256), dimension(:,:), allocatable cinfo
logical, dimension(:,:,:), allocatable linfo
real(dp), dimension(:,:,:), allocatable vtime
real(dp), dimension(:,:,:), allocatable tintrp
real(dp), dimension(:,:), allocatable fscale
character(len=46), dimension(0:nv) tname
real(dp), dimension(:,:,:), allocatable finfo
character(len=maxlen), dimension(6, 0:nv) vname
integer, dimension(:,:,:), allocatable iinfo
real(dp), parameter day2sec
integer, dimension(:), allocatable iic
type(t_clock) rclock
real(dp), dimension(:), allocatable time
integer, dimension(:), allocatable ntstart

References mod_ncparam::cinfo, mod_scalars::day2sec, mod_scalars::exit_flag, mod_ncparam::finfo, strings_mod::founderror(), mod_ncparam::fscale, get_hash_mod::get_hash(), mod_scalars::iic, mod_ncparam::iinfo, mod_ncparam::linfo, mod_parallel::master, mod_iounits::ncfile, mod_scalars::noerror, mod_scalars::ntstart, mod_scalars::rclock, mod_iounits::sourcefile, mod_iounits::stdout, mod_scalars::time, dateclock_mod::time_string(), mod_ncparam::tintrp, mod_ncparam::tname, mod_ncparam::vname, and mod_ncparam::vtime.

Referenced by get_ngfldr().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_ngfldr_pio()

subroutine get_ngfldr_pio ( integer, intent(in) ng,
integer, intent(in) model,
integer, intent(in) ifield,
type (file_desc_t), intent(inout) piofile,
integer, intent(in) nfiles,
type(t_io), dimension(nfiles), intent(inout) s,
logical, intent(in) recordless,
logical, intent(out) update,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) ubj,
integer, intent(in) ubk,
integer, intent(in) istr,
integer, intent(in) iend,
integer, intent(in) jrec,
real(r8), dimension(lbi:ubi,ubj,ubk), intent(inout) fout )

Definition at line 411 of file get_ngfldr.F.

415!***********************************************************************
416!
417 USE mod_param
418 USE mod_parallel
419 USE mod_iounits
420 USE mod_ncparam
422 USE mod_scalars
423!
424 USE dateclock_mod, ONLY : time_string
425# ifdef CHECKSUM
426 USE get_hash_mod, ONLY : get_hash
427# endif
428 USE inquiry_mod, ONLY : inquiry
429 USE strings_mod, ONLY : founderror
430!
431 implicit none
432!
433! Imported variable declarations.
434!
435 logical, intent(in) :: recordless
436 logical, intent(out) :: update
437!
438 integer, intent(in) :: ng, model, ifield, nfiles
439 integer, intent(in) :: LBi, UBi, UBj, UBk, Istr, Iend, Jrec
440!
441 TYPE (File_desc_t), intent(inout) :: pioFile
442 TYPE(T_IO), intent(inout) :: S(nfiles)
443!
444 real(r8), intent(inout) :: Fout(LBi:UBi,UBj,UBk)
445!
446! Local variable declarations.
447!
448 logical :: Linquire, Liocycle, Lmulti, Lonerec
449!
450 integer :: Nrec, Tindex, Trec, Vtype
451 integer :: i, ic, j, job, lend, lstr, npts, nvdim, status
452# ifdef CHECKSUM
453 integer(i8b) :: Fhash
454# endif
455!
456 real(r8) :: Aval, Fmax, Fmin
457
458 real(dp) :: Clength, Tdelta, Tend
459 real(dp) :: Tmax, Tmin, Tmono, Tscale, Tstr
460 real(dp) :: Tsec, Tval
461
462 real(r8), dimension((UBi-LBi+1)*UBj) :: Awrk
463!
464 character (len=22) :: t_code
465
466 character (len=*), parameter :: MyFile = &
467 & __FILE__//", get_ngfldr_pio"
468!
469 TYPE (My_VarDesc) :: TpioVar, VpioVar
470!
471 sourcefile=myfile
472!
473!-----------------------------------------------------------------------
474! Initialize.
475!-----------------------------------------------------------------------
476!
477 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
478!
479! Determine if inquiring about the field to process in input NetCDF
480! file(s). This usually happens on first call or when the field
481! time records are split (saved) in several multi-files.
482!
483 linquire=.false.
484 lmulti=.false.
485 IF (iic(ng).eq.0) linquire=.true.
486 IF (.not.linquire.and. &
487 & ((iinfo(10,ifield,ng).gt.1).and. &
488 & (linfo( 5,ifield,ng).or. &
489 & (finfo( 1,ifield,ng)*day2sec.gt.time(ng))))) THEN
490 linquire=.true.
491 lmulti=.true.
492 END IF
493!
494! If appropriate, inquire about the contents of input NetCDF file and
495! fill information arrays.
496!
497! Also, if appropriate, deactivate the Linfo(5,ifield,ng) switch after
498! the query for the LOWER snapshot interpolant from previous multifile
499! in the list. The switch was activated previously to indicate the
500! processing of the FIRST record of the file for the UPPER snapshot
501! interpolant.
502!
503 IF (linquire) THEN
504 job=-1
505 CALL inquiry (ng, model, job, ubk, iend, ubi, ifield, piofile, &
506 & lmulti, nfiles, s)
507 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
508 IF (linfo(5,ifield,ng)) THEN
509 linfo(5,ifield,ng)=.false.
510 END IF
511 IF (linfo(6,ifield,ng)) THEN
512 linfo(6,ifield,ng)=.false.
513 END IF
514 END IF
515!
516!-----------------------------------------------------------------------
517! If appropriate, read in new data.
518!-----------------------------------------------------------------------
519!
520 update=.false.
521 tmono=finfo(7,ifield,ng)
522!
523 IF ((tmono.gt.time(ng)).or.(iic(ng).eq.0).or. &
524 & (iic(ng).eq.ntstart(ng))) THEN
525!
526! Load properties for requested field from information arrays.
527!
528 liocycle=linfo( 2,ifield,ng)
529 lonerec =linfo( 3,ifield,ng)
530 vtype =iinfo( 1,ifield,ng)
531 vpiovar =dinfo( 1,ifield,ng)
532 tpiovar =dinfo( 2,ifield,ng)
533 nrec =iinfo( 4,ifield,ng)
534 tindex =iinfo( 8,ifield,ng)
535 trec =iinfo( 9,ifield,ng)
536 nvdim =iinfo(11,ifield,ng)
537 tmin =finfo( 1,ifield,ng)
538 tmax =finfo( 2,ifield,ng)
539 clength =finfo( 5,ifield,ng)
540 tscale =finfo( 6,ifield,ng)
541 ncfile =cinfo(ifield,ng)
542!
543 IF (liocycle) THEN
544 trec=mod(trec,nrec)-1
545 IF (trec.le.0) trec=nrec+trec
546 ELSE
547 trec=trec-1
548 END IF
549 iinfo(9,ifield,ng)=trec
550!
551 IF ((1.le.trec).and.(trec.le.nrec)) THEN
552!
553! Set rolling index for two-time record storage of input data. If
554! "UBk" is unity, input data is stored in recordless array by the
555! calling program.
556!
557 IF (ubk.eq.1) THEN
558 tindex=1
559 ELSE
560 tindex=3-tindex
561 END IF
562 iinfo(8,ifield,ng)=tindex
563!
564! Read in time coordinate.
565!
566 IF (.not.recordless.and.(tpiovar%vd%varID.ge.0).and. &
567 & (tpiovar%vd%varID.ne.vpiovar%vd%varID)) THEN
568 CALL pio_netcdf_get_time (ng, model, ncfile, tname(ifield), &
569 & rclock%DateNumber, tval, &
570 & piofile = piofile, &
571 & start = (/trec/), &
572 & total = (/1/))
573 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
574 IF (master) WRITE (stdout,50) trim(tname(ifield)), trec
575 RETURN
576 END IF
577 tval=tval*tscale
578 vtime(tindex,ifield,ng)=tval
579!
580! Activate switch Linfo(5,ifield,ng) if processing the FIRST record of
581! the file for the UPPER time snapshot. We need to get the LOWER time
582! snapshot from previous multifile.
583!
584 IF ((trec.eq.1).and.(tval*day2sec.ge.time(ng))) THEN
585 linfo(5,ifield,ng)=.true.
586 END IF
587 END IF
588!
589! Read in non-grided data. The conditional statement on Jrec is to
590! differentiate between reading a 3D and 2D array.
591!
592 IF (vpiovar%vd%varID.ge.0) THEN
593 fmin=0.0_r8
594 fmax=0.0_r8
595 IF (nvdim.eq.1) THEN
596 npts=iend-istr+1
597 CALL pio_netcdf_get_fvar (ng, model, ncfile, &
598 & vname(1,ifield), awrk, &
599 & piofile = piofile, &
600 & start = (/1/), &
601 & total = (/iend-istr+1/))
602 ELSE IF (nvdim.eq.2) THEN
603 IF (recordless) THEN
604 npts=(iend-istr+1)*jrec
605 CALL pio_netcdf_get_fvar (ng, model, ncfile, &
606 & vname(1,ifield), awrk, &
607 & piofile = piofile, &
608 & start = (/1,1/), &
609 & total = (/iend-istr+1,jrec/))
610 ELSE
611 npts=iend-istr+1
612 CALL pio_netcdf_get_fvar (ng, model, ncfile, &
613 & vname(1,ifield), awrk, &
614 & piofile = piofile, &
615 & start = (/1,trec/), &
616 & total = (/iend-istr+1,1/))
617 END IF
618 ELSE IF (nvdim.eq.3) THEN
619 npts=(iend-istr+1)*jrec
620 CALL pio_netcdf_get_fvar (ng, model, ncfile, &
621 & vname(1,ifield), awrk, &
622 & piofile = piofile, &
623 & start = (/1,1,trec/), &
624 & total = (/iend-istr+1,jrec,1/))
625 END IF
626 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
627 IF (master) WRITE (stdout,50) trim(vname(1,ifield)), trec
628 RETURN
629 END IF
630 fmin=awrk(1)*fscale(ifield,ng)
631 fmax=awrk(1)*fscale(ifield,ng)
632 ic=0
633 DO j=1,jrec
634 DO i=istr,iend
635 ic=ic+1
636 aval=awrk(ic)*fscale(ifield,ng)
637 fmin=min(fmin,aval)
638 fmax=max(fmax,aval)
639 fout(i,j,tindex)=aval
640 END DO
641 END DO
642 finfo(8,ifield,ng)=fmin
643 finfo(9,ifield,ng)=fmax
644# ifdef CHECKSUM
645 CALL get_hash (awrk, npts, fhash)
646# endif
647 IF (master) THEN
648 IF (ubk.eq.1) THEN
649 WRITE (stdout,60) trim(vname(2,ifield)), ng, fmin, fmax
650 ELSE
651 lstr=scan(ncfile,'/',back=.true.)+1
652 lend=len_trim(ncfile)
653 tsec=tval*day2sec
654 CALL time_string (tsec, t_code)
655 WRITE (stdout,70) trim(vname(2,ifield)), t_code, &
656 & ng, trec, tindex, ncfile(lstr:lend), &
657 & tmin, tmax, tval, fmin, fmax
658 END IF
659# ifdef CHECKSUM
660 WRITE (stdout,80) fhash
661# endif
662 END IF
663 update=.true.
664 END IF
665 END IF
666!
667! Decrease the local time variable "Tmono" by the interval between
668! snapshots. If the interval is negative, indicating cycling, add in
669! a cycle length. Load time value (sec) into "Tintrp" which used
670! during interpolation between snapshots.
671!
672 IF (.not.lonerec.and.(.not.recordless)) THEN
673 tdelta=vtime(3-tindex,ifield,ng)-vtime(tindex,ifield,ng)
674 IF (liocycle.and.(tdelta.lt.0.0_r8)) THEN
675 tdelta=tdelta+clength
676 END IF
677 tmono=tmono-tdelta*day2sec
678 finfo(7,ifield,ng)=tmono
679 tintrp(tindex,ifield,ng)=tmono
680 END IF
681 END IF
682!
683 10 FORMAT (/,' GET_NGFLDR_PIO - unable to find dimension ',a, &
684 & /,18x,'for variable: ',a,/,18x,'in file: ',a, &
685 & /,18x,'file is not CF compliant...')
686 20 FORMAT (/,' GET_NGFLDR_PIO - too small dimension for variable', &
687 & 1x,a,': ',i0,2x,i0)
688 30 FORMAT (/,' GET_NGFLDR_PIO - unable to find requested', &
689 & ' variable: ',a,/,18x,'in file: ',a)
690 40 FORMAT (/,' GET_NGFLDR_PIO - unable to open input NetCDF', &
691 & ' file: ',a)
692 50 FORMAT (/,' GET_NGFLDR_PIO - error while reading variable: ',a, &
693 & 2x,' at TIME index = ',i0)
694 60 FORMAT (2x,'GET_NGFLDR_PIO - ',a,/,22x,'(Grid = ',i2.2, &
695 & ', Min = ',1pe15.8,' Max = ',1pe15.8,')')
696 70 FORMAT (2x,'GET_NGFLDR_PIO - ',a,',',t75,a,/,22x, &
697 & '(Grid=',i2.2,', Rec=',i0,', Index=',i1, &
698 & ', File: ',a,')',/,22x, &
699 & '(Tmin= ', f15.4, ' Tmax= ', f15.4,')', &
700 & t71, 't = ', f15.4 ,/,22x, &
701 & '(Min = ', 1p,e15.8,0p,' Max = ',1p,e15.8,0p,')')
702# ifdef CHECKSUM
703 80 FORMAT (22x,'(CheckSum = ',i0,')')
704# endif
705!
706 RETURN
type(my_vardesc), dimension(:,:,:), pointer dinfo

References mod_ncparam::cinfo, mod_scalars::day2sec, mod_ncparam::dinfo, mod_scalars::exit_flag, mod_ncparam::finfo, strings_mod::founderror(), mod_ncparam::fscale, get_hash_mod::get_hash(), mod_scalars::iic, mod_ncparam::iinfo, mod_ncparam::linfo, mod_parallel::master, mod_iounits::ncfile, mod_scalars::noerror, mod_scalars::ntstart, mod_scalars::rclock, mod_iounits::sourcefile, mod_iounits::stdout, mod_scalars::time, dateclock_mod::time_string(), mod_ncparam::tintrp, mod_ncparam::tname, mod_ncparam::vname, and mod_ncparam::vtime.

Referenced by get_ngfldr().

Here is the call graph for this function:
Here is the caller graph for this function: