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.
zdfkpp.F90 in branches/2011/dev_r2787_LOCEAN3_TRA_TRP/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2011/dev_r2787_LOCEAN3_TRA_TRP/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfkpp.F90 @ 2789

Last change on this file since 2789 was 2789, checked in by cetlod, 13 years ago

Implementation of the merge of TRA/TRP : first guess, see ticket #842

  • Property svn:keywords set to Id
File size: 79.2 KB
Line 
1MODULE zdfkpp
2   !!======================================================================
3   !!                       ***  MODULE  zdfkpp  ***
4   !! Ocean physics:  vertical mixing coefficient compute from the KPP
5   !!                 turbulent closure parameterization
6   !!=====================================================================
7   !! History :  OPA  ! 2000-03 (W.G. Large, J. Chanut) Original code
8   !!            8.1  ! 2002-06 (J.M. Molines) for real case CLIPPER 
9   !!            8.2  ! 2003-10 (Chanut J.) re-writting
10   !!   NEMO     1.0  ! 2005-01 (C. Ethe, G. Madec) Free form, F90 + creation of tra_kpp routine
11   !!            3.3  ! 2010-10 (C. Ethe, G. Madec) reorganisation of initialisation phase + merge TRC-TRA
12   !!----------------------------------------------------------------------
13#if defined key_zdfkpp   ||   defined key_esopa
14   !!----------------------------------------------------------------------
15   !!   'key_zdfkpp'                                             KPP scheme
16   !!----------------------------------------------------------------------
17   !!   zdf_kpp      : update momentum and tracer Kz from a kpp scheme
18   !!   zdf_kpp_init : initialization, namelist read, and parameters control
19   !!   tra_kpp      : compute and add to the T & S trend the non-local flux
20   !!   trc_kpp      : compute and add to the passive tracer trend the non-local flux (lk_top=T)
21   !!----------------------------------------------------------------------
22   USE oce             ! ocean dynamics and active tracers
23   USE dom_oce         ! ocean space and time domain
24   USE zdf_oce         ! ocean vertical physics
25   USE sbc_oce         ! surface boundary condition: ocean
26   USE phycst          ! physical constants
27   USE eosbn2          ! equation of state
28   USE zdfddm          ! double diffusion mixing
29   USE in_out_manager  ! I/O manager
30   USE lib_mpp         ! MPP library
31   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
32   USE prtctl          ! Print control
33   USE trdmod_oce      ! ocean trends definition
34   USE trdtra          ! tracers trends
35
36   IMPLICIT NONE
37   PRIVATE
38
39   PUBLIC   zdf_kpp       ! routine called by step.F90
40   PUBLIC   zdf_kpp_init  ! routine called by opa.F90
41   PUBLIC   tra_kpp       ! routine called by step.F90
42   PUBLIC   trc_kpp       ! routine called by trcstp.F90
43
44   LOGICAL , PUBLIC, PARAMETER ::   lk_zdfkpp = .TRUE.    !: KPP vertical mixing flag
45
46   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   ghats    !: non-local scalar mixing term (gamma/<ws>o)
47   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   wt0      !: surface temperature flux for non local flux
48   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   ws0      !: surface salinity flux for non local flux
49   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   hkpp     !: boundary layer depth
50
51   !                                        !!* Namelist namzdf_kpp *
52   REAL(wp) ::   rn_difmiw  =  1.2e-04_wp    ! constant internal wave viscosity (m2/s)
53   REAL(wp) ::   rn_difsiw  =  1.2e-05_wp    ! constant internal wave diffusivity (m2/s)
54   REAL(wp) ::   rn_riinfty =  0.8_wp        ! local Richardson Number limit for shear instability
55   REAL(wp) ::   rn_difri   =  5.e-03_wp     ! maximum shear mixing at Rig = 0    (m2/s)
56   REAL(wp) ::   rn_bvsqcon = -1.e-09_wp     ! Brunt-Vaisala squared (1/s**2) for maximum convection
57   REAL(wp) ::   rn_difcon  =  1._wp         ! maximum mixing in interior convection (m2/s)
58   INTEGER  ::   nn_ave     =  1             ! = 0/1 flag for horizontal average on avt, avmu, avmv
59
60#if defined key_zdfddm
61   !                                        !!! ** Double diffusion Mixing
62   REAL(wp) ::   difssf  = 1.e-03_wp         ! maximum salt fingering mixing
63   REAL(wp) ::   Rrho0   = 1.9_wp            ! limit for salt  fingering mixing
64   REAL(wp) ::   difsdc  = 1.5e-06_wp        ! maximum diffusive convection mixing
65#endif
66   LOGICAL  ::   ln_kpprimix  = .TRUE.       ! Shear instability mixing
67
68   !                                        !!! ** General constants  **
69   REAL(wp) ::   epsln   = 1.0e-20_wp        ! a small positive number   
70   REAL(wp) ::   pthird  = 1._wp/3._wp       ! 1/3
71   REAL(wp) ::   pfourth = 1._wp/4._wp       ! 1/4
72
73   !                                        !!! ** Boundary Layer Turbulence Parameters  **
74   REAL(wp) ::   vonk     = 0.4_wp           ! von Karman's constant
75   REAL(wp) ::   epsilon  = 0.1_wp           ! nondimensional extent of the surface layer
76   REAL(wp) ::   rconc1   = 5.0_wp           ! standard flux profile function parmaeters
77   REAL(wp) ::   rconc2   = 16.0_wp          !         "        "
78   REAL(wp) ::   rconcm   = 8.38_wp          ! momentum flux profile fit
79   REAL(wp) ::   rconam   = 1.26_wp          !         "       "
80   REAL(wp) ::   rzetam   = -.20_wp          !         "       "       
81   REAL(wp) ::   rconcs   = 98.96_wp         !  scalar  flux profile fit
82   REAL(wp) ::   rconas   = -28.86_wp        !         "       "
83   REAL(wp) ::   rzetas   = -1.0_wp          !         "       " 
84   
85   !                                        !!! ** Boundary Layer Depth Diagnostic  **
86   REAL(wp) ::   Ricr     = 0.3_wp           ! critical bulk Richardson Number
87   REAL(wp) ::   rcekman  = 0.7_wp           ! coefficient for ekman depth 
88   REAL(wp) ::   rcmonob  = 1.0_wp           ! coefficient for Monin-Obukhov depth
89   REAL(wp) ::   rconcv   = 1.7_wp           ! ratio of interior buoyancy frequency to its value at entrainment depth
90   REAL(wp) ::   hbf      = 1.0_wp           ! fraction of bound. layer depth to which absorbed solar
91      !                                      ! rad. and contributes to surf. buo. forcing
92   REAL(wp) ::   Vtc                         ! function of rconcv,rconcs,epsilon,vonk,Ricr
93   
94   !                                        !!! ** Nonlocal Boundary Layer Mixing **
95   REAL(wp) ::   rcstar   = 5.0_wp           ! coefficient for convective nonlocal transport
96   REAL(wp) ::   rcs      = 1.0e-3_wp        ! conversion: mm/s ==> m/s   
97   REAL(wp) ::   rcg                         ! non-dimensional coefficient for nonlocal transport
98
99#if ! defined key_kppcustom
100   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   del     ! array for reference mean values of vertical integration
101#endif
102
103#if defined key_kpplktb
104   !                                         !!! ** Parameters for lookup table for turbulent velocity scales **
105   INTEGER, PARAMETER ::   nilktb   = 892     ! number of values for zehat in KPP lookup table
106   INTEGER, PARAMETER ::   njlktb   = 482     ! number of values for ustar in KPP lookup table
107   INTEGER, PARAMETER ::   nilktbm1 = nilktb-1   !
108   INTEGER, PARAMETER ::   njlktbm1 = njlktb-1   !
109
110   REAL(wp), DIMENSION(nilktb,njlktb) ::   wmlktb   ! lookup table for the turbulent vertical velocity scale (momentum)
111   REAL(wp), DIMENSION(nilktb,njlktb) ::   wslktb   ! lookup table for the turbulent vertical velocity scale (tracers)
112
113   REAL(wp) ::   dehatmin = -4.e-7_wp    ! minimum limit for zhat in lookup table (m3/s3)
114   REAL(wp) ::   dehatmax = 0._wp        ! maximum limit for zhat in lookup table (m3/s3)
115   REAL(wp) ::   ustmin   = 0._wp        ! minimum limit for ustar in lookup table (m/s)
116   REAL(wp) ::   ustmax   = 0.04_wp      ! maximum limit for ustar in lookup table (m/s)   
117   REAL(wp) ::   dezehat                 ! delta zhat in lookup table
118   REAL(wp) ::   deustar                 ! delta ustar in lookup table
119#endif
120   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:) ::   ratt   ! attenuation coef  (already defines in module traqsr,
121   !                                    ! but only if the solar radiation penetration is considered)
122   
123   !                                    !!! * penetrative solar radiation coefficient *
124   REAL(wp) ::   rabs = 0.58_wp          ! fraction associated with xsi1
125   REAL(wp) ::   xsi1 = 0.35_wp          ! first depth of extinction
126   REAL(wp) ::   xsi2 = 23.0_wp          ! second depth of extinction
127      !                           ! (default values: water type Ib)
128
129   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   etmean, eumean, evmean   ! coeff. used for hor. smoothing at t-, u- & v-points
130       
131#if defined key_c1d
132   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   rig    !: gradient Richardson number
133   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   rib    !: bulk Richardson number
134   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   buof   !: buoyancy forcing
135   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   mols   !: moning-Obukhov length scale
136   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   ekdp   !: Ekman depth
137#endif
138
139   INTEGER  ::   jip = 62 , jjp = 111
140
141   !! * Substitutions
142#  include "domzgr_substitute.h90"
143#  include "vectopt_loop_substitute.h90"
144#  include  "zdfddm_substitute.h90"
145   !!----------------------------------------------------------------------
146   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
147   !! $Id$
148   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
149   !!----------------------------------------------------------------------
150CONTAINS
151
152   INTEGER FUNCTION zdf_kpp_alloc()
153      !!----------------------------------------------------------------------
154      !!                 ***  FUNCTION zdf_kpp_alloc  ***
155      !!----------------------------------------------------------------------
156      ALLOCATE( ghats(jpi,jpj,jpk), wt0(jpi,jpj), ws0(jpi,jpj), hkpp(jpi,jpj), &
157#if ! defined key_kpplktb
158         &      del(jpk,jpk),                                                  &
159#endif
160         &      ratt(jpk),                                                     &
161         &      etmean(jpi,jpj,jpk), eumean(jpi,jpj,jpk), evmean(jpi,jpj,jpk), &
162#if defined key_c1d
163         &      rig (jpi,jpj,jpk), rib(jpi,jpj,jpk), buof(jpi,jpj,jpk),        &
164         &      mols(jpi,jpj,jpk), ekdp(jpi,jpj),                              &
165#endif
166         &      STAT= zdf_kpp_alloc )
167         !
168      IF( lk_mpp             )   CALL mpp_sum ( zdf_kpp_alloc )
169      IF( zdf_kpp_alloc /= 0 )   CALL ctl_warn('zdf_kpp_alloc: failed to allocate arrays')
170   END FUNCTION zdf_kpp_alloc
171
172
173   SUBROUTINE zdf_kpp( kt )
174      !!----------------------------------------------------------------------
175      !!                   ***  ROUTINE zdf_kpp  ***
176      !!
177      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
178      !!      coefficients and non local mixing using K-profile parameterization
179      !!
180      !! ** Method :   The boundary layer depth hkpp is diagnosed at tracer points
181      !!      from profiles of buoyancy, and shear, and the surface forcing.
182      !!      Above hbl (sigma=-z/hbl <1) the mixing coefficients are computed from
183      !!
184      !!                      Kx =  hkpp  Wx(sigma) G(sigma) 
185      !!
186      !!             and the non local term ghat = Cs / Ws(sigma) / hkpp
187      !!      Below hkpp  the coefficients are the sum of mixing due to internal waves
188      !!      shear instability and double diffusion.
189      !!
190      !!      -1- Compute the now interior vertical mixing coefficients at all depths.
191      !!      -2- Diagnose the boundary layer depth.
192      !!      -3- Compute the now boundary layer vertical mixing coefficients.
193      !!      -4- Compute the now vertical eddy vicosity and diffusivity.
194      !!      -5- Smoothing
195      !!
196      !!        N.B. The computation is done from jk=2 to jpkm1
197      !!             Surface value of avt avmu avmv are set once a time to zero
198      !!             in routine zdf_kpp_init.
199      !!
200      !! ** Action  :   update the non-local terms ghats
201      !!                update avt, avmu, avmv (before vertical eddy coef.)
202      !!
203      !! References : Large W.G., Mc Williams J.C. and Doney S.C.             
204      !!         Reviews of Geophysics, 32, 4, November 1994
205      !!         Comments in the code refer to this paper, particularly
206      !!         the equation number. (LMD94, here after)
207      !!----------------------------------------------------------------------
208      USE oce     , zviscos => ua   ! temp. array for viscosities use ua as workspace
209      USE oce     , zdiffut => va   ! temp. array for diffusivities use sa as workspace
210      USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released, wrk_in_use_xz, wrk_not_released_xz
211      USE wrk_nemo, ONLY: zBo    => wrk_2d_1, &  ! Surface buoyancy forcing,
212                          zBosol => wrk_2d_2, &  ! friction velocity
213                          zustar => wrk_2d_3
214      USE wrk_nemo, ONLY: zmask  => wrk_2d_4
215!gm      USE wrk_nemo, ONLY: wrk_2d_5, wrk_2d_6, wrk_2d_7, wrk_2d_8, wrk_2d_9, &
216      USE wrk_nemo, ONLY:           wrk_2d_6, wrk_2d_7, wrk_2d_8, wrk_2d_9, &
217                          wrk_2d_10,wrk_2d_11
218      USE wrk_nemo, ONLY: wrk_1d_1,  wrk_1d_2,  wrk_1d_3,  wrk_1d_4,  &
219                          wrk_1d_5,  wrk_1d_6,  wrk_1d_7,  wrk_1d_8,  &
220                          wrk_1d_9,  wrk_1d_10, wrk_1d_11, wrk_1d_12, &
221                          wrk_1d_13, wrk_1d_14
222      USE wrk_nemo, ONLY: zblcm => wrk_xz_1, &   ! Boundary layer
223                          zblct => wrk_xz_2      !  diffusivities/viscosities
224#if defined key_zdfddm
225      USE wrk_nemo, ONLY: zdiffus => wrk_3d_1
226      USE wrk_nemo, ONLY: zblcs   => wrk_xz_3
227#endif
228      !!
229      INTEGER, INTENT( in  ) ::   kt   ! ocean time step
230      !!
231      INTEGER ::   ji, jj, jk              ! dummy loop indices
232      INTEGER ::   ikbot, jkmax, jkm1, jkp2   !
233
234      REAL(wp) ::   ztx, zty, zflageos, zstabl, zbuofdep,zucube     !
235      REAL(wp) ::   zrhos, zalbet, zbeta, zthermal, zhalin, zatt1   !
236      REAL(wp) ::   zref, zt, zs, zh, zu, zv, zrh                   ! Bulk richardson number
237      REAL(wp) ::   zrib, zrinum, zdVsq, zVtsq                      !
238      REAL(wp) ::   zehat, zeta, zhrib, zsig, zscale, zwst, zws, zwm   ! Velocity scales
239#if defined key_kpplktb
240      INTEGER ::    il, jl                                          ! Lookup table or Analytical functions
241      REAL(wp) ::   ud, zfrac, ufrac, zwam, zwbm, zwas, zwbs        !
242#else
243      REAL(wp) ::   zwsun, zwmun, zcons, zconm, zwcons, zwconm      !
244#endif
245      REAL(wp) ::   zsr, zbw, ze, zb, zd, zc, zaw, za, zb1, za1, zkw, zk0, zcomp , zrhd,zrhdr,zbvzed   ! In situ density
246#if ! defined key_kppcustom     
247      INTEGER  ::   jm                          ! dummy loop indices
248      REAL(wp) ::   zr1, zr2, zr3, zr4, zrhop   ! Compression terms
249#endif
250      REAL(wp) ::   zflag, ztemp, zrn2, zdep21, zdep32, zdep43
251      REAL(wp) ::   zdku2, zdkv2, ze3sqr, zsh2, zri, zfri          ! Interior richardson mixing
252!gm      REAL(wp), POINTER, DIMENSION(:,:) ::   zmoek                 ! Moning-Obukov limitation
253      REAL(wp), DIMENSION(jpi,0:2) ::   zmoek                 ! Moning-Obukov limitation
254      REAL(wp), POINTER, DIMENSION(:)   ::   zmoa, zekman               
255      REAL(wp)                          ::   zmob, zek
256      REAL(wp), POINTER, DIMENSION(:,:) ::   zdepw, zdift, zvisc   ! The pipe
257      REAL(wp), POINTER, DIMENSION(:,:) ::   zdept
258      REAL(wp), POINTER, DIMENSION(:,:) ::   zriblk
259      REAL(wp), POINTER, DIMENSION(:)   ::   zhmax, zria, zhbl 
260      REAL(wp) ::   zflagri, zflagek, zflagmo, zflagh, zflagkb   !
261      REAL(wp), POINTER, DIMENSION(:)   ::   za2m, za3m, zkmpm, za2t, za3t, zkmpt   ! Shape function (G)
262      REAL(wp) ::   zdelta, zdelta2, zdzup, zdzdn, zdzh, zvath, zgat1, zdat1, zkm1m, zkm1t
263#if defined key_zdfddm
264      REAL(wp) ::   zrrau, zds, zavdds, zavddt,zinr   ! double diffusion mixing
265      REAL(wp), POINTER, DIMENSION(:,:) ::     zdifs
266      REAL(wp), POINTER, DIMENSION(:)   ::   za2s, za3s, zkmps
267      REAL(wp) ::                            zkm1s
268#endif
269      !!--------------------------------------------------------------------
270     
271      IF( wrk_in_use(1, 1,2,3,4,5,6,7,8,9,10,11,12,13,14) .OR.   &
272          wrk_in_use(2, 1,2,3,4,5,6,7,8,9,10,11)          .OR.   &
273          wrk_in_use(3, 1)                                .OR. &
274          wrk_in_use_xz(1,2,3)                              ) THEN
275         CALL ctl_stop('zdf_kpp : requested workspace arrays unavailable.')   ;   RETURN
276      ENDIF
277      ! Set-up pointers to 2D spaces
278!gm      zmoek(1:jpi,0:2) => wrk_2d_5(1:jpi,1:3)
279      zdepw => wrk_2d_6(:,1:4)
280      zdift => wrk_2d_7(:,1:4)
281      zvisc => wrk_2d_8(:,1:4)
282      zdept => wrk_2d_9(:,1:3)
283      zriblk => wrk_2d_10(:,1:2)
284      ! 1D spaces
285      zmoa   => wrk_1d_1(1:jpi)
286      zekman => wrk_1d_2(1:jpi)
287      zhmax  => wrk_1d_3(1:jpi)
288      zria   => wrk_1d_4(1:jpi)
289      zhbl   => wrk_1d_5(1:jpi)
290      za2m   => wrk_1d_6(1:jpi)
291      za3m   => wrk_1d_7(1:jpi)
292      zkmpm  => wrk_1d_8(1:jpi)
293      za2t   => wrk_1d_9(1:jpi)
294      za3t   => wrk_1d_10(1:jpi)
295      zkmpt  => wrk_1d_11(1:jpi)
296#if defined key_zdfddm
297      zdifs => wrk_2d_11(:,1:4)
298      za2s  => wrk_1d_12(1:jpi)
299      za3s  => wrk_1d_13(1:jpi)
300      zkmps => wrk_1d_14(1:jpi)
301#endif
302
303      zviscos(:,:,:) = 0.
304      zblcm  (:,:  ) = 0. 
305      zdiffut(:,:,:) = 0.
306      zblct  (:,:  ) = 0. 
307#if defined key_zdfddm
308      zdiffus(:,:,:) = 0.
309      zblcs  (:,:  ) = 0. 
310#endif
311      ghats(:,:,:) = 0.
312     
313      zBo   (:,:) = 0.
314      zBosol(:,:) = 0.
315      zustar(:,:) = 0.
316
317
318      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
319      ! I. Interior diffusivity and viscosity at w points ( T interfaces)
320      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
321      DO jk = 2, jpkm1
322         DO jj = 2, jpjm1
323            DO ji = fs_2, fs_jpim1 
324               ! Mixing due to internal waves breaking
325               ! -------------------------------------
326               avmu(ji,jj,jk)  = rn_difmiw 
327               avt (ji,jj,jk)  = rn_difsiw             
328               ! Mixing due to vertical shear instability
329               ! -------------------------------------               
330               IF( ln_kpprimix ) THEN         
331                  ! Compute the gradient Richardson  number at interfaces (zri):
332                  ! LMD94, eq. 27 (is vertical smoothing needed : Rig=N^2 / (dz(u))^2 + (dz(v))^2
333                  zdku2 =   ( un(ji - 1,jj,jk - 1) - un(ji - 1,jj,jk) ) &
334                     &    * ( un(ji - 1,jj,jk - 1) - un(ji - 1,jj,jk) ) &
335                     &    + ( un(ji    ,jj,jk - 1) - un(ji    ,jj,jk) ) &
336                     &    * ( un(ji    ,jj,jk - 1) - un(ji    ,jj,jk) ) 
337                 
338                  zdkv2 =   ( vn(ji,jj - 1,jk - 1) - vn(ji,jj - 1,jk) ) &
339                     &    * ( vn(ji,jj - 1,jk - 1) - vn(ji,jj - 1,jk) ) &
340                     &    + ( vn(ji,    jj,jk - 1) - vn(ji,    jj,jk) ) &
341                     &    * ( vn(ji,    jj,jk - 1) - vn(ji,    jj,jk) ) 
342
343                  ze3sqr = 1. / ( fse3w(ji,jj,jk) * fse3w(ji,jj,jk) )
344                  ! Square of vertical shear  at interfaces
345                  zsh2   = 0.5 * ( zdku2 + zdkv2 ) * ze3sqr
346                  zri    = MAX( rn2(ji,jj,jk), 0. ) / ( zsh2 + epsln ) 
347#if defined key_c1d
348                  ! save the gradient richardson number
349                  rig(ji,jj,jk) = zri * tmask(ji,jj,jk)
350#endif                 
351                  ! Evaluate f of Ri (zri) for shear instability store in zfri
352                  ! LMD94, eq. 28a,b,c, figure 3 ; Rem: p1 is 3, hard coded
353                  zfri  = MAX( zri , 0. )
354                  zfri  = MIN( zfri / rn_riinfty , 1.0 )
355                  zfri  = ( 1.0 - zfri * zfri )
356                  zfri  = zfri * zfri  * zfri
357                  ! add shear contribution to mixing coef.
358                  avmu(ji,jj,jk) =  avmu(ji,jj,jk) + rn_difri * zfri   
359                  avt (ji,jj,jk) =  avt (ji,jj,jk) + rn_difri * zfri   
360               ENDIF
361#if defined key_zdfddm 
362               avs (ji,jj,jk) =  avt (ji,jj,jk)             
363               !  Double diffusion mixing ; NOT IN ROUTINE ZDFDDM.F90
364               ! ------------------------------------------------------------------
365               ! only retains positive value of rrau
366               zrrau = MAX( rrau(ji,jj,jk), epsln )
367               zds   = tsn(ji,jj,jk-1,jp_sal) - tsn(ji,jj,jk,jp_sal)
368               IF( zrrau > 1. .AND. zds > 0.) THEN
369                  !
370                  ! Salt fingering case.
371                  !---------------------
372                  ! Compute interior diffusivity for double diffusive mixing of
373                  ! salinity. Upper bound "zrrau" by "Rrho0"; (Rrho0=1.9, difcoefnuf=0.001).
374                  ! After that set interior diffusivity for double diffusive mixing
375                  ! of temperature
376                  zavdds = MIN( zrrau, Rrho0 )
377                  zavdds = ( zavdds - 1.0 ) / ( Rrho0 - 1.0 )
378                  zavdds = 1.0 - zavdds * zavdds 
379                  zavdds = zavdds * zavdds * zavdds 
380                  zavdds = difssf * zavdds 
381                  zavddt = 0.7 * zavdds
382               ELSEIF( zrrau < 1. .AND. zrrau > 0. .AND. zds < 0.) THEN
383                  !
384                  ! Diffusive convection case.
385                  !---------------------------
386                  ! Compute interior diffusivity for double diffusive mixing of
387                  ! temperature (Marmorino and Caldwell, 1976);
388                  ! Compute interior diffusivity for double diffusive mixing of salinity
389                  zinr   = 1. / zrrau
390                  zavddt = 0.909 * EXP( 4.6 * EXP( -0.54* ( zinr - 1. ) ) ) 
391                  zavddt = difsdc * zavddt
392                  IF( zrrau < 0.5) THEN
393                     zavdds = zavddt * 0.15 * zrrau
394                  ELSE
395                     zavdds = zavddt * (1.85 * zrrau - 0.85 ) 
396                  ENDIF
397               ELSE
398                  zavddt = 0.
399                  zavdds = 0.
400               ENDIF 
401               ! Add double diffusion contribution to temperature and salinity  mixing coefficients.
402               avt (ji,jj,jk) =  avt (ji,jj,jk) +  zavddt 
403               avs (ji,jj,jk) =  avs (ji,jj,jk) +  zavdds         
404#endif                     
405            END DO
406         END DO
407      END DO
408
409
410      ! Radiative (zBosol) and non radiative (zBo) surface buoyancy
411      !JMM at the time zdfkpp is called, q still holds the sum q + qsr
412      !---------------------------------------------------------------------
413      DO jj = 2, jpjm1
414         DO ji = fs_2, fs_jpim1     
415            IF( nn_eos < 1) THEN   
416               zt     = tsn(ji,jj,1,jp_tem)
417               zs     = tsn(ji,jj,1,jp_sal) - 35.0
418               zh     = fsdept(ji,jj,1)
419               !  potential volumic mass
420               zrhos  = rhop(ji,jj,1)
421               zalbet = ( ( ( - 0.255019e-07 * zt + 0.298357e-05 ) * zt   &   ! ratio alpha/beta
422                  &                               - 0.203814e-03 ) * zt   &
423                  &                               + 0.170907e-01 ) * zt   &
424                  &   + 0.665157e-01                                      &
425                  &   +     ( - 0.678662e-05 * zs                         &
426                  &           - 0.846960e-04 * zt + 0.378110e-02 ) * zs   &
427                  &   +   ( ( - 0.302285e-13 * zh                         &
428                  &           - 0.251520e-11 * zs                         &
429                  &           + 0.512857e-12 * zt * zt           ) * zh   &
430                  &           - 0.164759e-06 * zs                         &
431                  &        +(   0.791325e-08 * zt - 0.933746e-06 ) * zt   &
432                  &                               + 0.380374e-04 ) * zh
433
434               zbeta  = ( ( -0.415613e-09 * zt + 0.555579e-07 ) * zt      &   ! beta
435                  &                            - 0.301985e-05 ) * zt      &
436                  &   + 0.785567e-03                                      &
437                  &   + (     0.515032e-08 * zs                           &
438                  &         + 0.788212e-08 * zt - 0.356603e-06 ) * zs     &
439                  &   +(  (   0.121551e-17 * zh                           &
440                  &         - 0.602281e-15 * zs                           &
441                  &         - 0.175379e-14 * zt + 0.176621e-12 ) * zh     &
442                  &                             + 0.408195e-10   * zs     &
443                  &     + ( - 0.213127e-11 * zt + 0.192867e-09 ) * zt     &
444                  &                             - 0.121555e-07 ) * zh
445
446               zthermal = zbeta * zalbet / ( rcp * zrhos + epsln )
447               zhalin   = zbeta * tsn(ji,jj,1,jp_sal) * rcs
448            ELSE
449               zrhos    = rhop(ji,jj,1) + rau0 * ( 1. - tmask(ji,jj,1) )
450               zthermal = rn_alpha / ( rcp * zrhos + epsln )
451               zhalin   = rn_beta * tsn(ji,jj,1,jp_sal) * rcs
452            ENDIF
453            ! Radiative surface buoyancy force
454            zBosol(ji,jj) = grav * zthermal * qsr(ji,jj)
455            ! Non radiative surface buoyancy force
456            zBo   (ji,jj) = grav * zthermal * qns(ji,jj) -  grav * zhalin * ( emps(ji,jj)-rnf(ji,jj) ) 
457            ! Surface Temperature flux for non-local term
458            wt0(ji,jj) = - ( qsr(ji,jj) + qns(ji,jj) )* ro0cpr * tmask(ji,jj,1)
459            ! Surface salinity flux for non-local term
460            ws0(ji,jj) = - ( ( emps(ji,jj)-rnf(ji,jj) ) * tsn(ji,jj,1,jp_sal) * rcs ) * tmask(ji,jj,1) 
461         ENDDO
462      ENDDO
463
464      zflageos = 0.5 + SIGN( 0.5, nn_eos - 1. ) 
465      !  Compute surface buoyancy forcing, Monin Obukhov and Ekman depths 
466      !------------------------------------------------------------------   
467      DO jj = 2, jpjm1
468         DO ji = fs_2, fs_jpim1
469            !  Reference surface density = density at first T point level   
470            zrhos         = rhop(ji,jj,1) + zflageos * rau0 * ( 1. - tmask(ji,jj,1) ) 
471            ! Friction velocity (zustar), at T-point : LMD94 eq. 2
472            zustar(ji,jj) = SQRT( taum(ji,jj) / ( zrhos +  epsln ) )
473         ENDDO
474      ENDDO
475
476!CDIR NOVERRCHK 
477      !                                               ! ===============
478      DO jj = 2, jpjm1                                 !  Vertical slab
479         !                                             ! ===============
480         
481         !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
482         ! II Compute Boundary layer mixing coef. and diagnose the new boundary layer depth
483         !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
484         
485         ! Initialization
486         jkmax       = 0
487         zdept (:,:) = 0.
488         zdepw (:,:) = 0.
489         zriblk(:,:) = 0.
490         zmoek (:,:) = 0.
491         zvisc (:,:) = 0.
492         zdift (:,:) = 0.
493#if defined key_zdfddm
494         zdifs (:,:) = 0.
495#endif
496         zmask (:,:) = 0.
497         DO ji = fs_2, fs_jpim1
498            zria(ji ) = 0.
499            ! Maximum boundary layer depth
500            ikbot     = mbkt(ji,jj)    ! ikbot is the last T point in the water
501            zhmax(ji) = fsdept(ji,jj,ikbot) - 0.001     
502            ! Compute Monin obukhov length scale at the surface and Ekman depth:
503            zbuofdep   = zBo(ji,jj) + zBosol(ji,jj) * ratt(1)
504            zekman(ji) = rcekman * zustar(ji,jj) / ( ABS( ff(ji,jj) ) + epsln )
505            zucube     = zustar(ji,jj) * zustar(ji,jj) * zustar(ji,jj) 
506            zmoa(ji)   = zucube / ( vonk * ( zbuofdep + epsln ) )   
507#if defined key_c1d
508            ! store the surface buoyancy forcing
509            zstabl        = 0.5 + SIGN( 0.5, zbuofdep )
510            buof(ji,jj,1) = zbuofdep * tmask(ji,jj,1)
511            ! store the moning-oboukov length scale at surface
512            zmob          = zstabl * zmoa(ji) + ( 1.0 - zstabl ) * fsdept(ji,jj,1)
513            mols(ji,jj,1) = MIN( zmob , zhmax(ji) ) * tmask(ji,jj,1)
514            ! store Ekman depth
515            zek           = zstabl * zekman(ji) + ( 1.0 - zstabl ) * fsdept(ji,jj,1) 
516            ekdp(ji,jj )  = MIN( zek , zhmax(ji) ) * tmask(ji,jj,1) 
517#endif
518         END DO     
519         ! Compute the pipe
520         ! ---------------------
521         DO jk = 2, jpkm1
522            DO ji = fs_2, fs_jpim1
523               ! Compute bfsfc = Bo + radiative contribution down to hbf*depht
524               zbuofdep = zBo(ji,jj) + zBosol(ji,jj) * ratt(jk)
525               ! Flag (zstabl  = 1) if positive forcing
526               zstabl   =  0.5 + SIGN(  0.5, zbuofdep)
527
528               !   Compute bulk richardson number zrib at depht
529               !-------------------------------------------------------
530               !                           [Br - B(d)] * d         zrinum
531               !             Rib(z) = ----------------------- = -------------
532               !                       |Vr - V(d)|^2 + Vt(d)^2   zdVsq + zVtsq
533               !
534               ! First compute zt,zs,zu,zv = means in the surface layer < epsilon*depht 
535               ! Else surface values are taken at the first T level.
536               ! For stability, resolved vertical shear is computed with "before velocities".
537               zref = epsilon * fsdept(ji,jj,jk)
538#if defined key_kppcustom
539               ! zref = gdept(1)
540               zref = fsdept(ji,jj,1)
541               zt   = tsn(ji,jj,1,jp_tem)
542               zs   = tsn(ji,jj,1,jp_sal)
543               zrh  = rhop(ji,jj,1)
544               zu   = ( ub(ji,jj,1) + ub(ji - 1,jj    ,1) ) / MAX( 1. , umask(ji,jj,1) + umask(ji - 1,jj   ,1) )
545               zv   = ( vb(ji,jj,1) + vb(ji    ,jj - 1,1) ) / MAX( 1. , vmask(ji,jj,1) + vmask(ji   ,jj - 1,1) )
546#else
547               zt   = 0.
548               zs   = 0.
549               zu   = 0.
550               zv   = 0.
551               zrh  = 0.
552               ! vertically integration over the upper epsilon*gdept(jk) ; del () array is computed once in zdf_kpp_init
553               DO jm = 1, jpkm1
554                  zt   = zt  + del(jk,jm) * tsn(ji,jj,jm,jp_tem)
555                  zs   = zs  + del(jk,jm) * tsn(ji,jj,jm,jp_sal)
556                  zu   = zu  + 0.5 * del(jk,jm) &
557                     &            * ( ub(ji,jj,jm) + ub(ji - 1,jj,jm) ) &
558                     &            / MAX( 1. , umask(ji,jj,jm) + umask(ji - 1,jj,jm) )
559                  zv   = zv  + 0.5 * del(jk,jm) &
560                     &            * ( vb(ji,jj,jm) + vb(ji,jj - 1,jm) ) &
561                     &            / MAX( 1. , vmask(ji,jj,jm) + vmask(ji,jj - 1,jm) )
562                  zrh  = zrh + del(jk,jm) * rhop(ji,jj,jm)
563               END DO
564#endif
565               zsr = SQRT( ABS( tsn(ji,jj,jk,jp_sal) ) )
566               ! depth
567               zh = fsdept(ji,jj,jk)
568               ! compute compression terms on density
569               ze  = ( -3.508914e-8*zt-1.248266e-8 ) *zt-2.595994e-6
570               zbw = (  1.296821e-6*zt-5.782165e-9 ) *zt+1.045941e-4
571               zb  = zbw + ze * zs
572               
573               zd  = -2.042967e-2
574               zc  =   (-7.267926e-5*zt+2.598241e-3 ) *zt+0.1571896
575               zaw = ( ( 5.939910e-6*zt+2.512549e-3 ) *zt-0.1028859 ) *zt - 4.721788
576               za  = ( zd*zsr + zc ) *zs + zaw
577               
578               zb1 =   (-0.1909078*zt+7.390729 ) *zt-55.87545
579               za1 = ( ( 2.326469e-3*zt+1.553190)*zt-65.00517 ) *zt+1044.077
580               zkw = ( ( (-1.361629e-4*zt-1.852732e-2 ) *zt-30.41638 ) *zt + 2098.925 ) *zt+190925.6
581               zk0 = ( zb1*zsr + za1 )*zs + zkw
582               zcomp =   1.0 - zh / ( zk0 - zh * ( za - zh * zb ) )
583               
584#if defined key_kppcustom
585               ! potential density of water(zrh = zt,zs at level jk):
586               zrhdr = zrh / zcomp
587#else
588               ! potential density of water(ztref,zsref at level jk):
589               ! compute volumic mass pure water at atm pressure
590               IF ( nn_eos < 1 ) THEN
591                  zr1= ( ( ( ( 6.536332e-9*zt-1.120083e-6 )*zt+1.001685e-4)*zt   &
592                     &       -9.095290e-3 )*zt+6.793952e-2 )*zt+999.842594
593                  ! seawater volumic mass atm pressure
594                  zr2= ( ( ( 5.3875e-9*zt-8.2467e-7 ) *zt+7.6438e-5 ) *zt   &
595                     &   -4.0899e-3 ) *zt+0.824493
596                  zr3= ( -1.6546e-6*zt+1.0227e-4 ) *zt-5.72466e-3
597                  zr4= 4.8314e-4             
598                  ! potential volumic mass (reference to the surface)
599                  zrhop= ( zr4*zs + zr3*zsr + zr2 ) *zs + zr1                 
600                  zrhdr = zrhop / zcomp
601               ELSE
602                  zrhdr = zrh / zcomp
603               ENDIF
604#endif
605               
606               ! potential density of ambiant water at level jk :
607               zrhd   = ( rhd(ji,jj,jk) * rau0 + rau0 ) 
608               
609               ! And now the Rib number numerator .
610               zrinum = grav * ( zrhd - zrhdr ) / rau0
611               zrinum = zrinum * ( fsdept(ji,jj,jk) - zref ) * tmask(ji,jj,jk)
612           
613               ! Resolved shear contribution to Rib at depth T-point (zdVsq)
614               ztx    =   ( ub( ji , jj ,jk)   +  ub(ji - 1, jj ,jk) ) &
615                  &     / MAX( 1. , umask( ji , jj ,jk) + umask(ji - 1, jj ,jk) )   
616               zty    =   ( vb( ji , jj ,jk)   +  vb(ji  ,jj - 1,jk) ) &
617                  &     / MAX( 1., vmask( ji , jj ,jk) + vmask(ji  ,jj - 1,jk) ) 
618               
619               zdVsq  = ( zu - ztx ) * ( zu - ztx ) + ( zv - zty ) * ( zv - zty )
620               
621               ! Scalar turbulent velocity scale zws for hbl=gdept
622               zscale = zstabl + ( 1.0 - zstabl ) * epsilon
623               zehat  = vonk * zscale * fsdept(ji,jj,jk) * zbuofdep
624               zucube = zustar(ji,jj) * zustar(ji,jj) * zustar(ji,jj)             
625               zeta   = zehat / ( zucube + epsln )
626               
627               IF( zehat > 0. ) THEN
628                  ! Stable case
629                  zws  = vonk * zustar(ji,jj) / ( 1.0 + rconc1 * zeta )
630               ELSE
631                  ! Unstable case
632#if defined key_kpplktb
633                  ! use lookup table
634                  zd     = zehat - dehatmin
635                  il     = INT( zd / dezehat )
636                  il     = MIN( il, nilktbm1 )
637                  il     = MAX( il, 1 )
638                 
639                  ud     = zustar(ji,jj) - ustmin
640                  jl     = INT( ud / deustar )
641                  jl     = MIN( jl, njlktbm1 )
642                  jl     = MAX( jl, 1 )
643                 
644                  zfrac  = zd / dezehat - FLOAT( il ) 
645                  ufrac  = ud / deustar - FLOAT( jl )
646                  zwas   = ( 1. - zfrac ) * wslktb(il,jl+1) + zfrac * wslktb(il+1,jl+1)
647                  zwbs   = ( 1. - zfrac ) * wslktb(il,jl  ) + zfrac * wslktb(il+1,jl  )
648                  !
649                  zws    = ( 1. - ufrac ) * zwbs + ufrac * zwas
650#else
651                  ! use analytical functions:
652                  zcons  = 0.5 + SIGN( 0.5 , ( rzetas - zeta ) )
653                  zwcons = vonk * zustar(ji,jj) * ( ( ABS( rconas - rconcs * zeta ) )**pthird ) 
654                  zwsun  = vonk * zustar(ji,jj) * SQRT( ABS ( 1.0 - rconc2 * zeta ) )
655                  !
656                  zws    = zcons * zwcons +  ( 1.0 - zcons) * zwsun
657#endif
658               ENDIF
659               
660               ! Turbulent shear contribution to Rib (zVtsq) bv frequency at levels  ( ie T-point jk)
661               zrn2   = 0.5 * ( rn2(ji,jj,jk) + rn2(ji,jj,jk+1) )   
662               zbvzed = SQRT( ABS( zrn2 ) ) 
663               zVtsq  = fsdept(ji,jj,jk) * zws * zbvzed  * Vtc
664               
665               ! Finally, the bulk Richardson number at depth fsdept(i,j,k)
666               zrib  = zrinum   / ( zdVsq + zVtsq + epsln )
667 
668               ! Find subscripts around the boundary layer depth, build the pipe
669               ! ----------------------------------------------------------------
670
671               ! Flag (zflagri = 1) if zrib < Ricr 
672               zflagri = 0.5 + SIGN( 0.5, ( Ricr - zrib ) )
673               !  Flag (zflagh  = 1) if still within overall boundary layer
674               zflagh  = 0.5 + SIGN( 0.5, ( fsdept(ji,jj,1) - zdept(ji,2) ) )
675               
676               ! Ekman layer depth
677               zek     = zstabl * zekman(ji) + ( 1.0 - zstabl ) * zhmax(ji)
678               zflag   = 0.5 + SIGN( 0.5, ( zek - fsdept(ji,jj,jk-1) ) )
679               zek     = zflag * zek + ( 1.0 - zflag ) * zhmax(ji)
680               zflagek = 0.5 + SIGN( 0.5, ( zek - fsdept(ji,jj,jk) ) )
681               ! Flag (zflagmo = 1) if still within stable Monin-Obukhov and in water
682               zmob    = zucube / ( vonk * ( zbuofdep + epsln ) ) 
683               ztemp   = zstabl * zmob + ( 1.0 - zstabl) * zhmax(ji)
684               ztemp   = MIN( ztemp , zhmax(ji) ) 
685               zflagmo = 0.5 + SIGN( 0.5, ( ztemp - fsdept(ji,jj,jk) ) )             
686
687               ! No limitation by Monin Obukhov or Ekman depths:
688!               zflagek = 1.0
689!               zflagmo = 0.5 + SIGN( 0.5, ( zhmax(ji) - fsdept(ji,jj,jk) ) )
690
691               ! Load  pipe via zflagkb  for later calculations
692               ! Flag (zflagkb = 1) if zflagh = 1 and (zflagri = 0 or zflagek = 0 or zflagmo = 0)
693               zflagkb = zflagh * ( 1.0 - ( zflagri * zflagek * zflagmo ) )
694
695               zmask(ji,jk) = zflagh
696               jkp2         = MIN( jk+2 , ikbot )
697               jkm1         = MAX( jk-1 , 2 )
698               jkmax        = MAX( jkmax, jk * INT( REAL( zflagh+epsln ) ) ) 
699
700               zdept(ji,1)  = zdept(ji,1) + zflagkb * fsdept(ji,jj,jk-1) 
701               zdept(ji,2)  = zdept(ji,2) + zflagkb * fsdept(ji,jj,jk  ) 
702               zdept(ji,3)  = zdept(ji,3) + zflagkb * fsdept(ji,jj,jk+1) 
703
704               zdepw(ji,1)  = zdepw(ji,1) + zflagkb * fsdepw(ji,jj,jk-1) 
705               zdepw(ji,2)  = zdepw(ji,2) + zflagkb * fsdepw(ji,jj,jk  ) 
706               zdepw(ji,3)  = zdepw(ji,3) + zflagkb * fsdepw(ji,jj,jk+1)
707               zdepw(ji,4)  = zdepw(ji,4) + zflagkb * fsdepw(ji,jj,jkp2) 
708
709               zriblk(ji,1) = zriblk(ji,1) + zflagkb * zria(ji)
710               zriblk(ji,2) = zriblk(ji,2) + zflagkb * zrib
711
712               zmoek (ji,0) = zmoek (ji,0) + zflagkb * zek
713               zmoek (ji,1) = zmoek (ji,1) + zflagkb * zmoa(ji)
714               zmoek (ji,2) = zmoek (ji,2) + zflagkb * ztemp 
715               ! Save Monin Obukhov depth
716               zmoa  (ji)   = zmob
717           
718               zvisc(ji,1) = zvisc(ji,1) + zflagkb * avmu(ji,jj,jkm1)
719               zvisc(ji,2) = zvisc(ji,2) + zflagkb * avmu(ji,jj,jk  )
720               zvisc(ji,3) = zvisc(ji,3) + zflagkb * avmu(ji,jj,jk+1)
721               zvisc(ji,4) = zvisc(ji,4) + zflagkb * avmu(ji,jj,jkp2)
722               
723               zdift(ji,1) = zdift(ji,1) + zflagkb * avt (ji,jj,jkm1)
724               zdift(ji,2) = zdift(ji,2) + zflagkb * avt (ji,jj,jk  )
725               zdift(ji,3) = zdift(ji,3) + zflagkb * avt (ji,jj,jk+1)
726               zdift(ji,4) = zdift(ji,4) + zflagkb * avt (ji,jj,jkp2)
727
728#if defined key_zdfddm
729               zdifs(ji,1) = zdifs(ji,1) + zflagkb * avs (ji,jj,jkm1)
730               zdifs(ji,2) = zdifs(ji,2) + zflagkb * avs (ji,jj,jk  )
731               zdifs(ji,3) = zdifs(ji,3) + zflagkb * avs (ji,jj,jk+1)
732               zdifs(ji,4) = zdifs(ji,4) + zflagkb * avs (ji,jj,jkp2)
733#endif               
734               ! Save the Richardson number
735               zria  (ji)   = zrib 
736#if defined key_c1d
737               ! store buoyancy length scale
738               buof(ji,jj,jk) = zbuofdep * tmask(ji,jj,jk) 
739               ! store Monin Obukhov
740               zmob           = zstabl * zmob + ( 1.0 - zstabl) * fsdept(ji,jj,1)
741               mols(ji,jj,jk) = MIN( zmob , zhmax(ji) ) * tmask(ji,jj,jk) 
742               ! Bulk Richardson number
743               rib(ji,jj,jk)  = zrib * tmask(ji,jj,jk)             
744#endif               
745            END DO
746         END DO
747         !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
748         ! III PROCESS THE PIPE
749         !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
750         
751         DO ji = fs_2, fs_jpim1 
752           
753            ! Find the boundary layer depth zhbl
754            ! ----------------------------------------
755           
756            ! Interpolate monin Obukhov and critical Ri mumber depths   
757            ztemp = zdept(ji,2) - zdept(ji,1)
758            zflag = ( Ricr - zriblk(ji,1) ) / ( zriblk(ji,2) - zriblk(ji,1)  + epsln )
759            zhrib = zdept(ji,1) + zflag * ztemp     
760
761            IF( zriblk(ji,2) < Ricr ) zhrib = zhmax(ji) 
762         
763            IF( zmoek(ji,2) < zdept(ji,2) ) THEN
764               IF ( zmoek(ji,1) < 0. ) THEN
765                  zmob = zdept(ji,2) - epsln
766               ELSE
767                  zmob = ztemp + zmoek(ji,1) - zmoek(ji,2)
768                  zmob = ( zmoek(ji,1) * zdept(ji,2) - zmoek(ji,2) * zdept(ji,1) ) / zmob
769                  zmob = MAX( zmob , zdept(ji,1) + epsln )               
770               ENDIF
771            ELSE           
772               zmob = zhmax(ji) 
773            ENDIF
774            ztemp   = MIN( zmob , zmoek(ji,0) )
775                         
776            ! Finally, the boundary layer depth, zhbl
777            zhbl(ji) = MAX( fsdept(ji,jj,1) + epsln, MIN( zhrib , ztemp ) )
778           
779            ! Save hkpp for further diagnostics (optional)
780            hkpp(ji,jj) = zhbl(ji) * tmask(ji,jj,1) 
781         
782            ! Correct mask if zhbl < fsdepw(ji,jj,2) for no viscosity/diffusivity enhancement at fsdepw(ji,jj,2)
783            !     zflag = 1 if zhbl(ji) > fsdepw(ji,jj,2)
784            IF( zhbl(ji) < fsdepw(ji,jj,2) ) zmask(ji,2) = 0.
785         
786           
787            !  Velocity scales at depth zhbl
788            ! -----------------------------------
789           
790            !  Compute bouyancy forcing down to zhbl
791            ztemp    = -hbf * zhbl(ji)
792            zatt1    = 1.0 - ( rabs * EXP( ztemp / xsi1 ) + ( 1.0 - rabs ) * EXP( ztemp / xsi2 ) )
793            zbuofdep = zBo(ji,jj) + zBosol(ji,jj) * zatt1
794            zstabl   = 0.5 + SIGN( 0.5 , zbuofdep ) 
795
796            zbuofdep = zbuofdep + zstabl * epsln
797
798            zscale = zstabl + ( 1.0 - zstabl ) * epsilon         
799            zehat  = vonk * zscale * zhbl(ji) * zbuofdep
800            zucube = zustar(ji,jj) * zustar(ji,jj) * zustar(ji,jj)             
801            zeta   = zehat / ( zucube + epsln )
802           
803            IF( zehat > 0. ) THEN
804               ! Stable case
805               zws  = vonk * zustar(ji,jj) / ( 1.0 + rconc1 * zeta )
806               zwm  = zws
807            ELSE
808               ! Unstable case
809#if defined key_kpplktb
810               ! use lookup table
811               zd     = zehat - dehatmin
812               il     = INT( zd / dezehat )
813               il     = MIN( il, nilktbm1 )
814               il     = MAX( il, 1 )
815               
816               ud     = zustar(ji,jj) - ustmin
817               jl     = INT( ud / deustar )
818               jl     = MIN( jl, njlktbm1 )
819               jl     = MAX( jl, 1 )
820               
821               zfrac  = zd / dezehat - FLOAT( il ) 
822               ufrac  = ud / deustar - FLOAT( jl )
823               zwas   = ( 1. - zfrac ) * wslktb(il,jl+1) + zfrac * wslktb(il+1,jl+1)
824               zwbs   = ( 1. - zfrac ) * wslktb(il,jl  ) + zfrac * wslktb(il+1,jl  )
825               zwam   = ( 1. - zfrac ) * wmlktb(il,jl+1) + zfrac * wmlktb(il+1,jl+1)
826               zwbm   = ( 1. - zfrac ) * wmlktb(il,jl  ) + zfrac * wmlktb(il+1,jl  )
827               !
828               zws    = ( 1. - ufrac ) * zwbs + ufrac * zwas
829               zwm    = ( 1. - ufrac ) * zwbm + ufrac * zwam
830#else
831               ! use analytical functions
832               zconm  = 0.5 + SIGN( 0.5, ( rzetam - zeta) )
833               zcons  = 0.5 + SIGN( 0.5, ( rzetas - zeta) )
834               
835               ! Momentum : zeta < rzetam (zconm = 1)
836               ! Scalars  : zeta < rzetas (zcons = 1)
837               zwconm = zustar(ji,jj) * vonk * ( ( ABS( rconam - rconcm * zeta) )**pthird )
838               zwcons = zustar(ji,jj) * vonk * ( ( ABS( rconas - rconcs * zeta) )**pthird )
839               
840               ! Momentum : rzetam <= zeta < 0 (zconm = 0)
841               ! Scalars  : rzetas <= zeta < 0 (zcons = 0) 
842               zwmun  = SQRT( ABS( 1.0 - rconc2 * zeta ) )
843               zwsun  = vonk * zustar(ji,jj) * zwmun
844               zwmun  = vonk * zustar(ji,jj) * SQRT(zwmun)
845               !
846               zwm    = zconm * zwconm + ( 1.0 - zconm ) * zwmun
847               zws    = zcons * zwcons + ( 1.0 - zcons ) * zwsun
848               
849#endif
850            ENDIF
851           
852           
853            ! Viscosity, diffusivity values and derivatives at h
854            ! --------------------------------------------------------
855           
856            ! check between at which interfaces is located zhbl(ji)
857            ! ztemp = 1, zdepw(ji,2) < zhbl <  zdepw(ji,3)
858            ! ztemp = 0, zdepw(ji,1) < zhbl <  zdepw(ji,2)
859            ztemp  =  0.5 + SIGN( 0.5, ( zhbl(ji) - zdepw(ji,2) ) ) 
860            zdep21 =   zdepw(ji,2) - zdepw(ji,1) + epsln
861            zdep32 =   zdepw(ji,3) - zdepw(ji,2) + epsln
862            zdep43 =   zdepw(ji,4) - zdepw(ji,3) + epsln 
863           
864            ! Compute R as in LMD94, eq D5b
865            zdelta =  ( zhbl(ji) - zdepw(ji,2) ) *         ztemp   / zdep32   &
866               &    + ( zhbl(ji) - zdepw(ji,1) ) * ( 1.0 - ztemp ) / zdep21 
867           
868            ! Compute the vertical derivative of viscosities (zdzh) at z=zhbl(ji)
869            zdzup  =  ( zvisc(ji,2) - zvisc(ji,3) ) *         ztemp   / zdep32 &
870               &    + ( zvisc(ji,1) - zvisc(ji,2) ) * ( 1.0 - ztemp ) / zdep21
871           
872            zdzdn  =  ( zvisc(ji,3) - zvisc(ji,4) ) *         ztemp   / zdep43 &
873               &    + ( zvisc(ji,2) - zvisc(ji,3) ) * ( 1.0 - ztemp ) / zdep32
874           
875            ! LMD94, eq D5b :         
876            zdzh   = ( 1.0 - zdelta ) * zdzup + zdelta * zdzdn
877            zdzh   = MAX( zdzh , 0. )
878           
879            ! Compute viscosities (zvath) at z=zhbl(ji), LMD94 eq D5a
880            zvath  =          ztemp   * ( zvisc(ji,3) + zdzh * ( zdepw(ji,3) - zhbl(ji) ) ) &
881               &    + ( 1.0 - ztemp ) * ( zvisc(ji,2) + zdzh * ( zdepw(ji,2) - zhbl(ji) ) )
882           
883            ! Compute G (zgat1) and its derivative (zdat1) at z=hbl(ji), LMD94 eq 18
884           
885            ! Vertical derivative of velocity scale divided by velocity scale squared at z=hbl(ji)
886            ! (non zero only in stable conditions)
887            zflag  =  -zstabl * rconc1 * zbuofdep / ( zucube * zustar(ji,jj) + epsln )
888           
889            ! G at its derivative at z=hbl:
890            zgat1  = zvath  / ( zhbl(ji) * ( zwm + epsln )  )
891            zdat1  = -zdzh  / ( zwm + epsln ) -  zflag * zvath / zhbl(ji)
892           
893            ! G coefficients, LMD94 eq 17
894            za2m(ji) = -2.0 + 3.0 * zgat1 - zdat1
895            za3m(ji) =  1.0 - 2.0 * zgat1 + zdat1
896
897           
898            ! Compute the vertical derivative of temperature diffusivities (zdzh) at z=zhbl(ji)
899            zdzup  =  ( zdift(ji,2) - zdift(ji,3) ) *         ztemp   / zdep32 &
900               &    + ( zdift(ji,1) - zdift(ji,2) ) * ( 1.0 - ztemp ) / zdep21
901           
902            zdzdn  =  ( zdift(ji,3) - zdift(ji,4) ) *         ztemp   / zdep43 &
903               &    + ( zdift(ji,2) - zdift(ji,3) ) * ( 1.0 - ztemp ) / zdep32
904           
905            ! LMD94, eq D5b :         
906            zdzh   = ( 1.0 - zdelta ) * zdzup + zdelta * zdzdn
907            zdzh   = MAX( zdzh , 0. )
908           
909           
910            ! Compute diffusivities (zvath) at z=zhbl(ji), LMD94 eq D5a
911            zvath  =          ztemp   * ( zdift(ji,3) + zdzh * ( zdepw(ji,3) - zhbl(ji) ) ) &
912               &    + ( 1.0 - ztemp ) * ( zdift(ji,2) + zdzh * ( zdepw(ji,2) - zhbl(ji) ) )
913                       
914            ! G at its derivative at z=hbl:
915            zgat1  = zvath  / ( zhbl(ji) * ( zws + epsln )  )
916            zdat1  = -zdzh  / ( zws + epsln ) -  zflag * zvath / zhbl(ji)
917           
918            ! G coefficients, LMD94 eq 17
919            za2t(ji) = -2.0 + 3.0 * zgat1 - zdat1
920            za3t(ji) =  1.0 - 2.0 * zgat1 + zdat1
921
922#if defined key_zdfddm
923            ! Compute the vertical derivative of salinities diffusivities (zdzh) at z=zhbl(ji)
924            zdzup  =  ( zdifs(ji,2) - zdifs(ji,3) ) *         ztemp   / zdep32 &
925               &    + ( zdifs(ji,1) - zdifs(ji,2) ) * ( 1.0 - ztemp ) / zdep21
926           
927            zdzdn  =  ( zdifs(ji,3) - zdifs(ji,4) ) *         ztemp   / zdep43 &
928               &    + ( zdifs(ji,2) - zdifs(ji,3) ) * ( 1.0 - ztemp ) / zdep32
929           
930            ! LMD94, eq D5b :         
931            zdzh   = ( 1.0 - zdelta ) * zdzup + zdelta * zdzdn
932            zdzh   = MAX( zdzh , 0. )           
933           
934            ! Compute diffusivities (zvath) at z=zhbl(ji), LMD94 eq D5a
935            zvath  =          ztemp   * ( zdifs(ji,3) + zdzh * ( zdepw(ji,3) - zhbl(ji) ) ) &
936               &    + ( 1.0 - ztemp ) * ( zdifs(ji,2) + zdzh * ( zdepw(ji,2) - zhbl(ji) ) )
937                       
938            ! G at its derivative at z=hbl:
939            zgat1  = zvath  / ( zhbl(ji) * ( zws + epsln )  )
940            zdat1  = -zdzh  / ( zws + epsln ) -  zflag * zvath / zhbl(ji)
941           
942            ! G coefficients, LMD94 eq 17
943            za2s(ji) = -2.0 + 3.0 * zgat1 - zdat1
944            za3s(ji) =  1.0 - 2.0 * zgat1 + zdat1
945#endif
946
947            !-------------------turn off interior matching here------
948            !          za2(ji,1) = -2.0
949            !          za3(ji,1) =  1.0
950            !          za2(ji,2) = -2.0
951            !          za3(ji,2) =  1.0
952            !--------------------------------------------------------
953           
954            !  Compute Enhanced Mixing Coefficients (LMD94,eq D6)
955            ! ---------------------------------------------------------------
956           
957            ! Delta
958            zdelta  = ( zhbl(ji)  - zdept(ji,1) ) / ( zdept(ji,2) - zdept(ji,1) + epsln )
959            zdelta2 = zdelta * zdelta
960           
961            !  Mixing coefficients at first level above h (zdept(ji,1))
962            ! and at first interface in the pipe (zdepw(ji,2))
963           
964            ! At first T level above h (zdept(ji,1)) (always in the boundary layer)
965            zsig    = zdept(ji,1) / zhbl(ji)
966            ztemp   = zstabl * zsig  + ( 1.0 - zstabl ) * MIN( zsig , epsilon )
967            zehat   = vonk * ztemp * zhbl(ji) * zbuofdep
968            zeta    = zehat / ( zucube + epsln)
969            zwst    = vonk * zustar(ji,jj) / ( ABS( 1.0 + rconc1 * zeta ) + epsln)
970            zwm     = zstabl * zwst + ( 1.0 - zstabl ) * zwm
971            zws     = zstabl * zwst + ( 1.0 - zstabl ) * zws
972
973            zkm1m  = zhbl(ji) * zwm * zsig * ( 1.0 + zsig * ( za2m(ji) + zsig * za3m(ji) ) )
974            zkm1t  = zhbl(ji) * zws * zsig * ( 1.0 + zsig * ( za2t(ji) + zsig * za3t(ji) ) )
975#if defined key_zdfddm
976            zkm1s  = zhbl(ji) * zws * zsig * ( 1.0 + zsig * ( za2s(ji) + zsig * za3s(ji) ) )
977#endif                       
978            ! At first W level in the pipe (zdepw(ji,2)) (not always in the boundary layer ):
979            zsig    = MIN( zdepw(ji,2) / zhbl(ji) , 1.0 )
980            ztemp   = zstabl * zsig + ( 1.0 - zstabl ) * MIN( zsig , epsilon )
981            zehat   = vonk * ztemp * zhbl(ji) * zbuofdep
982            zeta    = zehat / ( zucube + epsln )
983            zwst    = vonk * zustar(ji,jj) / ( ABS( 1.0 + rconc1 * zeta ) + epsln)
984            zws     = zstabl * zws + ( 1.0 - zstabl ) * zws
985            zwm     = zstabl * zws + ( 1.0 - zstabl ) * zwm
986
987            zkmpm(ji) = zhbl(ji) * zwm * zsig * ( 1.0 + zsig * ( za2m(ji) + zsig * za3m(ji) ) )
988            zkmpt(ji) = zhbl(ji) * zws * zsig * ( 1.0 + zsig * ( za2t(ji) + zsig * za3t(ji) ) )
989#if defined key_zdfddm
990            zkmps(ji) = zhbl(ji) * zws * zsig * ( 1.0 + zsig * ( za2s(ji) + zsig * za3s(ji) ) )
991#endif 
992     
993            ! check if this point is in the boundary layer,else take interior viscosity/diffusivity:
994            zflag       = 0.5 + SIGN( 0.5, ( zhbl(ji) - zdepw(ji,2) ) )
995            zkmpm(ji) = zkmpm(ji) * zflag + ( 1.0 - zflag ) * zvisc(ji,2)
996            zkmpt(ji) = zkmpt(ji) * zflag + ( 1.0 - zflag ) * zdift(ji,2)
997#if defined key_zdfddm
998            zkmps(ji) = zkmps(ji) * zflag + ( 1.0 - zflag ) * zdifs(ji,2)
999#endif
1000
1001            ! Enhanced viscosity/diffusivity at zdepw(ji,2)
1002            ztemp     = ( 1.0 - 2.0 * zdelta + zdelta2 ) * zkm1m + zdelta2 * zkmpm(ji)
1003            zkmpm(ji) = ( 1.0 - zdelta ) * zvisc(ji,2) + zdelta * ztemp
1004            ztemp     = ( 1.0 - 2.0 * zdelta + zdelta2 ) * zkm1t + zdelta2 * zkmpt(ji)
1005            zkmpt(ji) = ( 1.0 - zdelta ) * zdift(ji,2) + zdelta * ztemp
1006#if defined key_zdfddm
1007            ztemp     = ( 1.0 - 2.0 * zdelta + zdelta2 ) * zkm1s + zdelta2 * zkmps(ji)
1008            zkmps(ji) = ( 1.0 - zdelta ) * zdifs(ji,2) + zdelta * ztemp
1009#endif           
1010
1011         END DO
1012         !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
1013         ! IV. Compute vertical eddy viscosity and diffusivity coefficients
1014         !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1015         
1016         DO jk  = 2, jkmax
1017           
1018            ! Compute turbulent velocity scales on the interfaces
1019            ! --------------------------------------------------------
1020            DO  ji = fs_2, fs_jpim1
1021               zbuofdep = zBo(ji,jj) + zBosol(ji,jj) * zatt1
1022               zstabl   = 0.5 + SIGN( 0.5 , zbuofdep ) 
1023               zbuofdep = zbuofdep + zstabl * epsln         
1024               zsig    = fsdepw(ji,jj,jk) / zhbl(ji)
1025               ztemp   = zstabl * zsig + ( 1. - zstabl ) * MIN( zsig , epsilon )
1026               zehat   = vonk * ztemp * zhbl(ji) * zbuofdep
1027               zucube  = zustar(ji,jj) * zustar(ji,jj) * zustar(ji,jj)
1028               zeta    = zehat / ( zucube + epsln )
1029
1030               IF( zehat > 0. ) THEN
1031                  ! Stable case
1032                  zws  = vonk * zustar(ji,jj) / ( 1.0 + rconc1 * zeta )
1033                  zwm  = zws
1034               ELSE
1035                  ! Unstable case
1036#if defined key_kpplktb
1037                  ! use lookup table
1038                  zd     = zehat - dehatmin
1039                  il     = INT( zd / dezehat )
1040                  il     = MIN( il, nilktbm1 )
1041                  il     = MAX( il, 1 )
1042                 
1043                  ud     = zustar(ji,jj) - ustmin
1044                  jl     = INT( ud / deustar )
1045                  jl     = MIN( jl, njlktbm1 )
1046                  jl     = MAX( jl, 1 )
1047                 
1048                  zfrac  = zd / dezehat - FLOAT( il ) 
1049                  ufrac  = ud / deustar - FLOAT( jl )
1050                  zwas   = ( 1. - zfrac ) * wslktb(il,jl+1) + zfrac * wslktb(il+1,jl+1)
1051                  zwbs   = ( 1. - zfrac ) * wslktb(il,jl  ) + zfrac * wslktb(il+1,jl  )
1052                  zwam   = ( 1. - zfrac ) * wmlktb(il,jl+1) + zfrac * wmlktb(il+1,jl+1)
1053                  zwbm   = ( 1. - zfrac ) * wmlktb(il,jl  ) + zfrac * wmlktb(il+1,jl  )
1054                  !
1055                  zws    = ( 1. - ufrac ) * zwbs + ufrac * zwas
1056                  zwm    = ( 1. - ufrac ) * zwbm + ufrac * zwam
1057#else
1058                  ! use analytical functions
1059                  zconm  = 0.5 + SIGN( 0.5, ( rzetam - zeta) )
1060                  zcons  = 0.5 + SIGN( 0.5, ( rzetas - zeta) )
1061                 
1062                  ! Momentum : zeta < rzetam (zconm = 1)
1063                  ! Scalars  : zeta < rzetas (zcons = 1)
1064                  zwconm = zustar(ji,jj) * vonk * ( ( ABS( rconam - rconcm * zeta) )**pthird )
1065                  zwcons = zustar(ji,jj) * vonk * ( ( ABS( rconas - rconcs * zeta) )**pthird )
1066                 
1067                  ! Momentum : rzetam <= zeta < 0 (zconm = 0)
1068                  ! Scalars  : rzetas <= zeta < 0 (zcons = 0) 
1069                  zwmun  = SQRT( ABS( 1.0 - rconc2 * zeta ) )
1070                  zwsun  = vonk * zustar(ji,jj) * zwmun
1071                  zwmun  = vonk * zustar(ji,jj) * SQRT(zwmun)
1072                  !
1073                  zwm    = zconm * zwconm + ( 1.0 - zconm ) * zwmun
1074                  zws    = zcons * zwcons + ( 1.0 - zcons ) * zwsun
1075                 
1076#endif
1077               ENDIF
1078               
1079               zblcm(ji,jk) = zhbl(ji) * zwm * zsig  * ( 1.0 + zsig * ( za2m(ji) + zsig * za3m(ji) ) )
1080               zblct(ji,jk) = zhbl(ji) * zws * zsig  * ( 1.0 + zsig * ( za2t(ji) + zsig * za3t(ji) ) )
1081#if defined key_zdfddm
1082               zblcs(ji,jk) = zhbl(ji) * zws * zsig  * ( 1.0 + zsig * ( za2s(ji) + zsig * za3s(ji) ) )
1083#endif             
1084               !  Compute Nonlocal transport term = ghats * <ws>o
1085               ! ----------------------------------------------------
1086               ghats(ji,jj,jk-1) = ( 1. - zstabl ) * rcg / ( zws * zhbl(ji) + epsln ) * tmask(ji,jj,jk)
1087
1088            END DO
1089         END DO     
1090         !  Combine interior and boundary layer coefficients and nonlocal term
1091         ! -----------------------------------------------------------------------
1092         DO jk = 2, jpkm1   
1093            DO ji = fs_2, fs_jpim1
1094               zflag = zmask(ji,jk) * zmask(ji,jk+1)
1095               zviscos(ji,jj,jk) = ( 1.0 - zmask(ji,jk) )         * avmu (ji,jj,jk) & ! interior viscosities
1096                  &              +                        zflag   * zblcm(ji,jk    ) & ! boundary layer viscosities
1097                  &              + zmask(ji,jk) * ( 1.0 - zflag ) * zkmpm(ji       )   ! viscosity enhancement at W_level near zhbl
1098               
1099               zviscos(ji,jj,jk) = zviscos(ji,jj,jk) * tmask(ji,jj,jk)   
1100
1101           
1102               zdiffut(ji,jj,jk) = ( 1.0 - zmask(ji,jk) )          * avt (ji,jj,jk) & ! interior diffusivities
1103                  &              +                        zflag   * zblct(ji,jk   ) & ! boundary layer diffusivities
1104                  &              + zmask(ji,jk) * ( 1.0 - zflag ) * zkmpt(ji      )   ! diffusivity enhancement at W_level near zhbl
1105                       
1106               zdiffut(ji,jj,jk) = zdiffut(ji,jj,jk) * tmask(ji,jj,jk) 
1107#if defined key_zdfddm
1108               zdiffus(ji,jj,jk) = ( 1.0 - zmask(ji,jk) )          * avs (ji,jj,jk) & ! interior diffusivities
1109                  &              +                        zflag   * zblcs(ji,jk   ) & ! boundary layer diffusivities
1110                  &              + zmask(ji,jk) * ( 1.0 - zflag ) * zkmps(ji      )   ! diffusivity enhancement at W_level near zhbl
1111                       
1112               zdiffus(ji,jj,jk) = zdiffus(ji,jj,jk) * tmask(ji,jj,jk) 
1113#endif               
1114               ! Non local flux in the boundary layer only
1115               ghats(ji,jj,jk-1) = zmask(ji,jk) * ghats(ji,jj,jk-1)
1116
1117            ENDDO
1118         END DO
1119         !                                             ! ===============
1120      END DO                                           !   End of slab
1121      !                                                ! ===============
1122
1123      ! Lateral boundary conditions on zvicos and zdiffus  (sign unchanged)
1124      CALL lbc_lnk( zviscos(:,:,:), 'U', 1. )  ; CALL lbc_lnk( zdiffut(:,:,:), 'W', 1. ) 
1125#if defined key_zdfddm 
1126      CALL lbc_lnk( zdiffus(:,:,:), 'W', 1. ) 
1127#endif
1128
1129      SELECT CASE ( nn_ave )
1130         !
1131         CASE ( 0 )             ! no viscosity and diffusivity smoothing
1132
1133            DO jk = 2, jpkm1
1134               DO jj = 2, jpjm1
1135                  DO ji = fs_2, fs_jpim1
1136                     avmu(ji,jj,jk) = ( zviscos(ji,jj,jk) + zviscos(ji+1,jj,jk) ) &
1137                        &  / MAX( 1., tmask(ji,jj,jk) + tmask (ji + 1,jj,jk) ) * umask(ji,jj,jk)
1138                     
1139                     avmv(ji,jj,jk) = ( zviscos(ji,jj,jk) + zviscos(ji,jj+1,jk) ) &
1140                        &  / MAX( 1., tmask(ji,jj,jk) + tmask (ji,jj+1,jk) ) * vmask(ji,jj,jk)
1141                     
1142                     avt (ji,jj,jk) =  zdiffut(ji,jj,jk) * tmask(ji,jj,jk) 
1143#if defined key_zdfddm     
1144                     avs (ji,jj,jk) =  zdiffus(ji,jj,jk) * tmask(ji,jj,jk) 
1145#endif
1146                  END DO
1147               END DO
1148            END DO
1149           
1150         CASE ( 1 )                ! viscosity and diffusivity smoothing
1151            !                     
1152            !           ( 1/2  1  1/2 )              ( 1/2  1/2 )             ( 1/2  1  1/2 )
1153            ! avt = 1/8 ( 1    2  1   )   avmu = 1/4 ( 1    1   )   avmv= 1/4 ( 1/2  1  1/2 )
1154            !           ( 1/2  1  1/2 )              ( 1/2  1/2 )
1155 
1156            DO jk = 2, jpkm1
1157               DO jj = 2, jpjm1
1158                  DO ji = fs_2, fs_jpim1
1159
1160                     avmu(ji,jj,jk) = (      zviscos(ji  ,jj  ,jk) + zviscos(ji+1,jj  ,jk)   &
1161                        &              +.5*( zviscos(ji  ,jj-1,jk) + zviscos(ji+1,jj-1,jk)   &
1162                        &                   +zviscos(ji  ,jj+1,jk) + zviscos(ji+1,jj+1,jk) ) ) * eumean(ji,jj,jk)
1163                     
1164                     avmv(ji,jj,jk) = (      zviscos(ji  ,jj  ,jk) + zviscos(ji  ,jj+1,jk)   &
1165                        &              +.5*( zviscos(ji-1,jj  ,jk) + zviscos(ji-1,jj+1,jk)   &
1166                        &                   +zviscos(ji+1,jj  ,jk) + zviscos(ji+1,jj+1,jk) ) ) * evmean(ji,jj,jk)
1167 
1168                     avt (ji,jj,jk) = ( .5*( zdiffut(ji-1,jj+1,jk) + zdiffut(ji-1,jj-1,jk)    &
1169                        &                   +zdiffut(ji+1,jj+1,jk) + zdiffut(ji+1,jj-1,jk) )  &
1170                        &              +1.*( zdiffut(ji-1,jj  ,jk) + zdiffut(ji  ,jj+1,jk)    &
1171                        &                   +zdiffut(ji  ,jj-1,jk) + zdiffut(ji+1,jj  ,jk) )  &
1172                        &              +2.*  zdiffut(ji  ,jj  ,jk)                          ) * etmean(ji,jj,jk)
1173#if defined key_zdfddm   
1174                     avs (ji,jj,jk) = ( .5*( zdiffus(ji-1,jj+1,jk) + zdiffus(ji-1,jj-1,jk)    &
1175                        &                   +zdiffus(ji+1,jj+1,jk) + zdiffus(ji+1,jj-1,jk) )  &
1176                        &              +1.*( zdiffus(ji-1,jj  ,jk) + zdiffus(ji  ,jj+1,jk)    &
1177                        &                   +zdiffus(ji  ,jj-1,jk) + zdiffus(ji+1,jj  ,jk) )  &
1178                        &              +2.*  zdiffus(ji  ,jj  ,jk)                          ) * etmean(ji,jj,jk) 
1179#endif               
1180                  END DO
1181               END DO
1182            END DO
1183         
1184         END SELECT
1185
1186         DO jk = 2, jpkm1                       ! vertical slab
1187            !
1188            !  Minimum value on the eddy diffusivity
1189            ! ----------------------------------------
1190            DO jj = 2, jpjm1
1191               DO ji = fs_2, fs_jpim1   ! vector opt.
1192                  avt(ji,jj,jk) = MAX( avt(ji,jj,jk), avtb(jk) ) * tmask(ji,jj,jk)
1193#if defined key_zdfddm 
1194                  avs(ji,jj,jk) = MAX( avs(ji,jj,jk), avtb(jk) ) * tmask(ji,jj,jk)
1195#endif
1196               END DO
1197            END DO
1198
1199            !
1200            ! Minimum value on the eddy viscosity
1201            ! ----------------------------------------
1202            DO jj = 1, jpj
1203               DO ji = 1, jpi
1204                  avmu(ji,jj,jk) = MAX( avmu(ji,jj,jk), avmb(jk) ) * umask(ji,jj,jk)
1205                  avmv(ji,jj,jk) = MAX( avmv(ji,jj,jk), avmb(jk) ) * vmask(ji,jj,jk)
1206               END DO
1207            END DO
1208            !
1209         END DO
1210
1211         ! Lateral boundary conditions on avt  (sign unchanged)
1212         CALL lbc_lnk( hkpp(:,:), 'T', 1. )
1213
1214         ! Lateral boundary conditions on avt  (sign unchanged)
1215         CALL lbc_lnk( avt(:,:,:), 'W', 1. )
1216#if defined key_zdfddm 
1217         CALL lbc_lnk( avs(:,:,:), 'W', 1. ) 
1218#endif
1219         ! Lateral boundary conditions (avmu,avmv) (U- and V- points, sign unchanged)
1220         CALL lbc_lnk( avmu(:,:,:), 'U', 1. )   ;    CALL lbc_lnk( avmv(:,:,:), 'V', 1. ) 
1221 
1222         IF(ln_ctl) THEN
1223#if defined key_zdfddm
1224            CALL prt_ctl(tab3d_1=avt , clinfo1=' kpp - t: ', tab3d_2=avs , clinfo2=' s: ', ovlap=1, kdim=jpk)
1225#else
1226            CALL prt_ctl(tab3d_1=avt , clinfo1=' kpp - t: ', ovlap=1, kdim=jpk)
1227#endif
1228            CALL prt_ctl(tab3d_1=avmu, clinfo1=' kpp - u: ', mask1=umask,  &
1229               &         tab3d_2=avmv, clinfo2=      ' v: ', mask2=vmask, ovlap=1, kdim=jpk)
1230         ENDIF
1231
1232      IF( wrk_not_released(1, 1,2,3,4,5,6,7,8,9,10,11,12,13,14) .OR. &
1233          wrk_not_released(2, 1,2,3,4,5,6,7,8,9,10,11)          .OR. &
1234          wrk_not_released(3, 1)                                .OR. &
1235          wrk_not_released_xz(1,2,3)  )   CALL ctl_stop('zdf_kpp : failed to release workspace arrays')
1236      !
1237   END SUBROUTINE zdf_kpp
1238
1239
1240   SUBROUTINE tra_kpp( kt )
1241      !!----------------------------------------------------------------------
1242      !!                  ***  ROUTINE tra_kpp  ***
1243      !!
1244      !! ** Purpose :   compute and add to the tracer trend the non-local tracer flux
1245      !!
1246      !! ** Method  :   ???
1247      !!----------------------------------------------------------------------
1248      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   ztrdt, ztrds   ! 3D workspace
1249      !!----------------------------------------------------------------------
1250      INTEGER, INTENT(in) :: kt
1251      INTEGER :: ji, jj, jk
1252
1253      IF( kt == nit000 ) THEN
1254         IF(lwp) WRITE(numout,*) 
1255         IF(lwp) WRITE(numout,*) 'tra_kpp : KPP non-local tracer fluxes'
1256         IF(lwp) WRITE(numout,*) '~~~~~~~   '
1257      ENDIF
1258
1259      IF( l_trdtra )   THEN                    !* Save ta and sa trends
1260         ALLOCATE( ztrdt(jpi,jpj,jpk) )   ;    ztrdt(:,:,:) = tsa(:,:,:,jp_tem)
1261         ALLOCATE( ztrds(jpi,jpj,jpk) )   ;    ztrds(:,:,:) = tsa(:,:,:,jp_sal)
1262      ENDIF
1263
1264      ! add non-local temperature and salinity flux ( in convective case only)
1265      DO jk = 1, jpkm1
1266         DO jj = 2, jpjm1 
1267            DO ji = fs_2, fs_jpim1
1268               tsa(ji,jj,jk,jp_tem) =  tsa(ji,jj,jk,jp_tem)                      &
1269                  &                 - (  ghats(ji,jj,jk  ) * avt  (ji,jj,jk  )   & 
1270                  &                    - ghats(ji,jj,jk+1) * avt  (ji,jj,jk+1) ) * wt0(ji,jj) / fse3t(ji,jj,jk)
1271               tsa(ji,jj,jk,jp_sal) =  tsa(ji,jj,jk,jp_sal)                      &
1272                  &                 - (  ghats(ji,jj,jk  ) * fsavs(ji,jj,jk  )   & 
1273                  &                    - ghats(ji,jj,jk+1) * fsavs(ji,jj,jk+1) ) * ws0(ji,jj) / fse3t(ji,jj,jk)
1274            END DO
1275         END DO
1276      END DO
1277
1278      ! save the non-local tracer flux trends for diagnostic
1279      IF( l_trdtra )   THEN
1280         ztrdt(:,:,:) = tsa(:,:,:,jp_tem) - ztrdt(:,:,:)
1281         ztrds(:,:,:) = tsa(:,:,:,jp_sal) - ztrds(:,:,:)
1282!!bug gm jpttdzdf ==> jpttkpp
1283         CALL trd_tra( kt, 'TRA', jp_tem, jptra_trd_zdf, ztrdt )
1284         CALL trd_tra( kt, 'TRA', jp_sal, jptra_trd_zdf, ztrds )
1285         DEALLOCATE( ztrdt )      ;     DEALLOCATE( ztrds )
1286      ENDIF
1287
1288      IF(ln_ctl) THEN
1289         CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' kpp  - Ta: ', mask1=tmask,   &
1290         &             tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
1291      ENDIF
1292
1293   END SUBROUTINE tra_kpp
1294
1295#if defined key_top
1296   !!----------------------------------------------------------------------
1297   !!   'key_top'                                                TOP models
1298   !!----------------------------------------------------------------------
1299   SUBROUTINE trc_kpp( kt )
1300      !!----------------------------------------------------------------------
1301      !!                  ***  ROUTINE trc_kpp  ***
1302      !!
1303      !! ** Purpose :   compute and add to the tracer trend the non-local
1304      !!                tracer flux
1305      !!
1306      !! ** Method  :   ???
1307      !!
1308      !! history :
1309      !!            9.0  ! 2005-11 (G. Madec)  Original code
1310      !!       NEMO 3.3  ! 2010-06 (C. Ethe )  Adapted to passive tracers
1311      !!----------------------------------------------------------------------
1312      USE trc
1313      USE prtctl_trc          ! Print control
1314      !
1315      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
1316      !
1317      INTEGER  ::   ji, jj, jk, jn      ! Dummy loop indices
1318      REAL(wp) ::   ztra, zflx
1319      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   ztrtrd
1320      !!----------------------------------------------------------------------
1321
1322      IF( kt == nit000 ) THEN
1323         IF(lwp) WRITE(numout,*) 
1324         IF(lwp) WRITE(numout,*) 'trc_kpp : KPP non-local tracer fluxes'
1325         IF(lwp) WRITE(numout,*) '~~~~~~~   '
1326      ENDIF
1327
1328      IF( l_trdtrc )  ALLOCATE( ztrtrd(jpi,jpj,jpk) )
1329      !
1330      DO jn = 1, jptra
1331         !
1332         IF( l_trdtrc )  ztrtrd(:,:,:)  = tra(:,:,:,jn)
1333         ! add non-local on passive tracer flux ( in convective case only)
1334         DO jk = 1, jpkm1
1335            DO jj = 2, jpjm1 
1336               DO ji = fs_2, fs_jpim1
1337                  ! Surface tracer flux for non-local term
1338                  zflx = - ( emps(ji,jj) * tra(ji,jj,1,jn) * rcs ) * tmask(ji,jj,1)
1339                  ! compute the trend
1340                  ztra = - ( ghats(ji,jj,jk  ) * fsavs(ji,jj,jk  )   &
1341                  &        - ghats(ji,jj,jk+1) * fsavs(ji,jj,jk+1) ) * zflx / fse3t(ji,jj,jk)
1342                  ! add the trend to the general trend
1343                  tra(ji,jj,jk,jn) = tra(ji,jj,jk,jn)  + ztra
1344               END DO
1345            END DO
1346         END DO
1347         ! save the non-local tracer flux trends for diagnostic
1348         IF( l_trdtrc )  ztrtrd(:,:,:)  = tra(:,:,:,jn) - ztrtrd(:,:,:)
1349         CALL trd_tra( kt, 'TRC', jn, jptra_trd_zdf, ztrtrd(:,:,:,jn) )
1350         !
1351      END DO
1352      IF( l_trdtrc )  DEALLOCATE( ztrtrd )
1353      IF( ln_ctl )   THEN
1354         WRITE(charout, FMT="(' kpp')")  ;  CALL prt_ctl_trc_info(charout)
1355         CALL prt_ctl_trc( tab4d=tra, mask=tmask, clinfo=clname, clinfo2='trd' )
1356      ENDIF
1357      !
1358   END SUBROUTINE trc_kpp
1359
1360#else
1361   !!----------------------------------------------------------------------
1362   !!   NO 'key_top'           DUMMY routine                  No TOP models
1363   !!----------------------------------------------------------------------
1364   SUBROUTINE trc_kpp( kt )         ! Dummy routine
1365      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
1366      WRITE(*,*) 'tra_kpp: You should not have seen this print! error?', kt
1367   END SUBROUTINE trc_kpp
1368#endif
1369
1370   SUBROUTINE zdf_kpp_init
1371      !!----------------------------------------------------------------------
1372      !!                  ***  ROUTINE zdf_kpp_init  ***
1373      !!                     
1374      !! ** Purpose :   Initialization of the vertical eddy diffivity and
1375      !!      viscosity when using a kpp turbulent closure scheme
1376      !!
1377      !! ** Method  :   Read the namkpp namelist and check the parameters
1378      !!      called at the first timestep (nit000)
1379      !!
1380      !! ** input   :   Namlist namkpp
1381      !!----------------------------------------------------------------------
1382      INTEGER  ::   ji, jj, jk     ! dummy loop indices
1383#if ! defined key_kppcustom
1384      INTEGER  ::   jm             ! dummy loop indices     
1385      REAL(wp) ::   zref, zdist    ! tempory scalars
1386#endif
1387#if defined key_kpplktb
1388      REAL(wp) ::   zustar, zucube, zustvk, zeta, zehat   ! tempory scalars
1389#endif
1390      REAL(wp) ::   zhbf           ! tempory scalars
1391      LOGICAL  ::   ll_kppcustom   ! 1st ocean level taken as surface layer
1392      LOGICAL  ::   ll_kpplktb     ! Lookup table for turbul. velocity scales
1393      !!
1394      NAMELIST/namzdf_kpp/ ln_kpprimix, rn_difmiw, rn_difsiw, rn_riinfty, rn_difri, rn_bvsqcon, rn_difcon, nn_ave
1395      !!----------------------------------------------------------------------
1396
1397      REWIND ( numnam )               ! Read Namelist namkpp : K-Profile Parameterisation
1398      READ   ( numnam, namzdf_kpp )
1399
1400      IF(lwp) THEN                    ! Control print
1401         WRITE(numout,*)
1402         WRITE(numout,*) 'zdf_kpp_init : K-Profile Parameterisation'
1403         WRITE(numout,*) '~~~~~~~~~~~~'
1404         WRITE(numout,*) '   Namelist namzdf_kpp : set tke mixing parameters'
1405         WRITE(numout,*) '     Shear instability mixing                      ln_kpprimix = ', ln_kpprimix
1406         WRITE(numout,*) '     max. internal wave viscosity                  rn_difmiw   = ', rn_difmiw
1407         WRITE(numout,*) '     max. internal wave diffusivity                rn_difsiw   = ', rn_difsiw
1408         WRITE(numout,*) '     Richardson Number limit for shear instability rn_riinfty  = ', rn_riinfty
1409         WRITE(numout,*) '     max. shear mixing at Rig = 0                  rn_difri    = ', rn_difri
1410         WRITE(numout,*) '     Brunt-Vaisala squared for max. convection     rn_bvsqcon  = ', rn_bvsqcon
1411         WRITE(numout,*) '     max. mix. in interior convec.                 rn_difcon   = ', rn_difcon
1412         WRITE(numout,*) '     horizontal average flag                       nn_ave      = ', nn_ave
1413      ENDIF
1414
1415      !                              ! allocate zdfkpp arrays
1416      IF( zdf_kpp_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_kpp_init : unable to allocate arrays' )
1417
1418      ll_kppcustom = .FALSE.
1419      ll_kpplktb   = .FALSE.
1420
1421#if defined key_kppcustom
1422      ll_kppcustom = .TRUE.
1423#endif
1424#if defined key_kpplktb
1425      ll_kpplktb   = .TRUE.
1426#endif
1427      IF(lwp) THEN
1428         WRITE(numout,*) '     Lookup table for turbul. velocity scales ll_kpplktb   = ', ll_kpplktb
1429         WRITE(numout,*) '     1st ocean level taken as surface layer   ll_kppcustom = ', ll_kppcustom
1430      ENDIF
1431     
1432      IF( lk_zdfddm) THEN
1433         IF(lwp) THEN
1434            WRITE(numout,*)
1435            WRITE(numout,*) '    Double diffusion mixing on temperature and salinity '
1436            WRITE(numout,*) '    CAUTION : done in routine zdfkpp, not in routine zdfddm '
1437         ENDIF
1438      ENDIF
1439     
1440
1441      !set constants not in namelist
1442      !-----------------------------
1443      Vtc  = rconcv * SQRT( 0.2 / ( rconcs * epsilon ) ) / ( vonk * vonk * Ricr )
1444      rcg  = rcstar * vonk * ( rconcs * vonk * epsilon )**pthird
1445
1446      IF(lwp) THEN
1447         WRITE(numout,*)
1448         WRITE(numout,*) '     Constant value for unreso. turbul. velocity shear Vtc = ', Vtc
1449         WRITE(numout,*) '     Non-dimensional coef. for nonlocal transport      rcg = ', rcg
1450       ENDIF
1451
1452      ! ratt is the attenuation coefficient for solar flux
1453      ! Should be different is s_coordinate
1454      DO jk = 1, jpk
1455         zhbf     = - fsdept(1,1,jk) * hbf
1456         ratt(jk) = 1.0 - ( rabs * EXP( zhbf / xsi1 ) + ( 1.0 - rabs ) * EXP( zhbf / xsi2 ) )       
1457      ENDDO
1458
1459      ! Horizontal average : initialization of weighting arrays
1460      ! -------------------
1461     
1462      SELECT CASE ( nn_ave )
1463
1464      CASE ( 0 )                ! no horizontal average
1465         IF(lwp) WRITE(numout,*) '          no horizontal average on avt, avmu, avmv'
1466         IF(lwp) WRITE(numout,*) '          only in very high horizontal resolution !'
1467         ! weighting mean arrays etmean, eumean and evmean
1468         !           ( 1  1 )                                          ( 1 )
1469         ! avt = 1/4 ( 1  1 )     avmu = 1/2 ( 1  1 )       avmv=  1/2 ( 1 )
1470         !                         
1471         etmean(:,:,:) = 0.e0
1472         eumean(:,:,:) = 0.e0
1473         evmean(:,:,:) = 0.e0
1474         
1475         DO jk = 1, jpkm1
1476            DO jj = 2, jpjm1
1477               DO ji = 2, jpim1   ! vector opt.
1478                  etmean(ji,jj,jk) = tmask(ji,jj,jk)                     &
1479                  &  / MAX( 1.,  umask(ji-1,jj  ,jk) + umask(ji,jj,jk)   &
1480                  &            + vmask(ji  ,jj-1,jk) + vmask(ji,jj,jk)  )
1481                 
1482                  eumean(ji,jj,jk) = umask(ji,jj,jk)                     &
1483                  &  / MAX( 1.,  tmask(ji,jj,jk) + tmask(ji+1,jj  ,jk)  )
1484
1485                  evmean(ji,jj,jk) = vmask(ji,jj,jk)                     &
1486                  &  / MAX( 1.,  tmask(ji,jj,jk) + tmask(ji  ,jj+1,jk)  )
1487               END DO
1488            END DO
1489         END DO
1490
1491      CASE ( 1 )                ! horizontal average
1492         IF(lwp) WRITE(numout,*) '          horizontal average on avt, avmu, avmv'
1493         ! weighting mean arrays etmean, eumean and evmean
1494         !           ( 1/2  1  1/2 )              ( 1/2  1/2 )             ( 1/2  1  1/2 )
1495         ! avt = 1/8 ( 1    2  1   )   avmu = 1/4 ( 1    1   )   avmv= 1/4 ( 1/2  1  1/2 )
1496         !           ( 1/2  1  1/2 )              ( 1/2  1/2 )
1497         etmean(:,:,:) = 0.e0
1498         eumean(:,:,:) = 0.e0
1499         evmean(:,:,:) = 0.e0
1500         
1501         DO jk = 1, jpkm1
1502            DO jj = 2, jpjm1
1503               DO ji = fs_2, fs_jpim1   ! vector opt.
1504                  etmean(ji,jj,jk) = tmask(ji, jj,jk)                           &
1505                     & / MAX( 1., 2.* tmask(ji,jj,jk)                           &
1506                     &      +.5 * ( tmask(ji-1,jj+1,jk) + tmask(ji-1,jj-1,jk)   &
1507                     &             +tmask(ji+1,jj+1,jk) + tmask(ji+1,jj-1,jk) ) &
1508                     &      +1. * ( tmask(ji-1,jj  ,jk) + tmask(ji  ,jj+1,jk)   &
1509                     &             +tmask(ji  ,jj-1,jk) + tmask(ji+1,jj  ,jk) ) )
1510                 
1511                  eumean(ji,jj,jk) = umask(ji,jj,jk)                        &
1512                     &  / MAX( 1.,   tmask(ji,jj  ,jk) + tmask(ji+1,jj  ,jk)   &
1513                     &       +.5 * ( tmask(ji,jj-1,jk) + tmask(ji+1,jj-1,jk)   &
1514                     &              +tmask(ji,jj+1,jk) + tmask(ji+1,jj+1,jk) )  )
1515                 
1516                  evmean(ji,jj,jk) = vmask(ji,jj,jk)                        &
1517                     &  / MAX( 1.,   tmask(ji  ,jj,jk) + tmask(ji  ,jj+1,jk)   &
1518                     &       +.5 * ( tmask(ji-1,jj,jk) + tmask(ji-1,jj+1,jk)   &
1519                     &              +tmask(ji+1,jj,jk) + tmask(ji+1,jj+1,jk) )  )
1520               END DO
1521            END DO
1522         END DO
1523
1524      CASE DEFAULT
1525         WRITE(ctmp1,*) '          bad flag value for nn_ave = ', nn_ave
1526         CALL ctl_stop( ctmp1 )
1527
1528      END SELECT
1529 
1530      ! Initialization of vertical eddy coef. to the background value
1531      ! -------------------------------------------------------------
1532      DO jk = 1, jpk
1533         avt (:,:,jk) = avtb(jk) * tmask(:,:,jk)
1534         avmu(:,:,jk) = avmb(jk) * umask(:,:,jk)
1535         avmv(:,:,jk) = avmb(jk) * vmask(:,:,jk)
1536      END DO
1537
1538      ! zero the surface flux for non local term and kpp mixed layer depth
1539      ! ------------------------------------------------------------------
1540      ghats(:,:,:) = 0.
1541      wt0  (:,:  ) = 0.
1542      ws0  (:,:  ) = 0.
1543      hkpp (:,:  ) = 0. ! just a diagnostic (not essential)
1544
1545#if ! defined key_kppcustom
1546      ! compute arrays (del, wz) for reference mean values
1547      ! (increase speed for vectorization key_kppcustom not defined)
1548      del(1:jpk, 1:jpk) = 0.
1549      DO jk = 1, jpk
1550         zref = epsilon * fsdept(1,1,jk)   
1551         DO jm = 1 , jpk
1552            zdist = zref - fsdepw(1,1,jm)   
1553            IF( zdist > 0.  ) THEN
1554               del(jk,jm) = MIN( zdist, fse3t(1,1,jm) ) / zref   
1555            ELSE
1556               del(jk,jm) = 0.
1557            ENDIF
1558         ENDDO
1559      ENDDO
1560#endif
1561
1562#if defined key_kpplktb
1563      ! build lookup table for turbulent velocity scales
1564      dezehat = ( dehatmax - dehatmin ) / nilktbm1
1565      deustar = ( ustmax   - ustmin   ) / njlktbm1
1566 
1567      DO jj = 1, njlktb
1568         zustar = ( jj - 1) * deustar + ustmin
1569         zustvk = vonk * zustar 
1570         zucube = zustar * zustar * zustar 
1571         DO ji = 1 , nilktb
1572            zehat = ( ji - 1 ) * dezehat + dehatmin
1573            zeta   = zehat / ( zucube + epsln )
1574            IF( zehat >= 0 ) THEN             ! Stable case
1575               wmlktb(ji,jj) = zustvk / ABS( 1.0 + rconc1 * zeta + epsln )                       
1576               wslktb(ji,jj) = wmlktb(ji,jj)
1577            ELSE                                ! Unstable case
1578               IF( zeta > rzetam ) THEN
1579                  wmlktb(ji,jj) = zustvk * ABS( 1.0    - rconc2 * zeta )**pfourth
1580               ELSE
1581                  wmlktb(ji,jj) = zustvk * ABS( rconam - rconcm * zeta )**pthird
1582               ENDIF
1583               
1584               IF( zeta > rzetas ) THEN
1585                  wslktb(ji,jj) = zustvk * SQRT( ABS( 1.0 - rconc2 * zeta ) )
1586               ELSE
1587                  wslktb(ji,jj) = zustvk * ABS( rconas - rconcs * zeta )**pthird
1588               ENDIF
1589            ENDIF
1590         END DO
1591      END DO
1592#endif
1593      !
1594   END SUBROUTINE zdf_kpp_init
1595
1596#else
1597   !!----------------------------------------------------------------------
1598   !!   Dummy module :                                        NO KPP scheme
1599   !!----------------------------------------------------------------------
1600   LOGICAL, PUBLIC, PARAMETER ::   lk_zdfkpp = .FALSE.   !: KPP flag
1601CONTAINS
1602   SUBROUTINE zdf_kpp_init           ! Dummy routine
1603      WRITE(*,*) 'zdf_kpp_init: You should not have seen this print! error?'
1604   END SUBROUTINE zdf_kpp_init
1605   SUBROUTINE zdf_kpp( kt )          ! Dummy routine
1606      WRITE(*,*) 'zdf_kpp: You should not have seen this print! error?', kt
1607   END SUBROUTINE zdf_kpp
1608   SUBROUTINE tra_kpp( kt )          ! Dummy routine
1609      WRITE(*,*) 'tra_kpp: You should not have seen this print! error?', kt
1610   END SUBROUTINE tra_kpp
1611   SUBROUTINE trc_kpp( kt )          ! Dummy routine
1612      WRITE(*,*) 'trc_kpp: You should not have seen this print! error?', kt
1613   END SUBROUTINE trc_kpp
1614#endif
1615
1616   !!======================================================================
1617END MODULE zdfkpp
Note: See TracBrowser for help on using the repository browser.