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