43 subroutine mc_coag(coag_kernel_type, env_state, aero_data, aero_state, &
44 del_t, tot_n_samp, tot_n_coag)
47 integer,
intent(in) :: coag_kernel_type
55 real(kind=
dp) :: del_t
57 integer,
intent(out) :: tot_n_samp
59 integer,
intent(out) :: tot_n_coag
61 integer :: b1, b2, c1, c2, b2_start
64 if (.not. aero_state%aero_sorted%coag_kernel_bounds_valid)
then
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.
77 if (integer_varray_n_entry( &
78 aero_state%aero_sorted%size_class%inverse(b1, c1)) == 0) &
86 if (integer_varray_n_entry( &
87 aero_state%aero_sorted%size_class%inverse(b2, c2)) == 0) &
90 aero_state, del_t, tot_n_samp, tot_n_coag, b1, b2, c1, c2)
102 aero_state, del_t, tot_n_samp, tot_n_coag, b1, b2, c1, c2)
105 integer,
intent(in) :: coag_kernel_type
113 real(kind=
dp) :: del_t
115 integer,
intent(inout) :: tot_n_samp
117 integer,
intent(inout) :: tot_n_coag
119 integer,
intent(in) :: b1
121 integer,
intent(in) :: b2
123 integer,
intent(in) :: c1
125 integer,
intent(in) :: c2
127 logical :: per_particle_coag_succeeded
129 real(kind=
dp) :: f_max, k_max
132 aero_state%aero_sorted%bin_grid, b1, b2, c1, c2)
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
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
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)
151 aero_data, aero_state, del_t, tot_n_samp, tot_n_coag, b1, b2, c1, c2, &
152 cc, per_particle_coag_succeeded)
155 integer,
intent(in) :: coag_kernel_type
157 real(kind=
dp),
intent(in) :: k_max
165 real(kind=
dp) :: del_t
167 integer,
intent(inout) :: tot_n_samp
169 integer,
intent(inout) :: tot_n_coag
171 integer,
intent(in) :: b1
173 integer,
intent(in) :: b2
175 integer,
intent(in) :: c1
177 integer,
intent(in) :: c2
179 integer,
intent(in) :: cc
181 logical,
intent(inout) :: per_particle_coag_succeeded
183 logical :: correct_weight_ordering
184 integer :: target_unif_entry, target_part, n_samp, n_coag, n_remove, bt, bs
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
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.
200 aero_state%aero_sorted%size_class%inverse(bs, cs)), k_max, del_t, &
201 n_source_per_target, accept_factor)
203 per_particle_coag_succeeded = .false.
208 aero_state%aero_sorted%bin_grid%edges(bt))
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
213 max_target_growth_factor = max_new_target_vol / min_target_vol
215 per_particle_coag_succeeded = .false.
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)
226 target_particle = aero_state%apa%particle(target_part)
228 coag_kernel_type, bs, cs, target_particle, n_source_per_target, &
229 accept_factor, n_samp, n_coag, n_remove, source_particle)
232 target_unif_entry, source_particle, cc)
234 tot_n_samp = tot_n_samp + n_samp
235 tot_n_coag = tot_n_coag + n_coag
239 per_particle_coag_succeeded = .true.
248 b2, c1, c2, cc, bt, bs, ct, cs, correct_weight_ordering)
255 integer,
intent(in) :: b1
257 integer,
intent(in) :: b2
259 integer,
intent(in) :: c1
261 integer,
intent(in) :: c2
263 integer,
intent(in) :: cc
265 integer,
intent(out) :: bt
267 integer,
intent(out) :: bs
269 integer,
intent(out) :: ct
271 integer,
intent(out) :: cs
273 logical,
intent(out) :: correct_weight_ordering
275 real(kind=
dp) :: nc1_min, nc1_max, nc2_min, nc2_max
276 logical :: monotone_increasing, monotone_decreasing
279 monotone_increasing, monotone_decreasing)
280 if (.not. monotone_decreasing)
then
281 correct_weight_ordering = .false.
286 bin_grid%edges(b1), bin_grid%edges(b1 + 1), nc1_min, nc1_max)
288 bin_grid%edges(b2), bin_grid%edges(b2 + 1), nc2_min, nc2_max)
291 correct_weight_ordering = .false.
292 if ((nc1_max < nc2_min) .and. (cc == c1))
then
297 correct_weight_ordering = .true.
299 if ((nc2_max < nc1_min) .and. (cc == c2))
then
304 correct_weight_ordering = .true.
315 integer,
intent(in) :: n_part
317 real(kind=
dp),
intent(in) :: k_max
319 real(kind=
dp),
intent(in) :: del_t
321 real(kind=
dp),
intent(out) :: n_source_per_target
323 real(kind=
dp),
intent(out) :: accept_factor
325 n_source_per_target = k_max * del_t * real(n_part, kind=
dp)
326 accept_factor = 1d0 / k_max
334 coag_kernel_type, bs, cs, coag_particle, n_samp_mean, &
335 accept_factor, n_samp, n_coag, n_remove, source_particle)
344 integer,
intent(in) :: coag_kernel_type
346 integer,
intent(in) :: bs
348 integer,
intent(in) :: cs
350 type(aero_particle_t),
intent(in) :: coag_particle
352 real(kind=
dp),
intent(in) :: n_samp_mean
354 real(kind=
dp),
intent(in) :: accept_factor
356 integer,
intent(out) :: n_samp
358 integer,
intent(out) :: n_coag
360 integer,
intent(out) :: n_remove
362 type(aero_particle_t),
intent(inout) :: source_particle
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
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
379 if (integer_varray_n_entry( &
380 aero_state%aero_sorted%size_class%inverse(bs, cs)) == 0)
then
388 coag_particle, aero_data)
389 target_id = coag_particle%id
390 ct = coag_particle%weight_class
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)
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
407 call aero_particle_zero(source_particle, aero_data)
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) &
419 vol_cv_max = maxval(vol_cv)
420 mean_95_conf_cv = vol_cv_max / sqrt(real(n_avg, kind=
dp)) &
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)
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
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
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
464 prob_coag_mean = prob_coag_tot / real(n_samp, kind=
dp)
466 "kernel upper bound estimation is too tight")
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))
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))
480 new_aero_component(i) = source_particle%component(sample(i))
482 call move_alloc(new_aero_component, source_particle%component)
484 if (n_comp < max_aero_component_size)
then
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
492 new_aero_component(i_comp+i-1) = source_particle%component(i)
494 i_comp = i_comp + n_comp
497 sample_size = n_comp_new - i_comp + 1
500 new_aero_component(i_comp+i-1) = &
501 source_particle%component(sample(i))
503 i_comp = n_comp_new + 1
506 call move_alloc(new_aero_component, source_particle%component)
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 = ' &
514 //
' is less than n_components = ' &
523 target_unif_entry, source_particle, cc)
530 integer,
intent(in) :: bt
532 integer,
intent(in) :: ct
534 integer,
intent(in) :: target_unif_entry
536 type(aero_particle_t),
intent(in) :: source_particle
538 integer,
intent(in) :: cc
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
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
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
554 aero_particle_radius(aero_state%apa%particle(target_part), &
556 call aero_particle_set_weight(aero_state%apa%particle(target_part), &
560 aero_state%apa%particle(target_part), aero_data)
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")
567 new_bin, new_group, cc)
575 aero_state%apa%particle(target_part), aero_data)
577 old_num_conc_target / new_num_conc_target, random_weight_group=.true.)
579 call assert(654300924, aero_state%apa%particle(target_part)%id &
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)
591 integer,
intent(in) :: coag_kernel_type
593 real(kind=
dp),
intent(in) :: k_max
601 real(kind=
dp) :: del_t
603 integer,
intent(inout) :: tot_n_samp
605 integer,
intent(inout) :: tot_n_coag
607 integer,
intent(in) :: b1
609 integer,
intent(in) :: b2
611 integer,
intent(in) :: c1
613 integer,
intent(in) :: c2
615 integer,
intent(in) :: cc
617 real(kind=
dp) :: n_samp_mean, accept_factor
618 integer :: i_samp, n_samp, n1, n2
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
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)) &
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
649 n_samp, accept_factor)
652 integer,
intent(in) :: ni
654 integer,
intent(in) :: nj
656 logical,
intent(in) :: same_bin
658 real(kind=
dp),
intent(in) :: k_max
660 real(kind=
dp),
intent(in) :: del_t
662 real(kind=
dp),
intent(out) :: n_samp_mean
664 integer,
intent(out) :: n_samp
666 real(kind=
dp),
intent(out) :: accept_factor
668 real(kind=
dp) :: r_samp
669 real(kind=
dp) :: n_possible
677 n_possible = real(ni, kind=
dp) * (real(nj, kind=
dp) - 1d0) / 2d0
679 n_possible = real(ni, kind=
dp) * real(nj, kind=
dp)
682 r_samp = k_max * del_t
683 n_samp_mean = r_samp * n_possible
685 accept_factor = 1d0 / k_max
713 c1, c2, cc, coag_kernel_type, accept_factor, did_coag)
722 integer,
intent(in) :: b1
724 integer,
intent(in) :: b2
726 integer,
intent(in) :: c1
728 integer,
intent(in) :: c2
730 integer,
intent(in) :: cc
732 integer,
intent(in) :: coag_kernel_type
734 real(kind=
dp),
intent(in) :: accept_factor
736 logical,
intent(out) :: did_coag
740 real(kind=
dp) :: p, k
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)
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
750 "kernel upper bound estimation is too tight, " &
751 //
"could be caused by changing env_state")
755 call coagulate(aero_data, aero_state, p1, p2, c1, c2, cc)
771 integer,
intent(in) :: b1
773 integer,
intent(in) :: b2
775 integer,
intent(in) :: c1
777 integer,
intent(in) :: c2
779 integer,
intent(out) :: i1
781 integer,
intent(out) :: i2
784 integer_varray_n_entry(aero_sorted%size_class%inverse(b1, c1)) >= 1)
786 integer_varray_n_entry(aero_sorted%size_class%inverse(b1, c1)))
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)
792 integer_varray_n_entry(aero_sorted%size_class%inverse(b2, c2)) - 1)
794 i2 = integer_varray_n_entry(aero_sorted%size_class%inverse(b2, c2))
797 call assert(271635751, integer_varray_n_entry( &
798 aero_sorted%size_class%inverse(b2, c2)) >= 1)
800 integer_varray_n_entry(aero_sorted%size_class%inverse(b2, c2)))
811 aero_weight_array, remove_1, remove_2, create_new, id_1_lost, &
812 id_2_lost, aero_info_1, aero_info_2)
815 type(aero_particle_t),
intent(in) :: pt1
817 type(aero_particle_t),
intent(in) :: pt2
820 type(aero_particle_t),
intent(inout) :: ptc
822 integer,
intent(in) :: c1
824 integer,
intent(in) :: c2
826 integer,
intent(in) :: cc
832 logical,
intent(out) :: remove_1
834 logical,
intent(out) :: remove_2
836 logical,
intent(out) :: create_new
838 logical,
intent(out) :: id_1_lost
840 logical,
intent(out) :: id_2_lost
842 type(aero_info_t),
intent(inout) :: aero_info_1
844 type(aero_info_t),
intent(inout) :: aero_info_2
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
851 call assert(371947172, pt1%id /= pt2%id)
855 r1 = aero_particle_radius(pt1, aero_data)
856 r2 = aero_particle_radius(pt2, aero_data)
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
882 if (remove_1 .and. remove_2)
then
883 if (aero_particle_volume(pt1) > aero_particle_volume(pt2))
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
900 call assert(975059559, id_2_lost .eqv. remove_2)
901 elseif (remove_2 .and. (.not. id_2_lost))
then
903 call assert(246529753, id_1_lost .eqv. remove_1)
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)
909 info_other_id = ptc%id
914 aero_info_1%id = pt1%id
915 aero_info_1%action = aero_info_coag
916 aero_info_1%other_id = info_other_id
918 aero_info_2%id = pt2%id
919 aero_info_2%action = aero_info_coag
920 aero_info_2%other_id = info_other_id
928 subroutine coagulate(aero_data, aero_state, p1, p2, c1, c2, cc)
935 integer,
intent(in) :: p1
937 integer,
intent(in) :: p2
939 integer,
intent(in) :: c1
941 integer,
intent(in) :: c2
943 integer,
intent(in) :: cc
945 type(aero_particle_t) :: ptc
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
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)
962 id_2_lost, aero_info_2)
966 id_1_lost, aero_info_1)
971 id_1_lost, aero_info_1)
975 id_2_lost, aero_info_2)
982 allow_resort=.false.)