PartMC
2.6.1
|
The aero_weight_t structure and associated subroutines. More...
Data Types | |
type | aero_weight_t |
An aerosol size distribution weighting function. More... | |
Functions/Subroutines | |
elemental subroutine | aero_weight_normalize (aero_weight) |
Sets the aero_weight to a non-zero normalized value. More... | |
elemental subroutine | aero_weight_scale (aero_weight, factor) |
Scale the weight by the given fraction, so new_weight = old_weight * factor . More... | |
elemental subroutine | aero_weight_combine (aero_weight, aero_weight_delta) |
Combine aero_weight_delta into aero_weight with a harmonic mean. More... | |
elemental subroutine | aero_weight_shift (aero_weight_from, aero_weight_to, sample_prop, overwrite_to) |
Adjust source and destination weights to reflect moving sample_prop proportion of particles from aero_weight_from to aero_weight_to . More... | |
real(kind=dp) function | aero_weight_num_conc_at_radius (aero_weight, radius) |
Compute the number concentration at a given radius (m^{-3}). More... | |
real(kind=dp) function | aero_weight_num_conc (aero_weight, aero_particle, aero_data) |
Compute the number concentration for a particle (m^{-3}). More... | |
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(). More... | |
subroutine | aero_weight_check_valid_exponent (aero_weight) |
Ensures that a weight function exponent is valid. More... | |
subroutine | aero_weight_check_monotonicity (aero_weight, monotone_increasing, monotone_decreasing) |
Determine whether a weight function is monotone increasing, monotone decreasing, or neither. More... | |
integer function | pmc_mpi_pack_size_aero_weight (val) |
Determines the number of bytes required to pack the given value. More... | |
subroutine | pmc_mpi_pack_aero_weight (buffer, position, val) |
Packs the given value into the buffer, advancing position. More... | |
subroutine | pmc_mpi_unpack_aero_weight (buffer, position, val) |
Unpacks the given value from the buffer, advancing position. More... | |
Variables | |
integer, parameter | aero_weight_type_invalid = 0 |
Type code for an undefined or invalid weighting. More... | |
integer, parameter | aero_weight_type_none = 1 |
Type code for no (or flat) weighting. More... | |
integer, parameter | aero_weight_type_power = 2 |
Type code for power function weighting. More... | |
integer, parameter | aero_weight_type_mfa = 3 |
Type code for MFA weighting. More... | |
The aero_weight_t structure and associated subroutines.
subroutine pmc_aero_weight::aero_weight_check_monotonicity | ( | type(aero_weight_t), intent(in) | aero_weight, |
logical, intent(out) | monotone_increasing, | ||
logical, intent(out) | monotone_decreasing | ||
) |
Determine whether a weight function is monotone increasing, monotone decreasing, or neither.
[in] | aero_weight | Aerosol weight array. |
[out] | monotone_increasing | Whether all weights are monotone increasing. |
[out] | monotone_decreasing | Whether all weights are monotone decreasing. |
Definition at line 218 of file aero_weight.F90.
subroutine pmc_aero_weight::aero_weight_check_valid_exponent | ( | type(aero_weight_t), intent(in) | aero_weight | ) |
Ensures that a weight function exponent is valid.
[in] | aero_weight | Aerosol weight array. |
Definition at line 196 of file aero_weight.F90.
elemental subroutine pmc_aero_weight::aero_weight_combine | ( | type(aero_weight_t), intent(inout) | aero_weight, |
type(aero_weight_t), intent(in) | aero_weight_delta | ||
) |
Combine aero_weight_delta
into aero_weight
with a harmonic mean.
[in,out] | aero_weight | Aerosol weight to add into. |
[in] | aero_weight_delta | Aerosol weight to add from. |
Definition at line 74 of file aero_weight.F90.
elemental subroutine pmc_aero_weight::aero_weight_normalize | ( | type(aero_weight_t), intent(inout) | aero_weight | ) |
Sets the aero_weight
to a non-zero normalized value.
[in,out] | aero_weight | Aerosol weight. |
Definition at line 47 of file aero_weight.F90.
real(kind=dp) function pmc_aero_weight::aero_weight_num_conc | ( | type(aero_weight_t), intent(in) | aero_weight, |
type(aero_particle_t), intent(in) | aero_particle, | ||
type(aero_data_t), intent(in) | aero_data | ||
) |
Compute the number concentration for a particle (m^{-3}).
[in] | aero_weight | Aerosol weight. |
[in] | aero_particle | Aerosol particle to compute number concentration for. |
[in] | aero_data | Aerosol data. |
Definition at line 151 of file aero_weight.F90.
real(kind=dp) function pmc_aero_weight::aero_weight_num_conc_at_radius | ( | type(aero_weight_t), intent(in) | aero_weight, |
real(kind=dp), intent(in) | radius | ||
) |
Compute the number concentration at a given radius (m^{-3}).
[in] | aero_weight | Aerosol weight. |
[in] | radius | Radius to compute number concentration at (m). |
Definition at line 127 of file aero_weight.F90.
real(kind=dp) function pmc_aero_weight::aero_weight_radius_at_num_conc | ( | type(aero_weight_t), intent(in) | aero_weight, |
real(kind=dp), intent(in) | num_conc | ||
) |
Compute the radius at a given number concentration (m). Inverse of aero_weight_num_conc_at_radius().
[in] | aero_weight | Aerosol weight. |
[in] | num_conc | Number concentration to compute the radius at (m^{-3}). |
Definition at line 170 of file aero_weight.F90.
elemental subroutine pmc_aero_weight::aero_weight_scale | ( | type(aero_weight_t), intent(inout) | aero_weight, |
real(kind=dp), intent(in) | factor | ||
) |
Scale the weight by the given fraction, so new_weight = old_weight * factor
.
[in,out] | aero_weight | Aerosol weight to scale. |
[in] | factor | Factor to scale by. |
Definition at line 60 of file aero_weight.F90.
elemental subroutine pmc_aero_weight::aero_weight_shift | ( | type(aero_weight_t), intent(inout) | aero_weight_from, |
type(aero_weight_t), intent(inout) | aero_weight_to, | ||
real(kind=dp), intent(in) | sample_prop, | ||
logical, intent(in), optional | overwrite_to | ||
) |
Adjust source and destination weights to reflect moving sample_prop
proportion of particles from aero_weight_from
to aero_weight_to
.
[in,out] | aero_weight_from | Aerosol weight to shift from. |
[in,out] | aero_weight_to | Aerosol weight to shift to. |
[in] | sample_prop | Proportion of particles being transfered. |
[in] | overwrite_to | Whether to overwrite the destination weight (default: no). |
Definition at line 91 of file aero_weight.F90.
subroutine pmc_aero_weight::pmc_mpi_pack_aero_weight | ( | character, dimension(:), intent(inout) | buffer, |
integer, intent(inout) | position, | ||
type(aero_weight_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 327 of file aero_weight.F90.
integer function pmc_aero_weight::pmc_mpi_pack_size_aero_weight | ( | type(aero_weight_t), intent(in) | val | ) |
Determines the number of bytes required to pack the given value.
[in] | val | Value to pack. |
Definition at line 312 of file aero_weight.F90.
subroutine pmc_aero_weight::pmc_mpi_unpack_aero_weight | ( | character, dimension(:), intent(inout) | buffer, |
integer, intent(inout) | position, | ||
type(aero_weight_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 352 of file aero_weight.F90.
integer, parameter pmc_aero_weight::aero_weight_type_invalid = 0 |
Type code for an undefined or invalid weighting.
Definition at line 24 of file aero_weight.F90.
integer, parameter pmc_aero_weight::aero_weight_type_mfa = 3 |
Type code for MFA weighting.
Definition at line 30 of file aero_weight.F90.
integer, parameter pmc_aero_weight::aero_weight_type_none = 1 |
Type code for no (or flat) weighting.
Definition at line 26 of file aero_weight.F90.
integer, parameter pmc_aero_weight::aero_weight_type_power = 2 |
Type code for power function weighting.
Definition at line 28 of file aero_weight.F90.