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