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