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

Go to the source code of this file.

Functions/Subroutines

subroutine get_2dfld (ng, model, ifield, ncid, piofile, nfiles, s, update, lbi, ubi, lbj, ubj, iout, irec, fmask, fout)
 
subroutine get_2dfld_nf90 (ng, model, ifield, ncid, nfiles, s, update, lbi, ubi, lbj, ubj, iout, irec, fmask, fout)
 
subroutine get_2dfld_pio (ng, model, ifield, piofile, nfiles, s, update, lbi, ubi, lbj, ubj, iout, irec, fmask, fout)
 

Function/Subroutine Documentation

◆ get_2dfld()

subroutine get_2dfld ( 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(out) update,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) iout,
integer, intent(in) irec,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) fmask,
real(r8), dimension(lbi:ubi,lbj:ubj,iout), intent(inout) fout )

Definition at line 2 of file get_2dfld.F.

12!
13!git $Id$
14!================================================== Hernan G. Arango ===
15! Copyright (c) 2002-2025 The ROMS Group !
16! Licensed under a MIT/X style license !
17! See License_ROMS.md !
18!=======================================================================
19! !
20! This routine reads in requested 2D field (point or grided) from !
21! specified NetCDF file. Forward time processing. !
22! !
23! On Input: !
24! !
25! ng Nested grid number. !
26! model Calling model identifier. !
27! ifield Field ID. !
28! ncid NetCDF file ID. !
29#if defined PIO_LIB && defined DISTRIBUTE
30! pioFile PIO file descriptor structure, TYPE(file_desc_t) !
31! pioFile%fh file handler !
32! pioFile%iosystem IO system descriptor (struct) !
33#endif
34! nfiles Number of input NetCDF files. !
35! S I/O derived type structure, TYPE(T_IO). !
36! LBi I-dimension Lower bound. !
37! UBi I-dimension Upper bound. !
38! LBj J-dimension Lower bound. !
39! UBj J-dimension Upper bound. !
40! Iout Size of the outer dimension, if any. Otherwise, !
41! Iout must be set to one by the calling program. !
42! Irec Number of 2D field records to read. !
43! Fmask Land/Sea mask, if any. !
44! !
45! On Output: !
46! !
47! Fout Read field. !
48! update Switch indicating reading of the requested field !
49! the current time step. !
50! !
51!=======================================================================
52!
53 USE mod_param
54 USE mod_parallel
55 USE mod_iounits
56 USE mod_ncparam
57 USE mod_scalars
58!
59 USE strings_mod, ONLY : founderror
60!
61 implicit none
62!
63! Imported variable declarations.
64!
65 logical, intent(out) :: update
66!
67 integer, intent(in) :: ng, model, ifield, nfiles, Iout, Irec
68 integer, intent(in) :: LBi, UBi, LBj, UBj
69 integer, intent(inout) :: ncid
70!
71#if defined PIO_LIB && defined DISTRIBUTE
72 TYPE (File_desc_t), intent(inout) :: pioFile
73#endif
74 TYPE(T_IO), intent(inout) :: S(nfiles)
75!
76#ifdef MASKING
77 real(r8), intent(in) :: Fmask(LBi:UBi,LBj:UBj)
78#endif
79 real(r8), intent(inout) :: Fout(LBi:UBi,LBj:UBj,Iout)
80!
81! Local variable declarations.
82!
83 character (len=*), parameter :: MyFile = &
84 & __FILE__
85!
86!-----------------------------------------------------------------------
87! Read in requested 2D field according to IO type.
88!-----------------------------------------------------------------------
89!
90 SELECT CASE (s(1)%IOtype)
91 CASE (io_nf90)
92 CALL get_2dfld_nf90 (ng, model, ifield, ncid, nfiles, s, &
93 & update, lbi, ubi, lbj, ubj, iout, irec, &
94#ifdef MASKING
95 & fmask, &
96#endif
97 & fout)
98
99#if defined PIO_LIB && defined DISTRIBUTE
100 CASE (io_pio)
101 CALL get_2dfld_pio (ng, model, ifield, piofile, nfiles, s, &
102 & update, lbi, ubi, lbj, ubj, iout, irec, &
103# ifdef MASKING
104 & fmask, &
105# endif
106 & fout)
107#endif
108 CASE DEFAULT
109 IF (master) WRITE (stdout,10) s(1)%IOtype
110 exit_flag=2
111 END SELECT
112 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
113!
114 10 FORMAT (' GET_2DFLD - Illegal input file type, io_type = ',i0, &
115 & /,13x,'Check KeyWord ''INP_LIB'' in ''roms.in''.')
116!
117 RETURN
subroutine get_2dfld_nf90(ng, model, ifield, ncid, nfiles, s, update, lbi, ubi, lbj, ubj, iout, irec, fmask, fout)
Definition get_2dfld.F:128
subroutine get_2dfld_pio(ng, model, ifield, piofile, nfiles, s, update, lbi, ubi, lbj, ubj, iout, irec, fmask, fout)
Definition get_2dfld.F:438
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_2dfld_nf90(), get_2dfld_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_2dfld_nf90()

subroutine get_2dfld_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(out) update,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) iout,
integer, intent(in) irec,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) fmask,
real(r8), dimension(lbi:ubi,lbj:ubj,iout), intent(inout) fout )

Definition at line 121 of file get_2dfld.F.

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

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

◆ get_2dfld_pio()

subroutine get_2dfld_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(out) update,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) iout,
integer, intent(in) irec,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) fmask,
real(r8), dimension(lbi:ubi,lbj:ubj,iout), intent(inout) fout )

Definition at line 431 of file get_2dfld.F.

438!***********************************************************************
439!
440 USE mod_param
441 USE mod_parallel
442 USE mod_iounits
443 USE mod_ncparam
445 USE mod_scalars
446!
447 USE dateclock_mod, ONLY : time_string
448 USE inquiry_mod, ONLY : inquiry
449 USE nf_fread2d_mod, ONLY : nf_fread2d
450 USE nf_fread3d_mod, ONLY : nf_fread3d
451 USE strings_mod, ONLY : founderror
452!
453 implicit none
454!
455! Imported variable declarations.
456!
457 logical, intent(out) :: update
458!
459 integer, intent(in) :: ng, model, ifield, nfiles, Iout, Irec
460 integer, intent(in) :: LBi, UBi, LBj, UBj
461!
462 TYPE (File_desc_t), intent(inout) :: pioFile
463 TYPE(T_IO), intent(inout) :: S(nfiles)
464!
465# ifdef MASKING
466 real(r8), intent(in) :: Fmask(LBi:UBi,LBj:UBj)
467# endif
468 real(r8), intent(inout) :: Fout(LBi:UBi,LBj:UBj,Iout)
469!
470! Local variable declarations.
471!
472 logical :: Lgridded, Linquire, Liocycle, Lmulti, Lonerec, Lregrid
473 logical :: special
474!
475 integer :: Nrec, Tindex, Trec, Vtype
476 integer :: i, job, lend, lstr, lvar, status
477 integer :: Vsize(4)
478# ifdef CHECKSUM
479 integer(i8b) :: Fhash
480# endif
481!
482 real(r8) :: Fmax, Fmin, Fval
483
484 real(dp) :: Clength, Tdelta, Tend
485 real(dp) :: Tmax, Tmin, Tmono, Tscale, Tstr
486 real(dp) :: Tsec, Tval
487!
488 character (len= 1) :: Rswitch
489 character (len=22) :: t_code
490
491 character (len=*), parameter :: MyFile = &
492 & __FILE__//", get_2dfld_pio"
493!
494 TYPE (IO_Desc_t), pointer :: ioDesc
495 TYPE (My_VarDesc) :: TpioVar, VpioVar
496!
497 sourcefile=myfile
498!
499!-----------------------------------------------------------------------
500! Initialize.
501!-----------------------------------------------------------------------
502!
503 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
504!
505! Determine if inquiring about the field to process in input NetCDF
506! file(s). This usually happens on first call or when the field
507! time records are split (saved) in several multi-files.
508!
509 linquire=.false.
510 lmulti=.false.
511 IF (iic(ng).eq.0) linquire=.true.
512 IF (.not.linquire.and. &
513 & ((iinfo(10,ifield,ng).gt.1).and. &
514 & (linfo( 6,ifield,ng).or. &
515 & (finfo( 2,ifield,ng)*day2sec.lt.time(ng))))) THEN
516 linquire=.true.
517 lmulti=.true.
518 END IF
519!
520! If appropriate, inquire about the contents of input NetCDF file and
521! fill information arrays.
522!
523! Also, if appropriate, deactivate the Linfo(6,ifield,ng) switch after
524! the query for the UPPER snapshot interpolant from previous multifile
525! in the list. The switch was activated previously to indicate the
526! processing of the FIRST record of the file for the LOWER snapshot
527! interpolant.
528!
529 IF (linquire) THEN
530 job=2
531 CALL inquiry (ng, model, job, iout, irec, 1, ifield, piofile, &
532 & lmulti, nfiles, s)
533 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
534 IF (linfo(6,ifield,ng)) THEN
535 linfo(6,ifield,ng)=.false.
536 END IF
537 END IF
538!
539!-----------------------------------------------------------------------
540! If appropriate, read in new data.
541!-----------------------------------------------------------------------
542!
543 update=.false.
544 lregrid=.false.
545 tmono=finfo(7,ifield,ng)
546!
547 IF ((tmono.lt.time(ng)).or.(iic(ng).eq.0).or. &
548 & (iic(ng).eq.ntstart(ng))) THEN
549!
550! Load properties for requested field from information arrays.
551!
552 lgridded=linfo(1,ifield,ng)
553 liocycle=linfo(2,ifield,ng)
554 lonerec =linfo(3,ifield,ng)
555 special =linfo(4,ifield,ng)
556 vtype =iinfo(1,ifield,ng)
557 vpiovar =dinfo(1,ifield,ng)
558 tpiovar =dinfo(2,ifield,ng)
559 nrec =iinfo(4,ifield,ng)
560 vsize(1)=iinfo(5,ifield,ng)
561 vsize(2)=iinfo(6,ifield,ng)
562 tindex =iinfo(8,ifield,ng)
563 trec =iinfo(9,ifield,ng)
564 tmin =finfo(1,ifield,ng)
565 tmax =finfo(2,ifield,ng)
566 clength =finfo(5,ifield,ng)
567 tscale =finfo(6,ifield,ng)
568 ncfile =cinfo(ifield,ng)
569!
570 IF (liocycle) THEN
571 trec=mod(trec,nrec)+1
572 ELSE
573 trec=trec+1
574 END IF
575 iinfo(9,ifield,ng)=trec
576!
577 IF (trec.le.nrec) THEN
578!
579! Set rolling index for two-time record storage of input data. If
580! "Iout" is unity, input data is stored in timeless array by the
581! calling program. If Irec > 1, this routine is used to read a 2D
582! field varying in another non-time dimension like tidal constituent
583! component.
584!
585 IF (.not.special.and.(irec.eq.1)) THEN
586 IF (iout.eq.1) THEN
587 tindex=1
588 ELSE
589 tindex=3-tindex
590 END IF
591 iinfo(8,ifield,ng)=tindex
592 END IF
593!
594! Read in time coordinate and scale it to day units.
595!
596 IF (.not.special.and.(tpiovar%vd%varID.ge.0)) THEN
597 CALL pio_netcdf_get_time (ng, model, ncfile, tname(ifield), &
598 & rclock%DateNumber, tval, &
599 & piofile = piofile, &
600 & start = (/trec/), &
601 & total = (/1/))
602 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
603 IF (master) WRITE (stdout,40) trim(tname(ifield)), trec
604 RETURN
605 END IF
606 tval=tval*tscale
607 vtime(tindex,ifield,ng)=tval
608!
609! Activate switch Linfo(6,ifield,ng) if processing the LAST record of
610! the file for the LOWER time snapshot. We need to get the UPPER time
611! snapshot from NEXT multifile.
612!
613 IF ((trec.eq.nrec).and.(tval*day2sec.le.time(ng))) THEN
614 linfo(6,ifield,ng)=.true.
615 END IF
616 END IF
617!
618! Read in 2D-grided or point data. Notice for special 2D fields, Vtype
619! is augmented by four indicating reading a 3D field (tidal data with
620! tidal constituents in the third dimension. Currently, all input tidal
621! constituents data is at RHO-points.
622!
623 IF (special) vtype=vtype+4
624!
625 IF (vpiovar%vd%varID.ge.0) THEN
626 SELECT CASE (vtype)
627 CASE (p2dvar)
628 IF (kind(fout).eq.8) THEN
629 iodesc => iodesc_dp_p2dvar(ng)
630 ELSE
631 iodesc => iodesc_sp_p2dvar(ng)
632 END IF
633 CASE (r2dvar)
634 IF (kind(fout).eq.8) THEN
635 iodesc => iodesc_dp_r2dvar(ng)
636 ELSE
637 iodesc => iodesc_sp_r2dvar(ng)
638 END IF
639 CASE (u2dvar)
640 IF (kind(fout).eq.8) THEN
641 iodesc => iodesc_dp_u2dvar(ng)
642 ELSE
643 iodesc => iodesc_sp_u2dvar(ng)
644 END IF
645 CASE (v2dvar)
646 IF (kind(fout).eq.8) THEN
647 iodesc => iodesc_dp_v2dvar(ng)
648 ELSE
649 iodesc => iodesc_sp_v2dvar(ng)
650 END IF
651# if defined SSH_TIDES || defined UV_TIDES
652 CASE (r2dvar+4)
653 IF (kind(fout).eq.8) THEN
654 iodesc => iodesc_dp_rtides(ng)
655 ELSE
656 iodesc => iodesc_sp_rtides(ng)
657 END IF
658 vpiovar%gtype=r2dvar+4
659# endif
660 END SELECT
661!
662 IF (lgridded) THEN
663 IF (special) THEN
664 vsize(3)=irec
665 status=nf_fread3d(ng, model, ncfile, piofile, &
666 & vname(1,ifield), vpiovar, &
667 & 0, iodesc, vsize, &
668 & lbi, ubi, lbj, ubj, 1, irec, &
669 & fscale(ifield,ng), fmin, fmax, &
670# ifdef MASKING
671 & fmask, &
672# endif
673# ifdef CHECKSUM
674 & fout, &
675 & checksum = fhash)
676# else
677 & fout)
678# endif
679 ELSE
680 status=nf_fread2d(ng, model, ncfile, piofile, &
681 & vname(1,ifield), vpiovar, &
682 & trec, iodesc, vsize, &
683 & lbi, ubi, lbj, ubj, &
684 & fscale(ifield,ng), fmin, fmax, &
685# ifdef MASKING
686 & fmask, &
687# endif
688 & fout(:,:,tindex), &
689# ifdef CHECKSUM
690 & checksum = fhash, &
691# endif
692 & lregrid = lregrid)
693
694 END IF
695 ELSE
696 CALL pio_netcdf_get_fvar (ng, model, ncfile, &
697 & vname(1,ifield), fval, &
698 & piofile = piofile, &
699 & start = (/trec/), &
700 & total = (/1/))
701 fval=fval*fscale(ifield,ng)
702 fpoint(tindex,ifield,ng)=fval
703 fmin=fval
704 fmax=fval
705 END IF
706 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
707 IF (master) WRITE (stdout,40) trim(vname(1,ifield)), trec
708 RETURN
709 END IF
710 finfo(8,ifield,ng)=fmin
711 finfo(9,ifield,ng)=fmax
712 IF (master) THEN
713 IF (special) THEN
714 WRITE (stdout,50) trim(vname(2,ifield)), ng, fmin, fmax
715 ELSE
716 lstr=scan(ncfile,'/',back=.true.)+1
717 lend=len_trim(ncfile)
718 lvar=min(46,len_trim(vname(2,ifield)))
719 tsec=tval*day2sec
720 CALL time_string (tsec, t_code)
721 IF (lregrid) THEN
722 rswitch='T'
723 ELSE
724 rswitch='F'
725 END IF
726 WRITE (stdout,60) vname(2,ifield)(1:lvar), t_code, &
727 & ng, trec, tindex, ncfile(lstr:lend), &
728 & tmin, tmax, tval, fmin, fmax, rswitch
729 END IF
730# ifdef CHECKSUM
731 WRITE (stdout,70) fhash
732# endif
733 END IF
734 update=.true.
735 END IF
736 END IF
737!
738! Increment the local time variable "Tmono" by the interval between
739! snapshots. If the interval is negative, indicating cycling, add in
740! a cycle length. Load time value (sec) into "Tintrp" which used
741! during interpolation between snapshots.
742!
743 IF (.not.lonerec.and.(.not.special)) THEN
744 tdelta=vtime(tindex,ifield,ng)-vtime(3-tindex,ifield,ng)
745 IF (liocycle.and.(tdelta.lt.0.0_r8)) THEN
746 tdelta=tdelta+clength
747 END IF
748 tmono=tmono+tdelta*day2sec
749 finfo(7,ifield,ng)=tmono
750 tintrp(tindex,ifield,ng)=tmono
751 END IF
752 END IF
753!
754 10 FORMAT (/,' GET_2DFLD_PIO - unable to find dimension ',a, &
755 & /,17x,'for variable: ',a,/,17x,'in file: ',a, &
756 & /,17x,'file is not CF compliant...')
757 20 FORMAT (/,' GET_2DFLD_PIO - unable to find requested variable:', &
758 & 1x,a,/,17x,'in ',a)
759 30 FORMAT (/,' GET_2DFLD_PIO - unable to open input NetCDF', &
760 & ' file: ',a)
761 40 FORMAT (/,' GET_2DFLD_PIO - error while reading variable: ',a, &
762 & 2x,' at TIME index = ',i0)
763 50 FORMAT (2x,'GET_2DFLD_PIO - ',a,/,22x,'(Grid = ',i2.2, &
764 & ', Min = ',1p,e15.8,0p,' Max = ',1p,e15.8,0p,')')
765 60 FORMAT (2x,'GET_2DFLD_PIO - ',a,',',t75,a,/,22x, &
766 & '(Grid=',i2.2,', Rec=',i0,', Index=',i1, &
767 & ', File: ',a,')',/,22x, &
768 & '(Tmin= ', f15.4, ' Tmax= ', f15.4,')', &
769 & t71, 't = ', f15.4 ,/,22x, &
770 & '(Min = ', 1p,e15.8,0p,' Max = ',1p,e15.8,0p,')', &
771 & t71, 'regrid = ',a)
772# ifdef CHECKSUM
773 70 FORMAT (22x,'(CheckSum = ',i0,')')
774# endif
775!
776 RETURN
type(my_vardesc), dimension(:,:,:), pointer dinfo
integer, parameter u2dvar
Definition mod_param.F:718
integer, parameter p2dvar
Definition mod_param.F:716
integer, parameter r2dvar
Definition mod_param.F:717
integer, parameter v2dvar
Definition mod_param.F:719
type(io_desc_t), dimension(:), pointer iodesc_sp_u2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_r2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_rtides
type(io_desc_t), dimension(:), pointer iodesc_dp_u2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_p2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_rtides
type(io_desc_t), dimension(:), pointer iodesc_dp_v2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_r2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_p2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_v2dvar

References mod_ncparam::cinfo, mod_scalars::day2sec, mod_ncparam::dinfo, mod_scalars::exit_flag, mod_ncparam::finfo, strings_mod::founderror(), mod_ncparam::fpoint, mod_ncparam::fscale, mod_scalars::iic, mod_ncparam::iinfo, mod_pio_netcdf::iodesc_dp_p2dvar, mod_pio_netcdf::iodesc_dp_r2dvar, mod_pio_netcdf::iodesc_dp_rtides, mod_pio_netcdf::iodesc_dp_u2dvar, mod_pio_netcdf::iodesc_dp_v2dvar, mod_pio_netcdf::iodesc_sp_p2dvar, mod_pio_netcdf::iodesc_sp_r2dvar, mod_pio_netcdf::iodesc_sp_rtides, mod_pio_netcdf::iodesc_sp_u2dvar, mod_pio_netcdf::iodesc_sp_v2dvar, mod_ncparam::linfo, mod_parallel::master, mod_iounits::ncfile, mod_scalars::noerror, mod_scalars::ntstart, mod_param::p2dvar, mod_param::r2dvar, mod_scalars::rclock, mod_iounits::sourcefile, mod_iounits::stdout, mod_scalars::time, dateclock_mod::time_string(), mod_ncparam::tintrp, mod_ncparam::tname, mod_param::u2dvar, mod_param::v2dvar, mod_ncparam::vname, and mod_ncparam::vtime.

Referenced by get_2dfld().

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