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 branches/UKMO/dev_r5518_GO6_package/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/UKMO/dev_r5518_GO6_package/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 6491

Last change on this file since 6491 was 6491, checked in by davestorkey, 8 years ago

Commit remaining changes to UKMO/r5518_GO6_package branch from following branches:
UKMO/dev_r5021_nn_etau_revision@6238
UKMO/dev_r5107_mld_zint@5534
UKMO/dev_r5107_eorca025_closea@6390
UKMO/restart_datestamp@5539
UKMO/icebergs_latent_heat@5821
UKMO/icebergs_restart_single_file_corrected@6480
UKMO/product_diagnostics@5971
UKMO/antarctic_partial_slip@5961
Custom merge into /branches/UKMO/dev_r5518_GO6_package/NEMOGCM: r5961 cf. r5958 of /branches/UKMO/antarctic_partial_slip/NEMOGCM@6490

Custom merge into /branches/UKMO/dev_r5518_GO6_package/NEMOGCM: r6349 cf. r5962 of /branches/UKMO/product_diagnostics/NEMOGCM@6490

Custom merge into /branches/UKMO/dev_r5518_GO6_package/NEMOGCM: r6480 cf. r6479 of /branches/UKMO/icebergs_restart_single_file_corrected/NEMOGCM@6490

Custom merge into /branches/UKMO/dev_r5518_GO6_package/NEMOGCM: r5986 cf. r5852 of /branches/UKMO/icebergs_restart_single_file/NEMOGCM@6490

Custom merge into /branches/UKMO/dev_r5518_GO6_package/NEMOGCM: r5821 cf. r5808 of /branches/UKMO/icebergs_latent_heat/NEMOGCM@6490

File size: 50.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   !!----------------------------------------------------------------------
30#if defined key_zdftke   ||   defined key_esopa
31   !!----------------------------------------------------------------------
32   !!   'key_zdftke'                                   TKE vertical physics
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 zdf_oce        ! vertical physics: ocean variables
46   USE zdfmxl         ! vertical physics: mixed layer
47   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
48   USE prtctl         ! Print control
49   USE in_out_manager ! I/O manager
50   USE iom            ! I/O manager library
51   USE lib_mpp        ! MPP library
52   USE wrk_nemo       ! work arrays
53   USE timing         ! Timing
54   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
55#if defined key_agrif
56   USE agrif_opa_interp
57   USE agrif_opa_update
58#endif
59
60
61
62   IMPLICIT NONE
63   PRIVATE
64
65   PUBLIC   zdf_tke        ! routine called in step module
66   PUBLIC   zdf_tke_init   ! routine called in opa module
67   PUBLIC   tke_rst        ! routine called in step module
68
69   LOGICAL , PUBLIC, PARAMETER ::   lk_zdftke = .TRUE.  !: TKE vertical mixing flag
70
71   !                      !!** Namelist  namzdf_tke  **
72   LOGICAL  ::   ln_mxl0   ! mixing length scale surface value as function of wind stress or not
73   INTEGER  ::   nn_mxl    ! type of mixing length (=0/1/2/3)
74   REAL(wp) ::   rn_mxl0   ! surface  min value of mixing length (kappa*z_o=0.4*0.1 m)  [m]
75   INTEGER  ::   nn_pdl    ! Prandtl number or not (ratio avt/avm) (=0/1)
76   REAL(wp) ::   rn_ediff  ! coefficient for avt: avt=rn_ediff*mxl*sqrt(e)
77   REAL(wp) ::   rn_ediss  ! coefficient of the Kolmogoroff dissipation
78   REAL(wp) ::   rn_ebb    ! coefficient of the surface input of tke
79   REAL(wp) ::   rn_emin   ! minimum value of tke           [m2/s2]
80   REAL(wp) ::   rn_emin0  ! surface minimum value of tke   [m2/s2]
81   REAL(wp) ::   rn_bshear ! background shear (>0) currently a numerical threshold (do not change it)
82   INTEGER  ::   nn_etau   ! type of depth penetration of surface tke (=0/1/2/3)
83   INTEGER  ::   nn_htau   ! type of tke profile of penetration (=0/1)
84   REAL(wp) ::   rn_efr    ! fraction of TKE surface value which penetrates in the ocean
85   REAL(wp) ::   rn_c      ! fraction of TKE added within the mixed layer by nn_etau
86   LOGICAL  ::   ln_lc     ! Langmuir cells (LC) as a source term of TKE or not
87   REAL(wp) ::   rn_lc     ! coef to compute vertical velocity of Langmuir cells
88
89   REAL(wp) ::   ri_cri    ! critic Richardson number (deduced from rn_ediff and rn_ediss values)
90   REAL(wp) ::   rmxl_min  ! minimum mixing length value (deduced from rn_ediff and rn_emin values)  [m]
91   REAL(wp) ::   rhtau                     ! coefficient to relate MLD to htau when nn_htau == 2
92   REAL(wp) ::   rhftau_add = 1.e-3_wp     ! add offset   applied to HF part of taum  (nn_etau=3)
93   REAL(wp) ::   rhftau_scl = 1.0_wp       ! scale factor applied to HF part of taum  (nn_etau=3)
94
95   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   htau           ! depth of tke penetration (nn_htau)
96   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   e_niw          !: TKE budget- near-inertial waves term
97   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   efr            ! surface boundary condition for nn_etau = 4
98   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dissl          ! now mixing lenght of dissipation
99#if defined key_c1d
100   !                                                                        !!** 1D cfg only  **   ('key_c1d')
101   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   e_dis, e_mix   !: dissipation and mixing turbulent lengh scales
102   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   e_pdl, e_ric   !: prandl and local Richardson numbers
103#endif
104
105   !! * Substitutions
106#  include "domzgr_substitute.h90"
107#  include "vectopt_loop_substitute.h90"
108   !!----------------------------------------------------------------------
109   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
110   !! $Id$
111   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
112   !!----------------------------------------------------------------------
113CONTAINS
114
115   INTEGER FUNCTION zdf_tke_alloc()
116      !!----------------------------------------------------------------------
117      !!                ***  FUNCTION zdf_tke_alloc  ***
118      !!----------------------------------------------------------------------
119      ALLOCATE(                                                                    &
120         &      efr  (jpi,jpj)     , e_niw(jpi,jpj,jpk) ,                         &     
121#if defined key_c1d
122         &      e_dis(jpi,jpj,jpk) , e_mix(jpi,jpj,jpk) ,                          &
123         &      e_pdl(jpi,jpj,jpk) , e_ric(jpi,jpj,jpk) ,                          &
124#endif
125         &      htau  (jpi,jpj)    , dissl(jpi,jpj,jpk) , STAT= zdf_tke_alloc      )
126         !
127      IF( lk_mpp             )   CALL mpp_sum ( zdf_tke_alloc )
128      IF( zdf_tke_alloc /= 0 )   CALL ctl_warn('zdf_tke_alloc: failed to allocate arrays')
129      !
130   END FUNCTION zdf_tke_alloc
131
132
133   SUBROUTINE zdf_tke( kt )
134      !!----------------------------------------------------------------------
135      !!                   ***  ROUTINE zdf_tke  ***
136      !!
137      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
138      !!              coefficients using a turbulent closure scheme (TKE).
139      !!
140      !! ** Method  :   The time evolution of the turbulent kinetic energy (tke)
141      !!              is computed from a prognostic equation :
142      !!         d(en)/dt = avm (d(u)/dz)**2             ! shear production
143      !!                  + d( avm d(en)/dz )/dz         ! diffusion of tke
144      !!                  + avt N^2                      ! stratif. destruc.
145      !!                  - rn_ediss / emxl en**(2/3)    ! Kolmogoroff dissipation
146      !!      with the boundary conditions:
147      !!         surface: en = max( rn_emin0, rn_ebb * taum )
148      !!         bottom : en = rn_emin
149      !!      The associated critical Richardson number is: ri_cri = 2/(2+rn_ediss/rn_ediff)
150      !!
151      !!        The now Turbulent kinetic energy is computed using the following
152      !!      time stepping: implicit for vertical diffusion term, linearized semi
153      !!      implicit for kolmogoroff dissipation term, and explicit forward for
154      !!      both buoyancy and shear production terms. Therefore a tridiagonal
155      !!      linear system is solved. Note that buoyancy and shear terms are
156      !!      discretized in a energy conserving form (Bruchard 2002).
157      !!
158      !!        The dissipative and mixing length scale are computed from en and
159      !!      the stratification (see tke_avn)
160      !!
161      !!        The now vertical eddy vicosity and diffusivity coefficients are
162      !!      given by:
163      !!              avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
164      !!              avt = max( avmb, pdl * avm                 ) 
165      !!              eav = max( avmb, avm )
166      !!      where pdl, the inverse of the Prandtl number is 1 if nn_pdl=0 and
167      !!      given by an empirical funtion of the localRichardson number if nn_pdl=1
168      !!
169      !! ** Action  :   compute en (now turbulent kinetic energy)
170      !!                update avt, avmu, avmv (before vertical eddy coef.)
171      !!
172      !! References : Gaspar et al., JGR, 1990,
173      !!              Blanke and Delecluse, JPO, 1991
174      !!              Mellor and Blumberg, JPO 2004
175      !!              Axell, JGR, 2002
176      !!              Bruchard OM 2002
177      !!----------------------------------------------------------------------
178      INTEGER, INTENT(in) ::   kt   ! ocean time step
179      !!----------------------------------------------------------------------
180      !
181      IF( kt /= nit000 ) THEN   ! restore before value to compute tke
182         avt (:,:,:) = avt_k (:,:,:) 
183         avm (:,:,:) = avm_k (:,:,:) 
184         avmu(:,:,:) = avmu_k(:,:,:) 
185         avmv(:,:,:) = avmv_k(:,:,:) 
186      ENDIF 
187      !
188      CALL tke_tke      ! now tke (en)
189      !
190      CALL tke_avn      ! now avt, avm, avmu, avmv
191      !
192      avt_k (:,:,:) = avt (:,:,:) 
193      avm_k (:,:,:) = avm (:,:,:) 
194      avmu_k(:,:,:) = avmu(:,:,:) 
195      avmv_k(:,:,:) = avmv(:,:,:) 
196      !
197#if defined key_agrif
198      ! Update child grid f => parent grid
199      IF( .NOT.Agrif_Root() )   CALL Agrif_Update_Tke( kt )      ! children only
200#endif     
201     !
202   END SUBROUTINE zdf_tke
203
204
205   SUBROUTINE tke_tke
206      !!----------------------------------------------------------------------
207      !!                   ***  ROUTINE tke_tke  ***
208      !!
209      !! ** Purpose :   Compute the now Turbulente Kinetic Energy (TKE)
210      !!
211      !! ** Method  : - TKE surface boundary condition
212      !!              - source term due to Langmuir cells (Axell JGR 2002) (ln_lc=T)
213      !!              - source term due to shear (saved in avmu, avmv arrays)
214      !!              - Now TKE : resolution of the TKE equation by inverting
215      !!                a tridiagonal linear system by a "methode de chasse"
216      !!              - increase TKE due to surface and internal wave breaking
217      !!
218      !! ** Action  : - en : now turbulent kinetic energy)
219      !!              - avmu, avmv : production of TKE by shear at u and v-points
220      !!                (= Kz dz[Ub] * dz[Un] )
221      !! ---------------------------------------------------------------------
222      INTEGER  ::   ji, jj, jk                      ! dummy loop arguments
223!!bfr      INTEGER  ::   ikbu, ikbv, ikbum1, ikbvm1      ! temporary scalar
224!!bfr      INTEGER  ::   ikbt, ikbumm1, ikbvmm1          ! temporary scalar
225      REAL(wp) ::   zrhoa  = 1.22                   ! Air density kg/m3
226      REAL(wp) ::   zcdrag = 1.5e-3                 ! drag coefficient
227      REAL(wp) ::   zbbrau, zesh2                   ! temporary scalars
228      REAL(wp) ::   zfact1, zfact2, zfact3          !    -         -
229      REAL(wp) ::   ztx2  , zty2  , zcof            !    -         -
230      REAL(wp) ::   ztau  , zdif                    !    -         -
231      REAL(wp) ::   zus   , zwlc  , zind            !    -         -
232      REAL(wp) ::   zzd_up, zzd_lw                  !    -         -
233!!bfr      REAL(wp) ::   zebot                           !    -         -
234      INTEGER , POINTER, DIMENSION(:,:  ) :: imlc
235      REAL(wp), POINTER, DIMENSION(:,:  ) :: zhlc
236      REAL(wp), POINTER, DIMENSION(:,:,:) :: zpelc, zdiag, zd_up, zd_lw
237      !!--------------------------------------------------------------------
238      !
239      IF( nn_timing == 1 )  CALL timing_start('tke_tke')
240      !
241      CALL wrk_alloc( jpi,jpj, imlc )    ! integer
242      CALL wrk_alloc( jpi,jpj, zhlc ) 
243      CALL wrk_alloc( jpi,jpj,jpk, zpelc, zdiag, zd_up, zd_lw ) 
244      !
245      zbbrau = rn_ebb / rau0       ! Local constant initialisation
246      zfact1 = -.5_wp * rdt 
247      zfact2 = 1.5_wp * rdt * rn_ediss
248      zfact3 = 0.5_wp       * rn_ediss
249      !
250      !
251      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
252      !                     !  Surface boundary condition on tke
253      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
254      IF ( ln_isfcav ) THEN
255         DO jj = 2, jpjm1            ! en(mikt(ji,jj))   = rn_emin
256            DO ji = fs_2, fs_jpim1   ! vector opt.
257               en(ji,jj,mikt(ji,jj))=rn_emin * tmask(ji,jj,1)
258            END DO
259         END DO
260      END IF
261      DO jj = 2, jpjm1            ! en(1)   = rn_ebb taum / rau0  (min value rn_emin0)
262         DO ji = fs_2, fs_jpim1   ! vector opt.
263            en(ji,jj,1) = MAX( rn_emin0, zbbrau * taum(ji,jj) ) * tmask(ji,jj,1)
264         END DO
265      END DO
266     
267!!bfr   - start commented area
268      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
269      !                     !  Bottom boundary condition on tke
270      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
271      !
272      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
273      ! Tests to date have found the bottom boundary condition on tke to have very little effect.
274      ! The condition is coded here for completion but commented out until there is proof that the
275      ! computational cost is justified
276      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
277      !                     en(bot)   = (rn_ebb0/rau0)*0.5*sqrt(u_botfr^2+v_botfr^2) (min value rn_emin)
278!CDIR NOVERRCHK
279!!    DO jj = 2, jpjm1
280!CDIR NOVERRCHK
281!!       DO ji = fs_2, fs_jpim1   ! vector opt.
282!!          ztx2 = bfrua(ji-1,jj) * ub(ji-1,jj,mbku(ji-1,jj)) + &
283!!                 bfrua(ji  ,jj) * ub(ji  ,jj,mbku(ji  ,jj) )
284!!          zty2 = bfrva(ji,jj  ) * vb(ji,jj  ,mbkv(ji,jj  )) + &
285!!                 bfrva(ji,jj-1) * vb(ji,jj-1,mbkv(ji,jj-1) )
286!!          zebot = 0.001875_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 )   !  where 0.001875 = (rn_ebb0/rau0) * 0.5 = 3.75*0.5/1000.
287!!          en (ji,jj,mbkt(ji,jj)+1) = MAX( zebot, rn_emin ) * tmask(ji,jj,1)
288!!       END DO
289!!    END DO
290!!bfr   - end commented area
291      !
292      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
293      IF( ln_lc ) THEN      !  Langmuir circulation source term added to tke       (Axell JGR 2002)
294         !                  !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
295         !
296         !                        !* total energy produce by LC : cumulative sum over jk
297         zpelc(:,:,1) =  MAX( rn2b(:,:,1), 0._wp ) * fsdepw(:,:,1) * fse3w(:,:,1)
298         DO jk = 2, jpk
299            zpelc(:,:,jk)  = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0._wp ) * fsdepw(:,:,jk) * fse3w(:,:,jk)
300         END DO
301         !                        !* finite Langmuir Circulation depth
302         zcof = 0.5 * 0.016 * 0.016 / ( zrhoa * zcdrag )
303         imlc(:,:) = mbkt(:,:) + 1       ! Initialization to the number of w ocean point (=2 over land)
304         DO jk = jpkm1, 2, -1
305            DO jj = 1, jpj               ! Last w-level at which zpelc>=0.5*us*us
306               DO ji = 1, jpi            !      with us=0.016*wind(starting from jpk-1)
307                  zus  = zcof * taum(ji,jj)
308                  IF( zpelc(ji,jj,jk) > zus )   imlc(ji,jj) = jk
309               END DO
310            END DO
311         END DO
312         !                               ! finite LC depth
313         DO jj = 1, jpj 
314            DO ji = 1, jpi
315               zhlc(ji,jj) = fsdepw(ji,jj,imlc(ji,jj))
316            END DO
317         END DO
318         zcof = 0.016 / SQRT( zrhoa * zcdrag )
319!CDIR NOVERRCHK
320         DO jk = 2, jpkm1         !* TKE Langmuir circulation source term added to en
321!CDIR NOVERRCHK
322            DO jj = 2, jpjm1
323!CDIR NOVERRCHK
324               DO ji = fs_2, fs_jpim1   ! vector opt.
325                  zus  = zcof * SQRT( taum(ji,jj) )           ! Stokes drift
326                  !                                           ! vertical velocity due to LC
327                  zind = 0.5 - SIGN( 0.5, fsdepw(ji,jj,jk) - zhlc(ji,jj) )
328                  zwlc = zind * rn_lc * zus * SIN( rpi * fsdepw(ji,jj,jk) / zhlc(ji,jj) )
329                  !                                           ! TKE Langmuir circulation source term
330                  en(ji,jj,jk) = en(ji,jj,jk) + rdt * ( 1._wp - fr_i(ji,jj) ) * ( zwlc * zwlc * zwlc ) /   &
331                     &   zhlc(ji,jj) * wmask(ji,jj,jk) * tmask(ji,jj,1)
332               END DO
333            END DO
334         END DO
335         !
336      ENDIF
337      !
338      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
339      !                     !  Now Turbulent kinetic energy (output in en)
340      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
341      !                     ! Resolution of a tridiagonal linear system by a "methode de chasse"
342      !                     ! computation from level 2 to jpkm1  (e(1) already computed and e(jpk)=0 ).
343      !                     ! zdiag : diagonal zd_up : upper diagonal zd_lw : lower diagonal
344      !
345      DO jk = 2, jpkm1           !* Shear production at uw- and vw-points (energy conserving form)
346         DO jj = 1, jpj                 ! here avmu, avmv used as workspace
347            DO ji = 1, jpi
348               avmu(ji,jj,jk) = avmu(ji,jj,jk) * ( un(ji,jj,jk-1) - un(ji,jj,jk) )   &
349                  &                            * ( ub(ji,jj,jk-1) - ub(ji,jj,jk) )   & 
350                  &                            / (  fse3uw_n(ji,jj,jk)               &
351                  &                              *  fse3uw_b(ji,jj,jk)  )
352               avmv(ji,jj,jk) = avmv(ji,jj,jk) * ( vn(ji,jj,jk-1) - vn(ji,jj,jk) )   &
353                  &                            * ( vb(ji,jj,jk-1) - vb(ji,jj,jk) )   &
354                  &                            / (  fse3vw_n(ji,jj,jk)               &
355                  &                              *  fse3vw_b(ji,jj,jk)  )
356            END DO
357         END DO
358      END DO
359      !
360      DO jk = 2, jpkm1           !* Matrix and right hand side in en
361         DO jj = 2, jpjm1
362            DO ji = fs_2, fs_jpim1   ! vector opt.
363               zcof   = zfact1 * tmask(ji,jj,jk)
364               zzd_up = zcof * ( avm  (ji,jj,jk+1) + avm  (ji,jj,jk  ) )   &  ! upper diagonal
365                  &          / ( fse3t(ji,jj,jk  ) * fse3w(ji,jj,jk  ) )
366               zzd_lw = zcof * ( avm  (ji,jj,jk  ) + avm  (ji,jj,jk-1) )   &  ! lower diagonal
367                  &          / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk  ) )
368                  !                                                           ! shear prod. at w-point weightened by mask
369               zesh2  =  ( avmu(ji-1,jj,jk) + avmu(ji,jj,jk) ) / MAX( 1._wp , umask(ji-1,jj,jk) + umask(ji,jj,jk) )   &
370                  &    + ( avmv(ji,jj-1,jk) + avmv(ji,jj,jk) ) / MAX( 1._wp , vmask(ji,jj-1,jk) + vmask(ji,jj,jk) )   
371                  !
372               zd_up(ji,jj,jk) = zzd_up            ! Matrix (zdiag, zd_up, zd_lw)
373               zd_lw(ji,jj,jk) = zzd_lw
374               zdiag(ji,jj,jk) = 1._wp - zzd_lw - zzd_up + zfact2 * dissl(ji,jj,jk) * tmask(ji,jj,jk)
375               !
376               !                                   ! right hand side in en
377               en(ji,jj,jk) = en(ji,jj,jk) + rdt * (  zesh2  -   avt(ji,jj,jk) * rn2(ji,jj,jk)    &
378                  &                                 + zfact3 * dissl(ji,jj,jk) * en (ji,jj,jk)  ) &
379                  &                                 * wmask(ji,jj,jk)
380            END DO
381         END DO
382      END DO
383      !                          !* Matrix inversion from level 2 (tke prescribed at level 1)
384      DO jk = 3, jpkm1                             ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
385         DO jj = 2, jpjm1
386            DO ji = fs_2, fs_jpim1    ! vector opt.
387               zdiag(ji,jj,jk) = zdiag(ji,jj,jk) - zd_lw(ji,jj,jk) * zd_up(ji,jj,jk-1) / zdiag(ji,jj,jk-1)
388            END DO
389         END DO
390      END DO
391      !
392      ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1
393      DO jj = 2, jpjm1
394         DO ji = fs_2, fs_jpim1   ! vector opt.
395            zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1)    ! Surface boudary conditions on tke
396         END DO
397      END DO
398      DO jk = 3, jpkm1
399         DO jj = 2, jpjm1
400            DO ji = fs_2, fs_jpim1    ! vector opt.
401               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)
402            END DO
403         END DO
404      END DO
405      !
406      ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
407      DO jj = 2, jpjm1
408         DO ji = fs_2, fs_jpim1   ! vector opt.
409            en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1)
410         END DO
411      END DO
412      DO jk = jpk-2, 2, -1
413         DO jj = 2, jpjm1
414            DO ji = fs_2, fs_jpim1    ! vector opt.
415               en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk)
416            END DO
417         END DO
418      END DO
419      DO jk = 2, jpkm1                             ! set the minimum value of tke
420         DO jj = 2, jpjm1
421            DO ji = fs_2, fs_jpim1   ! vector opt.
422               en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * wmask(ji,jj,jk)
423            END DO
424         END DO
425      END DO
426
427      !                                 ! Save TKE prior to nn_etau addition 
428      e_niw(:,:,:) = en(:,:,:) 
429     
430      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
431      !                            !  TKE due to surface and internal wave breaking
432      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
433      IF( nn_htau == 2 ) THEN           !* mixed-layer depth dependant length scale
434         DO jj = 2, jpjm1
435            DO ji = fs_2, fs_jpim1   ! vector opt.
436               htau(ji,jj) = rhtau * hmlp(ji,jj)
437            END DO
438         END DO
439      ENDIF
440#if defined key_iomput
441      !
442      CALL iom_put( "htau", htau(:,:) )  ! Check htau (even if constant in time)
443#endif
444      !
445      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction)
446         DO jk = 2, jpkm1
447            DO jj = 2, jpjm1
448               DO ji = fs_2, fs_jpim1   ! vector opt.
449                  en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
450                     &                                 * ( 1._wp - fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
451               END DO
452            END DO
453         END DO
454      ELSEIF( nn_etau == 2 ) THEN       !* act only at the base of the mixed layer (jk=nmln)  (rn_efr fraction)
455         DO jj = 2, jpjm1
456            DO ji = fs_2, fs_jpim1   ! vector opt.
457               jk = nmln(ji,jj)
458               en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
459                  &                                 * ( 1._wp - fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
460            END DO
461         END DO
462      ELSEIF( nn_etau == 3 ) THEN       !* penetration belox the mixed layer (HF variability)
463!CDIR NOVERRCHK
464         DO jk = 2, jpkm1
465!CDIR NOVERRCHK
466            DO jj = 2, jpjm1
467!CDIR NOVERRCHK
468               DO ji = fs_2, fs_jpim1   ! vector opt.
469                  ztx2 = utau(ji-1,jj  ) + utau(ji,jj)
470                  zty2 = vtau(ji  ,jj-1) + vtau(ji,jj)
471                  ztau = 0.5_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1)    ! module of the mean stress
472                  zdif = taum(ji,jj) - ztau                            ! mean of modulus - modulus of the mean
473                  zdif = rhftau_scl * MAX( 0._wp, zdif + rhftau_add )  ! apply some modifications...
474                  en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
475                     &                        * ( 1._wp - fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1)
476               END DO
477            END DO
478         END DO
479      ELSEIF( nn_etau == 4 ) THEN       !* column integral independant of htau (rn_efr must be scaled up)
480         IF( nn_htau == 2 ) THEN        ! efr dependant on time-varying htau
481            DO jj = 2, jpjm1
482               DO ji = fs_2, fs_jpim1   ! vector opt.
483                  efr(ji,jj) = rn_efr / ( htau(ji,jj) * ( 1._wp - EXP( -bathy(ji,jj) / htau(ji,jj) ) ) )
484               END DO
485            END DO
486         ENDIF
487         DO jk = 2, jpkm1
488            DO jj = 2, jpjm1
489               DO ji = fs_2, fs_jpim1   ! vector opt.
490                  en(ji,jj,jk) = en(ji,jj,jk) + efr(ji,jj) * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
491                     &                                                   * ( 1._wp - fr_i(ji,jj) )  * tmask(ji,jj,jk)
492               END DO
493            END DO
494         END DO
495      ENDIF
496      CALL lbc_lnk( en, 'W', 1. )      ! Lateral boundary conditions (sign unchanged)
497      !
498      DO jk = 2, jpkm1                             ! TKE budget: near-inertial waves term 
499         DO jj = 2, jpjm1 
500            DO ji = fs_2, fs_jpim1   ! vector opt. 
501               e_niw(ji,jj,jk) = en(ji,jj,jk) - e_niw(ji,jj,jk) 
502            END DO 
503         END DO 
504      END DO 
505     
506      CALL lbc_lnk( e_niw, 'W', 1. ) 
507      !
508      CALL wrk_dealloc( jpi,jpj, imlc )    ! integer
509      CALL wrk_dealloc( jpi,jpj, zhlc ) 
510      CALL wrk_dealloc( jpi,jpj,jpk, zpelc, zdiag, zd_up, zd_lw ) 
511      !
512      IF( nn_timing == 1 )  CALL timing_stop('tke_tke')
513      !
514   END SUBROUTINE tke_tke
515
516
517   SUBROUTINE tke_avn
518      !!----------------------------------------------------------------------
519      !!                   ***  ROUTINE tke_avn  ***
520      !!
521      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
522      !!
523      !! ** Method  :   At this stage, en, the now TKE, is known (computed in
524      !!              the tke_tke routine). First, the now mixing lenth is
525      !!      computed from en and the strafification (N^2), then the mixings
526      !!      coefficients are computed.
527      !!              - Mixing length : a first evaluation of the mixing lengh
528      !!      scales is:
529      !!                      mxl = sqrt(2*en) / N 
530      !!      where N is the brunt-vaisala frequency, with a minimum value set
531      !!      to rmxl_min (rn_mxl0) in the interior (surface) ocean.
532      !!        The mixing and dissipative length scale are bound as follow :
533      !!         nn_mxl=0 : mxl bounded by the distance to surface and bottom.
534      !!                        zmxld = zmxlm = mxl
535      !!         nn_mxl=1 : mxl bounded by the e3w and zmxld = zmxlm = mxl
536      !!         nn_mxl=2 : mxl bounded such that the vertical derivative of mxl is
537      !!                    less than 1 (|d/dz(mxl)|<1) and zmxld = zmxlm = mxl
538      !!         nn_mxl=3 : mxl is bounded from the surface to the bottom usings
539      !!                    |d/dz(xml)|<1 to obtain lup, and from the bottom to
540      !!                    the surface to obtain ldown. the resulting length
541      !!                    scales are:
542      !!                        zmxld = sqrt( lup * ldown )
543      !!                        zmxlm = min ( lup , ldown )
544      !!              - Vertical eddy viscosity and diffusivity:
545      !!                      avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
546      !!                      avt = max( avmb, pdlr * avm ) 
547      !!      with pdlr=1 if nn_pdl=0, pdlr=1/pdl=F(Ri) otherwise.
548      !!
549      !! ** Action  : - avt : now vertical eddy diffusivity (w-point)
550      !!              - avmu, avmv : now vertical eddy viscosity at uw- and vw-points
551      !!----------------------------------------------------------------------
552      INTEGER  ::   ji, jj, jk   ! dummy loop indices
553      REAL(wp) ::   zrn2, zraug, zcoef, zav     ! local scalars
554      REAL(wp) ::   zdku, zpdlr, zri, zsqen     !   -      -
555      REAL(wp) ::   zdkv, zemxl, zemlm, zemlp   !   -      -
556      REAL(wp), POINTER, DIMENSION(:,:,:) :: zmpdl, zmxlm, zmxld
557      !!--------------------------------------------------------------------
558      !
559      IF( nn_timing == 1 )  CALL timing_start('tke_avn')
560
561      CALL wrk_alloc( jpi,jpj,jpk, zmpdl, zmxlm, zmxld ) 
562
563      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
564      !                     !  Mixing length
565      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
566      !
567      !                     !* Buoyancy length scale: l=sqrt(2*e/n**2)
568      !
569      ! initialisation of interior minimum value (avoid a 2d loop with mikt)
570      zmxlm(:,:,:)  = rmxl_min   
571      zmxld(:,:,:)  = rmxl_min
572      !
573      IF( ln_mxl0 ) THEN            ! surface mixing length = F(stress) : l=vkarmn*2.e5*taum/(rau0*g)
574         DO jj = 2, jpjm1
575            DO ji = fs_2, fs_jpim1
576               zraug = vkarmn * 2.e5_wp / ( rau0 * grav )
577               zmxlm(ji,jj,1) = MAX( rn_mxl0, zraug * taum(ji,jj) * tmask(ji,jj,1) )
578            END DO
579         END DO
580      ELSE
581         zmxlm(:,:,1) = rn_mxl0
582      ENDIF
583      !
584!CDIR NOVERRCHK
585      DO jk = 2, jpkm1              ! interior value : l=sqrt(2*e/n^2)
586!CDIR NOVERRCHK
587         DO jj = 2, jpjm1
588!CDIR NOVERRCHK
589            DO ji = fs_2, fs_jpim1   ! vector opt.
590               zrn2 = MAX( rn2(ji,jj,jk), rsmall )
591               zmxlm(ji,jj,jk) = MAX(  rmxl_min,  SQRT( 2._wp * en(ji,jj,jk) / zrn2 ) )
592            END DO
593         END DO
594      END DO
595      !
596      !                     !* Physical limits for the mixing length
597      !
598      zmxld(:,:,1  ) = zmxlm(:,:,1)   ! surface set to the minimum value
599      zmxld(:,:,jpk) = rmxl_min       ! last level  set to the minimum value
600      !
601      SELECT CASE ( nn_mxl )
602      !
603      ! where wmask = 0 set zmxlm == fse3w
604      CASE ( 0 )           ! bounded by the distance to surface and bottom
605         DO jk = 2, jpkm1
606            DO jj = 2, jpjm1
607               DO ji = fs_2, fs_jpim1   ! vector opt.
608                  zemxl = MIN( fsdepw(ji,jj,jk) - fsdepw(ji,jj,mikt(ji,jj)), zmxlm(ji,jj,jk),   &
609                  &            fsdepw(ji,jj,mbkt(ji,jj)+1) - fsdepw(ji,jj,jk) )
610                  ! wmask prevent zmxlm = 0 if jk = mikt(ji,jj)
611                  zmxlm(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN(zmxlm(ji,jj,jk),fse3w(ji,jj,jk)) * (1 - wmask(ji,jj,jk))
612                  zmxld(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN(zmxlm(ji,jj,jk),fse3w(ji,jj,jk)) * (1 - wmask(ji,jj,jk))
613               END DO
614            END DO
615         END DO
616         !
617      CASE ( 1 )           ! bounded by the vertical scale factor
618         DO jk = 2, jpkm1
619            DO jj = 2, jpjm1
620               DO ji = fs_2, fs_jpim1   ! vector opt.
621                  zemxl = MIN( fse3w(ji,jj,jk), zmxlm(ji,jj,jk) )
622                  zmxlm(ji,jj,jk) = zemxl
623                  zmxld(ji,jj,jk) = zemxl
624               END DO
625            END DO
626         END DO
627         !
628      CASE ( 2 )           ! |dk[xml]| bounded by e3t :
629         DO jk = 2, jpkm1         ! from the surface to the bottom :
630            DO jj = 2, jpjm1
631               DO ji = fs_2, fs_jpim1   ! vector opt.
632                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
633               END DO
634            END DO
635         END DO
636         DO jk = jpkm1, 2, -1     ! from the bottom to the surface :
637            DO jj = 2, jpjm1
638               DO ji = fs_2, fs_jpim1   ! vector opt.
639                  zemxl = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
640                  zmxlm(ji,jj,jk) = zemxl
641                  zmxld(ji,jj,jk) = zemxl
642               END DO
643            END DO
644         END DO
645         !
646      CASE ( 3 )           ! lup and ldown, |dk[xml]| bounded by e3t :
647         DO jk = 2, jpkm1         ! from the surface to the bottom : lup
648            DO jj = 2, jpjm1
649               DO ji = fs_2, fs_jpim1   ! vector opt.
650                  zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
651               END DO
652            END DO
653         END DO
654         DO jk = jpkm1, 2, -1     ! from the bottom to the surface : ldown
655            DO jj = 2, jpjm1
656               DO ji = fs_2, fs_jpim1   ! vector opt.
657                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
658               END DO
659            END DO
660         END DO
661!CDIR NOVERRCHK
662         DO jk = 2, jpkm1
663!CDIR NOVERRCHK
664            DO jj = 2, jpjm1
665!CDIR NOVERRCHK
666               DO ji = fs_2, fs_jpim1   ! vector opt.
667                  zemlm = MIN ( zmxld(ji,jj,jk),  zmxlm(ji,jj,jk) )
668                  zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) )
669                  zmxlm(ji,jj,jk) = zemlm
670                  zmxld(ji,jj,jk) = zemlp
671               END DO
672            END DO
673         END DO
674         !
675      END SELECT
676      !
677# if defined key_c1d
678      e_dis(:,:,:) = zmxld(:,:,:)      ! c1d configuration : save mixing and dissipation turbulent length scales
679      e_mix(:,:,:) = zmxlm(:,:,:)
680# endif
681
682      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
683      !                     !  Vertical eddy viscosity and diffusivity  (avmu, avmv, avt)
684      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
685!CDIR NOVERRCHK
686      DO jk = 1, jpkm1            !* vertical eddy viscosity & diffivity at w-points
687!CDIR NOVERRCHK
688         DO jj = 2, jpjm1
689!CDIR NOVERRCHK
690            DO ji = fs_2, fs_jpim1   ! vector opt.
691               zsqen = SQRT( en(ji,jj,jk) )
692               zav   = rn_ediff * zmxlm(ji,jj,jk) * zsqen
693               avm  (ji,jj,jk) = MAX( zav,                  avmb(jk) ) * wmask(ji,jj,jk)
694               avt  (ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
695               dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk)
696            END DO
697         END DO
698      END DO
699      CALL lbc_lnk( avm, 'W', 1. )      ! Lateral boundary conditions (sign unchanged)
700      !
701      DO jk = 2, jpkm1            !* vertical eddy viscosity at wu- and wv-points
702         DO jj = 2, jpjm1
703            DO ji = fs_2, fs_jpim1   ! vector opt.
704               avmu(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji+1,jj  ,jk) ) * wumask(ji,jj,jk)
705               avmv(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji  ,jj+1,jk) ) * wvmask(ji,jj,jk)
706            END DO
707         END DO
708      END DO
709      CALL lbc_lnk( avmu, 'U', 1. )   ;   CALL lbc_lnk( avmv, 'V', 1. )      ! Lateral boundary conditions
710      !
711      IF( nn_pdl == 1 ) THEN      !* Prandtl number case: update avt
712         DO jk = 2, jpkm1
713            DO jj = 2, jpjm1
714               DO ji = fs_2, fs_jpim1   ! vector opt.
715                  zcoef = avm(ji,jj,jk) * 2._wp * fse3w(ji,jj,jk) * fse3w(ji,jj,jk)
716                  !                                          ! shear
717                  zdku = avmu(ji-1,jj,jk) * ( un(ji-1,jj,jk-1) - un(ji-1,jj,jk) ) * ( ub(ji-1,jj,jk-1) - ub(ji-1,jj,jk) )   &
718                    &  + avmu(ji  ,jj,jk) * ( un(ji  ,jj,jk-1) - un(ji  ,jj,jk) ) * ( ub(ji  ,jj,jk-1) - ub(ji  ,jj,jk) )
719                  zdkv = avmv(ji,jj-1,jk) * ( vn(ji,jj-1,jk-1) - vn(ji,jj-1,jk) ) * ( vb(ji,jj-1,jk-1) - vb(ji,jj-1,jk) )   &
720                    &  + avmv(ji,jj  ,jk) * ( vn(ji,jj  ,jk-1) - vn(ji,jj  ,jk) ) * ( vb(ji,jj  ,jk-1) - vb(ji,jj  ,jk) )
721                  !                                          ! local Richardson number
722                  zri   = MAX( rn2b(ji,jj,jk), 0._wp ) * zcoef / (zdku + zdkv + rn_bshear )
723                  zpdlr = MAX(  0.1_wp,  0.2 / MAX( 0.2 , zri )  )
724!!gm and even better with the use of the "true" ri_crit=0.22222...  (this change the results!)
725!!gm              zpdlr = MAX(  0.1_wp,  ri_crit / MAX( ri_crit , zri )  )
726                  avt(ji,jj,jk)   = MAX( zpdlr * avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
727# if defined key_c1d
728                  e_pdl(ji,jj,jk) = zpdlr * wmask(ji,jj,jk)  ! c1d configuration : save masked Prandlt number
729                  e_ric(ji,jj,jk) = zri   * wmask(ji,jj,jk)  ! c1d config. : save Ri
730# endif
731              END DO
732            END DO
733         END DO
734      ENDIF
735      CALL lbc_lnk( avt, 'W', 1. )                      ! Lateral boundary conditions on avt  (sign unchanged)
736
737      IF(ln_ctl) THEN
738         CALL prt_ctl( tab3d_1=en  , clinfo1=' tke  - e: ', tab3d_2=avt, clinfo2=' t: ', ovlap=1, kdim=jpk)
739         CALL prt_ctl( tab3d_1=avmu, clinfo1=' tke  - u: ', mask1=umask,                   &
740            &          tab3d_2=avmv, clinfo2=       ' v: ', mask2=vmask, ovlap=1, kdim=jpk )
741      ENDIF
742      !
743      CALL wrk_dealloc( jpi,jpj,jpk, zmpdl, zmxlm, zmxld ) 
744      !
745      IF( nn_timing == 1 )  CALL timing_stop('tke_avn')
746      !
747   END SUBROUTINE tke_avn
748
749
750   SUBROUTINE zdf_tke_init
751      !!----------------------------------------------------------------------
752      !!                  ***  ROUTINE zdf_tke_init  ***
753      !!                     
754      !! ** Purpose :   Initialization of the vertical eddy diffivity and
755      !!              viscosity when using a tke turbulent closure scheme
756      !!
757      !! ** Method  :   Read the namzdf_tke namelist and check the parameters
758      !!              called at the first timestep (nit000)
759      !!
760      !! ** input   :   Namlist namzdf_tke
761      !!
762      !! ** Action  :   Increase by 1 the nstop flag is setting problem encounter
763      !!----------------------------------------------------------------------
764      INTEGER ::   ji, jj, jk   ! dummy loop indices
765      INTEGER ::   ios, ierr
766      !!
767      NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb , rn_emin  ,   &
768         &                 rn_emin0, rn_bshear, nn_mxl , ln_mxl0  ,   &
769         &                 rn_mxl0 , nn_pdl   , ln_lc  , rn_lc    ,   &
770         &                 nn_etau , nn_htau  , rn_efr , rn_c   
771      !!----------------------------------------------------------------------
772
773      REWIND( numnam_ref )              ! Namelist namzdf_tke in reference namelist : Turbulent Kinetic Energy
774      READ  ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901)
775901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in reference namelist', lwp )
776
777      REWIND( numnam_cfg )              ! Namelist namzdf_tke in configuration namelist : Turbulent Kinetic Energy
778      READ  ( numnam_cfg, namzdf_tke, IOSTAT = ios, ERR = 902 )
779902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in configuration namelist', lwp )
780      IF(lwm) WRITE ( numond, namzdf_tke )
781      !
782      ri_cri   = 2._wp    / ( 2._wp + rn_ediss / rn_ediff )   ! resulting critical Richardson number
783      rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) )    ! resulting minimum length to recover molecular viscosity
784      !
785      IF(lwp) THEN                    !* Control print
786         WRITE(numout,*)
787         WRITE(numout,*) 'zdf_tke_init : tke turbulent closure scheme - initialisation'
788         WRITE(numout,*) '~~~~~~~~~~~~'
789         WRITE(numout,*) '   Namelist namzdf_tke : set tke mixing parameters'
790         WRITE(numout,*) '      coef. to compute avt                        rn_ediff  = ', rn_ediff
791         WRITE(numout,*) '      Kolmogoroff dissipation coef.               rn_ediss  = ', rn_ediss
792         WRITE(numout,*) '      tke surface input coef.                     rn_ebb    = ', rn_ebb
793         WRITE(numout,*) '      minimum value of tke                        rn_emin   = ', rn_emin
794         WRITE(numout,*) '      surface minimum value of tke                rn_emin0  = ', rn_emin0
795         WRITE(numout,*) '      background shear (>0)                       rn_bshear = ', rn_bshear
796         WRITE(numout,*) '      mixing length type                          nn_mxl    = ', nn_mxl
797         WRITE(numout,*) '      prandl number flag                          nn_pdl    = ', nn_pdl
798         WRITE(numout,*) '      surface mixing length = F(stress) or not    ln_mxl0   = ', ln_mxl0
799         WRITE(numout,*) '      surface  mixing length minimum value        rn_mxl0   = ', rn_mxl0
800         WRITE(numout,*) '      flag to take into acc.  Langmuir circ.      ln_lc     = ', ln_lc
801         WRITE(numout,*) '      coef to compute verticla velocity of LC     rn_lc     = ', rn_lc
802         WRITE(numout,*) '      test param. to add tke induced by wind      nn_etau   = ', nn_etau
803         WRITE(numout,*) '      flag for computation of exp. tke profile    nn_htau   = ', nn_htau
804         WRITE(numout,*) '      fraction of en which pene. the thermocline  rn_efr    = ', rn_efr
805         WRITE(numout,*) '      fraction of TKE added within the mixed layer by nn_etau rn_c    = ', rn_c
806         WRITE(numout,*)
807         WRITE(numout,*) '      critical Richardson nb with your parameters  ri_cri = ', ri_cri
808      ENDIF
809      !
810      !                              ! allocate tke arrays
811      IF( zdf_tke_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_init : unable to allocate arrays' )
812      !
813      !                               !* Check of some namelist values
814      IF( nn_mxl  < 0   .OR.  nn_mxl  > 3 )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' )
815      IF( nn_pdl  < 0   .OR.  nn_pdl  > 1 )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' )
816      IF( nn_htau < 0  .OR.  nn_htau > 5 )   CALL ctl_stop( 'bad flag: nn_htau is 0 to 5    ' )
817      IF( nn_etau == 3 .AND. .NOT. ln_cpl )   CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' )
818
819      IF( ln_mxl0 ) THEN
820         IF(lwp) WRITE(numout,*) '   use a surface mixing length = F(stress) :   set rn_mxl0 = rmxl_min'
821         rn_mxl0 = rmxl_min
822      ENDIF
823     
824      IF( nn_etau == 2  ) THEN
825          ierr = zdf_mxl_alloc()
826          nmln(:,:) = nlb10           ! Initialization of nmln
827      ENDIF
828
829      IF( nn_etau /= 0 .and. nn_htau == 2 ) THEN
830          ierr = zdf_mxl_alloc()
831          nmln(:,:) = nlb10           ! Initialization of nmln
832      ENDIF
833
834      !                               !* depth of penetration of surface tke
835      IF( nn_etau /= 0 ) THEN     
836         htau(:,:) = 0._wp
837         SELECT CASE( nn_htau )             ! Choice of the depth of penetration
838         CASE( 0 )                                 ! constant depth penetration (here 10 meters)
839            htau(:,:) = 10._wp
840         CASE( 1 )                                 ! F(latitude) : 0.5m to 30m poleward of 40 degrees
841            htau(:,:) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )           
842         CASE( 2 )                                 ! fraction of depth-integrated TKE within mixed-layer
843            rhtau = -1._wp / LOG( 1._wp - rn_c )
844         CASE( 3 )                                 ! F(latitude) : 0.5m to 15m poleward of 20 degrees
845            htau(:,:) = MAX(  0.5_wp, MIN( 15._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )
846         CASE( 4 )                                 ! F(latitude) : 0.5m to 10m/30m poleward of 13/40 degrees north/south
847            DO jj = 2, jpjm1
848               DO ji = fs_2, fs_jpim1   ! vector opt.
849                  IF( gphit(ji,jj) <= 0._wp ) THEN
850                     htau(ji,jj) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   )
851                  ELSE
852                     htau(ji,jj) = MAX(  0.5_wp, MIN( 10._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   )
853                  ENDIF
854               END DO
855            END DO
856         CASE ( 5 )                                ! F(latitude) : 0.5m to 10m poleward of 13 degrees north/south,
857            DO jj = 2, jpjm1                       !               10m to 30m between 30/45 degrees south
858               DO ji = fs_2, fs_jpim1   ! vector opt.
859                  IF( gphit(ji,jj) <= -30._wp ) THEN
860                     htau(ji,jj) = MAX(  10._wp, MIN( 30._wp, 55._wp* ABS( SIN( rpi/120._wp * ( gphit(ji,jj) + 23._wp ) ) ) )   )
861                  ELSE
862                     htau(ji,jj) = MAX(  0.5_wp, MIN( 10._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   )
863                  ENDIF
864               END DO
865            END DO
866         END SELECT
867         !
868         IF( nn_etau == 4 .AND. nn_htau /= 2 ) THEN            ! efr dependant on constant htau
869            DO jj = 2, jpjm1
870               DO ji = fs_2, fs_jpim1   ! vector opt.
871                  efr(ji,jj) = rn_efr / ( htau(ji,jj) * ( 1._wp - EXP( -bathy(ji,jj) / htau(ji,jj) ) ) )
872               END DO
873            END DO
874         ENDIF
875      ENDIF
876      !                               !* set vertical eddy coef. to the background value
877      DO jk = 1, jpk
878         avt (:,:,jk) = avtb(jk) * wmask (:,:,jk)
879         avm (:,:,jk) = avmb(jk) * wmask (:,:,jk)
880         avmu(:,:,jk) = avmb(jk) * wumask(:,:,jk)
881         avmv(:,:,jk) = avmb(jk) * wvmask(:,:,jk)
882      END DO
883      dissl(:,:,:) = 1.e-12_wp
884      !                             
885      CALL tke_rst( nit000, 'READ' )  !* read or initialize all required files
886      !
887   END SUBROUTINE zdf_tke_init
888
889
890   SUBROUTINE tke_rst( kt, cdrw )
891     !!---------------------------------------------------------------------
892     !!                   ***  ROUTINE tke_rst  ***
893     !!                     
894     !! ** Purpose :   Read or write TKE file (en) in restart file
895     !!
896     !! ** Method  :   use of IOM library
897     !!                if the restart does not contain TKE, en is either
898     !!                set to rn_emin or recomputed
899     !!----------------------------------------------------------------------
900     INTEGER         , INTENT(in) ::   kt     ! ocean time-step
901     CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
902     !
903     INTEGER ::   jit, jk   ! dummy loop indices
904     INTEGER ::   id1, id2, id3, id4, id5, id6   ! local integers
905     !!----------------------------------------------------------------------
906     !
907     IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
908        !                                   ! ---------------
909        IF( ln_rstart ) THEN                   !* Read the restart file
910           id1 = iom_varid( numror, 'en'   , ldstop = .FALSE. )
911           id2 = iom_varid( numror, 'avt'  , ldstop = .FALSE. )
912           id3 = iom_varid( numror, 'avm'  , ldstop = .FALSE. )
913           id4 = iom_varid( numror, 'avmu' , ldstop = .FALSE. )
914           id5 = iom_varid( numror, 'avmv' , ldstop = .FALSE. )
915           id6 = iom_varid( numror, 'dissl', ldstop = .FALSE. )
916           !
917           IF( id1 > 0 ) THEN                       ! 'en' exists
918              CALL iom_get( numror, jpdom_autoglo, 'en', en )
919              IF( MIN( id2, id3, id4, id5, id6 ) > 0 ) THEN        ! all required arrays exist
920                 CALL iom_get( numror, jpdom_autoglo, 'avt'  , avt   )
921                 CALL iom_get( numror, jpdom_autoglo, 'avm'  , avm   )
922                 CALL iom_get( numror, jpdom_autoglo, 'avmu' , avmu  )
923                 CALL iom_get( numror, jpdom_autoglo, 'avmv' , avmv  )
924                 CALL iom_get( numror, jpdom_autoglo, 'dissl', dissl )
925              ELSE                                                 ! one at least array is missing
926                 CALL tke_avn                                          ! compute avt, avm, avmu, avmv and dissl (approximation)
927              ENDIF
928           ELSE                                     ! No TKE array found: initialisation
929              IF(lwp) WRITE(numout,*) ' ===>>>> : previous run without tke scheme, en computed by iterative loop'
930              en (:,:,:) = rn_emin * tmask(:,:,:)
931              CALL tke_avn                               ! recompute avt, avm, avmu, avmv and dissl (approximation)
932              !
933              avt_k (:,:,:) = avt (:,:,:)
934              avm_k (:,:,:) = avm (:,:,:)
935              avmu_k(:,:,:) = avmu(:,:,:)
936              avmv_k(:,:,:) = avmv(:,:,:)
937              !
938              DO jit = nit000 + 1, nit000 + 10   ;   CALL zdf_tke( jit )   ;   END DO
939           ENDIF
940        ELSE                                   !* Start from rest
941           en(:,:,:) = rn_emin * tmask(:,:,:)
942           DO jk = 1, jpk                           ! set the Kz to the background value
943              avt (:,:,jk) = avtb(jk) * wmask (:,:,jk)
944              avm (:,:,jk) = avmb(jk) * wmask (:,:,jk)
945              avmu(:,:,jk) = avmb(jk) * wumask(:,:,jk)
946              avmv(:,:,jk) = avmb(jk) * wvmask(:,:,jk)
947           END DO
948        ENDIF
949        !
950     ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
951        !                                   ! -------------------
952        IF(lwp) WRITE(numout,*) '---- tke-rst ----'
953        CALL iom_rstput( kt, nitrst, numrow, 'en'   , en     )
954        CALL iom_rstput( kt, nitrst, numrow, 'avt'  , avt_k  )
955        CALL iom_rstput( kt, nitrst, numrow, 'avm'  , avm_k  )
956        CALL iom_rstput( kt, nitrst, numrow, 'avmu' , avmu_k )
957        CALL iom_rstput( kt, nitrst, numrow, 'avmv' , avmv_k )
958        CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl  )
959        !
960     ENDIF
961     !
962   END SUBROUTINE tke_rst
963
964#else
965   !!----------------------------------------------------------------------
966   !!   Dummy module :                                        NO TKE scheme
967   !!----------------------------------------------------------------------
968   LOGICAL, PUBLIC, PARAMETER ::   lk_zdftke = .FALSE.   !: TKE flag
969CONTAINS
970   SUBROUTINE zdf_tke_init           ! Dummy routine
971   END SUBROUTINE zdf_tke_init
972   SUBROUTINE zdf_tke( kt )          ! Dummy routine
973      WRITE(*,*) 'zdf_tke: You should not have seen this print! error?', kt
974   END SUBROUTINE zdf_tke
975   SUBROUTINE tke_rst( kt, cdrw )
976     CHARACTER(len=*) ::   cdrw
977     WRITE(*,*) 'tke_rst: You should not have seen this print! error?', kt, cdwr
978   END SUBROUTINE tke_rst
979#endif
980
981   !!======================================================================
982END MODULE zdftke
Note: See TracBrowser for help on using the repository browser.