MODULE step2d_mod ! !svn $Id: step2d.F 294 2009-01-09 21:37:26Z arango $ !======================================================================= ! Copyright (c) 2002-2009 The ROMS/TOMS Group ! ! Licensed under a MIT/X style license Hernan G. Arango ! ! See License_ROMS.txt Alexander F. Shchepetkin ! !==================================================== John C. Warner === ! ! ! This subroutine performs a fast (predictor or corrector) time-step ! ! for the free-surface and 2D momentum nonlinear equations. ! It also calculates the time filtering variables over all fast-time ! ! steps to damp high frequency signals in 3D applications. ! ! ! !======================================================================= ! implicit none ! PRIVATE PUBLIC :: step2d CONTAINS SUBROUTINE step2d (ng, tile) ! !svn $Id: step2d_LF_AM3.h 381 2009-08-11 19:50:39Z arango $ !======================================================================= ! ! ! Nonlinear shallow-water primitive equations predictor (Leap-frog) ! ! and corrector (Adams-Moulton) time-stepping engine. ! ! ! !======================================================================= ! USE mod_param USE mod_coupling USE mod_forces USE mod_grid USE mod_ocean USE mod_sources USE mod_stepping ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile ! ! 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 wclock_on (ng, iNLM, 9) CALL step2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, N(ng), & & IminS, ImaxS, JminS, JmaxS, & & krhs(ng), kstp(ng), knew(ng), & & nstp(ng), nnew(ng), & & Msrc(ng), Nsrc(ng), & & SOURCES(ng) % Isrc, SOURCES(ng) % Jsrc, & & SOURCES(ng) % Dsrc, SOURCES(ng) % Qbar, & & GRID(ng) % pmask, GRID(ng) % rmask, & & GRID(ng) % umask, GRID(ng) % vmask, & & GRID(ng) % rmask_wet, GRID(ng) % rmask_full, & & GRID(ng) % umask_wet, GRID(ng) % umask_full, & & GRID(ng) % vmask_wet, GRID(ng) % vmask_full, & & GRID(ng) % rmask_wet_avg, & & GRID(ng) % fomn, GRID(ng) % h, & & GRID(ng) % om_u, GRID(ng) % om_v, & & GRID(ng) % on_u, GRID(ng) % on_v, & & GRID(ng) % omn, & & GRID(ng) % pm, GRID(ng) % pn, & & GRID(ng) % dndx, GRID(ng) % dmde, & & COUPLING(ng) % rhoA, COUPLING(ng) % rhoS, & & COUPLING(ng) % DU_avg1, COUPLING(ng) % DU_avg2, & & COUPLING(ng) % DV_avg1, COUPLING(ng) % DV_avg2, & & COUPLING(ng) % Zt_avg1, & & COUPLING(ng) % rufrc, COUPLING(ng) % rvfrc, & & OCEAN(ng) % ru, OCEAN(ng) % rv, & & OCEAN(ng) % rubar, OCEAN(ng) % rvbar, & & OCEAN(ng) % rzeta, & & OCEAN(ng) % ubar, OCEAN(ng) % vbar, & & OCEAN(ng) % zeta) CALL wclock_off (ng, iNLM, 9) RETURN END SUBROUTINE step2d ! !*********************************************************************** SUBROUTINE step2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, UBk, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, knew, & & nstp, nnew, & & Msrc, Nsrc, Isrc, Jsrc, Dsrc, Qbar, & & pmask, rmask, umask, vmask, & & rmask_wet, rmask_full, & & umask_wet, umask_full, & & vmask_wet, vmask_full, & & rmask_wet_avg, & & fomn, h, & & om_u, om_v, on_u, on_v, omn, pm, pn, & & dndx, dmde, & & rhoA, rhoS, & & DU_avg1, DU_avg2, & & DV_avg1, DV_avg2, & & Zt_avg1, & & rufrc, rvfrc, ru, rv, & & rubar, rvbar, rzeta, & & ubar, vbar, zeta) !*********************************************************************** ! USE mod_param USE mod_scalars ! USE mp_exchange_mod, ONLY : mp_exchange2d USE wetdry_mod, ONLY : wetdry_tile USE u2dbc_mod, ONLY : u2dbc_tile USE v2dbc_mod, ONLY : v2dbc_tile USE zetabc_mod, ONLY : zetabc_tile ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile integer, intent(in) :: LBi, UBi, LBj, UBj, UBk integer, intent(in) :: IminS, ImaxS, JminS, JmaxS integer, intent(in) :: krhs, kstp, knew integer, intent(in) :: nstp, nnew integer, intent(in) :: Msrc, Nsrc ! integer, intent(in) :: Isrc(:) integer, intent(in) :: Jsrc(:) real(r8), intent(in) :: Dsrc(:) real(r8), intent(in) :: Qbar(:) real(r8), intent(in) :: pmask(LBi:,LBj:) real(r8), intent(in) :: rmask(LBi:,LBj:) real(r8), intent(in) :: umask(LBi:,LBj:) real(r8), intent(in) :: vmask(LBi:,LBj:) real(r8), intent(in) :: fomn(LBi:,LBj:) real(r8), intent(in) :: om_u(LBi:,LBj:) real(r8), intent(in) :: om_v(LBi:,LBj:) real(r8), intent(in) :: on_u(LBi:,LBj:) real(r8), intent(in) :: on_v(LBi:,LBj:) real(r8), intent(in) :: omn(LBi:,LBj:) real(r8), intent(in) :: pm(LBi:,LBj:) real(r8), intent(in) :: pn(LBi:,LBj:) real(r8), intent(in) :: dndx(LBi:,LBj:) real(r8), intent(in) :: dmde(LBi:,LBj:) real(r8), intent(in) :: rhoA(LBi:,LBj:) real(r8), intent(in) :: rhoS(LBi:,LBj:) real(r8), intent(inout) :: DU_avg1(LBi:,LBj:) real(r8), intent(inout) :: DU_avg2(LBi:,LBj:) real(r8), intent(inout) :: DV_avg1(LBi:,LBj:) real(r8), intent(inout) :: DV_avg2(LBi:,LBj:) real(r8), intent(inout) :: Zt_avg1(LBi:,LBj:) real(r8), intent(inout) :: rufrc(LBi:,LBj:) real(r8), intent(inout) :: rvfrc(LBi:,LBj:) real(r8), intent(inout) :: ru(LBi:,LBj:,0:,:) real(r8), intent(inout) :: rv(LBi:,LBj:,0:,:) real(r8), intent(inout) :: rmask_full(LBi:,LBj:) real(r8), intent(inout) :: rmask_wet(LBi:,LBj:) real(r8), intent(inout) :: umask_full(LBi:,LBj:) real(r8), intent(inout) :: umask_wet(LBi:,LBj:) real(r8), intent(inout) :: vmask_full(LBi:,LBj:) real(r8), intent(inout) :: vmask_wet(LBi:,LBj:) real(r8), intent(inout) :: rmask_wet_avg(LBi:,LBj:) real(r8), intent(inout) :: h(LBi:,LBj:) real(r8), intent(inout) :: rubar(LBi:,LBj:,:) real(r8), intent(inout) :: rvbar(LBi:,LBj:,:) real(r8), intent(inout) :: rzeta(LBi:,LBj:,:) real(r8), intent(inout) :: ubar(LBi:,LBj:,:) real(r8), intent(inout) :: vbar(LBi:,LBj:,:) real(r8), intent(inout) :: zeta(LBi:,LBj:,:) ! ! Local variable declarations. ! logical :: CORRECTOR_2D_STEP logical :: EWperiodic=.FALSE. logical :: NSperiodic=.FALSE. integer :: i, j, ptsk integer :: is real(r8) :: cff, cff1, cff2, cff3, cff4, cff5, cff6, cff7 real(r8) :: fac, fac1, fac2, fac3 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dgrad real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dnew real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Drhs real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Drhs_p real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dstp real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: DUon real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: DVom real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: UFe real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: UFx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: VFe real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: VFx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: grad real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: gzeta real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: gzeta2 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: gzetaSA real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rhs_ubar real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rhs_vbar real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rhs_zeta real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: zeta_new real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: zwrk real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: wetdry ! !----------------------------------------------------------------------- ! 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) ! ptsk=3-kstp CORRECTOR_2D_STEP=.not.PREDICTOR_2D_STEP(ng) ! !----------------------------------------------------------------------- ! Compute total depth (m) and vertically integrated mass fluxes. !----------------------------------------------------------------------- ! DO j=-2+JstrV,MIN(Jend+1,Mm(ng))+1 DO i=-2+IstrU,MIN(Iend+1,Lm(ng))+1 Drhs(i,j)=zeta(i,j,krhs)+h(i,j) END DO END DO DO j=-2+JstrV,MIN(Jend+1,Mm(ng))+1 DO i=-1+IstrU,MIN(Iend+1,Lm(ng))+1 cff=0.5_r8*on_u(i,j) cff1=cff*(Drhs(i,j)+Drhs(i-1,j)) DUon(i,j)=ubar(i,j,krhs)*cff1 END DO END DO DO j=-1+JstrV,MIN(Jend+1,Mm(ng))+1 DO i=-2+IstrU,MIN(Iend+1,Lm(ng))+1 cff=0.5_r8*om_v(i,j) cff1=cff*(Drhs(i,j)+Drhs(i,j-1)) DVom(i,j)=vbar(i,j,krhs)*cff1 END DO END DO ! ! Do a special exchange to avoid having three ghost points for ! high order numerical stencil. Notice that a private array is ! passed to the exchange routine. It will also apply periodic ! boundary conditions if no partitions in I- or J-directions. ! CALL mp_exchange2d (ng, tile, iNLM, 2, & & IminS, ImaxS, JminS, JmaxS, & & NghostPoints, EWperiodic, NSperiodic, & & DUon, DVom) ! !----------------------------------------------------------------------- ! Compute time averaged fields over all short time-steps. !----------------------------------------------------------------------- ! IF (PREDICTOR_2D_STEP(ng)) THEN IF (iif(ng).eq.1) THEN ! ! Reset arrays for 2D fields averaged within the short time-steps. ! cff2=(-1.0_r8/12.0_r8)*weight(2,iif(ng)+1,ng) DO j=JstrR,JendR DO i=IstrR,IendR Zt_avg1(i,j)=0.0_r8 END DO DO i=Istr,IendR DU_avg1(i,j)=0.0_r8 DU_avg2(i,j)=cff2*DUon(i,j) END DO END DO DO j=Jstr,JendR DO i=IstrR,IendR DV_avg1(i,j)=0.0_r8 DV_avg2(i,j)=cff2*DVom(i,j) END DO END DO ELSE ! ! Accumulate field averages of previous time-step after they are ! computed in the previous corrector step, updated their boundaries, ! and synchronized. ! cff1=weight(1,iif(ng)-1,ng) cff2=(8.0_r8/12.0_r8)*weight(2,iif(ng) ,ng)- & & (1.0_r8/12.0_r8)*weight(2,iif(ng)+1,ng) DO j=JstrR,JendR DO i=IstrR,IendR Zt_avg1(i,j)=Zt_avg1(i,j)+cff1*zeta(i,j,krhs) END DO DO i=Istr,IendR DU_avg1(i,j)=DU_avg1(i,j)+cff1*DUon(i,j) DU_avg2(i,j)=DU_avg2(i,j)+cff2*DUon(i,j) END DO END DO DO j=Jstr,JendR DO i=IstrR,IendR DV_avg1(i,j)=DV_avg1(i,j)+cff1*DVom(i,j) DV_avg2(i,j)=DV_avg2(i,j)+cff2*DVom(i,j) END DO END DO END IF ELSE IF (iif(ng).eq.1) THEN cff2=weight(2,iif(ng),ng) ELSE cff2=(5.0_r8/12.0_r8)*weight(2,iif(ng),ng) END IF DO j=JstrR,JendR DO i=Istr,IendR DU_avg2(i,j)=DU_avg2(i,j)+cff2*DUon(i,j) END DO END DO DO j=Jstr,JendR DO i=IstrR,IendR DV_avg2(i,j)=DV_avg2(i,j)+cff2*DVom(i,j) END DO END DO END IF ! ! After all fast time steps are completed, recompute S-coordinate ! surfaces according to the new free surface field. Apply boundary ! conditions to time averaged fields. ! IF ((iif(ng).eq.(nfast(ng)+1)).and.PREDICTOR_2D_STEP(ng)) THEN CALL mp_exchange2d (ng, tile, iNLM, 3, & & LBi, UBi, LBj, UBj, & & NghostPoints, EWperiodic, NSperiodic, & & Zt_avg1, DU_avg1, DV_avg1) END IF ! !----------------------------------------------------------------------- ! Compute new wet/dry masks. !----------------------------------------------------------------------- ! CALL wetdry_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & rmask, umask, vmask, & & h, zeta(:,:,kstp), & & DU_avg1, DV_avg1, & & rmask_wet_avg, & & rmask_full, umask_full, vmask_full, & & rmask_wet, umask_wet, vmask_wet) ! ! Do not perform the actual time stepping during the auxiliary ! (nfast(ng)+1) time step. ! IF (iif(ng).gt.nfast(ng)) RETURN ! !======================================================================= ! Time step free-surface equation. !======================================================================= ! ! During the first time-step, the predictor step is Forward-Euler ! and the corrector step is Backward-Euler. Otherwise, the predictor ! step is Leap-frog and the corrector step is Adams-Moulton. ! fac=1000.0_r8/rho0 IF (iif(ng).eq.1) THEN cff1=dtfast(ng) DO j=JstrV-1,Jend DO i=IstrU-1,Iend rhs_zeta(i,j)=(DUon(i,j)-DUon(i+1,j))+ & & (DVom(i,j)-DVom(i,j+1)) zeta_new(i,j)=zeta(i,j,kstp)+ & & pm(i,j)*pn(i,j)*cff1*rhs_zeta(i,j) zeta_new(i,j)=zeta_new(i,j)*rmask(i,j) Dnew(i,j)=zeta_new(i,j)+h(i,j) ! zwrk(i,j)=0.5_r8*(zeta(i,j,kstp)+zeta_new(i,j)) gzeta(i,j)=(fac+rhoS(i,j))*zwrk(i,j) gzeta2(i,j)=gzeta(i,j)*zwrk(i,j) gzetaSA(i,j)=zwrk(i,j)*(rhoS(i,j)-rhoA(i,j)) END DO END DO ELSE IF (PREDICTOR_2D_STEP(ng)) THEN cff1=2.0_r8*dtfast(ng) cff4=4.0_r8/25.0_r8 cff5=1.0_r8-2.0_r8*cff4 DO j=JstrV-1,Jend DO i=IstrU-1,Iend rhs_zeta(i,j)=(DUon(i,j)-DUon(i+1,j))+ & & (DVom(i,j)-DVom(i,j+1)) zeta_new(i,j)=zeta(i,j,kstp)+ & & pm(i,j)*pn(i,j)*cff1*rhs_zeta(i,j) zeta_new(i,j)=zeta_new(i,j)*rmask(i,j) Dnew(i,j)=zeta_new(i,j)+h(i,j) ! zwrk(i,j)=cff5*zeta(i,j,krhs)+ & & cff4*(zeta(i,j,kstp)+zeta_new(i,j)) gzeta(i,j)=(fac+rhoS(i,j))*zwrk(i,j) gzeta2(i,j)=gzeta(i,j)*zwrk(i,j) gzetaSA(i,j)=zwrk(i,j)*(rhoS(i,j)-rhoA(i,j)) END DO END DO ELSE IF (CORRECTOR_2D_STEP) THEN cff1=dtfast(ng)*5.0_r8/12.0_r8 cff2=dtfast(ng)*8.0_r8/12.0_r8 cff3=dtfast(ng)*1.0_r8/12.0_r8 cff4=2.0_r8/5.0_r8 cff5=1.0_r8-cff4 DO j=JstrV-1,Jend DO i=IstrU-1,Iend cff=cff1*((DUon(i,j)-DUon(i+1,j))+ & & (DVom(i,j)-DVom(i,j+1))) zeta_new(i,j)=zeta(i,j,kstp)+ & & pm(i,j)*pn(i,j)*(cff+ & & cff2*rzeta(i,j,kstp)- & & cff3*rzeta(i,j,ptsk)) zeta_new(i,j)=zeta_new(i,j)*rmask(i,j) Dnew(i,j)=zeta_new(i,j)+h(i,j) ! zwrk(i,j)=cff5*zeta_new(i,j)+cff4*zeta(i,j,krhs) gzeta(i,j)=(fac+rhoS(i,j))*zwrk(i,j) gzeta2(i,j)=gzeta(i,j)*zwrk(i,j) gzetaSA(i,j)=zwrk(i,j)*(rhoS(i,j)-rhoA(i,j)) END DO END DO END IF ! ! Load new free-surface values into shared array at both predictor ! and corrector steps. ! Modify new free-surface to Ensure that depth is > Dcrit for masked ! cells. ! DO j=Jstr,Jend DO i=Istr,Iend zeta(i,j,knew)=zeta_new(i,j) zeta(i,j,knew)=zeta(i,j,knew)+ & & (Dcrit(ng)-h(i,j))*(1.0_r8-rmask(i,j)) END DO END DO ! ! If predictor step, load right-side-term into shared array. ! IF (PREDICTOR_2D_STEP(ng)) THEN DO j=Jstr,Jend DO i=Istr,Iend rzeta(i,j,krhs)=rhs_zeta(i,j) END DO END DO CALL mp_exchange2d (ng, tile, iNLM, 1, & & LBi, UBi, LBj, UBj, & & NghostPoints, EWperiodic, NSperiodic, & & rzeta(:,:,krhs)) END IF ! ! Set free-surface lateral boundary conditions. ! CALL zetabc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, knew, & & zeta) CALL mp_exchange2d (ng, tile, iNLM, 1, & & LBi, UBi, LBj, UBj, & & NghostPoints, EWperiodic, NSperiodic, & & zeta(:,:,knew)) ! !======================================================================= ! Compute right-hand-side for the 2D momentum equations. !======================================================================= ! !----------------------------------------------------------------------- ! Compute pressure gradient terms. !----------------------------------------------------------------------- ! cff1=0.5_r8*g cff2=1.0_r8/3.0_r8 DO j=Jstr,Jend DO i=IstrU,Iend rhs_ubar(i,j)=cff1*on_u(i,j)* & & ((h(i-1,j)+ & & h(i ,j))* & & (gzeta(i-1,j)- & & gzeta(i ,j))+ & & (h(i-1,j)- & & h(i ,j))* & & (gzetaSA(i-1,j)+ & & gzetaSA(i ,j)+ & & cff2*(rhoA(i-1,j)- & & rhoA(i ,j))* & & (zwrk(i-1,j)- & & zwrk(i ,j)))+ & & (gzeta2(i-1,j)- & & gzeta2(i ,j))) END DO IF (j.ge.JstrV) THEN DO i=Istr,Iend rhs_vbar(i,j)=cff1*om_v(i,j)* & & ((h(i,j-1)+ & & h(i,j ))* & & (gzeta(i,j-1)- & & gzeta(i,j ))+ & & (h(i,j-1)- & & h(i,j ))* & & (gzetaSA(i,j-1)+ & & gzetaSA(i,j )+ & & cff2*(rhoA(i,j-1)- & & rhoA(i,j ))* & & (zwrk(i,j-1)- & & zwrk(i,j )))+ & & (gzeta2(i,j-1)- & & gzeta2(i,j ))) END DO END IF END DO ! !----------------------------------------------------------------------- ! Add in horizontal advection of momentum. !----------------------------------------------------------------------- ! ! ! Fourth-order, centered differences advection. ! DO j=Jstr,Jend DO i=MAX(IstrU-1,2),MIN(Iend+1,Lm(ng)) grad (i,j)=ubar(i-1,j,krhs)-2.0_r8*ubar(i,j,krhs)+ & & ubar(i+1,j,krhs) Dgrad(i,j)=DUon(i-1,j)-2.0_r8*DUon(i,j)+DUon(i+1,j) END DO END DO IF (Istr.eq.1) THEN DO j=Jstr,Jend grad (Istr,j)=grad (Istr+1,j) Dgrad(Istr,j)=Dgrad(Istr+1,j) END DO END IF IF (Iend.eq.Lm(ng)) THEN DO j=Jstr,Jend grad (Iend+1,j)=grad (Iend,j) Dgrad(Iend+1,j)=Dgrad(Iend,j) END DO END IF cff=1.0_r8/6.0_r8 DO j=Jstr,Jend DO i=IstrU-1,Iend UFx(i,j)=0.25_r8*(ubar(i ,j,krhs)+ & & ubar(i+1,j,krhs)- & & cff*(grad (i,j)+grad (i+1,j)))* & & (DUon(i,j)+DUon(i+1,j)- & & cff*(Dgrad(i,j)+Dgrad(i+1,j))) END DO END DO DO j=MAX(Jstr-1,1),MIN(Jend+1,Mm(ng)) DO i=IstrU,Iend grad(i,j)=ubar(i,j-1,krhs)-2.0_r8*ubar(i,j,krhs)+ & & ubar(i,j+1,krhs) END DO END DO IF (Jstr.eq.1) THEN DO i=IstrU,Iend grad(i,Jstr-1)=grad(i,Jstr) END DO END IF IF (Jend.eq.Mm(ng)) THEN DO i=IstrU,Iend grad(i,Jend+1)=grad(i,Jend) END DO END IF DO j=Jstr,Jend+1 DO i=IstrU-1,Iend Dgrad(i,j)=DVom(i-1,j)-2.0_r8*DVom(i,j)+DVom(i+1,j) END DO END DO cff=1.0_r8/6.0_r8 DO j=Jstr,Jend+1 DO i=IstrU,Iend UFe(i,j)=0.25_r8*(ubar(i,j ,krhs)+ & & ubar(i,j-1,krhs)- & & cff*(grad (i,j)+grad (i,j-1)))* & & (DVom(i,j)+DVom(i-1,j)- & & cff*(Dgrad(i,j)+Dgrad(i-1,j))) END DO END DO DO j=JstrV,Jend DO i=MAX(Istr-1,1),MIN(Iend+1,Lm(ng)) grad(i,j)=vbar(i-1,j,krhs)-2.0_r8*vbar(i,j,krhs)+ & & vbar(i+1,j,krhs) END DO END DO IF (Istr.eq.1) THEN DO j=JstrV,Jend grad(Istr-1,j)=grad(Istr,j) END DO END IF IF (Iend.eq.Lm(ng)) THEN DO j=JstrV,Jend grad(Iend+1,j)=grad(Iend,j) END DO END IF DO j=JstrV-1,Jend DO i=Istr,Iend+1 Dgrad(i,j)=DUon(i,j-1)-2.0_r8*DUon(i,j)+DUon(i,j+1) END DO END DO cff=1.0_r8/6.0_r8 DO j=JstrV,Jend DO i=Istr,Iend+1 VFx(i,j)=0.25_r8*(vbar(i ,j,krhs)+ & & vbar(i-1,j,krhs)- & & cff*(grad (i,j)+grad (i-1,j)))* & & (DUon(i,j)+DUon(i,j-1)- & & cff*(Dgrad(i,j)+Dgrad(i,j-1))) END DO END DO DO j=MAX(JstrV-1,2),MIN(Jend+1,Mm(ng)) DO i=Istr,Iend grad(i,j)=vbar(i,j-1,krhs)-2.0_r8*vbar(i,j,krhs)+ & & vbar(i,j+1,krhs) Dgrad(i,j)=DVom(i,j-1)-2.0_r8*DVom(i,j)+DVom(i,j+1) END DO END DO IF (Jstr.eq.1) THEN DO i=Istr,Iend grad (i,Jstr)=grad (i,Jstr+1) Dgrad(i,Jstr)=Dgrad(i,Jstr+1) END DO END IF IF (Jend.eq.Mm(ng)) THEN DO i=Istr,Iend grad (i,Jend+1)=grad (i,Jend) Dgrad(i,Jend+1)=Dgrad(i,Jend) END DO END IF cff=1.0_r8/6.0_r8 DO j=JstrV-1,Jend DO i=Istr,Iend VFe(i,j)=0.25_r8*(vbar(i,j ,krhs)+ & & vbar(i,j+1,krhs)- & & cff*(grad (i,j)+grad (i,j+1)))* & & (DVom(i,j)+DVom(i,j+1)- & & cff*(Dgrad(i,j)+Dgrad(i,j+1))) END DO END DO DO j=Jstr,Jend DO i=IstrU,Iend fac=(UFx(i,j)-UFx(i-1,j))+ & & (UFe(i,j+1)-UFe(i,j)) rhs_ubar(i,j)=rhs_ubar(i,j)-fac END DO END DO DO j=JstrV,Jend DO i=Istr,Iend fac=(VFx(i+1,j)-VFx(i,j))+ & & (VFe(i,j)-VFe(i,j-1)) rhs_vbar(i,j)=rhs_vbar(i,j)-fac END DO END DO ! !----------------------------------------------------------------------- ! Add in Coriolis and curvilinear transformation terms, if any. !----------------------------------------------------------------------- ! DO j=JstrV-1,Jend DO i=IstrU-1,Iend cff=0.5_r8*Drhs(i,j)*( & & fomn(i,j) & & +0.5_r8*((vbar(i,j ,krhs)+ & & vbar(i,j+1,krhs))*dndx(i,j)- & & (ubar(i ,j,krhs)+ & & ubar(i+1,j,krhs))*dmde(i,j)) & & ) UFx(i,j)=cff*(vbar(i,j ,krhs)+ & & vbar(i,j+1,krhs)) VFe(i,j)=cff*(ubar(i ,j,krhs)+ & & ubar(i+1,j,krhs)) END DO END DO DO j=Jstr,Jend DO i=IstrU,Iend fac1=0.5_r8*(UFx(i,j)+UFx(i-1,j)) rhs_ubar(i,j)=rhs_ubar(i,j)+fac1 END DO END DO DO j=JstrV,Jend DO i=Istr,Iend fac1=0.5_r8*(VFe(i,j)+VFe(i,j-1)) rhs_vbar(i,j)=rhs_vbar(i,j)-fac1 END DO END DO ! !----------------------------------------------------------------------- ! Coupling between 2D and 3D equations. !----------------------------------------------------------------------- ! ! Before the predictor step of the first barotropic time-step, ! arrays "rufrc" and "rvfrc" contain the vertical integrals of ! the 3D right-hand-side terms for momentum equations (including ! surface and bottom stresses, if so prescribed). ! ! Convert them into forcing terms by subtracting the fast time ! "rhs_ubar" and "rhs_vbar" from them; Also, immediately apply ! these forcing terms "rhs_ubar" and "rhs_vbar". ! ! From now on, these newly computed forcing terms will remain ! constant during the fast time stepping and will added to ! "rhs_ubar" and "rhs_vbar" during all subsequent time steps. ! IF (iif(ng).eq.1.and.PREDICTOR_2D_STEP(ng)) THEN IF (iic(ng).eq.ntfirst(ng)) THEN DO j=Jstr,Jend DO i=IstrU,Iend rufrc(i,j)=rufrc(i,j)-rhs_ubar(i,j) rhs_ubar(i,j)=rhs_ubar(i,j)+rufrc(i,j) ru(i,j,0,nstp)=rufrc(i,j) END DO END DO DO j=JstrV,Jend DO i=Istr,Iend rvfrc(i,j)=rvfrc(i,j)-rhs_vbar(i,j) rhs_vbar(i,j)=rhs_vbar(i,j)+rvfrc(i,j) rv(i,j,0,nstp)=rvfrc(i,j) END DO END DO ELSE IF (iic(ng).eq.(ntfirst(ng)+1)) THEN DO j=Jstr,Jend DO i=IstrU,Iend rufrc(i,j)=rufrc(i,j)-rhs_ubar(i,j) rhs_ubar(i,j)=rhs_ubar(i,j)+ & & 1.5_r8*rufrc(i,j)-0.5_r8*ru(i,j,0,nnew) ru(i,j,0,nstp)=rufrc(i,j) END DO END DO DO j=JstrV,Jend DO i=Istr,Iend rvfrc(i,j)=rvfrc(i,j)-rhs_vbar(i,j) rhs_vbar(i,j)=rhs_vbar(i,j)+ & & 1.5_r8*rvfrc(i,j)-0.5_r8*rv(i,j,0,nnew) rv(i,j,0,nstp)=rvfrc(i,j) END DO END DO ELSE cff1=23.0_r8/12.0_r8 cff2=16.0_r8/12.0_r8 cff3= 5.0_r8/12.0_r8 DO j=Jstr,Jend DO i=IstrU,Iend rufrc(i,j)=rufrc(i,j)-rhs_ubar(i,j) rhs_ubar(i,j)=rhs_ubar(i,j)+ & & cff1*rufrc(i,j)- & & cff2*ru(i,j,0,nnew)+ & & cff3*ru(i,j,0,nstp) ru(i,j,0,nstp)=rufrc(i,j) END DO END DO DO j=JstrV,Jend DO i=Istr,Iend rvfrc(i,j)=rvfrc(i,j)-rhs_vbar(i,j) rhs_vbar(i,j)=rhs_vbar(i,j)+ & & cff1*rvfrc(i,j)- & & cff2*rv(i,j,0,nnew)+ & & cff3*rv(i,j,0,nstp) rv(i,j,0,nstp)=rvfrc(i,j) END DO END DO END IF ELSE DO j=Jstr,Jend DO i=IstrU,Iend rhs_ubar(i,j)=rhs_ubar(i,j)+rufrc(i,j) END DO END DO DO j=JstrV,Jend DO i=Istr,Iend rhs_vbar(i,j)=rhs_vbar(i,j)+rvfrc(i,j) END DO END DO END IF ! !======================================================================= ! Time step 2D momentum equations. !======================================================================= ! ! Compute total water column depth. ! DO j=JstrV-1,Jend DO i=IstrU-1,Iend Dstp(i,j)=zeta(i,j,kstp)+h(i,j) END DO END DO ! ! During the first time-step, the predictor step is Forward-Euler ! and the corrector step is Backward-Euler. Otherwise, the predictor ! step is Leap-frog and the corrector step is Adams-Moulton. ! IF (iif(ng).eq.1) THEN cff1=0.5_r8*dtfast(ng) cff2=1.0_r8/cff1 DO j=Jstr,Jend DO i=IstrU,Iend cff=(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j)) fac=1.0_r8/(Dnew(i,j)+Dnew(i-1,j)) ubar(i,j,knew)=(ubar(i,j,kstp)* & & (Dstp(i,j)+Dstp(i-1,j))+ & & cff*cff1*rhs_ubar(i,j))*fac ubar(i,j,knew)=ubar(i,j,knew)*umask(i,j) cff5=ABS(ABS(umask_wet(i,j))-1.0_r8) cff6=0.5_r8+DSIGN(0.5_r8,ubar(i,j,knew))*umask_wet(i,j) cff7=0.5_r8*umask_wet(i,j)*cff5+cff6*(1.0_r8-cff5) ubar(i,j,knew)=ubar(i,j,knew)*cff7 fac1=cff2/cff rhs_ubar(i,j)=(ubar(i,j,knew)*(Dnew(i,j)+Dnew(i-1,j))- & & ubar(i,j,kstp)*(Dstp(i,j)+Dstp(i-1,j)))* & & fac1 END DO END DO DO j=JstrV,Jend DO i=Istr,Iend cff=(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1)) fac=1.0_r8/(Dnew(i,j)+Dnew(i,j-1)) vbar(i,j,knew)=(vbar(i,j,kstp)* & & (Dstp(i,j)+Dstp(i,j-1))+ & & cff*cff1*rhs_vbar(i,j))*fac vbar(i,j,knew)=vbar(i,j,knew)*vmask(i,j) cff5=ABS(ABS(vmask_wet(i,j))-1.0_r8) cff6=0.5_r8+DSIGN(0.5_r8,vbar(i,j,knew))*vmask_wet(i,j) cff7=0.5_r8*vmask_wet(i,j)*cff5+cff6*(1.0_r8-cff5) vbar(i,j,knew)=vbar(i,j,knew)*cff7 fac1=cff2/cff rhs_vbar(i,j)=(vbar(i,j,knew)*(Dnew(i,j)+Dnew(i,j-1))- & & vbar(i,j,kstp)*(Dstp(i,j)+Dstp(i,j-1)))* & & fac1 END DO END DO ELSE IF (PREDICTOR_2D_STEP(ng)) THEN cff1=dtfast(ng) cff2=1.0_r8/cff1 DO j=Jstr,Jend DO i=IstrU,Iend cff=(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j)) fac=1.0_r8/(Dnew(i,j)+Dnew(i-1,j)) ubar(i,j,knew)=(ubar(i,j,kstp)* & & (Dstp(i,j)+Dstp(i-1,j))+ & & cff*cff1*rhs_ubar(i,j))*fac ubar(i,j,knew)=ubar(i,j,knew)*umask(i,j) cff5=ABS(ABS(umask_wet(i,j))-1.0_r8) cff6=0.5_r8+DSIGN(0.5_r8,ubar(i,j,knew))*umask_wet(i,j) cff7=0.5_r8*umask_wet(i,j)*cff5+cff6*(1.0_r8-cff5) ubar(i,j,knew)=ubar(i,j,knew)*cff7 fac1=cff2/cff rhs_ubar(i,j)=(ubar(i,j,knew)*(Dnew(i,j)+Dnew(i-1,j))- & & ubar(i,j,kstp)*(Dstp(i,j)+Dstp(i-1,j)))* & & fac1 END DO END DO DO j=JstrV,Jend DO i=Istr,Iend cff=(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1)) fac=1.0_r8/(Dnew(i,j)+Dnew(i,j-1)) vbar(i,j,knew)=(vbar(i,j,kstp)* & & (Dstp(i,j)+Dstp(i,j-1))+ & & cff*cff1*rhs_vbar(i,j))*fac vbar(i,j,knew)=vbar(i,j,knew)*vmask(i,j) cff5=ABS(ABS(vmask_wet(i,j))-1.0_r8) cff6=0.5_r8+DSIGN(0.5_r8,vbar(i,j,knew))*vmask_wet(i,j) cff7=0.5_r8*vmask_wet(i,j)*cff5+cff6*(1.0_r8-cff5) vbar(i,j,knew)=vbar(i,j,knew)*cff7 fac1=cff2/cff rhs_vbar(i,j)=(vbar(i,j,knew)*(Dnew(i,j)+Dnew(i,j-1))- & & vbar(i,j,kstp)*(Dstp(i,j)+Dstp(i,j-1)))* & & fac1 END DO END DO ELSE IF (CORRECTOR_2D_STEP) THEN cff1=0.5_r8*dtfast(ng)*5.0_r8/12.0_r8 cff2=0.5_r8*dtfast(ng)*8.0_r8/12.0_r8 cff3=0.5_r8*dtfast(ng)*1.0_r8/12.0_r8 cff4=1.0_r8/cff1 DO j=Jstr,Jend DO i=IstrU,Iend cff=(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j)) fac=1.0_r8/(Dnew(i,j)+Dnew(i-1,j)) ubar(i,j,knew)=(ubar(i,j,kstp)* & & (Dstp(i,j)+Dstp(i-1,j))+ & & cff*(cff1*rhs_ubar(i,j)+ & & cff2*rubar(i,j,kstp)- & & cff3*rubar(i,j,ptsk)))*fac ubar(i,j,knew)=ubar(i,j,knew)*umask(i,j) cff5=ABS(ABS(umask_wet(i,j))-1.0_r8) cff6=0.5_r8+DSIGN(0.5_r8,ubar(i,j,knew))*umask_wet(i,j) cff7=0.5_r8*umask_wet(i,j)*cff5+cff6*(1.0_r8-cff5) ubar(i,j,knew)=ubar(i,j,knew)*cff7 fac1=1.0_r8/cff rhs_ubar(i,j)=((ubar(i,j,knew)*(Dnew(i,j)+Dnew(i-1,j))- & & ubar(i,j,kstp)*(Dstp(i,j)+Dstp(i-1,j)))* & & fac1- & & cff2*rubar(i,j,kstp)+ & & cff3*rubar(i,j,ptsk))*cff4 END DO END DO DO j=JstrV,Jend DO i=Istr,Iend cff=(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1)) fac=1.0_r8/(Dnew(i,j)+Dnew(i,j-1)) vbar(i,j,knew)=(vbar(i,j,kstp)* & & (Dstp(i,j)+Dstp(i,j-1))+ & & cff*(cff1*rhs_vbar(i,j)+ & & cff2*rvbar(i,j,kstp)- & & cff3*rvbar(i,j,ptsk)))*fac vbar(i,j,knew)=vbar(i,j,knew)*vmask(i,j) cff5=ABS(ABS(vmask_wet(i,j))-1.0_r8) cff6=0.5_r8+DSIGN(0.5_r8,vbar(i,j,knew))*vmask_wet(i,j) cff7=0.5_r8*vmask_wet(i,j)*cff5+cff6*(1.0_r8-cff5) vbar(i,j,knew)=vbar(i,j,knew)*cff7 fac1=1.0_r8/cff rhs_vbar(i,j)=((vbar(i,j,knew)*(Dnew(i,j)+Dnew(i,j-1))- & & vbar(i,j,kstp)*(Dstp(i,j)+Dstp(i,j-1)))* & & fac1- & & cff2*rvbar(i,j,kstp)+ & & cff3*rvbar(i,j,ptsk))*cff4 END DO END DO END IF ! ! If predictor step, load right-side-term into shared arrays for ! future use during the subsequent corrector step. ! IF (PREDICTOR_2D_STEP(ng)) THEN DO j=Jstr,Jend DO i=IstrU,Iend rubar(i,j,krhs)=rhs_ubar(i,j) END DO END DO DO j=JstrV,Jend DO i=Istr,Iend rvbar(i,j,krhs)=rhs_vbar(i,j) END DO END DO END IF ! !----------------------------------------------------------------------- ! Apply lateral boundary conditions. !----------------------------------------------------------------------- ! CALL u2dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, knew, & & ubar, vbar, zeta) CALL v2dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & krhs, kstp, knew, & & ubar, vbar, zeta) ! !----------------------------------------------------------------------- ! Apply mass point sources. !----------------------------------------------------------------------- ! DO j=Jstr-1,Jend+1 DO i=Istr-1,Iend+1 Dnew(i,j)=zeta(i,j,knew)+h(i,j) END DO END DO DO is=1,Nsrc i=Isrc(is) j=Jsrc(is) IF (((IstrR.le.i).and.(i.le.IendR)).and. & & ((JstrR.le.j).and.(j.le.JendR))) THEN IF (INT(Dsrc(is)).eq.0) THEN cff=1.0_r8/(on_u(i,j)*0.5_r8*(Dnew(i-1,j)+Dnew(i,j))) ubar(i,j,knew)=Qbar(is)*cff DU_avg1(i,j)=Qbar(is) ELSE cff=1.0_r8/(om_v(i,j)*0.5_r8*(Dnew(i,j-1)+Dnew(i,j))) vbar(i,j,knew)=Qbar(is)*cff DV_avg1(i,j)=Qbar(is) END IF END IF END DO ! !----------------------------------------------------------------------- ! Exchange boundary information. !----------------------------------------------------------------------- ! CALL mp_exchange2d (ng, tile, iNLM, 2, & & LBi, UBi, LBj, UBj, & & NghostPoints, EWperiodic, NSperiodic, & & ubar(:,:,knew), & & vbar(:,:,knew)) RETURN END SUBROUTINE step2d_tile END MODULE step2d_mod