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