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