Line data Source code
1 : ! Include shortname defintions, so that the F77 code does not have to be modified to
2 : ! reference the CARMA structure.
3 : #include "carma_globaer.h"
4 :
5 : !! This routine calculates particle source terms due to evaporation <evappe>.
6 : !!
7 : !! @author Andy Ackerman
8 : !! @version Aug-2001
9 213271518 : subroutine evapp(carma, cstate, iz, rc)
10 :
11 : ! types
12 : use carma_precision_mod
13 : use carma_enums_mod
14 : use carma_constants_mod
15 : use carma_types_mod
16 : use carmastate_mod
17 : use carma_mod
18 :
19 : implicit none
20 :
21 : type(carma_type), intent(in) :: carma !! the carma object
22 : type(carmastate_type), intent(inout) :: cstate !! the carma state object
23 : integer, intent(in) :: iz !! z index
24 : integer, intent(inout) :: rc !! return code, negative indicates failure
25 :
26 : ! Local declarations
27 : integer :: ibin !! bin index
28 : integer :: ielem !! element index
29 : integer :: ig !! source group index
30 : integer :: ip !! source number concentration element
31 : integer :: ic
32 : integer :: ic1 !! element of first core mass in group
33 : integer :: iecore
34 : integer :: ieto
35 : integer :: igto
36 : integer :: iavg
37 : logical :: evap_total
38 : real(kind=f) :: sig_mono
39 : real(kind=f) :: coretot
40 : real(kind=f) :: coremom
41 : real(kind=f) :: smf
42 : integer :: nbin
43 :
44 :
45 : ! Define criterion for monodisperse core mass distributions
46 213271518 : sig_mono = sqrt( ALMOST_ZERO )
47 :
48 : ! Loop over source groups (from which evaporation is being treated)
49 639814554 : do ig = 1, NGROUP
50 :
51 426543036 : ip = ienconc(ig)
52 :
53 : ! No evaporation unless particles are volatile
54 639814554 : if( itype(ip) .eq. I_VOLATILE )then
55 :
56 : ! Make sure that these always get intializaed, since they can
57 : ! cause problems in other parts of the code if they aren't.
58 8957403756 : totevap(:,ig) = .false.
59 8957403756 : cmf(:,ig) = 0._f
60 :
61 426543036 : if (pconmax(iz, ig) > FEW_PC) then
62 :
63 420623935 : ic1 = icorelem(1,ig)
64 :
65 : ! Loop over source bins and calculate temporary evaporation source
66 : ! for droplets in next smaller bin assuming no total evaporation <evdrop>
67 8833102635 : do ibin = 1, NBIN
68 8412478700 : evdrop = pc(iz,ibin,ip)*evaplg(ibin,ig)
69 :
70 : ! Check for evaporation of a sufficient number of droplets
71 : ! if( evdrop .gt. 0._f .and. pc(iz,ibin,ip) .gt. SMALL_PC )then
72 8833102635 : if( evdrop .gt. 0._f )then
73 :
74 : ! No cores: transfer droplets within group
75 244226427 : if( ic1 .eq. 0 )then
76 244226427 : call evap_ingrp(carma,cstate,iz,ibin,ig,ip,rc)
77 : else
78 :
79 : ! First core is not involatile (therefore none are)
80 : ! -- this is a hack until enforced/checked in setupbins() --
81 : ! transfer droplets within group
82 : !
83 0 : if( itype(ic1) .ne. I_COREMASS )then
84 0 : call evap_ingrp(carma,cstate,iz,ibin,ig,ip,rc)
85 : else
86 :
87 : ! Have cores: calculate <evcore> the amount of the source term
88 : ! by number <evdrop> associated with total evaporation of secondary cores
89 0 : coretot = pc(iz,ibin,ic1)
90 0 : do ic = 2, ncore(ig)
91 0 : iecore = icorelem(ic,ig)
92 0 : if( itype(iecore) .eq. I_COREMASS )then
93 0 : coretot = coretot + pc(iz,ibin,iecore)
94 : endif
95 : enddo
96 0 : do ic = 2, ncore(ig)
97 0 : iecore = icorelem(ic,ig)
98 0 : if( itype(iecore) .eq. I_COREMASS )then
99 0 : evcore(ic) = evdrop*pc(iz,ibin,iecore)/coretot
100 : endif
101 : enddo
102 :
103 : ! Calculate average particle core mass and fraction
104 0 : coreavg = coretot / pc(iz,ibin,ip)
105 0 : coreavg = min( rmass(ibin,ig), coreavg )
106 0 : cmf(ibin,ig) = coreavg / rmass(ibin,ig)
107 : ! cmf(ibin,ig) = max( 0., min( ONE, cmf(ibin,ig) ) )
108 :
109 : ! Get target number concentration element and group for total evaporation
110 : ! and evaluate logical flags regarding position on CN bin and index of
111 : ! target CN bin
112 0 : ieto = ievp2elem(ic1)
113 :
114 : ! To treat internal mixtures, it is possible for the condensate to
115 : ! totally evaporate and have core mass, but for there not to be another
116 : ! group to which the core mass should go. So allow no evp2elem, but
117 : ! always use the in group evaporation.
118 0 : if (ieto == 0) then
119 0 : nuc_small = .false.
120 : else
121 0 : igto = igelem(ieto)
122 :
123 0 : too_small = coreavg .lt. rmass(1,igto)
124 0 : nbin = NBIN
125 0 : too_big = coreavg .gt. rmass(nbin,igto)
126 :
127 0 : if( .not. (too_small .or. too_big) )then
128 : iavg = log( coreavg / rmassmin(igto) ) / &
129 0 : log( rmrat(igto) ) + 2
130 0 : iavg = min( iavg, NBIN )
131 : endif
132 :
133 : ! Only consider size of evaporating cores relative to nuc_small
134 : ! when treating core second moment for this particle group
135 0 : if( if_sec_mom(ig) )then
136 0 : nuc_small = coreavg .lt. rmass(1,igto)
137 : else
138 0 : nuc_small = .false.
139 : endif
140 : end if
141 :
142 : ! Want total evaporation when
143 : ! cores smaller than smallest nucleated
144 : ! OR evaporating droplets are in bin 1
145 : ! OR droplets will be created with core mass fraction > 1
146 : evap_total = nuc_small .or. ibin .eq. 1 .or. &
147 0 : rmrat(ig)*cmf(ibin,ig) .gt. ONE
148 :
149 : ! No core second moment: evaporate to monodisperse CN cores or within group.!
150 0 : if( .not. if_sec_mom(ig) )then
151 :
152 0 : if( evap_total .and. (ieto /= 0) )then
153 0 : call evap_mono(carma,cstate,iz,ibin,ig,iavg,ieto,igto,rc)
154 : else
155 0 : call evap_ingrp(carma,cstate,iz,ibin,ig,ip,rc)
156 : endif
157 :
158 : ! Have core second moments: evaporate to mono- or polydisperse CN cores
159 : ! or within group. First calculate average core second moment <coremom>,
160 : ! second moment fraction <smf>, and square of the logarithm of the geometric
161 : ! standard deviation of the assumed core mass distribution <coresig>.
162 : else
163 :
164 0 : coremom = pc(iz,ibin,imomelem(ig)) / pc(iz,ibin,ip)
165 0 : smf = coremom / rmass(ibin,ig)**2
166 0 : coresig = log( smf / cmf(ibin,ig)**2 )
167 :
168 : ! Want total evaporation for above reasons
169 : ! OR droplets will be created with core moment fraction > 1
170 0 : evap_total = evap_total .or. rmrat(ig)**2*smf .gt. ONE
171 :
172 0 : if( evap_total .and. (ieto /= 0) )then
173 :
174 : ! Want monodisperse total evaporation when
175 : ! cores smaller than smallest nucleated
176 : ! OR evaporating core distribution is narrow
177 : ! Otherwise want polydisperse total evaporation
178 0 : if( nuc_small .or. coresig .le. sig_mono )then
179 0 : call evap_mono(carma,cstate,iz,ibin,ig,iavg,ieto,igto,rc)
180 : else
181 0 : call evap_poly(carma,cstate,iz,ibin,ig,iavg,ieto,igto,rc)
182 : endif
183 :
184 : ! Droplet evaporation within group
185 : else
186 0 : call evap_ingrp(carma,cstate,iz,ibin,ig,ip,rc)
187 : endif
188 : endif ! if_sec_mom(ig)
189 : endif ! itype(ic1)
190 : endif ! ic1=0
191 : endif ! evaplg > 0
192 : enddo ! ibin=1,NBIN
193 : endif ! enough particles
194 : endif ! volatile particles
195 : enddo ! ig=1,NGROUP
196 :
197 : ! Return to caller with evaporation production terms evaluated.
198 213271518 : return
199 213271518 : end
|