|
PartMC
2.8.0
|
Functions/Subroutines | |
| subroutine | ice_nucleation_immersion_freezing (aero_state, aero_data, env_state, del_t, immersion_freezing_scheme_type, freezing_rate, do_freezing_naive, INAS_a, INAS_b) |
| Main subroutine for immersion freezing simulation. More... | |
| subroutine | ice_nucleation_singular_initialize (aero_state, aero_data, INAS_a, INAS_b) |
| Initialization for the singular scheme, sampling the freezing temperature for each particle. More... | |
| subroutine | ice_nucleation_immersion_freezing_singular (aero_state, aero_data, env_state) |
| Simulation for singular scheme, deciding whether to freeze for each particle. Run in each time step. More... | |
| subroutine | ice_nucleation_immersion_freezing_time_dependent (aero_state, aero_data, env_state, del_t, immersion_freezing_scheme_type, freezing_rate) |
| Simulation for time-dependent scheme (e.g., ABIFM, constant rate), deciding whether to freeze for each particle. Run in each time step. This subroutine applies the binned-tau leaping algorithm for speeding up. More... | |
| subroutine | ice_nucleation_immersion_freezing_time_dependent_naive (aero_state, aero_data, env_state, del_t, immersion_freezing_scheme_type, freezing_rate) |
| Simulation for time-dependent scheme (e.g., ABIFM, constant rate), deciding whether to freeze for each particle. Run in each time step. This subroutine applies the naive algorithm that checks each particle. More... | |
| subroutine | ice_nucleation_melting (aero_state, aero_data, env_state) |
| Simulates melting: if the environmental temperature is above the freezing temperature of water, all particles are set to be unfrozen. More... | |
| real(kind=dp) function | abifm_pfrz_particle (aero_particle, aero_data, a_w_ice, del_t) |
| Calculating the freezing probability for the particle (i_part) using ABIFM method (Knopf et al.,2013) More... | |
| subroutine | abifm_max_spec (aero_data, a_w_ice, i_spec_max, j_het_max) |
| Finding the maximum heterogeneous ice nucleation rate coefficient. More... | |
| real(kind=dp) function | abifm_pfrz_max (diameter_max, aero_data, j_het_max, del_t) |
| Calculating the maximum freezing probability for particles in one bin using ABIFM method (Knopf et al.,2013). Only used by the binned-tau leaping algorithm. More... | |
Variables | |
| integer, parameter | immersion_freezing_scheme_invalid = 0 |
| Type code for an undefined or invalid immersion freezing scheme. More... | |
| integer, parameter | immersion_freezing_scheme_const = 1 |
| Type code for constant ice nucleation rate (J_het) immersion freezing scheme. More... | |
| integer, parameter | immersion_freezing_scheme_singular = 2 |
| Type code for the singular (INAS) immersion freezing scheme. More... | |
| integer, parameter | immersion_freezing_scheme_abifm = 3 |
| Type code for the ABIFM immersion freezing scheme. More... | |
| subroutine pmc_ice_nucleation::abifm_max_spec | ( | type(aero_data_t), intent(in) | aero_data, |
| real(kind=dp), intent(in) | a_w_ice, | ||
| integer, intent(out) | i_spec_max, | ||
| real(kind=dp), intent(out) | j_het_max | ||
| ) |
Finding the maximum heterogeneous ice nucleation rate coefficient.
| [in] | aero_data | Aerosol data. |
| [in] | a_w_ice | The water activity w.r.t. ice. |
| [out] | i_spec_max | The index of the maximum J_het species. |
| [out] | j_het_max | The maximum value of J_het among all species. |
Definition at line 426 of file ice_nucleation.F90.
| real(kind=dp) function pmc_ice_nucleation::abifm_pfrz_max | ( | real(kind=dp), intent(in) | diameter_max, |
| type(aero_data_t), intent(in) | aero_data, | ||
| real(kind=dp), intent(in) | j_het_max, | ||
| real(kind=dp), intent(in) | del_t | ||
| ) |
Calculating the maximum freezing probability for particles in one bin using ABIFM method (Knopf et al.,2013). Only used by the binned-tau leaping algorithm.
| [in] | aero_data | Aerosol data. |
| [in] | diameter_max | Maximum diameter. |
| [in] | del_t | Time interval. |
| [in] | j_het_max | Maximum J_het among all species. |
Definition at line 460 of file ice_nucleation.F90.
| real(kind=dp) function pmc_ice_nucleation::abifm_pfrz_particle | ( | type(aero_particle_t), intent(in) | aero_particle, |
| type(aero_data_t), intent(in) | aero_data, | ||
| real(kind=dp), intent(in) | a_w_ice, | ||
| real(kind=dp), intent(in) | del_t | ||
| ) |
Calculating the freezing probability for the particle (i_part) using ABIFM method (Knopf et al.,2013)
| [in] | aero_particle | Aerosol particle. |
| [in] | aero_data | Aerosol data. |
| [in] | a_w_ice | The water activity w.r.t. ice. |
| [in] | del_t | Time interval. |
Definition at line 383 of file ice_nucleation.F90.
| subroutine pmc_ice_nucleation::ice_nucleation_immersion_freezing | ( | type(aero_state_t), intent(inout) | aero_state, |
| type(aero_data_t), intent(in) | aero_data, | ||
| type(env_state_t), intent(inout) | env_state, | ||
| real(kind=dp), intent(in) | del_t, | ||
| integer, intent(in) | immersion_freezing_scheme_type, | ||
| real(kind=dp), intent(in) | freezing_rate, | ||
| logical, intent(in) | do_freezing_naive, | ||
| real(kind=dp), intent(in) | INAS_a, | ||
| real(kind=dp), intent(in) | INAS_b | ||
| ) |
Main subroutine for immersion freezing simulation.
| [in,out] | aero_state | Aerosol state. |
| [in] | aero_data | Aerosol data. |
| [in,out] | env_state | Environment state. |
| [in] | del_t | Total time to integrate. |
| [in] | immersion_freezing_scheme_type | Immersion freezing scheme type. |
| [in] | freezing_rate | Freezing rate (only used for the constant rate scheme). |
| [in] | do_freezing_naive | Whether to use the naive algorithm for time-dependent scheme. (If false, use the binned tau-leaping algorithm.) |
| [in] | inas_a | Slope parameter for the INAS parameterization (singular scheme only). |
| [in] | inas_b | Intercept parameter for the INAS parameterization (singular scheme only). |
Definition at line 34 of file ice_nucleation.F90.
| subroutine pmc_ice_nucleation::ice_nucleation_immersion_freezing_singular | ( | type(aero_state_t), intent(inout) | aero_state, |
| type(aero_data_t), intent(in) | aero_data, | ||
| type(env_state_t), intent(inout) | env_state | ||
| ) |
Simulation for singular scheme, deciding whether to freeze for each particle. Run in each time step.
| [in,out] | aero_state | Aerosol state. |
| [in] | aero_data | Aerosol data. |
| [in,out] | env_state | Environment state. |
Definition at line 126 of file ice_nucleation.F90.
| subroutine pmc_ice_nucleation::ice_nucleation_immersion_freezing_time_dependent | ( | type(aero_state_t), intent(inout) | aero_state, |
| type(aero_data_t), intent(in) | aero_data, | ||
| type(env_state_t), intent(inout) | env_state, | ||
| real(kind=dp), intent(in) | del_t, | ||
| integer, intent(in) | immersion_freezing_scheme_type, | ||
| real(kind=dp), intent(in) | freezing_rate | ||
| ) |
Simulation for time-dependent scheme (e.g., ABIFM, constant rate), deciding whether to freeze for each particle. Run in each time step. This subroutine applies the binned-tau leaping algorithm for speeding up.
| [in,out] | aero_state | Aerosol state. |
| [in] | aero_data | Aerosol data. |
| [in,out] | env_state | Environment state. |
| [in] | del_t | Total time to integrate. |
| [in] | freezing_rate | Freezing rate (only used for the constant rate scheme). |
| [in] | immersion_freezing_scheme_type | Immersion freezing scheme type. |
Definition at line 172 of file ice_nucleation.F90.
| subroutine pmc_ice_nucleation::ice_nucleation_immersion_freezing_time_dependent_naive | ( | type(aero_state_t), intent(inout) | aero_state, |
| type(aero_data_t), intent(in) | aero_data, | ||
| type(env_state_t), intent(inout) | env_state, | ||
| real(kind=dp), intent(in) | del_t, | ||
| integer, intent(in) | immersion_freezing_scheme_type, | ||
| real(kind=dp), intent(in) | freezing_rate | ||
| ) |
Simulation for time-dependent scheme (e.g., ABIFM, constant rate), deciding whether to freeze for each particle. Run in each time step. This subroutine applies the naive algorithm that checks each particle.
| [in,out] | aero_state | Aerosol state. |
| [in] | aero_data | Aerosol data. |
| [in,out] | env_state | Environment state. |
| [in] | del_t | Total time to integrate. |
| [in] | immersion_freezing_scheme_type | Type of the immersion freezing scheme |
| [in] | freezing_rate | Freezing rate (only used for the constant rate scheme). |
Definition at line 289 of file ice_nucleation.F90.
| subroutine pmc_ice_nucleation::ice_nucleation_melting | ( | type(aero_state_t), intent(inout) | aero_state, |
| type(aero_data_t), intent(in) | aero_data, | ||
| type(env_state_t), intent(inout) | env_state | ||
| ) |
Simulates melting: if the environmental temperature is above the freezing temperature of water, all particles are set to be unfrozen.
| [in,out] | aero_state | Aerosol state. |
| [in] | aero_data | Aerosol data. |
| [in,out] | env_state | Environment state. |
Definition at line 358 of file ice_nucleation.F90.
| subroutine pmc_ice_nucleation::ice_nucleation_singular_initialize | ( | type(aero_state_t), intent(inout) | aero_state, |
| type(aero_data_t), intent(in) | aero_data, | ||
| real(kind=dp), intent(in) | INAS_a, | ||
| real(kind=dp), intent(in) | INAS_b | ||
| ) |
Initialization for the singular scheme, sampling the freezing temperature for each particle.
| [in,out] | aero_state | Aerosol state. |
| [in] | aero_data | Aerosol data. |
| [in] | inas_a | Slope parameter for the INAS parameterization. |
| [in] | inas_b | Intercept parameter for the INAS parameterization. |
Definition at line 91 of file ice_nucleation.F90.
| integer, parameter pmc_ice_nucleation::immersion_freezing_scheme_abifm = 3 |
Type code for the ABIFM immersion freezing scheme.
Definition at line 27 of file ice_nucleation.F90.
| integer, parameter pmc_ice_nucleation::immersion_freezing_scheme_const = 1 |
Type code for constant ice nucleation rate (J_het) immersion freezing scheme.
Definition at line 23 of file ice_nucleation.F90.
| integer, parameter pmc_ice_nucleation::immersion_freezing_scheme_invalid = 0 |
Type code for an undefined or invalid immersion freezing scheme.
Definition at line 20 of file ice_nucleation.F90.
| integer, parameter pmc_ice_nucleation::immersion_freezing_scheme_singular = 2 |
Type code for the singular (INAS) immersion freezing scheme.
Definition at line 25 of file ice_nucleation.F90.