updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / wrftladj / module_bc_tl.F
blob17289a1fb7ef75a74d9edb55822355d774320d7e
1 !WRF+/TL:MODEL_LAYER:BOUNDARY
2 !Created by Ning Pan, 2010-08
5 MODULE g_module_bc
7    USE module_configure
8    USE module_wrf_error
9    USE module_model_constants
11    IMPLICIT NONE
13 !  set the bdyzone.  We are hardwiring this here and we'll
14 !  decide later where it should be set and stored
16    INTEGER, PARAMETER            :: bdyzone = 4
17    INTEGER, PARAMETER            :: bdyzone_x = bdyzone
18    INTEGER, PARAMETER            :: bdyzone_y = bdyzone
20 CONTAINS
22 !-----------------------------------
24    SUBROUTINE g_set_physical_bc2d( dat,g_dat, variable_in,  &
25                                  config_flags,           & 
26                                  ids,ide, jds,jde,   & ! domain dims
27                                  ims,ime, jms,jme,   & ! memory dims
28                                  ips,ipe, jps,jpe,   & ! patch  dims
29                                  its,ite, jts,jte   )      
31       IMPLICIT NONE
33       INTEGER,      INTENT(IN   )    :: ids,ide, jds,jde
34       INTEGER,      INTENT(IN   )    :: ims,ime, jms,jme
35       INTEGER,      INTENT(IN   )    :: ips,ipe, jps,jpe
36       INTEGER,      INTENT(IN   )    :: its,ite, jts,jte
37       CHARACTER,    INTENT(IN   )    :: variable_in
39       CHARACTER                      :: variable
41       REAL,  DIMENSION( ims:ime , jms:jme ) :: dat,g_dat
42       TYPE( grid_config_rec_type ) config_flags
44       INTEGER  :: i, j, istag, jstag, itime
46       LOGICAL  :: debug, open_bc_copy
48 !------------
50       debug = .false.
52       open_bc_copy = .false.
54       variable = variable_in
55       IF ( variable_in .ge. 'A' .and. variable_in .le. 'Z' ) THEN
56         variable = CHAR( ICHAR(variable_in) - ICHAR('A') + ICHAR('a') )
57       ENDIF
58       IF ((variable == 'u') .or. (variable == 'v') .or.  &
59           (variable == 'w') .or. (variable == 't') .or.  &
60           (variable == 'x') .or. (variable == 'y') .or.  &
61           (variable == 'r') .or. (variable == 'p') ) open_bc_copy = .true.
63 !  begin, first set a staggering variable
65       istag = -1
66       jstag = -1
68       IF ((variable == 'u') .or. (variable == 'x')) istag = 0
69       IF ((variable == 'v') .or. (variable == 'y')) jstag = 0
71       if(debug) then
72         write(6,*) ' in bc2d, var is ',variable, istag, jstag
73         write(6,*) ' b.cs are ',  &
74       config_flags%periodic_x,  &
75       config_flags%periodic_y
76       end if
77       
78       IF ( variable == 'd' ) then  !JDM
79          istag = 0
80          jstag = 0
81       ENDIF
82       IF ( variable == 'e' ) then  !JDM
83          istag = 0
84       ENDIF
85       IF ( variable == 'f' ) then  !JDM
86          jstag = 0
87       ENDIF
89       periodicity_x:  IF( ( config_flags%periodic_x ) ) THEN 
90         IF ( ( ids == ips ) .and.  ( ide == ipe ) ) THEN  ! test if east and west both on-processor 
91           IF ( its == ids ) THEN
93             DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
94             DO i = 0,-(bdyzone-1),-1
95               g_dat(ids+i-1,j) = g_dat(ide+i-1,j)
96               dat(ids+i-1,j) = dat(ide+i-1,j)
97             ENDDO
98             ENDDO
100           ENDIF
102           IF ( ite == ide ) THEN
104             DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
105 !!          DO i = 1 , bdyzone
106             DO i = -istag , bdyzone
107               g_dat(ide+i+istag,j) = g_dat(ids+i+istag,j)
108               dat(ide+i+istag,j) = dat(ids+i+istag,j)
109             ENDDO
110             ENDDO
112           ENDIF
113         ENDIF
115       ELSE 
117         symmetry_xs: IF( ( config_flags%symmetric_xs ) .and.  &
118                          ( its == ids )                  )  THEN
120           IF ( (variable /= 'u') .and. (variable /= 'x') ) THEN
122             DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
123             DO i = 1, bdyzone
124               g_dat(ids-i,j) = g_dat(ids+i-1,j)
125               dat(ids-i,j) = dat(ids+i-1,j)
126             ENDDO
127             ENDDO
129           ELSE
131             IF( variable == 'u' ) THEN
133               DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
134               DO i = 0, bdyzone-1
135                 g_dat(ids-i,j) = - g_dat(ids+i,j)
136                 dat(ids-i,j) = - dat(ids+i,j)
137               ENDDO
138               ENDDO
140             ELSE
142               DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
143               DO i = 0, bdyzone-1
144                 g_dat(ids-i,j) =   g_dat(ids+i,j)
145                 dat(ids-i,j) =   dat(ids+i,j)
146               ENDDO
147               ENDDO
149             END IF
151           ENDIF
153         ENDIF symmetry_xs
156 !  now the symmetry boundary at xe
158         symmetry_xe: IF( ( config_flags%symmetric_xe ) .and.  &
159                          ( ite == ide )                  )  THEN
161           IF ( (variable /= 'u') .and. (variable /= 'x') ) THEN
163             DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
164             DO i = 1, bdyzone
165               g_dat(ide+i-1,j) = g_dat(ide-i,j)
166               dat(ide+i-1,j) = dat(ide-i,j)
167             ENDDO
168             ENDDO
170           ELSE
172             IF (variable == 'u' ) THEN
174               DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
175               DO i = 0, bdyzone-1
176                 g_dat(ide+i,j) = - g_dat(ide-i,j)
177                 dat(ide+i,j) = - dat(ide-i,j)
178               ENDDO
179               ENDDO
182             ELSE
184               DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
185               DO i = 0, bdyzone-1
186                 g_dat(ide+i,j) = g_dat(ide-i,j)
187                 dat(ide+i,j) = dat(ide-i,j)
188               ENDDO
189               ENDDO
191             END IF
193           END IF 
195         END IF symmetry_xe
198 !  set open b.c in X copy into boundary zone here.  WCS, 19 March 2000
200         open_xs: IF( ( config_flags%open_xs   .or. &
201                        config_flags%specified .or. &
202                        config_flags%nested            ) .and.  &
203                          ( its == ids ) .and. open_bc_copy  )  THEN
205             DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
206               g_dat(ids-1,j) = g_dat(ids,j)
207               dat(ids-1,j) = dat(ids,j)
208               g_dat(ids-2,j) = g_dat(ids,j)
209               dat(ids-2,j) = dat(ids,j)
210               g_dat(ids-3,j) = g_dat(ids,j)
211               dat(ids-3,j) = dat(ids,j)
212             ENDDO
214         ENDIF open_xs
217 !  now the open boundary copy at xe
219         open_xe: IF( ( config_flags%open_xe   .or. &
220                        config_flags%specified .or. &
221                        config_flags%nested            ) .and.  &
222                           ( ite == ide ) .and. open_bc_copy  )  THEN
224           IF ( variable /= 'u' .and. variable /= 'x') THEN
226             DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
227               g_dat(ide  ,j) = g_dat(ide-1,j) 
228               dat(ide  ,j) = dat(ide-1,j) 
229               g_dat(ide+1,j) = g_dat(ide-1,j) 
230               dat(ide+1,j) = dat(ide-1,j) 
231               g_dat(ide+2,j) = g_dat(ide-1,j) 
232               dat(ide+2,j) = dat(ide-1,j) 
233             ENDDO
235           ELSE
237             DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
238               g_dat(ide+1,j) = g_dat(ide,j)
239               dat(ide+1,j) = dat(ide,j)
240               g_dat(ide+2,j) = g_dat(ide,j)
241               dat(ide+2,j) = dat(ide,j)
242               g_dat(ide+3,j) = g_dat(ide,j)
243               dat(ide+3,j) = dat(ide,j)
244             ENDDO
246           END IF 
248         END IF open_xe
250 !  end open b.c in X copy into boundary zone addition.  WCS, 19 March 2000
252       END IF periodicity_x
254 !  same procedure in y
256       periodicity_y:  IF( ( config_flags%periodic_y ) ) THEN
257         IF ( ( jds == jps ) .and. ( jde == jpe ) )  THEN    ! test of both north and south on processor
259           IF( jts == jds ) then
261             DO j = 0, -(bdyzone-1), -1
262             DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
263               g_dat(i,jds+j-1) = g_dat(i,jde+j-1)
264               dat(i,jds+j-1) = dat(i,jde+j-1)
265             ENDDO
266             ENDDO
268           END IF
270           IF( jte == jde ) then
272             DO j = -jstag, bdyzone
273             DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
274               g_dat(i,jde+j+jstag) = g_dat(i,jds+j+jstag)
275               dat(i,jde+j+jstag) = dat(i,jds+j+jstag)
276             ENDDO
277             ENDDO
279           END IF
281         END IF
283       ELSE
285         symmetry_ys: IF( ( config_flags%symmetric_ys ) .and.  &
286                          ( jts == jds)                   )  THEN
288           IF ( (variable /= 'v') .and. (variable /= 'y') ) THEN
290             DO j = 1, bdyzone
291             DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
292               g_dat(i,jds-j) = g_dat(i,jds+j-1) 
293               dat(i,jds-j) = dat(i,jds+j-1) 
294             ENDDO
295             ENDDO
297           ELSE
299             IF (variable == 'v') THEN
301               DO j = 1, bdyzone
302               DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
303                 g_dat(i,jds-j) = - g_dat(i,jds+j) 
304                 dat(i,jds-j) = - dat(i,jds+j) 
305               ENDDO              
306               ENDDO
308             ELSE
310               DO j = 1, bdyzone
311               DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
312                 g_dat(i,jds-j) = g_dat(i,jds+j) 
313                 dat(i,jds-j) = dat(i,jds+j) 
314               ENDDO              
315               ENDDO
317             END IF
319           ENDIF
321         ENDIF symmetry_ys
323 !  now the symmetry boundary at ye
325         symmetry_ye: IF( ( config_flags%symmetric_ye ) .and.  &
326                          ( jte == jde )                  )  THEN
328           IF ( (variable /= 'v') .and. (variable /= 'y') ) THEN
330             DO j = 1, bdyzone
331             DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
332               g_dat(i,jde+j-1) = g_dat(i,jde-j) 
333               dat(i,jde+j-1) = dat(i,jde-j) 
334             ENDDO                               
335             ENDDO
337           ELSE
339             IF (variable == 'v' ) THEN
341               DO j = 1, bdyzone
342               DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
343                 g_dat(i,jde+j) = - g_dat(i,jde-j)
344                 dat(i,jde+j) = - dat(i,jde-j)
345               ENDDO                               
346               ENDDO
348             ELSE
350               DO j = 1, bdyzone
351               DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
352                 g_dat(i,jde+j) = g_dat(i,jde-j)
353                 dat(i,jde+j) = dat(i,jde-j)
354               ENDDO                               
355               ENDDO
357             END IF
359           ENDIF
361         END IF symmetry_ye
363 !  set open b.c in Y copy into boundary zone here.  WCS, 19 March 2000
365         open_ys: IF( ( config_flags%open_ys   .or. &
366                        config_flags%polar     .or. &
367                        config_flags%specified .or. &
368                        config_flags%nested            ) .and.  &
369                          ( jts == jds) .and. open_bc_copy )  THEN
371             DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
372               g_dat(i,jds-1) = g_dat(i,jds) 
373               dat(i,jds-1) = dat(i,jds) 
374               g_dat(i,jds-2) = g_dat(i,jds) 
375               dat(i,jds-2) = dat(i,jds) 
376               g_dat(i,jds-3) = g_dat(i,jds) 
377               dat(i,jds-3) = dat(i,jds) 
378             ENDDO
380         ENDIF open_ys
382 !  now the open boundary copy at ye
384         open_ye: IF( ( config_flags%open_ye   .or. &
385                        config_flags%polar     .or. &
386                        config_flags%specified .or. &
387                        config_flags%nested            ) .and.  &
388                          ( jte == jde ) .and. open_bc_copy )  THEN
390           IF  (variable /= 'v' .and. variable /= 'y' ) THEN
392             DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
393               g_dat(i,jde  ) = g_dat(i,jde-1) 
394               dat(i,jde  ) = dat(i,jde-1) 
395               g_dat(i,jde+1) = g_dat(i,jde-1) 
396               dat(i,jde+1) = dat(i,jde-1) 
397               g_dat(i,jde+2) = g_dat(i,jde-1) 
398               dat(i,jde+2) = dat(i,jde-1) 
399             ENDDO                               
401           ELSE
403             DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
404               g_dat(i,jde+1) = g_dat(i,jde) 
405               dat(i,jde+1) = dat(i,jde) 
406               g_dat(i,jde+2) = g_dat(i,jde) 
407               dat(i,jde+2) = dat(i,jde) 
408               g_dat(i,jde+3) = g_dat(i,jde) 
409               dat(i,jde+3) = dat(i,jde) 
410             ENDDO                               
412           ENDIF
414         END IF open_ye
415       
416 !  end open b.c in Y copy into boundary zone addition.  WCS, 19 March 2000
418       END IF periodicity_y
420 !  fix corners for doubly periodic domains
422       IF ( config_flags%periodic_x .and. config_flags%periodic_y &
423            .and. (ids == ips) .and. (ide == ipe)                 &
424            .and. (jds == jps) .and. (jde == jpe)                   ) THEN
426          IF ( (its == ids) .and. (jts == jds) ) THEN  ! lower left corner fill
427            DO j = 0, -(bdyzone-1), -1
428            DO i = 0, -(bdyzone-1), -1
429              g_dat(ids+i-1,jds+j-1) = g_dat(ide+i-1,jde+j-1)
430              dat(ids+i-1,jds+j-1) = dat(ide+i-1,jde+j-1)
431            ENDDO
432            ENDDO
433          END IF
435          IF ( (ite == ide) .and. (jts == jds) ) THEN  ! lower right corner fill
436            DO j = 0, -(bdyzone-1), -1
437            DO i = 1, bdyzone
438              g_dat(ide+i+istag,jds+j-1) = g_dat(ids+i+istag,jde+j-1)
439              dat(ide+i+istag,jds+j-1) = dat(ids+i+istag,jde+j-1)
440            ENDDO
441            ENDDO
442          END IF
444          IF ( (ite == ide) .and. (jte == jde) ) THEN  ! upper right corner fill
445            DO j = 1, bdyzone
446            DO i = 1, bdyzone
447              g_dat(ide+i+istag,jde+j+jstag) = g_dat(ids+i+istag,jds+j+jstag)
448              dat(ide+i+istag,jde+j+jstag) = dat(ids+i+istag,jds+j+jstag)
449            ENDDO
450            ENDDO
451          END IF
453          IF ( (its == ids) .and. (jte == jde) ) THEN  ! upper left corner fill
454            DO j = 1, bdyzone
455            DO i = 0, -(bdyzone-1), -1
456              g_dat(ids+i-1,jde+j+jstag) = g_dat(ide+i-1,jds+j+jstag)
457              dat(ids+i-1,jde+j+jstag) = dat(ide+i-1,jds+j+jstag)
458            ENDDO
459            ENDDO
460          END IF
462        END IF
464    END SUBROUTINE g_set_physical_bc2d
466 !-----------------------------------
468    SUBROUTINE g_set_physical_bc3d( dat,g_dat, variable_in,        &
469                                config_flags,                   & 
470                                ids,ide, jds,jde, kds,kde,  & ! domain dims
471                                ims,ime, jms,jme, kms,kme,  & ! memory dims
472                                ips,ipe, jps,jpe, kps,kpe,  & ! patch  dims
473                                its,ite, jts,jte, kts,kte )
475       IMPLICIT NONE
477       INTEGER,      INTENT(IN   )    :: ids,ide, jds,jde, kds,kde
478       INTEGER,      INTENT(IN   )    :: ims,ime, jms,jme, kms,kme
479       INTEGER,      INTENT(IN   )    :: ips,ipe, jps,jpe, kps,kpe
480       INTEGER,      INTENT(IN   )    :: its,ite, jts,jte, kts,kte
481       CHARACTER,    INTENT(IN   )    :: variable_in
483       CHARACTER                      :: variable
485       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ) :: dat,g_dat
486       TYPE( grid_config_rec_type ) config_flags
488       INTEGER  :: i, j, k, istag, jstag, itime, k_end
490       LOGICAL  :: debug, open_bc_copy
492 !------------
494       debug = .false.
496       open_bc_copy = .false.
498       variable = variable_in
499       IF ( variable_in .ge. 'A' .and. variable_in .le. 'Z' ) THEN
500         variable = CHAR( ICHAR(variable_in) - ICHAR('A') + ICHAR('a') )
501       ENDIF
503       IF ((variable == 'u') .or. (variable == 'v') .or.     &
504           (variable == 'w') .or. (variable == 't') .or.     &
505           (variable == 'd') .or. (variable == 'e') .or. &
506           (variable == 'x') .or. (variable == 'y') .or. &
507           (variable == 'f') .or. (variable == 'r') .or. &
508           (variable == 'p')                        ) open_bc_copy = .true.
510 !  begin, first set a staggering variable
512       istag = -1
513       jstag = -1
514       k_end = max(1,min(kde-1,kte))
517       IF ((variable == 'u') .or. (variable == 'x')) istag = 0
518       IF ((variable == 'v') .or. (variable == 'y')) jstag = 0
519       IF ((variable == 'd') .or. (variable == 'xy')) then
520          istag = 0
521          jstag = 0
522       ENDIF
523       IF ((variable == 'e') ) then
524          istag = 0
525          k_end = min(kde,kte)
526       ENDIF
528       IF ((variable == 'f') ) then
529          jstag = 0
530          k_end = min(kde,kte)
531       ENDIF
533       IF ( variable == 'w')  k_end = min(kde,kte)
535 !      k_end = kte
537       if(debug) then
538         write(6,*) ' in bc, var is ',variable, istag, jstag, kte, k_end
539         write(6,*) ' b.cs are ',  &
540       config_flags%periodic_x,  &
541       config_flags%periodic_y
542       end if
544       periodicity_x:  IF( ( config_flags%periodic_x ) ) THEN
546         IF ( ( ids == ips ) .and. ( ide == ipe ) ) THEN  ! test if both east and west on-processor
547           IF ( its == ids ) THEN
549             DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
550             DO k = kts, k_end
551             DO i = 0,-(bdyzone-1),-1
552               g_dat(ids+i-1,k,j) = g_dat(ide+i-1,k,j)
553               dat(ids+i-1,k,j) = dat(ide+i-1,k,j)
554             ENDDO
555             ENDDO
556             ENDDO
558           ENDIF
561           IF ( ite == ide ) THEN
563             DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
564             DO k = kts, k_end
565             DO i = -istag , bdyzone
566               g_dat(ide+i+istag,k,j) = g_dat(ids+i+istag,k,j)
567               dat(ide+i+istag,k,j) = dat(ids+i+istag,k,j)
568             ENDDO
569             ENDDO
570             ENDDO
572           ENDIF
574         ENDIF
576       ELSE 
578         symmetry_xs: IF( ( config_flags%symmetric_xs ) .and.  &
579                          ( its == ids )                  )  THEN
581           IF ( (variable /= 'u') .and. (variable /= 'x') ) THEN
583             DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
584             DO k = kts, k_end
585             DO i = 1, bdyzone
586               g_dat(ids-i,k,j) = g_dat(ids+i-1,k,j)
587               dat(ids-i,k,j) = dat(ids+i-1,k,j)
588             ENDDO
589             ENDDO
590             ENDDO
592           ELSE
594             IF ( variable == 'u' ) THEN
596               DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
597               DO k = kts, k_end
598               DO i = 1, bdyzone
599                 g_dat(ids-i,k,j) = - g_dat(ids+i,k,j)
600                 dat(ids-i,k,j) = - dat(ids+i,k,j)
601               ENDDO
602               ENDDO
603               ENDDO
605             ELSE
607               DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
608               DO k = kts, k_end
609               DO i = 1, bdyzone
610                 g_dat(ids-i,k,j) = g_dat(ids+i,k,j)
611                 dat(ids-i,k,j) = dat(ids+i,k,j)
612               ENDDO
613               ENDDO
614               ENDDO
616             END IF
618           ENDIF
620         ENDIF symmetry_xs
623 !  now the symmetry boundary at xe
625         symmetry_xe: IF( ( config_flags%symmetric_xe ) .and.  &
626                          ( ite == ide )                  )  THEN
628           IF ( (variable /= 'u') .and. (variable /= 'x') ) THEN
630             DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
631             DO k = kts, k_end
632             DO i = 1, bdyzone
633               g_dat(ide+i-1,k,j) = g_dat(ide-i,k,j)
634               dat(ide+i-1,k,j) = dat(ide-i,k,j)
635             ENDDO
636             ENDDO
637             ENDDO
639           ELSE
641             IF (variable == 'u') THEN
643               DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
644               DO k = kts, k_end
645               DO i = 1, bdyzone
646                 g_dat(ide+i,k,j) = - g_dat(ide-i,k,j)
647                 dat(ide+i,k,j) = - dat(ide-i,k,j)
648               ENDDO
649               ENDDO
650               ENDDO
652             ELSE
654               DO j = MAX(jds,jts-1), MIN(jte+1,jde+jstag)
655               DO k = kts, k_end
656               DO i = 1, bdyzone
657                 g_dat(ide+i,k,j) = g_dat(ide-i,k,j)
658                 dat(ide+i,k,j) = dat(ide-i,k,j)
659               ENDDO
660               ENDDO
661               ENDDO
663              END IF
665           END IF 
667         END IF symmetry_xe
669 !  set open b.c in X copy into boundary zone here.  WCS, 19 March 2000
671         open_xs: IF( ( config_flags%open_xs   .or. &
672                        config_flags%specified .or. &
673                        config_flags%nested            ) .and.  &
674                          ( its == ids ) .and. open_bc_copy  )  THEN
676             DO j = jts-bdyzone, MIN(jte,jde+jstag)+bdyzone
677             DO k = kts, k_end
678               g_dat(ids-1,k,j) = g_dat(ids,k,j)
679               dat(ids-1,k,j) = dat(ids,k,j)
680               g_dat(ids-2,k,j) = g_dat(ids,k,j)
681               dat(ids-2,k,j) = dat(ids,k,j)
682               g_dat(ids-3,k,j) = g_dat(ids,k,j)
683               dat(ids-3,k,j) = dat(ids,k,j)
684             ENDDO
685             ENDDO
687         ENDIF open_xs
690 !  now the open_xe boundary copy 
692         open_xe: IF( ( config_flags%open_xe   .or. &
693                        config_flags%specified .or. &
694                        config_flags%nested            ) .and.  &
695                          ( ite == ide ) .and. open_bc_copy )  THEN
697           IF (variable /= 'u' .and. variable /= 'x' ) THEN
699             DO j = jts-bdyzone, MIN(jte,jde+jstag)+bdyzone
700             DO k = kts, k_end
701               g_dat(ide  ,k,j) = g_dat(ide-1,k,j)
702               dat(ide  ,k,j) = dat(ide-1,k,j)
703               g_dat(ide+1,k,j) = g_dat(ide-1,k,j)
704               dat(ide+1,k,j) = dat(ide-1,k,j)
705               g_dat(ide+2,k,j) = g_dat(ide-1,k,j)
706               dat(ide+2,k,j) = dat(ide-1,k,j)
707             ENDDO
708             ENDDO
710           ELSE
712 !!!!!!! I am not sure about this one!  JM 20020402
713             DO j = MAX(jds,jts-1)-bdyzone, MIN(jte+1,jde+jstag)+bdyzone
714             DO k = kts, k_end
715               g_dat(ide+1,k,j) = g_dat(ide,k,j)
716               dat(ide+1,k,j) = dat(ide,k,j)
717               g_dat(ide+2,k,j) = g_dat(ide,k,j)
718               dat(ide+2,k,j) = dat(ide,k,j)
719               g_dat(ide+3,k,j) = g_dat(ide,k,j)
720               dat(ide+3,k,j) = dat(ide,k,j)
721             ENDDO
722             ENDDO
724           END IF 
726         END IF open_xe
728 !  end open b.c in X copy into boundary zone addition.  WCS, 19 March 2000
730       END IF periodicity_x
732 !  same procedure in y
734       periodicity_y:  IF( ( config_flags%periodic_y ) ) THEN
735         IF ( ( jds == jps ) .and. ( jde == jpe ) )  THEN      ! test if both north and south on processor
736           IF( jts == jds ) then
738             DO j = 0, -(bdyzone-1), -1
739             DO k = kts, k_end
740             DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
741               g_dat(i,k,jds+j-1) = g_dat(i,k,jde+j-1)
742               dat(i,k,jds+j-1) = dat(i,k,jde+j-1)
743             ENDDO
744             ENDDO
745             ENDDO
747           END IF
749           IF( jte == jde ) then
751             DO j = -jstag, bdyzone
752             DO k = kts, k_end
753             DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
754               g_dat(i,k,jde+j+jstag) = g_dat(i,k,jds+j+jstag)
755               dat(i,k,jde+j+jstag) = dat(i,k,jds+j+jstag)
756             ENDDO
757             ENDDO
758             ENDDO
760           END IF
762         END IF
764       ELSE
766         symmetry_ys: IF( ( config_flags%symmetric_ys ) .and.  &
767                          ( jts == jds)                   )  THEN
769           IF ( (variable /= 'v') .and. (variable /= 'y') ) THEN
771             DO j = 1, bdyzone
772             DO k = kts, k_end
773             DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
774               g_dat(i,k,jds-j) = g_dat(i,k,jds+j-1) 
775               dat(i,k,jds-j) = dat(i,k,jds+j-1) 
776             ENDDO                               
777             ENDDO
778             ENDDO
780           ELSE
782             IF (variable == 'v') THEN
784               DO j = 1, bdyzone
785               DO k = kts, k_end
786               DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
787                 g_dat(i,k,jds-j) = - g_dat(i,k,jds+j) 
788                 dat(i,k,jds-j) = - dat(i,k,jds+j) 
789               ENDDO              
790               ENDDO
791               ENDDO
793             ELSE
795               DO j = 1, bdyzone
796               DO k = kts, k_end
797               DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
798                 g_dat(i,k,jds-j) = g_dat(i,k,jds+j) 
799                 dat(i,k,jds-j) = dat(i,k,jds+j) 
800               ENDDO              
801               ENDDO
802               ENDDO
804             END IF
806           ENDIF
808         ENDIF symmetry_ys
810 !  now the symmetry boundary at ye
812         symmetry_ye: IF( ( config_flags%symmetric_ye ) .and.  &
813                          ( jte == jde )                  )  THEN
815           IF ( (variable /= 'v') .and. (variable /= 'y') ) THEN
817             DO j = 1, bdyzone
818             DO k = kts, k_end
819             DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
820               g_dat(i,k,jde+j-1) = g_dat(i,k,jde-j) 
821               dat(i,k,jde+j-1) = dat(i,k,jde-j) 
822             ENDDO                               
823             ENDDO
824             ENDDO
826           ELSE
828             IF ( variable == 'v' ) THEN
830               DO j = 1, bdyzone
831               DO k = kts, k_end
832               DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
833                 g_dat(i,k,jde+j) = - g_dat(i,k,jde-j) 
834                 dat(i,k,jde+j) = - dat(i,k,jde-j) 
835               ENDDO                               
836               ENDDO
837               ENDDO
839             ELSE
841               DO j = 1, bdyzone
842               DO k = kts, k_end
843               DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
844                 g_dat(i,k,jde+j) = g_dat(i,k,jde-j) 
845                 dat(i,k,jde+j) = dat(i,k,jde-j) 
846               ENDDO                               
847               ENDDO
848               ENDDO
850             END IF
852           ENDIF
854         END IF symmetry_ye
855       
856 !  set open b.c in Y copy into boundary zone here.  WCS, 19 March 2000
858         open_ys: IF( ( config_flags%open_ys   .or. &
859                        config_flags%polar     .or. &
860                        config_flags%specified .or. &
861                        config_flags%nested            ) .and.  &
862                          ( jts == jds) .and. open_bc_copy )  THEN
864             DO k = kts, k_end
865             DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
866               g_dat(i,k,jds-1) = g_dat(i,k,jds) 
867               dat(i,k,jds-1) = dat(i,k,jds) 
868               g_dat(i,k,jds-2) = g_dat(i,k,jds) 
869               dat(i,k,jds-2) = dat(i,k,jds) 
870               g_dat(i,k,jds-3) = g_dat(i,k,jds) 
871               dat(i,k,jds-3) = dat(i,k,jds) 
872             ENDDO
873             ENDDO
875         ENDIF open_ys
877 !  now the open boundary copy at ye
879         open_ye: IF( ( config_flags%open_ye   .or. &
880                        config_flags%polar     .or. &
881                        config_flags%specified .or. &
882                        config_flags%nested            ) .and.  &
883                          ( jte == jde ) .and. open_bc_copy )  THEN
885           IF (variable /= 'v' .and. variable /= 'y' ) THEN
887             DO k = kts, k_end
888             DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
889               g_dat(i,k,jde  ) = g_dat(i,k,jde-1) 
890               dat(i,k,jde  ) = dat(i,k,jde-1) 
891               g_dat(i,k,jde+1) = g_dat(i,k,jde-1) 
892               dat(i,k,jde+1) = dat(i,k,jde-1) 
893               g_dat(i,k,jde+2) = g_dat(i,k,jde-1) 
894               dat(i,k,jde+2) = dat(i,k,jde-1) 
895             ENDDO                               
896             ENDDO
898           ELSE
900             DO k = kts, k_end
901             DO i = MAX(ids,its-1), MIN(ite+1,ide+istag)
902               g_dat(i,k,jde+1) = g_dat(i,k,jde) 
903               dat(i,k,jde+1) = dat(i,k,jde) 
904               g_dat(i,k,jde+2) = g_dat(i,k,jde) 
905               dat(i,k,jde+2) = dat(i,k,jde) 
906               g_dat(i,k,jde+3) = g_dat(i,k,jde) 
907               dat(i,k,jde+3) = dat(i,k,jde) 
908             ENDDO                               
909             ENDDO
911           ENDIF
913       END IF open_ye
915 !  end open b.c in Y copy into boundary zone addition.  WCS, 19 March 2000
917       END IF periodicity_y
919 !  fix corners for doubly periodic domains
921       IF ( config_flags%periodic_x .and. config_flags%periodic_y &
922            .and. (ids == ips) .and. (ide == ipe)                 &
923            .and. (jds == jps) .and. (jde == jpe)                   ) THEN
925          IF ( (its == ids) .and. (jts == jds) ) THEN  ! lower left corner fill
926            DO j = 0, -(bdyzone-1), -1
927            DO k = kts, k_end
928            DO i = 0, -(bdyzone-1), -1
929              g_dat(ids+i-1,k,jds+j-1) = g_dat(ide+i-1,k,jde+j-1)
930              dat(ids+i-1,k,jds+j-1) = dat(ide+i-1,k,jde+j-1)
931            ENDDO
932            ENDDO
933            ENDDO
934          END IF
936          IF ( (ite == ide) .and. (jts == jds) ) THEN  ! lower right corner fill
937            DO j = 0, -(bdyzone-1), -1
938            DO k = kts, k_end
939            DO i = 1, bdyzone
940              g_dat(ide+i+istag,k,jds+j-1) = g_dat(ids+i+istag,k,jde+j-1)
941              dat(ide+i+istag,k,jds+j-1) = dat(ids+i+istag,k,jde+j-1)
942            ENDDO
943            ENDDO
944            ENDDO
945          END IF
947          IF ( (ite == ide) .and. (jte == jde) ) THEN  ! upper right corner fill
948            DO j = 1, bdyzone
949            DO k = kts, k_end
950            DO i = 1, bdyzone
951              g_dat(ide+i+istag,k,jde+j+jstag) = g_dat(ids+i+istag,k,jds+j+jstag)
952              dat(ide+i+istag,k,jde+j+jstag) = dat(ids+i+istag,k,jds+j+jstag)
953            ENDDO
954            ENDDO
955            ENDDO
956          END IF
958          IF ( (its == ids) .and. (jte == jde) ) THEN  ! upper left corner fill
959            DO j = 1, bdyzone
960            DO k = kts, k_end
961            DO i = 0, -(bdyzone-1), -1
962              g_dat(ids+i-1,k,jde+j+jstag) = g_dat(ide+i-1,k,jds+j+jstag)
963              dat(ids+i-1,k,jde+j+jstag) = dat(ide+i-1,k,jds+j+jstag)
964            ENDDO
965            ENDDO
966            ENDDO
967          END IF
969        END IF
971    END SUBROUTINE g_set_physical_bc3d
973 !------------------------------------------------------------------------
975    SUBROUTINE g_init_module_bc
976    END SUBROUTINE g_init_module_bc
978 !------------------------------------------------------------------------
980    SUBROUTINE g_relax_bdytend ( field, g_field, field_tend, g_field_tend, &
981                                 field_bdy_xs, g_field_bdy_xs,  &
982                                 field_bdy_xe, g_field_bdy_xe,  &
983                                 field_bdy_ys, g_field_bdy_ys,  &
984                                 field_bdy_ye, g_field_bdy_ye,  &
985                                 field_bdy_tend_xs, g_field_bdy_tend_xs,  &
986                                 field_bdy_tend_xe, g_field_bdy_tend_xe,  &
987                                 field_bdy_tend_ys, g_field_bdy_tend_ys,  &
988                                 field_bdy_tend_ye, g_field_bdy_tend_ye,  &
989                                 variable_in, config_flags,             &
990                                 spec_bdy_width, spec_zone, relax_zone, &
991                                 dtbc, fcx, gcx,             &
992                                 ids,ide, jds,jde, kds,kde,  & ! domain dims
993                                 ims,ime, jms,jme, kms,kme,  & ! memory dims
994                                 ips,ipe, jps,jpe, kps,kpe,  & ! patch  dims
995                                 its,ite, jts,jte, kts,kte  )
997    IMPLICIT NONE
999    INTEGER,   INTENT(IN) :: ids,ide, jds,jde, kds,kde
1000    INTEGER,   INTENT(IN) :: ims,ime, jms,jme, kms,kme
1001    INTEGER,   INTENT(IN) :: ips,ipe, jps,jpe, kps,kpe
1002    INTEGER,   INTENT(IN) :: its,ite, jts,jte, kts,kte
1003    INTEGER,   INTENT(IN) :: spec_bdy_width, spec_zone, relax_zone
1004    REAL,      INTENT(IN) :: dtbc
1005    CHARACTER, INTENT(IN) :: variable_in
1007    REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(IN) :: g_field
1008    REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(IN) :: field
1009    REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: g_field_tend
1010    REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: field_tend
1011    REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: g_field_bdy_xs, g_field_bdy_xe
1012    REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: field_bdy_xs, field_bdy_xe
1013    REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: g_field_bdy_ys, g_field_bdy_ye
1014    REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: field_bdy_ys, field_bdy_ye
1015    REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: g_field_bdy_tend_xs, g_field_bdy_tend_xe
1016    REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: field_bdy_tend_xs, field_bdy_tend_xe
1017    REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: g_field_bdy_tend_ys, g_field_bdy_tend_ye
1018    REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: field_bdy_tend_ys, field_bdy_tend_ye
1019    REAL, DIMENSION( spec_bdy_width ), INTENT(IN) :: fcx, gcx
1020    TYPE(grid_config_rec_type),        INTENT(IN) :: config_flags
1023    CALL g_relax_bdytend_core ( field, g_field, field_tend, g_field_tend, &
1024                        field_bdy_xs, g_field_bdy_xs,  &
1025                        field_bdy_xe, g_field_bdy_xe,  &
1026                        field_bdy_ys, g_field_bdy_ys,  &
1027                        field_bdy_ye, g_field_bdy_ye,  &
1028                        field_bdy_tend_xs, g_field_bdy_tend_xs,  &
1029                        field_bdy_tend_xe, g_field_bdy_tend_xe,  &
1030                        field_bdy_tend_ys, g_field_bdy_tend_ys,  &
1031                        field_bdy_tend_ye, g_field_bdy_tend_ye,  &
1032                        variable_in, config_flags,             &
1033                        spec_bdy_width, spec_zone, relax_zone, &
1034                        dtbc, fcx, gcx,             &
1035                        ids,ide, jds,jde, kds,kde,  & ! domain dims
1036                        ims,ime, jms,jme, kms,kme,  & ! memory dims
1037                        ips,ipe, jps,jpe, kps,kpe,  & ! patch  dims
1038                        its,ite, jts,jte, kts,kte,  & ! patch  dims
1039                        ims,ime, jms,jme, kms,kme )  ! dimension of the field argument
1041    END SUBROUTINE g_relax_bdytend
1043 ! version that allows tile-sized version of field. Note, caller should define the
1044 ! field to be -+1 of tile size in each dimension because routine is going off onto halo
1045 ! for example, see relax_bdytend in dyn_em/module_bc_em.F 
1046    SUBROUTINE g_relax_bdytend_tile ( field, g_field, field_tend, g_field_tend, & 
1047                           field_bdy_xs, g_field_bdy_xs,  &
1048                           field_bdy_xe, g_field_bdy_xe,  &
1049                           field_bdy_ys, g_field_bdy_ys,  &
1050                           field_bdy_ye, g_field_bdy_ye,  &
1051                           field_bdy_tend_xs, g_field_bdy_tend_xs,  &
1052                           field_bdy_tend_xe, g_field_bdy_tend_xe,  &
1053                           field_bdy_tend_ys, g_field_bdy_tend_ys,  &
1054                           field_bdy_tend_ye, g_field_bdy_tend_ye,  &
1055                           variable_in, config_flags,             &
1056                           spec_bdy_width, spec_zone, relax_zone, &
1057                           dtbc, fcx, gcx,             &
1058                           ids,ide, jds,jde, kds,kde,  & ! domain dims
1059                           ims,ime, jms,jme, kms,kme,  & ! memory dims
1060                           ips,ipe, jps,jpe, kps,kpe,  & ! patch  dims
1061                           its,ite, jts,jte, kts,kte,  &
1062                           iXs,iXe, jXs,jXe, kXs,kXe   &  ! dims of first argument
1063                           )
1065    IMPLICIT NONE
1067    INTEGER,   INTENT(IN) :: ids,ide, jds,jde, kds,kde
1068    INTEGER,   INTENT(IN) :: ims,ime, jms,jme, kms,kme
1069    INTEGER,   INTENT(IN) :: ips,ipe, jps,jpe, kps,kpe
1070    INTEGER,   INTENT(IN) :: its,ite, jts,jte, kts,kte
1071    INTEGER,   INTENT(IN) :: iXs,iXe, jXs,jXe, kXs,kXe
1072    INTEGER,   INTENT(IN) :: spec_bdy_width, spec_zone, relax_zone
1073    REAL,      INTENT(IN) :: dtbc
1074    CHARACTER, INTENT(IN) :: variable_in
1076    REAL, DIMENSION(iXs:iXe, kXs:kXe, jXs:jXe ), INTENT(IN   ) :: g_field
1077    REAL, DIMENSION(iXs:iXe, kXs:kXe, jXs:jXe ), INTENT(IN   ) :: field
1078    REAL, DIMENSION(ims:ime, kms:kme, jms:jme ), INTENT(INOUT) :: g_field_tend
1079    REAL, DIMENSION(ims:ime, kms:kme, jms:jme ), INTENT(INOUT) :: field_tend
1080    REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: g_field_bdy_xs, g_field_bdy_xe
1081    REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: field_bdy_xs, field_bdy_xe
1082    REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: g_field_bdy_ys, g_field_bdy_ye
1083    REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: field_bdy_ys, field_bdy_ye
1084    REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: g_field_bdy_tend_xs, g_field_bdy_tend_xe
1085    REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: field_bdy_tend_xs, field_bdy_tend_xe
1086    REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: g_field_bdy_tend_ys, g_field_bdy_tend_ye
1087    REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: field_bdy_tend_ys, field_bdy_tend_ye
1088    REAL, DIMENSION(spec_bdy_width), INTENT(IN) :: fcx, gcx
1089    TYPE(grid_config_rec_type),      INTENT(IN) :: config_flags
1092    CALL g_relax_bdytend_core ( field, g_field, field_tend, g_field_tend, &
1093                        field_bdy_xs, g_field_bdy_xs,  &
1094                        field_bdy_xe, g_field_bdy_xe,  &
1095                        field_bdy_ys, g_field_bdy_ys,  &
1096                        field_bdy_ye, g_field_bdy_ye,  &
1097                        field_bdy_tend_xs, g_field_bdy_tend_xs,  &
1098                        field_bdy_tend_xe, g_field_bdy_tend_xe,  &
1099                        field_bdy_tend_ys, g_field_bdy_tend_ys,  &
1100                        field_bdy_tend_ye, g_field_bdy_tend_ye,  &
1101                        variable_in, config_flags,             &
1102                        spec_bdy_width, spec_zone, relax_zone, &
1103                        dtbc, fcx, gcx,             &
1104                        ids,ide, jds,jde, kds,kde,  & ! domain dims
1105                        ims,ime, jms,jme, kms,kme,  & ! memory dims
1106                        ips,ipe, jps,jpe, kps,kpe,  & ! patch  dims
1107                        its,ite, jts,jte, kts,kte,  &
1108                        iXs,iXe, jXs,jXe, kXs,kXe )  ! dimension of the field argument
1110    END SUBROUTINE g_relax_bdytend_tile
1112 !  Tapenade 3.5 (r3785) - 22 Mar 2011 18:35
1114 !  Differentiation of relax_bdytend_core in forward (tangent) mode:
1115 !   variations   of useful results: field_tend
1116 !   with respect to varying inputs: field field_bdy_xe field_bdy_tend_xe
1117 !                field_tend field_bdy_xs field_bdy_tend_xs field_bdy_ye
1118 !                field_bdy_tend_ye field_bdy_ys field_bdy_tend_ys
1119 !   RW status of diff variables: field:in field_bdy_xe:in field_bdy_tend_xe:in
1120 !                field_tend:in-out field_bdy_xs:in field_bdy_tend_xs:in
1121 !                field_bdy_ye:in field_bdy_tend_ye:in field_bdy_ys:in
1122 !                field_bdy_tend_ys:in
1123 ! domain dims
1124 ! memory dims
1125 ! patch  dims
1126 ! patch  dims
1127 ! field (1st arg) dims; might be tile or patch
1128 SUBROUTINE G_RELAX_BDYTEND_CORE(field, fieldd, field_tend, field_tendd, &
1129 &  field_bdy_xs, field_bdy_xsd, field_bdy_xe, field_bdy_xed, field_bdy_ys&
1130 &  , field_bdy_ysd, field_bdy_ye, field_bdy_yed, field_bdy_tend_xs, &
1131 &  field_bdy_tend_xsd, field_bdy_tend_xe, field_bdy_tend_xed, &
1132 &  field_bdy_tend_ys, field_bdy_tend_ysd, field_bdy_tend_ye, &
1133 &  field_bdy_tend_yed, variable_in, config_flags, spec_bdy_width, &
1134 &  spec_zone, relax_zone, dtbc, fcx, gcx, ids, ide, jds, jde, kds, kde, &
1135 &  ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe, its, ite, &
1136 &  jts, jte, kts, kte, ixs, ixe, jxs, jxe, kxs, kxe)
1137   IMPLICIT NONE
1138   INTEGER, INTENT(IN) :: ids, ide, jds, jde, kds, kde
1139   INTEGER, INTENT(IN) :: ims, ime, jms, jme, kms, kme
1140   INTEGER, INTENT(IN) :: ips, ipe, jps, jpe, kps, kpe
1141   INTEGER, INTENT(IN) :: its, ite, jts, jte, kts, kte
1142   INTEGER, INTENT(IN) :: ixs, ixe, jxs, jxe, kxs, kxe
1143   INTEGER, INTENT(IN) :: spec_bdy_width, spec_zone, relax_zone
1144   REAL, INTENT(IN) :: dtbc
1145   CHARACTER, INTENT(IN) :: variable_in
1146   REAL, DIMENSION(ixs:ixe, kxs:kxe, jxs:jxe), INTENT(IN) :: field
1147   REAL, DIMENSION(ixs:ixe, kxs:kxe, jxs:jxe), INTENT(IN) :: fieldd
1148   REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: &
1149 &  field_tend
1150   REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: &
1151 &  field_tendd
1152   REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: &
1153 &  field_bdy_xs, field_bdy_xe
1154   REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: &
1155 &  field_bdy_xsd, field_bdy_xed
1156   REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: &
1157 &  field_bdy_ys, field_bdy_ye
1158   REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: &
1159 &  field_bdy_ysd, field_bdy_yed
1160   REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: &
1161 &  field_bdy_tend_xs, field_bdy_tend_xe
1162   REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: &
1163 &  field_bdy_tend_xsd, field_bdy_tend_xed
1164   REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: &
1165 &  field_bdy_tend_ys, field_bdy_tend_ye
1166   REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: &
1167 &  field_bdy_tend_ysd, field_bdy_tend_yed
1168   REAL, DIMENSION(spec_bdy_width), INTENT(IN) :: fcx, gcx
1169   TYPE(GRID_CONFIG_REC_TYPE) :: config_flags
1170   CHARACTER :: variable
1171   INTEGER :: i, j, k, ibs, ibe, jbs, jbe, itf, jtf, ktf, im1, ip1
1172   INTEGER :: b_dist, b_limit
1173   REAL :: fls0, fls1, fls2, fls3, fls4
1174   REAL :: fls0d, fls1d, fls2d, fls3d, fls4d
1175   LOGICAL :: periodic_x
1176   INTEGER :: min8
1177   INTEGER :: min7
1178   INTEGER :: min6
1179   INTEGER :: min5
1180   INTEGER :: min4
1181   INTEGER :: min3
1182   INTEGER :: min2
1183   INTEGER :: min1
1184   INTEGER :: max8
1185   INTEGER :: max7
1186   INTEGER :: max6
1187   INTEGER :: max5
1188   INTEGER :: max4
1189   INTEGER :: max3
1190   INTEGER :: max2
1191   INTEGER :: max1
1192   periodic_x = config_flags%periodic_x
1193   variable = variable_in
1194   IF (variable .EQ. 'U') variable = 'u'
1195   IF (variable .EQ. 'V') variable = 'v'
1196   IF (variable .EQ. 'M') variable = 'm'
1197   IF (variable .EQ. 'H') variable = 'h'
1198   ibs = ids
1199   ibe = ide - 1
1200   IF (ite .GT. ide - 1) THEN
1201     itf = ide - 1
1202   ELSE
1203     itf = ite
1204   END IF
1205   jbs = jds
1206   jbe = jde - 1
1207   IF (jte .GT. jde - 1) THEN
1208     jtf = jde - 1
1209   ELSE
1210     jtf = jte
1211   END IF
1212   ktf = kde - 1
1213   IF (variable .EQ. 'u') ibe = ide
1214   IF (variable .EQ. 'u') THEN
1215     IF (ite .GT. ide) THEN
1216       itf = ide
1217     ELSE
1218       itf = ite
1219     END IF
1220   END IF
1221   IF (variable .EQ. 'v') jbe = jde
1222   IF (variable .EQ. 'v') THEN
1223     IF (jte .GT. jde) THEN
1224       jtf = jde
1225     ELSE
1226       jtf = jte
1227     END IF
1228   END IF
1229   IF (variable .EQ. 'm') ktf = kte
1230   IF (variable .EQ. 'h') ktf = kte
1231   IF (jts - jbs .LT. relax_zone) THEN
1232     IF (jts .LT. jbs + spec_zone) THEN
1233       max1 = jbs + spec_zone
1234     ELSE
1235       max1 = jts
1236     END IF
1237     IF (jtf .GT. jbs + relax_zone - 1) THEN
1238       min1 = jbs + relax_zone - 1
1239     ELSE
1240       min1 = jtf
1241     END IF
1242 ! Y-start boundary
1243     DO j=max1,min1
1244       b_dist = j - jbs
1245       b_limit = b_dist
1246       IF (periodic_x) b_limit = 0
1247       DO k=kts,ktf
1248         IF (its .LT. b_limit + ibs) THEN
1249           max2 = b_limit + ibs
1250         ELSE
1251           max2 = its
1252         END IF
1253         IF (itf .GT. ibe - b_limit) THEN
1254           min2 = ibe - b_limit
1255         ELSE
1256           min2 = itf
1257         END IF
1258         DO i=max2,min2
1259           IF (i - 1 .LT. ibs) THEN
1260             im1 = ibs
1261           ELSE
1262             im1 = i - 1
1263           END IF
1264           IF (i + 1 .GT. ibe) THEN
1265             ip1 = ibe
1266           ELSE
1267             ip1 = i + 1
1268           END IF
1269           fls0d = field_bdy_ysd(i, k, b_dist+1) + dtbc*&
1270 &            field_bdy_tend_ysd(i, k, b_dist+1) - fieldd(i, k, j)
1271           fls0 = field_bdy_ys(i, k, b_dist+1) + dtbc*field_bdy_tend_ys(i&
1272 &            , k, b_dist+1) - field(i, k, j)
1273           fls1d = field_bdy_ysd(im1, k, b_dist+1) + dtbc*&
1274 &            field_bdy_tend_ysd(im1, k, b_dist+1) - fieldd(im1, k, j)
1275           fls1 = field_bdy_ys(im1, k, b_dist+1) + dtbc*field_bdy_tend_ys&
1276 &            (im1, k, b_dist+1) - field(im1, k, j)
1277           fls2d = field_bdy_ysd(ip1, k, b_dist+1) + dtbc*&
1278 &            field_bdy_tend_ysd(ip1, k, b_dist+1) - fieldd(ip1, k, j)
1279           fls2 = field_bdy_ys(ip1, k, b_dist+1) + dtbc*field_bdy_tend_ys&
1280 &            (ip1, k, b_dist+1) - field(ip1, k, j)
1281           fls3d = field_bdy_ysd(i, k, b_dist) + dtbc*field_bdy_tend_ysd(&
1282 &            i, k, b_dist) - fieldd(i, k, j-1)
1283           fls3 = field_bdy_ys(i, k, b_dist) + dtbc*field_bdy_tend_ys(i, &
1284 &            k, b_dist) - field(i, k, j-1)
1285           fls4d = field_bdy_ysd(i, k, b_dist+2) + dtbc*&
1286 &            field_bdy_tend_ysd(i, k, b_dist+2) - fieldd(i, k, j+1)
1287           fls4 = field_bdy_ys(i, k, b_dist+2) + dtbc*field_bdy_tend_ys(i&
1288 &            , k, b_dist+2) - field(i, k, j+1)
1289           field_tendd(i, k, j) = field_tendd(i, k, j) + fcx(b_dist+1)*&
1290 &            fls0d - gcx(b_dist+1)*(fls1d+fls2d+fls3d+fls4d-4.*fls0d)
1291           field_tend(i, k, j) = field_tend(i, k, j) + fcx(b_dist+1)*fls0&
1292 &            - gcx(b_dist+1)*(fls1+fls2+fls3+fls4-4.*fls0)
1293         END DO
1294       END DO
1295     END DO
1296   END IF
1297   IF (jbe - jtf .LT. relax_zone) THEN
1298     IF (jts .LT. jbe - relax_zone + 1) THEN
1299       max3 = jbe - relax_zone + 1
1300     ELSE
1301       max3 = jts
1302     END IF
1303     IF (jtf .GT. jbe - spec_zone) THEN
1304       min3 = jbe - spec_zone
1305     ELSE
1306       min3 = jtf
1307     END IF
1308 ! Y-end boundary
1309     DO j=max3,min3
1310       b_dist = jbe - j
1311       b_limit = b_dist
1312       IF (periodic_x) b_limit = 0
1313       DO k=kts,ktf
1314         IF (its .LT. b_limit + ibs) THEN
1315           max4 = b_limit + ibs
1316         ELSE
1317           max4 = its
1318         END IF
1319         IF (itf .GT. ibe - b_limit) THEN
1320           min4 = ibe - b_limit
1321         ELSE
1322           min4 = itf
1323         END IF
1324         DO i=max4,min4
1325           IF (i - 1 .LT. ibs) THEN
1326             im1 = ibs
1327           ELSE
1328             im1 = i - 1
1329           END IF
1330           IF (i + 1 .GT. ibe) THEN
1331             ip1 = ibe
1332           ELSE
1333             ip1 = i + 1
1334           END IF
1335           fls0d = field_bdy_yed(i, k, b_dist+1) + dtbc*&
1336 &            field_bdy_tend_yed(i, k, b_dist+1) - fieldd(i, k, j)
1337           fls0 = field_bdy_ye(i, k, b_dist+1) + dtbc*field_bdy_tend_ye(i&
1338 &            , k, b_dist+1) - field(i, k, j)
1339           fls1d = field_bdy_yed(im1, k, b_dist+1) + dtbc*&
1340 &            field_bdy_tend_yed(im1, k, b_dist+1) - fieldd(im1, k, j)
1341           fls1 = field_bdy_ye(im1, k, b_dist+1) + dtbc*field_bdy_tend_ye&
1342 &            (im1, k, b_dist+1) - field(im1, k, j)
1343           fls2d = field_bdy_yed(ip1, k, b_dist+1) + dtbc*&
1344 &            field_bdy_tend_yed(ip1, k, b_dist+1) - fieldd(ip1, k, j)
1345           fls2 = field_bdy_ye(ip1, k, b_dist+1) + dtbc*field_bdy_tend_ye&
1346 &            (ip1, k, b_dist+1) - field(ip1, k, j)
1347           fls3d = field_bdy_yed(i, k, b_dist) + dtbc*field_bdy_tend_yed(&
1348 &            i, k, b_dist) - fieldd(i, k, j+1)
1349           fls3 = field_bdy_ye(i, k, b_dist) + dtbc*field_bdy_tend_ye(i, &
1350 &            k, b_dist) - field(i, k, j+1)
1351           fls4d = field_bdy_yed(i, k, b_dist+2) + dtbc*&
1352 &            field_bdy_tend_yed(i, k, b_dist+2) - fieldd(i, k, j-1)
1353           fls4 = field_bdy_ye(i, k, b_dist+2) + dtbc*field_bdy_tend_ye(i&
1354 &            , k, b_dist+2) - field(i, k, j-1)
1355           field_tendd(i, k, j) = field_tendd(i, k, j) + fcx(b_dist+1)*&
1356 &            fls0d - gcx(b_dist+1)*(fls1d+fls2d+fls3d+fls4d-4.*fls0d)
1357           field_tend(i, k, j) = field_tend(i, k, j) + fcx(b_dist+1)*fls0&
1358 &            - gcx(b_dist+1)*(fls1+fls2+fls3+fls4-4.*fls0)
1359         END DO
1360       END DO
1361     END DO
1362   END IF
1363   IF (.NOT.periodic_x) THEN
1364     IF (its - ibs .LT. relax_zone) THEN
1365       IF (its .LT. ibs + spec_zone) THEN
1366         max5 = ibs + spec_zone
1367       ELSE
1368         max5 = its
1369       END IF
1370       IF (itf .GT. ibs + relax_zone - 1) THEN
1371         min5 = ibs + relax_zone - 1
1372       ELSE
1373         min5 = itf
1374       END IF
1375 ! X-start boundary
1376       DO i=max5,min5
1377         b_dist = i - ibs
1378         DO k=kts,ktf
1379           IF (jts .LT. b_dist + jbs + 1) THEN
1380             max6 = b_dist + jbs + 1
1381           ELSE
1382             max6 = jts
1383           END IF
1384           IF (jtf .GT. jbe - b_dist - 1) THEN
1385             min6 = jbe - b_dist - 1
1386           ELSE
1387             min6 = jtf
1388           END IF
1389           DO j=max6,min6
1390             fls0d = field_bdy_xsd(j, k, b_dist+1) + dtbc*&
1391 &              field_bdy_tend_xsd(j, k, b_dist+1) - fieldd(i, k, j)
1392             fls0 = field_bdy_xs(j, k, b_dist+1) + dtbc*field_bdy_tend_xs&
1393 &              (j, k, b_dist+1) - field(i, k, j)
1394             fls1d = field_bdy_xsd(j-1, k, b_dist+1) + dtbc*&
1395 &              field_bdy_tend_xsd(j-1, k, b_dist+1) - fieldd(i, k, j-1)
1396             fls1 = field_bdy_xs(j-1, k, b_dist+1) + dtbc*&
1397 &              field_bdy_tend_xs(j-1, k, b_dist+1) - field(i, k, j-1)
1398             fls2d = field_bdy_xsd(j+1, k, b_dist+1) + dtbc*&
1399 &              field_bdy_tend_xsd(j+1, k, b_dist+1) - fieldd(i, k, j+1)
1400             fls2 = field_bdy_xs(j+1, k, b_dist+1) + dtbc*&
1401 &              field_bdy_tend_xs(j+1, k, b_dist+1) - field(i, k, j+1)
1402             fls3d = field_bdy_xsd(j, k, b_dist) + dtbc*&
1403 &              field_bdy_tend_xsd(j, k, b_dist) - fieldd(i-1, k, j)
1404             fls3 = field_bdy_xs(j, k, b_dist) + dtbc*field_bdy_tend_xs(j&
1405 &              , k, b_dist) - field(i-1, k, j)
1406             fls4d = field_bdy_xsd(j, k, b_dist+2) + dtbc*&
1407 &              field_bdy_tend_xsd(j, k, b_dist+2) - fieldd(i+1, k, j)
1408             fls4 = field_bdy_xs(j, k, b_dist+2) + dtbc*field_bdy_tend_xs&
1409 &              (j, k, b_dist+2) - field(i+1, k, j)
1410             field_tendd(i, k, j) = field_tendd(i, k, j) + fcx(b_dist+1)*&
1411 &              fls0d - gcx(b_dist+1)*(fls1d+fls2d+fls3d+fls4d-4.*fls0d)
1412             field_tend(i, k, j) = field_tend(i, k, j) + fcx(b_dist+1)*&
1413 &              fls0 - gcx(b_dist+1)*(fls1+fls2+fls3+fls4-4.*fls0)
1414           END DO
1415         END DO
1416       END DO
1417     END IF
1418     IF (ibe - itf .LT. relax_zone) THEN
1419       IF (its .LT. ibe - relax_zone + 1) THEN
1420         max7 = ibe - relax_zone + 1
1421       ELSE
1422         max7 = its
1423       END IF
1424       IF (itf .GT. ibe - spec_zone) THEN
1425         min7 = ibe - spec_zone
1426       ELSE
1427         min7 = itf
1428       END IF
1429 ! X-end boundary
1430       DO i=max7,min7
1431         b_dist = ibe - i
1432         DO k=kts,ktf
1433           IF (jts .LT. b_dist + jbs + 1) THEN
1434             max8 = b_dist + jbs + 1
1435           ELSE
1436             max8 = jts
1437           END IF
1438           IF (jtf .GT. jbe - b_dist - 1) THEN
1439             min8 = jbe - b_dist - 1
1440           ELSE
1441             min8 = jtf
1442           END IF
1443           DO j=max8,min8
1444             fls0d = field_bdy_xed(j, k, b_dist+1) + dtbc*&
1445 &              field_bdy_tend_xed(j, k, b_dist+1) - fieldd(i, k, j)
1446             fls0 = field_bdy_xe(j, k, b_dist+1) + dtbc*field_bdy_tend_xe&
1447 &              (j, k, b_dist+1) - field(i, k, j)
1448             fls1d = field_bdy_xed(j-1, k, b_dist+1) + dtbc*&
1449 &              field_bdy_tend_xed(j-1, k, b_dist+1) - fieldd(i, k, j-1)
1450             fls1 = field_bdy_xe(j-1, k, b_dist+1) + dtbc*&
1451 &              field_bdy_tend_xe(j-1, k, b_dist+1) - field(i, k, j-1)
1452             fls2d = field_bdy_xed(j+1, k, b_dist+1) + dtbc*&
1453 &              field_bdy_tend_xed(j+1, k, b_dist+1) - fieldd(i, k, j+1)
1454             fls2 = field_bdy_xe(j+1, k, b_dist+1) + dtbc*&
1455 &              field_bdy_tend_xe(j+1, k, b_dist+1) - field(i, k, j+1)
1456             fls3d = field_bdy_xed(j, k, b_dist) + dtbc*&
1457 &              field_bdy_tend_xed(j, k, b_dist) - fieldd(i+1, k, j)
1458             fls3 = field_bdy_xe(j, k, b_dist) + dtbc*field_bdy_tend_xe(j&
1459 &              , k, b_dist) - field(i+1, k, j)
1460             fls4d = field_bdy_xed(j, k, b_dist+2) + dtbc*&
1461 &              field_bdy_tend_xed(j, k, b_dist+2) - fieldd(i-1, k, j)
1462             fls4 = field_bdy_xe(j, k, b_dist+2) + dtbc*field_bdy_tend_xe&
1463 &              (j, k, b_dist+2) - field(i-1, k, j)
1464             field_tendd(i, k, j) = field_tendd(i, k, j) + fcx(b_dist+1)*&
1465 &              fls0d - gcx(b_dist+1)*(fls1d+fls2d+fls3d+fls4d-4.*fls0d)
1466             field_tend(i, k, j) = field_tend(i, k, j) + fcx(b_dist+1)*&
1467 &              fls0 - gcx(b_dist+1)*(fls1+fls2+fls3+fls4-4.*fls0)
1468           END DO
1469         END DO
1470       END DO
1471     END IF
1472   END IF
1473 END SUBROUTINE G_RELAX_BDYTEND_CORE
1474 !------------------------------------------------------------------------
1476    SUBROUTINE g_spec_bdytend ( field_tend, g_field_tend,     &
1477                                field_bdy_xs, g_field_bdy_xs, &
1478                                field_bdy_xe, g_field_bdy_xe, &
1479                                field_bdy_ys, g_field_bdy_ys, &
1480                                field_bdy_ye, g_field_bdy_ye, &
1481                                field_bdy_tend_xs, g_field_bdy_tend_xs, &
1482                                field_bdy_tend_xe, g_field_bdy_tend_xe, &
1483                                field_bdy_tend_ys, g_field_bdy_tend_ys, &
1484                                field_bdy_tend_ye, g_field_bdy_tend_ye, &
1485                                variable_in, config_flags, & 
1486                                spec_bdy_width, spec_zone, &
1487                                ids,ide, jds,jde, kds,kde,  & ! domain dims
1488                                ims,ime, jms,jme, kms,kme,  & ! memory dims
1489                                ips,ipe, jps,jpe, kps,kpe,  & ! patch  dims
1490                                its,ite, jts,jte, kts,kte )
1492 !  spec_bdy_width is only used to dimension the boundary arrays.
1493 !  spec_zone is the width of the outer specified b.c.s that are set here.
1495       IMPLICIT NONE
1497       INTEGER,   INTENT(IN) :: ids,ide, jds,jde, kds,kde
1498       INTEGER,   INTENT(IN) :: ims,ime, jms,jme, kms,kme
1499       INTEGER,   INTENT(IN) :: ips,ipe, jps,jpe, kps,kpe
1500       INTEGER,   INTENT(IN) :: its,ite, jts,jte, kts,kte
1501       INTEGER,   INTENT(IN) :: spec_bdy_width, spec_zone
1502       CHARACTER, INTENT(IN) :: variable_in
1504       REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(OUT) :: g_field_tend
1505       REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(OUT) :: field_tend
1507       REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: g_field_bdy_xs, g_field_bdy_xe
1508       REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: g_field_bdy_ys, g_field_bdy_ye 
1509       REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: g_field_bdy_tend_xs, g_field_bdy_tend_xe
1510       REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: g_field_bdy_tend_ys, g_field_bdy_tend_ye 
1511       REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: field_bdy_xs, field_bdy_xe
1512       REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: field_bdy_ys, field_bdy_ye 
1513       REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: field_bdy_tend_xs, field_bdy_tend_xe
1514       REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: field_bdy_tend_ys, field_bdy_tend_ye 
1515       TYPE(grid_config_rec_type), INTENT(IN) :: config_flags
1517       CHARACTER  :: variable
1518       INTEGER    :: i, j, k, ibs, ibe, jbs, jbe, itf, jtf, ktf
1519       INTEGER    :: b_dist, b_limit
1520       LOGICAL    :: periodic_x
1523       periodic_x = config_flags%periodic_x
1525       variable = variable_in
1527       IF (variable == 'U') variable = 'u'
1528       IF (variable == 'V') variable = 'v'
1529       IF (variable == 'M') variable = 'm'
1530       IF (variable == 'H') variable = 'h'
1532       ibs = ids
1533       ibe = ide-1
1534       itf = min(ite,ide-1)
1535       jbs = jds
1536       jbe = jde-1
1537       jtf = min(jte,jde-1)
1538       ktf = kde-1
1539       IF (variable == 'u') ibe = ide
1540       IF (variable == 'u') itf = min(ite,ide)
1541       IF (variable == 'v') jbe = jde
1542       IF (variable == 'v') jtf = min(jte,jde)
1543       IF (variable == 'm') ktf = kte
1544       IF (variable == 'h') ktf = kte
1546       IF (jts - jbs .lt. spec_zone) THEN
1547 ! Y-start boundary
1548         DO j = jts, min(jtf,jbs+spec_zone-1)
1549           b_dist = j - jbs
1550           b_limit = b_dist
1551           IF(periodic_x)b_limit = 0
1552           DO k = kts, ktf
1553             DO i = max(its,b_limit+ibs), min(itf,ibe-b_limit)
1554               g_field_tend(i,k,j) = g_field_bdy_tend_ys(i, k, b_dist+1)
1555               field_tend(i,k,j) = field_bdy_tend_ys(i, k, b_dist+1)
1556             ENDDO
1557           ENDDO
1558         ENDDO
1559       ENDIF 
1560       IF (jbe - jtf .lt. spec_zone) THEN 
1561 ! Y-end boundary 
1562         DO j = max(jts,jbe-spec_zone+1), jtf 
1563           b_dist = jbe - j 
1564           b_limit = b_dist
1565           IF(periodic_x)b_limit = 0
1566           DO k = kts, ktf 
1567             DO i = max(its,b_limit+ibs), min(itf,ibe-b_limit)
1568               g_field_tend(i,k,j) = g_field_bdy_tend_ye(i, k, b_dist+1)
1569               field_tend(i,k,j) = field_bdy_tend_ye(i, k, b_dist+1)
1570             ENDDO
1571           ENDDO
1572         ENDDO
1573       ENDIF 
1575     IF(.NOT.periodic_x)THEN
1576       IF (its - ibs .lt. spec_zone) THEN
1577 ! X-start boundary
1578         DO i = its, min(itf,ibs+spec_zone-1)
1579           b_dist = i - ibs
1580           DO k = kts, ktf
1581             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
1582               g_field_tend(i,k,j) = g_field_bdy_tend_xs(j, k, b_dist+1)
1583               field_tend(i,k,j) = field_bdy_tend_xs(j, k, b_dist+1)
1584             ENDDO
1585           ENDDO
1586         ENDDO
1587       ENDIF 
1589       IF (ibe - itf .lt. spec_zone) THEN
1590 ! X-end boundary
1591         DO i = max(its,ibe-spec_zone+1), itf
1592           b_dist = ibe - i
1593           DO k = kts, ktf
1594             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
1595               g_field_tend(i,k,j) = g_field_bdy_tend_xe(j, k, b_dist+1)
1596               field_tend(i,k,j) = field_bdy_tend_xe(j, k, b_dist+1)
1597             ENDDO
1598           ENDDO
1599         ENDDO
1600       ENDIF 
1601     ENDIF
1603    END SUBROUTINE g_spec_bdytend
1605 !------------------------------------------------------------------------
1607    SUBROUTINE g_couple_bdy ( field, g_field,  &
1608                              variable_in, config_flags, & 
1609                              spec_zone,       &
1610                              mu, g_mu, msf,   &
1611                              ids,ide, jds,jde, kds,kde,  & ! domain dims
1612                              ims,ime, jms,jme, kms,kme,  & ! memory dims
1613                              ips,ipe, jps,jpe, kps,kpe,  & ! patch  dims
1614                              its,ite, jts,jte, kts,kte )
1616 !  This subroutine adds the tendencies in the boundary specified region.
1617 !  spec_zone is the width of the outer specified b.c.s that are set here.
1618 !  (JD August 2000)
1620       IMPLICIT NONE
1622       INTEGER,      INTENT(IN   )    :: ids,ide, jds,jde, kds,kde
1623       INTEGER,      INTENT(IN   )    :: ims,ime, jms,jme, kms,kme
1624       INTEGER,      INTENT(IN   )    :: ips,ipe, jps,jpe, kps,kpe
1625       INTEGER,      INTENT(IN   )    :: its,ite, jts,jte, kts,kte
1626       INTEGER,      INTENT(IN   )    :: spec_zone
1627       CHARACTER,    INTENT(IN   )    :: variable_in
1628       REAL, DIMENSION( ims:ime , jms:jme ),    INTENT(IN   ) :: g_mu
1629       REAL, DIMENSION( ims:ime , jms:jme ),    INTENT(IN   ) :: mu
1630       REAL, DIMENSION( ims:ime , jms:jme ),    INTENT(IN   ) :: msf
1631       REAL, DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(INOUT) :: g_field
1632       REAL, DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(INOUT) :: field
1633       TYPE( grid_config_rec_type ) config_flags
1635       CHARACTER  :: variable
1636       INTEGER    :: i, j, k, ibs, ibe, jbs, jbe, itf, jtf, ktf
1637       INTEGER    :: b_dist, b_limit
1638       LOGICAL    :: periodic_x
1640       periodic_x = config_flags%periodic_x
1642       variable = variable_in
1644       IF (variable == 'U') variable = 'u'
1645       IF (variable == 'V') variable = 'v'
1646       IF (variable == 'T') variable = 't'
1647       IF (variable == 'H') variable = 'h'
1648       IF (variable == 'W') variable = 'w'
1650       ibs = ids
1651       ibe = ide-1
1652       itf = min(ite,ide-1)
1653       jbs = jds
1654       jbe = jde-1
1655       jtf = min(jte,jde-1)
1656       ktf = kde-1
1657       IF (variable == 'u') ibe = ide
1658       IF (variable == 'u') itf = min(ite,ide)
1659       IF (variable == 'v') jbe = jde
1660       IF (variable == 'v') jtf = min(jte,jde)
1661       IF (variable == 'h') ktf = kte
1662       IF (variable == 'w') ktf = kte
1664       IF (jts - jbs .lt. spec_zone) THEN
1665 ! Y-start boundary
1666         DO j = jts, min(jtf,jbs+spec_zone-1)
1667           b_dist = j - jbs
1668           b_limit = b_dist
1669           IF(periodic_x)b_limit = 0
1670           DO k = kts, ktf
1671             DO i = max(its,b_limit+ibs), min(itf,ibe-b_limit)
1672               if (variable == 't' .or. variable == 'h') then 
1673                  g_field(i,k,j) = g_field(i,k,j)*mu(i,j) + field(i,k,j)*g_mu(i,j)
1674                  field(i,k,j) = field(i,k,j)*mu(i,j)
1675               else 
1676                  g_field(i,k,j) = (g_field(i,k,j)*mu(i,j)+field(i,k,j)*g_mu(i,j)) / msf(i,j)
1677                  field(i,k,j) = field(i,k,j)*mu(i,j)/msf(i,j)
1678               end if
1679             ENDDO
1680           ENDDO
1681         ENDDO
1682       ENDIF 
1683       IF (jbe - jtf .lt. spec_zone) THEN 
1684 ! Y-end boundary 
1685         DO j = max(jts,jbe-spec_zone+1), jtf 
1686           b_dist = jbe - j 
1687           b_limit = b_dist
1688           IF(periodic_x)b_limit = 0
1689           DO k = kts, ktf 
1690             DO i = max(its,b_limit+ibs), min(itf,ibe-b_limit)
1691               if (variable == 't' .or. variable == 'h') then 
1692                  g_field(i,k,j) = g_field(i,k,j)*mu(i,j) + field(i,k,j)*g_mu(i,j)
1693                  field(i,k,j) = field(i,k,j)*mu(i,j)
1694               else 
1695                  g_field(i,k,j) = (g_field(i,k,j)*mu(i,j)+field(i,k,j)*g_mu(i,j)) / msf(i,j)
1696                  field(i,k,j) = field(i,k,j)*mu(i,j)/msf(i,j)
1697               end if
1698             ENDDO
1699           ENDDO
1700         ENDDO
1701       ENDIF 
1703     IF(.NOT.periodic_x)THEN
1704       IF (its - ibs .lt. spec_zone) THEN
1705 ! X-start boundary
1706         DO i = its, min(itf,ibs+spec_zone-1)
1707           b_dist = i - ibs
1708           DO k = kts, ktf
1709             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
1710               if (variable == 't' .or. variable == 'h') then 
1711                  g_field(i,k,j) = g_field(i,k,j)*mu(i,j) + field(i,k,j)*g_mu(i,j)
1712                  field(i,k,j) = field(i,k,j)*mu(i,j)
1713               else 
1714                  g_field(i,k,j) = (g_field(i,k,j)*mu(i,j)+field(i,k,j)*g_mu(i,j)) / msf(i,j)
1715                  field(i,k,j) = field(i,k,j)*mu(i,j)/msf(i,j)
1716               end if
1717             ENDDO
1718           ENDDO
1719         ENDDO
1720       ENDIF 
1722       IF (ibe - itf .lt. spec_zone) THEN
1723 ! X-end boundary
1724         DO i = max(its,ibe-spec_zone+1), itf
1725           b_dist = ibe - i
1726           DO k = kts, ktf
1727             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
1728               if (variable == 't' .or. variable == 'h') then 
1729                  g_field(i,k,j) = g_field(i,k,j)*mu(i,j) + field(i,k,j)*g_mu(i,j)
1730                  field(i,k,j) = field(i,k,j)*mu(i,j)
1731               else 
1732                  g_field(i,k,j) = (g_field(i,k,j)*mu(i,j)+field(i,k,j)*g_mu(i,j)) / msf(i,j)
1733                  field(i,k,j) = field(i,k,j)*mu(i,j)/msf(i,j)
1734               end if
1735             ENDDO
1736           ENDDO
1737         ENDDO
1738       ENDIF 
1739     ENDIF
1741    END SUBROUTINE g_couple_bdy
1743 !------------------------------------------------------------------------
1745    SUBROUTINE g_uncouple_bdy(  field, g_field,  &
1746                                variable_in, config_flags, & 
1747                                spec_zone,       &
1748                                mu, g_mu, msf,   &
1749                                ids,ide, jds,jde, kds,kde,  & ! domain dims
1750                                ims,ime, jms,jme, kms,kme,  & ! memory dims
1751                                ips,ipe, jps,jpe, kps,kpe,  & ! patch  dims
1752                                its,ite, jts,jte, kts,kte )
1754 !  This subroutine adds the tendencies in the boundary specified region.
1755 !  spec_zone is the width of the outer specified b.c.s that are set here.
1756 !  (JD August 2000)
1758       IMPLICIT NONE
1760       INTEGER,      INTENT(IN   )    :: ids,ide, jds,jde, kds,kde
1761       INTEGER,      INTENT(IN   )    :: ims,ime, jms,jme, kms,kme
1762       INTEGER,      INTENT(IN   )    :: ips,ipe, jps,jpe, kps,kpe
1763       INTEGER,      INTENT(IN   )    :: its,ite, jts,jte, kts,kte
1764       INTEGER,      INTENT(IN   )    :: spec_zone
1765       CHARACTER,    INTENT(IN   )    :: variable_in
1766       REAL, DIMENSION( ims:ime , jms:jme ),    INTENT(IN   ) :: g_mu
1767       REAL, DIMENSION( ims:ime , jms:jme ),    INTENT(IN   ) :: mu
1768       REAL, DIMENSION( ims:ime , jms:jme ),    INTENT(IN   ) :: msf
1769       REAL, DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(INOUT) :: g_field
1770       REAL, DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(INOUT) :: field
1771       TYPE( grid_config_rec_type ) config_flags
1773       CHARACTER  :: variable
1774       INTEGER    :: i, j, k, ibs, ibe, jbs, jbe, itf, jtf, ktf
1775       INTEGER    :: b_dist, b_limit
1776       LOGICAL    :: periodic_x
1778       periodic_x = config_flags%periodic_x
1780       variable = variable_in
1782       IF (variable == 'U') variable = 'u'
1783       IF (variable == 'V') variable = 'v'
1784       IF (variable == 'T') variable = 't'
1785       IF (variable == 'H') variable = 'h'
1786       IF (variable == 'W') variable = 'w'
1788       ibs = ids
1789       ibe = ide-1
1790       itf = min(ite,ide-1)
1791       jbs = jds
1792       jbe = jde-1
1793       jtf = min(jte,jde-1)
1794       ktf = kde-1
1795       IF (variable == 'u') ibe = ide
1796       IF (variable == 'u') itf = min(ite,ide)
1797       IF (variable == 'v') jbe = jde
1798       IF (variable == 'v') jtf = min(jte,jde)
1799       IF (variable == 'h') ktf = kte
1800       IF (variable == 'w') ktf = kte
1802       IF (jts - jbs .lt. spec_zone) THEN
1803 ! Y-start boundary
1804         DO j = jts, min(jtf,jbs+spec_zone-1)
1805           b_dist = j - jbs
1806           b_limit = b_dist
1807           IF(periodic_x)b_limit = 0
1808           DO k = kts, ktf
1809             DO i = max(its,b_limit+ibs), min(itf,ibe-b_limit)
1810               if (variable == 't' .or. variable == 'h') then 
1811                  g_field(i,k,j) = g_field(i,k,j)/mu(i,j) &
1812                                 - field(i,k,j)*g_mu(i,j)/(mu(i,j)*mu(i,j))
1813                  field(i,k,j) = field(i,k,j)/mu(i,j)
1814               else 
1815                  g_field(i,k,j) = ( g_field(i,k,j)/mu(i,j) &
1816                                   - field(i,k,j)*g_mu(i,j)/(mu(i,j)*mu(i,j)) ) &
1817                                 * msf(i,j)
1818                  field(i,k,j) = field(i,k,j)/mu(i,j)*msf(i,j)
1819               end if
1820             ENDDO
1821           ENDDO
1822         ENDDO
1823       ENDIF 
1824       IF (jbe - jtf .lt. spec_zone) THEN 
1825 ! Y-end boundary 
1826         DO j = max(jts,jbe-spec_zone+1), jtf 
1827           b_dist = jbe - j 
1828           b_limit = b_dist
1829           IF(periodic_x)b_limit = 0
1830           DO k = kts, ktf 
1831             DO i = max(its,b_limit+ibs), min(itf,ibe-b_limit)
1832               if (variable == 't' .or. variable == 'h') then 
1833                  g_field(i,k,j) = g_field(i,k,j)/mu(i,j) &
1834                                 - field(i,k,j)*g_mu(i,j)/(mu(i,j)*mu(i,j))
1835                  field(i,k,j) = field(i,k,j)/mu(i,j)
1836               else 
1837                  g_field(i,k,j) = ( g_field(i,k,j)/mu(i,j) &
1838                                   - field(i,k,j)*g_mu(i,j)/(mu(i,j)*mu(i,j)) ) &
1839                                 * msf(i,j)
1840                  field(i,k,j) = field(i,k,j)/mu(i,j)*msf(i,j)
1841               end if
1842             ENDDO
1843           ENDDO
1844         ENDDO
1845       ENDIF 
1847     IF(.NOT.periodic_x)THEN
1848       IF (its - ibs .lt. spec_zone) THEN
1849 ! X-start boundary
1850         DO i = its, min(itf,ibs+spec_zone-1)
1851           b_dist = i - ibs
1852           DO k = kts, ktf
1853             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
1854               if (variable == 't' .or. variable == 'h') then 
1855                  g_field(i,k,j) = g_field(i,k,j)/mu(i,j) &
1856                                 - field(i,k,j)*g_mu(i,j)/(mu(i,j)*mu(i,j))
1857                  field(i,k,j) = field(i,k,j)/mu(i,j)
1858               else 
1859                  g_field(i,k,j) = ( g_field(i,k,j)/mu(i,j) &
1860                                   - field(i,k,j)*g_mu(i,j)/(mu(i,j)*mu(i,j)) ) &
1861                                 * msf(i,j)
1862                  field(i,k,j) = field(i,k,j)/mu(i,j)*msf(i,j)
1863               end if
1864             ENDDO
1865           ENDDO
1866         ENDDO
1867       ENDIF 
1869       IF (ibe - itf .lt. spec_zone) THEN
1870 ! X-end boundary
1871         DO i = max(its,ibe-spec_zone+1), itf
1872           b_dist = ibe - i
1873           DO k = kts, ktf
1874             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
1875               if (variable == 't' .or. variable == 'h') then 
1876                  g_field(i,k,j) = g_field(i,k,j)/mu(i,j) &
1877                                 - field(i,k,j)*g_mu(i,j)/(mu(i,j)*mu(i,j))
1878                  field(i,k,j) = field(i,k,j)/mu(i,j)
1879               else 
1880                  g_field(i,k,j) = ( g_field(i,k,j)/mu(i,j) &
1881                                   - field(i,k,j)*g_mu(i,j)/(mu(i,j)*mu(i,j)) ) &
1882                                 * msf(i,j)
1883                  field(i,k,j) = field(i,k,j)/mu(i,j)*msf(i,j)
1884               end if
1885             ENDDO
1886           ENDDO
1887         ENDDO
1888       ENDIF 
1889     ENDIF
1891    END SUBROUTINE g_uncouple_bdy 
1892 !------------------------------------------------------------------------
1894    SUBROUTINE g_spec_bdyupdate(  field, g_field,  &
1895                                  field_tend, g_field_tend, dt,            &
1896                                  variable_in, config_flags, & 
1897                                  spec_zone,                  &
1898                                  ids,ide, jds,jde, kds,kde,  & ! domain dims
1899                                  ims,ime, jms,jme, kms,kme,  & ! memory dims
1900                                  ips,ipe, jps,jpe, kps,kpe,  & ! patch  dims
1901                                  its,ite, jts,jte, kts,kte )
1903 !  This subroutine adds the tendencies in the boundary specified region.
1904 !  spec_zone is the width of the outer specified b.c.s that are set here.
1905 !  (JD August 2000)
1907       IMPLICIT NONE
1909       INTEGER,      INTENT(IN   )    :: ids,ide, jds,jde, kds,kde
1910       INTEGER,      INTENT(IN   )    :: ims,ime, jms,jme, kms,kme
1911       INTEGER,      INTENT(IN   )    :: ips,ipe, jps,jpe, kps,kpe
1912       INTEGER,      INTENT(IN   )    :: its,ite, jts,jte, kts,kte
1913       INTEGER,      INTENT(IN   )    :: spec_zone
1914       CHARACTER,    INTENT(IN   )    :: variable_in
1915       REAL,         INTENT(IN   )    :: dt
1918       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(INOUT) :: g_field
1919       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(INOUT) :: field
1920       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(IN   ) :: g_field_tend
1921       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(IN   ) :: field_tend
1922       TYPE( grid_config_rec_type ) config_flags
1924       CHARACTER  :: variable
1925       INTEGER    :: i, j, k, ibs, ibe, jbs, jbe, itf, jtf, ktf
1926       INTEGER    :: b_dist, b_limit
1927       LOGICAL    :: periodic_x
1929       periodic_x = config_flags%periodic_x
1931       variable = variable_in
1933       IF (variable == 'U') variable = 'u'
1934       IF (variable == 'V') variable = 'v'
1935       IF (variable == 'M') variable = 'm'
1936       IF (variable == 'H') variable = 'h'
1938       ibs = ids
1939       ibe = ide-1
1940       itf = min(ite,ide-1)
1941       jbs = jds
1942       jbe = jde-1
1943       jtf = min(jte,jde-1)
1944       ktf = kde-1
1945       IF (variable == 'u') ibe = ide
1946       IF (variable == 'u') itf = min(ite,ide)
1947       IF (variable == 'v') jbe = jde
1948       IF (variable == 'v') jtf = min(jte,jde)
1949       IF (variable == 'm') ktf = kte
1950       IF (variable == 'h') ktf = kte
1952       IF (jts - jbs .lt. spec_zone) THEN
1953 ! Y-start boundary
1954         DO j = jts, min(jtf,jbs+spec_zone-1)
1955           b_dist = j - jbs
1956           b_limit = b_dist
1957           IF(periodic_x)b_limit = 0
1958           DO k = kts, ktf
1959             DO i = max(its,b_limit+ibs), min(itf,ibe-b_limit)
1960               g_field(i,k,j) = g_field(i,k,j) + dt*g_field_tend(i,k,j) 
1961               field(i,k,j) = field(i,k,j) + dt*field_tend(i,k,j) 
1962             ENDDO
1963           ENDDO
1964         ENDDO
1965       ENDIF 
1966       IF (jbe - jtf .lt. spec_zone) THEN 
1967 ! Y-end boundary 
1968         DO j = max(jts,jbe-spec_zone+1), jtf 
1969           b_dist = jbe - j 
1970           b_limit = b_dist
1971           IF(periodic_x)b_limit = 0
1972           DO k = kts, ktf 
1973             DO i = max(its,b_limit+ibs), min(itf,ibe-b_limit)
1974               g_field(i,k,j) = g_field(i,k,j) + dt*g_field_tend(i,k,j) 
1975               field(i,k,j) = field(i,k,j) + dt*field_tend(i,k,j) 
1976             ENDDO
1977           ENDDO
1978         ENDDO
1979       ENDIF 
1981     IF(.NOT.periodic_x)THEN
1982       IF (its - ibs .lt. spec_zone) THEN
1983 ! X-start boundary
1984         DO i = its, min(itf,ibs+spec_zone-1)
1985           b_dist = i - ibs
1986           DO k = kts, ktf
1987             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
1988               g_field(i,k,j) = g_field(i,k,j) + dt*g_field_tend(i,k,j) 
1989               field(i,k,j) = field(i,k,j) + dt*field_tend(i,k,j) 
1990             ENDDO
1991           ENDDO
1992         ENDDO
1993       ENDIF 
1995       IF (ibe - itf .lt. spec_zone) THEN
1996 ! X-end boundary
1997         DO i = max(its,ibe-spec_zone+1), itf
1998           b_dist = ibe - i
1999           DO k = kts, ktf
2000             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
2001               g_field(i,k,j) = g_field(i,k,j) + dt*g_field_tend(i,k,j) 
2002               field(i,k,j) = field(i,k,j) + dt*field_tend(i,k,j) 
2003             ENDDO
2004           ENDDO
2005         ENDDO
2006       ENDIF 
2007     ENDIF
2009    END SUBROUTINE g_spec_bdyupdate
2010 !------------------------------------------------------------------------
2011 !        Generated by TAPENADE     (INRIA, Tropics team)
2012 !  Tapenade 3.10 (r5498) - 20 Jan 2015 09:48
2014 !  Differentiation of spec_bdy_final in forward (tangent) mode:
2015 !   variations   of useful results: field
2016 !   with respect to varying inputs: field field_bdy_xe field_bdy_tend_xe
2017 !                field_bdy_xs field_bdy_tend_xs field_bdy_ye field_bdy_tend_ye
2018 !                field_bdy_ys field_bdy_tend_ys mu
2019 !   RW status of diff variables: field:in-out field_bdy_xe:in field_bdy_tend_xe:in
2020 !                field_bdy_xs:in field_bdy_tend_xs:in field_bdy_ye:in
2021 !                field_bdy_tend_ye:in field_bdy_ys:in field_bdy_tend_ys:in
2022 !                mu:in
2023 ! domain dims
2024 ! memory dims
2025 ! patch  dims
2026 SUBROUTINE g_SPEC_BDY_FINAL(field, fieldd, mu, mud, msf, field_bdy_xs, &
2027 & field_bdy_xsd, field_bdy_xe, field_bdy_xed, field_bdy_ys, &
2028 & field_bdy_ysd, field_bdy_ye, field_bdy_yed, field_bdy_tend_xs, &
2029 & field_bdy_tend_xsd, field_bdy_tend_xe, field_bdy_tend_xed, &
2030 & field_bdy_tend_ys, field_bdy_tend_ysd, field_bdy_tend_ye, &
2031 & field_bdy_tend_yed, variable_in, config_flags, spec_bdy_width, &
2032 & spec_zone, dtbc, ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms&
2033 & , kme, ips, ipe, jps, jpe, kps, kpe, its, ite, jts, jte, kts, kte)
2034   IMPLICIT NONE
2035   INTEGER, INTENT(IN) :: ids, ide, jds, jde, kds, kde
2036   INTEGER, INTENT(IN) :: ims, ime, jms, jme, kms, kme
2037   INTEGER, INTENT(IN) :: ips, ipe, jps, jpe, kps, kpe
2038   INTEGER, INTENT(IN) :: its, ite, jts, jte, kts, kte
2039   INTEGER, INTENT(IN) :: spec_bdy_width, spec_zone
2040   REAL, INTENT(IN) :: dtbc
2041   CHARACTER, INTENT(IN) :: variable_in
2042   REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: field
2043   REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: fieldd
2044   REAL, DIMENSION(ims:ime, jms:jme), INTENT(IN) :: mu, msf
2045   REAL, DIMENSION(ims:ime, jms:jme), INTENT(IN) :: mud
2046   REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: &
2047 & field_bdy_xs, field_bdy_xe
2048   REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: &
2049 & field_bdy_xsd, field_bdy_xed
2050   REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: &
2051 & field_bdy_ys, field_bdy_ye
2052   REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: &
2053 & field_bdy_ysd, field_bdy_yed
2054   REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: &
2055 & field_bdy_tend_xs, field_bdy_tend_xe
2056   REAL, DIMENSION(jms:jme, kds:kde, spec_bdy_width), INTENT(IN) :: &
2057 & field_bdy_tend_xsd, field_bdy_tend_xed
2058   REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: &
2059 & field_bdy_tend_ys, field_bdy_tend_ye
2060   REAL, DIMENSION(ims:ime, kds:kde, spec_bdy_width), INTENT(IN) :: &
2061 & field_bdy_tend_ysd, field_bdy_tend_yed
2062   TYPE(GRID_CONFIG_REC_TYPE) :: config_flags
2063   CHARACTER :: variable
2064   INTEGER :: i, j, k, ibs, ibe, jbs, jbe, itf, jtf, ktf, im1, ip1
2065   INTEGER :: b_dist, b_limit
2066   REAL :: bfield, xmsf, xmu
2067   REAL :: bfieldd, xmud
2068   LOGICAL :: periodic_x, msfcouple, mucouple
2069   INTEGER :: min6
2070   INTEGER :: min5
2071   INTEGER :: min4
2072   INTEGER :: min3
2073   INTEGER :: min2
2074   INTEGER :: min1
2075   INTEGER :: max6
2076   INTEGER :: max5
2077   INTEGER :: max4
2078   INTEGER :: max3
2079   INTEGER :: max2
2080   INTEGER :: max1
2081   periodic_x = config_flags%periodic_x
2082   variable = variable_in
2083   IF (variable .EQ. 'U') variable = 'u'
2084   IF (variable .EQ. 'V') variable = 'v'
2085   IF (variable .EQ. 'W') variable = 'w'
2086   IF (variable .EQ. 'M') variable = 'm'
2087   IF (variable .EQ. 'T') variable = 't'
2088   IF (variable .EQ. 'H') variable = 'h'
2089   ibs = ids
2090   ibe = ide - 1
2091   IF (ite .GT. ide - 1) THEN
2092     itf = ide - 1
2093   ELSE
2094     itf = ite
2095   END IF
2096   jbs = jds
2097   jbe = jde - 1
2098   IF (jte .GT. jde - 1) THEN
2099     jtf = jde - 1
2100   ELSE
2101     jtf = jte
2102   END IF
2103   ktf = kde - 1
2104   IF (variable .EQ. 'u') ibe = ide
2105   IF (variable .EQ. 'u') THEN
2106     IF (ite .GT. ide) THEN
2107       itf = ide
2108     ELSE
2109       itf = ite
2110     END IF
2111   END IF
2112   IF (variable .EQ. 'v') jbe = jde
2113   IF (variable .EQ. 'v') THEN
2114     IF (jte .GT. jde) THEN
2115       jtf = jde
2116     ELSE
2117       jtf = jte
2118     END IF
2119   END IF
2120   IF (variable .EQ. 'm') ktf = kde
2121   IF (variable .EQ. 'h') ktf = kde
2122   IF (variable .EQ. 'w') ktf = kde
2123   msfcouple = .false.
2124   mucouple = .true.
2125   IF ((variable .EQ. 'u' .OR. variable .EQ. 'v') .OR. variable .EQ. 'w'&
2126 & ) msfcouple = .true.
2127   IF (variable .EQ. 'm') mucouple = .false.
2128   xmsf = 1.
2129   xmu = 1.
2130   IF (jts - jbs .LT. spec_zone) THEN
2131     IF (jtf .GT. jbs + spec_zone - 1) THEN
2132       min1 = jbs + spec_zone - 1
2133       xmud = 0.0
2134     ELSE
2135       min1 = jtf
2136       xmud = 0.0
2137     END IF
2138 ! Y-start boundary
2139     DO j=jts,min1
2140       b_dist = j - jbs
2141       b_limit = b_dist
2142       IF (periodic_x) b_limit = 0
2143       DO k=kts,ktf
2144         IF (its .LT. b_limit + ibs) THEN
2145           max1 = b_limit + ibs
2146         ELSE
2147           max1 = its
2148         END IF
2149         IF (itf .GT. ibe - b_limit) THEN
2150           min3 = ibe - b_limit
2151         ELSE
2152           min3 = itf
2153         END IF
2154         DO i=max1,min3
2155           bfieldd = field_bdy_ysd(i, k, b_dist+1) + dtbc*&
2156 &           field_bdy_tend_ysd(i, k, b_dist+1)
2157           bfield = field_bdy_ys(i, k, b_dist+1) + dtbc*field_bdy_tend_ys&
2158 &           (i, k, b_dist+1)
2159           IF (msfcouple) xmsf = msf(i, j)
2160           IF (mucouple) THEN
2161             xmud = mud(i, j)
2162             xmu = mu(i, j)
2163           END IF
2164           fieldd(i, k, j) = (xmsf*bfieldd*xmu-xmsf*bfield*xmud)/xmu**2
2165           field(i, k, j) = xmsf*bfield/xmu
2166         END DO
2167       END DO
2168     END DO
2169   ELSE
2170     xmud = 0.0
2171   END IF
2172   IF (jbe - jtf .LT. spec_zone) THEN
2173     IF (jts .LT. jbe - spec_zone + 1) THEN
2174       max2 = jbe - spec_zone + 1
2175     ELSE
2176       max2 = jts
2177     END IF
2178 ! Y-end boundary
2179     DO j=max2,jtf
2180       b_dist = jbe - j
2181       b_limit = b_dist
2182       IF (periodic_x) b_limit = 0
2183       DO k=kts,ktf
2184         IF (its .LT. b_limit + ibs) THEN
2185           max3 = b_limit + ibs
2186         ELSE
2187           max3 = its
2188         END IF
2189         IF (itf .GT. ibe - b_limit) THEN
2190           min4 = ibe - b_limit
2191         ELSE
2192           min4 = itf
2193         END IF
2194         DO i=max3,min4
2195           bfieldd = field_bdy_yed(i, k, b_dist+1) + dtbc*&
2196 &           field_bdy_tend_yed(i, k, b_dist+1)
2197           bfield = field_bdy_ye(i, k, b_dist+1) + dtbc*field_bdy_tend_ye&
2198 &           (i, k, b_dist+1)
2199           IF (msfcouple) xmsf = msf(i, j)
2200           IF (mucouple) THEN
2201             xmud = mud(i, j)
2202             xmu = mu(i, j)
2203           END IF
2204           fieldd(i, k, j) = (xmsf*bfieldd*xmu-xmsf*bfield*xmud)/xmu**2
2205           field(i, k, j) = xmsf*bfield/xmu
2206         END DO
2207       END DO
2208     END DO
2209   END IF
2210   IF (.NOT.periodic_x) THEN
2211     IF (its - ibs .LT. spec_zone) THEN
2212       IF (itf .GT. ibs + spec_zone - 1) THEN
2213         min2 = ibs + spec_zone - 1
2214       ELSE
2215         min2 = itf
2216       END IF
2217 ! X-start boundary
2218       DO i=its,min2
2219         b_dist = i - ibs
2220         DO k=kts,ktf
2221           IF (jts .LT. b_dist + jbs + 1) THEN
2222             max4 = b_dist + jbs + 1
2223           ELSE
2224             max4 = jts
2225           END IF
2226           IF (jtf .GT. jbe - b_dist - 1) THEN
2227             min5 = jbe - b_dist - 1
2228           ELSE
2229             min5 = jtf
2230           END IF
2231           DO j=max4,min5
2232             bfieldd = field_bdy_xsd(j, k, b_dist+1) + dtbc*&
2233 &             field_bdy_tend_xsd(j, k, b_dist+1)
2234             bfield = field_bdy_xs(j, k, b_dist+1) + dtbc*&
2235 &             field_bdy_tend_xs(j, k, b_dist+1)
2236             IF (msfcouple) xmsf = msf(i, j)
2237             IF (mucouple) THEN
2238               xmud = mud(i, j)
2239               xmu = mu(i, j)
2240             END IF
2241             fieldd(i, k, j) = (xmsf*bfieldd*xmu-xmsf*bfield*xmud)/xmu**2
2242             field(i, k, j) = xmsf*bfield/xmu
2243           END DO
2244         END DO
2245       END DO
2246     END IF
2247     IF (ibe - itf .LT. spec_zone) THEN
2248       IF (its .LT. ibe - spec_zone + 1) THEN
2249         max5 = ibe - spec_zone + 1
2250       ELSE
2251         max5 = its
2252       END IF
2253 ! X-end boundary
2254       DO i=max5,itf
2255         b_dist = ibe - i
2256         DO k=kts,ktf
2257           IF (jts .LT. b_dist + jbs + 1) THEN
2258             max6 = b_dist + jbs + 1
2259           ELSE
2260             max6 = jts
2261           END IF
2262           IF (jtf .GT. jbe - b_dist - 1) THEN
2263             min6 = jbe - b_dist - 1
2264           ELSE
2265             min6 = jtf
2266           END IF
2267           DO j=max6,min6
2268             bfieldd = field_bdy_xed(j, k, b_dist+1) + dtbc*&
2269 &             field_bdy_tend_xed(j, k, b_dist+1)
2270             bfield = field_bdy_xe(j, k, b_dist+1) + dtbc*&
2271 &             field_bdy_tend_xe(j, k, b_dist+1)
2272             IF (msfcouple) xmsf = msf(i, j)
2273             IF (mucouple) THEN
2274               xmud = mud(i, j)
2275               xmu = mu(i, j)
2276             END IF
2277             fieldd(i, k, j) = (xmsf*bfieldd*xmu-xmsf*bfield*xmud)/xmu**2
2278             field(i, k, j) = xmsf*bfield/xmu
2279           END DO
2280         END DO
2281       END DO
2282     END IF
2283   END IF
2284 END SUBROUTINE g_SPEC_BDY_FINAL
2285 !------------------------------------------------------------------------
2287    SUBROUTINE g_zero_grad_bdy (  field, g_field,                    &
2288                                variable_in, config_flags, & 
2289                                spec_zone,                  &
2290                                ids,ide, jds,jde, kds,kde,  & ! domain dims
2291                                ims,ime, jms,jme, kms,kme,  & ! memory dims
2292                                ips,ipe, jps,jpe, kps,kpe,  & ! patch  dims
2293                                its,ite, jts,jte, kts,kte )
2295 !  This subroutine sets zero gradient conditions in the boundary specified region.
2296 !  spec_zone is the width of the outer specified b.c.s that are set here.
2297 !  (JD August 2000)
2299       IMPLICIT NONE
2301       INTEGER,      INTENT(IN   )    :: ids,ide, jds,jde, kds,kde
2302       INTEGER,      INTENT(IN   )    :: ims,ime, jms,jme, kms,kme
2303       INTEGER,      INTENT(IN   )    :: ips,ipe, jps,jpe, kps,kpe
2304       INTEGER,      INTENT(IN   )    :: its,ite, jts,jte, kts,kte
2305       INTEGER,      INTENT(IN   )    :: spec_zone
2306       CHARACTER,    INTENT(IN   )    :: variable_in
2309       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(INOUT) :: g_field
2310       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(INOUT) :: field
2311       TYPE( grid_config_rec_type ) config_flags
2313       CHARACTER  :: variable
2314       INTEGER    :: i, j, k, ibs, ibe, jbs, jbe, itf, jtf, ktf, i_inner, j_inner
2315       INTEGER    :: b_dist, b_limit
2316       LOGICAL    :: periodic_x
2318       periodic_x = config_flags%periodic_x
2320       variable = variable_in
2322       IF (variable == 'U') variable = 'u'
2323       IF (variable == 'V') variable = 'v'
2325       ibs = ids
2326       ibe = ide-1
2327       itf = min(ite,ide-1)
2328       jbs = jds
2329       jbe = jde-1
2330       jtf = min(jte,jde-1)
2331       ktf = kde-1
2332       IF (variable == 'u') ibe = ide
2333       IF (variable == 'u') itf = min(ite,ide)
2334       IF (variable == 'v') jbe = jde
2335       IF (variable == 'v') jtf = min(jte,jde)
2336       IF (variable == 'w') ktf = kde
2338       IF (jts - jbs .lt. spec_zone) THEN
2339 ! Y-start boundary
2340         DO j = jts, min(jtf,jbs+spec_zone-1)
2341           b_dist = j - jbs
2342           b_limit = b_dist
2343           IF(periodic_x)b_limit = 0
2344           DO k = kts, ktf
2345             DO i = max(its,b_limit+ibs), min(itf,ibe-b_limit)
2346               i_inner = max(i,ibs+spec_zone)
2347               i_inner = min(i_inner,ibe-spec_zone)
2348               IF(periodic_x)i_inner = i
2349               g_field(i,k,j) = g_field(i_inner,k,jbs+spec_zone)
2350               field(i,k,j) = field(i_inner,k,jbs+spec_zone)
2351             ENDDO
2352           ENDDO
2353         ENDDO
2354       ENDIF 
2355       IF (jbe - jtf .lt. spec_zone) THEN 
2356 ! Y-end boundary 
2357         DO j = max(jts,jbe-spec_zone+1), jtf 
2358           b_dist = jbe - j 
2359           b_limit = b_dist
2360           IF(periodic_x)b_limit = 0
2361           DO k = kts, ktf 
2362             DO i = max(its,b_limit+ibs), min(itf,ibe-b_limit)
2363               i_inner = max(i,ibs+spec_zone)
2364               i_inner = min(i_inner,ibe-spec_zone)
2365               IF(periodic_x)i_inner = i
2366               g_field(i,k,j) = g_field(i_inner,k,jbe-spec_zone)
2367               field(i,k,j) = field(i_inner,k,jbe-spec_zone)
2368             ENDDO
2369           ENDDO
2370         ENDDO
2371       ENDIF 
2373     IF(.NOT.periodic_x)THEN
2374       IF (its - ibs .lt. spec_zone) THEN
2375 ! X-start boundary
2376         DO i = its, min(itf,ibs+spec_zone-1)
2377           b_dist = i - ibs
2378           DO k = kts, ktf
2379             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
2380               j_inner = max(j,jbs+spec_zone)
2381               j_inner = min(j_inner,jbe-spec_zone)
2382               g_field(i,k,j) = g_field(ibs+spec_zone,k,j_inner)
2383               field(i,k,j) = field(ibs+spec_zone,k,j_inner)
2384             ENDDO
2385           ENDDO
2386         ENDDO
2387       ENDIF 
2389       IF (ibe - itf .lt. spec_zone) THEN
2390 ! X-end boundary
2391         DO i = max(its,ibe-spec_zone+1), itf
2392           b_dist = ibe - i
2393           DO k = kts, ktf
2394             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
2395               j_inner = max(j,jbs+spec_zone)
2396               j_inner = min(j_inner,jbe-spec_zone)
2397               g_field(i,k,j) = g_field(ibe-spec_zone,k,j_inner)
2398               field(i,k,j) = field(ibe-spec_zone,k,j_inner)
2399             ENDDO
2400           ENDDO
2401         ENDDO
2402       ENDIF 
2403     ENDIF
2405    END SUBROUTINE g_zero_grad_bdy
2406 !------------------------------------------------------------------------
2408    SUBROUTINE g_flow_dep_bdy ( field, g_field,     &
2409                                u, v, config_flags, & 
2410                                spec_zone,                  &
2411                                ids,ide, jds,jde, kds,kde,  & ! domain dims
2412                                ims,ime, jms,jme, kms,kme,  & ! memory dims
2413                                ips,ipe, jps,jpe, kps,kpe,  & ! patch  dims
2414                                its,ite, jts,jte, kts,kte )
2416       IMPLICIT NONE
2418       INTEGER, INTENT(IN) :: ids,ide, jds,jde, kds,kde
2419       INTEGER, INTENT(IN) :: ims,ime, jms,jme, kms,kme
2420       INTEGER, INTENT(IN) :: ips,ipe, jps,jpe, kps,kpe
2421       INTEGER, INTENT(IN) :: its,ite, jts,jte, kts,kte
2422       INTEGER, INTENT(IN) :: spec_zone
2424       REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: g_field
2425       REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: field
2426       REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(IN   ) :: u
2427       REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(IN   ) :: v
2428       TYPE(grid_config_rec_type), INTENT(IN) :: config_flags
2430       INTEGER    :: i, j, k, ibs, ibe, jbs, jbe, itf, jtf, ktf, i_inner, j_inner
2431       INTEGER    :: b_dist, b_limit
2432       LOGICAL    :: periodic_x
2435       periodic_x = config_flags%periodic_x
2437       ibs = ids
2438       ibe = ide-1
2439       itf = min(ite,ide-1)
2440       jbs = jds
2441       jbe = jde-1
2442       jtf = min(jte,jde-1)
2443       ktf = kde-1
2445       IF (jts - jbs .lt. spec_zone) THEN
2446 ! Y-start boundary
2447         DO j = jts, min(jtf,jbs+spec_zone-1)
2448           b_dist = j - jbs
2449           b_limit = b_dist
2450           IF(periodic_x)b_limit = 0
2451           DO k = kts, ktf
2452             DO i = max(its,b_limit+ibs), min(itf,ibe-b_limit)
2453               i_inner = max(i,ibs+spec_zone)
2454               i_inner = min(i_inner,ibe-spec_zone)
2455               IF(periodic_x)i_inner = i
2456               IF(v(i,k,j) .lt. 0.)THEN
2457                 g_field(i,k,j) = g_field(i_inner,k,jbs+spec_zone)
2458                 field(i,k,j) = field(i_inner,k,jbs+spec_zone)
2459               ELSE
2460                 g_field(i,k,j) = 0.
2461                 field(i,k,j) = 0.
2462               ENDIF
2463             ENDDO
2464           ENDDO
2465         ENDDO
2466       ENDIF 
2467       IF (jbe - jtf .lt. spec_zone) THEN 
2468 ! Y-end boundary 
2469         DO j = max(jts,jbe-spec_zone+1), jtf 
2470           b_dist = jbe - j 
2471           b_limit = b_dist
2472           IF(periodic_x)b_limit = 0
2473           DO k = kts, ktf 
2474             DO i = max(its,b_limit+ibs), min(itf,ibe-b_limit)
2475               i_inner = max(i,ibs+spec_zone)
2476               i_inner = min(i_inner,ibe-spec_zone)
2477               IF(periodic_x)i_inner = i
2478               IF(v(i,k,j+1) .gt. 0.)THEN
2479                 g_field(i,k,j) = g_field(i_inner,k,jbe-spec_zone)
2480                 field(i,k,j) = field(i_inner,k,jbe-spec_zone)
2481               ELSE
2482                 g_field(i,k,j) = 0.
2483                 field(i,k,j) = 0.
2484               ENDIF
2485             ENDDO
2486           ENDDO
2487         ENDDO
2488       ENDIF 
2490     IF(.NOT.periodic_x)THEN
2491       IF (its - ibs .lt. spec_zone) THEN
2492 ! X-start boundary
2493         DO i = its, min(itf,ibs+spec_zone-1)
2494           b_dist = i - ibs
2495           DO k = kts, ktf
2496             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
2497               j_inner = max(j,jbs+spec_zone)
2498               j_inner = min(j_inner,jbe-spec_zone)
2499               IF(u(i,k,j) .lt. 0.)THEN
2500                 g_field(i,k,j) = g_field(ibs+spec_zone,k,j_inner)
2501                 field(i,k,j) = field(ibs+spec_zone,k,j_inner)
2502               ELSE
2503                 g_field(i,k,j) = 0.
2504                 field(i,k,j) = 0.
2505               ENDIF
2506             ENDDO
2507           ENDDO
2508         ENDDO
2509       ENDIF 
2511       IF (ibe - itf .lt. spec_zone) THEN
2512 ! X-end boundary
2513         DO i = max(its,ibe-spec_zone+1), itf
2514           b_dist = ibe - i
2515           DO k = kts, ktf
2516             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
2517               j_inner = max(j,jbs+spec_zone)
2518               j_inner = min(j_inner,jbe-spec_zone)
2519               IF(u(i+1,k,j) .gt. 0.)THEN
2520                 g_field(i,k,j) = g_field(ibe-spec_zone,k,j_inner)
2521                 field(i,k,j) = field(ibe-spec_zone,k,j_inner)
2522               ELSE
2523                 g_field(i,k,j) = 0.
2524                 field(i,k,j) = 0.
2525               ENDIF
2526             ENDDO
2527           ENDDO
2528         ENDDO
2529       ENDIF 
2530     ENDIF
2532    END SUBROUTINE g_flow_dep_bdy
2534 !        Generated by TAPENADE     (INRIA, Tropics team)
2535 !  Tapenade 3.6 (r4165) - 21 sep 2011 20:54
2537 !  Differentiation of flow_dep_bdy_qnn in forward (tangent) mode:
2538 !   variations   of useful results: field
2539 !   with respect to varying inputs: field
2540 !   RW status of diff variables: field:in-out
2541 ! domain dims
2542 ! memory dims
2543 ! patch  dims
2544 SUBROUTINE G_FLOW_DEP_BDY_QNN(field, fieldd, u, v, config_flags, &
2545 &  spec_zone, ccn_conc, ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme&
2546 &  , ips, ipe, jps, jpe, kps, kpe, its, ite, jts, jte, kts, kte)
2547   IMPLICIT NONE
2548   INTEGER, INTENT(IN) :: ids, ide, jds, jde, kds, kde
2549   INTEGER, INTENT(IN) :: ims, ime, jms, jme, kms, kme
2550   INTEGER, INTENT(IN) :: ips, ipe, jps, jpe, kps, kpe
2551   INTEGER, INTENT(IN) :: its, ite, jts, jte, kts, kte
2552   INTEGER, INTENT(IN) :: spec_zone
2553   REAL,    INTENT(IN) :: ccn_conc ! RAS
2554   REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: field
2555   REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(INOUT) :: fieldd
2556   REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(IN) :: u
2557   REAL, DIMENSION(ims:ime, kms:kme, jms:jme), INTENT(IN) :: v
2558   TYPE(GRID_CONFIG_REC_TYPE) :: config_flags
2559   INTEGER :: i, j, k, ibs, ibe, jbs, jbe, itf, jtf, ktf, i_inner, &
2560 &  j_inner
2561   INTEGER :: b_dist, b_limit
2562   LOGICAL :: periodic_x
2563   INTEGER :: min6
2564   INTEGER :: min5
2565   INTEGER :: min4
2566   INTEGER :: min3
2567   INTEGER :: min2
2568   INTEGER :: min1
2569   INTRINSIC MAX
2570   INTRINSIC MIN
2571   INTEGER :: max6
2572   INTEGER :: max5
2573   INTEGER :: max4
2574   INTEGER :: max3
2575   INTEGER :: max2
2576   INTEGER :: max1
2577   periodic_x = config_flags%periodic_x
2578   ibs = ids
2579   ibe = ide - 1
2580   IF (ite .GT. ide - 1) THEN
2581     itf = ide - 1
2582   ELSE
2583     itf = ite
2584   END IF
2585   jbs = jds
2586   jbe = jde - 1
2587   IF (jte .GT. jde - 1) THEN
2588     jtf = jde - 1
2589   ELSE
2590     jtf = jte
2591   END IF
2592   ktf = kde - 1
2593   IF (jts - jbs .LT. spec_zone) THEN
2594     IF (jtf .GT. jbs + spec_zone - 1) THEN
2595       min1 = jbs + spec_zone - 1
2596     ELSE
2597       min1 = jtf
2598     END IF
2599 ! Y-start boundary
2600     DO j=jts,min1
2601       b_dist = j - jbs
2602       b_limit = b_dist
2603       IF (periodic_x) b_limit = 0
2604       DO k=kts,ktf
2605         IF (its .LT. b_limit + ibs) THEN
2606           max1 = b_limit + ibs
2607         ELSE
2608           max1 = its
2609         END IF
2610         IF (itf .GT. ibe - b_limit) THEN
2611           min3 = ibe - b_limit
2612         ELSE
2613           min3 = itf
2614         END IF
2615         DO i=max1,min3
2616           IF (i .LT. ibs + spec_zone) THEN
2617             i_inner = ibs + spec_zone
2618           ELSE
2619             i_inner = i
2620           END IF
2621           IF (i_inner .GT. ibe - spec_zone) THEN
2622             i_inner = ibe - spec_zone
2623           ELSE
2624             i_inner = i_inner
2625           END IF
2626           IF (periodic_x) i_inner = i
2627           IF (v(i, k, j) .LT. 0.) THEN
2628             fieldd(i, k, j) = fieldd(i_inner, k, jbs+spec_zone)
2629             field(i, k, j) = field(i_inner, k, jbs+spec_zone)
2630           ELSE
2631             fieldd(i, k, j) = 0.0
2632             field(i, k, j) = ccn_conc
2633           END IF
2634         END DO
2635       END DO
2636     END DO
2637   END IF
2638   IF (jbe - jtf .LT. spec_zone) THEN
2639     IF (jts .LT. jbe - spec_zone + 1) THEN
2640       max2 = jbe - spec_zone + 1
2641     ELSE
2642       max2 = jts
2643     END IF
2644 ! Y-end boundary
2645     DO j=max2,jtf
2646       b_dist = jbe - j
2647       b_limit = b_dist
2648       IF (periodic_x) b_limit = 0
2649       DO k=kts,ktf
2650         IF (its .LT. b_limit + ibs) THEN
2651           max3 = b_limit + ibs
2652         ELSE
2653           max3 = its
2654         END IF
2655         IF (itf .GT. ibe - b_limit) THEN
2656           min4 = ibe - b_limit
2657         ELSE
2658           min4 = itf
2659         END IF
2660         DO i=max3,min4
2661           IF (i .LT. ibs + spec_zone) THEN
2662             i_inner = ibs + spec_zone
2663           ELSE
2664             i_inner = i
2665           END IF
2666           IF (i_inner .GT. ibe - spec_zone) THEN
2667             i_inner = ibe - spec_zone
2668           ELSE
2669             i_inner = i_inner
2670           END IF
2671           IF (periodic_x) i_inner = i
2672           IF (v(i, k, j+1) .GT. 0.) THEN
2673             fieldd(i, k, j) = fieldd(i_inner, k, jbe-spec_zone)
2674             field(i, k, j) = field(i_inner, k, jbe-spec_zone)
2675           ELSE
2676             fieldd(i, k, j) = 0.0
2677             field(i, k, j) = ccn_conc
2678           END IF
2679         END DO
2680       END DO
2681     END DO
2682   END IF
2683   IF (.NOT.periodic_x) THEN
2684     IF (its - ibs .LT. spec_zone) THEN
2685       IF (itf .GT. ibs + spec_zone - 1) THEN
2686         min2 = ibs + spec_zone - 1
2687       ELSE
2688         min2 = itf
2689       END IF
2690 ! X-start boundary
2691       DO i=its,min2
2692         b_dist = i - ibs
2693         DO k=kts,ktf
2694           IF (jts .LT. b_dist + jbs + 1) THEN
2695             max4 = b_dist + jbs + 1
2696           ELSE
2697             max4 = jts
2698           END IF
2699           IF (jtf .GT. jbe - b_dist - 1) THEN
2700             min5 = jbe - b_dist - 1
2701           ELSE
2702             min5 = jtf
2703           END IF
2704           DO j=max4,min5
2705             IF (j .LT. jbs + spec_zone) THEN
2706               j_inner = jbs + spec_zone
2707             ELSE
2708               j_inner = j
2709             END IF
2710             IF (j_inner .GT. jbe - spec_zone) THEN
2711               j_inner = jbe - spec_zone
2712             ELSE
2713               j_inner = j_inner
2714             END IF
2715             IF (u(i, k, j) .LT. 0.) THEN
2716               fieldd(i, k, j) = fieldd(ibs+spec_zone, k, j_inner)
2717               field(i, k, j) = field(ibs+spec_zone, k, j_inner)
2718             ELSE
2719               fieldd(i, k, j) = 0.0
2720               field(i, k, j) = ccn_conc
2721             END IF
2722           END DO
2723         END DO
2724       END DO
2725     END IF
2726     IF (ibe - itf .LT. spec_zone) THEN
2727       IF (its .LT. ibe - spec_zone + 1) THEN
2728         max5 = ibe - spec_zone + 1
2729       ELSE
2730         max5 = its
2731       END IF
2732 ! X-end boundary
2733       DO i=max5,itf
2734         b_dist = ibe - i
2735         DO k=kts,ktf
2736           IF (jts .LT. b_dist + jbs + 1) THEN
2737             max6 = b_dist + jbs + 1
2738           ELSE
2739             max6 = jts
2740           END IF
2741           IF (jtf .GT. jbe - b_dist - 1) THEN
2742             min6 = jbe - b_dist - 1
2743           ELSE
2744             min6 = jtf
2745           END IF
2746           DO j=max6,min6
2747             IF (j .LT. jbs + spec_zone) THEN
2748               j_inner = jbs + spec_zone
2749             ELSE
2750               j_inner = j
2751             END IF
2752             IF (j_inner .GT. jbe - spec_zone) THEN
2753               j_inner = jbe - spec_zone
2754             ELSE
2755               j_inner = j_inner
2756             END IF
2757             IF (u(i+1, k, j) .GT. 0.) THEN
2758               fieldd(i, k, j) = fieldd(ibe-spec_zone, k, j_inner)
2759               field(i, k, j) = field(ibe-spec_zone, k, j_inner)
2760             ELSE
2761               fieldd(i, k, j) = 0.0
2762               field(i, k, j) = ccn_conc
2763             END IF
2764           END DO
2765         END DO
2766       END DO
2767     END IF
2768   END IF
2769 END SUBROUTINE G_FLOW_DEP_BDY_QNN
2771 !---------------------------------------------ntu3m---------------------------
2772    SUBROUTINE g_flow_dep_bdy_fixed_inflow(field,g_field,u,v,config_flags, &
2773                                  spec_zone,ids,ide,jds,jde,kds,kde,ims,   &
2774                                  ime,jms,jme,kms,kme,ips,ipe,jps,jpe,kps, &
2775                                  kpe,its,ite,jts,jte,kts,kte)
2776 !-----------------------------------------------------------------------------
2777       IMPLICIT NONE
2779       INTEGER, INTENT(IN) :: ids,ide,jds,jde,kds,kde,ims,ime,jms,jme,kms,  &
2780                              kme,ips,ipe,jps,jpe,kps,kpe,its,ite,jts,jte,  &
2781                              kts,kte,spec_zone
2782       REAL, DIMENSION(ims:ime,kms:kme,jms:jme), INTENT(INOUT) :: g_field,  &
2783                                                                  field
2784       REAL, DIMENSION(ims:ime,kms:kme,jms:jme), INTENT(IN) :: u,v
2785       TYPE(grid_config_rec_type), INTENT(IN) :: config_flags
2786       INTEGER :: i,j,k,ibs,ibe,jbs,jbe,itf,jtf,ktf,i_inner,j_inner,        &
2787                  b_dist,b_limit
2788       LOGICAL :: periodic_x
2790       periodic_x = config_flags%periodic_x
2791       ibs = ids
2792       ibe = ide-1
2793       itf = min(ite,ide-1)
2794       jbs = jds
2795       jbe = jde-1
2796       jtf = min(jte,jde-1)
2797       ktf = kde-1
2799       IF (jts - jbs .lt. spec_zone) THEN
2800 ! Y-start boundary
2801          DO j = jts, min(jtf,jbs+spec_zone-1)
2802             b_dist = j - jbs
2803             b_limit = b_dist
2804             IF (periodic_x) b_limit = 0
2805             DO k = kts, ktf
2806                DO i = max(its,b_limit+ibs), min(itf,ibe-b_limit)
2807                   i_inner = max(i,ibs+spec_zone)
2808                   i_inner = min(i_inner,ibe-spec_zone)
2809                   IF (periodic_x) i_inner = i
2810                   IF (v(i,k,j) .lt. 0.) THEN
2811                      g_field(i,k,j) = g_field(i_inner,k,jbs+spec_zone)
2812                      field(i,k,j) = field(i_inner,k,jbs+spec_zone)
2813                   ELSE
2814                      g_field(i,k,j) = g_field(i,k,j)
2815                      field(i,k,j) = field(i,k,j)
2816                   ENDIF
2817                ENDDO
2818             ENDDO
2819          ENDDO
2820       ENDIF
2821       IF (jbe - jtf .lt. spec_zone) THEN
2822 ! Y-end boundary
2823          DO j = max(jts,jbe-spec_zone+1), jtf
2824             b_dist = jbe - j
2825             b_limit = b_dist
2826             IF (periodic_x) b_limit = 0
2827             DO k = kts, ktf
2828                DO i = max(its,b_limit+ibs), min(itf,ibe-b_limit)
2829                   i_inner = max(i,ibs+spec_zone)
2830                   i_inner = min(i_inner,ibe-spec_zone)
2831                   IF (periodic_x) i_inner = i
2832                   IF (v(i,k,j+1) .gt. 0.) THEN
2833                      g_field(i,k,j) = g_field(i_inner,k,jbe-spec_zone)
2834                      field(i,k,j) = field(i_inner,k,jbe-spec_zone)
2835                   ELSE
2836                      g_field(i,k,j) = g_field(i,k,j)
2837                      field(i,k,j) = field(i,k,j)
2838                   ENDIF
2839                ENDDO
2840             ENDDO
2841          ENDDO
2842       ENDIF
2844    IF (.NOT.periodic_x) THEN
2845       IF (its - ibs .lt. spec_zone) THEN
2846 ! X-start boundary
2847          DO i = its, min(itf,ibs+spec_zone-1)
2848             b_dist = i - ibs
2849             DO k = kts, ktf
2850                DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
2851                   j_inner = max(j,jbs+spec_zone)
2852                   j_inner = min(j_inner,jbe-spec_zone)
2853                   IF (u(i,k,j) .lt. 0.) THEN
2854                      g_field(i,k,j) = g_field(ibs+spec_zone,k,j_inner)
2855                      field(i,k,j) = field(ibs+spec_zone,k,j_inner)
2856                   ELSE
2857                      g_field(i,k,j) = g_field(i,k,j)
2858                      field(i,k,j) = field(i,k,j)
2859                   ENDIF
2860                ENDDO
2861             ENDDO
2862          ENDDO
2863       ENDIF
2865       IF (ibe - itf .lt. spec_zone) THEN
2866 ! X-end boundary
2867          DO i = max(its,ibe-spec_zone+1), itf
2868             b_dist = ibe - i
2869             DO k = kts, ktf
2870                DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
2871                   j_inner = max(j,jbs+spec_zone)
2872                   j_inner = min(j_inner,jbe-spec_zone)
2873                   IF (u(i+1,k,j) .gt. 0.) THEN
2874                      g_field(i,k,j) = g_field(ibe-spec_zone,k,j_inner)
2875                      field(i,k,j) = field(ibe-spec_zone,k,j_inner)
2876                   ELSE
2877                      g_field(i,k,j) = g_field(i,k,j)
2878                      field(i,k,j) = field(i,k,j)
2879                   ENDIF
2880                ENDDO
2881             ENDDO
2882          ENDDO
2883       ENDIF
2884    ENDIF
2886    END SUBROUTINE g_flow_dep_bdy_fixed_inflow
2887 !-----------------------------------------------ntu3m----------------------------
2889 END MODULE g_module_bc