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