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