PartMC
2.2.1
|
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_coag_kernel module. 00007 00008 !> Generic coagulation kernel. 00009 module pmc_coag_kernel 00010 00011 use pmc_env_state 00012 use pmc_bin_grid 00013 use pmc_aero_particle 00014 use pmc_aero_data 00015 use pmc_aero_weight 00016 use pmc_coag_kernel_sedi 00017 use pmc_coag_kernel_additive 00018 use pmc_coag_kernel_constant 00019 use pmc_coag_kernel_brown 00020 use pmc_coag_kernel_zero 00021 00022 !> Maximum length of a mode type. 00023 integer, parameter :: COAG_KERNEL_TYPE_LEN = 20 00024 00025 !> Type code for an undefined or invalid kernel. 00026 integer, parameter :: COAG_KERNEL_TYPE_INVALID = 0 00027 !> Type code for a sedimentation kernel. 00028 integer, parameter :: COAG_KERNEL_TYPE_SEDI = 1 00029 !> Type code for an additive kernel. 00030 integer, parameter :: COAG_KERNEL_TYPE_ADDITIVE = 2 00031 !> Type code for a constant kernel. 00032 integer, parameter :: COAG_KERNEL_TYPE_CONSTANT = 3 00033 !> Type code for a Brownian kernel. 00034 integer, parameter :: COAG_KERNEL_TYPE_BROWN = 4 00035 !> Type code for a zero kernel. 00036 integer, parameter :: COAG_KERNEL_TYPE_ZERO = 5 00037 00038 contains 00039 00040 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00041 00042 !> Return a string representation of a kernel type. 00043 character(len=COAG_KERNEL_TYPE_LEN) function coag_kernel_type_to_string( & 00044 coag_kernel_type) 00045 00046 !> Coagulation kernel type. 00047 integer, intent(in) :: coag_kernel_type 00048 00049 if (coag_kernel_type == COAG_KERNEL_TYPE_INVALID) then 00050 coag_kernel_type_to_string = "invalid" 00051 elseif (coag_kernel_type == COAG_KERNEL_TYPE_SEDI) then 00052 coag_kernel_type_to_string = "sedi" 00053 elseif (coag_kernel_type == COAG_KERNEL_TYPE_ADDITIVE) then 00054 coag_kernel_type_to_string = "additive" 00055 elseif (coag_kernel_type == COAG_KERNEL_TYPE_CONSTANT) then 00056 coag_kernel_type_to_string = "constant" 00057 elseif (coag_kernel_type == COAG_KERNEL_TYPE_BROWN) then 00058 coag_kernel_type_to_string = "brown" 00059 elseif (coag_kernel_type == COAG_KERNEL_TYPE_ZERO) then 00060 coag_kernel_type_to_string = "zero" 00061 else 00062 coag_kernel_type_to_string = "unknown" 00063 end if 00064 00065 end function coag_kernel_type_to_string 00066 00067 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00068 00069 !> Evalulate a coagulation kernel function. 00070 subroutine kernel(coag_kernel_type, aero_particle_1, aero_particle_2, & 00071 aero_data, env_state, k) 00072 00073 !> Coagulation kernel type. 00074 integer, intent(in) :: coag_kernel_type 00075 !> First particle. 00076 type(aero_particle_t), intent(in) :: aero_particle_1 00077 !> Second particle. 00078 type(aero_particle_t), intent(in) :: aero_particle_2 00079 !> Aerosol data. 00080 type(aero_data_t), intent(in) :: aero_data 00081 !> Environment state. 00082 type(env_state_t), intent(in) :: env_state 00083 !> Kernel k(a,b) (m^3/s). 00084 real(kind=dp), intent(out) :: k 00085 00086 if (coag_kernel_type == COAG_KERNEL_TYPE_SEDI) then 00087 call kernel_sedi(aero_particle_1, aero_particle_2, & 00088 aero_data, env_state, k) 00089 elseif (coag_kernel_type == COAG_KERNEL_TYPE_ADDITIVE) then 00090 call kernel_additive(aero_particle_1, aero_particle_2, & 00091 aero_data, env_state, k) 00092 elseif (coag_kernel_type == COAG_KERNEL_TYPE_CONSTANT) then 00093 call kernel_constant(aero_particle_1, aero_particle_2, & 00094 aero_data, env_state, k) 00095 elseif (coag_kernel_type == COAG_KERNEL_TYPE_BROWN) then 00096 call kernel_brown(aero_particle_1, aero_particle_2, & 00097 aero_data, env_state, k) 00098 elseif (coag_kernel_type == COAG_KERNEL_TYPE_ZERO) then 00099 call kernel_zero(aero_particle_1, aero_particle_2, & 00100 aero_data, env_state, k) 00101 else 00102 call die_msg(200724934, "Unknown kernel type: " & 00103 // trim(integer_to_string(coag_kernel_type))) 00104 end if 00105 00106 end subroutine kernel 00107 00108 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00109 00110 !> Compute the minimum and maximum coagulation kernel. 00111 subroutine kernel_minmax(coag_kernel_type, v1, v2, aero_data, env_state, & 00112 k_min, k_max) 00113 00114 !> Coagulation kernel type. 00115 integer, intent(in) :: coag_kernel_type 00116 !> Volume of first particle (m^3). 00117 real(kind=dp), intent(in) :: v1 00118 !> Volume of second particle (m^3). 00119 real(kind=dp), intent(in) :: v2 00120 !> Aerosol data. 00121 type(aero_data_t), intent(in) :: aero_data 00122 !> Environment state. 00123 type(env_state_t), intent(in) :: env_state 00124 !> Minimum kernel value (m^3/s). 00125 real(kind=dp), intent(out) :: k_min 00126 !> Maximum kernel value (m^3/s). 00127 real(kind=dp), intent(out) :: k_max 00128 00129 if (coag_kernel_type == COAG_KERNEL_TYPE_SEDI) then 00130 call kernel_sedi_minmax(v1, v2, aero_data, env_state, k_min, k_max) 00131 elseif (coag_kernel_type == COAG_KERNEL_TYPE_ADDITIVE) then 00132 call kernel_additive_minmax(v1, v2, aero_data, env_state, k_min, k_max) 00133 elseif (coag_kernel_type == COAG_KERNEL_TYPE_CONSTANT) then 00134 call kernel_constant_minmax(v1, v2, aero_data, env_state, k_min, k_max) 00135 elseif (coag_kernel_type == COAG_KERNEL_TYPE_BROWN) then 00136 call kernel_brown_minmax(v1, v2, aero_data, env_state, k_min, k_max) 00137 elseif (coag_kernel_type == COAG_KERNEL_TYPE_ZERO) then 00138 call kernel_zero_minmax(v1, v2, aero_data, env_state, k_min, k_max) 00139 else 00140 call die_msg(330498208, "Unknown kernel type: " & 00141 // trim(integer_to_string(coag_kernel_type))) 00142 end if 00143 00144 end subroutine kernel_minmax 00145 00146 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00147 00148 !> Compute the kernel value with the given weight. 00149 subroutine weighted_kernel(coag_kernel_type, aero_particle_1, & 00150 aero_particle_2, aero_data, aero_weight, env_state, k) 00151 00152 !> Coagulation kernel type. 00153 integer, intent(in) :: coag_kernel_type 00154 !> First particle. 00155 type(aero_particle_t), intent(in) :: aero_particle_1 00156 !> Second particle. 00157 type(aero_particle_t), intent(in) :: aero_particle_2 00158 !> Aerosol data. 00159 type(aero_data_t), intent(in) :: aero_data 00160 !> Aerosol weight. 00161 type(aero_weight_t), intent(in) :: aero_weight 00162 !> Environment state. 00163 type(env_state_t), intent(in) :: env_state 00164 !> Coagulation kernel. 00165 real(kind=dp), intent(out) :: k 00166 00167 real(kind=dp) :: unweighted_k 00168 real(kind=dp) :: radius_1, radius_2, radius_1_plus_2 00169 real(kind=dp) :: weight_1, weight_2, weight_1_plus_2, weight_min 00170 00171 call kernel(coag_kernel_type, aero_particle_1, aero_particle_2, & 00172 aero_data, env_state, unweighted_k) 00173 radius_1 = aero_particle_radius(aero_particle_1) 00174 radius_2 = aero_particle_radius(aero_particle_2) 00175 radius_1_plus_2 = vol2rad(rad2vol(radius_1) + rad2vol(radius_2)) 00176 weight_1 = aero_weight_num_conc_at_radius(aero_weight, radius_1) & 00177 * aero_weight%comp_vol 00178 weight_2 = aero_weight_num_conc_at_radius(aero_weight, radius_2) & 00179 * aero_weight%comp_vol 00180 weight_1_plus_2 = aero_weight_num_conc_at_radius(aero_weight, & 00181 radius_1_plus_2) * aero_weight%comp_vol 00182 weight_min = min(weight_1, weight_2, weight_1_plus_2) 00183 k = unweighted_k * weight_1 * weight_2 / weight_min 00184 00185 end subroutine weighted_kernel 00186 00187 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00188 00189 !> Compute the kernel value with the given number concentration 00190 !> weighting. 00191 subroutine num_conc_weighted_kernel(coag_kernel_type, aero_particle_1, & 00192 aero_particle_2, aero_data, aero_weight_array, env_state, k) 00193 00194 !> Coagulation kernel type. 00195 integer, intent(in) :: coag_kernel_type 00196 !> First particle. 00197 type(aero_particle_t), intent(in) :: aero_particle_1 00198 !> Second particle. 00199 type(aero_particle_t), intent(in) :: aero_particle_2 00200 !> Aerosol data. 00201 type(aero_data_t), intent(in) :: aero_data 00202 !> Aerosol weight array. 00203 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00204 !> Environment state. 00205 type(env_state_t), intent(in) :: env_state 00206 !> Coagulation kernel. 00207 real(kind=dp), intent(out) :: k 00208 00209 real(kind=dp) :: unweighted_k, radius_1, radius_2 00210 00211 call kernel(coag_kernel_type, aero_particle_1, aero_particle_2, & 00212 aero_data, env_state, unweighted_k) 00213 radius_1 = aero_particle_radius(aero_particle_1) 00214 radius_2 = aero_particle_radius(aero_particle_2) 00215 k = unweighted_k * coag_num_conc_factor(aero_weight_array, & 00216 radius_1, radius_2) 00217 00218 end subroutine num_conc_weighted_kernel 00219 00220 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00221 00222 !> Compute the minimum and maximum kernel value with the given weight. 00223 subroutine weighted_kernel_minmax(coag_kernel_type, v1, v2, aero_data, & 00224 aero_weight, env_state, k_min, k_max) 00225 00226 !> Coagulation kernel type. 00227 integer, intent(in) :: coag_kernel_type 00228 !> Volume of first particle. 00229 real(kind=dp), intent(in) :: v1 00230 !> Volume of second particle. 00231 real(kind=dp), intent(in) :: v2 00232 !> Aerosol data. 00233 type(aero_data_t), intent(in) :: aero_data 00234 !> Aerosol weight. 00235 type(aero_weight_t), intent(in) :: aero_weight 00236 !> Environment state. 00237 type(env_state_t), intent(in) :: env_state 00238 !> Coagulation kernel minimum value. 00239 real(kind=dp), intent(out) :: k_min 00240 !> Coagulation kernel maximum value. 00241 real(kind=dp), intent(out) :: k_max 00242 00243 real(kind=dp) :: unweighted_k_min, unweighted_k_max 00244 real(kind=dp) :: weight_1, weight_2, weight_1_plus_2, weight_min 00245 00246 call kernel_minmax(coag_kernel_type, v1, v2, aero_data, env_state, & 00247 unweighted_k_min, unweighted_k_max) 00248 00249 weight_1 = aero_weight_num_conc_at_radius(aero_weight, vol2rad(v1)) & 00250 * aero_weight%comp_vol 00251 weight_2 = aero_weight_num_conc_at_radius(aero_weight, vol2rad(v2)) & 00252 * aero_weight%comp_vol 00253 weight_1_plus_2 = aero_weight_num_conc_at_radius(aero_weight, & 00254 vol2rad(v1 + v2)) * aero_weight%comp_vol 00255 weight_min = min(weight_1, weight_2, weight_1_plus_2) 00256 k_min = unweighted_k_min * weight_1 * weight_2 / weight_min 00257 k_max = unweighted_k_max * weight_1 * weight_2 / weight_min 00258 00259 end subroutine weighted_kernel_minmax 00260 00261 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00262 00263 !> Computes an array of kernel values for each bin pair. k(i,j) is 00264 !> the kernel value at the centers of bins i and j. This assumes the 00265 !> kernel is only a function of the particle volumes. 00266 subroutine bin_kernel(n_bin, bin_r, aero_data, coag_kernel_type, & 00267 env_state, k) 00268 00269 !> Number of bins. 00270 integer, intent(in) :: n_bin 00271 !> Radii of particles in bins (m). 00272 real(kind=dp), intent(in) :: bin_r(n_bin) 00273 !> Aerosol data. 00274 type(aero_data_t), intent(in) :: aero_data 00275 !> Coagulation kernel type. 00276 integer, intent(in) :: coag_kernel_type 00277 !> Environment state. 00278 type(env_state_t), intent(in) :: env_state 00279 !> Kernel values. 00280 real(kind=dp), intent(out) :: k(n_bin,n_bin) 00281 00282 integer :: i, j 00283 type(aero_particle_t) :: aero_particle_1, aero_particle_2 00284 00285 call aero_particle_allocate_size(aero_particle_1, aero_data%n_spec, & 00286 aero_data%n_source) 00287 call aero_particle_allocate_size(aero_particle_2, aero_data%n_spec, & 00288 aero_data%n_source) 00289 do i = 1,n_bin 00290 do j = 1,n_bin 00291 aero_particle_1%vol(1) = rad2vol(bin_r(i)) 00292 aero_particle_2%vol(1) = rad2vol(bin_r(j)) 00293 call kernel(coag_kernel_type, aero_particle_1, aero_particle_2, & 00294 aero_data, env_state, k(i,j)) 00295 end do 00296 end do 00297 call aero_particle_deallocate(aero_particle_1) 00298 call aero_particle_deallocate(aero_particle_2) 00299 00300 end subroutine bin_kernel 00301 00302 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00303 00304 !> Estimate an array of minimum and maximum kernel values. Given 00305 !> particles v1 in bin b1 and v2 in bin b2, it is probably true that 00306 !> <tt>k_min(b1,b2) <= kernel(v1,v2) <= k_max(b1,b2)</tt>. 00307 subroutine est_k_minmax_binned(bin_grid, coag_kernel_type, aero_data, & 00308 aero_weight, env_state, k_min, k_max) 00309 00310 !> Bin_grid. 00311 type(bin_grid_t), intent(in) :: bin_grid 00312 !> Coagulation kernel type. 00313 integer, intent(in) :: coag_kernel_type 00314 !> Aerosol data. 00315 type(aero_data_t), intent(in) :: aero_data 00316 !> Aerosol weight. 00317 type(aero_weight_t), intent(in) :: aero_weight 00318 !> Environment state. 00319 type(env_state_t), intent(in) :: env_state 00320 !> Minimum kernel vals. 00321 real(kind=dp), intent(out) :: k_min(bin_grid%n_bin,bin_grid%n_bin) 00322 !> Maximum kernel vals. 00323 real(kind=dp), intent(out) :: k_max(bin_grid%n_bin,bin_grid%n_bin) 00324 00325 integer i, j 00326 00327 do i = 1,bin_grid%n_bin 00328 do j = 1,bin_grid%n_bin 00329 call est_k_minmax_for_bin(bin_grid, coag_kernel_type, i, j, & 00330 aero_data, aero_weight, env_state, k_min(i,j), k_max(i,j)) 00331 end do 00332 end do 00333 00334 end subroutine est_k_minmax_binned 00335 00336 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00337 00338 !> Estimate an array of minimum and maximum kernel values. Given 00339 !> particles v1 in bin b1 and v2 in bin b2, it is probably true that 00340 !> <tt>k_min(b1,b2) <= kernel(v1,v2) <= k_max(b1,b2)</tt>. 00341 subroutine est_k_minmax_binned_unweighted(bin_grid, coag_kernel_type, & 00342 aero_data, env_state, k_min, k_max) 00343 00344 !> Bin_grid. 00345 type(bin_grid_t), intent(in) :: bin_grid 00346 !> Coagulation kernel type. 00347 integer, intent(in) :: coag_kernel_type 00348 !> Aerosol data. 00349 type(aero_data_t), intent(in) :: aero_data 00350 !> Environment state. 00351 type(env_state_t), intent(in) :: env_state 00352 !> Minimum kernel vals. 00353 real(kind=dp), intent(out) :: k_min(bin_grid%n_bin,bin_grid%n_bin) 00354 !> Maximum kernel vals. 00355 real(kind=dp), intent(out) :: k_max(bin_grid%n_bin,bin_grid%n_bin) 00356 00357 integer i, j 00358 00359 do i = 1,bin_grid%n_bin 00360 do j = 1,bin_grid%n_bin 00361 call est_k_minmax_for_bin_unweighted(bin_grid, coag_kernel_type, & 00362 i, j, aero_data, env_state, k_min(i,j), k_max(i,j)) 00363 end do 00364 end do 00365 00366 end subroutine est_k_minmax_binned_unweighted 00367 00368 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00369 00370 !> Samples within bins b1 and b2 to find the minimum and maximum 00371 !> value of the kernel between particles from the two bins. 00372 subroutine est_k_minmax_for_bin(bin_grid, coag_kernel_type, b1, b2, & 00373 aero_data, aero_weight, env_state, k_min, k_max) 00374 00375 !> Bin_grid. 00376 type(bin_grid_t), intent(in) :: bin_grid 00377 !> Coagulation kernel type. 00378 integer, intent(in) :: coag_kernel_type 00379 !> First bin. 00380 integer, intent(in) :: b1 00381 !> Second bin. 00382 integer, intent(in) :: b2 00383 !> Aerosol data. 00384 type(aero_data_t), intent(in) :: aero_data 00385 !> Aerosol weight. 00386 type(aero_weight_t), intent(in) :: aero_weight 00387 !> Environment state. 00388 type(env_state_t), intent(in) :: env_state 00389 !> Minimum kernel value. 00390 real(kind=dp), intent(out) :: k_min 00391 !> Maximum kernel value. 00392 real(kind=dp), intent(out) :: k_max 00393 00394 !> Number of sample points per bin. 00395 integer, parameter :: n_sample = 3 00396 !> Over-estimation scale factor parameter. 00397 real(kind=dp), parameter :: over_scale = 1.5d0 00398 00399 real(kind=dp) :: v1, v2, v1_high, v1_low, v2_high, v2_low 00400 real(kind=dp) :: new_k_min, new_k_max 00401 integer :: i, j 00402 00403 ! v1_low < bin_v(b1) < v1_high 00404 v1_low = rad2vol(bin_grid%edge_radius(b1)) 00405 v1_high = rad2vol(bin_grid%edge_radius(b1 + 1)) 00406 00407 ! v2_low < bin_v(b2) < v2_high 00408 v2_low = rad2vol(bin_grid%edge_radius(b2)) 00409 v2_high = rad2vol(bin_grid%edge_radius(b2 + 1)) 00410 00411 do i = 1,n_sample 00412 do j = 1,n_sample 00413 v1 = interp_linear_disc(v1_low, v1_high, n_sample, i) 00414 v2 = interp_linear_disc(v2_low, v2_high, n_sample, j) 00415 call weighted_kernel_minmax(coag_kernel_type, v1, v2, aero_data, & 00416 aero_weight, env_state, new_k_min, new_k_max) 00417 if ((i == 1) .and. (j == 1)) then 00418 k_min = new_k_min 00419 k_max = new_k_max 00420 else 00421 k_min = min(k_min, new_k_min) 00422 k_max = max(k_max, new_k_max) 00423 end if 00424 end do 00425 end do 00426 00427 k_max = k_max * over_scale 00428 00429 end subroutine est_k_minmax_for_bin 00430 00431 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00432 00433 !> Samples within bins b1 and b2 to find the minimum and maximum 00434 !> value of the kernel between particles from the two bins. 00435 subroutine est_k_minmax_for_bin_unweighted(bin_grid, coag_kernel_type, & 00436 b1, b2, aero_data, env_state, k_min, k_max) 00437 00438 !> Bin_grid. 00439 type(bin_grid_t), intent(in) :: bin_grid 00440 !> Coagulation kernel type. 00441 integer, intent(in) :: coag_kernel_type 00442 !> First bin. 00443 integer, intent(in) :: b1 00444 !> Second bin. 00445 integer, intent(in) :: b2 00446 !> Aerosol data. 00447 type(aero_data_t), intent(in) :: aero_data 00448 !> Environment state. 00449 type(env_state_t), intent(in) :: env_state 00450 !> Minimum kernel value. 00451 real(kind=dp), intent(out) :: k_min 00452 !> Maximum kernel value. 00453 real(kind=dp), intent(out) :: k_max 00454 00455 !> Number of sample points per bin. 00456 integer, parameter :: n_sample = 3 00457 !> Over-estimation scale factor parameter. 00458 real(kind=dp), parameter :: over_scale = 2d0 00459 00460 real(kind=dp) :: v1, v2, v1_high, v1_low, v2_high, v2_low 00461 real(kind=dp) :: new_k_min, new_k_max 00462 integer :: i, j 00463 00464 ! v1_low < bin_v(b1) < v1_high 00465 v1_low = rad2vol(bin_grid%edge_radius(b1)) 00466 v1_high = rad2vol(bin_grid%edge_radius(b1 + 1)) 00467 00468 ! v2_low < bin_v(b2) < v2_high 00469 v2_low = rad2vol(bin_grid%edge_radius(b2)) 00470 v2_high = rad2vol(bin_grid%edge_radius(b2 + 1)) 00471 00472 do i = 1,n_sample 00473 do j = 1,n_sample 00474 v1 = interp_linear_disc(v1_low, v1_high, n_sample, i) 00475 v2 = interp_linear_disc(v2_low, v2_high, n_sample, j) 00476 call kernel_minmax(coag_kernel_type, v1, v2, aero_data, & 00477 env_state, new_k_min, new_k_max) 00478 if ((i == 1) .and. (j == 1)) then 00479 k_min = new_k_min 00480 k_max = new_k_max 00481 else 00482 k_min = min(k_min, new_k_min) 00483 k_max = max(k_max, new_k_max) 00484 end if 00485 end do 00486 end do 00487 00488 k_max = k_max * over_scale 00489 00490 end subroutine est_k_minmax_for_bin_unweighted 00491 00492 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00493 00494 !> Coagulation scale factor due to number concentrations. 00495 real(kind=dp) function coag_num_conc_factor(aero_weight_array, & 00496 r_1, r_2) 00497 00498 !> Aerosol weight array. 00499 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00500 !> Radius of first particle. 00501 real(kind=dp), intent(in) :: r_1 00502 !> Radius of second particle. 00503 real(kind=dp), intent(in) :: r_2 00504 00505 real(kind=dp) :: r_12, nc_1, nc_2, nc_12, nc_min 00506 00507 r_12 = vol2rad(rad2vol(r_1) + rad2vol(r_2)) 00508 nc_1 = aero_weight_array_num_conc_at_radius(aero_weight_array, r_1) 00509 nc_2 = aero_weight_array_num_conc_at_radius(aero_weight_array, r_2) 00510 nc_12 = aero_weight_array_num_conc_at_radius(aero_weight_array, r_12) 00511 nc_min = min(nc_1, nc_2, nc_12) 00512 coag_num_conc_factor = nc_1 * nc_2 / nc_min 00513 00514 end function coag_num_conc_factor 00515 00516 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00517 00518 !> Determine the minimum and maximum number concentration factors 00519 !> for coagulation. 00520 subroutine max_coag_num_conc_factor(aero_weight_array, bin_grid, i_bin, & 00521 j_bin, f_max) 00522 00523 !> Aerosol weight array. 00524 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00525 !> Bin grid. 00526 type(bin_grid_t), intent(in) :: bin_grid 00527 !> First bin number. 00528 integer, intent(in) :: i_bin 00529 !> Second bin number. 00530 integer, intent(in) :: j_bin 00531 !> Maximum coagulation factor. 00532 real(kind=dp), intent(out) :: f_max 00533 00534 real(kind=dp) :: i_r_min, i_r_max, j_r_min, j_r_max, ij_r_min, ij_r_max 00535 real(kind=dp) :: nc_i_max, nc_j_max, nc_i_min, nc_j_min, nc_min 00536 logical :: monotone_increasing, monotone_decreasing 00537 00538 call aero_weight_array_check_monotonicity(aero_weight_array, & 00539 monotone_increasing, monotone_decreasing) 00540 call assert(121527417, monotone_increasing .or. monotone_decreasing) 00541 00542 i_r_min = bin_grid%edge_radius(i_bin) 00543 i_r_max = bin_grid%edge_radius(i_bin + 1) 00544 j_r_min = bin_grid%edge_radius(j_bin) 00545 j_r_max = bin_grid%edge_radius(j_bin + 1) 00546 ij_r_min = i_r_min + j_r_min 00547 ij_r_max = i_r_max + j_r_max 00548 00549 if (monotone_increasing) then 00550 nc_i_max = aero_weight_array_num_conc_at_radius(aero_weight_array, & 00551 i_r_max) 00552 nc_j_max = aero_weight_array_num_conc_at_radius(aero_weight_array, & 00553 i_r_max) 00554 nc_i_min = aero_weight_array_num_conc_at_radius(aero_weight_array, & 00555 i_r_min) 00556 nc_j_min = aero_weight_array_num_conc_at_radius(aero_weight_array, & 00557 i_r_min) 00558 nc_min = min(nc_i_min, nc_j_min) 00559 f_max = nc_i_max * nc_j_max / nc_min 00560 else 00561 call assert(990892385, monotone_decreasing) 00562 nc_i_max = aero_weight_array_num_conc_at_radius(aero_weight_array, & 00563 i_r_min) 00564 nc_j_max = aero_weight_array_num_conc_at_radius(aero_weight_array, & 00565 i_r_min) 00566 nc_min = aero_weight_array_num_conc_at_radius(aero_weight_array, & 00567 ij_r_max) 00568 f_max = nc_i_max * nc_j_max / nc_min 00569 end if 00570 00571 end subroutine max_coag_num_conc_factor 00572 00573 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00574 00575 !> Determine the minimum and maximum number concentration factors 00576 !> for coagulation. 00577 subroutine max_coag_num_conc_factor_better(aero_weight_array, bin_grid, i_bin, & 00578 j_bin, f_max) 00579 00580 !> Aerosol weight array. 00581 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00582 !> Bin grid. 00583 type(bin_grid_t), intent(in) :: bin_grid 00584 !> First bin number. 00585 integer, intent(in) :: i_bin 00586 !> Second bin number. 00587 integer, intent(in) :: j_bin 00588 !> Maximum coagulation factor. 00589 real(kind=dp), intent(out) :: f_max 00590 00591 integer, parameter :: n_sample = 5 00592 00593 real(kind=dp) :: i_r_min, i_r_max, j_r_min, j_r_max, i_r, j_r, f 00594 integer :: i_sample, j_sample 00595 00596 i_r_min = bin_grid%edge_radius(i_bin) 00597 i_r_max = bin_grid%edge_radius(i_bin + 1) 00598 j_r_min = bin_grid%edge_radius(j_bin) 00599 j_r_max = bin_grid%edge_radius(j_bin + 1) 00600 00601 f_max = 0d0 00602 do i_sample = 1,n_sample 00603 do j_sample = 1,n_sample 00604 i_r = interp_linear_disc(i_r_min, i_r_max, n_sample, i_sample) 00605 j_r = interp_linear_disc(j_r_min, j_r_max, n_sample, j_sample) 00606 f = coag_num_conc_factor(aero_weight_array, i_r, j_r) 00607 f_max = max(f_max, f) 00608 end do 00609 end do 00610 00611 end subroutine max_coag_num_conc_factor_better 00612 00613 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00614 00615 !> Read the specification for a kernel type from a spec file and 00616 !> generate it. 00617 subroutine spec_file_read_coag_kernel_type(file, coag_kernel_type) 00618 00619 !> Spec file. 00620 type(spec_file_t), intent(inout) :: file 00621 !> Kernel type. 00622 integer, intent(out) :: coag_kernel_type 00623 00624 character(len=SPEC_LINE_MAX_VAR_LEN) :: kernel_name 00625 00626 !> \page input_format_coag_kernel Input File Format: Coagulation Kernel 00627 !! 00628 !! The coagulation kernel is specified by the parameter: 00629 !! - \b coag_kernel (string): the type of coagulation kernel --- 00630 !! must be one of: \c sedi for the gravitational sedimentation 00631 !! kernel; \c additive for the additive kernel; \c constant 00632 !! for the constant kernel; \c brown for the Brownian kernel, 00633 !! or \c zero for no coagulation 00634 !! 00635 !! See also: 00636 !! - \ref spec_file_format --- the input file text format 00637 00638 call spec_file_read_string(file, 'coag_kernel', kernel_name) 00639 if (trim(kernel_name) == 'sedi') then 00640 coag_kernel_type = COAG_KERNEL_TYPE_SEDI 00641 elseif (trim(kernel_name) == 'additive') then 00642 coag_kernel_type = COAG_KERNEL_TYPE_ADDITIVE 00643 elseif (trim(kernel_name) == 'constant') then 00644 coag_kernel_type = COAG_KERNEL_TYPE_CONSTANT 00645 elseif (trim(kernel_name) == 'brown') then 00646 coag_kernel_type = COAG_KERNEL_TYPE_BROWN 00647 elseif (trim(kernel_name) == 'zero') then 00648 coag_kernel_type = COAG_KERNEL_TYPE_ZERO 00649 else 00650 call spec_file_die_msg(920761229, file, & 00651 "Unknown coagulation kernel type: " // trim(kernel_name)) 00652 end if 00653 00654 end subroutine spec_file_read_coag_kernel_type 00655 00656 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00657 00658 end module pmc_coag_kernel