xref: /honee/index.md (revision d783cc74b86284fc2d343457b9abcec1031cc945)
1*d783cc74SJed Brown(example-petsc-navier-stokes)=
2*d783cc74SJed Brown
3*d783cc74SJed Brown# Compressible Navier-Stokes mini-app
4*d783cc74SJed Brown
5*d783cc74SJed BrownThis example is located in the subdirectory {file}`examples/fluids`.
6*d783cc74SJed 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*d783cc74SJed 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*d783cc74SJed Brown
9*d783cc74SJed BrownThe mathematical formulation (from {cite}`giraldoetal2010`, cf. SE3) is given in what follows.
10*d783cc74SJed BrownThe compressible Navier-Stokes equations in conservative form are
11*d783cc74SJed Brown
12*d783cc74SJed Brown$$
13*d783cc74SJed Brown\begin{aligned}
14*d783cc74SJed Brown\frac{\partial \rho}{\partial t} + \nabla \cdot \bm{U} &= 0 \\
15*d783cc74SJed 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*d783cc74SJed 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*d783cc74SJed Brown\end{aligned}
18*d783cc74SJed Brown$$ (eq-ns)
19*d783cc74SJed Brown
20*d783cc74SJed 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*d783cc74SJed 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*d783cc74SJed Brown
23*d783cc74SJed Brown$$
24*d783cc74SJed BrownP = \left( {c_p}/{c_v} -1\right) \left( E - {\bm{U}\cdot\bm{U}}/{(2 \rho)} - \rho g z \right) \, ,
25*d783cc74SJed Brown$$ (eq-state)
26*d783cc74SJed Brown
27*d783cc74SJed 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*d783cc74SJed Brown
29*d783cc74SJed BrownThe system {math:numref}`eq-ns` can be rewritten in vector form
30*d783cc74SJed Brown
31*d783cc74SJed Brown$$
32*d783cc74SJed Brown\frac{\partial \bm{q}}{\partial t} + \nabla \cdot \bm{F}(\bm{q}) -S(\bm{q}) = 0 \, ,
33*d783cc74SJed Brown$$ (eq-vector-ns)
34*d783cc74SJed Brown
35*d783cc74SJed Brownfor the state variables 5-dimensional vector
36*d783cc74SJed Brown
37*d783cc74SJed Brown$$
38*d783cc74SJed 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*d783cc74SJed Brown$$
40*d783cc74SJed Brown
41*d783cc74SJed Brownwhere the flux and the source terms, respectively, are given by
42*d783cc74SJed Brown
43*d783cc74SJed Brown$$
44*d783cc74SJed Brown\begin{aligned}
45*d783cc74SJed Brown\bm{F}(\bm{q}) &=
46*d783cc74SJed Brown\begin{pmatrix}
47*d783cc74SJed Brown    \bm{U}\\
48*d783cc74SJed Brown    {(\bm{U} \otimes \bm{U})}/{\rho} + P \bm{I}_3 -  \bm{\sigma} \\
49*d783cc74SJed Brown    {(E + P)\bm{U}}/{\rho} - \bm{u}  \cdot \bm{\sigma} - k \nabla T
50*d783cc74SJed Brown\end{pmatrix} ,\\
51*d783cc74SJed BrownS(\bm{q}) &=
52*d783cc74SJed Brown- \begin{pmatrix}
53*d783cc74SJed Brown    0\\
54*d783cc74SJed Brown    \rho g \bm{\hat{k}}\\
55*d783cc74SJed Brown    0
56*d783cc74SJed Brown\end{pmatrix}.
57*d783cc74SJed Brown\end{aligned}
58*d783cc74SJed Brown$$
59*d783cc74SJed Brown
60*d783cc74SJed BrownLet the discrete solution be
61*d783cc74SJed Brown
62*d783cc74SJed Brown$$
63*d783cc74SJed Brown\bm{q}_N (\bm{x},t)^{(e)} = \sum_{k=1}^{P}\psi_k (\bm{x})\bm{q}_k^{(e)}
64*d783cc74SJed Brown$$
65*d783cc74SJed Brown
66*d783cc74SJed Brownwith $P=p+1$ the number of nodes in the element $e$.
67*d783cc74SJed BrownWe use tensor-product bases $\psi_{kji} = h_i(X_0)h_j(X_1)h_k(X_2)$.
68*d783cc74SJed Brown
69*d783cc74SJed BrownFor the time discretization, we use two types of time stepping schemes.
70*d783cc74SJed Brown
71*d783cc74SJed Brown- Explicit time-stepping method
72*d783cc74SJed Brown
73*d783cc74SJed 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*d783cc74SJed Brown
75*d783cc74SJed Brown  $$
76*d783cc74SJed Brown  \bm{q}_N^{n+1} = \bm{q}_N^n + \Delta t \sum_{i=1}^{s} b_i k_i \, ,
77*d783cc74SJed Brown  $$
78*d783cc74SJed Brown
79*d783cc74SJed Brown  where
80*d783cc74SJed Brown
81*d783cc74SJed Brown  $$
82*d783cc74SJed Brown  \begin{aligned}
83*d783cc74SJed Brown     k_1 &= f(t^n, \bm{q}_N^n)\\
84*d783cc74SJed Brown     k_2 &= f(t^n + c_2 \Delta t, \bm{q}_N^n + \Delta t (a_{21} k_1))\\
85*d783cc74SJed 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*d783cc74SJed Brown     \vdots&\\
87*d783cc74SJed 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*d783cc74SJed Brown  \end{aligned}
89*d783cc74SJed Brown  $$
90*d783cc74SJed Brown
91*d783cc74SJed Brown  and with
92*d783cc74SJed Brown
93*d783cc74SJed Brown  $$
94*d783cc74SJed Brown  f(t^n, \bm{q}_N^n) = - [\nabla \cdot \bm{F}(\bm{q}_N)]^n + [S(\bm{q}_N)]^n \, .
95*d783cc74SJed Brown  $$
96*d783cc74SJed Brown
97*d783cc74SJed Brown- Implicit time-stepping method
98*d783cc74SJed Brown
99*d783cc74SJed 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*d783cc74SJed Brown  The implicit formulation solves nonlinear systems for $\bm q_N$:
101*d783cc74SJed Brown
102*d783cc74SJed Brown  $$
103*d783cc74SJed Brown  \bm f(\bm q_N) \equiv \bm g(t^{n+1}, \bm{q}_N, \bm{\dot{q}}_N) = 0 \, ,
104*d783cc74SJed Brown  $$ (eq-ts-implicit-ns)
105*d783cc74SJed Brown
106*d783cc74SJed Brown  where the time derivative $\bm{\dot q}_N$ is defined by
107*d783cc74SJed Brown
108*d783cc74SJed Brown  $$
109*d783cc74SJed Brown  \bm{\dot{q}}_N(\bm q_N) = \alpha \bm q_N + \bm z_N
110*d783cc74SJed Brown  $$
111*d783cc74SJed Brown
112*d783cc74SJed 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*d783cc74SJed Brown  Each nonlinear system {math:numref}`eq-ts-implicit-ns` will correspond to a weak form, as explained below.
114*d783cc74SJed Brown  In determining how difficult a given problem is to solve, we consider the Jacobian of {math:numref}`eq-ts-implicit-ns`,
115*d783cc74SJed Brown
116*d783cc74SJed Brown  $$
117*d783cc74SJed 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*d783cc74SJed Brown  $$
119*d783cc74SJed Brown
120*d783cc74SJed 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*d783cc74SJed 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*d783cc74SJed 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*d783cc74SJed Brown
124*d783cc74SJed 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*d783cc74SJed Brown
126*d783cc74SJed Brown$$
127*d783cc74SJed 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*d783cc74SJed Brown$$
129*d783cc74SJed Brown
130*d783cc74SJed 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*d783cc74SJed Brown
132*d783cc74SJed BrownIntegrating by parts on the divergence term, we arrive at the weak form,
133*d783cc74SJed Brown
134*d783cc74SJed Brown$$
135*d783cc74SJed Brown\begin{aligned}
136*d783cc74SJed Brown\int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
137*d783cc74SJed Brown- \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
138*d783cc74SJed Brown+ \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm q_N) \cdot \widehat{\bm{n}} \,dS
139*d783cc74SJed Brown  &= 0 \, , \; \forall \bm v \in \mathcal{V}_p \,,
140*d783cc74SJed Brown\end{aligned}
141*d783cc74SJed Brown$$ (eq-weak-vector-ns)
142*d783cc74SJed Brown
143*d783cc74SJed Brownwhere $\bm{F}(\bm q_N) \cdot \widehat{\bm{n}}$ is typically replaced with a boundary condition.
144*d783cc74SJed Brown
145*d783cc74SJed Brown:::{note}
146*d783cc74SJed 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*d783cc74SJed Brown:::
148*d783cc74SJed Brown
149*d783cc74SJed 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*d783cc74SJed Brown
151*d783cc74SJed 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*d783cc74SJed 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*d783cc74SJed Brown
154*d783cc74SJed Brown- **SUPG** (streamline-upwind/Petrov-Galerkin)
155*d783cc74SJed Brown
156*d783cc74SJed 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*d783cc74SJed Brown  The weak form for this method is given as
158*d783cc74SJed Brown
159*d783cc74SJed Brown  $$
160*d783cc74SJed Brown  \begin{aligned}
161*d783cc74SJed Brown  \int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
162*d783cc74SJed Brown  - \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
163*d783cc74SJed Brown  + \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm{q}_N) \cdot \widehat{\bm{n}} \,dS & \\
164*d783cc74SJed Brown  + \int_{\Omega} \bm{P}(\bm v)^T \, \left( \frac{\partial \bm{q}_N}{\partial t} \, + \,
165*d783cc74SJed Brown  \nabla \cdot \bm{F} \, (\bm{q}_N) - \bm{S}(\bm{q}_N) \right) \,dV &= 0
166*d783cc74SJed Brown  \, , \; \forall \bm v \in \mathcal{V}_p
167*d783cc74SJed Brown  \end{aligned}
168*d783cc74SJed Brown  $$ (eq-weak-vector-ns-supg)
169*d783cc74SJed Brown
170*d783cc74SJed Brown  This stabilization technique can be selected using the option `-stab supg`.
171*d783cc74SJed Brown
172*d783cc74SJed Brown- **SU** (streamline-upwind)
173*d783cc74SJed Brown
174*d783cc74SJed 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*d783cc74SJed Brown
176*d783cc74SJed Brown  $$
177*d783cc74SJed Brown  \begin{aligned}
178*d783cc74SJed Brown  \int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
179*d783cc74SJed Brown  - \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
180*d783cc74SJed Brown  + \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm{q}_N) \cdot \widehat{\bm{n}} \,dS & \\
181*d783cc74SJed Brown  + \int_{\Omega} \bm{P}(\bm v)^T \, \nabla \cdot \bm{F} \, (\bm{q}_N) \,dV
182*d783cc74SJed Brown  & = 0 \, , \; \forall \bm v \in \mathcal{V}_p
183*d783cc74SJed Brown  \end{aligned}
184*d783cc74SJed Brown  $$ (eq-weak-vector-ns-su)
185*d783cc74SJed Brown
186*d783cc74SJed Brown  This stabilization technique can be selected using the option `-stab su`.
187*d783cc74SJed Brown
188*d783cc74SJed 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*d783cc74SJed BrownIt is defined as
190*d783cc74SJed Brown
191*d783cc74SJed Brown$$
192*d783cc74SJed 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*d783cc74SJed Brown$$
194*d783cc74SJed Brown
195*d783cc74SJed Brownwhere parameter $\bm{\tau} \in \mathbb R^{3\times 3}$ is an intrinsic time/space scale matrix.
196*d783cc74SJed Brown
197*d783cc74SJed BrownCurrently, this demo provides three types of problems/physical models that can be selected at run time via the option `-problem`.
198*d783cc74SJed 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*d783cc74SJed Brown
200*d783cc74SJed Brown(problem-advection)=
201*d783cc74SJed Brown
202*d783cc74SJed Brown## Advection
203*d783cc74SJed Brown
204*d783cc74SJed BrownA simplified version of system {math:numref}`eq-ns`, only accounting for the transport of total energy, is given by
205*d783cc74SJed Brown
206*d783cc74SJed Brown$$
207*d783cc74SJed Brown\frac{\partial E}{\partial t} + \nabla \cdot (\bm{u} E ) = 0 \, ,
208*d783cc74SJed Brown$$ (eq-advection)
209*d783cc74SJed Brown
210*d783cc74SJed 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*d783cc74SJed Brown
212*d783cc74SJed Brown- **Rotation**
213*d783cc74SJed Brown
214*d783cc74SJed Brown  In this case, a uniform circular velocity field transports the blob of total energy.
215*d783cc74SJed Brown  We have solved {math:numref}`eq-advection` applying zero energy density $E$, and no-flux for $\bm{u}$ on the boundaries.
216*d783cc74SJed Brown
217*d783cc74SJed Brown  The $3D$ version of this test case can be run with:
218*d783cc74SJed Brown
219*d783cc74SJed Brown  ```
220*d783cc74SJed Brown  ./navierstokes -problem advection -wind_type rotation
221*d783cc74SJed Brown  ```
222*d783cc74SJed Brown
223*d783cc74SJed Brown  while the $2D$ version with:
224*d783cc74SJed Brown
225*d783cc74SJed Brown  ```
226*d783cc74SJed Brown  ./navierstokes -problem advection2d -wind_type rotation
227*d783cc74SJed Brown  ```
228*d783cc74SJed Brown
229*d783cc74SJed Brown- **Translation**
230*d783cc74SJed Brown
231*d783cc74SJed 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*d783cc74SJed Brown
233*d783cc74SJed 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*d783cc74SJed Brown
235*d783cc74SJed Brown  $$
236*d783cc74SJed 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*d783cc74SJed Brown  $$
238*d783cc74SJed Brown
239*d783cc74SJed 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*d783cc74SJed Brown  The weak form boundary integral in {math:numref}`eq-weak-vector-ns` for outflow boundary conditions is defined as
241*d783cc74SJed Brown
242*d783cc74SJed Brown  $$
243*d783cc74SJed 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*d783cc74SJed Brown  $$
245*d783cc74SJed Brown
246*d783cc74SJed Brown  The $3D$ version of this test case problem can be run with:
247*d783cc74SJed Brown
248*d783cc74SJed Brown  ```
249*d783cc74SJed Brown  ./navierstokes -problem advection -wind_type translation -wind_translation .5,-1,0
250*d783cc74SJed Brown  ```
251*d783cc74SJed Brown
252*d783cc74SJed Brown  while the $2D$ version with:
253*d783cc74SJed Brown
254*d783cc74SJed Brown  ```
255*d783cc74SJed Brown  ./navierstokes -problem advection2d -wind_type translation -wind_translation 1,-.5
256*d783cc74SJed Brown  ```
257*d783cc74SJed Brown
258*d783cc74SJed Brown(problem-euler-vortex)=
259*d783cc74SJed Brown
260*d783cc74SJed Brown## Isentropic Vortex
261*d783cc74SJed Brown
262*d783cc74SJed 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*d783cc74SJed Brown
264*d783cc74SJed Brown$$
265*d783cc74SJed Brown\begin{aligned}
266*d783cc74SJed Brown\frac{\partial \rho}{\partial t} + \nabla \cdot \bm{U} &= 0 \\
267*d783cc74SJed Brown\frac{\partial \bm{U}}{\partial t} + \nabla \cdot \left( \frac{\bm{U} \otimes \bm{U}}{\rho} + P \bm{I}_3 \right) &= 0 \\
268*d783cc74SJed Brown\frac{\partial E}{\partial t} + \nabla \cdot \left( \frac{(E + P)\bm{U}}{\rho} \right) &= 0 \, , \\
269*d783cc74SJed Brown\end{aligned}
270*d783cc74SJed Brown$$ (eq-euler)
271*d783cc74SJed Brown
272*d783cc74SJed 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*d783cc74SJed Brown
274*d783cc74SJed Brown$$
275*d783cc74SJed 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*d783cc74SJed Brown$$
277*d783cc74SJed Brown
278*d783cc74SJed 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*d783cc74SJed BrownThere is no perturbation in the entropy $S=P/\rho^\gamma$ ($\delta S=0)$.
280*d783cc74SJed Brown
281*d783cc74SJed BrownThis problem can be run with:
282*d783cc74SJed Brown
283*d783cc74SJed Brown```
284*d783cc74SJed Brown./navierstokes -problem euler_vortex -mean_velocity .5,-.8,0.
285*d783cc74SJed Brown```
286*d783cc74SJed Brown
287*d783cc74SJed Brown(problem-density-current)=
288*d783cc74SJed Brown
289*d783cc74SJed Brown## Density Current
290*d783cc74SJed Brown
291*d783cc74SJed 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*d783cc74SJed 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*d783cc74SJed Brown
294*d783cc74SJed Brown$$
295*d783cc74SJed 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*d783cc74SJed Brown$$
297*d783cc74SJed Brown
298*d783cc74SJed Brownwhere $P_0$ is the atmospheric pressure.
299*d783cc74SJed 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*d783cc74SJed BrownThis problem can be run with:
301*d783cc74SJed Brown
302*d783cc74SJed Brown```
303*d783cc74SJed Brown./navierstokes -problem density_current
304*d783cc74SJed Brown```
305