Line data Source code
1 : module aerosol_properties_mod
2 : use shr_kind_mod, only: r8 => shr_kind_r8
3 :
4 : implicit none
5 :
6 : private
7 :
8 : public :: aerosol_properties
9 :
10 : !> aerosol_properties defines the configuration of any aerosol package (using
11 : !! any aerosol representation) based on user specification. These values are
12 : !! set during initialization and do not vary during the simulation.
13 : !!
14 : !! Each aerosol package (e.g., MAM, CARMA, etc) must extend the abstract
15 : !! aerosol_properties class to define the details of their configuration. Any
16 : !! package must implement each of the deferred procedures of the abstract
17 : !! aerosol_properties class, may include additional private data members and
18 : !! type-bound procedures, and may override functions of the abstract class.
19 : !!
20 : !! Please see the modal_aerosol_properties module for an example of how the
21 : !! aerosol_properties class can be extended for a specific aerosol package.
22 : type, abstract :: aerosol_properties
23 : private
24 : integer :: nbins_ = 0 ! number of aerosol bins
25 : integer :: ncnst_tot_ = 0 ! total number of constituents
26 : integer, allocatable :: nmasses_(:) ! number of species masses
27 : integer, allocatable :: nspecies_(:) ! number of species
28 : integer, allocatable :: indexer_(:,:) ! unique indices of the aerosol elements
29 : real(r8), allocatable :: alogsig_(:) ! natural log of geometric deviation of the number distribution for aerosol bin
30 : real(r8), allocatable :: f1_(:) ! eq 28 Abdul-Razzak et al 1998
31 : real(r8), allocatable :: f2_(:) ! eq 29 Abdul-Razzak et al 1998
32 : ! Abdul-Razzak, H., S.J. Ghan, and C. Rivera-Carpio, A parameterization of aerosol activation,
33 : ! 1, Singleaerosoltype. J. Geophys. Res., 103, 6123-6132, 1998.
34 : real(r8) :: soa_equivso4_factor_ = -huge(1._r8)
35 : real(r8) :: pom_equivso4_factor_ = -huge(1._r8)
36 : contains
37 : procedure :: initialize => aero_props_init
38 : procedure,private :: nbins_0list
39 : procedure(aero_nbins_rlist), deferred :: nbins_rlist
40 : generic :: nbins => nbins_0list,nbins_rlist
41 : procedure :: ncnst_tot
42 : procedure,private :: nspecies_per_bin
43 : procedure(aero_nspecies_rlist), deferred :: nspecies_per_bin_rlist
44 : procedure,private :: nspecies_all_bins
45 : generic :: nspecies => nspecies_all_bins,nspecies_per_bin,nspecies_per_bin_rlist
46 : procedure,private :: n_masses_all_bins
47 : procedure,private :: n_masses_per_bin
48 : generic :: nmasses => n_masses_all_bins,n_masses_per_bin
49 : procedure :: indexer
50 : procedure :: maxsat
51 : procedure(aero_amcube), deferred :: amcube
52 : procedure :: alogsig_0list
53 : procedure(aero_alogsig_rlist), deferred :: alogsig_rlist
54 : generic :: alogsig => alogsig_0list,alogsig_rlist
55 : procedure(aero_number_transported), deferred :: number_transported
56 : procedure(aero_props_get), deferred :: get
57 : procedure(aero_actfracs), deferred :: actfracs
58 : procedure(aero_num_names), deferred :: num_names
59 : procedure(aero_mmr_names), deferred :: mmr_names
60 : procedure(aero_amb_num_name), deferred :: amb_num_name
61 : procedure(aero_amb_mmr_name), deferred :: amb_mmr_name
62 : procedure(aero_species_type), deferred :: species_type
63 : procedure(aero_icenuc_updates_num), deferred :: icenuc_updates_num
64 : procedure(aero_icenuc_updates_mmr), deferred :: icenuc_updates_mmr
65 : procedure(aero_apply_num_limits), deferred :: apply_number_limits
66 : procedure(aero_hetfrz_species), deferred :: hetfrz_species
67 : procedure :: soa_equivso4_factor ! SOA Hygroscopicity / Sulfate Hygroscopicity
68 : procedure :: pom_equivso4_factor ! POM Hygroscopicity / Sulfate Hygroscopicity
69 : procedure(aero_soluble), deferred :: soluble
70 : procedure(aero_min_mass_mean_rad), deferred :: min_mass_mean_rad
71 : procedure(aero_optics_params), deferred :: optics_params
72 : procedure(aero_bin_name), deferred :: bin_name
73 : procedure(aero_scav_diam), deferred :: scav_diam
74 : procedure(aero_resuspension_resize), deferred :: resuspension_resize
75 : procedure(aero_rebin_bulk_fluxes), deferred :: rebin_bulk_fluxes
76 : procedure(aero_hydrophilic), deferred :: hydrophilic
77 :
78 : procedure :: final=>aero_props_final
79 : end type aerosol_properties
80 :
81 : integer,public, parameter :: aero_name_len = 32 ! common length of aersols names, species, etc
82 :
83 : abstract interface
84 :
85 : !------------------------------------------------------------------------------
86 : ! returns number of transported aerosol constituents
87 : !------------------------------------------------------------------------------
88 : integer function aero_number_transported(self)
89 : import :: aerosol_properties
90 : class(aerosol_properties), intent(in) :: self
91 : end function aero_number_transported
92 :
93 : !------------------------------------------------------------------------
94 : ! returns aerosol properties:
95 : ! density
96 : ! hygroscopicity
97 : ! species type
98 : ! species name
99 : ! short wave species refractive indices
100 : ! long wave species refractive indices
101 : ! species morphology
102 : !------------------------------------------------------------------------
103 : subroutine aero_props_get(self, bin_ndx, species_ndx, list_ndx, density, hygro, &
104 : spectype, specname, specmorph, refindex_sw, refindex_lw)
105 : import :: aerosol_properties, r8
106 : class(aerosol_properties), intent(in) :: self
107 : integer, intent(in) :: bin_ndx ! bin index
108 : integer, intent(in) :: species_ndx ! species index
109 : integer, optional, intent(in) :: list_ndx ! climate or a diagnostic list number
110 : real(r8), optional, intent(out) :: density ! density (kg/m3)
111 : real(r8), optional, intent(out) :: hygro ! hygroscopicity
112 : character(len=*), optional, intent(out) :: spectype ! species type
113 : character(len=*), optional, intent(out) :: specname ! species name
114 : character(len=*), optional, intent(out) :: specmorph ! species morphology
115 : complex(r8), pointer, optional, intent(out) :: refindex_sw(:) ! short wave species refractive indices
116 : complex(r8), pointer, optional, intent(out) :: refindex_lw(:) ! long wave species refractive indices
117 :
118 : end subroutine aero_props_get
119 :
120 : !------------------------------------------------------------------------
121 : ! returns optics type and table parameters
122 : !------------------------------------------------------------------------
123 : subroutine aero_optics_params(self, list_ndx, bin_ndx, opticstype, extpsw, abspsw, asmpsw, absplw, &
124 : refrtabsw, refitabsw, refrtablw, refitablw, ncoef, prefr, prefi, sw_hygro_ext_wtp, &
125 : sw_hygro_ssa_wtp, sw_hygro_asm_wtp, lw_hygro_ext_wtp, wgtpct, nwtp, &
126 : sw_hygro_coreshell_ext, sw_hygro_coreshell_ssa, sw_hygro_coreshell_asm, lw_hygro_coreshell_ext, &
127 : corefrac, bcdust, kap, relh, nfrac, nbcdust, nkap, nrelh )
128 :
129 : import :: aerosol_properties, r8
130 :
131 : class(aerosol_properties), intent(in) :: self
132 : integer, intent(in) :: bin_ndx ! bin index
133 : integer, intent(in) :: list_ndx ! rad climate/diags list
134 :
135 : character(len=*), optional, intent(out) :: opticstype
136 :
137 : ! refactive index table parameters
138 : real(r8), optional, pointer :: extpsw(:,:,:,:) ! short wave specific extinction
139 : real(r8), optional, pointer :: abspsw(:,:,:,:) ! short wave specific absorption
140 : real(r8), optional, pointer :: asmpsw(:,:,:,:) ! short wave asymmetry factor
141 : real(r8), optional, pointer :: absplw(:,:,:,:) ! long wave specific absorption
142 : real(r8), optional, pointer :: refrtabsw(:,:) ! table of short wave real refractive indices for aerosols
143 : real(r8), optional, pointer :: refitabsw(:,:) ! table of short wave imaginary refractive indices for aerosols
144 : real(r8), optional, pointer :: refrtablw(:,:) ! table of long wave real refractive indices for aerosols
145 : real(r8), optional, pointer :: refitablw(:,:) ! table of long wave imaginary refractive indices for aerosols
146 : integer, optional, intent(out) :: ncoef ! number of chebychev polynomials
147 : integer, optional, intent(out) :: prefr ! number of real refractive indices in table
148 : integer, optional, intent(out) :: prefi ! number of imaginary refractive indices in table
149 :
150 : ! hygrowghtpct table parameters
151 : real(r8), optional, pointer :: sw_hygro_ext_wtp(:,:) ! short wave extinction table
152 : real(r8), optional, pointer :: sw_hygro_ssa_wtp(:,:) ! short wave single-scatter albedo table
153 : real(r8), optional, pointer :: sw_hygro_asm_wtp(:,:) ! short wave asymmetry table
154 : real(r8), optional, pointer :: lw_hygro_ext_wtp(:,:) ! long wave absorption table
155 : real(r8), optional, pointer :: wgtpct(:) ! weight precent of H2SO4/H2O solution
156 : integer, optional, intent(out) :: nwtp ! number of weight precent values
157 :
158 : ! hygrocoreshell table parameters
159 : real(r8), optional, pointer :: sw_hygro_coreshell_ext(:,:,:,:,:) ! short wave extinction table
160 : real(r8), optional, pointer :: sw_hygro_coreshell_ssa(:,:,:,:,:) ! short wave single-scatter albedo table
161 : real(r8), optional, pointer :: sw_hygro_coreshell_asm(:,:,:,:,:) ! short wave asymmetry table
162 : real(r8), optional, pointer :: lw_hygro_coreshell_ext(:,:,:,:,:) ! long wave absorption table
163 : real(r8), optional, pointer :: corefrac(:) ! core fraction dimension values
164 : real(r8), optional, pointer :: bcdust(:) ! bc/(bc + dust) fraction dimension values
165 : real(r8), optional, pointer :: kap(:) ! hygroscopicity dimension values
166 : real(r8), optional, pointer :: relh(:) ! relative humidity dimension values
167 : integer, optional, intent(out) :: nfrac ! core fraction dimension size
168 : integer, optional, intent(out) :: nbcdust ! bc/(bc + dust) fraction dimension size
169 : integer, optional, intent(out) :: nkap ! hygroscopicity dimension size
170 : integer, optional, intent(out) :: nrelh ! relative humidity dimension size
171 :
172 : end subroutine aero_optics_params
173 :
174 : !------------------------------------------------------------------------
175 : ! returns species type
176 : !------------------------------------------------------------------------
177 : subroutine aero_species_type(self, bin_ndx, species_ndx, spectype)
178 : import :: aerosol_properties
179 : class(aerosol_properties), intent(in) :: self
180 : integer, intent(in) :: bin_ndx ! bin number
181 : integer, intent(in) :: species_ndx ! species number
182 : character(len=*), intent(out) :: spectype ! species type
183 :
184 : end subroutine aero_species_type
185 :
186 : !------------------------------------------------------------------------
187 : ! returns mass and number activation fractions
188 : !------------------------------------------------------------------------
189 : subroutine aero_actfracs(self, bin_ndx, smc, smax, fn, fm )
190 : import :: aerosol_properties, r8
191 : class(aerosol_properties), intent(in) :: self
192 : integer, intent(in) :: bin_ndx ! bin index
193 : real(r8),intent(in) :: smc ! critical supersaturation for particles of bin radius
194 : real(r8),intent(in) :: smax ! maximum supersaturation for multiple competing aerosols
195 : real(r8),intent(out) :: fn ! activation fraction for aerosol number
196 : real(r8),intent(out) :: fm ! activation fraction for aerosol mass
197 :
198 : end subroutine aero_actfracs
199 :
200 : !------------------------------------------------------------------------
201 : ! returns constituents names of aerosol number mixing ratios
202 : !------------------------------------------------------------------------
203 : subroutine aero_num_names(self, bin_ndx, name_a, name_c)
204 : import :: aerosol_properties
205 : class(aerosol_properties), intent(in) :: self
206 : integer, intent(in) :: bin_ndx ! bin number
207 : character(len=*), intent(out) :: name_a ! constituent name of ambient aerosol number dens
208 : character(len=*), intent(out) :: name_c ! constituent name of cloud-borne aerosol number dens
209 : end subroutine aero_num_names
210 :
211 : !------------------------------------------------------------------------
212 : ! returns constituents names of aerosol mass mixing ratios
213 : !------------------------------------------------------------------------
214 : subroutine aero_mmr_names(self, bin_ndx, species_ndx, name_a, name_c)
215 : import :: aerosol_properties
216 : class(aerosol_properties), intent(in) :: self
217 : integer, intent(in) :: bin_ndx ! bin number
218 : integer, intent(in) :: species_ndx ! species number
219 : character(len=*), intent(out) :: name_a ! constituent name of ambient aerosol MMR
220 : character(len=*), intent(out) :: name_c ! constituent name of cloud-borne aerosol MMR
221 : end subroutine aero_mmr_names
222 :
223 : !------------------------------------------------------------------------
224 : ! returns constituent name of ambient aerosol number mixing ratios
225 : !------------------------------------------------------------------------
226 : subroutine aero_amb_num_name(self, bin_ndx, name)
227 : import :: aerosol_properties
228 : class(aerosol_properties), intent(in) :: self
229 : integer, intent(in) :: bin_ndx ! bin number
230 : character(len=*), intent(out) :: name ! constituent name of ambient aerosol number dens
231 :
232 : end subroutine aero_amb_num_name
233 :
234 : !------------------------------------------------------------------------
235 : ! returns constituent name of ambient aerosol mass mixing ratios
236 : !------------------------------------------------------------------------
237 : subroutine aero_amb_mmr_name(self, bin_ndx, species_ndx, name)
238 : import :: aerosol_properties
239 : class(aerosol_properties), intent(in) :: self
240 : integer, intent(in) :: bin_ndx ! bin number
241 : integer, intent(in) :: species_ndx ! species number
242 : character(len=*), intent(out) :: name ! constituent name of ambient aerosol MMR
243 :
244 : end subroutine aero_amb_mmr_name
245 :
246 : !------------------------------------------------------------------------------
247 : ! returns radius^3 (m3) of a given bin number
248 : !------------------------------------------------------------------------------
249 : pure elemental real(r8) function aero_amcube(self, bin_ndx, volconc, numconc)
250 : import :: aerosol_properties, r8
251 :
252 : class(aerosol_properties), intent(in) :: self
253 : integer, intent(in) :: bin_ndx ! bin number
254 : real(r8), intent(in) :: volconc ! volume conc (m3/m3)
255 : real(r8), intent(in) :: numconc ! number conc (1/m3)
256 :
257 : end function aero_amcube
258 :
259 : !------------------------------------------------------------------------------
260 : ! returns TRUE if Ice Nucleation tendencies are applied to given aerosol bin number
261 : !------------------------------------------------------------------------------
262 : function aero_icenuc_updates_num(self, bin_ndx) result(res)
263 : import :: aerosol_properties
264 : class(aerosol_properties), intent(in) :: self
265 : integer, intent(in) :: bin_ndx ! bin number
266 :
267 : logical :: res
268 :
269 : end function aero_icenuc_updates_num
270 :
271 : !------------------------------------------------------------------------------
272 : ! returns TRUE if Ice Nucleation tendencies are applied to a given species within a bin
273 : !------------------------------------------------------------------------------
274 : function aero_icenuc_updates_mmr(self, bin_ndx, species_ndx) result(res)
275 : import :: aerosol_properties
276 : class(aerosol_properties), intent(in) :: self
277 : integer, intent(in) :: bin_ndx ! bin number
278 : integer, intent(in) :: species_ndx ! species number
279 :
280 : logical :: res
281 :
282 : end function aero_icenuc_updates_mmr
283 :
284 : !------------------------------------------------------------------------------
285 : ! apply max / min to number concentration
286 : !------------------------------------------------------------------------------
287 : subroutine aero_apply_num_limits( self, naerosol, vaerosol, istart, istop, m )
288 : import :: aerosol_properties, r8
289 : class(aerosol_properties), intent(in) :: self
290 : real(r8), intent(inout) :: naerosol(:) ! number conc (1/m3)
291 : real(r8), intent(in) :: vaerosol(:) ! volume conc (m3/m3)
292 : integer, intent(in) :: istart ! start column index (1 <= istart <= istop <= pcols)
293 : integer, intent(in) :: istop ! stop column index
294 : integer, intent(in) :: m ! mode or bin index
295 :
296 : end subroutine aero_apply_num_limits
297 :
298 : !------------------------------------------------------------------------------
299 : ! returns TRUE if species `spc_ndx` in aerosol subset `bin_ndx` contributes to
300 : ! the particles' ability to act as heterogeneous freezing nuclei
301 : !------------------------------------------------------------------------------
302 : function aero_hetfrz_species(self, bin_ndx, spc_ndx) result(res)
303 : import :: aerosol_properties
304 : class(aerosol_properties), intent(in) :: self
305 : integer, intent(in) :: bin_ndx ! bin number
306 : integer, intent(in) :: spc_ndx ! species number
307 :
308 : logical :: res
309 :
310 : end function aero_hetfrz_species
311 :
312 : !------------------------------------------------------------------------------
313 : ! returns minimum mass mean radius (meters)
314 : !------------------------------------------------------------------------------
315 : function aero_min_mass_mean_rad(self,bin_ndx,species_ndx) result(minrad)
316 : import :: aerosol_properties, r8
317 : class(aerosol_properties), intent(in) :: self
318 : integer, intent(in) :: bin_ndx ! bin number
319 : integer, intent(in) :: species_ndx ! species number
320 :
321 : real(r8) :: minrad ! meters
322 :
323 : end function aero_min_mass_mean_rad
324 :
325 : !------------------------------------------------------------------------------
326 : ! returns TRUE if soluble
327 : !------------------------------------------------------------------------------
328 : logical function aero_soluble(self,bin_ndx)
329 : import :: aerosol_properties
330 : class(aerosol_properties), intent(in) :: self
331 : integer, intent(in) :: bin_ndx ! bin number
332 :
333 : end function aero_soluble
334 :
335 : !------------------------------------------------------------------------------
336 : ! returns the total number of bins for a given radiation list index
337 : !------------------------------------------------------------------------------
338 : function aero_nbins_rlist(self, list_ndx) result(res)
339 : import :: aerosol_properties
340 : class(aerosol_properties), intent(in) :: self
341 : integer, intent(in) :: list_ndx ! radiation list number
342 :
343 : integer :: res
344 :
345 : end function aero_nbins_rlist
346 :
347 : !------------------------------------------------------------------------------
348 : ! returns number of species in a bin for a given radiation list index
349 : !------------------------------------------------------------------------------
350 : function aero_nspecies_rlist(self, list_ndx, bin_ndx) result(res)
351 : import :: aerosol_properties
352 : class(aerosol_properties), intent(in) :: self
353 : integer, intent(in) :: list_ndx ! radiation list number
354 : integer, intent(in) :: bin_ndx ! bin number
355 :
356 : integer :: res
357 :
358 : end function aero_nspecies_rlist
359 :
360 : !------------------------------------------------------------------------------
361 : ! returns the natural log of geometric standard deviation of the number
362 : ! distribution for radiation list number and aerosol bin
363 : !------------------------------------------------------------------------------
364 : function aero_alogsig_rlist(self, list_ndx, bin_ndx) result(res)
365 : import :: aerosol_properties, r8
366 : class(aerosol_properties), intent(in) :: self
367 : integer, intent(in) :: list_ndx ! radiation list number
368 : integer, intent(in) :: bin_ndx ! bin number
369 :
370 : real(r8) :: res
371 :
372 : end function aero_alogsig_rlist
373 :
374 : !------------------------------------------------------------------------------
375 : ! returns name for a given radiation list number and aerosol bin
376 : !------------------------------------------------------------------------------
377 : function aero_bin_name(self, list_ndx, bin_ndx) result(name)
378 : import :: aerosol_properties, r8
379 : class(aerosol_properties), intent(in) :: self
380 : integer, intent(in) :: list_ndx ! radiation list number
381 : integer, intent(in) :: bin_ndx ! bin number
382 :
383 : character(len=32) name
384 :
385 : end function aero_bin_name
386 :
387 : !------------------------------------------------------------------------------
388 : ! returns scavenging diameter for a given aerosol bin number
389 : !------------------------------------------------------------------------------
390 : function aero_scav_diam(self, bin_ndx) result(diam)
391 : import :: aerosol_properties, r8
392 : class(aerosol_properties), intent(in) :: self
393 : integer, intent(in) :: bin_ndx ! bin number
394 :
395 : real(r8) :: diam
396 :
397 : end function aero_scav_diam
398 :
399 : !------------------------------------------------------------------------------
400 : ! adjust aerosol concentration tendencies to create larger sizes of aerosols
401 : ! during resuspension
402 : !------------------------------------------------------------------------------
403 : subroutine aero_resuspension_resize(self, dcondt)
404 : import :: aerosol_properties, r8
405 :
406 : class(aerosol_properties), intent(in) :: self
407 : real(r8), intent(inout) :: dcondt(:)
408 :
409 : end subroutine aero_resuspension_resize
410 :
411 : !------------------------------------------------------------------------------
412 : ! returns bulk deposition fluxes of the specified species type
413 : ! rebinned to specified diameter limits
414 : !------------------------------------------------------------------------------
415 : subroutine aero_rebin_bulk_fluxes(self, bulk_type, dep_fluxes, diam_edges, bulk_fluxes, &
416 : error_code, error_string)
417 : import :: aerosol_properties, r8
418 : class(aerosol_properties), intent(in) :: self
419 : character(len=*),intent(in) :: bulk_type ! aerosol type to rebin
420 : real(r8), intent(in) :: dep_fluxes(:) ! kg/m2
421 : real(r8), intent(in) :: diam_edges(:) ! meters
422 : real(r8), intent(out) :: bulk_fluxes(:) ! kg/m2
423 : integer, intent(out) :: error_code ! error code (0 if no error)
424 : character(len=*), intent(out) :: error_string ! error string
425 :
426 : end subroutine aero_rebin_bulk_fluxes
427 :
428 : !------------------------------------------------------------------------------
429 : ! Returns TRUE if bin is hydrophilic, otherwise FALSE
430 : !------------------------------------------------------------------------------
431 : logical function aero_hydrophilic(self, bin_ndx)
432 : import :: aerosol_properties
433 : class(aerosol_properties), intent(in) :: self
434 : integer, intent(in) :: bin_ndx ! bin number
435 : end function aero_hydrophilic
436 :
437 : end interface
438 :
439 : contains
440 :
441 : !------------------------------------------------------------------------------
442 : ! object initializer
443 : !------------------------------------------------------------------------------
444 6144 : subroutine aero_props_init(self, nbin, ncnst, nspec, nmasses, alogsig, f1,f2, ierr )
445 : class(aerosol_properties), intent(inout) :: self
446 : integer, intent(in) :: nbin ! number of bins
447 : integer, intent(in) :: ncnst ! total number of constituents
448 : integer, intent(in) :: nspec(nbin) ! number of species in each bin
449 : integer, intent(in) :: nmasses(nbin) ! number of masses in each bin
450 : real(r8),intent(in) :: alogsig(nbin) ! natural log of the standard deviation (sigma) of the aerosol bins
451 : real(r8),intent(in) :: f1(nbin) ! eq 28 Abdul-Razzak et al 1998
452 : real(r8),intent(in) :: f2(nbin) ! eq 29 Abdul-Razzak et al 1998
453 : integer,intent(out) :: ierr
454 :
455 : integer :: imas,ibin,indx
456 : character(len=*),parameter :: prefix = 'aerosol_properties::aero_props_init: '
457 :
458 : real(r8), parameter :: spechygro_so4 = 0.507_r8 ! Sulfate hygroscopicity
459 : real(r8), parameter :: spechygro_soa = 0.14_r8 ! SOA hygroscopicity
460 : real(r8), parameter :: spechygro_pom = 0.1_r8 ! POM hygroscopicity
461 :
462 6144 : ierr = 0
463 :
464 18432 : allocate(self%nspecies_(nbin),stat=ierr)
465 6144 : if( ierr /= 0 ) then
466 : return
467 : end if
468 12288 : allocate(self%nmasses_(nbin),stat=ierr)
469 6144 : if( ierr /= 0 ) then
470 : return
471 : end if
472 18432 : allocate(self%alogsig_(nbin),stat=ierr)
473 6144 : if( ierr /= 0 ) then
474 : return
475 : end if
476 12288 : allocate(self%f1_(nbin),stat=ierr)
477 6144 : if( ierr /= 0 ) then
478 : return
479 : end if
480 12288 : allocate(self%f2_(nbin),stat=ierr)
481 6144 : if( ierr /= 0 ) then
482 : return
483 : end if
484 :
485 49152 : allocate( self%indexer_(nbin,0:maxval(nmasses)),stat=ierr )
486 6144 : if( ierr /= 0 ) then
487 : return
488 : end if
489 :
490 : ! Local indexing compresses the mode and number/mass indices into one index.
491 : ! This indexing is used by the pointer arrays used to reference state and pbuf
492 : ! fields. We add number = 0, total mass = 1 (if available), and mass from each
493 : ! constituency into mm.
494 :
495 221184 : self%indexer_ = -1
496 : indx = 0
497 :
498 30720 : do ibin=1,nbin
499 147456 : do imas = 0,nmasses(ibin)
500 116736 : indx = indx+1
501 141312 : self%indexer_(ibin,imas) = indx
502 : end do
503 : end do
504 :
505 6144 : self%nbins_ = nbin
506 6144 : self%ncnst_tot_ = ncnst
507 30720 : self%nmasses_(:) = nmasses(:)
508 30720 : self%nspecies_(:) = nspec(:)
509 30720 : self%alogsig_(:) = alogsig(:)
510 30720 : self%f1_(:) = f1(:)
511 30720 : self%f2_(:) = f2(:)
512 :
513 6144 : self%soa_equivso4_factor_ = spechygro_soa/spechygro_so4
514 6144 : self%pom_equivso4_factor_ = spechygro_pom/spechygro_so4
515 :
516 : end subroutine aero_props_init
517 :
518 : !------------------------------------------------------------------------------
519 : ! Object clean
520 : !------------------------------------------------------------------------------
521 1536 : subroutine aero_props_final(self)
522 : class(aerosol_properties), intent(inout) :: self
523 :
524 1536 : if (allocated(self%nspecies_)) then
525 1536 : deallocate(self%nspecies_)
526 : end if
527 1536 : if (allocated(self%nmasses_)) then
528 1536 : deallocate(self%nmasses_)
529 : end if
530 1536 : if (allocated(self%indexer_)) then
531 1536 : deallocate(self%indexer_)
532 : endif
533 1536 : if (allocated(self%alogsig_)) then
534 1536 : deallocate(self%alogsig_)
535 : endif
536 1536 : if (allocated(self%f1_)) then
537 1536 : deallocate(self%f1_)
538 : endif
539 1536 : if (allocated(self%f2_)) then
540 1536 : deallocate(self%f2_)
541 : endif
542 :
543 1536 : self%nbins_ = 0
544 1536 : self%ncnst_tot_ = 0
545 :
546 1536 : end subroutine aero_props_final
547 :
548 : !------------------------------------------------------------------------------
549 : ! returns number of species in a bin
550 : !------------------------------------------------------------------------------
551 12767931338 : pure function nspecies_per_bin(self,bin_ndx) result(val)
552 : class(aerosol_properties), intent(in) :: self
553 : integer, intent(in) :: bin_ndx ! bin number
554 : integer :: val
555 :
556 12767931338 : val = self%nspecies_(bin_ndx)
557 12767931338 : end function nspecies_per_bin
558 :
559 : !------------------------------------------------------------------------------
560 : ! returns number of species for all bins
561 : !------------------------------------------------------------------------------
562 3072 : pure function nspecies_all_bins(self) result(arr)
563 : class(aerosol_properties), intent(in) :: self
564 : integer :: arr(self%nbins_)
565 :
566 15360 : arr(:) = self%nspecies_(:)
567 :
568 3072 : end function nspecies_all_bins
569 :
570 : !------------------------------------------------------------------------------
571 : ! returns number of species masses in a given bin number
572 : !------------------------------------------------------------------------------
573 6192344848 : pure function n_masses_per_bin(self,bin_ndx) result(val)
574 : class(aerosol_properties), intent(in) :: self
575 : integer, intent(in) :: bin_ndx ! bin number
576 : integer :: val
577 :
578 6192344848 : val = self%nmasses_(bin_ndx)
579 6192344848 : end function n_masses_per_bin
580 :
581 : !------------------------------------------------------------------------------
582 : ! returns an array of number of species masses for all bins
583 : !------------------------------------------------------------------------------
584 4608 : pure function n_masses_all_bins(self) result(arr)
585 : class(aerosol_properties), intent(in) :: self
586 : integer :: arr(self%nbins_)
587 :
588 23040 : arr(:) = self%nmasses_(:)
589 4608 : end function n_masses_all_bins
590 :
591 : !------------------------------------------------------------------------------
592 : ! returns a single index for given bin and species
593 : !------------------------------------------------------------------------------
594 35090384453 : pure integer function indexer(self,bin_ndx,species_ndx)
595 : class(aerosol_properties), intent(in) :: self
596 : integer, intent(in) :: bin_ndx ! bin number
597 : integer, intent(in) :: species_ndx ! species number
598 :
599 35090384453 : indexer = self%indexer_(bin_ndx,species_ndx)
600 35090384453 : end function indexer
601 :
602 : !------------------------------------------------------------------------------
603 : ! returns the total number of bins
604 : !------------------------------------------------------------------------------
605 6708844464 : pure function nbins_0list(self) result(nbins)
606 : class(aerosol_properties), intent(in) :: self
607 : integer :: nbins
608 :
609 6708844464 : nbins = self%nbins_
610 6708844464 : end function nbins_0list
611 :
612 : !------------------------------------------------------------------------------
613 : ! returns number of constituents (or elements) totaled across all bins
614 : !------------------------------------------------------------------------------
615 4475208 : pure integer function ncnst_tot(self)
616 : class(aerosol_properties), intent(in) :: self
617 :
618 4475208 : ncnst_tot = self%ncnst_tot_
619 4475208 : end function ncnst_tot
620 :
621 : !------------------------------------------------------------------------------
622 : ! returns the natural log of geometric standard deviation of the number distribution for aerosol bin
623 : !------------------------------------------------------------------------------
624 3371132784 : pure real(r8) function alogsig_0list(self, bin_ndx)
625 : class(aerosol_properties), intent(in) :: self
626 : integer, intent(in) :: bin_ndx ! bin number
627 :
628 3371132784 : alogsig_0list = self%alogsig_(bin_ndx)
629 3371132784 : end function alogsig_0list
630 :
631 : !------------------------------------------------------------------------------
632 : ! returns maximum supersaturation
633 : !------------------------------------------------------------------------------
634 365272372 : function maxsat(self, zeta,eta,smc) result(smax)
635 :
636 : !-------------------------------------------------------------------------
637 : ! Calculates maximum supersaturation for multiple competing aerosols.
638 : !
639 : ! Abdul-Razzak and Ghan, A parameterization of aerosol activation.
640 : ! 2. Multiple aerosol types. J. Geophys. Res., 105, 6837-6844., 2000
641 : !-------------------------------------------------------------------------
642 :
643 : class(aerosol_properties), intent(in) :: self
644 : real(r8), intent(in) :: zeta(self%nbins_) ! Abdul-Razzak and Ghan eq 10
645 : real(r8), intent(in) :: eta(self%nbins_) ! Abdul-Razzak and Ghan eq 11
646 : real(r8), intent(in) :: smc(self%nbins_) ! critical supersaturation
647 :
648 : real(r8) :: smax ! maximum supersaturation
649 :
650 : integer :: m
651 : integer :: nbins
652 : real(r8) :: sum, g1, g2, g1sqrt, g2sqrt
653 :
654 : real(r8), parameter :: small_maxsat = 1.e-20_r8 ! for weak forcing
655 : real(r8), parameter :: large_maxsat = 1.e20_r8 ! for small eta
656 :
657 365272372 : smax=0.0_r8
658 365272372 : nbins = self%nbins_
659 :
660 365272372 : check_loop: do m=1,nbins
661 365272372 : if((zeta(m) > 1.e5_r8*eta(m)) .or. (smc(m)*smc(m) > 1.e5_r8*eta(m))) then
662 : ! weak forcing -- essentially none activated
663 0 : smax=small_maxsat
664 : else
665 : ! significant activation of this mode -- calc activation all modes
666 : exit check_loop
667 : endif
668 : ! No significant activation in any mode. Do nothing.
669 365272372 : if (m == nbins) return
670 : enddo check_loop
671 :
672 : sum=0.0_r8
673 :
674 1826361860 : do m=1,nbins
675 1826361860 : if(eta(m) > 1.e-20_r8)then
676 : ! from Abdul-Razzak and Ghan 2000
677 1461089488 : g1=zeta(m)/eta(m)
678 1461089488 : g1sqrt=sqrt(g1)
679 1461089488 : g1=g1sqrt*g1
680 1461089488 : g2=smc(m)/sqrt(eta(m)+3._r8*zeta(m))
681 1461089488 : g2sqrt=sqrt(g2)
682 1461089488 : g2=g2sqrt*g2
683 1461089488 : sum=sum+(self%f1_(m)*g1+self%f2_(m)*g2)/(smc(m)*smc(m))
684 : else
685 : sum=large_maxsat
686 : endif
687 : enddo
688 :
689 365272372 : smax=1._r8/sqrt(sum)
690 :
691 365272372 : end function maxsat
692 :
693 : !------------------------------------------------------------------------------
694 : ! returns the ratio of SOA Hygroscopicity / Sulfate Hygroscopicity
695 : !------------------------------------------------------------------------------
696 8935056 : pure real(r8) function soa_equivso4_factor(self)
697 : class(aerosol_properties), intent(in) :: self
698 :
699 8935056 : soa_equivso4_factor = self%soa_equivso4_factor_
700 :
701 8935056 : end function soa_equivso4_factor
702 :
703 : !------------------------------------------------------------------------------
704 : ! returns the ratio of POM Hygroscopicity / Sulfate Hygroscopicity
705 : !------------------------------------------------------------------------------
706 13402584 : pure real(r8) function pom_equivso4_factor(self)
707 : class(aerosol_properties), intent(in) :: self
708 :
709 13402584 : pom_equivso4_factor = self%pom_equivso4_factor_
710 :
711 13402584 : end function pom_equivso4_factor
712 :
713 0 : end module aerosol_properties_mod
|