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_r7963_nemo_v3_6_AGRIF-3_AGRIFVVL/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2017/dev_r7963_nemo_v3_6_AGRIF-3_AGRIFVVL/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfgls.F90 @ 7988

Last change on this file since 7988 was 7988, checked in by jchanut, 7 years ago

Add AGRIF proper AGRIF bcs to GLS and TKE + vvl update

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