MODULE ad_set_depth_mod ! !git $Id$ !svn $Id: ad_set_depth.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 routine computes the time evolving depths of the model grid ! ! and its associated vertical transformation metric (thickness). ! ! ! ! Currently, two vertical coordinate transformations are available ! ! with various possible vertical stretching, C(s), functions, (see ! ! routine "set_scoord.F" for details). ! ! ! ! BASIC STATE variables needed: NONE ! ! Independent Variables: ad_Hz, ad_z_r, ad_z_w ! ! ! !======================================================================= ! implicit none ! PRIVATE PUBLIC :: ad_set_depth, ad_set_depth_tile ! CONTAINS ! !*********************************************************************** SUBROUTINE ad_set_depth (ng, tile, model) !*********************************************************************** ! USE mod_param USE mod_coupling USE mod_grid USE mod_ocean USE mod_stepping ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model ! ! Local variable declarations. ! character (len=*), parameter :: MyFile = & & "ROMS/Adjoint/ad_set_depth.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 ! CALL wclock_on (ng, model, 12, 57, MyFile) CALL ad_set_depth_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nstp(ng), nnew(ng), & & GRID(ng) % h, & & GRID(ng) % ad_h, & & COUPLING(ng) % Zt_avg1, & & COUPLING(ng) % ad_Zt_avg1, & & GRID(ng) % ad_Hz, & & GRID(ng) % ad_z_r, & & GRID(ng) % ad_z_w) CALL wclock_off (ng, model, 12, 74, MyFile) ! RETURN END SUBROUTINE ad_set_depth ! !*********************************************************************** SUBROUTINE ad_set_depth_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nstp, nnew, & & h, ad_h, & & Zt_avg1, ad_Zt_avg1, & & ad_Hz, ad_z_r, ad_z_w) !*********************************************************************** ! USE mod_param USE mod_scalars ! USE ad_exchange_2d_mod USE ad_exchange_3d_mod USE mp_exchange_mod, ONLY : ad_mp_exchange2d, ad_mp_exchange3d ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model integer, intent(in) :: LBi, UBi, LBj, UBj integer, intent(in) :: IminS, ImaxS, JminS, JmaxS integer, intent(in) :: nstp, nnew ! real(r8), intent(in) :: h(LBi:,LBj:) real(r8), intent(in) :: Zt_avg1(LBi:,LBj:) real(r8), intent(inout) :: ad_h(LBi:,LBj:) real(r8), intent(inout) :: ad_Zt_avg1(LBi:,LBj:) real(r8), intent(inout) :: ad_Hz(LBi:,LBj:,:) real(r8), intent(inout) :: ad_z_r(LBi:,LBj:,:) real(r8), intent(inout) :: ad_z_w(LBi:,LBj:,0:) ! ! Local variable declarations. ! integer :: i, j, k real(r8) :: cff, cff_r, cff1_r, cff2_r, cff_w, cff1_w, cff2_w real(r8) :: hinv, hwater, z_r0, z_w0 real(r8) :: adfac, ad_hinv, ad_hwater, ad_z_r0, ad_z_w0 real(r8) :: ad_cff2_r, ad_cff2_w ! !----------------------------------------------------------------------- ! 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_cff2_r=0.0_r8 ad_cff2_w=0.0_r8 ad_z_r0=0.0_r8 ad_z_w0=0.0_r8 ad_hinv=0.0_r8 ad_hwater=0.0_r8 ! !----------------------------------------------------------------------- ! Compute time evolving adjoint depths and vertical thicknesses. !----------------------------------------------------------------------- ! ! Exchange boundary information. ! !^ CALL mp_exchange3d (ng, tile, model, 2, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), & !^ & NghostPoints, & !^ & EWperiodic(ng), NSperiodic(ng), & !^ & tl_z_r, tl_Hz) !^ CALL ad_mp_exchange3d (ng, tile, model, 2, & & LBi, UBi, LBj, UBj, 1, N(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_z_r, ad_Hz) !^ CALL mp_exchange3d (ng, tile, model, 1, & !^ & LBi, UBi, LBj, UBj, 0, N(ng), & !^ & NghostPoints, & !^ & EWperiodic(ng), NSperiodic(ng), & !^ & tl_z_w) !^ CALL ad_mp_exchange3d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, 0, N(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_z_w) !^ CALL mp_exchange2d (ng, tile, model, 1, & !^ & LBi, UBi, LBj, UBj, & !^ & NghostPoints, & !^ & EWperiodic(ng), NSperiodic(ng), & !^ & tl_h) !^ CALL ad_mp_exchange2d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_h) ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN !^ CALL exchange_r3d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), & !^ & tl_Hz) !^ CALL ad_exchange_r3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_Hz) !^ CALL exchange_r3d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), & !^ & tl_z_r) !^ CALL ad_exchange_r3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_z_r) !^ CALL exchange_w3d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, 0, N(ng), & !^ & tl_z_w) !^ CALL ad_exchange_w3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 0, N(ng), & & ad_z_w) !^ CALL exchange_r2d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & tl_h) !^ CALL ad_exchange_r2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_h) END IF ! !----------------------------------------------------------------------- ! Original formulation: Compute vertical depths (meters, negative) at ! RHO- and W-points, and vertical grid ! thicknesses. Various stretching functions are possible. ! ! z_w(x,y,s,t) = Zo_w + zeta(x,y,t) * [1.0 + Zo_w / h(x,y)] ! ! Zo_w = hc * [s(k) - C(k)] + C(k) * h(x,y) ! !----------------------------------------------------------------------- ! IF (Vtransform(ng).eq.1) THEN DO j=JstrT,JendT DO k=N(ng),1,-1 cff_w=hc(ng)*(SCALARS(ng)%sc_w(k)-SCALARS(ng)%Cs_w(k)) cff_r=hc(ng)*(SCALARS(ng)%sc_r(k)-SCALARS(ng)%Cs_r(k)) cff1_r=SCALARS(ng)%Cs_r(k) cff1_w=SCALARS(ng)%Cs_w(k) DO i=IstrT,IendT hwater=h(i,j) hinv=1.0_r8/hwater z_w0=cff_w+cff1_w*hwater z_r0=cff_r+cff1_r*hwater !^ tl_Hz(i,j,k)=tl_z_w(i,j,k)-tl_z_w(i,j,k-1) !^ ad_z_w(i,j,k )=ad_z_w(i,j,k )+ad_Hz(i,j,k) ad_z_w(i,j,k-1)=ad_z_w(i,j,k-1)-ad_Hz(i,j,k) ad_Hz(i,j,k)=0.0_r8 !^ tl_z_r(i,j,k)=tl_z_r0+ & !^ & tl_Zt_avg1(i,j)*(1.0_r8+z_r0*hinv)+ & !^ & Zt_avg1(i,j)*(tl_z_r0*hinv+z_r0*tl_hinv) !^ adfac=Zt_avg1(i,j)*ad_z_r(i,j,k) ad_z_r0=ad_z_r0+hinv*adfac+ad_z_r(i,j,k) ad_hinv=ad_hinv+z_r0*adfac ad_Zt_avg1(i,j)=ad_Zt_avg1(i,j)+ & & (1.0_r8+z_r0*hinv)*ad_z_r(i,j,k) ad_z_r(i,j,k)=0.0_r8 !^ tl_z_r0=cff1_r*tl_hwater !^ ad_hwater=ad_hwater+cff1_r*ad_z_r0 ad_z_r0=0.0_r8 !^ tl_z_w(i,j,k)=tl_z_w0+ & !^ & tl_Zt_avg1(i,j)*(1.0_r8+z_w0*hinv)+ & !^ & Zt_avg1(i,j)*(tl_z_w0*hinv+z_w0*tl_hinv) !^ adfac=Zt_avg1(i,j)*ad_z_w(i,j,k) ad_z_w0=ad_z_w0+hinv*adfac+ad_z_w(i,j,k) ad_hinv=ad_hinv+z_w0*adfac ad_Zt_avg1(i,j)=ad_Zt_avg1(i,j)+ & & (1.0_r8+z_w0*hinv)*ad_z_w(i,j,k) ad_z_w(i,j,k)=0.0_r8 !^ tl_z_w0=cff1_w*tl_hwater !^ ad_hwater=ad_hwater+cff1_w*ad_z_w0 ad_z_w0=0.0_r8 !^ tl_hinv=-hinv*hinv*tl_hwater !^ ad_hwater=ad_hwater-hinv*hinv*ad_hinv ad_hinv=0.0_r8 !^ tl_hwater=tl_h(i,j) !^ ad_h(i,j)=ad_h(i,j)+ad_hwater ad_hwater=0.0_r8 END DO END DO DO i=IstrT,IendT !^ tl_z_w(i,j,0)=-tl_h(i,j) !^ ad_h(i,j)=ad_h(i,j)-ad_z_w(i,j,0) ad_z_w(i,j,0)=0.0_r8 END DO END DO ! !----------------------------------------------------------------------- ! New formulation: Compute vertical depths (meters, negative) at ! RHO- and W-points, and vertical grid thicknesses. ! Various stretching functions are possible. ! ! z_w(x,y,s,t) = zeta(x,y,t) + [zeta(x,y,t)+ h(x,y)] * Zo_w ! ! Zo_w = [hc * s(k) + C(k) * h(x,y)] / [hc + h(x,y)] ! !----------------------------------------------------------------------- ! ELSE IF (Vtransform(ng).eq.2) THEN DO j=JstrT,JendT DO k=N(ng),1,-1 cff_r=hc(ng)*SCALARS(ng)%sc_r(k) cff_w=hc(ng)*SCALARS(ng)%sc_w(k) cff1_r=SCALARS(ng)%Cs_r(k) cff1_w=SCALARS(ng)%Cs_w(k) DO i=IstrT,IendT hwater=h(i,j) hinv=1.0_r8/(hc(ng)+hwater) cff2_r=(cff_r+cff1_r*hwater)*hinv cff2_w=(cff_w+cff1_w*hwater)*hinv !^ tl_Hz(i,j,k)=tl_z_w(i,j,k)-tl_z_w(i,j,k-1) !^ ad_z_w(i,j,k )=ad_z_w(i,j,k )+ad_Hz(i,j,k) ad_z_w(i,j,k-1)=ad_z_w(i,j,k-1)-ad_Hz(i,j,k) ad_Hz(i,j,k)=0.0_r8 !^ tl_z_r(i,j,k)=tl_Zt_avg1(i,j)+ & !^ & (tl_Zt_avg1(i,j)+tl_hwater)*cff2_r+ & !^ & (Zt_avg1(i,j)+hwater)*tl_cff2_r !^ adfac=cff2_r*ad_z_r(i,j,k) ad_cff2_r=ad_cff2_r+ & & (Zt_avg1(i,j)+hwater)*ad_z_r(i,j,k) ad_hwater=ad_hwater+adfac ad_Zt_avg1(i,j)=ad_Zt_avg1(i,j)+ad_z_r(i,j,k)+adfac ad_z_r(i,j,k)=0.0_r8 !^ tl_cff2_r=cff1_r*tl_hwater*hinv+ & !^ & (cff_r+cff1_r*hwater)*tl_hinv !^ ad_hinv=ad_hinv+ & & (cff_r+cff1_r*hwater)*ad_cff2_r ad_hwater=ad_hwater+ & & cff1_r*hinv*ad_cff2_r ad_cff2_r=0.0_r8 !^ tl_z_w(i,j,k)=tl_Zt_avg1(i,j)+ & !^ & (tl_Zt_avg1(i,j)+tl_hwater)*cff2_w+ & !^ & (Zt_avg1(i,j)+hwater)*tl_cff2_w !^ adfac=cff2_w*ad_z_w(i,j,k) ad_cff2_w=ad_cff2_w+ & & (Zt_avg1(i,j)+hwater)*ad_z_w(i,j,k) ad_hwater=ad_hwater+adfac ad_Zt_avg1(i,j)=ad_Zt_avg1(i,j)+ad_z_w(i,j,k)+adfac ad_z_w(i,j,k)=0.0_r8 !^ tl_cff2_w=cff1_w*tl_hwater*hinv+ & !^ & (cff_w+cff1_w*hwater)*tl_hinv !^ ad_hinv=ad_hinv+ & & (cff_w+cff1_w*hwater)*ad_cff2_w ad_hwater=ad_hwater+ & & cff1_w*hinv*ad_cff2_w ad_cff2_w=0.0_r8 !^ tl_hinv=-hinv*hinv*tl_hwater !^ ad_hwater=ad_hwater-hinv*hinv*ad_hinv ad_hinv=0.0_r8 !^ tl_hwater=tl_h(i,j) !^ ad_h(i,j)=ad_h(i,j)+ad_hwater ad_hwater=0.0_r8 END DO END DO DO i=IstrT,IendT !^ tl_z_w(i,j,0)=-tl_h(i,j) !^ ad_h(i,j)=ad_h(i,j)-ad_z_w(i,j,0) ad_z_w(i,j,0)=0.0_r8 END DO END DO END IF RETURN END SUBROUTINE ad_set_depth_tile END MODULE ad_set_depth_mod