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

source: branches/2016/dev_merge_2016/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_ts.F90 @ 7412

Last change on this file since 7412 was 7412, checked in by lovato, 8 years ago

Merge dev_NOC_CMCC_merge_2016 into branch

  • Property svn:keywords set to Id
File size: 62.8 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_oce   , ONLY: ln_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) ) >                &
378                     &    MAX( -bathy(ji,jj)               , -bathy(ji+1,jj) ) .AND.            &
379                     &    MAX(   sshn(ji,jj) + bathy(ji,jj),   sshn(ji+1,jj) + bathy(ji+1,jj) ) &
380                     &                                                         > rn_wdmin1 + rn_wdmin2
381                ll_tmp2 = MAX(   sshn(ji,jj)               ,   sshn(ji+1,jj) ) >                &
382                     &    MAX( -bathy(ji,jj)               , -bathy(ji+1,jj) ) + rn_wdmin1 + rn_wdmin2
383   
384                IF(ll_tmp1) THEN
385                  zcpx(ji,jj) = 1.0_wp
386                ELSE IF(ll_tmp2) THEN
387                  ! no worries about  sshn(ji+1,jj) -  sshn(ji  ,jj) = 0, it won't happen ! here
388                  zcpx(ji,jj) = ABS( (sshn(ji+1,jj) + bathy(ji+1,jj) - sshn(ji,jj) - bathy(ji,jj)) &
389                              &    / (sshn(ji+1,jj) -  sshn(ji  ,jj)) )
390                ELSE
391                  zcpx(ji,jj) = 0._wp
392                END IF
393         
394                ll_tmp1 = MIN(   sshn(ji,jj)               ,   sshn(ji,jj+1) ) >                &
395                     &    MAX( -bathy(ji,jj)               , -bathy(ji,jj+1) ) .AND.            &
396                     &    MAX(   sshn(ji,jj) + bathy(ji,jj),   sshn(ji,jj+1) + bathy(ji,jj+1) ) &
397                     &                                                         > rn_wdmin1 + rn_wdmin2
398                ll_tmp2 = MAX(   sshn(ji,jj)               ,   sshn(ji,jj+1) ) >                &
399                     &    MAX( -bathy(ji,jj)               , -bathy(ji,jj+1) ) + rn_wdmin1 + rn_wdmin2
400   
401                IF(ll_tmp1) THEN
402                  zcpy(ji,jj) = 1.0_wp
403                ELSE IF(ll_tmp2) THEN
404                  ! no worries about  sshn(ji,jj+1) -  sshn(ji,jj  ) = 0, it won't happen ! here
405                  zcpy(ji,jj) = ABS( (sshn(ji,jj+1) + bathy(ji,jj+1) - sshn(ji,jj) - bathy(ji,jj)) &
406                              &    / (sshn(ji,jj+1) -  sshn(ji,jj  )) )
407                ELSE
408                  zcpy(ji,jj) = 0._wp
409                END IF
410              END DO
411           END DO
412 
413           DO jj = 2, jpjm1
414              DO ji = 2, jpim1
415                 zu_trd(ji,jj) = zu_trd(ji,jj) - grav * ( sshn(ji+1,jj  ) - sshn(ji  ,jj ) )   &
416                        &                        * r1_e1u(ji,jj) * zcpx(ji,jj)
417                 zv_trd(ji,jj) = zv_trd(ji,jj) - grav * ( sshn(ji  ,jj+1) - sshn(ji  ,jj ) )   &
418                        &                        * r1_e2v(ji,jj) * zcpy(ji,jj)
419              END DO
420           END DO
421
422         ELSE
423
424           DO jj = 2, jpjm1
425              DO ji = fs_2, fs_jpim1   ! vector opt.
426                 zu_trd(ji,jj) = zu_trd(ji,jj) - grav * (  sshn(ji+1,jj  ) - sshn(ji  ,jj  )  ) * r1_e1u(ji,jj)
427                 zv_trd(ji,jj) = zv_trd(ji,jj) - grav * (  sshn(ji  ,jj+1) - sshn(ji  ,jj  )  ) * r1_e2v(ji,jj) 
428              END DO
429           END DO
430        ENDIF
431
432      ENDIF
433
434      DO jj = 2, jpjm1                          ! Remove coriolis term (and possibly spg) from barotropic trend
435         DO ji = fs_2, fs_jpim1
436             zu_frc(ji,jj) = zu_frc(ji,jj) - zu_trd(ji,jj) * ssumask(ji,jj)
437             zv_frc(ji,jj) = zv_frc(ji,jj) - zv_trd(ji,jj) * ssvmask(ji,jj)
438          END DO
439      END DO 
440      !
441      !                 ! Add bottom stress contribution from baroclinic velocities:     
442      IF (ln_bt_fw) THEN
443         DO jj = 2, jpjm1                         
444            DO ji = fs_2, fs_jpim1   ! vector opt.
445               ikbu = mbku(ji,jj)       
446               ikbv = mbkv(ji,jj)   
447               zwx(ji,jj) = un(ji,jj,ikbu) - un_b(ji,jj) ! NOW bottom baroclinic velocities
448               zwy(ji,jj) = vn(ji,jj,ikbv) - vn_b(ji,jj)
449            END DO
450         END DO
451      ELSE
452         DO jj = 2, jpjm1
453            DO ji = fs_2, fs_jpim1   ! vector opt.
454               ikbu = mbku(ji,jj)       
455               ikbv = mbkv(ji,jj)   
456               zwx(ji,jj) = ub(ji,jj,ikbu) - ub_b(ji,jj) ! BEFORE bottom baroclinic velocities
457               zwy(ji,jj) = vb(ji,jj,ikbv) - vb_b(ji,jj)
458            END DO
459         END DO
460      ENDIF
461      !
462      ! Note that the "unclipped" bottom friction parameter is used even with explicit drag
463      IF( ln_wd ) THEN
464        zu_frc(:,:) = zu_frc(:,:) + MAX(r1_hu_n(:,:) * bfrua(:,:),-1._wp / rdtbt) * zwx(:,:)
465        zv_frc(:,:) = zv_frc(:,:) + MAX(r1_hv_n(:,:) * bfrva(:,:),-1._wp / rdtbt) * zwy(:,:)
466      ELSE
467        zu_frc(:,:) = zu_frc(:,:) + r1_hu_n(:,:) * bfrua(:,:) * zwx(:,:)
468        zv_frc(:,:) = zv_frc(:,:) + r1_hv_n(:,:) * bfrva(:,:) * zwy(:,:)
469      END IF
470      !
471      !                                         ! Add top 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               iktu = miku(ji,jj)
476               iktv = mikv(ji,jj)
477               zwx(ji,jj) = un(ji,jj,iktu) - un_b(ji,jj) ! NOW top baroclinic velocities
478               zwy(ji,jj) = vn(ji,jj,iktv) - vn_b(ji,jj)
479            END DO
480         END DO
481      ELSE
482         DO jj = 2, jpjm1
483            DO ji = fs_2, fs_jpim1   ! vector opt.
484               iktu = miku(ji,jj)
485               iktv = mikv(ji,jj)
486               zwx(ji,jj) = ub(ji,jj,iktu) - ub_b(ji,jj) ! BEFORE top baroclinic velocities
487               zwy(ji,jj) = vb(ji,jj,iktv) - vb_b(ji,jj)
488            END DO
489         END DO
490      ENDIF
491      !
492      ! Note that the "unclipped" top friction parameter is used even with explicit drag
493      zu_frc(:,:) = zu_frc(:,:) + r1_hu_n(:,:) * tfrua(:,:) * zwx(:,:)
494      zv_frc(:,:) = zv_frc(:,:) + r1_hv_n(:,:) * tfrva(:,:) * zwy(:,:)
495      !       
496      IF (ln_bt_fw) THEN                        ! Add wind forcing
497         zu_frc(:,:) =  zu_frc(:,:) + zraur * utau(:,:) * r1_hu_n(:,:)
498         zv_frc(:,:) =  zv_frc(:,:) + zraur * vtau(:,:) * r1_hv_n(:,:)
499      ELSE
500         zu_frc(:,:) =  zu_frc(:,:) + zraur * z1_2 * ( utau_b(:,:) + utau(:,:) ) * r1_hu_n(:,:)
501         zv_frc(:,:) =  zv_frc(:,:) + zraur * z1_2 * ( vtau_b(:,:) + vtau(:,:) ) * r1_hv_n(:,:)
502      ENDIF 
503      !
504      IF ( ln_apr_dyn ) THEN                    ! Add atm pressure forcing
505         IF (ln_bt_fw) THEN
506            DO jj = 2, jpjm1             
507               DO ji = fs_2, fs_jpim1   ! vector opt.
508                  zu_spg =  grav * (  ssh_ib (ji+1,jj  ) - ssh_ib (ji,jj) ) * r1_e1u(ji,jj)
509                  zv_spg =  grav * (  ssh_ib (ji  ,jj+1) - ssh_ib (ji,jj) ) * r1_e2v(ji,jj)
510                  zu_frc(ji,jj) = zu_frc(ji,jj) + zu_spg
511                  zv_frc(ji,jj) = zv_frc(ji,jj) + zv_spg
512               END DO
513            END DO
514         ELSE
515            DO jj = 2, jpjm1             
516               DO ji = fs_2, fs_jpim1   ! vector opt.
517                  zu_spg =  grav * z1_2 * (  ssh_ib (ji+1,jj  ) - ssh_ib (ji,jj)    &
518                      &                    + ssh_ibb(ji+1,jj  ) - ssh_ibb(ji,jj)  ) * r1_e1u(ji,jj)
519                  zv_spg =  grav * z1_2 * (  ssh_ib (ji  ,jj+1) - ssh_ib (ji,jj)    &
520                      &                    + ssh_ibb(ji  ,jj+1) - ssh_ibb(ji,jj)  ) * r1_e2v(ji,jj)
521                  zu_frc(ji,jj) = zu_frc(ji,jj) + zu_spg
522                  zv_frc(ji,jj) = zv_frc(ji,jj) + zv_spg
523               END DO
524            END DO
525         ENDIF
526      ENDIF
527      !                                   !* Right-Hand-Side of the barotropic ssh equation
528      !                                   ! -----------------------------------------------
529      !                                         ! Surface net water flux and rivers
530      IF (ln_bt_fw) THEN
531         zssh_frc(:,:) = zraur * ( emp(:,:) - rnf(:,:) + fwfisf(:,:) )
532      ELSE
533         zssh_frc(:,:) = zraur * z1_2 * (  emp(:,:) + emp_b(:,:) - rnf(:,:) - rnf_b(:,:)   &
534                &                        + fwfisf(:,:) + fwfisf_b(:,:)                     )
535      ENDIF
536#if defined key_asminc
537      !                                         ! Include the IAU weighted SSH increment
538      IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN
539         zssh_frc(:,:) = zssh_frc(:,:) - ssh_iau(:,:)
540      ENDIF
541#endif
542      !                                   !* Fill boundary data arrays for AGRIF
543      !                                   ! ------------------------------------
544#if defined key_agrif
545         IF( .NOT.Agrif_Root() ) CALL agrif_dta_ts( kt )
546#endif
547      !
548      ! -----------------------------------------------------------------------
549      !  Phase 2 : Integration of the barotropic equations
550      ! -----------------------------------------------------------------------
551      !
552      !                                             ! ==================== !
553      !                                             !    Initialisations   !
554      !                                             ! ==================== ! 
555      ! Initialize barotropic variables:     
556      IF( ll_init )THEN
557         sshbb_e(:,:) = 0._wp
558         ubb_e  (:,:) = 0._wp
559         vbb_e  (:,:) = 0._wp
560         sshb_e (:,:) = 0._wp
561         ub_e   (:,:) = 0._wp
562         vb_e   (:,:) = 0._wp
563      ENDIF
564
565      !
566      IF (ln_bt_fw) THEN                  ! FORWARD integration: start from NOW fields                   
567         sshn_e(:,:) =    sshn(:,:)           
568         un_e  (:,:) =    un_b(:,:)           
569         vn_e  (:,:) =    vn_b(:,:)
570         !
571         hu_e  (:,:) =    hu_n(:,:)       
572         hv_e  (:,:) =    hv_n(:,:) 
573         hur_e (:,:) = r1_hu_n(:,:)   
574         hvr_e (:,:) = r1_hv_n(:,:)
575      ELSE                                ! CENTRED integration: start from BEFORE fields
576         sshn_e(:,:) =    sshb(:,:)
577         un_e  (:,:) =    ub_b(:,:)         
578         vn_e  (:,:) =    vb_b(:,:)
579         !
580         hu_e  (:,:) =    hu_b(:,:)       
581         hv_e  (:,:) =    hv_b(:,:) 
582         hur_e (:,:) = r1_hu_b(:,:)   
583         hvr_e (:,:) = r1_hv_b(:,:)
584      ENDIF
585      !
586      !
587      !
588      ! Initialize sums:
589      ua_b  (:,:) = 0._wp       ! After barotropic velocities (or transport if flux form)         
590      va_b  (:,:) = 0._wp
591      ssha  (:,:) = 0._wp       ! Sum for after averaged sea level
592      un_adv(:,:) = 0._wp       ! Sum for now transport issued from ts loop
593      vn_adv(:,:) = 0._wp
594      !                                             ! ==================== !
595      DO jn = 1, icycle                             !  sub-time-step loop  !
596         !                                          ! ==================== !
597         !                                                !* Update the forcing (BDY and tides)
598         !                                                !  ------------------
599         ! Update only tidal forcing at open boundaries
600         IF( ln_bdy      .AND. ln_tide )   CALL bdy_dta_tides( kt, kit=jn, time_offset= noffset+1 )
601         IF( ln_tide_pot .AND. ln_tide )   CALL upd_tide     ( kt, kit=jn, time_offset= noffset   )
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         ! Duplicate sea level across open boundaries (this is only cosmetic if linssh=T)
694         IF( ln_bdy )   CALL bdy_ssh( ssha_e )
695
696#if defined key_agrif
697         IF( .NOT.Agrif_Root() )   CALL agrif_ssh_ts( jn )
698#endif
699         
700         ! Sea Surface Height at u-,v-points (vvl case only)
701         IF( .NOT.ln_linssh ) THEN                               
702            DO jj = 2, jpjm1
703               DO ji = 2, jpim1      ! NO Vector Opt.
704                  zsshu_a(ji,jj) = z1_2 * ssumask(ji,jj) * r1_e1e2u(ji,jj)    &
705                     &              * ( e1e2t(ji  ,jj  )  * ssha_e(ji  ,jj  ) &
706                     &              +   e1e2t(ji+1,jj  )  * ssha_e(ji+1,jj  ) )
707                  zsshv_a(ji,jj) = z1_2 * ssvmask(ji,jj) * r1_e1e2v(ji,jj)    &
708                     &              * ( e1e2t(ji  ,jj  )  * ssha_e(ji  ,jj  ) &
709                     &              +   e1e2t(ji  ,jj+1)  * ssha_e(ji  ,jj+1) )
710               END DO
711            END DO
712            CALL lbc_lnk_multi( zsshu_a, 'U', 1._wp, zsshv_a, 'V', 1._wp )
713         ENDIF   
714         !                                 
715         ! Half-step back interpolation of SSH for surface pressure computation:
716         !----------------------------------------------------------------------
717         IF ((jn==1).AND.ll_init) THEN
718           za0=1._wp                        ! Forward-backward
719           za1=0._wp                           
720           za2=0._wp
721           za3=0._wp
722         ELSEIF ((jn==2).AND.ll_init) THEN  ! AB2-AM3 Coefficients; bet=0 ; gam=-1/6 ; eps=1/12
723           za0= 1.0833333333333_wp          ! za0 = 1-gam-eps
724           za1=-0.1666666666666_wp          ! za1 = gam
725           za2= 0.0833333333333_wp          ! za2 = eps
726           za3= 0._wp             
727         ELSE                               ! AB3-AM4 Coefficients; bet=0.281105 ; eps=0.013 ; gam=0.0880
728           za0=0.614_wp                     ! za0 = 1/2 +   gam + 2*eps   
729           za1=0.285_wp                     ! za1 = 1/2 - 2*gam - 3*eps
730           za2=0.088_wp                     ! za2 = gam
731           za3=0.013_wp                     ! za3 = eps
732         ENDIF
733         !
734         zsshp2_e(:,:) = za0 *  ssha_e(:,:) + za1 *  sshn_e (:,:) &
735          &            + za2 *  sshb_e(:,:) + za3 *  sshbb_e(:,:)
736
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) ) >                &
741                     &    MAX(   -bathy(ji,jj)               ,   -bathy(ji+1,jj) ) .AND.            &
742                     &    MAX( zsshp2_e(ji,jj) + bathy(ji,jj), zsshp2_e(ji+1,jj) + bathy(ji+1,jj) ) &
743                     &                                                             > rn_wdmin1 + rn_wdmin2
744                ll_tmp2 = MAX( zsshp2_e(ji,jj)               , zsshp2_e(ji+1,jj) ) >                &
745                     &    MAX(   -bathy(ji,jj)               ,   -bathy(ji+1,jj) ) + rn_wdmin1 + rn_wdmin2
746   
747                IF(ll_tmp1) THEN
748                  zcpx(ji,jj) = 1.0_wp
749                ELSE IF(ll_tmp2) THEN
750                  ! no worries about  zsshp2_e(ji+1,jj) - zsshp2_e(ji  ,jj) = 0, it won't happen ! here
751                  zcpx(ji,jj) = ABS( (zsshp2_e(ji+1,jj) +    bathy(ji+1,jj) - zsshp2_e(ji,jj) - bathy(ji,jj)) &
752                              &    / (zsshp2_e(ji+1,jj) - zsshp2_e(ji  ,jj)) )
753                ELSE
754                  zcpx(ji,jj) = 0._wp
755                END IF
756         
757                ll_tmp1 = MIN( zsshp2_e(ji,jj)               , zsshp2_e(ji,jj+1) ) >                &
758                     &    MAX(   -bathy(ji,jj)               ,   -bathy(ji,jj+1) ) .AND.            &
759                     &    MAX( zsshp2_e(ji,jj) + bathy(ji,jj), zsshp2_e(ji,jj+1) + bathy(ji,jj+1) ) &
760                     &                                                             > rn_wdmin1 + rn_wdmin2
761                ll_tmp2 = MAX( zsshp2_e(ji,jj)               , zsshp2_e(ji,jj+1) ) >                &
762                     &    MAX(   -bathy(ji,jj)               ,   -bathy(ji,jj+1) ) + rn_wdmin1 + rn_wdmin2
763   
764                IF(ll_tmp1) THEN
765                  zcpy(ji,jj) = 1.0_wp
766                ELSE IF(ll_tmp2) THEN
767                  ! no worries about  zsshp2_e(ji,jj+1) - zsshp2_e(ji,jj  ) = 0, it won't happen ! here
768                  zcpy(ji,jj) = ABS( (zsshp2_e(ji,jj+1) +    bathy(ji,jj+1) - zsshp2_e(ji,jj) - bathy(ji,jj)) &
769                              &    / (zsshp2_e(ji,jj+1) - zsshp2_e(ji,jj  )) )
770                ELSE
771                  zcpy(ji,jj) = 0._wp
772                END IF
773              END DO
774           END DO
775         END IF
776         !
777         ! Compute associated depths at U and V points:
778         IF( .NOT.ln_linssh  .AND. .NOT.ln_dynadv_vec ) THEN   !* Vector form
779            !                                       
780            DO jj = 2, jpjm1                           
781               DO ji = 2, jpim1
782                  zx1 = z1_2 * ssumask(ji  ,jj) *  r1_e1e2u(ji  ,jj)    &
783                     &      * ( e1e2t(ji  ,jj  ) * zsshp2_e(ji  ,jj)    &
784                     &      +   e1e2t(ji+1,jj  ) * zsshp2_e(ji+1,jj  ) )
785                  zy1 = z1_2 * ssvmask(ji  ,jj) *  r1_e1e2v(ji  ,jj  )  &
786                     &       * ( e1e2t(ji ,jj  ) * zsshp2_e(ji  ,jj  )  &
787                     &       +   e1e2t(ji ,jj+1) * zsshp2_e(ji  ,jj+1) )
788                  zhust_e(ji,jj) = hu_0(ji,jj) + zx1 
789                  zhvst_e(ji,jj) = hv_0(ji,jj) + zy1
790               END DO
791            END DO
792
793         ENDIF
794         !
795         ! Add Coriolis trend:
796         ! zwz array below or triads normally depend on sea level with ln_linssh=F and should be updated
797         ! at each time step. We however keep them constant here for optimization.
798         ! Recall that zwx and zwy arrays hold fluxes at this stage:
799         ! zwx(:,:) = e2u(:,:) * ua_e(:,:) * zhup2_e(:,:)   ! fluxes at jn+0.5
800         ! zwy(:,:) = e1v(:,:) * va_e(:,:) * zhvp2_e(:,:)
801         !
802         IF( ln_dynvor_ene .OR. ln_dynvor_mix ) THEN     !==  energy conserving or mixed scheme  ==!
803            DO jj = 2, jpjm1
804               DO ji = fs_2, fs_jpim1   ! vector opt.
805                  zy1 = ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) ) * r1_e1u(ji,jj)
806                  zy2 = ( zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) * r1_e1u(ji,jj)
807                  zx1 = ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) ) * r1_e2v(ji,jj)
808                  zx2 = ( zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) * r1_e2v(ji,jj)
809                  zu_trd(ji,jj) = z1_4 * ( zwz(ji  ,jj-1) * zy1 + zwz(ji,jj) * zy2 )
810                  zv_trd(ji,jj) =-z1_4 * ( zwz(ji-1,jj  ) * zx1 + zwz(ji,jj) * zx2 )
811               END DO
812            END DO
813            !
814         ELSEIF ( ln_dynvor_ens ) THEN                   !==  enstrophy conserving scheme  ==!
815            DO jj = 2, jpjm1
816               DO ji = fs_2, fs_jpim1   ! vector opt.
817                  zy1 =   z1_8 * ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) &
818                   &             + zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) * r1_e1u(ji,jj)
819                  zx1 = - z1_8 * ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) &
820                   &             + zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) * r1_e2v(ji,jj)
821                  zu_trd(ji,jj)  = zy1 * ( zwz(ji  ,jj-1) + zwz(ji,jj) )
822                  zv_trd(ji,jj)  = zx1 * ( zwz(ji-1,jj  ) + zwz(ji,jj) )
823               END DO
824            END DO
825            !
826         ELSEIF ( ln_dynvor_een ) THEN                   !==  energy and enstrophy conserving scheme  ==!
827            DO jj = 2, jpjm1
828               DO ji = fs_2, fs_jpim1   ! vector opt.
829                  zu_trd(ji,jj) = + z1_12 * r1_e1u(ji,jj) * (  ftne(ji,jj  ) * zwy(ji  ,jj  ) &
830                     &                                       + ftnw(ji+1,jj) * zwy(ji+1,jj  ) &
831                     &                                       + ftse(ji,jj  ) * zwy(ji  ,jj-1) & 
832                     &                                       + ftsw(ji+1,jj) * zwy(ji+1,jj-1) )
833                  zv_trd(ji,jj) = - z1_12 * r1_e2v(ji,jj) * (  ftsw(ji,jj+1) * zwx(ji-1,jj+1) & 
834                     &                                       + ftse(ji,jj+1) * zwx(ji  ,jj+1) &
835                     &                                       + ftnw(ji,jj  ) * zwx(ji-1,jj  ) & 
836                     &                                       + ftne(ji,jj  ) * zwx(ji  ,jj  ) )
837               END DO
838            END DO
839            !
840         ENDIF
841         !
842         ! Add tidal astronomical forcing if defined
843         IF ( ln_tide.AND.ln_tide_pot ) THEN
844            DO jj = 2, jpjm1
845               DO ji = fs_2, fs_jpim1   ! vector opt.
846                  zu_spg = grav * ( pot_astro(ji+1,jj) - pot_astro(ji,jj) ) * r1_e1u(ji,jj)
847                  zv_spg = grav * ( pot_astro(ji,jj+1) - pot_astro(ji,jj) ) * r1_e2v(ji,jj)
848                  zu_trd(ji,jj) = zu_trd(ji,jj) + zu_spg
849                  zv_trd(ji,jj) = zv_trd(ji,jj) + zv_spg
850               END DO
851            END DO
852         ENDIF
853         !
854         ! Add bottom stresses:
855         zu_trd(:,:) = zu_trd(:,:) + bfrua(:,:) * un_e(:,:) * hur_e(:,:)
856         zv_trd(:,:) = zv_trd(:,:) + bfrva(:,:) * vn_e(:,:) * hvr_e(:,:)
857         !
858         ! Add top stresses:
859         zu_trd(:,:) = zu_trd(:,:) + tfrua(:,:) * un_e(:,:) * hur_e(:,:)
860         zv_trd(:,:) = zv_trd(:,:) + tfrva(:,:) * vn_e(:,:) * hvr_e(:,:)
861         !
862         ! Surface pressure trend:
863
864         IF( ln_wd ) THEN
865           DO jj = 2, jpjm1
866              DO ji = 2, jpim1 
867                 ! Add surface pressure gradient
868                 zu_spg = - grav * ( zsshp2_e(ji+1,jj) - zsshp2_e(ji,jj) ) * r1_e1u(ji,jj)
869                 zv_spg = - grav * ( zsshp2_e(ji,jj+1) - zsshp2_e(ji,jj) ) * r1_e2v(ji,jj)
870                 zwx(ji,jj) = zu_spg * zcpx(ji,jj) * wdmask(ji,jj) * wdmask(ji+1, jj) 
871                 zwy(ji,jj) = zv_spg * zcpy(ji,jj) * wdmask(ji,jj) * wdmask(ji, jj+1)
872              END DO
873           END DO
874         ELSE
875           DO jj = 2, jpjm1
876              DO ji = fs_2, fs_jpim1   ! vector opt.
877                 ! Add surface pressure gradient
878                 zu_spg = - grav * ( zsshp2_e(ji+1,jj) - zsshp2_e(ji,jj) ) * r1_e1u(ji,jj)
879                 zv_spg = - grav * ( zsshp2_e(ji,jj+1) - zsshp2_e(ji,jj) ) * r1_e2v(ji,jj)
880                 zwx(ji,jj) = zu_spg
881                 zwy(ji,jj) = zv_spg
882              END DO
883           END DO
884         END IF
885
886         !
887         ! Set next velocities:
888         IF( ln_dynadv_vec .OR. ln_linssh ) THEN   !* Vector form
889            DO jj = 2, jpjm1
890               DO ji = fs_2, fs_jpim1   ! vector opt.
891                  ua_e(ji,jj) = (                                 un_e(ji,jj)   & 
892                            &     + rdtbt * (                      zwx(ji,jj)   &
893                            &                                 + zu_trd(ji,jj)   &
894                            &                                 + zu_frc(ji,jj) ) & 
895                            &   ) * ssumask(ji,jj)
896
897                  va_e(ji,jj) = (                                 vn_e(ji,jj)   &
898                            &     + rdtbt * (                      zwy(ji,jj)   &
899                            &                                 + zv_trd(ji,jj)   &
900                            &                                 + zv_frc(ji,jj) ) &
901                            &   ) * ssvmask(ji,jj)
902               END DO
903            END DO
904            !
905         ELSE                                      !* Flux form
906            DO jj = 2, jpjm1
907               DO ji = fs_2, fs_jpim1   ! vector opt.
908
909                  zhura = hu_0(ji,jj) + zsshu_a(ji,jj)
910                  zhvra = hv_0(ji,jj) + zsshv_a(ji,jj)
911                  zhura = ssumask(ji,jj)/(zhura + 1._wp - ssumask(ji,jj))
912                  zhvra = ssvmask(ji,jj)/(zhvra + 1._wp - ssvmask(ji,jj))
913
914                  ua_e(ji,jj) = (                hu_e(ji,jj)  *   un_e(ji,jj)   & 
915                            &     + rdtbt * ( zhust_e(ji,jj)  *    zwx(ji,jj)   & 
916                            &               + zhup2_e(ji,jj)  * zu_trd(ji,jj)   &
917                            &               +    hu_n(ji,jj)  * zu_frc(ji,jj) ) &
918                            &   ) * zhura
919
920                  va_e(ji,jj) = (                hv_e(ji,jj)  *   vn_e(ji,jj)   &
921                            &     + rdtbt * ( zhvst_e(ji,jj)  *    zwy(ji,jj)   &
922                            &               + zhvp2_e(ji,jj)  * zv_trd(ji,jj)   &
923                            &               +    hv_n(ji,jj)  * zv_frc(ji,jj) ) &
924                            &   ) * zhvra
925               END DO
926            END DO
927         ENDIF
928         !
929         IF( .NOT.ln_linssh ) THEN                     !* Update ocean depth (variable volume case only)
930            hu_e (:,:) = hu_0(:,:) + zsshu_a(:,:)
931            hv_e (:,:) = hv_0(:,:) + zsshv_a(:,:)
932            hur_e(:,:) = ssumask(:,:) / ( hu_e(:,:) + 1._wp - ssumask(:,:) )
933            hvr_e(:,:) = ssvmask(:,:) / ( hv_e(:,:) + 1._wp - ssvmask(:,:) )
934            !
935         ENDIF
936         !                                             !* domain lateral boundary
937         CALL lbc_lnk_multi( ua_e, 'U', -1._wp, va_e , 'V', -1._wp )
938         !
939         !                                                 ! open boundaries
940         IF( ln_bdy )   CALL bdy_dyn2d( jn, ua_e, va_e, un_e, vn_e, hur_e, hvr_e, ssha_e )
941
942#if defined key_agrif                                                           
943         IF( .NOT.Agrif_Root() )  CALL agrif_dyn_ts( jn )  ! Agrif
944#endif
945         !                                             !* Swap
946         !                                             !  ----
947         ubb_e  (:,:) = ub_e  (:,:)
948         ub_e   (:,:) = un_e  (:,:)
949         un_e   (:,:) = ua_e  (:,:)
950         !
951         vbb_e  (:,:) = vb_e  (:,:)
952         vb_e   (:,:) = vn_e  (:,:)
953         vn_e   (:,:) = va_e  (:,:)
954         !
955         sshbb_e(:,:) = sshb_e(:,:)
956         sshb_e (:,:) = sshn_e(:,:)
957         sshn_e (:,:) = ssha_e(:,:)
958
959         !                                             !* Sum over whole bt loop
960         !                                             !  ----------------------
961         za1 = wgtbtp1(jn)                                   
962         IF( ln_dynadv_vec .OR. ln_linssh ) THEN    ! Sum velocities
963            ua_b  (:,:) = ua_b  (:,:) + za1 * ua_e  (:,:) 
964            va_b  (:,:) = va_b  (:,:) + za1 * va_e  (:,:) 
965         ELSE                                              ! Sum transports
966            ua_b  (:,:) = ua_b  (:,:) + za1 * ua_e  (:,:) * hu_e (:,:)
967            va_b  (:,:) = va_b  (:,:) + za1 * va_e  (:,:) * hv_e (:,:)
968         ENDIF
969         !                                   ! Sum sea level
970         ssha(:,:) = ssha(:,:) + za1 * ssha_e(:,:)
971         !                                                 ! ==================== !
972      END DO                                               !        end loop      !
973      !                                                    ! ==================== !
974      ! -----------------------------------------------------------------------------
975      ! Phase 3. update the general trend with the barotropic trend
976      ! -----------------------------------------------------------------------------
977      !
978      ! Set advection velocity correction:
979      zwx(:,:) = un_adv(:,:)
980      zwy(:,:) = vn_adv(:,:)
981      IF( ( kt == nit000 .AND. neuler==0 ) .OR. .NOT.ln_bt_fw ) THEN     
982         un_adv(:,:) = zwx(:,:) * r1_hu_n(:,:)
983         vn_adv(:,:) = zwy(:,:) * r1_hv_n(:,:)
984      ELSE
985         un_adv(:,:) = z1_2 * ( ub2_b(:,:) + zwx(:,:) ) * r1_hu_n(:,:)
986         vn_adv(:,:) = z1_2 * ( vb2_b(:,:) + zwy(:,:) ) * r1_hv_n(:,:)
987      END IF
988
989      IF( ln_bt_fw ) THEN ! Save integrated transport for next computation
990         ub2_b(:,:) = zwx(:,:)
991         vb2_b(:,:) = zwy(:,:)
992      ENDIF
993      !
994      ! Update barotropic trend:
995      IF(ln_wd) THEN
996        IF( ln_dynadv_vec .OR. ln_linssh ) THEN
997           DO jk=1,jpkm1
998              ua(:,:,jk) = ua(:,:,jk) + ( ua_b(:,:) - ub_b(:,:) ) * z1_2dt_b * wdmask(:,:)
999              va(:,:,jk) = va(:,:,jk) + ( va_b(:,:) - vb_b(:,:) ) * z1_2dt_b * wdmask(:,:)
1000           END DO
1001        ELSE
1002           ! At this stage, ssha has been corrected: compute new depths at velocity points
1003           DO jj = 1, jpjm1
1004              DO ji = 1, jpim1      ! NO Vector Opt.
1005                 zsshu_a(ji,jj) = z1_2 * umask(ji,jj,1)  * r1_e1e2u(ji,jj) &
1006                    &              * ( e1e2t(ji  ,jj) * ssha(ji  ,jj)    &
1007                    &              +   e1e2t(ji+1,jj) * ssha(ji+1,jj) )
1008                 zsshv_a(ji,jj) = z1_2 * vmask(ji,jj,1)  * r1_e1e2v(ji,jj) &
1009                    &              * ( e1e2t(ji,jj  ) * ssha(ji,jj  )    &
1010                    &              +   e1e2t(ji,jj+1) * ssha(ji,jj+1) )
1011              END DO
1012           END DO
1013           CALL lbc_lnk_multi( zsshu_a, 'U', 1._wp, zsshv_a, 'V', 1._wp ) ! Boundary conditions
1014           !
1015           DO jk=1,jpkm1
1016              ua(:,:,jk) = ua(:,:,jk) + r1_hu_n(:,:) * ( ua_b(:,:) - ub_b(:,:) * hu_b(:,:) ) * z1_2dt_b * wdmask(:,:)
1017              va(:,:,jk) = va(:,:,jk) + r1_hv_n(:,:) * ( va_b(:,:) - vb_b(:,:) * hv_b(:,:) ) * z1_2dt_b * wdmask(:,:)
1018           END DO
1019           ! Save barotropic velocities not transport:
1020           ua_b(:,:) =  ua_b(:,:) / ( hu_0(:,:) + zsshu_a(:,:) + 1._wp - ssumask(:,:) )
1021           va_b(:,:) =  va_b(:,:) / ( hv_0(:,:) + zsshv_a(:,:) + 1._wp - ssvmask(:,:) )
1022        ENDIF
1023      ELSE
1024        IF( ln_dynadv_vec .OR. ln_linssh ) THEN
1025           DO jk=1,jpkm1
1026              ua(:,:,jk) = ua(:,:,jk) + ( ua_b(:,:) - ub_b(:,:) ) * z1_2dt_b
1027              va(:,:,jk) = va(:,:,jk) + ( va_b(:,:) - vb_b(:,:) ) * z1_2dt_b
1028           END DO
1029        ELSE
1030           ! At this stage, ssha has been corrected: compute new depths at velocity points
1031           DO jj = 1, jpjm1
1032              DO ji = 1, jpim1      ! NO Vector Opt.
1033                 zsshu_a(ji,jj) = z1_2 * umask(ji,jj,1)  * r1_e1e2u(ji,jj) &
1034                    &              * ( e1e2t(ji  ,jj) * ssha(ji  ,jj)    &
1035                    &              +   e1e2t(ji+1,jj) * ssha(ji+1,jj) )
1036                 zsshv_a(ji,jj) = z1_2 * vmask(ji,jj,1)  * r1_e1e2v(ji,jj) &
1037                    &              * ( e1e2t(ji,jj  ) * ssha(ji,jj  )    &
1038                    &              +   e1e2t(ji,jj+1) * ssha(ji,jj+1) )
1039              END DO
1040           END DO
1041           CALL lbc_lnk_multi( zsshu_a, 'U', 1._wp, zsshv_a, 'V', 1._wp ) ! Boundary conditions
1042           !
1043           DO jk=1,jpkm1
1044              ua(:,:,jk) = ua(:,:,jk) + r1_hu_n(:,:) * ( ua_b(:,:) - ub_b(:,:) * hu_b(:,:) ) * z1_2dt_b
1045              va(:,:,jk) = va(:,:,jk) + r1_hv_n(:,:) * ( va_b(:,:) - vb_b(:,:) * hv_b(:,:) ) * z1_2dt_b
1046           END DO
1047           ! Save barotropic velocities not transport:
1048           ua_b(:,:) =  ua_b(:,:) / ( hu_0(:,:) + zsshu_a(:,:) + 1._wp - ssumask(:,:) )
1049           va_b(:,:) =  va_b(:,:) / ( hv_0(:,:) + zsshv_a(:,:) + 1._wp - ssvmask(:,:) )
1050        ENDIF
1051
1052      END IF
1053      !
1054      DO jk = 1, jpkm1
1055         ! Correct velocities:
1056         un(:,:,jk) = ( un(:,:,jk) + un_adv(:,:) - un_b(:,:) ) * umask(:,:,jk)
1057         vn(:,:,jk) = ( vn(:,:,jk) + vn_adv(:,:) - vn_b(:,:) ) * vmask(:,:,jk)
1058         !
1059      END DO
1060      !
1061      CALL iom_put(  "ubar", un_adv(:,:)      )    ! barotropic i-current
1062      CALL iom_put(  "vbar", vn_adv(:,:)      )    ! barotropic i-current
1063      !
1064#if defined key_agrif
1065      ! Save time integrated fluxes during child grid integration
1066      ! (used to update coarse grid transports at next time step)
1067      !
1068      IF( .NOT.Agrif_Root() .AND. ln_bt_fw ) THEN
1069         IF( Agrif_NbStepint() == 0 ) THEN
1070            ub2_i_b(:,:) = 0._wp
1071            vb2_i_b(:,:) = 0._wp
1072         END IF
1073         !
1074         za1 = 1._wp / REAL(Agrif_rhot(), wp)
1075         ub2_i_b(:,:) = ub2_i_b(:,:) + za1 * ub2_b(:,:)
1076         vb2_i_b(:,:) = vb2_i_b(:,:) + za1 * vb2_b(:,:)
1077      ENDIF
1078#endif     
1079      !                                   !* write time-spliting arrays in the restart
1080      IF( lrst_oce .AND.ln_bt_fw )   CALL ts_rst( kt, 'WRITE' )
1081      !
1082      CALL wrk_dealloc( jpi,jpj,   zsshp2_e, zhdiv )
1083      CALL wrk_dealloc( jpi,jpj,   zu_trd, zv_trd )
1084      CALL wrk_dealloc( jpi,jpj,   zwx, zwy, zssh_frc, zu_frc, zv_frc )
1085      CALL wrk_dealloc( jpi,jpj,   zhup2_e, zhvp2_e, zhust_e, zhvst_e )
1086      CALL wrk_dealloc( jpi,jpj,   zsshu_a, zsshv_a                                   )
1087      CALL wrk_dealloc( jpi,jpj,   zhf )
1088      IF( ln_wd ) CALL wrk_dealloc( jpi, jpj, zcpx, zcpy )
1089      !
1090      IF ( ln_diatmb ) THEN
1091         CALL iom_put( "baro_u" , un_b*umask(:,:,1)+zmdi*(1-umask(:,:,1 ) ) )  ! Barotropic  U Velocity
1092         CALL iom_put( "baro_v" , vn_b*vmask(:,:,1)+zmdi*(1-vmask(:,:,1 ) ) )  ! Barotropic  V Velocity
1093      ENDIF
1094      IF( nn_timing == 1 )  CALL timing_stop('dyn_spg_ts')
1095      !
1096   END SUBROUTINE dyn_spg_ts
1097
1098
1099   SUBROUTINE ts_wgt( ll_av, ll_fw, jpit, zwgt1, zwgt2)
1100      !!---------------------------------------------------------------------
1101      !!                   ***  ROUTINE ts_wgt  ***
1102      !!
1103      !! ** Purpose : Set time-splitting weights for temporal averaging (or not)
1104      !!----------------------------------------------------------------------
1105      LOGICAL, INTENT(in) ::   ll_av      ! temporal averaging=.true.
1106      LOGICAL, INTENT(in) ::   ll_fw      ! forward time splitting =.true.
1107      INTEGER, INTENT(inout) :: jpit      ! cycle length   
1108      REAL(wp), DIMENSION(3*nn_baro), INTENT(inout) ::   zwgt1, & ! Primary weights
1109                                                         zwgt2    ! Secondary weights
1110     
1111      INTEGER ::  jic, jn, ji                      ! temporary integers
1112      REAL(wp) :: za1, za2
1113      !!----------------------------------------------------------------------
1114
1115      zwgt1(:) = 0._wp
1116      zwgt2(:) = 0._wp
1117
1118      ! Set time index when averaged value is requested
1119      IF (ll_fw) THEN
1120         jic = nn_baro
1121      ELSE
1122         jic = 2 * nn_baro
1123      ENDIF
1124
1125      ! Set primary weights:
1126      IF (ll_av) THEN
1127           ! Define simple boxcar window for primary weights
1128           ! (width = nn_baro, centered around jic)     
1129         SELECT CASE ( nn_bt_flt )
1130              CASE( 0 )  ! No averaging
1131                 zwgt1(jic) = 1._wp
1132                 jpit = jic
1133
1134              CASE( 1 )  ! Boxcar, width = nn_baro
1135                 DO jn = 1, 3*nn_baro
1136                    za1 = ABS(float(jn-jic))/float(nn_baro) 
1137                    IF (za1 < 0.5_wp) THEN
1138                      zwgt1(jn) = 1._wp
1139                      jpit = jn
1140                    ENDIF
1141                 ENDDO
1142
1143              CASE( 2 )  ! Boxcar, width = 2 * nn_baro
1144                 DO jn = 1, 3*nn_baro
1145                    za1 = ABS(float(jn-jic))/float(nn_baro) 
1146                    IF (za1 < 1._wp) THEN
1147                      zwgt1(jn) = 1._wp
1148                      jpit = jn
1149                    ENDIF
1150                 ENDDO
1151              CASE DEFAULT   ;   CALL ctl_stop( 'unrecognised value for nn_bt_flt' )
1152         END SELECT
1153
1154      ELSE ! No time averaging
1155         zwgt1(jic) = 1._wp
1156         jpit = jic
1157      ENDIF
1158   
1159      ! Set secondary weights
1160      DO jn = 1, jpit
1161        DO ji = jn, jpit
1162             zwgt2(jn) = zwgt2(jn) + zwgt1(ji)
1163        END DO
1164      END DO
1165
1166      ! Normalize weigths:
1167      za1 = 1._wp / SUM(zwgt1(1:jpit))
1168      za2 = 1._wp / SUM(zwgt2(1:jpit))
1169      DO jn = 1, jpit
1170        zwgt1(jn) = zwgt1(jn) * za1
1171        zwgt2(jn) = zwgt2(jn) * za2
1172      END DO
1173      !
1174   END SUBROUTINE ts_wgt
1175
1176
1177   SUBROUTINE ts_rst( kt, cdrw )
1178      !!---------------------------------------------------------------------
1179      !!                   ***  ROUTINE ts_rst  ***
1180      !!
1181      !! ** Purpose : Read or write time-splitting arrays in restart file
1182      !!----------------------------------------------------------------------
1183      INTEGER         , INTENT(in) ::   kt         ! ocean time-step
1184      CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
1185      !
1186      !!----------------------------------------------------------------------
1187      !
1188      IF( TRIM(cdrw) == 'READ' ) THEN
1189         CALL iom_get( numror, jpdom_autoglo, 'ub2_b'  , ub2_b  (:,:) )   
1190         CALL iom_get( numror, jpdom_autoglo, 'vb2_b'  , vb2_b  (:,:) ) 
1191         IF( .NOT.ln_bt_av ) THEN
1192            CALL iom_get( numror, jpdom_autoglo, 'sshbb_e'  , sshbb_e(:,:) )   
1193            CALL iom_get( numror, jpdom_autoglo, 'ubb_e'    ,   ubb_e(:,:) )   
1194            CALL iom_get( numror, jpdom_autoglo, 'vbb_e'    ,   vbb_e(:,:) )
1195            CALL iom_get( numror, jpdom_autoglo, 'sshb_e'   ,  sshb_e(:,:) ) 
1196            CALL iom_get( numror, jpdom_autoglo, 'ub_e'     ,    ub_e(:,:) )   
1197            CALL iom_get( numror, jpdom_autoglo, 'vb_e'     ,    vb_e(:,:) )
1198         ENDIF
1199#if defined key_agrif
1200         ! Read time integrated fluxes
1201         IF ( .NOT.Agrif_Root() ) THEN
1202            CALL iom_get( numror, jpdom_autoglo, 'ub2_i_b'  , ub2_i_b(:,:) )   
1203            CALL iom_get( numror, jpdom_autoglo, 'vb2_i_b'  , vb2_i_b(:,:) )
1204         ENDIF
1205#endif
1206      !
1207      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN
1208         CALL iom_rstput( kt, nitrst, numrow, 'ub2_b'   , ub2_b  (:,:) )
1209         CALL iom_rstput( kt, nitrst, numrow, 'vb2_b'   , vb2_b  (:,:) )
1210         !
1211         IF (.NOT.ln_bt_av) THEN
1212            CALL iom_rstput( kt, nitrst, numrow, 'sshbb_e'  , sshbb_e(:,:) ) 
1213            CALL iom_rstput( kt, nitrst, numrow, 'ubb_e'    ,   ubb_e(:,:) )
1214            CALL iom_rstput( kt, nitrst, numrow, 'vbb_e'    ,   vbb_e(:,:) )
1215            CALL iom_rstput( kt, nitrst, numrow, 'sshb_e'   ,  sshb_e(:,:) )
1216            CALL iom_rstput( kt, nitrst, numrow, 'ub_e'     ,    ub_e(:,:) )
1217            CALL iom_rstput( kt, nitrst, numrow, 'vb_e'     ,    vb_e(:,:) )
1218         ENDIF
1219#if defined key_agrif
1220         ! Save time integrated fluxes
1221         IF ( .NOT.Agrif_Root() ) THEN
1222            CALL iom_rstput( kt, nitrst, numrow, 'ub2_i_b'  , ub2_i_b(:,:) )
1223            CALL iom_rstput( kt, nitrst, numrow, 'vb2_i_b'  , vb2_i_b(:,:) )
1224         ENDIF
1225#endif
1226      ENDIF
1227      !
1228   END SUBROUTINE ts_rst
1229
1230
1231   SUBROUTINE dyn_spg_ts_init
1232      !!---------------------------------------------------------------------
1233      !!                   ***  ROUTINE dyn_spg_ts_init  ***
1234      !!
1235      !! ** Purpose : Set time splitting options
1236      !!----------------------------------------------------------------------
1237      INTEGER  ::   ji ,jj              ! dummy loop indices
1238      REAL(wp) ::   zxr2, zyr2, zcmax   ! local scalar
1239      REAL(wp), POINTER, DIMENSION(:,:) ::   zcu
1240      !!----------------------------------------------------------------------
1241      !
1242      ! Max courant number for ext. grav. waves
1243      !
1244      CALL wrk_alloc( jpi,jpj,   zcu )
1245      !
1246      DO jj = 1, jpj
1247         DO ji =1, jpi
1248            zxr2 = r1_e1t(ji,jj) * r1_e1t(ji,jj)
1249            zyr2 = r1_e2t(ji,jj) * r1_e2t(ji,jj)
1250            zcu(ji,jj) = SQRT( grav * ht_0(ji,jj) * (zxr2 + zyr2) )
1251         END DO
1252      END DO
1253      !
1254      zcmax = MAXVAL( zcu(:,:) )
1255      IF( lk_mpp )   CALL mpp_max( zcmax )
1256
1257      ! Estimate number of iterations to satisfy a max courant number= rn_bt_cmax
1258      IF( ln_bt_auto )   nn_baro = CEILING( rdt / rn_bt_cmax * zcmax)
1259     
1260      rdtbt = rdt / REAL( nn_baro , wp )
1261      zcmax = zcmax * rdtbt
1262                     ! Print results
1263      IF(lwp) WRITE(numout,*)
1264      IF(lwp) WRITE(numout,*) 'dyn_spg_ts : split-explicit free surface'
1265      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
1266      IF( ln_bt_auto ) THEN
1267         IF(lwp) WRITE(numout,*) '     ln_ts_auto=.true. Automatically set nn_baro '
1268         IF(lwp) WRITE(numout,*) '     Max. courant number allowed: ', rn_bt_cmax
1269      ELSE
1270         IF(lwp) WRITE(numout,*) '     ln_ts_auto=.false.: Use nn_baro in namelist '
1271      ENDIF
1272
1273      IF(ln_bt_av) THEN
1274         IF(lwp) WRITE(numout,*) '     ln_bt_av=.true.  => Time averaging over nn_baro time steps is on '
1275      ELSE
1276         IF(lwp) WRITE(numout,*) '     ln_bt_av=.false. => No time averaging of barotropic variables '
1277      ENDIF
1278      !
1279      !
1280      IF(ln_bt_fw) THEN
1281         IF(lwp) WRITE(numout,*) '     ln_bt_fw=.true.  => Forward integration of barotropic variables '
1282      ELSE
1283         IF(lwp) WRITE(numout,*) '     ln_bt_fw =.false.=> Centred integration of barotropic variables '
1284      ENDIF
1285      !
1286#if defined key_agrif
1287      ! Restrict the use of Agrif to the forward case only
1288      IF( .NOT.ln_bt_fw .AND. .NOT.Agrif_Root() )   CALL ctl_stop( 'AGRIF not implemented if ln_bt_fw=.FALSE.' )
1289#endif
1290      !
1291      IF(lwp) WRITE(numout,*)    '     Time filter choice, nn_bt_flt: ', nn_bt_flt
1292      SELECT CASE ( nn_bt_flt )
1293         CASE( 0 )      ;   IF(lwp) WRITE(numout,*) '           Dirac'
1294         CASE( 1 )      ;   IF(lwp) WRITE(numout,*) '           Boxcar: width = nn_baro'
1295         CASE( 2 )      ;   IF(lwp) WRITE(numout,*) '           Boxcar: width = 2*nn_baro' 
1296         CASE DEFAULT   ;   CALL ctl_stop( 'unrecognised value for nn_bt_flt: should 0,1,2' )
1297      END SELECT
1298      !
1299      IF(lwp) WRITE(numout,*) ' '
1300      IF(lwp) WRITE(numout,*) '     nn_baro = ', nn_baro
1301      IF(lwp) WRITE(numout,*) '     Barotropic time step [s] is :', rdtbt
1302      IF(lwp) WRITE(numout,*) '     Maximum Courant number is   :', zcmax
1303      !
1304      IF( .NOT.ln_bt_av .AND. .NOT.ln_bt_fw ) THEN
1305         CALL ctl_stop( 'dynspg_ts ERROR: No time averaging => only forward integration is possible' )
1306      ENDIF
1307      IF( zcmax>0.9_wp ) THEN
1308         CALL ctl_stop( 'dynspg_ts ERROR: Maximum Courant number is greater than 0.9: Inc. nn_baro !' )         
1309      ENDIF
1310      !
1311      CALL wrk_dealloc( jpi,jpj,   zcu )
1312      !
1313   END SUBROUTINE dyn_spg_ts_init
1314
1315   !!======================================================================
1316END MODULE dynspg_ts
Note: See TracBrowser for help on using the repository browser.