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

source: trunk/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfgls.F90 @ 6140

Last change on this file since 6140 was 6140, checked in by timgraham, 8 years ago

Merge of branches/2015/dev_merge_2015 back into trunk. Merge excludes NEMOGCM/TOOLS/OBSTOOLS/ for now due to issues with the change of file type. Will sort these manually with further commits.

Branch merged as follows:
In the working copy of branch ran:
svn merge svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/trunk@HEAD
Small conflicts due to bug fixes applied to trunk since the dev_merge_2015 was copied. Bug fixes were applied to the branch as well so these were easy to resolve.
Branch committed at this stage

In working copy run:
svn switch svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/trunk
to switch working copy

Run:
svn merge --reintegrate svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/branches/2015/dev_merge_2015
to merge the branch into the trunk and then commit - no conflicts at this stage.

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