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/UKMO_MERGE_2019/src/OCE/DYN – NEMO

source: NEMO/branches/2019/UKMO_MERGE_2019/src/OCE/DYN/sshwzv.F90 @ 12068

Last change on this file since 12068 was 12068, checked in by davestorkey, 4 years ago

2019/UKMO_MERGE_2019 : Merging in changes from ENHANCE-02_ISF_nemo.

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