xref: /libCEED/examples/fluids/index.md (revision bc7bbd5dca53df9bb3e78c945912120808c904be)
1bcb2dfaeSJed Brown(example-petsc-navier-stokes)=
2bcb2dfaeSJed Brown
3bcb2dfaeSJed Brown# Compressible Navier-Stokes mini-app
4bcb2dfaeSJed Brown
5bcb2dfaeSJed BrownThis example is located in the subdirectory {file}`examples/fluids`.
6bcb2dfaeSJed 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).
7bcb2dfaeSJed 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.
8bcb2dfaeSJed Brown
9*bc7bbd5dSLeila Ghaffari## Running the mini-app
10*bc7bbd5dSLeila Ghaffari
11*bc7bbd5dSLeila Ghaffari```{include} README.md
12*bc7bbd5dSLeila Ghaffari:start-after: inclusion-fluids-marker
13*bc7bbd5dSLeila Ghaffari```
14*bc7bbd5dSLeila Ghaffari## The Navier-Stokes equations
15*bc7bbd5dSLeila Ghaffari
16bcb2dfaeSJed BrownThe mathematical formulation (from {cite}`giraldoetal2010`, cf. SE3) is given in what follows.
17bcb2dfaeSJed BrownThe compressible Navier-Stokes equations in conservative form are
18bcb2dfaeSJed Brown
19bcb2dfaeSJed Brown$$
20bcb2dfaeSJed Brown\begin{aligned}
21bcb2dfaeSJed Brown\frac{\partial \rho}{\partial t} + \nabla \cdot \bm{U} &= 0 \\
22bcb2dfaeSJed 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 \\
23bcb2dfaeSJed 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 \, , \\
24bcb2dfaeSJed Brown\end{aligned}
25bcb2dfaeSJed Brown$$ (eq-ns)
26bcb2dfaeSJed Brown
27bcb2dfaeSJed 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.
288791656fSJed BrownIn equations {eq}`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
29bcb2dfaeSJed Brown
30bcb2dfaeSJed Brown$$
31bcb2dfaeSJed BrownP = \left( {c_p}/{c_v} -1\right) \left( E - {\bm{U}\cdot\bm{U}}/{(2 \rho)} - \rho g z \right) \, ,
32bcb2dfaeSJed Brown$$ (eq-state)
33bcb2dfaeSJed Brown
34bcb2dfaeSJed 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).
35bcb2dfaeSJed Brown
368791656fSJed BrownThe system {eq}`eq-ns` can be rewritten in vector form
37bcb2dfaeSJed Brown
38bcb2dfaeSJed Brown$$
39bcb2dfaeSJed Brown\frac{\partial \bm{q}}{\partial t} + \nabla \cdot \bm{F}(\bm{q}) -S(\bm{q}) = 0 \, ,
40bcb2dfaeSJed Brown$$ (eq-vector-ns)
41bcb2dfaeSJed Brown
42bcb2dfaeSJed Brownfor the state variables 5-dimensional vector
43bcb2dfaeSJed Brown
44bcb2dfaeSJed Brown$$
45bcb2dfaeSJed 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}
46bcb2dfaeSJed Brown$$
47bcb2dfaeSJed Brown
48bcb2dfaeSJed Brownwhere the flux and the source terms, respectively, are given by
49bcb2dfaeSJed Brown
50bcb2dfaeSJed Brown$$
51bcb2dfaeSJed Brown\begin{aligned}
52bcb2dfaeSJed Brown\bm{F}(\bm{q}) &=
53bcb2dfaeSJed Brown\begin{pmatrix}
54bcb2dfaeSJed Brown    \bm{U}\\
55bcb2dfaeSJed Brown    {(\bm{U} \otimes \bm{U})}/{\rho} + P \bm{I}_3 -  \bm{\sigma} \\
56bcb2dfaeSJed Brown    {(E + P)\bm{U}}/{\rho} - \bm{u}  \cdot \bm{\sigma} - k \nabla T
57bcb2dfaeSJed Brown\end{pmatrix} ,\\
58bcb2dfaeSJed BrownS(\bm{q}) &=
59bcb2dfaeSJed Brown- \begin{pmatrix}
60bcb2dfaeSJed Brown    0\\
61bcb2dfaeSJed Brown    \rho g \bm{\hat{k}}\\
62bcb2dfaeSJed Brown    0
63bcb2dfaeSJed Brown\end{pmatrix}.
64bcb2dfaeSJed Brown\end{aligned}
65bcb2dfaeSJed Brown$$
66bcb2dfaeSJed Brown
67bcb2dfaeSJed BrownLet the discrete solution be
68bcb2dfaeSJed Brown
69bcb2dfaeSJed Brown$$
70bcb2dfaeSJed Brown\bm{q}_N (\bm{x},t)^{(e)} = \sum_{k=1}^{P}\psi_k (\bm{x})\bm{q}_k^{(e)}
71bcb2dfaeSJed Brown$$
72bcb2dfaeSJed Brown
73bcb2dfaeSJed Brownwith $P=p+1$ the number of nodes in the element $e$.
74bcb2dfaeSJed BrownWe use tensor-product bases $\psi_{kji} = h_i(X_0)h_j(X_1)h_k(X_2)$.
75bcb2dfaeSJed Brown
76bcb2dfaeSJed BrownFor the time discretization, we use two types of time stepping schemes.
77bcb2dfaeSJed Brown
78bcb2dfaeSJed Brown- Explicit time-stepping method
79bcb2dfaeSJed Brown
80bcb2dfaeSJed 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)
81bcb2dfaeSJed Brown
82bcb2dfaeSJed Brown  $$
83bcb2dfaeSJed Brown  \bm{q}_N^{n+1} = \bm{q}_N^n + \Delta t \sum_{i=1}^{s} b_i k_i \, ,
84bcb2dfaeSJed Brown  $$
85bcb2dfaeSJed Brown
86bcb2dfaeSJed Brown  where
87bcb2dfaeSJed Brown
88bcb2dfaeSJed Brown  $$
89bcb2dfaeSJed Brown  \begin{aligned}
90bcb2dfaeSJed Brown     k_1 &= f(t^n, \bm{q}_N^n)\\
91bcb2dfaeSJed Brown     k_2 &= f(t^n + c_2 \Delta t, \bm{q}_N^n + \Delta t (a_{21} k_1))\\
92bcb2dfaeSJed Brown     k_3 &= f(t^n + c_3 \Delta t, \bm{q}_N^n + \Delta t (a_{31} k_1 + a_{32} k_2))\\
93bcb2dfaeSJed Brown     \vdots&\\
94bcb2dfaeSJed 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)\\
95bcb2dfaeSJed Brown  \end{aligned}
96bcb2dfaeSJed Brown  $$
97bcb2dfaeSJed Brown
98bcb2dfaeSJed Brown  and with
99bcb2dfaeSJed Brown
100bcb2dfaeSJed Brown  $$
101bcb2dfaeSJed Brown  f(t^n, \bm{q}_N^n) = - [\nabla \cdot \bm{F}(\bm{q}_N)]^n + [S(\bm{q}_N)]^n \, .
102bcb2dfaeSJed Brown  $$
103bcb2dfaeSJed Brown
104bcb2dfaeSJed Brown- Implicit time-stepping method
105bcb2dfaeSJed Brown
106bcb2dfaeSJed 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).
107bcb2dfaeSJed Brown  The implicit formulation solves nonlinear systems for $\bm q_N$:
108bcb2dfaeSJed Brown
109bcb2dfaeSJed Brown  $$
110bcb2dfaeSJed Brown  \bm f(\bm q_N) \equiv \bm g(t^{n+1}, \bm{q}_N, \bm{\dot{q}}_N) = 0 \, ,
111bcb2dfaeSJed Brown  $$ (eq-ts-implicit-ns)
112bcb2dfaeSJed Brown
113bcb2dfaeSJed Brown  where the time derivative $\bm{\dot q}_N$ is defined by
114bcb2dfaeSJed Brown
115bcb2dfaeSJed Brown  $$
116bcb2dfaeSJed Brown  \bm{\dot{q}}_N(\bm q_N) = \alpha \bm q_N + \bm z_N
117bcb2dfaeSJed Brown  $$
118bcb2dfaeSJed Brown
119bcb2dfaeSJed 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.).
1208791656fSJed Brown  Each nonlinear system {eq}`eq-ts-implicit-ns` will correspond to a weak form, as explained below.
1218791656fSJed Brown  In determining how difficult a given problem is to solve, we consider the Jacobian of {eq}`eq-ts-implicit-ns`,
122bcb2dfaeSJed Brown
123bcb2dfaeSJed Brown  $$
124bcb2dfaeSJed 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}.
125bcb2dfaeSJed Brown  $$
126bcb2dfaeSJed Brown
127bcb2dfaeSJed 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).
128bcb2dfaeSJed 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.
129bcb2dfaeSJed 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.
130bcb2dfaeSJed Brown
1318791656fSJed BrownTo obtain a finite element discretization, we first multiply the strong form {eq}`eq-vector-ns` by a test function $\bm v \in H^1(\Omega)$ and integrate,
132bcb2dfaeSJed Brown
133bcb2dfaeSJed Brown$$
134bcb2dfaeSJed 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\,,
135bcb2dfaeSJed Brown$$
136bcb2dfaeSJed Brown
137bcb2dfaeSJed 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).
138bcb2dfaeSJed Brown
139bcb2dfaeSJed BrownIntegrating by parts on the divergence term, we arrive at the weak form,
140bcb2dfaeSJed Brown
141bcb2dfaeSJed Brown$$
142bcb2dfaeSJed Brown\begin{aligned}
143bcb2dfaeSJed Brown\int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
144bcb2dfaeSJed Brown- \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
145bcb2dfaeSJed Brown+ \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm q_N) \cdot \widehat{\bm{n}} \,dS
146bcb2dfaeSJed Brown  &= 0 \, , \; \forall \bm v \in \mathcal{V}_p \,,
147bcb2dfaeSJed Brown\end{aligned}
148bcb2dfaeSJed Brown$$ (eq-weak-vector-ns)
149bcb2dfaeSJed Brown
150bcb2dfaeSJed Brownwhere $\bm{F}(\bm q_N) \cdot \widehat{\bm{n}}$ is typically replaced with a boundary condition.
151bcb2dfaeSJed Brown
152bcb2dfaeSJed Brown:::{note}
153bcb2dfaeSJed 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.
154bcb2dfaeSJed Brown:::
155bcb2dfaeSJed Brown
1568791656fSJed BrownWe solve {eq}`eq-weak-vector-ns` using a Galerkin discretization (default) or a stabilized method, as is necessary for most real-world flows.
157bcb2dfaeSJed Brown
158bcb2dfaeSJed 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.
159bcb2dfaeSJed BrownOur formulation follows {cite}`hughesetal2010`, which offers a comprehensive review of stabilization and shock-capturing methods for continuous finite element discretization of compressible flows.
160bcb2dfaeSJed Brown
161bcb2dfaeSJed Brown- **SUPG** (streamline-upwind/Petrov-Galerkin)
162bcb2dfaeSJed Brown
1638791656fSJed Brown  In this method, the weighted residual of the strong form {eq}`eq-vector-ns` is added to the Galerkin formulation {eq}`eq-weak-vector-ns`.
164bcb2dfaeSJed Brown  The weak form for this method is given as
165bcb2dfaeSJed Brown
166bcb2dfaeSJed Brown  $$
167bcb2dfaeSJed Brown  \begin{aligned}
168bcb2dfaeSJed Brown  \int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
169bcb2dfaeSJed Brown  - \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
170bcb2dfaeSJed Brown  + \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm{q}_N) \cdot \widehat{\bm{n}} \,dS & \\
171bcb2dfaeSJed Brown  + \int_{\Omega} \bm{P}(\bm v)^T \, \left( \frac{\partial \bm{q}_N}{\partial t} \, + \,
172bcb2dfaeSJed Brown  \nabla \cdot \bm{F} \, (\bm{q}_N) - \bm{S}(\bm{q}_N) \right) \,dV &= 0
173bcb2dfaeSJed Brown  \, , \; \forall \bm v \in \mathcal{V}_p
174bcb2dfaeSJed Brown  \end{aligned}
175bcb2dfaeSJed Brown  $$ (eq-weak-vector-ns-supg)
176bcb2dfaeSJed Brown
177bcb2dfaeSJed Brown  This stabilization technique can be selected using the option `-stab supg`.
178bcb2dfaeSJed Brown
179bcb2dfaeSJed Brown- **SU** (streamline-upwind)
180bcb2dfaeSJed Brown
1818791656fSJed Brown  This method is a simplified version of *SUPG* {eq}`eq-weak-vector-ns-supg` which is developed for debugging/comparison purposes. The weak form for this method is
182bcb2dfaeSJed Brown
183bcb2dfaeSJed Brown  $$
184bcb2dfaeSJed Brown  \begin{aligned}
185bcb2dfaeSJed Brown  \int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
186bcb2dfaeSJed Brown  - \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
187bcb2dfaeSJed Brown  + \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm{q}_N) \cdot \widehat{\bm{n}} \,dS & \\
188bcb2dfaeSJed Brown  + \int_{\Omega} \bm{P}(\bm v)^T \, \nabla \cdot \bm{F} \, (\bm{q}_N) \,dV
189bcb2dfaeSJed Brown  & = 0 \, , \; \forall \bm v \in \mathcal{V}_p
190bcb2dfaeSJed Brown  \end{aligned}
191bcb2dfaeSJed Brown  $$ (eq-weak-vector-ns-su)
192bcb2dfaeSJed Brown
193bcb2dfaeSJed Brown  This stabilization technique can be selected using the option `-stab su`.
194bcb2dfaeSJed Brown
1958791656fSJed BrownIn both {eq}`eq-weak-vector-ns-su` and {eq}`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.
196bcb2dfaeSJed BrownIt is defined as
197bcb2dfaeSJed Brown
198bcb2dfaeSJed Brown$$
199bcb2dfaeSJed 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\,,
200bcb2dfaeSJed Brown$$
201bcb2dfaeSJed Brown
202bcb2dfaeSJed Brownwhere parameter $\bm{\tau} \in \mathbb R^{3\times 3}$ is an intrinsic time/space scale matrix.
203bcb2dfaeSJed Brown
204bcb2dfaeSJed BrownCurrently, this demo provides three types of problems/physical models that can be selected at run time via the option `-problem`.
205bcb2dfaeSJed 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.
206bcb2dfaeSJed Brown
207bcb2dfaeSJed Brown(problem-advection)=
208bcb2dfaeSJed Brown
209bcb2dfaeSJed Brown## Advection
210bcb2dfaeSJed Brown
2118791656fSJed BrownA simplified version of system {eq}`eq-ns`, only accounting for the transport of total energy, is given by
212bcb2dfaeSJed Brown
213bcb2dfaeSJed Brown$$
214bcb2dfaeSJed Brown\frac{\partial E}{\partial t} + \nabla \cdot (\bm{u} E ) = 0 \, ,
215bcb2dfaeSJed Brown$$ (eq-advection)
216bcb2dfaeSJed Brown
217bcb2dfaeSJed 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.
218bcb2dfaeSJed Brown
219bcb2dfaeSJed Brown- **Rotation**
220bcb2dfaeSJed Brown
221bcb2dfaeSJed Brown  In this case, a uniform circular velocity field transports the blob of total energy.
2228791656fSJed Brown  We have solved {eq}`eq-advection` applying zero energy density $E$, and no-flux for $\bm{u}$ on the boundaries.
223bcb2dfaeSJed Brown
224bcb2dfaeSJed Brown- **Translation**
225bcb2dfaeSJed Brown
226bcb2dfaeSJed 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.
227bcb2dfaeSJed Brown
2288791656fSJed 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 {eq}`eq-weak-vector-ns` is defined as
229bcb2dfaeSJed Brown
230bcb2dfaeSJed Brown  $$
231bcb2dfaeSJed 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  \, ,
232bcb2dfaeSJed Brown  $$
233bcb2dfaeSJed Brown
234bcb2dfaeSJed 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.
2358791656fSJed Brown  The weak form boundary integral in {eq}`eq-weak-vector-ns` for outflow boundary conditions is defined as
236bcb2dfaeSJed Brown
237bcb2dfaeSJed Brown  $$
238bcb2dfaeSJed 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  \, ,
239bcb2dfaeSJed Brown  $$
240bcb2dfaeSJed Brown
241bcb2dfaeSJed Brown(problem-euler-vortex)=
242bcb2dfaeSJed Brown
243bcb2dfaeSJed Brown## Isentropic Vortex
244bcb2dfaeSJed Brown
245*bc7bbd5dSLeila GhaffariThree-dimensional Euler equations, which are simplified and nondimensionalized version of system {eq}`eq-ns` and account only for the convective fluxes, are given by
246bcb2dfaeSJed Brown
247bcb2dfaeSJed Brown$$
248bcb2dfaeSJed Brown\begin{aligned}
249bcb2dfaeSJed Brown\frac{\partial \rho}{\partial t} + \nabla \cdot \bm{U} &= 0 \\
250bcb2dfaeSJed Brown\frac{\partial \bm{U}}{\partial t} + \nabla \cdot \left( \frac{\bm{U} \otimes \bm{U}}{\rho} + P \bm{I}_3 \right) &= 0 \\
251bcb2dfaeSJed Brown\frac{\partial E}{\partial t} + \nabla \cdot \left( \frac{(E + P)\bm{U}}{\rho} \right) &= 0 \, , \\
252bcb2dfaeSJed Brown\end{aligned}
253bcb2dfaeSJed Brown$$ (eq-euler)
254bcb2dfaeSJed Brown
255*bc7bbd5dSLeila GhaffariFollowing the setup given in {cite}`zhang2011verification`, the mean flow for this problem is $\rho=1$, $P=1$, $T=P/\rho= 1$ (Specific Gas Constant, $R$, is 1), and $\bm{u}=(u_1,u_2,0)$ while the perturbation $\delta \bm{u}$, and $\delta T$ are defined as
256bcb2dfaeSJed Brown
257bcb2dfaeSJed Brown$$
258bcb2dfaeSJed 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}
259bcb2dfaeSJed Brown$$
260bcb2dfaeSJed Brown
261*bc7bbd5dSLeila Ghaffariwhere $(\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 ($\epsilon$ < 10).
262bcb2dfaeSJed BrownThere is no perturbation in the entropy $S=P/\rho^\gamma$ ($\delta S=0)$.
263bcb2dfaeSJed Brown
264bcb2dfaeSJed Brown(problem-density-current)=
265bcb2dfaeSJed Brown
266bcb2dfaeSJed Brown## Density Current
267bcb2dfaeSJed Brown
2688791656fSJed BrownFor this test problem (from {cite}`straka1993numerical`), we solve the full Navier-Stokes equations {eq}`eq-ns`, for which a cold air bubble (of radius $r_c$) drops by convection in a neutrally stratified atmosphere.
269bcb2dfaeSJed 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
270bcb2dfaeSJed Brown
271bcb2dfaeSJed Brown$$
272bcb2dfaeSJed 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}
273bcb2dfaeSJed Brown$$
274bcb2dfaeSJed Brown
275bcb2dfaeSJed Brownwhere $P_0$ is the atmospheric pressure.
276bcb2dfaeSJed BrownFor this problem, we have used no-slip and non-penetration boundary conditions for $\bm{u}$, and no-flux for mass and energy densities.
277