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

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

Suppress jki routines and associated key_mpp_omp

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