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