source: flexpart.git/src/wetdepo.f90 @ 54cbd6c

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

all f90 files for dry/wet backward mode

  • Property mode set to 100644
File size: 17.4 KB
RevLine 
[e200b7a]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
[54cbd6c]22subroutine wetdepo(itime,ltsample,loutnext,forreceptor)
23!                  i      i        i            i
[db712a8]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.                                        *
[54cbd6c]32! SEC: if forrecptor is true then the wetdeposition fraction is only applied *
33! on the xscav_frac and not on the xmass                                     *
[db712a8]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!*****************************************************************************
44!                                                                            *
45! Variables:                                                                 *
46! cc [0-1]           total cloud cover                                       *
47! convp [mm/h]       convective precipitation rate                           *
48! grfraction [0-1]   fraction of grid, for which precipitation occurs        *
49! ix,jy              indices of output grid cell for each particle           *
50! itime [s]          actual simulation time [s]                              *
51! jpart              particle index                                          *
52! ldeltat [s]        interval since radioactive decay was computed           *
53! lfr, cfr           area fraction covered by precipitation for large scale  *
54!                    and convective precipitation (dependent on prec. rate)  *
55! loutnext [s]       time for which gridded deposition is next output        *
56! loutstep [s]       interval at which gridded deposition is output          *
57! lsp [mm/h]         large scale precipitation rate                          *
58! ltsample [s]       interval over which mass is deposited                   *
59! prec [mm/h]        precipitation rate in subgrid, where precipitation occurs*
60! wetdeposit         mass that is wet deposited                              *
61! wetgrid            accumulated deposited mass on output grid               *
62! wetscav            scavenging coefficient                                  *
63!                                                                            *
64! Constants:                                                                 *
65!                                                                            *
66!*****************************************************************************
[e200b7a]67
68  use point_mod
69  use par_mod
70  use com_mod
71
72  implicit none
73
74  integer :: jpart,itime,ltsample,loutnext,ldeltat,i,j,ix,jy
[db712a8]75  integer :: ngrid,itage,nage,hz,il,interp_time, n
76  integer(kind=1) :: clouds_v
[4fbe7a5]77  integer :: ks, kp
[f75967d]78!  integer :: n1,n2, icbot,ictop, indcloud !TEST
[e200b7a]79  real :: S_i, act_temp, cl, cle ! in cloud scavenging
80  real :: clouds_h ! cloud height for the specific grid point
[8a65cb0]81  real :: xtn,ytn,lsp,convp,cc,grfraction(3),prec(3),wetscav,totprec
[e200b7a]82  real :: wetdeposit(maxspec),restmass
83  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
[db712a8]84!save lfr,cfr
[e200b7a]85
[5f9d14a]86  real, parameter :: lfr(5) = (/ 0.5,0.65,0.8,0.9,0.95/)
87  real, parameter :: cfr(5) = (/ 0.4,0.55,0.7,0.8,0.9 /)
[e200b7a]88
[8a65cb0]89!ZHG aerosol below-cloud scavenging removal polynomial constants for rain and snow
[5f9d14a]90  real, parameter :: bclr(6) = (/274.35758, 332839.59273, 226656.57259, 58005.91340, 6588.38582, 0.244984/) !rain (Laakso et al 2003)
91  real, parameter :: bcls(6) = (/22.7, 0.0, 0.0, 1321.0, 381.0, 0.0/) !now (Kyro et al 2009)
[8a65cb0]92  real :: frac_act, liq_frac, dquer_m
93
94  integer :: blc_count, inc_count
[d6a0977]95  real    :: Si_dummy, wetscav_dummy
[54cbd6c]96  logical :: readclouds_this_nest,forreceptor
[8a65cb0]97
98
[db712a8]99! Compute interval since radioactive decay of deposited mass was computed
100!************************************************************************
[e200b7a]101
102  if (itime.le.loutnext) then
103    ldeltat=itime-(loutnext-loutstep)
104  else                                  ! first half of next interval
105    ldeltat=itime-loutnext
106  endif
107
[db712a8]108! Loop over all particles
109!************************
[e200b7a]110
[8a65cb0]111  blc_count=0
112  inc_count=0
113
[e200b7a]114  do jpart=1,numpart
[4fbe7a5]115
[e200b7a]116    if (itra1(jpart).eq.-999999999) goto 20
117    if(ldirect.eq.1)then
118      if (itra1(jpart).gt.itime) goto 20
119    else
120      if (itra1(jpart).lt.itime) goto 20
121    endif
[8a65cb0]122
[db712a8]123! Determine age class of the particle
[e200b7a]124    itage=abs(itra1(jpart)-itramem(jpart))
125    do nage=1,nageclass
126      if (itage.lt.lage(nage)) goto 33
127    end do
[4fbe7a5]12833  continue
[e200b7a]129
130
[db712a8]131! Determine which nesting level to be used
132!*****************************************
[e200b7a]133
134    ngrid=0
135    do j=numbnests,1,-1
136      if ((xtra1(jpart).gt.xln(j)).and.(xtra1(jpart).lt.xrn(j)).and. &
[db712a8]137           (ytra1(jpart).gt.yln(j)).and.(ytra1(jpart).lt.yrn(j))) then
[e200b7a]138        ngrid=j
139        goto 23
140      endif
141    end do
[4fbe7a5]14223  continue
[e200b7a]143
144
[db712a8]145! Determine nested grid coordinates
146!**********************************
[1c3c778]147    readclouds_this_nest=.false.
[e200b7a]148
149    if (ngrid.gt.0) then
150      xtn=(xtra1(jpart)-xln(ngrid))*xresoln(ngrid)
151      ytn=(ytra1(jpart)-yln(ngrid))*yresoln(ngrid)
152      ix=int(xtn)
153      jy=int(ytn)
[1c3c778]154      if (readclouds_nest(ngrid)) readclouds_this_nest=.true.
[e200b7a]155    else
156      ix=int(xtra1(jpart))
157      jy=int(ytra1(jpart))
158    endif
159
160
[db712a8]161! Interpolate large scale precipitation, convective precipitation and
162! total cloud cover
163! Note that interpolated time refers to itime-0.5*ltsample [PS]
164!********************************************************************
[e200b7a]165    interp_time=nint(itime-0.5*ltsample)
[db712a8]166
[4fbe7a5]167    if (ngrid.eq.0) then
168      call interpol_rain(lsprec,convprec,tcc,nxmax,nymax, &
[db712a8]169           1,nx,ny,memind,real(xtra1(jpart)),real(ytra1(jpart)),1, &
170           memtime(1),memtime(2),interp_time,lsp,convp,cc)
[4fbe7a5]171    else
172      call interpol_rain_nests(lsprecn,convprecn,tccn, &
[db712a8]173           nxmaxn,nymaxn,1,maxnests,ngrid,nxn,nyn,memind,xtn,ytn,1, &
174           memtime(1),memtime(2),interp_time,lsp,convp,cc)
[4fbe7a5]175    endif
[e200b7a]176
[db712a8]177! get the level were the actual particle is in
[4fbe7a5]178    do il=2,nz
179      if (height(il).gt.ztra1(jpart)) then
180        hz=il-1
[db712a8]181!        goto 26
182        exit
[4fbe7a5]183      endif
184    end do
[db712a8]185!26  continue
[e200b7a]186
[4fbe7a5]187    n=memind(2)
188    if (abs(memtime(1)-interp_time).lt.abs(memtime(2)-interp_time)) &
[f75967d]189         n=memind(1)
[e200b7a]190
[4fbe7a5]191    if (ngrid.eq.0) then
192      clouds_v=clouds(ix,jy,hz,n)
193      clouds_h=cloudsh(ix,jy,n)
194    else
195      clouds_v=cloudsn(ix,jy,hz,n,ngrid)
[db712a8]196      clouds_h=cloudshn(ix,jy,n,ngrid)
[4fbe7a5]197    endif
[e200b7a]198
[db712a8]199! if there is no precipitation or the particle is above the clouds no
200! scavenging is done
[8a65cb0]201
[f75967d]202    if (clouds_v.le.1) goto 20
[db712a8]203
204! 1) Parameterization of the the area fraction of the grid cell where the
205!    precipitation occurs: the absolute limit is the total cloud cover, but
206!    for low precipitation rates, an even smaller fraction of the grid cell
207!    is used. Large scale precipitation occurs over larger areas than
208!    convective precipitation.
209!**************************************************************************
[e200b7a]210
211    if (lsp.gt.20.) then
212      i=5
213    else if (lsp.gt.8.) then
214      i=4
215    else if (lsp.gt.3.) then
216      i=3
217    else if (lsp.gt.1.) then
218      i=2
219    else
220      i=1
221    endif
222
223    if (convp.gt.20.) then
224      j=5
225    else if (convp.gt.8.) then
226      j=4
227    else if (convp.gt.3.) then
228      j=3
229    else if (convp.gt.1.) then
230      j=2
231    else
232      j=1
233    endif
234
[8a65cb0]235
[db712a8]236!ZHG oct 2014 : Calculated for 1) both 2) lsp 3) convp
237! Tentatively differentiate the grfraction for lsp and convp for treating differently the two forms
238! for now they are treated the same
[8a65cb0]239    grfraction(1)=max(0.05,cc*(lsp*lfr(i)+convp*cfr(j))/(lsp+convp))
240    grfraction(2)=max(0.05,cc*(lfr(i)))
241    grfraction(3)=max(0.05,cc*(cfr(j)))
242
[e200b7a]243
[db712a8]244! 2) Computation of precipitation rate in sub-grid cell
245!******************************************************
[8a65cb0]246    prec(1)=(lsp+convp)/grfraction(1)
247    prec(2)=(lsp)/grfraction(2)
248    prec(3)=(convp)/grfraction(3)
[e200b7a]249
250
[db712a8]251! 3) Computation of scavenging coefficients for all species
252!    Computation of wet deposition
253!**********************************************************
[e200b7a]254
[4fbe7a5]255    do ks=1,nspec      ! loop over species
256      wetdeposit(ks)=0.
257      wetscav=0.   
258
[5f9d14a]259      if (ngrid.gt.0) then
260        act_temp=ttn(ix,jy,hz,n,ngrid)
261      else
262        act_temp=tt(ix,jy,hz,n)
263      endif
[4fbe7a5]264
[db712a8]265
266!***********************
267! BELOW CLOUD SCAVENGING
268!*********************** 
[8a65cb0]269      if (clouds_v.ge.4) then !below cloud
270
[341f4b7]271! For gas: if positive below-cloud parameters (A or B), and dquer<=0
272!******************************************************************
273        if ((dquer(ks).le.0.).and.(weta_gas(ks).gt.0..or.wetb_gas(ks).gt.0.)) then
274          !        if (weta(ks).gt.0. .or. wetb(ks).gt.0.) then
[5f9d14a]275          blc_count=blc_count+1
[341f4b7]276          wetscav=weta_gas(ks)*prec(1)**wetb_gas(ks)
[8a65cb0]277
[341f4b7]278! For aerosols: if positive below-cloud parameters (Crain/Csnow or B), and dquer>0
279!*********************************************************************************
280        else if ((dquer(ks).gt.0.).and.(crain_aero(ks).gt.0..or.csnow_aero(ks).gt.0.)) then
281          blc_count=blc_count+1
[db712a8]282
[d6a0977]283!NIK 17.02.2015
[341f4b7]284! For the calculation here particle size needs to be in meter and not um as dquer is
285! changed in readreleases
286! For particles larger than 10 um use the largest size defined in the parameterizations (10um)
287          dquer_m=min(10.,dquer(ks))/1000000. !conversion from um to m
288
289! Rain:
290          if (act_temp .ge. 273. .and. crain_aero(ks).gt.0.)  then
291
[db712a8]292! ZHG 2014 : Particle RAIN scavenging coefficient based on Laakso et al 2003,
[341f4b7]293! the below-cloud scavenging (rain efficienty) parameter Crain (=crain_aero) from SPECIES file
294            wetscav=crain_aero(ks)*10**(bclr(1)+(bclr(2)*(log10(dquer_m))**(-4))+ &
295                 & (bclr(3)*(log10(dquer_m))**(-3))+ (bclr(4)*(log10(dquer_m))**(-2))+&
296                 &(bclr(5)*(log10(dquer_m))**(-1))+bclr(6)* (prec(1))**(0.5))
[8a65cb0]297
[341f4b7]298! Snow:
299          elseif (act_temp .lt. 273. .and. csnow_aero(ks).gt.0.)  then
[db712a8]300! ZHG 2014 : Particle SNOW scavenging coefficient based on Kyro et al 2009,
[341f4b7]301! the below-cloud scavenging (Snow efficiency) parameter Csnow (=csnow_aero) from SPECIES file
302            wetscav=csnow_aero(ks)*10**(bcls(1)+(bcls(2)*(log10(dquer_m))**(-4))+&
303                 &(bcls(3)*(log10(dquer_m))**(-3))+ (bcls(4)*(log10(dquer_m))**(-2))+&
304                 &(bcls(5)*(log10(dquer_m))**(-1))+ bcls(6)* (prec(1))**(0.5))
[8a65cb0]305
[341f4b7]306          endif
307         
[8a65cb0]308!             write(*,*) 'bl-cloud, act_temp=',act_temp, ',prec=',prec(1),',wetscav=', wetscav, ', jpart=',jpart
309
[341f4b7]310        endif ! gas or particle
311!      endif ! positive below-cloud scavenging parameters given in Species file
[d6a0977]312      endif !end BELOW
[4fbe7a5]313
[db712a8]314!********************
315! IN CLOUD SCAVENGING
316!********************
[d6a0977]317      if (clouds_v.lt.4) then ! In-cloud
[341f4b7]318! NIK 13 may 2015: only do incloud if positive in-cloud scavenging parameters are
319! given in species file, or if gas and positive Henry's constant
320        if ((ccn_aero(ks).gt.0. .or. in_aero(ks).gt.0.).or.(henry(ks).gt.0.and.dquer(ks).le.0)) then
[f75967d]321          inc_count=inc_count+1
[0f20c31]322! if negative coefficients (turned off) set to zero for use in equation
[341f4b7]323          if (ccn_aero(ks).lt.0.) ccn_aero(ks)=0.
324          if (in_aero(ks).lt.0.) in_aero(ks)=0.
[5f9d14a]325
[db712a8]326!ZHG 2015 Cloud liquid & ice water (CLWC+CIWC) from ECMWF
327! nested fields
328          if (ngrid.gt.0.and.readclouds_this_nest) then
[341f4b7]329            cl=ctwcn(ix,jy,n,ngrid)*(grfraction(1)/cc)
[db712a8]330          else if (ngrid.eq.0.and.readclouds) then
[341f4b7]331            cl=ctwc(ix,jy,n)*(grfraction(1)/cc)
[d6a0977]332          else                                  !parameterize cloudwater m2/m3
[db712a8]333!ZHG updated parameterization of cloud water to better reproduce the values coming from ECMWF
[d6a0977]334            cl=1.6E-6*prec(1)**0.36
[8a65cb0]335          endif
336
[db712a8]337!ZHG: Calculate the partition between liquid and water phase water.
338          if (act_temp .le. 253.) then
339            liq_frac=0
340          else if (act_temp .ge. 273.) then
341            liq_frac=1
342          else
343            liq_frac =((act_temp-273.)/(273.-253.))**2.
344          end if
345! ZHG: Calculate the aerosol partition based on cloud phase and Ai and Bi
[341f4b7]346          frac_act = liq_frac*ccn_aero(ks) +(1-liq_frac)*in_aero(ks)
[db712a8]347
348!ZHG Use the activated fraction and the liqid water to calculate the washout
[0f20c31]349
[db712a8]350! AEROSOL
[05cf28d]351!********
352          if (dquer(ks).gt.0.) then
[0f20c31]353            S_i= frac_act/cl
[8a65cb0]354
[db712a8]355! GAS
[05cf28d]356!****
357          else
[db712a8]358
[0f20c31]359            cle=(1-cl)/(henry(ks)*(r_air/3500.)*act_temp)+cl
[db712a8]360!REPLACE to switch old/ new scheme
[05cf28d]361          ! S_i=frac_act/cle
[d6a0977]362            S_i=1/cle
[0f20c31]363          endif ! gas or particle
[8a65cb0]364
[db712a8]365! scavenging coefficient based on Hertel et al 1995 - using the S_i for either gas or aerosol
366!OLD
367          if ((readclouds.and.ngrid.eq.0).or.(readclouds_this_nest.and.ngrid.gt.0)) then
[1c3c778]368            wetscav=incloud_ratio*S_i*(prec(1)/3.6E6)
[d6a0977]369          else
[1c3c778]370            wetscav=incloud_ratio*S_i*(prec(1)/3.6E6)/clouds_h
[d6a0977]371          endif
372
[8a65cb0]373
[0f20c31]374        endif ! positive in-cloud scavenging parameters given in Species file
[8a65cb0]375      endif !incloud
[d6a0977]376!END ZHG TEST
[db712a8]377
378!**************************************************
379! CALCULATE DEPOSITION
380!**************************************************
381!     if (wetscav.le.0) write (*,*) 'neg, or 0 wetscav!'
382!     +          ,wetscav,cle,cl,act_temp,prec,clouds_h,clouds_v
[e200b7a]383
[4fbe7a5]384      if (wetscav.gt.0.) then
[e200b7a]385        wetdeposit(ks)=xmass1(jpart,ks)* &
[5f9d14a]386             (1.-exp(-wetscav*abs(ltsample)))*grfraction(1)  ! wet deposition
[d6a0977]387!write(*,*) 'MASS DEPOSITED: PREC, WETSCAV, WETSCAVP', prec(1), wetdeposit(ks), xmass1(jpart,ks)* &
388!             (1.-exp(-wetscav_dummy*abs(ltsample)))*grfraction(1), clouds_v
389
390
[4fbe7a5]391      else ! if no scavenging
392        wetdeposit(ks)=0.
393      endif
[f13406c]394
[4fbe7a5]395      restmass = xmass1(jpart,ks)-wetdeposit(ks)
396      if (ioutputforeachrelease.eq.1) then
397        kp=npoint(jpart)
398      else
399        kp=1
400      endif
[54cbd6c]401      if (forreceptor .eqv. .false.) then
402         if (restmass .gt. smallnum) then
403           xmass1(jpart,ks)=restmass
[db712a8]404!   depostatistic
405!   wetdepo_sum(ks,kp)=wetdepo_sum(ks,kp)+wetdeposit(ks)
406!   depostatistic
[54cbd6c]407         else
408           xmass1(jpart,ks)=0.
409        endif
410     else ! for the backward deposition calculation
411         if (wetdeposit(ks).gt.0) then ! deposition occured
412                xscav_frac1(jpart,ks)=xscav_frac1(jpart,ks)*(-1.)* &
413                   wetdeposit(ks)/xmass1(jpart,ks)
414!                write (*,*) 'paricle kept: ',jpart,ks,wetdeposit(ks),xscav_frac1(jpart,ks)
415         else
416                xmass1(jpart,ks)=0.
417                xscav_frac1(jpart,ks)=0.
418         endif
419     endif
[db712a8]420!   Correct deposited mass to the last time step when radioactive decay of
421!   gridded deposited mass was calculated
[4fbe7a5]422      if (decay(ks).gt.0.) then
423        wetdeposit(ks)=wetdeposit(ks)*exp(abs(ldeltat)*decay(ks))
424      endif
[f13406c]425
[4fbe7a5]426    end do !all species
[e200b7a]427
[db712a8]428! Sabine Eckhardt, June 2008 create deposition runs only for forward runs
429! Add the wet deposition to accumulated amount on output grid and nested output grid
430!*****************************************************************************
[e200b7a]431
[4fbe7a5]432    if (ldirect.eq.1) then
433      call wetdepokernel(nclass(jpart),wetdeposit,real(xtra1(jpart)), &
[5f9d14a]434           real(ytra1(jpart)),nage,kp)
[4fbe7a5]435      if (nested_output.eq.1) call wetdepokernel_nest(nclass(jpart), &
[5f9d14a]436           wetdeposit,real(xtra1(jpart)),real(ytra1(jpart)),nage,kp)
[4fbe7a5]437    endif
[e200b7a]438
43920  continue
[4fbe7a5]440  end do ! all particles
[e200b7a]441
[db712a8]442! count the total number of below-cloud and in-cloud occurences:
[5f9d14a]443  tot_blc_count=tot_blc_count+blc_count
444  tot_inc_count=tot_inc_count+inc_count
[8a65cb0]445
[e200b7a]446end subroutine wetdepo
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG