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/2021/dev_r14116_HPC-10_mcastril_Mixed_Precision_implementation/src/OCE/DYN – NEMO

source: NEMO/branches/2021/dev_r14116_HPC-10_mcastril_Mixed_Precision_implementation/src/OCE/DYN/dynspg_ts.F90 @ 14944

Last change on this file since 14944 was 14944, checked in by sparonuz, 3 years ago

Added interface for function dyn_drg_init

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