| 1 | module turbdiff_mod |
|---|
| 2 | |
|---|
| 3 | implicit none |
|---|
| 4 | |
|---|
| 5 | contains |
|---|
| 6 | |
|---|
| 7 | subroutine turbdiff(ngrid,nlay,nq, & |
|---|
| 8 | ptimestep,pcapcal, & |
|---|
| 9 | pplay,pplev,pzlay,pzlev,pz0, & |
|---|
| 10 | pu,pv,pt,ppopsk,pq,ptsrf,pemis,pqsurf, & |
|---|
| 11 | pdtfi,pdqfi,pfluxsrf, & |
|---|
| 12 | Pdudif,pdvdif,pdtdif,pdtsrf,sensibFlux,pq2, & |
|---|
| 13 | pdqdif,pdqevap,pdqsdif,flux_u,flux_v) |
|---|
| 14 | |
|---|
| 15 | use radcommon_h, only : sigma, glat |
|---|
| 16 | use surfdat_h, only: dryness |
|---|
| 17 | use comcstfi_mod, only: rcp, g, r, cpp |
|---|
| 18 | use callkeys_mod, only: tracer,nosurf,kmixmin |
|---|
| 19 | use turb_mod, only : ustar |
|---|
| 20 | #ifdef MESOSCALE |
|---|
| 21 | use comm_wrf, only : comm_LATENT_HF |
|---|
| 22 | #endif |
|---|
| 23 | implicit none |
|---|
| 24 | |
|---|
| 25 | !================================================================== |
|---|
| 26 | ! |
|---|
| 27 | ! Purpose |
|---|
| 28 | ! ------- |
|---|
| 29 | ! Turbulent diffusion (mixing) for pot. T, U, V and tracers |
|---|
| 30 | ! |
|---|
| 31 | ! Implicit scheme |
|---|
| 32 | ! We start by adding to variables x the physical tendencies |
|---|
| 33 | ! already computed. We resolve the equation: |
|---|
| 34 | ! |
|---|
| 35 | ! x(t+1) = x(t) + dt * (dx/dt)phys(t) + dt * (dx/dt)difv(t+1) |
|---|
| 36 | ! |
|---|
| 37 | ! Authors |
|---|
| 38 | ! ------- |
|---|
| 39 | ! F. Hourdin, F. Forget, R. Fournier (199X) |
|---|
| 40 | ! R. Wordsworth, B. Charnay (2010) |
|---|
| 41 | ! J. Leconte (2012): To f90 |
|---|
| 42 | ! - Rewritten the diffusion scheme to conserve total enthalpy |
|---|
| 43 | ! by accounting for dissipation of turbulent kinetic energy. |
|---|
| 44 | ! - Accounting for lost mean flow kinetic energy should come soon. |
|---|
| 45 | ! |
|---|
| 46 | !================================================================== |
|---|
| 47 | |
|---|
| 48 | !----------------------------------------------------------------------- |
|---|
| 49 | ! declarations |
|---|
| 50 | ! ------------ |
|---|
| 51 | |
|---|
| 52 | ! arguments |
|---|
| 53 | ! --------- |
|---|
| 54 | INTEGER,INTENT(IN) :: ngrid |
|---|
| 55 | INTEGER,INTENT(IN) :: nlay |
|---|
| 56 | REAL,INTENT(IN) :: ptimestep |
|---|
| 57 | REAL,INTENT(IN) :: pplay(ngrid,nlay),pplev(ngrid,nlay+1) |
|---|
| 58 | REAL,INTENT(IN) :: pzlay(ngrid,nlay),pzlev(ngrid,nlay+1) |
|---|
| 59 | REAL,INTENT(IN) :: pu(ngrid,nlay),pv(ngrid,nlay) |
|---|
| 60 | REAL,INTENT(IN) :: pt(ngrid,nlay),ppopsk(ngrid,nlay) |
|---|
| 61 | REAL,INTENT(IN) :: ptsrf(ngrid) ! surface temperature (K) |
|---|
| 62 | REAL,INTENT(IN) :: pemis(ngrid) |
|---|
| 63 | REAL,INTENT(IN) :: pdtfi(ngrid,nlay) |
|---|
| 64 | REAL,INTENT(IN) :: pfluxsrf(ngrid) |
|---|
| 65 | REAL,INTENT(OUT) :: pdudif(ngrid,nlay),pdvdif(ngrid,nlay) |
|---|
| 66 | REAL,INTENT(OUT) :: pdtdif(ngrid,nlay) |
|---|
| 67 | REAL,INTENT(OUT) :: pdtsrf(ngrid) ! tendency (K/s) on surface temperature |
|---|
| 68 | REAL,INTENT(OUT) :: sensibFlux(ngrid) |
|---|
| 69 | REAL,INTENT(IN) :: pcapcal(ngrid) |
|---|
| 70 | REAL,INTENT(INOUT) :: pq2(ngrid,nlay+1) |
|---|
| 71 | REAL,INTENT(OUT) :: flux_u(ngrid),flux_v(ngrid) |
|---|
| 72 | |
|---|
| 73 | REAL,INTENT(IN) :: pz0 |
|---|
| 74 | |
|---|
| 75 | ! Tracers |
|---|
| 76 | ! -------- |
|---|
| 77 | integer,intent(in) :: nq |
|---|
| 78 | real,intent(in) :: pqsurf(ngrid,nq) |
|---|
| 79 | real,intent(in) :: pq(ngrid,nlay,nq), pdqfi(ngrid,nlay,nq) |
|---|
| 80 | real,intent(out) :: pdqdif(ngrid,nlay,nq) |
|---|
| 81 | real,intent(out) :: pdqsdif(ngrid,nq) |
|---|
| 82 | |
|---|
| 83 | ! local |
|---|
| 84 | ! ----- |
|---|
| 85 | integer ilev,ig,ilay,nlev |
|---|
| 86 | |
|---|
| 87 | REAL z4st,zdplanck(ngrid) |
|---|
| 88 | REAL zkv(ngrid,nlay+1),zkh(ngrid,nlay+1) |
|---|
| 89 | REAL zcdv(ngrid),zcdh(ngrid) |
|---|
| 90 | REAL zcdv_true(ngrid),zcdh_true(ngrid) |
|---|
| 91 | REAL zu(ngrid,nlay),zv(ngrid,nlay) |
|---|
| 92 | REAL zh(ngrid,nlay),zt(ngrid,nlay) |
|---|
| 93 | REAL ztsrf(ngrid) |
|---|
| 94 | REAL z1(ngrid),z2(ngrid) |
|---|
| 95 | REAL zmass(ngrid,nlay) |
|---|
| 96 | REAL zfluxv(ngrid,nlay),zfluxt(ngrid,nlay),zfluxq(ngrid,nlay) |
|---|
| 97 | REAL zb0(ngrid,nlay) |
|---|
| 98 | REAL zExner(ngrid,nlay),zovExner(ngrid,nlay) |
|---|
| 99 | REAL zcv(ngrid,nlay),zdv(ngrid,nlay) !inversion coefficient for winds |
|---|
| 100 | REAL zct(ngrid,nlay),zdt(ngrid,nlay) !inversion coefficient for temperature |
|---|
| 101 | REAL zcq(ngrid,nlay),zdq(ngrid,nlay) !inversion coefficient for tracers |
|---|
| 102 | REAL zcst1 |
|---|
| 103 | REAL zu2!, a |
|---|
| 104 | REAL zcq0(ngrid),zdq0(ngrid) |
|---|
| 105 | REAL zx_alf1(ngrid),zx_alf2(ngrid) |
|---|
| 106 | ! 1D eddy diffusion coefficient |
|---|
| 107 | REAL kzz_eddy(nlay) |
|---|
| 108 | REAL pmin_kzz |
|---|
| 109 | |
|---|
| 110 | LOGICAL,SAVE :: firstcall=.true. |
|---|
| 111 | !$OMP THREADPRIVATE(firstcall) |
|---|
| 112 | |
|---|
| 113 | ! Tracers |
|---|
| 114 | ! ------- |
|---|
| 115 | INTEGER iq |
|---|
| 116 | REAL zq(ngrid,nlay,nq) |
|---|
| 117 | REAL zqnoevap(ngrid,nlay) !special for water case to compute where evaporated water goes. |
|---|
| 118 | REAL pdqevap(ngrid,nlay) !special for water case to compute where evaporated water goes. |
|---|
| 119 | REAL zdmassevap(ngrid) |
|---|
| 120 | REAL rho(ngrid) ! near-surface air density |
|---|
| 121 | |
|---|
| 122 | ! Variables added for implicit latent heat inclusion |
|---|
| 123 | ! -------------------------------------------------- |
|---|
| 124 | real dqsat(ngrid),psat_temp,qsat(ngrid),psat(ngrid) |
|---|
| 125 | |
|---|
| 126 | integer, save :: ivap, iliq, iliq_surf,iice_surf ! also make liq for clarity on surface... |
|---|
| 127 | !$OMP THREADPRIVATE(ivap,iliq,iliq_surf,iice_surf) |
|---|
| 128 | |
|---|
| 129 | real, parameter :: karman=0.4 |
|---|
| 130 | real cd0, roughratio |
|---|
| 131 | |
|---|
| 132 | real dqsdif_total(ngrid) |
|---|
| 133 | real zq0(ngrid) |
|---|
| 134 | |
|---|
| 135 | |
|---|
| 136 | ! Coherence test |
|---|
| 137 | ! -------------- |
|---|
| 138 | |
|---|
| 139 | IF (firstcall) THEN |
|---|
| 140 | sensibFlux(:)=0. |
|---|
| 141 | firstcall=.false. |
|---|
| 142 | ENDIF |
|---|
| 143 | |
|---|
| 144 | !----------------------------------------------------------------------- |
|---|
| 145 | ! 1. Initialisation |
|---|
| 146 | ! ----------------- |
|---|
| 147 | |
|---|
| 148 | nlev=nlay+1 |
|---|
| 149 | |
|---|
| 150 | ! Calculate rho*dz, (P/Ps)**(R/cp) and dt*rho/dz=dt*rho**2 g/dp |
|---|
| 151 | ! with rho=p/RT=p/ (R Theta) (p/ps)**kappa |
|---|
| 152 | ! --------------------------------------------- |
|---|
| 153 | |
|---|
| 154 | DO ilay=1,nlay |
|---|
| 155 | DO ig=1,ngrid |
|---|
| 156 | zmass(ig,ilay)=(pplev(ig,ilay)-pplev(ig,ilay+1))/glat(ig) |
|---|
| 157 | zExner(ig,ilay)=(pplev(ig,ilay)/pplev(ig,1))**rcp |
|---|
| 158 | zovExner(ig,ilay)=1./ppopsk(ig,ilay) |
|---|
| 159 | ENDDO |
|---|
| 160 | ENDDO |
|---|
| 161 | |
|---|
| 162 | zcst1=4.*g*ptimestep/(R*R) |
|---|
| 163 | DO ilev=2,nlev-1 |
|---|
| 164 | DO ig=1,ngrid |
|---|
| 165 | zb0(ig,ilev)=pplev(ig,ilev)/(pt(ig,ilev-1)+pt(ig,ilev)) |
|---|
| 166 | zb0(ig,ilev)=zcst1*zb0(ig,ilev)*zb0(ig,ilev)/(pplay(ig,ilev-1)-pplay(ig,ilev)) |
|---|
| 167 | ENDDO |
|---|
| 168 | ENDDO |
|---|
| 169 | DO ig=1,ngrid |
|---|
| 170 | zb0(ig,1)=ptimestep*pplev(ig,1)/(R*ptsrf(ig)) |
|---|
| 171 | ENDDO |
|---|
| 172 | dqsdif_total(:)=0.0 |
|---|
| 173 | |
|---|
| 174 | !----------------------------------------------------------------------- |
|---|
| 175 | ! 2. Add the physical tendencies computed so far |
|---|
| 176 | ! ---------------------------------------------- |
|---|
| 177 | |
|---|
| 178 | DO ilev=1,nlay |
|---|
| 179 | DO ig=1,ngrid |
|---|
| 180 | zu(ig,ilev)=pu(ig,ilev) |
|---|
| 181 | zv(ig,ilev)=pv(ig,ilev) |
|---|
| 182 | zt(ig,ilev)=pt(ig,ilev)+pdtfi(ig,ilev)*ptimestep |
|---|
| 183 | zh(ig,ilev)=pt(ig,ilev)*zovExner(ig,ilev) !for call vdif_kc, but could be moved and computed there |
|---|
| 184 | ENDDO |
|---|
| 185 | ENDDO |
|---|
| 186 | if(tracer) then |
|---|
| 187 | DO iq =1, nq |
|---|
| 188 | DO ilev=1,nlay |
|---|
| 189 | DO ig=1,ngrid |
|---|
| 190 | zq(ig,ilev,iq)=pq(ig,ilev,iq) + pdqfi(ig,ilev,iq)*ptimestep |
|---|
| 191 | ENDDO |
|---|
| 192 | ENDDO |
|---|
| 193 | ENDDO |
|---|
| 194 | end if |
|---|
| 195 | |
|---|
| 196 | !----------------------------------------------------------------------- |
|---|
| 197 | ! 3. Turbulence scheme |
|---|
| 198 | ! -------------------- |
|---|
| 199 | ! |
|---|
| 200 | ! Source of turbulent kinetic energy at the surface |
|---|
| 201 | ! ------------------------------------------------- |
|---|
| 202 | ! Formula is Cd_0 = (karman / log[1+z1/z0])^2 |
|---|
| 203 | |
|---|
| 204 | DO ig=1,ngrid |
|---|
| 205 | roughratio = 1. + pzlay(ig,1)/pz0 |
|---|
| 206 | cd0 = karman/log(roughratio) |
|---|
| 207 | cd0 = cd0*cd0 |
|---|
| 208 | zcdv_true(ig) = cd0 |
|---|
| 209 | zcdh_true(ig) = cd0 |
|---|
| 210 | if(nosurf)then |
|---|
| 211 | zcdv_true(ig)=0.D+0 !JL12 disable atm/surface momentum flux |
|---|
| 212 | zcdh_true(ig)=0.D+0 !JL12 disable sensible heat flux |
|---|
| 213 | endif |
|---|
| 214 | ENDDO |
|---|
| 215 | |
|---|
| 216 | DO ig=1,ngrid |
|---|
| 217 | zu2=pu(ig,1)*pu(ig,1)+pv(ig,1)*pv(ig,1) |
|---|
| 218 | zcdv(ig)=zcdv_true(ig)*sqrt(zu2) |
|---|
| 219 | zcdh(ig)=zcdh_true(ig)*sqrt(zu2) |
|---|
| 220 | ENDDO |
|---|
| 221 | |
|---|
| 222 | ! Turbulent diffusion coefficients in the boundary layer |
|---|
| 223 | ! ------------------------------------------------------ |
|---|
| 224 | |
|---|
| 225 | call vdif_kc(ngrid,nlay,ptimestep,g,pzlev,pzlay,pu,pv,zh,zcdv_true,pq2,zkv,zkh) !JL12 why not call vdif_kc with updated winds and temperature |
|---|
| 226 | |
|---|
| 227 | ! Adding eddy mixing to mimic 3D general circulation in 1D |
|---|
| 228 | ! R. Wordsworth & F. Forget (2010) |
|---|
| 229 | if ((ngrid.eq.1)) then |
|---|
| 230 | ! kmixmin minimum eddy mix coeff in 1D |
|---|
| 231 | ! set up in inifis_mod.F90 - default value 1.0e-2 |
|---|
| 232 | do ilev=1,nlay |
|---|
| 233 | |
|---|
| 234 | ! Here to code your specific eddy mix coeff in 1D |
|---|
| 235 | ! Earth example that can be uncommented below |
|---|
| 236 | ! ------------------------------------------------- |
|---|
| 237 | ! *====== Earth kzz from Zahnle et al. 2006 ======* |
|---|
| 238 | ! ------------------------------------------------- |
|---|
| 239 | ! if(pzlev(1,ilev).le.11.0e3) then |
|---|
| 240 | ! kzz_eddy(ilev)=10.0 |
|---|
| 241 | ! pmin_kzz=pplev(1,ilev)*exp((pzlev(1,ilev)-11.0e3)*g/(r*zt(1,ilev))) |
|---|
| 242 | ! else |
|---|
| 243 | ! kzz_eddy(ilev)=0.1*(pplev(1,ilev)/pmin_kzz)**(-0.5) |
|---|
| 244 | ! kzz_eddy(ilev)=min(kzz_eddy(ilev),100.0) |
|---|
| 245 | ! endif |
|---|
| 246 | ! do ig=1,ngrid |
|---|
| 247 | ! zkh(ig,ilev) = max(kzz_eddy(ilev),zkh(ig,ilev)) |
|---|
| 248 | ! zkv(ig,ilev) = max(kzz_eddy(ilev),zkv(ig,ilev)) |
|---|
| 249 | ! end do |
|---|
| 250 | |
|---|
| 251 | do ig=1,ngrid |
|---|
| 252 | zkh(ig,ilev) = max(kmixmin,zkh(ig,ilev)) |
|---|
| 253 | zkv(ig,ilev) = max(kmixmin,zkv(ig,ilev)) |
|---|
| 254 | end do |
|---|
| 255 | end do |
|---|
| 256 | end if |
|---|
| 257 | |
|---|
| 258 | !JL12 change zkv at the surface by zcdv to calculate the surface momentum flux properly |
|---|
| 259 | DO ig=1,ngrid |
|---|
| 260 | zkv(ig,1)=zcdv(ig) |
|---|
| 261 | ENDDO |
|---|
| 262 | !we treat only winds, energy and tracers coefficients will be computed with upadted winds |
|---|
| 263 | !JL12 calculate the flux coefficients (tables multiplied elements by elements) |
|---|
| 264 | zfluxv(1:ngrid,1:nlay)=zkv(1:ngrid,1:nlay)*zb0(1:ngrid,1:nlay) |
|---|
| 265 | |
|---|
| 266 | !----------------------------------------------------------------------- |
|---|
| 267 | ! 4. Implicit inversion of u |
|---|
| 268 | ! -------------------------- |
|---|
| 269 | |
|---|
| 270 | ! u(t+1) = u(t) + dt * {(du/dt)phys}(t) + dt * {(du/dt)difv}(t+1) |
|---|
| 271 | ! avec |
|---|
| 272 | ! /zu/ = u(t) + dt * {(du/dt)phys}(t) (voir paragraphe 2.) |
|---|
| 273 | ! et |
|---|
| 274 | ! dt * {(du/dt)difv}(t+1) = dt * {(d/dz)[ Ku (du/dz) ]}(t+1) |
|---|
| 275 | ! donc les entrees sont /zcdv/ pour la condition a la limite sol |
|---|
| 276 | ! et /zkv/ = Ku |
|---|
| 277 | |
|---|
| 278 | DO ig=1,ngrid |
|---|
| 279 | z1(ig)=1./(zmass(ig,nlay)+zfluxv(ig,nlay)) |
|---|
| 280 | zcv(ig,nlay)=zmass(ig,nlay)*zu(ig,nlay)*z1(ig) |
|---|
| 281 | zdv(ig,nlay)=zfluxv(ig,nlay)*z1(ig) |
|---|
| 282 | ENDDO |
|---|
| 283 | |
|---|
| 284 | DO ilay=nlay-1,1,-1 |
|---|
| 285 | DO ig=1,ngrid |
|---|
| 286 | z1(ig)=1./(zmass(ig,ilay)+zfluxv(ig,ilay) + zfluxv(ig,ilay+1)*(1.-zdv(ig,ilay+1))) |
|---|
| 287 | zcv(ig,ilay)=(zmass(ig,ilay)*zu(ig,ilay)+zfluxv(ig,ilay+1)*zcv(ig,ilay+1))*z1(ig) |
|---|
| 288 | zdv(ig,ilay)=zfluxv(ig,ilay)*z1(ig) |
|---|
| 289 | ENDDO |
|---|
| 290 | ENDDO |
|---|
| 291 | |
|---|
| 292 | DO ig=1,ngrid |
|---|
| 293 | zu(ig,1)=zcv(ig,1) |
|---|
| 294 | ENDDO |
|---|
| 295 | DO ilay=2,nlay |
|---|
| 296 | DO ig=1,ngrid |
|---|
| 297 | zu(ig,ilay)=zcv(ig,ilay)+zdv(ig,ilay)*zu(ig,ilay-1) |
|---|
| 298 | ENDDO |
|---|
| 299 | ENDDO |
|---|
| 300 | |
|---|
| 301 | !----------------------------------------------------------------------- |
|---|
| 302 | ! 5. Implicit inversion of v |
|---|
| 303 | ! -------------------------- |
|---|
| 304 | |
|---|
| 305 | ! v(t+1) = v(t) + dt * {(dv/dt)phys}(t) + dt * {(dv/dt)difv}(t+1) |
|---|
| 306 | ! avec |
|---|
| 307 | ! /zv/ = v(t) + dt * {(dv/dt)phys}(t) (voir paragraphe 2.) |
|---|
| 308 | ! et |
|---|
| 309 | ! dt * {(dv/dt)difv}(t+1) = dt * {(d/dz)[ Kv (dv/dz) ]}(t+1) |
|---|
| 310 | ! donc les entrees sont /zcdv/ pour la condition a la limite sol |
|---|
| 311 | ! et /zkv/ = Kv |
|---|
| 312 | |
|---|
| 313 | DO ig=1,ngrid |
|---|
| 314 | z1(ig)=1./(zmass(ig,nlay)+zfluxv(ig,nlay)) |
|---|
| 315 | zcv(ig,nlay)=zmass(ig,nlay)*zv(ig,nlay)*z1(ig) |
|---|
| 316 | zdv(ig,nlay)=zfluxv(ig,nlay)*z1(ig) |
|---|
| 317 | ENDDO |
|---|
| 318 | |
|---|
| 319 | DO ilay=nlay-1,1,-1 |
|---|
| 320 | DO ig=1,ngrid |
|---|
| 321 | z1(ig)=1./(zmass(ig,ilay)+zfluxv(ig,ilay)+zfluxv(ig,ilay+1)*(1.-zdv(ig,ilay+1))) |
|---|
| 322 | zcv(ig,ilay)=(zmass(ig,ilay)*zv(ig,ilay)+zfluxv(ig,ilay+1)*zcv(ig,ilay+1))*z1(ig) |
|---|
| 323 | zdv(ig,ilay)=zfluxv(ig,ilay)*z1(ig) |
|---|
| 324 | ENDDO |
|---|
| 325 | ENDDO |
|---|
| 326 | |
|---|
| 327 | DO ig=1,ngrid |
|---|
| 328 | zv(ig,1)=zcv(ig,1) |
|---|
| 329 | ENDDO |
|---|
| 330 | DO ilay=2,nlay |
|---|
| 331 | DO ig=1,ngrid |
|---|
| 332 | zv(ig,ilay)=zcv(ig,ilay)+zdv(ig,ilay)*zv(ig,ilay-1) |
|---|
| 333 | ENDDO |
|---|
| 334 | ENDDO |
|---|
| 335 | |
|---|
| 336 | ! Calcul of wind stress |
|---|
| 337 | |
|---|
| 338 | DO ig=1,ngrid |
|---|
| 339 | flux_u(ig) = zfluxv(ig,1)/ptimestep*zu(ig,1) |
|---|
| 340 | flux_v(ig) = zfluxv(ig,1)/ptimestep*zv(ig,1) |
|---|
| 341 | ENDDO |
|---|
| 342 | |
|---|
| 343 | |
|---|
| 344 | !---------------------------------------------------------------------------- |
|---|
| 345 | ! 6. Implicit inversion of h, not forgetting the coupling with the ground |
|---|
| 346 | |
|---|
| 347 | ! h(t+1) = h(t) + dt * {(dh/dt)phys}(t) + dt * {(dh/dt)difv}(t+1) |
|---|
| 348 | ! avec |
|---|
| 349 | ! /zh/ = h(t) + dt * {(dh/dt)phys}(t) (voir paragraphe 2.) |
|---|
| 350 | ! et |
|---|
| 351 | ! dt * {(dh/dt)difv}(t+1) = dt * {(d/dz)[ Kh (dh/dz) ]}(t+1) |
|---|
| 352 | ! donc les entrees sont /zcdh/ pour la condition de raccord au sol |
|---|
| 353 | ! et /zkh/ = Kh |
|---|
| 354 | |
|---|
| 355 | ! Using the wind modified by friction for lifting and sublimation |
|---|
| 356 | ! --------------------------------------------------------------- |
|---|
| 357 | DO ig=1,ngrid |
|---|
| 358 | zu2 = zu(ig,1)*zu(ig,1)+zv(ig,1)*zv(ig,1) |
|---|
| 359 | zcdv(ig) = zcdv_true(ig)*sqrt(zu2) |
|---|
| 360 | zcdh(ig) = zcdh_true(ig)*sqrt(zu2) |
|---|
| 361 | zkh(ig,1)= zcdh(ig) |
|---|
| 362 | ustar(ig)=sqrt(zcdv_true(ig))*sqrt(zu2) |
|---|
| 363 | ENDDO |
|---|
| 364 | |
|---|
| 365 | |
|---|
| 366 | ! JL12 calculate the flux coefficients (tables multiplied elements by elements) |
|---|
| 367 | ! --------------------------------------------------------------- |
|---|
| 368 | zfluxq(1:ngrid,1:nlay)=zkh(1:ngrid,1:nlay)*zb0(1:ngrid,1:nlay) !JL12 we save zfluxq which doesn't need the Exner factor |
|---|
| 369 | zfluxt(1:ngrid,1:nlay)=zfluxq(1:ngrid,1:nlay)*zExner(1:ngrid,1:nlay) |
|---|
| 370 | |
|---|
| 371 | DO ig=1,ngrid |
|---|
| 372 | z1(ig)=1./(zmass(ig,nlay)+zfluxt(ig,nlay)*zovExner(ig,nlay)) |
|---|
| 373 | zct(ig,nlay)=zmass(ig,nlay)*zt(ig,nlay)*z1(ig) |
|---|
| 374 | zdt(ig,nlay)=zfluxt(ig,nlay)*zovExner(ig,nlay-1)*z1(ig) |
|---|
| 375 | ENDDO |
|---|
| 376 | |
|---|
| 377 | DO ilay=nlay-1,2,-1 |
|---|
| 378 | DO ig=1,ngrid |
|---|
| 379 | z1(ig)=1./(zmass(ig,ilay)+zfluxt(ig,ilay)*zovExner(ig,ilay)+ & |
|---|
| 380 | zfluxt(ig,ilay+1)*(zovExner(ig,ilay)-zdt(ig,ilay+1)*zovExner(ig,ilay+1))) |
|---|
| 381 | zct(ig,ilay)=(zmass(ig,ilay)*zt(ig,ilay)+zfluxt(ig,ilay+1)*zct(ig,ilay+1)*zovExner(ig,ilay+1))*z1(ig) |
|---|
| 382 | zdt(ig,ilay)=zfluxt(ig,ilay)*z1(ig)*zovExner(ig,ilay-1) |
|---|
| 383 | ENDDO |
|---|
| 384 | ENDDO |
|---|
| 385 | |
|---|
| 386 | !JL12 we treat last point afterward because zovExner(ig,ilay-1) does not exist there |
|---|
| 387 | DO ig=1,ngrid |
|---|
| 388 | z1(ig)=1./(zmass(ig,1)+zfluxt(ig,1)*zovExner(ig,1)+ & |
|---|
| 389 | zfluxt(ig,2)*(zovExner(ig,1)-zdt(ig,2)*zovExner(ig,2))) |
|---|
| 390 | zct(ig,1)=(zmass(ig,1)*zt(ig,1)+zfluxt(ig,2)*zct(ig,2)*zovExner(ig,2))*z1(ig) |
|---|
| 391 | zdt(ig,1)=zfluxt(ig,1)*z1(ig) |
|---|
| 392 | ENDDO |
|---|
| 393 | |
|---|
| 394 | |
|---|
| 395 | ! Calculate (d Planck / dT) at the interface temperature |
|---|
| 396 | ! ------------------------------------------------------ |
|---|
| 397 | |
|---|
| 398 | z4st=4.0*sigma*ptimestep |
|---|
| 399 | DO ig=1,ngrid |
|---|
| 400 | zdplanck(ig)=z4st*pemis(ig)*ptsrf(ig)*ptsrf(ig)*ptsrf(ig) |
|---|
| 401 | ENDDO |
|---|
| 402 | |
|---|
| 403 | ! Calculate temperature tendency at the interface (dry case) |
|---|
| 404 | ! ---------------------------------------------------------- |
|---|
| 405 | ! Sum of fluxes at interface at time t + \delta t gives change in T: |
|---|
| 406 | ! radiative fluxes |
|---|
| 407 | ! turbulent convective (sensible) heat flux |
|---|
| 408 | ! flux (if any) from subsurface |
|---|
| 409 | |
|---|
| 410 | ! if(.not.water) then |
|---|
| 411 | |
|---|
| 412 | DO ig=1,ngrid |
|---|
| 413 | z1(ig)=pcapcal(ig)*ptsrf(ig)+cpp*zfluxt(ig,1)*zct(ig,1)*zovExner(ig,1) & |
|---|
| 414 | + pfluxsrf(ig)*ptimestep + zdplanck(ig)*ptsrf(ig) |
|---|
| 415 | z2(ig) = pcapcal(ig)+zdplanck(ig)+cpp*zfluxt(ig,1)*(1.-zovExner(ig,1)*zdt(ig,1)) |
|---|
| 416 | ztsrf(ig) = z1(ig) / z2(ig) |
|---|
| 417 | pdtsrf(ig) = (ztsrf(ig) - ptsrf(ig))/ptimestep |
|---|
| 418 | zt(ig,1) = zct(ig,1) + zdt(ig,1)*ztsrf(ig) |
|---|
| 419 | ENDDO |
|---|
| 420 | ! JL12 note that the black body radiative flux emitted by the surface has been updated by the implicit scheme |
|---|
| 421 | |
|---|
| 422 | |
|---|
| 423 | ! Recalculate temperature to top of atmosphere, starting from ground |
|---|
| 424 | ! ------------------------------------------------------------------ |
|---|
| 425 | |
|---|
| 426 | DO ilay=2,nlay |
|---|
| 427 | DO ig=1,ngrid |
|---|
| 428 | zt(ig,ilay)=zct(ig,ilay)+zdt(ig,ilay)*zt(ig,ilay-1) |
|---|
| 429 | ENDDO |
|---|
| 430 | ENDDO |
|---|
| 431 | |
|---|
| 432 | ! endif ! not water |
|---|
| 433 | |
|---|
| 434 | !----------------------------------------------------------------------- |
|---|
| 435 | ! TRACERS (no vapour) |
|---|
| 436 | ! ------- |
|---|
| 437 | |
|---|
| 438 | if(tracer) then |
|---|
| 439 | |
|---|
| 440 | ! Calculate vertical flux from the bottom to the first layer (dust) |
|---|
| 441 | ! ----------------------------------------------------------------- |
|---|
| 442 | do ig=1,ngrid |
|---|
| 443 | rho(ig) = zb0(ig,1) /ptimestep |
|---|
| 444 | end do |
|---|
| 445 | |
|---|
| 446 | pdqsdif(:,:)=0. |
|---|
| 447 | |
|---|
| 448 | ! Implicit inversion of q |
|---|
| 449 | ! ----------------------- |
|---|
| 450 | do iq=1,nq |
|---|
| 451 | |
|---|
| 452 | if (iq.ne.ivap) then |
|---|
| 453 | |
|---|
| 454 | DO ig=1,ngrid |
|---|
| 455 | z1(ig)=1./(zmass(ig,nlay)+zfluxq(ig,nlay)) |
|---|
| 456 | zcq(ig,nlay)=zmass(ig,nlay)*zq(ig,nlay,iq)*z1(ig) |
|---|
| 457 | zdq(ig,nlay)=zfluxq(ig,nlay)*z1(ig) |
|---|
| 458 | ENDDO |
|---|
| 459 | |
|---|
| 460 | DO ilay=nlay-1,2,-1 |
|---|
| 461 | DO ig=1,ngrid |
|---|
| 462 | z1(ig)=1./(zmass(ig,ilay)+zfluxq(ig,ilay)+zfluxq(ig,ilay+1)*(1.-zdq(ig,ilay+1))) |
|---|
| 463 | zcq(ig,ilay)=(zmass(ig,ilay)*zq(ig,ilay,iq)+zfluxq(ig,ilay+1)*zcq(ig,ilay+1))*z1(ig) |
|---|
| 464 | zdq(ig,ilay)=zfluxq(ig,ilay)*z1(ig) |
|---|
| 465 | ENDDO |
|---|
| 466 | ENDDO |
|---|
| 467 | |
|---|
| 468 | ! Starting upward calculations for simple tracer mixing (e.g., dust) |
|---|
| 469 | do ig=1,ngrid |
|---|
| 470 | zq(ig,1,iq)=zcq(ig,1) |
|---|
| 471 | end do |
|---|
| 472 | |
|---|
| 473 | do ilay=2,nlay |
|---|
| 474 | do ig=1,ngrid |
|---|
| 475 | zq(ig,ilay,iq)=zcq(ig,ilay)+zdq(ig,ilay)*zq(ig,ilay-1,iq) |
|---|
| 476 | end do |
|---|
| 477 | end do |
|---|
| 478 | |
|---|
| 479 | endif ! if (iq.ne.ivap) |
|---|
| 480 | end do ! of do iq=1,nq |
|---|
| 481 | |
|---|
| 482 | ! Revmoed (AF24): Calculate temperature tendency including latent heat term |
|---|
| 483 | ! and assuming an infinite source of water on the ground |
|---|
| 484 | ! ------------------------------------------------------------------ |
|---|
| 485 | |
|---|
| 486 | endif ! tracer |
|---|
| 487 | |
|---|
| 488 | |
|---|
| 489 | !----------------------------------------------------------------------- |
|---|
| 490 | ! 8. Final calculation of the vertical diffusion tendencies |
|---|
| 491 | ! ----------------------------------------------------------------- |
|---|
| 492 | |
|---|
| 493 | do ilev = 1, nlay |
|---|
| 494 | do ig=1,ngrid |
|---|
| 495 | pdudif(ig,ilev)=(zu(ig,ilev)-(pu(ig,ilev)))/ptimestep |
|---|
| 496 | pdvdif(ig,ilev)=(zv(ig,ilev)-(pv(ig,ilev)))/ptimestep |
|---|
| 497 | pdtdif(ig,ilev)=( zt(ig,ilev)- pt(ig,ilev))/ptimestep-pdtfi(ig,ilev) |
|---|
| 498 | enddo |
|---|
| 499 | enddo |
|---|
| 500 | |
|---|
| 501 | DO ig=1,ngrid ! computing sensible heat flux (atm => surface) |
|---|
| 502 | sensibFlux(ig)=cpp*zfluxt(ig,1)/ptimestep*(zt(ig,1)*zovExner(ig,1)-ztsrf(ig)) |
|---|
| 503 | ENDDO |
|---|
| 504 | |
|---|
| 505 | if (tracer) then |
|---|
| 506 | do iq = 1, nq |
|---|
| 507 | do ilev = 1, nlay |
|---|
| 508 | do ig=1,ngrid |
|---|
| 509 | pdqdif(ig,ilev,iq)=(zq(ig,ilev,iq)-(pq(ig,ilev,iq)+pdqfi(ig,ilev,iq)*ptimestep))/ptimestep |
|---|
| 510 | enddo |
|---|
| 511 | enddo |
|---|
| 512 | enddo |
|---|
| 513 | endif |
|---|
| 514 | end subroutine turbdiff |
|---|
| 515 | |
|---|
| 516 | end module turbdiff_mod |
|---|