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.
zdfosm.F90 in NEMO/branches/NERC/dev_r11078_OSMOSIS_IMMERSE_Nurser_4.0/src/OCE/ZDF – NEMO

source: NEMO/branches/NERC/dev_r11078_OSMOSIS_IMMERSE_Nurser_4.0/src/OCE/ZDF/zdfosm.F90 @ 14515

Last change on this file since 14515 was 14515, checked in by agn, 3 years ago

slow down shoaling for both stabilizing surface forcing and FK

  • Property svn:keywords set to Id
File size: 162.2 KB
Line 
1MODULE zdfosm
2   !!======================================================================
3   !!                       ***  MODULE  zdfosm  ***
4   !! Ocean physics:  vertical mixing coefficient compute from the OSMOSIS
5   !!                 turbulent closure parameterization
6   !!=====================================================================
7   !!  History : NEMO 4.0  ! A. Grant, G. Nurser
8   !! 15/03/2017  Changed calculation of pycnocline thickness in unstable conditions and stable conditions AG
9   !! 15/03/2017  Calculation of pycnocline gradients for stable conditions changed. Pycnocline gradients now depend on stability of the OSBL. A.G
10   !! 06/06/2017  (1) Checks on sign of buoyancy jump in calculation of  OSBL depth.  A.G.
11   !!             (2) Removed variable zbrad0, zbradh and zbradav since they are not used.
12   !!             (3) Approximate treatment for shear turbulence.
13   !!                        Minimum values for zustar and zustke.
14   !!                        Add velocity scale, zvstr, that tends to zustar for large Langmuir numbers.
15   !!                        Limit maximum value for Langmuir number.
16   !!                        Use zvstr in definition of stability parameter zhol.
17   !!             (4) Modified parametrization of entrainment flux, changing original coefficient 0.0485 for Langmuir contribution to 0.135 * zla
18   !!             (5) For stable boundary layer add factor that depends on length of timestep to 'slow' collapse and growth. Make sure buoyancy jump not negative.
19   !!             (6) For unstable conditions when growth is over multiple levels, limit change to maximum of one level per cycle through loop.
20   !!             (7) Change lower limits for loops that calculate OSBL averages from 1 to 2. Large gradients between levels 1 and 2 can cause problems.
21   !!             (8) Change upper limits from ibld-1 to ibld.
22   !!             (9) Calculation of pycnocline thickness in unstable conditions. Check added to ensure that buoyancy jump is positive before calculating Ri.
23   !!            (10) Thickness of interface layer at base of the stable OSBL set by Richardson number. Gives continuity in transition from unstable OSBL.
24   !!            (11) Checks that buoyancy jump is poitive when calculating pycnocline profiles.
25   !!            (12) Replace zwstrl with zvstr in calculation of eddy viscosity.
26   !! 27/09/2017 (13) Calculate Stokes drift and Stokes penetration depth from wave information
27   !!            (14) Buoyancy flux due to entrainment changed to include contribution from shear turbulence.
28   !! 28/09/2017 (15) Calculation of Stokes drift moved into separate do-loops to allow for different options for the determining the Stokes drift to be added.
29   !!            (16) Calculation of Stokes drift from windspeed for PM spectrum (for testing, commented out)
30   !!            (17) Modification to Langmuir velocity scale to include effects due to the Stokes penetration depth (for testing, commented out)
31   !! ??/??/2018 (18) Revision to code structure, selected using key_osmldpth1. Inline code moved into subroutines. Changes to physics made,
32   !!                  (a) Pycnocline temperature and salinity profies changed for unstable layers
33   !!                  (b) The stable OSBL depth parametrization changed.
34   !! 16/05/2019 (19) Fox-Kemper parametrization of restratification through mixed layer eddies added to revised code.
35   !! 23/05/19   (20) Old code where key_osmldpth1` is *not* set removed, together with the key key_osmldpth1
36   !!----------------------------------------------------------------------
37
38   !!----------------------------------------------------------------------
39   !!   'ln_zdfosm'                                             OSMOSIS scheme
40   !!----------------------------------------------------------------------
41   !!   zdf_osm       : update momentum and tracer Kz from osm scheme
42   !!   zdf_osm_init  : initialization, namelist read, and parameters control
43   !!   osm_rst       : read (or initialize) and write osmosis restart fields
44   !!   tra_osm       : compute and add to the T & S trend the non-local flux
45   !!   trc_osm       : compute and add to the passive tracer trend the non-local flux (TBD)
46   !!   dyn_osm       : compute and add to u & v trensd the non-local flux
47   !!
48   !! Subroutines in revised code.
49   !!----------------------------------------------------------------------
50   USE oce            ! ocean dynamics and active tracers
51                      ! uses wn from previous time step (which is now wb) to calculate hbl
52   USE dom_oce        ! ocean space and time domain
53   USE zdf_oce        ! ocean vertical physics
54   USE sbc_oce        ! surface boundary condition: ocean
55   USE sbcwave        ! surface wave parameters
56   USE phycst         ! physical constants
57   USE eosbn2         ! equation of state
58   USE traqsr         ! details of solar radiation absorption
59   USE zdfddm         ! double diffusion mixing (avs array)
60!   USE zdfmxl         ! mixed layer depth
61   USE iom            ! I/O library
62   USE lib_mpp        ! MPP library
63   USE trd_oce        ! ocean trends definition
64   USE trdtra         ! tracers trends
65   !
66   USE in_out_manager ! I/O manager
67   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
68   USE prtctl         ! Print control
69   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined)
70
71   IMPLICIT NONE
72   PRIVATE
73
74   PUBLIC   zdf_osm       ! routine called by step.F90
75   PUBLIC   zdf_osm_init  ! routine called by nemogcm.F90
76   PUBLIC   osm_rst       ! routine called by step.F90
77   PUBLIC   tra_osm       ! routine called by step.F90
78   PUBLIC   trc_osm       ! routine called by trcstp.F90
79   PUBLIC   dyn_osm       ! routine called by step.F90
80
81   PUBLIC   ln_osm_mle    ! logical needed by tra_mle_init in tramle.F90
82
83   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   ghamu    !: non-local u-momentum flux
84   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   ghamv    !: non-local v-momentum flux
85   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   ghamt    !: non-local temperature flux (gamma/<ws>o)
86   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   ghams    !: non-local salinity flux (gamma/<ws>o)
87   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   etmean   !: averaging operator for avt
88   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   hbl      !: boundary layer depth
89   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   dh       ! depth of pycnocline
90   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   hml      ! ML depth
91   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   dstokes  !: penetration depth of the Stokes drift.
92
93   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)           ::   r1_ft    ! inverse of the modified Coriolis parameter at t-pts
94   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   hmle     ! Depth of layer affexted by mixed layer eddies in Fox-Kemper parametrization
95   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   dbdx_mle ! zonal buoyancy gradient in ML
96   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   dbdy_mle ! meridional buoyancy gradient in ML
97   INTEGER,  PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   mld_prof ! level of base of MLE layer.
98
99   !                      !!** Namelist  namzdf_osm  **
100   LOGICAL  ::   ln_use_osm_la      ! Use namelist  rn_osm_la
101
102   LOGICAL  ::   ln_osm_mle           !: flag to activate the Mixed Layer Eddy (MLE) parameterisation
103
104   REAL(wp) ::   rn_osm_la          ! Turbulent Langmuir number
105   REAL(wp) ::   rn_osm_dstokes     ! Depth scale of Stokes drift
106   REAL(wp) ::   rn_zdfosm_adjust_sd = 1.0 ! factor to reduce Stokes drift by
107   REAL(wp) ::   rn_osm_hblfrac = 0.1! for nn_osm_wave = 3/4 specify fraction in top of hbl
108   LOGICAL  ::   ln_zdfosm_ice_shelter      ! flag to activate ice sheltering
109   REAL(wp) ::   rn_osm_hbl0 = 10._wp       ! Initial value of hbl for 1D runs
110   INTEGER  ::   nn_ave             ! = 0/1 flag for horizontal average on avt
111   INTEGER  ::   nn_osm_wave = 0    ! = 0/1/2 flag for getting stokes drift from La# / PM wind-waves/Inputs into sbcwave
112   INTEGER  ::   nn_osm_SD_reduce   ! = 0/1/2 flag for getting effective stokes drift from surface value
113   LOGICAL  ::   ln_dia_osm         ! Use namelist  rn_osm_la
114
115   LOGICAL  ::   ln_kpprimix  = .true.  ! Shear instability mixing
116   REAL(wp) ::   rn_riinfty   = 0.7     ! local Richardson Number limit for shear instability
117   REAL(wp) ::   rn_difri    =  0.005   ! maximum shear mixing at Rig = 0    (m2/s)
118   LOGICAL  ::   ln_convmix  = .true.   ! Convective instability mixing
119   REAL(wp) ::   rn_difconv = 1._wp     ! diffusivity when unstable below BL  (m2/s)
120
121! OSMOSIS mixed layer eddy parametrization constants
122   INTEGER  ::   nn_osm_mle             ! = 0/1 flag for horizontal average on avt
123   REAL(wp) ::   rn_osm_mle_ce           ! MLE coefficient
124   !                                        ! parameters used in nn_osm_mle = 0 case
125   REAL(wp) ::   rn_osm_mle_lf               ! typical scale of mixed layer front
126   REAL(wp) ::   rn_osm_mle_time             ! time scale for mixing momentum across the mixed layer
127   !                                        ! parameters used in nn_osm_mle = 1 case
128   REAL(wp) ::   rn_osm_mle_lat              ! reference latitude for a 5 km scale of ML front
129   LOGICAL  ::   ln_osm_hmle_limit           ! If true arbitrarily restrict hmle to rn_osm_hmle_limit*zmld
130   REAL(wp) ::   rn_osm_hmle_limit           ! If ln_osm_hmle_limit true arbitrarily restrict hmle to rn_osm_hmle_limit*zmld
131   REAL(wp) ::   rn_osm_mle_rho_c        ! Density criterion for definition of MLD used by FK
132   REAL(wp) ::   r5_21 = 5.e0 / 21.e0   ! factor used in mle streamfunction computation
133   REAL(wp) ::   rb_c                   ! ML buoyancy criteria = g rho_c /rau0 where rho_c is defined in zdfmld
134   REAL(wp) ::   rc_f                   ! MLE coefficient (= rn_ce / (5 km * fo) ) in nn_osm_mle=1 case
135   REAL(wp) ::   rn_osm_mle_thresh          ! Threshold buoyancy for deepening of MLE layer below OSBL base.
136   REAL(wp) ::   rn_osm_bl_thresh          ! Threshold buoyancy for deepening of OSBL base.
137   REAL(wp) ::   rn_osm_mle_tau             ! Adjustment timescale for MLE.
138
139
140   !                                    !!! ** General constants  **
141   REAL(wp) ::   epsln   = 1.0e-20_wp   ! a small positive number to ensure no div by zero
142   REAL(wp) ::   depth_tol = 1.0e-6_wp  ! a small-ish positive number to give a hbl slightly shallower than gdepw
143   REAL(wp) ::   pthird  = 1._wp/3._wp  ! 1/3
144   REAL(wp) ::   p2third = 2._wp/3._wp  ! 2/3
145
146   INTEGER :: idebug = 236
147   INTEGER :: jdebug = 228
148   !!----------------------------------------------------------------------
149   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
150   !! $Id: zdfosm.F90 12317 2020-01-14 12:40:47Z agn $
151   !! Software governed by the CeCILL license (see ./LICENSE)
152   !!----------------------------------------------------------------------
153CONTAINS
154
155   INTEGER FUNCTION zdf_osm_alloc()
156      !!----------------------------------------------------------------------
157      !!                 ***  FUNCTION zdf_osm_alloc  ***
158      !!----------------------------------------------------------------------
159     ALLOCATE( ghamu(jpi,jpj,jpk), ghamv(jpi,jpj,jpk), ghamt(jpi,jpj,jpk),ghams(jpi,jpj,jpk), &
160          &       hbl(jpi,jpj), dh(jpi,jpj), hml(jpi,jpj), dstokes(jpi, jpj), &
161          &       etmean(jpi,jpj,jpk), STAT= zdf_osm_alloc )
162
163     ALLOCATE(  hmle(jpi,jpj), r1_ft(jpi,jpj), dbdx_mle(jpi,jpj), dbdy_mle(jpi,jpj), &
164          &       mld_prof(jpi,jpj), STAT= zdf_osm_alloc )
165
166!     ALLOCATE( ghamu(jpi,jpj,jpk), ghamv(jpi,jpj,jpk), ghamt(jpi,jpj,jpk),ghams(jpi,jpj,jpk), &    ! would ths be better ?
167!          &       hbl(jpi,jpj), dh(jpi,jpj), hml(jpi,jpj), dstokes(jpi, jpj), &
168!          &       etmean(jpi,jpj,jpk), STAT= zdf_osm_alloc )
169!     IF( zdf_osm_alloc /= 0 )   CALL ctl_warn('zdf_osm_alloc: failed to allocate zdf_osm arrays')
170!     IF ( ln_osm_mle ) THEN
171!        Allocate( hmle(jpi,jpj), r1_ft(jpi,jpj), STAT= zdf_osm_alloc )
172!        IF( zdf_osm_alloc /= 0 )   CALL ctl_warn('zdf_osm_alloc: failed to allocate zdf_osm mle arrays')
173!     ENDIF
174
175     IF( zdf_osm_alloc /= 0 )   CALL ctl_warn('zdf_osm_alloc: failed to allocate zdf_osm arrays')
176     CALL mpp_sum ( 'zdfosm', zdf_osm_alloc )
177   END FUNCTION zdf_osm_alloc
178
179
180   SUBROUTINE zdf_osm( kt, p_avm, p_avt )
181      !!----------------------------------------------------------------------
182      !!                   ***  ROUTINE zdf_osm  ***
183      !!
184      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
185      !!      coefficients and non local mixing using the OSMOSIS scheme
186      !!
187      !! ** Method :   The boundary layer depth hosm is diagnosed at tracer points
188      !!      from profiles of buoyancy, and shear, and the surface forcing.
189      !!      Above hbl (sigma=-z/hbl <1) the mixing coefficients are computed from
190      !!
191      !!                      Kx =  hosm  Wx(sigma) G(sigma)
192      !!
193      !!             and the non local term ghamt = Cs / Ws(sigma) / hosm
194      !!      Below hosm  the coefficients are the sum of mixing due to internal waves
195      !!      shear instability and double diffusion.
196      !!
197      !!      -1- Compute the now interior vertical mixing coefficients at all depths.
198      !!      -2- Diagnose the boundary layer depth.
199      !!      -3- Compute the now boundary layer vertical mixing coefficients.
200      !!      -4- Compute the now vertical eddy vicosity and diffusivity.
201      !!      -5- Smoothing
202      !!
203      !!        N.B. The computation is done from jk=2 to jpkm1
204      !!             Surface value of avt are set once a time to zero
205      !!             in routine zdf_osm_init.
206      !!
207      !! ** Action  :   update the non-local terms ghamts
208      !!                update avt (before vertical eddy coef.)
209      !!
210      !! References : Large W.G., Mc Williams J.C. and Doney S.C.
211      !!         Reviews of Geophysics, 32, 4, November 1994
212      !!         Comments in the code refer to this paper, particularly
213      !!         the equation number. (LMD94, here after)
214      !!----------------------------------------------------------------------
215      INTEGER                   , INTENT(in   ) ::   kt            ! ocean time step
216      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::  p_avm, p_avt   ! momentum and tracer Kz (w-points)
217      !!
218      INTEGER ::   ji, jj, jk                   ! dummy loop indices
219
220      INTEGER ::   jl                   ! dummy loop indices
221
222      INTEGER ::   ikbot, jkmax, jkm1, jkp2     !
223
224      REAL(wp) ::   ztx, zty, zflageos, zstabl, zbuofdep,zucube      !
225      REAL(wp) ::   zbeta, zthermal                                  !
226      REAL(wp) ::   zehat, zeta, zhrib, zsig, zscale, zwst, zws, zwm ! Velocity scales
227      REAL(wp) ::   zwsun, zwmun, zcons, zconm, zwcons, zwconm      !
228
229      REAL(wp) ::   zsr, zbw, ze, zb, zd, zc, zaw, za, zb1, za1, zkw, zk0, zcomp , zrhd,zrhdr,zbvzed   ! In situ density
230      INTEGER  ::   jm                          ! dummy loop indices
231      REAL(wp) ::   zr1, zr2, zr3, zr4, zrhop   ! Compression terms
232      REAL(wp) ::   zflag, zrn2, zdep21, zdep32, zdep43
233      REAL(wp) ::   zesh2, zri, zfri            ! Interior richardson mixing
234      REAL(wp) ::   zdelta, zdelta2, zdzup, zdzdn, zdzh, zvath, zgat1, zdat1, zkm1m, zkm1t
235      REAL(wp) :: zt,zs,zu,zv,zrh               ! variables used in constructing averages
236! Scales
237      REAL(wp), DIMENSION(jpi,jpj) :: zrad0     ! Surface solar temperature flux (deg m/s)
238      REAL(wp), DIMENSION(jpi,jpj) :: zradh     ! Radiative flux at bl base (Buoyancy units)
239      REAL(wp), DIMENSION(jpi,jpj) :: zradav    ! Radiative flux, bl average (Buoyancy Units)
240      REAL(wp), DIMENSION(jpi,jpj) :: zustar    ! friction velocity
241      REAL(wp), DIMENSION(jpi,jpj) :: zwstrl    ! Langmuir velocity scale
242      REAL(wp), DIMENSION(jpi,jpj) :: zvstr     ! Velocity scale that ends to zustar for large Langmuir number.
243      REAL(wp), DIMENSION(jpi,jpj) :: zwstrc    ! Convective velocity scale
244      REAL(wp), DIMENSION(jpi,jpj) :: zuw0      ! Surface u-momentum flux
245      REAL(wp), DIMENSION(jpi,jpj) :: zvw0      ! Surface v-momentum flux
246      REAL(wp), DIMENSION(jpi,jpj) :: zwth0     ! Surface heat flux (Kinematic)
247      REAL(wp), DIMENSION(jpi,jpj) :: zws0      ! Surface freshwater flux
248      REAL(wp), DIMENSION(jpi,jpj) :: zwb0      ! Surface buoyancy flux
249      REAL(wp), DIMENSION(jpi,jpj) :: zwb0tot   ! Total surface buoyancy flux including insolation
250      REAL(wp), DIMENSION(jpi,jpj) :: zwthav    ! Heat flux - bl average
251      REAL(wp), DIMENSION(jpi,jpj) :: zwsav     ! freshwater flux - bl average
252      REAL(wp), DIMENSION(jpi,jpj) :: zwbav     ! Buoyancy flux - bl average
253      REAL(wp), DIMENSION(jpi,jpj) :: zwb_ent   ! Buoyancy entrainment flux
254      REAL(wp), DIMENSION(jpi,jpj) :: zwb_min
255
256
257      REAL(wp), DIMENSION(jpi,jpj) :: zwb_fk_b  ! MLE buoyancy flux averaged over OSBL
258      REAL(wp), DIMENSION(jpi,jpj) :: zwb_fk    ! max MLE buoyancy flux
259      REAL(wp), DIMENSION(jpi,jpj) :: zdiff_mle ! extra MLE vertical diff
260      REAL(wp), DIMENSION(jpi,jpj) :: zvel_mle  ! velocity scale for dhdt with stable ML and FK
261
262      REAL(wp), DIMENSION(jpi,jpj) :: zustke    ! Surface Stokes drift
263      REAL(wp), DIMENSION(jpi,jpj) :: zla       ! Trubulent Langmuir number
264      REAL(wp), DIMENSION(jpi,jpj) :: zcos_wind ! Cos angle of surface stress
265      REAL(wp), DIMENSION(jpi,jpj) :: zsin_wind ! Sin angle of surface stress
266      REAL(wp), DIMENSION(jpi,jpj) :: zhol      ! Stability parameter for boundary layer
267      LOGICAL, DIMENSION(jpi,jpj)  :: lconv     ! unstable/stable bl
268      LOGICAL, DIMENSION(jpi,jpj)  :: lshear    ! Shear layers
269      LOGICAL, DIMENSION(jpi,jpj)  :: lpyc      ! OSBL pycnocline present
270      LOGICAL, DIMENSION(jpi,jpj)  :: lflux     ! surface flux extends below OSBL into MLE layer.
271      LOGICAL, DIMENSION(jpi,jpj)  :: lmle      ! MLE layer increases in hickness.
272
273      ! mixed-layer variables
274
275      INTEGER, DIMENSION(jpi,jpj) :: ibld ! level of boundary layer base
276      INTEGER, DIMENSION(jpi,jpj) :: imld ! level of mixed-layer depth (pycnocline top)
277      INTEGER, DIMENSION(jpi,jpj) :: jp_ext, jp_ext_mle ! offset for external level
278      INTEGER, DIMENSION(jpi, jpj) :: j_ddh ! Type of shear layer
279
280      REAL(wp) :: ztgrad,zsgrad,zbgrad ! Temporary variables used to calculate pycnocline gradients
281      REAL(wp) :: zugrad,zvgrad        ! temporary variables for calculating pycnocline shear
282
283      REAL(wp), DIMENSION(jpi,jpj) :: zhbl  ! bl depth - grid
284      REAL(wp), DIMENSION(jpi,jpj) :: zhml  ! ml depth - grid
285
286      REAL(wp), DIMENSION(jpi,jpj) :: zhmle ! MLE depth - grid
287      REAL(wp), DIMENSION(jpi,jpj) :: zmld  ! ML depth on grid
288
289      REAL(wp), DIMENSION(jpi,jpj) :: zdh   ! pycnocline depth - grid
290      REAL(wp), DIMENSION(jpi,jpj) :: zdhdt ! BL depth tendency
291      REAL(wp), DIMENSION(jpi,jpj) :: zdtdz_bl_ext,zdsdz_bl_ext,zdbdz_bl_ext              ! external temperature/salinity and buoyancy gradients
292      REAL(wp), DIMENSION(jpi,jpj) :: zdtdz_mle_ext,zdsdz_mle_ext,zdbdz_mle_ext              ! external temperature/salinity and buoyancy gradients
293      REAL(wp), DIMENSION(jpi,jpj) :: zdtdx, zdtdy, zdsdx, zdsdy      ! horizontal gradients for Fox-Kemper parametrization.
294
295      REAL(wp), DIMENSION(jpi,jpj) :: zt_bl,zs_bl,zu_bl,zv_bl,zb_bl  ! averages over the depth of the blayer
296      REAL(wp), DIMENSION(jpi,jpj) :: zt_ml,zs_ml,zu_ml,zv_ml,zb_ml  ! averages over the depth of the mixed layer
297      REAL(wp), DIMENSION(jpi,jpj) :: zt_mle,zs_mle,zu_mle,zv_mle,zb_mle  ! averages over the depth of the MLE layer
298      REAL(wp), DIMENSION(jpi,jpj) :: zdt_bl,zds_bl,zdu_bl,zdv_bl,zdb_bl ! difference between blayer average and parameter at base of blayer
299      REAL(wp), DIMENSION(jpi,jpj) :: zdt_ml,zds_ml,zdu_ml,zdv_ml,zdb_ml ! difference between mixed layer average and parameter at base of blayer
300      REAL(wp), DIMENSION(jpi,jpj) :: zdt_mle,zds_mle,zdu_mle,zdv_mle,zdb_mle ! difference between MLE layer average and parameter at base of blayer
301!      REAL(wp), DIMENSION(jpi,jpj) :: zwth_ent,zws_ent ! heat and salinity fluxes at the top of the pycnocline
302      REAL(wp) :: zwth_ent,zws_ent ! heat and salinity fluxes at the top of the pycnocline
303      REAL(wp) :: zuw_bse,zvw_bse  ! momentum fluxes at the top of the pycnocline
304      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zdtdz_pyc    ! parametrized gradient of temperature in pycnocline
305      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zdsdz_pyc    ! parametrised gradient of salinity in pycnocline
306      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zdbdz_pyc    ! parametrised gradient of buoyancy in the pycnocline
307      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zdudz_pyc    ! u-shear across the pycnocline
308      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zdvdz_pyc    ! v-shear across the pycnocline
309      REAL(wp), DIMENSION(jpi,jpj) :: zdbds_mle    ! Magnitude of horizontal buoyancy gradient.
310      ! Flux-gradient relationship variables
311      REAL(wp), DIMENSION(jpi, jpj) :: zshear ! Shear production.
312
313      REAL(wp) :: zl_c,zl_l,zl_eps  ! Used to calculate turbulence length scale.
314
315      REAL(wp) :: za_cubic, zb_cubic, zc_cubic, zd_cubic ! coefficients in cubic polynomial specifying diffusivity in pycnocline. 
316      REAL(wp), DIMENSION(jpi,jpj) :: zsc_wth_1,zsc_ws_1 ! Temporary scales used to calculate scalar non-gradient terms.
317      REAL(wp), DIMENSION(jpi,jpj) :: zsc_wth_pyc, zsc_ws_pyc ! Scales for pycnocline transport term/
318      REAL(wp), DIMENSION(jpi,jpj) :: zsc_uw_1,zsc_uw_2,zsc_vw_1,zsc_vw_2 ! Temporary scales for non-gradient momentum flux terms.
319      REAL(wp), DIMENSION(jpi,jpj) :: zhbl_t ! holds boundary layer depth updated by full timestep
320
321      ! For calculating Ri#-dependent mixing
322      REAL(wp), DIMENSION(jpi,jpj,jpk) :: z3du   ! u-shear^2
323      REAL(wp), DIMENSION(jpi,jpj,jpk) :: z3dv   ! v-shear^2
324      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zrimix ! spatial form of ri#-induced diffusion
325
326      ! Temporary variables
327      INTEGER :: inhml
328      REAL(wp) :: znd,znd_d,zznd_ml,zznd_pyc,zznd_d ! temporary non-dimensional depths used in various routines
329      REAL(wp) :: ztemp, zari, zpert, zzdhdt, zdb   ! temporary variables
330      REAL(wp) :: zthick, zz0, zz1 ! temporary variables
331      REAL(wp) :: zvel_max, zhbl_s ! temporary variables
332      REAL(wp) :: zfac, ztmp       ! temporary variable
333      REAL(wp) :: zus_x, zus_y     ! temporary Stokes drift
334      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zviscos ! viscosity
335      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zdiffut ! t-diffusivity
336      REAL(wp), DIMENSION(jpi,jpj) :: zalpha_pyc
337      REAL(wp), DIMENSION(jpi,jpj) :: ztau_sc_u ! dissipation timescale at baes of WML.
338      REAL(wp) :: zdelta_pyc, zwt_pyc_sc_1, zws_pyc_sc_1, zzeta_pyc
339      REAL(wp) :: zbuoy_pyc_sc, zomega, zvw_max
340      INTEGER :: ibld_ext=0                          ! does not have to be zero for modified scheme
341      REAL(wp) :: zgamma_b_nd, zgamma_b, zdhoh, ztau
342      REAL(wp) :: zzeta_s = 0._wp
343      REAL(wp) :: zzeta_v = 0.46
344      REAL(wp) :: zabsstke
345      REAL(wp) :: zsqrtpi, z_two_thirds, zproportion, ztransp, zthickness
346      REAL(wp) :: z2k_times_thickness, zsqrt_depth, zexp_depth, zdstokes0, zf, zexperfc
347
348      ! For debugging
349      INTEGER :: ikt
350      !!--------------------------------------------------------------------
351      !
352      ibld(:,:)   = 0     ; imld(:,:)  = 0
353      zrad0(:,:)  = 0._wp ; zradh(:,:) = 0._wp ; zradav(:,:)    = 0._wp ; zustar(:,:)    = 0._wp
354      zwstrl(:,:) = 0._wp ; zvstr(:,:) = 0._wp ; zwstrc(:,:)    = 0._wp ; zuw0(:,:)      = 0._wp
355      zvw0(:,:)   = 0._wp ; zwth0(:,:) = 0._wp ; zws0(:,:)      = 0._wp ; zwb0(:,:)      = 0._wp
356      zwthav(:,:) = 0._wp ; zwsav(:,:) = 0._wp ; zwbav(:,:)     = 0._wp ; zwb_ent(:,:)   = 0._wp
357      zustke(:,:) = 0._wp ; zla(:,:)   = 0._wp ; zcos_wind(:,:) = 0._wp ; zsin_wind(:,:) = 0._wp
358      zhol(:,:)   = 0._wp ; zwb0tot(:,:) = 0._wp
359      lconv(:,:)  = .FALSE.; lpyc(:,:) = .FALSE. ; lflux(:,:) = .FALSE. ;  lmle(:,:) = .FALSE.
360      ! mixed layer
361      ! no initialization of zhbl or zhml (or zdh?)
362      zhbl(:,:)    = 1._wp ; zhml(:,:)    = 1._wp ; zdh(:,:)      = 1._wp ; zdhdt(:,:)   = 0._wp
363      zt_bl(:,:)   = 0._wp ; zs_bl(:,:)   = 0._wp ; zu_bl(:,:)    = 0._wp
364      zv_bl(:,:)   = 0._wp ; zb_bl(:,:)  = 0._wp
365      zt_ml(:,:)   = 0._wp ; zs_ml(:,:)    = 0._wp ; zu_ml(:,:)   = 0._wp
366      zt_mle(:,:)   = 0._wp ; zs_mle(:,:)    = 0._wp ; zu_mle(:,:)   = 0._wp
367      zb_mle(:,:) = 0._wp
368      zv_ml(:,:)   = 0._wp ; zdt_bl(:,:)   = 0._wp ; zds_bl(:,:)  = 0._wp
369      zdu_bl(:,:)  = 0._wp ; zdv_bl(:,:)  = 0._wp ; zdb_bl(:,:)  = 0._wp
370      zdt_ml(:,:)  = 0._wp ; zds_ml(:,:)  = 0._wp ; zdu_ml(:,:)   = 0._wp ; zdv_ml(:,:)  = 0._wp
371      zdb_ml(:,:)  = 0._wp
372      zdt_mle(:,:)  = 0._wp ; zds_mle(:,:)  = 0._wp ; zdu_mle(:,:)   = 0._wp
373      zdv_mle(:,:)  = 0._wp ; zdb_mle(:,:)  = 0._wp
374      zwth_ent = 0._wp ; zws_ent = 0._wp
375      !
376      zdtdz_pyc(:,:,:) = 0._wp ; zdsdz_pyc(:,:,:) = 0._wp ; zdbdz_pyc(:,:,:) = 0._wp
377      zdudz_pyc(:,:,:) = 0._wp ; zdvdz_pyc(:,:,:) = 0._wp
378      !
379      zdtdz_bl_ext(:,:) = 0._wp ; zdsdz_bl_ext(:,:) = 0._wp ; zdbdz_bl_ext(:,:) = 0._wp
380
381      IF ( ln_osm_mle ) THEN  ! only initialise arrays if needed
382         zdtdx(:,:) = 0._wp ; zdtdy(:,:) = 0._wp ; zdsdx(:,:) = 0._wp
383         zdsdy(:,:) = 0._wp ; dbdx_mle(:,:) = 0._wp ; dbdy_mle(:,:) = 0._wp
384         zwb_fk(:,:) = 0._wp ; zvel_mle(:,:) = 0._wp; zdiff_mle(:,:) = 0._wp
385         zhmle(:,:) = 0._wp  ; zmld(:,:) = 0._wp
386      ENDIF
387      zwb_fk_b(:,:) = 0._wp   ! must be initialised even with ln_osm_mle=F as used in zdf_osm_calculate_dhdt
388
389      ! Flux-Gradient arrays.
390      zsc_wth_1(:,:)  = 0._wp ; zsc_ws_1(:,:)   = 0._wp ; zsc_uw_1(:,:)   = 0._wp
391      zsc_uw_2(:,:)   = 0._wp ; zsc_vw_1(:,:)   = 0._wp ; zsc_vw_2(:,:)   = 0._wp
392      zhbl_t(:,:)     = 0._wp ; zdhdt(:,:)      = 0._wp
393
394      zdiffut(:,:,:) = 0._wp ; zviscos(:,:,:) = 0._wp ; ghamt(:,:,:) = 0._wp
395      ghams(:,:,:)   = 0._wp ; ghamu(:,:,:)   = 0._wp ; ghamv(:,:,:) = 0._wp
396
397       ! hbl = MAX(hbl,epsln)
398      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
399      ! Calculate boundary layer scales
400      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
401
402      ! Assume two-band radiation model for depth of OSBL
403     zz0 =       rn_abs       ! surface equi-partition in 2-bands
404     zz1 =  1. - rn_abs
405     DO jj = 2, jpjm1
406        DO ji = 2, jpim1
407           ! Surface downward irradiance (so always +ve)
408           zrad0(ji,jj) = qsr(ji,jj) * r1_rau0_rcp
409           ! Downwards irradiance at base of boundary layer
410           zradh(ji,jj) = zrad0(ji,jj) * ( zz0 * EXP( -hbl(ji,jj)/rn_si0 ) + zz1 * EXP( -hbl(ji,jj)/rn_si1) )
411           ! Downwards irradiance averaged over depth of the OSBL
412           zradav(ji,jj) = zrad0(ji,jj) * ( zz0 * ( 1.0 - EXP( -hbl(ji,jj)/rn_si0 ) )*rn_si0 &
413                 &                         + zz1 * ( 1.0 - EXP( -hbl(ji,jj)/rn_si1 ) )*rn_si1 ) / hbl(ji,jj)
414        END DO
415     END DO
416     ! Turbulent surface fluxes and fluxes averaged over depth of the OSBL
417     DO jj = 2, jpjm1
418        DO ji = 2, jpim1
419           zthermal = rab_n(ji,jj,1,jp_tem)
420           zbeta    = rab_n(ji,jj,1,jp_sal)
421           ! Upwards surface Temperature flux for non-local term
422           zwth0(ji,jj) = - qns(ji,jj) * r1_rau0_rcp * tmask(ji,jj,1)
423           ! Upwards surface salinity flux for non-local term
424           zws0(ji,jj) = - ( ( emp(ji,jj)-rnf(ji,jj) ) * tsn(ji,jj,1,jp_sal)  + sfx(ji,jj) ) * r1_rau0 * tmask(ji,jj,1)
425           ! Non radiative upwards surface buoyancy flux
426           zwb0(ji,jj) = grav * zthermal * zwth0(ji,jj) -  grav * zbeta * zws0(ji,jj)
427           ! Total upwards surface buoyancy flux
428           zwb0tot(ji,jj) = zwb0(ji,jj) -  grav * zthermal * zrad0(ji,jj)
429           ! turbulent heat flux averaged over depth of OSBL
430           zwthav(ji,jj) = 0.5 * zwth0(ji,jj) - ( 0.5*( zrad0(ji,jj) + zradh(ji,jj) ) - zradav(ji,jj) )
431           ! turbulent salinity flux averaged over depth of the OBSL
432           zwsav(ji,jj) = 0.5 * zws0(ji,jj)
433           ! turbulent buoyancy flux averaged over the depth of the OBSBL
434           zwbav(ji,jj) = grav  * zthermal * zwthav(ji,jj) - grav  * zbeta * zwsav(ji,jj)
435           ! Surface upward velocity fluxes
436           zuw0(ji,jj) = - 0.5 * (utau(ji-1,jj) + utau(ji,jj)) * r1_rau0 * tmask(ji,jj,1)
437           zvw0(ji,jj) = - 0.5 * (vtau(ji,jj-1) + vtau(ji,jj)) * r1_rau0 * tmask(ji,jj,1)
438           ! Friction velocity (zustar), at T-point : LMD94 eq. 2
439           zustar(ji,jj) = MAX( SQRT( SQRT( zuw0(ji,jj) * zuw0(ji,jj) + zvw0(ji,jj) * zvw0(ji,jj) ) ), 1.0e-8 )
440           zcos_wind(ji,jj) = -zuw0(ji,jj) / ( zustar(ji,jj) * zustar(ji,jj) )
441           zsin_wind(ji,jj) = -zvw0(ji,jj) / ( zustar(ji,jj) * zustar(ji,jj) )
442        END DO
443     END DO
444     ! Calculate Stokes drift in direction of wind (zustke) and Stokes penetration depth (dstokes)
445     SELECT CASE (nn_osm_wave)
446     ! Assume constant La#=0.3
447     CASE(0)
448        DO jj = 2, jpjm1
449           DO ji = 2, jpim1
450              zus_x = zcos_wind(ji,jj) * zustar(ji,jj) / 0.3**2
451              zus_y = zsin_wind(ji,jj) * zustar(ji,jj) / 0.3**2
452              ! Linearly
453              zustke(ji,jj) = MAX ( SQRT( zus_x*zus_x + zus_y*zus_y), 1.0e-8 )
454              dstokes(ji,jj) = rn_osm_dstokes
455           END DO
456        END DO
457     ! Assume Pierson-Moskovitz wind-wave spectrum
458     CASE(1)
459        DO jj = 2, jpjm1
460           DO ji = 2, jpim1
461              ! Use wind speed wndm included in sbc_oce module
462              zustke(ji,jj) =  MAX ( 0.016 * wndm(ji,jj), 1.0e-8 )
463              dstokes(ji,jj) = MAX ( 0.12 * wndm(ji,jj)**2 / grav, 5.e-1)
464           END DO
465        END DO
466     ! Use ECMWF wave fields as output from SBCWAVE
467     CASE(2)
468        zfac =  2.0_wp * rpi / 16.0_wp
469
470        DO jj = 2, jpjm1
471           DO ji = 2, jpim1
472              IF (hsw(ji,jj) > 1.e-4) THEN
473                 ! Use  wave fields
474                 zabsstke = SQRT(ut0sd(ji,jj)**2 + vt0sd(ji,jj)**2)
475                 zustke(ji,jj) = MAX ( ( zcos_wind(ji,jj) * ut0sd(ji,jj) + zsin_wind(ji,jj)  * vt0sd(ji,jj) ), 1.0e-8)
476                 dstokes(ji,jj) = MAX (zfac * hsw(ji,jj)*hsw(ji,jj) / ( MAX(zabsstke * wmp(ji,jj), 1.0e-7 ) ), 5.0e-1)
477              ELSE
478                 ! Assume masking issue (e.g. ice in ECMWF reanalysis but not in model run)
479                 ! .. so default to Pierson-Moskowitz
480                 zustke(ji,jj) = MAX ( 0.016 * wndm(ji,jj), 1.0e-8 )
481                 dstokes(ji,jj) = MAX ( 0.12 * wndm(ji,jj)**2 / grav, 5.e-1)
482              END IF
483            END DO
484        END DO
485     END SELECT
486
487     IF (ln_zdfosm_ice_shelter) THEN
488        ! Reduce both Stokes drift and its depth scale by ocean fraction to represent sheltering by ice
489        DO jj = 2, jpjm1
490           DO ji = 2, jpim1
491              zustke(ji,jj) = zustke(ji,jj) * (1.0_wp - fr_i(ji,jj))
492              dstokes(ji,jj) = dstokes(ji,jj) * (1.0_wp - fr_i(ji,jj))
493           END DO
494        END DO
495     END IF
496
497     SELECT CASE (nn_osm_SD_reduce)
498     ! Reduce surface Stokes drift by a constant factor or following Breivik (2016) + van  Roekel (2012) or Grant (2020).
499     CASE(0)
500        ! The Langmur number from the ECMWF model (or from PM)  appears to give La<0.3 for wind-driven seas.
501        !    The coefficient rn_zdfosm_adjust_sd = 0.8 gives La=0.3  in this situation.
502        ! It could represent the effects of the spread of wave directions
503        ! around the mean wind. The effect of this adjustment needs to be tested.
504        IF(nn_osm_wave > 0) THEN
505           zustke(2:jpim1,2:jpjm1) = rn_zdfosm_adjust_sd * zustke(2:jpim1,2:jpjm1)
506        END IF
507     CASE(1)
508        ! van  Roekel (2012): consider average SD over top 10% of boundary layer
509        ! assumes approximate depth profile of SD from Breivik (2016)
510        zsqrtpi = SQRT(rpi)
511        z_two_thirds = 2.0_wp / 3.0_wp
512
513        DO jj = 2, jpjm1
514           DO ji = 2, jpim1
515              zthickness = rn_osm_hblfrac*hbl(ji,jj)
516              z2k_times_thickness =  zthickness * 2.0_wp / MAX( ABS( 5.97_wp * dstokes(ji,jj) ), 0.0000001_wp )
517              zsqrt_depth = SQRT(z2k_times_thickness)
518              zexp_depth  = EXP(-z2k_times_thickness)
519              zustke(ji,jj) = zustke(ji,jj) * (1.0_wp - zexp_depth  &
520                   &              - z_two_thirds * ( zsqrtpi*zsqrt_depth*z2k_times_thickness * ERFC(zsqrt_depth) &
521                   &              + 1.0_wp - (1.0_wp + z2k_times_thickness)*zexp_depth ) ) / z2k_times_thickness
522
523           END DO
524        END DO
525     CASE(2)
526        ! Grant (2020): Match to exponential with same SD and d/dz(Sd) at depth 10% of boundary layer
527        ! assumes approximate depth profile of SD from Breivik (2016)
528        zsqrtpi = SQRT(rpi)
529
530        DO jj = 2, jpjm1
531           DO ji = 2, jpim1
532              zthickness = rn_osm_hblfrac*hbl(ji,jj)
533              z2k_times_thickness =  zthickness * 2.0_wp / MAX( ABS( 5.97_wp * dstokes(ji,jj) ), 0.0000001_wp )
534
535              IF(z2k_times_thickness < 50._wp) THEN
536                 zsqrt_depth = SQRT(z2k_times_thickness)
537                 zexperfc = zsqrtpi * zsqrt_depth * ERFC(zsqrt_depth) * EXP(z2k_times_thickness)
538              ELSE
539                 ! asymptotic expansion of sqrt(pi)*zsqrt_depth*EXP(z2k_times_thickness)*ERFC(zsqrt_depth) for large z2k_times_thickness
540                 ! See Abramowitz and Stegun, Eq. 7.1.23
541                 ! zexperfc = 1._wp - (1/2)/(z2k_times_thickness)  + (3/4)/(z2k_times_thickness**2) - (15/8)/(z2k_times_thickness**3)
542                 zexperfc = ((- 1.875_wp/z2k_times_thickness + 0.75_wp)/z2k_times_thickness - 0.5_wp)/z2k_times_thickness + 1.0_wp
543              END IF
544              zf = z2k_times_thickness*(1.0_wp/zexperfc - 1.0_wp)
545              dstokes(ji,jj) = 5.97 * zf * dstokes(ji,jj)
546              zustke(ji,jj) = zustke(ji,jj) * EXP(z2k_times_thickness * ( 1.0_wp / (2. * zf) - 1.0_wp )) * ( 1.0_wp - zexperfc)
547           END DO
548        END DO
549     END SELECT
550
551     ! Langmuir velocity scale (zwstrl), La # (zla)
552     ! mixed scale (zvstr), convective velocity scale (zwstrc)
553     DO jj = 2, jpjm1
554        DO ji = 2, jpim1
555           ! Langmuir velocity scale (zwstrl), at T-point
556           zwstrl(ji,jj) = ( zustar(ji,jj) * zustar(ji,jj) * zustke(ji,jj) )**pthird
557           zla(ji,jj) = MAX(MIN(SQRT ( zustar(ji,jj) / ( zwstrl(ji,jj) + epsln ) )**3, 4.0), 0.2)
558           IF(zla(ji,jj) > 0.45) dstokes(ji,jj) = MIN(dstokes(ji,jj), 0.5_wp*hbl(ji,jj))
559           ! Velocity scale that tends to zustar for large Langmuir numbers
560           zvstr(ji,jj) = ( zwstrl(ji,jj)**3  + &
561                & ( 1.0 - EXP( -0.5 * zla(ji,jj)**2 ) ) * zustar(ji,jj) * zustar(ji,jj) * zustar(ji,jj) )**pthird
562
563           ! limit maximum value of Langmuir number as approximate treatment for shear turbulence.
564           ! Note zustke and zwstrl are not amended.
565           !
566           ! get convective velocity (zwstrc), stabilty scale (zhol) and logical conection flag lconv
567           IF ( zwbav(ji,jj) > 0.0) THEN
568              zwstrc(ji,jj) = ( 2.0 * zwbav(ji,jj) * 0.9 * hbl(ji,jj) )**pthird
569              zhol(ji,jj) = -0.9 * hbl(ji,jj) * 2.0 * zwbav(ji,jj) / (zvstr(ji,jj)**3 + epsln )
570            ELSE
571              zhol(ji,jj) = -hbl(ji,jj) *  2.0 * zwbav(ji,jj)/ (zvstr(ji,jj)**3  + epsln )
572           ENDIF
573        END DO
574     END DO
575
576     !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
577     ! Mixed-layer model - calculate averages over the boundary layer, and the change in the boundary layer depth
578     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
579     ! BL must be always 4 levels deep.
580     ! For calculation of lateral buoyancy gradients for FK in
581     ! zdf_osm_zmld_horizontal_gradients need halo values for ibld, so must
582     ! previously exist for hbl also.
583
584     ! agn 23/6/20: not clear all this is needed, as hbl checked after it is re-calculated anyway
585     ! ##########################################################################
586      hbl(:,:) = MAX(hbl(:,:), gdepw_n(:,:,4) )
587      ibld(:,:) = 4
588      DO jk = 5, jpkm1
589         DO jj = 1, jpj
590            DO ji = 1, jpi
591               IF ( hbl(ji,jj) >= gdepw_n(ji,jj,jk) ) THEN
592                  ibld(ji,jj) = MIN(mbkt(ji,jj), jk)
593               ENDIF
594            END DO
595         END DO
596      END DO
597     ! ##########################################################################
598
599      DO jj = 2, jpjm1
600         DO ji = 2, jpim1
601            zhbl(ji,jj) = gdepw_n(ji,jj,ibld(ji,jj))
602            imld(ji,jj) = MAX(3,ibld(ji,jj) - MAX( INT( dh(ji,jj) / e3t_n(ji, jj, ibld(ji,jj) )) , 1 ))
603            zhml(ji,jj) = gdepw_n(ji,jj,imld(ji,jj))
604            zdh(ji,jj) = zhbl(ji,jj) - zhml(ji,jj)
605         END DO
606      END DO
607      ! Averages over well-mixed and boundary layer
608      jp_ext(:,:) = 2
609      CALL zdf_osm_vertical_average(ibld, jp_ext, zt_bl, zs_bl, zb_bl, zu_bl, zv_bl, zdt_bl, zds_bl, zdb_bl, zdu_bl, zdv_bl)
610!      jp_ext(:,:) = ibld(:,:) - imld(:,:) + 1
611      CALL zdf_osm_vertical_average(imld-1, ibld-imld+1, zt_ml, zs_ml, zb_ml, zu_ml, zv_ml, zdt_ml, zds_ml, zdb_ml, zdu_ml, zdv_ml)
612! Velocity components in frame aligned with surface stress.
613      CALL zdf_osm_velocity_rotation( zcos_wind, zsin_wind, zu_ml, zv_ml, zdu_ml, zdv_ml )
614      CALL zdf_osm_velocity_rotation( zcos_wind, zsin_wind, zu_bl, zv_bl, zdu_bl, zdv_bl )
615! Determine the state of the OSBL, stable/unstable, shear/no shear
616      CALL zdf_osm_osbl_state( lconv, lshear, j_ddh, zwb_ent, zwb_min, zshear )
617
618      IF ( ln_osm_mle ) THEN
619! Fox-Kemper Scheme
620         mld_prof = 4
621         DO jk = 5, jpkm1
622           DO jj = 2, jpjm1
623             DO ji = 2, jpim1
624               IF ( hmle(ji,jj) >= gdepw_n(ji,jj,jk) ) mld_prof(ji,jj) = MIN(mbkt(ji,jj), jk)
625             END DO
626           END DO
627         END DO
628         jp_ext_mle(:,:) = 2
629        CALL zdf_osm_vertical_average(mld_prof, jp_ext_mle, zt_mle, zs_mle, zb_mle, zu_mle, zv_mle, zdt_mle, zds_mle, zdb_mle, zdu_mle, zdv_mle)
630
631         DO jj = 2, jpjm1
632           DO ji = 2, jpim1
633              zhmle(ji,jj) = gdepw_n(ji,jj,mld_prof(ji,jj))
634           END DO
635         END DO
636
637!! External gradient
638         CALL zdf_osm_external_gradients( ibld+2, zdtdz_bl_ext, zdsdz_bl_ext, zdbdz_bl_ext )
639         CALL zdf_osm_zmld_horizontal_gradients( zmld, zdtdx, zdtdy, zdsdx, zdsdy, dbdx_mle, dbdy_mle, zdbds_mle )
640         CALL zdf_osm_external_gradients( mld_prof, zdtdz_mle_ext, zdsdz_mle_ext, zdbdz_mle_ext )
641         CALL zdf_osm_osbl_state_fk( lpyc, lflux, lmle, zwb_fk )
642         CALL zdf_osm_mle_parameters( zmld, mld_prof, hmle, zhmle, zvel_mle, zdiff_mle )
643      ELSE    ! ln_osm_mle
644! FK not selected, Boundary Layer only.
645         lpyc(:,:) = .TRUE.
646         lflux(:,:) = .FALSE.
647         lmle(:,:) = .FALSE.
648         DO jj = 2, jpjm1
649           DO ji = 2, jpim1
650             IF ( lconv(ji,jj) .AND. zdb_bl(ji,jj) < rn_osm_bl_thresh ) lpyc(ji,jj) = .FALSE.
651           END DO
652         END DO
653      ENDIF   ! ln_osm_mle
654
655! Test if pycnocline well resolved
656      DO jj = 2, jpjm1
657        DO ji = 2,jpim1
658          IF (lconv(ji,jj) ) THEN
659             ztmp = 0.2 * zhbl(ji,jj) / e3w_n(ji,jj,ibld(ji,jj))
660             IF ( ztmp > 6 ) THEN
661      ! pycnocline well resolved
662               jp_ext(ji,jj) = 1
663             ELSE
664      ! pycnocline poorly resolved
665               jp_ext(ji,jj) = 0
666             ENDIF
667          ELSE
668      ! Stable conditions
669            jp_ext(ji,jj) = 0
670          ENDIF
671        END DO
672      END DO
673
674      CALL zdf_osm_vertical_average(ibld, jp_ext, zt_bl, zs_bl, zb_bl, zu_bl, zv_bl, zdt_bl, zds_bl, zdb_bl, zdu_bl, zdv_bl )
675!      jp_ext = ibld-imld+1
676      CALL zdf_osm_vertical_average(imld-1, ibld-imld+1, zt_ml, zs_ml, zb_ml, zu_ml, zv_ml, zdt_ml, zds_ml, zdb_ml, zdu_ml, zdv_ml)
677! Rate of change of hbl
678      CALL zdf_osm_calculate_dhdt( zdhdt )
679      DO jj = 2, jpjm1
680        DO ji = 2, jpim1
681          zhbl_t(ji,jj) = hbl(ji,jj) + (zdhdt(ji,jj) - wn(ji,jj,ibld(ji,jj)))* rn_rdt ! certainly need wn here, so subtract it
682               ! adjustment to represent limiting by ocean bottom
683          IF ( zhbl_t(ji,jj) >= gdepw_n(ji, jj, mbkt(ji,jj) + 1 ) ) THEN
684             zhbl_t(ji,jj) = MIN(zhbl_t(ji,jj), gdepw_n(ji,jj, mbkt(ji,jj) + 1) - depth_tol)! ht_n(:,:))
685             lpyc(ji,jj) = .FALSE.
686          ENDIF
687         END DO
688       END DO
689
690      imld(:,:) = ibld(:,:)           ! use imld to hold previous blayer index
691      ibld(:,:) = 4
692
693      DO jk = 4, jpkm1
694         DO jj = 2, jpjm1
695            DO ji = 2, jpim1
696               IF ( zhbl_t(ji,jj) >= gdepw_n(ji,jj,jk) ) THEN
697                  ibld(ji,jj) = jk
698               ENDIF
699            END DO
700         END DO
701      END DO
702
703!
704! Step through model levels taking account of buoyancy change to determine the effect on dhdt
705!
706      CALL zdf_osm_timestep_hbl( zdhdt )
707! is external level in bounds?
708
709      CALL zdf_osm_vertical_average( ibld, jp_ext, zt_bl, zs_bl, zb_bl, zu_bl, zv_bl, zdt_bl, zds_bl, zdb_bl, zdu_bl, zdv_bl )
710!
711!
712! Check to see if lpyc needs to be changed
713
714      CALL zdf_osm_pycnocline_thickness( dh, zdh )
715
716      DO jj = 2, jpjm1
717        DO ji = 2, jpim1
718          IF ( zdb_bl(ji,jj) < rn_osm_bl_thresh .or. ibld(ji,jj) + jp_ext(ji,jj) >= mbkt(ji,jj) .or. ibld(ji,jj)-imld(ji,jj) == 1 ) lpyc(ji,jj) = .FALSE. 
719        END DO
720      END DO       
721
722      dstokes(:,:) = MIN ( dstokes(:,:), hbl(:,:)/3. )  !  Limit delta for shallow boundary layers for calculating flux-gradient terms.
723!
724    ! Average over the depth of the mixed layer in the convective boundary layer
725!      jp_ext = ibld - imld +1
726      CALL zdf_osm_vertical_average( imld-1, ibld-imld+1, zt_ml, zs_ml, zb_ml, zu_ml, zv_ml, zdt_ml, zds_ml, zdb_ml, zdu_ml, zdv_ml )
727    ! rotate mean currents and changes onto wind align co-ordinates
728    !
729     CALL zdf_osm_velocity_rotation( zcos_wind, zsin_wind, zu_ml, zv_ml, zdu_ml, zdv_ml )
730     CALL zdf_osm_velocity_rotation( zcos_wind, zsin_wind, zu_bl, zv_bl, zdu_bl, zdv_bl )
731      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
732      !  Pycnocline gradients for scalars and velocity
733      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
734
735      CALL zdf_osm_external_gradients( ibld+2, zdtdz_bl_ext, zdsdz_bl_ext, zdbdz_bl_ext )
736      CALL zdf_osm_pycnocline_scalar_profiles( zdtdz_pyc, zdsdz_pyc, zdbdz_pyc, zalpha_pyc )
737      CALL zdf_osm_pycnocline_shear_profiles( zdudz_pyc, zdvdz_pyc )
738       !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
739       ! Eddy viscosity/diffusivity and non-gradient terms in the flux-gradient relationship
740       !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
741       CALL zdf_osm_diffusivity_viscosity( zdiffut, zviscos )
742
743       !
744       ! calculate non-gradient components of the flux-gradient relationships
745       !
746! Stokes term in scalar flux, flux-gradient relationship
747       WHERE ( lconv )
748          zsc_wth_1 = zwstrl**3 * zwth0 / ( zvstr**3 + 0.5 * zwstrc**3 + epsln)
749          !
750          zsc_ws_1 = zwstrl**3 * zws0 / ( zvstr**3 + 0.5 * zwstrc**3 + epsln )
751       ELSEWHERE
752          zsc_wth_1 = 2.0 * zwthav
753          !
754          zsc_ws_1 = 2.0 * zwsav
755       ENDWHERE
756
757
758       DO jj = 2, jpjm1
759          DO ji = 2, jpim1
760            IF ( lconv(ji,jj) ) THEN
761              DO jk = 2, imld(ji,jj)
762                 zznd_d = gdepw_n(ji,jj,jk) / dstokes(ji,jj)
763                 ghamt(ji,jj,jk) = ghamt(ji,jj,jk) + 1.35 * EXP ( -zznd_d ) * ( 1.0 - EXP ( -2.0 * zznd_d ) ) * zsc_wth_1(ji,jj)
764                 !
765                 ghams(ji,jj,jk) = ghams(ji,jj,jk) + 1.35 * EXP ( -zznd_d ) * ( 1.0 - EXP ( -2.0 * zznd_d ) ) *  zsc_ws_1(ji,jj)
766              END DO ! end jk loop
767            ELSE     ! else for if (lconv)
768 ! Stable conditions
769               DO jk = 2, ibld(ji,jj)
770                  zznd_d=gdepw_n(ji,jj,jk) / dstokes(ji,jj)
771                  ghamt(ji,jj,jk) = ghamt(ji,jj,jk) + 2.15 * EXP ( -0.85 * zznd_d ) &
772                       &          *                 ( 1.0 - EXP ( -4.0 * zznd_d ) ) * zsc_wth_1(ji,jj)
773                  !
774                  ghams(ji,jj,jk) = ghams(ji,jj,jk) + 2.15 * EXP ( -0.85 * zznd_d ) &
775                       &          *                 ( 1.0 - EXP ( -4.0 * zznd_d ) ) *  zsc_ws_1(ji,jj)
776               END DO
777            ENDIF               ! endif for check on lconv
778
779          END DO  ! end of ji loop
780       END DO     ! end of jj loop
781
782! Stokes term in flux-gradient relationship (note in zsc_uw_n don't use zvstr since term needs to go to zero as zwstrl goes to zero)
783       WHERE ( lconv )
784          zsc_uw_1 = ( zwstrl**3 + 0.5 * zwstrc**3 )**pthird * zustke / MAX( ( 1.0 - 1.0 * 6.5 * zla**(8.0/3.0) ), 0.2 )
785          zsc_uw_2 = ( zwstrl**3 + 0.5 * zwstrc**3 )**pthird * zustke / MIN( zla**(8.0/3.0) + epsln, 0.12 )
786          zsc_vw_1 = ff_t * zhml * zustke**3 * MIN( zla**(8.0/3.0), 0.12 ) / ( ( zvstr**3 + 0.5 * zwstrc**3 )**(2.0/3.0) + epsln )
787       ELSEWHERE
788          zsc_uw_1 = zustar**2
789          zsc_vw_1 = ff_t * zhbl * zustke**3 * MIN( zla**(8.0/3.0), 0.12 ) / (zvstr**2 + epsln)
790       ENDWHERE
791       IF(ln_dia_osm) THEN
792          IF ( iom_use("ghamu_00") ) CALL iom_put( "ghamu_00", wmask*ghamu )
793          IF ( iom_use("ghamv_00") ) CALL iom_put( "ghamv_00", wmask*ghamv )
794       END IF
795       DO jj = 2, jpjm1
796          DO ji = 2, jpim1
797             IF ( lconv(ji,jj) ) THEN
798                DO jk = 2, imld(ji,jj)
799                   zznd_d = gdepw_n(ji,jj,jk) / dstokes(ji,jj)
800                   ghamu(ji,jj,jk) = ghamu(ji,jj,jk) +      ( -0.05 * EXP ( -0.4 * zznd_d )   * zsc_uw_1(ji,jj)   &
801                        &          +                        0.00125 * EXP (      - zznd_d )   * zsc_uw_2(ji,jj) ) &
802                        &          *                          ( 1.0 - EXP ( -2.0 * zznd_d ) )
803!
804                   ghamv(ji,jj,jk) = ghamv(ji,jj,jk) - 0.65 *  0.15 * EXP (      - zznd_d )                       &
805                        &          *                          ( 1.0 - EXP ( -2.0 * zznd_d ) ) * zsc_vw_1(ji,jj)
806                END DO   ! end jk loop
807             ELSE
808! Stable conditions
809                DO jk = 2, ibld(ji,jj) ! corrected to ibld
810                   zznd_d = gdepw_n(ji,jj,jk) / dstokes(ji,jj)
811                   ghamu(ji,jj,jk) = ghamu(ji,jj,jk) - 0.75 *   1.3 * EXP ( -0.5 * zznd_d )                       &
812                        &                                   * ( 1.0 - EXP ( -4.0 * zznd_d ) ) * zsc_uw_1(ji,jj)
813                   ghamv(ji,jj,jk) = ghamv(ji,jj,jk) + 0._wp
814                END DO   ! end jk loop
815             ENDIF
816          END DO  ! ji loop
817       END DO     ! jj loo
818
819! Buoyancy term in flux-gradient relationship [note : includes ROI ratio (X0.3) and pressure (X0.5)]
820
821       WHERE ( lconv )
822          zsc_wth_1 = zwbav * zwth0 * ( 1.0 + EXP ( 0.2 * zhol ) ) * zhml / ( zvstr**3 + 0.5 * zwstrc**3 + epsln )
823          zsc_ws_1  = zwbav * zws0  * ( 1.0 + EXP ( 0.2 * zhol ) ) * zhml / ( zvstr**3 + 0.5 * zwstrc**3 + epsln )
824       ELSEWHERE
825          zsc_wth_1 = 0._wp
826          zsc_ws_1 = 0._wp
827       ENDWHERE
828
829       DO jj = 2, jpjm1
830          DO ji = 2, jpim1
831             IF (lconv(ji,jj) ) THEN
832                DO jk = 2, imld(ji,jj)
833                   zznd_ml = gdepw_n(ji,jj,jk) / zhml(ji,jj)
834                   ! calculate turbulent time scale
835                   zl_c = 0.9 * ( 1.0 - EXP ( - 5.0 * ( zznd_ml + zznd_ml**3 / 3.0 ) ) )                                           &
836                        &     * ( 1.0 - EXP ( -15.0 * (     1.2 - zznd_ml          ) ) )
837                   zl_l = 2.0 * ( 1.0 - EXP ( - 2.0 * ( zznd_ml + zznd_ml**3 / 3.0 ) ) )                                           &
838                        &     * ( 1.0 - EXP ( - 8.0 * (     1.15 - zznd_ml          ) ) ) * ( 1.0 + dstokes(ji,jj) / zhml (ji,jj) )
839                   zl_eps = zl_l + ( zl_c - zl_l ) / ( 1.0 + EXP ( -3.0 * LOG10 ( - zhol(ji,jj) ) ) ) ** (3.0 / 2.0)
840                   ! non-gradient buoyancy terms
841                   ghamt(ji,jj,jk) = ghamt(ji,jj,jk) + 0.3 * 0.4 * zsc_wth_1(ji,jj) * zl_eps / ( 0.15 + zznd_ml )
842                   ghams(ji,jj,jk) = ghams(ji,jj,jk) + 0.3 * 0.4 *  zsc_ws_1(ji,jj) * zl_eps / ( 0.15 + zznd_ml )
843                END DO
844               
845                IF ( lpyc(ji,jj) ) THEN
846                  ztau_sc_u(ji,jj) = zhml(ji,jj) / ( zvstr(ji,jj)**3 + zwstrc(ji,jj)**3 )**pthird
847                  ztau_sc_u(ji,jj) = ztau_sc_u(ji,jj) * ( 1.4 -0.4 / ( 1.0 + EXP( -3.5 * LOG10( -zhol(ji,jj) ) ) )**1.5 )
848                  zwth_ent =  -0.003 * ( 0.15 * zvstr(ji,jj)**3 + zwstrc(ji,jj)**3 )**pthird * ( 1.0 - zdh(ji,jj) /zhbl(ji,jj) ) * zdt_ml(ji,jj)                 
849                  zws_ent =  -0.003 * ( 0.15 * zvstr(ji,jj)**3 + zwstrc(ji,jj)**3 )**pthird * ( 1.0 - zdh(ji,jj) /zhbl(ji,jj) ) * zds_ml(ji,jj)
850! Cubic profile used for buoyancy term
851                  DO jk = 2, ibld(ji,jj)
852                    zznd_pyc = -( gdepw_n(ji,jj,jk) - zhbl(ji,jj) ) / zdh(ji,jj)
853                    ghamt(ji,jj,jk) = ghamt(ji,jj,jk) - 0.045 * ( ( zwth_ent * zdbdz_pyc(ji,jj,jk) ) * ztau_sc_u(ji,jj)**2 ) * MAX( ( 1.75 * zznd_pyc -0.15 * zznd_pyc**2 - 0.2 * zznd_pyc**3 ), 0.0 )
854
855                    ghams(ji,jj,jk) = ghams(ji,jj,jk) - 0.045 * ( ( zws_ent * zdbdz_pyc(ji,jj,jk) ) * ztau_sc_u(ji,jj)**2 ) * MAX( ( 1.75 * zznd_pyc -0.15 * zznd_pyc**2 - 0.2 * zznd_pyc**3 ), 0.0 )
856                  END DO
857                  !
858                  IF ( dh(ji,jj)  <  0.2*hbl(ji,jj) ) THEN
859                     zbuoy_pyc_sc = zalpha_pyc(ji,jj) * zdb_ml(ji,jj) / zdh(ji,jj) + zdbdz_bl_ext(ji,jj)
860                     zdelta_pyc = ( zvstr(ji,jj)**3 + zwstrc(ji,jj)**3 )**pthird / SQRT( MAX( zbuoy_pyc_sc, ( zvstr(ji,jj)**3 + zwstrc(ji,jj)**3 )**p2third / zdh(ji,jj)**2 ) )
861!
862                     zwt_pyc_sc_1 = 0.325 * ( zalpha_pyc(ji,jj) * zdt_ml(ji,jj) / zdh(ji,jj) + zdtdz_bl_ext(ji,jj) ) * zdelta_pyc**2 / zdh(ji,jj)
863!
864                     zws_pyc_sc_1 = 0.325 * ( zalpha_pyc(ji,jj) * zds_ml(ji,jj) / zdh(ji,jj) + zdsdz_bl_ext(ji,jj) ) * zdelta_pyc**2 / zdh(ji,jj)
865!
866                     zzeta_pyc = 0.15 - 0.175 / ( 1.0 + EXP( -3.5 * LOG10( -zhol(ji,jj) ) ) ) 
867                     DO jk = 2, ibld(ji,jj)
868                        zznd_pyc = -( gdepw_n(ji,jj,jk) - zhbl(ji,jj) ) / zdh(ji,jj)
869                        ghamt(ji,jj,jk) = ghamt(ji,jj,jk) + 0.05 * zwt_pyc_sc_1 * EXP( -0.25 * ( zznd_pyc / zzeta_pyc )**2 ) * zdh(ji,jj) / ( zvstr(ji,jj)**3 + zwstrc(ji,jj)**3 )**pthird
870!
871                        ghams(ji,jj,jk) = ghams(ji,jj,jk) + 0.05 * zws_pyc_sc_1 * EXP( -0.25 * ( zznd_pyc / zzeta_pyc )**2 ) * zdh(ji,jj) / ( zvstr(ji,jj)**3 + zwstrc(ji,jj)**3 )**pthird
872                     END DO
873                  END IF
874               ENDIF ! End of pycnocline                 
875             ELSE ! lconv test - stable conditions
876                DO jk = 2, ibld(ji,jj)
877                   ghamt(ji,jj,jk) = ghamt(ji,jj,jk) + zsc_wth_1(ji,jj)
878                   ghams(ji,jj,jk) = ghams(ji,jj,jk) +  zsc_ws_1(ji,jj)
879                END DO
880             ENDIF
881          END DO   ! ji loop
882       END DO      ! jj loop
883
884       WHERE ( lconv )
885          zsc_uw_1 = -zwb0 * zustar**2 * zhml / ( zvstr**3 + 0.5 * zwstrc**3 + epsln )
886          zsc_uw_2 =  zwb0 * zustke    * zhml / ( zvstr**3 + 0.5 * zwstrc**3 + epsln )**(2.0/3.0)
887          zsc_vw_1 = 0._wp
888       ELSEWHERE
889         zsc_uw_1 = 0._wp
890         zsc_vw_1 = 0._wp
891       ENDWHERE
892
893       DO jj = 2, jpjm1
894          DO ji = 2, jpim1
895             IF ( lconv(ji,jj) ) THEN
896                DO jk = 2 , imld(ji,jj)
897                   zznd_d = gdepw_n(ji,jj,jk) / dstokes(ji,jj)
898                   ghamu(ji,jj,jk) = ghamu(ji,jj,jk) + 0.3 * 0.5 * ( zsc_uw_1(ji,jj) +   0.125 * EXP( -0.5 * zznd_d )     &
899                        &                                                            * (   1.0 - EXP( -0.5 * zznd_d ) )   &
900                        &                                          * zsc_uw_2(ji,jj)                                    )
901                   ghamv(ji,jj,jk) = ghamv(ji,jj,jk) + zsc_vw_1(ji,jj)
902                END DO  ! jk loop
903             ELSE
904             ! stable conditions
905                DO jk = 2, ibld(ji,jj)
906                   ghamu(ji,jj,jk) = ghamu(ji,jj,jk) + zsc_uw_1(ji,jj)
907                   ghamv(ji,jj,jk) = ghamv(ji,jj,jk) + zsc_vw_1(ji,jj)
908                END DO
909             ENDIF
910          END DO        ! ji loop
911       END DO           ! jj loop
912
913       DO jj = 2, jpjm1
914         DO ji = 2, jpim1
915           IF( lconv(ji,jj) ) THEN
916             IF ( lpyc(ji,jj) ) THEN
917               IF ( j_ddh(ji,jj) == 0 ) THEN
918! Place holding code. Parametrization needs checking for these conditions.
919               zomega = ( 0.15 * zwstrl(ji,jj)**3 + zwstrc(ji,jj)**3 + 4.75 * ( zshear(ji,jj)* zhbl(ji,jj) ))**pthird
920               zuw_bse = -0.0035 * zomega * ( 1.0 - zdh(ji,jj) / zhbl(ji,jj) ) * zdu_ml(ji,jj)
921               zvw_bse = -0.0075 * zomega * ( 1.0 - zdh(ji,jj) / zhbl(ji,jj) ) * zdv_ml(ji,jj)
922             ELSE
923               zomega = ( 0.15 * zwstrl(ji,jj)**3 + zwstrc(ji,jj)**3 + 4.75 * ( zshear(ji,jj)* zhbl(ji,jj) ))**pthird
924               zuw_bse = -0.0035 * zomega * ( 1.0 - zdh(ji,jj) / zhbl(ji,jj) ) * zdu_ml(ji,jj)
925               zvw_bse = -0.0075 * zomega * ( 1.0 - zdh(ji,jj) / zhbl(ji,jj) ) * zdv_ml(ji,jj)
926             ENDIF
927             zd_cubic = zdh(ji,jj) / zhbl(ji,jj) * zuw0(ji,jj) - ( 2.0 + zdh(ji,jj) /zhml(ji,jj) ) * zuw_bse
928             zc_cubic = zuw_bse - zd_cubic
929! need ztau_sc_u to be available. Change to array.
930             DO jk = imld(ji,jj), ibld(ji,jj)
931                zznd_pyc = - ( gdepw_n(ji,jj,jk) - zhbl(ji,jj) ) / zdh(ji,jj)
932                ghamu(ji,jj,jk) = ghamu(ji,jj,jk) - 0.045 * ( ztau_sc_u(ji,jj)**2 ) * zuw_bse * &
933                     & ( zc_cubic * zznd_pyc**2 + zd_cubic * zznd_pyc**3 ) * ( 0.75 + 0.25 * zznd_pyc )**2 * zdbdz_pyc(ji,jj,jk)
934             END DO
935             zvw_max = 0.7 * ff_t(ji,jj) * ( zustke(ji,jj) * dstokes(ji,jj) + 0.75 * zustar(ji,jj) * zhml(ji,jj) )
936             zd_cubic = zvw_max * zdh(ji,jj) / zhml(ji,jj) - ( 2.0 + zdh(ji,jj) /zhml(ji,jj) ) * zvw_bse
937             zc_cubic = zvw_bse - zd_cubic
938             DO jk = imld(ji,jj), ibld(ji,jj)
939               zznd_pyc = -( gdepw_n(ji,jj,jk) -zhbl(ji,jj) ) / zdh(ji,jj)
940               ghamv(ji,jj,jk) = ghamv(ji,jj,jk) - 0.045 * ( ztau_sc_u(ji,jj)**2 ) * zvw_bse *  &
941                    & ( zc_cubic * zznd_pyc**2 + zd_cubic * zznd_pyc**3 ) * ( 0.75 + 0.25 * zznd_pyc )**2 * zdbdz_pyc(ji,jj,jk)
942             END DO
943             ENDIF  ! lpyc
944           ENDIF  ! lconv
945         END DO ! ji loop
946       END DO  ! jj loop
947
948       IF(ln_dia_osm) THEN
949          IF ( iom_use("ghamu_0") ) CALL iom_put( "ghamu_0", wmask*ghamu )
950          IF ( iom_use("zsc_uw_1_0") ) CALL iom_put( "zsc_uw_1_0", tmask(:,:,1)*zsc_uw_1 )
951       END IF
952! Transport term in flux-gradient relationship [note : includes ROI ratio (X0.3) ]
953
954       DO jj = 1, jpjm1
955          DO ji = 1, jpim1
956         
957            IF ( lconv(ji,jj) ) THEN
958              zsc_wth_1(ji,jj) = zwth0(ji,jj) / ( 1.0 - 0.56 * EXP( zhol(ji,jj) ) )
959              zsc_ws_1(ji,jj) = zws0(ji,jj) / (1.0 - 0.56 *EXP( zhol(ji,jj) ) )
960              IF ( lpyc(ji,jj) ) THEN
961! Pycnocline scales
962                 zsc_wth_pyc(ji,jj) = -0.003 * zwstrc(ji,jj) * ( 1.0 - zdh(ji,jj) /zhbl(ji,jj) ) * zdt_ml(ji,jj)
963                 zsc_ws_pyc(ji,jj) = -0.003 * zwstrc(ji,jj) * ( 1.0 - zdh(ji,jj) /zhbl(ji,jj) ) * zds_ml(ji,jj)
964               ENDIF
965            ELSE
966              zsc_wth_1(ji,jj) = 2.0 * zwthav(ji,jj)
967              zsc_ws_1(ji,jj) = zws0(ji,jj)
968            ENDIF
969          END DO
970        END DO
971
972       DO jj = 2, jpjm1
973          DO ji = 2, jpim1
974            IF ( lconv(ji,jj) ) THEN
975               DO jk = 2, imld(ji,jj)
976                  zznd_ml=gdepw_n(ji,jj,jk) / zhml(ji,jj)
977                  ghamt(ji,jj,jk) = ghamt(ji,jj,jk) + 0.3 * zsc_wth_1(ji,jj)                &
978                       &          * ( -2.0 + 2.75 * (       ( 1.0 + 0.6 * zznd_ml**4 )      &
979                       &                               - EXP(     - 6.0 * zznd_ml    ) ) )  &
980                       &          * ( 1.0 - EXP( - 15.0 * (         1.0 - zznd_ml    ) ) )
981                  !
982                  ghams(ji,jj,jk) = ghams(ji,jj,jk) + 0.3 * zsc_ws_1(ji,jj)  &
983                       &          * ( -2.0 + 2.75 * (       ( 1.0 + 0.6 * zznd_ml**4 )      &
984                       &                               - EXP(     - 6.0 * zznd_ml    ) ) )  &
985                       &          * ( 1.0 - EXP ( -15.0 * (         1.0 - zznd_ml    ) ) )
986               END DO
987               !
988               ! may need to comment out lpyc block
989               IF ( lpyc(ji,jj) ) THEN
990! pycnocline
991                 DO jk = imld(ji,jj), ibld(ji,jj)
992                   zznd_pyc = - ( gdepw_n(ji,jj,jk) - zhbl(ji,jj) ) / zdh(ji,jj)
993                   ghamt(ji,jj,jk) = ghamt(ji,jj,jk) + 4.0 * zsc_wth_pyc(ji,jj) * ( 0.48 - EXP( -1.5 * ( zznd_pyc -0.3)**2 ) ) 
994                   ghams(ji,jj,jk) = ghams(ji,jj,jk) + 4.0 * zsc_ws_pyc(ji,jj) * ( 0.48 - EXP( -1.5 * ( zznd_pyc -0.3)**2 ) ) 
995                 END DO
996              ENDIF
997            ELSE
998               IF( zdhdt(ji,jj) > 0. ) THEN
999                 DO jk = 2, ibld(ji,jj)
1000                    zznd_d = gdepw_n(ji,jj,jk) / dstokes(ji,jj)
1001                    znd = gdepw_n(ji,jj,jk) / zhbl(ji,jj)
1002                    ghamt(ji,jj,jk) = ghamt(ji,jj,jk) + 0.3 * ( -4.06 * EXP( -2.0 * zznd_d ) * (1.0 - EXP( -4.0 * zznd_d ) ) + &
1003                 &  7.5 * EXP ( -10.0 * ( 0.95 - znd )**2 ) * ( 1.0 - znd ) ) * zsc_wth_1(ji,jj)
1004                    ghams(ji,jj,jk) = ghams(ji,jj,jk) + 0.3 * ( -4.06 * EXP( -2.0 * zznd_d ) * (1.0 - EXP( -4.0 * zznd_d ) ) + &
1005                  &  7.5 * EXP ( -10.0 * ( 0.95 - znd )**2 ) * ( 1.0 - znd ) ) * zsc_ws_1(ji,jj)
1006                 END DO
1007               ENDIF
1008            ENDIF
1009          ENDDO    ! ji loop
1010       END DO      ! jj loop
1011
1012       WHERE ( lconv )
1013          zsc_uw_1 = zustar**2
1014          zsc_vw_1 = ff_t * zustke * zhml
1015       ELSEWHERE
1016          zsc_uw_1 = zustar**2
1017          zsc_uw_2 = (2.25 - 3.0 * ( 1.0 - EXP( -1.25 * 2.0 ) ) ) * ( 1.0 - EXP( -4.0 * 2.0 ) ) * zsc_uw_1
1018          zsc_vw_1 = ff_t * zustke * zhbl
1019          zsc_vw_2 = -0.11 * SIN( 3.14159 * ( 2.0 + 0.4 ) ) * EXP(-( 1.5 + 2.0 )**2 ) * zsc_vw_1
1020       ENDWHERE
1021
1022       DO jj = 2, jpjm1
1023          DO ji = 2, jpim1
1024             IF ( lconv(ji,jj) ) THEN
1025               DO jk = 2, imld(ji,jj)
1026                  zznd_ml = gdepw_n(ji,jj,jk) / zhml(ji,jj)
1027                  zznd_d = gdepw_n(ji,jj,jk) / dstokes(ji,jj)
1028                  ghamu(ji,jj,jk) = ghamu(ji,jj,jk)&
1029                       & + 0.3 * ( -2.0 + 2.5 * ( 1.0 + 0.1 * zznd_ml**4 ) - EXP ( -8.0 * zznd_ml ) ) * zsc_uw_1(ji,jj)
1030                  !
1031                  ghamv(ji,jj,jk) = ghamv(ji,jj,jk)&
1032                       & + 0.3 * 0.1 * ( EXP( -zznd_d ) + EXP( -5.0 * ( 1.0 - zznd_ml ) ) ) * zsc_vw_1(ji,jj)
1033               END DO
1034
1035             ELSE
1036               DO jk = 2, ibld(ji,jj)
1037                  znd = gdepw_n(ji,jj,jk) / zhbl(ji,jj)
1038                  zznd_d = gdepw_n(ji,jj,jk) / dstokes(ji,jj)
1039                  IF ( zznd_d <= 2.0 ) THEN
1040                     ghamu(ji,jj,jk) = ghamu(ji,jj,jk) + 0.5 * 0.3 &
1041                          &*  ( 2.25 - 3.0  * ( 1.0 - EXP( - 1.25 * zznd_d ) ) * ( 1.0 - EXP( -2.0 * zznd_d ) ) ) * zsc_uw_1(ji,jj)
1042                     !
1043                  ELSE
1044                     ghamu(ji,jj,jk) = ghamu(ji,jj,jk)&
1045                          & + 0.5 * 0.3 * ( 1.0 - EXP( -5.0 * ( 1.0 - znd ) ) ) * zsc_uw_2(ji,jj)
1046                     !
1047                  ENDIF
1048
1049                  ghamv(ji,jj,jk) = ghamv(ji,jj,jk)&
1050                       & + 0.3 * 0.15 * SIN( 3.14159 * ( 0.65 * zznd_d ) ) * EXP( -0.25 * zznd_d**2 ) * zsc_vw_1(ji,jj)
1051                  ghamv(ji,jj,jk) = ghamv(ji,jj,jk)&
1052                       & + 0.3 * 0.15 * EXP( -5.0 * ( 1.0 - znd ) ) * ( 1.0 - EXP( -20.0 * ( 1.0 - znd ) ) ) * zsc_vw_2(ji,jj)
1053               END DO
1054             ENDIF
1055          END DO
1056       END DO
1057
1058       IF(ln_dia_osm) THEN
1059          IF ( iom_use("ghamu_f") ) CALL iom_put( "ghamu_f", wmask*ghamu )
1060          IF ( iom_use("ghamv_f") ) CALL iom_put( "ghamv_f", wmask*ghamv )
1061          IF ( iom_use("zsc_uw_1_f") ) CALL iom_put( "zsc_uw_1_f", tmask(:,:,1)*zsc_uw_1 )
1062          IF ( iom_use("zsc_vw_1_f") ) CALL iom_put( "zsc_vw_1_f", tmask(:,:,1)*zsc_vw_1 )
1063          IF ( iom_use("zsc_uw_2_f") ) CALL iom_put( "zsc_uw_2_f", tmask(:,:,1)*zsc_uw_2 )
1064          IF ( iom_use("zsc_vw_2_f") ) CALL iom_put( "zsc_vw_2_f", tmask(:,:,1)*zsc_vw_2 )
1065       END IF
1066!
1067! Make surface forced velocity non-gradient terms go to zero at the base of the mixed layer.
1068
1069
1070 ! Make surface forced velocity non-gradient terms go to zero at the base of the boundary layer.
1071
1072      DO jj = 2, jpjm1
1073         DO ji = 2, jpim1
1074            IF ( .not. lconv(ji,jj) ) THEN
1075               DO jk = 2, ibld(ji,jj)
1076                  znd = -( gdepw_n(ji,jj,jk) - zhbl(ji,jj) ) / zhbl(ji,jj) !ALMG to think about
1077                  IF ( znd >= 0.0 ) THEN
1078                     ghamu(ji,jj,jk) = ghamu(ji,jj,jk) * ( 1.0 - EXP( -10.0 * znd**2 ) )
1079                     ghamv(ji,jj,jk) = ghamv(ji,jj,jk) * ( 1.0 - EXP( -10.0 * znd**2 ) )
1080                  ELSE
1081                     ghamu(ji,jj,jk) = 0._wp
1082                     ghamv(ji,jj,jk) = 0._wp
1083                  ENDIF
1084               END DO
1085            ENDIF
1086         END DO
1087      END DO
1088
1089      ! pynocline contributions
1090       DO jj = 2, jpjm1
1091          DO ji = 2, jpim1
1092            IF ( .not. lconv(ji,jj) ) THEN
1093             IF ( ibld(ji,jj) + jp_ext(ji,jj) < mbkt(ji,jj) ) THEN
1094                DO jk= 2, ibld(ji,jj)
1095                   ghamt(ji,jj,jk) = ghamt(ji,jj,jk) + zdiffut(ji,jj,jk) * zdtdz_pyc(ji,jj,jk)
1096                   ghams(ji,jj,jk) = ghams(ji,jj,jk) + zdiffut(ji,jj,jk) * zdsdz_pyc(ji,jj,jk)
1097                   ghamu(ji,jj,jk) = ghamu(ji,jj,jk) + zviscos(ji,jj,jk) * zdudz_pyc(ji,jj,jk)
1098                   ghamv(ji,jj,jk) = ghamv(ji,jj,jk) + zviscos(ji,jj,jk) * zdvdz_pyc(ji,jj,jk)
1099                END DO
1100             END IF
1101            END IF
1102          END DO
1103       END DO
1104      IF(ln_dia_osm) THEN
1105          IF ( iom_use("ghamu_b") ) CALL iom_put( "ghamu_b", wmask*ghamu )
1106          IF ( iom_use("ghamv_b") ) CALL iom_put( "ghamv_b", wmask*ghamv )
1107       END IF
1108
1109       DO jj=2, jpjm1
1110          DO ji = 2, jpim1
1111             ghamt(ji,jj,ibld(ji,jj)) = 0._wp
1112             ghams(ji,jj,ibld(ji,jj)) = 0._wp
1113             ghamu(ji,jj,ibld(ji,jj)) = 0._wp
1114             ghamv(ji,jj,ibld(ji,jj)) = 0._wp
1115          END DO       ! ji loop
1116       END DO          ! jj loop
1117
1118       IF(ln_dia_osm) THEN
1119          IF ( iom_use("ghamu_1") ) CALL iom_put( "ghamu_1", wmask*ghamu )
1120          IF ( iom_use("ghamv_1") ) CALL iom_put( "ghamv_1", wmask*ghamv )
1121          IF ( iom_use("zdudz_pyc") ) CALL iom_put( "zdudz_pyc", wmask*zdudz_pyc )
1122          IF ( iom_use("zdvdz_pyc") ) CALL iom_put( "zdvdz_pyc", wmask*zdvdz_pyc )
1123          IF ( iom_use("zviscos") ) CALL iom_put( "zviscos", wmask*zviscos )
1124       END IF
1125       !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
1126       ! Need to put in code for contributions that are applied explicitly to
1127       ! the prognostic variables
1128       !  1. Entrainment flux
1129       !
1130       !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1131
1132
1133
1134       ! rotate non-gradient velocity terms back to model reference frame
1135
1136       DO jj = 2, jpjm1
1137          DO ji = 2, jpim1
1138             DO jk = 2, ibld(ji,jj)
1139                ztemp = ghamu(ji,jj,jk)
1140                ghamu(ji,jj,jk) = ghamu(ji,jj,jk) * zcos_wind(ji,jj) - ghamv(ji,jj,jk) * zsin_wind(ji,jj)
1141                ghamv(ji,jj,jk) = ghamv(ji,jj,jk) * zcos_wind(ji,jj) + ztemp * zsin_wind(ji,jj)
1142             END DO
1143          END DO
1144       END DO
1145
1146       IF(ln_dia_osm) THEN
1147          IF ( iom_use("zdtdz_pyc") ) CALL iom_put( "zdtdz_pyc", wmask*zdtdz_pyc )
1148          IF ( iom_use("zdsdz_pyc") ) CALL iom_put( "zdsdz_pyc", wmask*zdsdz_pyc )
1149          IF ( iom_use("zdbdz_pyc") ) CALL iom_put( "zdbdz_pyc", wmask*zdbdz_pyc )
1150       END IF
1151
1152! KPP-style Ri# mixing
1153       IF( ln_kpprimix) THEN
1154          DO jk = 2, jpkm1           !* Shear production at uw- and vw-points (energy conserving form)
1155             DO jj = 1, jpjm1
1156                DO ji = 1, jpim1   ! vector opt.
1157                   z3du(ji,jj,jk) = 0.5 * (  un(ji,jj,jk-1) -  un(ji  ,jj,jk) )   &
1158                        &                 * (  ub(ji,jj,jk-1) -  ub(ji  ,jj,jk) ) * wumask(ji,jj,jk) &
1159                        &                 / (  e3uw_n(ji,jj,jk) * e3uw_b(ji,jj,jk) )
1160                   z3dv(ji,jj,jk) = 0.5 * (  vn(ji,jj,jk-1) -  vn(ji,jj  ,jk) )   &
1161                        &                 * (  vb(ji,jj,jk-1) -  vb(ji,jj  ,jk) ) * wvmask(ji,jj,jk) &
1162                        &                 / (  e3vw_n(ji,jj,jk) * e3vw_b(ji,jj,jk) )
1163                END DO
1164             END DO
1165          END DO
1166      !
1167         DO jk = 2, jpkm1
1168            DO jj = 2, jpjm1
1169               DO ji = 2, jpim1   ! vector opt.
1170                  !                                          ! shear prod. at w-point weightened by mask
1171                  zesh2  =  ( z3du(ji-1,jj,jk) + z3du(ji,jj,jk) ) / MAX( 1._wp , umask(ji-1,jj,jk) + umask(ji,jj,jk) )   &
1172                     &    + ( z3dv(ji,jj-1,jk) + z3dv(ji,jj,jk) ) / MAX( 1._wp , vmask(ji,jj-1,jk) + vmask(ji,jj,jk) )
1173                  !                                          ! local Richardson number
1174                  zri   = MAX( rn2b(ji,jj,jk), 0._wp ) / MAX(zesh2, epsln)
1175                  zfri =  MIN( zri / rn_riinfty , 1.0_wp )
1176                  zfri  = ( 1.0_wp - zfri * zfri )
1177                  zrimix(ji,jj,jk)  =  zfri * zfri  * zfri * wmask(ji, jj, jk)
1178                END DO
1179             END DO
1180          END DO
1181
1182          DO jj = 2, jpjm1
1183             DO ji = 2, jpim1
1184                DO jk = ibld(ji,jj) + 1, jpkm1
1185                   zdiffut(ji,jj,jk) = zrimix(ji,jj,jk)*rn_difri
1186                   zviscos(ji,jj,jk) = zrimix(ji,jj,jk)*rn_difri
1187                END DO
1188             END DO
1189          END DO
1190
1191       END IF ! ln_kpprimix = .true.
1192
1193! KPP-style set diffusivity large if unstable below BL
1194       IF( ln_convmix) THEN
1195          DO jj = 2, jpjm1
1196             DO ji = 2, jpim1
1197                DO jk = ibld(ji,jj) + 1, jpkm1
1198                  IF(  MIN( rn2(ji,jj,jk), rn2b(ji,jj,jk) ) <= -1.e-12 ) zdiffut(ji,jj,jk) = rn_difconv
1199                END DO
1200             END DO
1201          END DO
1202       END IF ! ln_convmix = .true.
1203
1204
1205
1206       IF ( ln_osm_mle ) THEN  ! set up diffusivity and non-gradient mixing
1207          DO jj = 2 , jpjm1
1208             DO ji = 2, jpim1
1209                 IF ( lflux(ji,jj) ) THEN ! MLE mixing extends below boundary layer
1210                ! Calculate MLE flux contribution from surface fluxes
1211                   DO jk = 1, ibld(ji,jj)
1212                     znd = gdepw_n(ji,jj,jk) / MAX(zhbl(ji,jj),epsln)
1213                     ghamt(ji,jj,jk) = ghamt(ji,jj,jk) - ( zwth0(ji,jj) - zrad0(ji,jj) ) * ( 1.0 - znd )
1214                     ghams(ji,jj,jk) = ghams(ji,jj,jk) - zws0(ji,jj) * ( 1.0 - znd )
1215                    END DO
1216                    DO jk = 1, mld_prof(ji,jj)
1217                      znd = gdepw_n(ji,jj,jk) / MAX(zhmle(ji,jj),epsln)
1218                      ghamt(ji,jj,jk) = ghamt(ji,jj,jk) +  ( zwth0(ji,jj) - zrad0(ji,jj) ) * ( 1.0 - znd )
1219                      ghams(ji,jj,jk) = ghams(ji,jj,jk) + zws0(ji,jj) * ( 1.0 -znd )
1220                    END DO
1221            ! Viscosity for MLEs
1222                    DO jk = 1, mld_prof(ji,jj)
1223                      znd = -gdepw_n(ji,jj,jk) / MAX(zhmle(ji,jj),epsln)
1224                      zdiffut(ji,jj,jk) = zdiffut(ji,jj,jk) + zdiff_mle(ji,jj) * ( 1.0 - ( 2.0 * znd + 1.0 )**2 ) * ( 1.0 + 5.0 / 21.0 * ( 2.0 * znd + 1.0 )** 2 )
1225                    END DO
1226                 ELSE
1227! Surface transports limited to OSBL.                 
1228            ! Viscosity for MLEs
1229                    DO jk = 1, mld_prof(ji,jj)
1230                      znd = -gdepw_n(ji,jj,jk) / MAX(zhmle(ji,jj),epsln)
1231                      zdiffut(ji,jj,jk) = zdiffut(ji,jj,jk) + zdiff_mle(ji,jj) * ( 1.0 - ( 2.0 * znd + 1.0 )**2 ) * ( 1.0 + 5.0 / 21.0 * ( 2.0 * znd + 1.0 )** 2 )
1232                    END DO
1233                 ENDIF
1234            END DO
1235          END DO
1236       ENDIF
1237
1238       IF(ln_dia_osm) THEN
1239          IF ( iom_use("zdtdz_pyc") ) CALL iom_put( "zdtdz_pyc", wmask*zdtdz_pyc )
1240          IF ( iom_use("zdsdz_pyc") ) CALL iom_put( "zdsdz_pyc", wmask*zdsdz_pyc )
1241          IF ( iom_use("zdbdz_pyc") ) CALL iom_put( "zdbdz_pyc", wmask*zdbdz_pyc )
1242       END IF
1243
1244
1245       ! Lateral boundary conditions on zvicos (sign unchanged), needed to caclulate viscosities on u and v grids
1246       !CALL lbc_lnk( zviscos(:,:,:), 'W', 1. )
1247
1248       ! GN 25/8: need to change tmask --> wmask
1249
1250     DO jk = 2, jpkm1
1251         DO jj = 2, jpjm1
1252             DO ji = 2, jpim1
1253                p_avt(ji,jj,jk) = MAX( zdiffut(ji,jj,jk), avtb(jk) ) * tmask(ji,jj,jk)
1254                p_avm(ji,jj,jk) = MAX( zviscos(ji,jj,jk), avmb(jk) ) * tmask(ji,jj,jk)
1255             END DO
1256         END DO
1257     END DO
1258      ! Lateral boundary conditions on ghamu and ghamv, currently on W-grid  (sign unchanged), needed to caclulate gham[uv] on u and v grids
1259     CALL lbc_lnk_multi( 'zdfosm', p_avt, 'W', 1. , p_avm, 'W', 1.,   &
1260      &                  ghamu, 'W', 1. , ghamv, 'W', 1. )
1261       DO jk = 2, jpkm1
1262           DO jj = 2, jpjm1
1263               DO ji = 2, jpim1
1264                  ghamu(ji,jj,jk) = ( ghamu(ji,jj,jk) + ghamu(ji+1,jj,jk) ) &
1265                     &  / MAX( 1., tmask(ji,jj,jk) + tmask (ji + 1,jj,jk) ) * umask(ji,jj,jk)
1266
1267                  ghamv(ji,jj,jk) = ( ghamv(ji,jj,jk) + ghamv(ji,jj+1,jk) ) &
1268                      &  / MAX( 1., tmask(ji,jj,jk) + tmask (ji,jj+1,jk) ) * vmask(ji,jj,jk)
1269
1270                  ghamt(ji,jj,jk) =  ghamt(ji,jj,jk) * tmask(ji,jj,jk)
1271                  ghams(ji,jj,jk) =  ghams(ji,jj,jk) * tmask(ji,jj,jk)
1272               END DO
1273           END DO
1274        END DO
1275        ! Lateral boundary conditions on final outputs for hbl,  on T-grid (sign unchanged)
1276        CALL lbc_lnk_multi( 'zdfosm', hbl, 'T', 1., dh, 'T', 1., hmle, 'T', 1. )
1277        ! Lateral boundary conditions on final outputs for gham[ts],  on W-grid  (sign unchanged)
1278        ! Lateral boundary conditions on final outputs for gham[uv],  on [UV]-grid  (sign unchanged)
1279        CALL lbc_lnk_multi( 'zdfosm', ghamt, 'W', 1. , ghams, 'W', 1.,   &
1280         &                  ghamu, 'U', -1. , ghamv, 'V', -1. )
1281
1282      IF(ln_dia_osm) THEN
1283         SELECT CASE (nn_osm_wave)
1284         ! Stokes drift set by assumimg onstant La#=0.3(=0)  or Pierson-Moskovitz spectrum (=1).
1285         CASE(0:1)
1286            IF ( iom_use("us_x") ) CALL iom_put( "us_x", tmask(:,:,1)*zustke*zcos_wind )   ! x surface Stokes drift
1287            IF ( iom_use("us_y") ) CALL iom_put( "us_y", tmask(:,:,1)*zustke*zsin_wind )  ! y surface Stokes drift
1288            IF ( iom_use("wind_wave_abs_power") ) CALL iom_put( "wind_wave_abs_power", 1000.*rau0*tmask(:,:,1)*zustar**2*zustke )
1289         ! Stokes drift read in from sbcwave  (=2).
1290         CASE(2:3)
1291            IF ( iom_use("us_x") ) CALL iom_put( "us_x", ut0sd*umask(:,:,1) )               ! x surface Stokes drift
1292            IF ( iom_use("us_y") ) CALL iom_put( "us_y", vt0sd*vmask(:,:,1) )               ! y surface Stokes drift
1293            IF ( iom_use("wmp") ) CALL iom_put( "wmp", wmp*tmask(:,:,1) )                   ! wave mean period
1294            IF ( iom_use("hsw") ) CALL iom_put( "hsw", hsw*tmask(:,:,1) )                   ! significant wave height
1295            IF ( iom_use("wmp_NP") ) CALL iom_put( "wmp_NP", (2.*rpi*1.026/(0.877*grav) )*wndm*tmask(:,:,1) )                  ! wave mean period from NP spectrum
1296            IF ( iom_use("hsw_NP") ) CALL iom_put( "hsw_NP", (0.22/grav)*wndm**2*tmask(:,:,1) )                   ! significant wave height from NP spectrum
1297            IF ( iom_use("wndm") ) CALL iom_put( "wndm", wndm*tmask(:,:,1) )                   ! U_10
1298            IF ( iom_use("wind_wave_abs_power") ) CALL iom_put( "wind_wave_abs_power", 1000.*rau0*tmask(:,:,1)*zustar**2* &
1299                 & SQRT(ut0sd**2 + vt0sd**2 ) )
1300         END SELECT
1301         IF ( iom_use("ghamt") ) CALL iom_put( "ghamt", tmask*ghamt )            ! <Tw_NL>
1302         IF ( iom_use("ghams") ) CALL iom_put( "ghams", tmask*ghams )            ! <Sw_NL>
1303         IF ( iom_use("ghamu") ) CALL iom_put( "ghamu", umask*ghamu )            ! <uw_NL>
1304         IF ( iom_use("ghamv") ) CALL iom_put( "ghamv", vmask*ghamv )            ! <vw_NL>
1305         IF ( iom_use("zwth0") ) CALL iom_put( "zwth0", tmask(:,:,1)*zwth0 )            ! <Tw_0>
1306         IF ( iom_use("zws0") ) CALL iom_put( "zws0", tmask(:,:,1)*zws0 )                ! <Sw_0>
1307         IF ( iom_use("zwb0") ) CALL iom_put( "zwb0", tmask(:,:,1)*zwb0 )                ! <Sw_0>
1308         IF ( iom_use("zwbav") ) CALL iom_put( "zwbav", tmask(:,:,1)*zwthav )         ! upward BL-avged turb buoyancy flux
1309         IF ( iom_use("hbl") ) CALL iom_put( "hbl", tmask(:,:,1)*hbl )                  ! boundary-layer depth
1310         IF ( iom_use("ibld") ) CALL iom_put( "ibld", tmask(:,:,1)*ibld )               ! boundary-layer max k
1311         IF ( iom_use("zdt_bl") ) CALL iom_put( "zdt_bl", tmask(:,:,1)*zdt_bl )           ! dt at ml base
1312         IF ( iom_use("zds_bl") ) CALL iom_put( "zds_bl", tmask(:,:,1)*zds_bl )           ! ds at ml base
1313         IF ( iom_use("zdb_bl") ) CALL iom_put( "zdb_bl", tmask(:,:,1)*zdb_bl )           ! db at ml base
1314         IF ( iom_use("zdu_bl") ) CALL iom_put( "zdu_bl", tmask(:,:,1)*zdu_bl )           ! du at ml base
1315         IF ( iom_use("zdv_bl") ) CALL iom_put( "zdv_bl", tmask(:,:,1)*zdv_bl )           ! dv at ml base
1316         IF ( iom_use("dh") ) CALL iom_put( "dh", tmask(:,:,1)*dh )               ! Initial boundary-layer depth
1317         IF ( iom_use("hml") ) CALL iom_put( "hml", tmask(:,:,1)*hml )               ! Initial boundary-layer depth
1318         IF ( iom_use("zdt_ml") ) CALL iom_put( "zdt_ml", tmask(:,:,1)*zdt_ml )           ! dt at ml base
1319         IF ( iom_use("zds_ml") ) CALL iom_put( "zds_ml", tmask(:,:,1)*zds_ml )           ! ds at ml base
1320         IF ( iom_use("zdb_ml") ) CALL iom_put( "zdb_ml", tmask(:,:,1)*zdb_ml )           ! db at ml base
1321         IF ( iom_use("dstokes") ) CALL iom_put( "dstokes", tmask(:,:,1)*dstokes )      ! Stokes drift penetration depth
1322         IF ( iom_use("zustke") ) CALL iom_put( "zustke", tmask(:,:,1)*zustke )            ! Stokes drift magnitude at T-points
1323         IF ( iom_use("zwstrc") ) CALL iom_put( "zwstrc", tmask(:,:,1)*zwstrc )         ! convective velocity scale
1324         IF ( iom_use("zwstrl") ) CALL iom_put( "zwstrl", tmask(:,:,1)*zwstrl )         ! Langmuir velocity scale
1325         IF ( iom_use("zustar") ) CALL iom_put( "zustar", tmask(:,:,1)*zustar )         ! friction velocity scale
1326         IF ( iom_use("zvstr") ) CALL iom_put( "zvstr", tmask(:,:,1)*zvstr )         ! mixed velocity scale
1327         IF ( iom_use("zla") ) CALL iom_put( "zla", tmask(:,:,1)*zla )         ! langmuir #
1328         IF ( iom_use("wind_power") ) CALL iom_put( "wind_power", 1000.*rau0*tmask(:,:,1)*zustar**3 ) ! BL depth internal to zdf_osm routine
1329         IF ( iom_use("wind_wave_power") ) CALL iom_put( "wind_wave_power", 1000.*rau0*tmask(:,:,1)*zustar**2*zustke )
1330         IF ( iom_use("zhbl") ) CALL iom_put( "zhbl", tmask(:,:,1)*zhbl )               ! BL depth internal to zdf_osm routine
1331         IF ( iom_use("zhml") ) CALL iom_put( "zhml", tmask(:,:,1)*zhml )               ! ML depth internal to zdf_osm routine
1332         IF ( iom_use("imld") ) CALL iom_put( "imld", tmask(:,:,1)*imld )               ! index for ML depth internal to zdf_osm routine
1333         IF ( iom_use("jp_ext") ) CALL iom_put( "jp_ext", tmask(:,:,1)*jp_ext )         ! =1 if pycnocline resolved internal to zdf_osm routine
1334         IF ( iom_use("j_ddh") ) CALL iom_put( "j_ddh", tmask(:,:,1)*j_ddh )            ! index forpyc thicknessh internal to zdf_osm routine
1335         IF ( iom_use("zshear") ) CALL iom_put( "zshear", tmask(:,:,1)*zshear )         ! shear production of TKE internal to zdf_osm routine
1336         IF ( iom_use("zdh") ) CALL iom_put( "zdh", tmask(:,:,1)*zdh )                  ! pyc thicknessh internal to zdf_osm routine
1337         IF ( iom_use("zhol") ) CALL iom_put( "zhol", tmask(:,:,1)*zhol )               ! ML depth internal to zdf_osm routine
1338         IF ( iom_use("zwth_ent") ) CALL iom_put( "zwth_ent", tmask(:,:,1)*zwth_ent )   ! upward turb temp entrainment flux
1339         IF ( iom_use("zwb_ent") ) CALL iom_put( "zwb_ent", tmask(:,:,1)*zwb_ent )      ! upward turb buoyancy entrainment flux
1340         IF ( iom_use("zws_ent") ) CALL iom_put( "zws_ent", tmask(:,:,1)*zws_ent )      ! upward turb salinity entrainment flux
1341         IF ( iom_use("zt_ml") ) CALL iom_put( "zt_ml", tmask(:,:,1)*zt_ml )            ! average T in ML
1342
1343         IF ( iom_use("hmle") ) CALL iom_put( "hmle", tmask(:,:,1)*hmle )               ! FK layer depth
1344         IF ( iom_use("zmld") ) CALL iom_put( "zmld", tmask(:,:,1)*zmld )               ! FK target layer depth
1345         IF ( iom_use("zwb_fk") ) CALL iom_put( "zwb_fk", tmask(:,:,1)*zwb_fk )         ! FK b flux
1346         IF ( iom_use("zwb_fk_b") ) CALL iom_put( "zwb_fk_b", tmask(:,:,1)*zwb_fk_b )   ! FK b flux averaged over ML
1347         IF ( iom_use("mld_prof") ) CALL iom_put( "mld_prof", tmask(:,:,1)*mld_prof )! FK layer max k
1348         IF ( iom_use("zdtdx") ) CALL iom_put( "zdtdx", umask(:,:,1)*zdtdx )            ! FK dtdx at u-pt
1349         IF ( iom_use("zdtdy") ) CALL iom_put( "zdtdy", vmask(:,:,1)*zdtdy )            ! FK dtdy at v-pt
1350         IF ( iom_use("zdsdx") ) CALL iom_put( "zdsdx", umask(:,:,1)*zdsdx )            ! FK dtdx at u-pt
1351         IF ( iom_use("zdsdy") ) CALL iom_put( "zdsdy", vmask(:,:,1)*zdsdy )            ! FK dsdy at v-pt
1352         IF ( iom_use("dbdx_mle") ) CALL iom_put( "dbdx_mle", umask(:,:,1)*dbdx_mle )            ! FK dbdx at u-pt
1353         IF ( iom_use("dbdy_mle") ) CALL iom_put( "dbdy_mle", vmask(:,:,1)*dbdy_mle )            ! FK dbdy at v-pt
1354         IF ( iom_use("zdiff_mle") ) CALL iom_put( "zdiff_mle", tmask(:,:,1)*zdiff_mle )! FK diff in MLE at t-pt
1355         IF ( iom_use("zvel_mle") ) CALL iom_put( "zvel_mle", tmask(:,:,1)*zdiff_mle )! FK diff in MLE at t-pt
1356
1357      END IF
1358
1359CONTAINS
1360! subroutine code changed, needs syntax checking.
1361  SUBROUTINE zdf_osm_diffusivity_viscosity( zdiffut, zviscos )
1362
1363!!---------------------------------------------------------------------
1364     !!                   ***  ROUTINE zdf_osm_diffusivity_viscosity  ***
1365     !!
1366     !! ** Purpose : Determines the eddy diffusivity and eddy viscosity profiles in the mixed layer and the pycnocline.
1367     !!
1368     !! ** Method  :
1369     !!
1370     !! !!----------------------------------------------------------------------
1371     REAL(wp), DIMENSION(:,:,:) :: zdiffut
1372     REAL(wp), DIMENSION(:,:,:) :: zviscos
1373! local
1374
1375! Scales used to calculate eddy diffusivity and viscosity profiles
1376      REAL(wp), DIMENSION(jpi,jpj) :: zdifml_sc, zvisml_sc
1377      REAL(wp), DIMENSION(jpi,jpj) :: zdifpyc_n_sc, zdifpyc_s_sc, zdifpyc_shr
1378      REAL(wp), DIMENSION(jpi,jpj) :: zvispyc_n_sc, zvispyc_s_sc,zvispyc_shr
1379      REAL(wp), DIMENSION(jpi,jpj) :: zbeta_d_sc, zbeta_v_sc
1380!
1381      REAL(wp) :: zvel_sc_pyc, zvel_sc_ml, zstab_fac
1382     
1383      REAL(wp), PARAMETER :: rn_dif_ml = 0.8, rn_vis_ml = 0.375
1384      REAL(wp), PARAMETER :: rn_dif_pyc = 0.15, rn_vis_pyc = 0.142
1385      REAL(wp), PARAMETER :: rn_vispyc_shr = 0.15
1386     
1387      DO jj = 2, jpjm1
1388          DO ji = 2, jpim1
1389             IF ( lconv(ji,jj) ) THEN
1390             
1391               zvel_sc_pyc = ( 0.15 * zvstr(ji,jj)**3 + zwstrc(ji,jj)**3 + 4.25 * zshear(ji,jj) * zhbl(ji,jj) )**pthird
1392               zvel_sc_ml = ( zvstr(ji,jj)**3 + 0.5 * zwstrc(ji,jj)**3 )**pthird
1393               zstab_fac = ( zhml(ji,jj) / zvel_sc_ml * ( 1.4 - 0.4 / ( 1.0 + EXP(-3.5 * LOG10(-zhol(ji,jj) ) ) )**1.25 ) )**2
1394
1395               zdifml_sc(ji,jj) = rn_dif_ml * zhml(ji,jj) * zvel_sc_ml
1396               zvisml_sc(ji,jj) = rn_vis_ml * zdifml_sc(ji,jj)
1397
1398               IF ( lpyc(ji,jj) ) THEN
1399                 zdifpyc_n_sc(ji,jj) =  rn_dif_pyc * zvel_sc_ml * zdh(ji,jj)
1400
1401                 IF ( lshear(ji,jj) .AND. j_ddh(ji,jj) /= 2 ) THEN
1402                   zdifpyc_n_sc(ji,jj) = zdifpyc_n_sc(ji,jj) + rn_vispyc_shr * ( zshear(ji,jj) * zhbl(ji,jj) )**pthird * zhbl(ji,jj)
1403                 ENDIF
1404               
1405                 zdifpyc_s_sc(ji,jj) = zwb_ent(ji,jj) + 0.0025 * zvel_sc_pyc * ( zhbl(ji,jj) / zdh(ji,jj) - 1.0 ) * ( zb_ml(ji,jj) - zb_bl(ji,jj) )
1406                 zdifpyc_s_sc(ji,jj) = 0.09 * zdifpyc_s_sc(ji,jj) * zstab_fac
1407                 zdifpyc_s_sc(ji,jj) = MAX( zdifpyc_s_sc(ji,jj), -0.5 * zdifpyc_n_sc(ji,jj) )
1408                 
1409                 zvispyc_n_sc(ji,jj) = 0.09 * zvel_sc_pyc * ( 1.0 - zhbl(ji,jj) / zdh(ji,jj) )**2 * ( 0.005 * ( zu_ml(ji,jj)-zu_bl(ji,jj) )**2 + 0.0075 * ( zv_ml(ji,jj)-zv_bl(ji,jj) )**2 ) / zdh(ji,jj)
1410                 zvispyc_n_sc(ji,jj) = rn_vis_pyc * zvel_sc_ml * zdh(ji,jj) + zvispyc_n_sc(ji,jj) * zstab_fac
1411                 IF ( lshear(ji,jj) .AND. j_ddh(ji,jj) /= 2 ) THEN
1412                   zvispyc_n_sc(ji,jj) = zvispyc_n_sc(ji,jj) + rn_vispyc_shr * ( zshear(ji,jj) * zhbl(ji,jj ) )**pthird * zhbl(ji,jj)
1413                 ENDIF
1414               
1415                 zvispyc_s_sc(ji,jj) = 0.09 * ( zwb_min(ji,jj) + 0.0025 * zvel_sc_pyc * ( zhbl(ji,jj) / zdh(ji,jj) - 1.0 ) * ( zb_ml(ji,jj) - zb_bl(ji,jj) ) )
1416                 zvispyc_s_sc(ji,jj) = zvispyc_s_sc(ji,jj) * zstab_fac
1417                 zvispyc_s_sc(ji,jj) = MAX( zvispyc_s_sc(ji,jj), -0.5 * zvispyc_n_sc(ji,jj) )
1418
1419                 zbeta_d_sc(ji,jj) = 1.0 - ( ( zdifpyc_n_sc(ji,jj) + 1.4 * zdifpyc_s_sc(ji,jj) ) / ( zdifml_sc(ji,jj) + epsln ) )**p2third
1420                 zbeta_v_sc(ji,jj) = 1.0 -  2.0 * ( zvispyc_n_sc(ji,jj) + zvispyc_s_sc(ji,jj) ) / ( zvisml_sc(ji,jj) + epsln )
1421               ELSE
1422                 zbeta_d_sc(ji,jj) = 1.0
1423                 zbeta_v_sc(ji,jj) = 1.0
1424               ENDIF
1425             ELSE
1426               zdifml_sc(ji,jj) = zvstr(ji,jj) * zhbl(ji,jj) * MAX( EXP ( -( zhol(ji,jj) / 0.6_wp )**2 ), 0.2_wp)
1427               zvisml_sc(ji,jj) = zvstr(ji,jj) * zhbl(ji,jj) * MAX( EXP ( -( zhol(ji,jj) / 0.6_wp )**2 ), 0.2_wp)
1428             END IF
1429          END DO
1430       END DO
1431!
1432       DO jj = 2, jpjm1
1433          DO ji = 2, jpim1
1434             IF ( lconv(ji,jj) ) THEN
1435                DO jk = 2, imld(ji,jj)   ! mixed layer diffusivity
1436                    zznd_ml = gdepw_n(ji,jj,jk) / zhml(ji,jj)
1437                    !
1438                    zdiffut(ji,jj,jk) = zdifml_sc(ji,jj) * zznd_ml * ( 1.0 - zbeta_d_sc(ji,jj) * zznd_ml )**1.5
1439                    !
1440                    zviscos(ji,jj,jk) = zvisml_sc(ji,jj) * zznd_ml * ( 1.0 - zbeta_v_sc(ji,jj) * zznd_ml ) &
1441      &            *                                      ( 1.0 - 0.5 * zznd_ml**2 )
1442                END DO
1443! pycnocline
1444                IF ( lpyc(ji,jj) ) THEN
1445! Diffusivity profile in the pycnocline given by cubic polynomial.
1446                   za_cubic = 0.5
1447                   zb_cubic = -1.75 * zdifpyc_s_sc(ji,jj) / zdifpyc_n_sc(ji,jj)
1448                   zd_cubic = ( zdh(ji,jj) * zdifml_sc(ji,jj) / zhml(ji,jj) * SQRT( 1.0 - zbeta_d_sc(ji,jj) ) * ( 2.5 * zbeta_d_sc(ji,jj) - 1.0 ) &
1449                        & - 0.85 * zdifpyc_s_sc(ji,jj) ) / MAX(zdifpyc_n_sc(ji,jj), 1.e-8)
1450                   zd_cubic = zd_cubic - zb_cubic - 2.0 * ( 1.0 - za_cubic  - zb_cubic )
1451                   zc_cubic = 1.0 - za_cubic - zb_cubic - zd_cubic
1452                   DO jk = imld(ji,jj) , ibld(ji,jj)
1453                     zznd_pyc = -( gdepw_n(ji,jj,jk) - zhbl(ji,jj) ) / MAX(zdh(ji,jj), 1.e-6)
1454                         !
1455                     zdiffut(ji,jj,jk) = zdifpyc_n_sc(ji,jj) * ( za_cubic + zb_cubic * zznd_pyc + zc_cubic * zznd_pyc**2 +   zd_cubic * zznd_pyc**3 )
1456
1457                     zdiffut(ji,jj,jk) = zdiffut(ji,jj,jk) + zdifpyc_s_sc(ji,jj) * ( 1.75 * zznd_pyc - 0.15 * zznd_pyc**2 - 0.2 * zznd_pyc**3 )
1458                   END DO
1459 ! viscosity profiles.
1460                   za_cubic = 0.5
1461                   zb_cubic = -1.75 * zvispyc_s_sc(ji,jj) / zvispyc_n_sc(ji,jj)
1462                   zd_cubic = ( 0.5 * zvisml_sc(ji,jj) * zdh(ji,jj) / zhml(ji,jj) - 0.85 * zvispyc_s_sc(ji,jj)  )  / MAX(zvispyc_n_sc(ji,jj), 1.e-8)
1463                   zd_cubic = zd_cubic - zb_cubic - 2.0 * ( 1.0 - za_cubic - zb_cubic )
1464                   zc_cubic = 1.0 - za_cubic - zb_cubic - zd_cubic
1465                   DO jk = imld(ji,jj) , ibld(ji,jj)
1466                      zznd_pyc = -( gdepw_n(ji,jj,jk) - zhbl(ji,jj) ) / MAX(zdh(ji,jj), 1.e-6)
1467                       zviscos(ji,jj,jk) = zvispyc_n_sc(ji,jj) * ( za_cubic + zb_cubic * zznd_pyc + zc_cubic * zznd_pyc**2 + zd_cubic * zznd_pyc**3 )
1468                       zviscos(ji,jj,jk) = zviscos(ji,jj,jk) + zvispyc_s_sc(ji,jj) * ( 1.75 * zznd_pyc - 0.15 * zznd_pyc**2 -0.2 * zznd_pyc**3 )
1469                   END DO
1470                   IF ( zdhdt(ji,jj) > 0._wp ) THEN
1471                    zdiffut(ji,jj,ibld(ji,jj)+1) = MAX( 0.5 * zdhdt(ji,jj) * e3w_n(ji,jj,ibld(ji,jj)+1), 1.0e-6 )
1472                    zviscos(ji,jj,ibld(ji,jj)+1) = MAX( 0.5 * zdhdt(ji,jj) * e3w_n(ji,jj,ibld(ji,jj)+1), 1.0e-6 )
1473                   ELSE
1474                     zdiffut(ji,jj,ibld(ji,jj)) = 0._wp
1475                     zviscos(ji,jj,ibld(ji,jj)) = 0._wp
1476                   ENDIF
1477                ENDIF
1478             ELSE
1479             ! stable conditions
1480                DO jk = 2, ibld(ji,jj)
1481                   zznd_ml = gdepw_n(ji,jj,jk) / zhbl(ji,jj)
1482                   zdiffut(ji,jj,jk) = 0.75 * zdifml_sc(ji,jj) * zznd_ml * ( 1.0 - zznd_ml )**1.5
1483                   zviscos(ji,jj,jk) = 0.375 * zvisml_sc(ji,jj) * zznd_ml * (1.0 - zznd_ml) * ( 1.0 - zznd_ml**2 )
1484                END DO
1485
1486                IF ( zdhdt(ji,jj) > 0._wp ) THEN
1487                   zdiffut(ji,jj,ibld(ji,jj)) = MAX(zdhdt(ji,jj), 1.0e-6) * e3w_n(ji, jj, ibld(ji,jj))
1488                   zviscos(ji,jj,ibld(ji,jj)) = MAX(zdhdt(ji,jj), 1.0e-6) * e3w_n(ji, jj, ibld(ji,jj))
1489                ENDIF
1490             ENDIF   ! end if ( lconv )
1491             !
1492          END DO  ! end of ji loop
1493       END DO     ! end of jj loop
1494       
1495  END SUBROUTINE zdf_osm_diffusivity_viscosity
1496 
1497  SUBROUTINE zdf_osm_osbl_state( lconv, lshear, j_ddh, zwb_ent, zwb_min, zshear )
1498
1499!!---------------------------------------------------------------------
1500     !!                   ***  ROUTINE zdf_osm_osbl_state  ***
1501     !!
1502     !! ** Purpose : Determines the state of the OSBL, stable/unstable, shear/ noshear. Also determines shear production, entrainment buoyancy flux and interfacial Richardson number
1503     !!
1504     !! ** Method  :
1505     !!
1506     !! !!----------------------------------------------------------------------
1507
1508     INTEGER, DIMENSION(jpi,jpj) :: j_ddh  ! j_ddh = 0, active shear layer; j_ddh=1, shear layer not active; j_ddh=2 shear production low.
1509     
1510     LOGICAL, DIMENSION(jpi,jpj) :: lconv, lshear
1511
1512     REAL(wp), DIMENSION(jpi,jpj) :: zwb_ent, zwb_min ! Buoyancy fluxes at base of well-mixed layer.
1513     REAL(wp), DIMENSION(jpi,jpj) :: zshear  ! production of TKE due to shear across the pycnocline
1514
1515! Local Variables
1516
1517     INTEGER :: jj, ji
1518     
1519     REAL(wp), DIMENSION(jpi,jpj) :: zekman
1520     REAL(wp), DIMENSION(jpi,jpj) :: zri_p, zri_b   ! Richardson numbers
1521     REAL(wp) :: zshear_u, zshear_v, zwb_shr
1522     REAL(wp) :: zwcor, zrf_conv, zrf_shear, zrf_langmuir, zr_stokes
1523
1524     REAL, PARAMETER :: za_shr = 0.4, zb_shr = 6.5, za_wb_s = 0.8
1525     REAL, PARAMETER :: zalpha_c = 0.2, zalpha_lc = 0.03     
1526     REAL, PARAMETER :: zalpha_ls = 0.06, zalpha_s = 0.15
1527     REAL, PARAMETER :: rn_ri_p_thresh = 27.0
1528     REAL, PARAMETER :: zri_c = 0.25
1529     REAL, PARAMETER :: zek = 4.0
1530     REAL, PARAMETER :: zrot=0._wp  ! dummy rotation rate of surface stress.
1531     
1532! Determins stability and set flag lconv
1533     DO jj = 2, jpjm1
1534       DO ji = 2, jpim1
1535         IF ( zhol(ji,jj) < 0._wp ) THEN
1536            lconv(ji,jj) = .TRUE.
1537          ELSE
1538             lconv(ji,jj) = .FALSE.
1539          ENDIF
1540       END DO
1541     END DO       
1542 
1543     zekman(:,:) = EXP( - zek * ABS( ff_t(:,:) ) * zhbl(:,:) / MAX(zustar(:,:), 1.e-8 ) )
1544     
1545     zshear(:,:) = 0._wp
1546     j_ddh(:,:) = 1     
1547 
1548     DO jj = 2, jpjm1
1549       DO ji = 2, jpim1
1550         IF ( lconv(ji,jj) ) THEN
1551            IF ( zdb_bl(ji,jj) > 0._wp ) THEN
1552              zri_p(ji,jj) = MAX (  SQRT( zdb_bl(ji,jj) * zdh(ji,jj) / MAX( zdu_bl(ji,jj)**2 + zdv_bl(ji,jj)**2, 1.e-8) )  *  ( zhbl(ji,jj) / zdh(ji,jj) ) * ( zvstr(ji,jj) / MAX( zustar(ji,jj), 1.e-6 ) )**2 &
1553                   & / MAX( zekman(ji,jj), 1.e-6 )  , 5._wp )
1554           
1555              IF ( ff_t(ji,jj) >= 0._wp ) THEN
1556                 !  Northern Hemisphere
1557                 zri_b(ji,jj) = zdb_ml(ji,jj) * zdh(ji,jj) / ( MAX( zdu_ml(ji,jj), 1.e-5 )**2 + MAX( -zdv_ml(ji,jj), 1.e-5)**2 )
1558              ELSE
1559                  !  Southern Hemisphere
1560                 zri_b(ji,jj) = zdb_ml(ji,jj) * zdh(ji,jj) / ( MAX( zdu_ml(ji,jj), 1.e-5 )**2 + MAX( zdv_ml(ji,jj), 1.e-5)**2 )
1561              ENDIF
1562              zshear(ji,jj) = za_shr * zekman(ji,jj) * ( MAX( zustar(ji,jj)**2 * zdu_ml(ji,jj) / zhbl(ji,jj), 0._wp ) + zb_shr * MAX( -ff_t(ji,jj) * zustke(ji,jj) * dstokes(ji,jj) * zdv_ml(ji,jj) / zhbl(ji,jj), 0._wp ) )
1563! Stability Dependence
1564              zshear(ji,jj) = zshear(ji,jj) * EXP( -0.75 * MAX(0._wp,( zri_b(ji,jj) - zri_c ) / zri_c ) )
1565!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1566! Test ensures j_ddh=0 is not selected. Change to zri_p<27 when  !
1567! full code available                                          !
1568!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1569              IF ( zshear(ji,jj) > 1.e-10 ) THEN
1570                IF ( zri_p(ji,jj) < rn_ri_p_thresh ) THEN
1571! Growing shear layer
1572                  j_ddh(ji,jj) = 0
1573                  lshear(ji,jj) = .TRUE.
1574                ELSE
1575                  j_ddh(ji,jj) = 1
1576!                IF ( zri_b <= 1.5 .and. zshear(ji,jj) > 0._wp ) THEN
1577! shear production large enough to determine layer charcteristics, but can't maintain a shear layer.
1578                lshear(ji,jj) = .TRUE.
1579!                ELSE
1580                ENDIF
1581              ELSE
1582                j_ddh(ji,jj) = 2
1583                lshear(ji,jj) = .FALSE.
1584              ENDIF
1585! Shear production may not be zero, but is small and doesn't determine characteristics of pycnocline.
1586!                  zshear(ji,jj) = 0.5 * zshear(ji,jj)
1587!                  lshear(ji,jj) = .FALSE.
1588!                ENDIF               
1589            ELSE                ! zdb_bl test, note zshear set to zero
1590              j_ddh(ji,jj) = 2
1591              lshear(ji,jj) = .FALSE.
1592            ENDIF
1593          ENDIF
1594       END DO
1595     END DO
1596 
1597! Calculate entrainment buoyancy flux due to surface fluxes.
1598
1599     DO jj = 2, jpjm1
1600       DO ji = 2, jpim1
1601         IF ( lconv(ji,jj) ) THEN
1602           zwcor = ABS(ff_t(ji,jj)) * zhbl(ji,jj) + epsln
1603           zrf_conv = TANH( ( zwstrc(ji,jj) / zwcor )**0.69 )
1604           zrf_shear = TANH( ( zustar(ji,jj) / zwcor )**0.69 )
1605           zrf_langmuir = TANH( ( zwstrl(ji,jj) / zwcor )**0.69 )
1606           IF (nn_osm_SD_reduce > 0 ) THEN
1607           ! Effective Stokes drift already reduced from surface value
1608              zr_stokes = 1.0_wp
1609           ELSE
1610            ! Effective Stokes drift only reduced by factor rn_zdfodm_adjust_sd,
1611             ! requires further reduction where BL is deep
1612              zr_stokes = 1.0 - EXP( -25.0 * dstokes(ji,jj) / hbl(ji,jj) &
1613            &                  * ( 1.0 + 4.0 * dstokes(ji,jj) / hbl(ji,jj) ) )
1614           END IF
1615           zwb_ent(ji,jj) = - 2.0 * zalpha_c * zrf_conv * zwbav(ji,jj) &
1616                  &                  - zalpha_s * zrf_shear * zustar(ji,jj)**3 /zhml(ji,jj) &
1617                  &         + zr_stokes * ( zalpha_s * EXP( -1.5 * zla(ji,jj) ) * zrf_shear * zustar(ji,jj)**3 &
1618                  &                                         - zrf_langmuir * zalpha_lc * zwstrl(ji,jj)**3 ) / zhml(ji,jj)
1619             !
1620         ENDIF
1621       END DO ! ji loop
1622     END DO   ! jj loop
1623
1624     zwb_min(:,:) = 0._wp
1625
1626     DO jj = 2, jpjm1
1627       DO ji = 2, jpim1
1628         IF ( lshear(ji,jj) ) THEN
1629           IF ( lconv(ji,jj) ) THEN
1630! Unstable OSBL
1631              zwb_shr = -za_wb_s * zri_b(ji,jj) * zshear(ji,jj)
1632              IF ( j_ddh(ji,jj) == 0 ) THEN
1633
1634! ! Developing shear layer, additional shear production possible.
1635
1636!                 zshear_u = MAX( zustar(ji,jj)**2 * MAX( zdu_ml(ji,jj), 0._wp ) /  zhbl(ji,jj), 0._wp )
1637!                 zshear(ji,jj) = zshear(ji,jj) + zshear_u * ( 1.0 - MIN( zri_p(ji,jj) / rn_ri_p_thresh, 1.d0 )**2 )
1638!                 zshear(ji,jj) = MIN( zshear(ji,jj), zshear_u )
1639               
1640!                 zwb_shr = zwb_shr - 0.25 * MAX ( zshear_u, 0._wp) * ( 1.0 - MIN( zri_p(ji,jj) / rn_ri_p_thresh, 1._wp )**2 )
1641!                 zwb_shr = MAX( zwb_shr, -0.25 * zshear_u )
1642               
1643              ENDIF               
1644              zwb_ent(ji,jj) = zwb_ent(ji,jj) + zwb_shr
1645!              zwb_min(ji,jj) = zwb_ent(ji,jj) + zdh(ji,jj) / zhbl(ji,jj) * zwb0(ji,jj)
1646           ELSE    ! IF ( lconv ) THEN - ENDIF
1647! Stable OSBL  - shear production not coded for first attempt.           
1648           ENDIF  ! lconv
1649         ENDIF  ! lshear
1650         IF ( lconv(ji,jj) ) THEN
1651! Unstable OSBL
1652            zwb_min(ji,jj) = zwb_ent(ji,jj) + zdh(ji,jj) / zhbl(ji,jj) * zwb0tot(ji,jj)
1653         ENDIF  ! lconv
1654       END DO   ! ji
1655     END DO     ! jj
1656   END SUBROUTINE zdf_osm_osbl_state
1657     
1658     
1659   SUBROUTINE zdf_osm_vertical_average( jnlev_av, jp_ext, zt, zs, zb, zu, zv, zdt, zds, zdb, zdu, zdv )
1660     !!---------------------------------------------------------------------
1661     !!                   ***  ROUTINE zdf_vertical_average  ***
1662     !!
1663     !! ** Purpose : Determines vertical averages from surface to jnlev.
1664     !!
1665     !! ** Method  : Averages are calculated from the surface to jnlev.
1666     !!              The external level used to calculate differences is ibld+ibld_ext
1667     !!
1668     !!----------------------------------------------------------------------
1669
1670        INTEGER, DIMENSION(jpi,jpj) :: jnlev_av  ! Number of levels to average over.
1671        INTEGER, DIMENSION(jpi,jpj) :: jp_ext
1672
1673        ! Alan: do we need zb?
1674        REAL(wp), DIMENSION(jpi,jpj) :: zt, zs, zb        ! Average temperature and salinity
1675        REAL(wp), DIMENSION(jpi,jpj) :: zu,zv         ! Average current components
1676        REAL(wp), DIMENSION(jpi,jpj) :: zdt, zds, zdb ! Difference between average and value at base of OSBL
1677        REAL(wp), DIMENSION(jpi,jpj) :: zdu, zdv      ! Difference for velocity components.
1678
1679        INTEGER :: jk, ji, jj, ibld_ext
1680        REAL(wp) :: zthick, zthermal, zbeta
1681
1682
1683        zt   = 0._wp
1684        zs   = 0._wp
1685        zu   = 0._wp
1686        zv   = 0._wp
1687        DO jj = 2, jpjm1                                 !  Vertical slab
1688         DO ji = 2, jpim1
1689            zthermal = rab_n(ji,jj,1,jp_tem) !ideally use ibld not 1??
1690            zbeta    = rab_n(ji,jj,1,jp_sal)
1691               ! average over depth of boundary layer
1692            zthick = epsln
1693            DO jk = 2, jnlev_av(ji,jj)
1694               zthick = zthick + e3t_n(ji,jj,jk)
1695               zt(ji,jj)   = zt(ji,jj)  + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_tem)
1696               zs(ji,jj)   = zs(ji,jj)  + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal)
1697               zu(ji,jj)   = zu(ji,jj)  + e3t_n(ji,jj,jk) &
1698                     &            * ( ub(ji,jj,jk) + ub(ji - 1,jj,jk) ) &
1699                     &            / MAX( 1. , umask(ji,jj,jk) + umask(ji - 1,jj,jk) )
1700               zv(ji,jj)   = zv(ji,jj)  + e3t_n(ji,jj,jk) &
1701                     &            * ( vb(ji,jj,jk) + vb(ji,jj - 1,jk) ) &
1702                     &            / MAX( 1. , vmask(ji,jj,jk) + vmask(ji,jj - 1,jk) )
1703            END DO
1704            zt(ji,jj) = zt(ji,jj) / zthick
1705            zs(ji,jj) = zs(ji,jj) / zthick
1706            zu(ji,jj) = zu(ji,jj) / zthick
1707            zv(ji,jj) = zv(ji,jj) / zthick
1708            zb(ji,jj) = grav * zthermal * zt(ji,jj) - grav * zbeta * zs(ji,jj)
1709            ibld_ext = jnlev_av(ji,jj) + jp_ext(ji,jj)
1710            IF ( ibld_ext < mbkt(ji,jj) ) THEN
1711              zdt(ji,jj) = zt(ji,jj) - tsn(ji,jj,ibld_ext,jp_tem)
1712              zds(ji,jj) = zs(ji,jj) - tsn(ji,jj,ibld_ext,jp_sal)
1713              zdu(ji,jj) = zu(ji,jj) - ( ub(ji,jj,ibld_ext) + ub(ji-1,jj,ibld_ext ) ) &
1714                     &    / MAX(1. , umask(ji,jj,ibld_ext ) + umask(ji-1,jj,ibld_ext ) )
1715              zdv(ji,jj) = zv(ji,jj) - ( vb(ji,jj,ibld_ext) + vb(ji,jj-1,ibld_ext ) ) &
1716                     &   / MAX(1. , vmask(ji,jj,ibld_ext ) + vmask(ji,jj-1,ibld_ext ) )
1717              zdb(ji,jj) = grav * zthermal * zdt(ji,jj) - grav * zbeta * zds(ji,jj)
1718            ELSE
1719              zdt(ji,jj) = 0._wp
1720              zds(ji,jj) = 0._wp
1721              zdu(ji,jj) = 0._wp
1722              zdv(ji,jj) = 0._wp
1723              zdb(ji,jj) = 0._wp
1724            ENDIF
1725         END DO
1726        END DO
1727   END SUBROUTINE zdf_osm_vertical_average
1728
1729   SUBROUTINE zdf_osm_velocity_rotation( zcos_w, zsin_w, zu, zv, zdu, zdv )
1730     !!---------------------------------------------------------------------
1731     !!                   ***  ROUTINE zdf_velocity_rotation  ***
1732     !!
1733     !! ** Purpose : Rotates frame of reference of averaged velocity components.
1734     !!
1735     !! ** Method  : The velocity components are rotated into frame specified by zcos_w and zsin_w
1736     !!
1737     !!----------------------------------------------------------------------
1738
1739        REAL(wp), DIMENSION(jpi,jpj) :: zcos_w, zsin_w       ! Cos and Sin of rotation angle
1740        REAL(wp), DIMENSION(jpi,jpj) :: zu, zv               ! Components of current
1741        REAL(wp), DIMENSION(jpi,jpj) :: zdu, zdv             ! Change in velocity components across pycnocline
1742
1743        INTEGER :: ji, jj
1744        REAL(wp) :: ztemp
1745
1746        DO jj = 2, jpjm1
1747           DO ji = 2, jpim1
1748              ztemp = zu(ji,jj)
1749              zu(ji,jj) = zu(ji,jj) * zcos_w(ji,jj) + zv(ji,jj) * zsin_w(ji,jj)
1750              zv(ji,jj) = zv(ji,jj) * zcos_w(ji,jj) - ztemp * zsin_w(ji,jj)
1751              ztemp = zdu(ji,jj)
1752              zdu(ji,jj) = zdu(ji,jj) * zcos_w(ji,jj) + zdv(ji,jj) * zsin_w(ji,jj)
1753              zdv(ji,jj) = zdv(ji,jj) * zcos_w(ji,jj) - ztemp * zsin_w(ji,jj)
1754           END DO
1755        END DO
1756    END SUBROUTINE zdf_osm_velocity_rotation
1757
1758    SUBROUTINE zdf_osm_osbl_state_fk( lpyc, lflux, lmle, zwb_fk )
1759     !!---------------------------------------------------------------------
1760     !!                   ***  ROUTINE zdf_osm_osbl_state_fk  ***
1761     !!
1762     !! ** Purpose : Determines the state of the OSBL and MLE layer. Info is returned in the logicals lpyc,lflux and lmle. Used with Fox-Kemper scheme.
1763     !!  lpyc :: determines whether pycnocline flux-grad relationship needs to be determined
1764     !!  lflux :: determines whether effects of surface flux extend below the base of the OSBL
1765     !!  lmle  :: determines whether the layer with MLE is increasing with time or if base is relaxing towards hbl.
1766     !!
1767     !! ** Method  :
1768     !!
1769     !!
1770     !!----------------------------------------------------------------------
1771     
1772! Outputs
1773      LOGICAL,  DIMENSION(jpi,jpj)  :: lpyc, lflux, lmle
1774      REAL(wp), DIMENSION(jpi,jpj)  :: zwb_fk
1775!
1776      REAL(wp), DIMENSION(jpi,jpj)  :: znd_param
1777      REAL(wp)                      :: zbuoy, ztmp, zpe_mle_layer
1778      REAL(wp)                      :: zpe_mle_ref, zdbdz_mle_int
1779     
1780      znd_param(:,:) = 0._wp
1781
1782        DO jj = 2, jpjm1
1783          DO ji = 2, jpim1         
1784             ztmp =  r1_ft(ji,jj) *  MIN( 111.e3_wp , e1u(ji,jj) ) / rn_osm_mle_lf
1785             zwb_fk(ji,jj) = rn_osm_mle_ce * hmle(ji,jj) * hmle(ji,jj) * ztmp * zdbds_mle(ji,jj) * zdbds_mle(ji,jj)
1786          END DO
1787        END DO       
1788        DO jj = 2, jpjm1
1789          DO ji = 2, jpim1
1790                    !
1791            IF ( lconv(ji,jj) ) THEN
1792              IF ( zhmle(ji,jj) > 1.2 * zhbl(ji,jj) ) THEN
1793                zt_mle(ji,jj) = ( zt_mle(ji,jj) * zhmle(ji,jj) - zt_bl(ji,jj) * zhbl(ji,jj) ) / ( zhmle(ji,jj) - zhbl(ji,jj) )
1794                zs_mle(ji,jj) = ( zs_mle(ji,jj) * zhmle(ji,jj) - zs_bl(ji,jj) * zhbl(ji,jj) ) / ( zhmle(ji,jj) - zhbl(ji,jj) )
1795                zb_mle(ji,jj) = ( zb_mle(ji,jj) * zhmle(ji,jj) - zb_bl(ji,jj) * zhbl(ji,jj) ) / ( zhmle(ji,jj) - zhbl(ji,jj) )
1796                zdbdz_mle_int = ( zb_bl(ji,jj) - ( 2.0 * zb_mle(ji,jj) -zb_bl(ji,jj) ) ) / ( zhmle(ji,jj) - zhbl(ji,jj) )
1797! Calculate potential energies of actual profile and reference profile.
1798                zpe_mle_layer = 0._wp
1799                zpe_mle_ref = 0._wp
1800                zthermal = rab_n(ji,jj,1,jp_tem)
1801                zbeta    = rab_n(ji,jj,1,jp_sal)
1802
1803                DO jk = ibld(ji,jj), mld_prof(ji,jj)
1804                  zbuoy = grav * ( zthermal * tsn(ji,jj,jk,jp_tem) - zbeta * tsn(ji,jj,jk,jp_sal) )
1805                  zpe_mle_layer = zpe_mle_layer + zbuoy * gdepw_n(ji,jj,jk) * e3w_n(ji,jj,jk)
1806                  zpe_mle_ref = zpe_mle_ref + ( zb_bl(ji,jj) - zdbdz_mle_int * ( gdepw_n(ji,jj,jk) - zhbl(ji,jj) ) ) * gdepw_n(ji,jj,jk) * e3w_n(ji,jj,jk)
1807                END DO
1808! Non-dimensional parameter to diagnose the presence of thermocline
1809                   
1810                znd_param(ji,jj) = ( zpe_mle_layer - zpe_mle_ref ) * ABS( ff_t(ji,jj) ) / ( MAX( zwb_fk(ji,jj), 1.0e-10 ) * zhmle(ji,jj) )
1811              ENDIF
1812            ENDIF
1813          END DO
1814        END DO
1815
1816! Diagnosis
1817        DO jj = 2, jpjm1
1818           DO ji = 2, jpim1
1819             IF ( lconv(ji,jj) ) THEN
1820               IF ( -2.0 * zwb_fk(ji,jj) / zwb_ent(ji,jj) > 0.5 ) THEN
1821                 IF ( zhmle(ji,jj) > 1.2 * zhbl(ji,jj) ) THEN
1822! MLE layer growing
1823                   IF ( znd_param (ji,jj) > 100. ) THEN
1824! Thermocline present
1825                     lflux(ji,jj) = .FALSE.
1826                     lmle(ji,jj) =.FALSE.
1827                   ELSE
1828! Thermocline not present
1829                     lflux(ji,jj) = .TRUE.
1830                     lmle(ji,jj) = .TRUE.
1831                   ENDIF  ! znd_param > 100
1832!
1833                   IF ( zdb_bl(ji,jj) < rn_osm_bl_thresh ) THEN
1834                     lpyc(ji,jj) = .FALSE.
1835                   ELSE
1836                      lpyc(ji,jj) = .TRUE.
1837                   ENDIF
1838                 ELSE
1839! MLE layer restricted to OSBL or just below.
1840                   IF ( zdb_bl(ji,jj) < rn_osm_bl_thresh ) THEN
1841! Weak stratification MLE layer can grow.
1842                     lpyc(ji,jj) = .FALSE.
1843                     lflux(ji,jj) = .TRUE.
1844                     lmle(ji,jj) = .TRUE.
1845                   ELSE
1846! Strong stratification
1847                     lpyc(ji,jj) = .TRUE.
1848                     lflux(ji,jj) = .FALSE.
1849                     lmle(ji,jj) = .FALSE.
1850                   ENDIF ! zdb_bl < rn_mle_thresh_bl and
1851                 ENDIF  ! zhmle > 1.2 zhbl
1852               ELSE
1853                 lpyc(ji,jj) = .TRUE.
1854                 lflux(ji,jj) = .FALSE.
1855                 lmle(ji,jj) = .FALSE.
1856                 IF ( zdb_bl(ji,jj) < rn_osm_bl_thresh ) lpyc(ji,jj) = .FALSE.
1857               ENDIF !  -2.0 * zwb_fk(ji,jj) / zwb_ent > 0.5
1858             ELSE
1859! Stable Boundary Layer
1860               lpyc(ji,jj) = .FALSE.
1861               lflux(ji,jj) = .FALSE.
1862               lmle(ji,jj) = .FALSE.
1863             ENDIF  ! lconv
1864           END DO
1865        END DO
1866    END SUBROUTINE zdf_osm_osbl_state_fk
1867
1868    SUBROUTINE zdf_osm_external_gradients(jbase, zdtdz, zdsdz, zdbdz )
1869     !!---------------------------------------------------------------------
1870     !!                   ***  ROUTINE zdf_osm_external_gradients  ***
1871     !!
1872     !! ** Purpose : Calculates the gradients below the OSBL
1873     !!
1874     !! ** Method  : Uses ibld and ibld_ext to determine levels to calculate the gradient.
1875     !!
1876     !!----------------------------------------------------------------------
1877
1878     INTEGER, DIMENSION(jpi,jpj)  :: jbase
1879     REAL(wp), DIMENSION(jpi,jpj) :: zdtdz, zdsdz, zdbdz   ! External gradients of temperature, salinity and buoyancy.
1880
1881     INTEGER :: jj, ji, jkb, jkb1
1882     REAL(wp) :: zthermal, zbeta
1883
1884
1885     DO jj = 2, jpjm1
1886        DO ji = 2, jpim1
1887           IF ( jbase(ji,jj)+1 < mbkt(ji,jj) ) THEN
1888              zthermal = rab_n(ji,jj,1,jp_tem) !ideally use ibld not 1??
1889              zbeta    = rab_n(ji,jj,1,jp_sal)
1890              jkb = jbase(ji,jj)
1891              jkb1 = MIN(jkb + 1, mbkt(ji,jj))
1892              zdtdz(ji,jj) = - ( tsn(ji,jj,jkb1,jp_tem) - tsn(ji,jj,jkb,jp_tem ) ) &
1893                   &   / e3t_n(ji,jj,jkb)
1894              zdsdz(ji,jj) = - ( tsn(ji,jj,jkb1,jp_sal) - tsn(ji,jj,jkb,jp_sal ) ) &
1895                   &   / e3t_n(ji,jj,jkb)
1896              zdbdz(ji,jj) = grav * zthermal * zdtdz(ji,jj) - grav * zbeta * zdsdz(ji,jj)
1897           ELSE
1898              zdtdz(ji,jj) = 0._wp
1899              zdsdz(ji,jj) = 0._wp
1900              zdbdz(ji,jj) = 0._wp
1901           END IF
1902        END DO
1903     END DO
1904    END SUBROUTINE zdf_osm_external_gradients
1905
1906    SUBROUTINE zdf_osm_pycnocline_scalar_profiles( zdtdz, zdsdz, zdbdz, zalpha )
1907
1908     REAL(wp), DIMENSION(jpi,jpj,jpk) :: zdtdz, zdsdz, zdbdz      ! gradients in the pycnocline
1909     REAL(wp), DIMENSION(jpi,jpj) :: zalpha
1910
1911     INTEGER :: jk, jj, ji
1912     REAL(wp) :: ztgrad, zsgrad, zbgrad
1913     REAL(wp) :: zgamma_b_nd, znd
1914     REAL(wp) :: zzeta_m, zzeta_en, zbuoy_pyc_sc
1915     REAL(wp), PARAMETER :: zgamma_b = 2.25, zzeta_sh = 0.15
1916
1917     DO jj = 2, jpjm1
1918        DO ji = 2, jpim1
1919           IF ( ibld(ji,jj) + jp_ext(ji,jj) < mbkt(ji,jj) ) THEN
1920              IF ( lconv(ji,jj) ) THEN  ! convective conditions
1921                IF ( lpyc(ji,jj) ) THEN
1922                   zzeta_m = 0.1 + 0.3 / ( 1.0 + EXP( -3.5 * LOG10( -zhol(ji,jj) ) ) )
1923                   zalpha(ji,jj) = 2.0 * ( 1.0 - ( 0.80 * zzeta_m + 0.5 * SQRT( 3.14159 / zgamma_b ) ) * zdbdz_bl_ext(ji,jj) * zdh(ji,jj) / zdb_ml(ji,jj) ) / ( 0.723 + SQRT( 3.14159 / zgamma_b ) )
1924                   zalpha(ji,jj) = MAX( zalpha(ji,jj), 0._wp )
1925
1926                   ztmp = 1._wp/MAX(zdh(ji,jj), epsln)
1927!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1928! Commented lines in this section are not needed in new code, once tested !
1929! can be removed                                                          !
1930!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1931!                   ztgrad = zalpha * zdt_ml(ji,jj) * ztmp + zdtdz_bl_ext(ji,jj)
1932!                   zsgrad = zalpha * zds_ml(ji,jj) * ztmp + zdsdz_bl_ext(ji,jj)
1933                   zbgrad = zalpha(ji,jj) * zdb_ml(ji,jj) * ztmp + zdbdz_bl_ext(ji,jj)
1934                   zgamma_b_nd = zdbdz_bl_ext(ji,jj) * zdh(ji,jj) / MAX(zdb_ml(ji,jj), epsln)
1935                   DO jk = 2, ibld(ji,jj)
1936                     znd = -( gdepw_n(ji,jj,jk) - zhbl(ji,jj) ) * ztmp
1937                     IF ( znd <= zzeta_m ) THEN
1938!                        zdtdz(ji,jj,jk) = zdtdz_bl_ext(ji,jj) + zalpha * zdt_ml(ji,jj) * ztmp * &
1939!                &        EXP( -6.0 * ( znd -zzeta_m )**2 )
1940!                        zdsdz(ji,jj,jk) = zdsdz_bl_ext(ji,jj) + zalpha * zds_ml(ji,jj) * ztmp * &
1941!                                  & EXP( -6.0 * ( znd -zzeta_m )**2 )
1942                        zdbdz(ji,jj,jk) = zdbdz_bl_ext(ji,jj) + zalpha(ji,jj) * zdb_ml(ji,jj) * ztmp * &
1943                                  & EXP( -6.0 * ( znd -zzeta_m )**2 )
1944                     ELSE
1945!                         zdtdz(ji,jj,jk) =  ztgrad * EXP( -zgamma_b * ( znd - zzeta_m )**2 )
1946!                         zdsdz(ji,jj,jk) =  zsgrad * EXP( -zgamma_b * ( znd - zzeta_m )**2 )
1947                         zdbdz(ji,jj,jk) =  zbgrad * EXP( -zgamma_b * ( znd - zzeta_m )**2 )
1948                     ENDIF
1949                  END DO
1950               ENDIF ! if no pycnocline pycnocline gradients set to zero
1951              ELSE
1952                 ! stable conditions
1953                 ! if pycnocline profile only defined when depth steady of increasing.
1954                 IF ( zdhdt(ji,jj) > 0.0 ) THEN        ! Depth increasing, or steady.
1955                    IF ( zdb_bl(ji,jj) > 0._wp ) THEN
1956                       IF ( zhol(ji,jj) >= 0.5 ) THEN      ! Very stable - 'thick' pycnocline
1957                          ztmp = 1._wp/MAX(zhbl(ji,jj), epsln)
1958                          ztgrad = zdt_bl(ji,jj) * ztmp
1959                          zsgrad = zds_bl(ji,jj) * ztmp
1960                          zbgrad = zdb_bl(ji,jj) * ztmp
1961                          DO jk = 2, ibld(ji,jj)
1962                             znd = gdepw_n(ji,jj,jk) * ztmp
1963                             zdtdz(ji,jj,jk) =  ztgrad * EXP( -15.0 * ( znd - 0.9 )**2 )
1964                             zdbdz(ji,jj,jk) =  zbgrad * EXP( -15.0 * ( znd - 0.9 )**2 )
1965                             zdsdz(ji,jj,jk) =  zsgrad * EXP( -15.0 * ( znd - 0.9 )**2 )
1966                          END DO
1967                       ELSE                                   ! Slightly stable - 'thin' pycnoline - needed when stable layer begins to form.
1968                          ztmp = 1._wp/MAX(zdh(ji,jj), epsln)
1969                          ztgrad = zdt_bl(ji,jj) * ztmp
1970                          zsgrad = zds_bl(ji,jj) * ztmp
1971                          zbgrad = zdb_bl(ji,jj) * ztmp
1972                          DO jk = 2, ibld(ji,jj)
1973                             znd = -( gdepw_n(ji,jj,jk) - zhml(ji,jj) ) * ztmp
1974                             zdtdz(ji,jj,jk) =  ztgrad * EXP( -1.75 * ( znd + 0.75 )**2 )
1975                             zdbdz(ji,jj,jk) =  zbgrad * EXP( -1.75 * ( znd + 0.75 )**2 )
1976                             zdsdz(ji,jj,jk) =  zsgrad * EXP( -1.75 * ( znd + 0.75 )**2 )
1977                          END DO
1978                       ENDIF ! IF (zhol >=0.5)
1979                    ENDIF    ! IF (zdb_bl> 0.)
1980                 ENDIF       ! IF (zdhdt >= 0) zdhdt < 0 not considered since pycnocline profile is zero and profile arrays are intialized to zero
1981              ENDIF          ! IF (lconv)
1982           ENDIF      ! IF ( ibld(ji,jj) < mbkt(ji,jj) )
1983        END DO
1984     END DO
1985
1986    END SUBROUTINE zdf_osm_pycnocline_scalar_profiles
1987
1988    SUBROUTINE zdf_osm_pycnocline_shear_profiles( zdudz, zdvdz )
1989      !!---------------------------------------------------------------------
1990      !!                   ***  ROUTINE zdf_osm_pycnocline_shear_profiles  ***
1991      !!
1992      !! ** Purpose : Calculates velocity shear in the pycnocline
1993      !!
1994      !! ** Method  :
1995      !!
1996      !!----------------------------------------------------------------------
1997
1998      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zdudz, zdvdz
1999
2000      INTEGER :: jk, jj, ji
2001      REAL(wp) :: zugrad, zvgrad, znd
2002      REAL(wp) :: zzeta_v = 0.45
2003      !
2004      DO jj = 2, jpjm1
2005         DO ji = 2, jpim1
2006            !
2007            IF ( ibld(ji,jj) + jp_ext(ji,jj) < mbkt(ji,jj) ) THEN
2008               IF ( lconv (ji,jj) ) THEN
2009                  ! Unstable conditions. Shouldn;t be needed with no pycnocline code.
2010!                  zugrad = 0.7 * zdu_ml(ji,jj) / zdh(ji,jj) + 0.3 * zustar(ji,jj)*zustar(ji,jj) / &
2011!                       &      ( ( ( zvstr(ji,jj)**3 + 0.5 * zwstrc(ji,jj)**3 )**pthird * zhml(ji,jj) ) * &
2012!                      &      MIN(zla(ji,jj)**(8.0/3.0) + epsln, 0.12 ))
2013                  !Alan is this right?
2014!                  zvgrad = ( 0.7 * zdv_ml(ji,jj) + &
2015!                       &    2.0 * ff_t(ji,jj) * zustke(ji,jj) * dstokes(ji,jj) / &
2016!                       &          ( ( zvstr(ji,jj)**3 + 0.5 * zwstrc(ji,jj)**3 )**pthird  + epsln ) &
2017!                       &      )/ (zdh(ji,jj)  + epsln )
2018!                  DO jk = 2, ibld(ji,jj) - 1 + ibld_ext
2019!                     znd = -( gdepw_n(ji,jj,jk) - zhbl(ji,jj) ) / (zdh(ji,jj) + epsln ) - zzeta_v
2020!                     IF ( znd <= 0.0 ) THEN
2021!                        zdudz(ji,jj,jk) = 1.25 * zugrad * EXP( 3.0 * znd )
2022!                        zdvdz(ji,jj,jk) = 1.25 * zvgrad * EXP( 3.0 * znd )
2023!                     ELSE
2024!                        zdudz(ji,jj,jk) = 1.25 * zugrad * EXP( -2.0 * znd )
2025!                        zdvdz(ji,jj,jk) = 1.25 * zvgrad * EXP( -2.0 * znd )
2026!                     ENDIF
2027!                  END DO
2028               ELSE
2029                  ! stable conditions
2030                  zugrad = 3.25 * zdu_bl(ji,jj) / zhbl(ji,jj)
2031                  zvgrad = 2.75 * zdv_bl(ji,jj) / zhbl(ji,jj)
2032                  DO jk = 2, ibld(ji,jj)
2033                     znd = gdepw_n(ji,jj,jk) / zhbl(ji,jj)
2034                     IF ( znd < 1.0 ) THEN
2035                        zdudz(ji,jj,jk) = zugrad * EXP( -40.0 * ( znd - 1.0 )**2 )
2036                     ELSE
2037                        zdudz(ji,jj,jk) = zugrad * EXP( -20.0 * ( znd - 1.0 )**2 )
2038                     ENDIF
2039                     zdvdz(ji,jj,jk) = zvgrad * EXP( -20.0 * ( znd - 0.85 )**2 )
2040                  END DO
2041               ENDIF
2042               !
2043            END IF      ! IF ( ibld(ji,jj) + ibld_ext < mbkt(ji,jj) )
2044         END DO
2045      END DO
2046    END SUBROUTINE zdf_osm_pycnocline_shear_profiles
2047
2048   SUBROUTINE zdf_osm_calculate_dhdt( zdhdt )
2049     !!---------------------------------------------------------------------
2050     !!                   ***  ROUTINE zdf_osm_calculate_dhdt  ***
2051     !!
2052     !! ** Purpose : Calculates the rate at which hbl changes.
2053     !!
2054     !! ** Method  :
2055     !!
2056     !!----------------------------------------------------------------------
2057
2058    REAL(wp), DIMENSION(jpi,jpj) :: zdhdt        ! Rate of change of hbl
2059
2060    INTEGER :: jj, ji
2061    REAL(wp) :: zgamma_b_nd, zgamma_dh_nd, zpert, zpsi
2062    REAL(wp) :: zvel_max,zddhdt
2063    REAL(wp), PARAMETER :: zzeta_m = 0.3
2064    REAL(wp), PARAMETER :: zgamma_c = 2.0
2065    REAL(wp), PARAMETER :: zdhoh = 0.1
2066    REAL(wp), PARAMETER :: zalpha_b = 0.3
2067    REAL, PARAMETER :: a_ddh = 2.5, a_ddh_2 = 3.5 ! also in pycnocline_depth
2068 
2069  DO jj = 2, jpjm1
2070     DO ji = 2, jpim1
2071       
2072       IF ( lshear(ji,jj) ) THEN
2073          IF ( lconv(ji,jj) ) THEN    ! Convective
2074
2075             IF ( ln_osm_mle ) THEN
2076
2077                IF ( hmle(ji,jj) > hbl(ji,jj) ) THEN
2078          ! Fox-Kemper buoyancy flux average over OSBL
2079                   zwb_fk_b(ji,jj) = zwb_fk(ji,jj) *  &
2080                        (1.0 + hmle(ji,jj) / ( 6.0 * hbl(ji,jj) ) * (-1.0 + ( 1.0 - 2.0 * hbl(ji,jj) / hmle(ji,jj))**3) )
2081                ELSE
2082                   zwb_fk_b(ji,jj) = 0.5 * zwb_fk(ji,jj) * hmle(ji,jj) / hbl(ji,jj)
2083                ENDIF
2084                zvel_max = ( zvstr(ji,jj)**3 + 0.5 * zwstrc(ji,jj)**3 )**p2third / hbl(ji,jj)
2085                IF ( ( zwb_ent(ji,jj) + 2.0 * zwb_fk_b(ji,jj) ) < 0.0 ) THEN
2086                   ! OSBL is deepening, entrainment > restratification
2087                   IF ( zdb_bl(ji,jj) > 1.0e-15 ) THEN
2088                      zgamma_b_nd = MAX( zdbdz_bl_ext(ji,jj), 0._wp ) * zdh(ji,jj) / zdb_ml(ji,jj)
2089                      zpsi = ( 1.0 - 0.5 * zdh(ji,jj) / zhbl(ji,jj) ) * ( zwb0(ji,jj) - MIN( ( zwb_min(ji,jj) + 2.0 * zwb_fk_b(ji,jj) ), 0._wp ) ) * zdh(ji,jj) / zhbl(ji,jj)
2090                      zpsi = zpsi + 1.75 * ( 1.0 - 0.5 * zdh(ji,jj) / zhbl(ji,jj) )*( zdh(ji,jj) / zhbl(ji,jj) + zgamma_b_nd ) * MIN( ( zwb_min(ji,jj) + 2.0 * zwb_fk_b(ji,jj) ), 0._wp )
2091                      zpsi = zalpha_b * MAX ( zpsi, 0._wp )
2092                      zdhdt(ji,jj) = -( zwb_ent(ji,jj) + 2.0 * zwb_fk_b(ji,jj) ) / ( zvel_max + MAX(zdb_bl(ji,jj), 1.0e-15) ) + zpsi / ( zvel_max + MAX( zdb_bl(ji,jj), 1.e-15 ) )
2093                      IF ( j_ddh(ji,jj) == 1 ) THEN
2094                        IF ( ( zwstrc(ji,jj) / zvstr(ji,jj) )**3 <= 0.5 ) THEN
2095                           zari = MIN( 1.5 * zdb_bl(ji,jj) / ( zhbl(ji,jj) * ( MAX(zdbdz_bl_ext(ji,jj),0._wp) + zdb_bl(ji,jj)**2 / MAX(4.5 * zvstr(ji,jj)**2 , 1.e-12 )) ), 0.2d0 )
2096                        ELSE
2097                           zari = MIN( 1.5 * zdb_bl(ji,jj) / ( zhbl(ji,jj) * ( MAX(zdbdz_bl_ext(ji,jj),0._wp) + zdb_bl(ji,jj)**2 / MAX(4.5 * zwstrc(ji,jj)**2 , 1.e-12 )) ), 0.2d0 )
2098                        ENDIF
2099! Relaxation to dh_ref = zari * hbl
2100                        zddhdt = -a_ddh_2 * ( 1.0 - zdh(ji,jj) / ( zari * zhbl(ji,jj) ) ) * zwb_ent(ji,jj) / zdb_bl(ji,jj)
2101                       
2102                      ELSE IF ( j_ddh(ji,jj) == 0 ) THEN
2103! Growing shear layer
2104                        zddhdt = -a_ddh * ( 1.0 - 1.6 * zdh(ji,jj) / zhbl(ji,jj) ) * zwb_ent(ji,jj) / zdb_bl(ji,jj)
2105                        zddhdt = EXP( - 4.0 * ABS( ff_t(ji,jj) ) * zhbl(ji,jj) / MAX(zustar(ji,jj), 1.e-8 ) ) * zddhdt
2106                      ELSE
2107                        zddhdt = 0._wp
2108                      ENDIF ! j_ddh
2109                      zdhdt(ji,jj) = zdhdt(ji,jj) + zalpha_b * ( 1.0 -0.5 * zdh(ji,jj) / zhbl(ji,jj) ) * zddhdt / ( zvel_max + MAX( zdb_bl(ji,jj), 1.0e-15 ) )
2110                   ELSE    ! zdb_bl >0
2111                       zdhdt(ji,jj) = -( zwb_ent(ji,jj) + 2.0 * zwb_fk_b(ji,jj) ) /  MAX( zvel_max, 1.0e-15)
2112                   ENDIF
2113                ELSE   ! zwb_min + 2*zwb_fk_b < 0
2114                   ! OSBL shoaling due to restratification flux. This is the velocity defined in Fox-Kemper et al (2008)
2115                   zdhdt(ji,jj) = - MIN(zvel_mle(ji,jj), hbl(ji,jj)/10800.)
2116
2117
2118                ENDIF
2119
2120             ELSE
2121                ! Fox-Kemper not used.
2122
2123                  zvel_max = - ( 1.0 + 1.0 * ( zvstr(ji,jj)**3 + 0.5 * zwstrc(ji,jj)**3 )**pthird * rn_rdt / hbl(ji,jj) ) * zwb_ent(ji,jj) / &
2124                  &        MAX((zvstr(ji,jj)**3 + 0.5 * zwstrc(ji,jj)**3 )**pthird, epsln)
2125                  zdhdt(ji,jj) = -zwb_ent(ji,jj) / ( zvel_max + MAX(zdb_bl(ji,jj), 1.0e-15) )
2126                ! added ajgn 23 July as temporay fix
2127
2128             ENDIF  ! ln_osm_mle
2129
2130            ELSE    ! lconv - Stable
2131                zdhdt(ji,jj) = ( 0.06 + 0.52 * zhol(ji,jj) / 2.0 ) * zvstr(ji,jj)**3 / hbl(ji,jj) + zwbav(ji,jj)
2132                IF ( zdhdt(ji,jj) < 0._wp ) THEN
2133                   ! For long timsteps factor in brackets slows the rapid collapse of the OSBL
2134                    zpert = 2.0 * ( 1.0 + 0.0 * 2.0 * zvstr(ji,jj) * rn_rdt / hbl(ji,jj) ) * zvstr(ji,jj)**2 / hbl(ji,jj)
2135                ELSE
2136                    zpert = MAX( zvstr(ji,jj)**2 / hbl(ji,jj), zdb_bl(ji,jj) )
2137                ENDIF
2138                zdhdt(ji,jj) = 2.0 * zdhdt(ji,jj) / MAX(zpert, epsln)
2139                zdhdt(ji,jj) = MAX(zdhdt(ji,jj), -hbl(ji,jj)/5400.)
2140            ENDIF  ! lconv
2141       ELSE ! lshear
2142         IF ( lconv(ji,jj) ) THEN    ! Convective
2143
2144             IF ( ln_osm_mle ) THEN
2145
2146                IF ( hmle(ji,jj) > hbl(ji,jj) ) THEN
2147          ! Fox-Kemper buoyancy flux average over OSBL
2148                   zwb_fk_b(ji,jj) = zwb_fk(ji,jj) *  &
2149                        (1.0 + hmle(ji,jj) / ( 6.0 * hbl(ji,jj) ) * (-1.0 + ( 1.0 - 2.0 * hbl(ji,jj) / hmle(ji,jj))**3) )
2150                ELSE
2151                   zwb_fk_b(ji,jj) = 0.5 * zwb_fk(ji,jj) * hmle(ji,jj) / hbl(ji,jj)
2152                ENDIF
2153                zvel_max = ( zwstrl(ji,jj)**3 + 0.5 * zwstrc(ji,jj)**3 )**p2third / hbl(ji,jj)
2154                IF ( ( zwb_ent(ji,jj) + 2.0 * zwb_fk_b(ji,jj) ) < 0.0 ) THEN
2155                   ! OSBL is deepening, entrainment > restratification
2156                   IF ( zdb_bl(ji,jj) > 0.0 .and. zdbdz_bl_ext(ji,jj) > 0.0 ) THEN
2157                      zdhdt(ji,jj) = -( zwb_ent(ji,jj) + 2.0 * zwb_fk_b(ji,jj) ) / ( zvel_max + MAX(zdb_bl(ji,jj), 1.0e-15) )
2158                   ELSE
2159                      zdhdt(ji,jj) = -( zwb_ent(ji,jj) + 2.0 * zwb_fk_b(ji,jj) ) /  MAX( zvel_max, 1.0e-15)
2160                   ENDIF
2161                ELSE
2162                   ! OSBL shoaling due to restratification flux. This is the velocity defined in Fox-Kemper et al (2008)
2163                   zdhdt(ji,jj) = - MIN(zvel_mle(ji,jj), hbl(ji,jj)/10800.)
2164
2165
2166                ENDIF
2167
2168             ELSE
2169                ! Fox-Kemper not used.
2170
2171                  zvel_max = -zwb_ent(ji,jj) / &
2172                  &        MAX((zvstr(ji,jj)**3 + 0.5 * zwstrc(ji,jj)**3 )**pthird, epsln)
2173                  zdhdt(ji,jj) = -zwb_ent(ji,jj) / ( zvel_max + MAX(zdb_bl(ji,jj), 1.0e-15) )
2174                ! added ajgn 23 July as temporay fix
2175
2176             ENDIF  ! ln_osm_mle
2177
2178            ELSE                        ! Stable
2179                zdhdt(ji,jj) = ( 0.06 + 0.52 * zhol(ji,jj) / 2.0 ) * zvstr(ji,jj)**3 / hbl(ji,jj) + zwbav(ji,jj)
2180                IF ( zdhdt(ji,jj) < 0._wp ) THEN
2181                   ! For long timsteps factor in brackets slows the rapid collapse of the OSBL
2182                    zpert = 2.0 * zvstr(ji,jj)**2 / hbl(ji,jj)
2183                ELSE
2184                    zpert = MAX( zvstr(ji,jj)**2 / hbl(ji,jj), zdb_bl(ji,jj) )
2185                ENDIF
2186                zdhdt(ji,jj) = 2.0 * zdhdt(ji,jj) / MAX(zpert, epsln)
2187                zdhdt(ji,jj) = MAX(zdhdt(ji,jj), -hbl(ji,jj)/5400.)
2188            ENDIF  ! lconv
2189         ENDIF ! lshear
2190       END DO
2191     END DO
2192    END SUBROUTINE zdf_osm_calculate_dhdt
2193
2194    SUBROUTINE zdf_osm_timestep_hbl( zdhdt )
2195     !!---------------------------------------------------------------------
2196     !!                   ***  ROUTINE zdf_osm_timestep_hbl  ***
2197     !!
2198     !! ** Purpose : Increments hbl.
2199     !!
2200     !! ** Method  : If thechange in hbl exceeds one model level the change is
2201     !!              is calculated by moving down the grid, changing the buoyancy
2202     !!              jump. This is to ensure that the change in hbl does not
2203     !!              overshoot a stable layer.
2204     !!
2205     !!----------------------------------------------------------------------
2206
2207
2208    REAL(wp), DIMENSION(jpi,jpj) :: zdhdt   ! rates of change of hbl.
2209
2210    INTEGER :: jk, jj, ji, jm
2211    REAL(wp) :: zhbl_s, zvel_max, zdb
2212    REAL(wp) :: zthermal, zbeta
2213
2214     DO jj = 2, jpjm1
2215         DO ji = 2, jpim1
2216           IF ( ibld(ji,jj) - imld(ji,jj) > 1 ) THEN
2217!
2218! If boundary layer changes by more than one level, need to check for stable layers between initial and final depths.
2219!
2220              zhbl_s = hbl(ji,jj)
2221              jm = imld(ji,jj)
2222              zthermal = rab_n(ji,jj,1,jp_tem)
2223              zbeta = rab_n(ji,jj,1,jp_sal)
2224
2225
2226              IF ( lconv(ji,jj) ) THEN
2227!unstable
2228
2229                 IF( ln_osm_mle ) THEN
2230                    zvel_max = ( zwstrl(ji,jj)**3 + zwstrc(ji,jj)**3 )**p2third / hbl(ji,jj)
2231                 ELSE
2232
2233                    zvel_max = -( 1.0 + 1.0 * ( zvstr(ji,jj)**3 + 0.5 * zwstrc(ji,jj)**3 )**pthird * rn_rdt / hbl(ji,jj) ) * zwb_ent(ji,jj) / &
2234                      &      ( zvstr(ji,jj)**3 + 0.5 * zwstrc(ji,jj)**3 )**pthird
2235
2236                 ENDIF
2237
2238                 DO jk = imld(ji,jj), ibld(ji,jj)
2239                    zdb = MAX( grav * ( zthermal * ( zt_bl(ji,jj) - tsn(ji,jj,jm,jp_tem) ) &
2240                         & - zbeta * ( zs_bl(ji,jj) - tsn(ji,jj,jm,jp_sal) ) ), &
2241                         &  0.0 ) + zvel_max
2242
2243
2244                    IF ( ln_osm_mle ) THEN
2245                       zhbl_s = zhbl_s + MIN( &
2246                        & rn_rdt * ( ( -zwb_ent(ji,jj) - 2.0 * zwb_fk_b(ji,jj) )/ zdb ) / FLOAT(ibld(ji,jj) - imld(ji,jj) ), &
2247                        & e3w_n(ji,jj,jm) )
2248                    ELSE
2249                      zhbl_s = zhbl_s + MIN( &
2250                        & rn_rdt * (  -zwb_ent(ji,jj) / zdb ) / FLOAT(ibld(ji,jj) - imld(ji,jj) ), &
2251                        & e3w_n(ji,jj,jm) )
2252                    ENDIF
2253
2254!                    zhbl_s = MIN(zhbl_s,  gdepw_n(ji,jj, mbkt(ji,jj) + 1) - depth_tol)
2255                    IF ( zhbl_s >= gdepw_n(ji,jj,mbkt(ji,jj) + 1) ) THEN
2256                      zhbl_s = MIN(zhbl_s,  gdepw_n(ji,jj, mbkt(ji,jj) + 1) - depth_tol)
2257                      lpyc(ji,jj) = .FALSE.
2258                    ENDIF
2259                    IF ( zhbl_s >= gdepw_n(ji,jj,jm+1) ) jm = jm + 1
2260                 END DO
2261                 hbl(ji,jj) = zhbl_s
2262                 ibld(ji,jj) = jm
2263             ELSE
2264! stable
2265                 DO jk = imld(ji,jj), ibld(ji,jj)
2266                    zdb = MAX( &
2267                         & grav * ( zthermal * ( zt_bl(ji,jj) - tsn(ji,jj,jm,jp_tem) )&
2268                         &           - zbeta * ( zs_bl(ji,jj) - tsn(ji,jj,jm,jp_sal) ) ),&
2269                         & 0.0 ) + &
2270             &       2.0 * zvstr(ji,jj)**2 / zhbl_s
2271
2272                    ! Alan is thuis right? I have simply changed hbli to hbl
2273                    zhol(ji,jj) = -zhbl_s / ( ( zvstr(ji,jj)**3 + epsln )/ zwbav(ji,jj) )
2274                    zdhdt(ji,jj) = -( zwbav(ji,jj) - 0.04 / 2.0 * zwstrl(ji,jj)**3 / zhbl_s - 0.15 / 2.0 * ( 1.0 - EXP( -1.5 * zla(ji,jj) ) ) * &
2275               &                  zustar(ji,jj)**3 / zhbl_s ) * ( 0.725 + 0.225 * EXP( -7.5 * zhol(ji,jj) ) )
2276                    zdhdt(ji,jj) = zdhdt(ji,jj) + zwbav(ji,jj)
2277                    zhbl_s = zhbl_s + MIN( zdhdt(ji,jj) / zdb * rn_rdt / FLOAT( ibld(ji,jj) - imld(ji,jj) ), e3w_n(ji,jj,jm) )
2278
2279!                    zhbl_s = MIN(zhbl_s, gdepw_n(ji,jj, mbkt(ji,jj) + 1) - depth_tol)
2280                    IF ( zhbl_s >= mbkt(ji,jj) + 1 ) THEN
2281                      zhbl_s = MIN(zhbl_s,  gdepw_n(ji,jj, mbkt(ji,jj) + 1) - depth_tol)
2282                      lpyc(ji,jj) = .FALSE.
2283                    ENDIF
2284                    IF ( zhbl_s >= gdepw_n(ji,jj,jm) ) jm = jm + 1
2285                 END DO
2286             ENDIF   ! IF ( lconv )
2287             hbl(ji,jj) = MAX(zhbl_s, gdepw_n(ji,jj,4) )
2288             ibld(ji,jj) = MAX(jm, 4 )
2289           ELSE
2290! change zero or one model level.
2291             hbl(ji,jj) = MAX(zhbl_t(ji,jj), gdepw_n(ji,jj,4) )
2292           ENDIF
2293           zhbl(ji,jj) = gdepw_n(ji,jj,ibld(ji,jj))
2294         END DO
2295      END DO
2296
2297    END SUBROUTINE zdf_osm_timestep_hbl
2298
2299    SUBROUTINE zdf_osm_pycnocline_thickness( dh, zdh )
2300      !!---------------------------------------------------------------------
2301      !!                   ***  ROUTINE zdf_osm_pycnocline_thickness  ***
2302      !!
2303      !! ** Purpose : Calculates thickness of the pycnocline
2304      !!
2305      !! ** Method  : The thickness is calculated from a prognostic equation
2306      !!              that relaxes the pycnocine thickness to a diagnostic
2307      !!              value. The time change is calculated assuming the
2308      !!              thickness relaxes exponentially. This is done to deal
2309      !!              with large timesteps.
2310      !!
2311      !!----------------------------------------------------------------------
2312
2313      REAL(wp), DIMENSION(jpi,jpj) :: dh, zdh     ! pycnocline thickness.
2314       !
2315      INTEGER :: jj, ji
2316      INTEGER :: inhml
2317      REAL(wp) :: zari, ztau, zdh_ref, zddhdt, zvel_max
2318      REAL, PARAMETER :: a_ddh = 2.5, a_ddh_2 = 3.5 ! also in pycnocline_depth
2319
2320    DO jj = 2, jpjm1
2321       DO ji = 2, jpim1
2322
2323         IF ( lshear(ji,jj) ) THEN
2324            IF ( lconv(ji,jj) ) THEN
2325             IF ( zdb_bl(ji,jj) > 1.0e-15) THEN
2326              IF ( j_ddh(ji,jj) == 0 ) THEN
2327                zvel_max = ( zvstr(ji,jj)**3 + 0.5 * zwstrc(ji,jj)**3 )**p2third / hbl(ji,jj)
2328! ddhdt for pycnocline determined in osm_calculate_dhdt
2329                zddhdt = -a_ddh * ( 1.0 - 1.6 * zdh(ji,jj) / zhbl(ji,jj) ) * zwb_ent(ji,jj) / ( zvel_max + MAX( zdb_bl(ji,jj), 1.0e-15 ) )
2330                zddhdt = EXP( - 4.0 * ABS( ff_t(ji,jj) ) * zhbl(ji,jj) / MAX(zustar(ji,jj), 1.e-8 ) ) * zddhdt
2331! maximum limit for how thick the shear layer can grow relative to the thickness of the boundary kayer
2332                dh(ji,jj) = MIN( dh(ji,jj) + zddhdt * rn_rdt, 0.625 * hbl(ji,jj) )
2333              ELSE
2334! Need to recalculate because hbl has been updated.
2335                IF ( ( zwstrc(ji,jj) / zvstr(ji,jj) )**3 <= 0.5 ) THEN
2336                  zari = MIN( 1.5 * zdb_bl(ji,jj) / ( zhbl(ji,jj) * ( MAX(zdbdz_bl_ext(ji,jj),0._wp) + zdb_bl(ji,jj)**2 / MAX(4.5 * zvstr(ji,jj)**2 , 1.e-12 )) ), 0.2d0 )
2337                ELSE
2338                  zari = MIN( 1.5 * zdb_bl(ji,jj) / ( zhbl(ji,jj) * ( MAX(zdbdz_bl_ext(ji,jj),0._wp) + zdb_bl(ji,jj)**2 / MAX(4.5 * zwstrc(ji,jj)**2 , 1.e-12 )) ), 0.2d0 )
2339                ENDIF
2340                ztau = MAX( zdb_bl(ji,jj) * ( zari * hbl(ji,jj) ) / ( a_ddh_2 * MAX(-zwb_ent(ji,jj), 1.e-12) ), 2.0 * rn_rdt )
2341                dh(ji,jj) = dh(ji,jj) * EXP( -rn_rdt / ztau ) + zari * zhbl(ji,jj) * ( 1.0 - EXP( -rn_rdt / ztau ) )
2342                IF ( dh(ji,jj) >= hbl(ji,jj) ) dh(ji,jj) = zari * zhbl(ji,jj)
2343              ENDIF
2344             ELSE
2345              ztau = MAX( MAX( hbl(ji,jj) / ( zvstr(ji,jj)**3 + 0.5 * zwstrc(ji,jj)**3 )**pthird, epsln), 2.0 * rn_rdt )
2346              dh(ji,jj) = dh(ji,jj) * EXP( -rn_rdt / ztau ) + 0.2 * zhbl(ji,jj) * ( 1.0 - EXP( -rn_rdt / ztau ) )
2347              IF ( dh(ji,jj) > hbl(ji,jj) ) dh(ji,jj) = 0.2 * hbl(ji,jj)
2348             ENDIF
2349            ELSE ! lconv
2350! Initially shear only for entraining OSBL. Stable code will be needed if extended to stable OSBL
2351
2352               ztau = hbl(ji,jj) / MAX(zvstr(ji,jj), epsln)
2353               IF ( zdhdt(ji,jj) >= 0.0 ) THEN    ! probably shouldn't include wm here
2354                  ! boundary layer deepening
2355                  IF ( zdb_bl(ji,jj) > 0._wp ) THEN
2356                     ! pycnocline thickness set by stratification - use same relationship as for neutral conditions.
2357                     zari = MIN( 4.5 * ( zvstr(ji,jj)**2 ) &
2358                          & /  MAX(zdb_bl(ji,jj) * zhbl(ji,jj), epsln ) + 0.01  , 0.2 )
2359                     zdh_ref = MIN( zari, 0.2 ) * hbl(ji,jj)
2360                  ELSE
2361                     zdh_ref = 0.2 * hbl(ji,jj)
2362                  ENDIF
2363               ELSE     ! IF(dhdt < 0)
2364                  zdh_ref = 0.2 * hbl(ji,jj)
2365               ENDIF    ! IF (dhdt >= 0)
2366               dh(ji,jj) = dh(ji,jj) * EXP( -rn_rdt / ztau )+ zdh_ref * ( 1.0 - EXP( -rn_rdt / ztau ) )
2367               IF ( zdhdt(ji,jj) < 0._wp .and. dh(ji,jj) >= hbl(ji,jj) ) dh(ji,jj) = zdh_ref       ! can be a problem with dh>hbl for rapid collapse
2368               ! Alan: this hml is never defined or used -- do we need it?
2369            ENDIF
2370             
2371         ELSE   ! lshear 
2372! for lshear = .FALSE. calculate ddhdt here
2373
2374             IF ( lconv(ji,jj) ) THEN
2375
2376               IF( ln_osm_mle ) THEN
2377                  IF ( ( zwb_ent(ji,jj) + 2.0 * zwb_fk_b(ji,jj) ) < 0._wp ) THEN
2378                     ! OSBL is deepening. Note wb_fk_b is zero if ln_osm_mle=F
2379                     IF ( zdb_bl(ji,jj) > 0._wp .and. zdbdz_bl_ext(ji,jj) > 0._wp)THEN
2380                        IF ( ( zwstrc(ji,jj) / MAX(zvstr(ji,jj), epsln) )**3 <= 0.5 ) THEN  ! near neutral stability
2381                           zari = MIN( 1.5 * zdb_bl(ji,jj) / ( zhbl(ji,jj) * ( MAX(zdbdz_bl_ext(ji,jj),0._wp) + zdb_bl(ji,jj)**2 / MAX(4.5 * zvstr(ji,jj)**2 , 1.e-12 )) ), 0.2d0 )
2382                        ELSE                                                     ! unstable
2383                           zari = MIN( 1.5 * zdb_bl(ji,jj) / ( zhbl(ji,jj) * ( MAX(zdbdz_bl_ext(ji,jj),0._wp) + zdb_bl(ji,jj)**2 / MAX(4.5 * zwstrc(ji,jj)**2 , 1.e-12 )) ), 0.2d0 )
2384                        ENDIF
2385                        ztau = 0.2 * hbl(ji,jj) / MAX(epsln, (zvstr(ji,jj)**3 + 0.5 *zwstrc(ji,jj)**3)**pthird)
2386                        zdh_ref = zari * hbl(ji,jj)
2387                     ELSE
2388                        ztau = 0.2 * hbl(ji,jj) / MAX(epsln, (zvstr(ji,jj)**3 + 0.5 *zwstrc(ji,jj)**3)**pthird)
2389                        zdh_ref = 0.2 * hbl(ji,jj)
2390                     ENDIF
2391                  ELSE
2392                     ztau = 0.2 * hbl(ji,jj) /  MAX(epsln, (zvstr(ji,jj)**3 + 0.5 *zwstrc(ji,jj)**3)**pthird)
2393                     zdh_ref = 0.2 * hbl(ji,jj)
2394                  ENDIF
2395               ELSE ! ln_osm_mle
2396                  IF ( zdb_bl(ji,jj) > 0._wp .and. zdbdz_bl_ext(ji,jj) > 0._wp)THEN
2397                     IF ( ( zwstrc(ji,jj) / MAX(zvstr(ji,jj), epsln) )**3 <= 0.5 ) THEN  ! near neutral stability
2398                        zari = MIN( 1.5 * zdb_bl(ji,jj) / ( zhbl(ji,jj) * ( MAX(zdbdz_bl_ext(ji,jj),0._wp) + zdb_bl(ji,jj)**2 / MAX(4.5 * zvstr(ji,jj)**2 , 1.e-12 )) ), 0.2d0 )
2399                     ELSE                                                     ! unstable
2400                        zari = MIN( 1.5 * zdb_bl(ji,jj) / ( zhbl(ji,jj) * ( MAX(zdbdz_bl_ext(ji,jj),0._wp) + zdb_bl(ji,jj)**2 / MAX(4.5 * zwstrc(ji,jj)**2 , 1.e-12 )) ), 0.2d0 )
2401                     ENDIF
2402                     ztau = hbl(ji,jj) / MAX(epsln, (zvstr(ji,jj)**3 + 0.5 *zwstrc(ji,jj)**3)**pthird)
2403                     zdh_ref = zari * hbl(ji,jj)
2404                  ELSE
2405                     ztau = hbl(ji,jj) / MAX(epsln, (zvstr(ji,jj)**3 + 0.5 *zwstrc(ji,jj)**3)**pthird)
2406                     zdh_ref = 0.2 * hbl(ji,jj)
2407                  ENDIF
2408
2409               END IF  ! ln_osm_mle
2410
2411               dh(ji,jj) = dh(ji,jj) * EXP( -rn_rdt / ztau ) + zdh_ref * ( 1.0 - EXP( -rn_rdt / ztau ) )
2412!               IF ( zdhdt(ji,jj) < 0._wp .and. dh(ji,jj) >= hbl(ji,jj) ) dh(ji,jj) = zdh_ref
2413               IF ( dh(ji,jj) >= hbl(ji,jj) ) dh(ji,jj) = zdh_ref
2414               ! Alan: this hml is never defined or used
2415            ELSE   ! IF (lconv)
2416               ztau = hbl(ji,jj) / MAX(zvstr(ji,jj), epsln)
2417               IF ( zdhdt(ji,jj) >= 0.0 ) THEN    ! probably shouldn't include wm here
2418                  ! boundary layer deepening
2419                  IF ( zdb_bl(ji,jj) > 0._wp ) THEN
2420                     ! pycnocline thickness set by stratification - use same relationship as for neutral conditions.
2421                     zari = MIN( 4.5 * ( zvstr(ji,jj)**2 ) &
2422                          & /  MAX(zdb_bl(ji,jj) * zhbl(ji,jj), epsln ) + 0.01  , 0.2 )
2423                     zdh_ref = MIN( zari, 0.2 ) * hbl(ji,jj)
2424                  ELSE
2425                     zdh_ref = 0.2 * hbl(ji,jj)
2426                  ENDIF
2427               ELSE     ! IF(dhdt < 0)
2428                  zdh_ref = 0.2 * hbl(ji,jj)
2429               ENDIF    ! IF (dhdt >= 0)
2430               dh(ji,jj) = dh(ji,jj) * EXP( -rn_rdt / ztau )+ zdh_ref * ( 1.0 - EXP( -rn_rdt / ztau ) )
2431               IF ( zdhdt(ji,jj) < 0._wp .and. dh(ji,jj) >= hbl(ji,jj) ) dh(ji,jj) = zdh_ref       ! can be a problem with dh>hbl for rapid collapse
2432            ENDIF       ! IF (lconv)
2433         ENDIF  ! lshear
2434 
2435         hml(ji,jj) = hbl(ji,jj) - dh(ji,jj)
2436         inhml = MAX( INT( dh(ji,jj) / MAX(e3t_n(ji,jj,ibld(ji,jj)), 1.e-3) ) , 1 )
2437         imld(ji,jj) = MAX( ibld(ji,jj) - inhml, 3)
2438         zhml(ji,jj) = gdepw_n(ji,jj,imld(ji,jj))
2439         zdh(ji,jj) = zhbl(ji,jj) - zhml(ji,jj)
2440       END DO
2441     END DO
2442
2443    END SUBROUTINE zdf_osm_pycnocline_thickness
2444
2445
2446   SUBROUTINE zdf_osm_zmld_horizontal_gradients( zmld, zdtdx, zdtdy, zdsdx, zdsdy, dbdx_mle, dbdy_mle, zdbds_mle )
2447      !!----------------------------------------------------------------------
2448      !!                  ***  ROUTINE zdf_osm_horizontal_gradients  ***
2449      !!
2450      !! ** Purpose :   Calculates horizontal gradients of buoyancy for use with Fox-Kemper parametrization.
2451      !!
2452      !! ** Method  :
2453      !!
2454      !! References: Fox-Kemper et al., JPO, 38, 1145-1165, 2008
2455      !!             Fox-Kemper and Ferrari, JPO, 38, 1166-1179, 2008
2456
2457
2458      REAL(wp), DIMENSION(jpi,jpj)     :: dbdx_mle, dbdy_mle ! MLE horiz gradients at u & v points
2459      REAL(wp), DIMENSION(jpi,jpj)     :: zdbds_mle ! Magnitude of horizontal buoyancy gradient.
2460      REAL(wp), DIMENSION(jpi,jpj)     :: zmld ! ==  estimated FK BLD used for MLE horiz gradients  == !
2461      REAL(wp), DIMENSION(jpi,jpj)     :: zdtdx, zdtdy, zdsdx, zdsdy
2462
2463      INTEGER  ::   ji, jj, jk          ! dummy loop indices
2464      INTEGER  ::   ii, ij, ik, ikmax   ! local integers
2465      REAL(wp)                         :: zc
2466      REAL(wp)                         :: zN2_c           ! local buoyancy difference from 10m value
2467      REAL(wp), DIMENSION(jpi,jpj)     :: ztm, zsm, zLf_NH, zLf_MH
2468      REAL(wp), DIMENSION(jpi,jpj,jpts):: ztsm_midu, ztsm_midv, zabu, zabv
2469      REAL(wp), DIMENSION(jpi,jpj)     :: zmld_midu, zmld_midv
2470!!----------------------------------------------------------------------
2471      !
2472      !                                      !==  MLD used for MLE  ==!
2473
2474      mld_prof(:,:)  = nlb10               ! Initialization to the number of w ocean point
2475      zmld(:,:)  = 0._wp               ! here hmlp used as a dummy variable, integrating vertically N^2
2476      zN2_c = grav * rn_osm_mle_rho_c * r1_rau0   ! convert density criteria into N^2 criteria
2477      DO jk = nlb10, jpkm1
2478         DO jj = 1, jpj                ! Mixed layer level: w-level
2479            DO ji = 1, jpi
2480               ikt = mbkt(ji,jj)
2481               zmld(ji,jj) = zmld(ji,jj) + MAX( rn2b(ji,jj,jk) , 0._wp ) * e3w_n(ji,jj,jk)
2482               IF( zmld(ji,jj) < zN2_c )   mld_prof(ji,jj) = MIN( jk , ikt ) + 1   ! Mixed layer level
2483            END DO
2484         END DO
2485      END DO
2486      DO jj = 1, jpj
2487         DO ji = 1, jpi
2488            mld_prof(ji,jj) = MAX(mld_prof(ji,jj),ibld(ji,jj))
2489            zmld(ji,jj) = gdepw_n(ji,jj,mld_prof(ji,jj))
2490         END DO
2491      END DO
2492      ! ensure mld_prof .ge. ibld
2493      !
2494      ikmax = MIN( MAXVAL( mld_prof(:,:) ), jpkm1 )                  ! max level of the computation
2495      !
2496      ztm(:,:) = 0._wp
2497      zsm(:,:) = 0._wp
2498      DO jk = 1, ikmax                                 ! MLD and mean buoyancy and N2 over the mixed layer
2499         DO jj = 1, jpj
2500            DO ji = 1, jpi
2501               zc = e3t_n(ji,jj,jk) * REAL( MIN( MAX( 0, mld_prof(ji,jj)-jk ) , 1  )  )    ! zc being 0 outside the ML t-points
2502               ztm(ji,jj) = ztm(ji,jj) + zc * tsn(ji,jj,jk,jp_tem)
2503               zsm(ji,jj) = zsm(ji,jj) + zc * tsn(ji,jj,jk,jp_sal)
2504            END DO
2505         END DO
2506      END DO
2507      ! average temperature and salinity.
2508      ztm(:,:) = ztm(:,:) / MAX( e3t_n(:,:,1), zmld(:,:) )
2509      zsm(:,:) = zsm(:,:) / MAX( e3t_n(:,:,1), zmld(:,:) )
2510      ! calculate horizontal gradients at u & v points
2511
2512      DO jj = 2, jpjm1
2513         DO ji = 1, jpim1
2514            zdtdx(ji,jj) = ( ztm(ji+1,jj) - ztm( ji,jj) )  * umask(ji,jj,1) / e1u(ji,jj)
2515            zdsdx(ji,jj) = ( zsm(ji+1,jj) - zsm( ji,jj) )  * umask(ji,jj,1) / e1u(ji,jj)
2516            zmld_midu(ji,jj) = 0.25_wp * (zmld(ji+1,jj) + zmld( ji,jj))
2517            ztsm_midu(ji,jj,jp_tem) = 0.5_wp * ( ztm(ji+1,jj) + ztm( ji,jj) )
2518            ztsm_midu(ji,jj,jp_sal) = 0.5_wp * ( zsm(ji+1,jj) + zsm( ji,jj) )
2519         END DO
2520      END DO
2521
2522      DO jj = 1, jpjm1
2523         DO ji = 2, jpim1
2524            zdtdy(ji,jj) = ( ztm(ji,jj+1) - ztm( ji,jj) ) * vmask(ji,jj,1) / e1v(ji,jj)
2525            zdsdy(ji,jj) = ( zsm(ji,jj+1) - zsm( ji,jj) ) * vmask(ji,jj,1) / e1v(ji,jj)
2526            zmld_midv(ji,jj) = 0.25_wp * (zmld(ji,jj+1) + zmld( ji,jj))
2527            ztsm_midv(ji,jj,jp_tem) = 0.5_wp * ( ztm(ji,jj+1) + ztm( ji,jj) )
2528            ztsm_midv(ji,jj,jp_sal) = 0.5_wp * ( zsm(ji,jj+1) + zsm( ji,jj) )
2529         END DO
2530      END DO
2531
2532      ! ensure salinity > 0 in unset values so EOS doesn't give FP error with fpe0 on
2533      ztsm_midu(:,jpj,jp_sal) = 10.
2534      ztsm_midv(jpi,:,jp_sal) = 10.
2535
2536      CALL eos_rab(ztsm_midu, zmld_midu, zabu)
2537      CALL eos_rab(ztsm_midv, zmld_midv, zabv)
2538
2539      DO jj = 2, jpjm1
2540         DO ji = 1, jpim1
2541            dbdx_mle(ji,jj) = grav*(zdtdx(ji,jj)*zabu(ji,jj,jp_tem) - zdsdx(ji,jj)*zabu(ji,jj,jp_sal))
2542         END DO
2543      END DO
2544      DO jj = 1, jpjm1
2545         DO ji = 2, jpim1
2546            dbdy_mle(ji,jj) = grav*(zdtdy(ji,jj)*zabv(ji,jj,jp_tem) - zdsdy(ji,jj)*zabv(ji,jj,jp_sal))
2547         END DO
2548      END DO
2549
2550      DO jj = 2, jpjm1
2551        DO ji = 2, jpim1
2552           ztmp =  r1_ft(ji,jj) *  MIN( 111.e3_wp , e1u(ji,jj) ) / rn_osm_mle_lf
2553           zdbds_mle(ji,jj) = SQRT( 0.5_wp * ( dbdx_mle(ji,jj) * dbdx_mle(ji,jj) + dbdy_mle(ji,jj) * dbdy_mle(ji,jj) &
2554                & + dbdx_mle(ji-1,jj) * dbdx_mle(ji-1,jj) + dbdy_mle(ji,jj-1) * dbdy_mle(ji,jj-1) ) )
2555        END DO
2556      END DO
2557     
2558 END SUBROUTINE zdf_osm_zmld_horizontal_gradients
2559  SUBROUTINE zdf_osm_mle_parameters( zmld,  mld_prof, hmle, zhmle, zvel_mle, zdiff_mle )
2560      !!----------------------------------------------------------------------
2561      !!                  ***  ROUTINE zdf_osm_mle_parameters  ***
2562      !!
2563      !! ** Purpose :   Timesteps the mixed layer eddy depth, hmle and calculates the mixed layer eddy fluxes for buoyancy, heat and salinity.
2564      !!
2565      !! ** Method  :
2566      !!
2567      !! References: Fox-Kemper et al., JPO, 38, 1145-1165, 2008
2568      !!             Fox-Kemper and Ferrari, JPO, 38, 1166-1179, 2008
2569
2570      REAL(wp), DIMENSION(jpi,jpj)     :: zmld ! ==  estimated FK BLD used for MLE horiz gradients  == !
2571      INTEGER, DIMENSION(jpi,jpj)      :: mld_prof
2572      REAL(wp), DIMENSION(jpi,jpj)     :: hmle, zhmle, zwb_fk, zvel_mle, zdiff_mle
2573      INTEGER  ::   ji, jj, jk          ! dummy loop indices
2574      INTEGER  ::   ii, ij, ik, jkb, jkb1  ! local integers
2575      INTEGER , DIMENSION(jpi,jpj)     :: inml_mle
2576      REAL(wp) ::  ztmp, zdbdz, zdtdz, zdsdz, zthermal,zbeta, zbuoy, zdb_mle
2577
2578   ! Calculate vertical buoyancy, heat and salinity fluxes due to MLE.
2579
2580      DO jj = 2, jpjm1
2581        DO ji = 2, jpim1
2582          IF ( lconv(ji,jj) ) THEN
2583             ztmp =  r1_ft(ji,jj) *  MIN( 111.e3_wp , e1u(ji,jj) ) / rn_osm_mle_lf
2584      ! This velocity scale, defined in Fox-Kemper et al (2008), is needed for calculating dhdt.
2585             zvel_mle(ji,jj) = zdbds_mle(ji,jj) * ztmp * hmle(ji,jj) * tmask(ji,jj,1)
2586             zdiff_mle(ji,jj) = 5.e-4_wp * rn_osm_mle_ce * ztmp * zdbds_mle(ji,jj) * zhmle(ji,jj)**2
2587          ENDIF
2588        END DO
2589      END DO
2590   ! Timestep mixed layer eddy depth.
2591      DO jj = 2, jpjm1
2592        DO ji = 2, jpim1
2593           IF ( lmle(ji,jj) ) THEN  ! MLE layer growing.
2594! Buoyancy gradient at base of MLE layer.
2595              zthermal = rab_n(ji,jj,1,jp_tem)
2596              zbeta    = rab_n(ji,jj,1,jp_sal)
2597              jkb = mld_prof(ji,jj)
2598              jkb1 = MIN(jkb + 1, mbkt(ji,jj))
2599!             
2600              zbuoy = grav * ( zthermal * tsn(ji,jj,mld_prof(ji,jj)+2,jp_tem) - zbeta * tsn(ji,jj,mld_prof(ji,jj)+2,jp_sal) )
2601              zdb_mle = zb_bl(ji,jj) - zbuoy 
2602! Timestep hmle.
2603              hmle(ji,jj) = hmle(ji,jj) + zwb0tot(ji,jj) * rn_rdt / zdb_mle
2604           ELSE
2605              IF ( zhmle(ji,jj) > zhbl(ji,jj) ) THEN
2606                 hmle(ji,jj) = hmle(ji,jj) - ( hmle(ji,jj) - hbl(ji,jj) ) * rn_rdt / rn_osm_mle_tau
2607              ELSE
2608                 hmle(ji,jj) = hmle(ji,jj) - 10.0 * ( hmle(ji,jj) - hbl(ji,jj) ) * rn_rdt /rn_osm_mle_tau
2609              ENDIF
2610           ENDIF
2611           hmle(ji,jj) = MAX(MIN(hmle(ji,jj), ht_n(ji,jj)),  gdepw_n(ji,jj,4))
2612           IF(ln_osm_hmle_limit) hmle(ji,jj) = MIN(hmle(ji,jj), rn_osm_hmle_limit*hbl(ji,jj) )
2613           ! For now try just set hmle to zmld
2614           hmle(ji,jj) = zmld(ji,jj)
2615        END DO
2616      END DO
2617
2618      mld_prof = 4
2619      DO jk = 5, jpkm1
2620        DO jj = 2, jpjm1
2621          DO ji = 2, jpim1
2622            IF ( hmle(ji,jj) >= gdepw_n(ji,jj,jk) ) mld_prof(ji,jj) = MIN(mbkt(ji,jj), jk)
2623          END DO
2624        END DO
2625      END DO
2626      DO jj = 2, jpjm1
2627         DO ji = 2, jpim1
2628            zhmle(ji,jj) = gdepw_n(ji,jj, mld_prof(ji,jj))
2629         END DO
2630       END DO
2631END SUBROUTINE zdf_osm_mle_parameters
2632
2633END SUBROUTINE zdf_osm
2634
2635
2636   SUBROUTINE zdf_osm_init
2637     !!----------------------------------------------------------------------
2638     !!                  ***  ROUTINE zdf_osm_init  ***
2639     !!
2640     !! ** Purpose :   Initialization of the vertical eddy diffivity and
2641     !!      viscosity when using a osm turbulent closure scheme
2642     !!
2643     !! ** Method  :   Read the namosm namelist and check the parameters
2644     !!      called at the first timestep (nit000)
2645     !!
2646     !! ** input   :   Namlist namosm
2647     !!----------------------------------------------------------------------
2648     INTEGER  ::   ios            ! local integer
2649     INTEGER  ::   ji, jj, jk     ! dummy loop indices
2650     REAL z1_t2
2651     !!
2652     NAMELIST/namzdf_osm/ ln_use_osm_la, rn_osm_la, rn_osm_dstokes, nn_ave &
2653          & ,nn_osm_wave, ln_dia_osm, rn_osm_hbl0, rn_zdfosm_adjust_sd &
2654          & ,ln_kpprimix, rn_riinfty, rn_difri, ln_convmix, rn_difconv, nn_osm_wave &
2655          & ,nn_osm_SD_reduce, ln_osm_mle, rn_osm_hblfrac, rn_osm_bl_thresh, ln_zdfosm_ice_shelter
2656! Namelist for Fox-Kemper parametrization.
2657      NAMELIST/namosm_mle/ nn_osm_mle, rn_osm_mle_ce, rn_osm_mle_lf, rn_osm_mle_time, rn_osm_mle_lat,&
2658           & rn_osm_mle_rho_c, rn_osm_mle_thresh, rn_osm_mle_tau, ln_osm_hmle_limit, rn_osm_hmle_limit
2659
2660     !!----------------------------------------------------------------------
2661     !
2662     REWIND( numnam_ref )              ! Namelist namzdf_osm in reference namelist : Osmosis ML model
2663     READ  ( numnam_ref, namzdf_osm, IOSTAT = ios, ERR = 901)
2664901  IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_osm in reference namelist' )
2665
2666     REWIND( numnam_cfg )              ! Namelist namzdf_tke in configuration namelist : Turbulent Kinetic Energy
2667     READ  ( numnam_cfg, namzdf_osm, IOSTAT = ios, ERR = 902 )
2668902  IF( ios >  0 ) CALL ctl_nam ( ios , 'namzdf_osm in configuration namelist' )
2669     IF(lwm) WRITE ( numond, namzdf_osm )
2670
2671     IF(lwp) THEN                    ! Control print
2672        WRITE(numout,*)
2673        WRITE(numout,*) 'zdf_osm_init : OSMOSIS Parameterisation'
2674        WRITE(numout,*) '~~~~~~~~~~~~'
2675        WRITE(numout,*) '   Namelist namzdf_osm : set osm mixing parameters'
2676        WRITE(numout,*) '     Use  rn_osm_la                                ln_use_osm_la = ', ln_use_osm_la
2677        WRITE(numout,*) '     Use  MLE in OBL, i.e. Fox-Kemper param        ln_osm_mle = ', ln_osm_mle
2678        WRITE(numout,*) '     Turbulent Langmuir number                     rn_osm_la   = ', rn_osm_la
2679        WRITE(numout,*) '     Stokes drift reduction factor                 rn_zdfosm_adjust_sd   = ', rn_zdfosm_adjust_sd
2680        WRITE(numout,*) '     Initial hbl for 1D runs                       rn_osm_hbl0   = ', rn_osm_hbl0
2681        WRITE(numout,*) '     Depth scale of Stokes drift                   rn_osm_dstokes = ', rn_osm_dstokes
2682        WRITE(numout,*) '     horizontal average flag                       nn_ave      = ', nn_ave
2683        WRITE(numout,*) '     Stokes drift                                  nn_osm_wave = ', nn_osm_wave
2684        SELECT CASE (nn_osm_wave)
2685        CASE(0)
2686           WRITE(numout,*) '     calculated assuming constant La#=0.3'
2687        CASE(1)
2688           WRITE(numout,*) '     calculated from Pierson Moskowitz wind-waves'
2689        CASE(2)
2690           WRITE(numout,*) '     calculated from ECMWF wave fields'
2691         END SELECT
2692        WRITE(numout,*) '     Stokes drift reduction                        nn_osm_SD_reduce', nn_osm_SD_reduce
2693        WRITE(numout,*) '     fraction of hbl to average SD over/fit'
2694        WRITE(numout,*) '     exponential with nn_osm_SD_reduce = 1 or 2    rn_osm_hblfrac =  ', rn_osm_hblfrac
2695        SELECT CASE (nn_osm_SD_reduce)
2696        CASE(0)
2697           WRITE(numout,*) '     No reduction'
2698        CASE(1)
2699           WRITE(numout,*) '     Average SD over upper rn_osm_hblfrac of BL'
2700        CASE(2)
2701           WRITE(numout,*) '     Fit exponential to slope rn_osm_hblfrac of BL'
2702        END SELECT
2703        WRITE(numout,*) '     reduce surface SD and depth scale under ice   ln_zdfosm_ice_shelter=', ln_zdfosm_ice_shelter
2704        WRITE(numout,*) '     Output osm diagnostics                       ln_dia_osm  = ',  ln_dia_osm
2705        WRITE(numout,*) '         Threshold used to define BL              rn_osm_bl_thresh  = ', rn_osm_bl_thresh, 'm^2/s'
2706        WRITE(numout,*) '     Use KPP-style shear instability mixing       ln_kpprimix = ', ln_kpprimix
2707        WRITE(numout,*) '     local Richardson Number limit for shear instability rn_riinfty = ', rn_riinfty
2708        WRITE(numout,*) '     maximum shear diffusivity at Rig = 0    (m2/s) rn_difri = ', rn_difri
2709        WRITE(numout,*) '     Use large mixing below BL when unstable       ln_convmix = ', ln_convmix
2710        WRITE(numout,*) '     diffusivity when unstable below BL     (m2/s) rn_difconv = ', rn_difconv
2711     ENDIF
2712
2713
2714     !                              ! Check wave coupling settings !
2715     !                         ! Further work needed - see ticket #2447 !
2716     IF( nn_osm_wave == 2 ) THEN
2717        IF (.NOT. ( ln_wave .AND. ln_sdw )) &
2718           & CALL ctl_stop( 'zdf_osm_init : ln_zdfosm and nn_osm_wave=2, ln_wave and ln_sdw must be true' )
2719     END IF
2720
2721     !                              ! allocate zdfosm arrays
2722     IF( zdf_osm_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_osm_init : unable to allocate arrays' )
2723
2724
2725     IF( ln_osm_mle ) THEN
2726! Initialise Fox-Kemper parametrization
2727         REWIND( numnam_ref )              ! Namelist namosm_mle in reference namelist : Tracer advection scheme
2728         READ  ( numnam_ref, namosm_mle, IOSTAT = ios, ERR = 903)
2729903      IF( ios /= 0 )   CALL ctl_nam ( ios , 'namosm_mle in reference namelist')
2730
2731         REWIND( numnam_cfg )              ! Namelist namosm_mle in configuration namelist : Tracer advection scheme
2732         READ  ( numnam_cfg, namosm_mle, IOSTAT = ios, ERR = 904 )
2733904      IF( ios >  0 )   CALL ctl_nam ( ios , 'namosm_mle in configuration namelist')
2734         IF(lwm) WRITE ( numond, namosm_mle )
2735
2736         IF(lwp) THEN                     ! Namelist print
2737            WRITE(numout,*)
2738            WRITE(numout,*) 'zdf_osm_init : initialise mixed layer eddy (MLE)'
2739            WRITE(numout,*) '~~~~~~~~~~~~~'
2740            WRITE(numout,*) '   Namelist namosm_mle : '
2741            WRITE(numout,*) '         MLE type: =0 standard Fox-Kemper ; =1 new formulation        nn_osm_mle    = ', nn_osm_mle
2742            WRITE(numout,*) '         magnitude of the MLE (typical value: 0.06 to 0.08)           rn_osm_mle_ce    = ', rn_osm_mle_ce
2743            WRITE(numout,*) '         scale of ML front (ML radius of deformation) (nn_osm_mle=0)      rn_osm_mle_lf     = ', rn_osm_mle_lf, 'm'
2744            WRITE(numout,*) '         maximum time scale of MLE                    (nn_osm_mle=0)      rn_osm_mle_time   = ', rn_osm_mle_time, 's'
2745            WRITE(numout,*) '         reference latitude (degrees) of MLE coef.    (nn_osm_mle=1)      rn_osm_mle_lat    = ', rn_osm_mle_lat, 'deg'
2746            WRITE(numout,*) '         Density difference used to define ML for FK              rn_osm_mle_rho_c  = ', rn_osm_mle_rho_c
2747            WRITE(numout,*) '         Threshold used to define MLE for FK                      rn_osm_mle_thresh  = ', rn_osm_mle_thresh, 'm^2/s'
2748            WRITE(numout,*) '         Timescale for OSM-FK                                         rn_osm_mle_tau  = ', rn_osm_mle_tau, 's'
2749            WRITE(numout,*) '         switch to limit hmle                                      ln_osm_hmle_limit  = ', ln_osm_hmle_limit
2750            WRITE(numout,*) '         fraction of zmld to limit hmle to if ln_osm_hmle_limit =.T.  rn_osm_hmle_limit = ', rn_osm_hmle_limit
2751         ENDIF         !
2752     ENDIF
2753      !
2754      IF(lwp) THEN
2755         WRITE(numout,*)
2756         IF( ln_osm_mle ) THEN
2757            WRITE(numout,*) '   ==>>>   Mixed Layer Eddy induced transport added to OSMOSIS BL calculation'
2758            IF( nn_osm_mle == 0 )   WRITE(numout,*) '              Fox-Kemper et al 2010 formulation'
2759            IF( nn_osm_mle == 1 )   WRITE(numout,*) '              New formulation'
2760         ELSE
2761            WRITE(numout,*) '   ==>>>   Mixed Layer induced transport NOT added to OSMOSIS BL calculation'
2762         ENDIF
2763      ENDIF
2764      !
2765      IF( ln_osm_mle ) THEN                ! MLE initialisation
2766         !
2767         rb_c = grav * rn_osm_mle_rho_c /rau0        ! Mixed Layer buoyancy criteria
2768         IF(lwp) WRITE(numout,*)
2769         IF(lwp) WRITE(numout,*) '      ML buoyancy criteria = ', rb_c, ' m/s2 '
2770         IF(lwp) WRITE(numout,*) '      associated ML density criteria defined in zdfmxl = ', rn_osm_mle_rho_c, 'kg/m3'
2771         !
2772         IF( nn_osm_mle == 0 ) THEN           ! MLE array allocation & initialisation            !
2773!
2774         ELSEIF( nn_osm_mle == 1 ) THEN           ! MLE array allocation & initialisation
2775            rc_f = rn_osm_mle_ce/ (  5.e3_wp * 2._wp * omega * SIN( rad * rn_osm_mle_lat )  )
2776            !
2777         ENDIF
2778         !                                ! 1/(f^2+tau^2)^1/2 at t-point (needed in both nn_osm_mle case)
2779         z1_t2 = 2.e-5
2780         do jj=1,jpj
2781            do ji = 1,jpi
2782               r1_ft(ji,jj) = MIN(1./( ABS(ff_t(ji,jj)) + epsln ), ABS(ff_t(ji,jj))/z1_t2**2)
2783            end do
2784         end do
2785         ! z1_t2 = 1._wp / ( rn_osm_mle_time * rn_osm_mle_timeji,jj )
2786         ! r1_ft(:,:) = 1._wp / SQRT(  ff_t(:,:) * ff_t(:,:) + z1_t2  )
2787         !
2788      ENDIF
2789
2790     call osm_rst( nit000, 'READ' ) !* read or initialize hbl, dh, hmle
2791
2792
2793     IF( ln_zdfddm) THEN
2794        IF(lwp) THEN
2795           WRITE(numout,*)
2796           WRITE(numout,*) '    Double diffusion mixing on temperature and salinity '
2797           WRITE(numout,*) '    CAUTION : done in routine zdfosm, not in routine zdfddm '
2798        ENDIF
2799     ENDIF
2800
2801
2802     !set constants not in namelist
2803     !-----------------------------
2804
2805     IF(lwp) THEN
2806        WRITE(numout,*)
2807     ENDIF
2808
2809     IF (nn_osm_wave == 0) THEN
2810        dstokes(:,:) = rn_osm_dstokes
2811     END IF
2812
2813     ! Horizontal average : initialization of weighting arrays
2814     ! -------------------
2815
2816     SELECT CASE ( nn_ave )
2817
2818     CASE ( 0 )                ! no horizontal average
2819        IF(lwp) WRITE(numout,*) '          no horizontal average on avt'
2820        IF(lwp) WRITE(numout,*) '          only in very high horizontal resolution !'
2821        ! weighting mean arrays etmean
2822        !           ( 1  1 )
2823        ! avt = 1/4 ( 1  1 )
2824        !
2825        etmean(:,:,:) = 0.e0
2826
2827        DO jk = 1, jpkm1
2828           DO jj = 2, jpjm1
2829              DO ji = 2, jpim1   ! vector opt.
2830                 etmean(ji,jj,jk) = tmask(ji,jj,jk)                     &
2831                      &  / MAX( 1.,  umask(ji-1,jj  ,jk) + umask(ji,jj,jk)   &
2832                      &            + vmask(ji  ,jj-1,jk) + vmask(ji,jj,jk)  )
2833              END DO
2834           END DO
2835        END DO
2836
2837     CASE ( 1 )                ! horizontal average
2838        IF(lwp) WRITE(numout,*) '          horizontal average on avt'
2839        ! weighting mean arrays etmean
2840        !           ( 1/2  1  1/2 )
2841        ! avt = 1/8 ( 1    2  1   )
2842        !           ( 1/2  1  1/2 )
2843        etmean(:,:,:) = 0.e0
2844
2845        DO jk = 1, jpkm1
2846           DO jj = 2, jpjm1
2847              DO ji = 2, jpim1   ! vector opt.
2848                 etmean(ji,jj,jk) = tmask(ji, jj,jk)                           &
2849                      & / MAX( 1., 2.* tmask(ji,jj,jk)                           &
2850                      &      +.5 * ( tmask(ji-1,jj+1,jk) + tmask(ji-1,jj-1,jk)   &
2851                      &             +tmask(ji+1,jj+1,jk) + tmask(ji+1,jj-1,jk) ) &
2852                      &      +1. * ( tmask(ji-1,jj  ,jk) + tmask(ji  ,jj+1,jk)   &
2853                      &             +tmask(ji  ,jj-1,jk) + tmask(ji+1,jj  ,jk) ) )
2854              END DO
2855           END DO
2856        END DO
2857
2858     CASE DEFAULT
2859        WRITE(ctmp1,*) '          bad flag value for nn_ave = ', nn_ave
2860        CALL ctl_stop( ctmp1 )
2861
2862     END SELECT
2863
2864     ! Initialization of vertical eddy coef. to the background value
2865     ! -------------------------------------------------------------
2866     DO jk = 1, jpk
2867        avt (:,:,jk) = avtb(jk) * tmask(:,:,jk)
2868     END DO
2869
2870     ! zero the surface flux for non local term and osm mixed layer depth
2871     ! ------------------------------------------------------------------
2872     ghamt(:,:,:) = 0.
2873     ghams(:,:,:) = 0.
2874     ghamu(:,:,:) = 0.
2875     ghamv(:,:,:) = 0.
2876     !
2877     IF( lwxios ) THEN
2878        CALL iom_set_rstw_var_active('wn')
2879        CALL iom_set_rstw_var_active('hbl')
2880        CALL iom_set_rstw_var_active('dh')
2881        IF( ln_osm_mle ) THEN
2882            CALL iom_set_rstw_var_active('hmle')
2883        END IF
2884     ENDIF
2885   END SUBROUTINE zdf_osm_init
2886
2887
2888   SUBROUTINE osm_rst( kt, cdrw )
2889     !!---------------------------------------------------------------------
2890     !!                   ***  ROUTINE osm_rst  ***
2891     !!
2892     !! ** Purpose :   Read or write BL fields in restart file
2893     !!
2894     !! ** Method  :   use of IOM library. If the restart does not contain
2895     !!                required fields, they are recomputed from stratification
2896     !!----------------------------------------------------------------------
2897
2898     INTEGER, INTENT(in) :: kt
2899     CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
2900
2901     INTEGER ::   id1, id2, id3   ! iom enquiry index
2902     INTEGER  ::   ji, jj, jk     ! dummy loop indices
2903     INTEGER  ::   iiki, ikt ! local integer
2904     REAL(wp) ::   zhbf           ! tempory scalars
2905     REAL(wp) ::   zN2_c           ! local scalar
2906     REAL(wp) ::   rho_c = 0.01_wp    !: density criterion for mixed layer depth
2907     INTEGER, DIMENSION(jpi,jpj) :: imld_rst ! level of mixed-layer depth (pycnocline top)
2908     !!----------------------------------------------------------------------
2909     !
2910     !!-----------------------------------------------------------------------------
2911     ! If READ/WRITE Flag is 'READ', try to get hbl from restart file. If successful then return
2912     !!-----------------------------------------------------------------------------
2913     IF( TRIM(cdrw) == 'READ'.AND. ln_rstart) THEN
2914        id1 = iom_varid( numror, 'wn'   , ldstop = .FALSE. )
2915        IF( id1 > 0 ) THEN                       ! 'wn' exists; read
2916           CALL iom_get( numror, jpdom_autoglo, 'wn', wn, ldxios = lrxios )
2917           WRITE(numout,*) ' ===>>>> :  wn read from restart file'
2918        ELSE
2919           wn(:,:,:) = 0._wp
2920           WRITE(numout,*) ' ===>>>> :  wn not in restart file, set to zero initially'
2921        END IF
2922
2923        id1 = iom_varid( numror, 'hbl'   , ldstop = .FALSE. )
2924        id2 = iom_varid( numror, 'dh'   , ldstop = .FALSE. )
2925        IF( id1 > 0 .AND. id2 > 0) THEN                       ! 'hbl' exists; read and return
2926           CALL iom_get( numror, jpdom_autoglo, 'hbl' , hbl , ldxios = lrxios )
2927           CALL iom_get( numror, jpdom_autoglo, 'dh', dh, ldxios = lrxios  )
2928           WRITE(numout,*) ' ===>>>> :  hbl & dh read from restart file'
2929           IF( ln_osm_mle ) THEN
2930              id3 = iom_varid( numror, 'hmle'   , ldstop = .FALSE. )
2931              IF( id3 > 0) THEN
2932                 CALL iom_get( numror, jpdom_autoglo, 'hmle' , hmle , ldxios = lrxios )
2933                 WRITE(numout,*) ' ===>>>> :  hmle read from restart file'
2934              ELSE
2935                 WRITE(numout,*) ' ===>>>> :  hmle not found, set to hbl'
2936                 hmle(:,:) = hbl(:,:)            ! Initialise MLE depth.
2937              END IF
2938           END IF
2939           RETURN
2940        ELSE                      ! 'hbl' & 'dh' not in restart file, recalculate
2941           WRITE(numout,*) ' ===>>>> : previous run without osmosis scheme, hbl computed from stratification'
2942        END IF
2943     END IF
2944
2945     !!-----------------------------------------------------------------------------
2946     ! If READ/WRITE Flag is 'WRITE', write hbl into the restart file, then return
2947     !!-----------------------------------------------------------------------------
2948     IF( TRIM(cdrw) == 'WRITE') THEN     !* Write hbli into the restart file, then return
2949        IF(lwp) WRITE(numout,*) '---- osm-rst ----'
2950         CALL iom_rstput( kt, nitrst, numrow, 'wn'     , wn,   ldxios = lwxios )
2951         CALL iom_rstput( kt, nitrst, numrow, 'hbl'    , hbl,  ldxios = lwxios )
2952         CALL iom_rstput( kt, nitrst, numrow, 'dh'     , dh,   ldxios = lwxios )
2953         IF( ln_osm_mle ) THEN
2954            CALL iom_rstput( kt, nitrst, numrow, 'hmle', hmle, ldxios = lwxios )
2955         END IF
2956        RETURN
2957     END IF
2958
2959     !!-----------------------------------------------------------------------------
2960     ! Getting hbl, no restart file with hbl, so calculate from surface stratification
2961     !!-----------------------------------------------------------------------------
2962     IF( lwp ) WRITE(numout,*) ' ===>>>> : calculating hbl computed from stratification'
2963     ! w-level of the mixing and mixed layers
2964     CALL eos_rab( tsn, rab_n )
2965     CALL bn2(tsn, rab_n, rn2)
2966     imld_rst(:,:)  = nlb10         ! Initialization to the number of w ocean point
2967     hbl(:,:)  = 0._wp              ! here hbl used as a dummy variable, integrating vertically N^2
2968     zN2_c = grav * rho_c * r1_rau0 ! convert density criteria into N^2 criteria
2969     !
2970     hbl(:,:)  = 0._wp              ! here hbl used as a dummy variable, integrating vertically N^2
2971     DO jk = 1, jpkm1
2972        DO jj = 1, jpj              ! Mixed layer level: w-level
2973           DO ji = 1, jpi
2974              ikt = mbkt(ji,jj)
2975              hbl(ji,jj) = hbl(ji,jj) + MAX( rn2(ji,jj,jk) , 0._wp ) * e3w_n(ji,jj,jk)
2976              IF( hbl(ji,jj) < zN2_c )   imld_rst(ji,jj) = MIN( jk , ikt ) + 1   ! Mixed layer level
2977           END DO
2978        END DO
2979     END DO
2980     !
2981     DO jj = 1, jpj
2982        DO ji = 1, jpi
2983           iiki = MAX(4,imld_rst(ji,jj))
2984           hbl (ji,jj) = gdepw_n(ji,jj,iiki  )    ! Turbocline depth
2985           dh (ji,jj) = e3t_n(ji,jj,iiki-1  )     ! Turbocline depth
2986        END DO
2987     END DO
2988
2989     WRITE(numout,*) ' ===>>>> : hbl computed from stratification'
2990
2991     IF( ln_osm_mle ) THEN
2992        hmle(:,:) = hbl(:,:)            ! Initialise MLE depth.
2993        WRITE(numout,*) ' ===>>>> : hmle set = to hbl'
2994     END IF
2995
2996     wn(:,:,:) = 0._wp
2997     WRITE(numout,*) ' ===>>>> :  wn not in restart file, set to zero initially'
2998   END SUBROUTINE osm_rst
2999
3000
3001   SUBROUTINE tra_osm( kt )
3002      !!----------------------------------------------------------------------
3003      !!                  ***  ROUTINE tra_osm  ***
3004      !!
3005      !! ** Purpose :   compute and add to the tracer trend the non-local tracer flux
3006      !!
3007      !! ** Method  :   ???
3008      !!----------------------------------------------------------------------
3009      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   ztrdt, ztrds   ! 3D workspace
3010      !!----------------------------------------------------------------------
3011      INTEGER, INTENT(in) :: kt
3012      INTEGER :: ji, jj, jk
3013      !
3014      IF( kt == nit000 ) THEN
3015         IF(lwp) WRITE(numout,*)
3016         IF(lwp) WRITE(numout,*) 'tra_osm : OSM non-local tracer fluxes'
3017         IF(lwp) WRITE(numout,*) '~~~~~~~   '
3018      ENDIF
3019
3020      IF( l_trdtra )   THEN                    !* Save ta and sa trends
3021         ALLOCATE( ztrdt(jpi,jpj,jpk) )   ;    ztrdt(:,:,:) = tsa(:,:,:,jp_tem)
3022         ALLOCATE( ztrds(jpi,jpj,jpk) )   ;    ztrds(:,:,:) = tsa(:,:,:,jp_sal)
3023      ENDIF
3024
3025      DO jk = 1, jpkm1
3026         DO jj = 2, jpjm1
3027            DO ji = 2, jpim1
3028               tsa(ji,jj,jk,jp_tem) =  tsa(ji,jj,jk,jp_tem)                      &
3029                  &                 - (  ghamt(ji,jj,jk  )  &
3030                  &                    - ghamt(ji,jj,jk+1) ) /e3t_n(ji,jj,jk)
3031               tsa(ji,jj,jk,jp_sal) =  tsa(ji,jj,jk,jp_sal)                      &
3032                  &                 - (  ghams(ji,jj,jk  )  &
3033                  &                    - ghams(ji,jj,jk+1) ) / e3t_n(ji,jj,jk)
3034            END DO
3035         END DO
3036      END DO
3037
3038      ! save the non-local tracer flux trends for diagnostics
3039      IF( l_trdtra )   THEN
3040         ztrdt(:,:,:) = tsa(:,:,:,jp_tem) - ztrdt(:,:,:)
3041         ztrds(:,:,:) = tsa(:,:,:,jp_sal) - ztrds(:,:,:)
3042
3043         CALL trd_tra( kt, 'TRA', jp_tem, jptra_osm, ztrdt )
3044         CALL trd_tra( kt, 'TRA', jp_sal, jptra_osm, ztrds )
3045         DEALLOCATE( ztrdt )      ;     DEALLOCATE( ztrds )
3046      ENDIF
3047
3048      IF(ln_ctl) THEN
3049         CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' osm  - Ta: ', mask1=tmask,   &
3050         &             tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
3051      ENDIF
3052      !
3053   END SUBROUTINE tra_osm
3054
3055
3056   SUBROUTINE trc_osm( kt )          ! Dummy routine
3057      !!----------------------------------------------------------------------
3058      !!                  ***  ROUTINE trc_osm  ***
3059      !!
3060      !! ** Purpose :   compute and add to the passive tracer trend the non-local
3061      !!                 passive tracer flux
3062      !!
3063      !!
3064      !! ** Method  :   ???
3065      !!----------------------------------------------------------------------
3066      !
3067      !!----------------------------------------------------------------------
3068      INTEGER, INTENT(in) :: kt
3069      WRITE(*,*) 'trc_osm: Not written yet', kt
3070   END SUBROUTINE trc_osm
3071
3072
3073   SUBROUTINE dyn_osm( kt )
3074      !!----------------------------------------------------------------------
3075      !!                  ***  ROUTINE dyn_osm  ***
3076      !!
3077      !! ** Purpose :   compute and add to the velocity trend the non-local flux
3078      !! copied/modified from tra_osm
3079      !!
3080      !! ** Method  :   ???
3081      !!----------------------------------------------------------------------
3082      INTEGER, INTENT(in) ::   kt   !
3083      !
3084      INTEGER :: ji, jj, jk   ! dummy loop indices
3085      !!----------------------------------------------------------------------
3086      !
3087      IF( kt == nit000 ) THEN
3088         IF(lwp) WRITE(numout,*)
3089         IF(lwp) WRITE(numout,*) 'dyn_osm : OSM non-local velocity'
3090         IF(lwp) WRITE(numout,*) '~~~~~~~   '
3091      ENDIF
3092      !code saving tracer trends removed, replace with trdmxl_oce
3093
3094      DO jk = 1, jpkm1           ! add non-local u and v fluxes
3095         DO jj = 2, jpjm1
3096            DO ji = 2, jpim1
3097               ua(ji,jj,jk) =  ua(ji,jj,jk)                      &
3098                  &                 - (  ghamu(ji,jj,jk  )  &
3099                  &                    - ghamu(ji,jj,jk+1) ) / e3u_n(ji,jj,jk)
3100               va(ji,jj,jk) =  va(ji,jj,jk)                      &
3101                  &                 - (  ghamv(ji,jj,jk  )  &
3102                  &                    - ghamv(ji,jj,jk+1) ) / e3v_n(ji,jj,jk)
3103            END DO
3104         END DO
3105      END DO
3106      !
3107      ! code for saving tracer trends removed
3108      !
3109   END SUBROUTINE dyn_osm
3110
3111   !!======================================================================
3112
3113END MODULE zdfosm
Note: See TracBrowser for help on using the repository browser.