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

Functions/Subroutines

subroutine u3dbc (ng, tile, nout)
 
subroutine, public u3dbc_tile (ng, tile, lbi, ubi, lbj, ubj, ubk, imins, imaxs, jmins, jmaxs, nstp, nout, u)
 

Function/Subroutine Documentation

◆ u3dbc()

subroutine u3dbc_mod::u3dbc ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) nout )
private

Definition at line 25 of file u3dbc_im.F.

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
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 u3dbc_tile().

Here is the call graph for this function:

◆ u3dbc_tile()

subroutine public u3dbc_mod::u3dbc_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) u )

Definition at line 50 of file u3dbc_im.F.

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, phi, 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,nout)-u(istr+2,j,k,nout)
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# ifdef IMPLICIT_NUDGING
128 IF (tau.gt.0.0_r8) tau=1.0_r8/tau
129# else
130 tau=tau*dt(ng)
131# endif
132 END IF
133
134 IF ((dudt*dudx).lt.0.0_r8) dudt=0.0_r8
135 IF ((dudt*(grad(istr+1,j )+ &
136 & grad(istr+1,j+1))).gt.0.0_r8) THEN
137 dude=grad(istr+1,j )
138 ELSE
139 dude=grad(istr+1,j+1)
140 END IF
141 cff=max(dudx*dudx+dude*dude,eps)
142 cx=dudt*dudx
143# ifdef RADIATION_2D
144 ce=min(cff,max(dudt*dude,-cff))
145# else
146 ce=0.0_r8
147# endif
148# if defined CELERITY_WRITE && defined FORWARD_WRITE
149 boundary(ng)%u_west_Cx(j,k)=cx
150 boundary(ng)%u_west_Ce(j,k)=ce
151 boundary(ng)%u_west_C2(j,k)=cff
152# endif
153 u(istr,j,k,nout)=(cff*u(istr ,j,k,nstp)+ &
154 & cx *u(istr+1,j,k,nout)- &
155 & max(ce,0.0_r8)*grad(istr,j )- &
156 & min(ce,0.0_r8)*grad(istr,j+1))/ &
157 & (cff+cx)
158
159 IF (lbc(iwest,isuvel,ng)%nudging) THEN
160# ifdef IMPLICIT_NUDGING
161 phi=dt(ng)/(tau+dt(ng))
162 u(istr,j,k,nout)=(1.0_r8-phi)*u(istr,j,k,nout)+ &
163 & phi*boundary(ng)%u_west(j,k)
164# else
165 u(istr,j,k,nout)=u(istr,j,k,nout)+ &
166 & tau*(boundary(ng)%u_west(j,k)- &
167 & u(istr,j,k,nstp))
168# endif
169 END IF
170# ifdef MASKING
171 u(istr,j,k,nout)=u(istr,j,k,nout)* &
172 & grid(ng)%umask(istr,j)
173# endif
174# ifdef WET_DRY
175 u(istr,j,k,nout)=u(istr,j,k,nout)* &
176 & grid(ng)%umask_wet(istr,j)
177# endif
178 END IF
179 END DO
180 END DO
181!
182! Western edge, clamped boundary condition.
183!
184 ELSE IF (lbc(iwest,isuvel,ng)%clamped) THEN
185 DO k=1,n(ng)
186 DO j=jstr,jend
187 IF (lbc_apply(ng)%west(j)) THEN
188 u(istr,j,k,nout)=boundary(ng)%u_west(j,k)
189# ifdef MASKING
190 u(istr,j,k,nout)=u(istr,j,k,nout)* &
191 & grid(ng)%umask(istr,j)
192# endif
193# ifdef WET_DRY
194 u(istr,j,k,nout)=u(istr,j,k,nout)* &
195 & grid(ng)%umask_wet(istr,j)
196# endif
197 END IF
198 END DO
199 END DO
200!
201! Western edge, gradient boundary condition.
202!
203 ELSE IF (lbc(iwest,isuvel,ng)%gradient) THEN
204 DO k=1,n(ng)
205 DO j=jstr,jend
206 IF (lbc_apply(ng)%west(j)) THEN
207 u(istr,j,k,nout)=u(istr+1,j,k,nout)
208# ifdef MASKING
209 u(istr,j,k,nout)=u(istr,j,k,nout)* &
210 & grid(ng)%umask(istr,j)
211# endif
212# ifdef WET_DRY
213 u(istr,j,k,nout)=u(istr,j,k,nout)* &
214 & grid(ng)%umask_wet(istr,j)
215# endif
216 END IF
217 END DO
218 END DO
219!
220! Western edge, closed boundary condition.
221!
222 ELSE IF (lbc(iwest,isuvel,ng)%closed) THEN
223 DO k=1,n(ng)
224 DO j=jstr,jend
225 IF (lbc_apply(ng)%west(j)) THEN
226 u(istr,j,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 eastern edge.
235!-----------------------------------------------------------------------
236!
237 IF (domain(ng)%Eastern_Edge(tile)) THEN
238!
239! Eastern edge, implicit upstream radiation condition.
240!
241 IF (lbc(ieast,isuvel,ng)%radiation) THEN
242 DO k=1,n(ng)
243 DO j=jstr,jend+1
244 grad(iend ,j)=u(iend ,j ,k,nstp)- &
245 & u(iend ,j-1,k,nstp)
246 grad(iend+1,j)=u(iend+1,j ,k,nstp)- &
247 & u(iend+1,j-1,k,nstp)
248 END DO
249 DO j=jstr,jend
250 IF (lbc_apply(ng)%east(j)) THEN
251 dudt=u(iend,j,k,nstp)-u(iend ,j,k,nout)
252 dudx=u(iend,j,k,nout)-u(iend-1,j,k,nout)
253
254 IF (lbc(ieast,isuvel,ng)%nudging) THEN
255 IF (lnudgem3clm(ng)) THEN
256 obc_out=0.5_r8* &
257 & (clima(ng)%M3nudgcof(iend ,j,k)+ &
258 & clima(ng)%M3nudgcof(iend+1,j,k))
259 obc_in =obcfac(ng)*obc_out
260 ELSE
261 obc_out=m3obc_out(ng,ieast)
262 obc_in =m3obc_in(ng,ieast)
263 END IF
264 IF ((dudt*dudx).lt.0.0_r8) THEN
265 tau=obc_in
266 ELSE
267 tau=obc_out
268 END IF
269# ifdef IMPLICIT_NUDGING
270 IF (tau.gt.0.0_r8) tau=1.0_r8/tau
271# else
272 tau=tau*dt(ng)
273# endif
274 END IF
275
276 IF ((dudt*dudx).lt.0.0_r8) dudt=0.0_r8
277 IF ((dudt*(grad(iend,j )+ &
278 & grad(iend,j+1))).gt.0.0_r8) THEN
279 dude=grad(iend,j )
280 ELSE
281 dude=grad(iend,j+1)
282 END IF
283 cff=max(dudx*dudx+dude*dude,eps)
284 cx=dudt*dudx
285# ifdef RADIATION_2D
286 ce=min(cff,max(dudt*dude,-cff))
287# else
288 ce=0.0_r8
289# endif
290# if defined CELERITY_WRITE && defined FORWARD_WRITE
291 boundary(ng)%u_east_Cx(j,k)=cx
292 boundary(ng)%u_east_Ce(j,k)=ce
293 boundary(ng)%u_east_C2(j,k)=cff
294# endif
295 u(iend+1,j,k,nout)=(cff*u(iend+1,j,k,nstp)+ &
296 & cx *u(iend ,j,k,nout)- &
297 & max(ce,0.0_r8)*grad(iend+1,j )- &
298 & min(ce,0.0_r8)*grad(iend+1,j+1))/ &
299 & (cff+cx)
300
301 IF (lbc(ieast,isuvel,ng)%nudging) THEN
302# ifdef IMPLICIT_NUDGING
303 phi=dt(ng)/(tau+dt(ng))
304 u(iend+1,j,k,nout)=(1.0_r8-phi)*u(iend+1,j,k,nout)+ &
305 & phi*boundary(ng)%u_east(j,k)
306# else
307 u(iend+1,j,k,nout)=u(iend+1,j,k,nout)+ &
308 & tau*(boundary(ng)%u_east(j,k)- &
309 & u(iend+1,j,k,nstp))
310# endif
311 END IF
312# ifdef MASKING
313 u(iend+1,j,k,nout)=u(iend+1,j,k,nout)* &
314 & grid(ng)%umask(iend+1,j)
315# endif
316# ifdef WET_DRY
317 u(iend+1,j,k,nout)=u(iend+1,j,k,nout)* &
318 & grid(ng)%umask_wet(iend+1,j)
319# endif
320 END IF
321 END DO
322 END DO
323!
324! Eastern edge, clamped boundary condition.
325!
326 ELSE IF (lbc(ieast,isuvel,ng)%clamped) THEN
327 DO k=1,n(ng)
328 DO j=jstr,jend
329 IF (lbc_apply(ng)%east(j)) THEN
330 u(iend+1,j,k,nout)=boundary(ng)%u_east(j,k)
331# ifdef MASKING
332 u(iend+1,j,k,nout)=u(iend+1,j,k,nout)* &
333 & grid(ng)%umask(iend+1,j)
334# endif
335# ifdef WET_DRY
336 u(iend+1,j,k,nout)=u(iend+1,j,k,nout)* &
337 & grid(ng)%umask_wet(iend+1,j)
338# endif
339 END IF
340 END DO
341 END DO
342!
343! Eastern edge, gradient boundary condition.
344!
345 ELSE IF (lbc(ieast,isuvel,ng)%gradient) THEN
346 DO k=1,n(ng)
347 DO j=jstr,jend
348 IF (lbc_apply(ng)%east(j)) THEN
349 u(iend+1,j,k,nout)=u(iend,j,k,nout)
350# ifdef MASKING
351 u(iend+1,j,k,nout)=u(iend+1,j,k,nout)* &
352 & grid(ng)%umask(iend+1,j)
353# endif
354# ifdef WET_DRY
355 u(iend+1,j,k,nout)=u(iend+1,j,k,nout)* &
356 & grid(ng)%umask_wet(iend+1,j)
357# endif
358 END IF
359 END DO
360 END DO
361!
362! Eastern edge, closed boundary condition.
363!
364 ELSE IF (lbc(ieast,isuvel,ng)%closed) THEN
365 DO k=1,n(ng)
366 DO j=jstr,jend
367 IF (lbc_apply(ng)%east(j)) THEN
368 u(iend+1,j,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 southern edge.
377!-----------------------------------------------------------------------
378!
379 IF (domain(ng)%Southern_Edge(tile)) THEN
380!
381! Southern edge, implicit upstream radiation condition.
382!
383 IF (lbc(isouth,isuvel,ng)%radiation) THEN
384 DO k=1,n(ng)
385 DO i=istru-1,iend
386 grad(i,jstr-1)=u(i+1,jstr-1,k,nstp)- &
387 & u(i ,jstr-1,k,nstp)
388 grad(i,jstr )=u(i+1,jstr ,k,nstp)- &
389 & u(i ,jstr ,k,nstp)
390 END DO
391 DO i=istru,iend
392 IF (lbc_apply(ng)%south(i)) THEN
393 dudt=u(i,jstr,k,nstp)-u(i,jstr ,k,nout)
394 dude=u(i,jstr,k,nout)-u(i,jstr+1,k,nout)
395
396 IF (lbc(isouth,isuvel,ng)%nudging) THEN
397 IF (lnudgem3clm(ng)) THEN
398 obc_out=0.5_r8* &
399 & (clima(ng)%M3nudgcof(i-1,jstr-1,k)+ &
400 & clima(ng)%M3nudgcof(i ,jstr-1,k))
401 obc_in =obcfac(ng)*obc_out
402 ELSE
403 obc_out=m3obc_out(ng,isouth)
404 obc_in =m3obc_in(ng,isouth)
405 END IF
406 IF ((dudt*dude).lt.0.0_r8) THEN
407 tau=obc_in
408 ELSE
409 tau=obc_out
410 END IF
411# ifdef IMPLICIT_NUDGING
412 IF (tau.gt.0.0_r8) tau=1.0_r8/tau
413# else
414 tau=tau*dt(ng)
415# endif
416 END IF
417
418 IF ((dudt*dude).lt.0.0_r8) dudt=0.0_r8
419 IF ((dudt*(grad(i-1,jstr)+ &
420 & grad(i ,jstr))).gt.0.0_r8) THEN
421 dudx=grad(i-1,jstr)
422 ELSE
423 dudx=grad(i ,jstr)
424 END IF
425 cff=max(dudx*dudx+dude*dude,eps)
426# ifdef RADIATION_2D
427 cx=min(cff,max(dudt*dudx,-cff))
428# else
429 cx=0.0_r8
430# endif
431 ce=dudt*dude
432# if defined CELERITY_WRITE && defined FORWARD_WRITE
433 boundary(ng)%u_south_Cx(i,k)=cx
434 boundary(ng)%u_south_Ce(i,k)=ce
435 boundary(ng)%u_south_C2(i,k)=cff
436# endif
437 u(i,jstr-1,k,nout)=(cff*u(i,jstr-1,k,nstp)+ &
438 & ce *u(i,jstr ,k,nout)- &
439 & max(cx,0.0_r8)*grad(i-1,jstr-1)- &
440 & min(cx,0.0_r8)*grad(i ,jstr-1))/ &
441 & (cff+ce)
442
443 IF (lbc(isouth,isuvel,ng)%nudging) THEN
444# ifdef IMPLICIT_NUDGING
445 phi=dt(ng)/(tau+dt(ng))
446 u(i,jstr-1,k,nout)=(1.0_r8-phi)*u(i,jstr-1,k,nout)+ &
447 & phi*boundary(ng)%u_south(i,k)
448# else
449 u(i,jstr-1,k,nout)=u(i,jstr-1,k,nout)+ &
450 & tau*(boundary(ng)%u_south(i,k)- &
451 & u(i,jstr-1,k,nstp))
452# endif
453 END IF
454# ifdef MASKING
455 u(i,jstr-1,k,nout)=u(i,jstr-1,k,nout)* &
456 & grid(ng)%umask(i,jstr-1)
457# endif
458# ifdef WET_DRY
459 u(i,jstr-1,k,nout)=u(i,jstr-1,k,nout)* &
460 & grid(ng)%umask_wet(i,jstr-1)
461# endif
462 END IF
463 END DO
464 END DO
465!
466! Southern edge, clamped boundary condition.
467!
468 ELSE IF (lbc(isouth,isuvel,ng)%clamped) THEN
469 DO k=1,n(ng)
470 DO i=istru,iend
471 IF (lbc_apply(ng)%south(i)) THEN
472 u(i,jstr-1,k,nout)=boundary(ng)%u_south(i,k)
473# ifdef MASKING
474 u(i,jstr-1,k,nout)=u(i,jstr-1,k,nout)* &
475 & grid(ng)%umask(i,jstr-1)
476# endif
477# ifdef WET_DRY
478 u(i,jstr-1,k,nout)=u(i,jstr-1,k,nout)* &
479 & grid(ng)%umask_wet(i,jstr-1)
480# endif
481 END IF
482 END DO
483 END DO
484!
485! Southern edge, gradient boundary condition.
486!
487 ELSE IF (lbc(isouth,isuvel,ng)%gradient) THEN
488 DO k=1,n(ng)
489 DO i=istru,iend
490 IF (lbc_apply(ng)%south(i)) THEN
491 u(i,jstr-1,k,nout)=u(i,jstr,k,nout)
492# ifdef MASKING
493 u(i,jstr-1,k,nout)=u(i,jstr-1,k,nout)* &
494 & grid(ng)%umask(i,jstr-1)
495# endif
496# ifdef WET_MASK
497 u(i,jstr-1,k,nout)=u(i,jstr-1,k,nout)* &
498 & grid(ng)%umask_wet(i,jstr-1)
499# endif
500 END IF
501 END DO
502 END DO
503!
504! Southern edge, closed boundary condition: free slip (gamma2=1) or
505! no slip (gamma2=-1).
506!
507 ELSE IF (lbc(isouth,isuvel,ng)%closed) THEN
508 IF (ewperiodic(ng)) THEN
509 imin=istru
510 imax=iend
511 ELSE
512 imin=istr
513 imax=iendr
514 END IF
515 DO k=1,n(ng)
516 DO i=imin,imax
517 IF (lbc_apply(ng)%south(i)) THEN
518 u(i,jstr-1,k,nout)=gamma2(ng)*u(i,jstr,k,nout)
519# ifdef MASKING
520 u(i,jstr-1,k,nout)=u(i,jstr-1,k,nout)* &
521 & grid(ng)%umask(i,jstr-1)
522# endif
523# ifdef WET_DRY
524 u(i,jstr-1,k,nout)=u(i,jstr-1,k,nout)* &
525 & grid(ng)%umask_wet(i,jstr-1)
526# endif
527 END IF
528 END DO
529 END DO
530 END IF
531 END IF
532!
533!-----------------------------------------------------------------------
534! Lateral boundary conditions at the northern edge.
535!-----------------------------------------------------------------------
536!
537 IF (domain(ng)%Northern_Edge(tile)) THEN
538!
539! Northern edge, implicit upstream radiation condition.
540!
541 IF (lbc(inorth,isuvel,ng)%radiation) THEN
542 DO k=1,n(ng)
543 DO i=istru-1,iend
544 grad(i,jend )=u(i+1,jend ,k,nstp)- &
545 & u(i ,jend ,k,nstp)
546 grad(i,jend+1)=u(i+1,jend+1,k,nstp)- &
547 & u(i ,jend+1,k,nstp)
548 END DO
549 DO i=istru,iend
550 IF (lbc_apply(ng)%north(i)) THEN
551 dudt=u(i,jend,k,nstp)-u(i,jend ,k,nout)
552 dude=u(i,jend,k,nout)-u(i,jend-1,k,nout)
553
554 IF (lbc(inorth,isuvel,ng)%nudging) THEN
555 IF (lnudgem3clm(ng)) THEN
556 obc_out=0.5_r8* &
557 & (clima(ng)%M3nudgcof(i-1,jend+1,k)+ &
558 & clima(ng)%M3nudgcof(i ,jend+1,k))
559 obc_in =obcfac(ng)*obc_out
560 ELSE
561 obc_out=m3obc_out(ng,inorth)
562 obc_in =m3obc_in(ng,inorth)
563 END IF
564 IF ((dudt*dude).lt.0.0_r8) THEN
565 tau=obc_in
566 ELSE
567 tau=obc_out
568 END IF
569# ifdef IMPLICIT_NUDGING
570 IF (tau.gt.0.0_r8) tau=1.0_r8/tau
571# else
572 tau=tau*dt(ng)
573# endif
574 END IF
575
576 IF ((dudt*dude).lt.0.0_r8) dudt=0.0_r8
577 IF ((dudt*(grad(i-1,jend)+ &
578 & grad(i ,jend))).gt.0.0_r8) THEN
579 dudx=grad(i-1,jend)
580 ELSE
581 dudx=grad(i ,jend)
582 END IF
583 cff=max(dudx*dudx+dude*dude,eps)
584# ifdef RADIATION_2D
585 cx=min(cff,max(dudt*dudx,-cff))
586# else
587 cx=0.0_r8
588# endif
589 ce=dudt*dude
590# if defined CELERITY_WRITE && defined FORWARD_WRITE
591 boundary(ng)%u_north_Cx(i,k)=cx
592 boundary(ng)%u_north_Ce(i,k)=ce
593 boundary(ng)%u_north_C2(i,k)=cff
594# endif
595 u(i,jend+1,k,nout)=(cff*u(i,jend+1,k,nstp)+ &
596 & ce *u(i,jend ,k,nout)- &
597 & max(cx,0.0_r8)*grad(i-1,jend+1)- &
598 & min(cx,0.0_r8)*grad(i ,jend+1))/ &
599 & (cff+ce)
600
601 IF (lbc(inorth,isuvel,ng)%nudging) THEN
602# ifdef IMPLICIT_NUDGING
603 phi=dt(ng)/(tau+dt(ng))
604 u(i,jend+1,k,nout)=(1.0_r8-phi)*u(i,jend+1,k,nout)+ &
605 & phi*boundary(ng)%u_north(i,k)
606# else
607 u(i,jend+1,k,nout)=u(i,jend+1,k,nout)+ &
608 & tau*(boundary(ng)%u_north(i,k)- &
609 & u(i,jend+1,k,nstp))
610# endif
611 END IF
612# ifdef MASKING
613 u(i,jend+1,k,nout)=u(i,jend+1,k,nout)* &
614 & grid(ng)%umask(i,jend+1)
615# endif
616# ifdef WET_DRY
617 u(i,jend+1,k,nout)=u(i,jend+1,k,nout)* &
618 & grid(ng)%umask_wet(i,jend+1)
619# endif
620 END IF
621 END DO
622 END DO
623!
624! Northern edge, clamped boundary condition.
625!
626 ELSE IF (lbc(inorth,isuvel,ng)%clamped) THEN
627 DO k=1,n(ng)
628 DO i=istru,iend
629 IF (lbc_apply(ng)%north(i)) THEN
630 u(i,jend+1,k,nout)=boundary(ng)%u_north(i,k)
631# ifdef MASKING
632 u(i,jend+1,k,nout)=u(i,jend+1,k,nout)* &
633 & grid(ng)%umask(i,jend+1)
634# endif
635# ifdef WET_DRY
636 u(i,jend+1,k,nout)=u(i,jend+1,k,nout)* &
637 & grid(ng)%umask_wet(i,jend+1)
638# endif
639 END IF
640 END DO
641 END DO
642!
643! Northern edge, gradient boundary condition.
644!
645 ELSE IF (lbc(inorth,isuvel,ng)%gradient) THEN
646 DO k=1,n(ng)
647 DO i=istru,iend
648 IF (lbc_apply(ng)%north(i)) THEN
649 u(i,jend+1,k,nout)=u(i,jend,k,nout)
650# ifdef MASKING
651 u(i,jend+1,k,nout)=u(i,jend+1,k,nout)* &
652 & grid(ng)%umask(i,jend+1)
653# endif
654# ifdef WET_DRY
655 u(i,jend+1,k,nout)=u(i,jend+1,k,nout)* &
656 & grid(ng)%umask_wet(i,jend+1)
657# endif
658 END IF
659 END DO
660 END DO
661!
662! Northern edge, closed boundary condition: free slip (gamma2=1) or
663! no slip (gamma2=-1).
664!
665 ELSE IF (lbc(inorth,isuvel,ng)%closed) THEN
666 IF (ewperiodic(ng)) THEN
667 imin=istru
668 imax=iend
669 ELSE
670 imin=istr
671 imax=iendr
672 END IF
673 DO k=1,n(ng)
674 DO i=imin,imax
675 IF (lbc_apply(ng)%north(i)) THEN
676 u(i,jend+1,k,nout)=gamma2(ng)*u(i,jend,k,nout)
677# ifdef MASKING
678 u(i,jend+1,k,nout)=u(i,jend+1,k,nout)* &
679 & grid(ng)%umask(i,jend+1)
680# endif
681# ifdef WET_DRY
682 u(i,jend+1,k,nout)=u(i,jend+1,k,nout)* &
683 & grid(ng)%umask_wet(i,jend+1)
684# endif
685 END IF
686 END DO
687 END DO
688 END IF
689 END IF
690!
691!-----------------------------------------------------------------------
692! Boundary corners.
693!-----------------------------------------------------------------------
694!
695 IF (.not.(ewperiodic(ng).or.nsperiodic(ng))) THEN
696 IF (domain(ng)%SouthWest_Corner(tile)) THEN
697 IF (lbc_apply(ng)%south(istr ).and. &
698 & lbc_apply(ng)%west (jstr-1)) THEN
699 DO k=1,n(ng)
700 u(istr,jstr-1,k,nout)=0.5_r8*(u(istr+1,jstr-1,k,nout)+ &
701 & u(istr ,jstr ,k,nout))
702 END DO
703 END IF
704 END IF
705 IF (domain(ng)%SouthEast_Corner(tile)) THEN
706 IF (lbc_apply(ng)%south(iend+1).and. &
707 & lbc_apply(ng)%east (jstr-1)) THEN
708 DO k=1,n(ng)
709 u(iend+1,jstr-1,k,nout)=0.5_r8*(u(iend ,jstr-1,k,nout)+ &
710 & u(iend+1,jstr ,k,nout))
711 END DO
712 END IF
713 END IF
714 IF (domain(ng)%NorthWest_Corner(tile)) THEN
715 IF (lbc_apply(ng)%north(istr ).and. &
716 & lbc_apply(ng)%west (jend+1)) THEN
717 DO k=1,n(ng)
718 u(istr,jend+1,k,nout)=0.5_r8*(u(istr ,jend ,k,nout)+ &
719 & u(istr+1,jend+1,k,nout))
720 END DO
721 END IF
722 END IF
723 IF (domain(ng)%NorthEast_Corner(tile)) THEN
724 IF (lbc_apply(ng)%north(iend+1).and. &
725 & lbc_apply(ng)%east (jend+1)) THEN
726 DO k=1,n(ng)
727 u(iend+1,jend+1,k,nout)=0.5_r8*(u(iend+1,jend ,k,nout)+ &
728 & u(iend ,jend+1,k,nout))
729 END DO
730 END IF
731 END IF
732 END IF
733
734 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 isuvel
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

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::inorth, mod_scalars::isouth, mod_ncparam::isuvel, mod_scalars::iwest, mod_param::lbc, mod_boundary::lbc_apply, mod_scalars::lnudgem3clm, mod_scalars::m3obc_in, mod_scalars::m3obc_out, mod_param::n, mod_scalars::nsperiodic, and mod_scalars::obcfac.

Referenced by ini_fields_mod::ini_fields_tile(), ini_adjust_mod::ini_perturb_tile(), step3d_uv_mod::step3d_uv_tile(), and u3dbc().

Here is the caller graph for this function: