ROMS
Loading...
Searching...
No Matches
tl_u2dbc_im.F
Go to the documentation of this file.
1#include "cppdefs.h"
3#ifdef TANGENT
4!
5!git $Id$
6!================================================== Hernan G. Arango ===
7! Copyright (c) 2002-2025 The ROMS Group Andrew M. Moore !
8! Licensed under a MIT/X style license !
9! See License_ROMS.md !
10!=======================================================================
11! !
12! This subroutine sets tangent linear lateral boundary conditions for !
13! vertically integrated U-velocity. It updates the specified "kout" !
14! index. !
15! !
16! BASIC STATE variables needed: ubar !
17! !
18!=======================================================================
19!
20 implicit none
21!
22 PRIVATE
23 PUBLIC :: tl_u2dbc, tl_u2dbc_tile
24!
25 CONTAINS
26!
27!***********************************************************************
28 SUBROUTINE tl_u2dbc (ng, tile, kout)
29!***********************************************************************
30!
31 USE mod_param
32 USE mod_ocean
33 USE mod_stepping
34!
35! Imported variable declarations.
36!
37 integer, intent(in) :: ng, tile, kout
38!
39! Local variable declarations.
40!
41# include "tile.h"
42!
43 CALL tl_u2dbc_tile (ng, tile, &
44 & lbi, ubi, lbj, ubj, &
45 & imins, imaxs, jmins, jmaxs, &
46 & krhs(ng), kstp(ng), kout, &
47 & ocean(ng) % ubar, &
48 & ocean(ng) % vbar, &
49 & ocean(ng) % zeta, &
50 & ocean(ng) % tl_ubar, &
51 & ocean(ng) % tl_vbar, &
52 & ocean(ng) % tl_zeta)
53
54 RETURN
55 END SUBROUTINE tl_u2dbc
56!
57!***********************************************************************
58 SUBROUTINE tl_u2dbc_tile (ng, tile, &
59 & LBi, UBi, LBj, UBj, &
60 & IminS, ImaxS, JminS, JmaxS, &
61 & krhs, kstp, kout, &
62 & ubar, vbar, zeta, &
63 & tl_ubar, tl_vbar, tl_zeta)
64!***********************************************************************
65!
66 USE mod_param
67 USE mod_boundary
68 USE mod_clima
69 USE mod_forces
70 USE mod_grid
71 USE mod_ncparam
72 USE mod_scalars
73!
74! Imported variable declarations.
75!
76 integer, intent(in) :: ng, tile
77 integer, intent(in) :: lbi, ubi, lbj, ubj
78 integer, intent(in) :: imins, imaxs, jmins, jmaxs
79 integer, intent(in) :: krhs, kstp, kout
80!
81# ifdef ASSUMED_SHAPE
82 real(r8), intent(in) :: ubar(lbi:,lbj:,:)
83 real(r8), intent(in) :: vbar(lbi:,lbj:,:)
84 real(r8), intent(in) :: zeta(lbi:,lbj:,:)
85 real(r8), intent(in) :: tl_vbar(lbi:,lbj:,:)
86 real(r8), intent(in) :: tl_zeta(lbi:,lbj:,:)
87
88 real(r8), intent(inout) :: tl_ubar(lbi:,lbj:,:)
89# else
90 real(r8), intent(in) :: ubar(lbi:ubi,lbj:ubj,:)
91 real(r8), intent(in) :: vbar(lbi:ubi,lbj:ubj,:)
92 real(r8), intent(in) :: zeta(lbi:ubi,lbj:ubj,:)
93 real(r8), intent(in) :: tl_vbar(lbi:ubi,lbj:ubj,:)
94 real(r8), intent(in) :: tl_zeta(lbi:ubi,lbj:ubj,:)
95
96 real(r8), intent(inout) :: tl_ubar(lbi:ubi,lbj:ubj,:)
97# endif
98!
99! Local variable declarations.
100!
101 integer :: imin, imax
102 integer :: i, j, know
103
104 real(r8) :: ce, cx, zx
105 real(r8) :: bry_pgr, bry_cor, bry_str
106 real(r8) :: cff, cff1, cff2, cff3, dt2d
107 real(r8) :: obc_in, obc_out, tau
108# if defined ATM_PRESS && defined PRESS_COMPENSATE
109 real(r8) :: oneatm, fac
110# endif
111
112 real(r8) :: tl_ce, tl_cx, tl_zx
113 real(r8) :: tl_bry_pgr, tl_bry_cor, tl_bry_str, tl_bry_val
114 real(r8) :: tl_cff, tl_cff1, tl_cff2, tl_cff3
115
116 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_grad
117
118# include "set_bounds.h"
119!
120!-----------------------------------------------------------------------
121! Set time-indices
122!-----------------------------------------------------------------------
123!
124 IF (first_2d_step) THEN
125 know=krhs
126 dt2d=dtfast(ng)
127 ELSE IF (predictor_2d_step(ng)) THEN
128 know=krhs
129 dt2d=2.0_r8*dtfast(ng)
130 ELSE
131 know=kstp
132 dt2d=dtfast(ng)
133 END IF
134# if defined ATM_PRESS && defined PRESS_COMPENSATE
135 oneatm=1013.25_r8 ! 1 atm = 1013.25 mb
136 fac=100.0_r8/(g*rho0)
137# endif
138!
139!-----------------------------------------------------------------------
140! Lateral boundary conditions at the western edge.
141!-----------------------------------------------------------------------
142!
143 IF (domain(ng)%Western_Edge(tile)) THEN
144!
145! Western edge, implicit upstream radiation condition.
146!
147 IF (tl_lbc(iwest,isubar,ng)%radiation) THEN
148 IF (iic(ng).ne.0) THEN
149 DO j=jstr,jend+1
150!^ grad(Istr,j)=ubar(Istr,j ,know)- &
151!^ & ubar(Istr,j-1,know)
152!^
153 tl_grad(istr,j)=0.0_r8
154 END DO
155 DO j=jstr,jend
156 IF (lbc_apply(ng)%west(j)) THEN
157# if defined CELERITY_READ && defined FORWARD_READ
158 IF (tl_lbc(iwest,isubar,ng)%nudging) THEN
159 IF (lnudgem2clm(ng)) THEN
160 obc_out=0.5_r8* &
161 & (clima(ng)%M2nudgcof(istr-1,j)+ &
162 & clima(ng)%M2nudgcof(istr ,j))
163 obc_in =obcfac(ng)*obc_out
164 ELSE
165 obc_out=m2obc_out(ng,iwest)
166 obc_in =m2obc_in(ng,iwest)
167 END IF
168 IF (boundary(ng)%ubar_west_Cx(j).lt.0.0_r8) THEN
169 tau=obc_in
170 ELSE
171 tau=obc_out
172 END IF
173 tau=tau*dt2d
174 END IF
175 cx=boundary(ng)%ubar_west_Cx(j)
176# ifdef RADIATION_2D
177 ce=boundary(ng)%ubar_west_Ce(j)
178# else
179 ce=0.0_r8
180# endif
181 cff=boundary(ng)%ubar_west_C2(j)
182# endif
183!^ ubar(Istr,j,kout)=(cff*ubar(Istr ,j,know)+ &
184!^ & Cx *ubar(Istr+1,j,kout)- &
185!^ & MAX(Ce,0.0_r8)*grad(Istr,j )- &
186!^ & MIN(Ce,0.0_r8)*grad(Istr,j+1))/ &
187!^ & (cff+Cx)
188!^
189 tl_ubar(istr,j,kout)=(cff*tl_ubar(istr ,j,know)+ &
190 & cx *tl_ubar(istr+1,j,kout)- &
191 & max(ce,0.0_r8)* &
192 & tl_grad(istr,j )- &
193 & min(ce,0.0_r8)* &
194 & tl_grad(istr,j+1))/ &
195 & (cff+cx)
196
197 IF (tl_lbc(iwest,isubar,ng)%nudging) THEN
198!^ ubar(Istr,j,kout)=ubar(Istr,j,kout)+ &
199!^ & tau*(BOUNDARY(ng)%ubar_west(j)- &
200!^ & ubar(Istr,j,know)) &
201!^
202 tl_ubar(istr,j,kout)=tl_ubar(istr,j,kout)- &
203 & tau*tl_ubar(istr,j,know)
204 END IF
205# ifdef MASKING
206!^ ubar(Istr,j,kout)=ubar(Istr,j,kout)* &
207!^ & GRID(ng)%umask(Istr,j)
208!^
209 tl_ubar(istr,j,kout)=tl_ubar(istr,j,kout)* &
210 & grid(ng)%umask(istr,j)
211# endif
212 END IF
213 END DO
214 END IF
215!
216! Western edge, Flather boundary condition.
217!
218 ELSE IF (tl_lbc(iwest,isubar,ng)%Flather) THEN
219 DO j=jstr,jend
220 IF (lbc_apply(ng)%west(j)) THEN
221# if defined SSH_TIDES_NOT_YET && !defined UV_TIDES_NOT_YET
222 IF (tl_lbc(iwest,isfsur,ng)%acquire) THEN
223 bry_pgr=-g*(zeta(istr,j,know)- &
224 & boundary(ng)%zeta_west(j))* &
225 & 0.5_r8*grid(ng)%pm(istr,j)
226 tl_bry_pgr=-g*tl_zeta(istr,j,know)* &
227 & 0.5_r8*grid(ng)%pm(istr,j)
228# ifdef ADJUST_BOUNDARY
229 IF (lobc(iwest,isubar,ng)) THEN
230 tl_bry_pgr=tl_bry_pgr+ &
231 & g*boundary(ng)%tl_zeta_west(j)* &
232 & 0.5_r8*grid(ng)%pm(istr,j)
233 END IF
234# endif
235 ELSE
236 bry_pgr=-g*(zeta(istr ,j,know)- &
237 & zeta(istr-1,j,know))* &
238 & 0.5_r8*(grid(ng)%pm(istr-1,j)+ &
239 & grid(ng)%pm(istr ,j))
240 tl_bry_pgr=-g*(tl_zeta(istr ,j,know)- &
241 & tl_zeta(istr-1,j,know))* &
242 & 0.5_r8*(grid(ng)%pm(istr-1,j)+ &
243 & grid(ng)%pm(istr ,j))
244 END IF
245# ifdef UV_COR
246 bry_cor=0.125_r8*(vbar(istr-1,j ,know)+ &
247 & vbar(istr-1,j+1,know)+ &
248 & vbar(istr ,j ,know)+ &
249 & vbar(istr ,j+1,know))* &
250 & (grid(ng)%f(istr-1,j)+ &
251 & grid(ng)%f(istr ,j))
252 tl_bry_cor=0.125_r8*(tl_vbar(istr-1,j ,know)+ &
253 & tl_vbar(istr-1,j+1,know)+ &
254 & tl_vbar(istr ,j ,know)+ &
255 & tl_vbar(istr ,j+1,know))* &
256 & (grid(ng)%f(istr-1,j)+ &
257 & grid(ng)%f(istr ,j))
258# else
259 bry_cor=0.0_r8
260 tl_bry_cor=0.0_r8
261# endif
262 cff1=1.0_r8/(0.5_r8*(grid(ng)%h(istr-1,j)+ &
263 & zeta(istr-1,j,know)+ &
264 & grid(ng)%h(istr ,j)+ &
265 & zeta(istr ,j,know)))
266 tl_cff1=-cff1*cff1*0.5_r8*(grid(ng)%tl_h(istr-1,j)+ &
267 & tl_zeta(istr-1,j,know)+ &
268 & grid(ng)%tl_h(istr ,j)+ &
269 & tl_zeta(istr ,j,know))
270 bry_str=cff1*(forces(ng)%sustr(istr,j)- &
271 & forces(ng)%bustr(istr,j))
272 tl_bry_str=tl_cff1*(forces(ng)%sustr(istr,j)- &
273 & forces(ng)%bustr(istr,j))+ &
274 & cff1*(forces(ng)%tl_sustr(istr,j)- &
275 & forces(ng)%tl_bustr(istr,j))
276 cx=1.0_r8/sqrt(g*0.5_r8*(grid(ng)%h(istr-1,j)+ &
277 & zeta(istr-1,j,know)+ &
278 & grid(ng)%h(istr ,j)+ &
279 & zeta(istr ,j,know)))
280 tl_cx=-cx*cx*cx*0.25_r8*g*(grid(ng)%tl_h(istr-1,j)+ &
281 & tl_zeta(istr-1,j,know)+ &
282 & grid(ng)%tl_h(istr ,j)+ &
283 & tl_zeta(istr ,j,know))
284 cff2=grid(ng)%om_u(istr,j)*cx
285 tl_cff2=grid(ng)%om_u(istr,j)*tl_cx
286!^ bry_val=ubar(Istr+1,j,know)+ &
287!^ & cff2*(bry_pgr+ &
288!^ & bry_cor+ &
289!^ & bry_str)
290!^
291 tl_bry_val=tl_ubar(istr+1,j,know)+ &
292 & tl_cff2*(bry_pgr+ &
293 & bry_cor+ &
294 & bry_str)+ &
295 & cff2*(tl_bry_pgr+ &
296 & tl_bry_cor+ &
297 & tl_bry_str)
298# else
299!^ bry_val=BOUNDARY(ng)%ubar_west(j)
300!^
301# ifdef ADJUST_BOUNDARY
302 IF (lobc(iwest,isubar,ng)) THEN
303 tl_bry_val=boundary(ng)%tl_ubar_west(j)
304 ELSE
305 tl_bry_val=0.0_r8
306 END IF
307# else
308 tl_bry_val=0.0_r8
309# endif
310# endif
311 cff=1.0_r8/(0.5_r8*(grid(ng)%h(istr-1,j)+ &
312 & zeta(istr-1,j,know)+ &
313 & grid(ng)%h(istr ,j)+ &
314 & zeta(istr ,j,know)))
315 tl_cff=-cff*cff*(0.5_r8*(grid(ng)%tl_h(istr-1,j)+ &
316 & tl_zeta(istr-1,j,know)+ &
317 & grid(ng)%tl_h(istr ,j)+ &
318 & tl_zeta(istr ,j,know)))
319 cx=sqrt(g*cff)
320 tl_cx=0.5_r8*g*tl_cff/cx
321# if defined ATM_PRESS && defined PRESS_COMPENSATE
322!^ ubar(Istr,j,kout)=bry_val- &
323!^ & Cx*(0.5_r8* &
324!^ & (zeta(Istr-1,j,know)+ &
325!^ & zeta(Istr ,j,know)+ &
326!^ & fac*(FORCES(ng)%Pair(Istr-1,j)+ &
327!^ & FORCES(ng)%Pair(Istr ,j)- &
328!^ & 2.0_r8*OneAtm))- &
329!^ & BOUNDARY(ng)%zeta_west(j))
330!^
331 tl_ubar(istr,j,kout)=tl_bry_val- &
332 & tl_cx* &
333 & (0.5_r8* &
334 & (zeta(istr-1,j,know)+ &
335 & zeta(istr ,j,know)+ &
336 & fac*(forces(ng)%Pair(istr-1,j)+ &
337 & forces(ng)%Pair(istr ,j)- &
338 & 2.0_r8*oneatm))- &
339 & boundary(ng)%zeta_west(j))- &
340 & cx* &
341 & (0.5_r8* &
342 & (tl_zeta(istr-1,j,know)+ &
343 & tl_zeta(istr ,j,know)))
344# else
345!^ ubar(Istr,j,kout)=bry_val- &
346!^ & Cx*(0.5_r8*(zeta(Istr-1,j,know)+ &
347!^ & zeta(Istr ,j,know))- &
348!^ & BOUNDARY(ng)%zeta_west(j))
349!^
350 tl_ubar(istr,j,kout)=tl_bry_val- &
351 & tl_cx* &
352 & (0.5_r8*(zeta(istr-1,j,know)+ &
353 & zeta(istr ,j,know))- &
354 & boundary(ng)%zeta_west(j))- &
355 & cx* &
356 & (0.5_r8*(tl_zeta(istr-1,j,know)+ &
357 & tl_zeta(istr ,j,know)))
358# endif
359# ifdef ADJUST_BOUNDARY
360 IF (lobc(iwest,isubar,ng)) THEN
361 tl_ubar(istr,j,kout)=tl_ubar(istr,j,kout)+ &
362 & cx*boundary(ng)%tl_zeta_west(j)
363 END IF
364# endif
365# ifdef MASKING
366!^ ubar(Istr,j,kout)=ubar(Istr,j,kout)* &
367!^ & GRID(ng)%umask(Istr,j)
368!^
369 tl_ubar(istr,j,kout)=tl_ubar(istr,j,kout)* &
370 & grid(ng)%umask(istr,j)
371# endif
372 END IF
373 END DO
374!
375! Western edge, Shchepetkin boundary condition (Maison et al., 2010).
376!
377 ELSE IF (tl_lbc(iwest,isubar,ng)%Shchepetkin) THEN
378 DO j=jstr,jend
379 IF (lbc_apply(ng)%west(j)) THEN
380# if defined SSH_TIDES_NOT_YET && !defined UV_TIDES_NOT_YET
381 IF (tl_lbc(iwest,isfsur,ng)%acquire) THEN
382 bry_pgr=-g*(zeta(istr,j,know)- &
383 & boundary(ng)%zeta_west(j))* &
384 & 0.5_r8*grid(ng)%pm(istr,j)
385 tl_bry_pgr=-g*tl_zeta(istr,j,know)* &
386 & 0.5_r8*grid(ng)%pm(istr,j)
387# ifdef ADJUST_BOUNDARY
388 IF (lobc(iwest,isubar,ng)) THEN
389 tl_bry_pgr=tl_bry_pgr+ &
390 & g*boundary(ng)%tl_zeta_west(j)* &
391 & 0.5_r8*grid(ng)%pm(istr,j)
392 END IF
393# endif
394 ELSE
395 bry_pgr=-g*(zeta(istr ,j,know)- &
396 & zeta(istr-1,j,know))* &
397 & 0.5_r8*(grid(ng)%pm(istr-1,j)+ &
398 & grid(ng)%pm(istr ,j))
399 tl_bry_pgr=-g*(tl_zeta(istr ,j,know)- &
400 & tl_zeta(istr-1,j,know))* &
401 & 0.5_r8*(grid(ng)%pm(istr-1,j)+ &
402 & grid(ng)%pm(istr ,j))
403 END IF
404# ifdef UV_COR
405 bry_cor=0.125_r8*(vbar(istr-1,j ,know)+ &
406 & vbar(istr-1,j+1,know)+ &
407 & vbar(istr ,j ,know)+ &
408 & vbar(istr ,j+1,know))* &
409 & (grid(ng)%f(istr-1,j)+ &
410 & grid(ng)%f(istr ,j))
411 tl_bry_cor=0.125_r8*(tl_vbar(istr-1,j ,know)+ &
412 & tl_vbar(istr-1,j+1,know)+ &
413 & tl_vbar(istr ,j ,know)+ &
414 & tl_vbar(istr ,j+1,know))* &
415 & (grid(ng)%f(istr-1,j)+ &
416 & grid(ng)%f(istr ,j))
417# else
418 bry_cor=0.0_r8
419 tl_bry_cor=0.0_r8
420# endif
421 cff1=1.0_r8/(0.5_r8*(grid(ng)%h(istr-1,j)+ &
422 & zeta(istr-1,j,know)+ &
423 & grid(ng)%h(istr ,j)+ &
424 & zeta(istr ,j,know)))
425 tl_cff1=-cff1*cff1*0.5_r8*(grid(ng)%tl_h(istr-1,j)+ &
426 & tl_zeta(istr-1,j,know)+ &
427 & grid(ng)%tl_h(istr ,j)+ &
428 & tl_zeta(istr ,j,know))
429 bry_str=cff1*(forces(ng)%sustr(istr,j)- &
430 & forces(ng)%bustr(istr,j))
431 tl_bry_str=tl_cff1*(forces(ng)%sustr(istr,j)- &
432 & forces(ng)%bustr(istr,j))+ &
433 & cff1*(forces(ng)%tl_sustr(istr,j)- &
434 & forces(ng)%tl_bustr(istr,j))
435 cx=1.0_r8/sqrt(g*0.5_r8*(grid(ng)%h(istr-1,j)+ &
436 & zeta(istr-1,j,know)+ &
437 & grid(ng)%h(istr ,j)+ &
438 & zeta(istr ,j,know)))
439 tl_cx=-cx*cx*cx*0.25_r8*g*(grid(ng)%tl_h(istr-1,j)+ &
440 & tl_zeta(istr-1,j,know)+ &
441 & grid(ng)%tl_h(istr ,j)+ &
442 & tl_zeta(istr ,j,know))
443 cff2=grid(ng)%om_u(istr,j)*cx
444 tl_cff2=grid(ng)%om_u(istr,j)*tl_cx
445!^ bry_val=ubar(Istr+1,j,know)+ &
446!^ & cff2*(bry_pgr+ &
447!^ & bry_cor+ &
448!^ & bry_str)
449!^
450 tl_bry_val=tl_ubar(istr+1,j,know)+ &
451 & tl_cff2*(bry_pgr+ &
452 & bry_cor+ &
453 & bry_str)+ &
454 & cff2*(tl_bry_pgr+ &
455 & tl_bry_cor+ &
456 & tl_bry_str)
457# else
458!^ bry_val=BOUNDARY(ng)%ubar_west(j)
459!^
460# ifdef ADJUST_BOUNDARY
461 IF (lobc(iwest,isubar,ng)) THEN
462 tl_bry_val=boundary(ng)%tl_ubar_west(j)
463 ELSE
464 tl_bry_val=0.0_r8
465 END IF
466# else
467 tl_bry_val=0.0_r8
468# endif
469# endif
470# ifdef WET_DRY_NOT_YET
471 cff=0.5_r8*(grid(ng)%h(istr-1,j)+ &
472 & zeta(istr-1,j,know)+ &
473 & grid(ng)%h(istr ,j)+ &
474 & zeta(istr ,j,know))
475 tl_cff=0.5_r8*(grid(ng)%tl_h(istr-1,j)+ &
476 & tl_zeta(istr-1,j,know)+ &
477 & grid(ng)%tl_h(istr ,j)+ &
478 & tl_zeta(istr ,j,know))
479# else
480 cff=0.5_r8*(grid(ng)%h(istr-1,j)+ &
481 & grid(ng)%h(istr ,j))
482 tl_cff=0.5_r8*(grid(ng)%tl_h(istr-1,j)+ &
483 & grid(ng)%tl_h(istr ,j))
484# endif
485 cff1=sqrt(g/cff)
486 tl_cff1=-0.5_r8*cff1*tl_cff/cff
487 cx=dt2d*cff1*cff*0.5_r8*(grid(ng)%pm(istr-1,j)+ &
488 & grid(ng)%pm(istr ,j))
489 tl_cx=dt2d*0.5_r8*(grid(ng)%pm(istr-1,j)+ &
490 & grid(ng)%pm(istr ,j))* &
491 & (cff1*tl_cff+ &
492 & tl_cff1*cff)
493 zx=(0.5_r8+cx)*zeta(istr ,j,know)+ &
494 & (0.5_r8-cx)*zeta(istr-1,j,know)
495 tl_zx=(0.5_r8+cx)*tl_zeta(istr ,j,know)+ &
496 & (0.5_r8-cx)*tl_zeta(istr-1,j,know)+ &
497 & tl_cx*(zeta(istr ,j,know)- &
498 & zeta(istr-1,j,know))
499 IF (cx.gt.co) THEN
500 cff2=(1.0_r8-co/cx)**2
501 tl_cff2=2.0_r8*cff2*co*tl_cx/(cx*cx)
502 cff3=zeta(istr,j,kout)+ &
503 & cx*zeta(istr-1,j,know)- &
504 & (1.0_r8+cx)*zeta(istr,j,know)
505 tl_cff3=tl_zeta(istr,j,kout)+ &
506 & cx*tl_zeta(istr-1,j,know)+ &
507 & tl_cx*(zeta(istr-1,j,know)+ &
508 & zeta(istr ,j,know))- &
509 & (1.0_r8+cx)*tl_zeta(istr,j,know)
510 zx=zx+cff2*cff3
511 tl_zx=tl_zx+cff2*tl_cff3+ &
512 & tl_cff2*cff3
513 END IF
514!^ ubar(Istr,j,kout)=0.5_r8* &
515!^ & ((1.0_r8-Cx)*ubar(Istr,j,know)+ &
516!^ & Cx*ubar(Istr+1,j,know)+ &
517!^ & bry_val- &
518!^ & cff1*(Zx-BOUNDARY(ng)%zeta_west(j)))
519!^
520 tl_ubar(istr,j,kout)=0.5_r8* &
521 & ((1.0_r8-cx)* &
522 & tl_ubar(istr,j,know)- &
523 & tl_cx*(ubar(istr ,j,know)- &
524 & ubar(istr+1,j,know))+ &
525 & cx*tl_ubar(istr+1,j,know)+ &
526 & tl_bry_val- &
527 & tl_cff1* &
528 & (zx-boundary(ng)%zeta_west(j))- &
529 & cff1*tl_zx)
530# ifdef ADJUST_BOUNDARY
531 IF (lobc(iwest,isubar,ng)) THEN
532 tl_ubar(istr,j,kout)=tl_ubar(istr,j,kout)+ &
533 & 0.5_r8*cff1* &
534 & boundary(ng)%tl_zeta_west(j)
535 END IF
536# endif
537# ifdef MASKING
538!^ ubar(Istr,j,kout)=ubar(Istr,j,kout)* &
539!^ & GRID(ng)%umask(Istr,j)
540!^
541 tl_ubar(istr,j,kout)=tl_ubar(istr,j,kout)* &
542 & grid(ng)%umask(istr,j)
543# endif
544 END IF
545 END DO
546!
547! Western edge, clamped boundary condition.
548!
549 ELSE IF (tl_lbc(iwest,isubar,ng)%clamped) THEN
550 DO j=jstr,jend
551 IF (lbc_apply(ng)%west(j)) THEN
552!^ ubar(Istr,j,kout)=BOUNDARY(ng)%ubar_west(j)
553!^
554# ifdef ADJUST_BOUNDARY
555 IF (lobc(iwest,isubar,ng)) THEN
556 tl_ubar(istr,j,kout)=boundary(ng)%tl_ubar_west(j)
557 ELSE
558 tl_ubar(istr,j,kout)=0.0_r8
559 END IF
560# else
561 tl_ubar(istr,j,kout)=0.0_r8
562# endif
563# ifdef MASKING
564!^ ubar(Istr,j,kout)=ubar(Istr,j,kout)* &
565!^ & GRID(ng)%umask(Istr,j)
566!^
567 tl_ubar(istr,j,kout)=tl_ubar(istr,j,kout)* &
568 & grid(ng)%umask(istr,j)
569# endif
570 END IF
571 END DO
572!
573! Western edge, gradient boundary condition.
574!
575 ELSE IF (tl_lbc(iwest,isubar,ng)%gradient) THEN
576 DO j=jstr,jend
577 IF (lbc_apply(ng)%west(j)) THEN
578!^ ubar(Istr,j,kout)=ubar(Istr+1,j,kout)
579!^
580 tl_ubar(istr,j,kout)=tl_ubar(istr+1,j,kout)
581# ifdef MASKING
582!^ ubar(Istr,j,kout)=ubar(Istr,j,kout)* &
583!^ & GRID(ng)%umask(Istr,j)
584!^
585 tl_ubar(istr,j,kout)=tl_ubar(istr,j,kout)* &
586 & grid(ng)%umask(istr,j)
587# endif
588 END IF
589 END DO
590!
591! Western edge, reduced-physics boundary condition.
592!
593 ELSE IF (tl_lbc(iwest,isubar,ng)%reduced) THEN
594 DO j=jstr,jend
595 IF (lbc_apply(ng)%west(j)) THEN
596 IF (tl_lbc(iwest,isfsur,ng)%acquire) THEN
597!^ bry_pgr=-g*(zeta(Istr,j,know)- &
598!^ & BOUNDARY(ng)%zeta_west(j))* &
599!^ & 0.5_r8*GRID(ng)%pm(Istr,j)
600!^
601 tl_bry_pgr=-g*tl_zeta(istr,j,know)* &
602 & 0.5_r8*grid(ng)%pm(istr,j)
603# ifdef ADJUST_BOUNDARY
604 IF (lobc(iwest,isubar,ng)) THEN
605 tl_bry_pgr=tl_bry_pgr+ &
606 & g*boundary(ng)%tl_zeta_west(j)* &
607 & 0.5_r8*grid(ng)%pm(istr,j)
608 END IF
609# endif
610 ELSE
611!^ bry_pgr=-g*(zeta(Istr,j,know)- &
612!^ & zeta(Istr-1,j,know))* &
613!^ & 0.5_r8*(GRID(ng)%pm(Istr-1,j)+ &
614!^ & GRID(ng)%pm(Istr ,j))
615!^
616 tl_bry_pgr=-g*(tl_zeta(istr ,j,know)- &
617 & tl_zeta(istr-1,j,know))* &
618 & 0.5_r8*(grid(ng)%pm(istr-1,j)+ &
619 & grid(ng)%pm(istr ,j))
620 END IF
621# ifdef UV_COR
622!^ bry_cor=0.125_r8*(vbar(Istr-1,j ,know)+ &
623!^ & vbar(Istr-1,j+1,know)+ &
624!^ & vbar(Istr ,j ,know)+ &
625!^ & vbar(Istr ,j+1,know))* &
626!^ & (GRID(ng)%f(Istr-1,j)+ &
627!^ & GRID(ng)%f(Istr ,j))
628!^
629 tl_bry_cor=0.125_r8*(tl_vbar(istr-1,j ,know)+ &
630 & tl_vbar(istr-1,j+1,know)+ &
631 & tl_vbar(istr ,j ,know)+ &
632 & tl_vbar(istr ,j+1,know))* &
633 & (grid(ng)%f(istr-1,j)+ &
634 & grid(ng)%f(istr ,j))
635# else
636!^ bry_cor=0.0_r8
637!^
638 tl_bry_cor=0.0_r8
639# endif
640 cff=1.0_r8/(0.5_r8*(grid(ng)%h(istr-1,j)+ &
641 & zeta(istr-1,j,know)+ &
642 & grid(ng)%h(istr ,j)+ &
643 & zeta(istr ,j,know)))
644 tl_cff=-cff*cff*0.5_r8*(grid(ng)%tl_h(istr-1,j)+ &
645 & tl_zeta(istr-1,j,know)+ &
646 & grid(ng)%tl_h(istr ,j)+ &
647 & tl_zeta(istr ,j,know))
648!^ bry_str=cff*(FORCES(ng)%sustr(Istr,j)- &
649!^ & FORCES(ng)%bustr(Istr,j))
650!^
651 tl_bry_str=tl_cff*(forces(ng)%sustr(istr,j)- &
652 & forces(ng)%bustr(istr,j))+ &
653 & cff*(forces(ng)%tl_sustr(istr,j)- &
654 & forces(ng)%tl_bustr(istr,j))
655!^ ubar(Istr,j,kout)=ubar(Istr,j,know)+ &
656!^ & dt2d*(bry_pgr+ &
657!^ & bry_cor+ &
658!^ & bry_str)
659!^
660 tl_ubar(istr,j,kout)=tl_ubar(istr,j,know)+ &
661 & dt2d*(tl_bry_pgr+ &
662 & tl_bry_cor+ &
663 & tl_bry_str)
664# ifdef MASKING
665!^ ubar(Istr,j,kout)=ubar(Istr,j,kout)* &
666!^ & GRID(ng)%umask(Istr,j)
667!^
668 tl_ubar(istr,j,kout)=tl_ubar(istr,j,kout)* &
669 & grid(ng)%umask(istr,j)
670# endif
671 END IF
672 END DO
673!
674! Western edge, closed boundary condition.
675!
676 ELSE IF (tl_lbc(iwest,isubar,ng)%closed) THEN
677 DO j=jstr,jend
678 IF (lbc_apply(ng)%west(j)) THEN
679!^ ubar(Istr,j,kout)=0.0_r8
680!^
681 tl_ubar(istr,j,kout)=0.0_r8
682 END IF
683 END DO
684 END IF
685 END IF
686!
687!-----------------------------------------------------------------------
688! Lateral boundary conditions at the eastern edge.
689!-----------------------------------------------------------------------
690!
691 IF (domain(ng)%Eastern_Edge(tile)) THEN
692!
693! Eastern edge, implicit upstream radiation condition.
694!
695 IF (tl_lbc(ieast,isubar,ng)%radiation) THEN
696 IF (iic(ng).ne.0) THEN
697 DO j=jstr,jend+1
698!^ grad(Iend+1,j)=ubar(Iend+1,j ,know)- &
699!^ & ubar(Iend+1,j-1,know)
700!^
701 tl_grad(iend+1,j)=0.0_r8
702 END DO
703 DO j=jstr,jend
704 IF (lbc_apply(ng)%east(j)) THEN
705# if defined CELERITY_READ && defined FORWARD_READ
706 IF (tl_lbc(ieast,isubar,ng)%nudging) THEN
707 IF (lnudgem2clm(ng)) THEN
708 obc_out=0.5_r8* &
709 & (clima(ng)%M2nudgcof(iend ,j)+ &
710 & clima(ng)%M2nudgcof(iend+1,j))
711 obc_in =obcfac(ng)*obc_out
712 ELSE
713 obc_out=m2obc_out(ng,ieast)
714 obc_in =m2obc_in(ng,ieast)
715 END IF
716 IF (boundary(ng)%ubar_east_Cx(j).lt.0.0_r8) THEN
717 tau=obc_in
718 ELSE
719 tau=obc_out
720 END IF
721 tau=tau*dt2d
722 END IF
723 cx=boundary(ng)%ubar_east_Cx(j)
724# ifdef RADIATION_2D
725 ce=boundary(ng)%ubar_east_Ce(j)
726# else
727 ce=0.0_r8
728# endif
729 cff=boundary(ng)%ubar_east_C2(j)
730# endif
731!^ ubar(Iend+1,j,kout)=(cff*ubar(Iend+1,j,know)+ &
732!^ & Cx *ubar(Iend ,j,kout)- &
733!^ & MAX(Ce,0.0_r8)*grad(Iend+1,j )- &
734!^ & MIN(Ce,0.0_r8)*grad(Iend+1,j+1))/ &
735!^ & (cff+Cx)
736!^
737 tl_ubar(iend+1,j,kout)=(cff*tl_ubar(iend+1,j,know)+ &
738 & cx *tl_ubar(iend ,j,kout)- &
739 & max(ce,0.0_r8)* &
740 & tl_grad(iend+1,j )- &
741 & min(ce,0.0_r8)* &
742 & tl_grad(iend+1,j+1))/ &
743 & (cff+cx)
744
745 IF (tl_lbc(ieast,isubar,ng)%nudging) THEN
746!^ ubar(Iend+1,j,kout)=ubar(Iend+1,j,kout)+ &
747!^ & tau*(BOUNDARY(ng)%ubar_east(j)- &
748!^ & ubar(Iend+1,j,know))
749!^
750 tl_ubar(iend+1,j,kout)=tl_ubar(iend+1,j,kout)- &
751 & tau*tl_ubar(iend+1,j,know)
752 END IF
753# ifdef MASKING
754!^ ubar(Iend+1,j,kout)=ubar(Iend+1,j,kout)* &
755!^ & GRID(ng)%umask(Iend+1,j)
756!^
757 tl_ubar(iend+1,j,kout)=tl_ubar(iend+1,j,kout)* &
758 & grid(ng)%umask(iend+1,j)
759# endif
760 END IF
761 END DO
762 END IF
763!
764! Eastern edge, Flather boundary condition.
765!
766 ELSE IF (tl_lbc(ieast,isubar,ng)%Flather) THEN
767 DO j=jstr,jend
768 IF (lbc_apply(ng)%east(j)) THEN
769# if defined SSH_TIDES_NOT_YET && !defined UV_TIDES_NOT_YET
770 IF (tl_lbc(ieast,isfsur,ng)%acquire) THEN
771 bry_pgr=-g*(boundary(ng)%zeta_east(j)- &
772 & zeta(iend,j,know))* &
773 & 0.5_r8*grid(ng)%pm(iend,j)
774 tl_bry_pgr=g*tl_zeta(iend,j,know)* &
775 & 0.5_r8*grid(ng)%pm(iend,j)
776# ifdef ADJUST_BOUNDARY
777 IF (lobc(ieast,isubar,ng)) THEN
778 tl_bry_pgr=tl_bry_pgr- &
779 & g*boundary(ng)%tl_zeta_east(j)* &
780 & 0.5_r8*grid(ng)%pm(iend,j)
781 END IF
782# endif
783 ELSE
784 bry_pgr=-g*(zeta(iend+1,j,know)- &
785 & zeta(iend ,j,know))* &
786 & 0.5_r8*(grid(ng)%pm(iend ,j)+ &
787 & grid(ng)%pm(iend+1,j))
788 tl_bry_pgr=-g*(tl_zeta(iend+1,j,know)- &
789 & tl_zeta(iend ,j,know))* &
790 & 0.5_r8*(grid(ng)%pm(iend ,j)+ &
791 & grid(ng)%pm(iend+1,j))
792 END IF
793# ifdef UV_COR
794 bry_cor=0.125_r8*(vbar(iend ,j ,know)+ &
795 & vbar(iend ,j+1,know)+ &
796 & vbar(iend+1,j ,know)+ &
797 & vbar(iend+1,j+1,know))* &
798 & (grid(ng)%f(iend ,j)+ &
799 & grid(ng)%f(iend+1,j))
800 tl_bry_cor=0.125_r8*(tl_vbar(iend ,j ,know)+ &
801 & tl_vbar(iend ,j+1,know)+ &
802 & tl_vbar(iend+1,j ,know)+ &
803 & tl_vbar(iend+1,j+1,know))* &
804 & (grid(ng)%f(iend ,j)+ &
805 & grid(ng)%f(iend+1,j))
806# else
807 bry_cor=0.0_r8
808 tl_bry_cor=0.0_r8
809# endif
810 cff1=1.0_r8/(0.5_r8*(grid(ng)%h(iend ,j)+ &
811 & zeta(iend ,j,know)+ &
812 & grid(ng)%h(iend+1,j)+ &
813 & zeta(iend+1,j,know)))
814 tl_cff1=-cff1*cff1*0.5_r8*(grid(ng)%tl_h(iend ,j)+ &
815 & tl_zeta(iend ,j,know)+ &
816 & grid(ng)%tl_h(iend+1,j)+ &
817 & tl_zeta(iend+1,j,know))
818 bry_str=cff1*(forces(ng)%sustr(iend+1,j)- &
819 & forces(ng)%bustr(iend+1,j))
820 tl_bry_str=tl_cff1*(forces(ng)%sustr(iend+1,j)- &
821 & forces(ng)%bustr(iend+1,j))+ &
822 & cff1*(forces(ng)%tl_sustr(iend+1,j)- &
823 & forces(ng)%tl_bustr(iend+1,j))
824 cx=1.0_r8/sqrt(g*0.5_r8*(grid(ng)%h(iend+1,j)+ &
825 & zeta(iend+1,j,know)+ &
826 & grid(ng)%h(iend ,j)+ &
827 & zeta(iend ,j,know)))
828 tl_cx=-cx*cx*cx*0.25_r8*g*(grid(ng)%tl_h(iend+1,j)+ &
829 & tl_zeta(iend+1,j,know)+ &
830 & grid(ng)%tl_h(iend ,j)+ &
831 & tl_zeta(iend ,j,know))
832 cff2=grid(ng)%om_u(iend+1,j)*cx
833 tl_cff2=grid(ng)%om_u(iend+1,j)*tl_cx
834!^ bry_val=ubar(Iend,j,know)+ &
835!^ & cff2*(bry_pgr+ &
836!^ & bry_cor+ &
837!^ & bry_str)
838!^
839 tl_bry_val=tl_ubar(iend,j,know)+ &
840 & tl_cff2*(bry_pgr+ &
841 & bry_cor+ &
842 & bry_str)+ &
843 & cff2*(tl_bry_pgr+ &
844 & tl_bry_cor+ &
845 & tl_bry_str)
846# else
847!^ bry_val=BOUNDARY(ng)%ubar_east(j)
848!^
849# ifdef ADJUST_BOUNDARY
850 IF (lobc(ieast,isubar,ng)) THEN
851 tl_bry_val=boundary(ng)%tl_ubar_east(j)
852 ELSE
853 tl_bry_val=0.0_r8
854 END IF
855# else
856 tl_bry_val=0.0_r8
857# endif
858# endif
859 cff=1.0_r8/(0.5_r8*(grid(ng)%h(iend ,j)+ &
860 & zeta(iend ,j,know)+ &
861 & grid(ng)%h(iend+1,j)+ &
862 & zeta(iend+1,j,know)))
863 tl_cff=-cff*cff*(0.5_r8*(grid(ng)%tl_h(iend ,j)+ &
864 & tl_zeta(iend ,j,know)+ &
865 & grid(ng)%tl_h(iend+1,j)+ &
866 & tl_zeta(iend+1,j,know)))
867 cx=sqrt(g*cff)
868 tl_cx=0.5_r8*g*tl_cff/cx
869# if defined ATM_PRESS && defined PRESS_COMPENSATE
870!^ ubar(Iend+1,j,kout)=bry_val+ &
871!^ & Cx*(0.5_r8* &
872!^ & (zeta(Iend ,j,know)+ &
873!^ & zeta(Iend+1,j,know)+ &
874!^ & fac*(FORCES(ng)%Pair(Iend ,j)+ &
875!^ & FORCES(ng)%Pair(Iend+1,j)- &
876!^ & 2.0_r8*OneAtm))- &
877!^ & BOUNDARY(ng)%zeta_east(j))
878!^
879 tl_ubar(iend+1,j,kout)=tl_bry_val+ &
880 & tl_cx* &
881 & (0.5_r8* &
882 & (zeta(iend ,j,know)+ &
883 & zeta(iend+1,j,know)+ &
884 & fac*(forces(ng)%Pair(iend ,j)+ &
885 & forces(ng)%Pair(iend+1,j)- &
886 & 2.0_r8*oneatm))- &
887 & boundary(ng)%zeta_east(j))+ &
888 & cx* &
889 & (0.5_r8*(tl_zeta(iend ,j,know)+ &
890 & tl_zeta(iend+1,j,know)))
891# else
892!^ ubar(Iend+1,j,kout)=bry_val+ &
893!^ & Cx*(0.5_r8*(zeta(Iend ,j,know)+ &
894!^ & zeta(Iend+1,j,know))- &
895!^ & BOUNDARY(ng)%zeta_east(j))
896!^
897 tl_ubar(iend+1,j,kout)=tl_bry_val+ &
898 & tl_cx* &
899 & (0.5_r8*(zeta(iend ,j,know)+ &
900 & zeta(iend+1,j,know))- &
901 & boundary(ng)%zeta_east(j))+ &
902 & cx* &
903 & (0.5_r8*(tl_zeta(iend ,j,know)+ &
904 & tl_zeta(iend+1,j,know)))
905# endif
906# ifdef ADJUST_BOUNDARY
907 IF (lobc(ieast,isubar,ng)) THEN
908 tl_ubar(iend+1,j,kout)=tl_ubar(iend+1,j,kout)- &
909 & cx*boundary(ng)%tl_zeta_east(j)
910 END IF
911# endif
912# ifdef MASKING
913!^ & ubar(Iend+1,j,kout)=ubar(Iend+1,j,kout)* &
914!^ & GRID(ng)%umask(Iend+1,j)
915!^
916 tl_ubar(iend+1,j,kout)=tl_ubar(iend+1,j,kout)* &
917 & grid(ng)%umask(iend+1,j)
918# endif
919 END IF
920 END DO
921!
922! Eastern edge, Shchepetkin boundary condition (Maison et al., 2010).
923!
924 ELSE IF (tl_lbc(ieast,isubar,ng)%Shchepetkin) THEN
925 DO j=jstr,jend
926 IF (lbc_apply(ng)%east(j)) THEN
927# if defined SSH_TIDES_NOT_YET && !defined UV_TIDES_NOT_YET
928 IF (tl_lbc(ieast,isfsur,ng)%acquire) THEN
929 bry_pgr=-g*(boundary(ng)%zeta_east(j)- &
930 & zeta(iend,j,know))* &
931 & 0.5_r8*grid(ng)%pm(iend,j)
932 tl_bry_pgr=g*tl_zeta(iend,j,know)* &
933 & 0.5_r8*grid(ng)%pm(iend,j)
934# ifdef ADJUST_BOUNDARY
935 IF (lobc(ieast,isubar,ng)) THEN
936 tl_bry_pgr=tl_bry_pgr- &
937 & g*boundary(ng)%tl_zeta_east(j)* &
938 & 0.5_r8*grid(ng)%pm(iend,j)
939 END IF
940# endif
941 ELSE
942 bry_pgr=-g*(zeta(iend+1,j,know)- &
943 & zeta(iend ,j,know))* &
944 & 0.5_r8*(grid(ng)%pm(iend ,j)+ &
945 & grid(ng)%pm(iend+1,j))
946 tl_bry_pgr=-g*(tl_zeta(iend+1,j,know)- &
947 & tl_zeta(iend ,j,know))* &
948 & 0.5_r8*(grid(ng)%pm(iend ,j)+ &
949 & grid(ng)%pm(iend+1,j))
950 END IF
951# ifdef UV_COR
952 bry_cor=0.125_r8*(vbar(iend ,j ,know)+ &
953 & vbar(iend ,j+1,know)+ &
954 & vbar(iend+1,j ,know)+ &
955 & vbar(iend+1,j+1,know))* &
956 & (grid(ng)%f(iend ,j)+ &
957 & grid(ng)%f(iend+1,j))
958 tl_bry_cor=0.125_r8*(tl_vbar(iend ,j ,know)+ &
959 & tl_vbar(iend ,j+1,know)+ &
960 & tl_vbar(iend+1,j ,know)+ &
961 & tl_vbar(iend+1,j+1,know))* &
962 & (grid(ng)%f(iend ,j)+ &
963 & grid(ng)%f(iend+1,j))
964# else
965 bry_cor=0.0_r8
966 tl_bry_cor=0.0_r8
967# endif
968 cff1=1.0_r8/(0.5_r8*(grid(ng)%h(iend ,j)+ &
969 & zeta(iend ,j,know)+ &
970 & grid(ng)%h(iend+1,j)+ &
971 & zeta(iend+1,j,know)))
972 tl_cff1=-cff1*cff1*0.5_r8*(grid(ng)%tl_h(iend ,j)+ &
973 & tl_zeta(iend ,j,know)+ &
974 & grid(ng)%tl_h(iend+1,j)+ &
975 & tl_zeta(iend+1,j,know))
976 bry_str=cff1*(forces(ng)%sustr(iend+1,j)- &
977 & forces(ng)%bustr(iend+1,j))
978 tl_bry_str=tl_cff1*(forces(ng)%sustr(iend+1,j)- &
979 & forces(ng)%bustr(iend+1,j))+ &
980 & cff1*(forces(ng)%tl_sustr(iend+1,j)- &
981 & forces(ng)%tl_bustr(iend+1,j))
982 cx=1.0_r8/sqrt(g*0.5_r8*(grid(ng)%h(iend+1,j)+ &
983 & zeta(iend+1,j,know)+ &
984 & grid(ng)%h(iend ,j)+ &
985 & zeta(iend ,j,know)))
986 tl_cx=-cx*cx*cx*0.25_r8*g*(grid(ng)%tl_h(iend+1,j)+ &
987 & tl_zeta(iend+1,j,know)+ &
988 & grid(ng)%tl_h(iend ,j)+ &
989 & tl_zeta(iend ,j,know))
990 cff2=grid(ng)%om_u(iend+1,j)*cx
991 tl_cff2=grid(ng)%om_u(iend+1,j)*tl_cx
992!^ bry_val=ubar(Iend,j,know)+ &
993!^ & cff2*(bry_pgr+ &
994!^ & bry_cor+ &
995!^ & bry_str)
996!^
997 tl_bry_val=tl_ubar(iend,j,know)+ &
998 & tl_cff2*(bry_pgr+ &
999 & bry_cor+ &
1000 & bry_str)+ &
1001 & cff2*(tl_bry_pgr+ &
1002 & tl_bry_cor+ &
1003 & tl_bry_str)
1004# else
1005!^ bry_val=BOUNDARY(ng)%ubar_east(j)
1006!^
1007# ifdef ADJUST_BOUNDARY
1008 IF (lobc(ieast,isubar,ng)) THEN
1009 tl_bry_val=boundary(ng)%tl_ubar_east(j)
1010 ELSE
1011 tl_bry_val=0.0_r8
1012 END IF
1013# else
1014 tl_bry_val=0.0_r8
1015# endif
1016# endif
1017# ifdef WET_DRY_NOT_YET
1018 cff=0.5_r8*(grid(ng)%h(iend ,j)+ &
1019 & zeta(iend ,j,know)+ &
1020 & grid(ng)%h(iend+1,j)+ &
1021 & zeta(iend+1,j,know))
1022 tl_cff=0.5_r8*(grid(ng)%tl_h(iend ,j)+ &
1023 & tl_zeta(iend ,j,know)+ &
1024 & grid(ng)%tl_h(iend+1,j)+ &
1025 & tl_zeta(iend+1,j,know))
1026# else
1027 cff=0.5_r8*(grid(ng)%h(iend ,j)+ &
1028 & grid(ng)%h(iend+1,j))
1029 tl_cff=0.5_r8*(grid(ng)%tl_h(iend ,j)+ &
1030 & grid(ng)%tl_h(iend+1,j))
1031# endif
1032 cff1=sqrt(g/cff)
1033 tl_cff1=-0.5_r8*cff1*tl_cff/cff
1034 cx=dt2d*cff1*cff*0.5_r8*(grid(ng)%pm(iend ,j)+ &
1035 & grid(ng)%pm(iend+1,j))
1036 tl_cx=dt2d*0.5_r8*(grid(ng)%pm(iend ,j)+ &
1037 & grid(ng)%pm(iend+1,j))* &
1038 & (cff1*tl_cff+ &
1039 & tl_cff1*cff)
1040 zx=(0.5_r8+cx)*zeta(iend ,j,know)+ &
1041 & (0.5_r8-cx)*zeta(iend+1,j,know)
1042 tl_zx=(0.5_r8+cx)*tl_zeta(iend ,j,know)+ &
1043 & (0.5_r8-cx)*tl_zeta(iend+1,j,know)+ &
1044 & tl_cx*(zeta(iend ,j,know)- &
1045 & zeta(iend+1,j,know))
1046 IF (cx.gt.co) THEN
1047 cff2=(1.0_r8-co/cx)**2
1048 tl_cff2=2.0_r8*cff2*co*tl_cx/(cx*cx)
1049 cff3=zeta(iend,j,kout)+ &
1050 & cx*zeta(iend+1,j,know)- &
1051 & (1.0_r8+cx)*zeta(iend,j,know)
1052 tl_cff3=tl_zeta(iend,j,kout)+ &
1053 & cx*tl_zeta(iend+1,j,know)+ &
1054 & tl_cx*(zeta(iend ,j,know)+ &
1055 & zeta(iend+1,j,know))- &
1056 & (1.0_r8+cx)*tl_zeta(iend,j,know)
1057 zx=zx+cff2*cff3
1058 tl_zx=tl_zx+cff2*tl_cff3+ &
1059 & tl_cff2*cff3
1060 END IF
1061!^ ubar(Iend+1,j,kout)=0.5_r8* &
1062!^ & ((1.0_r8-Cx)*ubar(Iend+1,j,know)+ &
1063!^ & Cx*ubar(Iend,j,know)+ &
1064!^ & bry_val+ &
1065!^ & cff1*(Zx-BOUNDARY(ng)%zeta_east(j)))
1066!^
1067 tl_ubar(iend+1,j,kout)=0.5_r8* &
1068 & ((1.0_r8-cx)* &
1069 & tl_ubar(iend+1,j,know)+ &
1070 & tl_cx*(ubar(iend ,j,know)- &
1071 & ubar(iend+1,j,know))+ &
1072 & cx*tl_ubar(iend,j,know)+ &
1073 & tl_bry_val+ &
1074 & tl_cff1* &
1075 & (zx-boundary(ng)%zeta_east(j))- &
1076 & cff1*tl_zx)
1077# ifdef ADJUST_BOUNDARY
1078 IF (lobc(ieast,isubar,ng)) THEN
1079 tl_ubar(iend+1,j,kout)=tl_ubar(iend+1,j,kout)- &
1080 & 0.5_r8*cff1* &
1081 & boundary(ng)%tl_zeta_east(j)
1082 END IF
1083# endif
1084# ifdef MASKING
1085!^ & ubar(Iend+1,j,kout)=ubar(Iend+1,j,kout)* &
1086!^ & GRID(ng)%umask(Iend+1,j)
1087!^
1088 tl_ubar(iend+1,j,kout)=tl_ubar(iend+1,j,kout)* &
1089 & grid(ng)%umask(iend+1,j)
1090# endif
1091 END IF
1092 END DO
1093!
1094! Eastern edge, clamped boundary condition.
1095!
1096 ELSE IF (tl_lbc(ieast,isubar,ng)%clamped) THEN
1097 DO j=jstr,jend
1098 IF (lbc_apply(ng)%east(j)) THEN
1099!^ ubar(Iend+1,j,kout)=BOUNDARY(ng)%ubar_east(j)
1100!^
1101# ifdef ADJUST_BOUNDARY
1102 IF (lobc(ieast,isubar,ng)) THEN
1103 tl_ubar(iend+1,j,kout)=boundary(ng)%tl_ubar_east(j)
1104 ELSE
1105 tl_ubar(iend+1,j,kout)=0.0_r8
1106 END IF
1107# else
1108 tl_ubar(iend+1,j,kout)=0.0_r8
1109# endif
1110# ifdef MASKING
1111!^ ubar(Iend+1,j,kout)=ubar(Iend+1,j,kout)* &
1112!^ & GRID(ng)%umask(Iend+1,j)
1113!^
1114 tl_ubar(iend+1,j,kout)=tl_ubar(iend+1,j,kout)* &
1115 & grid(ng)%umask(iend+1,j)
1116# endif
1117 END IF
1118 END DO
1119!
1120! Eastern edge, gradient boundary condition.
1121!
1122 ELSE IF (tl_lbc(ieast,isubar,ng)%gradient) THEN
1123 DO j=jstr,jend
1124 IF (lbc_apply(ng)%east(j)) THEN
1125!^ ubar(Iend+1,j,kout)=ubar(Iend,j,kout)
1126!^
1127 tl_ubar(iend+1,j,kout)=tl_ubar(iend,j,kout)
1128# ifdef MASKING
1129!^ ubar(Iend+1,j,kout)=ubar(Iend+1,j,kout)* &
1130!^ & GRID(ng)%umask(Iend+1,j)
1131!^
1132 tl_ubar(iend+1,j,kout)=tl_ubar(iend+1,j,kout)* &
1133 & grid(ng)%umask(iend+1,j)
1134# endif
1135 END IF
1136 END DO
1137!
1138! Eastern edge, reduced-physics boundary condition.
1139!
1140 ELSE IF (tl_lbc(ieast,isubar,ng)%reduced) THEN
1141 DO j=jstr,jend
1142 IF (lbc_apply(ng)%east(j)) THEN
1143 IF (tl_lbc(ieast,isfsur,ng)%acquire) THEN
1144!^ bry_pgr=-g*(BOUNDARY(ng)%zeta_east(j)- &
1145!^ & zeta(Iend,j,know))* &
1146!^ & 0.5_r8*GRID(ng)%pm(Iend,j)
1147!^
1148 tl_bry_pgr=g*tl_zeta(iend,j,know)* &
1149 & 0.5_r8*grid(ng)%pm(iend,j)
1150# ifdef ADJUST_BOUNDARY
1151 IF (lobc(ieast,isubar,ng)) THEN
1152 tl_bry_pgr=tl_bry_pgr- &
1153 & g*boundary(ng)%tl_zeta_east(j)* &
1154 & 0.5_r8*grid(ng)%pm(iend,j)
1155 END IF
1156# endif
1157 ELSE
1158!^ bry_pgr=-g*(zeta(Iend+1,j,know)- &
1159!^ & zeta(Iend ,j,know))* &
1160!^ & 0.5_r8*(GRID(ng)%pm(Iend ,j)+ &
1161!^ & GRID(ng)%pm(Iend+1,j))
1162!^
1163 tl_bry_pgr=-g*(tl_zeta(iend+1,j,know)- &
1164 & tl_zeta(iend ,j,know))* &
1165 & 0.5_r8*(grid(ng)%pm(iend ,j)+ &
1166 & grid(ng)%pm(iend+1,j))
1167 END IF
1168# ifdef UV_COR
1169!^ bry_cor=0.125_r8*(vbar(Iend ,j ,know)+ &
1170!^ & vbar(Iend ,j+1,know)+ &
1171!^ & vbar(Iend+1,j ,know)+ &
1172!^ & vbar(Iend+1,j+1,know))* &
1173!^ & (GRID(ng)%f(Iend ,j)+ &
1174!^ & GRID(ng)%f(Iend+1,j))
1175!^
1176 tl_bry_cor=0.125_r8*(tl_vbar(iend, j ,know)+ &
1177 & tl_vbar(iend ,j+1,know)+ &
1178 & tl_vbar(iend+1,j ,know)+ &
1179 & tl_vbar(iend+1,j+1,know))* &
1180 & (grid(ng)%f(iend ,j)+ &
1181 & grid(ng)%f(iend+1,j))
1182# else
1183!^ bry_cor=0.0_r8
1184!^
1185 tl_bry_cor=0.0_r8
1186# endif
1187 cff=1.0_r8/(0.5_r8*(grid(ng)%h(iend ,j)+ &
1188 & zeta(iend ,j,know)+ &
1189 & grid(ng)%h(iend+1,j)+ &
1190 & zeta(iend+1,j,know)))
1191 tl_cff=-cff*cff*0.5_r8*(grid(ng)%tl_h(iend ,j)+ &
1192 & tl_zeta(iend ,j,know)+ &
1193 & grid(ng)%tl_h(iend+1,j)+ &
1194 & tl_zeta(iend+1,j,know))
1195!^ bry_str=cff*(FORCES(ng)%sustr(Iend+1,j)- &
1196!^ & FORCES(ng)%bustr(Iend+1,j))
1197!^
1198 tl_bry_str=tl_cff*(forces(ng)%sustr(iend+1,j)- &
1199 & forces(ng)%bustr(iend+1,j))+ &
1200 & cff*(forces(ng)%tl_sustr(iend+1,j)- &
1201 & forces(ng)%tl_bustr(iend+1,j))
1202!^ ubar(Iend+1,j,kout)=ubar(Iend+1,j,know)+ &
1203!^ & dt2d*(bry_pgr+ &
1204!^ & bry_cor+ &
1205!^ & bry_str)
1206!^
1207 tl_ubar(iend+1,j,kout)=tl_ubar(iend+1,j,know)+ &
1208 & dt2d*(tl_bry_pgr+ &
1209 & tl_bry_cor+ &
1210 & tl_bry_str)
1211# ifdef MASKING
1212!^ ubar(Iend+1,j,kout)=ubar(Iend+1,j,kout)* &
1213!^ & GRID(ng)%umask(Iend+1,j)
1214!^
1215 tl_ubar(iend+1,j,kout)=tl_ubar(iend+1,j,kout)* &
1216 & grid(ng)%umask(iend+1,j)
1217# endif
1218 END IF
1219 END DO
1220!
1221! Eastern edge, closed boundary condition.
1222!
1223 ELSE IF (tl_lbc(ieast,isubar,ng)%closed) THEN
1224 DO j=jstr,jend
1225 IF (lbc_apply(ng)%east(j)) THEN
1226!^ ubar(Iend+1,j,kout)=0.0_r8
1227!^
1228 tl_ubar(iend+1,j,kout)=0.0_r8
1229 END IF
1230 END DO
1231 END IF
1232 END IF
1233!
1234!-----------------------------------------------------------------------
1235! Lateral boundary conditions at the southern edge.
1236!-----------------------------------------------------------------------
1237!
1238 IF (domain(ng)%Southern_Edge(tile)) THEN
1239!
1240! Southern edge, implicit upstream radiation condition.
1241!
1242 IF (tl_lbc(isouth,isubar,ng)%radiation) THEN
1243 IF (iic(ng).ne.0) THEN
1244 DO i=istru-1,iend
1245!^ grad(i,Jstr-1)=ubar(i+1,Jstr-1,know)- &
1246!^ & ubar(i ,Jstr-1,know)
1247!^
1248 tl_grad(i,jstr-1)=0.0_r8
1249 END DO
1250 DO i=istru,iend
1251 IF (lbc_apply(ng)%south(i)) THEN
1252# if defined CELERITY_READ && defined FORWARD_READ
1253 IF (tl_lbc(isouth,isubar,ng)%nudging) THEN
1254 IF (lnudgem2clm(ng)) THEN
1255 obc_out=0.5_r8* &
1256 & (clima(ng)%M2nudgcof(i-1,jstr-1)+ &
1257 & clima(ng)%M2nudgcof(i ,jstr-1))
1258 obc_in =obcfac(ng)*obc_out
1259 ELSE
1260 obc_out=m2obc_out(ng,isouth)
1261 obc_in =m2obc_in(ng,isouth)
1262 END IF
1263 IF (boundary(ng)%ubar_south_Ce(i).lt.0.0_r8) THEN
1264 tau=obc_in
1265 ELSE
1266 tau=obc_out
1267 END IF
1268 tau=tau*dt2d
1269 END IF
1270# ifdef RADIATION_2D
1271 cx=boundary(ng)%ubar_south_Cx(i)
1272# else
1273 cx=0.0_r8
1274# endif
1275 ce=boundary(ng)%ubar_south_Ce(i)
1276 cff=boundary(ng)%ubar_south_C2(i)
1277# endif
1278!^ ubar(i,Jstr-1,kout)=(cff*ubar(i,Jstr-1,know)+ &
1279!^ & Ce *ubar(i,Jstr ,kout)- &
1280!^ & MAX(Cx,0.0_r8)*grad(i-1,Jstr-1)- &
1281!^ & MIN(Cx,0.0_r8)*grad(i ,Jstr-1))/ &
1282!^ & (cff+Ce)
1283!^
1284 tl_ubar(i,jstr-1,kout)=(cff*tl_ubar(i,jstr-1,know)+ &
1285 & ce *tl_ubar(i,jstr ,kout)- &
1286 & max(cx,0.0_r8)* &
1287 & tl_grad(i-1,jstr-1)- &
1288 & min(cx,0.0_r8)* &
1289 & tl_grad(i ,jstr-1))/ &
1290 & (cff+ce)
1291
1292 IF (tl_lbc(isouth,isubar,ng)%nudging) THEN
1293!^ ubar(i,Jstr-1,kout)=ubar(i,Jstr-1,kout)+ &
1294!^ & tau*(BOUNDARY(ng)%ubar_south(i)- &
1295!^ & ubar(i,Jstr-1,know))
1296!^
1297 tl_ubar(i,jstr-1,kout)=tl_ubar(i,jstr-1,kout)- &
1298 & tau*tl_ubar(i,jstr-1,know)
1299 END IF
1300# ifdef MASKING
1301!^ ubar(i,Jstr-1,kout)=ubar(i,Jstr-1,kout)* &
1302!^ & GRID(ng)%umask(i,Jstr-1)
1303!^
1304 tl_ubar(i,jstr-1,kout)=tl_ubar(i,jstr-1,kout)* &
1305 & grid(ng)%umask(i,jstr-1)
1306# endif
1307 END IF
1308 END DO
1309 END IF
1310!
1311! Southern edge, Chapman boundary condition.
1312!
1313 ELSE IF (tl_lbc(isouth,isubar,ng)%Flather.or. &
1314 & tl_lbc(isouth,isubar,ng)%reduced.or. &
1315 & tl_lbc(isouth,isubar,ng)%Shchepetkin) THEN
1316 DO i=istru,iend
1317 IF (lbc_apply(ng)%south(i)) THEN
1318 cff=dt2d*0.5_r8*(grid(ng)%pn(i-1,jstr)+ &
1319 & grid(ng)%pn(i ,jstr))
1320 cff1=sqrt(g*0.5_r8*(grid(ng)%h(i-1,jstr)+ &
1321 & zeta(i-1,jstr,know)+ &
1322 & grid(ng)%h(i ,jstr)+ &
1323 & zeta(i ,jstr,know)))
1324 tl_cff1=0.25_r8*g*(grid(ng)%tl_h(i-1,jstr)+ &
1325 & tl_zeta(i-1,jstr,know)+ &
1326 & grid(ng)%tl_h(i ,jstr)+ &
1327 & tl_zeta(i ,jstr,know))/cff1
1328 ce=cff*cff1
1329 tl_ce=cff*tl_cff1
1330 cff2=1.0_r8/(1.0_r8+ce)
1331 tl_cff2=-cff2*cff2*tl_ce
1332!^ ubar(i,Jstr-1,kout)=cff2*(ubar(i,Jstr-1,know)+ &
1333!^ & Ce*ubar(i,Jstr,kout))
1334!^
1335 tl_ubar(i,jstr-1,kout)=tl_cff2*(ubar(i,jstr-1,know)+ &
1336 & ce*ubar(i,jstr,kout))+ &
1337 & cff2*(tl_ubar(i,jstr-1,know)+ &
1338 & tl_ce*ubar(i,jstr,kout)+ &
1339 & ce*tl_ubar(i,jstr,kout))
1340# ifdef MASKING
1341!^ ubar(i,Jstr-1,kout)=ubar(i,Jstr-1,kout)* &
1342!^ & GRID(ng)%umask(i,Jstr-1)
1343!^
1344 tl_ubar(i,jstr-1,kout)=tl_ubar(i,jstr-1,kout)* &
1345 & grid(ng)%umask(i,jstr-1)
1346# endif
1347 END IF
1348 END DO
1349!
1350! Southern edge, clamped boundary condition.
1351!
1352 ELSE IF (tl_lbc(isouth,isubar,ng)%clamped) THEN
1353 DO i=istru,iend
1354 IF (lbc_apply(ng)%south(i)) THEN
1355!^ ubar(i,Jstr-1,kout)=BOUNDARY(ng)%ubar_south(i)
1356!^
1357# ifdef ADJUST_BOUNDARY
1358 IF (lobc(isouth,isubar,ng)) THEN
1359 tl_ubar(i,jstr-1,kout)=boundary(ng)%tl_ubar_south(i)
1360 ELSE
1361 tl_ubar(i,jstr-1,kout)=0.0_r8
1362 END IF
1363# else
1364 tl_ubar(i,jstr-1,kout)=0.0_r8
1365# endif
1366# ifdef MASKING
1367!^ ubar(i,Jstr-1,kout)=ubar(i,Jstr-1,kout)* &
1368!^ & GRID(ng)%umask(i,Jstr-1)
1369!^
1370 tl_ubar(i,jstr-1,kout)=tl_ubar(i,jstr-1,kout)* &
1371 & grid(ng)%umask(i,jstr-1)
1372# endif
1373 END IF
1374 END DO
1375!
1376! Southern edge, gradient boundary condition.
1377!
1378 ELSE IF (tl_lbc(isouth,isubar,ng)%gradient) THEN
1379 DO i=istru,iend
1380 IF (lbc_apply(ng)%south(i)) THEN
1381!^ ubar(i,Jstr-1,kout)=ubar(i,Jstr,kout)
1382!^
1383 tl_ubar(i,jstr-1,kout)=tl_ubar(i,jstr,kout)
1384# ifdef MASKING
1385!^ ubar(i,Jstr-1,kout)=ubar(i,Jstr-1,kout)* &
1386!^ & GRID(ng)%umask(i,Jstr-1)
1387!^
1388 tl_ubar(i,jstr-1,kout)=tl_ubar(i,jstr-1,kout)* &
1389 & grid(ng)%umask(i,jstr-1)
1390# endif
1391 END IF
1392 END DO
1393!
1394! Southern edge, closed boundary condition: free slip (gamma2=1) or
1395! no slip (gamma2=-1).
1396!
1397 ELSE IF (tl_lbc(isouth,isubar,ng)%closed) THEN
1398 IF (ewperiodic(ng)) THEN
1399 imin=istru
1400 imax=iend
1401 ELSE
1402 imin=istr
1403 imax=iendr
1404 END IF
1405 DO i=imin,imax
1406 IF (lbc_apply(ng)%south(i)) THEN
1407!^ ubar(i,Jstr-1,kout)=gamma2(ng)*ubar(i,Jstr,kout)
1408!^
1409 tl_ubar(i,jstr-1,kout)=gamma2(ng)*tl_ubar(i,jstr,kout)
1410# ifdef MASKING
1411!^ ubar(i,Jstr-1,kout)=ubar(i,Jstr-1,kout)* &
1412!^ & GRID(ng)%umask(i,Jstr-1)
1413!^
1414 tl_ubar(i,jstr-1,kout)=tl_ubar(i,jstr-1,kout)* &
1415 & grid(ng)%umask(i,jstr-1)
1416# endif
1417 END IF
1418 END DO
1419 END IF
1420 END IF
1421!
1422!-----------------------------------------------------------------------
1423! Lateral boundary conditions at the northern edge.
1424!-----------------------------------------------------------------------
1425!
1426 IF (domain(ng)%Northern_Edge(tile)) THEN
1427!
1428! Northern edge, implicit upstream radiation condition.
1429!
1430 IF (tl_lbc(inorth,isubar,ng)%radiation) THEN
1431 IF (iic(ng).ne.0) THEN
1432 DO i=istru-1,iend
1433!^ grad(i,Jend+1)=ubar(i+1,Jend+1,know)- &
1434!^ & ubar(i ,Jend+1,know)
1435!^
1436 tl_grad(i,jend+1)=0.0_r8
1437 END DO
1438 DO i=istru,iend
1439 IF (lbc_apply(ng)%north(i)) THEN
1440# if defined CELERITY_READ && defined FORWARD_READ
1441 IF (tl_lbc(inorth,isubar,ng)%nudging) THEN
1442 IF (lnudgem2clm(ng)) THEN
1443 obc_out=0.5_r8* &
1444 & (clima(ng)%M2nudgcof(i-1,jend+1)+ &
1445 & clima(ng)%M2nudgcof(i ,jend+1))
1446 obc_in =obcfac(ng)*obc_out
1447 ELSE
1448 obc_out=m2obc_out(ng,inorth)
1449 obc_in =m2obc_in(ng,inorth)
1450 END IF
1451 IF (boundary(ng)%ubar_north_Ce(i).lt.0.0_r8) THEN
1452 tau=obc_in
1453 ELSE
1454 tau=obc_out
1455 END IF
1456 tau=tau*dt2d
1457 END IF
1458# ifdef RADIATION_2D
1459 cx=boundary(ng)%ubar_north_Cx(i)
1460# else
1461 cx=0.0_r8
1462# endif
1463 ce=boundary(ng)%ubar_north_Ce(i)
1464 cff=boundary(ng)%ubar_north_C2(i)
1465# endif
1466!^ ubar(i,Jend+1,kout)=(cff*ubar(i,Jend+1,know)+ &
1467!^ & Ce *ubar(i,Jend ,kout)- &
1468!^ & MAX(Cx,0.0_r8)*grad(i-1,Jend+1)- &
1469!^ & MIN(Cx,0.0_r8)*grad(i ,Jend+1))/ &
1470!^ & (cff+Ce)
1471!^
1472 tl_ubar(i,jend+1,kout)=(cff*tl_ubar(i,jend+1,know)+ &
1473 & ce *tl_ubar(i,jend ,kout)- &
1474 & max(cx,0.0_r8)* &
1475 & tl_grad(i-1,jend+1)- &
1476 & min(cx,0.0_r8)* &
1477 & tl_grad(i ,jend+1))/ &
1478 & (cff+ce)
1479
1480 IF (tl_lbc(inorth,isubar,ng)%nudging) THEN
1481!^ ubar(i,Jend+1,kout)=ubar(i,Jend+1,kout)+ &
1482!^ & tau*(BOUNDARY(ng)%ubar_north(i)- &
1483!^ & ubar(i,Jend+1,know))
1484!^
1485 tl_ubar(i,jend+1,kout)=tl_ubar(i,jend+1,kout)- &
1486 & tau*tl_ubar(i,jend+1,know)
1487 END IF
1488# ifdef MASKING
1489!^ ubar(i,Jend+1,kout)=ubar(i,Jend+1,kout)* &
1490!^ & GRID(ng)%umask(i,Jend+1)
1491!^
1492 tl_ubar(i,jend+1,kout)=tl_ubar(i,jend+1,kout)* &
1493 & grid(ng)%umask(i,jend+1)
1494# endif
1495 END IF
1496 END DO
1497 END IF
1498!
1499! Northern edge, Chapman boundary condition.
1500!
1501 ELSE IF (tl_lbc(inorth,isubar,ng)%Flather.or. &
1502 & tl_lbc(inorth,isubar,ng)%reduced.or. &
1503 & tl_lbc(inorth,isubar,ng)%Shchepetkin) THEN
1504 DO i=istru,iend
1505 IF (lbc_apply(ng)%north(i)) THEN
1506 cff=dt2d*0.5_r8*(grid(ng)%pn(i-1,jend)+ &
1507 & grid(ng)%pn(i ,jend))
1508 cff1=sqrt(g*0.5_r8*(grid(ng)%h(i-1,jend)+ &
1509 & zeta(i-1,jend,know)+ &
1510 & grid(ng)%h(i ,jend)+ &
1511 & zeta(i ,jend,know)))
1512 tl_cff1=0.25_r8*g*(grid(ng)%tl_h(i-1,jend)+ &
1513 & tl_zeta(i-1,jend,know)+ &
1514 & grid(ng)%tl_h(i ,jend)+ &
1515 & tl_zeta(i ,jend,know))/cff1
1516 ce=cff*cff1
1517 tl_ce=cff*tl_cff1
1518 cff2=1.0_r8/(1.0_r8+ce)
1519 tl_cff2=-cff2*cff2*tl_ce
1520!^ ubar(i,Jend+1,kout)=cff2*(ubar(i,Jend+1,know)+ &
1521!^ & Ce*ubar(i,Jend,kout))
1522!^
1523 tl_ubar(i,jend+1,kout)=tl_cff2*(ubar(i,jend+1,know)+ &
1524 & ce*ubar(i,jend,kout))+ &
1525 & cff2*(tl_ubar(i,jend+1,know)+ &
1526 & tl_ce*ubar(i,jend,kout)+ &
1527 & ce*tl_ubar(i,jend,kout))
1528# ifdef MASKING
1529!^ ubar(i,Jend+1,kout)=ubar(i,Jend+1,kout)* &
1530!^ & GRID(ng)%umask(i,Jend+1)
1531!^
1532 tl_ubar(i,jend+1,kout)=tl_ubar(i,jend+1,kout)* &
1533 & grid(ng)%umask(i,jend+1)
1534# endif
1535 END IF
1536 END DO
1537!
1538! Northern edge, clamped boundary condition.
1539!
1540 ELSE IF (tl_lbc(inorth,isubar,ng)%clamped) THEN
1541 DO i=istru,iend
1542 IF (lbc_apply(ng)%north(i)) THEN
1543!^ ubar(i,Jend+1,kout)=BOUNDARY(ng)%ubar_north(i)
1544!^
1545# ifdef ADJUST_BOUNDARY
1546 IF (lobc(inorth,isubar,ng)) THEN
1547 tl_ubar(i,jend+1,kout)=boundary(ng)%tl_ubar_north(i)
1548 ELSE
1549 tl_ubar(i,jend+1,kout)=0.0_r8
1550 END IF
1551# else
1552 tl_ubar(i,jend+1,kout)=0.0_r8
1553# endif
1554# ifdef MASKING
1555!^ ubar(i,Jend+1,kout)=ubar(i,Jend+1,kout)* &
1556!^ & GRID(ng)%umask(i,Jend+1)
1557!^
1558 tl_ubar(i,jend+1,kout)=tl_ubar(i,jend+1,kout)* &
1559 & grid(ng)%umask(i,jend+1)
1560# endif
1561 END IF
1562 END DO
1563!
1564! Northern edge, gradient boundary condition.
1565!
1566 ELSE IF (tl_lbc(inorth,isubar,ng)%gradient) THEN
1567 DO i=istru,iend
1568 IF (lbc_apply(ng)%north(i)) THEN
1569!^ ubar(i,Jend+1,kout)=ubar(i,Jend,kout)
1570!^
1571 tl_ubar(i,jend+1,kout)=tl_ubar(i,jend,kout)
1572# ifdef MASKING
1573!^ ubar(i,Jend+1,kout)=ubar(i,Jend+1,kout)* &
1574!^ & GRID(ng)%umask(i,Jend+1)
1575!^
1576 tl_ubar(i,jend+1,kout)=tl_ubar(i,jend+1,kout)* &
1577 & grid(ng)%umask(i,jend+1)
1578# endif
1579 END IF
1580 END DO
1581!
1582! Northern edge, closed boundary condition: free slip (gamma2=1) or
1583! no slip (gamma2=-1).
1584!
1585 ELSE IF (tl_lbc(inorth,isubar,ng)%closed) THEN
1586 IF (ewperiodic(ng)) THEN
1587 imin=istru
1588 imax=iend
1589 ELSE
1590 imin=istr
1591 imax=iendr
1592 END IF
1593 DO i=imin,imax
1594 IF (lbc_apply(ng)%north(i)) THEN
1595!^ ubar(i,Jend+1,kout)=gamma2(ng)*ubar(i,Jend,kout)
1596!^
1597 tl_ubar(i,jend+1,kout)=gamma2(ng)*tl_ubar(i,jend,kout)
1598
1599# ifdef MASKING
1600!^ ubar(i,Jend+1,kout)=ubar(i,Jend+1,kout)* &
1601!^ & GRID(ng)%GRID(ng)%umask(i,Jend+1)
1602!^
1603 tl_ubar(i,jend+1,kout)=tl_ubar(i,jend+1,kout)* &
1604 & grid(ng)%umask(i,jend+1)
1605# endif
1606 END IF
1607 END DO
1608 END IF
1609 END IF
1610!
1611!-----------------------------------------------------------------------
1612! Boundary corners.
1613!-----------------------------------------------------------------------
1614!
1615 IF (.not.(ewperiodic(ng).or.nsperiodic(ng))) THEN
1616 IF (domain(ng)%SouthWest_Corner(tile)) THEN
1617 IF (lbc_apply(ng)%south(istr ).and. &
1618 & lbc_apply(ng)%west (jstr-1)) THEN
1619!^ ubar(Istr,Jstr-1,kout)=0.5_r8*(ubar(Istr+1,Jstr-1,kout)+ &
1620!^ & ubar(Istr ,Jstr ,kout))
1621!^
1622 tl_ubar(istr,jstr-1,kout)=0.5_r8* &
1623 & (tl_ubar(istr+1,jstr-1,kout)+ &
1624 & tl_ubar(istr ,jstr ,kout))
1625 END IF
1626 END IF
1627 IF (domain(ng)%SouthEast_Corner(tile)) THEN
1628 IF (lbc_apply(ng)%south(iend+1).and. &
1629 & lbc_apply(ng)%east (jstr-1)) THEN
1630!^ ubar(Iend+1,Jstr-1,kout)=0.5_r8*(ubar(Iend ,Jstr-1,kout)+ &
1631!^ & ubar(Iend+1,Jstr ,kout))
1632!^
1633 tl_ubar(iend+1,jstr-1,kout)=0.5_r8* &
1634 & (tl_ubar(iend ,jstr-1,kout)+ &
1635 & tl_ubar(iend+1,jstr ,kout))
1636 END IF
1637 END IF
1638 IF (domain(ng)%NorthWest_Corner(tile)) THEN
1639 IF (lbc_apply(ng)%north(istr ).and. &
1640 & lbc_apply(ng)%west (jend+1)) THEN
1641!^ ubar(Istr,Jend+1,kout)=0.5_r8*(ubar(Istr ,Jend ,kout)+ &
1642!^ & ubar(Istr+1,Jend+1,kout))
1643!^
1644 tl_ubar(istr,jend+1,kout)=0.5_r8* &
1645 & (tl_ubar(istr ,jend ,kout)+ &
1646 & tl_ubar(istr+1,jend+1,kout))
1647 END IF
1648 END IF
1649 IF (domain(ng)%NorthEast_Corner(tile)) THEN
1650 IF (lbc_apply(ng)%north(iend+1).and. &
1651 & lbc_apply(ng)%east (jend+1)) THEN
1652!^ ubar(Iend+1,Jend+1,kout)=0.5_r8*(ubar(Iend+1,Jend ,kout)+ &
1653!^ & ubar(Iend ,Jend+1,kout))
1654!^
1655 tl_ubar(iend+1,jend+1,kout)=0.5_r8* &
1656 & (tl_ubar(iend+1,jend ,kout)+ &
1657 & tl_ubar(iend ,jend+1,kout))
1658 END IF
1659 END IF
1660 END IF
1661
1662# if defined WET_DRY_NOT_YET
1663!
1664!-----------------------------------------------------------------------
1665! Impose wetting and drying conditions.
1666!
1667! HGA: need TLM code here.
1668!-----------------------------------------------------------------------
1669!
1670 IF (.not.ewperiodic(ng)) THEN
1671 IF (domain(ng)%Western_Edge(tile)) THEN
1672 DO j=jstr,jend
1673 IF (lbc_apply(ng)%west(j).or. &
1674 & lbc(iwest,isubar,ng)%nested) THEN
1675!^ cff1=ABS(ABS(GRID(ng)%umask_wet(Istr,j))-1.0_r8)
1676!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(Istr,j,kout))* &
1677!^ & GRID(ng)%umask_wet(Istr,j)
1678!^ cff=0.5_r8*GRID(ng)%umask_wet(Istr,j)*cff1+ &
1679!^ & cff2*(1.0_r8-cff1)
1680!^ ubar(Istr,j,kout)=ubar(Istr,j,kout)*cff
1681 END IF
1682 END DO
1683 END IF
1684 IF (domain(ng)%Eastern_Edge(tile)) THEN
1685 DO j=jstr,jend
1686 IF (lbc_apply(ng)%east(j).or. &
1687 & lbc(ieast,isubar,ng)%nested) THEN
1688!^ cff1=ABS(ABS(GRID(ng)%umask_wet(Iend+1,j))-1.0_r8)
1689!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(Iend+1,j,kout))* &
1690!^ & GRID(ng)%umask_wet(Iend+1,j)
1691!^ cff=0.5_r8*GRID(ng)%umask_wet(Iend+1,j)*cff1+ &
1692!^ & cff2*(1.0_r8-cff1)
1693!^ ubar(Iend+1,j,kout)=ubar(Iend+1,j,kout)*cff
1694 END IF
1695 END DO
1696 END IF
1697 END IF
1698!
1699 IF (.not.nsperiodic(ng)) THEN
1700 IF (domain(ng)%Southern_Edge(tile)) THEN
1701 DO i=istru,iend
1702 IF (lbc_apply(ng)%south(i).or. &
1703 & lbc(isouth,isubar,ng)%nested) THEN
1704!^ cff1=ABS(ABS(GRID(ng)%umask_wet(i,Jstr-1))-1.0_r8)
1705!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(i,Jstr-1,kout))* &
1706!^ & GRID(ng)%umask_wet(i,Jstr-1)
1707!^ cff=0.5_r8*GRID(ng)%umask_wet(i,Jstr-1)*cff1+ &
1708!^ & cff2*(1.0_r8-cff1)
1709!^ ubar(i,Jstr-1,kout)=ubar(i,Jstr-1,kout)*cff
1710 END IF
1711 END DO
1712 END IF
1713 IF (domain(ng)%Northern_Edge(tile)) THEN
1714 DO i=istr,iend
1715 IF (lbc_apply(ng)%north(i).or. &
1716 & lbc(inorth,isubar,ng)%nested) THEN
1717!^ cff1=ABS(ABS(GRID(ng)%umask_wet(i,Jend+1))-1.0_r8)
1718!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(i,Jend+1,kout))* &
1719!^ & GRID(ng)%umask_wet(i,Jend+1)
1720!^ cff=0.5_r8*GRID(ng)%umask_wet(i,Jend+1)*cff1+ &
1721!^ & cff2*(1.0_r8-cff1)
1722!^ ubar(i,Jend+1,kout)=ubar(i,Jend+1,kout)*cff
1723 END IF
1724 END DO
1725 END IF
1726 END IF
1727!
1728 IF (.not.(ewperiodic(ng).or.nsperiodic(ng))) THEN
1729 IF (domain(ng)%SouthWest_Corner(tile)) THEN
1730 IF ((lbc_apply(ng)%south(istr ).and. &
1731 & lbc_apply(ng)%west (jstr-1)).or. &
1732 & (lbc(iwest,isubar,ng)%nested.and. &
1733 & lbc(isouth,isubar,ng)%nested)) THEN
1734!^ cff1=ABS(ABS(GRID(ng)%umask_wet(Istr,Jstr-1))-1.0_r8)
1735!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(Istr,Jstr-1,kout))* &
1736!^ & GRID(ng)%umask_wet(Istr,Jstr-1)
1737!^ cff=0.5_r8*GRID(ng)%umask_wet(Istr,Jstr-1)*cff1+ &
1738!^ & cff2*(1.0_r8-cff1)
1739!^ ubar(Istr,Jstr-1,kout)=ubar(Istr,Jstr-1,kout)*cff
1740 END IF
1741 END IF
1742 IF (domain(ng)%SouthEast_Corner(tile)) THEN
1743 IF ((lbc_apply(ng)%south(iend+1).and. &
1744 & lbc_apply(ng)%east (jstr-1)).or. &
1745 & (lbc(ieast,isubar,ng)%nested.and. &
1746 & lbc(isouth,isubar,ng)%nested)) THEN
1747!^ cff1=ABS(ABS(GRID(ng)%umask_wet(Iend+1,Jstr-1))-1.0_r8)
1748!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(Iend+1,Jstr-1,kout))* &
1749!^ & GRID(ng)%umask_wet(Iend+1,Jstr-1)
1750!^ cff=0.5_r8*GRID(ng)%umask_wet(Iend+1,Jstr-1)*cff1+ &
1751!^ & cff2*(1.0_r8-cff1)
1752!^ ubar(Iend+1,Jstr-1,kout)=ubar(Iend+1,Jstr-1,kout)*cff
1753 END IF
1754 END IF
1755 IF (domain(ng)%NorthWest_Corner(tile)) THEN
1756 IF ((lbc_apply(ng)%north(istr ).and. &
1757 & lbc_apply(ng)%west (jend+1)).or. &
1758 & (lbc(iwest,isubar,ng)%nested.and. &
1759 & lbc(inorth,isubar,ng)%nested)) THEN
1760!^ cff1=ABS(ABS(GRID(ng)%umask_wet(Istr,Jend+1))-1.0_r8)
1761!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(Istr,Jend+1,kout))* &
1762!^ & GRID(ng)%umask_wet(Istr,Jend+1)
1763!^ cff=0.5_r8*GRID(ng)%umask_wet(Istr,Jend+1)*cff1+ &
1764!^ & cff2*(1.0_r8-cff1)
1765!^ ubar(Istr,Jend+1,kout)=ubar(Istr,Jend+1,kout)*cff
1766 END IF
1767 END IF
1768 IF (domain(ng)%NorthEast_Corner(tile)) THEN
1769 IF ((lbc_apply(ng)%north(iend+1).and. &
1770 & lbc_apply(ng)%east (jend+1)).or. &
1771 & (lbc(ieast,isubar,ng)%nested.and. &
1772 & lbc(inorth,isubar,ng)%nested)) THEN
1773!^ cff1=ABS(ABS(GRID(ng)%umask_wet(Iend+1,Jend+1))-1.0_r8)
1774!^ cff2=0.5_r8+DSIGN(0.5_r8,ubar(Iend+1,Jend+1,kout))* &
1775!^ & GRID(ng)%umask_wet(Iend+1,Jend+1)
1776!^ cff=0.5_r8*GRID(ng)%umask_wet(Iend+1,Jend+1)*cff1+ &
1777!^ & cff2*(1.0_r8-cff1)
1778!^ ubar(Iend+1,Jend+1,kout)=ubar(Iend+1,Jend+1,kout)*cff
1779 END IF
1780 END IF
1781 END IF
1782# endif
1783!
1784 RETURN
1785 END SUBROUTINE tl_u2dbc_tile
1786#endif
1787 END MODULE tl_u2dbc_mod
type(t_boundary), dimension(:), allocatable boundary
type(t_apply), dimension(:), allocatable lbc_apply
type(t_clima), dimension(:), allocatable clima
Definition mod_clima.F:153
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
integer isfsur
integer isubar
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
type(t_lbc), dimension(:,:,:), allocatable tl_lbc
Definition mod_param.F:379
type(t_lbc), dimension(:,:,:), allocatable lbc
Definition mod_param.F:375
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
real(r8) co
logical, dimension(:), allocatable lnudgem2clm
integer, dimension(:), allocatable iic
logical, dimension(:,:,:), allocatable lobc
logical, dimension(:), allocatable ewperiodic
integer, parameter iwest
logical, dimension(:), allocatable nsperiodic
logical, dimension(:), allocatable predictor_2d_step
real(dp), dimension(:), allocatable obcfac
real(r8), dimension(:), allocatable gamma2
integer, parameter isouth
real(dp), dimension(:), allocatable dtfast
real(dp), dimension(:,:), allocatable m2obc_out
integer, parameter ieast
real(dp) g
real(dp) rho0
integer, parameter inorth
real(dp), dimension(:,:), allocatable m2obc_in
integer, dimension(:), allocatable kstp
integer, dimension(:), allocatable krhs
subroutine, public tl_u2dbc(ng, tile, kout)
Definition tl_u2dbc_im.F:29
subroutine, public tl_u2dbc_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, krhs, kstp, kout, ubar, vbar, zeta, tl_ubar, tl_vbar, tl_zeta)
Definition tl_u2dbc_im.F:64