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

Last change on this file since 7439 was 7439, checked in by lovato, 7 years ago

#1811 - Reinstate developments from dev_NOC_CMCC branch lost in merge

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