PartMC 2.1.4
|
00001 ! Copyright (C) 2005-2011 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 #ifdef PMC_USE_MPI 00020 use mpi 00021 #endif 00022 00023 contains 00024 00025 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00026 00027 !> Do coagulation for time del_t. 00028 subroutine mc_coag(coag_kernel_type, bin_grid, env_state, aero_data, & 00029 aero_weight, aero_state, del_t, k_max, tot_n_samp, tot_n_coag) 00030 00031 !> Coagulation kernel type. 00032 integer, intent(in) :: coag_kernel_type 00033 !> Bin grid. 00034 type(bin_grid_t), intent(in) :: bin_grid 00035 !> Environment state. 00036 type(env_state_t), intent(in) :: env_state 00037 !> Aerosol data. 00038 type(aero_data_t), intent(in) :: aero_data 00039 !> Aerosol weight. 00040 type(aero_weight_t), intent(in) :: aero_weight 00041 !> Aerosol state. 00042 type(aero_state_t), intent(inout) :: aero_state 00043 !> Timestep for coagulation. 00044 real(kind=dp) :: del_t 00045 !> Maximum kernel. 00046 real(kind=dp), intent(in) :: k_max(bin_grid%n_bin,bin_grid%n_bin) 00047 !> Total number of samples tested. 00048 integer, intent(out) :: tot_n_samp 00049 !> Number of coagulation events. 00050 integer, intent(out) :: tot_n_coag 00051 00052 logical :: did_coag 00053 integer :: i, j, n_samp, i_samp 00054 real(kind=dp) :: accept_factor 00055 00056 tot_n_samp = 0 00057 tot_n_coag = 0 00058 do i = 1,bin_grid%n_bin 00059 do j = 1,bin_grid%n_bin 00060 call compute_n_samp(aero_state%bin(i)%n_part, & 00061 aero_state%bin(j)%n_part, i == j, k_max(i,j), & 00062 aero_state%comp_vol, del_t, n_samp, accept_factor) 00063 tot_n_samp = tot_n_samp + n_samp 00064 do i_samp = 1,n_samp 00065 ! check we still have enough particles to coagulate 00066 if ((aero_state%bin(i)%n_part < 1) & 00067 .or. (aero_state%bin(j)%n_part < 1) & 00068 .or. ((i == j) .and. (aero_state%bin(i)%n_part < 2))) then 00069 exit 00070 end if 00071 call maybe_coag_pair(bin_grid, env_state, aero_data, & 00072 aero_weight, aero_state, i, j, coag_kernel_type, & 00073 accept_factor, did_coag) 00074 if (did_coag) tot_n_coag = tot_n_coag + 1 00075 enddo 00076 enddo 00077 enddo 00078 00079 end subroutine mc_coag 00080 00081 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00082 00083 !> Compute the number of samples required for the pair of bins. 00084 subroutine compute_n_samp(ni, nj, same_bin, k_max, comp_vol, & 00085 del_t, n_samp, accept_factor) 00086 00087 !> Number particles in first bin. 00088 integer, intent(in) :: ni 00089 !> Number particles in second bin. 00090 integer, intent(in) :: nj 00091 !> Whether first bin is second bin. 00092 logical, intent(in) :: same_bin 00093 !> Maximum kernel value (s^{-1}). 00094 real(kind=dp), intent(in) :: k_max 00095 !> Computational volume (m^3). 00096 real(kind=dp), intent(in) :: comp_vol 00097 !> Timestep (s). 00098 real(kind=dp), intent(in) :: del_t 00099 !> Number of samples per timestep. 00100 integer, intent(out) :: n_samp 00101 !> Scale factor for accept probability (1). 00102 real(kind=dp), intent(out) :: accept_factor 00103 00104 real(kind=dp) :: r_samp, n_samp_mean 00105 real(kind=dp) :: n_possible ! use real(kind=dp) to avoid integer overflow 00106 ! could use integer*8 or integer(kind = 8) 00107 ! or di = selected_int_kind(18), integer(kind=di) 00108 ! to represent 10^{-18} to 10^{18} 00109 00110 if (same_bin) then 00111 ! don't change this to ni * (ni - 1) as the ni/nj distinction 00112 ! is important for coagulation_dist, which also calls this 00113 n_possible = real(ni, kind=dp) * (real(nj, kind=dp) - 1d0) / 2d0 00114 else 00115 n_possible = real(ni, kind=dp) * real(nj, kind=dp) / 2d0 00116 endif 00117 00118 r_samp = k_max / comp_vol * del_t 00119 n_samp_mean = r_samp * n_possible 00120 n_samp = rand_poisson(n_samp_mean) 00121 accept_factor = 1d0 / k_max 00122 00123 ! possible variants: 00124 ! A: accept_factor = 1d0 / k_max 00125 ! B: accept_factor = del_t * n_possible & 00126 ! / (real(n_samp, kind=dp) * comp_vol) 00127 ! timings of test suite as of 2010-12-22T17:12:14-0600: 00128 ! A with n_samp = prob_round(n_samp_mean): 00129 ! 159.82 162.18 156.28 00130 ! A with n_samp = rand_poisson(n_samp_mean): 00131 ! 151.93 158.38 174.74 157.65 00132 ! B with n_samp = ceiling(n_samp_mean): 00133 ! 196.06 200.23 192.41 00134 ! B with n_samp = ceiling(n_samp_mean + 1 * sqrt(n_samp_mean)): 00135 ! 189.78 211.12 195.19 00136 ! B with n_samp = ceiling(n_samp_mean + 3 * sqrt(n_samp_mean)): 00137 ! 214.60 201.25 203.55 00138 00139 end subroutine compute_n_samp 00140 00141 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00142 00143 !> Choose a random pair for potential coagulation and test its 00144 !> probability of coagulation. If it happens, do the coagulation and 00145 !> update all structures. 00146 !! 00147 !! The probability of a coagulation will be taken as <tt>(kernel / 00148 !! k_max)</tt>. 00149 subroutine maybe_coag_pair(bin_grid, env_state, aero_data, aero_weight, & 00150 aero_state, b1, b2, coag_kernel_type, accept_factor, did_coag) 00151 00152 !> Bin grid. 00153 type(bin_grid_t), intent(in) :: bin_grid 00154 !> Environment state. 00155 type(env_state_t), intent(in) :: env_state 00156 !> Aerosol data. 00157 type(aero_data_t), intent(in) :: aero_data 00158 !> Aerosol weight. 00159 type(aero_weight_t), intent(in) :: aero_weight 00160 !> Aerosol state. 00161 type(aero_state_t), intent(inout) :: aero_state 00162 !> Bin of first particle. 00163 integer, intent(in) :: b1 00164 !> Bin of second particle. 00165 integer, intent(in) :: b2 00166 !> Coagulation kernel type. 00167 integer, intent(in) :: coag_kernel_type 00168 !> Scale factor for accept probability (1). 00169 real(kind=dp), intent(in) :: accept_factor 00170 !> Whether a coagulation occured. 00171 logical, intent(out) :: did_coag 00172 00173 integer :: s1, s2 00174 real(kind=dp) :: p, k 00175 00176 did_coag = .false. 00177 00178 call assert(210827476, aero_state%bin(b1)%n_part >= 1) 00179 call assert(368973460, aero_state%bin(b2)%n_part >= 1) 00180 if (b1 == b2) then 00181 call assert(528541565, aero_state%bin(b1)%n_part >= 2) 00182 end if 00183 00184 call find_rand_pair(aero_state, b1, b2, s1, s2) 00185 call weighted_kernel(coag_kernel_type, aero_state%bin(b1)%particle(s1), & 00186 aero_state%bin(b2)%particle(s2), aero_data, aero_weight, & 00187 env_state, k) 00188 p = k * accept_factor 00189 00190 if (pmc_random() .lt. p) then 00191 call coagulate(bin_grid, aero_data, aero_weight, aero_state, & 00192 b1, s1, b2, s2) 00193 did_coag = .true. 00194 end if 00195 00196 end subroutine maybe_coag_pair 00197 00198 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00199 00200 !> Given bins b1 and b2, find a random pair of particles (b1, s1) 00201 !> and (b2, s2) that are not the same particle particle as each 00202 !> other. 00203 subroutine find_rand_pair(aero_state, b1, b2, s1, s2) 00204 00205 !> Aerosol state. 00206 type(aero_state_t), intent(in) :: aero_state 00207 !> Bin number of first particle. 00208 integer, intent(in) :: b1 00209 !> Bin number of second particle. 00210 integer, intent(in) :: b2 00211 !> First rand particle. 00212 integer, intent(out) :: s1 00213 !> Second rand particle. 00214 integer, intent(out) :: s2 00215 00216 ! check we have enough particles to avoid being stuck in an 00217 ! infinite loop below 00218 call assert(362349482, aero_state%bin(b1)%n_part >= 1) 00219 call assert(479121681, aero_state%bin(b2)%n_part >= 1) 00220 if (b1 == b2) then 00221 call assert(161928491, aero_state%bin(b1)%n_part >= 2) 00222 end if 00223 00224 do 00225 s1 = pmc_rand_int(aero_state%bin(b1)%n_part) 00226 s2 = pmc_rand_int(aero_state%bin(b2)%n_part) 00227 if ((b1 /= b2) .or. (s1 /= s2)) then 00228 ! stop generating if we have two distinct particles 00229 exit 00230 end if 00231 end do 00232 00233 end subroutine find_rand_pair 00234 00235 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00236 00237 !> Actually coagulate particle_1 and particle_2 to form particle_new 00238 !> and compute weighting effects, including which particles should 00239 !> be lost and which gained. 00240 subroutine coagulate_weighting(particle_1, particle_2, particle_new, & 00241 aero_data, aero_weight, remove_1, remove_2, create_new, & 00242 id_1_lost, id_2_lost, aero_info_1, aero_info_2) 00243 00244 !> First coagulating aerosol particle. 00245 type(aero_particle_t), intent(in) :: particle_1 00246 !> Second coagulating aerosol particle. 00247 type(aero_particle_t), intent(in) :: particle_2 00248 !> Combined aerosol particle resulting from coagulation of particle_1 00249 !> and particle_2. 00250 type(aero_particle_t), intent(inout) :: particle_new 00251 !> Aerosol data. 00252 type(aero_data_t), intent(in) :: aero_data 00253 !> Aerosol weight. 00254 type(aero_weight_t), intent(in) :: aero_weight 00255 !> Whether to remove particle_1. 00256 logical, intent(out) :: remove_1 00257 !> Whether to remove particle_2. 00258 logical, intent(out) :: remove_2 00259 !> Whether to create particle_new. 00260 logical, intent(out) :: create_new 00261 !> Whether the ID of particle_1 will be lost due to coagulation. 00262 logical, intent(out) :: id_1_lost 00263 !> Whether the ID of particle_2 will be lost due to coagulation. 00264 logical, intent(out) :: id_2_lost 00265 !> The removal information associated with particle_1. 00266 type(aero_info_t), intent(inout) :: aero_info_1 00267 !> The removal information associated with particle_2. 00268 type(aero_info_t), intent(inout) :: aero_info_2 00269 00270 real(kind=dp) :: radius_1, radius_2, radius_new 00271 real(kind=dp) :: weight_1, weight_2, weight_new, weight_min 00272 real(kind=dp) :: prob_remove_1, prob_remove_2, prob_create_new 00273 integer :: info_other_id 00274 00275 call assert(371947172, particle_1%id /= particle_2%id) 00276 00277 ! decide which old particles are to be removed and whether to 00278 ! create the resulting coagulated particle 00279 if (aero_weight%type == AERO_WEIGHT_TYPE_NONE) then 00280 remove_1 = .true. 00281 remove_2 = .true. 00282 create_new = .true. 00283 elseif ((aero_weight%type == AERO_WEIGHT_TYPE_POWER) & 00284 .or. (aero_weight%type == AERO_WEIGHT_TYPE_MFA)) then 00285 radius_1 = aero_particle_radius(particle_1) 00286 radius_2 = aero_particle_radius(particle_2) 00287 radius_new = vol2rad(rad2vol(radius_1) + rad2vol(radius_2)) 00288 weight_1 = aero_weight_value(aero_weight, radius_1) 00289 weight_2 = aero_weight_value(aero_weight, radius_2) 00290 weight_new = aero_weight_value(aero_weight, radius_new) 00291 weight_min = min(weight_1, weight_2, weight_new) 00292 prob_remove_1 = weight_min / weight_1 00293 prob_remove_2 = weight_min / weight_2 00294 prob_create_new = weight_min / weight_new 00295 remove_1 = (pmc_random() < prob_remove_1) 00296 if (aero_weight%type == AERO_WEIGHT_TYPE_MFA) then 00297 remove_2 = .not. remove_1 00298 else 00299 remove_2 = (pmc_random() < prob_remove_2) 00300 end if 00301 create_new = (pmc_random() < prob_create_new) 00302 else 00303 call die_msg(886524113, "unknown aero_weight type: " & 00304 // trim(integer_to_string(aero_weight%type))) 00305 end if 00306 00307 ! figure out what to do about the ID numbers of the various 00308 ! particles --- we try to preserve particle IDs as much as 00309 ! possible 00310 if (create_new) then 00311 id_1_lost = .false. 00312 id_2_lost = .false. 00313 if (remove_1 .and. remove_2) then 00314 if (aero_particle_volume(particle_1) & 00315 > aero_particle_volume(particle_2)) then 00316 id_2_lost = .true. 00317 else 00318 id_1_lost = .true. 00319 end if 00320 end if 00321 else 00322 id_1_lost = remove_1 00323 id_2_lost = remove_2 00324 end if 00325 00326 ! create a new particle and set its ID 00327 if (create_new) then 00328 call aero_particle_deallocate(particle_new) 00329 call aero_particle_allocate_size(particle_new, aero_data%n_spec, & 00330 aero_data%n_source) 00331 call aero_particle_coagulate(particle_1, particle_2, particle_new) 00332 if (remove_1 .and. (.not. id_1_lost)) then 00333 particle_new%id = particle_1%id 00334 call assert(975059559, id_2_lost .eqv. remove_2) 00335 elseif (remove_2 .and. (.not. id_2_lost)) then 00336 particle_new%id = particle_2%id 00337 call assert(246529753, id_1_lost .eqv. remove_1) 00338 else 00339 call aero_particle_new_id(particle_new) 00340 call assert(852038606, id_1_lost .eqv. remove_1) 00341 call assert(254018921, id_2_lost .eqv. remove_2) 00342 end if 00343 info_other_id = particle_new%id 00344 else 00345 info_other_id = 0 00346 end if 00347 00348 aero_info_1%id = particle_1%id 00349 aero_info_1%action = AERO_INFO_COAG 00350 aero_info_1%other_id = info_other_id 00351 00352 aero_info_2%id = particle_2%id 00353 aero_info_2%action = AERO_INFO_COAG 00354 aero_info_2%other_id = info_other_id 00355 00356 end subroutine coagulate_weighting 00357 00358 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00359 00360 !> Join together particles (b1, s1) and (b2, s2), updating all 00361 !> particle and bin structures to reflect the change. 00362 subroutine coagulate(bin_grid, aero_data, aero_weight, aero_state, & 00363 b1, s1, b2, s2) 00364 00365 !> Bin grid. 00366 type(bin_grid_t), intent(in) :: bin_grid 00367 !> Aerosol data. 00368 type(aero_data_t), intent(in) :: aero_data 00369 !> Aerosol weight. 00370 type(aero_weight_t), intent(in) :: aero_weight 00371 !> Aerosol state. 00372 type(aero_state_t), intent(inout) :: aero_state 00373 !> First particle (bin number). 00374 integer, intent(in) :: b1 00375 !> First particle (number in bin). 00376 integer, intent(in) :: s1 00377 !> Second particle (bin number). 00378 integer, intent(in) :: b2 00379 !> Second particle (number in bin). 00380 integer, intent(in) :: s2 00381 00382 type(aero_particle_t), pointer :: particle_1, particle_2 00383 type(aero_particle_t) :: particle_new 00384 integer :: bn 00385 type(aero_info_t) :: aero_info_1, aero_info_2 00386 logical :: remove_1, remove_2, create_new, id_1_lost, id_2_lost 00387 00388 call aero_particle_allocate(particle_new) 00389 call aero_info_allocate(aero_info_1) 00390 call aero_info_allocate(aero_info_2) 00391 00392 particle_1 => aero_state%bin(b1)%particle(s1) 00393 particle_2 => aero_state%bin(b2)%particle(s2) 00394 00395 call coagulate_weighting(particle_1, particle_2, particle_new, & 00396 aero_data, aero_weight, remove_1, remove_2, create_new, & 00397 id_1_lost, id_2_lost, aero_info_1, aero_info_2) 00398 00399 ! remove old particles 00400 if ((b1 == b2) .and. (s2 > s1)) then 00401 ! handle a tricky corner case where we have to watch for s2 or 00402 ! s1 being the last entry in the array and being repacked when 00403 ! the other one is removed 00404 if (remove_2) then 00405 call aero_state_remove_particle(aero_state, b2, s2, & 00406 id_2_lost, aero_info_2) 00407 end if 00408 if (remove_1) then 00409 call aero_state_remove_particle(aero_state, b1, s1, & 00410 id_1_lost, aero_info_1) 00411 end if 00412 else 00413 if (remove_1) then 00414 call aero_state_remove_particle(aero_state, b1, s1, & 00415 id_1_lost, aero_info_1) 00416 end if 00417 if (remove_2) then 00418 call aero_state_remove_particle(aero_state, b2, s2, & 00419 id_2_lost, aero_info_2) 00420 end if 00421 end if 00422 00423 ! add new particle 00424 if (create_new) then 00425 bn = aero_particle_in_bin(particle_new, bin_grid) 00426 call aero_state_add_particle(aero_state, bn, particle_new) 00427 end if 00428 00429 call aero_info_deallocate(aero_info_1) 00430 call aero_info_deallocate(aero_info_2) 00431 call aero_particle_deallocate(particle_new) 00432 00433 end subroutine coagulate 00434 00435 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00436 00437 end module pmc_coagulation