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