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

Last change on this file since 1601 was 1601, checked in by ctlod, 15 years ago

Doctor naming of OPA namelist variables , see ticket: #526

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 73.3 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, ztau, 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            ztx           = 0.5 * ( utau(ji,jj) + utau(ji - 1, jj    ) )
472            zty           = 0.5 * ( vtau(ji,jj) + vtau(ji    , jj - 1) )
473            ztau          = SQRT( ztx * ztx + zty * zty )
474            zustar(ji,jj) = SQRT( ztau / ( zrhos +  epsln ) )
475         ENDDO
476      ENDDO
477
478!CDIR NOVERRCHK 
479      !                                               ! ===============
480      DO jj = 2, jpjm1                                 !  Vertical slab
481         !                                             ! ===============
482         
483         !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
484         ! II Compute Boundary layer mixing coef. and diagnose the new boundary layer depth
485         !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
486         
487         ! Initialization
488         jkmax       = 0
489         zdept (:,:) = 0.
490         zdepw (:,:) = 0.
491         zriblk(:,:) = 0.
492         zmoek (:,:) = 0.
493         zvisc (:,:) = 0.
494         zdift (:,:) = 0.
495#if defined key_zdfddm
496         zdifs (:,:) = 0.
497#endif
498         zmask (:,:) = 0.
499         DO ji = fs_2, fs_jpim1
500            zria(ji ) = 0.
501            ! Maximum boundary layer depth
502            ikbot     = mbathy(ji,jj) - 1 ! ikbot is the last T point in the water
503            zhmax(ji) = fsdept(ji,jj,ikbot) - 0.001     
504            ! Compute Monin obukhov length scale at the surface and Ekman depth:
505            zbuofdep   = zBo(ji,jj) + zBosol(ji,jj) * ratt(1)
506            zekman(ji) = rcekman * zustar(ji,jj) / ( ABS( ff(ji,jj) ) + epsln )
507            zucube     = zustar(ji,jj) * zustar(ji,jj) * zustar(ji,jj) 
508            zmoa(ji)   = zucube / ( vonk * ( zbuofdep + epsln ) )   
509#if defined key_c1d
510            ! store the surface buoyancy forcing
511            zstabl        = 0.5 + SIGN( 0.5, zbuofdep )
512            buof(ji,jj,1) = zbuofdep * tmask(ji,jj,1)
513            ! store the moning-oboukov length scale at surface
514            zmob          = zstabl * zmoa(ji) + ( 1.0 - zstabl ) * fsdept(ji,jj,1)
515            mols(ji,jj,1) = MIN( zmob , zhmax(ji) ) * tmask(ji,jj,1)
516            ! store Ekman depth
517            zek           = zstabl * zekman(ji) + ( 1.0 - zstabl ) * fsdept(ji,jj,1) 
518            ekdp(ji,jj )  = MIN( zek , zhmax(ji) ) * tmask(ji,jj,1) 
519#endif
520         END DO     
521         ! Compute the pipe
522         ! ---------------------
523         DO jk = 2, jpkm1
524            DO ji = fs_2, fs_jpim1
525               ! Compute bfsfc = Bo + radiative contribution down to hbf*depht
526               zbuofdep = zBo(ji,jj) + zBosol(ji,jj) * ratt(jk)
527               ! Flag (zstabl  = 1) if positive forcing
528               zstabl   =  0.5 + SIGN(  0.5, zbuofdep)
529
530               !   Compute bulk richardson number zrib at depht
531               !-------------------------------------------------------
532               !                           [Br - B(d)] * d         zrinum
533               !             Rib(z) = ----------------------- = -------------
534               !                       |Vr - V(d)|^2 + Vt(d)^2   zdVsq + zVtsq
535               !
536               ! First compute zt,zs,zu,zv = means in the surface layer < epsilon*depht 
537               ! Else surface values are taken at the first T level.
538               ! For stability, resolved vertical shear is computed with "before velocities".
539               zref = epsilon * fsdept(ji,jj,jk)
540#if defined key_kppcustom
541               ! zref = gdept(1)
542               zref = fsdept(ji,jj,1)
543               zt   = tn(ji,jj,1)
544               zs   = sn(ji,jj,1)
545               zrh  = rhop(ji,jj,1)
546               zu   = ( ub(ji,jj,1) + ub(ji - 1,jj    ,1) ) / MAX( 1. , umask(ji,jj,1) + umask(ji - 1,jj   ,1) )
547               zv   = ( vb(ji,jj,1) + vb(ji    ,jj - 1,1) ) / MAX( 1. , vmask(ji,jj,1) + vmask(ji   ,jj - 1,1) )
548#else
549               zt   = 0.
550               zs   = 0.
551               zu   = 0.
552               zv   = 0.
553               zrh  = 0.
554               ! vertically integration over the upper epsilon*gdept(jk) ; del () array is computed once in zdf_kpp_init
555               DO jm = 1, jpkm1
556                  zt   = zt  + del(jk,jm) * tn(ji,jj,jm)
557                  zs   = zs  + del(jk,jm) * sn(ji,jj,jm)
558                  zu   = zu  + 0.5 * del(jk,jm) &
559                     &            * ( ub(ji,jj,jm) + ub(ji - 1,jj,jm) ) &
560                     &            / MAX( 1. , umask(ji,jj,jm) + umask(ji - 1,jj,jm) )
561                  zv   = zv  + 0.5 * del(jk,jm) &
562                     &            * ( vb(ji,jj,jm) + vb(ji,jj - 1,jm) ) &
563                     &            / MAX( 1. , vmask(ji,jj,jm) + vmask(ji,jj - 1,jm) )
564                  zrh  = zrh + del(jk,jm) * rhop(ji,jj,jm)
565               END DO
566#endif
567               zsr = SQRT( ABS( sn(ji,jj,jk) ) )
568               ! depth
569               zh = fsdept(ji,jj,jk)
570               ! compute compression terms on density
571               ze  = ( -3.508914e-8*zt-1.248266e-8 ) *zt-2.595994e-6
572               zbw = (  1.296821e-6*zt-5.782165e-9 ) *zt+1.045941e-4
573               zb  = zbw + ze * zs
574               
575               zd  = -2.042967e-2
576               zc  =   (-7.267926e-5*zt+2.598241e-3 ) *zt+0.1571896
577               zaw = ( ( 5.939910e-6*zt+2.512549e-3 ) *zt-0.1028859 ) *zt - 4.721788
578               za  = ( zd*zsr + zc ) *zs + zaw
579               
580               zb1 =   (-0.1909078*zt+7.390729 ) *zt-55.87545
581               za1 = ( ( 2.326469e-3*zt+1.553190)*zt-65.00517 ) *zt+1044.077
582               zkw = ( ( (-1.361629e-4*zt-1.852732e-2 ) *zt-30.41638 ) *zt + 2098.925 ) *zt+190925.6
583               zk0 = ( zb1*zsr + za1 )*zs + zkw
584               zcomp =   1.0 - zh / ( zk0 - zh * ( za - zh * zb ) )
585               
586#if defined key_kppcustom
587               ! potential density of water(zrh = zt,zs at level jk):
588               zrhdr = zrh / zcomp
589#else
590               ! potential density of water(ztref,zsref at level jk):
591               ! compute volumic mass pure water at atm pressure
592               IF ( nn_eos < 1 ) THEN
593                  zr1= ( ( ( ( 6.536332e-9*zt-1.120083e-6 )*zt+1.001685e-4)*zt   &
594                     &       -9.095290e-3 )*zt+6.793952e-2 )*zt+999.842594
595                  ! seawater volumic mass atm pressure
596                  zr2= ( ( ( 5.3875e-9*zt-8.2467e-7 ) *zt+7.6438e-5 ) *zt   &
597                     &   -4.0899e-3 ) *zt+0.824493
598                  zr3= ( -1.6546e-6*zt+1.0227e-4 ) *zt-5.72466e-3
599                  zr4= 4.8314e-4             
600                  ! potential volumic mass (reference to the surface)
601                  zrhop= ( zr4*zs + zr3*zsr + zr2 ) *zs + zr1                 
602                  zrhdr = zrhop / zcomp
603               ELSE
604                  zrhdr = zrh / zcomp
605               ENDIF
606#endif
607               
608               ! potential density of ambiant water at level jk :
609               zrhd   = ( rhd(ji,jj,jk) * rau0 + rau0 ) 
610               
611               ! And now the Rib number numerator .
612               zrinum = grav * ( zrhd - zrhdr ) / rau0
613               zrinum = zrinum * ( fsdept(ji,jj,jk) - zref ) * tmask(ji,jj,jk)
614           
615               ! Resolved shear contribution to Rib at depth T-point (zdVsq)
616               ztx    =   ( ub( ji , jj ,jk)   +  ub(ji - 1, jj ,jk) ) &
617                  &     / MAX( 1. , umask( ji , jj ,jk) + umask(ji - 1, jj ,jk) )   
618               zty    =   ( vb( ji , jj ,jk)   +  vb(ji  ,jj - 1,jk) ) &
619                  &     / MAX( 1., vmask( ji , jj ,jk) + vmask(ji  ,jj - 1,jk) ) 
620               
621               zdVsq  = ( zu - ztx ) * ( zu - ztx ) + ( zv - zty ) * ( zv - zty )
622               
623               ! Scalar turbulent velocity scale zws for hbl=gdept
624               zscale = zstabl + ( 1.0 - zstabl ) * epsilon
625               zehat  = vonk * zscale * fsdept(ji,jj,jk) * zbuofdep
626               zucube = zustar(ji,jj) * zustar(ji,jj) * zustar(ji,jj)             
627               zeta   = zehat / ( zucube + epsln )
628               
629               IF( zehat > 0. ) THEN
630                  ! Stable case
631                  zws  = vonk * zustar(ji,jj) / ( 1.0 + rconc1 * zeta )
632               ELSE
633                  ! Unstable case
634#if defined key_kpplktb
635                  ! use lookup table
636                  zd     = zehat - dehatmin
637                  il     = INT( zd / dezehat )
638                  il     = MIN( il, nilktbm1 )
639                  il     = MAX( il, 1 )
640                 
641                  ud     = zustar(ji,jj) - ustmin
642                  jl     = INT( ud / deustar )
643                  jl     = MIN( jl, njlktbm1 )
644                  jl     = MAX( jl, 1 )
645                 
646                  zfrac  = zd / dezehat - FLOAT( il ) 
647                  ufrac  = ud / deustar - FLOAT( jl )
648                  zwas   = ( 1. - zfrac ) * wslktb(il,jl+1) + zfrac * wslktb(il+1,jl+1)
649                  zwbs   = ( 1. - zfrac ) * wslktb(il,jl  ) + zfrac * wslktb(il+1,jl  )
650                  !
651                  zws    = ( 1. - ufrac ) * zwbs + ufrac * zwas
652#else
653                  ! use analytical functions:
654                  zcons  = 0.5 + SIGN( 0.5 , ( rzetas - zeta ) )
655                  zwcons = vonk * zustar(ji,jj) * ( ( ABS( rconas - rconcs * zeta ) )**pthird ) 
656                  zwsun  = vonk * zustar(ji,jj) * SQRT( ABS ( 1.0 - rconc2 * zeta ) )
657                  !
658                  zws    = zcons * zwcons +  ( 1.0 - zcons) * zwsun
659#endif
660               ENDIF
661               
662               ! Turbulent shear contribution to Rib (zVtsq) bv frequency at levels  ( ie T-point jk)
663               zrn2   = 0.5 * ( rn2(ji,jj,jk) + rn2(ji,jj,jk+1) )   
664               zbvzed = SQRT( ABS( zrn2 ) ) 
665               zVtsq  = fsdept(ji,jj,jk) * zws * zbvzed  * Vtc
666               
667               ! Finally, the bulk Richardson number at depth fsdept(i,j,k)
668               zrib  = zrinum   / ( zdVsq + zVtsq + epsln )
669 
670               ! Find subscripts around the boundary layer depth, build the pipe
671               ! ----------------------------------------------------------------
672
673               ! Flag (zflagri = 1) if zrib < Ricr 
674               zflagri = 0.5 + SIGN( 0.5, ( Ricr - zrib ) )
675               !  Flag (zflagh  = 1) if still within overall boundary layer
676               zflagh  = 0.5 + SIGN( 0.5, ( fsdept(ji,jj,1) - zdept(ji,2) ) )
677               
678               ! Ekman layer depth
679               zek     = zstabl * zekman(ji) + ( 1.0 - zstabl ) * zhmax(ji)
680               zflag   = 0.5 + SIGN( 0.5, ( zek - fsdept(ji,jj,jk-1) ) )
681               zek     = zflag * zek + ( 1.0 - zflag ) * zhmax(ji)
682               zflagek = 0.5 + SIGN( 0.5, ( zek - fsdept(ji,jj,jk) ) )
683               ! Flag (zflagmo = 1) if still within stable Monin-Obukhov and in water
684               zmob    = zucube / ( vonk * ( zbuofdep + epsln ) ) 
685               ztemp   = zstabl * zmob + ( 1.0 - zstabl) * zhmax(ji)
686               ztemp   = MIN( ztemp , zhmax(ji) ) 
687               zflagmo = 0.5 + SIGN( 0.5, ( ztemp - fsdept(ji,jj,jk) ) )             
688
689               ! No limitation by Monin Obukhov or Ekman depths:
690!               zflagek = 1.0
691!               zflagmo = 0.5 + SIGN( 0.5, ( zhmax(ji) - fsdept(ji,jj,jk) ) )
692
693               ! Load  pipe via zflagkb  for later calculations
694               ! Flag (zflagkb = 1) if zflagh = 1 and (zflagri = 0 or zflagek = 0 or zflagmo = 0)
695               zflagkb = zflagh * ( 1.0 - ( zflagri * zflagek * zflagmo ) )
696
697               zmask(ji,jk) = zflagh
698               jkp2         = MIN( jk+2 , ikbot )
699               jkm1         = MAX( jk-1 , 2 )
700               jkmax        = MAX( jkmax, jk * INT( REAL( zflagh+epsln ) ) ) 
701
702               zdept(ji,1)  = zdept(ji,1) + zflagkb * fsdept(ji,jj,jk-1) 
703               zdept(ji,2)  = zdept(ji,2) + zflagkb * fsdept(ji,jj,jk  ) 
704               zdept(ji,3)  = zdept(ji,3) + zflagkb * fsdept(ji,jj,jk+1) 
705
706               zdepw(ji,1)  = zdepw(ji,1) + zflagkb * fsdepw(ji,jj,jk-1) 
707               zdepw(ji,2)  = zdepw(ji,2) + zflagkb * fsdepw(ji,jj,jk  ) 
708               zdepw(ji,3)  = zdepw(ji,3) + zflagkb * fsdepw(ji,jj,jk+1)
709               zdepw(ji,4)  = zdepw(ji,4) + zflagkb * fsdepw(ji,jj,jkp2) 
710
711               zriblk(ji,1) = zriblk(ji,1) + zflagkb * zria(ji)
712               zriblk(ji,2) = zriblk(ji,2) + zflagkb * zrib
713
714               zmoek (ji,0) = zmoek (ji,0) + zflagkb * zek
715               zmoek (ji,1) = zmoek (ji,1) + zflagkb * zmoa(ji)
716               zmoek (ji,2) = zmoek (ji,2) + zflagkb * ztemp 
717               ! Save Monin Obukhov depth
718               zmoa  (ji)   = zmob
719           
720               zvisc(ji,1) = zvisc(ji,1) + zflagkb * avmu(ji,jj,jkm1)
721               zvisc(ji,2) = zvisc(ji,2) + zflagkb * avmu(ji,jj,jk  )
722               zvisc(ji,3) = zvisc(ji,3) + zflagkb * avmu(ji,jj,jk+1)
723               zvisc(ji,4) = zvisc(ji,4) + zflagkb * avmu(ji,jj,jkp2)
724               
725               zdift(ji,1) = zdift(ji,1) + zflagkb * avt (ji,jj,jkm1)
726               zdift(ji,2) = zdift(ji,2) + zflagkb * avt (ji,jj,jk  )
727               zdift(ji,3) = zdift(ji,3) + zflagkb * avt (ji,jj,jk+1)
728               zdift(ji,4) = zdift(ji,4) + zflagkb * avt (ji,jj,jkp2)
729
730#if defined key_zdfddm
731               zdifs(ji,1) = zdifs(ji,1) + zflagkb * avs (ji,jj,jkm1)
732               zdifs(ji,2) = zdifs(ji,2) + zflagkb * avs (ji,jj,jk  )
733               zdifs(ji,3) = zdifs(ji,3) + zflagkb * avs (ji,jj,jk+1)
734               zdifs(ji,4) = zdifs(ji,4) + zflagkb * avs (ji,jj,jkp2)
735#endif               
736               ! Save the Richardson number
737               zria  (ji)   = zrib 
738#if defined key_c1d
739               ! store buoyancy length scale
740               buof(ji,jj,jk) = zbuofdep * tmask(ji,jj,jk) 
741               ! store Monin Obukhov
742               zmob           = zstabl * zmob + ( 1.0 - zstabl) * fsdept(ji,jj,1)
743               mols(ji,jj,jk) = MIN( zmob , zhmax(ji) ) * tmask(ji,jj,jk) 
744               ! Bulk Richardson number
745               rib(ji,jj,jk)  = zrib * tmask(ji,jj,jk)             
746#endif               
747            END DO
748         END DO
749         !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
750         ! III PROCESS THE PIPE
751         !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
752         
753         DO ji = fs_2, fs_jpim1 
754           
755            ! Find the boundary layer depth zhbl
756            ! ----------------------------------------
757           
758            ! Interpolate monin Obukhov and critical Ri mumber depths   
759            ztemp = zdept(ji,2) - zdept(ji,1)
760            zflag = ( Ricr - zriblk(ji,1) ) / ( zriblk(ji,2) - zriblk(ji,1)  + epsln )
761            zhrib = zdept(ji,1) + zflag * ztemp     
762
763            IF( zriblk(ji,2) < Ricr ) zhrib = zhmax(ji) 
764         
765            IF( zmoek(ji,2) < zdept(ji,2) ) THEN
766               IF ( zmoek(ji,1) < 0. ) THEN
767                  zmob = zdept(ji,2) - epsln
768               ELSE
769                  zmob = ztemp + zmoek(ji,1) - zmoek(ji,2)
770                  zmob = ( zmoek(ji,1) * zdept(ji,2) - zmoek(ji,2) * zdept(ji,1) ) / zmob
771                  zmob = MAX( zmob , zdept(ji,1) + epsln )               
772               ENDIF
773            ELSE           
774               zmob = zhmax(ji) 
775            ENDIF
776            ztemp   = MIN( zmob , zmoek(ji,0) )
777                         
778            ! Finally, the boundary layer depth, zhbl
779            zhbl(ji) = MAX( fsdept(ji,jj,1) + epsln, MIN( zhrib , ztemp ) )
780           
781            ! Save hkpp for further diagnostics (optional)
782            hkpp(ji,jj) = zhbl(ji) * tmask(ji,jj,1) 
783         
784            ! Correct mask if zhbl < fsdepw(ji,jj,2) for no viscosity/diffusivity enhancement at fsdepw(ji,jj,2)
785            !     zflag = 1 if zhbl(ji) > fsdepw(ji,jj,2)
786            IF( zhbl(ji) < fsdepw(ji,jj,2) ) zmask(ji,2) = 0.
787         
788           
789            !  Velocity scales at depth zhbl
790            ! -----------------------------------
791           
792            !  Compute bouyancy forcing down to zhbl
793            ztemp    = -hbf * zhbl(ji)
794            zatt1    = 1.0 - ( rabs * EXP( ztemp / xsi1 ) + ( 1.0 - rabs ) * EXP( ztemp / xsi2 ) )
795            zbuofdep = zBo(ji,jj) + zBosol(ji,jj) * zatt1
796            zstabl   = 0.5 + SIGN( 0.5 , zbuofdep ) 
797
798            zbuofdep = zbuofdep + zstabl * epsln
799
800            zscale = zstabl + ( 1.0 - zstabl ) * epsilon         
801            zehat  = vonk * zscale * zhbl(ji) * zbuofdep
802            zucube = zustar(ji,jj) * zustar(ji,jj) * zustar(ji,jj)             
803            zeta   = zehat / ( zucube + epsln )
804           
805            IF( zehat > 0. ) THEN
806               ! Stable case
807               zws  = vonk * zustar(ji,jj) / ( 1.0 + rconc1 * zeta )
808               zwm  = zws
809            ELSE
810               ! Unstable case
811#if defined key_kpplktb
812               ! use lookup table
813               zd     = zehat - dehatmin
814               il     = INT( zd / dezehat )
815               il     = MIN( il, nilktbm1 )
816               il     = MAX( il, 1 )
817               
818               ud     = zustar(ji,jj) - ustmin
819               jl     = INT( ud / deustar )
820               jl     = MIN( jl, njlktbm1 )
821               jl     = MAX( jl, 1 )
822               
823               zfrac  = zd / dezehat - FLOAT( il ) 
824               ufrac  = ud / deustar - FLOAT( jl )
825               zwas   = ( 1. - zfrac ) * wslktb(il,jl+1) + zfrac * wslktb(il+1,jl+1)
826               zwbs   = ( 1. - zfrac ) * wslktb(il,jl  ) + zfrac * wslktb(il+1,jl  )
827               zwam   = ( 1. - zfrac ) * wmlktb(il,jl+1) + zfrac * wmlktb(il+1,jl+1)
828               zwbm   = ( 1. - zfrac ) * wmlktb(il,jl  ) + zfrac * wmlktb(il+1,jl  )
829               !
830               zws    = ( 1. - ufrac ) * zwbs + ufrac * zwas
831               zwm    = ( 1. - ufrac ) * zwbm + ufrac * zwam
832#else
833               ! use analytical functions
834               zconm  = 0.5 + SIGN( 0.5, ( rzetam - zeta) )
835               zcons  = 0.5 + SIGN( 0.5, ( rzetas - zeta) )
836               
837               ! Momentum : zeta < rzetam (zconm = 1)
838               ! Scalars  : zeta < rzetas (zcons = 1)
839               zwconm = zustar(ji,jj) * vonk * ( ( ABS( rconam - rconcm * zeta) )**pthird )
840               zwcons = zustar(ji,jj) * vonk * ( ( ABS( rconas - rconcs * zeta) )**pthird )
841               
842               ! Momentum : rzetam <= zeta < 0 (zconm = 0)
843               ! Scalars  : rzetas <= zeta < 0 (zcons = 0) 
844               zwmun  = SQRT( ABS( 1.0 - rconc2 * zeta ) )
845               zwsun  = vonk * zustar(ji,jj) * zwmun
846               zwmun  = vonk * zustar(ji,jj) * SQRT(zwmun)
847               !
848               zwm    = zconm * zwconm + ( 1.0 - zconm ) * zwmun
849               zws    = zcons * zwcons + ( 1.0 - zcons ) * zwsun
850               
851#endif
852            ENDIF
853           
854           
855            ! Viscosity, diffusivity values and derivatives at h
856            ! --------------------------------------------------------
857           
858            ! check between at which interfaces is located zhbl(ji)
859            ! ztemp = 1, zdepw(ji,2) < zhbl <  zdepw(ji,3)
860            ! ztemp = 0, zdepw(ji,1) < zhbl <  zdepw(ji,2)
861            ztemp  =  0.5 + SIGN( 0.5, ( zhbl(ji) - zdepw(ji,2) ) ) 
862            zdep21 =   zdepw(ji,2) - zdepw(ji,1) + epsln
863            zdep32 =   zdepw(ji,3) - zdepw(ji,2) + epsln
864            zdep43 =   zdepw(ji,4) - zdepw(ji,3) + epsln 
865           
866            ! Compute R as in LMD94, eq D5b
867            zdelta =  ( zhbl(ji) - zdepw(ji,2) ) *         ztemp   / zdep32   &
868               &    + ( zhbl(ji) - zdepw(ji,1) ) * ( 1.0 - ztemp ) / zdep21 
869           
870            ! Compute the vertical derivative of viscosities (zdzh) at z=zhbl(ji)
871            zdzup  =  ( zvisc(ji,2) - zvisc(ji,3) ) *         ztemp   / zdep32 &
872               &    + ( zvisc(ji,1) - zvisc(ji,2) ) * ( 1.0 - ztemp ) / zdep21
873           
874            zdzdn  =  ( zvisc(ji,3) - zvisc(ji,4) ) *         ztemp   / zdep43 &
875               &    + ( zvisc(ji,2) - zvisc(ji,3) ) * ( 1.0 - ztemp ) / zdep32
876           
877            ! LMD94, eq D5b :         
878            zdzh   = ( 1.0 - zdelta ) * zdzup + zdelta * zdzdn
879            zdzh   = MAX( zdzh , 0. )
880           
881            ! Compute viscosities (zvath) at z=zhbl(ji), LMD94 eq D5a
882            zvath  =          ztemp   * ( zvisc(ji,3) + zdzh * ( zdepw(ji,3) - zhbl(ji) ) ) &
883               &    + ( 1.0 - ztemp ) * ( zvisc(ji,2) + zdzh * ( zdepw(ji,2) - zhbl(ji) ) )
884           
885            ! Compute G (zgat1) and its derivative (zdat1) at z=hbl(ji), LMD94 eq 18
886           
887            ! Vertical derivative of velocity scale divided by velocity scale squared at z=hbl(ji)
888            ! (non zero only in stable conditions)
889            zflag  =  -zstabl * rconc1 * zbuofdep / ( zucube * zustar(ji,jj) + epsln )
890           
891            ! G at its derivative at z=hbl:
892            zgat1  = zvath  / ( zhbl(ji) * ( zwm + epsln )  )
893            zdat1  = -zdzh  / ( zwm + epsln ) -  zflag * zvath / zhbl(ji)
894           
895            ! G coefficients, LMD94 eq 17
896            za2m(ji) = -2.0 + 3.0 * zgat1 - zdat1
897            za3m(ji) =  1.0 - 2.0 * zgat1 + zdat1
898
899           
900            ! Compute the vertical derivative of temperature diffusivities (zdzh) at z=zhbl(ji)
901            zdzup  =  ( zdift(ji,2) - zdift(ji,3) ) *         ztemp   / zdep32 &
902               &    + ( zdift(ji,1) - zdift(ji,2) ) * ( 1.0 - ztemp ) / zdep21
903           
904            zdzdn  =  ( zdift(ji,3) - zdift(ji,4) ) *         ztemp   / zdep43 &
905               &    + ( zdift(ji,2) - zdift(ji,3) ) * ( 1.0 - ztemp ) / zdep32
906           
907            ! LMD94, eq D5b :         
908            zdzh   = ( 1.0 - zdelta ) * zdzup + zdelta * zdzdn
909            zdzh   = MAX( zdzh , 0. )
910           
911           
912            ! Compute diffusivities (zvath) at z=zhbl(ji), LMD94 eq D5a
913            zvath  =          ztemp   * ( zdift(ji,3) + zdzh * ( zdepw(ji,3) - zhbl(ji) ) ) &
914               &    + ( 1.0 - ztemp ) * ( zdift(ji,2) + zdzh * ( zdepw(ji,2) - zhbl(ji) ) )
915                       
916            ! G at its derivative at z=hbl:
917            zgat1  = zvath  / ( zhbl(ji) * ( zws + epsln )  )
918            zdat1  = -zdzh  / ( zws + epsln ) -  zflag * zvath / zhbl(ji)
919           
920            ! G coefficients, LMD94 eq 17
921            za2t(ji) = -2.0 + 3.0 * zgat1 - zdat1
922            za3t(ji) =  1.0 - 2.0 * zgat1 + zdat1
923
924#if defined key_zdfddm
925            ! Compute the vertical derivative of salinities diffusivities (zdzh) at z=zhbl(ji)
926            zdzup  =  ( zdifs(ji,2) - zdifs(ji,3) ) *         ztemp   / zdep32 &
927               &    + ( zdifs(ji,1) - zdifs(ji,2) ) * ( 1.0 - ztemp ) / zdep21
928           
929            zdzdn  =  ( zdifs(ji,3) - zdifs(ji,4) ) *         ztemp   / zdep43 &
930               &    + ( zdifs(ji,2) - zdifs(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            ! Compute diffusivities (zvath) at z=zhbl(ji), LMD94 eq D5a
937            zvath  =          ztemp   * ( zdifs(ji,3) + zdzh * ( zdepw(ji,3) - zhbl(ji) ) ) &
938               &    + ( 1.0 - ztemp ) * ( zdifs(ji,2) + zdzh * ( zdepw(ji,2) - zhbl(ji) ) )
939                       
940            ! G at its derivative at z=hbl:
941            zgat1  = zvath  / ( zhbl(ji) * ( zws + epsln )  )
942            zdat1  = -zdzh  / ( zws + epsln ) -  zflag * zvath / zhbl(ji)
943           
944            ! G coefficients, LMD94 eq 17
945            za2s(ji) = -2.0 + 3.0 * zgat1 - zdat1
946            za3s(ji) =  1.0 - 2.0 * zgat1 + zdat1
947#endif
948
949            !-------------------turn off interior matching here------
950            !          za2(ji,1) = -2.0
951            !          za3(ji,1) =  1.0
952            !          za2(ji,2) = -2.0
953            !          za3(ji,2) =  1.0
954            !--------------------------------------------------------
955           
956            !  Compute Enhanced Mixing Coefficients (LMD94,eq D6)
957            ! ---------------------------------------------------------------
958           
959            ! Delta
960            zdelta  = ( zhbl(ji)  - zdept(ji,1) ) / ( zdept(ji,2) - zdept(ji,1) + epsln )
961            zdelta2 = zdelta * zdelta
962           
963            !  Mixing coefficients at first level above h (zdept(ji,1))
964            ! and at first interface in the pipe (zdepw(ji,2))
965           
966            ! At first T level above h (zdept(ji,1)) (always in the boundary layer)
967            zsig    = zdept(ji,1) / zhbl(ji)
968            ztemp   = zstabl * zsig  + ( 1.0 - zstabl ) * MIN( zsig , epsilon )
969            zehat   = vonk * ztemp * zhbl(ji) * zbuofdep
970            zeta    = zehat / ( zucube + epsln)
971            zwst    = vonk * zustar(ji,jj) / ( ABS( 1.0 + rconc1 * zeta ) + epsln)
972            zwm     = zstabl * zwst + ( 1.0 - zstabl ) * zwm
973            zws     = zstabl * zwst + ( 1.0 - zstabl ) * zws
974
975            zkm1m  = zhbl(ji) * zwm * zsig * ( 1.0 + zsig * ( za2m(ji) + zsig * za3m(ji) ) )
976            zkm1t  = zhbl(ji) * zws * zsig * ( 1.0 + zsig * ( za2t(ji) + zsig * za3t(ji) ) )
977#if defined key_zdfddm
978            zkm1s  = zhbl(ji) * zws * zsig * ( 1.0 + zsig * ( za2s(ji) + zsig * za3s(ji) ) )
979#endif                       
980            ! At first W level in the pipe (zdepw(ji,2)) (not always in the boundary layer ):
981            zsig    = MIN( zdepw(ji,2) / zhbl(ji) , 1.0 )
982            ztemp   = zstabl * zsig + ( 1.0 - zstabl ) * MIN( zsig , epsilon )
983            zehat   = vonk * ztemp * zhbl(ji) * zbuofdep
984            zeta    = zehat / ( zucube + epsln )
985            zwst    = vonk * zustar(ji,jj) / ( ABS( 1.0 + rconc1 * zeta ) + epsln)
986            zws     = zstabl * zws + ( 1.0 - zstabl ) * zws
987            zwm     = zstabl * zws + ( 1.0 - zstabl ) * zwm
988
989            zkmpm(ji) = zhbl(ji) * zwm * zsig * ( 1.0 + zsig * ( za2m(ji) + zsig * za3m(ji) ) )
990            zkmpt(ji) = zhbl(ji) * zws * zsig * ( 1.0 + zsig * ( za2t(ji) + zsig * za3t(ji) ) )
991#if defined key_zdfddm
992            zkmps(ji) = zhbl(ji) * zws * zsig * ( 1.0 + zsig * ( za2s(ji) + zsig * za3s(ji) ) )
993#endif 
994     
995            ! check if this point is in the boundary layer,else take interior viscosity/diffusivity:
996            zflag       = 0.5 + SIGN( 0.5, ( zhbl(ji) - zdepw(ji,2) ) )
997            zkmpm(ji) = zkmpm(ji) * zflag + ( 1.0 - zflag ) * zvisc(ji,2)
998            zkmpt(ji) = zkmpt(ji) * zflag + ( 1.0 - zflag ) * zdift(ji,2)
999#if defined key_zdfddm
1000            zkmps(ji) = zkmps(ji) * zflag + ( 1.0 - zflag ) * zdifs(ji,2)
1001#endif
1002
1003            ! Enhanced viscosity/diffusivity at zdepw(ji,2)
1004            ztemp     = ( 1.0 - 2.0 * zdelta + zdelta2 ) * zkm1m + zdelta2 * zkmpm(ji)
1005            zkmpm(ji) = ( 1.0 - zdelta ) * zvisc(ji,2) + zdelta * ztemp
1006            ztemp     = ( 1.0 - 2.0 * zdelta + zdelta2 ) * zkm1t + zdelta2 * zkmpt(ji)
1007            zkmpt(ji) = ( 1.0 - zdelta ) * zdift(ji,2) + zdelta * ztemp
1008#if defined key_zdfddm
1009            ztemp     = ( 1.0 - 2.0 * zdelta + zdelta2 ) * zkm1s + zdelta2 * zkmps(ji)
1010            zkmps(ji) = ( 1.0 - zdelta ) * zdifs(ji,2) + zdelta * ztemp
1011#endif           
1012
1013         END DO
1014         !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
1015         ! IV. Compute vertical eddy viscosity and diffusivity coefficients
1016         !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1017         
1018         DO jk  = 2, jkmax
1019           
1020            ! Compute turbulent velocity scales on the interfaces
1021            ! --------------------------------------------------------
1022            DO  ji = fs_2, fs_jpim1
1023               zbuofdep = zBo(ji,jj) + zBosol(ji,jj) * zatt1
1024               zstabl   = 0.5 + SIGN( 0.5 , zbuofdep ) 
1025               zbuofdep = zbuofdep + zstabl * epsln         
1026               zsig    = fsdepw(ji,jj,jk) / zhbl(ji)
1027               ztemp   = zstabl * zsig + ( 1. - zstabl ) * MIN( zsig , epsilon )
1028               zehat   = vonk * ztemp * zhbl(ji) * zbuofdep
1029               zucube  = zustar(ji,jj) * zustar(ji,jj) * zustar(ji,jj)
1030               zeta    = zehat / ( zucube + epsln )
1031
1032               IF( zehat > 0. ) THEN
1033                  ! Stable case
1034                  zws  = vonk * zustar(ji,jj) / ( 1.0 + rconc1 * zeta )
1035                  zwm  = zws
1036               ELSE
1037                  ! Unstable case
1038#if defined key_kpplktb
1039                  ! use lookup table
1040                  zd     = zehat - dehatmin
1041                  il     = INT( zd / dezehat )
1042                  il     = MIN( il, nilktbm1 )
1043                  il     = MAX( il, 1 )
1044                 
1045                  ud     = zustar(ji,jj) - ustmin
1046                  jl     = INT( ud / deustar )
1047                  jl     = MIN( jl, njlktbm1 )
1048                  jl     = MAX( jl, 1 )
1049                 
1050                  zfrac  = zd / dezehat - FLOAT( il ) 
1051                  ufrac  = ud / deustar - FLOAT( jl )
1052                  zwas   = ( 1. - zfrac ) * wslktb(il,jl+1) + zfrac * wslktb(il+1,jl+1)
1053                  zwbs   = ( 1. - zfrac ) * wslktb(il,jl  ) + zfrac * wslktb(il+1,jl  )
1054                  zwam   = ( 1. - zfrac ) * wmlktb(il,jl+1) + zfrac * wmlktb(il+1,jl+1)
1055                  zwbm   = ( 1. - zfrac ) * wmlktb(il,jl  ) + zfrac * wmlktb(il+1,jl  )
1056                  !
1057                  zws    = ( 1. - ufrac ) * zwbs + ufrac * zwas
1058                  zwm    = ( 1. - ufrac ) * zwbm + ufrac * zwam
1059#else
1060                  ! use analytical functions
1061                  zconm  = 0.5 + SIGN( 0.5, ( rzetam - zeta) )
1062                  zcons  = 0.5 + SIGN( 0.5, ( rzetas - zeta) )
1063                 
1064                  ! Momentum : zeta < rzetam (zconm = 1)
1065                  ! Scalars  : zeta < rzetas (zcons = 1)
1066                  zwconm = zustar(ji,jj) * vonk * ( ( ABS( rconam - rconcm * zeta) )**pthird )
1067                  zwcons = zustar(ji,jj) * vonk * ( ( ABS( rconas - rconcs * zeta) )**pthird )
1068                 
1069                  ! Momentum : rzetam <= zeta < 0 (zconm = 0)
1070                  ! Scalars  : rzetas <= zeta < 0 (zcons = 0) 
1071                  zwmun  = SQRT( ABS( 1.0 - rconc2 * zeta ) )
1072                  zwsun  = vonk * zustar(ji,jj) * zwmun
1073                  zwmun  = vonk * zustar(ji,jj) * SQRT(zwmun)
1074                  !
1075                  zwm    = zconm * zwconm + ( 1.0 - zconm ) * zwmun
1076                  zws    = zcons * zwcons + ( 1.0 - zcons ) * zwsun
1077                 
1078#endif
1079               ENDIF
1080               
1081               zblcm(ji,jk) = zhbl(ji) * zwm * zsig  * ( 1.0 + zsig * ( za2m(ji) + zsig * za3m(ji) ) )
1082               zblct(ji,jk) = zhbl(ji) * zws * zsig  * ( 1.0 + zsig * ( za2t(ji) + zsig * za3t(ji) ) )
1083#if defined key_zdfddm
1084               zblcs(ji,jk) = zhbl(ji) * zws * zsig  * ( 1.0 + zsig * ( za2s(ji) + zsig * za3s(ji) ) )
1085#endif             
1086               !  Compute Nonlocal transport term = ghats * <ws>o
1087               ! ----------------------------------------------------
1088               ghats(ji,jj,jk-1) = ( 1. - zstabl ) * rcg / ( zws * zhbl(ji) + epsln ) * tmask(ji,jj,jk)
1089
1090            END DO
1091         END DO     
1092         !  Combine interior and boundary layer coefficients and nonlocal term
1093         ! -----------------------------------------------------------------------
1094         DO jk = 2, jpkm1   
1095            DO ji = fs_2, fs_jpim1
1096               zflag = zmask(ji,jk) * zmask(ji,jk+1)
1097               zviscos(ji,jj,jk) = ( 1.0 - zmask(ji,jk) )         * avmu (ji,jj,jk) & ! interior viscosities
1098                  &              +                        zflag   * zblcm(ji,jk    ) & ! boundary layer viscosities
1099                  &              + zmask(ji,jk) * ( 1.0 - zflag ) * zkmpm(ji       )   ! viscosity enhancement at W_level near zhbl
1100               
1101               zviscos(ji,jj,jk) = zviscos(ji,jj,jk) * tmask(ji,jj,jk)   
1102
1103           
1104               zdiffut(ji,jj,jk) = ( 1.0 - zmask(ji,jk) )          * avt (ji,jj,jk) & ! interior diffusivities
1105                  &              +                        zflag   * zblct(ji,jk   ) & ! boundary layer diffusivities
1106                  &              + zmask(ji,jk) * ( 1.0 - zflag ) * zkmpt(ji      )   ! diffusivity enhancement at W_level near zhbl
1107                       
1108               zdiffut(ji,jj,jk) = zdiffut(ji,jj,jk) * tmask(ji,jj,jk) 
1109#if defined key_zdfddm
1110               zdiffus(ji,jj,jk) = ( 1.0 - zmask(ji,jk) )          * avs (ji,jj,jk) & ! interior diffusivities
1111                  &              +                        zflag   * zblcs(ji,jk   ) & ! boundary layer diffusivities
1112                  &              + zmask(ji,jk) * ( 1.0 - zflag ) * zkmps(ji      )   ! diffusivity enhancement at W_level near zhbl
1113                       
1114               zdiffus(ji,jj,jk) = zdiffus(ji,jj,jk) * tmask(ji,jj,jk) 
1115#endif               
1116               ! Non local flux in the boundary layer only
1117               ghats(ji,jj,jk-1) = zmask(ji,jk) * ghats(ji,jj,jk-1)
1118
1119            ENDDO
1120         END DO
1121         !                                             ! ===============
1122      END DO                                           !   End of slab
1123      !                                                ! ===============
1124
1125      ! Lateral boundary conditions on zvicos and zdiffus  (sign unchanged)
1126      CALL lbc_lnk( zviscos(:,:,:), 'U', 1. )  ; CALL lbc_lnk( zdiffut(:,:,:), 'W', 1. ) 
1127#if defined key_zdfddm 
1128      CALL lbc_lnk( zdiffus(:,:,:), 'W', 1. ) 
1129#endif
1130
1131      SELECT CASE ( nn_ave )
1132         !
1133         CASE ( 0 )             ! no viscosity and diffusivity smoothing
1134
1135            DO jk = 2, jpkm1
1136               DO jj = 2, jpjm1
1137                  DO ji = fs_2, fs_jpim1
1138                     avmu(ji,jj,jk) = ( zviscos(ji,jj,jk) + zviscos(ji+1,jj,jk) ) &
1139                        &  / MAX( 1., tmask(ji,jj,jk) + tmask (ji + 1,jj,jk) ) * umask(ji,jj,jk)
1140                     
1141                     avmv(ji,jj,jk) = ( zviscos(ji,jj,jk) + zviscos(ji,jj+1,jk) ) &
1142                        &  / MAX( 1., tmask(ji,jj,jk) + tmask (ji,jj+1,jk) ) * vmask(ji,jj,jk)
1143                     
1144                     avt (ji,jj,jk) =  zdiffut(ji,jj,jk) * tmask(ji,jj,jk) 
1145#if defined key_zdfddm     
1146                     avs (ji,jj,jk) =  zdiffus(ji,jj,jk) * tmask(ji,jj,jk) 
1147#endif
1148                  END DO
1149               END DO
1150            END DO
1151           
1152         CASE ( 1 )                ! viscosity and diffusivity smoothing
1153            !                     
1154            !           ( 1/2  1  1/2 )              ( 1/2  1/2 )             ( 1/2  1  1/2 )
1155            ! avt = 1/8 ( 1    2  1   )   avmu = 1/4 ( 1    1   )   avmv= 1/4 ( 1/2  1  1/2 )
1156            !           ( 1/2  1  1/2 )              ( 1/2  1/2 )
1157 
1158            DO jk = 2, jpkm1
1159               DO jj = 2, jpjm1
1160                  DO ji = fs_2, fs_jpim1
1161
1162                     avmu(ji,jj,jk) = (      zviscos(ji  ,jj  ,jk) + zviscos(ji+1,jj  ,jk)   &
1163                        &              +.5*( zviscos(ji  ,jj-1,jk) + zviscos(ji+1,jj-1,jk)   &
1164                        &                   +zviscos(ji  ,jj+1,jk) + zviscos(ji+1,jj+1,jk) ) ) * eumean(ji,jj,jk)
1165                     
1166                     avmv(ji,jj,jk) = (      zviscos(ji  ,jj  ,jk) + zviscos(ji  ,jj+1,jk)   &
1167                        &              +.5*( zviscos(ji-1,jj  ,jk) + zviscos(ji-1,jj+1,jk)   &
1168                        &                   +zviscos(ji+1,jj  ,jk) + zviscos(ji+1,jj+1,jk) ) ) * evmean(ji,jj,jk)
1169 
1170                     avt (ji,jj,jk) = ( .5*( zdiffut(ji-1,jj+1,jk) + zdiffut(ji-1,jj-1,jk)    &
1171                        &                   +zdiffut(ji+1,jj+1,jk) + zdiffut(ji+1,jj-1,jk) )  &
1172                        &              +1.*( zdiffut(ji-1,jj  ,jk) + zdiffut(ji  ,jj+1,jk)    &
1173                        &                   +zdiffut(ji  ,jj-1,jk) + zdiffut(ji+1,jj  ,jk) )  &
1174                        &              +2.*  zdiffut(ji  ,jj  ,jk)                          ) * etmean(ji,jj,jk)
1175#if defined key_zdfddm   
1176                     avs (ji,jj,jk) = ( .5*( zdiffus(ji-1,jj+1,jk) + zdiffus(ji-1,jj-1,jk)    &
1177                        &                   +zdiffus(ji+1,jj+1,jk) + zdiffus(ji+1,jj-1,jk) )  &
1178                        &              +1.*( zdiffus(ji-1,jj  ,jk) + zdiffus(ji  ,jj+1,jk)    &
1179                        &                   +zdiffus(ji  ,jj-1,jk) + zdiffus(ji+1,jj  ,jk) )  &
1180                        &              +2.*  zdiffus(ji  ,jj  ,jk)                          ) * etmean(ji,jj,jk) 
1181#endif               
1182                  END DO
1183               END DO
1184            END DO
1185         
1186         END SELECT
1187
1188         DO jk = 2, jpkm1                       ! vertical slab
1189            !
1190            !  Minimum value on the eddy diffusivity
1191            ! ----------------------------------------
1192            DO jj = 2, jpjm1
1193               DO ji = fs_2, fs_jpim1   ! vector opt.
1194                  avt(ji,jj,jk) = MAX( avt(ji,jj,jk), avtb(jk) ) * tmask(ji,jj,jk)
1195#if defined key_zdfddm 
1196                  avs(ji,jj,jk) = MAX( avs(ji,jj,jk), avtb(jk) ) * tmask(ji,jj,jk)
1197#endif
1198               END DO
1199            END DO
1200
1201            !
1202            ! Minimum value on the eddy viscosity
1203            ! ----------------------------------------
1204            DO jj = 1, jpj
1205               DO ji = 1, jpi
1206                  avmu(ji,jj,jk) = MAX( avmu(ji,jj,jk), avmb(jk) ) * umask(ji,jj,jk)
1207                  avmv(ji,jj,jk) = MAX( avmv(ji,jj,jk), avmb(jk) ) * vmask(ji,jj,jk)
1208               END DO
1209            END DO
1210            !
1211         END DO
1212
1213         ! Lateral boundary conditions on avt  (sign unchanged)
1214         CALL lbc_lnk( hkpp(:,:), 'T', 1. )
1215
1216         ! Lateral boundary conditions on avt  (sign unchanged)
1217         CALL lbc_lnk( avt(:,:,:), 'W', 1. )
1218#if defined key_zdfddm 
1219         CALL lbc_lnk( avs(:,:,:), 'W', 1. ) 
1220#endif
1221         ! Lateral boundary conditions (avmu,avmv) (U- and V- points, sign unchanged)
1222         CALL lbc_lnk( avmu(:,:,:), 'U', 1. )   ;    CALL lbc_lnk( avmv(:,:,:), 'V', 1. ) 
1223 
1224         IF(ln_ctl) THEN
1225#if defined key_zdfddm
1226            CALL prt_ctl(tab3d_1=avt , clinfo1=' kpp - t: ', tab3d_2=avs , clinfo2=' s: ', ovlap=1, kdim=jpk)
1227#else
1228            CALL prt_ctl(tab3d_1=avt , clinfo1=' kpp - t: ', ovlap=1, kdim=jpk)
1229#endif
1230            CALL prt_ctl(tab3d_1=avmu, clinfo1=' kpp - u: ', mask1=umask,  &
1231               &         tab3d_2=avmv, clinfo2=      ' v: ', mask2=vmask, ovlap=1, kdim=jpk)
1232         ENDIF
1233
1234   END SUBROUTINE zdf_kpp
1235
1236
1237   SUBROUTINE tra_kpp( kt )
1238      !!----------------------------------------------------------------------
1239      !!                  ***  ROUTINE tra_kpp  ***
1240      !!
1241      !! ** Purpose :   compute and add to the tracer trend the non-local
1242      !!      tracer flux
1243      !!
1244      !! ** Method  :   ???
1245      !!
1246      !! history :
1247      !!     9.0  ! 05-11 (G. Madec)  Original code
1248      !!----------------------------------------------------------------------
1249      !! * Modules used
1250      USE oce, ONLY :    ztrdt => ua,       & ! use ua as 3D workspace
1251                         ztrds => va          ! use va as 3D workspace
1252      !!----------------------------------------------------------------------
1253      INTEGER, INTENT(in) :: kt
1254      INTEGER :: ji, jj, jk
1255
1256      IF( kt == nit000 ) THEN
1257         IF(lwp) WRITE(numout,*)
1258         IF(lwp) WRITE(numout,*) 'tra_kpp : KPP non-local tracer fluxes'
1259         IF(lwp) WRITE(numout,*) '~~~~~~~   '
1260      ENDIF
1261
1262
1263      ! Save ta and sa trends
1264      IF( l_trdtra )   THEN
1265         ztrdt(:,:,:) = ta(:,:,:)
1266         ztrds(:,:,:) = sa(:,:,:)
1267      ENDIF
1268
1269      ! add non-local temperature and salinity flux ( in convective case only)
1270      DO jk = 1, jpkm1
1271         DO jj = 2, jpjm1
1272            DO ji = fs_2, fs_jpim1
1273               ta(ji,jj,jk) =  ta(ji,jj,jk)                           &
1274                  &         - ( ghats(ji,jj,jk  ) * avt(ji,jj,jk  )   &
1275                  &           - ghats(ji,jj,jk+1) * avt(ji,jj,jk+1) ) * wt0(ji,jj) / fse3t(ji,jj,jk)
1276               sa(ji,jj,jk) = sa(ji,jj,jk)                              &
1277                  &         - ( ghats(ji,jj,jk  ) * fsavs(ji,jj,jk  )   &
1278                  &           - ghats(ji,jj,jk+1) * fsavs(ji,jj,jk+1) ) * ws0(ji,jj) / fse3t(ji,jj,jk)
1279            END DO
1280         END DO
1281      END DO
1282
1283      ! save the non-local tracer flux trends for diagnostic
1284      IF( l_trdtra )   THEN
1285         ztrdt(:,:,:) = ta(:,:,:) - ztrdt(:,:,:)
1286         ztrds(:,:,:) = sa(:,:,:) - ztrds(:,:,:)
1287!!bug gm jpttdzdf ==> jpttkpp
1288         CALL trd_mod(ztrdt, ztrds, jptra_trd_zdf, 'TRA', kt)
1289      ENDIF
1290
1291      IF(ln_ctl) THEN 
1292         CALL prt_ctl( tab3d_1=ta, clinfo1=' kpp  - Ta: ', mask1=tmask,   &
1293         &             tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
1294      ENDIF
1295
1296   END SUBROUTINE tra_kpp
1297
1298
1299   SUBROUTINE zdf_kpp_init
1300      !!----------------------------------------------------------------------
1301      !!                  ***  ROUTINE zdf_kpp_init  ***
1302      !!                     
1303      !! ** Purpose :   Initialization of the vertical eddy diffivity and
1304      !!      viscosity when using a kpp turbulent closure scheme
1305      !!
1306      !! ** Method  :   Read the namkpp namelist and check the parameters
1307      !!      called at the first timestep (nit000)
1308      !!
1309      !! ** input   :   Namlist namkpp
1310      !!
1311      !!
1312      !! history :
1313      !!     8.1  ! 00-02 (J. Chanut) KPP Mixing
1314      !!     9.0  ! 05-01 (C. Ethe) F90 : free form
1315      !!----------------------------------------------------------------------
1316      !! * local declarations
1317
1318      INTEGER    ::   &
1319         ji, jj, jk             ! dummy loop indices
1320     
1321#if ! defined key_kppcustom
1322      INTEGER    ::   &
1323         jm                       ! dummy loop indices     
1324      REAL(wp)   ::              & !!! tempory scalars
1325         zref, zdist
1326#endif
1327
1328#if defined key_kpplktb
1329      REAL(wp)   ::              & !!! tempory scalars
1330         zustar,    &
1331         zucube, zustvk,         & 
1332         zeta, zehat
1333#endif
1334      REAL(wp)   ::             & !!! tempory scalars
1335         zhbf
1336      LOGICAL ::                &
1337         ll_kppcustom,          &  ! 1st ocean level taken as surface layer
1338         ll_kpplktb                ! Lookup table for turbul. velocity scales
1339      !!
1340      NAMELIST/namzdf_kpp/ ln_kpprimix, rn_difmiw, rn_difsiw, rn_riinfty, rn_difri, rn_bvsqcon, rn_difcon, nn_ave
1341      !!----------------------------------------------------------------------
1342
1343      REWIND ( numnam )               ! Read Namelist namkpp : K-Profile Parameterisation
1344      READ   ( numnam, namzdf_kpp )
1345
1346      IF(lwp) THEN                    ! Control print
1347         WRITE(numout,*)
1348         WRITE(numout,*) 'zdf_kpp_init : K-Profile Parameterisation'
1349         WRITE(numout,*) '~~~~~~~~~~~~'
1350         WRITE(numout,*) '   Namelist namzdf_kpp : set tke mixing parameters'
1351         WRITE(numout,*) '     Shear instability mixing                      ln_kpprimix = ', ln_kpprimix
1352         WRITE(numout,*) '     max. internal wave viscosity                  rn_difmiw   = ', rn_difmiw
1353         WRITE(numout,*) '     max. internal wave diffusivity                rn_difsiw   = ', rn_difsiw
1354         WRITE(numout,*) '     Richardson Number limit for shear instability rn_riinfty  = ', rn_riinfty
1355         WRITE(numout,*) '     max. shear mixing at Rig = 0                  rn_difri    = ', rn_difri
1356         WRITE(numout,*) '     Brunt-Vaisala squared for max. convection     rn_bvsqcon  = ', rn_bvsqcon
1357         WRITE(numout,*) '     max. mix. in interior convec.                 rn_difcon   = ', rn_difcon
1358         WRITE(numout,*) '     horizontal average flag                       nn_ave      = ', nn_ave
1359      ENDIF
1360
1361      ll_kppcustom = .FALSE.
1362      ll_kpplktb   = .FALSE.
1363
1364#if defined key_kppcustom
1365      ll_kppcustom = .TRUE.
1366#endif
1367#if defined key_kpplktb
1368      ll_kpplktb   = .TRUE.
1369#endif
1370      IF(lwp) THEN
1371         WRITE(numout,*) '     Lookup table for turbul. velocity scales ll_kpplktb   = ', ll_kpplktb
1372         WRITE(numout,*) '     1st ocean level taken as surface layer   ll_kppcustom = ', ll_kppcustom
1373      ENDIF
1374     
1375      IF( lk_zdfddm) THEN
1376         IF(lwp) THEN
1377            WRITE(numout,*)
1378            WRITE(numout,*) '    Double diffusion mixing on temperature and salinity '
1379            WRITE(numout,*) '    CAUTION : done in routine zdfkpp, not in routine zdfddm '
1380         ENDIF
1381      ENDIF
1382     
1383         
1384
1385      !set constants not in namelist
1386      !-----------------------------
1387      Vtc  = rconcv * SQRT( 0.2 / ( rconcs * epsilon ) ) / ( vonk * vonk * Ricr )
1388      rcg  = rcstar * vonk * ( rconcs * vonk * epsilon )**pthird
1389
1390      IF(lwp) THEN
1391         WRITE(numout,*)
1392         WRITE(numout,*) '     Constant value for unreso. turbul. velocity shear Vtc = ', Vtc
1393         WRITE(numout,*) '     Non-dimensional coef. for nonlocal transport      rcg = ', rcg
1394       ENDIF
1395
1396      ! ratt is the attenuation coefficient for solar flux
1397      ! Should be different is s_coordinate
1398      DO jk = 1, jpk
1399         zhbf     = - fsdept(1,1,jk) * hbf
1400         ratt(jk) = 1.0 - ( rabs * EXP( zhbf / xsi1 ) + ( 1.0 - rabs ) * EXP( zhbf / xsi2 ) )       
1401      ENDDO
1402
1403      ! Horizontal average : initialization of weighting arrays
1404      ! -------------------
1405     
1406      SELECT CASE ( nn_ave )
1407
1408      CASE ( 0 )                ! no horizontal average
1409         IF(lwp) WRITE(numout,*) '          no horizontal average on avt, avmu, avmv'
1410         IF(lwp) WRITE(numout,*) '          only in very high horizontal resolution !'
1411         ! weighting mean arrays etmean, eumean and evmean
1412         !           ( 1  1 )                                          ( 1 )
1413         ! avt = 1/4 ( 1  1 )     avmu = 1/2 ( 1  1 )       avmv=  1/2 ( 1 )
1414         !                         
1415         etmean(:,:,:) = 0.e0
1416         eumean(:,:,:) = 0.e0
1417         evmean(:,:,:) = 0.e0
1418         
1419         DO jk = 1, jpkm1
1420            DO jj = 2, jpjm1
1421               DO ji = 2, jpim1   ! vector opt.
1422                  etmean(ji,jj,jk) = tmask(ji,jj,jk)                     &
1423                  &  / MAX( 1.,  umask(ji-1,jj  ,jk) + umask(ji,jj,jk)   &
1424                  &            + vmask(ji  ,jj-1,jk) + vmask(ji,jj,jk)  )
1425                 
1426                  eumean(ji,jj,jk) = umask(ji,jj,jk)                     &
1427                  &  / MAX( 1.,  tmask(ji,jj,jk) + tmask(ji+1,jj  ,jk)  )
1428
1429                  evmean(ji,jj,jk) = vmask(ji,jj,jk)                     &
1430                  &  / MAX( 1.,  tmask(ji,jj,jk) + tmask(ji  ,jj+1,jk)  )
1431               END DO
1432            END DO
1433         END DO
1434
1435      CASE ( 1 )                ! horizontal average
1436         IF(lwp) WRITE(numout,*) '          horizontal average on avt, avmu, avmv'
1437         ! weighting mean arrays etmean, eumean and evmean
1438         !           ( 1/2  1  1/2 )              ( 1/2  1/2 )             ( 1/2  1  1/2 )
1439         ! avt = 1/8 ( 1    2  1   )   avmu = 1/4 ( 1    1   )   avmv= 1/4 ( 1/2  1  1/2 )
1440         !           ( 1/2  1  1/2 )              ( 1/2  1/2 )
1441         etmean(:,:,:) = 0.e0
1442         eumean(:,:,:) = 0.e0
1443         evmean(:,:,:) = 0.e0
1444         
1445         DO jk = 1, jpkm1
1446            DO jj = 2, jpjm1
1447               DO ji = fs_2, fs_jpim1   ! vector opt.
1448                  etmean(ji,jj,jk) = tmask(ji, jj,jk)                           &
1449                     & / MAX( 1., 2.* tmask(ji,jj,jk)                           &
1450                     &      +.5 * ( tmask(ji-1,jj+1,jk) + tmask(ji-1,jj-1,jk)   &
1451                     &             +tmask(ji+1,jj+1,jk) + tmask(ji+1,jj-1,jk) ) &
1452                     &      +1. * ( tmask(ji-1,jj  ,jk) + tmask(ji  ,jj+1,jk)   &
1453                     &             +tmask(ji  ,jj-1,jk) + tmask(ji+1,jj  ,jk) ) )
1454                 
1455                  eumean(ji,jj,jk) = umask(ji,jj,jk)                        &
1456                     &  / MAX( 1.,   tmask(ji,jj  ,jk) + tmask(ji+1,jj  ,jk)   &
1457                     &       +.5 * ( tmask(ji,jj-1,jk) + tmask(ji+1,jj-1,jk)   &
1458                     &              +tmask(ji,jj+1,jk) + tmask(ji+1,jj+1,jk) )  )
1459                 
1460                  evmean(ji,jj,jk) = vmask(ji,jj,jk)                        &
1461                     &  / MAX( 1.,   tmask(ji  ,jj,jk) + tmask(ji  ,jj+1,jk)   &
1462                     &       +.5 * ( tmask(ji-1,jj,jk) + tmask(ji-1,jj+1,jk)   &
1463                     &              +tmask(ji+1,jj,jk) + tmask(ji+1,jj+1,jk) )  )
1464               END DO
1465            END DO
1466         END DO
1467
1468      CASE DEFAULT
1469         WRITE(ctmp1,*) '          bad flag value for nn_ave = ', nn_ave
1470         CALL ctl_stop( ctmp1 )
1471
1472      END SELECT
1473 
1474      ! Initialization of vertical eddy coef. to the background value
1475      ! -------------------------------------------------------------
1476      DO jk = 1, jpk
1477         avt (:,:,jk) = avtb(jk) * tmask(:,:,jk)
1478         avmu(:,:,jk) = avmb(jk) * umask(:,:,jk)
1479         avmv(:,:,jk) = avmb(jk) * vmask(:,:,jk)
1480      END DO
1481
1482      ! zero the surface flux for non local term and kpp mixed layer depth
1483      ! ------------------------------------------------------------------
1484      ghats(:,:,:) = 0.
1485      wt0  (:,:  ) = 0.
1486      ws0  (:,:  ) = 0.
1487      hkpp (:,:  ) = 0. ! just a diagnostic (not essential)
1488
1489#if ! defined key_kppcustom
1490      ! compute arrays (del, wz) for reference mean values
1491      ! (increase speed for vectorization key_kppcustom not defined)
1492      del(1:jpk, 1:jpk) = 0.
1493      DO jk = 1, jpk
1494         zref = epsilon * fsdept(1,1,jk)   
1495         DO jm = 1 , jpk
1496            zdist = zref - fsdepw(1,1,jm)   
1497            IF( zdist > 0.  ) THEN
1498               del(jk,jm) = MIN( zdist, fse3t(1,1,jm) ) / zref   
1499            ELSE
1500               del(jk,jm) = 0.
1501            ENDIF
1502         ENDDO
1503      ENDDO
1504#endif
1505
1506#if defined key_kpplktb
1507      ! build lookup table for turbulent velocity scales
1508      dezehat = ( dehatmax - dehatmin ) / nilktbm1
1509      deustar = ( ustmax   - ustmin   ) / njlktbm1
1510 
1511      DO jj = 1, njlktb
1512         zustar = ( jj - 1) * deustar + ustmin
1513         zustvk = vonk * zustar 
1514         zucube = zustar * zustar * zustar 
1515         DO ji = 1 , nilktb
1516            zehat = ( ji - 1 ) * dezehat + dehatmin
1517            zeta   = zehat / ( zucube + epsln )
1518            IF( zehat >= 0 ) THEN             ! Stable case
1519               wmlktb(ji,jj) = zustvk / ABS( 1.0 + rconc1 * zeta + epsln )                       
1520               wslktb(ji,jj) = wmlktb(ji,jj)
1521            ELSE                                ! Unstable case
1522               IF( zeta > rzetam ) THEN
1523                  wmlktb(ji,jj) = zustvk * ABS( 1.0    - rconc2 * zeta )**pfourth
1524               ELSE
1525                  wmlktb(ji,jj) = zustvk * ABS( rconam - rconcm * zeta )**pthird
1526               ENDIF
1527               
1528               IF( zeta > rzetas ) THEN
1529                  wslktb(ji,jj) = zustvk * SQRT( ABS( 1.0 - rconc2 * zeta ) )
1530               ELSE
1531                  wslktb(ji,jj) = zustvk * ABS( rconas - rconcs * zeta )**pthird
1532               ENDIF
1533            ENDIF
1534         END DO
1535      END DO
1536#endif
1537   END SUBROUTINE zdf_kpp_init
1538
1539#else
1540   !!----------------------------------------------------------------------
1541   !!   Dummy module :                                        NO KPP scheme
1542   !!----------------------------------------------------------------------
1543   LOGICAL, PUBLIC, PARAMETER ::   lk_zdfkpp = .FALSE.   !: KPP flag
1544CONTAINS
1545   SUBROUTINE zdf_kpp( kt )          ! Empty routine
1546      WRITE(*,*) 'zdf_kpp: You should not have seen this print! error?', kt
1547   END SUBROUTINE zdf_kpp
1548   SUBROUTINE tra_kpp( kt )          ! Empty routine
1549      WRITE(*,*) 'tra_kpp: You should not have seen this print! error?', kt
1550   END SUBROUTINE tra_kpp
1551#endif
1552
1553   !!======================================================================
1554END MODULE zdfkpp
Note: See TracBrowser for help on using the repository browser.