PartMC  2.8.0
scenario.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2016 Nicole Riemer and Matthew West
2 ! Licensed under the GNU General Public License version 2 or (at your
3 ! option) any later version. See the file COPYING for details.
4 
5 !> \file
6 !> The pmc_scenario module.
7 
8 !> The scenario_t structure and associated subroutines.
10 
11  use pmc_gas_state
12  use pmc_aero_dist
13  use pmc_util
14  use pmc_env_state
15  use pmc_aero_state
16  use pmc_spec_file
17  use pmc_aero_data
18  use pmc_gas_data
19  use pmc_chamber
20  use pmc_mpi
21 #ifdef PMC_USE_MPI
22  use mpi
23 #endif
24 
25  !> Type code for an undefined or invalid loss function.
26  integer, parameter :: scenario_loss_function_invalid = 0
27  !> Type code for a zero loss function.
28  integer, parameter :: scenario_loss_function_none = 1
29  !> Type code for a constant loss function.
30  integer, parameter :: scenario_loss_function_constant = 2
31  !> Type code for a loss rate function proportional to particle volume.
32  integer, parameter :: scenario_loss_function_volume = 3
33  !> Type code for a loss rate function based on dry deposition
34  integer, parameter :: scenario_loss_function_drydep = 4
35  !> Type code for a loss rate function for chamber experiments.
36  integer, parameter :: scenario_loss_function_chamber = 5
37 
38  !> Parameter to switch between algorithms for particle loss.
39  !! A value of 0 will always use the naive algorithm, and
40  !! a value of 1 will always use the accept-reject algorithm.
41  real(kind=dp), parameter :: scenario_loss_alg_threshold = 1.0d0
42 
43  !> Scenario data.
44  !!
45  !! This is everything needed to drive the scenario being simulated.
46  !!
47  !! The temperature, pressure, emissions and background states are profiles
48  !! prescribed as functions of time by giving a number of times and
49  !! the corresponding data. Simple data such as temperature and pressure is
50  !! linearly interpolated between times, with constant interpolation
51  !! outside of the range of times. Gases and aerosols are
52  !! interpolated with gas_state_interp_1d() and
53  !! aero_dist_interp_1d(), respectively.
55  !> Temperature set-point times (s).
56  real(kind=dp), allocatable :: temp_time(:)
57  !> Temperatures at set-points (K).
58  real(kind=dp), allocatable :: temp(:)
59 
60  !> Pressure set-point times (s).
61  real(kind=dp), allocatable :: pressure_time(:)
62  !> Pressures at set-points (Pa).
63  real(kind=dp), allocatable :: pressure(:)
64 
65  !> Height set-point times (s).
66  real(kind=dp), allocatable :: height_time(:)
67  !> Heights at set-points (m).
68  real(kind=dp), allocatable :: height(:)
69 
70  !> Gas emission set-point times (s).
71  real(kind=dp), allocatable :: gas_emission_time(:)
72  !> Gas emisssion rate scales at set-points (1).
73  real(kind=dp), allocatable :: gas_emission_rate_scale(:)
74  !> Gas emission rates at set-points (mol m^{-2} s^{-1}).
75  type(gas_state_t), allocatable :: gas_emission(:)
76 
77  !> Gas-background dilution set-point times (s).
78  real(kind=dp), allocatable :: gas_dilution_time(:)
79  !> Gas-background dilution rates at set-points (s^{-1}).
80  real(kind=dp), allocatable :: gas_dilution_rate(:)
81  !> Background gas mixing ratios at set-points (ppb).
82  type(gas_state_t), allocatable :: gas_background(:)
83 
84  !> Aerosol emission set-points times (s).
85  real(kind=dp), allocatable :: aero_emission_time(:)
86  !> Aerosol emission rate scales at set-points (1).
87  real(kind=dp), allocatable :: aero_emission_rate_scale(:)
88  !> Aerosol emissions at set-points (# m^{-2} s^{-1}).
89  type(aero_dist_t), allocatable :: aero_emission(:)
90 
91  !> Aerosol-background dilution set-point times (s).
92  real(kind=dp), allocatable :: aero_dilution_time(:)
93  !> Aerosol-background dilution rates at set-points (s^{-1}).
94  real(kind=dp), allocatable :: aero_dilution_rate(:)
95  !> Aerosol background at set-points (# m^{-3}).
96  type(aero_dist_t), allocatable :: aero_background(:)
97 
98  !> Type of loss rate function.
99  integer :: loss_function_type
100  !> Chamber parameters for wall loss and sedimentation.
101  type(chamber_t) :: chamber
102  end type scenario_t
103 
104 contains
105 
106 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
107 
108  !> Initialize the time-dependent contents of the
109  !> environment. Thereafter scenario_update_env_state() should be used.
110  subroutine scenario_init_env_state(scenario, env_state, time)
111 
112  !> Scenario data.
113  type(scenario_t), intent(in) :: scenario
114  !> Environment state to update.
115  type(env_state_t), intent(inout) :: env_state
116  !> Current time (s).
117  real(kind=dp), intent(in) :: time
118 
119  env_state%temp = interp_1d(scenario%temp_time, scenario%temp, time)
120  env_state%pressure = interp_1d(scenario%pressure_time, &
121  scenario%pressure, time)
122  env_state%height = interp_1d(scenario%height_time, scenario%height, time)
123  env_state%elapsed_time = time
124  ! FIXME: should compute this at some point
125  env_state%solar_zenith_angle = 0d0
126 
127  end subroutine scenario_init_env_state
128 
129 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
130 
131  !> Update time-dependent contents of the environment.
132  !> scenario_init_env_state() should have been called at the start.
133  subroutine scenario_update_env_state(scenario, env_state, time)
134 
135  !> Scenario data.
136  type(scenario_t), intent(in) :: scenario
137  !> Environment state to update.
138  type(env_state_t), intent(inout) :: env_state
139  !> Current time (s).
140  real(kind=dp), intent(in) :: time
141 
142  !> Ambient water vapor pressure (Pa).
143  real(kind=dp) :: pmv_old, pmv_new
144  !> Ambient pressure (Pa)
145  real(kind=dp) :: pressure_old
146  !> Ambient temperature (K)
147  real(kind=dp) :: temp_old
148 
149  ! Update temperature and pressure and adjust relative humidity to maintain
150  ! water mixing ratio.
151 
152  pmv_old = env_state_sat_vapor_pressure(env_state) * env_state%rel_humid
153  pressure_old = env_state%pressure
154  temp_old = env_state%temp
155 
156  env_state%temp = interp_1d(scenario%temp_time, scenario%temp, time)
157  env_state%pressure = interp_1d(scenario%pressure_time, &
158  scenario%pressure, time)
159 
160  pmv_new = pmv_old * env_state%pressure / pressure_old
161  env_state%rel_humid = pmv_new / env_state_sat_vapor_pressure(env_state)
162 
163  env_state%height = interp_1d(scenario%height_time, scenario%height, time)
164  env_state%elapsed_time = time
165 
166  end subroutine scenario_update_env_state
167 
168 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
169 
170  !> Do gas emissions and background dilution.
171  !!
172  !! Emissions are given as an areal rate \f$e(t)\f$, and then divided
173  !! by the current box height \f$h(t)\f$ to obtain a volume
174  !! rate. There is also a dimensionless rate scaling \f$r(t)\f$. All
175  !! input functions are asusumed constant over the timestep, so the
176  !! concentration \f$c(t)\f$ change is given by
177  !! \f[
178  !! c(t) = c(0) + \frac{r t}{h} e.
179  !! \f]
180  !!
181  !! We model dilution by considering a gas concentration \f$c(t)\f$
182  !! in a box of height \f$h(t)\f$, subject to first-order dilution
183  !! with a rate \f$\lambda(t)\f$. Then the effective dilution rate is
184  !! \f[
185  !! \lambda_{\rm eff}(t) = \lambda(t) + \frac{\dot{h}(t)}{h(t)}
186  !! \f]
187  !! and the evolution of \f$c(t)\f$ is given by
188  !! \f[
189  !! \dot{c}(t) = - \lambda_{\rm eff}(t) c(t).
190  !! \f]
191  !! Solving this with separation of variables gives
192  !! \f[
193  !! \frac{c(t)}{c(0)} = \frac{h(0)}{h(t)}
194  !! \exp\left( - \int_0^t \lambda(t)\,dt\right).
195  !! \f]
196  !! If we define \f$p = c(t)/c(0)\f$ to be the remaining proportion
197  !! of the initial concentration, and \f$b\f$ to be the constant
198  !! background concentration, then we have
199  !! \f[
200  !! c(t) = p(t) c(0) + (1 - p(t)) b.
201  !! \f]
202  !! We assume constant \f$\lambda\f$ and we only do entrainment with
203  !! increasing height \f$h(t)\f$, so we have
204  !! \f[
205  !! p(t) = \min\left(1, \frac{h(0)}{h(t)}\right) \exp(-\lambda t).
206  !! \f]
207  subroutine scenario_update_gas_state(scenario, delta_t, env_state, &
208  old_env_state, gas_data, gas_state)
209 
210  !> Scenario data.
211  type(scenario_t), intent(in) :: scenario
212  !> Time increment to update over.
213  real(kind=dp), intent(in) :: delta_t
214  !> Current environment.
215  type(env_state_t), intent(in) :: env_state
216  !> Previous environment.
217  type(env_state_t), intent(in) :: old_env_state
218  !> Gas data values.
219  type(gas_data_t), intent(in) :: gas_data
220  !> Gas state to update.
221  type(gas_state_t), intent(inout) :: gas_state
222 
223  real(kind=dp) :: emission_rate_scale, dilution_rate, p
224  type(gas_state_t) :: emissions, background
225 
226  ! emissions
227  if (size(scenario%gas_emission) > 0) then
228  call gas_state_interp_1d(scenario%gas_emission, &
229  scenario%gas_emission_time, scenario%gas_emission_rate_scale, &
230  env_state%elapsed_time, emissions, emission_rate_scale)
231  call gas_state_molar_conc_to_ppb(emissions, env_state)
232  p = emission_rate_scale * delta_t / env_state%height
233  call gas_state_add_scaled(gas_state, emissions, p)
234  end if
235 #ifndef PMC_USE_WRF
236  ! dilution
237  call gas_state_interp_1d(scenario%gas_background, &
238  scenario%gas_dilution_time, scenario%gas_dilution_rate, &
239  env_state%elapsed_time, background, dilution_rate)
240  p = exp(- dilution_rate * delta_t)
241  if (env_state%height > old_env_state%height) then
242  p = p * old_env_state%height / env_state%height
243  end if
244  call gas_state_scale(gas_state, p)
245  call gas_state_add_scaled(gas_state, background, 1d0 - p)
246 #endif
247  call gas_state_ensure_nonnegative(gas_state)
248 
249  end subroutine scenario_update_gas_state
250 
251 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
252 
253  !> Do emissions and background dilution for a particle aerosol
254  !> distribution.
255  !!
256  !! See scenario_update_gas_state() for a description of the
257  !! model. We additionally scale the number concentration to account
258  !! for temperature changes.
259  subroutine scenario_update_aero_state(scenario, delta_t, env_state, &
260  old_env_state, aero_data, aero_state, n_emit, n_dil_in, n_dil_out, &
261  allow_doubling, allow_halving)
262 
263  !> Scenario data.
264  type(scenario_t), intent(in) :: scenario
265  !> Time increment to update over.
266  real(kind=dp), intent(in) :: delta_t
267  !> Current environment.
268  type(env_state_t), intent(in) :: env_state
269  !> Previous environment.
270  type(env_state_t), intent(in) :: old_env_state
271  !> Aero data values.
272  type(aero_data_t), intent(in) :: aero_data
273  !> Aero state to update.
274  type(aero_state_t), intent(inout) :: aero_state
275  !> Number of emitted particles.
276  integer, intent(out) :: n_emit
277  !> Number of diluted-in particles.
278  integer, intent(out) :: n_dil_in
279  !> Number of diluted-out particles.
280  integer, intent(out) :: n_dil_out
281  !> Whether to allow doubling of the population.
282  logical, intent(in) :: allow_doubling
283  !> Whether to allow halving of the population.
284  logical, intent(in) :: allow_halving
285 
286  real(kind=dp), parameter :: sample_timescale = 3600.0d0
287  real(kind=dp) :: characteristic_factor, sample_timescale_effective
288  real(kind=dp) :: emission_rate_scale, dilution_rate, p
289  type(aero_dist_t) :: emissions, background
290  type(aero_state_t) :: aero_state_delta
291 
292  sample_timescale_effective = max(1.0, min(sample_timescale, &
293  env_state%elapsed_time))
294  characteristic_factor = sample_timescale_effective / delta_t
295 
296  ! emissions
297  if (size(scenario%aero_emission) > 0) then
298  call aero_dist_interp_1d(scenario%aero_emission, &
299  scenario%aero_emission_time, scenario%aero_emission_rate_scale, &
300  env_state%elapsed_time, emissions, emission_rate_scale)
301  p = emission_rate_scale * delta_t / env_state%height
302  call aero_state_add_aero_dist_sample(aero_state, aero_data, &
303  emissions, p, characteristic_factor, env_state%elapsed_time, &
304  allow_doubling, allow_halving, n_emit)
305  end if
306 #ifndef PMC_USE_WRF
307  ! dilution
308  call aero_dist_interp_1d(scenario%aero_background, &
309  scenario%aero_dilution_time, scenario%aero_dilution_rate, &
310  env_state%elapsed_time, background, dilution_rate)
311  p = exp(- dilution_rate * delta_t)
312  if (env_state%height > old_env_state%height) then
313  p = p * old_env_state%height / env_state%height
314  end if
315  ! loss to background
316  call aero_state_zero(aero_state_delta)
317  call aero_state_copy_weight(aero_state, aero_state_delta)
318  call aero_state_sample_particles(aero_state, aero_state_delta, &
319  aero_data, 1d0 - p, aero_info_dilution)
320  n_dil_out = aero_state_total_particles(aero_state_delta)
321  ! addition from background
322  call aero_state_add_aero_dist_sample(aero_state, aero_data, &
323  background, 1d0 - p, characteristic_factor, env_state%elapsed_time, &
324  allow_doubling, allow_halving, n_dil_in)
325 
326  ! particle loss function
327  call scenario_particle_loss(scenario, delta_t, aero_data, aero_state, &
328  env_state)
329 #endif
330 
331 #ifdef PMC_USE_WRF
332  call aero_weight_array_scale(aero_state%awa, &
333  old_env_state%inverse_density * (1.0d0 / env_state%inverse_density))
334 #else
335  ! update computational volume
336  call aero_weight_array_scale(aero_state%awa, &
337  old_env_state%temp * env_state%pressure &
338  / (env_state%temp * old_env_state%pressure))
339 #endif
340 
341  end subroutine scenario_update_aero_state
342 
343 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
344 
345  !> Do emissions and background dilution from the environment for a
346  !> binned aerosol distribution.
347  !!
348  !! See scenario_update_gas_state() for a description of the model.
349  subroutine scenario_update_aero_binned(scenario, delta_t, env_state, &
350  old_env_state, bin_grid, aero_data, aero_binned)
351 
352  !> Scenario data.
353  type(scenario_t), intent(in) :: scenario
354  !> Time increment to update over.
355  real(kind=dp), intent(in) :: delta_t
356  !> Current environment.
357  type(env_state_t), intent(in) :: env_state
358  !> Previous environment.
359  type(env_state_t), intent(in) :: old_env_state
360  !> Bin grid.
361  type(bin_grid_t), intent(in) :: bin_grid
362  !> Aero data values.
363  type(aero_data_t), intent(in) :: aero_data
364  !> Aero binned to update.
365  type(aero_binned_t), intent(inout) :: aero_binned
366 
367  real(kind=dp) :: emission_rate_scale, dilution_rate, p
368  type(aero_dist_t) :: emissions, background
369  type(aero_binned_t) :: emissions_binned, background_binned
370 
371  ! emissions
372  call aero_dist_interp_1d(scenario%aero_emission, &
373  scenario%aero_emission_time, scenario%aero_emission_rate_scale, &
374  env_state%elapsed_time, emissions, emission_rate_scale)
375  call aero_binned_add_aero_dist(emissions_binned, bin_grid, aero_data, &
376  emissions)
377  p = emission_rate_scale * delta_t / env_state%height
378  call aero_binned_add_scaled(aero_binned, emissions_binned, p)
379 
380  ! dilution
381  call aero_dist_interp_1d(scenario%aero_background, &
382  scenario%aero_dilution_time, scenario%aero_dilution_rate, &
383  env_state%elapsed_time, background, dilution_rate)
384  call aero_binned_add_aero_dist(background_binned, bin_grid, aero_data, &
385  background)
386  p = exp(- dilution_rate * delta_t)
387  if (env_state%height > old_env_state%height) then
388  p = p * old_env_state%height / env_state%height
389  end if
390  call aero_binned_scale(aero_binned, p)
391  call aero_binned_add_scaled(aero_binned, background_binned, 1d0 - p)
392 
393  end subroutine scenario_update_aero_binned
394 
395 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
396 
397  !> Evaluate a loss rate function.
398  real(kind=dp) function scenario_loss_rate(scenario, vol, density, &
399  aero_data, env_state)
400 
401  !> Scenario data.
402  type(scenario_t), intent(in) :: scenario
403  !> Volume of particle (m^3).
404  real(kind=dp), intent(in) :: vol
405  !> Density of particle (kg/m^3).
406  real(kind=dp), intent(in) :: density
407  !> Aerosol data.
408  type(aero_data_t), intent(in) :: aero_data
409  !> Environment state.
410  type(env_state_t), intent(in) :: env_state
411 
412  scenario_loss_rate = 0d0
413  if (scenario%loss_function_type == scenario_loss_function_invalid) then
414  scenario_loss_rate = 0d0
415  else if (scenario%loss_function_type == scenario_loss_function_none) then
416  scenario_loss_rate = 0d0
417  else if (scenario%loss_function_type == scenario_loss_function_constant) &
418  then
419  scenario_loss_rate = 1d-3
420  else if (scenario%loss_function_type == scenario_loss_function_volume) then
421  scenario_loss_rate = 1d15*vol
422  else if (scenario%loss_function_type == scenario_loss_function_drydep) &
423  then
425  aero_data, env_state)
426  else if (scenario%loss_function_type == scenario_loss_function_chamber) &
427  then
428  scenario_loss_rate = chamber_loss_rate_wall(scenario%chamber, vol, &
429  aero_data, env_state) &
430  + chamber_loss_rate_sedi(scenario%chamber, vol, density, &
431  aero_data, env_state)
432  else
433  call die_msg(201594391, "Unknown loss function id: " &
434  // trim(integer_to_string(scenario%loss_function_type)))
435  end if
436 
437  end function scenario_loss_rate
438 
439 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
440 
441  !> Compute and return the dry deposition rate for a given particle.
442  !! All equations used here are written in detail in the file
443  !! \c doc/deposition/deposition.tex.
444  real(kind=dp) function scenario_loss_rate_dry_dep(vol, density, aero_data, &
445  env_state)
446 
447  !> Particle volume (m^3).
448  real(kind=dp), intent(in) :: vol
449  !> Particle density (kg m^-3).
450  real(kind=dp), intent(in) :: density
451  !> Aerosol data.
452  type(aero_data_t), intent(in) :: aero_data
453  !> Environment state.
454  type(env_state_t), intent(in) :: env_state
455 
456  real(kind=dp) :: v_d
457  real(kind=dp) :: v_s
458  real(kind=dp) :: d_p
459  real(kind=dp) :: density_air
460  real(kind=dp) :: visc_d, visc_k
461  real(kind=dp) :: gas_speed, gas_mean_free_path
462  real(kind=dp) :: knud, cunning
463  real(kind=dp) :: grav
464  real(kind=dp) :: r_s, r_a
465  real(kind=dp) :: alpha, beta, gamma, a, eps_0
466  real(kind=dp) :: diff_p
467  real(kind=dp) :: von_karman
468  real(kind=dp) :: st, sc, u_star
469  real(kind=dp) :: e_b, e_im, e_in, r1
470  real(kind=dp) :: u_mean, z_ref, z_rough
471 
472  ! User set variables
473  u_mean = 5.0d0 ! Mean wind speed at reference height
474  z_ref = 20.0d0 ! Reference height
475  ! Setting for LUC = 7, SC = 1 - See Table 3
476  z_rough = .1d0 ! According to land category
477  a = 2.0d0 / 1000.0d0 ! Dependent on land type
478  alpha = 1.2d0 ! From table
479  beta = 2.0d0 ! From text
480  gamma = .54d0 ! From table
481 
482  ! particle diameter
483  d_p = aero_data_vol2diam(aero_data, vol)
484  ! density of air
485  density_air = (const%air_molec_weight * env_state%pressure) &
486  / (const%univ_gas_const * env_state%temp)
487  ! dynamic viscosity
488  visc_d = 1.8325d-5 * (416.16 / (env_state%temp + 120.0d0)) &
489  * (env_state%temp / 296.16)**1.5d0
490  ! kinematic viscosity
491  visc_k = visc_d / density_air
492  ! gas speed
493  gas_speed = &
494  sqrt((8.0d0 * const%boltzmann * env_state%temp * const%avagadro) / &
495  (const%pi * const%air_molec_weight))
496  ! gas free path
497  gas_mean_free_path = (2.0d0 * visc_d) / (density_air * gas_speed)
498  ! knudson number
499  knud = (2.0d0 * gas_mean_free_path) / d_p
500  ! cunningham correction factor
501  cunning = 1.0d0 + knud * (1.257d0 + 0.4d0 * exp(-1.1d0 / knud))
502  ! gravity
503  grav = 9.81d0
504  ! Compute V_s
505  v_s = (density * d_p**2.0d0 * grav * cunning) / (18.0d0 * visc_d)
506 
507  ! Aerodynamic resistance
508  ! For neutral stability
509  u_star = .4d0 * u_mean / log(z_ref / z_rough)
510  r_a = (1.0d0 / (.4d0 * u_star)) * log(z_ref / z_rough)
511  ! Brownian diffusion efficiency
512  diff_p = (const%boltzmann * env_state%temp * cunning) / &
513  (3.d0 * const%pi * visc_d * d_p)
514  sc = visc_k / diff_p
515  e_b = sc**(-gamma)
516 
517  ! Interception efficiency
518  ! Characteristic radius of large collectors
519  e_in = .5d0 * (d_p / a)**2.0d0
520 
521  ! Impaction efficiency
522  st = (v_s * u_star) / (grav * a)
523  e_im = (st / (alpha + st))**beta
524 
525  ! Rebound correction
526  r1 = exp(-st**.5d0)
527 
528  ! Surface resistance
529  eps_0 = 3.0d0 ! Taken to be 3
530  r_s = 1.0d0 / (eps_0 * u_star * (e_b + e_in + e_im) * r1)
531 
532  ! Dry deposition
533  v_d = v_s + (1.0d0 / (r_a + r_s + r_a * r_s * v_s))
534 
535  ! The loss rate
536  scenario_loss_rate_dry_dep = v_d / env_state%height
537 
538  end function scenario_loss_rate_dry_dep
539 
540 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
541 
542  !> Compute and return the max loss rate function for a given volume.
543  real(kind=dp) function scenario_loss_rate_max(scenario, vol, aero_data, &
544  env_state)
545 
546  !> Scenario data.
547  type(scenario_t), intent(in) :: scenario
548  !> Particle volume (m^3).
549  real(kind=dp), intent(in) :: vol
550  !> Aerosol data.
551  type(aero_data_t), intent(in) :: aero_data
552  !> Environment state.
553  type(env_state_t), intent(in) :: env_state
554 
555  !> Number of density sample points.
556  integer, parameter :: n_sample = 3
557 
558  real(kind=dp) :: d, d_min, d_max, loss
559  integer :: i
560 
561  d_min = minval(aero_data%density)
562  d_max = maxval(aero_data%density)
563 
565  do i = 1,n_sample
566  d = interp_linear_disc(d_min, d_max, n_sample, i)
567  loss = scenario_loss_rate(scenario, vol, d, aero_data, env_state)
569  end do
570 
571  end function scenario_loss_rate_max
572 
573 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
574 
575  !> Compute an upper bound on the maximum kernel value for each
576  !> bin.
577  !! Value over_scale is multiplied to the maximum sampled value
578  !! to get the upper bound. A tighter bound may be reached if over_scale
579  !! is smaller, but that also risks falling below a kernel value.
580  subroutine scenario_loss_rate_bin_max(scenario, bin_grid, aero_data, &
581  env_state, loss_max)
582 
583  !> Scenario data.
584  type(scenario_t), intent(in) :: scenario
585  !> Bin_grid.
586  type(bin_grid_t), intent(in) :: bin_grid
587  !> Aerosol data.
588  type(aero_data_t), intent(in) :: aero_data
589  !> Environment state.
590  type(env_state_t), intent(in) :: env_state
591  !> Maximum loss vals.
592  real(kind=dp), intent(out) :: loss_max(bin_grid_size(bin_grid))
593 
594  !> Number of sample points per bin.
595  integer, parameter :: n_sample = 3
596  !> Over-estimation scale factor parameter.
597  real(kind=dp), parameter :: over_scale = 2d0
598 
599  real(kind=dp) :: v_low, v_high, vol, r, r_max
600  integer :: b, i
601 
602  do b = 1,bin_grid_size(bin_grid)
603  v_low = aero_data_rad2vol(aero_data, bin_grid%edges(b))
604  v_high = aero_data_rad2vol(aero_data, bin_grid%edges(b + 1))
605  r_max = 0d0
606  do i = 1,n_sample
607  vol = interp_linear_disc(v_low, v_high, n_sample, i)
608  r = scenario_loss_rate_max(scenario, vol, aero_data, env_state)
609  r_max = max(r_max, r)
610  end do
611  loss_max(b) = r_max * over_scale
612  end do
613 
614  end subroutine scenario_loss_rate_bin_max
615 
616 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
617 
618  !> Performs stochastic particle loss for one time-step.
619  !! If a particle \c i_part has a scenario_loss_rate() value of rate, then the
620  !! probability p will be removed by this function is
621  !! <tt>1 - exp(-delta_t*rate)</tt>.
622  !! Uses an accept-reject algorithm for efficiency, in which a particle
623  !! is first sampled with rate <tt>1 - exp(-delta_t*over_rate) </tt>
624  !! and then accepted with rate
625  !! <tt>(1 - exp(-delta_t*rate))/(1 - exp(-delta_t*over_rate))</tt>.
626  subroutine scenario_particle_loss(scenario, delta_t, aero_data, aero_state, &
627  env_state)
628 
629  !> Scenario data.
630  type(scenario_t), intent(in) :: scenario
631  !> Time increment to update over.
632  real(kind=dp), intent(in) :: delta_t
633  !> Aerosol data.
634  type(aero_data_t), intent(in) :: aero_data
635  !> Aerosol state.
636  type(aero_state_t), intent(inout) :: aero_state
637  !> Environment state.
638  type(env_state_t), intent(in) :: env_state
639 
640  integer :: c, b, s, i_part
641  real(kind=dp) :: over_rate, over_prob, rand_real, rand_geom
642 
643  if (scenario%loss_function_type == scenario_loss_function_none .or. &
644  scenario%loss_function_type == scenario_loss_function_invalid) return
645 
646  if (scenario_loss_alg_threshold <= 0d0) then
647  ! use naive algorithm for everything
648  do i_part = aero_state%apa%n_part, 1, -1
649  call scenario_try_single_particle_loss(scenario, delta_t, &
650  aero_data, aero_state, env_state, i_part, 1d0)
651  end do
652  return
653  end if
654 
655  call aero_state_sort(aero_state, aero_data)
656 
657  if (.not. aero_state%aero_sorted%removal_rate_bounds_valid) then
658  call scenario_loss_rate_bin_max(scenario, &
659  aero_state%aero_sorted%bin_grid, aero_data, env_state, &
660  aero_state%aero_sorted%removal_rate_max)
661  aero_state%aero_sorted%removal_rate_bounds_valid = .true.
662  end if
663 
664  do c = 1,aero_sorted_n_class(aero_state%aero_sorted)
665  do b = 1,aero_sorted_n_bin(aero_state%aero_sorted)
666  over_rate = aero_state%aero_sorted%removal_rate_max(b)
667  if (delta_t * over_rate <= 0d0) cycle
668  over_prob = 1d0 - exp(-delta_t * over_rate)
669  if (over_prob >= scenario_loss_alg_threshold) then
670  ! use naive algorithm over bin
671  do s = aero_state%aero_sorted%size_class%inverse(b, c)%n_entry, &
672  1,-1
673  i_part = &
674  aero_state%aero_sorted%size_class%inverse(b, c)%entry(s)
675  call scenario_try_single_particle_loss(scenario, delta_t, &
676  aero_data, aero_state, env_state, i_part, 1d0)
677  end do
678  else
679  ! use accept-reject algorithm over bin
680  s = aero_state%aero_sorted%size_class%inverse(b, c)%n_entry + 1
681  do while (.true.)
682  rand_real = pmc_random()
683  if (rand_real <= 0d0) exit
684  rand_geom = -log(rand_real) / (delta_t * over_rate) + 1d0
685  if (rand_geom >= real(s, kind=dp)) exit
686  s = s - floor(rand_geom)
687 
688  ! note: floor(rand_geom) is a random geometric variable
689  ! with accept probability 1 - exp(-delta_t*over_rate)
690 
691  i_part = &
692  aero_state%aero_sorted%size_class%inverse(b, c)%entry(s)
693  call scenario_try_single_particle_loss(scenario, delta_t, &
694  aero_data, aero_state, env_state, i_part, over_prob)
695  end do
696  end if
697  end do
698  end do
699 
700  !call aero_state_check_sort(aero_state)
701 
702  end subroutine scenario_particle_loss
703 
704 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
705 
706  !> Test a candidate particle to see if it should be removed,
707  !> and remove if necessary.
708  !! Particle is removed with probability
709  !! (1d0 - exp(-delta_t*rate))/over_prob, where rate is the loss function
710  !! evaluated for the given particle.
711  subroutine scenario_try_single_particle_loss(scenario, delta_t, &
712  aero_data, aero_state, env_state, i_part, over_prob)
713 
714  !> Scenario data.
715  type(scenario_t), intent(in) :: scenario
716  !> Time increment to update over.
717  real(kind=dp), intent(in) :: delta_t
718  !> Aerosol data.
719  type(aero_data_t), intent(in) :: aero_data
720  !> Aerosol state.
721  type(aero_state_t), intent(inout) :: aero_state
722  !> Environment state.
723  type(env_state_t), intent(in) :: env_state
724  !> Index of particle to attempt removal
725  integer, intent(in) :: i_part
726  !> Overestimated removal probability used previously
727  real(kind=dp), intent(in) :: over_prob
728 
729  real(kind=dp) :: prob, rate, vol, density
730  type(aero_info_t) :: aero_info
731 
732  vol = aero_particle_volume(aero_state%apa%particle(i_part))
733  density = aero_particle_density(aero_state%apa%particle(i_part), aero_data)
734  rate = scenario_loss_rate(scenario, vol, density, aero_data, env_state)
735  prob = 1d0 - exp(-delta_t * rate)
736  call warn_assert_msg(295846288, prob <= over_prob, &
737  "particle loss upper bound estimation is too tight: " &
738  // trim(real_to_string(prob)) // " > " &
739  // trim(real_to_string(over_prob)) )
740  if (pmc_random() * over_prob > prob) return
741 
742  aero_info%id = aero_state%apa%particle(i_part)%id
743  aero_info%action = aero_info_dilution
744  aero_info%other_id = 0
745  call aero_state_remove_particle_with_info(aero_state, i_part, aero_info)
746 
747  end subroutine scenario_try_single_particle_loss
748 
749 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
750 
751  !> Whether any of the contained aerosol modes are of the given type.
752  elemental logical function scenario_contains_aero_mode_type(scenario, &
753  aero_mode_type)
754 
755  !> Scenario data.
756  type(scenario_t), intent(in) :: scenario
757  !> Aerosol mode type to test for.
758  integer, intent(in) :: aero_mode_type
759 
761  = any(aero_dist_contains_aero_mode_type(scenario%aero_emission, &
762  aero_mode_type)) &
763  .or. any(aero_dist_contains_aero_mode_type(scenario%aero_background, &
764  aero_mode_type))
765 
767 
768 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
769 
770  !> Read environment data from an spec file.
771  subroutine spec_file_read_scenario(file, gas_data, aero_data, &
772  read_aero_weight_classes, scenario)
773 
774  !> Spec file.
775  type(spec_file_t), intent(inout) :: file
776  !> Gas data values.
777  type(gas_data_t), intent(in) :: gas_data
778  !> Aerosol data.
779  type(aero_data_t), intent(inout) :: aero_data
780  !> Whether the weight classes for each source are specified in inputs.
781  logical, intent(in) :: read_aero_weight_classes
782  !> Scenario data.
783  type(scenario_t), intent(inout) :: scenario
784 
785  character(len=PMC_MAX_FILENAME_LEN) :: sub_filename
786  type(spec_file_t) :: sub_file
787  character(len=SPEC_LINE_MAX_VAR_LEN) :: function_name
788 
789  ! note that we have to hard-code the list for doxygen below
790 
791  !> \page input_format_scenario Input File Format: Scenario
792  !!
793  !! The scenario parameters are:
794  !! <ul>
795  !! <li> \b temp_profile (string): the name of the file from which to
796  !! read the temperature profile --- the file format should be
797  !! \subpage input_format_temp_profile
798  !! <li> \b pressure_profile (string): the name of the file from which to
799  !! read the pressure profile --- the file format should be
800  !! \subpage input_format_pressure_profile
801  !! <li> \b height_profile (string): the name of the file from which
802  !! to read the mixing layer height profile --- the file format
803  !! should be \subpage input_format_height_profile
804  !! <li> \b gas_emissions (string): the name of the file from which to
805  !! read the gas emissions profile --- the file format should be
806  !! \subpage input_format_gas_profile
807  !! <li> \b gas_background (string): the name of the file from which
808  !! to read the gas background profile --- the file format should
809  !! be \subpage input_format_gas_profile
810  !! <li> \b aero_emissions (string): the name of the file from which
811  !! to read the aerosol emissions profile --- the file format
812  !! should be \subpage input_format_aero_dist_profile
813  !! <li> \b aero_background (string): the name of the file from which
814  !! to read the aerosol background profile --- the file format
815  !! should be \subpage input_format_aero_dist_profile
816  !! <li> \b loss_function (string): the type of loss function ---
817  !! must be one of: \c none for no particle loss, \c constant
818  !! for constant loss rate, \c volume for particle loss proportional
819  !! to particle volume, \c drydep for particle loss proportional
820  !! to dry deposition velocity, or \c chamber for a chamber model.
821  !! If \c loss_function is \c chamber, then the following
822  !! parameters must also be provided:
823  !! - \subpage input_format_chamber
824  !! </ul>
825  !!
826  !! See also:
827  !! - \ref spec_file_format --- the input file text format
828 
829  ! temperature profile
830  call spec_file_read_string(file, "temp_profile", sub_filename)
831  call spec_file_open(sub_filename, sub_file)
832  call spec_file_read_timed_real_array(sub_file, "temp", &
833  scenario%temp_time, scenario%temp)
834  call spec_file_close(sub_file)
835 
836  ! pressure profile
837  call spec_file_read_string(file, "pressure_profile", sub_filename)
838  call spec_file_open(sub_filename, sub_file)
839  call spec_file_read_timed_real_array(sub_file, "pressure", &
840  scenario%pressure_time, scenario%pressure)
841  call spec_file_close(sub_file)
842 
843  ! height profile
844  call spec_file_read_string(file, "height_profile", sub_filename)
845  call spec_file_open(sub_filename, sub_file)
846  call spec_file_read_timed_real_array(sub_file, "height", &
847  scenario%height_time, scenario%height)
848  call spec_file_close(sub_file)
849 
850  ! gas emissions profile
851  call spec_file_read_string(file, "gas_emissions", sub_filename)
852  call spec_file_open(sub_filename, sub_file)
853  call spec_file_read_gas_states_times_rates(sub_file, gas_data, &
854  scenario%gas_emission_time, scenario%gas_emission_rate_scale, &
855  scenario%gas_emission)
856  call spec_file_close(sub_file)
857 
858  ! gas background profile
859  call spec_file_read_string(file, "gas_background", sub_filename)
860  call spec_file_open(sub_filename, sub_file)
861  call spec_file_read_gas_states_times_rates(sub_file, gas_data, &
862  scenario%gas_dilution_time, scenario%gas_dilution_rate, &
863  scenario%gas_background)
864  call spec_file_close(sub_file)
865 
866  ! aerosol emissions profile
867  call spec_file_read_string(file, "aero_emissions", sub_filename)
868  call spec_file_open(sub_filename, sub_file)
869  call spec_file_read_aero_dists_times_rates(sub_file, aero_data, &
870  read_aero_weight_classes, scenario%aero_emission_time, &
871  scenario%aero_emission_rate_scale, scenario%aero_emission)
872  call spec_file_close(sub_file)
873 
874  ! aerosol background profile
875  call spec_file_read_string(file, "aero_background", sub_filename)
876  call spec_file_open(sub_filename, sub_file)
877  call spec_file_read_aero_dists_times_rates(sub_file, aero_data, &
878  read_aero_weight_classes, scenario%aero_dilution_time, &
879  scenario%aero_dilution_rate, scenario%aero_background)
880  call spec_file_close(sub_file)
881 
882  ! loss function
883  call spec_file_read_string(file, 'loss_function', function_name)
884  if (trim(function_name) == 'none') then
885  scenario%loss_function_type = scenario_loss_function_none
886  else if (trim(function_name) == 'constant') then
887  scenario%loss_function_type = scenario_loss_function_constant
888  else if (trim(function_name) == 'volume') then
889  scenario%loss_function_type = scenario_loss_function_volume
890  else if (trim(function_name) == 'drydep') then
891  scenario%loss_function_type = scenario_loss_function_drydep
892  else if (trim(function_name) == 'chamber') then
893  scenario%loss_function_type = scenario_loss_function_chamber
894  call spec_file_read_chamber(file, scenario%chamber)
895  else
896  call spec_file_die_msg(518248400, file, &
897  "Unknown loss function type: " // trim(function_name))
898  end if
899 
900  end subroutine spec_file_read_scenario
901 
902  ! the following blocks belong in the subroutine above, but they are
903  ! outside because Doxygen 1.8.7 doesn't resolve references when
904  ! multiple \page blocks are in one subroutine
905 
906  !> \page input_format_temp_profile Input File Format: Temperature Profile
907  !!
908  !! A temperature profile input file must consist of two lines:
909  !! - the first line must begin with \c time and should be followed
910  !! by \f$N\f$ space-separated real scalars, giving the times (in
911  !! s after the start of the simulation) of the temperature set
912  !! points --- the times must be in increasing order
913  !! - the second line must begin with \c temp and should be followed
914  !! by \f$N\f$ space-separated real scalars, giving the
915  !! temperatures (in K) at the corresponding times
916  !!
917  !! The temperature profile is linearly interpolated between the
918  !! specified times, while before the first time it takes the first
919  !! temperature value and after the last time it takes the last
920  !! temperature value.
921  !!
922  !! Example:
923  !! <pre>
924  !! time 0 600 1800 # time (in s) after simulation start
925  !! temp 270 290 280 # temperature (in K)
926  !! </pre>
927  !! Here the temperature starts at 270&nbsp;K at the start of the
928  !! simulation, rises to 290&nbsp;K after 10&nbsp;min, and then
929  !! falls again to 280&nbsp;K at 30&nbsp;min. Between these times
930  !! the temperature is linearly interpolated, while after
931  !! 30&nbsp;min it is held constant at 280&nbsp;K.
932  !!
933  !! See also:
934  !! - \ref spec_file_format --- the input file text format
935  !! - \ref input_format_scenario --- the environment data
936  !! containing the temperature profile
937 
938  !> \page input_format_pressure_profile Input File Format: Pressure Profile
939  !!
940  !! A pressure profile input file must consist of two lines:
941  !! - the first line must begin with \c time and should be followed
942  !! by \f$N\f$ space-separated real scalars, giving the times (in
943  !! s after the start of the simulation) of the pressure set
944  !! points --- the times must be in increasing order
945  !! - the second line must begin with \c pressure and should be followed
946  !! by \f$N\f$ space-separated real scalars, giving the
947  !! pressures (in Pa) at the corresponding times
948  !!
949  !! The pressure profile is linearly interpolated between the
950  !! specified times, while before the first time it takes the first
951  !! pressure value and after the last time it takes the last
952  !! pressure value.
953  !!
954  !! Example:
955  !! <pre>
956  !! time 0 600 1800 # time (in s) after simulation start
957  !! pressure 1e5 9e4 7.5e4 # pressure (in Pa)
958  !! </pre>
959  !! Here the pressure starts at 1e5&nbsp;Pa at the start of the
960  !! simulation, decreases to 9e4&nbsp;Pa after 10&nbsp;min, and then
961  !! decreases further to 7.5e4&nbsp;Pa at 30&nbsp;min. Between these times
962  !! the pressure is linearly interpolated, while after
963  !! 30&nbsp;min it is held constant at 7.5e4&nbsp;Pa.
964  !!
965  !! See also:
966  !! - \ref spec_file_format --- the input file text format
967  !! - \ref input_format_scenario --- the environment data
968  !! containing the pressure profile
969 
970  !> \page input_format_height_profile Input File Format: Mixing Layer Height Profile
971  !!
972  !! A mixing layer height profile input file must consist of two
973  !! lines:
974  !! - the first line must begin with \c time and should be followed
975  !! by \f$N\f$ space-separated real scalars, giving the times (in
976  !! s after the start of the simulation) of the height set
977  !! points --- the times must be in increasing order
978  !! - the second line must begin with \c height and should be
979  !! followed by \f$N\f$ space-separated real scalars, giving the
980  !! mixing layer heights (in m) at the corresponding times
981  !!
982  !! The mixing layer height profile is linearly interpolated
983  !! between the specified times, while before the first time it
984  !! takes the first height value and after the last time it takes
985  !! the last height value.
986  !!
987  !! Example:
988  !! <pre>
989  !! time 0 600 1800 # time (in s) after simulation start
990  !! height 500 1000 800 # mixing layer height (in m)
991  !! </pre>
992  !! Here the mixing layer height starts at 500&nbsp;m at the start
993  !! of the simulation, rises to 1000&nbsp;m after 10&nbsp;min, and
994  !! then falls again to 800&nbsp;m at 30&nbsp;min. Between these
995  !! times the mixing layer height is linearly interpolated, while
996  !! after 30&nbsp;min it is held constant at 800&nbsp;m.
997  !!
998  !! See also:
999  !! - \ref spec_file_format --- the input file text format
1000  !! - \ref input_format_scenario --- the environment data
1001  !! containing the mixing layer height profile
1002 
1003 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1004 
1005  !> Determines the number of bytes required to pack the given value.
1006  integer function pmc_mpi_pack_size_scenario(val)
1007 
1008  !> Value to pack.
1009  type(scenario_t), intent(in) :: val
1010 
1011  integer :: total_size, i, n
1012 
1013  total_size = &
1014  pmc_mpi_pack_size_real_array(val%temp_time) &
1015  + pmc_mpi_pack_size_real_array(val%temp) &
1016  + pmc_mpi_pack_size_real_array(val%pressure_time) &
1017  + pmc_mpi_pack_size_real_array(val%pressure) &
1018  + pmc_mpi_pack_size_real_array(val%height_time) &
1019  + pmc_mpi_pack_size_real_array(val%height) &
1020  + pmc_mpi_pack_size_real_array(val%gas_emission_time) &
1021  + pmc_mpi_pack_size_real_array(val%gas_emission_rate_scale) &
1022  + pmc_mpi_pack_size_real_array(val%gas_dilution_time) &
1023  + pmc_mpi_pack_size_real_array(val%gas_dilution_rate) &
1024  + pmc_mpi_pack_size_real_array(val%aero_emission_time) &
1025  + pmc_mpi_pack_size_real_array(val%aero_emission_rate_scale) &
1026  + pmc_mpi_pack_size_real_array(val%aero_dilution_time) &
1027  + pmc_mpi_pack_size_real_array(val%aero_dilution_rate) &
1028  + pmc_mpi_pack_size_integer(val%loss_function_type) &
1029  + pmc_mpi_pack_size_chamber(val%chamber)
1030  if (allocated(val%gas_emission_time)) then
1031  do i = 1,size(val%gas_emission)
1032  total_size = total_size &
1033  + pmc_mpi_pack_size_gas_state(val%gas_emission(i))
1034  end do
1035  end if
1036  if (allocated(val%gas_dilution_time)) then
1037  do i = 1,size(val%gas_background)
1038  total_size = total_size &
1039  + pmc_mpi_pack_size_gas_state(val%gas_background(i))
1040  end do
1041  end if
1042  if (allocated(val%aero_emission_time)) then
1043  do i = 1,size(val%aero_emission)
1044  total_size = total_size &
1045  + pmc_mpi_pack_size_aero_dist(val%aero_emission(i))
1046  end do
1047  end if
1048  if (allocated(val%aero_dilution_time)) then
1049  do i = 1,size(val%aero_background)
1050  total_size = total_size &
1051  + pmc_mpi_pack_size_aero_dist(val%aero_background(i))
1052  end do
1053  end if
1054 
1055  pmc_mpi_pack_size_scenario = total_size
1056 
1057  end function pmc_mpi_pack_size_scenario
1058 
1059 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1060 
1061  !> Packs the given value into the buffer, advancing position.
1062  subroutine pmc_mpi_pack_scenario(buffer, position, val)
1063 
1064  !> Memory buffer.
1065  character, intent(inout) :: buffer(:)
1066  !> Current buffer position.
1067  integer, intent(inout) :: position
1068  !> Value to pack.
1069  type(scenario_t), intent(in) :: val
1070 
1071 #ifdef PMC_USE_MPI
1072  integer :: prev_position, i
1073 
1074  prev_position = position
1075  call pmc_mpi_pack_real_array(buffer, position, val%temp_time)
1076  call pmc_mpi_pack_real_array(buffer, position, val%temp)
1077  call pmc_mpi_pack_real_array(buffer, position, val%pressure_time)
1078  call pmc_mpi_pack_real_array(buffer, position, val%pressure)
1079  call pmc_mpi_pack_real_array(buffer, position, val%height_time)
1080  call pmc_mpi_pack_real_array(buffer, position, val%height)
1081  call pmc_mpi_pack_real_array(buffer, position, val%gas_emission_time)
1082  call pmc_mpi_pack_real_array(buffer, position, val%gas_emission_rate_scale)
1083  call pmc_mpi_pack_real_array(buffer, position, val%gas_dilution_time)
1084  call pmc_mpi_pack_real_array(buffer, position, val%gas_dilution_rate)
1085  call pmc_mpi_pack_real_array(buffer, position, val%aero_emission_time)
1086  call pmc_mpi_pack_real_array(buffer, position, &
1087  val%aero_emission_rate_scale)
1088  call pmc_mpi_pack_real_array(buffer, position, val%aero_dilution_time)
1089  call pmc_mpi_pack_real_array(buffer, position, val%aero_dilution_rate)
1090  call pmc_mpi_pack_integer(buffer, position, val%loss_function_type)
1091  call pmc_mpi_pack_chamber(buffer, position, val%chamber)
1092  if (allocated(val%gas_emission_time)) then
1093  do i = 1,size(val%gas_emission)
1094  call pmc_mpi_pack_gas_state(buffer, position, val%gas_emission(i))
1095  end do
1096  end if
1097  if (allocated(val%gas_dilution_time)) then
1098  do i = 1,size(val%gas_background)
1099  call pmc_mpi_pack_gas_state(buffer, position, val%gas_background(i))
1100  end do
1101  end if
1102  if (allocated(val%aero_emission_time)) then
1103  do i = 1,size(val%aero_emission)
1104  call pmc_mpi_pack_aero_dist(buffer, position, val%aero_emission(i))
1105  end do
1106  end if
1107  if (allocated(val%aero_dilution_time)) then
1108  do i = 1,size(val%aero_background)
1109  call pmc_mpi_pack_aero_dist(buffer, position, val%aero_background(i))
1110  end do
1111  end if
1112  call assert(639466930, &
1113  position - prev_position <= pmc_mpi_pack_size_scenario(val))
1114 #endif
1115 
1116  end subroutine pmc_mpi_pack_scenario
1117 
1118 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1119 
1120  !> Unpacks the given value from the buffer, advancing position.
1121  subroutine pmc_mpi_unpack_scenario(buffer, position, val)
1122 
1123  !> Memory buffer.
1124  character, intent(inout) :: buffer(:)
1125  !> Current buffer position.
1126  integer, intent(inout) :: position
1127  !> Value to pack.
1128  type(scenario_t), intent(inout) :: val
1129 
1130 #ifdef PMC_USE_MPI
1131  integer :: prev_position, i
1132 
1133  prev_position = position
1134  call pmc_mpi_unpack_real_array(buffer, position, val%temp_time)
1135  call pmc_mpi_unpack_real_array(buffer, position, val%temp)
1136  call pmc_mpi_unpack_real_array(buffer, position, val%pressure_time)
1137  call pmc_mpi_unpack_real_array(buffer, position, val%pressure)
1138  call pmc_mpi_unpack_real_array(buffer, position, val%height_time)
1139  call pmc_mpi_unpack_real_array(buffer, position, val%height)
1140  call pmc_mpi_unpack_real_array(buffer, position, val%gas_emission_time)
1141  call pmc_mpi_unpack_real_array(buffer, position, &
1142  val%gas_emission_rate_scale)
1143  call pmc_mpi_unpack_real_array(buffer, position, val%gas_dilution_time)
1144  call pmc_mpi_unpack_real_array(buffer, position, val%gas_dilution_rate)
1145  call pmc_mpi_unpack_real_array(buffer, position, val%aero_emission_time)
1146  call pmc_mpi_unpack_real_array(buffer, position, &
1147  val%aero_emission_rate_scale)
1148  call pmc_mpi_unpack_real_array(buffer, position, val%aero_dilution_time)
1149  call pmc_mpi_unpack_real_array(buffer, position, val%aero_dilution_rate)
1150  call pmc_mpi_unpack_integer(buffer, position, val%loss_function_type)
1151  call pmc_mpi_unpack_chamber(buffer, position, val%chamber)
1152  if (allocated(val%gas_emission)) deallocate(val%gas_emission)
1153  if (allocated(val%gas_background)) deallocate(val%gas_background)
1154  if (allocated(val%aero_emission)) deallocate(val%aero_emission)
1155  if (allocated(val%aero_background)) deallocate(val%aero_background)
1156  if (allocated(val%gas_emission_time)) then
1157  allocate(val%gas_emission(size(val%gas_emission_time)))
1158  do i = 1,size(val%gas_emission)
1159  call pmc_mpi_unpack_gas_state(buffer, position, val%gas_emission(i))
1160  end do
1161  end if
1162  if (allocated(val%gas_dilution_time)) then
1163  allocate(val%gas_background(size(val%gas_dilution_time)))
1164  do i = 1,size(val%gas_background)
1165  call pmc_mpi_unpack_gas_state(buffer, position, &
1166  val%gas_background(i))
1167  end do
1168  end if
1169  if (allocated(val%aero_emission_time)) then
1170  allocate(val%aero_emission(size(val%aero_emission_time)))
1171  do i = 1,size(val%aero_emission)
1172  call pmc_mpi_unpack_aero_dist(buffer, position, val%aero_emission(i))
1173  end do
1174  end if
1175  if (allocated(val%aero_dilution_time)) then
1176  allocate(val%aero_background(size(val%aero_dilution_time)))
1177  do i = 1,size(val%aero_background)
1178  call pmc_mpi_unpack_aero_dist(buffer, position, &
1179  val%aero_background(i))
1180  end do
1181  end if
1182  call assert(611542570, &
1183  position - prev_position <= pmc_mpi_pack_size_scenario(val))
1184 #endif
1185 
1186  end subroutine pmc_mpi_unpack_scenario
1187 
1188 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1189 
1190 end module pmc_scenario
pmc_aero_state::aero_state_zero
subroutine aero_state_zero(aero_state)
Resets an aero_state to have zero particles per bin.
Definition: aero_state.F90:350
pmc_scenario::scenario_contains_aero_mode_type
elemental logical function scenario_contains_aero_mode_type(scenario, aero_mode_type)
Whether any of the contained aerosol modes are of the given type.
Definition: scenario.F90:754
pmc_scenario::scenario_update_aero_binned
subroutine scenario_update_aero_binned(scenario, delta_t, env_state, old_env_state, bin_grid, aero_data, aero_binned)
Do emissions and background dilution from the environment for a binned aerosol distribution.
Definition: scenario.F90:351
pmc_gas_data::gas_data_t
Constant gas data.
Definition: gas_data.F90:35
pmc_mpi
Wrapper functions for MPI.
Definition: mpi.F90:13
pmc_scenario::scenario_loss_rate
real(kind=dp) function scenario_loss_rate(scenario, vol, density, aero_data, env_state)
Evaluate a loss rate function.
Definition: scenario.F90:400
pmc_spec_file::spec_file_close
subroutine spec_file_close(file)
Close a spec file.
Definition: spec_file.F90:135
pmc_gas_state::gas_state_ensure_nonnegative
subroutine gas_state_ensure_nonnegative(gas_state)
Set any negative values to zero.
Definition: gas_state.F90:167
pmc_aero_state::aero_state_add_aero_dist_sample
subroutine aero_state_add_aero_dist_sample(aero_state, aero_data, aero_dist, sample_prop, characteristic_factor, create_time, allow_doubling, allow_halving, n_part_add)
Generates a Poisson sample of an aero_dist, adding to aero_state, with the given sample proportion.
Definition: aero_state.F90:761
pmc_aero_dist::pmc_mpi_pack_aero_dist
subroutine pmc_mpi_pack_aero_dist(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: aero_dist.F90:426
pmc_scenario::scenario_loss_function_volume
integer, parameter scenario_loss_function_volume
Type code for a loss rate function proportional to particle volume.
Definition: scenario.F90:32
pmc_scenario
The scenario_t structure and associated subroutines.
Definition: scenario.F90:9
pmc_scenario::scenario_t
Scenario data.
Definition: scenario.F90:54
pmc_gas_data
The gas_data_t structure and associated subroutines.
Definition: gas_data.F90:9
pmc_chamber::pmc_mpi_pack_chamber
subroutine pmc_mpi_pack_chamber(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: chamber.F90:206
pmc_scenario::pmc_mpi_pack_scenario
subroutine pmc_mpi_pack_scenario(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: scenario.F90:1063
pmc_util::warn_assert_msg
subroutine warn_assert_msg(code, condition_ok, warning_msg)
Prints a warning message if condition_ok is false.
Definition: util.F90:60
pmc_util::die_msg
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:135
pmc_scenario::scenario_loss_function_constant
integer, parameter scenario_loss_function_constant
Type code for a constant loss function.
Definition: scenario.F90:30
pmc_scenario::scenario_loss_function_none
integer, parameter scenario_loss_function_none
Type code for a zero loss function.
Definition: scenario.F90:28
pmc_constants::dp
integer, parameter dp
Kind of a double precision real number.
Definition: constants.F90:12
pmc_aero_dist::aero_dist_interp_1d
subroutine aero_dist_interp_1d(aero_dist_list, time_list, rate_list, time, aero_dist, rate)
Determine the current aero_dist and rate by interpolating at the current time with the lists of aero_...
Definition: aero_dist.F90:177
pmc_spec_file
Reading formatted text input.
Definition: spec_file.F90:43
pmc_gas_state::gas_state_add_scaled
subroutine gas_state_add_scaled(gas_state, gas_state_delta, alpha)
Adds the given gas_state_delta scaled by alpha.
Definition: gas_state.F90:124
pmc_util::interp_1d
real(kind=dp) function interp_1d(x_vals, y_vals, x)
1D linear interpolation.
Definition: util.F90:627
pmc_bin_grid::bin_grid_size
elemental integer function bin_grid_size(bin_grid)
Return the number of bins in the grid, or -1 if the bin grid is not allocated.
Definition: bin_grid.F90:51
pmc_chamber::pmc_mpi_pack_size_chamber
integer function pmc_mpi_pack_size_chamber(val)
Determines the number of bytes required to pack the given value.
Definition: chamber.F90:187
pmc_scenario::pmc_mpi_unpack_scenario
subroutine pmc_mpi_unpack_scenario(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: scenario.F90:1122
pmc_env_state::env_state_t
Current environment state.
Definition: env_state.F90:29
pmc_util::assert
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:104
pmc_scenario::scenario_update_gas_state
subroutine scenario_update_gas_state(scenario, delta_t, env_state, old_env_state, gas_data, gas_state)
Do gas emissions and background dilution.
Definition: scenario.F90:209
pmc_scenario::scenario_update_aero_state
subroutine scenario_update_aero_state(scenario, delta_t, env_state, old_env_state, aero_data, aero_state, n_emit, n_dil_in, n_dil_out, allow_doubling, allow_halving)
Do emissions and background dilution for a particle aerosol distribution.
Definition: scenario.F90:262
pmc_aero_state
The aero_state_t structure and assocated subroutines.
Definition: aero_state.F90:9
pmc_mpi::pmc_mpi_pack_size_real_array
integer function pmc_mpi_pack_size_real_array(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:527
pmc_scenario::scenario_loss_rate_max
real(kind=dp) function scenario_loss_rate_max(scenario, vol, aero_data, env_state)
Compute and return the max loss rate function for a given volume.
Definition: scenario.F90:545
pmc_util::interp_linear_disc
real(kind=dp) function interp_linear_disc(x_1, x_n, n, i)
Linear interpolation over discrete indices.
Definition: util.F90:664
pmc_scenario::scenario_update_env_state
subroutine scenario_update_env_state(scenario, env_state, time)
Update time-dependent contents of the environment. scenario_init_env_state() should have been called ...
Definition: scenario.F90:134
pmc_spec_file::spec_file_t
An input file with extra data for printing messages.
Definition: spec_file.F90:59
pmc_scenario::scenario_loss_rate_bin_max
subroutine scenario_loss_rate_bin_max(scenario, bin_grid, aero_data, env_state, loss_max)
Compute an upper bound on the maximum kernel value for each bin. Value over_scale is multiplied to th...
Definition: scenario.F90:582
pmc_scenario::scenario_loss_alg_threshold
real(kind=dp), parameter scenario_loss_alg_threshold
Parameter to switch between algorithms for particle loss. A value of 0 will always use the naive algo...
Definition: scenario.F90:41
pmc_gas_state
The gas_state_t structure and associated subroutines.
Definition: gas_state.F90:9
pmc_gas_state::pmc_mpi_unpack_gas_state
subroutine pmc_mpi_unpack_gas_state(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: gas_state.F90:547
pmc_util::integer_to_string
character(len=pmc_util_convert_string_len) function integer_to_string(val)
Convert an integer to a string format.
Definition: util.F90:767
pmc_util::real_to_string
character(len=pmc_util_convert_string_len) function real_to_string(val)
Convert a real to a string format.
Definition: util.F90:799
pmc_aero_dist::pmc_mpi_unpack_aero_dist
subroutine pmc_mpi_unpack_aero_dist(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: aero_dist.F90:456
pmc_aero_dist
The aero_dist_t structure and associated subroutines.
Definition: aero_dist.F90:18
pmc_gas_state::gas_state_scale
subroutine gas_state_scale(gas_state, alpha)
Scale a gas state.
Definition: gas_state.F90:86
pmc_chamber::pmc_mpi_unpack_chamber
subroutine pmc_mpi_unpack_chamber(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: chamber.F90:233
pmc_aero_state::aero_state_remove_particle_with_info
subroutine aero_state_remove_particle_with_info(aero_state, i_part, aero_info)
Remove the given particle and record the removal.
Definition: aero_state.F90:410
pmc_aero_state::aero_state_total_particles
integer function aero_state_total_particles(aero_state, i_group, i_class)
Returns the total number of particles in an aerosol distribution.
Definition: aero_state.F90:288
pmc_chamber::chamber_t
Definition: chamber.F90:18
pmc_constants::const
type(const_t), save const
Fixed variable for accessing the constant's values.
Definition: constants.F90:75
pmc_chamber::chamber_loss_rate_sedi
real(kind=dp) function chamber_loss_rate_sedi(chamber, vol, density, aero_data, env_state)
Calculate the loss rate due to sedimentation in chamber. Based on Eq. 37b in Naumann 2003 J....
Definition: chamber.F90:118
pmc_spec_file::spec_file_open
subroutine spec_file_open(filename, file)
Open a spec file for reading.
Definition: spec_file.F90:112
pmc_env_state::env_state_sat_vapor_pressure
real(kind=dp) function env_state_sat_vapor_pressure(env_state)
Computes the current saturation vapor pressure (Pa).
Definition: env_state.F90:166
pmc_spec_file::spec_file_read_timed_real_array
subroutine spec_file_read_timed_real_array(file, name, times, vals)
Read an a time-indexed array of real data.
Definition: spec_file.F90:688
pmc_scenario::pmc_mpi_pack_size_scenario
integer function pmc_mpi_pack_size_scenario(val)
Determines the number of bytes required to pack the given value.
Definition: scenario.F90:1007
pmc_env_state
The env_state_t structure and associated subroutines.
Definition: env_state.F90:9
pmc_gas_state::gas_state_t
Current state of the gas mixing ratios in the system.
Definition: gas_state.F90:33
pmc_aero_data::aero_data_t
Aerosol material properties and associated data.
Definition: aero_data.F90:55
pmc_mpi::pmc_mpi_unpack_integer
subroutine pmc_mpi_unpack_integer(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1139
pmc_aero_state::aero_state_copy_weight
subroutine aero_state_copy_weight(aero_state_from, aero_state_to)
Copies weighting information for an aero_state.
Definition: aero_state.F90:171
pmc_aero_dist::aero_dist_t
A complete aerosol distribution, consisting of several modes.
Definition: aero_dist.F90:33
pmc_scenario::scenario_loss_function_drydep
integer, parameter scenario_loss_function_drydep
Type code for a loss rate function based on dry deposition.
Definition: scenario.F90:34
pmc_util
Common utility subroutines.
Definition: util.F90:9
pmc_scenario::scenario_init_env_state
subroutine scenario_init_env_state(scenario, env_state, time)
Initialize the time-dependent contents of the environment. Thereafter scenario_update_env_state() sho...
Definition: scenario.F90:111
pmc_gas_state::pmc_mpi_pack_gas_state
subroutine pmc_mpi_pack_gas_state(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: gas_state.F90:524
pmc_spec_file::spec_file_die_msg
subroutine spec_file_die_msg(code, file, msg)
Exit with an error message containing filename and line number.
Definition: spec_file.F90:74
pmc_aero_dist::pmc_mpi_pack_size_aero_dist
integer function pmc_mpi_pack_size_aero_dist(val)
Determines the number of bytes required to pack the given value.
Definition: aero_dist.F90:404
pmc_aero_data::aero_data_vol2diam
real(kind=dp) elemental function aero_data_vol2diam(aero_data, v)
Convert mass-equivalent volume (m^3) to geometric diameter (m).
Definition: aero_data.F90:117
pmc_spec_file::spec_file_read_string
subroutine spec_file_read_string(file, name, var)
Read a string from a spec file that must have a given name.
Definition: spec_file.F90:605
pmc_scenario::scenario_try_single_particle_loss
subroutine scenario_try_single_particle_loss(scenario, delta_t, aero_data, aero_state, env_state, i_part, over_prob)
Test a candidate particle to see if it should be removed, and remove if necessary....
Definition: scenario.F90:713
pmc_chamber
Definition: chamber.F90:8
pmc_gas_state::gas_state_interp_1d
subroutine gas_state_interp_1d(gas_state_list, time_list, rate_list, time, gas_state, rate)
Determine the current gas_state and rate by interpolating at the current time with the lists of gas_s...
Definition: gas_state.F90:239
pmc_mpi::pmc_mpi_pack_size_integer
integer function pmc_mpi_pack_size_integer(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:345
pmc_mpi::pmc_mpi_pack_integer
subroutine pmc_mpi_pack_integer(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:711
pmc_aero_data::aero_data_rad2vol
real(kind=dp) elemental function aero_data_rad2vol(aero_data, r)
Convert geometric radius (m) to mass-equivalent volume (m^3).
Definition: aero_data.F90:132
pmc_aero_data
The aero_data_t structure and associated subroutines.
Definition: aero_data.F90:9
pmc_mpi::pmc_mpi_pack_real_array
subroutine pmc_mpi_pack_real_array(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:927
pmc_bin_grid::bin_grid_t
1D grid, either logarithmic or linear.
Definition: bin_grid.F90:33
pmc_scenario::scenario_loss_rate_dry_dep
real(kind=dp) function scenario_loss_rate_dry_dep(vol, density, aero_data, env_state)
Compute and return the dry deposition rate for a given particle. All equations used here are written ...
Definition: scenario.F90:446
pmc_chamber::chamber_loss_rate_wall
real(kind=dp) function chamber_loss_rate_wall(chamber, vol, aero_data, env_state)
Calculate the loss rate due to wall diffusion in chamber. Based on Eq. 37a in Naumann 2003 J....
Definition: chamber.F90:91
pmc_gas_state::gas_state_molar_conc_to_ppb
subroutine gas_state_molar_conc_to_ppb(gas_state, env_state)
Convert (mol m^{-3}) to (ppb).
Definition: gas_state.F90:181
pmc_scenario::scenario_particle_loss
subroutine scenario_particle_loss(scenario, delta_t, aero_data, aero_state, env_state)
Performs stochastic particle loss for one time-step. If a particle i_part has a scenario_loss_rate() ...
Definition: scenario.F90:628
pmc_scenario::scenario_loss_function_chamber
integer, parameter scenario_loss_function_chamber
Type code for a loss rate function for chamber experiments.
Definition: scenario.F90:36
pmc_aero_dist::aero_dist_contains_aero_mode_type
elemental logical function aero_dist_contains_aero_mode_type(aero_dist, aero_mode_type)
Whether any of the modes are of the given type.
Definition: aero_dist.F90:155
pmc_gas_state::pmc_mpi_pack_size_gas_state
integer function pmc_mpi_pack_size_gas_state(val)
Determines the number of bytes required to pack the given value.
Definition: gas_state.F90:511
pmc_mpi::pmc_mpi_unpack_real_array
subroutine pmc_mpi_unpack_real_array(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1356
pmc_scenario::scenario_loss_function_invalid
integer, parameter scenario_loss_function_invalid
Type code for an undefined or invalid loss function.
Definition: scenario.F90:26
pmc_aero_state::aero_state_t
The current collection of aerosol particles.
Definition: aero_state.F90:69
pmc_aero_state::aero_state_sort
subroutine aero_state_sort(aero_state, aero_data, bin_grid, all_procs_same)
Sorts the particles if necessary.
Definition: aero_state.F90:3219
pmc_rand::pmc_random
real(kind=dp) function pmc_random()
Returns a random number between 0 and 1.
Definition: rand.F90:139
pmc_aero_state::aero_state_sample_particles
subroutine aero_state_sample_particles(aero_state_from, aero_state_to, aero_data, sample_prob, removal_action)
Generates a random sample by removing particles from aero_state_from and adding them to aero_state_to...
Definition: aero_state.F90:858