xref: /honee/qfunctions/eulervortex.h (revision d83e1a91d1ed05ca0043f14430b9f3bc9bc33a34)
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   // Inputs
287   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],
288         (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
289   // Outputs
290   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
291 
292   EulerContext     context = (EulerContext)ctx;
293   const CeedScalar c_tau   = context->c_tau;
294   const CeedScalar gamma   = 1.4;
295 
296   CeedPragmaSIMD
297       // Quadrature Point Loop
298       for (CeedInt i = 0; i < Q; i++) {
299     // Setup
300     // -- Interp in
301     const CeedScalar rho      = q[0][i];
302     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
303     const CeedScalar E        = q[4][i];
304     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
305     const CeedScalar dU[3][3] = {
306         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
307         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
308         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
309     };
310     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
311     // -- Interp-to-Interp q_data
312     const CeedScalar wdetJ = q_data[0][i];
313     // -- Interp-to-Grad q_data
314     // ---- Inverse of change of coordinate matrix: X_i,j
315     const CeedScalar dXdx[3][3] = {
316         {q_data[1][i], q_data[2][i], q_data[3][i]},
317         {q_data[4][i], q_data[5][i], q_data[6][i]},
318         {q_data[7][i], q_data[8][i], q_data[9][i]}
319     };
320     // dU/dx
321     CeedScalar drhodx[3]       = {0.};
322     CeedScalar dEdx[3]         = {0.};
323     CeedScalar dUdx[3][3]      = {{0.}};
324     CeedScalar dXdxdXdxT[3][3] = {{0.}};
325     for (CeedInt j = 0; j < 3; j++) {
326       for (CeedInt k = 0; k < 3; k++) {
327         drhodx[j] += drho[k] * dXdx[k][j];
328         dEdx[j] += dE[k] * dXdx[k][j];
329         for (CeedInt l = 0; l < 3; l++) {
330           dUdx[j][k] += dU[j][l] * dXdx[l][k];
331           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
332         }
333       }
334     }
335     // Pressure
336     const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic,
337                      P = E_internal * (gamma - 1.);  // P = pressure
338 
339     // The Physics
340     // Zero v and dv so all future terms can safely sum into it
341     for (CeedInt j = 0; j < 5; j++) {
342       v[j][i] = 0.;
343       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
344     }
345 
346     // -- Density
347     // ---- u rho
348     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]);
349     // -- Momentum
350     // ---- rho (u x u) + P I3
351     for (CeedInt j = 0; j < 3; j++) {
352       for (CeedInt k = 0; k < 3; k++) {
353         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] +
354                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
355       }
356     }
357     // -- Total Energy Density
358     // ---- (E + P) u
359     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]);
360 
361     // --Stabilization terms
362     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
363     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
364     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
365 
366     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
367     CeedScalar dqdx[5][3];
368     for (CeedInt j = 0; j < 3; j++) {
369       dqdx[0][j] = drhodx[j];
370       dqdx[4][j] = dEdx[j];
371       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
372     }
373 
374     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
375     CeedScalar strong_conv[5] = {0.};
376     for (CeedInt j = 0; j < 3; j++) {
377       for (CeedInt k = 0; k < 5; k++) {
378         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
379       }
380     }
381 
382     // Stabilization
383     // -- Tau elements
384     const CeedScalar sound_speed = sqrt(gamma * P / rho);
385     CeedScalar       Tau_x[3]    = {0.};
386     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
387 
388     // -- Stabilization method: none or SU
389     CeedScalar stab[5][3] = {{0.}};
390     switch (context->stabilization) {
391       case 0:  // Galerkin
392         break;
393       case 1:  // SU
394         for (CeedInt j = 0; j < 3; j++) {
395           for (CeedInt k = 0; k < 5; k++) {
396             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
397           }
398         }
399 
400         for (CeedInt j = 0; j < 5; j++) {
401           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]);
402         }
403         break;
404       case 2:  // SUPG is not implemented for explicit scheme
405         break;
406     }
407 
408   }  // End Quadrature Point Loop
409 
410   // Return
411   return 0;
412 }
413 
414 // *****************************************************************************
415 // This QFunction implements the Euler equations with (mentioned above)
416 //   with implicit time stepping method
417 //
418 // *****************************************************************************
419 CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
420   // Inputs
421   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],
422         (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
423   // Outputs
424   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
425 
426   EulerContext     context = (EulerContext)ctx;
427   const CeedScalar c_tau   = context->c_tau;
428   const CeedScalar gamma   = 1.4;
429 
430   CeedPragmaSIMD
431       // Quadrature Point Loop
432       for (CeedInt i = 0; i < Q; i++) {
433     // Setup
434     // -- Interp in
435     const CeedScalar rho      = q[0][i];
436     const CeedScalar u[3]     = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
437     const CeedScalar E        = q[4][i];
438     const CeedScalar drho[3]  = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
439     const CeedScalar dU[3][3] = {
440         {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
441         {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
442         {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
443     };
444     const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
445     // -- Interp-to-Interp q_data
446     const CeedScalar wdetJ = q_data[0][i];
447     // -- Interp-to-Grad q_data
448     // ---- Inverse of change of coordinate matrix: X_i,j
449     const CeedScalar dXdx[3][3] = {
450         {q_data[1][i], q_data[2][i], q_data[3][i]},
451         {q_data[4][i], q_data[5][i], q_data[6][i]},
452         {q_data[7][i], q_data[8][i], q_data[9][i]}
453     };
454     // dU/dx
455     CeedScalar drhodx[3]       = {0.};
456     CeedScalar dEdx[3]         = {0.};
457     CeedScalar dUdx[3][3]      = {{0.}};
458     CeedScalar dXdxdXdxT[3][3] = {{0.}};
459     for (CeedInt j = 0; j < 3; j++) {
460       for (CeedInt k = 0; k < 3; k++) {
461         drhodx[j] += drho[k] * dXdx[k][j];
462         dEdx[j] += dE[k] * dXdx[k][j];
463         for (CeedInt l = 0; l < 3; l++) {
464           dUdx[j][k] += dU[j][l] * dXdx[l][k];
465           dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l];  // dXdx_j,k * dXdx_k,j
466         }
467       }
468     }
469     const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic,
470                      P = E_internal * (gamma - 1.);  // P = pressure
471 
472     // The Physics
473     // Zero v and dv so all future terms can safely sum into it
474     for (CeedInt j = 0; j < 5; j++) {
475       v[j][i] = 0.;
476       for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
477     }
478     //-----mass matrix
479     for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i];
480 
481     // -- Density
482     // ---- u rho
483     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]);
484     // -- Momentum
485     // ---- rho (u x u) + P I3
486     for (CeedInt j = 0; j < 3; j++) {
487       for (CeedInt k = 0; k < 3; k++) {
488         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] +
489                                     (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
490       }
491     }
492     // -- Total Energy Density
493     // ---- (E + P) u
494     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]);
495 
496     // -- Stabilization terms
497     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
498     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
499     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
500 
501     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
502     CeedScalar dqdx[5][3];
503     for (CeedInt j = 0; j < 3; j++) {
504       dqdx[0][j] = drhodx[j];
505       dqdx[4][j] = dEdx[j];
506       for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
507     }
508 
509     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
510     CeedScalar strong_conv[5] = {0.};
511     for (CeedInt j = 0; j < 3; j++) {
512       for (CeedInt k = 0; k < 5; k++) {
513         for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
514       }
515     }
516 
517     // ---- Strong residual
518     CeedScalar strong_res[5];
519     for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j];
520 
521     // Stabilization
522     // -- Tau elements
523     const CeedScalar sound_speed = sqrt(gamma * P / rho);
524     CeedScalar       Tau_x[3]    = {0.};
525     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
526 
527     // -- Stabilization method: none, SU, or SUPG
528     CeedScalar stab[5][3] = {{0.}};
529     switch (context->stabilization) {
530       case 0:  // Galerkin
531         break;
532       case 1:  // SU
533         for (CeedInt j = 0; j < 3; j++) {
534           for (CeedInt k = 0; k < 5; k++) {
535             for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
536           }
537         }
538 
539         for (CeedInt j = 0; j < 5; j++) {
540           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]);
541         }
542         break;
543       case 2:  // SUPG
544         for (CeedInt j = 0; j < 3; j++) {
545           for (CeedInt k = 0; k < 5; k++) {
546             for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l];
547           }
548         }
549 
550         for (CeedInt j = 0; j < 5; j++) {
551           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]);
552         }
553         break;
554     }
555   }  // End Quadrature Point Loop
556 
557   // Return
558   return 0;
559 }
560 // *****************************************************************************
561 // This QFunction sets the inflow boundary conditions for
562 //   the traveling vortex problem.
563 //
564 //  Prescribed T_inlet and P_inlet are converted to conservative variables
565 //      and applied weakly.
566 //
567 // *****************************************************************************
568 CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
569   // Inputs
570   const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
571   // Outputs
572   CeedScalar(*v)[CEED_Q_VLA]     = (CeedScalar(*)[CEED_Q_VLA])out[0];
573   EulerContext     context       = (EulerContext)ctx;
574   const int        euler_test    = context->euler_test;
575   const bool       implicit      = context->implicit;
576   CeedScalar      *mean_velocity = context->mean_velocity;
577   const CeedScalar cv            = 2.5;
578   const CeedScalar R             = 1.;
579   CeedScalar       T_inlet;
580   CeedScalar       P_inlet;
581 
582   // For test cases 1 and 3 the background velocity is zero
583   if (euler_test == 1 || euler_test == 3) {
584     for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.;
585   }
586 
587   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
588   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
589   else T_inlet = P_inlet = 1.;
590 
591   CeedPragmaSIMD
592       // Quadrature Point Loop
593       for (CeedInt i = 0; i < Q; i++) {
594     // Setup
595     // -- Interp-to-Interp q_data
596     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
597     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
598     // We can effect this by swapping the sign on this weight
599     const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i];
600     // ---- Normal vect
601     const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]};
602 
603     // face_normal = Normal vector of the face
604     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
605     // The Physics
606     // Zero v so all future terms can safely sum into it
607     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.;
608 
609     // Implementing in/outflow BCs
610     if (face_normal > 0) {
611     } else {  // inflow
612       const CeedScalar rho_inlet       = P_inlet / (R * T_inlet);
613       const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.;
614       // incoming total energy
615       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
616 
617       // The Physics
618       // -- Density
619       v[0][i] -= wdetJb * rho_inlet * face_normal;
620 
621       // -- Momentum
622       for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho_inlet * face_normal * mean_velocity[j] + norm[j] * P_inlet);
623 
624       // -- Total Energy Density
625       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
626     }
627 
628   }  // End Quadrature Point Loop
629   return 0;
630 }
631 
632 // *****************************************************************************
633 // This QFunction sets the outflow boundary conditions for
634 //   the Euler solver.
635 //
636 //  Outflow BCs:
637 //    The validity of the weak form of the governing equations is
638 //      extended to the outflow.
639 //
640 // *****************************************************************************
641 CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
642   // Inputs
643   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];
644   // Outputs
645   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
646   EulerContext context       = (EulerContext)ctx;
647   const bool   implicit      = context->implicit;
648   CeedScalar  *mean_velocity = context->mean_velocity;
649 
650   const CeedScalar gamma = 1.4;
651 
652   CeedPragmaSIMD
653       // Quadrature Point Loop
654       for (CeedInt i = 0; i < Q; i++) {
655     // Setup
656     // -- Interp in
657     const CeedScalar rho  = q[0][i];
658     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
659     const CeedScalar E    = q[4][i];
660 
661     // -- Interp-to-Interp q_data
662     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
663     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
664     // We can effect this by swapping the sign on this weight
665     const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i];
666     // ---- Normal vectors
667     const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]};
668 
669     // face_normal = Normal vector of the face
670     const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
671     // The Physics
672     // Zero v so all future terms can safely sum into it
673     for (CeedInt j = 0; j < 5; j++) v[j][i] = 0;
674 
675     // Implementing in/outflow BCs
676     if (face_normal > 0) {  // outflow
677       const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.;
678       const CeedScalar P         = (E - E_kinetic * rho) * (gamma - 1.);              // pressure
679       const CeedScalar u_normal  = norm[0] * u[0] + norm[1] * u[1] + norm[2] * u[2];  // Normal velocity
680       // The Physics
681       // -- Density
682       v[0][i] -= wdetJb * rho * u_normal;
683 
684       // -- Momentum
685       for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P);
686 
687       // -- Total Energy Density
688       v[4][i] -= wdetJb * u_normal * (E + P);
689     }
690   }  // End Quadrature Point Loop
691   return 0;
692 }
693 
694 // *****************************************************************************
695 
696 #endif  // eulervortex_h
697