[dced13c] | 1 | !********************************************************************** |
---|
| 2 | ! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010 * |
---|
| 3 | ! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa, * |
---|
| 4 | ! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann * |
---|
| 5 | ! * |
---|
| 6 | ! This file is part of FLEXPART. * |
---|
| 7 | ! * |
---|
| 8 | ! FLEXPART is free software: you can redistribute it and/or modify * |
---|
| 9 | ! it under the terms of the GNU General Public License as published by* |
---|
| 10 | ! the Free Software Foundation, either version 3 of the License, or * |
---|
| 11 | ! (at your option) any later version. * |
---|
| 12 | ! * |
---|
| 13 | ! FLEXPART is distributed in the hope that it will be useful, * |
---|
| 14 | ! but WITHOUT ANY WARRANTY; without even the implied warranty of * |
---|
| 15 | ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * |
---|
| 16 | ! GNU General Public License for more details. * |
---|
| 17 | ! * |
---|
| 18 | ! You should have received a copy of the GNU General Public License * |
---|
| 19 | ! along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. * |
---|
| 20 | !********************************************************************** |
---|
| 21 | |
---|
| 22 | subroutine advance_rec(itime,nrelpoint,ldt,up,vp,wp, & |
---|
| 23 | usigold,vsigold,wsigold,nstop,xt,yt,zt,prob,icbt) |
---|
| 24 | ! i i i/oi/oi/o |
---|
| 25 | ! i/o i/o i/o o i/oi/oi/o i/o i/o |
---|
| 26 | !***************************************************************************** |
---|
| 27 | ! * |
---|
| 28 | ! Calculation of turbulent particle trajectories utilizing a * |
---|
| 29 | ! zero-acceleration scheme, which is corrected by a numerically more * |
---|
| 30 | ! accurate Petterssen scheme whenever possible. * |
---|
| 31 | ! * |
---|
| 32 | ! Particle positions are read in, incremented, and returned to the calling * |
---|
| 33 | ! program. * |
---|
| 34 | ! * |
---|
| 35 | ! In different regions of the atmosphere (PBL vs. free troposphere), * |
---|
| 36 | ! different parameters are needed for advection, parameterizing turbulent * |
---|
| 37 | ! velocities, etc. For efficiency, different interpolation routines have * |
---|
| 38 | ! been written for these different cases, with the disadvantage that there * |
---|
| 39 | ! exist several routines doing almost the same. They all share the * |
---|
| 40 | ! included file 'interpol_mod'. The following * |
---|
| 41 | ! interpolation routines are used: * |
---|
| 42 | ! * |
---|
| 43 | ! interpol_all(_nests) interpolates everything (called inside the PBL) * |
---|
| 44 | ! interpol_misslev(_nests) if a particle moves vertically in the PBL, * |
---|
| 45 | ! additional parameters are interpolated if it * |
---|
| 46 | ! crosses a model level * |
---|
| 47 | ! interpol_wind(_nests) interpolates the wind and determines the * |
---|
| 48 | ! standard deviation of the wind (called outside * |
---|
| 49 | ! PBL) also interpolates potential vorticity * |
---|
| 50 | ! interpol_wind_short(_nests) only interpolates the wind (needed for the * |
---|
| 51 | ! Petterssen scheme) * |
---|
| 52 | ! interpol_vdep(_nests) interpolates deposition velocities * |
---|
| 53 | ! * |
---|
| 54 | ! * |
---|
| 55 | ! Author: A. Stohl * |
---|
| 56 | ! * |
---|
| 57 | ! 16 December 1997 * |
---|
| 58 | ! * |
---|
| 59 | ! Changes: * |
---|
| 60 | ! * |
---|
| 61 | ! 8 April 2000: Deep convection parameterization * |
---|
| 62 | ! * |
---|
| 63 | ! May 2002: Petterssen scheme introduced * |
---|
| 64 | ! * |
---|
| 65 | !***************************************************************************** |
---|
| 66 | ! * |
---|
| 67 | ! Variables: * |
---|
| 68 | ! icbt 1 if particle not transferred to forbidden state, * |
---|
| 69 | ! else -1 * |
---|
| 70 | ! dawsave accumulated displacement in along-wind direction * |
---|
| 71 | ! dcwsave accumulated displacement in cross-wind direction * |
---|
| 72 | ! dxsave accumulated displacement in longitude * |
---|
| 73 | ! dysave accumulated displacement in latitude * |
---|
| 74 | ! h [m] Mixing height * |
---|
| 75 | ! lwindinterv [s] time interval between two wind fields * |
---|
| 76 | ! itime [s] time at which this subroutine is entered * |
---|
| 77 | ! itimec [s] actual time, which is incremented in this subroutine * |
---|
| 78 | ! href [m] height for which dry deposition velocity is calculated * |
---|
| 79 | ! ladvance [s] Total integration time period * |
---|
| 80 | ! ldirect 1 forward, -1 backward * |
---|
| 81 | ! ldt [s] Time step for the next integration * |
---|
| 82 | ! lsynctime [s] Synchronisation interval of FLEXPART * |
---|
| 83 | ! ngrid index which grid is to be used * |
---|
| 84 | ! nrand index for a variable to be picked from rannumb * |
---|
| 85 | ! nstop if > 1 particle has left domain and must be stopped * |
---|
| 86 | ! prob probability of absorption due to dry deposition * |
---|
| 87 | ! rannumb(maxrand) normally distributed random variables * |
---|
| 88 | ! rhoa air density * |
---|
| 89 | ! rhograd vertical gradient of the air density * |
---|
| 90 | ! up,vp,wp random velocities due to turbulence (along wind, cross * |
---|
| 91 | ! wind, vertical wind * |
---|
| 92 | ! usig,vsig,wsig mesoscale wind fluctuations * |
---|
| 93 | ! usigold,vsigold,wsigold like usig, etc., but for the last time step * |
---|
| 94 | ! vdepo Deposition velocities for all species * |
---|
| 95 | ! xt,yt,zt Particle position * |
---|
| 96 | ! * |
---|
| 97 | !***************************************************************************** |
---|
| 98 | |
---|
| 99 | use point_mod |
---|
| 100 | use par_mod |
---|
| 101 | use com_mod |
---|
| 102 | use interpol_mod |
---|
| 103 | use hanna_mod |
---|
| 104 | use cmapf_mod |
---|
| 105 | use random_mod, only: ran3 |
---|
| 106 | |
---|
| 107 | implicit none |
---|
| 108 | |
---|
| 109 | real(kind=dp) :: xt,yt |
---|
| 110 | real :: zt,xts,yts,weight |
---|
| 111 | integer :: itime,itimec,nstop,ldt,i,j,k,nrand,loop,memindnext,mind |
---|
| 112 | integer :: ngr,nix,njy,ks,nsp,nrelpoint |
---|
| 113 | real :: dz,dz1,dz2,xlon,ylat,xpol,ypol,gridsize |
---|
| 114 | real :: ru,rv,rw,dt,ux,vy,cosfact,xtn,ytn,tropop |
---|
| 115 | real :: prob(maxspec),up,vp,wp,dxsave,dysave,dawsave |
---|
| 116 | real :: dcwsave |
---|
| 117 | real :: usigold,vsigold,wsigold,r,rs |
---|
| 118 | real :: uold,vold,wold,vdepo(maxspec) |
---|
| 119 | !real uprof(nzmax),vprof(nzmax),wprof(nzmax) |
---|
| 120 | !real usigprof(nzmax),vsigprof(nzmax),wsigprof(nzmax) |
---|
| 121 | !real rhoprof(nzmax),rhogradprof(nzmax) |
---|
| 122 | real :: rhoa,rhograd,delz,dtf,rhoaux,dtftlw,uxscale,wpscale |
---|
| 123 | integer(kind=2) :: icbt |
---|
| 124 | real,parameter :: eps=nxmax/3.e5,eps2=1.e-9 |
---|
| 125 | real :: ptot_lhh,Q_lhh,phi_lhh,ath,bth !modified by mc |
---|
| 126 | real :: old_wp_buf,dcas,dcas1,del_test !added by mc |
---|
| 127 | integer :: i_well,jj,flagrein !test well mixed: modified by mc |
---|
| 128 | |
---|
| 129 | |
---|
| 130 | integer :: idummy = -7 |
---|
| 131 | real :: settling = 0. |
---|
| 132 | |
---|
| 133 | |
---|
| 134 | nstop=0 |
---|
| 135 | do i=1,nmixz |
---|
| 136 | indzindicator(i)=.true. |
---|
| 137 | end do |
---|
| 138 | |
---|
| 139 | |
---|
| 140 | if (DRYDEP) then ! reset probability for deposition |
---|
| 141 | do ks=1,nspec |
---|
| 142 | depoindicator(ks)=.true. |
---|
| 143 | prob(ks)=0. |
---|
| 144 | end do |
---|
| 145 | endif |
---|
| 146 | |
---|
| 147 | dxsave=0. ! reset position displacements |
---|
| 148 | dysave=0. ! due to mean wind |
---|
| 149 | dawsave=0. ! and turbulent wind |
---|
| 150 | dcwsave=0. |
---|
| 151 | |
---|
| 152 | itimec=itime |
---|
| 153 | |
---|
| 154 | nrand=int(ran3(idummy)*real(maxrand-1))+1 |
---|
| 155 | |
---|
| 156 | |
---|
| 157 | ! Determine whether lat/long grid or polarstereographic projection |
---|
| 158 | ! is to be used |
---|
| 159 | ! Furthermore, determine which nesting level to be used |
---|
| 160 | !***************************************************************** |
---|
| 161 | |
---|
| 162 | if (nglobal.and.(yt.gt.switchnorthg)) then |
---|
| 163 | ngrid=-1 |
---|
| 164 | else if (sglobal.and.(yt.lt.switchsouthg)) then |
---|
| 165 | ngrid=-2 |
---|
| 166 | else |
---|
| 167 | ngrid=0 |
---|
| 168 | do j=numbnests,1,-1 |
---|
| 169 | if ((xt.gt.xln(j)+eps).and.(xt.lt.xrn(j)-eps).and. & |
---|
| 170 | (yt.gt.yln(j)+eps).and.(yt.lt.yrn(j)-eps)) then |
---|
| 171 | ngrid=j |
---|
| 172 | goto 23 |
---|
| 173 | endif |
---|
| 174 | end do |
---|
| 175 | 23 continue |
---|
| 176 | endif |
---|
| 177 | |
---|
| 178 | |
---|
| 179 | !*************************** |
---|
| 180 | ! Interpolate necessary data |
---|
| 181 | !*************************** |
---|
| 182 | |
---|
| 183 | if (abs(itime-memtime(1)).lt.abs(itime-memtime(2))) then |
---|
| 184 | memindnext=1 |
---|
| 185 | else |
---|
| 186 | memindnext=2 |
---|
| 187 | endif |
---|
| 188 | |
---|
| 189 | ! Determine nested grid coordinates |
---|
| 190 | !********************************** |
---|
| 191 | |
---|
| 192 | if (ngrid.gt.0) then |
---|
| 193 | xtn=(xt-xln(ngrid))*xresoln(ngrid) |
---|
| 194 | ytn=(yt-yln(ngrid))*yresoln(ngrid) |
---|
| 195 | ix=int(xtn) |
---|
| 196 | jy=int(ytn) |
---|
| 197 | nix=nint(xtn) |
---|
| 198 | njy=nint(ytn) |
---|
| 199 | else |
---|
| 200 | ix=int(xt) |
---|
| 201 | jy=int(yt) |
---|
| 202 | nix=nint(xt) |
---|
| 203 | njy=nint(yt) |
---|
| 204 | endif |
---|
| 205 | ixp=ix+1 |
---|
| 206 | jyp=jy+1 |
---|
| 207 | |
---|
| 208 | |
---|
| 209 | ! Compute maximum mixing height around particle position |
---|
| 210 | !******************************************************* |
---|
| 211 | |
---|
| 212 | h=0. |
---|
| 213 | if (ngrid.le.0) then |
---|
| 214 | do k=1,2 |
---|
| 215 | mind=memind(k) ! eso: compatibility with 3-field version |
---|
| 216 | do j=jy,jyp |
---|
| 217 | do i=ix,ixp |
---|
| 218 | if (hmix(i,j,1,mind).gt.h) h=hmix(i,j,1,mind) |
---|
| 219 | end do |
---|
| 220 | end do |
---|
| 221 | end do |
---|
| 222 | tropop=tropopause(nix,njy,1,1) |
---|
| 223 | else |
---|
| 224 | do k=1,2 |
---|
| 225 | mind=memind(k) |
---|
| 226 | do j=jy,jyp |
---|
| 227 | do i=ix,ixp |
---|
| 228 | if (hmixn(i,j,1,mind,ngrid).gt.h) h=hmixn(i,j,1,mind,ngrid) |
---|
| 229 | end do |
---|
| 230 | end do |
---|
| 231 | end do |
---|
| 232 | tropop=tropopausen(nix,njy,1,1,ngrid) |
---|
| 233 | endif |
---|
| 234 | |
---|
| 235 | zeta=zt/h |
---|
| 236 | |
---|
| 237 | |
---|
| 238 | |
---|
| 239 | !************************************************************* |
---|
| 240 | ! If particle is in the PBL, interpolate once and then make a |
---|
| 241 | ! time loop until end of interval is reached |
---|
| 242 | !************************************************************* |
---|
| 243 | |
---|
| 244 | if (zeta.le.1.) then |
---|
| 245 | |
---|
| 246 | ! BEGIN TIME LOOP |
---|
| 247 | !================ |
---|
| 248 | |
---|
| 249 | loop=0 |
---|
| 250 | 100 loop=loop+1 |
---|
| 251 | if (method.eq.1) then |
---|
| 252 | ldt=min(ldt,abs(lsynctime-itimec+itime)) |
---|
| 253 | itimec=itimec+ldt*ldirect |
---|
| 254 | else |
---|
| 255 | ldt=abs(lsynctime) |
---|
| 256 | itimec=itime+lsynctime |
---|
| 257 | endif |
---|
| 258 | dt=real(ldt) |
---|
| 259 | |
---|
| 260 | zeta=zt/h |
---|
| 261 | |
---|
| 262 | |
---|
| 263 | if (loop.eq.1) then |
---|
| 264 | if (ngrid.le.0) then |
---|
| 265 | xts=real(xt) |
---|
| 266 | yts=real(yt) |
---|
| 267 | call interpol_all(itime,xts,yts,zt) |
---|
| 268 | else |
---|
| 269 | call interpol_all_nests(itime,xtn,ytn,zt) |
---|
| 270 | endif |
---|
| 271 | |
---|
| 272 | else |
---|
| 273 | |
---|
| 274 | |
---|
| 275 | ! Determine the level below the current position for u,v,rho |
---|
| 276 | !*********************************************************** |
---|
| 277 | |
---|
| 278 | do i=2,nz |
---|
| 279 | if (height(i).gt.zt) then |
---|
| 280 | indz=i-1 |
---|
| 281 | indzp=i |
---|
| 282 | goto 6 |
---|
| 283 | endif |
---|
| 284 | end do |
---|
| 285 | 6 continue |
---|
| 286 | |
---|
| 287 | ! If one of the levels necessary is not yet available, |
---|
| 288 | ! calculate it |
---|
| 289 | !***************************************************** |
---|
| 290 | |
---|
| 291 | do i=indz,indzp |
---|
| 292 | if (indzindicator(i)) then |
---|
| 293 | if (ngrid.le.0) then |
---|
| 294 | call interpol_misslev(i) |
---|
| 295 | else |
---|
| 296 | call interpol_misslev_nests(i) |
---|
| 297 | endif |
---|
| 298 | endif |
---|
| 299 | end do |
---|
| 300 | endif |
---|
| 301 | |
---|
| 302 | |
---|
| 303 | ! Vertical interpolation of u,v,w,rho and drhodz |
---|
| 304 | !*********************************************** |
---|
| 305 | |
---|
| 306 | ! Vertical distance to the level below and above current position |
---|
| 307 | ! both in terms of (u,v) and (w) fields |
---|
| 308 | !**************************************************************** |
---|
| 309 | |
---|
| 310 | dz=1./(height(indzp)-height(indz)) |
---|
| 311 | dz1=(zt-height(indz))*dz |
---|
| 312 | dz2=(height(indzp)-zt)*dz |
---|
| 313 | |
---|
| 314 | u=dz1*uprof(indzp)+dz2*uprof(indz) |
---|
| 315 | v=dz1*vprof(indzp)+dz2*vprof(indz) |
---|
| 316 | w=dz1*wprof(indzp)+dz2*wprof(indz) |
---|
| 317 | rhoa=dz1*rhoprof(indzp)+dz2*rhoprof(indz) |
---|
| 318 | rhograd=dz1*rhogradprof(indzp)+dz2*rhogradprof(indz) |
---|
| 319 | |
---|
| 320 | |
---|
| 321 | ! Compute the turbulent disturbances |
---|
| 322 | ! Determine the sigmas and the timescales |
---|
| 323 | !**************************************** |
---|
| 324 | |
---|
| 325 | if (turbswitch) then |
---|
| 326 | call hanna(zt) |
---|
| 327 | else |
---|
| 328 | call hanna1(zt) |
---|
| 329 | endif |
---|
| 330 | |
---|
| 331 | |
---|
| 332 | !***************************************** |
---|
| 333 | ! Determine the new diffusivity velocities |
---|
| 334 | !***************************************** |
---|
| 335 | |
---|
| 336 | ! Horizontal components |
---|
| 337 | !********************** |
---|
| 338 | |
---|
| 339 | if (nrand+1.gt.maxrand) nrand=1 |
---|
| 340 | if (dt/tlu.lt..5) then |
---|
| 341 | up=(1.-dt/tlu)*up+rannumb(nrand)*sigu*sqrt(2.*dt/tlu) |
---|
| 342 | else |
---|
| 343 | ru=exp(-dt/tlu) |
---|
| 344 | up=ru*up+rannumb(nrand)*sigu*sqrt(1.-ru**2) |
---|
| 345 | endif |
---|
| 346 | if (dt/tlv.lt..5) then |
---|
| 347 | vp=(1.-dt/tlv)*vp+rannumb(nrand+1)*sigv*sqrt(2.*dt/tlv) |
---|
| 348 | else |
---|
| 349 | rv=exp(-dt/tlv) |
---|
| 350 | vp=rv*vp+rannumb(nrand+1)*sigv*sqrt(1.-rv**2) |
---|
| 351 | endif |
---|
| 352 | nrand=nrand+2 |
---|
| 353 | |
---|
| 354 | |
---|
| 355 | if (nrand+ifine.gt.maxrand) nrand=1 |
---|
| 356 | rhoaux=rhograd/rhoa |
---|
| 357 | dtf=dt*fine |
---|
| 358 | |
---|
| 359 | dtftlw=dtf/tlw |
---|
| 360 | |
---|
| 361 | ! Loop over ifine short time steps for vertical component |
---|
| 362 | !******************************************************** |
---|
| 363 | |
---|
| 364 | do i=1,ifine |
---|
| 365 | |
---|
| 366 | ! Determine the drift velocity and density correction velocity |
---|
| 367 | !************************************************************* |
---|
| 368 | |
---|
| 369 | if (turbswitch) then |
---|
| 370 | if (dtftlw.lt..5) then |
---|
| 371 | !************************************************************* |
---|
| 372 | !************** CBL options added by mc see routine cblf90 *** |
---|
| 373 | if (cblflag.eq.1) then !modified by mc |
---|
| 374 | if (-h/ol.gt.5) then !modified by mc |
---|
| 375 | !if (ol.lt.0.) then !modified by mc |
---|
| 376 | !if (ol.gt.0.) then !modified by mc : for test |
---|
| 377 | !print *,zt,wp,ath,bth,tlw,dtf,'prima' |
---|
| 378 | flagrein=0 |
---|
| 379 | nrand=nrand+1 |
---|
| 380 | old_wp_buf=wp |
---|
| 381 | call cbl(wp,zt,ust,wst,h,rhoa,rhograd,sigw,dsigwdz,tlw,ptot_lhh,Q_lhh,phi_lhh,ath,bth,ol,flagrein) !inside the routine for inverse time |
---|
| 382 | wp=(wp+ath*dtf+bth*rannumb(nrand)*sqrt(dtf))*real(icbt) |
---|
| 383 | ! wp=(wp+ath*dtf+bth*gasdev2(mydum)*sqrt(dtf))*real(icbt) |
---|
| 384 | delz=wp*dtf |
---|
| 385 | if (flagrein.eq.1) then |
---|
| 386 | call re_initialize_particle(zt,ust,wst,h,sigw,old_wp_buf,nrand,ol) |
---|
| 387 | wp=old_wp_buf |
---|
| 388 | delz=wp*dtf |
---|
| 389 | nan_count=nan_count+1 |
---|
| 390 | end if |
---|
| 391 | !print *,zt,wp,ath,bth,tlw,dtf,rannumb(nrand+i),icbt |
---|
| 392 | !pause |
---|
| 393 | else |
---|
| 394 | nrand=nrand+1 |
---|
| 395 | old_wp_buf=wp |
---|
| 396 | ath=-wp/tlw+sigw*dsigwdz+wp*wp/sigw*dsigwdz+sigw*sigw/rhoa*rhograd !1-note for inverse time should be -wp/tlw*ldirect+... calculated for wp=-wp |
---|
| 397 | !2-but since ldirect =-1 for inverse time and this must be calculated for (-wp) and |
---|
| 398 | !3-the gaussian pdf is symmetric (i.e. pdf(w)=pdf(-w) ldirect can be discarded |
---|
| 399 | bth=sigw*rannumb(nrand)*sqrt(2.*dtftlw) |
---|
| 400 | wp=(wp+ath*dtf+bth)*real(icbt) |
---|
| 401 | delz=wp*dtf |
---|
| 402 | del_test=(1.-wp)/wp !catch infinity value |
---|
| 403 | if (isnan(wp).or.isnan(del_test)) then |
---|
| 404 | nrand=nrand+1 |
---|
| 405 | wp=sigw*rannumb(nrand) |
---|
| 406 | delz=wp*dtf |
---|
| 407 | nan_count2=nan_count2+1 |
---|
| 408 | !print *,'NaN coutner equal to:', nan_count,'reduce ifine if this number became a non-negligible fraction of the particle number' |
---|
| 409 | end if |
---|
| 410 | end if |
---|
| 411 | !******************** END CBL option ******************************* |
---|
| 412 | !******************************************************************* |
---|
| 413 | else |
---|
| 414 | wp=((1.-dtftlw)*wp+rannumb(nrand+i)*sqrt(2.*dtftlw) & |
---|
| 415 | +dtf*(dsigwdz+rhoaux*sigw))*real(icbt) |
---|
| 416 | delz=wp*sigw*dtf |
---|
| 417 | end if |
---|
| 418 | else |
---|
| 419 | rw=exp(-dtftlw) |
---|
| 420 | wp=(rw*wp+rannumb(nrand+i)*sqrt(1.-rw**2) & |
---|
| 421 | +tlw*(1.-rw)*(dsigwdz+rhoaux*sigw))*real(icbt) |
---|
| 422 | delz=wp*sigw*dtf |
---|
| 423 | endif |
---|
| 424 | |
---|
| 425 | else |
---|
| 426 | rw=exp(-dtftlw) |
---|
| 427 | wp=(rw*wp+rannumb(nrand+i)*sqrt(1.-rw**2)*sigw & |
---|
| 428 | +tlw*(1.-rw)*(dsigw2dz+rhoaux*sigw**2))*real(icbt) |
---|
| 429 | delz=wp*dtf |
---|
| 430 | endif |
---|
| 431 | |
---|
| 432 | !**************************************************** |
---|
| 433 | ! Compute turbulent vertical displacement of particle |
---|
| 434 | !**************************************************** |
---|
| 435 | |
---|
| 436 | if (abs(delz).gt.h) delz=mod(delz,h) |
---|
| 437 | |
---|
| 438 | ! Determine if particle transfers to a "forbidden state" below the ground |
---|
| 439 | ! or above the mixing height |
---|
| 440 | !************************************************************************ |
---|
| 441 | |
---|
| 442 | if (delz.lt.-zt) then ! reflection at ground |
---|
| 443 | icbt=-1 |
---|
| 444 | zt=-zt-delz |
---|
| 445 | else if (delz.gt.(h-zt)) then ! reflection at h |
---|
| 446 | icbt=-1 |
---|
| 447 | zt=-zt-delz+2.*h |
---|
| 448 | else ! no reflection |
---|
| 449 | icbt=1 |
---|
| 450 | zt=zt+delz |
---|
| 451 | endif |
---|
| 452 | |
---|
| 453 | if (i.ne.ifine) then |
---|
| 454 | zeta=zt/h |
---|
| 455 | call hanna_short(zt) |
---|
| 456 | endif |
---|
| 457 | |
---|
| 458 | end do |
---|
| 459 | if (cblflag.ne.1) nrand=nrand+i |
---|
| 460 | |
---|
| 461 | ! Determine time step for next integration |
---|
| 462 | !***************************************** |
---|
| 463 | |
---|
| 464 | if (turbswitch) then |
---|
| 465 | ldt=int(min(tlw,h/max(2.*abs(wp*sigw),1.e-5), & |
---|
| 466 | 0.5/abs(dsigwdz))*ctl) |
---|
| 467 | else |
---|
| 468 | ldt=int(min(tlw,h/max(2.*abs(wp),1.e-5))*ctl) |
---|
| 469 | endif |
---|
| 470 | ldt=max(ldt,mintime) |
---|
| 471 | |
---|
| 472 | ! Determine probability of deposition |
---|
| 473 | !************************************ |
---|
| 474 | |
---|
| 475 | if ((DRYDEP).and.(zt.lt.2.*href)) then |
---|
| 476 | do ks=1,nspec |
---|
| 477 | if (DRYDEPSPEC(ks)) then |
---|
| 478 | if (depoindicator(ks)) then |
---|
| 479 | if (ngrid.le.0) then |
---|
| 480 | call interpol_vdep(ks,vdepo(ks)) |
---|
| 481 | else |
---|
| 482 | call interpol_vdep_nests(ks,vdepo(ks)) |
---|
| 483 | endif |
---|
| 484 | endif |
---|
| 485 | ! correction by Petra Seibert, 10 April 2001 |
---|
| 486 | ! this formulation means that prob(n) = 1 - f(0)*...*f(n) |
---|
| 487 | ! where f(n) is the exponential term |
---|
| 488 | prob(ks)=1.+(prob(ks)-1.)* & |
---|
| 489 | exp(-vdepo(ks)*abs(dt)/(2.*href)) |
---|
| 490 | endif |
---|
| 491 | end do |
---|
| 492 | endif |
---|
| 493 | |
---|
| 494 | endif |
---|
| 495 | |
---|
| 496 | end subroutine advance_rec |
---|
| 497 | |
---|