PartMC
2.2.0
|
00001 ! Copyright (C) 2005-2012 Nicole Riemer and Matthew West 00002 ! Licensed under the GNU General Public License version 2 or (at your 00003 ! option) any later version. See the file COPYING for details. 00004 00005 !> \file 00006 !> The pmc_coagulation module. 00007 00008 !> Aerosol particle coagulation. 00009 module pmc_coagulation 00010 00011 use pmc_bin_grid 00012 use pmc_aero_data 00013 use pmc_util 00014 use pmc_env_state 00015 use pmc_aero_state 00016 use pmc_aero_weight 00017 use pmc_mpi 00018 use pmc_coag_kernel 00019 use pmc_aero_sorted 00020 #ifdef PMC_USE_MPI 00021 use mpi 00022 #endif 00023 00024 !> Minimum number of coagulation events per large particle for which 00025 !> accelerated coagulation is used. 00026 real(kind=dp), parameter :: COAG_ACCEL_N_EVENT = 1d0 00027 !> Maximum allowed coefficient-of-variation due to undersampling in 00028 !> accelerated coagulation. 00029 real(kind=dp), parameter :: COAG_ACCEL_MAX_CV = 0.1d0 00030 00031 contains 00032 00033 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00034 00035 !> Do coagulation for time del_t. 00036 subroutine mc_coag(coag_kernel_type, env_state, aero_data, aero_state, & 00037 del_t, tot_n_samp, tot_n_coag) 00038 00039 !> Coagulation kernel type. 00040 integer, intent(in) :: coag_kernel_type 00041 !> Environment state. 00042 type(env_state_t), intent(in) :: env_state 00043 !> Aerosol data. 00044 type(aero_data_t), intent(in) :: aero_data 00045 !> Aerosol state. 00046 type(aero_state_t), intent(inout) :: aero_state 00047 !> Timestep for coagulation. 00048 real(kind=dp) :: del_t 00049 !> Total number of samples tested. 00050 integer, intent(out) :: tot_n_samp 00051 !> Number of coagulation events. 00052 integer, intent(out) :: tot_n_coag 00053 00054 integer :: i_bin, j_bin, i_group, j_group, j_bin_start 00055 00056 call aero_state_sort(aero_state) 00057 if (.not. aero_state%aero_sorted%coag_kernel_bounds_valid) then 00058 call est_k_minmax_binned_unweighted(aero_state%aero_sorted%bin_grid, & 00059 coag_kernel_type, aero_data, env_state, & 00060 aero_state%aero_sorted%coag_kernel_min, & 00061 aero_state%aero_sorted%coag_kernel_max) 00062 aero_state%aero_sorted%coag_kernel_bounds_valid = .true. 00063 end if 00064 00065 tot_n_samp = 0 00066 tot_n_coag = 0 00067 !do i_group = 1,size(aero_state%aero_sorted%group) 00068 !do j_group = 1,i_group 00069 do i_bin = 1,size(aero_state%aero_sorted%size%inverse) 00070 if (aero_state%aero_sorted%size%inverse(i_bin)%n_entry == 0) & 00071 cycle 00072 !if (i_group == j_group) then 00073 ! j_bin_start = i_bin 00074 !else 00075 ! j_bin_start = 1 00076 !end if 00077 do j_bin = i_bin,size(aero_state%aero_sorted%size%inverse) 00078 if (aero_state%aero_sorted%size%inverse(j_bin)%n_entry == 0) & 00079 cycle 00080 call mc_coag_for_bin(coag_kernel_type, env_state, aero_data, & 00081 aero_state, del_t, tot_n_samp, tot_n_coag, & 00082 i_bin, j_bin) 00083 end do 00084 end do 00085 !end do 00086 !end do 00087 00088 end subroutine mc_coag 00089 00090 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00091 00092 !> Do coagulation for time del_t for the given bins. 00093 subroutine mc_coag_for_bin(coag_kernel_type, env_state, aero_data, & 00094 aero_state, del_t, tot_n_samp, tot_n_coag, i_bin, j_bin) 00095 00096 !> Coagulation kernel type. 00097 integer, intent(in) :: coag_kernel_type 00098 !> Environment state. 00099 type(env_state_t), intent(in) :: env_state 00100 !> Aerosol data. 00101 type(aero_data_t), intent(in) :: aero_data 00102 !> Aerosol state. 00103 type(aero_state_t), intent(inout) :: aero_state 00104 !> Timestep for coagulation. 00105 real(kind=dp) :: del_t 00106 !> Total number of samples tested. 00107 integer, intent(inout) :: tot_n_samp 00108 !> Number of coagulation events. 00109 integer, intent(inout) :: tot_n_coag 00110 !> First bin number. 00111 integer, intent(in) :: i_bin 00112 !> First weight group number. 00113 !integer, intent(in) :: i_group 00114 !> Second bin number. 00115 integer, intent(in) :: j_bin 00116 !> Second weight group number. 00117 !integer, intent(in) :: j_group 00118 00119 logical :: per_particle_coag_succeeded 00120 real(kind=dp) :: f_max, k_max 00121 00122 call max_coag_num_conc_factor_better(aero_state%aero_weight, & 00123 aero_state%aero_sorted%bin_grid, i_bin, j_bin, f_max) 00124 k_max = aero_state%aero_sorted%coag_kernel_max(i_bin, j_bin) * f_max 00125 00126 call try_per_particle_coag(coag_kernel_type, k_max, env_state, aero_data, & 00127 aero_state, del_t, tot_n_samp, tot_n_coag, i_bin, j_bin, & 00128 per_particle_coag_succeeded) 00129 if (per_particle_coag_succeeded) return 00130 00131 call per_set_coag(coag_kernel_type, k_max, env_state, aero_data, & 00132 aero_state, del_t, tot_n_samp, tot_n_coag, i_bin, j_bin) 00133 00134 end subroutine mc_coag_for_bin 00135 00136 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00137 00138 !> Attempt per-particle coagulation. 00139 subroutine try_per_particle_coag(coag_kernel_type, k_max, env_state, & 00140 aero_data, aero_state, del_t, tot_n_samp, tot_n_coag, i_bin, j_bin, & 00141 per_particle_coag_succeeded) 00142 00143 !> Coagulation kernel type. 00144 integer, intent(in) :: coag_kernel_type 00145 !> Maximum coagulation kernel (s^{-1} m^3). 00146 real(kind=dp), intent(in) :: k_max 00147 !> Environment state. 00148 type(env_state_t), intent(in) :: env_state 00149 !> Aerosol data. 00150 type(aero_data_t), intent(in) :: aero_data 00151 !> Aerosol state. 00152 type(aero_state_t), intent(inout) :: aero_state 00153 !> Timestep for coagulation. 00154 real(kind=dp) :: del_t 00155 !> Total number of samples tested. 00156 integer, intent(inout) :: tot_n_samp 00157 !> Number of coagulation events. 00158 integer, intent(inout) :: tot_n_coag 00159 !> First bin number. 00160 integer, intent(in) :: i_bin 00161 !> First weight group number. 00162 !integer, intent(in) :: i_group 00163 !> Second bin number. 00164 integer, intent(in) :: j_bin 00165 !> Second weight group number. 00166 !integer, intent(in) :: j_group 00167 !> Whether we succeeded in doing per-particle coag. 00168 logical, intent(inout) :: per_particle_coag_succeeded 00169 00170 logical :: correct_weight_ordering 00171 integer :: target_unif_entry, target_part, n_samp, n_coag, n_remove 00172 integer :: target_bin, source_bin 00173 real(kind=dp) :: n_source_per_target, accept_factor 00174 type(aero_particle_t) :: target_particle, source_particle 00175 00176 call determine_target_and_source(aero_state%aero_weight, & 00177 aero_state%aero_sorted%bin_grid, i_bin, j_bin, target_bin, & 00178 source_bin, correct_weight_ordering) 00179 if (.not. correct_weight_ordering) then 00180 per_particle_coag_succeeded = .false. 00181 return 00182 end if 00183 00184 call compute_n_source( & 00185 aero_state%aero_sorted%size%inverse(source_bin)%n_entry, & 00186 k_max, del_t, n_source_per_target, accept_factor) 00187 if (n_source_per_target < COAG_ACCEL_N_EVENT) then 00188 per_particle_coag_succeeded = .false. 00189 return 00190 end if 00191 00192 call aero_particle_allocate(target_particle) 00193 call aero_particle_allocate(source_particle) 00194 00195 ! work backwards to avoid particle movement issues 00196 do target_unif_entry & 00197 = aero_state%aero_sorted%size%inverse(target_bin)%n_entry,1,-1 00198 target_part & 00199 = aero_state%aero_sorted%size%inverse(target_bin)%entry(target_unif_entry) 00200 ! need to copy coag_particle as the underlying storage may be 00201 ! rearranged due to removals 00202 call aero_particle_copy(aero_state%apa%particle(target_part), & 00203 target_particle) 00204 call sample_source_particle(aero_state, aero_data, env_state, & 00205 coag_kernel_type, source_bin, target_particle, & 00206 n_source_per_target, accept_factor, n_samp, n_coag, n_remove, & 00207 source_particle) 00208 if (n_coag > 0) then 00209 call coag_target_with_source(aero_state, target_bin, & 00210 target_unif_entry, source_particle) 00211 end if 00212 tot_n_samp = tot_n_samp + n_samp 00213 tot_n_coag = tot_n_coag + n_coag 00214 ! we discard n_remove information at present 00215 end do 00216 00217 call aero_particle_deallocate(target_particle) 00218 call aero_particle_deallocate(source_particle) 00219 00220 per_particle_coag_succeeded = .true. 00221 00222 end subroutine try_per_particle_coag 00223 00224 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00225 00226 !> Determine the source and target particle bin/group for 00227 !> per-particle coagulation, if possible. 00228 subroutine determine_target_and_source(aero_weight_array, bin_grid, i_bin, & 00229 j_bin, target_bin, source_bin, correct_weight_ordering) 00230 00231 !> Aero weight array. 00232 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00233 !> Bin grid. 00234 type(bin_grid_t), intent(in) :: bin_grid 00235 !> First bin number. 00236 integer, intent(in) :: i_bin 00237 !> Second bin number. 00238 integer, intent(in) :: j_bin 00239 !> Target bin number. 00240 integer, intent(out) :: target_bin 00241 !> Source bin number. 00242 integer, intent(out) :: source_bin 00243 !> Whether the weight ordering is correct for per-particle coagulation. 00244 logical, intent(out) :: correct_weight_ordering 00245 00246 real(kind=dp) :: i_nc_min, i_nc_max, j_nc_min, j_nc_max 00247 logical :: monotone_increasing, monotone_decreasing 00248 00249 call aero_weight_array_check_monotonicity(aero_weight_array, & 00250 monotone_increasing, monotone_decreasing) 00251 if (.not. monotone_decreasing) then 00252 correct_weight_ordering = .false. 00253 return 00254 end if 00255 00256 call aero_weight_array_minmax_num_conc(aero_weight_array, & 00257 bin_grid%edge_radius(i_bin), bin_grid%edge_radius(i_bin + 1), & 00258 i_nc_min, i_nc_max) 00259 call aero_weight_array_minmax_num_conc(aero_weight_array, & 00260 bin_grid%edge_radius(j_bin), bin_grid%edge_radius(j_bin + 1), & 00261 j_nc_min, j_nc_max) 00262 00263 ! we have already confirmed monotone_decreasing weights above 00264 correct_weight_ordering = .false. 00265 if (i_nc_max < j_nc_min) then 00266 target_bin = i_bin 00267 source_bin = j_bin 00268 correct_weight_ordering = .true. 00269 end if 00270 if (j_nc_max < i_nc_min) then 00271 target_bin = j_bin 00272 source_bin = i_bin 00273 correct_weight_ordering = .true. 00274 end if 00275 00276 end subroutine determine_target_and_source 00277 00278 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00279 00280 subroutine compute_n_source(n_part, k_max, del_t, n_source_per_target, & 00281 accept_factor) 00282 00283 !> Number of particles available as partners. 00284 integer, intent(in) :: n_part 00285 !> Maximum coagulation kernel (s^{-1} m^3). 00286 real(kind=dp), intent(in) :: k_max 00287 !> Timestep (s). 00288 real(kind=dp), intent(in) :: del_t 00289 !> Mean number of source particles per target particle. 00290 real(kind=dp), intent(out) :: n_source_per_target 00291 !> Accept factor for samples. 00292 real(kind=dp), intent(out) :: accept_factor 00293 00294 n_source_per_target = k_max * del_t * real(n_part, kind=dp) 00295 accept_factor = 1d0 / k_max 00296 00297 end subroutine compute_n_source 00298 00299 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00300 00301 !> Sample coagulation partners for a single coagulation event. 00302 subroutine sample_source_particle(aero_state, aero_data, env_state, & 00303 coag_kernel_type, source_bin, coag_particle, n_samp_mean, & 00304 accept_factor, n_samp, n_coag, n_remove, source_particle) 00305 00306 !> Aerosol state. 00307 type(aero_state_t), intent(inout) :: aero_state 00308 !> Aerosol data. 00309 type(aero_data_t), intent(in) :: aero_data 00310 !> Environment state. 00311 type(env_state_t), intent(in) :: env_state 00312 !> Coagulation kernel type. 00313 integer, intent(in) :: coag_kernel_type 00314 !> Bin to sample particles from. 00315 integer, intent(in) :: source_bin 00316 !> Aerosol particle that coagulation will be with. 00317 type(aero_particle_t), intent(in) :: coag_particle 00318 !> Mean number of samples to use. 00319 real(kind=dp), intent(in) :: n_samp_mean 00320 !> Probability factor of accepting samples. 00321 real(kind=dp), intent(in) :: accept_factor 00322 !> Number of samples used. 00323 integer, intent(out) :: n_samp 00324 !> Number of coagulations performed. 00325 integer, intent(out) :: n_coag 00326 !> Number of particles removed. 00327 integer, intent(out) :: n_remove 00328 !> Sampled average coagulation partner particle. 00329 type(aero_particle_t), intent(inout) :: source_particle 00330 00331 real(kind=dp) :: prob_remove_i, prob_remove_source_max 00332 real(kind=dp) :: prob_coag, prob_coag_tot, prob_coag_mean 00333 real(kind=dp) :: num_conc_i, num_conc_source_min, num_conc_target, k 00334 real(kind=dp) :: vol_sq(aero_data%n_spec), vol_mean(aero_data%n_spec) 00335 real(kind=dp) :: vol_cv(aero_data%n_spec), vol_cv_max, mean_95_conf_cv 00336 integer :: n_samp_remove, n_samp_extra, n_samp_total, n_avg, i_samp 00337 integer :: i_unif_entry, i_part, target_id, new_bin 00338 type(aero_particle_t), pointer :: i_particle 00339 type(aero_info_t) :: aero_info 00340 00341 if (aero_state%aero_sorted%size%inverse(source_bin)%n_entry == 0) then 00342 n_samp = 0 00343 n_remove = 0 00344 n_coag = 0 00345 return 00346 end if 00347 00348 num_conc_target = aero_weight_array_num_conc( & 00349 aero_state%aero_weight, coag_particle) 00350 target_id = coag_particle%id 00351 00352 num_conc_source_min = aero_weight_array_num_conc_at_radius( & 00353 aero_state%aero_weight, & 00354 aero_state%aero_sorted%bin_grid%edge_radius(source_bin)) 00355 prob_remove_source_max = num_conc_target / num_conc_source_min 00356 call assert(653606684, prob_remove_source_max <= 1d0) 00357 00358 n_samp_remove = rand_poisson(prob_remove_source_max * n_samp_mean) 00359 n_samp_extra = rand_poisson((1d0 - prob_remove_source_max) * n_samp_mean) 00360 n_samp_total = n_samp_remove + n_samp_extra 00361 00362 n_avg = 0 00363 n_samp = 0 00364 n_remove = 0 00365 prob_coag_tot = 0d0 00366 call aero_particle_deallocate(source_particle) 00367 call aero_particle_allocate_size(source_particle, & 00368 aero_data%n_spec, aero_data%n_source) 00369 vol_sq = 0d0 00370 00371 ! FIXME: Can't we just do n_samp = 1,n_samp_total and shift tests 00372 ! to the end? 00373 do i_samp = 1,n_samp_total 00374 if (aero_state%aero_sorted%size%inverse(source_bin)%n_entry == 0) & 00375 exit 00376 if ((n_samp > n_samp_remove) .and. (n_avg >= 2)) then 00377 vol_mean = source_particle%vol / real(n_avg, kind=dp) 00378 where(vol_mean > 0d0) & 00379 vol_cv = sqrt(vol_sq / real(n_avg, kind=dp) - (vol_mean)**2) 00380 / vol_mean 00381 vol_cv_max = maxval(vol_cv) 00382 mean_95_conf_cv = vol_cv_max / sqrt(real(n_avg, kind=dp)) 00383 * student_t_95_coeff(n_avg) 00384 ! FIXME: We are using just the max of the diagonal of the 00385 ! covariance matrix. Is this well-justified? 00386 if (mean_95_conf_cv < COAG_ACCEL_MAX_CV) exit 00387 end if 00388 n_samp = n_samp + 1 00389 ! FIXME: We are sampling with replacement. Is this a problem? 00390 i_unif_entry & 00391 = pmc_rand_int(aero_state%aero_sorted%size%inverse(source_bin)%n_entry) 00392 i_part = aero_state%aero_sorted%size%inverse(source_bin)%entry(i_unif_entry) 00393 i_particle => aero_state%apa%particle(i_part) 00394 ! re-get j_part as particle ordering may be changing 00395 call num_conc_weighted_kernel(coag_kernel_type, i_particle, & 00396 coag_particle, aero_data, aero_state%aero_weight, env_state, k) 00397 prob_coag = k * accept_factor 00398 prob_coag_tot = prob_coag_tot + prob_coag 00399 if (pmc_random() < prob_coag) then 00400 n_avg = n_avg + 1 00401 call aero_particle_coagulate(source_particle, i_particle, & 00402 source_particle) 00403 vol_sq = vol_sq + i_particle%vol**2 00404 if (i_samp <= n_samp_remove) then 00405 num_conc_i = aero_weight_array_num_conc(aero_state%aero_weight, & 00406 i_particle) 00407 prob_remove_i = num_conc_target / num_conc_i 00408 if (pmc_random() < prob_remove_i / prob_remove_source_max) then 00409 n_remove = n_remove + 1 00410 call aero_info_allocate(aero_info) 00411 aero_info%id = i_particle%id 00412 aero_info%action = AERO_INFO_COAG 00413 aero_info%other_id = target_id 00414 call aero_state_remove_particle_with_info(aero_state, & 00415 i_part, aero_info) 00416 call aero_info_deallocate(aero_info) 00417 end if 00418 end if 00419 end if 00420 end do 00421 00422 if (n_avg == 0) then 00423 n_coag = 0 00424 return 00425 end if 00426 00427 prob_coag_mean = prob_coag_tot / real(n_samp, kind=dp) ! note not n_avg 00428 call warn_assert_msg(383983415, prob_coag_mean <= 1d0, & 00429 "kernel upper bound estimation is too tight") 00430 n_coag = rand_binomial(n_samp_total, prob_coag_mean) 00431 source_particle%vol = source_particle%vol & 00432 * (real(n_coag, kind=dp) / real(n_avg, kind=dp)) 00433 00434 end subroutine sample_source_particle 00435 00436 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00437 00438 !> Coagulate a sampled source particle with a target particle. 00439 subroutine coag_target_with_source(aero_state, target_bin, & 00440 target_unif_entry, source_particle) 00441 00442 !> Aerosol state. 00443 type(aero_state_t), intent(inout) :: aero_state 00444 !> Bin of coagulating particle. 00445 integer, intent(in) :: target_bin 00446 !> Entry-in-bin of coagulating particle. 00447 integer, intent(in) :: target_unif_entry 00448 !> Sampled particle to coagulate with. 00449 type(aero_particle_t), intent(in) :: source_particle 00450 00451 integer :: target_part, target_id, new_bin, new_group 00452 real(kind=dp) :: old_num_conc_target, new_num_conc_target 00453 00454 target_part & 00455 = aero_state%aero_sorted%size%inverse(target_bin)%entry(target_unif_entry) 00456 target_id = aero_state%apa%particle(target_part)%id 00457 old_num_conc_target & 00458 = aero_weight_array_num_conc(aero_state%aero_weight, & 00459 aero_state%apa%particle(target_part)) 00460 call aero_particle_coagulate(aero_state%apa%particle(target_part), & 00461 source_particle, aero_state%apa%particle(target_part)) 00462 aero_state%apa%particle(target_part)%id = target_id 00463 ! assign to a randomly chosen group 00464 new_group = aero_weight_array_rand_group(aero_state%aero_weight, & 00465 aero_particle_radius(aero_state%apa%particle(target_part))) 00466 call aero_particle_set_group(aero_state%apa%particle(target_part), new_group) 00467 ! fix bin due to composition changes 00468 new_bin = aero_sorted_particle_in_bin(aero_state%aero_sorted, & 00469 aero_state%apa%particle(target_part)) 00470 if ((new_bin < 1) & 00471 .or. (new_bin > aero_state%aero_sorted%bin_grid%n_bin)) then 00472 call die_msg(765620746, "particle outside of bin_grid: " & 00473 // "try reducing the timestep del_t") 00474 end if 00475 call aero_sorted_move_particle(aero_state%aero_sorted, target_part, & 00476 new_bin, new_group) 00477 ! Adjust particle number to account for weight changes 00478 ! target_bin/target_group/target_entry are invalid, but 00479 ! target_part is still good. We are treating all particles in all 00480 ! groups together, and randomly reassigning between groups above, 00481 ! so here we can't use aero_state_reweight_particle(), as that 00482 ! assumes we are staying in the same weight group. 00483 new_num_conc_target & 00484 = aero_weight_array_num_conc(aero_state%aero_weight, & 00485 aero_state%apa%particle(target_part)) 00486 call aero_state_dup_particle(aero_state, target_part, & 00487 old_num_conc_target / new_num_conc_target, random_weight_group=.true.) 00488 ! we should only be doing this for decreasing weights 00489 call assert(654300924, aero_state%apa%particle(target_part)%id == target_id) 00490 00491 end subroutine coag_target_with_source 00492 00493 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00494 00495 !> Do set-wise coagulation. 00496 subroutine per_set_coag(coag_kernel_type, k_max, env_state, aero_data, & 00497 aero_state, del_t, tot_n_samp, tot_n_coag, i_bin, j_bin) 00498 00499 !> Coagulation kernel type. 00500 integer, intent(in) :: coag_kernel_type 00501 !> Maximum coagulation kernel (s^{-1} m^3). 00502 real(kind=dp), intent(in) :: k_max 00503 !> Environment state. 00504 type(env_state_t), intent(in) :: env_state 00505 !> Aerosol data. 00506 type(aero_data_t), intent(in) :: aero_data 00507 !> Aerosol state. 00508 type(aero_state_t), intent(inout) :: aero_state 00509 !> Timestep for coagulation. 00510 real(kind=dp) :: del_t 00511 !> Total number of samples tested. 00512 integer, intent(inout) :: tot_n_samp 00513 !> Number of coagulation events. 00514 integer, intent(inout) :: tot_n_coag 00515 !> First bin number. 00516 integer, intent(in) :: i_bin 00517 !> Second bin number. 00518 integer, intent(in) :: j_bin 00519 00520 real(kind=dp) :: n_samp_mean, accept_factor 00521 integer :: i_samp, n_samp 00522 logical :: did_coag 00523 00524 call compute_n_samp( & 00525 aero_state%aero_sorted%size%inverse(i_bin)%n_entry, & 00526 aero_state%aero_sorted%size%inverse(j_bin)%n_entry, & 00527 (i_bin == j_bin), k_max, del_t, n_samp_mean, n_samp, accept_factor) 00528 tot_n_samp = tot_n_samp + n_samp 00529 00530 do i_samp = 1,n_samp 00531 ! check we still have enough particles to coagulate 00532 if (((aero_state%aero_sorted%size%inverse(i_bin)%n_entry < 2) & 00533 .and. (i_bin == j_bin)) & 00534 .or. (aero_state%aero_sorted%size%inverse(i_bin)%n_entry < 1) & 00535 .or. (aero_state%aero_sorted%size%inverse(j_bin)%n_entry < 1)) & 00536 exit 00537 call maybe_coag_pair(env_state, aero_data, aero_state, i_bin, j_bin, & 00538 coag_kernel_type, accept_factor, did_coag) 00539 if (did_coag) tot_n_coag = tot_n_coag + 1 00540 end do 00541 00542 end subroutine per_set_coag 00543 00544 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00545 00546 !> Compute the number of samples required for the pair of bins. 00547 subroutine compute_n_samp(ni, nj, same_bin, k_max, del_t, n_samp_mean, & 00548 n_samp, accept_factor) 00549 00550 !> Number particles in first bin. 00551 integer, intent(in) :: ni 00552 !> Number particles in second bin. 00553 integer, intent(in) :: nj 00554 !> Whether first bin is second bin. 00555 logical, intent(in) :: same_bin 00556 !> Maximum kernel value (s^{-1} m^3). 00557 real(kind=dp), intent(in) :: k_max 00558 !> Timestep (s). 00559 real(kind=dp), intent(in) :: del_t 00560 !> Mean number of samples per timestep. 00561 real(kind=dp), intent(out) :: n_samp_mean 00562 !> Number of samples per timestep. 00563 integer, intent(out) :: n_samp 00564 !> Scale factor for accept probability (1). 00565 real(kind=dp), intent(out) :: accept_factor 00566 00567 real(kind=dp) :: r_samp 00568 real(kind=dp) :: n_possible ! use real(kind=dp) to avoid integer overflow 00569 ! could use integer*8 or integer(kind = 8) 00570 ! or di = selected_int_kind(18), integer(kind=di) 00571 ! to represent 10^{-18} to 10^{18} 00572 00573 if (same_bin) then 00574 ! don't change this to ni * (ni - 1) as the ni/nj distinction 00575 ! is important for coagulation_dist, which also calls this 00576 n_possible = real(ni, kind=dp) * (real(nj, kind=dp) - 1d0) / 2d0 00577 else 00578 n_possible = real(ni, kind=dp) * real(nj, kind=dp) 00579 end if 00580 00581 r_samp = k_max * del_t 00582 n_samp_mean = r_samp * n_possible 00583 n_samp = rand_poisson(n_samp_mean) 00584 accept_factor = 1d0 / k_max 00585 00586 ! possible variants: 00587 ! A: accept_factor = 1d0 / k_max 00588 ! B: accept_factor = del_t * n_possible & 00589 ! / (real(n_samp, kind=dp) * comp_vol) 00590 ! timings of test suite as of 2010-12-22T17:12:14-0600: 00591 ! A with n_samp = prob_round(n_samp_mean): 00592 ! 159.82 162.18 156.28 00593 ! A with n_samp = rand_poisson(n_samp_mean): 00594 ! 151.93 158.38 174.74 157.65 00595 ! B with n_samp = ceiling(n_samp_mean): 00596 ! 196.06 200.23 192.41 00597 ! B with n_samp = ceiling(n_samp_mean + 1 * sqrt(n_samp_mean)): 00598 ! 189.78 211.12 195.19 00599 ! B with n_samp = ceiling(n_samp_mean + 3 * sqrt(n_samp_mean)): 00600 ! 214.60 201.25 203.55 00601 00602 end subroutine compute_n_samp 00603 00604 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00605 00606 !> Choose a random pair for potential coagulation and test its 00607 !> probability of coagulation. If it happens, do the coagulation and 00608 !> update all structures. 00609 !! 00610 !! The probability of a coagulation will be taken as <tt>(kernel / 00611 !! k_max)</tt>. 00612 subroutine maybe_coag_pair(env_state, aero_data, aero_state, b1, b2, & 00613 coag_kernel_type, accept_factor, did_coag) 00614 00615 !> Environment state. 00616 type(env_state_t), intent(in) :: env_state 00617 !> Aerosol data. 00618 type(aero_data_t), intent(in) :: aero_data 00619 !> Aerosol state. 00620 type(aero_state_t), intent(inout) :: aero_state 00621 !> Bin of first particle. 00622 integer, intent(in) :: b1 00623 !> Bin of second particle. 00624 integer, intent(in) :: b2 00625 !> Coagulation kernel type. 00626 integer, intent(in) :: coag_kernel_type 00627 !> Scale factor for accept probability (1). 00628 real(kind=dp), intent(in) :: accept_factor 00629 !> Whether a coagulation occured. 00630 logical, intent(out) :: did_coag 00631 00632 integer :: s1, s2 00633 integer :: p1, p2 00634 real(kind=dp) :: p, k 00635 00636 call find_rand_pair(aero_state%aero_sorted, b1, b2, s1, s2) 00637 p1 = aero_state%aero_sorted%size%inverse(b1)%entry(s1) 00638 p2 = aero_state%aero_sorted%size%inverse(b2)%entry(s2) 00639 call num_conc_weighted_kernel(coag_kernel_type, & 00640 aero_state%apa%particle(p1), aero_state%apa%particle(p2), aero_data, & 00641 aero_state%aero_weight, env_state, k) 00642 p = k * accept_factor 00643 call warn_assert_msg(532446093, p <= 1d0, & 00644 "kernel upper bound estimation is too tight, " & 00645 // "could be caused by changing env_state") 00646 00647 did_coag = .false. 00648 if (pmc_random() .lt. p) then 00649 call coagulate(aero_data, aero_state, p1, p2) 00650 did_coag = .true. 00651 end if 00652 00653 end subroutine maybe_coag_pair 00654 00655 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00656 00657 !> Given bins b1 and b2, find a random pair of particles (b1, s1) 00658 !> and (b2, s2) that are not the same particle particle as each 00659 !> other. 00660 subroutine find_rand_pair(aero_sorted, b1, b2, s1, s2) 00661 00662 !> Aerosol sorted data. 00663 type(aero_sorted_t), intent(in) :: aero_sorted 00664 !> Bin number of first particle. 00665 integer, intent(in) :: b1 00666 !> Bin number of second particle. 00667 integer, intent(in) :: b2 00668 !> First rand particle. 00669 integer, intent(out) :: s1 00670 !> Second rand particle. 00671 integer, intent(out) :: s2 00672 00673 call assert(619608562, aero_sorted%size%inverse(b1)%n_entry >= 1) 00674 s1 = pmc_rand_int(aero_sorted%size%inverse(b1)%n_entry) 00675 00676 if (b1 == b2) then 00677 call assert(956184336, aero_sorted%size%inverse(b2)%n_entry >= 2) 00678 s2 = pmc_rand_int(aero_sorted%size%inverse(b2)%n_entry - 1) 00679 if (s2 == s1) then 00680 s2 = aero_sorted%size%inverse(b2)%n_entry 00681 end if 00682 else 00683 call assert(271635751, aero_sorted%size%inverse(b2)%n_entry >= 1) 00684 s2 = pmc_rand_int(aero_sorted%size%inverse(b2)%n_entry) 00685 end if 00686 00687 end subroutine find_rand_pair 00688 00689 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00690 00691 !> Actually coagulate particle_1 and particle_2 to form particle_new 00692 !> and compute weighting effects, including which particles should 00693 !> be lost and which gained. 00694 subroutine coagulate_weighting(particle_1, particle_2, particle_new, & 00695 aero_data, aero_weight_array, remove_1, remove_2, create_new, & 00696 id_1_lost, id_2_lost, aero_info_1, aero_info_2) 00697 00698 !> First coagulating aerosol particle. 00699 type(aero_particle_t), intent(in) :: particle_1 00700 !> Second coagulating aerosol particle. 00701 type(aero_particle_t), intent(in) :: particle_2 00702 !> Combined aerosol particle resulting from coagulation of particle_1 00703 !> and particle_2. 00704 type(aero_particle_t), intent(inout) :: particle_new 00705 !> Aerosol data. 00706 type(aero_data_t), intent(in) :: aero_data 00707 !> Aerosol weight array. 00708 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00709 !> Whether to remove particle_1. 00710 logical, intent(out) :: remove_1 00711 !> Whether to remove particle_2. 00712 logical, intent(out) :: remove_2 00713 !> Whether to create particle_new. 00714 logical, intent(out) :: create_new 00715 !> Whether the ID of particle_1 will be lost due to coagulation. 00716 logical, intent(out) :: id_1_lost 00717 !> Whether the ID of particle_2 will be lost due to coagulation. 00718 logical, intent(out) :: id_2_lost 00719 !> The removal information associated with particle_1. 00720 type(aero_info_t), intent(inout) :: aero_info_1 00721 !> The removal information associated with particle_2. 00722 type(aero_info_t), intent(inout) :: aero_info_2 00723 00724 real(kind=dp) :: radius_1, radius_2, radius_new 00725 real(kind=dp) :: num_conc_min, num_conc_1, num_conc_2 00726 real(kind=dp) :: num_conc_new 00727 real(kind=dp) :: prob_remove_1, prob_remove_2, prob_create_new 00728 integer :: info_other_id, new_group 00729 00730 call assert(371947172, particle_1%id /= particle_2%id) 00731 00732 ! decide which old particles are to be removed and whether to 00733 ! create the resulting coagulated particle 00734 radius_1 = aero_particle_radius(particle_1) 00735 radius_2 = aero_particle_radius(particle_2) 00736 radius_new = vol2rad(rad2vol(radius_1) + rad2vol(radius_2)) 00737 num_conc_1 = aero_weight_array_num_conc_at_radius(aero_weight_array, & 00738 radius_1) 00739 num_conc_2 = aero_weight_array_num_conc_at_radius(aero_weight_array, & 00740 radius_2) 00741 num_conc_new = aero_weight_array_num_conc_at_radius(aero_weight_array, & 00742 radius_new) 00743 new_group = aero_weight_array_rand_group(aero_weight_array, radius_new) 00744 num_conc_min = min(num_conc_1, num_conc_2, num_conc_new) 00745 prob_remove_1 = num_conc_min / num_conc_1 00746 prob_remove_2 = num_conc_min / num_conc_2 00747 prob_create_new = num_conc_min / num_conc_new 00748 remove_1 = (pmc_random() < prob_remove_1) 00749 ! FIXME 00750 !if (aero_weight%type == AERO_WEIGHT_TYPE_MFA) then 00751 ! remove_2 = .not. remove_1 00752 !else 00753 remove_2 = (pmc_random() < prob_remove_2) 00754 !end if 00755 create_new = (pmc_random() < prob_create_new) 00756 00757 ! figure out what to do about the ID numbers of the various 00758 ! particles --- we try to preserve particle IDs as much as 00759 ! possible 00760 if (create_new) then 00761 id_1_lost = .false. 00762 id_2_lost = .false. 00763 if (remove_1 .and. remove_2) then 00764 if (aero_particle_volume(particle_1) & 00765 > aero_particle_volume(particle_2)) then 00766 id_2_lost = .true. 00767 else 00768 id_1_lost = .true. 00769 end if 00770 end if 00771 else 00772 id_1_lost = remove_1 00773 id_2_lost = remove_2 00774 end if 00775 00776 ! create a new particle and set its ID 00777 if (create_new) then 00778 call aero_particle_deallocate(particle_new) 00779 call aero_particle_allocate_size(particle_new, aero_data%n_spec, & 00780 aero_data%n_source) 00781 call aero_particle_coagulate(particle_1, particle_2, particle_new) 00782 particle_new%weight_group = new_group 00783 if (remove_1 .and. (.not. id_1_lost)) then 00784 particle_new%id = particle_1%id 00785 call assert(975059559, id_2_lost .eqv. remove_2) 00786 elseif (remove_2 .and. (.not. id_2_lost)) then 00787 particle_new%id = particle_2%id 00788 call assert(246529753, id_1_lost .eqv. remove_1) 00789 else 00790 call aero_particle_new_id(particle_new) 00791 call assert(852038606, id_1_lost .eqv. remove_1) 00792 call assert(254018921, id_2_lost .eqv. remove_2) 00793 end if 00794 info_other_id = particle_new%id 00795 else 00796 info_other_id = 0 00797 end if 00798 00799 aero_info_1%id = particle_1%id 00800 aero_info_1%action = AERO_INFO_COAG 00801 aero_info_1%other_id = info_other_id 00802 00803 aero_info_2%id = particle_2%id 00804 aero_info_2%action = AERO_INFO_COAG 00805 aero_info_2%other_id = info_other_id 00806 00807 end subroutine coagulate_weighting 00808 00809 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00810 00811 !> Join together particles \c p1 and \c p2, updating all 00812 !> particle and bin structures to reflect the change. 00813 subroutine coagulate(aero_data, aero_state, p1, p2) 00814 00815 !> Aerosol data. 00816 type(aero_data_t), intent(in) :: aero_data 00817 !> Aerosol state. 00818 type(aero_state_t), intent(inout) :: aero_state 00819 !> First particle index. 00820 integer, intent(in) :: p1 00821 !> Second particle index. 00822 integer, intent(in) :: p2 00823 00824 type(aero_particle_t), pointer :: particle_1, particle_2 00825 type(aero_particle_t) :: particle_new 00826 integer :: bn 00827 type(aero_info_t) :: aero_info_1, aero_info_2 00828 logical :: remove_1, remove_2, create_new, id_1_lost, id_2_lost 00829 00830 call aero_particle_allocate(particle_new) 00831 call aero_info_allocate(aero_info_1) 00832 call aero_info_allocate(aero_info_2) 00833 00834 particle_1 => aero_state%apa%particle(p1) 00835 particle_2 => aero_state%apa%particle(p2) 00836 00837 call coagulate_weighting(particle_1, particle_2, particle_new, & 00838 aero_data, aero_state%aero_weight, remove_1, remove_2, & 00839 create_new, id_1_lost, id_2_lost, aero_info_1, aero_info_2) 00840 00841 ! remove old particles 00842 if (p2 > p1) then 00843 ! handle a tricky corner case where we have to watch for p2 or 00844 ! p1 being the last entry in the array and being repacked when 00845 ! the other one is removed 00846 if (remove_2) then 00847 call aero_state_remove_particle(aero_state, p2, & 00848 id_2_lost, aero_info_2) 00849 end if 00850 if (remove_1) then 00851 call aero_state_remove_particle(aero_state, p1, & 00852 id_1_lost, aero_info_1) 00853 end if 00854 else 00855 if (remove_1) then 00856 call aero_state_remove_particle(aero_state, p1, & 00857 id_1_lost, aero_info_1) 00858 end if 00859 if (remove_2) then 00860 call aero_state_remove_particle(aero_state, p2, & 00861 id_2_lost, aero_info_2) 00862 end if 00863 end if 00864 00865 ! add new particle 00866 if (create_new) then 00867 call aero_state_add_particle(aero_state, particle_new, & 00868 allow_resort=.false.) 00869 end if 00870 00871 call aero_info_deallocate(aero_info_1) 00872 call aero_info_deallocate(aero_info_2) 00873 call aero_particle_deallocate(particle_new) 00874 00875 end subroutine coagulate 00876 00877 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00878 00879 end module pmc_coagulation