PartMC
2.2.1
|
The aero_weight_t structure and associated subroutines. More...
Data Types | |
type | aero_weight_t |
An aerosol size distribution weighting function. More... | |
Public Member Functions | |
elemental subroutine | aero_weight_allocate (aero_weight) |
Allocates an aero_weight. | |
elemental subroutine | aero_weight_allocate_size (aero_weight) |
Allocates an aero_weight of the given size. | |
elemental subroutine | aero_weight_deallocate (aero_weight) |
Free all storage. | |
elemental subroutine | aero_weight_zero (aero_weight) |
Zeros the contents of the aero_weight . | |
elemental subroutine | aero_weight_zero_comp_vol (aero_weight) |
Zeros the aero_weight computational volume. | |
elemental subroutine | aero_weight_copy (aero_weight_from, aero_weight_to) |
Copy an aero_weight. | |
subroutine | aero_weight_array_copy (aero_weight_array_from, aero_weight_array_to) |
Copy an aero_weight array. | |
elemental subroutine | aero_weight_scale_comp_vol (aero_weight, fraction) |
Scale the computational volume by the given fraction, so new_comp_vol = old_comp_vol * fraction . | |
elemental subroutine | aero_weight_add_comp_vol (aero_weight, aero_weight_delta) |
Add the computational volume of aero_weight_delta to aero_weight . | |
elemental subroutine | aero_weight_transfer_comp_vol (aero_weight_from, aero_weight_to, sample_prop) |
Transfer the computational volume from aero_weight_from to aero_weight_to , weighted by sample_prop . | |
real(kind=dp) function | aero_weight_num_conc_at_radius (aero_weight, radius) |
Compute the number concentration at a given radius (m^{-3}). | |
real(kind=dp) function | aero_weight_num_conc (aero_weight, aero_particle) |
Compute the number concentration for a particle (m^{-3}). | |
real(kind=dp) function | aero_weight_array_single_num_conc (aero_weight_array, aero_particle) |
Compute the number concentration for a particle (m^{-3}). | |
real(kind=dp) function | aero_weight_array_num_conc_at_radius (aero_weight_array, radius) |
Compute the total number concentration at a given radius (m^3). | |
real(kind=dp) function | aero_weight_array_num_conc (aero_weight_array, aero_particle) |
Compute the number concentration for a particle (m^{-3}). | |
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(). | |
logical function | aero_weight_array_check_flat (aero_weight_array) |
Check whether a given aero_weight array is flat in total. | |
subroutine | aero_weight_array_check_monotonicity (aero_weight_array, monotone_increasing, monotone_decreasing) |
Determine whether all weight functions in an array are monotone increasing, monotone decreasing, or neither. | |
subroutine | aero_weight_array_minmax_num_conc (aero_weight_array, radius_1, radius_2, num_conc_min, num_conc_max) |
Compute the maximum and minimum number concentrations between the given radii. | |
integer function | aero_weight_array_rand_group (aero_weight_array, radius) |
Choose a random group at the given radius, with probability proportional to group volume at that radius. | |
subroutine | spec_file_read_aero_weight (file, aero_weight) |
Read an aero_weight from a spec file. | |
integer function | pmc_mpi_pack_size_aero_weight (val) |
Determines the number of bytes required to pack the given value. | |
subroutine | pmc_mpi_pack_aero_weight (buffer, position, val) |
Packs the given value into the buffer, advancing position. | |
subroutine | pmc_mpi_unpack_aero_weight (buffer, position, val) |
Unpacks the given value from the buffer, advancing position. | |
subroutine | aero_weight_netcdf_dim_aero_weight (aero_weight_array, ncid, dimid_aero_weight) |
Write the aero_weight dimension to the given NetCDF file if it is not already present and in any case return the associated dimid. | |
subroutine | aero_weight_array_output_netcdf (aero_weight_array, ncid) |
Write full state. | |
subroutine | aero_weight_array_input_netcdf (aero_weight_array, ncid) |
Read full state. | |
Public Attributes | |
integer, parameter | AERO_WEIGHT_TYPE_INVALID = 0 |
Type code for an undefined or invalid weighting. | |
integer, parameter | AERO_WEIGHT_TYPE_NONE = 1 |
Type code for no (or flat) weighting. | |
integer, parameter | AERO_WEIGHT_TYPE_POWER = 2 |
Type code for power function weighting. | |
integer, parameter | AERO_WEIGHT_TYPE_MFA = 3 |
Type code for MFA weighting. |
The aero_weight_t structure and associated subroutines.
Definition at line 9 of file aero_weight.F90.
elemental subroutine pmc_aero_weight::aero_weight_add_comp_vol | ( | type(aero_weight_t), intent(inout) | aero_weight, |
type(aero_weight_t), intent(in) | aero_weight_delta | ||
) |
Add the computational volume of aero_weight_delta
to aero_weight
.
[in,out] | aero_weight | Aerosol weight to add volume to. |
[in] | aero_weight_delta | Aerosol weight to add volume from. |
Definition at line 162 of file aero_weight.F90.
elemental subroutine pmc_aero_weight::aero_weight_allocate | ( | type(aero_weight_t), intent(out) | aero_weight | ) |
Allocates an aero_weight.
[out] | aero_weight | Aerosol weight. |
Definition at line 48 of file aero_weight.F90.
elemental subroutine pmc_aero_weight::aero_weight_allocate_size | ( | type(aero_weight_t), intent(out) | aero_weight | ) |
Allocates an aero_weight of the given size.
[out] | aero_weight | Aerosol weight. |
Definition at line 60 of file aero_weight.F90.
logical function pmc_aero_weight::aero_weight_array_check_flat | ( | type(aero_weight_t), dimension(:), intent(in) | aero_weight_array | ) |
Check whether a given aero_weight array is flat in total.
[in] | aero_weight_array | Aerosol weight array. |
Definition at line 321 of file aero_weight.F90.
subroutine pmc_aero_weight::aero_weight_array_check_monotonicity | ( | type(aero_weight_t), dimension(:), intent(in) | aero_weight_array, |
logical, intent(out) | monotone_increasing, | ||
logical, intent(out) | monotone_decreasing | ||
) |
Determine whether all weight functions in an array are monotone increasing, monotone decreasing, or neither.
[in] | aero_weight_array | Aerosol weight array. |
[out] | monotone_increasing | Whether all weights are monotone increasing. |
[out] | monotone_decreasing | Whether all weights are monotone decreasing. |
Definition at line 355 of file aero_weight.F90.
subroutine pmc_aero_weight::aero_weight_array_copy | ( | type(aero_weight_t), dimension(:), intent(in) | aero_weight_array_from, |
type(aero_weight_t), dimension(:), intent(inout), allocatable | aero_weight_array_to | ||
) |
Copy an aero_weight array.
[in] | aero_weight_array_from | Aerosol weight array original. |
[in,out] | aero_weight_array_to | Aerosol weight array copy. |
Definition at line 126 of file aero_weight.F90.
subroutine pmc_aero_weight::aero_weight_array_input_netcdf | ( | type(aero_weight_t), dimension(:), intent(inout), allocatable | aero_weight_array, |
integer, intent(in) | ncid | ||
) |
Read full state.
[in,out] | aero_weight_array | Environment state to read. |
[in] | ncid | NetCDF file ID, in data mode. |
Definition at line 675 of file aero_weight.F90.
subroutine pmc_aero_weight::aero_weight_array_minmax_num_conc | ( | type(aero_weight_t), dimension(:), intent(in) | aero_weight_array, |
real(kind=dp), intent(in) | radius_1, | ||
real(kind=dp), intent(in) | radius_2, | ||
real(kind=dp), intent(out) | num_conc_min, | ||
real(kind=dp), intent(out) | num_conc_max | ||
) |
Compute the maximum and minimum number concentrations between the given radii.
[in] | aero_weight_array | Aerosol weight. |
[in] | radius_1 | First radius. |
[in] | radius_2 | Second radius. |
[out] | num_conc_min | Minimum number concentration. |
[out] | num_conc_max | Maximum number concentration. |
Definition at line 394 of file aero_weight.F90.
real(kind=dp) function pmc_aero_weight::aero_weight_array_num_conc | ( | type(aero_weight_t), dimension(:), intent(in) | aero_weight_array, |
type(aero_particle_t), intent(in) | aero_particle | ||
) |
Compute the number concentration for a particle (m^{-3}).
[in] | aero_weight_array | Aerosol weight array. |
[in] | aero_particle | Aerosol particle to compute number concentration for. |
Definition at line 279 of file aero_weight.F90.
real(kind=dp) function pmc_aero_weight::aero_weight_array_num_conc_at_radius | ( | type(aero_weight_t), dimension(:), intent(in) | aero_weight_array, |
real(kind=dp), intent(in) | radius | ||
) |
Compute the total number concentration at a given radius (m^3).
[in] | aero_weight_array | Aerosol weight array. |
[in] | radius | Radius to compute number concentration at (m). |
Definition at line 255 of file aero_weight.F90.
subroutine pmc_aero_weight::aero_weight_array_output_netcdf | ( | type(aero_weight_t), dimension(:), intent(in) | aero_weight_array, |
integer, intent(in) | ncid | ||
) |
Write full state.
[in] | aero_weight_array | Aero weight to write. |
[in] | ncid | NetCDF file ID, in data mode. |
Definition at line 626 of file aero_weight.F90.
integer function pmc_aero_weight::aero_weight_array_rand_group | ( | type(aero_weight_t), dimension(:), intent(in) | aero_weight_array, |
real(kind=dp), intent(in) | radius | ||
) |
Choose a random group at the given radius, with probability proportional to group volume at that radius.
[in] | aero_weight_array | Aerosol weight array. |
[in] | radius | Radius to sample group at (m). |
Definition at line 428 of file aero_weight.F90.
real(kind=dp) function pmc_aero_weight::aero_weight_array_single_num_conc | ( | type(aero_weight_t), dimension(:), intent(in) | aero_weight_array, |
type(aero_particle_t), intent(in) | aero_particle | ||
) |
Compute the number concentration for a particle (m^{-3}).
[in] | aero_weight_array | Aerosol weight array. |
[in] | aero_particle | Aerosol particle to compute number concentration for. |
Definition at line 239 of file aero_weight.F90.
elemental subroutine pmc_aero_weight::aero_weight_copy | ( | type(aero_weight_t), intent(in) | aero_weight_from, |
type(aero_weight_t), intent(inout) | aero_weight_to | ||
) |
Copy an aero_weight.
[in] | aero_weight_from | Aerosol weight original. |
[in,out] | aero_weight_to | Aerosol weight copy. |
Definition at line 109 of file aero_weight.F90.
elemental subroutine pmc_aero_weight::aero_weight_deallocate | ( | type(aero_weight_t), intent(in) | aero_weight | ) |
Free all storage.
[in] | aero_weight | Aerosol weight. |
Definition at line 72 of file aero_weight.F90.
subroutine pmc_aero_weight::aero_weight_netcdf_dim_aero_weight | ( | type(aero_weight_t), dimension(:), intent(in) | aero_weight_array, |
integer, intent(in) | ncid, | ||
integer, intent(out) | dimid_aero_weight | ||
) |
Write the aero_weight
dimension to the given NetCDF file if it is not already present and in any case return the associated dimid.
[in] | aero_weight_array | Aero_weight structure array. |
[in] | ncid | NetCDF file ID, in data mode. |
[out] | dimid_aero_weight | Dimid of the species dimension. |
Definition at line 582 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 | ||
) |
Compute the number concentration for a particle (m^{-3}).
[in] | aero_weight | Aerosol weight. |
[in] | aero_particle | Aerosol particle to compute number concentration for. |
Definition at line 223 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 198 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 296 of file aero_weight.F90.
elemental subroutine pmc_aero_weight::aero_weight_scale_comp_vol | ( | type(aero_weight_t), intent(inout) | aero_weight, |
real(kind=dp), intent(in) | fraction | ||
) |
Scale the computational volume by the given fraction, so new_comp_vol = old_comp_vol * fraction
.
[in,out] | aero_weight | Aerosol weight to halve. |
[in] | fraction | Fraction to scale computational volume by. |
Definition at line 148 of file aero_weight.F90.
elemental subroutine pmc_aero_weight::aero_weight_transfer_comp_vol | ( | 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 | ||
) |
Transfer the computational volume from aero_weight_from
to aero_weight_to
, weighted by sample_prop
.
[in,out] | aero_weight_from | Aerosol weight to take volume from. |
[in,out] | aero_weight_to | Aerosol weight to add volume to. |
[in] | sample_prop | Proportion of from volume to transfer. |
Definition at line 177 of file aero_weight.F90.
elemental subroutine pmc_aero_weight::aero_weight_zero | ( | type(aero_weight_t), intent(inout) | aero_weight | ) |
Zeros the contents of the aero_weight
.
[in,out] | aero_weight | Aerosol weight. |
Definition at line 82 of file aero_weight.F90.
elemental subroutine pmc_aero_weight::aero_weight_zero_comp_vol | ( | type(aero_weight_t), intent(inout) | aero_weight | ) |
Zeros the aero_weight
computational volume.
[in,out] | aero_weight | Aerosol weight. |
Definition at line 97 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 528 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 512 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 554 of file aero_weight.F90.
subroutine pmc_aero_weight::spec_file_read_aero_weight | ( | type(spec_file_t), intent(inout) | file, |
type(aero_weight_t), intent(inout) | aero_weight | ||
) |
Read an aero_weight from a spec file.
[in,out] | file | Spec file. |
[in,out] | aero_weight | Aerosol weight. |
Definition at line 450 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 23 of file aero_weight.F90.
integer, parameter pmc_aero_weight::AERO_WEIGHT_TYPE_MFA = 3 |
Type code for MFA weighting.
Definition at line 29 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 25 of file aero_weight.F90.
integer, parameter pmc_aero_weight::AERO_WEIGHT_TYPE_POWER = 2 |
Type code for power function weighting.
Definition at line 27 of file aero_weight.F90.