Sigma Relative Velocity Ao2008¶
Particula Index / Particula / Dynamics / Coagulation / Turbulent Dns Kernel / Sigma Relative Velocity Ao2008
Auto-generated documentation for particula.dynamics.coagulation.turbulent_dns_kernel.sigma_relative_velocity_ao2008 module.
VelocityCorrelationTerms¶
Show source in sigma_relative_velocity_ao2008.py:37
Parameters from computing velocity correlation terms.
Signature¶
class VelocityCorrelationTerms(NamedTuple): ...
_compute_cross_correlation_velocity¶
Show source in sigma_relative_velocity_ao2008.py:256
Compute cross-correlation of fluctuating velocities for two droplets.
This function calculates the cross-correlation of the fluctuating velocities of two droplets using the following equation:
Where the equation is
- ⟨v'¹ v'²⟩ = (u'² f₂(R) / (τ_p1 τ_p2)) *
[b₁ d₁ Φ(c₁, e₁) - b₁ d₂ Φ(c₁, e₂) - b₂ d₁ Φ(c₂, e₁) + b₂ d₂ Φ(c₂, e₂)]
- v'¹, v'² are the fluctuating velocities for droplets 1 and 2.
- u' (fluid_rms_velocity) : Fluid RMS fluctuation velocity [m/s].
- τ_p1, τ_p2 : Inertia timescales of droplets 1 and 2 [s].
- f₂(R) : Longitudinal velocity correlation function.
- Φ(c, e) : Function Φ computed using
get_phi_ao2008
.
Arguments¶
- fluid_rms_velocity : Fluid RMS fluctuation velocity [m/s].
- collisional_radius : Distance between two colliding droplets [m].
- particle_inertia_time : Inertia timescale of droplet 1 [s].
- particle_velocity : Droplet velocity [m/s].
- taylor_microscale : Taylor microscale [m].
- eulerian_integral_length : Eulerian integral length scale [m].
- velocity_correlation_terms : Velocity correlation coefficients [-].
Returns¶
- Cross-correlation velocity [m²/s²].
Examples¶
import numpy as np
ccv = _compute_cross_correlation_velocity(
fluid_rms_velocity=0.3,
collisional_radius=np.array([1e-4, 2e-4]),
particle_inertia_time=np.array([1.0, 1.2]),
particle_velocity=np.array([0.1, 0.2]),
taylor_microscale=0.01,
eulerian_integral_length=0.1,
velocity_correlation_terms=VelocityCorrelationTerms(
b1=0.1, b2=0.2, d1=0.3, d2=0.4, c1=0.5, c2=0.6, e1=0.7, e2=0.8
)
)
# Output: array([...])
References¶
- Ayala, O. et al. (2008). Effects of turbulence on the geometric collision rate of sedimenting droplets. Part 2. Theory and parameterization. New Journal of Physics, 10. https://doi.org/10.1088/1367-2630/10/7/075016
Signature¶
@validate_inputs(
{
"fluid_rms_velocity": "positive",
"collisional_radius": "positive",
"particle_inertia_time": "positive",
"particle_velocity": "positive",
"taylor_microscale": "positive",
"eulerian_integral_length": "positive",
}
)
def _compute_cross_correlation_velocity(
fluid_rms_velocity: float,
collisional_radius: Union[float, NDArray[np.float64]],
particle_inertia_time: Union[float, NDArray[np.float64]],
particle_velocity: Union[float, NDArray[np.float64]],
taylor_microscale: float,
eulerian_integral_length: float,
velocity_correlation_terms: VelocityCorrelationTerms,
) -> Union[float, NDArray[np.float64]]: ...
See also¶
_compute_rms_fluctuation_velocity¶
Show source in sigma_relative_velocity_ao2008.py:151
Compute RMS fluctuation velocity for the k-th droplet.
This function calculates the square of the RMS fluctuation velocity for the k-th droplet using the following equation:
Where the equation is:
- ⟨(v'ᵏ)²⟩ = (u'² / τ_pk) * [b₁ d₁ Ψ(c₁, e₁) - b₁ d₂ Ψ(c₁, e₂)
- b₂ d₁ Ψ(c₂, e₁) + b₂ d₂ Ψ(c₂, e₂)]
- v'ᵏ is the fluctuating velocity for droplet k.
- u' (fluid_rms_velocity) : Fluid RMS fluctuation velocity [m/s].
- τ_pk (particle_inertia_time) : Inertia timescale of the droplet k [s].
- Ψ(c, e) : Function Ψ computed using
get_psi_ao2008
.
Arguments¶
- fluid_rms_velocity : Fluid RMS fluctuation velocity [m/s].
- particle_inertia_time : Inertia timescale of the droplet k [s].
- particle_velocity : Droplet velocity [m/s].
- velocity_correlation_terms : Velocity correlation coefficients [-].
Returns¶
- RMS fluctuation velocity [m²/s²].
Examples¶
import numpy as np
rms_fluct = _compute_rms_fluctuation_velocity(
fluid_rms_velocity=0.3,
particle_inertia_time=np.array([1.0, 1.2]),
particle_velocity=np.array([0.1, 0.2]),
velocity_correlation_terms=VelocityCorrelationTerms(
b1=0.1, b2=0.2, d1=0.3, d2=0.4, c1=0.5, c2=0.6, e1=0.7, e2=0.8
)
)
# Output: array([...])
References¶
- Ayala, O. et al. (2008). Effects of turbulence on the geometric collision rate of sedimenting droplets. Part 2. Theory and parameterization. New Journal of Physics, 10. https://doi.org/10.1088/1367-2630/10/7/075016
Signature¶
@validate_inputs({"fluid_rms_velocity": "positive", "particle_inertia_time": "positive"})
def _compute_rms_fluctuation_velocity(
fluid_rms_velocity: float,
particle_inertia_time: Union[float, NDArray[np.float64]],
particle_velocity: Union[float, NDArray[np.float64]],
velocity_correlation_terms: VelocityCorrelationTerms,
) -> Union[float, NDArray[np.float64]]: ...
See also¶
get_relative_velocity_variance¶
Show source in sigma_relative_velocity_ao2008.py:50
Compute the variance of particle relative-velocity fluctuations.
This function calculates the variance of particle relative-velocity fluctuations using the following equation:
Where the equation is:
- σ² = ⟨(v'²)⟩₁ + ⟨(v'²)⟩₂ - 2⟨v'¹ v'²⟩
- v'¹, v'² are the fluctuating velocities for droplets 1 and 2.
Arguments¶
- fluid_rms_velocity : Fluid RMS fluctuation velocity [m/s].
- collisional_radius : Distance between two colliding droplets [m].
- particle_inertia_time : Inertia timescale of droplet 1 [s].
- particle_velocity : Droplet velocity [m/s].
- taylor_microscale : Taylor microscale [m].
- eulerian_integral_length : Eulerian integral length scale [m].
- lagrangian_integral_time : Lagrangian integral time scale [s].
Returns¶
- σ² : Variance of the particle relative-velocity fluctuation [m²/s²].
Examples¶
import numpy as np
sigma_sq = get_relative_velocity_variance(
fluid_rms_velocity=0.3,
collisional_radius=np.array([1e-4, 2e-4]),
particle_inertia_time=np.array([1.0, 1.2]),
particle_velocity=np.array([0.1, 0.2]),
taylor_microscale=0.01,
eulerian_integral_length=0.1,
lagrangian_integral_time=0.5,
lagrangian_taylor_microscale_time=0.05
)
# Output: array([...])
References¶
- Ayala, O. et al. (2008). Effects of turbulence on the geometric collision rate of sedimenting droplets. Part 2. Theory and parameterization. New Journal of Physics, 10. https://doi.org/10.1088/1367-2630/10/7/075016
Signature¶
@validate_inputs(
{
"fluid_rms_velocity": "positive",
"collisional_radius": "positive",
"particle_inertia_time": "positive",
"particle_velocity": "positive",
}
)
def get_relative_velocity_variance(
fluid_rms_velocity: float,
collisional_radius: NDArray[np.float64],
particle_inertia_time: NDArray[np.float64],
particle_velocity: NDArray[np.float64],
taylor_microscale: float,
eulerian_integral_length: float,
lagrangian_integral_time: float,
lagrangian_taylor_microscale_time: float,
) -> Union[float, NDArray[np.float64]]: ...