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

Functions/Subroutines

subroutine v3dbc (ng, tile, nout)
 
subroutine, public v3dbc_tile (ng, tile, lbi, ubi, lbj, ubj, ubk, imins, imaxs, jmins, jmaxs, nstp, nout, v)
 

Function/Subroutine Documentation

◆ v3dbc()

subroutine v3dbc_mod::v3dbc ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) nout )
private

Definition at line 25 of file v3dbc_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 v3dbc_tile (ng, tile, &
41 & lbi, ubi, lbj, ubj, n(ng), &
42 & imins, imaxs, jmins, jmaxs, &
43 & nstp(ng), nout, &
44 & ocean(ng) % v)
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 v3dbc_tile().

Here is the call graph for this function:

◆ v3dbc_tile()

subroutine public v3dbc_mod::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) v )

Definition at line 50 of file v3dbc_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) :: v(LBi:,LBj:,:,:)
73# else
74 real(r8), intent(inout) :: v(LBi:UBi,LBj:UBj,UBk,2)
75# endif
76!
77! Local variable declarations.
78!
79 integer :: Jmin, Jmax
80 integer :: i, j, k
81
82 real(r8), parameter :: eps = 1.0e-20_r8
83
84 real(r8) :: Ce, Cx, cff, dVde, dVdt, dVdx
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 southern edge.
93!-----------------------------------------------------------------------
94!
95 IF (domain(ng)%Southern_Edge(tile)) THEN
96!
97! Southern edge, implicit upstream radiation condition.
98!
99 IF (lbc(isouth,isvvel,ng)%radiation) THEN
100 DO k=1,n(ng)
101 DO i=istr,iend+1
102 grad(i,jstr )=v(i ,jstr ,k,nstp)- &
103 & v(i-1,jstr ,k,nstp)
104 grad(i,jstr+1)=v(i ,jstr+1,k,nstp)- &
105 & v(i-1,jstr+1,k,nstp)
106 END DO
107 DO i=istr,iend
108 IF (lbc_apply(ng)%south(i)) THEN
109 dvdt=v(i,jstr+1,k,nstp)-v(i,jstr+1,k,nout)
110 dvde=v(i,jstr+1,k,nout)-v(i,jstr+2,k,nout)
111
112 IF (lbc(isouth,isvvel,ng)%nudging) THEN
113 IF (lnudgem3clm(ng)) THEN
114 obc_out=0.5_r8* &
115 & (clima(ng)%M3nudgcof(i,jstr-1,k)+ &
116 & clima(ng)%M3nudgcof(i,jstr ,k))
117 obc_in =obcfac(ng)*obc_out
118 ELSE
119 obc_out=m3obc_out(ng,isouth)
120 obc_in =m3obc_in(ng,isouth)
121 END IF
122 IF ((dvdt*dvde).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 ((dvdt*dvde).lt.0.0_r8) dvdt=0.0_r8
135 IF ((dvdt*(grad(i ,jstr+1)+ &
136 & grad(i+1,jstr+1))).gt.0.0_r8) THEN
137 dvdx=grad(i ,jstr+1)
138 ELSE
139 dvdx=grad(i+1,jstr+1)
140 END IF
141 cff=max(dvdx*dvdx+dvde*dvde,eps)
142# ifdef RADIATION_2D
143 cx=min(cff,max(dvdt*dvdx,-cff))
144# else
145 cx=0.0_r8
146# endif
147 ce=dvdt*dvde
148# if defined CELERITY_WRITE && defined FORWARD_WRITE
149 boundary(ng)%v_south_Cx(i,k)=cx
150 boundary(ng)%v_south_Ce(i,k)=ce
151 boundary(ng)%v_south_C2(i,k)=cff
152# endif
153 v(i,jstr,k,nout)=(cff*v(i,jstr ,k,nstp)+ &
154 & ce *v(i,jstr+1,k,nout)- &
155 & max(cx,0.0_r8)*grad(i ,jstr)- &
156 & min(cx,0.0_r8)*grad(i+1,jstr))/ &
157 & (cff+ce)
158
159 IF (lbc(isouth,isvvel,ng)%nudging) THEN
160# ifdef IMPLICIT_NUDGING
161 phi=dt(ng)/(tau+dt(ng))
162 v(i,jstr,k,nout)=(1.0_r8-phi)*v(i,jstr,k,nout)+ &
163 & phi*boundary(ng)%v_south(i,k)
164# else
165 v(i,jstr,k,nout)=v(i,jstr,k,nout)+ &
166 & tau*(boundary(ng)%v_south(i,k)- &
167 & v(i,jstr,k,nstp))
168# endif
169 END IF
170# ifdef MASKING
171 v(i,jstr,k,nout)=v(i,jstr,k,nout)* &
172 & grid(ng)%vmask(i,jstr)
173# endif
174# ifdef WET_DRY
175 v(i,jstr,k,nout)=v(i,jstr,k,nout)* &
176 & grid(ng)%vmask_wet(i,jstr)
177# endif
178 END IF
179 END DO
180 END DO
181!
182! Southern edge, clamped boundary condition.
183!
184 ELSE IF (lbc(isouth,isvvel,ng)%clamped) THEN
185 DO k=1,n(ng)
186 DO i=istr,iend
187 IF (lbc_apply(ng)%south(i)) THEN
188 v(i,jstr,k,nout)=boundary(ng)%v_south(i,k)
189# ifdef MASKING
190 v(i,jstr,k,nout)=v(i,jstr,k,nout)* &
191 & grid(ng)%vmask(i,jstr)
192# endif
193# ifdef WET_DRY
194 v(i,jstr,k,nout)=v(i,jstr,k,nout)* &
195 & grid(ng)%vmask_wet(i,jstr)
196# endif
197 END IF
198 END DO
199 END DO
200!
201! Southern edge, gradient boundary condition.
202!
203 ELSE IF (lbc(isouth,isvvel,ng)%gradient) THEN
204 DO k=1,n(ng)
205 DO i=istr,iend
206 IF (lbc_apply(ng)%south(i)) THEN
207 v(i,jstr,k,nout)=v(i,jstr+1,k,nout)
208# ifdef MASKING
209 v(i,jstr,k,nout)=v(i,jstr,k,nout)* &
210 & grid(ng)%vmask(i,jstr)
211# endif
212# ifdef WET_DRY
213 v(i,jstr,k,nout)=v(i,jstr,k,nout)* &
214 & grid(ng)%vmask_wet(i,jstr)
215# endif
216 END IF
217 END DO
218 END DO
219!
220! Southern edge, closed boundary condition.
221!
222 ELSE IF (lbc(isouth,isvvel,ng)%closed) THEN
223 DO k=1,n(ng)
224 DO i=istr,iend
225 IF (lbc_apply(ng)%south(i)) THEN
226 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 (lbc(inorth,isvvel,ng)%radiation) THEN
242 DO k=1,n(ng)
243 DO i=istr,iend+1
244 grad(i,jend )=v(i ,jend ,k,nstp)- &
245 & v(i-1,jend ,k,nstp)
246 grad(i,jend+1)=v(i ,jend+1,k,nstp)- &
247 & v(i-1,jend+1,k,nstp)
248 END DO
249 DO i=istr,iend
250 IF (lbc_apply(ng)%north(i)) THEN
251 dvdt=v(i,jend,k,nstp)-v(i,jend ,k,nout)
252 dvde=v(i,jend,k,nout)-v(i,jend-1,k,nout)
253
254 IF (lbc(inorth,isvvel,ng)%nudging) THEN
255 IF (lnudgem3clm(ng)) THEN
256 obc_out=0.5_r8* &
257 & (clima(ng)%M3nudgcof(i,jend ,k)+ &
258 & clima(ng)%M3nudgcof(i,jend+1,k))
259 obc_in =obcfac(ng)*obc_out
260 ELSE
261 obc_out=m3obc_out(ng,inorth)
262 obc_in =m3obc_in(ng,inorth)
263 END IF
264 IF ((dvdt*dvde).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 ((dvdt*dvde).lt.0.0_r8) dvdt=0.0_r8
277 IF ((dvdt*(grad(i ,jend)+ &
278 & grad(i+1,jend))).gt.0.0_r8) THEN
279 dvdx=grad(i ,jend)
280 ELSE
281 dvdx=grad(i+1,jend)
282 END IF
283 cff=max(dvdx*dvdx+dvde*dvde,eps)
284# ifdef RADIATION_2D
285 cx=min(cff,max(dvdt*dvdx,-cff))
286# else
287 cx=0.0_r8
288# endif
289 ce=dvdt*dvde
290# if defined CELERITY_WRITE && defined FORWARD_WRITE
291 boundary(ng)%v_north_Cx(i,k)=cx
292 boundary(ng)%v_north_Ce(i,k)=ce
293 boundary(ng)%v_north_C2(i,k)=cff
294# endif
295 v(i,jend+1,k,nout)=(cff*v(i,jend+1,k,nstp)+ &
296 & ce *v(i,jend ,k,nout)- &
297 & max(cx,0.0_r8)*grad(i ,jend+1)- &
298 & min(cx,0.0_r8)*grad(i+1,jend+1))/ &
299 & (cff+ce)
300
301 IF (lbc(inorth,isvvel,ng)%nudging) THEN
302# ifdef IMPLICIT_NUDGING
303 phi=dt(ng)/(tau+dt(ng))
304 v(i,jend+1,k,nout)=(1.0_r8-phi)*v(i,jend+1,k,nout)+ &
305 & phi*boundary(ng)%v_north(i,k)
306# else
307 v(i,jend+1,k,nout)=v(i,jend+1,k,nout)+ &
308 & tau*(boundary(ng)%v_north(i,k)- &
309 & v(i,jend+1,k,nstp))
310# endif
311 END IF
312# ifdef MASKING
313 v(i,jend+1,k,nout)=v(i,jend+1,k,nout)* &
314 & grid(ng)%vmask(i,jend+1)
315# endif
316# ifdef WET_DRY
317 v(i,jend+1,k,nout)=v(i,jend+1,k,nout)* &
318 & grid(ng)%vmask_wet(i,jend+1)
319# endif
320 END IF
321 END DO
322 END DO
323!
324! Northern edge, clamped boundary condition.
325!
326 ELSE IF (lbc(inorth,isvvel,ng)%clamped) THEN
327 DO k=1,n(ng)
328 DO i=istr,iend
329 IF (lbc_apply(ng)%north(i)) THEN
330 v(i,jend+1,k,nout)=boundary(ng)%v_north(i,k)
331# ifdef MASKING
332 v(i,jend+1,k,nout)=v(i,jend+1,k,nout)* &
333 & grid(ng)%vmask(i,jend+1)
334# endif
335# ifdef WET_DRY
336 v(i,jend+1,k,nout)=v(i,jend+1,k,nout)* &
337 & grid(ng)%vmask_wet(i,jend+1)
338# endif
339 END IF
340 END DO
341 END DO
342!
343! Northern edge, gradient boundary condition.
344!
345 ELSE IF (lbc(inorth,isvvel,ng)%gradient) THEN
346 DO k=1,n(ng)
347 DO i=istr,iend
348 IF (lbc_apply(ng)%north(i)) THEN
349 v(i,jend+1,k,nout)=v(i,jend,k,nout)
350# ifdef MASKING
351 v(i,jend+1,k,nout)=v(i,jend+1,k,nout)* &
352 & grid(ng)%vmask(i,jend+1)
353# endif
354# ifdef WET_DRY
355 v(i,jend+1,k,nout)=v(i,jend+1,k,nout)* &
356 & grid(ng)%vmask_wet(i,jend+1)
357# endif
358 END IF
359 END DO
360 END DO
361!
362! Northern edge, closed boundary condition.
363!
364 ELSE IF (lbc(inorth,isvvel,ng)%closed) THEN
365 DO k=1,n(ng)
366 DO i=istr,iend
367 IF (lbc_apply(ng)%north(i)) THEN
368 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 (lbc(iwest,isvvel,ng)%radiation) THEN
384 DO k=1,n(ng)
385 DO j=jstrv-1,jend
386 grad(istr-1,j)=v(istr-1,j+1,k,nstp)- &
387 & v(istr-1,j ,k,nstp)
388 grad(istr ,j)=v(istr ,j+1,k,nstp)- &
389 & v(istr ,j ,k,nstp)
390 END DO
391 DO j=jstrv,jend
392 IF (lbc_apply(ng)%west(j)) THEN
393 dvdt=v(istr,j,k,nstp)-v(istr ,j,k,nout)
394 dvdx=v(istr,j,k,nout)-v(istr+1,j,k,nout)
395
396 IF (lbc(iwest,isvvel,ng)%nudging) THEN
397 IF (lnudgem3clm(ng)) THEN
398 obc_out=0.5_r8* &
399 & (clima(ng)%M3nudgcof(istr-1,j-1,k)+ &
400 & clima(ng)%M3nudgcof(istr-1,j ,k))
401 obc_in =obcfac(ng)*obc_out
402 ELSE
403 obc_out=m3obc_out(ng,iwest)
404 obc_in =m3obc_in(ng,iwest)
405 END IF
406 IF ((dvdt*dvdx).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 ((dvdt*dvdx).lt.0.0_r8) dvdt=0.0_r8
419 IF ((dvdt*(grad(istr,j-1)+ &
420 & grad(istr,j ))).gt.0.0_r8) THEN
421 dvde=grad(istr,j-1)
422 ELSE
423 dvde=grad(istr,j )
424 END IF
425 cff=max(dvdx*dvdx+dvde*dvde,eps)
426 cx=dvdt*dvdx
427# ifdef RADIATION_2D
428 ce=min(cff,max(dvdt*dvde,-cff))
429# else
430 ce=0.0_r8
431# endif
432# if defined CELERITY_WRITE && defined FORWARD_WRITE
433 boundary(ng)%v_west_Cx(j,k)=cx
434 boundary(ng)%v_west_Ce(j,k)=ce
435 boundary(ng)%v_west_C2(j,k)=cff
436# endif
437 v(istr-1,j,k,nout)=(cff*v(istr-1,j,k,nstp)+ &
438 & cx *v(istr ,j,k,nout)- &
439 & max(ce,0.0_r8)*grad(istr-1,j-1)- &
440 & min(ce,0.0_r8)*grad(istr-1,j ))/ &
441 & (cff+cx)
442
443 IF (lbc(iwest,isvvel,ng)%nudging) THEN
444# ifdef IMPLICIT_NUDGING
445 phi=dt(ng)/(tau+dt(ng))
446 v(istr-1,j,k,nout)=(1.0_r8-phi)*v(istr-1,j,k,nout)+ &
447 & phi*boundary(ng)%v_west(j,k)
448# else
449 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)+ &
450 & tau*(boundary(ng)%v_west(j,k)- &
451 & v(istr-1,j,k,nstp))
452# endif
453 END IF
454# ifdef MASKING
455 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
456 & grid(ng)%vmask(istr-1,j)
457# endif
458# ifdef WET_DRY
459 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
460 & grid(ng)%vmask_wet(istr-1,j)
461# endif
462 END IF
463 END DO
464 END DO
465!
466! Western edge, clamped boundary condition.
467!
468 ELSE IF (lbc(iwest,isvvel,ng)%clamped) THEN
469 DO k=1,n(ng)
470 DO j=jstrv,jend
471 IF (lbc_apply(ng)%west(j)) THEN
472 v(istr-1,j,k,nout)=boundary(ng)%v_west(j,k)
473# ifdef MASKING
474 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
475 & grid(ng)%vmask(istr-1,j)
476# endif
477# ifdef WET_DRY
478 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
479 & grid(ng)%vmask_wet(istr-1,j)
480# endif
481 END IF
482 END DO
483 END DO
484!
485! Western edge, gradient boundary condition.
486!
487 ELSE IF (lbc(iwest,isvvel,ng)%gradient) THEN
488 DO k=1,n(ng)
489 DO j=jstrv,jend
490 IF (lbc_apply(ng)%west(j)) THEN
491 v(istr-1,j,k,nout)=v(istr,j,k,nout)
492# ifdef MASKING
493 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
494 & grid(ng)%vmask(istr-1,j)
495# endif
496# ifdef WET_DRY
497 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
498 & grid(ng)%vmask_wet(istr-1,j)
499# endif
500 END IF
501 END DO
502 END DO
503!
504! Western edge, closed boundary condition: free slip (gamma2=1) or
505! no slip (gamma2=-1).
506!
507 ELSE IF (lbc(iwest,isvvel,ng)%closed) THEN
508 IF (nsperiodic(ng)) THEN
509 jmin=jstrv
510 jmax=jend
511 ELSE
512 jmin=jstr
513 jmax=jendr
514 END IF
515 DO k=1,n(ng)
516 DO j=jmin,jmax
517 IF (lbc_apply(ng)%west(j)) THEN
518 v(istr-1,j,k,nout)=gamma2(ng)*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# endif
523# ifdef WET_DRY
524 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
525 & grid(ng)%vmask_wet(istr-1,j)
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 eastern edge.
535!-----------------------------------------------------------------------
536!
537 IF (domain(ng)%Eastern_Edge(tile)) THEN
538!
539! Eastern edge, implicit upstream radiation condition.
540!
541 IF (lbc(ieast,isvvel,ng)%radiation) THEN
542 DO k=1,n(ng)
543 DO j=jstrv-1,jend
544 grad(iend ,j)=v(iend ,j+1,k,nstp)- &
545 & v(iend ,j ,k,nstp)
546 grad(iend+1,j)=v(iend+1,j+1,k,nstp)- &
547 & v(iend+1,j ,k,nstp)
548 END DO
549 DO j=jstrv,jend
550 IF (lbc_apply(ng)%east(j)) THEN
551 dvdt=v(iend,j,k,nstp)-v(iend ,j,k,nout)
552 dvdx=v(iend,j,k,nout)-v(iend-1,j,k,nout)
553
554 IF (lbc(ieast,isvvel,ng)%nudging) THEN
555 IF (lnudgem3clm(ng)) THEN
556 obc_out=0.5_r8* &
557 & (clima(ng)%M3nudgcof(iend+1,j-1,k)+ &
558 & clima(ng)%M3nudgcof(iend+1,j ,k))
559 obc_in =obcfac(ng)*obc_out
560 ELSE
561 obc_out=m3obc_out(ng,ieast)
562 obc_in =m3obc_in(ng,ieast)
563 END IF
564 IF ((dvdt*dvdx).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 ((dvdt*dvdx).lt.0.0_r8) dvdt=0.0_r8
577 IF ((dvdt*(grad(iend,j-1)+ &
578 & grad(iend,j ))).gt.0.0_r8) THEN
579 dvde=grad(iend,j-1)
580 ELSE
581 dvde=grad(iend,j )
582 END IF
583 cff=max(dvdx*dvdx+dvde*dvde,eps)
584 cx=dvdt*dvdx
585# ifdef RADIATION_2D
586 ce=min(cff,max(dvdt*dvde,-cff))
587# else
588 ce=0.0_r8
589# endif
590# if defined CELERITY_WRITE && defined FORWARD_WRITE
591 boundary(ng)%v_east_Cx(j,k)=cx
592 boundary(ng)%v_east_Ce(j,k)=ce
593 boundary(ng)%v_east_C2(j,k)=cff
594# endif
595 v(iend+1,j,k,nout)=(cff*v(iend+1,j,k,nstp)+ &
596 & cx *v(iend ,j,k,nout)- &
597 & max(ce,0.0_r8)*grad(iend+1,j-1)- &
598 & min(ce,0.0_r8)*grad(iend+1,j ))/ &
599 & (cff+cx)
600
601 IF (lbc(ieast,isvvel,ng)%nudging) THEN
602# ifdef IMPLICIT_NUDGING
603 phi=dt(ng)/(tau+dt(ng))
604 v(iend+1,j,k,nout)=(1.0_r8-phi)*v(iend+1,j,k,nout)+ &
605 & phi*boundary(ng)%v_east(j,k)
606# else
607 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)+ &
608 & tau*(boundary(ng)%v_east(j,k)- &
609 & v(iend+1,j,k,nstp))
610# endif
611 END IF
612# ifdef MASKING
613 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
614 & grid(ng)%vmask(iend+1,j)
615# endif
616# ifdef WET_DRY
617 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
618 & grid(ng)%vmask_wet(iend+1,j)
619# endif
620 END IF
621 END DO
622 END DO
623!
624! Eastern edge, clamped boundary condition.
625!
626 ELSE IF (lbc(ieast,isvvel,ng)%clamped) THEN
627 DO k=1,n(ng)
628 DO j=jstrv,jend
629 IF (lbc_apply(ng)%east(j)) THEN
630 v(iend+1,j,k,nout)=boundary(ng)%v_east(j,k)
631# ifdef MASKING
632 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
633 & grid(ng)%vmask(iend+1,j)
634# endif
635# ifdef WET_DRY
636 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
637 & grid(ng)%vmask_wet(iend+1,j)
638# endif
639 END IF
640 END DO
641 END DO
642!
643! Eastern edge, gradient boundary condition.
644!
645 ELSE IF (lbc(ieast,isvvel,ng)%gradient) THEN
646 DO k=1,n(ng)
647 DO j=jstrv,jend
648 IF (lbc_apply(ng)%east(j)) THEN
649 v(iend+1,j,k,nout)=v(iend,j,k,nout)
650# ifdef MASKING
651 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
652 & grid(ng)%vmask(iend+1,j)
653# endif
654# ifdef WET_DRY
655 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
656 & grid(ng)%vmask_wet(iend+1,j)
657# endif
658 END IF
659 END DO
660 END DO
661!
662! Eastern edge, closed boundary condition: free slip (gamma2=1) or
663! no slip (gamma2=-1).
664!
665 ELSE IF (lbc(ieast,isvvel,ng)%closed) THEN
666 IF (nsperiodic(ng)) THEN
667 jmin=jstrv
668 jmax=jend
669 ELSE
670 jmin=jstr
671 jmax=jendr
672 END IF
673 DO k=1,n(ng)
674 DO j=jmin,jmax
675 IF (lbc_apply(ng)%east(j)) THEN
676 v(iend+1,j,k,nout)=gamma2(ng)*v(iend,j,k,nout)
677# ifdef MASKING
678 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
679 & grid(ng)%vmask(iend+1,j)
680# endif
681# ifdef WET_DRY
682 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
683 & grid(ng)%vmask_wet(iend+1,j)
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-1).and. &
698 & lbc_apply(ng)%west (jstr )) THEN
699 DO k=1,n(ng)
700 v(istr-1,jstr,k,nout)=0.5_r8*(v(istr ,jstr ,k,nout)+ &
701 & v(istr-1,jstr+1,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 )) THEN
708 DO k=1,n(ng)
709 v(iend+1,jstr,k,nout)=0.5_r8*(v(iend ,jstr ,k,nout)+ &
710 & v(iend+1,jstr+1,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-1).and. &
716 & lbc_apply(ng)%west (jend+1)) THEN
717 DO k=1,n(ng)
718 v(istr-1,jend+1,k,nout)=0.5_r8*(v(istr-1,jend ,k,nout)+ &
719 & v(istr ,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 v(iend+1,jend+1,k,nout)=0.5_r8*(v(iend+1,jend ,k,nout)+ &
728 & v(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 isvvel
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::isvvel, 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 v3dbc().

Here is the caller graph for this function: