source: branches/jerome/src_flexwrf_v3.1/wetdepo.f90 @ 16

Last change on this file since 16 was 16, checked in by jebri, 11 years ago

sources for flexwrf v3.1

File size: 15.9 KB
Line 
1!***********************************************************************
2!* Copyright 2012,2013                                                *
3!* Jerome Brioude, Delia Arnold, Andreas Stohl, Wayne Angevine,       *
4!* John Burkhart, Massimo Cassiani, Adam Dingwell, Richard C Easter, Sabine Eckhardt,*
5!* Stephanie Evan, Jerome D Fast, Don Morton, Ignacio Pisso,          *
6!* Petra Seibert, Gerard Wotawa, Caroline Forster, Harald Sodemann,   *
7!*                                                                     *
8!* This file is part of FLEXPART WRF           
9!                                                                     *
10! FLEXPART is free software: you can redistribute it and/or modify    *
11! it under the terms of the GNU General Public License as published by*
12! the Free Software Foundation, either version 3 of the License, or   *
13! (at your option) any later version.                                 *
14!                                                                     *
15! FLEXPART is distributed in the hope that it will be useful,         *
16! but WITHOUT ANY WARRANTY; without even the implied warranty of      *
17! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
18! GNU General Public License for more details.                        *
19!                                                                     *
20! You should have received a copy of the GNU General Public License   *
21! along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
22!**********************************************************************
23
24subroutine wetdepo(itime,ltsample,loutnext)
25  !                     i      i        i
26  !*****************************************************************************
27  !                                                                            *
28  ! Calculation of wet deposition using the concept of scavenging coefficients.*
29  ! For lack of detailed information, washout and rainout are jointly treated. *
30  ! It is assumed that precipitation does not occur uniformly within the whole *
31  ! grid cell, but that only a fraction of the grid cell experiences rainfall. *
32  ! This fraction is parameterized from total cloud cover and rates of large   *
33  ! scale and convective precipitation.                                        *
34  !                                                                            *
35  !    Author: A. Stohl                                                        *
36  !                                                                            *
37  !    1 December 1996                                                         *
38  !                                                                            *
39  ! Correction by Petra Seibert, Sept 2002:                                    *
40  ! use centred precipitation data for integration                             *
41  ! Code may not be correct for decay of deposition!                           *
42  !                                                                            *
43  !* D Arnold 2012 01 25 implement patch developed by P.Seibert to avoid       *
44  !* problems with cloud fraction                                              *
45  !* set cloud cover to 1  when there is precip & use diagnosed clouds in      *
46  !* vertransform* for in-cloud below-cloud scavenging                         *
47  !* Set fraction to max val from grid cells larger than 5 km.                 *
48  !* Missing: new fractions of grid cell affected by rain , coming from        *
49  !* modelling studies with 150 km grid-sizes!                                 *
50  !* precip interpolation is under check by developers                         *
51  !* modifications identified by "CDA"                                         *
52!*
53!* Petra Seibert, 2011/2012: Fixing some deficiencies in this modification
54  !********************************************************************************
55  !*****************************************************************************
56  !                                                                            *
57  ! Variables:                                                                 *
58  ! cc [0-1]           total cloud cover                                       *
59  ! convp [mm/h]       convective precipitation rate                           *
60  ! grfraction [0-1]   fraction of grid, for which precipitation occurs        *
61  ! ix,jy              indices of output grid cell for each particle           *
62  ! itime [s]          actual simulation time [s]                              *
63  ! jpart              particle index                                          *
64  ! ldeltat [s]        interval since radioactive decay was computed           *
65  ! lfr, cfr           area fraction covered by precipitation for large scale  *
66  !                    and convective precipitation (dependent on prec. rate)  *
67  ! loutnext [s]       time for which gridded deposition is next output        *
68  ! loutstep [s]       interval at which gridded deposition is output          *
69  ! lsp [mm/h]         large scale precipitation rate                          *
70  ! ltsample [s]       interval over which mass is deposited                   *
71  ! prec [mm/h]        precipitation rate in subgrid, where precipitation occurs*
72  ! wetdeposit         mass that is wet deposited                              *
73  ! wetgrid            accumulated deposited mass on output grid               *
74  ! wetscav            scavenging coefficient                                  *
75  !                                                                            *
76  ! Constants:                                                                 *
77  !                                                                            *
78  !*****************************************************************************
79
80  use point_mod
81  use par_mod
82  use com_mod
83! use ieee_arithmetic
84  implicit none
85
86  integer :: jpart,itime,ltsample,loutnext,ldeltat,i,j,ix,jy
87!  integer :: ngrid,itage,nage,hz,il,interp_time, n, clouds_v
88  integer :: ngrid,itage,nage,kz,il,interp_time, n
89!  integer :: ks, kp
90  integer :: ks, kp, icbot,ictop, indcloud
91
92  real :: S_i, act_temp, cl, cle ! in cloud scavenging
93!  real :: clouds_h ! cloud height for the specific grid point
94
95!  real :: xtn,ytn,lsp,convp,cc,grfraction,prec,wetscav
96  real :: xtn,ytn,lsp,convp,cc,fraction,prec,wetscav,precsub,f
97  real :: wetscavold
98!  real :: wetdeposit(maxspec),restmass
99  real :: wetdeposit(maxspec),restmass
100
101
102  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
103  save lfr,cfr
104
105
106  real :: lfr(5) = (/ 0.5,0.65,0.8,0.9,0.95/)
107  real :: cfr(5) = (/ 0.4,0.55,0.7,0.8,0.9 /)
108
109  ! Compute interval since radioactive decay of deposited mass was computed
110  !************************************************************************
111
112  if (itime.le.loutnext) then
113    ldeltat=itime-(loutnext-loutstep)
114  else                                  ! first half of next interval
115    ldeltat=itime-loutnext
116  endif
117
118
119  ! Loop over all particles
120  !************************
121
122  do jpart=1,numpart
123    if (itra1(jpart).eq.-999999999) goto 20
124    if(ldirect.eq.1)then
125      if (itra1(jpart).gt.itime) goto 20
126    else
127      if (itra1(jpart).lt.itime) goto 20
128    endif
129  ! Determine age class of the particle
130    itage=abs(itra1(jpart)-itramem(jpart))
131    do nage=1,nageclass
132      if (itage.lt.lage(nage)) goto 33
133    end do
13433   continue
135
136
137  ! Determine which nesting level to be used
138  !*****************************************
139
140    ngrid=0
141    do j=numbnests,1,-1
142      if ((xtra1(jpart).gt.xln(j)).and.(xtra1(jpart).lt.xrn(j)).and. &
143           (ytra1(jpart).gt.yln(j)).and.(ytra1(jpart).lt.yrn(j))) then
144        ngrid=j
145        goto 23
146      endif
147    end do
14823   continue
149
150
151  ! Determine nested grid coordinates
152  !**********************************
153
154    if (ngrid.gt.0) then
155      xtn=(xtra1(jpart)-xln(ngrid))*xresoln(ngrid)
156      ytn=(ytra1(jpart)-yln(ngrid))*yresoln(ngrid)
157      ix=int(xtn)
158      jy=int(ytn)
159    else
160      ix=int(xtra1(jpart))
161      jy=int(ytra1(jpart))
162    endif
163
164
165  ! Interpolate large scale precipitation, convective precipitation and
166  ! total cloud cover
167  ! Note that interpolated time refers to itime-0.5*ltsample [PS]
168  !********************************************************************
169    interp_time=nint(itime-0.5*ltsample)
170!   if (jpart.eq.103) print*,real(xtra1(jpart)),real(ytra1(jpart)),jpart
171
172! CDA part of the bug fix from P. Seiber includes modification of
173! CDA interpol_rain routines.
174!    if (ieee_is_nan(xtra1(jpart))) then
175      if (xtra1(jpart).ne.xtra1(jpart)) then
176     print*,jpart,xtra1(jpart),ytra1(jpart),itra1(jpart)
177       endif
178    if (ngrid.eq.0) then
179!CDA old code:
180!      call interpol_rain(lsprec,convprec,tcc,nxmax,nymax, &
181!           1,nx,ny,memind,real(xtra1(jpart)),real(ytra1(jpart)),1, &
182!           memtime(1),memtime(2),interp_time,lsp,convp,cc)
183!CDA new code:
184       
185          call interpol_rain(lsprec,convprec,tcc, &
186           icloudbot,icloudthck,nxmax,nymax,1,nx,ny, &
187!          memind,sngl(xtra1(jpart)),sngl(ytra1(jpart)),1,memtime(1), &
188           memind,real(xtra1(jpart)),real(ytra1(jpart)),1,memtime(1), &
189           memtime(2),interp_time,lsp,convp,cc,icbot,ictop,icmv)
190
191    else
192! CDA old code:
193!      call interpol_rain_nests(lsprecn,convprecn,tccn, &
194!           nxmaxn,nymaxn,1,maxnests,ngrid,nxn,nyn,memind,xtn,ytn,1, &
195!           memtime(1),memtime(2),interp_time,lsp,convp,cc)
196! CDA new code:
197      call interpol_rain_nests(lsprecn,convprecn,tccn, &
198           icloudbotn,icloudthckn,nxmaxn,nymaxn,1,&
199             maxnests,ngrid,nxn,nyn, memind,xtn,ytn, &
200           1,memtime(1), &
201           memtime(2),interp_time,lsp,convp,cc,icbot,ictop,icmv)
202
203    endif
204
205!    if ((lsp.lt.0.01).and.(convp.lt.0.01)) goto 20
206
207!CPS 2012: subtract a small value, eg 0.01 mm/h, to remove spurious precip
208        prec = lsp+convp
209        precsub = 0.01
210        if (prec .lt. precsub) then
211          goto 20
212        else
213          f = (prec-precsub)/prec
214          lsp = f*lsp
215          convp = f*convp
216        endif
217
218
219  ! get the level were the actual particle is in
220      do il=2,nz
221        if (height(il).gt.ztra1(jpart)) then
222          kz=il-1
223          goto 26
224        endif
225      end do
22626     continue
227
228  n=memind(2)
229  if (abs(memtime(1)-interp_time).lt.abs(memtime(2)-interp_time)) &
230       n=memind(1)
231
232  ! if there is no precipitation or the particle is above the clouds no
233  ! scavenging is done
234! CDA  PS part of fix
235
236        if (ztra1(jpart) .le. float(ictop)) then
237          if (ztra1(jpart) .gt. float(icbot)) then
238            indcloud = 2 ! in-cloud
239          else
240            indcloud = 1 ! below-cloud
241          endif
242        elseif (ictop .eq. icmv) then
243          indcloud = 0 ! no cloud found, use old scheme
244        else
245          goto 20 ! above cloud
246        endif
247
248!  if (ngrid.eq.0) then
249!     clouds_v=clouds(ix,jy,hz,n)
250!     clouds_h=cloudsh(ix,jy,n)
251!  else
252!     clouds_v=cloudsn(ix,jy,hz,n,ngrid)
253!     clouds_h=cloudsnh(ix,jy,n,ngrid)
254!  endif
255  !write(*,*) 'there is
256  !    + precipitation',(clouds(ix,jy,ihz,n),ihz=1,20),lsp,convp,hz
257!  if (clouds_v.le.1) goto 20
258  !write (*,*) 'there is scavenging'
259
260  ! 1) Parameterization of the the area fraction of the grid cell where the
261  !    precipitation occurs: the absolute limit is the total cloud cover, but
262  !    for low precipitation rates, an even smaller fraction of the grid cell
263  !    is used. Large scale precipitation occurs over larger areas than
264  !    convective precipitation.
265  !**************************************************************************
266
267    if (lsp.gt.20.) then
268      i=5
269    else if (lsp.gt.8.) then
270      i=4
271    else if (lsp.gt.3.) then
272      i=3
273    else if (lsp.gt.1.) then
274      i=2
275    else
276      i=1
277    endif
278
279    if (convp.gt.20.) then
280      j=5
281    else if (convp.gt.8.) then
282      j=4
283    else if (convp.gt.3.) then
284      j=3
285    else if (convp.gt.1.) then
286      j=2
287    else
288      j=1
289    endif
290
291  !CDA
292     !  if (dx.le.5000.) then
293     !    j=5
294     !    i=5
295     !  endif
296   ! for the moment i set all to 5
297       i=5
298       j=5
299  !CDA
300  !CDA cc (coming from CLDFRA works as a mask for the precipitation
301  ! in vertransform* we diagnose the clouds. We consider not that
302  ! when there is precip, there is cloud as diagnosed in vertransform
303  ! this is the same approx used in other models HYSPLIT, for insance
304  ! set cc to 1
305    cc=1.
306
307    fraction=max(0.05,cc*(lsp*lfr(i)+convp*cfr(j))/(lsp+convp))
308
309  ! 2) Computation of precipitation rate in sub-grid cell
310  !******************************************************
311
312    prec=(lsp+convp)/fraction
313
314  ! 3) Computation of scavenging coefficients for all species
315  !    Computation of wet deposition
316  !**********************************************************
317
318    do 10 ks=1,nspec                                  ! loop over species
319      wetdeposit(ks)=0.
320
321          if (weta(ks).gt.0.) then
322            if (indcloud .eq. 1) then ! BELOW CLOUD SCAVENGING
323!C               for aerosols and not highliy soluble substances weta=5E-6
324              wetscav=weta(ks)*prec**wetb(ks)                ! scavenging coeff
325!c             write(*,*) 'bel. wetscav: ',wetscav
326            elseif (indcloud .eq. 2) then !  IN CLOUD SCAVENGING
327              if (ngrid.gt.0) then
328                 act_temp=ttn(ix,jy,kz,n,ngrid)
329              else
330                 act_temp=tt(ix,jy,kz,n)
331              endif
332              cl=2E-7*prec**0.36
333              if (dquer(ks).gt.0) then ! is particle
334                S_i=0.9/cl
335              else ! is gas
336                cle=(1-cl)/(henry(ks)*(r_air/3500.)*act_temp)+cl
337                S_i=1/cle
338              endif
339              wetscav=S_i*prec/3.6E6/(ictop-icbot) ! 3.6e6 converts mm/h to m/s
340              wetscavold = 2.e-5*prec**0.8 ! 1.e-4*prec**0.62
341              wetscav = min(0.1*wetscav,wetscavold) ! 0.1 is ASt current setting
342!PS tihb new scheme
343!PS tihc for no cloud, 1.e-5*prec**0.8 ! 1.e-4*prec**0.62
344!PS tihd limit to wetscavold= 1.e-4*prec**0.8
345!PS tihe mulitiply wetscav with .1
346!PS tihf set the wetscavold to 1.e-5prec**0.8
347!PS tihg set the wetscav & wetscavold to 2.e-5*prec**0.8 + weta=2.e-6 (instead of 1.e-6)
348
349            else ! PS: no cloud diagnosed, old scheme,
350!CPS          using with fixed a,b for simplicity, one may wish to change!!
351!             wetscav = 1.e-4*prec**0.62
352              wetscav = 2.e-5*prec**0.8
353            endif
354
355            wetdeposit(ks)=xmass1(jpart,ks)*  &
356              (1.-exp(-wetscav*abs(ltsample)))*fraction  ! wet deposition
357! CDA test
358!        if (wetdeposit(ks).gt.0) then
359!           write(*,*) 'wetdepo: ',wetdeposit(ks),ks
360!        endif
361
362
363            restmass = xmass1(jpart,ks)-wetdeposit(ks)
364            if (ioutputforeachrelease.eq.1) then
365              kp=npoint(jpart)
366            else
367              kp=1
368            endif
369            if (restmass .gt. smallnum) then
370              xmass1(jpart,ks)=restmass
371!cccccccccccccccc depostatistic
372!c            wetdepo_sum(ks,kp)=wetdepo_sum(ks,kp)+wetdeposit(ks)
373!cccccccccccccccc depostatistic
374            else
375              xmass1(jpart,ks)=0.
376            endif
377!C Correct deposited mass to the last time step when radioactive decay of
378!C gridded deposited mass was calculated
379            if (decay(ks).gt.0.) then
380              wetdeposit(ks)=wetdeposit(ks)*exp(abs(ldeltat)*decay(ks))
381            endif
382          else  ! weta(k)<0
383             wetdeposit(ks)=0.
384          endif
38510      continue
386
387
388
389  ! Sabine Eckhard, June 2008 create deposition runs only for forward runs
390  ! Add the wet deposition to accumulated amount on output grid and nested output grid
391  !*****************************************************************************
392
393    if (ldirect.eq.1) then
394!   print*,'kp',kp,jpart,npoint(jpart)
395  ! CDA added itage for the kernel not to be applied during the first hours
396    call wetdepokernel(nclass(jpart),wetdeposit,real(xtra1(jpart)), &
397         real(ytra1(jpart)),itage,nage,kp)
398    if (nested_output.eq.1) call wetdepokernel_nest(nclass(jpart), &
399         wetdeposit,real(xtra1(jpart)),real(ytra1(jpart)),itage, &
400         nage,kp)
401    endif
402
40320  continue
404  end do
405
406end subroutine wetdepo
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG