Line data Source code
1 : module modal_aerosol_properties_mod
2 : use shr_kind_mod, only: r8 => shr_kind_r8
3 : use physconst, only: pi
4 : use aerosol_properties_mod, only: aerosol_properties, aero_name_len
5 : use rad_constituents, only: rad_cnst_get_info, rad_cnst_get_mode_props, rad_cnst_get_aer_props
6 :
7 : implicit none
8 :
9 : private
10 :
11 : public :: modal_aerosol_properties
12 :
13 : type, extends(aerosol_properties) :: modal_aerosol_properties
14 : private
15 : real(r8), allocatable :: exp45logsig_(:)
16 : real(r8), allocatable :: voltonumblo_(:)
17 : real(r8), allocatable :: voltonumbhi_(:)
18 : contains
19 : procedure :: number_transported
20 : procedure :: get
21 : procedure :: amcube
22 : procedure :: actfracs
23 : procedure :: num_names
24 : procedure :: mmr_names
25 : procedure :: amb_num_name
26 : procedure :: amb_mmr_name
27 : procedure :: species_type
28 : procedure :: icenuc_updates_num
29 : procedure :: icenuc_updates_mmr
30 : procedure :: apply_number_limits
31 : procedure :: hetfrz_species
32 : procedure :: optics_params
33 : procedure :: nbins_rlist
34 : procedure :: nspecies_per_bin_rlist
35 : procedure :: alogsig_rlist
36 : procedure :: soluble
37 : procedure :: min_mass_mean_rad
38 : procedure :: bin_name
39 :
40 : final :: destructor
41 : end type modal_aerosol_properties
42 :
43 : interface modal_aerosol_properties
44 : procedure :: constructor
45 : end interface modal_aerosol_properties
46 :
47 : contains
48 :
49 : !------------------------------------------------------------------------------
50 : !------------------------------------------------------------------------------
51 4608 : function constructor() result(newobj)
52 :
53 : type(modal_aerosol_properties), pointer :: newobj
54 :
55 : integer :: m, nmodes, ncnst_tot
56 : real(r8) :: dgnumlo
57 : real(r8) :: dgnumhi
58 4608 : integer,allocatable :: nspecies(:)
59 4608 : real(r8),allocatable :: sigmag(:)
60 4608 : real(r8),allocatable :: alogsig(:)
61 4608 : real(r8),allocatable :: f1(:)
62 4608 : real(r8),allocatable :: f2(:)
63 : integer :: ierr
64 :
65 4608 : allocate(newobj,stat=ierr)
66 4608 : if( ierr /= 0 ) then
67 4608 : nullify(newobj)
68 : return
69 : end if
70 :
71 4608 : call rad_cnst_get_info(0, nmodes=nmodes)
72 :
73 13824 : allocate(nspecies(nmodes),stat=ierr)
74 4608 : if( ierr /= 0 ) then
75 4608 : nullify(newobj)
76 : return
77 : end if
78 13824 : allocate(alogsig(nmodes),stat=ierr)
79 4608 : if( ierr /= 0 ) then
80 4608 : nullify(newobj)
81 : return
82 : end if
83 9216 : allocate( f1(nmodes),stat=ierr )
84 4608 : if( ierr /= 0 ) then
85 4608 : nullify(newobj)
86 : return
87 : end if
88 9216 : allocate( f2(nmodes),stat=ierr )
89 4608 : if( ierr /= 0 ) then
90 4608 : nullify(newobj)
91 : return
92 : end if
93 :
94 9216 : allocate(sigmag(nmodes),stat=ierr)
95 4608 : if( ierr /= 0 ) then
96 4608 : nullify(newobj)
97 : return
98 : end if
99 9216 : allocate(newobj%exp45logsig_(nmodes),stat=ierr)
100 4608 : if( ierr /= 0 ) then
101 4608 : nullify(newobj)
102 : return
103 : end if
104 9216 : allocate(newobj%voltonumblo_(nmodes),stat=ierr)
105 4608 : if( ierr /= 0 ) then
106 4608 : nullify(newobj)
107 : return
108 : end if
109 9216 : allocate(newobj%voltonumbhi_(nmodes),stat=ierr)
110 4608 : if( ierr /= 0 ) then
111 4608 : nullify(newobj)
112 : return
113 : end if
114 :
115 4608 : ncnst_tot = 0
116 :
117 23040 : do m = 1, nmodes
118 18432 : call rad_cnst_get_info(0, m, nspec=nspecies(m))
119 :
120 18432 : ncnst_tot = ncnst_tot + nspecies(m) + 1
121 :
122 : call rad_cnst_get_mode_props(0, m, sigmag=sigmag(m), &
123 18432 : dgnumhi=dgnumhi, dgnumlo=dgnumlo )
124 :
125 18432 : alogsig(m) = log(sigmag(m))
126 :
127 18432 : newobj%exp45logsig_(m) = exp(4.5_r8*alogsig(m)*alogsig(m))
128 :
129 18432 : f1(m) = 0.5_r8*exp(2.5_r8*alogsig(m)*alogsig(m))
130 18432 : f2(m) = 1._r8 + 0.25_r8*alogsig(m)
131 :
132 0 : newobj%voltonumblo_(m) = 1._r8 / ( (pi/6._r8)* &
133 18432 : (dgnumlo**3._r8)*exp(4.5_r8*alogsig(m)**2._r8) )
134 0 : newobj%voltonumbhi_(m) = 1._r8 / ( (pi/6._r8)* &
135 23040 : (dgnumhi**3._r8)*exp(4.5_r8*alogsig(m)**2._r8) )
136 :
137 : end do
138 :
139 4608 : call newobj%initialize(nmodes,ncnst_tot,nspecies,nspecies,alogsig,f1,f2,ierr)
140 4608 : if( ierr /= 0 ) then
141 4608 : nullify(newobj)
142 : return
143 : end if
144 4608 : deallocate(nspecies)
145 4608 : deallocate(alogsig)
146 4608 : deallocate(sigmag)
147 4608 : deallocate(f1)
148 4608 : deallocate(f2)
149 :
150 4608 : end function constructor
151 :
152 : !------------------------------------------------------------------------------
153 : !------------------------------------------------------------------------------
154 1536 : subroutine destructor(self)
155 : type(modal_aerosol_properties), intent(inout) :: self
156 :
157 1536 : if (allocated(self%exp45logsig_)) then
158 1536 : deallocate(self%exp45logsig_)
159 : end if
160 1536 : if (allocated(self%voltonumblo_)) then
161 1536 : deallocate(self%voltonumblo_)
162 : end if
163 1536 : if (allocated(self%voltonumbhi_)) then
164 1536 : deallocate(self%voltonumbhi_)
165 : end if
166 :
167 1536 : call self%final()
168 :
169 1536 : end subroutine destructor
170 :
171 : !------------------------------------------------------------------------------
172 : ! returns number of transported aerosol constituents
173 : !------------------------------------------------------------------------------
174 1536 : integer function number_transported(self)
175 : class(modal_aerosol_properties), intent(in) :: self
176 : ! to be implemented later
177 1536 : number_transported = -1
178 1536 : end function number_transported
179 :
180 : !------------------------------------------------------------------------
181 : ! returns aerosol properties:
182 : ! density
183 : ! hygroscopicity
184 : ! species type
185 : ! short wave species refractive indices
186 : ! long wave species refractive indices
187 : ! species morphology
188 : !------------------------------------------------------------------------
189 2491704957 : subroutine get(self, bin_ndx, species_ndx, list_ndx, density, hygro, &
190 : spectype, specmorph, refindex_sw, refindex_lw)
191 :
192 : class(modal_aerosol_properties), intent(in) :: self
193 : integer, intent(in) :: bin_ndx ! bin index
194 : integer, intent(in) :: species_ndx ! species index
195 : integer, optional, intent(in) :: list_ndx ! climate or a diagnostic list number
196 : real(r8), optional, intent(out) :: density ! density (kg/m3)
197 : real(r8), optional, intent(out) :: hygro ! hygroscopicity
198 : character(len=*), optional, intent(out) :: spectype ! species type
199 : character(len=*), optional, intent(out) :: specmorph ! species morphology
200 : complex(r8), pointer, optional, intent(out) :: refindex_sw(:) ! short wave species refractive indices
201 : complex(r8), pointer, optional, intent(out) :: refindex_lw(:) ! long wave species refractive indices
202 :
203 : integer :: ilist
204 :
205 2491704957 : if (present(list_ndx)) then
206 2060488800 : ilist = list_ndx
207 : else
208 431216157 : ilist = 0
209 : end if
210 :
211 : call rad_cnst_get_aer_props(ilist, bin_ndx, species_ndx, &
212 : density_aer=density, hygro_aer=hygro, spectype=spectype, &
213 7315203828 : refindex_aer_sw=refindex_sw, refindex_aer_lw=refindex_lw)
214 :
215 2491704957 : if (present(specmorph)) then
216 0 : specmorph = 'UNKNOWN'
217 : end if
218 :
219 2491704957 : end subroutine get
220 :
221 : !------------------------------------------------------------------------
222 : ! returns optics type and table parameters
223 : !------------------------------------------------------------------------
224 495360 : subroutine optics_params(self, list_ndx, bin_ndx, opticstype, extpsw, abspsw, asmpsw, absplw, &
225 : refrtabsw, refitabsw, refrtablw, refitablw, ncoef, prefr, prefi, sw_hygro_ext_wtp, &
226 : sw_hygro_ssa_wtp, sw_hygro_asm_wtp, lw_hygro_ext_wtp, wgtpct, nwtp, &
227 : sw_hygro_coreshell_ext, sw_hygro_coreshell_ssa, sw_hygro_coreshell_asm, lw_hygro_coreshell_ext, &
228 : corefrac, bcdust, kap, relh, nfrac, nbcdust, nkap, nrelh )
229 :
230 : class(modal_aerosol_properties), intent(in) :: self
231 : integer, intent(in) :: bin_ndx ! bin index
232 : integer, intent(in) :: list_ndx ! rad climate/diags list
233 :
234 : character(len=*), optional, intent(out) :: opticstype
235 :
236 : ! refactive index table parameters
237 : real(r8), optional, pointer :: extpsw(:,:,:,:) ! short wave specific extinction
238 : real(r8), optional, pointer :: abspsw(:,:,:,:) ! short wave specific absorption
239 : real(r8), optional, pointer :: asmpsw(:,:,:,:) ! short wave asymmetry factor
240 : real(r8), optional, pointer :: absplw(:,:,:,:) ! long wave specific absorption
241 : real(r8), optional, pointer :: refrtabsw(:,:) ! table of short wave real refractive indices for aerosols
242 : real(r8), optional, pointer :: refitabsw(:,:) ! table of short wave imaginary refractive indices for aerosols
243 : real(r8), optional, pointer :: refrtablw(:,:) ! table of long wave real refractive indices for aerosols
244 : real(r8), optional, pointer :: refitablw(:,:) ! table of long wave imaginary refractive indices for aerosols
245 : integer, optional, intent(out) :: ncoef ! number of chebychev polynomials
246 : integer, optional, intent(out) :: prefr ! number of real refractive indices in table
247 : integer, optional, intent(out) :: prefi ! number of imaginary refractive indices in table
248 :
249 : ! hygrowghtpct table parameters
250 : real(r8), optional, pointer :: sw_hygro_ext_wtp(:,:) ! short wave extinction table
251 : real(r8), optional, pointer :: sw_hygro_ssa_wtp(:,:) ! short wave single-scatter albedo table
252 : real(r8), optional, pointer :: sw_hygro_asm_wtp(:,:) ! short wave asymmetry table
253 : real(r8), optional, pointer :: lw_hygro_ext_wtp(:,:) ! long wave absorption table
254 : real(r8), optional, pointer :: wgtpct(:) ! weight precent of H2SO4/H2O solution
255 : integer, optional, intent(out) :: nwtp ! number of weight precent values
256 :
257 : ! hygrocoreshell table parameters
258 : real(r8), optional, pointer :: sw_hygro_coreshell_ext(:,:,:,:,:) ! short wave extinction table
259 : real(r8), optional, pointer :: sw_hygro_coreshell_ssa(:,:,:,:,:) ! short wave single-scatter albedo table
260 : real(r8), optional, pointer :: sw_hygro_coreshell_asm(:,:,:,:,:) ! short wave asymmetry table
261 : real(r8), optional, pointer :: lw_hygro_coreshell_ext(:,:,:,:,:) ! long wave absorption table
262 : real(r8), optional, pointer :: corefrac(:) ! core fraction dimension values
263 : real(r8), optional, pointer :: bcdust(:) ! bc/(bc + dust) fraction dimension values
264 : real(r8), optional, pointer :: kap(:) ! hygroscopicity dimension values
265 : real(r8), optional, pointer :: relh(:) ! relative humidity dimension values
266 : integer, optional, intent(out) :: nfrac ! core fraction dimension size
267 : integer, optional, intent(out) :: nbcdust ! bc/(bc + dust) fraction dimension size
268 : integer, optional, intent(out) :: nkap ! hygroscopicity dimension size
269 : integer, optional, intent(out) :: nrelh ! relative humidity dimension size
270 :
271 : ! refactive index table parameters
272 : call rad_cnst_get_mode_props(list_ndx, bin_ndx, &
273 : opticstype=opticstype, &
274 : extpsw=extpsw, &
275 : abspsw=abspsw, &
276 : asmpsw=asmpsw, &
277 : absplw=absplw, &
278 : refrtabsw=refrtabsw, &
279 : refitabsw=refitabsw, &
280 : refrtablw=refrtablw, &
281 : refitablw=refitablw, &
282 : ncoef=ncoef, &
283 : prefr=prefr, &
284 2724480 : prefi=prefi)
285 :
286 : ! hygrowghtpct table parameters
287 495360 : if (present(sw_hygro_ext_wtp)) then
288 0 : nullify(sw_hygro_ext_wtp)
289 : end if
290 495360 : if (present(sw_hygro_ssa_wtp)) then
291 0 : nullify(sw_hygro_ssa_wtp)
292 : end if
293 495360 : if (present(sw_hygro_asm_wtp)) then
294 0 : nullify(sw_hygro_asm_wtp)
295 : end if
296 495360 : if (present(lw_hygro_ext_wtp)) then
297 0 : nullify(lw_hygro_ext_wtp)
298 : end if
299 495360 : if (present(wgtpct)) then
300 0 : nullify(wgtpct)
301 : end if
302 495360 : if (present(nwtp)) then
303 0 : nwtp = -1
304 : end if
305 :
306 : ! hygrocoreshell table parameters
307 495360 : if (present(sw_hygro_coreshell_ext)) then
308 0 : nullify(sw_hygro_coreshell_ext)
309 : end if
310 495360 : if (present(sw_hygro_coreshell_ssa)) then
311 0 : nullify(sw_hygro_coreshell_ssa)
312 : end if
313 495360 : if (present(sw_hygro_coreshell_asm)) then
314 0 : nullify(sw_hygro_coreshell_asm)
315 : end if
316 495360 : if (present(lw_hygro_coreshell_ext)) then
317 0 : nullify(lw_hygro_coreshell_ext)
318 : end if
319 495360 : if (present(corefrac)) then
320 0 : nullify(corefrac)
321 : end if
322 495360 : if (present(bcdust)) then
323 0 : nullify(bcdust)
324 : end if
325 495360 : if (present(kap)) then
326 0 : nullify(kap)
327 : end if
328 495360 : if (present(relh)) then
329 0 : nullify(relh)
330 : end if
331 495360 : if (present(nfrac)) then
332 0 : nfrac = -1
333 : end if
334 495360 : if (present(nbcdust)) then
335 0 : nbcdust = -1
336 : end if
337 495360 : if (present(nkap)) then
338 0 : nkap = -1
339 : end if
340 495360 : if (present(nrelh)) then
341 0 : nrelh = -1
342 : end if
343 :
344 495360 : end subroutine optics_params
345 :
346 : !------------------------------------------------------------------------------
347 : ! returns radius^3 (m3) of a given bin number
348 : !------------------------------------------------------------------------------
349 900892578 : pure elemental real(r8) function amcube(self, bin_ndx, volconc, numconc)
350 :
351 : class(modal_aerosol_properties), intent(in) :: self
352 : integer, intent(in) :: bin_ndx ! bin number
353 : real(r8), intent(in) :: volconc ! volume conc (m3/m3)
354 : real(r8), intent(in) :: numconc ! number conc (1/m3)
355 :
356 900892578 : amcube = (3._r8*volconc/(4._r8*pi*self%exp45logsig_(bin_ndx)*numconc))
357 :
358 900892578 : end function amcube
359 :
360 : !------------------------------------------------------------------------------
361 : ! returns mass and number activation fractions
362 : !------------------------------------------------------------------------------
363 62887772 : subroutine actfracs(self, bin_ndx, smc, smax, fn, fm )
364 : use shr_spfn_mod, only: erf => shr_spfn_erf
365 : class(modal_aerosol_properties), intent(in) :: self
366 : integer, intent(in) :: bin_ndx ! bin index
367 : real(r8),intent(in) :: smc ! critical supersaturation for particles of bin radius
368 : real(r8),intent(in) :: smax ! maximum supersaturation for multiple competing aerosols
369 : real(r8),intent(out) :: fn ! activation fraction for aerosol number
370 : real(r8),intent(out) :: fm ! activation fraction for aerosol mass
371 :
372 : real(r8) :: x,y
373 : real(r8), parameter :: twothird = 2._r8/3._r8
374 : real(r8), parameter :: sq2 = sqrt(2._r8)
375 :
376 62887772 : x=twothird*(log(smc)-log(smax))/(sq2*self%alogsig(bin_ndx))
377 62887772 : y=x-1.5_r8*sq2*self%alogsig(bin_ndx)
378 :
379 62887772 : fn = 0.5_r8*(1._r8-erf(x))
380 62887772 : fm = 0.5_r8*(1._r8-erf(y))
381 :
382 62887772 : end subroutine actfracs
383 :
384 : !------------------------------------------------------------------------
385 : ! returns constituents names of aerosol number mixing ratios
386 : !------------------------------------------------------------------------
387 6144 : subroutine num_names(self, bin_ndx, name_a, name_c)
388 : class(modal_aerosol_properties), intent(in) :: self
389 : integer, intent(in) :: bin_ndx ! bin number
390 : character(len=*), intent(out) :: name_a ! constituent name of ambient aerosol number dens
391 : character(len=*), intent(out) :: name_c ! constituent name of cloud-borne aerosol number dens
392 :
393 6144 : call rad_cnst_get_info(0,bin_ndx, num_name=name_a, num_name_cw=name_c)
394 6144 : end subroutine num_names
395 :
396 : !------------------------------------------------------------------------
397 : ! returns constituents names of aerosol mass mixing ratios
398 : !------------------------------------------------------------------------
399 23040 : subroutine mmr_names(self, bin_ndx, species_ndx, name_a, name_c)
400 : class(modal_aerosol_properties), intent(in) :: self
401 : integer, intent(in) :: bin_ndx ! bin number
402 : integer, intent(in) :: species_ndx ! species number
403 : character(len=*), intent(out) :: name_a ! constituent name of ambient aerosol MMR
404 : character(len=*), intent(out) :: name_c ! constituent name of cloud-borne aerosol MMR
405 :
406 23040 : call rad_cnst_get_info(0, bin_ndx, species_ndx, spec_name=name_a, spec_name_cw=name_c)
407 23040 : end subroutine mmr_names
408 :
409 : !------------------------------------------------------------------------
410 : ! returns constituent name of ambient aerosol number mixing ratios
411 : !------------------------------------------------------------------------
412 1536 : subroutine amb_num_name(self, bin_ndx, name)
413 : class(modal_aerosol_properties), intent(in) :: self
414 : integer, intent(in) :: bin_ndx ! bin number
415 : character(len=*), intent(out) :: name ! constituent name of ambient aerosol number dens
416 :
417 1536 : call rad_cnst_get_info(0,bin_ndx, num_name=name)
418 :
419 1536 : end subroutine amb_num_name
420 :
421 : !------------------------------------------------------------------------
422 : ! returns constituent name of ambient aerosol mass mixing ratios
423 : !------------------------------------------------------------------------
424 1536 : subroutine amb_mmr_name(self, bin_ndx, species_ndx, name)
425 : class(modal_aerosol_properties), intent(in) :: self
426 : integer, intent(in) :: bin_ndx ! bin number
427 : integer, intent(in) :: species_ndx ! species number
428 : character(len=*), intent(out) :: name ! constituent name of ambient aerosol MMR
429 :
430 1536 : call rad_cnst_get_info(0, bin_ndx, species_ndx, spec_name=name)
431 :
432 1536 : end subroutine amb_mmr_name
433 :
434 : !------------------------------------------------------------------------
435 : ! returns species type
436 : !------------------------------------------------------------------------
437 1287742116 : subroutine species_type(self, bin_ndx, species_ndx, spectype)
438 : class(modal_aerosol_properties), intent(in) :: self
439 : integer, intent(in) :: bin_ndx ! bin number
440 : integer, intent(in) :: species_ndx ! species number
441 : character(len=*), intent(out) :: spectype ! species type
442 :
443 1287742116 : call rad_cnst_get_info(0, bin_ndx, species_ndx, spec_type=spectype)
444 :
445 1287742116 : end subroutine species_type
446 :
447 : !------------------------------------------------------------------------------
448 : ! returns TRUE if Ice Nucleation tendencies are applied to given aerosol bin number
449 : !------------------------------------------------------------------------------
450 721813776 : function icenuc_updates_num(self, bin_ndx) result(res)
451 : class(modal_aerosol_properties), intent(in) :: self
452 : integer, intent(in) :: bin_ndx ! bin number
453 :
454 : logical :: res
455 :
456 : character(len=aero_name_len) :: spectype
457 : character(len=aero_name_len) :: modetype
458 : integer :: spc_ndx
459 :
460 721813776 : res = .false.
461 :
462 721813776 : call rad_cnst_get_info(0, bin_ndx, mode_type=modetype)
463 721813776 : if (.not.(modetype=='coarse' .or. modetype=='coarse_dust')) then
464 : return
465 : end if
466 :
467 721813776 : do spc_ndx = 1, self%nspecies(bin_ndx)
468 541360332 : call self%species_type( bin_ndx, spc_ndx, spectype)
469 721813776 : if (spectype=='dust') res = .true.
470 : end do
471 :
472 180453444 : end function icenuc_updates_num
473 :
474 : !------------------------------------------------------------------------------
475 : ! returns TRUE if Ice Nucleation tendencies are applied to a given species within a bin
476 : !------------------------------------------------------------------------------
477 541361868 : function icenuc_updates_mmr(self, bin_ndx, species_ndx) result(res)
478 : class(modal_aerosol_properties), intent(in) :: self
479 : integer, intent(in) :: bin_ndx ! bin number
480 : integer, intent(in) :: species_ndx ! species number
481 :
482 : logical :: res
483 :
484 : character(len=32) :: spectype
485 : character(len=32) :: modetype
486 :
487 541361868 : res = .false.
488 :
489 541361868 : if (species_ndx>0) then
490 :
491 541360332 : call rad_cnst_get_info(0, bin_ndx, mode_type=modetype)
492 541360332 : if (.not.(modetype=='coarse' .or. modetype=='coarse_dust')) then
493 : return
494 : end if
495 :
496 541360332 : call self%species_type( bin_ndx, species_ndx, spectype)
497 541360332 : if (spectype=='dust') res = .true.
498 : end if
499 :
500 : end function icenuc_updates_mmr
501 :
502 : !------------------------------------------------------------------------------
503 : ! apply max / min to number concentration
504 : !------------------------------------------------------------------------------
505 114238028 : subroutine apply_number_limits( self, naerosol, vaerosol, istart, istop, m )
506 : class(modal_aerosol_properties), intent(in) :: self
507 : real(r8), intent(inout) :: naerosol(:) ! number conc (1/m3)
508 : real(r8), intent(in) :: vaerosol(:) ! volume conc (m3/m3)
509 : integer, intent(in) :: istart ! start column index (1 <= istart <= istop <= pcols)
510 : integer, intent(in) :: istop ! stop column index
511 : integer, intent(in) :: m ! mode or bin index
512 :
513 : integer :: i
514 :
515 : ! adjust number so that dgnumlo < dgnum < dgnumhi
516 : ! -- the diameter falls within the lower and upper limits which are
517 : ! represented by voltonumhi and voltonumblo values, respectively
518 1099968664 : do i = istart, istop
519 985730636 : naerosol(i) = max(naerosol(i), vaerosol(i)*self%voltonumbhi_(m))
520 1099968664 : naerosol(i) = min(naerosol(i), vaerosol(i)*self%voltonumblo_(m))
521 : end do
522 :
523 114238028 : end subroutine apply_number_limits
524 :
525 : !------------------------------------------------------------------------------
526 : ! returns TRUE if species `spc_ndx` in aerosol subset `bin_ndx` contributes to
527 : ! the particles' ability to act as heterogeneous freezing nuclei
528 : !------------------------------------------------------------------------------
529 46080 : function hetfrz_species(self, bin_ndx, spc_ndx) result(res)
530 : class(modal_aerosol_properties), intent(in) :: self
531 : integer, intent(in) :: bin_ndx ! bin number
532 : integer, intent(in) :: spc_ndx ! species number
533 :
534 : logical :: res
535 :
536 : character(len=aero_name_len) :: mode_name, species_type
537 :
538 46080 : res = .false.
539 :
540 46080 : call rad_cnst_get_info(0, bin_ndx, mode_type=mode_name)
541 :
542 46080 : if ((trim(mode_name)/='aitken')) then
543 :
544 33792 : call self%species_type(bin_ndx, spc_ndx, species_type)
545 :
546 33792 : if ((trim(species_type)=='black-c').or.(trim(species_type)=='dust')) then
547 :
548 12288 : res = .true.
549 :
550 : end if
551 :
552 : end if
553 :
554 46080 : end function hetfrz_species
555 :
556 : !------------------------------------------------------------------------------
557 : ! returns TRUE if soluble
558 : !------------------------------------------------------------------------------
559 705888 : logical function soluble(self,bin_ndx)
560 : class(modal_aerosol_properties), intent(in) :: self
561 : integer, intent(in) :: bin_ndx ! bin number
562 :
563 : character(len=aero_name_len) :: mode_name
564 :
565 705888 : call rad_cnst_get_info(0, bin_ndx, mode_type=mode_name)
566 :
567 705888 : soluble = trim(mode_name)/='primary_carbon'
568 :
569 705888 : end function soluble
570 :
571 : !------------------------------------------------------------------------------
572 : ! returns minimum mass mean radius (meters)
573 : !------------------------------------------------------------------------------
574 705888 : function min_mass_mean_rad(self,bin_ndx,species_ndx) result(minrad)
575 : class(modal_aerosol_properties), intent(in) :: self
576 : integer, intent(in) :: bin_ndx ! bin number
577 : integer, intent(in) :: species_ndx ! species number
578 :
579 : real(r8) :: minrad ! meters
580 :
581 : integer :: nmodes
582 : character(len=aero_name_len) :: species_type, mode_type
583 :
584 705888 : call self%species_type(bin_ndx, species_ndx, spectype=species_type)
585 1411776 : select case ( trim(species_type) )
586 : case('dust')
587 352944 : call rad_cnst_get_info(0, bin_ndx, mode_type=mode_type)
588 1058832 : select case ( trim(mode_type) )
589 : case ('accum','fine_dust')
590 176472 : minrad = 0.258e-6_r8
591 : case ('coarse','coarse_dust')
592 176472 : minrad = 1.576e-6_r8
593 : case default
594 705888 : minrad = -huge(1._r8)
595 : end select
596 : case('black-c')
597 352944 : call rad_cnst_get_info(0, nmodes=nmodes)
598 352944 : if (nmodes==3) then
599 : minrad = 0.04e-6_r8
600 : else
601 352944 : minrad = 0.067e-6_r8 ! from emission size
602 : endif
603 : case default
604 1411776 : minrad = -huge(1._r8)
605 : end select
606 :
607 705888 : end function min_mass_mean_rad
608 :
609 : !------------------------------------------------------------------------------
610 : ! returns the total number of bins for a given radiation list index
611 : !------------------------------------------------------------------------------
612 61920 : function nbins_rlist(self, list_ndx) result(res)
613 : class(modal_aerosol_properties), intent(in) :: self
614 : integer, intent(in) :: list_ndx ! radiation list number
615 :
616 : integer :: res
617 :
618 61920 : call rad_cnst_get_info(list_ndx, nmodes=res)
619 :
620 61920 : end function nbins_rlist
621 :
622 : !------------------------------------------------------------------------------
623 : ! returns number of species in a bin for a given radiation list index
624 : !------------------------------------------------------------------------------
625 526677120 : function nspecies_per_bin_rlist(self, list_ndx, bin_ndx) result(res)
626 : class(modal_aerosol_properties), intent(in) :: self
627 : integer, intent(in) :: list_ndx ! radiation list number
628 : integer, intent(in) :: bin_ndx ! bin number
629 :
630 : integer :: res
631 :
632 526677120 : call rad_cnst_get_info(list_ndx, bin_ndx, nspec=res)
633 :
634 526677120 : end function nspecies_per_bin_rlist
635 :
636 : !------------------------------------------------------------------------------
637 : ! returns the natural log of geometric standard deviation of the number
638 : ! distribution for radiation list number and aerosol bin
639 : !------------------------------------------------------------------------------
640 247680 : function alogsig_rlist(self, list_ndx, bin_ndx) result(res)
641 : class(modal_aerosol_properties), intent(in) :: self
642 : integer, intent(in) :: list_ndx ! radiation list number
643 : integer, intent(in) :: bin_ndx ! bin number
644 :
645 : real(r8) :: res
646 :
647 : real(r8) :: sig
648 :
649 247680 : call rad_cnst_get_mode_props(list_ndx, bin_ndx, sigmag=sig)
650 247680 : res = log(sig)
651 :
652 247680 : end function alogsig_rlist
653 :
654 : !------------------------------------------------------------------------------
655 : ! returns name for a given radiation list number and aerosol bin
656 : !------------------------------------------------------------------------------
657 24576 : function bin_name(self, list_ndx, bin_ndx) result(name)
658 : class(modal_aerosol_properties), intent(in) :: self
659 : integer, intent(in) :: list_ndx ! radiation list number
660 : integer, intent(in) :: bin_ndx ! bin number
661 :
662 : character(len=32) name
663 :
664 24576 : call rad_cnst_get_info(list_ndx, bin_ndx, mode_type=name)
665 :
666 24576 : end function bin_name
667 :
668 4608 : end module modal_aerosol_properties_mod
|