PartMC
2.6.1
|
Water condensation onto aerosol particles. More...
Data Types | |
type | condense_rates_inputs_t |
Internal-use structure for storing the inputs for the rate-calculation function. More... | |
type | condense_rates_outputs_t |
Internal-use structure for storing the outputs from the rate-calculation function. More... | |
Functions/Subroutines | |
subroutine | condense_particles (aero_state, aero_data, env_state_initial, env_state_final, del_t) |
Do condensation to all the particles for a given time interval, including updating the environment to account for the lost water vapor. More... | |
subroutine | condense_rates (inputs, outputs) |
Compute the rate of change of particle diameter and relative humidity for a single particle, together with the derivatives of the rates with respect to the input variables. More... | |
subroutine | condense_equilib_particle (env_state, aero_data, aero_particle) |
Determine the water equilibrium state of a single particle. More... | |
subroutine | condense_equilib_particles (env_state, aero_data, aero_state) |
Call condense_equilib_particle() on each particle in the aerosol to ensure that every particle has its water content in equilibrium. More... | |
Variables | |
logical, parameter | condense_do_test_jac_solve = .false. |
Whether to numerically test the Jacobian-solve function during execution (for debugging only). More... | |
logical, parameter | condense_do_test_counts = .false. |
Whether to print call-counts for helper routines during execution (for debugging only). More... | |
integer, parameter | pmc_condense_solver_success = 0 |
Result code indicating successful completion. More... | |
integer, parameter | pmc_condense_solver_init_y = 1 |
Result code indicating failure to allocate y vector. More... | |
integer, parameter | pmc_condense_solver_init_abstol = 2 |
Result code indicating failure to allocate abstol vector. More... | |
integer, parameter | pmc_condense_solver_init_cvode_mem = 3 |
Result code indicating failure to create the solver. More... | |
integer, parameter | pmc_condense_solver_init_cvode = 4 |
Result code indicating failure to initialize the solver. More... | |
integer, parameter | pmc_condense_solver_svtol = 5 |
Result code indicating failure to set tolerances. More... | |
integer, parameter | pmc_condense_solver_set_max_steps = 6 |
Result code indicating failure to set maximum steps. More... | |
integer, parameter | pmc_condense_solver_fail = 7 |
Result code indicating failure of the solver. More... | |
type(aero_data_t) | condense_saved_aero_data |
Internal-use variable for storing the aerosol data during calls to the ODE solver. More... | |
type(env_state_t) | condense_saved_env_state_initial |
Internal-use variable for storing the initial environment state during calls to the ODE solver. More... | |
real(kind=dp) | condense_saved_tdot |
Internal-use variable for storing the rate of change of the temperature during calls to the ODE solver. More... | |
real(kind=dp) | condense_saved_pdot |
Internal-use variable for storing the rate of change of the pressure during calls to the ODE solver. More... | |
real(kind=dp), dimension(:), allocatable | condense_saved_kappa |
Internal-use variable for storing the per-particle kappa values during calls to the ODE solver. More... | |
real(kind=dp), dimension(:), allocatable | condense_saved_d_dry |
Internal-use variable for storing the per-particle dry diameters during calls to the ODE solver. More... | |
real(kind=dp), dimension(:), allocatable | condense_saved_num_conc |
Internal-use variable for storing the per-particle number concentrations during calls to the ODE solver. More... | |
integer, save | condense_count_vf |
Internal-use variable for counting calls to the vector field subroutine. More... | |
integer, save | condense_count_solve |
Internal-use variable for counting calls to the Jacobian-solving subroutine. More... | |
Water condensation onto aerosol particles.
The model here assumes that the temperature and pressure are prescribed as functions of time, while water content per-particle and relative humidity are to be calculated by integrating their rates of change.
The state of the system is defined by the per-particle wet diameters and the relative humidity . The state vector stores these in the order . The time-derivative of the state vector and the Jacobian (derivative of the time-derivative with repsect to the state) all conform to this ordering.
The SUNDIALS ODE solver is used to compute the system evolution using an implicit method. The system Jacobian is explicitly inverted as its structure is very simple.
All equations used in this file are written in detail in the file doc/condense.tex
.
subroutine pmc_condense::condense_equilib_particle | ( | type(env_state_t), intent(in) | env_state, |
type(aero_data_t), intent(in) | aero_data, | ||
type(aero_particle_t), intent(inout) | aero_particle | ||
) |
Determine the water equilibrium state of a single particle.
[in] | env_state | Environment state. |
[in] | aero_data | Aerosol data. |
[in,out] | aero_particle | Particle. |
Definition at line 701 of file condense.F90.
subroutine pmc_condense::condense_equilib_particles | ( | type(env_state_t), intent(inout) | env_state, |
type(aero_data_t), intent(in) | aero_data, | ||
type(aero_state_t), intent(inout) | aero_state | ||
) |
Call condense_equilib_particle() on each particle in the aerosol to ensure that every particle has its water content in equilibrium.
[in,out] | env_state | Environment state. |
[in] | aero_data | Aerosol data. |
[in,out] | aero_state | Aerosol state. |
Definition at line 746 of file condense.F90.
subroutine pmc_condense::condense_particles | ( | type(aero_state_t), intent(inout) | aero_state, |
type(aero_data_t), intent(in) | aero_data, | ||
type(env_state_t), intent(in) | env_state_initial, | ||
type(env_state_t), intent(inout) | env_state_final, | ||
real(kind=dp), intent(in) | del_t | ||
) |
Do condensation to all the particles for a given time interval, including updating the environment to account for the lost water vapor.
[in,out] | aero_state | Aerosol state. |
[in] | aero_data | Aerosol data. |
[in] | env_state_initial | Environment state at the start of the timestep. |
[in,out] | env_state_final | Environment state at the end of the timestep. The rel_humid value will be ignored and overwritten with the condensation-computed value. |
[in] | del_t | Total time to integrate. |
Definition at line 147 of file condense.F90.
subroutine pmc_condense::condense_rates | ( | type(condense_rates_inputs_t), intent(in) | inputs, |
type(condense_rates_outputs_t), intent(out) | outputs | ||
) |
Compute the rate of change of particle diameter and relative humidity for a single particle, together with the derivatives of the rates with respect to the input variables.
[in] | inputs | Inputs to rates. |
[out] | outputs | Outputs rates. |
Definition at line 371 of file condense.F90.
integer, save pmc_condense::condense_count_solve |
Internal-use variable for counting calls to the Jacobian-solving subroutine.
Definition at line 138 of file condense.F90.
integer, save pmc_condense::condense_count_vf |
Internal-use variable for counting calls to the vector field subroutine.
Definition at line 135 of file condense.F90.
logical, parameter pmc_condense::condense_do_test_counts = .false. |
Whether to print call-counts for helper routines during execution (for debugging only).
Definition at line 46 of file condense.F90.
logical, parameter pmc_condense::condense_do_test_jac_solve = .false. |
Whether to numerically test the Jacobian-solve function during execution (for debugging only).
Definition at line 43 of file condense.F90.
type(aero_data_t) pmc_condense::condense_saved_aero_data |
Internal-use variable for storing the aerosol data during calls to the ODE solver.
Definition at line 113 of file condense.F90.
real(kind=dp), dimension(:), allocatable pmc_condense::condense_saved_d_dry |
Internal-use variable for storing the per-particle dry diameters during calls to the ODE solver.
Definition at line 128 of file condense.F90.
type(env_state_t) pmc_condense::condense_saved_env_state_initial |
Internal-use variable for storing the initial environment state during calls to the ODE solver.
Definition at line 116 of file condense.F90.
real(kind=dp), dimension(:), allocatable pmc_condense::condense_saved_kappa |
Internal-use variable for storing the per-particle kappa values during calls to the ODE solver.
Definition at line 125 of file condense.F90.
real(kind=dp), dimension(:), allocatable pmc_condense::condense_saved_num_conc |
Internal-use variable for storing the per-particle number concentrations during calls to the ODE solver.
Definition at line 131 of file condense.F90.
real(kind=dp) pmc_condense::condense_saved_pdot |
Internal-use variable for storing the rate of change of the pressure during calls to the ODE solver.
Definition at line 122 of file condense.F90.
real(kind=dp) pmc_condense::condense_saved_tdot |
Internal-use variable for storing the rate of change of the temperature during calls to the ODE solver.
Definition at line 119 of file condense.F90.
integer, parameter pmc_condense::pmc_condense_solver_fail = 7 |
Result code indicating failure of the solver.
Definition at line 63 of file condense.F90.
integer, parameter pmc_condense::pmc_condense_solver_init_abstol = 2 |
Result code indicating failure to allocate abstol
vector.
Definition at line 53 of file condense.F90.
integer, parameter pmc_condense::pmc_condense_solver_init_cvode = 4 |
Result code indicating failure to initialize the solver.
Definition at line 57 of file condense.F90.
integer, parameter pmc_condense::pmc_condense_solver_init_cvode_mem = 3 |
Result code indicating failure to create the solver.
Definition at line 55 of file condense.F90.
integer, parameter pmc_condense::pmc_condense_solver_init_y = 1 |
Result code indicating failure to allocate y
vector.
Definition at line 51 of file condense.F90.
integer, parameter pmc_condense::pmc_condense_solver_set_max_steps = 6 |
Result code indicating failure to set maximum steps.
Definition at line 61 of file condense.F90.
integer, parameter pmc_condense::pmc_condense_solver_success = 0 |
Result code indicating successful completion.
Definition at line 49 of file condense.F90.
integer, parameter pmc_condense::pmc_condense_solver_svtol = 5 |
Result code indicating failure to set tolerances.
Definition at line 59 of file condense.F90.