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 @ 3865

Last change on this file since 3865 was 3865, checked in by acc, 11 years ago

Branch 2013/dev_r3858_NOC_ZTC, #863. Nearly complete port of 2011/dev_r2739_LOCEAN8_ZTC development branch into v3.5aplha base. Compiles and runs but currently unstable after 8 timesteps with ORCA2_LIM reference configuration.

  • Property svn:keywords set to Id
File size: 45.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#if defined key_vvl
210               ub_b(ji,jj) = ub_b(ji,jj) + fse3u_b(ji,jj,jk)* ub(ji,jj,jk)   *umask(ji,jj,jk) 
211               vb_b(ji,jj) = vb_b(ji,jj) + fse3v_b(ji,jj,jk)* vb(ji,jj,jk)   *vmask(ji,jj,jk)
212#else
213               ub_b(ji,jj) = ub_b(ji,jj) + fse3u_0(ji,jj,jk) * ub(ji,jj,jk)  * umask(ji,jj,jk)
214               vb_b(ji,jj) = vb_b(ji,jj) + fse3v_0(ji,jj,jk) * vb(ji,jj,jk)  * vmask(ji,jj,jk)
215#endif
216            END DO
217         END DO
218      END DO
219
220      ! before inverse water column height at u- and v- points
221      IF( lk_vvl ) THEN
222         zhu_b(:,:) = 0.
223         zhv_b(:,:) = 0.
224         DO jk = 1, jpk
225            zhu_b(:,:) = zhu_b(:,:) + fse3u_b(:,:,jk) * umask(:,:,jk)
226            zhv_b(:,:) = zhv_b(:,:) + fse3v_b(:,:,jk) * vmask(:,:,jk)
227         END DO
228         zhu_b(:,:) = umask(:,:,1) / ( zhu_b(:,:) + 1. - umask(:,:,1) )
229         zhv_b(:,:) = vmask(:,:,1) / ( zhv_b(:,:) + 1. - vmask(:,:,1) )
230      ELSE
231         zhu_b(:,:) = hur(:,:)
232         zhv_b(:,:) = hvr(:,:)
233      ENDIF
234
235      !                                   !* baroclinic momentum trend (remove the vertical mean trend)
236      DO jk = 1, jpkm1                    ! --------------------------
237         DO jj = 2, jpjm1
238            DO ji = fs_2, fs_jpim1   ! vector opt.
239               ua(ji,jj,jk) = ua(ji,jj,jk) - zua(ji,jj) * hur(ji,jj)
240               va(ji,jj,jk) = va(ji,jj,jk) - zva(ji,jj) * hvr(ji,jj)
241            END DO
242         END DO
243      END DO
244
245      !                                   !* barotropic Coriolis trends * H (vorticity scheme dependent)
246      !                                   ! ---------------------------====
247      zwx(:,:) = zun(:,:) * e2u(:,:)                   ! now transport
248      zwy(:,:) = zvn(:,:) * e1v(:,:)
249      !
250      IF( ln_dynvor_ene .OR. ln_dynvor_mix ) THEN      ! energy conserving or mixed scheme
251         DO jj = 2, jpjm1
252            DO ji = fs_2, fs_jpim1   ! vector opt.
253               zy1 = ( zwy(ji,jj-1) + zwy(ji+1,jj-1) ) / e1u(ji,jj)
254               zy2 = ( zwy(ji,jj  ) + zwy(ji+1,jj  ) ) / e1u(ji,jj)
255               zx1 = ( zwx(ji-1,jj) + zwx(ji-1,jj+1) ) / e2v(ji,jj)
256               zx2 = ( zwx(ji  ,jj) + zwx(ji  ,jj+1) ) / e2v(ji,jj)
257               ! energy conserving formulation for planetary vorticity term
258               zcu(ji,jj) = z1_4 * ( ff(ji  ,jj-1) * zy1 + ff(ji,jj) * zy2 )
259               zcv(ji,jj) =-z1_4 * ( ff(ji-1,jj  ) * zx1 + ff(ji,jj) * zx2 )
260            END DO
261         END DO
262         !
263      ELSEIF ( ln_dynvor_ens ) THEN                    ! enstrophy conserving scheme
264         DO jj = 2, jpjm1
265            DO ji = fs_2, fs_jpim1   ! vector opt.
266               zy1 =   z1_8 * ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) + zwy(ji,jj) + zwy(ji+1,jj  ) ) / e1u(ji,jj)
267               zx1 = - z1_8 * ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) + zwx(ji,jj) + zwx(ji  ,jj+1) ) / e2v(ji,jj)
268               zcu(ji,jj)  = zy1 * ( ff(ji  ,jj-1) + ff(ji,jj) )
269               zcv(ji,jj)  = zx1 * ( ff(ji-1,jj  ) + ff(ji,jj) )
270            END DO
271         END DO
272         !
273      ELSEIF ( ln_dynvor_een ) THEN                    ! enstrophy and energy conserving scheme
274         DO jj = 2, jpjm1
275            DO ji = fs_2, fs_jpim1   ! vector opt.
276               zcu(ji,jj) = + z1_4 / e1u(ji,jj) * (  ftne(ji,jj  ) * zwy(ji  ,jj  ) + ftnw(ji+1,jj) * zwy(ji+1,jj  )   &
277                  &                                + ftse(ji,jj  ) * zwy(ji  ,jj-1) + ftsw(ji+1,jj) * zwy(ji+1,jj-1) )
278               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)   &
279                  &                                + ftnw(ji,jj  ) * zwx(ji-1,jj  ) + ftne(ji,jj  ) * zwx(ji  ,jj  ) )
280            END DO
281         END DO
282         !
283      ENDIF
284
285      !                                   !* Right-Hand-Side of the barotropic momentum equation
286      !                                   ! ----------------------------------------------------
287      IF( lk_vvl ) THEN                         ! Variable volume : remove both Coriolis and Surface pressure gradient
288         DO jj = 2, jpjm1 
289            DO ji = fs_2, fs_jpim1   ! vector opt.
290               zcu(ji,jj) = zcu(ji,jj) - grav * (  ( rhd(ji+1,jj  ,1) + 1 ) * sshn(ji+1,jj  )    &
291                  &                              - ( rhd(ji  ,jj  ,1) + 1 ) * sshn(ji  ,jj  )  ) * hu(ji,jj) / e1u(ji,jj)
292               zcv(ji,jj) = zcv(ji,jj) - grav * (  ( rhd(ji  ,jj+1,1) + 1 ) * sshn(ji  ,jj+1)    &
293                  &                              - ( rhd(ji  ,jj  ,1) + 1 ) * sshn(ji  ,jj  )  ) * hv(ji,jj) / e2v(ji,jj)
294            END DO
295         END DO
296      ENDIF
297
298      DO jj = 2, jpjm1                             ! Remove coriolis term (and possibly spg) from barotropic trend
299         DO ji = fs_2, fs_jpim1
300             zua(ji,jj) = zua(ji,jj) - zcu(ji,jj)
301             zva(ji,jj) = zva(ji,jj) - zcv(ji,jj)
302          END DO
303      END DO
304
305                   
306      !                                             ! Remove barotropic contribution of bottom friction
307      !                                             ! from the barotropic transport trend
308      zcoef = -1._wp * z1_2dt_b
309
310      IF(ln_bfrimp) THEN
311      !                                   ! Remove the bottom stress trend from 3-D sea surface level gradient
312      !                                   ! and Coriolis forcing in case of 3D semi-implicit bottom friction
313        DO jj = 2, jpjm1         
314           DO ji = fs_2, fs_jpim1
315              ikbu = mbku(ji,jj)
316              ikbv = mbkv(ji,jj)
317              ua_btm = zcu(ji,jj) * z2dt_bf * hur(ji,jj) * umask (ji,jj,ikbu)
318              va_btm = zcv(ji,jj) * z2dt_bf * hvr(ji,jj) * vmask (ji,jj,ikbv)
319
320              zua(ji,jj) = zua(ji,jj) - bfrua(ji,jj) * ua_btm
321              zva(ji,jj) = zva(ji,jj) - bfrva(ji,jj) * va_btm
322           END DO
323        END DO
324
325      ELSE
326
327# if defined key_vectopt_loop
328        DO jj = 1, 1
329           DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
330# else
331        DO jj = 2, jpjm1
332           DO ji = 2, jpim1
333# endif
334            ! Apply stability criteria for bottom friction
335            !RBbug for vvl and external mode we may need to use varying fse3
336            !!gm  Rq: the bottom e3 present the smallest variation, the use of e3u_0 is not a big approx.
337              zbfru(ji,jj) = MAX(  bfrua(ji,jj) , fse3u(ji,jj,mbku(ji,jj)) * zcoef  )
338              zbfrv(ji,jj) = MAX(  bfrva(ji,jj) , fse3v(ji,jj,mbkv(ji,jj)) * zcoef  )
339           END DO
340        END DO
341
342        IF( lk_vvl ) THEN
343           DO jj = 2, jpjm1
344              DO ji = fs_2, fs_jpim1   ! vector opt.
345                 zua(ji,jj) = zua(ji,jj) - zbfru(ji,jj) * ub_b(ji,jj)   &
346                    &       / ( hu_0(ji,jj) + sshu_b(ji,jj) + 1._wp - umask(ji,jj,1) )
347                 zva(ji,jj) = zva(ji,jj) - zbfrv(ji,jj) * vb_b(ji,jj)   &
348                    &       / ( hv_0(ji,jj) + sshv_b(ji,jj) + 1._wp - vmask(ji,jj,1) )
349              END DO
350           END DO
351        ELSE
352           DO jj = 2, jpjm1
353              DO ji = fs_2, fs_jpim1   ! vector opt.
354                 zua(ji,jj) = zua(ji,jj) - zbfru(ji,jj) * ub_b(ji,jj) * hur(ji,jj)
355                 zva(ji,jj) = zva(ji,jj) - zbfrv(ji,jj) * vb_b(ji,jj) * hvr(ji,jj)
356              END DO
357           END DO
358        ENDIF
359      END IF    ! end (ln_bfrimp)
360
361                   
362      !                                   !* d/dt(Ua), Ub, Vn (Vertical mean velocity)
363      !                                   ! --------------------------
364      zua(:,:) = zua(:,:) * hur(:,:)
365      zva(:,:) = zva(:,:) * hvr(:,:)
366      !
367      IF( lk_vvl ) THEN
368         ub_b(:,:) = ub_b(:,:) * umask(:,:,1) / ( hu_0(:,:) + sshu_b(:,:) + 1._wp - umask(:,:,1) )
369         vb_b(:,:) = vb_b(:,:) * vmask(:,:,1) / ( hv_0(:,:) + sshv_b(:,:) + 1._wp - vmask(:,:,1) )
370      ELSE
371         ub_b(:,:) = ub_b(:,:) * hur(:,:)
372         vb_b(:,:) = vb_b(:,:) * hvr(:,:)
373      ENDIF
374      ub_b(:,:) = ub_b(:,:) * zhu_b(:,:)
375      vb_b(:,:) = vb_b(:,:) * zhv_b(:,:)
376
377      ! -----------------------------------------------------------------------
378      !  Phase 2 : Integration of the barotropic equations with time splitting
379      ! -----------------------------------------------------------------------
380      !
381      !                                             ! ==================== !
382      !                                             !    Initialisations   !
383      !                                             ! ==================== !
384      icycle = 2  * nn_baro            ! Number of barotropic sub time-step
385     
386      !                                ! Start from NOW field
387      hu_e   (:,:) = hu   (:,:)            ! ocean depth at u- and v-points
388      hv_e   (:,:) = hv   (:,:) 
389      hur_e  (:,:) = hur  (:,:)            ! ocean depth inverted at u- and v-points
390      hvr_e  (:,:) = hvr  (:,:)
391!RBbug     zsshb_e(:,:) = sshn (:,:) 
392      zsshb_e(:,:) = sshn_b(:,:)           ! sea surface height (before and now)
393      sshn_e (:,:) = sshn (:,:)
394     
395      zun_e  (:,:) = un_b (:,:)            ! barotropic velocity (external)
396      zvn_e  (:,:) = vn_b (:,:)
397      zub_e  (:,:) = un_b (:,:)
398      zvb_e  (:,:) = vn_b (:,:)
399
400      zu_sum  (:,:) = un_b (:,:)           ! summation
401      zv_sum  (:,:) = vn_b (:,:)
402      zssh_sum(:,:) = sshn (:,:)
403
404#if defined key_obc
405      ! set ssh corrections to 0
406      ! ssh corrections are applied to normal velocities (Flather's algorithm) and averaged over the barotropic loop
407      IF( lp_obc_east  )   sshfoe_b(:,:) = 0._wp
408      IF( lp_obc_west  )   sshfow_b(:,:) = 0._wp
409      IF( lp_obc_south )   sshfos_b(:,:) = 0._wp
410      IF( lp_obc_north )   sshfon_b(:,:) = 0._wp
411#endif
412
413      !                                             ! ==================== !
414      DO jn = 1, icycle                             !  sub-time-step loop  ! (from NOW to AFTER+1)
415         !                                          ! ==================== !
416         z2dt_e = 2. * ( rdt / nn_baro )
417         IF( jn == 1 )   z2dt_e = rdt / nn_baro
418
419         !                                                !* Update the forcing (BDY and tides)
420         !                                                !  ------------------
421         IF( lk_obc )   CALL obc_dta_bt ( kt, jn   )
422         IF( lk_bdy )   CALL bdy_dta ( kt, jit=jn, time_offset=+1 )
423         IF ( ln_tide_pot .AND. lk_tide) CALL upd_tide( kt, jn )
424
425         !                                                !* after ssh_e
426         !                                                !  -----------
427         DO jj = 2, jpjm1                                 ! Horizontal divergence of barotropic transports
428            DO ji = fs_2, fs_jpim1   ! vector opt.
429               zhdiv(ji,jj) = (   e2u(ji  ,jj) * zun_e(ji  ,jj) * hu_e(ji  ,jj)     &
430                  &             - e2u(ji-1,jj) * zun_e(ji-1,jj) * hu_e(ji-1,jj)     &
431                  &             + e1v(ji,jj  ) * zvn_e(ji,jj  ) * hv_e(ji,jj  )     &
432                  &             - e1v(ji,jj-1) * zvn_e(ji,jj-1) * hv_e(ji,jj-1)   ) / ( e1t(ji,jj) * e2t(ji,jj) )
433            END DO
434         END DO
435         !
436#if defined key_obc
437         !                                                     ! OBC : zhdiv must be zero behind the open boundary
438!!  mpp remark: The zeroing of hdiv can probably be extended to 1->jpi/jpj for the correct row/column
439         IF( lp_obc_east  )   zhdiv(nie0p1:nie1p1,nje0  :nje1  ) = 0._wp      ! east
440         IF( lp_obc_west  )   zhdiv(niw0  :niw1  ,njw0  :njw1  ) = 0._wp      ! west
441         IF( lp_obc_north )   zhdiv(nin0  :nin1  ,njn0p1:njn1p1) = 0._wp      ! north
442         IF( lp_obc_south )   zhdiv(nis0  :nis1  ,njs0  :njs1  ) = 0._wp      ! south
443#endif
444#if defined key_bdy
445         zhdiv(:,:) = zhdiv(:,:) * bdytmask(:,:)               ! BDY mask
446#endif
447         !
448         DO jj = 2, jpjm1                                      ! leap-frog on ssh_e
449            DO ji = fs_2, fs_jpim1   ! vector opt.
450               ssha_e(ji,jj) = ( zsshb_e(ji,jj) - z2dt_e * ( zraur * ( emp(ji,jj)-rnf(ji,jj) ) + zhdiv(ji,jj) ) ) * tmask(ji,jj,1) 
451            END DO
452         END DO
453
454         !                                                !* after barotropic velocities (vorticity scheme dependent)
455         !                                                !  --------------------------- 
456         zwx(:,:) = e2u(:,:) * zun_e(:,:) * hu_e(:,:)     ! now_e transport
457         zwy(:,:) = e1v(:,:) * zvn_e(:,:) * hv_e(:,:)
458         !
459         IF( ln_dynvor_ene .OR. ln_dynvor_mix ) THEN      !==  energy conserving or mixed scheme  ==!
460            DO jj = 2, jpjm1
461               DO ji = fs_2, fs_jpim1   ! vector opt.
462                  ! surface pressure gradient
463                  IF( lk_vvl) THEN
464                     zu_spg = -grav * (  ( rhd(ji+1,jj ,1) + 1 ) * sshn_e(ji+1,jj  )    &
465                        &              - ( rhd(ji  ,jj ,1) + 1 ) * sshn_e(ji  ,jj  )  ) / e1u(ji,jj)
466                     zv_spg = -grav * (  ( rhd(ji ,jj+1,1) + 1 ) * sshn_e(ji  ,jj+1)    &
467                        &              - ( rhd(ji ,jj  ,1) + 1 ) * sshn_e(ji  ,jj  )  ) / e2v(ji,jj)
468                  ELSE
469                     zu_spg = -grav * ( sshn_e(ji+1,jj) - sshn_e(ji,jj) ) / e1u(ji,jj)
470                     zv_spg = -grav * ( sshn_e(ji,jj+1) - sshn_e(ji,jj) ) / e2v(ji,jj)
471                  ENDIF
472                  ! add tidal astronomical forcing
473                  IF ( ln_tide_pot .AND. lk_tide ) THEN
474                  zu_spg = zu_spg + grav * ( pot_astro(ji+1,jj) - pot_astro(ji,jj) ) / e1u(ji,jj)
475                  zv_spg = zv_spg + grav * ( pot_astro(ji,jj+1) - pot_astro(ji,jj) ) / e2v(ji,jj)
476                  ENDIF
477                  ! energy conserving formulation for planetary vorticity term
478                  zy1 = ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) ) / e1u(ji,jj)
479                  zy2 = ( zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) / e1u(ji,jj)
480                  zx1 = ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) ) / e2v(ji,jj)
481                  zx2 = ( zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) / e2v(ji,jj)
482                  zu_cor = z1_4 * ( ff(ji  ,jj-1) * zy1 + ff(ji,jj) * zy2 ) * hur_e(ji,jj)
483                  zv_cor =-z1_4 * ( ff(ji-1,jj  ) * zx1 + ff(ji,jj) * zx2 ) * hvr_e(ji,jj)
484                  ! after velocities with implicit bottom friction
485
486                  IF( ln_bfrimp ) THEN      ! implicit bottom friction
487                     !   A new method to implement the implicit bottom friction.
488                     !   H. Liu
489                     !   Sept 2011
490                     ua_e(ji,jj) = umask(ji,jj,1) * ( zub_e(ji,jj) +                                            &
491                      &                               z2dt_e * ( zu_cor + zu_spg + zu_sld + zu_asp )            &
492                      &                               / ( 1._wp      - z2dt_e * bfrua(ji,jj) * hur_e(ji,jj) ) )
493                     ua_e(ji,jj) = ( ua_e(ji,jj) + z2dt_e *   zua(ji,jj)  ) * umask(ji,jj,1)   
494                     !
495                     va_e(ji,jj) = vmask(ji,jj,1) * ( zvb_e(ji,jj) +                                            &
496                      &                               z2dt_e * ( zv_cor + zv_spg + zv_sld + zv_asp )            &
497                      &                               / ( 1._wp      - z2dt_e * bfrva(ji,jj) * hvr_e(ji,jj) ) )
498                     va_e(ji,jj) = ( va_e(ji,jj) + z2dt_e *   zva(ji,jj)  ) * vmask(ji,jj,1)   
499                     !
500                  ELSE
501                     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)   &
502                      &           / ( 1._wp         - z2dt_e * bfrua(ji,jj) * hur_e(ji,jj) )
503                     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)   &
504                      &           / ( 1._wp         - z2dt_e * bfrva(ji,jj) * hvr_e(ji,jj) )
505                  ENDIF
506               END DO
507            END DO
508            !
509         ELSEIF ( ln_dynvor_ens ) THEN                    !==  enstrophy conserving scheme  ==!
510            DO jj = 2, jpjm1
511               DO ji = fs_2, fs_jpim1   ! vector opt.
512                   ! surface pressure gradient
513                  IF( lk_vvl) THEN
514                     zu_spg = -grav * (  ( rhd(ji+1,jj ,1) + 1 ) * sshn_e(ji+1,jj  )    &
515                        &              - ( rhd(ji  ,jj ,1) + 1 ) * sshn_e(ji  ,jj  )  ) / e1u(ji,jj)
516                     zv_spg = -grav * (  ( rhd(ji ,jj+1,1) + 1 ) * sshn_e(ji  ,jj+1)    &
517                        &              - ( rhd(ji ,jj  ,1) + 1 ) * sshn_e(ji  ,jj  )  ) / e2v(ji,jj)
518                  ELSE
519                     zu_spg = -grav * ( sshn_e(ji+1,jj) - sshn_e(ji,jj) ) / e1u(ji,jj)
520                     zv_spg = -grav * ( sshn_e(ji,jj+1) - sshn_e(ji,jj) ) / e2v(ji,jj)
521                  ENDIF
522                  ! add tidal astronomical forcing
523                  IF ( ln_tide_pot .AND. lk_tide ) THEN
524                  zu_spg = zu_spg + grav * ( pot_astro(ji+1,jj) - pot_astro(ji,jj) ) / e1u(ji,jj)
525                  zv_spg = zv_spg + grav * ( pot_astro(ji,jj+1) - pot_astro(ji,jj) ) / e2v(ji,jj)
526                  ENDIF
527                  ! enstrophy conserving formulation for planetary vorticity term
528                  zy1 =   z1_8 * ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) + zwy(ji,jj) + zwy(ji+1,jj  ) ) / e1u(ji,jj)
529                  zx1 = - z1_8 * ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) + zwx(ji,jj) + zwx(ji  ,jj+1) ) / e2v(ji,jj)
530                  zu_cor  = zy1 * ( ff(ji  ,jj-1) + ff(ji,jj) ) * hur_e(ji,jj)
531                  zv_cor  = zx1 * ( ff(ji-1,jj  ) + ff(ji,jj) ) * hvr_e(ji,jj)
532                  ! after velocities with implicit bottom friction
533                  IF( ln_bfrimp ) THEN
534                     !   A new method to implement the implicit bottom friction.
535                     !   H. Liu
536                     !   Sept 2011
537                     ua_e(ji,jj) = umask(ji,jj,1) * ( zub_e(ji,jj) +                                            &
538                      &                               z2dt_e * ( zu_cor + zu_spg + zu_sld + zu_asp )            &
539                      &                               / ( 1._wp      - z2dt_e * bfrua(ji,jj) * hur_e(ji,jj) ) )
540                     ua_e(ji,jj) = ( ua_e(ji,jj) + z2dt_e *   zua(ji,jj)  ) * umask(ji,jj,1)   
541                     !
542                     va_e(ji,jj) = vmask(ji,jj,1) * ( zvb_e(ji,jj) +                                            &
543                      &                               z2dt_e * ( zv_cor + zv_spg + zv_sld + zv_asp )            &
544                      &                               / ( 1._wp      - z2dt_e * bfrva(ji,jj) * hvr_e(ji,jj) ) )
545                     va_e(ji,jj) = ( va_e(ji,jj) + z2dt_e *   zva(ji,jj)  ) * vmask(ji,jj,1)   
546                     !
547                  ELSE
548                     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)   &
549                     &            / ( 1._wp        - z2dt_e * bfrua(ji,jj) * hur_e(ji,jj) )
550                     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)   &
551                     &            / ( 1._wp        - z2dt_e * bfrva(ji,jj) * hvr_e(ji,jj) )
552                  ENDIF
553               END DO
554            END DO
555            !
556         ELSEIF ( ln_dynvor_een ) THEN                    !==  energy and enstrophy conserving scheme  ==!
557            DO jj = 2, jpjm1
558               DO ji = fs_2, fs_jpim1   ! vector opt.
559                  ! surface pressure gradient
560                  IF( lk_vvl) THEN
561                     zu_spg = -grav * (  ( rhd(ji+1,jj  ,1) + 1 ) * sshn_e(ji+1,jj  )    &
562                        &              - ( rhd(ji  ,jj  ,1) + 1 ) * sshn_e(ji  ,jj  )  ) / e1u(ji,jj)
563                     zv_spg = -grav * (  ( rhd(ji  ,jj+1,1) + 1 ) * sshn_e(ji  ,jj+1)    &
564                        &              - ( rhd(ji  ,jj  ,1) + 1 ) * sshn_e(ji  ,jj  )  ) / e2v(ji,jj)
565                  ELSE
566                     zu_spg = -grav * ( sshn_e(ji+1,jj) - sshn_e(ji,jj) ) / e1u(ji,jj)
567                     zv_spg = -grav * ( sshn_e(ji,jj+1) - sshn_e(ji,jj) ) / e2v(ji,jj)
568                  ENDIF
569                  ! add tidal astronomical forcing
570                  IF ( ln_tide_pot .AND. lk_tide ) THEN
571                  zu_spg = zu_spg + grav * ( pot_astro(ji+1,jj) - pot_astro(ji,jj) ) / e1u(ji,jj)
572                  zv_spg = zv_spg + grav * ( pot_astro(ji,jj+1) - pot_astro(ji,jj) ) / e2v(ji,jj)
573                  ENDIF
574                  ! energy/enstrophy conserving formulation for planetary vorticity term
575                  zu_cor = + z1_4 / e1u(ji,jj) * (  ftne(ji,jj  ) * zwy(ji  ,jj  ) + ftnw(ji+1,jj) * zwy(ji+1,jj  )   &
576                     &                           + ftse(ji,jj  ) * zwy(ji  ,jj-1) + ftsw(ji+1,jj) * zwy(ji+1,jj-1) ) * hur_e(ji,jj)
577                  zv_cor = - z1_4 / e2v(ji,jj) * (  ftsw(ji,jj+1) * zwx(ji-1,jj+1) + ftse(ji,jj+1) * zwx(ji  ,jj+1)   &
578                     &                           + ftnw(ji,jj  ) * zwx(ji-1,jj  ) + ftne(ji,jj  ) * zwx(ji  ,jj  ) ) * hvr_e(ji,jj)
579                  ! after velocities with implicit bottom friction
580                  IF( ln_bfrimp ) THEN
581                     !   A new method to implement the implicit bottom friction.
582                     !   H. Liu
583                     !   Sept 2011
584                     ua_e(ji,jj) = umask(ji,jj,1) * ( zub_e(ji,jj) +                                            &
585                      &                               z2dt_e * ( zu_cor + zu_spg + zu_sld + zu_asp )            &
586                      &                               / ( 1._wp      - z2dt_e * bfrua(ji,jj) * hur_e(ji,jj) ) )
587                     ua_e(ji,jj) = ( ua_e(ji,jj) + z2dt_e *   zua(ji,jj)  ) * umask(ji,jj,1)   
588                     !
589                     va_e(ji,jj) = vmask(ji,jj,1) * ( zvb_e(ji,jj) +                                            &
590                      &                               z2dt_e * ( zv_cor + zv_spg + zv_sld + zv_asp )            &
591                      &                               / ( 1._wp      - z2dt_e * bfrva(ji,jj) * hvr_e(ji,jj) ) )
592                     va_e(ji,jj) = ( va_e(ji,jj) + z2dt_e *   zva(ji,jj)  ) * vmask(ji,jj,1)   
593                     !
594                  ELSE
595                     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)   &
596                     &            / ( 1._wp        - z2dt_e * bfrua(ji,jj) * hur_e(ji,jj) )
597                     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)   &
598                     &            / ( 1._wp        - z2dt_e * bfrva(ji,jj) * hvr_e(ji,jj) )
599                  ENDIF
600               END DO
601            END DO
602            !
603         ENDIF
604         !                                                     !* domain lateral boundary
605         !                                                     !  -----------------------
606
607                                                               ! OBC open boundaries
608         IF( lk_obc               )   CALL obc_fla_ts ( ua_e, va_e, sshn_e, ssha_e )
609
610                                                               ! BDY open boundaries
611#if defined key_bdy
612         pssh => sshn_e
613         phur => hur_e
614         phvr => hvr_e
615         pu2d => ua_e
616         pv2d => va_e
617
618         IF( lk_bdy )   CALL bdy_dyn2d( kt ) 
619#endif
620
621         !
622         CALL lbc_lnk( ua_e  , 'U', -1. )                      ! local domain boundaries
623         CALL lbc_lnk( va_e  , 'V', -1. )
624         CALL lbc_lnk( ssha_e, 'T',  1. )
625
626         zu_sum  (:,:) = zu_sum  (:,:) + ua_e  (:,:)           ! Sum over sub-time-steps
627         zv_sum  (:,:) = zv_sum  (:,:) + va_e  (:,:) 
628         zssh_sum(:,:) = zssh_sum(:,:) + ssha_e(:,:) 
629
630         !                                                !* Time filter and swap
631         !                                                !  --------------------
632         IF( jn == 1 ) THEN                                     ! Swap only (1st Euler time step)
633            zsshb_e(:,:) = sshn_e(:,:)
634            zub_e  (:,:) = zun_e (:,:)
635            zvb_e  (:,:) = zvn_e (:,:)
636            sshn_e (:,:) = ssha_e(:,:)
637            zun_e  (:,:) = ua_e  (:,:)
638            zvn_e  (:,:) = va_e  (:,:)
639         ELSE                                                   ! Swap + Filter
640            zsshb_e(:,:) = atfp * ( zsshb_e(:,:) + ssha_e(:,:) ) + atfp1 * sshn_e(:,:)
641            zub_e  (:,:) = atfp * ( zub_e  (:,:) + ua_e  (:,:) ) + atfp1 * zun_e (:,:)
642            zvb_e  (:,:) = atfp * ( zvb_e  (:,:) + va_e  (:,:) ) + atfp1 * zvn_e (:,:)
643            sshn_e (:,:) = ssha_e(:,:)
644            zun_e  (:,:) = ua_e  (:,:)
645            zvn_e  (:,:) = va_e  (:,:)
646         ENDIF
647
648         IF( lk_vvl ) THEN                                !* Update ocean depth (variable volume case only)
649            !                                             !  ------------------
650            DO jj = 1, jpjm1                                    ! Sea Surface Height at u- & v-points
651               DO ji = 1, fs_jpim1   ! Vector opt.
652                  zsshun_e(ji,jj) = 0.5_wp * umask(ji,jj,1) / ( e1u(ji,jj) * e2u(ji,jj) )       &
653                     &                     * ( e1t(ji  ,jj) * e2t(ji  ,jj) * sshn_e(ji  ,jj)    &
654                     &                     +   e1t(ji+1,jj) * e2t(ji+1,jj) * sshn_e(ji+1,jj) )
655                  zsshvn_e(ji,jj) = 0.5_wp * vmask(ji,jj,1) / ( e1v(ji,jj) * e2v(ji,jj) )       &
656                     &                     * ( e1t(ji,jj  ) * e2t(ji,jj  ) * sshn_e(ji,jj  )    &
657                     &                     +   e1t(ji,jj+1) * e2t(ji,jj+1) * sshn_e(ji,jj+1) )
658               END DO
659            END DO
660            CALL lbc_lnk( zsshun_e, 'U', 1. )                   ! lateral boundaries conditions
661            CALL lbc_lnk( zsshvn_e, 'V', 1. ) 
662            !
663            hu_e (:,:) = hu_0(:,:) + zsshun_e(:,:)              ! Ocean depth at U- and V-points
664            hv_e (:,:) = hv_0(:,:) + zsshvn_e(:,:)
665            hur_e(:,:) = umask(:,:,1) / ( hu_e(:,:) + 1._wp - umask(:,:,1) )
666            hvr_e(:,:) = vmask(:,:,1) / ( hv_e(:,:) + 1._wp - vmask(:,:,1) )
667            !
668         ENDIF
669         !                                                 ! ==================== !
670      END DO                                               !        end loop      !
671      !                                                    ! ==================== !
672
673#if defined key_obc
674      IF( lp_obc_east  )   sshfoe_b(:,:) = zcoef * sshfoe_b(:,:)     !!gm totally useless ?????
675      IF( lp_obc_west  )   sshfow_b(:,:) = zcoef * sshfow_b(:,:)
676      IF( lp_obc_north )   sshfon_b(:,:) = zcoef * sshfon_b(:,:)
677      IF( lp_obc_south )   sshfos_b(:,:) = zcoef * sshfos_b(:,:)
678#endif
679
680      ! -----------------------------------------------------------------------------
681      ! Phase 3. update the general trend with the barotropic trend
682      ! -----------------------------------------------------------------------------
683      !
684      !                                   !* Time average ==> after barotropic u, v, ssh
685      zcoef =  1._wp / ( 2 * nn_baro  + 1 ) 
686      zu_sum(:,:) = zcoef * zu_sum  (:,:) 
687      zv_sum(:,:) = zcoef * zv_sum  (:,:) 
688      !
689      !                                   !* update the general momentum trend
690      DO jk=1,jpkm1
691         ua(:,:,jk) = ua(:,:,jk) + ( zu_sum(:,:) - ub_b(:,:) ) * z1_2dt_b
692         va(:,:,jk) = va(:,:,jk) + ( zv_sum(:,:) - vb_b(:,:) ) * z1_2dt_b
693      END DO
694      un_b  (:,:) =  zu_sum(:,:) 
695      vn_b  (:,:) =  zv_sum(:,:) 
696      sshn_b(:,:) = zcoef * zssh_sum(:,:) 
697      !
698      !                                   !* write time-spliting arrays in the restart
699      IF( lrst_oce )   CALL ts_rst( kt, 'WRITE' )
700      !
701      CALL wrk_dealloc( jpi, jpj, zsshun_e, zsshvn_e, zsshb_e, zssh_sum, zhdiv     )
702      CALL wrk_dealloc( jpi, jpj, zua, zva, zun, zvn, zun_e, zvn_e, zub_e, zvb_e   )
703      CALL wrk_dealloc( jpi, jpj, zcu, zcv, zwx, zwy, zbfru, zbfrv, zu_sum, zv_sum )
704      CALL wrk_dealloc( jpi, jpj, zhu_b, zhv_b                                     )
705      !
706      IF( nn_timing == 1 )  CALL timing_stop('dyn_spg_ts')
707      !
708   END SUBROUTINE dyn_spg_ts
709
710
711   SUBROUTINE ts_rst( kt, cdrw )
712      !!---------------------------------------------------------------------
713      !!                   ***  ROUTINE ts_rst  ***
714      !!
715      !! ** Purpose : Read or write time-splitting arrays in restart file
716      !!----------------------------------------------------------------------
717      INTEGER         , INTENT(in) ::   kt         ! ocean time-step
718      CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
719      !
720      REAL(wp), POINTER, DIMENSION(:,:) :: zzhu_b, zzhv_b
721      INTEGER ::  ji, jk        ! dummy loop indices
722      !!----------------------------------------------------------------------
723      !
724      IF( TRIM(cdrw) == 'READ' ) THEN
725         IF( iom_varid( numror, 'un_b', ldstop = .FALSE. ) > 0 ) THEN
726            CALL iom_get( numror, jpdom_autoglo, 'un_b'  , un_b  (:,:) )   ! external velocity issued
727            CALL iom_get( numror, jpdom_autoglo, 'vn_b'  , vn_b  (:,:) )   ! from barotropic loop
728         ELSE
729            CALL wrk_alloc( jpi, jpj, zzhu_b, zzhv_b )
730            un_b (:,:) = 0._wp
731            vn_b (:,:) = 0._wp
732            ! vertical sum
733            IF( lk_vopt_loop ) THEN          ! vector opt., forced unroll
734               DO jk = 1, jpkm1
735                  DO ji = 1, jpij
736                     un_b(ji,1) = un_b(ji,1) + fse3u_n(ji,1,jk) * un(ji,1,jk)
737                     vn_b(ji,1) = vn_b(ji,1) + fse3v_n(ji,1,jk) * vn(ji,1,jk)
738                  END DO
739               END DO
740            ELSE                             ! No  vector opt.
741               DO jk = 1, jpkm1
742                  un_b(:,:) = un_b(:,:) + fse3u_n(:,:,jk) * un(:,:,jk)
743                  vn_b(:,:) = vn_b(:,:) + fse3v_n(:,:,jk) * vn(:,:,jk)
744               END DO
745            ENDIF
746            un_b (:,:) = un_b(:,:) * hur(:,:)
747            vn_b (:,:) = vn_b(:,:) * hvr(:,:)
748         ENDIF
749
750         ! Vertically integrated velocity (before)
751         IF (neuler/=0) THEN
752            ub_b (:,:) = 0._wp
753            vb_b (:,:) = 0._wp
754
755            ! vertical sum
756            IF( lk_vopt_loop ) THEN          ! vector opt., forced unroll
757               DO jk = 1, jpkm1
758                  DO ji = 1, jpij
759                     ub_b(ji,1) = ub_b(ji,1) + fse3u_b(ji,1,jk) * ub(ji,1,jk)
760                     vb_b(ji,1) = vb_b(ji,1) + fse3v_b(ji,1,jk) * vb(ji,1,jk)
761                  END DO
762               END DO
763            ELSE                             ! No  vector opt.
764               DO jk = 1, jpkm1
765                  ub_b(:,:) = ub_b(:,:) + fse3u_b(:,:,jk) * ub(:,:,jk)
766                  vb_b(:,:) = vb_b(:,:) + fse3v_b(:,:,jk) * vb(:,:,jk)
767               END DO
768            ENDIF
769
770            IF( lk_vvl ) THEN
771               CALL wrk_alloc( jpi, jpj, zzhu_b, zzhv_b )
772               ub_b (:,:) = 0.
773               vb_b (:,:) = 0.
774               zzhu_b(:,:) = 0.
775               zzhv_b(:,:) = 0.
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                        zzhu_b(ji,1) = zhu_b(ji,1) + fse3u_b(ji,1,jk) * umask(ji,1,jk)
783                        zzhv_b(ji,1) = zhv_b(ji,1) + fse3v_b(ji,1,jk) * vmask(ji,1,jk)
784                     END DO
785                  END DO
786               ELSE                             ! No  vector opt.
787                  DO jk = 1, jpkm1
788                     ub_b  (:,:) = ub_b  (:,:) + fse3u_b(:,:,jk) * ub   (:,:,jk)
789                     vb_b  (:,:) = vb_b  (:,:) + fse3v_b(:,:,jk) * vb   (:,:,jk)
790                     zzhu_b(:,:) = zzhu_b(:,:) + fse3u_b(:,:,jk) * umask(:,:,jk)
791                     zzhv_b(:,:) = zzhv_b(:,:) + fse3v_b(:,:,jk) * vmask(:,:,jk)
792                  END DO
793               ENDIF
794               ub_b(:,:) = ub_b(:,:) / ( zzhu_b(:,:) + 1. - umask(:,:,1) )
795               vb_b(:,:) = vb_b(:,:) / ( zzhv_b(:,:) + 1. - vmask(:,:,1) )
796               CALL wrk_dealloc( jpi, jpj, zzhu_b, zzhv_b )
797            ELSE             
798               ub_b (:,:) = 0.e0
799               vb_b (:,:) = 0.e0
800               ! vertical sum
801               IF( lk_vopt_loop ) THEN          ! vector opt., forced unroll
802                  DO jk = 1, jpkm1
803                     DO ji = 1, jpij
804                        ub_b(ji,1) = ub_b(ji,1) + fse3u_b(ji,1,jk) * ub(ji,1,jk)
805                        vb_b(ji,1) = vb_b(ji,1) + fse3v_b(ji,1,jk) * vb(ji,1,jk)
806                     END DO
807                  END DO
808               ELSE                             ! No  vector opt.
809                  DO jk = 1, jpkm1
810                     ub_b(:,:) = ub_b(:,:) + fse3u_b(:,:,jk) * ub(:,:,jk)
811                     vb_b(:,:) = vb_b(:,:) + fse3v_b(:,:,jk) * vb(:,:,jk)
812                  END DO
813               ENDIF
814               ub_b(:,:) = ub_b(:,:) * hur(:,:)
815               vb_b(:,:) = vb_b(:,:) * hvr(:,:)
816            ENDIF
817         ELSE                                 ! neuler==0
818            ub_b (:,:) = un_b (:,:)
819            vb_b (:,:) = vn_b (:,:)
820         ENDIF
821
822         IF( iom_varid( numror, 'sshn_b', ldstop = .FALSE. ) > 0 ) THEN
823            CALL iom_get( numror, jpdom_autoglo, 'sshn_b' , sshn_b (:,:) )   ! filtered ssh
824         ELSE
825            sshn_b(:,:) = sshb(:,:)   ! if not in restart set previous time mean to current baroclinic before value   
826         ENDIF
827      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN
828         CALL iom_rstput( kt, nitrst, numrow, 'un_b'   , un_b  (:,:) )   ! external velocity and ssh
829         CALL iom_rstput( kt, nitrst, numrow, 'vn_b'   , vn_b  (:,:) )   ! issued from barotropic loop
830         CALL iom_rstput( kt, nitrst, numrow, 'sshn_b' , sshn_b(:,:) )   !
831      ENDIF
832      !
833   END SUBROUTINE ts_rst
834
835#else
836   !!----------------------------------------------------------------------
837   !!   Default case :   Empty module   No standart free surface cst volume
838   !!----------------------------------------------------------------------
839CONTAINS
840   INTEGER FUNCTION dyn_spg_ts_alloc()    ! Dummy function
841      dyn_spg_ts_alloc = 0
842   END FUNCTION dyn_spg_ts_alloc
843   SUBROUTINE dyn_spg_ts( kt )            ! Empty routine
844      INTEGER, INTENT(in) :: kt
845      WRITE(*,*) 'dyn_spg_ts: You should not have seen this print! error?', kt
846   END SUBROUTINE dyn_spg_ts
847   SUBROUTINE ts_rst( kt, cdrw )          ! Empty routine
848      INTEGER         , INTENT(in) ::   kt         ! ocean time-step
849      CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
850      WRITE(*,*) 'ts_rst    : You should not have seen this print! error?', kt, cdrw
851   END SUBROUTINE ts_rst   
852#endif
853   
854   !!======================================================================
855END MODULE dynspg_ts
Note: See TracBrowser for help on using the repository browser.