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 @ 1695

Last change on this file since 1695 was 1695, checked in by smasson, 14 years ago

wind stress module directly at T-point, see ticket:577

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