source: flexpart.git/src/get_wetscav.f90 @ c9cf570

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since c9cf570 was c9cf570, checked in by Sabine <sabine.eckhardt@…>, 7 years ago

first version with seperated get_wetscav

  • Property mode set to 100644
File size: 14.1 KB
Line 
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
22subroutine get_wetscav(itime,ltsample,loutnext,jpart,grfraction,inc_count,blc_count,wetscav)
23!                          i      i        i     i     o           o
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
10623  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
141
142! get the level were the actual particle is in
143    do il=2,nz
144      if (height(il).gt.ztra1(jpart)) then
145        hz=il-1
146!        goto 26
147        exit
148      endif
149    end do
150!26  continue
151
152    n=memind(2)
153    if (abs(memtime(1)-interp_time).lt.abs(memtime(2)-interp_time)) &
154         n=memind(1)
155
156    if (ngrid.eq.0) then
157      clouds_v=clouds(ix,jy,hz,n)
158      clouds_h=cloudsh(ix,jy,n)
159    else
160      clouds_v=cloudsn(ix,jy,hz,n,ngrid)
161      clouds_h=cloudshn(ix,jy,n,ngrid)
162    endif
163
164! if there is no precipitation or the particle is above the clouds no
165! scavenging is done
166
167    if (clouds_v.le.1) goto 20
168
169! 1) Parameterization of the the area fraction of the grid cell where the
170!    precipitation occurs: the absolute limit is the total cloud cover, but
171!    for low precipitation rates, an even smaller fraction of the grid cell
172!    is used. Large scale precipitation occurs over larger areas than
173!    convective precipitation.
174!**************************************************************************
175
176    if (lsp.gt.20.) then
177      i=5
178    else if (lsp.gt.8.) then
179      i=4
180    else if (lsp.gt.3.) then
181      i=3
182    else if (lsp.gt.1.) then
183      i=2
184    else
185      i=1
186    endif
187
188    if (convp.gt.20.) then
189      j=5
190    else if (convp.gt.8.) then
191      j=4
192    else if (convp.gt.3.) then
193      j=3
194    else if (convp.gt.1.) then
195      j=2
196    else
197      j=1
198    endif
199
200
201!ZHG oct 2014 : Calculated for 1) both 2) lsp 3) convp
202! Tentatively differentiate the grfraction for lsp and convp for treating differently the two forms
203! for now they are treated the same
204    grfraction(1)=max(0.05,cc*(lsp*lfr(i)+convp*cfr(j))/(lsp+convp))
205    grfraction(2)=max(0.05,cc*(lfr(i)))
206    grfraction(3)=max(0.05,cc*(cfr(j)))
207
208
209! 2) Computation of precipitation rate in sub-grid cell
210!******************************************************
211    prec(1)=(lsp+convp)/grfraction(1)
212    prec(2)=(lsp)/grfraction(2)
213    prec(3)=(convp)/grfraction(3)
214
215
216! 3) Computation of scavenging coefficients for all species
217!    Computation of wet deposition
218!**********************************************************
219
220    do ks=1,nspec      ! loop over species
221      wetdeposit(ks)=0.
222      wetscav=0.   
223
224      if (ngrid.gt.0) then
225        act_temp=ttn(ix,jy,hz,n,ngrid)
226      else
227        act_temp=tt(ix,jy,hz,n)
228      endif
229
230
231!***********************
232! BELOW CLOUD SCAVENGING
233!*********************** 
234      if (clouds_v.ge.4) then !below cloud
235
236! For gas: if positive below-cloud parameters (A or B), and dquer<=0
237!******************************************************************
238        if ((dquer(ks).le.0.).and.(weta_gas(ks).gt.0..or.wetb_gas(ks).gt.0.)) then
239          !        if (weta(ks).gt.0. .or. wetb(ks).gt.0.) then
240          blc_count=blc_count+1
241          wetscav=weta_gas(ks)*prec(1)**wetb_gas(ks)
242
243! For aerosols: if positive below-cloud parameters (Crain/Csnow or B), and dquer>0
244!*********************************************************************************
245        else if ((dquer(ks).gt.0.).and.(crain_aero(ks).gt.0..or.csnow_aero(ks).gt.0.)) then
246          blc_count=blc_count+1
247
248!NIK 17.02.2015
249! For the calculation here particle size needs to be in meter and not um as dquer is
250! changed in readreleases
251! For particles larger than 10 um use the largest size defined in the parameterizations (10um)
252          dquer_m=min(10.,dquer(ks))/1000000. !conversion from um to m
253
254! Rain:
255          if (act_temp .ge. 273. .and. crain_aero(ks).gt.0.)  then
256
257! ZHG 2014 : Particle RAIN scavenging coefficient based on Laakso et al 2003,
258! the below-cloud scavenging (rain efficienty) parameter Crain (=crain_aero) from SPECIES file
259            wetscav=crain_aero(ks)*10**(bclr(1)+(bclr(2)*(log10(dquer_m))**(-4))+ &
260                 & (bclr(3)*(log10(dquer_m))**(-3))+ (bclr(4)*(log10(dquer_m))**(-2))+&
261                 &(bclr(5)*(log10(dquer_m))**(-1))+bclr(6)* (prec(1))**(0.5))
262
263! Snow:
264          elseif (act_temp .lt. 273. .and. csnow_aero(ks).gt.0.)  then
265! ZHG 2014 : Particle SNOW scavenging coefficient based on Kyro et al 2009,
266! the below-cloud scavenging (Snow efficiency) parameter Csnow (=csnow_aero) from SPECIES file
267            wetscav=csnow_aero(ks)*10**(bcls(1)+(bcls(2)*(log10(dquer_m))**(-4))+&
268                 &(bcls(3)*(log10(dquer_m))**(-3))+ (bcls(4)*(log10(dquer_m))**(-2))+&
269                 &(bcls(5)*(log10(dquer_m))**(-1))+ bcls(6)* (prec(1))**(0.5))
270
271          endif
272         
273!             write(*,*) 'bl-cloud, act_temp=',act_temp, ',prec=',prec(1),',wetscav=', wetscav, ', jpart=',jpart
274
275        endif ! gas or particle
276!      endif ! positive below-cloud scavenging parameters given in Species file
277      endif !end BELOW
278
279!********************
280! IN CLOUD SCAVENGING
281!********************
282      if (clouds_v.lt.4) then ! In-cloud
283! NIK 13 may 2015: only do incloud if positive in-cloud scavenging parameters are
284! given in species file, or if gas and positive Henry's constant
285        if ((ccn_aero(ks).gt.0. .or. in_aero(ks).gt.0.).or.(henry(ks).gt.0.and.dquer(ks).le.0)) then
286          inc_count=inc_count+1
287! if negative coefficients (turned off) set to zero for use in equation
288          if (ccn_aero(ks).lt.0.) ccn_aero(ks)=0.
289          if (in_aero(ks).lt.0.) in_aero(ks)=0.
290
291!ZHG 2015 Cloud liquid & ice water (CLWC+CIWC) from ECMWF
292! nested fields
293          if (ngrid.gt.0.and.readclouds_this_nest) then
294            cl=ctwcn(ix,jy,n,ngrid)*(grfraction(1)/cc)
295          else if (ngrid.eq.0.and.readclouds) then
296            cl=ctwc(ix,jy,n)*(grfraction(1)/cc)
297          else                                  !parameterize cloudwater m2/m3
298!ZHG updated parameterization of cloud water to better reproduce the values coming from ECMWF
299            cl=1.6E-6*prec(1)**0.36
300          endif
301
302!ZHG: Calculate the partition between liquid and water phase water.
303          if (act_temp .le. 253.) then
304            liq_frac=0
305          else if (act_temp .ge. 273.) then
306            liq_frac=1
307          else
308            liq_frac =((act_temp-273.)/(273.-253.))**2.
309          end if
310! ZHG: Calculate the aerosol partition based on cloud phase and Ai and Bi
311          frac_act = liq_frac*ccn_aero(ks) +(1-liq_frac)*in_aero(ks)
312
313!ZHG Use the activated fraction and the liqid water to calculate the washout
314
315! AEROSOL
316!********
317          if (dquer(ks).gt.0.) then
318            S_i= frac_act/cl
319
320! GAS
321!****
322          else
323
324            cle=(1-cl)/(henry(ks)*(r_air/3500.)*act_temp)+cl
325!REPLACE to switch old/ new scheme
326          ! S_i=frac_act/cle
327            S_i=1/cle
328          endif ! gas or particle
329
330! scavenging coefficient based on Hertel et al 1995 - using the S_i for either gas or aerosol
331!OLD
332          if ((readclouds.and.ngrid.eq.0).or.(readclouds_this_nest.and.ngrid.gt.0)) then
333            wetscav=incloud_ratio*S_i*(prec(1)/3.6E6)
334          else
335            wetscav=incloud_ratio*S_i*(prec(1)/3.6E6)/clouds_h
336          endif
337        endif ! positive in-cloud scavenging parameters given in Species file
338      endif !incloud
339
340    end do ! loop over species
34120  continue
342
343end subroutine get_wetscav
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG