source: flexpart.git/src/get_wetscav.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 12.7 KB
Line 
1subroutine get_wetscav(itime,ltsample,loutnext,jpart,ks,grfraction,inc_count,blc_count,wetscav)
2!                          i      i        i     i   i    o           o          o       o
3!*****************************************************************************
4!                                                                            *
5! Calculation of wet deposition using the concept of scavenging coefficients.*
6! For lack of detailed information, washout and rainout are jointly treated. *
7! It is assumed that precipitation does not occur uniformly within the whole *
8! grid cell, but that only a fraction of the grid cell experiences rainfall. *
9! This fraction is parameterized from total cloud cover and rates of large   *
10! scale and convective precipitation.                                        *
11!                                                                            *
12!    Author: A. Stohl                                                        *
13!                                                                            *
14!    1 December 1996                                                         *
15!                                                                            *
16! Correction by Petra Seibert, Sept 2002:                                    *
17! use centred precipitation data for integration                             *
18! Code may not be correct for decay of deposition!                           *
19!                                                                            *
20!*****************************************************************************
21!                                                                            *
22! Variables:                                                                 *
23! cc [0-1]           total cloud cover                                       *
24! convp [mm/h]       convective precipitation rate                           *
25! grfraction [0-1]   fraction of grid, for which precipitation occurs        *
26! ix,jy              indices of output grid cell for each particle           *
27! itime [s]          actual simulation time [s]                              *
28! jpart              particle index                                          *
29! lfr, cfr           area fraction covered by precipitation for large scale  *
30!                    and convective precipitation (dependent on prec. rate)  *
31! loutnext [s]       time for which gridded deposition is next output        *
32! loutstep [s]       interval at which gridded deposition is output          *
33! lsp [mm/h]         large scale precipitation rate                          *
34! ltsample [s]       interval over which mass is deposited                   *
35! prec [mm/h]        precipitation rate in subgrid, where precipitation occurs*
36! wetgrid            accumulated deposited mass on output grid               *
37! wetscav            scavenging coefficient                                  *
38!                                                                            *
39! Constants:                                                                 *
40!                                                                            *
41!*****************************************************************************
42
43  use point_mod
44  use par_mod
45  use com_mod
46
47  implicit none
48
49  integer :: jpart,itime,ltsample,loutnext,i,j,ix,jy
50  integer :: ngrid,hz,il,interp_time, n
51  integer(kind=1) :: clouds_v
52  integer :: ks, kp
53  integer(selected_int_kind(16)), dimension(nspec) :: blc_count, inc_count
54
55!  integer :: n1,n2, icbot,ictop, indcloud !TEST
56  real :: S_i, act_temp, cl, cle ! in cloud scavenging
57  real :: clouds_h ! cloud height for the specific grid point
58  real :: xtn,ytn,lsp,convp,cc,grfraction(3),prec(3),wetscav,totprec
59  real :: restmass
60  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
61!save lfr,cfr
62
63  real, parameter :: lfr(5) = (/ 0.5,0.65,0.8,0.9,0.95/)
64  real, parameter :: cfr(5) = (/ 0.4,0.55,0.7,0.8,0.9 /)
65
66!ZHG aerosol below-cloud scavenging removal polynomial constants for rain and snow
67  real, parameter :: bclr(6) = (/274.35758, 332839.59273, 226656.57259, 58005.91340, 6588.38582, 0.244984/) !rain (Laakso et al 2003)
68  real, parameter :: bcls(6) = (/22.7, 0.0, 0.0, 1321.0, 381.0, 0.0/) !now (Kyro et al 2009)
69  real :: frac_act, liq_frac, ice_frac, dquer_m
70
71  real    :: Si_dummy, wetscav_dummy
72  logical :: readclouds_this_nest
73
74
75   wetscav=0.
76
77! Determine which nesting level to be used
78!*****************************************
79    ngrid=0
80    do j=numbnests,1,-1
81      if ((xtra1(jpart).gt.xln(j)).and.(xtra1(jpart).lt.xrn(j)).and. &
82           (ytra1(jpart).gt.yln(j)).and.(ytra1(jpart).lt.yrn(j))) then
83        ngrid=j
84        goto 23
85      endif
86    end do
8723  continue
88
89
90! Determine nested grid coordinates
91!**********************************
92    readclouds_this_nest=.false.
93
94    if (ngrid.gt.0) then
95      xtn=(xtra1(jpart)-xln(ngrid))*xresoln(ngrid)
96      ytn=(ytra1(jpart)-yln(ngrid))*yresoln(ngrid)
97      ix=int(xtn)
98      jy=int(ytn)
99      if (readclouds_nest(ngrid)) readclouds_this_nest=.true.
100    else
101      ix=int(xtra1(jpart))
102      jy=int(ytra1(jpart))
103    endif
104
105! Interpolate large scale precipitation, convective precipitation and
106! total cloud cover
107! Note that interpolated time refers to itime-0.5*ltsample [PS]
108!********************************************************************
109    interp_time=nint(itime-0.5*ltsample)
110
111    n=memind(2)
112    if (abs(memtime(1)-interp_time).lt.abs(memtime(2)-interp_time)) &
113         n=memind(1)
114
115    if (ngrid.eq.0) then
116      call interpol_rain(lsprec,convprec,tcc,nxmax,nymax, &
117           1,nx,ny,n,real(xtra1(jpart)),real(ytra1(jpart)),1, &
118           memtime(1),memtime(2),interp_time,lsp,convp,cc)
119    else
120      call interpol_rain_nests(lsprecn,convprecn,tccn, &
121           nxmaxn,nymaxn,1,maxnests,ngrid,nxn,nyn,n,xtn,ytn,1, &
122           memtime(1),memtime(2),interp_time,lsp,convp,cc)
123    endif
124
125!  If total precipitation is less than 0.01 mm/h - no scavenging occurs
126    if ((lsp.lt.0.01).and.(convp.lt.0.01)) goto 20
127
128! get the level were the actual particle is in
129    do il=2,nz
130      if (height(il).gt.ztra1(jpart)) then
131        hz=il-1
132        exit
133      endif
134    end do
135
136    if (ngrid.eq.0) then
137      clouds_v=clouds(ix,jy,hz,n)
138      clouds_h=cloudsh(ix,jy,n)
139    else
140      clouds_v=cloudsn(ix,jy,hz,n,ngrid)
141      clouds_h=cloudshn(ix,jy,n,ngrid)
142    endif
143
144! if there is no precipitation or the particle is above the clouds no
145! scavenging is done
146
147    if (clouds_v.le.1) goto 20
148
149! 1) Parameterization of the the area fraction of the grid cell where the
150!    precipitation occurs: the absolute limit is the total cloud cover, but
151!    for low precipitation rates, an even smaller fraction of the grid cell
152!    is used. Large scale precipitation occurs over larger areas than
153!    convective precipitation.
154!**************************************************************************
155
156    if (lsp.gt.20.) then
157      i=5
158    else if (lsp.gt.8.) then
159      i=4
160    else if (lsp.gt.3.) then
161      i=3
162    else if (lsp.gt.1.) then
163      i=2
164    else
165      i=1
166    endif
167
168    if (convp.gt.20.) then
169      j=5
170    else if (convp.gt.8.) then
171      j=4
172    else if (convp.gt.3.) then
173      j=3
174    else if (convp.gt.1.) then
175      j=2
176    else
177      j=1
178    endif
179
180
181!ZHG oct 2014 : Calculated for 1) both 2) lsp 3) convp - 2 and 3 not used removed by SE
182! Tentatively differentiate the grfraction for lsp and convp for treating differently the two forms
183! for now they are treated the same
184    grfraction(1)=max(0.05,cc*(lsp*lfr(i)+convp*cfr(j))/(lsp+convp))
185
186! 2) Computation of precipitation rate in sub-grid cell
187!******************************************************
188    prec(1)=(lsp+convp)/grfraction(1)
189
190! 3) Computation of scavenging coefficients for all species
191!    Computation of wet deposition
192!**********************************************************
193
194      if (ngrid.gt.0) then
195        act_temp=ttn(ix,jy,hz,n,ngrid)
196      else
197        act_temp=tt(ix,jy,hz,n)
198      endif
199
200!***********************
201! BELOW CLOUD SCAVENGING
202!*********************** 
203      if (clouds_v.ge.4) then !below cloud
204
205! For gas: if positive below-cloud parameters (A or B), and dquer<=0
206!******************************************************************
207        if ((dquer(ks).le.0.).and.(weta_gas(ks).gt.0..or.wetb_gas(ks).gt.0.)) then
208          blc_count(ks)=blc_count(ks)+1
209          wetscav=weta_gas(ks)*prec(1)**wetb_gas(ks)
210
211! For aerosols: if positive below-cloud parameters (Crain/Csnow or B), and dquer>0
212!*********************************************************************************
213        else if ((dquer(ks).gt.0.).and.(crain_aero(ks).gt.0..or.csnow_aero(ks).gt.0.)) then
214          blc_count(ks)=blc_count(ks)+1
215
216!NIK 17.02.2015
217! For the calculation here particle size needs to be in meter and not um as dquer is
218! changed in readreleases
219! For particles larger than 10 um use the largest size defined in the parameterizations (10um)
220          dquer_m=min(10.,dquer(ks))/1000000. !conversion from um to m
221
222! Rain:
223          if (act_temp .ge. 273. .and. crain_aero(ks).gt.0.)  then
224
225! ZHG 2014 : Particle RAIN scavenging coefficient based on Laakso et al 2003,
226! the below-cloud scavenging (rain efficienty) parameter Crain (=crain_aero) from SPECIES file
227            wetscav=crain_aero(ks)*10**(bclr(1)+(bclr(2)*(log10(dquer_m))**(-4))+ &
228                 & (bclr(3)*(log10(dquer_m))**(-3))+ (bclr(4)*(log10(dquer_m))**(-2))+&
229                 &(bclr(5)*(log10(dquer_m))**(-1))+bclr(6)* (prec(1))**(0.5))
230
231! Snow:
232          elseif (act_temp .lt. 273. .and. csnow_aero(ks).gt.0.)  then
233! ZHG 2014 : Particle SNOW scavenging coefficient based on Kyro et al 2009,
234! the below-cloud scavenging (Snow efficiency) parameter Csnow (=csnow_aero) from SPECIES file
235            wetscav=csnow_aero(ks)*10**(bcls(1)+(bcls(2)*(log10(dquer_m))**(-4))+&
236                 &(bcls(3)*(log10(dquer_m))**(-3))+ (bcls(4)*(log10(dquer_m))**(-2))+&
237                 &(bcls(5)*(log10(dquer_m))**(-1))+ bcls(6)* (prec(1))**(0.5))
238
239          endif
240         
241        endif ! gas or particle
242!      endif ! positive below-cloud scavenging parameters given in Species file
243      endif !end BELOW
244
245!********************
246! IN CLOUD SCAVENGING
247!********************
248      if (clouds_v.lt.4) then ! In-cloud
249! NIK 13 may 2015: only do incloud if positive in-cloud scavenging parameters are
250! given in species file, or if gas and positive Henry's constant
251        if ((ccn_aero(ks).gt.0. .or. in_aero(ks).gt.0.).or.(henry(ks).gt.0.and.dquer(ks).le.0)) then
252          inc_count(ks)=inc_count(ks)+1
253! if negative coefficients (turned off) set to zero for use in equation
254          if (ccn_aero(ks).lt.0.) ccn_aero(ks)=0.
255          if (in_aero(ks).lt.0.) in_aero(ks)=0.
256
257!ZHG 2015 Cloud liquid & ice water (CLWC+CIWC) from ECMWF
258! nested fields
259          if (ngrid.gt.0.and.readclouds_this_nest) then
260            cl=ctwcn(ix,jy,n,ngrid)*(grfraction(1)/cc)
261          else if (ngrid.eq.0.and.readclouds) then
262            cl=ctwc(ix,jy,n)*(grfraction(1)/cc)
263          else                                  !parameterize cloudwater m2/m3
264!ZHG updated parameterization of cloud water to better reproduce the values coming from ECMWF
265! sec test
266!           cl=1E6*1E-7*prec(1)**0.3 !Sec GFS new
267            cl=1E6*2E-7*prec(1)**0.36 !Sec ECMWF new, is also suitable for GFS
268!           cl=2E-7*prec(1)**0.36 !Andreas
269!           cl=1.6E-6*prec(1)**0.36 !Henrik
270          endif
271
272!ZHG: Calculate the partition between liquid and water phase water.
273          if (act_temp .le. 253.) then
274            liq_frac=0
275            ice_frac=1
276          else if (act_temp .ge. 273.) then
277            liq_frac=1
278            ice_frac=0
279          else
280! sec bugfix after FLEXPART paper review, liq_frac was 1-liq_frac
281! IP bugfix v10.4, calculate ice_frac and liq_frac
282            ice_frac= ((act_temp-273.)/(273.-253.))**2.
283            !liq_frac = 1-ice_frac   !((act_temp-253.)/(273.-253.))**2.
284            liq_frac=max(0.,1.-ice_frac)
285          end if
286! ZHG: Calculate the aerosol partition based on cloud phase and Ai and Bi
287!         frac_act = liq_frac*ccn_aero(ks) +(1-liq_frac)*in_aero(ks)
288! IP, use ice_frac and liq_frac
289          frac_act = liq_frac*ccn_aero(ks) + ice_frac*in_aero(ks)
290
291!ZHG Use the activated fraction and the liqid water to calculate the washout
292
293! AEROSOL
294!********
295          if (dquer(ks).gt.0.) then
296            S_i= frac_act/cl
297! GAS
298!****
299          else
300            cle=(1-cl)/(henry(ks)*(r_air/3500.)*act_temp)+cl
301            S_i=1/cle
302          endif ! gas or particle
303
304! scavenging coefficient based on Hertel et al 1995 - using the S_i for either gas or aerosol
305!SEC wetscav fix, the cloud height is no longer needed, it gives wrong results
306            wetscav=incloud_ratio*S_i*(prec(1)/3.6E6)
307        endif ! positive in-cloud scavenging parameters given in Species file
308      endif !incloud
309
310
31120  continue
312
313end subroutine get_wetscav
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG