PartMC  2.6.1
aero_dist.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2012 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_dist module.
7 
8 !> The aero_dist_t structure and associated subroutines.
9 !!
10 !! The initial size distributions are computed as number densities, so
11 !! they can be used for both sectional and particle-resolved
12 !! simulations. The routine dist_to_n() converts a number concentration
13 !! distribution to an actual number of particles ready for a
14 !! particle-resolved simulation.
15 !!
16 !! Initial distributions should be normalized so that <tt>sum(n_den) =
17 !! 1/log_width</tt>.
19 
20  use pmc_bin_grid
21  use pmc_util
22  use pmc_constants
23  use pmc_spec_file
24  use pmc_aero_data
25  use pmc_aero_mode
26  use pmc_mpi
27  use pmc_rand
28 #ifdef PMC_USE_MPI
29  use mpi
30 #endif
31 
32  !> A complete aerosol distribution, consisting of several modes.
34  !> Internally mixed modes.
35  type(aero_mode_t), allocatable :: mode(:)
36  end type aero_dist_t
37 
38 contains
39 
40 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
41 
42  !> Return the number of modes.
43  elemental integer function aero_dist_n_mode(aero_dist)
44 
45  !> Aero data structure.
46  type(aero_dist_t), intent(in) :: aero_dist
47 
48  if (allocated(aero_dist%mode)) then
49  aero_dist_n_mode = size(aero_dist%mode)
50  else
52  end if
53 
54  end function aero_dist_n_mode
55 
56 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
57 
58  !> Returns the total number concentration of a distribution. (#/m^3)
59  real(kind=dp) function aero_dist_total_num_conc(aero_dist)
60 
61  !> Aerosol distribution.
62  type(aero_dist_t), intent(in) :: aero_dist
63 
64  integer :: i_mode
65 
67  do i_mode = 1,aero_dist_n_mode(aero_dist)
69  + aero_mode_total_num_conc(aero_dist%mode(i_mode))
70  end do
71 
72  end function aero_dist_total_num_conc
73 
74 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
75 
76  !> Returns the total number of particles of a distribution.
77  real(kind=dp) function aero_dist_number(aero_dist, aero_weight)
78 
79  !> Aerosol distribution.
80  type(aero_dist_t), intent(in) :: aero_dist
81  !> Aerosol weight.
82  type(aero_weight_t), intent(in) :: aero_weight
83 
84  integer :: i_mode
85 
86  aero_dist_number = 0d0
87  do i_mode = 1,aero_dist_n_mode(aero_dist)
89  + aero_mode_number(aero_dist%mode(i_mode), aero_weight)
90  end do
91 
92  end function aero_dist_number
93 
94 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
95 
96  !> Return the binned number concentration for an aero_dist.
97  subroutine aero_dist_num_conc(aero_dist, bin_grid, aero_data, &
98  num_conc)
99 
100  !> Aero dist for which to compute number concentration.
101  type(aero_dist_t), intent(in) :: aero_dist
102  !> Bin grid.
103  type(bin_grid_t), intent(in) :: bin_grid
104  !> Aerosol data.
105  type(aero_data_t), intent(in) :: aero_data
106  !> Number concentration (#(ln(r))d(ln(r))).
107  real(kind=dp), intent(out) :: num_conc(bin_grid_size(bin_grid))
108 
109  integer :: i_mode
110  real(kind=dp) :: mode_num_conc(size(num_conc, 1))
111 
112  num_conc = 0d0
113  do i_mode = 1,aero_dist_n_mode(aero_dist)
114  call aero_mode_num_conc(aero_dist%mode(i_mode), bin_grid, &
115  aero_data, mode_num_conc)
116  num_conc = num_conc + mode_num_conc
117  end do
118 
119  end subroutine aero_dist_num_conc
120 
121 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
122 
123  !> Return the binned per-species volume concentration for an
124  !> aero_dist.
125  subroutine aero_dist_vol_conc(aero_dist, bin_grid, aero_data, &
126  vol_conc)
127 
128  !> Aero dist for which to compute volume concentration.
129  type(aero_dist_t), intent(in) :: aero_dist
130  !> Bin grid.
131  type(bin_grid_t), intent(in) :: bin_grid
132  !> Aerosol data.
133  type(aero_data_t), intent(in) :: aero_data
134  !> Volume concentration (V(ln(r))d(ln(r))).
135  real(kind=dp), intent(out) :: vol_conc(bin_grid_size(bin_grid), &
136  aero_data_n_spec(aero_data))
137 
138  integer :: i_mode
139  real(kind=dp) :: mode_vol_conc(size(vol_conc, 1), size(vol_conc, 2))
140 
141  vol_conc = 0d0
142  do i_mode = 1,aero_dist_n_mode(aero_dist)
143  call aero_mode_vol_conc(aero_dist%mode(i_mode), bin_grid, &
144  aero_data, mode_vol_conc)
145  vol_conc = vol_conc + mode_vol_conc
146  end do
147 
148  end subroutine aero_dist_vol_conc
149 
150 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
151 
152  !> Whether any of the modes are of the given type.
153  elemental logical function aero_dist_contains_aero_mode_type(aero_dist, &
154  aero_mode_type)
155 
156  !> Aerosol distribution.
157  type(aero_dist_t), intent(in) :: aero_dist
158  !> Aerosol mode type to test for.
159  integer, intent(in) :: aero_mode_type
160 
161  integer :: i_mode
162 
164  do i_mode = 1,aero_dist_n_mode(aero_dist)
166  .or. (aero_dist%mode(i_mode)%type == aero_mode_type)
167  end do
168 
170 
171 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
172 
173  !> Determine the current aero_dist and rate by interpolating at the
174  !> current time with the lists of aero_dists and rates.
175  subroutine aero_dist_interp_1d(aero_dist_list, time_list, &
176  rate_list, time, aero_dist, rate)
177 
178  !> Gas states.
179  type(aero_dist_t), intent(in) :: aero_dist_list(:)
180  !> Times (s).
181  real(kind=dp), intent(in) :: time_list(size(aero_dist_list))
182  !> Rates (s^{-1}).
183  real(kind=dp), intent(in) :: rate_list(size(aero_dist_list))
184  !> Current time (s).
185  real(kind=dp), intent(in) :: time
186  !> Current gas state.
187  type(aero_dist_t), intent(inout) :: aero_dist
188  !> Current rate (s^{-1}).
189  real(kind=dp), intent(out) :: rate
190 
191  integer :: n, p, n_bin, n_spec, i, i_new
192  real(kind=dp) :: y, alpha
193 
194  n = size(aero_dist_list)
195  p = find_1d(n, time_list, time)
196  if (p == 0) then
197  ! before the start, just use the first state and rate
198  aero_dist = aero_dist_list(1)
199  rate = rate_list(1)
200  elseif (p == n) then
201  ! after the end, just use the last state and rate
202  aero_dist = aero_dist_list(n)
203  rate = rate_list(n)
204  else
205  ! in the middle, use the previous dist
206  aero_dist = aero_dist_list(p)
207  rate = rate_list(p)
208  end if
209 
210  end subroutine aero_dist_interp_1d
211 
212 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
213 
214  !> Read continuous aerosol distribution composed of several modes.
215  subroutine spec_file_read_aero_dist(file, aero_data, aero_dist)
216 
217  !> Spec file to read data from.
218  type(spec_file_t), intent(inout) :: file
219  !> Aero_data data.
220  type(aero_data_t), intent(inout) :: aero_data
221  !> Aerosol dist.
222  type(aero_dist_t), intent(inout) :: aero_dist
223 
224  type(aero_mode_t) :: aero_mode
225  type(aero_mode_t), allocatable :: new_modes(:)
226  integer :: n, i
227  logical :: eof
228 
229  ! note that the <p> is needed below to force correct paragraph
230  ! breaking by doxygen
231 
232  !> \page input_format_aero_dist Input File Format: Aerosol Distribution
233  !!
234  !! <p>An aerosol distribution file consists of zero or more modes,
235  !! each in the format described by \subpage input_format_aero_mode
236  !!
237  !! See also:
238  !! - \ref spec_file_format --- the input file text format
239  !! - \ref input_format_aero_mode --- the format for each mode
240  !! of an aerosol distribution
241 
242  allocate(aero_dist%mode(0))
243  call spec_file_read_aero_mode(file, aero_data, aero_mode, eof)
244  do while (.not. eof)
245  n = size(aero_dist%mode)
246  allocate(new_modes(n + 1))
247  do i = 1,n
248  new_modes(i) = aero_dist%mode(i)
249  end do
250  new_modes(n + 1) = aero_mode
251  ! Next line is commented out due to a reported memory leak with valgrind
252  ! using both gcc 4.6.3 and 4.9.2.
253  !call move_alloc(new_modes, aero_dist%mode)
254  ! Next two lines are here to avoid issues with reallocation of lhs by
255  ! assignment.
256  deallocate(aero_dist%mode)
257  allocate(aero_dist%mode(n+1))
258  aero_dist%mode = new_modes
259  deallocate(new_modes)
260  ! Next line appears to work in gcc 4.9.2 but led to problems with older
261  ! versions of gcc (4.6.3). This will be ideal in the future as it avoids
262  ! the copying of the aero_dist%mode array each loop.
263  !aero_dist%mode = [aero_dist%mode, aero_mode]
264  call spec_file_read_aero_mode(file, aero_data, aero_mode, eof)
265  end do
266 
267  end subroutine spec_file_read_aero_dist
268 
269 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
270 
271  !> Read an array of aero_dists with associated times and rates from
272  !> the given file.
273  subroutine spec_file_read_aero_dists_times_rates(file, aero_data, &
274  times, rates, aero_dists)
275 
276  !> Spec file to read data from.
277  type(spec_file_t), intent(inout) :: file
278  !> Aero data.
279  type(aero_data_t), intent(inout) :: aero_data
280  !> Times (s).
281  real(kind=dp), allocatable :: times(:)
282  !> Rates (s^{-1}).
283  real(kind=dp), allocatable :: rates(:)
284  !> Aero dists.
285  type(aero_dist_t), allocatable :: aero_dists(:)
286 
287  type(spec_line_t) :: aero_dist_line
288  type(spec_file_t) :: aero_dist_file
289  integer :: n_time, i_time
290  character(len=SPEC_LINE_MAX_VAR_LEN), allocatable :: names(:)
291  real(kind=dp), allocatable :: data(:,:)
292 
293  !> \page input_format_aero_dist_profile Input File Format: Aerosol Distribution Profile
294  !!
295  !! An aerosol distribution profile input file must consist of
296  !! three lines:
297  !! - the first line must begin with \c time and should be followed
298  !! by \f$N\f$ space-separated real scalars, giving the times (in
299  !! s after the start of the simulation) of the aerosol
300  !! distrbution set points --- the times must be in increasing
301  !! order
302  !! - the second line must begin with \c rate and should be
303  !! followed by \f$N\f$ space-separated real scalars, giving the
304  !! values at the corresponding times
305  !! - the third line must begin with \c dist and should be followed
306  !! by \f$N\f$ space-separated filenames, each specifying an
307  !! aerosol distribution in the format \ref input_format_aero_dist
308  !! at the corresponding time
309  !!
310  !! The units of the \c rate line depend on the type of aerosol
311  !! distribution profile:
312  !! - Emissions aerosol profiles have rates with units m/s ---
313  !! the aerosol distribution number concentrations are multiplied
314  !! by the rate to give an emission rate with unit #/(m^2 s)
315  !! which is then divided by the current mixing layer height
316  !! to give a per-volume emission rate.
317  !! - Background aerosol profiles have rates with units \f$\rm
318  !! s^{-1}\f$, which is the dilution rate between the background
319  !! and the simulated air parcel. That is, if the simulated
320  !! number concentration is \f$N\f$ and the background number
321  !! concentration is \f$N_{\rm back}\f$, then dilution is modeled
322  !! as \f$\dot{N} = r N_{\rm back} - r N\f$, where \f$r\f$ is the
323  !! rate.
324  !!
325  !! Between the specified times the aerosol profile is interpolated
326  !! step-wise and kept constant at its last value. That is, if the
327  !! times are \f$t_i\f$, the rates are \f$r_i\f$, and the aerosol
328  !! distributions are \f$a_i\f$ (all with \f$i = 1,\ldots,n\f$),
329  !! then between times \f$t_i\f$ and \f$t_{i+1}\f$ the aerosol
330  !! state is constant at \f$r_i a_i\f$. Before time \f$t_1\f$ the
331  !! aerosol state is \f$r_1 a_1\f$, while after time \f$t_n\f$ it
332  !! is \f$r_n a_n\f$.
333  !!
334  !! Example: an emissions aerosol profile could be:
335  !! <pre>
336  !! time 0 600 1800 # time (in s) after sim start
337  !! rate 1 0.5 1 # scaling factor in m/s
338  !! dist dist1.dat dist2.dat dist3.dat # aerosol distribution files
339  !! </pre>
340  !! Here the emissions between 0&nbsp;min and 10&nbsp;min are given
341  !! by <tt>dist1.dat</tt> (with the number concentration
342  !! interpreted as having units 1/(m^2 s)), the emissions between
343  !! 10&nbsp;min and 30&nbsp;min are given by <tt>dist2.dat</tt>
344  !! (scaled by 0.5), while the emissions after 30&nbsp;min are
345  !! given by <tt>dist3.dat</tt>.
346  !!
347  !! See also:
348  !! - \ref spec_file_format --- the input file text format
349  !! - \ref input_format_aero_data --- the aerosol species list
350  !! and material data
351  !! - \ref input_format_aero_dist --- the format of the
352  !! instantaneous aerosol distribution files
353 
354  ! read the data from the file
355  call spec_file_read_real_named_array(file, 2, names, data)
356  call spec_file_read_line_no_eof(file, aero_dist_line)
357  call spec_file_check_line_name(file, aero_dist_line, "dist")
358  call spec_file_check_line_length(file, aero_dist_line, size(data, 2))
359 
360  ! check the data size
361  if (size(names) /= 2) then
362  call die_msg(530745700, &
363  trim(file%name) // ' must contain exactly two data lines')
364  end if
365  if (trim(names(1)) /= 'time') then
366  call die_msg(570205795, 'row 1 in ' // trim(file%name) &
367  // ' must start with: time not: ' // trim(names(1)))
368  end if
369  if (trim(names(2)) /= 'rate') then
370  call die_msg(221270915, 'row 2 in ' // trim(file%name) &
371  // ' must start with: rate not: ' // trim(names(1)))
372  end if
373  n_time = size(data, 2)
374  if (n_time < 1) then
375  call die_msg(457229710, 'each line in ' // trim(file%name) &
376  // ' must contain at least one data value')
377  end if
378 
379  ! copy over the data
380  times = data(1,:)
381  rates = data(2,:)
382  if (allocated(aero_dists)) deallocate(aero_dists)
383  allocate(aero_dists(n_time))
384  do i_time = 1,n_time
385  call spec_file_open(aero_dist_line%data(i_time), aero_dist_file)
386  call spec_file_read_aero_dist(aero_dist_file, &
387  aero_data, aero_dists(i_time))
388  call spec_file_close(aero_dist_file)
389  end do
390 
391  end subroutine spec_file_read_aero_dists_times_rates
392 
393 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
394 
395  !> Determines the number of bytes required to pack the given value.
396  integer function pmc_mpi_pack_size_aero_dist(val)
397 
398  !> Value to pack.
399  type(aero_dist_t), intent(in) :: val
400 
401  integer :: i, total_size
402 
403  if (allocated(val%mode)) then
405  do i = 1,size(val%mode)
406  total_size = total_size + pmc_mpi_pack_size_aero_mode(val%mode(i))
407  end do
408  else
409  total_size = pmc_mpi_pack_size_integer(-1)
410  end if
411  pmc_mpi_pack_size_aero_dist = total_size
412 
413  end function pmc_mpi_pack_size_aero_dist
414 
415 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
416 
417  !> Packs the given value into the buffer, advancing position.
418  subroutine pmc_mpi_pack_aero_dist(buffer, position, val)
419 
420  !> Memory buffer.
421  character, intent(inout) :: buffer(:)
422  !> Current buffer position.
423  integer, intent(inout) :: position
424  !> Value to pack.
425  type(aero_dist_t), intent(in) :: val
426 
427 #ifdef PMC_USE_MPI
428  integer :: prev_position, i
429 
430  prev_position = position
431  if (allocated(val%mode)) then
432  call pmc_mpi_pack_integer(buffer, position, aero_dist_n_mode(val))
433  do i = 1,size(val%mode)
434  call pmc_mpi_pack_aero_mode(buffer, position, val%mode(i))
435  end do
436  else
437  call pmc_mpi_pack_integer(buffer, position, -1)
438  end if
439  call assert(440557910, &
440  position - prev_position <= pmc_mpi_pack_size_aero_dist(val))
441 #endif
442 
443  end subroutine pmc_mpi_pack_aero_dist
444 
445 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
446 
447  !> Unpacks the given value from the buffer, advancing position.
448  subroutine pmc_mpi_unpack_aero_dist(buffer, position, val)
449 
450  !> Memory buffer.
451  character, intent(inout) :: buffer(:)
452  !> Current buffer position.
453  integer, intent(inout) :: position
454  !> Value to pack.
455  type(aero_dist_t), intent(inout) :: val
456 
457 #ifdef PMC_USE_MPI
458  integer :: prev_position, i, n
459 
460  prev_position = position
461  call pmc_mpi_unpack_integer(buffer, position, n)
462  if (allocated(val%mode)) deallocate(val%mode)
463  if (n >= 0) then
464  allocate(val%mode(n))
465  do i = 1,n
466  call pmc_mpi_unpack_aero_mode(buffer, position, val%mode(i))
467  end do
468  end if
469  call assert(742535268, &
470  position - prev_position <= pmc_mpi_pack_size_aero_dist(val))
471 #endif
472 
473  end subroutine pmc_mpi_unpack_aero_dist
474 
475 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
476 
477 end module pmc_aero_dist
pmc_aero_dist::aero_dist_n_mode
elemental integer function aero_dist_n_mode(aero_dist)
Return the number of modes.
Definition: aero_dist.F90:44
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:236
pmc_mpi
Wrapper functions for MPI.
Definition: mpi.F90:13
pmc_spec_file::spec_file_close
subroutine spec_file_close(file)
Close a spec file.
Definition: spec_file.F90:135
pmc_aero_dist::pmc_mpi_pack_aero_dist
subroutine pmc_mpi_pack_aero_dist(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: aero_dist.F90:419
pmc_util::die_msg
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:134
pmc_constants::dp
integer, parameter dp
Kind of a double precision real number.
Definition: constants.F90:12
pmc_aero_dist::aero_dist_interp_1d
subroutine aero_dist_interp_1d(aero_dist_list, time_list, rate_list, time, aero_dist, rate)
Determine the current aero_dist and rate by interpolating at the current time with the lists of aero_...
Definition: aero_dist.F90:177
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:1102
pmc_spec_file
Reading formatted text input.
Definition: spec_file.F90:43
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_util::assert
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:103
pmc_aero_dist::aero_dist_total_num_conc
real(kind=dp) function aero_dist_total_num_conc(aero_dist)
Returns the total number concentration of a distribution. (#/m^3)
Definition: aero_dist.F90:60
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:378
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:108
pmc_spec_file::spec_file_t
An input file with extra data for printing messages.
Definition: spec_file.F90:59
pmc_aero_dist::aero_dist_number
real(kind=dp) function aero_dist_number(aero_dist, aero_weight)
Returns the total number of particles of a distribution.
Definition: aero_dist.F90:78
pmc_aero_dist::aero_dist_vol_conc
subroutine aero_dist_vol_conc(aero_dist, bin_grid, aero_data, vol_conc)
Return the binned per-species volume concentration for an aero_dist.
Definition: aero_dist.F90:127
pmc_aero_dist::pmc_mpi_unpack_aero_dist
subroutine pmc_mpi_unpack_aero_dist(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: aero_dist.F90:449
pmc_aero_dist
The aero_dist_t structure and associated subroutines.
Definition: aero_dist.F90:18
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_mode::aero_mode_t
An aerosol size distribution mode.
Definition: aero_mode.F90:54
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:1156
pmc_spec_file::spec_file_open
subroutine spec_file_open(filename, file)
Open a spec file for reading.
Definition: spec_file.F90:112
pmc_spec_file::spec_file_read_line_no_eof
subroutine spec_file_read_line_no_eof(file, line)
Read a spec_line from the spec_file. This will always succeed or error out, so should only be called ...
Definition: spec_file.F90:298
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_data::aero_data_t
Aerosol material properties and associated data.
Definition: aero_data.F90:49
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:818
pmc_constants
Physical constants.
Definition: constants.F90:9
pmc_aero_dist::aero_dist_t
A complete aerosol distribution, consisting of several modes.
Definition: aero_dist.F90:33
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:1124
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:495
pmc_aero_mode
The aero_mode_t structure and associated subroutines.
Definition: aero_mode.F90:9
pmc_aero_dist::pmc_mpi_pack_size_aero_dist
integer function pmc_mpi_pack_size_aero_dist(val)
Determines the number of bytes required to pack the given value.
Definition: aero_dist.F90:397
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:561
pmc_aero_data
The aero_data_t structure and associated subroutines.
Definition: aero_data.F90:9
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_aero_dist::aero_dist_num_conc
subroutine aero_dist_num_conc(aero_dist, bin_grid, aero_data, num_conc)
Return the binned number concentration for an aero_dist.
Definition: aero_dist.F90:99
pmc_aero_dist::aero_dist_contains_aero_mode_type
elemental logical function aero_dist_contains_aero_mode_type(aero_dist, aero_mode_type)
Whether any of the modes are of the given type.
Definition: aero_dist.F90:155
pmc_util::find_1d
integer function find_1d(n, x_vals, x)
Find the position of a real number in an arbitrary 1D array.
Definition: util.F90:592
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:413