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