25 SUBROUTINE v3dbc (ng, tile, nout)
34 integer,
intent(in) :: ng, tile, nout
41 & lbi, ubi, lbj, ubj,
n(ng), &
42 & imins, imaxs, jmins, jmaxs, &
51 & LBi, UBi, LBj, UBj, UBk, &
52 & IminS, ImaxS, JminS, JmaxS, &
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
72 real(r8),
intent(inout) :: v(LBi:,LBj:,:,:)
74 real(r8),
intent(inout) :: v(LBi:UBi,LBj:UBj,UBk,2)
82 real(r8),
parameter :: eps = 1.0e-20_r8
84 real(r8) :: Ce, Cx, cff, dVde, dVdt, dVdx
85 real(r8) :: obc_in, obc_out, tau
87 real(r8),
dimension(IminS:ImaxS,JminS:JmaxS) :: grad
89# include "set_bounds.h"
95 IF (
domain(ng)%Southern_Edge(tile))
THEN
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)
109 dvdt=v(i,jstr+1,k,nstp)-v(i,jstr+1,k,nout)
110 dvde=v(i,jstr+1,k,nstp)-v(i,jstr+2,k,nstp)
115 & (
clima(ng)%M3nudgcof(i,jstr-1,k)+ &
116 &
clima(ng)%M3nudgcof(i,jstr ,k))
117 obc_in =
obcfac(ng)*obc_out
122 IF ((dvdt*dvde).lt.0.0_r8)
THEN
130 IF ((dvdt*dvde).lt.0.0_r8) dvdt=0.0_r8
131 IF ((dvdt*(grad(i ,jstr+1)+ &
132 & grad(i+1,jstr+1))).gt.0.0_r8)
THEN
135 dvdx=grad(i+1,jstr+1)
137 cff=dvdt/max(dvdx*dvdx+dvde*dvde,eps)
139 cx=min(1.0_r8,max(cff*dvdx,-1.0_r8))
143 ce=min(1.0_r8,cff*dvde)
144# if defined CELERITY_WRITE && defined FORWARD_WRITE
149 v(i,jstr,k,nout)=(1.0_r8-ce)*v(i,jstr,k,nstp)+ &
150 & ce*v(i,jstr+1,k,nstp)- &
151 & max(cx,0.0_r8)*grad(i ,jstr)- &
152 & min(cx,0.0_r8)*grad(i+1,jstr)
155 v(i,jstr,k,nout)=v(i,jstr,k,nout)+ &
160 v(i,jstr,k,nout)=v(i,jstr,k,nout)* &
161 &
grid(ng)%vmask(i,jstr)
164 v(i,jstr,k,nout)=v(i,jstr,k,nout)* &
165 &
grid(ng)%vmask_wet(i,jstr)
177 v(i,jstr,k,nout)=
boundary(ng)%v_south(i,k)
179 v(i,jstr,k,nout)=v(i,jstr,k,nout)* &
180 &
grid(ng)%vmask(i,jstr)
183 v(i,jstr,k,nout)=v(i,jstr,k,nout)* &
184 &
grid(ng)%vmask_wet(i,jstr)
196 v(i,jstr,k,nout)=v(i,jstr+1,k,nout)
198 v(i,jstr,k,nout)=v(i,jstr,k,nout)* &
199 &
grid(ng)%vmask(i,jstr)
202 v(i,jstr,k,nout)=v(i,jstr,k,nout)* &
203 &
grid(ng)%vmask_wet(i,jstr)
215 v(i,jstr,k,nout)=0.0_r8
226 IF (
domain(ng)%Northern_Edge(tile))
THEN
233 grad(i,jend )=v(i ,jend ,k,nstp)- &
234 & v(i-1,jend ,k,nstp)
235 grad(i,jend+1)=v(i ,jend+1,k,nstp)- &
236 & v(i-1,jend+1,k,nstp)
240 dvdt=v(i,jend,k,nstp)-v(i,jend ,k,nout)
241 dvde=v(i,jend,k,nstp)-v(i,jend-1,k,nstp)
246 & (
clima(ng)%M3nudgcof(i,jend ,k)+ &
247 &
clima(ng)%M3nudgcof(i,jend+1,k))
248 obc_in =
obcfac(ng)*obc_out
253 IF ((dvdt*dvde).lt.0.0_r8)
THEN
261 IF ((dvdt*dvde).lt.0.0_r8) dvdt=0.0_r8
262 IF ((dvdt*(grad(i ,jend)+ &
263 & grad(i+1,jend))).gt.0.0_r8)
THEN
268 cff=dvdt/max(dvdx*dvdx+dvde*dvde,eps)
270 cx=min(1.0_r8,max(cff*dvdx,-1.0_r8))
274 ce=min(1.0_r8,cff*dvde)
275# if defined CELERITY_WRITE && defined FORWARD_WRITE
280 v(i,jend+1,k,nout)=(1.0_r8-ce)*v(i,jend+1,k,nstp)+ &
281 & ce*v(i,jend,k,nstp)- &
282 & max(cx,0.0_r8)*grad(i ,jend+1)- &
283 & min(cx,0.0_r8)*grad(i+1,jend+1)
286 v(i,jend+1,k,nout)=v(i,jend+1,k,nout)+ &
288 & v(i,jend+1,k,nstp))
291 v(i,jend+1,k,nout)=v(i,jend+1,k,nout)* &
292 &
grid(ng)%vmask(i,jend+1)
295 v(i,jend+1,k,nout)=v(i,jend+1,k,nout)* &
296 &
grid(ng)%vmask_wet(i,jend+1)
308 v(i,jend+1,k,nout)=
boundary(ng)%v_north(i,k)
310 v(i,jend+1,k,nout)=v(i,jend+1,k,nout)* &
311 &
grid(ng)%vmask(i,jend+1)
314 v(i,jend+1,k,nout)=v(i,jend+1,k,nout)* &
315 &
grid(ng)%vmask_wet(i,jend+1)
327 v(i,jend+1,k,nout)=v(i,jend,k,nout)
329 v(i,jend+1,k,nout)=v(i,jend+1,k,nout)* &
330 &
grid(ng)%vmask(i,jend+1)
333 v(i,jend+1,k,nout)=v(i,jend+1,k,nout)* &
334 &
grid(ng)%vmask_wet(i,jend+1)
346 v(i,jend+1,k,nout)=0.0_r8
357 IF (
domain(ng)%Western_Edge(tile))
THEN
364 grad(istr-1,j)=v(istr-1,j+1,k,nstp)- &
365 & v(istr-1,j ,k,nstp)
366 grad(istr ,j)=v(istr ,j+1,k,nstp)- &
371 dvdt=v(istr,j,k,nstp)-v(istr ,j,k,nout)
372 dvdx=v(istr,j,k,nstp)-v(istr+1,j,k,nstp)
377 & (
clima(ng)%M3nudgcof(istr-1,j-1,k)+ &
378 &
clima(ng)%M3nudgcof(istr-1,j ,k))
379 obc_in =
obcfac(ng)*obc_out
384 IF ((dvdt*dvdx).lt.0.0_r8)
THEN
392 IF ((dvdt*dvdx).lt.0.0_r8) dvdt=0.0_r8
393 IF ((dvdt*(grad(istr,j-1)+ &
394 & grad(istr,j ))).gt.0.0_r8)
THEN
399 cff=dvdt/max(dvdx*dvdx+dvde*dvde,eps)
400 cx=min(1.0_r8,cff*dvdx)
402 ce=min(1.0_r8,max(cff*dvde,-1.0_r8))
406# if defined CELERITY_WRITE && defined FORWARD_WRITE
411 v(istr-1,j,k,nout)=(1.0_r8-cx)*v(istr-1,j,k,nstp)+ &
412 & cx*v(istr,j,k,nstp)- &
413 & max(ce,0.0_r8)*grad(istr-1,j-1)- &
414 & min(ce,0.0_r8)*grad(istr-1,j )
417 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)+ &
419 & v(istr-1,j,k,nstp))
422 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
423 &
grid(ng)%vmask(istr-1,j)
426 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
427 &
grid(ng)%vmask_wet(istr-1,j)
439 v(istr-1,j,k,nout)=
boundary(ng)%v_west(j,k)
441 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
442 &
grid(ng)%vmask(istr-1,j)
445 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
446 &
grid(ng)%vmask_wet(istr-1,j)
458 v(istr-1,j,k,nout)=v(istr,j,k,nout)
460 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
461 &
grid(ng)%vmask(istr-1,j)
464 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
465 &
grid(ng)%vmask_wet(istr-1,j)
485 v(istr-1,j,k,nout)=
gamma2(ng)*v(istr,j,k,nout)
487 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
488 &
grid(ng)%vmask(istr-1,j)
491 v(istr-1,j,k,nout)=v(istr-1,j,k,nout)* &
492 &
grid(ng)%vmask_wet(istr-1,j)
504 IF (
domain(ng)%Eastern_Edge(tile))
THEN
511 grad(iend ,j)=v(iend ,j+1,k,nstp)- &
513 grad(iend+1,j)=v(iend+1,j+1,k,nstp)- &
514 & v(iend+1,j ,k,nstp)
518 dvdt=v(iend,j,k,nstp)-v(iend ,j,k,nout)
519 dvdx=v(iend,j,k,nstp)-v(iend-1,j,k,nstp)
524 & (
clima(ng)%M3nudgcof(iend+1,j-1,k)+ &
525 &
clima(ng)%M3nudgcof(iend+1,j ,k))
526 obc_in =
obcfac(ng)*obc_out
531 IF ((dvdt*dvdx).lt.0.0_r8)
THEN
539 IF ((dvdt*dvdx).lt.0.0_r8) dvdt=0.0_r8
540 IF ((dvdt*(grad(iend,j-1)+ &
541 & grad(iend,j ))).gt.0.0_r8)
THEN
546 cff=dvdt/max(dvdx*dvdx+dvde*dvde,eps)
547 cx=min(1.0_r8,cff*dvdx)
549 ce=min(1.0_r8,max(cff*dvde,-1.0_r8))
553# if defined CELERITY_WRITE && defined FORWARD_WRITE
558 v(iend+1,j,k,nout)=(1.0_r8-cx)*v(iend+1,j,k,nstp)+ &
559 & cx*v(iend,j,k,nstp)- &
560 & max(ce,0.0_r8)*grad(iend+1,j-1)- &
561 & min(ce,0.0_r8)*grad(iend+1,j )
564 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)+ &
566 & v(iend+1,j,k,nstp))
569 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
570 &
grid(ng)%vmask(iend+1,j)
573 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
574 &
grid(ng)%vmask_wet(iend+1,j)
586 v(iend+1,j,k,nout)=
boundary(ng)%v_east(j,k)
588 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
589 &
grid(ng)%vmask(iend+1,j)
592 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
593 &
grid(ng)%vmask_wet(iend+1,j)
605 v(iend+1,j,k,nout)=v(iend,j,k,nout)
607 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
608 &
grid(ng)%vmask(iend+1,j)
611 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
612 &
grid(ng)%vmask_wet(iend+1,j)
632 v(iend+1,j,k,nout)=
gamma2(ng)*v(iend,j,k,nout)
634 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
635 &
grid(ng)%vmask(iend+1,j)
638 v(iend+1,j,k,nout)=v(iend+1,j,k,nout)* &
639 &
grid(ng)%vmask_wet(iend+1,j)
652 IF (
domain(ng)%SouthWest_Corner(tile))
THEN
656 v(istr-1,jstr,k,nout)=0.5_r8*(v(istr ,jstr ,k,nout)+ &
657 & v(istr-1,jstr+1,k,nout))
661 IF (
domain(ng)%SouthEast_Corner(tile))
THEN
665 v(iend+1,jstr,k,nout)=0.5_r8*(v(iend ,jstr ,k,nout)+ &
666 & v(iend+1,jstr+1,k,nout))
670 IF (
domain(ng)%NorthWest_Corner(tile))
THEN
674 v(istr-1,jend+1,k,nout)=0.5_r8*(v(istr-1,jend ,k,nout)+ &
675 & v(istr ,jend+1,k,nout))
679 IF (
domain(ng)%NorthEast_Corner(tile))
THEN
683 v(iend+1,jend+1,k,nout)=0.5_r8*(v(iend+1,jend ,k,nout)+ &
684 & v(iend ,jend+1,k,nout))
type(t_boundary), dimension(:), allocatable boundary
type(t_apply), dimension(:), allocatable lbc_apply
type(t_clima), dimension(:), allocatable clima
type(t_grid), dimension(:), allocatable grid
type(t_ocean), dimension(:), allocatable ocean
integer, dimension(:), allocatable n
type(t_lbc), dimension(:,:,:), allocatable lbc
type(t_domain), dimension(:), allocatable domain
real(dp), dimension(:,:), allocatable m3obc_out
real(dp), dimension(:), allocatable dt
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
real(dp), dimension(:), allocatable obcfac
real(r8), dimension(:), allocatable gamma2
logical, dimension(:), allocatable lnudgem3clm
integer, parameter isouth
integer, parameter inorth
real(dp), dimension(:,:), allocatable m3obc_in
integer, dimension(:), allocatable nstp
subroutine v3dbc(ng, tile, nout)
subroutine, public v3dbc_tile(ng, tile, lbi, ubi, lbj, ubj, ubk, imins, imaxs, jmins, jmaxs, nstp, nout, v)