ROMS
Loading...
Searching...
No Matches
extract_field.F
Go to the documentation of this file.
1#include "cppdefs.h"
3
4#ifdef GRID_EXTRACT
5!
6!git $Id$
7!=======================================================================
8! Copyright (c) 2002-2025 The ROMS Group !
9! Licensed under a MIT/X style license !
10! See License_ROMS.md Hernan G. Arango !
11!=======================================================================
12! !
13! These routines extract output fields at coaser resolution using !
14! decimation or horizontal interpolation. Both extraction strategies !
15! required a grid NetCDF (XTRNAME) provided at execution. !
16! !
17! If Extract_Flag > 1: Decimation !
18! !
19! It decimates the field solution at the prescribed integer factor, !
20! Extract_Flag . For example, if Extract_Flag=2 (recommended), the !
21! output fieldsare written at every other point, resulting in coarser !
22! data resolution. This strategy is advantageous in mixed resolution, !
23! split 4D-Var applications where the outer loop background (prior) !
24! trajectory may be computed at a higher resolution than in the inner !
25! loop minimization to accelerate the calculations. For decimation to !
26! work, the number of parent grid RHO-points (0: Lm+1, 0:Mm+1) must !
27! be multiples of dec_fator. That is, !
28! !
29! MOD(Lm+1, Extract_Flag) = 0 !
30! MOD(Mm+1, Extract_Flag) = 0 !
31! !
32! If Extract_Flag = 1: Interpolation !
33! !
34!=======================================================================
35!
36 USE mod_param
37 USE mod_parallel
38 USE mod_extract
39 USE mod_grid
40 USE mod_iounits
41 USE mod_ncparam
42 USE mod_scalars
43!
44# ifdef DISTRIBUTE
45 USE distribute_mod, ONLY : mp_gather2d, mp_gather2d_xtr
46# endif
47 USE roms_interpolate_mod, ONLY : cinterp2d, linterp2d, hindices
48 USE strings_mod, ONLY : founderror
49!
50 implicit none
51!
52! Interfaces for sam name routine overloading.
53!
54# ifdef ADJUST_BOUNDARY
55 INTERFACE extract_boundary
56 MODULE PROCEDURE extract_boundary2d
57 MODULE PROCEDURE extract_boundary3d
58 END INTERFACE extract_boundary
59!
60# endif
61 INTERFACE extract_field
62 MODULE PROCEDURE extract_field2d
63 MODULE PROCEDURE extract_field3d
64 MODULE PROCEDURE extract_field4d
65 END INTERFACE extract_field
66!
67 PUBLIC :: interp_coords
68!
69# ifdef ADJUST_BOUNDARY
70 PRIVATE :: average_boundary2d
71 PRIVATE :: average_boundary3d
72# endif
73 PRIVATE :: average_field2d
74 PRIVATE :: average_field3d
75# ifdef ADJUST_BOUNDARY
76 PRIVATE :: decimate_boundary2d
77 PRIVATE :: decimate_boundary3d
78# endif
79 PRIVATE :: decimate_field2d
80 PRIVATE :: decimate_field3d
81 PRIVATE :: decimate_field4d
82 PRIVATE :: interp_field2d
83 PRIVATE :: interp_field3d
84 PRIVATE :: interp_field4d
85 PRIVATE :: regrid_field2d
86 PRIVATE :: regrid_field3d
87!
88! Module parameters.
89!
90 logical, parameter :: ParallelOutput = .false.
91!
92 integer, parameter :: Bilinear = 0 ! bilinear interpolation
93 integer, parameter :: Bicubic = 1 ! bicubic interpolation
94!
95 CONTAINS
96!
97# ifdef ADJUST_BOUNDARY
98!
99 SUBROUTINE extract_boundary2d (ng, model, tile, &
100 & gtype, ncvname, tindex, &
101 & Extract_flag, &
102 & Imin, Imax, Jmin, Jmax, &
103 & Nrec, &
104 & start, total, &
105 & Npts, Bdat)
106!
107!=======================================================================
108! !
109! Extracts/samples boundary data from donor 2D array. !
110! !
111! On Input: !
112! !
113! ng Nested grid number (integer) !
114! model Calling model identifier (integer) !
115! gtype Staggered C-grid type (integer) !
116! ncvname NetCDF variable name (string) !
117! tindex Time record index to process (integer) !
118! Extract_Flag Extraction/decimation flag (integer) !
119! Imin Donor boundary starting data I-index (integer) !
120! Imax Donor boundary ending data I-index (integer) !
121! Jmin Donor boundary starting data J-index (integer) !
122! Jmax Donor boundary ending data J-index (integer) !
123! Nrec Number of boundary records (integer) !
124! Npts Length of packed boundary data (integer) !
125! Bdat Packed global donor 2D boundary data (real, 1D) !
126! !
127! On Output: !
128! !
129! start Start index where the first of the data values will !
130! be written along each dimension (integer) !
131! total Number of data values to be written along each !
132! dimension (integer) !
133! Bdat Extracted 2D boundary data (real, 1D) !
134! !
135!=======================================================================
136!
137! Imported variable declarations.
138!
139 integer, intent(in) :: ng, model, tile
140 integer, intent(in) :: gtype, tindex, Extract_flag
141 integer, intent(in) :: Imin, Imax, Jmin, Jmax
142 integer, intent(in) :: Nrec
143 integer, intent(in) :: Npts
144 integer, intent(out) :: start(:), total(:)
145!
146 real(r8), intent(inout) :: Bdat(:)
147!
148 character (len=*), intent(in) :: ncvname
149!
150! Local variable declarations.
151!
152 character (len=*), parameter :: MyFile = &
153 & __FILE__//", extract_boundary2d"
154!
155 sourcefile=myfile
156!
157!------------------------------------------------------------------------
158! Decimate or interpolate 2D boundary data.
159!------------------------------------------------------------------------
160!
161 IF (extract_flag.ge.2) THEN
162 CALL decimate_boundary2d (ng, model, tile, &
163 & gtype, ncvname, tindex, extract_flag, &
164 & imin, imax, jmin, jmax, &
165 & nrec, &
166 & start, total, &
167 & npts, bdat)
168 ELSE IF (extract_flag.eq.1) THEN
169 END IF
170!
171 RETURN
172 END SUBROUTINE extract_boundary2d
173!
174 SUBROUTINE extract_boundary3d (ng, model, tile, &
175 & gtype, ncvname, tindex, &
176 & Extract_flag, &
177 & Imin, Imax, Jmin, Jmax, LBk, UBk, &
178 & Nrec, &
179 & start, total, &
180 & Npts, Bdat)
181!
182!=======================================================================
183! !
184! Extracts/samples boundary data from donor 2D array. !
185! !
186! On Input: !
187! !
188! ng Nested grid number (integer) !
189! model Calling model identifier (integer) !
190! gtype Staggered C-grid type (integer) !
191! ncvname NetCDF variable name (string) !
192! tindex Time record index to process (integer) !
193! Extract_Flag Extraction/decimation flag (integer) !
194! Imin Donor boundary starting data I-index (integer) !
195! Imax Donor boundary ending data I-index (integer) !
196! Jmin Donor boundary starting data J-index (integer) !
197! Jmax Donor boundary ending data J-index (integer) !
198! LBk Donor boundary K-dimension lower bound (integer) !
199! UBk Donor boundary K-dimension upper bound (integer) !
200! Nrec Number of boundary records (integer) !
201! Npts Length of packed boundary data (integer) !
202! Bdat Packed global donor 2D boundary data (real, 1D) !
203! !
204! On Output: !
205! !
206! start Start index where the first of the data values will !
207! be written along each dimension (integer) !
208! total Number of data values to be written along each !
209! dimension (integer) !
210! Bdat Extracted 3D boundary data (real, 1D) !
211! !
212!=======================================================================
213!
214! Imported variable declarations.
215!
216 integer, intent(in) :: ng, model, tile
217 integer, intent(in) :: gtype, tindex, Extract_flag
218 integer, intent(in) :: Imin, Imax, Jmin, Jmax, LBk, UBk
219 integer, intent(in) :: Nrec
220 integer, intent(in) :: Npts
221 integer, intent(out) :: start(:), total(:)
222!
223 real(r8), intent(inout) :: Bdat(:)
224!
225 character (len=*), intent(in) :: ncvname
226!
227! Local variable declarations.
228!
229 integer :: ij
230!
231 character (len=*), parameter :: MyFile = &
232 & __FILE__//", extract_boundary3d"
233!
234!------------------------------------------------------------------------
235! Decimate or interpolate 2D boundary data.
236!------------------------------------------------------------------------
237!
238 IF (extract_flag.ge.2) THEN
239 CALL decimate_boundary3d (ng, model, tile, &
240 & gtype, ncvname, tindex, extract_flag, &
241 & imin, imax, jmin, jmax, lbk, ubk, &
242 & nrec, &
243 & start, total, &
244 & npts, bdat)
245 ELSE IF (extract_flag.eq.1) THEN
246 END IF
247!
248 RETURN
249 END SUBROUTINE extract_boundary3d
250# endif
251!
252 SUBROUTINE extract_field2d (ng, model, tile, &
253 & gtype, ifield, tindex, Extract_Flag, &
254 & Imin, Imax, Jmin, Jmax, &
255 & start, total, &
256 & Npts, Fdat)
257!
258!=======================================================================
259! !
260! It decimates or interpolates data from donor 2D fields, packed as !
261! a 1D global array, to specified extract grid geometry. !
262! !
263! On Input: !
264! !
265! ng Nested grid number (integer) !
266! model Calling model identifier (integer) !
267! tile Domain partition (integer) !
268! gtype Staggered C-grid type (integer) !
269! ifield Field metadata index (integer) !
270! tindex Time record index to process (integer) !
271! Extract_Flag Extraction/decimation flag (integer) !
272! Imin Donor field starting data I-index (integer) !
273! Imax Donor field ending data I-index (integer) !
274! Jmin Donor field starting data J-index (integer) !
275! Jmax Donor field ending data J-index (integer) !
276! Npts Length of packed field data (integer)
277! Fdat Packed global donor 2D field data (real 1D array) !
278! !
279! On Output: !
280! !
281! start Start index where the first of the data values will !
282! be written along each dimension (integer) !
283! total Number of data values to be written along each !
284! dimension (integer) !
285! Fdat Extracted 2D field data (real) !
286! !
287!=======================================================================
288!
289! Imported variable declarations.
290!
291 integer, intent(in) :: ng, model, tile
292 integer, intent(in) :: gtype, ifield, tindex, Extract_Flag
293 integer, intent(in) :: Imin, Imax, Jmin, Jmax
294 integer, intent(in) :: Npts
295 integer, intent(out) :: start(:), total(:)
296!
297 real(r8), intent(inout) :: Fdat(:)
298!
299! Local variable declarations.
300!
301 character (len=*), parameter :: MyFile = &
302 & __FILE__//", extract_field2d"
303!
304 sourcefile=myfile
305!
306!------------------------------------------------------------------------
307! Decimate or interpolate 2D field.
308!------------------------------------------------------------------------
309!
310 IF (extract_flag.ge.2) THEN
311 CALL decimate_field2d (ng, model, tile, &
312 & gtype, ifield, tindex, extract_flag, &
313 & imin, imax, jmin, jmax, &
314 & start, total, &
315 & npts, fdat)
316 ELSE IF (extract_flag.eq.1) THEN
317 IF (paralleloutput) THEN
318 CALL interp_field2d (ng, model, tile, &
319 & gtype, ifield, tindex, &
320 & imin, imax, jmin, jmax, &
321 & start, total, &
322 & npts, fdat)
323 ELSE
324 CALL interp_field2d_global (ng, model, tile, &
325 & gtype, ifield, tindex, &
326 & imin, imax, jmin, jmax, &
327 & start, total, &
328 & npts, fdat)
329 END IF
330 END IF
331!
332 RETURN
333 END SUBROUTINE extract_field2d
334!
335 SUBROUTINE extract_field3d (ng, model, tile, &
336 & gtype, ifield, tindex, Extract_Flag, &
337 & Imin, Imax, Jmin, Jmax, Kmin, Kmax, &
338 & start, total, &
339 & Npts, Fdat)
340!
341!=======================================================================
342! !
343! It decimates or interpolates data from donor 3D fields, packed as !
344! a 1D global array, to specified extract grid geometry. !
345! !
346! On Input: !
347! !
348! ng Nested grid number (integer) !
349! model Calling model identifier (integer) !
350! tile Domain partition (integer) !
351! gtype Staggered C-grid type (integer) !
352! ifield Field metadata index (integer) !
353! tindex Time record index to process (integer) !
354! Extract_Flag Extraction/decimation flag (integer) !
355! Imin Donor field starting data I-index (integer) !
356! Imax Donor field ending data I-index (integer) !
357! Jmin Donor field starting data J-index (integer) !
358! Jmax Donor field ending data J-index (integer) !
359! Kmin Donor field starting data K-index (integer) !
360! Kmax Donor field ending data K-index (integer) !
361! Fdat Packed global donor 2D field data (real 1D array) !
362! !
363! On Output: !
364! !
365! start Start index where the first of the data values will !
366! be written along each dimension (integer) !
367! total Number of data values to be written along each !
368! dimension (integer) !
369! Fdat Extracted 2D field data (real) !
370! !
371!=======================================================================
372!
373! Imported variable declarations.
374!
375 integer, intent(in) :: ng, model, tile
376 integer, intent(in) :: gtype, ifield, tindex, Extract_Flag
377 integer, intent(in) :: Imin, Imax, Jmin, Jmax, Kmin, Kmax
378 integer, intent(in) :: Npts
379 integer, intent(out) :: start(:), total(:)
380!
381 real(r8), intent(inout) :: Fdat(:)
382!
383! Local variable declarations.
384!
385 character (len=*), parameter :: MyFile = &
386 & __FILE__//", extract_field3d"
387!
388 sourcefile=myfile
389!
390!------------------------------------------------------------------------
391! Decimate or interpolate 3D field.
392!------------------------------------------------------------------------
393!
394 IF (extract_flag.ge.2) THEN
395 CALL decimate_field3d (ng, model, tile, &
396 & gtype, ifield, tindex, extract_flag, &
397 & imin, imax, jmin, jmax, &
398 & kmin, kmax, &
399 & start, total, &
400 & npts, fdat)
401 ELSE IF (extract_flag.eq.1) THEN
402 IF (paralleloutput) THEN
403 CALL interp_field3d (ng, model, tile, &
404 & gtype, ifield, tindex, &
405 & imin, imax, jmin, jmax, &
406 & kmin, kmax, &
407 & start, total, &
408 & npts, fdat)
409 ELSE
410 CALL interp_field3d_global (ng, model, tile, &
411 & gtype, ifield, tindex, &
412 & imin, imax, jmin, jmax, &
413 & kmin, kmax, &
414 & start, total, &
415 & npts, fdat)
416 END IF
417 END IF
418!
419 RETURN
420 END SUBROUTINE extract_field3d
421!
422 SUBROUTINE extract_field4d (ng, model, tile, &
423 & gtype, ifield, tindex, Extract_Flag, &
424 & Imin, Imax, Jmin, Jmax, Kmin, Kmax, &
425 & Fourth, Loff, &
426 & start, total, &
427 & Npts, Fdat)
428!
429!=======================================================================
430! !
431! It decimates or interpolates data from donor 4D fields, packed as !
432! a 1D global array, to specified extract grid geometry. The field !
433! data is processed by 3D slices to reduce memory requirements in !
434! the calling routine. !
435! !
436! On Input: !
437! !
438! ng Nested grid number (integer) !
439! model Calling model identifier (integer) !
440! tile Domain partition (integer) !
441! gtype Staggered C-grid type (integer) !
442! ifield Field metadata index (integer) !
443! tindex Time record index to process (integer) !
444! Extract_Flag Extraction/decimation flag (integer) !
445! Imin Donor field starting data I-index (integer) !
446! Imax Donor field ending data I-index (integer) !
447! Jmin Donor field starting data J-index (integer) !
448! Jmax Donor field ending data J-index (integer) !
449! Kmin Donor field starting data K-index (integer) !
450! Kmax Donor field ending data K-index (integer) !
451! fourth Donor Fouth dimension index to process (integer) !
452! Loff Fourth dimension couter offset (integer) !
453! Fdat Packed global donor 3D field data (real 1D array) !
454! !
455! On Output: !
456! !
457! start Start index where the first of the data values will !
458! be written along each dimension (integer) !
459! total Number of data values to be written along each !
460! dimension (integer) !
461! Npts Number of points processed in Fdat (integer) !
462! Fdat Extracted 3D field data (real) !
463! !
464!=======================================================================
465!
466! Imported variable declarations.
467!
468 integer, intent(in) :: ng, model, tile
469 integer, intent(in) :: gtype, ifield, tindex, Extract_Flag
470 integer, intent(in) :: Imin, Imax, Jmin, Jmax, Kmin, Kmax
471 integer, intent(in) :: fourth, Loff
472 integer, intent(in) :: Npts
473 integer, intent(out) :: start(:), total(:)
474!
475 real(r8), intent(inout) :: Fdat(:)
476!
477! Local variable declarations.
478!
479 character (len=*), parameter :: MyFile = &
480 & __FILE__//", extract_field4d"
481!
482 sourcefile=myfile
483!
484!------------------------------------------------------------------------
485! Decimate or interpolate 4D field.
486!------------------------------------------------------------------------
487!
488 IF (extract_flag.ge.2) THEN
489 CALL decimate_field4d (ng, model, tile, &
490 & gtype, ifield, tindex, extract_flag, &
491 & imin, imax, jmin, jmax, &
492 & kmin, kmax, fourth, loff, &
493 & start, total, &
494 & npts, fdat)
495 ELSE IF (extract_flag.eq.1) THEN
496 IF (paralleloutput) THEN
497 CALL interp_field4d (ng, model, tile, &
498 & gtype, ifield, tindex, &
499 & imin, imax, jmin, jmax, &
500 & kmin, kmax, fourth, loff, &
501 & start, total, &
502 & npts, fdat)
503 ELSE
504 CALL interp_field4d_global (ng, model, tile, &
505 & gtype, ifield, tindex, &
506 & imin, imax, jmin, jmax, &
507 & kmin, kmax, fourth, loff, &
508 & start, total, &
509 & npts, fdat)
510 END IF
511 END IF
512!
513 RETURN
514 END SUBROUTINE extract_field4d
515
516# ifdef ADJUST_BOUNDARY
517!
518 SUBROUTINE average_boundary2d (ng, model, tile, &
519 & gtype, ncvname, Extract_Flag, &
520 & Imin, Imax, Jmin, Jmax, Nrec, &
521 & Npts, Bdat, Bavg)
522!
523!=======================================================================
524! !
525! It decimates data from donor 2D boundary packed as 1D global array. !
526! !
527! On Input: !
528! !
529! ng Nested grid number (integer) !
530! model Calling model identifier (integer) !
531! tile Domain partition (integer) !
532! gtype Staggered C-grid type (integer) !
533! ncvname NetCDF variable name (string) !
534! Extract_Flag Extraction/decimation flag (integer) !
535! Imin Donor boundary starting data I-index (integer) !
536! Imax Donor boundary ending data I-index (integer) !
537! Jmin Donor boundary starting data J-index (integer) !
538! Jmax Donor boundary ending data J-index (integer) !
539! Nrec Number of boundary records (integer) !
540! Npts Length of packed boundary data (integer) !
541! Bdat Packed global donor 2D boundary data (real; 1D) !
542! !
543! On Output: !
544! !
545! Favg Averaged 2D boundary data (real 1D array) !
546! !
547!=======================================================================
548!
549! Imported variable declarations.
550!
551 integer, intent(in) :: ng, model, tile
552 integer, intent(in) :: gtype, Extract_Flag
553 integer, intent(in) :: Imin, Imax, Jmin, Jmax, Nrec
554 integer, intent(in) :: Npts
555!
556 real(r8), intent(in) :: Bdat(:)
557 real(r8), intent(out) :: Bavg(:)
558!
559 character (len=*), intent(in) :: ncvname
560!
561! Local variable declarations.
562!
563 logical, dimension(4) :: bounded
564!
565 integer :: bc, i, ib, ij, ir, j, rc
566 integer :: IorJ
567!
568 character (len=*), parameter :: MyFile = &
569 & __FILE__//", average_boundary2d"
570!
571 sourcefile=myfile
572!
573!------------------------------------------------------------------------
574! Average 2D U- and V-boundary data at the appropiate location.
575!------------------------------------------------------------------------
576!
577 iorj=iobounds(ng)%IorJ
578!
579! Set switch to process boundary data by their associated tiles.
580!
581 bounded(iwest )=domain(ng)%Western_Edge (tile)
582 bounded(ieast )=domain(ng)%Eastern_Edge (tile)
583 bounded(isouth)=domain(ng)%Southern_Edge(tile)
584 bounded(inorth)=domain(ng)%Northern_Edge(tile)
585!
586! The logic below works well for Extract_Flag=2.
587!
588 IF (extract_flag.ge.2) THEN
589!
590 SELECT CASE (gtype)
591!
592! ---------------------
593 CASE (u2dvar, u3dvar) ! U points field:
594! --------------------- Average data at RHO-points
595! for southern and northern
596 bavg=bdat ! boundaries only
597 DO ir=1,nrec
598 rc=(ir-1)*iorj*4
599 DO ib=1,4
600 bc=(ib-1)*iorj+rc
601 IF (bounded(ib).and. &
602 & ((ib.eq.isouth).or.(ib.eq.inorth))) THEN
603 DO i=imin,imax-1
604 ij=i+bc
605 bavg(ij)=0.5_r8*(bdat(ij)+bdat(ij+1))
606 END DO
607 END IF
608 END DO
609 END DO
610!
611! ---------------------
612 CASE (v2dvar, v3dvar) ! V points field:
613! --------------------- Average data at RHO-points
614! for western and eastern
615 bavg=bdat ! boundaries only
616 DO ir=1,nrec
617 rc=(ir-1)*iorj*4
618 DO ib=1,4
619 bc=(ib-1)*iorj+rc
620 IF (bounded(ib).and. &
621 & ((ib.eq.iwest).or.(ib.eq.ieast))) THEN
622 DO j=jmin,jmax-1
623 ij=j+bc
624 bavg(ij)=0.5_r8*(bdat(ij)+bdat(ij+1))
625 END DO
626 END IF
627 END DO
628 END DO
629!
630 END SELECT
631 END IF
632!
633 RETURN
634 END SUBROUTINE average_boundary2d
635!
636 SUBROUTINE average_boundary3d (ng, model, tile, &
637 & gtype, ncvname, Extract_Flag, &
638 & Imin, Imax, Jmin, Jmax, LBk, UBk, &
639 & Nrec, &
640 & Npts, Bdat, Bavg)
641!
642!=======================================================================
643! !
644! It decimates data from donor 3D boundary packed as 1D global array. !
645! !
646! On Input: !
647! !
648! ng Nested grid number (integer) !
649! model Calling model identifier (integer) !
650! tile Domain partition (integer) !
651! gtype Staggered C-grid type (integer) !
652! ncvname NetCDF variable name (string) !
653! Extract_Flag Extraction/decimation flag (integer) !
654! Imin Donor boundary starting data I-index (integer) !
655! Imax Donor boundary ending data I-index (integer) !
656! Jmin Donor boundary starting data J-index (integer) !
657! Jmax Donor boundary ending data J-index (integer) !
658! LBk Donor boundary K-dimension lower bound (integer) !
659! UBk Donor boundary K-dimension upper bound (integer) !
660! Nrec Number of boundary records (integer) !
661! Npts Length of packed boundary data (integer) !
662! Bdat Packed global donor 2D boundary data (real; 1D) !
663! !
664! On Output: !
665! !
666! Favg Averaged 3D boundary data (real 1D array) !
667! !
668!=======================================================================
669!
670! Imported variable declarations.
671!
672 integer, intent(in) :: ng, model, tile
673 integer, intent(in) :: gtype, Extract_Flag
674 integer, intent(in) :: Imin, Imax, Jmin, Jmax, LBk, UBk, Nrec
675 integer, intent(in) :: Npts
676!
677 real(r8), intent(in) :: Bdat(:)
678 real(r8), intent(out) :: Bavg(:)
679!
680 character (len=*), intent(in) :: ncvname
681!
682! Local variable declarations.
683!
684 logical, dimension(4) :: bounded
685!
686 integer :: bc, i, ib, ij, ir, j, k, kc, rc
687 integer :: IJKlen, IorJ, Klen
688!
689 character (len=*), parameter :: MyFile = &
690 & __FILE__//", average_boundary3d"
691!
692 sourcefile=myfile
693!
694!------------------------------------------------------------------------
695! Average 3D U- and V-boundary data at the appropiate location.
696!------------------------------------------------------------------------
697!
698 iorj=iobounds(ng)%IorJ
699!
700! Set switch to process boundary data by their associated tiles.
701!
702 bounded(iwest )=domain(ng)%Western_Edge (tile)
703 bounded(ieast )=domain(ng)%Eastern_Edge (tile)
704 bounded(isouth)=domain(ng)%Southern_Edge(tile)
705 bounded(inorth)=domain(ng)%Northern_Edge(tile)
706!
707! The logic below works well for Extract_Flag=2.
708!
709 IF (extract_flag.ge.2) THEN
710 klen=ubk-lbk+1
711 ijklen=iorj*klen
712!
713 SELECT CASE (gtype)
714!
715! ---------------------
716 CASE (u3dvar) ! U points field:
717! --------------------- Average data at RHO-points
718! for southern and northern
719 bavg=bdat ! boundaries only
720 DO ir=1,nrec
721 rc=(ir-1)*iorj*4
722 DO ib=1,4
723 bc=(ib-1)*ijklen+rc
724 IF (bounded(ib).and. &
725 & ((ib.eq.isouth).or.(ib.eq.inorth))) THEN
726 DO k=lbk,ubk
727 kc=(k-lbk)*iorj+bc
728 DO i=imin,imax-1
729 ij=i+kc
730 bavg(ij)=0.5_r8*(bdat(ij)+bdat(ij+1))
731 END DO
732 END DO
733 END IF
734 END DO
735 END DO
736!
737! ---------------------
738 CASE (v2dvar, v3dvar) ! V points field:
739! --------------------- Average data at RHO-points
740! for western and eastern
741 bavg=bdat ! boundaries only
742 DO ir=1,nrec
743 rc=(ir-1)*iorj*4
744 DO ib=1,4
745 bc=(ib-1)*ijklen+rc
746 IF (bounded(ib).and. &
747 & ((ib.eq.iwest).or.(ib.eq.ieast))) THEN
748 DO k=lbk,ubk
749 kc=(k-lbk)*iorj+bc
750 DO j=jmin,jmax-1
751 ij=j+kc
752 bavg(ij)=0.5_r8*(bdat(ij)+bdat(ij+1))
753 END DO
754 END DO
755 END IF
756 END DO
757 END DO
758!
759 END SELECT
760 END IF
761!
762 RETURN
763 END SUBROUTINE average_boundary3d
764#endif
765!
766 SUBROUTINE average_field2d (ng, model, tile, &
767 & gtype, ifield, Extract_Flag, &
768 & Imin, Imax, Jmin, Jmax, &
769 & Npts, Fdat, Favg)
770!
771!=======================================================================
772! !
773! It averages staggered 2D U- and V-fields at the appropiate location !
774! to facilitate decimation. !
775! !
776! On Input: !
777! !
778! ng Nested grid number (integer) !
779! model Calling model identifier (integer) !
780! tile Domain partition (integer) !
781! gtype Staggered C-grid type (integer) !
782! ifield Field metadata index (integer) !
783! Extract_Flag Extraction/decimation flag (integer) !
784! Imin Donor field starting data I-index (integer) !
785! Imax Donor field ending data I-index (integer) !
786! Jmin Donor field starting data J-index (integer) !
787! Jmax Donor field ending data J-index (integer) !
788! Fdat Packed global donor 2D field data (real 1D array) !
789! !
790! On Output: !
791! !
792! Favg Averaged 2D field data (real 1D array) !
793! !
794!=======================================================================
795!
796! Imported variable declarations.
797!
798 integer, intent(in) :: ng, model, tile
799 integer, intent(in) :: gtype, ifield, Extract_Flag
800 integer, intent(in) :: Imin, Imax, Jmin, Jmax
801 integer, intent(in) :: Npts
802!
803 real(r8), intent(in) :: Fdat(:)
804 real(r8), intent(out) :: Favg(:)
805!
806! Local variable declarations.
807!
808 integer :: i, j, ij
809 integer :: Ioff, Joff, Isize, Jsize
810!
811 character (len=*), parameter :: MyFile = &
812 & __FILE__//", average_field2d"
813!
814!------------------------------------------------------------------------
815! Average 2D U- and V-field data at the appropiate location.
816!------------------------------------------------------------------------
817!
818 IF (extract_flag.ge.2) THEN
819 isize=imax-imin+1
820 jsize=jmax-jmin+1
821!
822 SELECT CASE (gtype)
823!
824! ---------------------
825 CASE (u2dvar, u3dvar) ! U points field
826! ---------------------
827!
828 favg=fdat
829 ioff=0
830 joff=1
831 DO j=jmin,jmax ! Average data at RHO points
832 DO i=imin,imax-1
833 ij=i+ioff+(j-1+joff)*isize
834 favg(ij)=0.5_r8*(fdat(ij)+fdat(ij+1))
835 END DO
836 END DO
837!
838! ---------------------
839 CASE (v2dvar, v3dvar) ! V points field
840! ---------------------
841!
842 favg=fdat
843 ioff=1
844 joff=0
845 DO j=jmin,jmax-1 ! Average data at RHO points
846 DO i=imin,imax
847 ij=i+ioff+(j-1+joff)*isize
848 favg(ij)=0.5_r8*(fdat(ij)+fdat(ij+isize))
849 END DO
850 END DO
851!
852 END SELECT
853 END IF
854!
855 RETURN
856 END SUBROUTINE average_field2d
857!
858 SUBROUTINE average_field3d (ng, model, tile, &
859 & gtype, ifield, Extract_Flag, &
860 & Imin, Imax, Jmin, Jmax, Kmin, Kmax, &
861 & Npts, Fdat, Favg)
862!
863!=======================================================================
864! !
865! It averages staggered 3D U- and V-fields at the appropiate location !
866! to facilitate decimation. !
867! !
868! On Input: !
869! !
870! ng Nested grid number (integer) !
871! model Calling model identifier (integer) !
872! tile Domain partition (integer) !
873! gtype Staggered C-grid type (integer) !
874! ifield Field metadata index (integer) !
875! Extract_Flag Extraction/decimation flag (integer) !
876! Imin Donor field starting data I-index (integer) !
877! Imax Donor field ending data I-index (integer) !
878! Jmin Donor field starting data J-index (integer) !
879! Jmax Donor field ending data J-index (integer) !
880! Kmin Donor field starting data K-index (integer) !
881! Kmax Donor field ending data K-index (integer) !
882! Fdat Packed global donor 2D field data (real 1D array) !
883! !
884! On Output: !
885! !
886! Favg Averaged 2D field data (real 1D array) !
887! !
888!=======================================================================
889!
890! Imported variable declarations.
891!
892 integer, intent(in) :: ng, model, tile
893 integer, intent(in) :: gtype, ifield, Extract_Flag
894 integer, intent(in) :: Imin, Imax, Jmin, Jmax, Kmin, Kmax
895 integer, intent(in) :: Npts
896!
897 real(r8), intent(in) :: Fdat(:)
898 real(r8), intent(out) :: Favg(:)
899!
900! Local variable declarations.
901!
902 integer :: i, j, k, ij, ijk
903 integer :: Ioff, Joff, Koff
904 integer :: Isize, Jsize, Ksize, IJsize
905!
906 character (len=*), parameter :: MyFile = &
907 & __FILE__//", average_field3d"
908!
909!------------------------------------------------------------------------
910! Average 3D U- and V-field data at the appropiate location.
911!------------------------------------------------------------------------
912!
913 IF (extract_flag.ge.2) THEN
914 isize=imax-imin+1
915 jsize=jmax-jmin+1
916 ksize=kmax-kmin+1
917 ijsize=isize*jsize
918!
919 SELECT CASE (abs(gtype))
920!
921! -------------
922 CASE (u3dvar) ! U points field
923! -------------
924!
925 favg=fdat
926 ioff=0
927 joff=1
928 koff=0
929 DO k=kmin,kmax ! Average data at RHO points
930 DO j=jmin,jmax
931 DO i=imin,imax-1
932 ij=i+ioff+(j-1+joff)*isize
933 ijk=ij+(k-1+koff)*ijsize
934 favg(ijk)=0.5_r8*(fdat(ijk)+fdat(ijk+1))
935 END DO
936 END DO
937 END DO
938!
939! -------------
940 CASE (v3dvar) ! V points field
941! -------------
942!
943 favg=fdat
944 ioff=1
945 joff=0
946 koff=0
947 DO k=kmin,kmax ! Average data at RHO points
948 DO j=jmin,jmax-1
949 DO i=imin,imax
950 ij=i+ioff+(j-1+joff)*isize
951 ijk=ij+(k-1+koff)*ijsize
952 favg(ijk)=0.5_r8*(fdat(ijk)+fdat(ijk+isize))
953 END DO
954 END DO
955 END DO
956!
957 END SELECT
958 END IF
959!
960 RETURN
961 END SUBROUTINE average_field3d
962
963# ifdef ADJUST_BOUNDARY
964!
965 SUBROUTINE decimate_boundary2d (ng, model, tile, &
966 & gtype, ncvname, tindex, &
967 & Extract_Flag, &
968 & Imin, Imax, Jmin, Jmax, &
969 & Nrec, &
970 & start, total, &
971 & Npts, Bdat)
972!
973!=======================================================================
974! !
975! It decimates data from donor 2D boundary packed as 1D global array. !
976! !
977! On Input: !
978! !
979! ng Nested grid number (integer) !
980! model Calling model identifier (integer) !
981! tile Domain partition (integer) !
982! gtype Staggered C-grid type (integer) !
983! ncvname NetCDF variable name (string) !
984! tindex Time record index to process (integer) !
985! Extract_Flag Extraction/decimation flag (integer) !
986! Imin Donor boundary starting data I-index (integer) !
987! Imax Donor boundary ending data I-index (integer) !
988! Jmin Donor boundary starting data J-index (integer) !
989! Jmax Donor boundary ending data J-index (integer) !
990! Nrec Number of boundary records (integer) !
991! Npts Length of packed boundary data (integer) !
992! Bdat Packed global donor 2D boundary data (real; 1D) !
993! !
994! On Output: !
995! !
996! start Start index where the first of the data values will !
997! be written along each dimension (integer) !
998! total Number of data values to be written along each !
999! dimension (integer) !
1000! Bdat Decimated 2D boundary data (real; 1D) !
1001! !
1002!=======================================================================
1003!
1004! Imported variable declarations.
1005!
1006 integer, intent(in) :: ng, model, tile
1007 integer, intent(in) :: gtype, tindex, Extract_Flag
1008 integer, intent(in) :: Imin, Imax, Jmin, Jmax
1009 integer, intent(in) :: Nrec, Npts
1010 integer, intent(out) :: start(:), total(:)
1011!
1012 real(r8), intent(inout) :: Bdat(:)
1013!
1014 character (len=*), intent(in) :: ncvname
1015!
1016! Local variable declarations.
1017!
1018 integer :: bc, ib, ic, ij, ir, ifactor, rc
1019 integer :: Idim, IJdim, IorJ, Mpts
1020!
1021 real(r8) :: Bwrk(SIZE(Bdat))
1022!
1023 character (len=*), parameter :: MyFile = &
1024 & __FILE__//", decimate_boundary2d"
1025!
1026 sourcefile=myfile
1027!
1028!------------------------------------------------------------------------
1029! Decimate input packed 2D boundary by requested factor.
1030!------------------------------------------------------------------------
1031!
1032 ifactor=abs(extract_flag)
1033 iorj=iobounds(ng)%IorJ
1034!
1035 SELECT CASE (gtype)
1036!
1037! ---------------------
1038 CASE (r2dvar, r3dvar) ! RHO points field
1039! ---------------------
1040!
1041 bwrk=bdat
1042 bdat=0.0_r8
1043 ic=0
1044 DO ir=1,nrec
1045 rc=(ir-1)*iorj*4
1046 DO ib=1,4
1047 bc=(ib-1)*iorj+rc
1048 DO ij=1,iorj,ifactor
1049 ic=ic+1
1050 bdat(ic)=bwrk(ij)
1051 END DO
1052 IF ((ir.eq.1).and.(ib.eq.1)) ijdim=ic
1053 END DO
1054 END DO
1055 mpts=ic
1056!
1057! ---------------------
1058 CASE (u2dvar, u3dvar) ! U points field
1059! ---------------------
1060!
1061 CALL average_boundary2d (ng, model, tile, &
1062 & gtype, ncvname, extract_flag, &
1063 & imin, imax, jmin, jmax, nrec, &
1064 & npts, bdat, bwrk)
1065!
1066 bdat=0.0_r8
1067 ic=0
1068 DO ir=1,nrec
1069 rc=(ir-1)*iorj*4
1070 DO ib=1,4
1071 bc=(ib-1)*iorj+rc
1072 DO ij=1,iorj,ifactor
1073 ic=ic+1
1074 bdat(ic)=bwrk(ij)
1075 END DO
1076 IF ((ir.eq.1).and.(ib.eq.1)) ijdim=ic
1077 END DO
1078 END DO
1079 mpts=ic
1080!
1081! ---------------------
1082 CASE (v2dvar, v3dvar) ! V points field
1083! ---------------------
1084!
1085 CALL average_boundary2d (ng, model, tile, &
1086 & gtype, ncvname, extract_flag, &
1087 & imin, imax, jmin, jmax, nrec, &
1088 & npts, bdat, bwrk)
1089!
1090 bdat=0.0_r8
1091 ic=0
1092 DO ir=1,nrec
1093 rc=(ir-1)*iorj*4
1094 DO ib=1,4
1095 bc=(ib-1)*iorj+rc
1096 DO ij=1,iorj,ifactor
1097 ic=ic+1
1098 bdat(ic)=bwrk(ij)
1099 END DO
1100 IF ((ir.eq.1).and.(ib.eq.1)) ijdim=ic
1101 END DO
1102 END DO
1103 mpts=ic
1104!
1105! ------------
1106 CASE DEFAULT
1107! ------------
1108!
1109 IF (master) THEN
1110 WRITE (stdout,10) gtype, &
1111 & 'not supported for decimation:', &
1112 & trim(ncvname)
1113 END IF
1114 exit_flag=3
1115!
1116 END SELECT
1117 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1118!
1119! Set start and total vectors needed to write into output NetCDF file.
1120!
1121 start(1)=1
1122 total(1)=ijdim
1123 start(2)=1
1124 total(2)=4
1125 start(3)=1
1126 total(3)=nrec
1127 start(4)=tindex
1128 total(4)=1
1129!
1130 10 FORMAT (' DECIMATE_BOUNDARY2D - Staggered variable, gtype = ', &
1131 & i0,/,20x,a,2x,a)
1132!
1133 RETURN
1134 END SUBROUTINE decimate_boundary2d
1135!
1136 SUBROUTINE decimate_boundary3d (ng, model, tile, &
1137 & gtype, ncvname, tindex, &
1138 & Extract_Flag, &
1139 & Imin, Imax, Jmin, Jmax, LBk, UBk, &
1140 & Nrec, &
1141 & start, total, &
1142 & Npts, Bdat)
1143!
1144!=======================================================================
1145! !
1146! It decimates data from donor 3D boundary packed as 1D global array. !
1147! !
1148! On Input: !
1149! !
1150! ng Nested grid number (integer) !
1151! model Calling model identifier (integer) !
1152! tile Domain partition (integer) !
1153! gtype Staggered C-grid type (integer) !
1154! ncvname NetCDF variable name (string) !
1155! tindex Time record index to process (integer) !
1156! Extract_Flag Extraction/decimation flag (integer) !
1157! Imin Donor boundary starting data I-index (integer) !
1158! Imax Donor boundary ending data I-index (integer) !
1159! Jmin Donor boundary starting data J-index (integer) !
1160! Jmax Donor boundary ending data J-index (integer) !
1161! LBk Donor boundary K-dimension lower bound (integer) !
1162! UBk Donor boundary K-dimension upper bound (integer) !
1163! Nrec Number of boundary records (integer) !
1164! Npts Length of packed boundary data (integer) !
1165! Bdat Packed global donor 3D boundary data (real; 1D) !
1166! !
1167! On Output: !
1168! !
1169! start Start index where the first of the data values will !
1170! be written along each dimension (integer) !
1171! total Number of data values to be written along each !
1172! dimension (integer) !
1173! Bdat Decimated 3D boundary data (real; 1D) !
1174! !
1175!=======================================================================
1176!
1177! Imported variable declarations.
1178!
1179 integer, intent(in) :: ng, model, tile
1180 integer, intent(in) :: gtype, tindex, Extract_Flag
1181 integer, intent(in) :: Imin, Imax, Jmin, Jmax, LBk, UBk
1182 integer, intent(in) :: Nrec, Npts
1183 integer, intent(out) :: start(:), total(:)
1184!
1185 real(r8), intent(inout) :: Bdat(:)
1186!
1187 character (len=*), intent(in) :: ncvname
1188!
1189! Local variable declarations.
1190!
1191 integer :: bc, ib, ic, ij, ir, ifactor, k, rc
1192 integer :: Idim, IJdim, IJKlen, IorJ, Klen, Mpts
1193!
1194 real(r8) :: Bwrk(SIZE(Bdat))
1195!
1196 character (len=*), parameter :: MyFile = &
1197 & __FILE__//", decimate_boundary3d"
1198!
1199 sourcefile=myfile
1200!
1201!------------------------------------------------------------------------
1202! Decimate input packed 3D boundary by requested factor.
1203!------------------------------------------------------------------------
1204!
1205 ifactor=abs(extract_flag)
1206 iorj=iobounds(ng)%IorJ
1207 klen=ubk-lbk+1
1208 ijklen=iorj*klen
1209!
1210 SELECT CASE (gtype)
1211!
1212! ---------------------
1213 CASE (r2dvar, r3dvar) ! RHO points field
1214! ---------------------
1215!
1216 bwrk=bdat
1217 bdat=0.0_r8
1218 ic=0
1219 DO ir=1,nrec
1220 rc=(ir-1)*ijklen*4
1221 DO ib=1,4
1222 bc=(ib-1)*ijklen+rc
1223 DO k=lbk,ubk
1224 DO ij=1,iorj,ifactor
1225 ic=ic+1
1226 bwrk(ic)=bdat(ij)
1227 END DO
1228 IF ((ir.eq.1).and.(ib.eq.1).and.(k.eq.lbk)) ijdim=ic
1229 END DO
1230 END DO
1231 END DO
1232 mpts=ic
1233!
1234! ---------------------
1235 CASE (u2dvar, u3dvar) ! U points field
1236! ---------------------
1237!
1238 CALL average_boundary3d (ng, model, tile, &
1239 & gtype, ncvname, extract_flag, &
1240 & imin, imax, jmin, jmax, lbk, ubk, &
1241 & nrec, &
1242 & npts, bdat, bwrk)
1243!
1244 bdat=0.0_r8
1245 ic=0
1246 DO ir=1,nrec
1247 rc=(ir-1)*ijklen*4
1248 DO ib=1,4
1249 bc=(ib-1)*ijklen+rc
1250 DO k=lbk,ubk
1251 DO ij=1,iorj,ifactor
1252 ic=ic+1
1253 bwrk(ic)=bdat(ij)
1254 END DO
1255 IF ((ir.eq.1).and.(ib.eq.1).and.(k.eq.lbk)) ijdim=ic
1256 END DO
1257 END DO
1258 END DO
1259 mpts=ic
1260!
1261! ---------------------
1262 CASE (v2dvar, v3dvar) ! V points field
1263! ---------------------
1264!
1265 CALL average_boundary3d (ng, model, tile, &
1266 & gtype, ncvname, extract_flag, &
1267 & imin, imax, jmin, jmax, lbk, ubk, &
1268 & nrec, &
1269 & npts, bdat, bwrk)
1270!
1271 bdat=0.0_r8
1272 ic=0
1273 DO ir=1,nrec
1274 rc=(ir-1)*ijklen*4
1275 DO ib=1,4
1276 bc=(ib-1)*ijklen+rc
1277 DO k=lbk,ubk
1278 DO ij=1,iorj,ifactor
1279 ic=ic+1
1280 bwrk(ic)=bdat(ij)
1281 END DO
1282 IF ((ir.eq.1).and.(ib.eq.1).and.(k.eq.lbk)) ijdim=ic
1283 END DO
1284 END DO
1285 END DO
1286 mpts=ic
1287!
1288! ------------
1289 CASE DEFAULT
1290! ------------
1291!
1292 IF (master) THEN
1293 WRITE (stdout,10) gtype, &
1294 & 'not supported for decimation:', &
1295 & trim(ncvname)
1296 END IF
1297 exit_flag=3
1298!
1299 END SELECT
1300 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1301!
1302! Set start and total vectors needed to write into output NetCDF file.
1303!
1304 start(1)=1
1305 total(1)=ijdim
1306 start(2)=1
1307 total(2)=klen
1308 start(3)=1
1309 total(3)=4
1310 start(4)=1
1311 total(4)=nrec
1312 start(5)=tindex
1313 total(5)=1
1314!
1315 10 FORMAT (' DECIMATE_BOUNDARY3D - Staggered variable, gtype = ', &
1316 & i0,/,20x,a,2x,a)
1317!
1318 RETURN
1319 END SUBROUTINE decimate_boundary3d
1320# endif
1321!
1322 SUBROUTINE decimate_field2d (ng, model, tile, &
1323 & gtype, ifield, tindex, Extract_Flag, &
1324 & Imin, Imax, Jmin, Jmax, &
1325 & start, total, &
1326 & Npts, Fdat)
1327!
1328!=======================================================================
1329! !
1330! It decimates data from donor 2D fiels packed as 1D global array. !
1331! !
1332! On Input: !
1333! !
1334! ng Nested grid number (integer) !
1335! model Calling model identifier (integer) !
1336! tile Domain partition (integer) !
1337! gtype Staggered C-grid type (integer) !
1338! ifield Field metadata index (integer) !
1339! tindex Time record index to process (integer) !
1340! Extract_Flag Extraction/decimation flag (integer) !
1341! Imin Donor field starting data I-index (integer) !
1342! Imax Donor field ending data I-index (integer) !
1343! Jmin Donor field starting data J-index (integer) !
1344! Jmax Donor field ending data J-index (integer) !
1345! Fdat Packed global donor 2D field data (real 1D array) !
1346! !
1347! On Output: !
1348! !
1349! start Start index where the first of the data values will !
1350! be written along each dimension (integer) !
1351! total Number of data values to be written along each !
1352! dimension (integer) !
1353! Fdat Decimated 2D field data (real) !
1354! !
1355!=======================================================================
1356!
1357! Imported variable declarations.
1358!
1359 integer, intent(in) :: ng, model, tile
1360 integer, intent(in) :: gtype, ifield, tindex, Extract_Flag
1361 integer, intent(in) :: Imin, Imax, Jmin, Jmax
1362 integer, intent(in) :: Npts
1363 integer, intent(out) :: start(:), total(:)
1364!
1365 real(r8), intent(inout) :: Fdat(:)
1366!
1367! Local variable declarations.
1368!
1369 integer :: i, j, ij, ic, jc, ifactor
1370 integer :: Idim, Jdim, Ioff, Joff, Isize, Jsize, Mpts
1371!
1372 real(r8) :: Fwrk(SIZE(Fdat))
1373!
1374 character (len=*), parameter :: MyFile = &
1375 & __FILE__//", decimate_field2d"
1376!
1377 sourcefile=myfile
1378!
1379!------------------------------------------------------------------------
1380! Decimate input packed 2D field by requested factor.
1381!------------------------------------------------------------------------
1382!
1383 ifactor=abs(extract_flag)
1384 isize=imax-imin+1
1385 jsize=jmax-jmin+1
1386!
1387 SELECT CASE (gtype)
1388!
1389! ---------------------
1390 CASE (r2dvar, r3dvar) ! RHO points field
1391! ---------------------
1392!
1393 fwrk=fdat
1394 fdat=0.0_r8
1395 ioff=1
1396 joff=1
1397 ic=0
1398 jc=0
1399 DO j=jmin,jmax,ifactor ! decimate
1400 jc=jc+1
1401 DO i=imin,imax,ifactor
1402 ij=i+ioff+(j-1+joff)*isize
1403 ic=ic+1
1404 fdat(ic)=fwrk(ij)
1405 END DO
1406 IF (j.eq.jmin) idim=ic
1407 END DO
1408 jdim=jc
1409 mpts=ic
1410!
1411! ---------------------
1412 CASE (u2dvar, u3dvar) ! U points field
1413! ---------------------
1414!
1415 CALL average_field2d (ng, model, tile, &
1416 & gtype, ifield, extract_flag, &
1417 & imin, imax, jmin, jmax, &
1418 & npts, fdat, fwrk)
1419!
1420 fdat=0.0_r8
1421 ioff=0
1422 joff=1
1423 ic=0
1424 jc=0
1425 DO j=jmin,jmax,ifactor ! decimate
1426 jc=jc+1
1427 DO i=imin,imax,ifactor
1428 ij=i+ioff+(j-1+joff)*isize
1429 ic=ic+1
1430# ifdef MASKING
1431 fdat(ic)=fwrk(ij)*extract(ng)%Gmask_u(ic)
1432# else
1433 fdat(ic)=fwrk(ij)
1434# endif
1435 END DO
1436 IF (j.eq.jmin) idim=ic
1437 END DO
1438 jdim=jc
1439 mpts=ic
1440!
1441! ---------------------
1442 CASE (v2dvar, v3dvar) ! V points field
1443! ---------------------
1444!
1445 CALL average_field2d (ng, model, tile, &
1446 & gtype, ifield, extract_flag, &
1447 & imin, imax, jmin, jmax, &
1448 & npts, fdat, fwrk)
1449!
1450 fdat=0.0_r8
1451 ioff=1
1452 joff=0
1453 ic=0
1454 jc=0
1455 DO j=jmin,jmax,ifactor ! decimate
1456 jc=jc+1
1457 DO i=imin,imax,ifactor
1458 ij=i+ioff+(j-1+joff)*isize
1459 ic=ic+1
1460# ifdef MASKING
1461 fdat(ic)=fwrk(ij)*extract(ng)%Gmask_v(ic)
1462# else
1463 fdat(ic)=fwrk(ij)
1464# endif
1465 END DO
1466 IF (j.eq.jmin) idim=ic
1467 END DO
1468 jdim=jc
1469 mpts=ic
1470!
1471! ------------
1472 CASE DEFAULT
1473! ------------
1474!
1475 IF (master) THEN
1476 WRITE (stdout,10) gtype, &
1477 & 'not supported for decimation:', &
1478 & trim(vname(1,ifield))
1479 END IF
1480 exit_flag=3
1481!
1482 END SELECT
1483 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1484!
1485! Set start and total vectors needed to write into output NetCDF file.
1486!
1487 IF (gtype.gt.0) THEN
1488 start(1)=1
1489 total(1)=idim
1490 start(2)=1
1491 total(2)=jdim
1492 start(3)=tindex
1493 total(3)=1
1494# ifdef MASKING
1495 ELSE
1496 start(1)=1
1497 total(1)=mpts
1498 start(2)=tindex
1499 total(2)=1
1500# endif
1501 END IF
1502!
1503 10 FORMAT (' DECIMATE_FIELD2D - Staggered variable, gtype = ', i0, &
1504 & /,20x,a,2x,a)
1505!
1506 RETURN
1507 END SUBROUTINE decimate_field2d
1508!
1509 SUBROUTINE decimate_field3d (ng, model, tile, &
1510 & gtype, ifield, tindex, Extract_Flag, &
1511 & Imin, Imax, Jmin, Jmax, Kmin, Kmax, &
1512 & start, total, &
1513 & Npts, Fdat)
1514!
1515!=======================================================================
1516! !
1517! It decimates data from donor 2D fiels packed as 1D global array. !
1518! !
1519! On Input: !
1520! !
1521! ng Nested grid number (integer) !
1522! model Calling model identifier (integer) !
1523! tile Domain partition (integer) !
1524! gtype Staggered C-grid type (integer) !
1525! ifield Field metadata index (integer) !
1526! tindex Time record index to process (integer) !
1527! Extract_Flag Extraction/decimation flag (integer) !
1528! Imin Donor field starting data I-index (integer) !
1529! Imax Donor field ending data I-index (integer) !
1530! Jmin Donor field starting data J-index (integer) !
1531! Jmax Donor field ending data J-index (integer) !
1532! Kmin Donor field starting data K-index (integer) !
1533! Kmax Donor field ending data K-index (integer) !
1534! Fdat Packed global donor 2D field data (real 1D array) !
1535! !
1536! On Output: !
1537! !
1538! start Start index where the first of the data values will !
1539! be written along each dimension (integer) !
1540! total Number of data values to be written along each !
1541! dimension (integer) !
1542! Fdat Extracted 2D field data (real) !
1543! !
1544!=======================================================================
1545!
1546! Imported variable declarations.
1547!
1548 integer, intent(in) :: ng, model, tile
1549 integer, intent(in) :: gtype, ifield, tindex, Extract_Flag
1550 integer, intent(in) :: Imin, Imax, Jmin, Jmax, Kmin, Kmax
1551 integer, intent(in) :: Npts
1552 integer, intent(out) :: start(:), total(:)
1553!
1554 real(r8), intent(inout) :: Fdat(:)
1555!
1556! Local variable declarations.
1557!
1558 integer :: i, j, k, ij, ijk, ic, jc, ifactor, mc
1559 integer :: Idim, Jdim, Ioff, Joff, Koff
1560 integer :: Isize, IJsize, Jsize, Ksize, Mpts
1561!
1562 real(r8) :: Fwrk(SIZE(Fdat))
1563!
1564 character (len=*), parameter :: MyFile = &
1565 & __FILE__//", decimate_field3d"
1566!
1567 sourcefile=myfile
1568!
1569!------------------------------------------------------------------------
1570! Decimate input packed 3D field by requested factor.
1571!------------------------------------------------------------------------
1572!
1573 ifactor=abs(extract_flag)
1574 isize=imax-imin+1
1575 jsize=jmax-jmin+1
1576 ksize=kmax-kmin+1
1577 ijsize=isize*jsize
1578!
1579 SELECT CASE (gtype)
1580!
1581! -------------
1582 CASE (r3dvar) ! RHO points field
1583! -------------
1584!
1585 fwrk=fdat
1586 fdat=0.0_r8
1587 ioff=1
1588 joff=1
1589 koff=0
1590 ic=0
1591 jc=0
1592 DO k=kmin,kmax
1593 DO j=jmin,jmax,ifactor
1594 DO i=imin,imax,ifactor
1595 ij=i+ioff+(j-1+joff)*isize
1596 ijk=ij+(k-1+koff)*ijsize
1597 ic=ic+1
1598 fdat(ic)=fwrk(ijk)
1599 END DO
1600 IF ((j.eq.jmin).and.(k.eq.kmin)) idim=ic
1601 IF (k.eq.kmin) jc=jc+1
1602 END DO
1603 END DO
1604 jdim=jc
1605 mpts=ic
1606!
1607! -------------
1608 CASE (u3dvar) ! U points field
1609! -------------
1610!
1611 CALL average_field3d (ng, model, tile, &
1612 & gtype, ifield, extract_flag, &
1613 & imin, imax, jmin, jmax, kmin, kmax, &
1614 & npts, fdat, fwrk)
1615!
1616 fdat=0.0_r8
1617 ioff=0
1618 joff=1
1619 koff=0
1620 ic=0
1621 jc=0
1622 DO k=kmin,kmax
1623# ifdef MASKING
1624 mc=0
1625# endif
1626 DO j=jmin,jmax,ifactor
1627 DO i=imin,imax,ifactor
1628 ij=i+ioff+(j-1+joff)*isize
1629 ijk=ij+(k-1+koff)*ijsize
1630 ic=ic+1
1631# ifdef MASKING
1632 mc=mc+1
1633 fdat(ic)=fwrk(ijk)*extract(ng)%Gmask_u(mc)
1634# else
1635 fdat(ic)=fwrk(ijk)
1636# endif
1637 END DO
1638 IF ((j.eq.jmin).and.(k.eq.kmin)) idim=ic
1639 IF (k.eq.kmin) jc=jc+1
1640 END DO
1641 END DO
1642 jdim=jc
1643 mpts=ic
1644!
1645! -------------
1646 CASE (v3dvar) ! V points field
1647! -------------
1648!
1649 CALL average_field3d (ng, model, tile, &
1650 & gtype, ifield, extract_flag, &
1651 & imin, imax, jmin, jmax, kmin, kmax, &
1652 & npts, fdat, fwrk)
1653!
1654 fdat=0.0_r8
1655 ioff=1
1656 joff=0
1657 koff=0
1658 ic=0
1659 jc=0
1660 DO k=kmin,kmax
1661# ifdef MASKING
1662 mc=0
1663# endif
1664 DO j=jmin,jmax,ifactor
1665 DO i=imin,imax,ifactor
1666 ij=i+ioff+(j-1+joff)*isize
1667 ijk=ij+(k-1+koff)*ijsize
1668 ic=ic+1
1669# ifdef MASKING
1670 mc=mc+1
1671 fdat(ic)=fwrk(ijk)*extract(ng)%Gmask_v(mc)
1672# else
1673 fdat(ic)=fwrk(ijk)
1674# endif
1675 END DO
1676 IF ((j.eq.jmin).and.(k.eq.kmin)) idim=ic
1677 IF (k.eq.kmin) jc=jc+1
1678 END DO
1679 END DO
1680 jdim=jc
1681 mpts=ic
1682!
1683! -------------
1684 CASE (w3dvar) ! W points field
1685! -------------
1686!
1687 fwrk=fdat
1688 fdat=0.0_r8
1689 ioff=1
1690 joff=1
1691 koff=1
1692 ic=0
1693 jc=0
1694 DO k=kmin,kmax
1695 DO j=jmin,jmax,ifactor
1696 DO i=imin,imax,ifactor
1697 ij=i+ioff+(j-1+joff)*isize
1698 ijk=ij+(k-1+koff)*ijsize
1699 ic=ic+1
1700 fdat(ic)=fwrk(ijk)
1701 END DO
1702 IF ((j.eq.jmin).and.(k.eq.kmin)) idim=ic
1703 IF (k.eq.kmin) jc=jc+1
1704 END DO
1705 END DO
1706 jdim=jc
1707 mpts=ic
1708!
1709! -------------
1710 CASE DEFAULT
1711! -------------
1712!
1713 IF (master) THEN
1714 WRITE (stdout,10) gtype, &
1715 & 'not supported for decimation:', &
1716 & trim(vname(1,ifield))
1717 END IF
1718 exit_flag=3
1719!
1720 END SELECT
1721 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1722!
1723! Set start and total vectors needed to write into output NetCDF file.
1724!
1725 IF (gtype.gt.0) THEN
1726 start(1)=1
1727 total(1)=idim
1728 start(2)=1
1729 total(2)=jdim
1730 start(3)=1
1731 total(3)=ksize
1732 start(4)=tindex
1733 total(4)=1
1734# ifdef MASKING
1735 ELSE
1736 start(1)=1
1737 total(1)=mpts
1738 start(2)=tindex
1739 total(2)=1
1740# endif
1741 END IF
1742!
1743 10 FORMAT (' DECIMATE_FIELD3D - Staggered variable, gtype = ', i0, &
1744 & /,20x,a,2x,a)
1745!
1746 RETURN
1747 END SUBROUTINE decimate_field3d
1748!
1749 SUBROUTINE decimate_field4d (ng, model, tile, &
1750 & gtype, ifield, tindex, Extract_Flag, &
1751 & Imin, Imax, Jmin, Jmax, Kmin, Kmax, &
1752 & Fourth, Loff, &
1753 & start, total, &
1754 & Npts, Fdat)
1755!
1756!=======================================================================
1757! !
1758! It decimates data from donor 4D field packed as 1D global array. !
1759! The data is processed by 3D slices to reduce memory requirements !
1760! in the calling routine. !
1761! !
1762! On Input: !
1763! !
1764! ng Nested grid number (integer) !
1765! model Calling model identifier (integer) !
1766! tile Domain partition (integer) !
1767! gtype Staggered C-grid type (integer) !
1768! ifield Field metadata index (integer) !
1769! tindex Time record index to process (integer) !
1770! Extract_Flag Extraction/decimation flag (integer) !
1771! Imin Donor field starting data I-index (integer) !
1772! Imax Donor field ending data I-index (integer) !
1773! Jmin Donor field starting data J-index (integer) !
1774! Jmax Donor field ending data J-index (integer) !
1775! Kmin Donor field starting data K-index (integer) !
1776! Kmax Donor field ending data K-index (integer) !
1777! fourth Donor Fouth dimension index to process (integer) !
1778! Loff Fourth dimension couter offset (integer) !
1779! Fdat Packed global donor 3D field data (real 1D array) !
1780! !
1781! On Output: !
1782! !
1783! start Start index where the first of the data values will !
1784! be written along each dimension (integer) !
1785! total Number of data values to be written along each !
1786! dimension (integer) !
1787! Fdat Extracted 3D field data (real) !
1788! !
1789!=======================================================================
1790!
1791! Imported variable declarations.
1792!
1793 integer, intent(in) :: ng, model, tile
1794 integer, intent(in) :: gtype, ifield, tindex, Extract_Flag
1795 integer, intent(in) :: Imin, Imax, Jmin, Jmax, Kmin, Kmax
1796 integer, intent(in) :: fourth, Loff
1797 integer, intent(in) :: Npts
1798 integer, intent(out) :: start(:), total(:)
1799!
1800 real(r8), intent(inout) :: Fdat(:)
1801!
1802! Local variable declarations.
1803!
1804 integer :: i, j, k, ij, ijk, ic, jc, ifactor, mc
1805 integer :: Idim, Jdim, Ioff, Joff, Koff
1806 integer :: Isize, IJsize, Jsize, Ksize, Mpts
1807!
1808 real(r8) :: Fwrk(SIZE(Fdat))
1809!
1810 character (len=*), parameter :: MyFile = &
1811 & __FILE__//", decimate_field4d"
1812!
1813 sourcefile=myfile
1814!
1815!------------------------------------------------------------------------
1816! Decimate input packed 4D field by requested factor.
1817!------------------------------------------------------------------------
1818!
1819 ifactor=abs(extract_flag)
1820 isize=imax-imin+1
1821 jsize=jmax-jmin+1
1822 ksize=kmax-kmin+1
1823 ijsize=isize*jsize
1824!
1825 SELECT CASE (gtype)
1826!
1827! -------------
1828 CASE (r3dvar) ! RHO points field
1829! -------------
1830!
1831 fwrk=fdat
1832 fdat=0.0_r8
1833 ioff=1
1834 joff=1
1835 koff=0
1836 ic=0
1837 jc=0
1838 DO k=kmin,kmax
1839 DO j=jmin,jmax,ifactor
1840 DO i=imin,imax,ifactor
1841 ij=i+ioff+(j-1+joff)*isize
1842 ijk=ij+(k-1+koff)*ijsize
1843 ic=ic+1
1844 fdat(ic)=fwrk(ijk)
1845 END DO
1846 IF ((j.eq.jmin).and.(k.eq.kmin)) idim=ic
1847 IF (k.eq.kmin) jc=jc+1
1848 END DO
1849 END DO
1850 jdim=jc
1851 mpts=ic
1852!
1853! -------------
1854 CASE (u3dvar) ! U points field
1855! -------------
1856!
1857 CALL average_field3d (ng, model, tile, &
1858 & gtype, ifield, extract_flag, &
1859 & imin, imax, jmin, jmax, kmin, kmax, &
1860 & npts, fdat, fwrk)
1861!
1862 fdat=0.0_r8
1863 ioff=0
1864 joff=1
1865 koff=0
1866 ic=0
1867 jc=0
1868 DO k=kmin,kmax
1869# ifdef MASKING
1870 mc=0
1871# endif
1872 DO j=jmin,jmax,ifactor
1873 DO i=imin,imax,ifactor
1874 ij=i+ioff+(j-1+joff)*isize
1875 ijk=ij+(k-1+koff)*ijsize
1876 ic=ic+1
1877# ifdef MASKING
1878 mc=mc+1
1879 fdat(ic)=fwrk(ijk)*extract(ng)%Gmask_u(mc)
1880# else
1881 fdat(ic)=fwrk(ijk)
1882# endif
1883 END DO
1884 IF ((j.eq.jmin).and.(k.eq.kmin)) idim=ic
1885 IF (k.eq.kmin) jc=jc+1
1886 END DO
1887 END DO
1888 jdim=jc
1889 mpts=ic
1890!
1891! -------------
1892 CASE (v3dvar) ! V points field
1893! -------------
1894!
1895 CALL average_field3d (ng, model, tile, &
1896 & gtype, ifield, extract_flag, &
1897 & imin, imax, jmin, jmax, kmin, kmax, &
1898 & npts, fdat, fwrk)
1899!
1900 fdat=0.0_r8
1901 ioff=1
1902 joff=0
1903 koff=0
1904 ic=0
1905 jc=0
1906 DO k=kmin,kmax
1907# ifdef MASKING
1908 mc=0
1909# endif
1910 DO j=jmin,jmax,ifactor
1911 DO i=imin,imax,ifactor
1912 ij=i+ioff+(j-1+joff)*isize
1913 ijk=ij+(k-1+koff)*ijsize
1914 ic=ic+1
1915# ifdef MASKING
1916 mc=mc+1
1917 fdat(ic)=fwrk(ijk)*extract(ng)%Gmask_v(mc)
1918# else
1919 fdat(ic)=fwrk(ijk)
1920# endif
1921 END DO
1922 IF ((j.eq.jmin).and.(k.eq.kmin)) idim=ic
1923 IF (k.eq.kmin) jc=jc+1
1924 END DO
1925 END DO
1926 jdim=jc
1927 mpts=ic
1928!
1929! -------------
1930 CASE (w3dvar) ! W points field
1931! -------------
1932!
1933 fwrk=fdat
1934 fdat=0.0_r8
1935 ioff=1
1936 joff=1
1937 koff=1
1938 ic=0
1939 jc=0
1940 DO k=kmin,kmax
1941 DO j=jmin,jmax,ifactor
1942 DO i=imin,imax,ifactor
1943 ij=i+ioff+(j-1+joff)*isize
1944 ijk=ij+(k-1+koff)*ijsize
1945 ic=ic+1
1946 fdat(ic)=fwrk(ijk)
1947 END DO
1948 IF ((j.eq.jmin).and.(k.eq.kmin)) idim=ic
1949 IF (k.eq.kmin) jc=jc+1
1950 END DO
1951 END DO
1952 jdim=jc
1953 mpts=ic
1954!
1955! -------------
1956 CASE DEFAULT
1957! -------------
1958!
1959 IF (master) WRITE (stdout,10) trim(vname(1,ifield)), gtype
1960 10 FORMAT (' DECIMATE_FIELD4D - Staggered variable, gtype = ', &
1961 & i0,/,20x,'not supported for field: ',a)
1962 exit_flag=3
1963!
1964 END SELECT
1965 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1966!
1967! Set start and total vectors needed to write into output NetCDF file.
1968!
1969 IF (gtype.gt.0) THEN
1970 start(1)=1
1971 total(1)=idim
1972 start(2)=1
1973 total(2)=jdim
1974 start(3)=1
1975 total(3)=ksize
1976 start(4)=fourth+loff
1977 total(4)=1
1978 start(5)=tindex
1979 total(5)=1
1980# ifdef MASKING
1981 ELSE
1982 start(1)=1+(fourth+loff-1)*mpts
1983 total(1)=mpts
1984 start(2)=tindex
1985 total(2)=1
1986# endif
1987 END IF
1988!
1989 RETURN
1990 END SUBROUTINE decimate_field4d
1991!
1992 SUBROUTINE interp_coords (ng, tile, model, gtype, &
1993 & Ainp, &
1994# ifdef MASKING
1995 & Minp, &
1996# endif
1997 & Xinp, Yinp, &
1998# ifdef MASKING
1999 & Mout, &
2000# endif
2001 & Xout, Yout, &
2002 & Iout, Jout)
2003!
2004!=======================================================================
2005! !
2006! Suppose fields are extracted with the interpolation method. In that !
2007! case, this routine computes the fractional I- and J-coordinates !
2008! (Iout, Jout) of the extracted data in terms of the donor grid !
2009! application to facilitate the interpolation. Ipos and Jpos are the !
2010! donor interpolant left-bottom cell corner coordinates containing !
2011! each extracted field grid point. !
2012! !
2013! The 3D fields are interpolated level by level, so only horizontal !
2014! fractional coordinates are needed for each staggered location. !
2015! !
2016! Since the interpolation is done during the writing of the output !
2017! extraction fields, the strategy here is to set the geometry arrays !
2018! as 1D global arrays to facilitate the interpolation between the !
2019! donor grid-tiled data in parallel applications. The interpolation !
2020! stencil is trivial in fractional global coordinates for bilinear !
2021! or bicubic methodologies. !
2022! !
2023! On Input: !
2024! !
2025! ng Nested grid number (integer) !
2026! model Calling model identifier (integer) !
2027! tile Domain partition (integer) !
2028! gtype Staggered C-grid type (integer) !
2029! !
2030! On Output: !
2031! !
2032! Ainp Donor grid global curvilinear angle (1D real) !
2033# ifdef MASKING
2034! Minp Donor grid global land/sea mask (1D real) !
2035# endif
2036! Xinp Donor grid global X-coordinate (1D real) !
2037! Yinp Donor grid global Y-coordinate (1D real) !
2038# ifdef MASKING
2039! Mout Interpolated data global land/sea mask (1D real) !
2040# endif
2041! Xout Interpolated data global X-coordinates (1D real) !
2042! Yout Interpolated data global Y-coordinates (1D real) !
2043! Iout Global fractional I-coordinate (1D real) !
2044! Jout Global fractional J-coordinate (1D real) !
2045! !
2046!=======================================================================
2047!
2048! Imported variable declarations.
2049!
2050 integer, intent(in) :: ng, tile, model, gtype
2051!
2052 real(r8), intent(out) :: Ainp(:)
2053# ifdef MASKING
2054 real(r8), intent(out) :: Minp(:)
2055# endif
2056 real(r8), intent(out) :: Xinp(:)
2057 real(r8), intent(out) :: Yinp(:)
2058# ifdef MASKING
2059 real(r8), intent(out) :: Mout(:)
2060# endif
2061 real(r8), intent(out) :: Xout(:)
2062 real(r8), intent(out) :: Yout(:)
2063 real(r8), intent(out) :: Iout(:)
2064 real(r8), intent(out) :: Jout(:)
2065!
2066! Local variable declarations.
2067!
2068 integer :: LBi_inp, UBi_inp, LBj_inp, UBj_inp
2069 integer :: LBi_out, UBi_out, LBj_out, UBj_out
2070 integer :: Is_inp, Ie_inp, Js_inp, Je_inp
2071 integer :: Is_out, Ie_out, Js_out, Je_out
2072 integer :: Isize, Jsize, Msize, Nsize
2073 integer :: Cgrid, Npts, ghost, i, ic, j
2074!
2075 real(dp) :: scale
2076!
2077 real(r8), allocatable :: angle(:,:)
2078# ifdef MASKING
2079 real(r8), pointer :: mask_inp(:,:), mask_out(:,:)
2080# endif
2081 real(r8), pointer :: Xi(:,:), Yi(:,:)
2082 real(r8), pointer :: Xo(:,:), Yo(:,:)
2083!
2084 character (len=*), parameter :: MyFile = &
2085 & __FILE__//", interp_coord"
2086
2087# include "set_bounds.h"
2088!
2089!-----------------------------------------------------------------------
2090! Compute fractional coordinates (Ipos, Jpos).
2091!-----------------------------------------------------------------------
2092!
2093 scale=1.0_dp
2094!
2095! Get input donor and output extract grids arrays tiled bounds.
2096!
2097 lbi_inp=bounds(ng)%LBi(tile)
2098 ubi_inp=bounds(ng)%UBi(tile)
2099 lbj_inp=bounds(ng)%LBj(tile)
2100 ubj_inp=bounds(ng)%UBj(tile)
2101!
2102 lbi_out=xtr_bounds(ng)%LBi(tile)
2103 ubi_out=xtr_bounds(ng)%UBi(tile)
2104 lbj_out=xtr_bounds(ng)%LBj(tile)
2105 ubj_out=xtr_bounds(ng)%UBj(tile)
2106!
2107! Get input donor and output extract grids I/O bounds and coordinates
2108! according to staggered locations.
2109!
2110 SELECT CASE (gtype)
2111!
2112! -----------------------------
2113 CASE (p2dvar, p3dvar) ! PSI-points
2114! -----------------------------
2115!
2116 cgrid=1
2117!
2118 is_inp=iobounds(ng)%ILB_psi
2119 ie_inp=iobounds(ng)%IUB_psi
2120 js_inp=iobounds(ng)%JLB_psi
2121 je_inp=iobounds(ng)%JUB_psi
2122!
2123 is_out=xtr_iobounds(ng)%ILB_psi
2124 ie_out=xtr_iobounds(ng)%IUB_psi
2125 js_out=xtr_iobounds(ng)%JLB_psi
2126 je_out=xtr_iobounds(ng)%JUB_psi
2127!
2128 IF (spherical) THEN
2129 xi => grid(ng)%lonp
2130 yi => grid(ng)%latp
2131 xo => extract(ng)%lonp
2132 yo => extract(ng)%latp
2133 ELSE
2134 xi => grid(ng)%xp
2135 yi => grid(ng)%yp
2136 xo => extract(ng)%xp
2137 yo => extract(ng)%yp
2138 END IF
2139# ifdef MASKING
2140 mask_inp => grid(ng)%pmask
2141 mask_out => extract(ng)%pmask
2142# endif
2143!
2144! -----------------------------
2145 CASE (r2dvar, r3dvar, w3dvar) ! RHO-points
2146! -----------------------------
2147!
2148 cgrid=2
2149!
2150 is_inp=iobounds(ng)%ILB_rho
2151 ie_inp=iobounds(ng)%IUB_rho
2152 js_inp=iobounds(ng)%JLB_rho
2153 je_inp=iobounds(ng)%JUB_rho
2154!
2155 is_out=xtr_iobounds(ng)%ILB_rho
2156 ie_out=xtr_iobounds(ng)%IUB_rho
2157 js_out=xtr_iobounds(ng)%JLB_rho
2158 je_out=xtr_iobounds(ng)%JUB_rho
2159!
2160 IF (spherical) THEN
2161 xi => grid(ng)%lonr
2162 yi => grid(ng)%latr
2163 xo => extract(ng)%lonr
2164 yo => extract(ng)%latr
2165 ELSE
2166 xi => grid(ng)%xr
2167 yi => grid(ng)%yr
2168 xo => extract(ng)%xr
2169 yo => extract(ng)%yr
2170 END IF
2171# ifdef MASKING
2172 mask_inp => grid(ng)%rmask
2173 mask_out => extract(ng)%rmask
2174# endif
2175!
2176! -----------------------------
2177 CASE (u2dvar, u3dvar) ! U-points
2178! -----------------------------
2179!
2180 cgrid=3
2181!
2182 is_inp=iobounds(ng)%ILB_u
2183 ie_inp=iobounds(ng)%IUB_u
2184 js_inp=iobounds(ng)%JLB_u
2185 je_inp=iobounds(ng)%JUB_u
2186!
2187 is_out=xtr_iobounds(ng)%ILB_u
2188 ie_out=xtr_iobounds(ng)%IUB_u
2189 js_out=xtr_iobounds(ng)%JLB_u
2190 je_out=xtr_iobounds(ng)%JUB_u
2191!
2192 IF (spherical) THEN
2193 xi => grid(ng)%lonu
2194 yi => grid(ng)%latu
2195 xo => extract(ng)%lonu
2196 yo => extract(ng)%latu
2197 ELSE
2198 xi => grid(ng)%xu
2199 yi => grid(ng)%yu
2200 xo => extract(ng)%xu
2201 yo => extract(ng)%yu
2202 END IF
2203# ifdef MASKING
2204 mask_inp => grid(ng)%umask
2205 mask_out => extract(ng)%umask
2206
2207# endif
2208!
2209! -----------------------------
2210 CASE (v2dvar, v3dvar) ! V-points
2211! -----------------------------
2212!
2213 cgrid=4
2214!
2215 is_inp=iobounds(ng)%ILB_v
2216 ie_inp=iobounds(ng)%IUB_v
2217 js_inp=iobounds(ng)%JLB_v
2218 je_inp=iobounds(ng)%JUB_v
2219!
2220 is_out=xtr_iobounds(ng)%ILB_v
2221 ie_out=xtr_iobounds(ng)%IUB_v
2222 js_out=xtr_iobounds(ng)%JLB_v
2223 je_out=xtr_iobounds(ng)%JUB_v
2224!
2225 IF (spherical) THEN
2226 xi => grid(ng)%lonv
2227 yi => grid(ng)%latv
2228 xo => extract(ng)%lonv
2229 yo => extract(ng)%latv
2230 ELSE
2231 xi => grid(ng)%xv
2232 yi => grid(ng)%yv
2233 xo => extract(ng)%xv
2234 yo => extract(ng)%yv
2235 END IF
2236# ifdef MASKING
2237 mask_inp => grid(ng)%vmask
2238 mask_out => extract(ng)%vmask
2239# endif
2240!
2241 END SELECT
2242!
2243 isize=ie_inp-is_inp+1
2244 jsize=je_inp-js_inp+1
2245 msize=ie_out-is_out+1
2246 nsize=je_out-js_out+1
2247!
2248! Set donor grid curvilinear rotation angle.
2249!
2250 IF (.not.allocated(angle)) THEN
2251 allocate ( angle(lbi_inp:ubi_inp,lbj_inp:ubj_inp) )
2252 END IF
2253!
2254 IF (cgrid.eq.1) THEN ! PSI-points
2255 DO j=jstrp,jendt
2256 DO i=istrp,iendt
2257 angle(i,j)=0.25_r8*(grid(ng)%angler(i-1,j-1)+ &
2258 & grid(ng)%angler(i-1,j )+ &
2259 & grid(ng)%angler(i ,j-1)+ &
2260 & grid(ng)%angler(i ,j ))
2261 END DO
2262 END DO
2263 ELSE IF (cgrid.eq.2) THEN ! RHO-points
2264 DO j=jstrt,jendt
2265 DO i=istrt,iendt
2266 angle(i,j)=grid(ng)%angler(i,j)
2267 END DO
2268 END DO
2269 ELSE IF (cgrid.eq.3) THEN ! U-points
2270 DO j=jstrt,jendt
2271 DO i=istrp,iendt
2272 angle(i,j)=0.5_r8*(grid(ng)%angler(i-1,j)+ &
2273 & grid(ng)%angler(i ,j))
2274 END DO
2275 END DO
2276 ELSE IF (cgrid.eq.4) THEN ! V-points
2277 DO j=jstrp,jendt
2278 DO i=istrt,iendt
2279 angle(i,j)=0.5_r8*(grid(ng)%angler(i,j-1)+ &
2280 & grid(ng)%angler(i,j ))
2281 END DO
2282 END DO
2283 END IF
2284
2285# ifdef DISTRIBUTE
2286!
2287! If distributed-memory, collect donor grid location data from all
2288! spawned nodes and store it into a global scratch 1D array, packed
2289! in column-major order to facilitate interpolation.
2290
2291# ifdef MASKING
2292!
2293 CALL mp_gather2d (ng, model, &
2294 & lbi_inp, ubi_inp, lbj_inp, ubj_inp, &
2295 & 0, gtype, scale, &
2296 & mask_inp, &
2297 & mask_inp, npts, minp, .false.)
2298# endif
2299!
2300 CALL mp_gather2d (ng, model, &
2301 & lbi_inp, ubi_inp, lbj_inp, ubj_inp, &
2302 & 0, gtype, scale, &
2303# ifdef MASKING
2304 & mask_inp, &
2305# endif
2306 & angle, npts, ainp, .false.)
2307!
2308 CALL mp_gather2d (ng, model, &
2309 & lbi_inp, ubi_inp, lbj_inp, ubj_inp, &
2310 & 0, gtype, scale, &
2311# ifdef MASKING
2312 & mask_inp, &
2313# endif
2314 & xi, npts, xinp, .false.)
2315!
2316 CALL mp_gather2d (ng, model, &
2317 & lbi_inp, ubi_inp, lbj_inp, ubj_inp, &
2318 & 0, gtype, scale, &
2319# ifdef MASKING
2320 & mask_inp, &
2321# endif
2322 & yi, npts, yinp, .false.)
2323
2324!
2325 CALL mp_gather2d (ng, model, &
2326 & lbi_inp, ubi_inp, lbj_inp, ubj_inp, &
2327 & 0, gtype, scale, &
2328# ifdef MASKING
2329 & mask_inp, &
2330# endif
2331 & xi, npts, xinp, .false.)
2332!
2333 CALL mp_gather2d (ng, model, &
2334 & lbi_inp, ubi_inp, lbj_inp, ubj_inp, &
2335 & 0, gtype, scale, &
2336# ifdef MASKING
2337 & mask_inp, &
2338# endif
2339 & yi, npts, yinp, .false.)
2340
2341# ifdef MASKING
2342!
2343 CALL mp_gather2d_xtr (ng, model, &
2344 & lbi_out, ubi_out, lbj_out, ubj_out, &
2345 & 0, gtype, scale, &
2346 & mask_out, &
2347 & mask_out, npts, mout, .false.)
2348# endif
2349!
2350 CALL mp_gather2d_xtr (ng, model, &
2351 & lbi_out, ubi_out, lbj_out, ubj_out, &
2352 & 0, gtype, scale, &
2353# ifdef MASKING
2354 & mask_out, &
2355# endif
2356 & xo, npts, xout, .false.)
2357!
2358 CALL mp_gather2d_xtr (ng, model, &
2359 & lbi_out, ubi_out, lbj_out, ubj_out, &
2360 & 0, gtype, scale, &
2361# ifdef MASKING
2362 & mask_out, &
2363# endif
2364 & yo, npts, yout, .false.)
2365
2366# else
2367
2368!
2369! If serial or shared-memory applications, pack donor grid data into a
2370! global 1D array in column-major order.
2371!
2372 ainp=pack(angle, .true.)
2373# ifdef MASKING
2374 minp=pack(mask_inp, .true.)
2375# endif
2376 xinp=pack(xi, .true.)
2377 yinp=pack(yi, .true.)
2378!
2379# ifdef MASKING
2380 mout=pack(mask_out, .true.)
2381# endif
2382 xout=pack(xo, .true.)
2383 yout=pack(yo, .true.)
2384# endif
2385
2386# ifdef DISTRIBUTE
2387!
2388! Find the donor grid cell's fractional indices (Iout, Jout) containing
2389! the extracted field points. Thus, computing the interpolation weights
2390! is trivial. Here, each mpi-process calculates the global 1D-packed
2391! Iout and Jout indices once to avoid collecting communications when
2392! writing the extracted data.
2393!
2394 CALL hindices (ng, &
2395 & 1, isize, 1, jsize, &
2396 & 1, isize, 1, jsize, &
2397 & ainp, xinp, yinp, &
2398 & 1, msize, 1, nsize, &
2399 & 1, msize, 1, nsize, &
2400 & xout, yout, &
2401 & iout, jout, &
2402 & 0.0_r8, .false.)
2403# else
2404!
2405! Find the fractional indices (Iout,Jout) of the donor grid cells
2406! containing extracted field points.
2407!
2408 CALL hindices (ng, &
2409 & 1, isize, 1, jsize, &
2410 & 1, isize, 1, jsize, &
2411 & ainp, xinp, yinp, &
2412 & lbi_out, ubi_out, lbj_out, ubj_out, &
2413 & is_out, ie_out, js_out, je_out, &
2414 & xout, yout, &
2415 & iout, jout, &
2416 & 0.0_r8, .false.)
2417# endif
2418!
2419! Deallocate.
2420!
2421 IF (allocated(angle)) deallocate ( angle )
2422!
2423 RETURN
2424 END SUBROUTINE interp_coords
2425!
2426 SUBROUTINE interp_field2d (ng, model, tile, &
2427 & gtype, ifield, tindex, &
2428 & Imin, Imax, Jmin, Jmax, &
2429 & start, total, &
2430 & Npts, Fdat)
2431!
2432!=======================================================================
2433! !
2434! It interpolates data from donor 2D fields, packed as a 1D global, !
2435! array, to specified extract grid geometry. !
2436! !
2437! On Input: !
2438! !
2439! ng Nested grid number (integer) !
2440! model Calling model identifier (integer) !
2441! tile Domain partition (integer) !
2442! gtype Staggered C-grid type (integer) !
2443! ifield Field metadata index (integer) !
2444! tindex Time record index to process (integer) !
2445! Imin Donor field starting data I-index (integer) !
2446! Imax Donor field ending data I-index (integer) !
2447! Jmin Donor field starting data J-index (integer) !
2448! Jmax Donor field ending data J-index (integer) !
2449! Fdat Packed global donor 2D field data (real 1D array) !
2450! !
2451! On Output: !
2452! !
2453! start Start index where the first of the data values will !
2454! be written along each dimension (integer) !
2455! total Number of data values to be written along each !
2456! dimension (integer) !
2457! Fdat Extracted 2D field data (real) !
2458! !
2459!=======================================================================
2460!
2461! Imported variable declarations.
2462!
2463 integer, intent(in) :: ng, model, tile
2464 integer, intent(in) :: gtype, ifield, tindex
2465 integer, intent(in) :: Imin, Imax, Jmin, Jmax
2466 integer, intent(in) :: Npts
2467 integer, intent(out) :: start(:), total(:)
2468!
2469 real(r8), intent(inout) :: Fdat(:)
2470!
2471! Local variable declarations.
2472!
2473 logical :: LandFill
2474!
2475 integer :: ghost, i, j, ij, ic, jc, ifactor, method
2476 integer :: Cgrid, Ilen, Jlen, Ioff, Joff, Isize, Jsize
2477 integer :: Istr, Iend, Jstr, Jend
2478 integer :: LBi, UBi, LBj, UBj
2479!
2480 real(r8) :: Fmin, Fmax
2481
2482 real(r8) :: Fwrk(Npts)
2483!
2484 character (len=*), parameter :: MyFile = &
2485 & __FILE__//", extract_field2d"
2486!
2487 sourcefile=myfile
2488!
2489!------------------------------------------------------------------------
2490! Interpolate
2491!------------------------------------------------------------------------
2492!
2493 landfill=.true.
2494 ghost=0
2495 method=bilinear
2496!
2497 isize=imax-imin+1
2498 jsize=jmax-jmin+1
2499!
2500 lbi=xtr_bounds(ng)%LBi(tile)
2501 ubi=xtr_bounds(ng)%UBi(tile)
2502 lbj=xtr_bounds(ng)%LBj(tile)
2503 ubj=xtr_bounds(ng)%UBj(tile)
2504!
2505! At this point, the donor data is packed as global 1D arrays to
2506! facilitate interpolating to the extracted grid tiled data.
2507!
2508 SELECT CASE (gtype)
2509!
2510! ---------------------
2511 CASE (r2dvar, r3dvar) ! RHO points field
2512! ---------------------
2513!
2514 cgrid=2
2515 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
2516 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
2517 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
2518 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
2519!
2520 fwrk=fdat(1:npts)
2521!
2522 IF (spherical) THEN
2523 CALL regrid_field2d (ng, model, tile, gtype, ifield, &
2524 & method, landfill, &
2525 & 1, isize, 1, jsize, &
2526 & grid(ng) % Gx_rho, &
2527 & grid(ng) % Gy_rho, &
2528# ifdef MASKING
2529 & grid(ng) % Gmask_rho, &
2530# endif
2531 & fwrk, &
2532 & istr, iend, jstr, jend, &
2533 & lbi, ubi, lbj, ubj, &
2534 & extract(ng) % Iout_rho, &
2535 & extract(ng) % Jout_rho, &
2536 & extract(ng) % lonr, &
2537 & extract(ng) % latr, &
2538# ifdef MASKING
2539 & extract(ng) % rmask, &
2540# endif
2541 & npts, fdat, fmin, fmax)
2542 ELSE
2543 CALL regrid_field2d (ng, model, tile, gtype, ifield, &
2544 & method, landfill, &
2545 & 1, isize, 1, jsize, &
2546 & grid(ng) % Gx_rho, &
2547 & grid(ng) % Gy_rho, &
2548# ifdef MASKING
2549 & grid(ng) % Gmask_rho, &
2550# endif
2551 & fwrk, &
2552 & istr, iend, jstr, jend, &
2553 & lbi, ubi, lbj, ubj, &
2554 & extract(ng) % Iout_rho, &
2555 & extract(ng) % Jout_rho, &
2556 & extract(ng) % xr, &
2557 & extract(ng) % Yr, &
2558# ifdef MASKING
2559 & extract(ng) % rmask, &
2560# endif
2561 & npts, fdat, fmin, fmax)
2562 END IF
2563!
2564! ---------------------
2565 CASE (u2dvar, u3dvar) ! U points field
2566! ---------------------
2567!
2568 cgrid=3
2569 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
2570 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
2571 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
2572 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
2573!
2574 fwrk=fdat(1:npts)
2575!
2576 IF (spherical) THEN
2577 CALL regrid_field2d (ng, model, tile, gtype, ifield, &
2578 & method, landfill, &
2579 & 1, isize, 1, jsize, &
2580 & grid(ng) % Gx_u, &
2581 & grid(ng) % Gy_u, &
2582# ifdef MASKING
2583 & grid(ng) % Gmask_u, &
2584# endif
2585 & fwrk, &
2586 & istr, iend, jstr, jend, &
2587 & lbi, ubi, lbj, ubj, &
2588 & extract(ng) % Iout_u, &
2589 & extract(ng) % Jout_u, &
2590 & extract(ng) % lonu, &
2591 & extract(ng) % latu, &
2592# ifdef MASKING
2593 & extract(ng) % umask, &
2594# endif
2595 & npts, fdat, fmin, fmax)
2596 ELSE
2597 CALL regrid_field2d (ng, model, tile, gtype, ifield, &
2598 & method, landfill, &
2599 & 1, isize, 1, jsize, &
2600 & grid(ng) % Gx_u, &
2601 & grid(ng) % Gy_u, &
2602# ifdef MASKING
2603 & grid(ng) % Gmask_u, &
2604# endif
2605 & fwrk, &
2606 & istr, iend, jstr, jend, &
2607 & lbi, ubi, lbj, ubj, &
2608 & extract(ng) % Iout_u, &
2609 & extract(ng) % Jout_u, &
2610 & extract(ng) % xu, &
2611 & extract(ng) % yu, &
2612# ifdef MASKING
2613 & extract(ng) % umask, &
2614# endif
2615 & npts, fdat, fmin, fmax)
2616 END IF
2617!
2618! ---------------------
2619 CASE (v2dvar, v3dvar) ! V points field
2620! ---------------------
2621!
2622 cgrid=4
2623 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
2624 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
2625 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
2626 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
2627!
2628 fwrk=fdat(1:npts)
2629!
2630 IF (spherical) THEN
2631 CALL regrid_field2d (ng, model, tile, gtype, ifield, &
2632 & method, landfill, &
2633 & 1, isize, 1, jsize, &
2634 & grid(ng) % Gx_v, &
2635 & grid(ng) % Gy_v, &
2636# ifdef MASKING
2637 & grid(ng) % Gmask_v, &
2638# endif
2639 & fwrk, &
2640 & istr, iend, jstr, jend, &
2641 & lbi, ubi, lbj, ubj, &
2642 & extract(ng) % Iout_v, &
2643 & extract(ng) % Jout_v, &
2644 & extract(ng) % lonv, &
2645 & extract(ng) % latv, &
2646# ifdef MASKING
2647 & extract(ng) % vmask, &
2648# endif
2649 & npts, fdat, fmin, fmax)
2650 ELSE
2651 CALL regrid_field2d (ng, model, tile, gtype, ifield, &
2652 & method, landfill, &
2653 & 1, isize, 1, jsize, &
2654 & grid(ng) % Gx_v, &
2655 & grid(ng) % Gy_v, &
2656# ifdef MASKING
2657 & grid(ng) % Gmask_v, &
2658# endif
2659 & fwrk, &
2660 & istr, iend, jstr, jend, &
2661 & lbi, ubi, lbj, ubj, &
2662 & extract(ng) % Iout_v, &
2663 & extract(ng) % Jout_v, &
2664 & extract(ng) % lonv, &
2665 & extract(ng) % latv, &
2666# ifdef MASKING
2667 & extract(ng) % vmask, &
2668# endif
2669 & npts, fdat, fmin, fmax)
2670 END IF
2671!
2672! ------------
2673 CASE DEFAULT
2674! ------------
2675!
2676 IF (master) THEN
2677 WRITE (stdout,10) gtype, &
2678 & 'not supported for interpolation:', &
2679 & trim(vname(1,ifield))
2680 END IF
2681 exit_flag=3
2682!
2683 END SELECT
2684 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2685!
2686! Set start and total vectors needed to write into output NetCDF file.
2687!
2688 IF (gtype.gt.0) THEN
2689 start(1)=1
2690 total(1)=iend-istr+1
2691 start(2)=1
2692 total(2)=jend-jstr+1
2693 start(3)=tindex
2694 total(3)=1
2695# ifdef MASKING
2696 ELSE
2697 start(1)=1
2698 total(1)=npts
2699 start(2)=tindex
2700 total(2)=1
2701# endif
2702 END IF
2703!
2704 10 FORMAT (' INTERP_FIELD2D - Staggered variable, gtype = ', i0, &
2705 & /,18x,a,2x,a)
2706!
2707 RETURN
2708 END SUBROUTINE interp_field2d
2709!
2710 SUBROUTINE interp_field3d (ng, model, tile, &
2711 & gtype, ifield, tindex, &
2712 & Imin, Imax, Jmin, Jmax, Kmin, Kmax, &
2713 & start, total, &
2714 & Npts, Fdat)
2715!
2716!=======================================================================
2717! !
2718! It interpolates data from donor 3D fields, packed as a 1D global, !
2719! array, to specified extract grid geometry. !
2720! !
2721! On Input: !
2722! !
2723! ng Nested grid number (integer) !
2724! model Calling model identifier (integer) !
2725! tile Domain partition (integer) !
2726! gtype Staggered C-grid type (integer) !
2727! ifield Field metadata index (integer) !
2728! tindex Time record index to process (integer) !
2729! Extract_Flag Extraction/decimation flag (integer) !
2730! Imin Donor field starting data I-index (integer) !
2731! Imax Donor field ending data I-index (integer) !
2732! Jmin Donor field starting data J-index (integer) !
2733! Jmax Donor field ending data J-index (integer) !
2734! Kmin Donor field starting data K-index (integer) !
2735! Kmax Donor field ending data K-index (integer) !
2736! Fdat Packed global donor 2D field data (real 1D array) !
2737! !
2738! On Output: !
2739! !
2740! start Start index where the first of the data values will !
2741! be written along each dimension (integer) !
2742! total Number of data values to be written along each !
2743! dimension (integer) !
2744! Npts Number of points processed in Fwrk. !
2745! Fdat Extracted 2D field data (real) !
2746! !
2747!=======================================================================
2748!
2749! Imported variable declarations.
2750!
2751 integer, intent(in) :: ng, model, tile
2752 integer, intent(in) :: gtype, ifield, tindex
2753 integer, intent(in) :: Imin, Imax, Jmin, Jmax, Kmin, Kmax
2754 integer, intent(in) :: Npts
2755 integer, intent(out) :: start(:), total(:)
2756!
2757 real(r8), intent(inout) :: Fdat(:)
2758!
2759! Local variable declarations.
2760!
2761 logical :: LandFill
2762!
2763 integer :: ghost, i, j, k, ij, ijk, ic, jc, ifactor, method
2764 integer :: Ilen, Jlen, Ioff, Joff, Koff
2765 integer :: Cgrid, Isize, IJsize, Jsize, Ksize
2766 integer :: Istr, Iend, Jstr, Jend
2767 integer :: LBi, UBi, LBj, UBj
2768!
2769 real(r8) :: Fmin, Fmax
2770
2771 real(r8) :: Fwrk(Npts)
2772!
2773 character (len=*), parameter :: MyFile = &
2774 & __FILE__//", interp_field3d"
2775!
2776 sourcefile=myfile
2777!
2778!------------------------------------------------------------------------
2779! Interpolate
2780!------------------------------------------------------------------------
2781!
2782 landfill=.true.
2783 ghost=0
2784 method=bilinear
2785!
2786 isize=imax-imin+1
2787 jsize=jmax-jmin+1
2788 ksize=kmax-kmin+1
2789!
2790 lbi=xtr_bounds(ng)%LBi(tile)
2791 ubi=xtr_bounds(ng)%UBi(tile)
2792 lbj=xtr_bounds(ng)%LBj(tile)
2793 ubj=xtr_bounds(ng)%UBj(tile)
2794!
2795 SELECT CASE (gtype)
2796!
2797! ---------------------
2798 CASE (r3dvar, w3dvar) ! RHO and W points field
2799! ---------------------
2800!
2801 cgrid=2
2802 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
2803 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
2804 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
2805 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
2806!
2807 fwrk=fdat(1:npts)
2808!
2809 IF (spherical) THEN
2810 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
2811 & method, landfill, &
2812 & 1, isize, 1, jsize, kmin, kmax, &
2813 & grid(ng) % Gx_rho, &
2814 & grid(ng) % Gy_rho, &
2815# ifdef MASKING
2816 & grid(ng) % Gmask_rho, &
2817# endif
2818 & fwrk, &
2819 & istr, iend, jstr, jend, &
2820 & lbi, ubi, lbj, ubj, &
2821 & extract(ng) % Iout_rho, &
2822 & extract(ng) % Jout_rho, &
2823 & extract(ng) % lonr, &
2824 & extract(ng) % latr, &
2825# ifdef MASKING
2826 & extract(ng) % rmask, &
2827# endif
2828 & npts, fdat, fmin, fmax)
2829 ELSE
2830 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
2831 & method, landfill, &
2832 & 1, isize, 1, jsize, kmin, kmax, &
2833 & grid(ng) % Gx_rho, &
2834 & grid(ng) % Gy_rho, &
2835# ifdef MASKING
2836 & grid(ng) % Gmask_rho, &
2837# endif
2838 & fwrk, &
2839 & istr, iend, jstr, jend, &
2840 & lbi, ubi, lbj, ubj, &
2841 & extract(ng) % Iout_rho, &
2842 & extract(ng) % Jout_rho, &
2843 & extract(ng) % xr, &
2844 & extract(ng) % Yr, &
2845# ifdef MASKING
2846 & extract(ng) % rmask, &
2847# endif
2848 & npts, fdat, fmin, fmax)
2849 END IF
2850!
2851! ---------------------
2852 CASE (u3dvar) ! U points field
2853! ---------------------
2854!
2855 cgrid=3
2856 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
2857 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
2858 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
2859 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
2860!
2861 fwrk=fdat(1:npts)
2862!
2863 IF (spherical) THEN
2864 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
2865 & method, landfill, &
2866 & 1, isize, 1, jsize, kmin, kmax, &
2867 & grid(ng) % Gx_u, &
2868 & grid(ng) % Gy_u, &
2869# ifdef MASKING
2870 & grid(ng) % Gmask_u, &
2871# endif
2872 & fwrk, &
2873 & istr, iend, jstr, jend, &
2874 & lbi, ubi, lbj, ubj, &
2875 & extract(ng) % Iout_u, &
2876 & extract(ng) % Jout_u, &
2877 & extract(ng) % lonu, &
2878 & extract(ng) % latu, &
2879# ifdef MASKING
2880 & extract(ng) % umask, &
2881# endif
2882 & npts, fdat, fmin, fmax)
2883 ELSE
2884 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
2885 & method, landfill, &
2886 & 1, isize, 1, jsize, kmin, kmax, &
2887 & grid(ng) % Gx_u, &
2888 & grid(ng) % Gy_u, &
2889# ifdef MASKING
2890 & grid(ng) % Gmask_u, &
2891# endif
2892 & fwrk, &
2893 & istr, iend, jstr, jend, &
2894 & lbi, ubi, lbj, ubj, &
2895 & extract(ng) % Iout_u, &
2896 & extract(ng) % Jout_u, &
2897 & extract(ng) % xu, &
2898 & extract(ng) % yu, &
2899# ifdef MASKING
2900 & extract(ng) % umask, &
2901# endif
2902 & npts, fdat, fmin, fmax)
2903 END IF
2904!
2905! ---------------------
2906 CASE (v3dvar) ! V points field
2907! ---------------------
2908!
2909 cgrid=4
2910 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
2911 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
2912 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
2913 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
2914!
2915 fwrk=fdat(1:npts)
2916!
2917 IF (spherical) THEN
2918 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
2919 & method, landfill, &
2920 & 1, isize, 1, jsize, kmin, kmax, &
2921 & grid(ng) % Gx_v, &
2922 & grid(ng) % Gy_v, &
2923# ifdef MASKING
2924 & grid(ng) % Gmask_v, &
2925# endif
2926 & fwrk, &
2927 & istr, iend, jstr, jend, &
2928 & lbi, ubi, lbj, ubj, &
2929 & extract(ng) % Iout_v, &
2930 & extract(ng) % Jout_v, &
2931 & extract(ng) % lonv, &
2932 & extract(ng) % latv, &
2933# ifdef MASKING
2934 & extract(ng) % vmask, &
2935# endif
2936 & npts, fdat, fmin, fmax)
2937 ELSE
2938 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
2939 & method, landfill, &
2940 & 1, isize, 1, jsize, kmin, kmax, &
2941 & grid(ng) % Gx_v, &
2942 & grid(ng) % Gy_v, &
2943# ifdef MASKING
2944 & grid(ng) % Gmask_v, &
2945# endif
2946 & fwrk, &
2947 & istr, iend, jstr, jend, &
2948 & lbi, ubi, lbj, ubj, &
2949 & extract(ng) % Iout_v, &
2950 & extract(ng) % Jout_v, &
2951 & extract(ng) % lonv, &
2952 & extract(ng) % latv, &
2953# ifdef MASKING
2954 & extract(ng) % vmask, &
2955# endif
2956 & npts, fdat, fmin, fmax)
2957 END IF
2958!
2959! ------------
2960 CASE DEFAULT
2961! ------------
2962!
2963 IF (master) THEN
2964 WRITE (stdout,10) gtype, &
2965 & 'not supported for interpolation:', &
2966 & trim(vname(1,ifield))
2967 END IF
2968 exit_flag=3
2969!
2970 END SELECT
2971 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2972!
2973! Set start and total vectors needed to write into output NetCDF file.
2974!
2975 IF (gtype.gt.0) THEN
2976 start(1)=1
2977 total(1)=iend-istr+1
2978 start(2)=1
2979 total(2)=jend-jstr+1
2980 start(3)=1
2981 total(3)=ksize
2982 start(4)=tindex
2983 total(4)=1
2984# ifdef MASKING
2985 ELSE
2986 start(1)=1
2987 total(1)=npts
2988 start(2)=tindex
2989 total(2)=1
2990# endif
2991 END IF
2992!
2993 10 FORMAT (' INTERP_FIELD3D - Staggered variable, gtype = ', i0, &
2994 & /,18x,a,2x,a)
2995!
2996 RETURN
2997 END SUBROUTINE interp_field3d
2998!
2999 SUBROUTINE interp_field4d (ng, model, tile, &
3000 & gtype, ifield, tindex, &
3001 & Imin, Imax, Jmin, Jmax, Kmin, Kmax, &
3002 & Fourth, Loff, &
3003 & start, total, &
3004 & Npts, Fdat)
3005!
3006!=======================================================================
3007! !
3008! It interpolates data from donor 2D fields, packed as a 1D global, !
3009! array, to specified extract grid geometry. The data is processed !
3010! by 3D slices to reduce memory requirements in the calling routine. !
3011! !
3012! On Input: !
3013! !
3014! ng Nested grid number (integer) !
3015! model Calling model identifier (integer) !
3016! tile Domain partition (integer) !
3017! gtype Staggered C-grid type (integer) !
3018! ifield Field metadata index (integer) !
3019! tindex Time record index to process (integer) !
3020! Imin Donor field starting data I-index (integer) !
3021! Imax Donor field ending data I-index (integer) !
3022! Jmin Donor field starting data J-index (integer) !
3023! Jmax Donor field ending data J-index (integer) !
3024! Kmin Donor field starting data K-index (integer) !
3025! Kmax Donor field ending data K-index (integer) !
3026! fourth Donor Fouth dimension index to process (integer) !
3027! Loff Fourth dimension couter offset (integer) !
3028! Fdat Packed global donor 3D field data (real 1D array) !
3029! !
3030! On Output: !
3031! !
3032! start Start index where the first of the data values will !
3033! be written along each dimension (integer) !
3034! total Number of data values to be written along each !
3035! dimension (integer) !
3036! Fdat Extracted 3D field data (real) !
3037! !
3038!=======================================================================
3039!
3040! Imported variable declarations.
3041!
3042 integer, intent(in) :: ng, model, tile
3043 integer, intent(in) :: gtype, ifield, tindex
3044 integer, intent(in) :: Imin, Imax, Jmin, Jmax, Kmin, Kmax
3045 integer, intent(in) :: fourth, Loff
3046 integer, intent(in) :: Npts
3047 integer, intent(out) :: start(:), total(:)
3048!
3049 real(r8), intent(inout) :: Fdat(:)
3050!
3051! Local variable declarations.
3052!
3053 logical :: LandFill
3054!
3055 integer :: ghost, i, j, k, ij, ijk, ic, jc, ifactor, method
3056 integer :: Cgrid, Ilen, Jlen, Ioff, Joff, Koff
3057 integer :: Isize, IJsize, Jsize, Ksize
3058 integer :: Istr, Iend, Jstr, Jend
3059 integer :: LBi, UBi, LBj, UBj
3060!
3061 real(r8) :: Fmin, Fmax
3062
3063 real(r8) :: Fwrk(Npts)
3064!
3065 character (len=*), parameter :: MyFile = &
3066 & __FILE__//", interp_field4d"
3067!
3068 sourcefile=myfile
3069!
3070!------------------------------------------------------------------------
3071! Interpolate.
3072!------------------------------------------------------------------------
3073!
3074 landfill=.true.
3075 ghost=0
3076 method=bilinear
3077!
3078 isize=imax-imin+1
3079 jsize=jmax-jmin+1
3080!
3081 lbi=xtr_bounds(ng)%LBi(tile)
3082 ubi=xtr_bounds(ng)%UBi(tile)
3083 lbj=xtr_bounds(ng)%LBj(tile)
3084 ubj=xtr_bounds(ng)%UBj(tile)
3085!
3086 SELECT CASE (gtype)
3087!
3088! ---------------------
3089 CASE (r3dvar, w3dvar) ! RHO and W points field
3090! ---------------------
3091!
3092 cgrid=2
3093 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
3094 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
3095 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
3096 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
3097!
3098 fwrk=fdat(1:npts)
3099!
3100 IF (spherical) THEN
3101 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
3102 & method, landfill, &
3103 & 1, isize, 1, jsize, kmin, kmax, &
3104 & grid(ng) % Gx_rho, &
3105 & grid(ng) % Gy_rho, &
3106# ifdef MASKING
3107 & grid(ng) % Gmask_rho, &
3108# endif
3109 & fwrk, &
3110 & istr, iend, jstr, jend, &
3111 & lbi, ubi, lbj, ubj, &
3112 & extract(ng) % Iout_rho, &
3113 & extract(ng) % Jout_rho, &
3114 & extract(ng) % lonr, &
3115 & extract(ng) % latr, &
3116# ifdef MASKING
3117 & extract(ng) % rmask, &
3118# endif
3119 & npts, fdat, fmin, fmax)
3120 ELSE
3121 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
3122 & method, landfill, &
3123 & 1, isize, 1, jsize, kmin, kmax, &
3124 & grid(ng) % Gx_rho, &
3125 & grid(ng) % Gy_rho, &
3126# ifdef MASKING
3127 & grid(ng) % Gmask_rho, &
3128# endif
3129 & fwrk, &
3130 & istr, iend, jstr, jend, &
3131 & lbi, ubi, lbj, ubj, &
3132 & extract(ng) % Iout_rho, &
3133 & extract(ng) % Jout_rho, &
3134 & extract(ng) % xr, &
3135 & extract(ng) % Yr, &
3136# ifdef MASKING
3137 & extract(ng) % rmask, &
3138# endif
3139 & npts, fdat, fmin, fmax)
3140 END IF
3141!
3142! ---------------------
3143 CASE (u3dvar) ! U points field
3144! ---------------------
3145!
3146 cgrid=3
3147 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
3148 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
3149 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
3150 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
3151!
3152 fwrk=fdat(1:npts)
3153!
3154 IF (spherical) THEN
3155 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
3156 & method, landfill, &
3157 & 1, isize, 1, jsize, kmin, kmax, &
3158 & grid(ng) % Gx_u, &
3159 & grid(ng) % Gy_u, &
3160# ifdef MASKING
3161 & grid(ng) % Gmask_u, &
3162# endif
3163 & fwrk, &
3164 & istr, iend, jstr, jend, &
3165 & lbi, ubi, lbj, ubj, &
3166 & extract(ng) % Iout_u, &
3167 & extract(ng) % Jout_u, &
3168 & extract(ng) % lonu, &
3169 & extract(ng) % latu, &
3170# ifdef MASKING
3171 & extract(ng) % umask, &
3172# endif
3173 & npts, fdat, fmin, fmax)
3174 ELSE
3175 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
3176 & method, landfill, &
3177 & 1, isize, 1, jsize, kmin, kmax, &
3178 & grid(ng) % Gx_u, &
3179 & grid(ng) % Gy_u, &
3180# ifdef MASKING
3181 & grid(ng) % Gmask_u, &
3182# endif
3183 & fwrk, &
3184 & istr, iend, jstr, jend, &
3185 & lbi, ubi, lbj, ubj, &
3186 & extract(ng) % Iout_u, &
3187 & extract(ng) % Jout_u, &
3188 & extract(ng) % xu, &
3189 & extract(ng) % yu, &
3190# ifdef MASKING
3191 & extract(ng) % umask, &
3192# endif
3193 & npts, fdat, fmin, fmax)
3194 END IF
3195!
3196! ---------------------
3197 CASE (v3dvar) ! V points field
3198! ---------------------
3199!
3200 cgrid=4
3201 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
3202 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
3203 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
3204 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
3205!
3206 fwrk=fdat(1:npts)
3207!
3208 IF (spherical) THEN
3209 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
3210 & method, landfill, &
3211 & 1, isize, 1, jsize, kmin, kmax, &
3212 & grid(ng) % Gx_v, &
3213 & grid(ng) % Gy_v, &
3214# ifdef MASKING
3215 & grid(ng) % Gmask_v, &
3216# endif
3217 & fwrk, &
3218 & istr, iend, jstr, jend, &
3219 & lbi, ubi, lbj, ubj, &
3220 & extract(ng) % Iout_v, &
3221 & extract(ng) % Jout_v, &
3222 & extract(ng) % lonv, &
3223 & extract(ng) % latv, &
3224# ifdef MASKING
3225 & extract(ng) % vmask, &
3226# endif
3227 & npts, fdat, fmin, fmax)
3228 ELSE
3229 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
3230 & method, landfill, &
3231 & 1, isize, 1, jsize, kmin, kmax, &
3232 & grid(ng) % Gx_v, &
3233 & grid(ng) % Gy_v, &
3234# ifdef MASKING
3235 & grid(ng) % Gmask_v, &
3236# endif
3237 & fwrk, &
3238 & istr, iend, jstr, jend, &
3239 & lbi, ubi, lbj, ubj, &
3240 & extract(ng) % Iout_v, &
3241 & extract(ng) % Jout_v, &
3242 & extract(ng) % lonv, &
3243 & extract(ng) % latv, &
3244# ifdef MASKING
3245 & extract(ng) % vmask, &
3246# endif
3247 & npts, fdat, fmin, fmax)
3248 END IF
3249!
3250! ------------
3251 CASE DEFAULT
3252! ------------
3253!
3254 IF (master) THEN
3255 WRITE (stdout,10) gtype, &
3256 & 'not supported for interpolation:', &
3257 & trim(vname(1,ifield))
3258 END IF
3259 exit_flag=3
3260!
3261 END SELECT
3262 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
3263!
3264! Set start and total vectors needed to write into output NetCDF file.
3265!
3266 IF (gtype.gt.0) THEN
3267 start(1)=1
3268 total(1)=iend-istr+1
3269 start(2)=1
3270 total(2)=jend-jstr+1
3271 start(3)=1
3272 total(3)=ksize
3273 start(4)=fourth+loff
3274 total(4)=1
3275 start(5)=tindex
3276 total(5)=1
3277# ifdef MASKING
3278 ELSE
3279 start(1)=1+(fourth+loff-1)*npts
3280 total(1)=npts
3281 start(2)=tindex
3282 total(2)=1
3283# endif
3284 END IF
3285!
3286 10 FORMAT (' INTERP_FIELD4D - Staggered variable, gtype = ', i0, &
3287 & /,18x,a,2x,a)
3288!
3289 RETURN
3290 END SUBROUTINE interp_field4d
3291!
3292 SUBROUTINE interp_field2d_global (ng, model, tile, &
3293 & gtype, ifield, tindex, &
3294 & Imin, Imax, Jmin, Jmax, &
3295 & start, total, &
3296 & Npts, Fdat)
3297!
3298!=======================================================================
3299! !
3300! It interpolates data from donor 2D fields, packed as a 1D global, !
3301! array, to specified extract grid geometry. !
3302! !
3303! On Input: !
3304! !
3305! ng Nested grid number (integer) !
3306! model Calling model identifier (integer) !
3307! tile Domain partition (integer) !
3308! gtype Staggered C-grid type (integer) !
3309! ifield Field metadata index (integer) !
3310! tindex Time record index to process (integer) !
3311! Imin Donor field starting data I-index (integer) !
3312! Imax Donor field ending data I-index (integer) !
3313! Jmin Donor field starting data J-index (integer) !
3314! Jmax Donor field ending data J-index (integer) !
3315! Fdat Packed global donor 2D field data (real 1D array) !
3316! !
3317! On Output: !
3318! !
3319! start Start index where the first of the data values will !
3320! be written along each dimension (integer) !
3321! total Number of data values to be written along each !
3322! dimension (integer) !
3323! Fdat Extracted 2D field data (real) !
3324! !
3325!=======================================================================
3326!
3327! Imported variable declarations.
3328!
3329 integer, intent(in) :: ng, model, tile
3330 integer, intent(in) :: gtype, ifield, tindex
3331 integer, intent(in) :: Imin, Imax, Jmin, Jmax
3332 integer, intent(in) :: Npts
3333 integer, intent(out) :: start(:), total(:)
3334!
3335 real(r8), intent(inout) :: Fdat(:)
3336!
3337! Local variable declarations.
3338!
3339 logical :: LandFill
3340!
3341 integer :: ghost, i, j, ij, ic, jc, ifactor, method
3342 integer :: Cgrid, Ilen, Jlen, Ioff, Joff, Isize, Jsize
3343 integer :: Istr, Iend, Jstr, Jend
3344 integer :: LBi, UBi, LBj, UBj
3345!
3346 real(r8) :: Fmin, Fmax
3347
3348 real(r8) :: Fwrk(Npts)
3349!
3350 character (len=*), parameter :: MyFile = &
3351 & __FILE__//", extract_field2d_global"
3352!
3353 sourcefile=myfile
3354!
3355!------------------------------------------------------------------------
3356! Interpolate
3357!------------------------------------------------------------------------
3358!
3359 landfill=.true.
3360 ghost=0
3361 method=bilinear
3362!
3363 isize=imax-imin+1
3364 jsize=jmax-jmin+1
3365!
3366 lbi=xtr_bounds(ng)%LBi(tile)
3367 ubi=xtr_bounds(ng)%UBi(tile)
3368 lbj=xtr_bounds(ng)%LBj(tile)
3369 ubj=xtr_bounds(ng)%UBj(tile)
3370!
3371! At this point, the donor data is packed as global 1D arrays to
3372! facilitate interpolating to the extracted grid tiled data.
3373!
3374 SELECT CASE (gtype)
3375!
3376! ---------------------
3377 CASE (r2dvar, r3dvar) ! RHO points field
3378! ---------------------
3379!
3380 cgrid=2
3381 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
3382 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
3383 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
3384 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
3385!
3386 fwrk=fdat(1:npts)
3387!
3388 CALL regrid_field2d (ng, model, tile, gtype, ifield, &
3389 & method, landfill, &
3390 & imin, imax, jmin, jmax, &
3391 & grid(ng) % Gx_rho, &
3392 & grid(ng) % Gy_rho, &
3393# ifdef MASKING
3394 & grid(ng) % Gmask_rho, &
3395# endif
3396 & fwrk, &
3397 & istr, iend, jstr, jend, &
3398 & istr, iend, jstr, jend, &
3399 & extract(ng) % Iout_rho, &
3400 & extract(ng) % Jout_rho, &
3401 & extract(ng) % Gx_rho, &
3402 & extract(ng) % Gy_rho, &
3403# ifdef MASKING
3404 & extract(ng) % Gmask_rho, &
3405# endif
3406 & npts, fdat, fmin, fmax)
3407!
3408! ---------------------
3409 CASE (u2dvar, u3dvar) ! U points field
3410! ---------------------
3411!
3412 cgrid=3
3413 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
3414 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
3415 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
3416 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
3417!
3418 fwrk=fdat(1:npts)
3419!
3420 CALL regrid_field2d (ng, model, tile, gtype, ifield, &
3421 & method, landfill, &
3422 & imin, imax, jmin, jmax, &
3423 & grid(ng) % Gx_u, &
3424 & grid(ng) % Gy_u, &
3425# ifdef MASKING
3426 & grid(ng) % Gmask_u, &
3427# endif
3428 & fwrk, &
3429 & istr, iend, jstr, jend, &
3430 & istr, iend, jstr, jend, &
3431 & extract(ng) % Iout_u, &
3432 & extract(ng) % Jout_u, &
3433 & extract(ng) % Gx_u, &
3434 & extract(ng) % Gx_u, &
3435# ifdef MASKING
3436 & extract(ng) % Gmask_u, &
3437# endif
3438 & npts, fdat, fmin, fmax)
3439!
3440! ---------------------
3441 CASE (v2dvar, v3dvar) ! V points field
3442! ---------------------
3443!
3444 cgrid=4
3445 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
3446 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
3447 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
3448 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
3449!
3450 fwrk=fdat(1:npts)
3451!
3452 CALL regrid_field2d (ng, model, tile, gtype, ifield, &
3453 & method, landfill, &
3454 & imin, imax, jmin, jmax, &
3455 & grid(ng) % Gx_v, &
3456 & grid(ng) % Gy_v, &
3457# ifdef MASKING
3458 & grid(ng) % Gmask_v, &
3459# endif
3460 & fwrk, &
3461 & istr, iend, jstr, jend, &
3462 & istr, iend, jstr, jend, &
3463 & extract(ng) % Iout_v, &
3464 & extract(ng) % Jout_v, &
3465 & extract(ng) % Gx_v, &
3466 & extract(ng) % Gy_v, &
3467# ifdef MASKING
3468 & extract(ng) % Gmask_v, &
3469# endif
3470 & npts, fdat, fmin, fmax)
3471!
3472! ------------
3473 CASE DEFAULT
3474! ------------
3475!
3476 IF (master) THEN
3477 WRITE (stdout,10) gtype, &
3478 & 'not supported for interpolation:', &
3479 & trim(vname(1,ifield))
3480 END IF
3481 exit_flag=3
3482!
3483 END SELECT
3484 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
3485!
3486! Set start and total vectors needed to write into output NetCDF file.
3487!
3488 IF (gtype.gt.0) THEN
3489 start(1)=1
3490 total(1)=iend-istr+1
3491 start(2)=1
3492 total(2)=jend-jstr+1
3493 start(3)=tindex
3494 total(3)=1
3495# ifdef MASKING
3496 ELSE
3497 start(1)=1
3498 total(1)=npts
3499 start(2)=tindex
3500 total(2)=1
3501# endif
3502 END IF
3503!
3504 10 FORMAT (' INTERP_FIELD2D_GLOBAL - Staggered variable, gtype = ', &
3505 & i0,/,18x,a,2x,a)
3506!
3507 RETURN
3508 END SUBROUTINE interp_field2d_global
3509!
3510 SUBROUTINE interp_field3d_global (ng, model, tile, &
3511 & gtype, ifield, tindex, &
3512 & Imin, Imax, Jmin, Jmax, &
3513 & Kmin, Kmax, &
3514 & start, total, &
3515 & Npts, Fdat)
3516!
3517!=======================================================================
3518! !
3519! It interpolates data from donor 3D fields, packed as a 1D global, !
3520! array, to specified extract grid geometry. !
3521! !
3522! On Input: !
3523! !
3524! ng Nested grid number (integer) !
3525! model Calling model identifier (integer) !
3526! tile Domain partition (integer) !
3527! gtype Staggered C-grid type (integer) !
3528! ifield Field metadata index (integer) !
3529! tindex Time record index to process (integer) !
3530! Imin Donor field starting data I-index (integer) !
3531! Imax Donor field ending data I-index (integer) !
3532! Jmin Donor field starting data J-index (integer) !
3533! Jmax Donor field ending data J-index (integer) !
3534! Kmin Donor field starting data K-index (integer) !
3535! Kmax Donor field ending data K-index (integer) !
3536! Fdat Packed global donor 2D field data (real 1D array) !
3537! !
3538! On Output: !
3539! !
3540! start Start index where the first of the data values will !
3541! be written along each dimension (integer) !
3542! total Number of data values to be written along each !
3543! dimension (integer) !
3544! Npts Number of points processed in Fwrk. !
3545! Fdat Extracted 2D field data (real) !
3546! !
3547!=======================================================================
3548!
3549! Imported variable declarations.
3550!
3551 integer, intent(in) :: ng, model, tile
3552 integer, intent(in) :: gtype, ifield, tindex
3553 integer, intent(in) :: Imin, Imax, Jmin, Jmax, Kmin, Kmax
3554 integer, intent(in) :: Npts
3555 integer, intent(out) :: start(:), total(:)
3556!
3557 real(r8), intent(inout) :: Fdat(:)
3558!
3559! Local variable declarations.
3560!
3561 logical :: LandFill
3562!
3563 integer :: ghost, i, j, k, ij, ijk, ic, jc, ifactor, method
3564 integer :: Ilen, Jlen, Ioff, Joff, Koff
3565 integer :: Cgrid, Isize, IJsize, Jsize, Ksize
3566 integer :: Istr, Iend, Jstr, Jend
3567 integer :: LBi, UBi, LBj, UBj
3568!
3569 real(r8) :: Fmin, Fmax
3570
3571 real(r8) :: Fwrk(Npts)
3572!
3573 character (len=*), parameter :: MyFile = &
3574 & __FILE__//", interp_field3d_global"
3575!
3576 sourcefile=myfile
3577!
3578!------------------------------------------------------------------------
3579! Interpolate
3580!------------------------------------------------------------------------
3581!
3582 landfill=.true.
3583 ghost=0
3584 method=bilinear
3585!
3586 isize=imax-imin+1
3587 jsize=jmax-jmin+1
3588 ksize=kmax-kmin+1
3589!
3590 lbi=xtr_bounds(ng)%LBi(tile)
3591 ubi=xtr_bounds(ng)%UBi(tile)
3592 lbj=xtr_bounds(ng)%LBj(tile)
3593 ubj=xtr_bounds(ng)%UBj(tile)
3594!
3595 SELECT CASE (gtype)
3596!
3597! ---------------------
3598 CASE (r3dvar, w3dvar) ! RHO and W points field
3599! ---------------------
3600!
3601 cgrid=2
3602 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
3603 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
3604 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
3605 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
3606!
3607 fwrk=fdat(1:npts)
3608!
3609 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
3610 & method, landfill, &
3611 & imin, imax, jmin, jmax, kmin, kmax, &
3612 & grid(ng) % Gx_rho, &
3613 & grid(ng) % Gy_rho, &
3614# ifdef MASKING
3615 & grid(ng) % Gmask_rho, &
3616# endif
3617 & fwrk, &
3618 & istr, iend, jstr, jend, &
3619 & istr, iend, jstr, jend, &
3620 & extract(ng) % Iout_rho, &
3621 & extract(ng) % Jout_rho, &
3622 & extract(ng) % Gx_rho, &
3623 & extract(ng) % Gy_rho, &
3624# ifdef MASKING
3625 & extract(ng) % Gmask_rho, &
3626# endif
3627 & npts, fdat, fmin, fmax)
3628!
3629! ---------------------
3630 CASE (u3dvar) ! U points field
3631! ---------------------
3632!
3633 cgrid=3
3634 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
3635 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
3636 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
3637 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
3638!
3639 fwrk=fdat(1:npts)
3640!
3641 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
3642 & method, landfill, &
3643 & imin, imax, jmin, jmax, kmin, kmax, &
3644 & grid(ng) % Gx_u, &
3645 & grid(ng) % Gy_u, &
3646# ifdef MASKING
3647 & grid(ng) % Gmask_u, &
3648# endif
3649 & fwrk, &
3650 & istr, iend, jstr, jend, &
3651 & istr, iend, jstr, jend, &
3652 & extract(ng) % Iout_u, &
3653 & extract(ng) % Jout_u, &
3654 & extract(ng) % Gx_u, &
3655 & extract(ng) % Gy_u, &
3656# ifdef MASKING
3657 & extract(ng) % Gmask_u, &
3658# endif
3659 & npts, fdat, fmin, fmax)
3660!
3661! ---------------------
3662 CASE (v3dvar) ! V points field
3663! ---------------------
3664!
3665 cgrid=4
3666 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
3667 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
3668 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
3669 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
3670!
3671 fwrk=fdat(1:npts)
3672!
3673 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
3674 & method, landfill, &
3675 & imin, imax, jmin, jmax, kmin, kmax, &
3676 & grid(ng) % Gx_v, &
3677 & grid(ng) % Gy_v, &
3678# ifdef MASKING
3679 & grid(ng) % Gmask_v, &
3680# endif
3681 & fwrk, &
3682 & istr, iend, jstr, jend, &
3683 & istr, iend, jstr, jend, &
3684 & extract(ng) % Iout_v, &
3685 & extract(ng) % Jout_v, &
3686 & extract(ng) % Gx_v, &
3687 & extract(ng) % Gy_v, &
3688# ifdef MASKING
3689 & extract(ng) % Gmask_v, &
3690# endif
3691 & npts, fdat, fmin, fmax)
3692!
3693! ------------
3694 CASE DEFAULT
3695! ------------
3696!
3697 IF (master) THEN
3698 WRITE (stdout,10) gtype, &
3699 & 'not supported for interpolation:', &
3700 & trim(vname(1,ifield))
3701 END IF
3702 exit_flag=3
3703!
3704 END SELECT
3705 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
3706!
3707! Set start and total vectors needed to write into output NetCDF file.
3708!
3709 IF (gtype.gt.0) THEN
3710 start(1)=1
3711 total(1)=iend-istr+1
3712 start(2)=1
3713 total(2)=jend-jstr+1
3714 start(3)=1
3715 total(3)=ksize
3716 start(4)=tindex
3717 total(4)=1
3718# ifdef MASKING
3719 ELSE
3720 start(1)=1
3721 total(1)=npts
3722 start(2)=tindex
3723 total(2)=1
3724# endif
3725 END IF
3726!
3727 10 FORMAT (' INTERP_FIELD3D_GLOBAL - Staggered variable, gtype = ', &
3728 & i0,/,18x,a,2x,a)
3729!
3730 RETURN
3731 END SUBROUTINE interp_field3d_global
3732!
3733 SUBROUTINE interp_field4d_global (ng, model, tile, &
3734 & gtype, ifield, tindex, &
3735 & Imin, Imax, Jmin, Jmax, &
3736 & Kmin, Kmax, Fourth, Loff, &
3737 & start, total, &
3738 & Npts, Fdat)
3739!
3740!=======================================================================
3741! !
3742! It interpolates data from donor 2D fields, packed as a 1D global, !
3743! array, to specified extract grid geometry. The data is processed !
3744! by 3D slices to reduce memory requirements in the calling routine. !
3745! !
3746! On Input: !
3747! !
3748! ng Nested grid number (integer) !
3749! model Calling model identifier (integer) !
3750! tile Domain partition (integer) !
3751! gtype Staggered C-grid type (integer) !
3752! ifield Field metadata index (integer) !
3753! tindex Time record index to process (integer) !
3754! Imin Donor field starting data I-index (integer) !
3755! Imax Donor field ending data I-index (integer) !
3756! Jmin Donor field starting data J-index (integer) !
3757! Jmax Donor field ending data J-index (integer) !
3758! Kmin Donor field starting data K-index (integer) !
3759! Kmax Donor field ending data K-index (integer) !
3760! fourth Donor Fouth dimension index to process (integer) !
3761! Loff Fourth dimension couter offset (integer) !
3762! Fdat Packed global donor 3D field data (real 1D array) !
3763! !
3764! On Output: !
3765! !
3766! start Start index where the first of the data values will !
3767! be written along each dimension (integer) !
3768! total Number of data values to be written along each !
3769! dimension (integer) !
3770! Fdat Extracted 3D field data (real) !
3771! !
3772!=======================================================================
3773!
3774! Imported variable declarations.
3775!
3776 integer, intent(in) :: ng, model, tile
3777 integer, intent(in) :: gtype, ifield, tindex
3778 integer, intent(in) :: Imin, Imax, Jmin, Jmax, Kmin, Kmax
3779 integer, intent(in) :: fourth, Loff
3780 integer, intent(in) :: Npts
3781 integer, intent(out) :: start(:), total(:)
3782!
3783 real(r8), intent(inout) :: Fdat(:)
3784!
3785! Local variable declarations.
3786!
3787 logical :: LandFill
3788!
3789 integer :: ghost, i, j, k, ij, ijk, ic, jc, ifactor, method
3790 integer :: Cgrid, Ilen, Jlen, Ioff, Joff, Koff
3791 integer :: Isize, IJsize, Jsize, Ksize
3792 integer :: Istr, Iend, Jstr, Jend
3793 integer :: LBi, UBi, LBj, UBj
3794!
3795 real(r8) :: Fmin, Fmax
3796
3797 real(r8) :: Fwrk(Npts)
3798!
3799 character (len=*), parameter :: MyFile = &
3800 & __FILE__//", interp_field4d_global"
3801!
3802 sourcefile=myfile
3803!
3804!------------------------------------------------------------------------
3805! Interpolate.
3806!------------------------------------------------------------------------
3807!
3808 landfill=.true.
3809 ghost=0
3810 method=bilinear
3811!
3812 isize=imax-imin+1
3813 jsize=jmax-jmin+1
3814!
3815 lbi=xtr_bounds(ng)%LBi(tile)
3816 ubi=xtr_bounds(ng)%UBi(tile)
3817 lbj=xtr_bounds(ng)%LBj(tile)
3818 ubj=xtr_bounds(ng)%UBj(tile)
3819!
3820 SELECT CASE (gtype)
3821!
3822! ---------------------
3823 CASE (r3dvar, w3dvar) ! RHO and W points field
3824! ---------------------
3825!
3826 cgrid=2
3827 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
3828 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
3829 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
3830 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
3831!
3832 fwrk=fdat(1:npts)
3833!
3834 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
3835 & method, landfill, &
3836 & imin, imax, jmin, jmax, kmin, kmax, &
3837 & grid(ng) % Gx_rho, &
3838 & grid(ng) % Gy_rho, &
3839# ifdef MASKING
3840 & grid(ng) % Gmask_rho, &
3841# endif
3842 & fwrk, &
3843 & istr, iend, jstr, jend, &
3844 & istr, iend, jstr, jend, &
3845 & extract(ng) % Iout_rho, &
3846 & extract(ng) % Jout_rho, &
3847 & extract(ng) % Gx_rho, &
3848 & extract(ng) % Gy_rho, &
3849# ifdef MASKING
3850 & extract(ng) % Gmask_rho, &
3851# endif
3852 & npts, fdat, fmin, fmax)
3853!
3854! ---------------------
3855 CASE (u3dvar) ! U points field
3856! ---------------------
3857!
3858 cgrid=3
3859 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
3860 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
3861 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
3862 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
3863!
3864 fwrk=fdat(1:npts)
3865!
3866 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
3867 & method, landfill, &
3868 & imin, imax, jmin, jmax, kmin, kmax, &
3869 & grid(ng) % Gx_u, &
3870 & grid(ng) % Gy_u, &
3871# ifdef MASKING
3872 & grid(ng) % Gmask_u, &
3873# endif
3874 & fwrk, &
3875 & istr, iend, jstr, jend, &
3876 & istr, iend, jstr, jend, &
3877 & extract(ng) % Iout_u, &
3878 & extract(ng) % Jout_u, &
3879 & extract(ng) % Gx_u, &
3880 & extract(ng) % Gy_u, &
3881# ifdef MASKING
3882 & extract(ng) % Gmask_u, &
3883# endif
3884 & npts, fdat, fmin, fmax)
3885!
3886! ---------------------
3887 CASE (v3dvar) ! V points field
3888! ---------------------
3889!
3890 cgrid=4
3891 istr=xtr_bounds(ng)%Imin(cgrid,ghost,tile)
3892 iend=xtr_bounds(ng)%Imax(cgrid,ghost,tile)
3893 jstr=xtr_bounds(ng)%Jmin(cgrid,ghost,tile)
3894 jend=xtr_bounds(ng)%Jmax(cgrid,ghost,tile)
3895!
3896 fwrk=fdat(1:npts)
3897!
3898 CALL regrid_field3d (ng, model, tile, gtype, ifield, &
3899 & method, landfill, &
3900 & 1, isize, 1, jsize, kmin, kmax, &
3901 & grid(ng) % Gx_v, &
3902 & grid(ng) % Gy_v, &
3903# ifdef MASKING
3904 & grid(ng) % Gmask_v, &
3905# endif
3906 & fwrk, &
3907 & istr, iend, jstr, jend, &
3908 & istr, iend, jstr, jend, &
3909 & extract(ng) % Iout_v, &
3910 & extract(ng) % Jout_v, &
3911 & extract(ng) % Gx_v, &
3912 & extract(ng) % Gy_v, &
3913# ifdef MASKING
3914 & extract(ng) % Gmask_v, &
3915# endif
3916 & npts, fdat, fmin, fmax)
3917!
3918! ------------
3919 CASE DEFAULT
3920! ------------
3921!
3922 IF (master) THEN
3923 WRITE (stdout,10) gtype, &
3924 & 'not supported for interpolation:', &
3925 & trim(vname(1,ifield))
3926 END IF
3927 exit_flag=3
3928!
3929 END SELECT
3930 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
3931!
3932! Set start and total vectors needed to write into output NetCDF file.
3933!
3934 IF (gtype.gt.0) THEN
3935 start(1)=1
3936 total(1)=iend-istr+1
3937 start(2)=1
3938 total(2)=jend-jstr+1
3939 start(3)=1
3940 total(3)=ksize
3941 start(4)=fourth+loff
3942 total(4)=1
3943 start(5)=tindex
3944 total(5)=1
3945# ifdef MASKING
3946 ELSE
3947 start(1)=1+(fourth+loff-1)*npts
3948 total(1)=npts
3949 start(2)=tindex
3950 total(2)=1
3951# endif
3952 END IF
3953!
3954 10 FORMAT (' INTERP_FIELD4D_GLOBAL - Staggered variable, gtype = ', &
3955 & i0,/,18x,a,2x,a)
3956!
3957 RETURN
3958 END SUBROUTINE interp_field4d_global
3959!
3960 SUBROUTINE regrid_field2d (ng, model, tile, gtype, ifield, &
3961 & method, LandFill, &
3962 & LBx, UBx, LBy, UBy, &
3963 & Xinp, Yinp, &
3964# ifdef MASKING
3965 & maskInp, &
3966# endif
3967 & Finp, &
3968 & Istr, Iend, Jstr, Jend, &
3969 & LBi, UBi, LBj, UBj, &
3970 & Iout, Jout, &
3971 & Xout, Yout, &
3972# ifdef MASKING
3973 & maskOut, &
3974# endif
3975 & Npts, Fout, Fmin, Fmax)
3976!
3977!=======================================================================
3978! !
3979! It extracts 2D fields from donor grid by interpolation. !
3980! !
3981! On Input: !
3982! !
3983! ng Nested grid number (integer) !
3984! model Calling model identifier (integer) !
3985! tile Domain partition (integer) !
3986! gtype Staggered C-grid type (integer) !
3987! ifield Field metadata index (integer) !
3988! method Interpolation method (integer) !
3989! LandFill Switch to set fill value on land points (logical) !
3990! LBx Donor data I-dimension lower bound (integer) !
3991! UBx Donor data I-dimension upper bound (integer) !
3992! LBy Donor data J-dimension lower bound (integer) !
3993! UBy Donor data J-dimension upper bound (integer) !
3994! Xinp Donor data X-locations (real array) !
3995! Yinp Donor data Y-locations (real array) !
3996# ifdef MASKING
3997! maskInp Donor data land/sea masking (real array) !
3998# endif
3999! Finp Donor 2D field to interpolate from (real array) !
4000! Istr Interpolsted field Starting tile I-index (integer) !
4001! Iend Interpolated field ending tile I-index (integer) !
4002! Jstr Interpolated field starting tile J-index (integer) !
4003! Jend Interpolated field ending tile J-index (integer) !
4004! LBi Interpolated field I-dimension Lower bound (integer) !
4005! UBi Interpolated field I-dimension Upper bound (integer) !
4006! LBj Interpolated field J-dimension Lower bound (integer) !
4007! UBj Interpolated field J-dimension Upper bound (integer) !
4008! Iout I-fractional Xinp grid cell containing Xout (real) !
4009! Jout J-fractional Yinp grid cell containing Yout (real) !
4010! Xout Interpolated field X-locations (real array) !
4011! Yout Interpolated field Y-locations (real array) !
4012# ifdef MASKING
4013! maskOut Interpolated field land/sea masking (real array) !
4014# endif
4015! !
4016! On Output: !
4017! !
4018! Npts Number of points processed in Fout (integer) !
4019! Fout Interpolated 2D field data (real array) !
4020! Fmin Interpolated field minimum value (real) !
4021! Fmax Interpolated field maximum value (real) !
4022! !
4023!=======================================================================
4024!
4025! Imported variable declarations.
4026!
4027 logical, intent(in) :: LandFill
4028!
4029 integer, intent(in) :: ng, model, tile, gtype, ifield, method
4030 integer, intent(in) :: LBx, UBx, LBy, UBy
4031 integer, intent(in) :: LBi, UBi, LBj, UBj
4032 integer, intent(in) :: Istr, Iend, Jstr, Jend
4033 integer, intent(in) :: Npts
4034!
4035 real(r8), intent(in) :: Xinp(LBx:UBx,LBy:UBy)
4036 real(r8), intent(in) :: Yinp(LBx:UBx,LBy:UBy)
4037# ifdef MASKING
4038 real(r8), intent(in) :: maskInp(LBx:UBx,LBy:UBy)
4039# endif
4040 real(r8), intent(in) :: Finp(LBx:UBx,LBy:UBy)
4041 real(r8), intent(in) :: Iout(LBi:UBi,LBj:UBj)
4042 real(r8), intent(in) :: Jout(LBi:UBi,LBj:UBj)
4043 real(r8), intent(in) :: Xout(LBi:UBi,LBj:UBj)
4044 real(r8), intent(in) :: Yout(LBi:UBi,LBj:UBj)
4045# ifdef MASKING
4046 real(r8), intent(in) :: maskOut(LBi:UBi,LBj:UBj)
4047# endif
4048 real(r8), intent(out) :: Fout(LBi:UBi,LBj:UBj)
4049 real(r8), intent(out) :: Fmin, Fmax
4050!
4051! Local variable declarations.
4052!
4053 integer :: i, ic, j
4054!
4055 character (len=*), parameter :: MyFile = &
4056 & __FILE__//", regrid_field2d"
4057!
4058!-----------------------------------------------------------------------
4059! Interpolate 2D field.
4060!-----------------------------------------------------------------------
4061!
4062 SELECT CASE (abs(method))
4063 CASE (bilinear)
4064 CALL linterp2d (ng, &
4065 & lbx, ubx, lby, uby, &
4066 & xinp, yinp, finp, &
4067 & lbi, ubi, lbj, ubj, &
4068 & istr, iend, jstr, jend, &
4069 & iout, jout, &
4070 & xout, yout, &
4071 & fout)
4072 CASE (bicubic)
4073 CALL cinterp2d (ng, &
4074 & lbx, ubx, lby, uby, &
4075 & xinp, yinp, finp, &
4076 & lbi, ubi, lbj, ubj, &
4077 & istr, iend, jstr, jend, &
4078 & iout, jout, &
4079 & xout, yout, &
4080 & fout)
4081 CASE DEFAULT
4082 IF (master) WRITE(stdout,10) method
4083 exit_flag=5
4084 END SELECT
4085 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
4086!
4087!-----------------------------------------------------------------------
4088! If requested, set fill value in land areas.
4089!-----------------------------------------------------------------------
4090!
4091 fmin=spval
4092 fmax=-spval
4093 ic=0
4094 DO j=jstr,jend
4095 DO i=istr,iend
4096 ic=ic+1
4097# ifdef MASKING
4098 IF ((maskout(i,j).lt.1.0_r8).and.landfill) THEN
4099 fout(i,j)=spval
4100 ELSE
4101 fmin=min(fmin,fout(i,j))
4102 fmax=max(fmax,fout(i,j))
4103 END IF
4104# else
4105 fmin=min(fmin,fout(i,j))
4106 fmax=max(fmax,fout(i,j))
4107# endif
4108 END DO
4109 END DO
4110!
4111 10 FORMAT (' REGRID_FIEL2D - Illegal interpolation method =', i0)
4112!
4113 RETURN
4114 END SUBROUTINE regrid_field2d
4115!
4116 SUBROUTINE regrid_field3d (ng, model, tile, gtype, ifield, &
4117 & method, LandFill, &
4118 & LBx, UBx, LBy, UBy, Kmin, Kmax, &
4119 & Xinp, Yinp, &
4120# ifdef MASKING
4121 & maskInp, &
4122# endif
4123 & Finp, &
4124 & Istr, Iend, Jstr, Jend, &
4125 & LBi, UBi, LBj, UBj, &
4126 & Iout, Jout, &
4127 & Xout, Yout, &
4128# ifdef MASKING
4129 & maskOut, &
4130# endif
4131 & Npts, Fout, Fmin, Fmax)
4132!
4133!=======================================================================
4134! !
4135! It extracts 3D fields from donor grid by interpolation level-by- !
4136! level. !
4137! !
4138! On Input: !
4139! !
4140! ng Nested grid number (integer) !
4141! model Calling model identifier (integer) !
4142! tile Domain partition (integer) !
4143! gtype Staggered C-grid type (integer) !
4144! ifield Field metadata index (integer) !
4145! method Interpolation method (integer) !
4146! LandFill Switch to set fill value on land points (logical) !
4147! LBx Donor data I-dimension lower bound (integer) !
4148! UBx Donor data I-dimension upper bound (integer) !
4149! LBy Donor data J-dimension lower bound (integer) !
4150! UBy Donor data J-dimension upper bound (integer) !
4151! Kmin Donor data K-dimension lower bound (integer) !
4152! Kmax Donor data K-dimension upper bound (integer) !
4153! Xinp Donor data X-locations (real array) !
4154! Yinp Donor data Y-locations (real array) !
4155# ifdef MASKING
4156! maskInp Donor data land/sea masking (real array) !
4157# endif
4158! Finp Donor 3D field to interpolate from (real array) !
4159! Istr Interpolsted field Starting tile I-index (integer) !
4160! Iend Interpolated field ending tile I-index (integer) !
4161! Jstr Interpolated field starting tile J-index (integer) !
4162! Jend Interpolated field ending tile J-index (integer) !
4163! LBi Interpolated field I-dimension Lower bound (integer) !
4164! UBi Interpolated field I-dimension Upper bound (integer) !
4165! LBj Interpolated field J-dimension Lower bound (integer) !
4166! UBj Interpolated field J-dimension Upper bound (integer) !
4167! Iout I-fractional Xinp grid cell containing Xout (real) !
4168! Jout J-fractional Yinp grid cell containing Yout (real) !
4169! Xout Interpolated field X-locations (real array) !
4170! Yout Interpolated field Y-locations (real array) !
4171# ifdef MASKING
4172! maskOut Interpolated field land/sea masking (real array) !
4173# endif
4174! !
4175! On Output: !
4176! !
4177! Npts Number of points processed in Fout (integer) !
4178! Fout Interpolated 3D field data (real array) !
4179! Fmin Interpolated field minimum value (real) !
4180! Fmax Interpolated field maximum value (real) !
4181! !
4182!=======================================================================
4183!
4184! Imported variable declarations.
4185!
4186 logical, intent(in) :: landFill
4187!
4188 integer, intent(in) :: ng, model, tile, gtype, ifield, method
4189 integer, intent(in) :: LBx, UBx, LBy, UBy
4190 integer, intent(in) :: Kmin, Kmax
4191 integer, intent(in) :: Istr, Iend, Jstr, Jend
4192 integer, intent(in) :: LBi, UBi, LBj, UBj
4193 integer, intent(in) :: Npts
4194!
4195 real(r8), intent(in) :: Xinp(LBx:UBx,LBy:UBy)
4196 real(r8), intent(in) :: Yinp(LBx:UBx,LBy:UBy)
4197# ifdef MASKING
4198 real(r8), intent(in) :: maskInp(LBx:UBx,LBy:UBy)
4199# endif
4200 real(r8), intent(in) :: Finp(LBx:UBx,LBy:UBy,Kmin:Kmax)
4201 real(r8), intent(in) :: Iout(LBi:UBi,LBj:UBj)
4202 real(r8), intent(in) :: Jout(LBi:UBi,LBj:UBj)
4203 real(r8), intent(in) :: Xout(LBi:UBi,LBj:UBj)
4204 real(r8), intent(in) :: Yout(LBi:UBi,LBj:UBj)
4205# ifdef MASKING
4206 real(r8), intent(in) :: maskOut(LBi:UBi,LBj:UBj)
4207# endif
4208 real(r8), intent(out) :: Fout(LBi:UBi,LBj:UBj,Kmin:Kmax)
4209 real(r8), intent(out) :: Fmin, Fmax
4210!
4211! Local variable declarations.
4212!
4213 integer :: i, ic, j, k
4214!
4215 character (len=*), parameter :: MyFile = &
4216 & __FILE__//", regrid_field3d"
4217!
4218!-----------------------------------------------------------------------
4219! Interpolate 2D field.
4220!-----------------------------------------------------------------------
4221!
4222 SELECT CASE (abs(method))
4223 CASE (bilinear)
4224 DO k=kmin,kmax
4225 CALL linterp2d (ng, &
4226 & lbx, ubx, lby, uby, &
4227 & xinp, yinp, finp(:,:,k), &
4228 & lbi, ubi, lbj, ubj, &
4229 & istr, iend, jstr, jend, &
4230 & iout, jout, &
4231 & xout, yout, &
4232 & fout(lbi:,lbj:,k))
4233 END DO
4234 CASE (bicubic)
4235 DO k=kmin,kmax
4236 CALL cinterp2d (ng, &
4237 & lbx, ubx, lby, uby, &
4238 & xinp, yinp, finp(:,:,k), &
4239 & lbi, ubi, lbj, ubj, &
4240 & istr, iend, jstr, jend, &
4241 & iout, jout, &
4242 & xout, yout, &
4243 & fout(lbi:,lbj:,k))
4244 END DO
4245 CASE DEFAULT
4246 IF (master) WRITE(stdout,10) method
4247 exit_flag=5
4248 END SELECT
4249 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
4250!
4251!-----------------------------------------------------------------------
4252! If requested, set fill value in land areas.
4253!-----------------------------------------------------------------------
4254!
4255 fmin=spval
4256 fmax=-spval
4257 ic=0
4258 DO k=kmin,kmax
4259 DO j=jstr,jend
4260 DO i=istr,iend
4261 ic=ic+1
4262# ifdef MASKING
4263 IF ((maskout(i,j).lt.1.0_r8).and.landfill) THEN
4264 fout(i,j,k)=spval
4265 ELSE
4266 fmin=min(fmin,fout(i,j,k))
4267 fmax=max(fmax,fout(i,j,k))
4268 END IF
4269# else
4270 fmin=min(fmin,fout(i,j,k))
4271 fmax=max(fmax,fout(i,j,k))
4272# endif
4273 END DO
4274 END DO
4275 END DO
4276!
4277 10 FORMAT (' REGRID_FIEL3D - Illegal interpolation method =', i0)
4278!
4279 RETURN
4280 END SUBROUTINE regrid_field3d
4281!
4282#endif
4283 END MODULE extract_field_mod
subroutine hindices(ng, lbi, ubi, lbj, ubj, is, ie, js, je, angler, xgrd, ygrd, lbm, ubm, lbn, ubn, ms, me, ns, ne, xpos, ypos, ipos, jpos, ijspv, rectangular)
subroutine cinterp2d(ng, lbx, ubx, lby, uby, xinp, yinp, finp, lbi, ubi, lbj, ubj, istr, iend, jstr, jend, iout, jout, xout, yout, fout, minval, maxval)
subroutine linterp2d(ng, lbx, ubx, lby, uby, xinp, yinp, finp, lbi, ubi, lbj, ubj, istr, iend, jstr, jend, iout, jout, xout, yout, fout, minval, maxval)
subroutine mp_gather2d(ng, model, lbi, ubi, lbj, ubj, tindex, gtype, ascl, amask, a, npts, awrk, setfillval)
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
integer stdout
character(len=256) sourcefile
character(len=maxlen), dimension(6, 0:nv) vname
logical master
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
integer, parameter r3dvar
Definition mod_param.F:721
type(t_iobounds), dimension(:), allocatable iobounds
Definition mod_param.F:282
integer, parameter u3dvar
Definition mod_param.F:722
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, parameter u2dvar
Definition mod_param.F:718
integer, parameter w3dvar
Definition mod_param.F:724
integer, parameter p2dvar
Definition mod_param.F:716
integer, parameter r2dvar
Definition mod_param.F:717
integer, parameter v2dvar
Definition mod_param.F:719
integer, parameter p3dvar
Definition mod_param.F:720
integer, parameter v3dvar
Definition mod_param.F:723
logical spherical
real(dp), parameter spval
integer, parameter iwest
integer exit_flag
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth
integer noerror
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52