PartMC  2.3.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)
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 (aero_state%aero_sorted%size_class%inverse(b1, &
78  c1)%n_entry == 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 (aero_state%aero_sorted%size_class%inverse(b2, &
87  c2)%n_entry == 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_state%aero_sorted%bin_grid, b1, &
132  b2, c1, c2)
133  call max_coag_num_conc_factor(aero_state%awa, &
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(aero_state%aero_sorted%size_class%inverse(bs, &
200  cs)%n_entry, k_max, del_t, n_source_per_target, accept_factor)
201  if (n_source_per_target < coag_accel_n_event) then
202  per_particle_coag_succeeded = .false.
203  return
204  end if
205 
206  min_target_vol = rad2vol(aero_state%aero_sorted%bin_grid%edges(bt))
207  max_source_vol = rad2vol(aero_state%aero_sorted%bin_grid%edges(bs+1))
208  max_new_target_vol = min_target_vol + max_source_vol * n_source_per_target
209  ! check for unacceptably large volume change
210  max_target_growth_factor = max_new_target_vol / min_target_vol
211  if (max_target_growth_factor > max_allowable_growth_factor) then
212  per_particle_coag_succeeded = .false.
213  return
214  end if
215 
216  call aero_particle_allocate(target_particle)
217  call aero_particle_allocate(source_particle)
218 
219  ! work backwards to avoid particle movement issues
220  do target_unif_entry &
221  = aero_state%aero_sorted%size_class%inverse(bt, ct)%n_entry,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  call aero_particle_copy(aero_state%apa%particle(target_part), &
227  target_particle)
228  call sample_source_particle(aero_state, aero_data, env_state, &
229  coag_kernel_type, bs, cs, target_particle, n_source_per_target, &
230  accept_factor, n_samp, n_coag, n_remove, source_particle)
231  if (n_coag > 0) then
232  call coag_target_with_source(aero_state, bt, ct, target_unif_entry, &
233  source_particle, cc)
234  end if
235  tot_n_samp = tot_n_samp + n_samp
236  tot_n_coag = tot_n_coag + n_coag
237  ! we discard n_remove information at present
238  end do
239 
240  call aero_particle_deallocate(target_particle)
241  call aero_particle_deallocate(source_particle)
242 
243  per_particle_coag_succeeded = .true.
244 
245  end subroutine try_per_particle_coag
246 
247 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
248 
249  !> Determine the source and target particle bin/group for
250  !> per-particle coagulation, if possible.
251  subroutine determine_target_and_source(aero_weight_array, bin_grid, b1, &
252  b2, c1, c2, cc, bt, bs, ct, cs, correct_weight_ordering)
253 
254  !> Aero weight array.
255  type(aero_weight_array_t), intent(in) :: aero_weight_array
256  !> Bin grid.
257  type(bin_grid_t), intent(in) :: bin_grid
258  !> First bin number.
259  integer, intent(in) :: b1
260  !> Second bin number.
261  integer, intent(in) :: b2
262  !> First weight class.
263  integer, intent(in) :: c1
264  !> Second weight class.
265  integer, intent(in) :: c2
266  !> Coagulated weight class.
267  integer, intent(in) :: cc
268  !> Target bin number.
269  integer, intent(out) :: bt
270  !> Source bin number.
271  integer, intent(out) :: bs
272  !> Target weight class.
273  integer, intent(out) :: ct
274  !> Source weight class.
275  integer, intent(out) :: cs
276  !> Whether the weight ordering is correct for per-particle coagulation.
277  logical, intent(out) :: correct_weight_ordering
278 
279  real(kind=dp) :: nc1_min, nc1_max, nc2_min, nc2_max
280  logical :: monotone_increasing, monotone_decreasing
281 
282  call aero_weight_array_check_monotonicity(aero_weight_array, &
283  monotone_increasing, monotone_decreasing)
284  if (.not. monotone_decreasing) then
285  correct_weight_ordering = .false.
286  return
287  end if
288 
289  call aero_weight_array_minmax_num_conc(aero_weight_array, c1, &
290  bin_grid%edges(b1), bin_grid%edges(b1 + 1), nc1_min, nc1_max)
291  call aero_weight_array_minmax_num_conc(aero_weight_array, c2, &
292  bin_grid%edges(b2), bin_grid%edges(b2 + 1), nc2_min, nc2_max)
293 
294  ! we have already confirmed monotone_decreasing weights above
295  correct_weight_ordering = .false.
296  if ((nc1_max < nc2_min) .and. (cc == c1)) then
297  bt = b1
298  bs = b2
299  ct = c1
300  cs = c2
301  correct_weight_ordering = .true.
302  end if
303  if ((nc2_max < nc1_min) .and. (cc == c2)) then
304  bt = b2
305  bs = b1
306  ct = c2
307  cs = c1
308  correct_weight_ordering = .true.
309  end if
310 
311  end subroutine determine_target_and_source
312 
313 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
314 
315  subroutine compute_n_source(n_part, k_max, del_t, n_source_per_target, &
316  accept_factor)
317 
318  !> Number of particles available as partners.
319  integer, intent(in) :: n_part
320  !> Maximum coagulation kernel (s^{-1} m^3).
321  real(kind=dp), intent(in) :: k_max
322  !> Timestep (s).
323  real(kind=dp), intent(in) :: del_t
324  !> Mean number of source particles per target particle.
325  real(kind=dp), intent(out) :: n_source_per_target
326  !> Accept factor for samples.
327  real(kind=dp), intent(out) :: accept_factor
328 
329  n_source_per_target = k_max * del_t * real(n_part, kind=dp)
330  accept_factor = 1d0 / k_max
331 
332  end subroutine compute_n_source
333 
334 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
335 
336  !> Sample coagulation partners for a single coagulation event.
337  subroutine sample_source_particle(aero_state, aero_data, env_state, &
338  coag_kernel_type, bs, cs, coag_particle, n_samp_mean, &
339  accept_factor, n_samp, n_coag, n_remove, source_particle)
340 
341  !> Aerosol state.
342  type(aero_state_t), intent(inout) :: aero_state
343  !> Aerosol data.
344  type(aero_data_t), intent(in) :: aero_data
345  !> Environment state.
346  type(env_state_t), intent(in) :: env_state
347  !> Coagulation kernel type.
348  integer, intent(in) :: coag_kernel_type
349  !> Bin to sample particles from.
350  integer, intent(in) :: bs
351  !> Weight class to sample particles from.
352  integer, intent(in) :: cs
353  !> Aerosol particle that coagulation will be with.
354  type(aero_particle_t), intent(in) :: coag_particle
355  !> Mean number of samples to use.
356  real(kind=dp), intent(in) :: n_samp_mean
357  !> Probability factor of accepting samples.
358  real(kind=dp), intent(in) :: accept_factor
359  !> Number of samples used.
360  integer, intent(out) :: n_samp
361  !> Number of coagulations performed.
362  integer, intent(out) :: n_coag
363  !> Number of particles removed.
364  integer, intent(out) :: n_remove
365  !> Sampled average coagulation partner particle.
366  type(aero_particle_t), intent(inout) :: source_particle
367 
368  real(kind=dp) :: prob_remove_i, prob_remove_source_max
369  real(kind=dp) :: prob_coag, prob_coag_tot, prob_coag_mean
370  real(kind=dp) :: num_conc_i, num_conc_source_min, num_conc_target, k
371  real(kind=dp) :: vol_sq(aero_data%n_spec), vol_mean(aero_data%n_spec)
372  real(kind=dp) :: vol_cv(aero_data%n_spec), vol_cv_max, mean_95_conf_cv
373  integer :: n_samp_remove, n_samp_extra, n_samp_total, n_avg, i_samp
374  integer :: i_unif_entry, i_part, target_id, new_bin, ct
375  type(aero_particle_t), pointer :: i_particle
376  type(aero_info_t) :: aero_info
377 
378  if (aero_state%aero_sorted%size_class%inverse(bs, cs)%n_entry == 0) then
379  n_samp = 0
380  n_remove = 0
381  n_coag = 0
382  return
383  end if
384 
385  num_conc_target = aero_weight_array_num_conc(aero_state%awa, coag_particle)
386  target_id = coag_particle%id
387  ct = coag_particle%weight_class
388 
389  num_conc_source_min = aero_weight_array_num_conc_at_radius( &
390  aero_state%awa, cs, aero_state%aero_sorted%bin_grid%edges(bs))
391  prob_remove_source_max = num_conc_target / num_conc_source_min
392  call assert(653606684, prob_remove_source_max <= 1d0)
393 
394  n_samp_remove = rand_poisson(prob_remove_source_max * n_samp_mean)
395  n_samp_extra = rand_poisson((1d0 - prob_remove_source_max) * n_samp_mean)
396  n_samp_total = n_samp_remove + n_samp_extra
397 
398  n_avg = 0
399  n_samp = 0
400  n_remove = 0
401  prob_coag_tot = 0d0
402  call aero_particle_deallocate(source_particle)
403  call aero_particle_allocate_size(source_particle, &
404  aero_data%n_spec, aero_data%n_source)
405  vol_sq = 0d0
406 
407  ! FIXME: Can't we just do n_samp = 1,n_samp_total and shift tests
408  ! to the end?
409  do i_samp = 1,n_samp_total
410  if (aero_state%aero_sorted%size_class%inverse(bs, cs)%n_entry == 0) exit
411  if ((n_samp > n_samp_remove) .and. (n_avg >= 2)) then
412  vol_mean = source_particle%vol / real(n_avg, kind=dp)
413  where(vol_mean > 0d0) &
414  vol_cv = sqrt(vol_sq / real(n_avg, kind=dp) - (vol_mean)**2) &
415  / vol_mean
416  vol_cv_max = maxval(vol_cv)
417  mean_95_conf_cv = vol_cv_max / sqrt(real(n_avg, kind=dp)) &
418  * student_t_95_coeff(n_avg)
419  ! FIXME: We are using just the max of the diagonal of the
420  ! covariance matrix. Is this well-justified?
421  if (mean_95_conf_cv < coag_accel_max_cv) exit
422  end if
423  n_samp = n_samp + 1
424  ! FIXME: We are sampling with replacement. Is this a problem?
425  i_unif_entry &
426  = pmc_rand_int(aero_state%aero_sorted%size_class%inverse(bs, &
427  cs)%n_entry)
428  i_part = aero_state%aero_sorted%size_class%inverse(bs, &
429  cs)%entry(i_unif_entry)
430  i_particle => aero_state%apa%particle(i_part)
431  ! re-get j_part as particle ordering may be changing
432  call num_conc_weighted_kernel(coag_kernel_type, i_particle, &
433  coag_particle, cs, ct, ct, aero_data, aero_state%awa, env_state, k)
434  prob_coag = k * accept_factor
435  prob_coag_tot = prob_coag_tot + prob_coag
436  if (pmc_random() < prob_coag) then
437  n_avg = n_avg + 1
438  call aero_particle_coagulate(source_particle, i_particle, &
439  source_particle)
440  vol_sq = vol_sq + i_particle%vol**2
441  if (i_samp <= n_samp_remove) then
442  num_conc_i = aero_weight_array_num_conc(aero_state%awa, &
443  i_particle)
444  prob_remove_i = num_conc_target / num_conc_i
445  if (pmc_random() < prob_remove_i / prob_remove_source_max) then
446  n_remove = n_remove + 1
447  call aero_info_allocate(aero_info)
448  aero_info%id = i_particle%id
449  aero_info%action = aero_info_coag
450  aero_info%other_id = target_id
451  call aero_state_remove_particle_with_info(aero_state, &
452  i_part, aero_info)
453  call aero_info_deallocate(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 
471  end subroutine sample_source_particle
472 
473 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
474 
475  !> Coagulate a sampled source particle with a target particle.
476  subroutine coag_target_with_source(aero_state, bt, ct, target_unif_entry, &
477  source_particle, cc)
478 
479  !> Aerosol state.
480  type(aero_state_t), intent(inout) :: aero_state
481  !> Bin of coagulating particle.
482  integer, intent(in) :: bt
483  !> Weight class of coagulating particle.
484  integer, intent(in) :: ct
485  !> Entry-in-bin of coagulating particle.
486  integer, intent(in) :: target_unif_entry
487  !> Sampled particle to coagulate with.
488  type(aero_particle_t), intent(in) :: source_particle
489  !> Weight class for coagulated particle.
490  integer, intent(in) :: cc
491 
492  integer :: target_part, target_id, new_bin, new_group
493  real(kind=dp) :: old_num_conc_target, new_num_conc_target
494 
495  target_part = aero_state%aero_sorted%size_class%inverse(bt, &
496  ct)%entry(target_unif_entry)
497  target_id = aero_state%apa%particle(target_part)%id
498  old_num_conc_target = aero_weight_array_num_conc(aero_state%awa, &
499  aero_state%apa%particle(target_part))
500  call aero_particle_coagulate(aero_state%apa%particle(target_part), &
501  source_particle, aero_state%apa%particle(target_part))
502  aero_state%apa%particle(target_part)%id = target_id
503  ! assign to a randomly chosen group
504  new_group = aero_weight_array_rand_group(aero_state%awa, cc, &
505  aero_particle_radius(aero_state%apa%particle(target_part)))
506  call aero_particle_set_weight(aero_state%apa%particle(target_part), &
507  new_group, cc)
508  ! fix bin due to composition changes
509  new_bin = aero_sorted_particle_in_bin(aero_state%aero_sorted, &
510  aero_state%apa%particle(target_part))
511  if ((new_bin < 1) &
512  .or. (new_bin > aero_state%aero_sorted%bin_grid%n_bin)) then
513  call die_msg(765620746, "particle outside of bin_grid: " &
514  // "try reducing the timestep del_t")
515  end if
516  call aero_sorted_move_particle(aero_state%aero_sorted, target_part, &
517  new_bin, new_group, cc)
518  ! Adjust particle number to account for weight changes
519  ! bt/target_group/target_entry are invalid, but
520  ! target_part is still good. We are treating all particles in all
521  ! groups together, and randomly reassigning between groups above,
522  ! so here we can't use aero_state_reweight_particle(), as that
523  ! assumes we are staying in the same weight group.
524  new_num_conc_target = aero_weight_array_num_conc(aero_state%awa, &
525  aero_state%apa%particle(target_part))
526  call aero_state_dup_particle(aero_state, target_part, &
527  old_num_conc_target / new_num_conc_target, random_weight_group=.true.)
528  ! we should only be doing this for decreasing weights
529  call assert(654300924, aero_state%apa%particle(target_part)%id &
530  == target_id)
531 
532  end subroutine coag_target_with_source
533 
534 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
535 
536  !> Do set-wise coagulation.
537  subroutine per_set_coag(coag_kernel_type, k_max, env_state, aero_data, &
538  aero_state, del_t, tot_n_samp, tot_n_coag, b1, b2, c1, c2, cc)
539 
540  !> Coagulation kernel type.
541  integer, intent(in) :: coag_kernel_type
542  !> Maximum coagulation kernel (s^{-1} m^3).
543  real(kind=dp), intent(in) :: k_max
544  !> Environment state.
545  type(env_state_t), intent(in) :: env_state
546  !> Aerosol data.
547  type(aero_data_t), intent(in) :: aero_data
548  !> Aerosol state.
549  type(aero_state_t), intent(inout) :: aero_state
550  !> Timestep for coagulation.
551  real(kind=dp) :: del_t
552  !> Total number of samples tested.
553  integer, intent(inout) :: tot_n_samp
554  !> Number of coagulation events.
555  integer, intent(inout) :: tot_n_coag
556  !> First bin number.
557  integer, intent(in) :: b1
558  !> Second bin number.
559  integer, intent(in) :: b2
560  !> First weight class.
561  integer, intent(in) :: c1
562  !> Second weight class.
563  integer, intent(in) :: c2
564  !> Coagulated weight class.
565  integer, intent(in) :: cc
566 
567  real(kind=dp) :: n_samp_mean, accept_factor
568  integer :: i_samp, n_samp, n1, n2
569  logical :: did_coag
570 
571  n1 = aero_state%aero_sorted%size_class%inverse(b1, c1)%n_entry
572  n2 = aero_state%aero_sorted%size_class%inverse(b2, c2)%n_entry
573  call compute_n_samp(n1, n2, ((b1 == b2) .and. (c1 == c2)), k_max, del_t, &
574  n_samp_mean, n_samp, accept_factor)
575  tot_n_samp = tot_n_samp + n_samp
576 
577  do i_samp = 1,n_samp
578  ! check we still have enough particles to coagulate
579  n1 = aero_state%aero_sorted%size_class%inverse(b1, c1)%n_entry
580  n2 = aero_state%aero_sorted%size_class%inverse(b2, c2)%n_entry
581  if (((n1 < 2) .and. (b1 == b2) .and. (c1 == c2)) &
582  .or. (n1 < 1) .or. (n2 < 1)) &
583  exit
584  call maybe_coag_pair(env_state, aero_data, aero_state, b1, b2, c1, c2, &
585  cc, coag_kernel_type, accept_factor, did_coag)
586  if (did_coag) tot_n_coag = tot_n_coag + 1
587  end do
588 
589  end subroutine per_set_coag
590 
591 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
592 
593  !> Compute the number of samples required for the pair of bins.
594  subroutine compute_n_samp(ni, nj, same_bin, k_max, del_t, n_samp_mean, &
595  n_samp, accept_factor)
596 
597  !> Number particles in first bin.
598  integer, intent(in) :: ni
599  !> Number particles in second bin.
600  integer, intent(in) :: nj
601  !> Whether first bin is second bin.
602  logical, intent(in) :: same_bin
603  !> Maximum kernel value (s^{-1} m^3).
604  real(kind=dp), intent(in) :: k_max
605  !> Timestep (s).
606  real(kind=dp), intent(in) :: del_t
607  !> Mean number of samples per timestep.
608  real(kind=dp), intent(out) :: n_samp_mean
609  !> Number of samples per timestep.
610  integer, intent(out) :: n_samp
611  !> Scale factor for accept probability (1).
612  real(kind=dp), intent(out) :: accept_factor
613 
614  real(kind=dp) :: r_samp
615  real(kind=dp) :: n_possible ! use real(kind=dp) to avoid integer overflow
616  ! could use integer*8 or integer(kind = 8)
617  ! or di = selected_int_kind(18), integer(kind=di)
618  ! to represent 10^{-18} to 10^{18}
619 
620  if (same_bin) then
621  ! don't change this to ni * (ni - 1) as the ni/nj distinction
622  ! is important for coagulation_dist, which also calls this
623  n_possible = real(ni, kind=dp) * (real(nj, kind=dp) - 1d0) / 2d0
624  else
625  n_possible = real(ni, kind=dp) * real(nj, kind=dp)
626  end if
627 
628  r_samp = k_max * del_t
629  n_samp_mean = r_samp * n_possible
630  n_samp = rand_poisson(n_samp_mean)
631  accept_factor = 1d0 / k_max
632 
633  ! possible variants:
634  ! A: accept_factor = 1d0 / k_max
635  ! B: accept_factor = del_t * n_possible / (n_samp * V)
636  ! timings of test suite as of 2010-12-22T17:12:14-0600:
637  ! A with n_samp = prob_round(n_samp_mean):
638  ! 159.82 162.18 156.28
639  ! A with n_samp = rand_poisson(n_samp_mean):
640  ! 151.93 158.38 174.74 157.65
641  ! B with n_samp = ceiling(n_samp_mean):
642  ! 196.06 200.23 192.41
643  ! B with n_samp = ceiling(n_samp_mean + 1 * sqrt(n_samp_mean)):
644  ! 189.78 211.12 195.19
645  ! B with n_samp = ceiling(n_samp_mean + 3 * sqrt(n_samp_mean)):
646  ! 214.60 201.25 203.55
647 
648  end subroutine compute_n_samp
649 
650 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
651 
652  !> Choose a random pair for potential coagulation and test its
653  !> probability of coagulation. If it happens, do the coagulation and
654  !> update all structures.
655  !!
656  !! The probability of a coagulation will be taken as <tt>(kernel /
657  !! k_max)</tt>.
658  subroutine maybe_coag_pair(env_state, aero_data, aero_state, b1, b2, &
659  c1, c2, cc, coag_kernel_type, accept_factor, did_coag)
660 
661  !> Environment state.
662  type(env_state_t), intent(in) :: env_state
663  !> Aerosol data.
664  type(aero_data_t), intent(in) :: aero_data
665  !> Aerosol state.
666  type(aero_state_t), intent(inout) :: aero_state
667  !> Bin of first particle.
668  integer, intent(in) :: b1
669  !> Bin of second particle.
670  integer, intent(in) :: b2
671  !> First weight class.
672  integer, intent(in) :: c1
673  !> Second weight class.
674  integer, intent(in) :: c2
675  !> Coagulated weight class.
676  integer, intent(in) :: cc
677  !> Coagulation kernel type.
678  integer, intent(in) :: coag_kernel_type
679  !> Scale factor for accept probability (1).
680  real(kind=dp), intent(in) :: accept_factor
681  !> Whether a coagulation occured.
682  logical, intent(out) :: did_coag
683 
684  integer :: i1, i2
685  integer :: p1, p2
686  real(kind=dp) :: p, k
687 
688  call find_rand_pair(aero_state%aero_sorted, b1, b2, c1, c2, i1, i2)
689  p1 = aero_state%aero_sorted%size_class%inverse(b1, c1)%entry(i1)
690  p2 = aero_state%aero_sorted%size_class%inverse(b2, c2)%entry(i2)
691  call num_conc_weighted_kernel(coag_kernel_type, &
692  aero_state%apa%particle(p1), aero_state%apa%particle(p2), &
693  c1, c2, cc, aero_data, aero_state%awa, env_state, k)
694  p = k * accept_factor
695  call warn_assert_msg(532446093, p <= 1d0, &
696  "kernel upper bound estimation is too tight, " &
697  // "could be caused by changing env_state")
698 
699  did_coag = .false.
700  if (pmc_random() .lt. p) then
701  call coagulate(aero_data, aero_state, p1, p2, c1, c2, cc)
702  did_coag = .true.
703  end if
704 
705  end subroutine maybe_coag_pair
706 
707 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
708 
709  !> Given <tt>(b1, c1)</tt> and <tt>(b2, c2)</tt>, find a random pair
710  !> of particles <tt>(b1, c1, i1)</tt> and <tt>(b2, c2, i2)</tt> that
711  !> are not the same particle particle as each other.
712  subroutine find_rand_pair(aero_sorted, b1, b2, c1, c2, i1, i2)
713 
714  !> Aerosol sorted data.
715  type(aero_sorted_t), intent(in) :: aero_sorted
716  !> Bin number of first particle.
717  integer, intent(in) :: b1
718  !> Bin number of second particle.
719  integer, intent(in) :: b2
720  !> Weight class of first particle.
721  integer, intent(in) :: c1
722  !> Weight class of second particle.
723  integer, intent(in) :: c2
724  !> First rand particle.
725  integer, intent(out) :: i1
726  !> Second rand particle.
727  integer, intent(out) :: i2
728 
729  call assert(619608562, aero_sorted%size_class%inverse(b1, c1)%n_entry >= 1)
730  i1 = pmc_rand_int(aero_sorted%size_class%inverse(b1, c1)%n_entry)
731 
732  if ((b1 == b2) .and. (c1 == c2)) then
733  call assert(956184336, &
734  aero_sorted%size_class%inverse(b2, c2)%n_entry >= 2)
735  i2 = pmc_rand_int(aero_sorted%size_class%inverse(b2, c2)%n_entry - 1)
736  if (i2 == i1) then
737  i2 = aero_sorted%size_class%inverse(b2, c2)%n_entry
738  end if
739  else
740  call assert(271635751, &
741  aero_sorted%size_class%inverse(b2, c2)%n_entry >= 1)
742  i2 = pmc_rand_int(aero_sorted%size_class%inverse(b2, c2)%n_entry)
743  end if
744 
745  end subroutine find_rand_pair
746 
747 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
748 
749  !> Actually coagulate pt1 and pt2 to form ptc and compute weighting
750  !> effects, including which particles should be lost and which
751  !> gained.
752  subroutine coagulate_weighting(pt1, pt2, ptc, c1, c2, cc, aero_data, &
753  aero_weight_array, remove_1, remove_2, create_new, id_1_lost, &
754  id_2_lost, aero_info_1, aero_info_2)
755 
756  !> First coagulating aerosol particle.
757  type(aero_particle_t), intent(in) :: pt1
758  !> Second coagulating aerosol particle.
759  type(aero_particle_t), intent(in) :: pt2
760  !> Combined aerosol particle resulting from coagulation of pt1
761  !> and pt2.
762  type(aero_particle_t), intent(inout) :: ptc
763  !> First weight class.
764  integer, intent(in) :: c1
765  !> Second weight class.
766  integer, intent(in) :: c2
767  !> Coagulated weight class.
768  integer, intent(in) :: cc
769  !> Aerosol data.
770  type(aero_data_t), intent(in) :: aero_data
771  !> Aerosol weight array.
772  type(aero_weight_array_t), intent(in) :: aero_weight_array
773  !> Whether to remove pt1.
774  logical, intent(out) :: remove_1
775  !> Whether to remove pt2.
776  logical, intent(out) :: remove_2
777  !> Whether to create ptc.
778  logical, intent(out) :: create_new
779  !> Whether the ID of pt1 will be lost due to coagulation.
780  logical, intent(out) :: id_1_lost
781  !> Whether the ID of pt2 will be lost due to coagulation.
782  logical, intent(out) :: id_2_lost
783  !> The removal information associated with pt1.
784  type(aero_info_t), intent(inout) :: aero_info_1
785  !> The removal information associated with pt2.
786  type(aero_info_t), intent(inout) :: aero_info_2
787 
788  real(kind=dp) :: r1, r2, rc, nc_min, nc1, nc2, ncc
789  real(kind=dp) :: prob_remove_1, prob_remove_2, prob_create_new
790  integer :: info_other_id, new_group
791 
792  call assert(371947172, pt1%id /= pt2%id)
793 
794  ! decide which old particles are to be removed and whether to
795  ! create the resulting coagulated particle
796  r1 = aero_particle_radius(pt1)
797  r2 = aero_particle_radius(pt2)
798  rc = vol2rad(rad2vol(r1) + rad2vol(r2))
799  nc1 = aero_weight_array_num_conc_at_radius(aero_weight_array, c1, r1)
800  nc2 = aero_weight_array_num_conc_at_radius(aero_weight_array, c2, r2)
801  ncc = aero_weight_array_num_conc_at_radius(aero_weight_array, cc, rc)
802  new_group = aero_weight_array_rand_group(aero_weight_array, cc, rc)
803  nc_min = min(nc1, nc2, ncc)
804  prob_remove_1 = nc_min / nc1
805  prob_remove_2 = nc_min / nc2
806  prob_create_new = nc_min / ncc
807  remove_1 = (pmc_random() < prob_remove_1)
808  ! FIXME
809  !if (aero_weight%type == AERO_WEIGHT_TYPE_MFA) then
810  ! remove_2 = .not. remove_1
811  !else
812  remove_2 = (pmc_random() < prob_remove_2)
813  !end if
814  create_new = (pmc_random() < prob_create_new)
815 
816  ! figure out what to do about the ID numbers of the various
817  ! particles --- we try to preserve particle IDs as much as
818  ! possible
819  if (create_new) then
820  id_1_lost = .false.
821  id_2_lost = .false.
822  if (remove_1 .and. remove_2) then
823  if (aero_particle_volume(pt1) > aero_particle_volume(pt2)) then
824  id_2_lost = .true.
825  else
826  id_1_lost = .true.
827  end if
828  end if
829  else
830  id_1_lost = remove_1
831  id_2_lost = remove_2
832  end if
833 
834  ! create a new particle and set its ID
835  if (create_new) then
836  call aero_particle_deallocate(ptc)
837  call aero_particle_allocate_size(ptc, aero_data%n_spec, &
838  aero_data%n_source)
839  call aero_particle_coagulate(pt1, pt2, ptc)
840  call aero_particle_set_weight(ptc, new_group, cc)
841  if (remove_1 .and. (.not. id_1_lost)) then
842  ptc%id = pt1%id
843  call assert(975059559, id_2_lost .eqv. remove_2)
844  elseif (remove_2 .and. (.not. id_2_lost)) then
845  ptc%id = pt2%id
846  call assert(246529753, id_1_lost .eqv. remove_1)
847  else
848  call aero_particle_new_id(ptc)
849  call assert(852038606, id_1_lost .eqv. remove_1)
850  call assert(254018921, id_2_lost .eqv. remove_2)
851  end if
852  info_other_id = ptc%id
853  else
854  info_other_id = 0
855  end if
856 
857  aero_info_1%id = pt1%id
858  aero_info_1%action = aero_info_coag
859  aero_info_1%other_id = info_other_id
860 
861  aero_info_2%id = pt2%id
862  aero_info_2%action = aero_info_coag
863  aero_info_2%other_id = info_other_id
864 
865  end subroutine coagulate_weighting
866 
867 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
868 
869  !> Join together particles \c p1 and \c p2, updating all particle
870  !> and bin structures to reflect the change.
871  subroutine coagulate(aero_data, aero_state, p1, p2, c1, c2, cc)
872 
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  !> First particle index.
878  integer, intent(in) :: p1
879  !> Second particle index.
880  integer, intent(in) :: p2
881  !> First weight class.
882  integer, intent(in) :: c1
883  !> Second weight class.
884  integer, intent(in) :: c2
885  !> Coagulated weight class.
886  integer, intent(in) :: cc
887 
888  type(aero_particle_t), pointer :: pt1, pt2
889  type(aero_particle_t) :: ptc
890  integer :: bn
891  type(aero_info_t) :: aero_info_1, aero_info_2
892  logical :: remove_1, remove_2, create_new, id_1_lost, id_2_lost
893 
894  call aero_particle_allocate(ptc)
895  call aero_info_allocate(aero_info_1)
896  call aero_info_allocate(aero_info_2)
897 
898  pt1 => aero_state%apa%particle(p1)
899  pt2 => aero_state%apa%particle(p2)
900 
901  call coagulate_weighting(pt1, pt2, ptc, c1, c2, cc, aero_data, &
902  aero_state%awa, remove_1, remove_2, create_new, id_1_lost, &
903  id_2_lost, aero_info_1, aero_info_2)
904 
905  ! remove old particles
906  if (p2 > p1) then
907  ! handle a tricky corner case where we have to watch for p2 or
908  ! p1 being the last entry in the array and being repacked when
909  ! the other one is removed
910  if (remove_2) then
911  call aero_state_remove_particle(aero_state, p2, &
912  id_2_lost, aero_info_2)
913  end if
914  if (remove_1) then
915  call aero_state_remove_particle(aero_state, p1, &
916  id_1_lost, aero_info_1)
917  end if
918  else
919  if (remove_1) then
920  call aero_state_remove_particle(aero_state, p1, &
921  id_1_lost, aero_info_1)
922  end if
923  if (remove_2) then
924  call aero_state_remove_particle(aero_state, p2, &
925  id_2_lost, aero_info_2)
926  end if
927  end if
928 
929  ! add new particle
930  if (create_new) then
931  call aero_state_add_particle(aero_state, ptc, allow_resort=.false.)
932  end if
933 
934  call aero_info_deallocate(aero_info_1)
935  call aero_info_deallocate(aero_info_2)
936  call aero_particle_deallocate(ptc)
937 
938  end subroutine coagulate
939 
940 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
941 
942 end module pmc_coagulation
subroutine aero_particle_set_weight(aero_particle, i_group, i_class)
Sets the aerosol particle weight group.
The aero_sorted_t structure and assocated subroutines.
Definition: aero_sorted.F90:9
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.
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.
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_weight_array_check_monotonicity(aero_weight_array, monotone_increasing, monotone_decreasing)
Determine whether all weight functions in an array are monotone increasing, monotone decreasing...
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...
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
integer function rand_binomial(n, p)
Generate a Binomial-distributed random number with the given parameters.
Definition: rand.F90:319
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.
The aero_weight_t structure and associated subroutines.
Definition: aero_weight.F90:9
subroutine max_coag_num_conc_factor(aero_weight_array, bin_grid, i_bin, j_bin, i_class, j_class, ij_class, f_max)
Determine the minimum and maximum number concentration factors for coagulation.
subroutine aero_particle_deallocate(aero_particle)
Deallocates memory associated with an aero_particle_t.
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.
The env_state_t structure and associated subroutines.
Definition: env_state.F90:9
subroutine aero_particle_allocate_size(aero_particle, n_spec, n_source)
Allocates an aero_particle_t of the given size.
Aerosol particle coagulation.
Definition: coagulation.F90:11
subroutine aero_state_add_particle(aero_state, aero_particle, allow_resort)
Add the given particle.
Definition: aero_state.F90:330
subroutine coag_target_with_source(aero_state, bt, ct, target_unif_entry, source_particle, cc)
Coagulate a sampled source particle with a target particle.
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...
elemental real(kind=dp) function aero_particle_radius(aero_particle)
Total radius of the particle (m).
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...
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...
real(kind=dp) function aero_weight_array_num_conc(aero_weight_array, aero_particle)
Compute the number concentration for a particle (m^{-3}).
subroutine compute_n_source(n_part, k_max, del_t, n_source_per_target, accept_factor)
subroutine warn_assert_msg(code, condition_ok, warning_msg)
Prints a warning message if condition_ok is false.
Definition: util.F90:58
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.
Common utility subroutines.
Definition: util.F90:9
integer function aero_sorted_n_bin(aero_sorted)
Returns the number of size bins.
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, c2, i2) that are not the same particle particle as each other.
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_state_dup_particle(aero_state, 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:450
An array of aerosol size distribution weighting functions.
integer function pmc_rand_int(n)
Returns a random integer between 1 and n.
Definition: rand.F90:172
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:389
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:510
The aero_state_t structure and assocated subroutines.
Definition: aero_state.F90:9
Wrapper functions for MPI.
Definition: mpi.F90:13
integer function aero_sorted_n_class(aero_sorted)
Returns the number of weight classes.
subroutine aero_sorted_move_particle(aero_sorted, i_part, new_bin, new_group, new_class)
Move a particle to a different bin and group.
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.
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).
The aero_weight_array_t structure and associated subroutines.
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:43
subroutine aero_particle_allocate(aero_particle)
Allocates memory in an aero_particle_t.
integer function rand_poisson(mean)
Generate a Poisson-distributed random number with the given mean.
Definition: rand.F90:252
The current collection of aerosol particles.
Definition: aero_state.F90:63
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.
1D grid, either logarithmic or linear.
Definition: bin_grid.F90:33
The bin_grid_t structure and associated subroutines.
Definition: bin_grid.F90:9
subroutine aero_particle_new_id(aero_particle)
Assigns a globally-unique new ID number to the particle.
integer function coag_dest_class(aero_weight_array, bin_grid, i_bin, j_bin, i_class, j_class)
Determine the weight class in which coagulated particles will be placed.
Sorting of particles into bins.
Definition: aero_sorted.F90:46
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.
subroutine aero_particle_coagulate(aero_particle_1, aero_particle_2, aero_particle_new)
Coagulate two particles together to make a new one. The new particle will not have its ID set...
real(kind=dp) function pmc_random()
Returns a random number between 0 and 1.
Definition: rand.F90:138
real(kind=dp) elemental function vol2rad(v)
Convert volume (m^3) to radius (m).
Definition: util.F90:238
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. If it happens, do the coagulation and update all structures.
The stats_t type and associated subroutines.
Definition: stats.F90:9
elemental real(kind=dp) function aero_particle_volume(aero_particle)
Total volume of the particle (m^3).
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 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...
integer function aero_sorted_particle_in_bin(aero_sorted, aero_particle)
Find the bin number that contains a given particle.
Generic coagulation kernel.
Definition: coag_kernel.F90:9
subroutine aero_particle_copy(aero_particle_from, aero_particle_to)
Copies a particle.