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

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

source: trunk/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfkpp.F90 @ 3042

Last change on this file since 3042 was 2715, checked in by rblod, 13 years ago

First attempt to put dynamic allocation on the trunk

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