[c9cf570] | 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 | |
---|
[92a74b2] | 22 | subroutine get_wetscav(itime,ltsample,loutnext,jpart,ks,grfraction,inc_count,blc_count,wetscav) |
---|
[8ee24a5] | 23 | ! i i i i i o o o o |
---|
[c9cf570] | 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 | ! lfr, cfr area fraction covered by precipitation for large scale * |
---|
| 51 | ! and convective precipitation (dependent on prec. rate) * |
---|
| 52 | ! loutnext [s] time for which gridded deposition is next output * |
---|
| 53 | ! loutstep [s] interval at which gridded deposition is output * |
---|
| 54 | ! lsp [mm/h] large scale precipitation rate * |
---|
| 55 | ! ltsample [s] interval over which mass is deposited * |
---|
| 56 | ! prec [mm/h] precipitation rate in subgrid, where precipitation occurs* |
---|
| 57 | ! wetdeposit mass that is wet deposited * |
---|
| 58 | ! wetgrid accumulated deposited mass on output grid * |
---|
| 59 | ! wetscav scavenging coefficient * |
---|
| 60 | ! * |
---|
| 61 | ! Constants: * |
---|
| 62 | ! * |
---|
| 63 | !***************************************************************************** |
---|
| 64 | |
---|
| 65 | use point_mod |
---|
| 66 | use par_mod |
---|
| 67 | use com_mod |
---|
| 68 | |
---|
| 69 | implicit none |
---|
| 70 | |
---|
| 71 | integer :: jpart,itime,ltsample,loutnext,i,j,ix,jy |
---|
| 72 | integer :: ngrid,hz,il,interp_time, n |
---|
| 73 | integer(kind=1) :: clouds_v |
---|
| 74 | integer :: ks, kp |
---|
| 75 | integer :: inc_count, blc_count |
---|
| 76 | ! integer :: n1,n2, icbot,ictop, indcloud !TEST |
---|
| 77 | real :: S_i, act_temp, cl, cle ! in cloud scavenging |
---|
| 78 | real :: clouds_h ! cloud height for the specific grid point |
---|
| 79 | real :: xtn,ytn,lsp,convp,cc,grfraction(3),prec(3),wetscav,totprec |
---|
| 80 | real :: wetdeposit(maxspec),restmass |
---|
| 81 | real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled |
---|
| 82 | !save lfr,cfr |
---|
| 83 | |
---|
| 84 | real, parameter :: lfr(5) = (/ 0.5,0.65,0.8,0.9,0.95/) |
---|
| 85 | real, parameter :: cfr(5) = (/ 0.4,0.55,0.7,0.8,0.9 /) |
---|
| 86 | |
---|
| 87 | !ZHG aerosol below-cloud scavenging removal polynomial constants for rain and snow |
---|
| 88 | real, parameter :: bclr(6) = (/274.35758, 332839.59273, 226656.57259, 58005.91340, 6588.38582, 0.244984/) !rain (Laakso et al 2003) |
---|
| 89 | real, parameter :: bcls(6) = (/22.7, 0.0, 0.0, 1321.0, 381.0, 0.0/) !now (Kyro et al 2009) |
---|
| 90 | real :: frac_act, liq_frac, dquer_m |
---|
| 91 | |
---|
| 92 | real :: Si_dummy, wetscav_dummy |
---|
| 93 | logical :: readclouds_this_nest |
---|
| 94 | |
---|
| 95 | |
---|
| 96 | ! Determine which nesting level to be used |
---|
| 97 | !***************************************** |
---|
| 98 | ngrid=0 |
---|
| 99 | do j=numbnests,1,-1 |
---|
| 100 | if ((xtra1(jpart).gt.xln(j)).and.(xtra1(jpart).lt.xrn(j)).and. & |
---|
| 101 | (ytra1(jpart).gt.yln(j)).and.(ytra1(jpart).lt.yrn(j))) then |
---|
| 102 | ngrid=j |
---|
| 103 | goto 23 |
---|
| 104 | endif |
---|
| 105 | end do |
---|
| 106 | 23 continue |
---|
| 107 | |
---|
| 108 | |
---|
| 109 | ! Determine nested grid coordinates |
---|
| 110 | !********************************** |
---|
| 111 | readclouds_this_nest=.false. |
---|
| 112 | |
---|
| 113 | if (ngrid.gt.0) then |
---|
| 114 | xtn=(xtra1(jpart)-xln(ngrid))*xresoln(ngrid) |
---|
| 115 | ytn=(ytra1(jpart)-yln(ngrid))*yresoln(ngrid) |
---|
| 116 | ix=int(xtn) |
---|
| 117 | jy=int(ytn) |
---|
| 118 | if (readclouds_nest(ngrid)) readclouds_this_nest=.true. |
---|
| 119 | else |
---|
| 120 | ix=int(xtra1(jpart)) |
---|
| 121 | jy=int(ytra1(jpart)) |
---|
| 122 | endif |
---|
| 123 | |
---|
| 124 | |
---|
| 125 | ! Interpolate large scale precipitation, convective precipitation and |
---|
| 126 | ! total cloud cover |
---|
| 127 | ! Note that interpolated time refers to itime-0.5*ltsample [PS] |
---|
| 128 | !******************************************************************** |
---|
| 129 | interp_time=nint(itime-0.5*ltsample) |
---|
| 130 | |
---|
| 131 | if (ngrid.eq.0) then |
---|
| 132 | call interpol_rain(lsprec,convprec,tcc,nxmax,nymax, & |
---|
| 133 | 1,nx,ny,memind,real(xtra1(jpart)),real(ytra1(jpart)),1, & |
---|
| 134 | memtime(1),memtime(2),interp_time,lsp,convp,cc) |
---|
| 135 | else |
---|
| 136 | call interpol_rain_nests(lsprecn,convprecn,tccn, & |
---|
| 137 | nxmaxn,nymaxn,1,maxnests,ngrid,nxn,nyn,memind,xtn,ytn,1, & |
---|
| 138 | memtime(1),memtime(2),interp_time,lsp,convp,cc) |
---|
| 139 | endif |
---|
| 140 | |
---|
[8ee24a5] | 141 | ! If total precipitation is less than 0.01 mm/h - no scavenging occurs |
---|
| 142 | if ((lsp.lt.0.01).and.(convp.lt.0.01)) goto 20 |
---|
[c9cf570] | 143 | |
---|
| 144 | ! get the level were the actual particle is in |
---|
| 145 | do il=2,nz |
---|
| 146 | if (height(il).gt.ztra1(jpart)) then |
---|
| 147 | hz=il-1 |
---|
| 148 | ! goto 26 |
---|
| 149 | exit |
---|
| 150 | endif |
---|
| 151 | end do |
---|
| 152 | !26 continue |
---|
| 153 | |
---|
| 154 | n=memind(2) |
---|
| 155 | if (abs(memtime(1)-interp_time).lt.abs(memtime(2)-interp_time)) & |
---|
| 156 | n=memind(1) |
---|
| 157 | |
---|
| 158 | if (ngrid.eq.0) then |
---|
| 159 | clouds_v=clouds(ix,jy,hz,n) |
---|
| 160 | clouds_h=cloudsh(ix,jy,n) |
---|
| 161 | else |
---|
| 162 | clouds_v=cloudsn(ix,jy,hz,n,ngrid) |
---|
| 163 | clouds_h=cloudshn(ix,jy,n,ngrid) |
---|
| 164 | endif |
---|
| 165 | |
---|
| 166 | ! if there is no precipitation or the particle is above the clouds no |
---|
| 167 | ! scavenging is done |
---|
| 168 | |
---|
| 169 | if (clouds_v.le.1) goto 20 |
---|
| 170 | |
---|
| 171 | ! 1) Parameterization of the the area fraction of the grid cell where the |
---|
| 172 | ! precipitation occurs: the absolute limit is the total cloud cover, but |
---|
| 173 | ! for low precipitation rates, an even smaller fraction of the grid cell |
---|
| 174 | ! is used. Large scale precipitation occurs over larger areas than |
---|
| 175 | ! convective precipitation. |
---|
| 176 | !************************************************************************** |
---|
| 177 | |
---|
| 178 | if (lsp.gt.20.) then |
---|
| 179 | i=5 |
---|
| 180 | else if (lsp.gt.8.) then |
---|
| 181 | i=4 |
---|
| 182 | else if (lsp.gt.3.) then |
---|
| 183 | i=3 |
---|
| 184 | else if (lsp.gt.1.) then |
---|
| 185 | i=2 |
---|
| 186 | else |
---|
| 187 | i=1 |
---|
| 188 | endif |
---|
| 189 | |
---|
| 190 | if (convp.gt.20.) then |
---|
| 191 | j=5 |
---|
| 192 | else if (convp.gt.8.) then |
---|
| 193 | j=4 |
---|
| 194 | else if (convp.gt.3.) then |
---|
| 195 | j=3 |
---|
| 196 | else if (convp.gt.1.) then |
---|
| 197 | j=2 |
---|
| 198 | else |
---|
| 199 | j=1 |
---|
| 200 | endif |
---|
| 201 | |
---|
| 202 | |
---|
| 203 | !ZHG oct 2014 : Calculated for 1) both 2) lsp 3) convp |
---|
| 204 | ! Tentatively differentiate the grfraction for lsp and convp for treating differently the two forms |
---|
| 205 | ! for now they are treated the same |
---|
| 206 | grfraction(1)=max(0.05,cc*(lsp*lfr(i)+convp*cfr(j))/(lsp+convp)) |
---|
| 207 | grfraction(2)=max(0.05,cc*(lfr(i))) |
---|
| 208 | grfraction(3)=max(0.05,cc*(cfr(j))) |
---|
| 209 | |
---|
| 210 | |
---|
| 211 | ! 2) Computation of precipitation rate in sub-grid cell |
---|
| 212 | !****************************************************** |
---|
| 213 | prec(1)=(lsp+convp)/grfraction(1) |
---|
| 214 | prec(2)=(lsp)/grfraction(2) |
---|
| 215 | prec(3)=(convp)/grfraction(3) |
---|
| 216 | |
---|
| 217 | |
---|
| 218 | ! 3) Computation of scavenging coefficients for all species |
---|
| 219 | ! Computation of wet deposition |
---|
| 220 | !********************************************************** |
---|
| 221 | |
---|
| 222 | wetdeposit(ks)=0. |
---|
| 223 | wetscav=0. |
---|
| 224 | |
---|
| 225 | if (ngrid.gt.0) then |
---|
| 226 | act_temp=ttn(ix,jy,hz,n,ngrid) |
---|
| 227 | else |
---|
| 228 | act_temp=tt(ix,jy,hz,n) |
---|
| 229 | endif |
---|
| 230 | |
---|
| 231 | |
---|
| 232 | !*********************** |
---|
| 233 | ! BELOW CLOUD SCAVENGING |
---|
| 234 | !*********************** |
---|
| 235 | if (clouds_v.ge.4) then !below cloud |
---|
| 236 | |
---|
| 237 | ! For gas: if positive below-cloud parameters (A or B), and dquer<=0 |
---|
| 238 | !****************************************************************** |
---|
| 239 | if ((dquer(ks).le.0.).and.(weta_gas(ks).gt.0..or.wetb_gas(ks).gt.0.)) then |
---|
| 240 | ! if (weta(ks).gt.0. .or. wetb(ks).gt.0.) then |
---|
| 241 | blc_count=blc_count+1 |
---|
| 242 | wetscav=weta_gas(ks)*prec(1)**wetb_gas(ks) |
---|
| 243 | |
---|
| 244 | ! For aerosols: if positive below-cloud parameters (Crain/Csnow or B), and dquer>0 |
---|
| 245 | !********************************************************************************* |
---|
| 246 | else if ((dquer(ks).gt.0.).and.(crain_aero(ks).gt.0..or.csnow_aero(ks).gt.0.)) then |
---|
| 247 | blc_count=blc_count+1 |
---|
| 248 | |
---|
| 249 | !NIK 17.02.2015 |
---|
| 250 | ! For the calculation here particle size needs to be in meter and not um as dquer is |
---|
| 251 | ! changed in readreleases |
---|
| 252 | ! For particles larger than 10 um use the largest size defined in the parameterizations (10um) |
---|
| 253 | dquer_m=min(10.,dquer(ks))/1000000. !conversion from um to m |
---|
| 254 | |
---|
| 255 | ! Rain: |
---|
| 256 | if (act_temp .ge. 273. .and. crain_aero(ks).gt.0.) then |
---|
| 257 | |
---|
| 258 | ! ZHG 2014 : Particle RAIN scavenging coefficient based on Laakso et al 2003, |
---|
| 259 | ! the below-cloud scavenging (rain efficienty) parameter Crain (=crain_aero) from SPECIES file |
---|
| 260 | wetscav=crain_aero(ks)*10**(bclr(1)+(bclr(2)*(log10(dquer_m))**(-4))+ & |
---|
| 261 | & (bclr(3)*(log10(dquer_m))**(-3))+ (bclr(4)*(log10(dquer_m))**(-2))+& |
---|
| 262 | &(bclr(5)*(log10(dquer_m))**(-1))+bclr(6)* (prec(1))**(0.5)) |
---|
| 263 | |
---|
| 264 | ! Snow: |
---|
| 265 | elseif (act_temp .lt. 273. .and. csnow_aero(ks).gt.0.) then |
---|
| 266 | ! ZHG 2014 : Particle SNOW scavenging coefficient based on Kyro et al 2009, |
---|
| 267 | ! the below-cloud scavenging (Snow efficiency) parameter Csnow (=csnow_aero) from SPECIES file |
---|
| 268 | wetscav=csnow_aero(ks)*10**(bcls(1)+(bcls(2)*(log10(dquer_m))**(-4))+& |
---|
| 269 | &(bcls(3)*(log10(dquer_m))**(-3))+ (bcls(4)*(log10(dquer_m))**(-2))+& |
---|
| 270 | &(bcls(5)*(log10(dquer_m))**(-1))+ bcls(6)* (prec(1))**(0.5)) |
---|
| 271 | |
---|
| 272 | endif |
---|
| 273 | |
---|
| 274 | ! write(*,*) 'bl-cloud, act_temp=',act_temp, ',prec=',prec(1),',wetscav=', wetscav, ', jpart=',jpart |
---|
| 275 | |
---|
| 276 | endif ! gas or particle |
---|
| 277 | ! endif ! positive below-cloud scavenging parameters given in Species file |
---|
| 278 | endif !end BELOW |
---|
| 279 | |
---|
| 280 | !******************** |
---|
| 281 | ! IN CLOUD SCAVENGING |
---|
| 282 | !******************** |
---|
| 283 | if (clouds_v.lt.4) then ! In-cloud |
---|
| 284 | ! NIK 13 may 2015: only do incloud if positive in-cloud scavenging parameters are |
---|
| 285 | ! given in species file, or if gas and positive Henry's constant |
---|
| 286 | if ((ccn_aero(ks).gt.0. .or. in_aero(ks).gt.0.).or.(henry(ks).gt.0.and.dquer(ks).le.0)) then |
---|
| 287 | inc_count=inc_count+1 |
---|
| 288 | ! if negative coefficients (turned off) set to zero for use in equation |
---|
| 289 | if (ccn_aero(ks).lt.0.) ccn_aero(ks)=0. |
---|
| 290 | if (in_aero(ks).lt.0.) in_aero(ks)=0. |
---|
| 291 | |
---|
| 292 | !ZHG 2015 Cloud liquid & ice water (CLWC+CIWC) from ECMWF |
---|
| 293 | ! nested fields |
---|
| 294 | if (ngrid.gt.0.and.readclouds_this_nest) then |
---|
| 295 | cl=ctwcn(ix,jy,n,ngrid)*(grfraction(1)/cc) |
---|
| 296 | else if (ngrid.eq.0.and.readclouds) then |
---|
| 297 | cl=ctwc(ix,jy,n)*(grfraction(1)/cc) |
---|
| 298 | else !parameterize cloudwater m2/m3 |
---|
| 299 | !ZHG updated parameterization of cloud water to better reproduce the values coming from ECMWF |
---|
| 300 | cl=1.6E-6*prec(1)**0.36 |
---|
| 301 | endif |
---|
| 302 | |
---|
| 303 | !ZHG: Calculate the partition between liquid and water phase water. |
---|
| 304 | if (act_temp .le. 253.) then |
---|
| 305 | liq_frac=0 |
---|
| 306 | else if (act_temp .ge. 273.) then |
---|
| 307 | liq_frac=1 |
---|
| 308 | else |
---|
| 309 | liq_frac =((act_temp-273.)/(273.-253.))**2. |
---|
| 310 | end if |
---|
| 311 | ! ZHG: Calculate the aerosol partition based on cloud phase and Ai and Bi |
---|
| 312 | frac_act = liq_frac*ccn_aero(ks) +(1-liq_frac)*in_aero(ks) |
---|
| 313 | |
---|
| 314 | !ZHG Use the activated fraction and the liqid water to calculate the washout |
---|
| 315 | |
---|
| 316 | ! AEROSOL |
---|
| 317 | !******** |
---|
| 318 | if (dquer(ks).gt.0.) then |
---|
| 319 | S_i= frac_act/cl |
---|
| 320 | |
---|
| 321 | ! GAS |
---|
| 322 | !**** |
---|
| 323 | else |
---|
| 324 | |
---|
| 325 | cle=(1-cl)/(henry(ks)*(r_air/3500.)*act_temp)+cl |
---|
| 326 | !REPLACE to switch old/ new scheme |
---|
| 327 | ! S_i=frac_act/cle |
---|
| 328 | S_i=1/cle |
---|
| 329 | endif ! gas or particle |
---|
| 330 | |
---|
| 331 | ! scavenging coefficient based on Hertel et al 1995 - using the S_i for either gas or aerosol |
---|
| 332 | !OLD |
---|
| 333 | if ((readclouds.and.ngrid.eq.0).or.(readclouds_this_nest.and.ngrid.gt.0)) then |
---|
| 334 | wetscav=incloud_ratio*S_i*(prec(1)/3.6E6) |
---|
| 335 | else |
---|
| 336 | wetscav=incloud_ratio*S_i*(prec(1)/3.6E6)/clouds_h |
---|
| 337 | endif |
---|
| 338 | endif ! positive in-cloud scavenging parameters given in Species file |
---|
| 339 | endif !incloud |
---|
| 340 | |
---|
| 341 | 20 continue |
---|
| 342 | |
---|
| 343 | end subroutine get_wetscav |
---|