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/releases/r4.0/r4.0-HEAD/src/OCE/ZDF – NEMO

source: NEMO/releases/r4.0/r4.0-HEAD/src/OCE/ZDF/zdftke.F90 @ 12703

Last change on this file since 12703 was 12703, checked in by mathiot, 4 years ago

ticket #2406: too long line (> 132 cols)

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