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

source: branches/dev_r2586_dynamic_mem/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfkpp.F90 @ 2636

Last change on this file since 2636 was 2636, checked in by gm, 13 years ago

dynamic mem: #785 ; move ctl_stop & warn in lib_mpp to avoid a circular dependency + ctl_stop improvment

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