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 NEMO/branches/2019/ENHANCE-02_ISF_nemo/src/OCE/DYN – NEMO

source: NEMO/branches/2019/ENHANCE-02_ISF_nemo/src/OCE/DYN/sshwzv.F90 @ 11395

Last change on this file since 11395 was 11395, checked in by mathiot, 5 years ago

ENHANCE-02_ISF_nemo : Initial commit isf simplification (add ISF directory, moved isf routine in and split isf cavity and isf parametrisation, ...) (ticket #2142)

  • Property svn:keywords set to Id
File size: 17.3 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   !!            3.3  !  2011-10  (M. Leclair) split former ssh_wzv routine and remove all vvl related work
11   !!----------------------------------------------------------------------
12
13   !!----------------------------------------------------------------------
14   !!   ssh_nxt       : after ssh
15   !!   ssh_swp       : filter ans swap the ssh arrays
16   !!   wzv           : compute now vertical velocity
17   !!----------------------------------------------------------------------
18   USE oce            ! ocean dynamics and tracers variables
19   USE isf            ! ice shelf
20   USE dom_oce        ! ocean space and time domain variables
21   USE sbc_oce        ! surface boundary condition: ocean
22   USE domvvl         ! Variable volume
23   USE divhor         ! horizontal divergence
24   USE phycst         ! physical constants
25   USE bdy_oce , ONLY : ln_bdy, bdytmask   ! Open BounDarY
26   USE bdydyn2d       ! bdy_ssh routine
27#if defined key_agrif
28   USE agrif_oce_interp
29#endif
30   !
31   USE iom 
32   USE in_out_manager ! I/O manager
33   USE restart        ! only for lrst_oce
34   USE prtctl         ! Print control
35   USE lbclnk         ! ocean lateral boundary condition (or mpp link)
36   USE lib_mpp        ! MPP library
37   USE timing         ! Timing
38   USE wet_dry        ! Wetting/Drying flux limiting
39
40   IMPLICIT NONE
41   PRIVATE
42
43   PUBLIC   ssh_nxt    ! called by step.F90
44   PUBLIC   wzv        ! called by step.F90
45   PUBLIC   wAimp      ! called by step.F90
46   PUBLIC   ssh_swp    ! called by step.F90
47
48   !! * Substitutions
49#  include "vectopt_loop_substitute.h90"
50   !!----------------------------------------------------------------------
51   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
52   !! $Id$
53   !! Software governed by the CeCILL license (see ./LICENSE)
54   !!----------------------------------------------------------------------
55CONTAINS
56
57   SUBROUTINE ssh_nxt( kt )
58      !!----------------------------------------------------------------------
59      !!                ***  ROUTINE ssh_nxt  ***
60      !!                   
61      !! ** Purpose :   compute the after ssh (ssha)
62      !!
63      !! ** Method  : - Using the incompressibility hypothesis, the ssh increment
64      !!      is computed by integrating the horizontal divergence and multiply by
65      !!      by the time step.
66      !!
67      !! ** action  :   ssha, after sea surface height
68      !!
69      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
70      !!----------------------------------------------------------------------
71      INTEGER, INTENT(in) ::   kt   ! time step
72      !
73      INTEGER  ::   jk            ! dummy loop indice
74      REAL(wp) ::   z2dt, zcoef   ! local scalars
75      REAL(wp), DIMENSION(jpi,jpj) ::   zhdiv   ! 2D workspace
76      !!----------------------------------------------------------------------
77      !
78      IF( ln_timing )   CALL timing_start('ssh_nxt')
79      !
80      IF( kt == nit000 ) THEN
81         IF(lwp) WRITE(numout,*)
82         IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
83         IF(lwp) WRITE(numout,*) '~~~~~~~ '
84      ENDIF
85      !
86      z2dt = 2._wp * rdt                          ! set time step size (Euler/Leapfrog)
87      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
88      zcoef = 0.5_wp * r1_rau0
89
90      !                                           !------------------------------!
91      !                                           !   After Sea Surface Height   !
92      !                                           !------------------------------!
93      IF(ln_wd_il) THEN
94         CALL wad_lmt(sshb, zcoef * (emp_b(:,:) + emp(:,:)), z2dt)
95      ENDIF
96
97      CALL div_hor( kt )                               ! Horizontal divergence
98      !
99      zhdiv(:,:) = 0._wp
100      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
101        zhdiv(:,:) = zhdiv(:,:) + e3t_n(:,:,jk) * hdivn(:,:,jk)
102      END DO
103      !                                                ! Sea surface elevation time stepping
104      ! In time-split case we need a first guess of the ssh after (using the baroclinic timestep) in order to
105      ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations.
106      !
107      ssha(:,:) = (  sshb(:,:) - z2dt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:)
108      !
109#if defined key_agrif
110      CALL agrif_ssh( kt )
111#endif
112      !
113      IF ( .NOT.ln_dynspg_ts ) THEN
114         IF( ln_bdy ) THEN
115            CALL lbc_lnk( 'sshwzv', ssha, 'T', 1. )    ! Not sure that's necessary
116            CALL bdy_ssh( ssha )             ! Duplicate sea level across open boundaries
117         ENDIF
118      ENDIF
119      !                                           !------------------------------!
120      !                                           !           outputs            !
121      !                                           !------------------------------!
122      !
123      IF(ln_ctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask )
124      !
125      IF( ln_timing )   CALL timing_stop('ssh_nxt')
126      !
127   END SUBROUTINE ssh_nxt
128
129   
130   SUBROUTINE wzv( kt )
131      !!----------------------------------------------------------------------
132      !!                ***  ROUTINE wzv  ***
133      !!                   
134      !! ** Purpose :   compute the now vertical velocity
135      !!
136      !! ** Method  : - Using the incompressibility hypothesis, the vertical
137      !!      velocity is computed by integrating the horizontal divergence 
138      !!      from the bottom to the surface minus the scale factor evolution.
139      !!        The boundary conditions are w=0 at the bottom (no flux) and.
140      !!
141      !! ** action  :   wn      : now vertical velocity
142      !!
143      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
144      !!----------------------------------------------------------------------
145      INTEGER, INTENT(in) ::   kt   ! time step
146      !
147      INTEGER  ::   ji, jj, jk   ! dummy loop indices
148      REAL(wp) ::   z1_2dt       ! local scalars
149      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   zhdiv
150      !!----------------------------------------------------------------------
151      !
152      IF( ln_timing )   CALL timing_start('wzv')
153      !
154      IF( kt == nit000 ) THEN
155         IF(lwp) WRITE(numout,*)
156         IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
157         IF(lwp) WRITE(numout,*) '~~~~~ '
158         !
159         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
160      ENDIF
161      !                                           !------------------------------!
162      !                                           !     Now Vertical Velocity    !
163      !                                           !------------------------------!
164      z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
165      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt
166      !
167      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! z_tilde and layer cases
168         ALLOCATE( zhdiv(jpi,jpj,jpk) ) 
169         !
170         DO jk = 1, jpkm1
171            ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t)
172            ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way)
173            DO jj = 2, jpjm1
174               DO ji = fs_2, fs_jpim1   ! vector opt.
175                  zhdiv(ji,jj,jk) = r1_e1e2t(ji,jj) * ( un_td(ji,jj,jk) - un_td(ji-1,jj,jk) + vn_td(ji,jj,jk) - vn_td(ji,jj-1,jk) )
176               END DO
177            END DO
178         END DO
179         CALL lbc_lnk('sshwzv', zhdiv, 'T', 1.)  ! - ML - Perhaps not necessary: not used for horizontal "connexions"
180         !                             ! Is it problematic to have a wrong vertical velocity in boundary cells?
181         !                             ! Same question holds for hdivn. Perhaps just for security
182         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
183            ! computation of w
184            wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk) + zhdiv(:,:,jk)    &
185               &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )     ) * tmask(:,:,jk)
186         END DO
187         !          IF( ln_vvl_layer ) wn(:,:,:) = 0.e0
188         DEALLOCATE( zhdiv ) 
189      ELSE   ! z_star and linear free surface cases
190         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
191            ! computation of w
192            wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk)                 &
193               &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )  ) * tmask(:,:,jk)
194         END DO
195      ENDIF
196
197      IF( ln_bdy ) THEN
198         DO jk = 1, jpkm1
199            wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
200         END DO
201      ENDIF
202      !
203#if defined key_agrif 
204      IF( .NOT. AGRIF_Root() ) THEN
205         IF ((nbondi ==  1).OR.(nbondi == 2)) wn(nlci-1 , :     ,:) = 0.e0      ! east
206         IF ((nbondi == -1).OR.(nbondi == 2)) wn(2      , :     ,:) = 0.e0      ! west
207         IF ((nbondj ==  1).OR.(nbondj == 2)) wn(:      ,nlcj-1 ,:) = 0.e0      ! north
208         IF ((nbondj == -1).OR.(nbondj == 2)) wn(:      ,2      ,:) = 0.e0      ! south
209      ENDIF 
210#endif 
211      !
212      IF( ln_timing )   CALL timing_stop('wzv')
213      !
214   END SUBROUTINE wzv
215
216
217   SUBROUTINE ssh_swp( kt )
218      !!----------------------------------------------------------------------
219      !!                    ***  ROUTINE ssh_nxt  ***
220      !!
221      !! ** Purpose :   achieve the sea surface  height time stepping by
222      !!              applying Asselin time filter and swapping the arrays
223      !!              ssha  already computed in ssh_nxt 
224      !!
225      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
226      !!              from the filter, see Leclair and Madec 2010) and swap :
227      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
228      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
229      !!                sshn = ssha
230      !!
231      !! ** action  : - sshb, sshn   : before & now sea surface height
232      !!                               ready for the next time step
233      !!
234      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
235      !!----------------------------------------------------------------------
236      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
237      !
238      REAL(wp) ::   zcoef   ! local scalar
239      !!----------------------------------------------------------------------
240      !
241      IF( ln_timing )   CALL timing_start('ssh_swp')
242      !
243      IF( kt == nit000 ) THEN
244         IF(lwp) WRITE(numout,*)
245         IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
246         IF(lwp) WRITE(numout,*) '~~~~~~~ '
247      ENDIF
248      !              !==  Euler time-stepping: no filter, just swap  ==!
249      IF ( neuler == 0 .AND. kt == nit000 ) THEN
250         sshn(:,:) = ssha(:,:)                              ! now    <-- after  (before already = now)
251         !
252      ELSE           !==  Leap-Frog time-stepping: Asselin filter + swap  ==!
253         !                                                  ! before <-- now filtered
254         sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) )
255         IF( .NOT.ln_linssh ) THEN                          ! before <-- with forcing removed
256            zcoef = atfp * rdt * r1_rau0
257            sshb(:,:) = sshb(:,:) - zcoef * (     emp_b(:,:) - emp   (:,:)   &
258               &                             -    rnf_b(:,:) + rnf   (:,:)   &
259               &                             + fwfisf_cav_b(:,:) - fwfisf_cav(:,:)   &
260               &                             + fwfisf_par_b(:,:) - fwfisf_par(:,:)   ) * ssmask(:,:)
261         ENDIF
262         sshn(:,:) = ssha(:,:)                              ! now <-- after
263      ENDIF
264      !
265      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask )
266      !
267      IF( ln_timing )   CALL timing_stop('ssh_swp')
268      !
269   END SUBROUTINE ssh_swp
270
271   SUBROUTINE wAimp( kt )
272      !!----------------------------------------------------------------------
273      !!                ***  ROUTINE wAimp  ***
274      !!                   
275      !! ** Purpose :   compute the Courant number and partition vertical velocity
276      !!                if a proportion needs to be treated implicitly
277      !!
278      !! ** Method  : -
279      !!
280      !! ** action  :   wn      : now vertical velocity (to be handled explicitly)
281      !!            :   wi      : now vertical velocity (for implicit treatment)
282      !!
283      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
284      !!----------------------------------------------------------------------
285      INTEGER, INTENT(in) ::   kt   ! time step
286      !
287      INTEGER  ::   ji, jj, jk   ! dummy loop indices
288      REAL(wp)             ::   zCu, zcff, z1_e3t                     ! local scalars
289      REAL(wp) , PARAMETER ::   Cu_min = 0.15_wp                      ! local parameters
290      REAL(wp) , PARAMETER ::   Cu_max = 0.27                         ! local parameters
291      REAL(wp) , PARAMETER ::   Cu_cut = 2._wp*Cu_max - Cu_min        ! local parameters
292      REAL(wp) , PARAMETER ::   Fcu    = 4._wp*Cu_max*(Cu_max-Cu_min) ! local parameters
293      !!----------------------------------------------------------------------
294      !
295      IF( ln_timing )   CALL timing_start('wAimp')
296      !
297      IF( kt == nit000 ) THEN
298         IF(lwp) WRITE(numout,*)
299         IF(lwp) WRITE(numout,*) 'wAimp : Courant number-based partitioning of now vertical velocity '
300         IF(lwp) WRITE(numout,*) '~~~~~ '
301         wi(:,:,:) = 0._wp
302      ENDIF
303      !
304      !
305      DO jk = 1, jpkm1            ! calculate Courant numbers
306         DO jj = 2, jpjm1
307            DO ji = 2, fs_jpim1   ! vector opt.
308               z1_e3t = 1._wp / e3t_n(ji,jj,jk)
309               Cu_adv(ji,jj,jk) = 2._wp * rdt * ( ( MAX( wn(ji,jj,jk) , 0._wp ) - MIN( wn(ji,jj,jk+1) , 0._wp ) )   &  ! 2*rdt and not r2dt (for restartability)
310                  &                             + ( MAX( e2u(ji  ,jj)*e3u_n(ji  ,jj,jk)*un(ji  ,jj,jk), 0._wp ) -   &
311                  &                                 MIN( e2u(ji-1,jj)*e3u_n(ji-1,jj,jk)*un(ji-1,jj,jk), 0._wp ) )   &
312                  &                               * r1_e1e2t(ji,jj)                                                 &
313                  &                             + ( MAX( e1v(ji,jj  )*e3v_n(ji,jj  ,jk)*vn(ji,jj  ,jk), 0._wp ) -   &
314                  &                                 MIN( e1v(ji,jj-1)*e3v_n(ji,jj-1,jk)*vn(ji,jj-1,jk), 0._wp ) )   &
315                  &                               * r1_e1e2t(ji,jj)                                                 &
316                  &                             ) * z1_e3t
317            END DO
318         END DO
319      END DO
320      CALL lbc_lnk( 'sshwzv', Cu_adv, 'T', 1. )
321      !
322      CALL iom_put("Courant",Cu_adv)
323      !
324      IF( MAXVAL( Cu_adv(:,:,:) ) > Cu_min ) THEN       ! Quick check if any breaches anywhere
325         DO jk = jpkm1, 2, -1                           ! or scan Courant criterion and partition
326            DO jj = 1, jpj                              ! w where necessary
327               DO ji = 1, jpi
328                  !
329                  zCu = MAX( Cu_adv(ji,jj,jk) , Cu_adv(ji,jj,jk-1) )
330! alt:
331!                  IF ( wn(ji,jj,jk) > 0._wp ) THEN
332!                     zCu =  Cu_adv(ji,jj,jk)
333!                  ELSE
334!                     zCu =  Cu_adv(ji,jj,jk-1)
335!                  ENDIF
336                  !
337                  IF( zCu <= Cu_min ) THEN              !<-- Fully explicit
338                     zcff = 0._wp
339                  ELSEIF( zCu < Cu_cut ) THEN           !<-- Mixed explicit
340                     zcff = ( zCu - Cu_min )**2
341                     zcff = zcff / ( Fcu + zcff )
342                  ELSE                                  !<-- Mostly implicit
343                     zcff = ( zCu - Cu_max )/ zCu
344                  ENDIF
345                  zcff = MIN(1._wp, zcff)
346                  !
347                  wi(ji,jj,jk) =           zcff   * wn(ji,jj,jk)
348                  wn(ji,jj,jk) = ( 1._wp - zcff ) * wn(ji,jj,jk)
349                  !
350                  Cu_adv(ji,jj,jk) = zcff               ! Reuse array to output coefficient
351               END DO
352            END DO
353         END DO
354         Cu_adv(:,:,1) = 0._wp 
355      ELSE
356         ! Fully explicit everywhere
357         Cu_adv(:,:,:) = 0._wp                          ! Reuse array to output coefficient
358         wi    (:,:,:) = 0._wp
359      ENDIF
360      CALL iom_put("wimp",wi) 
361      CALL iom_put("wi_cff",Cu_adv)
362      CALL iom_put("wexp",wn)
363      !
364      IF( ln_timing )   CALL timing_stop('wAimp')
365      !
366   END SUBROUTINE wAimp
367   !!======================================================================
368END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.