#include "cppdefs.h"
      MODULE set_tides_mod
#if defined NONLINEAR && (defined SSH_TIDES || defined UV_TIDES)
!
!svn $Id: set_tides.F 645 2013-01-22 23:21:54Z arango $
!================================================== Hernan G. Arango ===
!  Copyright (c) 2002-2013 The ROMS/TOMS Group        Robert Hetland   !
!    Licensed under a MIT/X style license                              !
!    See License_ROMS.txt                                              !
!=======================================================================
!                                                                      !
!  This routine adds tidal elevation (m) and tidal currents (m/s) to   !
!  sea surface height and 2D momentum climatologies, respectively.     !
!                                                                      !
!=======================================================================
!
      implicit none
!
      PRIVATE
      PUBLIC  :: set_tides
!
      CONTAINS
!
!***********************************************************************
      SUBROUTINE set_tides (ng, tile)
!***********************************************************************
!
      USE mod_param
# ifdef CLIMATOLOGY
      USE mod_clima
# endif
      USE mod_grid
      USE mod_tides
      USE mod_stepping
!start nilsmk 15-05-2013
# ifdef ADD_FS_INV_BARO
      USE mod_forces
# endif
!stop nilsmk 15-05-2013
!
!  Imported variable declarations.
!
      integer, intent(in) :: ng, tile
!
!  Local variable declarations.
!
# include "tile.h"
!
# ifdef PROFILE
      CALL wclock_on (ng, iNLM, 11)
# endif
      CALL set_tides_tile (ng, tile,                                    &
     &                     LBi, UBi, LBj, UBj,                          &
     &                     IminS, ImaxS, JminS, JmaxS,                  &
     &                     NTC(ng),                                     &
     &                     GRID(ng) % angler,                           &
# ifdef MASKING
     &                     GRID(ng) % rmask,                            &
     &                     GRID(ng) % umask,                            &
     &                     GRID(ng) % vmask,                            &
# endif
# ifdef SSH_TIDES
     &                     TIDES(ng) % SSH_Tamp,                        &
     &                     TIDES(ng) % SSH_Tphase,                      &
# endif
!start nilsmk 15-05-2013
# ifdef ADD_FS_INV_BARO
     &                     FORCES(ng) % Pair,                             &
# endif
!stop nilsmk 15-05-2013
# ifdef UV_TIDES
     &                     TIDES(ng) % UV_Tangle,                       &
     &                     TIDES(ng) % UV_Tphase,                       &
     &                     TIDES(ng) % UV_Tmajor,                       &
     &                     TIDES(ng) % UV_Tminor,                       &
# endif
# if defined AVERAGES_DETIDE && (defined SSH_TIDES || defined UV_TIDES)
     &                     TIDES(ng) % SinOmega,                        &
     &                     TIDES(ng) % CosOmega,                        &
# endif
# ifdef ZCLIMATOLOGY
     &                     CLIMA(ng) % ssh,                             &
# endif
# ifdef M2CLIMATOLOGY
     &                     CLIMA(ng) % ubarclm,                         &
     &                     CLIMA(ng) % vbarclm,                         &
# endif
     &                     TIDES(ng) % Tperiod)
# ifdef PROFILE
      CALL wclock_off (ng, iNLM, 11)
# endif

      RETURN
      END SUBROUTINE set_tides
!
!***********************************************************************
      SUBROUTINE set_tides_tile (ng, tile,                              &
     &                           LBi, UBi, LBj, UBj,                    &
     &                           IminS, ImaxS, JminS, JmaxS,            &
     &                           NTC,                                   &
     &                           angler,                                &
# ifdef MASKING
     &                           rmask, umask, vmask,                   &
# endif
# ifdef SSH_TIDES
     &                           SSH_Tamp, SSH_Tphase,                  &
# endif
!start nilsmk 15-05-2013
# ifdef ADD_FS_INV_BARO
     &                           pair,                                  &
# endif
!stop nilsmk 15-05-2013
# ifdef UV_TIDES
     &                           UV_Tangle, UV_Tphase,                  &
     &                           UV_Tmajor, UV_Tminor,                  &
# endif
# if defined AVERAGES_DETIDE && (defined SSH_TIDES || defined UV_TIDES)
     &                           SinOmega, CosOmega,                    &
# endif
# ifdef ZCLIMATOLOGY
     &                           ssh,                                   &
# endif
# ifdef M2CLIMATOLOGY
     &                           ubarclm, vbarclm,                      &
# endif
     &                           Tperiod)
!***********************************************************************
!
      USE mod_param
      USE mod_boundary
      USE mod_ncparam
      USE mod_scalars
!
# ifdef DISTRIBUTE
      USE distribute_mod, ONLY : mp_boundary
# endif
      USE exchange_2d_mod
# ifdef DISTRIBUTE
      USE mp_exchange_mod, ONLY : mp_exchange2d
# endif
!
!  Imported variables declarations.
!
      integer, intent(in) :: ng, tile
      integer, intent(in) :: LBi, UBi, LBj, UBj
      integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
      integer, intent(in) :: NTC
!
# ifdef ASSUMED_SHAPE
      real(r8), intent(in) :: angler(LBi:,LBj:)
#  ifdef MASKING
      real(r8), intent(in) :: rmask(LBi:,LBj:)
      real(r8), intent(in) :: umask(LBi:,LBj:)
      real(r8), intent(in) :: vmask(LBi:,LBj:)
#  endif
      real(r8), intent(in) :: Tperiod(MTC)
#  ifdef SSH_TIDES
      real(r8), intent(in) :: SSH_Tamp(LBi:,LBj:,:)
      real(r8), intent(in) :: SSH_Tphase(LBi:,LBj:,:)
#  endif
!start nilsmk 15-05-2013
# ifdef ADD_FS_INV_BARO
     real(r8), intent(in) :: pair(LBi:,LBj:)
# endif
!stop nilsmk 15-05-2013
#  ifdef UV_TIDES
      real(r8), intent(in) :: UV_Tangle(LBi:,LBj:,:)
      real(r8), intent(in) :: UV_Tmajor(LBi:,LBj:,:)
      real(r8), intent(in) :: UV_Tminor(LBi:,LBj:,:)
      real(r8), intent(in) :: UV_Tphase(LBi:,LBj:,:)
#  endif
#  if defined AVERAGES_DETIDE && (defined SSH_TIDES || defined UV_TIDES)
      real(r8), intent(inout) :: SinOmega(:)
      real(r8), intent(inout) :: CosOmega(:)
#  endif
#  ifdef ZCLIMATOLOGY
      real(r8), intent(inout) :: ssh(LBi:,LBj:)
#  endif
#  ifdef M2CLIMATOLOGY
      real(r8), intent(inout) :: ubarclm(LBi:,LBj:)
      real(r8), intent(inout) :: vbarclm(LBi:,LBj:)
#  endif
# else
      real(r8), intent(in) :: angler(LBi:UBi,LBj:UBj)
#  ifdef MASKING
      real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
      real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
      real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
#  endif
      real(r8), intent(in) :: Tperiod(MTC)
#  ifdef SSH_TIDES
      real(r8), intent(in) :: SSH_Tamp(LBi:UBi,LBj:UBj,MTC)
      real(r8), intent(in) :: SSH_Tphase(LBi:UBi,LBj:UBj,MTC)
#  endif
!start nilsmk 15-05-2013
# ifdef ADD_FS_INV_BARO
     real(r8), intent(in) :: pair(LBi:UBi,LBj:UBj,MTC)
# endif
!stop nilsmk 15-05-2013
#  ifdef UV_TIDES
      real(r8), intent(in) :: UV_Tangle(LBi:UBi,LBj:UBj,MTC)
      real(r8), intent(in) :: UV_Tmajor(LBi:UBi,LBj:UBj,MTC)
      real(r8), intent(in) :: UV_Tminor(LBi:UBi,LBj:UBj,MTC)
      real(r8), intent(in) :: UV_Tphase(LBi:UBi,LBj:UBj,MTC)
#  endif
#  if defined AVERAGES_DETIDE && (defined SSH_TIDES || defined UV_TIDES)
      real(r8), intent(inout) :: SinOmega(MTC)
      real(r8), intent(inout) :: CosOmega(MTC)
#  endif
#  ifdef ZCLIMATOLOGY
      real(r8), intent(inout) :: ssh(LBi:UBi,LBj:UBj)
#  endif
#  ifdef M2CLIMATOLOGY
      real(r8), intent(inout) :: ubarclm(LBi:UBi,LBj:UBj)
      real(r8), intent(inout) :: vbarclm(LBi:UBi,LBj:UBj)
#  endif
# endif
!
!  Local variables declarations.
!
      logical :: update

# ifdef DISTRIBUTE
      integer :: ILB, IUB, JLB, JUB
# endif
      integer :: i, itide, j

      real(r8) :: Cangle, Cphase, Sangle, Sphase
      real(r8) :: angle, cff, phase, omega, ramp
      real(r8) :: bry_cor, bry_pgr, bry_str, bry_val

      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Etide
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Utide
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Uwrk
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Vtide
      real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: Vwrk

# include "set_bounds.h"

# ifdef DISTRIBUTE
!
!  Lower and upper bounds for nontiled (global values) boundary arrays.
!
      ILB=BOUNDS(ng)%LBi(-1)
      IUB=BOUNDS(ng)%UBi(-1)
      JLB=BOUNDS(ng)%LBj(-1)
      JUB=BOUNDS(ng)%UBj(-1)
# endif
!
!  Set time-ramping parameter.
!
# ifdef RAMP_TIDES
      ramp=TANH((tdays(ng)-dstart)/1.0_r8)
# else
      ramp=1.0_r8
# endif
# if defined AVERAGES_DETIDE && (defined SSH_TIDES || defined UV_TIDES)
!
!-----------------------------------------------------------------------
!  Compute harmonic used to detide output fields.
!-----------------------------------------------------------------------
!
      cff=2.0_r8*pi*(time(ng)-tide_start*day2sec)
      DO itide=1,NTC
        IF (Tperiod(itide).gt.0.0_r8) THEN
          omega=cff/Tperiod(itide)
          SinOmega(itide)=SIN(omega)
          CosOmega(itide)=COS(omega)
        ELSE
          SinOmega(itide)=0.0_r8
          CosOmega(itide)=0.0_r8
        END IF
      END DO
# endif
# ifdef SSH_TIDES
!
!-----------------------------------------------------------------------
!  Add tidal elevation (m) to sea surface height climatology.
!-----------------------------------------------------------------------
!
      Etide(:,:)=0.0_r8
      cff=2.0_r8*pi*(time(ng)-tide_start*day2sec)
      DO itide=1,NTC
        IF (Tperiod(itide).gt.0.0_r8) THEN
          omega=cff/Tperiod(itide)
          DO j=JstrR,JendR
            DO i=IstrR,IendR
              Etide(i,j)=Etide(i,j)+                                    &
     &                   ramp*SSH_Tamp(i,j,itide)*                      &
     &                   COS(omega-SSH_Tphase(i,j,itide))
#  ifdef MASKING
              Etide(i,j)=Etide(i,j)*rmask(i,j)
#  endif
            END DO
          END DO
        END IF
      END DO
#  if defined ZCLIMATOLOGY && defined ADD_FSOBC
!
!  Add sub-tidal forcing and adjust climatology to include tides.
!
      DO j=JstrR,JendR
        DO i=IstrR,IendR
          ssh(i,j)=ssh(i,j)+Etide(i,j)
!start nilsmk 15-05-2013
# ifdef ADD_FS_INV_BARO
          ssh(i,j)=ssh(i,j)-((pair(i,j)-1013.5_r8)/100.0_r8) !adds 0.01m per hPa
# endif
!stop nilsmk 15-05-2013
        END DO
      END DO
      IF (EWperiodic(ng).or.NSperiodic(ng)) THEN
        CALL exchange_r2d_tile (ng, tile,                               &
     &                          LBi, UBi, LBj, UBj,                     &
     &                          ssh)
      END IF
#   ifdef DISTRIBUTE
      CALL mp_exchange2d (ng, tile, iNLM, 1,                            &
     &                    LBi, UBi, LBj, UBj,                           &
     &                    NghostPoints,                                 &
     &                    EWperiodic(ng), NSperiodic(ng),               &
     &                    ssh)
#   endif
#  endif
!
!  If appropriate, load tidal forcing into boundary arrays.  The "zeta"
!  boundary arrays are important for the Flather or reduced physics
!  boundary conditions for 2D momentum. To avoid having two boundary
!  points for these arrays, the values of "zeta_west" and "zeta_east"
!  are averaged at u-points.  Similarly, the values of "zeta_south"
!  and "zeta_north" is averaged at v-points. Noticed that these
!  arrays are also used for the clamped conditions for free-surface.
!  This averaging is less important for that type ob boundary
!  conditions.
!
      IF (LBC(iwest,isFsur,ng)%acquire.or.                              &
     &    LBC(iwest,isUbar,ng)%acquire.or.                              &
     &    LBC(iwest,isVbar,ng)%acquire) THEN
        update=.FALSE.
        IF (DOMAIN(ng)%Western_Edge(tile)) THEN
          DO j=JstrR,JendR
#  ifdef ADD_FSOBC
            BOUNDARY(ng)%zeta_west(j)=BOUNDARY(ng)%zeta_west(j)+        &
     &                                0.5_r8*(Etide(Istr-1,j)+          &
     &                                        Etide(Istr  ,j))
#  else
            BOUNDARY(ng)%zeta_west(j)=0.5_r8*(Etide(Istr-1,j)+          &
     &                                        Etide(Istr  ,j))
#  endif
!start nilsmk 15-05-2013
# ifdef ADD_FS_INV_BARO
            BOUNDARY(ng)%zeta_west(j)=BOUNDARY(ng)%zeta_west(j)-        &
     &                 (((0.5_r8*(pair(Istr-1,j)+pair(Istr,j)))-        &
     &                 1013.5_r8)/100.0_r8) !adds 0.01m per hPa
# endif
!stop nilsmk 15-05-2013
          END DO
          update=.TRUE.
        END IF
#  ifdef DISTRIBUTE
        CALL mp_boundary (ng, iNLM, JstrR, JendR,                       &
     &                    JLB, JUB, 1, 1, update,                       &
     &                    BOUNDARY(ng)%zeta_west)
#  endif
      END IF
!
      IF (LBC(ieast,isFsur,ng)%acquire.or.                              &
     &    LBC(ieast,isUbar,ng)%acquire.or.                              &
     &    LBC(ieast,isVbar,ng)%acquire) THEN
        update=.FALSE.
        IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN
          DO j=JstrR,JendR
#  ifdef ADD_FSOBC
            BOUNDARY(ng)%zeta_east(j)=BOUNDARY(ng)%zeta_east(j)+        &
     &                                0.5_r8*(Etide(Iend  ,j)+          &
     &                                        Etide(Iend+1,j))
#  else
            BOUNDARY(ng)%zeta_east(j)=0.5_r8*(Etide(Iend  ,j)+          &
     &                                        Etide(Iend+1,j))
#  endif
!start nilsmk 15-05-2013
# ifdef ADD_FS_INV_BARO
          BOUNDARY(ng)%zeta_east(j)=BOUNDARY(ng)%zeta_east(j)-          &
     &                 (((0.5_r8*(pair(Iend,j)+pair(Iend+1,j)))-        &
     &                 1013.5_r8)/100.0_r8) !adds 0.01m per hPa
# endif
!stop nilsmk 15-05-2013
          END DO
          update=.TRUE.
        END IF
#  ifdef DISTRIBUTE
        CALL mp_boundary (ng, iNLM, JstrR, JendR,                       &
     &                    JLB, JUB, 1, 1, update,                       &
     &                    BOUNDARY(ng)%zeta_east)
#  endif
      END IF
!
      IF (LBC(isouth,isFsur,ng)%acquire.or.                             &
     &    LBC(isouth,isUbar,ng)%acquire.or.                             &
     &    LBC(isouth,isVbar,ng)%acquire) THEN
        update=.FALSE.
        IF (DOMAIN(ng)%Southern_Edge(tile)) THEN
          DO i=IstrR,IendR
#  ifdef ADD_FSOBC
            BOUNDARY(ng)%zeta_south(i)=BOUNDARY(ng)%zeta_south(i)+      &
     &                                 0.5_r8*(Etide(i,Jstr-1)+         &
     &                                         Etide(i,Jstr  ))
#  else
            BOUNDARY(ng)%zeta_south(i)=0.5_r8*(Etide(i,Jstr-1)+         &
     &                                         Etide(i,Jstr  ))
#  endif
!start nilsmk 15-05-2013
# ifdef ADD_FS_INV_BARO
          BOUNDARY(ng)%zeta_south(i)=BOUNDARY(ng)%zeta_south(i)-          &
     &                  (((0.5_r8*(pair(i,Jstr-1)+pair(i,Jstr)))-         &
     &                  1013.5_r8)/100.0_r8) !adds 0.01m per hPa
# endif
!stop nilsmk 15-05-2013
          END DO
          update=.TRUE.
        END IF
#  ifdef DISTRIBUTE
        CALL mp_boundary (ng, iNLM, IstrR, IendR,                       &
     &                    ILB, IUB, 1, 1, update,                       &
     &                    BOUNDARY(ng)%zeta_south)
#  endif
      END IF
!
      IF (LBC(inorth,isFsur,ng)%acquire.or.                             &
     &    LBC(inorth,isUbar,ng)%acquire.or.                             &
     &    LBC(inorth,isVbar,ng)%acquire) THEN
        update=.FALSE.
        IF (DOMAIN(ng)%Northern_Edge(tile)) THEN
          DO i=IstrR,IendR
#  ifdef ADD_FSOBC
            BOUNDARY(ng)%zeta_north(i)=BOUNDARY(ng)%zeta_north(i)+      &
     &                                 0.5_r8*(Etide(i,Jend  )+         &
     &                                         Etide(i,Jend+1))
#  else
            BOUNDARY(ng)%zeta_north(i)=0.5_r8*(Etide(i,Jend  )+         &
     &                                         Etide(i,Jend+1))
#  endif
!start nilsmk 15-05-2013
# ifdef ADD_FS_INV_BARO
          BOUNDARY(ng)%zeta_north(i)=BOUNDARY(ng)%zeta_north(i)-          &
     &                  (((0.5_r8*(pair(i,Jend)+pair(i,Jend+1)))-         &
     &                  1013.5_r8)/100.0_r8) !adds 0.01m per hPa
# endif
!stop nilsmk 15-05-2013
          END DO
          update=.TRUE.
        END IF
#  ifdef DISTRIBUTE
        CALL mp_boundary (ng, iNLM, IstrR, IendR,                       &
     &                    ILB, IUB, 1, 1, update,                       &
     &                    BOUNDARY(ng)%zeta_north)
#  endif
      END IF
# endif

# if defined UV_TIDES
!
!-----------------------------------------------------------------------
!  Add tidal currents (m/s) to 2D momentum climatologies.
!-----------------------------------------------------------------------
!
      Utide(:,:)=0.0_r8
      Vtide(:,:)=0.0_r8
      cff=2.0_r8*pi*(time(ng)-tide_start*day2sec)
      DO itide=1,NTC
        IF (Tperiod(itide).gt.0.0_r8) THEN
          omega=cff/Tperiod(itide)
          DO j=MIN(JstrR,Jstr-1),JendR
            DO i=MIN(IstrR,Istr-1),IendR
              angle=UV_Tangle(i,j,itide)-angler(i,j)
              Cangle=COS(angle)
              Sangle=SIN(angle)
              phase=omega-UV_Tphase(i,j,itide)
              Cphase=COS(phase)
              Sphase=SIN(phase)
              Uwrk(i,j)=UV_Tmajor(i,j,itide)*Cangle*Cphase-             &
     &                  UV_Tminor(i,j,itide)*Sangle*Sphase
              Vwrk(i,j)=UV_Tmajor(i,j,itide)*Sangle*Cphase+             &
     &                  UV_Tminor(i,j,itide)*Cangle*Sphase
            END DO
          END DO
          DO j=JstrR,JendR
            DO i=Istr,IendR
              Utide(i,j)=Utide(i,j)+                                    &
     &                   ramp*0.5_r8*(Uwrk(i-1,j)+Uwrk(i,j))
#  ifdef MASKING
              Utide(i,j)=Utide(i,j)*umask(i,j)
#  endif
            END DO
          END DO
          DO j=Jstr,JendR
            DO i=IstrR,IendR
              Vtide(i,j)=(Vtide(i,j)+                                   &
     &                    ramp*0.5_r8*(Vwrk(i,j-1)+Vwrk(i,j)))
#  ifdef MASKING
              Vtide(i,j)=Vtide(i,j)*vmask(i,j)
#  endif
            END DO
          END DO
        END IF
      END DO
#  if defined M2CLIMATOLOGY && defined ADD_M2OBC
!
!  Add sub-tidal forcing and adjust climatology to include tides.
!
      DO j=JstrR,JendR
        DO i=Istr,IendR
          ubarclm(i,j)=ubarclm(i,j)+Utide(i,j)
        END DO
      END DO
      DO j=Jstr,JendR
        DO i=IstrR,IendR
          vbarclm(i,j)=vbarclm(i,j)+Vtide(i,j)
        END DO
      END DO
      IF (EWperiodic(ng).or.NSperiodic(ng)) THEN
        CALL exchange_u2d_tile (ng, tile,                               &
     &                          LBi, UBi, LBj, UBj,                     &
     &                          ubarclm)
        CALL exchange_v2d_tile (ng, tile,                               &
     &                          LBi, UBi, LBj, UBj,                     &
     &                          vbarclm)
      END IF
#   ifdef DISTRIBUTE
      CALL mp_exchange2d (ng, tile, iNLM, 2,                            &
     &                    LBi, UBi, LBj, UBj,                           &
     &                    NghostPoints,                                 &
     &                    EWperiodic(ng), NSperiodic(ng),               &
     &                    ubarclm, vbarclm)
#   endif
#  endif
!
!  If appropriate, load tidal forcing into boundary arrays.
!
      IF (LBC(iwest,isUbar,ng)%acquire.and.                             &
     &    LBC(iwest,isVbar,ng)%acquire) THEN
        update=.FALSE.
        IF (DOMAIN(ng)%Western_Edge(tile)) THEN
          DO j=JstrR,JendR
#  ifdef ADD_M2OBC
            BOUNDARY(ng)%ubar_west(j)=BOUNDARY(ng)%ubar_west(j)+        &
     &                                Utide(Istr,j)
#  else
            BOUNDARY(ng)%ubar_west(j)=Utide(Istr,j)
#  endif
          END DO
          DO j=Jstr,JendR
#  ifdef ADD_M2OBC
            BOUNDARY(ng)%vbar_west(j)=BOUNDARY(ng)%vbar_west(j)+        &
     &                                Vtide(Istr-1,j)
#  else
            BOUNDARY(ng)%vbar_west(j)=Vtide(Istr-1,j)
#  endif
          END DO
          update=.TRUE.
        END IF
#  ifdef DISTRIBUTE
        CALL mp_boundary (ng, iNLM, JstrR, JendR,                       &
     &                    JLB, JUB, 1, 1, update,                       &
     &                    BOUNDARY(ng)%ubar_west)
        CALL mp_boundary (ng, iNLM, Jstr,  JendR,                       &
     &                    JLB, JUB, 1, 1, update,                       &
     &                    BOUNDARY(ng)%vbar_west)
#  endif
      END IF
!
      IF (LBC(ieast,isUbar,ng)%acquire.and.                             &
     &    LBC(ieast,isVbar,ng)%acquire) THEN
        update=.FALSE.
        IF (DOMAIN(ng)%Eastern_Edge(tile)) THEN
          DO j=JstrR,JendR
#  ifdef ADD_M2OBC
            BOUNDARY(ng)%ubar_east(j)=BOUNDARY(ng)%ubar_east(j)+        &
     &                                Utide(Iend+1,j)
#  else
            BOUNDARY(ng)%ubar_east(j)=Utide(Iend+1,j)
#  endif
          END DO
          DO j=Jstr,JendR
#  ifdef ADD_M2OBC
            BOUNDARY(ng)%vbar_east(j)=BOUNDARY(ng)%vbar_east(j)+        &
     &                                Vtide(Iend+1,j)
#  else
            BOUNDARY(ng)%vbar_east(j)=Vtide(Iend+1,j)
#  endif
          END DO
          update=.TRUE.
        END IF
#  ifdef DISTRIBUTE
        CALL mp_boundary (ng, iNLM, JstrR, JendR,                       &
     &                    JLB, JUB, 1, 1, update,                       &
     &                    BOUNDARY(ng)%ubar_east)
        CALL mp_boundary (ng, iNLM, Jstr,  JendR,                       &
     &                    JLB, JUB, 1, 1, update,                       &
     &                    BOUNDARY(ng)%vbar_east)
#  endif
      END IF
!
      IF (LBC(isouth,isUbar,ng)%acquire.and.                            &
     &    LBC(isouth,isVbar,ng)%acquire) THEN
        update=.FALSE.
        IF (DOMAIN(ng)%Southern_Edge(tile)) THEN
          DO i=Istr,IendR
#  ifdef ADD_M2OBC
            BOUNDARY(ng)%ubar_south(i)=BOUNDARY(ng)%ubar_south(i)+      &
     &                                 Utide(i,Jstr-1)
#  else
            BOUNDARY(ng)%ubar_south(i)=Utide(i,Jstr-1)
#  endif
          END DO
          DO i=IstrR,IendR
#  ifdef ADD_M2OBC
            BOUNDARY(ng)%vbar_south(i)=BOUNDARY(ng)%vbar_south(i)+      &
     &                                 Vtide(i,Jstr)
#  else
            BOUNDARY(ng)%vbar_south(i)=Vtide(i,Jstr)
#  endif
          END DO
          update=.TRUE.
        END IF
#  ifdef DISTRIBUTE
        CALL mp_boundary (ng, iNLM, Istr,  IendR,                       &
     &                    ILB, IUB, 1, 1, update,                       &
     &                    BOUNDARY(ng)%ubar_south)
        CALL mp_boundary (ng, iNLM, IstrR, IendR,                       &
     &                    ILB, IUB, 1, 1, update,                       &
     &                    BOUNDARY(ng)%vbar_south)
#  endif
      END IF
!
      IF (LBC(inorth,isUbar,ng)%acquire.and.                            &
     &    LBC(inorth,isVbar,ng)%acquire) THEN
        update=.FALSE.
        IF (DOMAIN(ng)%Northern_Edge(tile)) THEN
          DO i=Istr,IendR
#  ifdef ADD_M2OBC
            BOUNDARY(ng)%ubar_north(i)=BOUNDARY(ng)%ubar_north(i)+      &
     &                                 Utide(i,Jend+1)
#  else
            BOUNDARY(ng)%ubar_north(i)=Utide(i,Jend+1)
#  endif
          END DO
          DO i=IstrR,IendR
#  ifdef ADD_M2OBC
            BOUNDARY(ng)%vbar_north(i)=BOUNDARY(ng)%vbar_north(i)+      &
     &                                 Vtide(i,Jend+1)
#  else
            BOUNDARY(ng)%vbar_north(i)=Vtide(i,Jend+1)
#  endif
          END DO
          update=.TRUE.
        END IF
#  ifdef DISTRIBUTE
        CALL mp_boundary (ng, iNLM, Istr,  IendR,                       &
     &                    ILB, IUB, 1, 1, update,                       &
     &                    BOUNDARY(ng)%ubar_north)
        CALL mp_boundary (ng, iNLM, IstrR, IendR,                       &
     &                    ILB, IUB, 1, 1, update,                       &
     &                    BOUNDARY(ng)%vbar_north)
#  endif
      END IF
# endif

      RETURN
      END SUBROUTINE set_tides_tile
#endif
      END MODULE set_tides_mod
