New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
zdfkpp.F90 in branches/DEV_r2106_LOCEAN2010/NEMO/OPA_SRC/ZDF – NEMO

source: branches/DEV_r2106_LOCEAN2010/NEMO/OPA_SRC/ZDF/zdfkpp.F90 @ 2236

Last change on this file since 2236 was 2236, checked in by cetlod, 14 years ago

First guess of NEMO_v3.3

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