MODULE t3dmix_mod ! !svn $Id: t3dmix.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 routine computes horizontal mixing of tracers. ! ! ! !======================================================================= ! implicit none PRIVATE PUBLIC t3dmix2 CONTAINS SUBROUTINE t3dmix2 (ng, tile) ! !svn $Id: t3dmix2_geo.h 362 2009-07-02 21:48:48Z arango $ !************************************************** Hernan G. Arango *** ! Copyright (c) 2002-2009 The ROMS/TOMS Group ! ! Licensed under a MIT/X style license ! ! See License_ROMS.txt ! !*********************************************************************** ! ! ! This subroutine computes horizontal harmonic mixing of tracers ! ! along geopotential surfaces. ! ! ! !*********************************************************************** ! USE mod_param USE mod_grid USE mod_mixing USE mod_ocean 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, 25) CALL t3dmix2_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nrhs(ng), nstp(ng), nnew(ng), & & GRID(ng) % umask, & & GRID(ng) % vmask, & & GRID(ng) % om_v, & & GRID(ng) % on_u, & & GRID(ng) % pm, & & GRID(ng) % pn, & & GRID(ng) % Hz, & & GRID(ng) % z_r, & & MIXING(ng) % diff2, & & OCEAN(ng) % t) CALL wclock_off (ng, iNLM, 25) RETURN END SUBROUTINE t3dmix2 ! !*********************************************************************** SUBROUTINE t3dmix2_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nrhs, nstp, nnew, & & umask, vmask, & & om_v, on_u, pm, pn, & & Hz, z_r, & & diff2, & & t) !*********************************************************************** ! USE mod_param 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) :: nrhs, nstp, nnew real(r8), intent(in) :: umask(LBi:,LBj:) real(r8), intent(in) :: vmask(LBi:,LBj:) real(r8), intent(in) :: diff2(LBi:,LBj:,:) real(r8), intent(in) :: om_v(LBi:,LBj:) real(r8), intent(in) :: on_u(LBi:,LBj:) real(r8), intent(in) :: pm(LBi:,LBj:) real(r8), intent(in) :: pn(LBi:,LBj:) real(r8), intent(in) :: Hz(LBi:,LBj:,:) real(r8), intent(in) :: z_r(LBi:,LBj:,:) real(r8), intent(inout) :: t(LBi:,LBj:,:,:,:) ! ! Local variable declarations. ! integer :: i, itrc, j, k, k1, k2 real(r8) :: cff, cff1, cff2, cff3, cff4 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FE real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FX real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: FS real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTdz real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTdx real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTde real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dZdx real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dZde ! !----------------------------------------------------------------------- ! 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) ! !----------------------------------------------------------------------- ! Compute horizontal harmonic diffusion along geopotential surfaces. !----------------------------------------------------------------------- ! ! Compute horizontal and vertical gradients. Notice the recursive ! blocking sequence. The vertical placement of the gradients is: ! ! dTdx,dTde(:,:,k1) k rho-points ! dTdx,dTde(:,:,k2) k+1 rho-points ! FS,dTdz(:,:,k1) k-1/2 W-points ! FS,dTdz(:,:,k2) k+1/2 W-points ! T_LOOP : DO itrc=1,NT(ng) k2=1 K_LOOP : DO k=0,N(ng) k1=k2 k2=3-k1 IF (k.lt.N(ng)) THEN DO j=Jstr,Jend DO i=Istr,Iend+1 cff=0.5_r8*(pm(i,j)+pm(i-1,j)) cff=cff*umask(i,j) dZdx(i,j,k2)=cff*(z_r(i ,j,k+1)- & & z_r(i-1,j,k+1)) dTdx(i,j,k2)=cff*(t(i ,j,k+1,nrhs,itrc)- & & t(i-1,j,k+1,nrhs,itrc)) END DO END DO DO j=Jstr,Jend+1 DO i=Istr,Iend cff=0.5_r8*(pn(i,j)+pn(i,j-1)) cff=cff*vmask(i,j) dZde(i,j,k2)=cff*(z_r(i,j ,k+1)- & & z_r(i,j-1,k+1)) dTde(i,j,k2)=cff*(t(i,j ,k+1,nrhs,itrc)- & & t(i,j-1,k+1,nrhs,itrc)) END DO END DO END IF IF ((k.eq.0).or.(k.eq.N(ng))) THEN DO j=Jstr-1,Jend+1 DO i=Istr-1,Iend+1 dTdz(i,j,k2)=0.0_r8 FS(i,j,k2)=0.0_r8 END DO END DO ELSE DO j=Jstr-1,Jend+1 DO i=Istr-1,Iend+1 cff=1.0_r8/(z_r(i,j,k+1)-z_r(i,j,k)) dTdz(i,j,k2)=cff*(t(i,j,k+1,nrhs,itrc)- & & t(i,j,k ,nrhs,itrc)) END DO END DO END IF ! ! Compute components of the rotated tracer flux (T m3/s) along ! geopotential surfaces. ! IF (k.gt.0) THEN DO j=Jstr,Jend DO i=Istr,Iend+1 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i-1,j,itrc))* & & on_u(i,j) FX(i,j)=cff* & & (Hz(i,j,k)+Hz(i-1,j,k))* & & (dTdx(i,j,k1)- & & 0.5_r8*(MIN(dZdx(i,j,k1),0.0_r8)* & & (dTdz(i-1,j,k1)+ & & dTdz(i ,j,k2))+ & & MAX(dZdx(i,j,k1),0.0_r8)* & & (dTdz(i-1,j,k2)+ & & dTdz(i ,j,k1)))) END DO END DO DO j=Jstr,Jend+1 DO i=Istr,Iend cff=0.25_r8*(diff2(i,j,itrc)+diff2(i,j-1,itrc))* & & om_v(i,j) FE(i,j)=cff* & & (Hz(i,j,k)+Hz(i,j-1,k))* & & (dTde(i,j,k1)- & & 0.5_r8*(MIN(dZde(i,j,k1),0.0_r8)* & & (dTdz(i,j-1,k1)+ & & dTdz(i,j ,k2))+ & & MAX(dZde(i,j,k1),0.0_r8)* & & (dTdz(i,j-1,k2)+ & & dTdz(i,j ,k1)))) END DO END DO IF (k.lt.N(ng)) THEN DO j=Jstr,Jend DO i=Istr,Iend cff=0.5_r8*diff2(i,j,itrc) cff1=MIN(dZdx(i ,j,k1),0.0_r8) cff2=MIN(dZdx(i+1,j,k2),0.0_r8) cff3=MAX(dZdx(i ,j,k2),0.0_r8) cff4=MAX(dZdx(i+1,j,k1),0.0_r8) FS(i,j,k2)=cff* & & (cff1*(cff1*dTdz(i,j,k2)-dTdx(i ,j,k1))+ & & cff2*(cff2*dTdz(i,j,k2)-dTdx(i+1,j,k2))+ & & cff3*(cff3*dTdz(i,j,k2)-dTdx(i ,j,k2))+ & & cff4*(cff4*dTdz(i,j,k2)-dTdx(i+1,j,k1))) cff1=MIN(dZde(i,j ,k1),0.0_r8) cff2=MIN(dZde(i,j+1,k2),0.0_r8) cff3=MAX(dZde(i,j ,k2),0.0_r8) cff4=MAX(dZde(i,j+1,k1),0.0_r8) FS(i,j,k2)=FS(i,j,k2)+ & & cff* & & (cff1*(cff1*dTdz(i,j,k2)-dTde(i,j ,k1))+ & & cff2*(cff2*dTdz(i,j,k2)-dTde(i,j+1,k2))+ & & cff3*(cff3*dTdz(i,j,k2)-dTde(i,j ,k2))+ & & cff4*(cff4*dTdz(i,j,k2)-dTde(i,j+1,k1))) END DO END DO END IF ! ! Time-step harmonic, geopotential diffusion term (m Tunits). ! DO j=Jstr,Jend DO i=Istr,Iend cff=dt(ng)*pm(i,j)*pn(i,j)* & & (FX(i+1,j )-FX(i,j)+ & & FE(i ,j+1)-FE(i,j))+ & & dt(ng)*(FS(i,j,k2)-FS(i,j,k1)) t(i,j,k,nnew,itrc)=t(i,j,k,nnew,itrc)+cff cff1=1.0_r8/Hz(i,j,k) t(i,j,k,3,itrc)=cff1*t(i,j,k,nnew,itrc) END DO END DO END IF END DO K_LOOP END DO T_LOOP RETURN END SUBROUTINE t3dmix2_tile END MODULE t3dmix_mod