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.
sshwzv.F90 in branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 3211

Last change on this file since 3211 was 3211, checked in by spickles2, 12 years ago

Stephen Pickles, 11 Dec 2011

Commit to bring the rest of the DCSE NEMO development branch
in line with the latest development version. This includes
array index re-ordering of all OPA_SRC/.

  • Property svn:keywords set to Id
File size: 24.0 KB
Line 
1MODULE sshwzv   
2   !!==============================================================================
3   !!                       ***  MODULE  sshwzv  ***
4   !! Ocean dynamics : sea surface height and vertical velocity
5   !!==============================================================================
6   !! History :  3.1  !  2009-02  (G. Madec, M. Leclair)  Original code
7   !!            3.3  !  2010-04  (M. Leclair, G. Madec)  modified LF-RA
8   !!             -   !  2010-05  (K. Mogensen, A. Weaver, M. Martin, D. Lea) Assimilation interface
9   !!             -   !  2010-09  (D.Storkey and E.O'Dea) bug fixes for BDY module
10   !!----------------------------------------------------------------------
11
12   !!----------------------------------------------------------------------
13   !!   ssh_wzv        : after ssh & now vertical velocity
14   !!   ssh_nxt        : filter ans swap the ssh arrays
15   !!----------------------------------------------------------------------
16   USE oce             ! ocean dynamics and tracers variables
17   USE dom_oce         ! ocean space and time domain variables
18   USE sbc_oce         ! surface boundary condition: ocean
19   USE domvvl          ! Variable volume
20   USE divcur          ! hor. divergence and curl      (div & cur routines)
21   USE iom             ! I/O library
22   USE restart         ! only for lrst_oce
23   USE in_out_manager  ! I/O manager
24   USE prtctl          ! Print control
25   USE phycst
26   USE lbclnk          ! ocean lateral boundary condition (or mpp link)
27   USE lib_mpp         ! MPP library
28   USE obc_par         ! open boundary cond. parameter
29   USE obc_oce
30   USE bdy_oce
31   USE diaar5, ONLY:   lk_diaar5
32   USE iom
33   USE sbcrnf, ONLY: h_rnf, nk_rnf   ! River runoff
34#if defined key_agrif
35   USE agrif_opa_update
36   USE agrif_opa_interp
37#endif
38#if defined key_asminc   
39   USE asminc          ! Assimilation increment
40#endif
41
42   IMPLICIT NONE
43   PRIVATE
44
45   PUBLIC   ssh_wzv    ! called by step.F90
46   PUBLIC   ssh_nxt    ! called by step.F90
47
48   !! * Control permutation of array indices
49#  include "oce_ftrans.h90"
50#  include "dom_oce_ftrans.h90"
51#  include "sbc_oce_ftrans.h90"
52#  include "domvvl_ftrans.h90"
53#  include "obc_oce_ftrans.h90"
54#if defined key_asminc 
55#  include "asminc_ftrans.h90"
56#endif
57
58   !! * Substitutions
59#  include "domzgr_substitute.h90"
60#  include "vectopt_loop_substitute.h90"
61   !!----------------------------------------------------------------------
62   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
63   !! $Id$
64   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
65   !!----------------------------------------------------------------------
66CONTAINS
67
68   SUBROUTINE ssh_wzv( kt ) 
69      !!----------------------------------------------------------------------
70      !!                ***  ROUTINE ssh_wzv  ***
71      !!                   
72      !! ** Purpose :   compute the after ssh (ssha), the now vertical velocity
73      !!              and update the now vertical coordinate (lk_vvl=T).
74      !!
75      !! ** Method  : - Using the incompressibility hypothesis, the vertical
76      !!      velocity is computed by integrating the horizontal divergence 
77      !!      from the bottom to the surface minus the scale factor evolution.
78      !!        The boundary conditions are w=0 at the bottom (no flux) and.
79      !!
80      !! ** action  :   ssha    : after sea surface height
81      !!                wn      : now vertical velocity
82      !!                sshu_a, sshv_a, sshf_a  : after sea surface height (lk_vvl=T)
83      !!                hu, hv, hur, hvr        : ocean depth and its inverse at u-,v-points
84      !!
85      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
86      !!----------------------------------------------------------------------
87      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
88      USE oce     , ONLY:   z3d   => ta                           ! ta used as 3D workspace
89      USE wrk_nemo, ONLY:   zhdiv => wrk_2d_1 , z2d => wrk_2d_2   ! 2D workspace
90      !! DCSE_NEMO: need additional directives for renamed module variables
91!FTRANS z3d :I :I :z
92      !
93      INTEGER, INTENT(in) ::   kt   ! time step
94      !
95      INTEGER  ::   ji, jj, jk   ! dummy loop indices
96      REAL(wp) ::   zcoefu, zcoefv, zcoeff, z2dt, z1_2dt, z1_rau0   ! local scalars
97      !!----------------------------------------------------------------------
98
99      IF( wrk_in_use(2, 1,2) ) THEN
100         CALL ctl_stop('ssh_wzv: requested workspace arrays unavailable')   ;   RETURN
101      ENDIF
102
103      IF( kt == nit000 ) THEN
104         !
105         IF(lwp) WRITE(numout,*)
106         IF(lwp) WRITE(numout,*) 'ssh_wzv : after sea surface height and now vertical velocity '
107         IF(lwp) WRITE(numout,*) '~~~~~~~ '
108         !
109         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
110         !
111         IF( lk_vvl ) THEN                    ! before and now Sea SSH at u-, v-, f-points (vvl case only)
112            DO jj = 1, jpjm1
113               DO ji = 1, jpim1                    ! caution: use of Vector Opt. not possible
114#if defined key_z_first
115                  zcoefu = 0.5  * umask_1(ji,jj) / ( e1u(ji,jj) * e2u(ji,jj) )
116                  zcoefv = 0.5  * vmask_1(ji,jj) / ( e1v(ji,jj) * e2v(ji,jj) )
117                  zcoeff = 0.25 * umask_1(ji,jj) * umask_1(ji,jj+1)
118#else
119                  zcoefu = 0.5  * umask(ji,jj,1) / ( e1u(ji,jj) * e2u(ji,jj) )
120                  zcoefv = 0.5  * vmask(ji,jj,1) / ( e1v(ji,jj) * e2v(ji,jj) )
121                  zcoeff = 0.25 * umask(ji,jj,1) * umask(ji,jj+1,1)
122#endif
123                  sshu_b(ji,jj) = zcoefu * ( e1t(ji  ,jj) * e2t(ji  ,jj) * sshb(ji  ,jj)     &
124                     &                     + e1t(ji+1,jj) * e2t(ji+1,jj) * sshb(ji+1,jj) )
125                  sshv_b(ji,jj) = zcoefv * ( e1t(ji,jj  ) * e2t(ji,jj  ) * sshb(ji,jj  )     &
126                     &                     + e1t(ji,jj+1) * e2t(ji,jj+1) * sshb(ji,jj+1) )
127                  sshu_n(ji,jj) = zcoefu * ( e1t(ji  ,jj) * e2t(ji  ,jj) * sshn(ji  ,jj)     &
128                     &                     + e1t(ji+1,jj) * e2t(ji+1,jj) * sshn(ji+1,jj) )
129                  sshv_n(ji,jj) = zcoefv * ( e1t(ji,jj  ) * e2t(ji,jj  ) * sshn(ji,jj  )     &
130                     &                     + e1t(ji,jj+1) * e2t(ji,jj+1) * sshn(ji,jj+1) )
131               END DO
132            END DO
133            CALL lbc_lnk( sshu_b, 'U', 1. )   ;   CALL lbc_lnk( sshu_n, 'U', 1. )
134            CALL lbc_lnk( sshv_b, 'V', 1. )   ;   CALL lbc_lnk( sshv_n, 'V', 1. )
135            DO jj = 1, jpjm1
136               DO ji = 1, jpim1      ! NO Vector Opt.
137#if defined key_z_first
138                  sshf_n(ji,jj) = 0.5  * umask_1(ji,jj) * umask_1(ji,jj+1)                   &
139                       &               / ( e1f(ji,jj  ) * e2f(ji,jj  ) )                     &
140                       &               * ( e1u(ji,jj  ) * e2u(ji,jj  ) * sshu_n(ji,jj  )     &
141                       &                 + e1u(ji,jj+1) * e2u(ji,jj+1) * sshu_n(ji,jj+1) )
142#else
143                  sshf_n(ji,jj) = 0.5  * umask(ji,jj,1) * umask(ji,jj+1,1)                   &
144                       &               / ( e1f(ji,jj  ) * e2f(ji,jj  ) )                     &
145                       &               * ( e1u(ji,jj  ) * e2u(ji,jj  ) * sshu_n(ji,jj  )     &
146                       &                 + e1u(ji,jj+1) * e2u(ji,jj+1) * sshu_n(ji,jj+1) )
147#endif
148               END DO
149            END DO
150            CALL lbc_lnk( sshf_n, 'F', 1. )
151         ENDIF
152         !
153      ENDIF
154
155      !                                           !------------------------------------------!
156      IF( lk_vvl ) THEN                           !  Regridding: Update Now Vertical coord.  !   (only in vvl case)
157         !                                        !------------------------------------------!
158#if defined key_z_first
159         fsdept(:,:,1:jpkm1) = fsdept_n(:,:,1:jpkm1)   ! now local depths stored in fsdep. arrays
160         fsdepw(:,:,1:jpkm1) = fsdepw_n(:,:,1:jpkm1)
161         fsde3w(:,:,1:jpkm1) = fsde3w_n(:,:,1:jpkm1)
162         !
163         fse3t (:,:,1:jpkm1) = fse3t_n (:,:,1:jpkm1)   ! vertical scale factors stored in fse3. arrays
164         fse3u (:,:,1:jpkm1) = fse3u_n (:,:,1:jpkm1)
165         fse3v (:,:,1:jpkm1) = fse3v_n (:,:,1:jpkm1)
166         fse3f (:,:,1:jpkm1) = fse3f_n (:,:,1:jpkm1)
167         fse3w (:,:,1:jpkm1) = fse3w_n (:,:,1:jpkm1)
168         fse3uw(:,:,1:jpkm1) = fse3uw_n(:,:,1:jpkm1)
169         fse3vw(:,:,1:jpkm1) = fse3vw_n(:,:,1:jpkm1)
170#else
171         DO jk = 1, jpkm1
172            fsdept(:,:,jk) = fsdept_n(:,:,jk)         ! now local depths stored in fsdep. arrays
173            fsdepw(:,:,jk) = fsdepw_n(:,:,jk)
174            fsde3w(:,:,jk) = fsde3w_n(:,:,jk)
175            !
176            fse3t (:,:,jk) = fse3t_n (:,:,jk)         ! vertical scale factors stored in fse3. arrays
177            fse3u (:,:,jk) = fse3u_n (:,:,jk)
178            fse3v (:,:,jk) = fse3v_n (:,:,jk)
179            fse3f (:,:,jk) = fse3f_n (:,:,jk)
180            fse3w (:,:,jk) = fse3w_n (:,:,jk)
181            fse3uw(:,:,jk) = fse3uw_n(:,:,jk)
182            fse3vw(:,:,jk) = fse3vw_n(:,:,jk)
183         END DO
184#endif
185         !
186         hu(:,:) = hu_0(:,:) + sshu_n(:,:)            ! now ocean depth (at u- and v-points)
187         hv(:,:) = hv_0(:,:) + sshv_n(:,:)
188         !                                            ! now masked inverse of the ocean depth (at u- and v-points)
189#if defined key_z_first
190         hur(:,:) = umask_1(:,:) / ( hu(:,:) + 1._wp - umask_1(:,:) )
191         hvr(:,:) = vmask_1(:,:) / ( hv(:,:) + 1._wp - vmask_1(:,:) )
192#else
193         hur(:,:) = umask(:,:,1) / ( hu(:,:) + 1._wp - umask(:,:,1) )
194         hvr(:,:) = vmask(:,:,1) / ( hv(:,:) + 1._wp - vmask(:,:,1) )
195#endif
196         !
197      ENDIF
198      !
199      CALL div_cur( kt )                              ! Horizontal divergence & Relative vorticity
200      !
201      z2dt = 2._wp * rdt                              ! set time step size (Euler/Leapfrog)
202      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
203
204      !                                           !------------------------------!
205      !                                           !   After Sea Surface Height   !
206      !                                           !------------------------------!
207      zhdiv(:,:) = 0._wp
208#if defined key_z_first
209      DO jj = 1, jpj
210         DO ji = 1, jpi
211            DO jk = 1, jpkm1                           ! Horizontal divergence of barotropic transports
212               zhdiv(ji,jj) = zhdiv(ji,jj) + fse3t(ji,jj,jk) * hdivn(ji,jj,jk)
213            END DO
214         END DO
215      END DO
216#else
217      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
218        zhdiv(:,:) = zhdiv(:,:) + fse3t(:,:,jk) * hdivn(:,:,jk)
219      END DO
220#endif
221      !                                                ! Sea surface elevation time stepping
222      ! In forward Euler time stepping case, the same formulation as in the leap-frog case can be used
223      ! because emp_b field is initialized with the vlaues of emp field. Hence, 0.5 * ( emp + emp_b ) = emp
224      z1_rau0 = 0.5 / rau0
225#if defined key_z_first
226      ssha(:,:) = (  sshb(:,:) - z2dt * ( z1_rau0 * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * tmask_1(:,:)
227#else
228      ssha(:,:) = (  sshb(:,:) - z2dt * ( z1_rau0 * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * tmask(:,:,1)
229#endif
230
231#if defined key_agrif
232      CALL agrif_ssh( kt )
233#endif
234#if defined key_obc
235      IF( Agrif_Root() ) THEN
236         ssha(:,:) = ssha(:,:) * obctmsk(:,:)
237         CALL lbc_lnk( ssha, 'T', 1. )                 ! absolutly compulsory !! (jmm)
238      ENDIF
239#endif
240#if defined key_bdy
241      ssha(:,:) = ssha(:,:) * bdytmask(:,:)
242      CALL lbc_lnk( ssha, 'T', 1. ) 
243#endif
244
245      !                                                ! Sea Surface Height at u-,v- and f-points (vvl case only)
246      IF( lk_vvl ) THEN                                ! (required only in key_vvl case)
247         DO jj = 1, jpjm1
248            DO ji = 1, jpim1      ! NO Vector Opt.
249#if defined key_z_first
250               sshu_a(ji,jj) = 0.5  * umask_1(ji,jj) / ( e1u(ji  ,jj) * e2u(ji  ,jj) )                   &
251                  &                                  * ( e1t(ji  ,jj) * e2t(ji  ,jj) * ssha(ji  ,jj)     &
252                  &                                    + e1t(ji+1,jj) * e2t(ji+1,jj) * ssha(ji+1,jj) )
253               sshv_a(ji,jj) = 0.5  * vmask_1(ji,jj) / ( e1v(ji,jj  ) * e2v(ji,jj  ) )                   &
254                  &                                  * ( e1t(ji,jj  ) * e2t(ji,jj  ) * ssha(ji,jj  )     &
255                  &                                    + e1t(ji,jj+1) * e2t(ji,jj+1) * ssha(ji,jj+1) )
256#else
257               sshu_a(ji,jj) = 0.5  * umask(ji,jj,1) / ( e1u(ji  ,jj) * e2u(ji  ,jj) )                   &
258                  &                                  * ( e1t(ji  ,jj) * e2t(ji  ,jj) * ssha(ji  ,jj)     &
259                  &                                    + e1t(ji+1,jj) * e2t(ji+1,jj) * ssha(ji+1,jj) )
260               sshv_a(ji,jj) = 0.5  * vmask(ji,jj,1) / ( e1v(ji,jj  ) * e2v(ji,jj  ) )                   &
261                  &                                  * ( e1t(ji,jj  ) * e2t(ji,jj  ) * ssha(ji,jj  )     &
262                  &                                    + e1t(ji,jj+1) * e2t(ji,jj+1) * ssha(ji,jj+1) )
263#endif
264            END DO
265         END DO
266         CALL lbc_lnk( sshu_a, 'U', 1. )   ;   CALL lbc_lnk( sshv_a, 'V', 1. )      ! Boundaries conditions
267      ENDIF
268     
269#if defined key_asminc
270      !                                                ! Include the IAU weighted SSH increment
271      IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN
272         CALL ssh_asm_inc( kt )
273         ssha(:,:) = ssha(:,:) + z2dt * ssh_iau(:,:)
274      ENDIF
275#endif
276
277      !                                           !------------------------------!
278      !                                           !     Now Vertical Velocity    !
279      !                                           !------------------------------!
280      z1_2dt = 1.e0 / z2dt
281#if defined key_z_first
282      DO jj = 1, jpj
283         DO ji = 1, jpi
284            DO jk = jpkm1, 1, -1                      ! integrate from the bottom the hor. divergence
285                wn(ji,jj,jk) = wn(ji,jj,jk+1)                               &
286                   &         -   fse3t_n(ji,jj,jk) * hdivn(ji,jj,jk)        &
287                   &         - ( fse3t_a(ji,jj,jk) - fse3t_b(ji,jj,jk) )    &
288                   &            * tmask(ji,jj,jk) * z1_2dt
289#if defined key_bdy
290                wn(ji,jj,jk) = wn(ji,jj,jk) * bdytmask(ji,jj)
291#endif
292            END DO
293         END DO
294      END DO
295#else
296      DO jk = jpkm1, 1, -1                             ! integrate from the bottom the hor. divergence
297         ! - ML - need 3 lines here because replacement of fse3t by its expression yields too long lines otherwise
298         wn(:,:,jk) = wn(:,:,jk+1) -   fse3t_n(:,:,jk) * hdivn(:,:,jk)        &
299            &                      - ( fse3t_a(:,:,jk) - fse3t_b(:,:,jk) )    &
300            &                         * tmask(:,:,jk) * z1_2dt
301#if defined key_bdy
302         wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
303#endif
304      END DO
305#endif
306
307      !                                           !------------------------------!
308      !                                           !           outputs            !
309      !                                           !------------------------------!
310      CALL iom_put( "woce", wn                    )   ! vertical velocity
311      CALL iom_put( "ssh" , sshn                  )   ! sea surface height
312      CALL iom_put( "ssh2", sshn(:,:) * sshn(:,:) )   ! square of sea surface height
313      IF( lk_diaar5 ) THEN                            ! vertical mass transport & its square value
314         ! Caution: in the VVL case, it only correponds to the baroclinic mass transport.
315         z2d(:,:) = rau0 * e1t(:,:) * e2t(:,:)
316#if defined key_z_first
317         DO jj = 1, jpj
318            DO ji = 1, jpi
319               DO jk = 1, jpk
320                  z3d(ji,jj,jk) = wn(ji,jj,jk) * z2d(ji,jj)
321               END DO
322            END DO
323         END DO
324#else
325         DO jk = 1, jpk
326            z3d(:,:,jk) = wn(:,:,jk) * z2d(:,:)
327         END DO
328#endif
329         CALL iom_put( "w_masstr" , z3d                     ) 
330         CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) )
331      ENDIF
332      !
333      IF(ln_ctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask, ovlap=1 )
334      !
335      IF( wrk_not_released(2, 1,2) )   CALL ctl_stop('ssh_wzv: failed to release workspace arrays')
336      !
337   END SUBROUTINE ssh_wzv
338
339
340   SUBROUTINE ssh_nxt( kt )
341      !!----------------------------------------------------------------------
342      !!                    ***  ROUTINE ssh_nxt  ***
343      !!
344      !! ** Purpose :   achieve the sea surface  height time stepping by
345      !!              applying Asselin time filter and swapping the arrays
346      !!              ssha  already computed in ssh_wzv 
347      !!
348      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
349      !!              from the filter, see Leclair and Madec 2010) and swap :
350      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
351      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
352      !!                sshn = ssha
353      !!
354      !! ** action  : - sshb, sshn   : before & now sea surface height
355      !!                               ready for the next time step
356      !!
357      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
358      !!----------------------------------------------------------------------
359      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
360      !!
361      INTEGER  ::   ji, jj   ! dummy loop indices
362      REAL(wp) ::   zec      ! temporary scalar
363      !!----------------------------------------------------------------------
364
365      IF( kt == nit000 ) THEN
366         IF(lwp) WRITE(numout,*)
367         IF(lwp) WRITE(numout,*) 'ssh_nxt : next sea surface height (Asselin time filter + swap)'
368         IF(lwp) WRITE(numout,*) '~~~~~~~ '
369      ENDIF
370
371      !                       !--------------------------!
372      IF( lk_vvl ) THEN       !  Variable volume levels  !     (ssh at t-, u-, v, f-points)
373         !                    !--------------------------!
374         !
375         IF( neuler == 0 .AND. kt == nit000 ) THEN    !** Euler time-stepping at first time-step : no filter
376            sshn  (:,:) = ssha  (:,:)                       ! now <-- after  (before already = now)
377            sshu_n(:,:) = sshu_a(:,:)
378            sshv_n(:,:) = sshv_a(:,:)
379            DO jj = 1, jpjm1                                ! ssh now at f-point
380               DO ji = 1, jpim1      ! NO Vector Opt.
381#if defined key_z_first
382                  sshf_n(ji,jj) = 0.5  * umask_1(ji,jj) * umask_1(ji,jj+1)                 &
383                     &               / ( e1f(ji,jj  ) * e2f(ji,jj  ) )                     &
384                     &               * ( e1u(ji,jj  ) * e2u(ji,jj  ) * sshu_n(ji,jj  )     &
385                     &                 + e1u(ji,jj+1) * e2u(ji,jj+1) * sshu_n(ji,jj+1) )
386#else
387                  sshf_n(ji,jj) = 0.5  * umask(ji,jj,1) * umask(ji,jj+1,1)                 &
388                     &               / ( e1f(ji,jj  ) * e2f(ji,jj  ) )                     &
389                     &               * ( e1u(ji,jj  ) * e2u(ji,jj  ) * sshu_n(ji,jj  )     &
390                     &                 + e1u(ji,jj+1) * e2u(ji,jj+1) * sshu_n(ji,jj+1) )
391#endif
392               END DO
393            END DO
394            CALL lbc_lnk( sshf_n, 'F', 1. )                 ! Boundaries conditions
395            !
396         ELSE                                         !** Leap-Frog time-stepping: Asselin filter + swap
397            zec = atfp * rdt / rau0
398            DO jj = 1, jpj
399               DO ji = 1, jpi                               ! before <-- now filtered
400#if defined key_z_first
401                  sshb  (ji,jj) = sshn  (ji,jj) + atfp * ( sshb(ji,jj) - 2 * sshn(ji,jj) + ssha(ji,jj) )   &
402                     &                          - zec  * ( emp_b(ji,jj) - emp(ji,jj) ) * tmask_1(ji,jj)
403#else
404                  sshb  (ji,jj) = sshn  (ji,jj) + atfp * ( sshb(ji,jj) - 2 * sshn(ji,jj) + ssha(ji,jj) )   &
405                     &                          - zec  * ( emp_b(ji,jj) - emp(ji,jj) ) * tmask(ji,jj,1)
406#endif
407                  sshn  (ji,jj) = ssha  (ji,jj)             ! now <-- after
408                  sshu_n(ji,jj) = sshu_a(ji,jj)
409                  sshv_n(ji,jj) = sshv_a(ji,jj)
410               END DO
411            END DO
412            DO jj = 1, jpjm1                                ! ssh now at f-point
413               DO ji = 1, jpim1      ! NO Vector Opt.
414#if defined key_z_first
415                  sshf_n(ji,jj) = 0.5  * umask_1(ji,jj) * umask_1(ji,jj+1)                 &
416                     &               / ( e1f(ji,jj  ) * e2f(ji,jj  ) )                     &
417                     &               * ( e1u(ji,jj  ) * e2u(ji,jj  ) * sshu_n(ji,jj  )     &
418                     &                 + e1u(ji,jj+1) * e2u(ji,jj+1) * sshu_n(ji,jj+1) )
419#else
420                  sshf_n(ji,jj) = 0.5  * umask(ji,jj,1) * umask(ji,jj+1,1)                 &
421                     &               / ( e1f(ji,jj  ) * e2f(ji,jj  ) )                     &
422                     &               * ( e1u(ji,jj  ) * e2u(ji,jj  ) * sshu_n(ji,jj  )     &
423                     &                 + e1u(ji,jj+1) * e2u(ji,jj+1) * sshu_n(ji,jj+1) )
424#endif
425               END DO
426            END DO
427            CALL lbc_lnk( sshf_n, 'F', 1. )                 ! Boundaries conditions
428            !
429            DO jj = 1, jpjm1                                ! ssh before at u- & v-points
430               DO ji = 1, jpim1      ! NO Vector Opt.
431#if defined key_z_first
432                  sshu_b(ji,jj) = 0.5  * umask_1(ji,jj) / ( e1u(ji  ,jj) * e2u(ji  ,jj) )                   &
433                     &                                  * ( e1t(ji  ,jj) * e2t(ji  ,jj) * sshb(ji  ,jj)     &
434                     &                                    + e1t(ji+1,jj) * e2t(ji+1,jj) * sshb(ji+1,jj) )
435                  sshv_b(ji,jj) = 0.5  * vmask_1(ji,jj) / ( e1v(ji,jj  ) * e2v(ji,jj  ) )                   &
436                     &                                  * ( e1t(ji,jj  ) * e2t(ji,jj  ) * sshb(ji,jj  )     &
437                     &                                    + e1t(ji,jj+1) * e2t(ji,jj+1) * sshb(ji,jj+1) )
438#else
439                  sshu_b(ji,jj) = 0.5  * umask(ji,jj,1) / ( e1u(ji  ,jj) * e2u(ji  ,jj) )                   &
440                     &                                  * ( e1t(ji  ,jj) * e2t(ji  ,jj) * sshb(ji  ,jj)     &
441                     &                                    + e1t(ji+1,jj) * e2t(ji+1,jj) * sshb(ji+1,jj) )
442                  sshv_b(ji,jj) = 0.5  * vmask(ji,jj,1) / ( e1v(ji,jj  ) * e2v(ji,jj  ) )                   &
443                     &                                  * ( e1t(ji,jj  ) * e2t(ji,jj  ) * sshb(ji,jj  )     &
444                     &                                    + e1t(ji,jj+1) * e2t(ji,jj+1) * sshb(ji,jj+1) )
445#endif
446               END DO
447            END DO
448            CALL lbc_lnk( sshu_b, 'U', 1. )
449            CALL lbc_lnk( sshv_b, 'V', 1. )            !  Boundaries conditions
450            !
451         ENDIF
452         !                    !--------------------------!
453      ELSE                    !        fixed levels      !     (ssh at t-point only)
454         !                    !--------------------------!
455         !
456         IF( neuler == 0 .AND. kt == nit000 ) THEN    !** Euler time-stepping at first time-step : no filter
457            sshn(:,:) = ssha(:,:)                           ! now <-- after  (before already = now)
458            !
459         ELSE                                               ! Leap-Frog time-stepping: Asselin filter + swap
460            DO jj = 1, jpj
461               DO ji = 1, jpi                               ! before <-- now filtered
462                  sshb(ji,jj) = sshn(ji,jj) + atfp * ( sshb(ji,jj) - 2 * sshn(ji,jj) + ssha(ji,jj) )
463                  sshn(ji,jj) = ssha(ji,jj)                 ! now <-- after
464               END DO
465            END DO
466         ENDIF
467         !
468      ENDIF
469      !
470      ! Update velocity at AGRIF zoom boundaries
471#if defined key_agrif
472      IF ( .NOT.Agrif_Root() ) CALL Agrif_Update_Dyn( kt )
473#endif
474      !
475      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask, ovlap=1 )
476      !
477   END SUBROUTINE ssh_nxt
478
479   !!======================================================================
480END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.