PartMC  2.8.0
aero_mode.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2016 Nicole Riemer and Matthew West
2 ! Licensed under the GNU General Public License version 2 or (at your
3 ! option) any later version. See the file COPYING for details.
4 
5 !> \file
6 !> The pmc_aero_mode module.
7 
8 !> The aero_mode_t structure and associated subroutines.
10 
11  use pmc_bin_grid
12  use pmc_util
13  use pmc_constants
14  use pmc_spec_file
15  use pmc_aero_data
16  use pmc_aero_weight
17  use pmc_mpi
18  use pmc_rand
19 #ifdef PMC_USE_MPI
20  use mpi
21 #endif
22 
23  !> Maximum length of a mode name.
24  integer, parameter :: aero_mode_name_len = 300
25  !> Maximum length of a mode type.
26  integer, parameter :: aero_mode_type_len = 20
27 
28  !> Type code for an undefined or invalid mode.
29  integer, parameter :: aero_mode_type_invalid = 0
30  !> Type code for a log-normal mode.
31  integer, parameter :: aero_mode_type_log_normal = 1
32  !> Type code for an exponential mode.
33  integer, parameter :: aero_mode_type_exp = 2
34  !> Type code for a mono-disperse mode.
35  integer, parameter :: aero_mode_type_mono = 3
36  !> Type code for a sampled mode.
37  integer, parameter :: aero_mode_type_sampled = 4
38 
39  !> Type code for an undefined for invalid diameter type.
40  integer, parameter :: aero_mode_diam_type_invalid = 0
41  !> Type code for geometric diameter.
42  integer, parameter :: aero_mode_diam_type_geometric = 1
43  !> Type code for mobility equivalent diameter.
44  integer, parameter :: aero_mode_diam_type_mobility = 2
45 
46  !> An aerosol size distribution mode.
47  !!
48  !! Each mode is assumed to be fully internally mixed so that every
49  !! particle has the same composition. The composition is stored in
50  !! \c vol_frac, while the other parameters define the size
51  !! distribution (with \c type defining the type of size distribution
52  !! function). See \ref input_format_mass_frac for descriptions of
53  !! the parameters relevant to each mode type.
55  !> Mode name, used to track particle sources.
56  character(len=AERO_MODE_NAME_LEN) :: name
57  !> Mode type (given by module constants).
58  integer :: type
59  !> Characteristic radius, with meaning dependent on mode type (m).
60  real(kind=dp) :: char_radius
61  !> Log base 10 of geometric standard deviation of radius, (m).
62  real(kind=dp) :: log10_std_dev_radius
63  !> Sample bin radii [length <tt>(N + 1)</tt>] (m).
64  real(kind=dp), allocatable :: sample_radius(:)
65  !> Sample bin number concentrations [length <tt>N</tt>] (m^{-3}).
66  real(kind=dp), allocatable :: sample_num_conc(:)
67  !> Total number concentration of mode (#/m^3).
68  real(kind=dp) :: num_conc
69  !> Species fractions by volume [length \c aero_data_n_spec(aero_data)] (1).
70  real(kind=dp), allocatable :: vol_frac(:)
71  !> Species fraction standard deviation
72  !> [length \c aero_data_n_spec(aero_data)] (1).
73  real(kind=dp), allocatable :: vol_frac_std(:)
74  !> Source number.
75  integer :: source
76  !> Class number.
77  integer :: weight_class
78  end type aero_mode_t
79 
80 contains
81 
82 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
83 
84  !> Return a string representation of a kernel type.
85  character(len=AERO_MODE_TYPE_LEN) function aero_mode_type_to_string(type)
86 
87  !> Aero mode type.
88  integer, intent(in) :: type
89 
90  if (type == aero_mode_type_invalid) then
91  aero_mode_type_to_string = "invalid"
92  elseif (type == aero_mode_type_log_normal) then
93  aero_mode_type_to_string = "log_normal"
94  elseif (type == aero_mode_type_exp) then
96  elseif (type == aero_mode_type_mono) then
98  elseif (type == aero_mode_type_sampled) then
99  aero_mode_type_to_string = "sampled"
100  else
101  aero_mode_type_to_string = "unknown"
102  end if
103 
104  end function aero_mode_type_to_string
105 
106 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
107 
108  !> Returns the total number concentration of a mode. (#/m^3)
109  real(kind=dp) function aero_mode_total_num_conc(aero_mode)
110 
111  !> Aerosol mode.
112  type(aero_mode_t), intent(in) :: aero_mode
113 
115  if ((aero_mode%type == aero_mode_type_log_normal) &
116  .or. (aero_mode%type == aero_mode_type_exp) &
117  .or. (aero_mode%type == aero_mode_type_mono)) then
118  aero_mode_total_num_conc = aero_mode%num_conc
119  elseif (aero_mode%type == aero_mode_type_sampled) then
120  aero_mode_total_num_conc = sum(aero_mode%sample_num_conc)
121  else
122  call die_msg(719625922, "unknown aero_mode type: " &
123  // trim(integer_to_string(aero_mode%type)))
124  end if
125 
126  end function aero_mode_total_num_conc
127 
128 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
129 
130  !> Compute a log-normal distribution.
131  subroutine num_conc_log_normal(total_num_conc, geom_mean_radius, &
132  log10_sigma_g, bin_grid, num_conc)
133 
134  !> Total number concentration of the mode (m^{-3}).
135  real(kind=dp), intent(in) :: total_num_conc
136  !> Geometric mean radius (m).
137  real(kind=dp), intent(in) :: geom_mean_radius
138  !> log_10(geom. std. dev.) (1).
139  real(kind=dp), intent(in) :: log10_sigma_g
140  !> Bin grid.
141  type(bin_grid_t), intent(in) :: bin_grid
142  !> Number concentration (#(ln(r))d(ln(r))).
143  real(kind=dp), intent(out) :: num_conc(bin_grid_size(bin_grid))
144 
145  integer :: k
146 
147  do k = 1,bin_grid_size(bin_grid)
148  num_conc(k) = total_num_conc / (sqrt(2d0 * const%pi) &
149  * log10_sigma_g) * dexp(-(dlog10(bin_grid%centers(k)) &
150  - dlog10(geom_mean_radius))**2d0 &
151  / (2d0 * log10_sigma_g**2d0)) / dlog(10d0)
152  end do
153 
154  ! The formula above was originally for a distribution in
155  ! log_10(r), while we are using log_e(r) for our bin grid. The
156  ! division by dlog(10) at the end corrects for this.
157 
158  ! Remember that log_e(r) = log_10(r) * log_e(10).
159 
160  end subroutine num_conc_log_normal
161 
162 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
163 
164  !> Compute a log-normal distribution in volume.
165  subroutine vol_conc_log_normal(total_num_conc, &
166  geom_mean_radius, log10_sigma_g, bin_grid, aero_data, vol_conc)
167 
168  !> Total number concentration of the mode (m^{-3}).
169  real(kind=dp), intent(in) :: total_num_conc
170  !> Geometric mean radius (m).
171  real(kind=dp), intent(in) :: geom_mean_radius
172  !> log_10(geom. std. dev.) (1).
173  real(kind=dp), intent(in) :: log10_sigma_g
174  !> Bin grid.
175  type(bin_grid_t), intent(in) :: bin_grid
176  !> Aerosol data.
177  type(aero_data_t), intent(in) :: aero_data
178  !> Volume concentration (V(ln(r))d(ln(r))).
179  real(kind=dp), intent(out) :: vol_conc(bin_grid_size(bin_grid))
180 
181  real(kind=dp) :: num_conc(bin_grid_size(bin_grid))
182 
183  call num_conc_log_normal(total_num_conc, geom_mean_radius, &
184  log10_sigma_g, bin_grid, num_conc)
185 
186  vol_conc = num_conc * aero_data_rad2vol(aero_data, bin_grid%centers)
187 
188  end subroutine vol_conc_log_normal
189 
190 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
191 
192  !> Exponential distribution in volume.
193  !!
194  !! \f[ n(v) = \frac{1}{\rm mean-vol} \exp(- v / {\rm mean-vol}) \f]
195  !! Normalized so that sum(num_conc(k) * log_width) = 1
196  subroutine num_conc_exp(total_num_conc, radius_at_mean_vol, &
197  bin_grid, aero_data, num_conc)
198 
199  !> Total number concentration of the mode (m^{-3}).
200  real(kind=dp), intent(in) :: total_num_conc
201  !> Radius at mean volume (m).
202  real(kind=dp), intent(in) :: radius_at_mean_vol
203  !> Bin grid.
204  type(bin_grid_t), intent(in) :: bin_grid
205  !> Aerosol data.
206  type(aero_data_t), intent(in) :: aero_data
207  !> Number concentration (#(ln(r))d(ln(r))).
208  real(kind=dp), intent(out) :: num_conc(bin_grid_size(bin_grid))
209 
210  integer :: k
211  real(kind=dp) :: mean_vol, num_conc_vol
212 
213  mean_vol = aero_data_rad2vol(aero_data, radius_at_mean_vol)
214  do k = 1,bin_grid_size(bin_grid)
215  num_conc_vol = total_num_conc / mean_vol &
216  * exp(-(aero_data_rad2vol(aero_data, bin_grid%centers(k)) &
217  / mean_vol))
218  call vol_to_lnr(bin_grid%centers(k), num_conc_vol, num_conc(k))
219  end do
220 
221  end subroutine num_conc_exp
222 
223 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
224 
225  !> Exponential distribution in volume.
226  subroutine vol_conc_exp(total_num_conc, radius_at_mean_vol, &
227  bin_grid, aero_data, vol_conc)
228 
229  !> Total number concentration of the mode (m^{-3}).
230  real(kind=dp), intent(in) :: total_num_conc
231  !> Radius at mean volume (m).
232  real(kind=dp), intent(in) :: radius_at_mean_vol
233  !> Bin grid.
234  type(bin_grid_t), intent(in) :: bin_grid
235  !> Aerosol data.
236  type(aero_data_t), intent(in) :: aero_data
237  !> Volume concentration (V(ln(r))d(ln(r))).
238  real(kind=dp), intent(out) :: vol_conc(bin_grid_size(bin_grid))
239 
240  real(kind=dp) :: num_conc(bin_grid_size(bin_grid))
241 
242  call num_conc_exp(total_num_conc, radius_at_mean_vol, &
243  bin_grid, aero_data, num_conc)
244  vol_conc = num_conc * aero_data_rad2vol(aero_data, bin_grid%centers)
245 
246  end subroutine vol_conc_exp
247 
248 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
249 
250  !> Mono-disperse distribution.
251  !> Normalized so that sum(num_conc(k) * log_width) = 1
252  subroutine num_conc_mono(total_num_conc, radius, bin_grid, num_conc)
253 
254  !> Total number concentration of the mode (m^{-3}).
255  real(kind=dp), intent(in) :: total_num_conc
256  !> Radius of each particle (m^3).
257  real(kind=dp), intent(in) :: radius
258  !> Bin grid.
259  type(bin_grid_t), intent(in) :: bin_grid
260  !> Number concentration (#(ln(r))d(ln(r))).
261  real(kind=dp), intent(out) :: num_conc(bin_grid_size(bin_grid))
262 
263  integer :: k
264 
265  num_conc = 0d0
266  k = bin_grid_find(bin_grid, radius)
267  if ((k < 1) .or. (k > bin_grid_size(bin_grid))) then
268  call warn_msg(825666877, "monodisperse radius outside of bin_grid")
269  else
270  num_conc(k) = total_num_conc / bin_grid%widths(k)
271  end if
272 
273  end subroutine num_conc_mono
274 
275 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
276 
277  !> Mono-disperse distribution in volume.
278  subroutine vol_conc_mono(total_num_conc, radius, &
279  bin_grid, aero_data, vol_conc)
280 
281  !> Total number concentration of the mode (m^{-3}).
282  real(kind=dp), intent(in) :: total_num_conc
283  !> Radius of each particle (m^3).
284  real(kind=dp), intent(in) :: radius
285  !> Bin grid.
286  type(bin_grid_t), intent(in) :: bin_grid
287  !> Aerosol data.
288  type(aero_data_t), intent(in) :: aero_data
289  !> Volume concentration (V(ln(r))d(ln(r))).
290  real(kind=dp), intent(out) :: vol_conc(bin_grid_size(bin_grid))
291 
292  integer :: k
293 
294  vol_conc = 0d0
295  k = bin_grid_find(bin_grid, radius)
296  if ((k < 1) .or. (k > bin_grid_size(bin_grid))) then
297  call warn_msg(420930707, "monodisperse radius outside of bin_grid")
298  else
299  vol_conc(k) = total_num_conc / bin_grid%widths(k) &
300  * aero_data_rad2vol(aero_data, radius)
301  end if
302 
303  end subroutine vol_conc_mono
304 
305 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
306 
307  !> Sampled distribution, not normalized.
308  subroutine num_conc_sampled(sample_radius, sample_num_conc, bin_grid, &
309  num_conc)
310 
311  !> Sampled radius bin edges (m).
312  real(kind=dp), intent(in) :: sample_radius(:)
313  !> Sampled number concentrations (m^{-3}).
314  real(kind=dp), intent(in) :: sample_num_conc(:)
315  !> Bin grid.
316  type(bin_grid_t), intent(in) :: bin_grid
317  !> Number concentration (#(ln(r))d(ln(r))).
318  real(kind=dp), intent(out) :: num_conc(bin_grid_size(bin_grid))
319 
320  integer :: i_sample, n_sample, i_lower, i_upper, i_bin
321  real(kind=dp) :: r_lower, r_upper
322  real(kind=dp) :: r_bin_lower, r_bin_upper, r1, r2, ratio
323 
324  n_sample = size(sample_num_conc)
325  call assert(188766208, size(sample_radius) == n_sample + 1)
326  call assert(295384037, n_sample >= 1)
327 
328  num_conc = 0d0
329  do i_sample = 1,n_sample
330  r_lower = sample_radius(i_sample)
331  r_upper = sample_radius(i_sample + 1)
332  i_lower = bin_grid_find(bin_grid, r_lower)
333  i_upper = bin_grid_find(bin_grid, r_upper)
334  if (i_upper < 1) cycle
335  if (i_lower > bin_grid_size(bin_grid)) cycle
336  i_lower = max(1, i_lower)
337  i_upper = min(bin_grid_size(bin_grid), i_upper)
338  do i_bin = i_lower,i_upper
339  r_bin_lower = bin_grid%edges(i_bin)
340  r_bin_upper = bin_grid%edges(i_bin + 1)
341  r1 = max(r_lower, r_bin_lower)
342  r2 = min(r_upper, r_bin_upper)
343  ratio = (log(r2) - log(r1)) / (log(r_upper) - log(r_lower))
344  num_conc(i_bin) = num_conc(i_bin) + ratio &
345  * sample_num_conc(i_sample) / bin_grid%widths(i_bin)
346  end do
347  end do
348 
349  end subroutine num_conc_sampled
350 
351 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
352 
353  !> Sampled distribution in volume.
354  subroutine vol_conc_sampled(sample_radius, sample_num_conc, &
355  bin_grid, aero_data, vol_conc)
356 
357  !> Sampled radius bin edges (m).
358  real(kind=dp), intent(in) :: sample_radius(:)
359  !> Sampled number concentrations (m^{-3}).
360  real(kind=dp), intent(in) :: sample_num_conc(:)
361  !> Bin grid.
362  type(bin_grid_t), intent(in) :: bin_grid
363  !> Aerosol data.
364  type(aero_data_t), intent(in) :: aero_data
365  !> Volume concentration (V(ln(r))d(ln(r))).
366  real(kind=dp), intent(out) :: vol_conc(bin_grid_size(bin_grid))
367 
368  real(kind=dp) :: num_conc(bin_grid_size(bin_grid))
369 
370  call num_conc_sampled(sample_radius, sample_num_conc, bin_grid, num_conc)
371  vol_conc = num_conc * aero_data_rad2vol(aero_data, bin_grid%centers)
372 
373  end subroutine vol_conc_sampled
374 
375 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
376 
377  !> Return the binned number concentration for an aero_mode.
378  subroutine aero_mode_num_conc(aero_mode, bin_grid, aero_data, &
379  num_conc)
380 
381  !> Aero mode for which to compute number concentration.
382  type(aero_mode_t), intent(in) :: aero_mode
383  !> Bin grid.
384  type(bin_grid_t), intent(in) :: bin_grid
385  !> Aerosol data.
386  type(aero_data_t), intent(in) :: aero_data
387  !> Number concentration (#(ln(r))d(ln(r))).
388  real(kind=dp), intent(out) :: num_conc(bin_grid_size(bin_grid))
389 
390  if (aero_mode%type == aero_mode_type_log_normal) then
391  call num_conc_log_normal(aero_mode%num_conc, aero_mode%char_radius, &
392  aero_mode%log10_std_dev_radius, bin_grid, num_conc)
393  elseif (aero_mode%type == aero_mode_type_exp) then
394  call num_conc_exp(aero_mode%num_conc, &
395  aero_mode%char_radius, bin_grid, aero_data, num_conc)
396  elseif (aero_mode%type == aero_mode_type_mono) then
397  call num_conc_mono(aero_mode%num_conc, aero_mode%char_radius, &
398  bin_grid, num_conc)
399  elseif (aero_mode%type == aero_mode_type_sampled) then
400  call num_conc_sampled(aero_mode%sample_radius, &
401  aero_mode%sample_num_conc, bin_grid, num_conc)
402  else
403  call die_msg(223903246, "unknown aero_mode type: " &
404  // trim(integer_to_string(aero_mode%type)))
405  end if
406 
407  end subroutine aero_mode_num_conc
408 
409 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
410 
411  !> Return the binned per-species volume concentration for an
412  !> aero_mode.
413  subroutine aero_mode_vol_conc(aero_mode, bin_grid, aero_data, &
414  vol_conc)
415 
416  !> Aero mode for which to compute volume concentration.
417  type(aero_mode_t), intent(in) :: aero_mode
418  !> Bin grid.
419  type(bin_grid_t), intent(in) :: bin_grid
420  !> Aerosol data.
421  type(aero_data_t), intent(in) :: aero_data
422  !> Volume concentration (V(ln(r))d(ln(r))).
423  real(kind=dp), intent(out) :: vol_conc(bin_grid_size(bin_grid), &
424  aero_data_n_spec(aero_data))
425 
426  integer :: i_spec
427  real(kind=dp) :: vol_conc_total(bin_grid_size(bin_grid))
428 
429  if (aero_mode%type == aero_mode_type_log_normal) then
430  call vol_conc_log_normal(aero_mode%num_conc, &
431  aero_mode%char_radius, aero_mode%log10_std_dev_radius, &
432  bin_grid, aero_data, vol_conc_total)
433  elseif (aero_mode%type == aero_mode_type_exp) then
434  call vol_conc_exp(aero_mode%num_conc, &
435  aero_mode%char_radius, bin_grid, aero_data, vol_conc_total)
436  elseif (aero_mode%type == aero_mode_type_mono) then
437  call vol_conc_mono(aero_mode%num_conc, &
438  aero_mode%char_radius, bin_grid, aero_data, vol_conc_total)
439  elseif (aero_mode%type == aero_mode_type_sampled) then
440  call vol_conc_sampled(aero_mode%sample_radius, &
441  aero_mode%sample_num_conc, bin_grid, aero_data, vol_conc_total)
442  else
443  call die_msg(314169653, "Unknown aero_mode type: " &
444  // trim(integer_to_string(aero_mode%type)))
445  end if
446  call assert_msg(756593082, sum(aero_mode%vol_frac_std) == 0d0, &
447  "cannot convert species fractions with non-zero standard deviation " &
448  // "to binned distributions")
449  do i_spec = 1,aero_data_n_spec(aero_data)
450  vol_conc(:,i_spec) = vol_conc_total * aero_mode%vol_frac(i_spec)
451  end do
452 
453  end subroutine aero_mode_vol_conc
454 
455 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
456 
457  !> Compute weighted sampled number concentrations.
458  subroutine aero_mode_weighted_sampled_num_conc(aero_mode, aero_weight, &
459  weighted_num_conc)
460 
461  !> Aerosol mode.
462  type(aero_mode_t), intent(in) :: aero_mode
463  !> Aerosol weight.
464  type(aero_weight_t), intent(in) :: aero_weight
465  !> Weighted number concentration
466  real(kind=dp), intent(out) :: weighted_num_conc(:)
467 
468  integer :: i_sample
469  real(kind=dp) :: x0, x1
470 
471  call assert(256667423, aero_mode%type == aero_mode_type_sampled)
472  call assert(878731017, &
473  size(weighted_num_conc) == size(aero_mode%sample_num_conc))
474 
475  if (aero_weight%type == aero_weight_type_none) then
476  weighted_num_conc = aero_mode%sample_num_conc
477  elseif ((aero_weight%type == aero_weight_type_power) &
478  .or. (aero_weight%type == aero_weight_type_mfa)) then
479  do i_sample = 1,size(aero_mode%sample_num_conc)
480  x0 = log(aero_mode%sample_radius(i_sample))
481  x1 = log(aero_mode%sample_radius(i_sample + 1))
482  weighted_num_conc(i_sample) = aero_mode%sample_num_conc(i_sample) &
483  / aero_weight%exponent * (exp(- aero_weight%exponent * x0) &
484  - exp(- aero_weight%exponent * x1)) / (x1 - x0)
485  end do
486  else
487  call die_msg(576124393, "unknown aero_weight type: " &
488  // trim(integer_to_string(aero_weight%type)))
489  end if
490 
492 
493 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
494 
495  !> Return the total number of computational particles for an \c aero_mode.
496  real(kind=dp) function aero_mode_number(aero_mode, aero_weight)
497 
498  !> Aero_mode to sample radius from.
499  type(aero_mode_t), intent(in) :: aero_mode
500  !> Aero weight.
501  type(aero_weight_t), intent(in) :: aero_weight
502 
503  real(kind=dp) :: x_mean_prime
504  real(kind=dp), allocatable :: weighted_num_conc(:)
505 
506  aero_mode_number = 0d0
507  if (aero_mode%type == aero_mode_type_log_normal) then
508  if (aero_weight%type == aero_weight_type_none) then
509  aero_mode_number = aero_mode%num_conc / aero_weight%magnitude
510  elseif ((aero_weight%type == aero_weight_type_power) &
511  .or. (aero_weight%type == aero_weight_type_mfa)) then
512  x_mean_prime = log10(aero_mode%char_radius) &
513  - aero_weight%exponent * aero_mode%log10_std_dev_radius**2 &
514  * log(10d0)
515  aero_mode_number = aero_mode%num_conc / aero_weight%magnitude &
516  * exp((x_mean_prime**2 - log10(aero_mode%char_radius)**2) &
517  / (2d0 * aero_mode%log10_std_dev_radius**2))
518  else
519  call die_msg(466668240, "unknown aero_weight type: " &
520  // trim(integer_to_string(aero_weight%type)))
521  end if
522  elseif (aero_mode%type == aero_mode_type_exp) then
523  if (aero_weight%type == aero_weight_type_none) then
524  aero_mode_number = aero_mode%num_conc / aero_weight%magnitude
525  else
526  call die_msg(822252601, &
527  "cannot use exponential modes with weighting")
528  end if
529  elseif (aero_mode%type == aero_mode_type_mono) then
530  aero_mode_number = aero_mode%num_conc &
531  / aero_weight_num_conc_at_radius(aero_weight, &
532  aero_mode%char_radius)
533  elseif (aero_mode%type == aero_mode_type_sampled) then
534  allocate(weighted_num_conc(size(aero_mode%sample_num_conc)))
535  call aero_mode_weighted_sampled_num_conc(aero_mode, aero_weight, &
536  weighted_num_conc)
537  aero_mode_number = sum(weighted_num_conc) / aero_weight%magnitude
538  deallocate(weighted_num_conc)
539  else
540  call die_msg(901140225, "unknown aero_mode type: " &
541  // trim(integer_to_string(aero_mode%type)))
542  end if
543 
544  end function aero_mode_number
545 
546 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
547 
548  !> Return a radius randomly sampled from the mode distribution.
549  subroutine aero_mode_sample_radius(aero_mode, aero_data, aero_weight, &
550  radius, prob_threshold)
551 
552  !> Aero_mode to sample radius from.
553  type(aero_mode_t), intent(in) :: aero_mode
554  !> Aerosol data.
555  type(aero_data_t), intent(in) :: aero_data
556  !> Aero weight.
557  type(aero_weight_t), intent(in) :: aero_weight
558  !> Sampled radius (m).
559  real(kind=dp), intent(out) :: radius
560  !> Threshold for out of range sampling of normal distribution.
561  real(kind=dp), intent(in) :: prob_threshold
562 
563  real(kind=dp) :: x_mean_prime, x0, x1, x, r, inv_nc0, inv_nc1, inv_nc
564  integer :: i_sample
565  real(kind=dp), allocatable :: weighted_num_conc(:)
566 
567  if (aero_mode%type == aero_mode_type_log_normal) then
568  if (aero_weight%type == aero_weight_type_none) then
569  x_mean_prime = log10(aero_mode%char_radius)
570  elseif ((aero_weight%type == aero_weight_type_power) &
571  .or. (aero_weight%type == aero_weight_type_mfa)) then
572  x_mean_prime = log10(aero_mode%char_radius) &
573  - aero_weight%exponent * aero_mode%log10_std_dev_radius**2 &
574  * log(10d0)
575  else
576  call die_msg(517376844, "unknown aero_weight type: " &
577  // trim(integer_to_string(aero_weight%type)))
578  end if
579  radius = 10d0**rand_normal(x_mean_prime, &
580  aero_mode%log10_std_dev_radius, prob_threshold)
581  elseif (aero_mode%type == aero_mode_type_sampled) then
582  allocate(weighted_num_conc(size(aero_mode%sample_num_conc)))
583  call aero_mode_weighted_sampled_num_conc(aero_mode, aero_weight, &
584  weighted_num_conc)
585  i_sample = sample_cts_pdf(weighted_num_conc)
586  deallocate(weighted_num_conc)
587  if ((aero_weight%type == aero_weight_type_none) &
588  .or. (((aero_weight%type == aero_weight_type_power) &
589  .or. (aero_weight%type == aero_weight_type_mfa)) &
590  .and. (aero_weight%exponent == 0d0))) then
591  x0 = log(aero_mode%sample_radius(i_sample))
592  x1 = log(aero_mode%sample_radius(i_sample + 1))
593  r = pmc_random()
594  x = (1d0 - r) * x0 + r * x1
595  radius = exp(x)
596  elseif ((aero_weight%type == aero_weight_type_power) &
597  .or. (aero_weight%type == aero_weight_type_mfa)) then
598  inv_nc0 = 1d0 / aero_weight_num_conc_at_radius(aero_weight, &
599  aero_mode%sample_radius(i_sample))
600  inv_nc1 = 1d0 / aero_weight_num_conc_at_radius(aero_weight, &
601  aero_mode%sample_radius(i_sample + 1))
602  r = pmc_random()
603  inv_nc = (1d0 - r) * inv_nc0 + r * inv_nc1
604  radius = aero_weight_radius_at_num_conc(aero_weight, 1d0 / inv_nc)
605  else
606  call die_msg(769131141, "unknown aero_weight type: " &
607  // trim(integer_to_string(aero_weight%type)))
608  end if
609  elseif (aero_mode%type == aero_mode_type_exp) then
610  if (aero_weight%type == aero_weight_type_none) then
611  radius = aero_data_vol2rad(aero_data, -aero_data_rad2vol(aero_data, &
612  aero_mode%char_radius) * log(pmc_random()))
613  elseif ((aero_weight%type == aero_weight_type_power) &
614  .or. (aero_weight%type == aero_weight_type_mfa)) then
615  call die_msg(678481276, &
616  "cannot use exponential modes with weighting")
617  else
618  call die_msg(301787712, "unknown aero_weight type: " &
619  // trim(integer_to_string(aero_weight%type)))
620  end if
621  elseif (aero_mode%type == aero_mode_type_mono) then
622  radius = aero_mode%char_radius
623  else
624  call die_msg(749122931, "Unknown aero_mode type: " &
625  // trim(integer_to_string(aero_mode%type)))
626  end if
627 
628  end subroutine aero_mode_sample_radius
629 
630 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
631 
632  !> Return an array of volumes randomly sampled from the volume fractions.
633  subroutine aero_mode_sample_vols(aero_mode, total_vol, vols)
634 
635  !> Aero_mode to sample from.
636  type(aero_mode_t), intent(in) :: aero_mode
637  !> Total volume (m^3).
638  real(kind=dp), intent(in) :: total_vol
639  !> Sampled volumes (m^3).
640  real(kind=dp), intent(out) :: vols(size(aero_mode%vol_frac))
641 
642  integer, parameter :: aero_mode_max_sample_loops = 1000000
643 
644  integer :: i_sample
645  real(kind=dp) :: offset
646 
647  ! sampling of volume fractions, normalized to sum to 1 by
648  ! projecting out the mean direction, and accepted if all
649  ! non-negative, otherwise rejected and repeated (accept-reject)
650  do i_sample = 1,aero_mode_max_sample_loops
651  call rand_normal_array_1d(aero_mode%vol_frac, aero_mode%vol_frac_std, &
652  vols)
653  ! add the correct amount of (offset * vol_frac) to vols to ensure
654  ! that sum(vols) is 1
655  offset = 1d0 - sum(vols)
656  vols = vols + offset * aero_mode%vol_frac
657  if (minval(vols) >= 0d0) exit
658  end do
659  if (i_sample == aero_mode_max_sample_loops) then
660  call die_msg(549015143, "Unable to sample non-negative volumes for " &
661  // "mode: " // trim(aero_mode%name))
662  end if
663  vols = vols / sum(vols) * total_vol
664 
665  end subroutine aero_mode_sample_vols
666 
667 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
668 
669  !> Get the weight class for an aero_mode.
670  integer function aero_mode_get_weight_class(aero_mode)
671 
672  !> Aero_mode to get weight class for.
673  type(aero_mode_t) :: aero_mode
674 
675  aero_mode_get_weight_class = aero_mode%weight_class
676 
677  end function aero_mode_get_weight_class
678 
679 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
680 
681  !> Read volume fractions from a data file.
682  subroutine spec_file_read_vol_frac(file, aero_data, vol_frac, vol_frac_std)
683 
684  !> Spec file to read mass fractions from.
685  type(spec_file_t), intent(inout) :: file
686  !> Aero_data data.
687  type(aero_data_t), intent(in) :: aero_data
688  !> Aerosol species volume fractions.
689  real(kind=dp), allocatable, intent(inout) :: vol_frac(:)
690  !> Aerosol species volume fraction standard deviations.
691  real(kind=dp), allocatable, intent(inout) :: vol_frac_std(:)
692 
693  integer :: n_species, species, i
694  character(len=SPEC_LINE_MAX_VAR_LEN), allocatable :: species_name(:)
695  real(kind=dp), allocatable :: species_data(:,:)
696  real(kind=dp) :: tot_vol_frac
697 
698  !> \page input_format_mass_frac Input File Format: Aerosol Mass Fractions
699  !!
700  !! An aerosol mass fractions file must consist of one line per
701  !! aerosol species, with each line having the species name
702  !! followed by the species mass fraction in each aerosol
703  !! particle. The valid species names are those specfied by the
704  !! \ref input_format_aero_data file, but not all species have to
705  !! be listed. Any missing species will have proportions of
706  !! zero. If the proportions do not sum to 1 then they will be
707  !! normalized before use. For example, a mass fractions file file
708  !! could contain:
709  !! <pre>
710  !! # species proportion
711  !! OC 0.3
712  !! BC 0.7
713  !! </pre>
714  !! indicating that the particles are 30% organic carbon and 70%
715  !! black carbon.
716  !!
717  !! Optionally, the standard deviation can also be provided for
718  !! each species as a second number on each line. For example,
719  !! <pre>
720  !! # species proportion std_dev
721  !! OC 0.3 0.1
722  !! BC 0.7 0.2
723  !! </pre>
724  !! indicates that the particles are on average 30% OC and 70% BC,
725  !! but may vary to have particles with 20% OC and 80% BC, or 40%
726  !! OC and 60% BC, for example. The standard deviations will be
727  !! normalized by the sum of the proportions.
728  !!
729  !! Either all species in a given file must have standard
730  !! deviations or none of them can.
731  !!
732  !! See also:
733  !! - \ref spec_file_format --- the input file text format
734  !! - \ref input_format_aero_dist --- the format for a complete
735  !! aerosol distribution with several modes
736  !! - \ref input_format_aero_mode --- the format for each mode
737  !! of an aerosol distribution
738 
739  ! read the aerosol data from the specified file
740  call spec_file_read_real_named_array(file, 0, species_name, &
741  species_data)
742 
743  ! check the data size
744  n_species = size(species_data, 1)
745  if (n_species < 1) then
746  call die_msg(628123166, 'file ' // trim(file%name) &
747  // ' must contain at least one line of data')
748  end if
749  if ((size(species_data, 2) /= 1) .and. (size(species_data, 2) /= 2)) then
750  call die_msg(427666881, 'each line in file ' // trim(file%name) &
751  // ' must contain exactly one or two data values')
752  end if
753 
754  ! copy over the data
755  if (allocated(vol_frac)) deallocate(vol_frac)
756  if (allocated(vol_frac_std)) deallocate(vol_frac_std)
757  allocate(vol_frac(aero_data_n_spec(aero_data)))
758  allocate(vol_frac_std(aero_data_n_spec(aero_data)))
759  vol_frac = 0d0
760  vol_frac_std = 0d0
761  do i = 1,n_species
762  species = aero_data_spec_by_name(aero_data, species_name(i))
763  if (species == 0) then
764  call die_msg(775942501, 'unknown species ' // trim(species_name(i)) &
765  // ' in file ' // trim(file%name))
766  end if
767  vol_frac(species) = species_data(i, 1)
768  if (size(species_data, 2) == 2) then
769  vol_frac_std(species) = species_data(i, 2)
770  end if
771  end do
772 
773  ! convert mass fractions to volume fractions
774  vol_frac = vol_frac / aero_data%density
775  vol_frac_std = vol_frac_std / aero_data%density
776 
777  ! normalize
778  tot_vol_frac = sum(vol_frac)
779  if ((minval(vol_frac) < 0d0) .or. (tot_vol_frac <= 0d0)) then
780  call die_msg(356648030, 'fractions in ' // trim(file%name) &
781  // ' are not positive')
782  end if
783  if (minval(vol_frac_std) < 0d0) then
784  call die_msg(676576501, 'standard deviations in ' // trim(file%name) &
785  // ' are not positive')
786  end if
787  vol_frac = vol_frac / tot_vol_frac
788  vol_frac_std = vol_frac_std / tot_vol_frac
789 
790  end subroutine spec_file_read_vol_frac
791 
792 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
793 
794  !> Read a size distribution from a data file.
795  subroutine spec_file_read_size_dist(file, sample_radius, sample_num_conc)
796 
797  !> Spec file to read size distribution from.
798  type(spec_file_t), intent(inout) :: file
799  !> Sample radius values (m).
800  real(kind=dp), allocatable, intent(inout) :: sample_radius(:)
801  !> Sample number concentrations (m^{-3}).
802  real(kind=dp), allocatable, intent(inout) :: sample_num_conc(:)
803 
804  character(len=SPEC_LINE_MAX_VAR_LEN), allocatable :: names(:)
805  real(kind=dp), allocatable :: data(:,:)
806  integer :: n_sample, i_sample
807 
808  !> \page input_format_size_dist Input File Format: Size Distribution
809  !!
810  !! A size distribution file must consist of two lines:
811  !! - the first line must begin with \c diam and be followed by
812  !! \f$N + 1\f$ space-separated real scalars, giving the diameters
813  !! \f$D_1,\ldots,D_{N+1}\f$ of bin edges (m) --- these must be
814  !! in increasing order, so \f$D_i < D_{i+1}\f$
815  !! - the second line must begin with \c num_conc and be followed
816  !! by \f$N\f$ space-separated real scalars, giving the number
817  !! concenrations \f$C_1,\ldots,C_N\f$ in each bin (#/m^3) ---
818  !! \f$C_i\f$ is the total number concentrations of particles
819  !! with diameters in \f$[D_i, D_{i+1}]\f$
820  !!
821  !! The resulting size distribution is taken to be piecewise
822  !! constant in log-diameter coordinates.
823  !!
824  !! Example: a size distribution could be:
825  !! <pre>
826  !! diam 1e-7 1e-6 1e-5 # bin edge diameters (m)
827  !! num_conc 1e9 1e8 # bin number concentrations (m^{-3})
828  !! </pre>
829  !! This distribution has 1e9 particles per cubic meter with
830  !! diameters between 0.1 micron and 1 micron, and 1e8 particles
831  !! per cubic meter with diameters between 1 micron and 10 micron.
832  !!
833  !! See also:
834  !! - \ref spec_file_format --- the input file text format
835  !! - \ref input_format_aero_dist --- the format for a complete
836  !! aerosol distribution with several modes
837  !! - \ref input_format_aero_mode --- the format for each mode
838  !! of an aerosol distribution
839 
840  ! read the data from the file
841  call spec_file_read_real_named_array(file, 1, names, data)
842  call spec_file_assert_msg(311818741, file, size(names) == 1, &
843  'must contain a line starting with "diam"')
844  call spec_file_check_name(file, 'diam', names(1))
845  n_sample = size(data,2) - 1
846  call spec_file_assert_msg(669011124, file, n_sample >= 1, &
847  'must have at least two diam values')
848 
849  if (allocated(sample_radius)) deallocate(sample_radius)
850  allocate(sample_radius(n_sample + 1))
851  sample_radius = diam2rad(data(1,:))
852  do i_sample = 1,n_sample
853  call spec_file_assert_msg(528089871, file, &
854  sample_radius(i_sample) < sample_radius(i_sample + 1), &
855  'diam values must be strictly increasing')
856  end do
857 
858  call spec_file_read_real_named_array(file, 1, names, data)
859  call spec_file_assert_msg(801676496, file, size(names) == 1, &
860  'must contain a line starting with "num_conc"')
861  call spec_file_check_name(file, 'num_conc', names(1))
862 
863  call spec_file_assert_msg(721029144, file, size(data, 2) == n_sample, &
864  'must have one fewer num_conc than diam values')
865 
866  if (allocated(sample_num_conc)) deallocate(sample_num_conc)
867  allocate(sample_num_conc(n_sample))
868  sample_num_conc = data(1,:)
869  do i_sample = 1,n_sample
870  call spec_file_assert_msg(356490397, file, &
871  sample_num_conc(i_sample) >= 0d0, &
872  'num_conc values must be non-negative')
873  end do
874 
875  end subroutine spec_file_read_size_dist
876 
877 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
878 
879  !> Read one mode of an aerosol distribution (number concentration,
880  !> volume fractions, and mode shape).
881  subroutine spec_file_read_aero_mode(file, aero_data, &
882  read_aero_weight_classes, aero_mode, eof)
883 
884  !> Spec file.
885  type(spec_file_t), intent(inout) :: file
886  !> Aero_data data.
887  type(aero_data_t), intent(inout) :: aero_data
888  !> Aerosol mode.
889  type(aero_mode_t), intent(inout) :: aero_mode
890  !> If eof instead of reading data.
891  logical :: eof
892  !> Whether the weight classes for each source are specified in inputs.
893  logical, intent(in) :: read_aero_weight_classes
894 
895  character(len=SPEC_LINE_MAX_VAR_LEN) :: tmp_str, mode_type, diam_type_str
896  character(len=SPEC_LINE_MAX_VAR_LEN) :: mass_frac_filename
897  character(len=SPEC_LINE_MAX_VAR_LEN) :: size_dist_filename
898  character(len=SPEC_LINE_MAX_VAR_LEN) :: weight_class_name
899  type(spec_line_t) :: line
900  type(spec_file_t) :: mass_frac_file, size_dist_file
901  real(kind=dp) :: diam, temp, pressure
902  integer :: diam_type, i_radius
903  real(kind=dp) :: log10_r0_mob, log10_r1_mob, log10_r2_mob, &
904  r0_mob, r2_mob, r0_geom, r2_geom, log10_r0_geom, log10_r2_geom
905 
906  ! note that doxygen's automatic list creation breaks on the list
907  ! below for some reason, so we just use html format
908 
909  !> \page input_format_aero_mode Input File Format: Aerosol Distribution Mode
910  !!
911  !! An aerosol distribution mode has the parameters:
912  !! <ul>
913  !! <li> \b mode_name (string): the name of the mode (for
914  !! informational purposes only)
915  !! <li> \b mass_frac (string): name of file from which to read the
916  !! species mass fractions --- the file format should
917  !! be \subpage input_format_mass_frac
918  !! <li> \b diam_type (string): the type of diameter for the mode
919  !! --- must be one of: \c geometric for geometric diameter;
920  !! or \c mobility for mobility equivalent diameter
921  !! <li> if \c diam_type is \c mobility then the following
922  !! parameters are:
923  !! <ul>
924  !! <li> \b temp (real, unit K): the temperate at which the
925  !! mobility diameters were measured
926  !! <li> \b pressure (real, unit Pa): the pressure at which the
927  !! mobility diameters were measured
928  !! </ul>
929  !! <li> \b mode_type (string): the functional form of the mode ---
930  !! must be one of: \c log_normal for a log-normal
931  !! distribution; \c exp for an exponential distribution; \c
932  !! mono for a mono-disperse distribution; or \c sampled for a
933  !! sampled distribution
934  !! <li> if \c mode_type is \c log_normal then the mode distribution
935  !! shape is
936  !! \f[ n(\log D) {\rm d}\log D
937  !! = \frac{N_{\rm total}}{\sqrt{2\pi} \log \sigma_{\rm g}}
938  !! \exp\left(\frac{(\log D - \log D_{\rm gn})^2}{2 \log ^2
939  !! \sigma_{\rm g}}\right)
940  !! {\rm d}\log D \f]
941  !! and the following parameters are:
942  !! <ul>
943  !! <li> \b num_conc (real, unit 1/m^3): the total number
944  !! concentration \f$N_{\rm total}\f$ of the mode
945  !! <li> \b geom_mean_diam (real, unit m): the geometric mean
946  !! diameter \f$D_{\rm gn}\f$
947  !! <li> \b log10_geom_std_dev (real, dimensionless):
948  !! \f$\log_{10}\f$ of the geometric standard deviation
949  !! \f$\sigma_{\rm g}\f$ of the diameter
950  !! </ul>
951  !! <li> if \c mode_type is \c exp then the mode distribution shape is
952  !! \f[ n(v) {\rm d}v = \frac{N_{\rm total}}{v_{\rm \mu}}
953  !! \exp\left(- \frac{v}{v_{\rm \mu}}\right)
954  !! {\rm d}v \f]
955  !! and the following parameters are:
956  !! <ul>
957  !! <li> \b num_conc (real, unit 1/m^3): the total number
958  !! concentration \f$N_{\rm total}\f$ of the mode
959  !! <li> \b diam_at_mean_vol (real, unit m): the diameter
960  !! \f$D_{\rm \mu}\f$ such that \f$v_{\rm \mu}
961  !! = \frac{\pi}{6} D^3_{\rm \mu}\f$
962  !! </ul>
963  !! <li> if \c mode_type is \c mono then the mode distribution shape
964  !! is a delta distribution at diameter \f$D_0\f$ and the
965  !! following parameters are:
966  !! <ul>
967  !! <li> \b num_conc (real, unit 1/m^3): the total number
968  !! concentration \f$N_{\rm total}\f$ of the mode
969  !! <li> \b radius (real, unit m): the radius \f$R_0\f$ of the
970  !! particles, so that \f$D_0 = 2 R_0\f$
971  !! </ul>
972  !! <li> if \c mode_type is \c sampled then the mode distribution
973  !! shape is piecewise constant (in log-diameter coordinates)
974  !! and the following parameters are:
975  !! <ul>
976  !! <li> \b size_dist (string): name of file from which to
977  !! read the size distribution --- the file format should
978  !! be \subpage input_format_size_dist
979  !! </ul>
980  !! </ul>
981  !!
982  !! Example:
983  !! <pre>
984  !! mode_name diesel # mode name (descriptive only)
985  !! mass_frac comp_diesel.dat # mass fractions in each aerosol particle
986  !! mode_type log_normal # type of distribution
987  !! num_conc 1.6e8 # particle number density (#/m^3)
988  !! geom_mean_diam 2.5e-8 # geometric mean diameter (m)
989  !! log10_geom_std_dev 0.24 # log_10 of geometric standard deviation
990  !! </pre>
991  !!
992  !! See also:
993  !! - \ref spec_file_format --- the input file text format
994  !! - \ref input_format_aero_dist --- the format for a complete
995  !! aerosol distribution with several modes
996  !! - \ref input_format_mass_frac --- the format for the mass
997  !! fractions file
998 
999  call spec_file_read_line(file, line, eof)
1000  if (.not. eof) then
1001  call spec_file_check_line_name(file, line, "mode_name")
1002  call spec_file_check_line_length(file, line, 1)
1003  tmp_str = line%data(1) ! hack to avoid gfortran warning
1004  aero_mode%name = tmp_str(1:aero_mode_name_len)
1005  aero_mode%source = aero_data_source_by_name(aero_data, aero_mode%name)
1006  if (read_aero_weight_classes) then
1007  call spec_file_read_string(file, 'weight_class_name', &
1008  weight_class_name)
1009  else
1010  weight_class_name = aero_mode%name
1011  end if
1012  aero_mode%weight_class = aero_data_weight_class_by_name(aero_data, &
1013  weight_class_name)
1014  call spec_file_read_string(file, 'mass_frac', mass_frac_filename)
1015  call spec_file_open(mass_frac_filename, mass_frac_file)
1016  call spec_file_read_vol_frac(mass_frac_file, aero_data, &
1017  aero_mode%vol_frac, aero_mode%vol_frac_std)
1018  call spec_file_close(mass_frac_file)
1019 
1020  call spec_file_read_string(file, 'diam_type', diam_type_str)
1021  if (trim(diam_type_str) == 'geometric') then
1022  diam_type = aero_mode_diam_type_geometric
1023  elseif (trim(diam_type_str) == 'mobility') then
1024  diam_type = aero_mode_diam_type_mobility
1025  call spec_file_read_real(file, 'temp', temp)
1026  call spec_file_read_real(file, 'pressure', pressure)
1027  else
1028  call spec_file_die_msg(804343794, file, &
1029  "Unknown diam_type: " // trim(diam_type_str))
1030  end if
1031 
1032  call spec_file_read_string(file, 'mode_type', mode_type)
1033  aero_mode%sample_radius = [ real(kind=dp) :: ]
1034  aero_mode%sample_num_conc = [ real(kind=dp) :: ]
1035  if (trim(mode_type) == 'log_normal') then
1036  aero_mode%type = aero_mode_type_log_normal
1037  call spec_file_read_real(file, 'num_conc', aero_mode%num_conc)
1038  call spec_file_read_real(file, 'geom_mean_diam', diam)
1039  call spec_file_read_real(file, 'log10_geom_std_dev', &
1040  aero_mode%log10_std_dev_radius)
1041  if (diam_type == aero_mode_diam_type_geometric) then
1042  aero_mode%char_radius = diam2rad(diam)
1043  elseif (diam_type == aero_mode_diam_type_mobility) then
1044  aero_mode%char_radius &
1046  diam2rad(diam), temp, pressure)
1047 
1048  ! Convert log10_std_dev_radius from mobility to geometric radius.
1049  ! We do this by finding points +/- one std dev in mobility
1050  ! radius, converting them to geometric radius, and then using
1051  ! the distance between them as a measure of 2 * std_dev in
1052  ! geometric radius.
1053  log10_r1_mob = log10(diam2rad(diam))
1054  log10_r0_mob = log10_r1_mob - aero_mode%log10_std_dev_radius
1055  log10_r2_mob = log10_r1_mob + aero_mode%log10_std_dev_radius
1056  r0_mob = 10**log10_r0_mob
1057  r2_mob = 10**log10_r2_mob
1058  r0_geom = aero_data_mobility_rad_to_geometric_rad(aero_data, &
1059  r0_mob, temp, pressure)
1060  r2_geom = aero_data_mobility_rad_to_geometric_rad(aero_data, &
1061  r2_mob, temp, pressure)
1062  log10_r0_geom = log10(r0_geom)
1063  log10_r2_geom = log10(r2_geom)
1064  aero_mode%log10_std_dev_radius &
1065  = (log10_r2_geom - log10_r0_geom) / 2d0
1066  else
1067  call die_msg(532966100, "Diameter type not handled: " &
1068  // integer_to_string(diam_type))
1069  end if
1070  elseif (trim(mode_type) == 'exp') then
1071  aero_mode%type = aero_mode_type_exp
1072  call spec_file_read_real(file, 'num_conc', aero_mode%num_conc)
1073  call spec_file_read_real(file, 'diam_at_mean_vol', diam)
1074  if (diam_type == aero_mode_diam_type_geometric) then
1075  aero_mode%char_radius = diam2rad(diam)
1076  elseif (diam_type == aero_mode_diam_type_mobility) then
1077  aero_mode%char_radius &
1079  diam2rad(diam), temp, pressure)
1080  else
1081  call die_msg(585104460, "Diameter type not handled: " &
1082  // integer_to_string(diam_type))
1083  end if
1084  elseif (trim(mode_type) == 'mono') then
1085  aero_mode%type = aero_mode_type_mono
1086  call spec_file_read_real(file, 'num_conc', aero_mode%num_conc)
1087  call spec_file_read_real(file, 'diam', diam)
1088  if (diam_type == aero_mode_diam_type_geometric) then
1089  aero_mode%char_radius = diam2rad(diam)
1090  elseif (diam_type == aero_mode_diam_type_mobility) then
1091  aero_mode%char_radius &
1093  diam2rad(diam), temp, pressure)
1094  else
1095  call die_msg(902864269, "Diameter type not handled: " &
1096  // integer_to_string(diam_type))
1097  end if
1098  elseif (trim(mode_type) == 'sampled') then
1099  aero_mode%type = aero_mode_type_sampled
1100  call spec_file_read_string(file, 'size_dist', size_dist_filename)
1101  call spec_file_open(size_dist_filename, size_dist_file)
1102  call spec_file_read_size_dist(size_dist_file, &
1103  aero_mode%sample_radius, aero_mode%sample_num_conc)
1104  call spec_file_close(size_dist_file)
1105  if (diam_type == aero_mode_diam_type_geometric) then
1106  ! do nothing
1107  elseif (diam_type == aero_mode_diam_type_mobility) then
1108  do i_radius = 1,size(aero_mode%sample_radius)
1109  aero_mode%sample_radius(i_radius) &
1111  aero_mode%sample_radius(i_radius), temp, pressure)
1112  end do
1113  else
1114  call die_msg(239088838, "Diameter type not handled: " &
1115  // integer_to_string(diam_type))
1116  end if
1117  else
1118  call spec_file_die_msg(729472928, file, &
1119  "Unknown distribution mode type: " // trim(mode_type))
1120  end if
1121  end if
1122 
1123  end subroutine spec_file_read_aero_mode
1124 
1125 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1126 
1127  !> Determines the number of bytes required to pack the given value.
1128  integer function pmc_mpi_pack_size_aero_mode(val)
1129 
1130  !> Value to pack.
1131  type(aero_mode_t), intent(in) :: val
1132 
1134  pmc_mpi_pack_size_string(val%name) &
1135  + pmc_mpi_pack_size_integer(val%type) &
1136  + pmc_mpi_pack_size_real(val%char_radius) &
1137  + pmc_mpi_pack_size_real(val%log10_std_dev_radius) &
1138  + pmc_mpi_pack_size_real_array(val%sample_radius) &
1139  + pmc_mpi_pack_size_real_array(val%sample_num_conc) &
1140  + pmc_mpi_pack_size_real(val%num_conc) &
1141  + pmc_mpi_pack_size_real_array(val%vol_frac) &
1142  + pmc_mpi_pack_size_real_array(val%vol_frac_std) &
1143  + pmc_mpi_pack_size_integer(val%source) &
1144  + pmc_mpi_pack_size_integer(val%weight_class)
1145 
1146  end function pmc_mpi_pack_size_aero_mode
1147 
1148 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1149 
1150  !> Packs the given value into the buffer, advancing position.
1151  subroutine pmc_mpi_pack_aero_mode(buffer, position, val)
1152 
1153  !> Memory buffer.
1154  character, intent(inout) :: buffer(:)
1155  !> Current buffer position.
1156  integer, intent(inout) :: position
1157  !> Value to pack.
1158  type(aero_mode_t), intent(in) :: val
1159 
1160 #ifdef PMC_USE_MPI
1161  integer :: prev_position
1162 
1163  prev_position = position
1164  call pmc_mpi_pack_string(buffer, position, val%name)
1165  call pmc_mpi_pack_integer(buffer, position, val%type)
1166  call pmc_mpi_pack_real(buffer, position, val%char_radius)
1167  call pmc_mpi_pack_real(buffer, position, val%log10_std_dev_radius)
1168  call pmc_mpi_pack_real_array(buffer, position, val%sample_radius)
1169  call pmc_mpi_pack_real_array(buffer, position, val%sample_num_conc)
1170  call pmc_mpi_pack_real(buffer, position, val%num_conc)
1171  call pmc_mpi_pack_real_array(buffer, position, val%vol_frac)
1172  call pmc_mpi_pack_real_array(buffer, position, val%vol_frac_std)
1173  call pmc_mpi_pack_integer(buffer, position, val%source)
1174  call pmc_mpi_pack_integer(buffer, position, val%weight_class)
1175  call assert(497092471, &
1176  position - prev_position <= pmc_mpi_pack_size_aero_mode(val))
1177 #endif
1178 
1179  end subroutine pmc_mpi_pack_aero_mode
1180 
1181 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1182 
1183  !> Unpacks the given value from the buffer, advancing position.
1184  subroutine pmc_mpi_unpack_aero_mode(buffer, position, val)
1185 
1186  !> Memory buffer.
1187  character, intent(inout) :: buffer(:)
1188  !> Current buffer position.
1189  integer, intent(inout) :: position
1190  !> Value to pack.
1191  type(aero_mode_t), intent(inout) :: val
1192 
1193 #ifdef PMC_USE_MPI
1194  integer :: prev_position
1195 
1196  prev_position = position
1197  call pmc_mpi_unpack_string(buffer, position, val%name)
1198  call pmc_mpi_unpack_integer(buffer, position, val%type)
1199  call pmc_mpi_unpack_real(buffer, position, val%char_radius)
1200  call pmc_mpi_unpack_real(buffer, position, val%log10_std_dev_radius)
1201  call pmc_mpi_unpack_real_array(buffer, position, val%sample_radius)
1202  call pmc_mpi_unpack_real_array(buffer, position, val%sample_num_conc)
1203  call pmc_mpi_unpack_real(buffer, position, val%num_conc)
1204  call pmc_mpi_unpack_real_array(buffer, position, val%vol_frac)
1205  call pmc_mpi_unpack_real_array(buffer, position, val%vol_frac_std)
1206  call pmc_mpi_unpack_integer(buffer, position, val%source)
1207  call pmc_mpi_unpack_integer(buffer, position, val%weight_class)
1208  call assert(874467577, &
1209  position - prev_position <= pmc_mpi_pack_size_aero_mode(val))
1210 #endif
1211 
1212  end subroutine pmc_mpi_unpack_aero_mode
1213 
1214 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1215 
1216 end module pmc_aero_mode
pmc_aero_mode::num_conc_exp
subroutine num_conc_exp(total_num_conc, radius_at_mean_vol, bin_grid, aero_data, num_conc)
Exponential distribution in volume.
Definition: aero_mode.F90:198
pmc_aero_mode::aero_mode_type_len
integer, parameter aero_mode_type_len
Maximum length of a mode type.
Definition: aero_mode.F90:26
pmc_aero_mode::vol_conc_mono
subroutine vol_conc_mono(total_num_conc, radius, bin_grid, aero_data, vol_conc)
Mono-disperse distribution in volume.
Definition: aero_mode.F90:280
pmc_aero_weight::aero_weight_type_none
integer, parameter aero_weight_type_none
Type code for no (or flat) weighting.
Definition: aero_weight.F90:26
pmc_aero_data::aero_data_n_spec
elemental integer function aero_data_n_spec(aero_data)
Return the number of aerosol species, or -1 if uninitialized.
Definition: aero_data.F90:246
pmc_mpi
Wrapper functions for MPI.
Definition: mpi.F90:13
pmc_aero_mode::aero_mode_type_to_string
character(len=aero_mode_type_len) function aero_mode_type_to_string(type)
Return a string representation of a kernel type.
Definition: aero_mode.F90:86
pmc_aero_mode::aero_mode_diam_type_mobility
integer, parameter aero_mode_diam_type_mobility
Type code for mobility equivalent diameter.
Definition: aero_mode.F90:44
pmc_aero_weight
The aero_weight_t structure and associated subroutines.
Definition: aero_weight.F90:9
pmc_mpi::pmc_mpi_pack_size_real
integer function pmc_mpi_pack_size_real(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:385
pmc_aero_mode::num_conc_mono
subroutine num_conc_mono(total_num_conc, radius, bin_grid, num_conc)
Mono-disperse distribution. Normalized so that sum(num_conc(k) * log_width) = 1.
Definition: aero_mode.F90:253
pmc_aero_mode::aero_mode_sample_vols
subroutine aero_mode_sample_vols(aero_mode, total_vol, vols)
Return an array of volumes randomly sampled from the volume fractions.
Definition: aero_mode.F90:634
pmc_spec_file::spec_file_close
subroutine spec_file_close(file)
Close a spec file.
Definition: spec_file.F90:135
pmc_bin_grid::vol_to_lnr
subroutine vol_to_lnr(r, f_vol, f_lnr)
Convert a concentration f(vol)d(vol) to f(ln(r))d(ln(r)) where vol = 4/3 pi r^3.
Definition: bin_grid.F90:68
pmc_aero_mode::aero_mode_diam_type_invalid
integer, parameter aero_mode_diam_type_invalid
Type code for an undefined for invalid diameter type.
Definition: aero_mode.F90:40
pmc_aero_mode::vol_conc_exp
subroutine vol_conc_exp(total_num_conc, radius_at_mean_vol, bin_grid, aero_data, vol_conc)
Exponential distribution in volume.
Definition: aero_mode.F90:228
pmc_aero_data::aero_data_source_by_name
integer function aero_data_source_by_name(aero_data, name)
Returns the number of the source in aero_data with the given name, or adds the source if it doesn't e...
Definition: aero_data.F90:324
pmc_aero_mode::aero_mode_type_sampled
integer, parameter aero_mode_type_sampled
Type code for a sampled mode.
Definition: aero_mode.F90:37
pmc_aero_mode::aero_mode_type_exp
integer, parameter aero_mode_type_exp
Type code for an exponential mode.
Definition: aero_mode.F90:33
pmc_util::die_msg
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:135
pmc_aero_mode::vol_conc_sampled
subroutine vol_conc_sampled(sample_radius, sample_num_conc, bin_grid, aero_data, vol_conc)
Sampled distribution in volume.
Definition: aero_mode.F90:356
pmc_constants::dp
integer, parameter dp
Kind of a double precision real number.
Definition: constants.F90:12
pmc_aero_mode::pmc_mpi_pack_size_aero_mode
integer function pmc_mpi_pack_size_aero_mode(val)
Determines the number of bytes required to pack the given value.
Definition: aero_mode.F90:1129
pmc_mpi::pmc_mpi_pack_string
subroutine pmc_mpi_pack_string(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:786
pmc_spec_file
Reading formatted text input.
Definition: spec_file.F90:43
pmc_aero_data::aero_data_weight_class_by_name
integer function aero_data_weight_class_by_name(aero_data, name)
Returns the number of the weight class with the given name, or adds the weight class if it doesn't ex...
Definition: aero_data.F90:348
pmc_bin_grid::bin_grid_size
elemental integer function bin_grid_size(bin_grid)
Return the number of bins in the grid, or -1 if the bin grid is not allocated.
Definition: bin_grid.F90:51
pmc_mpi::pmc_mpi_pack_real
subroutine pmc_mpi_pack_real(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:761
pmc_util::assert
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:104
pmc_aero_mode::vol_conc_log_normal
subroutine vol_conc_log_normal(total_num_conc, geom_mean_radius, log10_sigma_g, bin_grid, aero_data, vol_conc)
Compute a log-normal distribution in volume.
Definition: aero_mode.F90:167
pmc_aero_mode::aero_mode_num_conc
subroutine aero_mode_num_conc(aero_mode, bin_grid, aero_data, num_conc)
Return the binned number concentration for an aero_mode.
Definition: aero_mode.F90:380
pmc_mpi::pmc_mpi_pack_size_real_array
integer function pmc_mpi_pack_size_real_array(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:527
pmc_aero_mode::aero_mode_total_num_conc
real(kind=dp) function aero_mode_total_num_conc(aero_mode)
Returns the total number concentration of a mode. (#/m^3)
Definition: aero_mode.F90:110
pmc_mpi::pmc_mpi_unpack_string
subroutine pmc_mpi_unpack_string(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1214
pmc_aero_mode::aero_mode_type_mono
integer, parameter aero_mode_type_mono
Type code for a mono-disperse mode.
Definition: aero_mode.F90:35
pmc_spec_file::spec_file_t
An input file with extra data for printing messages.
Definition: spec_file.F90:59
pmc_util::warn_msg
subroutine warn_msg(code, warning_msg, already_warned)
Prints a warning message.
Definition: util.F90:38
pmc_spec_file::spec_file_read_real
subroutine spec_file_read_real(file, name, var)
Read a real number from a spec file that must have the given name.
Definition: spec_file.F90:563
pmc_aero_weight::aero_weight_type_power
integer, parameter aero_weight_type_power
Type code for power function weighting.
Definition: aero_weight.F90:28
pmc_util::integer_to_string
character(len=pmc_util_convert_string_len) function integer_to_string(val)
Convert an integer to a string format.
Definition: util.F90:767
pmc_util::assert_msg
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
Definition: util.F90:78
pmc_aero_weight::aero_weight_type_mfa
integer, parameter aero_weight_type_mfa
Type code for MFA weighting.
Definition: aero_weight.F90:30
pmc_util::diam2rad
real(kind=dp) elemental function diam2rad(d)
Convert diameter (m) to radius (m).
Definition: util.F90:278
pmc_spec_file::spec_file_check_name
subroutine spec_file_check_name(file, name, read_name)
Checks that the read_name is the same as name.
Definition: spec_file.F90:410
pmc_spec_file::spec_file_check_line_length
subroutine spec_file_check_line_length(file, line, length)
Check that the length of the line data is as given.
Definition: spec_file.F90:432
pmc_aero_weight::aero_weight_num_conc_at_radius
real(kind=dp) function aero_weight_num_conc_at_radius(aero_weight, radius)
Compute the number concentration at a given radius (m^{-3}).
Definition: aero_weight.F90:128
pmc_aero_mode::aero_mode_t
An aerosol size distribution mode.
Definition: aero_mode.F90:54
pmc_aero_data::aero_data_mobility_rad_to_geometric_rad
real(kind=dp) function aero_data_mobility_rad_to_geometric_rad(aero_data, mobility_rad, temp, pressure)
Convert mobility equivalent radius (m) to geometric radius (m^3).
Definition: aero_data.F90:226
pmc_constants::const
type(const_t), save const
Fixed variable for accessing the constant's values.
Definition: constants.F90:75
pmc_rand::sample_cts_pdf
integer function sample_cts_pdf(pdf)
Sample the given continuous probability density function.
Definition: rand.F90:551
pmc_aero_mode::aero_mode_type_invalid
integer, parameter aero_mode_type_invalid
Type code for an undefined or invalid mode.
Definition: aero_mode.F90:29
pmc_aero_mode::pmc_mpi_unpack_aero_mode
subroutine pmc_mpi_unpack_aero_mode(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: aero_mode.F90:1185
pmc_spec_file::spec_file_open
subroutine spec_file_open(filename, file)
Open a spec file for reading.
Definition: spec_file.F90:112
pmc_aero_weight::aero_weight_t
An aerosol size distribution weighting function.
Definition: aero_weight.F90:33
pmc_aero_mode::aero_mode_name_len
integer, parameter aero_mode_name_len
Maximum length of a mode name.
Definition: aero_mode.F90:24
pmc_rand
Random number generators.
Definition: rand.F90:9
pmc_spec_file::spec_file_read_real_named_array
subroutine spec_file_read_real_named_array(file, max_lines, names, vals)
Read an array of named lines with real data. All lines must have the same number of data elements.
Definition: spec_file.F90:650
pmc_aero_mode::num_conc_sampled
subroutine num_conc_sampled(sample_radius, sample_num_conc, bin_grid, num_conc)
Sampled distribution, not normalized.
Definition: aero_mode.F90:310
pmc_aero_mode::aero_mode_diam_type_geometric
integer, parameter aero_mode_diam_type_geometric
Type code for geometric diameter.
Definition: aero_mode.F90:42
pmc_aero_data::aero_data_vol2rad
real(kind=dp) elemental function aero_data_vol2rad(aero_data, v)
Convert mass-equivalent volume (m^3) to geometric radius (m).
Definition: aero_data.F90:102
pmc_rand::rand_normal_array_1d
subroutine rand_normal_array_1d(mean, stddev, val)
Generates a vector of normally distributed random numbers with the given means and standard deviation...
Definition: rand.F90:495
pmc_aero_mode::aero_mode_sample_radius
subroutine aero_mode_sample_radius(aero_mode, aero_data, aero_weight, radius, prob_threshold)
Return a radius randomly sampled from the mode distribution.
Definition: aero_mode.F90:551
pmc_aero_data::aero_data_t
Aerosol material properties and associated data.
Definition: aero_data.F90:55
pmc_mpi::pmc_mpi_unpack_integer
subroutine pmc_mpi_unpack_integer(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1139
pmc_constants
Physical constants.
Definition: constants.F90:9
pmc_spec_file::spec_file_read_line
subroutine spec_file_read_line(file, line, eof)
Read a spec_line from the spec_file.
Definition: spec_file.F90:215
pmc_aero_mode::pmc_mpi_pack_aero_mode
subroutine pmc_mpi_pack_aero_mode(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: aero_mode.F90:1152
pmc_aero_mode::aero_mode_type_log_normal
integer, parameter aero_mode_type_log_normal
Type code for a log-normal mode.
Definition: aero_mode.F90:31
pmc_util
Common utility subroutines.
Definition: util.F90:9
pmc_aero_mode::aero_mode_number
real(kind=dp) function aero_mode_number(aero_mode, aero_weight)
Return the total number of computational particles for an aero_mode.
Definition: aero_mode.F90:497
pmc_spec_file::spec_file_die_msg
subroutine spec_file_die_msg(code, file, msg)
Exit with an error message containing filename and line number.
Definition: spec_file.F90:74
pmc_aero_mode
The aero_mode_t structure and associated subroutines.
Definition: aero_mode.F90:9
pmc_spec_file::spec_file_read_string
subroutine spec_file_read_string(file, name, var)
Read a string from a spec file that must have a given name.
Definition: spec_file.F90:605
pmc_aero_mode::aero_mode_weighted_sampled_num_conc
subroutine aero_mode_weighted_sampled_num_conc(aero_mode, aero_weight, weighted_num_conc)
Compute weighted sampled number concentrations.
Definition: aero_mode.F90:460
pmc_mpi::pmc_mpi_pack_size_integer
integer function pmc_mpi_pack_size_integer(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:345
pmc_bin_grid
The bin_grid_t structure and associated subroutines.
Definition: bin_grid.F90:9
pmc_mpi::pmc_mpi_pack_integer
subroutine pmc_mpi_pack_integer(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:711
pmc_aero_data::aero_data_rad2vol
real(kind=dp) elemental function aero_data_rad2vol(aero_data, r)
Convert geometric radius (m) to mass-equivalent volume (m^3).
Definition: aero_data.F90:132
pmc_aero_data
The aero_data_t structure and associated subroutines.
Definition: aero_data.F90:9
pmc_mpi::pmc_mpi_pack_real_array
subroutine pmc_mpi_pack_real_array(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:927
pmc_spec_file::spec_file_check_line_name
subroutine spec_file_check_line_name(file, line, name)
Check that the name of the line data is as given.
Definition: spec_file.F90:390
pmc_bin_grid::bin_grid_t
1D grid, either logarithmic or linear.
Definition: bin_grid.F90:33
pmc_spec_file::spec_file_assert_msg
subroutine spec_file_assert_msg(code, file, condition_ok, msg)
Exit with an error message containing filename and line number if condition_ok is ....
Definition: spec_file.F90:91
pmc_aero_mode::num_conc_log_normal
subroutine num_conc_log_normal(total_num_conc, geom_mean_radius, log10_sigma_g, bin_grid, num_conc)
Compute a log-normal distribution.
Definition: aero_mode.F90:133
pmc_mpi::pmc_mpi_unpack_real
subroutine pmc_mpi_unpack_real(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1189
pmc_mpi::pmc_mpi_unpack_real_array
subroutine pmc_mpi_unpack_real_array(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1356
pmc_rand::rand_normal
real(kind=dp) function rand_normal(mean, stddev, prob_threshold)
Generates a normally distributed random number with the given mean and standard deviation.
Definition: rand.F90:411
pmc_aero_data::aero_data_spec_by_name
integer function aero_data_spec_by_name(aero_data, name)
Returns the number of the species in aero_data with the given name, or returns 0 if there is no such ...
Definition: aero_data.F90:295
pmc_aero_mode::aero_mode_get_weight_class
integer function aero_mode_get_weight_class(aero_mode)
Get the weight class for an aero_mode.
Definition: aero_mode.F90:671
pmc_mpi::pmc_mpi_pack_size_string
integer function pmc_mpi_pack_size_string(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:405
pmc_rand::pmc_random
real(kind=dp) function pmc_random()
Returns a random number between 0 and 1.
Definition: rand.F90:139
pmc_bin_grid::bin_grid_find
integer function bin_grid_find(bin_grid, val)
Find the bin number that contains a given value.
Definition: bin_grid.F90:144
pmc_aero_weight::aero_weight_radius_at_num_conc
real(kind=dp) function aero_weight_radius_at_num_conc(aero_weight, num_conc)
Compute the radius at a given number concentration (m). Inverse of aero_weight_num_conc_at_radius().
Definition: aero_weight.F90:171
pmc_aero_mode::aero_mode_vol_conc
subroutine aero_mode_vol_conc(aero_mode, bin_grid, aero_data, vol_conc)
Return the binned per-species volume concentration for an aero_mode.
Definition: aero_mode.F90:415