MODULE step3d_t_mod ! !svn $Id: step3d_t.F 412 2009-12-03 20:46:03Z 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 time-steps tracer equations. Notice that advective ! ! and diffusive terms are time-stepped differently. ! ! ! !======================================================================= ! implicit none PRIVATE PUBLIC :: step3d_t CONTAINS ! !*********************************************************************** SUBROUTINE step3d_t (ng, tile) !*********************************************************************** ! USE mod_param USE mod_grid USE mod_mixing 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, 35) CALL step3d_t_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nrhs(ng), nstp(ng), nnew(ng), & & Msrc(ng), Nsrc(ng), & & SOURCES(ng) % Isrc, & & SOURCES(ng) % Jsrc, & & SOURCES(ng) % Tsrc, & & SOURCES(ng) % Dsrc, & & GRID(ng) % rmask, & & GRID(ng) % umask, & & GRID(ng) % vmask, & & GRID(ng) % rmask_wet, & & GRID(ng) % umask_wet, & & GRID(ng) % vmask_wet, & & GRID(ng) % omn, & & GRID(ng) % om_u, & & GRID(ng) % om_v, & & GRID(ng) % on_u, & & GRID(ng) % on_v, & & GRID(ng) % pm, & & GRID(ng) % pn, & & GRID(ng) % Hz, & & GRID(ng) % Huon, & & GRID(ng) % Hvom, & & GRID(ng) % z_r, & & MIXING(ng) % Akt, & & OCEAN(ng) % u, & & OCEAN(ng) % v, & & OCEAN(ng) % W, & & OCEAN(ng) % t) CALL wclock_off (ng, iNLM, 35) RETURN END SUBROUTINE step3d_t ! !*********************************************************************** SUBROUTINE step3d_t_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nrhs, nstp, nnew, & & Msrc, Nsrc, & & Isrc, Jsrc, Tsrc, & & Dsrc, & & rmask, umask, vmask, & & rmask_wet, umask_wet, vmask_wet, & & omn, om_u, om_v, on_u, on_v, & & pm, pn, & & Hz, Huon, Hvom, & & z_r, & & Akt, & & u, v, & & W, & & t) !*********************************************************************** ! USE mod_param USE mod_ncparam USE mod_scalars ! USE mp_exchange_mod, ONLY : mp_exchange4d USE mpdata_adiff_mod USE t3dbc_mod, ONLY : t3dbc_tile ! ! 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 integer, intent(in) :: Msrc, Nsrc ! integer, intent(in) :: Isrc(:) integer, intent(in) :: Jsrc(:) real(r8), intent(in) :: Tsrc(:,:,:) real(r8), intent(in) :: Dsrc(:) 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) :: rmask_wet(LBi:,LBj:) real(r8), intent(in) :: umask_wet(LBi:,LBj:) real(r8), intent(in) :: vmask_wet(LBi:,LBj:) real(r8), intent(in) :: omn(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) :: pm(LBi:,LBj:) real(r8), intent(in) :: pn(LBi:,LBj:) real(r8), intent(in) :: Hz(LBi:,LBj:,:) real(r8), intent(in) :: Huon(LBi:,LBj:,:) real(r8), intent(in) :: Hvom(LBi:,LBj:,:) real(r8), intent(in) :: z_r(LBi:,LBj:,:) real(r8), intent(in) :: Akt(LBi:,LBj:,0:,:) real(r8), intent(in) :: u(LBi:,LBj:,:,:) real(r8), intent(in) :: v(LBi:,LBj:,:,:) real(r8), intent(in) :: W(LBi:,LBj:,0:) real(r8), intent(inout) :: t(LBi:,LBj:,:,:,:) ! ! Local variable declarations. ! logical :: EWperiodic=.FALSE. logical :: NSperiodic=.FALSE. integer :: i, is, itrc, j, k, ltrc integer :: idiag real(r8), parameter :: eps = 1.0E-16_r8 real(r8) :: cff, cff1, cff2, cff3 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: CF real(r8), dimension(IminS:ImaxS,0:N(ng)) :: BC real(r8), dimension(IminS:ImaxS,0:N(ng)) :: DC real(r8), dimension(IminS:ImaxS,0:N(ng)) :: FC real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FE real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FX real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: curv real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: grad real(r8), dimension(IminS:ImaxS,JminS:JmaxS,N(ng)) :: oHz real(r8), dimension(IminS:ImaxS,JminS:JmaxS,N(ng),NT(ng)) :: Ta real(r8), dimension(IminS:ImaxS,JminS:JmaxS,N(ng)) :: Ua real(r8), dimension(IminS:ImaxS,JminS:JmaxS,N(ng)) :: Va real(r8), dimension(IminS:ImaxS,JminS:JmaxS,0:N(ng)) :: Wa ! !----------------------------------------------------------------------- ! 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) ! !----------------------------------------------------------------------- ! Time-step horizontal advection term. !----------------------------------------------------------------------- ! ! Compute inverse thickness. ! DO k=1,N(ng) DO j=MAX(Jstr-2,0),MIN(Jend+2,Mm(ng)+1) DO i=MAX(Istr-2,0),MIN(Iend+2,Lm(ng)+1) oHz(i,j,k)=1.0_r8/Hz(i,j,k) END DO END DO END DO ! ! Compute horizontal tracer advection fluxes. ! ! ! The MPDATA algorithm requires a three-point footprint, so exchange ! boundary data on t(:,:,:,3,:) so other processes computed earlier ! (horizontal diffusion, biology, or sediment) are accounted. ! CALL mp_exchange4d (ng, tile, iNLM, 1, & & LBi, UBi, LBj, UBj, 1, N(ng), 1, NT(ng), & & NghostPoints, EWperiodic, NSperiodic, & & t(:,:,:,3,:)) T_LOOP : DO itrc=1,NT(ng) K_LOOP : DO k=1,N(ng) ! ! First-order, upstream differences horizontal advective fluxes. ! DO j=MAX(Jstr-2,1),MIN(Jend+2,Mm(ng)) DO i=MAX(Istr-2,1),MIN(Iend+3,Lm(ng)+1) cff1=MAX(Huon(i,j,k),0.0_r8) cff2=MIN(Huon(i,j,k),0.0_r8) FX(i,j)=cff1*t(i-1,j,k,3,itrc)+ & & cff2*t(i ,j,k,3,itrc) END DO END DO DO j=MAX(Jstr-2,1),MIN(Jend+3,Mm(ng)+1) DO i=MAX(Istr-2,1),MIN(Iend+2,Lm(ng)) cff1=MAX(Hvom(i,j,k),0.0_r8) cff2=MIN(Hvom(i,j,k),0.0_r8) FE(i,j)=cff1*t(i,j-1,k,3,itrc)+ & & cff2*t(i,j ,k,3,itrc) END DO END DO ! ! Apply tracers point sources to the horizontal advection terms. ! DO is=1,Nsrc i=Isrc(is) j=Jsrc(is) IF (INT(Dsrc(is)).eq.0) THEN IF (((Istr-2.le.i).and.(i.le.Iend+2)).and. & & ((Jstr-2.le.j).and.(j.le.Jend+2))) THEN IF (LtracerSrc(itrc,ng)) THEN FX(i,j)=Huon(i,j,k)*Tsrc(is,k,itrc) ELSE IF ((rmask(i ,j).eq.0.0_r8).and. & & (rmask(i-1,j).eq.1.0_r8)) THEN FX(i,j)=Huon(i,j,k)*t(i-1,j,k,3,itrc) ELSE IF ((rmask(i ,j).eq.1.0_r8).and. & & (rmask(i-1,j).eq.0.0_r8)) THEN FX(i,j)=Huon(i,j,k)*t(i ,j,k,3,itrc) END IF END IF END IF ELSE IF (INT(Dsrc(is)).eq.1) THEN IF (((Istr-2.le.i).and.(i.le.Iend+2)).and. & & ((Jstr-2.le.j).and.(j.le.Jend+2))) THEN IF (LtracerSrc(itrc,ng)) THEN FE(i,j)=Hvom(i,j,k)*Tsrc(is,k,itrc) ELSE IF ((rmask(i,j ).eq.0.0_r8).and. & & (rmask(i,j-1).eq.1.0_r8)) THEN FE(i,j)=Hvom(i,j,k)*t(i,j-1,k,3,itrc) ELSE IF ((rmask(i,j ).eq.1.0_r8).and. & & (rmask(i,j-1).eq.0.0_r8)) THEN FE(i,j)=Hvom(i,j,k)*t(i,j ,k,3,itrc) END IF END IF END IF END IF END DO ! ! Time-step horizontal advection term. ! DO j=MAX(Jstr-2,1),MIN(Jend+2,Mm(ng)) DO i=MAX(Istr-2,1),MIN(Iend+2,Lm(ng)) cff=dt(ng)*pm(i,j)*pn(i,j) cff1=cff*(FX(i+1,j)-FX(i,j)+ & & FE(i,j+1)-FE(i,j)) cff2=Hz(i,j,k)+ & & cff*(Huon(i+1,j,k)-Huon(i,j,k)+ & & Hvom(i,j+1,k)-Hvom(i,j,k)+ & & (W(i,j,k)-W(i,j,k-1))) ! Hz_old Ta(i,j,k,itrc)=t(i,j,k,3,itrc)*cff2-cff1 END DO END DO END DO K_LOOP END DO T_LOOP ! !----------------------------------------------------------------------- ! Time-step vertical advection term. !----------------------------------------------------------------------- ! DO j=MAX(Jstr-2,1),MIN(Jend+2,Mm(ng)) DO itrc=1,NT(ng) ! ! First_order, upstream differences vertical advective flux. ! DO i=MAX(Istr-2,1),MIN(Iend+2,Lm(ng)) DO k=1,N(ng)-1 cff1=MAX(W(i,j,k),0.0_r8) cff2=MIN(W(i,j,k),0.0_r8) FC(i,k)=cff1*t(i,j,k ,3,itrc)+ & & cff2*t(i,j,k+1,3,itrc) END DO FC(i,0)=0.0_r8 FC(i,N(ng))=0.0_r8 END DO ! ! Time-step vertical advection term. ! DO i=MAX(Istr-2,1),MIN(Iend+2,Lm(ng)) CF(i,0)=dt(ng)*pm(i,j)*pn(i,j) END DO DO k=1,N(ng) DO i=MAX(Istr-2,1),MIN(Iend+2,Lm(ng)) cff1=CF(i,0)*(FC(i,k)-FC(i,k-1)) Ta(i,j,k,itrc)=(Ta(i,j,k,itrc)-cff1)*oHz(i,j,k) END DO END DO END DO END DO ! !----------------------------------------------------------------------- ! Compute anti-diffusive velocities to corrected advected tracers ! using MPDATA recursive method. Notice that pipelined J-loop ended. !----------------------------------------------------------------------- ! DO itrc=1,NT(ng) CALL mpdata_adiff_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nrhs, & & rmask, umask, vmask, & & rmask_wet, umask_wet, vmask_wet, & & pm, pn, omn, om_u, on_v, & & z_r, oHz, & & Huon, Hvom, W, & & t(:,:,:,3,itrc), & & Ta(:,:,:,itrc), Ua, Va, Wa) ! ! Compute anti-difussive corrected advection fluxes. ! DO k=1,N(ng) DO j=Jstr,Jend DO i=Istr,Iend+1 cff1=MAX(Ua(i,j,k),0.0_r8) cff2=MIN(Ua(i,j,k),0.0_r8) FX(i,j)=(cff1*Ta(i-1,j,k,itrc)+ & & cff2*Ta(i ,j,k,itrc))* & & 0.5_r8*(Hz(i,j,k)+Hz(i-1,j,k))*on_u(i,j) END DO END DO DO j=Jstr,Jend+1 DO i=Istr,Iend cff1=MAX(Va(i,j,k),0.0_r8) cff2=MIN(Va(i,j,k),0.0_r8) FE(i,j)=(cff1*Ta(i,j-1,k,itrc)+ & & cff2*Ta(i,j ,k,itrc))* & & 0.5_r8*(Hz(i,j,k)+Hz(i,j-1,k))*om_v(i,j) END DO END DO ! ! Time-step corrected horizontal advection (Tunits m). ! DO j=Jstr,Jend DO i=Istr,Iend cff1=dt(ng)*pm(i,j)*pn(i,j)* & & (FX(i+1,j)-FX(i,j)+ & & FE(i,j+1)-FE(i,j)) t(i,j,k,nnew,itrc)=Ta(i,j,k,itrc)*Hz(i,j,k)-cff1 END DO END DO END DO ! ! Compute anti-difussive corrected vertical advection flux. ! DO j=Jstr,Jend DO k=1,N(ng)-1 DO i=Istr,Iend cff1=MAX(Wa(i,j,k),0.0_r8) cff2=MIN(Wa(i,j,k),0.0_r8) FC(i,k)=cff1*Ta(i,j,k ,itrc)+ & & cff2*Ta(i,j,k+1,itrc) END DO END DO DO i=Istr,Iend FC(i,0)=0.0_r8 FC(i,N(ng))=0.0_r8 END DO ! ! Time-step corrected vertical advection (Tunits). ! DO i=Istr,Iend CF(i,0)=dt(ng)*pm(i,j)*pn(i,j) END DO DO k=1,N(ng) DO i=Istr,Iend cff1=CF(i,0)*(FC(i,k)-FC(i,k-1)) t(i,j,k,nnew,itrc)=t(i,j,k,nnew,itrc)-cff1 END DO END DO END DO END DO ! ! Start pipelined J-loop. ! DO j=Jstr,Jend ! !----------------------------------------------------------------------- ! Time-step vertical diffusion term. !----------------------------------------------------------------------- ! DO itrc=1,NT(ng) ltrc=MIN(NAT,itrc) ! ! Compute off-diagonal coefficients FC [lambda*dt*Akt/Hz] for the ! implicit vertical diffusion terms at future time step, located ! at horizontal RHO-points and vertical W-points. ! Also set FC at the top and bottom levels. ! cff=-dt(ng)*lambda DO k=1,N(ng)-1 DO i=Istr,Iend cff1=1.0_r8/(z_r(i,j,k+1)-z_r(i,j,k)) FC(i,k)=cff*cff1*Akt(i,j,k,ltrc) END DO END DO DO i=Istr,Iend FC(i,0)=0.0_r8 FC(i,N(ng))=0.0_r8 END DO ! ! Compute diagonal matrix coefficients BC and load right-hand-side ! terms for the tracer equation into DC. ! DO k=1,N(ng) DO i=Istr,Iend BC(i,k)=Hz(i,j,k)-FC(i,k)-FC(i,k-1) DC(i,k)=t(i,j,k,nnew,itrc) END DO END DO ! ! Solve the tridiagonal system. ! DO i=Istr,Iend cff=1.0_r8/BC(i,1) CF(i,1)=cff*FC(i,1) DC(i,1)=cff*DC(i,1) END DO DO k=2,N(ng)-1 DO i=Istr,Iend cff=1.0_r8/(BC(i,k)-FC(i,k-1)*CF(i,k-1)) CF(i,k)=cff*FC(i,k) DC(i,k)=cff*(DC(i,k)-FC(i,k-1)*DC(i,k-1)) END DO END DO ! ! Compute new solution by back substitution. ! DO i=Istr,Iend DC(i,N(ng))=(DC(i,N(ng))-FC(i,N(ng)-1)*DC(i,N(ng)-1))/ & & (BC(i,N(ng))-FC(i,N(ng)-1)*CF(i,N(ng)-1)) t(i,j,N(ng),nnew,itrc)=DC(i,N(ng)) END DO DO k=N(ng)-1,1,-1 DO i=Istr,Iend DC(i,k)=DC(i,k)-CF(i,k)*DC(i,k+1) t(i,j,k,nnew,itrc)=DC(i,k) END DO END DO END DO END DO ! !----------------------------------------------------------------------- ! Apply lateral boundary conditions and, if appropriate, nudge ! to tracer data and apply Land/Sea mask. !----------------------------------------------------------------------- ! DO itrc=1,NT(ng) ! ! Set lateral boundary conditions. ! CALL t3dbc_tile (ng, tile, itrc, & & LBi, UBi, LBj, UBj, N(ng), NT(ng), & & IminS, ImaxS, JminS, JmaxS, & & nstp, nnew, & & t) ! ! Apply Land/Sea mask. ! DO k=1,N(ng) DO j=JstrR,JendR DO i=IstrR,IendR t(i,j,k,nnew,itrc)=t(i,j,k,nnew,itrc)*rmask(i,j) END DO END DO END DO END DO ! ! Exchange boundary data. ! CALL mp_exchange4d (ng, tile, iNLM, 1, & & LBi, UBi, LBj, UBj, 1, N(ng), 1, NT(ng), & & NghostPoints, EWperiodic, NSperiodic, & & t(:,:,:,nnew,:)) RETURN END SUBROUTINE step3d_t_tile END MODULE step3d_t_mod