Line data Source code
1 : module aerosol_state_mod
2 : use shr_kind_mod, only: r8 => shr_kind_r8
3 : use aerosol_properties_mod, only: aerosol_properties, aero_name_len
4 : use physconst, only: pi
5 :
6 : implicit none
7 :
8 : private
9 :
10 : public :: aerosol_state
11 : public :: ptr2d_t
12 :
13 : !> aerosol_state defines the interface to the time-varying aerosol state
14 : !! variables (e.g., mixing ratios, number concentrations). This includes the
15 : !! aerosol portion of the overall model state.
16 : !!
17 : !! Each aerosol package (e.g., MAM, CARMA, etc) must extend the aerosol_state
18 : !! class to allow access to the state information (transported and not transported)
19 : !! of the aerosol package. Any package must implement each of the deferred
20 : !! procedures of the abstract aerosol_state class, may include additional private
21 : !! data members and type-bound procedures, and may override functions of the
22 : !! abstract class.
23 : !!
24 : !! Please see the modal_aerosol_state module for an example of how the aerosol_state
25 : !! class can be extended for a specific aerosol package.
26 : type, abstract :: aerosol_state
27 : contains
28 : procedure(aero_get_transported), deferred :: get_transported
29 : procedure(aero_set_transported), deferred :: set_transported
30 : procedure(aero_get_amb_total_bin_mmr), deferred :: ambient_total_bin_mmr
31 : procedure(aero_get_state_mmr), deferred :: get_ambient_mmr_0list
32 : procedure(aero_get_list_mmr), deferred :: get_ambient_mmr_rlist
33 : generic :: get_ambient_mmr=>get_ambient_mmr_0list,get_ambient_mmr_rlist
34 : procedure(aero_get_state_mmr), deferred :: get_cldbrne_mmr
35 : procedure(aero_get_state_num), deferred :: get_ambient_num
36 : procedure(aero_get_state_num), deferred :: get_cldbrne_num
37 : procedure(aero_get_states), deferred :: get_states
38 : procedure(aero_update_bin), deferred :: update_bin
39 : procedure :: loadaer
40 : procedure(aero_icenuc_size_wght_arr), deferred :: icenuc_size_wght_arr
41 : procedure(aero_icenuc_size_wght_val), deferred :: icenuc_size_wght_val
42 : generic :: icenuc_size_wght => icenuc_size_wght_arr,icenuc_size_wght_val
43 : procedure :: icenuc_type_wght_base
44 : procedure :: icenuc_type_wght => icenuc_type_wght_base
45 : procedure :: nuclice_get_numdens
46 : procedure :: get_amb_species_numdens
47 : procedure :: get_cld_species_numdens
48 : procedure :: coated_frac
49 : procedure :: mass_mean_radius
50 : procedure :: watact_mfactor
51 : procedure(aero_hetfrz_size_wght), deferred :: hetfrz_size_wght
52 : procedure(aero_hygroscopicity), deferred :: hygroscopicity
53 : procedure(aero_water_uptake), deferred :: water_uptake
54 : procedure(aero_wgtpct), deferred :: wgtpct
55 : procedure :: refractive_index_sw
56 : procedure :: refractive_index_lw
57 : procedure(aero_volume), deferred :: dry_volume
58 : procedure(aero_volume), deferred :: wet_volume
59 : procedure(aero_volume), deferred :: water_volume
60 : procedure(aero_wet_diam), deferred :: wet_diameter
61 : procedure :: convcld_actfrac
62 : procedure :: sol_factb_interstitial
63 : end type aerosol_state
64 :
65 : ! for state fields
66 : type ptr2d_t
67 : real(r8), pointer :: fld(:,:)
68 : end type ptr2d_t
69 :
70 : real(r8), parameter :: per_cm3 = 1.e-6_r8 ! factor for m-3 to cm-3 conversions
71 : real(r8), parameter :: per_m3 = 1.e6_r8 ! factor for cm-3 to m-3 conversions
72 : real(r8), parameter :: kg2mug = 1.e9_r8 ! factor for kg to micrograms (mug) conversions
73 :
74 : abstract interface
75 :
76 : !------------------------------------------------------------------------
77 : ! Total aerosol mass mixing ratio for a bin in a given grid box location (column and layer)
78 : !------------------------------------------------------------------------
79 : function aero_get_amb_total_bin_mmr(self, aero_props, bin_ndx, col_ndx, lyr_ndx) result(mmr_tot)
80 : import :: aerosol_state, aerosol_properties, r8
81 : class(aerosol_state), intent(in) :: self
82 : class(aerosol_properties), intent(in) :: aero_props ! aerosol properties object
83 : integer, intent(in) :: bin_ndx ! bin index
84 : integer, intent(in) :: col_ndx ! column index
85 : integer, intent(in) :: lyr_ndx ! vertical layer index
86 :
87 : real(r8) :: mmr_tot ! mass mixing ratios totaled for all species
88 :
89 : end function aero_get_amb_total_bin_mmr
90 :
91 : !------------------------------------------------------------------------
92 : ! returns aerosol mass mixing ratio for a given species index and bin index
93 : !------------------------------------------------------------------------
94 : subroutine aero_get_state_mmr(self, species_ndx, bin_ndx, mmr)
95 : import :: aerosol_state, r8
96 : class(aerosol_state), intent(in) :: self
97 : integer, intent(in) :: species_ndx ! species index
98 : integer, intent(in) :: bin_ndx ! bin index
99 : real(r8), pointer :: mmr(:,:) ! mass mixing ratios (ncol,nlev)
100 : end subroutine aero_get_state_mmr
101 :
102 : !------------------------------------------------------------------------
103 : ! returns aerosol mass mixing ratio for a given species index, bin index
104 : ! and raditaion climate or diagnsotic list number
105 : !------------------------------------------------------------------------
106 : subroutine aero_get_list_mmr(self, list_ndx, species_ndx, bin_ndx, mmr)
107 : import :: aerosol_state, r8
108 : class(aerosol_state), intent(in) :: self
109 : integer, intent(in) :: list_ndx ! rad climate/diagnostic list index
110 : integer, intent(in) :: species_ndx ! species index
111 : integer, intent(in) :: bin_ndx ! bin index
112 : real(r8), pointer :: mmr(:,:) ! mass mixing ratios (ncol,nlev)
113 : end subroutine aero_get_list_mmr
114 :
115 : !------------------------------------------------------------------------
116 : ! returns aerosol number mixing ratio for a given species index and bin index
117 : !------------------------------------------------------------------------
118 : subroutine aero_get_state_num(self, bin_ndx, num)
119 : import :: aerosol_state, r8
120 : class(aerosol_state), intent(in) :: self
121 : integer, intent(in) :: bin_ndx ! bin index
122 : real(r8), pointer :: num(:,:) ! number densities (ncol,nlev)
123 : end subroutine aero_get_state_num
124 :
125 : !------------------------------------------------------------------------
126 : ! returns interstitial and cloud-borne aerosol states
127 : !------------------------------------------------------------------------
128 : subroutine aero_get_states( self, aero_props, raer, qqcw )
129 : import :: aerosol_state, aerosol_properties, ptr2d_t
130 :
131 : class(aerosol_state), intent(in) :: self
132 : class(aerosol_properties), intent(in) :: aero_props ! properties of the aerosol model
133 : type(ptr2d_t), intent(out) :: raer(:) ! state of interstitial aerosols
134 : type(ptr2d_t), intent(out) :: qqcw(:) ! state of cloud-borne aerosols
135 :
136 : end subroutine aero_get_states
137 :
138 : !------------------------------------------------------------------------------
139 : ! sets transported components
140 : ! This updates the aerosol model state from the host transported aerosol constituents array.
141 : ! (mass mixing ratios or number mixing ratios)
142 : !------------------------------------------------------------------------------
143 : subroutine aero_set_transported( self, transported_array )
144 : import :: aerosol_state, r8
145 : class(aerosol_state), intent(inout) :: self
146 : real(r8), intent(in) :: transported_array(:,:,:)
147 : end subroutine aero_set_transported
148 :
149 : !------------------------------------------------------------------------------
150 : ! returns transported components
151 : ! This updates the transported aerosol constituent array to match the aerosol model state.
152 : ! (mass mixing ratios or number mixing ratios)
153 : !------------------------------------------------------------------------------
154 : subroutine aero_get_transported( self, transported_array )
155 : import :: aerosol_state, r8
156 : class(aerosol_state), intent(in) :: self
157 : real(r8), intent(out) :: transported_array(:,:,:)
158 : end subroutine aero_get_transported
159 :
160 : !------------------------------------------------------------------------------
161 : ! return aerosol bin size weights for a given bin
162 : !------------------------------------------------------------------------------
163 : subroutine aero_icenuc_size_wght_arr(self, bin_ndx, ncol, nlev, species_type, use_preexisting_ice, wght)
164 : import :: aerosol_state, r8
165 : class(aerosol_state), intent(in) :: self
166 : integer, intent(in) :: bin_ndx ! bin number
167 : integer, intent(in) :: ncol ! number of columns
168 : integer, intent(in) :: nlev ! number of vertical levels
169 : character(len=*), intent(in) :: species_type ! species type
170 : logical, intent(in) :: use_preexisting_ice ! pre-existing ice flag
171 : real(r8), intent(out) :: wght(:,:)
172 :
173 : end subroutine aero_icenuc_size_wght_arr
174 :
175 : !------------------------------------------------------------------------------
176 : ! return aerosol bin size weights for a given bin, column and vertical layer
177 : !------------------------------------------------------------------------------
178 : subroutine aero_icenuc_size_wght_val(self, bin_ndx, col_ndx, lyr_ndx, species_type, use_preexisting_ice, wght)
179 : import :: aerosol_state, r8
180 : class(aerosol_state), intent(in) :: self
181 : integer, intent(in) :: bin_ndx ! bin number
182 : integer, intent(in) :: col_ndx ! column index
183 : integer, intent(in) :: lyr_ndx ! vertical layer index
184 : character(len=*), intent(in) :: species_type ! species type
185 : logical, intent(in) :: use_preexisting_ice ! pre-existing ice flag
186 : real(r8), intent(out) :: wght
187 :
188 : end subroutine aero_icenuc_size_wght_val
189 :
190 : !------------------------------------------------------------------------------
191 : ! updates state and tendency
192 : !------------------------------------------------------------------------------
193 : subroutine aero_update_bin( self, bin_ndx, col_ndx, lyr_ndx, delmmr_sum, delnum_sum, tnd_ndx, dtime, tend )
194 : import :: aerosol_state, r8
195 : class(aerosol_state), intent(in) :: self
196 : integer, intent(in) :: bin_ndx ! bin number
197 : integer, intent(in) :: col_ndx ! column index
198 : integer, intent(in) :: lyr_ndx ! vertical layer index
199 : real(r8),intent(in) :: delmmr_sum ! mass mixing ratio change summed over all species in bin
200 : real(r8),intent(in) :: delnum_sum ! number mixing ratio change summed over all species in bin
201 : integer, intent(in) :: tnd_ndx ! tendency index
202 : real(r8),intent(in) :: dtime ! time step size (sec)
203 : real(r8),intent(inout) :: tend(:,:,:) ! tendency
204 :
205 : end subroutine aero_update_bin
206 :
207 : !------------------------------------------------------------------------------
208 : ! returns the volume-weighted fractions of aerosol subset `bin_ndx` that can act
209 : ! as heterogeneous freezing nuclei
210 : !------------------------------------------------------------------------------
211 : function aero_hetfrz_size_wght(self, bin_ndx, ncol, nlev) result(wght)
212 : import :: aerosol_state, r8
213 : class(aerosol_state), intent(in) :: self
214 : integer, intent(in) :: bin_ndx ! bin number
215 : integer, intent(in) :: ncol ! number of columns
216 : integer, intent(in) :: nlev ! number of vertical levels
217 :
218 : real(r8) :: wght(ncol,nlev)
219 :
220 : end function aero_hetfrz_size_wght
221 :
222 : !------------------------------------------------------------------------------
223 : ! returns hygroscopicity for a given radiation diagnostic list number and
224 : ! bin number
225 : !------------------------------------------------------------------------------
226 : subroutine aero_hygroscopicity(self, list_ndx, bin_ndx, kappa)
227 : import :: aerosol_state, r8
228 : class(aerosol_state), intent(in) :: self
229 : integer, intent(in) :: list_ndx ! rad climate/diagnostic list index
230 : integer, intent(in) :: bin_ndx ! bin number
231 :
232 : real(r8), intent(out) :: kappa(:,:) ! hygroscopicity (ncol,nlev)
233 :
234 : end subroutine aero_hygroscopicity
235 :
236 : !------------------------------------------------------------------------------
237 : ! returns aerosol wet diameter and aerosol water concentration for a given
238 : ! radiation diagnostic list number and bin number
239 : !------------------------------------------------------------------------------
240 : subroutine aero_water_uptake(self, aero_props, list_idx, bin_idx, ncol, nlev, dgnumwet, qaerwat)
241 : import :: aerosol_state, aerosol_properties, r8
242 :
243 : class(aerosol_state), intent(in) :: self
244 : class(aerosol_properties), intent(in) :: aero_props
245 : integer, intent(in) :: list_idx ! rad climate/diags list number
246 : integer, intent(in) :: bin_idx ! bin number
247 : integer, intent(in) :: ncol ! number of columns
248 : integer, intent(in) :: nlev ! number of levels
249 : real(r8),intent(out) :: dgnumwet(ncol,nlev) ! aerosol wet diameter (m)
250 : real(r8),intent(out) :: qaerwat(ncol,nlev) ! aerosol water concentration (g/g)
251 :
252 : end subroutine aero_water_uptake
253 :
254 : !------------------------------------------------------------------------------
255 : ! aerosol weight precent of H2SO4/H2O solution
256 : !------------------------------------------------------------------------------
257 : function aero_wgtpct(self, ncol, nlev) result(wtp)
258 : import :: aerosol_state, r8
259 : class(aerosol_state), intent(in) :: self
260 : integer, intent(in) :: ncol,nlev
261 : real(r8) :: wtp(ncol,nlev) ! weight precent of H2SO4/H2O solution for given icol, ilev
262 :
263 : end function aero_wgtpct
264 :
265 : !------------------------------------------------------------------------------
266 : ! aerosol volume interface
267 : !------------------------------------------------------------------------------
268 : function aero_volume(self, aero_props, list_idx, bin_idx, ncol, nlev) result(vol)
269 : import :: aerosol_state, aerosol_properties, r8
270 :
271 : class(aerosol_state), intent(in) :: self
272 : class(aerosol_properties), intent(in) :: aero_props
273 : integer, intent(in) :: list_idx ! rad climate/diags list number
274 : integer, intent(in) :: bin_idx ! bin number
275 : integer, intent(in) :: ncol ! number of columns
276 : integer, intent(in) :: nlev ! number of levels
277 :
278 : real(r8) :: vol(ncol,nlev) ! m3/kg
279 :
280 : end function aero_volume
281 :
282 : !------------------------------------------------------------------------------
283 : ! aerosol wet diameter
284 : !------------------------------------------------------------------------------
285 : function aero_wet_diam(self, bin_idx, ncol, nlev) result(diam)
286 : import :: aerosol_state, r8
287 :
288 : class(aerosol_state), intent(in) :: self
289 : integer, intent(in) :: bin_idx ! bin number
290 : integer, intent(in) :: ncol ! number of columns
291 : integer, intent(in) :: nlev ! number of levels
292 :
293 : real(r8) :: diam(ncol,nlev)
294 :
295 : end function aero_wet_diam
296 :
297 : end interface
298 :
299 : contains
300 :
301 : !------------------------------------------------------------------------------
302 : ! returns aerosol number, volume concentrations, and bulk hygroscopicity
303 : !------------------------------------------------------------------------------
304 2555442000 : subroutine loadaer( self, aero_props, istart, istop, k, m, cs, phase, &
305 851814000 : naerosol, vaerosol, hygro, errnum, errstr, pom_hygro)
306 :
307 : use aerosol_properties_mod, only: aerosol_properties
308 :
309 : ! input arguments
310 : class(aerosol_state), intent(in) :: self
311 : class(aerosol_properties), intent(in) :: aero_props
312 :
313 : integer, intent(in) :: istart ! start column index (1 <= istart <= istop <= pcols)
314 : integer, intent(in) :: istop ! stop column index
315 : integer, intent(in) :: k ! level index
316 : integer, intent(in) :: m ! mode or bin index
317 : real(r8), intent(in) :: cs(:,:) ! air density (kg/m3)
318 : integer, intent(in) :: phase ! phase of aerosol: 1 for interstitial, 2 for cloud-borne, 3 for sum
319 :
320 : ! output arguments
321 : real(r8), intent(out) :: naerosol(:) ! number conc (1/m3)
322 : real(r8), intent(out) :: vaerosol(:) ! volume conc (m3/m3)
323 : real(r8), intent(out) :: hygro(:) ! bulk hygroscopicity of mode
324 :
325 : integer , intent(out) :: errnum
326 : character(len=*), intent(out) :: errstr
327 :
328 : real(r8), optional, intent(in) :: pom_hygro ! POM hygroscopicity override
329 :
330 : ! internal
331 851814000 : real(r8), pointer :: raer(:,:) ! interstitial aerosol mass, number mixing ratios
332 851814000 : real(r8), pointer :: qqcw(:,:) ! cloud-borne aerosol mass, number mixing ratios
333 : real(r8) :: specdens, spechygro
334 : character(len=aero_name_len) :: spectype
335 :
336 851814000 : real(r8) :: vol(istart:istop) ! aerosol volume mixing ratio
337 : integer :: i, l
338 : !-------------------------------------------------------------------------------
339 851814000 : errnum = 0
340 :
341 7279400160 : do i = istart, istop
342 6427586160 : vaerosol(i) = 0._r8
343 7279400160 : hygro(i) = 0._r8
344 : end do
345 :
346 3833163000 : do l = 1, aero_props%nspecies(m)
347 :
348 2981349000 : call self%get_ambient_mmr(l,m, raer)
349 2981349000 : call self%get_cldbrne_mmr(l,m, qqcw)
350 2981349000 : call aero_props%get(m,l, density=specdens, hygro=spechygro, spectype=spectype)
351 2981349000 : if (present(pom_hygro)) then
352 0 : if (spectype=='p-organic'.and.pom_hygro>0._r8) then
353 0 : spechygro=pom_hygro
354 : endif
355 : endif
356 :
357 2981349000 : if (phase == 3) then
358 22427919360 : do i = istart, istop
359 22427919360 : vol(i) = max(raer(i,k) + qqcw(i,k), 0._r8)/specdens
360 : end do
361 1524990600 : else if (phase == 2) then
362 0 : do i = istart, istop
363 0 : vol(i) = max(qqcw(i,k), 0._r8)/specdens
364 : end do
365 1524990600 : else if (phase == 1) then
366 3049981200 : do i = istart, istop
367 3049981200 : vol(i) = max(raer(i,k), 0._r8)/specdens
368 : end do
369 : else
370 0 : errnum = -1
371 0 : write(errstr,*)'phase = ',phase,' in aerosol_state::loadaer not recognized'
372 : return
373 : end if
374 :
375 29311063560 : do i = istart, istop
376 22496551560 : vaerosol(i) = vaerosol(i) + vol(i)
377 25477900560 : hygro(i) = hygro(i) + vol(i)*spechygro
378 : end do
379 :
380 : end do
381 :
382 7279400160 : do i = istart, istop
383 7279400160 : if (vaerosol(i) > 1.0e-30_r8) then
384 2942080969 : hygro(i) = hygro(i)/(vaerosol(i))
385 2942080969 : vaerosol(i) = vaerosol(i)*cs(i,k)
386 : else
387 3485505191 : hygro(i) = 0.0_r8
388 3485505191 : vaerosol(i) = 0.0_r8
389 : end if
390 : end do
391 :
392 : ! aerosol number mixing ratios (#/kg)
393 851814000 : call self%get_ambient_num(m, raer)
394 851814000 : call self%get_cldbrne_num(m, qqcw)
395 851814000 : if (phase == 3) then
396 6407976960 : do i = istart, istop
397 6407976960 : naerosol(i) = (raer(i,k) + qqcw(i,k))*cs(i,k) ! #/kg -> #/m3
398 : end do
399 435711600 : else if (phase == 2) then
400 0 : do i = istart, istop
401 0 : naerosol(i) = qqcw(i,k)*cs(i,k)
402 : end do
403 : else
404 871423200 : do i = istart, istop
405 871423200 : naerosol(i) = raer(i,k)*cs(i,k)
406 : end do
407 : end if
408 :
409 : ! adjust number
410 851814000 : call aero_props%apply_number_limits( naerosol, vaerosol, istart, istop, m )
411 :
412 851814000 : end subroutine loadaer
413 :
414 : !------------------------------------------------------------------------------
415 : ! returns ambient aerosol number density for a given bin number and species type
416 : !------------------------------------------------------------------------------
417 29030400 : subroutine get_amb_species_numdens(self, bin_ndx, ncol, nlev, species_type, aero_props, rho, numdens)
418 : use aerosol_properties_mod, only: aerosol_properties
419 : class(aerosol_state), intent(in) :: self
420 : integer, intent(in) :: bin_ndx ! bin number
421 : integer, intent(in) :: ncol ! number of columns
422 : integer, intent(in) :: nlev ! number of vertical levels
423 : character(len=*), intent(in) :: species_type ! species type
424 : class(aerosol_properties), intent(in) :: aero_props ! aerosol properties object
425 : real(r8), intent(in) :: rho(:,:) ! air density (kg m-3)
426 : real(r8), intent(out) :: numdens(:,:) ! species number densities (#/cm^3)
427 :
428 29030400 : real(r8), pointer :: num(:,:)
429 58060800 : real(r8) :: type_wght(ncol,nlev)
430 58060800 : real(r8) :: size_wght(ncol,nlev)
431 :
432 31323801600 : size_wght = self%hetfrz_size_wght(bin_ndx, ncol, nlev)
433 :
434 29030400 : call self%icenuc_type_wght_base(bin_ndx, ncol, nlev, species_type, aero_props, rho, type_wght)
435 :
436 29030400 : call self%get_ambient_num(bin_ndx, num)
437 :
438 31323801600 : numdens(:ncol,:) = num(:ncol,:)*rho(:ncol,:)*type_wght(:ncol,:)*size_wght(:ncol,:)*per_cm3
439 :
440 29030400 : end subroutine get_amb_species_numdens
441 :
442 : !------------------------------------------------------------------------------
443 : ! returns cloud-borne aerosol number density for a given bin number and species type
444 : !------------------------------------------------------------------------------
445 9676800 : subroutine get_cld_species_numdens(self, bin_ndx, ncol, nlev, species_type, aero_props, rho, numdens)
446 : use aerosol_properties_mod, only: aerosol_properties
447 : class(aerosol_state), intent(in) :: self
448 : integer, intent(in) :: bin_ndx ! bin number
449 : integer, intent(in) :: ncol ! number of columns
450 : integer, intent(in) :: nlev ! number of vertical levels
451 : character(len=*), intent(in) :: species_type ! species type
452 : class(aerosol_properties), intent(in) :: aero_props ! aerosol properties object
453 : real(r8), intent(in) :: rho(:,:) ! air density (kg m-3)
454 : real(r8), intent(out) :: numdens(:,:) ! number densities (#/cm^3)
455 :
456 9676800 : real(r8), pointer :: num(:,:)
457 19353600 : real(r8) :: type_wght(ncol,nlev)
458 19353600 : real(r8) :: size_wght(ncol,nlev)
459 :
460 10441267200 : size_wght = self%hetfrz_size_wght(bin_ndx, ncol, nlev)
461 :
462 9676800 : call self%icenuc_type_wght_base(bin_ndx, ncol, nlev, species_type, aero_props, rho, type_wght, cloud_borne=.true.)
463 :
464 9676800 : call self%get_cldbrne_num(bin_ndx, num)
465 :
466 10441267200 : numdens(:ncol,:) = num(:ncol,:)*rho(:ncol,:)*type_wght(:ncol,:)*size_wght(:ncol,:)*per_cm3
467 :
468 9676800 : end subroutine get_cld_species_numdens
469 :
470 : !------------------------------------------------------------------------------
471 : ! returns aerosol type weights for a given aerosol type and bin
472 : !------------------------------------------------------------------------------
473 82252800 : subroutine icenuc_type_wght_base(self, bin_ndx, ncol, nlev, species_type, aero_props, rho, wght, cloud_borne)
474 :
475 : use aerosol_properties_mod, only: aerosol_properties
476 :
477 : class(aerosol_state), intent(in) :: self
478 : integer, intent(in) :: bin_ndx ! bin number
479 : integer, intent(in) :: ncol ! number of columns
480 : integer, intent(in) :: nlev ! number of vertical levels
481 : character(len=*), intent(in) :: species_type ! species type
482 : class(aerosol_properties), intent(in) :: aero_props ! aerosol properties object
483 : real(r8), intent(in) :: rho(:,:) ! air density (kg m-3)
484 : real(r8), intent(out) :: wght(:,:) ! type weights
485 : logical, optional, intent(in) :: cloud_borne ! if TRUE cloud-borne aerosols are used
486 : ! otherwise ambient aerosols are used
487 :
488 164505600 : real(r8) :: mass(ncol,nlev)
489 164505600 : real(r8) :: totalmass(ncol,nlev)
490 82252800 : real(r8), pointer :: aer_bin(:,:)
491 :
492 : character(len=aero_name_len) :: spectype, sptype
493 : integer :: ispc
494 : logical :: cldbrne
495 :
496 82252800 : if (present(cloud_borne)) then
497 9676800 : cldbrne = cloud_borne
498 : else
499 : cldbrne = .false.
500 : end if
501 :
502 88750771200 : wght(:,:) = 0._r8
503 88750771200 : totalmass(:,:) = 0._r8
504 88750771200 : mass(:,:) = 0._r8
505 :
506 82252800 : if (species_type=='sulfate_strat') then
507 9676800 : sptype = 'sulfate'
508 : else
509 72576000 : sptype = species_type
510 : end if
511 :
512 527385600 : do ispc = 1, aero_props%nspecies(bin_ndx)
513 :
514 445132800 : if (cldbrne) then
515 58060800 : call self%get_cldbrne_mmr(ispc, bin_ndx, aer_bin)
516 : else
517 387072000 : call self%get_ambient_mmr(ispc, bin_ndx, aer_bin)
518 : end if
519 445132800 : call aero_props%species_type(bin_ndx, ispc, spectype=spectype)
520 :
521 >48029*10^7 : totalmass(:ncol,:) = totalmass(:ncol,:) + aer_bin(:ncol,:)*rho(:ncol,:)
522 :
523 527385600 : if (trim(spectype) == trim(sptype)) then
524 88750771200 : mass(:ncol,:) = mass(:ncol,:) + aer_bin(:ncol,:)*rho(:ncol,:)
525 : end if
526 :
527 : end do
528 :
529 88750771200 : where (totalmass(:ncol,:) > 0._r8)
530 82252800 : wght(:ncol,:) = mass(:ncol,:)/totalmass(:ncol,:)
531 : end where
532 :
533 164505600 : end subroutine icenuc_type_wght_base
534 :
535 : !------------------------------------------------------------------------------
536 241920 : subroutine nuclice_get_numdens(self, aero_props, use_preexisting_ice, ncol, nlev, rho, dust_num_col, sulf_num_col, soot_num_col, sulf_num_tot_col )
537 :
538 : class(aerosol_state), intent(in) :: self
539 : class(aerosol_properties), intent(in) :: aero_props ! aerosol properties object
540 :
541 : logical, intent(in) :: use_preexisting_ice
542 : integer, intent(in) :: ncol ! number of columns
543 : integer, intent(in) :: nlev ! number of vertical levels
544 : real(r8), intent(in) :: rho(:,:) ! air density (kg m-3)
545 : real(r8), intent(out) :: dust_num_col(:,:) ! dust number densities (#/cm^3)
546 : real(r8), intent(out) :: sulf_num_col(:,:) ! sulfate number densities (#/cm^3)
547 : real(r8), intent(out) :: soot_num_col(:,:) ! soot number densities (#/cm^3)
548 : real(r8), intent(out) :: sulf_num_tot_col(:,:) ! stratopsheric sulfate number densities (#/cm^3)
549 :
550 : integer :: ibin,ispc
551 : character(len=aero_name_len) :: spectype
552 483840 : real(r8) :: size_wghts(ncol,nlev)
553 483840 : real(r8) :: type_wghts(ncol,nlev)
554 :
555 241920 : real(r8), pointer :: num_col(:,:)
556 :
557 288126720 : dust_num_col(:,:) = 0._r8
558 288126720 : sulf_num_col(:,:) = 0._r8
559 288126720 : soot_num_col(:,:) = 0._r8
560 288126720 : sulf_num_tot_col(:,:) = 0._r8
561 :
562 : ! collect number densities (#/cm^3) for dust, sulfate, and soot
563 9918720 : do ibin = 1,aero_props%nbins()
564 :
565 9676800 : call self%get_ambient_num(ibin, num_col)
566 :
567 43545600 : do ispc = 1,aero_props%nspecies(ibin)
568 :
569 33868800 : call aero_props%species_type(ibin, ispc, spectype)
570 :
571 33868800 : call self%icenuc_size_wght(ibin, ncol, nlev, spectype, use_preexisting_ice, size_wghts)
572 :
573 33868800 : call self%icenuc_type_wght(ibin, ncol, nlev, spectype, aero_props, rho, type_wghts)
574 :
575 77414400 : select case ( trim(spectype) )
576 : case('dust')
577 4838400 : dust_num_col(:ncol,:) = dust_num_col(:ncol,:) &
578 5225472000 : + size_wghts(:ncol,:)*type_wghts(:ncol,:)*num_col(:ncol,:)*rho(:ncol,:)*per_cm3
579 : case('sulfate')
580 : ! This order of ops gives bit-for-bit results for cam5 phys ( use_preexisting_ice = .false. )
581 9676800 : sulf_num_col(:ncol,:) = sulf_num_col(:ncol,:) &
582 10450944000 : + num_col(:ncol,:)*rho(:ncol,:)*per_cm3 * size_wghts(:ncol,:)*type_wghts(:ncol,:)
583 : case('black-c')
584 4838400 : soot_num_col(:ncol,:) = soot_num_col(:ncol,:) &
585 5288371200 : + size_wghts(:ncol,:)*type_wghts(:ncol,:)*num_col(:ncol,:)*rho(:ncol,:)*per_cm3
586 : end select
587 :
588 : enddo
589 :
590 : ! stratospheric sulfates -- special case not included in the species loop above
591 9676800 : call self%icenuc_size_wght(ibin, ncol, nlev, 'sulfate_strat', use_preexisting_ice, size_wghts)
592 9676800 : call self%icenuc_type_wght(ibin, ncol, nlev, 'sulfate_strat', aero_props, rho, type_wghts)
593 9676800 : sulf_num_tot_col(:ncol,:) = sulf_num_tot_col(:ncol,:) &
594 10451185920 : + size_wghts(:ncol,:)*type_wghts(:ncol,:)*num_col(:ncol,:)*rho(:ncol,:)*per_cm3
595 :
596 : enddo
597 :
598 241920 : end subroutine nuclice_get_numdens
599 :
600 : !------------------------------------------------------------------------------
601 : ! returns the fraction of particle surface area of aerosol subset `bin_ndx` covered
602 : ! by at least a monolayer of species `species_type` [0-1]
603 : !------------------------------------------------------------------------------
604 9676800 : function coated_frac(self, bin_ndx, species_type, ncol, nlev, aero_props, radius) result(frac)
605 :
606 : class(aerosol_state), intent(in) :: self
607 : integer, intent(in) :: bin_ndx ! bin number
608 : character(len=*), intent(in) :: species_type ! species type
609 : integer, intent(in) :: ncol ! number of columns
610 : integer, intent(in) :: nlev ! number of vertical levels
611 : class(aerosol_properties), intent(in) :: aero_props ! aerosol properties object
612 : real(r8), intent(in) :: radius(:,:) ! m
613 :
614 : real(r8) :: frac(ncol,nlev) ! coated fraction
615 :
616 : !------------coated variables--------------------
617 : real(r8), parameter :: n_so4_monolayers_dust = 1.0_r8 ! number of so4(+nh4) monolayers needed to coat a dust particle
618 : real(r8), parameter :: dr_so4_monolayers_dust = n_so4_monolayers_dust * 4.76e-10_r8
619 19353600 : real(r8) :: vol_shell(ncol,nlev)
620 19353600 : real(r8) :: vol_core(ncol,nlev)
621 : real(r8) :: alnsg, fac_volsfc
622 19353600 : real(r8) :: tmp1(ncol,nlev), tmp2(ncol,nlev)
623 9676800 : real(r8),pointer :: sulf_mmr(:,:)
624 9676800 : real(r8),pointer :: soa_mmr(:,:)
625 9676800 : real(r8),pointer :: pom_mmr(:,:)
626 9676800 : real(r8),pointer :: aer_mmr(:,:)
627 :
628 : integer :: sulf_ndx
629 : integer :: soa_ndx
630 : integer :: pom_ndx
631 : integer :: species_ndx
632 :
633 : real(r8) :: specdens_so4
634 : real(r8) :: specdens_pom
635 : real(r8) :: specdens_soa
636 : real(r8) :: specdens
637 :
638 : character(len=aero_name_len) :: spectype
639 : integer :: ispc
640 :
641 10441267200 : frac = -huge(1._r8)
642 :
643 9676800 : sulf_ndx = -1
644 9676800 : pom_ndx = -1
645 9676800 : soa_ndx = -1
646 9676800 : species_ndx = -1
647 :
648 67737600 : do ispc = 1, aero_props%nspecies(bin_ndx)
649 58060800 : call aero_props%species_type(bin_ndx, ispc, spectype)
650 :
651 116121600 : select case ( trim(spectype) )
652 : case('sulfate')
653 9676800 : sulf_ndx = ispc
654 : case('p-organic')
655 9676800 : pom_ndx = ispc
656 : case('s-organic')
657 116121600 : soa_ndx = ispc
658 : end select
659 67737600 : if (spectype==species_type) then
660 9676800 : species_ndx = ispc
661 : end if
662 : end do
663 :
664 10441267200 : vol_shell(:ncol,:) = 0._r8
665 :
666 9676800 : if (sulf_ndx>0) then
667 9676800 : call aero_props%get(bin_ndx, sulf_ndx, density=specdens_so4)
668 9676800 : call self%get_ambient_mmr(sulf_ndx, bin_ndx, sulf_mmr)
669 10441267200 : vol_shell(:ncol,:) = vol_shell(:ncol,:) + sulf_mmr(:ncol,:)/specdens_so4
670 : end if
671 9676800 : if (pom_ndx>0) then
672 9676800 : call aero_props%get(bin_ndx, pom_ndx, density=specdens_pom)
673 9676800 : call self%get_ambient_mmr(pom_ndx, bin_ndx, pom_mmr)
674 10441267200 : vol_shell(:ncol,:) = vol_shell(:ncol,:) + pom_mmr(:ncol,:)*aero_props%pom_equivso4_factor()/specdens_pom
675 : end if
676 9676800 : if (soa_ndx>0) then
677 9676800 : call aero_props%get(bin_ndx, soa_ndx, density=specdens_soa)
678 9676800 : call self%get_ambient_mmr(soa_ndx, bin_ndx, soa_mmr)
679 10441267200 : vol_shell(:ncol,:) = vol_shell(:ncol,:) + soa_mmr(:ncol,:)*aero_props%soa_equivso4_factor()/specdens_soa
680 : end if
681 :
682 9676800 : call aero_props%get(bin_ndx, species_ndx, density=specdens)
683 9676800 : call self%get_ambient_mmr(species_ndx, bin_ndx, aer_mmr)
684 10441267200 : vol_core(:ncol,:) = aer_mmr(:ncol,:)/specdens
685 :
686 9676800 : alnsg = aero_props%alogsig(bin_ndx)
687 9676800 : fac_volsfc = exp(2.5_r8*alnsg**2)
688 :
689 10441267200 : tmp1(:ncol,:) = vol_shell(:ncol,:)*(radius(:ncol,:)*2._r8)*fac_volsfc
690 10441267200 : tmp2(:ncol,:) = max(6.0_r8*dr_so4_monolayers_dust*vol_core(:ncol,:), 0.0_r8)
691 :
692 10441267200 : where(tmp1(:ncol,:)>0._r8 .and. tmp2(:ncol,:)>0._r8)
693 : frac(:ncol,:) = tmp1(:ncol,:)/tmp2(:ncol,:)
694 : elsewhere
695 : frac(:ncol,:) = 0.001_r8
696 : end where
697 :
698 10441267200 : where(frac(:ncol,:)>1._r8)
699 : frac(:ncol,:) = 1._r8
700 : end where
701 10441267200 : where(frac(:ncol,:) < 0.001_r8)
702 : frac(:ncol,:) = 0.001_r8
703 : end where
704 :
705 9676800 : end function coated_frac
706 :
707 : !------------------------------------------------------------------------------
708 : ! returns the radius [m] of particles in aerosol subset `bin_ndx` assuming all particles are
709 : ! the same size and only species `species_ndx` contributes to the particle volume
710 : !------------------------------------------------------------------------------
711 9676800 : function mass_mean_radius(self, bin_ndx, species_ndx, ncol, nlev, aero_props, rho) result(radius)
712 :
713 : class(aerosol_state), intent(in) :: self
714 : integer, intent(in) :: bin_ndx ! bin number
715 : integer, intent(in) :: species_ndx ! species number
716 : integer, intent(in) :: ncol ! number of columns
717 : integer, intent(in) :: nlev ! number of vertical levels
718 : class(aerosol_properties), intent(in) :: aero_props ! aerosol properties object
719 : real(r8), intent(in) :: rho(:,:) ! air density (kg m-3)
720 :
721 : real(r8) :: radius(ncol,nlev) ! m
722 :
723 : character(len=aero_name_len) :: species_type
724 19353600 : real(r8) :: aer_numdens(ncol,nlev) ! kg/m3
725 19353600 : real(r8) :: aer_massdens(ncol,nlev) ! kg/m3
726 9676800 : real(r8),pointer :: aer_mmr(:,:) ! kg/kg
727 :
728 : real(r8) :: specdens,minrad
729 19353600 : real(r8) :: wght(ncol,nlev)
730 : integer :: i,k
731 :
732 10441267200 : wght = self%hetfrz_size_wght(bin_ndx, ncol, nlev)
733 :
734 9676800 : call aero_props%species_type(bin_ndx, species_ndx, spectype=species_type)
735 :
736 9676800 : call aero_props%get(bin_ndx, species_ndx, density=specdens) ! kg/m3
737 9676800 : call self%get_ambient_mmr(species_ndx, bin_ndx, aer_mmr) ! kg/kg
738 9676800 : call self%get_amb_species_numdens(bin_ndx, ncol, nlev, species_type, aero_props, rho, aer_numdens) ! #/cm3
739 :
740 10441267200 : aer_massdens(:ncol,:) = aer_mmr(:ncol,:)*rho(:ncol,:)*wght(:ncol,:) ! kg/m3
741 :
742 9676800 : minrad = aero_props%min_mass_mean_rad(bin_ndx, species_ndx)
743 :
744 687052800 : do k = 1,nlev
745 10441267200 : do i = 1,ncol
746 10431590400 : if (aer_massdens(i,k)*1.0e-3_r8 > 1.0e-30_r8 .and. aer_numdens(i,k) > 1.0e-3_r8) then
747 507777076 : radius(i,k) = (3._r8/(4*pi*specdens)*aer_massdens(i,k)/(aer_numdens(i,k)*per_m3))**(1._r8/3._r8) ! m
748 : else
749 9246437324 : radius(i,k) = minrad
750 : end if
751 : end do
752 : end do
753 :
754 9676800 : end function mass_mean_radius
755 :
756 : !------------------------------------------------------------------------------
757 : ! calculates water activity mass factor -- density*(1.-(OC+BC)/(OC+BC+SO4)) [mug m-3]
758 : ! of species `species_type` in subset `bin_ndx`
759 : !------------------------------------------------------------------------------
760 9676800 : subroutine watact_mfactor(self, bin_ndx, species_type, ncol, nlev, aero_props, rho, wact_factor)
761 :
762 : class(aerosol_state), intent(in) :: self
763 : integer, intent(in) :: bin_ndx ! bin number
764 : character(len=*), intent(in) :: species_type ! species type
765 : integer, intent(in) :: ncol ! number of columns
766 : integer, intent(in) :: nlev ! number of vertical levels
767 : class(aerosol_properties), intent(in) :: aero_props ! aerosol properties object
768 : real(r8), intent(in) :: rho(:,:) ! air density (kg m-3)
769 : real(r8), intent(out) :: wact_factor(:,:) ! water activity factor -- density*(1.-(OC+BC)/(OC+BC+SO4)) [mug m-3]
770 :
771 9676800 : real(r8), pointer :: aer_mmr(:,:)
772 9676800 : real(r8), pointer :: bin_num(:,:)
773 19353600 : real(r8) :: tot2_mmr(ncol,nlev)
774 19353600 : real(r8) :: tot1_mmr(ncol,nlev)
775 19353600 : real(r8) :: aer_numdens(ncol,nlev)
776 : integer :: ispc
777 : character(len=aero_name_len) :: spectype
778 :
779 19353600 : real(r8) :: awcam(ncol,nlev) ! mass density [mug m-3]
780 19353600 : real(r8) :: awfacm(ncol,nlev) ! mass factor ! (OC+BC)/(OC+BC+SO4)
781 :
782 10441267200 : tot2_mmr = 0.0_r8
783 10441267200 : tot1_mmr = 0.0_r8
784 :
785 9676800 : if (aero_props%soluble(bin_ndx)) then
786 :
787 67737600 : do ispc = 1, aero_props%nspecies(bin_ndx)
788 :
789 58060800 : call aero_props%species_type(bin_ndx, ispc, spectype)
790 :
791 58060800 : if (trim(spectype)=='black-c' .or. trim(spectype)=='p-organic' .or. trim(spectype)=='s-organic') then
792 29030400 : call self%get_ambient_mmr(ispc, bin_ndx, aer_mmr)
793 31323801600 : tot2_mmr(:ncol,:) = tot2_mmr(:ncol,:) + aer_mmr(:ncol,:)
794 : end if
795 67737600 : if (trim(spectype)=='sulfate') then
796 9676800 : call self%get_ambient_mmr(ispc, bin_ndx, aer_mmr)
797 10441267200 : tot1_mmr(:ncol,:) = tot1_mmr(:ncol,:) + aer_mmr(:ncol,:)
798 : end if
799 : end do
800 :
801 : end if
802 :
803 10441267200 : tot1_mmr(:ncol,:) = tot1_mmr(:ncol,:) + tot2_mmr(:ncol,:)
804 :
805 9676800 : call self%get_amb_species_numdens(bin_ndx, ncol, nlev, species_type, aero_props, rho, aer_numdens) ! #/cm3
806 9676800 : call self%get_ambient_num(bin_ndx, bin_num) ! #/kg
807 :
808 10441267200 : where(bin_num(:ncol,:)>0._r8)
809 9676800 : awcam(:ncol,:) = ((aer_numdens(:ncol,:)*per_m3/bin_num(:ncol,:)) * tot1_mmr(:ncol,:)) * kg2mug ! [mug m-3]
810 : elsewhere
811 : awcam(:ncol,:) = 0._r8
812 : end where
813 :
814 10441267200 : where(tot1_mmr(:ncol,:)>0)
815 : awfacm(:ncol,:) = tot2_mmr(:ncol,:) / tot1_mmr(:ncol,:)
816 : elsewhere
817 : awfacm(:ncol,:) = 0._r8
818 : end where
819 :
820 10441267200 : wact_factor(:ncol,:) = awcam(:ncol,:)*(1._r8-awfacm(:ncol,:))
821 :
822 9676800 : end subroutine watact_mfactor
823 :
824 : !------------------------------------------------------------------------------
825 : ! aerosol short wave refactive index
826 : !------------------------------------------------------------------------------
827 0 : function refractive_index_sw(self, ncol, ilev, ilist, ibin, iwav, aero_props) result(crefin)
828 :
829 : class(aerosol_state), intent(in) :: self
830 : integer, intent(in) :: ncol ! number of columes
831 : integer, intent(in) :: ilev ! level index
832 : integer, intent(in) :: ilist ! radiation diagnostics list index
833 : integer, intent(in) :: ibin ! bin index
834 : integer, intent(in) :: iwav ! wave length index
835 : class(aerosol_properties), intent(in) :: aero_props ! aerosol properties object
836 :
837 : complex(r8) :: crefin(ncol) ! complex refractive index
838 :
839 0 : real(r8), pointer :: specmmr(:,:) ! species mass mixing ratio
840 0 : complex(r8), pointer :: specrefindex(:) ! species refractive index
841 : real(r8) :: specdens ! species density (kg/m3)
842 : integer :: ispec, icol
843 0 : real(r8) :: vol(ncol)
844 :
845 0 : crefin(:ncol) = (0._r8, 0._r8)
846 :
847 0 : do ispec = 1, aero_props%nspecies(ilist,ibin)
848 :
849 0 : call self%get_ambient_mmr(ilist,ispec,ibin,specmmr)
850 0 : call aero_props%get(ibin, ispec, list_ndx=ilist, density=specdens, refindex_sw=specrefindex)
851 :
852 0 : do icol = 1, ncol
853 0 : vol(icol) = specmmr(icol,ilev)/specdens
854 0 : crefin(icol) = crefin(icol) + vol(icol)*specrefindex(iwav)
855 : end do
856 : end do
857 :
858 0 : end function refractive_index_sw
859 :
860 : !------------------------------------------------------------------------------
861 : ! aerosol long wave refactive index
862 : !------------------------------------------------------------------------------
863 0 : function refractive_index_lw(self, ncol, ilev, ilist, ibin, iwav, aero_props) result(crefin)
864 :
865 : class(aerosol_state), intent(in) :: self
866 : integer, intent(in) :: ncol ! number of columes
867 : integer, intent(in) :: ilev ! level index
868 : integer, intent(in) :: ilist ! radiation diagnostics list index
869 : integer, intent(in) :: ibin ! bin index
870 : integer, intent(in) :: iwav ! wave length index
871 : class(aerosol_properties), intent(in) :: aero_props ! aerosol properties object
872 :
873 : complex(r8) :: crefin(ncol) ! complex refractive index
874 :
875 0 : real(r8), pointer :: specmmr(:,:) ! species mass mixing ratio
876 0 : complex(r8), pointer :: specrefindex(:) ! species refractive index
877 : real(r8) :: specdens ! species density (kg/m3)
878 : integer :: ispec, icol
879 0 : real(r8) :: vol(ncol)
880 :
881 0 : crefin(:ncol) = (0._r8, 0._r8)
882 :
883 0 : do ispec = 1, aero_props%nspecies(ilist,ibin)
884 :
885 0 : call self%get_ambient_mmr(ilist,ispec,ibin,specmmr)
886 0 : call aero_props%get(ibin, ispec, list_ndx=ilist, density=specdens, refindex_lw=specrefindex)
887 :
888 0 : do icol = 1, ncol
889 0 : vol(icol) = specmmr(icol,ilev)/specdens
890 0 : crefin(icol) = crefin(icol) + vol(icol)*specrefindex(iwav)
891 : end do
892 : end do
893 :
894 0 : end function refractive_index_lw
895 :
896 : !------------------------------------------------------------------------------
897 : ! prescribed aerosol activation fraction for convective cloud
898 : !------------------------------------------------------------------------------
899 22579200 : function convcld_actfrac(self, ibin, ispc, ncol, nlev) result(frac)
900 :
901 : class(aerosol_state), intent(in) :: self
902 : integer, intent(in) :: ibin ! bin index
903 : integer, intent(in) :: ispc ! species index
904 : integer, intent(in) :: ncol ! number of columns
905 : integer, intent(in) :: nlev ! number of vertical levels
906 :
907 : real(r8) :: frac(ncol,nlev)
908 :
909 12181478400 : frac = 0.8_r8 ! rce 2010/05/02
910 :
911 11289600 : end function convcld_actfrac
912 :
913 : !------------------------------------------------------------------------------
914 : ! below cloud solubility factor for interstitial aerosols
915 : !------------------------------------------------------------------------------
916 12902400 : function sol_factb_interstitial(self, bin_ndx, ncol, nlev, aero_props) result(sol_factb)
917 :
918 : class(aerosol_state), intent(in) :: self
919 : integer, intent(in) :: bin_ndx ! bin number
920 : integer, intent(in) :: ncol ! number of columns
921 : integer, intent(in) :: nlev ! number of vertical levels
922 : class(aerosol_properties), intent(in) :: aero_props ! aerosol properties object
923 :
924 : real(r8) :: sol_factb(ncol,nlev)
925 :
926 3225600 : real(r8), pointer :: aer_mmr(:,:)
927 6451200 : real(r8) :: totmmr(ncol,nlev)
928 6451200 : real(r8) :: solmmr(ncol,nlev)
929 : integer :: ispc
930 : character(len=aero_name_len) :: spectype
931 :
932 3480422400 : sol_factb(:,:) = 0.0_r8
933 :
934 3480422400 : totmmr(:,:) = 0._r8
935 3480422400 : solmmr(:,:) = 0._r8
936 :
937 14515200 : do ispc = 1, aero_props%nspecies(bin_ndx)
938 :
939 11289600 : call aero_props%species_type(bin_ndx, ispc, spectype)
940 11289600 : call self%get_ambient_mmr(ispc, bin_ndx, aer_mmr)
941 :
942 12181478400 : totmmr(:ncol,:) = totmmr(:ncol,:) + aer_mmr(:ncol,:)
943 :
944 11289600 : if (trim(spectype) == 'sulfate') then
945 3480422400 : solmmr(:ncol,:) = solmmr(:ncol,:) + aer_mmr(:ncol,:)*0.5_r8
946 : end if
947 11289600 : if (trim(spectype) == 'p-organic') then
948 1740211200 : solmmr(:ncol,:) = solmmr(:ncol,:) + aer_mmr(:ncol,:)*0.2_r8
949 : end if
950 11289600 : if (trim(spectype) == 's-organic') then
951 1740211200 : solmmr(:ncol,:) = solmmr(:ncol,:) + aer_mmr(:ncol,:)*0.2_r8
952 : end if
953 11289600 : if (trim(spectype) == 'dust') then
954 1740211200 : solmmr(:ncol,:) = solmmr(:ncol,:) + aer_mmr(:ncol,:)*0.1_r8
955 : end if
956 14515200 : if (trim(spectype) == 'seasalt') then
957 1740211200 : solmmr(:ncol,:) = solmmr(:ncol,:) + aer_mmr(:ncol,:)*0.8_r8
958 : end if
959 :
960 : end do !nspec
961 :
962 3480422400 : where ( totmmr > 0._r8 )
963 : sol_factb = solmmr/totmmr
964 : end where
965 :
966 3480422400 : where ( sol_factb > 0.8_r8 )
967 : sol_factb = 0.8_r8
968 : end where
969 3480422400 : where ( sol_factb < 0.1_r8 )
970 : sol_factb = 0.1_r8
971 : end where
972 :
973 3225600 : end function sol_factb_interstitial
974 :
975 :
976 0 : end module aerosol_state_mod
|