[1531] | 1 | MODULE zdftke |
---|
[1239] | 2 | !!====================================================================== |
---|
[1531] | 3 | !! *** MODULE zdftke *** |
---|
[1239] | 4 | !! Ocean physics: vertical mixing coefficient computed from the tke |
---|
| 5 | !! turbulent closure parameterization |
---|
| 6 | !!===================================================================== |
---|
[1492] | 7 | !! History : OPA ! 1991-03 (b. blanke) Original code |
---|
| 8 | !! 7.0 ! 1991-11 (G. Madec) bug fix |
---|
| 9 | !! 7.1 ! 1992-10 (G. Madec) new mixing length and eav |
---|
| 10 | !! 7.2 ! 1993-03 (M. Guyon) symetrical conditions |
---|
| 11 | !! 7.3 ! 1994-08 (G. Madec, M. Imbard) nn_pdl flag |
---|
| 12 | !! 7.5 ! 1996-01 (G. Madec) s-coordinates |
---|
| 13 | !! 8.0 ! 1997-07 (G. Madec) lbc |
---|
| 14 | !! 8.1 ! 1999-01 (E. Stretta) new option for the mixing length |
---|
| 15 | !! NEMO 1.0 ! 2002-06 (G. Madec) add tke_init routine |
---|
| 16 | !! - ! 2004-10 (C. Ethe ) 1D configuration |
---|
| 17 | !! 2.0 ! 2006-07 (S. Masson) distributed restart using iom |
---|
| 18 | !! 3.0 ! 2008-05 (C. Ethe, G.Madec) : update TKE physics: |
---|
| 19 | !! ! - tke penetration (wind steering) |
---|
| 20 | !! ! - suface condition for tke & mixing length |
---|
| 21 | !! ! - Langmuir cells |
---|
| 22 | !! - ! 2008-05 (J.-M. Molines, G. Madec) 2D form of avtb |
---|
| 23 | !! - ! 2008-06 (G. Madec) style + DOCTOR name for namelist parameters |
---|
| 24 | !! - ! 2008-12 (G. Reffray) stable discretization of the production term |
---|
| 25 | !! 3.2 ! 2009-06 (G. Madec, S. Masson) TKE restart compatible with key_cpl |
---|
| 26 | !! ! + cleaning of the parameters + bugs correction |
---|
[2528] | 27 | !! 3.3 ! 2010-10 (C. Ethe, G. Madec) reorganisation of initialisation phase |
---|
[5120] | 28 | !! 3.6 ! 2014-11 (P. Mathiot) add ice shelf capability |
---|
[1239] | 29 | !!---------------------------------------------------------------------- |
---|
[5836] | 30 | #if defined key_zdftke |
---|
[1239] | 31 | !!---------------------------------------------------------------------- |
---|
[1531] | 32 | !! 'key_zdftke' TKE vertical physics |
---|
[1239] | 33 | !!---------------------------------------------------------------------- |
---|
[3625] | 34 | !! zdf_tke : update momentum and tracer Kz from a tke scheme |
---|
| 35 | !! tke_tke : tke time stepping: update tke at now time step (en) |
---|
| 36 | !! tke_avn : compute mixing length scale and deduce avm and avt |
---|
| 37 | !! zdf_tke_init : initialization, namelist read, and parameters control |
---|
| 38 | !! tke_rst : read/write tke restart in ocean restart file |
---|
[1239] | 39 | !!---------------------------------------------------------------------- |
---|
[2528] | 40 | USE oce ! ocean: dynamics and active tracers variables |
---|
| 41 | USE phycst ! physical constants |
---|
| 42 | USE dom_oce ! domain: ocean |
---|
| 43 | USE domvvl ! domain: variable volume layer |
---|
[1492] | 44 | USE sbc_oce ! surface boundary condition: ocean |
---|
[2528] | 45 | USE zdf_oce ! vertical physics: ocean variables |
---|
| 46 | USE zdfmxl ! vertical physics: mixed layer |
---|
[1492] | 47 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
| 48 | USE prtctl ! Print control |
---|
| 49 | USE in_out_manager ! I/O manager |
---|
| 50 | USE iom ! I/O manager library |
---|
[2715] | 51 | USE lib_mpp ! MPP library |
---|
[3294] | 52 | USE wrk_nemo ! work arrays |
---|
| 53 | USE timing ! Timing |
---|
[3625] | 54 | USE lib_fortran ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) |
---|
[5656] | 55 | #if defined key_agrif |
---|
| 56 | USE agrif_opa_interp |
---|
| 57 | USE agrif_opa_update |
---|
| 58 | #endif |
---|
[1239] | 59 | |
---|
| 60 | IMPLICIT NONE |
---|
| 61 | PRIVATE |
---|
| 62 | |
---|
[2528] | 63 | PUBLIC zdf_tke ! routine called in step module |
---|
| 64 | PUBLIC zdf_tke_init ! routine called in opa module |
---|
| 65 | PUBLIC tke_rst ! routine called in step module |
---|
[1239] | 66 | |
---|
[2715] | 67 | LOGICAL , PUBLIC, PARAMETER :: lk_zdftke = .TRUE. !: TKE vertical mixing flag |
---|
[1239] | 68 | |
---|
[4147] | 69 | ! !!** Namelist namzdf_tke ** |
---|
| 70 | LOGICAL :: ln_mxl0 ! mixing length scale surface value as function of wind stress or not |
---|
| 71 | INTEGER :: nn_mxl ! type of mixing length (=0/1/2/3) |
---|
| 72 | REAL(wp) :: rn_mxl0 ! surface min value of mixing length (kappa*z_o=0.4*0.1 m) [m] |
---|
| 73 | INTEGER :: nn_pdl ! Prandtl number or not (ratio avt/avm) (=0/1) |
---|
| 74 | REAL(wp) :: rn_ediff ! coefficient for avt: avt=rn_ediff*mxl*sqrt(e) |
---|
| 75 | REAL(wp) :: rn_ediss ! coefficient of the Kolmogoroff dissipation |
---|
| 76 | REAL(wp) :: rn_ebb ! coefficient of the surface input of tke |
---|
| 77 | REAL(wp) :: rn_emin ! minimum value of tke [m2/s2] |
---|
| 78 | REAL(wp) :: rn_emin0 ! surface minimum value of tke [m2/s2] |
---|
| 79 | REAL(wp) :: rn_bshear ! background shear (>0) currently a numerical threshold (do not change it) |
---|
| 80 | INTEGER :: nn_etau ! type of depth penetration of surface tke (=0/1/2/3) |
---|
| 81 | INTEGER :: nn_htau ! type of tke profile of penetration (=0/1) |
---|
| 82 | REAL(wp) :: rn_efr ! fraction of TKE surface value which penetrates in the ocean |
---|
| 83 | LOGICAL :: ln_lc ! Langmuir cells (LC) as a source term of TKE or not |
---|
| 84 | REAL(wp) :: rn_lc ! coef to compute vertical velocity of Langmuir cells |
---|
[1239] | 85 | |
---|
[4147] | 86 | REAL(wp) :: ri_cri ! critic Richardson number (deduced from rn_ediff and rn_ediss values) |
---|
| 87 | REAL(wp) :: rmxl_min ! minimum mixing length value (deduced from rn_ediff and rn_emin values) [m] |
---|
[2528] | 88 | REAL(wp) :: rhftau_add = 1.e-3_wp ! add offset applied to HF part of taum (nn_etau=3) |
---|
| 89 | REAL(wp) :: rhftau_scl = 1.0_wp ! scale factor applied to HF part of taum (nn_etau=3) |
---|
[1239] | 90 | |
---|
[2715] | 91 | REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:) :: htau ! depth of tke penetration (nn_htau) |
---|
| 92 | REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: dissl ! now mixing lenght of dissipation |
---|
[5656] | 93 | REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: apdlr ! now mixing lenght of dissipation |
---|
[2715] | 94 | #if defined key_c1d |
---|
| 95 | ! !!** 1D cfg only ** ('key_c1d') |
---|
| 96 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: e_dis, e_mix !: dissipation and mixing turbulent lengh scales |
---|
| 97 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: e_pdl, e_ric !: prandl and local Richardson numbers |
---|
| 98 | #endif |
---|
[1492] | 99 | |
---|
[1239] | 100 | !! * Substitutions |
---|
| 101 | # include "vectopt_loop_substitute.h90" |
---|
| 102 | !!---------------------------------------------------------------------- |
---|
[5836] | 103 | !! NEMO/OPA 3.7 , NEMO Consortium (2015) |
---|
[2528] | 104 | !! $Id$ |
---|
| 105 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
[1239] | 106 | !!---------------------------------------------------------------------- |
---|
| 107 | CONTAINS |
---|
| 108 | |
---|
[2715] | 109 | INTEGER FUNCTION zdf_tke_alloc() |
---|
| 110 | !!---------------------------------------------------------------------- |
---|
| 111 | !! *** FUNCTION zdf_tke_alloc *** |
---|
| 112 | !!---------------------------------------------------------------------- |
---|
| 113 | ALLOCATE( & |
---|
| 114 | #if defined key_c1d |
---|
| 115 | & e_dis(jpi,jpj,jpk) , e_mix(jpi,jpj,jpk) , & |
---|
| 116 | & e_pdl(jpi,jpj,jpk) , e_ric(jpi,jpj,jpk) , & |
---|
| 117 | #endif |
---|
[5836] | 118 | & htau (jpi,jpj) , dissl(jpi,jpj,jpk) , & |
---|
| 119 | & apdlr(jpi,jpj,jpk) , STAT= zdf_tke_alloc ) |
---|
[2715] | 120 | ! |
---|
| 121 | IF( lk_mpp ) CALL mpp_sum ( zdf_tke_alloc ) |
---|
| 122 | IF( zdf_tke_alloc /= 0 ) CALL ctl_warn('zdf_tke_alloc: failed to allocate arrays') |
---|
| 123 | ! |
---|
| 124 | END FUNCTION zdf_tke_alloc |
---|
| 125 | |
---|
| 126 | |
---|
[1531] | 127 | SUBROUTINE zdf_tke( kt ) |
---|
[1239] | 128 | !!---------------------------------------------------------------------- |
---|
[1531] | 129 | !! *** ROUTINE zdf_tke *** |
---|
[1239] | 130 | !! |
---|
| 131 | !! ** Purpose : Compute the vertical eddy viscosity and diffusivity |
---|
[1492] | 132 | !! coefficients using a turbulent closure scheme (TKE). |
---|
[1239] | 133 | !! |
---|
[1492] | 134 | !! ** Method : The time evolution of the turbulent kinetic energy (tke) |
---|
| 135 | !! is computed from a prognostic equation : |
---|
| 136 | !! d(en)/dt = avm (d(u)/dz)**2 ! shear production |
---|
| 137 | !! + d( avm d(en)/dz )/dz ! diffusion of tke |
---|
| 138 | !! + avt N^2 ! stratif. destruc. |
---|
| 139 | !! - rn_ediss / emxl en**(2/3) ! Kolmogoroff dissipation |
---|
[1239] | 140 | !! with the boundary conditions: |
---|
[1695] | 141 | !! surface: en = max( rn_emin0, rn_ebb * taum ) |
---|
[1239] | 142 | !! bottom : en = rn_emin |
---|
[1492] | 143 | !! The associated critical Richardson number is: ri_cri = 2/(2+rn_ediss/rn_ediff) |
---|
| 144 | !! |
---|
| 145 | !! The now Turbulent kinetic energy is computed using the following |
---|
| 146 | !! time stepping: implicit for vertical diffusion term, linearized semi |
---|
| 147 | !! implicit for kolmogoroff dissipation term, and explicit forward for |
---|
| 148 | !! both buoyancy and shear production terms. Therefore a tridiagonal |
---|
| 149 | !! linear system is solved. Note that buoyancy and shear terms are |
---|
| 150 | !! discretized in a energy conserving form (Bruchard 2002). |
---|
| 151 | !! |
---|
| 152 | !! The dissipative and mixing length scale are computed from en and |
---|
| 153 | !! the stratification (see tke_avn) |
---|
| 154 | !! |
---|
| 155 | !! The now vertical eddy vicosity and diffusivity coefficients are |
---|
| 156 | !! given by: |
---|
| 157 | !! avm = max( avtb, rn_ediff * zmxlm * en^1/2 ) |
---|
| 158 | !! avt = max( avmb, pdl * avm ) |
---|
[1239] | 159 | !! eav = max( avmb, avm ) |
---|
[1492] | 160 | !! where pdl, the inverse of the Prandtl number is 1 if nn_pdl=0 and |
---|
| 161 | !! given by an empirical funtion of the localRichardson number if nn_pdl=1 |
---|
[1239] | 162 | !! |
---|
| 163 | !! ** Action : compute en (now turbulent kinetic energy) |
---|
| 164 | !! update avt, avmu, avmv (before vertical eddy coef.) |
---|
| 165 | !! |
---|
| 166 | !! References : Gaspar et al., JGR, 1990, |
---|
| 167 | !! Blanke and Delecluse, JPO, 1991 |
---|
| 168 | !! Mellor and Blumberg, JPO 2004 |
---|
| 169 | !! Axell, JGR, 2002 |
---|
[1492] | 170 | !! Bruchard OM 2002 |
---|
[1239] | 171 | !!---------------------------------------------------------------------- |
---|
[1492] | 172 | INTEGER, INTENT(in) :: kt ! ocean time step |
---|
| 173 | !!---------------------------------------------------------------------- |
---|
[1481] | 174 | ! |
---|
[5656] | 175 | #if defined key_agrif |
---|
| 176 | ! interpolation parent grid => child grid for avm_k ( ex : at west border: update column 1 and 2) |
---|
| 177 | IF( .NOT.Agrif_Root() ) CALL Agrif_Tke |
---|
| 178 | #endif |
---|
| 179 | ! |
---|
[3632] | 180 | IF( kt /= nit000 ) THEN ! restore before value to compute tke |
---|
| 181 | avt (:,:,:) = avt_k (:,:,:) |
---|
| 182 | avm (:,:,:) = avm_k (:,:,:) |
---|
| 183 | avmu(:,:,:) = avmu_k(:,:,:) |
---|
| 184 | avmv(:,:,:) = avmv_k(:,:,:) |
---|
| 185 | ENDIF |
---|
| 186 | ! |
---|
[2528] | 187 | CALL tke_tke ! now tke (en) |
---|
[1492] | 188 | ! |
---|
[2528] | 189 | CALL tke_avn ! now avt, avm, avmu, avmv |
---|
| 190 | ! |
---|
[3632] | 191 | avt_k (:,:,:) = avt (:,:,:) |
---|
| 192 | avm_k (:,:,:) = avm (:,:,:) |
---|
| 193 | avmu_k(:,:,:) = avmu(:,:,:) |
---|
| 194 | avmv_k(:,:,:) = avmv(:,:,:) |
---|
| 195 | ! |
---|
[5656] | 196 | #if defined key_agrif |
---|
| 197 | ! Update child grid f => parent grid |
---|
| 198 | IF( .NOT.Agrif_Root() ) CALL Agrif_Update_Tke( kt ) ! children only |
---|
| 199 | #endif |
---|
| 200 | ! |
---|
| 201 | END SUBROUTINE zdf_tke |
---|
[1239] | 202 | |
---|
[1492] | 203 | |
---|
[1481] | 204 | SUBROUTINE tke_tke |
---|
[1239] | 205 | !!---------------------------------------------------------------------- |
---|
[1492] | 206 | !! *** ROUTINE tke_tke *** |
---|
| 207 | !! |
---|
| 208 | !! ** Purpose : Compute the now Turbulente Kinetic Energy (TKE) |
---|
| 209 | !! |
---|
| 210 | !! ** Method : - TKE surface boundary condition |
---|
[2528] | 211 | !! - source term due to Langmuir cells (Axell JGR 2002) (ln_lc=T) |
---|
[1492] | 212 | !! - source term due to shear (saved in avmu, avmv arrays) |
---|
| 213 | !! - Now TKE : resolution of the TKE equation by inverting |
---|
| 214 | !! a tridiagonal linear system by a "methode de chasse" |
---|
| 215 | !! - increase TKE due to surface and internal wave breaking |
---|
| 216 | !! |
---|
| 217 | !! ** Action : - en : now turbulent kinetic energy) |
---|
| 218 | !! - avmu, avmv : production of TKE by shear at u and v-points |
---|
| 219 | !! (= Kz dz[Ub] * dz[Un] ) |
---|
[1239] | 220 | !! --------------------------------------------------------------------- |
---|
[1705] | 221 | INTEGER :: ji, jj, jk ! dummy loop arguments |
---|
[2528] | 222 | !!bfr INTEGER :: ikbu, ikbv, ikbum1, ikbvm1 ! temporary scalar |
---|
| 223 | !!bfr INTEGER :: ikbt, ikbumm1, ikbvmm1 ! temporary scalar |
---|
[1705] | 224 | REAL(wp) :: zrhoa = 1.22 ! Air density kg/m3 |
---|
| 225 | REAL(wp) :: zcdrag = 1.5e-3 ! drag coefficient |
---|
| 226 | REAL(wp) :: zbbrau, zesh2 ! temporary scalars |
---|
| 227 | REAL(wp) :: zfact1, zfact2, zfact3 ! - - |
---|
| 228 | REAL(wp) :: ztx2 , zty2 , zcof ! - - |
---|
| 229 | REAL(wp) :: ztau , zdif ! - - |
---|
| 230 | REAL(wp) :: zus , zwlc , zind ! - - |
---|
| 231 | REAL(wp) :: zzd_up, zzd_lw ! - - |
---|
[2528] | 232 | !!bfr REAL(wp) :: zebot ! - - |
---|
[5836] | 233 | INTEGER , POINTER, DIMENSION(:,: ) :: imlc |
---|
| 234 | REAL(wp), POINTER, DIMENSION(:,: ) :: zhlc |
---|
| 235 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zpelc, zdiag, zd_up, zd_lw, z3du, z3dv |
---|
[5656] | 236 | REAL(wp) :: zri ! local Richardson number |
---|
[1239] | 237 | !!-------------------------------------------------------------------- |
---|
[1492] | 238 | ! |
---|
[3294] | 239 | IF( nn_timing == 1 ) CALL timing_start('tke_tke') |
---|
| 240 | ! |
---|
[5836] | 241 | CALL wrk_alloc( jpi,jpj, imlc ) ! integer |
---|
| 242 | CALL wrk_alloc( jpi,jpj, zhlc ) |
---|
| 243 | CALL wrk_alloc( jpi,jpj,jpk, zpelc, zdiag, zd_up, zd_lw, z3du, z3dv ) |
---|
[3294] | 244 | ! |
---|
[1695] | 245 | zbbrau = rn_ebb / rau0 ! Local constant initialisation |
---|
[2528] | 246 | zfact1 = -.5_wp * rdt |
---|
| 247 | zfact2 = 1.5_wp * rdt * rn_ediss |
---|
| 248 | zfact3 = 0.5_wp * rn_ediss |
---|
[1492] | 249 | ! |
---|
[5120] | 250 | ! |
---|
[1492] | 251 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 252 | ! ! Surface boundary condition on tke |
---|
| 253 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[5120] | 254 | IF ( ln_isfcav ) THEN |
---|
| 255 | DO jj = 2, jpjm1 ! en(mikt(ji,jj)) = rn_emin |
---|
| 256 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[5836] | 257 | en(ji,jj,mikt(ji,jj)) = rn_emin * tmask(ji,jj,1) |
---|
[5120] | 258 | END DO |
---|
| 259 | END DO |
---|
| 260 | END IF |
---|
[1695] | 261 | DO jj = 2, jpjm1 ! en(1) = rn_ebb taum / rau0 (min value rn_emin0) |
---|
[1481] | 262 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[5120] | 263 | en(ji,jj,1) = MAX( rn_emin0, zbbrau * taum(ji,jj) ) * tmask(ji,jj,1) |
---|
[1481] | 264 | END DO |
---|
| 265 | END DO |
---|
[2528] | 266 | |
---|
| 267 | !!bfr - start commented area |
---|
[1492] | 268 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 269 | ! ! Bottom boundary condition on tke |
---|
| 270 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[1719] | 271 | ! |
---|
| 272 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 273 | ! Tests to date have found the bottom boundary condition on tke to have very little effect. |
---|
| 274 | ! The condition is coded here for completion but commented out until there is proof that the |
---|
| 275 | ! computational cost is justified |
---|
| 276 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 277 | ! en(bot) = (rn_ebb0/rau0)*0.5*sqrt(u_botfr^2+v_botfr^2) (min value rn_emin) |
---|
| 278 | !! DO jj = 2, jpjm1 |
---|
| 279 | !! DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[2528] | 280 | !! ztx2 = bfrua(ji-1,jj) * ub(ji-1,jj,mbku(ji-1,jj)) + & |
---|
| 281 | !! bfrua(ji ,jj) * ub(ji ,jj,mbku(ji ,jj) ) |
---|
| 282 | !! zty2 = bfrva(ji,jj ) * vb(ji,jj ,mbkv(ji,jj )) + & |
---|
| 283 | !! bfrva(ji,jj-1) * vb(ji,jj-1,mbkv(ji,jj-1) ) |
---|
[1719] | 284 | !! zebot = 0.001875_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 ) ! where 0.001875 = (rn_ebb0/rau0) * 0.5 = 3.75*0.5/1000. |
---|
[2528] | 285 | !! en (ji,jj,mbkt(ji,jj)+1) = MAX( zebot, rn_emin ) * tmask(ji,jj,1) |
---|
[1719] | 286 | !! END DO |
---|
| 287 | !! END DO |
---|
[2528] | 288 | !!bfr - end commented area |
---|
[1492] | 289 | ! |
---|
| 290 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[2528] | 291 | IF( ln_lc ) THEN ! Langmuir circulation source term added to tke (Axell JGR 2002) |
---|
[1492] | 292 | ! !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
[1239] | 293 | ! |
---|
[1492] | 294 | ! !* total energy produce by LC : cumulative sum over jk |
---|
[6140] | 295 | zpelc(:,:,1) = MAX( rn2b(:,:,1), 0._wp ) * gdepw_n(:,:,1) * e3w_n(:,:,1) |
---|
[1239] | 296 | DO jk = 2, jpk |
---|
[6140] | 297 | zpelc(:,:,jk) = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0._wp ) * gdepw_n(:,:,jk) * e3w_n(:,:,jk) |
---|
[1239] | 298 | END DO |
---|
[1492] | 299 | ! !* finite Langmuir Circulation depth |
---|
[1705] | 300 | zcof = 0.5 * 0.016 * 0.016 / ( zrhoa * zcdrag ) |
---|
[2528] | 301 | imlc(:,:) = mbkt(:,:) + 1 ! Initialization to the number of w ocean point (=2 over land) |
---|
[1239] | 302 | DO jk = jpkm1, 2, -1 |
---|
[1492] | 303 | DO jj = 1, jpj ! Last w-level at which zpelc>=0.5*us*us |
---|
| 304 | DO ji = 1, jpi ! with us=0.016*wind(starting from jpk-1) |
---|
[1705] | 305 | zus = zcof * taum(ji,jj) |
---|
[1239] | 306 | IF( zpelc(ji,jj,jk) > zus ) imlc(ji,jj) = jk |
---|
| 307 | END DO |
---|
| 308 | END DO |
---|
| 309 | END DO |
---|
[1492] | 310 | ! ! finite LC depth |
---|
| 311 | DO jj = 1, jpj |
---|
[1239] | 312 | DO ji = 1, jpi |
---|
[6140] | 313 | zhlc(ji,jj) = gdepw_n(ji,jj,imlc(ji,jj)) |
---|
[1239] | 314 | END DO |
---|
| 315 | END DO |
---|
[1705] | 316 | zcof = 0.016 / SQRT( zrhoa * zcdrag ) |
---|
[1492] | 317 | DO jk = 2, jpkm1 !* TKE Langmuir circulation source term added to en |
---|
[1239] | 318 | DO jj = 2, jpjm1 |
---|
| 319 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[1705] | 320 | zus = zcof * SQRT( taum(ji,jj) ) ! Stokes drift |
---|
[1492] | 321 | ! ! vertical velocity due to LC |
---|
[6140] | 322 | zind = 0.5 - SIGN( 0.5, gdepw_n(ji,jj,jk) - zhlc(ji,jj) ) |
---|
| 323 | zwlc = zind * rn_lc * zus * SIN( rpi * gdepw_n(ji,jj,jk) / zhlc(ji,jj) ) |
---|
[1492] | 324 | ! ! TKE Langmuir circulation source term |
---|
[6497] | 325 | en(ji,jj,jk) = en(ji,jj,jk) + rdt * (1._wp - fr_i(ji,jj) ) * ( zwlc * zwlc * zwlc ) & |
---|
| 326 | & / zhlc(ji,jj) * wmask(ji,jj,jk) * tmask(ji,jj,1) |
---|
[1239] | 327 | END DO |
---|
| 328 | END DO |
---|
| 329 | END DO |
---|
| 330 | ! |
---|
| 331 | ENDIF |
---|
[1492] | 332 | ! |
---|
| 333 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 334 | ! ! Now Turbulent kinetic energy (output in en) |
---|
| 335 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 336 | ! ! Resolution of a tridiagonal linear system by a "methode de chasse" |
---|
| 337 | ! ! computation from level 2 to jpkm1 (e(1) already computed and e(jpk)=0 ). |
---|
| 338 | ! ! zdiag : diagonal zd_up : upper diagonal zd_lw : lower diagonal |
---|
| 339 | ! |
---|
| 340 | DO jk = 2, jpkm1 !* Shear production at uw- and vw-points (energy conserving form) |
---|
[5656] | 341 | DO jj = 1, jpjm1 |
---|
| 342 | DO ji = 1, fs_jpim1 ! vector opt. |
---|
| 343 | z3du(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk ) + avm(ji+1,jj,jk) ) & |
---|
| 344 | & * ( un(ji,jj,jk-1) - un(ji ,jj,jk) ) & |
---|
[5803] | 345 | & * ( ub(ji,jj,jk-1) - ub(ji ,jj,jk) ) * wumask(ji,jj,jk) & |
---|
[6140] | 346 | & / ( e3uw_n(ji,jj,jk) * e3uw_b(ji,jj,jk) ) |
---|
[5656] | 347 | z3dv(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk ) + avm(ji,jj+1,jk) ) & |
---|
| 348 | & * ( vn(ji,jj,jk-1) - vn(ji,jj ,jk) ) & |
---|
[5803] | 349 | & * ( vb(ji,jj,jk-1) - vb(ji,jj ,jk) ) * wvmask(ji,jj,jk) & |
---|
[6140] | 350 | & / ( e3vw_n(ji,jj,jk) * e3vw_b(ji,jj,jk) ) |
---|
[1492] | 351 | END DO |
---|
| 352 | END DO |
---|
| 353 | END DO |
---|
| 354 | ! |
---|
[5656] | 355 | IF( nn_pdl == 1 ) THEN !* Prandtl number case: compute apdlr |
---|
| 356 | ! Note that zesh2 is also computed in the next loop. |
---|
| 357 | ! We decided to compute it twice to keep code readability and avoid an IF case in the DO loops |
---|
| 358 | DO jk = 2, jpkm1 |
---|
| 359 | DO jj = 2, jpjm1 |
---|
| 360 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 361 | ! ! shear prod. at w-point weightened by mask |
---|
| 362 | zesh2 = ( z3du(ji-1,jj,jk) + z3du(ji,jj,jk) ) / MAX( 1._wp , umask(ji-1,jj,jk) + umask(ji,jj,jk) ) & |
---|
| 363 | & + ( z3dv(ji,jj-1,jk) + z3dv(ji,jj,jk) ) / MAX( 1._wp , vmask(ji,jj-1,jk) + vmask(ji,jj,jk) ) |
---|
| 364 | ! ! local Richardson number |
---|
| 365 | zri = MAX( rn2b(ji,jj,jk), 0._wp ) * avm(ji,jj,jk) / ( zesh2 + rn_bshear ) |
---|
| 366 | apdlr(ji,jj,jk) = MAX( 0.1_wp, ri_cri / MAX( ri_cri , zri ) ) |
---|
| 367 | |
---|
| 368 | END DO |
---|
| 369 | END DO |
---|
| 370 | END DO |
---|
| 371 | ! |
---|
| 372 | ENDIF |
---|
[5836] | 373 | ! |
---|
[5120] | 374 | DO jk = 2, jpkm1 !* Matrix and right hand side in en |
---|
| 375 | DO jj = 2, jpjm1 |
---|
| 376 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[1492] | 377 | zcof = zfact1 * tmask(ji,jj,jk) |
---|
[6497] | 378 | # if defined key_zdftmx_new |
---|
| 379 | ! key_zdftmx_new: New internal wave-driven param: set a minimum value for Kz on TKE (ensure numerical stability) |
---|
| 380 | zzd_up = zcof * MAX( avm(ji,jj,jk+1) + avm(ji,jj,jk), 2.e-5_wp ) & ! upper diagonal |
---|
| 381 | & / ( e3t_n(ji,jj,jk ) * e3w_n(ji,jj,jk ) ) |
---|
| 382 | zzd_lw = zcof * MAX( avm(ji,jj,jk) + avm(ji,jj,jk-1), 2.e-5_wp ) & ! lower diagonal |
---|
| 383 | & / ( e3t_n(ji,jj,jk-1) * e3w_n(ji,jj,jk ) ) |
---|
| 384 | # else |
---|
[1492] | 385 | zzd_up = zcof * ( avm (ji,jj,jk+1) + avm (ji,jj,jk ) ) & ! upper diagonal |
---|
[6140] | 386 | & / ( e3t_n(ji,jj,jk ) * e3w_n(ji,jj,jk ) ) |
---|
[1492] | 387 | zzd_lw = zcof * ( avm (ji,jj,jk ) + avm (ji,jj,jk-1) ) & ! lower diagonal |
---|
[6140] | 388 | & / ( e3t_n(ji,jj,jk-1) * e3w_n(ji,jj,jk ) ) |
---|
[6497] | 389 | # endif |
---|
[5656] | 390 | ! ! shear prod. at w-point weightened by mask |
---|
| 391 | zesh2 = ( z3du(ji-1,jj,jk) + z3du(ji,jj,jk) ) / MAX( 1._wp , umask(ji-1,jj,jk) + umask(ji,jj,jk) ) & |
---|
| 392 | & + ( z3dv(ji,jj-1,jk) + z3dv(ji,jj,jk) ) / MAX( 1._wp , vmask(ji,jj-1,jk) + vmask(ji,jj,jk) ) |
---|
| 393 | ! |
---|
[1492] | 394 | zd_up(ji,jj,jk) = zzd_up ! Matrix (zdiag, zd_up, zd_lw) |
---|
| 395 | zd_lw(ji,jj,jk) = zzd_lw |
---|
[2528] | 396 | zdiag(ji,jj,jk) = 1._wp - zzd_lw - zzd_up + zfact2 * dissl(ji,jj,jk) * tmask(ji,jj,jk) |
---|
[1239] | 397 | ! |
---|
[1492] | 398 | ! ! right hand side in en |
---|
[1481] | 399 | en(ji,jj,jk) = en(ji,jj,jk) + rdt * ( zesh2 - avt(ji,jj,jk) * rn2(ji,jj,jk) & |
---|
[4990] | 400 | & + zfact3 * dissl(ji,jj,jk) * en (ji,jj,jk) ) & |
---|
[5120] | 401 | & * wmask(ji,jj,jk) |
---|
[1239] | 402 | END DO |
---|
[5120] | 403 | END DO |
---|
| 404 | END DO |
---|
| 405 | ! !* Matrix inversion from level 2 (tke prescribed at level 1) |
---|
| 406 | DO jk = 3, jpkm1 ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1 |
---|
| 407 | DO jj = 2, jpjm1 |
---|
| 408 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[1492] | 409 | zdiag(ji,jj,jk) = zdiag(ji,jj,jk) - zd_lw(ji,jj,jk) * zd_up(ji,jj,jk-1) / zdiag(ji,jj,jk-1) |
---|
[1239] | 410 | END DO |
---|
[5120] | 411 | END DO |
---|
| 412 | END DO |
---|
[5836] | 413 | DO jj = 2, jpjm1 ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1 |
---|
[5120] | 414 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 415 | zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1) ! Surface boudary conditions on tke |
---|
| 416 | END DO |
---|
| 417 | END DO |
---|
| 418 | DO jk = 3, jpkm1 |
---|
| 419 | DO jj = 2, jpjm1 |
---|
| 420 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[1492] | 421 | zd_lw(ji,jj,jk) = en(ji,jj,jk) - zd_lw(ji,jj,jk) / zdiag(ji,jj,jk-1) *zd_lw(ji,jj,jk-1) |
---|
[1239] | 422 | END DO |
---|
[5120] | 423 | END DO |
---|
| 424 | END DO |
---|
[5836] | 425 | DO jj = 2, jpjm1 ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk |
---|
[5120] | 426 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[1492] | 427 | en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1) |
---|
[5120] | 428 | END DO |
---|
| 429 | END DO |
---|
| 430 | DO jk = jpk-2, 2, -1 |
---|
| 431 | DO jj = 2, jpjm1 |
---|
| 432 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[1492] | 433 | en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk) |
---|
[1239] | 434 | END DO |
---|
[5120] | 435 | END DO |
---|
| 436 | END DO |
---|
| 437 | DO jk = 2, jpkm1 ! set the minimum value of tke |
---|
| 438 | DO jj = 2, jpjm1 |
---|
| 439 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 440 | en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * wmask(ji,jj,jk) |
---|
[1239] | 441 | END DO |
---|
| 442 | END DO |
---|
| 443 | END DO |
---|
| 444 | |
---|
[1492] | 445 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 446 | ! ! TKE due to surface and internal wave breaking |
---|
| 447 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[6140] | 448 | !!gm BUG : in the exp remove the depth of ssh !!! |
---|
| 449 | |
---|
| 450 | |
---|
[2528] | 451 | IF( nn_etau == 1 ) THEN !* penetration below the mixed layer (rn_efr fraction) |
---|
[1492] | 452 | DO jk = 2, jpkm1 |
---|
[1239] | 453 | DO jj = 2, jpjm1 |
---|
| 454 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[6140] | 455 | en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -gdepw_n(ji,jj,jk) / htau(ji,jj) ) & |
---|
[5120] | 456 | & * ( 1._wp - fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1) |
---|
[1239] | 457 | END DO |
---|
| 458 | END DO |
---|
[1492] | 459 | END DO |
---|
[2528] | 460 | ELSEIF( nn_etau == 2 ) THEN !* act only at the base of the mixed layer (jk=nmln) (rn_efr fraction) |
---|
[1492] | 461 | DO jj = 2, jpjm1 |
---|
| 462 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 463 | jk = nmln(ji,jj) |
---|
[6140] | 464 | en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -gdepw_n(ji,jj,jk) / htau(ji,jj) ) & |
---|
[5120] | 465 | & * ( 1._wp - fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1) |
---|
[1239] | 466 | END DO |
---|
[1492] | 467 | END DO |
---|
[2528] | 468 | ELSEIF( nn_etau == 3 ) THEN !* penetration belox the mixed layer (HF variability) |
---|
[1705] | 469 | DO jk = 2, jpkm1 |
---|
| 470 | DO jj = 2, jpjm1 |
---|
| 471 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 472 | ztx2 = utau(ji-1,jj ) + utau(ji,jj) |
---|
| 473 | zty2 = vtau(ji ,jj-1) + vtau(ji,jj) |
---|
[4990] | 474 | ztau = 0.5_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1) ! module of the mean stress |
---|
[2528] | 475 | zdif = taum(ji,jj) - ztau ! mean of modulus - modulus of the mean |
---|
| 476 | zdif = rhftau_scl * MAX( 0._wp, zdif + rhftau_add ) ! apply some modifications... |
---|
[6140] | 477 | en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -gdepw_n(ji,jj,jk) / htau(ji,jj) ) & |
---|
[5120] | 478 | & * ( 1._wp - fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1) |
---|
[1705] | 479 | END DO |
---|
| 480 | END DO |
---|
| 481 | END DO |
---|
[1239] | 482 | ENDIF |
---|
[1492] | 483 | CALL lbc_lnk( en, 'W', 1. ) ! Lateral boundary conditions (sign unchanged) |
---|
| 484 | ! |
---|
[5836] | 485 | CALL wrk_dealloc( jpi,jpj, imlc ) ! integer |
---|
| 486 | CALL wrk_dealloc( jpi,jpj, zhlc ) |
---|
| 487 | CALL wrk_dealloc( jpi,jpj,jpk, zpelc, zdiag, zd_up, zd_lw, z3du, z3dv ) |
---|
[2715] | 488 | ! |
---|
[3294] | 489 | IF( nn_timing == 1 ) CALL timing_stop('tke_tke') |
---|
| 490 | ! |
---|
[1239] | 491 | END SUBROUTINE tke_tke |
---|
| 492 | |
---|
[1492] | 493 | |
---|
| 494 | SUBROUTINE tke_avn |
---|
[1239] | 495 | !!---------------------------------------------------------------------- |
---|
[1492] | 496 | !! *** ROUTINE tke_avn *** |
---|
[1239] | 497 | !! |
---|
[1492] | 498 | !! ** Purpose : Compute the vertical eddy viscosity and diffusivity |
---|
| 499 | !! |
---|
| 500 | !! ** Method : At this stage, en, the now TKE, is known (computed in |
---|
| 501 | !! the tke_tke routine). First, the now mixing lenth is |
---|
| 502 | !! computed from en and the strafification (N^2), then the mixings |
---|
| 503 | !! coefficients are computed. |
---|
| 504 | !! - Mixing length : a first evaluation of the mixing lengh |
---|
| 505 | !! scales is: |
---|
| 506 | !! mxl = sqrt(2*en) / N |
---|
| 507 | !! where N is the brunt-vaisala frequency, with a minimum value set |
---|
[2528] | 508 | !! to rmxl_min (rn_mxl0) in the interior (surface) ocean. |
---|
[1492] | 509 | !! The mixing and dissipative length scale are bound as follow : |
---|
| 510 | !! nn_mxl=0 : mxl bounded by the distance to surface and bottom. |
---|
| 511 | !! zmxld = zmxlm = mxl |
---|
| 512 | !! nn_mxl=1 : mxl bounded by the e3w and zmxld = zmxlm = mxl |
---|
| 513 | !! nn_mxl=2 : mxl bounded such that the vertical derivative of mxl is |
---|
| 514 | !! less than 1 (|d/dz(mxl)|<1) and zmxld = zmxlm = mxl |
---|
| 515 | !! nn_mxl=3 : mxl is bounded from the surface to the bottom usings |
---|
| 516 | !! |d/dz(xml)|<1 to obtain lup, and from the bottom to |
---|
| 517 | !! the surface to obtain ldown. the resulting length |
---|
| 518 | !! scales are: |
---|
| 519 | !! zmxld = sqrt( lup * ldown ) |
---|
| 520 | !! zmxlm = min ( lup , ldown ) |
---|
| 521 | !! - Vertical eddy viscosity and diffusivity: |
---|
| 522 | !! avm = max( avtb, rn_ediff * zmxlm * en^1/2 ) |
---|
| 523 | !! avt = max( avmb, pdlr * avm ) |
---|
| 524 | !! with pdlr=1 if nn_pdl=0, pdlr=1/pdl=F(Ri) otherwise. |
---|
| 525 | !! |
---|
| 526 | !! ** Action : - avt : now vertical eddy diffusivity (w-point) |
---|
| 527 | !! - avmu, avmv : now vertical eddy viscosity at uw- and vw-points |
---|
[1239] | 528 | !!---------------------------------------------------------------------- |
---|
[2715] | 529 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 530 | REAL(wp) :: zrn2, zraug, zcoef, zav ! local scalars |
---|
[5836] | 531 | REAL(wp) :: zdku, zri, zsqen ! - - |
---|
[2715] | 532 | REAL(wp) :: zdkv, zemxl, zemlm, zemlp ! - - |
---|
[3294] | 533 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zmpdl, zmxlm, zmxld |
---|
[1239] | 534 | !!-------------------------------------------------------------------- |
---|
[3294] | 535 | ! |
---|
| 536 | IF( nn_timing == 1 ) CALL timing_start('tke_avn') |
---|
[1239] | 537 | |
---|
[3294] | 538 | CALL wrk_alloc( jpi,jpj,jpk, zmpdl, zmxlm, zmxld ) |
---|
| 539 | |
---|
[1492] | 540 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 541 | ! ! Mixing length |
---|
| 542 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 543 | ! |
---|
| 544 | ! !* Buoyancy length scale: l=sqrt(2*e/n**2) |
---|
| 545 | ! |
---|
[5120] | 546 | ! initialisation of interior minimum value (avoid a 2d loop with mikt) |
---|
| 547 | zmxlm(:,:,:) = rmxl_min |
---|
| 548 | zmxld(:,:,:) = rmxl_min |
---|
| 549 | ! |
---|
[2528] | 550 | IF( ln_mxl0 ) THEN ! surface mixing length = F(stress) : l=vkarmn*2.e5*taum/(rau0*g) |
---|
[4990] | 551 | DO jj = 2, jpjm1 |
---|
| 552 | DO ji = fs_2, fs_jpim1 |
---|
[5120] | 553 | zraug = vkarmn * 2.e5_wp / ( rau0 * grav ) |
---|
| 554 | zmxlm(ji,jj,1) = MAX( rn_mxl0, zraug * taum(ji,jj) * tmask(ji,jj,1) ) |
---|
[4990] | 555 | END DO |
---|
| 556 | END DO |
---|
| 557 | ELSE |
---|
[5120] | 558 | zmxlm(:,:,1) = rn_mxl0 |
---|
[1239] | 559 | ENDIF |
---|
| 560 | ! |
---|
[5120] | 561 | DO jk = 2, jpkm1 ! interior value : l=sqrt(2*e/n^2) |
---|
| 562 | DO jj = 2, jpjm1 |
---|
| 563 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[1239] | 564 | zrn2 = MAX( rn2(ji,jj,jk), rsmall ) |
---|
[5836] | 565 | zmxlm(ji,jj,jk) = MAX( rmxl_min, SQRT( 2._wp * en(ji,jj,jk) / zrn2 ) ) |
---|
[1239] | 566 | END DO |
---|
| 567 | END DO |
---|
| 568 | END DO |
---|
[1492] | 569 | ! |
---|
| 570 | ! !* Physical limits for the mixing length |
---|
| 571 | ! |
---|
[5836] | 572 | zmxld(:,:, 1 ) = zmxlm(:,:,1) ! surface set to the minimum value |
---|
[2528] | 573 | zmxld(:,:,jpk) = rmxl_min ! last level set to the minimum value |
---|
[1492] | 574 | ! |
---|
[1239] | 575 | SELECT CASE ( nn_mxl ) |
---|
| 576 | ! |
---|
[5836] | 577 | !!gm Not sure of that coding for ISF.... |
---|
[6140] | 578 | ! where wmask = 0 set zmxlm == e3w_n |
---|
[1239] | 579 | CASE ( 0 ) ! bounded by the distance to surface and bottom |
---|
[5120] | 580 | DO jk = 2, jpkm1 |
---|
| 581 | DO jj = 2, jpjm1 |
---|
| 582 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[6140] | 583 | zemxl = MIN( gdepw_n(ji,jj,jk) - gdepw_n(ji,jj,mikt(ji,jj)), zmxlm(ji,jj,jk), & |
---|
| 584 | & gdepw_n(ji,jj,mbkt(ji,jj)+1) - gdepw_n(ji,jj,jk) ) |
---|
[5120] | 585 | ! wmask prevent zmxlm = 0 if jk = mikt(ji,jj) |
---|
[6140] | 586 | zmxlm(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN(zmxlm(ji,jj,jk),e3w_n(ji,jj,jk)) * (1 - wmask(ji,jj,jk)) |
---|
| 587 | zmxld(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN(zmxlm(ji,jj,jk),e3w_n(ji,jj,jk)) * (1 - wmask(ji,jj,jk)) |
---|
[1239] | 588 | END DO |
---|
| 589 | END DO |
---|
| 590 | END DO |
---|
| 591 | ! |
---|
| 592 | CASE ( 1 ) ! bounded by the vertical scale factor |
---|
[5120] | 593 | DO jk = 2, jpkm1 |
---|
| 594 | DO jj = 2, jpjm1 |
---|
| 595 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[6140] | 596 | zemxl = MIN( e3w_n(ji,jj,jk), zmxlm(ji,jj,jk) ) |
---|
[1239] | 597 | zmxlm(ji,jj,jk) = zemxl |
---|
| 598 | zmxld(ji,jj,jk) = zemxl |
---|
| 599 | END DO |
---|
| 600 | END DO |
---|
| 601 | END DO |
---|
| 602 | ! |
---|
| 603 | CASE ( 2 ) ! |dk[xml]| bounded by e3t : |
---|
[5120] | 604 | DO jk = 2, jpkm1 ! from the surface to the bottom : |
---|
| 605 | DO jj = 2, jpjm1 |
---|
| 606 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[6140] | 607 | zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + e3t_n(ji,jj,jk-1), zmxlm(ji,jj,jk) ) |
---|
[1239] | 608 | END DO |
---|
[5120] | 609 | END DO |
---|
| 610 | END DO |
---|
| 611 | DO jk = jpkm1, 2, -1 ! from the bottom to the surface : |
---|
| 612 | DO jj = 2, jpjm1 |
---|
| 613 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[6140] | 614 | zemxl = MIN( zmxlm(ji,jj,jk+1) + e3t_n(ji,jj,jk+1), zmxlm(ji,jj,jk) ) |
---|
[1239] | 615 | zmxlm(ji,jj,jk) = zemxl |
---|
| 616 | zmxld(ji,jj,jk) = zemxl |
---|
| 617 | END DO |
---|
| 618 | END DO |
---|
| 619 | END DO |
---|
| 620 | ! |
---|
| 621 | CASE ( 3 ) ! lup and ldown, |dk[xml]| bounded by e3t : |
---|
[5120] | 622 | DO jk = 2, jpkm1 ! from the surface to the bottom : lup |
---|
| 623 | DO jj = 2, jpjm1 |
---|
| 624 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[6140] | 625 | zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + e3t_n(ji,jj,jk-1), zmxlm(ji,jj,jk) ) |
---|
[1239] | 626 | END DO |
---|
[5120] | 627 | END DO |
---|
| 628 | END DO |
---|
| 629 | DO jk = jpkm1, 2, -1 ! from the bottom to the surface : ldown |
---|
| 630 | DO jj = 2, jpjm1 |
---|
| 631 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[6140] | 632 | zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + e3t_n(ji,jj,jk+1), zmxlm(ji,jj,jk) ) |
---|
[1239] | 633 | END DO |
---|
| 634 | END DO |
---|
| 635 | END DO |
---|
| 636 | DO jk = 2, jpkm1 |
---|
| 637 | DO jj = 2, jpjm1 |
---|
| 638 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 639 | zemlm = MIN ( zmxld(ji,jj,jk), zmxlm(ji,jj,jk) ) |
---|
| 640 | zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) ) |
---|
| 641 | zmxlm(ji,jj,jk) = zemlm |
---|
| 642 | zmxld(ji,jj,jk) = zemlp |
---|
| 643 | END DO |
---|
| 644 | END DO |
---|
| 645 | END DO |
---|
| 646 | ! |
---|
| 647 | END SELECT |
---|
[1492] | 648 | ! |
---|
[1239] | 649 | # if defined key_c1d |
---|
[1492] | 650 | e_dis(:,:,:) = zmxld(:,:,:) ! c1d configuration : save mixing and dissipation turbulent length scales |
---|
[1239] | 651 | e_mix(:,:,:) = zmxlm(:,:,:) |
---|
| 652 | # endif |
---|
| 653 | |
---|
[1492] | 654 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 655 | ! ! Vertical eddy viscosity and diffusivity (avmu, avmv, avt) |
---|
| 656 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 657 | DO jk = 1, jpkm1 !* vertical eddy viscosity & diffivity at w-points |
---|
[1239] | 658 | DO jj = 2, jpjm1 |
---|
| 659 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 660 | zsqen = SQRT( en(ji,jj,jk) ) |
---|
| 661 | zav = rn_ediff * zmxlm(ji,jj,jk) * zsqen |
---|
[5120] | 662 | avm (ji,jj,jk) = MAX( zav, avmb(jk) ) * wmask(ji,jj,jk) |
---|
| 663 | avt (ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk) |
---|
[1239] | 664 | dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk) |
---|
| 665 | END DO |
---|
| 666 | END DO |
---|
| 667 | END DO |
---|
[1492] | 668 | CALL lbc_lnk( avm, 'W', 1. ) ! Lateral boundary conditions (sign unchanged) |
---|
| 669 | ! |
---|
[5120] | 670 | DO jk = 2, jpkm1 !* vertical eddy viscosity at wu- and wv-points |
---|
| 671 | DO jj = 2, jpjm1 |
---|
| 672 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 673 | avmu(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji+1,jj ,jk) ) * wumask(ji,jj,jk) |
---|
| 674 | avmv(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji ,jj+1,jk) ) * wvmask(ji,jj,jk) |
---|
[4990] | 675 | END DO |
---|
[1239] | 676 | END DO |
---|
| 677 | END DO |
---|
| 678 | CALL lbc_lnk( avmu, 'U', 1. ) ; CALL lbc_lnk( avmv, 'V', 1. ) ! Lateral boundary conditions |
---|
[1492] | 679 | ! |
---|
| 680 | IF( nn_pdl == 1 ) THEN !* Prandtl number case: update avt |
---|
[5120] | 681 | DO jk = 2, jpkm1 |
---|
| 682 | DO jj = 2, jpjm1 |
---|
| 683 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[5656] | 684 | avt(ji,jj,jk) = MAX( apdlr(ji,jj,jk) * avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * tmask(ji,jj,jk) |
---|
[1492] | 685 | # if defined key_c1d |
---|
[5656] | 686 | e_pdl(ji,jj,jk) = apdlr(ji,jj,jk) * wmask(ji,jj,jk) ! c1d configuration : save masked Prandlt number |
---|
[5836] | 687 | !!gm bug NO zri here.... |
---|
| 688 | !!gm remove the specific diag for c1d ! |
---|
[5656] | 689 | e_ric(ji,jj,jk) = zri * wmask(ji,jj,jk) ! c1d config. : save Ri |
---|
[1239] | 690 | # endif |
---|
| 691 | END DO |
---|
| 692 | END DO |
---|
| 693 | END DO |
---|
| 694 | ENDIF |
---|
| 695 | CALL lbc_lnk( avt, 'W', 1. ) ! Lateral boundary conditions on avt (sign unchanged) |
---|
| 696 | |
---|
| 697 | IF(ln_ctl) THEN |
---|
| 698 | CALL prt_ctl( tab3d_1=en , clinfo1=' tke - e: ', tab3d_2=avt, clinfo2=' t: ', ovlap=1, kdim=jpk) |
---|
| 699 | CALL prt_ctl( tab3d_1=avmu, clinfo1=' tke - u: ', mask1=umask, & |
---|
| 700 | & tab3d_2=avmv, clinfo2= ' v: ', mask2=vmask, ovlap=1, kdim=jpk ) |
---|
| 701 | ENDIF |
---|
| 702 | ! |
---|
[3294] | 703 | CALL wrk_dealloc( jpi,jpj,jpk, zmpdl, zmxlm, zmxld ) |
---|
| 704 | ! |
---|
| 705 | IF( nn_timing == 1 ) CALL timing_stop('tke_avn') |
---|
| 706 | ! |
---|
[1492] | 707 | END SUBROUTINE tke_avn |
---|
[1239] | 708 | |
---|
[1492] | 709 | |
---|
[2528] | 710 | SUBROUTINE zdf_tke_init |
---|
[1239] | 711 | !!---------------------------------------------------------------------- |
---|
[2528] | 712 | !! *** ROUTINE zdf_tke_init *** |
---|
[1239] | 713 | !! |
---|
| 714 | !! ** Purpose : Initialization of the vertical eddy diffivity and |
---|
[1492] | 715 | !! viscosity when using a tke turbulent closure scheme |
---|
[1239] | 716 | !! |
---|
[1601] | 717 | !! ** Method : Read the namzdf_tke namelist and check the parameters |
---|
[1492] | 718 | !! called at the first timestep (nit000) |
---|
[1239] | 719 | !! |
---|
[1601] | 720 | !! ** input : Namlist namzdf_tke |
---|
[1239] | 721 | !! |
---|
| 722 | !! ** Action : Increase by 1 the nstop flag is setting problem encounter |
---|
| 723 | !!---------------------------------------------------------------------- |
---|
| 724 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
[4147] | 725 | INTEGER :: ios |
---|
[1239] | 726 | !! |
---|
[2528] | 727 | NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb , rn_emin , & |
---|
| 728 | & rn_emin0, rn_bshear, nn_mxl , ln_mxl0 , & |
---|
| 729 | & rn_mxl0 , nn_pdl , ln_lc , rn_lc , & |
---|
| 730 | & nn_etau , nn_htau , rn_efr |
---|
[1239] | 731 | !!---------------------------------------------------------------------- |
---|
[2715] | 732 | ! |
---|
[4147] | 733 | REWIND( numnam_ref ) ! Namelist namzdf_tke in reference namelist : Turbulent Kinetic Energy |
---|
| 734 | READ ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901) |
---|
| 735 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in reference namelist', lwp ) |
---|
| 736 | |
---|
| 737 | REWIND( numnam_cfg ) ! Namelist namzdf_tke in configuration namelist : Turbulent Kinetic Energy |
---|
| 738 | READ ( numnam_cfg, namzdf_tke, IOSTAT = ios, ERR = 902 ) |
---|
| 739 | 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in configuration namelist', lwp ) |
---|
[4624] | 740 | IF(lwm) WRITE ( numond, namzdf_tke ) |
---|
[2715] | 741 | ! |
---|
[2528] | 742 | ri_cri = 2._wp / ( 2._wp + rn_ediss / rn_ediff ) ! resulting critical Richardson number |
---|
[6497] | 743 | # if defined key_zdftmx_new |
---|
| 744 | ! key_zdftmx_new: New internal wave-driven param: specified value of rn_emin & rmxl_min are used |
---|
| 745 | rn_emin = 1.e-10_wp |
---|
| 746 | rmxl_min = 1.e-03_wp |
---|
| 747 | IF(lwp) THEN ! Control print |
---|
| 748 | WRITE(numout,*) |
---|
| 749 | WRITE(numout,*) 'zdf_tke_init : New tidal mixing case: force rn_emin = 1.e-10 and rmxl_min = 1.e-3 ' |
---|
| 750 | WRITE(numout,*) '~~~~~~~~~~~~' |
---|
| 751 | ENDIF |
---|
| 752 | # else |
---|
[2528] | 753 | rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) ) ! resulting minimum length to recover molecular viscosity |
---|
[6497] | 754 | # endif |
---|
[2715] | 755 | ! |
---|
[1492] | 756 | IF(lwp) THEN !* Control print |
---|
[1239] | 757 | WRITE(numout,*) |
---|
[2528] | 758 | WRITE(numout,*) 'zdf_tke_init : tke turbulent closure scheme - initialisation' |
---|
| 759 | WRITE(numout,*) '~~~~~~~~~~~~' |
---|
[1601] | 760 | WRITE(numout,*) ' Namelist namzdf_tke : set tke mixing parameters' |
---|
[1705] | 761 | WRITE(numout,*) ' coef. to compute avt rn_ediff = ', rn_ediff |
---|
| 762 | WRITE(numout,*) ' Kolmogoroff dissipation coef. rn_ediss = ', rn_ediss |
---|
| 763 | WRITE(numout,*) ' tke surface input coef. rn_ebb = ', rn_ebb |
---|
| 764 | WRITE(numout,*) ' minimum value of tke rn_emin = ', rn_emin |
---|
| 765 | WRITE(numout,*) ' surface minimum value of tke rn_emin0 = ', rn_emin0 |
---|
| 766 | WRITE(numout,*) ' background shear (>0) rn_bshear = ', rn_bshear |
---|
| 767 | WRITE(numout,*) ' mixing length type nn_mxl = ', nn_mxl |
---|
| 768 | WRITE(numout,*) ' prandl number flag nn_pdl = ', nn_pdl |
---|
| 769 | WRITE(numout,*) ' surface mixing length = F(stress) or not ln_mxl0 = ', ln_mxl0 |
---|
[2528] | 770 | WRITE(numout,*) ' surface mixing length minimum value rn_mxl0 = ', rn_mxl0 |
---|
| 771 | WRITE(numout,*) ' flag to take into acc. Langmuir circ. ln_lc = ', ln_lc |
---|
| 772 | WRITE(numout,*) ' coef to compute verticla velocity of LC rn_lc = ', rn_lc |
---|
[1705] | 773 | WRITE(numout,*) ' test param. to add tke induced by wind nn_etau = ', nn_etau |
---|
| 774 | WRITE(numout,*) ' flag for computation of exp. tke profile nn_htau = ', nn_htau |
---|
| 775 | WRITE(numout,*) ' fraction of en which pene. the thermocline rn_efr = ', rn_efr |
---|
[1239] | 776 | WRITE(numout,*) |
---|
[1601] | 777 | WRITE(numout,*) ' critical Richardson nb with your parameters ri_cri = ', ri_cri |
---|
[1239] | 778 | ENDIF |
---|
[2715] | 779 | ! |
---|
| 780 | ! ! allocate tke arrays |
---|
| 781 | IF( zdf_tke_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'zdf_tke_init : unable to allocate arrays' ) |
---|
| 782 | ! |
---|
[1492] | 783 | ! !* Check of some namelist values |
---|
[4990] | 784 | IF( nn_mxl < 0 .OR. nn_mxl > 3 ) CALL ctl_stop( 'bad flag: nn_mxl is 0, 1 or 2 ' ) |
---|
| 785 | IF( nn_pdl < 0 .OR. nn_pdl > 1 ) CALL ctl_stop( 'bad flag: nn_pdl is 0 or 1 ' ) |
---|
| 786 | IF( nn_htau < 0 .OR. nn_htau > 1 ) CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2 ' ) |
---|
[5407] | 787 | IF( nn_etau == 3 .AND. .NOT. ln_cpl ) CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' ) |
---|
[1239] | 788 | |
---|
[2528] | 789 | IF( ln_mxl0 ) THEN |
---|
| 790 | IF(lwp) WRITE(numout,*) ' use a surface mixing length = F(stress) : set rn_mxl0 = rmxl_min' |
---|
| 791 | rn_mxl0 = rmxl_min |
---|
| 792 | ENDIF |
---|
| 793 | |
---|
[1492] | 794 | IF( nn_etau == 2 ) CALL zdf_mxl( nit000 ) ! Initialization of nmln |
---|
[1239] | 795 | |
---|
[1492] | 796 | ! !* depth of penetration of surface tke |
---|
| 797 | IF( nn_etau /= 0 ) THEN |
---|
[1601] | 798 | SELECT CASE( nn_htau ) ! Choice of the depth of penetration |
---|
[2528] | 799 | CASE( 0 ) ! constant depth penetration (here 10 meters) |
---|
| 800 | htau(:,:) = 10._wp |
---|
| 801 | CASE( 1 ) ! F(latitude) : 0.5m to 30m poleward of 40 degrees |
---|
| 802 | htau(:,:) = MAX( 0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) ) ) |
---|
[1492] | 803 | END SELECT |
---|
| 804 | ENDIF |
---|
| 805 | ! !* set vertical eddy coef. to the background value |
---|
[1239] | 806 | DO jk = 1, jpk |
---|
[5120] | 807 | avt (:,:,jk) = avtb(jk) * wmask (:,:,jk) |
---|
| 808 | avm (:,:,jk) = avmb(jk) * wmask (:,:,jk) |
---|
| 809 | avmu(:,:,jk) = avmb(jk) * wumask(:,:,jk) |
---|
| 810 | avmv(:,:,jk) = avmb(jk) * wvmask(:,:,jk) |
---|
[1239] | 811 | END DO |
---|
[2528] | 812 | dissl(:,:,:) = 1.e-12_wp |
---|
[2715] | 813 | ! |
---|
| 814 | CALL tke_rst( nit000, 'READ' ) !* read or initialize all required files |
---|
[1239] | 815 | ! |
---|
[2528] | 816 | END SUBROUTINE zdf_tke_init |
---|
[1239] | 817 | |
---|
| 818 | |
---|
[1531] | 819 | SUBROUTINE tke_rst( kt, cdrw ) |
---|
[1239] | 820 | !!--------------------------------------------------------------------- |
---|
[1531] | 821 | !! *** ROUTINE tke_rst *** |
---|
[1239] | 822 | !! |
---|
| 823 | !! ** Purpose : Read or write TKE file (en) in restart file |
---|
| 824 | !! |
---|
| 825 | !! ** Method : use of IOM library |
---|
| 826 | !! if the restart does not contain TKE, en is either |
---|
[1537] | 827 | !! set to rn_emin or recomputed |
---|
[1239] | 828 | !!---------------------------------------------------------------------- |
---|
[2715] | 829 | INTEGER , INTENT(in) :: kt ! ocean time-step |
---|
| 830 | CHARACTER(len=*), INTENT(in) :: cdrw ! "READ"/"WRITE" flag |
---|
[1239] | 831 | ! |
---|
[1481] | 832 | INTEGER :: jit, jk ! dummy loop indices |
---|
[2715] | 833 | INTEGER :: id1, id2, id3, id4, id5, id6 ! local integers |
---|
[1239] | 834 | !!---------------------------------------------------------------------- |
---|
| 835 | ! |
---|
[1481] | 836 | IF( TRIM(cdrw) == 'READ' ) THEN ! Read/initialise |
---|
| 837 | ! ! --------------- |
---|
| 838 | IF( ln_rstart ) THEN !* Read the restart file |
---|
| 839 | id1 = iom_varid( numror, 'en' , ldstop = .FALSE. ) |
---|
| 840 | id2 = iom_varid( numror, 'avt' , ldstop = .FALSE. ) |
---|
| 841 | id3 = iom_varid( numror, 'avm' , ldstop = .FALSE. ) |
---|
| 842 | id4 = iom_varid( numror, 'avmu' , ldstop = .FALSE. ) |
---|
| 843 | id5 = iom_varid( numror, 'avmv' , ldstop = .FALSE. ) |
---|
| 844 | id6 = iom_varid( numror, 'dissl', ldstop = .FALSE. ) |
---|
| 845 | ! |
---|
| 846 | IF( id1 > 0 ) THEN ! 'en' exists |
---|
[1239] | 847 | CALL iom_get( numror, jpdom_autoglo, 'en', en ) |
---|
[1481] | 848 | IF( MIN( id2, id3, id4, id5, id6 ) > 0 ) THEN ! all required arrays exist |
---|
| 849 | CALL iom_get( numror, jpdom_autoglo, 'avt' , avt ) |
---|
| 850 | CALL iom_get( numror, jpdom_autoglo, 'avm' , avm ) |
---|
| 851 | CALL iom_get( numror, jpdom_autoglo, 'avmu' , avmu ) |
---|
| 852 | CALL iom_get( numror, jpdom_autoglo, 'avmv' , avmv ) |
---|
| 853 | CALL iom_get( numror, jpdom_autoglo, 'dissl', dissl ) |
---|
[1492] | 854 | ELSE ! one at least array is missing |
---|
| 855 | CALL tke_avn ! compute avt, avm, avmu, avmv and dissl (approximation) |
---|
[1481] | 856 | ENDIF |
---|
| 857 | ELSE ! No TKE array found: initialisation |
---|
| 858 | IF(lwp) WRITE(numout,*) ' ===>>>> : previous run without tke scheme, en computed by iterative loop' |
---|
[1239] | 859 | en (:,:,:) = rn_emin * tmask(:,:,:) |
---|
[1492] | 860 | CALL tke_avn ! recompute avt, avm, avmu, avmv and dissl (approximation) |
---|
[5112] | 861 | ! |
---|
| 862 | avt_k (:,:,:) = avt (:,:,:) |
---|
| 863 | avm_k (:,:,:) = avm (:,:,:) |
---|
| 864 | avmu_k(:,:,:) = avmu(:,:,:) |
---|
| 865 | avmv_k(:,:,:) = avmv(:,:,:) |
---|
| 866 | ! |
---|
[1531] | 867 | DO jit = nit000 + 1, nit000 + 10 ; CALL zdf_tke( jit ) ; END DO |
---|
[1239] | 868 | ENDIF |
---|
[1481] | 869 | ELSE !* Start from rest |
---|
| 870 | en(:,:,:) = rn_emin * tmask(:,:,:) |
---|
| 871 | DO jk = 1, jpk ! set the Kz to the background value |
---|
[5120] | 872 | avt (:,:,jk) = avtb(jk) * wmask (:,:,jk) |
---|
| 873 | avm (:,:,jk) = avmb(jk) * wmask (:,:,jk) |
---|
| 874 | avmu(:,:,jk) = avmb(jk) * wumask(:,:,jk) |
---|
| 875 | avmv(:,:,jk) = avmb(jk) * wvmask(:,:,jk) |
---|
[1481] | 876 | END DO |
---|
[1239] | 877 | ENDIF |
---|
[1481] | 878 | ! |
---|
| 879 | ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN ! Create restart file |
---|
| 880 | ! ! ------------------- |
---|
[1531] | 881 | IF(lwp) WRITE(numout,*) '---- tke-rst ----' |
---|
[3632] | 882 | CALL iom_rstput( kt, nitrst, numrow, 'en' , en ) |
---|
| 883 | CALL iom_rstput( kt, nitrst, numrow, 'avt' , avt_k ) |
---|
| 884 | CALL iom_rstput( kt, nitrst, numrow, 'avm' , avm_k ) |
---|
| 885 | CALL iom_rstput( kt, nitrst, numrow, 'avmu' , avmu_k ) |
---|
| 886 | CALL iom_rstput( kt, nitrst, numrow, 'avmv' , avmv_k ) |
---|
| 887 | CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl ) |
---|
[1481] | 888 | ! |
---|
[1239] | 889 | ENDIF |
---|
| 890 | ! |
---|
[1531] | 891 | END SUBROUTINE tke_rst |
---|
[1239] | 892 | |
---|
| 893 | #else |
---|
| 894 | !!---------------------------------------------------------------------- |
---|
| 895 | !! Dummy module : NO TKE scheme |
---|
| 896 | !!---------------------------------------------------------------------- |
---|
[1531] | 897 | LOGICAL, PUBLIC, PARAMETER :: lk_zdftke = .FALSE. !: TKE flag |
---|
[1239] | 898 | CONTAINS |
---|
[2528] | 899 | SUBROUTINE zdf_tke_init ! Dummy routine |
---|
| 900 | END SUBROUTINE zdf_tke_init |
---|
| 901 | SUBROUTINE zdf_tke( kt ) ! Dummy routine |
---|
[1531] | 902 | WRITE(*,*) 'zdf_tke: You should not have seen this print! error?', kt |
---|
| 903 | END SUBROUTINE zdf_tke |
---|
| 904 | SUBROUTINE tke_rst( kt, cdrw ) |
---|
[1492] | 905 | CHARACTER(len=*) :: cdrw |
---|
[1531] | 906 | WRITE(*,*) 'tke_rst: You should not have seen this print! error?', kt, cdwr |
---|
| 907 | END SUBROUTINE tke_rst |
---|
[1239] | 908 | #endif |
---|
| 909 | |
---|
| 910 | !!====================================================================== |
---|
[1531] | 911 | END MODULE zdftke |
---|