PartMC  2.8.0
coag_kernel_additive.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2012 Nicole Riemer and Matthew West
2 ! Licensed under the GNU General Public License version 2 or (at your
3 ! option) any later version. See the file COPYING for details.
4 
5 !> \file
6 !> The pmc_coag_kernel_additive module.
7 
8 !> Additive coagulation kernel.
10 
11  use pmc_bin_grid
12  use pmc_env_state
13  use pmc_util
14  use pmc_constants
15  use pmc_constants
16  use pmc_aero_binned
17  use pmc_aero_data
18  use pmc_aero_dist
19  use pmc_aero_data
21 
22 contains
23 
24 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
25 
26  !> Additive coagulation kernel.
27  subroutine kernel_additive(aero_particle_1, aero_particle_2, &
28  aero_data, env_state, k)
29 
30  !> First particle.
31  type(aero_particle_t), intent(in) :: aero_particle_1
32  !> Second particle.
33  type(aero_particle_t), intent(in) :: aero_particle_2
34  !> Aerosol data.
35  type(aero_data_t), intent(in) :: aero_data
36  !> Environment state.
37  type(env_state_t), intent(in) :: env_state
38  !> Coagulation kernel.
39  real(kind=dp), intent(out) :: k
40 
41  k = env_state%additive_kernel_coefficient &
42  * (aero_particle_volume(aero_particle_1) &
43  + aero_particle_volume(aero_particle_2))
44 
45  end subroutine kernel_additive
46 
47 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
48 
49  !> Minimum and maximum values of the additive kernel.
50  subroutine kernel_additive_minmax(v1, v2, aero_data, env_state, k_min, k_max)
51 
52  !> Volume of first particle.
53  real(kind=dp), intent(in) :: v1
54  !> Volume of second particle.
55  real(kind=dp), intent(in) :: v2
56  !> Aerosol data.
57  type(aero_data_t), intent(in) :: aero_data
58  !> Environment state.
59  type(env_state_t), intent(in) :: env_state
60  !> Coagulation kernel minimum value.
61  real(kind=dp), intent(out) :: k_min
62  !> Coagulation kernel maximum value.
63  real(kind=dp), intent(out) :: k_max
64 
65  k_min = env_state%additive_kernel_coefficient * (v1 + v2)
66  k_max = k_min
67 
68  end subroutine kernel_additive_minmax
69 
70 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
71 
72  !> Exact solution with the additive coagulation kernel and
73  !> exponential initial condition.
74  !!
75  !! Given input paramaters \f$R\f$ and \f$N_0\f$ we let the mean
76  !! volume be \f$v_\mu = \frac{4\pi}{3} R^3\f$ and define the
77  !! rescaled times \f$\tau = N_0 v_\mu \beta_1 t\f$ and \f$T = 1 -
78  !! e^{-\tau}\f$, where \f$\beta_1\f$ is the fixed additive kernel scaling
79  !! coefficient. Then the solution is
80  !! \f[
81  !! n(D,t) \ {\rm d}\ln D
82  !! = \frac{\pi}{2} D^3
83  !! \frac{N_0}{v} \frac{1 - T}{\sqrt{T}}
84  !! \exp\left(-(1 + T) \frac{v}{v_\mu}\right)
85  !! I_1\left(2 \frac{v}{v_\mu} \sqrt{T}\right) {\rm d}\ln D
86  !! \f]
87  !! where \f$I_1(x)\f$ is the <a
88  !! href="http://en.wikipedia.org/wiki/Bessel_function">modified
89  !! Bessel function of the first kind</a> and \f$v = \frac{\pi}{6}
90  !! D^3\f$.
91  !!
92  !! For small \f$x\f$ we have \f$I_1(x) \approx \frac{x}{2}\f$, so
93  !! this solution has initial condition
94  !! \f[
95  !! n(D,t) \ {\rm d}\ln D
96  !! = \frac{\pi}{2} D^3 \frac{N_0}{v_\mu}
97  !! \exp\left(-\frac{v}{v_\mu}\right) {\rm d}\ln D
98  !! \f]
99  subroutine soln_additive_exp(bin_grid, aero_data, time, num_conc, &
100  radius_at_mean_vol, env_state, aero_binned)
101 
102  !> Bin grid.
103  type(bin_grid_t), intent(in) :: bin_grid
104  !> Aerosol data.
105  type(aero_data_t), intent(in) :: aero_data
106  !> Current time.
107  real(kind=dp), intent(in) :: time
108  !> Particle number concentration (#/m^3).
109  real(kind=dp), intent(in) :: num_conc
110  !> Mean init radius (m).
111  real(kind=dp), intent(in) :: radius_at_mean_vol
112  !> Environment state.
113  type(env_state_t), intent(in) :: env_state
114  !> Output state.
115  type(aero_binned_t), intent(inout) :: aero_binned
116 
117  real(kind=dp) :: tau, t, rat_v, nn, b, x, mean_vol
118  integer :: k
119 
120  call aero_binned_set_sizes(aero_binned, bin_grid_size(bin_grid), &
121  aero_data_n_spec(aero_data))
122 
123  mean_vol = aero_data_rad2vol(aero_data, radius_at_mean_vol)
124  if (time .eq. 0d0) then
125  do k = 1,bin_grid_size(bin_grid)
126  aero_binned%num_conc(k) = const%pi/2d0 &
127  * (2d0 * bin_grid%centers(k))**3 * num_conc / mean_vol &
128  * exp(-(aero_data_rad2vol(aero_data, bin_grid%centers(k)) &
129  / mean_vol))
130  end do
131  else
132  tau = num_conc * mean_vol * env_state%additive_kernel_coefficient * time
133  t = 1d0 - exp(-tau)
134  do k = 1,bin_grid_size(bin_grid)
135  rat_v = aero_data_rad2vol(aero_data, bin_grid%centers(k)) / mean_vol
136  x = 2d0 * rat_v * sqrt(t)
137  if (x .lt. 500d0) then
138  call bessi1(x, b)
139  nn = num_conc / aero_data_rad2vol(aero_data, &
140  bin_grid%centers(k)) &
141  * (1d0 - t) / sqrt(t) * exp(-((1d0 + t) * rat_v)) * b
142  else
143  ! For very large volumes we can use the asymptotic
144  ! approximation I_1(x) \approx e^x / sqrt(2 pi x) and
145  ! simplify the result to avoid the overflow from
146  ! multiplying a huge bessel function result by a very
147  ! tiny exponential.
148  nn = num_conc / aero_data_rad2vol(aero_data, &
149  bin_grid%centers(k)) &
150  * (1d0 - t) / sqrt(t) &
151  * exp((2d0*sqrt(t) - t - 1d0) * rat_v) &
152  / sqrt(4d0 * const%pi * rat_v * sqrt(t))
153  end if
154  aero_binned%num_conc(k) = const%pi/2d0 &
155  * (2d0 * bin_grid%centers(k))**3 * nn
156  end do
157  end if
158 
159  aero_binned%vol_conc = 0d0
160  do k = 1,bin_grid_size(bin_grid)
161  aero_binned%vol_conc(k,1) = aero_data_rad2vol(aero_data, &
162  bin_grid%centers(k)) * aero_binned%num_conc(k)
163  end do
164 
165  end subroutine soln_additive_exp
166 
167 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
168 
169  !> Modified Bessel function of the first kind \f$ I_1(x) \f$.
170  subroutine bessi1(x, r)
171 
172  !> Function argument.
173  real(kind=dp), intent(in) :: x
174  !> Function value.
175  real(kind=dp), intent(out) :: r
176 
177  call calci1 (x, r, 1 )
178 
179  end subroutine bessi1
180 
181 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
182 
183  !> Calculates modified Bessel functions of the first kind \f$ I_1(x) \f$.
184  subroutine calci1 ( arg, result, jint )
185 
186  !*************************************************************************
187  !
188  !! CALCI1 computes various I1 Bessel functions.
189  !
190  ! Discussion:
191  !
192  ! This routine computes modified Bessel functioons of the first kind
193  ! and order one, I1(X) and EXP(-ABS(X))*I1(X), for real
194  ! arguments X.
195  !
196  ! The main computation evaluates slightly modified forms of
197  ! minimax approximations generated by Blair and Edwards, Chalk
198  ! River (Atomic Energy of Canada Limited) Report AECL-4928,
199  ! October, 1974.
200  !
201  ! Licensing:
202  !
203  ! This code is distributed under the GNU LGPL license.
204  !
205  ! Modified:
206  !
207  ! 03 April 2007
208  !
209  ! Author:
210  !
211  ! Original FORTRAN77 version by William Cody, Laura Stoltz.
212  ! FORTRAN90 version by John Burkardt.
213  !
214  ! Parameters:
215  !
216  ! Input, real ( kind = 8 ) ARG, the argument. If JINT = 1, then
217  ! the argument must be less than XMAX.
218  !
219  ! Output, real ( kind = 8 ) RESULT, the value of the function,
220  ! which depends on the input value of JINT:
221  ! 1, RESULT = I1(x);
222  ! 2, RESULT = exp(-x) * I1(x);
223  !
224  ! Input, integer ( kind = 4 ) JINT, chooses the function to be computed.
225  ! 1, I1(x);
226  ! 2, exp(-x) * I1(x);
227  !
228 
229  real ( kind = 8 ) a
230  real ( kind = 8 ) arg
231  real ( kind = 8 ) b
232  real ( kind = 8 ) exp40
233  real ( kind = 8 ) forty
234  integer ( kind = 4 ) j
235  integer ( kind = 4 ) jint
236  real ( kind = 8 ) one5
237  real ( kind = 8 ) p(15)
238  real ( kind = 8 ) pbar
239  real ( kind = 8 ) pp(8)
240  real ( kind = 8 ) q(5)
241  real ( kind = 8 ) qq(6)
242  real ( kind = 8 ) rec15
243  real ( kind = 8 ) result
244  real ( kind = 8 ) sump
245  real ( kind = 8 ) sumq
246  real ( kind = 8 ) two25
247  real ( kind = 8 ) x
248  real ( kind = 8 ) xinf
249  real ( kind = 8 ) xmax
250  real ( kind = 8 ) xsmall
251  real ( kind = 8 ) xx
252  !
253  ! Mathematical constants
254  !
255  data one5 / 15.0d0 /
256  data exp40 / 2.353852668370199854d17 /
257  data forty / 40.0d0 /
258  data rec15 / 6.6666666666666666666d-2 /
259  data two25 / 225.0d0 /
260  !
261  ! Machine-dependent constants
262  !
263  data xsmall /5.55d-17/
264  data xinf /1.79d308/
265  data xmax /713.987d0/
266  !
267  ! Coefficients for XSMALL <= ABS(ARG) < 15.0
268  !
269  data p/-1.9705291802535139930d-19,-6.5245515583151902910d-16, &
270  -1.1928788903603238754d-12,-1.4831904935994647675d-09, &
271  -1.3466829827635152875d-06,-9.1746443287817501309d-04, &
272  -4.7207090827310162436d-01,-1.8225946631657315931d+02, &
273  -5.1894091982308017540d+04,-1.0588550724769347106d+07, &
274  -1.4828267606612366099d+09,-1.3357437682275493024d+11, &
275  -6.9876779648010090070d+12,-1.7732037840791591320d+14, &
276  -1.4577180278143463643d+15/
277  data q/-4.0076864679904189921d+03, 7.4810580356655069138d+06, &
278  -8.0059518998619764991d+09, 4.8544714258273622913d+12, &
279  -1.3218168307321442305d+15/
280  !
281  ! Coefficients for 15.0 <= ABS(ARG)
282  !
283  data pp/-6.0437159056137600000d-02, 4.5748122901933459000d-01, &
284  -4.2843766903304806403d-01, 9.7356000150886612134d-02, &
285  -3.2457723974465568321d-03,-3.6395264712121795296d-04, &
286  1.6258661867440836395d-05,-3.6347578404608223492d-07/
287  data qq/-3.8806586721556593450d+00, 3.2593714889036996297d+00, &
288  -8.5017476463217924408d-01, 7.4212010813186530069d-02, &
289  -2.2835624489492512649d-03, 3.7510433111922824643d-05/
290  data pbar/3.98437500d-01/
291 
292  x = abs( arg )
293  !
294  ! Return for ABS(ARG) < XSMALL.
295  !
296  if ( x < xsmall ) then
297 
298  result = 0.5d+00 * x
299  !
300  ! XSMALL <= ABS(ARG) < 15.0.
301  !
302  else if ( x < one5 ) then
303 
304  xx = x * x
305  sump = p(1)
306  do j = 2, 15
307  sump = sump * xx + p(j)
308  end do
309  xx = xx - two25
310 
311  sumq = (((( &
312  xx + q(1) ) &
313  * xx + q(2) ) &
314  * xx + q(3) ) &
315  * xx + q(4) ) &
316  * xx + q(5)
317 
318  result = ( sump / sumq ) * x
319 
320  if ( jint == 2 ) then
321  result = result * exp( -x )
322  end if
323 
324  else if ( jint == 1 .and. xmax < x ) then
325 
326  result = xinf
327 
328  else
329  !
330  ! 15.0 <= ABS(ARG).
331  !
332  xx = 1.0d+00 / x - rec15
333 
334  sump = (((((( &
335  pp(1) &
336  * xx + pp(2) ) &
337  * xx + pp(3) ) &
338  * xx + pp(4) ) &
339  * xx + pp(5) ) &
340  * xx + pp(6) ) &
341  * xx + pp(7) ) &
342  * xx + pp(8)
343 
344  sumq = ((((( &
345  xx + qq(1) ) &
346  * xx + qq(2) ) &
347  * xx + qq(3) ) &
348  * xx + qq(4) ) &
349  * xx + qq(5) ) &
350  * xx + qq(6)
351 
352  result = sump / sumq
353 
354  if ( jint /= 1 ) then
355  result = ( result + pbar ) / sqrt( x )
356  else
357  !
358  ! Calculation reformulated to avoid premature overflow.
359  !
360  if ( xmax - one5 < x ) then
361  a = exp( x - forty )
362  b = exp40
363  else
364  a = exp( x )
365  b = 1.0d+00
366  end if
367 
368  result = ( ( result * a + pbar * a ) / sqrt( x ) ) * b
369 
370  end if
371  end if
372 
373  if ( arg < 0.0d+00 ) then
374  result = -result
375  end if
376 
377  return
378 
379  end subroutine calci1
380 
381 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
382 
383 end module pmc_coag_kernel_additive
pmc_aero_particle::aero_particle_t
Single aerosol particle data structure.
Definition: aero_particle.F90:30
pmc_aero_data::aero_data_n_spec
elemental integer function aero_data_n_spec(aero_data)
Return the number of aerosol species, or -1 if uninitialized.
Definition: aero_data.F90:246
pmc_aero_binned::aero_binned_set_sizes
subroutine aero_binned_set_sizes(aero_binned, n_bin, n_spec)
Set the number of bins and species in an aero_binned_t.
Definition: aero_binned.F90:72
pmc_coag_kernel_additive
Additive coagulation kernel.
Definition: coag_kernel_additive.F90:9
pmc_aero_particle
The aero_particle_t structure and associated subroutines.
Definition: aero_particle.F90:9
pmc_coag_kernel_additive::soln_additive_exp
subroutine soln_additive_exp(bin_grid, aero_data, time, num_conc, radius_at_mean_vol, env_state, aero_binned)
Exact solution with the additive coagulation kernel and exponential initial condition.
Definition: coag_kernel_additive.F90:101
pmc_coag_kernel_additive::kernel_additive
subroutine kernel_additive(aero_particle_1, aero_particle_2, aero_data, env_state, k)
Additive coagulation kernel.
Definition: coag_kernel_additive.F90:29
pmc_constants::dp
integer, parameter dp
Kind of a double precision real number.
Definition: constants.F90:12
pmc_bin_grid::bin_grid_size
elemental integer function bin_grid_size(bin_grid)
Return the number of bins in the grid, or -1 if the bin grid is not allocated.
Definition: bin_grid.F90:51
pmc_env_state::env_state_t
Current environment state.
Definition: env_state.F90:29
pmc_aero_dist
The aero_dist_t structure and associated subroutines.
Definition: aero_dist.F90:18
pmc_constants::const
type(const_t), save const
Fixed variable for accessing the constant's values.
Definition: constants.F90:75
pmc_env_state
The env_state_t structure and associated subroutines.
Definition: env_state.F90:9
pmc_aero_data::aero_data_t
Aerosol material properties and associated data.
Definition: aero_data.F90:55
pmc_constants
Physical constants.
Definition: constants.F90:9
pmc_coag_kernel_additive::calci1
subroutine calci1(arg, result, jint)
Calculates modified Bessel functions of the first kind .
Definition: coag_kernel_additive.F90:185
pmc_util
Common utility subroutines.
Definition: util.F90:9
pmc_aero_binned
The aero_binned_t structure and associated subroutines.
Definition: aero_binned.F90:9
pmc_aero_binned::aero_binned_t
Aerosol number and volume distributions stored per bin.
Definition: aero_binned.F90:37
pmc_bin_grid
The bin_grid_t structure and associated subroutines.
Definition: bin_grid.F90:9
pmc_aero_data::aero_data_rad2vol
real(kind=dp) elemental function aero_data_rad2vol(aero_data, r)
Convert geometric radius (m) to mass-equivalent volume (m^3).
Definition: aero_data.F90:132
pmc_aero_data
The aero_data_t structure and associated subroutines.
Definition: aero_data.F90:9
pmc_coag_kernel_additive::bessi1
subroutine bessi1(x, r)
Modified Bessel function of the first kind .
Definition: coag_kernel_additive.F90:171
pmc_bin_grid::bin_grid_t
1D grid, either logarithmic or linear.
Definition: bin_grid.F90:33
pmc_aero_particle::aero_particle_volume
elemental real(kind=dp) function aero_particle_volume(aero_particle)
Total volume of the particle (m^3).
Definition: aero_particle.F90:318
pmc_coag_kernel_additive::kernel_additive_minmax
subroutine kernel_additive_minmax(v1, v2, aero_data, env_state, k_min, k_max)
Minimum and maximum values of the additive kernel.
Definition: coag_kernel_additive.F90:51