PartMC  2.8.0
coagulation.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2012 Nicole Riemer
2 ! Copyright (C) 2005-2015 Matthew West
3 ! Copyright (C) 2015 Jeffrey Curtis
4 ! Licensed under the GNU General Public License version 2 or (at your
5 ! option) any later version. See the file COPYING for details.
6 
7 !> \file
8 !> The pmc_coagulation module.
9 
10 !> Aerosol particle coagulation.
12 
13  use pmc_bin_grid
14  use pmc_aero_data
15  use pmc_util
16  use pmc_stats
17  use pmc_env_state
18  use pmc_aero_state
19  use pmc_aero_weight
21  use pmc_mpi
22  use pmc_coag_kernel
23  use pmc_aero_sorted
24 #ifdef PMC_USE_MPI
25  use mpi
26 #endif
27 
28  !> Minimum number of coagulation events per large particle for which
29  !> accelerated coagulation is used.
30  real(kind=dp), parameter :: coag_accel_n_event = 1d0
31  !> Maximum allowed coefficient-of-variation due to undersampling in
32  !> accelerated coagulation.
33  real(kind=dp), parameter :: coag_accel_max_cv = 0.1d0
34  !> Maximum allowable growth factor of a target particle volume within one
35  !> timestep when using accelerated coagulation.
36  real(kind=dp), parameter :: max_allowable_growth_factor = 1.5d0
37 
38 contains
39 
40 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
41 
42  !> Do coagulation for time del_t.
43  subroutine mc_coag(coag_kernel_type, env_state, aero_data, aero_state, &
44  del_t, tot_n_samp, tot_n_coag)
45 
46  !> Coagulation kernel type.
47  integer, intent(in) :: coag_kernel_type
48  !> Environment state.
49  type(env_state_t), intent(in) :: env_state
50  !> Aerosol data.
51  type(aero_data_t), intent(in) :: aero_data
52  !> Aerosol state.
53  type(aero_state_t), intent(inout) :: aero_state
54  !> Timestep for coagulation.
55  real(kind=dp) :: del_t
56  !> Total number of samples tested.
57  integer, intent(out) :: tot_n_samp
58  !> Number of coagulation events.
59  integer, intent(out) :: tot_n_coag
60 
61  integer :: b1, b2, c1, c2, b2_start
62 
63  call aero_state_sort(aero_state, aero_data)
64  if (.not. aero_state%aero_sorted%coag_kernel_bounds_valid) then
65  call est_k_minmax_binned_unweighted(aero_state%aero_sorted%bin_grid, &
66  coag_kernel_type, aero_data, env_state, &
67  aero_state%aero_sorted%coag_kernel_min, &
68  aero_state%aero_sorted%coag_kernel_max)
69  aero_state%aero_sorted%coag_kernel_bounds_valid = .true.
70  end if
71 
72  tot_n_samp = 0
73  tot_n_coag = 0
74  do c1 = 1,aero_sorted_n_class(aero_state%aero_sorted)
75  do c2 = 1,c1
76  do b1 = 1,aero_sorted_n_bin(aero_state%aero_sorted)
77  if (integer_varray_n_entry( &
78  aero_state%aero_sorted%size_class%inverse(b1, c1)) == 0) &
79  cycle
80  if (c1 == c2) then
81  b2_start = b1
82  else
83  b2_start = 1
84  end if
85  do b2 = b2_start,aero_sorted_n_bin(aero_state%aero_sorted)
86  if (integer_varray_n_entry( &
87  aero_state%aero_sorted%size_class%inverse(b2, c2)) == 0) &
88  cycle
89  call mc_coag_for_bin(coag_kernel_type, env_state, aero_data, &
90  aero_state, del_t, tot_n_samp, tot_n_coag, b1, b2, c1, c2)
91  end do
92  end do
93  end do
94  end do
95 
96  end subroutine mc_coag
97 
98 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
99 
100  !> Do coagulation for time del_t for the given bins.
101  subroutine mc_coag_for_bin(coag_kernel_type, env_state, aero_data, &
102  aero_state, del_t, tot_n_samp, tot_n_coag, b1, b2, c1, c2)
103 
104  !> Coagulation kernel type.
105  integer, intent(in) :: coag_kernel_type
106  !> Environment state.
107  type(env_state_t), intent(in) :: env_state
108  !> Aerosol data.
109  type(aero_data_t), intent(in) :: aero_data
110  !> Aerosol state.
111  type(aero_state_t), intent(inout) :: aero_state
112  !> Timestep for coagulation.
113  real(kind=dp) :: del_t
114  !> Total number of samples tested.
115  integer, intent(inout) :: tot_n_samp
116  !> Number of coagulation events.
117  integer, intent(inout) :: tot_n_coag
118  !> First bin number.
119  integer, intent(in) :: b1
120  !> Second bin number.
121  integer, intent(in) :: b2
122  !> First weight class.
123  integer, intent(in) :: c1
124  !> Second weight class.
125  integer, intent(in) :: c2
126 
127  logical :: per_particle_coag_succeeded
128  integer :: cc
129  real(kind=dp) :: f_max, k_max
130 
131  cc = coag_dest_class(aero_state%awa, aero_data, &
132  aero_state%aero_sorted%bin_grid, b1, b2, c1, c2)
133  call max_coag_num_conc_factor(aero_state%awa, aero_data, &
134  aero_state%aero_sorted%bin_grid, b1, b2, c1, c2, cc, f_max)
135  k_max = aero_state%aero_sorted%coag_kernel_max(b1, b2) * f_max
136 
137  call try_per_particle_coag(coag_kernel_type, k_max, env_state, aero_data, &
138  aero_state, del_t, tot_n_samp, tot_n_coag, b1, b2, c1, c2, cc, &
139  per_particle_coag_succeeded)
140  if (per_particle_coag_succeeded) return
141 
142  call per_set_coag(coag_kernel_type, k_max, env_state, aero_data, &
143  aero_state, del_t, tot_n_samp, tot_n_coag, b1, b2, c1, c2, cc)
144 
145  end subroutine mc_coag_for_bin
146 
147 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
148 
149  !> Attempt per-particle coagulation.
150  subroutine try_per_particle_coag(coag_kernel_type, k_max, env_state, &
151  aero_data, aero_state, del_t, tot_n_samp, tot_n_coag, b1, b2, c1, c2, &
152  cc, per_particle_coag_succeeded)
153 
154  !> Coagulation kernel type.
155  integer, intent(in) :: coag_kernel_type
156  !> Maximum coagulation kernel (s^{-1} m^3).
157  real(kind=dp), intent(in) :: k_max
158  !> Environment state.
159  type(env_state_t), intent(in) :: env_state
160  !> Aerosol data.
161  type(aero_data_t), intent(in) :: aero_data
162  !> Aerosol state.
163  type(aero_state_t), intent(inout) :: aero_state
164  !> Timestep for coagulation.
165  real(kind=dp) :: del_t
166  !> Total number of samples tested.
167  integer, intent(inout) :: tot_n_samp
168  !> Number of coagulation events.
169  integer, intent(inout) :: tot_n_coag
170  !> First bin number.
171  integer, intent(in) :: b1
172  !> Second bin number.
173  integer, intent(in) :: b2
174  !> First weight class.
175  integer, intent(in) :: c1
176  !> Second weight class.
177  integer, intent(in) :: c2
178  !> Coagulated weight class.
179  integer, intent(in) :: cc
180  !> Whether we succeeded in doing per-particle coag.
181  logical, intent(inout) :: per_particle_coag_succeeded
182 
183  logical :: correct_weight_ordering
184  integer :: target_unif_entry, target_part, n_samp, n_coag, n_remove, bt, bs
185  integer :: ct, cs
186  real(kind=dp) :: n_source_per_target, accept_factor
187  real(kind=dp) :: min_target_vol, max_source_vol, max_new_target_vol
188  real(kind=dp) :: max_target_growth_factor
189  type(aero_particle_t) :: target_particle, source_particle
190 
191  call determine_target_and_source(aero_state%awa, &
192  aero_state%aero_sorted%bin_grid, b1, b2, c1, c2, cc, bt, bs, ct, cs, &
193  correct_weight_ordering)
194  if (.not. correct_weight_ordering) then
195  per_particle_coag_succeeded = .false.
196  return
197  end if
198 
199  call compute_n_source(integer_varray_n_entry( &
200  aero_state%aero_sorted%size_class%inverse(bs, cs)), k_max, del_t, &
201  n_source_per_target, accept_factor)
202  if (n_source_per_target < coag_accel_n_event) then
203  per_particle_coag_succeeded = .false.
204  return
205  end if
206 
207  min_target_vol = aero_data_rad2vol(aero_data, &
208  aero_state%aero_sorted%bin_grid%edges(bt))
209  max_source_vol = aero_data_rad2vol(aero_data, &
210  aero_state%aero_sorted%bin_grid%edges(bs+1))
211  max_new_target_vol = min_target_vol + max_source_vol * n_source_per_target
212  ! check for unacceptably large volume change
213  max_target_growth_factor = max_new_target_vol / min_target_vol
214  if (max_target_growth_factor > max_allowable_growth_factor) then
215  per_particle_coag_succeeded = .false.
216  return
217  end if
218 
219  ! work backwards to avoid particle movement issues
220  do target_unif_entry = integer_varray_n_entry( &
221  aero_state%aero_sorted%size_class%inverse(bt, ct)),1,-1
222  target_part = aero_state%aero_sorted%size_class%inverse(bt, &
223  ct)%entry(target_unif_entry)
224  ! need to copy coag_particle as the underlying storage may be
225  ! rearranged due to removals
226  target_particle = aero_state%apa%particle(target_part)
227  call sample_source_particle(aero_state, aero_data, env_state, &
228  coag_kernel_type, bs, cs, target_particle, n_source_per_target, &
229  accept_factor, n_samp, n_coag, n_remove, source_particle)
230  if (n_coag > 0) then
231  call coag_target_with_source(aero_state, aero_data, bt, ct, &
232  target_unif_entry, source_particle, cc)
233  end if
234  tot_n_samp = tot_n_samp + n_samp
235  tot_n_coag = tot_n_coag + n_coag
236  ! we discard n_remove information at present
237  end do
238 
239  per_particle_coag_succeeded = .true.
240 
241  end subroutine try_per_particle_coag
242 
243 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
244 
245  !> Determine the source and target particle bin/group for
246  !> per-particle coagulation, if possible.
247  subroutine determine_target_and_source(aero_weight_array, bin_grid, b1, &
248  b2, c1, c2, cc, bt, bs, ct, cs, correct_weight_ordering)
249 
250  !> Aero weight array.
251  type(aero_weight_array_t), intent(in) :: aero_weight_array
252  !> Bin grid.
253  type(bin_grid_t), intent(in) :: bin_grid
254  !> First bin number.
255  integer, intent(in) :: b1
256  !> Second bin number.
257  integer, intent(in) :: b2
258  !> First weight class.
259  integer, intent(in) :: c1
260  !> Second weight class.
261  integer, intent(in) :: c2
262  !> Coagulated weight class.
263  integer, intent(in) :: cc
264  !> Target bin number.
265  integer, intent(out) :: bt
266  !> Source bin number.
267  integer, intent(out) :: bs
268  !> Target weight class.
269  integer, intent(out) :: ct
270  !> Source weight class.
271  integer, intent(out) :: cs
272  !> Whether the weight ordering is correct for per-particle coagulation.
273  logical, intent(out) :: correct_weight_ordering
274 
275  real(kind=dp) :: nc1_min, nc1_max, nc2_min, nc2_max
276  logical :: monotone_increasing, monotone_decreasing
277 
278  call aero_weight_array_check_monotonicity(aero_weight_array, &
279  monotone_increasing, monotone_decreasing)
280  if (.not. monotone_decreasing) then
281  correct_weight_ordering = .false.
282  return
283  end if
284 
285  call aero_weight_array_minmax_num_conc(aero_weight_array, c1, &
286  bin_grid%edges(b1), bin_grid%edges(b1 + 1), nc1_min, nc1_max)
287  call aero_weight_array_minmax_num_conc(aero_weight_array, c2, &
288  bin_grid%edges(b2), bin_grid%edges(b2 + 1), nc2_min, nc2_max)
289 
290  ! we have already confirmed monotone_decreasing weights above
291  correct_weight_ordering = .false.
292  if ((nc1_max < nc2_min) .and. (cc == c1)) then
293  bt = b1
294  bs = b2
295  ct = c1
296  cs = c2
297  correct_weight_ordering = .true.
298  end if
299  if ((nc2_max < nc1_min) .and. (cc == c2)) then
300  bt = b2
301  bs = b1
302  ct = c2
303  cs = c1
304  correct_weight_ordering = .true.
305  end if
306 
307  end subroutine determine_target_and_source
308 
309 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
310 
311  subroutine compute_n_source(n_part, k_max, del_t, n_source_per_target, &
312  accept_factor)
313 
314  !> Number of particles available as partners.
315  integer, intent(in) :: n_part
316  !> Maximum coagulation kernel (s^{-1} m^3).
317  real(kind=dp), intent(in) :: k_max
318  !> Timestep (s).
319  real(kind=dp), intent(in) :: del_t
320  !> Mean number of source particles per target particle.
321  real(kind=dp), intent(out) :: n_source_per_target
322  !> Accept factor for samples.
323  real(kind=dp), intent(out) :: accept_factor
324 
325  n_source_per_target = k_max * del_t * real(n_part, kind=dp)
326  accept_factor = 1d0 / k_max
327 
328  end subroutine compute_n_source
329 
330 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
331 
332  !> Sample coagulation partners for a single coagulation event.
333  subroutine sample_source_particle(aero_state, aero_data, env_state, &
334  coag_kernel_type, bs, cs, coag_particle, n_samp_mean, &
335  accept_factor, n_samp, n_coag, n_remove, source_particle)
336 
337  !> Aerosol state.
338  type(aero_state_t), intent(inout) :: aero_state
339  !> Aerosol data.
340  type(aero_data_t), intent(in) :: aero_data
341  !> Environment state.
342  type(env_state_t), intent(in) :: env_state
343  !> Coagulation kernel type.
344  integer, intent(in) :: coag_kernel_type
345  !> Bin to sample particles from.
346  integer, intent(in) :: bs
347  !> Weight class to sample particles from.
348  integer, intent(in) :: cs
349  !> Aerosol particle that coagulation will be with.
350  type(aero_particle_t), intent(in) :: coag_particle
351  !> Mean number of samples to use.
352  real(kind=dp), intent(in) :: n_samp_mean
353  !> Probability factor of accepting samples.
354  real(kind=dp), intent(in) :: accept_factor
355  !> Number of samples used.
356  integer, intent(out) :: n_samp
357  !> Number of coagulations performed.
358  integer, intent(out) :: n_coag
359  !> Number of particles removed.
360  integer, intent(out) :: n_remove
361  !> Sampled average coagulation partner particle.
362  type(aero_particle_t), intent(inout) :: source_particle
363 
364  real(kind=dp) :: prob_remove_i, prob_remove_source_max
365  real(kind=dp) :: prob_coag, prob_coag_tot, prob_coag_mean
366  real(kind=dp) :: num_conc_i, num_conc_source_min, num_conc_target, k
367  real(kind=dp) :: vol_sq(aero_data_n_spec(aero_data))
368  real(kind=dp) :: vol_mean(aero_data_n_spec(aero_data))
369  real(kind=dp) :: vol_cv(aero_data_n_spec(aero_data)), vol_cv_max
370  real(kind=dp) :: mean_95_conf_cv
371  integer :: n_samp_remove, n_samp_extra, n_samp_total, n_avg, i_samp
372  integer :: i_unif_entry, i_part, new_bin, ct, i_comp, i
373  integer(kind=8) :: target_id
374  integer :: n_comp_new, n_comp, sample_size
375  integer, allocatable :: sample(:)
376  type(aero_component_t), allocatable :: new_aero_component(:)
377  type(aero_info_t) :: aero_info
378 
379  if (integer_varray_n_entry( &
380  aero_state%aero_sorted%size_class%inverse(bs, cs)) == 0) then
381  n_samp = 0
382  n_remove = 0
383  n_coag = 0
384  return
385  end if
386 
387  num_conc_target = aero_weight_array_num_conc(aero_state%awa, &
388  coag_particle, aero_data)
389  target_id = coag_particle%id
390  ct = coag_particle%weight_class
391 
392  num_conc_source_min = aero_weight_array_num_conc_at_radius( &
393  aero_state%awa, cs, aero_state%aero_sorted%bin_grid%edges(bs))
394  prob_remove_source_max = num_conc_target / num_conc_source_min
395  call assert(653606684, prob_remove_source_max <= 1d0)
396 
397  n_samp_remove = rand_poisson(prob_remove_source_max * n_samp_mean)
398  n_samp_extra = rand_poisson((1d0 - prob_remove_source_max) * n_samp_mean)
399  n_samp_total = n_samp_remove + n_samp_extra
400 
401  n_avg = 0
402  n_samp = 0
403  n_remove = 0
404  prob_coag_tot = 0d0
405  vol_sq = 0d0
406 
407  call aero_particle_zero(source_particle, aero_data)
408 
409  ! FIXME: Can't we just do n_samp = 1,n_samp_total and shift tests
410  ! to the end?
411  do i_samp = 1,n_samp_total
412  if (integer_varray_n_entry( &
413  aero_state%aero_sorted%size_class%inverse(bs, cs)) == 0) exit
414  if ((n_samp > n_samp_remove) .and. (n_avg >= 2)) then
415  vol_mean = source_particle%vol / real(n_avg, kind=dp)
416  where(vol_mean > 0d0) &
417  vol_cv = sqrt(vol_sq / real(n_avg, kind=dp) - (vol_mean)**2) &
418  / vol_mean
419  vol_cv_max = maxval(vol_cv)
420  mean_95_conf_cv = vol_cv_max / sqrt(real(n_avg, kind=dp)) &
421  * student_t_95_coeff(n_avg)
422  ! FIXME: We are using just the max of the diagonal of the
423  ! covariance matrix. Is this well-justified?
424  if (mean_95_conf_cv < coag_accel_max_cv) exit
425  end if
426  n_samp = n_samp + 1
427  ! FIXME: We are sampling with replacement. Is this a problem?
428  i_unif_entry = pmc_rand_int(integer_varray_n_entry( &
429  aero_state%aero_sorted%size_class%inverse(bs, cs)))
430  i_part = aero_state%aero_sorted%size_class%inverse(bs, &
431  cs)%entry(i_unif_entry)
432  ! re-get j_part as particle ordering may be changing
433  call num_conc_weighted_kernel(coag_kernel_type, &
434  aero_state%apa%particle(i_part), coag_particle, cs, ct, ct, &
435  aero_data, aero_state%awa, env_state, k)
436  prob_coag = k * accept_factor
437  prob_coag_tot = prob_coag_tot + prob_coag
438  if (pmc_random() < prob_coag) then
439  n_avg = n_avg + 1
440  call aero_particle_coagulate(source_particle, &
441  aero_state%apa%particle(i_part), source_particle)
442  vol_sq = vol_sq + aero_state%apa%particle(i_part)%vol**2
443  if (i_samp <= n_samp_remove) then
444  num_conc_i = aero_weight_array_num_conc(aero_state%awa, &
445  aero_state%apa%particle(i_part), aero_data)
446  prob_remove_i = num_conc_target / num_conc_i
447  if (pmc_random() < prob_remove_i / prob_remove_source_max) then
448  n_remove = n_remove + 1
449  aero_info%id = aero_state%apa%particle(i_part)%id
450  aero_info%action = aero_info_coag
451  aero_info%other_id = target_id
452  call aero_state_remove_particle_with_info(aero_state, &
453  i_part, aero_info)
454  end if
455  end if
456  end if
457  end do
458 
459  if (n_avg == 0) then
460  n_coag = 0
461  return
462  end if
463 
464  prob_coag_mean = prob_coag_tot / real(n_samp, kind=dp) ! note not n_avg
465  call warn_assert_msg(383983415, prob_coag_mean <= 1d0, &
466  "kernel upper bound estimation is too tight")
467  n_coag = rand_binomial(n_samp_total, prob_coag_mean)
468  source_particle%vol = source_particle%vol &
469  * (real(n_coag, kind=dp) / real(n_avg, kind=dp))
470  source_particle%n_primary_parts = prob_round( &
471  source_particle%n_primary_parts * (real(n_coag, kind=dp) / n_avg))
472 
473  n_comp = aero_particle_n_components(source_particle)
474  n_comp_new = min(source_particle%n_primary_parts, &
475  max_aero_component_size)
476  if (n_comp_new < n_comp) then
477  allocate(new_aero_component(n_comp_new))
478  call sample_without_replacement(n_comp_new, n_comp, sample)
479  do i = 1,n_comp_new
480  new_aero_component(i) = source_particle%component(sample(i))
481  end do
482  call move_alloc(new_aero_component, source_particle%component)
483  else
484  if (n_comp < max_aero_component_size) then
485  ! make X copies plus sampling without replacment for leftovers
486  i_comp = 1
487  allocate(new_aero_component(n_comp_new))
488  do while (i_comp <= n_comp_new)
489  if (i_comp + n_comp <= n_comp_new) then
490  ! Append whole copy of source starting at i
491  do i = 1,n_comp
492  new_aero_component(i_comp+i-1) = source_particle%component(i)
493  end do
494  i_comp = i_comp + n_comp
495  else
496  ! sample a subset to fill in remaining entries
497  sample_size = n_comp_new - i_comp + 1
498  call sample_without_replacement(sample_size, n_comp, sample)
499  do i = 1,sample_size
500  new_aero_component(i_comp+i-1) = &
501  source_particle%component(sample(i))
502  end do
503  i_comp = n_comp_new + 1
504  end if
505  end do
506  call move_alloc(new_aero_component, source_particle%component)
507  end if
508  end if
509 
510  n_comp = aero_particle_n_components(source_particle)
511  call assert_msg(808144130, source_particle%n_primary_parts >= &
512  n_comp, 'n_primary_parts = ' &
513  // trim(integer_to_string(source_particle%n_primary_parts)) &
514  // ' is less than n_components = ' &
515  // trim(integer_to_string(n_comp)))
516 
517  end subroutine sample_source_particle
518 
519 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
520 
521  !> Coagulate a sampled source particle with a target particle.
522  subroutine coag_target_with_source(aero_state, aero_data, bt, ct, &
523  target_unif_entry, source_particle, cc)
524 
525  !> Aerosol state.
526  type(aero_state_t), intent(inout) :: aero_state
527  !> Aerosol data.
528  type(aero_data_t), intent(in) :: aero_data
529  !> Bin of coagulating particle.
530  integer, intent(in) :: bt
531  !> Weight class of coagulating particle.
532  integer, intent(in) :: ct
533  !> Entry-in-bin of coagulating particle.
534  integer, intent(in) :: target_unif_entry
535  !> Sampled particle to coagulate with.
536  type(aero_particle_t), intent(in) :: source_particle
537  !> Weight class for coagulated particle.
538  integer, intent(in) :: cc
539 
540  integer :: target_part, new_bin, new_group
541  integer(kind=8) :: target_id
542  real(kind=dp) :: old_num_conc_target, new_num_conc_target
543 
544  target_part = aero_state%aero_sorted%size_class%inverse(bt, &
545  ct)%entry(target_unif_entry)
546  target_id = aero_state%apa%particle(target_part)%id
547  old_num_conc_target = aero_weight_array_num_conc(aero_state%awa, &
548  aero_state%apa%particle(target_part), aero_data)
549  call aero_particle_coagulate(aero_state%apa%particle(target_part), &
550  source_particle, aero_state%apa%particle(target_part))
551  aero_state%apa%particle(target_part)%id = target_id
552  ! assign to a randomly chosen group
553  new_group = aero_weight_array_rand_group(aero_state%awa, cc, &
554  aero_particle_radius(aero_state%apa%particle(target_part), &
555  aero_data))
556  call aero_particle_set_weight(aero_state%apa%particle(target_part), &
557  new_group, cc)
558  ! fix bin due to composition changes
559  new_bin = aero_sorted_particle_in_bin(aero_state%aero_sorted, &
560  aero_state%apa%particle(target_part), aero_data)
561  if ((new_bin < 1) &
562  .or. (new_bin > bin_grid_size(aero_state%aero_sorted%bin_grid))) then
563  call die_msg(765620746, "particle outside of bin_grid: " &
564  // "try reducing the timestep del_t")
565  end if
566  call aero_sorted_move_particle(aero_state%aero_sorted, target_part, &
567  new_bin, new_group, cc)
568  ! Adjust particle number to account for weight changes
569  ! bt/target_group/target_entry are invalid, but
570  ! target_part is still good. We are treating all particles in all
571  ! groups together, and randomly reassigning between groups above,
572  ! so here we can't use aero_state_reweight_particle(), as that
573  ! assumes we are staying in the same weight group.
574  new_num_conc_target = aero_weight_array_num_conc(aero_state%awa, &
575  aero_state%apa%particle(target_part), aero_data)
576  call aero_state_dup_particle(aero_state, aero_data, target_part, &
577  old_num_conc_target / new_num_conc_target, random_weight_group=.true.)
578  ! we should only be doing this for decreasing weights
579  call assert(654300924, aero_state%apa%particle(target_part)%id &
580  == target_id)
581 
582  end subroutine coag_target_with_source
583 
584 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
585 
586  !> Do set-wise coagulation.
587  subroutine per_set_coag(coag_kernel_type, k_max, env_state, aero_data, &
588  aero_state, del_t, tot_n_samp, tot_n_coag, b1, b2, c1, c2, cc)
589 
590  !> Coagulation kernel type.
591  integer, intent(in) :: coag_kernel_type
592  !> Maximum coagulation kernel (s^{-1} m^3).
593  real(kind=dp), intent(in) :: k_max
594  !> Environment state.
595  type(env_state_t), intent(in) :: env_state
596  !> Aerosol data.
597  type(aero_data_t), intent(in) :: aero_data
598  !> Aerosol state.
599  type(aero_state_t), intent(inout) :: aero_state
600  !> Timestep for coagulation.
601  real(kind=dp) :: del_t
602  !> Total number of samples tested.
603  integer, intent(inout) :: tot_n_samp
604  !> Number of coagulation events.
605  integer, intent(inout) :: tot_n_coag
606  !> First bin number.
607  integer, intent(in) :: b1
608  !> Second bin number.
609  integer, intent(in) :: b2
610  !> First weight class.
611  integer, intent(in) :: c1
612  !> Second weight class.
613  integer, intent(in) :: c2
614  !> Coagulated weight class.
615  integer, intent(in) :: cc
616 
617  real(kind=dp) :: n_samp_mean, accept_factor
618  integer :: i_samp, n_samp, n1, n2
619  logical :: did_coag
620 
621  n1 = integer_varray_n_entry( &
622  aero_state%aero_sorted%size_class%inverse(b1, c1))
623  n2 = integer_varray_n_entry( &
624  aero_state%aero_sorted%size_class%inverse(b2, c2))
625  call compute_n_samp(n1, n2, ((b1 == b2) .and. (c1 == c2)), k_max, del_t, &
626  n_samp_mean, n_samp, accept_factor)
627  tot_n_samp = tot_n_samp + n_samp
628 
629  do i_samp = 1,n_samp
630  ! check we still have enough particles to coagulate
631  n1 = integer_varray_n_entry( &
632  aero_state%aero_sorted%size_class%inverse(b1, c1))
633  n2 = integer_varray_n_entry( &
634  aero_state%aero_sorted%size_class%inverse(b2, c2))
635  if (((n1 < 2) .and. (b1 == b2) .and. (c1 == c2)) &
636  .or. (n1 < 1) .or. (n2 < 1)) &
637  exit
638  call maybe_coag_pair(env_state, aero_data, aero_state, b1, b2, c1, c2, &
639  cc, coag_kernel_type, accept_factor, did_coag)
640  if (did_coag) tot_n_coag = tot_n_coag + 1
641  end do
642 
643  end subroutine per_set_coag
644 
645 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
646 
647  !> Compute the number of samples required for the pair of bins.
648  subroutine compute_n_samp(ni, nj, same_bin, k_max, del_t, n_samp_mean, &
649  n_samp, accept_factor)
650 
651  !> Number particles in first bin.
652  integer, intent(in) :: ni
653  !> Number particles in second bin.
654  integer, intent(in) :: nj
655  !> Whether first bin is second bin.
656  logical, intent(in) :: same_bin
657  !> Maximum kernel value (s^{-1} m^3).
658  real(kind=dp), intent(in) :: k_max
659  !> Timestep (s).
660  real(kind=dp), intent(in) :: del_t
661  !> Mean number of samples per timestep.
662  real(kind=dp), intent(out) :: n_samp_mean
663  !> Number of samples per timestep.
664  integer, intent(out) :: n_samp
665  !> Scale factor for accept probability (1).
666  real(kind=dp), intent(out) :: accept_factor
667 
668  real(kind=dp) :: r_samp
669  real(kind=dp) :: n_possible ! use real(kind=dp) to avoid integer overflow
670  ! could use integer*8 or integer(kind = 8)
671  ! or di = selected_int_kind(18), integer(kind=di)
672  ! to represent 10^{-18} to 10^{18}
673 
674  if (same_bin) then
675  ! don't change this to ni * (ni - 1) as the ni/nj distinction
676  ! is important for coagulation_dist, which also calls this
677  n_possible = real(ni, kind=dp) * (real(nj, kind=dp) - 1d0) / 2d0
678  else
679  n_possible = real(ni, kind=dp) * real(nj, kind=dp)
680  end if
681 
682  r_samp = k_max * del_t
683  n_samp_mean = r_samp * n_possible
684  n_samp = rand_poisson(n_samp_mean)
685  accept_factor = 1d0 / k_max
686 
687  ! possible variants:
688  ! A: accept_factor = 1d0 / k_max
689  ! B: accept_factor = del_t * n_possible / (n_samp * V)
690  ! timings of test suite as of 2010-12-22T17:12:14-0600:
691  ! A with n_samp = prob_round(n_samp_mean):
692  ! 159.82 162.18 156.28
693  ! A with n_samp = rand_poisson(n_samp_mean):
694  ! 151.93 158.38 174.74 157.65
695  ! B with n_samp = ceiling(n_samp_mean):
696  ! 196.06 200.23 192.41
697  ! B with n_samp = ceiling(n_samp_mean + 1 * sqrt(n_samp_mean)):
698  ! 189.78 211.12 195.19
699  ! B with n_samp = ceiling(n_samp_mean + 3 * sqrt(n_samp_mean)):
700  ! 214.60 201.25 203.55
701 
702  end subroutine compute_n_samp
703 
704 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
705 
706  !> Choose a random pair for potential coagulation and test its
707  !> probability of coagulation. If it happens, do the coagulation and
708  !> update all structures.
709  !!
710  !! The probability of a coagulation will be taken as <tt>(kernel /
711  !! k_max)</tt>.
712  subroutine maybe_coag_pair(env_state, aero_data, aero_state, b1, b2, &
713  c1, c2, cc, coag_kernel_type, accept_factor, did_coag)
714 
715  !> Environment state.
716  type(env_state_t), intent(in) :: env_state
717  !> Aerosol data.
718  type(aero_data_t), intent(in) :: aero_data
719  !> Aerosol state.
720  type(aero_state_t), intent(inout) :: aero_state
721  !> Bin of first particle.
722  integer, intent(in) :: b1
723  !> Bin of second particle.
724  integer, intent(in) :: b2
725  !> First weight class.
726  integer, intent(in) :: c1
727  !> Second weight class.
728  integer, intent(in) :: c2
729  !> Coagulated weight class.
730  integer, intent(in) :: cc
731  !> Coagulation kernel type.
732  integer, intent(in) :: coag_kernel_type
733  !> Scale factor for accept probability (1).
734  real(kind=dp), intent(in) :: accept_factor
735  !> Whether a coagulation occured.
736  logical, intent(out) :: did_coag
737 
738  integer :: i1, i2
739  integer :: p1, p2
740  real(kind=dp) :: p, k
741 
742  call find_rand_pair(aero_state%aero_sorted, b1, b2, c1, c2, i1, i2)
743  p1 = aero_state%aero_sorted%size_class%inverse(b1, c1)%entry(i1)
744  p2 = aero_state%aero_sorted%size_class%inverse(b2, c2)%entry(i2)
745  call num_conc_weighted_kernel(coag_kernel_type, &
746  aero_state%apa%particle(p1), aero_state%apa%particle(p2), &
747  c1, c2, cc, aero_data, aero_state%awa, env_state, k)
748  p = k * accept_factor
749  call warn_assert_msg(532446093, p <= 1d0, &
750  "kernel upper bound estimation is too tight, " &
751  // "could be caused by changing env_state")
752 
753  did_coag = .false.
754  if (pmc_random() .lt. p) then
755  call coagulate(aero_data, aero_state, p1, p2, c1, c2, cc)
756  did_coag = .true.
757  end if
758 
759  end subroutine maybe_coag_pair
760 
761 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
762 
763  !> Given <tt>(b1, c1)</tt> and <tt>(b2, c2)</tt>, find a random pair
764  !> of particles <tt>(b1, c1, i1)</tt> and <tt>(b2, c2, i2)</tt> that
765  !> are not the same particle particle as each other.
766  subroutine find_rand_pair(aero_sorted, b1, b2, c1, c2, i1, i2)
767 
768  !> Aerosol sorted data.
769  type(aero_sorted_t), intent(in) :: aero_sorted
770  !> Bin number of first particle.
771  integer, intent(in) :: b1
772  !> Bin number of second particle.
773  integer, intent(in) :: b2
774  !> Weight class of first particle.
775  integer, intent(in) :: c1
776  !> Weight class of second particle.
777  integer, intent(in) :: c2
778  !> First rand particle.
779  integer, intent(out) :: i1
780  !> Second rand particle.
781  integer, intent(out) :: i2
782 
783  call assert(619608562, &
784  integer_varray_n_entry(aero_sorted%size_class%inverse(b1, c1)) >= 1)
785  i1 = pmc_rand_int( &
786  integer_varray_n_entry(aero_sorted%size_class%inverse(b1, c1)))
787 
788  if ((b1 == b2) .and. (c1 == c2)) then
789  call assert(956184336, integer_varray_n_entry( &
790  aero_sorted%size_class%inverse(b2, c2)) >= 2)
791  i2 = pmc_rand_int( &
792  integer_varray_n_entry(aero_sorted%size_class%inverse(b2, c2)) - 1)
793  if (i2 == i1) then
794  i2 = integer_varray_n_entry(aero_sorted%size_class%inverse(b2, c2))
795  end if
796  else
797  call assert(271635751, integer_varray_n_entry( &
798  aero_sorted%size_class%inverse(b2, c2)) >= 1)
799  i2 = pmc_rand_int( &
800  integer_varray_n_entry(aero_sorted%size_class%inverse(b2, c2)))
801  end if
802 
803  end subroutine find_rand_pair
804 
805 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
806 
807  !> Actually coagulate pt1 and pt2 to form ptc and compute weighting
808  !> effects, including which particles should be lost and which
809  !> gained.
810  subroutine coagulate_weighting(pt1, pt2, ptc, c1, c2, cc, aero_data, &
811  aero_weight_array, remove_1, remove_2, create_new, id_1_lost, &
812  id_2_lost, aero_info_1, aero_info_2)
813 
814  !> First coagulating aerosol particle.
815  type(aero_particle_t), intent(in) :: pt1
816  !> Second coagulating aerosol particle.
817  type(aero_particle_t), intent(in) :: pt2
818  !> Combined aerosol particle resulting from coagulation of pt1
819  !> and pt2.
820  type(aero_particle_t), intent(inout) :: ptc
821  !> First weight class.
822  integer, intent(in) :: c1
823  !> Second weight class.
824  integer, intent(in) :: c2
825  !> Coagulated weight class.
826  integer, intent(in) :: cc
827  !> Aerosol data.
828  type(aero_data_t), intent(in) :: aero_data
829  !> Aerosol weight array.
830  type(aero_weight_array_t), intent(in) :: aero_weight_array
831  !> Whether to remove pt1.
832  logical, intent(out) :: remove_1
833  !> Whether to remove pt2.
834  logical, intent(out) :: remove_2
835  !> Whether to create ptc.
836  logical, intent(out) :: create_new
837  !> Whether the ID of pt1 will be lost due to coagulation.
838  logical, intent(out) :: id_1_lost
839  !> Whether the ID of pt2 will be lost due to coagulation.
840  logical, intent(out) :: id_2_lost
841  !> The removal information associated with pt1.
842  type(aero_info_t), intent(inout) :: aero_info_1
843  !> The removal information associated with pt2.
844  type(aero_info_t), intent(inout) :: aero_info_2
845 
846  real(kind=dp) :: r1, r2, rc, nc_min, nc1, nc2, ncc
847  real(kind=dp) :: prob_remove_1, prob_remove_2, prob_create_new
848  integer(kind=8) :: info_other_id
849  integer :: new_group
850 
851  call assert(371947172, pt1%id /= pt2%id)
852 
853  ! decide which old particles are to be removed and whether to
854  ! create the resulting coagulated particle
855  r1 = aero_particle_radius(pt1, aero_data)
856  r2 = aero_particle_radius(pt2, aero_data)
857  rc = aero_data_vol2rad(aero_data, aero_data_rad2vol(aero_data, r1) &
858  + aero_data_rad2vol(aero_data, r2))
859  nc1 = aero_weight_array_num_conc_at_radius(aero_weight_array, c1, r1)
860  nc2 = aero_weight_array_num_conc_at_radius(aero_weight_array, c2, r2)
861  ncc = aero_weight_array_num_conc_at_radius(aero_weight_array, cc, rc)
862  new_group = aero_weight_array_rand_group(aero_weight_array, cc, rc)
863  nc_min = min(nc1, nc2, ncc)
864  prob_remove_1 = nc_min / nc1
865  prob_remove_2 = nc_min / nc2
866  prob_create_new = nc_min / ncc
867  remove_1 = (pmc_random() < prob_remove_1)
868  ! FIXME
869  !if (aero_weight%type == AERO_WEIGHT_TYPE_MFA) then
870  ! remove_2 = .not. remove_1
871  !else
872  remove_2 = (pmc_random() < prob_remove_2)
873  !end if
874  create_new = (pmc_random() < prob_create_new)
875 
876  ! figure out what to do about the ID numbers of the various
877  ! particles --- we try to preserve particle IDs as much as
878  ! possible
879  if (create_new) then
880  id_1_lost = .false.
881  id_2_lost = .false.
882  if (remove_1 .and. remove_2) then
883  if (aero_particle_volume(pt1) > aero_particle_volume(pt2)) then
884  id_2_lost = .true.
885  else
886  id_1_lost = .true.
887  end if
888  end if
889  else
890  id_1_lost = remove_1
891  id_2_lost = remove_2
892  end if
893 
894  ! create a new particle and set its ID
895  if (create_new) then
896  call aero_particle_coagulate(pt1, pt2, ptc)
897  call aero_particle_set_weight(ptc, new_group, cc)
898  if (remove_1 .and. (.not. id_1_lost)) then
899  ptc%id = pt1%id
900  call assert(975059559, id_2_lost .eqv. remove_2)
901  elseif (remove_2 .and. (.not. id_2_lost)) then
902  ptc%id = pt2%id
903  call assert(246529753, id_1_lost .eqv. remove_1)
904  else
905  call aero_particle_new_id(ptc)
906  call assert(852038606, id_1_lost .eqv. remove_1)
907  call assert(254018921, id_2_lost .eqv. remove_2)
908  end if
909  info_other_id = ptc%id
910  else
911  info_other_id = 0
912  end if
913 
914  aero_info_1%id = pt1%id
915  aero_info_1%action = aero_info_coag
916  aero_info_1%other_id = info_other_id
917 
918  aero_info_2%id = pt2%id
919  aero_info_2%action = aero_info_coag
920  aero_info_2%other_id = info_other_id
921 
922  end subroutine coagulate_weighting
923 
924 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
925 
926  !> Join together particles \c p1 and \c p2, updating all particle
927  !> and bin structures to reflect the change.
928  subroutine coagulate(aero_data, aero_state, p1, p2, c1, c2, cc)
929 
930  !> Aerosol data.
931  type(aero_data_t), intent(in) :: aero_data
932  !> Aerosol state.
933  type(aero_state_t), intent(inout) :: aero_state
934  !> First particle index.
935  integer, intent(in) :: p1
936  !> Second particle index.
937  integer, intent(in) :: p2
938  !> First weight class.
939  integer, intent(in) :: c1
940  !> Second weight class.
941  integer, intent(in) :: c2
942  !> Coagulated weight class.
943  integer, intent(in) :: cc
944 
945  type(aero_particle_t) :: ptc
946  integer :: bn
947  type(aero_info_t) :: aero_info_1, aero_info_2
948  logical :: remove_1, remove_2, create_new, id_1_lost, id_2_lost
949 
950  call coagulate_weighting(aero_state%apa%particle(p1), &
951  aero_state%apa%particle(p2), ptc, c1, c2, cc, aero_data, &
952  aero_state%awa, remove_1, remove_2, create_new, id_1_lost, &
953  id_2_lost, aero_info_1, aero_info_2)
954 
955  ! remove old particles
956  if (p2 > p1) then
957  ! handle a tricky corner case where we have to watch for p2 or
958  ! p1 being the last entry in the array and being repacked when
959  ! the other one is removed
960  if (remove_2) then
961  call aero_state_remove_particle(aero_state, p2, &
962  id_2_lost, aero_info_2)
963  end if
964  if (remove_1) then
965  call aero_state_remove_particle(aero_state, p1, &
966  id_1_lost, aero_info_1)
967  end if
968  else
969  if (remove_1) then
970  call aero_state_remove_particle(aero_state, p1, &
971  id_1_lost, aero_info_1)
972  end if
973  if (remove_2) then
974  call aero_state_remove_particle(aero_state, p2, &
975  id_2_lost, aero_info_2)
976  end if
977  end if
978 
979  ! add new particle
980  if (create_new) then
981  call aero_state_add_particle(aero_state, ptc, aero_data, &
982  allow_resort=.false.)
983  end if
984 
985  end subroutine coagulate
986 
987 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
988 
989 end module pmc_coagulation
pmc_coagulation::sample_source_particle
subroutine sample_source_particle(aero_state, aero_data, env_state, coag_kernel_type, bs, cs, coag_particle, n_samp_mean, accept_factor, n_samp, n_coag, n_remove, source_particle)
Sample coagulation partners for a single coagulation event.
Definition: coagulation.F90:336
pmc_rand::prob_round
integer function prob_round(val)
Round val to floor(val) or ceiling(val) with probability proportional to the relative distance from v...
Definition: rand.F90:215
pmc_aero_state::aero_state_remove_particle
subroutine aero_state_remove_particle(aero_state, i_part, record_removal, aero_info)
Remove the given particle and possibly record the removal.
Definition: aero_state.F90:429
pmc_aero_data::aero_data_n_spec
elemental integer function aero_data_n_spec(aero_data)
Return the number of aerosol species, or -1 if uninitialized.
Definition: aero_data.F90:246
pmc_coagulation::coagulate_weighting
subroutine coagulate_weighting(pt1, pt2, ptc, c1, c2, cc, aero_data, aero_weight_array, remove_1, remove_2, create_new, id_1_lost, id_2_lost, aero_info_1, aero_info_2)
Actually coagulate pt1 and pt2 to form ptc and compute weighting effects, including which particles s...
Definition: coagulation.F90:813
pmc_mpi
Wrapper functions for MPI.
Definition: mpi.F90:13
pmc_aero_weight
The aero_weight_t structure and associated subroutines.
Definition: aero_weight.F90:9
pmc_coagulation::determine_target_and_source
subroutine determine_target_and_source(aero_weight_array, bin_grid, b1, b2, c1, c2, cc, bt, bs, ct, cs, correct_weight_ordering)
Determine the source and target particle bin/group for per-particle coagulation, if possible.
Definition: coagulation.F90:249
pmc_aero_sorted::aero_sorted_particle_in_bin
integer function aero_sorted_particle_in_bin(aero_sorted, aero_particle, aero_data)
Find the bin number that contains a given particle.
Definition: aero_sorted.F90:380
pmc_rand::rand_poisson
integer function rand_poisson(mean)
Generate a Poisson-distributed random number with the given mean.
Definition: rand.F90:253
pmc_aero_weight_array::aero_weight_array_check_monotonicity
subroutine aero_weight_array_check_monotonicity(aero_weight_array, monotone_increasing, monotone_decreasing)
Determine whether all weight functions in an array are monotone increasing, monotone decreasing,...
Definition: aero_weight_array.F90:306
pmc_coagulation::maybe_coag_pair
subroutine maybe_coag_pair(env_state, aero_data, aero_state, b1, b2, c1, c2, cc, coag_kernel_type, accept_factor, did_coag)
Choose a random pair for potential coagulation and test its probability of coagulation....
Definition: coagulation.F90:714
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_aero_sorted::aero_sorted_n_bin
integer function aero_sorted_n_bin(aero_sorted)
Returns the number of size bins.
Definition: aero_sorted.F90:81
pmc_constants::dp
integer, parameter dp
Kind of a double precision real number.
Definition: constants.F90:12
pmc_coag_kernel
Generic coagulation kernel.
Definition: coag_kernel.F90:9
pmc_aero_state::aero_state_add_particle
subroutine aero_state_add_particle(aero_state, aero_particle, aero_data, allow_resort)
Add the given particle.
Definition: aero_state.F90:367
pmc_aero_weight_array
The aero_weight_array_t structure and associated subroutines.
Definition: aero_weight_array.F90:10
pmc_coagulation::coag_accel_max_cv
real(kind=dp), parameter coag_accel_max_cv
Maximum allowed coefficient-of-variation due to undersampling in accelerated coagulation.
Definition: coagulation.F90:33
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_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_aero_state
The aero_state_t structure and assocated subroutines.
Definition: aero_state.F90:9
pmc_aero_weight_array::aero_weight_array_t
An array of aerosol size distribution weighting functions.
Definition: aero_weight_array.F90:35
pmc_coagulation::compute_n_source
subroutine compute_n_source(n_part, k_max, del_t, n_source_per_target, accept_factor)
Definition: coagulation.F90:313
pmc_rand::pmc_rand_int
integer function pmc_rand_int(n)
Returns a random integer between 1 and n.
Definition: rand.F90:173
pmc_coagulation::find_rand_pair
subroutine find_rand_pair(aero_sorted, b1, b2, c1, c2, i1, i2)
Given (b1, c1) and (b2, c2), find a random pair of particles (b1, c1, i1) and (b2,...
Definition: coagulation.F90:767
pmc_aero_weight_array::aero_weight_array_num_conc
real(kind=dp) function aero_weight_array_num_conc(aero_weight_array, aero_particle, aero_data)
Compute the number concentration for a particle (m^{-3}).
Definition: aero_weight_array.F90:259
pmc_aero_weight_array::aero_weight_array_num_conc_at_radius
real(kind=dp) function aero_weight_array_num_conc_at_radius(aero_weight_array, i_class, radius)
Compute the total number concentration at a given radius (m^3).
Definition: aero_weight_array.F90:234
pmc_coagulation::coag_target_with_source
subroutine coag_target_with_source(aero_state, aero_data, bt, ct, target_unif_entry, source_particle, cc)
Coagulate a sampled source particle with a target particle.
Definition: coagulation.F90:524
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_coagulation
Aerosol particle coagulation.
Definition: coagulation.F90:11
pmc_util::assert_msg
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
Definition: util.F90:78
pmc_coagulation::coag_accel_n_event
real(kind=dp), parameter coag_accel_n_event
Minimum number of coagulation events per large particle for which accelerated coagulation is used.
Definition: coagulation.F90:30
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_coagulation::per_set_coag
subroutine per_set_coag(coag_kernel_type, k_max, env_state, aero_data, aero_state, del_t, tot_n_samp, tot_n_coag, b1, b2, c1, c2, cc)
Do set-wise coagulation.
Definition: coagulation.F90:589
pmc_coagulation::compute_n_samp
subroutine compute_n_samp(ni, nj, same_bin, k_max, del_t, n_samp_mean, n_samp, accept_factor)
Compute the number of samples required for the pair of bins.
Definition: coagulation.F90:650
pmc_aero_state::aero_state_dup_particle
subroutine aero_state_dup_particle(aero_state, aero_data, i_part, n_part_mean, random_weight_group)
Add copies or remove a particle, with a given mean number of resulting particles.
Definition: aero_state.F90:490
pmc_coag_kernel::est_k_minmax_binned_unweighted
subroutine est_k_minmax_binned_unweighted(bin_grid, coag_kernel_type, aero_data, env_state, k_min, k_max)
Estimate an array of minimum and maximum kernel values. Given particles v1 in bin b1 and v2 in bin b2...
Definition: coag_kernel.F90:255
pmc_aero_weight_array::aero_weight_array_rand_group
integer function aero_weight_array_rand_group(aero_weight_array, i_class, radius)
Choose a random group at the given radius, with probability inversely proportional to group weight at...
Definition: aero_weight_array.F90:376
pmc_env_state
The env_state_t structure and associated subroutines.
Definition: env_state.F90:9
pmc_aero_sorted
The aero_sorted_t structure and assocated subroutines.
Definition: aero_sorted.F90:9
pmc_aero_data::aero_data_vol2rad
real(kind=dp) elemental function aero_data_vol2rad(aero_data, v)
Convert mass-equivalent volume (m^3) to geometric radius (m).
Definition: aero_data.F90:102
pmc_aero_weight_array::aero_weight_array_minmax_num_conc
subroutine aero_weight_array_minmax_num_conc(aero_weight_array, i_class, radius_1, radius_2, num_conc_min, num_conc_max)
Compute the maximum and minimum number concentrations between the given radii.
Definition: aero_weight_array.F90:340
pmc_coagulation::mc_coag
subroutine mc_coag(coag_kernel_type, env_state, aero_data, aero_state, del_t, tot_n_samp, tot_n_coag)
Do coagulation for time del_t.
Definition: coagulation.F90:45
pmc_aero_data::aero_data_t
Aerosol material properties and associated data.
Definition: aero_data.F90:55
pmc_rand::sample_without_replacement
subroutine sample_without_replacement(sample_size, population_size, sample)
Generates a vector of indices for a sample set from a population without replacement.
Definition: rand.F90:516
pmc_aero_sorted::aero_sorted_t
Sorting of particles into bins.
Definition: aero_sorted.F90:47
pmc_stats::student_t_95_coeff
real(kind=dp) function student_t_95_coeff(n_sample)
Return a fairly tight upper-bound on the Student's t coefficient for the 95% confidence interval.
Definition: stats.F90:511
pmc_util
Common utility subroutines.
Definition: util.F90:9
pmc_coagulation::mc_coag_for_bin
subroutine mc_coag_for_bin(coag_kernel_type, env_state, aero_data, aero_state, del_t, tot_n_samp, tot_n_coag, b1, b2, c1, c2)
Do coagulation for time del_t for the given bins.
Definition: coagulation.F90:103
pmc_coagulation::coagulate
subroutine coagulate(aero_data, aero_state, p1, p2, c1, c2, cc)
Join together particles p1 and p2, updating all particle and bin structures to reflect the change.
Definition: coagulation.F90:929
pmc_coag_kernel::coag_dest_class
integer function coag_dest_class(aero_weight_array, aero_data, bin_grid, i_bin, j_bin, i_class, j_class)
Determine the weight class in which coagulated particles will be placed.
Definition: coag_kernel.F90:384
pmc_bin_grid
The bin_grid_t structure and associated subroutines.
Definition: bin_grid.F90:9
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_coagulation::max_allowable_growth_factor
real(kind=dp), parameter max_allowable_growth_factor
Maximum allowable growth factor of a target particle volume within one timestep when using accelerate...
Definition: coagulation.F90:36
pmc_aero_sorted::aero_sorted_n_class
integer function aero_sorted_n_class(aero_sorted)
Returns the number of weight classes.
Definition: aero_sorted.F90:105
pmc_coag_kernel::max_coag_num_conc_factor
subroutine max_coag_num_conc_factor(aero_weight_array, aero_data, bin_grid, i_bin, j_bin, i_class, j_class, ij_class, f_max)
Determine the minimum and maximum number concentration factors for coagulation.
Definition: coag_kernel.F90:424
pmc_bin_grid::bin_grid_t
1D grid, either logarithmic or linear.
Definition: bin_grid.F90:33
pmc_stats
The stats_t type and associated subroutines.
Definition: stats.F90:9
pmc_coag_kernel::num_conc_weighted_kernel
subroutine num_conc_weighted_kernel(coag_kernel_type, aero_particle_1, aero_particle_2, i_class, j_class, ij_class, aero_data, aero_weight_array, env_state, k)
Compute the kernel value with the given number concentration weighting.
Definition: coag_kernel.F90:178
pmc_rand::rand_binomial
integer function rand_binomial(n, p)
Generate a Binomial-distributed random number with the given parameters.
Definition: rand.F90:320
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_sorted::aero_sorted_move_particle
subroutine aero_sorted_move_particle(aero_sorted, i_part, new_bin, new_group, new_class)
Move a particle to a different bin and group.
Definition: aero_sorted.F90:477
pmc_coagulation::try_per_particle_coag
subroutine try_per_particle_coag(coag_kernel_type, k_max, env_state, aero_data, aero_state, del_t, tot_n_samp, tot_n_coag, b1, b2, c1, c2, cc, per_particle_coag_succeeded)
Attempt per-particle coagulation.
Definition: coagulation.F90:153