ROMS
Loading...
Searching...
No Matches
nf_fread3d_bry_mod::nf_fread3d_bry Interface Reference

Public Member Functions

integer function nf90_fread3d_bry (ng, model, ncname, ncid, ncvname, ncvarid, tindex, gtype, lbij, ubij, lbk, ubk, nrec, ascl, amin, amax, abry, checksum)
 
integer function pio_fread3d_bry (ng, model, ncname, piofile, ncvname, piovar, tindex, piodesc, lbij, ubij, lbk, ubk, nrec, ascl, amin, amax, abry, checksum)
 

Detailed Description

Definition at line 71 of file nf_fread3d_bry.F.

Member Function/Subroutine Documentation

◆ nf90_fread3d_bry()

integer function nf_fread3d_bry_mod::nf_fread3d_bry::nf90_fread3d_bry ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
integer, intent(in) ncid,
character (len=*), intent(in) ncvname,
integer, intent(in) ncvarid,
integer, intent(in) tindex,
integer, intent(in) gtype,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) lbk,
integer, intent(in) ubk,
integer, intent(in) nrec,
real(dp), intent(in) ascl,
real(r8), intent(out) amin,
real(r8), intent(out) amax,
real(r8), dimension(lbij:,:,:,:), intent(out) abry,
integer(i8b), intent(out), optional checksum )

Definition at line 81 of file nf_fread3d_bry.F.

87!***********************************************************************
88!
89 USE mod_netcdf
90!
91! Imported variable declarations.
92!
93 integer, intent(in) :: ng, model, ncid, ncvarid, tindex, gtype
94 integer, intent(in) :: LBij, UBij, LBk, UBk, Nrec
95!
96 integer(i8b), intent(out), optional :: checksum
97!
98 real(dp), intent(in) :: Ascl
99 real(r8), intent(out) :: Amin
100 real(r8), intent(out) :: Amax
101!
102 character (len=*), intent(in) :: ncname
103 character (len=*), intent(in) :: ncvname
104!
105#ifdef ASSUMED_SHAPE
106 real(r8), intent(out) :: Abry(LBij:,:,:,:)
107#else
108 real(r8), intent(out) :: Abry(LBij:UBij,LBk:UBk,4,Nrec)
109#endif
110!
111! Local variable declarations.
112!
113 logical :: Lchecksum
114 logical, dimension(3) :: foundit
115 logical, dimension(4) :: bounded
116!
117 integer :: ghost, i, ib, ij, ir, j, k, tile
118 integer :: IorJ, IJKlen, Imin, Imax, Jmin, Jmax, Klen, Npts
119 integer :: Cgrid, Istr, Iend, Jstr, Jend
120 integer, dimension(5) :: start, total
121
122 integer :: status
123!
124 real(r8) :: Afactor, Aoffset, Aspval
125!
126 real(r8), allocatable :: Cwrk(:) ! used for checksum
127
128 real(r8), dimension(3) :: AttValue
129
130#if !defined PARALLEL_IO && defined DISTRIBUTE
131 real(r8), dimension(3) :: rbuffer
132#endif
133 real(r8), dimension(LBij:UBij,LBk:UBk,4,Nrec) :: wrk
134!
135 character (len=12), dimension(3) :: AttName
136
137 character (len=*), parameter :: MyFile = &
138 & __FILE__//", nf90_fread3d_bry"
139!
140!-----------------------------------------------------------------------
141! Set starting and ending indices to process.
142!-----------------------------------------------------------------------
143!
144 status=nf90_noerr
145!
146! Set first and last grid point according to staggered C-grid
147! classification.
148!
149! Notice that (Imin,Jmin) and (Imax,Jmax) are the corner of the
150! computational tile. If ghost=0, ghost points are not processed.
151! They will be processed elsewhere by the appropriate call to any
152! of the routines in "mp_exchange.F". If ghost=1, the ghost points
153! are read.
154!
155 IF (model.eq.iadm) THEN
156 ghost=0 ! non-overlapping, no ghost points
157 ELSE
158 ghost=1 ! overlapping, read ghost points
159 END IF
160
161 SELECT CASE (gtype)
162 CASE (p2dvar, p3dvar)
163 cgrid=1
164 CASE (r2dvar, r3dvar)
165 cgrid=2
166 CASE (u2dvar, u3dvar)
167 cgrid=3
168 CASE (v2dvar, v3dvar)
169 cgrid=4
170 CASE DEFAULT
171 cgrid=2
172 END SELECT
173
174#ifdef DISTRIBUTE
175 tile=myrank
176 imin=bounds(ng)%Imin(cgrid,ghost,tile)
177 imax=bounds(ng)%Imax(cgrid,ghost,tile)
178 jmin=bounds(ng)%Jmin(cgrid,ghost,tile)
179 jmax=bounds(ng)%Jmax(cgrid,ghost,tile)
180#else
181 tile=-1
182 imin=lbij
183 imax=ubij
184 jmin=lbij
185 jmax=ubij
186#endif
187 iorj=iobounds(ng)%IorJ
188 klen=ubk-lbk+1
189 ijklen=iorj*klen
190 npts=ijklen*4*nrec
191!
192! Get tile bounds.
193!
194 istr=bounds(ng)%Istr(tile)
195 iend=bounds(ng)%Iend(tile)
196 jstr=bounds(ng)%Jstr(tile)
197 jend=bounds(ng)%Jend(tile)
198!
199! Set switch to process boundary data by their associated tiles.
200!
201 bounded(iwest )=domain(ng)%Western_Edge (tile)
202 bounded(ieast )=domain(ng)%Eastern_Edge (tile)
203 bounded(isouth)=domain(ng)%Southern_Edge(tile)
204 bounded(inorth)=domain(ng)%Northern_Edge(tile)
205!
206! Set NetCDF dimension counters for processing requested field.
207!
208 start(1)=1
209 total(1)=iorj
210 start(2)=1
211 total(2)=klen
212 start(3)=1
213 total(3)=4
214 start(4)=1
215 total(4)=nrec
216 start(5)=tindex
217 total(5)=1
218!
219! Check if the following attributes: "scale_factor", "add_offset", and
220! "_FillValue" are present in the input NetCDF variable:
221!
222! If the "scale_value" attribute is present, the data is multiplied by
223! this factor after reading.
224! If the "add_offset" attribute is present, this value is added to the
225! data after reading.
226! If both "scale_factor" and "add_offset" attributes are present, the
227! data are first scaled before the offset is added.
228! If the "_FillValue" attribute is present, the data having this value
229! is treated a missing and it is replaced with zero. This feature it is
230! usually related with the land/sea masking.
231!
232 attname(1)='scale_factor'
233 attname(2)='add_offset '
234 attname(3)='_FillValue '
235
236 CALL netcdf_get_fatt (ng, model, ncname, ncvarid, attname, &
237 & attvalue, foundit, &
238 & ncid = ncid)
239 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
240 status=ioerror
241 RETURN
242 END IF
243
244 IF (.not.foundit(1)) THEN
245 afactor=1.0_r8
246 ELSE
247 afactor=attvalue(1)
248 END IF
249
250 IF (.not.foundit(2)) THEN
251 aoffset=0.0_r8
252 ELSE
253 aoffset=attvalue(2)
254 END IF
255
256 IF (.not.foundit(3)) THEN
257 aspval=spval_check
258 ELSE
259 aspval=attvalue(3)
260 END IF
261!
262! Initialize checsum value.
263!
264 IF (PRESENT(checksum)) THEN
265 lchecksum=.true.
266 checksum=0_i8b
267 ELSE
268 lchecksum=.false.
269 END IF
270!
271!-----------------------------------------------------------------------
272! Read in requested data and scale it.
273!-----------------------------------------------------------------------
274!
275 wrk=0.0_r8
276
277 IF (inpthread) THEN
278 status=nf90_get_var(ncid, ncvarid, wrk(lbij:,lbk:,:,:), &
279 & start, total)
280 IF (status.eq.nf90_noerr) THEN
281 amin=spval
282 amax=-spval
283 DO ir=1,nrec
284 DO ib=1,4
285 DO k=lbk,ubk
286 DO ij=lbij,ubij
287 IF (abs(wrk(ij,k,ib,ir)).ge.abs(aspval)) THEN
288 wrk(ij,k,ib,ir)=0.0_r8
289 ELSE
290 wrk(ij,k,ib,ir)=ascl* &
291 & (afactor*wrk(ij,k,ib,ir)+aoffset)
292 amin=min(amin,wrk(ij,k,ib,ir))
293 amax=max(amax,wrk(ij,k,ib,ir))
294 END IF
295 END DO
296 END DO
297 END DO
298 END DO
299 IF ((abs(amin).ge.abs(aspval)).and. &
300 & (abs(amax).ge.abs(aspval))) THEN
301 amin=0.0_r8 ! the entire data is all
302 amax=0.0_r8 ! field value, _FillValue
303 END IF
304!
305 IF (lchecksum) THEN
306 npts=(ubij-lbij+1)*(ubk-lbk+1)*nrec*4
307 IF (.not.allocated(cwrk)) allocate ( cwrk(npts) )
308 cwrk=pack(wrk(lbij:ubij, lbk:ubk, 1:4, 1:nrec), .true.)
309 CALL get_hash (cwrk, npts, checksum)
310 IF (allocated(cwrk)) deallocate (cwrk)
311 END IF
312 END IF
313 END IF
314
315#if !defined PARALLEL_IO && defined DISTRIBUTE
316!
317 rbuffer(1)=real(status,r8)
318 rbuffer(2)=amin
319 rbuffer(3)=amax
320 CALL mp_bcastf (ng, model, rbuffer)
321 status=int(rbuffer(1))
322 amin=rbuffer(2)
323 amax=rbuffer(3)
324#endif
325!
326 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
327 exit_flag=2
328 ioerror=status
329 RETURN
330 END IF
331
332#if !defined PARALLEL_IO && defined DISTRIBUTE
333!
334! Broadcast data to all spawned nodes.
335!
336 CALL mp_bcastf (ng, model, wrk)
337#endif
338!
339!-----------------------------------------------------------------------
340! Unpack read data.
341!-----------------------------------------------------------------------
342!
343 abry=0.0_r8
344
345 DO ir=1,nrec
346 DO ib=1,4
347 IF (bounded(ib)) THEN
348 IF ((ib.eq.iwest).or.(ib.eq.ieast)) THEN
349 DO k=lbk,ubk
350 DO j=jmin,jmax
351 abry(j,k,ib,ir)=wrk(j,k,ib,ir)
352 END DO
353 END DO
354 ELSE IF ((ib.eq.isouth).or.(ib.eq.inorth)) THEN
355 DO k=lbk,ubk
356 DO i=imin,imax
357 abry(i,k,ib,ir)=wrk(i,k,ib,ir)
358 END DO
359 END DO
360 END IF
361 END IF
362 END DO
363 END DO
364!
365 RETURN

References mod_param::bounds, mod_param::domain, mod_scalars::exit_flag, strings_mod::founderror(), get_hash_mod::get_hash(), mod_param::iadm, mod_scalars::ieast, mod_scalars::inorth, mod_parallel::inpthread, mod_param::iobounds, mod_iounits::ioerror, mod_scalars::isouth, mod_scalars::iwest, mod_parallel::myrank, mod_scalars::noerror, mod_param::p2dvar, mod_param::p3dvar, mod_param::r2dvar, mod_param::r3dvar, mod_scalars::spval, mod_scalars::spval_check, mod_param::u2dvar, mod_param::u3dvar, mod_param::v2dvar, and mod_param::v3dvar.

Here is the call graph for this function:

◆ pio_fread3d_bry()

integer function nf_fread3d_bry_mod::nf_fread3d_bry::pio_fread3d_bry ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
type (file_desc_t), intent(inout) piofile,
character (len=*), intent(in) ncvname,
type (my_vardesc), intent(inout) piovar,
integer, intent(in) tindex,
type (io_desc_t), intent(inout) piodesc,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) lbk,
integer, intent(in) ubk,
integer, intent(in) nrec,
real(dp), intent(in) ascl,
real(r8), intent(out) amin,
real(r8), intent(out) amax,
real(r8), dimension(lbij:,:,:,:), intent(out) abry,
integer(i8b), intent(out), optional checksum )

Definition at line 371 of file nf_fread3d_bry.F.

377!***********************************************************************
378!
380!
381! Imported variable declarations.
382!
383 integer, intent(in) :: ng, model, tindex
384 integer, intent(in) :: LBij, UBij, LBk, UBk, Nrec
385!
386 integer(i8b), intent(out), optional :: checksum
387!
388 real(dp), intent(in) :: Ascl
389 real(r8), intent(out) :: Amin
390 real(r8), intent(out) :: Amax
391!
392 character (len=*), intent(in) :: ncname
393 character (len=*), intent(in) :: ncvname
394!
395# ifdef ASSUMED_SHAPE
396 real(r8), intent(out) :: Abry(LBij:,:,:,:)
397# else
398 real(r8), intent(out) :: Abry(LBij:UBij,LBk:UBk,4,Nrec)
399# endif
400!
401 TYPE (File_desc_t), intent(inout) :: pioFile
402 TYPE (IO_Desc_t), intent(inout) :: pioDesc
403 TYPE (My_VarDesc), intent(inout) :: pioVar
404!
405! Local variable declarations.
406!
407 logical :: Lchecksum
408 logical, dimension(3) :: foundit
409 logical, dimension(4) :: bounded
410!
411 integer :: i, ib, ij, ir, j, k
412 integer :: Cgrid, dkind, ghost, gtype, tile
413 integer :: IorJ, IJKlen, Imin, Imax, Jmin, Jmax, Klen, Npts
414 integer :: Istr, Iend, Jstr, Jend
415 integer, dimension(5) :: start, total
416
417 integer :: status
418!
419 real(r8) :: Afactor, Aoffset, Aspval
420!
421 real(r8), allocatable :: Cwrk(:) ! used for checksum
422
423 real(r8), dimension(3) :: AttValue
424
425 real(r8), allocatable :: wrk(:,:,:,:)
426!
427 character (len=12), dimension(3) :: AttName
428
429 character (len=*), parameter :: MyFile = &
430 & __FILE__//", pio_fread3d_bry"
431!
432!-----------------------------------------------------------------------
433! Set starting and ending indices to process.
434!-----------------------------------------------------------------------
435!
436 status=pio_noerr
437!
438! Set first and last grid point according to staggered C-grid
439! classification.
440!
441! Notice that (Imin,Jmin) and (Imax,Jmax) are the corner of the
442! computational tile. If ghost=0, ghost points are not processed.
443! They will be processed elsewhere by the appropriate call to any
444! of the routines in "mp_exchange.F". If ghost=1, the ghost points
445! are read.
446!
447 IF (model.eq.iadm) THEN
448 ghost=0 ! non-overlapping, no ghost points
449 ELSE
450 ghost=1 ! overlapping, read ghost points
451 END IF
452 dkind=piovar%dkind
453 gtype=piovar%gtype
454!
455 SELECT CASE (gtype)
456 CASE (p2dvar, p3dvar)
457 cgrid=1
458 CASE (r2dvar, r3dvar)
459 cgrid=2
460 CASE (u2dvar, u3dvar)
461 cgrid=3
462 CASE (v2dvar, v3dvar)
463 cgrid=4
464 CASE DEFAULT
465 cgrid=2
466 END SELECT
467!
468 tile=myrank
469 imin=bounds(ng)%Imin(cgrid,ghost,tile)
470 imax=bounds(ng)%Imax(cgrid,ghost,tile)
471 jmin=bounds(ng)%Jmin(cgrid,ghost,tile)
472 jmax=bounds(ng)%Jmax(cgrid,ghost,tile)
473!
474 iorj=iobounds(ng)%IorJ
475 klen=ubk-lbk+1
476 ijklen=iorj*klen
477 npts=ijklen*4*nrec
478!
479! Get tile bounds.
480!
481 istr=bounds(ng)%Istr(tile)
482 iend=bounds(ng)%Iend(tile)
483 jstr=bounds(ng)%Jstr(tile)
484 jend=bounds(ng)%Jend(tile)
485!
486! Set switch to process boundary data by their associated tiles.
487!
488 bounded(iwest )=domain(ng)%Western_Edge (tile)
489 bounded(ieast )=domain(ng)%Eastern_Edge (tile)
490 bounded(isouth)=domain(ng)%Southern_Edge(tile)
491 bounded(inorth)=domain(ng)%Northern_Edge(tile)
492!
493! Set NetCDF dimension counters for processing requested field.
494!
495 start(1)=1
496 total(1)=iorj
497 start(2)=1
498 total(2)=klen
499 start(3)=1
500 total(3)=4
501 start(4)=1
502 total(4)=nrec
503 start(5)=tindex
504 total(5)=1
505!
506! Check if the following attributes: "scale_factor", "add_offset", and
507! "_FillValue" are present in the input NetCDF variable:
508!
509! If the "scale_value" attribute is present, the data is multiplied by
510! this factor after reading.
511! If the "add_offset" attribute is present, this value is added to the
512! data after reading.
513! If both "scale_factor" and "add_offset" attributes are present, the
514! data are first scaled before the offset is added.
515! If the "_FillValue" attribute is present, the data having this value
516! is treated a missing and it is replaced with zero. This feature it is
517! usually related with the land/sea masking.
518!
519 attname(1)='scale_factor'
520 attname(2)='add_offset '
521 attname(3)='_FillValue '
522
523 CALL pio_netcdf_get_fatt (ng, model, ncname, piovar%vd, attname, &
524 & attvalue, foundit, &
525 & piofile = piofile)
526 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
527 status=ioerror
528 RETURN
529 END IF
530
531 IF (.not.foundit(1)) THEN
532 afactor=1.0_r8
533 ELSE
534 afactor=attvalue(1)
535 END IF
536
537 IF (.not.foundit(2)) THEN
538 aoffset=0.0_r8
539 ELSE
540 aoffset=attvalue(2)
541 END IF
542
543 IF (.not.foundit(3)) THEN
544 aspval=spval_check
545 ELSE
546 aspval=attvalue(3)
547 END IF
548!
549! Initialize checsum value.
550!
551 IF (PRESENT(checksum)) THEN
552 lchecksum=.true.
553 checksum=0_i8b
554 ELSE
555 lchecksum=.false.
556 END IF
557!
558!-----------------------------------------------------------------------
559! Read in requested data and scale it.
560!-----------------------------------------------------------------------
561!
562 IF (.not.allocated(wrk)) THEN
563 allocate ( wrk(0:iorj-1,lbk:ubk,4,nrec) )
564 wrk=0.0_r8
565 END IF
566!
567 status=pio_get_var(piofile, piovar%vd, start, total, &
568 & wrk(0:,lbk:,:,:))
569 IF (status.eq.pio_noerr) THEN
570 amin=spval
571 amax=-spval
572 DO ir=1,nrec
573 DO ib=1,4
574 DO k=lbk,ubk
575 DO ij=0,iorj-1
576 IF (abs(wrk(ij,k,ib,ir)).ge.abs(aspval)) THEN
577 wrk(ij,k,ib,ir)=0.0_r8
578 ELSE
579 wrk(ij,k,ib,ir)=ascl*(afactor*wrk(ij,k,ib,ir)+aoffset)
580 amin=min(amin,wrk(ij,k,ib,ir))
581 amax=max(amax,wrk(ij,k,ib,ir))
582 END IF
583 END DO
584 END DO
585 END DO
586 END DO
587 IF ((abs(amin).ge.abs(aspval)).and. &
588 & (abs(amax).ge.abs(aspval))) THEN
589 amin=0.0_r8 ! the entire data is all
590 amax=0.0_r8 ! field value, _FillValue
591 END IF
592!
593 IF (lchecksum) THEN
594 npts=(ubij-lbij+1)*(ubk-lbk+1)*nrec*4
595 IF (.not.allocated(cwrk)) allocate ( cwrk(npts) )
596 cwrk=pack(wrk(lbij:ubij, lbk:ubk, 1:4, 1:nrec), .true.)
597 CALL get_hash (cwrk, npts, checksum)
598 IF (allocated(cwrk)) deallocate (cwrk)
599 END IF
600 END IF
601!
602 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
603 exit_flag=2
604 ioerror=status
605 RETURN
606 END IF
607!
608!-----------------------------------------------------------------------
609! Unpack read data.
610!-----------------------------------------------------------------------
611!
612 abry=0.0_r8
613
614 DO ir=1,nrec
615 DO ib=1,4
616 IF (bounded(ib)) THEN
617 IF ((ib.eq.iwest).or.(ib.eq.ieast)) THEN
618 DO k=lbk,ubk
619 DO j=jmin,jmax
620 abry(j,k,ib,ir)=wrk(j,k,ib,ir)
621 END DO
622 END DO
623 ELSE IF ((ib.eq.isouth).or.(ib.eq.inorth)) THEN
624 DO k=lbk,ubk
625 DO i=imin,imax
626 abry(i,k,ib,ir)=wrk(i,k,ib,ir)
627 END DO
628 END DO
629 END IF
630 END IF
631 END DO
632 END DO
633!
634! Deallocate local array.
635!
636 IF (allocated(wrk)) deallocate (wrk)
637!
638 RETURN

References mod_param::bounds, mod_param::domain, mod_scalars::exit_flag, strings_mod::founderror(), get_hash_mod::get_hash(), mod_param::iadm, mod_scalars::ieast, mod_scalars::inorth, mod_param::iobounds, mod_iounits::ioerror, mod_scalars::isouth, mod_scalars::iwest, mod_parallel::myrank, mod_scalars::noerror, mod_param::p2dvar, mod_param::p3dvar, mod_param::r2dvar, mod_param::r3dvar, mod_scalars::spval, mod_scalars::spval_check, mod_param::u2dvar, mod_param::u3dvar, mod_param::v2dvar, and mod_param::v3dvar.

Here is the call graph for this function:

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