xref: /honee/qfunctions/advection.h (revision 5a275bb67bb8312ddcb117417b4aaae9b609ecba)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 
4 /// @file
5 /// Advection initial condition and operator for Navier-Stokes example using PETSc
6 #include <ceed.h>
7 #include <math.h>
8 
9 #include "advection_types.h"
10 #include "newtonian_state.h"
11 #include "newtonian_types.h"
12 #include "stabilization_types.h"
13 #include "utils.h"
14 
15 // *****************************************************************************
16 // This QFunction sets the initial conditions and the boundary conditions
17 //   for two test cases: ROTATION and TRANSLATION
18 //
19 // -- ROTATION (default)
20 //      Initial Conditions:
21 //        Mass Density:
22 //          Constant mass density of 1.0
23 //        Momentum Density:
24 //          Rotational field in x,y
25 //        Energy Density:
26 //          Maximum of 1. x0 decreasing linearly to 0. as radial distance
27 //            increases to (1.-r/rc), then 0. everywhere else
28 //
29 //      Boundary Conditions:
30 //        Mass Density:
31 //          0.0 flux
32 //        Momentum Density:
33 //          0.0
34 //        Energy Density:
35 //          0.0 flux
36 //
37 // -- TRANSLATION
38 //      Initial Conditions:
39 //        Mass Density:
40 //          Constant mass density of 1.0
41 //        Momentum Density:
42 //           Constant rectilinear field in x,y
43 //        Energy Density:
44 //          Maximum of 1. x0 decreasing linearly to 0. as radial distance
45 //            increases to (1.-r/rc), then 0. everywhere else
46 //
47 //      Boundary Conditions:
48 //        Mass Density:
49 //          0.0 flux
50 //        Momentum Density:
51 //          0.0
52 //        Energy Density:
53 //          Inflow BCs:
54 //            E = E_wind
55 //          Outflow BCs:
56 //            E = E(boundary)
57 //          Both In/Outflow BCs for E are applied weakly in the
58 //            QFunction "Advection2d_Sur"
59 //
60 // *****************************************************************************
61 
62 // *****************************************************************************
63 // This helper function provides the exact, time-dependent solution and IC formulation for 2D advection
64 // *****************************************************************************
65 CEED_QFUNCTION_HELPER int Exact_AdvectionGeneric(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
66   const SetupContextAdv context = (SetupContextAdv)ctx;
67   const CeedScalar      rc      = context->rc;
68   const CeedScalar      lx      = context->lx;
69   const CeedScalar      ly      = context->ly;
70   const CeedScalar      lz      = dim == 2 ? 0. : context->lz;
71   const CeedScalar     *wind    = context->wind;
72 
73   const CeedScalar center[3] = {0.5 * lx, 0.5 * ly, 0.5 * lz};
74   const CeedScalar theta     = dim == 2 ? M_PI / 3 : M_PI;
75   const CeedScalar x0[3]     = {center[0] + .25 * lx * cos(theta + time), center[1] + .25 * ly * sin(theta + time), 0.5 * lz};
76 
77   const CeedScalar x = X[0], y = X[1], z = dim == 2 ? 0. : X[2];
78 
79   switch (context->wind_type) {
80     case ADVDIF_WIND_ROTATION:
81       q[0] = 1.;
82       q[1] = -(y - center[1]);
83       q[2] = (x - center[0]);
84       q[3] = 0;
85       break;
86     case ADVDIF_WIND_TRANSLATION:
87       q[0] = 1.;
88       q[1] = wind[0];
89       q[2] = wind[1];
90       q[3] = dim == 2 ? 0. : wind[2];
91       break;
92     default:
93       return 1;
94   }
95 
96   switch (context->initial_condition_type) {
97     case ADVDIF_IC_BUBBLE_SPHERE:
98     case ADVDIF_IC_BUBBLE_CYLINDER: {
99       CeedScalar r = sqrt(Square(x - x0[0]) + Square(y - x0[1]) + Square(z - x0[2]));
100 
101       switch (context->bubble_continuity_type) {
102         // original continuous, smooth shape
103         case ADVDIF_BUBBLE_CONTINUITY_SMOOTH:
104           q[4] = r <= rc ? (1. - r / rc) : 0.;
105           break;
106         // discontinuous, sharp back half shape
107         case ADVDIF_BUBBLE_CONTINUITY_BACK_SHARP:
108           q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) : 0.;
109           break;
110         // attempt to define a finite thickness that will get resolved under grid refinement
111         case ADVDIF_BUBBLE_CONTINUITY_THICK:
112           q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) * fmin(1.0, (center[1] - y) / 1.25) : 0.;
113           break;
114         case ADVDIF_BUBBLE_CONTINUITY_COSINE:
115           q[4] = r <= rc ? .5 + .5 * cos(r * M_PI / rc) : 0;
116           break;
117       }
118       break;
119     }
120 
121     case ADVDIF_IC_COSINE_HILL: {
122       CeedScalar r          = sqrt(Square(x - center[0]) + Square(y - center[1]));
123       CeedScalar half_width = context->lx / 2;
124       q[4]                  = r > half_width ? 0. : cos(2 * M_PI * r / half_width + M_PI) + 1.;
125     } break;
126 
127     case ADVDIF_IC_SKEW: {
128       CeedScalar       skewed_barrier[3]  = {wind[0], wind[1], 0};
129       CeedScalar       inflow_to_point[3] = {x - context->lx / 2, y, 0};
130       CeedScalar       cross_product[3]   = {0};
131       const CeedScalar boundary_threshold = 20 * CEED_EPSILON;
132       Cross3(skewed_barrier, inflow_to_point, cross_product);
133 
134       q[4] = cross_product[2] > boundary_threshold ? 0 : 1;
135       if ((x < boundary_threshold && wind[0] < boundary_threshold) ||                // outflow at -x boundary
136           (y < boundary_threshold && wind[1] < boundary_threshold) ||                // outflow at -y boundary
137           (x > context->lx - boundary_threshold && wind[0] > boundary_threshold) ||  // outflow at +x boundary
138           (y > context->ly - boundary_threshold && wind[1] > boundary_threshold)     // outflow at +y boundary
139       ) {
140         q[4] = 0;
141       }
142     } break;
143 
144     case ADVDIF_IC_WAVE: {
145       CeedScalar theta = context->wave_frequency * DotN(X, wind, dim) + context->wave_phase;
146       switch (context->wave_type) {
147         case ADVDIF_WAVE_SINE:
148           q[4] = sin(theta);
149           break;
150         case ADVDIF_WAVE_SQUARE:
151           q[4] = sin(theta) > 100 * CEED_EPSILON ? 1 : -1;
152           break;
153       }
154     }
155   }
156   return 0;
157 }
158 
159 // *****************************************************************************
160 // This QFunction sets the initial conditions for 3D advection
161 // *****************************************************************************
162 CEED_QFUNCTION(ICsAdvection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
163   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
164   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
165 
166   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
167     const CeedScalar x[]  = {X[0][i], X[1][i], X[2][i]};
168     CeedScalar       q[5] = {0.};
169 
170     Exact_AdvectionGeneric(3, 0., x, 5, q, ctx);
171     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
172   }
173   return 0;
174 }
175 
176 // *****************************************************************************
177 // This QFunction sets the initial conditions for 2D advection
178 // *****************************************************************************
179 CEED_QFUNCTION(ICsAdvection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
180   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
181   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
182   const SetupContextAdv context    = (SetupContextAdv)ctx;
183 
184   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
185     const CeedScalar x[]  = {X[0][i], X[1][i]};
186     CeedScalar       q[5] = {0.};
187 
188     Exact_AdvectionGeneric(2, context->time, x, 5, q, ctx);
189     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
190   }
191   return 0;
192 }
193 
194 CEED_QFUNCTION_HELPER void StatePhysicalGradientFromReference_ND(CeedInt N, CeedInt Q, CeedInt i, NewtonianIdealGasContext gas, State s,
195                                                                  StateVariable state_var, const CeedScalar *grad_q, const CeedScalar *dXdx,
196                                                                  State *grad_s) {
197   switch (N) {
198     case 2: {
199       for (CeedInt k = 0; k < 2; k++) {
200         CeedScalar dqi[5];
201         for (CeedInt j = 0; j < 5; j++) {
202           dqi[j] = grad_q[(Q * 5) * 0 + Q * j + i] * dXdx[0 * N + k] + grad_q[(Q * 5) * 1 + Q * j + i] * dXdx[1 * N + k];
203         }
204         grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var);
205       }
206       CeedScalar U[5] = {0.};
207       grad_s[2]       = StateFromU(gas, U);
208     } break;
209     case 3:
210       // Cannot directly use StatePhysicalGradientFromReference helper functions due to SYCL online compiler incompatabilities
211       for (CeedInt k = 0; k < 3; k++) {
212         CeedScalar dqi[5];
213         for (CeedInt j = 0; j < 5; j++) {
214           dqi[j] = grad_q[(Q * 5) * 0 + Q * j + i] * dXdx[0 * N + k] + grad_q[(Q * 5) * 1 + Q * j + i] * dXdx[1 * N + k] +
215                    grad_q[(Q * 5) * 2 + Q * j + i] * dXdx[2 * N + k];
216         }
217         grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var);
218       }
219       break;
220   }
221 }
222 
223 // @brief Calculate the stabilization constant \tau
224 CEED_QFUNCTION_HELPER CeedScalar Tau(AdvectionContext context, const State s, const CeedScalar *dXdx, CeedInt dim) {
225   switch (context->stabilization_tau) {
226     case STAB_TAU_CTAU: {
227       CeedScalar uX[3] = {0.};
228 
229       MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
230       return context->CtauS / sqrt(DotN(uX, uX, dim));
231     } break;
232     case STAB_TAU_ADVDIFF_SHAKIB: {
233       CeedScalar gijd_mat[9] = {0.}, gij_uj[3] = {0.};
234 
235       MatMatN(dXdx, dXdx, dim, CEED_TRANSPOSE, CEED_NOTRANSPOSE, gijd_mat);
236       MatVecNM(gijd_mat, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, gij_uj);
237       return 1 / sqrt(Square(2 * context->Ctau_t / context->dt) + DotN(s.Y.velocity, gij_uj, dim) * context->Ctau_a +
238                       Square(context->diffusion_coeff) * DotN(gijd_mat, gijd_mat, dim * dim) * context->Ctau_d);
239     } break;
240     default:
241       return 0.;
242   }
243 }
244 
245 // *****************************************************************************
246 // This QFunction implements Advection for implicit time stepping method
247 // *****************************************************************************
248 CEED_QFUNCTION_HELPER int IFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
249   AdvectionContext context = (AdvectionContext)ctx;
250 
251   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
252   const CeedScalar(*grad_q)            = in[1];
253   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
254   const CeedScalar(*q_data)            = in[3];
255   const CeedScalar(*divFdiff)          = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? in[5] : NULL;
256 
257   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
258   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
259 
260   NewtonianIdealGasContext         gas;
261   struct NewtonianIdealGasContext_ gas_struct = {0};
262   gas                                         = &gas_struct;
263 
264   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
265     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
266     const State      s     = StateFromU(gas, qi);
267 
268     CeedScalar wdetJ, dXdx[9];
269     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
270     State grad_s[3];
271     StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
272 
273     const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
274 
275     for (CeedInt f = 0; f < 4; f++) {
276       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
277       v[f][i] = wdetJ * q_dot[f][i];                          // K Mass/transient term
278     }
279 
280     CeedScalar div_u = 0;
281     for (CeedInt j = 0; j < dim; j++) {
282       for (CeedInt k = 0; k < dim; k++) {
283         div_u += grad_s[k].Y.velocity[j];
284       }
285     }
286     CeedScalar uX[3] = {0.};
287     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
288     CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
289 
290     v[4][i] = wdetJ * q_dot[4][i];  // transient part (ALWAYS)
291     if (context->strong_form) {     // Strong Galerkin convection term: v div(E u)
292       v[4][i] += wdetJ * strong_conv;
293     } else {  // Weak Galerkin convection term: -dv \cdot (E u)
294       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = -wdetJ * s.U.E_total * uX[j];
295     }
296 
297     {  // Diffusion
298       CeedScalar Fe[3], Fe_dXdx[3] = {0.};
299 
300       for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total;
301       MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx);
302       for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] -= wdetJ * Fe_dXdx[k];
303     }
304 
305     const CeedScalar TauS = Tau(context, s, dXdx, dim);
306     for (CeedInt j = 0; j < dim; j++) {
307       switch (context->stabilization) {
308         case STAB_NONE:
309           break;
310         case STAB_SU:
311           grad_v[j][4][i] += wdetJ * TauS * uX[j] * strong_conv;
312           break;
313         case STAB_SUPG: {
314           CeedScalar divFdiff_i = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? divFdiff[i] : 0.;
315           grad_v[j][4][i] += wdetJ * TauS * uX[j] * (q_dot[4][i] + strong_conv + divFdiff_i);
316         } break;
317       }
318     }
319   }
320   return 0;
321 }
322 
323 CEED_QFUNCTION(IFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
324   return IFunction_AdvectionGeneric(ctx, Q, in, out, 3);
325 }
326 
327 CEED_QFUNCTION(IFunction_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
328   return IFunction_AdvectionGeneric(ctx, Q, in, out, 2);
329 }
330 
331 CEED_QFUNCTION_HELPER int MassFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
332   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
333   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[1];
334   const CeedScalar(*q_data)            = in[2];
335 
336   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
337   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
338 
339   AdvectionContext                 context    = (AdvectionContext)ctx;
340   struct NewtonianIdealGasContext_ gas_struct = {0};
341   NewtonianIdealGasContext         gas        = &gas_struct;
342 
343   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
344     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
345     const State      s     = StateFromU(gas, qi);
346     CeedScalar       wdetJ, dXdx[9];
347     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
348 
349     for (CeedInt f = 0; f < 4; f++) {
350       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
351       v[f][i] = wdetJ * q_dot[f][i];                          // K Mass/transient term
352     }
353 
354     // Unstabilized mass term
355     v[4][i] = wdetJ * q_dot[4][i];
356 
357     // Stabilized mass term
358     CeedScalar uX[3] = {0.};
359     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
360     const CeedScalar TauS = Tau(context, s, dXdx, dim);
361     for (CeedInt j = 0; j < dim; j++) {
362       switch (context->stabilization) {
363         case STAB_NONE:
364         case STAB_SU:
365           grad_v[j][4][i] = 0;
366           break;  // These should be run with the unstabilized mass matrix anyways
367         case STAB_SUPG:
368           grad_v[j][4][i] = wdetJ * TauS * q_dot[4][i] * uX[j];
369           break;
370       }
371     }
372   }
373   return 0;
374 }
375 
376 CEED_QFUNCTION(MassFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
377   return MassFunction_AdvectionGeneric(ctx, Q, in, out, 3);
378 }
379 
380 CEED_QFUNCTION(MassFunction_Advection2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
381   return MassFunction_AdvectionGeneric(ctx, Q, in, out, 2);
382 }
383 
384 // *****************************************************************************
385 // This QFunction implements Advection for explicit time stepping method
386 // *****************************************************************************
387 CEED_QFUNCTION_HELPER int RHSFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
388   AdvectionContext context = (AdvectionContext)ctx;
389 
390   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
391   const CeedScalar(*grad_q)        = in[1];
392   const CeedScalar(*q_data)        = in[2];
393   const CeedScalar(*divFdiff)      = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? in[4] : NULL;
394 
395   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
396   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
397 
398   struct NewtonianIdealGasContext_ gas_struct = {0};
399   NewtonianIdealGasContext         gas        = &gas_struct;
400 
401   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
402     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
403     const State      s     = StateFromU(gas, qi);
404 
405     CeedScalar wdetJ, dXdx[9];
406     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
407     State grad_s[3];
408     StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
409 
410     const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
411 
412     for (CeedInt f = 0; f < 4; f++) {
413       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
414       v[f][i] = 0.;
415     }
416 
417     CeedScalar div_u = 0;
418     for (CeedInt j = 0; j < dim; j++) {
419       for (CeedInt k = 0; k < dim; k++) {
420         div_u += grad_s[k].Y.velocity[j];
421       }
422     }
423     CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
424 
425     CeedScalar uX[3] = {0.};
426     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
427 
428     if (context->strong_form) {  // Strong Galerkin convection term: v div(E u)
429       v[4][i] = -wdetJ * strong_conv;
430       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = 0;
431     } else {  // Weak Galerkin convection term: -dv \cdot (E u)
432       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = wdetJ * s.U.E_total * uX[j];
433       v[4][i] = 0.;
434     }
435 
436     {  // Diffusion
437       CeedScalar Fe[3], Fe_dXdx[3] = {0.};
438 
439       for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total;
440       MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx);
441       for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] += wdetJ * Fe_dXdx[k];
442     }
443 
444     const CeedScalar TauS = Tau(context, s, dXdx, dim);
445     for (CeedInt j = 0; j < dim; j++) {
446       switch (context->stabilization) {
447         case STAB_NONE:
448           break;
449         case STAB_SU:
450         case STAB_SUPG: {
451           CeedScalar divFdiff_i = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? divFdiff[i] : 0.;
452           grad_v[j][4][i] -= wdetJ * TauS * (strong_conv + divFdiff_i) * uX[j];
453         } break;
454       }
455     }
456   }
457   return 0;
458 }
459 
460 CEED_QFUNCTION(RHS_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
461   return RHSFunction_AdvectionGeneric(ctx, Q, in, out, 3);
462 }
463 
464 CEED_QFUNCTION(RHS_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
465   return RHSFunction_AdvectionGeneric(ctx, Q, in, out, 2);
466 }
467 
468 // *****************************************************************************
469 // This QFunction implements consistent outflow and inflow BCs
470 //      for advection
471 //
472 //  Inflow and outflow faces are determined based on sign(dot(wind, normal)):
473 //    sign(dot(wind, normal)) > 0 : outflow BCs
474 //    sign(dot(wind, normal)) < 0 : inflow BCs
475 //
476 //  Outflow BCs:
477 //    The validity of the weak form of the governing equations is extended to the outflow and the current values of E are applied.
478 //
479 //  Inflow BCs:
480 //    A prescribed Total Energy (E_wind) is applied weakly.
481 // *****************************************************************************
482 CEED_QFUNCTION_HELPER int Advection_InOutFlowGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
483   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
484   const CeedScalar(*q_data_sur)    = in[2];
485 
486   CeedScalar(*v)[CEED_Q_VLA]   = (CeedScalar(*)[CEED_Q_VLA])out[0];
487   AdvectionContext context     = (AdvectionContext)ctx;
488   const CeedScalar E_wind      = context->E_wind;
489   const CeedScalar strong_form = context->strong_form;
490   const bool       is_implicit = context->implicit;
491 
492   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
493     const CeedScalar rho  = q[0][i];
494     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
495     const CeedScalar E    = q[4][i];
496 
497     CeedScalar wdetJb, normal[3];
498     QdataBoundaryUnpack_ND(dim, Q, i, q_data_sur, &wdetJb, NULL, normal);
499     wdetJb *= is_implicit ? -1. : 1.;
500 
501     const CeedScalar u_normal = DotN(normal, u, dim);
502 
503     // No Change in density or momentum
504     for (CeedInt j = 0; j < 4; j++) {
505       v[j][i] = 0;
506     }
507     // Implementing in/outflow BCs
508     if (u_normal > 0) {  // outflow
509       v[4][i] = -(1 - strong_form) * wdetJb * E * u_normal;
510     } else {  // inflow
511       v[4][i] = -(1 - strong_form) * wdetJb * E_wind * u_normal;
512     }
513   }
514   return 0;
515 }
516 
517 CEED_QFUNCTION(Advection_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
518   return Advection_InOutFlowGeneric(ctx, Q, in, out, 3);
519 }
520 
521 CEED_QFUNCTION(Advection2d_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
522   return Advection_InOutFlowGeneric(ctx, Q, in, out, 2);
523 }
524 
525 // @brief Volume integral for RHS of divergence of diffusive flux direct projection
526 CEED_QFUNCTION_HELPER int DivDiffusiveFluxVolumeRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out,
527                                                                    const CeedInt dim) {
528   const CeedScalar(*Grad_q)       = in[0];
529   const CeedScalar(*q_data)       = in[1];
530   CeedScalar(*Grad_v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
531 
532   AdvectionContext context = (AdvectionContext)ctx;
533 
534   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
535     CeedScalar wdetJ, dXdx[9], F_diff[3] = {0.};
536 
537     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
538     {  // Get physical diffusive flux
539       CeedScalar Grad_qn[15], grad_E_ref[3];
540 
541       GradUnpackN(Q, i, 5, dim, Grad_q, Grad_qn);
542       CopyN(&Grad_qn[4 * dim], grad_E_ref, dim);
543       MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff);
544       ScaleN(F_diff, -context->diffusion_coeff, dim);
545     }
546 
547     CeedScalar F_diff_dXdx[3] = {0.};
548     MatVecNM(dXdx, F_diff, dim, dim, CEED_NOTRANSPOSE, F_diff_dXdx);
549     for (CeedInt k = 0; k < dim; k++) Grad_v[k][i] = -wdetJ * F_diff_dXdx[k];
550   }
551   return 0;
552 }
553 
554 CEED_QFUNCTION(DivDiffusiveFluxVolumeRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
555   return DivDiffusiveFluxVolumeRHS_AdvDif_Generic(ctx, Q, in, out, 2);
556 }
557 
558 CEED_QFUNCTION(DivDiffusiveFluxVolumeRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
559   return DivDiffusiveFluxVolumeRHS_AdvDif_Generic(ctx, Q, in, out, 3);
560 }
561 
562 // @brief Boundary integral for RHS of divergence of diffusive flux direct projection
563 CEED_QFUNCTION_HELPER int DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out,
564                                                                      const CeedInt dim) {
565   const CeedScalar(*Grad_q) = in[0];
566   const CeedScalar(*q_data) = in[1];
567   CeedScalar(*v)            = out[0];
568 
569   AdvectionContext context = (AdvectionContext)ctx;
570 
571   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
572     CeedScalar wdetJ, normal[3], dXdx[9], F_diff[3] = {0.};
573 
574     QdataBoundaryGradientUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx, normal);
575     {  // Get physical diffusive flux
576       CeedScalar Grad_qn[15], grad_E_ref[3];
577 
578       GradUnpackN(Q, i, 5, dim, Grad_q, Grad_qn);
579       CopyN(&Grad_qn[4 * dim], grad_E_ref, dim);
580       MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff);
581       ScaleN(F_diff, -context->diffusion_coeff, dim);
582     }
583 
584     v[i] = wdetJ * DotN(F_diff, normal, dim);
585   }
586   return 0;
587 }
588 
589 CEED_QFUNCTION(DivDiffusiveFluxBoundaryRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
590   return DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(ctx, Q, in, out, 2);
591 }
592 
593 CEED_QFUNCTION(DivDiffusiveFluxBoundaryRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
594   return DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(ctx, Q, in, out, 3);
595 }
596 
597 // @brief Volume integral for RHS of diffusive flux indirect projection
598 CEED_QFUNCTION_HELPER int DiffusiveFluxRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out,
599                                                           const CeedInt dim) {
600   const CeedScalar(*Grad_q)  = in[0];
601   const CeedScalar(*q_data)  = in[1];
602   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
603 
604   AdvectionContext context = (AdvectionContext)ctx;
605 
606   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
607     CeedScalar wdetJ, dXdx[9], F_diff[3] = {0.};
608 
609     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
610     {  // Get physical diffusive flux
611       CeedScalar Grad_qn[15], grad_E_ref[3];
612 
613       GradUnpackN(Q, i, 5, dim, Grad_q, Grad_qn);
614       CopyN(&Grad_qn[4 * dim], grad_E_ref, dim);
615       MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff);
616       ScaleN(F_diff, -context->diffusion_coeff, dim);
617     }
618     for (CeedInt k = 0; k < dim; k++) v[k][i] = wdetJ * F_diff[k];
619   }
620   return 0;
621 }
622 
623 CEED_QFUNCTION(DiffusiveFluxRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
624   return DiffusiveFluxRHS_AdvDif_Generic(ctx, Q, in, out, 2);
625 }
626 
627 CEED_QFUNCTION(DiffusiveFluxRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
628   return DiffusiveFluxRHS_AdvDif_Generic(ctx, Q, in, out, 3);
629 }
630