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

Last change on this file since 472 was 463, checked in by opalod, 18 years ago

nemo_v1_update_054:RB: take into account tracers and dynamics reorganization, change atsk to jki

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