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_additive module. 00007 00008 !> Additive coagulation kernel. 00009 module pmc_coag_kernel_additive 00010 00011 use pmc_bin_grid 00012 use pmc_env_state 00013 use pmc_util 00014 use pmc_constants 00015 use pmc_constants 00016 use pmc_aero_binned 00017 use pmc_aero_data 00018 use pmc_aero_dist 00019 use pmc_aero_data 00020 use pmc_aero_particle 00021 00022 !> Scaling coefficient for constant kernel. 00023 real(kind=dp), parameter :: beta_1 = 1000d0 00024 00025 contains 00026 00027 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00028 00029 !> Additive coagulation kernel. 00030 subroutine kernel_additive(aero_particle_1, aero_particle_2, & 00031 aero_data, env_state, k) 00032 00033 !> First particle. 00034 type(aero_particle_t), intent(in) :: aero_particle_1 00035 !> Second particle. 00036 type(aero_particle_t), intent(in) :: aero_particle_2 00037 !> Aerosol data. 00038 type(aero_data_t), intent(in) :: aero_data 00039 !> Environment state. 00040 type(env_state_t), intent(in) :: env_state 00041 !> Coagulation kernel. 00042 real(kind=dp), intent(out) :: k 00043 00044 k = beta_1 * (aero_particle_volume(aero_particle_1) & 00045 + aero_particle_volume(aero_particle_2)) 00046 00047 end subroutine kernel_additive 00048 00049 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00050 00051 !> Minimum and maximum values of the additive kernel. 00052 subroutine kernel_additive_minmax(v1, v2, aero_data, env_state, k_min, k_max) 00053 00054 !> Volume of first particle. 00055 real(kind=dp), intent(in) :: v1 00056 !> Volume of second particle. 00057 real(kind=dp), intent(in) :: v2 00058 !> Aerosol data. 00059 type(aero_data_t), intent(in) :: aero_data 00060 !> Environment state. 00061 type(env_state_t), intent(in) :: env_state 00062 !> Coagulation kernel minimum value. 00063 real(kind=dp), intent(out) :: k_min 00064 !> Coagulation kernel maximum value. 00065 real(kind=dp), intent(out) :: k_max 00066 00067 k_min = beta_1 * (v1 + v2) 00068 k_max = k_min 00069 00070 end subroutine kernel_additive_minmax 00071 00072 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00073 00074 !> Exact solution with the additive coagulation kernel and 00075 !> exponential initial condition. 00076 !! 00077 !! Given input paramaters \f$R\f$ and \f$N_0\f$ we let the mean 00078 !! volume be \f$v_\mu = \frac{4\pi}{3} R^3\f$ and define the 00079 !! rescaled times \f$\tau = N_0 v_\mu \beta_1 t\f$ and \f$T = 1 - 00080 !! e^{-\tau}\f$, where \f$\beta_1\f$ is the fixed kernel scaling 00081 !! parameter. Then the solution is 00082 !! \f[ 00083 !! n(D,t) \ {\rm d}\ln D 00084 !! = \frac{\pi}{2} D^3 00085 !! \frac{N_0}{v} \frac{1 - T}{\sqrt{T}} 00086 !! \exp\left(-(1 + T) \frac{v}{v_\mu}\right) 00087 !! I_1\left(2 \frac{v}{v_\mu} \sqrt{T}\right) {\rm d}\ln D 00088 !! \f] 00089 !! where \f$I_1(x)\f$ is the <a 00090 !! href="http://en.wikipedia.org/wiki/Bessel_function">modified 00091 !! Bessel function of the first kind</a> and \f$v = \frac{\pi}{6} 00092 !! D^3\f$. 00093 !! 00094 !! For small \f$x\f$ we have \f$I_1(x) \approx \frac{x}{2}\f$, so 00095 !! this solution has initial condition 00096 !! \f[ 00097 !! n(D,t) \ {\rm d}\ln D 00098 !! = \frac{\pi}{2} D^3 \frac{N_0}{v_\mu} 00099 !! \exp\left(-\frac{v}{v_\mu}\right) {\rm d}\ln D 00100 !! \f] 00101 subroutine soln_additive_exp(bin_grid, aero_data, time, num_conc, & 00102 radius_at_mean_vol, env_state, aero_binned) 00103 00104 !> Bin grid. 00105 type(bin_grid_t), intent(in) :: bin_grid 00106 !> Aerosol data. 00107 type(aero_data_t), intent(in) :: aero_data 00108 !> Current time. 00109 real(kind=dp), intent(in) :: time 00110 !> Particle number concentration (#/m^3). 00111 real(kind=dp), intent(in) :: num_conc 00112 !> Mean init radius (m). 00113 real(kind=dp), intent(in) :: radius_at_mean_vol 00114 !> Environment state. 00115 type(env_state_t), intent(in) :: env_state 00116 !> Output state. 00117 type(aero_binned_t), intent(inout) :: aero_binned 00118 00119 real(kind=dp) :: tau, T, rat_v, nn, b, x, mean_vol 00120 integer :: k 00121 00122 mean_vol = rad2vol(radius_at_mean_vol) 00123 if (time .eq. 0d0) then 00124 do k = 1,bin_grid%n_bin 00125 aero_binned%num_conc(k) = const%pi/2d0 & 00126 * (2d0 * bin_grid%center_radius(k))**3 * num_conc / mean_vol & 00127 * exp(-(rad2vol(bin_grid%center_radius(k)) / mean_vol)) 00128 end do 00129 else 00130 tau = num_conc * mean_vol * beta_1 * time 00131 T = 1d0 - exp(-tau) 00132 do k = 1,bin_grid%n_bin 00133 rat_v = rad2vol(bin_grid%center_radius(k)) / mean_vol 00134 x = 2d0 * rat_v * sqrt(T) 00135 if (x .lt. 500d0) then 00136 call bessi1(x, b) 00137 nn = num_conc / rad2vol(bin_grid%center_radius(k)) & 00138 * (1d0 - T) / sqrt(T) * exp(-((1d0 + T) * rat_v)) * b 00139 else 00140 ! For very large volumes we can use the asymptotic 00141 ! approximation I_1(x) \approx e^x / sqrt(2 pi x) and 00142 ! simplify the result to avoid the overflow from 00143 ! multiplying a huge bessel function result by a very 00144 ! tiny exponential. 00145 nn = num_conc / rad2vol(bin_grid%center_radius(k)) & 00146 * (1d0 - T) / sqrt(T) & 00147 * exp((2d0*sqrt(T) - T - 1d0) * rat_v) & 00148 / sqrt(4d0 * const%pi * rat_v * sqrt(T)) 00149 end if 00150 aero_binned%num_conc(k) = const%pi/2d0 & 00151 * (2d0 * bin_grid%center_radius(k))**3 * nn 00152 end do 00153 end if 00154 00155 aero_binned%vol_conc = 0d0 00156 do k = 1,bin_grid%n_bin 00157 aero_binned%vol_conc(k,1) = rad2vol(bin_grid%center_radius(k)) & 00158 * aero_binned%num_conc(k) 00159 end do 00160 00161 end subroutine soln_additive_exp 00162 00163 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00164 00165 !> Modified Bessel function of the first kind \f$ I_1(x) \f$. 00166 !! 00167 !! This looks like it was taken from Numerical Recipes. 00168 !! 00169 !! FIXME: Where did this code come from? What license does it have? 00170 subroutine bessi1(x, r) 00171 00172 !> Function argument. 00173 real(kind=dp), intent(in) :: x 00174 !> Function value. 00175 real(kind=dp), intent(out) :: r 00176 00177 real(kind=dp) ax 00178 real(kind=dp) p1,p2,p3,p4,p5,p6,p7,q1,q2,q3,q4,q5,q6,q7,q8,q9,y 00179 data p1,p2,p3,p4,p5,p6,p7/0.5d0,0.87890594d0,0.51498869d0, & 00180 0.15084934d0,0.2658733d-1,0.301532d-2,0.32411d-3/ 00181 data q1,q2,q3,q4,q5,q6,q7,q8,q9/0.39894228d0,-0.3988024d-1, & 00182 -0.362018d-2,0.163801d-2,-0.1031555d-1,0.2282967d-1, & 00183 -0.2895312d-1,0.1787654d-1,-0.420059d-2/ 00184 00185 if (abs(x) .lt. 3.75d0) then 00186 y = (x / 3.75d0)**2 00187 r = x*(p1+y*(p2+y*(p3+y*(p4+y*(p5+y*(p6+y*p7)))))) 00188 elseif (abs(x) .lt. 500d0) then 00189 ax = abs(x) 00190 y = 3.75d0 / ax 00191 r = (exp(ax)/sqrt(ax))*(q1+y*(q2+y*(q3+y*(q4+ & 00192 y*(q5+y*(q6+y*(q7+y*(q8+y*q9)))))))) 00193 if (x .lt. 0d0) r = -r 00194 else 00195 00196 end if 00197 00198 end subroutine bessi1 00199 00200 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00201 00202 end module pmc_coag_kernel_additive