MODULE tl_uv3dmix2_mod ! !git $Id$ !svn $Id: tl_uv3dmix2_s.h 1151 2023-02-09 03:08:53Z 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 routine computes tangent linear harmonic mixing of momentum, ! ! along constant S-surfaces, from the horizontal divergence of the ! ! stress tensor. A transverse isotropy is assumed so the stress ! ! tensor is splitted into vertical and horizontal subtensors. ! ! ! ! Reference: ! ! ! ! Wajsowicz, R.C, 1993: A consistent formulation of the ! ! anisotropic stress tensor for use in models of the ! ! large-scale ocean circulation, JCP, 105, 333-338. ! ! ! ! Sadourny, R. and K. Maynard, 1997: Formulations of ! ! lateral diffusion in geophysical fluid dynamics ! ! models, In Numerical Methods of Atmospheric and ! ! Oceanic Modelling. Lin, Laprise, and Ritchie, ! ! Eds., NRC Research Press, 547-556. ! ! ! ! Griffies, S.M. and R.W. Hallberg, 2000: Biharmonic ! ! friction with a Smagorinsky-like viscosity for ! ! use in large-scale eddy-permitting ocean models, ! ! Monthly Weather Rev., 128, 8, 2935-2946. ! ! ! ! BASIC STATE variables needed: visc2, u, v, Hz ! ! ! !======================================================================= ! implicit none ! PRIVATE PUBLIC tl_uv3dmix2 ! CONTAINS ! !*********************************************************************** SUBROUTINE tl_uv3dmix2 (ng, tile) !*********************************************************************** ! USE mod_param USE mod_coupling USE mod_grid USE mod_mixing USE mod_ocean USE mod_stepping ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile ! ! Local variable declarations. ! character (len=*), parameter :: MyFile = & & "ROMS/Tangent/tl_uv3dmix2_s.h" ! 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 ! CALL wclock_on (ng, iTLM, 30, 70, MyFile) CALL tl_uv3dmix2_s_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nrhs(ng), nnew(ng), & & GRID(ng) % pmask, & & GRID(ng) % Hz, & & GRID(ng) % tl_Hz, & & GRID(ng) % om_p, & & GRID(ng) % om_r, & & GRID(ng) % on_p, & & GRID(ng) % on_r, & & GRID(ng) % pm, & & GRID(ng) % pmon_p, & & GRID(ng) % pmon_r, & & GRID(ng) % pn, & & GRID(ng) % pnom_p, & & GRID(ng) % pnom_r, & & MIXING(ng) % visc2_p, & & MIXING(ng) % visc2_r, & & OCEAN(ng) % u, & & OCEAN(ng) % v, & & COUPLING(ng) % tl_rufrc, & & COUPLING(ng) % tl_rvfrc, & & OCEAN(ng) % tl_u, & & OCEAN(ng) % tl_v) CALL wclock_off (ng, iTLM, 30, 106, MyFile) ! RETURN END SUBROUTINE tl_uv3dmix2 ! !*********************************************************************** SUBROUTINE tl_uv3dmix2_s_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nrhs, nnew, & & pmask, & & Hz, tl_Hz, & & om_p, om_r, on_p, on_r, & & pm, pmon_p, pmon_r, & & pn, pnom_p, pnom_r, & & visc2_p, visc2_r, & & u, v, & & tl_rufrc, tl_rvfrc, & & tl_u, tl_v) !*********************************************************************** ! 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, nnew real(r8), intent(in) :: pmask(LBi:,LBj:) real(r8), intent(in) :: Hz(LBi:,LBj:,:) real(r8), intent(in) :: tl_Hz(LBi:,LBj:,:) real(r8), intent(in) :: om_p(LBi:,LBj:) real(r8), intent(in) :: om_r(LBi:,LBj:) real(r8), intent(in) :: on_p(LBi:,LBj:) real(r8), intent(in) :: on_r(LBi:,LBj:) real(r8), intent(in) :: pm(LBi:,LBj:) real(r8), intent(in) :: pmon_p(LBi:,LBj:) real(r8), intent(in) :: pmon_r(LBi:,LBj:) real(r8), intent(in) :: pn(LBi:,LBj:) real(r8), intent(in) :: pnom_p(LBi:,LBj:) real(r8), intent(in) :: pnom_r(LBi:,LBj:) real(r8), intent(in) :: visc2_p(LBi:,LBj:) real(r8), intent(in) :: visc2_r(LBi:,LBj:) real(r8), intent(in) :: u(LBi:,LBj:,:,:) real(r8), intent(in) :: v(LBi:,LBj:,:,:) real(r8), intent(inout) :: tl_rufrc(LBi:,LBj:) real(r8), intent(inout) :: tl_rvfrc(LBi:,LBj:) real(r8), intent(inout) :: tl_u(LBi:,LBj:,:,:) real(r8), intent(inout) :: tl_v(LBi:,LBj:,:,:) ! ! Local variable declarations. ! integer :: i, j, k real(r8) :: cff, tl_cff, tl_cff1, tl_cff2 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_UFe real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_VFe real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_UFx real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_VFx ! !----------------------------------------------------------------------- ! 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 ! !----------------------------------------------------------------------- ! Compute horizontal harmonic viscosity along constant S-surfaces. !----------------------------------------------------------------------- ! K_LOOP : DO k=1,N(ng) ! ! Compute flux-components of the horizontal divergence of the stress ! tensor (m5/s2) in XI- and ETA-directions. ! DO j=JstrV-1,Jend DO i=IstrU-1,Iend !^ cff=visc2_r(i,j)*Hz(i,j,k)*0.5_r8* & !^ & (pmon_r(i,j)* & !^ & ((pn(i ,j)+pn(i+1,j))*u(i+1,j,k,nrhs)- & !^ & (pn(i-1,j)+pn(i ,j))*u(i ,j,k,nrhs))- & !^ & pnom_r(i,j)* & !^ & ((pm(i,j )+pm(i,j+1))*v(i,j+1,k,nrhs)- & !^ & (pm(i,j-1)+pm(i,j ))*v(i,j ,k,nrhs))) !^ tl_cff=visc2_r(i,j)*0.5_r8* & & (tl_Hz(i,j,k)* & & (pmon_r(i,j)* & & ((pn(i ,j)+pn(i+1,j))*u(i+1,j,k,nrhs)- & & (pn(i-1,j)+pn(i ,j))*u(i ,j,k,nrhs))- & & pnom_r(i,j)* & & ((pm(i,j )+pm(i,j+1))*v(i,j+1,k,nrhs)- & & (pm(i,j-1)+pm(i,j ))*v(i,j ,k,nrhs)))+ & & Hz(i,j,k)* & & (pmon_r(i,j)* & & ((pn(i ,j)+pn(i+1,j))*tl_u(i+1,j,k,nrhs)- & & (pn(i-1,j)+pn(i ,j))*tl_u(i ,j,k,nrhs))- & & pnom_r(i,j)* & & ((pm(i,j )+pm(i,j+1))*tl_v(i,j+1,k,nrhs)- & & (pm(i,j-1)+pm(i,j ))*tl_v(i,j ,k,nrhs)))) !^ UFx(i,j)=on_r(i,j)*on_r(i,j)*cff !^ tl_UFx(i,j)=on_r(i,j)*on_r(i,j)*tl_cff !^ VFe(i,j)=om_r(i,j)*om_r(i,j)*cff !^ tl_VFe(i,j)=om_r(i,j)*om_r(i,j)*tl_cff END DO END DO DO j=Jstr,Jend+1 DO i=Istr,Iend+1 !^ cff=visc2_p(i,j)*0.125_r8*(Hz(i-1,j ,k)+Hz(i,j ,k)+ & !^ & Hz(i-1,j-1,k)+Hz(i,j-1,k))* & !^ & (pmon_p(i,j)* & !^ & ((pn(i ,j-1)+pn(i ,j))*v(i ,j,k,nrhs)- & !^ & (pn(i-1,j-1)+pn(i-1,j))*v(i-1,j,k,nrhs))+ & !^ & pnom_p(i,j)* & !^ & ((pm(i-1,j )+pm(i,j ))*u(i,j ,k,nrhs)- & !^ & (pm(i-1,j-1)+pm(i,j-1))*u(i,j-1,k,nrhs))) !^ tl_cff=visc2_p(i,j)*0.125_r8* & & ((tl_Hz(i-1,j ,k)+tl_Hz(i,j ,k)+ & & tl_Hz(i-1,j-1,k)+tl_Hz(i,j-1,k))* & & (pmon_p(i,j)* & & ((pn(i ,j-1)+pn(i ,j))*v(i ,j,k,nrhs)- & & (pn(i-1,j-1)+pn(i-1,j))*v(i-1,j,k,nrhs))+ & & pnom_p(i,j)* & & ((pm(i-1,j )+pm(i,j ))*u(i,j ,k,nrhs)- & & (pm(i-1,j-1)+pm(i,j-1))*u(i,j-1,k,nrhs)))+ & & (Hz(i-1,j ,k)+Hz(i,j ,k)+ & & Hz(i-1,j-1,k)+Hz(i,j-1,k))* & & (pmon_p(i,j)* & & ((pn(i ,j-1)+pn(i ,j))*tl_v(i ,j,k,nrhs)- & & (pn(i-1,j-1)+pn(i-1,j))*tl_v(i-1,j,k,nrhs))+ & & pnom_p(i,j)* & & ((pm(i-1,j )+pm(i,j ))*tl_u(i,j ,k,nrhs)- & & (pm(i-1,j-1)+pm(i,j-1))*tl_u(i,j-1,k,nrhs)))) !^ cff=cff*pmask(i,j) !^ tl_cff=tl_cff*pmask(i,j) !^ UFe(i,j)=om_p(i,j)*om_p(i,j)*cff !^ tl_UFe(i,j)=om_p(i,j)*om_p(i,j)*tl_cff !^ VFx(i,j)=on_p(i,j)*on_p(i,j)*cff !^ tl_VFx(i,j)=on_p(i,j)*on_p(i,j)*tl_cff END DO END DO ! ! Time-step harmonic, S-surfaces viscosity term. Notice that momentum ! at this stage is HzU and HzV and has m2/s units. Add contribution for ! barotropic forcing terms. ! DO j=Jstr,Jend DO i=IstrU,Iend cff=0.25_r8*(pm(i-1,j)+pm(i,j))*(pn(i-1,j)+pn(i,j)) !^ cff1=0.5_r8*((pn(i-1,j)+pn(i,j))* & !^ & (UFx(i,j )-UFx(i-1,j))+ & !^ & (pm(i-1,j)+pm(i,j))* & !^ & (UFe(i,j+1)-UFe(i ,j))) !^ tl_cff1=0.5_r8*((pn(i-1,j)+pn(i,j))* & & (tl_UFx(i,j )-tl_UFx(i-1,j))+ & & (pm(i-1,j)+pm(i,j))* & & (tl_UFe(i,j+1)-tl_UFe(i ,j))) !^ cff2=dt(ng)*cff*cff1 !^ tl_cff2=dt(ng)*cff*tl_cff1 !^ rufrc(i,j)=rufrc(i,j)+cff1 !^ tl_rufrc(i,j)=tl_rufrc(i,j)+tl_cff1 !^ u(i,j,k,nnew)=u(i,j,k,nnew)+cff2 !^ tl_u(i,j,k,nnew)=tl_u(i,j,k,nnew)+tl_cff2 END DO END DO DO j=JstrV,Jend DO i=Istr,Iend cff=0.25_r8*(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1)) !^ cff1=0.5_r8*((pn(i,j-1)+pn(i,j))* & !^ & (VFx(i+1,j)-VFx(i,j ))- & !^ & (pm(i,j-1)+pm(i,j))* & !^ & (VFe(i ,j)-VFe(i,j-1))) !^ tl_cff1=0.5_r8*((pn(i,j-1)+pn(i,j))* & & (tl_VFx(i+1,j)-tl_VFx(i,j ))- & & (pm(i,j-1)+pm(i,j))* & & (tl_VFe(i ,j)-tl_VFe(i,j-1))) !^ cff2=dt(ng)*cff*cff1 !^ tl_cff2=dt(ng)*cff*tl_cff1 !^ rvfrc(i,j)=rvfrc(i,j)+cff1 !^ tl_rvfrc(i,j)=tl_rvfrc(i,j)+tl_cff1 !^ v(i,j,k,nnew)=v(i,j,k,nnew)+cff2 !^ tl_v(i,j,k,nnew)=tl_v(i,j,k,nnew)+tl_cff2 END DO END DO END DO K_LOOP ! RETURN END SUBROUTINE tl_uv3dmix2_s_tile END MODULE tl_uv3dmix2_mod