source: flexpart.git/src/wetdepo.f90 @ f75967d

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since f75967d was f75967d, checked in by Espen Sollum ATMOS <eso@…>, 8 years ago

Fixed a bug with number of occurances of wet scavenging written to stdout

  • Property mode set to 100644
File size: 18.0 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 wetdepo(itime,ltsample,loutnext)
23  !                  i      i        i
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
73  integer :: ngrid,itage,nage,hz,il,interp_time, n, clouds_v
74  integer :: ks, kp
75!  integer :: n1,n2, icbot,ictop, indcloud !TEST
76  real :: S_i, act_temp, cl, cle ! in cloud scavenging
77  real :: clouds_h ! cloud height for the specific grid point
78  real :: xtn,ytn,lsp,convp,cc,grfraction(3),prec(3),wetscav,totprec
79  real :: wetdeposit(maxspec),restmass
80  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
81  !save lfr,cfr
82
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 /)
85
86!ZHG aerosol below-cloud scavenging removal polynomial constants for rain and snow
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)
89  real :: frac_act, liq_frac, dquer_m
90
91  integer :: blc_count, inc_count
92  real    :: Si_dummy, wetscav_dummy
93
94
95
96  ! Compute interval since radioactive decay of deposited mass was computed
97  !************************************************************************
98
99  if (itime.le.loutnext) then
100    ldeltat=itime-(loutnext-loutstep)
101  else                                  ! first half of next interval
102    ldeltat=itime-loutnext
103  endif
104
105  ! Loop over all particles
106  !************************
107
108  blc_count=0
109  inc_count=0
110
111  do jpart=1,numpart
112
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
119
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
12533  continue
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. &
134      (ytra1(jpart).gt.yln(j)).and.(ytra1(jpart).lt.yrn(j))) then
135        ngrid=j
136        goto 23
137      endif
138    end do
13923  continue
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)
161   
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
171
172!  If total precipitation is less than 0.01 mm/h - no scavenging occurs
173    if ((lsp.lt.0.01).and.(convp.lt.0.01)) goto 20
174
175  ! get the level were the actual particle is in
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
18226  continue
183
184    n=memind(2)
185    if (abs(memtime(1)-interp_time).lt.abs(memtime(2)-interp_time)) &
186         n=memind(1)
187
188    if (ngrid.eq.0) then
189      clouds_v=clouds(ix,jy,hz,n)
190      clouds_h=cloudsh(ix,jy,n)
191    else
192      ! new removal not implemented for nests yet
193      clouds_v=cloudsn(ix,jy,hz,n,ngrid)
194      clouds_h=cloudsnh(ix,jy,n,ngrid)
195    endif
196
197  ! if there is no precipitation or the particle is above the clouds no
198  ! scavenging is done
199
200    if (clouds_v.le.1) goto 20
201 
202  ! 1) Parameterization of the the area fraction of the grid cell where the
203  !    precipitation occurs: the absolute limit is the total cloud cover, but
204  !    for low precipitation rates, an even smaller fraction of the grid cell
205  !    is used. Large scale precipitation occurs over larger areas than
206  !    convective precipitation.
207  !**************************************************************************
208
209    if (lsp.gt.20.) then
210      i=5
211    else if (lsp.gt.8.) then
212      i=4
213    else if (lsp.gt.3.) then
214      i=3
215    else if (lsp.gt.1.) then
216      i=2
217    else
218      i=1
219    endif
220
221    if (convp.gt.20.) then
222      j=5
223    else if (convp.gt.8.) then
224      j=4
225    else if (convp.gt.3.) then
226      j=3
227    else if (convp.gt.1.) then
228      j=2
229    else
230      j=1
231    endif
232
233
234  !ZHG oct 2014 : Calculated for 1) both 2) lsp 3) convp
235  ! Tentatively differentiate the grfraction for lsp and convp for treating differently the two forms
236  ! for now they are treated the same
237    grfraction(1)=max(0.05,cc*(lsp*lfr(i)+convp*cfr(j))/(lsp+convp))
238    grfraction(2)=max(0.05,cc*(lfr(i)))
239    grfraction(3)=max(0.05,cc*(cfr(j)))
240
241
242  ! 2) Computation of precipitation rate in sub-grid cell
243  !******************************************************
244    prec(1)=(lsp+convp)/grfraction(1)
245    prec(2)=(lsp)/grfraction(2)
246    prec(3)=(convp)/grfraction(3)
247
248
249  ! 3) Computation of scavenging coefficients for all species
250  !    Computation of wet deposition
251  !**********************************************************
252
253    do ks=1,nspec      ! loop over species
254      wetdeposit(ks)=0.
255      wetscav=0.   
256
257      !ZHG test if it nested?
258      if (ngrid.gt.0) then
259        act_temp=ttn(ix,jy,hz,n,ngrid)
260      else
261        act_temp=tt(ix,jy,hz,n)
262      endif
263     
264
265  !***********************
266  ! BELOW CLOUD SCAVENGING
267  !*********************** 
268      if (clouds_v.ge.4) then !below cloud
269
270        if (weta(ks).gt.0. .or. wetb(ks).gt.0.) then !if positive below-cloud parameters given in SPECIES file (either A or B)
271          blc_count=blc_count+1
272
273!GAS
274          if (dquer(ks) .le. 0.) then  !is gas
275            wetscav=weta(ks)*prec(1)**wetb(ks)
276           
277!AEROSOL
278          else !is particle
279!NIK 17.02.2015
280! For the calculation here particle size needs to be in meter and not um as dquer is changed to in readreleases
281! for particles larger than 10 um use the largest size defined in the parameterizations (10um)
282            dquer_m=min(10.,dquer(ks))/1000000. !conversion from um to m
283            if (act_temp .ge. 273 .and. weta(ks).gt.0.)  then !Rain
284              ! ZHG 2014 : Particle RAIN scavenging coefficient based on Laakso et al 2003,
285              ! the below-cloud scavenging (rain efficienty)
286              ! parameter A (=weta) from SPECIES file
287              wetscav= weta(ks)*10**(bclr(1)+ (bclr(2)*(log10(dquer_m))**(-4))+(bclr(3)*(log10(dquer_m))**(-3))+ (bclr(4)* &
288                   (log10(dquer_m))**(-2))+ (bclr(5)*(log10(dquer_m))**(-1))+ bclr(6)* (prec(1))**(0.5))
289
290            elseif (act_temp .lt. 273 .and. wetb(ks).gt.0.)  then ! Snow
291              ! ZHG 2014 : Particle SNOW scavenging coefficient based on Kyro et al 2009,
292              ! the below-cloud scavenging (Snow efficiency)
293              ! parameter B (=wetb) from SPECIES file
294              wetscav= wetb(ks)*10**(bcls(1)+ (bcls(2)*(log10(dquer_m))**(-4))+(bcls(3)*(log10(dquer_m))**(-3))+ (bcls(4)* &
295                   (log10(dquer_m))**(-2))+ (bcls(5)*(log10(dquer_m))**(-1))+ bcls(6)* (prec(1))**(0.5))
296
297            endif
298
299!             write(*,*) 'bl-cloud, act_temp=',act_temp, ',prec=',prec(1),',wetscav=', wetscav, ', jpart=',jpart
300
301          endif !gas or particle
302        endif ! positive below-cloud scavenging parameters given in Species file
303      endif !end BELOW
304
305
306  !********************
307  ! IN CLOUD SCAVENGING
308      !******************************************************
309      if (clouds_v.lt.4) then ! In-cloud
310! NIK 13 may 2015: only do incloud if positive in-cloud scavenging parameters are given in species file
311        if (weta_in(ks).gt.0. .or. wetb_in(ks).gt.0.) then
312          inc_count=inc_count+1
313! if negative coefficients (turned off) set to zero for use in equation
314          if (weta_in(ks).lt.0.) weta_in(ks)=0.
315          if (wetb_in(ks).lt.0.) wetb_in(ks)=0.
316
317          !ZHG 2015 Cloud liquid & ice water (CLWC+CIWC) from ECMWF
318          if (readclouds) then                  !icloud_stats(ix,jy,4,n) has units kg/m2
319            cl =icloud_stats(ix,jy,4,n)*(grfraction(1)/cc)
320          else                                  !parameterize cloudwater m2/m3
321            !ZHG updated parameterization of cloud water to better reproduce the values coming from ECMWF
322            cl=1.6E-6*prec(1)**0.36
323          endif
324
325            !ZHG: Calculate the partition between liquid and water phase water.
326            if (act_temp .le. 253) then
327              liq_frac=0
328            else if (act_temp .ge. 273) then
329              liq_frac=1
330            else
331              liq_frac =((act_temp-273)/(273-253))**2
332            endif
333           ! ZHG: Calculate the aerosol partition based on cloud phase and Ai and Bi
334            frac_act = liq_frac*weta_in(ks) +(1-liq_frac)*wetb_in(ks)
335 
336  !ZHG Use the activated fraction and the liqid water to calculate the washout
337
338  ! AEROSOL
339          !**************************************************
340          if (dquer(ks).gt. 0.) then ! is particle
341
342            S_i= frac_act/cl
343
344          !*********************
345  ! GAS
346          else ! is gas
347               
348            cle=(1-cl)/(henry(ks)*(r_air/3500.)*act_temp)+cl
349            !REPLACE to switch old/ new scheme
350            ! S_i=frac_act/cle
351            S_i=1/cle
352          endif ! gas or particle
353
354  ! scavenging coefficient based on Hertel et al 1995 - using the S_i for either gas or aerosol
355           !OLD
356          if (readclouds) then
357            wetscav=S_i*(prec(1)/3.6E6)
358          else
359            wetscav=S_i*(prec(1)/3.6E6)/clouds_h
360          endif
361
362!ZHG 2015 TEST
363!          Si_dummy=frac_act/2E-7*prec(1)**0.36
364!           wetscav_dummy=Si_dummy*(prec(1)/3.6E6)/clouds_h
365!           if (clouds_v.lt.4) then
366!           talltest=talltest+1
367!if(talltest .eq. 1) OPEN(UNIT=199, FILE=utfil,FORM='FORMATTED',STATUS = 'UNKNOWN')
368!if(talltest .lt. 100001)  write(199,*) prec(1)/3.6E6, cl, clouds_h*2E-7*prec(1)**0.36,clouds_v,ytra1(jpart)-90
369!if(talltest .lt. 100001)  write(199,*) wetscav, wetscav_dummy,prec(1),ytra1(jpart)-90,clouds_v,cl
370!if(talltest .eq. 100001) CLOSE(199)
371!if(talltest .eq. 100001) STOP
372!
373!write(*,*)  'PREC kg/m2s CLOUD kg/m2', (prec(1)/3.6E6), cl !, '2E-7*prec(1)**0.36',  2E-7*prec(1)**0.36,'2E-7*prec(1)**0.36*clouds_h',2E-7*prec(1)**0.36*clouds_h
374!write(*,*)  'PREC kg/m2s LSP+convp kg/m2', prec(1), convp+lsp
375!write(*,*)  wetscav, wetscav_dummy
376!write(*,*) cc, grfraction(1), cc/grfraction(1)
377
378!write(*,*)  'Lmbda_old', (prec(1)/3.6E6)/(clouds_h*2E-7*prec(1)**0.36)
379
380
381!write(*,*) '**************************************************'
382!write(*,*)  'clouds_h', clouds_h, 'clouds_v',clouds_v,'abs(ltsample)', abs(ltsample)
383!write(*,*)  'readclouds', readclouds, 'wetscav',wetscav, 'wetscav_dummy', wetscav_dummy
384!write(*,*)  'S_i', S_i , 'Si_dummy', Si_dummy, 'prec(1)', prec(1)
385
386
387!           write(*,*) 'PRECIPITATION ,cl  ECMWF , cl PARAMETIZED, clouds_v, lat' &
388!                      ,prec(1)/3.6E6, cl, clouds_h*2E-7*prec(1)**0.36,clouds_v,ytra1(jpart)-90
389
390!endif
391
392        endif ! positive in-cloud scavenging parameters given in Species file
393      endif !incloud
394!END ZHG TEST
395     
396  !**************************************************
397  ! CALCULATE DEPOSITION
398  !**************************************************
399  !     if (wetscav.le.0) write (*,*) 'neg, or 0 wetscav!'
400  !     +          ,wetscav,cle,cl,act_temp,prec,clouds_h,clouds_v
401
402      if (wetscav.gt.0.) then
403        wetdeposit(ks)=xmass1(jpart,ks)* &
404             (1.-exp(-wetscav*abs(ltsample)))*grfraction(1)  ! wet deposition
405!write(*,*) 'MASS DEPOSITED: PREC, WETSCAV, WETSCAVP', prec(1), wetdeposit(ks), xmass1(jpart,ks)* &
406!             (1.-exp(-wetscav_dummy*abs(ltsample)))*grfraction(1), clouds_v
407
408
409      else ! if no scavenging
410        wetdeposit(ks)=0.
411      endif
412
413      restmass = xmass1(jpart,ks)-wetdeposit(ks)
414      if (ioutputforeachrelease.eq.1) then
415        kp=npoint(jpart)
416      else
417        kp=1
418      endif
419      if (restmass .gt. smallnum) then
420        xmass1(jpart,ks)=restmass
421  !   depostatistic
422  !   wetdepo_sum(ks,kp)=wetdepo_sum(ks,kp)+wetdeposit(ks)
423  !   depostatistic
424      else
425        xmass1(jpart,ks)=0.
426      endif
427  !   Correct deposited mass to the last time step when radioactive decay of
428  !   gridded deposited mass was calculated
429      if (decay(ks).gt.0.) then
430        wetdeposit(ks)=wetdeposit(ks)*exp(abs(ldeltat)*decay(ks))
431      endif
432
433
434    end do !all species
435
436  ! Sabine Eckhardt, June 2008 create deposition runs only for forward runs
437  ! Add the wet deposition to accumulated amount on output grid and nested output grid
438  !*****************************************************************************
439
440    if (ldirect.eq.1) then
441      call wetdepokernel(nclass(jpart),wetdeposit,real(xtra1(jpart)), &
442           real(ytra1(jpart)),nage,kp)
443      if (nested_output.eq.1) call wetdepokernel_nest(nclass(jpart), &
444           wetdeposit,real(xtra1(jpart)),real(ytra1(jpart)),nage,kp)
445    endif
446
44720  continue
448  end do ! all particles
449
450  ! count the total number of below-cloud and in-cloud occurences:
451  tot_blc_count=tot_blc_count+blc_count
452  tot_inc_count=tot_inc_count+inc_count
453
454end subroutine wetdepo
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG