xref: /libCEED/examples/fluids/index.md (revision bcb2dfae4c301ddfdddf58806f08f6e7d17f4ea5)
1*bcb2dfaeSJed Brown(example-petsc-navier-stokes)=
2*bcb2dfaeSJed Brown
3*bcb2dfaeSJed Brown# Compressible Navier-Stokes mini-app
4*bcb2dfaeSJed Brown
5*bcb2dfaeSJed BrownThis example is located in the subdirectory {file}`examples/fluids`.
6*bcb2dfaeSJed BrownIt solves the time-dependent Navier-Stokes equations of compressible gas dynamics in a static Eulerian three-dimensional frame using unstructured high-order finite/spectral element spatial discretizations and explicit or implicit high-order time-stepping (available in PETSc).
7*bcb2dfaeSJed BrownMoreover, the Navier-Stokes example has been developed using PETSc, so that the pointwise physics (defined at quadrature points) is separated from the parallelization and meshing concerns.
8*bcb2dfaeSJed Brown
9*bcb2dfaeSJed BrownThe mathematical formulation (from {cite}`giraldoetal2010`, cf. SE3) is given in what follows.
10*bcb2dfaeSJed BrownThe compressible Navier-Stokes equations in conservative form are
11*bcb2dfaeSJed Brown
12*bcb2dfaeSJed Brown$$
13*bcb2dfaeSJed Brown\begin{aligned}
14*bcb2dfaeSJed Brown\frac{\partial \rho}{\partial t} + \nabla \cdot \bm{U} &= 0 \\
15*bcb2dfaeSJed Brown\frac{\partial \bm{U}}{\partial t} + \nabla \cdot \left( \frac{\bm{U} \otimes \bm{U}}{\rho} + P \bm{I}_3 -\bm\sigma \right) + \rho g \bm{\hat k} &= 0 \\
16*bcb2dfaeSJed Brown\frac{\partial E}{\partial t} + \nabla \cdot \left( \frac{(E + P)\bm{U}}{\rho} -\bm{u} \cdot \bm{\sigma} - k \nabla T \right) &= 0 \, , \\
17*bcb2dfaeSJed Brown\end{aligned}
18*bcb2dfaeSJed Brown$$ (eq-ns)
19*bcb2dfaeSJed Brown
20*bcb2dfaeSJed Brownwhere $\bm{\sigma} = \mu(\nabla \bm{u} + (\nabla \bm{u})^T + \lambda (\nabla \cdot \bm{u})\bm{I}_3)$ is the Cauchy (symmetric) stress tensor, with $\mu$ the dynamic viscosity coefficient, and $\lambda = - 2/3$ the Stokes hypothesis constant.
21*bcb2dfaeSJed BrownIn equations {math:numref}`eq-ns`, $\rho$ represents the volume mass density, $U$ the momentum density (defined as $\bm{U}=\rho \bm{u}$, where $\bm{u}$ is the vector velocity field), $E$ the total energy density (defined as $E = \rho e$, where $e$ is the total energy), $\bm{I}_3$ represents the $3 \times 3$ identity matrix, $g$ the gravitational acceleration constant, $\bm{\hat{k}}$ the unit vector in the $z$ direction, $k$ the thermal conductivity constant, $T$ represents the temperature, and $P$ the pressure, given by the following equation of state
22*bcb2dfaeSJed Brown
23*bcb2dfaeSJed Brown$$
24*bcb2dfaeSJed BrownP = \left( {c_p}/{c_v} -1\right) \left( E - {\bm{U}\cdot\bm{U}}/{(2 \rho)} - \rho g z \right) \, ,
25*bcb2dfaeSJed Brown$$ (eq-state)
26*bcb2dfaeSJed Brown
27*bcb2dfaeSJed Brownwhere $c_p$ is the specific heat at constant pressure and $c_v$ is the specific heat at constant volume (that define $\gamma = c_p / c_v$, the specific heat ratio).
28*bcb2dfaeSJed Brown
29*bcb2dfaeSJed BrownThe system {math:numref}`eq-ns` can be rewritten in vector form
30*bcb2dfaeSJed Brown
31*bcb2dfaeSJed Brown$$
32*bcb2dfaeSJed Brown\frac{\partial \bm{q}}{\partial t} + \nabla \cdot \bm{F}(\bm{q}) -S(\bm{q}) = 0 \, ,
33*bcb2dfaeSJed Brown$$ (eq-vector-ns)
34*bcb2dfaeSJed Brown
35*bcb2dfaeSJed Brownfor the state variables 5-dimensional vector
36*bcb2dfaeSJed Brown
37*bcb2dfaeSJed Brown$$
38*bcb2dfaeSJed Brown\bm{q} =        \begin{pmatrix}            \rho \\            \bm{U} \equiv \rho \bm{ u }\\            E \equiv \rho e        \end{pmatrix}        \begin{array}{l}            \leftarrow\textrm{ volume mass density}\\            \leftarrow\textrm{ momentum density}\\            \leftarrow\textrm{ energy density}        \end{array}
39*bcb2dfaeSJed Brown$$
40*bcb2dfaeSJed Brown
41*bcb2dfaeSJed Brownwhere the flux and the source terms, respectively, are given by
42*bcb2dfaeSJed Brown
43*bcb2dfaeSJed Brown$$
44*bcb2dfaeSJed Brown\begin{aligned}
45*bcb2dfaeSJed Brown\bm{F}(\bm{q}) &=
46*bcb2dfaeSJed Brown\begin{pmatrix}
47*bcb2dfaeSJed Brown    \bm{U}\\
48*bcb2dfaeSJed Brown    {(\bm{U} \otimes \bm{U})}/{\rho} + P \bm{I}_3 -  \bm{\sigma} \\
49*bcb2dfaeSJed Brown    {(E + P)\bm{U}}/{\rho} - \bm{u}  \cdot \bm{\sigma} - k \nabla T
50*bcb2dfaeSJed Brown\end{pmatrix} ,\\
51*bcb2dfaeSJed BrownS(\bm{q}) &=
52*bcb2dfaeSJed Brown- \begin{pmatrix}
53*bcb2dfaeSJed Brown    0\\
54*bcb2dfaeSJed Brown    \rho g \bm{\hat{k}}\\
55*bcb2dfaeSJed Brown    0
56*bcb2dfaeSJed Brown\end{pmatrix}.
57*bcb2dfaeSJed Brown\end{aligned}
58*bcb2dfaeSJed Brown$$
59*bcb2dfaeSJed Brown
60*bcb2dfaeSJed BrownLet the discrete solution be
61*bcb2dfaeSJed Brown
62*bcb2dfaeSJed Brown$$
63*bcb2dfaeSJed Brown\bm{q}_N (\bm{x},t)^{(e)} = \sum_{k=1}^{P}\psi_k (\bm{x})\bm{q}_k^{(e)}
64*bcb2dfaeSJed Brown$$
65*bcb2dfaeSJed Brown
66*bcb2dfaeSJed Brownwith $P=p+1$ the number of nodes in the element $e$.
67*bcb2dfaeSJed BrownWe use tensor-product bases $\psi_{kji} = h_i(X_0)h_j(X_1)h_k(X_2)$.
68*bcb2dfaeSJed Brown
69*bcb2dfaeSJed BrownFor the time discretization, we use two types of time stepping schemes.
70*bcb2dfaeSJed Brown
71*bcb2dfaeSJed Brown- Explicit time-stepping method
72*bcb2dfaeSJed Brown
73*bcb2dfaeSJed Brown  The following explicit formulation is solved with the adaptive Runge-Kutta-Fehlberg (RKF4-5) method by default (any explicit time-stepping scheme available in PETSc can be chosen at runtime)
74*bcb2dfaeSJed Brown
75*bcb2dfaeSJed Brown  $$
76*bcb2dfaeSJed Brown  \bm{q}_N^{n+1} = \bm{q}_N^n + \Delta t \sum_{i=1}^{s} b_i k_i \, ,
77*bcb2dfaeSJed Brown  $$
78*bcb2dfaeSJed Brown
79*bcb2dfaeSJed Brown  where
80*bcb2dfaeSJed Brown
81*bcb2dfaeSJed Brown  $$
82*bcb2dfaeSJed Brown  \begin{aligned}
83*bcb2dfaeSJed Brown     k_1 &= f(t^n, \bm{q}_N^n)\\
84*bcb2dfaeSJed Brown     k_2 &= f(t^n + c_2 \Delta t, \bm{q}_N^n + \Delta t (a_{21} k_1))\\
85*bcb2dfaeSJed Brown     k_3 &= f(t^n + c_3 \Delta t, \bm{q}_N^n + \Delta t (a_{31} k_1 + a_{32} k_2))\\
86*bcb2dfaeSJed Brown     \vdots&\\
87*bcb2dfaeSJed Brown     k_i &= f\left(t^n + c_i \Delta t, \bm{q}_N^n + \Delta t \sum_{j=1}^s a_{ij} k_j \right)\\
88*bcb2dfaeSJed Brown  \end{aligned}
89*bcb2dfaeSJed Brown  $$
90*bcb2dfaeSJed Brown
91*bcb2dfaeSJed Brown  and with
92*bcb2dfaeSJed Brown
93*bcb2dfaeSJed Brown  $$
94*bcb2dfaeSJed Brown  f(t^n, \bm{q}_N^n) = - [\nabla \cdot \bm{F}(\bm{q}_N)]^n + [S(\bm{q}_N)]^n \, .
95*bcb2dfaeSJed Brown  $$
96*bcb2dfaeSJed Brown
97*bcb2dfaeSJed Brown- Implicit time-stepping method
98*bcb2dfaeSJed Brown
99*bcb2dfaeSJed Brown  This time stepping method which can be selected using the option `-implicit` is solved with Backward Differentiation Formula (BDF) method by default (similarly, any implicit time-stepping scheme available in PETSc can be chosen at runtime).
100*bcb2dfaeSJed Brown  The implicit formulation solves nonlinear systems for $\bm q_N$:
101*bcb2dfaeSJed Brown
102*bcb2dfaeSJed Brown  $$
103*bcb2dfaeSJed Brown  \bm f(\bm q_N) \equiv \bm g(t^{n+1}, \bm{q}_N, \bm{\dot{q}}_N) = 0 \, ,
104*bcb2dfaeSJed Brown  $$ (eq-ts-implicit-ns)
105*bcb2dfaeSJed Brown
106*bcb2dfaeSJed Brown  where the time derivative $\bm{\dot q}_N$ is defined by
107*bcb2dfaeSJed Brown
108*bcb2dfaeSJed Brown  $$
109*bcb2dfaeSJed Brown  \bm{\dot{q}}_N(\bm q_N) = \alpha \bm q_N + \bm z_N
110*bcb2dfaeSJed Brown  $$
111*bcb2dfaeSJed Brown
112*bcb2dfaeSJed Brown  in terms of $\bm z_N$ from prior state and $\alpha > 0$, both of which depend on the specific time integration scheme (backward difference formulas, generalized alpha, implicit Runge-Kutta, etc.).
113*bcb2dfaeSJed Brown  Each nonlinear system {math:numref}`eq-ts-implicit-ns` will correspond to a weak form, as explained below.
114*bcb2dfaeSJed Brown  In determining how difficult a given problem is to solve, we consider the Jacobian of {math:numref}`eq-ts-implicit-ns`,
115*bcb2dfaeSJed Brown
116*bcb2dfaeSJed Brown  $$
117*bcb2dfaeSJed Brown  \frac{\partial \bm f}{\partial \bm q_N} = \frac{\partial \bm g}{\partial \bm q_N} + \alpha \frac{\partial \bm g}{\partial \bm{\dot q}_N}.
118*bcb2dfaeSJed Brown  $$
119*bcb2dfaeSJed Brown
120*bcb2dfaeSJed Brown  The scalar "shift" $\alpha$ scales inversely with the time step $\Delta t$, so small time steps result in the Jacobian being dominated by the second term, which is a sort of "mass matrix", and typically well-conditioned independent of grid resolution with a simple preconditioner (such as Jacobi).
121*bcb2dfaeSJed Brown  In contrast, the first term dominates for large time steps, with a condition number that grows with the diameter of the domain and polynomial degree of the approximation space.
122*bcb2dfaeSJed Brown  Both terms are significant for time-accurate simulation and the setup costs of strong preconditioners must be balanced with the convergence rate of Krylov methods using weak preconditioners.
123*bcb2dfaeSJed Brown
124*bcb2dfaeSJed BrownTo obtain a finite element discretization, we first multiply the strong form {math:numref}`eq-vector-ns` by a test function $\bm v \in H^1(\Omega)$ and integrate,
125*bcb2dfaeSJed Brown
126*bcb2dfaeSJed Brown$$
127*bcb2dfaeSJed Brown\int_{\Omega} \bm v \cdot \left(\frac{\partial \bm{q}_N}{\partial t} + \nabla \cdot \bm{F}(\bm{q}_N) - \bm{S}(\bm{q}_N) \right) \,dV = 0 \, , \; \forall \bm v \in \mathcal{V}_p\,,
128*bcb2dfaeSJed Brown$$
129*bcb2dfaeSJed Brown
130*bcb2dfaeSJed Brownwith $\mathcal{V}_p = \{ \bm v(\bm x) \in H^{1}(\Omega_e) \,|\, \bm v(\bm x_e(\bm X)) \in P_p(\bm{I}), e=1,\ldots,N_e \}$ a mapped space of polynomials containing at least polynomials of degree $p$ (with or without the higher mixed terms that appear in tensor product spaces).
131*bcb2dfaeSJed Brown
132*bcb2dfaeSJed BrownIntegrating by parts on the divergence term, we arrive at the weak form,
133*bcb2dfaeSJed Brown
134*bcb2dfaeSJed Brown$$
135*bcb2dfaeSJed Brown\begin{aligned}
136*bcb2dfaeSJed Brown\int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
137*bcb2dfaeSJed Brown- \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
138*bcb2dfaeSJed Brown+ \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm q_N) \cdot \widehat{\bm{n}} \,dS
139*bcb2dfaeSJed Brown  &= 0 \, , \; \forall \bm v \in \mathcal{V}_p \,,
140*bcb2dfaeSJed Brown\end{aligned}
141*bcb2dfaeSJed Brown$$ (eq-weak-vector-ns)
142*bcb2dfaeSJed Brown
143*bcb2dfaeSJed Brownwhere $\bm{F}(\bm q_N) \cdot \widehat{\bm{n}}$ is typically replaced with a boundary condition.
144*bcb2dfaeSJed Brown
145*bcb2dfaeSJed Brown:::{note}
146*bcb2dfaeSJed BrownThe notation $\nabla \bm v \!:\! \bm F$ represents contraction over both fields and spatial dimensions while a single dot represents contraction in just one, which should be clear from context, e.g., $\bm v \cdot \bm S$ contracts over fields while $\bm F \cdot \widehat{\bm n}$ contracts over spatial dimensions.
147*bcb2dfaeSJed Brown:::
148*bcb2dfaeSJed Brown
149*bcb2dfaeSJed BrownWe solve {math:numref}`eq-weak-vector-ns` using a Galerkin discretization (default) or a stabilized method, as is necessary for most real-world flows.
150*bcb2dfaeSJed Brown
151*bcb2dfaeSJed BrownGalerkin methods produce oscillations for transport-dominated problems (any time the cell Péclet number is larger than 1), and those tend to blow up for nonlinear problems such as the Euler equations and (low-viscosity/poorly resolved) Navier-Stokes, in which case stabilization is necessary.
152*bcb2dfaeSJed BrownOur formulation follows {cite}`hughesetal2010`, which offers a comprehensive review of stabilization and shock-capturing methods for continuous finite element discretization of compressible flows.
153*bcb2dfaeSJed Brown
154*bcb2dfaeSJed Brown- **SUPG** (streamline-upwind/Petrov-Galerkin)
155*bcb2dfaeSJed Brown
156*bcb2dfaeSJed Brown  In this method, the weighted residual of the strong form {math:numref}`eq-vector-ns` is added to the Galerkin formulation {math:numref}`eq-weak-vector-ns`.
157*bcb2dfaeSJed Brown  The weak form for this method is given as
158*bcb2dfaeSJed Brown
159*bcb2dfaeSJed Brown  $$
160*bcb2dfaeSJed Brown  \begin{aligned}
161*bcb2dfaeSJed Brown  \int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
162*bcb2dfaeSJed Brown  - \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
163*bcb2dfaeSJed Brown  + \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm{q}_N) \cdot \widehat{\bm{n}} \,dS & \\
164*bcb2dfaeSJed Brown  + \int_{\Omega} \bm{P}(\bm v)^T \, \left( \frac{\partial \bm{q}_N}{\partial t} \, + \,
165*bcb2dfaeSJed Brown  \nabla \cdot \bm{F} \, (\bm{q}_N) - \bm{S}(\bm{q}_N) \right) \,dV &= 0
166*bcb2dfaeSJed Brown  \, , \; \forall \bm v \in \mathcal{V}_p
167*bcb2dfaeSJed Brown  \end{aligned}
168*bcb2dfaeSJed Brown  $$ (eq-weak-vector-ns-supg)
169*bcb2dfaeSJed Brown
170*bcb2dfaeSJed Brown  This stabilization technique can be selected using the option `-stab supg`.
171*bcb2dfaeSJed Brown
172*bcb2dfaeSJed Brown- **SU** (streamline-upwind)
173*bcb2dfaeSJed Brown
174*bcb2dfaeSJed Brown  This method is a simplified version of *SUPG* {math:numref}`eq-weak-vector-ns-supg` which is developed for debugging/comparison purposes. The weak form for this method is
175*bcb2dfaeSJed Brown
176*bcb2dfaeSJed Brown  $$
177*bcb2dfaeSJed Brown  \begin{aligned}
178*bcb2dfaeSJed Brown  \int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
179*bcb2dfaeSJed Brown  - \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
180*bcb2dfaeSJed Brown  + \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm{q}_N) \cdot \widehat{\bm{n}} \,dS & \\
181*bcb2dfaeSJed Brown  + \int_{\Omega} \bm{P}(\bm v)^T \, \nabla \cdot \bm{F} \, (\bm{q}_N) \,dV
182*bcb2dfaeSJed Brown  & = 0 \, , \; \forall \bm v \in \mathcal{V}_p
183*bcb2dfaeSJed Brown  \end{aligned}
184*bcb2dfaeSJed Brown  $$ (eq-weak-vector-ns-su)
185*bcb2dfaeSJed Brown
186*bcb2dfaeSJed Brown  This stabilization technique can be selected using the option `-stab su`.
187*bcb2dfaeSJed Brown
188*bcb2dfaeSJed BrownIn both {math:numref}`eq-weak-vector-ns-su` and {math:numref}`eq-weak-vector-ns-supg`, $\bm{P} \,$ is called the *perturbation to the test-function space*, since it modifies the original Galerkin method into *SUPG* or *SU* schemes.
189*bcb2dfaeSJed BrownIt is defined as
190*bcb2dfaeSJed Brown
191*bcb2dfaeSJed Brown$$
192*bcb2dfaeSJed Brown\bm{P}(\bm v) \equiv \left(\bm{\tau} \cdot \frac{\partial \bm{F} \, (\bm{q}_N)}{\partial \bm{q}_N} \right)^T \, \nabla \bm v\,,
193*bcb2dfaeSJed Brown$$
194*bcb2dfaeSJed Brown
195*bcb2dfaeSJed Brownwhere parameter $\bm{\tau} \in \mathbb R^{3\times 3}$ is an intrinsic time/space scale matrix.
196*bcb2dfaeSJed Brown
197*bcb2dfaeSJed BrownCurrently, this demo provides three types of problems/physical models that can be selected at run time via the option `-problem`.
198*bcb2dfaeSJed Brown{ref}`problem-advection`, the problem of the transport of energy in a uniform vector velocity field, {ref}`problem-euler-vortex`, the exact solution to the Euler equations, and the so called {ref}`problem-density-current` problem.
199*bcb2dfaeSJed Brown
200*bcb2dfaeSJed Brown(problem-advection)=
201*bcb2dfaeSJed Brown
202*bcb2dfaeSJed Brown## Advection
203*bcb2dfaeSJed Brown
204*bcb2dfaeSJed BrownA simplified version of system {math:numref}`eq-ns`, only accounting for the transport of total energy, is given by
205*bcb2dfaeSJed Brown
206*bcb2dfaeSJed Brown$$
207*bcb2dfaeSJed Brown\frac{\partial E}{\partial t} + \nabla \cdot (\bm{u} E ) = 0 \, ,
208*bcb2dfaeSJed Brown$$ (eq-advection)
209*bcb2dfaeSJed Brown
210*bcb2dfaeSJed Brownwith $\bm{u}$ the vector velocity field. In this particular test case, a blob of total energy (defined by a characteristic radius $r_c$) is transported by two different wind types.
211*bcb2dfaeSJed Brown
212*bcb2dfaeSJed Brown- **Rotation**
213*bcb2dfaeSJed Brown
214*bcb2dfaeSJed Brown  In this case, a uniform circular velocity field transports the blob of total energy.
215*bcb2dfaeSJed Brown  We have solved {math:numref}`eq-advection` applying zero energy density $E$, and no-flux for $\bm{u}$ on the boundaries.
216*bcb2dfaeSJed Brown
217*bcb2dfaeSJed Brown  The $3D$ version of this test case can be run with:
218*bcb2dfaeSJed Brown
219*bcb2dfaeSJed Brown  ```
220*bcb2dfaeSJed Brown  ./navierstokes -problem advection -wind_type rotation
221*bcb2dfaeSJed Brown  ```
222*bcb2dfaeSJed Brown
223*bcb2dfaeSJed Brown  while the $2D$ version with:
224*bcb2dfaeSJed Brown
225*bcb2dfaeSJed Brown  ```
226*bcb2dfaeSJed Brown  ./navierstokes -problem advection2d -wind_type rotation
227*bcb2dfaeSJed Brown  ```
228*bcb2dfaeSJed Brown
229*bcb2dfaeSJed Brown- **Translation**
230*bcb2dfaeSJed Brown
231*bcb2dfaeSJed Brown  In this case, a background wind with a constant rectilinear velocity field, enters the domain and transports the blob of total energy out of the domain.
232*bcb2dfaeSJed Brown
233*bcb2dfaeSJed Brown  For the inflow boundary conditions, a prescribed $E_{wind}$ is applied weakly on the inflow boundaries such that the weak form boundary integral in {math:numref}`eq-weak-vector-ns` is defined as
234*bcb2dfaeSJed Brown
235*bcb2dfaeSJed Brown  $$
236*bcb2dfaeSJed Brown  \int_{\partial \Omega_{inflow}} \bm v \cdot \bm{F}(\bm q_N) \cdot \widehat{\bm{n}} \,dS = \int_{\partial \Omega_{inflow}} \bm v \, E_{wind} \, \bm u \cdot \widehat{\bm{n}} \,dS  \, ,
237*bcb2dfaeSJed Brown  $$
238*bcb2dfaeSJed Brown
239*bcb2dfaeSJed Brown  For the outflow boundary conditions, we have used the current values of $E$, following {cite}`papanastasiou1992outflow` which extends the validity of the weak form of the governing equations to the outflow instead of replacing them with unknown essential or natural boundary conditions.
240*bcb2dfaeSJed Brown  The weak form boundary integral in {math:numref}`eq-weak-vector-ns` for outflow boundary conditions is defined as
241*bcb2dfaeSJed Brown
242*bcb2dfaeSJed Brown  $$
243*bcb2dfaeSJed Brown  \int_{\partial \Omega_{outflow}} \bm v \cdot \bm{F}(\bm q_N) \cdot \widehat{\bm{n}} \,dS = \int_{\partial \Omega_{outflow}} \bm v \, E \, \bm u \cdot \widehat{\bm{n}} \,dS  \, ,
244*bcb2dfaeSJed Brown  $$
245*bcb2dfaeSJed Brown
246*bcb2dfaeSJed Brown  The $3D$ version of this test case problem can be run with:
247*bcb2dfaeSJed Brown
248*bcb2dfaeSJed Brown  ```
249*bcb2dfaeSJed Brown  ./navierstokes -problem advection -wind_type translation -wind_translation .5,-1,0
250*bcb2dfaeSJed Brown  ```
251*bcb2dfaeSJed Brown
252*bcb2dfaeSJed Brown  while the $2D$ version with:
253*bcb2dfaeSJed Brown
254*bcb2dfaeSJed Brown  ```
255*bcb2dfaeSJed Brown  ./navierstokes -problem advection2d -wind_type translation -wind_translation 1,-.5
256*bcb2dfaeSJed Brown  ```
257*bcb2dfaeSJed Brown
258*bcb2dfaeSJed Brown(problem-euler-vortex)=
259*bcb2dfaeSJed Brown
260*bcb2dfaeSJed Brown## Isentropic Vortex
261*bcb2dfaeSJed Brown
262*bcb2dfaeSJed BrownThree-dimensional Euler equations, which are simplified version of system {math:numref}`eq-ns` and account only for the convective fluxes, are given by
263*bcb2dfaeSJed Brown
264*bcb2dfaeSJed Brown$$
265*bcb2dfaeSJed Brown\begin{aligned}
266*bcb2dfaeSJed Brown\frac{\partial \rho}{\partial t} + \nabla \cdot \bm{U} &= 0 \\
267*bcb2dfaeSJed Brown\frac{\partial \bm{U}}{\partial t} + \nabla \cdot \left( \frac{\bm{U} \otimes \bm{U}}{\rho} + P \bm{I}_3 \right) &= 0 \\
268*bcb2dfaeSJed Brown\frac{\partial E}{\partial t} + \nabla \cdot \left( \frac{(E + P)\bm{U}}{\rho} \right) &= 0 \, , \\
269*bcb2dfaeSJed Brown\end{aligned}
270*bcb2dfaeSJed Brown$$ (eq-euler)
271*bcb2dfaeSJed Brown
272*bcb2dfaeSJed BrownFollowing the setup given in {cite}`zhang2011verification`, the mean flow for this problem is $\rho=1$, $P=1$, $T=P/\rho= 1$, and $\bm{u}=(u_1,u_2,0)$ while the perturbation $\delta \bm{u}$, and $\delta T$ are defined as
273*bcb2dfaeSJed Brown
274*bcb2dfaeSJed Brown$$
275*bcb2dfaeSJed Brown\begin{aligned} (\delta u_1, \, \delta u_2) &= \frac{\epsilon}{2 \pi} \, e^{0.5(1-r^2)} \, (-\bar{y}, \, \bar{x}) \, , \\ \delta T &= - \frac{(\gamma-1) \, \epsilon^2}{8 \, \gamma \, \pi^2} \, e^{1-r^2} \, , \\ \end{aligned}
276*bcb2dfaeSJed Brown$$
277*bcb2dfaeSJed Brown
278*bcb2dfaeSJed Brownwhere $(\bar{x}, \, \bar{y}) = (x-x_c, \, y-y_c)$, $(x_c, \, y_c)$ represents the center of the domain, $r^2=\bar{x}^2 + \bar{y}^2$, and $\epsilon$ is the vortex strength.
279*bcb2dfaeSJed BrownThere is no perturbation in the entropy $S=P/\rho^\gamma$ ($\delta S=0)$.
280*bcb2dfaeSJed Brown
281*bcb2dfaeSJed BrownThis problem can be run with:
282*bcb2dfaeSJed Brown
283*bcb2dfaeSJed Brown```
284*bcb2dfaeSJed Brown./navierstokes -problem euler_vortex -mean_velocity .5,-.8,0.
285*bcb2dfaeSJed Brown```
286*bcb2dfaeSJed Brown
287*bcb2dfaeSJed Brown(problem-density-current)=
288*bcb2dfaeSJed Brown
289*bcb2dfaeSJed Brown## Density Current
290*bcb2dfaeSJed Brown
291*bcb2dfaeSJed BrownFor this test problem (from {cite}`straka1993numerical`), we solve the full Navier-Stokes equations {math:numref}`eq-ns`, for which a cold air bubble (of radius $r_c$) drops by convection in a neutrally stratified atmosphere.
292*bcb2dfaeSJed BrownIts initial condition is defined in terms of the Exner pressure, $\pi(\bm{x},t)$, and potential temperature, $\theta(\bm{x},t)$, that relate to the state variables via
293*bcb2dfaeSJed Brown
294*bcb2dfaeSJed Brown$$
295*bcb2dfaeSJed Brown\begin{aligned} \rho &= \frac{P_0}{( c_p - c_v)\theta(\bm{x},t)} \pi(\bm{x},t)^{\frac{c_v}{ c_p - c_v}} \, , \\ e &= c_v \theta(\bm{x},t) \pi(\bm{x},t) + \bm{u}\cdot \bm{u} /2 + g z \, , \end{aligned}
296*bcb2dfaeSJed Brown$$
297*bcb2dfaeSJed Brown
298*bcb2dfaeSJed Brownwhere $P_0$ is the atmospheric pressure.
299*bcb2dfaeSJed BrownFor this problem, we have used no-slip and non-penetration boundary conditions for $\bm{u}$, and no-flux for mass and energy densities.
300*bcb2dfaeSJed BrownThis problem can be run with:
301*bcb2dfaeSJed Brown
302*bcb2dfaeSJed Brown```
303*bcb2dfaeSJed Brown./navierstokes -problem density_current
304*bcb2dfaeSJed Brown```
305