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

Go to the source code of this file.

Functions/Subroutines

subroutine get_ngfld (ng, model, ifield, ncid, piofile, nfiles, s, recordless, update, lbi, ubi, ubj, ubk, istr, iend, jrec, fout)
 
subroutine get_ngfld_nf90 (ng, model, ifield, ncid, nfiles, s, recordless, update, lbi, ubi, ubj, ubk, istr, iend, jrec, fout)
 
subroutine get_ngfld_pio (ng, model, ifield, piofile, nfiles, s, recordless, update, lbi, ubi, ubj, ubk, istr, iend, jrec, fout)
 

Function/Subroutine Documentation

◆ get_ngfld()

subroutine get_ngfld ( integer, intent(in) ng,
integer, intent(in) model,
integer, intent(in) ifield,
integer 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 2 of file get_ngfld.F.

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

Referenced by get_data(), get_idata(), rp_get_data(), rp_get_idata(), tl_get_data(), and tl_get_idata().

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

◆ get_ngfld_nf90()

subroutine get_ngfld_nf90 ( integer, intent(in) ng,
integer, intent(in) model,
integer, intent(in) ifield,
integer 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 116 of file get_ngfld.F.

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

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

◆ get_ngfld_pio()

subroutine get_ngfld_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 410 of file get_ngfld.F.

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

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