source: branches/FLEXPART_9.1.3/src/wetdepo.f90-F9.1.2 @ 13

Last change on this file since 13 was 13, checked in by saeck, 11 years ago

update to wetdepo.f90

File size: 12.5 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  ! Modification by Sabine Eckhart to introduce a new in-/below-cloud
42  ! scheme, not dated
43  ! Petra Seibert, 2011/2012: Fixing some deficiencies in this modification
44  ! note that also other subroutines are affected by the fix
45  !*****************************************************************************
46  !                                                                            *
47  ! Variables:                                                                 *
48  ! cc [0-1]           total cloud cover                                       *
49  ! convp [mm/h]       convective precipitation rate                           *
50  ! grfraction [0-1]   fraction of grid, for which precipitation occurs        *
51  ! ix,jy              indices of output grid cell for each particle           *
52  ! itime [s]          actual simulation time [s]                              *
53  ! jpart              particle index                                          *
54  ! ldeltat [s]        interval since radioactive decay was computed           *
55  ! lfr, cfr           area fraction covered by precipitation for large scale  *
56  !                    and convective precipitation (dependent on prec. rate)  *
57  ! loutnext [s]       time for which gridded deposition is next output        *
58  ! loutstep [s]       interval at which gridded deposition is output          *
59  ! lsp [mm/h]         large scale precipitation rate                          *
60  ! ltsample [s]       interval over which mass is deposited                   *
61  ! prec [mm/h]        precipitation rate in subgrid, where precipitation occurs*
62  ! wetdeposit         mass that is wet deposited                              *
63  ! wetgrid            accumulated deposited mass on output grid               *
64  ! wetscav            scavenging coefficient                                  *
65  !                                                                            *
66  ! Constants:                                                                 *
67  !                                                                            *
68  !*****************************************************************************
69
70  use point_mod
71  use par_mod
72  use com_mod
73
74  implicit none
75
76  integer :: jpart,itime,ltsample,loutnext,ldeltat,i,j,ix,jy
77  integer :: ngrid,itage,nage,hz,il,interp_time, n, clouds_v
78  integer :: ks, kp
79  real :: S_i, act_temp, cl, cle ! in cloud scavenging
80  real :: clouds_h ! cloud height for the specific grid point
81  real :: xtn,ytn,lsp,convp,cc,grfraction,prec,wetscav
82  real :: wetdeposit(maxspec),restmass
83  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
84  save lfr,cfr
85
86
87  real :: lfr(5) = (/ 0.5,0.65,0.8,0.9,0.95/)
88  real :: cfr(5) = (/ 0.4,0.55,0.7,0.8,0.9 /)
89
90  ! Compute interval since radioactive decay of deposited mass was computed
91  !************************************************************************
92
93  if (itime.le.loutnext) then
94    ldeltat=itime-(loutnext-loutstep)
95  else                                  ! first half of next interval
96    ldeltat=itime-loutnext
97  endif
98
99
100  ! Loop over all particles
101  !************************
102
103  do jpart=1,numpart
104    if (itra1(jpart).eq.-999999999) goto 20
105    if(ldirect.eq.1)then
106      if (itra1(jpart).gt.itime) goto 20
107    else
108      if (itra1(jpart).lt.itime) goto 20
109    endif
110  ! Determine age class of the particle
111    itage=abs(itra1(jpart)-itramem(jpart))
112    do nage=1,nageclass
113      if (itage.lt.lage(nage)) goto 33
114    end do
11533   continue
116
117
118  ! Determine which nesting level to be used
119  !*****************************************
120
121    ngrid=0
122    do j=numbnests,1,-1
123      if ((xtra1(jpart).gt.xln(j)).and.(xtra1(jpart).lt.xrn(j)).and. &
124           (ytra1(jpart).gt.yln(j)).and.(ytra1(jpart).lt.yrn(j))) then
125        ngrid=j
126        goto 23
127      endif
128    end do
12923   continue
130
131
132  ! Determine nested grid coordinates
133  !**********************************
134
135    if (ngrid.gt.0) then
136      xtn=(xtra1(jpart)-xln(ngrid))*xresoln(ngrid)
137      ytn=(ytra1(jpart)-yln(ngrid))*yresoln(ngrid)
138      ix=int(xtn)
139      jy=int(ytn)
140    else
141      ix=int(xtra1(jpart))
142      jy=int(ytra1(jpart))
143    endif
144
145
146  ! Interpolate large scale precipitation, convective precipitation and
147  ! total cloud cover
148  ! Note that interpolated time refers to itime-0.5*ltsample [PS]
149  !********************************************************************
150    interp_time=nint(itime-0.5*ltsample)
151
152    if (ngrid.eq.0) then
153      call interpol_rain(lsprec,convprec,tcc,nxmax,nymax, &
154           1,nx,ny,memind,real(xtra1(jpart)),real(ytra1(jpart)),1, &
155           memtime(1),memtime(2),interp_time,lsp,convp,cc)
156    else
157      call interpol_rain_nests(lsprecn,convprecn,tccn, &
158           nxmaxn,nymaxn,1,maxnests,ngrid,nxn,nyn,memind,xtn,ytn,1, &
159           memtime(1),memtime(2),interp_time,lsp,convp,cc)
160    endif
161
162    if ((lsp.lt.0.01).and.(convp.lt.0.01)) goto 20
163
164  ! get the level were the actual particle is in
165      do il=2,nz
166        if (height(il).gt.ztra1(jpart)) then
167          hz=il-1
168          goto 26
169        endif
170      end do
17126     continue
172
173  n=memind(2)
174  if (abs(memtime(1)-interp_time).lt.abs(memtime(2)-interp_time)) &
175       n=memind(1)
176
177  ! if there is no precipitation or the particle is above the clouds no
178  ! scavenging is done
179
180  if (ngrid.eq.0) then
181     clouds_v=clouds(ix,jy,hz,n)
182     clouds_h=cloudsh(ix,jy,n)
183  else
184     clouds_v=cloudsn(ix,jy,hz,n,ngrid)
185     clouds_h=cloudsnh(ix,jy,n,ngrid)
186  endif
187  !write(*,*) 'there is
188  !    + precipitation',(clouds(ix,jy,ihz,n),ihz=1,20),lsp,convp,hz
189  if (clouds_v.le.1) goto 20
190  !write (*,*) 'there is scavenging'
191
192  ! 1) Parameterization of the the area fraction of the grid cell where the
193  !    precipitation occurs: the absolute limit is the total cloud cover, but
194  !    for low precipitation rates, an even smaller fraction of the grid cell
195  !    is used. Large scale precipitation occurs over larger areas than
196  !    convective precipitation.
197  !**************************************************************************
198
199    if (lsp.gt.20.) then
200      i=5
201    else if (lsp.gt.8.) then
202      i=4
203    else if (lsp.gt.3.) then
204      i=3
205    else if (lsp.gt.1.) then
206      i=2
207    else
208      i=1
209    endif
210
211    if (convp.gt.20.) then
212      j=5
213    else if (convp.gt.8.) then
214      j=4
215    else if (convp.gt.3.) then
216      j=3
217    else if (convp.gt.1.) then
218      j=2
219    else
220      j=1
221    endif
222
223    grfraction=max(0.05,cc*(lsp*lfr(i)+convp*cfr(j))/(lsp+convp))
224
225  ! 2) Computation of precipitation rate in sub-grid cell
226  !******************************************************
227
228    prec=(lsp+convp)/grfraction
229
230  ! 3) Computation of scavenging coefficients for all species
231  !    Computation of wet deposition
232  !**********************************************************
233
234    do ks=1,nspec                                  ! loop over species
235      wetdeposit(ks)=0.
236      if (weta(ks).gt.0.) then
237        if (clouds_v.ge.4) then
238  !          BELOW CLOUD SCAVENGING
239  !          for aerosols and not highliy soluble substances weta=5E-6
240          wetscav=weta(ks)*prec**wetb(ks)                ! scavenging coeff.
241  !        write(*,*) 'bel. wetscav: ',wetscav
242        else ! below_cloud clouds_v is lt 4 and gt 1 -> in cloud scavenging
243  !        IN CLOUD SCAVENGING
244  ! BUGFIX tt for nested fields should be ttn
245  ! sec may 2008
246          if (ngrid.gt.0) then
247             act_temp=ttn(ix,jy,hz,n,ngrid)
248          else
249             act_temp=tt(ix,jy,hz,n)
250          endif
251
252! NIK 31.01.2013: SPECIES defined parameters for the in-cloud scavening
253! weta_in=2.0E-07 (default)
254! wetb_in=0.36 (default)
255! wetc_in=0.9 (default)
256! wetd_in: Scaling factor for the total in-cloud scavenging (default 1.0-no scaling)
257          cl=weta_in(ks)*prec**wetb_in(ks)
258          if (dquer(ks).gt.0) then ! is particle
259            S_i=wetc_in(ks)/cl
260           else ! is gas
261            cle=(1-cl)/(henry(ks)*(r_air/3500.)*act_temp)+cl
262            S_i=1/cle
263           endif
264           wetscav=S_i*prec/3.6E6/clouds_h/wetd_in(ks)
265  !         write(*,*) 'in. wetscav:'
266  !    +          ,wetscav,cle,cl,act_temp,prec,clouds_h
267        endif
268
269
270  !      if (wetscav.le.0) write (*,*) 'neg, or 0 wetscav!'
271  !    +          ,wetscav,cle,cl,act_temp,prec,clouds_h,clouds_v
272        wetdeposit(ks)=xmass1(jpart,ks)* &
273             (1.-exp(-wetscav*abs(ltsample)))*grfraction  ! wet deposition
274  !      new particle mass:
275  !      if (wetdeposit(ks).gt.0) then
276  !         write(*,*) 'wetdepo: ',wetdeposit(ks),ks
277  !      endif
278        restmass = xmass1(jpart,ks)-wetdeposit(ks)
279         if (ioutputforeachrelease.eq.1) then
280           kp=npoint(jpart)
281         else
282           kp=1
283         endif
284        if (restmass .gt. smallnum) then
285          xmass1(jpart,ks)=restmass
286  !ccccccccccccccc depostatistic
287  !       wetdepo_sum(ks,kp)=wetdepo_sum(ks,kp)+wetdeposit(ks)
288  !ccccccccccccccc depostatistic
289        else
290          xmass1(jpart,ks)=0.
291        endif
292  ! Correct deposited mass to the last time step when radioactive decay of
293  ! gridded deposited mass was calculated
294       if (decay(ks).gt.0.) then
295          wetdeposit(ks)=wetdeposit(ks) &
296               *exp(abs(ldeltat)*decay(ks))
297       endif
298      else  ! weta(k)
299         wetdeposit(ks)=0.
300      endif ! weta(k)
301    end do
302
303  ! Sabine Eckhard, June 2008 create deposition runs only for forward runs
304  ! Add the wet deposition to accumulated amount on output grid and nested output grid
305  !*****************************************************************************
306
307    if (ldirect.eq.1) then
308    call wetdepokernel(nclass(jpart),wetdeposit,real(xtra1(jpart)), &
309         real(ytra1(jpart)),nage,kp)
310    if (nested_output.eq.1) call wetdepokernel_nest(nclass(jpart), &
311         wetdeposit,real(xtra1(jpart)),real(ytra1(jpart)), &
312         nage,kp)
313    endif
314
31520  continue
316  end do
317
318end subroutine wetdepo
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG