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/branches/MOI/NEMO_4.03_IODRAG/src/OCE/ZDF – NEMO

source: NEMO/branches/MOI/NEMO_4.03_IODRAG/src/OCE/ZDF/zdftke.F90 @ 12966

Last change on this file since 12966 was 12966, checked in by jchanut, 4 years ago

Switch to log-law bcs below ice in tke and gls - hard coded 3cm ice-ocean roughness: #2468

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