PartMC
2.2.0
|
The aero_dist_t structure and associated subroutines. More...
Data Types | |
type | aero_dist_t |
A complete aerosol distribution, consisting of several modes. More... | |
Public Member Functions | |
subroutine | aero_dist_allocate (aero_dist) |
Allocates an aero_dist. | |
subroutine | aero_dist_allocate_size (aero_dist, n_mode, n_spec) |
Allocates an aero_dist of the given size. | |
subroutine | aero_dist_deallocate (aero_dist) |
Free all storage. | |
subroutine | aero_dist_copy (aero_dist_from, aero_dist_to) |
Copy an aero_dist. | |
real(kind=dp) function | aero_dist_total_num_conc (aero_dist) |
Returns the total number concentration of a distribution. (#/m^3) | |
real(kind=dp) function | aero_dist_number (aero_dist, aero_weight) |
Returns the total number of particles of a distribution. | |
subroutine | aero_dist_num_conc (aero_dist, bin_grid, aero_data, num_conc) |
Return the binned number concentration for an aero_dist. | |
subroutine | aero_dist_vol_conc (aero_dist, bin_grid, aero_data, vol_conc) |
Return the binned per-species volume 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. | |
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_dists and rates. | |
subroutine | spec_file_read_aero_dist (file, aero_data, aero_dist) |
Read continuous aerosol distribution composed of several modes. | |
subroutine | spec_file_read_aero_dists_times_rates (file, aero_data, times, rates, aero_dists) |
Read an array of aero_dists with associated times and rates from the given file. | |
integer function | pmc_mpi_pack_size_aero_dist (val) |
Determines the number of bytes required to pack the given value. | |
subroutine | pmc_mpi_pack_aero_dist (buffer, position, val) |
Packs the given value into the buffer, advancing position. | |
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.
The initial size distributions are computed as number densities, so they can be used for both sectional and particle-resolved simulations. The routine dist_to_n() converts a number concentration distribution to an actual number of particles ready for a particle-resolved simulation.
Initial distributions should be normalized so that sum(n_den) = 1/log_width
.
Definition at line 18 of file aero_dist.F90.
subroutine pmc_aero_dist::aero_dist_allocate | ( | type(aero_dist_t), intent(out) | aero_dist | ) |
Allocates an aero_dist.
[out] | aero_dist | Aerosol distribution. |
Definition at line 45 of file aero_dist.F90.
subroutine pmc_aero_dist::aero_dist_allocate_size | ( | type(aero_dist_t), intent(out) | aero_dist, |
integer, intent(in) | n_mode, | ||
integer, intent(in) | n_spec | ||
) |
Allocates an aero_dist of the given size.
[out] | aero_dist | Aerosol distribution. |
[in] | n_mode | Number of modes. |
[in] | n_spec | Number of species. |
Definition at line 60 of file aero_dist.F90.
elemental logical function pmc_aero_dist::aero_dist_contains_aero_mode_type | ( | type(aero_dist_t), intent(in) | aero_dist, |
integer, intent(in) | aero_mode_type | ||
) |
Whether any of the modes are of the given type.
[in] | aero_dist | Aerosol distribution. |
[in] | aero_mode_type | Aerosol mode type to test for. |
Definition at line 217 of file aero_dist.F90.
subroutine pmc_aero_dist::aero_dist_copy | ( | type(aero_dist_t), intent(in) | aero_dist_from, |
type(aero_dist_t), intent(inout) | aero_dist_to | ||
) |
Copy an aero_dist.
[in] | aero_dist_from | Aero_dist original. |
[in,out] | aero_dist_to | Aero_dist copy. |
Definition at line 99 of file aero_dist.F90.
subroutine pmc_aero_dist::aero_dist_deallocate | ( | type(aero_dist_t), intent(inout) | aero_dist | ) |
Free all storage.
[in,out] | aero_dist | Aerosol distribution. |
Definition at line 82 of file aero_dist.F90.
subroutine pmc_aero_dist::aero_dist_interp_1d | ( | type(aero_dist_t), dimension(:), intent(in) | aero_dist_list, |
real(kind=dp), dimension(size(aero_dist_list)), intent(in) | time_list, | ||
real(kind=dp), dimension(size(aero_dist_list)), intent(in) | rate_list, | ||
real(kind=dp), intent(in) | time, | ||
type(aero_dist_t), intent(inout) | aero_dist, | ||
real(kind=dp), intent(out) | rate | ||
) |
Determine the current aero_dist and rate by interpolating at the current time with the lists of aero_dists and rates.
[in] | aero_dist_list | Gas states. |
[in] | time_list | Times (s). |
[in] | rate_list | Rates (s^{-1}). |
[in] | time | Current time (s). |
[in,out] | aero_dist | Current gas state. |
[out] | rate | Current rate (s^{-1}). |
Definition at line 234 of file aero_dist.F90.
subroutine pmc_aero_dist::aero_dist_num_conc | ( | type(aero_dist_t), intent(in) | aero_dist, |
type(bin_grid_t), intent(in) | bin_grid, | ||
type(aero_data_t), intent(in) | aero_data, | ||
real(kind=dp), dimension(bin_grid%n_bin), intent(out) | num_conc | ||
) |
Return the binned number concentration for an aero_dist.
[in] | aero_dist | Aero dist for which to compute number concentration. |
[in] | bin_grid | Bin grid. |
[in] | aero_data | Aerosol data. |
[out] | num_conc | Number concentration (#(ln(r))d(ln(r))). |
Definition at line 162 of file aero_dist.F90.
real(kind=dp) function pmc_aero_dist::aero_dist_number | ( | type(aero_dist_t), intent(in) | aero_dist, |
type(aero_weight_t), intent(in) | aero_weight | ||
) |
Returns the total number of particles of a distribution.
[in] | aero_dist | Aerosol distribution. |
[in] | aero_weight | Aerosol weight. |
Definition at line 142 of file aero_dist.F90.
real(kind=dp) function pmc_aero_dist::aero_dist_total_num_conc | ( | type(aero_dist_t), intent(in) | aero_dist | ) |
Returns the total number concentration of a distribution. (#/m^3)
[in] | aero_dist | Aerosol distribution. |
Definition at line 124 of file aero_dist.F90.
subroutine pmc_aero_dist::aero_dist_vol_conc | ( | type(aero_dist_t), intent(in) | aero_dist, |
type(bin_grid_t), intent(in) | bin_grid, | ||
type(aero_data_t), intent(in) | aero_data, | ||
real(kind=dp), dimension(bin_grid%n_bin, aero_data%n_spec), intent(out) | vol_conc | ||
) |
Return the binned per-species volume concentration for an aero_dist.
[in] | aero_dist | Aero dist for which to compute volume concentration. |
[in] | bin_grid | Bin grid. |
[in] | aero_data | Aerosol data. |
[out] | vol_conc | Volume concentration (V(ln(r))d(ln(r))). |
Definition at line 190 of file aero_dist.F90.
subroutine pmc_aero_dist::pmc_mpi_pack_aero_dist | ( | character, dimension(:), intent(inout) | buffer, |
integer, intent(inout) | position, | ||
type(aero_dist_t), intent(in) | val | ||
) |
Packs the given value into the buffer, advancing position.
[in,out] | buffer | Memory buffer. |
[in,out] | position | Current buffer position. |
[in] | val | Value to pack. |
Definition at line 482 of file aero_dist.F90.
integer function pmc_aero_dist::pmc_mpi_pack_size_aero_dist | ( | type(aero_dist_t), intent(in) | val | ) |
Determines the number of bytes required to pack the given value.
[in] | val | Value to pack. |
Definition at line 464 of file aero_dist.F90.
subroutine pmc_aero_dist::pmc_mpi_unpack_aero_dist | ( | character, dimension(:), intent(inout) | buffer, |
integer, intent(inout) | position, | ||
type(aero_dist_t), intent(inout) | val | ||
) |
Unpacks the given value from the buffer, advancing position.
[in,out] | buffer | Memory buffer. |
[in,out] | position | Current buffer position. |
[in,out] | val | Value to pack. |
Definition at line 508 of file aero_dist.F90.
subroutine pmc_aero_dist::spec_file_read_aero_dist | ( | type(spec_file_t), intent(inout) | file, |
type(aero_data_t), intent(inout) | aero_data, | ||
type(aero_dist_t), intent(inout) | aero_dist | ||
) |
Read continuous aerosol distribution composed of several modes.
[in,out] | file | Spec file to read data from. |
[in,out] | aero_data | Aero_data data. |
[in,out] | aero_dist | Aerosol dist. |
Definition at line 274 of file aero_dist.F90.
subroutine pmc_aero_dist::spec_file_read_aero_dists_times_rates | ( | type(spec_file_t), intent(inout) | file, |
type(aero_data_t), intent(inout) | aero_data, | ||
real(kind=dp), dimension(:), pointer | times, | ||
real(kind=dp), dimension(:), pointer | rates, | ||
type(aero_dist_t), dimension(:), pointer | aero_dists | ||
) |
Read an array of aero_dists with associated times and rates from the given file.
[in,out] | file | Spec file to read data from. |
[in,out] | aero_data | Aero data. |
times | Times (s). | |
rates | Rates (s^{-1}). | |
aero_dists | Aero dists. |
Definition at line 331 of file aero_dist.F90.