PartMC
2.1.5
|
00001 ! Copyright (C) 2005-2011 Nicole Riemer and Matthew West 00002 ! Copyright (C) 2009 Joseph Ching 00003 ! Licensed under the GNU General Public License version 2 or (at your 00004 ! option) any later version. See the file COPYING for details. 00005 00006 !> \file 00007 !> The pmc_condense module. 00008 00009 !> Water condensation onto aerosol particles. 00010 !! 00011 !! The model here assumes that the temperature \f$ T \f$ and pressure 00012 !! \f$ p \f$ are prescribed as functions of time, while water content 00013 !! per-particle and relative humidity are to be calculated by 00014 !! integrating their rates of change. 00015 !! 00016 !! The state of the system is defined by the per-particle wet 00017 !! diameters \f$ D_i \f$ and the relative humidity \f$ H \f$. The 00018 !! state vector stores these in the order \f$ (D_1,\ldots,D_n,H) 00019 !! \f$. The time-derivative of the state vector and the Jacobian 00020 !! (derivative of the time-derivative with repsect to the state) all 00021 !! conform to this ordering. 00022 !! 00023 !! The SUNDIALS ODE solver is used to compute the system evolution 00024 !! using an implicit method. The system Jacobian is explicitly 00025 !! inverveted as its structure is very simple. 00026 !! 00027 !! All equations used in this file are written in detail in the file 00028 !! \c doc/condense.tex. 00029 module pmc_condense 00030 00031 use pmc_aero_state 00032 use pmc_bin_grid 00033 use pmc_env_data 00034 use pmc_env_state 00035 use pmc_aero_data 00036 use pmc_util 00037 use pmc_aero_particle 00038 use pmc_constants 00039 #ifdef PMC_USE_SUNDIALS 00040 use iso_c_binding 00041 #endif 00042 00043 !> Whether to numerically test the Jacobian-solve function during 00044 !> execution (for debugging only). 00045 logical, parameter :: CONDENSE_DO_TEST_JAC_SOLVE = .false. 00046 !> Whether to print call-counts for helper routines during execution 00047 !> (for debugging only). 00048 logical, parameter :: CONDENSE_DO_TEST_COUNTS = .false. 00049 00050 !> Result code indicating successful completion. 00051 integer, parameter :: PMC_CONDENSE_SOLVER_SUCCESS = 0 00052 !> Result code indicating failure to allocate \c y vector. 00053 integer, parameter :: PMC_CONDENSE_SOLVER_INIT_Y = 1 00054 !> Result code indicating failure to allocate \c abstol vector. 00055 integer, parameter :: PMC_CONDENSE_SOLVER_INIT_ABSTOL = 2 00056 !> Result code indicating failure to create the solver. 00057 integer, parameter :: PMC_CONDENSE_SOLVER_INIT_CVODE_MEM = 3 00058 !> Result code indicating failure to initialize the solver. 00059 integer, parameter :: PMC_CONDENSE_SOLVER_INIT_CVODE = 4 00060 !> Result code indicating failure to set tolerances. 00061 integer, parameter :: PMC_CONDENSE_SOLVER_SVTOL = 5 00062 !> Result code indicating failure to set maximum steps. 00063 integer, parameter :: PMC_CONDENSE_SOLVER_SET_MAX_STEPS = 6 00064 !> Result code indicating failure of the solver. 00065 integer, parameter :: PMC_CONDENSE_SOLVER_FAIL = 7 00066 00067 !> Internal-use structure for storing the inputs for the 00068 !> rate-calculation function. 00069 type condense_rates_inputs_t 00070 !> Temperature (K). 00071 real(kind=dp) :: T 00072 !> Rate of change of temperature (K s^{-1}). 00073 real(kind=dp) :: Tdot 00074 !> Relative humidity (1). 00075 real(kind=dp) :: H 00076 !> Pressure (Pa). 00077 real(kind=dp) :: p 00078 !> Computational volume (m^3). 00079 real(kind=dp) :: V_comp 00080 !> Particle diameter (m). 00081 real(kind=dp) :: D 00082 !> Particle dry diameter (m). 00083 real(kind=dp) :: D_dry 00084 !> Kappa parameter (1). 00085 real(kind=dp) :: kappa 00086 end type condense_rates_inputs_t 00087 00088 !> Internal-use structure for storing the outputs from the 00089 !> rate-calculation function. 00090 type condense_rates_outputs_t 00091 !> Change rate of diameter (m s^{-1}). 00092 real(kind=dp) :: Ddot 00093 !> Change rate of relative humidity due to this particle (s^{-1}). 00094 real(kind=dp) :: Hdot_i 00095 !> Change rate of relative humidity due to environment changes (s^{-1}). 00096 real(kind=dp) :: Hdot_env 00097 !> Sensitivity of \c Ddot to input \c D (m s^{-1} m^{-1}). 00098 real(kind=dp) :: dDdot_dD 00099 !> Sensitivity of \c Ddot to input \c H (m s^{-1}). 00100 real(kind=dp) :: dDdot_dH 00101 !> Sensitivity of \c Hdot_i to input \c D (s^{-1} m^{-1}). 00102 real(kind=dp) :: dHdoti_dD 00103 !> Sensitivity of \c Hdot_i to input \c D (s^{-1}). 00104 real(kind=dp) :: dHdoti_dH 00105 !> Sensitivity of \c Hdot_env to input \c D (s^{-1} m^{-1}). 00106 real(kind=dp) :: dHdotenv_dD 00107 !> Sensitivity of \c Hdot_env to input \c D (s^{-1}). 00108 real(kind=dp) :: dHdotenv_dH 00109 end type condense_rates_outputs_t 00110 00111 !> Internal-use variable for storing the aerosol data during calls 00112 !> to the ODE solver. 00113 type(aero_data_t) :: condense_saved_aero_data 00114 !> Internal-use variable for storing the environment data during 00115 !> calls to the ODE solver. 00116 type(env_data_t) :: condense_saved_env_data 00117 !> Internal-use variable for storing the initial environment state 00118 !> during calls to the ODE solver. 00119 type(env_state_t) :: condense_saved_env_state_initial 00120 !> Internal-use variable for storing the inital computational volume 00121 !> during calls to the ODE solver. 00122 real(kind=dp) :: condense_saved_V_comp_initial 00123 !> Internal-use variable for storing the rate of change of the 00124 !> temperature during calls to the ODE solver. 00125 real(kind=dp) :: condense_saved_Tdot 00126 !> Internal-use variable for storing the per-particle kappa values 00127 !> during calls to the ODE solver. 00128 real(kind=dp), allocatable :: condense_saved_kappa(:) 00129 !> Internal-use variable for storing the per-particle dry diameters 00130 !> during calls to the ODE solver. 00131 real(kind=dp), allocatable :: condense_saved_D_dry(:) 00132 !> Internal-use variable for storing the per-particle weights during 00133 !> calls to the ODE solver. 00134 real(kind=dp), allocatable :: condense_saved_weight(:) 00135 00136 !> Internal-use variable for counting calls to the vector field 00137 !> subroutine. 00138 integer, save :: condense_count_vf 00139 !> Internal-use variable for counting calls to the Jacobian-solving 00140 !> subroutine. 00141 integer, save :: condense_count_solve 00142 00143 contains 00144 00145 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00146 00147 !> Do condensation to all the particles for a given time interval, 00148 !> including updating the environment to account for the lost 00149 !> water vapor. 00150 subroutine condense_particles(bin_grid, env_state, env_data, & 00151 aero_data, aero_weight, aero_state, del_t) 00152 00153 !> Bin grid. 00154 type(bin_grid_t), intent(in) :: bin_grid 00155 !> Environment state. 00156 type(env_state_t), intent(inout) :: env_state 00157 !> Environment data. 00158 type(env_data_t), intent(in) :: env_data 00159 !> Aerosol data. 00160 type(aero_data_t), intent(in) :: aero_data 00161 !> Aerosol weight. 00162 type(aero_weight_t), intent(in) :: aero_weight 00163 !> Aerosol state. 00164 type(aero_state_t), intent(inout) :: aero_state 00165 !> Total time to integrate. 00166 real(kind=dp), intent(in) :: del_t 00167 00168 integer :: i_bin, i_part, n_eqn, i_eqn 00169 type(aero_particle_t), pointer :: aero_particle 00170 real(kind=dp) :: state(aero_state%n_part + 1), init_time, final_time 00171 real(kind=dp) :: abs_tol_vector(aero_state%n_part + 1) 00172 real(kind=dp) :: weight, old_weight, new_weight 00173 type(env_state_t) :: env_state_final 00174 real(kind=dp) :: water_vol_initial, water_vol_final, d_water_vol 00175 real(kind=dp) :: vapor_vol_initial, vapor_vol_final, d_vapor_vol 00176 real(kind=dp) :: V_comp_final, water_rel_error 00177 #ifdef PMC_USE_SUNDIALS 00178 real(kind=c_double), target :: state_f(aero_state%n_part + 1) 00179 real(kind=c_double), target :: abstol_f(aero_state%n_part + 1) 00180 type(c_ptr) :: state_f_p, abstol_f_p 00181 integer(kind=c_int) :: n_eqn_f, solver_stat 00182 real(kind=c_double) :: reltol_f, t_initial_f, t_final_f 00183 #endif 00184 type(aero_particle_t) :: new_aero_particle 00185 type(aero_info_t) :: aero_info 00186 integer :: n_copies, i_dup 00187 00188 #ifdef PMC_USE_SUNDIALS 00189 #ifndef DOXYGEN_SKIP_DOC 00190 interface 00191 integer(kind=c_int) function condense_solver(neq, x_f, abstol_f, & 00192 reltol_f, t_initial_f, t_final_f) bind(c) 00193 use iso_c_binding 00194 integer(kind=c_int), value :: neq 00195 type(c_ptr), value :: x_f 00196 type(c_ptr), value :: abstol_f 00197 real(kind=c_double), value :: reltol_f 00198 real(kind=c_double), value :: t_initial_f 00199 real(kind=c_double), value :: t_final_f 00200 end function condense_solver 00201 end interface 00202 #endif 00203 #endif 00204 00205 ! initial water volume in the aerosol particles in volume V_comp 00206 water_vol_initial = 0d0 00207 do i_bin = 1,bin_grid%n_bin 00208 do i_part = 1,aero_state%bin(i_bin)%n_part 00209 aero_particle => aero_state%bin(i_bin)%particle(i_part) 00210 weight = aero_weight_value(aero_weight, & 00211 aero_particle_radius(aero_particle)) 00212 water_vol_initial = water_vol_initial & 00213 + aero_particle%vol(aero_data%i_water) * weight 00214 end do 00215 end do 00216 00217 ! save data for use within the timestepper 00218 call aero_data_allocate(condense_saved_aero_data) 00219 call env_data_allocate(condense_saved_env_data) 00220 call env_state_allocate(condense_saved_env_state_initial) 00221 00222 call aero_data_copy(aero_data, condense_saved_aero_data) 00223 call env_data_copy(env_data, condense_saved_env_data) 00224 call env_state_copy(env_state, condense_saved_env_state_initial) 00225 00226 condense_saved_V_comp_initial = aero_state%comp_vol 00227 00228 call env_state_allocate(env_state_final) 00229 call env_state_copy(env_state, env_state_final) 00230 call env_data_update_state(env_data, env_state_final, & 00231 env_state_final%elapsed_time + del_t, update_rel_humid = .false.) 00232 condense_saved_Tdot = (env_state_final%temp - env_state%temp) / del_t 00233 00234 ! construct initial state vector from aero_state and env_state 00235 allocate(condense_saved_kappa(aero_state%n_part)) 00236 allocate(condense_saved_D_dry(aero_state%n_part)) 00237 allocate(condense_saved_weight(aero_state%n_part)) 00238 i_eqn = 0 00239 do i_bin = 1,bin_grid%n_bin 00240 ! work backwards for consistency with the later weight 00241 ! adjustment, which has specific ordering requirements 00242 do i_part = aero_state%bin(i_bin)%n_part,1,-1 00243 i_eqn = i_eqn + 1 00244 aero_particle => aero_state%bin(i_bin)%particle(i_part) 00245 condense_saved_kappa(i_eqn) & 00246 = aero_particle_solute_kappa(aero_particle, aero_data) 00247 condense_saved_D_dry(i_eqn) = vol2diam(& 00248 aero_particle_solute_volume(aero_particle, aero_data)) 00249 condense_saved_weight(i_eqn) = aero_weight_value(aero_weight, & 00250 aero_particle_radius(aero_particle)) 00251 state(i_eqn) = aero_particle_diameter(aero_particle) 00252 abs_tol_vector(i_eqn) = max(1d-30, & 00253 1d-8 * (state(i_eqn) - condense_saved_D_dry(i_eqn))) 00254 end do 00255 end do 00256 state(aero_state%n_part + 1) = env_state%rel_humid 00257 abs_tol_vector(aero_state%n_part + 1) = 1d-10 00258 00259 #ifdef PMC_USE_SUNDIALS 00260 ! call SUNDIALS solver 00261 n_eqn = aero_state%n_part + 1 00262 n_eqn_f = int(n_eqn, kind=c_int) 00263 reltol_f = real(1d-8, kind=c_double) 00264 t_initial_f = real(0, kind=c_double) 00265 t_final_f = real(del_t, kind=c_double) 00266 do i_eqn = 1,n_eqn 00267 state_f(i_eqn) = real(state(i_eqn), kind=c_double) 00268 abstol_f(i_eqn) = real(abs_tol_vector(i_eqn), kind=c_double) 00269 end do 00270 state_f_p = c_loc(state_f) 00271 abstol_f_p = c_loc(abstol_f) 00272 condense_count_vf = 0 00273 condense_count_solve = 0 00274 solver_stat = condense_solver(n_eqn_f, state_f_p, abstol_f_p, reltol_f, & 00275 t_initial_f, t_final_f) 00276 call condense_check_solve(solver_stat) 00277 if (CONDENSE_DO_TEST_COUNTS) then 00278 write(0,*) 'condense_count_vf ', condense_count_vf 00279 write(0,*) 'condense_count_solve ', condense_count_solve 00280 end if 00281 do i_eqn = 1,n_eqn 00282 state(i_eqn) = real(state_f(i_eqn), kind=dp) 00283 end do 00284 #endif 00285 00286 ! unpack result state vector into env_state 00287 env_state%rel_humid = state(aero_state%n_part + 1) 00288 00289 ! unpack result state vector into aero_state, compute the final 00290 ! water volume in the aerosol particles in volume V_comp, and 00291 ! adjust particle number to account for weight changes 00292 call aero_particle_allocate(new_aero_particle) 00293 call aero_info_allocate(aero_info) 00294 water_vol_final = 0d0 00295 i_eqn = 0 00296 do i_bin = 1,bin_grid%n_bin 00297 ! work backwards so any additions and removals will only affect 00298 ! particles that we've already dealt with 00299 do i_part = aero_state%bin(i_bin)%n_part,1,-1 00300 i_eqn = i_eqn + 1 00301 aero_particle => aero_state%bin(i_bin)%particle(i_part) 00302 old_weight = aero_weight_value(aero_weight, & 00303 aero_particle_radius(aero_particle)) 00304 00305 ! translate output back to particle 00306 aero_particle%vol(aero_data%i_water) = diam2vol(state(i_eqn)) & 00307 - aero_particle_solute_volume(aero_particle, aero_data) 00308 00309 ! ensure volumes stay positive 00310 aero_particle%vol(aero_data%i_water) = max(0d0, & 00311 aero_particle%vol(aero_data%i_water)) 00312 00313 ! add up total water volume, using old weights 00314 water_vol_final = water_vol_final & 00315 + aero_particle%vol(aero_data%i_water) * old_weight 00316 00317 ! adjust particle number to account for weight changes 00318 if (aero_weight%type /= AERO_WEIGHT_TYPE_NONE) then 00319 new_weight = aero_weight_value(aero_weight, & 00320 aero_particle_radius(aero_particle)) 00321 n_copies = prob_round(old_weight / new_weight) 00322 if (n_copies == 0) then 00323 aero_info%id = aero_particle%id 00324 aero_info%action = AERO_INFO_WEIGHT 00325 aero_info%other_id = 0 00326 call aero_state_remove_particle_with_info(aero_state, & 00327 i_bin, i_part, aero_info) 00328 elseif (n_copies > 1) then 00329 do i_dup = 1,(n_copies - 1) 00330 call aero_particle_copy(aero_particle, new_aero_particle) 00331 call aero_particle_new_id(new_aero_particle) 00332 ! this might be adding into the wrong bin, but 00333 ! that's necessary as we might not have processed 00334 ! the correct bin yet. 00335 call aero_state_add_particle(aero_state, i_bin, & 00336 new_aero_particle) 00337 ! re-get the particle pointer, which may have 00338 ! changed due to reallocations caused by adding 00339 aero_particle => aero_state%bin(i_bin)%particle(i_part) 00340 end do 00341 end if 00342 end if 00343 end do 00344 end do 00345 ! We've modified particle diameters, so we need to update which 00346 ! bins they are in. 00347 call aero_state_resort(bin_grid, aero_state) 00348 call aero_particle_deallocate(new_aero_particle) 00349 call aero_info_deallocate(aero_info) 00350 00351 ! Check that water removed from particles equals water added to 00352 ! vapor. Note that water concentration is not conserved (due to 00353 ! V_comp changes), and we need to consider particle weightings 00354 ! correctly. 00355 V_comp_final = condense_saved_V_comp_initial & 00356 * env_state_final%temp / condense_saved_env_state_initial%temp 00357 vapor_vol_initial = aero_data%molec_weight(aero_data%i_water) & 00358 / (const%univ_gas_const * condense_saved_env_state_initial%temp) & 00359 * env_state_sat_vapor_pressure(condense_saved_env_state_initial) & 00360 * condense_saved_env_state_initial%rel_humid & 00361 * condense_saved_V_comp_initial & 00362 / aero_particle_water_density(aero_data) 00363 vapor_vol_final = aero_data%molec_weight(aero_data%i_water) & 00364 / (const%univ_gas_const * env_state_final%temp) & 00365 * env_state_sat_vapor_pressure(env_state_final) & 00366 * env_state%rel_humid & 00367 * V_comp_final / aero_particle_water_density(aero_data) 00368 d_vapor_vol = vapor_vol_final - vapor_vol_initial 00369 d_water_vol = water_vol_final - water_vol_initial 00370 water_rel_error = (d_vapor_vol + d_water_vol) & 00371 / (vapor_vol_final + water_vol_final) 00372 call warn_assert_msg(477865387, abs(water_rel_error) < 1d-6, & 00373 "condensation water imbalance too high: " & 00374 // trim(real_to_string(water_rel_error))) 00375 00376 deallocate(condense_saved_kappa) 00377 deallocate(condense_saved_D_dry) 00378 deallocate(condense_saved_weight) 00379 call env_state_deallocate(env_state_final) 00380 call aero_data_deallocate(condense_saved_aero_data) 00381 call env_data_deallocate(condense_saved_env_data) 00382 call env_state_deallocate(condense_saved_env_state_initial) 00383 00384 end subroutine condense_particles 00385 00386 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00387 00388 #ifdef PMC_USE_SUNDIALS 00389 !> Check the return code from the condense_solver() function. 00390 subroutine condense_check_solve(value) 00391 00392 !> Return code to check. 00393 integer(kind=c_int), intent(in) :: value 00394 00395 if (value == PMC_CONDENSE_SOLVER_SUCCESS) then 00396 return 00397 elseif (value == PMC_CONDENSE_SOLVER_INIT_Y) then 00398 call die_msg(123749472, "condense_solver: " & 00399 // "failed to allocate y vector") 00400 elseif (value == PMC_CONDENSE_SOLVER_INIT_ABSTOL) then 00401 call die_msg(563665949, "condense_solver: " & 00402 // "failed to allocate abstol vector") 00403 elseif (value == PMC_CONDENSE_SOLVER_INIT_CVODE_MEM) then 00404 call die_msg(700541443, "condense_solver: " & 00405 // "failed to create the solver") 00406 elseif (value == PMC_CONDENSE_SOLVER_INIT_CVODE) then 00407 call die_msg(297559183, "condense_solver: " & 00408 // "failure to initialize the solver") 00409 elseif (value == PMC_CONDENSE_SOLVER_SVTOL) then 00410 call die_msg(848342417, "condense_solver: " & 00411 // "failed to set tolerances") 00412 elseif (value == PMC_CONDENSE_SOLVER_SET_MAX_STEPS) then 00413 call die_msg(275591501, "condense_solver: " & 00414 // "failed to set maximum steps") 00415 elseif (value == PMC_CONDENSE_SOLVER_FAIL) then 00416 call die_msg(862254233, "condense_solver: solver failed") 00417 else 00418 call die_msg(635697577, "condense_solver: unknown return code: " & 00419 // trim(integer_to_string(value))) 00420 end if 00421 00422 end subroutine condense_check_solve 00423 #endif 00424 00425 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00426 00427 !> Fills in the \c env_state with the current environment state, 00428 !> taken from the \c state vector and from global variables. 00429 subroutine condense_current_env_state(n_eqn, time, state, env_state) 00430 00431 !> Length of state vector. 00432 integer, intent(in) :: n_eqn 00433 !> Current time (s). 00434 real(kind=dp), intent(in) :: time 00435 !> Current state vector. 00436 real(kind=dp), intent(in) :: state(n_eqn) 00437 !> Current environment state. 00438 type(env_state_t), intent(inout) :: env_state 00439 00440 call env_state_copy(condense_saved_env_state_initial, env_state) 00441 call env_data_update_state(condense_saved_env_data, & 00442 env_state, env_state%elapsed_time + time, & 00443 update_rel_humid = .false.) 00444 env_state%rel_humid = state(n_eqn) 00445 00446 end subroutine condense_current_env_state 00447 00448 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00449 00450 !> Compute the rate of change of particle diameter and relative 00451 !> humidity for a single particle, together with the derivatives of 00452 !> the rates with respect to the input variables. 00453 subroutine condense_rates(inputs, outputs) 00454 00455 !> Inputs to rates. 00456 type(condense_rates_inputs_t), intent(in) :: inputs 00457 !> Outputs rates. 00458 type(condense_rates_outputs_t), intent(out) :: outputs 00459 00460 real(kind=dp) :: rho_w, M_w, P_0, dP0_dT_div_P0, rho_air, k_a, D_v, U 00461 real(kind=dp) :: V, W, X, Y, Z, k_ap, dkap_dD, D_vp, dDvp_dD 00462 real(kind=dp) :: a_w, daw_dD, delta_star, h, dh_ddelta, dh_dD 00463 real(kind=dp) :: dh_dH, ddeltastar_dD, ddeltastar_dH 00464 integer :: newton_step 00465 00466 rho_w = const%water_density 00467 M_w = const%water_molec_weight 00468 P_0 = const%water_eq_vap_press & 00469 * 10d0**(7.45d0 * (inputs%T - const%water_freeze_temp) & 00470 / (inputs%T - 38d0)) 00471 dP0_dT_div_P0 = 7.45d0 * log(10d0) * (const%water_freeze_temp - 38d0) & 00472 / (inputs%T - 38d0)**2 00473 rho_air = const%air_molec_weight * inputs%p & 00474 / (const%univ_gas_const * inputs%T) 00475 00476 k_a = 1d-3 * (4.39d0 + 0.071d0 * inputs%T) 00477 D_v = 0.211d-4 / (inputs%p / const%air_std_press) & 00478 * (inputs%T / 273d0)**1.94d0 00479 U = const%water_latent_heat * rho_w / (4d0 * inputs%T) 00480 V = 4d0 * M_w * P_0 / (rho_w * const%univ_gas_const * inputs%T) 00481 W = const%water_latent_heat * M_w / (const%univ_gas_const * inputs%T) 00482 X = 4d0 * M_w * const%water_surf_eng & 00483 / (const%univ_gas_const * inputs%T * rho_w) 00484 Y = 2d0 * k_a / (const%accom_coeff * rho_air & 00485 * const%air_spec_heat) & 00486 * sqrt(2d0 * const%pi * const%air_molec_weight & 00487 / (const%univ_gas_const * inputs%T)) 00488 Z = 2d0 * D_v / const%accom_coeff * sqrt(2d0 * const%pi * M_w & 00489 / (const%univ_gas_const * inputs%T)) 00490 00491 outputs%Hdot_env = - dP0_dT_div_P0 * inputs%Tdot * inputs%H 00492 outputs%dHdotenv_dD = 0d0 00493 outputs%dHdotenv_dH = - dP0_dT_div_P0 * inputs%Tdot 00494 00495 if (inputs%D <= inputs%D_dry) then 00496 k_ap = k_a / (1d0 + Y / inputs%D_dry) 00497 dkap_dD = 0d0 00498 D_vp = D_v / (1d0 + Z / inputs%D_dry) 00499 dDvp_dD = 0d0 00500 a_w = 0d0 00501 daw_dD = 0d0 00502 00503 delta_star = U * V * D_vp * inputs%H / k_ap 00504 00505 outputs%Ddot = k_ap * delta_star / (U * inputs%D_dry) 00506 outputs%Hdot_i = - 2d0 * const%pi / (V * inputs%V_comp) & 00507 * inputs%D_dry**2 * outputs%Ddot 00508 00509 dh_ddelta = k_ap 00510 dh_dD = 0d0 00511 dh_dH = - U * V * D_vp 00512 00513 ddeltastar_dD = - dh_dD / dh_ddelta 00514 ddeltastar_dH = - dh_dH / dh_ddelta 00515 00516 outputs%dDdot_dD = 0d0 00517 outputs%dDdot_dH = k_ap / (U * inputs%D_dry) * ddeltastar_dH 00518 outputs%dHdoti_dD = - 2d0 * const%pi / (V * inputs%V_comp) & 00519 * inputs%D_dry**2 * outputs%dDdot_dD 00520 outputs%dHdoti_dH = - 2d0 * const%pi / (V * inputs%V_comp) & 00521 * inputs%D_dry**2 * outputs%dDdot_dH 00522 00523 return 00524 end if 00525 00526 k_ap = k_a / (1d0 + Y / inputs%D) 00527 dkap_dD = k_a * Y / (inputs%D + Y)**2 00528 D_vp = D_v / (1d0 + Z / inputs%D) 00529 dDvp_dD = D_v * Z / (inputs%D + Z)**2 00530 a_w = (inputs%D**3 - inputs%D_dry**3) & 00531 / (inputs%D**3 + (inputs%kappa - 1d0) * inputs%D_dry**3) 00532 daw_dD = 3d0 * inputs%D**2 * inputs%kappa * inputs%D_dry**3 & 00533 / (inputs%D**3 + (inputs%kappa - 1d0) * inputs%D_dry**3)**2 00534 00535 delta_star = 0d0 00536 h = 0d0 00537 dh_ddelta = 1d0 00538 do newton_step = 1,5 00539 ! update delta_star first so when the newton loop ends we have 00540 ! h and dh_ddelta evaluated at the final delta_star value 00541 delta_star = delta_star - h / dh_ddelta 00542 h = k_ap * delta_star - U * V * D_vp & 00543 * (inputs%H - a_w / (1d0 + delta_star) & 00544 * exp(W * delta_star / (1d0 + delta_star) & 00545 + (X / inputs%D) / (1d0 + delta_star))) 00546 dh_ddelta = & 00547 k_ap - U * V * D_vp * a_w / (1d0 + delta_star)**2 & 00548 * (1d0 - W / (1d0 + delta_star) & 00549 + (X / inputs%D) / (1d0 + delta_star)) & 00550 * exp(W * delta_star / (1d0 + delta_star) & 00551 + (X / inputs%D) / (1d0 + delta_star)) 00552 end do 00553 call warn_assert_msg(387362320, & 00554 abs(h) < 1d3 * epsilon(1d0) * abs(U * V * D_vp * inputs%H), & 00555 "condensation newton loop did not satisfy convergence tolerance") 00556 00557 outputs%Ddot = k_ap * delta_star / (U * inputs%D) 00558 outputs%Hdot_i = - 2d0 * const%pi / (V * inputs%V_comp) & 00559 * inputs%D**2 * outputs%Ddot 00560 00561 dh_dD = dkap_dD * delta_star & 00562 - U * V * dDvp_dD * inputs%H + U * V & 00563 * (a_w * dDvp_dD + D_vp * daw_dD & 00564 - D_vp * a_w * (X / inputs%D**2) / (1d0 + delta_star)) & 00565 * (1d0 / (1d0 + delta_star)) & 00566 * exp((W * delta_star) / (1d0 + delta_star) & 00567 + (X / inputs%D) / (1d0 + delta_star)) 00568 dh_dH = - U * V * D_vp 00569 00570 ddeltastar_dD = - dh_dD / dh_ddelta 00571 ddeltastar_dH = - dh_dH / dh_ddelta 00572 00573 outputs%dDdot_dD = dkap_dD * delta_star / (U * inputs%D) & 00574 + k_ap * ddeltastar_dD / (U * inputs%D) & 00575 - k_ap * delta_star / (U * inputs%D**2) 00576 outputs%dDdot_dH = k_ap / (U * inputs%D) * ddeltastar_dH 00577 outputs%dHdoti_dD = - 2d0 * const%pi / (V * inputs%V_comp) & 00578 * (2d0 * inputs%D * outputs%Ddot + inputs%D**2 * outputs%dDdot_dD) 00579 outputs%dHdoti_dH = - 2d0 * const%pi / (V * inputs%V_comp) & 00580 * inputs%D**2 * outputs%dDdot_dH 00581 00582 end subroutine condense_rates 00583 00584 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00585 00586 #ifdef PMC_USE_SUNDIALS 00587 !> Compute the condensation rates (Ddot and Hdot) at the current 00588 !> value of the state (D and H). 00589 subroutine condense_vf_f(n_eqn, time, state_p, state_dot_p) bind(c) 00590 00591 !> Length of state vector. 00592 integer(kind=c_int), value, intent(in) :: n_eqn 00593 !> Current time (s). 00594 real(kind=c_double), value, intent(in) :: time 00595 !> Pointer to state data. 00596 type(c_ptr), value, intent(in) :: state_p 00597 !> Pointer to state_dot data. 00598 type(c_ptr), value, intent(in) :: state_dot_p 00599 00600 real(kind=c_double), pointer :: state(:) 00601 real(kind=c_double), pointer :: state_dot(:) 00602 real(kind=dp) :: Hdot 00603 integer :: i_part 00604 type(env_state_t) :: env_state 00605 type(condense_rates_inputs_t) :: inputs 00606 type(condense_rates_outputs_t) :: outputs 00607 00608 condense_count_vf = condense_count_vf + 1 00609 00610 call c_f_pointer(state_p, state, (/ n_eqn /)) 00611 call c_f_pointer(state_dot_p, state_dot, (/ n_eqn /)) 00612 00613 call env_state_allocate(env_state) 00614 call condense_current_env_state(n_eqn, time, state, env_state) 00615 00616 inputs%T = env_state%temp 00617 inputs%Tdot = condense_saved_Tdot 00618 inputs%H = env_state%rel_humid 00619 inputs%p = env_state%pressure 00620 00621 Hdot = 0d0 00622 do i_part = 1,(n_eqn - 1) 00623 inputs%D = state(i_part) 00624 inputs%D_dry = condense_saved_D_dry(i_part) 00625 inputs%V_comp = condense_saved_V_comp_initial & 00626 * env_state%temp / condense_saved_env_state_initial%temp & 00627 / condense_saved_weight(i_part) 00628 inputs%kappa = condense_saved_kappa(i_part) 00629 call condense_rates(inputs, outputs) 00630 state_dot(i_part) = outputs%Ddot 00631 Hdot = Hdot + outputs%Hdot_i 00632 end do 00633 Hdot = Hdot + outputs%Hdot_env 00634 00635 state_dot(n_eqn) = Hdot 00636 00637 call env_state_deallocate(env_state) 00638 00639 end subroutine condense_vf_f 00640 #endif 00641 00642 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00643 00644 #ifdef PMC_USE_SUNDIALS 00645 !> Compute the Jacobian given by the derivatives of the condensation 00646 !> rates (Ddot and Hdot) with respect to the input variables (D and 00647 !> H). 00648 subroutine condense_jac(n_eqn, time, state_p, dDdot_dD, dDdot_dH, & 00649 dHdot_dD, dHdot_dH) 00650 00651 !> Length of state vector. 00652 integer(kind=c_int), intent(in) :: n_eqn 00653 !> Current time (s). 00654 real(kind=c_double), intent(in) :: time 00655 !> Pointer to current state vector. 00656 type(c_ptr), intent(in) :: state_p 00657 !> Derivative of Ddot with respect to D. 00658 real(kind=dp), intent(out) :: dDdot_dD(n_eqn - 1) 00659 !> Derivative of Ddot with respect to H. 00660 real(kind=dp), intent(out) :: dDdot_dH(n_eqn - 1) 00661 !> Derivative of Hdot with respect to D. 00662 real(kind=dp), intent(out) :: dHdot_dD(n_eqn - 1) 00663 !> Derivative of Hdot with respect to H. 00664 real(kind=dp), intent(out) :: dHdot_dH 00665 00666 real(kind=c_double), pointer :: state(:) 00667 integer :: i_part 00668 type(env_state_t) :: env_state 00669 type(condense_rates_inputs_t) :: inputs 00670 type(condense_rates_outputs_t) :: outputs 00671 00672 call c_f_pointer(state_p, state, (/ n_eqn /)) 00673 00674 call env_state_allocate(env_state) 00675 call condense_current_env_state(n_eqn, time, state, env_state) 00676 00677 inputs%T = env_state%temp 00678 inputs%Tdot = condense_saved_Tdot 00679 inputs%H = env_state%rel_humid 00680 inputs%p = env_state%pressure 00681 00682 dHdot_dH = 0d0 00683 do i_part = 1,(n_eqn - 1) 00684 inputs%D = state(i_part) 00685 inputs%D_dry = condense_saved_D_dry(i_part) 00686 inputs%V_comp = condense_saved_V_comp_initial & 00687 * env_state%temp / condense_saved_env_state_initial%temp & 00688 / condense_saved_weight(i_part) 00689 inputs%kappa = condense_saved_kappa(i_part) 00690 call condense_rates(inputs, outputs) 00691 dDdot_dD(i_part) = outputs%dDdot_dD 00692 dDdot_dH(i_part) = outputs%dDdot_dH 00693 dHdot_dD(i_part) = outputs%dHdoti_dD + outputs%dHdotenv_dD 00694 dHdot_dH = dHdot_dH + outputs%dHdoti_dH 00695 end do 00696 dHdot_dH = dHdot_dH + outputs%dHdotenv_dH 00697 00698 call env_state_deallocate(env_state) 00699 00700 end subroutine condense_jac 00701 #endif 00702 00703 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00704 00705 #ifdef PMC_USE_SUNDIALS 00706 !> Solve the system \f$ Pz = r \f$ where \f$ P = I - \gamma J \f$ 00707 !> and \f$ J = \partial f / \partial y \f$. The solution is returned 00708 !> in the \f$ r \f$ vector. 00709 subroutine condense_jac_solve_f(n_eqn, time, state_p, state_dot_p, & 00710 rhs_p, gamma) bind(c) 00711 00712 !> Length of state vector. 00713 integer(kind=c_int), value, intent(in) :: n_eqn 00714 !> Current time (s). 00715 real(kind=c_double), value, intent(in) :: time 00716 !> Pointer to current state vector. 00717 type(c_ptr), value, intent(in) :: state_p 00718 !> Pointer to current state derivative vector. 00719 type(c_ptr), value, intent(in) :: state_dot_p 00720 !> Pointer to right-hand-side vector. 00721 type(c_ptr), value, intent(in) :: rhs_p 00722 !> Value of \c gamma scalar parameter. 00723 real(kind=c_double), value, intent(in) :: gamma 00724 00725 real(kind=c_double), pointer :: state(:), state_dot(:), rhs(:) 00726 real(kind=c_double) :: soln(n_eqn) 00727 real(kind=dp) :: dDdot_dD(n_eqn - 1), dDdot_dH(n_eqn - 1) 00728 real(kind=dp) :: dHdot_dD(n_eqn - 1), dHdot_dH 00729 real(kind=dp) :: lhs_n, rhs_n 00730 real(kind=c_double) :: residual(n_eqn) 00731 real(kind=dp) :: rhs_norm, soln_norm, residual_norm 00732 integer :: i_part 00733 00734 condense_count_solve = condense_count_solve + 1 00735 00736 call condense_jac(n_eqn, time, state_p, dDdot_dD, dDdot_dH, & 00737 dHdot_dD, dHdot_dH) 00738 00739 call c_f_pointer(state_p, state, (/ n_eqn /)) 00740 call c_f_pointer(state_dot_p, state_dot, (/ n_eqn /)) 00741 call c_f_pointer(rhs_p, rhs, (/ n_eqn /)) 00742 00743 !FIXME: write this all in matrix-vector notation, no i_part looping 00744 lhs_n = 1d0 - gamma * dHdot_dH 00745 rhs_n = rhs(n_eqn) 00746 do i_part = 1,(n_eqn - 1) 00747 lhs_n = lhs_n - (- gamma * dDdot_dH(i_part)) & 00748 * (- gamma * dHdot_dD(i_part)) / (1d0 - gamma * dDdot_dD(i_part)) 00749 rhs_n = rhs_n - (- gamma * dHdot_dD(i_part)) * rhs(i_part) & 00750 / (1d0 - gamma * dDdot_dD(i_part)) 00751 end do 00752 soln(n_eqn) = rhs_n / lhs_n 00753 00754 do i_part = 1,(n_eqn - 1) 00755 soln(i_part) = (rhs(i_part) & 00756 - (- gamma * dDdot_dH(i_part)) * soln(n_eqn)) & 00757 / (1d0 - gamma * dDdot_dD(i_part)) 00758 end do 00759 00760 if (CONDENSE_DO_TEST_JAC_SOLVE) then 00761 ! (I - g J) soln = rhs 00762 00763 ! residual = J soln 00764 residual(n_eqn) = sum(dHdot_dD * soln(1:(n_eqn-1))) & 00765 + dHdot_dH * soln(n_eqn) 00766 residual(1:(n_eqn-1)) = dDdot_dD * soln(1:(n_eqn-1)) & 00767 + dDdot_dH * soln(n_eqn) 00768 00769 residual = rhs - (soln - gamma * residual) 00770 rhs_norm = sqrt(sum(rhs**2)) 00771 soln_norm = sqrt(sum(soln**2)) 00772 residual_norm = sqrt(sum(residual**2)) 00773 write(0,*) 'rhs, soln, residual, residual/rhs = ', & 00774 rhs_norm, soln_norm, residual_norm, residual_norm / rhs_norm 00775 end if 00776 00777 rhs = soln 00778 00779 end subroutine condense_jac_solve_f 00780 #endif 00781 00782 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00783 00784 !> Determine the water equilibrium state of a single particle. 00785 subroutine condense_equilib_particle(env_state, aero_data, & 00786 aero_particle) 00787 00788 !> Environment state. 00789 type(env_state_t), intent(in) :: env_state 00790 !> Aerosol data. 00791 type(aero_data_t), intent(in) :: aero_data 00792 !> Particle. 00793 type(aero_particle_t), intent(inout) :: aero_particle 00794 00795 real(kind=dp) :: X, kappa, D_dry, D, g, dg_dD, a_w, daw_dD 00796 integer :: newton_step 00797 00798 X = 4d0 * const%water_molec_weight * const%water_surf_eng & 00799 / (const%univ_gas_const * env_state%temp & 00800 * const%water_density) 00801 kappa = aero_particle_solute_kappa(aero_particle, aero_data) 00802 D_dry = vol2diam(aero_particle_solute_volume(aero_particle, aero_data)) 00803 00804 D = D_dry 00805 g = 0d0 00806 dg_dD = 1d0 00807 do newton_step = 1,20 00808 D = D - g / dg_dD 00809 a_w = (D**3 - D_dry**3) / (D**3 + (kappa - 1d0) * D_dry**3) 00810 daw_dD = 3d0 * D**2 * kappa * D_dry**3 & 00811 / (D**3 + (kappa - 1d0) * D_dry**3)**2 00812 g = env_state%rel_humid - a_w * exp(X / D) 00813 dg_dD = - daw_dD * exp(X / D) + a_w * exp(X / D) * (X / D**2) 00814 end do 00815 call warn_assert_msg(426620001, abs(g) < 1d3 * epsilon(1d0), & 00816 "convergence problem in equilibriation") 00817 00818 aero_particle%vol(aero_data%i_water) = diam2vol(D) - diam2vol(D_dry) 00819 00820 end subroutine condense_equilib_particle 00821 00822 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00823 00824 !> Call condense_equilib_particle() on each particle in the aerosol 00825 !> to ensure that every particle has its water content in 00826 !> equilibrium. 00827 subroutine condense_equilib_particles(bin_grid, env_state, aero_data, & 00828 aero_weight, aero_state) 00829 00830 !> Bin grid. 00831 type(bin_grid_t), intent(in) :: bin_grid 00832 !> Environment state. 00833 type(env_state_t), intent(inout) :: env_state 00834 !> Aerosol data. 00835 type(aero_data_t), intent(in) :: aero_data 00836 !> Aerosol weight. 00837 type(aero_weight_t), intent(in) :: aero_weight 00838 !> Aerosol state. 00839 type(aero_state_t), intent(inout) :: aero_state 00840 00841 integer :: i_bin, i_part 00842 type(aero_particle_t), pointer :: aero_particle 00843 real(kind=dp) :: old_weight, new_weight 00844 type(aero_particle_t) :: new_aero_particle 00845 type(aero_info_t) :: aero_info 00846 integer :: n_copies, i_dup 00847 00848 call aero_particle_allocate(new_aero_particle) 00849 call aero_info_allocate(aero_info) 00850 do i_bin = 1,bin_grid%n_bin 00851 ! work backwards so any additions and removals will only affect 00852 ! particles that we've already dealt with 00853 do i_part = aero_state%bin(i_bin)%n_part,1,-1 00854 aero_particle => aero_state%bin(i_bin)%particle(i_part) 00855 old_weight = aero_weight_value(aero_weight, & 00856 aero_particle_radius(aero_particle)) 00857 00858 ! equilibriate the particle by adjusting its water content 00859 call condense_equilib_particle(env_state, aero_data, & 00860 aero_state%bin(i_bin)%particle(i_part)) 00861 00862 ! adjust particle number to account for weight changes 00863 if (aero_weight%type /= AERO_WEIGHT_TYPE_NONE) then 00864 new_weight = aero_weight_value(aero_weight, & 00865 aero_particle_radius(aero_particle)) 00866 n_copies = prob_round(old_weight / new_weight) 00867 if (n_copies == 0) then 00868 aero_info%id = aero_particle%id 00869 aero_info%action = AERO_INFO_WEIGHT 00870 aero_info%other_id = 0 00871 call aero_state_remove_particle_with_info(aero_state, & 00872 i_bin, i_part, aero_info) 00873 elseif (n_copies > 1) then 00874 do i_dup = 1,(n_copies - 1) 00875 call aero_particle_copy(aero_particle, new_aero_particle) 00876 call aero_particle_new_id(new_aero_particle) 00877 ! this might be adding into the wrong bin, but 00878 ! that's necessary as we might not have processed 00879 ! the correct bin yet. 00880 call aero_state_add_particle(aero_state, i_bin, & 00881 new_aero_particle) 00882 ! re-get the particle pointer, which may have 00883 ! changed due to reallocations caused by adding 00884 aero_particle => aero_state%bin(i_bin)%particle(i_part) 00885 end do 00886 end if 00887 end if 00888 end do 00889 end do 00890 ! We've modified particle diameters, so we need to update which 00891 ! bins they are in. 00892 call aero_state_resort(bin_grid, aero_state) 00893 00894 call aero_particle_deallocate(new_aero_particle) 00895 call aero_info_deallocate(aero_info) 00896 00897 end subroutine condense_equilib_particles 00898 00899 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00900 00901 end module pmc_condense