Go to the documentation of this file.
48 if (
allocated(aero_dist%mode))
then
82 type(aero_weight_t),
intent(in) :: aero_weight
110 real(kind=
dp) :: mode_num_conc(
size(num_conc, 1))
115 aero_data, mode_num_conc)
116 num_conc = num_conc + mode_num_conc
139 real(kind=
dp) :: mode_vol_conc(
size(vol_conc, 1),
size(vol_conc, 2))
144 aero_data, mode_vol_conc)
145 vol_conc = vol_conc + mode_vol_conc
159 integer,
intent(in) :: aero_mode_type
166 .or. (aero_dist%mode(i_mode)%type == aero_mode_type)
176 rate_list, time, aero_dist, rate)
181 real(kind=
dp),
intent(in) :: time_list(
size(aero_dist_list))
183 real(kind=
dp),
intent(in) :: rate_list(
size(aero_dist_list))
185 real(kind=
dp),
intent(in) :: time
189 real(kind=
dp),
intent(out) :: rate
191 integer :: n, p, n_bin, n_spec, i, i_new
192 real(kind=
dp) :: y, alpha
194 n =
size(aero_dist_list)
195 p =
find_1d(n, time_list, time)
198 aero_dist = aero_dist_list(1)
202 aero_dist = aero_dist_list(n)
206 aero_dist = aero_dist_list(p)
215 subroutine spec_file_read_aero_dist(file, aero_data, &
216 read_aero_weight_classes, aero_dist)
223 logical,
intent(in) :: read_aero_weight_classes
245 allocate(aero_dist%mode(0))
246 call spec_file_read_aero_mode(file, aero_data, read_aero_weight_classes, &
249 n =
size(aero_dist%mode)
250 allocate(new_modes(n + 1))
252 new_modes(i) = aero_dist%mode(i)
254 new_modes(n + 1) = aero_mode
260 deallocate(aero_dist%mode)
261 allocate(aero_dist%mode(n+1))
262 aero_dist%mode = new_modes
263 deallocate(new_modes)
268 call spec_file_read_aero_mode(file, aero_data, &
269 read_aero_weight_classes, aero_mode, eof)
272 end subroutine spec_file_read_aero_dist
278 subroutine spec_file_read_aero_dists_times_rates(file, aero_data, &
279 read_aero_weight_classes, times, rates, aero_dists)
286 logical,
intent(in) :: read_aero_weight_classes
288 real(kind=
dp),
allocatable :: times(:)
290 real(kind=
dp),
allocatable :: rates(:)
294 type(spec_line_t) :: aero_dist_line
296 integer :: n_time, i_time
297 character(len=SPEC_LINE_MAX_VAR_LEN),
allocatable :: names(:)
298 real(kind=
dp),
allocatable ::
data(:,:)
368 if (
size(names) /= 2)
then
370 trim(file%name) //
' must contain exactly two data lines')
372 if (trim(names(1)) /=
'time')
then
373 call die_msg(570205795,
'row 1 in ' // trim(file%name) &
374 //
' must start with: time not: ' // trim(names(1)))
376 if (trim(names(2)) /=
'rate')
then
377 call die_msg(221270915,
'row 2 in ' // trim(file%name) &
378 //
' must start with: rate not: ' // trim(names(1)))
380 n_time =
size(
data, 2)
382 call die_msg(457229710,
'each line in ' // trim(file%name) &
383 //
' must contain at least one data value')
389 if (
allocated(aero_dists))
deallocate(aero_dists)
390 allocate(aero_dists(n_time))
393 call spec_file_read_aero_dist(aero_dist_file, aero_data, &
394 read_aero_weight_classes, aero_dists(i_time))
398 end subroutine spec_file_read_aero_dists_times_rates
408 integer :: i, total_size
410 if (
allocated(val%mode))
then
412 do i = 1,
size(val%mode)
428 character,
intent(inout) :: buffer(:)
430 integer,
intent(inout) :: position
435 integer :: prev_position, i
437 prev_position = position
438 if (
allocated(val%mode))
then
440 do i = 1,
size(val%mode)
458 character,
intent(inout) :: buffer(:)
460 integer,
intent(inout) :: position
465 integer :: prev_position, i, n
467 prev_position = position
469 if (
allocated(val%mode))
deallocate(val%mode)
471 allocate(val%mode(n))
elemental integer function aero_dist_n_mode(aero_dist)
Return the number of modes.
elemental integer function aero_data_n_spec(aero_data)
Return the number of aerosol species, or -1 if uninitialized.
Wrapper functions for MPI.
subroutine spec_file_close(file)
Close a spec file.
subroutine pmc_mpi_pack_aero_dist(buffer, position, val)
Packs the given value into the buffer, advancing position.
subroutine die_msg(code, error_msg)
Error immediately.
integer, parameter dp
Kind of a double precision real number.
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_...
integer function pmc_mpi_pack_size_aero_mode(val)
Determines the number of bytes required to pack the given value.
Reading formatted text input.
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.
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
real(kind=dp) function aero_dist_total_num_conc(aero_dist)
Returns the total number concentration of a distribution. (#/m^3)
subroutine aero_mode_num_conc(aero_mode, bin_grid, aero_data, num_conc)
Return the binned number concentration for an aero_mode.
real(kind=dp) function aero_mode_total_num_conc(aero_mode)
Returns the total number concentration of a mode. (#/m^3)
An input file with extra data for printing messages.
real(kind=dp) function aero_dist_number(aero_dist, aero_weight)
Returns the total number of particles of a distribution.
subroutine aero_dist_vol_conc(aero_dist, bin_grid, aero_data, vol_conc)
Return the binned per-species volume concentration for an aero_dist.
subroutine pmc_mpi_unpack_aero_dist(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
The aero_dist_t structure and associated subroutines.
subroutine spec_file_check_line_length(file, line, length)
Check that the length of the line data is as given.
An aerosol size distribution mode.
subroutine pmc_mpi_unpack_aero_mode(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
subroutine spec_file_open(filename, file)
Open a spec file for reading.
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 ...
Random number generators.
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.
Aerosol material properties and associated data.
subroutine pmc_mpi_unpack_integer(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
A complete aerosol distribution, consisting of several modes.
subroutine pmc_mpi_pack_aero_mode(buffer, position, val)
Packs the given value into the buffer, advancing position.
Common utility subroutines.
real(kind=dp) function aero_mode_number(aero_mode, aero_weight)
Return the total number of computational particles for an aero_mode.
The aero_mode_t structure and associated subroutines.
integer function pmc_mpi_pack_size_aero_dist(val)
Determines the number of bytes required to pack the given value.
integer function pmc_mpi_pack_size_integer(val)
Determines the number of bytes required to pack the given value.
The bin_grid_t structure and associated subroutines.
subroutine pmc_mpi_pack_integer(buffer, position, val)
Packs the given value into the buffer, advancing position.
The aero_data_t structure and associated subroutines.
subroutine spec_file_check_line_name(file, line, name)
Check that the name of the line data is as given.
1D grid, either logarithmic or linear.
subroutine aero_dist_num_conc(aero_dist, bin_grid, aero_data, num_conc)
Return the binned number concentration for an aero_dist.
elemental logical function aero_dist_contains_aero_mode_type(aero_dist, aero_mode_type)
Whether any of the modes are of the given type.
integer function find_1d(n, x_vals, x)
Find the position of a real number in an arbitrary 1D array.
subroutine aero_mode_vol_conc(aero_mode, bin_grid, aero_data, vol_conc)
Return the binned per-species volume concentration for an aero_mode.