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

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

source: NEMO/branches/UKMO/dev_r12866_HPC-02_Daley_Tiling_trial_extra_halo/src/OCE/ZDF/zdftke.F90 @ 12958

Last change on this file since 12958 was 12958, checked in by hadcv, 4 years ago

Merge in trunk changes to r12933

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