xref: /libCEED/examples/fluids/index.md (revision 135921ecb6efac3644642fa3818aa472a1a2755d)
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
9## Running the mini-app
10
11```{include} README.md
12:start-after: inclusion-fluids-marker
13```
14## The Navier-Stokes equations
15
16The mathematical formulation (from {cite}`giraldoetal2010`, cf. SE3) is given in what follows.
17The compressible Navier-Stokes equations in conservative form are
18
19$$
20\begin{aligned}
21\frac{\partial \rho}{\partial t} + \nabla \cdot \bm{U} &= 0 \\
22\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 \\
23\frac{\partial E}{\partial t} + \nabla \cdot \left( \frac{(E + P)\bm{U}}{\rho} -\bm{u} \cdot \bm{\sigma} - k \nabla T \right) &= 0 \, , \\
24\end{aligned}
25$$ (eq-ns)
26
27where $\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.
28In 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
29
30$$
31P = \left( {c_p}/{c_v} -1\right) \left( E - {\bm{U}\cdot\bm{U}}/{(2 \rho)} - \rho g z \right) \, ,
32$$ (eq-state)
33
34where $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).
35
36The system {eq}`eq-ns` can be rewritten in vector form
37
38$$
39\frac{\partial \bm{q}}{\partial t} + \nabla \cdot \bm{F}(\bm{q}) -S(\bm{q}) = 0 \, ,
40$$ (eq-vector-ns)
41
42for the state variables 5-dimensional vector
43
44$$
45\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}
46$$
47
48where the flux and the source terms, respectively, are given by
49
50$$
51\begin{aligned}
52\bm{F}(\bm{q}) &=
53\underbrace{\begin{pmatrix}
54    \bm{U}\\
55    {(\bm{U} \otimes \bm{U})}/{\rho} + P \bm{I}_3 \\
56    {(E + P)\bm{U}}/{\rho}
57\end{pmatrix}}_{\bm F_{\text{adv}}} +
58\underbrace{\begin{pmatrix}
590 \\
60-  \bm{\sigma} \\
61 - \bm{u}  \cdot \bm{\sigma} - k \nabla T
62\end{pmatrix}}_{\bm F_{\text{diff}}},\\
63S(\bm{q}) &=
64- \begin{pmatrix}
65    0\\
66    \rho g \bm{\hat{k}}\\
67    0
68\end{pmatrix}.
69\end{aligned}
70$$ (eq-ns-flux)
71
72### Finite Element Formulation (Spatial Discretization)
73
74Let the discrete solution be
75
76$$
77\bm{q}_N (\bm{x},t)^{(e)} = \sum_{k=1}^{P}\psi_k (\bm{x})\bm{q}_k^{(e)}
78$$
79
80with $P=p+1$ the number of nodes in the element $e$.
81We use tensor-product bases $\psi_{kji} = h_i(X_0)h_j(X_1)h_k(X_2)$.
82
83To 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,
84
85$$
86\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\,,
87$$
88
89with $\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).
90
91Integrating by parts on the divergence term, we arrive at the weak form,
92
93$$
94\begin{aligned}
95\int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
96- \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
97+ \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm q_N) \cdot \widehat{\bm{n}} \,dS
98  &= 0 \, , \; \forall \bm v \in \mathcal{V}_p \,,
99\end{aligned}
100$$ (eq-weak-vector-ns)
101
102where $\bm{F}(\bm q_N) \cdot \widehat{\bm{n}}$ is typically replaced with a boundary condition.
103
104:::{note}
105The 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.
106:::
107
108### Time Discretization
109For the time discretization, we use two types of time stepping schemes through PETSc.
110
111#### Explicit time-stepping method
112
113  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)
114
115  $$
116  \bm{q}_N^{n+1} = \bm{q}_N^n + \Delta t \sum_{i=1}^{s} b_i k_i \, ,
117  $$
118
119  where
120
121  $$
122  \begin{aligned}
123     k_1 &= f(t^n, \bm{q}_N^n)\\
124     k_2 &= f(t^n + c_2 \Delta t, \bm{q}_N^n + \Delta t (a_{21} k_1))\\
125     k_3 &= f(t^n + c_3 \Delta t, \bm{q}_N^n + \Delta t (a_{31} k_1 + a_{32} k_2))\\
126     \vdots&\\
127     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)\\
128  \end{aligned}
129  $$
130
131  and with
132
133  $$
134  f(t^n, \bm{q}_N^n) = - [\nabla \cdot \bm{F}(\bm{q}_N)]^n + [S(\bm{q}_N)]^n \, .
135  $$
136
137#### Implicit time-stepping method
138
139  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).
140  The implicit formulation solves nonlinear systems for $\bm q_N$:
141
142  $$
143  \bm f(\bm q_N) \equiv \bm g(t^{n+1}, \bm{q}_N, \bm{\dot{q}}_N) = 0 \, ,
144  $$ (eq-ts-implicit-ns)
145
146  where the time derivative $\bm{\dot q}_N$ is defined by
147
148  $$
149  \bm{\dot{q}}_N(\bm q_N) = \alpha \bm q_N + \bm z_N
150  $$
151
152  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.).
153  Each nonlinear system {eq}`eq-ts-implicit-ns` will correspond to a weak form, as explained below.
154  In determining how difficult a given problem is to solve, we consider the Jacobian of {eq}`eq-ts-implicit-ns`,
155
156  $$
157  \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}.
158  $$
159
160  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).
161  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.
162  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.
163
164More details of PETSc's time stepping solvers can be found in the [TS User Guide](https://petsc.org/release/docs/manual/ts/).
165
166### Stabilization
167We solve {eq}`eq-weak-vector-ns` using a Galerkin discretization (default) or a stabilized method, as is necessary for most real-world flows.
168
169Galerkin 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.
170Our formulation follows {cite}`hughesetal2010`, which offers a comprehensive review of stabilization and shock-capturing methods for continuous finite element discretization of compressible flows.
171
172- **SUPG** (streamline-upwind/Petrov-Galerkin)
173
174  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`.
175  The weak form for this method is given as
176
177  $$
178  \begin{aligned}
179  \int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
180  - \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
181  + \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm{q}_N) \cdot \widehat{\bm{n}} \,dS & \\
182  + \int_{\Omega} \nabla\bm v \tcolon\left(\frac{\partial \bm F_{\text{adv}}}{\partial \bm q}\right) \bm\tau \left( \frac{\partial \bm{q}_N}{\partial t} \, + \,
183  \nabla \cdot \bm{F} \, (\bm{q}_N) - \bm{S}(\bm{q}_N) \right) \,dV &= 0
184  \, , \; \forall \bm v \in \mathcal{V}_p
185  \end{aligned}
186  $$ (eq-weak-vector-ns-supg)
187
188  This stabilization technique can be selected using the option `-stab supg`.
189
190- **SU** (streamline-upwind)
191
192  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
193
194  $$
195  \begin{aligned}
196  \int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
197  - \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
198  + \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm{q}_N) \cdot \widehat{\bm{n}} \,dS & \\
199  + \int_{\Omega} \nabla\bm v \tcolon\left(\frac{\partial \bm F_{\text{adv}}}{\partial \bm q}\right) \bm\tau \nabla \cdot \bm{F} \, (\bm{q}_N) \,dV
200  & = 0 \, , \; \forall \bm v \in \mathcal{V}_p
201  \end{aligned}
202  $$ (eq-weak-vector-ns-su)
203
204  This stabilization technique can be selected using the option `-stab su`.
205
206In both {eq}`eq-weak-vector-ns-su` and {eq}`eq-weak-vector-ns-supg`, $\bm\tau \in \mathbb R^{5\times 5}$ (field indices) is an intrinsic time scale matrix.
207The SUPG technique and the operator $\frac{\partial \bm F_{\text{adv}}}{\partial \bm q}$ (rather than its transpose) can be explained via an ansatz for subgrid state fluctuations $\tilde{\bm q} = -\bm\tau \bm r$ where $\bm r$ is a strong form residual.
208The forward variational form can be readily expressed by differentiating $\bm F_{\text{adv}}$ of {eq}`eq-ns-flux`
209
210$$
211\begin{aligned}
212\diff\bm F_{\text{adv}}(\diff\bm q; \bm q) &= \frac{\partial \bm F_{\text{adv}}}{\partial \bm q} \diff\bm q \\
213&= \begin{pmatrix}
214\diff\bm U \\
215(\diff\bm U \otimes \bm U + \bm U \otimes \diff\bm U)/\rho - (\bm U \otimes \bm U)/\rho^2 \diff\rho + \diff P \bm I_3 \\
216(E + P)\diff\bm U/\rho + (\diff E + \diff P)\bm U/\rho - (E + P) \bm U/\rho^2 \diff\rho
217\end{pmatrix},
218\end{aligned}
219$$
220
221where $\diff P$ is defined by differentiating {eq}`eq-state`.
222
223:::{dropdown} Stabilization scale $\bm\tau$
224A velocity vector $\bm u$ can be pulled back to the reference element as $\bm u_{\bm X} = \nabla_{\bm x}\bm X \cdot \bm u$, with units of reference length (non-dimensional) per second.
225To build intuition, consider a boundary layer element of dimension $(1, \epsilon)$, for which $\nabla_{\bm x} \bm X = \bigl(\begin{smallmatrix} 2 & \\ & 2/\epsilon \end{smallmatrix}\bigr)$.
226So a small normal component of velocity will be amplified (by a factor of the aspect ratio $1/\epsilon$) in this transformation.
227The ratio $\lVert \bm u \rVert / \lVert \bm u_{\bm X} \rVert$ is a covariant measure of (half) the element length in the direction of the velocity.
228A contravariant measure of element length in the direction of a unit vector $\hat{\bm n}$ is given by $\lVert \bigl(\nabla_{\bm X} \bm x\bigr)^T \hat{\bm n} \rVert$.
229While $\nabla_{\bm X} \bm x$ is readily computable, its inverse $\nabla_{\bm x} \bm X$ is needed directly in finite element methods and thus more convenient for our use.
230If we consider a parallelogram, the covariant measure is larger than the contravariant measure for vectors pointing between acute corners and the opposite holds for vectors between oblique corners.
231
232The cell Péclet number is classically defined by $\mathrm{Pe}_h = \lVert \bm u \rVert h / (2 \kappa)$ where $\kappa$ is the diffusivity (units of $m^2/s$).
233This can be generalized to arbitrary grids by defining the local Péclet number
234
235$$
236\mathrm{Pe} = \frac{\lVert \bm u \rVert^2}{\lVert \bm u_{\bm X} \rVert \kappa}.
237$$ (eq-peclet)
238
239For scalar advection-diffusion, the stabilization is a scalar
240
241$$
242\tau = \frac{\xi(\mathrm{Pe})}{\lVert \bm u_{\bm X} \rVert},
243$$ (eq-tau-advdiff)
244
245where $\xi(\mathrm{Pe}) = \coth \mathrm{Pe} - 1/\mathrm{Pe}$ approaches 1 at large local Péclet number.
246Note that $\tau$ has units of time and, in the transport-dominated limit, is proportional to element transit time in the direction of the propagating wave.
247For advection-diffusion, $\bm F(q) = \bm u q$, and thus the SU stabilization term is
248
249$$
250\nabla v \cdot \bm u \tau \bm u \cdot \nabla q = \nabla_{\bm X} v \cdot (\bm u_{\bm X} \tau \bm u_{\bm X}) \cdot \nabla_{\bm X} q .
251$$ (eq-su-stabilize-advdiff)
252
253where the term in parentheses is a rank-1 diffusivity tensor that has been pulled back to the reference element.
254See {cite}`hughesetal2010` equations 15-17 and 34-36 for further discussion of this formulation.
255
256For the Navier-Stokes and Euler equations, {cite}`whiting2003hierarchical` defines a $5\times 5$ diagonal stabilization $\mathrm{diag}(\tau_c, \tau_m, \tau_m, \tau_m, \tau_E)$ consisting of
2571. continuity stabilization $\tau_c$
2582. momentum stabilization $\tau_m$
2593. energy stabilization $\tau_E$
260
261The Navier-Stokes code in this example uses the following formulation for $\tau_c$, $\tau_m$, $\tau_E$:
262
263$$
264\begin{aligned}
265
266\tau_c &= \frac{C_c \mathcal{F}}{8\rho \trace(\bm g)} \\
267\tau_m &= \frac{C_m}{\mathcal{F}} \\
268\tau_E &= \frac{C_E}{\mathcal{F} c_v} \\
269\end{aligned}
270$$
271
272$$
273\mathcal{F} = \sqrt{ \rho^2 \left [ \left(\frac{2C_t}{\Delta t}\right)^2
274+ \bm u \cdot (\bm u \cdot  \bm g)
275+ C_v \mu^2 \Vert \bm g \Vert_F ^2\right]}
276$$
277
278where $\bm g = \nabla_{\bm x} \bm{X} \cdot \nabla_{\bm x} \bm{X}$ is the metric tensor and $\Vert \cdot \Vert_F$ is the Frobenius norm.
279This formulation is currently not available in the Euler code.
280
281In the Euler code, we follow {cite}`hughesetal2010` in defining a $3\times 3$ diagonal stabilization according to spatial criterion 2 (equation 27) as follows.
282
283$$
284\tau_{ii} = c_{\tau} \frac{2 \xi(\mathrm{Pe})}{(\lambda_{\max \text{abs}})_i \lVert \nabla_{x_i} \bm X \rVert}
285$$ (eq-tau-conservative)
286
287where $c_{\tau}$ is a multiplicative constant reported to be optimal at 0.5 for linear elements, $\hat{\bm n}_i$ is a unit vector in direction $i$, and $\nabla_{x_i} = \hat{\bm n}_i \cdot \nabla_{\bm x}$ is the derivative in direction $i$.
288The flux Jacobian $\frac{\partial \bm F_{\text{adv}}}{\partial \bm q} \cdot \hat{\bm n}_i$ in each direction $i$ is a $5\times 5$ matrix with spectral radius $(\lambda_{\max \text{abs}})_i$ equal to the fastest wave speed.
289The complete set of eigenvalues of the Euler flux Jacobian in direction $i$ are (e.g., {cite}`toro2009`)
290
291$$
292\Lambda_i = [u_i - a, u_i, u_i, u_i, u_i+a],
293$$ (eq-eigval-advdiff)
294
295where $u_i = \bm u \cdot \hat{\bm n}_i$ is the velocity component in direction $i$ and $a = \sqrt{\gamma P/\rho}$ is the sound speed for ideal gasses.
296Note that the first and last eigenvalues represent nonlinear acoustic waves while the middle three are linearly degenerate, carrying a contact wave (temperature) and transverse components of momentum.
297The fastest wave speed in direction $i$ is thus
298
299$$
300\lambda_{\max \text{abs}} \Bigl( \frac{\partial \bm F_{\text{adv}}}{\partial \bm q} \cdot \hat{\bm n}_i \Bigr) = |u_i| + a
301$$ (eq-wavespeed)
302
303Note that this wave speed is specific to ideal gases as $\gamma$ is an ideal gas parameter; other equations of state will yield a different acoustic wave speed.
304
305:::
306
307Currently, this demo provides three types of problems/physical models that can be selected at run time via the option `-problem`.
308{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.
309
310### Subgrid Stress Modeling
311
312When a fluid simulation is under-resolved (the smallest length scale resolved by the grid is much larger than the smallest physical scale, the [Kolmogorov length scale](https://en.wikipedia.org/wiki/Kolmogorov_microscales)), this is mathematically interpreted as filtering the Navier-Stokes equations.
313This is known as large-eddy simulation (LES), as only the "large" scales of turbulence are resolved.
314This filtering operation results in an extra stress-like term, $\bm{\tau}^r$, representing the effect of unresolved (or "subgrid" scale) structures in the flow.
315Denoting the filtering operation by $\overline \cdot$, the LES governing equations are:
316
317$$
318\frac{\partial \bm{\overline q}}{\partial t} + \nabla \cdot \bm{\overline F}(\bm{\overline q}) -S(\bm{\overline q}) = 0 \, ,
319$$ (eq-vector-les)
320
321where
322
323$$
324\bm{\overline F}(\bm{\overline q}) =
325\bm{F} (\bm{\overline q}) +
326\begin{pmatrix}
327    0\\
328     \bm{\tau}^r \\
329     \bm{u}  \cdot \bm{\tau}^r
330\end{pmatrix}
331$$ (eq-les-flux)
332
333More details on deriving the above expression, filtering, and large eddy simulation can be found in {cite}`popeTurbulentFlows2000`.
334To close the problem, the subgrid stress must be defined.
335For implicit LES, the subgrid stress is set to zero and the numerical properties of the discretized system are assumed to account for the effect of subgrid scale structures on the filtered solution field.
336For explicit LES, it is defined by a subgrid stress model.
337
338#### Data-driven SGS Model
339
340The data-driven SGS model implemented here uses a small neural network to compute the SGS term.
341The SGS tensor is calculated at nodes using an $L^2$ projection of the velocity gradient and grid anisotropy tensor, and then interpolated onto quadrature points.
342More details regarding the theoretical background of the model can be found in {cite}`prakashDDSGS2022` and {cite}`prakashDDSGSAnisotropic2022`.
343
344The neural network itself consists of 1 hidden layer and 20 neurons, using Leaky ReLU as its activation function.
345The slope parameter for the Leaky ReLU function is set via `-sgs_model_dd_leakyrelu_alpha`.
346The outputs of the network are assumed to be normalized on a min-max scale, so they must be rescaled by the original min-max bounds.
347Parameters for the neural network are put into files in a directory found in `-sgs_model_dd_parameter_dir`.
348These files store the network weights (`w1.dat` and `w2.dat`), biases (`b1.dat` and `b2.dat`), and scaling parameters (`OutScaling.dat`).
349The first row of each files stores the number of columns and rows in each file.
350Note that the weight coefficients are assumed to be in column-major order.
351This is done to keep consistent with legacy file compatibility.
352
353(problem-advection)=
354
355## Advection
356
357A simplified version of system {eq}`eq-ns`, only accounting for the transport of total energy, is given by
358
359$$
360\frac{\partial E}{\partial t} + \nabla \cdot (\bm{u} E ) = 0 \, ,
361$$ (eq-advection)
362
363with $\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.
364
365- **Rotation**
366
367  In this case, a uniform circular velocity field transports the blob of total energy.
368  We have solved {eq}`eq-advection` applying zero energy density $E$, and no-flux for $\bm{u}$ on the boundaries.
369
370- **Translation**
371
372  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.
373
374  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
375
376  $$
377  \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  \, ,
378  $$
379
380  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.
381  The weak form boundary integral in {eq}`eq-weak-vector-ns` for outflow boundary conditions is defined as
382
383  $$
384  \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  \, ,
385  $$
386
387(problem-euler-vortex)=
388
389## Isentropic Vortex
390
391Three-dimensional Euler equations, which are simplified and nondimensionalized version of system {eq}`eq-ns` and account only for the convective fluxes, are given by
392
393$$
394\begin{aligned}
395\frac{\partial \rho}{\partial t} + \nabla \cdot \bm{U} &= 0 \\
396\frac{\partial \bm{U}}{\partial t} + \nabla \cdot \left( \frac{\bm{U} \otimes \bm{U}}{\rho} + P \bm{I}_3 \right) &= 0 \\
397\frac{\partial E}{\partial t} + \nabla \cdot \left( \frac{(E + P)\bm{U}}{\rho} \right) &= 0 \, , \\
398\end{aligned}
399$$ (eq-euler)
400
401Following 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
402
403$$
404\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}
405$$
406
407where $(\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).
408There is no perturbation in the entropy $S=P/\rho^\gamma$ ($\delta S=0)$.
409
410(problem-shock-tube)=
411
412## Shock Tube
413
414This test problem is based on Sod's Shock Tube (from{cite}`sodshocktubewiki`), a canonical test case for discontinuity capturing in one dimension. For this problem, the three-dimensional Euler equations are formulated exactly as in the Isentropic Vortex problem. The default initial conditions are $P=1$, $\rho=1$ for the driver section and $P=0.1$, $\rho=0.125$ for the driven section. The initial velocity is zero in both sections. Slip boundary conditions are applied to the side walls and wall boundary conditions are applied at the end walls.
415
416SU upwinding and discontinuity capturing have been implemented into the explicit timestepping operator for this problem. Discontinuity capturing is accomplished using a modified version of the $YZ\beta$ operator described in {cite}`tezduyar2007yzb`. This discontinuity capturing scheme involves the introduction of a dissipation term of the form
417
418$$
419\int_{\Omega} \nu_{SHOCK} \nabla \bm v \!:\! \nabla \bm q dV
420$$
421
422The shock capturing viscosity is implemented following the first formulation described in {cite}`tezduyar2007yzb`. The characteristic velocity $u_{cha}$ is taken to be the acoustic speed while the reference density $\rho_{ref}$ is just the local density. Shock capturing viscosity is defined by the following
423
424$$
425\nu_{SHOCK} = \tau_{SHOCK} u_{cha}^2
426$$
427
428where,
429
430$$
431\tau_{SHOCK} = \frac{h_{SHOCK}}{2u_{cha}} \left( \frac{ \,|\, \nabla \rho \,|\, h_{SHOCK}}{\rho_{ref}} \right)^{\beta}
432$$
433
434$\beta$ is a tuning parameter set between 1 (smoother shocks) and 2 (sharper shocks. The parameter $h_{SHOCK}$ is a length scale that is proportional to the element length in the direction of the density gradient unit vector. This density gradient unit vector is defined as $\hat{\bm j} = \frac{\nabla \rho}{|\nabla \rho|}$. The original formulation of Tezduyar and Senga relies on the shape function gradient to define the element length scale, but this gradient is not available to qFunctions in libCEED. To avoid this problem, $h_{SHOCK}$ is defined in the current implementation as
435
436$$
437h_{SHOCK} = 2 \left( C_{YZB} \,|\, \bm p \,|\, \right)^{-1}
438$$
439
440where
441
442$$
443p_k = \hat{j}_i \frac{\partial \xi_i}{x_k}
444$$
445
446The constant $C_{YZB}$ is set to 0.1 for piecewise linear elements in the current implementation. Larger values approaching unity are expected with more robust stabilization and implicit timestepping.
447
448(problem-density-current)=
449
450## Gaussian Wave
451This test case is taken/inspired by that presented in {cite}`mengaldoCompressibleBC2014`. It is intended to test non-reflecting/Riemann boundary conditions. It's primarily intended for Euler equations, but has been implemented for the Navier-Stokes equations here for flexibility.
452
453The problem has a perturbed initial condition and lets it evolve in time. The initial condition contains a Gaussian perturbation in the pressure field:
454
455$$
456\begin{aligned}
457\rho &= \rho_\infty\left(1+A\exp\left(\frac{-(\bar{x}^2 + \bar{y}^2)}{2\sigma^2}\right)\right) \\
458\bm{U} &= \bm U_\infty \\
459E &= \frac{p_\infty}{\gamma -1}\left(1+A\exp\left(\frac{-(\bar{x}^2 + \bar{y}^2)}{2\sigma^2}\right)\right) + \frac{\bm U_\infty \cdot \bm U_\infty}{2\rho_\infty},
460\end{aligned}
461$$
462
463where $A$ and $\sigma$ are the amplitude and width of the perturbation, respectively, and $(\bar{x}, \bar{y}) = (x-x_e, y-y_e)$ is the distance to the epicenter of the perturbation, $(x_e, y_e)$.
464The simulation produces a strong acoustic wave and leaves behind a cold thermal bubble that advects at the fluid velocity.
465
466The boundary conditions are freestream in the x and y directions. When using an HLL (Harten, Lax, van Leer) Riemann solver {cite}`toro2009` (option `-freestream_riemann hll`), the acoustic waves exit the domain cleanly, but when the thermal bubble reaches the boundary, it produces strong thermal oscillations that become acoustic waves reflecting into the domain.
467This problem can be fixed using a more sophisticated Riemann solver such as HLLC {cite}`toro2009` (option `-freestream_riemann hllc`, which is default), which is a linear constant-pressure wave that transports temperature and transverse momentum at the fluid velocity.
468
469## Vortex Shedding - Flow past Cylinder
470This test case, based on {cite}`shakib1991femcfd`, is an example of using an externally provided mesh from Gmsh.
471A cylinder with diameter $D=1$ is centered at $(0,0)$ in a computational domain $-4.5 \leq x \leq 15.5$, $-4.5 \leq y \leq 4.5$.
472We solve this as a 3D problem with (default) one element in the $z$ direction.
473The domain is filled with an ideal gas at rest (zero velocity) with temperature 24.92 and pressure 7143.
474The viscosity is 0.01 and thermal conductivity is 14.34 to maintain a Prandtl number of 0.71, which is typical for air.
475At time $t=0$, this domain is subjected to freestream boundary conditions at the inflow (left) and Riemann-type outflow on the right, with exterior reference state at velocity $(1, 0, 0)$ giving Reynolds number $100$ and Mach number $0.01$.
476A symmetry (adiabatic free slip) condition is imposed at the top and bottom boundaries $(y = \pm 4.5)$ (zero normal velocity component, zero heat-flux).
477The cylinder wall is an adiabatic (no heat flux) no-slip boundary condition.
478As we evolve in time, eddies appear past the cylinder leading to a vortex shedding known as the vortex street, with shedding period of about 6.
479
480The Gmsh input file, `examples/fluids/meshes/cylinder.geo` is parametrized to facilitate experimenting with similar configurations.
481The Strouhal number (nondimensional shedding frequency) is sensitive to the size of the computational domain and boundary conditions.
482
483Forces on the cylinder walls are computed using the "reaction force" method, which is variationally consistent with the volume operator.
484Given the force components $\bm F = (F_x, F_y, F_z)$ and surface area $S = \pi D L_z$ where $L_z$ is the spanwise extent of the domain, we define the coefficients of lift and drag as
485
486$$
487\begin{aligned}
488C_L &= \frac{2 F_y}{\rho_\infty u_\infty^2 S} \\
489C_D &= \frac{2 F_x}{\rho_\infty u_\infty^2 S} \\
490\end{aligned}
491$$
492
493where $\rho_\infty, u_\infty$ are the freestream (inflow) density and velocity respectively.
494
495## Density Current
496
497For 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.
498Its 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
499
500$$
501\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}
502$$
503
504where $P_0$ is the atmospheric pressure.
505For this problem, we have used no-slip and non-penetration boundary conditions for $\bm{u}$, and no-flux for mass and energy densities.
506
507## Channel
508
509A compressible channel flow. Analytical solution given in
510{cite}`whitingStabilizedFEM1999`:
511
512$$ u_1 = u_{\max} \left [ 1 - \left ( \frac{x_2}{H}\right)^2 \right] \quad \quad u_2 = u_3 = 0$$
513$$T = T_w \left [ 1 + \frac{Pr \hat{E}c}{3} \left \{1 - \left(\frac{x_2}{H}\right)^4  \right \} \right]$$
514$$p = p_0 - \frac{2\rho_0 u_{\max}^2 x_1}{Re_H H}$$
515
516where $H$ is the channel half-height, $u_{\max}$ is the center velocity, $T_w$ is the temperature at the wall, $Pr=\frac{\mu}{c_p \kappa}$ is the Prandlt number, $\hat E_c = \frac{u_{\max}^2}{c_p T_w}$ is the modified Eckert number, and $Re_h = \frac{u_{\max}H}{\nu}$ is the Reynolds number.
517
518Boundary conditions are periodic in the streamwise direction, and no-slip and non-penetration boundary conditions at the walls.
519The flow is driven by a body force determined analytically from the fluid properties and setup parameters $H$ and $u_{\max}$.
520
521## Flat Plate Boundary Layer
522
523### Laminar Boundary Layer - Blasius
524
525Simulation of a laminar boundary layer flow, with the inflow being prescribed
526by a [Blasius similarity
527solution](https://en.wikipedia.org/wiki/Blasius_boundary_layer). At the inflow,
528the velocity is prescribed by the Blasius soution profile, density is set
529constant, and temperature is allowed to float. Using `weakT: true`, density is
530allowed to float and temperature is set constant. At the outlet, a user-set
531pressure is used for pressure in the inviscid flux terms (all other inviscid
532flux terms use interior solution values). The wall is a no-slip,
533no-penetration, no-heat flux condition. The top of the domain is treated as an
534outflow and is tilted at a downward angle to ensure that flow is always exiting
535it.
536
537### Turbulent Boundary Layer
538
539Simulating a turbulent boundary layer without modeling the turbulence requires
540resolving the turbulent flow structures. These structures may be introduced
541into the simulations either by allowing a laminar boundary layer naturally
542transition to turbulence, or imposing turbulent structures at the inflow. The
543latter approach has been taken here, specifically using a *synthetic turbulence
544generation* (STG) method.
545
546#### Synthetic Turbulence Generation (STG) Boundary Condition
547
548We use the STG method described in
549{cite}`shurSTG2014`. Below follows a re-description of the formulation to match
550the present notation, and then a description of the implementation and usage.
551
552##### Equation Formulation
553
554$$
555\bm{u}(\bm{x}, t) = \bm{\overline{u}}(\bm{x}) + \bm{C}(\bm{x}) \cdot \bm{v}'
556$$
557
558$$
559\begin{aligned}
560\bm{v}' &= 2 \sqrt{3/2} \sum^N_{n=1} \sqrt{q^n(\bm{x})} \bm{\sigma}^n \cos(\kappa^n \bm{d}^n \cdot \bm{\hat{x}}^n(\bm{x}, t) + \phi^n ) \\
561\bm{\hat{x}}^n &= \left[(x - U_0 t)\max(2\kappa_{\min}/\kappa^n, 0.1) , y, z  \right]^T
562\end{aligned}
563$$
564
565Here, we define the number of wavemodes $N$, set of random numbers $ \{\bm{\sigma}^n,
566\bm{d}^n, \phi^n\}_{n=1}^N$, the Cholesky decomposition of the Reynolds stress
567tensor $\bm{C}$ (such that $\bm{R} = \bm{CC}^T$ ), bulk velocity $U_0$,
568wavemode amplitude $q^n$, wavemode frequency $\kappa^n$, and $\kappa_{\min} =
5690.5 \min_{\bm{x}} (\kappa_e)$.
570
571$$
572\kappa_e = \frac{2\pi}{\min(2d_w, 3.0 l_t)}
573$$
574
575where $l_t$ is the turbulence length scale, and $d_w$ is the distance to the
576nearest wall.
577
578
579The set of wavemode frequencies is defined by a geometric distribution:
580
581$$
582\kappa^n = \kappa_{\min} (1 + \alpha)^{n-1} \ , \quad \forall n=1, 2, ... , N
583$$
584
585The wavemode amplitudes $q^n$ are defined by a model energy spectrum $E(\kappa)$:
586
587$$
588q^n = \frac{E(\kappa^n) \Delta \kappa^n}{\sum^N_{n=1} E(\kappa^n)\Delta \kappa^n} \ ,\quad \Delta \kappa^n = \kappa^n - \kappa^{n-1}
589$$
590
591$$ E(\kappa) = \frac{(\kappa/\kappa_e)^4}{[1 + 2.4(\kappa/\kappa_e)^2]^{17/6}} f_\eta f_{\mathrm{cut}} $$
592
593$$
594f_\eta = \exp \left[-(12\kappa /\kappa_\eta)^2 \right], \quad
595f_\mathrm{cut} = \exp \left( - \left [ \frac{4\max(\kappa-0.9\kappa_\mathrm{cut}, 0)}{\kappa_\mathrm{cut}} \right]^3 \right)
596$$
597
598$\kappa_\eta$ represents turbulent dissipation frequency, and is given as $2\pi
599(\nu^3/\varepsilon)^{-1/4}$ with $\nu$ the kinematic viscosity and
600$\varepsilon$ the turbulent dissipation. $\kappa_\mathrm{cut}$ approximates the
601effective cutoff frequency of the mesh (viewing the mesh as a filter on
602solution over $\Omega$) and is given by:
603
604$$
605\kappa_\mathrm{cut} = \frac{2\pi}{ 2\min\{ [\max(h_y, h_z, 0.3h_{\max}) + 0.1 d_w], h_{\max} \} }
606$$
607
608The enforcement of the boundary condition is identical to the blasius inflow;
609it weakly enforces velocity, with the option of weakly enforcing either density
610or temperature using the the `-weakT` flag.
611
612##### Initialization Data Flow
613
614Data flow for initializing function (which creates the context data struct) is
615given below:
616```{mermaid}
617flowchart LR
618    subgraph STGInflow.dat
619    y
620    lt[l_t]
621    eps
622    Rij[R_ij]
623    ubar
624    end
625
626    subgraph STGRand.dat
627    rand[RN Set];
628    end
629
630    subgraph User Input
631    u0[U0];
632    end
633
634    subgraph init[Create Context Function]
635    ke[k_e]
636    N;
637    end
638    lt --Calc-->ke --Calc-->kn
639    y --Calc-->ke
640
641    subgraph context[Context Data]
642    yC[y]
643    randC[RN Set]
644    Cij[C_ij]
645    u0 --Copy--> u0C[U0]
646    kn[k^n];
647    ubarC[ubar]
648    ltC[l_t]
649    epsC[eps]
650    end
651    ubar --Copy--> ubarC;
652    y --Copy--> yC;
653    lt --Copy--> ltC;
654    eps --Copy--> epsC;
655
656    rand --Copy--> randC;
657    rand --> N --Calc--> kn;
658    Rij --Calc--> Cij[C_ij]
659```
660
661This is done once at runtime. The spatially-varying terms are then evaluated at
662each quadrature point on-the-fly, either by interpolation (for $l_t$,
663$\varepsilon$, $C_{ij}$, and $\overline{\bm u}$) or by calculation (for $q^n$).
664
665The `STGInflow.dat` file is a table of values at given distances from the wall.
666These values are then interpolated to a physical location (node or quadrature
667point). It has the following format:
668```
669[Total number of locations] 14
670[d_w] [u_1] [u_2] [u_3] [R_11] [R_22] [R_33] [R_12] [R_13] [R_23] [sclr_1] [sclr_2] [l_t] [eps]
671```
672where each `[  ]` item is a number in scientific notation (ie. `3.1415E0`), and `sclr_1` and
673`sclr_2` are reserved for turbulence modeling variables. They are not used in
674this example.
675
676The `STGRand.dat` file is the table of the random number set, $\{\bm{\sigma}^n,
677\bm{d}^n, \phi^n\}_{n=1}^N$. It has the format:
678```
679[Number of wavemodes] 7
680[d_1] [d_2] [d_3] [phi] [sigma_1] [sigma_2] [sigma_3]
681```
682
683The following table is presented to help clarify the dimensionality of the
684numerous terms in the STG formulation.
685
686| Math                                           | Label    | $f(\bm{x})$?   | $f(n)$?   |
687| -----------------                              | -------- | -------------- | --------- |
688| $ \{\bm{\sigma}^n, \bm{d}^n, \phi^n\}_{n=1}^N$ | RN Set   | No             | Yes       |
689| $\bm{\overline{u}}$                            | ubar     | Yes            | No        |
690| $U_0$                                          | U0       | No             | No        |
691| $l_t$                                          | l_t      | Yes            | No        |
692| $\varepsilon$                                  | eps      | Yes            | No        |
693| $\bm{R}$                                       | R_ij     | Yes            | No        |
694| $\bm{C}$                                       | C_ij     | Yes            | No        |
695| $q^n$                                          | q^n      | Yes            | Yes       |
696| $\{\kappa^n\}_{n=1}^N$                         | k^n      | No             | Yes       |
697| $h_i$                                          | h_i      | Yes            | No        |
698| $d_w$                                          | d_w      | Yes            | No        |
699
700#### Internal Damping Layer (IDL)
701The STG inflow boundary condition creates large amplitude acoustic waves.
702We use an internal damping layer (IDL) to damp them out without disrupting the synthetic structures developing into natural turbulent structures. This implementation was inspired from
703{cite}`shurSTG2014`, but is implemented here as a ramped volumetric forcing
704term, similar to a sponge layer (see 8.4.2.4 in {cite}`colonius2023turbBC` for example). It takes the following form:
705
706$$
707S(\bm{q}) = -\sigma(\bm{x})\left.\frac{\partial \bm{q}}{\partial \bm{Y}}\right\rvert_{\bm{q}} \bm{Y}'
708$$
709
710where $\bm{Y}' = [P - P_\mathrm{ref}, \bm{0}, 0]^T$, and $\sigma(\bm{x})$ is a
711linear ramp starting at `-idl_start` with length `-idl_length` and an amplitude
712of inverse `-idl_decay_rate`. The damping is defined in terms of a pressure-primitive
713anomaly $\bm Y'$ converted to conservative source using $\partial
714\bm{q}/\partial \bm{Y}\rvert_{\bm{q}}$, which is linearized about the current
715flow state. $P_\mathrm{ref}$ is defined via the `-reference_pressure` flag.
716
717### Meshing
718
719The flat plate boundary layer example has custom meshing features to better
720resolve the flow. One of those is tilting the top of the domain, allowing for
721it to be a outflow boundary condition. The angle of this tilt is controlled by
722`-platemesh_top_angle`
723
724The primary meshing feature is the ability to grade the mesh, providing better
725resolution near the wall. There are two methods to do this; algorithmically, or
726specifying the node locations via a file. Algorithmically, a base node
727distribution is defined at the inlet (assumed to be $\min(x)$) and then
728linearly stretched/squeezed to match the slanted top boundary condition. Nodes
729are placed such that `-platemesh_Ndelta` elements are within
730`-platemesh_refine_height` of the wall. They are placed such that the element
731height matches a geometric growth ratio defined by `-platemesh_growth`. The
732remaining elements are then distributed from `-platemesh_refine_height` to the
733top of the domain linearly in logarithmic space.
734
735Alternatively, a file may be specified containing the locations of each node.
736The file should be newline delimited, with the first line specifying the number
737of points and the rest being the locations of the nodes. The node locations
738used exactly at the inlet (assumed to be $\min(x)$) and linearly
739stretched/squeezed to match the slanted top boundary condition. The file is
740specified via `-platemesh_y_node_locs_path`. If this flag is given an empty
741string, then the algorithmic approach will be performed.
742