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@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 6.2 KB
Line 
1subroutine wetdepo(itime,ltsample,loutnext)
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!*****************************************************************************
37
38  use point_mod
39  use par_mod
40  use com_mod
41
42  implicit none
43
44  integer :: jpart,itime,ltsample,loutnext,ldeltat
45  integer :: itage,nage
46  integer :: ks, kp
47  integer(selected_int_kind(16)), dimension(nspec) :: blc_count, inc_count
48  real :: grfraction(3),wetscav
49  real :: wetdeposit(maxspec),restmass
50  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
51
52! Compute interval since radioactive decay of deposited mass was computed
53!************************************************************************
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
61! Loop over all particles
62!************************
63
64  blc_count(:)=0
65  inc_count(:)=0
66
67  do jpart=1,numpart
68
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
75
76! Determine age class of the particle - nage is used for the kernel
77!******************************************************************
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
83
84    do ks=1,nspec      ! loop over species
85
86      if (WETDEPSPEC(ks).eqv..false.) cycle
87
88!**************************************************
89! CALCULATE DEPOSITION
90!**************************************************
91!       wetscav=0.
92       
93!        write(*,*) ks,dquer(ks), crain_aero(ks),csnow_aero(ks)
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
98
99      call get_wetscav(itime,ltsample,loutnext,jpart,ks,grfraction,inc_count,blc_count,wetscav)
100     
101
102      if (wetscav.gt.0.) then
103        wetdeposit(ks)=xmass1(jpart,ks)* &
104             (1.-exp(-wetscav*abs(ltsample)))*grfraction(1)  ! wet deposition
105      else ! if no scavenging
106        wetdeposit(ks)=0.
107      endif
108 
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
117!   depostatistic
118!   wetdepo_sum(ks,kp)=wetdepo_sum(ks,kp)+wetdeposit(ks)
119!   depostatistic
120      else
121        xmass1(jpart,ks)=0.
122      endif
123!   Correct deposited mass to the last time step when radioactive decay of
124!   gridded deposited mass was calculated
125      if (decay(ks).gt.0.) then
126        wetdeposit(ks)=wetdeposit(ks)*exp(abs(ldeltat)*decay(ks))
127      endif
128
129!    endif ! no deposition
130    end do ! loop over species
131
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!*****************************************************************************
135
136    if (ldirect.eq.1) then
137      call wetdepokernel(nclass(jpart),wetdeposit,real(xtra1(jpart)), &
138           real(ytra1(jpart)),nage,kp)
139      if (nested_output.eq.1) call wetdepokernel_nest(nclass(jpart), &
140           wetdeposit,real(xtra1(jpart)),real(ytra1(jpart)),nage,kp)
141    endif
142
14320  continue
144  end do ! all particles
145
146! count the total number of below-cloud and in-cloud occurences:
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)
149
150end subroutine wetdepo
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG