Particles module. More...
Public Member Functions | |
subroutine | allocate_particles |
Particles variables inizialization. More... | |
subroutine | initialize_particles |
Particles variables inizialization. More... | |
real *8 function | particles_settling_velocity (i_part, diam_in) |
Settling velocity. More... | |
real *8 function | particles_heat_capacity (i_part, diam_in) |
Heat capacity. More... | |
real *8 function | particles_density (i_part, diam_in) |
Particle density. More... | |
real *8 function | particles_beta (i_part, diam_i, diam_j) |
Brownian aggregation. More... | |
real *8 function | aggregation_kernel (i_part, diam1, diam2) |
Aggregation kernel. More... | |
real *8 function | collision_kernel (i_part, diam1, diam2) |
Collision kernel. More... | |
real *8 function | collision_efficiency (i_part, diam1, diam2) |
Collision efficiency. More... | |
real *8 function | coalescence_efficiency (i_part, diam1, diam2) |
Coalescence efficiency. More... | |
subroutine | eval_particles_moments (xi, w) |
Particles moments computation. More... | |
Public Attributes | |
integer | n_part |
real *8, dimension(:), allocatable | solid_partial_mass_fraction |
mass fraction of the particle phases with respect to the total solid More... | |
real *8, dimension(:), allocatable | solid_partial_volume_fraction |
volume fraction of the particle phases with respect to the total solid More... | |
real *8, dimension(:), allocatable | solid_mass_fraction |
mass fraction of the particle phases with respect to the mixture More... | |
real *8, dimension(:), allocatable | solid_volume_fraction |
volume fraction of the particle phases with respect to the mixture More... | |
real *8, dimension(:,:), allocatable | mom |
Moments of the particles diameter. More... | |
real *8, dimension(:,:), allocatable | set_mom |
Moments of the settling velocities. More... | |
real *8, dimension(:,:), allocatable | rhop_mom |
Moments of the densities. More... | |
real *8, dimension(:,:), allocatable | set_rhop_mom |
Moments of the densities times the settling velocities. More... | |
real *8, dimension(:,:), allocatable | cp_rhop_mom |
Moments of the heat capacities times the densities. More... | |
real *8, dimension(:,:), allocatable | set_cp_rhop_mom |
Moments of the settling velocities times the heat cap times the densities. More... | |
real *8, dimension(:,:), allocatable | set_cp_mom |
Moments of the settling velocities times the heat capacity. More... | |
real *8, dimension(:,:), allocatable | birth_term |
Term accounting for the birth of aggregates in the moments equations. More... | |
real *8, dimension(:,:), allocatable | death_term |
Term accounting for the loss of particles because of aggregation. More... | |
real *8 | shape_factor |
shape factor for settling velocity (Pfeiffer) More... | |
real *8, dimension(:), allocatable | diam1 |
First diameter for the density function. More... | |
real *8, dimension(:), allocatable | rho1 |
Density at d=diam1. More... | |
real *8, dimension(:), allocatable | diam2 |
Second diameter for the density function. More... | |
real *8, dimension(:), allocatable | rho2 |
Density at d=diam1. More... | |
real *8, dimension(:), allocatable | cp_part |
real *8, dimension(1:50, 0:100) | mom0 |
Initial (at the base) moments of the particles diameter. More... | |
character *10 | settling_model |
Settling model: . More... | |
character(len=20) | distribution |
Ditribution of the particles: . More... | |
character(len=20) | distribution_variable |
logical | aggregation |
Flag for the aggregation: . More... | |
Particles module.
This module contains the procedures and the variables related to the solid particles. In particular, the statistical moments of the properties of the particles are defined and evaluated in this module.
Definition at line 11 of file particles.f90.
real*8 function particles_module::aggregation_kernel | ( | integer, intent(in) | i_part, |
real*8, intent(in) | diam1, | ||
real*8, intent(in) | diam2 | ||
) |
Aggregation kernel.
This function evaluates the aggregation kernel, using the expression given in Textor et al. 2006.
[in] | i_part | particle phase index |
[in] | diam1 | particle diameter (m) |
[in] | diam2 | particle diameter (m) |
Definition at line 546 of file particles.f90.
subroutine particles_module::allocate_particles | ( | ) |
Particles variables inizialization.
This subroutine allocate and evaluate the variables defining the moments for the particles. The moments are then corrected, if needed, and then the abscissas and weights for the quadrature formulas are computed.
Definition at line 113 of file particles.f90.
real*8 function particles_module::coalescence_efficiency | ( | integer | i_part, |
real*8, intent(in) | diam1, | ||
real*8, intent(in) | diam2 | ||
) |
Coalescence efficiency.
[in] | i_part | particle phase index |
[in] | diam1 | particle diameter (m) |
[in] | diam2 | particle diameter (m) |
particle Stokes number
Critical Stokes number
Efficiency exponent
Partciles settling velocities
Definition at line 724 of file particles.f90.
real*8 function particles_module::collision_efficiency | ( | integer, intent(in) | i_part, |
real*8, intent(in) | diam1, | ||
real*8, intent(in) | diam2 | ||
) |
Collision efficiency.
[in] | i_part | particle phase index |
[in] | diam1 | particle diameter (m) |
[in] | diam2 | particle diameter (m) |
Partciles settling velocities
Definition at line 652 of file particles.f90.
real*8 function particles_module::collision_kernel | ( | integer, intent(in) | i_part, |
real*8, intent(in) | diam1, | ||
real*8, intent(in) | diam2 | ||
) |
Collision kernel.
[in] | i_part | particle phase index |
[in] | diam1 | particle diameter (m) |
[in] | diam2 | particle diameter (m) |
Brownian motion collisions kernel
Laminar and turbulent fluid shear collisions kernel
Differential sedimentation kernel
Boltzmann constant
Partciles settling velocities
Gravitational collision efficiency
Rate of dissipation of turbulent kinetic energy
Fluid shear
Air kinematic viscosity
Definition at line 578 of file particles.f90.
subroutine particles_module::eval_particles_moments | ( | real*8, dimension(n_part,n_nodes), intent(in) | xi, |
real*8, dimension(n_part,n_nodes), intent(in) | w | ||
) |
Particles moments computation.
This subroutine compute the moments of the particles properties (density, heat capacity and settling velocity) using the quadrature formulas.
[in] | xi | abscissas for the quadrature |
[out] | w | weights for the quadrature |
Definition at line 794 of file particles.f90.
subroutine particles_module::initialize_particles | ( | ) |
Particles variables inizialization.
This subroutine allocate and evaluate the variables defining the moments for the particles. The moments are then corrected, if needed, and then the abscissas and weights for the quadrature formulas are computed.
Definition at line 154 of file particles.f90.
real*8 function particles_module::particles_beta | ( | integer, intent(in) | i_part, |
real*8, intent(in) | diam_i, | ||
real*8, intent(in) | diam_j | ||
) |
Brownian aggregation.
This function evaluates the aggregation kernel using a Brownian formulation given in Marchisio et al., 2003.
[in] | i_part | particle phase index |
[in] | diam_i | first particle diameter (m or phi) |
[in] | diam_j | second particle diameter (m or phi) |
Definition at line 500 of file particles.f90.
real*8 function particles_module::particles_density | ( | integer, intent(in) | i_part, |
real*8, intent(in) | diam_in | ||
) |
Particle density.
This function evaluates the density of a particle given the size (diameter), using the expression given in Bonadonna and Phillips, 2003.
[in] | i_part | particle phase index |
[in] | diam_in | particle diameter (m or phi) |
Definition at line 440 of file particles.f90.
real*8 function particles_module::particles_heat_capacity | ( | integer, intent(in) | i_part, |
real*8, intent(in) | diam_in | ||
) |
Heat capacity.
This function evaluates the heat capacity of the particles given the size (diameter).
[in] | i_part | particle phase index |
[in] | diam_in | particle diameter (m or phi) |
Definition at line 405 of file particles.f90.
real*8 function particles_module::particles_settling_velocity | ( | integer, intent(in) | i_part, |
real*8, intent(in) | diam_in | ||
) |
Settling velocity.
This function evaluates the settling velocity of a particle given the size (diameter), using the expression given in Textor et al. 2006 or in Pfeiffer et al 2005, accordingly with the variable SETTLING MODEL specified in the input file.
[in] | i_part | particle phase index |
[in] | diam_in | particle diameter (m or phi) |
cross sectional area
Drag coefficients at Rey=100,1000
Drag coefficent for intermediate values of Re
Settling velocity at Rey=100,1000
Mass of the particle
Settling velocities
Reynolds numbers for the two solutions of the settling equation
Coefficients of the settling equation
Square root of the discriminant of the settling equation
Definition at line 218 of file particles.f90.
logical particles_module::aggregation |
Flag for the aggregation:
.
Definition at line 96 of file particles.f90.
real*8, dimension(:,:), allocatable particles_module::birth_term |
Term accounting for the birth of aggregates in the moments equations.
Definition at line 53 of file particles.f90.
real*8, dimension(:), allocatable particles_module::cp_part |
Definition at line 73 of file particles.f90.
real*8, dimension(:,:), allocatable particles_module::cp_rhop_mom |
Moments of the heat capacities times the densities.
Definition at line 44 of file particles.f90.
real*8, dimension(:,:), allocatable particles_module::death_term |
Term accounting for the loss of particles because of aggregation.
Definition at line 56 of file particles.f90.
real*8, dimension(:), allocatable particles_module::diam1 |
First diameter for the density function.
Definition at line 62 of file particles.f90.
real*8, dimension(:), allocatable particles_module::diam2 |
Second diameter for the density function.
Definition at line 68 of file particles.f90.
character(len=20) particles_module::distribution |
Ditribution of the particles:
.
Definition at line 89 of file particles.f90.
character(len=20) particles_module::distribution_variable |
Definition at line 91 of file particles.f90.
real*8, dimension(:,:), allocatable particles_module::mom |
Moments of the particles diameter.
Definition at line 32 of file particles.f90.
real*8, dimension(1:50,0:100) particles_module::mom0 |
Initial (at the base) moments of the particles diameter.
Definition at line 76 of file particles.f90.
integer particles_module::n_part |
Definition at line 17 of file particles.f90.
real*8, dimension(:), allocatable particles_module::rho1 |
Density at d=diam1.
Definition at line 65 of file particles.f90.
real*8, dimension(:), allocatable particles_module::rho2 |
Density at d=diam1.
Definition at line 71 of file particles.f90.
real*8, dimension(:,:), allocatable particles_module::rhop_mom |
Moments of the densities.
Definition at line 38 of file particles.f90.
real*8, dimension(:,:), allocatable particles_module::set_cp_mom |
Moments of the settling velocities times the heat capacity.
Definition at line 50 of file particles.f90.
real*8, dimension(:,:), allocatable particles_module::set_cp_rhop_mom |
Moments of the settling velocities times the heat cap times the densities.
Definition at line 47 of file particles.f90.
real*8, dimension(:,:), allocatable particles_module::set_mom |
Moments of the settling velocities.
Definition at line 35 of file particles.f90.
real*8, dimension(:,:), allocatable particles_module::set_rhop_mom |
Moments of the densities times the settling velocities.
Definition at line 41 of file particles.f90.
character*10 particles_module::settling_model |
Settling model:
.
Definition at line 82 of file particles.f90.
real*8 particles_module::shape_factor |
shape factor for settling velocity (Pfeiffer)
Definition at line 59 of file particles.f90.
real*8, dimension(:), allocatable particles_module::solid_mass_fraction |
mass fraction of the particle phases with respect to the mixture
Definition at line 26 of file particles.f90.
real*8, dimension(:), allocatable particles_module::solid_partial_mass_fraction |
mass fraction of the particle phases with respect to the total solid
Definition at line 20 of file particles.f90.
real*8, dimension(:), allocatable particles_module::solid_partial_volume_fraction |
volume fraction of the particle phases with respect to the total solid
Definition at line 23 of file particles.f90.
real*8, dimension(:), allocatable particles_module::solid_volume_fraction |
volume fraction of the particle phases with respect to the mixture
Definition at line 29 of file particles.f90.