source: flexpart.git/src/wetdepo.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 5 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 6.2 KB
RevLine 
[e200b7a]1subroutine wetdepo(itime,ltsample,loutnext)
[db712a8]2!                  i      i        i
3!*****************************************************************************
4!                                                                            *
5! Calculation of wet deposition using the concept of scavenging coefficients.*
6! For lack of detailed information, washout and rainout are jointly treated. *
7! It is assumed that precipitation does not occur uniformly within the whole *
8! grid cell, but that only a fraction of the grid cell experiences rainfall. *
9! This fraction is parameterized from total cloud cover and rates of large   *
10! scale and convective precipitation.                                        *
11!                                                                            *
12!    Author: A. Stohl                                                        *
13!                                                                            *
14!    1 December 1996                                                         *
15!                                                                            *
16! Correction by Petra Seibert, Sept 2002:                                    *
17! use centred precipitation data for integration                             *
18! Code may not be correct for decay of deposition!                           *
19!                                                                            *
20!*****************************************************************************
21!                                                                            *
22! Variables:                                                                 *
23! ix,jy              indices of output grid cell for each particle           *
24! itime [s]          actual simulation time [s]                              *
25! jpart              particle index                                          *
26! ldeltat [s]        interval since radioactive decay was computed           *
27! loutnext [s]       time for which gridded deposition is next output        *
28! loutstep [s]       interval at which gridded deposition is output          *
29! ltsample [s]       interval over which mass is deposited                   *
30! wetdeposit         mass that is wet deposited                              *
31! wetgrid            accumulated deposited mass on output grid               *
32! wetscav            scavenging coefficient                                  *
33!                                                                            *
34! Constants:                                                                 *
35!                                                                            *
36!*****************************************************************************
[e200b7a]37
38  use point_mod
39  use par_mod
40  use com_mod
41
42  implicit none
43
[db375cc]44  integer :: jpart,itime,ltsample,loutnext,ldeltat
45  integer :: itage,nage
[4fbe7a5]46  integer :: ks, kp
[6985a98]47  integer(selected_int_kind(16)), dimension(nspec) :: blc_count, inc_count
[c9cf570]48  real :: grfraction(3),wetscav
[e200b7a]49  real :: wetdeposit(maxspec),restmass
50  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
[8a65cb0]51
[db712a8]52! Compute interval since radioactive decay of deposited mass was computed
53!************************************************************************
[e200b7a]54
55  if (itime.le.loutnext) then
56    ldeltat=itime-(loutnext-loutstep)
57  else                                  ! first half of next interval
58    ldeltat=itime-loutnext
59  endif
60
[db712a8]61! Loop over all particles
62!************************
[e200b7a]63
[c8fc724]64  blc_count(:)=0
65  inc_count(:)=0
[8a65cb0]66
[e200b7a]67  do jpart=1,numpart
[4fbe7a5]68
[e200b7a]69    if (itra1(jpart).eq.-999999999) goto 20
70    if(ldirect.eq.1)then
71      if (itra1(jpart).gt.itime) goto 20
72    else
73      if (itra1(jpart).lt.itime) goto 20
74    endif
[8a65cb0]75
[c9cf570]76! Determine age class of the particle - nage is used for the kernel
[341f4b7]77!******************************************************************
[c9cf570]78     itage=abs(itra1(jpart)-itramem(jpart))
79     do nage=1,nageclass
80       if (itage.lt.lage(nage)) goto 33
81     end do
82 33  continue
[d6a0977]83
[92a74b2]84    do ks=1,nspec      ! loop over species
[8a65cb0]85
[6985a98]86      if (WETDEPSPEC(ks).eqv..false.) cycle
[db712a8]87
88!**************************************************
89! CALCULATE DEPOSITION
90!**************************************************
[0539b8f]91!       wetscav=0.
[0e085a8]92       
[7919911]93!        write(*,*) ks,dquer(ks), crain_aero(ks),csnow_aero(ks)
[0e085a8]94!       if (((dquer(ks).le.0.).and.(weta_gas(ks).gt.0..or.wetb_gas(ks).gt.0.)) &
95!          .or. &
96!          ((dquer(ks).gt.0.).and.(crain_aero(ks).gt.0..or.csnow_aero(ks).gt.0.).or. &
97!            (ccn_aero(ks).gt0) .or. (in_aero(ks).gt.0) .or. (henry(ks).gt.0)))  then
[585a533]98
[92a74b2]99      call get_wetscav(itime,ltsample,loutnext,jpart,ks,grfraction,inc_count,blc_count,wetscav)
[0e085a8]100     
[e200b7a]101
[4fbe7a5]102      if (wetscav.gt.0.) then
[e200b7a]103        wetdeposit(ks)=xmass1(jpart,ks)* &
[5f9d14a]104             (1.-exp(-wetscav*abs(ltsample)))*grfraction(1)  ! wet deposition
[4fbe7a5]105      else ! if no scavenging
106        wetdeposit(ks)=0.
107      endif
[c9cf570]108 
[4fbe7a5]109      restmass = xmass1(jpart,ks)-wetdeposit(ks)
110      if (ioutputforeachrelease.eq.1) then
111        kp=npoint(jpart)
112      else
113        kp=1
114      endif
115      if (restmass .gt. smallnum) then
116        xmass1(jpart,ks)=restmass
[db712a8]117!   depostatistic
118!   wetdepo_sum(ks,kp)=wetdepo_sum(ks,kp)+wetdeposit(ks)
119!   depostatistic
[4fbe7a5]120      else
121        xmass1(jpart,ks)=0.
122      endif
[db712a8]123!   Correct deposited mass to the last time step when radioactive decay of
124!   gridded deposited mass was calculated
[4fbe7a5]125      if (decay(ks).gt.0.) then
126        wetdeposit(ks)=wetdeposit(ks)*exp(abs(ldeltat)*decay(ks))
127      endif
[f13406c]128
[7919911]129!    endif ! no deposition
[92a74b2]130    end do ! loop over species
[e200b7a]131
[db712a8]132! Sabine Eckhardt, June 2008 create deposition runs only for forward runs
133! Add the wet deposition to accumulated amount on output grid and nested output grid
134!*****************************************************************************
[e200b7a]135
[4fbe7a5]136    if (ldirect.eq.1) then
137      call wetdepokernel(nclass(jpart),wetdeposit,real(xtra1(jpart)), &
[5f9d14a]138           real(ytra1(jpart)),nage,kp)
[4fbe7a5]139      if (nested_output.eq.1) call wetdepokernel_nest(nclass(jpart), &
[5f9d14a]140           wetdeposit,real(xtra1(jpart)),real(ytra1(jpart)),nage,kp)
[4fbe7a5]141    endif
[e200b7a]142
14320  continue
[4fbe7a5]144  end do ! all particles
[e200b7a]145
[db712a8]146! count the total number of below-cloud and in-cloud occurences:
[c8fc724]147  tot_blc_count(1:nspec)=tot_blc_count(1:nspec)+blc_count(1:nspec)
148  tot_inc_count(1:nspec)=tot_inc_count(1:nspec)+inc_count(1:nspec)
[8a65cb0]149
[e200b7a]150end subroutine wetdepo
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG