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

Functions/Subroutines

subroutine, public tl_v3dbc (ng, tile, nout)
 
subroutine, public tl_v3dbc_tile (ng, tile, lbi, ubi, lbj, ubj, ubk, imins, imaxs, jmins, jmaxs, nstp, nout, tl_v)
 

Function/Subroutine Documentation

◆ tl_v3dbc()

subroutine, public tl_v3dbc_mod::tl_v3dbc ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) nout )

Definition at line 27 of file tl_v3dbc_im.F.

28!***********************************************************************
29!
30 USE mod_param
31 USE mod_ocean
32 USE mod_stepping
33!
34! Imported variable declarations.
35!
36 integer, intent(in) :: ng, tile, nout
37!
38! Local variable declarations.
39!
40# include "tile.h"
41!
42 CALL tl_v3dbc_tile (ng, tile, &
43 & lbi, ubi, lbj, ubj, n(ng), &
44 & imins, imaxs, jmins, jmaxs, &
45 & nstp(ng), nout, &
46 & ocean(ng) % tl_v)
47
48 RETURN
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer, dimension(:), allocatable nstp

References mod_param::n, mod_stepping::nstp, mod_ocean::ocean, and tl_v3dbc_tile().

Here is the call graph for this function:

◆ tl_v3dbc_tile()

subroutine, public tl_v3dbc_mod::tl_v3dbc_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) ubk,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) nstp,
integer, intent(in) nout,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) tl_v )

Definition at line 52 of file tl_v3dbc_im.F.

57!***********************************************************************
58!
59 USE mod_param
60 USE mod_boundary
61 USE mod_clima
62 USE mod_grid
63 USE mod_ncparam
64 USE mod_scalars
65!
66! Imported variable declarations.
67!
68 integer, intent(in) :: ng, tile
69 integer, intent(in) :: LBi, UBi, LBj, UBj, UBk
70 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
71 integer, intent(in) :: nstp, nout
72!
73# ifdef ASSUMED_SHAPE
74 real(r8), intent(inout) :: tl_v(LBi:,LBj:,:,:)
75# else
76 real(r8), intent(inout) :: tl_v(LBi:UBi,LBj:UBj,UBk,2)
77# endif
78!
79! Local variable declarations.
80!
81 integer :: Jmin, Jmax
82 integer :: i, j, k
83
84 real(r8) :: Ce, Cx, cff
85 real(r8) :: obc_in, obc_out, tau
86
87 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_grad
88
89# include "set_bounds.h"
90!
91!-----------------------------------------------------------------------
92! Lateral boundary conditions at the southern edge.
93!-----------------------------------------------------------------------
94!
95 IF (domain(ng)%Southern_Edge(tile)) THEN
96!
97! Southern edge, implicit upstream radiation condition.
98!
99 IF (tl_lbc(isouth,isvvel,ng)%radiation) THEN
100 IF (iic(ng).ne.0) THEN
101 DO k=1,n(ng)
102 DO i=istr,iend+1
103!^ grad(i,Jstr)=v(i ,Jstr,k,nstp)- &
104!^ & v(i-1,Jstr,k,nstp)
105!^
106 tl_grad(i,jstr)=0.0_r8
107 END DO
108 DO i=istr,iend
109 IF (lbc_apply(ng)%south(i)) THEN
110# if defined CELERITY_READ && defined FORWARD_READ
111 IF (tl_lbc(isouth,isvvel,ng)%nudging) THEN
112 IF (lnudgem3clm(ng)) THEN
113 obc_out=0.5_r8* &
114 & (clima(ng)%M3nudgcof(i,jstr-1,k)+ &
115 & clima(ng)%M3nudgcof(i,jstr ,k))
116 obc_in =obcfac(ng)*obc_out
117 ELSE
118 obc_out=m3obc_out(ng,isouth)
119 obc_in =m3obc_in(ng,isouth)
120 END IF
121 IF (boundary(ng)%v_south_Ce(i,k).lt.0.0_r8) THEN
122 tau=obc_in
123 ELSE
124 tau=obc_out
125 END IF
126 tau=tau*dt(ng)
127 END IF
128# ifdef RADIATION_2D
129 cx=boundary(ng)%v_south_Cx(i,k)
130# else
131 cx=0.0_r8
132# endif
133 ce=boundary(ng)%v_south_Ce(i,k)
134 cff=boundary(ng)%v_south_C2(i,k)
135# endif
136!^ v(i,Jstr,k,nout)=(cff*v(i,Jstr ,k,nstp)+ &
137!^ & Ce *v(i,Jstr+1,k,nout)- &
138!^ & MAX(Cx,0.0_r8)*grad(i ,Jstr)- &
139!^ & MIN(Cx,0.0_r8)*grad(i+1,Jstr))/ &
140!^ & (cff+Ce)
141!^
142 tl_v(i,jstr,k,nout)=(cff*tl_v(i,jstr ,k,nstp)+ &
143 & ce *tl_v(i,jstr+1,k,nout)- &
144 & max(cx,0.0_r8)* &
145 & tl_grad(i ,jstr)- &
146 & min(cx,0.0_r8)* &
147 & tl_grad(i+1,jstr))/ &
148 & (cff+ce)
149
150 IF (tl_lbc(isouth,isvvel,ng)%nudging) THEN
151!^ v(i,Jstr,k,nout)=v(i,Jstr,k,nout)+ &
152!^ & tau*(BOUNDARY(ng)%v_south(i,k)- &
153!^ & v(i,Jstr,k,nstp))
154!^
155 tl_v(i,jstr,k,nout)=tl_v(i,jstr,k,nout)- &
156 & tau*tl_v(i,jstr,k,nstp)
157 END IF
158# ifdef MASKING
159!^ v(i,Jstr,k,nout)=v(i,Jstr,k,nout)* &
160!^ & GRID(ng)%vmask(i,Jstr)
161!^
162 tl_v(i,jstr,k,nout)=tl_v(i,jstr,k,nout)* &
163 & grid(ng)%vmask(i,jstr)
164# endif
165 END IF
166 END DO
167 END DO
168 END IF
169!
170! Southern edge, clamped boundary condition.
171!
172 ELSE IF (tl_lbc(isouth,isvvel,ng)%clamped) THEN
173 DO k=1,n(ng)
174 DO i=istr,iend
175 IF (lbc_apply(ng)%south(i)) THEN
176!^ v(i,Jstr,k,nout)=BOUNDARY(ng)%v_south(i,k)
177!^
178# ifdef ADJUST_BOUNDARY
179 IF (lobc(isouth,isvvel,ng)) THEN
180 tl_v(i,jstr,k,nout)=boundary(ng)%tl_v_south(i,k)
181 ELSE
182 tl_v(i,jstr,k,nout)=0.0_r8
183 END IF
184# else
185 tl_v(i,jstr,k,nout)=0.0_r8
186# endif
187# ifdef MASKING
188!^ v(i,Jstr,k,nout)=v(i,Jstr,k,nout)* &
189!^ & GRID(ng)%vmask(i,Jstr)
190!^
191 tl_v(i,jstr,k,nout)=tl_v(i,jstr,k,nout)* &
192 & grid(ng)%vmask(i,jstr)
193# endif
194 END IF
195 END DO
196 END DO
197!
198! Southern edge, gradient boundary condition.
199!
200 ELSE IF (tl_lbc(isouth,isvvel,ng)%gradient) THEN
201 DO k=1,n(ng)
202 DO i=istr,iend
203 IF (lbc_apply(ng)%south(i)) THEN
204!^ v(i,Jstr,k,nout)=v(i,Jstr+1,k,nout)
205!^
206 tl_v(i,jstr,k,nout)=tl_v(i,jstr+1,k,nout)
207# ifdef MASKING
208!^ v(i,Jstr,k,nout)=v(i,Jstr,k,nout)* &
209!^ & GRID(ng)%vmask(i,Jstr)
210!^
211 tl_v(i,jstr,k,nout)=tl_v(i,jstr,k,nout)* &
212 & grid(ng)%vmask(i,jstr)
213# endif
214 END IF
215 END DO
216 END DO
217!
218! Southern edge, closed boundary condition.
219!
220 ELSE IF (tl_lbc(isouth,isvvel,ng)%closed) THEN
221 DO k=1,n(ng)
222 DO i=istr,iend
223 IF (lbc_apply(ng)%south(i)) THEN
224!^ v(i,Jstr,k,nout)=0.0_r8
225!^
226 tl_v(i,jstr,k,nout)=0.0_r8
227 END IF
228 END DO
229 END DO
230 END IF
231 END IF
232!
233!-----------------------------------------------------------------------
234! Lateral boundary conditions at the northern edge.
235!-----------------------------------------------------------------------
236!
237 IF (domain(ng)%Northern_Edge(tile)) THEN
238!
239! Northern edge, implicit upstream radiation condition.
240!
241 IF (tl_lbc(inorth,isvvel,ng)%radiation) THEN
242 IF (iic(ng).ne.0) THEN
243 DO k=1,n(ng)
244 DO i=istr,iend+1
245!^ grad(i,Jend+1)=v(i ,Jend+1,k,nstp)- &
246!^ & v(i-1,Jend+1,k,nstp)
247!^
248 tl_grad(i,jend+1)=0.0_r8
249 END DO
250 DO i=istr,iend
251 IF (lbc_apply(ng)%north(i)) THEN
252# if defined CELERITY_READ && defined FORWARD_READ
253 IF (tl_lbc(inorth,isvvel,ng)%nudging) THEN
254 IF (lnudgem3clm(ng)) THEN
255 obc_out=0.5_r8* &
256 & (clima(ng)%M3nudgcof(i,jend ,k)+ &
257 & clima(ng)%M3nudgcof(i,jend+1,k))
258 obc_in =obcfac(ng)*obc_out
259 ELSE
260 obc_out=m3obc_out(ng,inorth)
261 obc_in =m3obc_in(ng,inorth)
262 END IF
263 IF (boundary(ng)%v_south_Ce(i,k).lt.0.0_r8) THEN
264 tau=obc_in
265 ELSE
266 tau=obc_out
267 END IF
268 tau=tau*dt(ng)
269 END IF
270# ifdef RADIATION_2D
271 cx=boundary(ng)%v_south_Cx(i,k)
272# else
273 cx=0.0_r8
274# endif
275 ce=boundary(ng)%v_south_Ce(i,k)
276 cff=boundary(ng)%v_south_C2(i,k)
277# endif
278!^ v(i,Jend+1,k,nout)=(cff*v(i,Jend+1,k,nstp)+ &
279!^ & Ce *v(i,Jend ,k,nout)- &
280!^ & MAX(Cx,0.0_r8)*grad(i ,Jend+1)- &
281!^ & MIN(Cx,0.0_r8)*grad(i+1,Jend+1))/ &
282!^ & (cff+Ce)
283!^
284 tl_v(i,jend+1,k,nout)=(cff*tl_v(i,jend+1,k,nstp)+ &
285 & ce *tl_v(i,jend ,k,nout)- &
286 & max(cx,0.0_r8)* &
287 & tl_grad(i ,jend+1)- &
288 & min(cx,0.0_r8)* &
289 & tl_grad(i+1,jend+1))/ &
290 & (cff+ce)
291
292 IF (tl_lbc(inorth,isvvel,ng)%nudging) THEN
293!^ v(i,Jend+1,k,nout)=v(i,Jend+1,k,nout)+ &
294!^ & tau*(BOUNDARY(ng)%v_north(i,k)- &
295!^ & v(i,Jend+1,k,nstp))
296!^
297 tl_v(i,jend+1,k,nout)=tl_v(i,jend+1,k,nout)- &
298 & tau*tl_v(i,jend+1,k,nstp)
299 END IF
300# ifdef MASKING
301!^ v(i,Jend+1,k,nout)=v(i,Jend+1,k,nout)* &
302!^ & GRID(ng)%vmask(i,Jend+1)
303!^
304 tl_v(i,jend+1,k,nout)=tl_v(i,jend+1,k,nout)* &
305 & grid(ng)%vmask(i,jend+1)
306# endif
307 END IF
308 END DO
309 END DO
310 END IF
311!
312! Northern edge, clamped boundary condition.
313!
314 ELSE IF (tl_lbc(inorth,isvvel,ng)%clamped) THEN
315 DO k=1,n(ng)
316 DO i=istr,iend
317 IF (lbc_apply(ng)%north(i)) THEN
318!^ v(i,Jend+1,k,nout)=BOUNDARY(ng)%v_north(i,k)
319!^
320# ifdef ADJUST_BOUNDARY
321 IF (lobc(inorth,isvvel,ng)) THEN
322 tl_v(i,jend+1,k,nout)=boundary(ng)%tl_v_north(i,k)
323 ELSE
324 tl_v(i,jend+1,k,nout)=0.0_r8
325 END IF
326# else
327 tl_v(i,jend+1,k,nout)=0.0_r8
328# endif
329# ifdef MASKING
330!^ v(i,Jend+1,k,nout)=v(i,Jend+1,k,nout)* &
331!^ & GRID(ng)%vmask(i,Jend+1)
332!^
333 tl_v(i,jend+1,k,nout)=tl_v(i,jend+1,k,nout)* &
334 & grid(ng)%vmask(i,jend+1)
335# endif
336 END IF
337 END DO
338 END DO
339!
340! Northern edge, gradient boundary condition.
341!
342 ELSE IF (tl_lbc(inorth,isvvel,ng)%gradient) THEN
343 DO k=1,n(ng)
344 DO i=istr,iend
345 IF (lbc_apply(ng)%north(i)) THEN
346!^ v(i,Jend+1,k,nout)=v(i,Jend,k,nout)
347!^
348 tl_v(i,jend+1,k,nout)=tl_v(i,jend,k,nout)
349# ifdef MASKING
350!^ v(i,Jend+1,k,nout)=v(i,Jend+1,k,nout)* &
351!^ & GRID(ng)%vmask(i,Jend+1)
352!^
353 tl_v(i,jend+1,k,nout)=tl_v(i,jend+1,k,nout)* &
354 & grid(ng)%vmask(i,jend+1)
355# endif
356 END IF
357 END DO
358 END DO
359!
360! Northern edge, closed boundary condition.
361!
362 ELSE IF (tl_lbc(inorth,isvvel,ng)%closed) THEN
363 DO k=1,n(ng)
364 DO i=istr,iend
365 IF (lbc_apply(ng)%north(i)) THEN
366!^ v(i,Jend+1,k,nout)=0.0_r8
367!^
368 tl_v(i,jend+1,k,nout)=0.0_r8
369 END IF
370 END DO
371 END DO
372 END IF
373 END IF
374!
375!-----------------------------------------------------------------------
376! Lateral boundary conditions at the western edge.
377!-----------------------------------------------------------------------
378!
379 IF (domain(ng)%Western_Edge(tile)) THEN
380!
381! Western edge, implicit upstream radiation condition.
382!
383 IF (tl_lbc(iwest,isvvel,ng)%radiation) THEN
384 IF (iic(ng).ne.0) THEN
385 DO k=1,n(ng)
386 DO j=jstrv-1,jend
387!^ grad(Istr-1,j)=v(Istr-1,j+1,k,nstp)- &
388!^ & v(Istr-1,j ,k,nstp)
389!^
390 tl_grad(istr-1,j)=0.0_r8
391 END DO
392 DO j=jstrv,jend
393 IF (lbc_apply(ng)%west(j)) THEN
394# if defined CELERITY_READ && defined FORWARD_READ
395 IF (tl_lbc(iwest,isvvel,ng)%nudging) THEN
396 IF (lnudgem3clm(ng)) THEN
397 obc_out=0.5_r8* &
398 & (clima(ng)%M3nudgcof(istr-1,j-1,k)+ &
399 & clima(ng)%M3nudgcof(istr-1,j ,k))
400 obc_in =obcfac(ng)*obc_out
401 ELSE
402 obc_out=m3obc_out(ng,iwest)
403 obc_in =m3obc_in(ng,iwest)
404 END IF
405 IF (boundary(ng)%v_west_Cx(j,k).lt.0.0_r8) THEN
406 tau=obc_in
407 ELSE
408 tau=obc_out
409 END IF
410 tau=tau*dt(ng)
411 END IF
412 cx=boundary(ng)%v_west_Cx(j,k)
413# ifdef RADIATION_2D
414 ce=boundary(ng)%v_west_Ce(j,k)
415# else
416 ce=0.0_r8
417# endif
418 cff=boundary(ng)%v_west_C2(j,k)
419# endif
420!^ v(Istr-1,j,k,nout)=(cff*v(Istr-1,j,k,nstp)+ &
421!^ & Cx *v(Istr ,j,k,nout)- &
422!^ & MAX(Ce,0.0_r8)*grad(Istr-1,j-1)- &
423!^ & MIN(Ce,0.0_r8)*grad(Istr-1,j ))/ &
424!^ & (cff+Cx)
425!^
426 tl_v(istr-1,j,k,nout)=(cff*tl_v(istr-1,j,k,nstp)+ &
427 & cx *tl_v(istr ,j,k,nout)- &
428 & max(ce,0.0_r8)* &
429 & tl_grad(istr-1,j-1)- &
430 & min(ce,0.0_r8)* &
431 & tl_grad(istr-1,j ))/ &
432 & (cff+cx)
433
434 IF (tl_lbc(iwest,isvvel,ng)%nudging) THEN
435!^ v(Istr-1,j,k,nout)=v(Istr-1,j,k,nout)+ &
436!^ & tau*(BOUNDARY(ng)%v_west(j,k)- &
437!^ & v(Istr-1,j,k,nstp))
438!^
439 tl_v(istr-1,j,k,nout)=tl_v(istr-1,j,k,nout)- &
440 & tau*tl_v(istr-1,j,k,nstp)
441 END IF
442# ifdef MASKING
443!^ v(Istr-1,j,k,nout)=v(Istr-1,j,k,nout)* &
444!^ & GRID(ng)%vmask(Istr-1,j)
445!^
446 tl_v(istr-1,j,k,nout)=tl_v(istr-1,j,k,nout)* &
447 & grid(ng)%vmask(istr-1,j)
448# endif
449 END IF
450 END DO
451 END DO
452 END IF
453!
454! Western edge, clamped boundary condition.
455!
456 ELSE IF (tl_lbc(iwest,isvvel,ng)%clamped) THEN
457 DO k=1,n(ng)
458 DO j=jstrv,jend
459 IF (lbc_apply(ng)%west(j)) THEN
460!^ v(Istr-1,j,k,nout)=BOUNDARY(ng)%v_west(j,k)
461!^
462# ifdef ADJUST_BOUNDARY
463 IF (lobc(iwest,isvvel,ng)) THEN
464 tl_v(istr-1,j,k,nout)=boundary(ng)%tl_v_west(j,k)
465 ELSE
466 tl_v(istr-1,j,k,nout)=0.0_r8
467 END IF
468# else
469 tl_v(istr-1,j,k,nout)=0.0_r8
470# endif
471# ifdef MASKING
472!^ v(Istr-1,j,k,nout)=v(Istr-1,j,k,nout)* &
473!^ & GRID(ng)%vmask(Istr-1,j)
474!^
475 tl_v(istr-1,j,k,nout)=tl_v(istr-1,j,k,nout)* &
476 & grid(ng)%vmask(istr-1,j)
477# endif
478 END IF
479 END DO
480 END DO
481!
482! Western edge, gradient boundary condition.
483!
484 ELSE IF (tl_lbc(iwest,isvvel,ng)%gradient) THEN
485 DO k=1,n(ng)
486 DO j=jstrv,jend
487 IF (lbc_apply(ng)%west(j)) THEN
488!^ v(Istr-1,j,k,nout)=v(Istr,j,k,nout)
489!^
490 tl_v(istr-1,j,k,nout)=tl_v(istr,j,k,nout)
491# ifdef MASKING
492!^ v(Istr-1,j,k,nout)=v(Istr-1,j,k,nout)* &
493!^ & GRID(ng)%vmask(Istr-1,j)
494!^
495 tl_v(istr-1,j,k,nout)=tl_v(istr-1,j,k,nout)* &
496 & grid(ng)%vmask(istr-1,j)
497# endif
498 END IF
499 END DO
500 END DO
501!
502! Western edge, closed boundary condition: free slip (gamma2=1) or
503! no slip (gamma2=-1).
504!
505 ELSE IF (tl_lbc(iwest,isvvel,ng)%closed) THEN
506 IF (nsperiodic(ng)) THEN
507 jmin=jstrv
508 jmax=jend
509 ELSE
510 jmin=jstr
511 jmax=jendr
512 END IF
513 DO k=1,n(ng)
514 DO j=jmin,jmax
515 IF (lbc_apply(ng)%west(j)) THEN
516!^ v(Istr-1,j,k,nout)=gamma2(ng)*v(Istr,j,k,nout)
517!^
518 tl_v(istr-1,j,k,nout)=gamma2(ng)*tl_v(istr,j,k,nout)
519# ifdef MASKING
520!^ v(Istr-1,j,k,nout)=v(Istr-1,j,k,nout)* &
521!^ & GRID(ng)%vmask(Istr-1,j)
522!^
523 tl_v(istr-1,j,k,nout)=tl_v(istr-1,j,k,nout)* &
524 & grid(ng)%vmask(istr-1,j)
525# endif
526 END IF
527 END DO
528 END DO
529 END IF
530 END IF
531!
532!-----------------------------------------------------------------------
533! Lateral boundary conditions at the eastern edge.
534!-----------------------------------------------------------------------
535!
536 IF (domain(ng)%Eastern_Edge(tile)) THEN
537!
538! Eastern edge, implicit upstream radiation condition.
539!
540 IF (tl_lbc(ieast,isvvel,ng)%radiation) THEN
541 IF (iic(ng).ne.0) THEN
542 DO k=1,n(ng)
543 DO j=jstrv-1,jend
544!^ grad(Iend+1,j)=v(Iend+1,j+1,k,nstp)- &
545!^ & v(Iend+1,j ,k,nstp)
546!^
547 tl_grad(iend+1,j)=0.0_r8
548 END DO
549 DO j=jstrv,jend
550 IF (lbc_apply(ng)%east(j)) THEN
551# if defined CELERITY_READ && defined FORWARD_READ
552 IF (tl_lbc(ieast,isvvel,ng)%nudging) THEN
553 IF (lnudgem3clm(ng)) THEN
554 obc_out=0.5_r8* &
555 & (clima(ng)%M3nudgcof(iend+1,j-1,k)+ &
556 & clima(ng)%M3nudgcof(iend+1,j ,k))
557 obc_in =obcfac(ng)*obc_out
558 ELSE
559 obc_out=m3obc_out(ng,ieast)
560 obc_in =m3obc_in(ng,ieast)
561 END IF
562 IF (boundary(ng)%v_east_Cx(j,k).lt.0.0_r8) THEN
563 tau=obc_in
564 ELSE
565 tau=obc_out
566 END IF
567 tau=tau*dt(ng)
568 END IF
569 cx=boundary(ng)%v_east_Cx(j,k)
570# ifdef RADIATION_2D
571 ce=boundary(ng)%v_east_Ce(j,k)
572# else
573 ce=0.0_r8
574# endif
575 cff=boundary(ng)%v_east_C2(j,k)
576# endif
577!^ v(Iend+1,j,k,nout)=(cff*v(Iend+1,j,k,nstp)+ &
578!^ & Cx *v(Iend ,j,k,nout)- &
579!^ & MAX(Ce,0.0_r8)*grad(Iend+1,j-1)- &
580!^ & MIN(Ce,0.0_r8)*grad(Iend+1,j ))/ &
581!^ & (cff+Cx)
582!^
583 tl_v(iend+1,j,k,nout)=(cff*tl_v(iend+1,j,k,nstp)+ &
584 & cx *tl_v(iend ,j,k,nout)- &
585 & max(ce,0.0_r8)* &
586 & tl_grad(iend+1,j-1)- &
587 & min(ce,0.0_r8)* &
588 & tl_grad(iend+1,j ))/ &
589 & (cff+cx)
590
591 IF (tl_lbc(ieast,isvvel,ng)%nudging) THEN
592!^ v(Iend+1,j,k,nout)=v(Iend+1,j,k,nout)+ &
593!^ & tau*(BOUNDARY(ng)%v_east(j,k)- &
594!^ & v(Iend+1,j,k,nstp))
595!^
596 tl_v(iend+1,j,k,nout)=tl_v(iend+1,j,k,nout)- &
597 & tau*tl_v(iend+1,j,k,nstp)
598 END IF
599# ifdef MASKING
600!^ v(Iend+1,j,k,nout)=v(Iend+1,j,k,nout)* &
601!^ & GRID(ng)%vmask(Iend+1,j)
602!^
603 tl_v(iend+1,j,k,nout)=tl_v(iend+1,j,k,nout)* &
604 & grid(ng)%vmask(iend+1,j)
605# endif
606 END IF
607 END DO
608 END DO
609 END IF
610!
611! Eastern edge, clamped boundary condition.
612!
613 ELSE IF (tl_lbc(ieast,isvvel,ng)%clamped) THEN
614 DO k=1,n(ng)
615 DO j=jstrv,jend
616 IF (lbc_apply(ng)%east(j)) THEN
617!^ v(Iend+1,j,k,nout)=BOUNDARY(ng)%v_east(j,k)
618!^
619# ifdef ADJUST_BOUNDARY
620 IF (lobc(ieast,isvvel,ng)) THEN
621 tl_v(iend+1,j,k,nout)=boundary(ng)%tl_v_east(j,k)
622 ELSE
623 tl_v(iend+1,j,k,nout)=0.0_r8
624 END IF
625# else
626 tl_v(iend+1,j,k,nout)=0.0_r8
627# endif
628# ifdef MASKING
629!^ v(Iend+1,j,k,nout)=v(Iend+1,j,k,nout)* &
630!^ & GRID(ng)%vmask(Iend+1,j)
631!^
632 tl_v(iend+1,j,k,nout)=tl_v(iend+1,j,k,nout)* &
633 & grid(ng)%vmask(iend+1,j)
634# endif
635 END IF
636 END DO
637 END DO
638!
639! Eastern edge, gradient boundary condition.
640!
641 ELSE IF (tl_lbc(ieast,isvvel,ng)%gradient) THEN
642 DO k=1,n(ng)
643 DO j=jstrv,jend
644 IF (lbc_apply(ng)%east(j)) THEN
645!^ v(Iend+1,j,k,nout)=v(Iend,j,k,nout)
646!^
647 tl_v(iend+1,j,k,nout)=tl_v(iend,j,k,nout)
648# ifdef MASKING
649!^ v(Iend+1,j,k,nout)=v(Iend+1,j,k,nout)* &
650!^ & GRID(ng)%vmask(Iend+1,j)
651!^
652 tl_v(iend+1,j,k,nout)=tl_v(iend+1,j,k,nout)* &
653 & grid(ng)%vmask(iend+1,j)
654# endif
655 END IF
656 END DO
657 END DO
658!
659! Eastern edge, closed boundary condition: free slip (gamma2=1) or
660! no slip (gamma2=-1).
661!
662 ELSE IF (tl_lbc(ieast,isvvel,ng)%closed) THEN
663 IF (nsperiodic(ng)) THEN
664 jmin=jstrv
665 jmax=jend
666 ELSE
667 jmin=jstr
668 jmax=jendr
669 END IF
670 DO k=1,n(ng)
671 DO j=jmin,jmax
672 IF (lbc_apply(ng)%east(j)) THEN
673!^ v(Iend+1,j,k,nout)=gamma2(ng)*v(Iend,j,k,nout)
674!^
675 tl_v(iend+1,j,k,nout)=gamma2(ng)*tl_v(iend,j,k,nout)
676# ifdef MASKING
677!^ v(Iend+1,j,k,nout)=v(Iend+1,j,k,nout)* &
678!^ & GRID(ng)%vmask(Iend+1,j)
679!^
680 tl_v(iend+1,j,k,nout)=tl_v(iend+1,j,k,nout)* &
681 & grid(ng)%vmask(iend+1,j)
682# endif
683 END IF
684 END DO
685 END DO
686 END IF
687 END IF
688!
689!-----------------------------------------------------------------------
690! Boundary corners.
691!-----------------------------------------------------------------------
692!
693 IF (.not.(ewperiodic(ng).or.nsperiodic(ng))) THEN
694 IF (domain(ng)%SouthWest_Corner(tile)) THEN
695 IF (lbc_apply(ng)%south(istr-1).and. &
696 & lbc_apply(ng)%west (jstr )) THEN
697 DO k=1,n(ng)
698!^ v(Istr-1,Jstr,k,nout)=0.5_r8*(v(Istr ,Jstr ,k,nout)+ &
699!^ & v(Istr-1,Jstr+1,k,nout))
700!^
701 tl_v(istr-1,jstr,k,nout)=0.5_r8* &
702 & (tl_v(istr ,jstr ,k,nout)+ &
703 & tl_v(istr-1,jstr+1,k,nout))
704 END DO
705 END IF
706 END IF
707 IF (domain(ng)%SouthEast_Corner(tile)) THEN
708 IF (lbc_apply(ng)%south(iend+1).and. &
709 & lbc_apply(ng)%east (jstr )) THEN
710 DO k=1,n(ng)
711!^ v(Iend+1,Jstr,k,nout)=0.5_r8*(v(Iend ,Jstr ,k,nout)+ &
712!^ & v(Iend+1,Jstr+1,k,nout))
713!^
714 tl_v(iend+1,jstr,k,nout)=0.5_r8* &
715 & (tl_v(iend ,jstr ,k,nout)+ &
716 & tl_v(iend+1,jstr+1,k,nout))
717 END DO
718 END IF
719 END IF
720 IF (domain(ng)%NorthWest_Corner(tile)) THEN
721 IF (lbc_apply(ng)%north(istr-1).and. &
722 & lbc_apply(ng)%west (jend+1)) THEN
723 DO k=1,n(ng)
724!^ v(Istr-1,Jend+1,k,nout)=0.5_r8*(v(Istr-1,Jend ,k,nout)+ &
725!^ & v(Istr ,Jend+1,k,nout))
726!^
727 tl_v(istr-1,jend+1,k,nout)=0.5_r8* &
728 & (tl_v(istr-1,jend ,k,nout)+ &
729 & tl_v(istr ,jend+1,k,nout))
730 END DO
731 END IF
732 END IF
733 IF (domain(ng)%NorthEast_Corner(tile)) THEN
734 IF (lbc_apply(ng)%north(iend+1).and. &
735 & lbc_apply(ng)%east (jend+1)) THEN
736 DO k=1,n(ng)
737!^ v(Iend+1,Jend+1,k,nout)=0.5_r8*(v(Iend+1,Jend ,k,nout)+ &
738!^ & v(Iend ,Jend+1,k,nout))
739!^
740 tl_v(iend+1,jend+1,k,nout)=0.5_r8* &
741 & (tl_v(iend+1,jend ,k,nout)+ &
742 & tl_v(iend ,jend+1,k,nout))
743 END DO
744 END IF
745 END IF
746 END IF
747
748 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_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
integer isvvel
type(t_lbc), dimension(:,:,:), allocatable tl_lbc
Definition mod_param.F:379
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
real(dp), dimension(:,:), allocatable m3obc_out
integer, dimension(:), allocatable iic
real(dp), dimension(:), allocatable dt
logical, dimension(:,:,:), allocatable lobc
logical, dimension(:), allocatable ewperiodic
integer, parameter iwest
logical, dimension(:), allocatable nsperiodic
real(dp), dimension(:), allocatable obcfac
real(r8), dimension(:), allocatable gamma2
logical, dimension(:), allocatable lnudgem3clm
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth
real(dp), dimension(:,:), allocatable m3obc_in

References mod_boundary::boundary, mod_clima::clima, mod_param::domain, mod_scalars::dt, mod_scalars::ewperiodic, mod_scalars::gamma2, mod_grid::grid, mod_scalars::ieast, mod_scalars::iic, mod_scalars::inorth, mod_scalars::isouth, mod_ncparam::isvvel, mod_scalars::iwest, mod_boundary::lbc_apply, mod_scalars::lnudgem3clm, mod_scalars::lobc, mod_scalars::m3obc_in, mod_scalars::m3obc_out, mod_param::n, mod_scalars::nsperiodic, mod_scalars::obcfac, and mod_param::tl_lbc.

Referenced by tl_ini_fields_mod::tl_ini_fields_tile(), ini_adjust_mod::tl_ini_perturb_tile(), tl_step3d_uv_mod::tl_step3d_uv_tile(), and tl_v3dbc().

Here is the caller graph for this function: