This is the first of two notebooks on the Euler equations. In this notebook, we discuss the equations and the structure of the exact solution to the Riemann problem. In Euler_approximate.ipynb, we will investigate approximate Riemann solvers.
In this chapter we study the system of hyperbolic PDEs that governs the motions of a compressible gas in the absence of viscosity. These consist of conservation laws for mass, momentum, and energy. Together, they are referred to as the compressible Euler equations, or simply the Euler equations. Our discussion here is fairly brief; for much more detail see (LeVeque, 2002) or (Toro, 2013).
We will use $\rho(x,t)$ to denote the fluid density and $u(x,t)$ for its velocity. Then the equation for conservation of mass is just the familiar continuity equation:
$$\rho_t + (\rho u)_x = 0.$$We discussed the conservation of momentum in a fluid already in Acoustics.ipynb. For convenience, we review the ideas here. The momentum density is given by the product of mass density and velocity, $\rho u$. The momentum flux has two components. First, the momentum is transported in the same way that the density is; this flux is given by the momentum density times the velocity: $\rho u^2$.
To understand the second term in the momentum flux, we must realize that a fluid is made up of many tiny molecules. The density and velocity we are modeling are average values over some small region of space. The individual molecules in that region are not all moving with exactly velocity $u$; that's just their average. Each molecule also has some additional random velocity component. These random velocities are what accounts for the pressure of the fluid, which we'll denote by $p$. These velocity components also lead to a net flux of momentum. Thus the momentum conservation equation is
$$(\rho u)_t + (\rho u^2 + p)_x = 0.$$This is very similar to the conservation of momentum equation in the shallow water equations, as discussed in Shallow_water.ipynb, in which case $hu$ is the momentum density and $\frac 1 2 gh^2$ is the hydrostatic pressure. For gas dynamics, a different expression must be used to compute the pressure $p$ from the conserved quantities. This relation is called the equation of state of the gas, as discussed further below.
The energy has two components: internal energy density $\rho e$ and kinetic energy density $\rho u^2/2$:
$$E = \rho e + \frac{1}{2}\rho u^2.$$Like the momentum flux, the energy flux involves both bulk transport ($Eu$) and transport due to pressure ($pu$):
$$E_t + (u(E+p)) = 0.$$You may have noticed that we have 4 unknowns (density, momentum, energy, and pressure) but only 3 conservation laws. We need one more relation to close the system. That relation, known as the equation of state, expresses how the pressure is related to the other quantities. We'll focus on the case of a polytropic ideal gas, for which
$$p = \rho e (\gamma-1).$$Here $\gamma$ is the ratio of specific heats, which for air is approximately 1.4.
We can write the three conservation laws as a single system $q_t + f(q)_x = 0$ by defining
\begin{align}
q & = \begin{pmatrix} \rho \\ \rho u \\ E\end{pmatrix}, &
f(q) & = \begin{pmatrix} \rho u \\ \rho u^2 + p \\ u(E+p)\end{pmatrix}.
\label{euler_conserved}
\end{align}
This is the 1D Euler system. In three dimensions, the equations are similar. We have two additional velocity components $v, w$, and their corresponding fluxes. Additionally, we have to account for fluxes in the $y$ and $z$ directions. We can write the full system as
$$ q_t + f(q)_x + g(q)_y + h(q)_z = 0$$
with
\begin{align}
q & = \begin{pmatrix} \rho \\ \rho u \\ \rho v \\ \rho w \\ E\end{pmatrix}, &
f(q) & = \begin{pmatrix} \rho u \\ \rho u^2 + p \\ \rho u v \\ \rho u w \\ u(E+p)\end{pmatrix} &
g(q) & = \begin{pmatrix} \rho v \\ \rho uv \\ \rho v^2 + p \\ \rho v w \\ v(E+p)\end{pmatrix} &
h(q) & = \begin{pmatrix} \rho w \\ \rho uw \\ \rho vw \\ \rho w^2 + p \\ w(E+p)\end{pmatrix}.
\end{align}
In the rest of the chapter we focus on the 1D system.
In our discussion of the structure of these equations, it is convenient to work with the primitive variables $(\rho, u, p)$ rather than the conserved variables. The quasilinear form is particularly simple in the primitive variables:
\begin{align} \label{euler_primitive} \begin{bmatrix} \rho \\ u \\ p \end{bmatrix}_t + \begin{bmatrix} u & \rho & 0 \\ 0 & u & 1/\rho \\ 0 & \gamma \rho & u \end{bmatrix} \begin{bmatrix} \rho \\ u \\ p \end{bmatrix}_x & = 0. \end{align}The eigenvalues of the flux Jacobian $f'(q)$ for the 1D Euler equations are:
\begin{align} \lambda_1 & = u-c & \lambda_2 & = u & \lambda_3 & = u+c \end{align}Here $c$ is the sound speed:
$$ c = \sqrt{\frac{\gamma p}{\rho}}.$$These are also the eigenvalues of the coefficient matrix appearing in (\ref{euler_primitive}), and show that acoustic waves propagate at speeds $\pm c$ relative to the fluid velocity $u$. There is also a characteristic speed $\lambda_2 =u$ corresponding to the transport of entropy at the fluid velocity, as discussed further below.
The eigenvectors of the coefficient matrix appearing in (\ref{euler_primitive}) are:
\begin{align}\label{euler_evecs} r_1 & = \begin{bmatrix} -\rho/c \\ 1 \\ - \rho c \end{bmatrix} & r_2 & = \begin{bmatrix} 1 \\ 0 \\ 0 \end{bmatrix} & r_3 & = \begin{bmatrix} \rho/c \\ 1 \\ \rho c \end{bmatrix}. \end{align}These vectors show the relation between jumps in the primitive variables across waves in each family. The eigenvectors of the flux Jacobian $f'(q)$ arising from the conservative form (\ref{euler_conserved}) would be different, and would give the relation between jumps in the conserved variables across each wave.
Notice that the second characteristic speed, $\lambda_2$, depends only on $u$ and that $u$ does not change as we move in the direction of $r_2$. In other words, the 2-characteristic velocity is constant on 2-integral curves. This is similar to the wave that carries changes in the tracer that we considered in Shallow_tracer.ipynb. We say this characteristic field is linearly degenerate; it admits neither shocks nor rarefactions. In a simple 2-wave, all characteristics are parallel. A jump in this family carries a change only in the density, and is referred to as a contact discontinuity.
The other two fields have characteristic velocities that do vary along the corresponding integral curves; thus the 1-wave and the 3-wave in any Riemann solution will be either a shock or a rarefaction. We say these characteristic fields are genuinely nonlinear.
Mathematically, the $p$th field is linearly degenerate if
$$\nabla \lambda_p(q) \cdot r_p(q) = 0$$and genuinely nonlinear if
$$\nabla \lambda_p(q) \cdot r_p(q) \ne 0.$$Another important quantity in gas dynamics is the specific entropy:
$$ s = c_v \log(p/\rho^\gamma) + C,$$where $c_v$ and $C$ are constants. From the expression (\ref{euler_evecs}) for the eigenvector $r_2$, we see that the pressure and velocity are constant across a 2-wave. A simple 2-wave is also called an entropy wave because a variation in density while the pressure remains constant requires a variation in the entropy of the gas as well. On the other hand a simple acoustic wave (a continuously varying pure 1-wave or 3-wave) has constant entropy throughout the wave; the specific entropy is a Riemann invariant for these families.
A shock wave (either a 1-wave or 3-wave) satisfies the Rankine-Hugoniot conditions and exhibits a jump in entropy. To be physically correct, the entropy of the gas must increase as gas molecules pass through the shock, leading to the entropy condition for selecting shock waves. We have already seen this term used in the context of shallow water flow even though the entropy condition in that case did not involve the physical entropy.
The level sets of these Riemann invariants are two-dimensional surfaces; the intersection of two appropriate level sets defines an integral curve.
The 2-integral curves, of course, are simply lines of constant pressure and velocity (with varying density). Since the field is linearly degenerate, these coincide with the Hugoniot loci. We can determine the form of the 1- and 3-integral curves using the Riemann invariants above. For a curve passing through $(\rho_0,u_0,p_0)$, we find
\begin{align} \rho(p) &= (p/p_0)^{1/\gamma} \rho_0,\\ u(p) & = u_0 \pm \frac{2c_0}{\gamma-1}\left(1-(p/p_0)^{(\gamma-1)/(2\gamma)}\right). \end{align}Here the plus sign is for 1-waves and the minus sign is for 3-waves.
Below we plot the projection of some integral curves on the $p-u$ plane.
In [ ]:
%matplotlib inline
In [ ]:
from exact_solvers import euler
from exact_solvers import euler_demos
from ipywidgets import widgets
from ipywidgets import interact
State = euler.Primitive_State
gamma = 1.4
If you wish to examine the Python code for this chapter, see:
In [ ]:
interact(euler.plot_integral_curves,
gamma=widgets.FloatSlider(min=1.1,max=3,value=1.4),
rho_0=widgets.FloatSlider(min=0.1,max=3.,value=1.,
description=r'$\rho_0$'));
The Hugoniot loci for 1- and 3-shocks are \begin{align} \rho(p) &= \left(\frac{1 + \beta p/p_0}{p/p_l + \beta} \right),\\ u(p) & = u_0 \pm \frac{2c_0}{\sqrt{2\gamma(\gamma-1)}} \left(\frac{1-p/p_0}{\sqrt{1+\beta p/p_0}}\right), \\ \end{align} where $\beta = (\gamma+1)/(\gamma-1)$. Here the plus sign is for 1-shocks and the minus sign is for 3-shocks.
Below we plot the projection of some Hugoniot loci on the $p-u$ plane.
In [ ]:
interact(euler.plot_hugoniot_loci,
gamma=widgets.FloatSlider(min=1.1,max=3,value=1.4),
rho_0=widgets.FloatSlider(min=0.1,max=3.,value=1.,
description=r'$\rho_0$'));
As mentioned above, a shock wave is physically relevant only if the entropy of the gas increases as the gas particles move through the shock. A discontinuity satisfying the Rankine-Hugoniot jump conditions that violates this entropy condition (an "entropy-violating shock") is not physically correct and should be replaced by a rarefaction wave in the Riemann solution.
This physical entropy condition is equivalent to the mathematical condition that for a 1-shock to be physically relevant the 1-characteristics must impinge on the shock. If the entropy condition is violated, the 1-characteristics would spread out, allowing the insertion of an expansion fan (rarefaction wave).
The general Riemann solution is found following the steps listed below. This is essentially the same procedure used to determine the correct solution to the Riemann problem for the shallow water equations in Shallow_water.ipynb, where more details are given.
The Euler equations are a system of three equations and the general Riemann solution consists of three waves, so we must determine two intermediate states rather than the one intermediate state in the shallow water equations. However, it is nearly as simple because of the fact that we know the pressure and velocity are constant across the 2-wave, and so there is a single intermediate pressure $p_m$ and velocity $u_m$ in both intermediate states, and it is only the density that takes different values $\rho_{m1}$ and $\rho_{m2}$. Moreover any jump in density is allowed across the 2-wave, and we have expressions given above for how $u(p)$ varies along any integral curve or Hugoniot locus, expressions that do not explicitly involve $\rho$. So we can determine the intermediate $p_m$ by finding the intersection point of two relevant curves, in step 3 of this general algorithm:
Step 4 above requires finding the structure of rarefaction waves. This can be done using the the fact that the Riemann invariants are constant through the rarefaction wave. See Chapter 14 of (LeVeque, 2002) for more details.
Here we present some representative examples of Riemann problems and solutions. The examples chosen are closely related to the examples used in Shallow_water.ipynb and you might want to refer back to that notebook and compare the results.
First we consider the classic shock tube problem. The initial condition consists of high density and pressure on the left, low density and pressure on the right and zero velocity on both sides. The solution is composed of a shock propagating to the right (3-shock), while a left-going rarefaction forms (1-rarefaction). In between these two waves, there is a jump in the density, which is the contact discontinuity (2-wave) in the linearly degenerate characteristic field.
Note that this set of initial conditions is analogous to the "dam break" problem for shallow water quations, and the resulting structure of the solution is very similar to that obtained when those equations are solved with the addition of a scalar tracer. However, in the Euler equations the entropy jump across a 2-wave does affect the fluid dynamics on either side, so this is not a passive tracer and solving the Riemann problem is slightly more complex.
In [ ]:
left_state = State(Density = 3.,
Velocity = 0.,
Pressure = 3.)
right_state = State(Density = 1.,
Velocity = 0.,
Pressure = 1.)
euler.riemann_solution(left_state,right_state)
Here is a plot of the solution in the phase plane, showing the integral curve connecting the left and middle states, and the Hugoniot locus connecting the middle and right states.
In [ ]:
euler.phase_plane_plot(left_state, right_state)
In [ ]:
left_state = State(Density = 1.,
Velocity = -3.,
Pressure = 1.)
right_state = State(Density = 1.,
Velocity = 3.,
Pressure = 1.)
euler.riemann_solution(left_state,right_state);
In [ ]:
euler.phase_plane_plot(left_state, right_state)
In [ ]:
left_state = State(Density = 1.,
Velocity = 3.,
Pressure = 1.)
right_state = State(Density = 1.,
Velocity = -3.,
Pressure = 1.)
euler.riemann_solution(left_state,right_state)
In [ ]:
euler.phase_plane_plot(left_state, right_state)
In the next plot of the Riemann solution in the $x$-$t$ plane, we also plot the trajectories of a set of particles initially distributed along the $x$ axis at $t=0$, with the spacing inversely proportional to the density. The evolution of the distance between particles gives an indication of how the density changes.
In [ ]:
left_state = State(Density = 3.,
Velocity = 0.,
Pressure = 3.)
right_state = State(Density = 1.,
Velocity = 0.,
Pressure = 1.)
euler.plot_riemann_trajectories(left_state, right_state)
Recall that the evolution of the distance between particles gives an indication of how the density changes. Note that it increases across the shock wave and decreases through the rarefaction wave, and that in general there is a jump in density across the contact discontinuity.
Next we plot the Riemann solution with the density plot also showing an advected color to help visualize the flow better. The fluid initially to the left of $x=0$ is colored red and that initially to the right of $x=0$ is colored blue, with stripes of different shades of these colors to help visualize the motion of the fluid.
Let's plot the Sod shock tube data with this colored tracer:
In [ ]:
def plot_with_stripes_t_slider(t):
euler_demos.plot_with_stripes(rho_l=3.,u_l=0.,p_l=3.,
rho_r=1.,u_r=0.,p_r=1.,
gamma=gamma,t=t)
interact(plot_with_stripes_t_slider,
t=widgets.FloatSlider(min=0.,max=1.,step=0.1,value=0.5));
Note the following in the figure above:
In [ ]:
euler_demos.euler_demo1()
A vacuum state (with zero pressure and density) in the Euler equations is similar to a dry state (with depth $h=0$) in the shallow water equations. It can arise in the solution of the Riemann problem in two ways:
In [ ]:
left_state = State(Density =0.,
Velocity = 0.,
Pressure = 0.)
right_state = State(Density = 1.,
Velocity = -3.,
Pressure = 1.)
euler.riemann_solution(left_state,right_state)
In [ ]:
euler.phase_plane_plot(left_state, right_state)
In [ ]:
left_state = State(Density =1.,
Velocity = -10.,
Pressure = 1.)
right_state = State(Density = 1.,
Velocity = 10.,
Pressure = 1.)
euler.riemann_solution(left_state,right_state)
In [ ]:
euler.phase_plane_plot(left_state, right_state)