[2048] | 1 | MODULE zdfgls |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE zdfgls *** |
---|
| 4 | !! Ocean physics: vertical mixing coefficient computed from the gls |
---|
| 5 | !! turbulent closure parameterization |
---|
| 6 | !!====================================================================== |
---|
[2397] | 7 | !! History : 3.0 ! 2009-09 (G. Reffray) Original code |
---|
| 8 | !! 3.3 ! 2010-10 (C. Bricaud) Add in the reference |
---|
[2048] | 9 | !!---------------------------------------------------------------------- |
---|
| 10 | #if defined key_zdfgls || defined key_esopa |
---|
| 11 | !!---------------------------------------------------------------------- |
---|
| 12 | !! 'key_zdfgls' Generic Length Scale vertical physics |
---|
| 13 | !!---------------------------------------------------------------------- |
---|
[3625] | 14 | !! zdf_gls : update momentum and tracer Kz from a gls scheme |
---|
| 15 | !! zdf_gls_init : initialization, namelist read, and parameters control |
---|
| 16 | !! gls_rst : read/write gls restart in ocean restart file |
---|
[2048] | 17 | !!---------------------------------------------------------------------- |
---|
| 18 | USE oce ! ocean dynamics and active tracers |
---|
| 19 | USE dom_oce ! ocean space and time domain |
---|
| 20 | USE domvvl ! ocean space and time domain : variable volume layer |
---|
| 21 | USE zdf_oce ! ocean vertical physics |
---|
[5109] | 22 | USE zdfbfr ! bottom friction (only for rn_bfrz0) |
---|
[2048] | 23 | USE sbc_oce ! surface boundary condition: ocean |
---|
| 24 | USE phycst ! physical constants |
---|
| 25 | USE zdfmxl ! mixed layer |
---|
[11277] | 26 | USE sbcwave |
---|
| 27 | ! |
---|
[2048] | 28 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
[2715] | 29 | USE lib_mpp ! MPP manager |
---|
[3294] | 30 | USE wrk_nemo ! work arrays |
---|
[2048] | 31 | USE prtctl ! Print control |
---|
| 32 | USE in_out_manager ! I/O manager |
---|
| 33 | USE iom ! I/O manager library |
---|
[3294] | 34 | USE timing ! Timing |
---|
[3625] | 35 | USE lib_fortran ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) |
---|
[2048] | 36 | |
---|
| 37 | IMPLICIT NONE |
---|
| 38 | PRIVATE |
---|
| 39 | |
---|
[2329] | 40 | PUBLIC zdf_gls ! routine called in step module |
---|
[2397] | 41 | PUBLIC zdf_gls_init ! routine called in opa module |
---|
[2329] | 42 | PUBLIC gls_rst ! routine called in step module |
---|
[2048] | 43 | |
---|
[2715] | 44 | LOGICAL , PUBLIC, PARAMETER :: lk_zdfgls = .TRUE. !: TKE vertical mixing flag |
---|
| 45 | ! |
---|
| 46 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: mxln !: now mixing length |
---|
| 47 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: zwall !: wall function |
---|
| 48 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ustars2 !: Squared surface velocity scale at T-points |
---|
| 49 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ustarb2 !: Squared bottom velocity scale at T-points |
---|
[2048] | 50 | |
---|
[11277] | 51 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: rsbc_tke1 |
---|
| 52 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: rsbc_tke3 |
---|
| 53 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: rsbc_psi1 |
---|
| 54 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: rtrans |
---|
| 55 | |
---|
[4147] | 56 | ! !! ** Namelist namzdf_gls ** |
---|
| 57 | LOGICAL :: ln_length_lim ! use limit on the dissipation rate under stable stratification (Galperin et al. 1988) |
---|
| 58 | LOGICAL :: ln_sigpsi ! Activate Burchard (2003) modification for k-eps closure & wave breaking mixing |
---|
[5109] | 59 | INTEGER :: nn_bc_surf ! surface boundary condition (=0/1) |
---|
| 60 | INTEGER :: nn_bc_bot ! bottom boundary condition (=0/1) |
---|
| 61 | INTEGER :: nn_z0_met ! Method for surface roughness computation |
---|
[4147] | 62 | INTEGER :: nn_stab_func ! stability functions G88, KC or Canuto (=0/1/2) |
---|
| 63 | INTEGER :: nn_clos ! closure 0/1/2/3 MY82/k-eps/k-w/gen |
---|
| 64 | REAL(wp) :: rn_clim_galp ! Holt 2008 value for k-eps: 0.267 |
---|
| 65 | REAL(wp) :: rn_epsmin ! minimum value of dissipation (m2/s3) |
---|
| 66 | REAL(wp) :: rn_emin ! minimum value of TKE (m2/s2) |
---|
| 67 | REAL(wp) :: rn_charn ! Charnock constant for surface breaking waves mixing : 1400. (standard) or 2.e5 (Stacey value) |
---|
[11277] | 68 | REAL(wp) :: rn_crban_default ! Craig and Banner constant for surface breaking waves mixing |
---|
[5109] | 69 | REAL(wp) :: rn_hsro ! Minimum surface roughness |
---|
[11277] | 70 | REAL(wp) :: rn_frac_hs ! Fraction of wave height as surface roughness (if nn_z0_met = 1, 3) |
---|
[2048] | 71 | |
---|
[2397] | 72 | REAL(wp) :: rcm_sf = 0.73_wp ! Shear free turbulence parameters |
---|
| 73 | REAL(wp) :: ra_sf = -2.0_wp ! Must be negative -2 < ra_sf < -1 |
---|
| 74 | REAL(wp) :: rl_sf = 0.2_wp ! 0 <rl_sf<vkarmn |
---|
| 75 | REAL(wp) :: rghmin = -0.28_wp |
---|
| 76 | REAL(wp) :: rgh0 = 0.0329_wp |
---|
| 77 | REAL(wp) :: rghcri = 0.03_wp |
---|
[2299] | 78 | REAL(wp) :: ra1 = 0.92_wp |
---|
| 79 | REAL(wp) :: ra2 = 0.74_wp |
---|
| 80 | REAL(wp) :: rb1 = 16.60_wp |
---|
| 81 | REAL(wp) :: rb2 = 10.10_wp |
---|
| 82 | REAL(wp) :: re2 = 1.33_wp |
---|
| 83 | REAL(wp) :: rl1 = 0.107_wp |
---|
| 84 | REAL(wp) :: rl2 = 0.0032_wp |
---|
| 85 | REAL(wp) :: rl3 = 0.0864_wp |
---|
| 86 | REAL(wp) :: rl4 = 0.12_wp |
---|
| 87 | REAL(wp) :: rl5 = 11.9_wp |
---|
| 88 | REAL(wp) :: rl6 = 0.4_wp |
---|
| 89 | REAL(wp) :: rl7 = 0.0_wp |
---|
| 90 | REAL(wp) :: rl8 = 0.48_wp |
---|
| 91 | REAL(wp) :: rm1 = 0.127_wp |
---|
| 92 | REAL(wp) :: rm2 = 0.00336_wp |
---|
| 93 | REAL(wp) :: rm3 = 0.0906_wp |
---|
| 94 | REAL(wp) :: rm4 = 0.101_wp |
---|
| 95 | REAL(wp) :: rm5 = 11.2_wp |
---|
| 96 | REAL(wp) :: rm6 = 0.4_wp |
---|
| 97 | REAL(wp) :: rm7 = 0.0_wp |
---|
| 98 | REAL(wp) :: rm8 = 0.318_wp |
---|
[11277] | 99 | REAL(wp) :: rtrans_default = 0.1_wp |
---|
[2397] | 100 | REAL(wp) :: rc02, rc02r, rc03, rc04 ! coefficients deduced from above parameters |
---|
[11277] | 101 | REAL(wp) :: rsbc_tke1_default, rsbc_tke2, rfact_tke ! - - - - |
---|
| 102 | REAL(wp) :: rsbc_psi2, rsbc_psi3, rfact_psi ! - - - - |
---|
[5109] | 103 | REAL(wp) :: rsbc_zs1, rsbc_zs2 ! - - - - |
---|
[2397] | 104 | REAL(wp) :: rc0, rc2, rc3, rf6, rcff, rc_diff ! - - - - |
---|
| 105 | REAL(wp) :: rs0, rs1, rs2, rs4, rs5, rs6 ! - - - - |
---|
| 106 | REAL(wp) :: rd0, rd1, rd2, rd3, rd4, rd5 ! - - - - |
---|
| 107 | REAL(wp) :: rsc_tke, rsc_psi, rpsi1, rpsi2, rpsi3, rsc_psi0 ! - - - - |
---|
| 108 | REAL(wp) :: rpsi3m, rpsi3p, rpp, rmm, rnn ! - - - - |
---|
[2299] | 109 | |
---|
[2048] | 110 | !! * Substitutions |
---|
| 111 | # include "domzgr_substitute.h90" |
---|
| 112 | # include "vectopt_loop_substitute.h90" |
---|
| 113 | !!---------------------------------------------------------------------- |
---|
[2287] | 114 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
[2715] | 115 | !! $Id$ |
---|
[2329] | 116 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
[2048] | 117 | !!---------------------------------------------------------------------- |
---|
| 118 | CONTAINS |
---|
| 119 | |
---|
[2715] | 120 | INTEGER FUNCTION zdf_gls_alloc() |
---|
| 121 | !!---------------------------------------------------------------------- |
---|
| 122 | !! *** FUNCTION zdf_gls_alloc *** |
---|
| 123 | !!---------------------------------------------------------------------- |
---|
[6204] | 124 | ALLOCATE( mxln(jpi,jpj,jpk), zwall(jpi,jpj,jpk) , & |
---|
[11277] | 125 | & rsbc_tke1(jpi,jpj), rsbc_tke3(jpi,jpj) , & |
---|
| 126 | & rsbc_psi1(jpi,jpj), rtrans(jpi,jpj), & |
---|
| 127 | & ustars2(jpi,jpj), ustarb2(jpi,jpj) , STAT= zdf_gls_alloc ) |
---|
[2715] | 128 | ! |
---|
| 129 | IF( lk_mpp ) CALL mpp_sum ( zdf_gls_alloc ) |
---|
| 130 | IF( zdf_gls_alloc /= 0 ) CALL ctl_warn('zdf_gls_alloc: failed to allocate arrays') |
---|
| 131 | END FUNCTION zdf_gls_alloc |
---|
| 132 | |
---|
| 133 | |
---|
[2048] | 134 | SUBROUTINE zdf_gls( kt ) |
---|
| 135 | !!---------------------------------------------------------------------- |
---|
| 136 | !! *** ROUTINE zdf_gls *** |
---|
| 137 | !! |
---|
| 138 | !! ** Purpose : Compute the vertical eddy viscosity and diffusivity |
---|
[2397] | 139 | !! coefficients using the GLS turbulent closure scheme. |
---|
[2048] | 140 | !!---------------------------------------------------------------------- |
---|
| 141 | INTEGER, INTENT(in) :: kt ! ocean time step |
---|
| 142 | INTEGER :: ji, jj, jk, ibot, ibotm1, dir ! dummy loop arguments |
---|
[2397] | 143 | REAL(wp) :: zesh2, zsigpsi, zcoef, zex1, zex2 ! local scalars |
---|
| 144 | REAL(wp) :: ztx2, zty2, zup, zdown, zcof ! - - |
---|
| 145 | REAL(wp) :: zratio, zrn2, zflxb, sh ! - - |
---|
| 146 | REAL(wp) :: prod, buoy, diss, zdiss, sm ! - - |
---|
| 147 | REAL(wp) :: gh, gm, shr, dif, zsqen, zav ! - - |
---|
[3294] | 148 | REAL(wp), POINTER, DIMENSION(:,: ) :: zdep |
---|
[5109] | 149 | REAL(wp), POINTER, DIMENSION(:,: ) :: zkar |
---|
[3294] | 150 | REAL(wp), POINTER, DIMENSION(:,: ) :: zflxs ! Turbulence fluxed induced by internal waves |
---|
| 151 | REAL(wp), POINTER, DIMENSION(:,: ) :: zhsro ! Surface roughness (surface waves) |
---|
| 152 | REAL(wp), POINTER, DIMENSION(:,:,:) :: eb ! tke at time before |
---|
| 153 | REAL(wp), POINTER, DIMENSION(:,:,:) :: mxlb ! mixing length at time before |
---|
| 154 | REAL(wp), POINTER, DIMENSION(:,:,:) :: shear ! vertical shear |
---|
| 155 | REAL(wp), POINTER, DIMENSION(:,:,:) :: eps ! dissipation rate |
---|
[5109] | 156 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zwall_psi ! Wall function use in the wb case (ln_sigpsi) |
---|
| 157 | REAL(wp), POINTER, DIMENSION(:,:,:) :: psi ! psi at time now |
---|
| 158 | REAL(wp), POINTER, DIMENSION(:,:,:) :: z_elem_a ! element of the first matrix diagonal |
---|
| 159 | REAL(wp), POINTER, DIMENSION(:,:,:) :: z_elem_b ! element of the second matrix diagonal |
---|
| 160 | REAL(wp), POINTER, DIMENSION(:,:,:) :: z_elem_c ! element of the third matrix diagonal |
---|
[2048] | 161 | !!-------------------------------------------------------------------- |
---|
[3294] | 162 | ! |
---|
| 163 | IF( nn_timing == 1 ) CALL timing_start('zdf_gls') |
---|
| 164 | ! |
---|
[5109] | 165 | CALL wrk_alloc( jpi,jpj, zdep, zkar, zflxs, zhsro ) |
---|
| 166 | CALL wrk_alloc( jpi,jpj,jpk, eb, mxlb, shear, eps, zwall_psi, z_elem_a, z_elem_b, z_elem_c, psi ) |
---|
| 167 | |
---|
[2048] | 168 | ! Preliminary computing |
---|
| 169 | |
---|
[2397] | 170 | ustars2 = 0._wp ; ustarb2 = 0._wp ; psi = 0._wp ; zwall_psi = 0._wp |
---|
[2048] | 171 | |
---|
[11277] | 172 | ! variable initialization |
---|
| 173 | IF( ln_phioc ) THEN |
---|
| 174 | IF( nn_wmix==jp_janssen ) THEN |
---|
| 175 | rsbc_tke1(:,:) = (-rsc_tke*rn_crban(:,:)/(rcm_sf*ra_sf*rl_sf))**(2._wp/3._wp) ! k_eps = 53.Dirichlet + Wave breaking |
---|
| 176 | rsbc_tke3(:,:) = rdt * rn_crban(:,:) ! Neumann + Wave breaking |
---|
| 177 | rsbc_psi1(:,:) = rc0**rpp * rsbc_tke1(:,:)**rmm * rl_sf**rnn ! Dirichlet + Wave breaking |
---|
| 178 | ELSE IF( nn_wmix==jp_craigbanner ) THEN |
---|
| 179 | rsbc_tke1(:,:) = -3._wp/2._wp*rn_crban(:,:)*ra_sf*rl_sf |
---|
| 180 | rsbc_tke3(:,:) = rdt * rn_crban(:,:) / rl_sf |
---|
| 181 | rtrans(:,:) = 0.2_wp / MAX( rsmall, rsbc_tke1(:,:)**(1./(-ra_sf*3._wp/2._wp))-1._wp ) |
---|
| 182 | ENDIF |
---|
| 183 | ENDIF |
---|
| 184 | |
---|
[3798] | 185 | IF( kt /= nit000 ) THEN ! restore before value to compute tke |
---|
| 186 | avt (:,:,:) = avt_k (:,:,:) |
---|
| 187 | avm (:,:,:) = avm_k (:,:,:) |
---|
| 188 | avmu(:,:,:) = avmu_k(:,:,:) |
---|
| 189 | avmv(:,:,:) = avmv_k(:,:,:) |
---|
| 190 | ENDIF |
---|
| 191 | |
---|
[2497] | 192 | ! Compute surface and bottom friction at T-points |
---|
[5109] | 193 | !CDIR NOVERRCHK |
---|
| 194 | DO jj = 2, jpjm1 |
---|
| 195 | !CDIR NOVERRCHK |
---|
| 196 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 197 | ! |
---|
| 198 | ! surface friction |
---|
[3625] | 199 | ustars2(ji,jj) = r1_rau0 * taum(ji,jj) * tmask(ji,jj,1) |
---|
[5109] | 200 | ! |
---|
| 201 | ! bottom friction (explicit before friction) |
---|
| 202 | ! Note that we chose here not to bound the friction as in dynbfr) |
---|
| 203 | ztx2 = ( bfrua(ji,jj) * ub(ji,jj,mbku(ji,jj)) + bfrua(ji-1,jj) * ub(ji-1,jj,mbku(ji-1,jj)) ) & |
---|
| 204 | & * ( 1._wp - 0.5_wp * umask(ji,jj,1) * umask(ji-1,jj,1) ) |
---|
| 205 | zty2 = ( bfrva(ji,jj) * vb(ji,jj,mbkv(ji,jj)) + bfrva(ji,jj-1) * vb(ji,jj-1,mbkv(ji,jj-1)) ) & |
---|
| 206 | & * ( 1._wp - 0.5_wp * vmask(ji,jj,1) * vmask(ji,jj-1,1) ) |
---|
| 207 | ustarb2(ji,jj) = SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1) |
---|
| 208 | END DO |
---|
| 209 | END DO |
---|
[2048] | 210 | |
---|
[5109] | 211 | ! Set surface roughness length |
---|
| 212 | SELECT CASE ( nn_z0_met ) |
---|
| 213 | ! |
---|
| 214 | CASE ( 0 ) ! Constant roughness |
---|
| 215 | zhsro(:,:) = rn_hsro |
---|
| 216 | CASE ( 1 ) ! Standard Charnock formula |
---|
| 217 | zhsro(:,:) = MAX(rsbc_zs1 * ustars2(:,:), rn_hsro) |
---|
| 218 | CASE ( 2 ) ! Roughness formulae according to Rascle et al., Ocean Modelling (2008) |
---|
| 219 | zdep(:,:) = 30.*TANH(2.*0.3/(28.*SQRT(MAX(ustars2(:,:),rsmall)))) ! Wave age (eq. 10) |
---|
| 220 | zhsro(:,:) = MAX(rsbc_zs2 * ustars2(:,:) * zdep(:,:)**1.5, rn_hsro) ! zhsro = rn_frac_hs * Hsw (eq. 11) |
---|
[11277] | 221 | CASE ( 3 ) ! Roughness given by the wave model (coupled or read in file) |
---|
| 222 | zhsro = MAX(rn_frac_hs * hsw, rn_hsro) |
---|
[5109] | 223 | END SELECT |
---|
[2048] | 224 | |
---|
| 225 | ! Compute shear and dissipation rate |
---|
| 226 | DO jk = 2, jpkm1 |
---|
| 227 | DO jj = 2, jpjm1 |
---|
| 228 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 229 | avmu(ji,jj,jk) = avmu(ji,jj,jk) * ( un(ji,jj,jk-1) - un(ji,jj,jk) ) & |
---|
| 230 | & * ( ub(ji,jj,jk-1) - ub(ji,jj,jk) ) & |
---|
| 231 | & / ( fse3uw_n(ji,jj,jk) & |
---|
| 232 | & * fse3uw_b(ji,jj,jk) ) |
---|
| 233 | avmv(ji,jj,jk) = avmv(ji,jj,jk) * ( vn(ji,jj,jk-1) - vn(ji,jj,jk) ) & |
---|
| 234 | & * ( vb(ji,jj,jk-1) - vb(ji,jj,jk) ) & |
---|
| 235 | & / ( fse3vw_n(ji,jj,jk) & |
---|
| 236 | & * fse3vw_b(ji,jj,jk) ) |
---|
[2299] | 237 | eps(ji,jj,jk) = rc03 * en(ji,jj,jk) * SQRT(en(ji,jj,jk)) / mxln(ji,jj,jk) |
---|
[2397] | 238 | END DO |
---|
| 239 | END DO |
---|
| 240 | END DO |
---|
[2048] | 241 | ! |
---|
| 242 | ! Lateral boundary conditions (avmu,avmv) (sign unchanged) |
---|
[2397] | 243 | CALL lbc_lnk( avmu, 'U', 1. ) ; CALL lbc_lnk( avmv, 'V', 1. ) |
---|
[2048] | 244 | |
---|
| 245 | ! Save tke at before time step |
---|
| 246 | eb (:,:,:) = en (:,:,:) |
---|
| 247 | mxlb(:,:,:) = mxln(:,:,:) |
---|
| 248 | |
---|
[2397] | 249 | IF( nn_clos == 0 ) THEN ! Mellor-Yamada |
---|
[2048] | 250 | DO jk = 2, jpkm1 |
---|
| 251 | DO jj = 2, jpjm1 |
---|
| 252 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[2450] | 253 | zup = mxln(ji,jj,jk) * fsdepw(ji,jj,mbkt(ji,jj)+1) |
---|
| 254 | zdown = vkarmn * fsdepw(ji,jj,jk) * ( -fsdepw(ji,jj,jk) + fsdepw(ji,jj,mbkt(ji,jj)+1) ) |
---|
[2397] | 255 | zcoef = ( zup / MAX( zdown, rsmall ) ) |
---|
| 256 | zwall (ji,jj,jk) = ( 1._wp + re2 * zcoef*zcoef ) * tmask(ji,jj,jk) |
---|
| 257 | END DO |
---|
| 258 | END DO |
---|
| 259 | END DO |
---|
[2048] | 260 | ENDIF |
---|
| 261 | |
---|
| 262 | !!---------------------------------!! |
---|
| 263 | !! Equation to prognostic k !! |
---|
| 264 | !!---------------------------------!! |
---|
| 265 | ! |
---|
| 266 | ! Now Turbulent kinetic energy (output in en) |
---|
| 267 | ! ------------------------------- |
---|
| 268 | ! Resolution of a tridiagonal linear system by a "methode de chasse" |
---|
| 269 | ! computation from level 2 to jpkm1 (e(1) computed after and e(jpk)=0 ). |
---|
| 270 | ! The surface boundary condition are set after |
---|
| 271 | ! The bottom boundary condition are also set after. In standard e(bottom)=0. |
---|
| 272 | ! z_elem_b : diagonal z_elem_c : upper diagonal z_elem_a : lower diagonal |
---|
| 273 | ! Warning : after this step, en : right hand side of the matrix |
---|
| 274 | |
---|
| 275 | DO jk = 2, jpkm1 |
---|
| 276 | DO jj = 2, jpjm1 |
---|
| 277 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 278 | ! |
---|
| 279 | ! shear prod. at w-point weightened by mask |
---|
| 280 | shear(ji,jj,jk) = ( avmu(ji-1,jj,jk) + avmu(ji,jj,jk) ) / MAX( 1.e0 , umask(ji-1,jj,jk) + umask(ji,jj,jk) ) & |
---|
| 281 | & + ( avmv(ji,jj-1,jk) + avmv(ji,jj,jk) ) / MAX( 1.e0 , vmask(ji,jj-1,jk) + vmask(ji,jj,jk) ) |
---|
| 282 | ! |
---|
| 283 | ! stratif. destruction |
---|
| 284 | buoy = - avt(ji,jj,jk) * rn2(ji,jj,jk) |
---|
| 285 | ! |
---|
| 286 | ! shear prod. - stratif. destruction |
---|
| 287 | diss = eps(ji,jj,jk) |
---|
| 288 | ! |
---|
[2397] | 289 | dir = 0.5_wp + SIGN( 0.5_wp, shear(ji,jj,jk) + buoy ) ! dir =1(=0) if shear(ji,jj,jk)+buoy >0(<0) |
---|
[2048] | 290 | ! |
---|
[2397] | 291 | zesh2 = dir*(shear(ji,jj,jk)+buoy)+(1._wp-dir)*shear(ji,jj,jk) ! production term |
---|
| 292 | zdiss = dir*(diss/en(ji,jj,jk)) +(1._wp-dir)*(diss-buoy)/en(ji,jj,jk) ! dissipation term |
---|
[2048] | 293 | ! |
---|
[2299] | 294 | ! Compute a wall function from 1. to rsc_psi*zwall/rsc_psi0 |
---|
[2048] | 295 | ! Note that as long that Dirichlet boundary conditions are NOT set at the first and last levels (GOTM style) |
---|
| 296 | ! there is no need to set a boundary condition for zwall_psi at the top and bottom boundaries. |
---|
[2299] | 297 | ! Otherwise, this should be rsc_psi/rsc_psi0 |
---|
[2397] | 298 | IF( ln_sigpsi ) THEN |
---|
| 299 | zsigpsi = MIN( 1._wp, zesh2 / eps(ji,jj,jk) ) ! 0. <= zsigpsi <= 1. |
---|
[3294] | 300 | zwall_psi(ji,jj,jk) = rsc_psi / & |
---|
| 301 | & ( zsigpsi * rsc_psi + (1._wp-zsigpsi) * rsc_psi0 / MAX( zwall(ji,jj,jk), 1._wp ) ) |
---|
[2048] | 302 | ELSE |
---|
[2397] | 303 | zwall_psi(ji,jj,jk) = 1._wp |
---|
[2048] | 304 | ENDIF |
---|
| 305 | ! |
---|
| 306 | ! building the matrix |
---|
[2299] | 307 | zcof = rfact_tke * tmask(ji,jj,jk) |
---|
[2048] | 308 | ! |
---|
| 309 | ! lower diagonal |
---|
| 310 | z_elem_a(ji,jj,jk) = zcof * ( avm (ji,jj,jk ) + avm (ji,jj,jk-1) ) & |
---|
| 311 | & / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk ) ) |
---|
| 312 | ! |
---|
| 313 | ! upper diagonal |
---|
| 314 | z_elem_c(ji,jj,jk) = zcof * ( avm (ji,jj,jk+1) + avm (ji,jj,jk ) ) & |
---|
| 315 | & / ( fse3t(ji,jj,jk ) * fse3w(ji,jj,jk) ) |
---|
| 316 | ! |
---|
| 317 | ! diagonal |
---|
[2397] | 318 | z_elem_b(ji,jj,jk) = 1._wp - z_elem_a(ji,jj,jk) - z_elem_c(ji,jj,jk) & |
---|
| 319 | & + rdt * zdiss * tmask(ji,jj,jk) |
---|
[2048] | 320 | ! |
---|
| 321 | ! right hand side in en |
---|
| 322 | en(ji,jj,jk) = en(ji,jj,jk) + rdt * zesh2 * tmask(ji,jj,jk) |
---|
| 323 | END DO |
---|
| 324 | END DO |
---|
| 325 | END DO |
---|
| 326 | ! |
---|
[2397] | 327 | z_elem_b(:,:,jpk) = 1._wp |
---|
[2048] | 328 | ! |
---|
| 329 | ! Set surface condition on zwall_psi (1 at the bottom) |
---|
[5109] | 330 | zwall_psi(:,:,1) = zwall_psi(:,:,2) |
---|
| 331 | zwall_psi(:,:,jpk) = 1. |
---|
| 332 | ! |
---|
[2048] | 333 | ! Surface boundary condition on tke |
---|
| 334 | ! --------------------------------- |
---|
| 335 | ! |
---|
[11277] | 336 | IF( ln_phioc ) THEN |
---|
| 337 | SELECT CASE ( nn_bc_surf ) |
---|
| 338 | ! |
---|
| 339 | CASE ( 0 ) ! Dirichlet case |
---|
| 340 | IF( nn_wmix==jp_janssen ) THEN |
---|
| 341 | ! First level |
---|
| 342 | en(:,:,1) = MAX( rsbc_tke1(:,:) * ustars2(:,:), rn_emin ) |
---|
| 343 | z_elem_a(:,:,1) = en(:,:,1) |
---|
| 344 | z_elem_c(:,:,1) = 0._wp |
---|
| 345 | z_elem_b(:,:,1) = 1._wp |
---|
| 346 | |
---|
| 347 | ! One level below |
---|
| 348 | en(:,:,2) = MAX( rsbc_tke1(:,:) * ustars2(:,:) * ((zhsro(:,:)+fsdepw(:,:,2))/zhsro(:,:) )**ra_sf, rn_emin ) |
---|
| 349 | z_elem_a(:,:,2) = 0._wp |
---|
| 350 | z_elem_c(:,:,2) = 0._wp |
---|
| 351 | z_elem_b(:,:,2) = 1._wp |
---|
| 352 | ELSE IF( nn_wmix==jp_craigbanner ) THEN |
---|
| 353 | en(:,:,1) = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1(:,:))**(2._wp/3._wp) |
---|
| 354 | en(:,:,1) = MAX(en(:,:,1), rn_emin) |
---|
| 355 | z_elem_a(:,:,1) = en(:,:,1) |
---|
| 356 | z_elem_c(:,:,1) = 0._wp |
---|
| 357 | z_elem_b(:,:,1) = 1._wp |
---|
| 358 | ! |
---|
| 359 | ! One level below |
---|
| 360 | en(:,:,2) = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1(:,:) * ((zhsro(:,:)+fsdepw(:,:,2)) & |
---|
| 361 | & / zhsro(:,:) )**(1.5_wp*ra_sf))**(2._wp/3._wp) |
---|
| 362 | en(:,:,2) = MAX(en(:,:,2), rn_emin ) |
---|
| 363 | z_elem_a(:,:,2) = 0._wp |
---|
| 364 | z_elem_c(:,:,2) = 0._wp |
---|
| 365 | z_elem_b(:,:,2) = 1._wp |
---|
| 366 | ! |
---|
| 367 | ENDIF |
---|
| 368 | CASE ( 1 ) ! Neumann boundary condition on d(e)/dz |
---|
| 369 | IF( nn_wmix==jp_janssen ) THEN |
---|
| 370 | ! Dirichlet conditions at k=1 |
---|
| 371 | en(:,:,1) = MAX( rsbc_tke1(:,:) * ustars2(:,:), rn_emin ) |
---|
| 372 | z_elem_a(:,:,1) = en(:,:,1) |
---|
| 373 | z_elem_c(:,:,1) = 0._wp |
---|
| 374 | z_elem_b(:,:,1) = 1._wp |
---|
| 375 | ! at k=2, set de/dz=Fw |
---|
| 376 | z_elem_b(:,:,2) = z_elem_b(:,:,2) + z_elem_a(:,:,2) ! Remove z_elem_a from z_elem_b |
---|
| 377 | z_elem_a(:,:,2) = 0._wp |
---|
| 378 | zflxs(:,:) = rsbc_tke3(:,:) * ustars2(:,:)**1.5_wp * ((zhsro(:,:)+fsdept(:,:,1) ) / zhsro(:,:) )**(1.5*ra_sf) |
---|
| 379 | en(:,:,2) = en(:,:,2) + zflxs(:,:) / fse3w(:,:,2) |
---|
| 380 | ELSE IF( nn_wmix==jp_craigbanner ) THEN |
---|
| 381 | ! Dirichlet conditions at k=1 |
---|
| 382 | en(:,:,1) = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1(:,:))**(2._wp/3._wp) |
---|
| 383 | en(:,:,1) = MAX(en(:,:,1), rn_emin) |
---|
| 384 | z_elem_a(:,:,1) = en(:,:,1) |
---|
| 385 | z_elem_c(:,:,1) = 0._wp |
---|
| 386 | z_elem_b(:,:,1) = 1._wp |
---|
| 387 | ! |
---|
| 388 | ! at k=2, set de/dz=Fw |
---|
| 389 | !cbr |
---|
| 390 | z_elem_b(:,:,2) = z_elem_b(:,:,2) + z_elem_a(:,:,2) ! Remove z_elem_a from z_elem_b |
---|
| 391 | z_elem_a(:,:,2) = 0._wp |
---|
| 392 | zkar(:,:) = (rl_sf + (vkarmn-rl_sf)*(1.-exp(-rtrans(:,:)*fsdept(:,:,1)/zhsro(:,:)) )) |
---|
| 393 | zflxs(:,:) = rsbc_tke3(:,:) * ustars2(:,:)**1.5_wp * zkar(:,:) & |
---|
| 394 | & * ((zhsro(:,:)+fsdept(:,:,1))/zhsro(:,:) )**(1.5_wp*ra_sf) |
---|
| 395 | |
---|
| 396 | en(:,:,2) = en(:,:,2) + zflxs(:,:)/fse3w(:,:,2) |
---|
| 397 | ! |
---|
| 398 | ENDIF |
---|
| 399 | END SELECT |
---|
| 400 | ELSE |
---|
[5109] | 401 | SELECT CASE ( nn_bc_surf ) |
---|
[2048] | 402 | ! |
---|
| 403 | CASE ( 0 ) ! Dirichlet case |
---|
[5109] | 404 | ! First level |
---|
[11277] | 405 | en(:,:,1) = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1_default)**(2._wp/3._wp) |
---|
[5109] | 406 | en(:,:,1) = MAX(en(:,:,1), rn_emin) |
---|
| 407 | z_elem_a(:,:,1) = en(:,:,1) |
---|
| 408 | z_elem_c(:,:,1) = 0._wp |
---|
| 409 | z_elem_b(:,:,1) = 1._wp |
---|
| 410 | ! |
---|
| 411 | ! One level below |
---|
[11277] | 412 | en(:,:,2) = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1_default * ((zhsro(:,:)+fsdepw(:,:,2))/zhsro(:,:) )**(1.5_wp*ra_sf))**(2._wp/3._wp) |
---|
[5109] | 413 | en(:,:,2) = MAX(en(:,:,2), rn_emin ) |
---|
| 414 | z_elem_a(:,:,2) = 0._wp |
---|
| 415 | z_elem_c(:,:,2) = 0._wp |
---|
| 416 | z_elem_b(:,:,2) = 1._wp |
---|
| 417 | ! |
---|
| 418 | ! |
---|
[2048] | 419 | CASE ( 1 ) ! Neumann boundary condition on d(e)/dz |
---|
[5109] | 420 | ! Dirichlet conditions at k=1 |
---|
[11277] | 421 | en(:,:,1) = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1_default)**(2._wp/3._wp) |
---|
[5109] | 422 | en(:,:,1) = MAX(en(:,:,1), rn_emin) |
---|
| 423 | z_elem_a(:,:,1) = en(:,:,1) |
---|
| 424 | z_elem_c(:,:,1) = 0._wp |
---|
| 425 | z_elem_b(:,:,1) = 1._wp |
---|
| 426 | ! |
---|
| 427 | ! at k=2, set de/dz=Fw |
---|
| 428 | !cbr |
---|
| 429 | z_elem_b(:,:,2) = z_elem_b(:,:,2) + z_elem_a(:,:,2) ! Remove z_elem_a from z_elem_b |
---|
| 430 | z_elem_a(:,:,2) = 0._wp |
---|
[11277] | 431 | zkar(:,:) = (rl_sf + (vkarmn-rl_sf)*(1.-exp(-rtrans_default*fsdept(:,:,1)/zhsro(:,:)) )) |
---|
| 432 | zflxs(:,:) = rsbc_tke2 * ustars2(:,:)**1.5_wp * zkar(:,:) * ((zhsro(:,:)+fsdept(:,:,1))/zhsro(:,:) )**(1.5_wp*ra_sf) |
---|
[5109] | 433 | |
---|
| 434 | en(:,:,2) = en(:,:,2) + zflxs(:,:)/fse3w(:,:,2) |
---|
| 435 | ! |
---|
| 436 | ! |
---|
[2048] | 437 | END SELECT |
---|
[11277] | 438 | ENDIF ! ln_phioc |
---|
[2048] | 439 | |
---|
| 440 | ! Bottom boundary condition on tke |
---|
| 441 | ! -------------------------------- |
---|
| 442 | ! |
---|
[5109] | 443 | SELECT CASE ( nn_bc_bot ) |
---|
[2048] | 444 | ! |
---|
| 445 | CASE ( 0 ) ! Dirichlet |
---|
[2397] | 446 | ! ! en(ibot) = u*^2 / Co2 and mxln(ibot) = rn_lmin |
---|
| 447 | ! ! Balance between the production and the dissipation terms |
---|
[2048] | 448 | !CDIR NOVERRCHK |
---|
[2397] | 449 | DO jj = 2, jpjm1 |
---|
[2048] | 450 | !CDIR NOVERRCHK |
---|
[2397] | 451 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[2450] | 452 | ibot = mbkt(ji,jj) + 1 ! k bottom level of w-point |
---|
| 453 | ibotm1 = mbkt(ji,jj) ! k-1 bottom level of w-point but >=1 |
---|
[2397] | 454 | ! |
---|
| 455 | ! Bottom level Dirichlet condition: |
---|
| 456 | z_elem_a(ji,jj,ibot ) = 0._wp |
---|
| 457 | z_elem_c(ji,jj,ibot ) = 0._wp |
---|
| 458 | z_elem_b(ji,jj,ibot ) = 1._wp |
---|
| 459 | en(ji,jj,ibot ) = MAX( rc02r * ustarb2(ji,jj), rn_emin ) |
---|
| 460 | ! |
---|
| 461 | ! Just above last level, Dirichlet condition again |
---|
| 462 | z_elem_a(ji,jj,ibotm1) = 0._wp |
---|
| 463 | z_elem_c(ji,jj,ibotm1) = 0._wp |
---|
| 464 | z_elem_b(ji,jj,ibotm1) = 1._wp |
---|
| 465 | en(ji,jj,ibotm1) = MAX( rc02r * ustarb2(ji,jj), rn_emin ) |
---|
| 466 | END DO |
---|
[2048] | 467 | END DO |
---|
[2397] | 468 | ! |
---|
[2048] | 469 | CASE ( 1 ) ! Neumman boundary condition |
---|
[2397] | 470 | ! |
---|
[2048] | 471 | !CDIR NOVERRCHK |
---|
[2397] | 472 | DO jj = 2, jpjm1 |
---|
[2048] | 473 | !CDIR NOVERRCHK |
---|
[2397] | 474 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[2450] | 475 | ibot = mbkt(ji,jj) + 1 ! k bottom level of w-point |
---|
| 476 | ibotm1 = mbkt(ji,jj) ! k-1 bottom level of w-point but >=1 |
---|
[2397] | 477 | ! |
---|
| 478 | ! Bottom level Dirichlet condition: |
---|
| 479 | z_elem_a(ji,jj,ibot) = 0._wp |
---|
| 480 | z_elem_c(ji,jj,ibot) = 0._wp |
---|
| 481 | z_elem_b(ji,jj,ibot) = 1._wp |
---|
| 482 | en(ji,jj,ibot) = MAX( rc02r * ustarb2(ji,jj), rn_emin ) |
---|
| 483 | ! |
---|
| 484 | ! Just above last level: Neumann condition |
---|
| 485 | z_elem_b(ji,jj,ibotm1) = z_elem_b(ji,jj,ibotm1) + z_elem_c(ji,jj,ibotm1) ! Remove z_elem_c from z_elem_b |
---|
| 486 | z_elem_c(ji,jj,ibotm1) = 0._wp |
---|
| 487 | END DO |
---|
[2048] | 488 | END DO |
---|
[2397] | 489 | ! |
---|
[2048] | 490 | END SELECT |
---|
| 491 | |
---|
| 492 | ! Matrix inversion (en prescribed at surface and the bottom) |
---|
| 493 | ! ---------------------------------------------------------- |
---|
| 494 | ! |
---|
| 495 | DO jk = 2, jpkm1 ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1 |
---|
| 496 | DO jj = 2, jpjm1 |
---|
| 497 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 498 | z_elem_b(ji,jj,jk) = z_elem_b(ji,jj,jk) - z_elem_a(ji,jj,jk) * z_elem_c(ji,jj,jk-1) / z_elem_b(ji,jj,jk-1) |
---|
| 499 | END DO |
---|
| 500 | END DO |
---|
| 501 | END DO |
---|
| 502 | DO jk = 2, jpk ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1 |
---|
| 503 | DO jj = 2, jpjm1 |
---|
| 504 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 505 | z_elem_a(ji,jj,jk) = en(ji,jj,jk) - z_elem_a(ji,jj,jk) / z_elem_b(ji,jj,jk-1) * z_elem_a(ji,jj,jk-1) |
---|
| 506 | END DO |
---|
| 507 | END DO |
---|
| 508 | END DO |
---|
| 509 | DO jk = jpk-1, 2, -1 ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk |
---|
| 510 | DO jj = 2, jpjm1 |
---|
| 511 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 512 | en(ji,jj,jk) = ( z_elem_a(ji,jj,jk) - z_elem_c(ji,jj,jk) * en(ji,jj,jk+1) ) / z_elem_b(ji,jj,jk) |
---|
| 513 | END DO |
---|
| 514 | END DO |
---|
| 515 | END DO |
---|
[2397] | 516 | ! ! set the minimum value of tke |
---|
[2048] | 517 | en(:,:,:) = MAX( en(:,:,:), rn_emin ) |
---|
[5109] | 518 | |
---|
[2048] | 519 | !!----------------------------------------!! |
---|
| 520 | !! Solve prognostic equation for psi !! |
---|
| 521 | !!----------------------------------------!! |
---|
| 522 | |
---|
| 523 | ! Set psi to previous time step value |
---|
| 524 | ! |
---|
| 525 | SELECT CASE ( nn_clos ) |
---|
| 526 | ! |
---|
| 527 | CASE( 0 ) ! k-kl (Mellor-Yamada) |
---|
[2397] | 528 | DO jk = 2, jpkm1 |
---|
| 529 | DO jj = 2, jpjm1 |
---|
| 530 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[3294] | 531 | psi(ji,jj,jk) = eb(ji,jj,jk) * mxlb(ji,jj,jk) |
---|
[2397] | 532 | END DO |
---|
| 533 | END DO |
---|
| 534 | END DO |
---|
| 535 | ! |
---|
[2048] | 536 | CASE( 1 ) ! k-eps |
---|
[2397] | 537 | DO jk = 2, jpkm1 |
---|
| 538 | DO jj = 2, jpjm1 |
---|
| 539 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 540 | psi(ji,jj,jk) = eps(ji,jj,jk) |
---|
| 541 | END DO |
---|
| 542 | END DO |
---|
| 543 | END DO |
---|
| 544 | ! |
---|
[2048] | 545 | CASE( 2 ) ! k-w |
---|
[2397] | 546 | DO jk = 2, jpkm1 |
---|
| 547 | DO jj = 2, jpjm1 |
---|
| 548 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[3294] | 549 | psi(ji,jj,jk) = SQRT( eb(ji,jj,jk) ) / ( rc0 * mxlb(ji,jj,jk) ) |
---|
[2397] | 550 | END DO |
---|
| 551 | END DO |
---|
| 552 | END DO |
---|
| 553 | ! |
---|
| 554 | CASE( 3 ) ! generic |
---|
| 555 | DO jk = 2, jpkm1 |
---|
| 556 | DO jj = 2, jpjm1 |
---|
| 557 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[3294] | 558 | psi(ji,jj,jk) = rc02 * eb(ji,jj,jk) * mxlb(ji,jj,jk)**rnn |
---|
[2397] | 559 | END DO |
---|
| 560 | END DO |
---|
| 561 | END DO |
---|
| 562 | ! |
---|
[2048] | 563 | END SELECT |
---|
| 564 | ! |
---|
| 565 | ! Now gls (output in psi) |
---|
| 566 | ! ------------------------------- |
---|
| 567 | ! Resolution of a tridiagonal linear system by a "methode de chasse" |
---|
| 568 | ! computation from level 2 to jpkm1 (e(1) already computed and e(jpk)=0 ). |
---|
| 569 | ! z_elem_b : diagonal z_elem_c : upper diagonal z_elem_a : lower diagonal |
---|
| 570 | ! Warning : after this step, en : right hand side of the matrix |
---|
| 571 | |
---|
| 572 | DO jk = 2, jpkm1 |
---|
| 573 | DO jj = 2, jpjm1 |
---|
| 574 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 575 | ! |
---|
| 576 | ! psi / k |
---|
| 577 | zratio = psi(ji,jj,jk) / eb(ji,jj,jk) |
---|
| 578 | ! |
---|
| 579 | ! psi3+ : stable : B=-KhN²<0 => N²>0 if rn2>0 dir = 1 (stable) otherwise dir = 0 (unstable) |
---|
[2397] | 580 | dir = 0.5_wp + SIGN( 0.5_wp, rn2(ji,jj,jk) ) |
---|
[2048] | 581 | ! |
---|
[2397] | 582 | rpsi3 = dir * rpsi3m + ( 1._wp - dir ) * rpsi3p |
---|
[2048] | 583 | ! |
---|
| 584 | ! shear prod. - stratif. destruction |
---|
[2299] | 585 | prod = rpsi1 * zratio * shear(ji,jj,jk) |
---|
[2048] | 586 | ! |
---|
| 587 | ! stratif. destruction |
---|
[2397] | 588 | buoy = rpsi3 * zratio * (- avt(ji,jj,jk) * rn2(ji,jj,jk) ) |
---|
[2048] | 589 | ! |
---|
| 590 | ! shear prod. - stratif. destruction |
---|
[2299] | 591 | diss = rpsi2 * zratio * zwall(ji,jj,jk) * eps(ji,jj,jk) |
---|
[2048] | 592 | ! |
---|
[2397] | 593 | dir = 0.5_wp + SIGN( 0.5_wp, prod + buoy ) ! dir =1(=0) if shear(ji,jj,jk)+buoy >0(<0) |
---|
[2048] | 594 | ! |
---|
[2397] | 595 | zesh2 = dir * ( prod + buoy ) + (1._wp - dir ) * prod ! production term |
---|
| 596 | zdiss = dir * ( diss / psi(ji,jj,jk) ) + (1._wp - dir ) * (diss-buoy) / psi(ji,jj,jk) ! dissipation term |
---|
[2048] | 597 | ! |
---|
| 598 | ! building the matrix |
---|
[2299] | 599 | zcof = rfact_psi * zwall_psi(ji,jj,jk) * tmask(ji,jj,jk) |
---|
[2048] | 600 | ! lower diagonal |
---|
| 601 | z_elem_a(ji,jj,jk) = zcof * ( avm (ji,jj,jk ) + avm (ji,jj,jk-1) ) & |
---|
| 602 | & / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk ) ) |
---|
| 603 | ! upper diagonal |
---|
| 604 | z_elem_c(ji,jj,jk) = zcof * ( avm (ji,jj,jk+1) + avm (ji,jj,jk ) ) & |
---|
| 605 | & / ( fse3t(ji,jj,jk ) * fse3w(ji,jj,jk) ) |
---|
| 606 | ! diagonal |
---|
[2397] | 607 | z_elem_b(ji,jj,jk) = 1._wp - z_elem_a(ji,jj,jk) - z_elem_c(ji,jj,jk) & |
---|
| 608 | & + rdt * zdiss * tmask(ji,jj,jk) |
---|
[2048] | 609 | ! |
---|
| 610 | ! right hand side in psi |
---|
| 611 | psi(ji,jj,jk) = psi(ji,jj,jk) + rdt * zesh2 * tmask(ji,jj,jk) |
---|
| 612 | END DO |
---|
| 613 | END DO |
---|
| 614 | END DO |
---|
| 615 | ! |
---|
[2397] | 616 | z_elem_b(:,:,jpk) = 1._wp |
---|
[2048] | 617 | |
---|
| 618 | ! Surface boundary condition on psi |
---|
| 619 | ! --------------------------------- |
---|
| 620 | ! |
---|
[5109] | 621 | SELECT CASE ( nn_bc_surf ) |
---|
[2048] | 622 | ! |
---|
| 623 | CASE ( 0 ) ! Dirichlet boundary conditions |
---|
[11277] | 624 | IF( ln_phioc ) THEN ! Wave induced mixing case |
---|
| 625 | ! en(1) = q2(1) = 0.5 * (15.8 * Ccb)^(2/3) * u*^2 |
---|
| 626 | ! balance between the production and the |
---|
| 627 | ! dissipation terms including the wave effect |
---|
| 628 | IF( nn_wmix==jp_janssen ) THEN |
---|
| 629 | ! Surface value |
---|
| 630 | zdep(:,:) = zhsro(:,:) * rl_sf ! Cosmetic |
---|
| 631 | psi (:,:,1) = rc0**rpp * en(:,:,1)**rmm * zdep(:,:)**rnn * tmask(:,:,1) |
---|
| 632 | z_elem_a(:,:,1) = psi(:,:,1) |
---|
| 633 | z_elem_c(:,:,1) = 0._wp |
---|
| 634 | z_elem_b(:,:,1) = 1._wp |
---|
| 635 | ! |
---|
| 636 | ! One level below |
---|
| 637 | zex1 = (rmm*ra_sf+rnn) |
---|
| 638 | zex2 = (rmm*ra_sf) |
---|
| 639 | zdep(:,:) = ( (zhsro(:,:) + fsdepw(:,:,2))**zex1 ) / zhsro(:,:)**zex2 |
---|
| 640 | psi (:,:,2) = rsbc_psi1(:,:) * ustars2(:,:)**rmm * zdep(:,:) * tmask(:,:,1) |
---|
| 641 | z_elem_a(:,:,2) = 0._wp |
---|
| 642 | z_elem_c(:,:,2) = 0._wp |
---|
| 643 | z_elem_b(:,:,2) = 1._wp |
---|
| 644 | ! |
---|
| 645 | ! |
---|
| 646 | ELSE IF( nn_wmix==jp_craigbanner ) THEN |
---|
| 647 | ! Surface value |
---|
| 648 | zdep(:,:) = zhsro(:,:) * rl_sf ! Cosmetic |
---|
| 649 | psi (:,:,1) = rc0**rpp * en(:,:,1)**rmm * zdep(:,:)**rnn * tmask(:,:,1) |
---|
| 650 | z_elem_a(:,:,1) = psi(:,:,1) |
---|
| 651 | z_elem_c(:,:,1) = 0._wp |
---|
| 652 | z_elem_b(:,:,1) = 1._wp |
---|
| 653 | ! |
---|
| 654 | ! One level below |
---|
| 655 | zkar(:,:) = (rl_sf + (vkarmn-rl_sf)*(1._wp-exp(-rtrans(:,:)*fsdepw(:,:,2)/zhsro(:,:) ))) |
---|
| 656 | zdep(:,:) = (zhsro(:,:) + fsdepw(:,:,2)) * zkar(:,:) |
---|
| 657 | psi (:,:,2) = rc0**rpp * en(:,:,2)**rmm * zdep(:,:)**rnn * tmask(:,:,1) |
---|
| 658 | z_elem_a(:,:,2) = 0._wp |
---|
| 659 | z_elem_c(:,:,2) = 0._wp |
---|
| 660 | z_elem_b(:,:,2) = 1._wp |
---|
| 661 | ! |
---|
| 662 | ! |
---|
| 663 | ENDIF |
---|
| 664 | ELSE ! Wave induced mixing case with default values |
---|
| 665 | ! Surface value |
---|
| 666 | zdep(:,:) = zhsro(:,:) * rl_sf ! Cosmetic |
---|
| 667 | psi (:,:,1) = rc0**rpp * en(:,:,1)**rmm * zdep(:,:)**rnn * tmask(:,:,1) |
---|
| 668 | z_elem_a(:,:,1) = psi(:,:,1) |
---|
| 669 | z_elem_c(:,:,1) = 0._wp |
---|
| 670 | z_elem_b(:,:,1) = 1._wp |
---|
| 671 | ! |
---|
| 672 | ! One level below |
---|
| 673 | zkar(:,:) = (rl_sf + (vkarmn-rl_sf)*(1._wp-exp(-rtrans_default*fsdepw(:,:,2)/zhsro(:,:) ))) |
---|
| 674 | zdep(:,:) = (zhsro(:,:) + fsdepw(:,:,2)) * zkar(:,:) |
---|
| 675 | psi (:,:,2) = rc0**rpp * en(:,:,2)**rmm * zdep(:,:)**rnn * tmask(:,:,1) |
---|
| 676 | z_elem_a(:,:,2) = 0._wp |
---|
| 677 | z_elem_c(:,:,2) = 0._wp |
---|
| 678 | z_elem_b(:,:,2) = 1._wp |
---|
| 679 | ! |
---|
| 680 | ! |
---|
| 681 | ENDIF |
---|
[2048] | 682 | CASE ( 1 ) ! Neumann boundary condition on d(psi)/dz |
---|
[11277] | 683 | IF( ln_phioc ) THEN ! Wave induced mixing case with forced/coupled fields |
---|
| 684 | IF( nn_wmix==jp_janssen ) THEN |
---|
| 685 | ! |
---|
| 686 | zdep(:,:) = rl_sf * zhsro(:,:) |
---|
| 687 | psi (:,:,1) = rc0**rpp * en(:,:,1)**rmm * zdep(:,:)**rnn * tmask(:,:,1) |
---|
| 688 | z_elem_a(:,:,1) = psi(:,:,1) |
---|
| 689 | z_elem_c(:,:,1) = 0._wp |
---|
| 690 | z_elem_b(:,:,1) = 1._wp |
---|
| 691 | ! |
---|
| 692 | ! Neumann condition at k=2 |
---|
| 693 | z_elem_b(:,:,2) = z_elem_b(:,:,2) + z_elem_a(:,:,2) ! Remove z_elem_a from z_elem_b |
---|
| 694 | z_elem_a(:,:,2) = 0._wp |
---|
| 695 | ! |
---|
| 696 | ! Set psi vertical flux at the surface: |
---|
| 697 | zdep(:,:) = (zhsro(:,:) + fsdept(:,:,1))**(rmm*ra_sf+rnn-1._wp) / zhsro(:,:)**(rmm*ra_sf) |
---|
| 698 | zflxs(:,:) = rsbc_psi3 * ( zwall_psi(:,:,1)*avm(:,:,1) + zwall_psi(:,:,2)*avm(:,:,2) ) & |
---|
| 699 | & * en(:,:,1)**rmm * zdep |
---|
| 700 | psi(:,:,2) = psi(:,:,2) + zflxs(:,:) / fse3w(:,:,2) |
---|
| 701 | ! |
---|
| 702 | ELSE IF( nn_wmix==jp_craigbanner ) THEN |
---|
| 703 | ! Surface value: Dirichlet |
---|
| 704 | zdep(:,:) = zhsro(:,:) * rl_sf |
---|
| 705 | psi (:,:,1) = rc0**rpp * en(:,:,1)**rmm * zdep(:,:)**rnn * tmask(:,:,1) |
---|
| 706 | z_elem_a(:,:,1) = psi(:,:,1) |
---|
| 707 | z_elem_c(:,:,1) = 0._wp |
---|
| 708 | z_elem_b(:,:,1) = 1._wp |
---|
| 709 | ! |
---|
| 710 | ! Neumann condition at k=2 |
---|
| 711 | z_elem_b(:,:,2) = z_elem_b(:,:,2) + z_elem_a(:,:,2) ! Remove z_elem_a from z_elem_b |
---|
| 712 | z_elem_a(:,:,2) = 0._wp |
---|
| 713 | ! |
---|
| 714 | ! Set psi vertical flux at the surface: |
---|
| 715 | zkar(:,:) = rl_sf + (vkarmn-rl_sf)*(1._wp-exp(-rtrans(:,:)*fsdept(:,:,1)/zhsro(:,:) )) ! Lengh scale slope |
---|
| 716 | zdep(:,:) = ((zhsro(:,:) + fsdept(:,:,1)) / zhsro(:,:))**(rmm*ra_sf) |
---|
| 717 | zflxs(:,:) = (rnn + rsbc_tke1(:,:) * (rnn + rmm*ra_sf) * zdep(:,:)) * & |
---|
| 718 | (1._wp + rsbc_tke1(:,:) * zdep(:,:))**(2._wp*rmm/3._wp-1_wp) |
---|
| 719 | zdep(:,:) = rsbc_psi1(:,:) * (zwall_psi(:,:,1)*avm(:,:,1)+zwall_psi(:,:,2)*avm(:,:,2)) * & |
---|
| 720 | & ustars2(:,:)**rmm * zkar(:,:)**rnn * (zhsro(:,:) + fsdept(:,:,1))**(rnn-1.) |
---|
| 721 | zflxs(:,:) = zdep(:,:) * zflxs(:,:) |
---|
| 722 | psi(:,:,2) = psi(:,:,2) + zflxs(:,:) / fse3w(:,:,2) |
---|
| 723 | ! |
---|
| 724 | ENDIF |
---|
| 725 | ELSE ! Wave induced mixing case with default values |
---|
| 726 | ! Surface value: Dirichlet |
---|
| 727 | zdep(:,:) = zhsro(:,:) * rl_sf |
---|
| 728 | psi (:,:,1) = rc0**rpp * en(:,:,1)**rmm * zdep(:,:)**rnn * tmask(:,:,1) |
---|
| 729 | z_elem_a(:,:,1) = psi(:,:,1) |
---|
| 730 | z_elem_c(:,:,1) = 0._wp |
---|
| 731 | z_elem_b(:,:,1) = 1._wp |
---|
| 732 | ! |
---|
| 733 | ! Neumann condition at k=2 |
---|
| 734 | z_elem_b(:,:,2) = z_elem_b(:,:,2) + z_elem_a(:,:,2) ! Remove z_elem_a from z_elem_b |
---|
| 735 | z_elem_a(:,:,2) = 0._wp |
---|
| 736 | ! |
---|
| 737 | ! Set psi vertical flux at the surface: |
---|
| 738 | zkar(:,:) = rl_sf + (vkarmn-rl_sf)*(1._wp-exp(-rtrans_default*fsdept(:,:,1)/zhsro(:,:) )) ! Lengh scale slope |
---|
| 739 | zdep(:,:) = ((zhsro(:,:) + fsdept(:,:,1)) / zhsro(:,:))**(rmm*ra_sf) |
---|
| 740 | zflxs(:,:) = (rnn + rsbc_tke1_default * (rnn + rmm*ra_sf) * zdep(:,:)) * & |
---|
| 741 | (1._wp + rsbc_tke1_default * zdep(:,:))**(2._wp*rmm/3._wp-1_wp) |
---|
| 742 | zdep(:,:) = rsbc_psi1(:,:) * (zwall_psi(:,:,1)*avm(:,:,1)+zwall_psi(:,:,2)*avm(:,:,2)) * & |
---|
| 743 | & ustars2(:,:)**rmm * zkar(:,:)**rnn * (zhsro(:,:) + fsdept(:,:,1))**(rnn-1.) |
---|
| 744 | zflxs(:,:) = zdep(:,:) * zflxs(:,:) |
---|
| 745 | psi(:,:,2) = psi(:,:,2) + zflxs(:,:) / fse3w(:,:,2) |
---|
| 746 | ! |
---|
| 747 | ! |
---|
| 748 | ENDIF |
---|
[2048] | 749 | END SELECT |
---|
| 750 | |
---|
| 751 | ! Bottom boundary condition on psi |
---|
| 752 | ! -------------------------------- |
---|
| 753 | ! |
---|
[5109] | 754 | SELECT CASE ( nn_bc_bot ) |
---|
[2048] | 755 | ! |
---|
[5109] | 756 | ! |
---|
[2048] | 757 | CASE ( 0 ) ! Dirichlet |
---|
[5109] | 758 | ! ! en(ibot) = u*^2 / Co2 and mxln(ibot) = vkarmn * rn_bfrz0 |
---|
[2397] | 759 | ! ! Balance between the production and the dissipation terms |
---|
[2048] | 760 | !CDIR NOVERRCHK |
---|
[2397] | 761 | DO jj = 2, jpjm1 |
---|
[2048] | 762 | !CDIR NOVERRCHK |
---|
[2397] | 763 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[2450] | 764 | ibot = mbkt(ji,jj) + 1 ! k bottom level of w-point |
---|
| 765 | ibotm1 = mbkt(ji,jj) ! k-1 bottom level of w-point but >=1 |
---|
[5109] | 766 | zdep(ji,jj) = vkarmn * rn_bfrz0 |
---|
[2397] | 767 | psi (ji,jj,ibot) = rc0**rpp * en(ji,jj,ibot)**rmm * zdep(ji,jj)**rnn |
---|
| 768 | z_elem_a(ji,jj,ibot) = 0._wp |
---|
| 769 | z_elem_c(ji,jj,ibot) = 0._wp |
---|
| 770 | z_elem_b(ji,jj,ibot) = 1._wp |
---|
| 771 | ! |
---|
| 772 | ! Just above last level, Dirichlet condition again (GOTM like) |
---|
[5109] | 773 | zdep(ji,jj) = vkarmn * ( rn_bfrz0 + fse3t(ji,jj,ibotm1) ) |
---|
[2397] | 774 | psi (ji,jj,ibotm1) = rc0**rpp * en(ji,jj,ibot )**rmm * zdep(ji,jj)**rnn |
---|
| 775 | z_elem_a(ji,jj,ibotm1) = 0._wp |
---|
| 776 | z_elem_c(ji,jj,ibotm1) = 0._wp |
---|
| 777 | z_elem_b(ji,jj,ibotm1) = 1._wp |
---|
| 778 | END DO |
---|
[2048] | 779 | END DO |
---|
[2397] | 780 | ! |
---|
[2048] | 781 | CASE ( 1 ) ! Neumman boundary condition |
---|
[2397] | 782 | ! |
---|
[2048] | 783 | !CDIR NOVERRCHK |
---|
[2397] | 784 | DO jj = 2, jpjm1 |
---|
[2048] | 785 | !CDIR NOVERRCHK |
---|
[2397] | 786 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[2450] | 787 | ibot = mbkt(ji,jj) + 1 ! k bottom level of w-point |
---|
| 788 | ibotm1 = mbkt(ji,jj) ! k-1 bottom level of w-point but >=1 |
---|
[2397] | 789 | ! |
---|
| 790 | ! Bottom level Dirichlet condition: |
---|
[5109] | 791 | zdep(ji,jj) = vkarmn * rn_bfrz0 |
---|
[2397] | 792 | psi (ji,jj,ibot) = rc0**rpp * en(ji,jj,ibot)**rmm * zdep(ji,jj)**rnn |
---|
| 793 | ! |
---|
| 794 | z_elem_a(ji,jj,ibot) = 0._wp |
---|
| 795 | z_elem_c(ji,jj,ibot) = 0._wp |
---|
| 796 | z_elem_b(ji,jj,ibot) = 1._wp |
---|
| 797 | ! |
---|
| 798 | ! Just above last level: Neumann condition with flux injection |
---|
| 799 | z_elem_b(ji,jj,ibotm1) = z_elem_b(ji,jj,ibotm1) + z_elem_c(ji,jj,ibotm1) ! Remove z_elem_c from z_elem_b |
---|
| 800 | z_elem_c(ji,jj,ibotm1) = 0. |
---|
| 801 | ! |
---|
| 802 | ! Set psi vertical flux at the bottom: |
---|
[5109] | 803 | zdep(ji,jj) = rn_bfrz0 + 0.5_wp*fse3t(ji,jj,ibotm1) |
---|
[2397] | 804 | zflxb = rsbc_psi2 * ( avm(ji,jj,ibot) + avm(ji,jj,ibotm1) ) & |
---|
| 805 | & * (0.5_wp*(en(ji,jj,ibot)+en(ji,jj,ibotm1)))**rmm * zdep(ji,jj)**(rnn-1._wp) |
---|
| 806 | psi(ji,jj,ibotm1) = psi(ji,jj,ibotm1) + zflxb / fse3w(ji,jj,ibotm1) |
---|
| 807 | END DO |
---|
[2048] | 808 | END DO |
---|
[2397] | 809 | ! |
---|
[2048] | 810 | END SELECT |
---|
| 811 | |
---|
| 812 | ! Matrix inversion |
---|
| 813 | ! ---------------- |
---|
| 814 | ! |
---|
| 815 | DO jk = 2, jpkm1 ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1 |
---|
| 816 | DO jj = 2, jpjm1 |
---|
| 817 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 818 | z_elem_b(ji,jj,jk) = z_elem_b(ji,jj,jk) - z_elem_a(ji,jj,jk) * z_elem_c(ji,jj,jk-1) / z_elem_b(ji,jj,jk-1) |
---|
| 819 | END DO |
---|
| 820 | END DO |
---|
| 821 | END DO |
---|
| 822 | DO jk = 2, jpk ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1 |
---|
| 823 | DO jj = 2, jpjm1 |
---|
| 824 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 825 | z_elem_a(ji,jj,jk) = psi(ji,jj,jk) - z_elem_a(ji,jj,jk) / z_elem_b(ji,jj,jk-1) * z_elem_a(ji,jj,jk-1) |
---|
| 826 | END DO |
---|
| 827 | END DO |
---|
| 828 | END DO |
---|
| 829 | DO jk = jpk-1, 2, -1 ! Third recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk |
---|
| 830 | DO jj = 2, jpjm1 |
---|
| 831 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 832 | psi(ji,jj,jk) = ( z_elem_a(ji,jj,jk) - z_elem_c(ji,jj,jk) * psi(ji,jj,jk+1) ) / z_elem_b(ji,jj,jk) |
---|
| 833 | END DO |
---|
| 834 | END DO |
---|
| 835 | END DO |
---|
| 836 | |
---|
| 837 | ! Set dissipation |
---|
| 838 | !---------------- |
---|
| 839 | |
---|
| 840 | SELECT CASE ( nn_clos ) |
---|
| 841 | ! |
---|
| 842 | CASE( 0 ) ! k-kl (Mellor-Yamada) |
---|
[2397] | 843 | DO jk = 1, jpkm1 |
---|
| 844 | DO jj = 2, jpjm1 |
---|
| 845 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[5109] | 846 | eps(ji,jj,jk) = rc03 * en(ji,jj,jk) * en(ji,jj,jk) * SQRT( en(ji,jj,jk) ) / MAX( psi(ji,jj,jk), rn_epsmin) |
---|
[2397] | 847 | END DO |
---|
| 848 | END DO |
---|
| 849 | END DO |
---|
| 850 | ! |
---|
[2048] | 851 | CASE( 1 ) ! k-eps |
---|
[2397] | 852 | DO jk = 1, jpkm1 |
---|
| 853 | DO jj = 2, jpjm1 |
---|
| 854 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 855 | eps(ji,jj,jk) = psi(ji,jj,jk) |
---|
| 856 | END DO |
---|
| 857 | END DO |
---|
| 858 | END DO |
---|
| 859 | ! |
---|
[2048] | 860 | CASE( 2 ) ! k-w |
---|
[2397] | 861 | DO jk = 1, jpkm1 |
---|
| 862 | DO jj = 2, jpjm1 |
---|
| 863 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 864 | eps(ji,jj,jk) = rc04 * en(ji,jj,jk) * psi(ji,jj,jk) |
---|
| 865 | END DO |
---|
| 866 | END DO |
---|
| 867 | END DO |
---|
| 868 | ! |
---|
| 869 | CASE( 3 ) ! generic |
---|
| 870 | zcoef = rc0**( 3._wp + rpp/rnn ) |
---|
| 871 | zex1 = ( 1.5_wp + rmm/rnn ) |
---|
| 872 | zex2 = -1._wp / rnn |
---|
| 873 | DO jk = 1, jpkm1 |
---|
| 874 | DO jj = 2, jpjm1 |
---|
| 875 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 876 | eps(ji,jj,jk) = zcoef * en(ji,jj,jk)**zex1 * psi(ji,jj,jk)**zex2 |
---|
| 877 | END DO |
---|
| 878 | END DO |
---|
| 879 | END DO |
---|
| 880 | ! |
---|
[2048] | 881 | END SELECT |
---|
| 882 | |
---|
| 883 | ! Limit dissipation rate under stable stratification |
---|
| 884 | ! -------------------------------------------------- |
---|
| 885 | DO jk = 1, jpkm1 ! Note that this set boundary conditions on mxln at the same time |
---|
| 886 | DO jj = 2, jpjm1 |
---|
| 887 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 888 | ! limitation |
---|
| 889 | eps(ji,jj,jk) = MAX( eps(ji,jj,jk), rn_epsmin ) |
---|
[2397] | 890 | mxln(ji,jj,jk) = rc03 * en(ji,jj,jk) * SQRT( en(ji,jj,jk) ) / eps(ji,jj,jk) |
---|
[2048] | 891 | ! Galperin criterium (NOTE : Not required if the proper value of C3 in stable cases is calculated) |
---|
| 892 | zrn2 = MAX( rn2(ji,jj,jk), rsmall ) |
---|
[5109] | 893 | IF (ln_length_lim) mxln(ji,jj,jk) = MIN( rn_clim_galp * SQRT( 2._wp * en(ji,jj,jk) / zrn2 ), mxln(ji,jj,jk) ) |
---|
[2048] | 894 | END DO |
---|
| 895 | END DO |
---|
| 896 | END DO |
---|
| 897 | |
---|
| 898 | ! |
---|
| 899 | ! Stability function and vertical viscosity and diffusivity |
---|
| 900 | ! --------------------------------------------------------- |
---|
| 901 | ! |
---|
| 902 | SELECT CASE ( nn_stab_func ) |
---|
| 903 | ! |
---|
| 904 | CASE ( 0 , 1 ) ! Galperin or Kantha-Clayson stability functions |
---|
[2397] | 905 | DO jk = 2, jpkm1 |
---|
| 906 | DO jj = 2, jpjm1 |
---|
| 907 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 908 | ! zcof = l²/q² |
---|
| 909 | zcof = mxlb(ji,jj,jk) * mxlb(ji,jj,jk) / ( 2._wp*eb(ji,jj,jk) ) |
---|
| 910 | ! Gh = -N²l²/q² |
---|
| 911 | gh = - rn2(ji,jj,jk) * zcof |
---|
| 912 | gh = MIN( gh, rgh0 ) |
---|
| 913 | gh = MAX( gh, rghmin ) |
---|
| 914 | ! Stability functions from Kantha and Clayson (if C2=C3=0 => Galperin) |
---|
| 915 | sh = ra2*( 1._wp-6._wp*ra1/rb1 ) / ( 1.-3.*ra2*gh*(6.*ra1+rb2*( 1._wp-rc3 ) ) ) |
---|
| 916 | sm = ( rb1**(-1._wp/3._wp) + ( 18._wp*ra1*ra1 + 9._wp*ra1*ra2*(1._wp-rc2) )*sh*gh ) / (1._wp-9._wp*ra1*ra2*gh) |
---|
| 917 | ! |
---|
| 918 | ! Store stability function in avmu and avmv |
---|
| 919 | avmu(ji,jj,jk) = rc_diff * sh * tmask(ji,jj,jk) |
---|
| 920 | avmv(ji,jj,jk) = rc_diff * sm * tmask(ji,jj,jk) |
---|
| 921 | END DO |
---|
[2048] | 922 | END DO |
---|
| 923 | END DO |
---|
[2397] | 924 | ! |
---|
[2048] | 925 | CASE ( 2, 3 ) ! Canuto stability functions |
---|
[2397] | 926 | DO jk = 2, jpkm1 |
---|
| 927 | DO jj = 2, jpjm1 |
---|
| 928 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 929 | ! zcof = l²/q² |
---|
| 930 | zcof = mxlb(ji,jj,jk)*mxlb(ji,jj,jk) / ( 2._wp * eb(ji,jj,jk) ) |
---|
| 931 | ! Gh = -N²l²/q² |
---|
| 932 | gh = - rn2(ji,jj,jk) * zcof |
---|
| 933 | gh = MIN( gh, rgh0 ) |
---|
| 934 | gh = MAX( gh, rghmin ) |
---|
| 935 | gh = gh * rf6 |
---|
| 936 | ! Gm = M²l²/q² Shear number |
---|
| 937 | shr = shear(ji,jj,jk) / MAX( avm(ji,jj,jk), rsmall ) |
---|
| 938 | gm = MAX( shr * zcof , 1.e-10 ) |
---|
| 939 | gm = gm * rf6 |
---|
| 940 | gm = MIN ( (rd0 - rd1*gh + rd3*gh*gh) / (rd2-rd4*gh) , gm ) |
---|
| 941 | ! Stability functions from Canuto |
---|
| 942 | rcff = rd0 - rd1*gh +rd2*gm + rd3*gh*gh - rd4*gh*gm + rd5*gm*gm |
---|
| 943 | sm = (rs0 - rs1*gh + rs2*gm) / rcff |
---|
| 944 | sh = (rs4 - rs5*gh + rs6*gm) / rcff |
---|
| 945 | ! |
---|
| 946 | ! Store stability function in avmu and avmv |
---|
| 947 | avmu(ji,jj,jk) = rc_diff * sh * tmask(ji,jj,jk) |
---|
| 948 | avmv(ji,jj,jk) = rc_diff * sm * tmask(ji,jj,jk) |
---|
| 949 | END DO |
---|
[2048] | 950 | END DO |
---|
| 951 | END DO |
---|
[2397] | 952 | ! |
---|
[2048] | 953 | END SELECT |
---|
| 954 | |
---|
| 955 | ! Boundary conditions on stability functions for momentum (Neumann): |
---|
| 956 | ! Lines below are useless if GOTM style Dirichlet conditions are used |
---|
[5109] | 957 | |
---|
| 958 | avmv(:,:,1) = avmv(:,:,2) |
---|
| 959 | |
---|
[2048] | 960 | DO jj = 2, jpjm1 |
---|
| 961 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[5109] | 962 | avmv(ji,jj,mbkt(ji,jj)+1) = avmv(ji,jj,mbkt(ji,jj)) |
---|
[2048] | 963 | END DO |
---|
| 964 | END DO |
---|
| 965 | |
---|
| 966 | ! Compute diffusivities/viscosities |
---|
| 967 | ! The computation below could be restrained to jk=2 to jpkm1 if GOTM style Dirichlet conditions are used |
---|
| 968 | DO jk = 1, jpk |
---|
| 969 | DO jj = 2, jpjm1 |
---|
| 970 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[2397] | 971 | zsqen = SQRT( 2._wp * en(ji,jj,jk) ) * mxln(ji,jj,jk) |
---|
| 972 | zav = zsqen * avmu(ji,jj,jk) |
---|
| 973 | avt(ji,jj,jk) = MAX( zav, avtb(jk) )*tmask(ji,jj,jk) ! apply mask for zdfmxl routine |
---|
| 974 | zav = zsqen * avmv(ji,jj,jk) |
---|
| 975 | avm(ji,jj,jk) = MAX( zav, avmb(jk) ) ! Note that avm is not masked at the surface and the bottom |
---|
[2048] | 976 | END DO |
---|
| 977 | END DO |
---|
| 978 | END DO |
---|
| 979 | ! |
---|
| 980 | ! Lateral boundary conditions (sign unchanged) |
---|
[2397] | 981 | avt(:,:,1) = 0._wp |
---|
[2048] | 982 | CALL lbc_lnk( avm, 'W', 1. ) ; CALL lbc_lnk( avt, 'W', 1. ) |
---|
| 983 | |
---|
| 984 | DO jk = 2, jpkm1 !* vertical eddy viscosity at u- and v-points |
---|
| 985 | DO jj = 2, jpjm1 |
---|
| 986 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[2397] | 987 | avmu(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji+1,jj ,jk) ) * umask(ji,jj,jk) |
---|
| 988 | avmv(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji ,jj+1,jk) ) * vmask(ji,jj,jk) |
---|
[2048] | 989 | END DO |
---|
| 990 | END DO |
---|
| 991 | END DO |
---|
[2397] | 992 | avmu(:,:,1) = 0._wp ; avmv(:,:,1) = 0._wp ! set surface to zero |
---|
| 993 | CALL lbc_lnk( avmu, 'U', 1. ) ; CALL lbc_lnk( avmv, 'V', 1. ) ! Lateral boundary conditions |
---|
[2048] | 994 | |
---|
| 995 | IF(ln_ctl) THEN |
---|
| 996 | CALL prt_ctl( tab3d_1=en , clinfo1=' gls - e: ', tab3d_2=avt, clinfo2=' t: ', ovlap=1, kdim=jpk) |
---|
| 997 | CALL prt_ctl( tab3d_1=avmu, clinfo1=' gls - u: ', mask1=umask, & |
---|
| 998 | & tab3d_2=avmv, clinfo2= ' v: ', mask2=vmask, ovlap=1, kdim=jpk ) |
---|
| 999 | ENDIF |
---|
| 1000 | ! |
---|
[3798] | 1001 | avt_k (:,:,:) = avt (:,:,:) |
---|
| 1002 | avm_k (:,:,:) = avm (:,:,:) |
---|
| 1003 | avmu_k(:,:,:) = avmu(:,:,:) |
---|
| 1004 | avmv_k(:,:,:) = avmv(:,:,:) |
---|
| 1005 | ! |
---|
[5109] | 1006 | CALL wrk_dealloc( jpi,jpj, zdep, zkar, zflxs, zhsro ) |
---|
[3294] | 1007 | CALL wrk_dealloc( jpi,jpj,jpk, eb, mxlb, shear, eps, zwall_psi, z_elem_a, z_elem_b, z_elem_c, psi ) |
---|
[2715] | 1008 | ! |
---|
[3294] | 1009 | IF( nn_timing == 1 ) CALL timing_stop('zdf_gls') |
---|
| 1010 | ! |
---|
| 1011 | ! |
---|
[2048] | 1012 | END SUBROUTINE zdf_gls |
---|
| 1013 | |
---|
[2329] | 1014 | |
---|
[2048] | 1015 | SUBROUTINE zdf_gls_init |
---|
| 1016 | !!---------------------------------------------------------------------- |
---|
| 1017 | !! *** ROUTINE zdf_gls_init *** |
---|
| 1018 | !! |
---|
| 1019 | !! ** Purpose : Initialization of the vertical eddy diffivity and |
---|
| 1020 | !! viscosity when using a gls turbulent closure scheme |
---|
| 1021 | !! |
---|
| 1022 | !! ** Method : Read the namzdf_gls namelist and check the parameters |
---|
| 1023 | !! called at the first timestep (nit000) |
---|
| 1024 | !! |
---|
| 1025 | !! ** input : Namlist namzdf_gls |
---|
| 1026 | !! |
---|
| 1027 | !! ** Action : Increase by 1 the nstop flag is setting problem encounter |
---|
| 1028 | !! |
---|
| 1029 | !!---------------------------------------------------------------------- |
---|
[2397] | 1030 | USE dynzdf_exp |
---|
| 1031 | USE trazdf_exp |
---|
| 1032 | ! |
---|
[2329] | 1033 | INTEGER :: jk ! dummy loop indices |
---|
[4147] | 1034 | INTEGER :: ios ! Local integer output status for namelist read |
---|
[2329] | 1035 | REAL(wp):: zcr ! local scalar |
---|
[2048] | 1036 | !! |
---|
| 1037 | NAMELIST/namzdf_gls/rn_emin, rn_epsmin, ln_length_lim, & |
---|
[5109] | 1038 | & rn_clim_galp, ln_sigpsi, rn_hsro, & |
---|
[11277] | 1039 | & rn_crban_default, rn_charn, rn_frac_hs,& |
---|
[5109] | 1040 | & nn_bc_surf, nn_bc_bot, nn_z0_met, & |
---|
[2048] | 1041 | & nn_stab_func, nn_clos |
---|
| 1042 | !!---------------------------------------------------------- |
---|
[3294] | 1043 | ! |
---|
| 1044 | IF( nn_timing == 1 ) CALL timing_start('zdf_gls_init') |
---|
| 1045 | ! |
---|
[4147] | 1046 | REWIND( numnam_ref ) ! Namelist namzdf_gls in reference namelist : Vertical eddy diffivity and viscosity using gls turbulent closure scheme |
---|
| 1047 | READ ( numnam_ref, namzdf_gls, IOSTAT = ios, ERR = 901) |
---|
| 1048 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_gls in reference namelist', lwp ) |
---|
[2048] | 1049 | |
---|
[4147] | 1050 | REWIND( numnam_cfg ) ! Namelist namzdf_gls in configuration namelist : Vertical eddy diffivity and viscosity using gls turbulent closure scheme |
---|
| 1051 | READ ( numnam_cfg, namzdf_gls, IOSTAT = ios, ERR = 902 ) |
---|
| 1052 | 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_gls in configuration namelist', lwp ) |
---|
[4624] | 1053 | IF(lwm) WRITE ( numond, namzdf_gls ) |
---|
[4147] | 1054 | |
---|
[2397] | 1055 | IF(lwp) THEN !* Control print |
---|
[2048] | 1056 | WRITE(numout,*) |
---|
| 1057 | WRITE(numout,*) 'zdf_gls_init : gls turbulent closure scheme' |
---|
| 1058 | WRITE(numout,*) '~~~~~~~~~~~~' |
---|
[2397] | 1059 | WRITE(numout,*) ' Namelist namzdf_gls : set gls mixing parameters' |
---|
[5109] | 1060 | WRITE(numout,*) ' minimum value of en rn_emin = ', rn_emin |
---|
| 1061 | WRITE(numout,*) ' minimum value of eps rn_epsmin = ', rn_epsmin |
---|
| 1062 | WRITE(numout,*) ' Limit dissipation rate under stable stratif. ln_length_lim = ', ln_length_lim |
---|
| 1063 | WRITE(numout,*) ' Galperin limit (Standard: 0.53, Holt: 0.26) rn_clim_galp = ', rn_clim_galp |
---|
| 1064 | WRITE(numout,*) ' TKE Surface boundary condition nn_bc_surf = ', nn_bc_surf |
---|
| 1065 | WRITE(numout,*) ' TKE Bottom boundary condition nn_bc_bot = ', nn_bc_bot |
---|
| 1066 | WRITE(numout,*) ' Modify psi Schmidt number (wb case) ln_sigpsi = ', ln_sigpsi |
---|
[11277] | 1067 | WRITE(numout,*) ' Craig and Banner coefficient (default) rn_crban = ', rn_crban_default |
---|
[2397] | 1068 | WRITE(numout,*) ' Charnock coefficient rn_charn = ', rn_charn |
---|
[5109] | 1069 | WRITE(numout,*) ' Surface roughness formula nn_z0_met = ', nn_z0_met |
---|
| 1070 | WRITE(numout,*) ' Wave height frac. (used if nn_z0_met=2) rn_frac_hs = ', rn_frac_hs |
---|
[2397] | 1071 | WRITE(numout,*) ' Stability functions nn_stab_func = ', nn_stab_func |
---|
| 1072 | WRITE(numout,*) ' Type of closure nn_clos = ', nn_clos |
---|
[5109] | 1073 | WRITE(numout,*) ' Surface roughness (m) rn_hsro = ', rn_hsro |
---|
| 1074 | WRITE(numout,*) ' Bottom roughness (m) (nambfr namelist) rn_bfrz0 = ', rn_bfrz0 |
---|
[2048] | 1075 | ENDIF |
---|
| 1076 | |
---|
[2715] | 1077 | ! !* allocate gls arrays |
---|
| 1078 | IF( zdf_gls_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'zdf_gls_init : unable to allocate arrays' ) |
---|
| 1079 | |
---|
[2397] | 1080 | ! !* Check of some namelist values |
---|
[11277] | 1081 | IF( nn_bc_surf < 0 .OR. nn_bc_surf > 1 ) CALL ctl_stop( 'zdf_gls_init: bad flag: nn_bc_surf is 0 or 1' ) |
---|
| 1082 | IF( nn_bc_surf < 0 .OR. nn_bc_surf > 1 ) CALL ctl_stop( 'zdf_gls_init: bad flag: nn_bc_surf is 0 or 1' ) |
---|
| 1083 | IF( nn_z0_met < 0 .OR. nn_z0_met > 3 ) CALL ctl_stop( 'zdf_gls_init: bad flag: nn_z0_met is 0, 1, 2 or 3' ) |
---|
| 1084 | IF( nn_z0_met == 3 .AND. .NOT.ln_rough ) CALL ctl_stop( 'zdf_gls_init: nn_z0_met=3 requires ln_rough=T' ) |
---|
| 1085 | IF( nn_z0_met .NE. 3 .AND. ln_rough ) THEN |
---|
| 1086 | CALL ctl_warn('W A R N I N G: ln_rough=.TRUE. but nn_z0_met is not 3 - resetting nn_z0_met to 3') |
---|
| 1087 | nn_z0_met = 3 |
---|
| 1088 | ENDIF |
---|
| 1089 | IF( nn_stab_func < 0 .OR. nn_stab_func > 3 ) CALL ctl_stop( 'zdf_gls_init: bad flag: nn_stab_func is 0, 1, 2 and 3' ) |
---|
| 1090 | IF( nn_clos < 0 .OR. nn_clos > 3 ) CALL ctl_stop( 'zdf_gls_init: bad flag: nn_clos is 0, 1, 2 or 3' ) |
---|
[2048] | 1091 | |
---|
[2715] | 1092 | SELECT CASE ( nn_clos ) !* set the parameters for the chosen closure |
---|
[2048] | 1093 | ! |
---|
[2715] | 1094 | CASE( 0 ) ! k-kl (Mellor-Yamada) |
---|
[2397] | 1095 | ! |
---|
| 1096 | IF(lwp) WRITE(numout,*) 'The choosen closure is k-kl closed to the classical Mellor-Yamada' |
---|
| 1097 | rpp = 0._wp |
---|
| 1098 | rmm = 1._wp |
---|
| 1099 | rnn = 1._wp |
---|
| 1100 | rsc_tke = 1.96_wp |
---|
| 1101 | rsc_psi = 1.96_wp |
---|
| 1102 | rpsi1 = 0.9_wp |
---|
| 1103 | rpsi3p = 1._wp |
---|
| 1104 | rpsi2 = 0.5_wp |
---|
| 1105 | ! |
---|
[2048] | 1106 | SELECT CASE ( nn_stab_func ) |
---|
[2397] | 1107 | CASE( 0, 1 ) ; rpsi3m = 2.53_wp ! G88 or KC stability functions |
---|
[5109] | 1108 | CASE( 2 ) ; rpsi3m = 2.62_wp ! Canuto A stability functions |
---|
[2397] | 1109 | CASE( 3 ) ; rpsi3m = 2.38 ! Canuto B stability functions (caution : constant not identified) |
---|
| 1110 | END SELECT |
---|
[2048] | 1111 | ! |
---|
[2715] | 1112 | CASE( 1 ) ! k-eps |
---|
[2397] | 1113 | ! |
---|
| 1114 | IF(lwp) WRITE(numout,*) 'The choosen closure is k-eps' |
---|
| 1115 | rpp = 3._wp |
---|
| 1116 | rmm = 1.5_wp |
---|
| 1117 | rnn = -1._wp |
---|
| 1118 | rsc_tke = 1._wp |
---|
[5109] | 1119 | rsc_psi = 1.2_wp ! Schmidt number for psi |
---|
[2397] | 1120 | rpsi1 = 1.44_wp |
---|
| 1121 | rpsi3p = 1._wp |
---|
| 1122 | rpsi2 = 1.92_wp |
---|
| 1123 | ! |
---|
| 1124 | SELECT CASE ( nn_stab_func ) |
---|
| 1125 | CASE( 0, 1 ) ; rpsi3m = -0.52_wp ! G88 or KC stability functions |
---|
| 1126 | CASE( 2 ) ; rpsi3m = -0.629_wp ! Canuto A stability functions |
---|
| 1127 | CASE( 3 ) ; rpsi3m = -0.566 ! Canuto B stability functions |
---|
[2048] | 1128 | END SELECT |
---|
[2397] | 1129 | ! |
---|
[2715] | 1130 | CASE( 2 ) ! k-omega |
---|
[2397] | 1131 | ! |
---|
| 1132 | IF(lwp) WRITE(numout,*) 'The choosen closure is k-omega' |
---|
| 1133 | rpp = -1._wp |
---|
| 1134 | rmm = 0.5_wp |
---|
| 1135 | rnn = -1._wp |
---|
| 1136 | rsc_tke = 2._wp |
---|
| 1137 | rsc_psi = 2._wp |
---|
| 1138 | rpsi1 = 0.555_wp |
---|
| 1139 | rpsi3p = 1._wp |
---|
| 1140 | rpsi2 = 0.833_wp |
---|
| 1141 | ! |
---|
| 1142 | SELECT CASE ( nn_stab_func ) |
---|
| 1143 | CASE( 0, 1 ) ; rpsi3m = -0.58_wp ! G88 or KC stability functions |
---|
| 1144 | CASE( 2 ) ; rpsi3m = -0.64_wp ! Canuto A stability functions |
---|
| 1145 | CASE( 3 ) ; rpsi3m = -0.64_wp ! Canuto B stability functions caution : constant not identified) |
---|
| 1146 | END SELECT |
---|
| 1147 | ! |
---|
[2715] | 1148 | CASE( 3 ) ! generic |
---|
[2397] | 1149 | ! |
---|
| 1150 | IF(lwp) WRITE(numout,*) 'The choosen closure is generic' |
---|
| 1151 | rpp = 2._wp |
---|
| 1152 | rmm = 1._wp |
---|
| 1153 | rnn = -0.67_wp |
---|
| 1154 | rsc_tke = 0.8_wp |
---|
| 1155 | rsc_psi = 1.07_wp |
---|
| 1156 | rpsi1 = 1._wp |
---|
| 1157 | rpsi3p = 1._wp |
---|
| 1158 | rpsi2 = 1.22_wp |
---|
| 1159 | ! |
---|
| 1160 | SELECT CASE ( nn_stab_func ) |
---|
| 1161 | CASE( 0, 1 ) ; rpsi3m = 0.1_wp ! G88 or KC stability functions |
---|
| 1162 | CASE( 2 ) ; rpsi3m = 0.05_wp ! Canuto A stability functions |
---|
| 1163 | CASE( 3 ) ; rpsi3m = 0.05_wp ! Canuto B stability functions caution : constant not identified) |
---|
| 1164 | END SELECT |
---|
| 1165 | ! |
---|
[2048] | 1166 | END SELECT |
---|
| 1167 | |
---|
| 1168 | ! |
---|
[2715] | 1169 | SELECT CASE ( nn_stab_func ) !* set the parameters of the stability functions |
---|
[2048] | 1170 | ! |
---|
[2715] | 1171 | CASE ( 0 ) ! Galperin stability functions |
---|
[2397] | 1172 | ! |
---|
| 1173 | IF(lwp) WRITE(numout,*) 'Stability functions from Galperin' |
---|
| 1174 | rc2 = 0._wp |
---|
| 1175 | rc3 = 0._wp |
---|
| 1176 | rc_diff = 1._wp |
---|
| 1177 | rc0 = 0.5544_wp |
---|
| 1178 | rcm_sf = 0.9884_wp |
---|
| 1179 | rghmin = -0.28_wp |
---|
| 1180 | rgh0 = 0.0233_wp |
---|
| 1181 | rghcri = 0.02_wp |
---|
| 1182 | ! |
---|
[2715] | 1183 | CASE ( 1 ) ! Kantha-Clayson stability functions |
---|
[2397] | 1184 | ! |
---|
| 1185 | IF(lwp) WRITE(numout,*) 'Stability functions from Kantha-Clayson' |
---|
| 1186 | rc2 = 0.7_wp |
---|
| 1187 | rc3 = 0.2_wp |
---|
| 1188 | rc_diff = 1._wp |
---|
| 1189 | rc0 = 0.5544_wp |
---|
| 1190 | rcm_sf = 0.9884_wp |
---|
| 1191 | rghmin = -0.28_wp |
---|
| 1192 | rgh0 = 0.0233_wp |
---|
| 1193 | rghcri = 0.02_wp |
---|
| 1194 | ! |
---|
[2715] | 1195 | CASE ( 2 ) ! Canuto A stability functions |
---|
[2397] | 1196 | ! |
---|
| 1197 | IF(lwp) WRITE(numout,*) 'Stability functions from Canuto A' |
---|
| 1198 | rs0 = 1.5_wp * rl1 * rl5*rl5 |
---|
| 1199 | rs1 = -rl4*(rl6+rl7) + 2._wp*rl4*rl5*(rl1-(1._wp/3._wp)*rl2-rl3) + 1.5_wp*rl1*rl5*rl8 |
---|
| 1200 | rs2 = -(3._wp/8._wp) * rl1*(rl6*rl6-rl7*rl7) |
---|
| 1201 | rs4 = 2._wp * rl5 |
---|
| 1202 | rs5 = 2._wp * rl4 |
---|
| 1203 | rs6 = (2._wp/3._wp) * rl5 * ( 3._wp*rl3*rl3 - rl2*rl2 ) - 0.5_wp * rl5*rl1 * (3._wp*rl3-rl2) & |
---|
| 1204 | & + 0.75_wp * rl1 * ( rl6 - rl7 ) |
---|
| 1205 | rd0 = 3._wp * rl5*rl5 |
---|
| 1206 | rd1 = rl5 * ( 7._wp*rl4 + 3._wp*rl8 ) |
---|
| 1207 | rd2 = rl5*rl5 * ( 3._wp*rl3*rl3 - rl2*rl2 ) - 0.75_wp*(rl6*rl6 - rl7*rl7 ) |
---|
| 1208 | rd3 = rl4 * ( 4._wp*rl4 + 3._wp*rl8) |
---|
| 1209 | rd4 = rl4 * ( rl2 * rl6 - 3._wp*rl3*rl7 - rl5*(rl2*rl2 - rl3*rl3 ) ) + rl5*rl8 * ( 3._wp*rl3*rl3 - rl2*rl2 ) |
---|
| 1210 | rd5 = 0.25_wp * ( rl2*rl2 - 3._wp *rl3*rl3 ) * ( rl6*rl6 - rl7*rl7 ) |
---|
| 1211 | rc0 = 0.5268_wp |
---|
| 1212 | rf6 = 8._wp / (rc0**6._wp) |
---|
| 1213 | rc_diff = SQRT(2._wp) / (rc0**3._wp) |
---|
| 1214 | rcm_sf = 0.7310_wp |
---|
| 1215 | rghmin = -0.28_wp |
---|
| 1216 | rgh0 = 0.0329_wp |
---|
| 1217 | rghcri = 0.03_wp |
---|
| 1218 | ! |
---|
[2715] | 1219 | CASE ( 3 ) ! Canuto B stability functions |
---|
[2397] | 1220 | ! |
---|
| 1221 | IF(lwp) WRITE(numout,*) 'Stability functions from Canuto B' |
---|
| 1222 | rs0 = 1.5_wp * rm1 * rm5*rm5 |
---|
| 1223 | rs1 = -rm4 * (rm6+rm7) + 2._wp * rm4*rm5*(rm1-(1._wp/3._wp)*rm2-rm3) + 1.5_wp * rm1*rm5*rm8 |
---|
| 1224 | rs2 = -(3._wp/8._wp) * rm1 * (rm6*rm6-rm7*rm7 ) |
---|
| 1225 | rs4 = 2._wp * rm5 |
---|
| 1226 | rs5 = 2._wp * rm4 |
---|
| 1227 | rs6 = (2._wp/3._wp) * rm5 * (3._wp*rm3*rm3-rm2*rm2) - 0.5_wp * rm5*rm1*(3._wp*rm3-rm2) + 0.75_wp * rm1*(rm6-rm7) |
---|
| 1228 | rd0 = 3._wp * rm5*rm5 |
---|
| 1229 | rd1 = rm5 * (7._wp*rm4 + 3._wp*rm8) |
---|
| 1230 | rd2 = rm5*rm5 * (3._wp*rm3*rm3 - rm2*rm2) - 0.75_wp * (rm6*rm6 - rm7*rm7) |
---|
| 1231 | rd3 = rm4 * ( 4._wp*rm4 + 3._wp*rm8 ) |
---|
| 1232 | rd4 = rm4 * ( rm2*rm6 -3._wp*rm3*rm7 - rm5*(rm2*rm2 - rm3*rm3) ) + rm5 * rm8 * ( 3._wp*rm3*rm3 - rm2*rm2 ) |
---|
| 1233 | rd5 = 0.25_wp * ( rm2*rm2 - 3._wp*rm3*rm3 ) * ( rm6*rm6 - rm7*rm7 ) |
---|
| 1234 | rc0 = 0.5268_wp !! rc0 = 0.5540_wp (Warner ...) to verify ! |
---|
| 1235 | rf6 = 8._wp / ( rc0**6._wp ) |
---|
| 1236 | rc_diff = SQRT(2._wp)/(rc0**3.) |
---|
| 1237 | rcm_sf = 0.7470_wp |
---|
| 1238 | rghmin = -0.28_wp |
---|
| 1239 | rgh0 = 0.0444_wp |
---|
| 1240 | rghcri = 0.0414_wp |
---|
| 1241 | ! |
---|
[2048] | 1242 | END SELECT |
---|
| 1243 | |
---|
[2715] | 1244 | ! !* Set Schmidt number for psi diffusion in the wave breaking case |
---|
| 1245 | ! ! See Eq. (13) of Carniel et al, OM, 30, 225-239, 2009 |
---|
| 1246 | ! ! or Eq. (17) of Burchard, JPO, 31, 3133-3145, 2001 |
---|
[5109] | 1247 | IF( ln_sigpsi ) THEN |
---|
| 1248 | ra_sf = -1.5 ! Set kinetic energy slope, then deduce rsc_psi and rl_sf |
---|
| 1249 | ! Verification: retrieve Burchard (2001) results by uncomenting the line below: |
---|
| 1250 | ! Note that the results depend on the value of rn_cm_sf which is constant (=rc0) in his work |
---|
| 1251 | ! ra_sf = -SQRT(2./3.*rc0**3./rn_cm_sf*rn_sc_tke)/vkarmn |
---|
| 1252 | rsc_psi0 = rsc_tke/(24.*rpsi2)*(-1.+(4.*rnn + ra_sf*(1.+4.*rmm))**2./(ra_sf**2.)) |
---|
[2048] | 1253 | ELSE |
---|
[2299] | 1254 | rsc_psi0 = rsc_psi |
---|
[2048] | 1255 | ENDIF |
---|
| 1256 | |
---|
[2715] | 1257 | ! !* Shear free turbulence parameters |
---|
[2048] | 1258 | ! |
---|
[5109] | 1259 | ra_sf = -4._wp*rnn*SQRT(rsc_tke) / ( (1._wp+4._wp*rmm)*SQRT(rsc_tke) & |
---|
| 1260 | & - SQRT(rsc_tke + 24._wp*rsc_psi0*rpsi2 ) ) |
---|
[2048] | 1261 | |
---|
[11277] | 1262 | IF( .NOT. ln_phioc .AND. rn_crban_default==0._wp ) THEN |
---|
[5109] | 1263 | rl_sf = vkarmn |
---|
| 1264 | ELSE |
---|
| 1265 | rl_sf = rc0 * SQRT(rc0/rcm_sf) * SQRT( ( (1._wp + 4._wp*rmm + 8._wp*rmm**2_wp)*rsc_tke & |
---|
| 1266 | & + 12._wp * rsc_psi0*rpsi2 - (1._wp + 4._wp*rmm) & |
---|
| 1267 | & *SQRT(rsc_tke*(rsc_tke & |
---|
| 1268 | & + 24._wp*rsc_psi0*rpsi2)) ) & |
---|
| 1269 | & /(12._wp*rnn**2.) & |
---|
| 1270 | & ) |
---|
| 1271 | ENDIF |
---|
| 1272 | |
---|
[2048] | 1273 | ! |
---|
[2715] | 1274 | IF(lwp) THEN !* Control print |
---|
[2048] | 1275 | WRITE(numout,*) |
---|
| 1276 | WRITE(numout,*) 'Limit values' |
---|
| 1277 | WRITE(numout,*) '~~~~~~~~~~~~' |
---|
[2299] | 1278 | WRITE(numout,*) 'Parameter m = ',rmm |
---|
| 1279 | WRITE(numout,*) 'Parameter n = ',rnn |
---|
| 1280 | WRITE(numout,*) 'Parameter p = ',rpp |
---|
| 1281 | WRITE(numout,*) 'rpsi1 = ',rpsi1 |
---|
| 1282 | WRITE(numout,*) 'rpsi2 = ',rpsi2 |
---|
| 1283 | WRITE(numout,*) 'rpsi3m = ',rpsi3m |
---|
| 1284 | WRITE(numout,*) 'rpsi3p = ',rpsi3p |
---|
| 1285 | WRITE(numout,*) 'rsc_tke = ',rsc_tke |
---|
| 1286 | WRITE(numout,*) 'rsc_psi = ',rsc_psi |
---|
| 1287 | WRITE(numout,*) 'rsc_psi0 = ',rsc_psi0 |
---|
| 1288 | WRITE(numout,*) 'rc0 = ',rc0 |
---|
[2048] | 1289 | WRITE(numout,*) |
---|
| 1290 | WRITE(numout,*) 'Shear free turbulence parameters:' |
---|
[2299] | 1291 | WRITE(numout,*) 'rcm_sf = ',rcm_sf |
---|
| 1292 | WRITE(numout,*) 'ra_sf = ',ra_sf |
---|
| 1293 | WRITE(numout,*) 'rl_sf = ',rl_sf |
---|
[2048] | 1294 | WRITE(numout,*) |
---|
| 1295 | ENDIF |
---|
| 1296 | |
---|
[2715] | 1297 | ! !* Constants initialization |
---|
[2397] | 1298 | rc02 = rc0 * rc0 ; rc02r = 1. / rc02 |
---|
| 1299 | rc03 = rc02 * rc0 |
---|
| 1300 | rc04 = rc03 * rc0 |
---|
[11277] | 1301 | rsbc_tke1_default = -3._wp/2._wp*rn_crban_default*ra_sf*rl_sf |
---|
| 1302 | rsbc_tke2 = rdt * rn_crban_default / rl_sf |
---|
| 1303 | zcr = MAX( rsmall, rsbc_tke1_default**(1./(-ra_sf*3._wp/2._wp))-1._wp ) |
---|
| 1304 | rtrans_default = 0.2_wp / zcr |
---|
[5109] | 1305 | rsbc_zs1 = rn_charn/grav ! Charnock formula for surface roughness |
---|
| 1306 | rsbc_zs2 = rn_frac_hs / 0.85_wp / grav * 665._wp ! Rascle formula for surface roughness |
---|
| 1307 | rsbc_psi1 = -0.5_wp * rdt * rc0**(rpp-2._wp*rmm) / rsc_psi |
---|
| 1308 | rsbc_psi2 = -0.5_wp * rdt * rc0**rpp * rnn * vkarmn**rnn / rsc_psi ! Neumann + NO Wave breaking |
---|
[2048] | 1309 | |
---|
[5109] | 1310 | rfact_tke = -0.5_wp / rsc_tke * rdt ! Cst used for the Diffusion term of tke |
---|
| 1311 | rfact_psi = -0.5_wp / rsc_psi * rdt ! Cst used for the Diffusion term of tke |
---|
| 1312 | |
---|
[2397] | 1313 | ! !* Wall proximity function |
---|
[2048] | 1314 | zwall (:,:,:) = 1._wp * tmask(:,:,:) |
---|
| 1315 | |
---|
[2397] | 1316 | ! !* set vertical eddy coef. to the background value |
---|
[2048] | 1317 | DO jk = 1, jpk |
---|
| 1318 | avt (:,:,jk) = avtb(jk) * tmask(:,:,jk) |
---|
| 1319 | avm (:,:,jk) = avmb(jk) * tmask(:,:,jk) |
---|
| 1320 | avmu(:,:,jk) = avmb(jk) * umask(:,:,jk) |
---|
| 1321 | avmv(:,:,jk) = avmb(jk) * vmask(:,:,jk) |
---|
| 1322 | END DO |
---|
[2715] | 1323 | ! |
---|
| 1324 | CALL gls_rst( nit000, 'READ' ) !* read or initialize all required files |
---|
[2048] | 1325 | ! |
---|
[3294] | 1326 | IF( nn_timing == 1 ) CALL timing_stop('zdf_gls_init') |
---|
| 1327 | ! |
---|
[2048] | 1328 | END SUBROUTINE zdf_gls_init |
---|
| 1329 | |
---|
[2329] | 1330 | |
---|
[2048] | 1331 | SUBROUTINE gls_rst( kt, cdrw ) |
---|
[2452] | 1332 | !!--------------------------------------------------------------------- |
---|
| 1333 | !! *** ROUTINE ts_rst *** |
---|
| 1334 | !! |
---|
| 1335 | !! ** Purpose : Read or write TKE file (en) in restart file |
---|
| 1336 | !! |
---|
| 1337 | !! ** Method : use of IOM library |
---|
| 1338 | !! if the restart does not contain TKE, en is either |
---|
| 1339 | !! set to rn_emin or recomputed (nn_igls/=0) |
---|
| 1340 | !!---------------------------------------------------------------------- |
---|
| 1341 | INTEGER , INTENT(in) :: kt ! ocean time-step |
---|
| 1342 | CHARACTER(len=*), INTENT(in) :: cdrw ! "READ"/"WRITE" flag |
---|
| 1343 | ! |
---|
| 1344 | INTEGER :: jit, jk ! dummy loop indices |
---|
[3294] | 1345 | INTEGER :: id1, id2, id3, id4, id5, id6 |
---|
[2452] | 1346 | INTEGER :: ji, jj, ikbu, ikbv |
---|
| 1347 | REAL(wp):: cbx, cby |
---|
| 1348 | !!---------------------------------------------------------------------- |
---|
| 1349 | ! |
---|
| 1350 | IF( TRIM(cdrw) == 'READ' ) THEN ! Read/initialise |
---|
| 1351 | ! ! --------------- |
---|
| 1352 | IF( ln_rstart ) THEN !* Read the restart file |
---|
| 1353 | id1 = iom_varid( numror, 'en' , ldstop = .FALSE. ) |
---|
| 1354 | id2 = iom_varid( numror, 'avt' , ldstop = .FALSE. ) |
---|
| 1355 | id3 = iom_varid( numror, 'avm' , ldstop = .FALSE. ) |
---|
| 1356 | id4 = iom_varid( numror, 'avmu' , ldstop = .FALSE. ) |
---|
| 1357 | id5 = iom_varid( numror, 'avmv' , ldstop = .FALSE. ) |
---|
| 1358 | id6 = iom_varid( numror, 'mxln' , ldstop = .FALSE. ) |
---|
| 1359 | ! |
---|
[3294] | 1360 | IF( MIN( id1, id2, id3, id4, id5, id6 ) > 0 ) THEN ! all required arrays exist |
---|
[2452] | 1361 | CALL iom_get( numror, jpdom_autoglo, 'en' , en ) |
---|
| 1362 | CALL iom_get( numror, jpdom_autoglo, 'avt' , avt ) |
---|
| 1363 | CALL iom_get( numror, jpdom_autoglo, 'avm' , avm ) |
---|
| 1364 | CALL iom_get( numror, jpdom_autoglo, 'avmu' , avmu ) |
---|
| 1365 | CALL iom_get( numror, jpdom_autoglo, 'avmv' , avmv ) |
---|
| 1366 | CALL iom_get( numror, jpdom_autoglo, 'mxln' , mxln ) |
---|
| 1367 | ELSE |
---|
| 1368 | IF(lwp) WRITE(numout,*) ' ===>>>> : previous run without gls scheme, en and mxln computed by iterative loop' |
---|
| 1369 | en (:,:,:) = rn_emin |
---|
[5109] | 1370 | mxln(:,:,:) = 0.05 |
---|
[4839] | 1371 | avt_k (:,:,:) = avt (:,:,:) |
---|
| 1372 | avm_k (:,:,:) = avm (:,:,:) |
---|
| 1373 | avmu_k(:,:,:) = avmu(:,:,:) |
---|
| 1374 | avmv_k(:,:,:) = avmv(:,:,:) |
---|
[2452] | 1375 | DO jit = nit000 + 1, nit000 + 10 ; CALL zdf_gls( jit ) ; END DO |
---|
| 1376 | ENDIF |
---|
| 1377 | ELSE !* Start from rest |
---|
| 1378 | IF(lwp) WRITE(numout,*) ' ===>>>> : Initialisation of en and mxln by background values' |
---|
| 1379 | en (:,:,:) = rn_emin |
---|
[5109] | 1380 | mxln(:,:,:) = 0.05 |
---|
[2452] | 1381 | ENDIF |
---|
| 1382 | ! |
---|
| 1383 | ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN ! Create restart file |
---|
| 1384 | ! ! ------------------- |
---|
| 1385 | IF(lwp) WRITE(numout,*) '---- gls-rst ----' |
---|
[5109] | 1386 | CALL iom_rstput( kt, nitrst, numrow, 'en' , en ) |
---|
[3798] | 1387 | CALL iom_rstput( kt, nitrst, numrow, 'avt' , avt_k ) |
---|
| 1388 | CALL iom_rstput( kt, nitrst, numrow, 'avm' , avm_k ) |
---|
[5109] | 1389 | CALL iom_rstput( kt, nitrst, numrow, 'avmu' , avmu_k ) |
---|
[3798] | 1390 | CALL iom_rstput( kt, nitrst, numrow, 'avmv' , avmv_k ) |
---|
| 1391 | CALL iom_rstput( kt, nitrst, numrow, 'mxln' , mxln ) |
---|
[2452] | 1392 | ! |
---|
| 1393 | ENDIF |
---|
| 1394 | ! |
---|
[2048] | 1395 | END SUBROUTINE gls_rst |
---|
| 1396 | |
---|
| 1397 | #else |
---|
| 1398 | !!---------------------------------------------------------------------- |
---|
| 1399 | !! Dummy module : NO TKE scheme |
---|
| 1400 | !!---------------------------------------------------------------------- |
---|
| 1401 | LOGICAL, PUBLIC, PARAMETER :: lk_zdfgls = .FALSE. !: TKE flag |
---|
| 1402 | CONTAINS |
---|
[2409] | 1403 | SUBROUTINE zdf_gls_init ! Empty routine |
---|
| 1404 | WRITE(*,*) 'zdf_gls_init: You should not have seen this print! error?' |
---|
| 1405 | END SUBROUTINE zdf_gls_init |
---|
[2048] | 1406 | SUBROUTINE zdf_gls( kt ) ! Empty routine |
---|
| 1407 | WRITE(*,*) 'zdf_gls: You should not have seen this print! error?', kt |
---|
| 1408 | END SUBROUTINE zdf_gls |
---|
[2397] | 1409 | SUBROUTINE gls_rst( kt, cdrw ) ! Empty routine |
---|
[2264] | 1410 | INTEGER , INTENT(in) :: kt ! ocean time-step |
---|
| 1411 | CHARACTER(len=*), INTENT(in) :: cdrw ! "READ"/"WRITE" flag |
---|
[2397] | 1412 | WRITE(*,*) 'gls_rst: You should not have seen this print! error?', kt, cdrw |
---|
[2264] | 1413 | END SUBROUTINE gls_rst |
---|
[2048] | 1414 | #endif |
---|
| 1415 | |
---|
| 1416 | !!====================================================================== |
---|
| 1417 | END MODULE zdfgls |
---|
[2397] | 1418 | |
---|