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 @ 7437

Last change on this file since 7437 was 7437, checked in by timgraham, 8 years ago

Remove wd[uv]flt variables

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