xref: /libCEED/examples/fluids/index.md (revision 8a94a473032dc6ed59a2cf0afe1d886fbdb591f4)
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
9bc7bbd5dSLeila Ghaffari## Running the mini-app
10bc7bbd5dSLeila Ghaffari
11bc7bbd5dSLeila Ghaffari```{include} README.md
12bc7bbd5dSLeila Ghaffari:start-after: inclusion-fluids-marker
13bc7bbd5dSLeila Ghaffari```
14bc7bbd5dSLeila Ghaffari## The Navier-Stokes equations
15bc7bbd5dSLeila 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}) &=
5311dee7daSJed Brown\underbrace{\begin{pmatrix}
54bcb2dfaeSJed Brown    \bm{U}\\
5511dee7daSJed Brown    {(\bm{U} \otimes \bm{U})}/{\rho} + P \bm{I}_3 \\
5611dee7daSJed Brown    {(E + P)\bm{U}}/{\rho}
5711dee7daSJed Brown\end{pmatrix}}_{\bm F_{\text{adv}}} +
5811dee7daSJed Brown\underbrace{\begin{pmatrix}
5911dee7daSJed Brown0 \\
6011dee7daSJed Brown-  \bm{\sigma} \\
6111dee7daSJed Brown - \bm{u}  \cdot \bm{\sigma} - k \nabla T
6211dee7daSJed Brown\end{pmatrix}}_{\bm F_{\text{diff}}},\\
63bcb2dfaeSJed BrownS(\bm{q}) &=
64bcb2dfaeSJed Brown- \begin{pmatrix}
65bcb2dfaeSJed Brown    0\\
66bcb2dfaeSJed Brown    \rho g \bm{\hat{k}}\\
67bcb2dfaeSJed Brown    0
68bcb2dfaeSJed Brown\end{pmatrix}.
69bcb2dfaeSJed Brown\end{aligned}
7011dee7daSJed Brown$$ (eq-ns-flux)
71bcb2dfaeSJed Brown
72bcb2dfaeSJed BrownLet the discrete solution be
73bcb2dfaeSJed Brown
74bcb2dfaeSJed Brown$$
75bcb2dfaeSJed Brown\bm{q}_N (\bm{x},t)^{(e)} = \sum_{k=1}^{P}\psi_k (\bm{x})\bm{q}_k^{(e)}
76bcb2dfaeSJed Brown$$
77bcb2dfaeSJed Brown
78bcb2dfaeSJed Brownwith $P=p+1$ the number of nodes in the element $e$.
79bcb2dfaeSJed BrownWe use tensor-product bases $\psi_{kji} = h_i(X_0)h_j(X_1)h_k(X_2)$.
80bcb2dfaeSJed Brown
81bcb2dfaeSJed BrownFor the time discretization, we use two types of time stepping schemes.
82bcb2dfaeSJed Brown
83bcb2dfaeSJed Brown- Explicit time-stepping method
84bcb2dfaeSJed Brown
85bcb2dfaeSJed 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)
86bcb2dfaeSJed Brown
87bcb2dfaeSJed Brown  $$
88bcb2dfaeSJed Brown  \bm{q}_N^{n+1} = \bm{q}_N^n + \Delta t \sum_{i=1}^{s} b_i k_i \, ,
89bcb2dfaeSJed Brown  $$
90bcb2dfaeSJed Brown
91bcb2dfaeSJed Brown  where
92bcb2dfaeSJed Brown
93bcb2dfaeSJed Brown  $$
94bcb2dfaeSJed Brown  \begin{aligned}
95bcb2dfaeSJed Brown     k_1 &= f(t^n, \bm{q}_N^n)\\
96bcb2dfaeSJed Brown     k_2 &= f(t^n + c_2 \Delta t, \bm{q}_N^n + \Delta t (a_{21} k_1))\\
97bcb2dfaeSJed Brown     k_3 &= f(t^n + c_3 \Delta t, \bm{q}_N^n + \Delta t (a_{31} k_1 + a_{32} k_2))\\
98bcb2dfaeSJed Brown     \vdots&\\
99bcb2dfaeSJed 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)\\
100bcb2dfaeSJed Brown  \end{aligned}
101bcb2dfaeSJed Brown  $$
102bcb2dfaeSJed Brown
103bcb2dfaeSJed Brown  and with
104bcb2dfaeSJed Brown
105bcb2dfaeSJed Brown  $$
106bcb2dfaeSJed Brown  f(t^n, \bm{q}_N^n) = - [\nabla \cdot \bm{F}(\bm{q}_N)]^n + [S(\bm{q}_N)]^n \, .
107bcb2dfaeSJed Brown  $$
108bcb2dfaeSJed Brown
109bcb2dfaeSJed Brown- Implicit time-stepping method
110bcb2dfaeSJed Brown
111bcb2dfaeSJed 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).
112bcb2dfaeSJed Brown  The implicit formulation solves nonlinear systems for $\bm q_N$:
113bcb2dfaeSJed Brown
114bcb2dfaeSJed Brown  $$
115bcb2dfaeSJed Brown  \bm f(\bm q_N) \equiv \bm g(t^{n+1}, \bm{q}_N, \bm{\dot{q}}_N) = 0 \, ,
116bcb2dfaeSJed Brown  $$ (eq-ts-implicit-ns)
117bcb2dfaeSJed Brown
118bcb2dfaeSJed Brown  where the time derivative $\bm{\dot q}_N$ is defined by
119bcb2dfaeSJed Brown
120bcb2dfaeSJed Brown  $$
121bcb2dfaeSJed Brown  \bm{\dot{q}}_N(\bm q_N) = \alpha \bm q_N + \bm z_N
122bcb2dfaeSJed Brown  $$
123bcb2dfaeSJed Brown
124bcb2dfaeSJed 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.).
1258791656fSJed Brown  Each nonlinear system {eq}`eq-ts-implicit-ns` will correspond to a weak form, as explained below.
1268791656fSJed Brown  In determining how difficult a given problem is to solve, we consider the Jacobian of {eq}`eq-ts-implicit-ns`,
127bcb2dfaeSJed Brown
128bcb2dfaeSJed Brown  $$
129bcb2dfaeSJed 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}.
130bcb2dfaeSJed Brown  $$
131bcb2dfaeSJed Brown
132bcb2dfaeSJed 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).
133bcb2dfaeSJed 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.
134bcb2dfaeSJed 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.
135bcb2dfaeSJed Brown
1368791656fSJed 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,
137bcb2dfaeSJed Brown
138bcb2dfaeSJed Brown$$
139bcb2dfaeSJed 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\,,
140bcb2dfaeSJed Brown$$
141bcb2dfaeSJed Brown
142bcb2dfaeSJed 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).
143bcb2dfaeSJed Brown
144bcb2dfaeSJed BrownIntegrating by parts on the divergence term, we arrive at the weak form,
145bcb2dfaeSJed Brown
146bcb2dfaeSJed Brown$$
147bcb2dfaeSJed Brown\begin{aligned}
148bcb2dfaeSJed Brown\int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
149bcb2dfaeSJed Brown- \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
150bcb2dfaeSJed Brown+ \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm q_N) \cdot \widehat{\bm{n}} \,dS
151bcb2dfaeSJed Brown  &= 0 \, , \; \forall \bm v \in \mathcal{V}_p \,,
152bcb2dfaeSJed Brown\end{aligned}
153bcb2dfaeSJed Brown$$ (eq-weak-vector-ns)
154bcb2dfaeSJed Brown
155bcb2dfaeSJed Brownwhere $\bm{F}(\bm q_N) \cdot \widehat{\bm{n}}$ is typically replaced with a boundary condition.
156bcb2dfaeSJed Brown
157bcb2dfaeSJed Brown:::{note}
158bcb2dfaeSJed 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.
159bcb2dfaeSJed Brown:::
160bcb2dfaeSJed Brown
1618791656fSJed BrownWe solve {eq}`eq-weak-vector-ns` using a Galerkin discretization (default) or a stabilized method, as is necessary for most real-world flows.
162bcb2dfaeSJed Brown
163bcb2dfaeSJed 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.
164bcb2dfaeSJed BrownOur formulation follows {cite}`hughesetal2010`, which offers a comprehensive review of stabilization and shock-capturing methods for continuous finite element discretization of compressible flows.
165bcb2dfaeSJed Brown
166bcb2dfaeSJed Brown- **SUPG** (streamline-upwind/Petrov-Galerkin)
167bcb2dfaeSJed Brown
1688791656fSJed 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`.
169bcb2dfaeSJed Brown  The weak form for this method is given as
170bcb2dfaeSJed Brown
171bcb2dfaeSJed Brown  $$
172bcb2dfaeSJed Brown  \begin{aligned}
173bcb2dfaeSJed Brown  \int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
174bcb2dfaeSJed Brown  - \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
175bcb2dfaeSJed Brown  + \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm{q}_N) \cdot \widehat{\bm{n}} \,dS & \\
17693844253SJed Brown  + \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} \, + \,
177bcb2dfaeSJed Brown  \nabla \cdot \bm{F} \, (\bm{q}_N) - \bm{S}(\bm{q}_N) \right) \,dV &= 0
178bcb2dfaeSJed Brown  \, , \; \forall \bm v \in \mathcal{V}_p
179bcb2dfaeSJed Brown  \end{aligned}
180bcb2dfaeSJed Brown  $$ (eq-weak-vector-ns-supg)
181bcb2dfaeSJed Brown
182bcb2dfaeSJed Brown  This stabilization technique can be selected using the option `-stab supg`.
183bcb2dfaeSJed Brown
184bcb2dfaeSJed Brown- **SU** (streamline-upwind)
185bcb2dfaeSJed Brown
1868791656fSJed 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
187bcb2dfaeSJed Brown
188bcb2dfaeSJed Brown  $$
189bcb2dfaeSJed Brown  \begin{aligned}
190bcb2dfaeSJed Brown  \int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right)  \,dV
191bcb2dfaeSJed Brown  - \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\
192bcb2dfaeSJed Brown  + \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm{q}_N) \cdot \widehat{\bm{n}} \,dS & \\
19393844253SJed Brown  + \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
194bcb2dfaeSJed Brown  & = 0 \, , \; \forall \bm v \in \mathcal{V}_p
195bcb2dfaeSJed Brown  \end{aligned}
196bcb2dfaeSJed Brown  $$ (eq-weak-vector-ns-su)
197bcb2dfaeSJed Brown
198bcb2dfaeSJed Brown  This stabilization technique can be selected using the option `-stab su`.
199bcb2dfaeSJed Brown
20093844253SJed BrownIn 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.
20193844253SJed BrownThe 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.
20288626eedSJames WrightThe forward variational form can be readily expressed by differentiating $\bm F_{\text{adv}}$ of {eq}`eq-ns-flux`
20311dee7daSJed Brown
20411dee7daSJed Brown$$
20511dee7daSJed Brown\begin{aligned}
20611dee7daSJed Brown\diff\bm F_{\text{adv}}(\diff\bm q; \bm q) &= \frac{\partial \bm F_{\text{adv}}}{\partial \bm q} \diff\bm q \\
20711dee7daSJed Brown&= \begin{pmatrix}
20811dee7daSJed Brown\diff\bm U \\
20911dee7daSJed Brown(\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 \\
21011dee7daSJed Brown(E + P)\diff\bm U/\rho + (\diff E + \diff P)\bm U/\rho - (E + P) \bm U/\rho^2 \diff\rho
21111dee7daSJed Brown\end{pmatrix},
21211dee7daSJed Brown\end{aligned}
21311dee7daSJed Brown$$
21411dee7daSJed Brown
21511dee7daSJed Brownwhere $\diff P$ is defined by differentiating {eq}`eq-state`.
21611dee7daSJed Brown
21711dee7daSJed Brown:::{dropdown} Stabilization scale $\bm\tau$
21811dee7daSJed BrownA 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.
21911dee7daSJed BrownTo 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)$.
22011dee7daSJed BrownSo a small normal component of velocity will be amplified (by a factor of the aspect ratio $1/\epsilon$) in this transformation.
221679c4372SJed BrownThe 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.
222d4f43295SJames WrightA 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$.
223679c4372SJed BrownWhile $\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.
224679c4372SJed BrownIf 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.
22511dee7daSJed Brown
22611dee7daSJed BrownThe 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$).
22711dee7daSJed BrownThis can be generalized to arbitrary grids by defining the local Péclet number
22811dee7daSJed Brown
22911dee7daSJed Brown$$
23011dee7daSJed Brown\mathrm{Pe} = \frac{\lVert \bm u \rVert^2}{\lVert \bm u_{\bm X} \rVert \kappa}.
23111dee7daSJed Brown$$ (eq-peclet)
23211dee7daSJed Brown
23311dee7daSJed BrownFor scalar advection-diffusion, the stabilization is a scalar
23411dee7daSJed Brown
23511dee7daSJed Brown$$
23611dee7daSJed Brown\tau = \frac{\xi(\mathrm{Pe})}{\lVert \bm u_{\bm X} \rVert},
23711dee7daSJed Brown$$ (eq-tau-advdiff)
23811dee7daSJed Brown
23911dee7daSJed Brownwhere $\xi(\mathrm{Pe}) = \coth \mathrm{Pe} - 1/\mathrm{Pe}$ approaches 1 at large local Péclet number.
24011dee7daSJed BrownNote 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.
24193844253SJed BrownFor advection-diffusion, $\bm F(q) = \bm u q$, and thus the SU stabilization term is
24211dee7daSJed Brown
24311dee7daSJed Brown$$
24493844253SJed Brown\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 .
24593844253SJed Brown$$ (eq-su-stabilize-advdiff)
24611dee7daSJed Brown
24793844253SJed Brownwhere the term in parentheses is a rank-1 diffusivity tensor that has been pulled back to the reference element.
24811dee7daSJed BrownSee {cite}`hughesetal2010` equations 15-17 and 34-36 for further discussion of this formulation.
24911dee7daSJed Brown
25088626eedSJames WrightFor 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
25111dee7daSJed Brown1. continuity stabilization $\tau_c$
25211dee7daSJed Brown2. momentum stabilization $\tau_m$
25311dee7daSJed Brown3. energy stabilization $\tau_E$
25411dee7daSJed Brown
25588626eedSJames WrightThe Navier-Stokes code in this example uses the following formulation for $\tau_c$, $\tau_m$, $\tau_E$:
25688626eedSJames Wright
25788626eedSJames Wright$$
25888626eedSJames Wright\begin{aligned}
25988626eedSJames Wright
26088626eedSJames Wright\tau_c &= \frac{C_c \mathcal{F}}{8\rho \trace(\bm g)} \\
26188626eedSJames Wright\tau_m &= \frac{C_m}{\mathcal{F}} \\
26288626eedSJames Wright\tau_E &= \frac{C_E}{\mathcal{F} c_v} \\
26388626eedSJames Wright\end{aligned}
26488626eedSJames Wright$$
26588626eedSJames Wright
26688626eedSJames Wright$$
26788626eedSJames Wright\mathcal{F} = \sqrt{ \rho^2 \left [ \left(\frac{2C_t}{\Delta t}\right)^2
26888626eedSJames Wright+ \bm u \cdot (\bm u \cdot  \bm g)
26988626eedSJames Wright+ C_v \mu^2 \Vert \bm g \Vert_F ^2\right]}
27088626eedSJames Wright$$
27188626eedSJames Wright
27288626eedSJames Wrightwhere $\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.
27388626eedSJames WrightThis formulation is currently not available in the Euler code.
27488626eedSJames Wright
27588626eedSJames WrightIn the Euler code, we follow {cite}`hughesetal2010` in defining a $3\times 3$ diagonal stabilization according to spatial criterion 2 (equation 27) as follows.
276c94bf672SLeila Ghaffari
277c94bf672SLeila Ghaffari$$
278679c4372SJed Brown\tau_{ii} = c_{\tau} \frac{2 \xi(\mathrm{Pe})}{(\lambda_{\max \text{abs}})_i \lVert \nabla_{x_i} \bm X \rVert}
279c94bf672SLeila Ghaffari$$ (eq-tau-conservative)
280c94bf672SLeila Ghaffari
281679c4372SJed Brownwhere $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$.
282679c4372SJed BrownThe 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.
283679c4372SJed BrownThe complete set of eigenvalues of the Euler flux Jacobian in direction $i$ are (e.g., {cite}`toro2009`)
284c94bf672SLeila Ghaffari
285c94bf672SLeila Ghaffari$$
286679c4372SJed Brown\Lambda_i = [u_i - a, u_i, u_i, u_i, u_i+a],
287c94bf672SLeila Ghaffari$$ (eq-eigval-advdiff)
288c94bf672SLeila Ghaffari
289679c4372SJed Brownwhere $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.
290679c4372SJed BrownNote 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.
291679c4372SJed BrownThe fastest wave speed in direction $i$ is thus
292c94bf672SLeila Ghaffari
293c94bf672SLeila Ghaffari$$
294679c4372SJed Brown\lambda_{\max \text{abs}} \Bigl( \frac{\partial \bm F_{\text{adv}}}{\partial \bm q} \cdot \hat{\bm n}_i \Bigr) = |u_i| + a
295c94bf672SLeila Ghaffari$$ (eq-wavespeed)
296c94bf672SLeila Ghaffari
297679c4372SJed BrownNote 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.
298c94bf672SLeila Ghaffari
29911dee7daSJed Brown:::
300bcb2dfaeSJed Brown
301bcb2dfaeSJed BrownCurrently, this demo provides three types of problems/physical models that can be selected at run time via the option `-problem`.
302bcb2dfaeSJed 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.
303bcb2dfaeSJed Brown
304bcb2dfaeSJed Brown(problem-advection)=
305bcb2dfaeSJed Brown
306bcb2dfaeSJed Brown## Advection
307bcb2dfaeSJed Brown
3088791656fSJed BrownA simplified version of system {eq}`eq-ns`, only accounting for the transport of total energy, is given by
309bcb2dfaeSJed Brown
310bcb2dfaeSJed Brown$$
311bcb2dfaeSJed Brown\frac{\partial E}{\partial t} + \nabla \cdot (\bm{u} E ) = 0 \, ,
312bcb2dfaeSJed Brown$$ (eq-advection)
313bcb2dfaeSJed Brown
314bcb2dfaeSJed 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.
315bcb2dfaeSJed Brown
316bcb2dfaeSJed Brown- **Rotation**
317bcb2dfaeSJed Brown
318bcb2dfaeSJed Brown  In this case, a uniform circular velocity field transports the blob of total energy.
3198791656fSJed Brown  We have solved {eq}`eq-advection` applying zero energy density $E$, and no-flux for $\bm{u}$ on the boundaries.
320bcb2dfaeSJed Brown
321bcb2dfaeSJed Brown- **Translation**
322bcb2dfaeSJed Brown
323bcb2dfaeSJed 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.
324bcb2dfaeSJed Brown
3258791656fSJed 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
326bcb2dfaeSJed Brown
327bcb2dfaeSJed Brown  $$
328bcb2dfaeSJed 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  \, ,
329bcb2dfaeSJed Brown  $$
330bcb2dfaeSJed Brown
331bcb2dfaeSJed 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.
3328791656fSJed Brown  The weak form boundary integral in {eq}`eq-weak-vector-ns` for outflow boundary conditions is defined as
333bcb2dfaeSJed Brown
334bcb2dfaeSJed Brown  $$
335bcb2dfaeSJed 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  \, ,
336bcb2dfaeSJed Brown  $$
337bcb2dfaeSJed Brown
338bcb2dfaeSJed Brown(problem-euler-vortex)=
339bcb2dfaeSJed Brown
340bcb2dfaeSJed Brown## Isentropic Vortex
341bcb2dfaeSJed Brown
342bc7bbd5dSLeila 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
343bcb2dfaeSJed Brown
344bcb2dfaeSJed Brown$$
345bcb2dfaeSJed Brown\begin{aligned}
346bcb2dfaeSJed Brown\frac{\partial \rho}{\partial t} + \nabla \cdot \bm{U} &= 0 \\
347bcb2dfaeSJed Brown\frac{\partial \bm{U}}{\partial t} + \nabla \cdot \left( \frac{\bm{U} \otimes \bm{U}}{\rho} + P \bm{I}_3 \right) &= 0 \\
348bcb2dfaeSJed Brown\frac{\partial E}{\partial t} + \nabla \cdot \left( \frac{(E + P)\bm{U}}{\rho} \right) &= 0 \, , \\
349bcb2dfaeSJed Brown\end{aligned}
350bcb2dfaeSJed Brown$$ (eq-euler)
351bcb2dfaeSJed Brown
352bc7bbd5dSLeila 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
353bcb2dfaeSJed Brown
354bcb2dfaeSJed Brown$$
355bcb2dfaeSJed 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}
356bcb2dfaeSJed Brown$$
357bcb2dfaeSJed Brown
358bc7bbd5dSLeila 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).
359bcb2dfaeSJed BrownThere is no perturbation in the entropy $S=P/\rho^\gamma$ ($\delta S=0)$.
360bcb2dfaeSJed Brown
361019b7682STimothy Aiken(problem-shock-tube)=
362019b7682STimothy Aiken
363019b7682STimothy Aiken## Shock Tube
364019b7682STimothy Aiken
365019b7682STimothy AikenThis 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.
366019b7682STimothy Aiken
367019b7682STimothy AikenSU 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
368019b7682STimothy Aiken
369019b7682STimothy Aiken$$
370019b7682STimothy Aiken\int_{\Omega} \nu_{SHOCK} \nabla \bm v \!:\! \nabla \bm q dV
371019b7682STimothy Aiken$$
372019b7682STimothy Aiken
373019b7682STimothy AikenThe 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
374019b7682STimothy Aiken
375019b7682STimothy Aiken$$
376019b7682STimothy Aiken\nu_{SHOCK} = \tau_{SHOCK} u_{cha}^2
377019b7682STimothy Aiken$$
378ba6664aeSJames Wright
379019b7682STimothy Aikenwhere,
380ba6664aeSJames Wright
381019b7682STimothy Aiken$$
382019b7682STimothy Aiken\tau_{SHOCK} = \frac{h_{SHOCK}}{2u_{cha}} \left( \frac{ \,|\, \nabla \rho \,|\, h_{SHOCK}}{\rho_{ref}} \right)^{\beta}
383019b7682STimothy Aiken$$
384019b7682STimothy Aiken
385ba6664aeSJames Wright$\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
386019b7682STimothy Aiken
387019b7682STimothy Aiken$$
388019b7682STimothy Aikenh_{SHOCK} = 2 \left( C_{YZB} \,|\, \bm p \,|\, \right)^{-1}
389019b7682STimothy Aiken$$
390ba6664aeSJames Wright
391019b7682STimothy Aikenwhere
392ba6664aeSJames Wright
393019b7682STimothy Aiken$$
394019b7682STimothy Aikenp_k = \hat{j}_i \frac{\partial \xi_i}{x_k}
395019b7682STimothy Aiken$$
396019b7682STimothy Aiken
397019b7682STimothy AikenThe 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.
398019b7682STimothy Aiken
399bcb2dfaeSJed Brown(problem-density-current)=
4007ec884f8SJames Wright
4017ec884f8SJames Wright## Newtonian Wave
4027ec884f8SJames WrightThis 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.
4037ec884f8SJames Wright
4047ec884f8SJames WrightThe problem has a perturbed initial condition and lets it evolve in time. The initial condition contains a Gaussian perturbation in the pressure field:
4057ec884f8SJames Wright
4067ec884f8SJames Wright$$
4077ec884f8SJames Wright\begin{aligned}
4087ec884f8SJames Wright\rho &= \rho_\infty\left(1+A\exp\left(\frac{-(\bar{x}^2 + \bar{y}^2)}{2\sigma^2}\right)\right) \\
4097ec884f8SJames Wright\bm{U} &= \bm U_\infty \\
4107ec884f8SJames WrightE &= \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},
4117ec884f8SJames Wright\end{aligned}
4127ec884f8SJames Wright$$
4137ec884f8SJames Wright
4147ec884f8SJames Wrightwhere $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)$.
415f1e435c9SJed BrownThe simulation produces a strong acoustic wave and leaves behind a cold thermal bubble that advects at the fluid velocity.
4167ec884f8SJames Wright
417f1e435c9SJed BrownThe 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.
418f1e435c9SJed BrownThis 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.
419d310b3d3SAdeleke O. Bankole
420d310b3d3SAdeleke O. Bankole## Vortex Shedding - Flow past Cylinder
421d310b3d3SAdeleke O. BankoleThis test case, based on {cite}`shakib1991femcfd`, provides an example of using an externally provide mesh from Gmsh.
422d310b3d3SAdeleke O. BankoleA 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$. We solve this as a 3D problem problem with (default) one element in the $z$ direction. The domain is filled with an ideal gas at rest (zero velocity) with temperature 24.92, pressure 7143, and viscosity 0.01. At time $t=0$, it is subjected to freestream boundary conditions at the inflow (left) and outflow (right), with exterior reference state at velocity $(1, 0, 0)$ giving Reynolds number $100$ and Mach number $0.01$. A symmetry (adiabatic free slip) condition is imposed at the top and bottom boundaries $(y = \pm 4.5)$ (zero normal velocity component, zero heat-flux). On the cylinder wall, a no-slip boundary condition, no heat-flux condition are imposed. As we evolve in time, eddies appear past the cylinder leading to a vortex shedding which is known as the vortex street.
423d310b3d3SAdeleke O. Bankole
424d310b3d3SAdeleke O. BankoleThe Gmsh input file, `examples/fluids/meshes/cylinder.geo` is parametrized to facilitate experimenting with similar configurations. The Strouhal number (nondimensional shedding frequency) is sensitive to the size of the computational domain and boundary conditions.
425bcb2dfaeSJed Brown
426bcb2dfaeSJed Brown## Density Current
427bcb2dfaeSJed Brown
4288791656fSJed 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.
429bcb2dfaeSJed 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
430bcb2dfaeSJed Brown
431bcb2dfaeSJed Brown$$
432bcb2dfaeSJed 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}
433bcb2dfaeSJed Brown$$
434bcb2dfaeSJed Brown
435bcb2dfaeSJed Brownwhere $P_0$ is the atmospheric pressure.
436bcb2dfaeSJed BrownFor this problem, we have used no-slip and non-penetration boundary conditions for $\bm{u}$, and no-flux for mass and energy densities.
43788626eedSJames Wright
43888626eedSJames Wright## Channel
43988626eedSJames Wright
44088626eedSJames WrightA compressible channel flow. Analytical solution given in
44188626eedSJames Wright{cite}`whitingStabilizedFEM1999`:
44288626eedSJames Wright
44388626eedSJames Wright$$ u_1 = u_{\max} \left [ 1 - \left ( \frac{x_2}{H}\right)^2 \right] \quad \quad u_2 = u_3 = 0$$
44488626eedSJames Wright$$T = T_w \left [ 1 + \frac{Pr \hat{E}c}{3} \left \{1 - \left(\frac{x_2}{H}\right)^4  \right \} \right]$$
44588626eedSJames Wright$$p = p_0 - \frac{2\rho_0 u_{\max}^2 x_1}{Re_H H}$$
44688626eedSJames Wright
44788626eedSJames Wrightwhere $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.
44888626eedSJames Wright
44988626eedSJames WrightBoundary conditions are periodic in the streamwise direction, and no-slip and non-penetration boundary conditions at the walls.
450a1df05f8SJed BrownThe flow is driven by a body force determined analytically from the fluid properties and setup parameters $H$ and $u_{\max}$.
45188626eedSJames Wright
452ba6664aeSJames Wright## Flat Plate Boundary Layer
453ba6664aeSJames Wright
454ba6664aeSJames Wright### Laminar Boundary Layer - Blasius
45588626eedSJames Wright
45688626eedSJames WrightSimulation of a laminar boundary layer flow, with the inflow being prescribed
45788626eedSJames Wrightby a [Blasius similarity
45888626eedSJames Wrightsolution](https://en.wikipedia.org/wiki/Blasius_boundary_layer). At the inflow,
459ba6664aeSJames Wrightthe velocity is prescribed by the Blasius soution profile, density is set
460ba6664aeSJames Wrightconstant, and temperature is allowed to float. Using `weakT: true`, density is
461ba6664aeSJames Wrightallowed to float and temperature is set constant. At the outlet, a user-set
462ba6664aeSJames Wrightpressure is used for pressure in the inviscid flux terms (all other inviscid
463520dae65SJames Wrightflux terms use interior solution values). The wall is a no-slip,
464520dae65SJames Wrightno-penetration, no-heat flux condition. The top of the domain is treated as an
465520dae65SJames Wrightoutflow and is tilted at a downward angle to ensure that flow is always exiting
466520dae65SJames Wrightit.
46788626eedSJames Wright
468ba6664aeSJames Wright### Turbulent Boundary Layer
469ba6664aeSJames Wright
470ba6664aeSJames WrightSimulating a turbulent boundary layer without modeling the turbulence requires
471ba6664aeSJames Wrightresolving the turbulent flow structures. These structures may be introduced
472ba6664aeSJames Wrightinto the simulations either by allowing a laminar boundary layer naturally
473ba6664aeSJames Wrighttransition to turbulence, or imposing turbulent structures at the inflow. The
474ba6664aeSJames Wrightlatter approach has been taken here, specifically using a *synthetic turbulence
475ba6664aeSJames Wrightgeneration* (STG) method.
476ba6664aeSJames Wright
477ba6664aeSJames Wright#### Synthetic Turbulence Generation (STG) Boundary Condition
478ba6664aeSJames Wright
479ba6664aeSJames WrightWe use the STG method described in
480ba6664aeSJames Wright{cite}`shurSTG2014`. Below follows a re-description of the formulation to match
481ba6664aeSJames Wrightthe present notation, and then a description of the implementation and usage.
482ba6664aeSJames Wright
483ba6664aeSJames Wright##### Equation Formulation
484ba6664aeSJames Wright
485ba6664aeSJames Wright$$
486ba6664aeSJames Wright\bm{u}(\bm{x}, t) = \bm{\overline{u}}(\bm{x}) + \bm{C}(\bm{x}) \cdot \bm{v}'
487ba6664aeSJames Wright$$
488ba6664aeSJames Wright
489ba6664aeSJames Wright$$
490ba6664aeSJames Wright\begin{aligned}
491ba6664aeSJames Wright\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 ) \\
492ba6664aeSJames Wright\bm{\hat{x}}^n &= \left[(x - U_0 t)\max(2\kappa_{\min}/\kappa^n, 0.1) , y, z  \right]^T
493ba6664aeSJames Wright\end{aligned}
494ba6664aeSJames Wright$$
495ba6664aeSJames Wright
496ba6664aeSJames WrightHere, we define the number of wavemodes $N$, set of random numbers $ \{\bm{\sigma}^n,
497ba6664aeSJames Wright\bm{d}^n, \phi^n\}_{n=1}^N$, the Cholesky decomposition of the Reynolds stress
498ba6664aeSJames Wrighttensor $\bm{C}$ (such that $\bm{R} = \bm{CC}^T$ ), bulk velocity $U_0$,
499ba6664aeSJames Wrightwavemode amplitude $q^n$, wavemode frequency $\kappa^n$, and $\kappa_{\min} =
500ba6664aeSJames Wright0.5 \min_{\bm{x}} (\kappa_e)$.
501ba6664aeSJames Wright
502ba6664aeSJames Wright$$
503ba6664aeSJames Wright\kappa_e = \frac{2\pi}{\min(2d_w, 3.0 l_t)}
504ba6664aeSJames Wright$$
505ba6664aeSJames Wright
506ba6664aeSJames Wrightwhere $l_t$ is the turbulence length scale, and $d_w$ is the distance to the
507ba6664aeSJames Wrightnearest wall.
508ba6664aeSJames Wright
509ba6664aeSJames Wright
510ba6664aeSJames WrightThe set of wavemode frequencies is defined by a geometric distribution:
511ba6664aeSJames Wright
512ba6664aeSJames Wright$$
513ba6664aeSJames Wright\kappa^n = \kappa_{\min} (1 + \alpha)^{n-1} \ , \quad \forall n=1, 2, ... , N
514ba6664aeSJames Wright$$
515ba6664aeSJames Wright
516ba6664aeSJames WrightThe wavemode amplitudes $q^n$ are defined by a model energy spectrum $E(\kappa)$:
517ba6664aeSJames Wright
518ba6664aeSJames Wright$$
519ba6664aeSJames Wrightq^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}
520ba6664aeSJames Wright$$
521ba6664aeSJames Wright
522ba6664aeSJames Wright$$ E(\kappa) = \frac{(\kappa/\kappa_e)^4}{[1 + 2.4(\kappa/\kappa_e)^2]^{17/6}} f_\eta f_{\mathrm{cut}} $$
523ba6664aeSJames Wright
524ba6664aeSJames Wright$$
525ba6664aeSJames Wrightf_\eta = \exp \left[-(12\kappa /\kappa_\eta)^2 \right], \quad
526ba6664aeSJames Wrightf_\mathrm{cut} = \exp \left( - \left [ \frac{4\max(\kappa-0.9\kappa_\mathrm{cut}, 0)}{\kappa_\mathrm{cut}} \right]^3 \right)
527ba6664aeSJames Wright$$
528ba6664aeSJames Wright
529ba6664aeSJames Wright$\kappa_\eta$ represents turbulent dissipation frequency, and is given as $2\pi
530ba6664aeSJames Wright(\nu^3/\varepsilon)^{-1/4}$ with $\nu$ the kinematic viscosity and
531ba6664aeSJames Wright$\varepsilon$ the turbulent dissipation. $\kappa_\mathrm{cut}$ approximates the
532ba6664aeSJames Wrighteffective cutoff frequency of the mesh (viewing the mesh as a filter on
533ba6664aeSJames Wrightsolution over $\Omega$) and is given by:
534ba6664aeSJames Wright
535ba6664aeSJames Wright$$
536ba6664aeSJames Wright\kappa_\mathrm{cut} = \frac{2\pi}{ 2\min\{ [\max(h_y, h_z, 0.3h_{\max}) + 0.1 d_w], h_{\max} \} }
537ba6664aeSJames Wright$$
538ba6664aeSJames Wright
539ba6664aeSJames WrightThe enforcement of the boundary condition is identical to the blasius inflow;
540ba6664aeSJames Wrightit weakly enforces velocity, with the option of weakly enforcing either density
541ba6664aeSJames Wrightor temperature using the the `-weakT` flag.
542ba6664aeSJames Wright
543ba6664aeSJames Wright##### Initialization Data Flow
544ba6664aeSJames Wright
545ba6664aeSJames WrightData flow for initializing function (which creates the context data struct) is
546ba6664aeSJames Wrightgiven below:
547ba6664aeSJames Wright```{mermaid}
548ba6664aeSJames Wrightflowchart LR
549ba6664aeSJames Wright    subgraph STGInflow.dat
550ba6664aeSJames Wright    y
551ba6664aeSJames Wright    lt[l_t]
552ba6664aeSJames Wright    eps
553ba6664aeSJames Wright    Rij[R_ij]
554ba6664aeSJames Wright    ubar
555ba6664aeSJames Wright    end
556ba6664aeSJames Wright
557ba6664aeSJames Wright    subgraph STGRand.dat
558ba6664aeSJames Wright    rand[RN Set];
559ba6664aeSJames Wright    end
560ba6664aeSJames Wright
561ba6664aeSJames Wright    subgraph User Input
562ba6664aeSJames Wright    u0[U0];
563ba6664aeSJames Wright    end
564ba6664aeSJames Wright
565ba6664aeSJames Wright    subgraph init[Create Context Function]
566ba6664aeSJames Wright    ke[k_e]
567ba6664aeSJames Wright    N;
568ba6664aeSJames Wright    end
569ba6664aeSJames Wright    lt --Calc-->ke --Calc-->kn
570ba6664aeSJames Wright    y --Calc-->ke
571ba6664aeSJames Wright
572ba6664aeSJames Wright    subgraph context[Context Data]
573ba6664aeSJames Wright    yC[y]
574ba6664aeSJames Wright    randC[RN Set]
575ba6664aeSJames Wright    Cij[C_ij]
576ba6664aeSJames Wright    u0 --Copy--> u0C[U0]
577ba6664aeSJames Wright    kn[k^n];
578ba6664aeSJames Wright    ubarC[ubar]
579ba6664aeSJames Wright    ltC[l_t]
580ba6664aeSJames Wright    epsC[eps]
581ba6664aeSJames Wright    end
582ba6664aeSJames Wright    ubar --Copy--> ubarC;
583ba6664aeSJames Wright    y --Copy--> yC;
584ba6664aeSJames Wright    lt --Copy--> ltC;
585ba6664aeSJames Wright    eps --Copy--> epsC;
586ba6664aeSJames Wright
587ba6664aeSJames Wright    rand --Copy--> randC;
588ba6664aeSJames Wright    rand --> N --Calc--> kn;
589ba6664aeSJames Wright    Rij --Calc--> Cij[C_ij]
590ba6664aeSJames Wright```
591ba6664aeSJames Wright
592ba6664aeSJames WrightThis is done once at runtime. The spatially-varying terms are then evaluated at
593ba6664aeSJames Wrighteach quadrature point on-the-fly, either by interpolation (for $l_t$,
594ba6664aeSJames Wright$\varepsilon$, $C_{ij}$, and $\overline{\bm u}$) or by calculation (for $q^n$).
595ba6664aeSJames Wright
596ba6664aeSJames WrightThe `STGInflow.dat` file is a table of values at given distances from the wall.
597ba6664aeSJames WrightThese values are then interpolated to a physical location (node or quadrature
598ba6664aeSJames Wrightpoint). It has the following format:
599ba6664aeSJames Wright```
600ba6664aeSJames Wright[Total number of locations] 14
601ba6664aeSJames Wright[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]
602ba6664aeSJames Wright```
603ba6664aeSJames Wrightwhere each `[  ]` item is a number in scientific notation (ie. `3.1415E0`), and `sclr_1` and
604ba6664aeSJames Wright`sclr_2` are reserved for turbulence modeling variables. They are not used in
605ba6664aeSJames Wrightthis example.
606ba6664aeSJames Wright
607ba6664aeSJames WrightThe `STGRand.dat` file is the table of the random number set, $\{\bm{\sigma}^n,
608ba6664aeSJames Wright\bm{d}^n, \phi^n\}_{n=1}^N$. It has the format:
609ba6664aeSJames Wright```
610ba6664aeSJames Wright[Number of wavemodes] 7
611ba6664aeSJames Wright[d_1] [d_2] [d_3] [phi] [sigma_1] [sigma_2] [sigma_3]
612ba6664aeSJames Wright```
613ba6664aeSJames Wright
614ba6664aeSJames WrightThe following table is presented to help clarify the dimensionality of the
615ba6664aeSJames Wrightnumerous terms in the STG formulation.
616ba6664aeSJames Wright
617ba6664aeSJames Wright| Math            | Label  | $f(\bm{x})$? | $f(n)$? |
618ba6664aeSJames Wright|-----------------|--------|--------------|---------|
619ba6664aeSJames Wright| $ \{\bm{\sigma}^n, \bm{d}^n, \phi^n\}_{n=1}^N$        | RN Set | No           | Yes     |
620ba6664aeSJames Wright| $\bm{\overline{u}}$ | ubar | Yes | No |
621ba6664aeSJames Wright| $U_0$           | U0     | No           | No      |
622ba6664aeSJames Wright| $l_t$           | l_t    | Yes          | No   |
623ba6664aeSJames Wright| $\varepsilon$   | eps    | Yes          | No   |
624ba6664aeSJames Wright| $\bm{R}$        | R_ij   | Yes          | No      |
625ba6664aeSJames Wright| $\bm{C}$        | C_ij   | Yes          | No      |
626ba6664aeSJames Wright| $q^n$           | q^n    | Yes           | Yes     |
627ba6664aeSJames Wright| $\{\kappa^n\}_{n=1}^N$ | k^n  | No           | Yes      |
628ba6664aeSJames Wright| $h_i$           | h_i    | Yes          | No   |
629ba6664aeSJames Wright| $d_w$           | d_w    | Yes          | No   |
63091eaef80SJames Wright
63191eaef80SJames Wright### Meshing
63291eaef80SJames Wright
63391eaef80SJames WrightThe flat plate boundary layer example has custom meshing features to better
63491eaef80SJames Wrightresolve the flow. One of those is tilting the top of the domain, allowing for
635*8a94a473SJed Brownit to be a outflow boundary condition. The angle of this tilt is controlled by
63691eaef80SJames Wright`-platemesh_top_angle`
63791eaef80SJames Wright
63891eaef80SJames WrightThe primary meshing feature is the ability to grade the mesh, providing better
63991eaef80SJames Wrightresolution near the wall. There are two methods to do this; algorithmically, or
64091eaef80SJames Wrightspecifying the node locations via a file. Algorithmically, a base node
64191eaef80SJames Wrightdistribution is defined at the inlet (assumed to be $\min(x)$) and then
64291eaef80SJames Wrightlinearly stretched/squeezed to match the slanted top boundary condition. Nodes
64391eaef80SJames Wrightare placed such that `-platemesh_Ndelta` elements are within
64491eaef80SJames Wright`-platemesh_refine_height` of the wall. They are placed such that the element
64591eaef80SJames Wrightheight matches a geometric growth ratio defined by `-platemesh_growth`. The
64691eaef80SJames Wrightremaining elements are then distributed from `-platemesh_refine_height` to the
64791eaef80SJames Wrighttop of the domain linearly in logarithmic space.
64891eaef80SJames Wright
64991eaef80SJames WrightAlternatively, a file may be specified containing the locations of each node.
65091eaef80SJames WrightThe file should be newline delimited, with the first line specifying the number
65191eaef80SJames Wrightof points and the rest being the locations of the nodes. The node locations
65291eaef80SJames Wrightused exactly at the inlet (assumed to be $\min(x)$) and linearly
65391eaef80SJames Wrightstretched/squeezed to match the slanted top boundary condition. The file is
65491eaef80SJames Wrightspecified via `-platemesh_y_node_locs_path`. If this flag is given an empty
65591eaef80SJames Wrightstring, then the algorithmic approach will be performed.
656