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 trunk/NEMO/OPA_SRC/DYN – NEMO

source: trunk/NEMO/OPA_SRC/DYN/dynspg_ts.F90 @ 1146

Last change on this file since 1146 was 1146, checked in by rblod, 16 years ago

Add svn Id (first try), see ticket #210

  • Property svn:eol-style set to native
  • Property svn:executable set to *
  • Property svn:keywords set to Id
File size: 34.5 KB
Line 
1MODULE dynspg_ts
2   !!======================================================================
3   !!                   ***  MODULE  dynspg_ts  ***
4   !! Ocean dynamics:  surface pressure gradient trend
5   !!======================================================================
6   !! History :   9.0  !  04-12  (L. Bessieres, G. Madec)  Original code
7   !!             " "  !  05-11  (V. Garnier, G. Madec)  optimization
8   !!             9.0  !  06-08  (S. Masson)  distributed restart using iom
9   !!              "   !  08-01  (R. Benshila)  change averaging method
10   !!              "   !  07-07  (D. Storkey) calls to BDY routines
11   !!---------------------------------------------------------------------
12#if defined key_dynspg_ts   ||   defined key_esopa
13   !!----------------------------------------------------------------------
14   !!   'key_dynspg_ts'         free surface cst volume with time splitting
15   !!----------------------------------------------------------------------
16   !!----------------------------------------------------------------------
17   !!   dyn_spg_ts  : compute surface pressure gradient trend using a time-
18   !!                 splitting scheme and add to the general trend
19   !!   ts_rst      : read/write the time-splitting restart fields in the ocean restart file
20   !!----------------------------------------------------------------------
21   !! * Modules used
22   USE oce             ! ocean dynamics and tracers
23   USE dom_oce         ! ocean space and time domain
24   USE sbc_oce         ! surface boundary condition: ocean
25   USE dynspg_oce      ! surface pressure gradient variables
26   USE phycst          ! physical constants
27   USE domvvl          ! variable volume
28   USE obcdta          ! open boundary condition data     
29   USE obcfla          ! Flather open boundary condition 
30   USE dynvor          ! vorticity term
31   USE obc_oce         ! Lateral open boundary condition
32   USE obc_par         ! open boundary condition parameters
33   USE bdy_oce         ! unstructured open boundaries
34   USE bdy_par         ! unstructured open boundaries
35   USE bdydta          ! unstructured open boundaries
36   USE bdydyn          ! unstructured open boundaries
37   USE bdytides        ! tidal forcing at unstructured open boundaries.
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
43   USE restart         ! only for lrst_oce
44
45   IMPLICIT NONE
46   PRIVATE
47
48   PUBLIC dyn_spg_ts  ! routine called by step.F90
49   PUBLIC ts_rst      ! routine called by istate.F90
50
51   REAL(wp), DIMENSION(jpi,jpj) ::  ftnw, ftne,   &  ! triad of coriolis parameter
52      &                             ftsw, ftse       ! (only used with een vorticity scheme)
53
54
55   !! * Substitutions
56#  include "domzgr_substitute.h90"
57#  include "vectopt_loop_substitute.h90"
58   !!----------------------------------------------------------------------
59   !!  OPA 9.0 , LOCEAN-IPSL (2005)
60   !! $Id$
61   !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt
62   !!----------------------------------------------------------------------
63
64CONTAINS
65
66   SUBROUTINE dyn_spg_ts( kt )
67      !!----------------------------------------------------------------------
68      !!                  ***  routine dyn_spg_ts  ***
69      !!
70      !! ** Purpose :   Compute the now trend due to the surface pressure
71      !!      gradient in case of free surface formulation with time-splitting.
72      !!      Add it to the general trend of momentum equation.
73      !!      Compute the free surface.
74      !!
75      !! ** Method  :   Free surface formulation with time-splitting
76      !!      -1- Save the vertically integrated trend. This general trend is
77      !!          held constant over the barotropic integration.
78      !!          The Coriolis force is removed from the general trend as the
79      !!          surface gradient and the Coriolis force are updated within
80      !!          the barotropic integration.
81      !!      -2- Barotropic loop : updates of sea surface height (ssha_e) and
82      !!          barotropic transports (ua_e and va_e) through barotropic
83      !!          momentum and continuity integration. Barotropic former
84      !!          variables are time averaging over the full barotropic cycle
85      !!          (= 2 * baroclinic time step) and saved in zsshX_b, zuX_b
86      !!          and zvX_b (X specifying after, now or before).
87      !!      -3- Update of sea surface height from time averaged barotropic
88      !!          variables.
89      !!        - apply lateral boundary conditions on sshn.
90      !!      -4- The new general trend becomes :
91      !!          ua = ua - sum_k(ua)/H + ( zua_b - sum_k(ub) )/H
92      !!
93      !! ** Action : - Update (ua,va) with the surf. pressure gradient trend
94      !!
95      !! References : Griffies et al., (2003): A technical guide to MOM4. NOAA/GFDL
96      !!---------------------------------------------------------------------
97      INTEGER, INTENT( in )  ::   kt           ! ocean time-step index
98
99      !! * Local declarations
100      INTEGER  ::  ji, jj, jk, jit             ! dummy loop indices
101      INTEGER  ::  icycle, ibaro               ! temporary scalar
102      REAL(wp) ::                           &
103         zraur, zcoef, z2dt_e, z2dt_b, zfac25,   &  ! temporary scalars
104         zfact1, zspgu, zcubt, zx1, zy1,    &  !     "        "
105         zfact2, zspgv, zcvbt, zx2, zy2        !     "        "
106      REAL(wp), DIMENSION(jpi,jpj) ::       &
107         zcu, zcv, zwx, zwy, zhdiv,         &  ! temporary arrays
108         zua, zva, zub, zvb,                &  !     "        "
109         zssha_b, zua_b, zva_b,             &  !     "        "
110         zub_e, zvb_e,                      &  !     "        "
111         zun_e, zvn_e                          !     "        "
112      !! Variable volume
113      REAL(wp), DIMENSION(jpi,jpj)     ::   &
114         zspgu_1, zspgv_1, zsshun_e, zsshvn_e                     ! 2D workspace
115      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zfse3un_e, zfse3vn_e  ! 3D workspace
116      !!----------------------------------------------------------------------
117
118      ! Arrays initialization
119      ! ---------------------
120      zua_b(:,:) = 0.e0   ;   zub_e(:,:) = 0.e0   ;   zun_e(:,:) = 0.e0
121      zva_b(:,:) = 0.e0   ;   zvb_e(:,:) = 0.e0   ;   zvn_e(:,:) = 0.e0
122      zhdiv(:,:) = 0.e0
123
124
125      IF( kt == nit000 ) THEN
126         !
127         IF(lwp) WRITE(numout,*)
128         IF(lwp) WRITE(numout,*) 'dyn_spg_ts : surface pressure gradient trend'
129         IF(lwp) WRITE(numout,*) '~~~~~~~~~~   free surface with time splitting'
130         IF(lwp) WRITE(numout,*) ' Number of sub cycle in 1 time-step (2 rdt) : icycle = ', FLOOR( 2*rdt/rdtbt )
131
132         CALL ts_rst( nit000, 'READ' )   ! read or initialize the following fields:
133         !                               ! sshb, sshn, sshb_b, sshn_b, un_b, vn_b
134
135         ssha_e(:,:) = sshn(:,:)
136         ua_e(:,:)   = un_b(:,:)
137         va_e(:,:)   = vn_b(:,:)
138         hu_e(:,:)   = hu(:,:)
139         hv_e(:,:)   = hv(:,:)
140
141         IF( ln_dynvor_een ) THEN
142            ftne(1,:) = 0.e0   ;   ftnw(1,:) = 0.e0   ;   ftse(1,:) = 0.e0   ;   ftsw(1,:) = 0.e0
143            DO jj = 2, jpj
144               DO ji = fs_2, jpi   ! vector opt.
145                  ftne(ji,jj) = ( ff(ji-1,jj  ) + ff(ji  ,jj  ) + ff(ji  ,jj-1) ) / 3.
146                  ftnw(ji,jj) = ( ff(ji-1,jj-1) + ff(ji-1,jj  ) + ff(ji  ,jj  ) ) / 3.
147                  ftse(ji,jj) = ( ff(ji  ,jj  ) + ff(ji  ,jj-1) + ff(ji-1,jj-1) ) / 3.
148                  ftsw(ji,jj) = ( ff(ji  ,jj-1) + ff(ji-1,jj-1) + ff(ji-1,jj  ) ) / 3.
149               END DO
150            END DO
151         ENDIF
152         !
153      ENDIF
154
155      ! Substract the surface pressure gradient from the momentum
156      ! ---------------------------------------------------------
157      IF( lk_vvl ) THEN    ! Variable volume
158
159         ! 0. Local initialization
160         ! -----------------------
161         zspgu_1(:,:) = 0.e0
162         zspgv_1(:,:) = 0.e0
163
164         ! 1. Surface pressure gradient (now)
165         ! ----------------------------
166         DO jj = 2, jpjm1
167            DO ji = fs_2, fs_jpim1   ! vector opt.
168               zspgu_1(ji,jj) = -grav * ( ( rhd(ji+1,jj  ,1) + 1 ) * sshn(ji+1,jj  ) &
169                  &                     - ( rhd(ji  ,jj  ,1) + 1 ) * sshn(ji  ,jj  ) ) / e1u(ji,jj)
170               zspgv_1(ji,jj) = -grav * ( ( rhd(ji  ,jj+1,1) + 1 ) * sshn(ji  ,jj+1) &
171                  &                     - ( rhd(ji  ,jj  ,1) + 1 ) * sshn(ji  ,jj  ) ) / e2v(ji,jj)
172            END DO
173         END DO
174
175         ! 2. Substract the surface pressure trend from the general trend
176         ! ------------------------------------------------------
177         DO jk = 1, jpkm1
178            DO jj = 2, jpjm1
179               DO ji = fs_2, fs_jpim1   ! vector opt.
180                  ua(ji,jj,jk) = ua(ji,jj,jk) - zspgu_1(ji,jj)
181                  va(ji,jj,jk) = va(ji,jj,jk) - zspgv_1(ji,jj)
182               END DO
183            END DO
184         END DO
185
186      ENDIF
187
188      ! Local constant initialization
189      ! --------------------------------
190      z2dt_b = 2.0 * rdt                                    ! baroclinic time step
191      IF ( neuler == 0 .AND. kt == nit000 ) z2dt_b = rdt
192      zfact1 = 0.5 * 0.25                                   ! coefficient for vorticity estimates
193      zfact2 = 0.5 * 0.5
194      zraur  = 1. / rauw                                    ! 1 / volumic mass of pure water
195     
196      ! -----------------------------------------------------------------------------
197      !  Phase 1 : Coupling between general trend and barotropic estimates (1st step)
198      ! -----------------------------------------------------------------------------
199
200      ! Vertically integrated quantities
201      ! --------------------------------
202      zua(:,:) = 0.e0
203      zva(:,:) = 0.e0
204      zub(:,:) = 0.e0
205      zvb(:,:) = 0.e0
206      zwx(:,:) = 0.e0
207      zwy(:,:) = 0.e0
208
209      ! vertical sum
210      IF( lk_vopt_loop ) THEN          ! vector opt., forced unroll
211         DO jk = 1, jpkm1
212            DO ji = 1, jpij
213               !                                                           ! Vertically integrated momentum trends
214               zua(ji,1) = zua(ji,1) + fse3u(ji,1,jk) * umask(ji,1,jk) * ua(ji,1,jk)
215               zva(ji,1) = zva(ji,1) + fse3v(ji,1,jk) * vmask(ji,1,jk) * va(ji,1,jk)
216               !                                                           ! Vertically integrated transports (before)
217               zub(ji,1) = zub(ji,1) + fse3u(ji,1,jk) * ub(ji,1,jk)
218               zvb(ji,1) = zvb(ji,1) + fse3v(ji,1,jk) * vb(ji,1,jk)
219               !                                                           ! Planetary vorticity transport fluxes (now)
220               zwx(ji,1) = zwx(ji,1) + e2u(ji,1) * fse3u(ji,1,jk) * un(ji,1,jk)
221               zwy(ji,1) = zwy(ji,1) + e1v(ji,1) * fse3v(ji,1,jk) * vn(ji,1,jk)
222            END DO
223         END DO
224      ELSE                             ! No  vector opt.
225         DO jk = 1, jpkm1
226            !                                                           ! Vertically integrated momentum trends
227            zua(:,:) = zua(:,:) + fse3u(:,:,jk) * umask(:,:,jk) * ua(:,:,jk)
228            zva(:,:) = zva(:,:) + fse3v(:,:,jk) * vmask(:,:,jk) * va(:,:,jk)
229            !                                                           ! Vertically integrated transports (before)
230            zub(:,:) = zub(:,:) + fse3u(:,:,jk) * ub(:,:,jk)
231            zvb(:,:) = zvb(:,:) + fse3v(:,:,jk) * vb(:,:,jk)
232            !                                                           ! Planetary vorticity (now)
233            zwx(:,:) = zwx(:,:) + e2u(:,:) * fse3u(:,:,jk) * un(:,:,jk)
234            zwy(:,:) = zwy(:,:) + e1v(:,:) * fse3v(:,:,jk) * vn(:,:,jk)
235         END DO
236      ENDIF
237
238      IF( ln_dynvor_ene .OR. ln_dynvor_mix ) THEN      ! energy conserving or mixed scheme
239         DO jj = 2, jpjm1
240            DO ji = fs_2, fs_jpim1   ! vector opt.
241               zy1 = ( zwy(ji,jj-1) + zwy(ji+1,jj-1) ) / e1u(ji,jj)
242               zy2 = ( zwy(ji,jj  ) + zwy(ji+1,jj  ) ) / e1u(ji,jj)
243               zx1 = ( zwx(ji-1,jj) + zwx(ji-1,jj+1) ) / e2v(ji,jj)
244               zx2 = ( zwx(ji  ,jj) + zwx(ji  ,jj+1) ) / e2v(ji,jj)
245               ! energy conserving formulation for planetary vorticity term
246               zcu(ji,jj) = zfact2 * ( ff(ji  ,jj-1) * zy1 + ff(ji,jj) * zy2 )
247               zcv(ji,jj) =-zfact2 * ( ff(ji-1,jj  ) * zx1 + ff(ji,jj) * zx2 )
248            END DO
249         END DO
250         !
251      ELSEIF ( ln_dynvor_ens ) THEN                    ! enstrophy conserving scheme
252         DO jj = 2, jpjm1
253            DO ji = fs_2, fs_jpim1   ! vector opt.
254               zy1 = zfact1 * ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1)   &
255                              + zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) / e1u(ji,jj)
256               zx1 =-zfact1 * ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1)   &
257                              + zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) / e2v(ji,jj)
258               zcu(ji,jj)  = zy1 * ( ff(ji  ,jj-1) + ff(ji,jj) )
259               zcv(ji,jj)  = zx1 * ( ff(ji-1,jj  ) + ff(ji,jj) )
260            END DO
261         END DO
262         !
263      ELSEIF ( ln_dynvor_een ) THEN                    ! enstrophy and energy conserving scheme
264         zfac25 = 0.25
265         DO jj = 2, jpjm1
266            DO ji = fs_2, fs_jpim1   ! vector opt.
267               zcu(ji,jj) = + zfac25 / e1u(ji,jj)   &
268                  &       * (  ftne(ji,jj  ) * zwy(ji  ,jj  ) + ftnw(ji+1,jj) * zwy(ji+1,jj  )   &
269                  &          + ftse(ji,jj  ) * zwy(ji  ,jj-1) + ftsw(ji+1,jj) * zwy(ji+1,jj-1) )
270               zcv(ji,jj) = - zfac25 / e2v(ji,jj)   &
271                  &       * (  ftsw(ji,jj+1) * zwx(ji-1,jj+1) + ftse(ji,jj+1) * zwx(ji  ,jj+1)   &
272                  &          + ftnw(ji,jj  ) * zwx(ji-1,jj  ) + ftne(ji,jj  ) * zwx(ji  ,jj  ) )
273            END DO
274         END DO
275         !
276      ENDIF
277
278
279      ! Remove barotropic trend from general momentum trend
280      ! ---------------------------------------------------
281      DO jk = 1 , jpkm1
282         DO jj = 2, jpjm1
283            DO ji = fs_2, fs_jpim1   ! vector opt.
284               ua(ji,jj,jk) = ua(ji,jj,jk) - zua(ji,jj) * hur(ji,jj)
285               va(ji,jj,jk) = va(ji,jj,jk) - zva(ji,jj) * hvr(ji,jj)
286            END DO
287         END DO
288      END DO
289
290      ! Remove coriolis term from barotropic trend
291      ! ------------------------------------------
292      DO jj = 2, jpjm1
293         DO ji = fs_2, fs_jpim1
294            zua(ji,jj) = zua(ji,jj) - zcu(ji,jj)
295            zva(ji,jj) = zva(ji,jj) - zcv(ji,jj)
296         END DO
297      END DO
298
299      ! -----------------------------------------------------------------------
300      !  Phase 2 : Integration of the barotropic equations with time splitting
301      ! -----------------------------------------------------------------------
302
303      ! Initialisations
304      !----------------
305      ! Number of iteration of the barotropic loop
306      ibaro = FLOOR( rdt / rdtbt )
307      icycle = 3 /2 * ibaro 
308
309      ! variables for the barotropic equations
310      sshb_e (:,:) = sshn_b(:,:)       ! (barotropic) sea surface height (before and now)
311      sshn_e (:,:) = sshn_b(:,:)
312      zub_e  (:,:) = un_b  (:,:)       ! barotropic transports issued from the barotropic equations (before and now)
313      zvb_e  (:,:) = vn_b  (:,:)
314      zun_e  (:,:) = un_b  (:,:)
315      zvn_e  (:,:) = vn_b  (:,:)
316      zssha_b(:,:) = 0.e0
317      zua_b  (:,:) = 0.e0
318      zva_b  (:,:) = 0.e0
319      hu_e   (:,:) = hu   (:,:)     ! (barotropic) ocean depth at u-point
320      hv_e   (:,:) = hv   (:,:)     ! (barotropic) ocean depth at v-point
321      IF( lk_vvl ) THEN
322         zsshun_e(:,:)    = sshu (:,:)     ! (barotropic) sea surface height (now) at u-point
323         zsshvn_e(:,:)    = sshv (:,:)     ! (barotropic) sea surface height (now) at v-point
324         zfse3un_e(:,:,:) = fse3u(:,:,:)   ! (barotropic) scale factors  at u-point
325         zfse3un_e(:,:,:) = fse3v(:,:,:)   ! (barotropic) scale factors  at v-point
326      ENDIF
327
328      ! set ssh corrections to 0
329      ! ssh corrections are applied to normal velocities (Flather's algorithm) and averaged over the barotropic loop
330#if defined key_obc
331      IF( lp_obc_east  )   sshfoe_b(:,:) = 0.e0
332      IF( lp_obc_west  )   sshfow_b(:,:) = 0.e0
333      IF( lp_obc_south )   sshfos_b(:,:) = 0.e0
334      IF( lp_obc_north )   sshfon_b(:,:) = 0.e0
335#endif
336
337      ! Barotropic integration over 2 baroclinic time steps
338      ! ---------------------------------------------------
339
340      !                                                    ! ==================== !
341      DO jit = 1, icycle                                   !  sub-time-step loop  !
342         !                                                 ! ==================== !
343         z2dt_e = 2. * rdtbt
344         IF ( jit == 1 )   z2dt_e = rdtbt
345
346         ! Time interpolation of open boundary condition data
347         IF( lk_obc )   CALL obc_dta_bt( kt, jit )
348         IF( lk_bdy .or. lk_bdy_tides)   CALL bdy_dta_bt( kt, jit+1 )
349
350         ! Horizontal divergence of barotropic transports
351         !--------------------------------------------------
352         zhdiv(:,:) = 0.e0
353         DO jj = 2, jpjm1
354            DO ji = fs_2, fs_jpim1   ! vector opt.
355               zhdiv(ji,jj) = ( e2u(ji  ,jj  ) * zun_e(ji  ,jj)              &
356                  &            -e2u(ji-1,jj  ) * zun_e(ji-1,jj)              &
357                  &            +e1v(ji  ,jj  ) * zvn_e(ji  ,jj)              &
358                  &            -e1v(ji  ,jj-1) * zvn_e(ji  ,jj-1) )          &
359                  &           / (e1t(ji,jj)*e2t(ji,jj))
360            END DO
361         END DO
362
363#if defined key_obc
364         ! open boundaries (div must be zero behind the open boundary)
365         !  mpp remark: The zeroing of hdiv can probably be extended to 1->jpi/jpj for the correct row/column
366         IF( lp_obc_east  )   zhdiv(nie0p1:nie1p1,nje0  :nje1)   = 0.e0      ! east
367         IF( lp_obc_west  )   zhdiv(niw0  :niw1  ,njw0  :njw1)   = 0.e0      ! west
368         IF( lp_obc_north )   zhdiv(nin0  :nin1  ,njn0p1:njn1p1) = 0.e0      ! north
369         IF( lp_obc_south )   zhdiv(nis0  :nis1  ,njs0  :njs1)   = 0.e0      ! south
370#endif
371
372         IF( lk_bdy .or. lk_bdy_tides ) THEN
373            DO jj = 1, jpj
374               DO ji = 1, jpi
375                  zhdiv(ji,jj) = zhdiv(ji,jj)*bdytmask(ji,jj)
376               END DO
377            END DO
378         ENDIF
379
380         ! Sea surface height from the barotropic system
381         !----------------------------------------------
382         DO jj = 2, jpjm1
383            DO ji = fs_2, fs_jpim1   ! vector opt.
384               ssha_e(ji,jj) = ( sshb_e(ji,jj) - z2dt_e *  ( zraur * emp(ji,jj)  &
385            &            +  zhdiv(ji,jj) ) ) * tmask(ji,jj,1)
386            END DO
387         END DO
388
389         ! evolution of the barotropic transport ( following the vorticity scheme used)
390         ! ----------------------------------------------------------------------------
391         zwx(:,:) = e2u(:,:) * zun_e(:,:)
392         zwy(:,:) = e1v(:,:) * zvn_e(:,:)
393
394         IF( ln_dynvor_ene .OR. ln_dynvor_mix ) THEN      ! energy conserving or mixed scheme
395            DO jj = 2, jpjm1
396               DO ji = fs_2, fs_jpim1   ! vector opt.
397                  ! surface pressure gradient
398                  IF( lk_vvl) THEN
399                     zspgu = -grav * ( ( rhd(ji+1,jj,1) + 1 ) * sshn_e(ji+1,jj) &
400                        &  - ( rhd(ji,jj,1) + 1 ) * sshn_e(ji,jj) ) * hu_e(ji,jj) / e1u(ji,jj)
401                     zspgv = -grav * ( ( rhd(ji,jj+1,1) + 1 ) * sshn_e(ji,jj+1) &
402                        &  - ( rhd(ji,jj,1) + 1 ) * sshn_e(ji,jj) ) * hv_e(ji,jj) / e2v(ji,jj)
403                  ELSE
404                     zspgu = -grav * ( sshn_e(ji+1,jj) - sshn_e(ji,jj) ) * hu(ji,jj) / e1u(ji,jj)
405                     zspgv = -grav * ( sshn_e(ji,jj+1) - sshn_e(ji,jj) ) * hv(ji,jj) / e2v(ji,jj)
406                  ENDIF
407                  ! energy conserving formulation for planetary vorticity term
408                  zy1 = ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) ) / e1u(ji,jj)
409                  zy2 = ( zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) / e1u(ji,jj)
410                  zx1 = ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) ) / e2v(ji,jj)
411                  zx2 = ( zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) / e2v(ji,jj)
412                  zcubt = zfact2 * ( ff(ji  ,jj-1) * zy1 + ff(ji,jj) * zy2 )
413                  zcvbt =-zfact2 * ( ff(ji-1,jj  ) * zx1 + ff(ji,jj) * zx2 )
414                  ! after transports
415                  ua_e(ji,jj) = ( zub_e(ji,jj) + z2dt_e * ( zcubt + zspgu + zua(ji,jj) ) ) * umask(ji,jj,1)
416                  va_e(ji,jj) = ( zvb_e(ji,jj) + z2dt_e * ( zcvbt + zspgv + zva(ji,jj) ) ) * vmask(ji,jj,1)
417               END DO
418            END DO
419            !
420         ELSEIF ( ln_dynvor_ens ) THEN                    ! enstrophy conserving scheme
421            DO jj = 2, jpjm1
422               DO ji = fs_2, fs_jpim1   ! vector opt.
423                  ! surface pressure gradient
424                  IF( lk_vvl) THEN
425                     zspgu = -grav * ( ( rhd(ji+1,jj,1) + 1 ) * sshn_e(ji+1,jj) &
426                        &  - ( rhd(ji,jj,1) + 1 ) * sshn_e(ji,jj) ) * hu_e(ji,jj) / e1u(ji,jj)
427                     zspgv = -grav * ( ( rhd(ji,jj+1,1) + 1 ) * sshn_e(ji,jj+1) &
428                        &  - ( rhd(ji,jj,1) + 1 ) * sshn_e(ji,jj) ) * hv_e(ji,jj) / e2v(ji,jj)
429                  ELSE
430                     zspgu = -grav * ( sshn_e(ji+1,jj) - sshn_e(ji,jj) ) * hu(ji,jj) / e1u(ji,jj)
431                     zspgv = -grav * ( sshn_e(ji,jj+1) - sshn_e(ji,jj) ) * hv(ji,jj) / e2v(ji,jj)
432                  ENDIF
433                  ! enstrophy conserving formulation for planetary vorticity term
434                  zy1 = zfact1 * ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1)   &
435                                 + zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) / e1u(ji,jj)
436                  zx1 =-zfact1 * ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1)   &
437                                 + zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) / e2v(ji,jj)
438                  zcubt  = zy1 * ( ff(ji  ,jj-1) + ff(ji,jj) )
439                  zcvbt  = zx1 * ( ff(ji-1,jj  ) + ff(ji,jj) )
440                  ! after transports
441                  ua_e(ji,jj) = ( zub_e(ji,jj) + z2dt_e * ( zcubt + zspgu + zua(ji,jj) ) ) * umask(ji,jj,1)
442                  va_e(ji,jj) = ( zvb_e(ji,jj) + z2dt_e * ( zcvbt + zspgv + zva(ji,jj) ) ) * vmask(ji,jj,1)
443               END DO
444            END DO
445            !
446         ELSEIF ( ln_dynvor_een ) THEN                    ! energy and enstrophy conserving scheme
447            zfac25 = 0.25
448            DO jj = 2, jpjm1
449               DO ji = fs_2, fs_jpim1   ! vector opt.
450                  ! surface pressure gradient
451                  IF( lk_vvl) THEN
452                     zspgu = -grav * ( ( rhd(ji+1,jj,1) + 1 ) * sshn_e(ji+1,jj) &
453                        &  - ( rhd(ji,jj,1) + 1 ) * sshn_e(ji,jj) ) * hu_e(ji,jj) / e1u(ji,jj)
454                     zspgv = -grav * ( ( rhd(ji,jj+1,1) + 1 ) * sshn_e(ji,jj+1) &
455                        &  - ( rhd(ji,jj,1) + 1 ) * sshn_e(ji,jj) ) * hv_e(ji,jj) / e2v(ji,jj)
456                  ELSE
457                     zspgu = -grav * ( sshn_e(ji+1,jj) - sshn_e(ji,jj) ) * hu(ji,jj) / e1u(ji,jj)
458                     zspgv = -grav * ( sshn_e(ji,jj+1) - sshn_e(ji,jj) ) * hv(ji,jj) / e2v(ji,jj)
459                  ENDIF
460                  ! energy/enstrophy conserving formulation for planetary vorticity term
461                  zcubt = + zfac25 / e1u(ji,jj) * (  ftne(ji,jj  ) * zwy(ji  ,jj  ) + ftnw(ji+1,jj) * zwy(ji+1,jj  )   &
462                     &                             + ftse(ji,jj  ) * zwy(ji  ,jj-1) + ftsw(ji+1,jj) * zwy(ji+1,jj-1) )
463                  zcvbt = - zfac25 / e2v(ji,jj) * (  ftsw(ji,jj+1) * zwx(ji-1,jj+1) + ftse(ji,jj+1) * zwx(ji  ,jj+1)   &
464                     &                             + ftnw(ji,jj  ) * zwx(ji-1,jj  ) + ftne(ji,jj  ) * zwx(ji  ,jj  ) )
465                  ! after transports
466                  ua_e(ji,jj) = ( zub_e(ji,jj) + z2dt_e * ( zcubt + zspgu + zua(ji,jj) ) ) * umask(ji,jj,1)
467                  va_e(ji,jj) = ( zvb_e(ji,jj) + z2dt_e * ( zcvbt + zspgv + zva(ji,jj) ) ) * vmask(ji,jj,1)
468               END DO
469            END DO
470            !
471         ENDIF
472
473         ! Flather's boundary condition for the barotropic loop :
474         !         - Update sea surface height on each open boundary
475         !         - Correct the barotropic transports
476         IF( lk_obc )   CALL obc_fla_ts
477         IF( lk_bdy .or. lk_bdy_tides )   CALL bdy_dyn_fla
478
479         ! ... Boundary conditions on ua_e, va_e, ssha_e
480         CALL lbc_lnk( ua_e  , 'U', -1. )
481         CALL lbc_lnk( va_e  , 'V', -1. )
482         CALL lbc_lnk( ssha_e, 'T',  1. )
483
484         ! temporal sum
485         !-------------
486         IF( jit >= ibaro/2 ) THEN
487            zssha_b(:,:) = zssha_b(:,:) + ssha_e(:,:)
488            zua_b  (:,:) = zua_b  (:,:) + ua_e  (:,:)
489            zva_b  (:,:) = zva_b  (:,:) + va_e  (:,:) 
490         ENDIF
491
492         ! Time filter and swap of dynamics arrays
493         ! ---------------------------------------
494         IF( neuler == 0 .AND. jit == 1 ) THEN   ! Euler (forward) time stepping
495            sshb_e (:,:) = sshn_e(:,:)
496            zub_e  (:,:) = zun_e (:,:)
497            zvb_e  (:,:) = zvn_e (:,:)
498            sshn_e (:,:) = ssha_e(:,:)
499            zun_e  (:,:) = ua_e  (:,:)
500            zvn_e  (:,:) = va_e  (:,:)
501         ELSE                                        ! Asselin filtering
502            sshb_e (:,:) = atfp * ( sshb_e(:,:) + ssha_e(:,:) ) + atfp1 * sshn_e(:,:)
503            zub_e  (:,:) = atfp * ( zub_e  (:,:) + ua_e  (:,:) ) + atfp1 * zun_e  (:,:)
504            zvb_e  (:,:) = atfp * ( zvb_e  (:,:) + va_e  (:,:) ) + atfp1 * zvn_e  (:,:)
505            sshn_e (:,:) = ssha_e(:,:)
506            zun_e  (:,:) = ua_e  (:,:)
507            zvn_e  (:,:) = va_e  (:,:)
508         ENDIF
509
510         IF( lk_vvl ) THEN ! Variable volume
511
512            ! Sea surface elevation
513            ! ---------------------
514            DO jj = 1, jpjm1
515               DO ji = 1,jpim1
516
517                  ! Sea Surface Height at u-point before
518                  zsshun_e(ji,jj) = 0.5 * umask(ji,jj,1) / ( e1u(ji,jj) * e2u(ji,jj) )   &
519                     &              * ( e1t(ji  ,jj) * e2t(ji  ,jj) * sshn_e(ji  ,jj)    &
520                     &                + e1t(ji+1,jj) * e2t(ji+1,jj) * sshn_e(ji+1,jj) )
521
522                  ! Sea Surface Height at v-point before
523                  zsshvn_e(ji,jj) = 0.5 * vmask(ji,jj,1) / ( e1v(ji,jj) * e2v(ji,jj) )   &
524                     &              * ( e1t(ji,jj  ) * e2t(ji,jj  ) * sshn_e(ji,jj  )    &
525                     &                + e1t(ji,jj+1) * e2t(ji,jj+1) * sshn_e(ji,jj+1) )
526
527               END DO
528            END DO
529
530            ! Boundaries conditions
531            CALL lbc_lnk( zsshun_e, 'U', 1. )    ;    CALL lbc_lnk( zsshvn_e, 'V', 1. )
532
533            ! Scale factors at before and after time step
534            ! -------------------------------------------
535            CALL dom_vvl_sf( zsshun_e, 'U', zfse3un_e ) ;   CALL dom_vvl_sf( zsshvn_e, 'V', zfse3vn_e )
536
537            ! Ocean depth at U- and V-points
538            hu_e(:,:) = 0.e0
539            hv_e(:,:) = 0.e0
540
541            DO jk = 1, jpk
542               hu_e(:,:) = hu_e(:,:) + zfse3un_e(:,:,jk) * umask(:,:,jk)
543               hv_e(:,:) = hv_e(:,:) + zfse3vn_e(:,:,jk) * vmask(:,:,jk)
544            END DO
545
546         ENDIF ! End variable volume case
547
548         !                                                 ! ==================== !
549      END DO                                               !        end loop      !
550      !                                                    ! ==================== !
551
552
553      ! Time average of after barotropic variables
554      zcoef =  1.e0 / ( ibaro + 1 )
555      zssha_b(:,:) = zcoef * zssha_b(:,:) 
556      zua_b  (:,:) = zcoef *  zua_b (:,:) 
557      zva_b  (:,:) = zcoef *  zva_b (:,:) 
558#if defined key_obc
559         IF( lp_obc_east  )   sshfoe_b(:,:) = zcoef * sshfoe_b(:,:)
560         IF( lp_obc_west  )   sshfow_b(:,:) = zcoef * sshfow_b(:,:)
561         IF( lp_obc_north )   sshfon_b(:,:) = zcoef * sshfon_b(:,:)
562         IF( lp_obc_south )   sshfos_b(:,:) = zcoef * sshfos_b(:,:)
563#endif
564     
565
566      ! ---------------------------------------------------------------------------
567      ! Phase 3 : Update sea surface height from time averaged barotropic variables
568      ! ---------------------------------------------------------------------------
569
570 
571      ! Horizontal divergence of time averaged barotropic transports
572      !-------------------------------------------------------------
573      IF( .NOT. lk_vvl ) THEN
574         DO jj = 2, jpjm1
575            DO ji = fs_2, fs_jpim1   ! vector opt.
576               zhdiv(ji,jj) = ( e2u(ji,jj) * un_b(ji,jj) - e2u(ji-1,jj  ) * un_b(ji-1,jj  )     &
577                  &            +e1v(ji,jj) * vn_b(ji,jj) - e1v(ji  ,jj-1) * vn_b(ji  ,jj-1) )   &
578                  &           / ( e1t(ji,jj) * e2t(ji,jj) )
579            END DO
580         END DO
581      ENDIF
582
583#if defined key_obc && ! defined key_vvl
584      ! open boundaries (div must be zero behind the open boundary)
585      !  mpp remark: The zeroing of hdiv can probably be extended to 1->jpi/jpj for the correct row/column
586      IF( lp_obc_east  )   zhdiv(nie0p1:nie1p1,nje0  :nje1)   = 0.e0    ! east
587      IF( lp_obc_west  )   zhdiv(niw0  :niw1  ,njw0  :njw1)   = 0.e0    ! west
588      IF( lp_obc_north )   zhdiv(nin0  :nin1  ,njn0p1:njn1p1) = 0.e0    ! north
589      IF( lp_obc_south )   zhdiv(nis0  :nis1  ,njs0  :njs1)   = 0.e0    ! south
590#endif
591
592      IF ( lk_bdy .or. lk_bdy_tides ) THEN
593         DO jj = 1, jpj
594           DO ji = 1, jpi
595             zhdiv(ji,jj) = zhdiv(ji,jj)*bdytmask(ji,jj)
596           END DO
597         END DO
598      ENDIF
599
600      ! sea surface height
601      !-------------------
602      IF( lk_vvl ) THEN
603         sshbb(:,:) = sshb(:,:)
604         sshb (:,:) = sshn(:,:)
605         sshn (:,:) = ssha(:,:)
606      ELSE
607         sshb (:,:) = sshn(:,:)
608         sshn(:,:) = (  sshb_b(:,:) - z2dt_b * ( zraur * emp(:,:) + zhdiv(:,:) )  ) * tmask(:,:,1)
609      ENDIF
610
611      ! ... Boundary conditions on sshn
612      IF( .NOT. lk_obc ) CALL lbc_lnk( sshn, 'T', 1. )
613
614
615      ! -----------------------------------------------------------------------------
616      ! Phase 4. Coupling between general trend and barotropic estimates - (2nd step)
617      ! -----------------------------------------------------------------------------
618
619      ! Swap on time averaged barotropic variables
620      ! ------------------------------------------
621      sshb_b(:,:) = sshn_b (:,:)
622      IF ( kt == nit000 ) zssha_b(:,:) = sshn(:,:)
623      sshn_b(:,:) = zssha_b(:,:)
624      un_b  (:,:) = zua_b  (:,:) 
625      vn_b  (:,:) = zva_b  (:,:) 
626   
627      ! add time averaged barotropic coriolis and surface pressure gradient
628      ! terms to the general momentum trend
629      ! --------------------------------------------------------------------
630      DO jk=1,jpkm1
631         ua(:,:,jk) = ua(:,:,jk) + hur(:,:) * ( zua_b(:,:) - zub(:,:) ) / z2dt_b
632         va(:,:,jk) = va(:,:,jk) + hvr(:,:) * ( zva_b(:,:) - zvb(:,:) ) / z2dt_b
633      END DO
634
635      ! write filtered free surface arrays in restart file
636      ! --------------------------------------------------
637      IF( lrst_oce )   CALL ts_rst( kt, 'WRITE' )
638
639      ! print sum trends (used for debugging)
640      IF( ln_ctl )     CALL prt_ctl( tab2d_1=sshn, clinfo1=' ssh      : ', mask1=tmask )
641      !
642   END SUBROUTINE dyn_spg_ts
643
644
645   SUBROUTINE ts_rst( kt, cdrw )
646      !!---------------------------------------------------------------------
647      !!                   ***  ROUTINE ts_rst  ***
648      !!
649      !! ** Purpose : Read or write time-splitting arrays in restart file
650      !!----------------------------------------------------------------------
651      INTEGER         , INTENT(in) ::   kt         ! ocean time-step
652      CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
653      !
654      INTEGER ::  ji, jk        ! dummy loop indices
655      !!----------------------------------------------------------------------
656      !
657      IF( TRIM(cdrw) == 'READ' ) THEN
658         IF( iom_varid( numror, 'sshn', ldstop = .FALSE. ) > 0 ) THEN
659            CALL iom_get( numror, jpdom_autoglo, 'sshb'  , sshb(:,:)   )
660            CALL iom_get( numror, jpdom_autoglo, 'sshn'  , sshn(:,:)   )
661            IF( neuler == 0 ) sshb(:,:) = sshn(:,:)
662         ELSE
663            IF( nn_rstssh == 1 ) THEN 
664               sshb(:,:) = 0.e0
665               sshn(:,:) = 0.e0
666            ENDIF
667         ENDIF
668         IF( iom_varid( numror, 'sshn_b', ldstop = .FALSE. ) > 0 ) THEN
669            CALL iom_get( numror, jpdom_autoglo, 'sshb_b', sshb_b(:,:) )   ! free surface issued
670            CALL iom_get( numror, jpdom_autoglo, 'sshn_b', sshn_b(:,:) )   ! from time-splitting loop
671            CALL iom_get( numror, jpdom_autoglo, 'un_b'  , un_b  (:,:) )   ! horizontal transports issued
672            CALL iom_get( numror, jpdom_autoglo, 'vn_b'  , vn_b  (:,:) )   ! from barotropic loop
673            IF( neuler == 0 ) sshb_b(:,:) = sshn_b(:,:)
674         ELSE
675            sshb_b(:,:) = sshb(:,:)
676            sshn_b(:,:) = sshn(:,:)
677            un_b  (:,:) = 0.e0
678            vn_b  (:,:) = 0.e0
679            ! vertical sum
680            IF( lk_vopt_loop ) THEN          ! vector opt., forced unroll
681               DO jk = 1, jpkm1
682                  DO ji = 1, jpij
683                     un_b(ji,1) = un_b(ji,1) + fse3u(ji,1,jk) * un(ji,1,jk)
684                     vn_b(ji,1) = vn_b(ji,1) + fse3v(ji,1,jk) * vn(ji,1,jk)
685                  END DO
686               END DO
687            ELSE                             ! No  vector opt.
688               DO jk = 1, jpkm1
689                  un_b(:,:) = un_b(:,:) + fse3u(:,:,jk) * un(:,:,jk)
690                  vn_b(:,:) = vn_b(:,:) + fse3v(:,:,jk) * vn(:,:,jk)
691               END DO
692            ENDIF
693         ENDIF
694      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN
695         CALL iom_rstput( kt, nitrst, numrow, 'sshb'  , sshb  (:,:) )
696         CALL iom_rstput( kt, nitrst, numrow, 'sshn'  , sshn  (:,:) )
697         CALL iom_rstput( kt, nitrst, numrow, 'sshb_b', sshb_b(:,:) )   ! free surface issued
698         CALL iom_rstput( kt, nitrst, numrow, 'sshn_b', sshn_b(:,:) )   ! from barotropic loop
699         CALL iom_rstput( kt, nitrst, numrow, 'un_b'  , un_b  (:,:) )   ! horizontal transports issued
700         CALL iom_rstput( kt, nitrst, numrow, 'vn_b'  , vn_b  (:,:) )   ! from barotropic loop
701      ENDIF
702      !
703   END SUBROUTINE ts_rst
704
705#else
706   !!----------------------------------------------------------------------
707   !!   Default case :   Empty module   No standart free surface cst volume
708   !!----------------------------------------------------------------------
709CONTAINS
710   SUBROUTINE dyn_spg_ts( kt )       ! Empty routine
711      WRITE(*,*) 'dyn_spg_ts: You should not have seen this print! error?', kt
712   END SUBROUTINE dyn_spg_ts
713   SUBROUTINE ts_rst( kt, cdrw )     ! Empty routine
714      INTEGER         , INTENT(in) ::   kt         ! ocean time-step
715      CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
716      WRITE(*,*) 'ts_rst    : You should not have seen this print! error?', kt, cdrw
717   END SUBROUTINE ts_rst   
718#endif
719   
720   !!======================================================================
721END MODULE dynspg_ts
Note: See TracBrowser for help on using the repository browser.