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

source: branches/NERC/dev_r3874_FASTNEt/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfkpp.F90 @ 6736

Last change on this file since 6736 was 6736, checked in by jamesharle, 8 years ago

FASTNEt code modifications

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