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.
zdftke.F90 in NEMO/trunk/src/OCE/ZDF – NEMO

source: NEMO/trunk/src/OCE/ZDF/zdftke.F90 @ 13226

Last change on this file since 13226 was 13226, checked in by orioltp, 4 years ago

Merging dev_r12512_HPC-04_mcastril_Mixed_Precision_implementation into the trunk

  • Property svn:keywords set to Id
File size: 41.5 KB
Line 
1MODULE zdftke
2   !!======================================================================
3   !!                       ***  MODULE  zdftke  ***
4   !! Ocean physics:  vertical mixing coefficient computed from the tke
5   !!                 turbulent closure parameterization
6   !!=====================================================================
7   !! History :  OPA  !  1991-03  (b. blanke)  Original code
8   !!            7.0  !  1991-11  (G. Madec)   bug fix
9   !!            7.1  !  1992-10  (G. Madec)   new mixing length and eav
10   !!            7.2  !  1993-03  (M. Guyon)   symetrical conditions
11   !!            7.3  !  1994-08  (G. Madec, M. Imbard)  nn_pdl flag
12   !!            7.5  !  1996-01  (G. Madec)   s-coordinates
13   !!            8.0  !  1997-07  (G. Madec)   lbc
14   !!            8.1  !  1999-01  (E. Stretta) new option for the mixing length
15   !!  NEMO      1.0  !  2002-06  (G. Madec) add tke_init routine
16   !!             -   !  2004-10  (C. Ethe )  1D configuration
17   !!            2.0  !  2006-07  (S. Masson)  distributed restart using iom
18   !!            3.0  !  2008-05  (C. Ethe,  G.Madec) : update TKE physics:
19   !!                 !           - tke penetration (wind steering)
20   !!                 !           - suface condition for tke & mixing length
21   !!                 !           - Langmuir cells
22   !!             -   !  2008-05  (J.-M. Molines, G. Madec)  2D form of avtb
23   !!             -   !  2008-06  (G. Madec)  style + DOCTOR name for namelist parameters
24   !!             -   !  2008-12  (G. Reffray) stable discretization of the production term
25   !!            3.2  !  2009-06  (G. Madec, S. Masson) TKE restart compatible with key_cpl
26   !!                 !                                + cleaning of the parameters + bugs correction
27   !!            3.3  !  2010-10  (C. Ethe, G. Madec) reorganisation of initialisation phase
28   !!            3.6  !  2014-11  (P. Mathiot) add ice shelf capability
29   !!            4.0  !  2017-04  (G. Madec)  remove CPP ddm key & avm at t-point only
30   !!             -   !  2017-05  (G. Madec)  add top/bottom friction as boundary condition (ln_drg)
31   !!----------------------------------------------------------------------
32
33   !!----------------------------------------------------------------------
34   !!   zdf_tke       : update momentum and tracer Kz from a tke scheme
35   !!   tke_tke       : tke time stepping: update tke at now time step (en)
36   !!   tke_avn       : compute mixing length scale and deduce avm and avt
37   !!   zdf_tke_init  : initialization, namelist read, and parameters control
38   !!   tke_rst       : read/write tke restart in ocean restart file
39   !!----------------------------------------------------------------------
40   USE oce            ! ocean: dynamics and active tracers variables
41   USE phycst         ! physical constants
42   USE dom_oce        ! domain: ocean
43   USE domvvl         ! domain: variable volume layer
44   USE sbc_oce        ! surface boundary condition: ocean
45   USE zdfdrg         ! vertical physics: top/bottom drag coef.
46   USE zdfmxl         ! vertical physics: mixed layer
47#if defined key_si3
48   USE ice, ONLY: hm_i, h_i
49#endif
50#if defined key_cice
51   USE sbc_ice, ONLY: h_i
52#endif
53   !
54   USE in_out_manager ! I/O manager
55   USE iom            ! I/O manager library
56   USE lib_mpp        ! MPP library
57   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
58   USE prtctl         ! Print control
59   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
60
61   IMPLICIT NONE
62   PRIVATE
63
64   PUBLIC   zdf_tke        ! routine called in step module
65   PUBLIC   zdf_tke_init   ! routine called in opa module
66   PUBLIC   tke_rst        ! routine called in step module
67
68   !                      !!** Namelist  namzdf_tke  **
69   LOGICAL  ::   ln_mxl0   ! mixing length scale surface value as function of wind stress or not
70   INTEGER  ::   nn_mxl    ! type of mixing length (=0/1/2/3)
71   REAL(wp) ::   rn_mxl0   ! surface  min value of mixing length (kappa*z_o=0.4*0.1 m)  [m]
72   INTEGER  ::      nn_mxlice ! type of scaling under sea-ice
73   REAL(wp) ::      rn_mxlice ! max constant ice thickness value when scaling under sea-ice ( nn_mxlice=1)
74   INTEGER  ::   nn_pdl    ! Prandtl number or not (ratio avt/avm) (=0/1)
75   REAL(wp) ::   rn_ediff  ! coefficient for avt: avt=rn_ediff*mxl*sqrt(e)
76   REAL(wp) ::   rn_ediss  ! coefficient of the Kolmogoroff dissipation
77   REAL(wp) ::   rn_ebb    ! coefficient of the surface input of tke
78   REAL(wp) ::   rn_emin   ! minimum value of tke           [m2/s2]
79   REAL(wp) ::   rn_emin0  ! surface minimum value of tke   [m2/s2]
80   REAL(wp) ::   rn_bshear ! background shear (>0) currently a numerical threshold (do not change it)
81   LOGICAL  ::   ln_drg    ! top/bottom friction forcing flag
82   INTEGER  ::   nn_etau   ! type of depth penetration of surface tke (=0/1/2/3)
83   INTEGER  ::      nn_htau   ! type of tke profile of penetration (=0/1)
84   REAL(wp) ::      rn_efr    ! fraction of TKE surface value which penetrates in the ocean
85   REAL(wp) ::      rn_eice   ! =0 ON below sea-ice, =4 OFF when ice fraction > 1/4   
86   LOGICAL  ::   ln_lc     ! Langmuir cells (LC) as a source term of TKE or not
87   REAL(wp) ::      rn_lc     ! coef to compute vertical velocity of Langmuir cells
88
89   REAL(wp) ::   ri_cri    ! critic Richardson number (deduced from rn_ediff and rn_ediss values)
90   REAL(wp) ::   rmxl_min  ! minimum mixing length value (deduced from rn_ediff and rn_emin values)  [m]
91   REAL(wp) ::   rhftau_add = 1.e-3_wp     ! add offset   applied to HF part of taum  (nn_etau=3)
92   REAL(wp) ::   rhftau_scl = 1.0_wp       ! scale factor applied to HF part of taum  (nn_etau=3)
93
94   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   htau    ! depth of tke penetration (nn_htau)
95   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dissl   ! now mixing lenght of dissipation
96   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   apdlr   ! now mixing lenght of dissipation
97
98   !! * Substitutions
99#  include "do_loop_substitute.h90"
100   !!----------------------------------------------------------------------
101   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
102   !! $Id$
103   !! Software governed by the CeCILL license (see ./LICENSE)
104   !!----------------------------------------------------------------------
105CONTAINS
106
107   INTEGER FUNCTION zdf_tke_alloc()
108      !!----------------------------------------------------------------------
109      !!                ***  FUNCTION zdf_tke_alloc  ***
110      !!----------------------------------------------------------------------
111      ALLOCATE( htau(jpi,jpj) , dissl(jpi,jpj,jpk) , apdlr(jpi,jpj,jpk) ,   STAT= zdf_tke_alloc )
112      !
113      CALL mpp_sum ( 'zdftke', zdf_tke_alloc )
114      IF( zdf_tke_alloc /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_alloc: failed to allocate arrays' )
115      !
116   END FUNCTION zdf_tke_alloc
117
118
119   SUBROUTINE zdf_tke( kt, Kbb, Kmm, p_sh2, p_avm, p_avt )
120      !!----------------------------------------------------------------------
121      !!                   ***  ROUTINE zdf_tke  ***
122      !!
123      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
124      !!              coefficients using a turbulent closure scheme (TKE).
125      !!
126      !! ** Method  :   The time evolution of the turbulent kinetic energy (tke)
127      !!              is computed from a prognostic equation :
128      !!         d(en)/dt = avm (d(u)/dz)**2             ! shear production
129      !!                  + d( avm d(en)/dz )/dz         ! diffusion of tke
130      !!                  + avt N^2                      ! stratif. destruc.
131      !!                  - rn_ediss / emxl en**(2/3)    ! Kolmogoroff dissipation
132      !!      with the boundary conditions:
133      !!         surface: en = max( rn_emin0, rn_ebb * taum )
134      !!         bottom : en = rn_emin
135      !!      The associated critical Richardson number is: ri_cri = 2/(2+rn_ediss/rn_ediff)
136      !!
137      !!        The now Turbulent kinetic energy is computed using the following
138      !!      time stepping: implicit for vertical diffusion term, linearized semi
139      !!      implicit for kolmogoroff dissipation term, and explicit forward for
140      !!      both buoyancy and shear production terms. Therefore a tridiagonal
141      !!      linear system is solved. Note that buoyancy and shear terms are
142      !!      discretized in a energy conserving form (Bruchard 2002).
143      !!
144      !!        The dissipative and mixing length scale are computed from en and
145      !!      the stratification (see tke_avn)
146      !!
147      !!        The now vertical eddy vicosity and diffusivity coefficients are
148      !!      given by:
149      !!              avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
150      !!              avt = max( avmb, pdl * avm                 ) 
151      !!              eav = max( avmb, avm )
152      !!      where pdl, the inverse of the Prandtl number is 1 if nn_pdl=0 and
153      !!      given by an empirical funtion of the localRichardson number if nn_pdl=1
154      !!
155      !! ** Action  :   compute en (now turbulent kinetic energy)
156      !!                update avt, avm (before vertical eddy coef.)
157      !!
158      !! References : Gaspar et al., JGR, 1990,
159      !!              Blanke and Delecluse, JPO, 1991
160      !!              Mellor and Blumberg, JPO 2004
161      !!              Axell, JGR, 2002
162      !!              Bruchard OM 2002
163      !!----------------------------------------------------------------------
164      INTEGER                   , INTENT(in   ) ::   kt             ! ocean time step
165      INTEGER                   , INTENT(in   ) ::   Kbb, Kmm       ! ocean time level indices
166      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   p_sh2          ! shear production term
167      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   p_avm, p_avt   !  momentum and tracer Kz (w-points)
168      !!----------------------------------------------------------------------
169      !
170      CALL tke_tke( Kbb, Kmm, p_sh2, p_avm, p_avt )   ! now tke (en)
171      !
172      CALL tke_avn( Kbb, Kmm,        p_avm, p_avt )   ! now avt, avm, dissl
173      !
174  END SUBROUTINE zdf_tke
175
176
177   SUBROUTINE tke_tke( Kbb, Kmm, p_sh2, p_avm, p_avt )
178      !!----------------------------------------------------------------------
179      !!                   ***  ROUTINE tke_tke  ***
180      !!
181      !! ** Purpose :   Compute the now Turbulente Kinetic Energy (TKE)
182      !!
183      !! ** Method  : - TKE surface boundary condition
184      !!              - source term due to Langmuir cells (Axell JGR 2002) (ln_lc=T)
185      !!              - source term due to shear (= Kz dz[Ub] * dz[Un] )
186      !!              - Now TKE : resolution of the TKE equation by inverting
187      !!                a tridiagonal linear system by a "methode de chasse"
188      !!              - increase TKE due to surface and internal wave breaking
189      !!             NB: when sea-ice is present, both LC parameterization
190      !!                 and TKE penetration are turned off when the ice fraction
191      !!                 is smaller than 0.25
192      !!
193      !! ** Action  : - en : now turbulent kinetic energy)
194      !! ---------------------------------------------------------------------
195      USE zdf_oce , ONLY : en   ! ocean vertical physics
196      !!
197      INTEGER                    , INTENT(in   ) ::   Kbb, Kmm       ! ocean time level indices
198      REAL(wp), DIMENSION(:,:,:) , INTENT(in   ) ::   p_sh2          ! shear production term
199      REAL(wp), DIMENSION(:,:,:) , INTENT(in   ) ::   p_avm, p_avt   ! vertical eddy viscosity & diffusivity (w-points)
200      !
201      INTEGER ::   ji, jj, jk              ! dummy loop arguments
202      REAL(wp) ::   zetop, zebot, zmsku, zmskv ! local scalars
203      REAL(wp) ::   zrhoa  = 1.22              ! Air density kg/m3
204      REAL(wp) ::   zcdrag = 1.5e-3            ! drag coefficient
205      REAL(wp) ::   zbbrau, zri                ! local scalars
206      REAL(wp) ::   zfact1, zfact2, zfact3     !   -         -
207      REAL(wp) ::   ztx2  , zty2  , zcof       !   -         -
208      REAL(wp) ::   ztau  , zdif               !   -         -
209      REAL(wp) ::   zus   , zwlc  , zind       !   -         -
210      REAL(wp) ::   zzd_up, zzd_lw             !   -         -
211      INTEGER , DIMENSION(jpi,jpj)     ::   imlc
212      REAL(wp), DIMENSION(jpi,jpj)     ::   zhlc, zfr_i
213      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zpelc, zdiag, zd_up, zd_lw
214      !!--------------------------------------------------------------------
215      !
216      zbbrau = rn_ebb / rho0       ! Local constant initialisation
217      zfact1 = -.5_wp * rn_Dt 
218      zfact2 = 1.5_wp * rn_Dt * rn_ediss
219      zfact3 = 0.5_wp       * rn_ediss
220      !
221      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
222      !                     !  Surface/top/bottom boundary condition on tke
223      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
224      !
225      DO_2D_00_00
226         en(ji,jj,1) = MAX( rn_emin0, zbbrau * taum(ji,jj) ) * tmask(ji,jj,1)
227      END_2D
228      !
229      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
230      !                     !  Bottom boundary condition on tke
231      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
232      !
233      !   en(bot)   = (ebb0/rho0)*0.5*sqrt(u_botfr^2+v_botfr^2) (min value rn_emin)
234      ! where ebb0 does not includes surface wave enhancement (i.e. ebb0=3.75)
235      ! Note that stress averaged is done using an wet-only calculation of u and v at t-point like in zdfsh2
236      !
237      IF( ln_drg ) THEN       !== friction used as top/bottom boundary condition on TKE
238         !
239         DO_2D_00_00
240            zmsku = ( 2. - umask(ji-1,jj,mbkt(ji,jj)) * umask(ji,jj,mbkt(ji,jj)) )
241            zmskv = ( 2. - vmask(ji,jj-1,mbkt(ji,jj)) * vmask(ji,jj,mbkt(ji,jj)) )
242            !                       ! where 0.001875 = (rn_ebb0/rho0) * 0.5 = 3.75*0.5/1000. (CAUTION CdU<0)
243            zebot = - 0.001875_wp * rCdU_bot(ji,jj) * SQRT(  ( zmsku*( uu(ji,jj,mbkt(ji,jj),Kbb)+uu(ji-1,jj,mbkt(ji,jj),Kbb) ) )**2  &
244               &                                           + ( zmskv*( vv(ji,jj,mbkt(ji,jj),Kbb)+vv(ji,jj-1,mbkt(ji,jj),Kbb) ) )**2  )
245            en(ji,jj,mbkt(ji,jj)+1) = MAX( zebot, rn_emin ) * ssmask(ji,jj)
246         END_2D
247         IF( ln_isfcav ) THEN       ! top friction
248            DO_2D_00_00
249               zmsku = ( 2. - umask(ji-1,jj,mikt(ji,jj)) * umask(ji,jj,mikt(ji,jj)) )
250               zmskv = ( 2. - vmask(ji,jj-1,mikt(ji,jj)) * vmask(ji,jj,mikt(ji,jj)) )
251               !                             ! where 0.001875 = (rn_ebb0/rho0) * 0.5 = 3.75*0.5/1000.  (CAUTION CdU<0)
252               zetop = - 0.001875_wp * rCdU_top(ji,jj) * SQRT(  ( zmsku*( uu(ji,jj,mikt(ji,jj),Kbb)+uu(ji-1,jj,mikt(ji,jj),Kbb) ) )**2  &
253                  &                                           + ( zmskv*( vv(ji,jj,mikt(ji,jj),Kbb)+vv(ji,jj-1,mikt(ji,jj),Kbb) ) )**2  )
254               ! (1._wp - tmask(ji,jj,1)) * ssmask(ji,jj) = 1 where ice shelves are present
255               en(ji,jj,mikt(ji,jj)) = en(ji,jj,1)           * tmask(ji,jj,1) &
256                  &                  + MAX( zetop, rn_emin ) * (1._wp - tmask(ji,jj,1)) * ssmask(ji,jj)
257            END_2D
258         ENDIF
259         !
260      ENDIF
261      !
262      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
263      IF( ln_lc ) THEN      !  Langmuir circulation source term added to tke   !   (Axell JGR 2002)
264         !                  !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
265         !
266         !                        !* total energy produce by LC : cumulative sum over jk
267         zpelc(:,:,1) =  MAX( rn2b(:,:,1), 0._wp ) * gdepw(:,:,1,Kmm) * e3w(:,:,1,Kmm)
268         DO jk = 2, jpk
269            zpelc(:,:,jk)  = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0._wp ) * gdepw(:,:,jk,Kmm) * e3w(:,:,jk,Kmm)
270         END DO
271         !                        !* finite Langmuir Circulation depth
272         zcof = 0.5 * 0.016 * 0.016 / ( zrhoa * zcdrag )
273         imlc(:,:) = mbkt(:,:) + 1       ! Initialization to the number of w ocean point (=2 over land)
274         DO_3DS_11_11( jpkm1, 2, -1 )
275            zus  = zcof * taum(ji,jj)
276            IF( zpelc(ji,jj,jk) > zus )   imlc(ji,jj) = jk
277         END_3D
278         !                               ! finite LC depth
279         DO_2D_11_11
280            zhlc(ji,jj) = gdepw(ji,jj,imlc(ji,jj),Kmm)
281         END_2D
282         zcof = 0.016 / SQRT( zrhoa * zcdrag )
283         DO_2D_00_00
284            zus  = zcof * SQRT( taum(ji,jj) )           ! Stokes drift
285            zfr_i(ji,jj) = ( 1._wp - 4._wp * fr_i(ji,jj) ) * zus * zus * zus * tmask(ji,jj,1) ! zus > 0. ok
286            IF (zfr_i(ji,jj) < 0. ) zfr_i(ji,jj) = 0.
287         END_2D
288         DO_3D_00_00( 2, jpkm1 )
289            IF ( zfr_i(ji,jj) /= 0. ) THEN               
290               ! vertical velocity due to LC   
291               IF ( gdepw(ji,jj,jk,Kmm) - zhlc(ji,jj) < 0 .AND. wmask(ji,jj,jk) /= 0. ) THEN
292                  !                                           ! vertical velocity due to LC
293                  zwlc = rn_lc * SIN( rpi * gdepw(ji,jj,jk,Kmm) / zhlc(ji,jj) )   ! warning: optimization: zus^3 is in zfr_i
294                  !                                           ! TKE Langmuir circulation source term
295                  en(ji,jj,jk) = en(ji,jj,jk) + rn_Dt * zfr_i(ji,jj) * ( zwlc * zwlc * zwlc ) / zhlc(ji,jj)
296               ENDIF
297            ENDIF
298         END_3D
299         !
300      ENDIF
301      !
302      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
303      !                     !  Now Turbulent kinetic energy (output in en)
304      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
305      !                     ! Resolution of a tridiagonal linear system by a "methode de chasse"
306      !                     ! computation from level 2 to jpkm1  (e(1) already computed and e(jpk)=0 ).
307      !                     ! zdiag : diagonal zd_up : upper diagonal zd_lw : lower diagonal
308      !
309      IF( nn_pdl == 1 ) THEN      !* Prandtl number = F( Ri )
310         DO_3D_00_00( 2, jpkm1 )
311            !                             ! local Richardson number
312            IF (rn2b(ji,jj,jk) <= 0.0_wp) then
313                zri = 0.0_wp
314            ELSE
315                zri = rn2b(ji,jj,jk) * p_avm(ji,jj,jk) / ( p_sh2(ji,jj,jk) + rn_bshear )
316            ENDIF
317            !                             ! inverse of Prandtl number
318            apdlr(ji,jj,jk) = MAX(  0.1_wp,  ri_cri / MAX( ri_cri , zri )  )
319         END_3D
320      ENDIF
321      !         
322      DO_3D_00_00( 2, jpkm1 )
323         zcof   = zfact1 * tmask(ji,jj,jk)
324         !                                   ! A minimum of 2.e-5 m2/s is imposed on TKE vertical
325         !                                   ! eddy coefficient (ensure numerical stability)
326         zzd_up = zcof * MAX(  p_avm(ji,jj,jk+1) + p_avm(ji,jj,jk  ) , 2.e-5_wp  )   &  ! upper diagonal
327            &          /    (  e3t(ji,jj,jk  ,Kmm) * e3w(ji,jj,jk  ,Kmm)  )
328         zzd_lw = zcof * MAX(  p_avm(ji,jj,jk  ) + p_avm(ji,jj,jk-1) , 2.e-5_wp  )   &  ! lower diagonal
329            &          /    (  e3t(ji,jj,jk-1,Kmm) * e3w(ji,jj,jk  ,Kmm)  )
330         !
331         zd_up(ji,jj,jk) = zzd_up            ! Matrix (zdiag, zd_up, zd_lw)
332         zd_lw(ji,jj,jk) = zzd_lw
333         zdiag(ji,jj,jk) = 1._wp - zzd_lw - zzd_up + zfact2 * dissl(ji,jj,jk) * wmask(ji,jj,jk)
334         !
335         !                                   ! right hand side in en
336         en(ji,jj,jk) = en(ji,jj,jk) + rn_Dt * (  p_sh2(ji,jj,jk)                        &   ! shear
337            &                                 - p_avt(ji,jj,jk) * rn2(ji,jj,jk)          &   ! stratification
338            &                                 + zfact3 * dissl(ji,jj,jk) * en(ji,jj,jk)  &   ! dissipation
339            &                                ) * wmask(ji,jj,jk)
340      END_3D
341      !                          !* Matrix inversion from level 2 (tke prescribed at level 1)
342      DO_3D_00_00( 3, jpkm1 )
343         zdiag(ji,jj,jk) = zdiag(ji,jj,jk) - zd_lw(ji,jj,jk) * zd_up(ji,jj,jk-1) / zdiag(ji,jj,jk-1)
344      END_3D
345      DO_2D_00_00
346         zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1)    ! Surface boudary conditions on tke
347      END_2D
348      DO_3D_00_00( 3, jpkm1 )
349         zd_lw(ji,jj,jk) = en(ji,jj,jk) - zd_lw(ji,jj,jk) / zdiag(ji,jj,jk-1) *zd_lw(ji,jj,jk-1)
350      END_3D
351      DO_2D_00_00
352         en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1)
353      END_2D
354      DO_3DS_00_00( jpk-2, 2, -1 )
355         en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk)
356      END_3D
357      DO_3D_00_00( 2, jpkm1 )
358         en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * wmask(ji,jj,jk)
359      END_3D
360      !
361      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
362      !                            !  TKE due to surface and internal wave breaking
363      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
364!!gm BUG : in the exp  remove the depth of ssh !!!
365!!gm       i.e. use gde3w in argument (gdepw(:,:,:,Kmm))
366     
367     
368      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction)
369         DO_3D_00_00( 2, jpkm1 )
370            en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -gdepw(ji,jj,jk,Kmm) / htau(ji,jj) )   &
371               &                                 * MAX(0.,1._wp - rn_eice *fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
372         END_3D
373      ELSEIF( nn_etau == 2 ) THEN       !* act only at the base of the mixed layer (jk=nmln)  (rn_efr fraction)
374         DO_2D_00_00
375            jk = nmln(ji,jj)
376            en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -gdepw(ji,jj,jk,Kmm) / htau(ji,jj) )   &
377               &                                 * MAX(0.,1._wp - rn_eice *fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
378         END_2D
379      ELSEIF( nn_etau == 3 ) THEN       !* penetration belox the mixed layer (HF variability)
380         DO_3D_00_00( 2, jpkm1 )
381            ztx2 = utau(ji-1,jj  ) + utau(ji,jj)
382            zty2 = vtau(ji  ,jj-1) + vtau(ji,jj)
383            ztau = 0.5_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1)    ! module of the mean stress
384            zdif = taum(ji,jj) - ztau                            ! mean of modulus - modulus of the mean
385            zdif = rhftau_scl * MAX( 0._wp, zdif + rhftau_add )  ! apply some modifications...
386            en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -gdepw(ji,jj,jk,Kmm) / htau(ji,jj) )   &
387               &                        * MAX(0.,1._wp - rn_eice *fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1)
388         END_3D
389      ENDIF
390      !
391   END SUBROUTINE tke_tke
392
393
394   SUBROUTINE tke_avn( Kbb, Kmm, p_avm, p_avt )
395      !!----------------------------------------------------------------------
396      !!                   ***  ROUTINE tke_avn  ***
397      !!
398      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
399      !!
400      !! ** Method  :   At this stage, en, the now TKE, is known (computed in
401      !!              the tke_tke routine). First, the now mixing lenth is
402      !!      computed from en and the strafification (N^2), then the mixings
403      !!      coefficients are computed.
404      !!              - Mixing length : a first evaluation of the mixing lengh
405      !!      scales is:
406      !!                      mxl = sqrt(2*en) / N 
407      !!      where N is the brunt-vaisala frequency, with a minimum value set
408      !!      to rmxl_min (rn_mxl0) in the interior (surface) ocean.
409      !!        The mixing and dissipative length scale are bound as follow :
410      !!         nn_mxl=0 : mxl bounded by the distance to surface and bottom.
411      !!                        zmxld = zmxlm = mxl
412      !!         nn_mxl=1 : mxl bounded by the e3w and zmxld = zmxlm = mxl
413      !!         nn_mxl=2 : mxl bounded such that the vertical derivative of mxl is
414      !!                    less than 1 (|d/dz(mxl)|<1) and zmxld = zmxlm = mxl
415      !!         nn_mxl=3 : mxl is bounded from the surface to the bottom usings
416      !!                    |d/dz(xml)|<1 to obtain lup, and from the bottom to
417      !!                    the surface to obtain ldown. the resulting length
418      !!                    scales are:
419      !!                        zmxld = sqrt( lup * ldown )
420      !!                        zmxlm = min ( lup , ldown )
421      !!              - Vertical eddy viscosity and diffusivity:
422      !!                      avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
423      !!                      avt = max( avmb, pdlr * avm ) 
424      !!      with pdlr=1 if nn_pdl=0, pdlr=1/pdl=F(Ri) otherwise.
425      !!
426      !! ** Action  : - avt, avm : now vertical eddy diffusivity and viscosity (w-point)
427      !!----------------------------------------------------------------------
428      USE zdf_oce , ONLY : en, avtb, avmb, avtb_2d   ! ocean vertical physics
429      !!
430      INTEGER                   , INTENT(in   ) ::   Kbb, Kmm       ! ocean time level indices
431      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_avm, p_avt   ! vertical eddy viscosity & diffusivity (w-points)
432      !
433      INTEGER  ::   ji, jj, jk   ! dummy loop indices
434      REAL(wp) ::   zrn2, zraug, zcoef, zav   ! local scalars
435      REAL(wp) ::   zdku,   zdkv, zsqen       !   -      -
436      REAL(wp) ::   zemxl, zemlm, zemlp, zmaxice       !   -      -
437      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zmxlm, zmxld   ! 3D workspace
438      !!--------------------------------------------------------------------
439      !
440      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
441      !                     !  Mixing length
442      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
443      !
444      !                     !* Buoyancy length scale: l=sqrt(2*e/n**2)
445      !
446      ! initialisation of interior minimum value (avoid a 2d loop with mikt)
447      zmxlm(:,:,:)  = rmxl_min   
448      zmxld(:,:,:)  = rmxl_min
449      !
450     IF( ln_mxl0 ) THEN            ! surface mixing length = F(stress) : l=vkarmn*2.e5*taum/(rho0*g)
451         !
452         zraug = vkarmn * 2.e5_wp / ( rho0 * grav )
453#if ! defined key_si3 && ! defined key_cice
454         DO_2D_00_00
455            zmxlm(ji,jj,1) =  zraug * taum(ji,jj) * tmask(ji,jj,1)
456         END_2D
457#else
458         SELECT CASE( nn_mxlice )             ! Type of scaling under sea-ice
459         !
460         CASE( 0 )                      ! No scaling under sea-ice
461            DO_2D_00_00
462               zmxlm(ji,jj,1) = zraug * taum(ji,jj) * tmask(ji,jj,1)
463            END_2D
464            !
465         CASE( 1 )                           ! scaling with constant sea-ice thickness
466            DO_2D_00_00
467               zmxlm(ji,jj,1) =  ( ( 1. - fr_i(ji,jj) ) * zraug * taum(ji,jj) + fr_i(ji,jj) * rn_mxlice ) * tmask(ji,jj,1)
468            END_2D
469            !
470         CASE( 2 )                                 ! scaling with mean sea-ice thickness
471            DO_2D_00_00
472#if defined key_si3
473               zmxlm(ji,jj,1) = ( ( 1. - fr_i(ji,jj) ) * zraug * taum(ji,jj) + fr_i(ji,jj) * hm_i(ji,jj) * 2. ) * tmask(ji,jj,1)
474#elif defined key_cice
475               zmaxice = MAXVAL( h_i(ji,jj,:) )
476               zmxlm(ji,jj,1) = ( ( 1. - fr_i(ji,jj) ) * zraug * taum(ji,jj) + fr_i(ji,jj) * zmaxice ) * tmask(ji,jj,1)
477#endif
478            END_2D
479            !
480         CASE( 3 )                                 ! scaling with max sea-ice thickness
481            DO_2D_00_00
482               zmaxice = MAXVAL( h_i(ji,jj,:) )
483               zmxlm(ji,jj,1) = ( ( 1. - fr_i(ji,jj) ) * zraug * taum(ji,jj) + fr_i(ji,jj) * zmaxice ) * tmask(ji,jj,1)
484            END_2D
485            !
486         END SELECT
487#endif
488         !
489         DO_2D_00_00
490            zmxlm(ji,jj,1) = MAX( rn_mxl0, zmxlm(ji,jj,1) )
491         END_2D
492         !
493      ELSE
494         zmxlm(:,:,1) = rn_mxl0
495      ENDIF
496
497      !
498      DO_3D_00_00( 2, jpkm1 )
499         zrn2 = MAX( rn2(ji,jj,jk), rsmall )
500         zmxlm(ji,jj,jk) = MAX(  rmxl_min,  SQRT( 2._wp * en(ji,jj,jk) / zrn2 )  )
501      END_3D
502      !
503      !                     !* Physical limits for the mixing length
504      !
505      zmxld(:,:, 1 ) = zmxlm(:,:,1)   ! surface set to the minimum value
506      zmxld(:,:,jpk) = rmxl_min       ! last level  set to the minimum value
507      !
508      SELECT CASE ( nn_mxl )
509      !
510 !!gm Not sure of that coding for ISF....
511      ! where wmask = 0 set zmxlm == e3w(:,:,:,Kmm)
512      CASE ( 0 )           ! bounded by the distance to surface and bottom
513         DO_3D_00_00( 2, jpkm1 )
514            zemxl = MIN( gdepw(ji,jj,jk,Kmm) - gdepw(ji,jj,mikt(ji,jj),Kmm), zmxlm(ji,jj,jk),   &
515            &            gdepw(ji,jj,mbkt(ji,jj)+1,Kmm) - gdepw(ji,jj,jk,Kmm) )
516            ! wmask prevent zmxlm = 0 if jk = mikt(ji,jj)
517            zmxlm(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN( zmxlm(ji,jj,jk) , e3w(ji,jj,jk,Kmm) ) * (1 - wmask(ji,jj,jk))
518            zmxld(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN( zmxlm(ji,jj,jk) , e3w(ji,jj,jk,Kmm) ) * (1 - wmask(ji,jj,jk))
519         END_3D
520         !
521      CASE ( 1 )           ! bounded by the vertical scale factor
522         DO_3D_00_00( 2, jpkm1 )
523            zemxl = MIN( e3w(ji,jj,jk,Kmm), zmxlm(ji,jj,jk) )
524            zmxlm(ji,jj,jk) = zemxl
525            zmxld(ji,jj,jk) = zemxl
526         END_3D
527         !
528      CASE ( 2 )           ! |dk[xml]| bounded by e3t :
529         DO_3D_00_00( 2, jpkm1 )
530            zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + e3t(ji,jj,jk-1,Kmm), zmxlm(ji,jj,jk) )
531         END_3D
532         DO_3DS_00_00( jpkm1, 2, -1 )
533            zemxl = MIN( zmxlm(ji,jj,jk+1) + e3t(ji,jj,jk+1,Kmm), zmxlm(ji,jj,jk) )
534            zmxlm(ji,jj,jk) = zemxl
535            zmxld(ji,jj,jk) = zemxl
536         END_3D
537         !
538      CASE ( 3 )           ! lup and ldown, |dk[xml]| bounded by e3t :
539         DO_3D_00_00( 2, jpkm1 )
540            zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + e3t(ji,jj,jk-1,Kmm), zmxlm(ji,jj,jk) )
541         END_3D
542         DO_3DS_00_00( jpkm1, 2, -1 )
543            zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + e3t(ji,jj,jk+1,Kmm), zmxlm(ji,jj,jk) )
544         END_3D
545         DO_3D_00_00( 2, jpkm1 )
546            zemlm = MIN ( zmxld(ji,jj,jk),  zmxlm(ji,jj,jk) )
547            zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) )
548            zmxlm(ji,jj,jk) = zemlm
549            zmxld(ji,jj,jk) = zemlp
550         END_3D
551         !
552      END SELECT
553      !
554      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
555      !                     !  Vertical eddy viscosity and diffusivity  (avm and avt)
556      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
557      DO_3D_00_00( 1, jpkm1 )
558         zsqen = SQRT( en(ji,jj,jk) )
559         zav   = rn_ediff * zmxlm(ji,jj,jk) * zsqen
560         p_avm(ji,jj,jk) = MAX( zav,                  avmb(jk) ) * wmask(ji,jj,jk)
561         p_avt(ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
562         dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk)
563      END_3D
564      !
565      !
566      IF( nn_pdl == 1 ) THEN      !* Prandtl number case: update avt
567         DO_3D_00_00( 2, jpkm1 )
568            p_avt(ji,jj,jk)   = MAX( apdlr(ji,jj,jk) * p_avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
569         END_3D
570      ENDIF
571      !
572      IF(sn_cfctl%l_prtctl) THEN
573         CALL prt_ctl( tab3d_1=en   , clinfo1=' tke  - e: ', tab3d_2=p_avt, clinfo2=' t: ', kdim=jpk)
574         CALL prt_ctl( tab3d_1=p_avm, clinfo1=' tke  - m: ', kdim=jpk )
575      ENDIF
576      !
577   END SUBROUTINE tke_avn
578
579
580   SUBROUTINE zdf_tke_init( Kmm )
581      !!----------------------------------------------------------------------
582      !!                  ***  ROUTINE zdf_tke_init  ***
583      !!                     
584      !! ** Purpose :   Initialization of the vertical eddy diffivity and
585      !!              viscosity when using a tke turbulent closure scheme
586      !!
587      !! ** Method  :   Read the namzdf_tke namelist and check the parameters
588      !!              called at the first timestep (nit000)
589      !!
590      !! ** input   :   Namlist namzdf_tke
591      !!
592      !! ** Action  :   Increase by 1 the nstop flag is setting problem encounter
593      !!----------------------------------------------------------------------
594      USE zdf_oce , ONLY : ln_zdfiwm   ! Internal Wave Mixing flag
595      !!
596      INTEGER, INTENT(in) ::   Kmm          ! time level index
597      INTEGER             ::   ji, jj, jk   ! dummy loop indices
598      INTEGER             ::   ios
599      !!
600      NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb   , rn_emin  ,  &
601         &                 rn_emin0, rn_bshear, nn_mxl   , ln_mxl0  ,  &
602         &                 rn_mxl0 , nn_mxlice, rn_mxlice,             &
603         &                 nn_pdl  , ln_drg   , ln_lc    , rn_lc,      &
604         &                 nn_etau , nn_htau  , rn_efr   , rn_eice 
605      !!----------------------------------------------------------------------
606      !
607      READ  ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901)
608901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in reference namelist' )
609
610      READ  ( numnam_cfg, namzdf_tke, IOSTAT = ios, ERR = 902 )
611902   IF( ios >  0 ) CALL ctl_nam ( ios , 'namzdf_tke in configuration namelist' )
612      IF(lwm) WRITE ( numond, namzdf_tke )
613      !
614      ri_cri   = 2._wp    / ( 2._wp + rn_ediss / rn_ediff )   ! resulting critical Richardson number
615      !
616      IF(lwp) THEN                    !* Control print
617         WRITE(numout,*)
618         WRITE(numout,*) 'zdf_tke_init : tke turbulent closure scheme - initialisation'
619         WRITE(numout,*) '~~~~~~~~~~~~'
620         WRITE(numout,*) '   Namelist namzdf_tke : set tke mixing parameters'
621         WRITE(numout,*) '      coef. to compute avt                        rn_ediff  = ', rn_ediff
622         WRITE(numout,*) '      Kolmogoroff dissipation coef.               rn_ediss  = ', rn_ediss
623         WRITE(numout,*) '      tke surface input coef.                     rn_ebb    = ', rn_ebb
624         WRITE(numout,*) '      minimum value of tke                        rn_emin   = ', rn_emin
625         WRITE(numout,*) '      surface minimum value of tke                rn_emin0  = ', rn_emin0
626         WRITE(numout,*) '      prandl number flag                          nn_pdl    = ', nn_pdl
627         WRITE(numout,*) '      background shear (>0)                       rn_bshear = ', rn_bshear
628         WRITE(numout,*) '      mixing length type                          nn_mxl    = ', nn_mxl
629         WRITE(numout,*) '         surface mixing length = F(stress) or not    ln_mxl0   = ', ln_mxl0
630         IF( ln_mxl0 ) THEN
631            WRITE(numout,*) '      type of scaling under sea-ice               nn_mxlice = ', nn_mxlice
632            IF( nn_mxlice == 1 ) &
633            WRITE(numout,*) '      ice thickness when scaling under sea-ice    rn_mxlice = ', rn_mxlice
634         ENDIF         
635         WRITE(numout,*) '         surface  mixing length minimum value        rn_mxl0   = ', rn_mxl0
636         WRITE(numout,*) '      top/bottom friction forcing flag            ln_drg    = ', ln_drg
637         WRITE(numout,*) '      Langmuir cells parametrization              ln_lc     = ', ln_lc
638         WRITE(numout,*) '         coef to compute vertical velocity of LC     rn_lc  = ', rn_lc
639         WRITE(numout,*) '      test param. to add tke induced by wind      nn_etau   = ', nn_etau
640         WRITE(numout,*) '          type of tke penetration profile            nn_htau   = ', nn_htau
641         WRITE(numout,*) '          fraction of TKE that penetrates            rn_efr    = ', rn_efr
642         WRITE(numout,*) '          below sea-ice:  =0 ON                      rn_eice   = ', rn_eice
643         WRITE(numout,*) '          =4 OFF when ice fraction > 1/4   '
644         IF( ln_drg ) THEN
645            WRITE(numout,*)
646            WRITE(numout,*) '   Namelist namdrg_top/_bot:   used values:'
647            WRITE(numout,*) '      top    ocean cavity roughness (m)          rn_z0(_top)= ', r_z0_top
648            WRITE(numout,*) '      Bottom seafloor     roughness (m)          rn_z0(_bot)= ', r_z0_bot
649         ENDIF
650         WRITE(numout,*)
651         WRITE(numout,*) '   ==>>>   critical Richardson nb with your parameters  ri_cri = ', ri_cri
652         WRITE(numout,*)
653      ENDIF
654      !
655      IF( ln_zdfiwm ) THEN          ! Internal wave-driven mixing
656         rn_emin  = 1.e-10_wp             ! specific values of rn_emin & rmxl_min are used
657         rmxl_min = 1.e-03_wp             ! associated avt minimum = molecular salt diffusivity (10^-9 m2/s)
658         IF(lwp) WRITE(numout,*) '   ==>>>   Internal wave-driven mixing case:   force   rn_emin = 1.e-10 and rmxl_min = 1.e-3'
659      ELSE                          ! standard case : associated avt minimum = molecular viscosity (10^-6 m2/s)
660         rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) )    ! resulting minimum length to recover molecular viscosity
661         IF(lwp) WRITE(numout,*) '   ==>>>   minimum mixing length with your parameters rmxl_min = ', rmxl_min
662      ENDIF
663      !
664      !                              ! allocate tke arrays
665      IF( zdf_tke_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_init : unable to allocate arrays' )
666      !
667      !                               !* Check of some namelist values
668      IF( nn_mxl  < 0   .OR.  nn_mxl  > 3 )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' )
669      IF( nn_pdl  < 0   .OR.  nn_pdl  > 1 )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' )
670      IF( nn_htau < 0   .OR.  nn_htau > 1 )   CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2 ' )
671      IF( nn_etau == 3 .AND. .NOT. ln_cpl )   CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' )
672      !
673      IF( ln_mxl0 ) THEN
674         IF(lwp) WRITE(numout,*)
675         IF(lwp) WRITE(numout,*) '   ==>>>   use a surface mixing length = F(stress) :   set rn_mxl0 = rmxl_min'
676         rn_mxl0 = rmxl_min
677      ENDIF
678     
679      IF( nn_etau == 2  )   CALL zdf_mxl( nit000, Kmm )      ! Initialization of nmln
680
681      !                               !* depth of penetration of surface tke
682      IF( nn_etau /= 0 ) THEN     
683         SELECT CASE( nn_htau )             ! Choice of the depth of penetration
684         CASE( 0 )                                 ! constant depth penetration (here 10 meters)
685            htau(:,:) = 10._wp
686         CASE( 1 )                                 ! F(latitude) : 0.5m to 30m poleward of 40 degrees
687            htau(:,:) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )           
688         END SELECT
689      ENDIF
690      !                                !* read or initialize all required files
691      CALL tke_rst( nit000, 'READ' )      ! (en, avt_k, avm_k, dissl)
692      !
693      IF( lwxios ) THEN
694         CALL iom_set_rstw_var_active('en')
695         CALL iom_set_rstw_var_active('avt_k')
696         CALL iom_set_rstw_var_active('avm_k')
697         CALL iom_set_rstw_var_active('dissl')
698      ENDIF
699   END SUBROUTINE zdf_tke_init
700
701
702   SUBROUTINE tke_rst( kt, cdrw )
703      !!---------------------------------------------------------------------
704      !!                   ***  ROUTINE tke_rst  ***
705      !!                     
706      !! ** Purpose :   Read or write TKE file (en) in restart file
707      !!
708      !! ** Method  :   use of IOM library
709      !!                if the restart does not contain TKE, en is either
710      !!                set to rn_emin or recomputed
711      !!----------------------------------------------------------------------
712      USE zdf_oce , ONLY : en, avt_k, avm_k   ! ocean vertical physics
713      !!
714      INTEGER         , INTENT(in) ::   kt     ! ocean time-step
715      CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
716      !
717      INTEGER ::   jit, jk              ! dummy loop indices
718      INTEGER ::   id1, id2, id3, id4   ! local integers
719      !!----------------------------------------------------------------------
720      !
721      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
722         !                                   ! ---------------
723         IF( ln_rstart ) THEN                   !* Read the restart file
724            id1 = iom_varid( numror, 'en'   , ldstop = .FALSE. )
725            id2 = iom_varid( numror, 'avt_k', ldstop = .FALSE. )
726            id3 = iom_varid( numror, 'avm_k', ldstop = .FALSE. )
727            id4 = iom_varid( numror, 'dissl', ldstop = .FALSE. )
728            !
729            IF( MIN( id1, id2, id3, id4 ) > 0 ) THEN      ! fields exist
730               CALL iom_get( numror, jpdom_autoglo, 'en'   , en   , ldxios = lrxios )
731               CALL iom_get( numror, jpdom_autoglo, 'avt_k', avt_k, ldxios = lrxios )
732               CALL iom_get( numror, jpdom_autoglo, 'avm_k', avm_k, ldxios = lrxios )
733               CALL iom_get( numror, jpdom_autoglo, 'dissl', dissl, ldxios = lrxios )
734            ELSE                                          ! start TKE from rest
735               IF(lwp) WRITE(numout,*)
736               IF(lwp) WRITE(numout,*) '   ==>>>   previous run without TKE scheme, set en to background values'
737               en   (:,:,:) = rn_emin * wmask(:,:,:)
738               dissl(:,:,:) = 1.e-12_wp
739               ! avt_k, avm_k already set to the background value in zdf_phy_init
740            ENDIF
741         ELSE                                   !* Start from rest
742            IF(lwp) WRITE(numout,*)
743            IF(lwp) WRITE(numout,*) '   ==>>>   start from rest: set en to the background value'
744            en   (:,:,:) = rn_emin * wmask(:,:,:)
745            dissl(:,:,:) = 1.e-12_wp
746            ! avt_k, avm_k already set to the background value in zdf_phy_init
747         ENDIF
748         !
749      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
750         !                                   ! -------------------
751         IF(lwp) WRITE(numout,*) '---- tke_rst ----'
752         IF( lwxios ) CALL iom_swap(      cwxios_context          ) 
753         CALL iom_rstput( kt, nitrst, numrow, 'en'   , en   , ldxios = lwxios )
754         CALL iom_rstput( kt, nitrst, numrow, 'avt_k', avt_k, ldxios = lwxios )
755         CALL iom_rstput( kt, nitrst, numrow, 'avm_k', avm_k, ldxios = lwxios )
756         CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl, ldxios = lwxios )
757         IF( lwxios ) CALL iom_swap(      cxios_context          )
758         !
759      ENDIF
760      !
761   END SUBROUTINE tke_rst
762
763   !!======================================================================
764END MODULE zdftke
Note: See TracBrowser for help on using the repository browser.