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/UKMO_MERGE_2019/src/OCE/DYN – NEMO

source: NEMO/branches/2019/UKMO_MERGE_2019/src/OCE/DYN/dynspg_ts.F90 @ 12068

Last change on this file since 12068 was 12068, checked in by davestorkey, 5 years ago

2019/UKMO_MERGE_2019 : Merging in changes from ENHANCE-02_ISF_nemo.

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