PartMC  2.2.0
coag_kernel_sedi.F90
Go to the documentation of this file.
00001 ! Copyright (C) 2005-2011 Nicole Riemer and Matthew West
00002 ! Copyright (C) Andreas Bott
00003 ! Licensed under the GNU General Public License version 2 or (at your
00004 ! option) any later version. See the file COPYING for details.
00005 
00006 !> \file
00007 !> The pmc_coag_kernel_sedi module.
00008 !!
00009 !! Contains code based on \c coad1d.f by Andreas Bott
00010 !!     - http://www.meteo.uni-bonn.de/mitarbeiter/ABott/
00011 !!     - Released under the GPL to Nicole Riemer (personal communication)
00012 !!     - A. Bott, A flux method for the numerical solution of the
00013 !!       stochastic collection equation, J. Atmos. Sci. 55, 2284-2293,
00014 !!       1998.
00015 
00016 !> Gravitational sedimentation coagulation kernel.
00017 module pmc_coag_kernel_sedi
00018 
00019   use pmc_env_state
00020   use pmc_constants
00021   use pmc_aero_data
00022   use pmc_aero_particle
00023   
00024 contains
00025   
00026 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00027 
00028   !> Sedimentation coagulation kernel.
00029   subroutine kernel_sedi(aero_particle_1, aero_particle_2, &
00030        aero_data, env_state, k)
00031 
00032     !> First particle.
00033     type(aero_particle_t), intent(in) :: aero_particle_1
00034     !> Second particle.
00035     type(aero_particle_t), intent(in) :: aero_particle_2
00036     !> Aerosol data.
00037     type(aero_data_t), intent(in) :: aero_data
00038     !> Environment state.
00039     type(env_state_t), intent(in) :: env_state
00040     !> Kernel \c k(a,b) (m^3/s).
00041     real(kind=dp), intent(out) :: k
00042 
00043     call kernel_sedi_helper(aero_particle_volume(aero_particle_1), &
00044          aero_particle_volume(aero_particle_2), k)
00045 
00046   end subroutine kernel_sedi
00047   
00048 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00049 
00050   !> Minimum and maximum values of the sedimentation coagulation.
00051   subroutine kernel_sedi_minmax(v1, v2, aero_data, env_state, k_min, k_max)
00052 
00053     !> Volume of first particle (m^3).
00054     real(kind=dp), intent(in) :: v1
00055     !> Volume of second particle (m^3).
00056     real(kind=dp), intent(in) :: v2
00057     !> Aerosol data.
00058     type(aero_data_t), intent(in) :: aero_data
00059     !> Environment state.
00060     type(env_state_t), intent(in) :: env_state
00061     !> Minimum kernel \c k(a,b) (m^3/s).
00062     real(kind=dp), intent(out) :: k_min
00063     !> Maximum kernel \c k(a,b) (m^3/s).
00064     real(kind=dp), intent(out) :: k_max
00065 
00066     call kernel_sedi_helper(v1, v2, k_min)
00067     k_max = k_min
00068 
00069   end subroutine kernel_sedi_minmax
00070 
00071 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00072 
00073   !> Helper function that does the actual sedimentation kernel computation.
00074   !!
00075   !! Helper function. Do not call directly. Instead use kernel_sedi().
00076   subroutine kernel_sedi_helper(v1, v2, k)
00077 
00078     !> Volume of first particle (m^3).
00079     real(kind=dp), intent(in) :: v1
00080     !> Volume of second particle (m^3).
00081     real(kind=dp), intent(in) :: v2
00082     !> Kernel k(a,b) (m^3/s).
00083     real(kind=dp), intent(out) :: k
00084 
00085     real(kind=dp) r1, r2, winf1, winf2, ec
00086     
00087     r1 = vol2rad(v1) ! m
00088     r2 = vol2rad(v2) ! m
00089     call fall_g(r1, winf1) ! winf1 in m/s
00090     call fall_g(r2, winf2) ! winf2 in m/s
00091     call effic(r1, r2, ec) ! ec is dimensionless
00092     k = ec * const%pi * (r1 + r2)**2 * abs(winf1 - winf2)
00093 
00094   end subroutine kernel_sedi_helper
00095 
00096 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00097 
00098   !> Finds the terminal velocity of a particle based on its size.
00099   subroutine fall_g(r, w_inf)
00100     
00101     !> Particle radius (m).
00102     real(kind=dp), intent(in) :: r
00103     !> Terminal velocity (m/s).
00104     real(kind=dp), intent(out) :: w_inf
00105     
00106     ! terminal velocity of falling drops
00107     real(kind=dp) eta, xlamb, rhow, rhoa, grav, cunh, t0, sigma
00108     real(kind=dp) stok, stb, phy, py, rr, x, y, xrey, bond
00109     integer i
00110     real(kind=dp) b(7),c(6)
00111     data b /-0.318657d1,0.992696d0,-0.153193d-2,-0.987059d-3, &
00112          -0.578878d-3,0.855176d-4,-0.327815d-5/
00113     data c /-0.500015d1,0.523778d1,-0.204914d1,0.475294d0, &
00114          -0.542819d-1,0.238449d-2/
00115     
00116     eta = 1.818d-4
00117     xlamb = 6.62d-6
00118     rhow = 1d0
00119     rhoa = 1.225d-3
00120     grav = 980.665d0
00121     cunh = 1.257d0 * xlamb
00122     t0 = 273.15d0
00123     sigma = 76.1d0 - 0.155d0 * (293.15d0 - t0)
00124     stok = 2d0 * grav * (rhow - rhoa) / (9d0 * eta)
00125     stb = 32d0 * rhoa * (rhow - rhoa) * grav / (3d0 * eta * eta)
00126     phy = sigma * sigma * sigma * rhoa * rhoa  &
00127          / (eta**4 * grav * (rhow - rhoa))
00128     py = phy**(1d0/6d0)
00129     
00130     ! rr: radius in cm-units
00131     rr = r * 1d2
00132     
00133     if (rr .le. 1d-3) then
00134        w_inf = stok * (rr * rr + cunh * rr)
00135     elseif (rr .gt. 1d-3 .and. rr .le. 5.35d-2) then
00136        x = log(stb * rr * rr * rr)
00137        y = 0d0
00138        do i = 1,7
00139           y = y + b(i) * (x**(i - 1))
00140        end do
00141        xrey = (1d0 + cunh/rr) * exp(y)
00142        w_inf = xrey * eta / (2d0 * rhoa * rr)
00143     elseif (rr .gt. 5.35d-2) then
00144        bond = grav * (rhow - rhoa) * rr**2 / sigma
00145        if (rr .gt. 0.35d0) then
00146           bond = grav * (rhow - rhoa) * 0.35d0**2 / sigma
00147        end if
00148        x = log(16d0 * bond * py / 3d0)
00149        y = 0d0
00150        do i = 1,6
00151           y = y + c(i) * (x**(i - 1))
00152        end do
00153        xrey = py * exp(y)
00154        w_inf = xrey * eta / (2d0 * rhoa * rr)
00155        if (rr .gt. 0.35d0) then
00156           w_inf = xrey * eta / (2d0 * rhoa * 0.35d0)
00157        end if
00158     end if
00159     w_inf = w_inf / 100d0
00160     
00161   end subroutine fall_g
00162   
00163 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00164 
00165   !> Coagulation efficiency.
00166   !!
00167   !! Determines the chance that two particles will actually coagulate,
00168   !! given that they approach close enough to do so.
00169   subroutine effic(r1, r2, ec)
00170 
00171     !> Radius of first particle (m).
00172     real(kind=dp), intent(in) :: r1
00173     !> Radius of second particle (m).
00174     real(kind=dp), intent(in) :: r2
00175     !> Collision efficiency (dimensionless).
00176     real(kind=dp), intent(out) :: ec
00177     
00178     real(kind=dp) :: r_small, r_big, rq, p, q, ek
00179     integer :: k, ir, kk, iq
00180     ! collision efficiencies of hall kernel
00181     real(kind=dp) :: rat(21),r0(15),ecoll(15,21)
00182 
00183     data r0 /6.0d0,8.0d0,10.0d0,15.0d0,20.0d0,25.0d0,30.0d0,40.0d0 &
00184          ,50.0d0,60.0d0,70.0d0,100.0d0,150.0d0,200.0d0,300.0d0/
00185     data rat /0.0d0,0.05d0,0.1d0,0.15d0,0.2d0,0.25d0,0.3d0,0.35d0 &
00186          ,0.4d0,0.45d0,0.5d0,0.55d0,0.6d0,0.65d0,0.7d0,0.75d0,0.8d0 &
00187          ,0.85d0,0.9d0,0.95d0,1.0d0/
00188     ! two-dimensional linear interpolation of the collision efficiency
00189     data ecoll /0.001d0,0.001d0,0.001d0,0.001d0,0.001d0,0.001d0 &
00190          ,0.001d0,0.001d0,0.001d0,0.001d0 ,0.001d0,0.001d0,0.001d0 &
00191          ,0.001d0,0.001d0,0.003d0,0.003d0,0.003d0,0.004d0,0.005d0 &
00192          ,0.005d0,0.005d0,0.010d0,0.100d0,0.050d0,0.200d0,0.500d0 &
00193          ,0.770d0,0.870d0,0.970d0 ,0.007d0,0.007d0,0.007d0,0.008d0 &
00194          ,0.009d0,0.010d0,0.010d0,0.070d0,0.400d0,0.430d0 ,0.580d0 &
00195          ,0.790d0,0.930d0,0.960d0,1.000d0,0.009d0,0.009d0,0.009d0 &
00196          ,0.012d0,0.015d0 ,0.010d0,0.020d0,0.280d0,0.600d0,0.640d0 &
00197          ,0.750d0,0.910d0,0.970d0,0.980d0,1.000d0 ,0.014d0,0.014d0 &
00198          ,0.014d0,0.015d0,0.016d0,0.030d0,0.060d0,0.500d0,0.700d0 &
00199          ,0.770d0 ,0.840d0,0.950d0,0.970d0,1.000d0,1.000d0,0.017d0 &
00200          ,0.017d0,0.017d0,0.020d0,0.022d0 ,0.060d0,0.100d0,0.620d0 &
00201          ,0.780d0,0.840d0,0.880d0,0.950d0,1.000d0,1.000d0,1.000d0 &
00202          ,0.030d0,0.030d0,0.024d0,0.022d0,0.032d0,0.062d0,0.200d0 &
00203          ,0.680d0,0.830d0,0.870d0 ,0.900d0,0.950d0,1.000d0,1.000d0 &
00204          ,1.000d0,0.025d0,0.025d0,0.025d0,0.036d0,0.043d0 ,0.130d0 &
00205          ,0.270d0,0.740d0,0.860d0,0.890d0,0.920d0,1.000d0,1.000d0 &
00206          ,1.000d0,1.000d0 ,0.027d0,0.027d0,0.027d0,0.040d0,0.052d0 &
00207          ,0.200d0,0.400d0,0.780d0,0.880d0,0.900d0 ,0.940d0,1.000d0 &
00208          ,1.000d0,1.000d0,1.000d0,0.030d0,0.030d0,0.030d0,0.047d0 &
00209          ,0.064d0 ,0.250d0,0.500d0,0.800d0,0.900d0,0.910d0,0.950d0 &
00210          ,1.000d0,1.000d0,1.000d0,1.000d0 ,0.040d0,0.040d0,0.033d0 &
00211          ,0.037d0,0.068d0,0.240d0,0.550d0,0.800d0,0.900d0,0.910d0 &
00212          ,0.950d0,1.000d0,1.000d0,1.000d0,1.000d0,0.035d0,0.035d0 &
00213          ,0.035d0,0.055d0,0.079d0 ,0.290d0,0.580d0,0.800d0,0.900d0 &
00214          ,0.910d0,0.950d0,1.000d0,1.000d0,1.000d0,1.000d0 ,0.037d0 &
00215          ,0.037d0,0.037d0,0.062d0,0.082d0,0.290d0,0.590d0,0.780d0 &
00216          ,0.900d0,0.910d0 ,0.950d0,1.000d0,1.000d0,1.000d0,1.000d0 &
00217          ,0.037d0,0.037d0,0.037d0,0.060d0,0.080d0 ,0.290d0,0.580d0 &
00218          ,0.770d0,0.890d0,0.910d0,0.950d0,1.000d0,1.000d0,1.000d0 &
00219          ,1.000d0 ,0.037d0,0.037d0,0.037d0,0.041d0,0.075d0,0.250d0 &
00220          ,0.540d0,0.760d0,0.880d0,0.920d0 ,0.950d0,1.000d0,1.000d0 &
00221          ,1.000d0,1.000d0,0.037d0,0.037d0,0.037d0,0.052d0,0.067d0 &
00222          ,0.250d0,0.510d0,0.770d0,0.880d0,0.930d0,0.970d0,1.000d0 &
00223          ,1.000d0,1.000d0,1.000d0 ,0.037d0,0.037d0,0.037d0,0.047d0 &
00224          ,0.057d0,0.250d0,0.490d0,0.770d0,0.890d0,0.950d0 ,1.000d0 &
00225          ,1.000d0,1.000d0,1.000d0,1.000d0,0.036d0,0.036d0,0.036d0 &
00226          ,0.042d0,0.048d0 ,0.230d0,0.470d0,0.780d0,0.920d0,1.000d0 &
00227          ,1.020d0,1.020d0,1.020d0,1.020d0,1.020d0 ,0.040d0,0.040d0 &
00228          ,0.035d0,0.033d0,0.040d0,0.112d0,0.450d0,0.790d0,1.010d0 &
00229          ,1.030d0 ,1.040d0,1.040d0,1.040d0,1.040d0,1.040d0,0.033d0 &
00230          ,0.033d0,0.033d0,0.033d0,0.033d0 ,0.119d0,0.470d0,0.950d0 &
00231          ,1.300d0,1.700d0,2.300d0,2.300d0,2.300d0,2.300d0,2.300d0 &
00232          ,0.027d0,0.027d0,0.027d0,0.027d0,0.027d0,0.125d0,0.520d0 &
00233          ,1.400d0,2.300d0,3.000d0 ,4.000d0,4.000d0,4.000d0,4.000d0 &
00234          ,4.000d0/
00235     
00236     r_small = min(r1 * 1d6, r2 * 1d6) ! um
00237     r_big = max(r1 * 1d6, r2 * 1d6) ! um
00238     rq = r_small / r_big
00239     
00240     ir = 1
00241     do k = 1, 15
00242        if (r_big .gt. r0(k)) then
00243           ir = k + 1
00244        end if
00245     end do
00246     
00247     iq = 1
00248     do kk = 1,21
00249        if (rq .gt. rat(kk)) then
00250           iq = kk + 1
00251        end if
00252     end do
00253     
00254     if (ir .lt. 16) then
00255        if (ir .ge. 2) then
00256           p = (r_big - r0(ir - 1)) / (r0(ir) - r0(ir - 1))
00257           q = (rq - rat(iq - 1)) / (rat(iq) - rat(iq - 1))
00258           ec = (1d0 - p) * (1d0 - q) * ecoll(ir - 1, iq - 1) &
00259                + p * (1d0 - q) * ecoll(ir, iq - 1) &
00260                + q * (1d0 - p) * ecoll(ir - 1, iq) &
00261                + p * q * ecoll(ir, iq)
00262        else
00263           q = (rq - rat(iq - 1)) / (rat(iq) - rat(iq - 1))
00264           ec = (1d0 - q) * ecoll(1, iq - 1) + q * ecoll(1, iq)
00265        end if
00266     else
00267        q = (rq - rat(iq - 1)) / (rat(iq) - rat(iq - 1))
00268        ek = (1d0 - q) * ecoll(15, iq - 1) + q * ecoll(15, iq)
00269        ec = min(ek, 1d0)
00270     end if
00271     
00272     if (ec .lt. 1d-20) stop 99
00273     
00274   end subroutine effic
00275   
00276 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00277 
00278 end module pmc_coag_kernel_sedi