source: flexpart.git/src/get_wetscav.f90

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

add SPDX-License-Identifier to all .f90 files

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