!
!git $Id$
!svn $Id: ad_step2d.F 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 performs a fast (predictor or corrector) time-step  !
!  for the free-surface and 2D momentum adjoint equations.             !
!  It also calculates the time filtering variables over all fast-time  !
!  steps to damp high frequency signals in 3D applications.            !
!                                                                      !
!=======================================================================
!
      MODULE ad_step2d_mod
!
!git $Id$
!svn $Id: ad_step2d_LF_AM3.h 1188 2023-08-03 19:26:47Z arango $
!=======================================================================
!                                                                      !
!  Adjoint shallow-water primitive equations predictor (Leap-frog)     !
!  and corrector (Adams-Moulton) time-stepping engine.                 !
!                                                                      !
!=======================================================================
!
      implicit none
!
      PRIVATE
      PUBLIC  :: ad_step2d
!
      CONTAINS
!
!***********************************************************************
      SUBROUTINE ad_step2d (ng, tile)
!***********************************************************************
!
      USE mod_param
      USE mod_ncparam
      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_step2d_LF_AM3.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, 9, 56, MyFile)
      CALL ad_step2d_tile (ng, tile,                                    &
     &                     LBi, UBi, LBj, UBj, N(ng),                   &
     &                     IminS, ImaxS, JminS, JmaxS,                  &
     &                     krhs(ng), kstp(ng), knew(ng),                &
     &                     nstp(ng), nnew(ng),                          &
     &                     GRID(ng) % pmask,       GRID(ng) % rmask,    &
     &                     GRID(ng) % umask,       GRID(ng) % vmask,    &
     &                     GRID(ng) % fomn,                             &
     &                     GRID(ng) % h,           GRID(ng) % ad_h,     &
     &                     GRID(ng) % om_u,        GRID(ng) % om_v,     &
     &                     GRID(ng) % on_u,        GRID(ng) % on_v,     &
     &                     GRID(ng) % omn,                              &
     &                     GRID(ng) % pm,          GRID(ng) % pn,       &
     &                     GRID(ng) % dndx,        GRID(ng) % dmde,     &
     &                     GRID(ng) % pmon_r,      GRID(ng) % pnom_r,   &
     &                     GRID(ng) % pmon_p,      GRID(ng) % pnom_p,   &
     &                     GRID(ng) % om_r,        GRID(ng) % on_r,     &
     &                     GRID(ng) % om_p,        GRID(ng) % on_p,     &
     &                     MIXING(ng) % visc2_p,   MIXING(ng) % visc2_r,&
     &                     COUPLING(ng) % ad_DU_avg1,                   &
     &                     COUPLING(ng) % ad_DU_avg2,                   &
     &                     COUPLING(ng) % ad_DV_avg1,                   &
     &                     COUPLING(ng) % ad_DV_avg2,                   &
     &                     COUPLING(ng) % Zt_avg1,                      &
     &                     COUPLING(ng) % ad_Zt_avg1,                   &
     &                     COUPLING(ng) % ad_rufrc,                     &
     &                     COUPLING(ng) % ad_rvfrc,                     &
     &                     OCEAN(ng) % ad_ru,                           &
     &                     OCEAN(ng) % ad_rv,                           &
     &                     OCEAN(ng) % rubar,      OCEAN(ng) % ad_rubar,&
     &                     OCEAN(ng) % rvbar,      OCEAN(ng) % ad_rvbar,&
     &                     OCEAN(ng) % rzeta,      OCEAN(ng) % ad_rzeta,&
     &                     OCEAN(ng) % ubar,       OCEAN(ng) % ad_ubar, &
     &                     OCEAN(ng) % vbar,       OCEAN(ng) % ad_vbar, &
     &                     OCEAN(ng) % zeta,       OCEAN(ng) % ad_zeta)
      CALL wclock_off (ng, iADM, 9, 162, MyFile)
!
      RETURN
      END SUBROUTINE ad_step2d
!
!***********************************************************************
      SUBROUTINE ad_step2d_tile (ng, tile,                              &
     &                           LBi, UBi, LBj, UBj, UBk,               &
     &                           IminS, ImaxS, JminS, JmaxS,            &
     &                           krhs, kstp, knew,                      &
     &                           nstp, nnew,                            &
     &                           pmask, rmask, umask, vmask,            &
     &                           fomn,                                  &
     &                           h, ad_h,                               &
     &                           om_u, om_v, on_u, on_v, omn, pm, pn,   &
     &                           dndx, dmde,                            &
     &                           pmon_r, pnom_r, pmon_p, pnom_p,        &
     &                           om_r, on_r, om_p, on_p,                &
     &                           visc2_p, visc2_r,                      &
     &                           ad_DU_avg1, ad_DU_avg2,                &
     &                           ad_DV_avg1, ad_DV_avg2,                &
     &                           Zt_avg1, ad_Zt_avg1,                   &
     &                           ad_rufrc, ad_rvfrc,                    &
     &                           ad_ru, ad_rv,                          &
     &                           rubar, ad_rubar,                       &
     &                           rvbar, ad_rvbar,                       &
     &                           rzeta, ad_rzeta,                       &
     &                           ubar, ad_ubar,                         &
     &                           vbar, ad_vbar,                         &
     &                           zeta, ad_zeta)
!***********************************************************************
!
      USE mod_param
      USE mod_clima
      USE mod_ncparam
      USE mod_scalars
      USE mod_sources
!
      USE ad_exchange_2d_mod
      USE exchange_2d_mod
      USE mp_exchange_mod,   ONLY : ad_mp_exchange2d
      USE mp_exchange_mod,   ONLY : mp_exchange2d
      USE obc_volcons_mod
      USE ad_obc_volcons_mod
      USE ad_u2dbc_mod,      ONLY : ad_u2dbc_tile
      USE ad_v2dbc_mod,      ONLY : ad_v2dbc_tile
      USE ad_zetabc_mod,     ONLY : ad_zetabc_tile
!
!  Imported variable declarations.
!
      integer, intent(in) :: ng, tile
      integer, intent(in) :: LBi, UBi, LBj, UBj, UBk
      integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
      integer, intent(in) :: krhs, kstp, knew
      integer, intent(in) :: nstp, nnew
!
      real(r8), intent(in) :: pmask(LBi:,LBj:)
      real(r8), intent(in) :: rmask(LBi:,LBj:)
      real(r8), intent(in) :: umask(LBi:,LBj:)
      real(r8), intent(in) :: vmask(LBi:,LBj:)
      real(r8), intent(in) :: fomn(LBi:,LBj:)
      real(r8), intent(in) :: h(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) :: omn(LBi:,LBj:)
      real(r8), intent(in) :: pm(LBi:,LBj:)
      real(r8), intent(in) :: pn(LBi:,LBj:)
      real(r8), intent(in) :: dndx(LBi:,LBj:)
      real(r8), intent(in) :: dmde(LBi:,LBj:)
      real(r8), intent(in) :: pmon_r(LBi:,LBj:)
      real(r8), intent(in) :: pnom_r(LBi:,LBj:)
      real(r8), intent(in) :: pmon_p(LBi:,LBj:)
      real(r8), intent(in) :: pnom_p(LBi:,LBj:)
      real(r8), intent(in) :: om_r(LBi:,LBj:)
      real(r8), intent(in) :: on_r(LBi:,LBj:)
      real(r8), intent(in) :: om_p(LBi:,LBj:)
      real(r8), intent(in) :: on_p(LBi:,LBj:)
      real(r8), intent(in) :: visc2_p(LBi:,LBj:)
      real(r8), intent(in) :: visc2_r(LBi:,LBj:)
      real(r8), intent(in) :: rubar(LBi:,LBj:,:)
      real(r8), intent(in) :: rvbar(LBi:,LBj:,:)
      real(r8), intent(in) :: rzeta(LBi:,LBj:,:)
      real(r8), intent(in) :: ubar(LBi:,LBj:,:)
      real(r8), intent(in) :: vbar(LBi:,LBj:,:)
      real(r8), intent(in) :: zeta(LBi:,LBj:,:)
      real(r8), intent(in) :: Zt_avg1(LBi:,LBj:)
      real(r8), intent(inout) :: ad_DU_avg1(LBi:,LBj:)
      real(r8), intent(inout) :: ad_DU_avg2(LBi:,LBj:)
      real(r8), intent(inout) :: ad_DV_avg1(LBi:,LBj:)
      real(r8), intent(inout) :: ad_DV_avg2(LBi:,LBj:)
      real(r8), intent(inout) :: ad_Zt_avg1(LBi:,LBj:)
      real(r8), intent(inout) :: ad_rufrc(LBi:,LBj:)
      real(r8), intent(inout) :: ad_rvfrc(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_h(LBi:,LBj:)
      real(r8), intent(inout) :: ad_rubar(LBi:,LBj:,:)
      real(r8), intent(inout) :: ad_rvbar(LBi:,LBj:,:)
      real(r8), intent(inout) :: ad_rzeta(LBi:,LBj:,:)
      real(r8), intent(inout) :: ad_ubar(LBi:,LBj:,:)
      real(r8), intent(inout) :: ad_vbar(LBi:,LBj:,:)
      real(r8), intent(inout) :: ad_zeta(LBi:,LBj:,:)
!
!  Local variable declarations.
!
      logical :: CORRECTOR_2D_STEP
!
      integer :: i, is, j, ptsk
!
      real(r8) :: cff, cff1, cff2, cff3, cff4, cff5, cff6, cff7
      real(r8) :: fac, fac1, fac2, fac3
      real(r8) :: ad_cff, ad_cff1, ad_cff2, ad_cff3, ad_cff4
      real(r8) :: ad_fac, ad_fac1
      real(r8) :: adfac, adfac1, adfac2, adfac3, adfac4
!
      real(r8), parameter :: IniVal = 0.0_r8
!
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dgrad
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dnew
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Drhs
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Drhs_p
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Dstp
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: DUon
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: DVom
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: grad
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: gzeta
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: gzeta2
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rhs_ubar
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rhs_vbar
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: rhs_zeta
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: zeta_new
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: zwrk
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_Dgrad
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_Dnew
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_Drhs
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_Drhs_p
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_Dstp
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_DUon
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_DVom
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_UFe
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_UFx
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_VFe
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_VFx
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_grad
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_gzeta
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_gzeta2
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_rhs_ubar
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_rhs_vbar
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_rhs_zeta
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_zeta_new
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_zwrk
!
!-----------------------------------------------------------------------
!  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
!
      ptsk=3-kstp
      CORRECTOR_2D_STEP=.not.PREDICTOR_2D_STEP(ng)
!
!-----------------------------------------------------------------------
!  Initialize adjoint private variables.
!-----------------------------------------------------------------------
!
      ad_cff=IniVal
      ad_cff1=IniVal
      ad_cff2=IniVal
      ad_cff3=IniVal
      ad_cff4=IniVal
      ad_fac=IniVal
      ad_fac1=IniVal
      DO j=JminS,JmaxS
        DO i=IminS,ImaxS
          ad_Dgrad(i,j)=IniVal
          ad_Dnew(i,j)=IniVal
          ad_Drhs(i,j)=IniVal
          ad_Drhs_p(i,j)=IniVal
          ad_Dstp(i,j)=IniVal
          ad_DUon(i,j)=IniVal
          ad_DVom(i,j)=IniVal
          ad_UFe(i,j)=IniVal
          ad_UFx(i,j)=IniVal
          ad_VFe(i,j)=IniVal
          ad_VFx(i,j)=IniVal
          ad_grad(i,j)=IniVal
          ad_gzeta(i,j)=IniVal
          ad_gzeta2(i,j)=IniVal
          ad_rhs_ubar(i,j)=IniVal
          ad_rhs_vbar(i,j)=IniVal
          ad_rhs_zeta(i,j)=IniVal
          ad_zeta_new(i,j)=IniVal
          ad_zwrk(i,j)=IniVal
          ad_DUon(i,j)=IniVal
          ad_DVom(i,j)=IniVal
        END DO
      END DO
!
!-----------------------------------------------------------------------
!  Compute BASIC STATE total depth (m) arrays and vertically
!  integerated mass fluxes.
!-----------------------------------------------------------------------
!
!  In distributed-memory, the I- and J-ranges are different and a
!  special exchange is done to avoid having three ghost points for
!  high order numerical stencils. Notice that a private array is
!  passed below to the exchange routine. It also applies periodic
!  boundary conditions, if appropriate and no partitions in I- or
!  J-directions.
!
      DO j=JstrV-2,Jendp2
        DO i=IstrU-2,Iendp2
          Dnew(i,j)=zeta(i,j,knew)+h(i,j)
          Drhs(i,j)=zeta(i,j,krhs)+h(i,j)
          Dstp(i,j)=zeta(i,j,kstp)+h(i,j)
        END DO
      END DO
      DO j=JstrV-2,Jendp2
        DO i=IstrU-1,Iendp2
          cff=0.5_r8*on_u(i,j)
          cff1=cff*(Drhs(i,j)+Drhs(i-1,j))
          DUon(i,j)=ubar(i,j,krhs)*cff1
        END DO
      END DO
      DO j=JstrV-1,Jendp2
        DO i=IstrU-2,Iendp2
          cff=0.5_r8*om_v(i,j)
          cff1=cff*(Drhs(i,j)+Drhs(i,j-1))
          DVom(i,j)=vbar(i,j,krhs)*cff1
        END DO
      END DO
      IF (EWperiodic(ng).or.NSperiodic(ng)) THEN
        CALL exchange_u2d_tile (ng, tile,                               &
     &                          IminS, ImaxS, JminS, JmaxS,             &
     &                          DUon)
        CALL exchange_v2d_tile (ng, tile,                               &
     &                          IminS, ImaxS, JminS, JmaxS,             &
     &                          DVom)
      END IF
      CALL mp_exchange2d (ng, tile, iADM, 2,                            &
     &                    IminS, ImaxS, JminS, JmaxS,                   &
     &                    NghostPoints,                                 &
     &                    EWperiodic(ng), NSperiodic(ng),               &
     &                    DUon, DVom)
!
!  Compute integral mass flux across open boundaries and adjust
!  for volume conservation. Compute BASIC STATE value.
!  This must be computed here instead of below.
!
      IF (ANY(ad_VolCons(:,ng))) THEN
        CALL obc_flux_tile (ng, tile,                                   &
     &                      LBi, UBi, LBj, UBj,                         &
     &                      IminS, ImaxS, JminS, JmaxS,                 &
     &                      knew,                                       &
     &                      umask, vmask,                               &
     &                      h, om_v, on_u,                              &
     &                      ubar, vbar, zeta)
!
!  Set vertically integrated mass fluxes DUon and DVom along the open
!  boundaries in such a way that the integral volume is conserved.
!
        CALL set_DUV_bc_tile (ng, tile,                                 &
     &                        LBi, UBi, LBj, UBj,                       &
     &                        IminS, ImaxS, JminS, JmaxS,               &
     &                        krhs,                                     &
     &                        umask, vmask,                             &
     &                        om_v, on_u,                               &
     &                        ubar, vbar,                               &
     &                        Drhs, DUon, DVom)
      END IF
!
!  Compute BASIC state depths at PSI-points for viscosity.
!
      DO j=Jstr,Jend+1
        DO i=Istr,Iend+1
          Drhs_p(i,j)=0.25_r8*(Drhs(i,j  )+Drhs(i-1,j  )+               &
     &                         Drhs(i,j-1)+Drhs(i-1,j-1))
        END DO
      END DO
!
!  Initialize BASIC STATE right-hand-side terms.
!
      DO j=Jstr,Jend
        DO i=Istr,Iend
          rhs_ubar(i,j)=0.0_r8
          rhs_vbar(i,j)=0.0_r8
        END DO
      END DO
!
!  Do not perform the actual time stepping during the auxiliary
!  (nfast(ng)+1) time step.
!
      STEP_LOOP : IF (iif(ng).le.nfast(ng)) THEN
!
!-----------------------------------------------------------------------
!  Exchange boundary information.
!-----------------------------------------------------------------------
!
!^      CALL mp_exchange2d (ng, tile, iTLM, 2,                          &
!^   &                      LBi, UBi, LBj, UBj,                         &
!^   &                      NghostPoints,                               &
!^   &                      EWperiodic(ng), NSperiodic(ng),             &
!^   &                      tl_ubar(:,:,knew),                          &
!^   &                      tl_vbar(:,:,knew))
!^
        CALL ad_mp_exchange2d (ng, tile, iADM, 2,                       &
     &                         LBi, UBi, LBj, UBj,                      &
     &                         NghostPoints,                            &
     &                         EWperiodic(ng), NSperiodic(ng),          &
     &                         ad_ubar(:,:,knew),                       &
     &                         ad_vbar(:,:,knew))
!
        IF (EWperiodic(ng).or.NSperiodic(ng)) THEN
!^        CALL exchange_v2d_tile (ng, tile,                             &
!^   &                            LBi, UBi, LBj, UBj,                   &
!^   &                            tl_vbar(:,:,knew))
!^
          CALL ad_exchange_v2d_tile (ng, tile,                          &
     &                               LBi, UBi, LBj, UBj,                &
     &                               ad_vbar(:,:,knew))
!^        CALL exchange_u2d_tile (ng, tile,                             &
!^   &                            LBi, UBi, LBj, UBj,                   &
!^   &                            tl_ubar(:,:,knew))
!^
          CALL ad_exchange_u2d_tile (ng, tile,                          &
     &                               LBi, UBi, LBj, UBj,                &
     &                               ad_ubar(:,:,knew))
        END IF
!
!-----------------------------------------------------------------------
!  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
                cff=1.0_r8/(on_u(i,j)*                                  &
     &                      0.5_r8*(zeta(i-1,j,knew)+h(i-1,j)+          &
     &                              zeta(i  ,j,knew)+h(i  ,j)))
!^              tl_ubar(i,j,knew)=SOURCES(ng)%tl_Qbar(is)*cff+          &
!^   &                            SOURCES(ng)%Qbar(is)*tl_cff
!^
                SOURCES(ng)%ad_Qbar(is)=SOURCES(ng)%ad_Qbar(is)+        &
     &                                  cff*ad_ubar(i,j,knew)
                ad_cff=ad_cff+                                          &
     &                 SOURCES(ng)%Qbar(is)*ad_ubar(i,j,knew)
                ad_ubar(i,j,knew)=0.0_r8
!^              tl_cff=-cff*cff*on_u(i,j)*                              &
!^   &                 0.5_r8*(tl_zeta(i-1,j,knew)+tl_h(i-1,j)+         &
!^   &                         tl_zeta(i  ,j,knew)+tl_h(i  ,j))
!^
                adfac=-cff*cff*on_u(i,j)*0.5_r8*ad_cff
                ad_h(i-1,j)=ad_h(i-1,j)+adfac
                ad_h(i  ,j)=ad_h(i  ,j)+adfac
                ad_zeta(i-1,j,knew)=ad_zeta(i-1,j,knew)+adfac
                ad_zeta(i  ,j,knew)=ad_zeta(i  ,j,knew)+adfac
                ad_cff=0.0_r8
              ELSE IF (INT(SOURCES(ng)%Dsrc(is)).eq.1) THEN
                cff=1.0_r8/(om_v(i,j)*                                  &
     &                      0.5_r8*(zeta(i,j-1,knew)+h(i,j-1)+          &
     &                              zeta(i,j  ,knew)+h(i,j  )))
!^              tl_vbar(i,j,knew)=SOURCES(ng)%tl_Qbar(is)*cff+          &
!^   &                            SOURCES(ng)%Qbar(is)*tl_cff
!^
                SOURCES(ng)%ad_Qbar(is)=SOURCES(ng)%ad_Qbar(is)+        &
     &                                  cff*ad_vbar(i,j,knew)
                ad_cff=ad_cff+                                          &
     &                 SOURCES(ng)%Qbar(is)*ad_vbar(i,j,knew)
                ad_vbar(i,j,knew)=0.0_r8
!^              tl_cff=-cff*cff*om_v(i,j)*                              &
!^   &                 0.5_r8*(tl_zeta(i,j-1,knew)+tl_h(i,j-1)+         &
!^   &                         tl_zeta(i,j  ,knew)+tl_h(i,j  ))
!^
                adfac=-cff*cff*om_v(i,j)*0.5_r8*ad_cff
                ad_h(i,j-1)=ad_h(i,j-1)+adfac
                ad_h(i,j  )=ad_h(i,j  )+adfac
                ad_zeta(i,j-1,knew)=ad_zeta(i,j-1,knew)+adfac
                ad_zeta(i,j  ,knew)=ad_zeta(i,j  ,knew)+adfac
                ad_cff=0.0_r8
              END IF
            END IF
          END DO
        END IF
!
!-----------------------------------------------------------------------
!  Apply adjoint lateral boundary conditions.
!-----------------------------------------------------------------------
!
!  Compute integral mass flux across open boundaries and adjust
!  for adjoint volume conservation.
!
        IF (ANY(ad_VolCons(:,ng))) THEN
!^        CALL tl_obc_flux_tile (ng, tile,                              &
!^   &                           LBi, UBi, LBj, UBj,                    &
!^   &                           IminS, ImaxS, JminS, JmaxS,            &
!^   &                           knew,                                  &
!^   &                           umask, vmask,                          &
!^   &                           h, tl_h, om_v, on_u,                   &
!^   &                           ubar, vbar, zeta,                      &
!^   &                           tl_ubar, tl_vbar, tl_zeta)
!^
          CALL ad_obc_flux_tile (ng, tile,                              &
     &                           LBi, UBi, LBj, UBj,                    &
     &                           IminS, ImaxS, JminS, JmaxS,            &
     &                           knew,                                  &
     &                           umask, vmask,                          &
     &                           h, ad_h, om_v, on_u,                   &
     &                           ubar, vbar, zeta,                      &
     &                           ad_ubar, ad_vbar, ad_zeta)
        END IF
!
!  Adjoint lateral boundary conditons.
!
!^      CALL tl_v2dbc_tile (ng, tile,                                   &
!^   &                      LBi, UBi, LBj, UBj,                         &
!^   &                      IminS, ImaxS, JminS, JmaxS,                 &
!^   &                      krhs, kstp, knew,                           &
!^   &                      ubar, vbar, zeta,                           &
!^   &                      tl_ubar, tl_vbar, tl_zeta)
!^
        CALL ad_v2dbc_tile (ng, tile,                                   &
     &                      LBi, UBi, LBj, UBj,                         &
     &                      IminS, ImaxS, JminS, JmaxS,                 &
     &                      krhs, kstp, knew,                           &
     &                      ubar, vbar, zeta,                           &
     &                      ad_ubar, ad_vbar, ad_zeta)
!^      CALL tl_u2dbc_tile (ng, tile,                                   &
!^   &                      LBi, UBi, LBj, UBj,                         &
!^   &                      IminS, ImaxS, JminS, JmaxS,                 &
!^   &                      krhs, kstp, knew,                           &
!^   &                      ubar, vbar, zeta,                           &
!^   &                      tl_ubar, tl_vbar, tl_zeta)
!^
        CALL ad_u2dbc_tile (ng, tile,                                   &
     &                      LBi, UBi, LBj, UBj,                         &
     &                      IminS, ImaxS, JminS, JmaxS,                 &
     &                      krhs, kstp, knew,                           &
     &                      ubar, vbar, zeta,                           &
     &                      ad_ubar, ad_vbar, ad_zeta)
!
!  If predictor step, load right-side-term into shared arrays for
!  future use during the subsequent corrector step.
!
        IF (PREDICTOR_2D_STEP(ng)) THEN
          DO j=JstrV,Jend
            DO i=Istr,Iend
!^            tl_rvbar(i,j,krhs)=tl_rhs_vbar(i,j)
!^
              ad_rhs_vbar(i,j)=ad_rhs_vbar(i,j)+ad_rvbar(i,j,krhs)
              ad_rvbar(i,j,krhs)=0.0_r8
            END DO
          END DO
          DO j=Jstr,Jend
            DO i=IstrU,Iend
!^            tl_rubar(i,j,krhs)=tl_rhs_ubar(i,j)
!^
              ad_rhs_ubar(i,j)=ad_rhs_ubar(i,j)+ad_rubar(i,j,krhs)
              ad_rubar(i,j,krhs)=0.0_r8
            END DO
          END DO
        END IF
!
!=======================================================================
!  Time step adjoint 2D momentum equations.
!=======================================================================
!
!  During the first time-step, the predictor step is Forward-Euler
!  and the corrector step is Backward-Euler. Otherwise, the predictor
!  step is Leap-frog and the corrector step is Adams-Moulton.
!
        IF (iif(ng).eq.1) THEN
          cff1=0.5_r8*dtfast(ng)
          DO j=JstrV,Jend
            DO i=Istr,Iend
              cff=(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))
              fac=1.0_r8/(Dnew(i,j)+Dnew(i,j-1))
!^            tl_vbar(i,j,knew)=tl_vbar(i,j,knew)*vmask(i,j)
!^
              ad_vbar(i,j,knew)=ad_vbar(i,j,knew)*vmask(i,j)
!^            tl_vbar(i,j,knew)=(tl_vbar(i,j,kstp)*                     &
!^   &                           (Dstp(i,j)+Dstp(i,j-1))+               &
!^   &                           vbar(i,j,kstp)*                        &
!^   &                           (tl_Dstp(i,j)+tl_Dstp(i,j-1))+         &
!^   &                           cff*cff1*tl_rhs_vbar(i,j))*fac+        &
!^   &                          (vbar(i,j,kstp)*                        &
!^   &                           (Dstp(i,j)+Dstp(i,j-1))+               &
!^   &                           cff*cff1*rhs_vbar(i,j))*tl_fac
!^
              adfac=fac*ad_vbar(i,j,knew)
              adfac1=adfac*(Dstp(i,j)+Dstp(i,j-1))
              adfac2=adfac*cff*cff1
              adfac3=adfac*vbar(i,j,kstp)
              ad_vbar(i,j,kstp)=ad_vbar(i,j,kstp)+adfac1
              ad_rhs_vbar(i,j)=ad_rhs_vbar(i,j)+adfac2
              ad_Dstp(i,j-1)=ad_Dstp(i,j-1)+adfac3
              ad_Dstp(i,j  )=ad_Dstp(i,j  )+adfac3
              ad_fac=ad_fac+                                            &
     &               (vbar(i,j,kstp)*(Dstp(i,j)+Dstp(i,j-1))+           &
     &                cff*cff1*rhs_vbar(i,j))*ad_vbar(i,j,knew)
              ad_vbar(i,j,knew)=0.0_r8
!^            tl_fac=-fac*fac*(tl_Dnew(i,j)+tl_Dnew(i,j-1))
!^
              adfac=-fac*fac*ad_fac
              ad_Dnew(i,j-1)=ad_Dnew(i,j-1)+adfac
              ad_Dnew(i,j  )=ad_Dnew(i,j  )+adfac
              ad_fac=0.0_r8
            END DO
          END DO
          DO j=Jstr,Jend
            DO i=IstrU,Iend
              cff=(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j))
              fac=1.0_r8/(Dnew(i,j)+Dnew(i-1,j))
!^            tl_ubar(i,j,knew)=tl_ubar(i,j,knew)*umask(i,j)
!^
              ad_ubar(i,j,knew)=ad_ubar(i,j,knew)*umask(i,j)
!^            tl_ubar(i,j,knew)=(tl_ubar(i,j,kstp)*                     &
!^   &                           (Dstp(i,j)+Dstp(i-1,j))+               &
!^   &                           ubar(i,j,kstp)*                        &
!^   &                           (tl_Dstp(i,j)+tl_Dstp(i-1,j))+         &
!^   &                           cff*cff1*tl_rhs_ubar(i,j))*fac+        &
!^   &                          (ubar(i,j,kstp)*                        &
!^   &                           (Dstp(i,j)+Dstp(i-1,j))+               &
!^   &                           cff*cff1*rhs_ubar(i,j))*tl_fac
!^
              adfac=fac*ad_ubar(i,j,knew)
              adfac1=adfac*(Dstp(i,j)+Dstp(i-1,j))
              adfac2=adfac*cff*cff1
              adfac3=adfac*ubar(i,j,kstp)
              ad_ubar(i,j,kstp)=ad_ubar(i,j,kstp)+adfac1
              ad_rhs_ubar(i,j)=ad_rhs_ubar(i,j)+adfac2
              ad_Dstp(i-1,j)=ad_Dstp(i-1,j)+adfac3
              ad_Dstp(i  ,j)=ad_Dstp(i  ,j)+adfac3
              ad_fac=ad_fac+                                            &
     &               (ubar(i,j,kstp)*(Dstp(i,j)+Dstp(i-1,j))+           &
     &                cff*cff1*rhs_ubar(i,j))*ad_ubar(i,j,knew)
              ad_ubar(i,j,knew)=0.0_r8
!^            tl_fac=-fac*fac*(tl_Dnew(i,j)+tl_Dnew(i-1,j))
!^
              adfac=-fac*fac*ad_fac
              ad_Dnew(i-1,j)=ad_Dnew(i-1,j)+adfac
              ad_Dnew(i  ,j)=ad_Dnew(i  ,j)+adfac
              ad_fac=0.0_r8
            END DO
          END DO
        ELSE IF (PREDICTOR_2D_STEP(ng)) THEN
          cff1=dtfast(ng)
          DO j=JstrV,Jend
            DO i=Istr,Iend
              cff=(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))
              fac=1.0_r8/(Dnew(i,j)+Dnew(i,j-1))
!^            tl_vbar(i,j,knew)=tl_vbar(i,j,knew)*vmask(i,j)
!^
              ad_vbar(i,j,knew)=ad_vbar(i,j,knew)*vmask(i,j)
!^            tl_vbar(i,j,knew)=(tl_vbar(i,j,kstp)*                     &
!^   &                           (Dstp(i,j)+Dstp(i,j-1))+               &
!^   &                           vbar(i,j,kstp)*                        &
!^   &                           (tl_Dstp(i,j)+tl_Dstp(i,j-1))+         &
!^   &                           cff*cff1*tl_rhs_vbar(i,j))*fac+        &
!^   &                          (vbar(i,j,kstp)*                        &
!^   &                           (Dstp(i,j)+Dstp(i,j-1))+               &
!^   &                           cff*cff1*rhs_vbar(i,j))*tl_fac
!^
              adfac=fac*ad_vbar(i,j,knew)
              adfac1=adfac*(Dstp(i,j)+Dstp(i,j-1))
              adfac2=adfac*cff*cff1
              adfac3=adfac*vbar(i,j,kstp)
              ad_vbar(i,j,kstp)=ad_vbar(i,j,kstp)+adfac1
              ad_rhs_vbar(i,j)=ad_rhs_vbar(i,j)+adfac2
              ad_Dstp(i,j-1)=ad_Dstp(i,j-1)+adfac3
              ad_Dstp(i,j  )=ad_Dstp(i,j  )+adfac3
              ad_fac=ad_fac+                                            &
     &               (vbar(i,j,kstp)*(Dstp(i,j)+Dstp(i,j-1))+           &
     &                cff*cff1*rhs_vbar(i,j))*ad_vbar(i,j,knew)
              ad_vbar(i,j,knew)=0.0_r8
!^            tl_fac=-fac*fac*(tl_Dnew(i,j)+tl_Dnew(i,j-1))
!^
              adfac=-fac*fac*ad_fac
              ad_Dnew(i,j-1)=ad_Dnew(i,j-1)+adfac
              ad_Dnew(i,j  )=ad_Dnew(i,j  )+adfac
              ad_fac=0.0_r8
            END DO
          END DO
          DO j=Jstr,Jend
            DO i=IstrU,Iend
              cff=(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j))
              fac=1.0_r8/(Dnew(i,j)+Dnew(i-1,j))
!^            tl_ubar(i,j,knew)=tl_ubar(i,j,knew)*umask(i,j)
!^
              ad_ubar(i,j,knew)=ad_ubar(i,j,knew)*umask(i,j)
!^            tl_ubar(i,j,knew)=(tl_ubar(i,j,kstp)*                     &
!^   &                           (Dstp(i,j)+Dstp(i-1,j))+               &
!^   &                           ubar(i,j,kstp)*                        &
!^   &                           (tl_Dstp(i,j)+tl_Dstp(i-1,j))+         &
!^   &                           cff*cff1*tl_rhs_ubar(i,j))*fac+        &
!^   &                          (ubar(i,j,kstp)*                        &
!^   &                           (Dstp(i,j)+Dstp(i-1,j))+               &
!^   &                           cff*cff1*rhs_ubar(i,j))*tl_fac
!^
              adfac=fac*ad_ubar(i,j,knew)
              adfac1=adfac*(Dstp(i,j)+Dstp(i-1,j))
              adfac2=adfac*cff*cff1
              adfac3=adfac*ubar(i,j,kstp)
              ad_ubar(i,j,kstp)=ad_ubar(i,j,kstp)+adfac1
              ad_rhs_ubar(i,j)=ad_rhs_ubar(i,j)+adfac2
              ad_Dstp(i-1,j)=ad_Dstp(i-1,j)+adfac3
              ad_Dstp(i  ,j)=ad_Dstp(i  ,j)+adfac3
              ad_fac=ad_fac+                                            &
     &               (ubar(i,j,kstp)*(Dstp(i,j)+Dstp(i-1,j))+           &
     &                cff*cff1*rhs_ubar(i,j))*ad_ubar(i,j,knew)
              ad_ubar(i,j,knew)=0.0_r8
!^            tl_fac=-fac*fac*(tl_Dnew(i,j)+tl_Dnew(i-1,j))
!^
              adfac=-fac*fac*ad_fac
              ad_Dnew(i-1,j)=ad_Dnew(i-1,j)+adfac
              ad_Dnew(i  ,j)=ad_Dnew(i  ,j)+adfac
              ad_fac=0.0_r8
            END DO
          END DO
        ELSE IF (CORRECTOR_2D_STEP) THEN
          cff1=0.5_r8*dtfast(ng)*5.0_r8/12.0_r8
          cff2=0.5_r8*dtfast(ng)*8.0_r8/12.0_r8
          cff3=0.5_r8*dtfast(ng)*1.0_r8/12.0_r8
          DO j=JstrV,Jend
            DO i=Istr,Iend
              cff=(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))
              fac=1.0_r8/(Dnew(i,j)+Dnew(i,j-1))
!^            tl_vbar(i,j,knew)=tl_vbar(i,j,knew)*vmask(i,j)
!^
              ad_vbar(i,j,knew)=ad_vbar(i,j,knew)*vmask(i,j)
!^            tl_vbar(i,j,knew)=(tl_vbar(i,j,kstp)*                     &
!^   &                           (Dstp(i,j)+Dstp(i,j-1))+               &
!^   &                           vbar(i,j,kstp)*                        &
!^   &                           (tl_Dstp(i,j)+tl_Dstp(i,j-1))+         &
!^   &                           cff*(cff1*tl_rhs_vbar(i,j)+            &
!^   &                                cff2*tl_rvbar(i,j,kstp)-          &
!^   &                                cff3*tl_rvbar(i,j,ptsk)))*fac+    &
!^   &                          (vbar(i,j,kstp)*                        &
!^   &                           (Dstp(i,j)+Dstp(i,j-1))+               &
!^   &                           cff*(cff1*rhs_vbar(i,j)+               &
!^   &                                cff2*rvbar(i,j,kstp)-             &
!^   &                                cff3*rvbar(i,j,ptsk)))*tl_fac
!^
              adfac=fac*ad_vbar(i,j,knew)
              adfac1=adfac*(Dstp(i,j)+Dstp(i,j-1))
              adfac2=adfac*cff
              adfac3=adfac*vbar(i,j,kstp)
              ad_vbar(i,j,kstp)=ad_vbar(i,j,kstp)+adfac1
              ad_rhs_vbar(i,j)=ad_rhs_vbar(i,j)+cff1*adfac2
              ad_rvbar(i,j,kstp)=ad_rvbar(i,j,kstp)+cff2*adfac2
              ad_rvbar(i,j,ptsk)=-cff3*adfac2
              ad_Dstp(i,j-1)=ad_Dstp(i,j-1)+adfac3
              ad_Dstp(i,j  )=ad_Dstp(i,j  )+adfac3
              ad_fac=ad_fac+                                            &
     &               (vbar(i,j,kstp)*(Dstp(i,j)+Dstp(i,j-1))+           &
     &                cff*(cff1*rhs_vbar(i,j)+                          &
     &                     cff2*rvbar(i,j,kstp)-                        &
     &                     cff3*rvbar(i,j,ptsk)))*ad_vbar(i,j,knew)
              ad_vbar(i,j,knew)=0.0_r8
!^            tl_fac=-fac*fac*(tl_Dnew(i,j)+tl_Dnew(i,j-1))
!^
              adfac=-fac*fac*ad_fac
              ad_Dnew(i,j-1)=ad_Dnew(i,j-1)+adfac
              ad_Dnew(i,j  )=ad_Dnew(i,j  )+adfac
              ad_fac=0.0_r8
            END DO
          END DO
          DO j=Jstr,Jend
            DO i=IstrU,Iend
              cff=(pm(i,j)+pm(i-1,j))*(pn(i,j)+pn(i-1,j))
              fac=1.0_r8/(Dnew(i,j)+Dnew(i-1,j))
!^            tl_ubar(i,j,knew)=tl_ubar(i,j,knew)*umask(i,j)
!^
              ad_ubar(i,j,knew)=ad_ubar(i,j,knew)*umask(i,j)
!^            tl_ubar(i,j,knew)=(tl_ubar(i,j,kstp)*                     &
!^   &                           (Dstp(i,j)+Dstp(i-1,j))+               &
!^   &                           ubar(i,j,kstp)*                        &
!^   &                           (tl_Dstp(i,j)+tl_Dstp(i-1,j))+         &
!^   &                           cff*(cff1*tl_rhs_ubar(i,j)+            &
!^   &                                cff2*tl_rubar(i,j,kstp)-          &
!^   &                                cff3*tl_rubar(i,j,ptsk)))*fac+    &
!^   &                          (ubar(i,j,kstp)*                        &
!^   &                           (Dstp(i,j)+Dstp(i-1,j))+               &
!^   &                           cff*(cff1*rhs_ubar(i,j)+               &
!^   &                                cff2*rubar(i,j,kstp)-             &
!^   &                                cff3*rubar(i,j,ptsk)))*tl_fac
!^
              adfac=fac*ad_ubar(i,j,knew)
              adfac1=adfac*(Dstp(i,j)+Dstp(i-1,j))
              adfac2=adfac*cff
              adfac3=adfac*ubar(i,j,kstp)
              ad_ubar(i,j,kstp)=ad_ubar(i,j,kstp)+adfac1
              ad_rhs_ubar(i,j)=ad_rhs_ubar(i,j)+cff1*adfac2
              ad_rubar(i,j,kstp)=ad_rubar(i,j,kstp)+cff2*adfac2
              ad_rubar(i,j,ptsk)=-cff3*adfac2
              ad_Dstp(i-1,j)=ad_Dstp(i-1,j)+adfac3
              ad_Dstp(i  ,j)=ad_Dstp(i  ,j)+adfac3
              ad_fac=ad_fac+                                            &
     &               (ubar(i,j,kstp)*(Dstp(i,j)+Dstp(i-1,j))+           &
     &                cff*(cff1*rhs_ubar(i,j)+                          &
     &                     cff2*rubar(i,j,kstp)-                        &
     &                     cff3*rubar(i,j,ptsk)))*ad_ubar(i,j,knew)
              ad_ubar(i,j,knew)=0.0_r8
!^            tl_fac=-fac*fac*(tl_Dnew(i,j)+tl_Dnew(i-1,j))
!^
              adfac=-fac*fac*ad_fac
              ad_Dnew(i-1,j)=ad_Dnew(i-1,j)+adfac
              ad_Dnew(i  ,j)=ad_Dnew(i  ,j)+adfac
              ad_fac=0.0_r8
            END DO
          END DO
        END IF
!
!  Compute adjoint total water column depth.
!
        DO j=JstrV-1,Jend
          DO i=IstrU-1,Iend
!^          tl_Dstp(i,j)=tl_zeta(i,j,kstp)+tl_h(i,j)
!^
            ad_zeta(i,j,kstp)=ad_zeta(i,j,kstp)+ad_Dstp(i,j)
            ad_h(i,j)=ad_h(i,j)+ad_Dstp(i,j)
            ad_Dstp(i,j)=0.0_r8
          END DO
        END DO
!
!=======================================================================
!  Compute right-hand-side for the 2D momentum equations.
!=======================================================================
!
!-----------------------------------------------------------------------
!  Adjoint Coupling between 2D and 3D equations.
!-----------------------------------------------------------------------
!
!  Before the predictor step of the first barotropic time-step,
!  arrays "rufrc" and "rvfrc" contain the vertical integrals of
!  the 3D right-hand-side terms for momentum equations (including
!  surface and bottom stresses, if so prescribed).
!
!  Convert them into forcing terms by subtracting the fast time
!  "rhs_ubar" and "rhs_vbar" from them; Also, immediately apply
!  these forcing terms "rhs_ubar" and "rhs_vbar".
!
!  From now on, these newly computed forcing terms will remain
!  constant during the fast time stepping and will added to
!  "rhs_ubar" and "rhs_vbar" during all subsequent time steps.
!
        IF (iif(ng).eq.1.and.PREDICTOR_2D_STEP(ng)) THEN
          IF (iic(ng).eq.ntfirst(ng)) THEN
            DO j=JstrV,Jend
              DO i=Istr,Iend
!^              tl_rv(i,j,0,nstp)=tl_rvfrc(i,j)
!^
                ad_rvfrc(i,j)=ad_rvfrc(i,j)+ad_rv(i,j,0,nstp)
                ad_rv(i,j,0,nstp)=0.0_r8
!^              tl_rhs_vbar(i,j)=tl_rhs_vbar(i,j)+tl_rvfrc(i,j)
!^
                ad_rvfrc(i,j)=ad_rvfrc(i,j)+ad_rhs_vbar(i,j)
!^              tl_rvfrc(i,j)=tl_rvfrc(i,j)-tl_rhs_vbar(i,j)
!^
                ad_rhs_vbar(i,j)=ad_rhs_vbar(i,j)-ad_rvfrc(i,j)
              END DO
            END DO
            DO j=Jstr,Jend
              DO i=IstrU,Iend
!^              tl_ru(i,j,0,nstp)=tl_rufrc(i,j)
!^
                ad_rufrc(i,j)=ad_rufrc(i,j)+ad_ru(i,j,0,nstp)
                ad_ru(i,j,0,nstp)=0.0_r8
!^              tl_rhs_ubar(i,j)=tl_rhs_ubar(i,j)+tl_rufrc(i,j)
!^
                ad_rufrc(i,j)=ad_rufrc(i,j)+ad_rhs_ubar(i,j)
!^              tl_rufrc(i,j)=tl_rufrc(i,j)-tl_rhs_ubar(i,j)
!^
                ad_rhs_ubar(i,j)=ad_rhs_ubar(i,j)-ad_rufrc(i,j)
              END DO
            END DO
          ELSE IF (iic(ng).eq.(ntfirst(ng)+1)) THEN
            DO j=JstrV,Jend
              DO i=Istr,Iend
!^              tl_rv(i,j,0,nstp)=tl_rvfrc(i,j)
!^
                ad_rvfrc(i,j)=ad_rvfrc(i,j)+ad_rv(i,j,0,nstp)
                ad_rv(i,j,0,nstp)=0.0_r8
!^              tl_rhs_vbar(i,j)=tl_rhs_vbar(i,j)+                      &
!^   &                           1.5_r8*tl_rvfrc(i,j)-                  &
!^   &                           0.5_r8*tl_rv(i,j,0,nnew)
!^
                ad_rvfrc(i,j)=ad_rvfrc(i,j)+1.5_r8*ad_rhs_vbar(i,j)
                ad_rv(i,j,0,nnew)=ad_rv(i,j,0,nnew)-                    &
     &                            0.5_r8*ad_rhs_vbar(i,j)
!^              tl_rvfrc(i,j)=tl_rvfrc(i,j)-tl_rhs_vbar(i,j)
!^
                ad_rhs_vbar(i,j)=ad_rhs_vbar(i,j)-ad_rvfrc(i,j)
              END DO
            END DO
            DO j=Jstr,Jend
              DO i=IstrU,Iend
!^              tl_ru(i,j,0,nstp)=tl_rufrc(i,j)
!^
                ad_rufrc(i,j)=ad_rufrc(i,j)+ad_ru(i,j,0,nstp)
                ad_ru(i,j,0,nstp)=0.0_r8
!^              tl_rhs_ubar(i,j)=tl_rhs_ubar(i,j)+                      &
!^   &                           1.5_r8*tl_rufrc(i,j)-                  &
!^   &                           0.5_r8*tl_ru(i,j,0,nnew)
!^
                ad_rufrc(i,j)=ad_rufrc(i,j)+1.5_r8*ad_rhs_ubar(i,j)
                ad_ru(i,j,0,nnew)=ad_ru(i,j,0,nnew)-                    &
     &                            0.5_r8*ad_rhs_ubar(i,j)
!^              tl_rufrc(i,j)=tl_rufrc(i,j)-tl_rhs_ubar(i,j)
!^
                ad_rhs_ubar(i,j)=ad_rhs_ubar(i,j)-ad_rufrc(i,j)
              END DO
            END DO
          ELSE
            cff1=23.0_r8/12.0_r8
            cff2=16.0_r8/12.0_r8
            cff3= 5.0_r8/12.0_r8
            DO j=JstrV,Jend
              DO i=Istr,Iend
!^              tl_rv(i,j,0,nstp)=tl_rvfrc(i,j)
!^
                ad_rvfrc(i,j)=ad_rvfrc(i,j)+ad_rv(i,j,0,nstp)
                ad_rv(i,j,0,nstp)=0.0_r8
!^              tl_rhs_vbar(i,j)=tl_rhs_vbar(i,j)+                      &
!^   &                           cff1*tl_rvfrc(i,j)-                    &
!^   &                           cff2*tl_rv(i,j,0,nnew)+                &
!^   &                           cff3*tl_rv(i,j,0,nstp)
!^
                ad_rvfrc(i,j)=ad_rvfrc(i,j)+cff1*ad_rhs_vbar(i,j)
                ad_rv(i,j,0,nnew)=ad_rv(i,j,0,nnew)-                    &
     &                            cff2*ad_rhs_vbar(i,j)
                ad_rv(i,j,0,nstp)=ad_rv(i,j,0,nstp)+                    &
     &                            cff3*ad_rhs_vbar(i,j)
!^              tl_rvfrc(i,j)=tl_rvfrc(i,j)-tl_rhs_vbar(i,j)
!^
                ad_rhs_vbar(i,j)=ad_rhs_vbar(i,j)-ad_rvfrc(i,j)
              END DO
            END DO
            DO j=Jstr,Jend
              DO i=IstrU,Iend
!^              tl_ru(i,j,0,nstp)=tl_rufrc(i,j)
!^
                ad_rufrc(i,j)=ad_rufrc(i,j)+ad_ru(i,j,0,nstp)
                ad_ru(i,j,0,nstp)=0.0_r8
!^              tl_rhs_ubar(i,j)=tl_rhs_ubar(i,j)+                      &
!^   &                           cff1*tl_rufrc(i,j)-                    &
!^   &                           cff2*tl_ru(i,j,0,nnew)+                &
!^   &                           cff3*tl_ru(i,j,0,nstp)
!^
                ad_rufrc(i,j)=ad_rufrc(i,j)+cff1*ad_rhs_ubar(i,j)
                ad_ru(i,j,0,nnew)=ad_ru(i,j,0,nnew)-                    &
     &                            cff2*ad_rhs_ubar(i,j)
                ad_ru(i,j,0,nstp)=ad_ru(i,j,0,nstp)+                    &
     &                            cff3*ad_rhs_ubar(i,j)
!^              tl_rufrc(i,j)=tl_rufrc(i,j)-tl_rhs_ubar(i,j)
!^
                ad_rhs_ubar(i,j)=ad_rhs_ubar(i,j)-ad_rufrc(i,j)
              END DO
            END DO
          END IF
        ELSE
          DO j=JstrV,Jend
            DO i=Istr,Iend
!^            tl_rhs_vbar(i,j)=tl_rhs_vbar(i,j)+tl_rvfrc(i,j)
!^
              ad_rvfrc(i,j)=ad_rvfrc(i,j)+ad_rhs_vbar(i,j)
            END DO
          END DO
          DO j=Jstr,Jend
            DO i=IstrU,Iend
!^            tl_rhs_ubar(i,j)=tl_rhs_ubar(i,j)+tl_rufrc(i,j)
!^
              ad_rufrc(i,j)=ad_rufrc(i,j)+ad_rhs_ubar(i,j)
            END DO
          END DO
        END IF
!
!-----------------------------------------------------------------------
!  Add in adjoint nudging of 2D momentum climatology.
!-----------------------------------------------------------------------
!
        IF (LnudgeM2CLM(ng)) THEN
          DO j=JstrV,Jend
            DO i=Istr,Iend
              cff=0.25_r8*(CLIMA(ng)%M2nudgcof(i,j-1)+                  &
     &                     CLIMA(ng)%M2nudgcof(i,j  ))*                 &
     &            om_v(i,j)*on_v(i,j)
!^            tl_rhs_vbar(i,j)=tl_rhs_vbar(i,j)+                        &
!^   &                         cff*((Drhs(i,j-1)+Drhs(i,j))*            &
!^   &                              (-tl_vbar(i,j,krhs))+               &
!^   &                              (tl_Drhs(i,j-1)+tl_Drhs(i,j))*      &
!^   &                              (CLIMA(ng)%vbarclm(i,j)-
!^   &                               vbar(i,j,krhs)))
!^
              adfac=cff*ad_rhs_vbar(i,j)
              adfac1=adfac*(Drhs(i,j-1)+Drhs(i,j))
              adfac2=adfac*(CLIMA(ng)%vbarclm(i,j)-vbar(i,j,krhs))
              ad_vbar(i,j,krhs)=ad_vbar(i,j,krhs)-adfac1
              ad_Drhs(i,j-1)=ad_Drhs(i,j-1)+adfac2
              ad_Drhs(i,j  )=ad_Drhs(i,j  )+adfac2
            END DO
          END DO
          DO j=Jstr,Jend
            DO i=IstrU,Iend
              cff=0.25_r8*(CLIMA(ng)%M2nudgcof(i-1,j)+                  &
     &                     CLIMA(ng)%M2nudgcof(i  ,j))*                 &
     &            om_u(i,j)*on_u(i,j)
!^            tl_rhs_ubar(i,j)=tl_rhs_ubar(i,j)+                        &
!^   &                         cff*((Drhs(i-1,j)+Drhs(i,j))*            &
!^   &                              (-tl_ubar(i,j,krhs))+               &
!^   &                              (tl_Drhs(i-1,j)+tl_Drhs(i,j))*      &
!^   &                              (CLIMA(ng)%ubarclm(i,j)-
!^   &                               ubar(i,j,krhs)))
!^
              adfac=cff*ad_rhs_ubar(i,j)
              adfac1=adfac*(Drhs(i-1,j)+Drhs(i,j))
              adfac2=adfac*(CLIMA(ng)%ubarclm(i,j)-ubar(i,j,krhs))
              ad_ubar(i,j,krhs)=ad_ubar(i,j,krhs)-adfac1
              ad_Drhs(i-1,j)=ad_Drhs(i-1,j)+adfac2
              ad_Drhs(i  ,j)=ad_Drhs(i  ,j)+adfac2
            END DO
          END DO
        END IF
!
!-----------------------------------------------------------------------
!  Compute basic state total depth at PSI-points.
!-----------------------------------------------------------------------
!
        DO j=Jstr,Jend+1
          DO i=Istr,Iend+1
            Drhs_p(i,j)=0.25_r8*(Drhs(i,j  )+Drhs(i-1,j  )+             &
     &                           Drhs(i,j-1)+Drhs(i-1,j-1))
          END DO
        END DO
!
!-----------------------------------------------------------------------
!  Add in adjoint horizontal harmonic viscosity.
!-----------------------------------------------------------------------
!
!  Add in harmonic viscosity.
!
        DO j=JstrV,Jend
          DO i=Istr,Iend
!^          tl_rhs_vbar(i,j)=tl_rhs_vbar(i,j)+tl_fac
!^
            ad_fac=ad_fac+ad_rhs_vbar(i,j)
!^          tl_fac=tl_cff1-tl_cff2
!^
            ad_cff1=ad_cff1+ad_fac
            ad_cff2=ad_cff2-ad_fac
            ad_fac=0.0_r8
!^          tl_cff2=0.5_r8*(pm(i,j-1)+pm(i,j))*                         &
!^   &              (tl_VFe(i  ,j)-tl_VFe(i,j-1))
!^
            adfac=0.5_r8*(pm(i,j-1)+pm(i,j))*ad_cff2
            ad_VFe(i,j-1)=ad_VFe(i,j-1)-adfac
            ad_VFe(i  ,j)=ad_VFe(i  ,j)+adfac
            ad_cff2=0.0_r8
!^          tl_cff1=0.5_r8*(pn(i,j-1)+pn(i,j))*                         &
!^   &              (tl_VFx(i+1,j)-tl_VFx(i,j  ))
!^
            adfac=0.5_r8*(pn(i,j-1)+pn(i,j))*ad_cff1
            ad_VFx(i  ,j)=ad_VFx(i  ,j)-adfac
            ad_VFx(i+1,j)=ad_VFx(i+1,j)+adfac
            ad_cff1=0.0_r8
          END DO
        END DO
        DO j=Jstr,Jend
          DO i=IstrU,Iend
!^          tl_rhs_ubar(i,j)=tl_rhs_ubar(i,j)+tl_fac
!^
            ad_fac=ad_fac+ad_rhs_ubar(i,j)
!^          tl_fac=tl_cff1+tl_cff2
!^
            ad_cff1=ad_cff1+ad_fac
            ad_cff2=ad_cff2+ad_fac
            ad_fac=0.0_r8
!^          tl_cff2=0.5_r8*(pm(i-1,j)+pm(i,j))*                         &
!^   &              (tl_UFe(i,j+1)-tl_UFe(i  ,j))
!^
            adfac=0.5_r8*(pm(i-1,j)+pm(i,j))*ad_cff2
            ad_UFe(i,j  )=ad_UFe(i,j  )-adfac
            ad_UFe(i,j+1)=ad_UFe(i,j+1)+adfac
            ad_cff2=0.0_r8
!^          tl_cff1=0.5_r8*(pn(i-1,j)+pn(i,j))*                         &
!^   &              (tl_UFx(i,j  )-tl_UFx(i-1,j))
!^
            adfac=0.5_r8*(pn(i-1,j)+pn(i,j))*ad_cff1
            ad_UFx(i-1,j)=ad_UFx(i-1,j)-adfac
            ad_UFx(i  ,j)=ad_UFx(i  ,j)+adfac
            ad_cff1=0.0_r8
          END DO
        END DO
!
!  Compute flux-components of the adjoint horizontal divergence of the
!  stress tensor (m5/s2) in XI- and ETA-directions.
!
        DO j=Jstr,Jend+1
          DO i=Istr,Iend+1
!^          tl_VFx(i,j)=on_p(i,j)*on_p(i,j)*tl_cff
!^          tl_UFe(i,j)=om_p(i,j)*om_p(i,j)*tl_cff
!^
            ad_cff=ad_cff+                                              &
     &             on_p(i,j)*on_p(i,j)*ad_VFx(i,j)+                     &
     &             om_p(i,j)*om_p(i,j)*ad_UFe(i,j)
            ad_VFx(i,j)=0.0_r8
            ad_UFe(i,j)=0.0_r8
!^          tl_cff=tl_cff*pmask(i,j)
!^
            ad_cff=ad_cff*pmask(i,j)
!^          tl_cff=visc2_p(i,j)*0.5_r8*                                 &
!^   &             (tl_Drhs_p(i,j)*                                     &
!^   &              (pmon_p(i,j)*                                       &
!^   &               ((pn(i  ,j-1)+pn(i  ,j))*vbar(i  ,j,krhs)-         &
!^   &                (pn(i-1,j-1)+pn(i-1,j))*vbar(i-1,j,krhs))+        &
!^   &               pnom_p(i,j)*                                       &
!^   &               ((pm(i-1,j  )+pm(i,j  ))*ubar(i,j  ,krhs)-         &
!^   &                (pm(i-1,j-1)+pm(i,j-1))*ubar(i,j-1,krhs)))+       &
!^   &              Drhs_p(i,j)*                                        &
!^   &              (pmon_p(i,j)*                                       &
!^   &               ((pn(i  ,j-1)+pn(i  ,j))*tl_vbar(i  ,j,krhs)-      &
!^   &                (pn(i-1,j-1)+pn(i-1,j))*tl_vbar(i-1,j,krhs))+     &
!^   &               pnom_p(i,j)*                                       &
!^   &               ((pm(i-1,j  )+pm(i,j  ))*tl_ubar(i,j  ,krhs)-      &
!^   &                (pm(i-1,j-1)+pm(i,j-1))*tl_ubar(i,j-1,krhs))))
!^
            adfac=visc2_p(i,j)*0.5_r8*ad_cff
            adfac1=adfac*Drhs_p(i,j)
            adfac2=adfac1*pmon_p(i,j)
            adfac3=adfac1*pnom_p(i,j)
            ad_Drhs_p(i,j)=ad_Drhs_p(i,j)+                              &
     &                     (pmon_p(i,j)*                                &
     &                      ((pn(i  ,j-1)+pn(i  ,j))*vbar(i  ,j,krhs)-  &
     &                       (pn(i-1,j-1)+pn(i-1,j))*vbar(i-1,j,krhs))+ &
     &                      pnom_p(i,j)*                                &
     &                      ((pm(i-1,j  )+pm(i,j  ))*ubar(i,j  ,krhs)-  &
     &                       (pm(i-1,j-1)+pm(i,j-1))*ubar(i,j-1,krhs)))*&
     &                     adfac
            ad_vbar(i-1,j,krhs)=ad_vbar(i-1,j,krhs)-                    &
     &                          (pn(i-1,j-1)+pn(i-1,j))*adfac2
            ad_vbar(i  ,j,krhs)=ad_vbar(i  ,j,krhs)+                    &
     &                          (pn(i  ,j-1)+pn(i  ,j))*adfac2
            ad_ubar(i,j-1,krhs)=ad_ubar(i,j-1,krhs)-                    &
     &                          (pm(i-1,j-1)+pm(i,j-1))*adfac3
            ad_ubar(i,j  ,krhs)=ad_ubar(i,j  ,krhs)+                    &
     &                          (pm(i-1,j  )+pm(i,j  ))*adfac3
            ad_cff=0.0_r8
          END DO
        END DO
        DO j=JstrV-1,Jend
          DO i=IstrU-1,Iend
!^          tl_VFe(i,j)=om_r(i,j)*om_r(i,j)*tl_cff
!^          tl_UFx(i,j)=on_r(i,j)*on_r(i,j)*tl_cff
!^
            ad_cff=ad_cff+                                              &
     &             om_r(i,j)*om_r(i,j)*ad_VFe(i,j)+                     &
     &             on_r(i,j)*on_r(i,j)*ad_UFx(i,j)
            ad_VFe(i,j)=0.0_r8
            ad_UFx(i,j)=0.0_r8
!^          tl_cff=visc2_r(i,j)*0.5_r8*                                 &
!^   &             (tl_Drhs(i,j)*                                       &
!^   &              (pmon_r(i,j)*                                       &
!^   &               ((pn(i  ,j)+pn(i+1,j))*ubar(i+1,j,krhs)-           &
!^   &                (pn(i-1,j)+pn(i  ,j))*ubar(i  ,j,krhs))-          &
!^   &               pnom_r(i,j)*                                       &
!^   &               ((pm(i,j  )+pm(i,j+1))*vbar(i,j+1,krhs)-           &
!^   &                (pm(i,j-1)+pm(i,j  ))*vbar(i,j  ,krhs)))+         &
!^   &              Drhs(i,j)*                                          &
!^   &              (pmon_r(i,j)*                                       &
!^   &               ((pn(i  ,j)+pn(i+1,j))*tl_ubar(i+1,j,krhs)-        &
!^   &                (pn(i-1,j)+pn(i  ,j))*tl_ubar(i  ,j,krhs))-       &
!^   &               pnom_r(i,j)*                                       &
!^   &               ((pm(i,j  )+pm(i,j+1))*tl_vbar(i,j+1,krhs)-        &
!^   &                (pm(i,j-1)+pm(i,j  ))*tl_vbar(i,j  ,krhs))))
!^
            adfac=visc2_r(i,j)*0.5_r8*ad_cff
            adfac1=adfac*Drhs(i,j)
            adfac2=adfac1*pmon_r(i,j)
            adfac3=adfac1*pnom_r(i,j)
            ad_Drhs(i,j)=ad_Drhs(i,j)+                                  &
     &                   (pmon_r(i,j)*                                  &
     &                    ((pn(i  ,j)+pn(i+1,j))*ubar(i+1,j,krhs)-      &
     &                     (pn(i-1,j)+pn(i  ,j))*ubar(i  ,j,krhs))-     &
     &                    pnom_r(i,j)*                                  &
     &                    ((pm(i,j  )+pm(i,j+1))*vbar(i,j+1,krhs)-      &
     &                     (pm(i,j-1)+pm(i,j  ))*vbar(i,j  ,krhs)))*    &
     &                   adfac
            ad_ubar(i  ,j,krhs)=ad_ubar(i  ,j,krhs)-                    &
     &                          (pn(i-1,j)+pn(i  ,j))*adfac2
            ad_ubar(i+1,j,krhs)=ad_ubar(i+1,j,krhs)+                    &
     &                          (pn(i  ,j)+pn(i+1,j))*adfac2
            ad_vbar(i,j  ,krhs)=ad_vbar(i,j  ,krhs)+                    &
     &                          (pm(i,j-1)+pm(i,j  ))*adfac3
            ad_vbar(i,j+1,krhs)=ad_vbar(i,j+1,krhs)-                    &
     &                          (pm(i,j  )+pm(i,j+1))*adfac3
            ad_cff=0.0_r8
          END DO
        END DO
!
!-----------------------------------------------------------------------
!  If horizontal mixing, compute adjoint total depth at PSI-points.
!-----------------------------------------------------------------------
!
        DO j=Jstr,Jend+1
          DO i=Istr,Iend+1
            Drhs_p(i,j)=0.25_r8*(Drhs(i,j  )+Drhs(i-1,j  )+             &
     &                           Drhs(i,j-1)+Drhs(i-1,j-1))
!^          tl_Drhs_p(i,j)=0.25_r8*(tl_Drhs(i,j  )+tl_Drhs(i-1,j  )+    &
!^   &                              tl_Drhs(i,j-1)+tl_Drhs(i-1,j-1))
!^
            adfac=0.25_r8*ad_Drhs_p(i,j)
            ad_Drhs(i-1,j  )=ad_Drhs(i-1,j  )+adfac
            ad_Drhs(i  ,j  )=ad_Drhs(i  ,j  )+adfac
            ad_Drhs(i-1,j-1)=ad_Drhs(i-1,j-1)+adfac
            ad_Drhs(i  ,j-1)=ad_Drhs(i  ,j-1)+adfac
            ad_Drhs_p(i,j)=0.0_r8
          END DO
        END DO
!
!-----------------------------------------------------------------------
!  Add in curvilinear transformation terms.
!-----------------------------------------------------------------------
!
      DO j=JstrV,Jend
        DO i=Istr,Iend
!^        tl_rhs_vbar(i,j)=tl_rhs_vbar(i,j)-tl_fac1
!^
          ad_fac1=ad_fac1-ad_rhs_vbar(i,j)
!^        tl_fac1=0.5_r8*(tl_VFe(i,j)+tl_VFe(i,j-1))
!^
          adfac=0.5_r8*ad_fac1
          ad_VFe(i,j-1)=ad_VFe(i,j-1)+adfac
          ad_VFe(i,j  )=ad_VFe(i,j  )+adfac
          ad_fac1=0.0_r8
        END DO
      END DO
      DO j=Jstr,Jend
        DO i=IstrU,Iend
!^        tl_rhs_ubar(i,j)=tl_rhs_ubar(i,j)+tl_fac1
!^
          ad_fac1=ad_fac1+ad_rhs_ubar(i,j)
!^        tl_fac1=0.5_r8*(tl_UFx(i,j)+tl_UFx(i-1,j))
!^
          adfac=0.5_r8*ad_fac1
          ad_UFx(i-1,j)=ad_UFx(i-1,j)+adfac
          ad_UFx(i  ,j)=ad_UFx(i  ,j)+adfac
          ad_fac1=0.0_r8
        END DO
      END DO
      DO j=JstrV-1,Jend
        DO i=IstrU-1,Iend
          cff1=0.5_r8*(vbar(i,j  ,krhs)+                                &
     &                 vbar(i,j+1,krhs))
          cff2=0.5_r8*(ubar(i  ,j,krhs)+                                &
     &                 ubar(i+1,j,krhs))
          cff3=cff1*dndx(i,j)
          cff4=cff2*dmde(i,j)
          cff=Drhs(i,j)*(cff3-cff4)
!^        tl_VFe(i,j)=tl_cff*cff2+cff*tl_cff2
!^        tl_UFx(i,j)=tl_cff*cff1+cff*tl_cff1
!^
          ad_cff=ad_cff+                                                &
     &           cff1*ad_UFx(i,j)+                                      &
     &           cff2*ad_VFe(i,j)
          ad_cff1=ad_cff1+cff*ad_UFx(i,j)
          ad_cff2=ad_cff2+cff*ad_VFe(i,j)
          ad_UFx(i,j)=0.0_r8
          ad_VFe(i,j)=0.0_r8
!^        tl_cff=tl_Drhs(i,j)*(cff3-cff4)+                              &
!^   &           Drhs(i,j)*(tl_cff3-tl_cff4)
!^
          adfac=Drhs(i,j)*ad_cff
          ad_cff4=ad_cff4-adfac
          ad_cff3=ad_cff3+adfac
          ad_Drhs(i,j)=ad_Drhs(i,j)+(cff3-cff4)*ad_cff
          ad_cff=0.0_r8
!^        tl_cff4=tl_cff2*dmde(i,j)
!^
          ad_cff2=ad_cff2+dmde(i,j)*ad_cff4
          ad_cff4=0.0_r8
!^        tl_cff3=tl_cff1*dndx(i,j)
!^
          ad_cff1=ad_cff1+dndx(i,j)*ad_cff3
          ad_cff3=0.0_r8
!^        tl_cff2=0.5_r8*(tl_ubar(i  ,j,krhs)+                          &
!^   &                    tl_ubar(i+1,j,krhs))
!^
          adfac=0.5_r8*ad_cff2
          ad_ubar(i  ,j,krhs)=ad_ubar(i  ,j,krhs)+adfac
          ad_ubar(i+1,j,krhs)=ad_ubar(i+1,j,krhs)+adfac
          ad_cff2=0.0_r8
!^        tl_cff1=0.5_r8*(tl_vbar(i,j  ,krhs)+                          &
!^   &                    tl_vbar(i,j+1,krhs))
!^
          adfac=0.5_r8*ad_cff1
          ad_vbar(i,j  ,krhs)=ad_vbar(i,j  ,krhs)+adfac
          ad_vbar(i,j+1,krhs)=ad_vbar(i,j+1,krhs)+adfac
          ad_cff1=0.0_r8
        END DO
      END DO
!
!-----------------------------------------------------------------------
!  Add in Coriolis term.
!-----------------------------------------------------------------------
!
      DO j=JstrV,Jend
        DO i=Istr,Iend
!^        tl_rhs_vbar(i,j)=tl_rhs_vbar(i,j)-tl_fac1
!^
          ad_fac1=ad_fac1-ad_rhs_vbar(i,j)
!^        tl_fac1=0.5_r8*(tl_VFe(i,j)+tl_VFe(i,j-1))
!^
          adfac=0.5_r8*ad_fac1
          ad_VFe(i,j-1)=ad_VFe(i,j-1)+adfac
          ad_VFe(i,j  )=ad_VFe(i,j  )+adfac
          ad_fac1=0.0_r8
        END DO
      END DO
      DO j=Jstr,Jend
        DO i=IstrU,Iend
!^        tl_rhs_ubar(i,j)=tl_rhs_ubar(i,j)+tl_fac1
!^
          ad_fac1=ad_fac1+ad_rhs_ubar(i,j)
!^        tl_fac1=0.5_r8*(tl_UFx(i,j)+tl_UFx(i-1,j))
!^
          adfac=0.5_r8*ad_fac1
          ad_UFx(i-1,j)=ad_UFx(i-1,j)+adfac
          ad_UFx(i  ,j)=ad_UFx(i  ,j)+adfac
          ad_fac1=0.0_r8
        END DO
      END DO
      DO j=JstrV-1,Jend
        DO i=IstrU-1,Iend
          cff=0.5_r8*Drhs(i,j)*fomn(i,j)
!^        tl_VFe(i,j)=tl_cff*(ubar(i  ,j,krhs)+                         &
!^   &                        ubar(i+1,j,krhs))+                        &
!^   &                cff*(tl_ubar(i  ,j,krhs)+                         &
!^   &                     tl_ubar(i+1,j,krhs))
!^
          adfac=cff*ad_VFe(i,j)
          ad_ubar(i  ,j,krhs)=ad_ubar(i  ,j,krhs)+adfac
          ad_ubar(i+1,j,krhs)=ad_ubar(i+1,j,krhs)+adfac
          ad_cff=ad_cff+(ubar(i  ,j,krhs)+                              &
     &                   ubar(i+1,j,krhs))*ad_VFe(i,j)
          ad_VFe(i,j)=0.0_r8
!^        tl_UFx(i,j)=tl_cff*(vbar(i,j  ,krhs)+                         &
!^   &                        vbar(i,j+1,krhs))+                        &
!^   &                cff*(tl_vbar(i,j  ,krhs)+                         &
!^   &                     tl_vbar(i,j+1,krhs))
!^
          adfac=cff*ad_UFx(i,j)
          ad_vbar(i,j  ,krhs)=ad_vbar(i,j  ,krhs)+adfac
          ad_vbar(i,j+1,krhs)=ad_vbar(i,j+1,krhs)+adfac
          ad_cff=ad_cff+(vbar(i,j  ,krhs)+                              &
     &                   vbar(i,j+1,krhs))*ad_UFx(i,j)
          ad_UFx(i,j)=0.0_r8
!^        tl_cff=0.5_r8*tl_Drhs(i,j)*fomn(i,j)
!^
          ad_Drhs(i,j)=ad_Drhs(i,j)+0.5_r8*fomn(i,j)*ad_cff
          ad_cff=0.0_r8
        END DO
      END DO
!
!-----------------------------------------------------------------------
!  Add in adjoint horizontal advection of momentum.
!-----------------------------------------------------------------------
!
        DO j=JstrV,Jend
          DO i=Istr,Iend
!^          tl_rhs_vbar(i,j)=tl_rhs_vbar(i,j)-tl_fac
!^
            ad_fac=ad_fac-ad_rhs_vbar(i,j)
!^          tl_fac=tl_cff1+tl_cff2
!^
            ad_cff1=ad_cff1+ad_fac
            ad_cff2=ad_cff2+ad_fac
            ad_fac=0.0_r8
!^          tl_cff2=tl_VFe(i,j)-tl_VFe(i,j-1)
!^
            ad_VFe(i,j-1)=ad_VFe(i,j-1)-ad_cff2
            ad_VFe(i,j  )=ad_VFe(i,j  )+ad_cff2
            ad_cff2=0.0_r8
!^          tl_cff1=tl_VFx(i+1,j)-tl_VFx(i,j)
!^
            ad_VFx(i  ,j)=ad_VFx(i  ,j)-ad_cff1
            ad_VFx(i+1,j)=ad_VFx(i+1,j)+ad_cff1
            ad_cff1=0.0_r8
          END DO
        END DO
        DO j=Jstr,Jend
          DO i=IstrU,Iend
!^          tl_rhs_ubar(i,j)=tl_rhs_ubar(i,j)-tl_fac
!^
            ad_fac=ad_fac-ad_rhs_ubar(i,j)
!^          tl_fac=tl_cff1+tl_cff2
!^
            ad_cff1=ad_cff1+ad_fac
            ad_cff2=ad_cff2+ad_fac
            ad_fac=0.0_r8
!^          tl_cff2=tl_UFe(i,j+1)-tl_UFe(i,j)
!^
            ad_UFe(i,j  )=ad_UFe(i,j  )-ad_cff2
            ad_UFe(i,j+1)=ad_UFe(i,j+1)+ad_cff2
            ad_cff2=0.0_r8
!^          tl_cff1=tl_UFx(i,j)-tl_UFx(i-1,j)
!^
            ad_UFx(i-1,j)=ad_UFx(i-1,j)-ad_cff1
            ad_UFx(i  ,j)=ad_UFx(i  ,j)+ad_cff1
            ad_cff1=0.0_r8
          END DO
        END DO
!
!  Fourth-order, centered differences advection.
!
        DO j=JstrVm1,Jendp1
          DO i=Istr,Iend
            grad (i,j)=vbar(i,j-1,krhs)-2.0_r8*vbar(i,j,krhs)+          &
     &                 vbar(i,j+1,krhs)
            Dgrad(i,j)=DVom(i,j-1)-2.0_r8*DVom(i,j)+DVom(i,j+1)
          END DO
        END DO
        IF (.not.(CompositeGrid(inorth,ng).or.NSperiodic(ng))) THEN
          IF (DOMAIN(ng)%Northern_Edge(tile)) THEN
            DO i=Istr,Iend
              grad (i,Jend+1)=grad (i,Jend)
              Dgrad(i,Jend+1)=Dgrad(i,Jend)
            END DO
          END IF
        END IF
        IF (.not.(CompositeGrid(isouth,ng).or.NSperiodic(ng))) THEN
          IF (DOMAIN(ng)%Southern_Edge(tile)) THEN
            DO i=Istr,Iend
              grad (i,Jstr)=grad (i,Jstr+1)
              Dgrad(i,Jstr)=Dgrad(i,Jstr+1)
            END DO
          END IF
        END IF
        cff=1.0_r8/6.0_r8
        DO j=JstrV-1,Jend
          DO i=Istr,Iend
!^          tl_VFe(i,j)=0.25_r8*                                        &
!^   &                  ((tl_vbar(i,j  ,krhs)+                          &
!^   &                    tl_vbar(i,j+1,krhs)-                          &
!^   &                    cff*(tl_grad (i,j)+tl_grad (i,j+1)))*         &
!^   &                   (DVom(i,j)+DVom(i,j+1)-                        &
!^   &                    cff*(Dgrad(i,j)+Dgrad(i,j+1)))+               &
!^   &                   (vbar(i,j  ,krhs)+                             &
!^   &                    vbar(i,j+1,krhs)-                             &
!^   &                    cff*(grad (i,j)+grad (i,j+1)))*               &
!^   &                   (tl_DVom(i,j)+tl_DVom(i,j+1)-                  &
!^   &                    cff*(tl_Dgrad(i,j)+tl_Dgrad(i,j+1))))
!^
            adfac=0.25_r8*ad_VFe(i,j)
            adfac1=adfac*(DVom(i,j)+DVom(i,j+1)-                        &
     &                    cff*(Dgrad(i,j)+Dgrad(i,j+1)))
            adfac2=adfac1*cff
            adfac3=adfac*(vbar(i,j  ,krhs)+                             &
     &                    vbar(i,j+1,krhs)-                             &
     &                    cff*(grad (i,j)+grad (i,j+1)))
            adfac4=adfac3*cff
            ad_vbar(i,j  ,krhs)=ad_vbar(i,j  ,krhs)+adfac1
            ad_vbar(i,j+1,krhs)=ad_vbar(i,j+1,krhs)+adfac1
            ad_grad (i,j  )=ad_grad (i,j  )-adfac2
            ad_grad (i,j+1)=ad_grad (i,j+1)-adfac2
            ad_DVom(i,j  )=ad_DVom(i,j  )+adfac3
            ad_DVom(i,j+1)=ad_DVom(i,j+1)+adfac3
            ad_Dgrad(i,j  )=ad_Dgrad(i,j  )-adfac4
            ad_Dgrad(i,j+1)=ad_Dgrad(i,j+1)-adfac4
            ad_VFe(i,j)=0.0_r8
          END DO
        END DO
        IF (.not.(CompositeGrid(inorth,ng).or.NSperiodic(ng))) THEN
          IF (DOMAIN(ng)%Northern_Edge(tile)) THEN
            DO i=Istr,Iend
!^            tl_Dgrad(i,Jend+1)=tl_Dgrad(i,Jend)
!^
              ad_Dgrad(i,Jend)=ad_Dgrad(i,Jend)+ad_Dgrad(i,Jend+1)
              ad_Dgrad(i,Jend+1)=0.0_r8
!^            tl_grad (i,Jend+1)=tl_grad (i,Jend)
!^
              ad_grad (i,Jend)=ad_grad (i,Jend)+ad_grad (i,Jend+1)
              ad_grad (i,Jend+1)=0.0_r8
            END DO
          END IF
        END IF
        IF (.not.(CompositeGrid(isouth,ng).or.NSperiodic(ng))) THEN
          IF (DOMAIN(ng)%Southern_Edge(tile)) THEN
            DO i=Istr,Iend
!^            tl_Dgrad(i,Jstr)=tl_Dgrad(i,Jstr+1)
!^
              ad_Dgrad(i,Jstr+1)=ad_Dgrad(i,Jstr+1)+ad_Dgrad(i,Jstr)
              ad_Dgrad(i,Jstr)=0.0_r8
!^            tl_grad (i,Jstr)=tl_grad (i,Jstr+1)
!^
              ad_grad (i,Jstr+1)=ad_grad (i,Jstr+1)+ad_grad (i,Jstr)
              ad_grad (i,Jstr)=0.0_r8
            END DO
          END IF
        END IF
        DO j=JstrVm1,Jendp1
          DO i=Istr,Iend
!^          tl_Dgrad(i,j)=tl_DVom(i,j-1)-2.0_r8*tl_DVom(i,j)+           &
!^   &                    tl_DVom(i,j+1)
!^
            ad_DVom(i,j-1)=ad_DVom(i,j-1)+ad_Dgrad(i,j)
            ad_DVom(i,j  )=ad_DVom(i,j  )-2.0_r8*ad_Dgrad(i,j)
            ad_DVom(i,j+1)=ad_DVom(i,j+1)+ad_Dgrad(i,j)
            ad_Dgrad(i,j)=0.0_r8
!^          tl_grad (i,j)=tl_vbar(i,j-1,krhs)-2.0_r8*tl_vbar(i,j,krhs)+ &
!^   &                    tl_vbar(i,j+1,krhs)
!^
            ad_vbar(i,j-1,krhs)=ad_vbar(i,j-1,krhs)+ad_grad(i,j)
            ad_vbar(i,j  ,krhs)=ad_vbar(i,j  ,krhs)-                    &
     &                          2.0_r8*ad_grad(i,j)
            ad_vbar(i,j+1,krhs)=ad_vbar(i,j+1,krhs)+ad_grad(i,j)
            ad_grad(i,j)=0.0_r8
          END DO
        END DO
        DO j=JstrV,Jend
          DO i=Istrm1,Iendp1
            grad(i,j)=vbar(i-1,j,krhs)-2.0_r8*vbar(i,j,krhs)+           &
     &                vbar(i+1,j,krhs)
          END DO
        END DO
        IF (.not.(CompositeGrid(iwest,ng).or.EWperiodic(ng))) THEN
          IF (DOMAIN(ng)%Western_Edge(tile)) THEN
            DO j=JstrV,Jend
              grad(Istr-1,j)=grad(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
              grad(Iend+1,j)=grad(Iend,j)
            END DO
          END IF
        END IF
        DO j=JstrV-1,Jend
          DO i=Istr,Iend+1
            Dgrad(i,j)=DUon(i,j-1)-2.0_r8*DUon(i,j)+DUon(i,j+1)
          END DO
        END DO
        cff=1.0_r8/6.0_r8
        DO j=JstrV,Jend
          DO i=Istr,Iend+1
!^          tl_VFx(i,j)=0.25_r8*                                        &
!^   &                  ((tl_vbar(i  ,j,krhs)+                          &
!^   &                    tl_vbar(i-1,j,krhs)-                          &
!^   &                    cff*(tl_grad (i,j)+tl_grad (i-1,j)))*         &
!^   &                   (DUon(i,j)+DUon(i,j-1)-                        &
!^   &                    cff*(Dgrad(i,j)+Dgrad(i,j-1)))+               &
!^   &                   (vbar(i  ,j,krhs)+                             &
!^   &                    vbar(i-1,j,krhs)-                             &
!^   &                    cff*(grad (i,j)+grad (i-1,j)))*               &
!^   &                   (tl_DUon(i,j)+tl_DUon(i,j-1)-                  &
!^   &                    cff*(tl_Dgrad(i,j)+tl_Dgrad(i,j-1))))
!^
            adfac=0.25_r8*ad_VFx(i,j)
            adfac1=adfac*(DUon(i,j)+DUon(i,j-1)-                        &
     &                    cff*(Dgrad(i,j)+Dgrad(i,j-1)))
            adfac2=adfac1*cff
            adfac3=adfac*(vbar(i  ,j,krhs)+                             &
     &                    vbar(i-1,j,krhs)-                             &
     &                    cff*(grad (i,j)+grad (i-1,j)))
            adfac4=adfac3*cff
            ad_vbar(i-1,j,krhs)=ad_vbar(i-1,j,krhs)+adfac1
            ad_vbar(i  ,j,krhs)=ad_vbar(i  ,j,krhs)+adfac1
            ad_grad (i-1,j)=ad_grad (i-1,j)-adfac2
            ad_grad (i  ,j)=ad_grad (i  ,j)-adfac2
            ad_DUon(i,j-1)=ad_DUon(i,j-1)+adfac3
            ad_DUon(i,j  )=ad_DUon(i,j  )+adfac3
            ad_Dgrad(i,j-1)=ad_Dgrad(i,j-1)-adfac4
            ad_Dgrad(i,j  )=ad_Dgrad(i,j  )-adfac4
            ad_VFx(i,j)=0.0_r8
          END DO
        END DO
        DO j=JstrV-1,Jend
          DO i=Istr,Iend+1
!^          tl_Dgrad(i,j)=tl_DUon(i,j-1)-2.0_r8*tl_DUon(i,j)+           &
!^   &                    tl_DUon(i,j+1)
!^
            ad_DUon(i,j-1)=ad_DUon(i,j-1)+ad_Dgrad(i,j)
            ad_DUon(i,j  )=ad_DUon(i,j  )-2.0_r8*ad_Dgrad(i,j)
            ad_DUon(i,j+1)=ad_DUon(i,j+1)+ad_Dgrad(i,j)
            ad_Dgrad(i,j)=0.0_r8
          END DO
        END DO
        IF (.not.(CompositeGrid(ieast,ng).or.EWperiodic(ng))) THEN
          IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN
            DO j=JstrV,Jend
!^            tl_grad(Iend+1,j)=tl_grad(Iend,j)
!^
              ad_grad(Iend,j)=ad_grad(Iend,j)+ad_grad(Iend+1,j)
              ad_grad(Iend+1,j)=0.0_r8
            END DO
          END IF
        END IF
        IF (.not.(CompositeGrid(iwest,ng).or.EWperiodic(ng))) THEN
          IF (DOMAIN(ng)%Western_Edge(tile)) THEN
            DO j=JstrV,Jend
!^            tl_grad(Istr-1,j)=tl_grad(Istr,j)
!^
              ad_grad(Istr,j)=ad_grad(Istr,j)+ad_grad(Istr-1,j)
              ad_grad(Istr-1,j)=0.0_r8
            END DO
           END IF
        END IF
        DO j=JstrV,Jend
          DO i=Istrm1,Iendp1
!^          tl_grad(i,j)=tl_vbar(i-1,j,krhs)-2.0_r8*tl_vbar(i,j,krhs)+  &
!^   &                   tl_vbar(i+1,j,krhs)
!^
            ad_vbar(i-1,j,krhs)=ad_vbar(i-1,j,krhs)+ad_grad(i,j)
            ad_vbar(i  ,j,krhs)=ad_vbar(i  ,j,krhs)-                    &
     &                          2.0_r8*ad_grad(i,j)
            ad_vbar(i+1,j,krhs)=ad_vbar(i+1,j,krhs)+ad_grad(i,j)
            ad_grad(i,j)=0.0_r8
          END DO
        END DO
        DO j=Jstrm1,Jendp1
          DO i=IstrU,Iend
            grad(i,j)=ubar(i,j-1,krhs)-2.0_r8*ubar(i,j,krhs)+           &
     &                ubar(i,j+1,krhs)
          END DO
        END DO
        IF (.not.(CompositeGrid(isouth,ng).or.NSperiodic(ng))) THEN
          IF (DOMAIN(ng)%Southern_Edge(tile)) THEN
            DO i=IstrU,Iend
              grad(i,Jstr-1)=grad(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
              grad(i,Jend+1)=grad(i,Jend)
            END DO
          END IF
        END IF
        DO j=Jstr,Jend+1
          DO i=IstrU-1,Iend
            Dgrad(i,j)=DVom(i-1,j)-2.0_r8*DVom(i,j)+DVom(i+1,j)
          END DO
        END DO
        cff=1.0_r8/6.0_r8
        DO j=Jstr,Jend+1
          DO i=IstrU,Iend
!^          tl_UFe(i,j)=0.25_r8*                                        &
!^   &                  ((tl_ubar(i,j  ,krhs)+                          &
!^   &                    tl_ubar(i,j-1,krhs)-                          &
!^   &                    cff*(tl_grad (i,j)+tl_grad (i,j-1)))*         &
!^   &                   (DVom(i,j)+DVom(i-1,j)-                        &
!^   &                    cff*(Dgrad(i,j)+Dgrad(i-1,j)))+               &
!^   &                   (ubar(i,j  ,krhs)+                             &
!^   &                    ubar(i,j-1,krhs)-                             &
!^   &                    cff*(grad (i,j)+grad (i,j-1)))*               &
!^   &                   (tl_DVom(i,j)+tl_DVom(i-1,j)-                  &
!^   &                    cff*(tl_Dgrad(i,j)+tl_Dgrad(i-1,j))))
!^
            adfac=0.25_r8*ad_UFe(i,j)
            adfac1=adfac*(DVom(i,j)+DVom(i-1,j)-                        &
     &                    cff*(Dgrad(i,j)+Dgrad(i-1,j)))
            adfac2=adfac1*cff
            adfac3=adfac*(ubar(i,j  ,krhs)+                             &
     &                    ubar(i,j-1,krhs)-                             &
     &                    cff*(grad (i,j)+grad (i,j-1)))
            adfac4=adfac3*cff
            ad_ubar(i,j-1,krhs)=ad_ubar(i,j-1,krhs)+adfac1
            ad_ubar(i,j  ,krhs)=ad_ubar(i,j  ,krhs)+adfac1
            ad_grad (i,j-1)=ad_grad (i,j-1)-adfac2
            ad_grad (i,j  )=ad_grad (i,j  )-adfac2
            ad_DVom(i-1,j)=ad_DVom(i-1,j)+adfac3
            ad_DVom(i  ,j)=ad_DVom(i  ,j)+adfac3
            ad_Dgrad(i-1,j)=ad_Dgrad(i-1,j)-adfac4
            ad_Dgrad(i  ,j)=ad_Dgrad(i  ,j)-adfac4
            ad_UFe(i,j)=0.0_r8
          END DO
        END DO
        DO j=Jstr,Jend+1
          DO i=IstrU-1,Iend
!^          tl_Dgrad(i,j)=tl_DVom(i-1,j)-2.0_r8*tl_DVom(i,j)+           &
!^   &                    tl_DVom(i+1,j)
!^
            ad_DVom(i-1,j)=ad_DVom(i-1,j)+ad_Dgrad(i,j)
            ad_DVom(i  ,j)=ad_DVom(i  ,j)-2.0_r8*ad_Dgrad(i,j)
            ad_DVom(i+1,j)=ad_DVom(i+1,j)+ad_Dgrad(i,j)
            ad_Dgrad(i,j)=0.0_r8
          END DO
        END DO
        IF (.not.(CompositeGrid(inorth,ng).or.NSperiodic(ng))) THEN
          IF (DOMAIN(ng)%Northern_Edge(tile)) THEN
            DO i=IstrU,Iend
!^            tl_grad(i,Jend+1)=tl_grad(i,Jend)
!^
              ad_grad(i,Jend)=ad_grad(i,Jend)+ad_grad(i,Jend+1)
              ad_grad(i,Jend+1)=0.0_r8
            END DO
          END IF
        END IF
        IF (.not.(CompositeGrid(isouth,ng).or.NSperiodic(ng))) THEN
          IF (DOMAIN(ng)%Southern_Edge(tile)) THEN
            DO i=IstrU,Iend
!^            tl_grad(i,Jstr-1)=tl_grad(i,Jstr)
!^
              ad_grad(i,Jstr)=ad_grad(i,Jstr)+ad_grad(i,Jstr-1)
              ad_grad(i,Jstr-1)=0.0_r8
            END DO
          END IF
        END IF
        DO j=Jstrm1,Jendp1
          DO i=IstrU,Iend
!^          tl_grad(i,j)=tl_ubar(i,j-1,krhs)-2.0_r8*tl_ubar(i,j,krhs)+  &
!^   &                   tl_ubar(i,j+1,krhs)
!^
            ad_ubar(i,j-1,krhs)=ad_ubar(i,j-1,krhs)+ad_grad(i,j)
            ad_ubar(i,j  ,krhs)=ad_ubar(i,j  ,krhs)-                    &
     &                          2.0_r8*ad_grad(i,j)
            ad_ubar(i,j+1,krhs)=ad_ubar(i,j+1,krhs)+ad_grad(i,j)
            ad_grad(i,j)=0.0_r8
          END DO
        END DO
        DO j=Jstr,Jend
          DO i=IstrUm1,Iendp1
            grad (i,j)=ubar(i-1,j,krhs)-2.0_r8*ubar(i,j,krhs)+          &
     &                 ubar(i+1,j,krhs)
            Dgrad(i,j)=DUon(i-1,j)-2.0_r8*DUon(i,j)+DUon(i+1,j)
          END DO
        END DO
        IF (.not.(CompositeGrid(iwest,ng).or.EWperiodic(ng))) THEN
          IF (DOMAIN(ng)%Western_Edge(tile)) THEN
            DO j=Jstr,Jend
              grad (Istr,j)=grad (Istr+1,j)
              Dgrad(Istr,j)=Dgrad(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
              grad (Iend+1,j)=grad (Iend,j)
              Dgrad(Iend+1,j)=Dgrad(Iend,j)
            END DO
          END IF
        END IF
        cff=1.0_r8/6.0_r8
        DO j=Jstr,Jend
          DO i=IstrU-1,Iend
!^          tl_UFx(i,j)=0.25_r8*                                        &
!^   &                  ((ubar(i  ,j,krhs)+                             &
!^   &                    ubar(i+1,j,krhs)-                             &
!^   &                    cff*(grad (i,j)+grad (i+1,j)))*               &
!^   &                   (tl_DUon(i,j)+tl_DUon(i+1,j)-                  &
!^   &                    cff*(tl_Dgrad(i,j)+tl_Dgrad(i+1,j)))+         &
!^   &                   (tl_ubar(i  ,j,krhs)+                          &
!^   &                    tl_ubar(i+1,j,krhs)-                          &
!^   &                    cff*(tl_grad (i,j)+tl_grad (i+1,j)))*         &
!^   &                   (DUon(i,j)+DUon(i+1,j)-                        &
!^   &                    cff*(Dgrad(i,j)+Dgrad(i+1,j))))
!^
            adfac=0.25_r8*ad_UFx(i,j)
            adfac1=adfac*(DUon(i,j)+DUon(i+1,j)-                        &
     &                    cff*(Dgrad(i,j)+Dgrad(i+1,j)))
            adfac2=adfac1*cff
            adfac3=adfac*(ubar(i  ,j,krhs)+                             &
     &                    ubar(i+1,j,krhs)-                             &
     &                    cff*(grad (i,j)+grad (i+1,j)))
            adfac4=adfac3*cff
            ad_ubar(i  ,j,krhs)=ad_ubar(i  ,j,krhs)+adfac1
            ad_ubar(i+1,j,krhs)=ad_ubar(i+1,j,krhs)+adfac1
            ad_grad (i  ,j)=ad_grad (i  ,j)-adfac2
            ad_grad (i+1,j)=ad_grad (i+1,j)-adfac2
            ad_DUon(i  ,j)=ad_DUon(i  ,j)+adfac3
            ad_DUon(i+1,j)=ad_DUon(i+1,j)+adfac3
            ad_Dgrad(i  ,j)=ad_Dgrad(i  ,j)-adfac4
            ad_Dgrad(i+1,j)=ad_Dgrad(i+1,j)-adfac4
            ad_UFx(i,j)=0.0_r8
          END DO
        END DO
        IF (.not.(CompositeGrid(ieast,ng).or.EWperiodic(ng))) THEN
          IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN
            DO j=Jstr,Jend
!^            tl_Dgrad(Iend+1,j)=tl_Dgrad(Iend,j)
!^
              ad_Dgrad(Iend,j)=ad_Dgrad(Iend,j)+ad_Dgrad(Iend+1,j)
              ad_Dgrad(Iend+1,j)=0.0_r8
!^            tl_grad (Iend+1,j)=tl_grad (Iend,j)
!^
              ad_grad (Iend,j)=ad_grad (Iend,j)+ad_grad (Iend+1,j)
              ad_grad (Iend+1,j)=0.0_r8
            END DO
          END IF
        END IF
        IF (.not.(CompositeGrid(iwest,ng).or.EWperiodic(ng))) THEN
          IF (DOMAIN(ng)%Western_Edge(tile)) THEN
            DO j=Jstr,Jend
!^            tl_Dgrad(Istr,j)=tl_Dgrad(Istr+1,j)
!^
              ad_Dgrad(Istr+1,j)=ad_Dgrad(Istr+1,j)+ad_Dgrad(Istr,j)
              ad_Dgrad(Istr,j)=0.0_r8
!^            tl_grad (Istr,j)=tl_grad (Istr+1,j)
!^
              ad_grad (Istr+1,j)=ad_grad (Istr+1,j)+ad_grad (Istr,j)
              ad_grad (Istr,j)=0.0_r8
            END DO
          END IF
        END IF
        DO j=Jstr,Jend
          DO i=IstrUm1,Iendp1
!^          tl_Dgrad(i,j)=tl_DUon(i-1,j)-2.0_r8*tl_DUon(i,j)+           &
!^   &                    tl_DUon(i+1,j)
!^
            ad_DUon(i-1,j)=ad_DUon(i-1,j)+ad_Dgrad(i,j)
            ad_DUon(i  ,j)=ad_DUon(i  ,j)-2.0_r8*ad_Dgrad(i,j)
            ad_DUon(i+1,j)=ad_DUon(i+1,j)+ad_Dgrad(i,j)
            ad_Dgrad(i,j)=0.0_r8
!^          tl_grad (i,j)=tl_ubar(i-1,j,krhs)-2.0_r8*tl_ubar(i,j,krhs)+ &
!^   &                    tl_ubar(i+1,j,krhs)
!^
            ad_ubar(i-1,j,krhs)=ad_ubar(i-1,j,krhs)+ad_grad (i,j)
            ad_ubar(i  ,j,krhs)=ad_ubar(i  ,j,krhs)-                    &
     &                          2.0_r8*ad_grad (i,j)
            ad_ubar(i+1,j,krhs)=ad_ubar(i+1,j,krhs)+ad_grad (i,j)
            ad_grad(i,j)=0.0_r8
          END DO
        END DO
!
!-----------------------------------------------------------------------
!  Compute adjoint pressure gradient terms.
!-----------------------------------------------------------------------
!
!  Compute BASIC STATE fields associated with pressure gradient and
!  time-stepping of adjoint free-surface.
!
        fac=1000.0_r8/rho0
        IF (iif(ng).eq.1) THEN
          cff1=dtfast(ng)
          DO j=JstrV-1,Jend
            DO i=IstrU-1,Iend
!^            rhs_zeta(i,j)=(DUon(i,j)-DUon(i+1,j))+                    &
!^   &                      (DVom(i,j)-DVom(i,j+1))
!^            zeta_new(i,j)=zeta(i,j,kstp)+                             &
!^   &                      pm(i,j)*pn(i,j)*cff1*rhs_zeta(i,j)
!^
!^                                                use background instead
              zeta_new(i,j)=zeta(i,j,knew)
              zeta_new(i,j)=zeta_new(i,j)*rmask(i,j)
              zwrk(i,j)=0.5_r8*(zeta(i,j,kstp)+zeta_new(i,j))
              gzeta(i,j)=zwrk(i,j)
              gzeta2(i,j)=zwrk(i,j)*zwrk(i,j)
            END DO
          END DO
        ELSE IF (PREDICTOR_2D_STEP(ng)) THEN
          cff1=2.0_r8*dtfast(ng)
          cff4=4.0_r8/25.0_r8
          cff5=1.0_r8-2.0_r8*cff4
          DO j=JstrV-1,Jend
            DO i=IstrU-1,Iend
!^            rhs_zeta(i,j)=(DUon(i,j)-DUon(i+1,j))+                    &
!^   &                      (DVom(i,j)-DVom(i,j+1))
!^            zeta_new(i,j)=zeta(i,j,kstp)+                             &
!^   &                      pm(i,j)*pn(i,j)*cff1*rhs_zeta(i,j)
!^
!^                                                use background instead
              zeta_new(i,j)=zeta(i,j,knew)
              zeta_new(i,j)=zeta_new(i,j)*rmask(i,j)
              zwrk(i,j)=cff5*zeta(i,j,krhs)+                            &
     &                  cff4*(zeta(i,j,kstp)+zeta_new(i,j))
              gzeta(i,j)=zwrk(i,j)
              gzeta2(i,j)=zwrk(i,j)*zwrk(i,j)
            END DO
          END DO
        ELSE IF (CORRECTOR_2D_STEP) THEN
          cff1=dtfast(ng)*5.0_r8/12.0_r8
          cff2=dtfast(ng)*8.0_r8/12.0_r8
          cff3=dtfast(ng)*1.0_r8/12.0_r8
          cff4=2.0_r8/5.0_r8
          cff5=1.0_r8-cff4
          DO j=JstrV-1,Jend
            DO i=IstrU-1,Iend
!^            cff=cff1*((DUon(i,j)-DUon(i+1,j))+                        &
!^   &                  (DVom(i,j)-DVom(i,j+1)))
!^            zeta_new(i,j)=zeta(i,j,kstp)+                             &
!^   &                      pm(i,j)*pn(i,j)*(cff+                       &
!^   &                                       cff2*rzeta(i,j,kstp)-      &
!^   &                                       cff3*rzeta(i,j,ptsk))
!^
!^                                                use background instead
              zeta_new(i,j)=zeta(i,j,knew)
              zeta_new(i,j)=zeta_new(i,j)*rmask(i,j)
              zwrk(i,j)=cff5*zeta_new(i,j)+cff4*zeta(i,j,krhs)
              gzeta(i,j)=zwrk(i,j)
              gzeta2(i,j)=zwrk(i,j)*zwrk(i,j)
            END DO
          END DO
        END IF
!
!  Compute adjoint pressure gradient.
!
        cff1=0.5_r8*g
        cff2=1.0_r8/3.0_r8
        DO j=Jstr,Jend
          IF (j.ge.JstrV) THEN
            DO i=Istr,Iend
!^            tl_rhs_vbar(i,j)=cff1*om_v(i,j)*                          &
!^   &                         ((tl_h(i,j-1)+                           &
!^   &                           tl_h(i,j  ))*                          &
!^   &                          (gzeta(i,j-1)-                          &
!^   &                           gzeta(i,j  ))+                         &
!^   &                          (h(i,j-1)+                              &
!^   &                           h(i,j  ))*                             &
!^   &                          (tl_gzeta(i,j-1)-                       &
!^   &                           tl_gzeta(i,j  ))+                      &
!^   &                          (tl_gzeta2(i,j-1)-                      &
!^   &                           tl_gzeta2(i,j  )
!^
              adfac=cff1*om_v(i,j)*ad_rhs_vbar(i,j)
              adfac1=adfac*(gzeta(i,j-1)-gzeta(i,j  ))
              adfac2=adfac*(h(i,j-1)+h(i,j  ))
              ad_h(i,j-1)=ad_h(i,j-1)+adfac1
              ad_h(i,j  )=ad_h(i,j  )+adfac1
              ad_gzeta(i,j-1)=ad_gzeta(i,j-1)+adfac2
              ad_gzeta(i,j  )=ad_gzeta(i,j  )-adfac2
              ad_gzeta2(i,j-1)=ad_gzeta2(i,j-1)+adfac
              ad_gzeta2(i,j  )=ad_gzeta2(i,j  )-adfac
              ad_rhs_vbar(i,j)=0.0_r8
            END DO
          END IF
          DO i=IstrU,Iend
!^          tl_rhs_ubar(i,j)=cff1*on_u(i,j)*                            &
!^   &                       ((tl_h(i-1,j)+                             &
!^   &                         tl_h(i ,j))*                             &
!^   &                        (gzeta(i-1,j)-                            &
!^   &                         gzeta(i  ,j))+                           &
!^   &                        (h(i-1,j)+                                &
!^   &                         h(i  ,j))*                               &
!^   &                        (tl_gzeta(i-1,j)-                         &
!^   &                         tl_gzeta(i  ,j))+                        &
!^   &                        (tl_gzeta2(i-1,j)-                        &
!^   &                         tl_gzeta2(i  ,j)))
!^
            adfac=cff1*on_u(i,j)*ad_rhs_ubar(i,j)
            adfac1=adfac*(gzeta(i-1,j)-gzeta(i  ,j))
            adfac2=adfac*(h(i-1,j)+h(i  ,j))
            ad_h(i-1,j)=ad_h(i-1,j)+adfac1
            ad_h(i  ,j)=ad_h(i  ,j)+adfac1
            ad_gzeta(i-1,j)=ad_gzeta(i-1,j)+adfac2
            ad_gzeta(i  ,j)=ad_gzeta(i  ,j)-adfac2
            ad_gzeta2(i-1,j)=ad_gzeta2(i-1,j)+adfac
            ad_gzeta2(i  ,j)=ad_gzeta2(i  ,j)-adfac
            ad_rhs_ubar(i,j)=0.0_r8
          END DO
        END DO
!
!  Set adjoint free-surface lateral boundary conditions.
!
!^      CALL mp_exchange2d (ng, tile, iTLM, 1,                          &
!^   &                      LBi, UBi, LBj, UBj,                         &
!^   &                      NghostPoints,                               &
!^   &                      EWperiodic(ng), NSperiodic(ng),             &
!^   &                      tl_zeta(:,:,knew))
!^
        CALL ad_mp_exchange2d (ng, tile, iADM, 1,                       &
     &                         LBi, UBi, LBj, UBj,                      &
     &                         NghostPoints,                            &
     &                         EWperiodic(ng), NSperiodic(ng),          &
     &                         ad_zeta(:,:,knew))
        IF (EWperiodic(ng).or.NSperiodic(ng)) THEN
!^        CALL exchange_r2d_tile (ng, tile,                             &
!^   &                            LBi, UBi, LBj, UBj,                   &
!^   &                            tl_zeta(:,:,knew))
!^
          CALL ad_exchange_r2d_tile (ng, tile,                          &
     &                               LBi, UBi, LBj, UBj,                &
     &                               ad_zeta(:,:,knew))
        END IF
!^      CALL tl_zetabc_tile (ng, tile,                                  &
!^   &                       LBi, UBi, LBj, UBj,                        &
!^   &                       IminS, ImaxS, JminS, JmaxS,                &
!^   &                       krhs, kstp, knew,                          &
!^   &                       zeta, tl_zeta)
!^
        CALL ad_zetabc_tile (ng, tile,                                  &
     &                       LBi, UBi, LBj, UBj,                        &
     &                       IminS, ImaxS, JminS, JmaxS,                &
     &                       krhs, kstp, knew,                          &
     &                       zeta, ad_zeta)
!
!  Apply adjoint mass point sources (volume vertical influx), if any.
!
!    Dsrc(is) = 2,  flow across grid cell w-face (positive or negative)
!
        IF (LwSrc(ng)) THEN
          DO is=1,Nsrc(ng)
            IF (INT(SOURCES(ng)%Dsrc(is)).eq.2) THEN
              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
!^              tl_zeta(i,j,knew)=tl_zeta(i,j,knew)+0.0_r8
!^
              END IF
            END IF
          END DO
        END IF
!
!  If adjoint predictor step, load right-side-term into shared array.
!
        IF (PREDICTOR_2D_STEP(ng)) THEN
!^        CALL mp_exchange2d (ng, tile, iTLM, 1,                        &
!^   &                        LBi, UBi, LBj, UBj,                       &
!^   &                        NghostPoints,                             &
!^   &                        EWperiodic(ng), NSperiodic(ng),           &
!^   &                        tl_rzeta(:,:,krhs))
!^
          CALL ad_mp_exchange2d (ng, tile, iADM, 1,                     &
     &                           LBi, UBi, LBj, UBj,                    &
     &                           NghostPoints,                          &
     &                           EWperiodic(ng), NSperiodic(ng),        &
     &                           ad_rzeta(:,:,krhs))
          IF (EWperiodic(ng).or.NSperiodic(ng)) THEN
!^          CALL exchange_r2d_tile (ng, tile,                           &
!^   &                              LBi, UBi, LBj, UBj,                 &
!^   &                              tl_rzeta(:,:,krhs))
!^
            CALL ad_exchange_r2d_tile (ng, tile,                        &
     &                                 LBi, UBi, LBj, UBj,              &
     &                                 ad_rzeta(:,:,krhs))
          END IF
          DO j=Jstr,Jend
            DO i=Istr,Iend
!^            tl_rzeta(i,j,krhs)=tl_rhs_zeta(i,j)
!^
              ad_rhs_zeta(i,j)=ad_rhs_zeta(i,j)+ad_rzeta(i,j,krhs)
              ad_rzeta(i,j,krhs)=0.0
            END DO
          END DO
        END IF
!
!  Load new adjoint free-surface values into shared array at both
!  predictor and corrector steps.
!
        DO j=Jstr,Jend
          DO i=Istr,Iend
!^          tl_zeta(i,j,knew)=tl_zeta_new(i,j)
!^
            ad_zeta_new(i,j)=ad_zeta_new(i,j)+ad_zeta(i,j,knew)
            ad_zeta(i,j,knew)=0.0_r8
          END DO
        END DO
!
!=======================================================================
!  Time step adjoint free-surface equation.
!=======================================================================
!
!  During the first time-step, the predictor step is Forward-Euler
!  and the corrector step is Backward-Euler. Otherwise, the predictor
!  step is Leap-frog and the corrector step is Adams-Moulton.
!
        IF (iif(ng).eq.1) THEN
          cff1=dtfast(ng)
          DO j=JstrV-1,Jend
            DO i=IstrU-1,Iend
!^            tl_gzeta2(i,j)=2.0_r8*tl_zwrk(i,j)*zwrk(i,j)
!^            tl_gzeta(i,j)=tl_zwrk(i,j)
!^
              ad_zwrk(i,j)=ad_zwrk(i,j)+                                &
     &                     2.0_r8*zwrk(i,j)*ad_gzeta2(i,j)+             &
     &                     ad_gzeta(i,j)
              ad_gzeta2(i,j)=0.0_r8
              ad_gzeta(i,j)=0.0_r8
!^            tl_zwrk(i,j)=0.5_r8*(tl_zeta(i,j,kstp)+tl_zeta_new(i,j))
!^
              adfac=0.5_r8*ad_zwrk(i,j)
              ad_zeta(i,j,kstp)=ad_zeta(i,j,kstp)+adfac
              ad_zeta_new(i,j)=ad_zeta_new(i,j)+adfac
              ad_zwrk(i,j)=0.0_r8
!^            tl_Dnew(i,j)=tl_zeta_new(i,j)+tl_h(i,j)
!^
              ad_zeta_new(i,j)=ad_zeta_new(i,j)+ad_Dnew(i,j)
              ad_h(i,j)=ad_h(i,j)+ad_Dnew(i,j)
              ad_Dnew(i,j)=0.0_r8
!^            tl_zeta_new(i,j)=tl_zeta_new(i,j)*rmask(i,j)
!^
              ad_zeta_new(i,j)=ad_zeta_new(i,j)*rmask(i,j)
!^            tl_zeta_new(i,j)=tl_zeta(i,j,kstp)+                       &
!^   &                         pm(i,j)*pn(i,j)*cff1*tl_rhs_zeta(i,j)
!^
              ad_zeta(i,j,kstp)=ad_zeta(i,j,kstp)+ad_zeta_new(i,j)
              ad_rhs_zeta(i,j)=ad_rhs_zeta(i,j)+                        &
     &                         pm(i,j)*pn(i,j)*cff1*ad_zeta_new(i,j)
              ad_zeta_new(i,j)=0.0_r8
!^            tl_rhs_zeta(i,j)=(tl_DUon(i,j)-tl_DUon(i+1,j))+           &
!^   &                         (tl_DVom(i,j)-tl_DVom(i,j+1))
!^
              ad_DUon(i  ,j  )=ad_DUon(i  ,j  )+ad_rhs_zeta(i,j)
              ad_DUon(i+1,j  )=ad_DUon(i+1,j  )-ad_rhs_zeta(i,j)
              ad_DVom(i  ,j  )=ad_DVom(i  ,j  )+ad_rhs_zeta(i,j)
              ad_DVom(i  ,j+1)=ad_DVom(i  ,j+1)-ad_rhs_zeta(i,j)
              ad_rhs_zeta(i,j)=0.0_r8
            END DO
          END DO
        ELSE IF (PREDICTOR_2D_STEP(ng)) THEN
          cff1=2.0_r8*dtfast(ng)
          cff4=4.0_r8/25.0_r8
          cff5=1.0_r8-2.0_r8*cff4
          DO j=JstrV-1,Jend
            DO i=IstrU-1,Iend
!^            tl_gzeta2(i,j)=2.0_r8*tl_zwrk(i,j)*zwrk(i,j)
!^            tl_gzeta(i,j)=tl_zwrk(i,j)
!^
              ad_zwrk(i,j)=ad_zwrk(i,j)+                                &
     &                     2.0_r8*zwrk(i,j)*ad_gzeta2(i,j)+             &
     &                     ad_gzeta(i,j)
              ad_gzeta2(i,j)=0.0_r8
              ad_gzeta(i,j)=0.0_r8
!^            tl_zwrk(i,j)=cff5*tl_zeta(i,j,krhs)+                      &
!^   &                     cff4*(tl_zeta(i,j,kstp)+tl_zeta_new(i,j))
!^
              adfac=cff4*ad_zwrk(i,j)
              ad_zeta(i,j,krhs)=ad_zeta(i,j,krhs)+cff5*ad_zwrk(i,j)
              ad_zeta(i,j,kstp)=ad_zeta(i,j,kstp)+adfac
              ad_zeta_new(i,j)=ad_zeta_new(i,j)+adfac
              ad_zwrk(i,j)=0.0_r8
!^            tl_Dnew(i,j)=tl_zeta_new(i,j)+tl_h(i,j)
!^
              ad_zeta_new(i,j)=ad_zeta_new(i,j)+ad_Dnew(i,j)
              ad_h(i,j)=ad_h(i,j)+ad_Dnew(i,j)
              ad_Dnew(i,j)=0.0_r8
!^            tl_zeta_new(i,j)=tl_zeta_new(i,j)*rmask(i,j)
!^
              ad_zeta_new(i,j)=ad_zeta_new(i,j)*rmask(i,j)
!^            tl_zeta_new(i,j)=tl_zeta(i,j,kstp)+                       &
!^   &                         pm(i,j)*pn(i,j)*cff1*tl_rhs_zeta(i,j)
!^
              ad_zeta(i,j,kstp)=ad_zeta(i,j,kstp)+ad_zeta_new(i,j)
              ad_rhs_zeta(i,j)=ad_rhs_zeta(i,j)+                        &
     &                         pm(i,j)*pn(i,j)*cff1*ad_zeta_new(i,j)
              ad_zeta_new(i,j)=0.0_r8
!^            tl_rhs_zeta(i,j)=(tl_DUon(i,j)-tl_DUon(i+1,j))+           &
!^   &                         (tl_DVom(i,j)-tl_DVom(i,j+1))
!^
              ad_DUon(i  ,j  )=ad_DUon(i  ,j  )+ad_rhs_zeta(i,j)
              ad_DUon(i+1,j  )=ad_DUon(i+1,j  )-ad_rhs_zeta(i,j)
              ad_DVom(i  ,j  )=ad_DVom(i  ,j  )+ad_rhs_zeta(i,j)
              ad_DVom(i  ,j+1)=ad_DVom(i  ,j+1)-ad_rhs_zeta(i,j)
              ad_rhs_zeta(i,j)=0.0_r8
            END DO
          END DO
        ELSE IF (CORRECTOR_2D_STEP) THEN
          cff1=dtfast(ng)*5.0_r8/12.0_r8
          cff2=dtfast(ng)*8.0_r8/12.0_r8
          cff3=dtfast(ng)*1.0_r8/12.0_r8
          cff4=2.0_r8/5.0_r8
          cff5=1.0_r8-cff4
          DO j=JstrV-1,Jend
            DO i=IstrU-1,Iend
!^            tl_gzeta(i,j)=tl_zwrk(i,j)
!^            tl_gzeta2(i,j)=2.0_r8*tl_zwrk(i,j)*zwrk(i,j)
!^
              ad_zwrk(i,j)=ad_zwrk(i,j)+                                &
     &                     2.0_r8*zwrk(i,j)*ad_gzeta2(i,j)+             &
     &                     ad_gzeta(i,j)
              ad_gzeta2(i,j)=0.0_r8
              ad_gzeta(i,j)=0.0_r8
!^            tl_zwrk(i,j)=cff5*tl_zeta_new(i,j)+cff4*tl_zeta(i,j,krhs)
!^
              ad_zeta_new(i,j)=ad_zeta_new(i,j)+cff5*ad_zwrk(i,j)
              ad_zeta(i,j,krhs)=ad_zeta(i,j,krhs)+cff4*ad_zwrk(i,j)
              ad_zwrk(i,j)=0.0_r8
!^            tl_Dnew(i,j)=tl_zeta_new(i,j)+tl_h(i,j)
!^
              ad_zeta_new(i,j)=ad_zeta_new(i,j)+ad_Dnew(i,j)
              ad_h(i,j)=ad_h(i,j)+ad_Dnew(i,j)
              ad_Dnew(i,j)=0.0_r8
!^            tl_zeta_new(i,j)=tl_zeta_new(i,j)*rmask(i,j)
!^
              ad_zeta_new(i,j)=ad_zeta_new(i,j)*rmask(i,j)
!^            tl_zeta_new(i,j)=tl_zeta(i,j,kstp)+                       &
!^   &                         pm(i,j)*pn(i,j)*(tl_cff+                 &
!^   &                                          cff2*tl_rzeta(i,j,kstp)-&
!^   &                                          cff3*tl_rzeta(i,j,ptsk))
!^
              adfac=pm(i,j)*pn(i,j)*ad_zeta_new(i,j)
              ad_zeta(i,j,kstp)=ad_zeta(i,j,kstp)+ad_zeta_new(i,j)
              ad_cff=ad_cff+adfac
              ad_rzeta(i,j,kstp)=ad_rzeta(i,j,kstp)+adfac*cff2
              ad_rzeta(i,j,ptsk)=-adfac*cff3
              ad_zeta_new(i,j)=0.0_r8
!^            tl_cff=cff1*((tl_DUon(i,j)-tl_DUon(i+1,j))+               &
!^   &                     (tl_DVom(i,j)-tl_DVom(i,j+1)))
!^
              adfac=cff1*ad_cff
              ad_DUon(i  ,j  )=ad_DUon(i  ,j  )+adfac
              ad_DUon(i+1,j  )=ad_DUon(i+1,j  )-adfac
              ad_DVom(i  ,j  )=ad_DVom(i  ,j  )+adfac
              ad_DVom(i  ,j+1)=ad_DVom(i  ,j+1)-adfac
              ad_cff=0.0_r8
            END DO
          END DO
        END IF
      END IF STEP_LOOP
!
!-----------------------------------------------------------------------
!  Compute adjoint time averaged fields over all short timesteps.
!-----------------------------------------------------------------------
!
!  After all fast time steps are completed, recompute S-coordinate
!  surfaces according to the new free surface field. Apply boundary
!  conditions to time averaged fields.
!
      IF ((iif(ng).eq.(nfast(ng)+1)).and.PREDICTOR_2D_STEP(ng)) THEN
!^      CALL mp_exchange2d (ng, tile, iTLM, 3,                          &
!^   &                      LBi, UBi, LBj, UBj,                         &
!^   &                      NghostPoints,                               &
!^   &                      EWperiodic(ng), NSperiodic(ng),             &
!^   &                      tl_Zt_avg1, tl_DU_avg1, tl_DV_avg1)
!^
        CALL ad_mp_exchange2d (ng, tile, iADM, 3,                       &
     &                         LBi, UBi, LBj, UBj,                      &
     &                         NghostPoints,                            &
     &                         EWperiodic(ng), NSperiodic(ng),          &
     &                         ad_Zt_avg1, ad_DU_avg1, ad_DV_avg1)
        IF (EWperiodic(ng).or.NSperiodic(ng)) THEN
!^        CALL exchange_v2d_tile (ng, tile,                             &
!^   &                            LBi, UBi, LBj, UBj,                   &
!^   &                            tl_DV_avg1)
!^
          CALL ad_exchange_v2d_tile (ng, tile,                          &
     &                               LBi, UBi, LBj, UBj,                &
     &                               ad_DV_avg1)
!^        CALL exchange_u2d_tile (ng, tile,                             &
!^   &                            LBi, UBi, LBj, UBj,                   &
!^   &                            tl_DU_avg1)
!^
          CALL ad_exchange_u2d_tile (ng, tile,                          &
     &                               LBi, UBi, LBj, UBj,                &
     &                               ad_DU_avg1)
!^        CALL exchange_r2d_tile (ng, tile,                             &
!^   &                            LBi, UBi, LBj, UBj,                   &
!^   &                            tl_Zt_avg1)
!^
          CALL ad_exchange_r2d_tile (ng, tile,                          &
     &                               LBi, UBi, LBj, UBj,                &
     &                               ad_Zt_avg1)
        END IF
      END IF
!
!  Compute time-averaged fields.
!
      IF (PREDICTOR_2D_STEP(ng)) THEN
        IF (iif(ng).eq.1) THEN
!
!  Reset arrays for 2D fields averaged within the short time-steps.
!
          cff2=(-1.0_r8/12.0_r8)*weight(2,iif(ng)+1,ng)
          DO j=JstrR,JendR
            DO i=IstrR,IendR
!^            tl_Zt_avg1(i,j)=0.0_r8
!^
              ad_Zt_avg1(i,j)=0.0_r8
            END DO
            DO i=Istr,IendR
!^            tl_DU_avg2(i,j)=cff2*tl_DUon(i,j)
!^
              ad_DUon(i,j)=ad_DUon(i,j)+cff2*ad_DU_avg2(i,j)
              ad_DU_avg2(i,j)=0.0_r8
!^            tl_DU_avg1(i,j)=0.0_r8
!^
              ad_DU_avg1(i,j)=0.0_r8
            END DO
          END DO
          DO j=Jstr,JendR
            DO i=IstrR,IendR
!^            tl_DV_avg2(i,j)=cff2*tl_DVom(i,j)
!^
              ad_DVom(i,j)=ad_DVom(i,j)+cff2*ad_DV_avg2(i,j)
              ad_DV_avg2(i,j)=0.0_r8
!^            tl_DV_avg1(i,j)=0.0_r8
!^
              ad_DV_avg1(i,j)=0.0_r8
            END DO
          END DO
        ELSE
!
!  Accumulate field averages of previous time-step after they are
!  computed in the previous corrector step, updated their boundaries,
!  and synchronized.
!
          cff1=weight(1,iif(ng)-1,ng)
          cff2=(8.0_r8/12.0_r8)*weight(2,iif(ng)  ,ng)-                 &
     &         (1.0_r8/12.0_r8)*weight(2,iif(ng)+1,ng)
          DO j=JstrR,JendR
            DO i=IstrR,IendR
!^            tl_Zt_avg1(i,j)=tl_Zt_avg1(i,j)+cff1*tl_zeta(i,j,krhs)
!^
              ad_zeta(i,j,krhs)=ad_zeta(i,j,krhs)+cff1*ad_Zt_avg1(i,j)
            END DO
            DO i=Istr,IendR
!^            tl_DU_avg2(i,j)=tl_DU_avg2(i,j)+cff2*tl_DUon(i,j)
!^
              ad_DUon(i,j)=ad_DUon(i,j)+                                &
     &                     cff2*ad_DU_avg2(i,j)
!^            tl_DU_avg1(i,j)=tl_DU_avg1(i,j)+cff1*tl_DUon(i,j)
!^
              ad_DUon(i,j)=ad_DUon(i,j)+                                &
     &                     cff1*ad_DU_avg1(i,j)
            END DO
          END DO
          DO j=Jstr,JendR
            DO i=IstrR,IendR
!^            tl_DV_avg2(i,j)=tl_DV_avg2(i,j)+cff2*tl_DVom(i,j)
!^
              ad_DVom(i,j)=ad_DVom(i,j)+                                &
     &                     cff2*ad_DV_avg2(i,j)
!^            tl_DV_avg1(i,j)=tl_DV_avg1(i,j)+cff1*tl_DVom(i,j)
!^
              ad_DVom(i,j)=ad_DVom(i,j)+                                &
     &                     cff1*ad_DV_avg1(i,j)
            END DO
          END DO
        END IF
      ELSE
        IF (iif(ng).eq.1) THEN
          cff2=weight(2,iif(ng),ng)
        ELSE
          cff2=(5.0_r8/12.0_r8)*weight(2,iif(ng),ng)
        END IF
        DO j=JstrR,JendR
          DO i=Istr,IendR
!^          tl_DV_avg2(i,j)=tl_DV_avg2(i,j)+cff2*tl_DVom(i,j)
!^
            ad_DVom(i,j)=ad_DVom(i,j)+cff2*ad_DV_avg2(i,j)
          END DO
        END DO
        DO j=Jstr,JendR
          DO i=IstrR,IendR
!^          tl_DU_avg2(i,j)=tl_DU_avg2(i,j)+cff2*tl_DUon(i,j)
!^
            ad_DUon(i,j)=ad_DUon(i,j)+cff2*ad_DU_avg2(i,j)
          END DO
        END DO
      END IF
!
!-----------------------------------------------------------------------
!  Compute total depth (m) and vertically integrated mass fluxes.
!-----------------------------------------------------------------------
!
!  Set vertically integrated mass fluxes DUon and DVom along the open
!  boundaries in such a way that the integral volume is conserved.
!
      IF (ANY(ad_VolCons(:,ng))) THEN
!^      CALL tl_set_DUV_bc_tile (ng, tile,                              &
!^   &                           LBi, UBi, LBj, UBj,                    &
!^   &                           IminS, ImaxS, JminS, JmaxS,            &
!^   &                           krhs,                                  &
!^   &                           umask, vmask,                          &
!^   &                           om_v, on_u, ubar, vbar,                &
!^   &                           tl_ubar, tl_vbar,                      &
!^   &                           Drhs, DUon, DVom,                      &
!^   &                           tl_Drhs, tl_DUon, tl_DVom)
!^
        CALL ad_set_DUV_bc_tile (ng, tile,                              &
     &                           LBi, UBi, LBj, UBj,                    &
     &                           IminS, ImaxS, JminS, JmaxS,            &
     &                           krhs,                                  &
     &                           umask, vmask,                          &
     &                           om_v, on_u, ubar, vbar,                &
     &                           ad_ubar, ad_vbar,                      &
     &                           Drhs, DUon, DVom,                      &
     &                           ad_Drhs, ad_DUon, ad_DVom)
      END IF
!
!  In distributed-memory, the I- and J-ranges are different and a
!  special exchange is done to avoid having three ghost points for
!  high order numerical stencils. Notice that a private array is
!  passed below to the exchange routine. It also applies periodic
!  boundary conditions, if appropriate and no partitions in I- or
!  J-directions.
!
!^    CALL mp_exchange2d (ng, tile, iTLM, 2,                            &
!^   &                    IminS, ImaxS, JminS, JmaxS,                   &
!^   &                    NghostPoints,                                 &
!^   &                    EWperiodic(ng), NSperiodic(ng),               &
!^   &                    tl_DUon, tl_DVom)
!^
      CALL ad_mp_exchange2d (ng, tile, iADM, 2,                         &
     &                       IminS, ImaxS, JminS, JmaxS,                &
     &                       NghostPoints,                              &
     &                       EWperiodic(ng), NSperiodic(ng),            &
     &                       ad_DUon, ad_DVom)
      IF (EWperiodic(ng).or.NSperiodic(ng)) THEN
!^      CALL exchange_u2d_tile (ng, tile,                               &
!^   &                          IminS, ImaxS, JminS, JmaxS,             &
!^   &                          tl_DUon)
!^
        CALL ad_exchange_u2d_tile (ng, tile,                            &
     &                             IminS, ImaxS, JminS, JmaxS,          &
     &                             ad_DUon)
!^      CALL exchange_v2d_tile (ng, tile,                               &
!^   &                          IminS, ImaxS, JminS, JmaxS,             &
!^   &                          tl_DVom)
!^
        CALL ad_exchange_v2d_tile (ng, tile,                            &
     &                             IminS, ImaxS, JminS, JmaxS,          &
     &                             ad_DVom)
      END IF
!
!  Compute adjoint adjoint vertically integrated mass fluxes.
!
      DO j=JstrV-1,Jendp2
        DO i=IstrU-2,Iendp2
          cff=0.5_r8*om_v(i,j)
          cff1=cff*(Drhs(i,j)+Drhs(i,j-1))
!^        tl_DVom(i,j)=tl_vbar(i,j,krhs)*cff1+                          &
!^   &                 vbar(i,j,krhs)*tl_cff1
!^
          ad_cff1=ad_cff1+vbar(i,j,krhs)*ad_DVom(i,j)
          ad_vbar(i,j,krhs)=ad_vbar(i,j,krhs)+cff1*ad_DVom(i,j)
          ad_DVom(i,j)=0.0_r8
!^        tl_cff1=cff*(tl_Drhs(i,j)+tl_Drhs(i,j-1))
!^
          adfac=cff*ad_cff1
          ad_Drhs(i,j-1)=ad_Drhs(i,j-1)+adfac
          ad_Drhs(i,j  )=ad_Drhs(i,j  )+adfac
          ad_cff1=0.0_r8
        END DO
      END DO
      DO j=JstrV-2,Jendp2
        DO i=IstrU-1,Iendp2
          cff=0.5_r8*on_u(i,j)
          cff1=cff*(Drhs(i,j)+Drhs(i-1,j))
!^        tl_DUon(i,j)=tl_ubar(i,j,krhs)*cff1+                          &
!^   &                 ubar(i,j,krhs)*tl_cff1
!^
          ad_cff1=ad_cff1+ubar(i,j,krhs)*ad_DUon(i,j)
          ad_ubar(i,j,krhs)=ad_ubar(i,j,krhs)+cff1*ad_DUon(i,j)
          ad_DUon(i,j)=0.0_r8
!^        tl_cff1=cff*(tl_Drhs(i,j)+tl_Drhs(i-1,j))
!^
          adfac=cff*ad_cff1
          ad_Drhs(i-1,j)=ad_Drhs(i-1,j)+adfac
          ad_Drhs(i  ,j)=ad_Drhs(i  ,j)+adfac
          ad_cff1=0.0_r8
        END DO
      END DO
!
!  Compute adjoint total depth.
!
      DO j=JstrV-2,Jendp2
        DO i=IstrU-2,Iendp2
!^        tl_Drhs(i,j)=tl_zeta(i,j,krhs)+tl_h(i,j)
!^
          ad_zeta(i,j,krhs)=ad_zeta(i,j,krhs)+ad_Drhs(i,j)
          ad_h(i,j)=ad_h(i,j)+ad_Drhs(i,j)
          ad_Drhs(i,j)=0.0_r8
        END DO
      END DO
!
      RETURN
      END SUBROUTINE ad_step2d_tile
      END MODULE ad_step2d_mod