particula.dynamics.coagulation.turbulent_dns_kernel.sigma_relative_velocity_ao2008¶
sigma_relative_velocity_ao2008
¶
Compute RMS fluctuation velocities and cross-correlations.
This module provides functions to compute the square of the RMS fluctuation velocity and the cross-correlation of the fluctuating velocities for colliding droplets, based on the theory of turbulent DNS kernels.
VelocityCorrelationTerms
¶
Bases: NamedTuple
Parameters from computing velocity correlation terms.
get_relative_velocity_variance
¶
get_relative_velocity_variance(fluid_rms_velocity: float, collisional_radius: NDArray[float64], particle_inertia_time: NDArray[float64], particle_velocity: NDArray[float64], taylor_microscale: float, eulerian_integral_length: float, lagrangian_integral_time: float, lagrangian_taylor_microscale_time: float) -> Union[float, NDArray[np.float64]]
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.
Parameters:
-
- 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:
-
Union[float, NDArray[float64]]–- σ² : 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
Source code in particula/dynamics/coagulation/turbulent_dns_kernel/sigma_relative_velocity_ao2008.py
47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 | |