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 NEMO/branches/2019/dev_r12072_MERGE_OPTION2_2019/src/OCE/DYN – NEMO

source: NEMO/branches/2019/dev_r12072_MERGE_OPTION2_2019/src/OCE/DYN/dynspg_ts.F90 @ 12246

Last change on this file since 12246 was 12246, checked in by smasson, 4 years ago

rev12232_dev_r12072_MERGE_OPTION2_2019: add modifications from dev_r12114_ticket_2263, results unchanged except SPITZ12 as explained in #2263

  • Property svn:keywords set to Id
File size: 80.9 KB
Line 
1MODULE dynspg_ts
2
3   !! Includes ROMS wd scheme with diagnostic outputs ; un and ua updates are commented out !
4
5   !!======================================================================
6   !!                   ***  MODULE  dynspg_ts  ***
7   !! Ocean dynamics:  surface pressure gradient trend, split-explicit scheme
8   !!======================================================================
9   !! History :   1.0  ! 2004-12  (L. Bessieres, G. Madec)  Original code
10   !!              -   ! 2005-11  (V. Garnier, G. Madec)  optimization
11   !!              -   ! 2006-08  (S. Masson)  distributed restart using iom
12   !!             2.0  ! 2007-07  (D. Storkey) calls to BDY routines
13   !!              -   ! 2008-01  (R. Benshila)  change averaging method
14   !!             3.2  ! 2009-07  (R. Benshila, G. Madec) Complete revisit associated to vvl reactivation
15   !!             3.3  ! 2010-09  (D. Storkey, E. O'Dea) update for BDY for Shelf configurations
16   !!             3.3  ! 2011-03  (R. Benshila, R. Hordoir, P. Oddo) update calculation of ub_b
17   !!             3.5  ! 2013-07  (J. Chanut) Switch to Forward-backward time stepping
18   !!             3.6  ! 2013-11  (A. Coward) Update for z-tilde compatibility
19   !!             3.7  ! 2015-11  (J. Chanut) free surface simplification
20   !!              -   ! 2016-12  (G. Madec, E. Clementi) update for Stoke-Drift divergence
21   !!             4.0  ! 2017-05  (G. Madec)  drag coef. defined at t-point (zdfdrg.F90)
22   !!---------------------------------------------------------------------
23
24   !!----------------------------------------------------------------------
25   !!   dyn_spg_ts     : compute surface pressure gradient trend using a time-splitting scheme
26   !!   dyn_spg_ts_init: initialisation of the time-splitting scheme
27   !!   ts_wgt         : set time-splitting weights for temporal averaging (or not)
28   !!   ts_rst         : read/write time-splitting fields in restart file
29   !!----------------------------------------------------------------------
30   USE oce             ! ocean dynamics and tracers
31   USE dom_oce         ! ocean space and time domain
32   USE sbc_oce         ! surface boundary condition: ocean
33   USE isf_oce         ! ice shelf variable (fwfisf)
34   USE zdf_oce         ! vertical physics: variables
35   USE zdfdrg          ! vertical physics: top/bottom drag coef.
36   USE sbcapr          ! surface boundary condition: atmospheric pressure
37   USE dynadv    , ONLY: ln_dynadv_vec
38   USE dynvor          ! vortivity scheme indicators
39   USE phycst          ! physical constants
40   USE dynvor          ! vorticity term
41   USE wet_dry         ! wetting/drying flux limter
42   USE bdy_oce         ! open boundary
43   USE bdyvol          ! open boundary volume conservation
44   USE bdytides        ! open boundary condition data
45   USE bdydyn2d        ! open boundary conditions on barotropic variables
46   USE tide_mod        !
47   USE sbcwave         ! surface wave
48   USE diatmb          ! Top,middle,bottom output
49#if defined key_agrif
50   USE agrif_oce_interp ! agrif
51   USE agrif_oce
52#endif
53#if defined key_asminc   
54   USE asminc          ! Assimilation increment
55#endif
56   !
57   USE in_out_manager  ! I/O manager
58   USE lib_mpp         ! distributed memory computing library
59   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
60   USE prtctl          ! Print control
61   USE iom             ! IOM library
62   USE restart         ! only for lrst_oce
63   USE diatmb          ! Top,middle,bottom output
64
65   USE iom   ! to remove
66
67   IMPLICIT NONE
68   PRIVATE
69
70   PUBLIC dyn_spg_ts        ! called by dyn_spg
71   PUBLIC dyn_spg_ts_init   !    -    - dyn_spg_init
72
73   !! Time filtered arrays at baroclinic time step:
74   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   un_adv , vn_adv   !: Advection vel. at "now" barocl. step
75   !
76   INTEGER, SAVE :: icycle      ! Number of barotropic sub-steps for each internal step nn_baro <= 2.5 nn_baro
77   REAL(wp),SAVE :: rdtbt       ! Barotropic time step
78   !
79   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:)   ::   wgtbtp1, wgtbtp2   ! 1st & 2nd weights used in time filtering of barotropic fields
80   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   zwz                ! ff_f/h at F points
81   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   ftnw, ftne         ! triad of coriolis parameter
82   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   ftsw, ftse         ! (only used with een vorticity scheme)
83
84   REAL(wp) ::   r1_12 = 1._wp / 12._wp   ! local ratios
85   REAL(wp) ::   r1_8  = 0.125_wp         !
86   REAL(wp) ::   r1_4  = 0.25_wp          !
87   REAL(wp) ::   r1_2  = 0.5_wp           !
88
89   !! * Substitutions
90#  include "vectopt_loop_substitute.h90"
91   !!----------------------------------------------------------------------
92   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
93   !! $Id$
94   !! Software governed by the CeCILL license (see ./LICENSE)
95   !!----------------------------------------------------------------------
96CONTAINS
97
98   INTEGER FUNCTION dyn_spg_ts_alloc()
99      !!----------------------------------------------------------------------
100      !!                  ***  routine dyn_spg_ts_alloc  ***
101      !!----------------------------------------------------------------------
102      INTEGER :: ierr(3)
103      !!----------------------------------------------------------------------
104      ierr(:) = 0
105      !
106      ALLOCATE( wgtbtp1(3*nn_baro), wgtbtp2(3*nn_baro), zwz(jpi,jpj), STAT=ierr(1) )
107      IF( ln_dynvor_een .OR. ln_dynvor_eeT )   &
108         &     ALLOCATE( ftnw(jpi,jpj) , ftne(jpi,jpj) , ftsw(jpi,jpj) , ftse(jpi,jpj), STAT=ierr(2)   )
109         !
110      ALLOCATE( un_adv(jpi,jpj), vn_adv(jpi,jpj)                    , STAT=ierr(3) )
111      !
112      dyn_spg_ts_alloc = MAXVAL( ierr(:) )
113      !
114      CALL mpp_sum( 'dynspg_ts', dyn_spg_ts_alloc )
115      IF( dyn_spg_ts_alloc /= 0 )   CALL ctl_stop( 'STOP', 'dyn_spg_ts_alloc: failed to allocate arrays' )
116      !
117   END FUNCTION dyn_spg_ts_alloc
118
119
120   SUBROUTINE dyn_spg_ts( kt )
121      !!----------------------------------------------------------------------
122      !!
123      !! ** Purpose : - Compute the now trend due to the explicit time stepping
124      !!              of the quasi-linear barotropic system, and add it to the
125      !!              general momentum trend.
126      !!
127      !! ** Method  : - split-explicit schem (time splitting) :
128      !!      Barotropic variables are advanced from internal time steps
129      !!      "n"   to "n+1" if ln_bt_fw=T
130      !!      or from
131      !!      "n-1" to "n+1" if ln_bt_fw=F
132      !!      thanks to a generalized forward-backward time stepping (see ref. below).
133      !!
134      !! ** Action :
135      !!      -Update the filtered free surface at step "n+1"      : ssha
136      !!      -Update filtered barotropic velocities at step "n+1" : ua_b, va_b
137      !!      -Compute barotropic advective fluxes at step "n"     : un_adv, vn_adv
138      !!      These are used to advect tracers and are compliant with discrete
139      !!      continuity equation taken at the baroclinic time steps. This
140      !!      ensures tracers conservation.
141      !!      - (ua, va) momentum trend updated with barotropic component.
142      !!
143      !! References : Shchepetkin and McWilliams, Ocean Modelling, 2005.
144      !!---------------------------------------------------------------------
145      INTEGER, INTENT(in)  ::   kt   ! ocean time-step index
146      !
147      INTEGER  ::   ji, jj, jk, jn        ! dummy loop indices
148      LOGICAL  ::   ll_fw_start           ! =T : forward integration
149      LOGICAL  ::   ll_init               ! =T : special startup of 2d equations
150      INTEGER  ::   noffset               ! local integers  : time offset for bdy update
151      REAL(wp) ::   r1_2dt_b, z1_hu, z1_hv          ! local scalars
152      REAL(wp) ::   za0, za1, za2, za3              !   -      -
153      REAL(wp) ::   zmdi, zztmp, zldg               !   -      -
154      REAL(wp) ::   zhu_bck, zhv_bck, zhdiv         !   -      -
155      REAL(wp) ::   zun_save, zvn_save              !   -      -
156      REAL(wp), DIMENSION(jpi,jpj) :: zu_trd, zu_frc, zu_spg, zssh_frc
157      REAL(wp), DIMENSION(jpi,jpj) :: zv_trd, zv_frc, zv_spg
158      REAL(wp), DIMENSION(jpi,jpj) :: zsshu_a, zhup2_e, zhtp2_e
159      REAL(wp), DIMENSION(jpi,jpj) :: zsshv_a, zhvp2_e, zsshp2_e
160      REAL(wp), DIMENSION(jpi,jpj) :: zCdU_u, zCdU_v   ! top/bottom stress at u- & v-points
161      REAL(wp), DIMENSION(jpi,jpj) :: zhU, zhV         ! fluxes
162      !
163      REAL(wp) ::   zwdramp                     ! local scalar - only used if ln_wd_dl = .True.
164
165      INTEGER  :: iwdg, jwdg, kwdg   ! short-hand values for the indices of the output point
166
167      REAL(wp) ::   zepsilon, zgamma            !   -      -
168      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zcpx, zcpy   ! Wetting/Dying gravity filter coef.
169      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: ztwdmask, zuwdmask, zvwdmask ! ROMS wetting and drying masks at t,u,v points
170      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zuwdav2, zvwdav2    ! averages over the sub-steps of zuwdmask and zvwdmask
171      REAL(wp) ::   zt0substep !   Time of day at the beginning of the time substep
172      !!----------------------------------------------------------------------
173      !
174      IF( ln_wd_il ) ALLOCATE( zcpx(jpi,jpj), zcpy(jpi,jpj) )
175      !                                         !* Allocate temporary arrays
176      IF( ln_wd_dl ) ALLOCATE( ztwdmask(jpi,jpj), zuwdmask(jpi,jpj), zvwdmask(jpi,jpj), zuwdav2(jpi,jpj), zvwdav2(jpi,jpj))
177      !
178      zmdi=1.e+20                               !  missing data indicator for masking
179      !
180      zwdramp = r_rn_wdmin1               ! simplest ramp
181!     zwdramp = 1._wp / (rn_wdmin2 - rn_wdmin1) ! more general ramp
182      !                                         ! inverse of baroclinic time step
183      IF( kt == nit000 .AND. neuler == 0 ) THEN   ;   r1_2dt_b = 1._wp / (         rdt )
184      ELSE                                        ;   r1_2dt_b = 1._wp / ( 2._wp * rdt )
185      ENDIF
186      !
187      ll_init     = ln_bt_av                    ! if no time averaging, then no specific restart
188      ll_fw_start = .FALSE.
189      !                                         ! time offset in steps for bdy data update
190      IF( .NOT.ln_bt_fw ) THEN   ;   noffset = - nn_baro
191      ELSE                       ;   noffset =   0 
192      ENDIF
193      !
194      IF( kt == nit000 ) THEN                   !* initialisation
195         !
196         IF(lwp) WRITE(numout,*)
197         IF(lwp) WRITE(numout,*) 'dyn_spg_ts : surface pressure gradient trend'
198         IF(lwp) WRITE(numout,*) '~~~~~~~~~~   free surface with time splitting'
199         IF(lwp) WRITE(numout,*)
200         !
201         IF( neuler == 0 )   ll_init=.TRUE.
202         !
203         IF( ln_bt_fw .OR. neuler == 0 ) THEN
204            ll_fw_start =.TRUE.
205            noffset     = 0
206         ELSE
207            ll_fw_start =.FALSE.
208         ENDIF
209         !                    ! Set averaging weights and cycle length:
210         CALL ts_wgt( ln_bt_av, ll_fw_start, icycle, wgtbtp1, wgtbtp2 )
211         !
212      ENDIF
213      !
214      ! If forward start at previous time step, and centered integration,
215      ! then update averaging weights:
216      IF (.NOT.ln_bt_fw .AND.( neuler==0 .AND. kt==nit000+1 ) ) THEN
217         ll_fw_start=.FALSE.
218         CALL ts_wgt( ln_bt_av, ll_fw_start, icycle, wgtbtp1, wgtbtp2 )
219      ENDIF
220      !
221                         
222      ! -----------------------------------------------------------------------------
223      !  Phase 1 : Coupling between general trend and barotropic estimates (1st step)
224      ! -----------------------------------------------------------------------------
225      !     
226      !
227      !                                   !=  zu_frc =  1/H e3*d/dt(Ua)  =!  (Vertical mean of Ua, the 3D trends)
228      !                                   !  ---------------------------  !
229      zu_frc(:,:) = SUM( e3u_n(:,:,:) * ua(:,:,:) * umask(:,:,:) , DIM=3 ) * r1_hu_n(:,:)
230      zv_frc(:,:) = SUM( e3v_n(:,:,:) * va(:,:,:) * vmask(:,:,:) , DIM=3 ) * r1_hv_n(:,:)
231      !
232      !
233      !                                   !=  Ua => baroclinic trend  =!   (remove its vertical mean)
234      DO jk = 1, jpkm1                    !  ------------------------  !
235         ua(:,:,jk) = ( ua(:,:,jk) - zu_frc(:,:) ) * umask(:,:,jk)
236         va(:,:,jk) = ( va(:,:,jk) - zv_frc(:,:) ) * vmask(:,:,jk)
237      END DO
238     
239!!gm  Question here when removing the Vertically integrated trends, we remove the vertically integrated NL trends on momentum....
240!!gm  Is it correct to do so ?   I think so...
241     
242      !                                   !=  remove 2D Coriolis and pressure gradient trends  =!
243      !                                   !  -------------------------------------------------  !
244      !
245      IF( kt == nit000 .OR. .NOT. ln_linssh )   CALL dyn_cor_2D_init   ! Set zwz, the barotropic Coriolis force coefficient
246      !       ! recompute zwz = f/depth  at every time step for (.NOT.ln_linssh) as the water colomn height changes
247      !
248      !                                         !* 2D Coriolis trends
249      zhU(:,:) = un_b(:,:) * hu_n(:,:) * e2u(:,:)        ! now fluxes
250      zhV(:,:) = vn_b(:,:) * hv_n(:,:) * e1v(:,:)        ! NB: FULL domain : put a value in last row and column
251      !
252      CALL dyn_cor_2d( hu_n, hv_n, un_b, vn_b, zhU, zhV,  &   ! <<== in
253         &                               zu_trd, zv_trd   )   ! ==>> out
254      !
255      IF( .NOT.ln_linssh ) THEN                 !* surface pressure gradient   (variable volume only)
256         !
257         IF( ln_wd_il ) THEN                       ! W/D : limiter applied to spgspg
258            CALL wad_spg( sshn, zcpx, zcpy )          ! Calculating W/D gravity filters, zcpx and zcpy
259            DO jj = 2, jpjm1
260               DO ji = 2, jpim1                ! SPG with the application of W/D gravity filters
261                  zu_trd(ji,jj) = zu_trd(ji,jj) - grav * ( sshn(ji+1,jj  ) - sshn(ji  ,jj ) )   &
262                     &                          * r1_e1u(ji,jj) * zcpx(ji,jj)  * wdrampu(ji,jj)  !jth
263                  zv_trd(ji,jj) = zv_trd(ji,jj) - grav * ( sshn(ji  ,jj+1) - sshn(ji  ,jj ) )   &
264                     &                          * r1_e2v(ji,jj) * zcpy(ji,jj)  * wdrampv(ji,jj)  !jth
265               END DO
266            END DO
267         ELSE                                      ! now suface pressure gradient
268            DO jj = 2, jpjm1
269               DO ji = fs_2, fs_jpim1   ! vector opt.
270                  zu_trd(ji,jj) = zu_trd(ji,jj) - grav * (  sshn(ji+1,jj  ) - sshn(ji  ,jj  )  ) * r1_e1u(ji,jj)
271                  zv_trd(ji,jj) = zv_trd(ji,jj) - grav * (  sshn(ji  ,jj+1) - sshn(ji  ,jj  )  ) * r1_e2v(ji,jj) 
272               END DO
273            END DO
274         ENDIF
275         !
276      ENDIF
277      !
278      DO jj = 2, jpjm1                          ! Remove coriolis term (and possibly spg) from barotropic trend
279         DO ji = fs_2, fs_jpim1
280             zu_frc(ji,jj) = zu_frc(ji,jj) - zu_trd(ji,jj) * ssumask(ji,jj)
281             zv_frc(ji,jj) = zv_frc(ji,jj) - zv_trd(ji,jj) * ssvmask(ji,jj)
282          END DO
283      END DO 
284      !
285      !                                   !=  Add bottom stress contribution from baroclinic velocities  =!
286      !                                   !  -----------------------------------------------------------  !
287      CALL dyn_drg_init( zu_frc, zv_frc,  zCdU_u, zCdU_v )      ! also provide the barotropic drag coefficients
288      !
289      !                                   !=  Add atmospheric pressure forcing  =!
290      !                                   !  ----------------------------------  !
291      IF( ln_apr_dyn ) THEN
292         IF( ln_bt_fw ) THEN                          ! FORWARD integration: use kt+1/2 pressure (NOW+1/2)
293            DO jj = 2, jpjm1             
294               DO ji = fs_2, fs_jpim1   ! vector opt.
295                  zu_frc(ji,jj) = zu_frc(ji,jj) + grav * (  ssh_ib (ji+1,jj  ) - ssh_ib (ji,jj) ) * r1_e1u(ji,jj)
296                  zv_frc(ji,jj) = zv_frc(ji,jj) + grav * (  ssh_ib (ji  ,jj+1) - ssh_ib (ji,jj) ) * r1_e2v(ji,jj)
297               END DO
298            END DO
299         ELSE                                         ! CENTRED integration: use kt-1/2 + kt+1/2 pressure (NOW)
300            zztmp = grav * r1_2
301            DO jj = 2, jpjm1             
302               DO ji = fs_2, fs_jpim1   ! vector opt.
303                  zu_frc(ji,jj) = zu_frc(ji,jj) + zztmp * (  ssh_ib (ji+1,jj  ) - ssh_ib (ji,jj)  &
304                       &                                   + ssh_ibb(ji+1,jj  ) - ssh_ibb(ji,jj)  ) * r1_e1u(ji,jj)
305                  zv_frc(ji,jj) = zv_frc(ji,jj) + zztmp * (  ssh_ib (ji  ,jj+1) - ssh_ib (ji,jj)  &
306                       &                                   + ssh_ibb(ji  ,jj+1) - ssh_ibb(ji,jj)  ) * r1_e2v(ji,jj)
307               END DO
308            END DO
309         ENDIF
310      ENDIF
311      !
312      !                                   !=  Add atmospheric pressure forcing  =!
313      !                                   !  ----------------------------------  !
314      IF( ln_bt_fw ) THEN                        ! Add wind forcing
315         DO jj = 2, jpjm1
316            DO ji = fs_2, fs_jpim1   ! vector opt.
317               zu_frc(ji,jj) =  zu_frc(ji,jj) + r1_rau0 * utau(ji,jj) * r1_hu_n(ji,jj)
318               zv_frc(ji,jj) =  zv_frc(ji,jj) + r1_rau0 * vtau(ji,jj) * r1_hv_n(ji,jj)
319            END DO
320         END DO
321      ELSE
322         zztmp = r1_rau0 * r1_2
323         DO jj = 2, jpjm1
324            DO ji = fs_2, fs_jpim1   ! vector opt.
325               zu_frc(ji,jj) =  zu_frc(ji,jj) + zztmp * ( utau_b(ji,jj) + utau(ji,jj) ) * r1_hu_n(ji,jj)
326               zv_frc(ji,jj) =  zv_frc(ji,jj) + zztmp * ( vtau_b(ji,jj) + vtau(ji,jj) ) * r1_hv_n(ji,jj)
327            END DO
328         END DO
329      ENDIF 
330      !
331      !              !----------------!
332      !              !==  sssh_frc  ==!   Right-Hand-Side of the barotropic ssh equation   (over the FULL domain)
333      !              !----------------!
334      !                                   !=  Net water flux forcing applied to a water column  =!
335      !                                   ! ---------------------------------------------------  !
336      IF (ln_bt_fw) THEN                          ! FORWARD integration: use kt+1/2 fluxes (NOW+1/2)
337         zssh_frc(:,:) = r1_rau0 * ( emp(:,:) - rnf(:,:) + fwfisf_cav(:,:) + fwfisf_par(:,:) )
338      ELSE                                        ! CENTRED integration: use kt-1/2 + kt+1/2 fluxes (NOW)
339         zztmp = r1_rau0 * r1_2
340         zssh_frc(:,:) = zztmp * (  emp(:,:)        + emp_b(:,:)                    &
341                &                 - rnf(:,:)        - rnf_b(:,:)                    &
342                &                 + fwfisf_cav(:,:) + fwfisf_cav_b(:,:)             &
343                &                 + fwfisf_par(:,:) + fwfisf_par_b(:,:)             )
344      ENDIF
345      !                                   !=  Add Stokes drift divergence  =!   (if exist)
346      IF( ln_sdw ) THEN                   !  -----------------------------  !
347         zssh_frc(:,:) = zssh_frc(:,:) + div_sd(:,:)
348      ENDIF
349      !
350      !                                         ! ice sheet coupling
351      IF ( ln_isf .AND. ln_isfcpl ) THEN
352         !
353         ! ice sheet coupling
354         IF( ln_rstart .AND. kt == nit000 ) THEN
355            zssh_frc(:,:) = zssh_frc(:,:) + risfcpl_ssh(:,:)
356         END IF
357         !
358         ! conservation option
359         IF( ln_isfcpl_cons ) THEN
360            zssh_frc(:,:) = zssh_frc(:,:) + risfcpl_cons_ssh(:,:)
361         END IF
362         !
363      END IF
364      !
365#if defined key_asminc
366      !                                   !=  Add the IAU weighted SSH increment  =!
367      !                                   !  ------------------------------------  !
368      IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN
369         zssh_frc(:,:) = zssh_frc(:,:) - ssh_iau(:,:)
370      ENDIF
371#endif
372      !                                   != Fill boundary data arrays for AGRIF
373      !                                   ! ------------------------------------
374#if defined key_agrif
375         IF( .NOT.Agrif_Root() ) CALL agrif_dta_ts( kt )
376#endif
377      !
378      ! -----------------------------------------------------------------------
379      !  Phase 2 : Integration of the barotropic equations
380      ! -----------------------------------------------------------------------
381      !
382      !                                             ! ==================== !
383      !                                             !    Initialisations   !
384      !                                             ! ==================== ! 
385      ! Initialize barotropic variables:     
386      IF( ll_init )THEN
387         sshbb_e(:,:) = 0._wp
388         ubb_e  (:,:) = 0._wp
389         vbb_e  (:,:) = 0._wp
390         sshb_e (:,:) = 0._wp
391         ub_e   (:,:) = 0._wp
392         vb_e   (:,:) = 0._wp
393      ENDIF
394      !
395      IF( ln_linssh ) THEN    ! mid-step ocean depth is fixed (hup2_e=hu_n=hu_0)
396         zhup2_e(:,:) = hu_n(:,:)
397         zhvp2_e(:,:) = hv_n(:,:)
398         zhtp2_e(:,:) = ht_n(:,:)
399      ENDIF
400      !
401      IF (ln_bt_fw) THEN                  ! FORWARD integration: start from NOW fields                   
402         sshn_e(:,:) =    sshn(:,:)           
403         un_e  (:,:) =    un_b(:,:)           
404         vn_e  (:,:) =    vn_b(:,:)
405         !
406         hu_e  (:,:) =    hu_n(:,:)       
407         hv_e  (:,:) =    hv_n(:,:) 
408         hur_e (:,:) = r1_hu_n(:,:)   
409         hvr_e (:,:) = r1_hv_n(:,:)
410      ELSE                                ! CENTRED integration: start from BEFORE fields
411         sshn_e(:,:) =    sshb(:,:)
412         un_e  (:,:) =    ub_b(:,:)         
413         vn_e  (:,:) =    vb_b(:,:)
414         !
415         hu_e  (:,:) =    hu_b(:,:)       
416         hv_e  (:,:) =    hv_b(:,:) 
417         hur_e (:,:) = r1_hu_b(:,:)   
418         hvr_e (:,:) = r1_hv_b(:,:)
419      ENDIF
420      !
421      ! Initialize sums:
422      ua_b  (:,:) = 0._wp       ! After barotropic velocities (or transport if flux form)         
423      va_b  (:,:) = 0._wp
424      ssha  (:,:) = 0._wp       ! Sum for after averaged sea level
425      un_adv(:,:) = 0._wp       ! Sum for now transport issued from ts loop
426      vn_adv(:,:) = 0._wp
427      !
428      IF( ln_wd_dl ) THEN
429         zuwdmask(:,:) = 0._wp  ! set to zero for definiteness (not sure this is necessary)
430         zvwdmask(:,:) = 0._wp  !
431         zuwdav2 (:,:) = 0._wp 
432         zvwdav2 (:,:) = 0._wp   
433      END IF 
434
435      !                                             ! ==================== !
436      DO jn = 1, icycle                             !  sub-time-step loop  !
437         !                                          ! ==================== !
438         !
439         l_full_nf_update = jn == icycle   ! false: disable full North fold update (performances) for jn = 1 to icycle-1
440         !
441         !                    !==  Update the forcing ==! (BDY and tides)
442         !
443         IF( ln_bdy      .AND. ln_tide )   CALL bdy_dta_tides( kt, kit=jn, pt_offset= REAL(noffset+1,wp) )
444         ! Update tide potential at the beginning of current time substep
445         IF( ln_tide_pot .AND. ln_tide ) THEN
446            zt0substep = REAL(nsec_day, wp) - 0.5_wp*rdt + (jn + noffset - 1) * rdt / REAL(nn_baro, wp)
447            CALL upd_tide(zt0substep)
448         END IF
449         !
450         !                    !==  extrapolation at mid-step  ==!   (jn+1/2)
451         !
452         !                       !* Set extrapolation coefficients for predictor step:
453         IF ((jn<3).AND.ll_init) THEN      ! Forward           
454           za1 = 1._wp                                         
455           za2 = 0._wp                       
456           za3 = 0._wp                       
457         ELSE                              ! AB3-AM4 Coefficients: bet=0.281105
458           za1 =  1.781105_wp              ! za1 =   3/2 +   bet
459           za2 = -1.06221_wp               ! za2 = -(1/2 + 2*bet)
460           za3 =  0.281105_wp              ! za3 = bet
461         ENDIF
462         !
463         !                       !* Extrapolate barotropic velocities at mid-step (jn+1/2)
464         !--        m+1/2               m                m-1           m-2       --!
465         !--       u      = (3/2+beta) u   -(1/2+2beta) u      + beta u          --!
466         !-------------------------------------------------------------------------!
467         ua_e(:,:) = za1 * un_e(:,:) + za2 * ub_e(:,:) + za3 * ubb_e(:,:)
468         va_e(:,:) = za1 * vn_e(:,:) + za2 * vb_e(:,:) + za3 * vbb_e(:,:)
469
470         IF( .NOT.ln_linssh ) THEN                        !* Update ocean depth (variable volume case only)
471            !                                             !  ------------------
472            ! Extrapolate Sea Level at step jit+0.5:
473            !--         m+1/2                 m                  m-1             m-2       --!
474            !--      ssh      = (3/2+beta) ssh   -(1/2+2beta) ssh      + beta ssh          --!
475            !--------------------------------------------------------------------------------!
476            zsshp2_e(:,:) = za1 * sshn_e(:,:)  + za2 * sshb_e(:,:) + za3 * sshbb_e(:,:)
477           
478            ! set wetting & drying mask at tracer points for this barotropic mid-step
479            IF( ln_wd_dl )   CALL wad_tmsk( zsshp2_e, ztwdmask )
480            !
481            !                          ! ocean t-depth at mid-step
482            zhtp2_e(:,:) = ht_0(:,:) + zsshp2_e(:,:)
483            !
484            !                          ! ocean u- and v-depth at mid-step   (separate DO-loops remove the need of a lbc_lnk)
485            DO jj = 1, jpj
486               DO ji = 1, jpim1   ! not jpi-column
487                  zhup2_e(ji,jj) = hu_0(ji,jj) + r1_2 * r1_e1e2u(ji,jj)                        &
488                       &                              * (  e1e2t(ji  ,jj) * zsshp2_e(ji  ,jj)  &
489                       &                                 + e1e2t(ji+1,jj) * zsshp2_e(ji+1,jj)  ) * ssumask(ji,jj)
490               END DO
491            END DO
492            DO jj = 1, jpjm1        ! not jpj-row
493               DO ji = 1, jpi
494                  zhvp2_e(ji,jj) = hv_0(ji,jj) + r1_2 * r1_e1e2v(ji,jj)                        &
495                       &                              * (  e1e2t(ji,jj  ) * zsshp2_e(ji,jj  )  &
496                       &                                 + e1e2t(ji,jj+1) * zsshp2_e(ji,jj+1)  ) * ssvmask(ji,jj)
497               END DO
498            END DO
499            !
500         ENDIF
501         !
502         !                    !==  after SSH  ==!   (jn+1)
503         !
504         !                             ! update (ua_e,va_e) to enforce volume conservation at open boundaries
505         !                             ! values of zhup2_e and zhvp2_e on the halo are not needed in bdy_vol2d
506         IF( ln_bdy .AND. ln_vol ) CALL bdy_vol2d( kt, jn, ua_e, va_e, zhup2_e, zhvp2_e )
507         !     
508         !                             ! resulting flux at mid-step (not over the full domain)
509         zhU(1:jpim1,1:jpj  ) = e2u(1:jpim1,1:jpj  ) * ua_e(1:jpim1,1:jpj  ) * zhup2_e(1:jpim1,1:jpj  )   ! not jpi-column
510         zhV(1:jpi  ,1:jpjm1) = e1v(1:jpi  ,1:jpjm1) * va_e(1:jpi  ,1:jpjm1) * zhvp2_e(1:jpi  ,1:jpjm1)   ! not jpj-row
511         !
512#if defined key_agrif
513         ! Set fluxes during predictor step to ensure volume conservation
514         IF( .NOT.Agrif_Root() .AND. ln_bt_fw ) CALL agrif_dyn_ts_flux( jn, zhU, zhV )
515#endif
516         IF( ln_wd_il )   CALL wad_lmt_bt(zhU, zhV, sshn_e, zssh_frc, rdtbt)    !!gm wad_lmt_bt use of lbc_lnk on zhU, zhV
517
518         IF( ln_wd_dl ) THEN           ! un_e and vn_e are set to zero at faces where
519            !                          ! the direction of the flow is from dry cells
520            CALL wad_Umsk( ztwdmask, zhU, zhV, un_e, vn_e, zuwdmask, zvwdmask )   ! not jpi colomn for U, not jpj row for V
521            !
522         ENDIF   
523         !
524         !
525         !     Compute Sea Level at step jit+1
526         !--           m+1        m                               m+1/2          --!
527         !--        ssh    =  ssh   - delta_t' * [ frc + div( flux      ) ]      --!
528         !-------------------------------------------------------------------------!
529         DO jj = 2, jpjm1        ! INNER domain                             
530            DO ji = 2, jpim1
531               zhdiv = (   zhU(ji,jj) - zhU(ji-1,jj) + zhV(ji,jj) - zhV(ji,jj-1)   ) * r1_e1e2t(ji,jj)
532               ssha_e(ji,jj) = (  sshn_e(ji,jj) - rdtbt * ( zssh_frc(ji,jj) + zhdiv )  ) * ssmask(ji,jj)
533            END DO
534         END DO
535         !
536         CALL lbc_lnk_multi( 'dynspg_ts', ssha_e, 'T', 1._wp,  zhU, 'U', -1._wp,  zhV, 'V', -1._wp )
537         !
538         !                             ! Sum over sub-time-steps to compute advective velocities
539         za2 = wgtbtp2(jn)             ! zhU, zhV hold fluxes extrapolated at jn+0.5
540         un_adv(:,:) = un_adv(:,:) + za2 * zhU(:,:) * r1_e2u(:,:)
541         vn_adv(:,:) = vn_adv(:,:) + za2 * zhV(:,:) * r1_e1v(:,:)
542         ! sum over sub-time-steps to decide which baroclinic velocities to set to zero (zuwdav2 is only used when ln_wd_dl_bc=True)
543         IF ( ln_wd_dl_bc ) THEN
544            zuwdav2(1:jpim1,1:jpj  ) = zuwdav2(1:jpim1,1:jpj  ) + za2 * zuwdmask(1:jpim1,1:jpj  )   ! not jpi-column
545            zvwdav2(1:jpi  ,1:jpjm1) = zvwdav2(1:jpi  ,1:jpjm1) + za2 * zvwdmask(1:jpi  ,1:jpjm1)   ! not jpj-row
546         END IF
547         !
548         ! Duplicate sea level across open boundaries (this is only cosmetic if linssh=T)
549         IF( ln_bdy )   CALL bdy_ssh( ssha_e )
550#if defined key_agrif
551         IF( .NOT.Agrif_Root() )   CALL agrif_ssh_ts( jn )
552#endif
553         
554         ! Sea Surface Height at u-,v-points (vvl case only)
555         IF( .NOT.ln_linssh ) THEN                               
556            DO jj = 2, jpjm1   ! INNER domain, will be extended to whole domain later
557               DO ji = 2, jpim1      ! NO Vector Opt.
558                  zsshu_a(ji,jj) = r1_2 * ssumask(ji,jj) * r1_e1e2u(ji,jj)    &
559                     &              * ( e1e2t(ji  ,jj  )  * ssha_e(ji  ,jj  ) &
560                     &              +   e1e2t(ji+1,jj  )  * ssha_e(ji+1,jj  ) )
561                  zsshv_a(ji,jj) = r1_2 * ssvmask(ji,jj) * r1_e1e2v(ji,jj)    &
562                     &              * ( e1e2t(ji  ,jj  )  * ssha_e(ji  ,jj  ) &
563                     &              +   e1e2t(ji  ,jj+1)  * ssha_e(ji  ,jj+1) )
564               END DO
565            END DO
566         ENDIF   
567         !         
568         ! Half-step back interpolation of SSH for surface pressure computation at step jit+1/2
569         !--            m+1/2           m+1              m               m-1              m-2     --!
570         !--        ssh'    =  za0 * ssh     +  za1 * ssh   +  za2 * ssh      +  za3 * ssh        --!
571         !------------------------------------------------------------------------------------------!
572         CALL ts_bck_interp( jn, ll_init, za0, za1, za2, za3 )   ! coeficients of the interpolation
573         zsshp2_e(:,:) = za0 *  ssha_e(:,:) + za1 *  sshn_e (:,:)   &
574            &          + za2 *  sshb_e(:,:) + za3 *  sshbb_e(:,:)
575         !
576         !                             ! Surface pressure gradient
577         zldg = ( 1._wp - rn_scal_load ) * grav    ! local factor
578         DO jj = 2, jpjm1                           
579            DO ji = 2, jpim1
580               zu_spg(ji,jj) = - zldg * ( zsshp2_e(ji+1,jj) - zsshp2_e(ji,jj) ) * r1_e1u(ji,jj)
581               zv_spg(ji,jj) = - zldg * ( zsshp2_e(ji,jj+1) - zsshp2_e(ji,jj) ) * r1_e2v(ji,jj)
582            END DO
583         END DO
584         IF( ln_wd_il ) THEN        ! W/D : gravity filters applied on pressure gradient
585            CALL wad_spg( zsshp2_e, zcpx, zcpy )   ! Calculating W/D gravity filters
586            zu_spg(2:jpim1,2:jpjm1) = zu_spg(2:jpim1,2:jpjm1) * zcpx(2:jpim1,2:jpjm1)
587            zv_spg(2:jpim1,2:jpjm1) = zv_spg(2:jpim1,2:jpjm1) * zcpy(2:jpim1,2:jpjm1)
588         ENDIF
589         !
590         ! Add Coriolis trend:
591         ! zwz array below or triads normally depend on sea level with ln_linssh=F and should be updated
592         ! at each time step. We however keep them constant here for optimization.
593         ! Recall that zhU and zhV hold fluxes at jn+0.5 (extrapolated not backward interpolated)
594         CALL dyn_cor_2d( zhup2_e, zhvp2_e, ua_e, va_e, zhU, zhV,    zu_trd, zv_trd   )
595         !
596         ! Add tidal astronomical forcing if defined
597         IF ( ln_tide .AND. ln_tide_pot ) THEN
598            DO jj = 2, jpjm1
599               DO ji = fs_2, fs_jpim1   ! vector opt.
600                  zu_trd(ji,jj) = zu_trd(ji,jj) + grav * ( pot_astro(ji+1,jj) - pot_astro(ji,jj) ) * r1_e1u(ji,jj)
601                  zv_trd(ji,jj) = zv_trd(ji,jj) + grav * ( pot_astro(ji,jj+1) - pot_astro(ji,jj) ) * r1_e2v(ji,jj)
602               END DO
603            END DO
604         ENDIF
605         !
606         ! Add bottom stresses:
607!jth do implicitly instead
608         IF ( .NOT. ll_wd ) THEN ! Revert to explicit for bit comparison tests in non wad runs
609            DO jj = 2, jpjm1
610               DO ji = fs_2, fs_jpim1   ! vector opt.
611                  zu_trd(ji,jj) = zu_trd(ji,jj) + zCdU_u(ji,jj) * un_e(ji,jj) * hur_e(ji,jj)
612                  zv_trd(ji,jj) = zv_trd(ji,jj) + zCdU_v(ji,jj) * vn_e(ji,jj) * hvr_e(ji,jj)
613               END DO
614            END DO
615         ENDIF
616         !
617         ! Set next velocities:
618         !     Compute barotropic speeds at step jit+1    (h : total height of the water colomn)
619         !--                              VECTOR FORM
620         !--   m+1                 m               /                                                       m+1/2           \    --!
621         !--  u     =             u   + delta_t' * \         (1-r)*g * grad_x( ssh') -         f * k vect u      +     frc /    --!
622         !--                                                                                                                    --!
623         !--                             FLUX FORM                                                                              --!
624         !--  m+1   __1__  /  m    m               /  m+1/2                             m+1/2              m+1/2    n      \ \  --!
625         !-- u    =   m+1 |  h  * u   + delta_t' * \ h     * (1-r)*g * grad_x( ssh') - h     * f * k vect u      + h * frc /  | --!
626         !--         h     \                                                                                                 /  --!
627         !------------------------------------------------------------------------------------------------------------------------!
628         IF( ln_dynadv_vec .OR. ln_linssh ) THEN      !* Vector form
629            DO jj = 2, jpjm1
630               DO ji = fs_2, fs_jpim1   ! vector opt.
631                  ua_e(ji,jj) = (                                 un_e(ji,jj)   & 
632                            &     + rdtbt * (                   zu_spg(ji,jj)   &
633                            &                                 + zu_trd(ji,jj)   &
634                            &                                 + zu_frc(ji,jj) ) & 
635                            &   ) * ssumask(ji,jj)
636
637                  va_e(ji,jj) = (                                 vn_e(ji,jj)   &
638                            &     + rdtbt * (                   zv_spg(ji,jj)   &
639                            &                                 + zv_trd(ji,jj)   &
640                            &                                 + zv_frc(ji,jj) ) &
641                            &   ) * ssvmask(ji,jj)
642               END DO
643            END DO
644            !
645         ELSE                           !* Flux form
646            DO jj = 2, jpjm1
647               DO ji = 2, jpim1
648                  !                    ! hu_e, hv_e hold depth at jn,  zhup2_e, zhvp2_e hold extrapolated depth at jn+1/2
649                  !                    ! backward interpolated depth used in spg terms at jn+1/2
650                  zhu_bck = hu_0(ji,jj) + r1_2*r1_e1e2u(ji,jj) * (  e1e2t(ji  ,jj) * zsshp2_e(ji  ,jj)    &
651                       &                                          + e1e2t(ji+1,jj) * zsshp2_e(ji+1,jj)  ) * ssumask(ji,jj)
652                  zhv_bck = hv_0(ji,jj) + r1_2*r1_e1e2v(ji,jj) * (  e1e2t(ji,jj  ) * zsshp2_e(ji,jj  )    &
653                       &                                          + e1e2t(ji,jj+1) * zsshp2_e(ji,jj+1)  ) * ssvmask(ji,jj)
654                  !                    ! inverse depth at jn+1
655                  z1_hu = ssumask(ji,jj) / ( hu_0(ji,jj) + zsshu_a(ji,jj) + 1._wp - ssumask(ji,jj) )
656                  z1_hv = ssvmask(ji,jj) / ( hv_0(ji,jj) + zsshv_a(ji,jj) + 1._wp - ssvmask(ji,jj) )
657                  !
658                  ua_e(ji,jj) = (               hu_e  (ji,jj) *   un_e (ji,jj)      & 
659                       &            + rdtbt * (  zhu_bck        * zu_spg (ji,jj)  &   !
660                       &                       + zhup2_e(ji,jj) * zu_trd (ji,jj)  &   !
661                       &                       +  hu_n  (ji,jj) * zu_frc (ji,jj)  )   ) * z1_hu
662                  !
663                  va_e(ji,jj) = (               hv_e  (ji,jj) *   vn_e (ji,jj)      &
664                       &            + rdtbt * (  zhv_bck        * zv_spg (ji,jj)  &   !
665                       &                       + zhvp2_e(ji,jj) * zv_trd (ji,jj)  &   !
666                       &                       +  hv_n  (ji,jj) * zv_frc (ji,jj)  )   ) * z1_hv
667               END DO
668            END DO
669         ENDIF
670!jth implicit bottom friction:
671         IF ( ll_wd ) THEN ! revert to explicit for bit comparison tests in non wad runs
672            DO jj = 2, jpjm1
673               DO ji = fs_2, fs_jpim1   ! vector opt.
674                     ua_e(ji,jj) =  ua_e(ji,jj) /(1.0 -   rdtbt * zCdU_u(ji,jj) * hur_e(ji,jj))
675                     va_e(ji,jj) =  va_e(ji,jj) /(1.0 -   rdtbt * zCdU_v(ji,jj) * hvr_e(ji,jj))
676               END DO
677            END DO
678         ENDIF
679       
680         IF( .NOT.ln_linssh ) THEN   !* Update ocean depth (variable volume case only)
681            hu_e (2:jpim1,2:jpjm1) = hu_0(2:jpim1,2:jpjm1) + zsshu_a(2:jpim1,2:jpjm1)
682            hur_e(2:jpim1,2:jpjm1) = ssumask(2:jpim1,2:jpjm1) / ( hu_e(2:jpim1,2:jpjm1) + 1._wp - ssumask(2:jpim1,2:jpjm1) )
683            hv_e (2:jpim1,2:jpjm1) = hv_0(2:jpim1,2:jpjm1) + zsshv_a(2:jpim1,2:jpjm1)
684            hvr_e(2:jpim1,2:jpjm1) = ssvmask(2:jpim1,2:jpjm1) / ( hv_e(2:jpim1,2:jpjm1) + 1._wp - ssvmask(2:jpim1,2:jpjm1) )
685            CALL lbc_lnk_multi( 'dynspg_ts', ua_e , 'U', -1._wp, va_e , 'V', -1._wp  &
686                 &                         , hu_e , 'U',  1._wp, hv_e , 'V',  1._wp  &
687                 &                         , hur_e, 'U',  1._wp, hvr_e, 'V',  1._wp  )
688         ELSE
689            CALL lbc_lnk_multi( 'dynspg_ts', ua_e , 'U', -1._wp, va_e , 'V', -1._wp  )
690         ENDIF
691         !
692         !
693         !                                                 ! open boundaries
694         IF( ln_bdy )   CALL bdy_dyn2d( jn, ua_e, va_e, un_e, vn_e, hur_e, hvr_e, ssha_e )
695#if defined key_agrif                                                           
696         IF( .NOT.Agrif_Root() )  CALL agrif_dyn_ts( jn )  ! Agrif
697#endif
698         !                                             !* Swap
699         !                                             !  ----
700         ubb_e  (:,:) = ub_e  (:,:)
701         ub_e   (:,:) = un_e  (:,:)
702         un_e   (:,:) = ua_e  (:,:)
703         !
704         vbb_e  (:,:) = vb_e  (:,:)
705         vb_e   (:,:) = vn_e  (:,:)
706         vn_e   (:,:) = va_e  (:,:)
707         !
708         sshbb_e(:,:) = sshb_e(:,:)
709         sshb_e (:,:) = sshn_e(:,:)
710         sshn_e (:,:) = ssha_e(:,:)
711
712         !                                             !* Sum over whole bt loop
713         !                                             !  ----------------------
714         za1 = wgtbtp1(jn)                                   
715         IF( ln_dynadv_vec .OR. ln_linssh ) THEN    ! Sum velocities
716            ua_b  (:,:) = ua_b  (:,:) + za1 * ua_e  (:,:) 
717            va_b  (:,:) = va_b  (:,:) + za1 * va_e  (:,:) 
718         ELSE                                       ! Sum transports
719            IF ( .NOT.ln_wd_dl ) THEN 
720               ua_b  (:,:) = ua_b  (:,:) + za1 * ua_e  (:,:) * hu_e (:,:)
721               va_b  (:,:) = va_b  (:,:) + za1 * va_e  (:,:) * hv_e (:,:)
722            ELSE
723               ua_b  (:,:) = ua_b  (:,:) + za1 * ua_e  (:,:) * hu_e (:,:) * zuwdmask(:,:)
724               va_b  (:,:) = va_b  (:,:) + za1 * va_e  (:,:) * hv_e (:,:) * zvwdmask(:,:)
725            END IF
726         ENDIF
727         !                                          ! Sum sea level
728         ssha(:,:) = ssha(:,:) + za1 * ssha_e(:,:)
729
730         !                                                 ! ==================== !
731      END DO                                               !        end loop      !
732      !                                                    ! ==================== !
733      ! -----------------------------------------------------------------------------
734      ! Phase 3. update the general trend with the barotropic trend
735      ! -----------------------------------------------------------------------------
736      !
737      ! Set advection velocity correction:
738      IF (ln_bt_fw) THEN
739         IF( .NOT.( kt == nit000 .AND. neuler==0 ) ) THEN
740            DO jj = 1, jpj
741               DO ji = 1, jpi
742                  zun_save = un_adv(ji,jj)
743                  zvn_save = vn_adv(ji,jj)
744                  !                          ! apply the previously computed correction
745                  un_adv(ji,jj) = r1_2 * ( ub2_b(ji,jj) + zun_save - atfp * un_bf(ji,jj) )
746                  vn_adv(ji,jj) = r1_2 * ( vb2_b(ji,jj) + zvn_save - atfp * vn_bf(ji,jj) )
747                  !                          ! Update corrective fluxes for next time step
748                  un_bf(ji,jj)  = atfp * un_bf(ji,jj) + ( zun_save - ub2_b(ji,jj) )
749                  vn_bf(ji,jj)  = atfp * vn_bf(ji,jj) + ( zvn_save - vb2_b(ji,jj) )
750                  !                          ! Save integrated transport for next computation
751                  ub2_b(ji,jj) = zun_save
752                  vb2_b(ji,jj) = zvn_save
753               END DO
754            END DO
755         ELSE
756            un_bf(:,:) = 0._wp            ! corrective fluxes for next time step set to zero
757            vn_bf(:,:) = 0._wp
758            ub2_b(:,:) = un_adv(:,:)      ! Save integrated transport for next computation
759            vb2_b(:,:) = vn_adv(:,:)
760         END IF
761      ENDIF
762
763
764      !
765      ! Update barotropic trend:
766      IF( ln_dynadv_vec .OR. ln_linssh ) THEN
767         DO jk=1,jpkm1
768            ua(:,:,jk) = ua(:,:,jk) + ( ua_b(:,:) - ub_b(:,:) ) * r1_2dt_b
769            va(:,:,jk) = va(:,:,jk) + ( va_b(:,:) - vb_b(:,:) ) * r1_2dt_b
770         END DO
771      ELSE
772         ! At this stage, ssha has been corrected: compute new depths at velocity points
773         DO jj = 1, jpjm1
774            DO ji = 1, jpim1      ! NO Vector Opt.
775               zsshu_a(ji,jj) = r1_2 * ssumask(ji,jj)  * r1_e1e2u(ji,jj) &
776                  &              * ( e1e2t(ji  ,jj) * ssha(ji  ,jj)      &
777                  &              +   e1e2t(ji+1,jj) * ssha(ji+1,jj) )
778               zsshv_a(ji,jj) = r1_2 * ssvmask(ji,jj)  * r1_e1e2v(ji,jj) &
779                  &              * ( e1e2t(ji,jj  ) * ssha(ji,jj  )      &
780                  &              +   e1e2t(ji,jj+1) * ssha(ji,jj+1) )
781            END DO
782         END DO
783         CALL lbc_lnk_multi( 'dynspg_ts', zsshu_a, 'U', 1._wp, zsshv_a, 'V', 1._wp ) ! Boundary conditions
784         !
785         DO jk=1,jpkm1
786            ua(:,:,jk) = ua(:,:,jk) + r1_hu_n(:,:) * ( ua_b(:,:) - ub_b(:,:) * hu_b(:,:) ) * r1_2dt_b
787            va(:,:,jk) = va(:,:,jk) + r1_hv_n(:,:) * ( va_b(:,:) - vb_b(:,:) * hv_b(:,:) ) * r1_2dt_b
788         END DO
789         ! Save barotropic velocities not transport:
790         ua_b(:,:) =  ua_b(:,:) / ( hu_0(:,:) + zsshu_a(:,:) + 1._wp - ssumask(:,:) )
791         va_b(:,:) =  va_b(:,:) / ( hv_0(:,:) + zsshv_a(:,:) + 1._wp - ssvmask(:,:) )
792      ENDIF
793
794
795      ! Correct velocities so that the barotropic velocity equals (un_adv, vn_adv) (in all cases) 
796      DO jk = 1, jpkm1
797         un(:,:,jk) = ( un(:,:,jk) + un_adv(:,:)*r1_hu_n(:,:) - un_b(:,:) ) * umask(:,:,jk)
798         vn(:,:,jk) = ( vn(:,:,jk) + vn_adv(:,:)*r1_hv_n(:,:) - vn_b(:,:) ) * vmask(:,:,jk)
799      END DO
800
801      IF ( ln_wd_dl .and. ln_wd_dl_bc) THEN 
802         ! need to set lbc here because not done prior time averaging
803         CALL lbc_lnk_multi( 'dynspg_ts', zuwdav2, 'U', 1._wp, zvwdav2, 'V', 1._wp)
804         DO jk = 1, jpkm1
805            un(:,:,jk) = ( un_adv(:,:)*r1_hu_n(:,:) &
806                       & + zuwdav2(:,:)*(un(:,:,jk) - un_adv(:,:)*r1_hu_n(:,:)) ) * umask(:,:,jk) 
807            vn(:,:,jk) = ( vn_adv(:,:)*r1_hv_n(:,:) & 
808                       & + zvwdav2(:,:)*(vn(:,:,jk) - vn_adv(:,:)*r1_hv_n(:,:)) ) * vmask(:,:,jk) 
809         END DO
810      END IF
811
812     
813      CALL iom_put(  "ubar", un_adv(:,:)*r1_hu_n(:,:) )    ! barotropic i-current
814      CALL iom_put(  "vbar", vn_adv(:,:)*r1_hv_n(:,:) )    ! barotropic i-current
815      !
816#if defined key_agrif
817      ! Save time integrated fluxes during child grid integration
818      ! (used to update coarse grid transports at next time step)
819      !
820      IF( .NOT.Agrif_Root() .AND. ln_bt_fw ) THEN
821         IF( Agrif_NbStepint() == 0 ) THEN
822            ub2_i_b(:,:) = 0._wp
823            vb2_i_b(:,:) = 0._wp
824         END IF
825         !
826         za1 = 1._wp / REAL(Agrif_rhot(), wp)
827         ub2_i_b(:,:) = ub2_i_b(:,:) + za1 * ub2_b(:,:)
828         vb2_i_b(:,:) = vb2_i_b(:,:) + za1 * vb2_b(:,:)
829      ENDIF
830#endif     
831      !                                   !* write time-spliting arrays in the restart
832      IF( lrst_oce .AND.ln_bt_fw )   CALL ts_rst( kt, 'WRITE' )
833      !
834      IF( ln_wd_il )   DEALLOCATE( zcpx, zcpy )
835      IF( ln_wd_dl )   DEALLOCATE( ztwdmask, zuwdmask, zvwdmask, zuwdav2, zvwdav2 )
836      !
837      IF( ln_diatmb ) THEN
838         CALL iom_put( "baro_u" , un_b*ssumask(:,:)+zmdi*(1.-ssumask(:,:) ) )  ! Barotropic  U Velocity
839         CALL iom_put( "baro_v" , vn_b*ssvmask(:,:)+zmdi*(1.-ssvmask(:,:) ) )  ! Barotropic  V Velocity
840      ENDIF
841      !
842   END SUBROUTINE dyn_spg_ts
843
844
845   SUBROUTINE ts_wgt( ll_av, ll_fw, jpit, zwgt1, zwgt2)
846      !!---------------------------------------------------------------------
847      !!                   ***  ROUTINE ts_wgt  ***
848      !!
849      !! ** Purpose : Set time-splitting weights for temporal averaging (or not)
850      !!----------------------------------------------------------------------
851      LOGICAL, INTENT(in) ::   ll_av      ! temporal averaging=.true.
852      LOGICAL, INTENT(in) ::   ll_fw      ! forward time splitting =.true.
853      INTEGER, INTENT(inout) :: jpit      ! cycle length   
854      REAL(wp), DIMENSION(3*nn_baro), INTENT(inout) ::   zwgt1, & ! Primary weights
855                                                         zwgt2    ! Secondary weights
856     
857      INTEGER ::  jic, jn, ji                      ! temporary integers
858      REAL(wp) :: za1, za2
859      !!----------------------------------------------------------------------
860
861      zwgt1(:) = 0._wp
862      zwgt2(:) = 0._wp
863
864      ! Set time index when averaged value is requested
865      IF (ll_fw) THEN
866         jic = nn_baro
867      ELSE
868         jic = 2 * nn_baro
869      ENDIF
870
871      ! Set primary weights:
872      IF (ll_av) THEN
873           ! Define simple boxcar window for primary weights
874           ! (width = nn_baro, centered around jic)     
875         SELECT CASE ( nn_bt_flt )
876              CASE( 0 )  ! No averaging
877                 zwgt1(jic) = 1._wp
878                 jpit = jic
879
880              CASE( 1 )  ! Boxcar, width = nn_baro
881                 DO jn = 1, 3*nn_baro
882                    za1 = ABS(float(jn-jic))/float(nn_baro) 
883                    IF (za1 < 0.5_wp) THEN
884                      zwgt1(jn) = 1._wp
885                      jpit = jn
886                    ENDIF
887                 ENDDO
888
889              CASE( 2 )  ! Boxcar, width = 2 * nn_baro
890                 DO jn = 1, 3*nn_baro
891                    za1 = ABS(float(jn-jic))/float(nn_baro) 
892                    IF (za1 < 1._wp) THEN
893                      zwgt1(jn) = 1._wp
894                      jpit = jn
895                    ENDIF
896                 ENDDO
897              CASE DEFAULT   ;   CALL ctl_stop( 'unrecognised value for nn_bt_flt' )
898         END SELECT
899
900      ELSE ! No time averaging
901         zwgt1(jic) = 1._wp
902         jpit = jic
903      ENDIF
904   
905      ! Set secondary weights
906      DO jn = 1, jpit
907        DO ji = jn, jpit
908             zwgt2(jn) = zwgt2(jn) + zwgt1(ji)
909        END DO
910      END DO
911
912      ! Normalize weigths:
913      za1 = 1._wp / SUM(zwgt1(1:jpit))
914      za2 = 1._wp / SUM(zwgt2(1:jpit))
915      DO jn = 1, jpit
916        zwgt1(jn) = zwgt1(jn) * za1
917        zwgt2(jn) = zwgt2(jn) * za2
918      END DO
919      !
920   END SUBROUTINE ts_wgt
921
922
923   SUBROUTINE ts_rst( kt, cdrw )
924      !!---------------------------------------------------------------------
925      !!                   ***  ROUTINE ts_rst  ***
926      !!
927      !! ** Purpose : Read or write time-splitting arrays in restart file
928      !!----------------------------------------------------------------------
929      INTEGER         , INTENT(in) ::   kt     ! ocean time-step
930      CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
931      !!----------------------------------------------------------------------
932      !
933      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
934         !                                   ! ---------------
935         IF( ln_rstart .AND. ln_bt_fw .AND. (neuler/=0) ) THEN    !* Read the restart file
936            CALL iom_get( numror, jpdom_autoglo, 'ub2_b'  , ub2_b  (:,:), ldxios = lrxios )   
937            CALL iom_get( numror, jpdom_autoglo, 'vb2_b'  , vb2_b  (:,:), ldxios = lrxios ) 
938            CALL iom_get( numror, jpdom_autoglo, 'un_bf'  , un_bf  (:,:), ldxios = lrxios )   
939            CALL iom_get( numror, jpdom_autoglo, 'vn_bf'  , vn_bf  (:,:), ldxios = lrxios ) 
940            IF( .NOT.ln_bt_av ) THEN
941               CALL iom_get( numror, jpdom_autoglo, 'sshbb_e'  , sshbb_e(:,:), ldxios = lrxios )   
942               CALL iom_get( numror, jpdom_autoglo, 'ubb_e'    ,   ubb_e(:,:), ldxios = lrxios )   
943               CALL iom_get( numror, jpdom_autoglo, 'vbb_e'    ,   vbb_e(:,:), ldxios = lrxios )
944               CALL iom_get( numror, jpdom_autoglo, 'sshb_e'   ,  sshb_e(:,:), ldxios = lrxios ) 
945               CALL iom_get( numror, jpdom_autoglo, 'ub_e'     ,    ub_e(:,:), ldxios = lrxios )   
946               CALL iom_get( numror, jpdom_autoglo, 'vb_e'     ,    vb_e(:,:), ldxios = lrxios )
947            ENDIF
948#if defined key_agrif
949            ! Read time integrated fluxes
950            IF ( .NOT.Agrif_Root() ) THEN
951               CALL iom_get( numror, jpdom_autoglo, 'ub2_i_b'  , ub2_i_b(:,:), ldxios = lrxios )   
952               CALL iom_get( numror, jpdom_autoglo, 'vb2_i_b'  , vb2_i_b(:,:), ldxios = lrxios )
953            ENDIF
954#endif
955         ELSE                                   !* Start from rest
956            IF(lwp) WRITE(numout,*)
957            IF(lwp) WRITE(numout,*) '   ==>>>   start from rest: set barotropic values to 0'
958            ub2_b (:,:) = 0._wp   ;   vb2_b (:,:) = 0._wp   ! used in the 1st interpol of agrif
959            un_adv(:,:) = 0._wp   ;   vn_adv(:,:) = 0._wp   ! used in the 1st interpol of agrif
960            un_bf (:,:) = 0._wp   ;   vn_bf (:,:) = 0._wp   ! used in the 1st update   of agrif
961#if defined key_agrif
962            IF ( .NOT.Agrif_Root() ) THEN
963               ub2_i_b(:,:) = 0._wp   ;   vb2_i_b(:,:) = 0._wp   ! used in the 1st update of agrif
964            ENDIF
965#endif
966         ENDIF
967         !
968      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
969         !                                   ! -------------------
970         IF(lwp) WRITE(numout,*) '---- ts_rst ----'
971         IF( lwxios ) CALL iom_swap(      cwxios_context          )
972         CALL iom_rstput( kt, nitrst, numrow, 'ub2_b'   , ub2_b  (:,:), ldxios = lwxios )
973         CALL iom_rstput( kt, nitrst, numrow, 'vb2_b'   , vb2_b  (:,:), ldxios = lwxios )
974         CALL iom_rstput( kt, nitrst, numrow, 'un_bf'   , un_bf  (:,:), ldxios = lwxios )
975         CALL iom_rstput( kt, nitrst, numrow, 'vn_bf'   , vn_bf  (:,:), ldxios = lwxios )
976         !
977         IF (.NOT.ln_bt_av) THEN
978            CALL iom_rstput( kt, nitrst, numrow, 'sshbb_e'  , sshbb_e(:,:), ldxios = lwxios ) 
979            CALL iom_rstput( kt, nitrst, numrow, 'ubb_e'    ,   ubb_e(:,:), ldxios = lwxios )
980            CALL iom_rstput( kt, nitrst, numrow, 'vbb_e'    ,   vbb_e(:,:), ldxios = lwxios )
981            CALL iom_rstput( kt, nitrst, numrow, 'sshb_e'   ,  sshb_e(:,:), ldxios = lwxios )
982            CALL iom_rstput( kt, nitrst, numrow, 'ub_e'     ,    ub_e(:,:), ldxios = lwxios )
983            CALL iom_rstput( kt, nitrst, numrow, 'vb_e'     ,    vb_e(:,:), ldxios = lwxios )
984         ENDIF
985#if defined key_agrif
986         ! Save time integrated fluxes
987         IF ( .NOT.Agrif_Root() ) THEN
988            CALL iom_rstput( kt, nitrst, numrow, 'ub2_i_b'  , ub2_i_b(:,:), ldxios = lwxios )
989            CALL iom_rstput( kt, nitrst, numrow, 'vb2_i_b'  , vb2_i_b(:,:), ldxios = lwxios )
990         ENDIF
991#endif
992         IF( lwxios ) CALL iom_swap(      cxios_context          )
993      ENDIF
994      !
995   END SUBROUTINE ts_rst
996
997
998   SUBROUTINE dyn_spg_ts_init
999      !!---------------------------------------------------------------------
1000      !!                   ***  ROUTINE dyn_spg_ts_init  ***
1001      !!
1002      !! ** Purpose : Set time splitting options
1003      !!----------------------------------------------------------------------
1004      INTEGER  ::   ji ,jj              ! dummy loop indices
1005      REAL(wp) ::   zxr2, zyr2, zcmax   ! local scalar
1006      REAL(wp), DIMENSION(jpi,jpj) ::   zcu
1007      INTEGER  :: inum
1008      !!----------------------------------------------------------------------
1009      !
1010      ! Max courant number for ext. grav. waves
1011      !
1012      DO jj = 1, jpj
1013         DO ji =1, jpi
1014            zxr2 = r1_e1t(ji,jj) * r1_e1t(ji,jj)
1015            zyr2 = r1_e2t(ji,jj) * r1_e2t(ji,jj)
1016            zcu(ji,jj) = SQRT( grav * MAX(ht_0(ji,jj),0._wp) * (zxr2 + zyr2) )
1017         END DO
1018      END DO
1019      !
1020      zcmax = MAXVAL( zcu(:,:) )
1021      CALL mpp_max( 'dynspg_ts', zcmax )
1022
1023      ! Estimate number of iterations to satisfy a max courant number= rn_bt_cmax
1024      IF( ln_bt_auto )   nn_baro = CEILING( rdt / rn_bt_cmax * zcmax)
1025     
1026      rdtbt = rdt / REAL( nn_baro , wp )
1027      zcmax = zcmax * rdtbt
1028      ! Print results
1029      IF(lwp) WRITE(numout,*)
1030      IF(lwp) WRITE(numout,*) 'dyn_spg_ts_init : split-explicit free surface'
1031      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~~'
1032      IF( ln_bt_auto ) THEN
1033         IF(lwp) WRITE(numout,*) '     ln_ts_auto =.true. Automatically set nn_baro '
1034         IF(lwp) WRITE(numout,*) '     Max. courant number allowed: ', rn_bt_cmax
1035      ELSE
1036         IF(lwp) WRITE(numout,*) '     ln_ts_auto=.false.: Use nn_baro in namelist   nn_baro = ', nn_baro
1037      ENDIF
1038
1039      IF(ln_bt_av) THEN
1040         IF(lwp) WRITE(numout,*) '     ln_bt_av =.true.  ==> Time averaging over nn_baro time steps is on '
1041      ELSE
1042         IF(lwp) WRITE(numout,*) '     ln_bt_av =.false. => No time averaging of barotropic variables '
1043      ENDIF
1044      !
1045      !
1046      IF(ln_bt_fw) THEN
1047         IF(lwp) WRITE(numout,*) '     ln_bt_fw=.true.  => Forward integration of barotropic variables '
1048      ELSE
1049         IF(lwp) WRITE(numout,*) '     ln_bt_fw =.false.=> Centred integration of barotropic variables '
1050      ENDIF
1051      !
1052#if defined key_agrif
1053      ! Restrict the use of Agrif to the forward case only
1054!!!      IF( .NOT.ln_bt_fw .AND. .NOT.Agrif_Root() )   CALL ctl_stop( 'AGRIF not implemented if ln_bt_fw=.FALSE.' )
1055#endif
1056      !
1057      IF(lwp) WRITE(numout,*)    '     Time filter choice, nn_bt_flt: ', nn_bt_flt
1058      SELECT CASE ( nn_bt_flt )
1059         CASE( 0 )      ;   IF(lwp) WRITE(numout,*) '           Dirac'
1060         CASE( 1 )      ;   IF(lwp) WRITE(numout,*) '           Boxcar: width = nn_baro'
1061         CASE( 2 )      ;   IF(lwp) WRITE(numout,*) '           Boxcar: width = 2*nn_baro' 
1062         CASE DEFAULT   ;   CALL ctl_stop( 'unrecognised value for nn_bt_flt: should 0,1, or 2' )
1063      END SELECT
1064      !
1065      IF(lwp) WRITE(numout,*) ' '
1066      IF(lwp) WRITE(numout,*) '     nn_baro = ', nn_baro
1067      IF(lwp) WRITE(numout,*) '     Barotropic time step [s] is :', rdtbt
1068      IF(lwp) WRITE(numout,*) '     Maximum Courant number is   :', zcmax
1069      !
1070      IF(lwp) WRITE(numout,*)    '     Time diffusion parameter rn_bt_alpha: ', rn_bt_alpha
1071      IF ((ln_bt_av.AND.nn_bt_flt/=0).AND.(rn_bt_alpha>0._wp)) THEN
1072         CALL ctl_stop( 'dynspg_ts ERROR: if rn_bt_alpha > 0, remove temporal averaging' )
1073      ENDIF
1074      !
1075      IF( .NOT.ln_bt_av .AND. .NOT.ln_bt_fw ) THEN
1076         CALL ctl_stop( 'dynspg_ts ERROR: No time averaging => only forward integration is possible' )
1077      ENDIF
1078      IF( zcmax>0.9_wp ) THEN
1079         CALL ctl_stop( 'dynspg_ts ERROR: Maximum Courant number is greater than 0.9: Inc. nn_baro !' )         
1080      ENDIF
1081      !
1082      !                             ! Allocate time-splitting arrays
1083      IF( dyn_spg_ts_alloc() /= 0    )   CALL ctl_stop('STOP', 'dyn_spg_init: failed to allocate dynspg_ts  arrays' )
1084      !
1085      !                             ! read restart when needed
1086      CALL ts_rst( nit000, 'READ' )
1087      !
1088      IF( lwxios ) THEN
1089! define variables in restart file when writing with XIOS
1090         CALL iom_set_rstw_var_active('ub2_b')
1091         CALL iom_set_rstw_var_active('vb2_b')
1092         CALL iom_set_rstw_var_active('un_bf')
1093         CALL iom_set_rstw_var_active('vn_bf')
1094         !
1095         IF (.NOT.ln_bt_av) THEN
1096            CALL iom_set_rstw_var_active('sshbb_e')
1097            CALL iom_set_rstw_var_active('ubb_e')
1098            CALL iom_set_rstw_var_active('vbb_e')
1099            CALL iom_set_rstw_var_active('sshb_e')
1100            CALL iom_set_rstw_var_active('ub_e')
1101            CALL iom_set_rstw_var_active('vb_e')
1102         ENDIF
1103#if defined key_agrif
1104         ! Save time integrated fluxes
1105         IF ( .NOT.Agrif_Root() ) THEN
1106            CALL iom_set_rstw_var_active('ub2_i_b')
1107            CALL iom_set_rstw_var_active('vb2_i_b')
1108         ENDIF
1109#endif
1110      ENDIF
1111      !
1112   END SUBROUTINE dyn_spg_ts_init
1113
1114   
1115   SUBROUTINE dyn_cor_2d_init
1116      !!---------------------------------------------------------------------
1117      !!                   ***  ROUTINE dyn_cor_2d_init  ***
1118      !!
1119      !! ** Purpose : Set time splitting options
1120      !! Set arrays to remove/compute coriolis trend.
1121      !! Do it once during initialization if volume is fixed, else at each long time step.
1122      !! Note that these arrays are also used during barotropic loop. These are however frozen
1123      !! although they should be updated in the variable volume case. Not a big approximation.
1124      !! To remove this approximation, copy lines below inside barotropic loop
1125      !! and update depths at T-F points (ht and zhf resp.) at each barotropic time step
1126      !!
1127      !! Compute zwz = f / ( height of the water colomn )
1128      !!----------------------------------------------------------------------
1129      INTEGER  ::   ji ,jj, jk              ! dummy loop indices
1130      REAL(wp) ::   z1_ht
1131      REAL(wp), DIMENSION(jpi,jpj) :: zhf
1132      !!----------------------------------------------------------------------
1133      !
1134      SELECT CASE( nvor_scheme )
1135      CASE( np_EEN )                != EEN scheme using e3f (energy & enstrophy scheme)
1136         SELECT CASE( nn_een_e3f )              !* ff_f/e3 at F-point
1137         CASE ( 0 )                                   ! original formulation  (masked averaging of e3t divided by 4)
1138            DO jj = 1, jpjm1
1139               DO ji = 1, jpim1
1140                  zwz(ji,jj) =   ( ht_n(ji  ,jj+1) + ht_n(ji+1,jj+1) +                    &
1141                       &             ht_n(ji  ,jj  ) + ht_n(ji+1,jj  )   ) * 0.25_wp 
1142                  IF( zwz(ji,jj) /= 0._wp )   zwz(ji,jj) = ff_f(ji,jj) / zwz(ji,jj)
1143               END DO
1144            END DO
1145         CASE ( 1 )                                   ! new formulation  (masked averaging of e3t divided by the sum of mask)
1146            DO jj = 1, jpjm1
1147               DO ji = 1, jpim1
1148                  zwz(ji,jj) =             (  ht_n  (ji  ,jj+1) + ht_n  (ji+1,jj+1)      &
1149                       &                      + ht_n  (ji  ,jj  ) + ht_n  (ji+1,jj  )  )   &
1150                       &       / ( MAX( 1._wp,  ssmask(ji  ,jj+1) + ssmask(ji+1,jj+1)      &
1151                       &                      + ssmask(ji  ,jj  ) + ssmask(ji+1,jj  )  )   )
1152                  IF( zwz(ji,jj) /= 0._wp )   zwz(ji,jj) = ff_f(ji,jj) / zwz(ji,jj)
1153               END DO
1154            END DO
1155         END SELECT
1156         CALL lbc_lnk( 'dynspg_ts', zwz, 'F', 1._wp )
1157         !
1158         ftne(1,:) = 0._wp ; ftnw(1,:) = 0._wp ; ftse(1,:) = 0._wp ; ftsw(1,:) = 0._wp
1159         DO jj = 2, jpj
1160            DO ji = 2, jpi
1161               ftne(ji,jj) = zwz(ji-1,jj  ) + zwz(ji  ,jj  ) + zwz(ji  ,jj-1)
1162               ftnw(ji,jj) = zwz(ji-1,jj-1) + zwz(ji-1,jj  ) + zwz(ji  ,jj  )
1163               ftse(ji,jj) = zwz(ji  ,jj  ) + zwz(ji  ,jj-1) + zwz(ji-1,jj-1)
1164               ftsw(ji,jj) = zwz(ji  ,jj-1) + zwz(ji-1,jj-1) + zwz(ji-1,jj  )
1165            END DO
1166         END DO
1167         !
1168      CASE( np_EET )                  != EEN scheme using e3t (energy conserving scheme)
1169         ftne(1,:) = 0._wp ; ftnw(1,:) = 0._wp ; ftse(1,:) = 0._wp ; ftsw(1,:) = 0._wp
1170         DO jj = 2, jpj
1171            DO ji = 2, jpi
1172               z1_ht = ssmask(ji,jj) / ( ht_n(ji,jj) + 1._wp - ssmask(ji,jj) )
1173               ftne(ji,jj) = ( ff_f(ji-1,jj  ) + ff_f(ji  ,jj  ) + ff_f(ji  ,jj-1) ) * z1_ht
1174               ftnw(ji,jj) = ( ff_f(ji-1,jj-1) + ff_f(ji-1,jj  ) + ff_f(ji  ,jj  ) ) * z1_ht
1175               ftse(ji,jj) = ( ff_f(ji  ,jj  ) + ff_f(ji  ,jj-1) + ff_f(ji-1,jj-1) ) * z1_ht
1176               ftsw(ji,jj) = ( ff_f(ji  ,jj-1) + ff_f(ji-1,jj-1) + ff_f(ji-1,jj  ) ) * z1_ht
1177            END DO
1178         END DO
1179         !
1180      CASE( np_ENE, np_ENS , np_MIX )  != all other schemes (ENE, ENS, MIX) except ENT !
1181         !
1182         zwz(:,:) = 0._wp
1183         zhf(:,:) = 0._wp
1184         
1185         !!gm  assume 0 in both cases (which is almost surely WRONG ! ) as hvatf has been removed
1186!!gm    A priori a better value should be something like :
1187!!gm          zhf(i,j) = masked sum of  ht(i,j) , ht(i+1,j) , ht(i,j+1) , (i+1,j+1)
1188!!gm                     divided by the sum of the corresponding mask
1189!!gm
1190!!           
1191         IF( .NOT.ln_sco ) THEN
1192 
1193   !!gm  agree the JC comment  : this should be done in a much clear way
1194 
1195   ! JC: It not clear yet what should be the depth at f-points over land in z-coordinate case
1196   !     Set it to zero for the time being
1197   !              IF( rn_hmin < 0._wp ) THEN    ;   jk = - INT( rn_hmin )                                      ! from a nb of level
1198   !              ELSE                          ;   jk = MINLOC( gdepw_0, mask = gdepw_0 > rn_hmin, dim = 1 )  ! from a depth
1199   !              ENDIF
1200   !              zhf(:,:) = gdepw_0(:,:,jk+1)
1201            !
1202         ELSE
1203            !
1204            !zhf(:,:) = hbatf(:,:)
1205            DO jj = 1, jpjm1
1206               DO ji = 1, jpim1
1207                  zhf(ji,jj) =    (   ht_0  (ji,jj  ) + ht_0  (ji+1,jj  )          &
1208                       &            + ht_0  (ji,jj+1) + ht_0  (ji+1,jj+1)   )      &
1209                       &     / MAX(   ssmask(ji,jj  ) + ssmask(ji+1,jj  )          &
1210                       &            + ssmask(ji,jj+1) + ssmask(ji+1,jj+1) , 1._wp  )
1211               END DO
1212            END DO
1213         ENDIF
1214         !
1215         DO jj = 1, jpjm1
1216            zhf(:,jj) = zhf(:,jj) * (1._wp- umask(:,jj,1) * umask(:,jj+1,1))
1217         END DO
1218         !
1219         DO jk = 1, jpkm1
1220            DO jj = 1, jpjm1
1221               zhf(:,jj) = zhf(:,jj) + e3f_n(:,jj,jk) * umask(:,jj,jk) * umask(:,jj+1,jk)
1222            END DO
1223         END DO
1224         CALL lbc_lnk( 'dynspg_ts', zhf, 'F', 1._wp )
1225         ! JC: TBC. hf should be greater than 0
1226         DO jj = 1, jpj
1227            DO ji = 1, jpi
1228               IF( zhf(ji,jj) /= 0._wp )   zwz(ji,jj) = 1._wp / zhf(ji,jj)
1229            END DO
1230         END DO
1231         zwz(:,:) = ff_f(:,:) * zwz(:,:)
1232      END SELECT
1233     
1234   END SUBROUTINE dyn_cor_2d_init
1235
1236
1237
1238   SUBROUTINE dyn_cor_2d( hu_n, hv_n, un_b, vn_b, zhU, zhV,    zu_trd, zv_trd   )
1239      !!---------------------------------------------------------------------
1240      !!                   ***  ROUTINE dyn_cor_2d  ***
1241      !!
1242      !! ** Purpose : Compute u and v coriolis trends
1243      !!----------------------------------------------------------------------
1244      INTEGER  ::   ji ,jj                             ! dummy loop indices
1245      REAL(wp) ::   zx1, zx2, zy1, zy2, z1_hu, z1_hv   !   -      -
1246      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ) :: hu_n, hv_n, un_b, vn_b, zhU, zhV
1247      REAL(wp), DIMENSION(jpi,jpj), INTENT(  out) :: zu_trd, zv_trd
1248      !!----------------------------------------------------------------------
1249      SELECT CASE( nvor_scheme )
1250      CASE( np_ENT )                ! enstrophy conserving scheme (f-point)
1251         DO jj = 2, jpjm1
1252            DO ji = 2, jpim1
1253               z1_hu = ssumask(ji,jj) / ( hu_n(ji,jj) + 1._wp - ssumask(ji,jj) )
1254               z1_hv = ssvmask(ji,jj) / ( hv_n(ji,jj) + 1._wp - ssvmask(ji,jj) )
1255               zu_trd(ji,jj) = + r1_4 * r1_e1e2u(ji,jj) * z1_hu                    &
1256                  &               * (  e1e2t(ji+1,jj)*ht_n(ji+1,jj)*ff_t(ji+1,jj) * ( vn_b(ji+1,jj) + vn_b(ji+1,jj-1) )   &
1257                  &                  + e1e2t(ji  ,jj)*ht_n(ji  ,jj)*ff_t(ji  ,jj) * ( vn_b(ji  ,jj) + vn_b(ji  ,jj-1) )   )
1258                  !
1259               zv_trd(ji,jj) = - r1_4 * r1_e1e2v(ji,jj) * z1_hv                    &
1260                  &               * (  e1e2t(ji,jj+1)*ht_n(ji,jj+1)*ff_t(ji,jj+1) * ( un_b(ji,jj+1) + un_b(ji-1,jj+1) )   & 
1261                  &                  + e1e2t(ji,jj  )*ht_n(ji,jj  )*ff_t(ji,jj  ) * ( un_b(ji,jj  ) + un_b(ji-1,jj  ) )   ) 
1262            END DO 
1263         END DO 
1264         !         
1265      CASE( np_ENE , np_MIX )        ! energy conserving scheme (t-point) ENE or MIX
1266         DO jj = 2, jpjm1
1267            DO ji = fs_2, fs_jpim1   ! vector opt.
1268               zy1 = ( zhV(ji,jj-1) + zhV(ji+1,jj-1) ) * r1_e1u(ji,jj)
1269               zy2 = ( zhV(ji,jj  ) + zhV(ji+1,jj  ) ) * r1_e1u(ji,jj)
1270               zx1 = ( zhU(ji-1,jj) + zhU(ji-1,jj+1) ) * r1_e2v(ji,jj)
1271               zx2 = ( zhU(ji  ,jj) + zhU(ji  ,jj+1) ) * r1_e2v(ji,jj)
1272               ! energy conserving formulation for planetary vorticity term
1273               zu_trd(ji,jj) =   r1_4 * ( zwz(ji  ,jj-1) * zy1 + zwz(ji,jj) * zy2 )
1274               zv_trd(ji,jj) = - r1_4 * ( zwz(ji-1,jj  ) * zx1 + zwz(ji,jj) * zx2 )
1275            END DO
1276         END DO
1277         !
1278      CASE( np_ENS )                ! enstrophy conserving scheme (f-point)
1279         DO jj = 2, jpjm1
1280            DO ji = fs_2, fs_jpim1   ! vector opt.
1281               zy1 =   r1_8 * ( zhV(ji  ,jj-1) + zhV(ji+1,jj-1) &
1282                 &            + zhV(ji  ,jj  ) + zhV(ji+1,jj  ) ) * r1_e1u(ji,jj)
1283               zx1 = - r1_8 * ( zhU(ji-1,jj  ) + zhU(ji-1,jj+1) &
1284                 &            + zhU(ji  ,jj  ) + zhU(ji  ,jj+1) ) * r1_e2v(ji,jj)
1285               zu_trd(ji,jj)  = zy1 * ( zwz(ji  ,jj-1) + zwz(ji,jj) )
1286               zv_trd(ji,jj)  = zx1 * ( zwz(ji-1,jj  ) + zwz(ji,jj) )
1287            END DO
1288         END DO
1289         !
1290      CASE( np_EET , np_EEN )      ! energy & enstrophy scheme (using e3t or e3f)         
1291         DO jj = 2, jpjm1
1292            DO ji = fs_2, fs_jpim1   ! vector opt.
1293               zu_trd(ji,jj) = + r1_12 * r1_e1u(ji,jj) * (  ftne(ji,jj  ) * zhV(ji  ,jj  ) &
1294                &                                         + ftnw(ji+1,jj) * zhV(ji+1,jj  ) &
1295                &                                         + ftse(ji,jj  ) * zhV(ji  ,jj-1) &
1296                &                                         + ftsw(ji+1,jj) * zhV(ji+1,jj-1) )
1297               zv_trd(ji,jj) = - r1_12 * r1_e2v(ji,jj) * (  ftsw(ji,jj+1) * zhU(ji-1,jj+1) &
1298                &                                         + ftse(ji,jj+1) * zhU(ji  ,jj+1) &
1299                &                                         + ftnw(ji,jj  ) * zhU(ji-1,jj  ) &
1300                &                                         + ftne(ji,jj  ) * zhU(ji  ,jj  ) )
1301            END DO
1302         END DO
1303         !
1304      END SELECT
1305      !
1306   END SUBROUTINE dyn_cor_2D
1307
1308
1309   SUBROUTINE wad_tmsk( pssh, ptmsk )
1310      !!----------------------------------------------------------------------
1311      !!                  ***  ROUTINE wad_lmt  ***
1312      !!                   
1313      !! ** Purpose :   set wetting & drying mask at tracer points
1314      !!              for the current barotropic sub-step
1315      !!
1316      !! ** Method  :   ???
1317      !!
1318      !! ** Action  :  ptmsk : wetting & drying t-mask
1319      !!----------------------------------------------------------------------
1320      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ) ::   pssh    !
1321      REAL(wp), DIMENSION(jpi,jpj), INTENT(  out) ::   ptmsk   !
1322      !
1323      INTEGER  ::   ji, jj   ! dummy loop indices
1324      !!----------------------------------------------------------------------
1325      !
1326      IF( ln_wd_dl_rmp ) THEN     
1327         DO jj = 1, jpj
1328            DO ji = 1, jpi                   
1329               IF    ( pssh(ji,jj) + ht_0(ji,jj) >  2._wp * rn_wdmin1 ) THEN 
1330                  !           IF    ( pssh(ji,jj) + ht_0(ji,jj) >          rn_wdmin2 ) THEN
1331                  ptmsk(ji,jj) = 1._wp
1332               ELSEIF( pssh(ji,jj) + ht_0(ji,jj) >          rn_wdmin1 ) THEN
1333                  ptmsk(ji,jj) = TANH( 50._wp*( ( pssh(ji,jj) + ht_0(ji,jj) -  rn_wdmin1 )*r_rn_wdmin1) )
1334               ELSE
1335                  ptmsk(ji,jj) = 0._wp
1336               ENDIF
1337            END DO
1338         END DO
1339      ELSE 
1340         DO jj = 1, jpj
1341            DO ji = 1, jpi                             
1342               IF ( pssh(ji,jj) + ht_0(ji,jj) >  rn_wdmin1 ) THEN   ;   ptmsk(ji,jj) = 1._wp
1343               ELSE                                                 ;   ptmsk(ji,jj) = 0._wp
1344               ENDIF
1345            END DO
1346         END DO
1347      ENDIF
1348      !
1349   END SUBROUTINE wad_tmsk
1350
1351
1352   SUBROUTINE wad_Umsk( pTmsk, phU, phV, pu, pv, pUmsk, pVmsk )
1353      !!----------------------------------------------------------------------
1354      !!                  ***  ROUTINE wad_lmt  ***
1355      !!                   
1356      !! ** Purpose :   set wetting & drying mask at tracer points
1357      !!              for the current barotropic sub-step
1358      !!
1359      !! ** Method  :   ???
1360      !!
1361      !! ** Action  :  ptmsk : wetting & drying t-mask
1362      !!----------------------------------------------------------------------
1363      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ) ::   pTmsk              ! W & D t-mask
1364      REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) ::   phU, phV, pu, pv   ! ocean velocities and transports
1365      REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) ::   pUmsk, pVmsk       ! W & D u- and v-mask
1366      !
1367      INTEGER  ::   ji, jj   ! dummy loop indices
1368      !!----------------------------------------------------------------------
1369      !
1370      DO jj = 1, jpj
1371         DO ji = 1, jpim1   ! not jpi-column
1372            IF ( phU(ji,jj) > 0._wp ) THEN   ;   pUmsk(ji,jj) = pTmsk(ji  ,jj) 
1373            ELSE                             ;   pUmsk(ji,jj) = pTmsk(ji+1,jj) 
1374            ENDIF
1375            phU(ji,jj) = pUmsk(ji,jj)*phU(ji,jj)
1376            pu (ji,jj) = pUmsk(ji,jj)*pu (ji,jj)
1377         END DO
1378      END DO
1379      !
1380      DO jj = 1, jpjm1   ! not jpj-row
1381         DO ji = 1, jpi
1382            IF ( phV(ji,jj) > 0._wp ) THEN   ;   pVmsk(ji,jj) = pTmsk(ji,jj  )
1383            ELSE                             ;   pVmsk(ji,jj) = pTmsk(ji,jj+1) 
1384            ENDIF
1385            phV(ji,jj) = pVmsk(ji,jj)*phV(ji,jj) 
1386            pv (ji,jj) = pVmsk(ji,jj)*pv (ji,jj)
1387         END DO
1388      END DO
1389      !
1390   END SUBROUTINE wad_Umsk
1391
1392
1393   SUBROUTINE wad_spg( sshn, zcpx, zcpy )
1394      !!---------------------------------------------------------------------
1395      !!                   ***  ROUTINE  wad_sp  ***
1396      !!
1397      !! ** Purpose :
1398      !!----------------------------------------------------------------------
1399      INTEGER  ::   ji ,jj               ! dummy loop indices
1400      LOGICAL  ::   ll_tmp1, ll_tmp2
1401      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ) :: sshn
1402      REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) :: zcpx, zcpy
1403      !!----------------------------------------------------------------------
1404      DO jj = 2, jpjm1
1405         DO ji = 2, jpim1 
1406            ll_tmp1 = MIN(  sshn(ji,jj)               ,  sshn(ji+1,jj) ) >                &
1407                 &      MAX( -ht_0(ji,jj)               , -ht_0(ji+1,jj) ) .AND.            &
1408                 &      MAX(  sshn(ji,jj) + ht_0(ji,jj) ,  sshn(ji+1,jj) + ht_0(ji+1,jj) )  &
1409                 &                                                         > rn_wdmin1 + rn_wdmin2
1410            ll_tmp2 = ( ABS( sshn(ji+1,jj)            -  sshn(ji  ,jj))  > 1.E-12 ).AND.( &
1411                 &      MAX(   sshn(ji,jj)              ,  sshn(ji+1,jj) ) >                &
1412                 &      MAX(  -ht_0(ji,jj)              , -ht_0(ji+1,jj) ) + rn_wdmin1 + rn_wdmin2 )
1413            IF(ll_tmp1) THEN
1414               zcpx(ji,jj) = 1.0_wp
1415            ELSEIF(ll_tmp2) THEN
1416               ! no worries about  sshn(ji+1,jj) -  sshn(ji  ,jj) = 0, it won't happen ! here
1417               zcpx(ji,jj) = ABS( (sshn(ji+1,jj) + ht_0(ji+1,jj) - sshn(ji,jj) - ht_0(ji,jj)) &
1418                    &    / (sshn(ji+1,jj) - sshn(ji  ,jj)) )
1419               zcpx(ji,jj) = max(min( zcpx(ji,jj) , 1.0_wp),0.0_wp)
1420            ELSE
1421               zcpx(ji,jj) = 0._wp
1422            ENDIF
1423            !
1424            ll_tmp1 = MIN(  sshn(ji,jj)               ,  sshn(ji,jj+1) ) >                &
1425                 &      MAX( -ht_0(ji,jj)               , -ht_0(ji,jj+1) ) .AND.            &
1426                 &      MAX(  sshn(ji,jj) + ht_0(ji,jj) ,  sshn(ji,jj+1) + ht_0(ji,jj+1) )  &
1427                 &                                                       > rn_wdmin1 + rn_wdmin2
1428            ll_tmp2 = ( ABS( sshn(ji,jj)              -  sshn(ji,jj+1))  > 1.E-12 ).AND.( &
1429                 &      MAX(   sshn(ji,jj)              ,  sshn(ji,jj+1) ) >                &
1430                 &      MAX(  -ht_0(ji,jj)              , -ht_0(ji,jj+1) ) + rn_wdmin1 + rn_wdmin2 )
1431           
1432            IF(ll_tmp1) THEN
1433               zcpy(ji,jj) = 1.0_wp
1434            ELSE IF(ll_tmp2) THEN
1435               ! no worries about  sshn(ji,jj+1) -  sshn(ji,jj  ) = 0, it won't happen ! here
1436               zcpy(ji,jj) = ABS( (sshn(ji,jj+1) + ht_0(ji,jj+1) - sshn(ji,jj) - ht_0(ji,jj)) &
1437                    &             / (sshn(ji,jj+1) - sshn(ji,jj  )) )
1438               zcpy(ji,jj) = MAX(  0._wp , MIN( zcpy(ji,jj) , 1.0_wp )  )
1439            ELSE
1440               zcpy(ji,jj) = 0._wp
1441            ENDIF
1442         END DO
1443      END DO
1444           
1445   END SUBROUTINE wad_spg
1446     
1447
1448
1449   SUBROUTINE dyn_drg_init( pu_RHSi, pv_RHSi, pCdU_u, pCdU_v )
1450      !!----------------------------------------------------------------------
1451      !!                  ***  ROUTINE dyn_drg_init  ***
1452      !!                   
1453      !! ** Purpose : - add the baroclinic top/bottom drag contribution to
1454      !!              the baroclinic part of the barotropic RHS
1455      !!              - compute the barotropic drag coefficients
1456      !!
1457      !! ** Method  :   computation done over the INNER domain only
1458      !!----------------------------------------------------------------------
1459      REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) ::   pu_RHSi, pv_RHSi   ! baroclinic part of the barotropic RHS
1460      REAL(wp), DIMENSION(jpi,jpj), INTENT(  out) ::   pCdU_u , pCdU_v    ! barotropic drag coefficients
1461      !
1462      INTEGER  ::   ji, jj   ! dummy loop indices
1463      INTEGER  ::   ikbu, ikbv, iktu, iktv
1464      REAL(wp) ::   zztmp
1465      REAL(wp), DIMENSION(jpi,jpj) ::   zu_i, zv_i
1466      !!----------------------------------------------------------------------
1467      !
1468      !                    !==  Set the barotropic drag coef.  ==!
1469      !
1470      IF( ln_isfcav ) THEN          ! top+bottom friction (ocean cavities)
1471         
1472         DO jj = 2, jpjm1
1473            DO ji = 2, jpim1     ! INNER domain
1474               pCdU_u(ji,jj) = r1_2*( rCdU_bot(ji+1,jj)+rCdU_bot(ji,jj) + rCdU_top(ji+1,jj)+rCdU_top(ji,jj) )
1475               pCdU_v(ji,jj) = r1_2*( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj) + rCdU_top(ji,jj+1)+rCdU_top(ji,jj) )
1476            END DO
1477         END DO
1478      ELSE                          ! bottom friction only
1479         DO jj = 2, jpjm1
1480            DO ji = 2, jpim1  ! INNER domain
1481               pCdU_u(ji,jj) = r1_2*( rCdU_bot(ji+1,jj)+rCdU_bot(ji,jj) )
1482               pCdU_v(ji,jj) = r1_2*( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj) )
1483            END DO
1484         END DO
1485      ENDIF
1486      !
1487      !                    !==  BOTTOM stress contribution from baroclinic velocities  ==!
1488      !
1489      IF( ln_bt_fw ) THEN                 ! FORWARD integration: use NOW bottom baroclinic velocities
1490         
1491         DO jj = 2, jpjm1
1492            DO ji = 2, jpim1  ! INNER domain
1493               ikbu = mbku(ji,jj)       
1494               ikbv = mbkv(ji,jj)   
1495               zu_i(ji,jj) = un(ji,jj,ikbu) - un_b(ji,jj)
1496               zv_i(ji,jj) = vn(ji,jj,ikbv) - vn_b(ji,jj)
1497            END DO
1498         END DO
1499      ELSE                                ! CENTRED integration: use BEFORE bottom baroclinic velocities
1500         
1501         DO jj = 2, jpjm1
1502            DO ji = 2, jpim1   ! INNER domain
1503               ikbu = mbku(ji,jj)       
1504               ikbv = mbkv(ji,jj)   
1505               zu_i(ji,jj) = ub(ji,jj,ikbu) - ub_b(ji,jj)
1506               zv_i(ji,jj) = vb(ji,jj,ikbv) - vb_b(ji,jj)
1507            END DO
1508         END DO
1509      ENDIF
1510      !
1511      IF( ln_wd_il ) THEN      ! W/D : use the "clipped" bottom friction   !!gm   explain WHY, please !
1512         zztmp = -1._wp / rdtbt
1513         DO jj = 2, jpjm1
1514            DO ji = 2, jpim1    ! INNER domain
1515               pu_RHSi(ji,jj) = pu_RHSi(ji,jj) + zu_i(ji,jj) *  wdrampu(ji,jj) * MAX(                                 & 
1516                    &                              r1_hu_n(ji,jj) * r1_2*( rCdU_bot(ji+1,jj)+rCdU_bot(ji,jj) ) , zztmp  )
1517               pv_RHSi(ji,jj) = pv_RHSi(ji,jj) + zv_i(ji,jj) *  wdrampv(ji,jj) * MAX(                                 & 
1518                    &                              r1_hv_n(ji,jj) * r1_2*( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj) ) , zztmp  )
1519            END DO
1520         END DO
1521      ELSE                    ! use "unclipped" drag (even if explicit friction is used in 3D calculation)
1522         
1523         DO jj = 2, jpjm1
1524            DO ji = 2, jpim1    ! INNER domain
1525               pu_RHSi(ji,jj) = pu_RHSi(ji,jj) + r1_hu_n(ji,jj) * r1_2*( rCdU_bot(ji+1,jj)+rCdU_bot(ji,jj) ) * zu_i(ji,jj)
1526               pv_RHSi(ji,jj) = pv_RHSi(ji,jj) + r1_hv_n(ji,jj) * r1_2*( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj) ) * zv_i(ji,jj)
1527            END DO
1528         END DO
1529      END IF
1530      !
1531      !                    !==  TOP stress contribution from baroclinic velocities  ==!   (no W/D case)
1532      !
1533      IF( ln_isfcav ) THEN
1534         !
1535         IF( ln_bt_fw ) THEN                ! FORWARD integration: use NOW top baroclinic velocity
1536           
1537            DO jj = 2, jpjm1
1538               DO ji = 2, jpim1   ! INNER domain
1539                  iktu = miku(ji,jj)
1540                  iktv = mikv(ji,jj)
1541                  zu_i(ji,jj) = un(ji,jj,iktu) - un_b(ji,jj)
1542                  zv_i(ji,jj) = vn(ji,jj,iktv) - vn_b(ji,jj)
1543               END DO
1544            END DO
1545         ELSE                                ! CENTRED integration: use BEFORE top baroclinic velocity
1546           
1547            DO jj = 2, jpjm1
1548               DO ji = 2, jpim1      ! INNER domain
1549                  iktu = miku(ji,jj)
1550                  iktv = mikv(ji,jj)
1551                  zu_i(ji,jj) = ub(ji,jj,iktu) - ub_b(ji,jj)
1552                  zv_i(ji,jj) = vb(ji,jj,iktv) - vb_b(ji,jj)
1553               END DO
1554            END DO
1555         ENDIF
1556         !
1557         !                    ! use "unclipped" top drag (even if explicit friction is used in 3D calculation)
1558         
1559         DO jj = 2, jpjm1
1560            DO ji = 2, jpim1    ! INNER domain
1561               pu_RHSi(ji,jj) = pu_RHSi(ji,jj) + r1_hu_n(ji,jj) * r1_2*( rCdU_top(ji+1,jj)+rCdU_top(ji,jj) ) * zu_i(ji,jj)
1562               pv_RHSi(ji,jj) = pv_RHSi(ji,jj) + r1_hv_n(ji,jj) * r1_2*( rCdU_top(ji,jj+1)+rCdU_top(ji,jj) ) * zv_i(ji,jj)
1563            END DO
1564         END DO
1565         !
1566      ENDIF
1567      !
1568   END SUBROUTINE dyn_drg_init
1569
1570   SUBROUTINE ts_bck_interp( jn, ll_init,       &   ! <== in
1571      &                      za0, za1, za2, za3 )   ! ==> out
1572      !!----------------------------------------------------------------------
1573      INTEGER ,INTENT(in   ) ::   jn                   ! index of sub time step
1574      LOGICAL ,INTENT(in   ) ::   ll_init              !
1575      REAL(wp),INTENT(  out) ::   za0, za1, za2, za3   ! Half-step back interpolation coefficient
1576      !
1577      REAL(wp) ::   zepsilon, zgamma                   !   -      -
1578      !!----------------------------------------------------------------------
1579      !                             ! set Half-step back interpolation coefficient
1580      IF    ( jn==1 .AND. ll_init ) THEN   !* Forward-backward
1581         za0 = 1._wp                       
1582         za1 = 0._wp                           
1583         za2 = 0._wp
1584         za3 = 0._wp
1585      ELSEIF( jn==2 .AND. ll_init ) THEN   !* AB2-AM3 Coefficients; bet=0 ; gam=-1/6 ; eps=1/12
1586         za0 = 1.0833333333333_wp                 ! za0 = 1-gam-eps
1587         za1 =-0.1666666666666_wp                 ! za1 = gam
1588         za2 = 0.0833333333333_wp                 ! za2 = eps
1589         za3 = 0._wp             
1590      ELSE                                 !* AB3-AM4 Coefficients; bet=0.281105 ; eps=0.013 ; gam=0.0880
1591         IF( rn_bt_alpha == 0._wp ) THEN      ! Time diffusion 
1592            za0 = 0.614_wp                        ! za0 = 1/2 +   gam + 2*eps
1593            za1 = 0.285_wp                        ! za1 = 1/2 - 2*gam - 3*eps
1594            za2 = 0.088_wp                        ! za2 = gam
1595            za3 = 0.013_wp                        ! za3 = eps
1596         ELSE                                 ! no time diffusion
1597            zepsilon = 0.00976186_wp - 0.13451357_wp * rn_bt_alpha
1598            zgamma   = 0.08344500_wp - 0.51358400_wp * rn_bt_alpha
1599            za0 = 0.5_wp + zgamma + 2._wp * rn_bt_alpha + 2._wp * zepsilon
1600            za1 = 1._wp - za0 - zgamma - zepsilon
1601            za2 = zgamma
1602            za3 = zepsilon
1603         ENDIF
1604      ENDIF
1605   END SUBROUTINE ts_bck_interp
1606
1607
1608   !!======================================================================
1609END MODULE dynspg_ts
Note: See TracBrowser for help on using the repository browser.