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/2019/dev_r11085_ASINTER-05_Brodeau_Advanced_Bulk/src/OCE/ZDF – NEMO

source: NEMO/branches/2019/dev_r11085_ASINTER-05_Brodeau_Advanced_Bulk/src/OCE/ZDF/zdftke.F90 @ 11831

Last change on this file since 11831 was 11831, checked in by laurent, 4 years ago

Update the branch to r11830 of the trunk!

  • Property svn:keywords set to Id
File size: 43.7 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, 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)     ::   zhlc, zfr_i
205      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zpelc, zdiag, zd_up, zd_lw
206      !!--------------------------------------------------------------------
207      !
208      zbbrau = rn_ebb / rau0       ! Local constant initialisation
209      zfact1 = -.5_wp * rdt 
210      zfact2 = 1.5_wp * rdt * rn_ediss
211      zfact3 = 0.5_wp       * rn_ediss
212      !
213      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
214      !                     !  Surface/top/bottom boundary condition on tke
215      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
216     
217      DO jj = 2, jpjm1            ! en(1)   = rn_ebb taum / rau0  (min value rn_emin0)
218         DO ji = fs_2, fs_jpim1   ! vector opt.
219            en(ji,jj,1) = MAX( rn_emin0, zbbrau * taum(ji,jj) ) * tmask(ji,jj,1)
220         END DO
221      END DO
222      IF ( ln_isfcav ) THEN
223         DO jj = 2, jpjm1            ! en(mikt(ji,jj))   = rn_emin
224            DO ji = fs_2, fs_jpim1   ! vector opt.
225               en(ji,jj,mikt(ji,jj)) = rn_emin * tmask(ji,jj,1)
226            END DO
227         END DO
228      ENDIF
229      !
230      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
231      !                     !  Bottom boundary condition on tke
232      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
233      !
234      !   en(bot)   = (ebb0/rau0)*0.5*sqrt(u_botfr^2+v_botfr^2) (min value rn_emin)
235      ! where ebb0 does not includes surface wave enhancement (i.e. ebb0=3.75)
236      ! Note that stress averaged is done using an wet-only calculation of u and v at t-point like in zdfsh2
237      !
238      IF( ln_drg ) THEN       !== friction used as top/bottom boundary condition on TKE
239         !
240         DO jj = 2, jpjm1           ! bottom friction
241            DO ji = fs_2, fs_jpim1     ! vector opt.
242               zmsku = ( 2. - umask(ji-1,jj,mbkt(ji,jj)) * umask(ji,jj,mbkt(ji,jj)) )
243               zmskv = ( 2. - vmask(ji,jj-1,mbkt(ji,jj)) * vmask(ji,jj,mbkt(ji,jj)) )
244               !                       ! where 0.001875 = (rn_ebb0/rau0) * 0.5 = 3.75*0.5/1000. (CAUTION CdU<0)
245               zebot = - 0.001875_wp * rCdU_bot(ji,jj) * SQRT(  ( zmsku*( ub(ji,jj,mbkt(ji,jj))+ub(ji-1,jj,mbkt(ji,jj)) ) )**2  &
246                  &                                           + ( zmskv*( vb(ji,jj,mbkt(ji,jj))+vb(ji,jj-1,mbkt(ji,jj)) ) )**2  )
247               en(ji,jj,mbkt(ji,jj)+1) = MAX( zebot, rn_emin ) * ssmask(ji,jj)
248            END DO
249         END DO
250         IF( ln_isfcav ) THEN       ! top friction
251            DO jj = 2, jpjm1
252               DO ji = fs_2, fs_jpim1   ! vector opt.
253                  zmsku = ( 2. - umask(ji-1,jj,mikt(ji,jj)) * umask(ji,jj,mikt(ji,jj)) )
254                  zmskv = ( 2. - vmask(ji,jj-1,mikt(ji,jj)) * vmask(ji,jj,mikt(ji,jj)) )
255                  !                             ! where 0.001875 = (rn_ebb0/rau0) * 0.5 = 3.75*0.5/1000.  (CAUTION CdU<0)
256                  zetop = - 0.001875_wp * rCdU_top(ji,jj) * SQRT(  ( zmsku*( ub(ji,jj,mikt(ji,jj))+ub(ji-1,jj,mikt(ji,jj)) ) )**2  &
257                     &                                           + ( zmskv*( vb(ji,jj,mikt(ji,jj))+vb(ji,jj-1,mikt(ji,jj)) ) )**2  )
258                  en(ji,jj,mikt(ji,jj)) = MAX( zetop, rn_emin ) * (1._wp - tmask(ji,jj,1))   ! masked at ocean surface
259               END DO
260            END DO
261         ENDIF
262         !
263      ENDIF
264      !
265      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
266      IF( ln_lc ) THEN      !  Langmuir circulation source term added to tke   !   (Axell JGR 2002)
267         !                  !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
268         !
269         !                        !* total energy produce by LC : cumulative sum over jk
270         zpelc(:,:,1) =  MAX( rn2b(:,:,1), 0._wp ) * pdepw(:,:,1) * p_e3w(:,:,1)
271         DO jk = 2, jpk
272            zpelc(:,:,jk)  = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0._wp ) * pdepw(:,:,jk) * p_e3w(:,:,jk)
273         END DO
274         !                        !* finite Langmuir Circulation depth
275         zcof = 0.5 * 0.016 * 0.016 / ( zrhoa * zcdrag )
276         imlc(:,:) = mbkt(:,:) + 1       ! Initialization to the number of w ocean point (=2 over land)
277         DO jk = jpkm1, 2, -1
278            DO jj = 1, jpj               ! Last w-level at which zpelc>=0.5*us*us
279               DO ji = 1, jpi            !      with us=0.016*wind(starting from jpk-1)
280                  zus  = zcof * taum(ji,jj)
281                  IF( zpelc(ji,jj,jk) > zus )   imlc(ji,jj) = jk
282               END DO
283            END DO
284         END DO
285         !                               ! finite LC depth
286         DO jj = 1, jpj 
287            DO ji = 1, jpi
288               zhlc(ji,jj) = pdepw(ji,jj,imlc(ji,jj))
289            END DO
290         END DO
291         zcof = 0.016 / SQRT( zrhoa * zcdrag )
292         DO jj = 2, jpjm1
293            DO ji = fs_2, fs_jpim1   ! vector opt.
294               zus  = zcof * SQRT( taum(ji,jj) )           ! Stokes drift
295               zfr_i(ji,jj) = ( 1._wp - 4._wp * fr_i(ji,jj) ) * zus * zus * zus * tmask(ji,jj,1) ! zus > 0. ok
296               IF (zfr_i(ji,jj) < 0. ) zfr_i(ji,jj) = 0.
297            END DO
298         END DO         
299         DO jk = 2, jpkm1         !* TKE Langmuir circulation source term added to en
300            DO jj = 2, jpjm1
301               DO ji = fs_2, fs_jpim1   ! vector opt.
302                  IF ( zfr_i(ji,jj) /= 0. ) THEN               
303                     ! vertical velocity due to LC   
304                     IF ( pdepw(ji,jj,jk) - zhlc(ji,jj) < 0 .AND. wmask(ji,jj,jk) /= 0. ) THEN
305                        !                                           ! vertical velocity due to LC
306                        zwlc = rn_lc * SIN( rpi * pdepw(ji,jj,jk) / zhlc(ji,jj) )   ! warning: optimization: zus^3 is in zfr_i
307                        !                                           ! TKE Langmuir circulation source term
308                        en(ji,jj,jk) = en(ji,jj,jk) + rdt * zfr_i(ji,jj) * ( zwlc * zwlc * zwlc ) / zhlc(ji,jj)
309                     ENDIF
310                  ENDIF
311               END DO
312            END DO
313         END DO
314         !
315      ENDIF
316      !
317      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
318      !                     !  Now Turbulent kinetic energy (output in en)
319      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
320      !                     ! Resolution of a tridiagonal linear system by a "methode de chasse"
321      !                     ! computation from level 2 to jpkm1  (e(1) already computed and e(jpk)=0 ).
322      !                     ! zdiag : diagonal zd_up : upper diagonal zd_lw : lower diagonal
323      !
324      IF( nn_pdl == 1 ) THEN      !* Prandtl number = F( Ri )
325         DO jk = 2, jpkm1
326            DO jj = 2, jpjm1
327               DO ji = 2, jpim1
328                  !                             ! local Richardson number
329                  zri = MAX( rn2b(ji,jj,jk), 0._wp ) * p_avm(ji,jj,jk) / ( p_sh2(ji,jj,jk) + rn_bshear )
330                  !                             ! inverse of Prandtl number
331                  apdlr(ji,jj,jk) = MAX(  0.1_wp,  ri_cri / MAX( ri_cri , zri )  )
332               END DO
333            END DO
334         END DO
335      ENDIF
336      !         
337      DO jk = 2, jpkm1           !* Matrix and right hand side in en
338         DO jj = 2, jpjm1
339            DO ji = fs_2, fs_jpim1   ! vector opt.
340               zcof   = zfact1 * tmask(ji,jj,jk)
341               !                                   ! A minimum of 2.e-5 m2/s is imposed on TKE vertical
342               !                                   ! eddy coefficient (ensure numerical stability)
343               zzd_up = zcof * MAX(  p_avm(ji,jj,jk+1) + p_avm(ji,jj,jk  ) , 2.e-5_wp  )   &  ! upper diagonal
344                  &          /    (  p_e3t(ji,jj,jk  ) * p_e3w(ji,jj,jk  )  )
345               zzd_lw = zcof * MAX(  p_avm(ji,jj,jk  ) + p_avm(ji,jj,jk-1) , 2.e-5_wp  )   &  ! lower diagonal
346                  &          /    (  p_e3t(ji,jj,jk-1) * p_e3w(ji,jj,jk  )  )
347               !
348               zd_up(ji,jj,jk) = zzd_up            ! Matrix (zdiag, zd_up, zd_lw)
349               zd_lw(ji,jj,jk) = zzd_lw
350               zdiag(ji,jj,jk) = 1._wp - zzd_lw - zzd_up + zfact2 * dissl(ji,jj,jk) * wmask(ji,jj,jk)
351               !
352               !                                   ! right hand side in en
353               en(ji,jj,jk) = en(ji,jj,jk) + rdt * (  p_sh2(ji,jj,jk)                          &   ! shear
354                  &                                 - p_avt(ji,jj,jk) * rn2(ji,jj,jk)          &   ! stratification
355                  &                                 + zfact3 * dissl(ji,jj,jk) * en(ji,jj,jk)  &   ! dissipation
356                  &                                ) * wmask(ji,jj,jk)
357            END DO
358         END DO
359      END DO
360      !                          !* Matrix inversion from level 2 (tke prescribed at level 1)
361      DO jk = 3, jpkm1                             ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
362         DO jj = 2, jpjm1
363            DO ji = fs_2, fs_jpim1    ! vector opt.
364               zdiag(ji,jj,jk) = zdiag(ji,jj,jk) - zd_lw(ji,jj,jk) * zd_up(ji,jj,jk-1) / zdiag(ji,jj,jk-1)
365            END DO
366         END DO
367      END DO
368      DO jj = 2, jpjm1                             ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1
369         DO ji = fs_2, fs_jpim1   ! vector opt.
370            zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1)    ! Surface boudary conditions on tke
371         END DO
372      END DO
373      DO jk = 3, jpkm1
374         DO jj = 2, jpjm1
375            DO ji = fs_2, fs_jpim1    ! vector opt.
376               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)
377            END DO
378         END DO
379      END DO
380      DO jj = 2, jpjm1                             ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
381         DO ji = fs_2, fs_jpim1   ! vector opt.
382            en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1)
383         END DO
384      END DO
385      DO jk = jpk-2, 2, -1
386         DO jj = 2, jpjm1
387            DO ji = fs_2, fs_jpim1    ! vector opt.
388               en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk)
389            END DO
390         END DO
391      END DO
392      DO jk = 2, jpkm1                             ! set the minimum value of tke
393         DO jj = 2, jpjm1
394            DO ji = fs_2, fs_jpim1   ! vector opt.
395               en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * wmask(ji,jj,jk)
396            END DO
397         END DO
398      END DO
399      !
400      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
401      !                            !  TKE due to surface and internal wave breaking
402      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
403!!gm BUG : in the exp  remove the depth of ssh !!!
404!!gm       i.e. use gde3w in argument (pdepw)
405     
406     
407      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction)
408         DO jk = 2, jpkm1                       ! rn_eice =0 ON below sea-ice, =4 OFF when ice fraction > 0.25
409            DO jj = 2, jpjm1
410               DO ji = fs_2, fs_jpim1   ! vector opt.
411                  en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -pdepw(ji,jj,jk) / htau(ji,jj) )   &
412                     &                                 * MAX(0.,1._wp - rn_eice *fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
413               END DO
414            END DO
415         END DO
416      ELSEIF( nn_etau == 2 ) THEN       !* act only at the base of the mixed layer (jk=nmln)  (rn_efr fraction)
417         DO jj = 2, jpjm1
418            DO ji = fs_2, fs_jpim1   ! vector opt.
419               jk = nmln(ji,jj)
420               en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -pdepw(ji,jj,jk) / htau(ji,jj) )   &
421                  &                                 * MAX(0.,1._wp - rn_eice *fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
422            END DO
423         END DO
424      ELSEIF( nn_etau == 3 ) THEN       !* penetration belox the mixed layer (HF variability)
425         DO jk = 2, jpkm1
426            DO jj = 2, jpjm1
427               DO ji = fs_2, fs_jpim1   ! vector opt.
428                  ztx2 = utau(ji-1,jj  ) + utau(ji,jj)
429                  zty2 = vtau(ji  ,jj-1) + vtau(ji,jj)
430                  ztau = 0.5_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1)    ! module of the mean stress
431                  zdif = taum(ji,jj) - ztau                            ! mean of modulus - modulus of the mean
432                  zdif = rhftau_scl * MAX( 0._wp, zdif + rhftau_add )  ! apply some modifications...
433                  en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -pdepw(ji,jj,jk) / htau(ji,jj) )   &
434                     &                        * MAX(0.,1._wp - rn_eice *fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1)
435               END DO
436            END DO
437         END DO
438      ENDIF
439      !
440   END SUBROUTINE tke_tke
441
442
443   SUBROUTINE tke_avn( pdepw, p_e3t, p_e3w, p_avm, p_avt )
444      !!----------------------------------------------------------------------
445      !!                   ***  ROUTINE tke_avn  ***
446      !!
447      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
448      !!
449      !! ** Method  :   At this stage, en, the now TKE, is known (computed in
450      !!              the tke_tke routine). First, the now mixing lenth is
451      !!      computed from en and the strafification (N^2), then the mixings
452      !!      coefficients are computed.
453      !!              - Mixing length : a first evaluation of the mixing lengh
454      !!      scales is:
455      !!                      mxl = sqrt(2*en) / N 
456      !!      where N is the brunt-vaisala frequency, with a minimum value set
457      !!      to rmxl_min (rn_mxl0) in the interior (surface) ocean.
458      !!        The mixing and dissipative length scale are bound as follow :
459      !!         nn_mxl=0 : mxl bounded by the distance to surface and bottom.
460      !!                        zmxld = zmxlm = mxl
461      !!         nn_mxl=1 : mxl bounded by the e3w and zmxld = zmxlm = mxl
462      !!         nn_mxl=2 : mxl bounded such that the vertical derivative of mxl is
463      !!                    less than 1 (|d/dz(mxl)|<1) and zmxld = zmxlm = mxl
464      !!         nn_mxl=3 : mxl is bounded from the surface to the bottom usings
465      !!                    |d/dz(xml)|<1 to obtain lup, and from the bottom to
466      !!                    the surface to obtain ldown. the resulting length
467      !!                    scales are:
468      !!                        zmxld = sqrt( lup * ldown )
469      !!                        zmxlm = min ( lup , ldown )
470      !!              - Vertical eddy viscosity and diffusivity:
471      !!                      avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
472      !!                      avt = max( avmb, pdlr * avm ) 
473      !!      with pdlr=1 if nn_pdl=0, pdlr=1/pdl=F(Ri) otherwise.
474      !!
475      !! ** Action  : - avt, avm : now vertical eddy diffusivity and viscosity (w-point)
476      !!----------------------------------------------------------------------
477      USE zdf_oce , ONLY : en, avtb, avmb, avtb_2d   ! ocean vertical physics
478      !!
479      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pdepw          ! depth (w-points)
480      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   p_e3t, p_e3w   ! level thickness (t- & w-points)
481      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_avm, p_avt   ! vertical eddy viscosity & diffusivity (w-points)
482      !
483      INTEGER  ::   ji, jj, jk   ! dummy loop indices
484      REAL(wp) ::   zrn2, zraug, zcoef, zav   ! local scalars
485      REAL(wp) ::   zdku,   zdkv, zsqen       !   -      -
486      REAL(wp) ::   zemxl, zemlm, zemlp       !   -      -
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      !
511      DO jk = 2, jpkm1              ! interior value : l=sqrt(2*e/n^2)
512         DO jj = 2, jpjm1
513            DO ji = fs_2, fs_jpim1   ! vector opt.
514               zrn2 = MAX( rn2(ji,jj,jk), rsmall )
515               zmxlm(ji,jj,jk) = MAX(  rmxl_min,  SQRT( 2._wp * en(ji,jj,jk) / zrn2 )  )
516            END DO
517         END DO
518      END DO
519      !
520      !                     !* Physical limits for the mixing length
521      !
522      zmxld(:,:, 1 ) = zmxlm(:,:,1)   ! surface set to the minimum value
523      zmxld(:,:,jpk) = rmxl_min       ! last level  set to the minimum value
524      !
525      SELECT CASE ( nn_mxl )
526      !
527 !!gm Not sure of that coding for ISF....
528      ! where wmask = 0 set zmxlm == p_e3w
529      CASE ( 0 )           ! bounded by the distance to surface and bottom
530         DO jk = 2, jpkm1
531            DO jj = 2, jpjm1
532               DO ji = fs_2, fs_jpim1   ! vector opt.
533                  zemxl = MIN( pdepw(ji,jj,jk) - pdepw(ji,jj,mikt(ji,jj)), zmxlm(ji,jj,jk),   &
534                  &            pdepw(ji,jj,mbkt(ji,jj)+1) - pdepw(ji,jj,jk) )
535                  ! wmask prevent zmxlm = 0 if jk = mikt(ji,jj)
536                  zmxlm(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN( zmxlm(ji,jj,jk) , p_e3w(ji,jj,jk) ) * (1 - wmask(ji,jj,jk))
537                  zmxld(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN( zmxlm(ji,jj,jk) , p_e3w(ji,jj,jk) ) * (1 - wmask(ji,jj,jk))
538               END DO
539            END DO
540         END DO
541         !
542      CASE ( 1 )           ! bounded by the vertical scale factor
543         DO jk = 2, jpkm1
544            DO jj = 2, jpjm1
545               DO ji = fs_2, fs_jpim1   ! vector opt.
546                  zemxl = MIN( p_e3w(ji,jj,jk), zmxlm(ji,jj,jk) )
547                  zmxlm(ji,jj,jk) = zemxl
548                  zmxld(ji,jj,jk) = zemxl
549               END DO
550            END DO
551         END DO
552         !
553      CASE ( 2 )           ! |dk[xml]| bounded by e3t :
554         DO jk = 2, jpkm1         ! from the surface to the bottom :
555            DO jj = 2, jpjm1
556               DO ji = fs_2, fs_jpim1   ! vector opt.
557                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + p_e3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
558               END DO
559            END DO
560         END DO
561         DO jk = jpkm1, 2, -1     ! from the bottom to the surface :
562            DO jj = 2, jpjm1
563               DO ji = fs_2, fs_jpim1   ! vector opt.
564                  zemxl = MIN( zmxlm(ji,jj,jk+1) + p_e3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
565                  zmxlm(ji,jj,jk) = zemxl
566                  zmxld(ji,jj,jk) = zemxl
567               END DO
568            END DO
569         END DO
570         !
571      CASE ( 3 )           ! lup and ldown, |dk[xml]| bounded by e3t :
572         DO jk = 2, jpkm1         ! from the surface to the bottom : lup
573            DO jj = 2, jpjm1
574               DO ji = fs_2, fs_jpim1   ! vector opt.
575                  zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + p_e3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
576               END DO
577            END DO
578         END DO
579         DO jk = jpkm1, 2, -1     ! from the bottom to the surface : ldown
580            DO jj = 2, jpjm1
581               DO ji = fs_2, fs_jpim1   ! vector opt.
582                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + p_e3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
583               END DO
584            END DO
585         END DO
586         DO jk = 2, jpkm1
587            DO jj = 2, jpjm1
588               DO ji = fs_2, fs_jpim1   ! vector opt.
589                  zemlm = MIN ( zmxld(ji,jj,jk),  zmxlm(ji,jj,jk) )
590                  zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) )
591                  zmxlm(ji,jj,jk) = zemlm
592                  zmxld(ji,jj,jk) = zemlp
593               END DO
594            END DO
595         END DO
596         !
597      END SELECT
598      !
599      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
600      !                     !  Vertical eddy viscosity and diffusivity  (avm and avt)
601      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
602      DO jk = 1, jpkm1            !* vertical eddy viscosity & diffivity at w-points
603         DO jj = 2, jpjm1
604            DO ji = fs_2, fs_jpim1   ! vector opt.
605               zsqen = SQRT( en(ji,jj,jk) )
606               zav   = rn_ediff * zmxlm(ji,jj,jk) * zsqen
607               p_avm(ji,jj,jk) = MAX( zav,                  avmb(jk) ) * wmask(ji,jj,jk)
608               p_avt(ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
609               dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk)
610            END DO
611         END DO
612      END DO
613      !
614      !
615      IF( nn_pdl == 1 ) THEN      !* Prandtl number case: update avt
616         DO jk = 2, jpkm1
617            DO jj = 2, jpjm1
618               DO ji = fs_2, fs_jpim1   ! vector opt.
619                  p_avt(ji,jj,jk)   = MAX( apdlr(ji,jj,jk) * p_avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * tmask(ji,jj,jk)
620              END DO
621            END DO
622         END DO
623      ENDIF
624      !
625      IF(ln_ctl) THEN
626         CALL prt_ctl( tab3d_1=en   , clinfo1=' tke  - e: ', tab3d_2=p_avt, clinfo2=' t: ', kdim=jpk)
627         CALL prt_ctl( tab3d_1=p_avm, clinfo1=' tke  - m: ', kdim=jpk )
628      ENDIF
629      !
630   END SUBROUTINE tke_avn
631
632
633   SUBROUTINE zdf_tke_init
634      !!----------------------------------------------------------------------
635      !!                  ***  ROUTINE zdf_tke_init  ***
636      !!                     
637      !! ** Purpose :   Initialization of the vertical eddy diffivity and
638      !!              viscosity when using a tke turbulent closure scheme
639      !!
640      !! ** Method  :   Read the namzdf_tke namelist and check the parameters
641      !!              called at the first timestep (nit000)
642      !!
643      !! ** input   :   Namlist namzdf_tke
644      !!
645      !! ** Action  :   Increase by 1 the nstop flag is setting problem encounter
646      !!----------------------------------------------------------------------
647      USE zdf_oce , ONLY : ln_zdfiwm   ! Internal Wave Mixing flag
648      !!
649      INTEGER ::   ji, jj, jk   ! dummy loop indices
650      INTEGER ::   ios
651      !!
652      NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb , rn_emin  ,          &
653         &                 rn_emin0, rn_bshear, nn_mxl , ln_mxl0  ,          &
654         &                 rn_mxl0 , nn_pdl   , ln_drg , ln_lc    , rn_lc,   &
655         &                 nn_etau , nn_htau  , rn_efr , rn_eice 
656      !!----------------------------------------------------------------------
657      !
658      REWIND( numnam_ref )              ! Namelist namzdf_tke in reference namelist : Turbulent Kinetic Energy
659      READ  ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901)
660901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in reference namelist' )
661
662      REWIND( numnam_cfg )              ! Namelist namzdf_tke in configuration namelist : Turbulent Kinetic Energy
663      READ  ( numnam_cfg, namzdf_tke, IOSTAT = ios, ERR = 902 )
664902   IF( ios >  0 ) CALL ctl_nam ( ios , 'namzdf_tke in configuration namelist' )
665      IF(lwm) WRITE ( numond, namzdf_tke )
666      !
667      ri_cri   = 2._wp    / ( 2._wp + rn_ediss / rn_ediff )   ! resulting critical Richardson number
668      !
669      IF(lwp) THEN                    !* Control print
670         WRITE(numout,*)
671         WRITE(numout,*) 'zdf_tke_init : tke turbulent closure scheme - initialisation'
672         WRITE(numout,*) '~~~~~~~~~~~~'
673         WRITE(numout,*) '   Namelist namzdf_tke : set tke mixing parameters'
674         WRITE(numout,*) '      coef. to compute avt                        rn_ediff  = ', rn_ediff
675         WRITE(numout,*) '      Kolmogoroff dissipation coef.               rn_ediss  = ', rn_ediss
676         WRITE(numout,*) '      tke surface input coef.                     rn_ebb    = ', rn_ebb
677         WRITE(numout,*) '      minimum value of tke                        rn_emin   = ', rn_emin
678         WRITE(numout,*) '      surface minimum value of tke                rn_emin0  = ', rn_emin0
679         WRITE(numout,*) '      prandl number flag                          nn_pdl    = ', nn_pdl
680         WRITE(numout,*) '      background shear (>0)                       rn_bshear = ', rn_bshear
681         WRITE(numout,*) '      mixing length type                          nn_mxl    = ', nn_mxl
682         WRITE(numout,*) '         surface mixing length = F(stress) or not    ln_mxl0   = ', ln_mxl0
683         WRITE(numout,*) '         surface  mixing length minimum value        rn_mxl0   = ', rn_mxl0
684         WRITE(numout,*) '      top/bottom friction forcing flag            ln_drg    = ', ln_drg
685         WRITE(numout,*) '      Langmuir cells parametrization              ln_lc     = ', ln_lc
686         WRITE(numout,*) '         coef to compute vertical velocity of LC     rn_lc  = ', rn_lc
687         WRITE(numout,*) '      test param. to add tke induced by wind      nn_etau   = ', nn_etau
688         WRITE(numout,*) '          type of tke penetration profile            nn_htau   = ', nn_htau
689         WRITE(numout,*) '          fraction of TKE that penetrates            rn_efr    = ', rn_efr
690         WRITE(numout,*) '          below sea-ice:  =0 ON                      rn_eice   = ', rn_eice
691         WRITE(numout,*) '          =4 OFF when ice fraction > 1/4   '
692         IF( ln_drg ) THEN
693            WRITE(numout,*)
694            WRITE(numout,*) '   Namelist namdrg_top/_bot:   used values:'
695            WRITE(numout,*) '      top    ocean cavity roughness (m)          rn_z0(_top)= ', r_z0_top
696            WRITE(numout,*) '      Bottom seafloor     roughness (m)          rn_z0(_bot)= ', r_z0_bot
697         ENDIF
698         WRITE(numout,*)
699         WRITE(numout,*) '   ==>>>   critical Richardson nb with your parameters  ri_cri = ', ri_cri
700         WRITE(numout,*)
701      ENDIF
702      !
703      IF( ln_zdfiwm ) THEN          ! Internal wave-driven mixing
704         rn_emin  = 1.e-10_wp             ! specific values of rn_emin & rmxl_min are used
705         rmxl_min = 1.e-03_wp             ! associated avt minimum = molecular salt diffusivity (10^-9 m2/s)
706         IF(lwp) WRITE(numout,*) '   ==>>>   Internal wave-driven mixing case:   force   rn_emin = 1.e-10 and rmxl_min = 1.e-3'
707      ELSE                          ! standard case : associated avt minimum = molecular viscosity (10^-6 m2/s)
708         rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) )    ! resulting minimum length to recover molecular viscosity
709         IF(lwp) WRITE(numout,*) '   ==>>>   minimum mixing length with your parameters rmxl_min = ', rmxl_min
710      ENDIF
711      !
712      !                              ! allocate tke arrays
713      IF( zdf_tke_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_init : unable to allocate arrays' )
714      !
715      !                               !* Check of some namelist values
716      IF( nn_mxl  < 0   .OR.  nn_mxl  > 3 )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' )
717      IF( nn_pdl  < 0   .OR.  nn_pdl  > 1 )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' )
718      IF( nn_htau < 0   .OR.  nn_htau > 1 )   CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2 ' )
719      IF( nn_etau == 3 .AND. .NOT. ln_cpl )   CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' )
720      !
721      IF( ln_mxl0 ) THEN
722         IF(lwp) WRITE(numout,*)
723         IF(lwp) WRITE(numout,*) '   ==>>>   use a surface mixing length = F(stress) :   set rn_mxl0 = rmxl_min'
724         rn_mxl0 = rmxl_min
725      ENDIF
726     
727      IF( nn_etau == 2  )   CALL zdf_mxl( nit000 )      ! Initialization of nmln
728
729      !                               !* depth of penetration of surface tke
730      IF( nn_etau /= 0 ) THEN     
731         SELECT CASE( nn_htau )             ! Choice of the depth of penetration
732         CASE( 0 )                                 ! constant depth penetration (here 10 meters)
733            htau(:,:) = 10._wp
734         CASE( 1 )                                 ! F(latitude) : 0.5m to 30m poleward of 40 degrees
735            htau(:,:) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )           
736         END SELECT
737      ENDIF
738      !                                !* read or initialize all required files
739      CALL tke_rst( nit000, 'READ' )      ! (en, avt_k, avm_k, dissl)
740      !
741      IF( lwxios ) THEN
742         CALL iom_set_rstw_var_active('en')
743         CALL iom_set_rstw_var_active('avt_k')
744         CALL iom_set_rstw_var_active('avm_k')
745         CALL iom_set_rstw_var_active('dissl')
746      ENDIF
747   END SUBROUTINE zdf_tke_init
748
749
750   SUBROUTINE tke_rst( kt, cdrw )
751      !!---------------------------------------------------------------------
752      !!                   ***  ROUTINE tke_rst  ***
753      !!                     
754      !! ** Purpose :   Read or write TKE file (en) in restart file
755      !!
756      !! ** Method  :   use of IOM library
757      !!                if the restart does not contain TKE, en is either
758      !!                set to rn_emin or recomputed
759      !!----------------------------------------------------------------------
760      USE zdf_oce , ONLY : en, avt_k, avm_k   ! ocean vertical physics
761      !!
762      INTEGER         , INTENT(in) ::   kt     ! ocean time-step
763      CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
764      !
765      INTEGER ::   jit, jk              ! dummy loop indices
766      INTEGER ::   id1, id2, id3, id4   ! local integers
767      !!----------------------------------------------------------------------
768      !
769      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
770         !                                   ! ---------------
771         IF( ln_rstart ) THEN                   !* Read the restart file
772            id1 = iom_varid( numror, 'en'   , ldstop = .FALSE. )
773            id2 = iom_varid( numror, 'avt_k', ldstop = .FALSE. )
774            id3 = iom_varid( numror, 'avm_k', ldstop = .FALSE. )
775            id4 = iom_varid( numror, 'dissl', ldstop = .FALSE. )
776            !
777            IF( MIN( id1, id2, id3, id4 ) > 0 ) THEN      ! fields exist
778               CALL iom_get( numror, jpdom_autoglo, 'en'   , en   , ldxios = lrxios )
779               CALL iom_get( numror, jpdom_autoglo, 'avt_k', avt_k, ldxios = lrxios )
780               CALL iom_get( numror, jpdom_autoglo, 'avm_k', avm_k, ldxios = lrxios )
781               CALL iom_get( numror, jpdom_autoglo, 'dissl', dissl, ldxios = lrxios )
782            ELSE                                          ! start TKE from rest
783               IF(lwp) WRITE(numout,*)
784               IF(lwp) WRITE(numout,*) '   ==>>>   previous run without TKE scheme, set en to background values'
785               en   (:,:,:) = rn_emin * wmask(:,:,:)
786               dissl(:,:,:) = 1.e-12_wp
787               ! avt_k, avm_k already set to the background value in zdf_phy_init
788            ENDIF
789         ELSE                                   !* Start from rest
790            IF(lwp) WRITE(numout,*)
791            IF(lwp) WRITE(numout,*) '   ==>>>   start from rest: set en to the background value'
792            en   (:,:,:) = rn_emin * wmask(:,:,:)
793            dissl(:,:,:) = 1.e-12_wp
794            ! avt_k, avm_k already set to the background value in zdf_phy_init
795         ENDIF
796         !
797      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
798         !                                   ! -------------------
799         IF(lwp) WRITE(numout,*) '---- tke_rst ----'
800         IF( lwxios ) CALL iom_swap(      cwxios_context          ) 
801         CALL iom_rstput( kt, nitrst, numrow, 'en'   , en   , ldxios = lwxios )
802         CALL iom_rstput( kt, nitrst, numrow, 'avt_k', avt_k, ldxios = lwxios )
803         CALL iom_rstput( kt, nitrst, numrow, 'avm_k', avm_k, ldxios = lwxios )
804         CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl, ldxios = lwxios )
805         IF( lwxios ) CALL iom_swap(      cxios_context          )
806         !
807      ENDIF
808      !
809   END SUBROUTINE tke_rst
810
811   !!======================================================================
812END MODULE zdftke
Note: See TracBrowser for help on using the repository browser.