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

source: branches/2016/dev_r6522_SIMPLIF_3/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_ts.F90 @ 6864

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

#1729 - trunk: removed key_tide from the code and set usage of ln_tide. Tested with AMM12 and ORCA2_LIM_PISCES.

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