xref: /libCEED/examples/fluids/qfunctions/eulervortex.h (revision 2b730f8b5a9c809740a0b3b302db43a719c636b1)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Euler traveling vortex initial condition and operator for Navier-Stokes
10 /// example using PETSc
11 
12 // Model from:
13 //   On the Order of Accuracy and Numerical Performance of Two Classes of
14 //   Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011).
15 
16 #ifndef eulervortex_h
17 #define eulervortex_h
18 
19 #include <ceed.h>
20 #include <math.h>
21 
22 #include "utils.h"
23 
24 typedef struct EulerContext_ *EulerContext;
25 struct EulerContext_ {
26   CeedScalar center[3];
27   CeedScalar curr_time;
28   CeedScalar vortex_strength;
29   CeedScalar c_tau;
30   CeedScalar mean_velocity[3];
31   bool       implicit;
32   int        euler_test;
33   int        stabilization;  // See StabilizationType: 0=none, 1=SU, 2=SUPG
34 };
35 
36 // *****************************************************************************
37 // This function sets the initial conditions
38 //
39 //   Temperature:
40 //     T   = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2)
41 //   Density:
42 //     rho = (T/S_vortex)^(1 / (gamma - 1))
43 //   Pressure:
44 //     P   = rho * T
45 //   Velocity:
46 //     ui  = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi)
47 //     r   = sqrt( (x - xc)**2 + (y - yc)**2 )
48 //   Velocity/Momentum Density:
49 //     Ui  = rho ui
50 //   Total Energy:
51 //     E   = P / (gamma - 1) + rho (u u)/2
52 //
53 // Constants:
54 //   cv              ,  Specific heat, constant volume
55 //   cp              ,  Specific heat, constant pressure
56 //   vortex_strength ,  Strength of vortex
57 //   center          ,  Location of bubble center
58 //   gamma  = cp / cv,  Specific heat ratio
59 //
60 // *****************************************************************************
61 
62 // *****************************************************************************
63 // This helper function provides support for the exact, time-dependent solution
64 //   (currently not implemented) and IC formulation for Euler traveling vortex
65 // *****************************************************************************
66 CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
67   // Context
68   const EulerContext context         = (EulerContext)ctx;
69   const CeedScalar   vortex_strength = context->vortex_strength;
70   const CeedScalar  *center          = context->center;  // Center of the domain
71   const CeedScalar  *mean_velocity   = context->mean_velocity;
72 
73   // Setup
74   const CeedScalar gamma = 1.4;
75   const CeedScalar cv    = 2.5;
76   const CeedScalar R     = 1.;
77   const CeedScalar x = X[0], y = X[1];  // Coordinates
78   // Vortex center
79   const CeedScalar xc = center[0] + mean_velocity[0] * time;
80   const CeedScalar yc = center[1] + mean_velocity[1] * time;
81 
82   const CeedScalar x0       = x - xc;
83   const CeedScalar y0       = y - yc;
84   const CeedScalar r        = sqrt(x0 * x0 + y0 * y0);
85   const CeedScalar C        = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI);
86   const CeedScalar delta_T  = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI);
87   const CeedScalar S_vortex = 1;  // no perturbation in the entropy P / rho^gamma
88   const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI);
89   CeedScalar       rho, P, T, E, u[3] = {0.};
90 
91   // Initial Conditions
92   switch (context->euler_test) {
93     case 0:  // Traveling vortex
94       T = 1 + delta_T;
95       // P = rho * T
96       // P = S * rho^gamma
97       // Solve for rho, then substitute for P
98       rho  = pow(T / S_vortex, 1 / (gamma - 1.));
99       P    = rho * T;
100       u[0] = mean_velocity[0] - C * y0;
101       u[1] = mean_velocity[1] + C * x0;
102 
103       // Assign exact solution
104       q[0] = rho;
105       q[1] = rho * u[0];
106       q[2] = rho * u[1];
107       q[3] = rho * u[2];
108       q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.;
109       break;
110     case 1:  // Constant zero velocity, density constant, total energy constant
111       rho = 1.;
112       E   = 2.;
113 
114       // Assign exact solution
115       q[0] = rho;
116       q[1] = rho * u[0];
117       q[2] = rho * u[1];
118       q[3] = rho * u[2];
119       q[4] = E;
120       break;
121     case 2:  // Constant nonzero velocity, density constant, total energy constant
122       rho  = 1.;
123       E    = 2.;
124       u[0] = mean_velocity[0];
125       u[1] = mean_velocity[1];
126 
127       // Assign exact solution
128       q[0] = rho;
129       q[1] = rho * u[0];
130       q[2] = rho * u[1];
131       q[3] = rho * u[2];
132       q[4] = E;
133       break;
134     case 3:  // Velocity zero, pressure constant
135       // (so density and internal energy will be non-constant),
136       // but the velocity should stay zero and the bubble won't diffuse
137       // (for Euler, where there is no thermal conductivity)
138       P   = 1.;
139       T   = 1. - S_bubble * exp(1. - r * r);
140       rho = P / (R * T);
141 
142       // Assign exact solution
143       q[0] = rho;
144       q[1] = rho * u[0];
145       q[2] = rho * u[1];
146       q[3] = rho * u[2];
147       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
148       break;
149     case 4:  // Constant nonzero velocity, pressure constant
150       // (so density and internal energy will be non-constant),
151       // it should be transported across the domain, but velocity stays constant
152       P    = 1.;
153       T    = 1. - S_bubble * exp(1. - r * r);
154       rho  = P / (R * T);
155       u[0] = mean_velocity[0];
156       u[1] = mean_velocity[1];
157 
158       // Assign exact solution
159       q[0] = rho;
160       q[1] = rho * u[0];
161       q[2] = rho * u[1];
162       q[3] = rho * u[2];
163       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
164       break;
165     case 5:  // non-smooth thermal bubble - cylinder
166       P    = 1.;
167       T    = 1. - (r < 1. ? S_bubble : 0.);
168       rho  = P / (R * T);
169       u[0] = mean_velocity[0];
170       u[1] = mean_velocity[1];
171 
172       // Assign exact solution
173       q[0] = rho;
174       q[1] = rho * u[0];
175       q[2] = rho * u[1];
176       q[3] = rho * u[2];
177       q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
178       break;
179   }
180   // Return
181   return 0;
182 }
183 
184 // *****************************************************************************
185 // Helper function for computing flux Jacobian
186 // *****************************************************************************
187 CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
188                                                         const CeedScalar gamma) {
189   CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2];  // Velocity square
190   for (CeedInt i = 0; i < 3; i++) {                           // Jacobian matrices for 3 directions
191     for (CeedInt j = 0; j < 3; j++) {                         // Rows of each Jacobian matrix
192       dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j];
193       for (CeedInt k = 0; k < 3; k++) {  // Columns of each Jacobian matrix
194         dF[i][0][k + 1]     = ((i == k) ? 1. : 0.);
195         dF[i][j + 1][k + 1] = ((j == k) ? u[i] : 0.) + ((i == k) ? u[j] : 0.) - ((i == j) ? u[k] : 0.) * (gamma - 1.);
196         dF[i][4][k + 1]     = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k];
197       }
198       dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.);
199     }
200     dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho);
201     dF[i][4][4] = u[i] * gamma;
202   }
203 }
204 
205 // *****************************************************************************
206 // Helper function for computing Tau elements (stabilization constant)
207 //   Model from:
208 //     Stabilized Methods for Compressible Flows, Hughes et al 2010
209 //
210 //   Spatial criterion #2 - Tau is a 3x3 diagonal matrix
211 //   Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
212 //
213 // Where
214 //   c_tau     = stabilization constant (0.5 is reported as "optimal")
215 //   h[i]      = 2 length(dxdX[i])
216 //   Pe        = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
217 //   Xi(Pe)    = coth Pe - 1. / Pe (1. at large local Peclet number )
218 //   rho(A[i]) = spectral radius of the convective flux Jacobian i,
219 //               wave speed in direction i
220 // *****************************************************************************
221 CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], const CeedScalar dXdx[3][3], const CeedScalar u[3], const CeedScalar sound_speed,
222                                        const CeedScalar c_tau) {
223   for (CeedInt i = 0; i < 3; i++) {
224     // length of element in direction i
225     CeedScalar h = 2 / sqrt(dXdx[0][i] * dXdx[0][i] + dXdx[1][i] * dXdx[1][i] + dXdx[2][i] * dXdx[2][i]);
226     // fastest wave in direction i
227     CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
228     Tau_x[i]                = c_tau * h / fastest_wave;
229   }
230 }
231 
232 // *****************************************************************************
233 // This QFunction sets the initial conditions for Euler traveling vortex
234 // *****************************************************************************
235 CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
236   // Inputs
237   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
238 
239   // Outputs
240   CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
241   const EulerContext context  = (EulerContext)ctx;
242 
243   CeedPragmaSIMD
244       // Quadrature Point Loop
245       for (CeedInt i = 0; i < Q; i++) {
246     const CeedScalar x[]  = {X[0][i], X[1][i], X[2][i]};
247     CeedScalar       q[5] = {0.};
248 
249     Exact_Euler(3, context->curr_time, x, 5, q, ctx);
250 
251     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
252   }  // End of Quadrature Point Loop
253 
254   // Return
255   return 0;
256 }
257 
258 // *****************************************************************************
259 // This QFunction implements the following formulation of Euler equations
260 //   with explicit time stepping method
261 //
262 // This is 3D Euler for compressible gas dynamics in conservation
263 //   form with state variables of density, momentum density, and total
264 //   energy density.
265 //
266 // State Variables: q = ( rho, U1, U2, U3, E )
267 //   rho - Mass Density
268 //   Ui  - Momentum Density,      Ui = rho ui
269 //   E   - Total Energy Density,  E  = P / (gamma - 1) + rho (u u)/2
270 //
271 // Euler Equations:
272 //   drho/dt + div( U )                   = 0
273 //   dU/dt   + div( rho (u x u) + P I3 )  = 0
274 //   dE/dt   + div( (E + P) u )           = 0
275 //
276 // Equation of State:
277 //   P = (gamma - 1) (E - rho (u u) / 2)
278 //
279 // Constants:
280 //   cv              ,  Specific heat, constant volume
281 //   cp              ,  Specific heat, constant pressure
282 //   g               ,  Gravity
283 //   gamma  = cp / cv,  Specific heat ratio
284 // *****************************************************************************
285 CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
286   // *INDENT-OFF*
287   // Inputs
288   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
289         (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
290   // Outputs
291   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
292 
293   EulerContext     context = (EulerContext)ctx;
294   const CeedScalar c_tau   = context->c_tau;
295   const CeedScalar gamma   = 1.4;
296 
297   CeedPragmaSIMD
298       // Quadrature Point Loop
299       for (CeedInt i = 0; i < Q; i++) {
300     // *INDENT-OFF*
301     // Setup
302     // -- Interp in
303     const CeedScalar rho      = q[0][i];
304     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
305     const CeedScalar E        = q[4][i];
306     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
307     const CeedScalar dU[3][3] = {
308         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
309         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
310         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
311     };
312     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
313     // -- Interp-to-Interp q_data
314     const CeedScalar wdetJ = q_data[0][i];
315     // -- Interp-to-Grad q_data
316     // ---- Inverse of change of coordinate matrix: X_i,j
317     // *INDENT-OFF*
318     const CeedScalar dXdx[3][3] = {
319         {q_data[1][i], q_data[2][i], q_data[3][i]},
320         {q_data[4][i], q_data[5][i], q_data[6][i]},
321         {q_data[7][i], q_data[8][i], q_data[9][i]}
322     };
323     // *INDENT-ON*
324     // dU/dx
325     CeedScalar drhodx[3]       = {0.};
326     CeedScalar dEdx[3]         = {0.};
327     CeedScalar dUdx[3][3]      = {{0.}};
328     CeedScalar dXdxdXdxT[3][3] = {{0.}};
329     for (CeedInt j = 0; j < 3; j++) {
330       for (CeedInt k = 0; k < 3; k++) {
331         drhodx[j] += drho[k] * dXdx[k][j];
332         dEdx[j] += dE[k] * dXdx[k][j];
333         for (CeedInt l = 0; l < 3; l++) {
334           dUdx[j][k] += dU[j][l] * dXdx[l][k];
335           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
336         }
337       }
338     }
339     // Pressure
340     const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic,
341                      P = E_internal * (gamma - 1.);  // P = pressure
342 
343     // The Physics
344     // Zero v and dv so all future terms can safely sum into it
345     for (CeedInt j = 0; j < 5; j++) {
346       v[j][i] = 0.;
347       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
348     }
349 
350     // -- Density
351     // ---- u rho
352     for (CeedInt j = 0; j < 3; j++) dv[j][0][i] += wdetJ * (rho * u[0] * dXdx[j][0] + rho * u[1] * dXdx[j][1] + rho * u[2] * dXdx[j][2]);
353     // -- Momentum
354     // ---- rho (u x u) + P I3
355     for (CeedInt j = 0; j < 3; j++) {
356       for (CeedInt k = 0; k < 3; k++) {
357         dv[k][j + 1][i] += wdetJ * ((rho * u[j] * u[0] + (j == 0 ? P : 0.)) * dXdx[k][0] + (rho * u[j] * u[1] + (j == 1 ? P : 0.)) * dXdx[k][1] +
358                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
359       }
360     }
361     // -- Total Energy Density
362     // ---- (E + P) u
363     for (CeedInt j = 0; j < 3; j++) dv[j][4][i] += wdetJ * (E + P) * (u[0] * dXdx[j][0] + u[1] * dXdx[j][1] + u[2] * dXdx[j][2]);
364 
365     // --Stabilization terms
366     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
367     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
368     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
369 
370     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
371     CeedScalar dqdx[5][3];
372     for (CeedInt j = 0; j < 3; j++) {
373       dqdx[0][j] = drhodx[j];
374       dqdx[4][j] = dEdx[j];
375       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
376     }
377 
378     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
379     CeedScalar strong_conv[5] = {0.};
380     for (CeedInt j = 0; j < 3; j++) {
381       for (CeedInt k = 0; k < 5; k++) {
382         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
383       }
384     }
385 
386     // Stabilization
387     // -- Tau elements
388     const CeedScalar sound_speed = sqrt(gamma * P / rho);
389     CeedScalar       Tau_x[3]    = {0.};
390     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
391 
392     // -- Stabilization method: none or SU
393     CeedScalar stab[5][3] = {{0.}};
394     switch (context->stabilization) {
395       case 0:  // Galerkin
396         break;
397       case 1:  // SU
398         for (CeedInt j = 0; j < 3; j++) {
399           for (CeedInt k = 0; k < 5; k++) {
400             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
401           }
402         }
403 
404         for (CeedInt j = 0; j < 5; j++) {
405           for (CeedInt k = 0; k < 3; k++) dv[k][j][i] -= wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]);
406         }
407         break;
408       case 2:  // SUPG is not implemented for explicit scheme
409         break;
410     }
411 
412   }  // End Quadrature Point Loop
413 
414   // Return
415   return 0;
416 }
417 
418 // *****************************************************************************
419 // This QFunction implements the Euler equations with (mentioned above)
420 //   with implicit time stepping method
421 //
422 // *****************************************************************************
423 CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
424   // *INDENT-OFF*
425   // Inputs
426   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
427         (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
428   // Outputs
429   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
430 
431   EulerContext     context = (EulerContext)ctx;
432   const CeedScalar c_tau   = context->c_tau;
433   const CeedScalar gamma   = 1.4;
434 
435   CeedPragmaSIMD
436       // Quadrature Point Loop
437       for (CeedInt i = 0; i < Q; i++) {
438     // *INDENT-OFF*
439     // Setup
440     // -- Interp in
441     const CeedScalar rho      = q[0][i];
442     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
443     const CeedScalar E        = q[4][i];
444     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
445     const CeedScalar dU[3][3] = {
446         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
447         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
448         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
449     };
450     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
451     // -- Interp-to-Interp q_data
452     const CeedScalar wdetJ = q_data[0][i];
453     // -- Interp-to-Grad q_data
454     // ---- Inverse of change of coordinate matrix: X_i,j
455     // *INDENT-OFF*
456     const CeedScalar dXdx[3][3] = {
457         {q_data[1][i], q_data[2][i], q_data[3][i]},
458         {q_data[4][i], q_data[5][i], q_data[6][i]},
459         {q_data[7][i], q_data[8][i], q_data[9][i]}
460     };
461     // *INDENT-ON*
462     // dU/dx
463     CeedScalar drhodx[3]       = {0.};
464     CeedScalar dEdx[3]         = {0.};
465     CeedScalar dUdx[3][3]      = {{0.}};
466     CeedScalar dXdxdXdxT[3][3] = {{0.}};
467     for (CeedInt j = 0; j < 3; j++) {
468       for (CeedInt k = 0; k < 3; k++) {
469         drhodx[j] += drho[k] * dXdx[k][j];
470         dEdx[j] += dE[k] * dXdx[k][j];
471         for (CeedInt l = 0; l < 3; l++) {
472           dUdx[j][k] += dU[j][l] * dXdx[l][k];
473           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
474         }
475       }
476     }
477     const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic,
478                      P = E_internal * (gamma - 1.);  // P = pressure
479 
480     // The Physics
481     // Zero v and dv so all future terms can safely sum into it
482     for (CeedInt j = 0; j < 5; j++) {
483       v[j][i] = 0.;
484       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
485     }
486     //-----mass matrix
487     for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i];
488 
489     // -- Density
490     // ---- u rho
491     for (CeedInt j = 0; j < 3; j++) dv[j][0][i] -= wdetJ * (rho * u[0] * dXdx[j][0] + rho * u[1] * dXdx[j][1] + rho * u[2] * dXdx[j][2]);
492     // -- Momentum
493     // ---- rho (u x u) + P I3
494     for (CeedInt j = 0; j < 3; j++) {
495       for (CeedInt k = 0; k < 3; k++) {
496         dv[k][j + 1][i] -= wdetJ * ((rho * u[j] * u[0] + (j == 0 ? P : 0.)) * dXdx[k][0] + (rho * u[j] * u[1] + (j == 1 ? P : 0.)) * dXdx[k][1] +
497                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
498       }
499     }
500     // -- Total Energy Density
501     // ---- (E + P) u
502     for (CeedInt j = 0; j < 3; j++) dv[j][4][i] -= wdetJ * (E + P) * (u[0] * dXdx[j][0] + u[1] * dXdx[j][1] + u[2] * dXdx[j][2]);
503 
504     // -- Stabilization terms
505     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
506     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
507     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
508 
509     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
510     CeedScalar dqdx[5][3];
511     for (CeedInt j = 0; j < 3; j++) {
512       dqdx[0][j] = drhodx[j];
513       dqdx[4][j] = dEdx[j];
514       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
515     }
516 
517     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
518     CeedScalar strong_conv[5] = {0.};
519     for (CeedInt j = 0; j < 3; j++) {
520       for (CeedInt k = 0; k < 5; k++) {
521         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
522       }
523     }
524 
525     // ---- Strong residual
526     CeedScalar strong_res[5];
527     for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j];
528 
529     // Stabilization
530     // -- Tau elements
531     const CeedScalar sound_speed = sqrt(gamma * P / rho);
532     CeedScalar       Tau_x[3]    = {0.};
533     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
534 
535     // -- Stabilization method: none, SU, or SUPG
536     CeedScalar stab[5][3] = {{0.}};
537     switch (context->stabilization) {
538       case 0:  // Galerkin
539         break;
540       case 1:  // SU
541         for (CeedInt j = 0; j < 3; j++) {
542           for (CeedInt k = 0; k < 5; k++) {
543             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
544           }
545         }
546 
547         for (CeedInt j = 0; j < 5; j++) {
548           for (CeedInt k = 0; k < 3; k++) dv[k][j][i] += wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]);
549         }
550         break;
551       case 2:  // SUPG
552         for (CeedInt j = 0; j < 3; j++) {
553           for (CeedInt k = 0; k < 5; k++) {
554             for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l];
555           }
556         }
557 
558         for (CeedInt j = 0; j < 5; j++) {
559           for (CeedInt k = 0; k < 3; k++) dv[k][j][i] += wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]);
560         }
561         break;
562     }
563   }  // End Quadrature Point Loop
564 
565   // Return
566   return 0;
567 }
568 // *****************************************************************************
569 // This QFunction sets the inflow boundary conditions for
570 //   the traveling vortex problem.
571 //
572 //  Prescribed T_inlet and P_inlet are converted to conservative variables
573 //      and applied weakly.
574 //
575 // *****************************************************************************
576 CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
577   // *INDENT-OFF*
578   // Inputs
579   const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
580   // Outputs
581   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
582   // *INDENT-ON*
583   EulerContext     context       = (EulerContext)ctx;
584   const int        euler_test    = context->euler_test;
585   const bool       implicit      = context->implicit;
586   CeedScalar      *mean_velocity = context->mean_velocity;
587   const CeedScalar cv            = 2.5;
588   const CeedScalar R             = 1.;
589   CeedScalar       T_inlet;
590   CeedScalar       P_inlet;
591 
592   // For test cases 1 and 3 the background velocity is zero
593   if (euler_test == 1 || euler_test == 3) {
594     for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.;
595   }
596 
597   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
598   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
599   else T_inlet = P_inlet = 1.;
600 
601   CeedPragmaSIMD
602       // Quadrature Point Loop
603       for (CeedInt i = 0; i < Q; i++) {
604     // Setup
605     // -- Interp-to-Interp q_data
606     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
607     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
608     // We can effect this by swapping the sign on this weight
609     const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i];
610     // ---- Normal vect
611     const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]};
612 
613     // face_normal = Normal vector of the face
614     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
615     // The Physics
616     // Zero v so all future terms can safely sum into it
617     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.;
618 
619     // Implementing in/outflow BCs
620     if (face_normal > 0) {
621     } else {  // inflow
622       const CeedScalar rho_inlet       = P_inlet / (R * T_inlet);
623       const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.;
624       // incoming total energy
625       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
626 
627       // The Physics
628       // -- Density
629       v[0][i] -= wdetJb * rho_inlet * face_normal;
630 
631       // -- Momentum
632       for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho_inlet * face_normal * mean_velocity[j] + norm[j] * P_inlet);
633 
634       // -- Total Energy Density
635       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
636     }
637 
638   }  // End Quadrature Point Loop
639   return 0;
640 }
641 
642 // *****************************************************************************
643 // This QFunction sets the outflow boundary conditions for
644 //   the Euler solver.
645 //
646 //  Outflow BCs:
647 //    The validity of the weak form of the governing equations is
648 //      extended to the outflow.
649 //
650 // *****************************************************************************
651 CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
652   // *INDENT-OFF*
653   // Inputs
654   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
655   // Outputs
656   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
657   // *INDENT-ON*
658   EulerContext context       = (EulerContext)ctx;
659   const bool   implicit      = context->implicit;
660   CeedScalar  *mean_velocity = context->mean_velocity;
661 
662   const CeedScalar gamma = 1.4;
663 
664   CeedPragmaSIMD
665       // Quadrature Point Loop
666       for (CeedInt i = 0; i < Q; i++) {
667     // Setup
668     // -- Interp in
669     const CeedScalar rho  = q[0][i];
670     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
671     const CeedScalar E    = q[4][i];
672 
673     // -- Interp-to-Interp q_data
674     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
675     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
676     // We can effect this by swapping the sign on this weight
677     const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i];
678     // ---- Normal vectors
679     const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]};
680 
681     // face_normal = Normal vector of the face
682     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
683     // The Physics
684     // Zero v so all future terms can safely sum into it
685     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0;
686 
687     // Implementing in/outflow BCs
688     if (face_normal > 0) {  // outflow
689       const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.;
690       const CeedScalar P         = (E - E_kinetic * rho) * (gamma - 1.);              // pressure
691       const CeedScalar u_normal  = norm[0] * u[0] + norm[1] * u[1] + norm[2] * u[2];  // Normal velocity
692       // The Physics
693       // -- Density
694       v[0][i] -= wdetJb * rho * u_normal;
695 
696       // -- Momentum
697       for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P);
698 
699       // -- Total Energy Density
700       v[4][i] -= wdetJb * u_normal * (E + P);
701     }
702   }  // End Quadrature Point Loop
703   return 0;
704 }
705 
706 // *****************************************************************************
707 
708 #endif  // eulervortex_h
709