MODULE ad_t3dmix2_mod ! !git $Id$ !svn $Id: ad_t3dmix2_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 subroutine computes adjoint horizontal harmonic mixing of ! ! tracers along S-coordinate levels surfaces. ! ! ! ! BASIC STATE variables needed: t, Hz ! ! ! !======================================================================= ! implicit none ! PRIVATE PUBLIC ad_t3dmix2 ! CONTAINS ! !*********************************************************************** SUBROUTINE ad_t3dmix2 (ng, tile) !*********************************************************************** ! 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. ! character (len=*), parameter :: MyFile = & & "ROMS/Adjoint/ad_t3dmix2_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, iADM, 24, 53, MyFile) CALL ad_t3dmix2_s_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) % Hz, & & GRID(ng) % ad_Hz, & & GRID(ng) % pmon_u, & & GRID(ng) % pnom_v, & & GRID(ng) % pm, & & GRID(ng) % pn, & & MIXING(ng) % diff2, & & OCEAN(ng) % t, & & OCEAN(ng) % ad_t) CALL wclock_off (ng, iADM, 24, 87, MyFile) ! RETURN END SUBROUTINE ad_t3dmix2 ! !*********************************************************************** SUBROUTINE ad_t3dmix2_s_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nrhs, nstp, nnew, & & umask, vmask, & & Hz, ad_Hz, & & pmon_u, pnom_v, pm, pn, & & diff2, & & t, ad_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) :: Hz(LBi:,LBj:,:) real(r8), intent(in) :: pmon_u(LBi:,LBj:) real(r8), intent(in) :: pnom_v(LBi:,LBj:) real(r8), intent(in) :: pm(LBi:,LBj:) real(r8), intent(in) :: pn(LBi:,LBj:) real(r8), intent(in) :: t(LBi:,LBj:,:,:,:) real(r8), intent(inout) :: ad_Hz(LBi:,LBj:,:) real(r8), intent(inout) :: ad_t(LBi:,LBj:,:,:,:) ! ! Local variable declarations. ! integer :: i, itrc, j, k real(r8) :: cff, cff1 real(r8) :: adfac, adfac1, adfac2, adfac3, adfac4, ad_cff, ad_cff1 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_FE real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_FX ! !----------------------------------------------------------------------- ! 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 ! !----------------------------------------------------------------------- ! Initialize adjoint private variables. !----------------------------------------------------------------------- ! ad_cff=0.0_r8 ad_cff1=0.0_r8 ad_FE(IminS:ImaxS,JminS:JmaxS)=0.0_r8 ad_FX(IminS:ImaxS,JminS:JmaxS)=0.0_r8 ! !---------------------------------------------------------------------- ! Compute adjoint horizontal harmonic diffusion along constant ! S-surfaces. !---------------------------------------------------------------------- ! DO itrc=1,NT(ng) DO k=1,N(ng) ! ! Time-step harmonic, S-surfaces diffusion term (m Tunits). ! DO j=Jstr,Jend DO i=Istr,Iend !^ tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)+tl_cff !^ ad_cff=ad_cff+ad_t(i,j,k,nnew,itrc) !^ tl_cff=dt(ng)*pm(i,j)*pn(i,j)* & !^ & (tl_FX(i+1,j)-tl_FX(i,j)+ & !^ & tl_FE(i,j+1)-tl_FE(i,j)) !^ adfac=dt(ng)*pm(i,j)*pn(i,j)*ad_cff ad_FX(i ,j)=ad_FX(i ,j)-adfac ad_FX(i+1,j)=ad_FX(i+1,j)+adfac ad_FE(i,j )=ad_FE(i,j )-adfac ad_FE(i,j+1)=ad_FE(i,j+1)+adfac ad_cff=0.0_r8 END DO END DO ! ! Compute XI- and ETA-components of diffusive tracer flux (T m3/s). ! DO j=Jstr,Jend+1 DO i=Istr,Iend cff=0.25_r8*(diff2(i,j,itrc)+diff2(i,j-1,itrc))* & & pnom_v(i,j) !^ tl_FE(i,j)=tl_FE(i,j)*vmask(i,j) !^ ad_FE(i,j)=ad_FE(i,j)*vmask(i,j) !^ tl_FE(i,j)=cff* & !^ & ((tl_Hz(i,j,k)+tl_Hz(i,j-1,k))* & !^ & (t(i,j ,k,nrhs,itrc)- & !^ & t(i,j-1,k,nrhs,itrc))+ & !^ & (Hz(i,j,k)+Hz(i,j-1,k))* & !^ & (tl_t(i,j ,k,nrhs,itrc)- & !^ & tl_t(i,j-1,k,nrhs,itrc))) !^ adfac=cff*ad_FE(i,j) adfac1=adfac*(t(i,j,k,nrhs,itrc)-t(i,j-1,k,nrhs,itrc)) adfac2=adfac*(Hz(i,j,k)+Hz(i,j-1,k)) ad_Hz(i,j-1,k)=ad_Hz(i,j-1,k)+adfac1 ad_Hz(i,j ,k)=ad_Hz(i,j ,k)+adfac1 ad_t(i,j-1,k,nrhs,itrc)=ad_t(i,j-1,k,nrhs,itrc)-adfac2 ad_t(i,j ,k,nrhs,itrc)=ad_t(i,j ,k,nrhs,itrc)+adfac2 ad_FE(i,j)=0.0_r8 END DO END DO DO j=Jstr,Jend DO i=Istr,Iend+1 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i-1,j,itrc))* & & pmon_u(i,j) !^ tl_FX(i,j)=tl_FX(i,j)*umask(i,j) !^ ad_FX(i,j)=ad_FX(i,j)*umask(i,j) !^ tl_FX(i,j)=cff* & !^ & ((tl_Hz(i,j,k)+tl_Hz(i-1,j,k))* & !^ & (t(i ,j,k,nrhs,itrc)- & !^ & t(i-1,j,k,nrhs,itrc))+ & !^ & (Hz(i,j,k)+Hz(i-1,j,k))* & !^ & (tl_t(i ,j,k,nrhs,itrc)- & !^ & tl_t(i-1,j,k,nrhs,itrc))) !^ adfac=cff*ad_FX(i,j) adfac1=adfac*(t(i,j,k,nrhs,itrc)-t(i-1,j,k,nrhs,itrc)) adfac2=adfac*(Hz(i,j,k)+Hz(i-1,j,k)) ad_Hz(i-1,j,k)=ad_Hz(i-1,j,k)+adfac1 ad_Hz(i ,j,k)=ad_Hz(i ,j,k)+adfac1 ad_t(i-1,j,k,nrhs,itrc)=ad_t(i-1,j,k,nrhs,itrc)-adfac2 ad_t(i ,j,k,nrhs,itrc)=ad_t(i ,j,k,nrhs,itrc)+adfac2 ad_FX(i,j)=0.0_r8 END DO END DO END DO END DO ! RETURN END SUBROUTINE ad_t3dmix2_s_tile END MODULE ad_t3dmix2_mod