PartMC  2.3.0
condense.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2012 Nicole Riemer and Matthew West
2 ! Copyright (C) 2009 Joseph Ching
3 ! Licensed under the GNU General Public License version 2 or (at your
4 ! option) any later version. See the file COPYING for details.
5 
6 !> \file
7 !> The pmc_condense module.
8 
9 !> Water condensation onto aerosol particles.
10 !!
11 !! The model here assumes that the temperature \f$ T \f$ and pressure
12 !! \f$ p \f$ are prescribed as functions of time, while water content
13 !! per-particle and relative humidity are to be calculated by
14 !! integrating their rates of change.
15 !!
16 !! The state of the system is defined by the per-particle wet
17 !! diameters \f$ D_i \f$ and the relative humidity \f$ H \f$. The
18 !! state vector stores these in the order \f$ (D_1,\ldots,D_n,H)
19 !! \f$. The time-derivative of the state vector and the Jacobian
20 !! (derivative of the time-derivative with repsect to the state) all
21 !! conform to this ordering.
22 !!
23 !! The SUNDIALS ODE solver is used to compute the system evolution
24 !! using an implicit method. The system Jacobian is explicitly
25 !! inverveted as its structure is very simple.
26 !!
27 !! All equations used in this file are written in detail in the file
28 !! \c doc/condense.tex.
30 
31  use pmc_aero_state
32  use pmc_env_state
33  use pmc_aero_data
34  use pmc_util
36  use pmc_constants
37 #ifdef PMC_USE_SUNDIALS
38  use iso_c_binding
39 #endif
40 
41  !> Whether to numerically test the Jacobian-solve function during
42  !> execution (for debugging only).
43  logical, parameter :: CONDENSE_DO_TEST_JAC_SOLVE = .false.
44  !> Whether to print call-counts for helper routines during execution
45  !> (for debugging only).
46  logical, parameter :: CONDENSE_DO_TEST_COUNTS = .false.
47 
48  !> Result code indicating successful completion.
49  integer, parameter :: PMC_CONDENSE_SOLVER_SUCCESS = 0
50  !> Result code indicating failure to allocate \c y vector.
51  integer, parameter :: PMC_CONDENSE_SOLVER_INIT_Y = 1
52  !> Result code indicating failure to allocate \c abstol vector.
53  integer, parameter :: PMC_CONDENSE_SOLVER_INIT_ABSTOL = 2
54  !> Result code indicating failure to create the solver.
55  integer, parameter :: PMC_CONDENSE_SOLVER_INIT_CVODE_MEM = 3
56  !> Result code indicating failure to initialize the solver.
57  integer, parameter :: PMC_CONDENSE_SOLVER_INIT_CVODE = 4
58  !> Result code indicating failure to set tolerances.
59  integer, parameter :: PMC_CONDENSE_SOLVER_SVTOL = 5
60  !> Result code indicating failure to set maximum steps.
61  integer, parameter :: PMC_CONDENSE_SOLVER_SET_MAX_STEPS = 6
62  !> Result code indicating failure of the solver.
63  integer, parameter :: PMC_CONDENSE_SOLVER_FAIL = 7
64 
65  !> Internal-use structure for storing the inputs for the
66  !> rate-calculation function.
68  !> Temperature (K).
69  real(kind=dp) :: T
70  !> Rate of change of temperature (K s^{-1}).
71  real(kind=dp) :: Tdot
72  !> Relative humidity (1).
73  real(kind=dp) :: H
74  !> Pressure (Pa).
75  real(kind=dp) :: p
76  !> Rate of change of pressure (Pa s^{-1}).
77  real(kind=dp) :: pdot
78  !> Computational volume (m^3).
79  real(kind=dp) :: V_comp
80  !> Particle diameter (m).
81  real(kind=dp) :: D
82  !> Particle dry diameter (m).
83  real(kind=dp) :: D_dry
84  !> Kappa parameter (1).
85  real(kind=dp) :: kappa
87 
88  !> Internal-use structure for storing the outputs from the
89  !> rate-calculation function.
91  !> Change rate of diameter (m s^{-1}).
92  real(kind=dp) :: Ddot
93  !> Change rate of relative humidity due to this particle (s^{-1}).
94  real(kind=dp) :: Hdot_i
95  !> Change rate of relative humidity due to environment changes (s^{-1}).
96  real(kind=dp) :: Hdot_env
97  !> Sensitivity of \c Ddot to input \c D (m s^{-1} m^{-1}).
98  real(kind=dp) :: dDdot_dD
99  !> Sensitivity of \c Ddot to input \c H (m s^{-1}).
100  real(kind=dp) :: dDdot_dH
101  !> Sensitivity of \c Hdot_i to input \c D (s^{-1} m^{-1}).
102  real(kind=dp) :: dHdoti_dD
103  !> Sensitivity of \c Hdot_i to input \c D (s^{-1}).
104  real(kind=dp) :: dHdoti_dH
105  !> Sensitivity of \c Hdot_env to input \c D (s^{-1} m^{-1}).
106  real(kind=dp) :: dHdotenv_dD
107  !> Sensitivity of \c Hdot_env to input \c D (s^{-1}).
108  real(kind=dp) :: dHdotenv_dH
109  end type condense_rates_outputs_t
110 
111  !> Internal-use variable for storing the aerosol data during calls
112  !> to the ODE solver.
113  type(aero_data_t) :: condense_saved_aero_data
114  !> Internal-use variable for storing the initial environment state
115  !> during calls to the ODE solver.
116  type(env_state_t) :: condense_saved_env_state_initial
117  !> Internal-use variable for storing the rate of change of the
118  !> temperature during calls to the ODE solver.
119  real(kind=dp) :: condense_saved_Tdot
120  !> Internal-use variable for storing the rate of change of the
121  !> pressure during calls to the ODE solver.
122  real(kind=dp) :: condense_saved_pdot
123  !> Internal-use variable for storing the per-particle kappa values
124  !> during calls to the ODE solver.
125  real(kind=dp), allocatable :: condense_saved_kappa(:)
126  !> Internal-use variable for storing the per-particle dry diameters
127  !> during calls to the ODE solver.
128  real(kind=dp), allocatable :: condense_saved_D_dry(:)
129  !> Internal-use variable for storing the per-particle number
130  !> concentrations during calls to the ODE solver.
131  real(kind=dp), allocatable :: condense_saved_num_conc(:)
132 
133  !> Internal-use variable for counting calls to the vector field
134  !> subroutine.
135  integer, save :: condense_count_vf
136  !> Internal-use variable for counting calls to the Jacobian-solving
137  !> subroutine.
138  integer, save :: condense_count_solve
139 
140 contains
141 
142 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
143 
144  !> Do condensation to all the particles for a given time interval,
145  !> including updating the environment to account for the lost
146  !> water vapor.
147  subroutine condense_particles(aero_state, aero_data, env_state_initial, &
148  env_state_final, del_t)
149 
150  !> Aerosol state.
151  type(aero_state_t), intent(inout) :: aero_state
152  !> Aerosol data.
153  type(aero_data_t), intent(in) :: aero_data
154  !> Environment state at the start of the timestep.
155  type(env_state_t), intent(in) :: env_state_initial
156  !> Environment state at the end of the timestep. The rel_humid
157  !> value will be ignored and overwritten with the
158  !> condensation-computed value.
159  type(env_state_t), intent(inout) :: env_state_final
160  !> Total time to integrate.
161  real(kind=dp), intent(in) :: del_t
162 
163  integer :: i_part, n_eqn, i_eqn
164  type(aero_particle_t), pointer :: aero_particle
165  real(kind=dp) :: state(aero_state%apa%n_part + 1), init_time, final_time
166  real(kind=dp) :: abs_tol_vector(aero_state%apa%n_part + 1)
167  real(kind=dp) :: num_conc
168  real(kind=dp) :: reweight_num_conc(aero_state%apa%n_part)
169  real(kind=dp) :: water_vol_conc_initial, water_vol_conc_final
170  real(kind=dp) :: vapor_vol_conc_initial, vapor_vol_conc_final
171  real(kind=dp) :: d_water_vol_conc, d_vapor_vol_conc
172  real(kind=dp) :: v_comp_ratio, water_rel_error
173 #ifdef PMC_USE_SUNDIALS
174  real(kind=c_double), target :: state_f(aero_state%apa%n_part + 1)
175  real(kind=c_double), target :: abstol_f(aero_state%apa%n_part + 1)
176  type(c_ptr) :: state_f_p, abstol_f_p
177  integer(kind=c_int) :: n_eqn_f, solver_stat
178  real(kind=c_double) :: reltol_f, t_initial_f, t_final_f
179 #endif
180 
181 #ifdef PMC_USE_SUNDIALS
182 #ifndef DOXYGEN_SKIP_DOC
183  interface
184  integer(kind=c_int) function condense_solver(neq, x_f, abstol_f, &
185  reltol_f, t_initial_f, t_final_f) bind(c)
186  use iso_c_binding
187  integer(kind=c_int), value :: neq
188  type(c_ptr), value :: x_f
189  type(c_ptr), value :: abstol_f
190  real(kind=c_double), value :: reltol_f
191  real(kind=c_double), value :: t_initial_f
192  real(kind=c_double), value :: t_final_f
193  end function condense_solver
194  end interface
195 #endif
196 #endif
197 
198  ! initial water concentration in the aerosol particles
199  water_vol_conc_initial = 0d0
200  do i_part = 1,aero_state%apa%n_part
201  aero_particle => aero_state%apa%particle(i_part)
202  num_conc = aero_weight_array_num_conc(aero_state%awa, aero_particle)
203  water_vol_conc_initial = water_vol_conc_initial &
204  + aero_particle%vol(aero_data%i_water) * num_conc
205  end do
206 
207  ! save data for use within the timestepper
208  call aero_data_allocate(condense_saved_aero_data)
209  call env_state_allocate(condense_saved_env_state_initial)
210  call aero_data_copy(aero_data, condense_saved_aero_data)
211  call env_state_copy(env_state_initial, condense_saved_env_state_initial)
212  condense_saved_tdot &
213  = (env_state_final%temp - env_state_initial%temp) / del_t
214  condense_saved_pdot &
215  = (env_state_final%pressure - env_state_initial%pressure) / del_t
216 
217  ! construct initial state vector from aero_state and env_state
218  allocate(condense_saved_kappa(aero_state%apa%n_part))
219  allocate(condense_saved_d_dry(aero_state%apa%n_part))
220  allocate(condense_saved_num_conc(aero_state%apa%n_part))
221  ! work backwards for consistency with the later number
222  ! concentration adjustment, which has specific ordering
223  ! requirements
224  do i_part = aero_state%apa%n_part,1,-1
225  aero_particle => aero_state%apa%particle(i_part)
226  condense_saved_kappa(i_part) &
227  = aero_particle_solute_kappa(aero_particle, aero_data)
228  condense_saved_d_dry(i_part) = vol2diam(&
229  aero_particle_solute_volume(aero_particle, aero_data))
230  condense_saved_num_conc(i_part) &
231  = aero_weight_array_num_conc(aero_state%awa, aero_particle)
232  state(i_part) = aero_particle_diameter(aero_particle)
233  abs_tol_vector(i_part) = max(1d-30, &
234  1d-8 * (state(i_part) - condense_saved_d_dry(i_part)))
235  end do
236  state(aero_state%apa%n_part + 1) = env_state_initial%rel_humid
237  abs_tol_vector(aero_state%apa%n_part + 1) = 1d-10
238 
239 #ifdef PMC_USE_SUNDIALS
240  ! call SUNDIALS solver
241  n_eqn = aero_state%apa%n_part + 1
242  n_eqn_f = int(n_eqn, kind=c_int)
243  reltol_f = real(1d-8, kind=c_double)
244  t_initial_f = real(0, kind=c_double)
245  t_final_f = real(del_t, kind=c_double)
246  do i_eqn = 1,n_eqn
247  state_f(i_eqn) = real(state(i_eqn), kind=c_double)
248  abstol_f(i_eqn) = real(abs_tol_vector(i_eqn), kind=c_double)
249  end do
250  state_f_p = c_loc(state_f)
251  abstol_f_p = c_loc(abstol_f)
252  condense_count_vf = 0
253  condense_count_solve = 0
254  solver_stat = condense_solver(n_eqn_f, state_f_p, abstol_f_p, reltol_f, &
255  t_initial_f, t_final_f)
256  call condense_check_solve(solver_stat)
257  if (condense_do_test_counts) then
258  write(0,*) 'condense_count_vf ', condense_count_vf
259  write(0,*) 'condense_count_solve ', condense_count_solve
260  end if
261  do i_eqn = 1,n_eqn
262  state(i_eqn) = real(state_f(i_eqn), kind=dp)
263  end do
264 #endif
265 
266  ! unpack result state vector into env_state_final
267  env_state_final%rel_humid = state(aero_state%apa%n_part + 1)
268 
269  ! unpack result state vector into aero_state, compute the final
270  ! water volume concentration, and adjust particle number to
271  ! account for number concentration changes
272  water_vol_conc_final = 0d0
273  call aero_state_num_conc_for_reweight(aero_state, reweight_num_conc)
274  do i_part = 1,aero_state%apa%n_part
275  aero_particle => aero_state%apa%particle(i_part)
276  num_conc = aero_weight_array_num_conc(aero_state%awa, aero_particle)
277 
278  ! translate output back to particle
279  aero_particle%vol(aero_data%i_water) = diam2vol(state(i_part)) &
280  - aero_particle_solute_volume(aero_particle, aero_data)
281 
282  ! ensure volumes stay positive
283  aero_particle%vol(aero_data%i_water) = max(0d0, &
284  aero_particle%vol(aero_data%i_water))
285 
286  ! add up total water volume, using old number concentrations
287  water_vol_conc_final = water_vol_conc_final &
288  + aero_particle%vol(aero_data%i_water) * num_conc
289  end do
290  ! adjust particles to account for weight changes
291  call aero_state_reweight(aero_state, reweight_num_conc)
292 
293  ! Check that water removed from particles equals water added to
294  ! vapor. Note that water concentration is not conserved (due to
295  ! computational volume changes), and we need to consider particle
296  ! weightings correctly.
297  v_comp_ratio = env_state_final%temp * env_state_initial%pressure &
298  / (env_state_initial%temp * env_state_final%pressure)
299  vapor_vol_conc_initial = aero_data%molec_weight(aero_data%i_water) &
300  / (const%univ_gas_const * env_state_initial%temp) &
301  * env_state_sat_vapor_pressure(env_state_initial) &
302  * env_state_initial%rel_humid &
303  / aero_particle_water_density(aero_data)
304  vapor_vol_conc_final = aero_data%molec_weight(aero_data%i_water) &
305  / (const%univ_gas_const * env_state_final%temp) &
306  * env_state_sat_vapor_pressure(env_state_final) &
307  * env_state_final%rel_humid &
308  * v_comp_ratio / aero_particle_water_density(aero_data)
309  d_vapor_vol_conc = vapor_vol_conc_final - vapor_vol_conc_initial
310  d_water_vol_conc = water_vol_conc_final - water_vol_conc_initial
311  water_rel_error = (d_vapor_vol_conc + d_water_vol_conc) &
312  / (vapor_vol_conc_final + water_vol_conc_final)
313  call warn_assert_msg(477865387, abs(water_rel_error) < 1d-6, &
314  "condensation water imbalance too high: " &
315  // trim(real_to_string(water_rel_error)))
316 
317  deallocate(condense_saved_kappa)
318  deallocate(condense_saved_d_dry)
319  deallocate(condense_saved_num_conc)
320  call aero_data_deallocate(condense_saved_aero_data)
321  call env_state_deallocate(condense_saved_env_state_initial)
322 
323  end subroutine condense_particles
324 
325 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
326 
327 #ifdef PMC_USE_SUNDIALS
328  !> Check the return code from the condense_solver() function.
329  subroutine condense_check_solve(value)
330 
331  !> Return code to check.
332  integer(kind=c_int), intent(in) :: value
333 
334  if (value == pmc_condense_solver_success) then
335  return
336  elseif (value == pmc_condense_solver_init_y) then
337  call die_msg(123749472, "condense_solver: " &
338  // "failed to allocate y vector")
339  elseif (value == pmc_condense_solver_init_abstol) then
340  call die_msg(563665949, "condense_solver: " &
341  // "failed to allocate abstol vector")
342  elseif (value == pmc_condense_solver_init_cvode_mem) then
343  call die_msg(700541443, "condense_solver: " &
344  // "failed to create the solver")
345  elseif (value == pmc_condense_solver_init_cvode) then
346  call die_msg(297559183, "condense_solver: " &
347  // "failure to initialize the solver")
348  elseif (value == pmc_condense_solver_svtol) then
349  call die_msg(848342417, "condense_solver: " &
350  // "failed to set tolerances")
351  elseif (value == pmc_condense_solver_set_max_steps) then
352  call die_msg(275591501, "condense_solver: " &
353  // "failed to set maximum steps")
354  elseif (value == pmc_condense_solver_fail) then
355  call die_msg(862254233, "condense_solver: solver failed")
356  else
357  call die_msg(635697577, "condense_solver: unknown return code: " &
358  // trim(integer_to_string(value)))
359  end if
360 
361  end subroutine condense_check_solve
362 #endif
363 
364 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
365 
366  !> Compute the rate of change of particle diameter and relative
367  !> humidity for a single particle, together with the derivatives of
368  !> the rates with respect to the input variables.
369  subroutine condense_rates(inputs, outputs)
370 
371  !> Inputs to rates.
372  type(condense_rates_inputs_t), intent(in) :: inputs
373  !> Outputs rates.
374  type(condense_rates_outputs_t), intent(out) :: outputs
375 
376  real(kind=dp) :: rho_w, m_w, p_0, dp0_dt_div_p0, rho_air, k_a, d_v, u
377  real(kind=dp) :: v, w, x, y, z, k_ap, dkap_dd, d_vp, ddvp_dd
378  real(kind=dp) :: a_w, daw_dd, delta_star, h, dh_ddelta, dh_dd
379  real(kind=dp) :: dh_dh, ddeltastar_dd, ddeltastar_dh
380  integer :: newton_step
381 
382  rho_w = const%water_density
383  m_w = const%water_molec_weight
384  p_0 = const%water_eq_vap_press &
385  * 10d0**(7.45d0 * (inputs%T - const%water_freeze_temp) &
386  / (inputs%T - 38d0))
387  dp0_dt_div_p0 = 7.45d0 * log(10d0) * (const%water_freeze_temp - 38d0) &
388  / (inputs%T - 38d0)**2
389  rho_air = const%air_molec_weight * inputs%p &
390  / (const%univ_gas_const * inputs%T)
391 
392  k_a = 1d-3 * (4.39d0 + 0.071d0 * inputs%T)
393  d_v = 0.211d-4 / (inputs%p / const%air_std_press) &
394  * (inputs%T / 273d0)**1.94d0
395  u = const%water_latent_heat * rho_w / (4d0 * inputs%T)
396  v = 4d0 * m_w * p_0 / (rho_w * const%univ_gas_const * inputs%T)
397  w = const%water_latent_heat * m_w / (const%univ_gas_const * inputs%T)
398  x = 4d0 * m_w * const%water_surf_eng &
399  / (const%univ_gas_const * inputs%T * rho_w)
400  y = 2d0 * k_a / (const%accom_coeff * rho_air &
401  * const%air_spec_heat) &
402  * sqrt(2d0 * const%pi * const%air_molec_weight &
403  / (const%univ_gas_const * inputs%T))
404  z = 2d0 * d_v / const%accom_coeff * sqrt(2d0 * const%pi * m_w &
405  / (const%univ_gas_const * inputs%T))
406 
407  outputs%Hdot_env = - dp0_dt_div_p0 * inputs%Tdot * inputs%H &
408  + inputs%H * inputs%pdot / inputs%p
409  outputs%dHdotenv_dD = 0d0
410  outputs%dHdotenv_dH = - dp0_dt_div_p0 * inputs%Tdot &
411  + inputs%pdot / inputs%p
412 
413  if (inputs%D <= inputs%D_dry) then
414  k_ap = k_a / (1d0 + y / inputs%D_dry)
415  dkap_dd = 0d0
416  d_vp = d_v / (1d0 + z / inputs%D_dry)
417  ddvp_dd = 0d0
418  a_w = 0d0
419  daw_dd = 0d0
420 
421  delta_star = u * v * d_vp * inputs%H / k_ap
422 
423  outputs%Ddot = k_ap * delta_star / (u * inputs%D_dry)
424  outputs%Hdot_i = - 2d0 * const%pi / (v * inputs%V_comp) &
425  * inputs%D_dry**2 * outputs%Ddot
426 
427  dh_ddelta = k_ap
428  dh_dd = 0d0
429  dh_dh = - u * v * d_vp
430 
431  ddeltastar_dd = - dh_dd / dh_ddelta
432  ddeltastar_dh = - dh_dh / dh_ddelta
433 
434  outputs%dDdot_dD = 0d0
435  outputs%dDdot_dH = k_ap / (u * inputs%D_dry) * ddeltastar_dh
436  outputs%dHdoti_dD = - 2d0 * const%pi / (v * inputs%V_comp) &
437  * inputs%D_dry**2 * outputs%dDdot_dD
438  outputs%dHdoti_dH = - 2d0 * const%pi / (v * inputs%V_comp) &
439  * inputs%D_dry**2 * outputs%dDdot_dH
440 
441  return
442  end if
443 
444  k_ap = k_a / (1d0 + y / inputs%D)
445  dkap_dd = k_a * y / (inputs%D + y)**2
446  d_vp = d_v / (1d0 + z / inputs%D)
447  ddvp_dd = d_v * z / (inputs%D + z)**2
448  a_w = (inputs%D**3 - inputs%D_dry**3) &
449  / (inputs%D**3 + (inputs%kappa - 1d0) * inputs%D_dry**3)
450  daw_dd = 3d0 * inputs%D**2 * inputs%kappa * inputs%D_dry**3 &
451  / (inputs%D**3 + (inputs%kappa - 1d0) * inputs%D_dry**3)**2
452 
453  delta_star = 0d0
454  h = 0d0
455  dh_ddelta = 1d0
456  do newton_step = 1,5
457  ! update delta_star first so when the newton loop ends we have
458  ! h and dh_ddelta evaluated at the final delta_star value
459  delta_star = delta_star - h / dh_ddelta
460  h = k_ap * delta_star - u * v * d_vp &
461  * (inputs%H - a_w / (1d0 + delta_star) &
462  * exp(w * delta_star / (1d0 + delta_star) &
463  + (x / inputs%D) / (1d0 + delta_star)))
464  dh_ddelta = &
465  k_ap - u * v * d_vp * a_w / (1d0 + delta_star)**2 &
466  * (1d0 - w / (1d0 + delta_star) &
467  + (x / inputs%D) / (1d0 + delta_star)) &
468  * exp(w * delta_star / (1d0 + delta_star) &
469  + (x / inputs%D) / (1d0 + delta_star))
470  end do
471  call warn_assert_msg(387362320, &
472  abs(h) < 1d3 * epsilon(1d0) * abs(u * v * d_vp * inputs%H), &
473  "condensation newton loop did not satisfy convergence tolerance")
474 
475  outputs%Ddot = k_ap * delta_star / (u * inputs%D)
476  outputs%Hdot_i = - 2d0 * const%pi / (v * inputs%V_comp) &
477  * inputs%D**2 * outputs%Ddot
478 
479  dh_dd = dkap_dd * delta_star &
480  - u * v * ddvp_dd * inputs%H + u * v &
481  * (a_w * ddvp_dd + d_vp * daw_dd &
482  - d_vp * a_w * (x / inputs%D**2) / (1d0 + delta_star)) &
483  * (1d0 / (1d0 + delta_star)) &
484  * exp((w * delta_star) / (1d0 + delta_star) &
485  + (x / inputs%D) / (1d0 + delta_star))
486  dh_dh = - u * v * d_vp
487 
488  ddeltastar_dd = - dh_dd / dh_ddelta
489  ddeltastar_dh = - dh_dh / dh_ddelta
490 
491  outputs%dDdot_dD = dkap_dd * delta_star / (u * inputs%D) &
492  + k_ap * ddeltastar_dd / (u * inputs%D) &
493  - k_ap * delta_star / (u * inputs%D**2)
494  outputs%dDdot_dH = k_ap / (u * inputs%D) * ddeltastar_dh
495  outputs%dHdoti_dD = - 2d0 * const%pi / (v * inputs%V_comp) &
496  * (2d0 * inputs%D * outputs%Ddot + inputs%D**2 * outputs%dDdot_dD)
497  outputs%dHdoti_dH = - 2d0 * const%pi / (v * inputs%V_comp) &
498  * inputs%D**2 * outputs%dDdot_dH
499 
500  end subroutine condense_rates
501 
502 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
503 
504 #ifdef PMC_USE_SUNDIALS
505  !> Compute the condensation rates (Ddot and Hdot) at the current
506  !> value of the state (D and H).
507  subroutine condense_vf_f(n_eqn, time, state_p, state_dot_p) bind(c)
508 
509  !> Length of state vector.
510  integer(kind=c_int), value, intent(in) :: n_eqn
511  !> Current time (s).
512  real(kind=c_double), value, intent(in) :: time
513  !> Pointer to state data.
514  type(c_ptr), value, intent(in) :: state_p
515  !> Pointer to state_dot data.
516  type(c_ptr), value, intent(in) :: state_dot_p
517 
518  real(kind=c_double), pointer :: state(:)
519  real(kind=c_double), pointer :: state_dot(:)
520  real(kind=dp) :: hdot
521  integer :: i_part
522  type(condense_rates_inputs_t) :: inputs
523  type(condense_rates_outputs_t) :: outputs
524 
525  condense_count_vf = condense_count_vf + 1
526 
527  call c_f_pointer(state_p, state, (/ n_eqn /))
528  call c_f_pointer(state_dot_p, state_dot, (/ n_eqn /))
529 
530  inputs%T = condense_saved_env_state_initial%temp &
531  + time * condense_saved_tdot
532  inputs%p = condense_saved_env_state_initial%pressure &
533  + time * condense_saved_pdot
534  inputs%Tdot = condense_saved_tdot
535  inputs%pdot = condense_saved_pdot
536  inputs%H = state(n_eqn)
537 
538  hdot = 0d0
539  do i_part = 1,(n_eqn - 1)
540  inputs%D = state(i_part)
541  inputs%D_dry = condense_saved_d_dry(i_part)
542  inputs%V_comp = (inputs%T &
543  * condense_saved_env_state_initial%pressure) &
544  / (condense_saved_env_state_initial%temp * inputs%p) &
545  / condense_saved_num_conc(i_part)
546  inputs%kappa = condense_saved_kappa(i_part)
547  call condense_rates(inputs, outputs)
548  state_dot(i_part) = outputs%Ddot
549  hdot = hdot + outputs%Hdot_i
550  end do
551  hdot = hdot + outputs%Hdot_env
552 
553  state_dot(n_eqn) = hdot
554 
555  end subroutine condense_vf_f
556 #endif
557 
558 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
559 
560 #ifdef PMC_USE_SUNDIALS
561  !> Compute the Jacobian given by the derivatives of the condensation
562  !> rates (Ddot and Hdot) with respect to the input variables (D and
563  !> H).
564  subroutine condense_jac(n_eqn, time, state_p, dDdot_dD, dDdot_dH, &
565  dhdot_dd, dhdot_dh)
566 
567  !> Length of state vector.
568  integer(kind=c_int), intent(in) :: n_eqn
569  !> Current time (s).
570  real(kind=c_double), intent(in) :: time
571  !> Pointer to current state vector.
572  type(c_ptr), intent(in) :: state_p
573  !> Derivative of Ddot with respect to D.
574  real(kind=dp), intent(out) :: dddot_dd(n_eqn - 1)
575  !> Derivative of Ddot with respect to H.
576  real(kind=dp), intent(out) :: dddot_dh(n_eqn - 1)
577  !> Derivative of Hdot with respect to D.
578  real(kind=dp), intent(out) :: dhdot_dd(n_eqn - 1)
579  !> Derivative of Hdot with respect to H.
580  real(kind=dp), intent(out) :: dhdot_dh
581 
582  real(kind=c_double), pointer :: state(:)
583  integer :: i_part
584  type(condense_rates_inputs_t) :: inputs
585  type(condense_rates_outputs_t) :: outputs
586 
587  call c_f_pointer(state_p, state, (/ n_eqn /))
588 
589  inputs%T = condense_saved_env_state_initial%temp &
590  + time * condense_saved_tdot
591  inputs%p = condense_saved_env_state_initial%pressure &
592  + time * condense_saved_pdot
593  inputs%Tdot = condense_saved_tdot
594  inputs%pdot = condense_saved_pdot
595  inputs%H = state(n_eqn)
596 
597  dhdot_dh = 0d0
598  do i_part = 1,(n_eqn - 1)
599  inputs%D = state(i_part)
600  inputs%D_dry = condense_saved_d_dry(i_part)
601  inputs%V_comp = (inputs%T &
602  * condense_saved_env_state_initial%pressure) &
603  / (condense_saved_env_state_initial%temp * inputs%p) &
604  / condense_saved_num_conc(i_part)
605  inputs%kappa = condense_saved_kappa(i_part)
606  call condense_rates(inputs, outputs)
607  dddot_dd(i_part) = outputs%dDdot_dD
608  dddot_dh(i_part) = outputs%dDdot_dH
609  dhdot_dd(i_part) = outputs%dHdoti_dD + outputs%dHdotenv_dD
610  dhdot_dh = dhdot_dh + outputs%dHdoti_dH
611  end do
612  dhdot_dh = dhdot_dh + outputs%dHdotenv_dH
613 
614  end subroutine condense_jac
615 #endif
616 
617 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
618 
619 #ifdef PMC_USE_SUNDIALS
620  !> Solve the system \f$ Pz = r \f$ where \f$ P = I - \gamma J \f$
621  !> and \f$ J = \partial f / \partial y \f$. The solution is returned
622  !> in the \f$ r \f$ vector.
623  subroutine condense_jac_solve_f(n_eqn, time, state_p, state_dot_p, &
624  rhs_p, gamma) bind(c)
625 
626  !> Length of state vector.
627  integer(kind=c_int), value, intent(in) :: n_eqn
628  !> Current time (s).
629  real(kind=c_double), value, intent(in) :: time
630  !> Pointer to current state vector.
631  type(c_ptr), value, intent(in) :: state_p
632  !> Pointer to current state derivative vector.
633  type(c_ptr), value, intent(in) :: state_dot_p
634  !> Pointer to right-hand-side vector.
635  type(c_ptr), value, intent(in) :: rhs_p
636  !> Value of \c gamma scalar parameter.
637  real(kind=c_double), value, intent(in) :: gamma
638 
639  real(kind=c_double), pointer :: state(:), state_dot(:), rhs(:)
640  real(kind=c_double) :: soln(n_eqn)
641  real(kind=dp) :: dddot_dd(n_eqn - 1), dddot_dh(n_eqn - 1)
642  real(kind=dp) :: dhdot_dd(n_eqn - 1), dhdot_dh
643  real(kind=dp) :: lhs_n, rhs_n
644  real(kind=c_double) :: residual(n_eqn)
645  real(kind=dp) :: rhs_norm, soln_norm, residual_norm
646  integer :: i_part
647 
648  condense_count_solve = condense_count_solve + 1
649 
650  call condense_jac(n_eqn, time, state_p, dddot_dd, dddot_dh, &
651  dhdot_dd, dhdot_dh)
652 
653  call c_f_pointer(state_p, state, (/ n_eqn /))
654  call c_f_pointer(state_dot_p, state_dot, (/ n_eqn /))
655  call c_f_pointer(rhs_p, rhs, (/ n_eqn /))
656 
657  !FIXME: write this all in matrix-vector notation, no i_part looping
658  lhs_n = 1d0 - gamma * dhdot_dh
659  rhs_n = rhs(n_eqn)
660  do i_part = 1,(n_eqn - 1)
661  lhs_n = lhs_n - (- gamma * dddot_dh(i_part)) &
662  * (- gamma * dhdot_dd(i_part)) / (1d0 - gamma * dddot_dd(i_part))
663  rhs_n = rhs_n - (- gamma * dhdot_dd(i_part)) * rhs(i_part) &
664  / (1d0 - gamma * dddot_dd(i_part))
665  end do
666  soln(n_eqn) = rhs_n / lhs_n
667 
668  do i_part = 1,(n_eqn - 1)
669  soln(i_part) = (rhs(i_part) &
670  - (- gamma * dddot_dh(i_part)) * soln(n_eqn)) &
671  / (1d0 - gamma * dddot_dd(i_part))
672  end do
673 
674  if (condense_do_test_jac_solve) then
675  ! (I - g J) soln = rhs
676 
677  ! residual = J soln
678  residual(n_eqn) = sum(dhdot_dd * soln(1:(n_eqn-1))) &
679  + dhdot_dh * soln(n_eqn)
680  residual(1:(n_eqn-1)) = dddot_dd * soln(1:(n_eqn-1)) &
681  + dddot_dh * soln(n_eqn)
682 
683  residual = rhs - (soln - gamma * residual)
684  rhs_norm = sqrt(sum(rhs**2))
685  soln_norm = sqrt(sum(soln**2))
686  residual_norm = sqrt(sum(residual**2))
687  write(0,*) 'rhs, soln, residual, residual/rhs = ', &
688  rhs_norm, soln_norm, residual_norm, residual_norm / rhs_norm
689  end if
690 
691  rhs = soln
692 
693  end subroutine condense_jac_solve_f
694 #endif
695 
696 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
697 
698  !> Determine the water equilibrium state of a single particle.
699  subroutine condense_equilib_particle(env_state, aero_data, &
700  aero_particle)
701 
702  !> Environment state.
703  type(env_state_t), intent(in) :: env_state
704  !> Aerosol data.
705  type(aero_data_t), intent(in) :: aero_data
706  !> Particle.
707  type(aero_particle_t), intent(inout) :: aero_particle
708 
709  real(kind=dp) :: x, kappa, d_dry, d, g, dg_dd, a_w, daw_dd
710  integer :: newton_step
711 
712  x = 4d0 * const%water_molec_weight * const%water_surf_eng &
713  / (const%univ_gas_const * env_state%temp &
714  * const%water_density)
715  kappa = aero_particle_solute_kappa(aero_particle, aero_data)
716  d_dry = vol2diam(aero_particle_solute_volume(aero_particle, aero_data))
717 
718  d = d_dry
719  g = 0d0
720  dg_dd = 1d0
721  do newton_step = 1,20
722  d = d - g / dg_dd
723  a_w = (d**3 - d_dry**3) / (d**3 + (kappa - 1d0) * d_dry**3)
724  daw_dd = 3d0 * d**2 * kappa * d_dry**3 &
725  / (d**3 + (kappa - 1d0) * d_dry**3)**2
726  g = env_state%rel_humid - a_w * exp(x / d)
727  dg_dd = - daw_dd * exp(x / d) + a_w * exp(x / d) * (x / d**2)
728  end do
729  call warn_assert_msg(426620001, abs(g) < 1d3 * epsilon(1d0), &
730  "convergence problem in equilibriation")
731 
732  aero_particle%vol(aero_data%i_water) = diam2vol(d) - diam2vol(d_dry)
733 
734  end subroutine condense_equilib_particle
735 
736 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
737 
738  !> Call condense_equilib_particle() on each particle in the aerosol
739  !> to ensure that every particle has its water content in
740  !> equilibrium.
741  subroutine condense_equilib_particles(env_state, aero_data, aero_state)
742 
743  !> Environment state.
744  type(env_state_t), intent(inout) :: env_state
745  !> Aerosol data.
746  type(aero_data_t), intent(in) :: aero_data
747  !> Aerosol state.
748  type(aero_state_t), intent(inout) :: aero_state
749 
750  integer :: i_part
751  real(kind=dp) :: reweight_num_conc(aero_state%apa%n_part)
752 
753  ! We're modifying particle diameters, so bin sorting is now invalid
754  aero_state%valid_sort = .false.
755 
756  call aero_state_num_conc_for_reweight(aero_state, reweight_num_conc)
757  do i_part = aero_state%apa%n_part,1,-1
758  call condense_equilib_particle(env_state, aero_data, &
759  aero_state%apa%particle(i_part))
760  end do
761  ! adjust particles to account for weight changes
762  call aero_state_reweight(aero_state, reweight_num_conc)
763 
764  end subroutine condense_equilib_particles
765 
766 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
767 
768 end module pmc_condense
real(kind=dp) elemental function vol2diam(v)
Convert volume (m^3) to diameter (m).
Definition: util.F90:250
The aero_data_t structure and associated subroutines.
Definition: aero_data.F90:9
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:133
real(kind=dp) function aero_particle_solute_kappa(aero_particle, aero_data)
Returns the average of the solute kappas (1).
real(kind=dp) elemental function diam2vol(d)
Convert diameter (m) to volume (m^3).
Definition: util.F90:298
Water condensation onto aerosol particles.
Definition: condense.F90:29
Physical constants.
Definition: constants.F90:9
The aero_particle_t structure and associated subroutines.
The env_state_t structure and associated subroutines.
Definition: env_state.F90:9
subroutine aero_state_num_conc_for_reweight(aero_state, reweight_num_conc)
Save the correct number concentrations for later use by aero_state_reweight().
Definition: aero_state.F90:522
elemental real(kind=dp) function aero_particle_diameter(aero_particle)
Total diameter of the particle (m).
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 it...
Definition: condense.F90:741
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_particle_water_density(aero_data)
Returns the water density (kg/m^3).
subroutine warn_assert_msg(code, condition_ok, warning_msg)
Prints a warning message if condition_ok is false.
Definition: util.F90:58
subroutine aero_data_deallocate(aero_data)
Frees all storage.
Definition: aero_data.F90:116
subroutine condense_equilib_particle(env_state, aero_data, aero_particle)
Determine the water equilibrium state of a single particle.
Definition: condense.F90:699
Common utility subroutines.
Definition: util.F90:9
Current environment state.
Definition: env_state.F90:26
subroutine env_state_deallocate(env_state)
Free all storage.
Definition: env_state.F90:78
int condense_solver(int neq, double *x_f, double *abstol_f, double reltol_f, double t_initial_f, double t_final_f)
Call the ODE solver.
The aero_state_t structure and assocated subroutines.
Definition: aero_state.F90:9
subroutine condense_rates(inputs, outputs)
Compute the rate of change of particle diameter and relative humidity for a single particle...
Definition: condense.F90:369
character(len=pmc_util_convert_string_len) function integer_to_string(val)
Convert an integer to a string format.
Definition: util.F90:743
subroutine aero_state_reweight(aero_state, reweight_num_conc)
Reweight all particles after their constituent volumes have been altered.
Definition: aero_state.F90:550
The current collection of aerosol particles.
Definition: aero_state.F90:63
Single aerosol particle data structure.
subroutine env_state_copy(env_from, env_to)
env_to = env_from
Definition: env_state.F90:138
Internal-use structure for storing the inputs for the rate-calculation function.
Definition: condense.F90:67
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...
Definition: condense.F90:147
subroutine aero_data_copy(aero_data_from, aero_data_to)
Copy structure.
Definition: aero_data.F90:134
real(kind=dp) function aero_particle_solute_volume(aero_particle, aero_data)
Returns the total solute volume (m^3).
Internal-use structure for storing the outputs from the rate-calculation function.
Definition: condense.F90:90
subroutine aero_data_allocate(aero_data)
Allocate storage for aero_data.
Definition: aero_data.F90:69
subroutine env_state_allocate(env_state)
Allocate an empty environment.
Definition: env_state.F90:56
real(kind=dp) function env_state_sat_vapor_pressure(env_state)
Computes the current saturation vapor pressure (Pa).
Definition: env_state.F90:192
character(len=pmc_util_convert_string_len) function real_to_string(val)
Convert a real to a string format.
Definition: util.F90:759
Aerosol material properties and associated data.
Definition: aero_data.F90:40