[e200b7a] | 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 wetdepo(itime,ltsample,loutnext) |
---|
[4fbe7a5] | 23 | ! i i i |
---|
[e200b7a] | 24 | !***************************************************************************** |
---|
| 25 | ! * |
---|
| 26 | ! Calculation of wet deposition using the concept of scavenging coefficients.* |
---|
| 27 | ! For lack of detailed information, washout and rainout are jointly treated. * |
---|
| 28 | ! It is assumed that precipitation does not occur uniformly within the whole * |
---|
| 29 | ! grid cell, but that only a fraction of the grid cell experiences rainfall. * |
---|
| 30 | ! This fraction is parameterized from total cloud cover and rates of large * |
---|
| 31 | ! scale and convective precipitation. * |
---|
| 32 | ! * |
---|
| 33 | ! Author: A. Stohl * |
---|
| 34 | ! * |
---|
| 35 | ! 1 December 1996 * |
---|
| 36 | ! * |
---|
| 37 | ! Correction by Petra Seibert, Sept 2002: * |
---|
| 38 | ! use centred precipitation data for integration * |
---|
| 39 | ! Code may not be correct for decay of deposition! * |
---|
| 40 | ! * |
---|
| 41 | !***************************************************************************** |
---|
| 42 | ! * |
---|
| 43 | ! Variables: * |
---|
| 44 | ! cc [0-1] total cloud cover * |
---|
| 45 | ! convp [mm/h] convective precipitation rate * |
---|
| 46 | ! grfraction [0-1] fraction of grid, for which precipitation occurs * |
---|
| 47 | ! ix,jy indices of output grid cell for each particle * |
---|
| 48 | ! itime [s] actual simulation time [s] * |
---|
| 49 | ! jpart particle index * |
---|
| 50 | ! ldeltat [s] interval since radioactive decay was computed * |
---|
| 51 | ! lfr, cfr area fraction covered by precipitation for large scale * |
---|
| 52 | ! and convective precipitation (dependent on prec. rate) * |
---|
| 53 | ! loutnext [s] time for which gridded deposition is next output * |
---|
| 54 | ! loutstep [s] interval at which gridded deposition is output * |
---|
| 55 | ! lsp [mm/h] large scale precipitation rate * |
---|
| 56 | ! ltsample [s] interval over which mass is deposited * |
---|
| 57 | ! prec [mm/h] precipitation rate in subgrid, where precipitation occurs* |
---|
| 58 | ! wetdeposit mass that is wet deposited * |
---|
| 59 | ! wetgrid accumulated deposited mass on output grid * |
---|
| 60 | ! wetscav scavenging coefficient * |
---|
| 61 | ! * |
---|
| 62 | ! Constants: * |
---|
| 63 | ! * |
---|
| 64 | !***************************************************************************** |
---|
| 65 | |
---|
| 66 | use point_mod |
---|
| 67 | use par_mod |
---|
| 68 | use com_mod |
---|
| 69 | |
---|
| 70 | implicit none |
---|
| 71 | |
---|
| 72 | integer :: jpart,itime,ltsample,loutnext,ldeltat,i,j,ix,jy |
---|
[4fbe7a5] | 73 | integer :: ngrid,itage,nage,hz,il,interp_time, n, clouds_v |
---|
| 74 | integer :: ks, kp |
---|
[8a65cb0] | 75 | integer :: n1,n2, icbot,ictop, indcloud !TEST |
---|
[e200b7a] | 76 | real :: S_i, act_temp, cl, cle ! in cloud scavenging |
---|
| 77 | real :: clouds_h ! cloud height for the specific grid point |
---|
[8a65cb0] | 78 | real :: xtn,ytn,lsp,convp,cc,grfraction(3),prec(3),wetscav,totprec |
---|
[e200b7a] | 79 | real :: wetdeposit(maxspec),restmass |
---|
| 80 | real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled |
---|
[5f9d14a] | 81 | !save lfr,cfr |
---|
[e200b7a] | 82 | |
---|
[5f9d14a] | 83 | real, parameter :: lfr(5) = (/ 0.5,0.65,0.8,0.9,0.95/) |
---|
| 84 | real, parameter :: cfr(5) = (/ 0.4,0.55,0.7,0.8,0.9 /) |
---|
[e200b7a] | 85 | |
---|
[8a65cb0] | 86 | !ZHG aerosol below-cloud scavenging removal polynomial constants for rain and snow |
---|
[5f9d14a] | 87 | real, parameter :: bclr(6) = (/274.35758, 332839.59273, 226656.57259, 58005.91340, 6588.38582, 0.244984/) !rain (Laakso et al 2003) |
---|
| 88 | real, parameter :: bcls(6) = (/22.7, 0.0, 0.0, 1321.0, 381.0, 0.0/) !now (Kyro et al 2009) |
---|
[8a65cb0] | 89 | real :: frac_act, liq_frac, dquer_m |
---|
| 90 | |
---|
| 91 | integer :: blc_count, inc_count |
---|
| 92 | |
---|
| 93 | |
---|
| 94 | |
---|
[e200b7a] | 95 | ! Compute interval since radioactive decay of deposited mass was computed |
---|
| 96 | !************************************************************************ |
---|
| 97 | |
---|
| 98 | if (itime.le.loutnext) then |
---|
| 99 | ldeltat=itime-(loutnext-loutstep) |
---|
| 100 | else ! first half of next interval |
---|
| 101 | ldeltat=itime-loutnext |
---|
| 102 | endif |
---|
| 103 | |
---|
| 104 | ! Loop over all particles |
---|
| 105 | !************************ |
---|
| 106 | |
---|
[8a65cb0] | 107 | |
---|
| 108 | blc_count=0 |
---|
| 109 | inc_count=0 |
---|
| 110 | |
---|
[e200b7a] | 111 | do jpart=1,numpart |
---|
[4fbe7a5] | 112 | |
---|
[e200b7a] | 113 | if (itra1(jpart).eq.-999999999) goto 20 |
---|
| 114 | if(ldirect.eq.1)then |
---|
| 115 | if (itra1(jpart).gt.itime) goto 20 |
---|
| 116 | else |
---|
| 117 | if (itra1(jpart).lt.itime) goto 20 |
---|
| 118 | endif |
---|
[8a65cb0] | 119 | |
---|
[e200b7a] | 120 | ! Determine age class of the particle |
---|
| 121 | itage=abs(itra1(jpart)-itramem(jpart)) |
---|
| 122 | do nage=1,nageclass |
---|
| 123 | if (itage.lt.lage(nage)) goto 33 |
---|
| 124 | end do |
---|
[4fbe7a5] | 125 | 33 continue |
---|
[e200b7a] | 126 | |
---|
| 127 | |
---|
| 128 | ! Determine which nesting level to be used |
---|
| 129 | !***************************************** |
---|
| 130 | |
---|
| 131 | ngrid=0 |
---|
| 132 | do j=numbnests,1,-1 |
---|
| 133 | if ((xtra1(jpart).gt.xln(j)).and.(xtra1(jpart).lt.xrn(j)).and. & |
---|
[4fbe7a5] | 134 | (ytra1(jpart).gt.yln(j)).and.(ytra1(jpart).lt.yrn(j))) then |
---|
[e200b7a] | 135 | ngrid=j |
---|
| 136 | goto 23 |
---|
| 137 | endif |
---|
| 138 | end do |
---|
[4fbe7a5] | 139 | 23 continue |
---|
[e200b7a] | 140 | |
---|
| 141 | |
---|
| 142 | ! Determine nested grid coordinates |
---|
| 143 | !********************************** |
---|
| 144 | |
---|
| 145 | if (ngrid.gt.0) then |
---|
| 146 | xtn=(xtra1(jpart)-xln(ngrid))*xresoln(ngrid) |
---|
| 147 | ytn=(ytra1(jpart)-yln(ngrid))*yresoln(ngrid) |
---|
| 148 | ix=int(xtn) |
---|
| 149 | jy=int(ytn) |
---|
| 150 | else |
---|
| 151 | ix=int(xtra1(jpart)) |
---|
| 152 | jy=int(ytra1(jpart)) |
---|
| 153 | endif |
---|
| 154 | |
---|
| 155 | |
---|
| 156 | ! Interpolate large scale precipitation, convective precipitation and |
---|
| 157 | ! total cloud cover |
---|
| 158 | ! Note that interpolated time refers to itime-0.5*ltsample [PS] |
---|
| 159 | !******************************************************************** |
---|
| 160 | interp_time=nint(itime-0.5*ltsample) |
---|
[5f9d14a] | 161 | |
---|
[4fbe7a5] | 162 | if (ngrid.eq.0) then |
---|
| 163 | call interpol_rain(lsprec,convprec,tcc,nxmax,nymax, & |
---|
| 164 | 1,nx,ny,memind,real(xtra1(jpart)),real(ytra1(jpart)),1, & |
---|
| 165 | memtime(1),memtime(2),interp_time,lsp,convp,cc) |
---|
| 166 | else |
---|
| 167 | call interpol_rain_nests(lsprecn,convprecn,tccn, & |
---|
| 168 | nxmaxn,nymaxn,1,maxnests,ngrid,nxn,nyn,memind,xtn,ytn,1, & |
---|
| 169 | memtime(1),memtime(2),interp_time,lsp,convp,cc) |
---|
| 170 | endif |
---|
[e200b7a] | 171 | |
---|
[8a65cb0] | 172 | ! If total precipitation is less than 0.01 mm/h - no scavenging occurs |
---|
[5f9d14a] | 173 | if ((lsp.lt.0.01).and.(convp.lt.0.01)) goto 20 |
---|
[8a65cb0] | 174 | |
---|
[e200b7a] | 175 | ! get the level were the actual particle is in |
---|
[4fbe7a5] | 176 | do il=2,nz |
---|
| 177 | if (height(il).gt.ztra1(jpart)) then |
---|
| 178 | hz=il-1 |
---|
| 179 | goto 26 |
---|
| 180 | endif |
---|
| 181 | end do |
---|
| 182 | 26 continue |
---|
[e200b7a] | 183 | |
---|
[4fbe7a5] | 184 | n=memind(2) |
---|
| 185 | if (abs(memtime(1)-interp_time).lt.abs(memtime(2)-interp_time)) & |
---|
| 186 | n=memind(1) |
---|
[e200b7a] | 187 | |
---|
[4fbe7a5] | 188 | if (ngrid.eq.0) then |
---|
| 189 | clouds_v=clouds(ix,jy,hz,n) |
---|
| 190 | clouds_h=cloudsh(ix,jy,n) |
---|
| 191 | else |
---|
| 192 | clouds_v=cloudsn(ix,jy,hz,n,ngrid) |
---|
| 193 | clouds_h=cloudsnh(ix,jy,n,ngrid) |
---|
| 194 | endif |
---|
[e200b7a] | 195 | |
---|
[8a65cb0] | 196 | ! if there is no precipitation or the particle is above the clouds no |
---|
| 197 | ! scavenging is done |
---|
| 198 | |
---|
| 199 | if (clouds_v.le.1) goto 20 |
---|
| 200 | |
---|
[e200b7a] | 201 | ! 1) Parameterization of the the area fraction of the grid cell where the |
---|
| 202 | ! precipitation occurs: the absolute limit is the total cloud cover, but |
---|
| 203 | ! for low precipitation rates, an even smaller fraction of the grid cell |
---|
| 204 | ! is used. Large scale precipitation occurs over larger areas than |
---|
| 205 | ! convective precipitation. |
---|
| 206 | !************************************************************************** |
---|
| 207 | |
---|
| 208 | if (lsp.gt.20.) then |
---|
| 209 | i=5 |
---|
| 210 | else if (lsp.gt.8.) then |
---|
| 211 | i=4 |
---|
| 212 | else if (lsp.gt.3.) then |
---|
| 213 | i=3 |
---|
| 214 | else if (lsp.gt.1.) then |
---|
| 215 | i=2 |
---|
| 216 | else |
---|
| 217 | i=1 |
---|
| 218 | endif |
---|
| 219 | |
---|
| 220 | if (convp.gt.20.) then |
---|
| 221 | j=5 |
---|
| 222 | else if (convp.gt.8.) then |
---|
| 223 | j=4 |
---|
| 224 | else if (convp.gt.3.) then |
---|
| 225 | j=3 |
---|
| 226 | else if (convp.gt.1.) then |
---|
| 227 | j=2 |
---|
| 228 | else |
---|
| 229 | j=1 |
---|
| 230 | endif |
---|
| 231 | |
---|
[8a65cb0] | 232 | |
---|
[5f9d14a] | 233 | !ZHG calculated for 1) both 2) lsp 3) convp |
---|
[8a65cb0] | 234 | grfraction(1)=max(0.05,cc*(lsp*lfr(i)+convp*cfr(j))/(lsp+convp)) |
---|
| 235 | grfraction(2)=max(0.05,cc*(lfr(i))) |
---|
| 236 | grfraction(3)=max(0.05,cc*(cfr(j))) |
---|
| 237 | |
---|
[e200b7a] | 238 | |
---|
| 239 | ! 2) Computation of precipitation rate in sub-grid cell |
---|
| 240 | !****************************************************** |
---|
[8a65cb0] | 241 | prec(1)=(lsp+convp)/grfraction(1) |
---|
| 242 | prec(2)=(lsp)/grfraction(2) |
---|
| 243 | prec(3)=(convp)/grfraction(3) |
---|
[e200b7a] | 244 | |
---|
| 245 | |
---|
| 246 | ! 3) Computation of scavenging coefficients for all species |
---|
| 247 | ! Computation of wet deposition |
---|
| 248 | !********************************************************** |
---|
| 249 | |
---|
[4fbe7a5] | 250 | do ks=1,nspec ! loop over species |
---|
| 251 | wetdeposit(ks)=0. |
---|
| 252 | wetscav=0. |
---|
| 253 | |
---|
[5f9d14a] | 254 | if (ngrid.gt.0) then |
---|
| 255 | act_temp=ttn(ix,jy,hz,n,ngrid) |
---|
| 256 | else |
---|
| 257 | act_temp=tt(ix,jy,hz,n) |
---|
| 258 | endif |
---|
[8a65cb0] | 259 | |
---|
[4fbe7a5] | 260 | |
---|
| 261 | !****i******************* |
---|
| 262 | ! BELOW CLOUD SCAVENGING |
---|
| 263 | !*********************** |
---|
[8a65cb0] | 264 | if (clouds_v.ge.4) then !below cloud |
---|
| 265 | |
---|
[0f20c31] | 266 | if (weta(ks).gt.0. .or. wetb(ks).gt.0.) then !if positive below-cloud parameters given in SPECIES file (either A or B) |
---|
[5f9d14a] | 267 | blc_count=blc_count+1 |
---|
[8a65cb0] | 268 | |
---|
[0f20c31] | 269 | |
---|
[8a65cb0] | 270 | !GAS |
---|
[5f9d14a] | 271 | if (dquer(ks) .le. 0.) then !is gas |
---|
| 272 | ! Gas scavenging coefficient based on Hertel et al 1995 using the below-cloud scavenging parameters A (=weta) and B (=wetb) from SPECIES file |
---|
[8a65cb0] | 273 | wetscav=weta(ks)*prec(1)**wetb(ks) |
---|
| 274 | |
---|
[5f9d14a] | 275 | !AEROSOL |
---|
| 276 | else !is particle |
---|
| 277 | !NIK 17.02.2015 |
---|
| 278 | ! For the calculation here particle size needs to be in meter and not um as dquer is changed to in readreleases |
---|
| 279 | dquer_m=dquer(ks)/1000000. !conversion from um to m |
---|
| 280 | |
---|
| 281 | !ZHG snow or rain removal is applied based on the temperature. |
---|
[0f20c31] | 282 | if (act_temp .ge. 273 .and. weta(ks).gt.0.) then !Rain |
---|
[8a65cb0] | 283 | |
---|
[5f9d14a] | 284 | !Particle RAIN scavenging coefficient based on Laakso et al 2003, the below-cloud scavenging (rain efficienty) parameter A (=weta) from SPECIES file |
---|
| 285 | wetscav= weta(ks)*10**(bclr(1)+ (bclr(2)*(log10(dquer_m))**(-4))+(bclr(3)*(log10(dquer_m))**(-3))+ (bclr(4)* & |
---|
| 286 | (log10(dquer_m))**(-2))+ (bclr(5)*(log10(dquer_m))**(-1))+ bclr(6)* (prec(1))**(0.5)) |
---|
[8a65cb0] | 287 | |
---|
[0f20c31] | 288 | elseif (act_temp .lt. 273 .and. wetb(ks).gt.0.) then !snow |
---|
[5f9d14a] | 289 | |
---|
| 290 | !Particle SNOW scavenging coefficient based on Kyro et al 2009, the below-cloud scavenging (Snow efficiency) parameter B (=wetb) from SPECIES file |
---|
| 291 | wetscav= wetb(ks)*10**(bcls(1)+ (bcls(2)*(log10(dquer_m))**(-4))+(bcls(3)*(log10(dquer_m))**(-3))+ (bcls(4)* & |
---|
| 292 | (log10(dquer_m))**(-2))+ (bcls(5)*(log10(dquer_m))**(-1))+ bcls(6)* (prec(1))**(0.5)) |
---|
[8a65cb0] | 293 | |
---|
[5f9d14a] | 294 | endif |
---|
[8a65cb0] | 295 | |
---|
| 296 | ! write(*,*) 'bl-cloud, act_temp=',act_temp, ',prec=',prec(1),',wetscav=', wetscav, ', jpart=',jpart |
---|
| 297 | |
---|
[5f9d14a] | 298 | endif !gas or particle |
---|
| 299 | endif ! positive below-cloud scavenging parameters given in Species file |
---|
[4fbe7a5] | 300 | endif !end below-cloud |
---|
| 301 | |
---|
| 302 | |
---|
| 303 | !******************** |
---|
| 304 | ! IN CLOUD SCAVENGING |
---|
| 305 | !******************** |
---|
[8a65cb0] | 306 | if (clouds_v.lt.4) then !in-cloud |
---|
| 307 | |
---|
[0f20c31] | 308 | ! NIK 13 may 2015: only do incloud if positive in-cloud scavenging parameters are given in species file |
---|
| 309 | if (weta_in(ks).gt.0. .or. wetb_in(ks).gt.0.) then !if positive in-cloud parameters given in SPECIES file (either Ai or Bi) |
---|
[5f9d14a] | 310 | |
---|
[0f20c31] | 311 | ! if negative coefficients (turned off) set to zero for use in equation |
---|
| 312 | if (weta_in(ks).lt.0.) weta_in(ks)=0. |
---|
| 313 | if (wetb_in(ks).lt.0.) wetb_in(ks)=0. |
---|
[5f9d14a] | 314 | |
---|
[0f20c31] | 315 | inc_count=inc_count+1 |
---|
| 316 | |
---|
| 317 | !ZHG liquid water parameterization (CLWC+CIWC) |
---|
| 318 | if (readclouds) then !get cloud water clwc & ciwc units Kg/Kg |
---|
| 319 | cl=clwc(ix,jy,hz,n)+ciwc(ix,jy,hz,n) |
---|
| 320 | else !parameterize cloudwater |
---|
| 321 | cl=2E-7*prec(1)**0.36 |
---|
[8a65cb0] | 322 | endif |
---|
| 323 | |
---|
[0f20c31] | 324 | if (act_temp .le. 253) then |
---|
| 325 | liq_frac=0 |
---|
| 326 | else if (act_temp .ge. 273) then |
---|
| 327 | liq_frac=1 |
---|
| 328 | else |
---|
| 329 | liq_frac =((act_temp-273)/(273-253))**2 |
---|
| 330 | endif |
---|
| 331 | |
---|
[5f9d14a] | 332 | ! ZHG calculate the activated fraction based on the In-cloud scavenging parameters Ai (=weta_in) and Bi (=wetb_in) from SPECIES file |
---|
| 333 | ! frac_act is the combined IN and CCN efficiency |
---|
| 334 | ! The default values are 0.9 for CCN and 0.1 IN |
---|
| 335 | ! This parameterization is based on Verheggen et al. (2007) & Cozich et al. (2006) |
---|
[0f20c31] | 336 | frac_act = liq_frac*weta_in(ks) +(1-liq_frac)*wetb_in(ks) |
---|
[8a65cb0] | 337 | |
---|
[5f9d14a] | 338 | !ZHG Use the activated fraction and the liqid water to calculate the washout |
---|
[0f20c31] | 339 | |
---|
| 340 | ! AEROSOL |
---|
| 341 | if (dquer(ks).gt. 0.) then ! is particle |
---|
| 342 | |
---|
| 343 | S_i= frac_act/cl |
---|
[8a65cb0] | 344 | |
---|
[5f9d14a] | 345 | ! GAS |
---|
[0f20c31] | 346 | else ! is gas |
---|
| 347 | |
---|
| 348 | cle=(1-cl)/(henry(ks)*(r_air/3500.)*act_temp)+cl |
---|
| 349 | S_i=frac_act/cle |
---|
| 350 | |
---|
| 351 | endif ! gas or particle |
---|
[8a65cb0] | 352 | |
---|
[5f9d14a] | 353 | ! scavenging coefficient based on Hertel et al 1995 - using the S_i for either gas or aerosol |
---|
[0f20c31] | 354 | wetscav=S_i*(prec(1)/3.6E6)/clouds_h |
---|
[8a65cb0] | 355 | |
---|
| 356 | ! write(*,*) 'in-cloud, act_temp=',act_temp,',prec=',prec(1),',wetscav=',wetscav,',jpart=',jpart,',clouds_h=,', & |
---|
| 357 | ! clouds_h,',cl=',cl, 'diff to old scheme=', cl-2E-7*prec(1)**0.36 |
---|
[0f20c31] | 358 | endif ! positive in-cloud scavenging parameters given in Species file |
---|
[8a65cb0] | 359 | endif !incloud |
---|
[5f9d14a] | 360 | |
---|
[4fbe7a5] | 361 | !************************************************** |
---|
| 362 | ! CALCULATE DEPOSITION |
---|
| 363 | !************************************************** |
---|
| 364 | ! if (wetscav.le.0) write (*,*) 'neg, or 0 wetscav!' |
---|
| 365 | ! + ,wetscav,cle,cl,act_temp,prec,clouds_h,clouds_v |
---|
[e200b7a] | 366 | |
---|
[4fbe7a5] | 367 | if (wetscav.gt.0.) then |
---|
[e200b7a] | 368 | wetdeposit(ks)=xmass1(jpart,ks)* & |
---|
[5f9d14a] | 369 | (1.-exp(-wetscav*abs(ltsample)))*grfraction(1) ! wet deposition |
---|
[4fbe7a5] | 370 | else ! if no scavenging |
---|
| 371 | wetdeposit(ks)=0. |
---|
| 372 | endif |
---|
[f13406c] | 373 | |
---|
[4fbe7a5] | 374 | restmass = xmass1(jpart,ks)-wetdeposit(ks) |
---|
| 375 | if (ioutputforeachrelease.eq.1) then |
---|
| 376 | kp=npoint(jpart) |
---|
| 377 | else |
---|
| 378 | kp=1 |
---|
| 379 | endif |
---|
| 380 | if (restmass .gt. smallnum) then |
---|
| 381 | xmass1(jpart,ks)=restmass |
---|
| 382 | ! depostatistic |
---|
| 383 | ! wetdepo_sum(ks,kp)=wetdepo_sum(ks,kp)+wetdeposit(ks) |
---|
| 384 | ! depostatistic |
---|
| 385 | else |
---|
| 386 | xmass1(jpart,ks)=0. |
---|
| 387 | endif |
---|
| 388 | ! Correct deposited mass to the last time step when radioactive decay of |
---|
| 389 | ! gridded deposited mass was calculated |
---|
| 390 | if (decay(ks).gt.0.) then |
---|
| 391 | wetdeposit(ks)=wetdeposit(ks)*exp(abs(ldeltat)*decay(ks)) |
---|
| 392 | endif |
---|
[f13406c] | 393 | |
---|
| 394 | |
---|
[4fbe7a5] | 395 | end do !all species |
---|
[e200b7a] | 396 | |
---|
[4fbe7a5] | 397 | ! Sabine Eckhardt, June 2008 create deposition runs only for forward runs |
---|
[e200b7a] | 398 | ! Add the wet deposition to accumulated amount on output grid and nested output grid |
---|
| 399 | !***************************************************************************** |
---|
| 400 | |
---|
[4fbe7a5] | 401 | if (ldirect.eq.1) then |
---|
| 402 | call wetdepokernel(nclass(jpart),wetdeposit,real(xtra1(jpart)), & |
---|
[5f9d14a] | 403 | real(ytra1(jpart)),nage,kp) |
---|
[4fbe7a5] | 404 | if (nested_output.eq.1) call wetdepokernel_nest(nclass(jpart), & |
---|
[5f9d14a] | 405 | wetdeposit,real(xtra1(jpart)),real(ytra1(jpart)),nage,kp) |
---|
[4fbe7a5] | 406 | endif |
---|
[e200b7a] | 407 | |
---|
| 408 | 20 continue |
---|
[4fbe7a5] | 409 | end do ! all particles |
---|
[e200b7a] | 410 | |
---|
[5f9d14a] | 411 | ! count the total number of below-cloud and in-cloud occurences: |
---|
| 412 | tot_blc_count=tot_blc_count+blc_count |
---|
| 413 | tot_inc_count=tot_inc_count+inc_count |
---|
[8a65cb0] | 414 | |
---|
[e200b7a] | 415 | end subroutine wetdepo |
---|