MODULE tl_rhs3d_mod ! !git $Id$ !svn $Id: tl_rhs3d.F 1180 2023-07-13 02:42:10Z arango $ !================================================== Hernan G. Arango === ! Copyright (c) 2002-2023 The ROMS/TOMS Group Andrew M. Moore ! ! Licensed under a MIT/X style license ! ! See License_ROMS.md ! !======================================================================= ! ! ! This subroutine evaluates tangent linear right-hand-side terms ! ! for 3D momentum and tracers equations ! ! ! ! BASIC STATE variables needed: Hz, Huon, HVom, u, v, W, uclm, vclm, ! ! sustr, svstr, bustr, bvstr ! ! ! !======================================================================= ! implicit none ! PRIVATE PUBLIC :: tl_rhs3d ! CONTAINS ! !*********************************************************************** SUBROUTINE tl_rhs3d (ng, tile) !*********************************************************************** ! USE mod_param USE mod_coupling USE mod_forces USE mod_grid USE mod_ocean USE mod_stepping ! USE tl_pre_step3d_mod, ONLY : tl_pre_step3d USE tl_prsgrd_mod, ONLY : tl_prsgrd USE tl_t3dmix2_mod, ONLY : tl_t3dmix2 USE tl_uv3dmix2_mod, ONLY : tl_uv3dmix2 ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile ! ! Local variable declarations. ! character (len=*), parameter :: MyFile = & & "ROMS/Tangent/tl_rhs3d.F" ! 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) directions and ! MAX(I,J) directions. ! LBij=BOUNDS(ng)%LBij UBij=BOUNDS(ng)%UBij ! !----------------------------------------------------------------------- ! Initialize computations for new time step of the 3D primitive ! variables. !----------------------------------------------------------------------- ! CALL tl_pre_step3d (ng, tile) ! !----------------------------------------------------------------------- ! Compute baroclinic pressure gradient. !----------------------------------------------------------------------- ! CALL tl_prsgrd (ng, tile) ! !----------------------------------------------------------------------- ! Compute horizontal harmonic mixing of tracer type variables. !----------------------------------------------------------------------- ! CALL tl_t3dmix2 (ng, tile) ! !----------------------------------------------------------------------- ! Compute right-hand-side terms for the 3D momentum equations. !----------------------------------------------------------------------- ! CALL wclock_on (ng, iTLM, 21, 128, MyFile) CALL tl_rhs3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nrhs(ng), & & GRID(ng) % Hz, & & GRID(ng) % tl_Hz, & & GRID(ng) % Huon, & & GRID(ng) % tl_Huon, & & GRID(ng) % Hvom, & & GRID(ng) % tl_Hvom, & & GRID(ng) % dmde, & & GRID(ng) % dndx, & & GRID(ng) % fomn, & & GRID(ng) % om_u, & & GRID(ng) % om_v, & & GRID(ng) % on_u, & & GRID(ng) % on_v, & & GRID(ng) % pm, & & GRID(ng) % pn, & & FORCES(ng) % bustr, & & FORCES(ng) % tl_bustr, & & FORCES(ng) % bvstr, & & FORCES(ng) % tl_bvstr, & & FORCES(ng) % sustr, & & FORCES(ng) % tl_sustr, & & FORCES(ng) % svstr, & & FORCES(ng) % tl_svstr, & & OCEAN(ng) % u, & & OCEAN(ng) % tl_u, & & OCEAN(ng) % v, & & OCEAN(ng) % tl_v, & & OCEAN(ng) % W, & & OCEAN(ng) % tl_W, & & COUPLING(ng) % tl_rufrc, & & COUPLING(ng) % tl_rvfrc, & & OCEAN(ng) % tl_ru, & & OCEAN(ng) % tl_rv) CALL wclock_off (ng, iTLM, 21, 190, MyFile) ! !----------------------------------------------------------------------- ! Compute horizontal, harmonic mixing of momentum. !----------------------------------------------------------------------- ! CALL tl_uv3dmix2 (ng, tile) RETURN END SUBROUTINE tl_rhs3d ! !*********************************************************************** SUBROUTINE tl_rhs3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nrhs, & & Hz, tl_Hz, & & Huon, tl_Huon, & & Hvom, tl_Hvom, & & dmde, dndx, & & fomn, & & om_u, om_v, on_u, on_v, pm, pn, & & bustr, tl_bustr, & & bvstr, tl_bvstr, & & sustr, tl_sustr, & & svstr, tl_svstr, & & u, tl_u, & & v, tl_v, & & W, tl_W, & & tl_rufrc, & & tl_rvfrc, & & tl_ru, tl_rv) !*********************************************************************** ! USE mod_param USE mod_clima 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 ! 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) :: dmde(LBi:,LBj:) real(r8), intent(in) :: dndx(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) :: pm(LBi:,LBj:) real(r8), intent(in) :: pn(LBi:,LBj:) real(r8), intent(in) :: bustr(LBi:,LBj:) real(r8), intent(in) :: bvstr(LBi:,LBj:) real(r8), intent(in) :: sustr(LBi:,LBj:) real(r8), intent(in) :: svstr(LBi:,LBj:) 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(in) :: tl_Hz(LBi:,LBj:,:) real(r8), intent(in) :: tl_Huon(LBi:,LBj:,:) real(r8), intent(in) :: tl_Hvom(LBi:,LBj:,:) real(r8), intent(in) :: tl_bustr(LBi:,LBj:) real(r8), intent(in) :: tl_bvstr(LBi:,LBj:) real(r8), intent(in) :: tl_sustr(LBi:,LBj:) real(r8), intent(in) :: tl_svstr(LBi:,LBj:) real(r8), intent(in) :: tl_u(LBi:,LBj:,:,:) real(r8), intent(in) :: tl_v(LBi:,LBj:,:,:) real(r8), intent(in) :: tl_W(LBi:,LBj:,0:) real(r8), intent(inout) :: tl_ru(LBi:,LBj:,0:,:) real(r8), intent(inout) :: tl_rv(LBi:,LBj:,0:,:) real(r8), intent(out) :: tl_rufrc(LBi:,LBj:) real(r8), intent(out) :: tl_rvfrc(LBi:,LBj:) ! ! Local variable declarations. ! integer :: i, j, k real(r8), parameter :: Gadv = -0.25_r8 real(r8) :: cff, cff1, cff2, cff3, cff4 real(r8) :: tl_cff, tl_cff1, tl_cff2, tl_cff3, tl_cff4 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: CF real(r8), dimension(IminS:ImaxS,0:N(ng)) :: DC real(r8), dimension(IminS:ImaxS,0:N(ng)) :: FC real(r8), dimension(IminS:ImaxS,0:N(ng)) :: tl_CF real(r8), dimension(IminS:ImaxS,0:N(ng)) :: tl_DC real(r8), dimension(IminS:ImaxS,0:N(ng)) :: tl_FC real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Huee real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Huxx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Hvee real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Hvxx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: UFx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: UFe real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Uwrk real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: VFx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: VFe real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Vwrk real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: uee real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: uxx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: vee real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: vxx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: wrk real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_Huee real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_Huxx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_Hvee real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_Hvxx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_UFx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_UFe real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_Uwrk real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_VFx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_VFe real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_Vwrk real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_uee real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_uxx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_vee real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_vxx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_wrk ! !----------------------------------------------------------------------- ! 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, IstrB, IstrP, IstrR, IstrT, IstrM, IstrU integer :: Iend, IendB, IendP, IendR, IendT integer :: Jstr, JstrB, JstrP, JstrR, JstrT, JstrM, JstrV integer :: Jend, JendB, JendP, JendR, JendT integer :: Istrm3, Istrm2, Istrm1, IstrUm2, IstrUm1 integer :: Iendp1, Iendp2, Iendp2i, Iendp3 integer :: Jstrm3, Jstrm2, Jstrm1, JstrVm2, JstrVm1 integer :: Jendp1, Jendp2, Jendp2i, Jendp3 ! Istr =BOUNDS(ng) % Istr (tile) IstrB =BOUNDS(ng) % IstrB (tile) IstrM =BOUNDS(ng) % IstrM (tile) IstrP =BOUNDS(ng) % IstrP (tile) IstrR =BOUNDS(ng) % IstrR (tile) IstrT =BOUNDS(ng) % IstrT (tile) IstrU =BOUNDS(ng) % IstrU (tile) Iend =BOUNDS(ng) % Iend (tile) IendB =BOUNDS(ng) % IendB (tile) IendP =BOUNDS(ng) % IendP (tile) IendR =BOUNDS(ng) % IendR (tile) IendT =BOUNDS(ng) % IendT (tile) Jstr =BOUNDS(ng) % Jstr (tile) JstrB =BOUNDS(ng) % JstrB (tile) JstrM =BOUNDS(ng) % JstrM (tile) JstrP =BOUNDS(ng) % JstrP (tile) JstrR =BOUNDS(ng) % JstrR (tile) JstrT =BOUNDS(ng) % JstrT (tile) JstrV =BOUNDS(ng) % JstrV (tile) Jend =BOUNDS(ng) % Jend (tile) JendB =BOUNDS(ng) % JendB (tile) JendP =BOUNDS(ng) % JendP (tile) JendR =BOUNDS(ng) % JendR (tile) JendT =BOUNDS(ng) % JendT (tile) ! Istrm3 =BOUNDS(ng) % Istrm3 (tile) ! Istr-3 Istrm2 =BOUNDS(ng) % Istrm2 (tile) ! Istr-2 Istrm1 =BOUNDS(ng) % Istrm1 (tile) ! Istr-1 IstrUm2=BOUNDS(ng) % IstrUm2(tile) ! IstrU-2 IstrUm1=BOUNDS(ng) % IstrUm1(tile) ! IstrU-1 Iendp1 =BOUNDS(ng) % Iendp1 (tile) ! Iend+1 Iendp2 =BOUNDS(ng) % Iendp2 (tile) ! Iend+2 Iendp2i=BOUNDS(ng) % Iendp2i(tile) ! Iend+2 interior Iendp3 =BOUNDS(ng) % Iendp3 (tile) ! Iend+3 Jstrm3 =BOUNDS(ng) % Jstrm3 (tile) ! Jstr-3 Jstrm2 =BOUNDS(ng) % Jstrm2 (tile) ! Jstr-2 Jstrm1 =BOUNDS(ng) % Jstrm1 (tile) ! Jstr-1 JstrVm2=BOUNDS(ng) % JstrVm2(tile) ! JstrV-2 JstrVm1=BOUNDS(ng) % JstrVm1(tile) ! JstrV-1 Jendp1 =BOUNDS(ng) % Jendp1 (tile) ! Jend+1 Jendp2 =BOUNDS(ng) % Jendp2 (tile) ! Jend+2 Jendp2i=BOUNDS(ng) % Jendp2i(tile) ! Jend+2 interior Jendp3 =BOUNDS(ng) % Jendp3 (tile) ! Jend+3 ! K_LOOP : DO k=1,N(ng) ! !----------------------------------------------------------------------- ! Add in Coriolis terms. !----------------------------------------------------------------------- ! DO j=JstrV-1,Jend DO i=IstrU-1,Iend cff=0.5_r8*Hz(i,j,k)*fomn(i,j) tl_cff=0.5_r8*tl_Hz(i,j,k)*fomn(i,j) !^ UFx(i,j)=cff*(v(i,j ,k,nrhs)+ & !^ & v(i,j+1,k,nrhs)) !^ tl_UFx(i,j)=tl_cff*(v(i,j ,k,nrhs)+ & & v(i,j+1,k,nrhs))+ & & cff*(tl_v(i,j ,k,nrhs)+ & & tl_v(i,j+1,k,nrhs)) !^ VFe(i,j)=cff*(u(i ,j,k,nrhs)+ & !^ & u(i+1,j,k,nrhs)) !^ tl_VFe(i,j)=tl_cff*(u(i ,j,k,nrhs)+ & & u(i+1,j,k,nrhs))+ & & cff*(tl_u(i ,j,k,nrhs)+ & & tl_u(i+1,j,k,nrhs)) END DO END DO DO j=Jstr,Jend DO i=IstrU,Iend !^ cff1=0.5_r8*(UFx(i,j)+UFx(i-1,j)) !^ tl_cff1=0.5_r8*(tl_UFx(i,j)+tl_UFx(i-1,j)) !^ ru(i,j,k,nrhs)=ru(i,j,k,nrhs)+cff1 !^ tl_ru(i,j,k,nrhs)=tl_ru(i,j,k,nrhs)+tl_cff1 END DO END DO DO j=JstrV,Jend DO i=Istr,Iend !^ cff1=0.5_r8*(VFe(i,j)+VFe(i,j-1)) !^ tl_cff1=0.5_r8*(tl_VFe(i,j)+tl_VFe(i,j-1)) !^ rv(i,j,k,nrhs)=rv(i,j,k,nrhs)-cff1 !^ tl_rv(i,j,k,nrhs)=tl_rv(i,j,k,nrhs)-tl_cff1 END DO END DO ! !----------------------------------------------------------------------- ! Add in curvilinear transformation terms. !----------------------------------------------------------------------- ! DO j=JstrV-1,Jend DO i=IstrU-1,Iend cff1=0.5_r8*(v(i,j ,k,nrhs)+ & & v(i,j+1,k,nrhs)) tl_cff1=0.5_r8*(tl_v(i,j ,k,nrhs)+ & & tl_v(i,j+1,k,nrhs)) cff2=0.5_r8*(u(i ,j,k,nrhs)+ & & u(i+1,j,k,nrhs)) tl_cff2=0.5_r8*(tl_u(i ,j,k,nrhs)+ & & tl_u(i+1,j,k,nrhs)) cff3=cff1*dndx(i,j) tl_cff3=tl_cff1*dndx(i,j) cff4=cff2*dmde(i,j) tl_cff4=tl_cff2*dmde(i,j) cff=Hz(i,j,k)*(cff3-cff4) tl_cff=tl_Hz(i,j,k)*(cff3-cff4)+ & & Hz(i,j,k)*(tl_cff3-tl_cff4) !^ UFx(i,j)=cff*cff1 !^ tl_UFx(i,j)=tl_cff*cff1+cff*tl_cff1 !^ VFe(i,j)=cff*cff2 !^ tl_VFe(i,j)=tl_cff*cff2+cff*tl_cff2 END DO END DO DO j=Jstr,Jend DO i=IstrU,Iend !^ cff1=0.5_r8*(UFx(i,j)+UFx(i-1,j)) !^ tl_cff1=0.5_r8*(tl_UFx(i,j)+tl_UFx(i-1,j)) !^ ru(i,j,k,nrhs)=ru(i,j,k,nrhs)+cff1 !^ tl_ru(i,j,k,nrhs)=tl_ru(i,j,k,nrhs)+tl_cff1 END DO END DO DO j=JstrV,Jend DO i=Istr,Iend !^ cff1=0.5_r8*(VFe(i,j)+VFe(i,j-1)) !^ tl_cff1=0.5_r8*(tl_VFe(i,j)+tl_VFe(i,j-1)) !^ rv(i,j,k,nrhs)=rv(i,j,k,nrhs)-cff1 !^ tl_rv(i,j,k,nrhs)=tl_rv(i,j,k,nrhs)-tl_cff1 END DO END DO ! !----------------------------------------------------------------------- ! Add in nudging of 3D momentum climatology. !----------------------------------------------------------------------- ! IF (LnudgeM3CLM(ng)) THEN DO j=Jstr,Jend DO i=IstrU,Iend cff=0.25_r8*(CLIMA(ng)%M3nudgcof(i-1,j,k)+ & & CLIMA(ng)%M3nudgcof(i ,j,k))* & & om_u(i,j)*on_u(i,j) !^ ru(i,j,k,nrhs)=ru(i,j,k,nrhs)+ & !^ & cff*(Hz(i-1,j,k)+Hz(i,j,k))* & !^ & (CLIMA(ng)%uclm(i,j,k)- & !^ & u(i,j,k,nrhs)) !^ tl_ru(i,j,k,nrhs)=tl_ru(i,j,k,nrhs)+ & & cff*((Hz(i-1,j,k)+Hz(i,j,k))* & & (-tl_u(i,j,k,nrhs))+ & & (tl_Hz(i-1,j,k)+tl_Hz(i,j,k))* & & (CLIMA(ng)%uclm(i,j,k)- & & u(i,j,k,nrhs))) END DO END DO DO j=JstrV,Jend DO i=Istr,Iend cff=0.25_r8*(CLIMA(ng)%M3nudgcof(i,j-1,k)+ & & CLIMA(ng)%M3nudgcof(i,j ,k))* & & om_v(i,j)*on_v(i,j) !^ rv(i,j,k,nrhs)=rv(i,j,k,nrhs)+ & !^ & cff*(Hz(i,j-1,k)+Hz(i,j,k))* & !^ & (CLIMA(ng)%vclm(i,j,k)- & !^ & v(i,j,k,nrhs)) !^ tl_rv(i,j,k,nrhs)=tl_rv(i,j,k,nrhs)+ & & cff*((Hz(i,j-1,k)+Hz(i,j,k))* & & (-tl_v(i,j,k,nrhs))+ & & (tl_Hz(i,j-1,k)+tl_Hz(i,j,k))* & & (CLIMA(ng)%vclm(i,j,k)- & & v(i,j,k,nrhs))) END DO END DO END IF ! !----------------------------------------------------------------------- ! Add in horizontal advection of momentum. !----------------------------------------------------------------------- ! ! Compute diagonal [UFx,VFe] and off-diagonal [UFe,VFx] components ! of tensor of momentum flux due to horizontal advection. ! DO j=Jstr,Jend DO i=IstrUm1,Iendp1 uxx(i,j)=u(i-1,j,k,nrhs)-2.0_r8*u(i,j,k,nrhs)+ & & u(i+1,j,k,nrhs) tl_uxx(i,j)=tl_u(i-1,j,k,nrhs)-2.0_r8*tl_u(i,j,k,nrhs)+ & & tl_u(i+1,j,k,nrhs) Huxx(i,j)=Huon(i-1,j,k)-2.0_r8*Huon(i,j,k)+Huon(i+1,j,k) tl_Huxx(i,j)=tl_Huon(i-1,j,k)-2.0_r8*tl_Huon(i,j,k)+ & & tl_Huon(i+1,j,k) END DO END DO IF (.not.(CompositeGrid(iwest,ng).or.EWperiodic(ng))) THEN IF (DOMAIN(ng)%Western_Edge(tile)) THEN DO j=Jstr,Jend uxx (Istr,j)=uxx (Istr+1,j) tl_uxx (Istr,j)=tl_uxx (Istr+1,j) Huxx(Istr,j)=Huxx(Istr+1,j) tl_Huxx(Istr,j)=tl_Huxx(Istr+1,j) END DO END IF END IF IF (.not.(CompositeGrid(ieast,ng).or.EWperiodic(ng))) THEN IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN DO j=Jstr,Jend uxx (Iend+1,j)=uxx (Iend,j) tl_uxx (Iend+1,j)=tl_uxx (Iend,j) Huxx(Iend+1,j)=Huxx(Iend,j) tl_Huxx(Iend+1,j)=tl_Huxx(Iend,j) END DO END IF END IF ! ! Third-order, upstream bias u-momentum advection with velocity ! dependent hyperdiffusion. ! DO j=Jstr,Jend DO i=IstrU-1,Iend cff1=u(i ,j,k,nrhs)+ & & u(i+1,j,k,nrhs) tl_cff1=tl_u(i ,j,k,nrhs)+ & & tl_u(i+1,j,k,nrhs) IF (cff1.gt.0.0_r8) THEN cff=uxx(i,j) tl_cff=tl_uxx(i,j) ELSE cff=uxx(i+1,j) tl_cff=tl_uxx(i+1,j) END IF !^ UFx(i,j)=0.25_r8*(cff1+Gadv*cff)* & !^ & (Huon(i ,j,k)+ & !^ & Huon(i+1,j,k)+ & !^ & Gadv*0.5_r8*(Huxx(i ,j)+ & !^ & Huxx(i+1,j))) !^ tl_UFx(i,j)=0.25_r8* & & ((tl_cff1+Gadv*tl_cff)* & & (Huon(i ,j,k)+ & & Huon(i+1,j,k)+ & & Gadv*0.5_r8*(Huxx(i ,j)+ & & Huxx(i+1,j)))+ & & (cff1+Gadv*cff)* & & (tl_Huon(i ,j,k)+ & & tl_Huon(i+1,j,k)+ & & Gadv*0.5_r8*(tl_Huxx(i ,j)+ & & tl_Huxx(i+1,j)))) END DO END DO DO j=Jstrm1,Jendp1 DO i=IstrU,Iend uee(i,j)=u(i,j-1,k,nrhs)-2.0_r8*u(i,j,k,nrhs)+ & & u(i,j+1,k,nrhs) tl_uee(i,j)=tl_u(i,j-1,k,nrhs)-2.0_r8*tl_u(i,j,k,nrhs)+ & & tl_u(i,j+1,k,nrhs) END DO END DO IF (.not.(CompositeGrid(isouth,ng).or.NSperiodic(ng))) THEN IF (DOMAIN(ng)%Southern_Edge(tile)) THEN DO i=IstrU,Iend uee(i,Jstr-1)=uee(i,Jstr) tl_uee(i,Jstr-1)=tl_uee(i,Jstr) END DO END IF END IF IF (.not.(CompositeGrid(inorth,ng).or.NSperiodic(ng))) THEN IF (DOMAIN(ng)%Northern_Edge(tile)) THEN DO i=IstrU,Iend uee(i,Jend+1)=uee(i,Jend) tl_uee(i,Jend+1)=tl_uee(i,Jend) END DO END IF END IF DO j=Jstr,Jend+1 DO i=IstrU-1,Iend Hvxx(i,j)=Hvom(i-1,j,k)-2.0_r8*Hvom(i,j,k)+Hvom(i+1,j,k) tl_Hvxx(i,j)=tl_Hvom(i-1,j,k)-2.0_r8*tl_Hvom(i,j,k)+ & & tl_Hvom(i+1,j,k) END DO END DO DO j=Jstr,Jend+1 DO i=IstrU,Iend cff1=u(i,j ,k,nrhs)+ & & u(i,j-1,k,nrhs) tl_cff1=tl_u(i,j,k,nrhs)+ & & tl_u(i,j-1,k,nrhs) cff2=Hvom(i,j,k)+Hvom(i-1,j,k) tl_cff2=tl_Hvom(i,j,k)+tl_Hvom(i-1,j,k) IF (cff2.gt.0.0_r8) THEN cff=uee(i,j-1) tl_cff=tl_uee(i,j-1) ELSE cff=uee(i,j) tl_cff=tl_uee(i,j) END IF !^ UFe(i,j)=0.25_r8*(cff1+Gadv*cff)* & !^ & (cff2+Gadv*0.5_r8*(Hvxx(i ,j)+ & !^ & Hvxx(i-1,j))) !^ tl_UFe(i,j)=0.25_r8* & & ((tl_cff1+Gadv*tl_cff)* & & (cff2+Gadv*0.5_r8*(Hvxx(i ,j)+ & & Hvxx(i-1,j)))+ & & (cff1+Gadv*cff)* & & (tl_cff2+Gadv*0.5_r8*(tl_Hvxx(i ,j)+ & & tl_Hvxx(i-1,j)))) END DO END DO DO j=JstrV,Jend DO i=Istrm1,Iendp1 vxx(i,j)=v(i-1,j,k,nrhs)-2.0_r8*v(i,j,k,nrhs)+ & & v(i+1,j,k,nrhs) tl_vxx(i,j)=tl_v(i-1,j,k,nrhs)-2.0_r8*tl_v(i,j,k,nrhs)+ & & tl_v(i+1,j,k,nrhs) END DO END DO IF (.not.(CompositeGrid(iwest,ng).or.EWperiodic(ng))) THEN IF (DOMAIN(ng)%Western_Edge(tile)) THEN DO j=JstrV,Jend vxx(Istr-1,j)=vxx(Istr,j) tl_vxx(Istr-1,j)=tl_vxx(Istr,j) END DO END IF END IF IF (.not.(CompositeGrid(ieast,ng).or.EWperiodic(ng))) THEN IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN DO j=JstrV,Jend vxx(Iend+1,j)=vxx(Iend,j) tl_vxx(Iend+1,j)=tl_vxx(Iend,j) END DO END IF END IF DO j=JstrV-1,Jend DO i=Istr,Iend+1 Huee(i,j)=Huon(i,j-1,k)-2.0_r8*Huon(i,j,k)+Huon(i,j+1,k) tl_Huee(i,j)=tl_Huon(i,j-1,k)-2.0_r8*tl_Huon(i,j,k)+ & & tl_Huon(i,j+1,k) END DO END DO ! ! Third-order, upstream bias v-momentum advection with velocity ! dependent hyperdiffusion. ! DO j=JstrV,Jend DO i=Istr,Iend+1 cff1=v(i ,j,k,nrhs)+ & & v(i-1,j,k,nrhs) tl_cff1=tl_v(i ,j,k,nrhs)+ & & tl_v(i-1,j,k,nrhs) cff2=Huon(i,j,k)+Huon(i,j-1,k) tl_cff2=tl_Huon(i,j,k)+tl_Huon(i,j-1,k) IF (cff2.gt.0.0_r8) THEN cff=vxx(i-1,j) tl_cff=tl_vxx(i-1,j) ELSE cff=vxx(i,j) tl_cff=tl_vxx(i,j) END IF !^ VFx(i,j)=0.25_r8*(cff1+Gadv*cff)* & !^ & (cff2+Gadv*0.5_r8*(Huee(i,j )+ & !^ & Huee(i,j-1))) !^ tl_VFx(i,j)=0.25_r8* & & ((tl_cff1+Gadv*tl_cff)* & & (cff2+Gadv*0.5_r8*(Huee(i,j )+ & & Huee(i,j-1)))+ & & (cff1+Gadv*cff)* & & (tl_cff2+Gadv*0.5_r8*(tl_Huee(i,j )+ & & tl_Huee(i,j-1)))) END DO END DO DO j=JstrVm1,Jendp1 DO i=Istr,Iend vee(i,j)=v(i,j-1,k,nrhs)-2.0_r8*v(i,j,k,nrhs)+ & & v(i,j+1,k,nrhs) tl_vee(i,j)=tl_v(i,j-1,k,nrhs)-2.0_r8*tl_v(i,j,k,nrhs)+ & & tl_v(i,j+1,k,nrhs) Hvee(i,j)=Hvom(i,j-1,k)-2.0_r8*Hvom(i,j,k)+Hvom(i,j+1,k) tl_Hvee(i,j)=tl_Hvom(i,j-1,k)-2.0_r8*tl_Hvom(i,j,k)+ & & tl_Hvom(i,j+1,k) END DO END DO IF (.not.(CompositeGrid(isouth,ng).or.NSperiodic(ng))) THEN IF (DOMAIN(ng)%Southern_Edge(tile)) THEN DO i=Istr,Iend vee (i,Jstr)=vee (i,Jstr+1) tl_vee (i,Jstr)=tl_vee (i,Jstr+1) Hvee(i,Jstr)=Hvee(i,Jstr+1) tl_Hvee(i,Jstr)=tl_Hvee(i,Jstr+1) END DO END IF END IF IF (.not.(CompositeGrid(inorth,ng).or.NSperiodic(ng))) THEN IF (DOMAIN(ng)%Northern_Edge(tile)) THEN DO i=Istr,Iend vee (i,Jend+1)=vee (i,Jend) tl_vee (i,Jend+1)=tl_vee (i,Jend) Hvee(i,Jend+1)=Hvee(i,Jend) tl_Hvee(i,Jend+1)=tl_Hvee(i,Jend) END DO END IF END IF DO j=JstrV-1,Jend DO i=Istr,Iend cff1=v(i,j ,k,nrhs)+ & & v(i,j+1,k,nrhs) tl_cff1=tl_v(i,j ,k,nrhs)+ & & tl_v(i,j+1,k,nrhs) IF (cff1.gt.0.0_r8) THEN cff=vee(i,j) tl_cff=tl_vee(i,j) ELSE cff=vee(i,j+1) tl_cff=tl_vee(i,j+1) END IF !^ VFe(i,j)=0.25_r8*(cff1+Gadv*cff)* & !^ & (Hvom(i,j ,k)+ & !^ & Hvom(i,j+1,k)+ & !^ & Gadv*0.5_r8*(Hvee(i,j )+ & !^ & Hvee(i,j+1))) !^ tl_VFe(i,j)=0.25_r8* & & ((tl_cff1+Gadv*tl_cff)* & & (Hvom(i,j ,k)+ & & Hvom(i,j+1,k)+ & & Gadv*0.5_r8*(Hvee(i,j )+ & & Hvee(i,j+1)))+ & & (cff1+Gadv*cff)* & & (tl_Hvom(i,j ,k)+ & & tl_Hvom(i,j+1,k)+ & & Gadv*0.5_r8*(tl_Hvee(i,j )+ & & tl_Hvee(i,j+1)))) END DO END DO ! ! Add in horizontal advection. ! DO j=Jstr,Jend DO i=IstrU,Iend !^ cff1=UFx(i,j)-UFx(i-1,j) !^ tl_cff1=tl_UFx(i,j)-tl_UFx(i-1,j) !^ cff2=UFe(i,j+1)-UFe(i,j) !^ tl_cff2=tl_UFe(i,j+1)-tl_UFe(i,j) !^ cff=cff1+cff2 !^ tl_cff=tl_cff1+tl_cff2 !^ ru(i,j,k,nrhs)=ru(i,j,k,nrhs)-cff !^ tl_ru(i,j,k,nrhs)=tl_ru(i,j,k,nrhs)-tl_cff END DO END DO DO j=JstrV,Jend DO i=Istr,Iend !^ cff1=VFx(i+1,j)-VFx(i,j) !^ tl_cff1=tl_VFx(i+1,j)-tl_VFx(i,j) !^ cff2=VFe(i,j)-VFe(i,j-1) !^ tl_cff2=tl_VFe(i,j)-tl_VFe(i,j-1) !^ cff=cff1+cff2 !^ tl_cff=tl_cff1+tl_cff2 !^ rv(i,j,k,nrhs)=rv(i,j,k,nrhs)-cff !^ tl_rv(i,j,k,nrhs)=tl_rv(i,j,k,nrhs)-tl_cff END DO END DO END DO K_LOOP ! J_LOOP : DO j=Jstr,Jend ! !----------------------------------------------------------------------- ! Add in vertical advection. !----------------------------------------------------------------------- ! cff1=9.0_r8/16.0_r8 cff2=1.0_r8/16.0_r8 DO k=2,N(ng)-2 DO i=IstrU,Iend !^ FC(i,k)=(cff1*(u(i,j,k ,nrhs)+ & !^ & u(i,j,k+1,nrhs))- & !^ & cff2*(u(i,j,k-1,nrhs)+ & !^ & u(i,j,k+2,nrhs)))* & !^ & (cff1*(W(i ,j,k)+ & !^ & W(i-1,j,k))- & !^ & cff2*(W(i+1,j,k)+ & !^ & W(i-2,j,k))) !^ tl_FC(i,k)=(cff1*(tl_u(i,j,k ,nrhs)+ & & tl_u(i,j,k+1,nrhs))- & & cff2*(tl_u(i,j,k-1,nrhs)+ & & tl_u(i,j,k+2,nrhs)))* & & (cff1*(W(i ,j,k)+ & & W(i-1,j,k))- & & cff2*(W(i+1,j,k)+ & & W(i-2,j,k)))+ & & (cff1*(u(i,j,k ,nrhs)+ & & u(i,j,k+1,nrhs))- & & cff2*(u(i,j,k-1,nrhs)+ & & u(i,j,k+2,nrhs)))* & & (cff1*(tl_W(i ,j,k)+ & & tl_W(i-1,j,k))- & & cff2*(tl_W(i+1,j,k)+ & & tl_W(i-2,j,k))) END DO END DO DO i=IstrU,Iend !^ FC(i,N(ng))=0.0_r8 !^ tl_FC(i,N(ng))=0.0_r8 !^ FC(i,N(ng)-1)=(cff1*(u(i,j,N(ng)-1,nrhs)+ & !^ & u(i,j,N(ng) ,nrhs))- & !^ & cff2*(u(i,j,N(ng)-2,nrhs)+ & !^ & u(i,j,N(ng) ,nrhs)))* & !^ & (cff1*(W(i ,j,N(ng)-1)+ & !^ & W(i-1,j,N(ng)-1))- & !^ & cff2*(W(i+1,j,N(ng)-1)+ & !^ & W(i-2,j,N(ng)-1))) !^ tl_FC(i,N(ng)-1)=(cff1*(tl_u(i,j,N(ng)-1,nrhs)+ & & tl_u(i,j,N(ng) ,nrhs))- & & cff2*(tl_u(i,j,N(ng)-2,nrhs)+ & & tl_u(i,j,N(ng) ,nrhs)))* & & (cff1*(W(i ,j,N(ng)-1)+ & & W(i-1,j,N(ng)-1))- & & cff2*(W(i+1,j,N(ng)-1)+ & & W(i-2,j,N(ng)-1)))+ & & (cff1*(u(i,j,N(ng)-1,nrhs)+ & & u(i,j,N(ng) ,nrhs))- & & cff2*(u(i,j,N(ng)-2,nrhs)+ & & u(i,j,N(ng) ,nrhs)))* & & (cff1*(tl_W(i ,j,N(ng)-1)+ & & tl_W(i-1,j,N(ng)-1))- & & cff2*(tl_W(i+1,j,N(ng)-1)+ & & tl_W(i-2,j,N(ng)-1))) !^ FC(i,1)=(cff1*(u(i,j,1,nrhs)+ & !^ & u(i,j,2,nrhs))- & !^ & cff2*(u(i,j,1,nrhs)+ & !^ & u(i,j,3,nrhs)))* & !^ & (cff1*(W(i ,j,1)+ & !^ & W(i-1,j,1))- & !^ & cff2*(W(i+1,j,1)+ & !^ & W(i-2,j,1))) !^ tl_FC(i,1)=(cff1*(tl_u(i,j,1,nrhs)+ & & tl_u(i,j,2,nrhs))- & & cff2*(tl_u(i,j,1,nrhs)+ & & tl_u(i,j,3,nrhs)))* & & (cff1*(W(i ,j,1)+ & & W(i-1,j,1))- & & cff2*(W(i+1,j,1)+ & & W(i-2,j,1)))+ & & (cff1*(u(i,j,1,nrhs)+ & & u(i,j,2,nrhs))- & & cff2*(u(i,j,1,nrhs)+ & & u(i,j,3,nrhs)))* & & (cff1*(tl_W(i ,j,1)+ & & tl_W(i-1,j,1))- & & cff2*(tl_W(i+1,j,1)+ & & tl_W(i-2,j,1))) !^ FC(i,0)=0.0_r8 !^ tl_FC(i,0)=0.0_r8 END DO DO k=1,N(ng) DO i=IstrU,Iend !^ cff=FC(i,k)-FC(i,k-1) !^ tl_cff=tl_FC(i,k)-tl_FC(i,k-1) !^ ru(i,j,k,nrhs)=ru(i,j,k,nrhs)-cff !^ tl_ru(i,j,k,nrhs)=tl_ru(i,j,k,nrhs)-tl_cff END DO END DO IF (j.ge.JstrV) THEN cff1=9.0_r8/16.0_r8 cff2=1.0_r8/16.0_r8 DO k=2,N(ng)-2 DO i=Istr,Iend !^ FC(i,k)=(cff1*(v(i,j,k ,nrhs)+ & !^ & v(i,j,k+1,nrhs))- & !^ & cff2*(v(i,j,k-1,nrhs)+ & !^ & v(i,j,k+2,nrhs)))* & !^ & (cff1*(W(i,j ,k)+ & !^ & W(i,j-1,k))- & !^ & cff2*(W(i,j+1,k)+ & !^ & W(i,j-2,k))) !^ tl_FC(i,k)=(cff1*(tl_v(i,j,k ,nrhs)+ & & tl_v(i,j,k+1,nrhs))- & & cff2*(tl_v(i,j,k-1,nrhs)+ & & tl_v(i,j,k+2,nrhs)))* & & (cff1*(W(i,j ,k)+ & & W(i,j-1,k))- & & cff2*(W(i,j+1,k)+ & & W(i,j-2,k)))+ & & (cff1*(v(i,j,k ,nrhs)+ & & v(i,j,k+1,nrhs))- & & cff2*(v(i,j,k-1,nrhs)+ & & v(i,j,k+2,nrhs)))* & & (cff1*(tl_W(i,j ,k)+ & & tl_W(i,j-1,k))- & & cff2*(tl_W(i,j+1,k)+ & & tl_W(i,j-2,k))) END DO END DO DO i=Istr,Iend !^ FC(i,N(ng))=0.0_r8 !^ tl_FC(i,N(ng))=0.0_r8 !^ FC(i,N(ng)-1)=(cff1*(v(i,j,N(ng)-1,nrhs)+ & !^ & v(i,j,N(ng) ,nrhs))- & !^ & cff2*(v(i,j,N(ng)-2,nrhs)+ & !^ & v(i,j,N(ng) ,nrhs)))* & !^ & (cff1*(W(i,j ,N(ng)-1)+ & !^ & W(i,j-1,N(ng)-1))- & !^ & cff2*(W(i,j+1,N(ng)-1)+ & !^ & W(i,j-2,N(ng)-1))) !^ tl_FC(i,N(ng)-1)=(cff1*(tl_v(i,j,N(ng)-1,nrhs)+ & & tl_v(i,j,N(ng) ,nrhs))- & & cff2*(tl_v(i,j,N(ng)-2,nrhs)+ & & tl_v(i,j,N(ng) ,nrhs)))* & & (cff1*(W(i,j ,N(ng)-1)+ & & W(i,j-1,N(ng)-1))- & & cff2*(W(i,j+1,N(ng)-1)+ & & W(i,j-2,N(ng)-1)))+ & & (cff1*(v(i,j,N(ng)-1,nrhs)+ & & v(i,j,N(ng) ,nrhs))- & & cff2*(v(i,j,N(ng)-2,nrhs)+ & & v(i,j,N(ng) ,nrhs)))* & & (cff1*(tl_W(i,j ,N(ng)-1)+ & & tl_W(i,j-1,N(ng)-1))- & & cff2*(tl_W(i,j+1,N(ng)-1)+ & & tl_W(i,j-2,N(ng)-1))) !^ FC(i,1)=(cff1*(v(i,j,1,nrhs)+ & !^ & v(i,j,2,nrhs))- & !^ & cff2*(v(i,j,1,nrhs)+ & !^ & v(i,j,3,nrhs)))* & !^ & (cff1*(W(i,j ,1)+ & !^ & W(i,j-1,1))- & !^ & cff2*(W(i,j+1,1)+ & !^ & W(i,j-2,1))) !^ tl_FC(i,1)=(cff1*(tl_v(i,j,1,nrhs)+ & & tl_v(i,j,2,nrhs))- & & cff2*(tl_v(i,j,1,nrhs)+ & & tl_v(i,j,3,nrhs)))* & & (cff1*(W(i,j ,1)+ & & W(i,j-1,1))- & & cff2*(W(i,j+1,1)+ & & W(i,j-2,1)))+ & & (cff1*(v(i,j,1,nrhs)+ & & v(i,j,2,nrhs))- & & cff2*(v(i,j,1,nrhs)+ & & v(i,j,3,nrhs)))* & & (cff1*(tl_W(i,j ,1)+ & & tl_W(i,j-1,1))- & & cff2*(tl_W(i,j+1,1)+ & & tl_W(i,j-2,1))) !^ FC(i,0)=0.0_r8 !^ tl_FC(i,0)=0.0_r8 END DO DO k=1,N(ng) DO i=Istr,Iend !^ cff=FC(i,k)-FC(i,k-1) !^ tl_cff=tl_FC(i,k)-tl_FC(i,k-1) !^ rv(i,j,k,nrhs)=rv(i,j,k,nrhs)-cff !^ tl_rv(i,j,k,nrhs)=tl_rv(i,j,k,nrhs)-tl_cff END DO END DO END IF ! !----------------------------------------------------------------------- ! Compute forcing term for the 2D momentum equations. !----------------------------------------------------------------------- ! ! Vertically integrate baroclinic right-hand-side terms. If not ! body force stresses, add in the difference between surface and ! bottom stresses. ! DO i=IstrU,Iend !^ rufrc(i,j)=ru(i,j,1,nrhs) !^ tl_rufrc(i,j)=tl_ru(i,j,1,nrhs) END DO DO k=2,N(ng) DO i=IstrU,Iend !^ rufrc(i,j)=rufrc(i,j)+ru(i,j,k,nrhs) !^ tl_rufrc(i,j)=tl_rufrc(i,j)+tl_ru(i,j,k,nrhs) END DO END DO DO i=IstrU,Iend cff=om_u(i,j)*on_u(i,j) !^ cff1= sustr(i,j)*cff !^ tl_cff1= tl_sustr(i,j)*cff !^ cff2=-bustr(i,j)*cff !^ tl_cff2=-tl_bustr(i,j)*cff !^ rufrc(i,j)=rufrc(i,j)+cff1+cff2 !^ tl_rufrc(i,j)=tl_rufrc(i,j)+tl_cff1+tl_cff2 END DO IF (j.ge.JstrV) THEN DO i=Istr,Iend !^ rvfrc(i,j)=rv(i,j,1,nrhs) !^ tl_rvfrc(i,j)=tl_rv(i,j,1,nrhs) END DO DO k=2,N(ng) DO i=Istr,Iend !^ rvfrc(i,j)=rvfrc(i,j)+rv(i,j,k,nrhs) !^ tl_rvfrc(i,j)=tl_rvfrc(i,j)+tl_rv(i,j,k,nrhs) END DO END DO DO i=Istr,Iend cff=om_v(i,j)*on_v(i,j) !^ cff1= svstr(i,j)*cff !^ tl_cff1= tl_svstr(i,j)*cff !^ cff2=-bvstr(i,j)*cff !^ tl_cff2=-tl_bvstr(i,j)*cff !^ rvfrc(i,j)=rvfrc(i,j)+cff1+cff2 !^ tl_rvfrc(i,j)=tl_rvfrc(i,j)+tl_cff1+tl_cff2 END DO END IF END DO J_LOOP ! RETURN END SUBROUTINE tl_rhs3d_tile END MODULE tl_rhs3d_mod