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/2017/dev_r8624_AGRIF3_VVL/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2017/dev_r8624_AGRIF3_VVL/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfgls.F90 @ 8741

Last change on this file since 8741 was 8741, checked in by jchanut, 6 years ago

AGRIF + vvl Main changes - #1965

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