MODULE v2dbc_mod ! !svn $Id: v2dbc_im.F 349 2009-04-17 19:56:13Z arango $ !======================================================================= ! Copyright (c) 2002-2009 The ROMS/TOMS Group ! ! Licensed under a MIT/X style license ! ! See License_ROMS.txt Hernan G. Arango ! !========================================== Alexander F. Shchepetkin === ! ! ! This subroutine sets lateral boundary conditions for vertically ! ! integrated V-velocity. ! ! ! !======================================================================= ! implicit none PRIVATE PUBLIC :: v2dbc, v2dbc_tile CONTAINS ! !*********************************************************************** SUBROUTINE v2dbc (ng, tile, kout) !*********************************************************************** ! USE mod_param USE mod_ocean USE mod_stepping ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, kout ! ! Local variable declarations. ! integer :: IminS, ImaxS, JminS, JmaxS integer :: LBi, UBi, LBj, UBj, LBij, UBij ! ! Set horizontal starting and ending indices for automatic private storage ! arrays. ! IminS=BOUNDS(ng)%Istr(tile)-3 ImaxS=BOUNDS(ng)%Iend(tile)+3 JminS=BOUNDS(ng)%Jstr(tile)-3 JmaxS=BOUNDS(ng)%Jend(tile)+3 ! ! Determine array lower and upper bounds in the I- and J-directions. ! LBi=BOUNDS(ng)%LBi(tile) UBi=BOUNDS(ng)%UBi(tile) LBj=BOUNDS(ng)%LBj(tile) UBj=BOUNDS(ng)%UBj(tile) ! ! Set array lower and upper bounds for MIN(I,J)- and MAX(I,J)-directions. ! LBij=BOUNDS(ng)%LBij UBij=BOUNDS(ng)%UBij ! CALL v2dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs(ng), kstp(ng), kout, & & OCEAN(ng) % ubar, & & OCEAN(ng) % vbar, & & OCEAN(ng) % zeta) RETURN END SUBROUTINE v2dbc ! !*********************************************************************** SUBROUTINE v2dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, kout, & & ubar, vbar, zeta) !*********************************************************************** ! USE mod_param USE mod_boundary USE mod_forces USE mod_grid USE mod_scalars ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile integer, intent(in) :: LBi, UBi, LBj, UBj integer, intent(in) :: IminS, ImaxS, JminS, JmaxS integer, intent(in) :: krhs, kstp, kout real(r8), intent(in) :: ubar(LBi:,LBj:,:) real(r8), intent(in) :: zeta(LBi:,LBj:,:) real(r8), intent(inout) :: vbar(LBi:,LBj:,:) ! ! Local variable declarations. ! integer :: i, j, know real(r8), parameter :: eps = 1.0E-20_r8 real(r8) :: Ce, Cx real(r8) :: bry_pgr, bry_cor, bry_str, bry_val real(r8):: cff, cff1, cff2, dt2d, dVde, dVdt, dVdx, tau real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: grad ! !----------------------------------------------------------------------- ! Set lower and upper tile bounds and staggered variables bounds for ! this horizontal domain partition. Notice that if tile=-1, it will ! set the values for the global grid. !----------------------------------------------------------------------- ! integer :: Istr, IstrR, IstrT, IstrU, Iend, IendR, IendT integer :: Jstr, JstrR, JstrT, JstrV, Jend, JendR, JendT ! Istr =BOUNDS(ng)%Istr (tile) IstrR=BOUNDS(ng)%IstrR(tile) IstrT=BOUNDS(ng)%IstrT(tile) IstrU=BOUNDS(ng)%IstrU(tile) Iend =BOUNDS(ng)%Iend (tile) IendR=BOUNDS(ng)%IendR(tile) IendT=BOUNDS(ng)%IendT(tile) Jstr =BOUNDS(ng)%Jstr (tile) JstrR=BOUNDS(ng)%JstrR(tile) JstrT=BOUNDS(ng)%JstrT(tile) JstrV=BOUNDS(ng)%JstrV(tile) Jend =BOUNDS(ng)%Jend (tile) JendR=BOUNDS(ng)%JendR(tile) JendT=BOUNDS(ng)%JendT(tile) ! !----------------------------------------------------------------------- ! Set time-indices !----------------------------------------------------------------------- ! IF (iif(ng).eq.1) THEN know=krhs dt2d=dtfast(ng) ELSE IF (PREDICTOR_2D_STEP(ng)) THEN know=krhs dt2d=2.0_r8*dtfast(ng) ELSE know=kstp dt2d=dtfast(ng) END IF ! !----------------------------------------------------------------------- ! Lateral boundary conditions at the southern edge. !----------------------------------------------------------------------- ! IF (Jstr.eq.1) THEN ! ! Southern edge, Flather boundary condition. ! DO i=Istr,Iend bry_val=BOUNDARY(ng)%vbar_south(i) cff=1.0_r8/(0.5_r8*(GRID(ng)%h(i,Jstr-1)+ & & zeta(i,Jstr-1,know)+ & & GRID(ng)%h(i,Jstr )+ & & zeta(i,Jstr ,know))) Ce=SQRT(g*cff) vbar(i,Jstr,kout)=bry_val- & & Ce*(0.5_r8*(zeta(i,Jstr-1,know)+ & & zeta(i,Jstr ,know))- & & BOUNDARY(ng)%zeta_south(i)) vbar(i,Jstr,kout)=vbar(i,Jstr,kout)* & & GRID(ng)%vmask(i,Jstr) END DO END IF ! !----------------------------------------------------------------------- ! Lateral boundary conditions at the northern edge. !----------------------------------------------------------------------- ! IF (Jend.eq.Mm(ng)) THEN ! ! Northern edge, closed boundary condition. ! DO i=Istr,Iend vbar(i,Jend+1,kout)=0.0_r8 END DO END IF ! !----------------------------------------------------------------------- ! Lateral boundary conditions at the western edge. !----------------------------------------------------------------------- ! IF (Istr.eq.1) THEN ! ! Western edge, Chapman boundary condition. ! DO j=JstrV,Jend cff=dt2d*0.5_r8*(GRID(ng)%pm(Istr,j-1)+ & & GRID(ng)%pm(Istr,j )) cff1=SQRT(g*0.5_r8*(GRID(ng)%h(Istr,j-1)+ & & zeta(Istr,j-1,know)+ & & GRID(ng)%h(Istr,j )+ & & zeta(Istr,j ,know))) Cx=cff*cff1 cff2=1.0_r8/(1.0_r8+Cx) vbar(Istr-1,j,kout)=cff2*(vbar(Istr-1,j,know)+ & & Cx*vbar(Istr,j,kout)) vbar(Istr-1,j,kout)=vbar(Istr-1,j,kout)* & & GRID(ng)%vmask(Istr-1,j) END DO END IF ! !----------------------------------------------------------------------- ! Lateral boundary conditions at the eastern edge. !----------------------------------------------------------------------- ! IF (Iend.eq.Lm(ng)) THEN ! ! Eastern edge, Chapman boundary condition. ! DO j=JstrV,Jend cff=dt2d*0.5_r8*(GRID(ng)%pm(Iend,j-1)+ & & GRID(ng)%pm(Iend,j )) cff1=SQRT(g*0.5_r8*(GRID(ng)%h(Iend,j-1)+ & & zeta(Iend,j-1,know)+ & & GRID(ng)%h(Iend,j )+ & & zeta(Iend,j ,know))) Cx=cff*cff1 cff2=1.0_r8/(1.0_r8+Cx) vbar(Iend+1,j,kout)=cff2*(vbar(Iend+1,j,know)+ & & Cx*vbar(Iend,j,kout)) vbar(Iend+1,j,kout)=vbar(Iend+1,j,kout)* & & GRID(ng)%vmask(Iend+1,j) END DO END IF ! !----------------------------------------------------------------------- ! Boundary corners. !----------------------------------------------------------------------- ! IF ((Jstr.eq.1).and.(Istr.eq.1)) THEN vbar(Istr-1,Jstr,kout)=0.5_r8*(vbar(Istr ,Jstr ,kout)+ & & vbar(Istr-1,Jstr+1,kout)) END IF IF ((Jstr.eq.1).and.(Iend.eq.Lm(ng))) THEN vbar(Iend+1,Jstr,kout)=0.5_r8*(vbar(Iend ,Jstr ,kout)+ & & vbar(Iend+1,Jstr+1,kout)) END IF IF ((Jend.eq.Mm(ng)).and.(Istr.eq.1)) THEN vbar(Istr-1,Jend+1,kout)=0.5_r8*(vbar(Istr-1,Jend ,kout)+ & & vbar(Istr ,Jend+1,kout)) END IF IF ((Jend.eq.Mm(ng)).and.(Iend.eq.Lm(ng))) THEN vbar(Iend+1,Jend+1,kout)=0.5_r8*(vbar(Iend+1,Jend ,kout)+ & & vbar(Iend ,Jend+1,kout)) END IF ! !----------------------------------------------------------------------- ! Impose wetting and drying conditions. !----------------------------------------------------------------------- ! IF (Istr.eq.1) THEN DO j=JstrV,Jend cff1=ABS(ABS(GRID(ng)%vmask_wet(Istr-1,j))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,vbar(Istr-1,j,kout))* & & GRID(ng)%vmask_wet(Istr-1,j) cff=0.5_r8*GRID(ng)%vmask_wet(Istr-1,j)*cff1+ & & cff2*(1.0_r8-cff1) vbar(Istr,j,kout)=vbar(Istr,j,kout)*cff END DO END IF IF (Iend.eq.Lm(ng)) THEN DO j=JstrV,Jend cff1=ABS(ABS(GRID(ng)%vmask_wet(Iend+1,j))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,vbar(Iend+1,j,kout))* & & GRID(ng)%vmask_wet(Iend+1,j) cff=0.5_r8*GRID(ng)%vmask_wet(Iend+1,j)*cff1+ & & cff2*(1.0_r8-cff1) vbar(Iend+1,j,kout)=vbar(Iend+1,j,kout)*cff END DO END IF IF (Jstr.eq.1) THEN DO i=Istr,Iend cff1=ABS(ABS(GRID(ng)%vmask_wet(i,Jstr))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,vbar(i,Jstr,kout))* & & GRID(ng)%vmask_wet(i,Jstr) cff=0.5_r8*GRID(ng)%vmask_wet(i,Jstr)*cff1+ & & cff2*(1.0_r8-cff1) vbar(i,Jstr,kout)=vbar(i,Jstr,kout)*cff END DO END IF IF (Jend.eq.Mm(ng)) THEN DO i=Istr,Iend cff1=ABS(ABS(GRID(ng)%vmask_wet(i,Jend+1))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,vbar(i,Jend+1,kout))* & & GRID(ng)%vmask_wet(i,Jend+1) cff=0.5_r8*GRID(ng)%vmask_wet(i,Jend+1)*cff1+ & & cff2*(1.0_r8-cff1) vbar(i,Jend+1,kout)=vbar(i,Jend+1,kout)*cff END DO END IF IF ((Jstr.eq.1).and.(Istr.eq.1)) THEN cff1=ABS(ABS(GRID(ng)%vmask_wet(Istr-1,Jstr))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,vbar(Istr-1,Jstr,kout))* & & GRID(ng)%vmask_wet(Istr-1,Jstr) cff=0.5_r8*GRID(ng)%vmask_wet(Istr-1,Jstr)*cff1+ & & cff2*(1.0_r8-cff1) vbar(Istr-1,Jstr,kout)=vbar(Istr-1,Jstr,kout)*cff END IF IF ((Jstr.eq.1).and.(Iend.eq.Lm(ng))) THEN cff1=ABS(ABS(GRID(ng)%vmask_wet(Iend+1,Jstr))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,vbar(Iend+1,Jstr,kout))* & & GRID(ng)%vmask_wet(Iend+1,Jstr) cff=0.5_r8*GRID(ng)%vmask_wet(Iend+1,Jstr)*cff1+ & & cff2*(1.0_r8-cff1) vbar(Iend+1,Jstr,kout)=vbar(Iend+1,Jstr,kout)*cff END IF IF ((Jend.eq.Mm(ng)).and.(Istr.eq.1)) THEN cff1=ABS(ABS(GRID(ng)%vmask_wet(Istr-1,Jend+1))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,vbar(Istr-1,Jend+1,kout))* & & GRID(ng)%vmask_wet(Istr-1,Jend+1) cff=0.5_r8*GRID(ng)%vmask_wet(Istr-1,Jend+1)*cff1+ & & cff2*(1.0_r8-cff1) vbar(Istr-1,Jend+1,kout)=vbar(Istr-1,Jend+1,kout)*cff END IF IF ((Jend.eq.Mm(ng)).and.(Iend.eq.Lm(ng))) THEN cff1=ABS(ABS(GRID(ng)%vmask_wet(Iend+1,Jend+1))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,vbar(Iend+1,Jend+1,kout))* & & GRID(ng)%vmask_wet(Iend+1,Jend+1) cff=0.5_r8*GRID(ng)%vmask_wet(Iend+1,Jend+1)*cff1+ & & cff2*(1.0_r8-cff1) vbar(Iend+1,Jend+1,kout)=vbar(Iend+1,Jend+1,kout)*cff END IF RETURN END SUBROUTINE v2dbc_tile END MODULE v2dbc_mod