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