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.
zdfgls.F90 in branches/UKMO/AMM15_v3_6_STABLE_package_collate_coupling/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/UKMO/AMM15_v3_6_STABLE_package_collate_coupling/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfgls.F90 @ 10478

Last change on this file since 10478 was 10478, checked in by jcastill, 5 years ago

Merged branch UKMO/r6232_HZG_WAVE-coupling@9868

File size: 67.4 KB
Line 
1MODULE zdfgls
2   !!======================================================================
3   !!                       ***  MODULE  zdfgls  ***
4   !! Ocean physics:  vertical mixing coefficient computed from the gls
5   !!                 turbulent closure parameterization
6   !!======================================================================
7   !! History :   3.0  !  2009-09  (G. Reffray)  Original code
8   !!             3.3  !  2010-10  (C. Bricaud)  Add in the reference
9   !!----------------------------------------------------------------------
10#if defined key_zdfgls   ||   defined key_esopa
11   !!----------------------------------------------------------------------
12   !!   'key_zdfgls'                 Generic Length Scale vertical physics
13   !!----------------------------------------------------------------------
14   !!   zdf_gls       : update momentum and tracer Kz from a gls scheme
15   !!   zdf_gls_init  : initialization, namelist read, and parameters control
16   !!   gls_rst       : read/write gls restart in ocean restart file
17   !!----------------------------------------------------------------------
18   USE oce            ! ocean dynamics and active tracers
19   USE dom_oce        ! ocean space and time domain
20   USE domvvl         ! ocean space and time domain : variable volume layer
21   USE zdf_oce        ! ocean vertical physics
22   USE zdfbfr         ! bottom friction (only for rn_bfrz0)
23   USE sbc_oce        ! surface boundary condition: ocean
24   USE phycst         ! physical constants
25   USE zdfmxl         ! mixed layer
26   USE sbcwave
27   !
28   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
29   USE lib_mpp        ! MPP manager
30   USE wrk_nemo       ! work arrays
31   USE prtctl         ! Print control
32   USE in_out_manager ! I/O manager
33   USE iom            ! I/O manager library
34   USE timing         ! Timing
35   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
36
37   IMPLICIT NONE
38   PRIVATE
39
40   PUBLIC   zdf_gls        ! routine called in step module
41   PUBLIC   zdf_gls_init   ! routine called in opa module
42   PUBLIC   gls_rst        ! routine called in step module
43
44   LOGICAL , PUBLIC, PARAMETER ::   lk_zdfgls = .TRUE.   !: TKE vertical mixing flag
45   !
46   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   mxln    !: now mixing length
47   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   zwall   !: wall function
48   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   ustars2 !: Squared surface velocity scale at T-points
49   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   ustarb2 !: Squared bottom  velocity scale at T-points
50
51   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   :: rsbc_tke1 
52   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   :: rsbc_tke3 
53   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   :: rsbc_psi1 
54   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   :: rtrans
55
56   !                              !! ** Namelist  namzdf_gls  **
57   LOGICAL  ::   ln_length_lim     ! use limit on the dissipation rate under stable stratification (Galperin et al. 1988)
58   LOGICAL  ::   ln_sigpsi         ! Activate Burchard (2003) modification for k-eps closure & wave breaking mixing
59   INTEGER  ::   nn_bc_surf        ! surface boundary condition (=0/1)
60   INTEGER  ::   nn_bc_bot         ! bottom boundary condition (=0/1)
61   INTEGER  ::   nn_z0_met         ! Method for surface roughness computation
62   INTEGER  ::   nn_stab_func      ! stability functions G88, KC or Canuto (=0/1/2)
63   INTEGER  ::   nn_clos           ! closure 0/1/2/3 MY82/k-eps/k-w/gen
64   REAL(wp) ::   rn_clim_galp      ! Holt 2008 value for k-eps: 0.267
65   REAL(wp) ::   rn_epsmin         ! minimum value of dissipation (m2/s3)
66   REAL(wp) ::   rn_emin           ! minimum value of TKE (m2/s2)
67   REAL(wp) ::   rn_charn          ! Charnock constant for surface breaking waves mixing : 1400. (standard) or 2.e5 (Stacey value)
68   REAL(wp) ::   rn_crban_default  ! Craig and Banner constant for surface breaking waves mixing
69   REAL(wp) ::   rn_hsro           ! Minimum surface roughness
70   REAL(wp) ::   rn_frac_hs        ! Fraction of wave height as surface roughness (if nn_z0_met = 1, 3)
71
72   REAL(wp) ::   rcm_sf        =  0.73_wp     ! Shear free turbulence parameters
73   REAL(wp) ::   ra_sf         = -2.0_wp      ! Must be negative -2 < ra_sf < -1
74   REAL(wp) ::   rl_sf         =  0.2_wp      ! 0 <rl_sf<vkarmn   
75   REAL(wp) ::   rghmin        = -0.28_wp
76   REAL(wp) ::   rgh0          =  0.0329_wp
77   REAL(wp) ::   rghcri        =  0.03_wp
78   REAL(wp) ::   ra1           =  0.92_wp
79   REAL(wp) ::   ra2           =  0.74_wp
80   REAL(wp) ::   rb1           = 16.60_wp
81   REAL(wp) ::   rb2           = 10.10_wp         
82   REAL(wp) ::   re2           =  1.33_wp         
83   REAL(wp) ::   rl1           =  0.107_wp
84   REAL(wp) ::   rl2           =  0.0032_wp
85   REAL(wp) ::   rl3           =  0.0864_wp
86   REAL(wp) ::   rl4           =  0.12_wp
87   REAL(wp) ::   rl5           = 11.9_wp
88   REAL(wp) ::   rl6           =  0.4_wp
89   REAL(wp) ::   rl7           =  0.0_wp
90   REAL(wp) ::   rl8           =  0.48_wp
91   REAL(wp) ::   rm1           =  0.127_wp
92   REAL(wp) ::   rm2           =  0.00336_wp
93   REAL(wp) ::   rm3           =  0.0906_wp
94   REAL(wp) ::   rm4           =  0.101_wp
95   REAL(wp) ::   rm5           = 11.2_wp
96   REAL(wp) ::   rm6           =  0.4_wp
97   REAL(wp) ::   rm7           =  0.0_wp
98   REAL(wp) ::   rm8           =  0.318_wp
99   REAL(wp) ::   rtrans_default =  0.1_wp
100   REAL(wp) ::   rc02, rc02r, rc03, rc04                          ! coefficients deduced from above parameters
101   REAL(wp) ::   rsbc_tke1_default, rsbc_tke2, rfact_tke          !     -           -           -        -
102   REAL(wp) ::   rsbc_psi2, rsbc_psi3, rfact_psi                  !     -           -           -        -
103   REAL(wp) ::   rsbc_zs1, rsbc_zs2                               !     -           -           -        -
104   REAL(wp) ::   rc0, rc2, rc3, rf6, rcff, rc_diff                !     -           -           -        -
105   REAL(wp) ::   rs0, rs1, rs2, rs4, rs5, rs6                     !     -           -           -        -
106   REAL(wp) ::   rd0, rd1, rd2, rd3, rd4, rd5                     !     -           -           -        -
107   REAL(wp) ::   rsc_tke, rsc_psi, rpsi1, rpsi2, rpsi3, rsc_psi0  !     -           -           -        -
108   REAL(wp) ::   rpsi3m, rpsi3p, rpp, rmm, rnn                    !     -           -           -        -
109
110   !! * Substitutions
111#  include "domzgr_substitute.h90"
112#  include "vectopt_loop_substitute.h90"
113   !!----------------------------------------------------------------------
114   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
115   !! $Id$
116   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
117   !!----------------------------------------------------------------------
118CONTAINS
119
120   INTEGER FUNCTION zdf_gls_alloc()
121      !!----------------------------------------------------------------------
122      !!                ***  FUNCTION zdf_gls_alloc  ***
123      !!----------------------------------------------------------------------
124      ALLOCATE( mxln(jpi,jpj,jpk) , zwall(jpi,jpj,jpk) ,     & 
125         &      ustars2(jpi,jpj)  , ustarb2(jpi,jpj)   ,     & 
126         &      rsbc_tke1(jpi,jpj), rsbc_tke3(jpi,jpj) ,     & 
127         &      rsbc_psi1(jpi,jpj), rtrans(jpi,jpj), STAT= zdf_gls_alloc )
128         !
129      IF( lk_mpp             )   CALL mpp_sum ( zdf_gls_alloc )
130      IF( zdf_gls_alloc /= 0 )   CALL ctl_warn('zdf_gls_alloc: failed to allocate arrays')
131   END FUNCTION zdf_gls_alloc
132
133
134   SUBROUTINE zdf_gls( kt )
135      !!----------------------------------------------------------------------
136      !!                   ***  ROUTINE zdf_gls  ***
137      !!
138      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
139      !!              coefficients using the GLS turbulent closure scheme.
140      !!----------------------------------------------------------------------
141      INTEGER, INTENT(in) ::   kt ! ocean time step
142      INTEGER  ::   ji, jj, jk, ibot, ibotm1, dir  ! dummy loop arguments
143      REAL(wp) ::   zesh2, zsigpsi, zcoef, zex1, zex2   ! local scalars
144      REAL(wp) ::   ztx2, zty2, zup, zdown, zcof        !   -      -
145      REAL(wp) ::   zratio, zrn2, zflxb, sh             !   -      -
146      REAL(wp) ::   prod, buoy, diss, zdiss, sm         !   -      -
147      REAL(wp) ::   gh, gm, shr, dif, zsqen, zav        !   -      -
148      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zdep
149      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zkar
150      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zflxs       ! Turbulence fluxed induced by internal waves
151      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zhsro       ! Surface roughness (surface waves)
152      REAL(wp), POINTER, DIMENSION(:,:,:) ::   eb          ! tke at time before
153      REAL(wp), POINTER, DIMENSION(:,:,:) ::   mxlb        ! mixing length at time before
154      REAL(wp), POINTER, DIMENSION(:,:,:) ::   shear       ! vertical shear
155      REAL(wp), POINTER, DIMENSION(:,:,:) ::   eps         ! dissipation rate
156      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zwall_psi   ! Wall function use in the wb case (ln_sigpsi)
157      REAL(wp), POINTER, DIMENSION(:,:,:) ::   psi         ! psi at time now
158      REAL(wp), POINTER, DIMENSION(:,:,:) ::   z_elem_a    ! element of the first  matrix diagonal
159      REAL(wp), POINTER, DIMENSION(:,:,:) ::   z_elem_b    ! element of the second matrix diagonal
160      REAL(wp), POINTER, DIMENSION(:,:,:) ::   z_elem_c    ! element of the third  matrix diagonal
161      !!--------------------------------------------------------------------
162      !
163      IF( nn_timing == 1 )  CALL timing_start('zdf_gls')
164      !
165      CALL wrk_alloc( jpi,jpj, zdep, zkar, zflxs, zhsro )
166      CALL wrk_alloc( jpi,jpj,jpk, eb, mxlb, shear, eps, zwall_psi, z_elem_a, z_elem_b, z_elem_c, psi  )
167     
168      ! Preliminary computing
169
170      ustars2 = 0._wp   ;   ustarb2 = 0._wp   ;   psi  = 0._wp   ;   zwall_psi = 0._wp
171
172      ! variable initialization
173      IF( ln_phioc ) THEN
174         IF( nn_wmix==jp_janssen ) THEN
175            rsbc_tke1(:,:) = (-rsc_tke*rn_crban(:,:)/(rcm_sf*ra_sf*rl_sf))**(2._wp/3._wp)  ! k_eps = 53.Dirichlet + Wave breaking 
176            rsbc_tke3(:,:) = rdt * rn_crban(:,:)                                           ! Neumann + Wave breaking
177            rsbc_psi1(:,:) = rc0**rpp * rsbc_tke1(:,:)**rmm * rl_sf**rnn                   ! Dirichlet + Wave breaking
178         ELSE IF( nn_wmix==jp_craigbanner ) THEN
179            rsbc_tke1(:,:) = -3._wp/2._wp*rn_crban(:,:)*ra_sf*rl_sf 
180            rsbc_tke3(:,:) = rdt * rn_crban(:,:) / rl_sf 
181            rtrans(:,:) = 0.2_wp / MAX( rsmall, rsbc_tke1(:,:)**(1./(-ra_sf*3._wp/2._wp))-1._wp ) 
182         ENDIF
183      ENDIF
184
185      IF( kt /= nit000 ) THEN   ! restore before value to compute tke
186         avt (:,:,:) = avt_k (:,:,:)
187         avm (:,:,:) = avm_k (:,:,:)
188         avmu(:,:,:) = avmu_k(:,:,:)
189         avmv(:,:,:) = avmv_k(:,:,:) 
190      ENDIF
191
192      ! Compute surface and bottom friction at T-points
193!CDIR NOVERRCHK         
194      DO jj = 2, jpjm1         
195!CDIR NOVERRCHK         
196         DO ji = fs_2, fs_jpim1   ! vector opt.         
197            !
198            ! surface friction
199            ustars2(ji,jj) = r1_rau0 * taum(ji,jj) * tmask(ji,jj,1)
200            !   
201            ! bottom friction (explicit before friction)       
202            ! Note that we chose here not to bound the friction as in dynbfr)   
203            ztx2 = (  bfrua(ji,jj)  * ub(ji,jj,mbku(ji,jj)) + bfrua(ji-1,jj) * ub(ji-1,jj,mbku(ji-1,jj))  )   &         
204               & * ( 1._wp - 0.5_wp * umask(ji,jj,1) * umask(ji-1,jj,1)  )     
205            zty2 = (  bfrva(ji,jj)  * vb(ji,jj,mbkv(ji,jj)) + bfrva(ji,jj-1) * vb(ji,jj-1,mbkv(ji,jj-1))  )   &         
206               & * ( 1._wp - 0.5_wp * vmask(ji,jj,1) * vmask(ji,jj-1,1)  )     
207            ustarb2(ji,jj) = SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1)         
208         END DO         
209      END DO   
210
211      ! Set surface roughness length
212      SELECT CASE ( nn_z0_met )
213      !
214      CASE ( 0 )             ! Constant roughness         
215         zhsro(:,:) = rn_hsro
216      CASE ( 1 )             ! Standard Charnock formula
217         zhsro(:,:) = MAX(rsbc_zs1 * ustars2(:,:), rn_hsro)
218      CASE ( 2 )             ! Roughness formulae according to Rascle et al., Ocean Modelling (2008)
219         zdep(:,:)  = 30.*TANH(2.*0.3/(28.*SQRT(MAX(ustars2(:,:),rsmall))))             ! Wave age (eq. 10)
220         zhsro(:,:) = MAX(rsbc_zs2 * ustars2(:,:) * zdep(:,:)**1.5, rn_hsro) ! zhsro = rn_frac_hs * Hsw (eq. 11)
221      CASE ( 3 )             ! Roughness given by the wave model (coupled or read in file) 
222         zhsro = MAX(rn_frac_hs * hsw, rn_hsro)
223      END SELECT
224
225      ! Compute shear and dissipation rate
226      DO jk = 2, jpkm1
227         DO jj = 2, jpjm1
228            DO ji = fs_2, fs_jpim1   ! vector opt.
229               avmu(ji,jj,jk) = avmu(ji,jj,jk) * ( un(ji,jj,jk-1) - un(ji,jj,jk) )   &
230                  &                            * ( ub(ji,jj,jk-1) - ub(ji,jj,jk) )   &
231                  &                            / (  fse3uw_n(ji,jj,jk)               &
232                  &                            *    fse3uw_b(ji,jj,jk) )
233               avmv(ji,jj,jk) = avmv(ji,jj,jk) * ( vn(ji,jj,jk-1) - vn(ji,jj,jk) )   &
234                  &                            * ( vb(ji,jj,jk-1) - vb(ji,jj,jk) )   &
235                  &                            / (  fse3vw_n(ji,jj,jk)               &
236                  &                            *    fse3vw_b(ji,jj,jk) )
237               eps(ji,jj,jk)  = rc03 * en(ji,jj,jk) * SQRT(en(ji,jj,jk)) / mxln(ji,jj,jk)
238            END DO
239         END DO
240      END DO
241      !
242      ! Lateral boundary conditions (avmu,avmv) (sign unchanged)
243      CALL lbc_lnk( avmu, 'U', 1. )   ;   CALL lbc_lnk( avmv, 'V', 1. )
244
245      ! Save tke at before time step
246      eb  (:,:,:) = en  (:,:,:)
247      mxlb(:,:,:) = mxln(:,:,:)
248
249      IF( nn_clos == 0 ) THEN    ! Mellor-Yamada
250         DO jk = 2, jpkm1
251            DO jj = 2, jpjm1 
252               DO ji = fs_2, fs_jpim1   ! vector opt.
253                  zup   = mxln(ji,jj,jk) * fsdepw(ji,jj,mbkt(ji,jj)+1)
254                  zdown = vkarmn * fsdepw(ji,jj,jk) * ( -fsdepw(ji,jj,jk) + fsdepw(ji,jj,mbkt(ji,jj)+1) )
255                  zcoef = ( zup / MAX( zdown, rsmall ) )
256                  zwall (ji,jj,jk) = ( 1._wp + re2 * zcoef*zcoef ) * tmask(ji,jj,jk)
257               END DO
258            END DO
259         END DO
260      ENDIF
261
262      !!---------------------------------!!
263      !!   Equation to prognostic k      !!
264      !!---------------------------------!!
265      !
266      ! Now Turbulent kinetic energy (output in en)
267      ! -------------------------------
268      ! Resolution of a tridiagonal linear system by a "methode de chasse"
269      ! computation from level 2 to jpkm1  (e(1) computed after and e(jpk)=0 ).
270      ! The surface boundary condition are set after
271      ! The bottom boundary condition are also set after. In standard e(bottom)=0.
272      ! z_elem_b : diagonal z_elem_c : upper diagonal z_elem_a : lower diagonal
273      ! Warning : after this step, en : right hand side of the matrix
274
275      DO jk = 2, jpkm1
276         DO jj = 2, jpjm1
277            DO ji = fs_2, fs_jpim1   ! vector opt.
278               !
279               ! shear prod. at w-point weightened by mask
280               shear(ji,jj,jk) =  ( avmu(ji-1,jj,jk) + avmu(ji,jj,jk) ) / MAX( 1.e0 , umask(ji-1,jj,jk) + umask(ji,jj,jk) )   &
281                  &             + ( avmv(ji,jj-1,jk) + avmv(ji,jj,jk) ) / MAX( 1.e0 , vmask(ji,jj-1,jk) + vmask(ji,jj,jk) )
282               !
283               ! stratif. destruction
284               buoy = - avt(ji,jj,jk) * rn2(ji,jj,jk)
285               !
286               ! shear prod. - stratif. destruction
287               diss = eps(ji,jj,jk)
288               !
289               dir = 0.5_wp + SIGN( 0.5_wp, shear(ji,jj,jk) + buoy )   ! dir =1(=0) if shear(ji,jj,jk)+buoy >0(<0)
290               !
291               zesh2 = dir*(shear(ji,jj,jk)+buoy)+(1._wp-dir)*shear(ji,jj,jk)          ! production term
292               zdiss = dir*(diss/en(ji,jj,jk))   +(1._wp-dir)*(diss-buoy)/en(ji,jj,jk) ! dissipation term
293               !
294               ! Compute a wall function from 1. to rsc_psi*zwall/rsc_psi0
295               ! Note that as long that Dirichlet boundary conditions are NOT set at the first and last levels (GOTM style)
296               ! there is no need to set a boundary condition for zwall_psi at the top and bottom boundaries.
297               ! Otherwise, this should be rsc_psi/rsc_psi0
298               IF( ln_sigpsi ) THEN
299                  zsigpsi = MIN( 1._wp, zesh2 / eps(ji,jj,jk) )     ! 0. <= zsigpsi <= 1.
300                  zwall_psi(ji,jj,jk) = rsc_psi /   & 
301                     &     (  zsigpsi * rsc_psi + (1._wp-zsigpsi) * rsc_psi0 / MAX( zwall(ji,jj,jk), 1._wp )  )
302               ELSE
303                  zwall_psi(ji,jj,jk) = 1._wp
304               ENDIF
305               !
306               ! building the matrix
307               zcof = rfact_tke * tmask(ji,jj,jk)
308               !
309               ! lower diagonal
310               z_elem_a(ji,jj,jk) = zcof * ( avm  (ji,jj,jk  ) + avm  (ji,jj,jk-1) )   &
311                  &                      / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk  ) )
312               !
313               ! upper diagonal
314               z_elem_c(ji,jj,jk) = zcof * ( avm  (ji,jj,jk+1) + avm  (ji,jj,jk  ) )   &
315                  &                      / ( fse3t(ji,jj,jk  ) * fse3w(ji,jj,jk) )
316               !
317               ! diagonal
318               z_elem_b(ji,jj,jk) = 1._wp - z_elem_a(ji,jj,jk) - z_elem_c(ji,jj,jk)  &
319                  &                       + rdt * zdiss * tmask(ji,jj,jk) 
320               !
321               ! right hand side in en
322               en(ji,jj,jk) = en(ji,jj,jk) + rdt * zesh2 * tmask(ji,jj,jk)
323            END DO
324         END DO
325      END DO
326      !
327      z_elem_b(:,:,jpk) = 1._wp
328      !
329      ! Set surface condition on zwall_psi (1 at the bottom)
330      zwall_psi(:,:,1) = zwall_psi(:,:,2)
331      zwall_psi(:,:,jpk) = 1.
332      !
333      ! Surface boundary condition on tke
334      ! ---------------------------------
335      !
336      IF( ln_phioc ) THEN
337         SELECT CASE ( nn_bc_surf ) 
338         !
339         CASE ( 0 )             ! Dirichlet case
340            IF( nn_wmix==jp_janssen ) THEN 
341               ! First level
342               en(:,:,1) = MAX( rsbc_tke1(:,:) * ustars2(:,:), rn_emin ) 
343               z_elem_a(:,:,1) = en(:,:,1) 
344               z_elem_c(:,:,1) = 0._wp 
345               z_elem_b(:,:,1) = 1._wp 
346 
347               ! One level below
348               en(:,:,2) = MAX( rsbc_tke1(:,:) * ustars2(:,:) * ((zhsro(:,:)+fsdepw(:,:,2))/zhsro(:,:) )**ra_sf, rn_emin ) 
349               z_elem_a(:,:,2) = 0._wp 
350               z_elem_c(:,:,2) = 0._wp 
351               z_elem_b(:,:,2) = 1._wp 
352            ELSE IF( nn_wmix==jp_craigbanner ) THEN
353               en(:,:,1) = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1(:,:))**(2._wp/3._wp) 
354               en(:,:,1) = MAX(en(:,:,1), rn_emin) 
355               z_elem_a(:,:,1) = en(:,:,1) 
356               z_elem_c(:,:,1) = 0._wp 
357               z_elem_b(:,:,1) = 1._wp 
358               
359               ! One level below
360               en(:,:,2) = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1(:,:) * ((zhsro(:,:)+fsdepw(:,:,2)) & 
361                   &            / zhsro(:,:) )**(1.5_wp*ra_sf))**(2._wp/3._wp) 
362               en(:,:,2) = MAX(en(:,:,2), rn_emin ) 
363               z_elem_a(:,:,2) = 0._wp 
364               z_elem_c(:,:,2) = 0._wp 
365               z_elem_b(:,:,2) = 1._wp 
366               !
367            ENDIF
368         CASE ( 1 )             ! Neumann boundary condition on d(e)/dz
369            IF( nn_wmix==jp_janssen ) THEN 
370               ! Dirichlet conditions at k=1
371               en(:,:,1) = MAX( rsbc_tke1(:,:) * ustars2(:,:), rn_emin ) 
372               z_elem_a(:,:,1) = en(:,:,1) 
373               z_elem_c(:,:,1) = 0._wp 
374               z_elem_b(:,:,1) = 1._wp 
375               ! at k=2, set de/dz=Fw
376               z_elem_b(:,:,2) = z_elem_b(:,:,2) + z_elem_a(:,:,2) ! Remove z_elem_a from z_elem_b
377               z_elem_a(:,:,2) = 0._wp 
378               zflxs(:,:) = rsbc_tke3(:,:) * ustars2(:,:)**1.5_wp * ((zhsro(:,:)+fsdept(:,:,1) ) / zhsro(:,:) )**(1.5*ra_sf) 
379               en(:,:,2) = en(:,:,2) + zflxs(:,:) / fse3w(:,:,2) 
380            ELSE IF( nn_wmix==jp_craigbanner ) THEN 
381               ! Dirichlet conditions at k=1
382               en(:,:,1)       = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1(:,:))**(2._wp/3._wp) 
383               en(:,:,1)       = MAX(en(:,:,1), rn_emin)       
384               z_elem_a(:,:,1) = en(:,:,1) 
385               z_elem_c(:,:,1) = 0._wp 
386               z_elem_b(:,:,1) = 1._wp 
387               !
388               ! at k=2, set de/dz=Fw
389               !cbr
390               z_elem_b(:,:,2) = z_elem_b(:,:,2) + z_elem_a(:,:,2) ! Remove z_elem_a from z_elem_b
391               z_elem_a(:,:,2) = 0._wp 
392               zkar(:,:)       = (rl_sf + (vkarmn-rl_sf)*(1.-exp(-rtrans(:,:)*fsdept(:,:,1)/zhsro(:,:)) )) 
393               zflxs(:,:)      = rsbc_tke3(:,:) * ustars2(:,:)**1.5_wp * zkar(:,:) & 
394                    &                           * ((zhsro(:,:)+fsdept(:,:,1))/zhsro(:,:) )**(1.5_wp*ra_sf) 
395
396               en(:,:,2) = en(:,:,2) + zflxs(:,:)/fse3w(:,:,2) 
397               !
398            ENDIF
399         END SELECT
400      ELSE
401      SELECT CASE ( nn_bc_surf )
402      !
403      CASE ( 0 )             ! Dirichlet case
404      ! First level
405      en(:,:,1) = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1_default)**(2._wp/3._wp)
406      en(:,:,1) = MAX(en(:,:,1), rn_emin) 
407      z_elem_a(:,:,1) = en(:,:,1)
408      z_elem_c(:,:,1) = 0._wp
409      z_elem_b(:,:,1) = 1._wp
410      !
411      ! One level below
412      en(:,:,2) = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1_default * ((zhsro(:,:)+fsdepw(:,:,2)) &
413          &            / zhsro(:,:) )**(1.5_wp*ra_sf))**(2._wp/3._wp)
414      en(:,:,2) = MAX(en(:,:,2), rn_emin )
415      z_elem_a(:,:,2) = 0._wp 
416      z_elem_c(:,:,2) = 0._wp
417      z_elem_b(:,:,2) = 1._wp
418      !
419      !
420      CASE ( 1 )             ! Neumann boundary condition on d(e)/dz
421      ! Dirichlet conditions at k=1
422      en(:,:,1)       = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1_default)**(2._wp/3._wp)
423      en(:,:,1)       = MAX(en(:,:,1), rn_emin)     
424      z_elem_a(:,:,1) = en(:,:,1)
425      z_elem_c(:,:,1) = 0._wp
426      z_elem_b(:,:,1) = 1._wp
427      !
428      ! at k=2, set de/dz=Fw
429      !cbr
430      z_elem_b(:,:,2) = z_elem_b(:,:,2) +  z_elem_a(:,:,2) ! Remove z_elem_a from z_elem_b
431      z_elem_a(:,:,2) = 0._wp
432      zkar(:,:)       = (rl_sf + (vkarmn-rl_sf)*(1.-exp(-rtrans_default*fsdept(:,:,1)/zhsro(:,:)) ))
433      zflxs(:,:)      = rsbc_tke2 * ustars2(:,:)**1.5_wp * zkar(:,:) &
434           &                      * ((zhsro(:,:)+fsdept(:,:,1))/zhsro(:,:) )**(1.5_wp*ra_sf)
435
436      en(:,:,2) = en(:,:,2) + zflxs(:,:)/fse3w(:,:,2)
437      !
438      !
439      END SELECT
440      ENDIF ! ln_phioc
441
442      ! Bottom boundary condition on tke
443      ! --------------------------------
444      !
445      SELECT CASE ( nn_bc_bot )
446      !
447      CASE ( 0 )             ! Dirichlet
448         !                      ! en(ibot) = u*^2 / Co2 and mxln(ibot) = rn_lmin
449         !                      ! Balance between the production and the dissipation terms
450!CDIR NOVERRCHK
451         DO jj = 2, jpjm1
452!CDIR NOVERRCHK
453            DO ji = fs_2, fs_jpim1   ! vector opt.
454               ibot   = mbkt(ji,jj) + 1      ! k   bottom level of w-point
455               ibotm1 = mbkt(ji,jj)          ! k-1 bottom level of w-point but >=1
456               !
457               ! Bottom level Dirichlet condition:
458               z_elem_a(ji,jj,ibot  ) = 0._wp
459               z_elem_c(ji,jj,ibot  ) = 0._wp
460               z_elem_b(ji,jj,ibot  ) = 1._wp
461               en(ji,jj,ibot  ) = MAX( rc02r * ustarb2(ji,jj), rn_emin )
462               !
463               ! Just above last level, Dirichlet condition again
464               z_elem_a(ji,jj,ibotm1) = 0._wp
465               z_elem_c(ji,jj,ibotm1) = 0._wp
466               z_elem_b(ji,jj,ibotm1) = 1._wp
467               en(ji,jj,ibotm1) = MAX( rc02r * ustarb2(ji,jj), rn_emin ) 
468            END DO
469         END DO
470         !
471      CASE ( 1 )             ! Neumman boundary condition
472         !                     
473!CDIR NOVERRCHK
474         DO jj = 2, jpjm1
475!CDIR NOVERRCHK
476            DO ji = fs_2, fs_jpim1   ! vector opt.
477               ibot   = mbkt(ji,jj) + 1      ! k   bottom level of w-point
478               ibotm1 = mbkt(ji,jj)          ! k-1 bottom level of w-point but >=1
479               !
480               ! Bottom level Dirichlet condition:
481               z_elem_a(ji,jj,ibot) = 0._wp
482               z_elem_c(ji,jj,ibot) = 0._wp
483               z_elem_b(ji,jj,ibot) = 1._wp
484               en(ji,jj,ibot) = MAX( rc02r * ustarb2(ji,jj), rn_emin )
485               !
486               ! Just above last level: Neumann condition
487               z_elem_b(ji,jj,ibotm1) = z_elem_b(ji,jj,ibotm1) + z_elem_c(ji,jj,ibotm1)   ! Remove z_elem_c from z_elem_b
488               z_elem_c(ji,jj,ibotm1) = 0._wp
489            END DO
490         END DO
491         !
492      END SELECT
493
494      ! Matrix inversion (en prescribed at surface and the bottom)
495      ! ----------------------------------------------------------
496      !
497      DO jk = 2, jpkm1                             ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
498         DO jj = 2, jpjm1
499            DO ji = fs_2, fs_jpim1    ! vector opt.
500               z_elem_b(ji,jj,jk) = z_elem_b(ji,jj,jk) - z_elem_a(ji,jj,jk) * z_elem_c(ji,jj,jk-1) / z_elem_b(ji,jj,jk-1)
501            END DO
502         END DO
503      END DO
504      DO jk = 2, jpk                               ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1
505         DO jj = 2, jpjm1
506            DO ji = fs_2, fs_jpim1    ! vector opt.
507               z_elem_a(ji,jj,jk) = en(ji,jj,jk) - z_elem_a(ji,jj,jk) / z_elem_b(ji,jj,jk-1) * z_elem_a(ji,jj,jk-1)
508            END DO
509         END DO
510      END DO
511      DO jk = jpk-1, 2, -1                         ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
512         DO jj = 2, jpjm1
513            DO ji = fs_2, fs_jpim1    ! vector opt.
514               en(ji,jj,jk) = ( z_elem_a(ji,jj,jk) - z_elem_c(ji,jj,jk) * en(ji,jj,jk+1) ) / z_elem_b(ji,jj,jk)
515            END DO
516         END DO
517      END DO
518      !                                            ! set the minimum value of tke
519      en(:,:,:) = MAX( en(:,:,:), rn_emin )
520
521      !!----------------------------------------!!
522      !!   Solve prognostic equation for psi    !!
523      !!----------------------------------------!!
524
525      ! Set psi to previous time step value
526      !
527      SELECT CASE ( nn_clos )
528      !
529      CASE( 0 )               ! k-kl  (Mellor-Yamada)
530         DO jk = 2, jpkm1
531            DO jj = 2, jpjm1
532               DO ji = fs_2, fs_jpim1   ! vector opt.
533                  psi(ji,jj,jk)  = eb(ji,jj,jk) * mxlb(ji,jj,jk)
534               END DO
535            END DO
536         END DO
537         !
538      CASE( 1 )               ! k-eps
539         DO jk = 2, jpkm1
540            DO jj = 2, jpjm1
541               DO ji = fs_2, fs_jpim1   ! vector opt.
542                  psi(ji,jj,jk)  = eps(ji,jj,jk)
543               END DO
544            END DO
545         END DO
546         !
547      CASE( 2 )               ! k-w
548         DO jk = 2, jpkm1
549            DO jj = 2, jpjm1
550               DO ji = fs_2, fs_jpim1   ! vector opt.
551                  psi(ji,jj,jk)  = SQRT( eb(ji,jj,jk) ) / ( rc0 * mxlb(ji,jj,jk) )
552               END DO
553            END DO
554         END DO
555         !
556      CASE( 3 )               ! generic
557         DO jk = 2, jpkm1
558            DO jj = 2, jpjm1
559               DO ji = fs_2, fs_jpim1   ! vector opt.
560                  psi(ji,jj,jk)  = rc02 * eb(ji,jj,jk) * mxlb(ji,jj,jk)**rnn 
561               END DO
562            END DO
563         END DO
564         !
565      END SELECT
566      !
567      ! Now gls (output in psi)
568      ! -------------------------------
569      ! Resolution of a tridiagonal linear system by a "methode de chasse"
570      ! computation from level 2 to jpkm1  (e(1) already computed and e(jpk)=0 ).
571      ! z_elem_b : diagonal z_elem_c : upper diagonal z_elem_a : lower diagonal
572      ! Warning : after this step, en : right hand side of the matrix
573
574      DO jk = 2, jpkm1
575         DO jj = 2, jpjm1
576            DO ji = fs_2, fs_jpim1   ! vector opt.
577               !
578               ! psi / k
579               zratio = psi(ji,jj,jk) / eb(ji,jj,jk) 
580               !
581               ! psi3+ : stable : B=-KhN²<0 => N²>0 if rn2>0 dir = 1 (stable) otherwise dir = 0 (unstable)
582               dir = 0.5_wp + SIGN( 0.5_wp, rn2(ji,jj,jk) )
583               !
584               rpsi3 = dir * rpsi3m + ( 1._wp - dir ) * rpsi3p
585               !
586               ! shear prod. - stratif. destruction
587               prod = rpsi1 * zratio * shear(ji,jj,jk)
588               !
589               ! stratif. destruction
590               buoy = rpsi3 * zratio * (- avt(ji,jj,jk) * rn2(ji,jj,jk) )
591               !
592               ! shear prod. - stratif. destruction
593               diss = rpsi2 * zratio * zwall(ji,jj,jk) * eps(ji,jj,jk)
594               !
595               dir = 0.5_wp + SIGN( 0.5_wp, prod + buoy )   ! dir =1(=0) if shear(ji,jj,jk)+buoy >0(<0)
596               !
597               zesh2 = dir * ( prod + buoy )          + (1._wp - dir ) * prod                        ! production term
598               zdiss = dir * ( diss / psi(ji,jj,jk) ) + (1._wp - dir ) * (diss-buoy) / psi(ji,jj,jk) ! dissipation term
599               !                                                       
600               ! building the matrix
601               zcof = rfact_psi * zwall_psi(ji,jj,jk) * tmask(ji,jj,jk)
602               ! lower diagonal
603               z_elem_a(ji,jj,jk) = zcof * ( avm  (ji,jj,jk  ) + avm  (ji,jj,jk-1) )   &
604                  &                      / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk  ) )
605               ! upper diagonal
606               z_elem_c(ji,jj,jk) = zcof * ( avm  (ji,jj,jk+1) + avm  (ji,jj,jk  ) )   &
607                  &                      / ( fse3t(ji,jj,jk  ) * fse3w(ji,jj,jk) )
608               ! diagonal
609               z_elem_b(ji,jj,jk) = 1._wp - z_elem_a(ji,jj,jk) - z_elem_c(ji,jj,jk)  &
610                  &                       + rdt * zdiss * tmask(ji,jj,jk)
611               !
612               ! right hand side in psi
613               psi(ji,jj,jk) = psi(ji,jj,jk) + rdt * zesh2 * tmask(ji,jj,jk)
614            END DO
615         END DO
616      END DO
617      !
618      z_elem_b(:,:,jpk) = 1._wp
619
620      ! Surface boundary condition on psi
621      ! ---------------------------------
622      !
623      SELECT CASE ( nn_bc_surf )
624      !
625      CASE ( 0 )             ! Dirichlet boundary conditions
626         IF( ln_phioc ) THEN   ! Wave induced mixing case
627                               ! en(1) = q2(1) = 0.5 * (15.8 * Ccb)^(2/3) * u*^2
628                               ! balance between the production and the
629                               ! dissipation terms including the wave effect
630            IF( nn_wmix==jp_janssen ) THEN 
631               ! Surface value
632               zdep(:,:)       = zhsro(:,:) * rl_sf ! Cosmetic
633               psi (:,:,1)     = rc0**rpp * en(:,:,1)**rmm * zdep(:,:)**rnn * tmask(:,:,1) 
634               z_elem_a(:,:,1) = psi(:,:,1)       
635               z_elem_c(:,:,1) = 0._wp 
636               z_elem_b(:,:,1) = 1._wp 
637               !
638               ! One level below
639               zex1 = (rmm*ra_sf+rnn) 
640               zex2 = (rmm*ra_sf) 
641               zdep(:,:) = ( (zhsro(:,:) + fsdepw(:,:,2))**zex1 ) / zhsro(:,:)**zex2 
642               psi (:,:,2) = rsbc_psi1(:,:) * ustars2(:,:)**rmm * zdep(:,:) * tmask(:,:,1) 
643               z_elem_a(:,:,2) = 0._wp 
644               z_elem_c(:,:,2) = 0._wp 
645               z_elem_b(:,:,2) = 1._wp 
646               
647               !
648            ELSE IF( nn_wmix==jp_craigbanner ) THEN 
649               ! Surface value
650               zdep(:,:)       = zhsro(:,:) * rl_sf ! Cosmetic
651               psi (:,:,1)     = rc0**rpp * en(:,:,1)**rmm * zdep(:,:)**rnn * tmask(:,:,1) 
652               z_elem_a(:,:,1) = psi(:,:,1) 
653               z_elem_c(:,:,1) = 0._wp 
654               z_elem_b(:,:,1) = 1._wp 
655               !
656               ! One level below
657               zkar(:,:)       = (rl_sf + (vkarmn-rl_sf)*(1._wp-exp(-rtrans(:,:)*fsdepw(:,:,2)/zhsro(:,:) ))) 
658               zdep(:,:)       = (zhsro(:,:) + fsdepw(:,:,2)) * zkar(:,:) 
659               psi (:,:,2)     = rc0**rpp * en(:,:,2)**rmm * zdep(:,:)**rnn * tmask(:,:,1) 
660               z_elem_a(:,:,2) = 0._wp 
661               z_elem_c(:,:,2) = 0._wp 
662               z_elem_b(:,:,2) = 1._wp 
663               
664               !
665            ENDIF
666         ELSE                  ! Wave induced mixing case with default values
667            ! Surface value
668            zdep(:,:)       = zhsro(:,:) * rl_sf ! Cosmetic
669            psi (:,:,1)     = rc0**rpp * en(:,:,1)**rmm * zdep(:,:)**rnn * tmask(:,:,1) 
670            z_elem_a(:,:,1) = psi(:,:,1) 
671            z_elem_c(:,:,1) = 0._wp 
672            z_elem_b(:,:,1) = 1._wp 
673            !
674            ! One level below
675            zkar(:,:)       = (rl_sf + (vkarmn-rl_sf)*(1._wp-exp(-rtrans_default*fsdepw(:,:,2)/zhsro(:,:) ))) 
676            zdep(:,:)       = (zhsro(:,:) + fsdepw(:,:,2)) * zkar(:,:) 
677            psi (:,:,2)     = rc0**rpp * en(:,:,2)**rmm * zdep(:,:)**rnn * tmask(:,:,1) 
678            z_elem_a(:,:,2) = 0._wp 
679            z_elem_c(:,:,2) = 0._wp 
680            z_elem_b(:,:,2) = 1._wp 
681           
682            !
683         ENDIF
684      CASE ( 1 )             ! Neumann boundary condition on d(psi)/dz
685         IF( ln_phioc ) THEN  ! Wave induced mixing case with forced/coupled fields
686            IF( nn_wmix==jp_janssen ) THEN 
687               !
688               zdep(:,:) = rl_sf * zhsro(:,:) 
689               psi (:,:,1) = rc0**rpp * en(:,:,1)**rmm * zdep(:,:)**rnn * tmask(:,:,1) 
690               z_elem_a(:,:,1) = psi(:,:,1) 
691               z_elem_c(:,:,1) = 0._wp 
692               z_elem_b(:,:,1) = 1._wp 
693               !
694               ! Neumann condition at k=2
695               z_elem_b(:,:,2) = z_elem_b(:,:,2) + z_elem_a(:,:,2) ! Remove z_elem_a from z_elem_b
696               z_elem_a(:,:,2) = 0._wp 
697               !
698               ! Set psi vertical flux at the surface:
699               zdep(:,:) = (zhsro(:,:) + fsdept(:,:,1))**(rmm*ra_sf+rnn-1._wp) / zhsro(:,:)**(rmm*ra_sf) 
700               zflxs(:,:) = rsbc_psi3 * ( zwall_psi(:,:,1)*avm(:,:,1) + zwall_psi(:,:,2)*avm(:,:,2) ) & 
701                  &                   * en(:,:,1)**rmm * zdep 
702               psi(:,:,2) = psi(:,:,2) + zflxs(:,:) / fse3w(:,:,2) 
703               !
704            ELSE IF( nn_wmix==jp_craigbanner ) THEN 
705               ! Surface value: Dirichlet
706               zdep(:,:)       = zhsro(:,:) * rl_sf 
707               psi (:,:,1)     = rc0**rpp * en(:,:,1)**rmm * zdep(:,:)**rnn * tmask(:,:,1) 
708               z_elem_a(:,:,1) = psi(:,:,1) 
709               z_elem_c(:,:,1) = 0._wp 
710               z_elem_b(:,:,1) = 1._wp 
711               !
712               ! Neumann condition at k=2
713               z_elem_b(:,:,2) = z_elem_b(:,:,2) + z_elem_a(:,:,2) ! Remove z_elem_a from z_elem_b
714               z_elem_a(:,:,2) = 0._wp 
715               !
716               ! Set psi vertical flux at the surface:
717               zkar(:,:) = rl_sf + (vkarmn-rl_sf)*(1._wp-exp(-rtrans(:,:)*fsdept(:,:,1)/zhsro(:,:) )) ! Lengh scale slope
718               zdep(:,:) = ((zhsro(:,:) + fsdept(:,:,1)) / zhsro(:,:))**(rmm*ra_sf) 
719               zflxs(:,:) = (rnn + rsbc_tke1(:,:) * (rnn + rmm*ra_sf) * zdep(:,:)) * & 
720                          (1._wp + rsbc_tke1(:,:) * zdep(:,:))**(2._wp*rmm/3._wp-1_wp) 
721               zdep(:,:) =  rsbc_psi1(:,:) * (zwall_psi(:,:,1)*avm(:,:,1)+zwall_psi(:,:,2)*avm(:,:,2)) * & 
722                      & ustars2(:,:)**rmm * zkar(:,:)**rnn * (zhsro(:,:) + fsdept(:,:,1))**(rnn-1.) 
723               zflxs(:,:) = zdep(:,:) * zflxs(:,:) 
724               psi(:,:,2) = psi(:,:,2) + zflxs(:,:) / fse3w(:,:,2) 
725               !   
726            ENDIF
727         ELSE                 ! Wave induced mixing case with default values
728            ! Surface value: Dirichlet
729            zdep(:,:)       = zhsro(:,:) * rl_sf 
730            psi (:,:,1)     = rc0**rpp * en(:,:,1)**rmm * zdep(:,:)**rnn * tmask(:,:,1) 
731            z_elem_a(:,:,1) = psi(:,:,1) 
732            z_elem_c(:,:,1) = 0._wp 
733            z_elem_b(:,:,1) = 1._wp 
734            !
735            ! Neumann condition at k=2
736            z_elem_b(:,:,2) = z_elem_b(:,:,2) +  z_elem_a(:,:,2) ! Remove z_elem_a from z_elem_b
737            z_elem_a(:,:,2) = 0._wp 
738            !
739            ! Set psi vertical flux at the surface:
740            zkar(:,:) = rl_sf + (vkarmn-rl_sf)*(1._wp-exp(-rtrans_default*fsdept(:,:,1)/zhsro(:,:) )) ! Lengh scale slope
741            zdep(:,:) = ((zhsro(:,:) + fsdept(:,:,1)) / zhsro(:,:))**(rmm*ra_sf) 
742            zflxs(:,:) = (rnn + rsbc_tke1_default * (rnn + rmm*ra_sf) * zdep(:,:)) * & 
743                       (1._wp + rsbc_tke1_default * zdep(:,:))**(2._wp*rmm/3._wp-1_wp) 
744            zdep(:,:) =  rsbc_psi1(:,:) * (zwall_psi(:,:,1)*avm(:,:,1)+zwall_psi(:,:,2)*avm(:,:,2)) * & 
745                   & ustars2(:,:)**rmm * zkar(:,:)**rnn * (zhsro(:,:) + fsdept(:,:,1))**(rnn-1.) 
746            zflxs(:,:) = zdep(:,:) * zflxs(:,:) 
747            psi(:,:,2) = psi(:,:,2) + zflxs(:,:) / fse3w(:,:,2) 
748            !   
749            !
750         ENDIF
751      END SELECT
752
753      ! Bottom boundary condition on psi
754      ! --------------------------------
755      !
756      SELECT CASE ( nn_bc_bot )
757      !
758      !
759      CASE ( 0 )             ! Dirichlet
760         !                      ! en(ibot) = u*^2 / Co2 and mxln(ibot) = vkarmn * rn_bfrz0
761         !                      ! Balance between the production and the dissipation terms
762!CDIR NOVERRCHK
763         DO jj = 2, jpjm1
764!CDIR NOVERRCHK
765            DO ji = fs_2, fs_jpim1   ! vector opt.
766               ibot   = mbkt(ji,jj) + 1      ! k   bottom level of w-point
767               ibotm1 = mbkt(ji,jj)          ! k-1 bottom level of w-point but >=1
768               zdep(ji,jj) = vkarmn * rn_bfrz0
769               psi (ji,jj,ibot) = rc0**rpp * en(ji,jj,ibot)**rmm * zdep(ji,jj)**rnn
770               z_elem_a(ji,jj,ibot) = 0._wp
771               z_elem_c(ji,jj,ibot) = 0._wp
772               z_elem_b(ji,jj,ibot) = 1._wp
773               !
774               ! Just above last level, Dirichlet condition again (GOTM like)
775               zdep(ji,jj) = vkarmn * ( rn_bfrz0 + fse3t(ji,jj,ibotm1) )
776               psi (ji,jj,ibotm1) = rc0**rpp * en(ji,jj,ibot  )**rmm * zdep(ji,jj)**rnn
777               z_elem_a(ji,jj,ibotm1) = 0._wp
778               z_elem_c(ji,jj,ibotm1) = 0._wp
779               z_elem_b(ji,jj,ibotm1) = 1._wp
780            END DO
781         END DO
782         !
783      CASE ( 1 )             ! Neumman boundary condition
784         !                     
785!CDIR NOVERRCHK
786         DO jj = 2, jpjm1
787!CDIR NOVERRCHK
788            DO ji = fs_2, fs_jpim1   ! vector opt.
789               ibot   = mbkt(ji,jj) + 1      ! k   bottom level of w-point
790               ibotm1 = mbkt(ji,jj)          ! k-1 bottom level of w-point but >=1
791               !
792               ! Bottom level Dirichlet condition:
793               zdep(ji,jj) = vkarmn * rn_bfrz0
794               psi (ji,jj,ibot) = rc0**rpp * en(ji,jj,ibot)**rmm * zdep(ji,jj)**rnn
795               !
796               z_elem_a(ji,jj,ibot) = 0._wp
797               z_elem_c(ji,jj,ibot) = 0._wp
798               z_elem_b(ji,jj,ibot) = 1._wp
799               !
800               ! Just above last level: Neumann condition with flux injection
801               z_elem_b(ji,jj,ibotm1) = z_elem_b(ji,jj,ibotm1) + z_elem_c(ji,jj,ibotm1) ! Remove z_elem_c from z_elem_b
802               z_elem_c(ji,jj,ibotm1) = 0.
803               !
804               ! Set psi vertical flux at the bottom:
805               zdep(ji,jj) = rn_bfrz0 + 0.5_wp*fse3t(ji,jj,ibotm1)
806               zflxb = rsbc_psi2 * ( avm(ji,jj,ibot) + avm(ji,jj,ibotm1) )   &
807                  &  * (0.5_wp*(en(ji,jj,ibot)+en(ji,jj,ibotm1)))**rmm * zdep(ji,jj)**(rnn-1._wp)
808               psi(ji,jj,ibotm1) = psi(ji,jj,ibotm1) + zflxb / fse3w(ji,jj,ibotm1)
809            END DO
810         END DO
811         !
812      END SELECT
813
814      ! Matrix inversion
815      ! ----------------
816      !
817      DO jk = 2, jpkm1                             ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
818         DO jj = 2, jpjm1
819            DO ji = fs_2, fs_jpim1    ! vector opt.
820               z_elem_b(ji,jj,jk) = z_elem_b(ji,jj,jk) - z_elem_a(ji,jj,jk) * z_elem_c(ji,jj,jk-1) / z_elem_b(ji,jj,jk-1)
821            END DO
822         END DO
823      END DO
824      DO jk = 2, jpk                               ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1
825         DO jj = 2, jpjm1
826            DO ji = fs_2, fs_jpim1    ! vector opt.
827               z_elem_a(ji,jj,jk) = psi(ji,jj,jk) - z_elem_a(ji,jj,jk) / z_elem_b(ji,jj,jk-1) * z_elem_a(ji,jj,jk-1)
828            END DO
829         END DO
830      END DO
831      DO jk = jpk-1, 2, -1                         ! Third recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
832         DO jj = 2, jpjm1
833            DO ji = fs_2, fs_jpim1    ! vector opt.
834               psi(ji,jj,jk) = ( z_elem_a(ji,jj,jk) - z_elem_c(ji,jj,jk) * psi(ji,jj,jk+1) ) / z_elem_b(ji,jj,jk)
835            END DO
836         END DO
837      END DO
838
839      ! Set dissipation
840      !----------------
841
842      SELECT CASE ( nn_clos )
843      !
844      CASE( 0 )               ! k-kl  (Mellor-Yamada)
845         DO jk = 1, jpkm1
846            DO jj = 2, jpjm1
847               DO ji = fs_2, fs_jpim1   ! vector opt.
848                  eps(ji,jj,jk) = rc03 * en(ji,jj,jk) * en(ji,jj,jk) * SQRT( en(ji,jj,jk) ) / MAX( psi(ji,jj,jk), rn_epsmin)
849               END DO
850            END DO
851         END DO
852         !
853      CASE( 1 )               ! k-eps
854         DO jk = 1, jpkm1
855            DO jj = 2, jpjm1
856               DO ji = fs_2, fs_jpim1   ! vector opt.
857                  eps(ji,jj,jk) = psi(ji,jj,jk)
858               END DO
859            END DO
860         END DO
861         !
862      CASE( 2 )               ! k-w
863         DO jk = 1, jpkm1
864            DO jj = 2, jpjm1
865               DO ji = fs_2, fs_jpim1   ! vector opt.
866                  eps(ji,jj,jk) = rc04 * en(ji,jj,jk) * psi(ji,jj,jk) 
867               END DO
868            END DO
869         END DO
870         !
871      CASE( 3 )               ! generic
872         zcoef = rc0**( 3._wp  + rpp/rnn )
873         zex1  =      ( 1.5_wp + rmm/rnn )
874         zex2  = -1._wp / rnn
875         DO jk = 1, jpkm1
876            DO jj = 2, jpjm1
877               DO ji = fs_2, fs_jpim1   ! vector opt.
878                  eps(ji,jj,jk) = zcoef * en(ji,jj,jk)**zex1 * psi(ji,jj,jk)**zex2
879               END DO
880            END DO
881         END DO
882         !
883      END SELECT
884
885      ! Limit dissipation rate under stable stratification
886      ! --------------------------------------------------
887      DO jk = 1, jpkm1 ! Note that this set boundary conditions on mxln at the same time
888         DO jj = 2, jpjm1
889            DO ji = fs_2, fs_jpim1    ! vector opt.
890               ! limitation
891               eps(ji,jj,jk)  = MAX( eps(ji,jj,jk), rn_epsmin )
892               mxln(ji,jj,jk)  = rc03 * en(ji,jj,jk) * SQRT( en(ji,jj,jk) ) / eps(ji,jj,jk)
893               ! Galperin criterium (NOTE : Not required if the proper value of C3 in stable cases is calculated)
894               zrn2 = MAX( rn2(ji,jj,jk), rsmall )
895               IF (ln_length_lim) mxln(ji,jj,jk) = MIN(  rn_clim_galp * SQRT( 2._wp * en(ji,jj,jk) / zrn2 ), mxln(ji,jj,jk) )
896            END DO
897         END DO
898      END DO 
899
900      !
901      ! Stability function and vertical viscosity and diffusivity
902      ! ---------------------------------------------------------
903      !
904      SELECT CASE ( nn_stab_func )
905      !
906      CASE ( 0 , 1 )             ! Galperin or Kantha-Clayson stability functions
907         DO jk = 2, jpkm1
908            DO jj = 2, jpjm1
909               DO ji = fs_2, fs_jpim1   ! vector opt.
910                  ! zcof =  l²/q²
911                  zcof = mxlb(ji,jj,jk) * mxlb(ji,jj,jk) / ( 2._wp*eb(ji,jj,jk) )
912                  ! Gh = -N²l²/q²
913                  gh = - rn2(ji,jj,jk) * zcof
914                  gh = MIN( gh, rgh0   )
915                  gh = MAX( gh, rghmin )
916                  ! Stability functions from Kantha and Clayson (if C2=C3=0 => Galperin)
917                  sh = ra2*( 1._wp-6._wp*ra1/rb1 ) / ( 1.-3.*ra2*gh*(6.*ra1+rb2*( 1._wp-rc3 ) ) )
918                  sm = ( rb1**(-1._wp/3._wp) + ( 18._wp*ra1*ra1 + 9._wp*ra1*ra2*(1._wp-rc2) )*sh*gh ) / (1._wp-9._wp*ra1*ra2*gh)
919                  !
920                  ! Store stability function in avmu and avmv
921                  avmu(ji,jj,jk) = rc_diff * sh * tmask(ji,jj,jk)
922                  avmv(ji,jj,jk) = rc_diff * sm * tmask(ji,jj,jk)
923               END DO
924            END DO
925         END DO
926         !
927      CASE ( 2, 3 )               ! Canuto stability functions
928         DO jk = 2, jpkm1
929            DO jj = 2, jpjm1
930               DO ji = fs_2, fs_jpim1   ! vector opt.
931                  ! zcof =  l²/q²
932                  zcof = mxlb(ji,jj,jk)*mxlb(ji,jj,jk) / ( 2._wp * eb(ji,jj,jk) )
933                  ! Gh = -N²l²/q²
934                  gh = - rn2(ji,jj,jk) * zcof
935                  gh = MIN( gh, rgh0   )
936                  gh = MAX( gh, rghmin )
937                  gh = gh * rf6
938                  ! Gm =  M²l²/q² Shear number
939                  shr = shear(ji,jj,jk) / MAX( avm(ji,jj,jk), rsmall )
940                  gm = MAX( shr * zcof , 1.e-10 )
941                  gm = gm * rf6
942                  gm = MIN ( (rd0 - rd1*gh + rd3*gh*gh) / (rd2-rd4*gh) , gm )
943                  ! Stability functions from Canuto
944                  rcff = rd0 - rd1*gh +rd2*gm + rd3*gh*gh - rd4*gh*gm + rd5*gm*gm
945                  sm = (rs0 - rs1*gh + rs2*gm) / rcff
946                  sh = (rs4 - rs5*gh + rs6*gm) / rcff
947                  !
948                  ! Store stability function in avmu and avmv
949                  avmu(ji,jj,jk) = rc_diff * sh * tmask(ji,jj,jk)
950                  avmv(ji,jj,jk) = rc_diff * sm * tmask(ji,jj,jk)
951               END DO
952            END DO
953         END DO
954         !
955      END SELECT
956
957      ! Boundary conditions on stability functions for momentum (Neumann):
958      ! Lines below are useless if GOTM style Dirichlet conditions are used
959
960      avmv(:,:,1) = avmv(:,:,2)
961
962      DO jj = 2, jpjm1
963         DO ji = fs_2, fs_jpim1   ! vector opt.
964            avmv(ji,jj,mbkt(ji,jj)+1) = avmv(ji,jj,mbkt(ji,jj))
965         END DO
966      END DO
967
968      ! Compute diffusivities/viscosities
969      ! The computation below could be restrained to jk=2 to jpkm1 if GOTM style Dirichlet conditions are used
970      DO jk = 1, jpk
971         DO jj = 2, jpjm1
972            DO ji = fs_2, fs_jpim1   ! vector opt.
973               zsqen         = SQRT( 2._wp * en(ji,jj,jk) ) * mxln(ji,jj,jk)
974               zav           = zsqen * avmu(ji,jj,jk)
975               avt(ji,jj,jk) = MAX( zav, avtb(jk) )*tmask(ji,jj,jk) ! apply mask for zdfmxl routine
976               zav           = zsqen * avmv(ji,jj,jk)
977               avm(ji,jj,jk) = MAX( zav, avmb(jk) ) ! Note that avm is not masked at the surface and the bottom
978            END DO
979         END DO
980      END DO
981      !
982      ! Lateral boundary conditions (sign unchanged)
983      avt(:,:,1)  = 0._wp
984      CALL lbc_lnk( avm, 'W', 1. )   ;   CALL lbc_lnk( avt, 'W', 1. )
985
986      DO jk = 2, jpkm1            !* vertical eddy viscosity at u- and v-points
987         DO jj = 2, jpjm1
988            DO ji = fs_2, fs_jpim1   ! vector opt.
989               avmu(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji+1,jj  ,jk) ) * umask(ji,jj,jk)
990               avmv(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji  ,jj+1,jk) ) * vmask(ji,jj,jk)
991            END DO
992         END DO
993      END DO
994      avmu(:,:,1) = 0._wp             ;   avmv(:,:,1) = 0._wp                 ! set surface to zero
995      CALL lbc_lnk( avmu, 'U', 1. )   ;   CALL lbc_lnk( avmv, 'V', 1. )       ! Lateral boundary conditions
996
997      IF(ln_ctl) THEN
998         CALL prt_ctl( tab3d_1=en  , clinfo1=' gls  - e: ', tab3d_2=avt, clinfo2=' t: ', ovlap=1, kdim=jpk)
999         CALL prt_ctl( tab3d_1=avmu, clinfo1=' gls  - u: ', mask1=umask,                   &
1000            &          tab3d_2=avmv, clinfo2=       ' v: ', mask2=vmask, ovlap=1, kdim=jpk )
1001      ENDIF
1002      !
1003      avt_k (:,:,:) = avt (:,:,:)
1004      avm_k (:,:,:) = avm (:,:,:)
1005      avmu_k(:,:,:) = avmu(:,:,:)
1006      avmv_k(:,:,:) = avmv(:,:,:)
1007      !
1008      CALL wrk_dealloc( jpi,jpj, zdep, zkar, zflxs, zhsro )
1009      CALL wrk_dealloc( jpi,jpj,jpk, eb, mxlb, shear, eps, zwall_psi, z_elem_a, z_elem_b, z_elem_c, psi )
1010      !
1011      IF( nn_timing == 1 )  CALL timing_stop('zdf_gls')
1012      !
1013      !
1014   END SUBROUTINE zdf_gls
1015
1016
1017   SUBROUTINE zdf_gls_init
1018      !!----------------------------------------------------------------------
1019      !!                  ***  ROUTINE zdf_gls_init  ***
1020      !!                     
1021      !! ** Purpose :   Initialization of the vertical eddy diffivity and
1022      !!      viscosity when using a gls turbulent closure scheme
1023      !!
1024      !! ** Method  :   Read the namzdf_gls namelist and check the parameters
1025      !!      called at the first timestep (nit000)
1026      !!
1027      !! ** input   :   Namlist namzdf_gls
1028      !!
1029      !! ** Action  :   Increase by 1 the nstop flag is setting problem encounter
1030      !!
1031      !!----------------------------------------------------------------------
1032      USE dynzdf_exp
1033      USE trazdf_exp
1034      !
1035      INTEGER ::   jk    ! dummy loop indices
1036      INTEGER ::   ios   ! Local integer output status for namelist read
1037      REAL(wp)::   zcr   ! local scalar
1038      !!
1039      NAMELIST/namzdf_gls/rn_emin, rn_epsmin, ln_length_lim, &
1040         &            rn_clim_galp, ln_sigpsi, rn_hsro,      &
1041         &            rn_crban_default, rn_charn, rn_frac_hs,&
1042         &            nn_bc_surf, nn_bc_bot, nn_z0_met,      &
1043         &            nn_stab_func, nn_clos
1044      !!----------------------------------------------------------
1045      !
1046      IF( nn_timing == 1 )  CALL timing_start('zdf_gls_init')
1047      !
1048      REWIND( numnam_ref )              ! Namelist namzdf_gls in reference namelist : Vertical eddy diffivity and viscosity using gls turbulent closure scheme
1049      READ  ( numnam_ref, namzdf_gls, IOSTAT = ios, ERR = 901)
1050901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_gls in reference namelist', lwp )
1051
1052      REWIND( numnam_cfg )              ! Namelist namzdf_gls in configuration namelist : Vertical eddy diffivity and viscosity using gls turbulent closure scheme
1053      READ  ( numnam_cfg, namzdf_gls, IOSTAT = ios, ERR = 902 )
1054902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_gls in configuration namelist', lwp )
1055      IF(lwm) WRITE ( numond, namzdf_gls )
1056
1057      IF(lwp) THEN                     !* Control print
1058         WRITE(numout,*)
1059         WRITE(numout,*) 'zdf_gls_init : gls turbulent closure scheme'
1060         WRITE(numout,*) '~~~~~~~~~~~~'
1061         WRITE(numout,*) '   Namelist namzdf_gls : set gls mixing parameters'
1062         WRITE(numout,*) '      minimum value of en                           rn_emin        = ', rn_emin
1063         WRITE(numout,*) '      minimum value of eps                          rn_epsmin      = ', rn_epsmin
1064         WRITE(numout,*) '      Limit dissipation rate under stable stratif.  ln_length_lim  = ', ln_length_lim
1065         WRITE(numout,*) '      Galperin limit (Standard: 0.53, Holt: 0.26)   rn_clim_galp   = ', rn_clim_galp
1066         WRITE(numout,*) '      TKE Surface boundary condition                nn_bc_surf     = ', nn_bc_surf
1067         WRITE(numout,*) '      TKE Bottom boundary condition                 nn_bc_bot      = ', nn_bc_bot
1068         WRITE(numout,*) '      Modify psi Schmidt number (wb case)           ln_sigpsi      = ', ln_sigpsi
1069         WRITE(numout,*) '      Craig and Banner coefficient (default)        rn_crban       = ', rn_crban_default
1070         WRITE(numout,*) '      Charnock coefficient                          rn_charn       = ', rn_charn
1071         WRITE(numout,*) '      Surface roughness formula                     nn_z0_met      = ', nn_z0_met
1072         WRITE(numout,*) '      Wave height frac. (used if nn_z0_met=2)       rn_frac_hs     = ', rn_frac_hs
1073         WRITE(numout,*) '      Stability functions                           nn_stab_func   = ', nn_stab_func
1074         WRITE(numout,*) '      Type of closure                               nn_clos        = ', nn_clos
1075         WRITE(numout,*) '      Surface roughness (m)                         rn_hsro        = ', rn_hsro
1076         WRITE(numout,*) '      Bottom roughness (m) (nambfr namelist)        rn_bfrz0       = ', rn_bfrz0
1077      ENDIF
1078
1079      !                                !* allocate gls arrays
1080      IF( zdf_gls_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_gls_init : unable to allocate arrays' )
1081
1082      !                                !* Check of some namelist values
1083      IF( nn_bc_surf < 0 .OR. nn_bc_surf > 1 ) CALL ctl_stop( 'zdf_gls_init: bad flag: nn_bc_surf is 0 or 1' ) 
1084      IF( nn_bc_surf < 0 .OR. nn_bc_surf > 1 ) CALL ctl_stop( 'zdf_gls_init: bad flag: nn_bc_surf is 0 or 1' ) 
1085      IF( nn_z0_met < 0 .OR. nn_z0_met > 3 ) CALL ctl_stop( 'zdf_gls_init: bad flag: nn_z0_met is 0, 1, 2 or 3' ) 
1086      IF( nn_z0_met == 3 .AND. .NOT.ln_rough ) CALL ctl_stop( 'zdf_gls_init: nn_z0_met=3 requires ln_rough=T' ) 
1087      IF( nn_z0_met .NE. 3 .AND. ln_rough ) THEN
1088         CALL ctl_warn('W A R N I N G:  ln_rough=.TRUE. but nn_z0_met is not 3 - resetting nn_z0_met to 3') 
1089         nn_z0_met = 3 
1090      ENDIF
1091      IF( nn_stab_func  < 0 .OR. nn_stab_func  > 3 ) CALL ctl_stop( 'zdf_gls_init: bad flag: nn_stab_func is 0, 1, 2 and 3' ) 
1092      IF( nn_clos       < 0 .OR. nn_clos       > 3 ) CALL ctl_stop( 'zdf_gls_init: bad flag: nn_clos is 0, 1, 2 or 3' )
1093
1094      SELECT CASE ( nn_clos )          !* set the parameters for the chosen closure
1095      !
1096      CASE( 0 )                              ! k-kl  (Mellor-Yamada)
1097         !
1098         IF(lwp) WRITE(numout,*) 'The choosen closure is k-kl closed to the classical Mellor-Yamada'
1099         rpp     = 0._wp
1100         rmm     = 1._wp
1101         rnn     = 1._wp
1102         rsc_tke = 1.96_wp
1103         rsc_psi = 1.96_wp
1104         rpsi1   = 0.9_wp
1105         rpsi3p  = 1._wp
1106         rpsi2   = 0.5_wp
1107         !
1108         SELECT CASE ( nn_stab_func )
1109         CASE( 0, 1 )   ;   rpsi3m = 2.53_wp       ! G88 or KC stability functions
1110         CASE( 2 )      ;   rpsi3m = 2.62_wp       ! Canuto A stability functions
1111         CASE( 3 )      ;   rpsi3m = 2.38          ! Canuto B stability functions (caution : constant not identified)
1112         END SELECT
1113         !
1114      CASE( 1 )                              ! k-eps
1115         !
1116         IF(lwp) WRITE(numout,*) 'The choosen closure is k-eps'
1117         rpp     =  3._wp
1118         rmm     =  1.5_wp
1119         rnn     = -1._wp
1120         rsc_tke =  1._wp
1121         rsc_psi =  1.2_wp  ! Schmidt number for psi
1122         rpsi1   =  1.44_wp
1123         rpsi3p  =  1._wp
1124         rpsi2   =  1.92_wp
1125         !
1126         SELECT CASE ( nn_stab_func )
1127         CASE( 0, 1 )   ;   rpsi3m = -0.52_wp      ! G88 or KC stability functions
1128         CASE( 2 )      ;   rpsi3m = -0.629_wp     ! Canuto A stability functions
1129         CASE( 3 )      ;   rpsi3m = -0.566        ! Canuto B stability functions
1130         END SELECT
1131         !
1132      CASE( 2 )                              ! k-omega
1133         !
1134         IF(lwp) WRITE(numout,*) 'The choosen closure is k-omega'
1135         rpp     = -1._wp
1136         rmm     =  0.5_wp
1137         rnn     = -1._wp
1138         rsc_tke =  2._wp
1139         rsc_psi =  2._wp
1140         rpsi1   =  0.555_wp
1141         rpsi3p  =  1._wp
1142         rpsi2   =  0.833_wp
1143         !
1144         SELECT CASE ( nn_stab_func )
1145         CASE( 0, 1 )   ;   rpsi3m = -0.58_wp       ! G88 or KC stability functions
1146         CASE( 2 )      ;   rpsi3m = -0.64_wp       ! Canuto A stability functions
1147         CASE( 3 )      ;   rpsi3m = -0.64_wp       ! Canuto B stability functions caution : constant not identified)
1148         END SELECT
1149         !
1150      CASE( 3 )                              ! generic
1151         !
1152         IF(lwp) WRITE(numout,*) 'The choosen closure is generic'
1153         rpp     = 2._wp
1154         rmm     = 1._wp
1155         rnn     = -0.67_wp
1156         rsc_tke = 0.8_wp
1157         rsc_psi = 1.07_wp
1158         rpsi1   = 1._wp
1159         rpsi3p  = 1._wp
1160         rpsi2   = 1.22_wp
1161         !
1162         SELECT CASE ( nn_stab_func )
1163         CASE( 0, 1 )   ;   rpsi3m = 0.1_wp         ! G88 or KC stability functions
1164         CASE( 2 )      ;   rpsi3m = 0.05_wp        ! Canuto A stability functions
1165         CASE( 3 )      ;   rpsi3m = 0.05_wp        ! Canuto B stability functions caution : constant not identified)
1166         END SELECT
1167         !
1168      END SELECT
1169
1170      !
1171      SELECT CASE ( nn_stab_func )     !* set the parameters of the stability functions
1172      !
1173      CASE ( 0 )                             ! Galperin stability functions
1174         !
1175         IF(lwp) WRITE(numout,*) 'Stability functions from Galperin'
1176         rc2     =  0._wp
1177         rc3     =  0._wp
1178         rc_diff =  1._wp
1179         rc0     =  0.5544_wp
1180         rcm_sf  =  0.9884_wp
1181         rghmin  = -0.28_wp
1182         rgh0    =  0.0233_wp
1183         rghcri  =  0.02_wp
1184         !
1185      CASE ( 1 )                             ! Kantha-Clayson stability functions
1186         !
1187         IF(lwp) WRITE(numout,*) 'Stability functions from Kantha-Clayson'
1188         rc2     =  0.7_wp
1189         rc3     =  0.2_wp
1190         rc_diff =  1._wp
1191         rc0     =  0.5544_wp
1192         rcm_sf  =  0.9884_wp
1193         rghmin  = -0.28_wp
1194         rgh0    =  0.0233_wp
1195         rghcri  =  0.02_wp
1196         !
1197      CASE ( 2 )                             ! Canuto A stability functions
1198         !
1199         IF(lwp) WRITE(numout,*) 'Stability functions from Canuto A'
1200         rs0 = 1.5_wp * rl1 * rl5*rl5
1201         rs1 = -rl4*(rl6+rl7) + 2._wp*rl4*rl5*(rl1-(1._wp/3._wp)*rl2-rl3) + 1.5_wp*rl1*rl5*rl8
1202         rs2 = -(3._wp/8._wp) * rl1*(rl6*rl6-rl7*rl7)
1203         rs4 = 2._wp * rl5
1204         rs5 = 2._wp * rl4
1205         rs6 = (2._wp/3._wp) * rl5 * ( 3._wp*rl3*rl3 - rl2*rl2 ) - 0.5_wp * rl5*rl1 * (3._wp*rl3-rl2)   &
1206            &                                                    + 0.75_wp * rl1 * ( rl6 - rl7 )
1207         rd0 = 3._wp * rl5*rl5
1208         rd1 = rl5 * ( 7._wp*rl4 + 3._wp*rl8 )
1209         rd2 = rl5*rl5 * ( 3._wp*rl3*rl3 - rl2*rl2 ) - 0.75_wp*(rl6*rl6 - rl7*rl7 )
1210         rd3 = rl4 * ( 4._wp*rl4 + 3._wp*rl8)
1211         rd4 = rl4 * ( rl2 * rl6 - 3._wp*rl3*rl7 - rl5*(rl2*rl2 - rl3*rl3 ) ) + rl5*rl8 * ( 3._wp*rl3*rl3 - rl2*rl2 )
1212         rd5 = 0.25_wp * ( rl2*rl2 - 3._wp *rl3*rl3 ) * ( rl6*rl6 - rl7*rl7 )
1213         rc0 = 0.5268_wp
1214         rf6 = 8._wp / (rc0**6._wp)
1215         rc_diff = SQRT(2._wp) / (rc0**3._wp)
1216         rcm_sf  =  0.7310_wp
1217         rghmin  = -0.28_wp
1218         rgh0    =  0.0329_wp
1219         rghcri  =  0.03_wp
1220         !
1221      CASE ( 3 )                             ! Canuto B stability functions
1222         !
1223         IF(lwp) WRITE(numout,*) 'Stability functions from Canuto B'
1224         rs0 = 1.5_wp * rm1 * rm5*rm5
1225         rs1 = -rm4 * (rm6+rm7) + 2._wp * rm4*rm5*(rm1-(1._wp/3._wp)*rm2-rm3) + 1.5_wp * rm1*rm5*rm8
1226         rs2 = -(3._wp/8._wp) * rm1 * (rm6*rm6-rm7*rm7 )
1227         rs4 = 2._wp * rm5
1228         rs5 = 2._wp * rm4
1229         rs6 = (2._wp/3._wp) * rm5 * (3._wp*rm3*rm3-rm2*rm2) - 0.5_wp * rm5*rm1*(3._wp*rm3-rm2) + 0.75_wp * rm1*(rm6-rm7)
1230         rd0 = 3._wp * rm5*rm5
1231         rd1 = rm5 * (7._wp*rm4 + 3._wp*rm8)
1232         rd2 = rm5*rm5 * (3._wp*rm3*rm3 - rm2*rm2) - 0.75_wp * (rm6*rm6 - rm7*rm7)
1233         rd3 = rm4 * ( 4._wp*rm4 + 3._wp*rm8 )
1234         rd4 = rm4 * ( rm2*rm6 -3._wp*rm3*rm7 - rm5*(rm2*rm2 - rm3*rm3) ) + rm5 * rm8 * ( 3._wp*rm3*rm3 - rm2*rm2 )
1235         rd5 = 0.25_wp * ( rm2*rm2 - 3._wp*rm3*rm3 ) * ( rm6*rm6 - rm7*rm7 )
1236         rc0 = 0.5268_wp            !!       rc0 = 0.5540_wp (Warner ...) to verify !
1237         rf6 = 8._wp / ( rc0**6._wp )
1238         rc_diff = SQRT(2._wp)/(rc0**3.)
1239         rcm_sf  =  0.7470_wp
1240         rghmin  = -0.28_wp
1241         rgh0    =  0.0444_wp
1242         rghcri  =  0.0414_wp
1243         !
1244      END SELECT
1245   
1246      !                                !* Set Schmidt number for psi diffusion in the wave breaking case
1247      !                                     ! See Eq. (13) of Carniel et al, OM, 30, 225-239, 2009
1248      !                                     !  or Eq. (17) of Burchard, JPO, 31, 3133-3145, 2001
1249      IF( ln_sigpsi ) THEN
1250         ra_sf = -1.5 ! Set kinetic energy slope, then deduce rsc_psi and rl_sf
1251         ! Verification: retrieve Burchard (2001) results by uncomenting the line below:
1252         ! Note that the results depend on the value of rn_cm_sf which is constant (=rc0) in his work
1253         ! ra_sf = -SQRT(2./3.*rc0**3./rn_cm_sf*rn_sc_tke)/vkarmn
1254         rsc_psi0 = rsc_tke/(24.*rpsi2)*(-1.+(4.*rnn + ra_sf*(1.+4.*rmm))**2./(ra_sf**2.))
1255      ELSE
1256         rsc_psi0 = rsc_psi
1257      ENDIF
1258 
1259      !                                !* Shear free turbulence parameters
1260      !
1261      ra_sf  = -4._wp*rnn*SQRT(rsc_tke) / ( (1._wp+4._wp*rmm)*SQRT(rsc_tke) &
1262               &                              - SQRT(rsc_tke + 24._wp*rsc_psi0*rpsi2 ) )
1263
1264      IF( .NOT. ln_phioc .AND. rn_crban_default==0._wp ) THEN
1265         rl_sf = vkarmn
1266      ELSE
1267         rl_sf = rc0 * SQRT(rc0/rcm_sf) * SQRT( ( (1._wp + 4._wp*rmm + 8._wp*rmm**2_wp)*rsc_tke          &
1268                 &                                       + 12._wp * rsc_psi0*rpsi2 - (1._wp + 4._wp*rmm) &
1269                 &                                                *SQRT(rsc_tke*(rsc_tke                 &
1270                 &                                                   + 24._wp*rsc_psi0*rpsi2)) )         &
1271                 &                                         /(12._wp*rnn**2.)                             &
1272                 &                                       )
1273      ENDIF
1274
1275      !
1276      IF(lwp) THEN                     !* Control print
1277         WRITE(numout,*)
1278         WRITE(numout,*) 'Limit values'
1279         WRITE(numout,*) '~~~~~~~~~~~~'
1280         WRITE(numout,*) 'Parameter  m = ',rmm
1281         WRITE(numout,*) 'Parameter  n = ',rnn
1282         WRITE(numout,*) 'Parameter  p = ',rpp
1283         WRITE(numout,*) 'rpsi1   = ',rpsi1
1284         WRITE(numout,*) 'rpsi2   = ',rpsi2
1285         WRITE(numout,*) 'rpsi3m  = ',rpsi3m
1286         WRITE(numout,*) 'rpsi3p  = ',rpsi3p
1287         WRITE(numout,*) 'rsc_tke = ',rsc_tke
1288         WRITE(numout,*) 'rsc_psi = ',rsc_psi
1289         WRITE(numout,*) 'rsc_psi0 = ',rsc_psi0
1290         WRITE(numout,*) 'rc0     = ',rc0
1291         WRITE(numout,*)
1292         WRITE(numout,*) 'Shear free turbulence parameters:'
1293         WRITE(numout,*) 'rcm_sf  = ',rcm_sf
1294         WRITE(numout,*) 'ra_sf   = ',ra_sf
1295         WRITE(numout,*) 'rl_sf   = ',rl_sf
1296         WRITE(numout,*)
1297      ENDIF
1298
1299      !                                !* Constants initialization
1300      rc02  = rc0  * rc0   ;   rc02r = 1. / rc02
1301      rc03  = rc02 * rc0
1302      rc04  = rc03 * rc0
1303      rsbc_tke1_default = -3._wp/2._wp*rn_crban_default*ra_sf*rl_sf 
1304      rsbc_tke2         = rdt * rn_crban_default / rl_sf 
1305      zcr               = MAX( rsmall, rsbc_tke1_default**(1./(-ra_sf*3._wp/2._wp))-1._wp ) 
1306      rtrans_default    = 0.2_wp / zcr
1307      rsbc_zs1  = rn_charn/grav                                          ! Charnock formula for surface roughness
1308      rsbc_zs2  = rn_frac_hs / 0.85_wp / grav * 665._wp                  ! Rascle formula for surface roughness
1309      rsbc_psi1 = -0.5_wp * rdt * rc0**(rpp-2._wp*rmm) / rsc_psi
1310      rsbc_psi2 = -0.5_wp * rdt * rc0**rpp * rnn * vkarmn**rnn / rsc_psi ! Neumann + NO Wave breaking
1311
1312      rfact_tke = -0.5_wp / rsc_tke * rdt                                ! Cst used for the Diffusion term of tke
1313      rfact_psi = -0.5_wp / rsc_psi * rdt                                ! Cst used for the Diffusion term of tke
1314
1315      !                                !* Wall proximity function
1316      zwall (:,:,:) = 1._wp * tmask(:,:,:)
1317
1318      !                                !* set vertical eddy coef. to the background value
1319      DO jk = 1, jpk
1320         avt (:,:,jk) = avtb(jk) * tmask(:,:,jk)
1321         avm (:,:,jk) = avmb(jk) * tmask(:,:,jk)
1322         avmu(:,:,jk) = avmb(jk) * umask(:,:,jk)
1323         avmv(:,:,jk) = avmb(jk) * vmask(:,:,jk)
1324      END DO
1325      !                             
1326      CALL gls_rst( nit000, 'READ' )   !* read or initialize all required files
1327      !
1328      IF( nn_timing == 1 )  CALL timing_stop('zdf_gls_init')
1329      !
1330   END SUBROUTINE zdf_gls_init
1331
1332
1333   SUBROUTINE gls_rst( kt, cdrw )
1334      !!---------------------------------------------------------------------
1335      !!                   ***  ROUTINE ts_rst  ***
1336      !!                     
1337      !! ** Purpose :   Read or write TKE file (en) in restart file
1338      !!
1339      !! ** Method  :   use of IOM library
1340      !!                if the restart does not contain TKE, en is either
1341      !!                set to rn_emin or recomputed (nn_igls/=0)
1342      !!----------------------------------------------------------------------
1343      INTEGER         , INTENT(in) ::   kt         ! ocean time-step
1344      CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
1345      !
1346      INTEGER ::   jit, jk   ! dummy loop indices
1347      INTEGER ::   id1, id2, id3, id4, id5, id6
1348      INTEGER ::   ji, jj, ikbu, ikbv
1349      REAL(wp)::   cbx, cby
1350      !!----------------------------------------------------------------------
1351      !
1352      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
1353         !                                   ! ---------------
1354         IF( ln_rstart ) THEN                   !* Read the restart file
1355            id1 = iom_varid( numror, 'en'   , ldstop = .FALSE. )
1356            id2 = iom_varid( numror, 'avt'  , ldstop = .FALSE. )
1357            id3 = iom_varid( numror, 'avm'  , ldstop = .FALSE. )
1358            id4 = iom_varid( numror, 'avmu' , ldstop = .FALSE. )
1359            id5 = iom_varid( numror, 'avmv' , ldstop = .FALSE. )
1360            id6 = iom_varid( numror, 'mxln' , ldstop = .FALSE. )
1361            !
1362            IF( MIN( id1, id2, id3, id4, id5, id6 ) > 0 ) THEN        ! all required arrays exist
1363               CALL iom_get( numror, jpdom_autoglo, 'en'    , en     )
1364               CALL iom_get( numror, jpdom_autoglo, 'avt'   , avt    )
1365               CALL iom_get( numror, jpdom_autoglo, 'avm'   , avm    )
1366               CALL iom_get( numror, jpdom_autoglo, 'avmu'  , avmu   )
1367               CALL iom_get( numror, jpdom_autoglo, 'avmv'  , avmv   )
1368               CALL iom_get( numror, jpdom_autoglo, 'mxln'  , mxln   )
1369            ELSE                       
1370               IF(lwp) WRITE(numout,*) ' ===>>>> : previous run without gls scheme, en and mxln computed by iterative loop'
1371               en  (:,:,:) = rn_emin
1372               mxln(:,:,:) = 0.05       
1373               avt_k (:,:,:) = avt (:,:,:)
1374               avm_k (:,:,:) = avm (:,:,:)
1375               avmu_k(:,:,:) = avmu(:,:,:)
1376               avmv_k(:,:,:) = avmv(:,:,:)
1377               DO jit = nit000 + 1, nit000 + 10   ;   CALL zdf_gls( jit )   ;   END DO
1378            ENDIF
1379         ELSE                                   !* Start from rest
1380            IF(lwp) WRITE(numout,*) ' ===>>>> : Initialisation of en and mxln by background values'
1381            en  (:,:,:) = rn_emin
1382            mxln(:,:,:) = 0.05       
1383         ENDIF
1384         !
1385      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
1386         !                                   ! -------------------
1387         IF(lwp) WRITE(numout,*) '---- gls-rst ----'
1388         CALL iom_rstput( kt, nitrst, numrow, 'en'   , en     ) 
1389         CALL iom_rstput( kt, nitrst, numrow, 'avt'  , avt_k  )
1390         CALL iom_rstput( kt, nitrst, numrow, 'avm'  , avm_k  )
1391         CALL iom_rstput( kt, nitrst, numrow, 'avmu' , avmu_k ) 
1392         CALL iom_rstput( kt, nitrst, numrow, 'avmv' , avmv_k )
1393         CALL iom_rstput( kt, nitrst, numrow, 'mxln' , mxln   )
1394         !
1395      ENDIF
1396      !
1397   END SUBROUTINE gls_rst
1398
1399#else
1400   !!----------------------------------------------------------------------
1401   !!   Dummy module :                                        NO TKE scheme
1402   !!----------------------------------------------------------------------
1403   LOGICAL, PUBLIC, PARAMETER ::   lk_zdfgls = .FALSE.   !: TKE flag
1404CONTAINS
1405   SUBROUTINE zdf_gls_init           ! Empty routine
1406      WRITE(*,*) 'zdf_gls_init: You should not have seen this print! error?'
1407   END SUBROUTINE zdf_gls_init
1408   SUBROUTINE zdf_gls( kt )          ! Empty routine
1409      WRITE(*,*) 'zdf_gls: You should not have seen this print! error?', kt
1410   END SUBROUTINE zdf_gls
1411   SUBROUTINE gls_rst( kt, cdrw )          ! Empty routine
1412      INTEGER         , INTENT(in) ::   kt         ! ocean time-step
1413      CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
1414      WRITE(*,*) 'gls_rst: You should not have seen this print! error?', kt, cdrw
1415   END SUBROUTINE gls_rst
1416#endif
1417
1418   !!======================================================================
1419END MODULE zdfgls
1420
Note: See TracBrowser for help on using the repository browser.