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

Functions/Subroutines

subroutine, public tl_v2dbc (ng, tile, kout)
 
subroutine, public tl_v2dbc_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, krhs, kstp, kout, ubar, vbar, zeta, tl_ubar, tl_vbar, tl_zeta)
 

Function/Subroutine Documentation

◆ tl_v2dbc()

subroutine, public tl_v2dbc_mod::tl_v2dbc ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) kout )

Definition at line 28 of file tl_v2dbc_im.F.

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_v2dbc_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
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, dimension(:), allocatable kstp
integer, dimension(:), allocatable krhs

References mod_stepping::krhs, mod_stepping::kstp, mod_ocean::ocean, and tl_v2dbc_tile().

Here is the call graph for this function:

◆ tl_v2dbc_tile()

subroutine, public tl_v2dbc_mod::tl_v2dbc_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) krhs,
integer, intent(in) kstp,
integer, intent(in) kout,
real(r8), dimension(lbi:,lbj:,:), intent(in) ubar,
real(r8), dimension(lbi:,lbj:,:), intent(in) vbar,
real(r8), dimension(lbi:,lbj:,:), intent(in) zeta,
real(r8), dimension(lbi:,lbj:,:), intent(in) tl_ubar,
real(r8), dimension(lbi:,lbj:,:), intent(inout) tl_vbar,
real(r8), dimension(lbi:,lbj:,:), intent(in) tl_zeta )

Definition at line 58 of file tl_v2dbc_im.F.

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_ubar(LBi:,LBj:,:)
86 real(r8), intent(in) :: tl_zeta(LBi:,LBj:,:)
87
88 real(r8), intent(inout) :: tl_vbar(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_ubar(LBi:UBi,LBj:UBj,:)
94 real(r8), intent(in) :: tl_zeta(LBi:UBi,LBj:UBj,:)
95
96 real(r8), intent(inout) :: tl_vbar(LBi:UBi,LBj:UBj,:)
97# endif
98!
99! Local variable declarations.
100!
101 integer :: Jmin, Jmax
102 integer :: i, j, know
103
104 real(r8) :: Ce, Cx, Ze
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_Ze
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 southern edge.
141!-----------------------------------------------------------------------
142!
143 IF (domain(ng)%Southern_Edge(tile)) THEN
144!
145! Southern edge, implicit upstream radiation condition.
146!
147 IF (tl_lbc(isouth,isvbar,ng)%radiation) THEN
148 IF (iic(ng).ne.0) THEN
149 DO i=istr,iend+1
150!^ grad(i,Jstr)=vbar(i ,Jstr,know)- &
151!^ & vbar(i-1,Jstr,know)
152!^
153 tl_grad(i,jstr)=0.0_r8
154 END DO
155 DO i=istr,iend
156 IF (lbc_apply(ng)%south(i)) THEN
157# if defined CELERITY_READ && defined FORWARD_READ
158 IF (tl_lbc(isouth,isvbar,ng)%nudging) THEN
159 IF (lnudgem2clm(ng)) THEN
160 obc_out=0.5_r8* &
161 & (clima(ng)%M2nudgcof(i,jstr-1)+ &
162 & clima(ng)%M2nudgcof(i,jstr ))
163 obc_in =obcfac(ng)*obc_out
164 ELSE
165 obc_out=m2obc_out(ng,isouth)
166 obc_in =m2obc_in(ng,isouth)
167 END IF
168 IF (boundary(ng)%vbar_south_Ce(i).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# ifdef RADIATION_2D
176 cx=boundary(ng)%vbar_south_Cx(i)
177# else
178 cx=0.0_r8
179# endif
180 ce=boundary(ng)%vbar_south_Ce(i)
181 cff=boundary(ng)%vbar_south_C2(i)
182# endif
183!^ vbar(i,Jstr,kout)=(cff*vbar(i,Jstr ,know)+ &
184!^ & Ce *vbar(i,Jstr+1,kout)- &
185!^ & MAX(Cx,0.0_r8)*grad(i ,Jstr)- &
186!^ & MIN(Cx,0.0_r8)*grad(i+1,Jstr))/ &
187!^ & (cff+Ce)
188!^
189 tl_vbar(i,jstr,kout)=(cff*tl_vbar(i,jstr ,know)+ &
190 & ce *tl_vbar(i,jstr+1,kout)- &
191 & max(cx,0.0_r8)* &
192 & tl_grad(i ,jstr)- &
193 & min(cx,0.0_r8)* &
194 & tl_grad(i+1,jstr))/ &
195 & (cff+ce)
196
197 IF (tl_lbc(isouth,isvbar,ng)%nudging) THEN
198!^ vbar(i,Jstr,kout)=vbar(i,Jstr,kout)+ &
199!^ & tau*(BOUNDARY(ng)%vbar_south(i)- &
200!^ & vbar(i,Jstr,know))
201!^
202 tl_vbar(i,jstr,kout)=tl_vbar(i,jstr,kout)- &
203 & tau*tl_vbar(i,jstr,know)
204 END IF
205# ifdef MASKING
206!^ vbar(i,Jstr,kout)=vbar(i,Jstr,kout)* &
207!^ & GRID(ng)%vmask(i,Jstr)
208!^
209 tl_vbar(i,jstr,kout)=tl_vbar(i,jstr,kout)* &
210 & grid(ng)%vmask(i,jstr)
211# endif
212 END IF
213 END DO
214 END IF
215!
216! Southern edge, Flather boundary condition.
217!
218 ELSE IF (tl_lbc(isouth,isvbar,ng)%Flather) THEN
219 DO i=istr,iend
220 IF (lbc_apply(ng)%south(i)) THEN
221# if defined SSH_TIDES_NOT_YET && !defined UV_TIDES_NOT_YET
222 IF (tl_lbc(isouth,isfsur,ng)%acquire) THEN
223 bry_pgr=-g*(zeta(i,jstr,know)- &
224 & boundary(ng)%zeta_south(i))* &
225 & 0.5_r8*grid(ng)%pn(i,jstr)
226 tl_bry_pgr=-g*tl_zeta(i,jstr,know)* &
227 & 0.5_r8*grid(ng)%pn(i,jstr)
228# ifdef ADJUST_BOUNDARY
229 IF (lobc(isouth,isvbar,ng)) THEN
230 tl_bry_pgr=tl_bry_pgr+ &
231 & g*boundary(ng)%tl_zeta_south(i)* &
232 & 0.5_r8*grid(ng)%pn(i,jstr)
233 END IF
234# endif
235 ELSE
236 bry_pgr=-g*(zeta(i,jstr ,know)- &
237 & zeta(i,jstr-1,know))* &
238 & 0.5_r8*(grid(ng)%pn(i,jstr-1)+ &
239 & grid(ng)%pn(i,jstr ))
240 tl_bry_pgr=-g*(tl_zeta(i,jstr ,know)- &
241 & tl_zeta(i,jstr-1,know))* &
242 & 0.5_r8*(grid(ng)%pn(i,jstr-1)+ &
243 & grid(ng)%pn(i,jstr ))
244 END IF
245# ifdef UV_COR
246 bry_cor=-0.125_r8*(ubar(i ,jstr-1,know)+ &
247 & ubar(i+1,jstr-1,know)+ &
248 & ubar(i ,jstr ,know)+ &
249 & ubar(i+1,jstr ,know))* &
250 & (grid(ng)%f(i,jstr-1)+ &
251 & grid(ng)%f(i,jstr ))
252 tl_bry_cor=-0.125_r8*(tl_ubar(i ,jstr-1,know)+ &
253 & tl_ubar(i+1,jstr-1,know)+ &
254 & tl_ubar(i ,jstr ,know)+ &
255 & tl_ubar(i+1,jstr ,know))* &
256 & (grid(ng)%f(i,jstr-1)+ &
257 & grid(ng)%f(i,jstr ))
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(i,jstr-1)+ &
263 & zeta(i,jstr-1,know)+ &
264 & grid(ng)%h(i,jstr )+ &
265 & zeta(i,jstr ,know)))
266 tl_cff1=-cff1*cff1*(0.5_r8*(grid(ng)%tl_h(i,jstr-1)+ &
267 & tl_zeta(i,jstr-1,know)+ &
268 & grid(ng)%tl_h(i,jstr )+ &
269 & tl_zeta(i,jstr ,know)))
270 bry_str=cff1*(forces(ng)%svstr(i,jstr)- &
271 & forces(ng)%bvstr(i,jstr))
272 tl_bry_str=tl_cff1*(forces(ng)%svstr(i,jstr)- &
273 & forces(ng)%bvstr(i,jstr))+ &
274 & cff1*(forces(ng)%tl_svstr(i,jstr)- &
275 & forces(ng)%tl_bvstr(i,jstr))
276 ce=1.0_r8/sqrt(g*0.5_r8*(grid(ng)%h(i,jstr-1)+ &
277 & zeta(i,jstr-1,know)+ &
278 & grid(ng)%h(i,jstr )+ &
279 & zeta(i,jstr ,know)))
280 tl_ce=-ce*ce*ce*0.25_r8*g*(grid(ng)%tl_h(i,jstr-1)+ &
281 & tl_zeta(i,jstr-1,know)+ &
282 & grid(ng)%tl_h(i,jstr )+ &
283 & tl_zeta(i,jstr ,know))
284 cff2=grid(ng)%on_v(i,jstr)*ce
285 tl_cff2=grid(ng)%on_v(i,jstr)*tl_ce
286!^ bry_val=vbar(i,Jstr+1,know)+ &
287!^ & cff2*(bry_pgr+ &
288!^ & bry_cor+ &
289!^ & bry_str)
290!^
291 tl_bry_val=tl_vbar(i,jstr+1,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)%vbar_south(i)
300!^
301# ifdef ADJUST_BOUNDARY
302 IF (lobc(isouth,isvbar,ng)) THEN
303 tl_bry_val=boundary(ng)%tl_vbar_south(i)
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(i,jstr-1)+ &
312 & zeta(i,jstr-1,know)+ &
313 & grid(ng)%h(i,jstr )+ &
314 & zeta(i,jstr ,know)))
315 tl_cff=-cff*cff*(0.5_r8*(grid(ng)%tl_h(i,jstr-1)+ &
316 & tl_zeta(i,jstr-1,know)+ &
317 & grid(ng)%tl_h(i,jstr )+ &
318 & tl_zeta(i,jstr ,know)))
319 ce=sqrt(g*cff)
320 tl_ce=0.5_r8*g*tl_cff/ce
321# if defined ATM_PRESS && defined PRESS_COMPENSATE
322!^ vbar(i,Jstr,kout)=bry_val- &
323!^ & Ce*(0.5_r8* &
324!^ & (zeta(i,Jstr-1,know)+ &
325!^ & zeta(i,Jstr ,know)+ &
326!^ & fac*(FORCES(ng)%Pair(i,Jstr-1)+ &
327!^ & FORCES(ng)%Pair(i,Jstr )- &
328!^ & 2.0_r8*OneAtm))- &
329!^ & BOUNDARY(ng)%zeta_south(i))
330!^
331 tl_vbar(i,jstr,kout)=tl_bry_val- &
332 & tl_ce* &
333 & (0.5_r8* &
334 & (zeta(i,jstr-1,know)+ &
335 & zeta(i,jstr ,know)+ &
336 & fac*(forces(ng)%Pair(i,jstr-1)+ &
337 & forces(ng)%Pair(i,jstr )- &
338 & 2.0_r8*oneatm))- &
339 & boundary(ng)%zeta_south(i))- &
340 & ce* &
341 & (0.5_r8*(tl_zeta(i,jstr-1,know)+ &
342 & tl_zeta(i,jstr ,know)))
343# else
344!^ vbar(i,Jstr,kout)=bry_val- &
345!^ & Ce*(0.5_r8*(zeta(i,Jstr-1,know)+ &
346!^ & zeta(i,Jstr ,know))- &
347!^ & BOUNDARY(ng)%zeta_south(i))
348!^
349 tl_vbar(i,jstr,kout)=tl_bry_val- &
350 & tl_ce* &
351 & (0.5_r8*(zeta(i,jstr-1,know)+ &
352 & zeta(i,jstr ,know))- &
353 & boundary(ng)%zeta_south(i))- &
354 & ce* &
355 & (0.5_r8*(tl_zeta(i,jstr-1,know)+ &
356 & tl_zeta(i,jstr ,know)))
357# endif
358# ifdef ADJUST_BOUNDARY
359 IF (lobc(isouth,isvbar,ng)) THEN
360 tl_vbar(i,jstr,kout)=tl_vbar(i,jstr,kout)+ &
361 & ce*boundary(ng)%tl_zeta_south(i)
362 END IF
363# endif
364# ifdef MASKING
365!^ vbar(i,Jstr,kout)=vbar(i,Jstr,kout)* &
366!^ & GRID(ng)%vmask(i,Jstr)
367!^
368 tl_vbar(i,jstr,kout)=tl_vbar(i,jstr,kout)* &
369 & grid(ng)%vmask(i,jstr)
370# endif
371 END IF
372 END DO
373!
374! Southern edge, Shchepetkin boundary condition (Maison et al., 2010).
375!
376 ELSE IF (tl_lbc(isouth,isvbar,ng)%Shchepetkin) THEN
377 DO i=istr,iend
378 IF (lbc_apply(ng)%south(i)) THEN
379# if defined SSH_TIDES_NOT_YET && !defined UV_TIDES_NOT_YET
380 IF (tl_lbc(isouth,isfsur,ng)%acquire) THEN
381 bry_pgr=-g*(zeta(i,jstr,know)- &
382 & boundary(ng)%zeta_south(i))* &
383 & 0.5_r8*grid(ng)%pn(i,jstr)
384 tl_bry_pgr=-g*tl_zeta(i,jstr,know)* &
385 & 0.5_r8*grid(ng)%pn(i,jstr)
386# ifdef ADJUST_BOUNDARY
387 IF (lobc(isouth,isvbar,ng)) THEN
388 tl_bry_pgr=tl_bry_pgr+ &
389 & g*boundary(ng)%tl_zeta_south(i)* &
390 & 0.5_r8*grid(ng)%pn(i,jstr)
391 END IF
392# endif
393 ELSE
394 bry_pgr=-g*(zeta(i,jstr ,know)- &
395 & zeta(i,jstr-1,know))* &
396 & 0.5_r8*(grid(ng)%pn(i,jstr-1)+ &
397 & grid(ng)%pn(i,jstr ))
398 tl_bry_pgr=-g*(tl_zeta(i,jstr ,know)- &
399 & tl_zeta(i,jstr-1,know))* &
400 & 0.5_r8*(grid(ng)%pn(i,jstr-1)+ &
401 & grid(ng)%pn(i,jstr ))
402 END IF
403# ifdef UV_COR
404 bry_cor=-0.125_r8*(ubar(i ,jstr-1,know)+ &
405 & ubar(i+1,jstr-1,know)+ &
406 & ubar(i ,jstr ,know)+ &
407 & ubar(i+1,jstr ,know))* &
408 & (grid(ng)%f(i,jstr-1)+ &
409 & grid(ng)%f(i,jstr ))
410 tl_bry_cor=-0.125_r8*(tl_ubar(i ,jstr-1,know)+ &
411 & tl_ubar(i+1,jstr-1,know)+ &
412 & tl_ubar(i ,jstr ,know)+ &
413 & tl_ubar(i+1,jstr ,know))* &
414 & (grid(ng)%f(i,jstr-1)+ &
415 & grid(ng)%f(i,jstr ))
416# else
417 bry_cor=0.0_r8
418 tl_bry_cor=0.0_r8
419# endif
420 cff1=1.0_r8/(0.5_r8*(grid(ng)%h(i,jstr-1)+ &
421 & zeta(i,jstr-1,know)+ &
422 & grid(ng)%h(i,jstr )+ &
423 & zeta(i,jstr ,know)))
424 tl_cff1=-cff1*cff1*(0.5_r8*(grid(ng)%tl_h(i,jstr-1)+ &
425 & tl_zeta(i,jstr-1,know)+ &
426 & grid(ng)%tl_h(i,jstr )+ &
427 & tl_zeta(i,jstr ,know)))
428 bry_str=cff1*(forces(ng)%svstr(i,jstr)- &
429 & forces(ng)%bvstr(i,jstr))
430 tl_bry_str=tl_cff1*(forces(ng)%svstr(i,jstr)- &
431 & forces(ng)%bvstr(i,jstr))+ &
432 & cff1*(forces(ng)%tl_svstr(i,jstr)- &
433 & forces(ng)%tl_bvstr(i,jstr))
434 ce=1.0_r8/sqrt(g*0.5_r8*(grid(ng)%h(i,jstr-1)+ &
435 & zeta(i,jstr-1,know)+ &
436 & grid(ng)%h(i,jstr )+ &
437 & zeta(i,jstr ,know)))
438 tl_ce=-ce*ce*ce*0.25_r8*g*(grid(ng)%tl_h(i,jstr-1)+ &
439 & tl_zeta(i,jstr-1,know)+ &
440 & grid(ng)%tl_h(i,jstr )+ &
441 & tl_zeta(i,jstr ,know))
442 cff2=grid(ng)%on_v(i,jstr)*ce
443 tl_cff2=grid(ng)%on_v(i,jstr)*tl_ce
444!^ bry_val=vbar(i,Jstr+1,know)+ &
445!^ & cff2*(bry_pgr+ &
446!^ & bry_cor+ &
447!^ & bry_str)
448!^
449 tl_bry_val=tl_vbar(i,jstr+1,know)+ &
450 & tl_cff2*(bry_pgr+ &
451 & bry_cor+ &
452 & bry_str)+ &
453 & cff2*(tl_bry_pgr+ &
454 & tl_bry_cor+ &
455 & tl_bry_str)
456# else
457!^ bry_val=BOUNDARY(ng)%vbar_south(i)
458!^
459# ifdef ADJUST_BOUNDARY
460 IF (lobc(isouth,isvbar,ng)) THEN
461 tl_bry_val=boundary(ng)%tl_vbar_south(i)
462 ELSE
463 tl_bry_val=0.0_r8
464 END IF
465# else
466 tl_bry_val=0.0_r8
467# endif
468# endif
469# ifdef WET_DRY_NOT_YET
470 cff=0.5_r8*(grid(ng)%h(i,jstr-1)+ &
471 & zeta(i,jstr-1,know)+ &
472 & grid(ng)%h(i,jstr )+ &
473 & zeta(i,jstr ,know))
474 tl_cff=0.5_r8*(grid(ng)%tl_h(i,jstr-1)+ &
475 & tl_zeta(i,jstr-1,know)+ &
476 & grid(ng)%tl_h(i,jstr )+ &
477 & tl_zeta(i,jstr ,know))
478# else
479 cff=0.5_r8*(grid(ng)%h(i,jstr-1)+ &
480 & grid(ng)%h(i,jstr ))
481 tl_cff=0.5_r8*(grid(ng)%tl_h(i,jstr-1)+ &
482 & grid(ng)%tl_h(i,jstr ))
483# endif
484 cff1=sqrt(g/cff)
485 tl_cff1=-0.5_r8*cff1*tl_cff/cff
486 ce=dt2d*cff1*cff*0.5_r8*(grid(ng)%pn(i,jstr-1)+ &
487 & grid(ng)%pn(i,jstr ))
488 tl_ce=dt2d*0.5_r8*(grid(ng)%pn(i,jstr-1)+ &
489 & grid(ng)%pn(i,jstr ))* &
490 & (cff1*tl_cff+ &
491 & tl_cff1*cff)
492 ze=(0.5_r8+ce)*zeta(i,jstr ,know)+ &
493 & (0.5_r8-ce)*zeta(i,jstr-1,know)
494 tl_ze=(0.5_r8+ce)*tl_zeta(i,jstr ,know)+ &
495 & (0.5_r8-ce)*tl_zeta(i,jstr-1,know)+ &
496 & tl_ce*(zeta(i,jstr ,know)- &
497 & zeta(i,jstr-1,know))
498 IF (ce.gt.co) THEN
499 cff2=(1.0_r8-co/ce)**2
500 tl_cff2=2.0_r8*cff2*co*tl_ce/(ce*ce)
501 cff3=zeta(i,jstr,kout)+ &
502 & ce*zeta(i,jstr-1,know)- &
503 & (1.0_r8+ce)*zeta(i,jstr,know)
504 tl_cff3=tl_zeta(i,jstr,kout)+ &
505 & ce*tl_zeta(i,jstr-1,know)+ &
506 & tl_ce*(zeta(i,jstr-1,know)+ &
507 & zeta(i,jstr ,know))- &
508 & (1.0_r8+ce)*tl_zeta(i,jstr,know)
509 ze=ze+cff2*cff3
510 tl_ze=tl_ze+cff2*tl_cff3+ &
511 & tl_cff2*cff3
512 END IF
513!^ vbar(i,Jstr,kout)=0.5_r8* &
514!^ & ((1.0_r8-Ce)*vbar(i,Jstr,know)+ &
515!^ & Ce*vbar(i,Jstr+1,know)+ &
516!^ & bry_val- &
517!^ & cff1*(Ze-BOUNDARY(ng)%zeta_south(i)))
518!^
519 tl_vbar(i,jstr,kout)=0.5_r8* &
520 & ((1.0_r8-ce)* &
521 & tl_vbar(i,jstr,know)- &
522 & tl_ce*(vbar(i,jstr ,know)- &
523 & vbar(i,jstr+1,know))+ &
524 & ce*tl_vbar(i,jstr+1,know)+ &
525 & tl_bry_val- &
526 & tl_cff1* &
527 & (ze-boundary(ng)%zeta_south(i))- &
528 & cff1*tl_ze)
529# ifdef ADJUST_BOUNDARY
530 IF (lobc(isouth,isvbar,ng)) THEN
531 tl_vbar(i,jstr,kout)=tl_vbar(i,jstr,kout)+ &
532 & 0.5_r8*cff1* &
533 & boundary(ng)%tl_zeta_south(i)
534 END IF
535# endif
536# ifdef MASKING
537!^ vbar(i,Jstr,kout)=vbar(i,Jstr,kout)* &
538!^ & GRID(ng)%vmask(i,Jstr)
539!^
540 tl_vbar(i,jstr,kout)=tl_vbar(i,jstr,kout)* &
541 & grid(ng)%vmask(i,jstr)
542# endif
543 END IF
544 END DO
545!
546! Southern edge, clamped boundary condition.
547!
548 ELSE IF (tl_lbc(isouth,isvbar,ng)%clamped) THEN
549 DO i=istr,iend
550 IF (lbc_apply(ng)%south(i)) THEN
551!^ vbar(i,Jstr,kout)=BOUNDARY(ng)%vbar_south(i)
552!^
553# ifdef ADJUST_BOUNDARY
554 IF (lobc(isouth,isvbar,ng)) THEN
555 tl_vbar(i,jstr,kout)=boundary(ng)%tl_vbar_south(i)
556 ELSE
557 tl_vbar(i,jstr,kout)=0.0_r8
558 END IF
559# else
560 tl_vbar(i,jstr,kout)=0.0_r8
561# endif
562# ifdef MASKING
563!^ vbar(i,Jstr,kout)=vbar(i,Jstr,kout)* &
564!^ & GRID(ng)%vmask(i,Jstr)
565!^
566 tl_vbar(i,jstr,kout)=tl_vbar(i,jstr,kout)* &
567 & grid(ng)%vmask(i,jstr)
568# endif
569 END IF
570 END DO
571!
572! Southern edge, gradient boundary condition.
573!
574 ELSE IF (tl_lbc(isouth,isvbar,ng)%gradient) THEN
575 DO i=istr,iend
576 IF (lbc_apply(ng)%south(i)) THEN
577!^ vbar(i,Jstr,kout)=vbar(i,Jstr+1,kout)
578!^
579 tl_vbar(i,jstr,kout)=tl_vbar(i,jstr+1,kout)
580# ifdef MASKING
581!^ vbar(i,Jstr,kout)=vbar(i,Jstr,kout)* &
582!^ & GRID(ng)%vmask(i,Jstr)
583!^
584 tl_vbar(i,jstr,kout)=tl_vbar(i,jstr,kout)* &
585 & grid(ng)%vmask(i,jstr)
586# endif
587 END IF
588 END DO
589!
590! Southern edge, reduced-physics boundary condition.
591!
592 ELSE IF (tl_lbc(isouth,isvbar,ng)%reduced) THEN
593 DO i=istr,iend
594 IF (lbc_apply(ng)%south(i)) THEN
595 IF (tl_lbc(isouth,isfsur,ng)%acquire) THEN
596!^ bry_pgr=-g*(zeta(i,Jstr,know)- &
597!^ & BOUNDARY(ng)%zeta_south(i))* &
598!^ & 0.5_r8*GRID(ng)%pn(i,Jstr)
599!^
600 tl_bry_pgr=-g*tl_zeta(i,jstr,know)* &
601 & 0.5_r8*grid(ng)%pn(i,jstr)
602# ifdef ADJUST_BOUNDARY
603 IF (lobc(inorth,isvbar,ng)) THEN
604 tl_bry_pgr=tl_bry_pgr+ &
605 & g*boundary(ng)%tl_zeta_south(i)* &
606 & 0.5_r8*grid(ng)%pn(i,jstr)
607 END IF
608# endif
609 ELSE
610!^ bry_pgr=-g*(zeta(i,Jstr,know)- &
611!^ & zeta(i,Jstr-1,know))* &
612!^ & 0.5_r8*(GRID(ng)%pn(i,Jstr-1)+ &
613!^ & GRID(ng)%pn(i,Jstr ))
614!^
615 tl_bry_pgr=-g*(tl_zeta(i,jstr ,know)- &
616 & tl_zeta(i,jstr-1,know))* &
617 & 0.5_r8*(grid(ng)%pn(i,jstr-1)+ &
618 & grid(ng)%pn(i,jstr ))
619 END IF
620# ifdef UV_COR
621!^ bry_cor=-0.125_r8*(ubar(i ,Jstr-1,know)+ &
622!^ & ubar(i+1,Jstr-1,know)+ &
623!^ & ubar(i ,Jstr ,know)+ &
624!^ & ubar(i+1,Jstr ,know))* &
625!^ & (GRID(ng)%f(i,Jstr-1)+ &
626!^ & GRID(ng)%f(i,Jstr ))
627!^
628 tl_bry_cor=-0.125_r8*(tl_ubar(i ,jstr-1,know)+ &
629 & tl_ubar(i+1,jstr-1,know)+ &
630 & tl_ubar(i ,jstr ,know)+ &
631 & tl_ubar(i+1,jstr ,know))* &
632 & (grid(ng)%f(i,jstr-1)+ &
633 & grid(ng)%f(i,jstr ))
634# else
635!^ bry_cor=0.0_r8
636!^
637 tl_bry_cor=0.0_r8
638# endif
639 cff=1.0_r8/(0.5_r8*(grid(ng)%h(i,jstr-1)+ &
640 & zeta(i,jstr-1,know)+ &
641 & grid(ng)%h(i,jstr )+ &
642 & zeta(i,jstr ,know)))
643 tl_cff=-cff*cff*0.5_r8*(grid(ng)%tl_h(i,jstr-1)+ &
644 & tl_zeta(i,jstr-1,know)+ &
645 & grid(ng)%tl_h(i,jstr )+ &
646 & tl_zeta(i,jstr ,know))
647!^ bry_str=cff*(FORCES(ng)%svstr(i,Jstr)- &
648!^ & FORCES(ng)%bvstr(i,Jstr))
649!^
650 tl_bry_str=tl_cff*(forces(ng)%svstr(i,jstr)- &
651 & forces(ng)%bvstr(i,jstr))+ &
652 & cff*(forces(ng)%tl_svstr(i,jstr)- &
653 & forces(ng)%tl_bvstr(i,jstr))
654
655!^ vbar(i,Jstr,kout)=vbar(i,Jstr,know)+ &
656!^ & dt2d*(bry_pgr+ &
657!^ & bry_cor+ &
658!^ & bry_str)
659!^
660 tl_vbar(i,jstr,kout)=tl_vbar(i,jstr,know)+ &
661 & dt2d*(tl_bry_pgr+ &
662 & tl_bry_cor+ &
663 & tl_bry_str)
664# ifdef MASKING
665!^ vbar(i,Jstr,kout)=vbar(i,Jstr,kout)* &
666!^ & GRID(ng)%vmask(i,Jstr)
667!^
668 tl_vbar(i,jstr,kout)=tl_vbar(i,jstr,kout)* &
669 & grid(ng)%vmask(i,jstr)
670# endif
671 END IF
672 END DO
673!
674! Southern edge, closed boundary condition.
675!
676 ELSE IF (tl_lbc(isouth,isvbar,ng)%closed) THEN
677 DO i=istr,iend
678 IF (lbc_apply(ng)%south(i)) THEN
679!^ vbar(i,Jstr,kout)=0.0_r8
680!^
681 tl_vbar(i,jstr,kout)=0.0_r8
682 END IF
683 END DO
684 END IF
685 END IF
686!
687!-----------------------------------------------------------------------
688! Lateral boundary conditions at the northern edge.
689!-----------------------------------------------------------------------
690!
691 IF (domain(ng)%Northern_Edge(tile)) THEN
692!
693! Northern edge, implicit upstream radiation condition.
694!
695 IF (tl_lbc(inorth,isvbar,ng)%radiation) THEN
696 IF (iic(ng).ne.0) THEN
697 DO i=istr,iend+1
698!^ grad(i,Jend+1)=vbar(i ,Jend+1,know)- &
699!^ & vbar(i-1,Jend+1,know)
700!^
701 tl_grad(i,jend+1)=0.0_r8
702 END DO
703 DO i=istr,iend
704 IF (lbc_apply(ng)%north(i)) THEN
705# if defined CELERITY_READ && defined FORWARD_READ
706 IF (tl_lbc(inorth,isvbar,ng)%nudging) THEN
707 IF (lnudgem2clm(ng)) THEN
708 obc_out=0.5_r8* &
709 & (clima(ng)%M2nudgcof(i,jend )+ &
710 & clima(ng)%M2nudgcof(i,jend+1))
711 obc_in =obcfac(ng)*obc_out
712 ELSE
713 obc_out=m2obc_out(ng,inorth)
714 obc_in =m2obc_in(ng,inorth)
715 END IF
716 IF (boundary(ng)%vbar_north_Ce(i).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# ifdef RADIATION_2D
724 cx=boundary(ng)%vbar_north_Cx(i)
725# else
726 cx=0.0_r8
727# endif
728 ce=boundary(ng)%vbar_north_Ce(i)
729 cff=boundary(ng)%vbar_north_C2(i)
730# endif
731!^ vbar(i,Jend+1,kout)=(cff*vbar(i,Jend+1,know)+ &
732!^ & Ce *vbar(i,Jend ,kout)- &
733!^ & MAX(Cx,0.0_r8)*grad(i ,Jend+1)- &
734!^ & MIN(Cx,0.0_r8)*grad(i+1,Jend+1))/ &
735!^ & (cff+Ce)
736!^
737 tl_vbar(i,jend+1,kout)=(cff*tl_vbar(i,jend+1,know)+ &
738 & ce *tl_vbar(i,jend ,kout)- &
739 & max(cx,0.0_r8)* &
740 & tl_grad(i ,jend+1)- &
741 & min(cx,0.0_r8)* &
742 & tl_grad(i+1,jend+1))/ &
743 & (cff+ce)
744
745 IF (tl_lbc(inorth,isvbar,ng)%nudging) THEN
746!^ vbar(i,Jend+1,kout)=vbar(i,Jend+1,kout)+ &
747!^ & tau*(BOUNDARY(ng)%vbar_north(i)- &
748!^ & vbar(i,Jend+1,know))
749!^
750 tl_vbar(i,jend+1,kout)=tl_vbar(i,jend+1,kout)- &
751 & tau*tl_vbar(i,jend+1,know)
752 END IF
753# ifdef MASKING
754!^ vbar(i,Jend+1,kout)=vbar(i,Jend+1,kout)* &
755!^ & GRID(ng)%vmask(i,Jend+1)
756!^
757 tl_vbar(i,jend+1,kout)=tl_vbar(i,jend+1,kout)* &
758 & grid(ng)%vmask(i,jend+1)
759# endif
760 END IF
761 END DO
762 END IF
763!
764! Northern edge, Flather boundary condition.
765!
766 ELSE IF (tl_lbc(inorth,isvbar,ng)%Flather) THEN
767 DO i=istr,iend
768 IF (lbc_apply(ng)%north(i)) THEN
769# if defined SSH_TIDES_NOT_YET && !defined UV_TIDES_NOT_YET
770 IF (tl_lbc(inorth,isfsur,ng)%acquire) THEN
771 bry_pgr=-g*(boundary(ng)%zeta_north(i)- &
772 & zeta(i,jend,know))* &
773 & 0.5_r8*grid(ng)%pn(i,jend)
774 tl_bry_pgr=g*tl_zeta(i,jend,know)* &
775 & 0.5_r8*grid(ng)%pn(i,jend)
776# ifdef ADJUST_BOUNDARY
777 IF (lobc(inorth,isvbar,ng)) THEN
778 tl_bry_pgr=tl_bry_pgr- &
779 & g*boundary(ng)%tl_zeta_north(i)* &
780 & 0.5_r8*grid(ng)%pn(i,jend)
781 END IF
782# endif
783 ELSE
784 bry_pgr=-g*(zeta(i,jend+1,know)- &
785 & zeta(i,jend ,know))* &
786 & 0.5_r8*(grid(ng)%pn(i,jend )+ &
787 & grid(ng)%pn(i,jend+1))
788 tl_bry_pgr=-g*(tl_zeta(i,jend+1,know)- &
789 & tl_zeta(i,jend ,know))* &
790 & 0.5_r8*(grid(ng)%pn(i,jend )+ &
791 & grid(ng)%pn(i,jend+1))
792 END IF
793# ifdef UV_COR
794 bry_cor=-0.125_r8*(ubar(i ,jend ,know)+ &
795 & ubar(i+1,jend ,know)+ &
796 & ubar(i ,jend+1,know)+ &
797 & ubar(i+1,jend+1,know))* &
798 & (grid(ng)%f(i,jend )+ &
799 & grid(ng)%f(i,jend+1))
800 tl_bry_cor=-0.125_r8*(tl_ubar(i ,jend ,know)+ &
801 & tl_ubar(i+1,jend ,know)+ &
802 & tl_ubar(i ,jend+1,know)+ &
803 & tl_ubar(i+1,jend+1,know))* &
804 & (grid(ng)%f(i,jend )+ &
805 & grid(ng)%f(i,jend+1))
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(i,jend )+ &
811 & zeta(i,jend ,know)+ &
812 & grid(ng)%h(i,jend+1)+ &
813 & zeta(i,jend+1,know)))
814 tl_cff1=-cff1*cff1*0.5_r8*(grid(ng)%tl_h(i,jend )+ &
815 & tl_zeta(i,jend ,know)+ &
816 & grid(ng)%tl_h(i,jend+1)+ &
817 & tl_zeta(i,jend+1,know))
818 bry_str=cff1*(forces(ng)%svstr(i,jend+1)- &
819 & forces(ng)%bvstr(i,jend+1))
820 tl_bry_str=tl_cff1*(forces(ng)%svstr(i,jend+1)- &
821 & forces(ng)%bvstr(i,jend+1))+ &
822 & cff1*(forces(ng)%tl_svstr(i,jend+1)- &
823 & forces(ng)%tl_bvstr(i,jend+1))
824 ce=1.0_r8/sqrt(g*0.5_r8*(grid(ng)%h(i,jend+1)+ &
825 & zeta(i,jend+1,know)+ &
826 & grid(ng)%h(i,jend )+ &
827 & zeta(i,jend ,know)))
828 tl_ce=-ce*ce*ce*0.25_r8*g*(grid(ng)%tl_h(i,jend+1)+ &
829 & tl_zeta(i,jend+1,know)+ &
830 & grid(ng)%tl_h(i,jend )+ &
831 & tl_zeta(i,jend ,know))
832 cff2=grid(ng)%on_v(i,jend+1)*ce
833 tl_cff2=grid(ng)%on_v(i,jend+1)*tl_ce
834!^ bry_val=vbar(i,Jend,know)+ &
835!^ & cff2*(bry_pgr+ &
836!^ & bry_cor+ &
837!^ & bry_str)
838!^
839 tl_bry_val=tl_vbar(i,jend,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)%vbar_north(i)
848!^
849# ifdef ADJUST_BOUNDARY
850 IF (lobc(inorth,isvbar,ng)) THEN
851 tl_bry_val=boundary(ng)%tl_vbar_north(i)
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(i,jend )+ &
860 & zeta(i,jend ,know)+ &
861 & grid(ng)%h(i,jend+1)+ &
862 & zeta(i,jend+1,know)))
863 tl_cff=-cff*cff*(0.5_r8*(grid(ng)%tl_h(i,jend )+ &
864 & tl_zeta(i,jend ,know)+ &
865 & grid(ng)%tl_h(i,jend+1)+ &
866 & tl_zeta(i,jend+1,know)))
867 ce=sqrt(g*cff)
868 tl_ce=0.5_r8*g*tl_cff/ce
869# if defined ATM_PRESS && defined PRESS_COMPENSATE
870!^ vbar(i,Jend+1,kout)=bry_val+ &
871!^ & Ce*(0.5_r8* &
872!^ & (zeta(i,Jend ,know)+ &
873!^ & zeta(i,Jend+1,know)+ &
874!^ & fac*(FORCES(ng)%Pair(i,Jend )+ &
875!^ & FORCES(ng)%Pair(i,Jend+1)- &
876!^ & 2.0_r8*OneAtm))- &
877!^ & BOUNDARY(ng)%zeta_north(i))
878!^
879 tl_vbar(i,jend+1,kout)=tl_bry_val+ &
880 & tl_ce* &
881 & (0.5_r8* &
882 & (zeta(i,jend ,know)+ &
883 & zeta(i,jend+1,know)+ &
884 & fac*(forces(ng)%Pair(i,jend )+ &
885 & forces(ng)%Pair(i,jend+1)- &
886 & 2.0_r8*oneatm))- &
887 & boundary(ng)%zeta_north(i))+ &
888 & ce* &
889 & (0.5_r8*(tl_zeta(i,jend ,know)+ &
890 & tl_zeta(i,jend+1,know)))
891# else
892!^ vbar(i,Jend+1,kout)=bry_val+ &
893!^ & Ce*(0.5_r8*(zeta(i,Jend ,know)+ &
894!^ & zeta(i,Jend+1,know))- &
895!^ & BOUNDARY(ng)%zeta_north(i))
896!^
897 tl_vbar(i,jend+1,kout)=tl_bry_val+ &
898 & tl_ce* &
899 & (0.5_r8*(zeta(i,jend ,know)+ &
900 & zeta(i,jend+1,know))- &
901 & boundary(ng)%zeta_north(i))+ &
902 & ce* &
903 & (0.5_r8*(tl_zeta(i,jend ,know)+ &
904 & tl_zeta(i,jend+1,know)))
905# endif
906# ifdef ADJUST_BOUNDARY
907 IF (lobc(inorth,isvbar,ng)) THEN
908 tl_vbar(i,jend+1,kout)=tl_vbar(i,jend+1,kout)- &
909 & ce*boundary(ng)%tl_zeta_north(i)
910 END IF
911# endif
912# ifdef MASKING
913!^ vbar(i,Jend+1,kout)=vbar(i,Jend+1,kout)* &
914!^ & GRID(ng)%vmask(i,Jend+1)
915!^
916 tl_vbar(i,jend+1,kout)=tl_vbar(i,jend+1,kout)* &
917 & grid(ng)%vmask(i,jend+1)
918# endif
919 END IF
920 END DO
921!
922! Northern edge, Shchepetkin boundary condition (Maison et al., 2010).
923!
924 ELSE IF (lbc(inorth,isvbar,ng)%Shchepetkin) THEN
925 DO i=istr,iend
926 IF (lbc_apply(ng)%north(i)) THEN
927# if defined SSH_TIDES_NOT_YET && !defined UV_TIDES_NOT_YET
928 IF (tl_lbc(inorth,isfsur,ng)%acquire) THEN
929 bry_pgr=-g*(boundary(ng)%zeta_north(i)- &
930 & zeta(i,jend,know))* &
931 & 0.5_r8*grid(ng)%pn(i,jend)
932 tl_bry_pgr=g*tl_zeta(i,jend,know)* &
933 & 0.5_r8*grid(ng)%pn(i,jend)
934# ifdef ADJUST_BOUNDARY
935 IF (lobc(inorth,isvbar,ng)) THEN
936 tl_bry_pgr=tl_bry_pgr- &
937 & g*boundary(ng)%tl_zeta_north(i)* &
938 & 0.5_r8*grid(ng)%pn(i,jend)
939 END IF
940# endif
941 ELSE
942 bry_pgr=-g*(zeta(i,jend+1,know)- &
943 & zeta(i,jend ,know))* &
944 & 0.5_r8*(grid(ng)%pn(i,jend )+ &
945 & grid(ng)%pn(i,jend+1))
946 tl_bry_pgr=-g*(tl_zeta(i,jend+1,know)- &
947 & tl_zeta(i,jend ,know))* &
948 & 0.5_r8*(grid(ng)%pn(i,jend )+ &
949 & grid(ng)%pn(i,jend+1))
950 END IF
951# ifdef UV_COR
952 bry_cor=-0.125_r8*(ubar(i ,jend ,know)+ &
953 & ubar(i+1,jend ,know)+ &
954 & ubar(i ,jend+1,know)+ &
955 & ubar(i+1,jend+1,know))* &
956 & (grid(ng)%f(i,jend )+ &
957 & grid(ng)%f(i,jend+1))
958 tl_bry_cor=-0.125_r8*(tl_ubar(i ,jend ,know)+ &
959 & tl_ubar(i+1,jend ,know)+ &
960 & tl_ubar(i ,jend+1,know)+ &
961 & tl_ubar(i+1,jend+1,know))* &
962 & (grid(ng)%f(i,jend )+ &
963 & grid(ng)%f(i,jend+1))
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(i,jend )+ &
969 & zeta(i,jend ,know)+ &
970 & grid(ng)%h(i,jend+1)+ &
971 & zeta(i,jend+1,know)))
972 tl_cff1=-cff1*cff1*0.5_r8*(grid(ng)%tl_h(i,jend )+ &
973 & tl_zeta(i,jend ,know)+ &
974 & grid(ng)%tl_h(i,jend+1)+ &
975 & tl_zeta(i,jend+1,know))
976 bry_str=cff1*(forces(ng)%svstr(i,jend+1)- &
977 & forces(ng)%bvstr(i,jend+1))
978 tl_bry_str=tl_cff1*(forces(ng)%svstr(i,jend+1)- &
979 & forces(ng)%bvstr(i,jend+1))+ &
980 & cff1*(forces(ng)%tl_svstr(i,jend+1)- &
981 & forces(ng)%tl_bvstr(i,jend+1))
982 ce=1.0_r8/sqrt(g*0.5_r8*(grid(ng)%h(i,jend+1)+ &
983 & zeta(i,jend+1,know)+ &
984 & grid(ng)%h(i,jend )+ &
985 & zeta(i,jend ,know)))
986 tl_ce=-ce*ce*ce*0.25_r8*g*(grid(ng)%tl_h(i,jend+1)+ &
987 & tl_zeta(i,jend+1,know)+ &
988 & grid(ng)%tl_h(i,jend )+ &
989 & tl_zeta(i,jend ,know))
990 cff2=grid(ng)%on_v(i,jend+1)*ce
991 tl_cff2=grid(ng)%on_v(i,jend+1)*tl_ce
992!^ bry_val=vbar(i,Jend,know)+ &
993!^ & cff2*(bry_pgr+ &
994!^ & bry_cor+ &
995!^ & bry_str)
996!^
997 tl_bry_val=tl_vbar(i,jend,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)%vbar_north(i)
1006!^
1007# ifdef ADJUST_BOUNDARY
1008 IF (lobc(inorth,isvbar,ng)) THEN
1009 tl_bry_val=boundary(ng)%tl_vbar_north(i)
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(i,jend )+ &
1019 & zeta(i,jend ,know)+ &
1020 & grid(ng)%h(i,jend+1)+ &
1021 & zeta(i,jend+1,know))
1022 tl_cff=0.5_r8*(grid(ng)%tl_h(i,jend )+ &
1023 & tl_zeta(i,jend ,know)+ &
1024 & grid(ng)%tl_h(i,jend+1)+ &
1025 & tl_zeta(i,jend+1,know))
1026# else
1027
1028 cff=0.5_r8*(grid(ng)%h(i,jend )+ &
1029 & grid(ng)%h(i,jend+1))
1030 tl_cff=0.5_r8*(grid(ng)%tl_h(i,jend )+ &
1031 & grid(ng)%tl_h(i,jend+1))
1032# endif
1033 cff1=sqrt(g/cff)
1034 tl_cff1=-0.5_r8*cff1*tl_cff/cff
1035 ce=dt2d*cff1*cff*0.5_r8*(grid(ng)%pn(i,jend )+ &
1036 & grid(ng)%pn(i,jend+1))
1037 tl_ce=dt2d*0.5_r8*(grid(ng)%pn(i,jend )+ &
1038 & grid(ng)%pn(i,jend+1))* &
1039 & (cff1*tl_cff+ &
1040 & tl_cff1*cff)
1041 ze=(0.5_r8+ce)*zeta(i,jend ,know)+ &
1042 & (0.5_r8-ce)*zeta(i,jend+1,know)
1043 tl_ze=(0.5_r8+ce)*tl_zeta(i,jend ,know)+ &
1044 & (0.5_r8-ce)*tl_zeta(i,jend+1,know)+ &
1045 & tl_ce*(zeta(i,jend ,know)- &
1046 & zeta(i,jend+1,know))
1047 IF (ce.gt.co) THEN
1048 cff2=(1.0_r8-co/ce)**2
1049 tl_cff2=2.0_r8*cff2*co*tl_ce/(ce*ce)
1050 cff3=zeta(i,jend,kout)+ &
1051 & ce*zeta(i,jend+1,know)- &
1052 & (1.0_r8+ce)*zeta(i,jend,know)
1053 tl_cff3=tl_zeta(i,jend,kout)+ &
1054 & ce*tl_zeta(i,jend+1,know)+ &
1055 & tl_ce*(zeta(i,jend ,know)+ &
1056 & zeta(i,jend+1,know))- &
1057 & (1.0_r8+ce)*tl_zeta(i,jend,know)
1058 ze=ze+cff2*cff3
1059 tl_ze=tl_ze+cff2*tl_cff3+ &
1060 & tl_cff2*cff3
1061 END IF
1062!^ vbar(i,Jend+1,kout)=0.5_r8* &
1063!^ & ((1.0_r8-Ce)*vbar(i,Jend+1,know)+ &
1064!^ & Ce*vbar(i,Jend,know)+ &
1065!^ & bry_val+ &
1066!^ & cff1*(Ze-BOUNDARY(ng)%zeta_north(i)))
1067!^
1068 tl_vbar(i,jend+1,kout)=0.5_r8* &
1069 & ((1.0_r8-ce)* &
1070 & tl_vbar(i,jend+1,know)+ &
1071 & tl_ce*(vbar(i,jend ,know)- &
1072 & vbar(i,jend+1,know))+ &
1073 & ce*tl_vbar(i,jend,know)+ &
1074 & tl_bry_val+ &
1075 & tl_cff1* &
1076 & (ze-boundary(ng)%zeta_north(i))- &
1077 & cff1*tl_ze)
1078# ifdef ADJUST_BOUNDARY
1079 IF (lobc(inorth,isvbar,ng)) THEN
1080 tl_vbar(i,jend+1,kout)=tl_vbar(i,jend+1,kout)- &
1081 & 0.5_r8*cff1* &
1082 & ce*boundary(ng)%tl_zeta_north(i)
1083 END IF
1084# endif
1085# ifdef MASKING
1086!^ vbar(i,Jend+1,kout)=vbar(i,Jend+1,kout)* &
1087!^ & GRID(ng)%vmask(i,Jend+1)
1088!^
1089 tl_vbar(i,jend+1,kout)=tl_vbar(i,jend+1,kout)* &
1090 & grid(ng)%vmask(i,jend+1)
1091# endif
1092 END IF
1093 END DO
1094!
1095! Northern edge, clamped boundary condition.
1096!
1097 ELSE IF (tl_lbc(inorth,isvbar,ng)%clamped) THEN
1098 DO i=istr,iend
1099 IF (lbc_apply(ng)%north(i)) THEN
1100!^ vbar(i,Jend+1,kout)=BOUNDARY(ng)%vbar_north(i)
1101!^
1102# ifdef ADJUST_BOUNDARY
1103 IF (lobc(inorth,isvbar,ng)) THEN
1104 tl_vbar(i,jend+1,kout)=boundary(ng)%tl_vbar_north(i)
1105 ELSE
1106 tl_vbar(i,jend+1,kout)=0.0_r8
1107 END IF
1108# else
1109 tl_vbar(i,jend+1,kout)=0.0_r8
1110# endif
1111# ifdef MASKING
1112!^ vbar(i,Jend+1,kout)=vbar(i,Jend+1,kout)* &
1113!^ & GRID(ng)%vmask(i,Jend+1)
1114!^
1115 tl_vbar(i,jend+1,kout)=tl_vbar(i,jend+1,kout)* &
1116 & grid(ng)%vmask(i,jend+1)
1117# endif
1118 END IF
1119 END DO
1120!
1121! Northern edge, gradient boundary condition.
1122!
1123 ELSE IF (tl_lbc(inorth,isvbar,ng)%gradient) THEN
1124 DO i=istr,iend
1125 IF (lbc_apply(ng)%north(i)) THEN
1126!^ vbar(i,Jend+1,kout)=vbar(i,Jend,kout)
1127!^
1128 tl_vbar(i,jend+1,kout)=tl_vbar(i,jend,kout)
1129# ifdef MASKING
1130!^ vbar(i,Jend+1,kout)=vbar(i,Jend+1,kout)* &
1131!^ & GRID(ng)%vmask(i,Jend+1)
1132!^
1133 tl_vbar(i,jend+1,kout)=tl_vbar(i,jend+1,kout)* &
1134 & grid(ng)%vmask(i,jend+1)
1135# endif
1136 END IF
1137 END DO
1138!
1139! Northern edge, reduced-physics boundary condition.
1140!
1141 ELSE IF (tl_lbc(inorth,isvbar,ng)%reduced) THEN
1142 DO i=istr,iend
1143 IF (lbc_apply(ng)%north(i)) THEN
1144 IF (tl_lbc(inorth,isfsur,ng)%acquire) THEN
1145!^ bry_pgr=-g*(BOUNDARY(ng)%zeta_north(i)- &
1146!^ & zeta(i,Jend,know))* &
1147!^ & 0.5_r8*GRID(ng)%pn(i,Jend)
1148!^
1149 tl_bry_pgr=g*tl_zeta(i,jend,know)* &
1150 & 0.5_r8*grid(ng)%pn(i,jend)
1151# ifdef ADJUST_BOUNDARY
1152 IF (lobc(inorth,isvbar,ng)) THEN
1153 tl_bry_pgr=tl_bry_pgr- &
1154 & g*boundary(ng)%tl_zeta_north(i)* &
1155 & 0.5_r8*grid(ng)%pn(i,jend)
1156 END IF
1157# endif
1158 ELSE
1159!^ bry_pgr=-g*(zeta(i,Jend+1,know)- &
1160!^ & zeta(i,Jend ,know))* &
1161!^ & 0.5_r8*(GRID(ng)%pn(i,Jend )+ &
1162!^ & GRID(ng)%pn(i,Jend+1))
1163!^
1164 tl_bry_pgr=-g*(tl_zeta(i,jend+1,know)- &
1165 & tl_zeta(i,jend ,know))* &
1166 & 0.5_r8*(grid(ng)%pn(i,jend )+ &
1167 & grid(ng)%pn(i,jend+1))
1168 END IF
1169# ifdef UV_COR
1170!^ bry_cor=-0.125_r8*(ubar(i ,Jend ,know)+ &
1171!^ & ubar(i+1,Jend ,know)+ &
1172!^ & ubar(i ,Jend+1,know)+ &
1173!^ & ubar(i+1,Jend+1,know))* &
1174!^ & (GRID(ng)%f(i,Jend )+ &
1175!^ & GRID(ng)%f(i,Jend+1))
1176!^
1177 tl_bry_cor=-0.125_r8*(tl_ubar(i ,jend ,know)+ &
1178 & tl_ubar(i+1,jend ,know)+ &
1179 & tl_ubar(i ,jend+1,know)+ &
1180 & tl_ubar(i+1,jend+1,know))* &
1181 & (grid(ng)%f(i,jend )+ &
1182 & grid(ng)%f(i,jend+1))
1183# else
1184!^ bry_cor=0.0_r8
1185!^
1186 tl_bry_cor=0.0_r8
1187# endif
1188 cff=1.0_r8/(0.5_r8*(grid(ng)%h(i,jend )+ &
1189 & zeta(i,jend ,know)+ &
1190 & grid(ng)%h(i,jend+1)+ &
1191 & zeta(i,jend+1,know)))
1192 tl_cff=-cff*cff*0.5_r8*(grid(ng)%tl_h(i,jend )+ &
1193 & tl_zeta(i,jend ,know)+ &
1194 & grid(ng)%tl_h(i,jend+1)+ &
1195 & tl_zeta(i,jend+1,know))
1196!^ bry_str=cff*(FORCES(ng)%svstr(i,Jend+1)- &
1197!^ & FORCES(ng)%bvstr(i,Jend+1))
1198!^
1199 tl_bry_str=tl_cff*(forces(ng)%svstr(i,jend+1)- &
1200 & forces(ng)%bvstr(i,jend+1))+ &
1201 & cff*(forces(ng)%tl_svstr(i,jend+1)- &
1202 & forces(ng)%tl_bvstr(i,jend+1))
1203!^ vbar(i,Jend+1,kout)=vbar(i,Jend+1,know)+ &
1204!^ & dt2d*(bry_pgr+ &
1205!^ & bry_cor+ &
1206!^ & bry_str)
1207!^
1208 tl_vbar(i,jend+1,kout)=tl_vbar(i,jend+1,know)+ &
1209 & dt2d*(tl_bry_pgr+ &
1210 & tl_bry_cor+ &
1211 & tl_bry_str)
1212# ifdef MASKING
1213!^ vbar(i,Jend+1,kout)=vbar(i,Jend+1,kout)* &
1214!^ & GRID(ng)%vmask(i,Jend+1)
1215!^
1216 tl_vbar(i,jend+1,kout)=tl_vbar(i,jend+1,kout)* &
1217 & grid(ng)%vmask(i,jend+1)
1218# endif
1219 END IF
1220 END DO
1221!
1222! Northern edge, closed boundary condition.
1223!
1224 ELSE IF (tl_lbc(inorth,isvbar,ng)%closed) THEN
1225 DO i=istr,iend
1226 IF (lbc_apply(ng)%north(i)) THEN
1227!^ vbar(i,Jend+1,kout)=0.0_r8
1228!^
1229 tl_vbar(i,jend+1,kout)=0.0_r8
1230 END IF
1231 END DO
1232 END IF
1233 END IF
1234!
1235!-----------------------------------------------------------------------
1236! Lateral boundary conditions at the western edge.
1237!-----------------------------------------------------------------------
1238!
1239 IF (domain(ng)%Western_Edge(tile)) THEN
1240!
1241! Western edge, implicit upstream radiation condition.
1242!
1243 IF (tl_lbc(iwest,isvbar,ng)%radiation) THEN
1244 IF (iic(ng).ne.0) THEN
1245 DO j=jstrv-1,jend
1246!^ grad(Istr-1,j)=vbar(Istr-1,j+1,know)- &
1247!^ & vbar(Istr-1,j ,know)
1248!^
1249 tl_grad(istr-1,j)=0.0_r8
1250 END DO
1251 DO j=jstrv,jend
1252 IF (lbc_apply(ng)%west(j)) THEN
1253# if defined CELERITY_READ && defined FORWARD_READ
1254 IF (tl_lbc(iwest,isvbar,ng)%nudging) THEN
1255 IF (lnudgem2clm(ng)) THEN
1256 obc_out=0.5_r8* &
1257 & (clima(ng)%M2nudgcof(istr-1,j-1)+ &
1258 & clima(ng)%M2nudgcof(istr-1,j ))
1259 obc_in =obcfac(ng)*obc_out
1260 ELSE
1261 obc_out=m2obc_out(ng,iwest)
1262 obc_in =m2obc_in(ng,iwest)
1263 END IF
1264 IF (boundary(ng)%vbar_west_Cx(j).lt.0.0_r8) THEN
1265 tau=obc_in
1266 ELSE
1267 tau=obc_out
1268 END IF
1269 tau=tau*dt2d
1270 END IF
1271 cx=boundary(ng)%vbar_west_Cx(j)
1272# ifdef RADIATION_2D
1273 ce=boundary(ng)%vbar_west_Ce(j)
1274# else
1275 ce=0.0_r8
1276# endif
1277 cff=boundary(ng)%vbar_west_C2(j)
1278# endif
1279!^ vbar(Istr-1,j,kout)=(cff*vbar(Istr-1,j,know)+ &
1280!^ & Cx *vbar(Istr ,j,kout)- &
1281!^ & MAX(Ce,0.0_r8)*grad(Istr-1,j-1)- &
1282!^ & MIN(Ce,0.0_r8)*grad(Istr-1,j ))/ &
1283!^ & (cff+Cx)
1284!^
1285 tl_vbar(istr-1,j,kout)=(cff*tl_vbar(istr-1,j,know)+ &
1286 & cx *tl_vbar(1,j,kout)- &
1287 & max(ce,0.0_r8)* &
1288 & tl_grad(istr-1,j-1)- &
1289 & min(ce,0.0_r8)* &
1290 & tl_grad(istr-1,j ))/ &
1291 & (cff+cx)
1292
1293 IF (tl_lbc(iwest,isvbar,ng)%nudging) THEN
1294!^ vbar(Istr-1,j,kout)=vbar(Istr-1,j,kout)+ &
1295!^ & tau*(BOUNDARY(ng)%vbar_west(j)- &
1296!^ & vbar(Istr-1,j,know))
1297!^
1298 tl_vbar(istr-1,j,kout)=tl_vbar(istr-1,j,kout)- &
1299 & tau*tl_vbar(1,j,know)
1300 END IF
1301# ifdef MASKING
1302!^ vbar(Istr-1,j,kout)=vbar(Istr-1,j,kout)* &
1303!^ & GRID(ng)%vmask(Istr-1,j)
1304!^
1305 tl_vbar(istr-1,j,kout)=tl_vbar(istr-1,j,kout)* &
1306 & grid(ng)%vmask(istr-1,j)
1307# endif
1308 END IF
1309 END DO
1310 END IF
1311!
1312! Western edge, Chapman boundary condition.
1313!
1314 ELSE IF (tl_lbc(iwest,isvbar,ng)%Flather.or. &
1315 & tl_lbc(iwest,isvbar,ng)%reduced.or. &
1316 & tl_lbc(iwest,isvbar,ng)%Shchepetkin) THEN
1317 DO j=jstrv,jend
1318 IF (lbc_apply(ng)%west(j)) THEN
1319 cff=dt2d*0.5_r8*(grid(ng)%pm(istr,j-1)+ &
1320 & grid(ng)%pm(istr,j ))
1321 cff1=sqrt(g*0.5_r8*(grid(ng)%h(istr,j-1)+ &
1322 & zeta(istr,j-1,know)+ &
1323 & grid(ng)%h(istr,j )+ &
1324 & zeta(istr,j ,know)))
1325 tl_cff1=0.25_r8*g*(grid(ng)%tl_h(istr,j-1)+ &
1326 & tl_zeta(istr,j-1,know)+ &
1327 & grid(ng)%tl_h(istr,j )+ &
1328 & tl_zeta(istr,j ,know))/cff1
1329 cx=cff*cff1
1330 tl_cx=cff*tl_cff1
1331 cff2=1.0_r8/(1.0_r8+cx)
1332 tl_cff2=-cff2*cff2*tl_cx
1333!^ vbar(Istr-1,j,kout)=cff2*(vbar(Istr-1,j,know)+ &
1334!^ & Cx*vbar(Istr,j,kout))
1335!^
1336 tl_vbar(istr-1,j,kout)=tl_cff2*(vbar(istr-1,j,know)+ &
1337 & cx*vbar(istr,j,kout))+ &
1338 & cff2*(tl_vbar(istr-1,j,know)+ &
1339 & tl_cx*vbar(istr,j,kout)+ &
1340 & cx*tl_vbar(istr,j,kout))
1341# ifdef MASKING
1342!^ vbar(Istr-1,j,kout)=vbar(Istr-1,j,kout)* &
1343!^ & GRID(ng)%vmask(Istr-1,j)
1344!^
1345 tl_vbar(istr-1,j,kout)=tl_vbar(istr-1,j,kout)* &
1346 & grid(ng)%vmask(istr-1,j)
1347# endif
1348 END IF
1349 END DO
1350!
1351! Western edge, clamped boundary condition.
1352!
1353 ELSE IF (tl_lbc(iwest,isvbar,ng)%clamped) THEN
1354 DO j=jstrv,jend
1355 IF (lbc_apply(ng)%west(j)) THEN
1356!^ vbar(Istr-1,j,kout)=BOUNDARY(ng)%vbar_west(j)
1357!^
1358# ifdef ADJUST_BOUNDARY
1359 IF (lobc(iwest,isvbar,ng)) THEN
1360 tl_vbar(istr-1,j,kout)=boundary(ng)%tl_vbar_west(j)
1361 ELSE
1362 tl_vbar(istr-1,j,kout)=0.0_r8
1363 END IF
1364# else
1365 tl_vbar(istr-1,j,kout)=0.0_r8
1366# endif
1367# ifdef MASKING
1368!^ vbar(Istr-1,j,kout)=vbar(Istr-1,j,kout)* &
1369!^ & GRID(ng)%vmask(Istr-1,j)
1370!^
1371 tl_vbar(istr-1,j,kout)=tl_vbar(istr-1,j,kout)* &
1372 & grid(ng)%vmask(istr-1,j)
1373# endif
1374 END IF
1375 END DO
1376!
1377! Western edge, gradient boundary condition.
1378!
1379 ELSE IF (tl_lbc(iwest,isvbar,ng)%gradient) THEN
1380 DO j=jstrv,jend
1381 IF (lbc_apply(ng)%west(j)) THEN
1382!^ vbar(Istr-1,j,kout)=vbar(Istr,j,kout)
1383!^
1384 tl_vbar(istr-1,j,kout)=tl_vbar(istr,j,kout)
1385# ifdef MASKING
1386!^ vbar(Istr-1,j,kout)=vbar(Istr-1,j,kout)* &
1387!^ & GRID(ng)%vmask(Istr-1,j)
1388!^
1389 tl_vbar(istr-1,j,kout)=tl_vbar(istr-1,j,kout)* &
1390 & grid(ng)%vmask(istr-1,j)
1391# endif
1392 END IF
1393 END DO
1394!
1395! Western edge, closed boundary condition: free slip (gamma2=1) or
1396! no slip (gamma2=-1).
1397!
1398 ELSE IF (tl_lbc(iwest,isvbar,ng)%closed) THEN
1399 IF (nsperiodic(ng)) THEN
1400 jmin=jstrv
1401 jmax=jend
1402 ELSE
1403 jmin=jstr
1404 jmax=jendr
1405 END IF
1406 DO j=jmin,jmax
1407 IF (lbc_apply(ng)%west(j)) THEN
1408!^ vbar(Istr-1,j,kout)=gamma2(ng)*vbar(Istr,j,kout)
1409!^
1410 tl_vbar(istr-1,j,kout)=gamma2(ng)*tl_vbar(istr,j,kout)
1411# ifdef MASKING
1412!^ vbar(Istr-1,j,kout)=vbar(Istr-1,j,kout)* &
1413!^ & GRID(ng)%vmask(Istr-1,j)
1414!^
1415 tl_vbar(istr-1,j,kout)=tl_vbar(istr-1,j,kout)* &
1416 & grid(ng)%vmask(istr-1,j)
1417# endif
1418 END IF
1419 END DO
1420 END IF
1421 END IF
1422!
1423!-----------------------------------------------------------------------
1424! Lateral boundary conditions at the eastern edge.
1425!-----------------------------------------------------------------------
1426!
1427 IF (domain(ng)%Eastern_Edge(tile)) THEN
1428!
1429! Eastern edge, implicit upstream radiation condition.
1430!
1431 IF (tl_lbc(ieast,isvbar,ng)%radiation) THEN
1432 IF (iic(ng).ne.0) THEN
1433 DO j=jstrv-1,jend
1434!^ grad(Iend+1,j)=vbar(Iend+1,j+1,know)- &
1435!^ & vbar(Iend+1,j ,know)
1436!^
1437 tl_grad(iend+1,j)=0.0_r8
1438 END DO
1439 DO j=jstrv,jend
1440 IF (lbc_apply(ng)%east(j)) THEN
1441# if defined CELERITY_READ && defined FORWARD_READ
1442 IF (tl_lbc(ieast,isvbar,ng)%nudging) THEN
1443 IF (lnudgem2clm(ng)) THEN
1444 obc_out=0.5_r8* &
1445 & (clima(ng)%M2nudgcof(iend+1,j-1)+ &
1446 & clima(ng)%M2nudgcof(iend+1,j ))
1447 obc_in =obcfac(ng)*obc_out
1448 ELSE
1449 obc_out=m2obc_out(ng,ieast)
1450 obc_in =m2obc_in(ng,ieast)
1451 END IF
1452 IF (boundary(ng)%vbar_east_Cx(j).lt.0.0_r8) THEN
1453 tau=obc_in
1454 ELSE
1455 tau=obc_out
1456 END IF
1457 tau=tau*dt2d
1458 END IF
1459 cx=boundary(ng)%vbar_east_Cx(j)
1460# ifdef RADIATION_2D
1461 ce=boundary(ng)%vbar_east_Ce(j)
1462# else
1463 ce=0.0_r8
1464# endif
1465 cff=boundary(ng)%vbar_east_C2(j)
1466# endif
1467!^ vbar(Iend+1,j,kout)=(cff*vbar(Iend+1,j,know)+ &
1468!^ & Cx *vbar(Iend ,j,kout)- &
1469!^ & MAX(Ce,0.0_r8)*grad(Iend+1,j-1)- &
1470!^ & MIN(Ce,0.0_r8)*grad(Iend+1,j ))/ &
1471!^ & (cff+Cx)
1472!^
1473 tl_vbar(iend+1,j,kout)=(cff*tl_vbar(iend+1,j,know)+ &
1474 & cx *tl_vbar(iend ,j,kout)- &
1475 & max(ce,0.0_r8)* &
1476 & tl_grad(iend+1,j-1)- &
1477 & min(ce,0.0_r8)* &
1478 & tl_grad(iend+1,j ))/ &
1479 & (cff+cx)
1480
1481 IF (tl_lbc(ieast,isvbar,ng)%nudging) THEN
1482!^ vbar(Iend+1,j,kout)=vbar(Iend+1,j,kout)+ &
1483!^ & tau*(BOUNDARY(ng)%vbar_east(j)- &
1484!^ & vbar(Iend+1,j,know))
1485!^
1486 tl_vbar(iend+1,j,kout)=tl_vbar(iend+1,j,kout)- &
1487 & tau*tl_vbar(iend+1,j,know)
1488 END IF
1489# ifdef MASKING
1490!^ vbar(Iend+1,j,kout)=vbar(Iend+1,j,kout)* &
1491!^ & GRID(ng)%vmask(Iend+1,j)
1492!^
1493 tl_vbar(iend+1,j,kout)=tl_vbar(iend+1,j,kout)* &
1494 & grid(ng)%vmask(iend+1,j)
1495# endif
1496 END IF
1497 END DO
1498 END IF
1499!
1500! Eastern edge, Chapman boundary condition.
1501!
1502 ELSE IF (tl_lbc(ieast,isvbar,ng)%Flather.or. &
1503 & tl_lbc(ieast,isvbar,ng)%reduced.or. &
1504 & tl_lbc(ieast,isvbar,ng)%Shchepetkin) THEN
1505 DO j=jstrv,jend
1506 IF (lbc_apply(ng)%east(j)) THEN
1507 cff=dt2d*0.5_r8*(grid(ng)%pm(iend,j-1)+ &
1508 & grid(ng)%pm(iend,j ))
1509 cff1=sqrt(g*0.5_r8*(grid(ng)%h(iend,j-1)+ &
1510 & zeta(iend,j-1,know)+ &
1511 & grid(ng)%h(iend,j )+ &
1512 & zeta(iend,j ,know)))
1513 tl_cff1=0.25_r8*g*(grid(ng)%tl_h(iend,j-1)+ &
1514 & tl_zeta(iend,j-1,know)+ &
1515 & grid(ng)%tl_h(iend,j )+ &
1516 & tl_zeta(iend,j ,know))/cff1
1517 cx=cff*cff1
1518 tl_cx=cff*tl_cff1
1519 cff2=1.0_r8/(1.0_r8+cx)
1520 tl_cff2=-cff2*cff2*tl_cx
1521!^ vbar(Iend+1,j,kout)=cff2*(vbar(Iend+1,j,know)+ &
1522!^ & Cx*vbar(Iend,j,kout))
1523!^
1524 tl_vbar(iend+1,j,kout)=tl_cff2*(vbar(iend+1,j,know)+ &
1525 & cx*vbar(iend,j,kout))+ &
1526 & cff2*(tl_vbar(iend+1,j,know)+ &
1527 & tl_cx*vbar(iend,j,kout)+ &
1528 & cx*tl_vbar(iend,j,kout))
1529# ifdef MASKING
1530!^ vbar(Iend+1,j,kout)=vbar(Iend+1,j,kout)* &
1531!^ & GRID(ng)%vmask(Iend+1,j)
1532!^
1533 tl_vbar(iend+1,j,kout)=tl_vbar(iend+1,j,kout)* &
1534 & grid(ng)%vmask(iend+1,j)
1535# endif
1536 END IF
1537 END DO
1538!
1539! Eastern edge, clamped boundary condition.
1540!
1541 ELSE IF (tl_lbc(ieast,isvbar,ng)%clamped) THEN
1542 DO j=jstrv,jend
1543 IF (lbc_apply(ng)%east(j)) THEN
1544!^ vbar(Iend+1,j,kout)=BOUNDARY(ng)%vbar_east(j)
1545!^
1546# ifdef ADJUST_BOUNDARY
1547 IF (lobc(ieast,isvbar,ng)) THEN
1548 tl_vbar(iend+1,j,kout)=boundary(ng)%tl_vbar_east(j)
1549 ELSE
1550 tl_vbar(iend+1,j,kout)=0.0_r8
1551 END IF
1552# else
1553 tl_vbar(iend+1,j,kout)=0.0_r8
1554# endif
1555# ifdef MASKING
1556!^ vbar(Iend+1,j,kout)=vbar(Iend+1,j,kout)* &
1557!^ & GRID(ng)%vmask(Iend+1,j)
1558!^
1559 tl_vbar(iend+1,j,kout)=tl_vbar(iend+1,j,kout)* &
1560 & grid(ng)%vmask(iend+1,j)
1561# endif
1562 END IF
1563 END DO
1564!
1565! Eastern edge, gradient boundary condition.
1566!
1567 ELSE IF (tl_lbc(ieast,isvbar,ng)%gradient) THEN
1568 DO j=jstrv,jend
1569 IF (lbc_apply(ng)%east(j)) THEN
1570!^ vbar(Iend+1,j,kout)=vbar(Iend,j,kout)
1571!^
1572 tl_vbar(iend+1,j,kout)=tl_vbar(iend,j,kout)
1573# ifdef MASKING
1574!^ vbar(Iend+1,j,kout)=vbar(Iend+1,j,kout)* &
1575!^ & GRID(ng)%vmask(Iend+1,j)
1576!^
1577 tl_vbar(iend+1,j,kout)=tl_vbar(iend+1,j,kout)* &
1578 & grid(ng)%vmask(iend+1,j)
1579# endif
1580 END IF
1581 END DO
1582!
1583! Eastern edge, closed boundary condition: free slip (gamma2=1) or
1584! no slip (gamma2=-1).
1585!
1586 ELSE IF (tl_lbc(ieast,isvbar,ng)%closed) THEN
1587 IF (nsperiodic(ng)) THEN
1588 jmin=jstrv
1589 jmax=jend
1590 ELSE
1591 jmin=jstr
1592 jmax=jendr
1593 END IF
1594 DO j=jmin,jmax
1595 IF (lbc_apply(ng)%east(j)) THEN
1596!^ vbar(Iend+1,j,kout)=gamma2(ng)*vbar(Iend,j,kout)
1597!^
1598 tl_vbar(iend+1,j,kout)=gamma2(ng)*tl_vbar(iend,j,kout)
1599# ifdef MASKING
1600!^ vbar(Iend+1,j,kout)=vbar(Iend+1,j,kout)* &
1601!^ & GRID(ng)%vmask(Iend+1,j)
1602!^
1603 tl_vbar(iend+1,j,kout)=tl_vbar(iend+1,j,kout)* &
1604 & grid(ng)%vmask(iend+1,j)
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-1).and. &
1618 & lbc_apply(ng)%west (jstr )) THEN
1619!^ vbar(Istr-1,Jstr,kout)=0.5_r8*(vbar(Istr ,Jstr ,kout)+ &
1620!^ & vbar(Istr-1,Jstr+1,kout))
1621!^
1622 tl_vbar(istr-1,jstr,kout)=0.5_r8* &
1623 & (tl_vbar(istr ,jstr ,kout)+ &
1624 & tl_vbar(istr-1,jstr+1,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 )) THEN
1630!^ vbar(Iend+1,Jstr,kout)=0.5_r8*(vbar(Iend ,Jstr ,kout)+ &
1631!^ & vbar(Iend+1,Jstr+1,kout))
1632!^
1633 tl_vbar(iend+1,jstr,kout)=0.5_r8* &
1634 & (tl_vbar(iend ,jstr ,kout)+ &
1635 & tl_vbar(iend+1,jstr+1,kout))
1636 END IF
1637 END IF
1638 IF (domain(ng)%NorthWest_Corner(tile)) THEN
1639 IF (lbc_apply(ng)%north(istr-1).and. &
1640 & lbc_apply(ng)%west (jend+1)) THEN
1641!^ vbar(Istr-1,Jend+1,kout)=0.5_r8*(vbar(Istr-1,Jend ,kout)+ &
1642!^ & vbar(Istr ,Jend+1,kout))
1643!^
1644 tl_vbar(istr-1,jend+1,kout)=0.5_r8* &
1645 & (tl_vbar(istr-1,jend ,kout)+ &
1646 & tl_vbar(istr ,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!^ vbar(Iend+1,Jend+1,kout)=0.5_r8*(vbar(Iend+1,Jend ,kout)+ &
1653!^ & vbar(Iend ,Jend+1,kout))
1654!^
1655 tl_vbar(iend+1,jend+1,kout)=0.5_r8* &
1656 & (tl_vbar(iend+1,jend ,kout)+ &
1657 & tl_vbar(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=jstrv,jend
1673 IF (lbc_apply(ng)%west(j).or. &
1674 & lbc(iwest,isvbar,ng)%nested) THEN
1675!^ cff1=ABS(ABS(GRID(ng)%vmask_wet(Istr-1,j))-1.0_r8)
1676!^ cff2=0.5_r8+DSIGN(0.5_r8,vbar(Istr-1,j,kout))* &
1677!^ & GRID(ng)%vmask_wet(Istr-1,j)
1678!^ cff=0.5_r8*GRID(ng)%vmask_wet(Istr-1,j)*cff1+ &
1679!^ & cff2*(1.0_r8-cff1)
1680!^ vbar(Istr,j,kout)=vbar(Istr,j,kout)*cff
1681 END IF
1682 END DO
1683 END IF
1684 IF (domain(ng)%Eastern_Edge(tile)) THEN
1685 DO j=jstrv,jend
1686 IF (lbc_apply(ng)%east(j).or. &
1687 & lbc(ieast,isvbar,ng)%nested) THEN
1688!^ cff1=ABS(ABS(GRID(ng)%vmask_wet(Iend+1,j))-1.0_r8)
1689!^ cff2=0.5_r8+DSIGN(0.5_r8,vbar(Iend+1,j,kout))* &
1690!^ & GRID(ng)%vmask_wet(Iend+1,j)
1691!^ cff=0.5_r8*GRID(ng)%vmask_wet(Iend+1,j)*cff1+ &
1692!^ & cff2*(1.0_r8-cff1)
1693!^ vbar(Iend+1,j,kout)=vbar(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=istr,iend
1702 IF (lbc_apply(ng)%south(i).or. &
1703 & lbc(isouth,isvbar,ng)%nested) THEN
1704!^ cff1=ABS(ABS(GRID(ng)%vmask_wet(i,Jstr))-1.0_r8)
1705!^ cff2=0.5_r8+DSIGN(0.5_r8,vbar(i,Jstr,kout))* &
1706!^ & GRID(ng)%vmask_wet(i,Jstr)
1707!^ cff=0.5_r8*GRID(ng)%vmask_wet(i,Jstr)*cff1+ &
1708!^ & cff2*(1.0_r8-cff1)
1709!^ vbar(i,Jstr,kout)=vbar(i,Jstr,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,isvbar,ng)%nested) THEN
1717!^ cff1=ABS(ABS(GRID(ng)%vmask_wet(i,Jend+1))-1.0_r8)
1718!^ cff2=0.5_r8+DSIGN(0.5_r8,vbar(i,Jend+1,kout))* &
1719!^ & GRID(ng)%vmask_wet(i,Jend+1)
1720!^ cff=0.5_r8*GRID(ng)%vmask_wet(i,Jend+1)*cff1+ &
1721!^ & cff2*(1.0_r8-cff1)
1722!^ vbar(i,Jend+1,kout)=vbar(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-1).and. &
1731 & lbc_apply(ng)%west (jstr )).or. &
1732 & (lbc(iwest,isvbar,ng)%nested.and. &
1733 & lbc(isouth,isvbar,ng)%nested)) THEN
1734!^ cff1=ABS(ABS(GRID(ng)%vmask_wet(Istr-1,Jstr))-1.0_r8)
1735!^ cff2=0.5_r8+DSIGN(0.5_r8,vbar(Istr-1,Jstr,kout))* &
1736!^ & GRID(ng)%vmask_wet(Istr-1,Jstr)
1737!^ cff=0.5_r8*GRID(ng)%vmask_wet(Istr-1,Jstr)*cff1+ &
1738!^ & cff2*(1.0_r8-cff1)
1739!^ vbar(Istr-1,Jstr,kout)=vbar(Istr-1,Jstr,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 )).or. &
1745 & (lbc(ieast,isvbar,ng)%nested.and. &
1746 & lbc(isouth,isvbar,ng)%nested)) THEN
1747!^ cff1=ABS(ABS(GRID(ng)%vmask_wet(Iend+1,Jstr))-1.0_r8)
1748!^ cff2=0.5_r8+DSIGN(0.5_r8,vbar(Iend+1,Jstr,kout))* &
1749!^ & GRID(ng)%vmask_wet(Iend+1,Jstr)
1750!^ cff=0.5_r8*GRID(ng)%vmask_wet(Iend+1,Jstr)*cff1+ &
1751!^ & cff2*(1.0_r8-cff1)
1752!^ vbar(Iend+1,Jstr,kout)=vbar(Iend+1,Jstr,kout)*cff
1753 END IF
1754 END IF
1755 IF (domain(ng)%NorthWest_Corner(tile)) THEN
1756 IF ((lbc_apply(ng)%north(istr-1).and. &
1757 & lbc_apply(ng)%west (jend+1)).or. &
1758 & (lbc(iwest,isvbar,ng)%nested.and. &
1759 & lbc(inorth,isvbar,ng)%nested)) THEN
1760!^ cff1=ABS(ABS(GRID(ng)%vmask_wet(Istr-1,Jend+1))-1.0_r8)
1761!^ cff2=0.5_r8+DSIGN(0.5_r8,vbar(Istr-1,Jend+1,kout))* &
1762!^ & GRID(ng)%vmask_wet(Istr-1,Jend+1)
1763!^ cff=0.5_r8*GRID(ng)%vmask_wet(Istr-1,Jend+1)*cff1+ &
1764!^ & cff2*(1.0_r8-cff1)
1765!^ vbar(Istr-1,Jend+1,kout)=vbar(Istr-1,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,isvbar,ng)%nested.and. &
1772 & lbc(inorth,isvbar,ng)%nested)) THEN
1773!^ cff1=ABS(ABS(GRID(ng)%vmask_wet(Iend+1,Jend+1))-1.0_r8)
1774!^ cff2=0.5_r8+DSIGN(0.5_r8,vbar(Iend+1,Jend+1,kout))* &
1775!^ & GRID(ng)%vmask_wet(Iend+1,Jend+1)
1776!^ cff=0.5_r8*GRID(ng)%vmask_wet(Iend+1,Jend+1)*cff1+ &
1777!^ & cff2*(1.0_r8-cff1)
1778!^ vbar(Iend+1,Jend+1,kout)=vbar(Iend+1,Jend+1,kout)*cff
1779 END IF
1780 END IF
1781 END IF
1782# endif
1783!
1784 RETURN
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 isvbar
integer isfsur
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

References mod_boundary::boundary, mod_clima::clima, mod_scalars::co, mod_param::domain, mod_scalars::dtfast, mod_scalars::ewperiodic, mod_forces::forces, mod_scalars::g, mod_scalars::gamma2, mod_grid::grid, mod_scalars::ieast, mod_scalars::iic, mod_scalars::inorth, mod_ncparam::isfsur, mod_scalars::isouth, mod_ncparam::isvbar, mod_scalars::iwest, mod_param::lbc, mod_boundary::lbc_apply, mod_scalars::lnudgem2clm, mod_scalars::lobc, mod_scalars::m2obc_in, mod_scalars::m2obc_out, mod_scalars::nsperiodic, mod_scalars::obcfac, mod_scalars::predictor_2d_step, mod_scalars::rho0, and mod_param::tl_lbc.

Referenced by tl_ini_fields_mod::tl_ini_fields_tile(), ini_adjust_mod::tl_ini_perturb_tile(), tl_step2d_mod::tl_step2d_tile(), tl_step2d_mod::tl_step2d_tile(), tl_step2d_mod::tl_step2d_tile(), and tl_v2dbc().

Here is the caller graph for this function: