ROMS
Loading...
Searching...
No Matches
pack_field_mod Module Reference

Functions/Subroutines

subroutine, public pack_boundary2d (ng, model, tile, gtype, ncvname, tindex, extract_flag, lbij, ubij, nrec, bscl, bdat, start, total, npts, bwrk)
 
subroutine, public pack_boundary3d (ng, model, tile, gtype, ncvname, tindex, extract_flag, lbij, ubij, lbk, ubk, nrec, bscl, bdat, start, total, npts, bwrk)
 
subroutine, public pack_field2d (ng, model, tile, gtype, ifield, tindex, landfill, extract_flag, lbi, ubi, lbj, ubj, amask, ascl, adat, start, total, npts, awrk)
 
subroutine, public pack_field3d (ng, model, tile, gtype, ifield, tindex, landfill, extract_flag, lbi, ubi, lbj, ubj, lbk, ubk, amask, ascl, adat, start, total, npts, awrk)
 
subroutine, public pack_field4d (ng, model, tile, gtype, ifield, tindex, landfill, extract_flag, lbi, ubi, lbj, ubj, lbk, ubk, lbt, ubt, fourth, amask, ascl, adat, start, total, npts, awrk)
 

Function/Subroutine Documentation

◆ pack_boundary2d()

subroutine, public pack_field_mod::pack_boundary2d ( integer, intent(in) ng,
integer, intent(in) model,
integer, intent(in) tile,
integer, intent(in) gtype,
character (len=*), intent(in) ncvname,
integer, intent(in) tindex,
integer, intent(in) extract_flag,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) nrec,
real(dp), intent(in) bscl,
real(r8), dimension(lbij:,:,:), intent(in) bdat,
integer, dimension(:), intent(out) start,
integer, dimension(:), intent(out) total,
integer, intent(in) npts,
real(r8), dimension(:), intent(out) bwrk )

Definition at line 52 of file pack_field.F.

58!
59!=======================================================================
60! !
61! Packs 2D boundary data into 1D array to be written into output !
62! NetCDF file elsewhere. !
63! !
64! On Input: !
65! !
66! ng Nested grid number (integer) !
67! model Calling model identifier (integer) !
68! tile Domain partition (integer) !
69! gtype Staggered C-grid type (integer) !
70! ncvname NetCDF variable name (string) !
71! tindex Time record index to process (integer) !
72! Extract_Flag Extraction flag interpolation/decimation (integer) !
73! LBij Donor IJ-dimension Lower bound (integer) !
74! UBij Donor IJ-dimension Upper bound (integer) !
75! Nrec Number of boundary records (integer) !
76! Bscl Factor to scale boundary data before writing (real) !
77! Bdat 2D boundary data to process (real) !
78! !
79! On Output: !
80! !
81! start Start index where the first of the data values will !
82! be written along each dimension (integer) !
83! total Number of data values to be written along each !
84! dimension (integer) !
85! Npts Number of points processed in Bwrk (integer) !
86! Bwrk Packed 2D boundary data (real) !
87! !
88!=======================================================================
89!
90! Imported variable declarations.
91!
92 integer, intent(in) :: ng, model, tile
93 integer, intent(in) :: gtype, tindex
94 integer, intent(in) :: Extract_Flag
95 integer, intent(in) :: LBij, UBij, Nrec
96 integer, intent(in) :: Npts
97 integer, intent(out) :: start(:), total(:)
98!
99 real(dp), intent(in) :: Bscl
100!
101# ifdef ASSUMED_SHAPE
102 real(r8), intent(in) :: Bdat(LBij:,:,:)
103# else
104 real(r8), intent(in) :: Bdat(LBij:LBij,4,Nrec)
105# endif
106 real(r8), intent(out) :: Bwrk(:)
107!
108 character (len=*), intent(in) :: ncvname
109!
110! Local variable declarations.
111!
112 logical, dimension(4) :: bounded
113!
114 integer :: bc, i, ib, ic, ir, j, rc
115 integer :: IorJ, Imin, Imax, Jmin, Jmax
116!
117 real(r8), parameter :: Bspv = 0.0_r8
118!
119 character (len=*), parameter :: MyFile = &
120 & __FILE__//", pack_boundary2d"
121!
122!------------------------------------------------------------------------
123! Pack 2D boundary data.
124!------------------------------------------------------------------------
125!
126! Set first and last grid point according to staggered C-grid
127! classification.
128!
129# ifdef DISTRIBUTE
130 SELECT CASE (gtype)
131!
132! ---------------------
133 CASE (p2dvar, p3dvar) ! PSI points field
134! ---------------------
135!
136 imin=bounds(ng)%Istr (tile)
137 imax=bounds(ng)%Iend (tile)
138 jmin=bounds(ng)%Jstr (tile)
139 jmax=bounds(ng)%Jend (tile)
140!
141! ---------------------
142 CASE (r2dvar, r3dvar) ! RHO points field
143! ---------------------
144!
145 imin=bounds(ng)%IstrR(tile)
146 imax=bounds(ng)%IendR(tile)
147 jmin=bounds(ng)%JstrR(tile)
148 jmax=bounds(ng)%JendR(tile)
149!
150! ---------------------
151 CASE (u2dvar, u3dvar) ! U points field
152! ---------------------
153!
154 imin=bounds(ng)%Istr (tile)
155 imax=bounds(ng)%IendR(tile)
156 jmin=bounds(ng)%JstrR(tile)
157 jmax=bounds(ng)%JendR(tile)
158!
159! ---------------------
160 CASE (v2dvar, v3dvar) ! V points field
161! ---------------------
162!
163 imin=bounds(ng)%IstrR(tile)
164 imax=bounds(ng)%IendR(tile)
165 jmin=bounds(ng)%Jstr (tile)
166 jmax=bounds(ng)%JendR(tile)
167!
168! ---------------------
169 CASE DEFAULT ! RHO points field
170! ---------------------
171!
172 imin=bounds(ng)%IstrR(tile)
173 imax=bounds(ng)%IendR(tile)
174 jmin=bounds(ng)%JstrR(tile)
175 jmax=bounds(ng)%JendR(tile)
176 END SELECT
177# else
178 imin=lbij
179 imax=ubij
180 jmin=lbij
181 jmax=ubij
182# endif
183!
184 iorj=iobounds(ng)%IorJ
185!
186! Set switch to process boundary data by their associated tiles.
187!
188 bounded(iwest )=domain(ng)%Western_Edge (tile)
189 bounded(ieast )=domain(ng)%Eastern_Edge (tile)
190 bounded(isouth)=domain(ng)%Southern_Edge(tile)
191 bounded(inorth)=domain(ng)%Northern_Edge(tile)
192!
193!-----------------------------------------------------------------------
194! Pack and scale output data.
195!-----------------------------------------------------------------------
196!
197 bwrk=bspv
198!
199 DO ir=1,nrec
200 rc=(ir-1)*iorj*4
201 DO ib=1,4
202 IF (bounded(ib)) THEN
203 bc=(ib-1)*iorj+rc
204 IF ((ib.eq.iwest).or.(ib.eq.ieast)) THEN
205 DO j=jmin,jmax
206 ic=1+(j-lbij)+bc
207 bwrk(ic)=bdat(j,ib,ir)*bscl
208# ifdef POSITIVE_ZERO
209 IF (abs(bwrk(ic)).eq.0.0_r8) THEN
210 bwrk(ic)=0.0_r8 ! impose positive zero
211 END IF
212# endif
213 END DO
214 ELSE IF ((ib.eq.isouth).or.(ib.eq.inorth)) THEN
215 DO i=imin,imax
216 ic=1+(i-lbij)+bc
217 bwrk(ic)=bdat(i,ib,ir)*bscl
218# ifdef POSITIVE_ZERO
219 IF (abs(bwrk(ic)).eq.0.0_r8) THEN
220 bwrk(ic)=0.0_r8 ! impose positive zero
221 END IF
222# endif
223 END DO
224 END IF
225 END IF
226 END DO
227 END DO
228
229# ifdef DISTRIBUTE
230!
231! If distributed-memory set-up, collect data from all spawned
232! processes.
233!
234 CALL mp_collect (ng, model, npts, bspv, bwrk)
235# endif
236!
237!-----------------------------------------------------------------------
238! If there is no extracting data, set the start and total vectors
239! needed by the NetCDF library for writing.
240!-----------------------------------------------------------------------
241!
242 IF (extract_flag.le.0) THEN
243 start(1)=1
244 total(1)=iorj
245 start(2)=1
246 total(2)=4
247 start(3)=1
248 total(3)=nrec
249 start(4)=tindex
250 total(4)=1
251
252# ifdef GRID_EXTRACT
253!
254!-----------------------------------------------------------------------
255! Otherwise, extract field by decimation or interpolation.
256!-----------------------------------------------------------------------
257!
258 ELSE IF (extract_flag.ge.1) THEN
259 CALL extract_boundary (ng, model, tile, &
260 & gtype, ncvname, tindex, &
261 & extract_flag, &
262 & imin, imax, jmin, jmax, nrec, &
263 & start, total, &
264 & npts, bwrk)
265# endif
266 END IF
267!
268 RETURN

References mod_param::bounds, mod_param::domain, mod_scalars::ieast, mod_scalars::inorth, mod_param::iobounds, mod_scalars::isouth, mod_scalars::iwest, mod_param::p2dvar, mod_param::p3dvar, mod_param::r2dvar, mod_param::r3dvar, mod_param::u2dvar, mod_param::u3dvar, mod_param::v2dvar, and mod_param::v3dvar.

Referenced by nf_fwrite2d_bry_mod::nf_fwrite2d_bry::nf90_fwrite2d_bry().

Here is the caller graph for this function:

◆ pack_boundary3d()

subroutine, public pack_field_mod::pack_boundary3d ( integer, intent(in) ng,
integer, intent(in) model,
integer, intent(in) tile,
integer, intent(in) gtype,
character (len=*), intent(in) ncvname,
integer, intent(in) tindex,
integer, intent(in) extract_flag,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) lbk,
integer, intent(in) ubk,
integer, intent(in) nrec,
real(dp), intent(in) bscl,
real(r8), dimension(lbij:,lbk:,:,:), intent(in) bdat,
integer, dimension(:), intent(out) start,
integer, dimension(:), intent(out) total,
integer, intent(in) npts,
real(r8), dimension(:), intent(out) bwrk )

Definition at line 271 of file pack_field.F.

277!
278!=======================================================================
279! !
280! Packs 3D boundary data into 1D array to be written into output !
281! NetCDF file elsewhere. !
282! !
283! On Input: !
284! !
285! ng Nested grid number (integer) !
286! model Calling model identifier (integer) !
287! tile Domain partition (integer) !
288! gtype Staggered C-grid type (integer) !
289! ncvname NetCDF variable name (string) !
290! tindex Time record index to process (integer) !
291! Extract_Flag Extraction flag interpolation/decimation (integer) !
292! LBij Donor IJ-dimension Lower bound (integer) !
293! UBij Donor IJ-dimension Upper bound (integer) !
294! LBk Donor K-dimension lower bound (integer) !
295! UBk Donor K-dimension upper bound (integer) !
296! Nrec Number of boundary records (integer) !
297! Bscl Factor to scale boundary data before writing (real) !
298! Bdat 3D boundary data to process (real) !
299! !
300! On Output: !
301! !
302! start Start index where the first of the data values will !
303! be written along each dimension (integer) !
304! total Number of data values to be written along each !
305! dimension (integer) !
306! Npts Number of points processed in Bwrk (integer) !
307! Bwrk Packed 3D boundary data (real) !
308! !
309!=======================================================================
310!
311! Imported variable declarations.
312!
313 integer, intent(in) :: ng, model, tile
314 integer, intent(in) :: gtype, tindex
315 integer, intent(in) :: Extract_Flag
316 integer, intent(in) :: LBij, UBij, LBk, UBk, Nrec
317 integer, intent(in) :: Npts
318 integer, intent(out) :: start(:), total(:)
319!
320 real(dp), intent(in) :: Bscl
321!
322# ifdef ASSUMED_SHAPE
323 real(r8), intent(in) :: Bdat(LBij:,LBk:,:,:)
324# else
325 real(r8), intent(in) :: Bdat(LBij:LBij,LBk:UBk,4,Nrec)
326# endif
327 real(r8), intent(out) :: Bwrk(:)
328!
329 character (len=*), intent(in) :: ncvname
330!
331! Local variable declarations.
332!
333 logical, dimension(4) :: bounded
334!
335 integer :: Imin, Imax, Jmin, Jmax
336 integer :: IorJ, IJKlen, Klen
337 integer :: bc, i, ib, ic, ir, j, k, kc, rc
338!
339 real(r8), parameter :: Bspv = 0.0_r8
340!
341 character (len=*), parameter :: MyFile = &
342 & __FILE__//", pack_boundary3d"
343!
344!------------------------------------------------------------------------
345! Pack 3D boundary data.
346!------------------------------------------------------------------------
347!
348! Set first and last grid point according to staggered C-grid
349! classification.
350!
351# ifdef DISTRIBUTE
352 SELECT CASE (gtype)
353!
354! ---------------------
355 CASE (p2dvar, p3dvar) ! PSI points field
356! ---------------------
357!
358 imin=bounds(ng)%Istr (tile)
359 imax=bounds(ng)%Iend (tile)
360 jmin=bounds(ng)%Jstr (tile)
361 jmax=bounds(ng)%Jend (tile)
362!
363! ---------------------
364 CASE (r2dvar, r3dvar) ! RHO points field
365! ---------------------
366!
367 imin=bounds(ng)%IstrR(tile)
368 imax=bounds(ng)%IendR(tile)
369 jmin=bounds(ng)%JstrR(tile)
370 jmax=bounds(ng)%JendR(tile)
371!
372! ---------------------
373 CASE (u2dvar, u3dvar) ! U points field
374! ---------------------
375!
376 imin=bounds(ng)%Istr (tile)
377 imax=bounds(ng)%IendR(tile)
378 jmin=bounds(ng)%JstrR(tile)
379 jmax=bounds(ng)%JendR(tile)
380!
381! ---------------------
382 CASE (v2dvar, v3dvar) ! V points field
383! ---------------------
384!
385 imin=bounds(ng)%IstrR(tile)
386 imax=bounds(ng)%IendR(tile)
387 jmin=bounds(ng)%Jstr (tile)
388 jmax=bounds(ng)%JendR(tile)
389!
390! ---------------------
391 CASE DEFAULT ! RHO points field
392! ---------------------
393!
394 imin=bounds(ng)%IstrR(tile)
395 imax=bounds(ng)%IendR(tile)
396 jmin=bounds(ng)%JstrR(tile)
397 jmax=bounds(ng)%JendR(tile)
398 END SELECT
399# else
400 imin=lbij
401 imax=ubij
402 jmin=lbij
403 jmax=ubij
404# endif
405!
406 iorj=iobounds(ng)%IorJ
407 klen=ubk-lbk+1
408 ijklen=iorj*klen
409!
410! Set switch to process boundary data by their associated tiles.
411!
412 bounded(iwest )=domain(ng)%Western_Edge (tile)
413 bounded(ieast )=domain(ng)%Eastern_Edge (tile)
414 bounded(isouth)=domain(ng)%Southern_Edge(tile)
415 bounded(inorth)=domain(ng)%Northern_Edge(tile)
416!
417!-----------------------------------------------------------------------
418! Pack and scale output data.
419!-----------------------------------------------------------------------
420!
421 bwrk=bspv
422!
423 DO ir=1,nrec
424 rc=(ir-1)*ijklen*4
425 DO ib=1,4
426 IF (bounded(ib)) THEN
427 bc=(ib-1)*ijklen+rc
428 IF ((ib.eq.iwest).or.(ib.eq.ieast)) THEN
429 DO k=lbk,ubk
430 kc=(k-lbk)*iorj+bc
431 DO j=jmin,jmax
432 ic=1+(j-lbij)+kc
433 bwrk(ic)=bdat(j,k,ib,ir)*bscl
434# ifdef POSITIVE_ZERO
435 IF (abs(bwrk(ic)).eq.0.0_r8) THEN
436 bwrk(ic)=0.0_r8 ! impose positive zero
437 END IF
438# endif
439 END DO
440 END DO
441 ELSE IF ((ib.eq.isouth).or.(ib.eq.inorth)) THEN
442 DO k=lbk,ubk
443 kc=(k-lbk)*iorj+bc
444 DO i=imin,imax
445 ic=1+(i-lbij)+kc
446 bwrk(ic)=bdat(i,k,ib,ir)*bscl
447# ifdef POSITIVE_ZERO
448 IF (abs(bwrk(ic)).eq.0.0_r8) THEN
449 bwrk(ic)=0.0_r8 ! impose positive zero
450 END IF
451# endif
452 END DO
453 END DO
454 END IF
455 END IF
456 END DO
457 END DO
458
459# ifdef DISTRIBUTE
460!
461! If distributed-memory set-up, collect data from all spawned
462! processes.
463!
464 CALL mp_collect (ng, model, npts, bspv, bwrk)
465# endif
466!
467!-----------------------------------------------------------------------
468! If there is no extracting data, set the start and total vectors
469! needed by the NetCDF library for writing.
470!-----------------------------------------------------------------------
471!
472 IF (extract_flag.le.0) THEN
473 start(1)=1
474 total(1)=iorj
475 start(2)=1
476 total(2)=klen
477 start(3)=1
478 total(3)=4
479 start(4)=1
480 total(4)=nrec
481 start(5)=tindex
482 total(5)=1
483
484# ifdef GRID_EXTRACT
485!
486!-----------------------------------------------------------------------
487! Otherwise, extract field by decimation or interpolation.
488!-----------------------------------------------------------------------
489!
490 ELSE IF (extract_flag.ge.1) THEN
491 CALL extract_boundary (ng, model, tile, &
492 & gtype, ncvname, tindex, &
493 & extract_flag, &
494 & imin, imax, jmin, jmax, lbk, ubk, nrec, &
495 & start, total, &
496 & npts, bwrk)
497# endif
498 END IF
499!
500 RETURN

References mod_param::bounds, mod_param::domain, mod_scalars::ieast, mod_scalars::inorth, mod_param::iobounds, mod_scalars::isouth, mod_scalars::iwest, mod_param::p2dvar, mod_param::p3dvar, mod_param::r2dvar, mod_param::r3dvar, mod_param::u2dvar, mod_param::u3dvar, mod_param::v2dvar, and mod_param::v3dvar.

Referenced by nf_fwrite3d_bry_mod::nf_fwrite3d_bry::nf90_fwrite3d_bry().

Here is the caller graph for this function:

◆ pack_field2d()

subroutine, public pack_field_mod::pack_field2d ( integer, intent(in) ng,
integer, intent(in) model,
integer, intent(in) tile,
integer, intent(in) gtype,
integer, intent(in) ifield,
integer, intent(in) tindex,
logical, intent(in) landfill,
integer, intent(in) extract_flag,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
real(r8), dimension(lbi:,lbj:), intent(in) amask,
real(dp), intent(in) ascl,
real(r8), dimension(lbi:,lbj:), intent(in) adat,
integer, dimension(:), intent(out) start,
integer, dimension(:), intent(out) total,
integer, intent(out) npts,
real(r8), dimension(:), intent(out) awrk )

Definition at line 504 of file pack_field.F.

513!
514!=======================================================================
515! !
516! Packs 2D field data into 1D array to be written into output NetCDF !
517! file elsewhere. !
518! !
519! On Input: !
520! !
521! ng Nested grid number (integer) !
522! model Calling model identifier (integer) !
523! tile Domain partition (integer) !
524! gtype Staggered C-grid type (integer) !
525! ifield Field metadata index (integer) !
526! tindex Time record index to process (integer) !
527! LandFill Switch to set fill value in land areas (logical) ! !
528! Extract_Flag Extraction flag interpolation/decimation (integer) !
529! LBi Donor field I-dimension Lower bound (integer) !
530! UBi Donor field I-dimension Upper bound (integer) !
531! LBj Donor field J-dimension Lower bound (integer) !
532! UBj Donor field J-dimension Upper bound (integer) !
533#ifdef MASKING
534! Amask Land/Sea mask, if any (real 2D array) !
535#endif
536! Ascl Factor to scale field after reading (real) !
537! Adat 2D field to process (real) !
538! !
539! On Output: !
540! !
541! start Start index where the first of the data values will !
542! be written along each dimension (integer) !
543! total Number of data values to be written along each !
544! dimension (integer) !
545! Npts Number of points processed in Awrk. !
546! Awrk Packed 2D field data (real) !
547! !
548!=======================================================================
549!
550! Imported variable declarations.
551!
552 logical, intent(in) :: LandFill
553!
554 integer, intent(in) :: ng, model, tile
555 integer, intent(in) :: gtype, ifield, tindex
556 integer, intent(in) :: Extract_Flag
557 integer, intent(in) :: LBi, UBi, LBj, UBj
558 integer, intent(out) :: Npts
559 integer, intent(out) :: start(:), total(:)
560!
561 real(dp), intent(in) :: Ascl
562!
563#ifdef ASSUMED_SHAPE
564# ifdef MASKING
565 real(r8), intent(in) :: Amask(LBi:,LBj:)
566# endif
567 real(r8), intent(in) :: Adat(LBi:,LBj:)
568#else
569# ifdef MASKING
570 real(r8), intent(in) :: Amask(LBi:UBi,LBj:UBj)
571# endif
572 real(r8), intent(in) :: Adat(LBi:UBi,LBj:UBj)
573#endif
574 real(r8), intent(out) :: Awrk(:)
575!
576! Local variable declarations.
577!
578 integer :: i, ic, j
579 integer :: Imin, Imax, Jmin, Jmax
580 integer :: Isize, Jsize, IJsize, MyType
581!
582 character (len=*), parameter :: MyFile = &
583 & __FILE__//", pack_field2d"
584
585!-----------------------------------------------------------------------
586! Set starting and ending indices to process.
587!-----------------------------------------------------------------------
588!
589! Set first and last grid point according to staggered C-grid
590! classification. Set loops offsets.
591#ifdef GRID_EXTRACT
592!
593! Here, the ELSE statements (Extract_Flag < 0) indicate the processing
594! of data already on the extracted grid, like geometry variables
595! requiring no extraction since they are read during initialization.
596#endif
597!
598 mytype=gtype
599!
600 SELECT CASE (abs(mytype))
601!
602! ---------------------
603 CASE (p2dvar, p3dvar) ! PSI points field
604! ---------------------
605!
606 IF (extract_flag.ge.0) THEN
607 imin=iobounds(ng)%ILB_psi
608 imax=iobounds(ng)%IUB_psi
609 jmin=iobounds(ng)%JLB_psi
610 jmax=iobounds(ng)%JUB_psi
611#ifdef GRID_EXTRACT
612 ELSE ! Already on extract grid
613 imin=xtr_iobounds(ng)%ILB_psi
614 imax=xtr_iobounds(ng)%IUB_psi
615 jmin=xtr_iobounds(ng)%JLB_psi
616 jmax=xtr_iobounds(ng)%JUB_psi
617#endif
618 END IF
619!
620! ---------------------
621 CASE (r2dvar, r3dvar) ! RHO points field
622! ---------------------
623!
624 IF (extract_flag.ge.0) THEN
625 imin=iobounds(ng)%ILB_rho
626 imax=iobounds(ng)%IUB_rho
627 jmin=iobounds(ng)%JLB_rho
628 jmax=iobounds(ng)%JUB_rho
629#ifdef GRID_EXTRACT
630 ELSE ! Already on extract grid
631 imin=xtr_iobounds(ng)%ILB_rho
632 imax=xtr_iobounds(ng)%IUB_rho
633 jmin=xtr_iobounds(ng)%JLB_rho
634 jmax=xtr_iobounds(ng)%JUB_rho
635#endif
636 END IF
637!
638! ---------------------
639 CASE (u2dvar, u3dvar) ! U points field
640! ---------------------
641!
642 IF (extract_flag.ge.0) THEN
643 imin=iobounds(ng)%ILB_u
644 imax=iobounds(ng)%IUB_u
645 jmin=iobounds(ng)%JLB_u
646 jmax=iobounds(ng)%JUB_u
647#ifdef GRID_EXTRACT
648 ELSE ! Already on extract grid
649 imin=xtr_iobounds(ng)%ILB_u
650 imax=xtr_iobounds(ng)%IUB_u
651 jmin=xtr_iobounds(ng)%JLB_u
652 jmax=xtr_iobounds(ng)%JUB_u
653#endif
654 END IF
655!
656! ---------------------
657 CASE (v2dvar, v3dvar) ! V points field
658! ---------------------
659!
660 IF (extract_flag.ge.0) THEN
661 imin=iobounds(ng)%ILB_v
662 imax=iobounds(ng)%IUB_v
663 jmin=iobounds(ng)%JLB_v
664 jmax=iobounds(ng)%JUB_v
665#ifdef GRID_EXTRACT
666 ELSE ! Already on extract grid
667 imin=xtr_iobounds(ng)%ILB_v
668 imax=xtr_iobounds(ng)%IUB_v
669 jmin=xtr_iobounds(ng)%JLB_v
670 jmax=xtr_iobounds(ng)%JUB_v
671#endif
672 END IF
673!
674! ---------------------
675 CASE DEFAULT ! RHO points field
676! ---------------------
677!
678 IF (extract_flag.ge.0) THEN
679 imin=iobounds(ng)%ILB_rho
680 imax=iobounds(ng)%IUB_rho
681 jmin=iobounds(ng)%JLB_rho
682 jmax=iobounds(ng)%JUB_rho
683#ifdef GRID_EXTRACT
684 ELSE ! Already on extract grid
685 imin=xtr_iobounds(ng)%ILB_rho
686 imax=xtr_iobounds(ng)%IUB_rho
687 jmin=xtr_iobounds(ng)%JLB_rho
688 jmax=xtr_iobounds(ng)%JUB_rho
689#endif
690 END IF
691 END SELECT
692!
693 isize=imax-imin+1
694 jsize=jmax-jmin+1
695 ijsize=isize*jsize
696
697#ifdef DISTRIBUTE
698!
699!-----------------------------------------------------------------------
700! If distributed-memory set-up, collect tile data from all spawned
701! nodes and store it into a global scratch 1D array, packed in column-
702! major order.
703# ifdef MASKING
704# ifdef WRITE_WATER
705! Remove land points and pack water points into 1D-array.
706# else
707! Overwrite masked points with special value.
708# endif
709# endif
710!-----------------------------------------------------------------------
711!
712 IF (extract_flag.ge.0) THEN
713 CALL mp_gather2d (ng, model, lbi, ubi, lbj, ubj, &
714 & tindex, gtype, ascl, &
715# ifdef MASKING
716 & amask, &
717# endif
718 & adat, npts, awrk, landfill)
719# ifdef GRID_EXTRACT
720 ELSE ! Already on extract grid
721 CALL mp_gather2d_xtr (ng, model, lbi, ubi, lbj, ubj, &
722 & tindex, gtype, ascl, &
723# ifdef MASKING
724 & amask, &
725# endif
726 & adat, npts, awrk, landfill)
727# endif
728 END IF
729
730#else
731!
732!-----------------------------------------------------------------------
733! If serial or shared-memory applications and serial output, pack data
734! into a global 1D array in column-major order.
735# ifdef MASKING
736# ifdef WRITE_WATER
737! Remove land points and pack water points into 1D-array.
738# else
739! Overwrite masked points with special value.
740# endif
741# endif
742!-----------------------------------------------------------------------
743!
744 IF (gtype.gt.0) THEN
745 ic=0
746 DO j=jmin,jmax
747 DO i=imin,imax
748 ic=ic+1
749 awrk(ic)=adat(i,j)*ascl
750# ifdef MASKING
751 IF ((amask(i,j).eq.0.0_r8).and.landfill) THEN
752 awrk(ic)=spval
753 END IF
754# endif
755 END DO
756 END DO
757 npts=ic
758# ifdef MASKING
759 ELSE
760 ic=0
761 DO j=jmin,jmax
762 DO i=imin,imax
763 IF (amask(i,j).gt.0.0_r8) THEN
764 ic=ic+1
765 awrk(ic)=adat(i,j)*ascl
766 END IF
767 END DO
768 END DO
769 npts=ic
770# endif
771 END IF
772#endif
773!
774!-----------------------------------------------------------------------
775! If there is no extracting data, set the start and total vectors
776! needed by the NetCDF library for writing.
777!-----------------------------------------------------------------------
778!
779 IF (extract_flag.le.0) THEN
780 IF (gtype.gt.0) THEN
781 start(1)=1
782 total(1)=isize
783 start(2)=1
784 total(2)=jsize
785 start(3)=tindex
786 total(3)=1
787#ifdef MASKING
788 ELSE
789 start(1)=1
790 total(1)=npts
791 start(2)=tindex
792 total(2)=1
793#endif
794 END IF
795
796#ifdef GRID_EXTRACT
797!
798!-----------------------------------------------------------------------
799! Otherwise, extract field by decimation or interpolation.
800!-----------------------------------------------------------------------
801!
802 ELSE IF (extract_flag.ge.1) THEN
803 CALL extract_field (ng, model, tile, &
804 & gtype, ifield, tindex, &
805 & extract_flag, &
806 & imin, imax, jmin, jmax, &
807 & start, total, &
808 & npts, awrk)
809#endif
810 END IF
811!
812 RETURN

References mod_param::iobounds, distribute_mod::mp_gather2d(), mod_param::p2dvar, mod_param::p3dvar, mod_param::r2dvar, mod_param::r3dvar, mod_scalars::spval, mod_param::u2dvar, mod_param::u3dvar, mod_param::v2dvar, and mod_param::v3dvar.

Referenced by nf_fwrite2d_mod::nf_fwrite2d::nf90_fwrite2d().

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

◆ pack_field3d()

subroutine, public pack_field_mod::pack_field3d ( integer, intent(in) ng,
integer, intent(in) model,
integer, intent(in) tile,
integer, intent(in) gtype,
integer, intent(in) ifield,
integer, intent(in) tindex,
logical, intent(in) landfill,
integer, intent(in) extract_flag,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) lbk,
integer, intent(in) ubk,
real(r8), dimension(lbi:,lbj:), intent(in) amask,
real(dp), intent(in) ascl,
real(r8), dimension(lbi:,lbj:,lbk:), intent(in) adat,
integer, dimension(:), intent(out) start,
integer, dimension(:), intent(out) total,
integer, intent(out) npts,
real(r8), dimension(:), intent(out) awrk )

Definition at line 815 of file pack_field.F.

824!
825!=======================================================================
826! !
827! Packs 3D field data into 1D array to be written into output NetCDF !
828! file elsewhere. !
829! !
830! On Input: !
831! !
832! ng Nested grid number (integer) !
833! model Calling model identifier (integer) !
834! tile Domain partition (integer) !
835! gtype Staggered C-grid type (integer) !
836! ifield Field metadata index (integer) !
837! tindex Time record index to process (integer) !
838! LandFill Switch to set fill value in land areas (logical) !
839! Extract_Flag Extraction flag interpolation/decimation (integer) !
840! LBi Field I-dimension Lower bound (integer) !
841! UBi Field I-dimension Upper bound (integer) !
842! LBj Field J-dimension Lower bound (integer) !
843! UBj Field J-dimension Upper bound (integer) !
844! LBk Field K-dimension Lower bound (integer) !
845! UBk Field K-dimension Upper bound (integer) !
846#ifdef MASKING
847! Amask Land/Sea mask, if any (real 2D array) !
848#endif
849! Ascl Factor to scale field after reading (real) !
850! Adat 3D field to process (real) !
851! !
852! On Output: !
853! !
854! start Start index where the first of the data values will !
855! be written along each dimension (integer) !
856! total Number of data values to be written along each !
857! dimension (integer) !
858! Npts Number of points processed in Awrk. !
859! Awrk Packed 3D field data (real) !
860! !
861!=======================================================================
862!
863! Imported variable declarations.
864!
865 logical, intent(in) :: LandFill
866!
867 integer, intent(in) :: ng, model, tile
868 integer, intent(in) :: gtype, ifield, tindex
869 integer, intent(in) :: Extract_Flag
870 integer, intent(in) :: LBi, UBi, LBj, UBj, LBk, UBk
871 integer, intent(out) :: Npts
872 integer, intent(out) :: start(:), total(:)
873!
874 real(dp), intent(in) :: Ascl
875!
876#ifdef ASSUMED_SHAPE
877# ifdef MASKING
878 real(r8), intent(in) :: Amask(LBi:,LBj:)
879# endif
880 real(r8), intent(in) :: Adat(LBi:,LBj:,LBk:)
881#else
882# ifdef MASKING
883 real(r8), intent(in) :: Amask(LBi:UBi,LBj:UBj)
884# endif
885 real(r8), intent(in) :: Adat(LBi:UBi,LBj:UBj,LBk:UBk)
886#endif
887 real(r8), intent(out) :: Awrk(:)
888!
889! Local variable declarations.
890!
891 integer :: i, j, k, ic
892 integer :: Imin, Imax, Jmin, Jmax, Koff
893 integer :: Isize, Jsize, Ksize, IJsize, MyType
894!
895 character (len=*), parameter :: MyFile = &
896 & __FILE__//", pack_field3d"
897!
898!------------------------------------------------------------------------
899! Pack 3D field data.
900!------------------------------------------------------------------------
901!
902! Set first and last grid point according to staggered C-grid
903! classification. Set loops offsets.
904!
905 mytype=gtype
906!
907 SELECT CASE (abs(mytype))
908!
909! ---------------------
910 CASE (p3dvar) ! PSI points field
911! ---------------------
912!
913 imin=iobounds(ng)%ILB_psi
914 imax=iobounds(ng)%IUB_psi
915 jmin=iobounds(ng)%JLB_psi
916 jmax=iobounds(ng)%JUB_psi
917!
918! ---------------------
919 CASE (r3dvar) ! RHO points field
920! ---------------------
921!
922 imin=iobounds(ng)%ILB_rho
923 imax=iobounds(ng)%IUB_rho
924 jmin=iobounds(ng)%JLB_rho
925 jmax=iobounds(ng)%JUB_rho
926!
927! ---------------------
928 CASE (u3dvar) ! U points field
929! ---------------------
930!
931 imin=iobounds(ng)%ILB_u
932 imax=iobounds(ng)%IUB_u
933 jmin=iobounds(ng)%JLB_u
934 jmax=iobounds(ng)%JUB_u
935!
936! ---------------------
937 CASE (v3dvar) ! V points field
938! ---------------------
939!
940 imin=iobounds(ng)%ILB_v
941 imax=iobounds(ng)%IUB_v
942 jmin=iobounds(ng)%JLB_v
943 jmax=iobounds(ng)%JUB_v
944!
945! ---------------------
946 CASE DEFAULT ! RHO points field
947! ---------------------
948!
949 imin=iobounds(ng)%ILB_rho
950 imax=iobounds(ng)%IUB_rho
951 jmin=iobounds(ng)%JLB_rho
952 jmax=iobounds(ng)%JUB_rho
953 END SELECT
954!
955 isize=imax-imin+1
956 jsize=jmax-jmin+1
957 ksize=ubk-lbk+1
958 ijsize=isize*jsize
959!
960 IF (lbk.eq.0) THEN
961 koff=0
962 ELSE
963 koff=1
964 END IF
965
966#ifdef DISTRIBUTE
967!
968!-----------------------------------------------------------------------
969! If distributed-memory set-up, collect tile data from all spawned
970! nodes and store it into a global scratch 1D array, packed in column-
971! major order.
972# ifdef MASKING
973# ifdef WRITE_WATER
974! Remove land points and pack water points into 1D-array.
975# else
976! Overwrite masked points with special value.
977# endif
978# endif
979!-----------------------------------------------------------------------
980!
981 CALL mp_gather3d (ng, model, lbi, ubi, lbj, ubj, lbk, ubk, &
982 & tindex, gtype, ascl, &
983# ifdef MASKING
984 & amask, &
985# endif
986 & adat, npts, awrk, landfill)
987#else
988!
989!-----------------------------------------------------------------------
990! If serial or shared-memory applications and serial output, pack data
991! into a global 1D array in column-major order.
992# ifdef MASKING
993# ifdef WRITE_WATER
994! Remove land points and pack water points into 1D-array.
995# else
996! Overwrite masked points with special value.
997# endif
998# endif
999!-----------------------------------------------------------------------
1000!
1001 IF (gtype.gt.0) THEN
1002 ic=0
1003 DO k=lbk,ubk
1004 DO j=jmin,jmax
1005 DO i=imin,imax
1006 ic=ic+1
1007 awrk(ic)=adat(i,j,k)*ascl
1008# ifdef MASKING
1009 IF ((amask(i,j).eq.0.0_r8).and.landfill) THEN
1010 awrk(ic)=spval
1011 END IF
1012# endif
1013 END DO
1014 END DO
1015 END DO
1016 npts=ic
1017# ifdef MASKING
1018 ELSE
1019 ic=0
1020 DO k=lbk,ubk
1021 DO j=jmin,jmax
1022 DO i=imin,imax
1023 IF (amask(i,j).gt.0.0_r8) THEN
1024 ic=ic+1
1025 awrk(ic)=adat(i,j,k)*ascl
1026 END IF
1027 END DO
1028 END DO
1029 END DO
1030 npts=ic
1031# endif
1032 END IF
1033#endif
1034!
1035!-----------------------------------------------------------------------
1036! If there is no extracting data, set the start and total vectors
1037! needed by the NetCDF library for writing.
1038!-----------------------------------------------------------------------
1039!
1040 IF (extract_flag.le.0) THEN
1041 IF (gtype.gt.0) THEN
1042 start(1)=1
1043 total(1)=isize
1044 start(2)=1
1045 total(2)=jsize
1046 start(3)=1
1047 total(3)=ksize
1048 start(4)=tindex
1049 total(4)=1
1050#ifdef MASKING
1051 ELSE
1052 start(1)=1
1053 total(1)=npts
1054 start(2)=tindex
1055 total(2)=1
1056#endif
1057 END IF
1058
1059#ifdef GRID_EXTRACT
1060!
1061!-----------------------------------------------------------------------
1062! Otherwise, extract field by decimation or interpolation.
1063!-----------------------------------------------------------------------
1064!
1065 ELSE IF (extract_flag.ge.1) THEN
1066 CALL extract_field (ng, model, tile, &
1067 & gtype, ifield, tindex, &
1068 & extract_flag, &
1069 & imin, imax, jmin, jmax, lbk, ubk, &
1070 & start, total, &
1071 & npts, awrk)
1072#endif
1073 END IF
1074
1075#ifdef POSITIVE_ZERO
1076!
1077!-----------------------------------------------------------------------
1078! Impose positive zero. Recall that IEEE allows both negative and
1079! positive zeros, which affects bit-by-bit comparison betwee different
1080! tile partitions.
1081!-----------------------------------------------------------------------
1082!
1083 DO ic=1,npts
1084 IF (abs(awrk(ic)).eq.0.0_r8) THEN
1085 awrk(ic)=0.0_r8
1086 END IF
1087 END DO
1088#endif
1089!
1090 RETURN

References mod_param::iobounds, distribute_mod::mp_gather3d(), mod_param::p3dvar, mod_param::r3dvar, mod_scalars::spval, mod_param::u3dvar, and mod_param::v3dvar.

Referenced by nf_fwrite3d_mod::nf_fwrite3d::nf90_fwrite3d().

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

◆ pack_field4d()

subroutine, public pack_field_mod::pack_field4d ( integer, intent(in) ng,
integer, intent(in) model,
integer, intent(in) tile,
integer, intent(in) gtype,
integer, intent(in) ifield,
integer, intent(in) tindex,
logical, intent(in) landfill,
integer, intent(in) extract_flag,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) lbk,
integer, intent(in) ubk,
integer, intent(in) lbt,
integer, intent(in) ubt,
integer, intent(in) fourth,
real(r8), dimension(lbi:,lbj:), intent(in) amask,
real(dp), intent(in) ascl,
real(r8), dimension(lbi:,lbj:,lbk:,lbt:), intent(in) adat,
integer, dimension(:), intent(out) start,
integer, dimension(:), intent(out) total,
integer, intent(out) npts,
real(r8), dimension(:), intent(out) awrk )

Definition at line 1093 of file pack_field.F.

1103!
1104!=======================================================================
1105! !
1106! Packs 4D field data into 1D array to be written into output NetCDF !
1107! file elsewhere. The field data is processed by 3D slices to reduce !
1108! memory requirements. !
1109! !
1110! On Input: !
1111! !
1112! ng Nested grid number (integer) !
1113! model Calling model identifier (integer) !
1114! tile Domain partition (integer) !
1115! gtype Staggered C-grid type (integer) !
1116! ifield Field metadata index (integer) !
1117! tindex Time record index to process (integer) !
1118! LandFill Switch to set fill value in land areas (logical) !
1119! Extract_Flag Extraction flag interpolation/decimation (integer) !
1120! LBi Field I-dimension Lower bound (integer) !
1121! UBi Field I-dimension Upper bound (integer) !
1122! LBj Field J-dimension Lower bound (integer) !
1123! UBj Field J-dimension Upper bound (integer) !
1124! LBk Field K-dimension Lower bound (integer) !
1125! UBk Field K-dimension Upper bound (integer) !
1126! LBt Field fourth-dimension Lower bound (integer) !
1127! UBt Field fourth-dimension Upper bound (integer) !
1128! fourth Fourth dimension index to process (integer) !
1129#ifdef MASKING
1130! Amask Land/Sea mask, if any (real 2D array) !
1131#endif
1132! Ascl Factor to scale field after reading (real) !
1133! Adat 4D field to process (real) !
1134! !
1135! On Output: !
1136! !
1137! start Start index where the first of the data values will !
1138! be written along each dimension (integer) !
1139! total Number of data values to be written along each !
1140! dimension (integer) !
1141! Npts Number of points processed in Awrk. !
1142! Awrk Packed 4D field data (real) !
1143! !
1144!=======================================================================
1145!
1146! Imported variable declarations.
1147!
1148 logical, intent(in) :: LandFill
1149!
1150 integer, intent(in) :: ng, model, tile
1151 integer, intent(in) :: gtype, ifield, tindex
1152 integer, intent(in) :: Extract_Flag
1153 integer, intent(in) :: LBi, UBi, LBj, UBj, LBk, UBk, LBt, UBt
1154 integer, intent(in) :: fourth
1155 integer, intent(out) :: Npts
1156 integer, intent(out) :: start(:), total(:)
1157!
1158 real(dp), intent(in) :: Ascl
1159!
1160#ifdef ASSUMED_SHAPE
1161# ifdef MASKING
1162 real(r8), intent(in) :: Amask(LBi:,LBj:)
1163# endif
1164 real(r8), intent(in) :: Adat(LBi:,LBj:,LBk:,LBt:)
1165#else
1166# ifdef MASKING
1167 real(r8), intent(in) :: Amask(LBi:UBi,LBj:UBj)
1168# endif
1169 real(r8), intent(in) :: Adat(LBi:UBi,LBj:UBj,LBk:UBk,LBt:UBt)
1170#endif
1171 real(r8), intent(out) :: Awrk(:)
1172!
1173! Local variable declarations.
1174!
1175 integer :: i, j, k, ic
1176 integer :: Imin, Imax, Jmin, Jmax, Kmin, Kmax, Koff, Loff
1177 integer :: Isize, Jsize, Ksize, IJsize, MyType
1178!
1179 character (len=*), parameter :: MyFile = &
1180 & __FILE__//", pack_field4d"
1181!
1182!------------------------------------------------------------------------
1183! Pack 4D field data by 3D slices.
1184!------------------------------------------------------------------------
1185!
1186! Set first and last grid point according to staggered C-grid
1187! classification. Set loops offsets.
1188!
1189 mytype=gtype
1190!
1191 SELECT CASE (abs(mytype))
1192!
1193! ---------------------
1194 CASE (p2dvar, p3dvar) ! PSI points field
1195! ---------------------
1196!
1197 imin=iobounds(ng)%ILB_psi
1198 imax=iobounds(ng)%IUB_psi
1199 jmin=iobounds(ng)%JLB_psi
1200 jmax=iobounds(ng)%JUB_psi
1201!
1202! ---------------------
1203 CASE (r2dvar, r3dvar) ! RHO points field
1204! ---------------------
1205!
1206 imin=iobounds(ng)%ILB_rho
1207 imax=iobounds(ng)%IUB_rho
1208 jmin=iobounds(ng)%JLB_rho
1209 jmax=iobounds(ng)%JUB_rho
1210!
1211! ---------------------
1212 CASE (u2dvar, u3dvar) ! U points field
1213! ---------------------
1214!
1215 imin=iobounds(ng)%ILB_u
1216 imax=iobounds(ng)%IUB_u
1217 jmin=iobounds(ng)%JLB_u
1218 jmax=iobounds(ng)%JUB_u
1219!
1220! ---------------------
1221 CASE (v2dvar, v3dvar) ! V points field
1222! ---------------------
1223!
1224 imin=iobounds(ng)%ILB_v
1225 imax=iobounds(ng)%IUB_v
1226 jmin=iobounds(ng)%JLB_v
1227 jmax=iobounds(ng)%JUB_v
1228!
1229! ---------------------
1230 CASE DEFAULT ! RHO points field
1231! ---------------------
1232!
1233 imin=iobounds(ng)%ILB_rho
1234 imax=iobounds(ng)%IUB_rho
1235 jmin=iobounds(ng)%JLB_rho
1236 jmax=iobounds(ng)%JUB_rho
1237 END SELECT
1238!
1239 isize=imax-imin+1
1240 jsize=jmax-jmin+1
1241 ksize=ubk-lbk+1
1242 ijsize=isize*jsize
1243!
1244 IF (lbk.eq.0) THEN
1245 koff=0
1246 ELSE
1247 koff=1
1248 END IF
1249!
1250 IF (lbt.eq.0) THEN
1251 loff=1
1252 ELSE
1253 loff=0
1254 END IF
1255
1256#ifdef DISTRIBUTE
1257!
1258!-----------------------------------------------------------------------
1259! If distributed-memory set-up, collect tile data from all spawned
1260! nodes and store it into a global scratch 1D array, packed in column-
1261! major order.
1262# ifdef MASKING
1263# ifdef WRITE_WATER
1264! Remove land points and pack water points into 1D-array.
1265# else
1266! Overwrite masked points with special value.
1267# endif
1268# endif
1269!-----------------------------------------------------------------------
1270!
1271! Process the data by 3D slices.
1272!
1273 CALL mp_gather3d (ng, model, lbi, ubi, lbj, ubj, lbk, ubk, &
1274 & tindex, gtype, ascl, &
1275# ifdef MASKING
1276 & amask, &
1277# endif
1278 & adat(:,:,:,fourth), npts, awrk, landfill)
1279#else
1280!
1281!-----------------------------------------------------------------------
1282! If serial or shared-memory applications and serial output, pack data
1283! into a global 1D array in column-major order.
1284# ifdef MASKING
1285# ifdef WRITE_WATER
1286! Remove land points and pack water points into 1D-array.
1287# else
1288! Overwrite masked points with special value.
1289# endif
1290# endif
1291!-----------------------------------------------------------------------
1292!
1293! The data is processed in 3D slices.
1294!
1295 IF (gtype.gt.0) THEN
1296 ic=0
1297 npts=ijsize*ksize
1298 DO k=lbk,ubk
1299 DO j=jmin,jmax
1300 DO i=imin,imax
1301 ic=ic+1
1302 awrk(ic)=adat(i,j,k,fourth)*ascl
1303# ifdef MASKING
1304 IF ((amask(i,j).eq.0.0_r8).and.landfill) THEN
1305 awrk(ic)=spval
1306 END IF
1307# endif
1308 END DO
1309 END DO
1310 END DO
1311# ifdef MASKING
1312 ELSE
1313 npts=0
1314 DO k=lbk,ubk
1315 DO j=jmin,jmax
1316 DO i=imin,imax
1317 IF (amask(i,j).gt.0.0_r8) THEN
1318 npts=npts+1
1319 awrk(npts)=adat(i,j,k,fourth)*ascl
1320 END IF
1321 END DO
1322 END DO
1323 END DO
1324# endif
1325 END IF
1326#endif
1327!
1328!-----------------------------------------------------------------------
1329! If there is no extracting data, set the start and total vectors
1330! needed by the NetCDF library for writing.
1331!-----------------------------------------------------------------------
1332!
1333 IF (extract_flag.le.0) THEN
1334 IF (gtype.gt.0) THEN
1335 start(1)=1
1336 total(1)=isize
1337 start(2)=1
1338 total(2)=jsize
1339 start(3)=1
1340 total(3)=ksize
1341 start(4)=fourth+loff
1342 total(4)=1
1343 start(5)=tindex
1344 total(5)=1
1345#ifdef MASKING
1346 ELSE
1347 start(1)=1+(fourth+loff-1)*npts
1348 total(1)=npts
1349 start(2)=tindex
1350 total(2)=1
1351#endif
1352 END IF
1353
1354#ifdef GRID_EXTRACT
1355!
1356!-----------------------------------------------------------------------
1357! Otherwise, extract field by decimation or interpolation.
1358!-----------------------------------------------------------------------
1359!
1360 ELSE IF (extract_flag.ge.1) THEN
1361 CALL extract_field (ng, model, tile, &
1362 & gtype, ifield, tindex, &
1363 & extract_flag, &
1364 & imin, imax, jmin, jmax, lbk, ubk, &
1365 & start, total, &
1366 & npts, awrk)
1367#endif
1368 END IF
1369
1370#ifdef POSITIVE_ZERO
1371!
1372!-----------------------------------------------------------------------
1373! Impose positive zero. Recall that IEEE allows both negative and
1374! positive zeros, which affects bit-by-bit comparison betwee different
1375! tile partitions.
1376!-----------------------------------------------------------------------
1377!
1378 DO ic=1,npts
1379 IF (abs(awrk(ic)).eq.0.0_r8) THEN
1380 awrk(ic)=0.0_r8
1381 END IF
1382 END DO
1383#endif
1384!
1385 RETURN

References mod_param::iobounds, distribute_mod::mp_gather3d(), mod_param::p2dvar, mod_param::p3dvar, mod_param::r2dvar, mod_param::r3dvar, mod_scalars::spval, mod_param::u2dvar, mod_param::u3dvar, mod_param::v2dvar, and mod_param::v3dvar.

Referenced by nf_fwrite4d_mod::nf_fwrite4d::nf90_fwrite4d().

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