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

source: branches/dev_r2586_dynamic_mem/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 2690

Last change on this file since 2690 was 2690, checked in by gm, 13 years ago

dynamic mem: #785 ; homogeneization of the coding style associated with dyn allocation

  • Property svn:keywords set to Id
File size: 18.5 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   !! * Substitutions
49#  include "domzgr_substitute.h90"
50#  include "vectopt_loop_substitute.h90"
51   !!----------------------------------------------------------------------
52   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
53   !! $Id$
54   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
55   !!----------------------------------------------------------------------
56CONTAINS
57
58   SUBROUTINE ssh_wzv( kt ) 
59      !!----------------------------------------------------------------------
60      !!                ***  ROUTINE ssh_wzv  ***
61      !!                   
62      !! ** Purpose :   compute the after ssh (ssha), the now vertical velocity
63      !!              and update the now vertical coordinate (lk_vvl=T).
64      !!
65      !! ** Method  : - Using the incompressibility hypothesis, the vertical
66      !!      velocity is computed by integrating the horizontal divergence 
67      !!      from the bottom to the surface minus the scale factor evolution.
68      !!        The boundary conditions are w=0 at the bottom (no flux) and.
69      !!
70      !! ** action  :   ssha    : after sea surface height
71      !!                wn      : now vertical velocity
72      !!                sshu_a, sshv_a, sshf_a  : after sea surface height (lk_vvl=T)
73      !!                hu, hv, hur, hvr        : ocean depth and its inverse at u-,v-points
74      !!
75      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
76      !!----------------------------------------------------------------------
77      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
78      USE oce     , ONLY:   z3d   => ta                           ! ta used as 3D workspace
79      USE wrk_nemo, ONLY:   zhdiv => wrk_2d_1 , z2d => wrk_2d_2   ! 2D workspace
80      !
81      INTEGER, INTENT(in) ::   kt   ! time step
82      !
83      INTEGER  ::   ji, jj, jk   ! dummy loop indices
84      REAL(wp) ::   zcoefu, zcoefv, zcoeff, z2dt, z1_2dt, z1_rau0   ! local scalars
85      !!----------------------------------------------------------------------
86
87      IF( wrk_in_use(2, 1,2) ) THEN
88         CALL ctl_stop('ssh_wzv: requested workspace arrays unavailable')   ;   RETURN
89      ENDIF
90
91      IF( kt == nit000 ) THEN
92         !
93         IF(lwp) WRITE(numout,*)
94         IF(lwp) WRITE(numout,*) 'ssh_wzv : after sea surface height and now vertical velocity '
95         IF(lwp) WRITE(numout,*) '~~~~~~~ '
96         !
97         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
98         !
99         IF( lk_vvl ) THEN                    ! before and now Sea SSH at u-, v-, f-points (vvl case only)
100            DO jj = 1, jpjm1
101               DO ji = 1, jpim1                    ! caution: use of Vector Opt. not possible
102                  zcoefu = 0.5  * umask(ji,jj,1) / ( e1u(ji,jj) * e2u(ji,jj) )
103                  zcoefv = 0.5  * vmask(ji,jj,1) / ( e1v(ji,jj) * e2v(ji,jj) )
104                  zcoeff = 0.25 * umask(ji,jj,1) * umask(ji,jj+1,1)
105                  sshu_b(ji,jj) = zcoefu * ( e1t(ji  ,jj) * e2t(ji  ,jj) * sshb(ji  ,jj)     &
106                     &                     + e1t(ji+1,jj) * e2t(ji+1,jj) * sshb(ji+1,jj) )
107                  sshv_b(ji,jj) = zcoefv * ( e1t(ji,jj  ) * e2t(ji,jj  ) * sshb(ji,jj  )     &
108                     &                     + e1t(ji,jj+1) * e2t(ji,jj+1) * sshb(ji,jj+1) )
109                  sshu_n(ji,jj) = zcoefu * ( e1t(ji  ,jj) * e2t(ji  ,jj) * sshn(ji  ,jj)     &
110                     &                     + e1t(ji+1,jj) * e2t(ji+1,jj) * sshn(ji+1,jj) )
111                  sshv_n(ji,jj) = zcoefv * ( e1t(ji,jj  ) * e2t(ji,jj  ) * sshn(ji,jj  )     &
112                     &                     + e1t(ji,jj+1) * e2t(ji,jj+1) * sshn(ji,jj+1) )
113               END DO
114            END DO
115            CALL lbc_lnk( sshu_b, 'U', 1. )   ;   CALL lbc_lnk( sshu_n, 'U', 1. )
116            CALL lbc_lnk( sshv_b, 'V', 1. )   ;   CALL lbc_lnk( sshv_n, 'V', 1. )
117            DO jj = 1, jpjm1
118               DO ji = 1, jpim1      ! NO Vector Opt.
119                  sshf_n(ji,jj) = 0.5  * umask(ji,jj,1) * umask(ji,jj+1,1)                   &
120                       &               / ( e1f(ji,jj  ) * e2f(ji,jj  ) )                     &
121                       &               * ( e1u(ji,jj  ) * e2u(ji,jj  ) * sshu_n(ji,jj  )     &
122                       &                 + e1u(ji,jj+1) * e2u(ji,jj+1) * sshu_n(ji,jj+1) )
123               END DO
124            END DO
125            CALL lbc_lnk( sshf_n, 'F', 1. )
126         ENDIF
127         !
128      ENDIF
129
130      !                                           !------------------------------------------!
131      IF( lk_vvl ) THEN                           !  Regridding: Update Now Vertical coord.  !   (only in vvl case)
132         !                                        !------------------------------------------!
133         DO jk = 1, jpkm1
134            fsdept(:,:,jk) = fsdept_n(:,:,jk)         ! now local depths stored in fsdep. arrays
135            fsdepw(:,:,jk) = fsdepw_n(:,:,jk)
136            fsde3w(:,:,jk) = fsde3w_n(:,:,jk)
137            !
138            fse3t (:,:,jk) = fse3t_n (:,:,jk)         ! vertical scale factors stored in fse3. arrays
139            fse3u (:,:,jk) = fse3u_n (:,:,jk)
140            fse3v (:,:,jk) = fse3v_n (:,:,jk)
141            fse3f (:,:,jk) = fse3f_n (:,:,jk)
142            fse3w (:,:,jk) = fse3w_n (:,:,jk)
143            fse3uw(:,:,jk) = fse3uw_n(:,:,jk)
144            fse3vw(:,:,jk) = fse3vw_n(:,:,jk)
145         END DO
146         !
147         hu(:,:) = hu_0(:,:) + sshu_n(:,:)            ! now ocean depth (at u- and v-points)
148         hv(:,:) = hv_0(:,:) + sshv_n(:,:)
149         !                                            ! now masked inverse of the ocean depth (at u- and v-points)
150         hur(:,:) = umask(:,:,1) / ( hu(:,:) + 1._wp - umask(:,:,1) )
151         hvr(:,:) = vmask(:,:,1) / ( hv(:,:) + 1._wp - vmask(:,:,1) )
152         !
153      ENDIF
154      !
155      CALL div_cur( kt )                              ! Horizontal divergence & Relative vorticity
156      !
157      z2dt = 2._wp * rdt                              ! set time step size (Euler/Leapfrog)
158      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
159
160      !                                           !------------------------------!
161      !                                           !   After Sea Surface Height   !
162      !                                           !------------------------------!
163      zhdiv(:,:) = 0._wp
164      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
165        zhdiv(:,:) = zhdiv(:,:) + fse3t(:,:,jk) * hdivn(:,:,jk)
166      END DO
167      !                                                ! Sea surface elevation time stepping
168      ! In forward Euler time stepping case, the same formulation as in the leap-frog case can be used
169      ! because emp_b field is initialized with the vlaues of emp field. Hence, 0.5 * ( emp + emp_b ) = emp
170      z1_rau0 = 0.5 / rau0
171      ssha(:,:) = (  sshb(:,:) - z2dt * ( z1_rau0 * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * tmask(:,:,1)
172
173#if defined key_agrif
174      CALL agrif_ssh( kt )
175#endif
176#if defined key_obc
177      IF( Agrif_Root() ) THEN
178         ssha(:,:) = ssha(:,:) * obctmsk(:,:)
179         CALL lbc_lnk( ssha, 'T', 1. )                 ! absolutly compulsory !! (jmm)
180      ENDIF
181#endif
182#if defined key_bdy
183      ssha(:,:) = ssha(:,:) * bdytmask(:,:)
184      CALL lbc_lnk( ssha, 'T', 1. ) 
185#endif
186
187      !                                                ! Sea Surface Height at u-,v- and f-points (vvl case only)
188      IF( lk_vvl ) THEN                                ! (required only in key_vvl case)
189         DO jj = 1, jpjm1
190            DO ji = 1, jpim1      ! NO Vector Opt.
191               sshu_a(ji,jj) = 0.5  * umask(ji,jj,1) / ( e1u(ji  ,jj) * e2u(ji  ,jj) )                   &
192                  &                                  * ( e1t(ji  ,jj) * e2t(ji  ,jj) * ssha(ji  ,jj)     &
193                  &                                    + e1t(ji+1,jj) * e2t(ji+1,jj) * ssha(ji+1,jj) )
194               sshv_a(ji,jj) = 0.5  * vmask(ji,jj,1) / ( e1v(ji,jj  ) * e2v(ji,jj  ) )                   &
195                  &                                  * ( e1t(ji,jj  ) * e2t(ji,jj  ) * ssha(ji,jj  )     &
196                  &                                    + e1t(ji,jj+1) * e2t(ji,jj+1) * ssha(ji,jj+1) )
197            END DO
198         END DO
199         CALL lbc_lnk( sshu_a, 'U', 1. )   ;   CALL lbc_lnk( sshv_a, 'V', 1. )      ! Boundaries conditions
200      ENDIF
201     
202#if defined key_asminc
203      !                                                ! Include the IAU weighted SSH increment
204      IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN
205         CALL ssh_asm_inc( kt )
206         ssha(:,:) = ssha(:,:) + z2dt * ssh_iau(:,:)
207      ENDIF
208#endif
209
210      !                                           !------------------------------!
211      !                                           !     Now Vertical Velocity    !
212      !                                           !------------------------------!
213      z1_2dt = 1.e0 / z2dt
214      DO jk = jpkm1, 1, -1                             ! integrate from the bottom the hor. divergence
215         ! - ML - need 3 lines here because replacement of fse3t by its expression yields too long lines otherwise
216         wn(:,:,jk) = wn(:,:,jk+1) -   fse3t_n(:,:,jk) * hdivn(:,:,jk)        &
217            &                      - ( fse3t_a(:,:,jk) - fse3t_b(:,:,jk) )    &
218            &                         * tmask(:,:,jk) * z1_2dt
219#if defined key_bdy
220         wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
221#endif
222      END DO
223
224      !                                           !------------------------------!
225      !                                           !           outputs            !
226      !                                           !------------------------------!
227      CALL iom_put( "woce", wn                    )   ! vertical velocity
228      CALL iom_put( "ssh" , sshn                  )   ! sea surface height
229      CALL iom_put( "ssh2", sshn(:,:) * sshn(:,:) )   ! square of sea surface height
230      IF( lk_diaar5 ) THEN                            ! vertical mass transport & its square value
231         ! Caution: in the VVL case, it only correponds to the baroclinic mass transport.
232         z2d(:,:) = rau0 * e1t(:,:) * e2t(:,:)
233         DO jk = 1, jpk
234            z3d(:,:,jk) = wn(:,:,jk) * z2d(:,:)
235         END DO
236         CALL iom_put( "w_masstr" , z3d                     ) 
237         CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) )
238      ENDIF
239      !
240      IF(ln_ctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask, ovlap=1 )
241      !
242      IF( wrk_not_released(2, 1,2) )   CALL ctl_stop('ssh_wzv: failed to release workspace arrays')
243      !
244   END SUBROUTINE ssh_wzv
245
246
247   SUBROUTINE ssh_nxt( kt )
248      !!----------------------------------------------------------------------
249      !!                    ***  ROUTINE ssh_nxt  ***
250      !!
251      !! ** Purpose :   achieve the sea surface  height time stepping by
252      !!              applying Asselin time filter and swapping the arrays
253      !!              ssha  already computed in ssh_wzv 
254      !!
255      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
256      !!              from the filter, see Leclair and Madec 2010) and swap :
257      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
258      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
259      !!                sshn = ssha
260      !!
261      !! ** action  : - sshb, sshn   : before & now sea surface height
262      !!                               ready for the next time step
263      !!
264      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
265      !!----------------------------------------------------------------------
266      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
267      !!
268      INTEGER  ::   ji, jj   ! dummy loop indices
269      REAL(wp) ::   zec      ! temporary scalar
270      !!----------------------------------------------------------------------
271
272      IF( kt == nit000 ) THEN
273         IF(lwp) WRITE(numout,*)
274         IF(lwp) WRITE(numout,*) 'ssh_nxt : next sea surface height (Asselin time filter + swap)'
275         IF(lwp) WRITE(numout,*) '~~~~~~~ '
276      ENDIF
277
278      !                       !--------------------------!
279      IF( lk_vvl ) THEN       !  Variable volume levels  !     (ssh at t-, u-, v, f-points)
280         !                    !--------------------------!
281         !
282         IF( neuler == 0 .AND. kt == nit000 ) THEN    !** Euler time-stepping at first time-step : no filter
283            sshn  (:,:) = ssha  (:,:)                       ! now <-- after  (before already = now)
284            sshu_n(:,:) = sshu_a(:,:)
285            sshv_n(:,:) = sshv_a(:,:)
286            DO jj = 1, jpjm1                                ! ssh now at f-point
287               DO ji = 1, jpim1      ! NO Vector Opt.
288                  sshf_n(ji,jj) = 0.5  * umask(ji,jj,1) * umask(ji,jj+1,1)                 &
289                     &               / ( e1f(ji,jj  ) * e2f(ji,jj  ) )                     &
290                     &               * ( e1u(ji,jj  ) * e2u(ji,jj  ) * sshu_n(ji,jj  )     &
291                     &                 + e1u(ji,jj+1) * e2u(ji,jj+1) * sshu_n(ji,jj+1) )
292               END DO
293            END DO
294            CALL lbc_lnk( sshf_n, 'F', 1. )                 ! Boundaries conditions
295            !
296         ELSE                                         !** Leap-Frog time-stepping: Asselin filter + swap
297            zec = atfp * rdt / rau0
298            DO jj = 1, jpj
299               DO ji = 1, jpi                               ! before <-- now filtered
300                  sshb  (ji,jj) = sshn  (ji,jj) + atfp * ( sshb(ji,jj) - 2 * sshn(ji,jj) + ssha(ji,jj) )   &
301                     &                          - zec  * ( emp_b(ji,jj) - emp(ji,jj) ) * tmask(ji,jj,1)
302                  sshn  (ji,jj) = ssha  (ji,jj)             ! now <-- after
303                  sshu_n(ji,jj) = sshu_a(ji,jj)
304                  sshv_n(ji,jj) = sshv_a(ji,jj)
305               END DO
306            END DO
307            DO jj = 1, jpjm1                                ! ssh now at f-point
308               DO ji = 1, jpim1      ! NO Vector Opt.
309                  sshf_n(ji,jj) = 0.5  * umask(ji,jj,1) * umask(ji,jj+1,1)                 &
310                     &               / ( e1f(ji,jj  ) * e2f(ji,jj  ) )                     &
311                     &               * ( e1u(ji,jj  ) * e2u(ji,jj  ) * sshu_n(ji,jj  )     &
312                     &                 + e1u(ji,jj+1) * e2u(ji,jj+1) * sshu_n(ji,jj+1) )
313               END DO
314            END DO
315            CALL lbc_lnk( sshf_n, 'F', 1. )                 ! Boundaries conditions
316            !
317            DO jj = 1, jpjm1                                ! ssh before at u- & v-points
318               DO ji = 1, jpim1      ! NO Vector Opt.
319                  sshu_b(ji,jj) = 0.5  * umask(ji,jj,1) / ( e1u(ji  ,jj) * e2u(ji  ,jj) )                   &
320                     &                                  * ( e1t(ji  ,jj) * e2t(ji  ,jj) * sshb(ji  ,jj)     &
321                     &                                    + e1t(ji+1,jj) * e2t(ji+1,jj) * sshb(ji+1,jj) )
322                  sshv_b(ji,jj) = 0.5  * vmask(ji,jj,1) / ( e1v(ji,jj  ) * e2v(ji,jj  ) )                   &
323                     &                                  * ( e1t(ji,jj  ) * e2t(ji,jj  ) * sshb(ji,jj  )     &
324                     &                                    + e1t(ji,jj+1) * e2t(ji,jj+1) * sshb(ji,jj+1) )
325               END DO
326            END DO
327            CALL lbc_lnk( sshu_b, 'U', 1. )
328            CALL lbc_lnk( sshv_b, 'V', 1. )            !  Boundaries conditions
329            !
330         ENDIF
331         !                    !--------------------------!
332      ELSE                    !        fixed levels      !     (ssh at t-point only)
333         !                    !--------------------------!
334         !
335         IF( neuler == 0 .AND. kt == nit000 ) THEN    !** Euler time-stepping at first time-step : no filter
336            sshn(:,:) = ssha(:,:)                           ! now <-- after  (before already = now)
337            !
338         ELSE                                               ! Leap-Frog time-stepping: Asselin filter + swap
339            DO jj = 1, jpj
340               DO ji = 1, jpi                               ! before <-- now filtered
341                  sshb(ji,jj) = sshn(ji,jj) + atfp * ( sshb(ji,jj) - 2 * sshn(ji,jj) + ssha(ji,jj) )
342                  sshn(ji,jj) = ssha(ji,jj)                 ! now <-- after
343               END DO
344            END DO
345         ENDIF
346         !
347      ENDIF
348      !
349      ! Update velocity at AGRIF zoom boundaries
350#if defined key_agrif
351      IF ( .NOT.Agrif_Root() ) CALL Agrif_Update_Dyn( kt )
352#endif
353      !
354      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask, ovlap=1 )
355      !
356   END SUBROUTINE ssh_nxt
357
358   !!======================================================================
359END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.