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