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

source: trunk/NEMO/OPA_SRC/ZDF/zdfkpp.F90 @ 255

Last change on this file since 255 was 255, checked in by opalod, 19 years ago

nemo_v1_update_002 : CT : Integration of the KPP turbulent closure scheme

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