source: branches/petra/src/wetdepo.f90 @ 37

Last change on this file since 37 was 37, checked in by pesei, 9 years ago

Wet dep quick fix and other small changes. Wet depo quick fix not final yet.

File size: 12.9 KB
Line 
1!**********************************************************************
2! Copyright 1998-2015                                                 *
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  ! PS, 2/2015: implement wet depo quick fix from 2011/12 in this version
42  !  it implements interpolated and improved clouds
43  !  Also, certain deficiencies for thin clouds fixed also
44  !  Pass itage to wetdepokernel
45  !                                                                            *
46  !*****************************************************************************
47  !                                                                            *
48  ! Variables:                                                                 *
49  ! cc [0-1]           total cloud cover                                       *
50  ! convp [mm/h]       convective precipitation rate                           *
51  ! grfraction [0-1]   fraction of grid, for which precipitation occurs        *
52  ! ix,jy              indices of output grid cell for each particle           *
53  ! itime [s]          actual simulation time [s]                              *
54  ! jpart              particle index                                          *
55  ! ldeltat [s]        interval since radioactive decay was computed           *
56  ! lfr, cfr           area fraction covered by precipitation for large scale  *
57  !                    and convective precipitation (dependent on prec. rate)  *
58  ! loutnext [s]       time for which gridded deposition is next output        *
59  ! loutstep [s]       interval at which gridded deposition is output          *
60  ! lsp [mm/h]         large scale precipitation rate                          *
61  ! ltsample [s]       interval over which mass is deposited                   *
62  ! prec [mm/h]        precipitation rate in subgrid, where precipitation occurs*
63  ! wetdeposit         mass that is wet deposited                              *
64  ! wetgrid            accumulated deposited mass on output grid               *
65  ! wetscav            scavenging coefficient                                  *
66  !                                                                            *
67  ! Constants:                                                                 *
68  !                                                                            *
69  !*****************************************************************************
70
71  use point_mod
72  use par_mod
73  use com_mod
74
75  implicit none
76
77  integer :: jpart,itime,ltsample,loutnext,ldeltat,i,j,ix,jy
78  integer :: ngrid,itage,nage,kz,il,interp_time,n
79  integer :: ks, kp, n1,n2, icbot,ictop, indcloud
80  real :: S_i, act_temp, cl, cle ! in cloud scavenging
81  real :: xtn,ytn,lsp,convp,cc,grfraction,prec,wetscav,wetscavold,precsub,f
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  particle_loop: &
104  do jpart=1,numpart
105
106    if (itra1(jpart).eq.-999999999) goto 20
107    if (ldirect.eq.1) then
108      if (itra1(jpart).gt.itime) goto 20
109    else
110      if (itra1(jpart).lt.itime) goto 20
111    endif
112   
113  ! Determine age class of the particle
114    itage=abs(itra1(jpart)-itramem(jpart))
115    do nage=1,nageclass
116      if (itage.lt.lage(nage)) goto 33
117    end do
11833  continue
119
120
121  ! Determine which nesting level to be used
122  !*****************************************
123
124    ngrid=0
125    do j=numbnests,1,-1
126      if ((xtra1(jpart).gt.xln(j)).and.(xtra1(jpart).lt.xrn(j)).and. &
127      (ytra1(jpart).gt.yln(j)).and.(ytra1(jpart).lt.yrn(j))) then
128        ngrid=j
129        goto 23
130      endif
131    enddo
13223  continue
133
134
135  ! Determine nested grid coordinates
136  !**********************************
137
138    if (ngrid.gt.0) then
139      xtn=real(xtra1(jpart)-xln(ngrid))*xresoln(ngrid)
140      ytn=real(ytra1(jpart)-yln(ngrid))*yresoln(ngrid)
141      ix=int(xtn)
142      jy=int(ytn)
143    else
144      ix=int(xtra1(jpart))
145      jy=int(ytra1(jpart))
146    endif
147
148
149  ! Interpolate large scale precipitation, convective precipitation and
150  ! total cloud cover
151  ! Note that interpolated time refers to itime-0.5*ltsample [PS]
152  !********************************************************************
153    interp_time=nint(itime-0.5*ltsample)
154
155    if (ngrid.eq.0) then
156      call interpol_rain(lsprec,convprec,tcc,nxmax,nymax,1,nx,ny,memind,&
157      real(xtra1(jpart)),real(ytra1(jpart)), &
158      1,memtime(1),memtime(2),interp_time,lsp,convp,cc,icbot,ictop)
159    else
160      call interpol_rain_nests(lsprecn,convprecn,tccn, &
161      nxmaxn,nymaxn,1,maxnests,ngrid,nxn,nyn,memind,xtn,ytn, &
162      1,memtime(1),memtime(2),interp_time,lsp,convp,cc,icbot,ictop)
163    endif
164
165! PS 2012/2015: subtract a small value, eg 0.01 mm/h,
166! to remove spurious precip; replaces previous code
167        prec = lsp+convp
168        precsub = 0.01
169        if (prec .lt. precsub) then
170          goto 20
171        else
172          f = (prec-precsub)/prec
173          lsp = f*lsp
174          convp = f*convp
175        endif
176
177  ! get the level were the actual particle is in
178    do il=2,nz
179      if (height(il).gt.ztra1(jpart)) then
180        kz=il-1
181        goto 26
182      endif
183    enddo
18426  continue
185
186    n=memind(2)
187    if (abs(memtime(1)-interp_time).lt.abs(memtime(2)-interp_time)) &
188      n=memind(1)
189
190  ! if there is no precipitation or the particle is above the clouds no
191  ! scavenging is done
192  ! PS: part of 2011/2012/2015 fix, replaces previous code
193 
194    if (ztra1(jpart) .le. float(ictop)) then
195      if (ztra1(jpart) .gt. float(icbot)) then
196        indcloud = 2 ! in-cloud
197      else
198        indcloud = 1 ! below-cloud
199      endif
200    elseif (ictop .eq. icmv) then
201      indcloud = 0 ! no cloud found, use old scheme
202    else
203      goto 20 ! above cloud
204    endif
205
206
207  ! 1) Parameterization of the the area fraction of the grid cell where the
208  !    precipitation occurs: the absolute limit is the total cloud cover, but
209  !    for low precipitation rates, an even smaller fraction of the grid cell
210  !    is used. Large scale precipitation occurs over larger areas than
211  !    convective precipitation.
212  !**************************************************************************
213
214    if (lsp.gt.20.) then
215      i=5
216    else if (lsp.gt.8.) then
217      i=4
218    else if (lsp.gt.3.) then
219      i=3
220    else if (lsp.gt.1.) then
221      i=2
222    else
223      i=1
224    endif
225
226    if (convp.gt.20.) then
227      j=5
228    else if (convp.gt.8.) then
229      j=4
230    else if (convp.gt.3.) then
231      j=3
232    else if (convp.gt.1.) then
233      j=2
234    else
235      j=1
236    endif
237
238    grfraction=max(0.05,cc*(lsp*lfr(i)+convp*cfr(j))/(lsp+convp))
239
240  ! 2) Computation of precipitation rate in sub-grid cell
241  !******************************************************
242
243    prec=(lsp+convp)/grfraction
244
245  ! 3) Computation of scavenging coefficients for all species
246  !    Computation of wet deposition
247  !**********************************************************
248
249    species_loop: &
250    do ks=1,nspec      ! loop over species
251   
252      wetdeposit(ks)=0.
253      wetscav=0.   
254
255  ! NIK 09.12.13: allowed to turn off either below cloud or in-cloud by including TWO separate blocks of if tests
256      if (weta(ks).gt.0.) then
257        ! positive below-cloud coefficient from SPECIES file
258       
259        if (indcloud .eq. 1) then ! below-cloud scavenging
260
261         ! for aerosols and not highly soluble substances weta=5E-6
262           wetscav=weta(ks)*prec**wetb(ks)  ! scavenging coefficient
263
264        elseif (indcloud .eq. 2) then ! in-cloud scavenging
265
266          if (ngrid.gt.0) then
267            act_temp=ttn(ix,jy,kz,n,ngrid)
268          else
269            act_temp=tt(ix,jy,kz,n)
270          endif
271
272! NIK 31.01.2013: SPECIES defined parameters for the in-cloud scavening
273! weta_in=2.0E-07 (default)
274! wetb_in=0.36 (default)
275! wetc_in=0.9 (default)
276! wetd_in: Scaling factor for the total in-cloud scavenging (default 1.0 - no scaling)
277        cl=weta_in(ks)*prec**wetb_in(ks)
278        if (dquer(ks).gt.0) then ! is particle
279          S_i=wetc_in(ks)/cl
280        else ! is gas
281          cle=(1.-cl)/(henry(ks)*(r_air/3500.)*act_temp)+cl
282          S_i=1./cle
283        endif
284        wetscav=S_i*prec/3.6E6/wetd_in(ks)/(ictop-icbot) ! 3.6e6 converts mm/h to m/s
285! PS - prevent extrem values of wetscav:
286        wetscavold = 2.e-5*prec**0.8     
287        wetscav = min( 0.1*wetscav, wetscavold ) ! here we sacle the above wetscav. This is now twice, also in wetd_in.
288
289      else ! PS: no cloud diagnosed, old scheme,
290
291! PS    using with fixed a,b for simplicity, one may wish to change!!
292        wetscav = 2.e-5*prec**0.8 ! was before: 1.e-4*prec**0.62
293
294      endif ! end regime cases
295     
296  ! calculate deposition
297      wetdeposit(ks)=xmass1(jpart,ks)* &
298        (1.-exp(-wetscav*abs(ltsample)))*grfraction
299      restmass = xmass1(jpart,ks)-wetdeposit(ks)
300      if (ioutputforeachrelease.eq.1) then
301        kp=npoint(jpart)
302      else
303        kp=1
304      endif
305      if (restmass .gt. smallnum) then
306        xmass1(jpart,ks)=restmass
307    !   depostatistic:
308   !!   wetdepo_sum(ks,kp)=wetdepo_sum(ks,kp)+wetdeposit(ks)
309      else
310        xmass1(jpart,ks)=0.
311      endif
312  !   Correct deposited mass to the last time step when radioactive decay of
313  !   gridded deposited mass was calculated
314      if (decay(ks).gt.0.) &
315        wetdeposit(ks)=wetdeposit(ks)*exp(abs(ldeltat)*decay(ks))
316
317    else ! weta(k) <= 0
318
319      wetdeposit(ks)=0.
320   
321    endif
322
323  enddo species_loop
324
325 ! Sabine Eckhardt, June 2008: write deposition only in forward runs
326 ! add the wet deposition from this step to accumulated amount
327 ! on output grid and nested output grid
328
329    if (ldirect.eq.1) then
330      call wetdepokernel(nclass(jpart),wetdeposit, &
331        real(xtra1(jpart)),real(ytra1(jpart)),itage,nage,kp)
332      if (nested_output.eq.1) call wetdepokernel_nest(nclass(jpart), &
333        wetdeposit,real(xtra1(jpart)),real(ytra1(jpart)),itage,nage,kp)
334    endif
335
33620  continue ! jump here for particles not to be treated
337  enddo particle_loop
338
339end subroutine wetdepo
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG