New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
dynspg_ts.F90 in branches/2017/dev_r8624_AGRIF3_VVL/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/2017/dev_r8624_AGRIF3_VVL/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_ts.F90 @ 8762

Last change on this file since 8762 was 8762, checked in by jchanut, 6 years ago

AGRIF + vvl: Final changes, update SETTE tests (these are ok except SAS) - #1965

  • Property svn:keywords set to Id
File size: 63.5 KB
RevLine 
[358]1MODULE dynspg_ts
2   !!======================================================================
[6140]3   !!                   ***  MODULE  dynspg_ts  ***
4   !! Ocean dynamics:  surface pressure gradient trend, split-explicit scheme
5   !!======================================================================
[1502]6   !! History :   1.0  ! 2004-12  (L. Bessieres, G. Madec)  Original code
7   !!              -   ! 2005-11  (V. Garnier, G. Madec)  optimization
8   !!              -   ! 2006-08  (S. Masson)  distributed restart using iom
9   !!             2.0  ! 2007-07  (D. Storkey) calls to BDY routines
10   !!              -   ! 2008-01  (R. Benshila)  change averaging method
11   !!             3.2  ! 2009-07  (R. Benshila, G. Madec) Complete revisit associated to vvl reactivation
[2528]12   !!             3.3  ! 2010-09  (D. Storkey, E. O'Dea) update for BDY for Shelf configurations
[2724]13   !!             3.3  ! 2011-03  (R. Benshila, R. Hordoir, P. Oddo) update calculation of ub_b
[4292]14   !!             3.5  ! 2013-07  (J. Chanut) Switch to Forward-backward time stepping
15   !!             3.6  ! 2013-11  (A. Coward) Update for z-tilde compatibility
[5930]16   !!             3.7  ! 2015-11  (J. Chanut) free surface simplification
[7646]17   !!              -   ! 2016-12  (G. Madec, E. Clementi) update for Stoke-Drift divergence
[2724]18   !!---------------------------------------------------------------------
[6140]19
[358]20   !!----------------------------------------------------------------------
[6140]21   !!   dyn_spg_ts     : compute surface pressure gradient trend using a time-splitting scheme
22   !!   dyn_spg_ts_init: initialisation of the time-splitting scheme
23   !!   ts_wgt         : set time-splitting weights for temporal averaging (or not)
24   !!   ts_rst         : read/write time-splitting fields in restart file
[358]25   !!----------------------------------------------------------------------
26   USE oce             ! ocean dynamics and tracers
27   USE dom_oce         ! ocean space and time domain
[888]28   USE sbc_oce         ! surface boundary condition: ocean
[6140]29   USE zdf_oce         ! Bottom friction coefts
[5120]30   USE sbcisf          ! ice shelf variable (fwfisf)
[6140]31   USE sbcapr          ! surface boundary condition: atmospheric pressure
32   USE dynadv    , ONLY: ln_dynadv_vec
[358]33   USE phycst          ! physical constants
34   USE dynvor          ! vorticity term
[6152]35   USE wet_dry         ! wetting/drying flux limter
[7646]36   USE bdy_oce         ! open boundary
[5930]37   USE bdytides        ! open boundary condition data
[3294]38   USE bdydyn2d        ! open boundary conditions on barotropic variables
[4292]39   USE sbctide         ! tides
40   USE updtide         ! tide potential
[7646]41   USE sbcwave         ! surface wave
[6140]42   !
43   USE in_out_manager  ! I/O manager
[358]44   USE lib_mpp         ! distributed memory computing library
45   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
46   USE prtctl          ! Print control
[2715]47   USE iom             ! IOM library
[4292]48   USE restart         ! only for lrst_oce
[3294]49   USE wrk_nemo        ! Memory Allocation
[4292]50   USE timing          ! Timing   
[6140]51   USE diatmb          ! Top,middle,bottom output
[4292]52#if defined key_agrif
53   USE agrif_opa_interp ! agrif
[8741]54   USE agrif_oce
[4292]55#endif
[4757]56#if defined key_asminc   
57   USE asminc          ! Assimilation increment
58#endif
[358]59
[6140]60
[358]61   IMPLICIT NONE
62   PRIVATE
63
[4292]64   PUBLIC dyn_spg_ts        ! routine called in dynspg.F90
65   PUBLIC dyn_spg_ts_alloc  !    "      "     "    "
66   PUBLIC dyn_spg_ts_init   !    "      "     "    "
[4496]67   PUBLIC ts_rst            !    "      "     "    "
[358]68
[4292]69   INTEGER, SAVE :: icycle  ! Number of barotropic sub-steps for each internal step nn_baro <= 2.5 nn_baro
70   REAL(wp),SAVE :: rdtbt   ! Barotropic time step
71
[6140]72   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:) ::   wgtbtp1, wgtbtp2   !: 1st & 2nd weights used in time filtering of barotropic fields
[4292]73
[7646]74   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::  zwz          !: ff_f/h at F points
[6140]75   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::  ftnw, ftne   !: triad of coriolis parameter
76   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::  ftsw, ftse   !: (only used with een vorticity scheme)
[508]77
[6140]78   !! Time filtered arrays at baroclinic time step:
[8741]79   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   un_adv , vn_adv     !: Advection fluxes at "now" barocl. step
[6140]80
[358]81   !! * Substitutions
82#  include "vectopt_loop_substitute.h90"
[2715]83   !!----------------------------------------------------------------------
[4292]84   !! NEMO/OPA 3.5 , NEMO Consortium (2013)
[5217]85   !! $Id$
[2715]86   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
87   !!----------------------------------------------------------------------
[358]88CONTAINS
89
[2715]90   INTEGER FUNCTION dyn_spg_ts_alloc()
91      !!----------------------------------------------------------------------
92      !!                  ***  routine dyn_spg_ts_alloc  ***
93      !!----------------------------------------------------------------------
[6140]94      INTEGER :: ierr(3)
[4292]95      !!----------------------------------------------------------------------
96      ierr(:) = 0
[6140]97      !
98      ALLOCATE( wgtbtp1(3*nn_baro), wgtbtp2(3*nn_baro), zwz(jpi,jpj), STAT=ierr(1) )
99      !
100      IF( ln_dynvor_een )   ALLOCATE( ftnw(jpi,jpj) , ftne(jpi,jpj) , & 
101         &                            ftsw(jpi,jpj) , ftse(jpi,jpj) , STAT=ierr(2) )
102         !
103      ALLOCATE( un_adv(jpi,jpj), vn_adv(jpi,jpj)                    , STAT=ierr(3) )
104      !
105      dyn_spg_ts_alloc = MAXVAL( ierr(:) )
106      !
[2715]107      IF( lk_mpp                )   CALL mpp_sum( dyn_spg_ts_alloc )
[5930]108      IF( dyn_spg_ts_alloc /= 0 )   CALL ctl_warn('dyn_spg_ts_alloc: failed to allocate arrays')
[2715]109      !
110   END FUNCTION dyn_spg_ts_alloc
111
[5836]112
[358]113   SUBROUTINE dyn_spg_ts( kt )
114      !!----------------------------------------------------------------------
115      !!
[6140]116      !! ** Purpose : - Compute the now trend due to the explicit time stepping
117      !!              of the quasi-linear barotropic system, and add it to the
118      !!              general momentum trend.
[358]119      !!
[6140]120      !! ** Method  : - split-explicit schem (time splitting) :
[4374]121      !!      Barotropic variables are advanced from internal time steps
122      !!      "n"   to "n+1" if ln_bt_fw=T
123      !!      or from
124      !!      "n-1" to "n+1" if ln_bt_fw=F
125      !!      thanks to a generalized forward-backward time stepping (see ref. below).
[358]126      !!
[4374]127      !! ** Action :
128      !!      -Update the filtered free surface at step "n+1"      : ssha
129      !!      -Update filtered barotropic velocities at step "n+1" : ua_b, va_b
[8741]130      !!      -Compute barotropic advective fluxes at step "n"     : un_adv, vn_adv
[4374]131      !!      These are used to advect tracers and are compliant with discrete
132      !!      continuity equation taken at the baroclinic time steps. This
133      !!      ensures tracers conservation.
[6140]134      !!      - (ua, va) momentum trend updated with barotropic component.
[358]135      !!
[6140]136      !! References : Shchepetkin and McWilliams, Ocean Modelling, 2005.
[358]137      !!---------------------------------------------------------------------
[1502]138      INTEGER, INTENT(in)  ::   kt   ! ocean time-step index
[2715]139      !
[4292]140      LOGICAL  ::   ll_fw_start        ! if true, forward integration
[4374]141      LOGICAL  ::   ll_init             ! if true, special startup of 2d equations
[6152]142      LOGICAL  ::   ll_tmp1, ll_tmp2            ! local logical variables used in W/D
[4292]143      INTEGER  ::   ji, jj, jk, jn        ! dummy loop indices
144      INTEGER  ::   ikbu, ikbv, noffset      ! local integers
[6140]145      INTEGER  ::   iktu, iktv               ! local integers
146      REAL(wp) ::   zmdi
[4292]147      REAL(wp) ::   zraur, z1_2dt_b, z2dt_bf    ! local scalars
[5930]148      REAL(wp) ::   zx1, zy1, zx2, zy2          !   -      -
149      REAL(wp) ::   z1_12, z1_8, z1_4, z1_2  !   -      -
150      REAL(wp) ::   zu_spg, zv_spg              !   -      -
151      REAL(wp) ::   zhura, zhvra          !   -      -
152      REAL(wp) ::   za0, za1, za2, za3    !   -      -
[3294]153      !
[5930]154      REAL(wp), POINTER, DIMENSION(:,:) :: zsshp2_e
[4292]155      REAL(wp), POINTER, DIMENSION(:,:) :: zu_trd, zv_trd, zu_frc, zv_frc, zssh_frc
[5930]156      REAL(wp), POINTER, DIMENSION(:,:) :: zwx, zwy, zhdiv
[4292]157      REAL(wp), POINTER, DIMENSION(:,:) :: zhup2_e, zhvp2_e, zhust_e, zhvst_e
158      REAL(wp), POINTER, DIMENSION(:,:) :: zsshu_a, zsshv_a
[4370]159      REAL(wp), POINTER, DIMENSION(:,:) :: zhf
[6152]160      REAL(wp), POINTER, DIMENSION(:,:) :: zcpx, zcpy                 ! Wetting/Dying gravity filter coef.
[358]161      !!----------------------------------------------------------------------
[3294]162      !
163      IF( nn_timing == 1 )  CALL timing_start('dyn_spg_ts')
164      !
[4374]165      !                                         !* Allocate temporary arrays
[6140]166      CALL wrk_alloc( jpi,jpj,   zsshp2_e, zhdiv )
167      CALL wrk_alloc( jpi,jpj,   zu_trd, zv_trd)
168      CALL wrk_alloc( jpi,jpj,   zwx, zwy, zssh_frc, zu_frc, zv_frc)
169      CALL wrk_alloc( jpi,jpj,   zhup2_e, zhvp2_e, zhust_e, zhvst_e)
[6152]170      CALL wrk_alloc( jpi,jpj,   zsshu_a, zsshv_a                  )
[6140]171      CALL wrk_alloc( jpi,jpj,   zhf )
[7646]172      IF( ln_wd ) CALL wrk_alloc( jpi, jpj, zcpx, zcpy )
[3294]173      !
[6140]174      zmdi=1.e+20                               !  missing data indicator for masking
[4292]175      !                                         !* Local constant initialization
176      z1_12 = 1._wp / 12._wp 
177      z1_8  = 0.125_wp                                   
178      z1_4  = 0.25_wp
179      z1_2  = 0.5_wp     
180      zraur = 1._wp / rau0
[6140]181      !                                            ! reciprocal of baroclinic time step
182      IF( kt == nit000 .AND. neuler == 0 ) THEN   ;   z2dt_bf =          rdt
183      ELSE                                        ;   z2dt_bf = 2.0_wp * rdt
[4292]184      ENDIF
185      z1_2dt_b = 1.0_wp / z2dt_bf 
186      !
[6140]187      ll_init     = ln_bt_av                       ! if no time averaging, then no specific restart
[4292]188      ll_fw_start = .FALSE.
[6140]189      !                                            ! time offset in steps for bdy data update
190      IF( .NOT.ln_bt_fw ) THEN   ;   noffset = - nn_baro
191      ELSE                       ;   noffset =   0 
192      ENDIF
[4292]193      !
194      IF( kt == nit000 ) THEN                !* initialisation
[508]195         !
[358]196         IF(lwp) WRITE(numout,*)
197         IF(lwp) WRITE(numout,*) 'dyn_spg_ts : surface pressure gradient trend'
198         IF(lwp) WRITE(numout,*) '~~~~~~~~~~   free surface with time splitting'
[4354]199         IF(lwp) WRITE(numout,*)
[1502]200         !
[6140]201         IF( neuler == 0 )   ll_init=.TRUE.
[1502]202         !
[6140]203         IF( ln_bt_fw .OR. neuler == 0 ) THEN
204            ll_fw_start =.TRUE.
205            noffset     = 0
[4292]206         ELSE
[6140]207            ll_fw_start =.FALSE.
[4292]208         ENDIF
209         !
210         ! Set averaging weights and cycle length:
[6140]211         CALL ts_wgt( ln_bt_av, ll_fw_start, icycle, wgtbtp1, wgtbtp2 )
[4292]212         !
213      ENDIF
214      !
215      ! Set arrays to remove/compute coriolis trend.
216      ! Do it once at kt=nit000 if volume is fixed, else at each long time step.
217      ! Note that these arrays are also used during barotropic loop. These are however frozen
[4374]218      ! although they should be updated in the variable volume case. Not a big approximation.
[4292]219      ! To remove this approximation, copy lines below inside barotropic loop
[4374]220      ! and update depths at T-F points (ht and zhf resp.) at each barotropic time step
[4292]221      !
[6140]222      IF( kt == nit000 .OR. .NOT.ln_linssh ) THEN
223         IF( ln_dynvor_een ) THEN               !==  EEN scheme  ==!
[7646]224            SELECT CASE( nn_een_e3f )              !* ff_f/e3 at F-point
[5836]225            CASE ( 0 )                                   ! original formulation  (masked averaging of e3t divided by 4)
226               DO jj = 1, jpjm1
227                  DO ji = 1, jpim1
[6140]228                     zwz(ji,jj) =   ( ht_n(ji  ,jj+1) + ht_n(ji+1,jj+1) +                    &
229                        &             ht_n(ji  ,jj  ) + ht_n(ji+1,jj  )   ) * 0.25_wp 
[7646]230                     IF( zwz(ji,jj) /= 0._wp )   zwz(ji,jj) = ff_f(ji,jj) / zwz(ji,jj)
[5836]231                  END DO
[5032]232               END DO
[5836]233            CASE ( 1 )                                   ! new formulation  (masked averaging of e3t divided by the sum of mask)
234               DO jj = 1, jpjm1
235                  DO ji = 1, jpim1
[6140]236                     zwz(ji,jj) =   ( ht_n(ji  ,jj+1) + ht_n(ji+1,jj+1) +                     &
237                        &             ht_n(ji  ,jj  ) + ht_n(ji+1,jj  )   )                   &
[5836]238                        &       / ( MAX( 1._wp, tmask(ji  ,jj+1, 1) + tmask(ji+1,jj+1, 1) +    &
[4292]239                        &                       tmask(ji  ,jj  , 1) + tmask(ji+1,jj  , 1) ) )
[7646]240                     IF( zwz(ji,jj) /= 0._wp )   zwz(ji,jj) = ff_f(ji,jj) / zwz(ji,jj)
[5836]241                  END DO
[4292]242               END DO
[5836]243            END SELECT
[4292]244            CALL lbc_lnk( zwz, 'F', 1._wp )
[5836]245            !
[7753]246            ftne(1,:) = 0._wp ; ftnw(1,:) = 0._wp ; ftse(1,:) = 0._wp ; ftsw(1,:) = 0._wp
[358]247            DO jj = 2, jpj
[5836]248               DO ji = 2, jpi
[4292]249                  ftne(ji,jj) = zwz(ji-1,jj  ) + zwz(ji  ,jj  ) + zwz(ji  ,jj-1)
250                  ftnw(ji,jj) = zwz(ji-1,jj-1) + zwz(ji-1,jj  ) + zwz(ji  ,jj  )
251                  ftse(ji,jj) = zwz(ji  ,jj  ) + zwz(ji  ,jj-1) + zwz(ji-1,jj-1)
252                  ftsw(ji,jj) = zwz(ji  ,jj-1) + zwz(ji-1,jj-1) + zwz(ji-1,jj  )
[358]253               END DO
254            END DO
[5836]255            !
256         ELSE                                !== all other schemes (ENE, ENS, MIX)
[7753]257            zwz(:,:) = 0._wp
258            zhf(:,:) = 0._wp
[7646]259           
260!!gm  assume 0 in both cases (xhich is almost surely WRONG ! ) as hvatf has been removed
261!!gm    A priori a better value should be something like :
262!!gm          zhf(i,j) = masked sum of  ht(i,j) , ht(i+1,j) , ht(i,j+1) , (i+1,j+1)
263!!gm                     divided by the sum of the corresponding mask
264!!gm
265!!           
266              IF ( .not. ln_sco ) THEN
267 
268   !!gm  agree the JC comment  : this should be done in a much clear way
269 
270   ! JC: It not clear yet what should be the depth at f-points over land in z-coordinate case
271   !     Set it to zero for the time being
272   !              IF( rn_hmin < 0._wp ) THEN    ;   jk = - INT( rn_hmin )                                      ! from a nb of level
273   !              ELSE                          ;   jk = MINLOC( gdepw_0, mask = gdepw_0 > rn_hmin, dim = 1 )  ! from a depth
274   !              ENDIF
275   !              zhf(:,:) = gdepw_0(:,:,jk+1)
276               ELSE
277                 !zhf(:,:) = hbatf(:,:)
278                 DO jj = 1, jpjm1
279                   DO ji = 1, jpim1
280                     zhf(ji,jj) = MAX( 0._wp,                                &
281                                & ( ht_0(ji  ,jj  )*tmask(ji  ,jj  ,1) +     &
282                                &   ht_0(ji+1,jj  )*tmask(ji+1,jj  ,1) +     &
283                                &   ht_0(ji  ,jj+1)*tmask(ji  ,jj+1,1) +     &
284                                &   ht_0(ji+1,jj+1)*tmask(ji+1,jj+1,1) ) /   &
285                                & ( tmask(ji  ,jj  ,1) + tmask(ji+1,jj  ,1) +&
286                                &   tmask(ji  ,jj+1,1) + tmask(ji+1,jj+1,1) +&
287                                &   rsmall  ) )
288                   END DO
289                 END DO
290              END IF
291 
292              DO jj = 1, jpjm1
[7753]293                 zhf(:,jj) = zhf(:,jj) * (1._wp- umask(:,jj,1) * umask(:,jj+1,1))
[7646]294              END DO
295!!gm end
[5836]296
[4292]297            DO jk = 1, jpkm1
298               DO jj = 1, jpjm1
[7753]299                  zhf(:,jj) = zhf(:,jj) + e3f_n(:,jj,jk) * umask(:,jj,jk) * umask(:,jj+1,jk)
[4292]300               END DO
301            END DO
[4370]302            CALL lbc_lnk( zhf, 'F', 1._wp )
[4292]303            ! JC: TBC. hf should be greater than 0
304            DO jj = 1, jpj
305               DO ji = 1, jpi
[4370]306                  IF( zhf(ji,jj) /= 0._wp )   zwz(ji,jj) = 1._wp / zhf(ji,jj) ! zhf is actually hf here but it saves an array
[4292]307               END DO
308            END DO
[7753]309            zwz(:,:) = ff_f(:,:) * zwz(:,:)
[358]310         ENDIF
[508]311      ENDIF
[1502]312      !
[4292]313      ! If forward start at previous time step, and centered integration,
314      ! then update averaging weights:
[5836]315      IF (.NOT.ln_bt_fw .AND.( neuler==0 .AND. kt==nit000+1 ) ) THEN
[4292]316         ll_fw_start=.FALSE.
317         CALL ts_wgt(ln_bt_av, ll_fw_start, icycle, wgtbtp1, wgtbtp2)
318      ENDIF
319                         
[358]320      ! -----------------------------------------------------------------------------
321      !  Phase 1 : Coupling between general trend and barotropic estimates (1st step)
322      ! -----------------------------------------------------------------------------
[1502]323      !     
[4292]324      !
[4354]325      !                                   !* e3*d/dt(Ua) (Vertically integrated)
[4292]326      !                                   ! --------------------------------------------------
[7753]327      zu_frc(:,:) = 0._wp
328      zv_frc(:,:) = 0._wp
[1502]329      !
330      DO jk = 1, jpkm1
[7753]331         zu_frc(:,:) = zu_frc(:,:) + e3u_n(:,:,jk) * ua(:,:,jk) * umask(:,:,jk)
332         zv_frc(:,:) = zv_frc(:,:) + e3v_n(:,:,jk) * va(:,:,jk) * vmask(:,:,jk)         
[1502]333      END DO
[4292]334      !
[7753]335      zu_frc(:,:) = zu_frc(:,:) * r1_hu_n(:,:)
336      zv_frc(:,:) = zv_frc(:,:) * r1_hv_n(:,:)
[4292]337      !
[7753]338      !
[1502]339      !                                   !* baroclinic momentum trend (remove the vertical mean trend)
[4292]340      DO jk = 1, jpkm1                    ! -----------------------------------------------------------
[1502]341         DO jj = 2, jpjm1
342            DO ji = fs_2, fs_jpim1   ! vector opt.
[4292]343               ua(ji,jj,jk) = ua(ji,jj,jk) - zu_frc(ji,jj) * umask(ji,jj,jk)
344               va(ji,jj,jk) = va(ji,jj,jk) - zv_frc(ji,jj) * vmask(ji,jj,jk)
[1502]345            END DO
[358]346         END DO
[1502]347      END DO
[7646]348     
349!!gm  Question here when removing the Vertically integrated trends, we remove the vertically integrated NL trends on momentum....
350!!gm  Is it correct to do so ?   I think so...
351     
352     
[4292]353      !                                   !* barotropic Coriolis trends (vorticity scheme dependent)
354      !                                   ! --------------------------------------------------------
[7753]355      zwx(:,:) = un_b(:,:) * hu_n(:,:) * e2u(:,:)        ! now fluxes
356      zwy(:,:) = vn_b(:,:) * hv_n(:,:) * e1v(:,:)
[1502]357      !
[358]358      IF( ln_dynvor_ene .OR. ln_dynvor_mix ) THEN      ! energy conserving or mixed scheme
359         DO jj = 2, jpjm1
360            DO ji = fs_2, fs_jpim1   ! vector opt.
[5836]361               zy1 = ( zwy(ji,jj-1) + zwy(ji+1,jj-1) ) * r1_e1u(ji,jj)
362               zy2 = ( zwy(ji,jj  ) + zwy(ji+1,jj  ) ) * r1_e1u(ji,jj)
363               zx1 = ( zwx(ji-1,jj) + zwx(ji-1,jj+1) ) * r1_e2v(ji,jj)
364               zx2 = ( zwx(ji  ,jj) + zwx(ji  ,jj+1) ) * r1_e2v(ji,jj)
[358]365               ! energy conserving formulation for planetary vorticity term
[4292]366               zu_trd(ji,jj) = z1_4 * ( zwz(ji  ,jj-1) * zy1 + zwz(ji,jj) * zy2 )
367               zv_trd(ji,jj) =-z1_4 * ( zwz(ji-1,jj  ) * zx1 + zwz(ji,jj) * zx2 )
[358]368            END DO
369         END DO
[508]370         !
[4374]371      ELSEIF ( ln_dynvor_ens ) THEN                    ! enstrophy conserving scheme
[358]372         DO jj = 2, jpjm1
373            DO ji = fs_2, fs_jpim1   ! vector opt.
[4292]374               zy1 =   z1_8 * ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) &
[5836]375                 &            + zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) * r1_e1u(ji,jj)
[4292]376               zx1 = - z1_8 * ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) &
[5836]377                 &            + zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) * r1_e2v(ji,jj)
[4292]378               zu_trd(ji,jj)  = zy1 * ( zwz(ji  ,jj-1) + zwz(ji,jj) )
379               zv_trd(ji,jj)  = zx1 * ( zwz(ji-1,jj  ) + zwz(ji,jj) )
[358]380            END DO
381         END DO
[508]382         !
[5836]383      ELSEIF ( ln_dynvor_een ) THEN  ! enstrophy and energy conserving scheme
[358]384         DO jj = 2, jpjm1
385            DO ji = fs_2, fs_jpim1   ! vector opt.
[5836]386               zu_trd(ji,jj) = + z1_12 * r1_e1u(ji,jj) * (  ftne(ji,jj  ) * zwy(ji  ,jj  ) &
387                &                                         + ftnw(ji+1,jj) * zwy(ji+1,jj  ) &
388                &                                         + ftse(ji,jj  ) * zwy(ji  ,jj-1) &
389                &                                         + ftsw(ji+1,jj) * zwy(ji+1,jj-1) )
390               zv_trd(ji,jj) = - z1_12 * r1_e2v(ji,jj) * (  ftsw(ji,jj+1) * zwx(ji-1,jj+1) &
391                &                                         + ftse(ji,jj+1) * zwx(ji  ,jj+1) &
392                &                                         + ftnw(ji,jj  ) * zwx(ji-1,jj  ) &
393                &                                         + ftne(ji,jj  ) * zwx(ji  ,jj  ) )
[358]394            END DO
395         END DO
[508]396         !
[4292]397      ENDIF 
398      !
[1502]399      !                                   !* Right-Hand-Side of the barotropic momentum equation
400      !                                   ! ----------------------------------------------------
[6140]401      IF( .NOT.ln_linssh ) THEN                 ! Variable volume : remove surface pressure gradient
[6152]402        IF( ln_wd ) THEN                        ! Calculating and applying W/D gravity filters
[7646]403           DO jj = 2, jpjm1
404              DO ji = 2, jpim1 
405                ll_tmp1 = MIN(   sshn(ji,jj)               ,   sshn(ji+1,jj) ) >                &
406                     &    MAX( -ht_wd(ji,jj)               , -ht_wd(ji+1,jj) ) .AND.            &
407                     &    MAX(   sshn(ji,jj) + ht_wd(ji,jj),   sshn(ji+1,jj) + ht_wd(ji+1,jj) ) &
408                     &                                                         > rn_wdmin1 + rn_wdmin2
409                ll_tmp2 = ( ABS( sshn(ji+1,jj)             -   sshn(ji  ,jj))  > 1.E-12 ).AND.( &
410                     &    MAX(   sshn(ji,jj)               ,   sshn(ji+1,jj) ) >                &
411                     &    MAX( -ht_wd(ji,jj)               , -ht_wd(ji+1,jj) ) + rn_wdmin1 + rn_wdmin2 )
412   
[6152]413                IF(ll_tmp1) THEN
[7646]414                  zcpx(ji,jj) = 1.0_wp
415                ELSE IF(ll_tmp2) THEN
416                  ! no worries about  sshn(ji+1,jj) -  sshn(ji  ,jj) = 0, it won't happen ! here
417                  zcpx(ji,jj) = ABS( (sshn(ji+1,jj) + ht_wd(ji+1,jj) - sshn(ji,jj) - ht_wd(ji,jj)) &
418                              &    / (sshn(ji+1,jj) -  sshn(ji  ,jj)) )
[6152]419                ELSE
[7646]420                  zcpx(ji,jj) = 0._wp
[6152]421                END IF
[7646]422         
423                ll_tmp1 = MIN(   sshn(ji,jj)               ,   sshn(ji,jj+1) ) >                &
424                     &    MAX( -ht_wd(ji,jj)               , -ht_wd(ji,jj+1) ) .AND.            &
425                     &    MAX(   sshn(ji,jj) + ht_wd(ji,jj),   sshn(ji,jj+1) + ht_wd(ji,jj+1) ) &
426                     &                                                         > rn_wdmin1 + rn_wdmin2
427                ll_tmp2 = ( ABS( sshn(ji,jj)               -   sshn(ji,jj+1))  > 1.E-12 ).AND.( &
428                     &    MAX(   sshn(ji,jj)               ,   sshn(ji,jj+1) ) >                &
429                     &    MAX( -ht_wd(ji,jj)               , -ht_wd(ji,jj+1) ) + rn_wdmin1 + rn_wdmin2 )
430   
[6152]431                IF(ll_tmp1) THEN
[7646]432                  zcpy(ji,jj) = 1.0_wp
433                ELSE IF(ll_tmp2) THEN
434                  ! no worries about  sshn(ji,jj+1) -  sshn(ji,jj  ) = 0, it won't happen ! here
435                  zcpy(ji,jj) = ABS( (sshn(ji,jj+1) + ht_wd(ji,jj+1) - sshn(ji,jj) - ht_wd(ji,jj)) &
436                              &    / (sshn(ji,jj+1) -  sshn(ji,jj  )) )
[6152]437                ELSE
[7646]438                  zcpy(ji,jj) = 0._wp
439                END IF
440              END DO
[6152]441           END DO
[7646]442 
[6152]443           DO jj = 2, jpjm1
444              DO ji = 2, jpim1
[7646]445                 zu_trd(ji,jj) = zu_trd(ji,jj) - grav * ( sshn(ji+1,jj  ) - sshn(ji  ,jj ) )   &
446                        &                        * r1_e1u(ji,jj) * zcpx(ji,jj)
447                 zv_trd(ji,jj) = zv_trd(ji,jj) - grav * ( sshn(ji  ,jj+1) - sshn(ji  ,jj ) )   &
448                        &                        * r1_e2v(ji,jj) * zcpy(ji,jj)
[6152]449              END DO
450           END DO
451
452         ELSE
453
454           DO jj = 2, jpjm1
455              DO ji = fs_2, fs_jpim1   ! vector opt.
456                 zu_trd(ji,jj) = zu_trd(ji,jj) - grav * (  sshn(ji+1,jj  ) - sshn(ji  ,jj  )  ) * r1_e1u(ji,jj)
457                 zv_trd(ji,jj) = zv_trd(ji,jj) - grav * (  sshn(ji  ,jj+1) - sshn(ji  ,jj  )  ) * r1_e2v(ji,jj) 
458              END DO
459           END DO
460        ENDIF
461
[1502]462      ENDIF
[358]463
[4292]464      DO jj = 2, jpjm1                          ! Remove coriolis term (and possibly spg) from barotropic trend
[358]465         DO ji = fs_2, fs_jpim1
[6140]466             zu_frc(ji,jj) = zu_frc(ji,jj) - zu_trd(ji,jj) * ssumask(ji,jj)
467             zv_frc(ji,jj) = zv_frc(ji,jj) - zv_trd(ji,jj) * ssvmask(ji,jj)
[3294]468          END DO
[4292]469      END DO 
470      !
471      !                 ! Add bottom stress contribution from baroclinic velocities:     
472      IF (ln_bt_fw) THEN
473         DO jj = 2, jpjm1                         
474            DO ji = fs_2, fs_jpim1   ! vector opt.
475               ikbu = mbku(ji,jj)       
476               ikbv = mbkv(ji,jj)   
477               zwx(ji,jj) = un(ji,jj,ikbu) - un_b(ji,jj) ! NOW bottom baroclinic velocities
478               zwy(ji,jj) = vn(ji,jj,ikbv) - vn_b(ji,jj)
479            END DO
480         END DO
[3294]481      ELSE
[4292]482         DO jj = 2, jpjm1
483            DO ji = fs_2, fs_jpim1   ! vector opt.
484               ikbu = mbku(ji,jj)       
485               ikbv = mbkv(ji,jj)   
486               zwx(ji,jj) = ub(ji,jj,ikbu) - ub_b(ji,jj) ! BEFORE bottom baroclinic velocities
487               zwy(ji,jj) = vb(ji,jj,ikbv) - vb_b(ji,jj)
488            END DO
489         END DO
490      ENDIF
[1502]491      !
[4292]492      ! Note that the "unclipped" bottom friction parameter is used even with explicit drag
[6152]493      IF( ln_wd ) THEN
[7753]494        zu_frc(:,:) = zu_frc(:,:) + MAX(r1_hu_n(:,:) * bfrua(:,:),-1._wp / rdtbt) * zwx(:,:)
495        zv_frc(:,:) = zv_frc(:,:) + MAX(r1_hv_n(:,:) * bfrva(:,:),-1._wp / rdtbt) * zwy(:,:)
[6152]496      ELSE
[7753]497        zu_frc(:,:) = zu_frc(:,:) + r1_hu_n(:,:) * bfrua(:,:) * zwx(:,:)
498        zv_frc(:,:) = zv_frc(:,:) + r1_hv_n(:,:) * bfrva(:,:) * zwy(:,:)
[6152]499      END IF
500      !
[6140]501      !                                         ! Add top stress contribution from baroclinic velocities:     
[7646]502      IF( ln_bt_fw ) THEN
[6140]503         DO jj = 2, jpjm1
504            DO ji = fs_2, fs_jpim1   ! vector opt.
505               iktu = miku(ji,jj)
506               iktv = mikv(ji,jj)
507               zwx(ji,jj) = un(ji,jj,iktu) - un_b(ji,jj) ! NOW top baroclinic velocities
508               zwy(ji,jj) = vn(ji,jj,iktv) - vn_b(ji,jj)
509            END DO
510         END DO
511      ELSE
512         DO jj = 2, jpjm1
513            DO ji = fs_2, fs_jpim1   ! vector opt.
514               iktu = miku(ji,jj)
515               iktv = mikv(ji,jj)
516               zwx(ji,jj) = ub(ji,jj,iktu) - ub_b(ji,jj) ! BEFORE top baroclinic velocities
517               zwy(ji,jj) = vb(ji,jj,iktv) - vb_b(ji,jj)
518            END DO
519         END DO
520      ENDIF
521      !
522      ! Note that the "unclipped" top friction parameter is used even with explicit drag
[7753]523      zu_frc(:,:) = zu_frc(:,:) + r1_hu_n(:,:) * tfrua(:,:) * zwx(:,:)
524      zv_frc(:,:) = zv_frc(:,:) + r1_hv_n(:,:) * tfrva(:,:) * zwy(:,:)
[6140]525      !       
[4292]526      IF (ln_bt_fw) THEN                        ! Add wind forcing
[7753]527         zu_frc(:,:) =  zu_frc(:,:) + zraur * utau(:,:) * r1_hu_n(:,:)
528         zv_frc(:,:) =  zv_frc(:,:) + zraur * vtau(:,:) * r1_hv_n(:,:)
[2724]529      ELSE
[7753]530         zu_frc(:,:) =  zu_frc(:,:) + zraur * z1_2 * ( utau_b(:,:) + utau(:,:) ) * r1_hu_n(:,:)
531         zv_frc(:,:) =  zv_frc(:,:) + zraur * z1_2 * ( vtau_b(:,:) + vtau(:,:) ) * r1_hv_n(:,:)
[4292]532      ENDIF 
533      !
534      IF ( ln_apr_dyn ) THEN                    ! Add atm pressure forcing
535         IF (ln_bt_fw) THEN
536            DO jj = 2, jpjm1             
537               DO ji = fs_2, fs_jpim1   ! vector opt.
[5836]538                  zu_spg =  grav * (  ssh_ib (ji+1,jj  ) - ssh_ib (ji,jj) ) * r1_e1u(ji,jj)
539                  zv_spg =  grav * (  ssh_ib (ji  ,jj+1) - ssh_ib (ji,jj) ) * r1_e2v(ji,jj)
[4292]540                  zu_frc(ji,jj) = zu_frc(ji,jj) + zu_spg
541                  zv_frc(ji,jj) = zv_frc(ji,jj) + zv_spg
542               END DO
543            END DO
544         ELSE
545            DO jj = 2, jpjm1             
546               DO ji = fs_2, fs_jpim1   ! vector opt.
547                  zu_spg =  grav * z1_2 * (  ssh_ib (ji+1,jj  ) - ssh_ib (ji,jj)    &
[5836]548                      &                    + ssh_ibb(ji+1,jj  ) - ssh_ibb(ji,jj)  ) * r1_e1u(ji,jj)
[4292]549                  zv_spg =  grav * z1_2 * (  ssh_ib (ji  ,jj+1) - ssh_ib (ji,jj)    &
[5836]550                      &                    + ssh_ibb(ji  ,jj+1) - ssh_ibb(ji,jj)  ) * r1_e2v(ji,jj)
[4292]551                  zu_frc(ji,jj) = zu_frc(ji,jj) + zu_spg
552                  zv_frc(ji,jj) = zv_frc(ji,jj) + zv_spg
553               END DO
554            END DO
555         ENDIF
[2724]556      ENDIF
[4292]557      !                                   !* Right-Hand-Side of the barotropic ssh equation
558      !                                   ! -----------------------------------------------
559      !                                         ! Surface net water flux and rivers
560      IF (ln_bt_fw) THEN
[7753]561         zssh_frc(:,:) = zraur * ( emp(:,:) - rnf(:,:) + fwfisf(:,:) )
[4292]562      ELSE
[7753]563         zssh_frc(:,:) = zraur * z1_2 * (  emp(:,:) + emp_b(:,:) - rnf(:,:) - rnf_b(:,:)   &
564                &                        + fwfisf(:,:) + fwfisf_b(:,:)                     )
[4292]565      ENDIF
[7646]566      !
567      IF( ln_sdw ) THEN                         ! Stokes drift divergence added if necessary
[7753]568         zssh_frc(:,:) = zssh_frc(:,:) + div_sd(:,:)
[7646]569      ENDIF
570      !
[4292]571#if defined key_asminc
572      !                                         ! Include the IAU weighted SSH increment
573      IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN
[7753]574         zssh_frc(:,:) = zssh_frc(:,:) - ssh_iau(:,:)
[4292]575      ENDIF
576#endif
[5656]577      !                                   !* Fill boundary data arrays for AGRIF
578      !                                   ! ------------------------------------
[4486]579#if defined key_agrif
580         IF( .NOT.Agrif_Root() ) CALL agrif_dta_ts( kt )
581#endif
[4292]582      !
[358]583      ! -----------------------------------------------------------------------
[4292]584      !  Phase 2 : Integration of the barotropic equations
[358]585      ! -----------------------------------------------------------------------
[1502]586      !
587      !                                             ! ==================== !
588      !                                             !    Initialisations   !
[4292]589      !                                             ! ==================== ! 
[4370]590      ! Initialize barotropic variables:     
[4770]591      IF( ll_init )THEN
[7753]592         sshbb_e(:,:) = 0._wp
593         ubb_e  (:,:) = 0._wp
594         vbb_e  (:,:) = 0._wp
595         sshb_e (:,:) = 0._wp
596         ub_e   (:,:) = 0._wp
597         vb_e   (:,:) = 0._wp
[4700]598      ENDIF
[6152]599
[4700]600      !
[4370]601      IF (ln_bt_fw) THEN                  ! FORWARD integration: start from NOW fields                   
[7753]602         sshn_e(:,:) =    sshn(:,:)           
603         un_e  (:,:) =    un_b(:,:)           
604         vn_e  (:,:) =    vn_b(:,:)
605         !
606         hu_e  (:,:) =    hu_n(:,:)       
607         hv_e  (:,:) =    hv_n(:,:) 
608         hur_e (:,:) = r1_hu_n(:,:)   
609         hvr_e (:,:) = r1_hv_n(:,:)
[4370]610      ELSE                                ! CENTRED integration: start from BEFORE fields
[7753]611         sshn_e(:,:) =    sshb(:,:)
612         un_e  (:,:) =    ub_b(:,:)         
613         vn_e  (:,:) =    vb_b(:,:)
614         !
615         hu_e  (:,:) =    hu_b(:,:)       
616         hv_e  (:,:) =    hv_b(:,:) 
617         hur_e (:,:) = r1_hu_b(:,:)   
618         hvr_e (:,:) = r1_hv_b(:,:)
[4292]619      ENDIF
620      !
621      !
[4370]622      !
[4292]623      ! Initialize sums:
[7753]624      ua_b  (:,:) = 0._wp       ! After barotropic velocities (or transport if flux form)         
625      va_b  (:,:) = 0._wp
626      ssha  (:,:) = 0._wp       ! Sum for after averaged sea level
627      un_adv(:,:) = 0._wp       ! Sum for now transport issued from ts loop
628      vn_adv(:,:) = 0._wp
[1502]629      !                                             ! ==================== !
[4292]630      DO jn = 1, icycle                             !  sub-time-step loop  !
[1502]631         !                                          ! ==================== !
[3294]632         !                                                !* Update the forcing (BDY and tides)
[1502]633         !                                                !  ------------------
[4292]634         ! Update only tidal forcing at open boundaries
[7646]635         IF( ln_bdy      .AND. ln_tide )   CALL bdy_dta_tides( kt, kit=jn, time_offset= noffset+1 )
636         IF( ln_tide_pot .AND. ln_tide )   CALL upd_tide     ( kt, kit=jn, time_offset= noffset   )
[4292]637         !
638         ! Set extrapolation coefficients for predictor step:
639         IF ((jn<3).AND.ll_init) THEN      ! Forward           
640           za1 = 1._wp                                         
641           za2 = 0._wp                       
642           za3 = 0._wp                       
643         ELSE                              ! AB3-AM4 Coefficients: bet=0.281105
644           za1 =  1.781105_wp              ! za1 =   3/2 +   bet
645           za2 = -1.06221_wp               ! za2 = -(1/2 + 2*bet)
646           za3 =  0.281105_wp              ! za3 = bet
647         ENDIF
[367]648
[4292]649         ! Extrapolate barotropic velocities at step jit+0.5:
[7753]650         ua_e(:,:) = za1 * un_e(:,:) + za2 * ub_e(:,:) + za3 * ubb_e(:,:)
651         va_e(:,:) = za1 * vn_e(:,:) + za2 * vb_e(:,:) + za3 * vbb_e(:,:)
[4292]652
[6140]653         IF( .NOT.ln_linssh ) THEN                        !* Update ocean depth (variable volume case only)
[4292]654            !                                             !  ------------------
655            ! Extrapolate Sea Level at step jit+0.5:
[7753]656            zsshp2_e(:,:) = za1 * sshn_e(:,:)  + za2 * sshb_e(:,:) + za3 * sshbb_e(:,:)
[4292]657            !
658            DO jj = 2, jpjm1                                    ! Sea Surface Height at u- & v-points
659               DO ji = 2, fs_jpim1   ! Vector opt.
[6140]660                  zwx(ji,jj) = z1_2 * ssumask(ji,jj)  * r1_e1e2u(ji,jj)     &
[5836]661                     &              * ( e1e2t(ji  ,jj) * zsshp2_e(ji  ,jj)  &
662                     &              +   e1e2t(ji+1,jj) * zsshp2_e(ji+1,jj) )
[6140]663                  zwy(ji,jj) = z1_2 * ssvmask(ji,jj)  * r1_e1e2v(ji,jj)     &
[5836]664                     &              * ( e1e2t(ji,jj  ) * zsshp2_e(ji,jj  )  &
665                     &              +   e1e2t(ji,jj+1) * zsshp2_e(ji,jj+1) )
[4292]666               END DO
667            END DO
[5429]668            CALL lbc_lnk_multi( zwx, 'U', 1._wp, zwy, 'V', 1._wp )
[4292]669            !
[7753]670            zhup2_e (:,:) = hu_0(:,:) + zwx(:,:)                ! Ocean depth at U- and V-points
671            zhvp2_e (:,:) = hv_0(:,:) + zwy(:,:)
[4370]672         ELSE
[7753]673            zhup2_e (:,:) = hu_n(:,:)
674            zhvp2_e (:,:) = hv_n(:,:)
[4292]675         ENDIF
676         !                                                !* after ssh
[1502]677         !                                                !  -----------
[4292]678         ! One should enforce volume conservation at open boundaries here
679         ! considering fluxes below:
680         !
[7753]681         zwx(:,:) = e2u(:,:) * ua_e(:,:) * zhup2_e(:,:)         ! fluxes at jn+0.5
682         zwy(:,:) = e1v(:,:) * va_e(:,:) * zhvp2_e(:,:)
[4486]683         !
684#if defined key_agrif
[6140]685         ! Set fluxes during predictor step to ensure volume conservation
686         IF( .NOT.Agrif_Root() .AND. ln_bt_fw ) THEN
[4486]687            IF((nbondi == -1).OR.(nbondi == 2)) THEN
688               DO jj=1,jpj
689                  zwx(2,jj) = ubdy_w(jj) * e2u(2,jj)
690               END DO
691            ENDIF
692            IF((nbondi ==  1).OR.(nbondi == 2)) THEN
693               DO jj=1,jpj
694                  zwx(nlci-2,jj) = ubdy_e(jj) * e2u(nlci-2,jj)
695               END DO
696            ENDIF
697            IF((nbondj == -1).OR.(nbondj == 2)) THEN
698               DO ji=1,jpi
699                  zwy(ji,2) = vbdy_s(ji) * e1v(ji,2)
700               END DO
701            ENDIF
702            IF((nbondj ==  1).OR.(nbondj == 2)) THEN
703               DO ji=1,jpi
704                  zwy(ji,nlcj-2) = vbdy_n(ji) * e1v(ji,nlcj-2)
705               END DO
706            ENDIF
707         ENDIF
708#endif
[6152]709         IF( ln_wd ) CALL wad_lmt_bt(zwx, zwy, sshn_e, zssh_frc, rdtbt)
[4486]710         !
711         ! Sum over sub-time-steps to compute advective velocities
712         za2 = wgtbtp2(jn)
[7753]713         un_adv(:,:) = un_adv(:,:) + za2 * zwx(:,:) * r1_e2u(:,:)
714         vn_adv(:,:) = vn_adv(:,:) + za2 * zwy(:,:) * r1_e1v(:,:)
[4486]715         !
716         ! Set next sea level:
[4292]717         DO jj = 2, jpjm1                                 
[358]718            DO ji = fs_2, fs_jpim1   ! vector opt.
[4292]719               zhdiv(ji,jj) = (   zwx(ji,jj) - zwx(ji-1,jj)   &
[5836]720                  &             + zwy(ji,jj) - zwy(ji,jj-1)   ) * r1_e1e2t(ji,jj)
[358]721            END DO
722         END DO
[7753]723         ssha_e(:,:) = (  sshn_e(:,:) - rdtbt * ( zssh_frc(:,:) + zhdiv(:,:) )  ) * ssmask(:,:)
724         
[4292]725         CALL lbc_lnk( ssha_e, 'T',  1._wp )
726
[6140]727         ! Duplicate sea level across open boundaries (this is only cosmetic if linssh=T)
[7646]728         IF( ln_bdy )   CALL bdy_ssh( ssha_e )
[4292]729#if defined key_agrif
[6140]730         IF( .NOT.Agrif_Root() )   CALL agrif_ssh_ts( jn )
[4292]731#endif
732         
733         ! Sea Surface Height at u-,v-points (vvl case only)
[6140]734         IF( .NOT.ln_linssh ) THEN                               
[4292]735            DO jj = 2, jpjm1
736               DO ji = 2, jpim1      ! NO Vector Opt.
[6140]737                  zsshu_a(ji,jj) = z1_2 * ssumask(ji,jj) * r1_e1e2u(ji,jj)    &
738                     &              * ( e1e2t(ji  ,jj  )  * ssha_e(ji  ,jj  ) &
739                     &              +   e1e2t(ji+1,jj  )  * ssha_e(ji+1,jj  ) )
740                  zsshv_a(ji,jj) = z1_2 * ssvmask(ji,jj) * r1_e1e2v(ji,jj)    &
741                     &              * ( e1e2t(ji  ,jj  )  * ssha_e(ji  ,jj  ) &
742                     &              +   e1e2t(ji  ,jj+1)  * ssha_e(ji  ,jj+1) )
[4292]743               END DO
[358]744            END DO
[5429]745            CALL lbc_lnk_multi( zsshu_a, 'U', 1._wp, zsshv_a, 'V', 1._wp )
[4292]746         ENDIF   
747         !                                 
748         ! Half-step back interpolation of SSH for surface pressure computation:
749         !----------------------------------------------------------------------
750         IF ((jn==1).AND.ll_init) THEN
751           za0=1._wp                        ! Forward-backward
752           za1=0._wp                           
753           za2=0._wp
754           za3=0._wp
755         ELSEIF ((jn==2).AND.ll_init) THEN  ! AB2-AM3 Coefficients; bet=0 ; gam=-1/6 ; eps=1/12
756           za0= 1.0833333333333_wp          ! za0 = 1-gam-eps
757           za1=-0.1666666666666_wp          ! za1 = gam
758           za2= 0.0833333333333_wp          ! za2 = eps
759           za3= 0._wp             
760         ELSE                               ! AB3-AM4 Coefficients; bet=0.281105 ; eps=0.013 ; gam=0.0880
761           za0=0.614_wp                     ! za0 = 1/2 +   gam + 2*eps   
762           za1=0.285_wp                     ! za1 = 1/2 - 2*gam - 3*eps
763           za2=0.088_wp                     ! za2 = gam
764           za3=0.013_wp                     ! za3 = eps
765         ENDIF
[6140]766         !
[7753]767         zsshp2_e(:,:) = za0 *  ssha_e(:,:) + za1 *  sshn_e (:,:) &
768          &            + za2 *  sshb_e(:,:) + za3 *  sshbb_e(:,:)
[6152]769         IF( ln_wd ) THEN                   ! Calculating and applying W/D gravity filters
770           DO jj = 2, jpjm1
[7646]771              DO ji = 2, jpim1 
772                ll_tmp1 = MIN( zsshp2_e(ji,jj)               , zsshp2_e(ji+1,jj) ) >                &
773                     &    MAX(   -ht_wd(ji,jj)               ,   -ht_wd(ji+1,jj) ) .AND.            &
774                     &    MAX( zsshp2_e(ji,jj) + ht_wd(ji,jj), zsshp2_e(ji+1,jj) + ht_wd(ji+1,jj) ) &
775                     &                                                             > rn_wdmin1 + rn_wdmin2
776                ll_tmp2 = (ABS(zsshp2_e(ji,jj)               - zsshp2_e(ji+1,jj))  > 1.E-12 ).AND.( &
777                     &    MAX( zsshp2_e(ji,jj)               , zsshp2_e(ji+1,jj) ) >                &
778                     &    MAX(   -ht_wd(ji,jj)               ,   -ht_wd(ji+1,jj) ) + rn_wdmin1 + rn_wdmin2 )
779   
780                IF(ll_tmp1) THEN
781                  zcpx(ji,jj) = 1.0_wp
782                ELSE IF(ll_tmp2) THEN
783                  ! no worries about  zsshp2_e(ji+1,jj) - zsshp2_e(ji  ,jj) = 0, it won't happen ! here
784                  zcpx(ji,jj) = ABS( (zsshp2_e(ji+1,jj) +    ht_wd(ji+1,jj) - zsshp2_e(ji,jj) - ht_wd(ji,jj)) &
785                              &    / (zsshp2_e(ji+1,jj) - zsshp2_e(ji  ,jj)) )
786                ELSE
787                  zcpx(ji,jj) = 0._wp
788                END IF
789         
790                ll_tmp1 = MIN( zsshp2_e(ji,jj)               , zsshp2_e(ji,jj+1) ) >                &
791                     &    MAX(   -ht_wd(ji,jj)               ,   -ht_wd(ji,jj+1) ) .AND.            &
792                     &    MAX( zsshp2_e(ji,jj) + ht_wd(ji,jj), zsshp2_e(ji,jj+1) + ht_wd(ji,jj+1) ) &
793                     &                                                             > rn_wdmin1 + rn_wdmin2
794                ll_tmp2 = (ABS(zsshp2_e(ji,jj)               - zsshp2_e(ji,jj+1))  > 1.E-12 ).AND.( &
795                     &    MAX( zsshp2_e(ji,jj)               , zsshp2_e(ji,jj+1) ) >                &
796                     &    MAX(   -ht_wd(ji,jj)               ,   -ht_wd(ji,jj+1) ) + rn_wdmin1 + rn_wdmin2 )
797   
798                IF(ll_tmp1) THEN
799                  zcpy(ji,jj) = 1.0_wp
800                ELSE IF(ll_tmp2) THEN
801                  ! no worries about  zsshp2_e(ji,jj+1) - zsshp2_e(ji,jj  ) = 0, it won't happen ! here
802                  zcpy(ji,jj) = ABS( (zsshp2_e(ji,jj+1) +    ht_wd(ji,jj+1) - zsshp2_e(ji,jj) - ht_wd(ji,jj)) &
803                              &    / (zsshp2_e(ji,jj+1) - zsshp2_e(ji,jj  )) )
804                ELSE
805                  zcpy(ji,jj) = 0._wp
806                END IF
[6152]807              END DO
[7646]808           END DO
809         END IF
[1502]810         !
[4292]811         ! Compute associated depths at U and V points:
[6140]812         IF( .NOT.ln_linssh  .AND. .NOT.ln_dynadv_vec ) THEN   !* Vector form
[4292]813            !                                       
814            DO jj = 2, jpjm1                           
815               DO ji = 2, jpim1
[6140]816                  zx1 = z1_2 * ssumask(ji  ,jj) *  r1_e1e2u(ji  ,jj)    &
[5836]817                     &      * ( e1e2t(ji  ,jj  ) * zsshp2_e(ji  ,jj)    &
818                     &      +   e1e2t(ji+1,jj  ) * zsshp2_e(ji+1,jj  ) )
[6140]819                  zy1 = z1_2 * ssvmask(ji  ,jj) *  r1_e1e2v(ji  ,jj  )  &
[5836]820                     &       * ( e1e2t(ji ,jj  ) * zsshp2_e(ji  ,jj  )  &
821                     &       +   e1e2t(ji ,jj+1) * zsshp2_e(ji  ,jj+1) )
[4292]822                  zhust_e(ji,jj) = hu_0(ji,jj) + zx1 
823                  zhvst_e(ji,jj) = hv_0(ji,jj) + zy1
824               END DO
825            END DO
[6152]826
[4292]827         ENDIF
828         !
829         ! Add Coriolis trend:
[6140]830         ! zwz array below or triads normally depend on sea level with ln_linssh=F and should be updated
[4292]831         ! at each time step. We however keep them constant here for optimization.
832         ! Recall that zwx and zwy arrays hold fluxes at this stage:
833         ! zwx(:,:) = e2u(:,:) * ua_e(:,:) * zhup2_e(:,:)   ! fluxes at jn+0.5
834         ! zwy(:,:) = e1v(:,:) * va_e(:,:) * zhvp2_e(:,:)
835         !
[6140]836         IF( ln_dynvor_ene .OR. ln_dynvor_mix ) THEN     !==  energy conserving or mixed scheme  ==!
[358]837            DO jj = 2, jpjm1
838               DO ji = fs_2, fs_jpim1   ! vector opt.
[5836]839                  zy1 = ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) ) * r1_e1u(ji,jj)
840                  zy2 = ( zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) * r1_e1u(ji,jj)
841                  zx1 = ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) ) * r1_e2v(ji,jj)
842                  zx2 = ( zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) * r1_e2v(ji,jj)
[4292]843                  zu_trd(ji,jj) = z1_4 * ( zwz(ji  ,jj-1) * zy1 + zwz(ji,jj) * zy2 )
844                  zv_trd(ji,jj) =-z1_4 * ( zwz(ji-1,jj  ) * zx1 + zwz(ji,jj) * zx2 )
[358]845               END DO
846            END DO
[508]847            !
[6140]848         ELSEIF ( ln_dynvor_ens ) THEN                   !==  enstrophy conserving scheme  ==!
[358]849            DO jj = 2, jpjm1
850               DO ji = fs_2, fs_jpim1   ! vector opt.
[4292]851                  zy1 =   z1_8 * ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) &
[5836]852                   &             + zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) * r1_e1u(ji,jj)
[4292]853                  zx1 = - z1_8 * ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) &
[5836]854                   &             + zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) * r1_e2v(ji,jj)
[4292]855                  zu_trd(ji,jj)  = zy1 * ( zwz(ji  ,jj-1) + zwz(ji,jj) )
856                  zv_trd(ji,jj)  = zx1 * ( zwz(ji-1,jj  ) + zwz(ji,jj) )
[358]857               END DO
858            END DO
[508]859            !
[6140]860         ELSEIF ( ln_dynvor_een ) THEN                   !==  energy and enstrophy conserving scheme  ==!
[358]861            DO jj = 2, jpjm1
862               DO ji = fs_2, fs_jpim1   ! vector opt.
[5836]863                  zu_trd(ji,jj) = + z1_12 * r1_e1u(ji,jj) * (  ftne(ji,jj  ) * zwy(ji  ,jj  ) &
864                     &                                       + ftnw(ji+1,jj) * zwy(ji+1,jj  ) &
865                     &                                       + ftse(ji,jj  ) * zwy(ji  ,jj-1) & 
866                     &                                       + ftsw(ji+1,jj) * zwy(ji+1,jj-1) )
867                  zv_trd(ji,jj) = - z1_12 * r1_e2v(ji,jj) * (  ftsw(ji,jj+1) * zwx(ji-1,jj+1) & 
868                     &                                       + ftse(ji,jj+1) * zwx(ji  ,jj+1) &
869                     &                                       + ftnw(ji,jj  ) * zwx(ji-1,jj  ) & 
870                     &                                       + ftne(ji,jj  ) * zwx(ji  ,jj  ) )
[358]871               END DO
872            END DO
[508]873            !
[358]874         ENDIF
[4292]875         !
876         ! Add tidal astronomical forcing if defined
[7646]877         IF ( ln_tide .AND. ln_tide_pot ) THEN
[4292]878            DO jj = 2, jpjm1
879               DO ji = fs_2, fs_jpim1   ! vector opt.
[5836]880                  zu_spg = grav * ( pot_astro(ji+1,jj) - pot_astro(ji,jj) ) * r1_e1u(ji,jj)
881                  zv_spg = grav * ( pot_astro(ji,jj+1) - pot_astro(ji,jj) ) * r1_e2v(ji,jj)
[4292]882                  zu_trd(ji,jj) = zu_trd(ji,jj) + zu_spg
883                  zv_trd(ji,jj) = zv_trd(ji,jj) + zv_spg
884               END DO
885            END DO
886         ENDIF
887         !
888         ! Add bottom stresses:
[7753]889         zu_trd(:,:) = zu_trd(:,:) + bfrua(:,:) * un_e(:,:) * hur_e(:,:)
890         zv_trd(:,:) = zv_trd(:,:) + bfrva(:,:) * vn_e(:,:) * hvr_e(:,:)
[4292]891         !
[7753]892         ! Add top stresses:
893         zu_trd(:,:) = zu_trd(:,:) + tfrua(:,:) * un_e(:,:) * hur_e(:,:)
894         zv_trd(:,:) = zv_trd(:,:) + tfrva(:,:) * vn_e(:,:) * hvr_e(:,:)
895         !
[4292]896         ! Surface pressure trend:
[6152]897
898         IF( ln_wd ) THEN
899           DO jj = 2, jpjm1
900              DO ji = 2, jpim1 
901                 ! Add surface pressure gradient
902                 zu_spg = - grav * ( zsshp2_e(ji+1,jj) - zsshp2_e(ji,jj) ) * r1_e1u(ji,jj)
903                 zv_spg = - grav * ( zsshp2_e(ji,jj+1) - zsshp2_e(ji,jj) ) * r1_e2v(ji,jj)
904                 zwx(ji,jj) = zu_spg * zcpx(ji,jj) 
905                 zwy(ji,jj) = zv_spg * zcpy(ji,jj)
906              END DO
907           END DO
908         ELSE
909           DO jj = 2, jpjm1
910              DO ji = fs_2, fs_jpim1   ! vector opt.
911                 ! Add surface pressure gradient
912                 zu_spg = - grav * ( zsshp2_e(ji+1,jj) - zsshp2_e(ji,jj) ) * r1_e1u(ji,jj)
913                 zv_spg = - grav * ( zsshp2_e(ji,jj+1) - zsshp2_e(ji,jj) ) * r1_e2v(ji,jj)
914                 zwx(ji,jj) = zu_spg
915                 zwy(ji,jj) = zv_spg
916              END DO
917           END DO
918         END IF
919
[4292]920         !
921         ! Set next velocities:
[6140]922         IF( ln_dynadv_vec .OR. ln_linssh ) THEN   !* Vector form
[4292]923            DO jj = 2, jpjm1
924               DO ji = fs_2, fs_jpim1   ! vector opt.
[5930]925                  ua_e(ji,jj) = (                                 un_e(ji,jj)   & 
[4292]926                            &     + rdtbt * (                      zwx(ji,jj)   &
927                            &                                 + zu_trd(ji,jj)   &
928                            &                                 + zu_frc(ji,jj) ) & 
[6140]929                            &   ) * ssumask(ji,jj)
[358]930
[5930]931                  va_e(ji,jj) = (                                 vn_e(ji,jj)   &
[4292]932                            &     + rdtbt * (                      zwy(ji,jj)   &
933                            &                                 + zv_trd(ji,jj)   &
934                            &                                 + zv_frc(ji,jj) ) &
[6140]935                            &   ) * ssvmask(ji,jj)
[4292]936               END DO
937            END DO
[6140]938            !
939         ELSE                                      !* Flux form
[4292]940            DO jj = 2, jpjm1
941               DO ji = fs_2, fs_jpim1   ! vector opt.
[3294]942
[6152]943                  IF( ln_wd ) THEN
944                    zhura = MAX(hu_0(ji,jj) + zsshu_a(ji,jj), rn_wdmin1)
945                    zhvra = MAX(hv_0(ji,jj) + zsshv_a(ji,jj), rn_wdmin1)
946                  ELSE
947                    zhura = hu_0(ji,jj) + zsshu_a(ji,jj)
948                    zhvra = hv_0(ji,jj) + zsshv_a(ji,jj)
949                  END IF
950                  zhura = ssumask(ji,jj)/(zhura + 1._wp - ssumask(ji,jj))
951                  zhvra = ssvmask(ji,jj)/(zhvra + 1._wp - ssvmask(ji,jj))
952
[5930]953                  ua_e(ji,jj) = (                hu_e(ji,jj)  *   un_e(ji,jj)   & 
[4292]954                            &     + rdtbt * ( zhust_e(ji,jj)  *    zwx(ji,jj)   & 
955                            &               + zhup2_e(ji,jj)  * zu_trd(ji,jj)   &
[6140]956                            &               +    hu_n(ji,jj)  * zu_frc(ji,jj) ) &
[4292]957                            &   ) * zhura
[358]958
[5930]959                  va_e(ji,jj) = (                hv_e(ji,jj)  *   vn_e(ji,jj)   &
[4292]960                            &     + rdtbt * ( zhvst_e(ji,jj)  *    zwy(ji,jj)   &
961                            &               + zhvp2_e(ji,jj)  * zv_trd(ji,jj)   &
[6140]962                            &               +    hv_n(ji,jj)  * zv_frc(ji,jj) ) &
[4292]963                            &   ) * zhvra
[592]964               END DO
965            END DO
[4292]966         ENDIF
967         !
[6140]968         IF( .NOT.ln_linssh ) THEN                     !* Update ocean depth (variable volume case only)
[6152]969            IF( ln_wd ) THEN
[7753]970              hu_e (:,:) = MAX(hu_0(:,:) + zsshu_a(:,:), rn_wdmin1)
971              hv_e (:,:) = MAX(hv_0(:,:) + zsshv_a(:,:), rn_wdmin1)
[6152]972            ELSE
[7753]973              hu_e (:,:) = hu_0(:,:) + zsshu_a(:,:)
974              hv_e (:,:) = hv_0(:,:) + zsshv_a(:,:)
[6152]975            END IF
[7753]976            hur_e(:,:) = ssumask(:,:) / ( hu_e(:,:) + 1._wp - ssumask(:,:) )
977            hvr_e(:,:) = ssvmask(:,:) / ( hv_e(:,:) + 1._wp - ssvmask(:,:) )
[1502]978            !
[1438]979         ENDIF
[6140]980         !                                             !* domain lateral boundary
981         CALL lbc_lnk_multi( ua_e, 'U', -1._wp, va_e , 'V', -1._wp )
[4292]982         !
[6140]983         !                                                 ! open boundaries
[7646]984         IF( ln_bdy )   CALL bdy_dyn2d( jn, ua_e, va_e, un_e, vn_e, hur_e, hvr_e, ssha_e )
[4486]985#if defined key_agrif                                                           
986         IF( .NOT.Agrif_Root() )  CALL agrif_dyn_ts( jn )  ! Agrif
[4292]987#endif
988         !                                             !* Swap
989         !                                             !  ----
[7753]990         ubb_e  (:,:) = ub_e  (:,:)
991         ub_e   (:,:) = un_e  (:,:)
992         un_e   (:,:) = ua_e  (:,:)
993         !
994         vbb_e  (:,:) = vb_e  (:,:)
995         vb_e   (:,:) = vn_e  (:,:)
996         vn_e   (:,:) = va_e  (:,:)
997         !
998         sshbb_e(:,:) = sshb_e(:,:)
999         sshb_e (:,:) = sshn_e(:,:)
1000         sshn_e (:,:) = ssha_e(:,:)
[4292]1001
1002         !                                             !* Sum over whole bt loop
1003         !                                             !  ----------------------
1004         za1 = wgtbtp1(jn)                                   
[6140]1005         IF( ln_dynadv_vec .OR. ln_linssh ) THEN    ! Sum velocities
[7753]1006            ua_b  (:,:) = ua_b  (:,:) + za1 * ua_e  (:,:) 
1007            va_b  (:,:) = va_b  (:,:) + za1 * va_e  (:,:) 
[6140]1008         ELSE                                              ! Sum transports
[7753]1009            ua_b  (:,:) = ua_b  (:,:) + za1 * ua_e  (:,:) * hu_e (:,:)
1010            va_b  (:,:) = va_b  (:,:) + za1 * va_e  (:,:) * hv_e (:,:)
[4292]1011         ENDIF
1012         !                                   ! Sum sea level
[7753]1013         ssha(:,:) = ssha(:,:) + za1 * ssha_e(:,:)
[358]1014         !                                                 ! ==================== !
1015      END DO                                               !        end loop      !
1016      !                                                    ! ==================== !
[1438]1017      ! -----------------------------------------------------------------------------
[1502]1018      ! Phase 3. update the general trend with the barotropic trend
[1438]1019      ! -----------------------------------------------------------------------------
[1502]1020      !
[4292]1021      ! Set advection velocity correction:
[8741]1022      IF (ln_bt_fw) THEN
1023         zwx(:,:) = un_adv(:,:)
1024         zwy(:,:) = vn_adv(:,:)
1025         IF( .NOT.( kt == nit000 .AND. neuler==0 ) ) THEN
1026            un_adv(:,:) = z1_2 * ( ub2_b(:,:) + zwx(:,:) ) 
1027            vn_adv(:,:) = z1_2 * ( vb2_b(:,:) + zwy(:,:) )
1028         END IF
1029         ! Save integrated transport for next computation
[7753]1030         ub2_b(:,:) = zwx(:,:)
1031         vb2_b(:,:) = zwy(:,:)
[4292]1032      ENDIF
1033      !
1034      ! Update barotropic trend:
[6140]1035      IF( ln_dynadv_vec .OR. ln_linssh ) THEN
[4292]1036         DO jk=1,jpkm1
[7753]1037            ua(:,:,jk) = ua(:,:,jk) + ( ua_b(:,:) - ub_b(:,:) ) * z1_2dt_b
1038            va(:,:,jk) = va(:,:,jk) + ( va_b(:,:) - vb_b(:,:) ) * z1_2dt_b
[4292]1039         END DO
1040      ELSE
[5930]1041         ! At this stage, ssha has been corrected: compute new depths at velocity points
1042         DO jj = 1, jpjm1
1043            DO ji = 1, jpim1      ! NO Vector Opt.
1044               zsshu_a(ji,jj) = z1_2 * umask(ji,jj,1)  * r1_e1e2u(ji,jj) &
1045                  &              * ( e1e2t(ji  ,jj) * ssha(ji  ,jj)    &
1046                  &              +   e1e2t(ji+1,jj) * ssha(ji+1,jj) )
1047               zsshv_a(ji,jj) = z1_2 * vmask(ji,jj,1)  * r1_e1e2v(ji,jj) &
1048                  &              * ( e1e2t(ji,jj  ) * ssha(ji,jj  )    &
1049                  &              +   e1e2t(ji,jj+1) * ssha(ji,jj+1) )
1050            END DO
1051         END DO
1052         CALL lbc_lnk_multi( zsshu_a, 'U', 1._wp, zsshv_a, 'V', 1._wp ) ! Boundary conditions
1053         !
[4292]1054         DO jk=1,jpkm1
[7753]1055            ua(:,:,jk) = ua(:,:,jk) + r1_hu_n(:,:) * ( ua_b(:,:) - ub_b(:,:) * hu_b(:,:) ) * z1_2dt_b
1056            va(:,:,jk) = va(:,:,jk) + r1_hv_n(:,:) * ( va_b(:,:) - vb_b(:,:) * hv_b(:,:) ) * z1_2dt_b
[4292]1057         END DO
1058         ! Save barotropic velocities not transport:
[7753]1059         ua_b(:,:) =  ua_b(:,:) / ( hu_0(:,:) + zsshu_a(:,:) + 1._wp - ssumask(:,:) )
1060         va_b(:,:) =  va_b(:,:) / ( hv_0(:,:) + zsshv_a(:,:) + 1._wp - ssvmask(:,:) )
[4292]1061      ENDIF
1062      !
1063      DO jk = 1, jpkm1
[7753]1064         ! Correct velocities:
[8741]1065         un(:,:,jk) = ( un(:,:,jk) + un_adv(:,:)*r1_hu_n(:,:) - un_b(:,:) ) * umask(:,:,jk)
1066         vn(:,:,jk) = ( vn(:,:,jk) + vn_adv(:,:)*r1_hv_n(:,:) - vn_b(:,:) ) * vmask(:,:,jk)
[7753]1067         !
[358]1068      END DO
[1502]1069      !
[8741]1070      CALL iom_put(  "ubar", un_adv(:,:)*r1_hu_n(:,:) )    ! barotropic i-current
1071      CALL iom_put(  "vbar", vn_adv(:,:)*r1_hv_n(:,:) )    ! barotropic i-current
[6140]1072      !
[4486]1073#if defined key_agrif
1074      ! Save time integrated fluxes during child grid integration
[5656]1075      ! (used to update coarse grid transports at next time step)
[4486]1076      !
[6140]1077      IF( .NOT.Agrif_Root() .AND. ln_bt_fw ) THEN
1078         IF( Agrif_NbStepint() == 0 ) THEN
[7753]1079            ub2_i_b(:,:) = 0._wp
1080            vb2_i_b(:,:) = 0._wp
[4486]1081         END IF
1082         !
1083         za1 = 1._wp / REAL(Agrif_rhot(), wp)
[7753]1084         ub2_i_b(:,:) = ub2_i_b(:,:) + za1 * ub2_b(:,:)
1085         vb2_i_b(:,:) = vb2_i_b(:,:) + za1 * vb2_b(:,:)
[4486]1086      ENDIF
1087#endif     
[1502]1088      !                                   !* write time-spliting arrays in the restart
[6140]1089      IF( lrst_oce .AND.ln_bt_fw )   CALL ts_rst( kt, 'WRITE' )
[508]1090      !
[6140]1091      CALL wrk_dealloc( jpi,jpj,   zsshp2_e, zhdiv )
1092      CALL wrk_dealloc( jpi,jpj,   zu_trd, zv_trd )
1093      CALL wrk_dealloc( jpi,jpj,   zwx, zwy, zssh_frc, zu_frc, zv_frc )
1094      CALL wrk_dealloc( jpi,jpj,   zhup2_e, zhvp2_e, zhust_e, zhvst_e )
1095      CALL wrk_dealloc( jpi,jpj,   zsshu_a, zsshv_a                                   )
1096      CALL wrk_dealloc( jpi,jpj,   zhf )
[7646]1097      IF( ln_wd ) CALL wrk_dealloc( jpi, jpj, zcpx, zcpy )
[1662]1098      !
[6140]1099      IF ( ln_diatmb ) THEN
1100         CALL iom_put( "baro_u" , un_b*umask(:,:,1)+zmdi*(1-umask(:,:,1 ) ) )  ! Barotropic  U Velocity
1101         CALL iom_put( "baro_v" , vn_b*vmask(:,:,1)+zmdi*(1-vmask(:,:,1 ) ) )  ! Barotropic  V Velocity
1102      ENDIF
[3294]1103      IF( nn_timing == 1 )  CALL timing_stop('dyn_spg_ts')
[2715]1104      !
[508]1105   END SUBROUTINE dyn_spg_ts
1106
[6140]1107
[4292]1108   SUBROUTINE ts_wgt( ll_av, ll_fw, jpit, zwgt1, zwgt2)
1109      !!---------------------------------------------------------------------
1110      !!                   ***  ROUTINE ts_wgt  ***
1111      !!
1112      !! ** Purpose : Set time-splitting weights for temporal averaging (or not)
1113      !!----------------------------------------------------------------------
1114      LOGICAL, INTENT(in) ::   ll_av      ! temporal averaging=.true.
1115      LOGICAL, INTENT(in) ::   ll_fw      ! forward time splitting =.true.
1116      INTEGER, INTENT(inout) :: jpit      ! cycle length   
1117      REAL(wp), DIMENSION(3*nn_baro), INTENT(inout) ::   zwgt1, & ! Primary weights
1118                                                         zwgt2    ! Secondary weights
1119     
1120      INTEGER ::  jic, jn, ji                      ! temporary integers
1121      REAL(wp) :: za1, za2
1122      !!----------------------------------------------------------------------
[508]1123
[4292]1124      zwgt1(:) = 0._wp
1125      zwgt2(:) = 0._wp
1126
1127      ! Set time index when averaged value is requested
1128      IF (ll_fw) THEN
1129         jic = nn_baro
1130      ELSE
1131         jic = 2 * nn_baro
1132      ENDIF
1133
1134      ! Set primary weights:
1135      IF (ll_av) THEN
1136           ! Define simple boxcar window for primary weights
1137           ! (width = nn_baro, centered around jic)     
1138         SELECT CASE ( nn_bt_flt )
1139              CASE( 0 )  ! No averaging
1140                 zwgt1(jic) = 1._wp
1141                 jpit = jic
1142
1143              CASE( 1 )  ! Boxcar, width = nn_baro
1144                 DO jn = 1, 3*nn_baro
1145                    za1 = ABS(float(jn-jic))/float(nn_baro) 
1146                    IF (za1 < 0.5_wp) THEN
1147                      zwgt1(jn) = 1._wp
1148                      jpit = jn
1149                    ENDIF
1150                 ENDDO
1151
1152              CASE( 2 )  ! Boxcar, width = 2 * nn_baro
1153                 DO jn = 1, 3*nn_baro
1154                    za1 = ABS(float(jn-jic))/float(nn_baro) 
1155                    IF (za1 < 1._wp) THEN
1156                      zwgt1(jn) = 1._wp
1157                      jpit = jn
1158                    ENDIF
1159                 ENDDO
1160              CASE DEFAULT   ;   CALL ctl_stop( 'unrecognised value for nn_bt_flt' )
1161         END SELECT
1162
1163      ELSE ! No time averaging
1164         zwgt1(jic) = 1._wp
1165         jpit = jic
1166      ENDIF
1167   
1168      ! Set secondary weights
1169      DO jn = 1, jpit
1170        DO ji = jn, jpit
1171             zwgt2(jn) = zwgt2(jn) + zwgt1(ji)
1172        END DO
1173      END DO
1174
1175      ! Normalize weigths:
1176      za1 = 1._wp / SUM(zwgt1(1:jpit))
1177      za2 = 1._wp / SUM(zwgt2(1:jpit))
1178      DO jn = 1, jpit
1179        zwgt1(jn) = zwgt1(jn) * za1
1180        zwgt2(jn) = zwgt2(jn) * za2
1181      END DO
1182      !
1183   END SUBROUTINE ts_wgt
1184
[6140]1185
[508]1186   SUBROUTINE ts_rst( kt, cdrw )
1187      !!---------------------------------------------------------------------
1188      !!                   ***  ROUTINE ts_rst  ***
1189      !!
1190      !! ** Purpose : Read or write time-splitting arrays in restart file
1191      !!----------------------------------------------------------------------
1192      INTEGER         , INTENT(in) ::   kt         ! ocean time-step
1193      CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
1194      !
1195      !!----------------------------------------------------------------------
1196      !
1197      IF( TRIM(cdrw) == 'READ' ) THEN
[4292]1198         CALL iom_get( numror, jpdom_autoglo, 'ub2_b'  , ub2_b  (:,:) )   
1199         CALL iom_get( numror, jpdom_autoglo, 'vb2_b'  , vb2_b  (:,:) ) 
[4370]1200         IF( .NOT.ln_bt_av ) THEN
[4292]1201            CALL iom_get( numror, jpdom_autoglo, 'sshbb_e'  , sshbb_e(:,:) )   
1202            CALL iom_get( numror, jpdom_autoglo, 'ubb_e'    ,   ubb_e(:,:) )   
1203            CALL iom_get( numror, jpdom_autoglo, 'vbb_e'    ,   vbb_e(:,:) )
1204            CALL iom_get( numror, jpdom_autoglo, 'sshb_e'   ,  sshb_e(:,:) ) 
1205            CALL iom_get( numror, jpdom_autoglo, 'ub_e'     ,    ub_e(:,:) )   
1206            CALL iom_get( numror, jpdom_autoglo, 'vb_e'     ,    vb_e(:,:) )
[508]1207         ENDIF
[4486]1208#if defined key_agrif
1209         ! Read time integrated fluxes
1210         IF ( .NOT.Agrif_Root() ) THEN
1211            CALL iom_get( numror, jpdom_autoglo, 'ub2_i_b'  , ub2_i_b(:,:) )   
1212            CALL iom_get( numror, jpdom_autoglo, 'vb2_i_b'  , vb2_i_b(:,:) )
1213         ENDIF
1214#endif
[4292]1215      !
1216      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN
1217         CALL iom_rstput( kt, nitrst, numrow, 'ub2_b'   , ub2_b  (:,:) )
1218         CALL iom_rstput( kt, nitrst, numrow, 'vb2_b'   , vb2_b  (:,:) )
1219         !
1220         IF (.NOT.ln_bt_av) THEN
1221            CALL iom_rstput( kt, nitrst, numrow, 'sshbb_e'  , sshbb_e(:,:) ) 
1222            CALL iom_rstput( kt, nitrst, numrow, 'ubb_e'    ,   ubb_e(:,:) )
1223            CALL iom_rstput( kt, nitrst, numrow, 'vbb_e'    ,   vbb_e(:,:) )
1224            CALL iom_rstput( kt, nitrst, numrow, 'sshb_e'   ,  sshb_e(:,:) )
1225            CALL iom_rstput( kt, nitrst, numrow, 'ub_e'     ,    ub_e(:,:) )
1226            CALL iom_rstput( kt, nitrst, numrow, 'vb_e'     ,    vb_e(:,:) )
1227         ENDIF
[4486]1228#if defined key_agrif
1229         ! Save time integrated fluxes
1230         IF ( .NOT.Agrif_Root() ) THEN
1231            CALL iom_rstput( kt, nitrst, numrow, 'ub2_i_b'  , ub2_i_b(:,:) )
1232            CALL iom_rstput( kt, nitrst, numrow, 'vb2_i_b'  , vb2_i_b(:,:) )
1233         ENDIF
1234#endif
[4292]1235      ENDIF
1236      !
1237   END SUBROUTINE ts_rst
[2528]1238
[6140]1239
1240   SUBROUTINE dyn_spg_ts_init
[4292]1241      !!---------------------------------------------------------------------
1242      !!                   ***  ROUTINE dyn_spg_ts_init  ***
1243      !!
1244      !! ** Purpose : Set time splitting options
1245      !!----------------------------------------------------------------------
[6140]1246      INTEGER  ::   ji ,jj              ! dummy loop indices
1247      REAL(wp) ::   zxr2, zyr2, zcmax   ! local scalar
1248      REAL(wp), POINTER, DIMENSION(:,:) ::   zcu
[4292]1249      !!----------------------------------------------------------------------
[4370]1250      !
[5930]1251      ! Max courant number for ext. grav. waves
[4370]1252      !
[6140]1253      CALL wrk_alloc( jpi,jpj,   zcu )
[4292]1254      !
[5930]1255      DO jj = 1, jpj
1256         DO ji =1, jpi
1257            zxr2 = r1_e1t(ji,jj) * r1_e1t(ji,jj)
1258            zyr2 = r1_e2t(ji,jj) * r1_e2t(ji,jj)
[7646]1259            zcu(ji,jj) = SQRT( grav * MAX(ht_0(ji,jj),0._wp) * (zxr2 + zyr2) )
[4370]1260         END DO
[5930]1261      END DO
1262      !
[5836]1263      zcmax = MAXVAL( zcu(:,:) )
[4292]1264      IF( lk_mpp )   CALL mpp_max( zcmax )
[2528]1265
[4370]1266      ! Estimate number of iterations to satisfy a max courant number= rn_bt_cmax
[6140]1267      IF( ln_bt_auto )   nn_baro = CEILING( rdt / rn_bt_cmax * zcmax)
[4292]1268     
[5836]1269      rdtbt = rdt / REAL( nn_baro , wp )
[4292]1270      zcmax = zcmax * rdtbt
1271                     ! Print results
1272      IF(lwp) WRITE(numout,*)
1273      IF(lwp) WRITE(numout,*) 'dyn_spg_ts : split-explicit free surface'
1274      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
[5930]1275      IF( ln_bt_auto ) THEN
1276         IF(lwp) WRITE(numout,*) '     ln_ts_auto=.true. Automatically set nn_baro '
[4370]1277         IF(lwp) WRITE(numout,*) '     Max. courant number allowed: ', rn_bt_cmax
[4292]1278      ELSE
[5930]1279         IF(lwp) WRITE(numout,*) '     ln_ts_auto=.false.: Use nn_baro in namelist '
[358]1280      ENDIF
[4292]1281
1282      IF(ln_bt_av) THEN
[4370]1283         IF(lwp) WRITE(numout,*) '     ln_bt_av=.true.  => Time averaging over nn_baro time steps is on '
[4292]1284      ELSE
[4370]1285         IF(lwp) WRITE(numout,*) '     ln_bt_av=.false. => No time averaging of barotropic variables '
[4292]1286      ENDIF
[508]1287      !
[4292]1288      !
1289      IF(ln_bt_fw) THEN
[4370]1290         IF(lwp) WRITE(numout,*) '     ln_bt_fw=.true.  => Forward integration of barotropic variables '
[4292]1291      ELSE
[4370]1292         IF(lwp) WRITE(numout,*) '     ln_bt_fw =.false.=> Centred integration of barotropic variables '
[4292]1293      ENDIF
1294      !
[4486]1295#if defined key_agrif
1296      ! Restrict the use of Agrif to the forward case only
[8762]1297      IF( .NOT.ln_bt_fw .AND. .NOT.Agrif_Root() )   CALL ctl_stop( 'AGRIF not implemented if ln_bt_fw=.FALSE.' )
[4486]1298#endif
1299      !
[4370]1300      IF(lwp) WRITE(numout,*)    '     Time filter choice, nn_bt_flt: ', nn_bt_flt
[4292]1301      SELECT CASE ( nn_bt_flt )
[6140]1302         CASE( 0 )      ;   IF(lwp) WRITE(numout,*) '           Dirac'
1303         CASE( 1 )      ;   IF(lwp) WRITE(numout,*) '           Boxcar: width = nn_baro'
1304         CASE( 2 )      ;   IF(lwp) WRITE(numout,*) '           Boxcar: width = 2*nn_baro' 
1305         CASE DEFAULT   ;   CALL ctl_stop( 'unrecognised value for nn_bt_flt: should 0,1,2' )
[4292]1306      END SELECT
1307      !
[4370]1308      IF(lwp) WRITE(numout,*) ' '
1309      IF(lwp) WRITE(numout,*) '     nn_baro = ', nn_baro
1310      IF(lwp) WRITE(numout,*) '     Barotropic time step [s] is :', rdtbt
1311      IF(lwp) WRITE(numout,*) '     Maximum Courant number is   :', zcmax
1312      !
[6140]1313      IF( .NOT.ln_bt_av .AND. .NOT.ln_bt_fw ) THEN
[4292]1314         CALL ctl_stop( 'dynspg_ts ERROR: No time averaging => only forward integration is possible' )
1315      ENDIF
[6140]1316      IF( zcmax>0.9_wp ) THEN
[4292]1317         CALL ctl_stop( 'dynspg_ts ERROR: Maximum Courant number is greater than 0.9: Inc. nn_baro !' )         
1318      ENDIF
1319      !
[6140]1320      CALL wrk_dealloc( jpi,jpj,   zcu )
[4292]1321      !
1322   END SUBROUTINE dyn_spg_ts_init
[508]1323
[358]1324   !!======================================================================
1325END MODULE dynspg_ts
Note: See TracBrowser for help on using the repository browser.