New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
zdftke.F90 in NEMO/trunk/src/OCE/ZDF – NEMO

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

Last change on this file since 13237 was 13237, checked in by smasson, 4 years ago

trunk: Mid-year merge, merge back KERNEL-06_techene_e3

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