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

Functions/Subroutines

subroutine, public rp_u3dbc (ng, tile, nout)
 
subroutine, public rp_u3dbc_tile (ng, tile, lbi, ubi, lbj, ubj, ubk, imins, imaxs, jmins, jmaxs, nstp, nout, tl_u)
 

Function/Subroutine Documentation

◆ rp_u3dbc()

subroutine, public rp_u3dbc_mod::rp_u3dbc ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) nout )

Definition at line 28 of file rp_u3dbc_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, nout
38!
39! Local variable declarations.
40!
41# include "tile.h"
42!
43 CALL rp_u3dbc_tile (ng, tile, &
44 & lbi, ubi, lbj, ubj, n(ng), &
45 & imins, imaxs, jmins, jmaxs, &
46 & nstp(ng), nout, &
47 & ocean(ng) % tl_u)
48
49 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 rp_u3dbc_tile().

Here is the call graph for this function:

◆ rp_u3dbc_tile()

subroutine, public rp_u3dbc_mod::rp_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) tl_u )

Definition at line 53 of file rp_u3dbc_im.F.

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

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

Referenced by rp_ini_fields_mod::rp_ini_fields_tile(), rp_step3d_uv_mod::rp_step3d_uv_tile(), and rp_u3dbc().

Here is the caller graph for this function: