PartMC  2.6.1
coag_kernel_constant.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_constant module.
7 
8 !> Constant coagulation kernel.
10 
11  use pmc_env_state
12  use pmc_bin_grid
13  use pmc_util
14  use pmc_constants
15  use pmc_aero_binned
16  use pmc_aero_data
17  use pmc_aero_dist
18  use pmc_aero_data
20 
21  !> Coefficient for constant kernel.
22  real(kind=dp), parameter :: beta_0 = 0.25d0 / (60d0 * 2d8)
23 
24 contains
25 
26 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
27 
28  !> Constant coagulation kernel.
29  subroutine kernel_constant(aero_particle_1, aero_particle_2, &
30  aero_data, env_state, k)
31 
32  !> First particle.
33  type(aero_particle_t), intent(in) :: aero_particle_1
34  !> Second particle.
35  type(aero_particle_t), intent(in) :: aero_particle_2
36  !> Aerosol data.
37  type(aero_data_t), intent(in) :: aero_data
38  !> Environment state.
39  type(env_state_t), intent(in) :: env_state
40  !> Coagulation kernel.
41  real(kind=dp), intent(out) :: k
42 
43  k = beta_0
44 
45  end subroutine kernel_constant
46 
47 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
48 
49  !> Minimum and maximum values of the constant coagulation kernel.
50  subroutine kernel_constant_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 = beta_0
66  k_max = beta_0
67 
68  end subroutine kernel_constant_minmax
69 
70 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
71 
72  !> Exact solution with a constant coagulation kernel and an
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 time \f$\tau = N_0 \beta_0 t\f$, where \f$\beta_0\f$ is
78  !! the fixed constant kernel value. We also set the parameter
79  !! \f$\lambda = 1\f$. Then the solution is
80  !! \f[
81  !! n(D,t) \ {\rm d}\ln D
82  !! = \frac{\pi}{2} D^3 \frac{N_0}{v_\mu} \frac{4}{(\tau + 2)^2}
83  !! \exp\left(-\frac{v}{v_\mu} \frac{2}{\tau + 2}
84  !! \exp(-\lambda \tau) - \lambda \tau\right) {\rm d}\ln D
85  !! \f]
86  !! This thus has initial condition
87  !! \f[
88  !! n(D,t) \ {\rm d}\ln D
89  !! = \frac{\pi}{2} D^3 \frac{N_0}{v_\mu}
90  !! \exp\left(-\frac{v}{v_\mu}\right) {\rm d}\ln D
91  !! \f]
92  subroutine soln_constant_exp(bin_grid, aero_data, time, num_conc, &
93  radius_at_mean_vol, env_state, aero_binned)
94 
95  !> Bin grid.
96  type(bin_grid_t), intent(in) :: bin_grid
97  !> Aerosol data.
98  type(aero_data_t), intent(in) :: aero_data
99  !> Current time.
100  real(kind=dp), intent(in) :: time
101  !> Particle number concentration (#/m^3).
102  real(kind=dp), intent(in) :: num_conc
103  !> Mean init radius (m).
104  real(kind=dp), intent(in) :: radius_at_mean_vol
105  !> Environment state.
106  type(env_state_t), intent(in) :: env_state
107  !> Output state.
108  type(aero_binned_t), intent(inout) :: aero_binned
109 
110  real(kind=dp) :: tau, t, rat_v, nn, b, sigma, mean_vol
111  integer :: k
112 
113  real(kind=dp), parameter :: lambda = 1d0
114 
115  call aero_binned_set_sizes(aero_binned, bin_grid_size(bin_grid), &
116  aero_data_n_spec(aero_data))
117 
118  mean_vol = aero_data_rad2vol(aero_data, radius_at_mean_vol)
119  if (time .eq. 0d0) then
120  do k = 1,bin_grid_size(bin_grid)
121  aero_binned%num_conc(k) = const%pi / 2d0 &
122  * (2d0 * bin_grid%centers(k))**3 * num_conc / mean_vol &
123  * exp(-(aero_data_rad2vol(aero_data, bin_grid%centers(k)) &
124  / mean_vol))
125  end do
126  else
127  tau = num_conc * beta_0 * time
128  do k = 1,bin_grid_size(bin_grid)
129  rat_v = aero_data_rad2vol(aero_data, bin_grid%centers(k)) / mean_vol
130  nn = 4d0 * num_conc / (mean_vol * ( tau + 2d0 ) ** 2d0) &
131  * exp(-2d0*rat_v/(tau+2d0)*exp(-lambda*tau)-lambda*tau)
132  aero_binned%num_conc(k) = const%pi / 2d0 &
133  * (2d0 * bin_grid%centers(k))**3d0 * nn
134  end do
135  end if
136 
137  aero_binned%vol_conc = 0d0
138  do k = 1,bin_grid_size(bin_grid)
139  aero_binned%vol_conc(k,1) = aero_data_rad2vol(aero_data, &
140  bin_grid%centers(k)) * aero_binned%num_conc(k)
141  end do
142 
143  end subroutine soln_constant_exp
144 
145 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
146 
147 end module pmc_coag_kernel_constant
pmc_aero_particle::aero_particle_t
Single aerosol particle data structure.
Definition: aero_particle.F90:26
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:236
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_aero_particle
The aero_particle_t structure and associated subroutines.
Definition: aero_particle.F90:9
pmc_coag_kernel_constant
Constant coagulation kernel.
Definition: coag_kernel_constant.F90:9
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_coag_kernel_constant::kernel_constant_minmax
subroutine kernel_constant_minmax(v1, v2, aero_data, env_state, k_min, k_max)
Minimum and maximum values of the constant coagulation kernel.
Definition: coag_kernel_constant.F90:51
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:73
pmc_coag_kernel_constant::kernel_constant
subroutine kernel_constant(aero_particle_1, aero_particle_2, aero_data, env_state, k)
Constant coagulation kernel.
Definition: coag_kernel_constant.F90:31
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:49
pmc_constants
Physical constants.
Definition: constants.F90:9
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_coag_kernel_constant::beta_0
real(kind=dp), parameter beta_0
Coefficient for constant kernel.
Definition: coag_kernel_constant.F90:22
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:122
pmc_aero_data
The aero_data_t structure and associated subroutines.
Definition: aero_data.F90:9
pmc_bin_grid::bin_grid_t
1D grid, either logarithmic or linear.
Definition: bin_grid.F90:33
pmc_coag_kernel_constant::soln_constant_exp
subroutine soln_constant_exp(bin_grid, aero_data, time, num_conc, radius_at_mean_vol, env_state, aero_binned)
Exact solution with a constant coagulation kernel and an exponential initial condition.
Definition: coag_kernel_constant.F90:94