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_OMP/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/UKMO/dev_r5518_GO6_package_OMP/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 9616

Last change on this file since 9616 was 9176, checked in by andmirek, 6 years ago

#2001: OMP directives

File size: 52.1 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, zcofa     !    -         -
230      REAL(wp) ::   ztau  , zdif, zdifa             !    -         -
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!$OMP PARALLEL DO
256         DO jj = 2, jpjm1            ! en(mikt(ji,jj))   = rn_emin
257            DO ji = fs_2, fs_jpim1   ! vector opt.
258               en(ji,jj,mikt(ji,jj))=rn_emin * tmask(ji,jj,1)
259            END DO
260         END DO
261      END IF
262!$OMP PARALLEL DO
263      DO jj = 2, jpjm1            ! en(1)   = rn_ebb taum / rau0  (min value rn_emin0)
264         DO ji = fs_2, fs_jpim1   ! vector opt.
265            en(ji,jj,1) = MAX( rn_emin0, zbbrau * taum(ji,jj) ) * tmask(ji,jj,1)
266         END DO
267      END DO
268     
269!!bfr   - start commented area
270      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
271      !                     !  Bottom boundary condition on tke
272      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
273      !
274      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
275      ! Tests to date have found the bottom boundary condition on tke to have very little effect.
276      ! The condition is coded here for completion but commented out until there is proof that the
277      ! computational cost is justified
278      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
279      !                     en(bot)   = (rn_ebb0/rau0)*0.5*sqrt(u_botfr^2+v_botfr^2) (min value rn_emin)
280!CDIR NOVERRCHK
281!!    DO jj = 2, jpjm1
282!CDIR NOVERRCHK
283!!       DO ji = fs_2, fs_jpim1   ! vector opt.
284!!          ztx2 = bfrua(ji-1,jj) * ub(ji-1,jj,mbku(ji-1,jj)) + &
285!!                 bfrua(ji  ,jj) * ub(ji  ,jj,mbku(ji  ,jj) )
286!!          zty2 = bfrva(ji,jj  ) * vb(ji,jj  ,mbkv(ji,jj  )) + &
287!!                 bfrva(ji,jj-1) * vb(ji,jj-1,mbkv(ji,jj-1) )
288!!          zebot = 0.001875_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 )   !  where 0.001875 = (rn_ebb0/rau0) * 0.5 = 3.75*0.5/1000.
289!!          en (ji,jj,mbkt(ji,jj)+1) = MAX( zebot, rn_emin ) * tmask(ji,jj,1)
290!!       END DO
291!!    END DO
292!!bfr   - end commented area
293      !
294      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
295      IF( ln_lc ) THEN      !  Langmuir circulation source term added to tke       (Axell JGR 2002)
296         !                  !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
297         !
298         !                        !* total energy produce by LC : cumulative sum over jk
299         zpelc(:,:,1) =  MAX( rn2b(:,:,1), 0._wp ) * fsdepw(:,:,1) * fse3w(:,:,1)
300!$OMP PARALLEL
301         DO jk = 2, jpk
302!$OMP DO
303            DO jj = 1, jpj
304                  zpelc(:,jj,jk)  = zpelc(:,jj,jk-1) + MAX( rn2b(:,jj,jk), 0._wp ) * fsdepw(:,jj,jk) * fse3w(:,jj,jk)
305            END DO
306!$OMP END DO
307         END DO
308!$OMP END PARALLEL
309         !                        !* finite Langmuir Circulation depth
310         zcof = 0.5 * 0.016 * 0.016 / ( zrhoa * zcdrag )
311         zcofa = 0.016 / SQRT( zrhoa * zcdrag )
312         imlc(:,:) = mbkt(:,:) + 1       ! Initialization to the number of w ocean point (=2 over land)
313!$OMP PARALLEL SHARED(imlc)
314         DO jk = jpkm1, 2, -1
315!$OMP DO PRIVATE(zus)
316            DO jj = 1, jpj               ! Last w-level at which zpelc>=0.5*us*us
317               DO ji = 1, jpi            !      with us=0.016*wind(starting from jpk-1)
318                  zus  = zcof * taum(ji,jj)
319                  IF( zpelc(ji,jj,jk) > zus )   imlc(ji,jj) = jk
320               END DO
321            END DO
322!$OMP END DO
323         END DO
324         !                               ! finite LC depth
325!$OMP DO
326         DO jj = 1, jpj 
327            DO ji = 1, jpi
328               zhlc(ji,jj) = fsdepw(ji,jj,imlc(ji,jj))
329            END DO
330         END DO
331!$OMP END DO
332!         zcof = 0.016 / SQRT( zrhoa * zcdrag )
333!$OMP DO PRIVATE(zus, zind, zwlc)
334         DO jk = 2, jpkm1         !* TKE Langmuir circulation source term added to en
335            DO jj = 2, jpjm1
336               DO ji = fs_2, fs_jpim1   ! vector opt.
337                  zus  = zcofa * SQRT( taum(ji,jj) )           ! Stokes drift
338                  !                                           ! vertical velocity due to LC
339                  zind = 0.5 - SIGN( 0.5, fsdepw(ji,jj,jk) - zhlc(ji,jj) )
340                  zwlc = zind * rn_lc * zus * SIN( rpi * fsdepw(ji,jj,jk) / zhlc(ji,jj) )
341                  !                                           ! TKE Langmuir circulation source term
342                  en(ji,jj,jk) = en(ji,jj,jk) + rdt * ( 1._wp - fr_i(ji,jj) ) * ( zwlc * zwlc * zwlc ) /   &
343                     &   zhlc(ji,jj) * wmask(ji,jj,jk) * tmask(ji,jj,1)
344               END DO
345            END DO
346         END DO
347!$OMP END DO
348!$OMP END PARALLEL
349         !
350      ENDIF
351      !
352      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
353      !                     !  Now Turbulent kinetic energy (output in en)
354      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
355      !                     ! Resolution of a tridiagonal linear system by a "methode de chasse"
356      !                     ! computation from level 2 to jpkm1  (e(1) already computed and e(jpk)=0 ).
357      !                     ! zdiag : diagonal zd_up : upper diagonal zd_lw : lower diagonal
358      !
359!$OMP PARALLEL DO
360      DO jk = 2, jpkm1           !* Shear production at uw- and vw-points (energy conserving form)
361         DO jj = 1, jpj                 ! here avmu, avmv used as workspace
362            DO ji = 1, jpi
363               avmu(ji,jj,jk) = avmu(ji,jj,jk) * ( un(ji,jj,jk-1) - un(ji,jj,jk) )   &
364                  &                            * ( ub(ji,jj,jk-1) - ub(ji,jj,jk) )   & 
365                  &                            / (  fse3uw_n(ji,jj,jk)               &
366                  &                              *  fse3uw_b(ji,jj,jk)  )
367               avmv(ji,jj,jk) = avmv(ji,jj,jk) * ( vn(ji,jj,jk-1) - vn(ji,jj,jk) )   &
368                  &                            * ( vb(ji,jj,jk-1) - vb(ji,jj,jk) )   &
369                  &                            / (  fse3vw_n(ji,jj,jk)               &
370                  &                              *  fse3vw_b(ji,jj,jk)  )
371            END DO
372         END DO
373      END DO
374      !
375!$OMP PARALLEL DO PRIVATE(zcof, zzd_up, zzd_lw, zesh2)
376      DO jk = 2, jpkm1           !* Matrix and right hand side in en
377         DO jj = 2, jpjm1
378            DO ji = fs_2, fs_jpim1   ! vector opt.
379               zcof   = zfact1 * tmask(ji,jj,jk)
380# if defined key_zdftmx_new
381               ! key_zdftmx_new: New internal wave-driven param: set a minimum value for Kz on TKE (ensure numerical stability)
382               zzd_up = zcof * ( MAX( avm(ji,jj,jk+1) + avm(ji,jj,jk), 2.e-5_wp ) )   &  ! upper diagonal
383                  &          / ( fse3t(ji,jj,jk  ) * fse3w(ji,jj,jk  ) )
384               zzd_lw = zcof * ( MAX( avm(ji,jj,jk) + avm(ji,jj,jk-1), 2.e-5_wp ) )   &  ! lower diagonal
385                  &          / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk  ) )
386# else
387               zzd_up = zcof * ( avm  (ji,jj,jk+1) + avm  (ji,jj,jk  ) )   &  ! upper diagonal
388                  &          / ( fse3t(ji,jj,jk  ) * fse3w(ji,jj,jk  ) )
389               zzd_lw = zcof * ( avm  (ji,jj,jk  ) + avm  (ji,jj,jk-1) )   &  ! lower diagonal
390                  &          / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk  ) )
391# endif
392                  !                                                           ! shear prod. at w-point weightened by mask
393               zesh2  =  ( avmu(ji-1,jj,jk) + avmu(ji,jj,jk) ) / MAX( 1._wp , umask(ji-1,jj,jk) + umask(ji,jj,jk) )   &
394                  &    + ( avmv(ji,jj-1,jk) + avmv(ji,jj,jk) ) / MAX( 1._wp , vmask(ji,jj-1,jk) + vmask(ji,jj,jk) )   
395                  !
396               zd_up(ji,jj,jk) = zzd_up            ! Matrix (zdiag, zd_up, zd_lw)
397               zd_lw(ji,jj,jk) = zzd_lw
398               zdiag(ji,jj,jk) = 1._wp - zzd_lw - zzd_up + zfact2 * dissl(ji,jj,jk) * tmask(ji,jj,jk)
399               !
400               !                                   ! right hand side in en
401               en(ji,jj,jk) = en(ji,jj,jk) + rdt * (  zesh2  -   avt(ji,jj,jk) * rn2(ji,jj,jk)    &
402                  &                                 + zfact3 * dissl(ji,jj,jk) * en (ji,jj,jk)  ) &
403                  &                                 * wmask(ji,jj,jk)
404            END DO
405         END DO
406      END DO
407      !                          !* Matrix inversion from level 2 (tke prescribed at level 1)
408!$OMP PARALLEL
409      DO jk = 3, jpkm1                             ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
410!$OMP DO
411         DO jj = 2, jpjm1
412            DO ji = fs_2, fs_jpim1    ! vector opt.
413               zdiag(ji,jj,jk) = zdiag(ji,jj,jk) - zd_lw(ji,jj,jk) * zd_up(ji,jj,jk-1) / zdiag(ji,jj,jk-1)
414            END DO
415         END DO
416!$OMP END DO
417      END DO
418!$OMP END PARALLEL
419      !
420      ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1
421!$OMP PARALLEL DO
422      DO jj = 2, jpjm1
423         DO ji = fs_2, fs_jpim1   ! vector opt.
424            zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1)    ! Surface boudary conditions on tke
425         END DO
426      END DO
427!$OMP PARALLEL
428      DO jk = 3, jpkm1
429!$OMP DO
430         DO jj = 2, jpjm1
431            DO ji = fs_2, fs_jpim1    ! vector opt.
432               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)
433            END DO
434         END DO
435!$OMP END DO
436      END DO
437!$OMP END PARALLEL
438      !
439      ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
440!$OMP PARALLEL DO
441      DO jj = 2, jpjm1
442         DO ji = fs_2, fs_jpim1   ! vector opt.
443            en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1)
444         END DO
445      END DO
446!$OMP PARALLEL
447      DO jk = jpk-2, 2, -1
448!$OMP DO
449         DO jj = 2, jpjm1
450            DO ji = fs_2, fs_jpim1    ! vector opt.
451               en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk)
452            END DO
453         END DO
454!$OMP END DO
455      END DO
456!$OMP END PARALLEL
457!$OMP PARALLEL DO
458      DO jk = 2, jpkm1                             ! set the minimum value of tke
459         DO jj = 2, jpjm1
460            DO ji = fs_2, fs_jpim1   ! vector opt.
461               en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * wmask(ji,jj,jk)
462            END DO
463         END DO
464      END DO
465
466      !                                 ! Save TKE prior to nn_etau addition 
467      e_niw(:,:,:) = en(:,:,:) 
468     
469      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
470      !                            !  TKE due to surface and internal wave breaking
471      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
472      IF( nn_htau == 2 ) THEN           !* mixed-layer depth dependant length scale
473!$OMP PARALLEL DO
474         DO jj = 2, jpjm1
475            DO ji = fs_2, fs_jpim1   ! vector opt.
476               htau(ji,jj) = rhtau * hmlp(ji,jj)
477            END DO
478         END DO
479      ENDIF
480#if defined key_iomput
481      !
482      CALL iom_put( "htau", htau(:,:) )  ! Check htau (even if constant in time)
483#endif
484      !
485      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction)
486!$OMP PARALLEL DO
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) + rn_efr * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
491                     &                                 * ( 1._wp - fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
492               END DO
493            END DO
494         END DO
495      ELSEIF( nn_etau == 2 ) THEN       !* act only at the base of the mixed layer (jk=nmln)  (rn_efr fraction)
496!$OMP PARALLEL DO PRIVATE(jk)
497         DO jj = 2, jpjm1
498            DO ji = fs_2, fs_jpim1   ! vector opt.
499               jk = nmln(ji,jj)
500               en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
501                  &                                 * ( 1._wp - fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
502            END DO
503         END DO
504      ELSEIF( nn_etau == 3 ) THEN       !* penetration belox the mixed layer (HF variability)
505!$OMP PARALLEL DO PRIVATE(ztx2, zty2, ztau, zdif, zdifa)
506         DO jk = 2, jpkm1
507            DO jj = 2, jpjm1
508               DO ji = fs_2, fs_jpim1   ! vector opt.
509                  ztx2 = utau(ji-1,jj  ) + utau(ji,jj)
510                  zty2 = vtau(ji  ,jj-1) + vtau(ji,jj)
511                  ztau = 0.5_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1)    ! module of the mean stress
512                  zdifa = taum(ji,jj) - ztau                            ! mean of modulus - modulus of the mean
513                  zdif = rhftau_scl * MAX( 0._wp, zdifa + rhftau_add )  ! apply some modifications...
514                  en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
515                     &                        * ( 1._wp - fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1)
516               END DO
517            END DO
518         END DO
519      ELSEIF( nn_etau == 4 ) THEN       !* column integral independant of htau (rn_efr must be scaled up)
520         IF( nn_htau == 2 ) THEN        ! efr dependant on time-varying htau
521!$OMP PARALLEL DO
522            DO jj = 2, jpjm1
523               DO ji = fs_2, fs_jpim1   ! vector opt.
524                  efr(ji,jj) = rn_efr / ( htau(ji,jj) * ( 1._wp - EXP( -bathy(ji,jj) / htau(ji,jj) ) ) )
525               END DO
526            END DO
527         ENDIF
528!$OMP PARALLEL DO
529         DO jk = 2, jpkm1
530            DO jj = 2, jpjm1
531               DO ji = fs_2, fs_jpim1   ! vector opt.
532                  en(ji,jj,jk) = en(ji,jj,jk) + efr(ji,jj) * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
533                     &                                                   * ( 1._wp - fr_i(ji,jj) )  * tmask(ji,jj,jk)
534               END DO
535            END DO
536         END DO
537      ENDIF
538      CALL lbc_lnk( en, 'W', 1. )      ! Lateral boundary conditions (sign unchanged)
539      !
540!$OMP PARALLEL DO
541      DO jk = 2, jpkm1                             ! TKE budget: near-inertial waves term 
542         DO jj = 2, jpjm1 
543            DO ji = fs_2, fs_jpim1   ! vector opt. 
544               e_niw(ji,jj,jk) = en(ji,jj,jk) - e_niw(ji,jj,jk) 
545            END DO 
546         END DO 
547      END DO 
548     
549      CALL lbc_lnk( e_niw, 'W', 1. ) 
550      !
551      CALL wrk_dealloc( jpi,jpj, imlc )    ! integer
552      CALL wrk_dealloc( jpi,jpj, zhlc ) 
553      CALL wrk_dealloc( jpi,jpj,jpk, zpelc, zdiag, zd_up, zd_lw ) 
554      !
555      IF( nn_timing == 1 )  CALL timing_stop('tke_tke')
556      !
557   END SUBROUTINE tke_tke
558
559
560   SUBROUTINE tke_avn
561      !!----------------------------------------------------------------------
562      !!                   ***  ROUTINE tke_avn  ***
563      !!
564      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
565      !!
566      !! ** Method  :   At this stage, en, the now TKE, is known (computed in
567      !!              the tke_tke routine). First, the now mixing lenth is
568      !!      computed from en and the strafification (N^2), then the mixings
569      !!      coefficients are computed.
570      !!              - Mixing length : a first evaluation of the mixing lengh
571      !!      scales is:
572      !!                      mxl = sqrt(2*en) / N 
573      !!      where N is the brunt-vaisala frequency, with a minimum value set
574      !!      to rmxl_min (rn_mxl0) in the interior (surface) ocean.
575      !!        The mixing and dissipative length scale are bound as follow :
576      !!         nn_mxl=0 : mxl bounded by the distance to surface and bottom.
577      !!                        zmxld = zmxlm = mxl
578      !!         nn_mxl=1 : mxl bounded by the e3w and zmxld = zmxlm = mxl
579      !!         nn_mxl=2 : mxl bounded such that the vertical derivative of mxl is
580      !!                    less than 1 (|d/dz(mxl)|<1) and zmxld = zmxlm = mxl
581      !!         nn_mxl=3 : mxl is bounded from the surface to the bottom usings
582      !!                    |d/dz(xml)|<1 to obtain lup, and from the bottom to
583      !!                    the surface to obtain ldown. the resulting length
584      !!                    scales are:
585      !!                        zmxld = sqrt( lup * ldown )
586      !!                        zmxlm = min ( lup , ldown )
587      !!              - Vertical eddy viscosity and diffusivity:
588      !!                      avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
589      !!                      avt = max( avmb, pdlr * avm ) 
590      !!      with pdlr=1 if nn_pdl=0, pdlr=1/pdl=F(Ri) otherwise.
591      !!
592      !! ** Action  : - avt : now vertical eddy diffusivity (w-point)
593      !!              - avmu, avmv : now vertical eddy viscosity at uw- and vw-points
594      !!----------------------------------------------------------------------
595      INTEGER  ::   ji, jj, jk   ! dummy loop indices
596      REAL(wp) ::   zrn2, zraug, zcoef, zav     ! local scalars
597      REAL(wp) ::   zdku, zpdlr, zri, zsqen     !   -      -
598      REAL(wp) ::   zdkv, zemxl, zemlm, zemlp   !   -      -
599      REAL(wp), POINTER, DIMENSION(:,:,:) :: zmpdl, zmxlm, zmxld
600      !!--------------------------------------------------------------------
601      !
602      IF( nn_timing == 1 )  CALL timing_start('tke_avn')
603
604      CALL wrk_alloc( jpi,jpj,jpk, zmpdl, zmxlm, zmxld ) 
605
606      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
607      !                     !  Mixing length
608      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
609      !
610      !                     !* Buoyancy length scale: l=sqrt(2*e/n**2)
611      !
612      ! initialisation of interior minimum value (avoid a 2d loop with mikt)
613      zmxlm(:,:,:)  = rmxl_min   
614      zmxld(:,:,:)  = rmxl_min
615      !
616      IF( ln_mxl0 ) THEN            ! surface mixing length = F(stress) : l=vkarmn*2.e5*taum/(rau0*g)
617!$OMP PARALLEL DO PRIVATE(zraug)
618         DO jj = 2, jpjm1
619            DO ji = fs_2, fs_jpim1
620               zraug = vkarmn * 2.e5_wp / ( rau0 * grav )
621               zmxlm(ji,jj,1) = MAX( rn_mxl0, zraug * taum(ji,jj) * tmask(ji,jj,1) )
622            END DO
623         END DO
624      ELSE
625         zmxlm(:,:,1) = rn_mxl0
626      ENDIF
627      !
628!$OMP PARALLEL DO PRIVATE(zrn2)
629      DO jk = 2, jpkm1              ! interior value : l=sqrt(2*e/n^2)
630         DO jj = 2, jpjm1
631            DO ji = fs_2, fs_jpim1   ! vector opt.
632               zrn2 = MAX( rn2(ji,jj,jk), rsmall )
633               zmxlm(ji,jj,jk) = MAX(  rmxl_min,  SQRT( 2._wp * en(ji,jj,jk) / zrn2 ) )
634            END DO
635         END DO
636      END DO
637      !
638      !                     !* Physical limits for the mixing length
639      !
640      zmxld(:,:,1  ) = zmxlm(:,:,1)   ! surface set to the minimum value
641      zmxld(:,:,jpk) = rmxl_min       ! last level  set to the minimum value
642      !
643      SELECT CASE ( nn_mxl )
644      !
645      ! where wmask = 0 set zmxlm == fse3w
646      CASE ( 0 )           ! bounded by the distance to surface and bottom
647!$OMP PARALLEL DO PRIVATE(zemxl)
648         DO jk = 2, jpkm1
649            DO jj = 2, jpjm1
650               DO ji = fs_2, fs_jpim1   ! vector opt.
651                  zemxl = MIN( fsdepw(ji,jj,jk) - fsdepw(ji,jj,mikt(ji,jj)), zmxlm(ji,jj,jk),   &
652                  &            fsdepw(ji,jj,mbkt(ji,jj)+1) - fsdepw(ji,jj,jk) )
653                  ! wmask prevent zmxlm = 0 if jk = mikt(ji,jj)
654                  zmxlm(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN(zmxlm(ji,jj,jk),fse3w(ji,jj,jk)) * (1 - wmask(ji,jj,jk))
655                  zmxld(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN(zmxlm(ji,jj,jk),fse3w(ji,jj,jk)) * (1 - wmask(ji,jj,jk))
656               END DO
657            END DO
658         END DO
659         !
660      CASE ( 1 )           ! bounded by the vertical scale factor
661!$OMP PARALLEL DO PRIVATE(zemxl)
662         DO jk = 2, jpkm1
663            DO jj = 2, jpjm1
664               DO ji = fs_2, fs_jpim1   ! vector opt.
665                  zemxl = MIN( fse3w(ji,jj,jk), zmxlm(ji,jj,jk) )
666                  zmxlm(ji,jj,jk) = zemxl
667                  zmxld(ji,jj,jk) = zemxl
668               END DO
669            END DO
670         END DO
671         !
672      CASE ( 2 )           ! |dk[xml]| bounded by e3t :
673!$OMP PARALLEL
674         DO jk = 2, jpkm1         ! from the surface to the bottom :
675!$OMP DO
676            DO jj = 2, jpjm1
677               DO ji = fs_2, fs_jpim1   ! vector opt.
678                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
679               END DO
680            END DO
681         END DO
682         DO jk = jpkm1, 2, -1     ! from the bottom to the surface :
683!$OMP DO PRIVATE(zemxl)
684            DO jj = 2, jpjm1
685               DO ji = fs_2, fs_jpim1   ! vector opt.
686                  zemxl = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
687                  zmxlm(ji,jj,jk) = zemxl
688                  zmxld(ji,jj,jk) = zemxl
689               END DO
690            END DO
691         END DO
692!$OMP END PARALLEL
693         !
694      CASE ( 3 )           ! lup and ldown, |dk[xml]| bounded by e3t :
695!$OMP PARALLEL
696         DO jk = 2, jpkm1         ! from the surface to the bottom : lup
697!$OMP DO
698            DO jj = 2, jpjm1
699               DO ji = fs_2, fs_jpim1   ! vector opt.
700                  zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
701               END DO
702            END DO
703         END DO
704         DO jk = jpkm1, 2, -1     ! from the bottom to the surface : ldown
705!$OMP DO
706            DO jj = 2, jpjm1
707               DO ji = fs_2, fs_jpim1   ! vector opt.
708                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
709               END DO
710            END DO
711         END DO
712!$OMP DO PRIVATE(zemlm, zemlp)
713         DO jk = 2, jpkm1
714            DO jj = 2, jpjm1
715               DO ji = fs_2, fs_jpim1   ! vector opt.
716                  zemlm = MIN ( zmxld(ji,jj,jk),  zmxlm(ji,jj,jk) )
717                  zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) )
718                  zmxlm(ji,jj,jk) = zemlm
719                  zmxld(ji,jj,jk) = zemlp
720               END DO
721            END DO
722         END DO
723!$OMP END PARALLEL
724         !
725      END SELECT
726      !
727# if defined key_c1d
728      e_dis(:,:,:) = zmxld(:,:,:)      ! c1d configuration : save mixing and dissipation turbulent length scales
729      e_mix(:,:,:) = zmxlm(:,:,:)
730# endif
731
732      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
733      !                     !  Vertical eddy viscosity and diffusivity  (avmu, avmv, avt)
734      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
735      DO jk = 1, jpkm1            !* vertical eddy viscosity & diffivity at w-points
736         DO jj = 2, jpjm1
737            DO ji = fs_2, fs_jpim1   ! vector opt.
738               zsqen = SQRT( en(ji,jj,jk) )
739               zav   = rn_ediff * zmxlm(ji,jj,jk) * zsqen
740               avm  (ji,jj,jk) = MAX( zav,                  avmb(jk) ) * wmask(ji,jj,jk)
741               avt  (ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
742               dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk)
743            END DO
744         END DO
745      END DO
746      CALL lbc_lnk( avm, 'W', 1. )      ! Lateral boundary conditions (sign unchanged)
747      !
748      DO jk = 2, jpkm1            !* vertical eddy viscosity at wu- and wv-points
749         DO jj = 2, jpjm1
750            DO ji = fs_2, fs_jpim1   ! vector opt.
751               avmu(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji+1,jj  ,jk) ) * wumask(ji,jj,jk)
752               avmv(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji  ,jj+1,jk) ) * wvmask(ji,jj,jk)
753            END DO
754         END DO
755      END DO
756      CALL lbc_lnk( avmu, 'U', 1. )   ;   CALL lbc_lnk( avmv, 'V', 1. )      ! Lateral boundary conditions
757      !
758      IF( nn_pdl == 1 ) THEN      !* Prandtl number case: update avt
759         DO jk = 2, jpkm1
760            DO jj = 2, jpjm1
761               DO ji = fs_2, fs_jpim1   ! vector opt.
762                  zcoef = avm(ji,jj,jk) * 2._wp * fse3w(ji,jj,jk) * fse3w(ji,jj,jk)
763                  !                                          ! shear
764                  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) )   &
765                    &  + avmu(ji  ,jj,jk) * ( un(ji  ,jj,jk-1) - un(ji  ,jj,jk) ) * ( ub(ji  ,jj,jk-1) - ub(ji  ,jj,jk) )
766                  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) )   &
767                    &  + avmv(ji,jj  ,jk) * ( vn(ji,jj  ,jk-1) - vn(ji,jj  ,jk) ) * ( vb(ji,jj  ,jk-1) - vb(ji,jj  ,jk) )
768                  !                                          ! local Richardson number
769                  zri   = MAX( rn2b(ji,jj,jk), 0._wp ) * zcoef / (zdku + zdkv + rn_bshear )
770                  zpdlr = MAX(  0.1_wp,  0.2 / MAX( 0.2 , zri )  )
771!!gm and even better with the use of the "true" ri_crit=0.22222...  (this change the results!)
772!!gm              zpdlr = MAX(  0.1_wp,  ri_crit / MAX( ri_crit , zri )  )
773                  avt(ji,jj,jk)   = MAX( zpdlr * avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
774# if defined key_c1d
775                  e_pdl(ji,jj,jk) = zpdlr * wmask(ji,jj,jk)  ! c1d configuration : save masked Prandlt number
776                  e_ric(ji,jj,jk) = zri   * wmask(ji,jj,jk)  ! c1d config. : save Ri
777# endif
778              END DO
779            END DO
780         END DO
781      ENDIF
782      CALL lbc_lnk( avt, 'W', 1. )                      ! Lateral boundary conditions on avt  (sign unchanged)
783
784      IF(ln_ctl) THEN
785         CALL prt_ctl( tab3d_1=en  , clinfo1=' tke  - e: ', tab3d_2=avt, clinfo2=' t: ', ovlap=1, kdim=jpk)
786         CALL prt_ctl( tab3d_1=avmu, clinfo1=' tke  - u: ', mask1=umask,                   &
787            &          tab3d_2=avmv, clinfo2=       ' v: ', mask2=vmask, ovlap=1, kdim=jpk )
788      ENDIF
789      !
790      CALL wrk_dealloc( jpi,jpj,jpk, zmpdl, zmxlm, zmxld ) 
791      !
792      IF( nn_timing == 1 )  CALL timing_stop('tke_avn')
793      !
794   END SUBROUTINE tke_avn
795
796
797   SUBROUTINE zdf_tke_init
798      !!----------------------------------------------------------------------
799      !!                  ***  ROUTINE zdf_tke_init  ***
800      !!                     
801      !! ** Purpose :   Initialization of the vertical eddy diffivity and
802      !!              viscosity when using a tke turbulent closure scheme
803      !!
804      !! ** Method  :   Read the namzdf_tke namelist and check the parameters
805      !!              called at the first timestep (nit000)
806      !!
807      !! ** input   :   Namlist namzdf_tke
808      !!
809      !! ** Action  :   Increase by 1 the nstop flag is setting problem encounter
810      !!----------------------------------------------------------------------
811      INTEGER ::   ji, jj, jk   ! dummy loop indices
812      INTEGER ::   ios, ierr
813      !!
814      NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb , rn_emin  ,   &
815         &                 rn_emin0, rn_bshear, nn_mxl , ln_mxl0  ,   &
816         &                 rn_mxl0 , nn_pdl   , ln_lc  , rn_lc    ,   &
817         &                 nn_etau , nn_htau  , rn_efr , rn_c   
818      !!----------------------------------------------------------------------
819
820      REWIND( numnam_ref )              ! Namelist namzdf_tke in reference namelist : Turbulent Kinetic Energy
821      READ  ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901)
822901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in reference namelist', lwp )
823
824      REWIND( numnam_cfg )              ! Namelist namzdf_tke in configuration namelist : Turbulent Kinetic Energy
825      READ  ( numnam_cfg, namzdf_tke, IOSTAT = ios, ERR = 902 )
826902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in configuration namelist', lwp )
827      IF(lwm) WRITE ( numond, namzdf_tke )
828      !
829      ri_cri   = 2._wp    / ( 2._wp + rn_ediss / rn_ediff )   ! resulting critical Richardson number
830# if defined key_zdftmx_new
831      ! key_zdftmx_new: New internal wave-driven param: specified value of rn_emin & rmxl_min are used
832      rn_emin  = 1.e-10_wp
833      rmxl_min = 1.e-03_wp
834      IF(lwp) THEN                  ! Control print
835         WRITE(numout,*)
836         WRITE(numout,*) 'zdf_tke_init :  New tidal mixing case: force rn_emin = 1.e-10 and rmxl_min = 1.e-3 '
837         WRITE(numout,*) '~~~~~~~~~~~~'
838      ENDIF
839# else
840      rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) )    ! resulting minimum length to recover molecular viscosity
841# endif
842      !
843      IF(lwp) THEN                    !* Control print
844         WRITE(numout,*)
845         WRITE(numout,*) 'zdf_tke_init : tke turbulent closure scheme - initialisation'
846         WRITE(numout,*) '~~~~~~~~~~~~'
847         WRITE(numout,*) '   Namelist namzdf_tke : set tke mixing parameters'
848         WRITE(numout,*) '      coef. to compute avt                        rn_ediff  = ', rn_ediff
849         WRITE(numout,*) '      Kolmogoroff dissipation coef.               rn_ediss  = ', rn_ediss
850         WRITE(numout,*) '      tke surface input coef.                     rn_ebb    = ', rn_ebb
851         WRITE(numout,*) '      minimum value of tke                        rn_emin   = ', rn_emin
852         WRITE(numout,*) '      surface minimum value of tke                rn_emin0  = ', rn_emin0
853         WRITE(numout,*) '      background shear (>0)                       rn_bshear = ', rn_bshear
854         WRITE(numout,*) '      mixing length type                          nn_mxl    = ', nn_mxl
855         WRITE(numout,*) '      prandl number flag                          nn_pdl    = ', nn_pdl
856         WRITE(numout,*) '      surface mixing length = F(stress) or not    ln_mxl0   = ', ln_mxl0
857         WRITE(numout,*) '      surface  mixing length minimum value        rn_mxl0   = ', rn_mxl0
858         WRITE(numout,*) '      flag to take into acc.  Langmuir circ.      ln_lc     = ', ln_lc
859         WRITE(numout,*) '      coef to compute verticla velocity of LC     rn_lc     = ', rn_lc
860         WRITE(numout,*) '      test param. to add tke induced by wind      nn_etau   = ', nn_etau
861         WRITE(numout,*) '      flag for computation of exp. tke profile    nn_htau   = ', nn_htau
862         WRITE(numout,*) '      fraction of en which pene. the thermocline  rn_efr    = ', rn_efr
863         WRITE(numout,*) '      fraction of TKE added within the mixed layer by nn_etau rn_c    = ', rn_c
864         WRITE(numout,*)
865         WRITE(numout,*) '      critical Richardson nb with your parameters  ri_cri = ', ri_cri
866      ENDIF
867      !
868      !                              ! allocate tke arrays
869      IF( zdf_tke_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_init : unable to allocate arrays' )
870      !
871      !                               !* Check of some namelist values
872      IF( nn_mxl  < 0   .OR.  nn_mxl  > 3 )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' )
873      IF( nn_pdl  < 0   .OR.  nn_pdl  > 1 )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' )
874      IF( nn_htau < 0  .OR.  nn_htau > 5 )   CALL ctl_stop( 'bad flag: nn_htau is 0 to 5    ' )
875      IF( nn_etau == 3 .AND. .NOT. ln_cpl )   CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' )
876
877      IF( ln_mxl0 ) THEN
878         IF(lwp) WRITE(numout,*) '   use a surface mixing length = F(stress) :   set rn_mxl0 = rmxl_min'
879         rn_mxl0 = rmxl_min
880      ENDIF
881     
882      IF( nn_etau == 2  ) THEN
883          ierr = zdf_mxl_alloc()
884          nmln(:,:) = nlb10           ! Initialization of nmln
885      ENDIF
886
887      IF( nn_etau /= 0 .and. nn_htau == 2 ) THEN
888          ierr = zdf_mxl_alloc()
889          nmln(:,:) = nlb10           ! Initialization of nmln
890      ENDIF
891
892      !                               !* depth of penetration of surface tke
893      IF( nn_etau /= 0 ) THEN     
894         htau(:,:) = 0._wp
895         SELECT CASE( nn_htau )             ! Choice of the depth of penetration
896         CASE( 0 )                                 ! constant depth penetration (here 10 meters)
897            htau(:,:) = 10._wp
898         CASE( 1 )                                 ! F(latitude) : 0.5m to 30m poleward of 40 degrees
899            htau(:,:) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )           
900         CASE( 2 )                                 ! fraction of depth-integrated TKE within mixed-layer
901            rhtau = -1._wp / LOG( 1._wp - rn_c )
902         CASE( 3 )                                 ! F(latitude) : 0.5m to 15m poleward of 20 degrees
903            htau(:,:) = MAX(  0.5_wp, MIN( 15._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )
904         CASE( 4 )                                 ! F(latitude) : 0.5m to 10m/30m poleward of 13/40 degrees north/south
905            DO jj = 2, jpjm1
906               DO ji = fs_2, fs_jpim1   ! vector opt.
907                  IF( gphit(ji,jj) <= 0._wp ) THEN
908                     htau(ji,jj) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   )
909                  ELSE
910                     htau(ji,jj) = MAX(  0.5_wp, MIN( 10._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   )
911                  ENDIF
912               END DO
913            END DO
914         CASE ( 5 )                                ! F(latitude) : 0.5m to 10m poleward of 13 degrees north/south,
915            DO jj = 2, jpjm1                       !               10m to 30m between 30/45 degrees south
916               DO ji = fs_2, fs_jpim1   ! vector opt.
917                  IF( gphit(ji,jj) <= -30._wp ) THEN
918                     htau(ji,jj) = MAX(  10._wp, MIN( 30._wp, 55._wp* ABS( SIN( rpi/120._wp * ( gphit(ji,jj) + 23._wp ) ) ) )   )
919                  ELSE
920                     htau(ji,jj) = MAX(  0.5_wp, MIN( 10._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   )
921                  ENDIF
922               END DO
923            END DO
924         END SELECT
925         !
926         IF( nn_etau == 4 .AND. nn_htau /= 2 ) THEN            ! efr dependant on constant htau
927            DO jj = 2, jpjm1
928               DO ji = fs_2, fs_jpim1   ! vector opt.
929                  efr(ji,jj) = rn_efr / ( htau(ji,jj) * ( 1._wp - EXP( -bathy(ji,jj) / htau(ji,jj) ) ) )
930               END DO
931            END DO
932         ENDIF
933      ENDIF
934      !                               !* set vertical eddy coef. to the background value
935!$OMP PARALLEL DO
936      DO jk = 1, jpk
937         avt (:,:,jk) = avtb(jk) * wmask (:,:,jk)
938         avm (:,:,jk) = avmb(jk) * wmask (:,:,jk)
939         avmu(:,:,jk) = avmb(jk) * wumask(:,:,jk)
940         avmv(:,:,jk) = avmb(jk) * wvmask(:,:,jk)
941      END DO
942      dissl(:,:,:) = 1.e-12_wp
943      !                             
944      CALL tke_rst( nit000, 'READ' )  !* read or initialize all required files
945      !
946   END SUBROUTINE zdf_tke_init
947
948
949   SUBROUTINE tke_rst( kt, cdrw )
950     !!---------------------------------------------------------------------
951     !!                   ***  ROUTINE tke_rst  ***
952     !!                     
953     !! ** Purpose :   Read or write TKE file (en) in restart file
954     !!
955     !! ** Method  :   use of IOM library
956     !!                if the restart does not contain TKE, en is either
957     !!                set to rn_emin or recomputed
958     !!----------------------------------------------------------------------
959     INTEGER         , INTENT(in) ::   kt     ! ocean time-step
960     CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
961     !
962     INTEGER ::   jit, jk   ! dummy loop indices
963     INTEGER ::   id1, id2, id3, id4, id5, id6   ! local integers
964     !!----------------------------------------------------------------------
965     !
966     IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
967        !                                   ! ---------------
968        IF( ln_rstart ) THEN                   !* Read the restart file
969           id1 = iom_varid( numror, 'en'   , ldstop = .FALSE. )
970           id2 = iom_varid( numror, 'avt'  , ldstop = .FALSE. )
971           id3 = iom_varid( numror, 'avm'  , ldstop = .FALSE. )
972           id4 = iom_varid( numror, 'avmu' , ldstop = .FALSE. )
973           id5 = iom_varid( numror, 'avmv' , ldstop = .FALSE. )
974           id6 = iom_varid( numror, 'dissl', ldstop = .FALSE. )
975           !
976           IF( id1 > 0 ) THEN                       ! 'en' exists
977              CALL iom_get( numror, jpdom_autoglo, 'en', en )
978              IF( MIN( id2, id3, id4, id5, id6 ) > 0 ) THEN        ! all required arrays exist
979                 CALL iom_get( numror, jpdom_autoglo, 'avt'  , avt   )
980                 CALL iom_get( numror, jpdom_autoglo, 'avm'  , avm   )
981                 CALL iom_get( numror, jpdom_autoglo, 'avmu' , avmu  )
982                 CALL iom_get( numror, jpdom_autoglo, 'avmv' , avmv  )
983                 CALL iom_get( numror, jpdom_autoglo, 'dissl', dissl )
984              ELSE                                                 ! one at least array is missing
985                 CALL tke_avn                                          ! compute avt, avm, avmu, avmv and dissl (approximation)
986              ENDIF
987           ELSE                                     ! No TKE array found: initialisation
988              IF(lwp) WRITE(numout,*) ' ===>>>> : previous run without tke scheme, en computed by iterative loop'
989              en (:,:,:) = rn_emin * tmask(:,:,:)
990              CALL tke_avn                               ! recompute avt, avm, avmu, avmv and dissl (approximation)
991              !
992              avt_k (:,:,:) = avt (:,:,:)
993              avm_k (:,:,:) = avm (:,:,:)
994              avmu_k(:,:,:) = avmu(:,:,:)
995              avmv_k(:,:,:) = avmv(:,:,:)
996              !
997              DO jit = nit000 + 1, nit000 + 10   ;   CALL zdf_tke( jit )   ;   END DO
998           ENDIF
999        ELSE                                   !* Start from rest
1000           en(:,:,:) = rn_emin * tmask(:,:,:)
1001!$OMP PARALLEL DO
1002           DO jk = 1, jpk                           ! set the Kz to the background value
1003              avt (:,:,jk) = avtb(jk) * wmask (:,:,jk)
1004              avm (:,:,jk) = avmb(jk) * wmask (:,:,jk)
1005              avmu(:,:,jk) = avmb(jk) * wumask(:,:,jk)
1006              avmv(:,:,jk) = avmb(jk) * wvmask(:,:,jk)
1007           END DO
1008        ENDIF
1009        !
1010     ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
1011        !                                   ! -------------------
1012        IF(lwp) WRITE(numout,*) '---- tke-rst ----'
1013        CALL iom_rstput( kt, nitrst, numrow, 'en'   , en     )
1014        CALL iom_rstput( kt, nitrst, numrow, 'avt'  , avt_k  )
1015        CALL iom_rstput( kt, nitrst, numrow, 'avm'  , avm_k  )
1016        CALL iom_rstput( kt, nitrst, numrow, 'avmu' , avmu_k )
1017        CALL iom_rstput( kt, nitrst, numrow, 'avmv' , avmv_k )
1018        CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl  )
1019        !
1020     ENDIF
1021     !
1022   END SUBROUTINE tke_rst
1023
1024#else
1025   !!----------------------------------------------------------------------
1026   !!   Dummy module :                                        NO TKE scheme
1027   !!----------------------------------------------------------------------
1028   LOGICAL, PUBLIC, PARAMETER ::   lk_zdftke = .FALSE.   !: TKE flag
1029CONTAINS
1030   SUBROUTINE zdf_tke_init           ! Dummy routine
1031   END SUBROUTINE zdf_tke_init
1032   SUBROUTINE zdf_tke( kt )          ! Dummy routine
1033      WRITE(*,*) 'zdf_tke: You should not have seen this print! error?', kt
1034   END SUBROUTINE zdf_tke
1035   SUBROUTINE tke_rst( kt, cdrw )
1036     CHARACTER(len=*) ::   cdrw
1037     WRITE(*,*) 'tke_rst: You should not have seen this print! error?', kt, cdwr
1038   END SUBROUTINE tke_rst
1039#endif
1040
1041   !!======================================================================
1042END MODULE zdftke
Note: See TracBrowser for help on using the repository browser.