New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
dynspg_ts.F90 in branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_ts.F90 @ 3953

Last change on this file since 3953 was 3953, checked in by gm, 11 years ago

dev_r3858_NOC_ZTC, #863 : activate tide potential in filtered ssh case + style in tide modules

  • Property svn:keywords set to Id
File size: 44.2 KB
Line 
1MODULE dynspg_ts
2   !!======================================================================
3   !! History :   1.0  ! 2004-12  (L. Bessieres, G. Madec)  Original code
4   !!              -   ! 2005-11  (V. Garnier, G. Madec)  optimization
5   !!              -   ! 2006-08  (S. Masson)  distributed restart using iom
6   !!             2.0  ! 2007-07  (D. Storkey) calls to BDY routines
7   !!              -   ! 2008-01  (R. Benshila)  change averaging method
8   !!             3.2  ! 2009-07  (R. Benshila, G. Madec) Complete revisit associated to vvl reactivation
9   !!             3.3  ! 2010-09  (D. Storkey, E. O'Dea) update for BDY for Shelf configurations
10   !!             3.3  ! 2011-03  (R. Benshila, R. Hordoir, P. Oddo) update calculation of ub_b
11   !!---------------------------------------------------------------------
12#if defined key_dynspg_ts   ||   defined key_esopa
13   !!----------------------------------------------------------------------
14   !!   'key_dynspg_ts'         free surface cst volume with time splitting
15   !!----------------------------------------------------------------------
16   !!   dyn_spg_ts  : compute surface pressure gradient trend using a time-
17   !!                 splitting scheme and add to the general trend
18   !!   ts_rst      : read/write the time-splitting restart fields in the ocean restart file
19   !!----------------------------------------------------------------------
20   USE oce             ! ocean dynamics and tracers
21   USE dom_oce         ! ocean space and time domain
22   USE sbc_oce         ! surface boundary condition: ocean
23   USE dynspg_oce      ! surface pressure gradient variables
24   USE phycst          ! physical constants
25   USE domvvl          ! variable volume
26   USE zdfbfr          ! bottom friction
27   USE dynvor          ! vorticity term
28   USE obc_oce         ! Lateral open boundary condition
29   USE obc_par         ! open boundary condition parameters
30   USE obcdta          ! open boundary condition data     
31   USE obcfla          ! Flather open boundary condition 
32   USE bdy_par         ! for lk_bdy
33   USE bdy_oce         ! Lateral open boundary condition
34   USE bdydta          ! open boundary condition data     
35   USE bdydyn2d        ! open boundary conditions on barotropic variables
36   USE sbctide
37   USE updtide
38   USE lib_mpp         ! distributed memory computing library
39   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
40   USE prtctl          ! Print control
41   USE in_out_manager  ! I/O manager
42   USE iom             ! IOM library
43   USE zdf_oce         ! Vertical diffusion
44   USE wrk_nemo        ! Memory Allocation
45   USE timing          ! Timing
46
47
48   IMPLICIT NONE
49   PRIVATE
50
51   PUBLIC dyn_spg_ts        ! routine called by step.F90
52   PUBLIC ts_rst            ! routine called by istate.F90
53   PUBLIC dyn_spg_ts_alloc  ! routine called by dynspg.F90
54
55
56   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::  ftnw, ftne   ! triad of coriolis parameter
57   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::  ftsw, ftse   ! (only used with een vorticity scheme)
58
59   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   un_b, vn_b   ! now    averaged velocity
60   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   ub_b, vb_b   ! before averaged velocity
61
62   !! * Substitutions
63#  include "domzgr_substitute.h90"
64#  include "vectopt_loop_substitute.h90"
65   !!----------------------------------------------------------------------
66   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
67   !! $Id$
68   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
69   !!----------------------------------------------------------------------
70CONTAINS
71
72   INTEGER FUNCTION dyn_spg_ts_alloc()
73      !!----------------------------------------------------------------------
74      !!                  ***  routine dyn_spg_ts_alloc  ***
75      !!----------------------------------------------------------------------
76      ALLOCATE( ftnw  (jpi,jpj) , ftne(jpi,jpj) , un_b(jpi,jpj) , vn_b(jpi,jpj) ,     &
77         &      ftsw  (jpi,jpj) , ftse(jpi,jpj) , ub_b(jpi,jpj) , vb_b(jpi,jpj) , STAT= dyn_spg_ts_alloc )
78         !
79      IF( lk_mpp                )   CALL mpp_sum( dyn_spg_ts_alloc )
80      IF( dyn_spg_ts_alloc /= 0 )   CALL ctl_warn('dynspg_oce_alloc: failed to allocate arrays')
81      !
82   END FUNCTION dyn_spg_ts_alloc
83
84
85   SUBROUTINE dyn_spg_ts( kt )
86      !!----------------------------------------------------------------------
87      !!                  ***  routine dyn_spg_ts  ***
88      !!
89      !! ** Purpose :   Compute the now trend due to the surface pressure
90      !!      gradient in case of free surface formulation with time-splitting.
91      !!      Add it to the general trend of momentum equation.
92      !!
93      !! ** Method  :   Free surface formulation with time-splitting
94      !!      -1- Save the vertically integrated trend. This general trend is
95      !!          held constant over the barotropic integration.
96      !!          The Coriolis force is removed from the general trend as the
97      !!          surface gradient and the Coriolis force are updated within
98      !!          the barotropic integration.
99      !!      -2- Barotropic loop : updates of sea surface height (ssha_e) and
100      !!          barotropic velocity (ua_e and va_e) through barotropic
101      !!          momentum and continuity integration. Barotropic former
102      !!          variables are time averaging over the full barotropic cycle
103      !!          (= 2 * baroclinic time step) and saved in uX_b
104      !!          and vX_b (X specifying after, now or before).
105      !!      -3- The new general trend becomes :
106      !!          ua = ua - sum_k(ua)/H + ( un_b - ub_b )
107      !!
108      !! ** Action : - Update (ua,va) with the surf. pressure gradient trend
109      !!
110      !! References : Griffies et al., (2003): A technical guide to MOM4. NOAA/GFDL
111      !!---------------------------------------------------------------------
112      !
113      INTEGER, INTENT(in)  ::   kt   ! ocean time-step index
114      !
115      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices
116      INTEGER  ::   icycle           ! local scalar
117      INTEGER  ::   ikbu, ikbv       ! local scalar
118      REAL(wp) ::   zraur, zcoef, z2dt_e, z1_2dt_b, z2dt_bf   ! local scalars
119      REAL(wp) ::   z1_8, zx1, zy1                            !   -      -
120      REAL(wp) ::   z1_4, zx2, zy2                            !   -      -
121      REAL(wp) ::   zu_spg, zu_cor, zu_sld, zu_asp            !   -      -
122      REAL(wp) ::   zv_spg, zv_cor, zv_sld, zv_asp            !   -      -
123      REAL(wp) ::   ua_btm, va_btm                            !   -      -
124      !
125      REAL(wp), POINTER, DIMENSION(:,:) :: zsshun_e, zsshvn_e, zsshb_e, zssh_sum, zhdiv 
126      REAL(wp), POINTER, DIMENSION(:,:) :: zua, zva, zun, zvn, zun_e, zvn_e, zub_e, zvb_e 
127      REAL(wp), POINTER, DIMENSION(:,:) :: zcu, zcv, zwx, zwy, zbfru, zbfrv, zu_sum, zv_sum
128      REAL(wp), POINTER, DIMENSION(:,:) :: zhu_b, zhv_b
129      !!----------------------------------------------------------------------
130      !
131      IF( nn_timing == 1 )  CALL timing_start('dyn_spg_ts')
132      !
133      CALL wrk_alloc( jpi, jpj, zsshun_e, zsshvn_e, zsshb_e, zssh_sum, zhdiv     )
134      CALL wrk_alloc( jpi, jpj, zua, zva, zun, zvn, zun_e, zvn_e, zub_e, zvb_e   )
135      CALL wrk_alloc( jpi, jpj, zcu, zcv, zwx, zwy, zbfru, zbfrv, zu_sum, zv_sum )
136      CALL wrk_alloc( jpi, jpj, zhu_b, zhv_b                                     )
137      !
138      IF( kt == nit000 ) THEN             !* initialisation
139         !
140         IF(lwp) WRITE(numout,*)
141         IF(lwp) WRITE(numout,*) 'dyn_spg_ts : surface pressure gradient trend'
142         IF(lwp) WRITE(numout,*) '~~~~~~~~~~   free surface with time splitting'
143         IF(lwp) WRITE(numout,*) ' Number of sub cycle in 1 time-step (2 rdt) : icycle = ',  2*nn_baro
144         !
145         CALL ts_rst( nit000, 'READ' )   ! read or initialize the following fields: un_b, vn_b 
146         !
147         ua_e  (:,:) = un_b (:,:)
148         va_e  (:,:) = vn_b (:,:)
149         hu_e  (:,:) = hu   (:,:)
150         hv_e  (:,:) = hv   (:,:)
151         hur_e (:,:) = hur  (:,:)
152         hvr_e (:,:) = hvr  (:,:)
153         IF( ln_dynvor_een ) THEN
154            ftne(1,:) = 0._wp   ;   ftnw(1,:) = 0._wp   ;   ftse(1,:) = 0._wp   ;   ftsw(1,:) = 0._wp
155            DO jj = 2, jpj
156               DO ji = fs_2, jpi   ! vector opt.
157                  ftne(ji,jj) = ( ff(ji-1,jj  ) + ff(ji  ,jj  ) + ff(ji  ,jj-1) ) / 3._wp
158                  ftnw(ji,jj) = ( ff(ji-1,jj-1) + ff(ji-1,jj  ) + ff(ji  ,jj  ) ) / 3._wp
159                  ftse(ji,jj) = ( ff(ji  ,jj  ) + ff(ji  ,jj-1) + ff(ji-1,jj-1) ) / 3._wp
160                  ftsw(ji,jj) = ( ff(ji  ,jj-1) + ff(ji-1,jj-1) + ff(ji-1,jj  ) ) / 3._wp
161               END DO
162            END DO
163         ENDIF
164         !
165      ENDIF
166
167      !                                                     !* Local constant initialization
168      z1_2dt_b = 1._wp / ( 2.0_wp * rdt )                   ! reciprocal of baroclinic time step
169      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt_b = 1.0_wp / rdt    ! reciprocal of baroclinic
170                                                                        ! time step (euler timestep)
171      z1_8     = 0.125_wp                                   ! coefficient for vorticity estimates
172      z1_4     = 0.25_wp       
173      zraur    = 1._wp / rau0                               ! 1 / volumic mass
174      !
175      zhdiv(:,:) = 0._wp                                    ! barotropic divergence
176      zu_sld = 0._wp   ;   zu_asp = 0._wp                   ! tides trends (lk_tide=F)
177      zv_sld = 0._wp   ;   zv_asp = 0._wp
178
179      IF( kt == nit000 .AND. neuler == 0) THEN              ! for implicit bottom friction
180        z2dt_bf = rdt
181      ELSE
182        z2dt_bf = 2.0_wp * rdt
183      ENDIF
184
185      ! -----------------------------------------------------------------------------
186      !  Phase 1 : Coupling between general trend and barotropic estimates (1st step)
187      ! -----------------------------------------------------------------------------
188      !     
189      !                                   !* e3*d/dt(Ua), e3*Ub, e3*Vn (Vertically integrated)
190      !                                   ! --------------------------
191      zua(:,:) = 0._wp   ;   zun(:,:) = 0._wp   ;   ub_b(:,:) = 0._wp
192      zva(:,:) = 0._wp   ;   zvn(:,:) = 0._wp   ;   vb_b(:,:) = 0._wp
193      !
194      DO jk = 1, jpkm1
195#if defined key_vectopt_loop
196         DO jj = 1, 1         !Vector opt. => forced unrolling
197            DO ji = 1, jpij
198#else
199         DO jj = 1, jpj
200            DO ji = 1, jpi
201#endif
202               !                                                                              ! now trend
203               zua(ji,jj) = zua(ji,jj) + fse3u_n(ji,jj,jk) * ua(ji,jj,jk) * umask(ji,jj,jk)
204               zva(ji,jj) = zva(ji,jj) + fse3v_n(ji,jj,jk) * va(ji,jj,jk) * vmask(ji,jj,jk)
205               !                                                                              ! now velocity
206               zun(ji,jj) = zun(ji,jj) + fse3u_n(ji,jj,jk) * un(ji,jj,jk)
207               zvn(ji,jj) = zvn(ji,jj) + fse3v_n(ji,jj,jk) * vn(ji,jj,jk)               
208               !
209               ub_b(ji,jj) = ub_b(ji,jj) + fse3u_b(ji,jj,jk) * ub(ji,jj,jk)  * umask(ji,jj,jk)
210               vb_b(ji,jj) = vb_b(ji,jj) + fse3v_b(ji,jj,jk) * vb(ji,jj,jk)  * vmask(ji,jj,jk)
211            END DO
212         END DO
213      END DO
214
215      ! before inverse water column height at u- and v- points
216      IF( lk_vvl ) THEN
217         zhu_b(:,:) = 0.
218         zhv_b(:,:) = 0.
219         DO jk = 1, jpk
220            zhu_b(:,:) = zhu_b(:,:) + fse3u_b(:,:,jk) * umask(:,:,jk)
221            zhv_b(:,:) = zhv_b(:,:) + fse3v_b(:,:,jk) * vmask(:,:,jk)
222         END DO
223         zhu_b(:,:) = umask(:,:,1) / ( zhu_b(:,:) + 1. - umask(:,:,1) )
224         zhv_b(:,:) = vmask(:,:,1) / ( zhv_b(:,:) + 1. - vmask(:,:,1) )
225      ELSE
226         zhu_b(:,:) = hur(:,:)
227         zhv_b(:,:) = hvr(:,:)
228      ENDIF
229
230      !                                   !* baroclinic momentum trend (remove the vertical mean trend)
231      DO jk = 1, jpkm1                    ! --------------------------
232         DO jj = 2, jpjm1
233            DO ji = fs_2, fs_jpim1   ! vector opt.
234               ua(ji,jj,jk) = ua(ji,jj,jk) - zua(ji,jj) * hur(ji,jj)
235               va(ji,jj,jk) = va(ji,jj,jk) - zva(ji,jj) * hvr(ji,jj)
236            END DO
237         END DO
238      END DO
239
240      !                                   !* barotropic Coriolis trends * H (vorticity scheme dependent)
241      !                                   ! ---------------------------====
242      zwx(:,:) = zun(:,:) * e2u(:,:)                   ! now transport
243      zwy(:,:) = zvn(:,:) * e1v(:,:)
244      !
245      IF( ln_dynvor_ene .OR. ln_dynvor_mix ) THEN      ! energy conserving or mixed scheme
246         DO jj = 2, jpjm1
247            DO ji = fs_2, fs_jpim1   ! vector opt.
248               zy1 = ( zwy(ji,jj-1) + zwy(ji+1,jj-1) ) / e1u(ji,jj)
249               zy2 = ( zwy(ji,jj  ) + zwy(ji+1,jj  ) ) / e1u(ji,jj)
250               zx1 = ( zwx(ji-1,jj) + zwx(ji-1,jj+1) ) / e2v(ji,jj)
251               zx2 = ( zwx(ji  ,jj) + zwx(ji  ,jj+1) ) / e2v(ji,jj)
252               ! energy conserving formulation for planetary vorticity term
253               zcu(ji,jj) = z1_4 * ( ff(ji  ,jj-1) * zy1 + ff(ji,jj) * zy2 )
254               zcv(ji,jj) =-z1_4 * ( ff(ji-1,jj  ) * zx1 + ff(ji,jj) * zx2 )
255            END DO
256         END DO
257         !
258      ELSEIF ( ln_dynvor_ens ) THEN                    ! enstrophy conserving scheme
259         DO jj = 2, jpjm1
260            DO ji = fs_2, fs_jpim1   ! vector opt.
261               zy1 =   z1_8 * ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) + zwy(ji,jj) + zwy(ji+1,jj  ) ) / e1u(ji,jj)
262               zx1 = - z1_8 * ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) + zwx(ji,jj) + zwx(ji  ,jj+1) ) / e2v(ji,jj)
263               zcu(ji,jj)  = zy1 * ( ff(ji  ,jj-1) + ff(ji,jj) )
264               zcv(ji,jj)  = zx1 * ( ff(ji-1,jj  ) + ff(ji,jj) )
265            END DO
266         END DO
267         !
268      ELSEIF ( ln_dynvor_een ) THEN                    ! enstrophy and energy conserving scheme
269         DO jj = 2, jpjm1
270            DO ji = fs_2, fs_jpim1   ! vector opt.
271               zcu(ji,jj) = + z1_4 / e1u(ji,jj) * (  ftne(ji,jj  ) * zwy(ji  ,jj  ) + ftnw(ji+1,jj) * zwy(ji+1,jj  )   &
272                  &                                + ftse(ji,jj  ) * zwy(ji  ,jj-1) + ftsw(ji+1,jj) * zwy(ji+1,jj-1) )
273               zcv(ji,jj) = - z1_4 / e2v(ji,jj) * (  ftsw(ji,jj+1) * zwx(ji-1,jj+1) + ftse(ji,jj+1) * zwx(ji  ,jj+1)   &
274                  &                                + ftnw(ji,jj  ) * zwx(ji-1,jj  ) + ftne(ji,jj  ) * zwx(ji  ,jj  ) )
275            END DO
276         END DO
277         !
278      ENDIF
279
280      !                                   !* Right-Hand-Side of the barotropic momentum equation
281      !                                   ! ----------------------------------------------------
282      IF( lk_vvl ) THEN                         ! Variable volume : remove both Coriolis and Surface pressure gradient
283         DO jj = 2, jpjm1 
284            DO ji = fs_2, fs_jpim1   ! vector opt.
285               zcu(ji,jj) = zcu(ji,jj) - grav * (  ( rhd(ji+1,jj  ,1) + 1 ) * sshn(ji+1,jj  )    &
286                  &                              - ( rhd(ji  ,jj  ,1) + 1 ) * sshn(ji  ,jj  )  ) * hu(ji,jj) / e1u(ji,jj)
287               zcv(ji,jj) = zcv(ji,jj) - grav * (  ( rhd(ji  ,jj+1,1) + 1 ) * sshn(ji  ,jj+1)    &
288                  &                              - ( rhd(ji  ,jj  ,1) + 1 ) * sshn(ji  ,jj  )  ) * hv(ji,jj) / e2v(ji,jj)
289            END DO
290         END DO
291      ENDIF
292
293      DO jj = 2, jpjm1                             ! Remove coriolis term (and possibly spg) from barotropic trend
294         DO ji = fs_2, fs_jpim1
295             zua(ji,jj) = zua(ji,jj) - zcu(ji,jj)
296             zva(ji,jj) = zva(ji,jj) - zcv(ji,jj)
297          END DO
298      END DO
299
300                   
301      !                                             ! Remove barotropic contribution of bottom friction
302      !                                             ! from the barotropic transport trend
303      zcoef = -1._wp * z1_2dt_b
304
305      IF(ln_bfrimp) THEN
306      !                                   ! Remove the bottom stress trend from 3-D sea surface level gradient
307      !                                   ! and Coriolis forcing in case of 3D semi-implicit bottom friction
308        DO jj = 2, jpjm1         
309           DO ji = fs_2, fs_jpim1
310              ikbu = mbku(ji,jj)
311              ikbv = mbkv(ji,jj)
312              ua_btm = zcu(ji,jj) * z2dt_bf * hur(ji,jj) * umask (ji,jj,ikbu)
313              va_btm = zcv(ji,jj) * z2dt_bf * hvr(ji,jj) * vmask (ji,jj,ikbv)
314
315              zua(ji,jj) = zua(ji,jj) - bfrua(ji,jj) * ua_btm
316              zva(ji,jj) = zva(ji,jj) - bfrva(ji,jj) * va_btm
317           END DO
318        END DO
319
320      ELSE
321
322# if defined key_vectopt_loop
323        DO jj = 1, 1
324           DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
325# else
326        DO jj = 2, jpjm1
327           DO ji = 2, jpim1
328# endif
329            ! Apply stability criteria for bottom friction
330            !RBbug for vvl and external mode we may need to use varying fse3
331            !!gm  Rq: the bottom e3 present the smallest variation, the use of e3u_0 is not a big approx.
332              zbfru(ji,jj) = MAX(  bfrua(ji,jj) , fse3u(ji,jj,mbku(ji,jj)) * zcoef  )
333              zbfrv(ji,jj) = MAX(  bfrva(ji,jj) , fse3v(ji,jj,mbkv(ji,jj)) * zcoef  )
334           END DO
335        END DO
336
337        DO jj = 2, jpjm1
338           DO ji = fs_2, fs_jpim1   ! vector opt.
339              zua(ji,jj) = zua(ji,jj) - zbfru(ji,jj) * ub_b(ji,jj) * zhu_b(ji,jj)
340              zva(ji,jj) = zva(ji,jj) - zbfrv(ji,jj) * vb_b(ji,jj) * zhv_b(ji,jj)
341           END DO
342        END DO
343      END IF    ! end (ln_bfrimp)
344
345                   
346      !                                   !* d/dt(Ua), Ub, Vn (Vertical mean velocity)
347      !                                   ! --------------------------
348      zua(:,:) = zua(:,:) * hur(:,:)
349      zva(:,:) = zva(:,:) * hvr(:,:)
350      !
351      ub_b(:,:) = ub_b(:,:) * zhu_b(:,:)
352      vb_b(:,:) = vb_b(:,:) * zhv_b(:,:)
353
354      ! -----------------------------------------------------------------------
355      !  Phase 2 : Integration of the barotropic equations with time splitting
356      ! -----------------------------------------------------------------------
357      !
358      !                                             ! ==================== !
359      !                                             !    Initialisations   !
360      !                                             ! ==================== !
361      icycle = 2  * nn_baro            ! Number of barotropic sub time-step
362     
363      !                                ! Start from NOW field
364      hu_e   (:,:) = hu   (:,:)            ! ocean depth at u- and v-points
365      hv_e   (:,:) = hv   (:,:) 
366      hur_e  (:,:) = hur  (:,:)            ! ocean depth inverted at u- and v-points
367      hvr_e  (:,:) = hvr  (:,:)
368!RBbug     zsshb_e(:,:) = sshn (:,:) 
369      zsshb_e(:,:) = sshn_b(:,:)           ! sea surface height (before and now)
370      sshn_e (:,:) = sshn (:,:)
371     
372      zun_e  (:,:) = un_b (:,:)            ! barotropic velocity (external)
373      zvn_e  (:,:) = vn_b (:,:)
374      zub_e  (:,:) = un_b (:,:)
375      zvb_e  (:,:) = vn_b (:,:)
376
377      zu_sum  (:,:) = un_b (:,:)           ! summation
378      zv_sum  (:,:) = vn_b (:,:)
379      zssh_sum(:,:) = sshn (:,:)
380
381#if defined key_obc
382      ! set ssh corrections to 0
383      ! ssh corrections are applied to normal velocities (Flather's algorithm) and averaged over the barotropic loop
384      IF( lp_obc_east  )   sshfoe_b(:,:) = 0._wp
385      IF( lp_obc_west  )   sshfow_b(:,:) = 0._wp
386      IF( lp_obc_south )   sshfos_b(:,:) = 0._wp
387      IF( lp_obc_north )   sshfon_b(:,:) = 0._wp
388#endif
389
390      !                                             ! ==================== !
391      DO jn = 1, icycle                             !  sub-time-step loop  ! (from NOW to AFTER+1)
392         !                                          ! ==================== !
393         z2dt_e = 2. * ( rdt / nn_baro )
394         IF( jn == 1 )   z2dt_e = rdt / nn_baro
395
396         !                                                !* Update the forcing (BDY and tides)
397         !                                                !  ------------------
398         IF( lk_obc                    )   CALL obc_dta_bt( kt, jn   )
399         IF( lk_bdy                    )   CALL bdy_dta   ( kt, jit=jn, time_offset=1 )
400         IF( ln_tide_pot .AND. lk_tide )   CALL upd_tide  ( kt, kit=jn, kbaro=nn_baro )
401
402         !                                                !* after ssh_e
403         !                                                !  -----------
404         DO jj = 2, jpjm1                                 ! Horizontal divergence of barotropic transports
405            DO ji = fs_2, fs_jpim1   ! vector opt.
406               zhdiv(ji,jj) = (   e2u(ji  ,jj) * zun_e(ji  ,jj) * hu_e(ji  ,jj)     &
407                  &             - e2u(ji-1,jj) * zun_e(ji-1,jj) * hu_e(ji-1,jj)     &
408                  &             + e1v(ji,jj  ) * zvn_e(ji,jj  ) * hv_e(ji,jj  )     &
409                  &             - e1v(ji,jj-1) * zvn_e(ji,jj-1) * hv_e(ji,jj-1)   ) / ( e1t(ji,jj) * e2t(ji,jj) )
410            END DO
411         END DO
412         !
413#if defined key_obc
414         !                                                     ! OBC : zhdiv must be zero behind the open boundary
415!!  mpp remark: The zeroing of hdiv can probably be extended to 1->jpi/jpj for the correct row/column
416         IF( lp_obc_east  )   zhdiv(nie0p1:nie1p1,nje0  :nje1  ) = 0._wp      ! east
417         IF( lp_obc_west  )   zhdiv(niw0  :niw1  ,njw0  :njw1  ) = 0._wp      ! west
418         IF( lp_obc_north )   zhdiv(nin0  :nin1  ,njn0p1:njn1p1) = 0._wp      ! north
419         IF( lp_obc_south )   zhdiv(nis0  :nis1  ,njs0  :njs1  ) = 0._wp      ! south
420#endif
421#if defined key_bdy
422         zhdiv(:,:) = zhdiv(:,:) * bdytmask(:,:)               ! BDY mask
423#endif
424         !
425         DO jj = 2, jpjm1                                      ! leap-frog on ssh_e
426            DO ji = fs_2, fs_jpim1   ! vector opt.
427               ssha_e(ji,jj) = ( zsshb_e(ji,jj) - z2dt_e * ( zraur * ( emp(ji,jj)-rnf(ji,jj) ) + zhdiv(ji,jj) ) ) * tmask(ji,jj,1) 
428            END DO
429         END DO
430
431         !                                                !* after barotropic velocities (vorticity scheme dependent)
432         !                                                !  --------------------------- 
433         zwx(:,:) = e2u(:,:) * zun_e(:,:) * hu_e(:,:)     ! now_e transport
434         zwy(:,:) = e1v(:,:) * zvn_e(:,:) * hv_e(:,:)
435         !
436         IF( ln_dynvor_ene .OR. ln_dynvor_mix ) THEN      !==  energy conserving or mixed scheme  ==!
437            DO jj = 2, jpjm1
438               DO ji = fs_2, fs_jpim1   ! vector opt.
439                  ! surface pressure gradient
440                  IF( lk_vvl) THEN
441                     zu_spg = -grav * (  ( rhd(ji+1,jj ,1) + 1 ) * sshn_e(ji+1,jj  )    &
442                        &              - ( rhd(ji  ,jj ,1) + 1 ) * sshn_e(ji  ,jj  )  ) / e1u(ji,jj)
443                     zv_spg = -grav * (  ( rhd(ji ,jj+1,1) + 1 ) * sshn_e(ji  ,jj+1)    &
444                        &              - ( rhd(ji ,jj  ,1) + 1 ) * sshn_e(ji  ,jj  )  ) / e2v(ji,jj)
445                  ELSE
446                     zu_spg = -grav * ( sshn_e(ji+1,jj) - sshn_e(ji,jj) ) / e1u(ji,jj)
447                     zv_spg = -grav * ( sshn_e(ji,jj+1) - sshn_e(ji,jj) ) / e2v(ji,jj)
448                  ENDIF
449                  ! add tidal astronomical forcing
450                  IF ( ln_tide_pot .AND. lk_tide ) THEN
451                  zu_spg = zu_spg + grav * ( pot_astro(ji+1,jj) - pot_astro(ji,jj) ) / e1u(ji,jj)
452                  zv_spg = zv_spg + grav * ( pot_astro(ji,jj+1) - pot_astro(ji,jj) ) / e2v(ji,jj)
453                  ENDIF
454                  ! energy conserving formulation for planetary vorticity term
455                  zy1 = ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) ) / e1u(ji,jj)
456                  zy2 = ( zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) / e1u(ji,jj)
457                  zx1 = ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) ) / e2v(ji,jj)
458                  zx2 = ( zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) / e2v(ji,jj)
459                  zu_cor = z1_4 * ( ff(ji  ,jj-1) * zy1 + ff(ji,jj) * zy2 ) * hur_e(ji,jj)
460                  zv_cor =-z1_4 * ( ff(ji-1,jj  ) * zx1 + ff(ji,jj) * zx2 ) * hvr_e(ji,jj)
461                  ! after velocities with implicit bottom friction
462
463                  IF( ln_bfrimp ) THEN      ! implicit bottom friction
464                     !   A new method to implement the implicit bottom friction.
465                     !   H. Liu
466                     !   Sept 2011
467                     ua_e(ji,jj) = umask(ji,jj,1) * ( zub_e(ji,jj) +                                            &
468                      &                               z2dt_e * ( zu_cor + zu_spg + zu_sld + zu_asp )            &
469                      &                               / ( 1._wp      - z2dt_e * bfrua(ji,jj) * hur_e(ji,jj) ) )
470                     ua_e(ji,jj) = ( ua_e(ji,jj) + z2dt_e *   zua(ji,jj)  ) * umask(ji,jj,1)   
471                     !
472                     va_e(ji,jj) = vmask(ji,jj,1) * ( zvb_e(ji,jj) +                                            &
473                      &                               z2dt_e * ( zv_cor + zv_spg + zv_sld + zv_asp )            &
474                      &                               / ( 1._wp      - z2dt_e * bfrva(ji,jj) * hvr_e(ji,jj) ) )
475                     va_e(ji,jj) = ( va_e(ji,jj) + z2dt_e *   zva(ji,jj)  ) * vmask(ji,jj,1)   
476                     !
477                  ELSE
478                     ua_e(ji,jj) = ( zub_e(ji,jj) + z2dt_e * ( zu_cor + zu_spg + zu_sld + zu_asp + zua(ji,jj))) * umask(ji,jj,1)   &
479                      &           / ( 1._wp         - z2dt_e * bfrua(ji,jj) * hur_e(ji,jj) )
480                     va_e(ji,jj) = ( zvb_e(ji,jj) + z2dt_e * ( zv_cor + zv_spg + zv_sld + zv_asp + zva(ji,jj))) * vmask(ji,jj,1)   &
481                      &           / ( 1._wp         - z2dt_e * bfrva(ji,jj) * hvr_e(ji,jj) )
482                  ENDIF
483               END DO
484            END DO
485            !
486         ELSEIF ( ln_dynvor_ens ) THEN                    !==  enstrophy conserving scheme  ==!
487            DO jj = 2, jpjm1
488               DO ji = fs_2, fs_jpim1   ! vector opt.
489                   ! surface pressure gradient
490                  IF( lk_vvl) THEN
491                     zu_spg = -grav * (  ( rhd(ji+1,jj ,1) + 1 ) * sshn_e(ji+1,jj  )    &
492                        &              - ( rhd(ji  ,jj ,1) + 1 ) * sshn_e(ji  ,jj  )  ) / e1u(ji,jj)
493                     zv_spg = -grav * (  ( rhd(ji ,jj+1,1) + 1 ) * sshn_e(ji  ,jj+1)    &
494                        &              - ( rhd(ji ,jj  ,1) + 1 ) * sshn_e(ji  ,jj  )  ) / e2v(ji,jj)
495                  ELSE
496                     zu_spg = -grav * ( sshn_e(ji+1,jj) - sshn_e(ji,jj) ) / e1u(ji,jj)
497                     zv_spg = -grav * ( sshn_e(ji,jj+1) - sshn_e(ji,jj) ) / e2v(ji,jj)
498                  ENDIF
499                  ! add tidal astronomical forcing
500                  IF ( ln_tide_pot .AND. lk_tide ) THEN
501                  zu_spg = zu_spg + grav * ( pot_astro(ji+1,jj) - pot_astro(ji,jj) ) / e1u(ji,jj)
502                  zv_spg = zv_spg + grav * ( pot_astro(ji,jj+1) - pot_astro(ji,jj) ) / e2v(ji,jj)
503                  ENDIF
504                  ! enstrophy conserving formulation for planetary vorticity term
505                  zy1 =   z1_8 * ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) + zwy(ji,jj) + zwy(ji+1,jj  ) ) / e1u(ji,jj)
506                  zx1 = - z1_8 * ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) + zwx(ji,jj) + zwx(ji  ,jj+1) ) / e2v(ji,jj)
507                  zu_cor  = zy1 * ( ff(ji  ,jj-1) + ff(ji,jj) ) * hur_e(ji,jj)
508                  zv_cor  = zx1 * ( ff(ji-1,jj  ) + ff(ji,jj) ) * hvr_e(ji,jj)
509                  ! after velocities with implicit bottom friction
510                  IF( ln_bfrimp ) THEN
511                     !   A new method to implement the implicit bottom friction.
512                     !   H. Liu
513                     !   Sept 2011
514                     ua_e(ji,jj) = umask(ji,jj,1) * ( zub_e(ji,jj) +                                            &
515                      &                               z2dt_e * ( zu_cor + zu_spg + zu_sld + zu_asp )            &
516                      &                               / ( 1._wp      - z2dt_e * bfrua(ji,jj) * hur_e(ji,jj) ) )
517                     ua_e(ji,jj) = ( ua_e(ji,jj) + z2dt_e *   zua(ji,jj)  ) * umask(ji,jj,1)   
518                     !
519                     va_e(ji,jj) = vmask(ji,jj,1) * ( zvb_e(ji,jj) +                                            &
520                      &                               z2dt_e * ( zv_cor + zv_spg + zv_sld + zv_asp )            &
521                      &                               / ( 1._wp      - z2dt_e * bfrva(ji,jj) * hvr_e(ji,jj) ) )
522                     va_e(ji,jj) = ( va_e(ji,jj) + z2dt_e *   zva(ji,jj)  ) * vmask(ji,jj,1)   
523                     !
524                  ELSE
525                     ua_e(ji,jj) = ( zub_e(ji,jj) + z2dt_e * ( zu_cor + zu_spg + zu_sld + zu_asp + zua(ji,jj))) * umask(ji,jj,1)   &
526                     &            / ( 1._wp        - z2dt_e * bfrua(ji,jj) * hur_e(ji,jj) )
527                     va_e(ji,jj) = ( zvb_e(ji,jj) + z2dt_e * ( zv_cor + zv_spg + zv_sld + zv_asp + zva(ji,jj))) * vmask(ji,jj,1)   &
528                     &            / ( 1._wp        - z2dt_e * bfrva(ji,jj) * hvr_e(ji,jj) )
529                  ENDIF
530               END DO
531            END DO
532            !
533         ELSEIF ( ln_dynvor_een ) THEN                    !==  energy and enstrophy conserving scheme  ==!
534            DO jj = 2, jpjm1
535               DO ji = fs_2, fs_jpim1   ! vector opt.
536                  ! surface pressure gradient
537                  IF( lk_vvl) THEN
538                     zu_spg = -grav * (  ( rhd(ji+1,jj  ,1) + 1 ) * sshn_e(ji+1,jj  )    &
539                        &              - ( rhd(ji  ,jj  ,1) + 1 ) * sshn_e(ji  ,jj  )  ) / e1u(ji,jj)
540                     zv_spg = -grav * (  ( rhd(ji  ,jj+1,1) + 1 ) * sshn_e(ji  ,jj+1)    &
541                        &              - ( rhd(ji  ,jj  ,1) + 1 ) * sshn_e(ji  ,jj  )  ) / e2v(ji,jj)
542                  ELSE
543                     zu_spg = -grav * ( sshn_e(ji+1,jj) - sshn_e(ji,jj) ) / e1u(ji,jj)
544                     zv_spg = -grav * ( sshn_e(ji,jj+1) - sshn_e(ji,jj) ) / e2v(ji,jj)
545                  ENDIF
546                  ! add tidal astronomical forcing
547                  IF ( ln_tide_pot .AND. lk_tide ) THEN
548                  zu_spg = zu_spg + grav * ( pot_astro(ji+1,jj) - pot_astro(ji,jj) ) / e1u(ji,jj)
549                  zv_spg = zv_spg + grav * ( pot_astro(ji,jj+1) - pot_astro(ji,jj) ) / e2v(ji,jj)
550                  ENDIF
551                  ! energy/enstrophy conserving formulation for planetary vorticity term
552                  zu_cor = + z1_4 / e1u(ji,jj) * (  ftne(ji,jj  ) * zwy(ji  ,jj  ) + ftnw(ji+1,jj) * zwy(ji+1,jj  )   &
553                     &                           + ftse(ji,jj  ) * zwy(ji  ,jj-1) + ftsw(ji+1,jj) * zwy(ji+1,jj-1) ) * hur_e(ji,jj)
554                  zv_cor = - z1_4 / e2v(ji,jj) * (  ftsw(ji,jj+1) * zwx(ji-1,jj+1) + ftse(ji,jj+1) * zwx(ji  ,jj+1)   &
555                     &                           + ftnw(ji,jj  ) * zwx(ji-1,jj  ) + ftne(ji,jj  ) * zwx(ji  ,jj  ) ) * hvr_e(ji,jj)
556                  ! after velocities with implicit bottom friction
557                  IF( ln_bfrimp ) THEN
558                     !   A new method to implement the implicit bottom friction.
559                     !   H. Liu
560                     !   Sept 2011
561                     ua_e(ji,jj) = umask(ji,jj,1) * ( zub_e(ji,jj) +                                            &
562                      &                               z2dt_e * ( zu_cor + zu_spg + zu_sld + zu_asp )            &
563                      &                               / ( 1._wp      - z2dt_e * bfrua(ji,jj) * hur_e(ji,jj) ) )
564                     ua_e(ji,jj) = ( ua_e(ji,jj) + z2dt_e *   zua(ji,jj)  ) * umask(ji,jj,1)   
565                     !
566                     va_e(ji,jj) = vmask(ji,jj,1) * ( zvb_e(ji,jj) +                                            &
567                      &                               z2dt_e * ( zv_cor + zv_spg + zv_sld + zv_asp )            &
568                      &                               / ( 1._wp      - z2dt_e * bfrva(ji,jj) * hvr_e(ji,jj) ) )
569                     va_e(ji,jj) = ( va_e(ji,jj) + z2dt_e *   zva(ji,jj)  ) * vmask(ji,jj,1)   
570                     !
571                  ELSE
572                     ua_e(ji,jj) = ( zub_e(ji,jj) + z2dt_e * ( zu_cor + zu_spg + zu_sld + zu_asp + zua(ji,jj))) * umask(ji,jj,1)   &
573                     &            / ( 1._wp        - z2dt_e * bfrua(ji,jj) * hur_e(ji,jj) )
574                     va_e(ji,jj) = ( zvb_e(ji,jj) + z2dt_e * ( zv_cor + zv_spg + zv_sld + zv_asp + zva(ji,jj))) * vmask(ji,jj,1)   &
575                     &            / ( 1._wp        - z2dt_e * bfrva(ji,jj) * hvr_e(ji,jj) )
576                  ENDIF
577               END DO
578            END DO
579            !
580         ENDIF
581         !                                                     !* domain lateral boundary
582         !                                                     !  -----------------------
583
584                                                               ! OBC open boundaries
585         IF( lk_obc               )   CALL obc_fla_ts ( ua_e, va_e, sshn_e, ssha_e )
586
587                                                               ! BDY open boundaries
588#if defined key_bdy
589         pssh => sshn_e
590         phur => hur_e
591         phvr => hvr_e
592         pu2d => ua_e
593         pv2d => va_e
594
595         IF( lk_bdy )   CALL bdy_dyn2d( kt ) 
596#endif
597
598         !
599         CALL lbc_lnk( ua_e  , 'U', -1. )                      ! local domain boundaries
600         CALL lbc_lnk( va_e  , 'V', -1. )
601         CALL lbc_lnk( ssha_e, 'T',  1. )
602
603         zu_sum  (:,:) = zu_sum  (:,:) + ua_e  (:,:)           ! Sum over sub-time-steps
604         zv_sum  (:,:) = zv_sum  (:,:) + va_e  (:,:) 
605         zssh_sum(:,:) = zssh_sum(:,:) + ssha_e(:,:) 
606
607         !                                                !* Time filter and swap
608         !                                                !  --------------------
609         IF( jn == 1 ) THEN                                     ! Swap only (1st Euler time step)
610            zsshb_e(:,:) = sshn_e(:,:)
611            zub_e  (:,:) = zun_e (:,:)
612            zvb_e  (:,:) = zvn_e (:,:)
613            sshn_e (:,:) = ssha_e(:,:)
614            zun_e  (:,:) = ua_e  (:,:)
615            zvn_e  (:,:) = va_e  (:,:)
616         ELSE                                                   ! Swap + Filter
617            zsshb_e(:,:) = atfp * ( zsshb_e(:,:) + ssha_e(:,:) ) + atfp1 * sshn_e(:,:)
618            zub_e  (:,:) = atfp * ( zub_e  (:,:) + ua_e  (:,:) ) + atfp1 * zun_e (:,:)
619            zvb_e  (:,:) = atfp * ( zvb_e  (:,:) + va_e  (:,:) ) + atfp1 * zvn_e (:,:)
620            sshn_e (:,:) = ssha_e(:,:)
621            zun_e  (:,:) = ua_e  (:,:)
622            zvn_e  (:,:) = va_e  (:,:)
623         ENDIF
624
625         IF( lk_vvl ) THEN                                !* Update ocean depth (variable volume case only)
626            !                                             !  ------------------
627            DO jj = 1, jpjm1                                    ! Sea Surface Height at u- & v-points
628               DO ji = 1, fs_jpim1   ! Vector opt.
629                  zsshun_e(ji,jj) = 0.5_wp * umask(ji,jj,1) / ( e1u(ji,jj) * e2u(ji,jj) )       &
630                     &                     * ( e1t(ji  ,jj) * e2t(ji  ,jj) * sshn_e(ji  ,jj)    &
631                     &                     +   e1t(ji+1,jj) * e2t(ji+1,jj) * sshn_e(ji+1,jj) )
632                  zsshvn_e(ji,jj) = 0.5_wp * vmask(ji,jj,1) / ( e1v(ji,jj) * e2v(ji,jj) )       &
633                     &                     * ( e1t(ji,jj  ) * e2t(ji,jj  ) * sshn_e(ji,jj  )    &
634                     &                     +   e1t(ji,jj+1) * e2t(ji,jj+1) * sshn_e(ji,jj+1) )
635               END DO
636            END DO
637            CALL lbc_lnk( zsshun_e, 'U', 1. )                   ! lateral boundaries conditions
638            CALL lbc_lnk( zsshvn_e, 'V', 1. ) 
639            !
640            hu_e (:,:) = hu_0(:,:) + zsshun_e(:,:)              ! Ocean depth at U- and V-points
641            hv_e (:,:) = hv_0(:,:) + zsshvn_e(:,:)
642            hur_e(:,:) = umask(:,:,1) / ( hu_e(:,:) + 1._wp - umask(:,:,1) )
643            hvr_e(:,:) = vmask(:,:,1) / ( hv_e(:,:) + 1._wp - vmask(:,:,1) )
644            !
645         ENDIF
646         !                                                 ! ==================== !
647      END DO                                               !        end loop      !
648      !                                                    ! ==================== !
649
650#if defined key_obc
651      IF( lp_obc_east  )   sshfoe_b(:,:) = zcoef * sshfoe_b(:,:)     !!gm totally useless ?????
652      IF( lp_obc_west  )   sshfow_b(:,:) = zcoef * sshfow_b(:,:)
653      IF( lp_obc_north )   sshfon_b(:,:) = zcoef * sshfon_b(:,:)
654      IF( lp_obc_south )   sshfos_b(:,:) = zcoef * sshfos_b(:,:)
655#endif
656
657      ! -----------------------------------------------------------------------------
658      ! Phase 3. update the general trend with the barotropic trend
659      ! -----------------------------------------------------------------------------
660      !
661      !                                   !* Time average ==> after barotropic u, v, ssh
662      zcoef =  1._wp / ( 2 * nn_baro  + 1 ) 
663      zu_sum(:,:) = zcoef * zu_sum  (:,:) 
664      zv_sum(:,:) = zcoef * zv_sum  (:,:) 
665      !
666      !                                   !* update the general momentum trend
667      DO jk=1,jpkm1
668         ua(:,:,jk) = ua(:,:,jk) + ( zu_sum(:,:) - ub_b(:,:) ) * z1_2dt_b
669         va(:,:,jk) = va(:,:,jk) + ( zv_sum(:,:) - vb_b(:,:) ) * z1_2dt_b
670      END DO
671      un_b  (:,:) =  zu_sum(:,:) 
672      vn_b  (:,:) =  zv_sum(:,:) 
673      sshn_b(:,:) = zcoef * zssh_sum(:,:) 
674      !
675      !                                   !* write time-spliting arrays in the restart
676      IF( lrst_oce )   CALL ts_rst( kt, 'WRITE' )
677      !
678      CALL wrk_dealloc( jpi, jpj, zsshun_e, zsshvn_e, zsshb_e, zssh_sum, zhdiv     )
679      CALL wrk_dealloc( jpi, jpj, zua, zva, zun, zvn, zun_e, zvn_e, zub_e, zvb_e   )
680      CALL wrk_dealloc( jpi, jpj, zcu, zcv, zwx, zwy, zbfru, zbfrv, zu_sum, zv_sum )
681      CALL wrk_dealloc( jpi, jpj, zhu_b, zhv_b                                     )
682      !
683      IF( nn_timing == 1 )  CALL timing_stop('dyn_spg_ts')
684      !
685   END SUBROUTINE dyn_spg_ts
686
687
688   SUBROUTINE ts_rst( kt, cdrw )
689      !!---------------------------------------------------------------------
690      !!                   ***  ROUTINE ts_rst  ***
691      !!
692      !! ** Purpose : Read or write time-splitting arrays in restart file
693      !!----------------------------------------------------------------------
694      INTEGER         , INTENT(in) ::   kt         ! ocean time-step
695      CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
696      !
697      REAL(wp), POINTER, DIMENSION(:,:) :: zzhu_b, zzhv_b
698      INTEGER ::  ji, jk        ! dummy loop indices
699      !!----------------------------------------------------------------------
700      !
701      IF( TRIM(cdrw) == 'READ' ) THEN
702         IF( iom_varid( numror, 'un_b', ldstop = .FALSE. ) > 0 ) THEN
703            CALL iom_get( numror, jpdom_autoglo, 'un_b'  , un_b  (:,:) )   ! external velocity issued
704            CALL iom_get( numror, jpdom_autoglo, 'vn_b'  , vn_b  (:,:) )   ! from barotropic loop
705         ELSE
706            un_b (:,:) = 0._wp
707            vn_b (:,:) = 0._wp
708            ! vertical sum
709            IF( lk_vopt_loop ) THEN          ! vector opt., forced unroll
710               DO jk = 1, jpkm1
711                  DO ji = 1, jpij
712                     un_b(ji,1) = un_b(ji,1) + fse3u_n(ji,1,jk) * un(ji,1,jk)
713                     vn_b(ji,1) = vn_b(ji,1) + fse3v_n(ji,1,jk) * vn(ji,1,jk)
714                  END DO
715               END DO
716            ELSE                             ! No  vector opt.
717               DO jk = 1, jpkm1
718                  un_b(:,:) = un_b(:,:) + fse3u_n(:,:,jk) * un(:,:,jk)
719                  vn_b(:,:) = vn_b(:,:) + fse3v_n(:,:,jk) * vn(:,:,jk)
720               END DO
721            ENDIF
722            un_b (:,:) = un_b(:,:) * hur(:,:)
723            vn_b (:,:) = vn_b(:,:) * hvr(:,:)
724         ENDIF
725
726         ! Vertically integrated velocity (before)
727         IF (neuler/=0) THEN
728            ub_b (:,:) = 0._wp
729            vb_b (:,:) = 0._wp
730
731            ! vertical sum
732            IF( lk_vopt_loop ) THEN          ! vector opt., forced unroll
733               DO jk = 1, jpkm1
734                  DO ji = 1, jpij
735                     ub_b(ji,1) = ub_b(ji,1) + fse3u_b(ji,1,jk) * ub(ji,1,jk)
736                     vb_b(ji,1) = vb_b(ji,1) + fse3v_b(ji,1,jk) * vb(ji,1,jk)
737                  END DO
738               END DO
739            ELSE                             ! No  vector opt.
740               DO jk = 1, jpkm1
741                  ub_b(:,:) = ub_b(:,:) + fse3u_b(:,:,jk) * ub(:,:,jk)
742                  vb_b(:,:) = vb_b(:,:) + fse3v_b(:,:,jk) * vb(:,:,jk)
743               END DO
744            ENDIF
745
746            IF( lk_vvl ) THEN
747               CALL wrk_alloc( jpi, jpj, zzhu_b, zzhv_b )
748               ub_b (:,:) = 0.
749               vb_b (:,:) = 0.
750               zzhu_b(:,:) = 0.
751               zzhv_b(:,:) = 0.
752               ! vertical sum
753               IF( lk_vopt_loop ) THEN          ! vector opt., forced unroll
754                  DO jk = 1, jpkm1
755                     DO ji = 1, jpij
756                        ub_b  (ji,1) = ub_b (ji,1)  + fse3u_b(ji,1,jk) * ub   (ji,1,jk)
757                        vb_b  (ji,1) = vb_b (ji,1)  + fse3v_b(ji,1,jk) * vb   (ji,1,jk)
758                        zzhu_b(ji,1) = zzhu_b(ji,1) + fse3u_b(ji,1,jk) * umask(ji,1,jk)
759                        zzhv_b(ji,1) = zzhv_b(ji,1) + fse3v_b(ji,1,jk) * vmask(ji,1,jk)
760                     END DO
761                  END DO
762               ELSE                             ! No  vector opt.
763                  DO jk = 1, jpkm1
764                     ub_b  (:,:) = ub_b  (:,:) + fse3u_b(:,:,jk) * ub   (:,:,jk)
765                     vb_b  (:,:) = vb_b  (:,:) + fse3v_b(:,:,jk) * vb   (:,:,jk)
766                     zzhu_b(:,:) = zzhu_b(:,:) + fse3u_b(:,:,jk) * umask(:,:,jk)
767                     zzhv_b(:,:) = zzhv_b(:,:) + fse3v_b(:,:,jk) * vmask(:,:,jk)
768                  END DO
769               ENDIF
770               ub_b(:,:) = ub_b(:,:) / ( zzhu_b(:,:) + 1. - umask(:,:,1) )
771               vb_b(:,:) = vb_b(:,:) / ( zzhv_b(:,:) + 1. - vmask(:,:,1) )
772               CALL wrk_dealloc( jpi, jpj, zzhu_b, zzhv_b )
773            ELSE             
774               ub_b (:,:) = 0.e0
775               vb_b (:,:) = 0.e0
776               ! vertical sum
777               IF( lk_vopt_loop ) THEN          ! vector opt., forced unroll
778                  DO jk = 1, jpkm1
779                     DO ji = 1, jpij
780                        ub_b(ji,1) = ub_b(ji,1) + fse3u_b(ji,1,jk) * ub(ji,1,jk)
781                        vb_b(ji,1) = vb_b(ji,1) + fse3v_b(ji,1,jk) * vb(ji,1,jk)
782                     END DO
783                  END DO
784               ELSE                             ! No  vector opt.
785                  DO jk = 1, jpkm1
786                     ub_b(:,:) = ub_b(:,:) + fse3u_b(:,:,jk) * ub(:,:,jk)
787                     vb_b(:,:) = vb_b(:,:) + fse3v_b(:,:,jk) * vb(:,:,jk)
788                  END DO
789               ENDIF
790               ub_b(:,:) = ub_b(:,:) * hur(:,:)
791               vb_b(:,:) = vb_b(:,:) * hvr(:,:)
792            ENDIF
793         ELSE                                 ! neuler==0
794            ub_b (:,:) = un_b (:,:)
795            vb_b (:,:) = vn_b (:,:)
796         ENDIF
797
798         IF( iom_varid( numror, 'sshn_b', ldstop = .FALSE. ) > 0 ) THEN
799            CALL iom_get( numror, jpdom_autoglo, 'sshn_b' , sshn_b (:,:) )   ! filtered ssh
800         ELSE
801            sshn_b(:,:) = sshb(:,:)   ! if not in restart set previous time mean to current baroclinic before value   
802         ENDIF
803      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN
804         CALL iom_rstput( kt, nitrst, numrow, 'un_b'   , un_b  (:,:) )   ! external velocity and ssh
805         CALL iom_rstput( kt, nitrst, numrow, 'vn_b'   , vn_b  (:,:) )   ! issued from barotropic loop
806         CALL iom_rstput( kt, nitrst, numrow, 'sshn_b' , sshn_b(:,:) )   !
807      ENDIF
808      !
809   END SUBROUTINE ts_rst
810
811#else
812   !!----------------------------------------------------------------------
813   !!   Default case :   Empty module   No standart free surface cst volume
814   !!----------------------------------------------------------------------
815CONTAINS
816   INTEGER FUNCTION dyn_spg_ts_alloc()    ! Dummy function
817      dyn_spg_ts_alloc = 0
818   END FUNCTION dyn_spg_ts_alloc
819   SUBROUTINE dyn_spg_ts( kt )            ! Empty routine
820      INTEGER, INTENT(in) :: kt
821      WRITE(*,*) 'dyn_spg_ts: You should not have seen this print! error?', kt
822   END SUBROUTINE dyn_spg_ts
823   SUBROUTINE ts_rst( kt, cdrw )          ! Empty routine
824      INTEGER         , INTENT(in) ::   kt         ! ocean time-step
825      CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
826      WRITE(*,*) 'ts_rst    : You should not have seen this print! error?', kt, cdrw
827   END SUBROUTINE ts_rst   
828#endif
829   
830   !!======================================================================
831END MODULE dynspg_ts
Note: See TracBrowser for help on using the repository browser.