MODULE u2dbc_mod ! !svn $Id: u2dbc_im.F 294 2009-01-09 21:37:26Z 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 U-velocity. ! ! ! !======================================================================= ! implicit none PRIVATE PUBLIC :: u2dbc, u2dbc_tile CONTAINS ! !*********************************************************************** SUBROUTINE u2dbc (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 u2dbc_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 u2dbc ! !*********************************************************************** SUBROUTINE u2dbc_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) :: vbar(LBi:,LBj:,:) real(r8), intent(in) :: zeta(LBi:,LBj:,:) real(r8), intent(inout) :: ubar(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, dUde, dUdt, dUdx, 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 western edge. !----------------------------------------------------------------------- ! IF (Istr.eq.1) THEN ! ! Western edge, Flather boundary condition. ! DO j=Jstr,Jend bry_val=BOUNDARY(ng)%ubar_west(j) cff=1.0_r8/(0.5_r8*(GRID(ng)%h(Istr-1,j)+ & & zeta(Istr-1,j,know)+ & & GRID(ng)%h(Istr ,j)+ & & zeta(Istr ,j,know))) Cx=SQRT(g*cff) ubar(Istr,j,kout)=bry_val- & & Cx*(0.5_r8*(zeta(Istr-1,j,know)+ & & zeta(Istr ,j,know))- & & BOUNDARY(ng)%zeta_west(j)) ubar(Istr,j,kout)=ubar(Istr,j,kout)* & & GRID(ng)%umask(Istr,j) END DO END IF ! !----------------------------------------------------------------------- ! Lateral boundary conditions at the eastern edge. !----------------------------------------------------------------------- ! IF (Iend.eq.Lm(ng)) THEN ! ! Eastern edge, Flather boundary condition. ! DO j=Jstr,Jend bry_val=BOUNDARY(ng)%ubar_east(j) cff=1.0_r8/(0.5_r8*(GRID(ng)%h(Iend ,j)+ & & zeta(Iend ,j,know)+ & & GRID(ng)%h(Iend+1,j)+ & & zeta(Iend+1,j,know))) Cx=SQRT(g*cff) ubar(Iend+1,j,kout)=bry_val+ & & Cx*(0.5_r8*(zeta(Iend ,j,know)+ & & zeta(Iend+1,j,know))- & & BOUNDARY(ng)%zeta_east(j)) ubar(Iend+1,j,kout)=ubar(Iend+1,j,kout)* & & GRID(ng)%umask(Iend+1,j) END DO END IF ! !----------------------------------------------------------------------- ! Lateral boundary conditions at the southern edge. !----------------------------------------------------------------------- ! IF (Jstr.eq.1) THEN ! ! Southern edge, Chapman boundary condition. ! DO i=IstrU,Iend cff=dt2d*0.5_r8*(GRID(ng)%pn(i-1,Jstr)+ & & GRID(ng)%pn(i ,Jstr)) cff1=SQRT(g*0.5_r8*(GRID(ng)%h(i-1,Jstr)+ & & zeta(i-1,Jstr,know)+ & & GRID(ng)%h(i ,Jstr)+ & & zeta(i ,Jstr,know))) Ce=cff*cff1 cff2=1.0_r8/(1.0_r8+Ce) ubar(i,Jstr-1,kout)=cff2*(ubar(i,Jstr-1,know)+ & & Ce*ubar(i,Jstr,kout)) ubar(i,Jstr-1,kout)=ubar(i,Jstr-1,kout)* & & GRID(ng)%umask(i,Jstr-1) END DO END IF ! !----------------------------------------------------------------------- ! Lateral boundary conditions at the northern edge. !----------------------------------------------------------------------- ! IF (Jend.eq.Mm(ng)) THEN ! ! Northern edge, closed boundary condition: free slip (gamma2=1) or ! no slip (gamma2=-1). ! DO i=Istr,IendR ubar(i,Jend+1,kout)=gamma2(ng)*ubar(i,Jend,kout) ubar(i,Jend+1,kout)=ubar(i,Jend+1,kout)* & & GRID(ng)%umask(i,Jend+1) END DO END IF ! !----------------------------------------------------------------------- ! Boundary corners. !----------------------------------------------------------------------- ! IF ((Jstr.eq.1).and.(Istr.eq.1)) THEN ubar(Istr,Jstr-1,kout)=0.5_r8*(ubar(Istr+1,Jstr-1,kout)+ & & ubar(Istr ,Jstr ,kout)) END IF IF ((Jstr.eq.1).and.(Iend.eq.Lm(ng))) THEN ubar(Iend+1,Jstr-1,kout)=0.5_r8*(ubar(Iend ,Jstr-1,kout)+ & & ubar(Iend+1,Jstr ,kout)) END IF IF ((Jend.eq.Mm(ng)).and.(Istr.eq.1)) THEN ubar(Istr,Jend+1,kout)=0.5_r8*(ubar(Istr ,Jend ,kout)+ & & ubar(Istr+1,Jend+1,kout)) END IF IF ((Jend.eq.Mm(ng)).and.(Iend.eq.Lm(ng))) THEN ubar(Iend+1,Jend+1,kout)=0.5_r8*(ubar(Iend+1,Jend ,kout)+ & & ubar(Iend ,Jend+1,kout)) END IF ! !----------------------------------------------------------------------- ! Impose wetting and drying conditions. !----------------------------------------------------------------------- ! IF (Istr.eq.1) THEN DO j=Jstr,Jend cff1=ABS(ABS(GRID(ng)%umask_wet(Istr,j))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,ubar(Istr,j,kout))* & & GRID(ng)%umask_wet(Istr,j) cff=0.5_r8*GRID(ng)%umask_wet(Istr,j)*cff1+ & & cff2*(1.0_r8-cff1) ubar(Istr,j,kout)=ubar(Istr,j,kout)*cff END DO END IF IF (Iend.eq.Lm(ng)) THEN DO j=Jstr,Jend cff1=ABS(ABS(GRID(ng)%umask_wet(Iend+1,j))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,ubar(Iend+1,j,kout))* & & GRID(ng)%umask_wet(Iend+1,j) cff=0.5_r8*GRID(ng)%umask_wet(Iend+1,j)*cff1+ & & cff2*(1.0_r8-cff1) ubar(Iend+1,j,kout)=ubar(Iend+1,j,kout)*cff END DO END IF IF (Jstr.eq.1) THEN DO i=IstrU,Iend cff1=ABS(ABS(GRID(ng)%umask_wet(i,Jstr-1))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,ubar(i,Jstr-1,kout))* & & GRID(ng)%umask_wet(i,Jstr-1) cff=0.5_r8*GRID(ng)%umask_wet(i,Jstr-1)*cff1+ & & cff2*(1.0_r8-cff1) ubar(i,Jstr-1,kout)=ubar(i,Jstr-1,kout)*cff END DO END IF IF (Jend.eq.Mm(ng)) THEN DO i=Istr,Iend cff1=ABS(ABS(GRID(ng)%umask_wet(i,Jend+1))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,ubar(i,Jend+1,kout))* & & GRID(ng)%umask_wet(i,Jend+1) cff=0.5_r8*GRID(ng)%umask_wet(i,Jend+1)*cff1+ & & cff2*(1.0_r8-cff1) ubar(i,Jend+1,kout)=ubar(i,Jend+1,kout)*cff END DO END IF IF ((Jstr.eq.1).and.(Istr.eq.1)) THEN cff1=ABS(ABS(GRID(ng)%umask_wet(Istr,Jstr-1))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,ubar(Istr,Jstr-1,kout))* & & GRID(ng)%umask_wet(Istr,Jstr-1) cff=0.5_r8*GRID(ng)%umask_wet(Istr,Jstr-1)*cff1+ & & cff2*(1.0_r8-cff1) ubar(Istr,Jstr-1,kout)=ubar(Istr,Jstr-1,kout)*cff END IF IF ((Jstr.eq.1).and.(Iend.eq.Lm(ng))) THEN cff1=ABS(ABS(GRID(ng)%umask_wet(Iend+1,Jstr-1))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,ubar(Iend+1,Jstr-1,kout))* & & GRID(ng)%umask_wet(Iend+1,Jstr-1) cff=0.5_r8*GRID(ng)%umask_wet(Iend+1,Jstr-1)*cff1+ & & cff2*(1.0_r8-cff1) ubar(Iend+1,Jstr-1,kout)=ubar(Iend+1,Jstr-1,kout)*cff END IF IF ((Jend.eq.Mm(ng)).and.(Istr.eq.1)) THEN cff1=ABS(ABS(GRID(ng)%umask_wet(Istr,Jend+1))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,ubar(Istr,Jend+1,kout))* & & GRID(ng)%umask_wet(Istr,Jend+1) cff=0.5_r8*GRID(ng)%umask_wet(Istr,Jend+1)*cff1+ & & cff2*(1.0_r8-cff1) ubar(Istr,Jend+1,kout)=ubar(Istr,Jend+1,kout)*cff END IF IF ((Jend.eq.Mm(ng)).and.(Iend.eq.Lm(ng))) THEN cff1=ABS(ABS(GRID(ng)%umask_wet(Iend+1,Jend+1))-1.0_r8) cff2=0.5_r8+DSIGN(0.5_r8,ubar(Iend+1,Jend+1,kout))* & & GRID(ng)%umask_wet(Iend+1,Jend+1) cff=0.5_r8*GRID(ng)%umask_wet(Iend+1,Jend+1)+cff1+ & & cff2*(1.0_r8-cff1) ubar(Iend+1,Jend+1,kout)=ubar(Iend+1,Jend+1,kout)*cff END IF RETURN END SUBROUTINE u2dbc_tile END MODULE u2dbc_mod