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

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

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/ZDF/zdftke.F90 @ 12340

Last change on this file since 12340 was 12340, checked in by acc, 4 years ago

Branch 2019/dev_r11943_MERGE_2019. This commit introduces basic do loop macro
substitution to the 2019 option 1, merge branch. These changes have been SETTE
tested. The only addition is the do_loop_substitute.h90 file in the OCE directory but
the macros defined therein are used throughout the code to replace identifiable, 2D-
and 3D- nested loop opening and closing statements with single-line alternatives. Code
indents are also adjusted accordingly.

The following explanation is taken from comments in the new header file:

This header file contains preprocessor definitions and macros used in the do-loop
substitutions introduced between version 4.0 and 4.2. The primary aim of these macros
is to assist in future applications of tiling to improve performance. This is expected
to be achieved by alternative versions of these macros in selected locations. The
initial introduction of these macros simply replaces all identifiable nested 2D- and
3D-loops with single line statements (and adjusts indenting accordingly). Do loops
are identifiable if they comform to either:

DO jk = ....

DO jj = .... DO jj = ...

DO ji = .... DO ji = ...
. OR .
. .

END DO END DO

END DO END DO

END DO

and white-space variants thereof.

Additionally, only loops with recognised jj and ji loops limits are treated; these are:
Lower limits of 1, 2 or fs_2
Upper limits of jpi, jpim1 or fs_jpim1 (for ji) or jpj, jpjm1 or fs_jpjm1 (for jj)

The macro naming convention takes the form: DO_2D_BT_LR where:

B is the Bottom offset from the PE's inner domain;
T is the Top offset from the PE's inner domain;
L is the Left offset from the PE's inner domain;
R is the Right offset from the PE's inner domain

So, given an inner domain of 2,jpim1 and 2,jpjm1, a typical example would replace:

DO jj = 2, jpj

DO ji = 1, jpim1
.
.

END DO

END DO

with:

DO_2D_01_10
.
.
END_2D

similar conventions apply to the 3D loops macros. jk loop limits are retained
through macro arguments and are not restricted. This includes the possibility of
strides for which an extra set of DO_3DS macros are defined.

In the example definition below the inner PE domain is defined by start indices of
(kIs, kJs) and end indices of (kIe, KJe)

#define DO_2D_00_00 DO jj = kJs, kJe ; DO ji = kIs, kIe
#define END_2D END DO ; END DO

TO DO:


Only conventional nested loops have been identified and replaced by this step. There are constructs such as:

DO jk = 2, jpkm1

z2d(:,:) = z2d(:,:) + e3w(:,:,jk,Kmm) * z3d(:,:,jk) * wmask(:,:,jk)

END DO

which may need to be considered.

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