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 72Let the discrete solution be 73 74$$ 75\bm{q}_N (\bm{x},t)^{(e)} = \sum_{k=1}^{P}\psi_k (\bm{x})\bm{q}_k^{(e)} 76$$ 77 78with $P=p+1$ the number of nodes in the element $e$. 79We use tensor-product bases $\psi_{kji} = h_i(X_0)h_j(X_1)h_k(X_2)$. 80 81For the time discretization, we use two types of time stepping schemes. 82 83- Explicit time-stepping method 84 85 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) 86 87 $$ 88 \bm{q}_N^{n+1} = \bm{q}_N^n + \Delta t \sum_{i=1}^{s} b_i k_i \, , 89 $$ 90 91 where 92 93 $$ 94 \begin{aligned} 95 k_1 &= f(t^n, \bm{q}_N^n)\\ 96 k_2 &= f(t^n + c_2 \Delta t, \bm{q}_N^n + \Delta t (a_{21} k_1))\\ 97 k_3 &= f(t^n + c_3 \Delta t, \bm{q}_N^n + \Delta t (a_{31} k_1 + a_{32} k_2))\\ 98 \vdots&\\ 99 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)\\ 100 \end{aligned} 101 $$ 102 103 and with 104 105 $$ 106 f(t^n, \bm{q}_N^n) = - [\nabla \cdot \bm{F}(\bm{q}_N)]^n + [S(\bm{q}_N)]^n \, . 107 $$ 108 109- Implicit time-stepping method 110 111 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). 112 The implicit formulation solves nonlinear systems for $\bm q_N$: 113 114 $$ 115 \bm f(\bm q_N) \equiv \bm g(t^{n+1}, \bm{q}_N, \bm{\dot{q}}_N) = 0 \, , 116 $$ (eq-ts-implicit-ns) 117 118 where the time derivative $\bm{\dot q}_N$ is defined by 119 120 $$ 121 \bm{\dot{q}}_N(\bm q_N) = \alpha \bm q_N + \bm z_N 122 $$ 123 124 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.). 125 Each nonlinear system {eq}`eq-ts-implicit-ns` will correspond to a weak form, as explained below. 126 In determining how difficult a given problem is to solve, we consider the Jacobian of {eq}`eq-ts-implicit-ns`, 127 128 $$ 129 \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}. 130 $$ 131 132 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). 133 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. 134 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. 135 136To 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, 137 138$$ 139\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\,, 140$$ 141 142with $\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). 143 144Integrating by parts on the divergence term, we arrive at the weak form, 145 146$$ 147\begin{aligned} 148\int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right) \,dV 149- \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\ 150+ \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm q_N) \cdot \widehat{\bm{n}} \,dS 151 &= 0 \, , \; \forall \bm v \in \mathcal{V}_p \,, 152\end{aligned} 153$$ (eq-weak-vector-ns) 154 155where $\bm{F}(\bm q_N) \cdot \widehat{\bm{n}}$ is typically replaced with a boundary condition. 156 157:::{note} 158The 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. 159::: 160 161We solve {eq}`eq-weak-vector-ns` using a Galerkin discretization (default) or a stabilized method, as is necessary for most real-world flows. 162 163Galerkin 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. 164Our formulation follows {cite}`hughesetal2010`, which offers a comprehensive review of stabilization and shock-capturing methods for continuous finite element discretization of compressible flows. 165 166- **SUPG** (streamline-upwind/Petrov-Galerkin) 167 168 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`. 169 The weak form for this method is given as 170 171 $$ 172 \begin{aligned} 173 \int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right) \,dV 174 - \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\ 175 + \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm{q}_N) \cdot \widehat{\bm{n}} \,dS & \\ 176 + \int_{\Omega} \mathcal{P}(\bm v)^T \, \left( \frac{\partial \bm{q}_N}{\partial t} \, + \, 177 \nabla \cdot \bm{F} \, (\bm{q}_N) - \bm{S}(\bm{q}_N) \right) \,dV &= 0 178 \, , \; \forall \bm v \in \mathcal{V}_p 179 \end{aligned} 180 $$ (eq-weak-vector-ns-supg) 181 182 This stabilization technique can be selected using the option `-stab supg`. 183 184- **SU** (streamline-upwind) 185 186 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 187 188 $$ 189 \begin{aligned} 190 \int_{\Omega} \bm v \cdot \left( \frac{\partial \bm{q}_N}{\partial t} - \bm{S}(\bm{q}_N) \right) \,dV 191 - \int_{\Omega} \nabla \bm v \!:\! \bm{F}(\bm{q}_N)\,dV & \\ 192 + \int_{\partial \Omega} \bm v \cdot \bm{F}(\bm{q}_N) \cdot \widehat{\bm{n}} \,dS & \\ 193 + \int_{\Omega} \mathcal{P}(\bm v)^T \, \nabla \cdot \bm{F} \, (\bm{q}_N) \,dV 194 & = 0 \, , \; \forall \bm v \in \mathcal{V}_p 195 \end{aligned} 196 $$ (eq-weak-vector-ns-su) 197 198 This stabilization technique can be selected using the option `-stab su`. 199 200In both {eq}`eq-weak-vector-ns-su` and {eq}`eq-weak-vector-ns-supg`, $\mathcal P$ is called the *perturbation to the test-function space*, since it modifies the original Galerkin method into *SUPG* or *SU* schemes. 201It is defined as 202 203$$ 204\mathcal P(\bm v) \equiv \bm{\tau} \left(\frac{\partial \bm{F}_{\text{adv}} (\bm{q}_N)}{\partial \bm{q}_N} \right) \, \nabla \bm v\,, 205$$ (eq-streamline-P) 206 207where parameter $\bm{\tau} \in \mathbb R^{3}$ (spatial index) or $\bm \tau \in \mathbb R^{5\times 5}$ (field indices) is an intrinsic time scale matrix. 208Most generally, we consider $\bm\tau \in \mathbb R^{3,5,5}$. 209This expression contains the advective flux Jacobian, which may be thought of as mapping from a 5-vector (state) to a $(5,3)$ tensor (flux) or from a $(5,3)$ tensor (gradient of state) to a 5-vector (time derivative of state); the latter is used in {eq}`eq-streamline-P` because it's applied to $\nabla\bm v$. 210The forward variational form can be readily expressed by differentiating $\bm F_{\text{adv}}$ of {eq}`eq-ns-flux` 211 212$$ 213\begin{aligned} 214\diff\bm F_{\text{adv}}(\diff\bm q; \bm q) &= \frac{\partial \bm F_{\text{adv}}}{\partial \bm q} \diff\bm q \\ 215&= \begin{pmatrix} 216\diff\bm U \\ 217(\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 \\ 218(E + P)\diff\bm U/\rho + (\diff E + \diff P)\bm U/\rho - (E + P) \bm U/\rho^2 \diff\rho 219\end{pmatrix}, 220\end{aligned} 221$$ 222 223where $\diff P$ is defined by differentiating {eq}`eq-state`. 224This action is also readily computed by forward-mode AD, but since $\bm v$ is a test function, we actually need the action of the adjoint to use {eq}`eq-streamline-P` in finite element computation; that can be computed by reverse-mode AD. 225We may equivalently write the stabilization term as 226 227$$ 228\mathcal P(\bm v)^T \bm r = \nabla \bm v \tcolon \left(\frac{\partial \bm F_{\text{adv}}}{\partial \bm q}\right)^T \, \bm\tau \bm r, 229$$ 230 231where $\bm r$ is the strong form residual and $\bm\tau$ is a $5\times 5$ matrix. 232 233:::{dropdown} Stabilization scale $\bm\tau$ 234A 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. 235To 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)$. 236So a small normal component of velocity will be amplified (by a factor of the aspect ratio $1/\epsilon$) in this transformation. 237The 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. 238A 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$. 239While $\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. 240If 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. 241 242The 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$). 243This can be generalized to arbitrary grids by defining the local Péclet number 244 245$$ 246\mathrm{Pe} = \frac{\lVert \bm u \rVert^2}{\lVert \bm u_{\bm X} \rVert \kappa}. 247$$ (eq-peclet) 248 249For scalar advection-diffusion, the stabilization is a scalar 250 251$$ 252\tau = \frac{\xi(\mathrm{Pe})}{\lVert \bm u_{\bm X} \rVert}, 253$$ (eq-tau-advdiff) 254 255where $\xi(\mathrm{Pe}) = \coth \mathrm{Pe} - 1/\mathrm{Pe}$ approaches 1 at large local Péclet number. 256Note 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. 257For advection-diffusion, $\bm F(q) = \bm u q$, and thus the perturbed test function is 258 259$$ 260\mathcal P(v) = \tau \bm u \cdot \nabla v = \tau \bm u_{\bm X} \nabla_{\bm X} v. 261$$ (eq-test-perturbation-advdiff) 262 263See {cite}`hughesetal2010` equations 15-17 and 34-36 for further discussion of this formulation. 264 265For 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 2661. continuity stabilization $\tau_c$ 2672. momentum stabilization $\tau_m$ 2683. energy stabilization $\tau_E$ 269 270The Navier-Stokes code in this example uses the following formulation for $\tau_c$, $\tau_m$, $\tau_E$: 271 272$$ 273\begin{aligned} 274 275\tau_c &= \frac{C_c \mathcal{F}}{8\rho \trace(\bm g)} \\ 276\tau_m &= \frac{C_m}{\mathcal{F}} \\ 277\tau_E &= \frac{C_E}{\mathcal{F} c_v} \\ 278\end{aligned} 279$$ 280 281$$ 282\mathcal{F} = \sqrt{ \rho^2 \left [ \left(\frac{2C_t}{\Delta t}\right)^2 283+ \bm u \cdot (\bm u \cdot \bm g) 284+ C_v \mu^2 \Vert \bm g \Vert_F ^2\right]} 285$$ 286 287where $\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. 288This formulation is currently not available in the Euler code. 289 290In the Euler code, we follow {cite}`hughesetal2010` in defining a $3\times 3$ diagonal stabilization according to spatial criterion 2 (equation 27) as follows. 291 292$$ 293\tau_{ii} = c_{\tau} \frac{2 \xi(\mathrm{Pe})}{(\lambda_{\max \text{abs}})_i \lVert \nabla_{x_i} \bm X \rVert} 294$$ (eq-tau-conservative) 295 296where $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$. 297The 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. 298The complete set of eigenvalues of the Euler flux Jacobian in direction $i$ are (e.g., {cite}`toro2009`) 299 300$$ 301\Lambda_i = [u_i - a, u_i, u_i, u_i, u_i+a], 302$$ (eq-eigval-advdiff) 303 304where $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. 305Note 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. 306The fastest wave speed in direction $i$ is thus 307 308$$ 309\lambda_{\max \text{abs}} \Bigl( \frac{\partial \bm F_{\text{adv}}}{\partial \bm q} \cdot \hat{\bm n}_i \Bigr) = |u_i| + a 310$$ (eq-wavespeed) 311 312Note 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. 313 314::: 315 316Currently, this demo provides three types of problems/physical models that can be selected at run time via the option `-problem`. 317{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. 318 319(problem-advection)= 320 321## Advection 322 323A simplified version of system {eq}`eq-ns`, only accounting for the transport of total energy, is given by 324 325$$ 326\frac{\partial E}{\partial t} + \nabla \cdot (\bm{u} E ) = 0 \, , 327$$ (eq-advection) 328 329with $\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. 330 331- **Rotation** 332 333 In this case, a uniform circular velocity field transports the blob of total energy. 334 We have solved {eq}`eq-advection` applying zero energy density $E$, and no-flux for $\bm{u}$ on the boundaries. 335 336- **Translation** 337 338 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. 339 340 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 341 342 $$ 343 \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 \, , 344 $$ 345 346 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. 347 The weak form boundary integral in {eq}`eq-weak-vector-ns` for outflow boundary conditions is defined as 348 349 $$ 350 \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 \, , 351 $$ 352 353(problem-euler-vortex)= 354 355## Isentropic Vortex 356 357Three-dimensional Euler equations, which are simplified and nondimensionalized version of system {eq}`eq-ns` and account only for the convective fluxes, are given by 358 359$$ 360\begin{aligned} 361\frac{\partial \rho}{\partial t} + \nabla \cdot \bm{U} &= 0 \\ 362\frac{\partial \bm{U}}{\partial t} + \nabla \cdot \left( \frac{\bm{U} \otimes \bm{U}}{\rho} + P \bm{I}_3 \right) &= 0 \\ 363\frac{\partial E}{\partial t} + \nabla \cdot \left( \frac{(E + P)\bm{U}}{\rho} \right) &= 0 \, , \\ 364\end{aligned} 365$$ (eq-euler) 366 367Following 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 368 369$$ 370\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} 371$$ 372 373where $(\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). 374There is no perturbation in the entropy $S=P/\rho^\gamma$ ($\delta S=0)$. 375 376(problem-shock-tube)= 377 378## Shock Tube 379 380This 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. 381 382SU 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 383 384$$ 385\int_{\Omega} \nu_{SHOCK} \nabla \bm v \!:\! \nabla \bm q dV 386$$ 387 388The 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 389 390$$ 391\nu_{SHOCK} = \tau_{SHOCK} u_{cha}^2 392$$ 393 394where, 395 396$$ 397\tau_{SHOCK} = \frac{h_{SHOCK}}{2u_{cha}} \left( \frac{ \,|\, \nabla \rho \,|\, h_{SHOCK}}{\rho_{ref}} \right)^{\beta} 398$$ 399 400$\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 401 402$$ 403h_{SHOCK} = 2 \left( C_{YZB} \,|\, \bm p \,|\, \right)^{-1} 404$$ 405 406where 407 408$$ 409p_k = \hat{j}_i \frac{\partial \xi_i}{x_k} 410$$ 411 412The 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. 413 414(problem-density-current)= 415 416## Newtonian Wave 417This 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. 418 419The problem has a perturbed initial condition and lets it evolve in time. The initial condition contains a Gaussian perturbation in the pressure field: 420 421$$ 422\begin{aligned} 423\rho &= \rho_\infty\left(1+A\exp\left(\frac{-(\bar{x}^2 + \bar{y}^2)}{2\sigma^2}\right)\right) \\ 424\bm{U} &= \bm U_\infty \\ 425E &= \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}, 426\end{aligned} 427$$ 428 429where $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)$. 430The simulation produces a strong acoustic wave and leaves behind a cold thermal bubble that advects at the fluid velocity. 431 432The 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. 433This 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. 434 435## Density Current 436 437For 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. 438Its 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 439 440$$ 441\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} 442$$ 443 444where $P_0$ is the atmospheric pressure. 445For this problem, we have used no-slip and non-penetration boundary conditions for $\bm{u}$, and no-flux for mass and energy densities. 446 447## Channel 448 449A compressible channel flow. Analytical solution given in 450{cite}`whitingStabilizedFEM1999`: 451 452$$ u_1 = u_{\max} \left [ 1 - \left ( \frac{x_2}{H}\right)^2 \right] \quad \quad u_2 = u_3 = 0$$ 453$$T = T_w \left [ 1 + \frac{Pr \hat{E}c}{3} \left \{1 - \left(\frac{x_2}{H}\right)^4 \right \} \right]$$ 454$$p = p_0 - \frac{2\rho_0 u_{\max}^2 x_1}{Re_H H}$$ 455 456where $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. 457 458Boundary conditions are periodic in the streamwise direction, and no-slip and non-penetration boundary conditions at the walls. 459The flow is driven by a body force determined analytically from the fluid properties and setup parameters $H$ and $u_{\max}$. 460 461## Flat Plate Boundary Layer 462 463### Laminar Boundary Layer - Blasius 464 465Simulation of a laminar boundary layer flow, with the inflow being prescribed 466by a [Blasius similarity 467solution](https://en.wikipedia.org/wiki/Blasius_boundary_layer). At the inflow, 468the velocity is prescribed by the Blasius soution profile, density is set 469constant, and temperature is allowed to float. Using `weakT: true`, density is 470allowed to float and temperature is set constant. At the outlet, a user-set 471pressure is used for pressure in the inviscid flux terms (all other inviscid 472flux terms use interior solution values). The wall is a no-slip, 473no-penetration, no-heat flux condition. The top of the domain is treated as an 474outflow and is tilted at a downward angle to ensure that flow is always exiting 475it. 476 477### Turbulent Boundary Layer 478 479Simulating a turbulent boundary layer without modeling the turbulence requires 480resolving the turbulent flow structures. These structures may be introduced 481into the simulations either by allowing a laminar boundary layer naturally 482transition to turbulence, or imposing turbulent structures at the inflow. The 483latter approach has been taken here, specifically using a *synthetic turbulence 484generation* (STG) method. 485 486#### Synthetic Turbulence Generation (STG) Boundary Condition 487 488We use the STG method described in 489{cite}`shurSTG2014`. Below follows a re-description of the formulation to match 490the present notation, and then a description of the implementation and usage. 491 492##### Equation Formulation 493 494$$ 495\bm{u}(\bm{x}, t) = \bm{\overline{u}}(\bm{x}) + \bm{C}(\bm{x}) \cdot \bm{v}' 496$$ 497 498$$ 499\begin{aligned} 500\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 ) \\ 501\bm{\hat{x}}^n &= \left[(x - U_0 t)\max(2\kappa_{\min}/\kappa^n, 0.1) , y, z \right]^T 502\end{aligned} 503$$ 504 505Here, we define the number of wavemodes $N$, set of random numbers $ \{\bm{\sigma}^n, 506\bm{d}^n, \phi^n\}_{n=1}^N$, the Cholesky decomposition of the Reynolds stress 507tensor $\bm{C}$ (such that $\bm{R} = \bm{CC}^T$ ), bulk velocity $U_0$, 508wavemode amplitude $q^n$, wavemode frequency $\kappa^n$, and $\kappa_{\min} = 5090.5 \min_{\bm{x}} (\kappa_e)$. 510 511$$ 512\kappa_e = \frac{2\pi}{\min(2d_w, 3.0 l_t)} 513$$ 514 515where $l_t$ is the turbulence length scale, and $d_w$ is the distance to the 516nearest wall. 517 518 519The set of wavemode frequencies is defined by a geometric distribution: 520 521$$ 522\kappa^n = \kappa_{\min} (1 + \alpha)^{n-1} \ , \quad \forall n=1, 2, ... , N 523$$ 524 525The wavemode amplitudes $q^n$ are defined by a model energy spectrum $E(\kappa)$: 526 527$$ 528q^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} 529$$ 530 531$$ E(\kappa) = \frac{(\kappa/\kappa_e)^4}{[1 + 2.4(\kappa/\kappa_e)^2]^{17/6}} f_\eta f_{\mathrm{cut}} $$ 532 533$$ 534f_\eta = \exp \left[-(12\kappa /\kappa_\eta)^2 \right], \quad 535f_\mathrm{cut} = \exp \left( - \left [ \frac{4\max(\kappa-0.9\kappa_\mathrm{cut}, 0)}{\kappa_\mathrm{cut}} \right]^3 \right) 536$$ 537 538$\kappa_\eta$ represents turbulent dissipation frequency, and is given as $2\pi 539(\nu^3/\varepsilon)^{-1/4}$ with $\nu$ the kinematic viscosity and 540$\varepsilon$ the turbulent dissipation. $\kappa_\mathrm{cut}$ approximates the 541effective cutoff frequency of the mesh (viewing the mesh as a filter on 542solution over $\Omega$) and is given by: 543 544$$ 545\kappa_\mathrm{cut} = \frac{2\pi}{ 2\min\{ [\max(h_y, h_z, 0.3h_{\max}) + 0.1 d_w], h_{\max} \} } 546$$ 547 548The enforcement of the boundary condition is identical to the blasius inflow; 549it weakly enforces velocity, with the option of weakly enforcing either density 550or temperature using the the `-weakT` flag. 551 552##### Initialization Data Flow 553 554Data flow for initializing function (which creates the context data struct) is 555given below: 556```{mermaid} 557flowchart LR 558 subgraph STGInflow.dat 559 y 560 lt[l_t] 561 eps 562 Rij[R_ij] 563 ubar 564 end 565 566 subgraph STGRand.dat 567 rand[RN Set]; 568 end 569 570 subgraph User Input 571 u0[U0]; 572 end 573 574 subgraph init[Create Context Function] 575 ke[k_e] 576 N; 577 end 578 lt --Calc-->ke --Calc-->kn 579 y --Calc-->ke 580 581 subgraph context[Context Data] 582 yC[y] 583 randC[RN Set] 584 Cij[C_ij] 585 u0 --Copy--> u0C[U0] 586 kn[k^n]; 587 ubarC[ubar] 588 ltC[l_t] 589 epsC[eps] 590 end 591 ubar --Copy--> ubarC; 592 y --Copy--> yC; 593 lt --Copy--> ltC; 594 eps --Copy--> epsC; 595 596 rand --Copy--> randC; 597 rand --> N --Calc--> kn; 598 Rij --Calc--> Cij[C_ij] 599``` 600 601This is done once at runtime. The spatially-varying terms are then evaluated at 602each quadrature point on-the-fly, either by interpolation (for $l_t$, 603$\varepsilon$, $C_{ij}$, and $\overline{\bm u}$) or by calculation (for $q^n$). 604 605The `STGInflow.dat` file is a table of values at given distances from the wall. 606These values are then interpolated to a physical location (node or quadrature 607point). It has the following format: 608``` 609[Total number of locations] 14 610[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] 611``` 612where each `[ ]` item is a number in scientific notation (ie. `3.1415E0`), and `sclr_1` and 613`sclr_2` are reserved for turbulence modeling variables. They are not used in 614this example. 615 616The `STGRand.dat` file is the table of the random number set, $\{\bm{\sigma}^n, 617\bm{d}^n, \phi^n\}_{n=1}^N$. It has the format: 618``` 619[Number of wavemodes] 7 620[d_1] [d_2] [d_3] [phi] [sigma_1] [sigma_2] [sigma_3] 621``` 622 623The following table is presented to help clarify the dimensionality of the 624numerous terms in the STG formulation. 625 626| Math | Label | $f(\bm{x})$? | $f(n)$? | 627|-----------------|--------|--------------|---------| 628| $ \{\bm{\sigma}^n, \bm{d}^n, \phi^n\}_{n=1}^N$ | RN Set | No | Yes | 629| $\bm{\overline{u}}$ | ubar | Yes | No | 630| $U_0$ | U0 | No | No | 631| $l_t$ | l_t | Yes | No | 632| $\varepsilon$ | eps | Yes | No | 633| $\bm{R}$ | R_ij | Yes | No | 634| $\bm{C}$ | C_ij | Yes | No | 635| $q^n$ | q^n | Yes | Yes | 636| $\{\kappa^n\}_{n=1}^N$ | k^n | No | Yes | 637| $h_i$ | h_i | Yes | No | 638| $d_w$ | d_w | Yes | No | 639 640### Meshing 641 642The flat plate boundary layer example has custom meshing features to better 643resolve the flow. One of those is tilting the top of the domain, allowing for 644it to be a outflow boundary condition. The angle of this tilt is controled by 645`-platemesh_top_angle` 646 647The primary meshing feature is the ability to grade the mesh, providing better 648resolution near the wall. There are two methods to do this; algorithmically, or 649specifying the node locations via a file. Algorithmically, a base node 650distribution is defined at the inlet (assumed to be $\min(x)$) and then 651linearly stretched/squeezed to match the slanted top boundary condition. Nodes 652are placed such that `-platemesh_Ndelta` elements are within 653`-platemesh_refine_height` of the wall. They are placed such that the element 654height matches a geometric growth ratio defined by `-platemesh_growth`. The 655remaining elements are then distributed from `-platemesh_refine_height` to the 656top of the domain linearly in logarithmic space. 657 658Alternatively, a file may be specified containing the locations of each node. 659The file should be newline delimited, with the first line specifying the number 660of points and the rest being the locations of the nodes. The node locations 661used exactly at the inlet (assumed to be $\min(x)$) and linearly 662stretched/squeezed to match the slanted top boundary condition. The file is 663specified via `-platemesh_y_node_locs_path`. If this flag is given an empty 664string, then the algorithmic approach will be performed. 665