PLUME-MoM-TSM
1.0
VolcanicPlumeModel
|
Particles module. More...
Functions/Subroutines | |
subroutine | allocate_particles |
Particles variables allocation. More... | |
subroutine | deallocate_particles |
real(wp) function | particles_settling_velocity (diam, rhop, shape_fact) |
Settling velocity. More... | |
real(wp) function | particles_heat_capacity (i_part, diam) |
Heat capacity. More... | |
real(wp) function | particles_density (i_part, phi) |
Particle density. More... | |
real(wp) function | particles_shape (i_part, phi) |
Particle shape factor. More... | |
real(wp) function | particles_beta (temp, visc, diam_i, diam_j, rho_i, rho_j, Vs_i, Vs_j, lw_vf, ice_vf, solid_mf) |
Particles aggregation. More... | |
real(wp) function | aggregation_kernel (diam_i, rho_i, Vs_i, diam_j, rho_j, Vs_j, lw_vf, ice_vf, solid_mf, temp, visc) |
Aggregation kernel. More... | |
real(wp) function | collision_kernel (diam_i, rho_i, Vs_i, diam_j, rho_j, Vs_j, temp, visc) |
Collision kernel. More... | |
real(wp) function | coalescence_efficiency (diam_i, rho_i, Vs_i, diam_j, rho_j, Vs_j, lw_vf, ice_vf, solid_mf) |
Coalescence efficiency. More... | |
subroutine | eval_particles_moments |
Particles moments computation. More... | |
subroutine | eval_quad_values |
Quadrature values computation. More... | |
subroutine | init_quadrature_points |
Quadrature initialization. More... | |
subroutine | phifromm |
Compute size from mass. More... | |
subroutine | update_aggregation (temp, visc, lw_vf, ice_vf, solid_mf) |
Aggregation terms. More... | |
subroutine | init_aggregation |
Aggregation initialization. More... | |
Variables | |
integer | n_part |
number of particle phases More... | |
real(wp), dimension(:), allocatable | solid_partial_mass_fraction |
mass fraction of the particle phases with respect to the total solid More... | |
real(wp), dimension(:), allocatable | solid_partial_mass_fraction0 |
init mass fraction of the particle phases with respect to the total solid More... | |
real(wp), dimension(:), allocatable | solid_partial_volume_fraction |
volume fraction of the particle phases with respect to the total solid More... | |
real(wp), dimension(:), allocatable | solid_mass_fraction |
mass fraction of the particle phases with respect to the mixture More... | |
real(wp), dimension(:), allocatable | solid_mass_fraction0 |
initial mass fraction of the particle phases with respect to the mixture More... | |
real(wp), dimension(:), allocatable | solid_volume_fraction |
volume fraction of the particle phases with respect to the mixture More... | |
real(wp), dimension(:,:), allocatable | bin_partial_mass_fraction |
mass fraction of the bins of particle with respect to the total solid More... | |
real(wp), dimension(:,:), allocatable | particle_loss_rate |
rate of particles lost from the plume in the integration steps ( kg s^-1) More... | |
real(wp), dimension(:,:), allocatable | cum_particle_loss_rate |
cumulative rate of particles lost up to the integration height ( kg s^-1) More... | |
real(wp), dimension(:,:,:), allocatable | mom |
Moments of the particles diameter. More... | |
real(wp), dimension(:,:,:), allocatable | mom0 |
Initial moments of the particles diameter. More... | |
real(wp), dimension(:,:,:), allocatable | set_mom |
Moments of the settling velocities. More... | |
real(wp), dimension(:,:,:), allocatable | set_cp_mom |
Moments of the settling velocities times the heat capacity. More... | |
real(wp), dimension(:,:,:), allocatable | birth_mom |
Term accounting for the birth of aggregates in the moments equations. More... | |
real(wp), dimension(:,:,:), allocatable | death_mom |
Term accounting for the loss of particles because of aggregation. More... | |
real(wp), dimension(:), allocatable | shape_factor |
shape factor for settling velocity (Pfeiffer) More... | |
real(wp), dimension(:), allocatable | phi1 |
First diameter for the density function. More... | |
real(wp), dimension(:), allocatable | rho1 |
Density at phi=phi1. More... | |
real(wp), dimension(:), allocatable | shape1 |
Shape factor at phi=phi1. More... | |
real(wp), dimension(:), allocatable | phi2 |
Second diameter for the density function. More... | |
real(wp), dimension(:), allocatable | rho2 |
Density at phi=phi2. More... | |
real(wp), dimension(:), allocatable | shape2 |
Shape factor at phi=phi2. More... | |
real(wp), dimension(:), allocatable | cp_part |
Heat capacity of particle phases. More... | |
real(wp) | cpsolid |
Average heat capacity of particles. More... | |
real(wp) | particles_beta0 |
character *10 | settling_model |
Settling model: . More... | |
character(len=20) | distribution |
Ditribution of the particles: . More... | |
logical, dimension(:), allocatable | aggregation_array |
Flag for the aggregation: . More... | |
real(wp), dimension(:), allocatable | aggregate_porosity |
Array for porosity volume fraction of aggregates. More... | |
character(len=20) | aggregation_model |
Aggregation kernel model: . More... | |
integer, dimension(:), allocatable | aggr_idx |
Index defining the couple aggregated-non aggregated. More... | |
real(wp), dimension(:,:,:), allocatable | m_quad |
Abscissa of quadrature formulas. More... | |
real(wp), dimension(:,:,:), allocatable | w_quad |
Weights of quadrature formulas. More... | |
real(wp), dimension(:,:,:), allocatable | phi_quad |
Particle size (phi-scale) at quadrature points. More... | |
real(wp), dimension(:,:,:), allocatable | diam_quad |
Particle diameters (meters) at quadrature points. More... | |
real(wp), dimension(:,:,:), allocatable | vol_quad |
Particle volumes at quadrature points. More... | |
real(wp), dimension(:,:,:), allocatable | rho_quad |
Particle densities at quadrature points. More... | |
real(wp), dimension(:,:,:), allocatable | shape_quad |
Particle densities at quadrature points. More... | |
real(wp), dimension(:,:,:), allocatable | set_vel_quad |
Particle settling velocities at quadrature points. More... | |
real(wp), dimension(:,:,:), allocatable | cp_quad |
Particle heat capacities at quadrature points. More... | |
real(wp), dimension(:,:,:), allocatable | f_quad |
Values of linear reconstructions at quadrature points. More... | |
real(wp) | t_part |
Particle temperature for aggregation (Costa model) More... | |
real(wp), dimension(:), allocatable | phil |
left boundaries of the sections in phi-scale More... | |
real(wp), dimension(:), allocatable | phir |
right boundaries of the sections in phi-scale More... | |
real(wp), dimension(:,:), allocatable | m |
boundaries of the sections in mass scale (kg) More... | |
logical, dimension(:,:,:,:,:), allocatable | q_flag |
logical defining if particles ip/is+jp/js aggregates on section ks More... | |
real(wp), dimension(:,:,:,:,:,:), allocatable | kernel_aggr |
aggregation kernel computed for ip/is+jp/js More... | |
real(wp), dimension(:,:,:,:,:,:,:), allocatable | a |
real(wp), dimension(:,:,:,:,:,:,:), allocatable | wij |
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.
real(wp) function particles_module::aggregation_kernel | ( | real(wp), intent(in) | diam_i, |
real(wp), intent(in) | rho_i, | ||
real(wp), intent(in) | Vs_i, | ||
real(wp), intent(in) | diam_j, | ||
real(wp), intent(in) | rho_j, | ||
real(wp), intent(in) | Vs_j, | ||
real(wp), intent(in) | lw_vf, | ||
real(wp), intent(in) | ice_vf, | ||
real(wp), intent(in) | solid_mf, | ||
real(wp), intent(in) | temp, | ||
real(wp), intent(in) | visc | ||
) |
Aggregation kernel.
This function evaluates the aggregation kernel, using the expression given in Textor et al. 2006.
[in] | diam_i | first particle diameter (m) |
[in] | rho_i | first particle density (kg/m3) |
[in] | Vs_i | first particle settling velocity (m/s) |
[in] | diam_j | second particle diameter (m) |
[in] | rho_j | second particle density (kg/m3) |
[in] | Vs_j | second particle settling velocity (m/s) |
[in] | lw_vf | liquid water mass fraction |
[in] | ice_vf | ice mass fraction |
[in] | solid_mf | ice mass fraction |
[in] | temp | mixture temperature (K) |
[in] | visc | air kinematic viscosity (m2s-1) |
Definition at line 821 of file particles.f90.
subroutine particles_module::allocate_particles | ( | ) |
Particles variables allocation.
This subroutine allocate 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 198 of file particles.f90.
real(wp) function particles_module::coalescence_efficiency | ( | real(wp), intent(in) | diam_i, |
real(wp), intent(in) | rho_i, | ||
real(wp), intent(in) | Vs_i, | ||
real(wp), intent(in) | diam_j, | ||
real(wp), intent(in) | rho_j, | ||
real(wp), intent(in) | Vs_j, | ||
real(wp), intent(in) | lw_vf, | ||
real(wp), intent(in) | ice_vf, | ||
real(wp), intent(in) | solid_mf | ||
) |
Coalescence efficiency.
[in] | diam_i | first particle diameter (m) |
[in] | rho_i | first particle density (kg/m3) |
[in] | Vs_i | first particle settling velocity (m/s) |
[in] | diam_j | second particle diameter (m) |
[in] | rho_j | second particle density (kg/m3) |
[in] | Vs_j | second particle settling velocity (m/s) |
[in] | lw_vf | liquid water mass fraction |
[in] | ice_vf | ice mass fraction |
[in] | solid_mf | ice mass fraction |
particle Stokes number
Critical Stokes number
Efficiency exponent
Definition at line 995 of file particles.f90.
real(wp) function particles_module::collision_kernel | ( | real(wp), intent(in) | diam_i, |
real(wp), intent(in) | rho_i, | ||
real(wp), intent(in) | Vs_i, | ||
real(wp), intent(in) | diam_j, | ||
real(wp), intent(in) | rho_j, | ||
real(wp), intent(in) | Vs_j, | ||
real(wp), intent(in) | temp, | ||
real(wp), intent(in) | visc | ||
) |
Collision kernel.
[in] | diam_i | first particle diameter (m) |
[in] | rho_i | first particle density (kg/m3) |
[in] | Vs_i | first particle settling velocity (m/s) |
[in] | diam_j | second particle diameter (m) |
[in] | rho_j | second particle density (kg/m3) |
[in] | Vs_j | second particle settling velocity (m/s) |
[in] | temp | mixture temperature (K) |
[in] | visc | air kinematic viscosity (m2s-1) |
Brownian motion collisions kernel (m3 s-1)
Laminar and turbulent fluid shear collisions kernel (m3 s-1)
Differential sedimentation kernel (m3 s-1)
Gravitational collision efficiency
Rate of dissipation of turbulent kinetic energy
Fluid shear
Definition at line 898 of file particles.f90.
subroutine particles_module::deallocate_particles | ( | ) |
subroutine particles_module::eval_particles_moments | ( | ) |
Particles moments computation.
This subroutine compute the moments of the particles properties (density, heat capacity and settling velocity) using the quadrature formulas.
Definition at line 1100 of file particles.f90.
subroutine particles_module::eval_quad_values | ( | ) |
Quadrature values computation.
This subroutine compute the values of the linear reconstructions at the quadrature points. These values changes with the moments and are used in the quadrature schemes to compute the integral for the additional moments.
Definition at line 1179 of file particles.f90.
subroutine particles_module::init_aggregation | ( | ) |
Aggregation initialization.
This subroutine compute the terms needed in the quadrature formulas used for the integral defining the aggregation.
Definition at line 1638 of file particles.f90.
subroutine particles_module::init_quadrature_points | ( | ) |
Quadrature initialization.
This subroutine compute the abscissas and weights for the quadrature schemes used by the method of moments. The values computed here are constant through the entire simulation and depend only on the sections.
Definition at line 1227 of file particles.f90.
real(wp) function particles_module::particles_beta | ( | real(wp), intent(in) | temp, |
real(wp), intent(in) | visc, | ||
real(wp), intent(in) | diam_i, | ||
real(wp), intent(in) | diam_j, | ||
real(wp), intent(in), optional | rho_i, | ||
real(wp), intent(in), optional | rho_j, | ||
real(wp), intent(in), optional | Vs_i, | ||
real(wp), intent(in), optional | Vs_j, | ||
real(wp), intent(in), optional | lw_vf, | ||
real(wp), intent(in), optional | ice_vf, | ||
real(wp), intent(in), optional | solid_mf | ||
) |
Particles aggregation.
This function evaluates the aggregation kernel using several forumation:
[in] | diam_i | first particle diameter (m) |
[in] | diam_j | second particle diameter (m) |
[in] | rho_i | first particle density (kg/m3) |
[in] | rho_j | second particle density (kg/m3) |
[in] | Vs_i | first particle settling velocity |
[in] | Vs_j | second particle settling velocity |
[in] | lw_vf | liquid water mass fraction |
[in] | ice_vf | ice mass fraction |
[in] | solid_mf | ice mass fraction |
Definition at line 731 of file particles.f90.
real(wp) function particles_module::particles_density | ( | integer, intent(in) | i_part, |
real(wp), intent(in) | phi | ||
) |
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] | phi | particle size (Krumbein phi scale) |
Definition at line 607 of file particles.f90.
real(wp) function particles_module::particles_heat_capacity | ( | integer, intent(in) | i_part, |
real(wp), intent(in) | diam | ||
) |
Heat capacity.
This function evaluates the heat capacity of the particles given the size (diameter).
[in] | i_part | particle phase index |
[in] | diam | particle diameter (m) |
Definition at line 583 of file particles.f90.
real(wp) function particles_module::particles_settling_velocity | ( | real(wp), intent(in) | diam, |
real(wp), intent(in) | rhop, | ||
real(wp), intent(in) | shape_fact | ||
) |
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] | diam | particle diameter (m) |
[in] | rhop | particle density (kg/m3) |
[in] | shape_fact | shape factor |
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 355 of file particles.f90.
real(wp) function particles_module::particles_shape | ( | integer, intent(in) | i_part, |
real(wp), intent(in) | phi | ||
) |
Particle shape factor.
This function evaluates the shape factor of a particle given the size (phi).
[in] | i_part | particle phase index |
[in] | phi | particle size (Krumbein phi scale) |
Definition at line 662 of file particles.f90.
subroutine particles_module::phifromm | ( | ) |
Compute size from mass.
This subroutine compute the size in the Krumbein phi-scale from the mass at the quadrature points. A bisection method is used to invert the size-density relationships and find the density from the mass. Once the mass is computed, additional variables are computed at the quadrature points: diam (m), vol (m3) , heat capacity (J*K-1*kg-1), phi.
Definition at line 1291 of file particles.f90.
subroutine particles_module::update_aggregation | ( | real(wp), intent(in) | temp, |
real(wp), intent(in) | visc, | ||
real(wp), intent(in) | lw_vf, | ||
real(wp), intent(in) | ice_vf, | ||
real(wp), intent(in) | solid_mf | ||
) |
Aggregation terms.
This subroutine compute the particles birth and death terms appearing in the momentum transport equations. The aggregation kernel is called inside from this subroutine.
[in] | temp | temperature (K) |
[in] | visc | viscosity |
[in] | lw_vf | liquid water mass fraction |
[in] | ice_vf | ice mass fraction |
[in] | solid_mf | ice mass fraction |
Definition at line 1478 of file particles.f90.
real(wp), dimension(:,:,:,:,:,:,:), allocatable particles_module::a |
Definition at line 178 of file particles.f90.
integer, dimension(:), allocatable particles_module::aggr_idx |
Index defining the couple aggregated-non aggregated.
Definition at line 128 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::aggregate_porosity |
Array for porosity volume fraction of aggregates.
Definition at line 118 of file particles.f90.
logical, dimension(:), allocatable particles_module::aggregation_array |
Flag for the aggregation:
.
Definition at line 115 of file particles.f90.
character(len=20) particles_module::aggregation_model |
Aggregation kernel model:
.
Definition at line 125 of file particles.f90.
real(wp), dimension(:,:), allocatable particles_module::bin_partial_mass_fraction |
mass fraction of the bins of particle with respect to the total solid
Definition at line 44 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::birth_mom |
Term accounting for the birth of aggregates in the moments equations.
Definition at line 65 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::cp_part |
Heat capacity of particle phases.
Definition at line 92 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::cp_quad |
Particle heat capacities at quadrature points.
Definition at line 155 of file particles.f90.
real(wp) particles_module::cpsolid |
Average heat capacity of particles.
Definition at line 95 of file particles.f90.
real(wp), dimension(:,:), allocatable particles_module::cum_particle_loss_rate |
cumulative rate of particles lost up to the integration height ( kg s^-1)
Definition at line 50 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::death_mom |
Term accounting for the loss of particles because of aggregation.
Definition at line 68 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::diam_quad |
Particle diameters (meters) at quadrature points.
Definition at line 140 of file particles.f90.
character(len=20) particles_module::distribution |
Ditribution of the particles:
.
Definition at line 110 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::f_quad |
Values of linear reconstructions at quadrature points.
Definition at line 158 of file particles.f90.
real(wp), dimension(:,:,:,:,:,:), allocatable particles_module::kernel_aggr |
aggregation kernel computed for ip/is+jp/js
Definition at line 176 of file particles.f90.
real(wp), dimension(:,:), allocatable particles_module::m |
boundaries of the sections in mass scale (kg)
Definition at line 170 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::m_quad |
Abscissa of quadrature formulas.
Definition at line 131 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::mom |
Moments of the particles diameter.
Definition at line 53 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::mom0 |
Initial moments of the particles diameter.
Definition at line 56 of file particles.f90.
integer particles_module::n_part |
number of particle phases
Definition at line 23 of file particles.f90.
real(wp), dimension(:,:), allocatable particles_module::particle_loss_rate |
rate of particles lost from the plume in the integration steps ( kg s^-1)
Definition at line 47 of file particles.f90.
real(wp) particles_module::particles_beta0 |
Definition at line 98 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::phi1 |
First diameter for the density function.
Definition at line 74 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::phi2 |
Second diameter for the density function.
Definition at line 83 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::phi_quad |
Particle size (phi-scale) at quadrature points.
Definition at line 137 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::phil |
left boundaries of the sections in phi-scale
Definition at line 164 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::phir |
right boundaries of the sections in phi-scale
Definition at line 167 of file particles.f90.
logical, dimension(:,:,:,:,:), allocatable particles_module::q_flag |
logical defining if particles ip/is+jp/js aggregates on section ks
Definition at line 173 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::rho1 |
Density at phi=phi1.
Definition at line 77 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::rho2 |
Density at phi=phi2.
Definition at line 86 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::rho_quad |
Particle densities at quadrature points.
Definition at line 146 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::set_cp_mom |
Moments of the settling velocities times the heat capacity.
Definition at line 62 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::set_mom |
Moments of the settling velocities.
Definition at line 59 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::set_vel_quad |
Particle settling velocities at quadrature points.
Definition at line 152 of file particles.f90.
character*10 particles_module::settling_model |
Settling model:
.
Definition at line 104 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::shape1 |
Shape factor at phi=phi1.
Definition at line 80 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::shape2 |
Shape factor at phi=phi2.
Definition at line 89 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::shape_factor |
shape factor for settling velocity (Pfeiffer)
Definition at line 71 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::shape_quad |
Particle densities at quadrature points.
Definition at line 149 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::solid_mass_fraction |
mass fraction of the particle phases with respect to the mixture
Definition at line 35 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::solid_mass_fraction0 |
initial mass fraction of the particle phases with respect to the mixture
Definition at line 38 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::solid_partial_mass_fraction |
mass fraction of the particle phases with respect to the total solid
Definition at line 26 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::solid_partial_mass_fraction0 |
init mass fraction of the particle phases with respect to the total solid
Definition at line 29 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::solid_partial_volume_fraction |
volume fraction of the particle phases with respect to the total solid
Definition at line 32 of file particles.f90.
real(wp), dimension(:), allocatable particles_module::solid_volume_fraction |
volume fraction of the particle phases with respect to the mixture
Definition at line 41 of file particles.f90.
real(wp) particles_module::t_part |
Particle temperature for aggregation (Costa model)
Definition at line 161 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::vol_quad |
Particle volumes at quadrature points.
Definition at line 143 of file particles.f90.
real(wp), dimension(:,:,:), allocatable particles_module::w_quad |
Weights of quadrature formulas.
Definition at line 134 of file particles.f90.
real(wp), dimension(:,:,:,:,:,:,:), allocatable particles_module::wij |
Definition at line 180 of file particles.f90.