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/2016/dev_r6393_NOC_WAD/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/2016/dev_r6393_NOC_WAD/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_ts.F90 @ 7016

Last change on this file since 7016 was 7016, checked in by acc, 8 years ago

Branch 2016/dev_r6393_NOC_WAD. Reinstated and fixed limiter in barotropic loop to resolve conservation issues with wetting and drying. Now maintains constant salinity in WAD_TEST_CASES 1,2 and 3.

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