PartMC  2.6.1
Functions/Subroutines | Variables
pmc_coag_kernel_additive Module Reference

Additive coagulation kernel. More...

Functions/Subroutines

subroutine kernel_additive (aero_particle_1, aero_particle_2, aero_data, env_state, k)
 Additive coagulation kernel. More...
 
subroutine kernel_additive_minmax (v1, v2, aero_data, env_state, k_min, k_max)
 Minimum and maximum values of the additive kernel. More...
 
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. More...
 
subroutine bessi1 (x, r)
 Modified Bessel function of the first kind $ I_1(x) $. More...
 
subroutine calci1 (arg, result, jint)
 Calculates modified Bessel functions of the first kind $ I_1(x) $. More...
 

Variables

real(kind=dp), parameter beta_1 = 1000d0
 Scaling coefficient for constant kernel. More...
 

Detailed Description

Additive coagulation kernel.

Function/Subroutine Documentation

◆ bessi1()

subroutine pmc_coag_kernel_additive::bessi1 ( real(kind=dp), intent(in)  x,
real(kind=dp), intent(out)  r 
)

Modified Bessel function of the first kind $ I_1(x) $.

Parameters
[in]xFunction argument.
[out]rFunction value.

Definition at line 172 of file coag_kernel_additive.F90.

◆ calci1()

subroutine pmc_coag_kernel_additive::calci1 ( real ( kind = 8 )  arg,
real ( kind = 8 )  result,
integer ( kind = 4 )  jint 
)

Calculates modified Bessel functions of the first kind $ I_1(x) $.

Definition at line 186 of file coag_kernel_additive.F90.

◆ kernel_additive()

subroutine pmc_coag_kernel_additive::kernel_additive ( type(aero_particle_t), intent(in)  aero_particle_1,
type(aero_particle_t), intent(in)  aero_particle_2,
type(aero_data_t), intent(in)  aero_data,
type(env_state_t), intent(in)  env_state,
real(kind=dp), intent(out)  k 
)

Additive coagulation kernel.

Parameters
[in]aero_particle_1First particle.
[in]aero_particle_2Second particle.
[in]aero_dataAerosol data.
[in]env_stateEnvironment state.
[out]kCoagulation kernel.

Definition at line 30 of file coag_kernel_additive.F90.

◆ kernel_additive_minmax()

subroutine pmc_coag_kernel_additive::kernel_additive_minmax ( real(kind=dp), intent(in)  v1,
real(kind=dp), intent(in)  v2,
type(aero_data_t), intent(in)  aero_data,
type(env_state_t), intent(in)  env_state,
real(kind=dp), intent(out)  k_min,
real(kind=dp), intent(out)  k_max 
)

Minimum and maximum values of the additive kernel.

Parameters
[in]v1Volume of first particle.
[in]v2Volume of second particle.
[in]aero_dataAerosol data.
[in]env_stateEnvironment state.
[out]k_minCoagulation kernel minimum value.
[out]k_maxCoagulation kernel maximum value.

Definition at line 52 of file coag_kernel_additive.F90.

◆ soln_additive_exp()

subroutine pmc_coag_kernel_additive::soln_additive_exp ( type(bin_grid_t), intent(in)  bin_grid,
type(aero_data_t), intent(in)  aero_data,
real(kind=dp), intent(in)  time,
real(kind=dp), intent(in)  num_conc,
real(kind=dp), intent(in)  radius_at_mean_vol,
type(env_state_t), intent(in)  env_state,
type(aero_binned_t), intent(inout)  aero_binned 
)

Exact solution with the additive coagulation kernel and exponential initial condition.

Given input paramaters $R$ and $N_0$ we let the mean volume be $v_\mu = \frac{4\pi}{3} R^3$ and define the rescaled times $\tau = N_0 v_\mu \beta_1 t$ and $T = 1 - e^{-\tau}$, where $\beta_1$ is the fixed kernel scaling parameter. Then the solution is

\[ n(D,t) \ {\rm d}\ln D = \frac{\pi}{2} D^3 \frac{N_0}{v} \frac{1 - T}{\sqrt{T}} \exp\left(-(1 + T) \frac{v}{v_\mu}\right) I_1\left(2 \frac{v}{v_\mu} \sqrt{T}\right) {\rm d}\ln D \]

where $I_1(x)$ is the modified Bessel function of the first kind and $v = \frac{\pi}{6} D^3$.

For small $x$ we have $I_1(x) \approx \frac{x}{2}$, so this solution has initial condition

\[ n(D,t) \ {\rm d}\ln D = \frac{\pi}{2} D^3 \frac{N_0}{v_\mu} \exp\left(-\frac{v}{v_\mu}\right) {\rm d}\ln D \]

Parameters
[in]bin_gridBin grid.
[in]aero_dataAerosol data.
[in]timeCurrent time.
[in]num_concParticle number concentration (#/m^3).
[in]radius_at_mean_volMean init radius (m).
[in]env_stateEnvironment state.
[in,out]aero_binnedOutput state.

Definition at line 101 of file coag_kernel_additive.F90.

Variable Documentation

◆ beta_1

real(kind=dp), parameter pmc_coag_kernel_additive::beta_1 = 1000d0

Scaling coefficient for constant kernel.

Definition at line 23 of file coag_kernel_additive.F90.