Types and Functions

ThermofluidQuantities.DensityMethod
Density(T::Temperature,p::Pressure[;gas=Air])

Compute the density from the temperature T and pressure p, using

$\rho = p/(\rho T)$

Can also switch the order of the arguments.

source
ThermofluidQuantities.DensityMethod
Density(ρ0::StagnationDensity,M::MachNumber,Isentropic[;gas=Air])

Compute the density corresponding to stagnation density ρ0 and Mach number M, using

$\rho = \rho_0 \left( 1 + \dfrac{\gamma-1}{2} M^2\right)^{-1/(\gamma-1)}$

source
ThermofluidQuantities.DensityRatioMethod
DensityRatio(M1::MachNumber,NormalShock[;gas=Air])

Return the ratio of densities after and before a normal shock, given the Mach number M1 before the shock, using the equation

$\dfrac{\rho_2}{\rho_1} = \dfrac{(\gamma+1)M_1^2}{2+(\gamma-1)M_1^2}$

source
ThermofluidQuantities.DensityRatioMethod
DensityRatio(pr::PressureRatio,Isentropic[;gas=Air])

Compute the ratio of densities between two points connected isentropically, given the ratio of pressures pr between those two points, using the isentropic relation

$\dfrac{\rho_2}{\rho_1} = \left(\dfrac{p_2}{p_1}\right)^{1/\gamma}$

source
ThermofluidQuantities.DensityRatioMethod
DensityRatio(Tr::TemperatureRatio,Isentropic[;gas=Air])

Compute the ratio of densities between two points connected isentropically, given the ratio of temperatures Tr between those two points, using the isentropic relation

$\dfrac{\rho_2}{\rho_1} = \left(\dfrac{T_2}{T_1}\right)^{1/(\gamma-1)}$

source
ThermofluidQuantities.EntropyMethod
Entropy(s1::Entropy,M1::MachNumber,NormalShock[;gas=Air])

Return the entropy after a normal shock, given the entropy before the shock s1 and the Mach number before the shock M1, using the equation

$s_2 - s_1 = c_v \ln \left(1 + \dfrac{2\gamma}{\gamma+1}(M_1^2-1)\right) - c_p \ln \left( \dfrac{(\gamma+1)M_1^2}{2+(\gamma-1)M_1^2}\right)$

source
ThermofluidQuantities.HeatFluxMethod
HeatFlux(T01::StagnationTemperature,T02::StagnationTemperature[;gas=Air])

Compute the heat flux required to change stagnation temperature T01 to stagnation temperature T02, using

$q = c_p(T_{02} - T_{01})$

source
ThermofluidQuantities.MachNumberMethod
MachNumber(A_over_Astar::AreaRatio,Isentropic[;gas=Air]) -> MachNumber, MachNumber

Compute the subsonic and supersonic Mach numbers for the given ratio of area to sonic area A_over_Astar by solving the equation

$\dfrac{A}{A^*} = \dfrac{1}{M} \left( \dfrac{2}{\gamma+1} \left( 1 + \dfrac{\gamma-1}{2}M^2\right)\right)^{(\gamma+1)/2/(\gamma-1)}$

Can omit the Isentropic argument.

source
ThermofluidQuantities.MachNumberMethod
MachNumber(ρ_over_ρ0::DensityRatio,Isentropic[;gas=Air])

Compute the Mach number corresponding to the ratio of density to stagnation density ρ_over_ρ0, using

$M = \left( \dfrac{2}{\gamma-1}\left(\dfrac{1}{(\rho/\rho_0)^{(\gamma-1)}} - 1 \right)\right)^{1/2}$

Can also supply the arguments for ρ and ρ0 separately, MachNumber(ρ,ρ0)

source
ThermofluidQuantities.MachNumberMethod
MachNumber(fL_over_D::FLOverD,FannoFlow[;gas=Air]) -> MachNumber, MachNumber

Compute the subsonic and supersonic Mach numbers, given the ratio friction factor times distance to sonic point, divided by diameter, fL_over_D, by solving the equation

$\dfrac{fL^*}{D} = \dfrac{1-M^2}{\gamma M^2} + \dfrac{\gamma+1}{2\gamma} \ln \left(\dfrac{(\gamma+1)M^2}{2+(\gamma-1)M^2}\right)$

source
ThermofluidQuantities.MachNumberMethod
MachNumber(M1::MachNumber,A1::AreaRatio,A2::AreaRatio,Isentropic[;gas=Air]) -> MachNumber, MachNumber

Compute the subsonic and supersonic Mach numbers at location 2 with area A2, when location 1 has Mach number M1 and area A1.

source
ThermofluidQuantities.MachNumberMethod
MachNumber(M1::MachNumber,NormalShock[;gas=Air])

Return the Mach number after a normal shock, given the Mach number M1 before the shock, using the equation

$M_2^2 = \dfrac{2 + (\gamma-1)M_1^2}{2\gamma M_1^2 - (\gamma-1)}$

source
ThermofluidQuantities.MachNumberMethod
MachNumber(p_over_p0::PressureRatio,Isentropic[;gas=Air])

Compute the Mach number corresponding to the ratio of pressure to stagnation pressure p_over_p0, using

$M = \left( \dfrac{2}{\gamma-1}\left(\dfrac{1}{(p/p_0)^{(\gamma-1)/\gamma}} - 1 \right)\right)^{1/2}$

Can also supply the arguments for p and p0 separately, MachNumber(p,p0)

source
ThermofluidQuantities.MachNumberMethod
MachNumber(T_over_Tstar::TemperatureRatio,FannoFlow[;gas=PerfectGas])

Compute the Mach number for a given ratio of temperature to the sonic temperature T_over_Tstar in Fanno flow, from the solution of

$\dfrac{T}{T^*} = \dfrac{\gamma+1}{2+(\gamma-1)M^2}$

source
ThermofluidQuantities.MachNumberMethod
MachNumber(T_over_T0::TemperatureRatio,Isentropic[;gas=Air])

Compute the Mach number corresponding to the ratio of temperature to stagnation temperature T_over_T0, using

$M = \left( \dfrac{2}{\gamma-1}\left(\dfrac{1}{T/T_0} - 1 \right)\right)^{1/2}$

Can also omit the Isentropic argument, MachNumber(T_over_T0). Can also supply the arguments for T and T0 separately, MachNumber(T,T0)

source
ThermofluidQuantities.PressureMethod
Pressure(T::Temperature,ρ::Density[;gas=Air])

Compute the pressure from the temperature T and density ρ, using

$p = \rho R T$

Can also switch the order of the arguments.

source
ThermofluidQuantities.PressureMethod
Pressure(p0::StagnationPressure,M::MachNumber,Isentropic[;gas=Air])

Compute the pressure corresponding to stagnation pressure p0 and Mach number M, using

$p = p_0 \left( 1 + \dfrac{\gamma-1}{2} M^2\right)^{-\gamma/(\gamma-1)}$

Can also omit the Isentropic argument, Temperature(T0,M).

source
ThermofluidQuantities.PressureRatioMethod
PressureRatio(dr::DensityRatio,Isentropic[;gas=Air])

Compute the ratio of pressures between two points connected isentropically, given the ratio of densities dr between those two points, using the isentropic relation

$\dfrac{p_2}{p_1} = \left(\dfrac{\rho_2}{\rho_1}\right)^{\gamma}$

source
ThermofluidQuantities.PressureRatioMethod
PressureRatio(M1::MachNumber,NormalShock[;gas=Air])

Return the ratio of pressures after and before a normal shock, given the Mach number M1 before the shock, using the equation

$\dfrac{p_2}{p_1} = 1 + 2\dfrac{\gamma}{\gamma+1}(M_1^2-1)$

source
ThermofluidQuantities.PressureRatioMethod
PressureRatio(Tr::TemperatureRatio,Isentropic[;gas=Air])

Compute the ratio of temperatures between two points connected isentropically, given the ratio of densities dr between those two points, using the isentropic relation

$\dfrac{p_2}{p_1} = \left(\dfrac{T_2}{T_1}\right)^{\gamma/(\gamma-1)}$

source
ThermofluidQuantities.SoundSpeedMethod
SoundSpeed(T::Temperature[,γ::SpecificHeatRatio=1.4][,R::GasConstant=287])

Compute the speed of sound of a perfect gas, based on the absolute temperature T. The ratio of specific heats γ and R can also be supplied, but default to the values for air at standard conditions (1.4 and 287 J/kg.K, respectively).

source
ThermofluidQuantities.StagnationDensityMethod
StagnationDensity(ρ::Density,M::MachNumber,Isentropic[;gas=Air])

Compute the stagnation density corresponding to density ρ and Mach number M, using

$\rho_0 = \rho \left( 1 + \dfrac{\gamma-1}{2} M^2\right)^{1/(\gamma-1)}$

source
ThermofluidQuantities.StagnationDensityMethod
StagnationDensity(T0::StagnationTemperature,p0::StagnationPressure[;gas=Air])

Compute the stagnation density from the stagnation temperature T0 and stagnation pressure p0, using

$\rho_0 = p_0/(\rho T_0)$

Can also switch the order of the arguments.

source
ThermofluidQuantities.StagnationPressureMethod
StagnationPressure(p::Pressure,M::MachNumber,Isentropic[;gas=Air])

Compute the stagnation pressure corresponding to pressure p and Mach number M, using

$p_0 = p \left( 1 + \dfrac{\gamma-1}{2} M^2\right)^{\gamma/(\gamma-1)}$

Can also omit the Isentropic argument, StagnationTemperature(T,M).

source
ThermofluidQuantities.StagnationPressureMethod
StagnationPressure(T0::StagnationTemperature,ρ0::StagnationDensity[;gas=Air])

Compute the stagnation pressure from the stagnation temperature T0 and stagnation density ρ0, using

$p_0 = \rho R T_0$

Can also switch the order of the arguments.

source
ThermofluidQuantities.StagnationPressureRatioMethod
StagnationPressureRatio(M1::MachNumber,NormalShock[;gas=Air])

Return the ratio of stagnation pressures after and before a normal shock, given the Mach number M1 before the shock, using the equation

$\dfrac{p_{02}}{p_{01}} = \left(\dfrac{(\gamma+1)M_1^2}{2+(\gamma-1)M_1^2}\right)^{\gamma/(\gamma-1)} \left(\dfrac{\gamma+1}{2\gamma M_1^2 - (\gamma-1)}\right)^{1/(\gamma-1)}$

source
ThermofluidQuantities.StagnationTemperatureMethod
StagnationTemperature(ρ0::StagnationDensity,p0::StagnationPressure[;gas=Air])

Compute the stagnation temperature from the stagnation density ρ0 and stagnation pressure p0, using

$T_0 = p_0/(\rho_0 R)$

Can also switch the order of the arguments.

source
ThermofluidQuantities.StagnationTemperatureMethod
StagnationTemperature(T::Temperature,M::MachNumber,Isentropic[;gas=Air])

Compute the stagnation temperature corresponding to temperature T and Mach number M, using

$T_0 = T \left( 1 + \dfrac{\gamma-1}{2} M^2\right)$

Can also omit the Isentropic argument, StagnationTemperature(T,M).

source
ThermofluidQuantities.TemperatureMethod
Temperature(ρ::Density,p::Pressure[;gas=Air])

Compute the temperature from the density ρ and pressure p, using

$T = p/(\rho R)$

Can also switch the order of the arguments.

source
ThermofluidQuantities.TemperatureMethod
Temperature(T0::StagnationTemperature,M::MachNumber,Isentropic[;gas=Air])

Compute the temperature corresponding to stagnation temperature T0 and Mach number M, using

$T = T_0 \left( 1 + \dfrac{\gamma-1}{2} M^2\right)^{-1}$

Can also omit the Isentropic argument, Temperature(T0,M).

source
ThermofluidQuantities.TemperatureRatioMethod
TemperatureRatio(dr::DensityRatio,Isentropic[;gas=Air])

Compute the ratio of temperatures between two points connected isentropically, given the ratio of densities dr between those two points, using the isentropic relation

$\dfrac{T_2}{T_1} = \left(\dfrac{\rho_2}{\rho_1}\right)^{\gamma-1}$

source
ThermofluidQuantities.TemperatureRatioMethod
TemperatureRatio(M1::MachNumber,NormalShock[;gas=Air])

Return the ratio of temperatures after and before a normal shock, given the Mach number M1 before the shock, using the equation

$\dfrac{T_2}{T_1} = 1 + 2\dfrac{(\gamma-1)}{(\gamma+1)^2}\dfrac{(1+\gamma M_1^2)(M_1^2-1)}{M_1^2}$

source
ThermofluidQuantities.TemperatureRatioMethod
TemperatureRatio(pr::PressureRatio,Isentropic[;gas=Air])

Compute the ratio of temperatures between two points connected isentropically, given the ratio of pressures pr between those two points, using the isentropic relation

$\dfrac{T_2}{T_1} = \left(\dfrac{p_2}{p_1}\right)^{(\gamma-1)/\gamma}$

source
Gasdynamics1D.AOverAStarMethod
AOverAStar(M::MachNumber,Isentropic[;gas=Air]) -> AreaRatio

Compute the ratio of local area to the sonic area for the given Mach number M, using

$\dfrac{A}{A^*} = \dfrac{1}{M} \left( \dfrac{2}{\gamma+1} \left( 1 + \dfrac{\gamma-1}{2}M^2\right)\right)^{(\gamma+1)/2/(\gamma-1)}$

Can omit the Isentropic argument.

source
Gasdynamics1D.AStarMethod
AStar(A::Area,M::MachNumber,Isentropic[;gas=Air]) -> Area

Compute the sonic area for the given area A and Mach number M, using

$A^* = A M \left( \dfrac{2}{\gamma+1} \left( 1 + \dfrac{\gamma-1}{2}M^2\right)\right)^{-(\gamma+1)/2/(\gamma-1)}$

Can omit the Isentropic argument.

source
Gasdynamics1D.FLStarOverDMethod
FLStarOverD(M::MachNumber,FannoFlow[;gas=Air])

Compute the friction factor times distance to sonic point, divided by diameter, for the given Mach number M, using the equation

$\dfrac{fL^*}{D} = \dfrac{1-M^2}{\gamma M^2} + \dfrac{\gamma+1}{2\gamma} \ln \left(\dfrac{(\gamma+1)M^2}{2+(\gamma-1)M^2}\right)$

source
Gasdynamics1D.P0OverPMethod
    P0OverP(M::MachNumber,Isentropic[;gas=PerfectGas]) -> PressureRatio

Compute the ratio of stagnation pressure to pressure for Mach number M.

$\dfrac{p_0}{p} = \left(1 + \dfrac{\gamma-1}{2}M^2 \right)^{\gamma/(\gamma-1)}$

source
Gasdynamics1D.P0OverP0StarMethod
P0OverP0Star(M::MachNumber,FannoFlow[;gas=PerfectGas]) -> StagnationPressureRatio

Compute the ratio of stagnation pressure to the sonic stagnation pressure in Fanno flow, given Mach number M, from

$\dfrac{p_0}{p_0^*} = \dfrac{1}{M} \left(\dfrac{2+(\gamma-1)M^2}{\gamma+1}\right)^{(\gamma+1)/2/(\gamma-1)}$

source
Gasdynamics1D.P0OverP0StarMethod
P0OverP0Star(M::MachNumber,RayleighFlow[;gas=Air]) -> PressureRatio

Compute the ratio of stagnation pressure to sonic stagnation pressure for the given Mach number for Rayleigh flow, using the equation

$\dfrac{p_0}{p_0^*} = \dfrac{(\gamma+1)}{1+\gamma M^2} \left[ \left( \dfrac{1}{\gamma+1} \right)\left(2 + (\gamma-1) M^2\right)\right]^{\gamma/(\gamma-1)}$

source
Gasdynamics1D.POverPStarMethod
POverPStar(M::MachNumber,FannoFlow[;gas=PerfectGas]) -> PressureRatio

Compute the ratio of pressure to the sonic pressure in Fanno flow, given Mach number M, from

$\dfrac{p}{p^*} = \dfrac{1}{M} \left(\dfrac{1+\gamma}{2+(\gamma-1)M^2}\right)^{1/2}$

source
Gasdynamics1D.POverPStarMethod
POverPStar(M::MachNumber,RayleighFlow[;gas=Air]) -> PressureRatio

Compute the ratio of pressure to sonic pressure for the given Mach number for Rayleigh flow, using the equation

$\dfrac{p}{p^*} = \dfrac{\gamma+1}{1+\gamma M^2}$

source
Gasdynamics1D.SubsonicMachNumberMethod
SubsonicMachNumber(A_over_Astar::AreaRatio,Isentropic[;gas=Air]) -> MachNumber

Compute the subsonic Mach number for the given ratio of area to sonic area A_over_Astar by solving the equation

$\dfrac{A}{A^*} = \dfrac{1}{M} \left( \dfrac{2}{\gamma+1} \left( 1 + \dfrac{\gamma-1}{2}M^2\right)\right)^{(\gamma+1)/2/(\gamma-1)}$

source
Gasdynamics1D.SubsonicMachNumberMethod
SubsonicMachNumber(fL_over_D::FLOverD,FannoFlow[;gas=Air]) -> MachNumber

Compute the subsonic Mach number, given the ratio friction factor times distance to sonic point, divided by diameter, fL_over_D, by solving the equation

$\dfrac{fL^*}{D} = \dfrac{1-M^2}{\gamma M^2} + \dfrac{\gamma+1}{2\gamma} \ln \left(\dfrac{(\gamma+1)M^2}{2+(\gamma-1)M^2}\right)$

source
Gasdynamics1D.SupersonicMachNumberMethod
SupersonicMachNumber(A_over_Astar::AreaRatio,Isentropic[;gas=Air]) -> MachNumber

Compute the supersonic Mach number for the given ratio of area to sonic area A_over_Astar by solving the equation

$\dfrac{A}{A^*} = \dfrac{1}{M} \left( \dfrac{2}{\gamma+1} \left( 1 + \dfrac{\gamma-1}{2}M^2\right)\right)^{(\gamma+1)/2/(\gamma-1)}$

source
Gasdynamics1D.SupersonicMachNumberMethod
SupersonicMachNumber(fL_over_D::FLOverD,FannoFlow[;gas=Air]) -> MachNumber

Compute the supersonic Mach number, given the ratio friction factor times distance to sonic point, divided by diameter, fL_over_D, by solving the equation

$\dfrac{fL^*}{D} = \dfrac{1-M^2}{\gamma M^2} + \dfrac{\gamma+1}{2\gamma} \ln \left(\dfrac{(\gamma+1)M^2}{2+(\gamma-1)M^2}\right)$

source
Gasdynamics1D.T0OverTMethod
    T0OverT(M::MachNumber,Isentropic[;gas=PerfectGas]) -> TemperatureRatio

Compute the ratio of stagnation temperature to temperature for Mach number M.

$\dfrac{T_0}{T} = 1 + \dfrac{\gamma-1}{2}M^2$

Can also omit the Isentropic argument, T0OverT(M)

source
Gasdynamics1D.T0OverT0StarMethod
T0OverT0Star(M::MachNumber,RayleighFlow[;gas=Air]) -> TemperatureRatio

Compute the ratio of stagnation temperature to sonic stagnation temperature for the given Mach number for Rayleigh flow, using the equation

$\dfrac{T_0}{T_0^*} = \dfrac{(\gamma+1)M^2 (2 + (\gamma-1)M^2)}{(1+\gamma M^2)^2}$

source
Gasdynamics1D.TOverTStarMethod
TOverTStar(M::MachNumber,FannoFlow[;gas=PerfectGas]) -> TemperatureRatio

Compute the ratio of temperature to the sonic temperature in Fanno flow, given Mach number M, from

$\dfrac{T}{T^*} =\dfrac{\gamma+1}{2+(\gamma-1)M^2}$

source
Gasdynamics1D.TOverTStarMethod
TOverTStar(M::MachNumber,RayleighFlow[;gas=Air]) -> TemperatureRatio

Compute the ratio of temperature to sonic temperature for the given Mach number for Rayleigh flow, using the equation

$\dfrac{T}{T^*} = \dfrac{(\gamma+1)^2 M^2}{(1+\gamma M^2)^2}$

source
Gasdynamics1D.VOverVStarMethod
VOverVStar(M::MachNumber,RayleighFlow[;gas=Air]) -> VelocityRatio

Compute the ratio of velocity to sonic velocity for the given Mach number for Rayleigh flow, using the equation

$\dfrac{V}{V^*} = \dfrac{(\gamma+1)M^2}{1+\gamma M^2}$

source
Gasdynamics1D.ρ0OverρMethod
    ρ0Overρ(M::MachNumber,Isentropic[;gas=PerfectGas]) -> DensityRatio

Compute the ratio of stagnation density to density for Mach number M.

$\dfrac{\rho_0}{\rho} = \left(1 + \dfrac{\gamma-1}{2}M^2 \right)^{1/(\gamma-1)}$

source
Gasdynamics1D.ρOverρStarMethod
ρOverρStar(M::MachNumber,FannoFlow[;gas=PerfectGas]) -> DensityRatio

Compute the ratio of density to the sonic density in Fanno flow, given Mach number M, from

$\dfrac{\rho}{\rho^*} = \dfrac{1}{M} \left(\dfrac{2+(\gamma-1)M^2}{\gamma+1}\right)^{1/2}$

source
Gasdynamics1D.ρOverρStarMethod
ρOverρStar(M::MachNumber,RayleighFlow[;gas=Air]) -> DensityRatio

Compute the ratio of density to sonic density for the given Mach number for Rayleigh flow, using the equation

$\dfrac{\rho}{\rho^*} = \dfrac{1+\gamma M^2}{(\gamma+1)M^2}$

source