MODULE ad_step3d_uv_mod ! !git $Id$ !svn $Id: ad_step3d_uv.F 1183 2023-07-27 04:20:57Z 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 time-steps the adjoint horizontal momentum ! ! equations. The vertical viscosity terms are time-stepped using ! ! an implicit algorithm. ! ! ! ! BASIC STATE variables needed: u, v, ru, rv, Akv, ! ! Hz, Huon, Hvom, z_r, z_w, ! ! DU_avg1, DU_avg2, DV_avg1, DV_avg2, ! ! Qsrc ! ! ! !======================================================================= ! implicit none ! PRIVATE PUBLIC :: ad_step3d_uv ! CONTAINS ! !*********************************************************************** SUBROUTINE ad_step3d_uv (ng, tile) !*********************************************************************** ! USE mod_param USE mod_coupling USE mod_forces 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_step3d_uv.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, iADM, 34, 59, MyFile) CALL ad_step3d_uv_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & knew(ng), nrhs(ng), nstp(ng), nnew(ng), & & GRID(ng) % umask, & & GRID(ng) % vmask, & & GRID(ng) % om_v, & & GRID(ng) % on_u, & & GRID(ng) % pm, & & GRID(ng) % pn, & & GRID(ng) % Hz, & & GRID(ng) % ad_Hz, & & GRID(ng) % z_r, & & GRID(ng) % ad_z_r, & & GRID(ng) % z_w, & & GRID(ng) % ad_z_w, & & MIXING(ng) % Akv, & & MIXING(ng) % ad_Akv, & & COUPLING(ng) % DU_avg1, & & COUPLING(ng) % ad_DU_avg1, & & COUPLING(ng) % DV_avg1, & & COUPLING(ng) % ad_DV_avg1, & & COUPLING(ng) % DU_avg2, & & COUPLING(ng) % ad_DU_avg2, & & COUPLING(ng) % DV_avg2, & & COUPLING(ng) % ad_DV_avg2, & & OCEAN(ng) % ad_ru, & & OCEAN(ng) % ad_rv, & & OCEAN(ng) % u, & & OCEAN(ng) % ad_u, & & OCEAN(ng) % v, & & OCEAN(ng) % ad_v, & & OCEAN(ng) % ad_ubar, & & OCEAN(ng) % ad_vbar, & & OCEAN(ng) % ad_ubar_sol, & & OCEAN(ng) % ad_vbar_sol, & & GRID(ng) % Huon, & & GRID(ng) % ad_Huon, & & GRID(ng) % Hvom, & & GRID(ng) % ad_Hvom) CALL wclock_off (ng, iADM, 34, 132, MyFile) ! RETURN END SUBROUTINE ad_step3d_uv ! !*********************************************************************** SUBROUTINE ad_step3d_uv_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & knew, nrhs, nstp, nnew, & & umask, vmask, & & om_v, on_u, pm, pn, & & Hz, ad_Hz, & & z_r, ad_z_r, & & z_w, ad_z_w, & & Akv, ad_Akv, & & DU_avg1, ad_DU_avg1, & & DV_avg1, ad_DV_avg1, & & DU_avg2, ad_DU_avg2, & & DV_avg2, ad_DV_avg2, & & ad_ru, ad_rv, & & u, ad_u, & & v, ad_v, & & ad_ubar, ad_vbar, & & ad_ubar_sol, ad_vbar_sol, & & Huon, ad_Huon, & & Hvom, ad_Hvom) !*********************************************************************** ! USE mod_param USE mod_scalars USE mod_sources ! USE ad_exchange_2d_mod USE ad_exchange_3d_mod USE exchange_3d_mod USE mp_exchange_mod, ONLY : ad_mp_exchange2d, ad_mp_exchange3d USE ad_u3dbc_mod, ONLY : ad_u3dbc_tile USE ad_v3dbc_mod, ONLY : ad_v3dbc_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) :: knew, nrhs, nstp, nnew ! real(r8), intent(in) :: umask(LBi:,LBj:) real(r8), intent(in) :: vmask(LBi:,LBj:) real(r8), intent(in) :: om_v(LBi:,LBj:) real(r8), intent(in) :: on_u(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) :: z_r(LBi:,LBj:,:) real(r8), intent(in) :: z_w(LBi:,LBj:,0:) real(r8), intent(in) :: Akv(LBi:,LBj:,0:) real(r8), intent(in) :: DU_avg1(LBi:,LBj:) real(r8), intent(in) :: DV_avg1(LBi:,LBj:) real(r8), intent(in) :: DU_avg2(LBi:,LBj:) real(r8), intent(in) :: DV_avg2(LBi:,LBj:) real(r8), intent(in) :: u(LBi:,LBj:,:,:) real(r8), intent(in) :: v(LBi:,LBj:,:,:) real(r8), intent(inout) :: Huon(LBi:,LBj:,:) real(r8), intent(inout) :: Hvom(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:) real(r8), intent(inout) :: ad_Akv(LBi:,LBj:,0:) real(r8), intent(inout) :: ad_DU_avg1(LBi:,LBj:) real(r8), intent(inout) :: ad_DV_avg1(LBi:,LBj:) real(r8), intent(inout) :: ad_DU_avg2(LBi:,LBj:) real(r8), intent(inout) :: ad_DV_avg2(LBi:,LBj:) real(r8), intent(inout) :: ad_ru(LBi:,LBj:,0:,:) real(r8), intent(inout) :: ad_rv(LBi:,LBj:,0:,:) real(r8), intent(inout) :: ad_u(LBi:,LBj:,:,:) real(r8), intent(inout) :: ad_v(LBi:,LBj:,:,:) real(r8), intent(inout) :: ad_ubar(LBi:,LBj:,:) real(r8), intent(inout) :: ad_vbar(LBi:,LBj:,:) real(r8), intent(inout) :: ad_Huon(LBi:,LBj:,:) real(r8), intent(inout) :: ad_Hvom(LBi:,LBj:,:) real(r8), intent(out) :: ad_ubar_sol(LBi:,LBj:) real(r8), intent(out) :: ad_vbar_sol(LBi:,LBj:) ! ! Local variable declarations. ! integer :: i, idiag, is, j, k ! real(r8) :: cff, cff1, cff2 real(r8) :: ad_cff, ad_cff1, ad_cff2 real(r8) :: adfac, adfac1, adfac2 ! real(r8), dimension(IminS:ImaxS) :: CF1 real(r8), dimension(IminS:ImaxS) :: FC1 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: AK real(r8), dimension(IminS:ImaxS,0:N(ng)) :: BC 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)) :: DC1 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: FC real(r8), dimension(IminS:ImaxS,N(ng)) :: Hzk real(r8), dimension(IminS:ImaxS,N(ng)) :: oHz real(r8), dimension(IminS:ImaxS,0:N(ng)) :: ad_AK real(r8), dimension(IminS:ImaxS,0:N(ng)) :: ad_BC real(r8), dimension(IminS:ImaxS,0:N(ng)) :: ad_CF real(r8), dimension(IminS:ImaxS,0:N(ng)) :: ad_DC real(r8), dimension(IminS:ImaxS,0:N(ng)) :: ad_FC real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_Hzk real(r8), dimension(IminS:ImaxS,N(ng)) :: ad_oHz ! ! !----------------------------------------------------------------------- ! 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 DO k=0,N(ng) DO i=IminS,ImaxS ad_AK(i,k)=0.0_r8 ad_BC(i,k)=0.0_r8 ad_CF(i,k)=0.0_r8 ad_DC(i,k)=0.0_r8 ad_FC(i,k)=0.0_r8 END DO END DO DO k=1,N(ng) DO i=IminS,ImaxS ad_Hzk(i,k)=0.0_r8 ad_oHz(i,k)=0.0_r8 END DO END DO ! !----------------------------------------------------------------------- ! Exchange boundary data. !----------------------------------------------------------------------- ! !^ CALL mp_exchange2d (ng, tile, iTLM, 4, & !^ & LBi, UBi, LBj, UBj, & !^ & NghostPoints, & !^ & EWperiodic(ng), NSperiodic(ng), & !^ & tl_ubar(:,:,1), tl_vbar(:,:,1), & !^ & tl_ubar(:,:,2), tl_vbar(:,:,2)) !^ CALL ad_mp_exchange2d (ng, tile, iADM, 4, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_ubar(:,:,1), ad_vbar(:,:,1), & & ad_ubar(:,:,2), ad_vbar(:,:,2)) !^ CALL mp_exchange3d (ng, tile, iTLM, 4, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), & !^ & NghostPoints, & !^ & EWperiodic(ng), NSperiodic(ng), & !^ & tl_u(:,:,:,nnew), tl_v(:,:,:,nnew), & !^ & tl_Huon, tl_Hvom) !^ CALL ad_mp_exchange3d (ng, tile, iADM, 4, & & LBi, UBi, LBj, UBj, 1, N(ng), & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & ad_u(:,:,:,nnew), ad_v(:,:,:,nnew), & & ad_Huon, ad_Hvom) ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN DO k=1,2 !^ CALL exchange_v2d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & tl_vbar(:,:,k)) !^ CALL ad_exchange_v2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_vbar(:,:,k)) !^ CALL exchange_u2d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, & !^ & tl_ubar(:,:,k)) !^ CALL ad_exchange_u2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & ad_ubar(:,:,k)) END DO ! CALL exchange_v3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & Hvom) !^ CALL exchange_v3d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), & !^ & tl_Hvom) !^ CALL ad_exchange_v3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_Hvom) CALL exchange_u3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & Huon) !^ CALL exchange_u3d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), & !^ & tl_Huon) !^ CALL ad_exchange_u3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_Huon) !^ CALL exchange_v3d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), & !^ & tl_v(:,:,:,nnew)) !^ CALL ad_exchange_v3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_v(:,:,:,nnew)) !^ CALL exchange_u3d_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, 1, N(ng), & !^ & tl_u(:,:,:,nnew)) !^ CALL ad_exchange_u3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, N(ng), & & ad_u(:,:,:,nnew)) END IF ! !----------------------------------------------------------------------- ! Couple 2D and 3D momentum equations. !----------------------------------------------------------------------- ! ! Couple adjoint velocity component in the ETA-direction. ! J_LOOP1 : DO j=JstrT,JendT IF (j.ge.Jstr) THEN ! ! First, compute BASIC STATE quantities. This section was computed ! three times in the original code. This can avoided by saving the ! intermediate values in scratch arrays DC1, CF1, and FC1 (HGA). ! DO i=IstrT,IendT DC(i,0)=0.0_r8 CF(i,0)=0.0_r8 FC(i,0)=0.0_r8 END DO DO k=1,N(ng) DO i=IstrT,IendT cff=0.5_r8*om_v(i,j) DC(i,k)=cff*(Hz(i,j,k)+Hz(i,j-1,k)) DC(i,0)=DC(i,0)+DC(i,k) CF(i,0)=CF(i,0)+ & & DC(i,k)*v(i,j,k,nnew) END DO END DO DO i=IstrT,IendT DC1(i,0)=DC(i,0) ! intermediate DC(i,0)=1.0_r8/DC(i,0) ! recursive CF1(i)=CF(i,0) ! intermediate CF(i,0)=DC(i,0)*(CF(i,0)-DV_avg1(i,j)) ! recursive END DO DO k=N(ng),1,-1 DO i=IstrT,IendT !^ Hvom(i,j,k)=0.5_r8*(Hvom(i,j,k)+v(i,j,k,nnew)*DC(i,k)) !^ FC(i,0)=FC(i,0)+Hvom(i,j,k) !^ FC(i,0)=FC(i,0)+ & & 0.5_r8*(Hvom(i,j,k)+v(i,j,k,nnew)*DC(i,k)) END DO END DO DO i=IstrT,IendT FC1(i)=FC(i,0) ! intermediate FC(i,0)=DC(i,0)*(FC(i,0)-DV_avg2(i,j)) ! recursive END DO ! ! Compute correct mass flux, Hz*v/m. ! DO k=1,N(ng) DO i=IstrT,IendT !^ tl_Hvom(i,j,k)=tl_Hvom(i,j,k)- & !^ & tl_DC(i,k)*FC(i,0)- & !^ & DC(i,k)*tl_FC(i,0) !^ ad_DC(i,k)=ad_DC(i,k)-ad_Hvom(i,j,k)*FC(i,0) ad_FC(i,0)=ad_FC(i,0)-ad_Hvom(i,j,k)*DC(i,k) END DO END DO DO i=IstrT,IendT !^ tl_FC(i,0)=tl_DC(i,0)*(FC1(i)-DV_avg2(i,j))+ & !^ & DC(i,0)*(tl_FC(i,0)-tl_DV_avg2(i,j)) !^ adfac=DC(i,0)*ad_FC(i,0) ad_DC(i,0)=ad_DC(i,0)+ad_FC(i,0)*(FC1(i)-DV_avg2(i,j)) ad_DV_avg2(i,j)=ad_DV_avg2(i,j)-adfac ad_FC(i,0)=adfac END DO !^ DO k=N(ng),1,-1 !^ DO k=1,N(ng) DO i=IstrT,IendT !^ tl_FC(i,0)=tl_FC(i,0)+tl_Hvom(i,j,k) !^ ad_Hvom(i,j,k)=ad_Hvom(i,j,k)+ad_FC(i,0) !^ tl_Hvom(i,j,k)=0.5_r8*(tl_Hvom(i,j,k)+ & !^ & tl_v(i,j,k,nnew)*DC(i,k)+ & !^ & v(i,j,k,nnew)*tl_DC(i,k)) !^ adfac=0.5_r8*ad_Hvom(i,j,k) ad_DC(i,k)=ad_DC(i,k)+v(i,j,k,nnew)*adfac ad_v(i,j,k,nnew)=ad_v(i,j,k,nnew)+DC(i,k)*adfac ad_Hvom(i,j,k)=adfac END DO END DO ! ! Replace only BOUNDARY POINTS incorrect vertical mean with more ! accurate barotropic component, vbar=DV_avg1/(D*om_v). Recall that, ! D=CF(:,0). The replacement is avoided if the boundary edge is ! periodic. The J-loop is pipelined, so we need to do a special ! test for the southern and northern domain edges. ! IF (.not.(CompositeGrid(inorth,ng).or.NSperiodic(ng))) THEN IF (j.eq.Mm(ng)+1) THEN ! northern boundary DO k=1,N(ng) ! J-loop pipelined DO i=Istr,Iend !^ tl_v(i,j,k,nnew)=tl_v(i,j,k,nnew)* & !^ & vmask(i,j) !^ ad_v(i,j,k,nnew)=ad_v(i,j,k,nnew)* & & vmask(i,j) !^ tl_v(i,j,k,nnew)=tl_v(i,j,k,nnew)- & !^ & tl_CF(i,0) !^ ad_CF(i,0)=ad_CF(i,0)-ad_v(i,j,k,nnew) END DO END DO END IF END IF ! IF (.not.(CompositeGrid(isouth,ng).or.NSperiodic(ng))) THEN IF (j.eq.1) THEN ! southern boundary DO k=1,N(ng) ! J-loop pipelined DO i=Istr,Iend !^ tl_v(i,j,k,nnew)=tl_v(i,j,k,nnew)* & !^ & vmask(i,j) !^ ad_v(i,j,k,nnew)=ad_v(i,j,k,nnew)* & & vmask(i,j) !^ tl_v(i,j,k,nnew)=tl_v(i,j,k,nnew)- & !^ & tl_CF(i,0) !^ ad_CF(i,0)=ad_CF(i,0)-ad_v(i,j,k,nnew) END DO END DO END IF END IF ! IF (.not.(CompositeGrid(ieast,ng).or.EWperiodic(ng))) THEN IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN DO k=1,N(ng) !^ tl_v(Iend+1,j,k,nnew)=tl_v(Iend+1,j,k,nnew)* & !^ & vmask(Iend+1,j) !^ ad_v(Iend+1,j,k,nnew)=ad_v(Iend+1,j,k,nnew)* & & vmask(Iend+1,j) !^ tl_v(Iend+1,j,k,nnew)=tl_v(Iend+1,j,k,nnew)- & !^ & tl_CF(Iend+1,0) !^ ad_CF(Iend+1,0)=ad_CF(Iend+1,0)- & & ad_v(Iend+1,j,k,nnew) END DO END IF END IF ! IF (.not.(CompositeGrid(iwest,ng).or.EWperiodic(ng))) THEN IF (DOMAIN(ng)%Western_Edge(tile)) THEN DO k=1,N(ng) !^ tl_v(Istr-1,j,k,nnew)=tl_v(Istr-1,j,k,nnew)* & !^ & vmask(Istr-1,j) !^ ad_v(Istr-1,j,k,nnew)=ad_v(Istr-1,j,k,nnew)* & & vmask(Istr-1,j) !^ tl_v(Istr-1,j,k,nnew)=tl_v(Istr-1,j,k,nnew)- & !^ & tl_CF(Istr-1,0) !^ ad_CF(Istr-1,0)=ad_CF(Istr-1,0)- & & ad_v(Istr-1,j,k,nnew) END DO END IF END IF ! ! Save adjoint solution for time-step iic(ng)-1 as the sum of time ! records 1 and 2. ! DO i=IstrT,IendT ad_vbar_sol(i,j)=ad_vbar(i,j,1)+ad_vbar(i,j,2) END DO ! ! Compute adjoint thicknesses of V-boxes DC(i,1:N), total depth of the ! water column DC(i,0), and incorrect vertical mean CF(i,0). Notice ! that barotropic component is replaced with its fast-time averaged ! values. ! DO i=IstrT,IendT !^ tl_vbar(i,j,2)=tl_vbar(i,j,1) !^ ad_vbar(i,j,1)=ad_vbar(i,j,1)+ad_vbar(i,j,2) ad_vbar(i,j,2)=0.0_r8 !^ tl_vbar(i,j,1)=tl_DC(i,0)*DV_avg1(i,j)+ & !^ & DC(i,0)*tl_DV_avg1(i,j) !^ ad_DC(i,0)=ad_DC(i,0)+ad_vbar(i,j,1)*DV_avg1(i,j) ad_DV_avg1(i,j)=ad_DV_avg1(i,j)+ad_vbar(i,j,1)*DC(i,0) ad_vbar(i,j,1)=0.0_r8 !^ tl_CF(i,0)=tl_DC(i,0)*(CF1(i)-DV_avg1(i,j))+ & !^ & DC(i,0)*(tl_CF(i,0)-tl_DV_avg1(i,j)) !^ adfac=DC(i,0)*ad_CF(i,0) ad_DC(i,0)=ad_DC(i,0)+ad_CF(i,0)*(CF1(i)-DV_avg1(i,j)) ad_DV_avg1(i,j)=ad_DV_avg1(i,j)-adfac ad_CF(i,0)=adfac !^ tl_DC(i,0)=-tl_DC(i,0/(DC1(i,0)*DC1(i,0)) !^ ad_DC(i,0)=-ad_DC(i,0)/(DC1(i,0)*DC1(i,0)) END DO DO i=IstrT,IendT DC(i,0)=0.0_r8 CF(i,0)=0.0_r8 FC(i,0)=0.0_r8 END DO DO k=1,N(ng) DO i=IstrT,IendT cff=0.5_r8*om_v(i,j) DC(i,k)=cff*(Hz(i,j,k)+Hz(i,j-1,k)) DC(i,0)=DC(i,0)+DC(i,k) CF(i,0)=CF(i,0)+ & & DC(i,k)*v(i,j,k,nnew) !^ tl_CF(i,0)=tl_CF(i,0)+ & !^ & tl_DC(i,k)*v(i,j,k,nnew)+ & !^ & DC(i,k)*tl_v(i,j,k,nnew) !^ ad_DC(i,k)=ad_DC(i,k)+ad_CF(i,0)*v(i,j,k,nnew) ad_v(i,j,k,nnew)=ad_v(i,j,k,nnew)+ad_CF(i,0)*DC(i,k) !^ tl_DC(i,0)=tl_DC(i,0)+tl_DC(i,k) !^ ad_DC(i,k)=ad_DC(i,k)+ad_DC(i,0) !^ tl_DC(i,k)=cff*(tl_Hz(i,j,k)+tl_Hz(i,j-1,k)) !^ adfac=cff*ad_DC(i,k) ad_Hz(i,j-1,k)=ad_Hz(i,j-1,k)+adfac ad_Hz(i,j ,k)=ad_Hz(i,j ,k)+adfac ad_DC(i,k)=0.0_r8 END DO END DO DO i=IstrT,IendT !^ tl_FC(i,0)=0.0_r8 !^ ad_FC(i,0)=0.0_r8 !^ tl_CF(i,0)=0.0_r8 !^ ad_CF(i,0)=0.0_r8 !^ tl_DC(i,0)=0.0_r8 !^ ad_DC(i,0)=0.0_r8 END DO END IF ! ! Couple adjoint velocity component in the XI-direction. First, ! compute BASIC STATE quantities. This section was computed three ! times in the original code. This can avoided by saving the ! intermediate values in scratch arrays DC1, CF1, and FC1 (HGA). ! DO i=IstrP,IendT DC(i,0)=0.0_r8 CF(i,0)=0.0_r8 FC(i,0)=0.0_r8 END DO DO k=1,N(ng) DO i=IstrP,IendT cff=0.5_r8*on_u(i,j) DC(i,k)=cff*(Hz(i,j,k)+Hz(i-1,j,k)) DC(i,0)=DC(i,0)+DC(i,k) CF(i,0)=CF(i,0)+ & & DC(i,k)*u(i,j,k,nnew) END DO END DO DO i=IstrP,IendT DC1(i,0)=DC(i,0) ! intermediate DC(i,0)=1.0_r8/DC(i,0) ! recursive CF1(i)=CF(i,0) ! intermediate CF(i,0)=DC(i,0)*(CF(i,0)-DU_avg1(i,j)) ! recursive END DO DO k=N(ng),1,-1 DO i=IstrP,IendT !^ Huon(i,j,k)=0.5_r8*(Huon(i,j,k)+u(i,j,k,nnew)*DC(i,k)) !^ FC(i,0)=FC(i,0)+Huon(i,j,k) !^ FC(i,0)=FC(i,0)+ & & 0.5_r8*(Huon(i,j,k)+u(i,j,k,nnew)*DC(i,k)) END DO END DO DO i=IstrP,IendT FC1(i)=FC(i,0) ! intermediate FC(i,0)=DC(i,0)*(FC(i,0)-DU_avg2(i,j)) ! recursive END DO ! ! Compute correct adjoint mass flux, Hz*u/n. ! DO k=1,N(ng) DO i=IstrP,IendT !^ tl_Huon(i,j,k)=tl_Huon(i,j,k)- & !^ & tl_DC(i,k)*FC(i,0)- & !^ & DC(i,k)*tl_FC(i,0) !^ ad_DC(i,k)=ad_DC(i,k)-ad_Huon(i,j,k)*FC(i,0) ad_FC(i,0)=ad_FC(i,0)-ad_Huon(i,j,k)*DC(i,k) END DO END DO DO i=IstrP,IendT !^ tl_FC(i,0)=tl_DC(i,0)*(FC1(i)-DU_avg2(i,j))+ & !^ & DC(i,0)*(tl_FC(i,0)-tl_DU_avg2(i,j)) !^ adfac=DC(i,0)*ad_FC(i,0) ad_DC(i,0)=ad_DC(i,0)+ad_FC(i,0)*(FC1(i)-DU_avg2(i,j)) ad_DU_avg2(i,j)=ad_DU_avg2(i,j)-adfac ad_FC(i,0)=adfac END DO !^ DO k=N(ng),1,-1 !^ DO k=1,N(ng) DO i=IstrP,IendT !^ tl_FC(i,0)=tl_FC(i,0)+tl_Huon(i,j,k) !^ ad_Huon(i,j,k)=ad_Huon(i,j,k)+ad_FC(i,0) !^ tl_Huon(i,j,k)=0.5_r8*(tl_Huon(i,j,k)+ & !^ & tl_u(i,j,k,nnew)*DC(i,k)+ & !^ & u(i,j,k,nnew)*tl_DC(i,k)) !^ adfac=0.5_r8*ad_Huon(i,j,k) ad_DC(i,k)=ad_DC(i,k)+u(i,j,k,nnew)*adfac ad_u(i,j,k,nnew)=ad_u(i,j,k,nnew)+DC(i,k)*adfac ad_Huon(i,j,k)=adfac END DO END DO ! ! Replace only BOUNDARY POINTS incorrect vertical mean with more ! accurate barotropic component, ubar=DU_avg1/(D*on_u). Recall that, ! D=CF(:,0). The replacement is avoided if the boundary edge is ! periodic. The J-loop is pipelined, so we need to do a special ! test for the southern and northern domain edges. ! IF (.not.(CompositeGrid(inorth,ng).or.NSperiodic(ng))) THEN IF (j.eq.Mm(ng)+1) THEN ! northern boundary DO k=1,N(ng) ! J-loop piplined DO i=IstrU,Iend !^ tl_u(i,j,k,nnew)=tl_u(i,j,k,nnew)* & !^ & umask(i,j) !^ ad_u(i,j,k,nnew)=ad_u(i,j,k,nnew)* & & umask(i,j) !^ tl_u(i,j,k,nnew)=tl_u(i,j,k,nnew)-tl_CF(i,0) !^ ad_CF(i,0)=ad_CF(i,0)-ad_u(i,j,k,nnew) END DO END DO END IF END IF ! IF (.not.(CompositeGrid(isouth,ng).or.NSperiodic(ng))) THEN IF (j.eq.0) THEN ! southern boundary DO k=1,N(ng) ! J-loop pipelined DO i=IstrU,Iend !^ tl_u(i,j,k,nnew)=tl_u(i,j,k,nnew)* & !^ & umask(i,j) !^ ad_u(i,j,k,nnew)=ad_u(i,j,k,nnew)* & & umask(i,j) !^ tl_u(i,j,k,nnew)=tl_u(i,j,k,nnew)-tl_CF(i,0) !^ ad_CF(i,0)=ad_CF(i,0)-ad_u(i,j,k,nnew) END DO END DO END IF END IF ! IF (.not.(CompositeGrid(ieast,ng).or.EWperiodic(ng))) THEN IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN DO k=1,N(ng) !^ tl_u(Iend+1,j,k,nnew)=tl_u(Iend+1,j,k,nnew)* & !^ & umask(Iend+1,j) !^ ad_u(Iend+1,j,k,nnew)=ad_u(Iend+1,j,k,nnew)* & & umask(Iend+1,j) !^ tl_u(Iend+1,j,k,nnew)=tl_u(Iend+1,j,k,nnew)- & !^ & tl_CF(Iend+1,0) !^ ad_CF(Iend+1,0)=ad_CF(Iend+1,0)- & & ad_u(Iend+1,j,k,nnew) END DO END IF END IF ! IF (.not.(CompositeGrid(iwest,ng).or.EWperiodic(ng))) THEN IF (DOMAIN(ng)%Western_Edge(tile)) THEN DO k=1,N(ng) !^ tl_u(Istr,j,k,nnew)=tl_u(Istr,j,k,nnew)* & !^ & umask(Istr,j) !^ ad_u(Istr,j,k,nnew)=ad_u(Istr,j,k,nnew)* & & umask(Istr,j) !^ tl_u(Istr,j,k,nnew)=tl_u(Istr,j,k,nnew)- & !^ & tl_CF(Istr,0) !^ ad_CF(Istr,0)=ad_CF(Istr,0)- & & ad_u(Istr,j,k,nnew) END DO END IF END IF ! ! Save adjoint solution for time-step iic(ng)-1 as the sum of time ! records 1 and 2. ! DO i=IstrP,IendT ad_ubar_sol(i,j)=ad_ubar(i,j,1)+ad_ubar(i,j,2) END DO ! ! Compute thicknesses of U-boxes DC(i,1:N), total depth of the water ! column DC(i,0), and incorrect vertical mean CF(i,0). Notice that ! barotropic component is replaced with its fast-time averaged ! values. ! DO i=IstrP,IendT !^ tl_ubar(i,j,2)=tl_ubar(i,j,1) !^ ad_ubar(i,j,1)=ad_ubar(i,j,1)+ad_ubar(i,j,2) ad_ubar(i,j,2)=0.0_r8 !^ tl_ubar(i,j,1)=tl_DC(i,0)*DU_avg1(i,j)+ & !^ & DC(i,0)*tl_DU_avg1(i,j) !^ ad_DC(i,0)=ad_DC(i,0)+ad_ubar(i,j,1)*DU_avg1(i,j) ad_DU_avg1(i,j)=ad_DU_avg1(i,j)+ad_ubar(i,j,1)*DC(i,0) ad_ubar(i,j,1)=0.0_r8 !^ tl_CF(i,0)=tl_DC(i,0)*(CF1(i)-DU_avg1(i,j))+ & !^ & DC(i,0)*(tl_CF(i,0)-tl_DU_avg1(i,j)) !^ adfac=DC(i,0)*ad_CF(i,0) ad_DC(i,0)=ad_DC(i,0)+ad_CF(i,0)*(CF1(i)-DU_avg1(i,j)) ad_DU_avg1(i,j)=ad_DU_avg1(i,j)-adfac ad_CF(i,0)=adfac !^ tl_DC(i,0)=-tl_DC(i,0)/(DC1(i,0)*DC1(i,0)) !^ ad_DC(i,0)=-ad_DC(i,0)/(DC1(i,0)*DC1(i,0)) END DO DO i=IstrP,IendT DC(i,0)=0.0_r8 CF(i,0)=0.0_r8 FC(i,0)=0.0_r8 END DO DO k=1,N(ng) DO i=IstrP,IendT cff=0.5_r8*on_u(i,j) DC(i,k)=0.5_r8*(Hz(i,j,k)+Hz(i-1,j,k))*on_u(i,j) DC(i,0)=DC(i,0)+DC(i,k) CF(i,0)=CF(i,0)+ & & DC(i,k)*u(i,j,k,nnew) !^ tl_CF(i,0)=tl_CF(i,0)+ & !^ & tl_DC(i,k)*u(i,j,k,nnew)+ & !^ & DC(i,k)*tl_u(i,j,k,nnew) !^ ad_DC(i,k)=ad_DC(i,k)+ad_CF(i,0)*u(i,j,k,nnew) ad_u(i,j,k,nnew)=ad_u(i,j,k,nnew)+ad_CF(i,0)*DC(i,k) !^ tl_DC(i,0)=tl_DC(i,0)+tl_DC(i,k) !^ ad_DC(i,k)=ad_DC(i,k)+ad_DC(i,0) !^ tl_DC(i,k)=cff*(tl_Hz(i,j,k)+tl_Hz(i-1,j,k)) !^ adfac=cff*ad_DC(i,k) ad_Hz(i-1,j,k)=ad_Hz(i-1,j,k)+adfac ad_Hz(i ,j,k)=ad_Hz(i ,j,k)+adfac ad_DC(i,k)=0.0_r8 END DO END DO DO i=IstrP,IendT !^ tl_FC(i,0)=0.0_r8 !^ ad_FC(i,0)=0.0_r8 !^ tl_CF(i,0)=0.0_r8 !^ ad_CF(i,0)=0.0_r8 !^ tl_DC(i,0)=0.0_r8 !^ ad_DC(i,0)=0.0_r8 END DO END DO J_LOOP1 ! !----------------------------------------------------------------------- ! Apply adjoint momentum transport point sources (like river runoff), ! if any. !----------------------------------------------------------------------- ! IF (LuvSrc(ng)) THEN DO is=1,Nsrc(ng) i=SOURCES(ng)%Isrc(is) j=SOURCES(ng)%Jsrc(is) IF (((IstrR.le.i).and.(i.le.IendR)).and. & & ((JstrR.le.j).and.(j.le.JendR))) THEN IF (INT(SOURCES(ng)%Dsrc(is)).eq.0) THEN DO k=1,N(ng) cff1=1.0_r8/(on_u(i,j)* & & 0.5_r8*(z_w(i-1,j,k)-z_w(i-1,j,k-1)+ & & z_w(i ,j,k)-z_w(i ,j,k-1))) !^ tl_u(i,j,k,nnew)=SOURCES(ng)%tl_Qsrc(is,k)*cff1+ & !^ & SOURCES(ng)%Qsrc(is,k)*tl_cff1 !^ SOURCES(ng)%ad_Qsrc(is,k)=SOURCES(ng)%ad_Qsrc(is,k)+ & & cff1*ad_u(i,j,k,nnew) ad_cff1=ad_cff1+ & & SOURCES(ng)%Qsrc(is,k)*ad_u(i,j,k,nnew) ad_u(i,j,k,nnew)=0.0_r8 !^ tl_cff1=-cff1*cff1*on_u(i,j)* & !^ & 0.5_r8*(tl_z_w(i-1,j,k)-tl_z_w(i-1,j,k-1)+ & !^ & tl_z_w(i ,j,k)-tl_z_w(i ,j,k-1)) !^ adfac=-0.5_r8*cff1*cff1*on_u(i,j)*ad_cff1 ad_z_w(i-1,j,k-1)=ad_z_w(i-1,j,k-1)-adfac ad_z_w(i ,j,k-1)=ad_z_w(i ,j,k-1)-adfac ad_z_w(i-1,j,k )=ad_z_w(i-1,j,k )+adfac ad_z_w(i ,j,k )=ad_z_w(i ,j,k )+adfac ad_cff1=0.0_r8 END DO ELSE IF (INT(SOURCES(ng)%Dsrc(is)).eq.1) THEN DO k=1,N(ng) cff1=1.0_r8/(om_v(i,j)* & & 0.5_r8*(z_w(i,j-1,k)-z_w(i,j-1,k-1)+ & & z_w(i,j ,k)-z_w(i,j ,k-1))) !^ tl_v(i,j,k,nnew)=SOURCES(ng)%tl_Qsrc(is,k)*cff1+ & !^ & SOURCES(ng)%Qsrc(is,k)*tl_cff1 !^ SOURCES(ng)%ad_Qsrc(is,k)=SOURCES(ng)%ad_Qsrc(is,k)+ & & cff1*ad_v(i,j,k,nnew) ad_cff1=ad_cff1+ & & SOURCES(ng)%Qsrc(is,k)*ad_v(i,j,k,nnew) ad_v(i,j,k,nnew)=0.0_r8 !^ tl_cff1=-cff1*cff1*om_v(i,j)* & !^ & 0.5_r8*(tl_z_w(i,j-1,k)-tl_z_w(i,j-1,k-1)+ & !^ & tl_z_w(i,j ,k)-tl_z_w(i,j ,k-1)) !^ adfac=-0.5_r8*cff1*cff1*om_v(i,j)*ad_cff1 ad_z_w(i,j-1,k-1)=ad_z_w(i,j-1,k-1)-adfac ad_z_w(i,j ,k-1)=ad_z_w(i,j ,k-1)-adfac ad_z_w(i,j-1,k )=ad_z_w(i,j-1,k )+adfac ad_z_w(i,j ,k )=ad_z_w(i,j ,k )+adfac ad_cff1=0.0_r8 END DO END IF END IF END DO END IF ! !----------------------------------------------------------------------- ! Set lateral adjoint boundary conditions. !----------------------------------------------------------------------- ! !^ CALL tl_v3dbc_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, N(ng), & !^ & IminS, ImaxS, JminS, JmaxS, & !^ & nstp, nnew, & !^ & v) !^ CALL ad_v3dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, N(ng), & & IminS, ImaxS, JminS, JmaxS, & & nstp, nnew, & & ad_v) !^ CALL tl_u3dbc_tile (ng, tile, & !^ & LBi, UBi, LBj, UBj, N(ng), & !^ & IminS, ImaxS, JminS, JmaxS, & !^ & nstp, nnew, & !^ & tl_u) !^ CALL ad_u3dbc_tile (ng, tile, & & LBi, UBi, LBj, UBj, N(ng), & & IminS, ImaxS, JminS, JmaxS, & & nstp, nnew, & & ad_u) ! !----------------------------------------------------------------------- ! Time step momentum equation in the ETA-direction. !----------------------------------------------------------------------- ! J_LOOP2 : DO j=Jstr,Jend IF (j.ge.JstrV) THEN ! ! Compute BASIC STATE quantities. ! DO i=Istr,Iend AK(i,0)=0.5_r8*(Akv(i,j-1,0)+ & & Akv(i,j ,0)) DO k=1,N(ng) AK(i,k)=0.5_r8*(Akv(i,j-1,k)+ & & Akv(i,j ,k)) Hzk(i,k)=0.5_r8*(Hz(i,j-1,k)+ & & Hz(i,j ,k)) END DO END DO ! ! Couple and update new adjoint solution. ! DO k=1,N(ng) DO i=Istr,Iend !^ tl_v(i,j,k,nnew)=tl_v(i,j,k,nnew)*vmask(i,j) !^ ad_v(i,j,k,nnew)=ad_v(i,j,k,nnew)*vmask(i,j) !^ tl_v(i,j,k,nnew)=tl_v(i,j,k,nnew)-tl_DC(i,0) !^ ad_DC(i,0)=ad_DC(i,0)-ad_v(i,j,k,nnew) END DO END DO ! ! Replace INTERIOR POINTS incorrect vertical mean with more accurate ! barotropic component, vbar=DV_avg1/(D*om_v). Recall that, D=CF(:,0). ! DO i=Istr,Iend CF(i,0)=Hzk(i,1) DC(i,0)=v(i,j,1,nnew)*Hzk(i,1) END DO DO k=2,N(ng) DO i=Istr,Iend CF(i,0)=CF(i,0)+Hzk(i,k) DC(i,0)=DC(i,0)+v(i,j,k,nnew)*Hzk(i,k) END DO END DO DO i=Istr,Iend DC1(i,0)=DC(i,0) ! intermediate cff1=1.0_r8/(CF(i,0)*om_v(i,j)) DC(i,0)=(DC(i,0)*om_v(i,j)-DV_avg1(i,j))*cff1 ! recursive !^ tl_DC(i,0)=(tl_DC(i,0)*om_v(i,j)-tl_DV_avg1(i,j))*cff1+ & !^ & (DC1(i,0)*om_v(i,j)-DV_avg1(i,j))*tl_cff1 !^ adfac=cff1*ad_DC(i,0) ad_DV_avg1(i,j)=ad_DV_avg1(i,j)-adfac ad_cff1=ad_cff1+ & & (DC1(i,0)*om_v(i,j)-DV_avg1(i,j))*ad_DC(i,0) ad_DC(i,0)=om_v(i,j)*adfac !^ tl_cff1=-cff1*cff1*tl_CF(i,0)*om_v(i,j) !^ ad_CF(i,0)=ad_CF(i,0)-cff1*cff1*om_v(i,j)*ad_cff1 ad_cff1=0.0_r8 END DO DO i=Istr,Iend CF(i,0)=Hzk(i,1) DC(i,0)=v(i,j,1,nnew)*Hzk(i,1) END DO DO k=2,N(ng) DO i=Istr,Iend CF(i,0)=CF(i,0)+Hzk(i,k) DC(i,0)=DC(i,0)+v(i,j,k,nnew)*Hzk(i,k) !^ tl_DC(i,0)=tl_DC(i,0)+ & !^ & tl_v(i,j,k,nnew)*Hzk(i,k)+ & !^ & v(i,j,k,nnew)*tl_Hzk(i,k) !^ ad_v(i,j,k,nnew)=ad_v(i,j,k,nnew)+ad_DC(i,0)*Hzk(i,k) ad_Hzk(i,k)=ad_Hzk(i,k)+v(i,j,k,nnew)*ad_DC(i,0) !^ tl_CF(i,0)=tl_CF(i,0)+tl_Hzk(i,k) !^ ad_Hzk(i,k)=ad_Hzk(i,k)+ad_CF(i,0) END DO END DO DO i=Istr,Iend CF(i,0)=Hzk(i,1) DC(i,0)=v(i,j,1,nnew)*Hzk(i,1) !^ tl_DC(i,0)=tl_v(i,j,1,nnew)*Hzk(i,1)+ & !^ & v(i,j,1,nnew)*tl_Hzk(i,1) !^ ad_v(i,j,1,nnew)=ad_v(i,j,1,nnew)+ad_DC(i,0)*Hzk(i,1) ad_Hzk(i,1)=ad_Hzk(i,1)+v(i,j,1,nnew)*ad_DC(i,0) ad_DC(i,0)=0.0_r8 !^ tl_CF(i,0)=tl_Hzk(i,1) !^ ad_Hzk(i,1)=ad_Hzk(i,1)+ad_CF(i,0) ad_CF(i,0)=0.0_r8 END DO ! ! Compute BASIC STATE off-diagonal coefficients [lambda*dt*Akv/Hz] ! for the implicit vertical viscosity term at horizontal V-points ! and vertical W-points. ! ! NOTE: The original code solves the tridiagonal system A*v=r where ! A is a matrix and v and r are vectors. We need to solve the ! tangent linear of this system which is A*tl_v+tl_A*v=tl_r. ! Here, tl_A*v and tl_r are known, so we must solve for tl_v ! by inverting A*tl_v=tl_r-tl_A*v. ! cff=-lambda*dt(ng)/0.5_r8 DO k=1,N(ng)-1 DO i=Istr,Iend cff1=1.0_r8/(z_r(i,j,k+1)+z_r(i,j-1,k+1)- & & z_r(i,j,k )-z_r(i,j-1,k )) FC(i,k)=cff*cff1*AK(i,k) END DO END DO DO i=Istr,Iend FC(i,0)=0.0_r8 FC(i,N(ng))=0.0_r8 END DO DO k=1,N(ng) DO i=Istr,Iend BC(i,k)=Hzk(i,k)-FC(i,k)-FC(i,k-1) END DO END DO DO i=Istr,Iend cff=1.0_r8/BC(i,1) CF(i,1)=cff*FC(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) END DO END DO ! ! Compute new adjoint solution by back substitution. ! (DC is a tangent linear variable here). ! !^ DO k=N(ng)-1,1,-1 !^ DO k=1,N(ng)-1 DO i=Istr,Iend !^ tl_v(i,j,k,nnew)=DC(i,k) !^ ad_DC(i,k)=ad_DC(i,k)+ad_v(i,j,k,nnew) ad_v(i,j,k,nnew)=0.0_r8 !^ DC(i,k)=DC(i,k)-CF(i,k)*DC(i,k+1) !^ ad_DC(i,k+1)=-CF(i,k)*ad_DC(i,k) END DO END DO DO i=Istr,Iend !^ tl_v(i,j,N(ng),nnew)=DC(i,N(ng)) !^ ad_DC(i,N(ng))=ad_DC(i,N(ng))+ad_v(i,j,N(ng),nnew) ad_v(i,j,N(ng),nnew)=0.0_r8 !^ 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)) !^ adfac=ad_DC(i,N(ng))/ & & (BC(i,N(ng))-FC(i,N(ng)-1)*CF(i,N(ng)-1)) ad_DC(i,N(ng)-1)=ad_DC(i,N(ng)-1)-FC(i,N(ng)-1)*adfac ad_DC(i,N(ng) )=adfac END DO ! ! Solve the adjoint tridiagonal system. ! (DC is a tangent linear variable here). ! !^ DO k=2,N(ng)-1 !^ DO k=N(ng)-1,2,-1 DO i=Istr,Iend cff=1.0_r8/(BC(i,k)-FC(i,k-1)*CF(i,k-1)) !^ DC(i,k)=cff*(DC(i,k)-FC(i,k-1)*DC(i,k-1)) !^ adfac=cff*ad_DC(i,k) ad_DC(i,k-1)=ad_DC(i,k-1)-FC(i,k-1)*adfac ad_DC(i,k )=adfac END DO END DO DO i=Istr,Iend cff=1.0_r8/BC(i,1) !^ DC(i,1)=cff*DC(i,1) !^ ad_DC(i,1)=cff*ad_DC(i,1) END DO DO i=Istr,Iend !^ DC(i,N(ng))=tl_v(i,j,N(ng),nnew)- & !^ & (tl_FC(i,N(ng)-1)*v(i,j,N(ng)-1,nnew)+ & !^ & tl_BC(i,N(ng) )*v(i,j,N(ng) ,nnew)) !^ ad_BC(i,N(ng) )=-v(i,j,N(ng) ,nnew)*ad_DC(i,N(ng)) ad_FC(i,N(ng)-1)=-v(i,j,N(ng)-1,nnew)*ad_DC(i,N(ng)) ad_v(i,j,N(ng),nnew)=ad_DC(i,N(ng)) ad_DC(i,N(ng))=0.0_r8 !^ DC(i,1)=tl_v(i,j,1,nnew)- & !^ & (tl_BC(i,1)*v(i,j,1,nnew)+ & !^ & tl_FC(i,1)*v(i,j,2,nnew)) !^ ad_FC(i,1)=-v(i,j,2,nnew)*ad_DC(i,1) ad_BC(i,1)=-v(i,j,1,nnew)*ad_DC(i,1) ad_v(i,j,1,nnew)=ad_DC(i,1) ad_DC(i,1)=0.0_r8 END DO DO k=2,N(ng)-1 DO i=Istr,Iend !^ DC(i,k)=tl_v(i,j,k,nnew)- & !^ & (tl_FC(i,k-1)*v(i,j,k-1,nnew)+ & !^ & tl_BC(i,k )*v(i,j,k ,nnew)+ & !^ & tl_FC(i,k )*v(i,j,k+1,nnew)) !^ ad_FC(i,k-1)=ad_FC(i,k-1)-v(i,j,k-1,nnew)*ad_DC(i,k) ad_FC(i,k )=ad_FC(i,k )-v(i,j,k+1,nnew)*ad_DC(i,k) ad_BC(i,k )=ad_BC(i,k )-v(i,j,k ,nnew)*ad_DC(i,k) ad_v(i,j,k,nnew)=ad_v(i,j,k,nnew)+ad_DC(i,k) ad_DC(i,k)=0.0_r8 END DO END DO DO k=1,N(ng) DO i=Istr,Iend !^ tl_BC(i,k)=tl_Hzk(i,k)-tl_FC(i,k)-tl_FC(i,k-1) !^ ad_Hzk(i,k)=ad_Hzk(i,k)+ad_BC(i,k) ad_FC(i,k-1)=ad_FC(i,k-1)-ad_BC(i,k) ad_FC(i,k )=ad_FC(i,k )-ad_BC(i,k) ad_BC(i,k)=0.0_r8 END DO END DO ! ! Compute adjoint off-diagonal coefficients [lambda*dt*Akv/Hz] for ! the implicit vertical viscosity term at horizontal V-points and ! vertical W-points. ! DO i=Istr,Iend !^ tl_FC(i,N(ng))=0.0_r8 !^ ad_FC(i,N(ng))=0.0_r8 !^ tl_FC(i,0)=0.0_r8 !^ ad_FC(i,0)=0.0_r8 END DO cff=-lambda*dt(ng)/0.5_r8 DO k=1,N(ng)-1 DO i=Istr,Iend cff1=1.0_r8/(z_r(i,j,k+1)+z_r(i,j-1,k+1)- & & z_r(i,j,k )-z_r(i,j-1,k )) !^ tl_FC(i,k)=cff*(tl_cff1*AK(i,k)+cff1*tl_AK(i,k)) !^ adfac=cff*ad_FC(i,k) ad_AK(i,k)=ad_AK(i,k)+cff1*adfac ad_cff1=ad_cff1+AK(i,k)*adfac ad_FC(i,k)=0.0_r8 !^ tl_cff1=-cff1*cff1*(tl_z_r(i,j,k+1)+tl_z_r(i,j-1,k+1)- & !^ & tl_z_r(i,j,k )-tl_z_r(i,j-1,k )) !^ adfac=-cff1*cff1*ad_cff1 ad_z_r(i,j-1,k )=ad_z_r(i,j-1,k )-adfac ad_z_r(i,j ,k )=ad_z_r(i,j ,k )-adfac ad_z_r(i,j-1,k+1)=ad_z_r(i,j-1,k+1)+adfac ad_z_r(i,j ,k+1)=ad_z_r(i,j ,k+1)+adfac ad_cff1=0.0_r8 END DO END DO ! ! Time step adjoint right-hand-side terms. ! IF (iic(ng).eq.ntfirst(ng)) THEN cff=0.25_r8*dt(ng) ELSE IF (iic(ng).eq.(ntfirst(ng)+1)) THEN cff=0.25_r8*dt(ng)*3.0_r8/2.0_r8 ELSE cff=0.25_r8*dt(ng)*23.0_r8/12.0_r8 END IF DO i=Istr,Iend DC(i,0)=cff*(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1)) END DO ! ! The BASIC STATE "v" used below must be in transport units, but "v" ! retrived is in m/s so we multiply by "Hzk". ! DO k=1,N(ng) DO i=Istr,Iend !^ tl_v(i,j,k,nnew)=tl_v(i,j,k,nnew)+ & !^ & DC(i,0)*tl_rv(i,j,k,nrhs) !^ ad_rv(i,j,k,nrhs)=ad_rv(i,j,k,nrhs)+ & & DC(i,0)*ad_v(i,j,k,nnew) END DO END DO DO i=Istr,Iend DO k=1,N(ng) !^ tl_Hzk(i,k)=0.5_r8*(tl_Hz(i,j-1,k)+ & !^ & tl_Hz(i,j ,k)) !^ adfac=0.5_r8*ad_Hzk(i,k) ad_Hz(i,j-1,k)=ad_Hz(i,j-1,k)+adfac ad_Hz(i,j ,k)=ad_Hz(i,j ,k)+adfac ad_Hzk(i,k)=0.0_r8 !^ tl_AK(i,k)=0.5_r8*(tl_Akv(i,j-1,k)+ & !^ & tl_Akv(i,j ,k)) !^ adfac=0.5_r8*ad_AK(i,k) ad_Akv(i,j-1,k)=ad_Akv(i,j-1,k)+adfac ad_Akv(i,j ,k)=ad_Akv(i,j ,k)+adfac ad_AK(i,k)=0.0_r8 END DO !^ tl_AK(i,0)=0.5_r8*(tl_Akv(i,j-1,0)+ & !^ & tl_Akv(i,j ,0)) !^ adfac=0.5_r8*ad_AK(i,0) ad_Akv(i,j-1,0)=ad_Akv(i,j-1,0)+adfac ad_Akv(i,j ,0)=ad_Akv(i,j ,0)+adfac ad_AK(i,0)=0.0_r8 END DO END IF ! !----------------------------------------------------------------------- ! Time step momentum equation in the XI-direction. !----------------------------------------------------------------------- ! ! Compute BASIC STATE quatities. ! DO i=IstrU,Iend AK(i,0)=0.5_r8*(Akv(i-1,j,0)+ & & Akv(i ,j,0)) DO k=1,N(ng) AK(i,k)=0.5_r8*(Akv(i-1,j,k)+ & & Akv(i ,j,k)) Hzk(i,k)=0.5_r8*(Hz(i-1,j,k)+ & & Hz(i ,j,k)) END DO END DO ! ! Couple and update new adjoint solution. ! DO k=1,N(ng) DO i=IstrU,Iend !^ tl_u(i,j,k,nnew)=tl_u(i,j,k,nnew)*umask(i,j) !^ ad_u(i,j,k,nnew)=ad_u(i,j,k,nnew)*umask(i,j) !^ tl_u(i,j,k,nnew)=tl_u(i,j,k,nnew)-tl_DC(i,0) !^ ad_DC(i,0)=ad_DC(i,0)-ad_u(i,j,k,nnew) END DO END DO ! ! Replace INTERIOR POINTS incorrect vertical mean with more accurate ! barotropic component, ubar=DU_avg1/(D*on_u). Recall that, D=CF(:,0). ! DO i=IstrU,Iend CF(i,0)=Hzk(i,1) DC(i,0)=u(i,j,1,nnew)*Hzk(i,1) END DO DO k=2,N(ng) DO i=IstrU,Iend CF(i,0)=CF(i,0)+Hzk(i,k) DC(i,0)=DC(i,0)+u(i,j,k,nnew)*Hzk(i,k) END DO END DO DO i=IstrU,Iend DC1(i,0)=DC(i,0) ! intermediate cff1=1.0/(CF(i,0)*on_u(i,j)) DC(i,0)=(DC(i,0)*on_u(i,j)-DU_avg1(i,j))*cff1 ! recursive !^ tl_DC(i,0)=(tl_DC(i,0)*on_u(i,j)-tl_DU_avg1(i,j))*cff1+ & !^ & (DC1(i,0)*on_u(i,j)-DU_avg1(i,j))*tl_cff1 !^ adfac=cff1*ad_DC(i,0) ad_DU_avg1(i,j)=ad_DU_avg1(i,j)-adfac ad_cff1=ad_cff1+ & & (DC1(i,0)*on_u(i,j)-DU_avg1(i,j))*ad_DC(i,0) ad_DC(i,0)=on_u(i,j)*adfac !^ tl_cff1=-cff1*cff1*tl_CF(i,0)*on_u(i,j) !^ ad_CF(i,0)=ad_CF(i,0)-cff1*cff1*on_u(i,j)*ad_cff1 ad_cff1=0.0_r8 END DO DO i=IstrU,Iend CF(i,0)=Hzk(i,1) DC(i,0)=u(i,j,1,nnew)*Hzk(i,1) END DO DO k=2,N(ng) DO i=IstrU,Iend CF(i,0)=CF(i,0)+Hzk(i,k) DC(i,0)=DC(i,0)+u(i,j,k,nnew)*Hzk(i,k) !^ tl_DC(i,0)=tl_DC(i,0)+ & !^ & tl_u(i,j,k,nnew)*Hzk(i,k)+ & !^ & u(i,j,k,nnew)*tl_Hzk(i,k) !^ ad_u(i,j,k,nnew)=ad_u(i,j,k,nnew)+ad_DC(i,0)*Hzk(i,k) ad_Hzk(i,k)=ad_Hzk(i,k)+u(i,j,k,nnew)*ad_DC(i,0) !^ tl_CF(i,0)=tl_CF(i,0)+tl_Hzk(i,k) !^ ad_Hzk(i,k)=ad_Hzk(i,k)+ad_CF(i,0) END DO END DO DO i=IstrU,Iend CF(i,0)=Hzk(i,1) DC(i,0)=u(i,j,1,nnew)*Hzk(i,1) !^ tl_DC(i,0)=tl_u(i,j,1,nnew)*Hzk(i,1)+ & !^ & u(i,j,1,nnew)*tl_Hzk(i,1) !^ ad_u(i,j,1,nnew)=ad_u(i,j,1,nnew)+ad_DC(i,0)*Hzk(i,1) ad_Hzk(i,1)=ad_Hzk(i,1)+u(i,j,1,nnew)*ad_DC(i,0) ad_DC(i,0)=0.0_r8 !^ tl_CF(i,0)=tl_Hzk(i,1) !^ ad_Hzk(i,1)=ad_Hzk(i,1)+ad_CF(i,0) ad_CF(i,0)=0.0_r8 END DO ! ! Compute BASIC STATE off-diagonal coefficients [lambda*dt*Akv/Hz] ! for the implicit vertical viscosity term at horizontal U-points ! and vertical W-points. ! ! NOTE: The original code solves the tridiagonal system A*u=r where ! A is a matrix and u and r are vectors. We need to solve the ! tangent linear of this system which is A*tl_u+tl_A*u=tl_r. ! Here, tl_A*u and tl_r are known, so we must solve for tl_u ! by inverting A*tl_u=tl_r-tl_A*u. ! cff=-lambda*dt(ng)/0.5_r8 DO k=1,N(ng)-1 DO i=IstrU,Iend cff1=1.0_r8/(z_r(i,j,k+1)+z_r(i-1,j,k+1)- & & z_r(i,j,k )-z_r(i-1,j,k )) FC(i,k)=cff*cff1*AK(i,k) END DO END DO DO i=IstrU,Iend FC(i,0)=0.0_r8 FC(i,N(ng))=0.0_r8 END DO DO k=1,N(ng) DO i=IstrU,Iend BC(i,k)=Hzk(i,k)-FC(i,k)-FC(i,k-1) END DO END DO DO i=IstrU,Iend cff=1.0_r8/BC(i,1) CF(i,1)=cff*FC(i,1) END DO DO k=2,N(ng)-1 DO i=IstrU,Iend cff=1.0_r8/(BC(i,k)-FC(i,k-1)*CF(i,k-1)) CF(i,k)=cff*FC(i,k) END DO END DO ! ! Compute new adjoint solution by back substitution. ! (DC is a tangent linear variable here). ! !^ DO k=N(ng)-1,1,-1 !^ DO k=1,N(ng)-1 DO i=IstrU,Iend !^ tl_u(i,j,k,nnew)=DC(i,k) !^ ad_DC(i,k)=ad_DC(i,k)+ad_u(i,j,k,nnew) ad_u(i,j,k,nnew)=0.0_r8 !^ DC(i,k)=DC(i,k)-CF(i,k)*DC(i,k+1) !^ ad_DC(i,k+1)=-CF(i,k)*ad_DC(i,k) END DO END DO DO i=IstrU,Iend !^ tl_u(i,j,N(ng),nnew)=DC(i,N(ng)) !^ ad_DC(i,N(ng))=ad_DC(i,N(ng))+ad_u(i,j,N(ng),nnew) ad_u(i,j,N(ng),nnew)=0.0_r8 !^ 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)) !^ adfac=ad_DC(i,N(ng))/ & & (BC(i,N(ng))-FC(i,N(ng)-1)*CF(i,N(ng)-1)) ad_DC(i,N(ng)-1)=ad_DC(i,N(ng)-1)-FC(i,N(ng)-1)*adfac ad_DC(i,N(ng) )=adfac END DO ! ! Solve the adjoint tridiagonal system. ! (DC is a tangent linear variable here). ! !^ DO k=2,N(ng)-1 !^ DO k=N(ng)-1,2,-1 DO i=IstrU,Iend cff=1.0_r8/(BC(i,k)-FC(i,k-1)*CF(i,k-1)) !^ DC(i,k)=cff*(DC(i,k)-FC(i,k-1)*DC(i,k-1)) !^ adfac=cff*ad_DC(i,k) ad_DC(i,k-1)=ad_DC(i,k-1)-FC(i,k-1)*adfac ad_DC(i,k )=adfac END DO END DO DO i=IstrU,Iend cff=1.0_r8/BC(i,1) !^ DC(i,1)=cff*DC(i,1) !^ ad_DC(i,1)=cff*ad_DC(i,1) END DO DO i=IstrU,Iend !^ DC(i,N(ng))=tl_u(i,j,N(ng),nnew)- & !^ & (tl_FC(i,N(ng)-1)*u(i,j,N(ng)-1,nnew)+ & !^ & tl_BC(i,N(ng) )*u(i,j,N(ng) ,nnew)) !^ ad_BC(i,N(ng) )=-u(i,j,N(ng) ,nnew)*ad_DC(i,N(ng)) ad_FC(i,N(ng)-1)=-u(i,j,N(ng)-1,nnew)*ad_DC(i,N(ng)) ad_u(i,j,N(ng),nnew)=ad_DC(i,N(ng)) ad_DC(i,N(ng))=0.0_r8 !^ DC(i,1)=tl_u(i,j,1,nnew)- & !^ & (tl_BC(i,1)*u(i,j,1,nnew)+ & !^ & tl_FC(i,1)*u(i,j,2,nnew)) !^ ad_FC(i,1)=-u(i,j,2,nnew)*ad_DC(i,1) ad_BC(i,1)=-u(i,j,1,nnew)*ad_DC(i,1) ad_u(i,j,1,nnew)=ad_DC(i,1) ad_DC(i,1)=0.0_r8 END DO DO k=2,N(ng)-1 DO i=IstrU,Iend !^ DC(i,k)=tl_u(i,j,k,nnew)- & !^ & (tl_FC(i,k-1)*u(i,j,k-1,nnew)+ & !^ & tl_BC(i,k )*u(i,j,k ,nnew)+ & !^ & tl_FC(i,k )*u(i,j,k+1,nnew)) !^ ad_FC(i,k-1)=ad_FC(i,k-1)-u(i,j,k-1,nnew)*ad_DC(i,k) ad_FC(i,k )=ad_FC(i,k )-u(i,j,k+1,nnew)*ad_DC(i,k) ad_BC(i,k )=ad_BC(i,k )-u(i,j,k ,nnew)*ad_DC(i,k) ad_u(i,j,k,nnew)=ad_u(i,j,k,nnew)+ad_DC(i,k) ad_DC(i,k)=0.0_r8 END DO END DO DO k=1,N(ng) DO i=IstrU,Iend !^ tl_BC(i,k)=tl_Hzk(i,k)-tl_FC(i,k)-tl_FC(i,k-1) !^ ad_Hzk(i,k)=ad_Hzk(i,k)+ad_BC(i,k) ad_FC(i,k-1)=ad_FC(i,k-1)-ad_BC(i,k) ad_FC(i,k )=ad_FC(i,k )-ad_BC(i,k) ad_BC(i,k)=0.0_r8 END DO END DO ! ! Compute adjoint off-diagonal coefficients [lambda*dt*Akv/Hz] for ! the implicit vertical viscosity term at horizontal U-points and ! vertical W-points. ! DO i=IstrU,Iend !^ tl_FC(i,N(ng))=0.0_r8 !^ ad_FC(i,N)=0.0_r8 !^ tl_FC(i,0)=0.0_r8 !^ ad_FC(i,0)=0.0_r8 END DO cff=-lambda*dt(ng)/0.5_r8 DO k=1,N(ng)-1 DO i=IstrU,Iend cff1=1.0_r8/(z_r(i,j,k+1)+z_r(i-1,j,k+1)- & & z_r(i,j,k )-z_r(i-1,j,k )) !^ tl_FC(i,k)=cff*(tl_cff1*AK(i,k)+cff1*tl_AK(i,k)) !^ adfac=cff*ad_FC(i,k) ad_AK(i,k)=ad_AK(i,k)+cff1*adfac ad_cff1=ad_cff1+AK(i,k)*adfac ad_FC(i,k)=0.0_r8 !^ tl_cff1=-cff1*cff1*(tl_z_r(i,j,k+1)+tl_z_r(i-1,j,k+1)- & !^ & tl_z_r(i,j,k )-tl_z_r(i-1,j,k )) !^ adfac=-cff1*cff1*ad_cff1 ad_z_r(i-1,j,k )=ad_z_r(i-1,j,k )-adfac ad_z_r(i ,j,k )=ad_z_r(i ,j,k )-adfac ad_z_r(i-1,j,k+1)=ad_z_r(i-1,j,k+1)+adfac ad_z_r(i ,j,k+1)=ad_z_r(i ,j,k+1)+adfac ad_cff1=0.0_r8 END DO END DO ! ! Time step adjoint right-hand-side terms. ! IF (iic(ng).eq.ntfirst(ng)) THEN cff=0.25_r8*dt(ng) ELSE IF (iic(ng).eq.(ntfirst(ng)+1)) THEN cff=0.25_r8*dt(ng)*3.0_r8/2.0_r8 ELSE cff=0.25_r8*dt(ng)*23.0_r8/12.0_r8 END IF DO i=IstrU,Iend DC(i,0)=cff*(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j)) END DO ! ! The BASIC STATE "u" used below must be in transport units, but "u" ! retrived is in m/s so we multiply by "Hzk". ! DO k=1,N(ng) DO i=IstrU,Iend !^ tl_u(i,j,k,nnew)=tl_u(i,j,k,nnew)+ & !^ & DC(i,0)*tl_ru(i,j,k,nrhs) !^ ad_ru(i,j,k,nrhs)=ad_ru(i,j,k,nrhs)+ & & DC(i,0)*ad_u(i,j,k,nnew) END DO END DO DO i=IstrU,Iend DO k=1,N(ng) !^ tl_Hzk(i,k)=0.5_r8*(tl_Hz(i-1,j,k)+ & !^ & tl_Hz(i ,j,k)) !^ adfac=0.5_r8*ad_Hzk(i,k) ad_Hz(i-1,j,k)=ad_Hz(i-1,j,k)+adfac ad_Hz(i ,j,k)=ad_Hz(i ,j,k)+adfac ad_Hzk(i,k)=0.0_r8 !^ tl_AK(i,k)=0.5_r8*(tl_Akv(i-1,j,k)+ & !^ & tl_Akv(i ,j,k)) !^ adfac=0.5_r8*ad_AK(i,k) ad_Akv(i-1,j,k)=ad_Akv(i-1,j,k)+adfac ad_Akv(i ,j,k)=ad_Akv(i ,j,k)+adfac ad_AK(i,k)=0.0_r8 END DO !^ tl_AK(i,0)=0.5_r8*(tl_Akv(i-1,j,0)+ & !^ & tl_Akv(i ,j,0)) !^ adfac=0.5_r8*ad_AK(i,0) ad_Akv(i-1,j,0)=ad_Akv(i-1,j,0)+adfac ad_Akv(i ,j,0)=ad_Akv(i ,j,0)+adfac ad_AK(i,0)=0.0_r8 END DO END DO J_LOOP2 ! RETURN END SUBROUTINE ad_step3d_uv_tile END MODULE ad_step3d_uv_mod