xref: /honee/qfunctions/advection.h (revision da02a6e7947cab8e173fd83bf2b4eda896c06b14)
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 CeedInt 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 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 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 ADVECTIONIC_BUBBLE_SPHERE:
98     case ADVECTIONIC_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 BUBBLE_CONTINUITY_SMOOTH:
104           q[4] = r <= rc ? (1. - r / rc) : 0.;
105           break;
106         // discontinuous, sharp back half shape
107         case 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 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 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 ADVECTIONIC_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 ADVECTIONIC_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 ADVECTIONIC_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     } break;
239     default:
240       return 0.;
241   }
242 }
243 
244 // *****************************************************************************
245 // This QFunction implements Advection for implicit time stepping method
246 // *****************************************************************************
247 CEED_QFUNCTION_HELPER void IFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
248   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
249   const CeedScalar(*grad_q)            = in[1];
250   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
251   const CeedScalar(*q_data)            = in[3];
252 
253   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
254   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
255 
256   AdvectionContext                 context = (AdvectionContext)ctx;
257   NewtonianIdealGasContext         gas;
258   struct NewtonianIdealGasContext_ gas_struct = {0};
259   gas                                         = &gas_struct;
260 
261   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
262     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
263     const State      s     = StateFromU(gas, qi);
264 
265     CeedScalar wdetJ, dXdx[9];
266     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
267     State grad_s[3];
268     StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
269 
270     const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
271 
272     for (CeedInt f = 0; f < 4; f++) {
273       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
274       v[f][i] = wdetJ * q_dot[f][i];                          // K Mass/transient term
275     }
276 
277     CeedScalar div_u = 0;
278     for (CeedInt j = 0; j < dim; j++) {
279       for (CeedInt k = 0; k < dim; k++) {
280         div_u += grad_s[k].Y.velocity[j];
281       }
282     }
283     CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
284     CeedScalar strong_res  = q_dot[4][i] + strong_conv;
285 
286     v[4][i] = wdetJ * q_dot[4][i];  // transient part (ALWAYS)
287 
288     CeedScalar uX[3] = {0.};
289     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
290 
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++) switch (context->stabilization) {
307         case STAB_NONE:
308           break;
309         case STAB_SU:
310           grad_v[j][4][i] += wdetJ * TauS * strong_conv * uX[j];
311           break;
312         case STAB_SUPG:
313           grad_v[j][4][i] += wdetJ * TauS * strong_res * uX[j];
314           break;
315       }
316   }
317 }
318 
319 CEED_QFUNCTION(IFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
320   IFunction_AdvectionGeneric(ctx, Q, in, out, 3);
321   return 0;
322 }
323 
324 CEED_QFUNCTION(IFunction_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
325   IFunction_AdvectionGeneric(ctx, Q, in, out, 2);
326   return 0;
327 }
328 
329 CEED_QFUNCTION_HELPER void MassFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
330   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
331   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[1];
332   const CeedScalar(*q_data)            = in[2];
333 
334   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
335   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
336 
337   AdvectionContext                 context    = (AdvectionContext)ctx;
338   struct NewtonianIdealGasContext_ gas_struct = {0};
339   NewtonianIdealGasContext         gas        = &gas_struct;
340 
341   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
342     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
343     const State      s     = StateFromU(gas, qi);
344     CeedScalar       wdetJ, dXdx[9];
345     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
346 
347     for (CeedInt f = 0; f < 4; f++) {
348       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
349       v[f][i] = wdetJ * q_dot[f][i];                          // K Mass/transient term
350     }
351 
352     // Unstabilized mass term
353     v[4][i] = wdetJ * q_dot[4][i];
354 
355     // Stabilized mass term
356     CeedScalar uX[3] = {0.};
357     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
358     const CeedScalar TauS = Tau(context, s, dXdx, dim);
359     for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) {
360         case STAB_NONE:
361         case STAB_SU:
362           grad_v[j][4][i] = 0;
363           break;  // These should be run with the unstabilized mass matrix anyways
364         case STAB_SUPG:
365           grad_v[j][4][i] = wdetJ * TauS * q_dot[4][i] * uX[j];
366           break;
367       }
368   }
369 }
370 
371 CEED_QFUNCTION(MassFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
372   MassFunction_AdvectionGeneric(ctx, Q, in, out, 3);
373   return 0;
374 }
375 
376 CEED_QFUNCTION(MassFunction_Advection2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
377   MassFunction_AdvectionGeneric(ctx, Q, in, out, 2);
378   return 0;
379 }
380 
381 // *****************************************************************************
382 // This QFunction implements Advection for explicit time stepping method
383 // *****************************************************************************
384 CEED_QFUNCTION_HELPER void RHSFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
385   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
386   const CeedScalar(*grad_q)        = in[1];
387   const CeedScalar(*q_data)        = in[2];
388 
389   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
390   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
391 
392   AdvectionContext                 context    = (AdvectionContext)ctx;
393   struct NewtonianIdealGasContext_ gas_struct = {0};
394   NewtonianIdealGasContext         gas        = &gas_struct;
395 
396   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
397     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
398     const State      s     = StateFromU(gas, qi);
399 
400     CeedScalar wdetJ, dXdx[9];
401     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
402     State grad_s[3];
403     StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
404 
405     const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
406 
407     for (CeedInt f = 0; f < 4; f++) {
408       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
409       v[f][i] = 0.;
410     }
411 
412     CeedScalar div_u = 0;
413     for (CeedInt j = 0; j < dim; j++) {
414       for (CeedInt k = 0; k < dim; k++) {
415         div_u += grad_s[k].Y.velocity[j];
416       }
417     }
418     CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
419 
420     CeedScalar uX[3] = {0.};
421     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
422 
423     if (context->strong_form) {  // Strong Galerkin convection term: v div(E u)
424       v[4][i] = -wdetJ * strong_conv;
425       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = 0;
426     } else {  // Weak Galerkin convection term: -dv \cdot (E u)
427       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = wdetJ * s.U.E_total * uX[j];
428       v[4][i] = 0.;
429     }
430 
431     {  // Diffusion
432       CeedScalar Fe[3], Fe_dXdx[3] = {0.};
433 
434       for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total;
435       MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx);
436       for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] += wdetJ * Fe_dXdx[k];
437     }
438 
439     const CeedScalar TauS = Tau(context, s, dXdx, dim);
440     for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) {
441         case STAB_NONE:
442           break;
443         case STAB_SU:
444         case STAB_SUPG:
445           grad_v[j][4][i] -= wdetJ * TauS * strong_conv * uX[j];
446           break;
447       }
448   }
449 }
450 
451 CEED_QFUNCTION(RHS_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
452   RHSFunction_AdvectionGeneric(ctx, Q, in, out, 3);
453   return 0;
454 }
455 
456 CEED_QFUNCTION(RHS_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
457   RHSFunction_AdvectionGeneric(ctx, Q, in, out, 2);
458   return 0;
459 }
460 
461 // *****************************************************************************
462 // This QFunction implements consistent outflow and inflow BCs
463 //      for advection
464 //
465 //  Inflow and outflow faces are determined based on sign(dot(wind, normal)):
466 //    sign(dot(wind, normal)) > 0 : outflow BCs
467 //    sign(dot(wind, normal)) < 0 : inflow BCs
468 //
469 //  Outflow BCs:
470 //    The validity of the weak form of the governing equations is extended to the outflow and the current values of E are applied.
471 //
472 //  Inflow BCs:
473 //    A prescribed Total Energy (E_wind) is applied weakly.
474 // *****************************************************************************
475 CEED_QFUNCTION(Advection_InOutFlowGeneric)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
476   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
477   const CeedScalar(*q_data_sur)    = in[2];
478 
479   CeedScalar(*v)[CEED_Q_VLA]   = (CeedScalar(*)[CEED_Q_VLA])out[0];
480   AdvectionContext context     = (AdvectionContext)ctx;
481   const CeedScalar E_wind      = context->E_wind;
482   const CeedScalar strong_form = context->strong_form;
483   const bool       is_implicit = context->implicit;
484 
485   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
486     const CeedScalar rho  = q[0][i];
487     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
488     const CeedScalar E    = q[4][i];
489 
490     CeedScalar wdetJb, normal[3];
491     QdataBoundaryUnpack_ND(dim, Q, i, q_data_sur, &wdetJb, NULL, normal);
492     wdetJb *= is_implicit ? -1. : 1.;
493 
494     const CeedScalar u_normal = DotN(normal, u, dim);
495 
496     // No Change in density or momentum
497     for (CeedInt j = 0; j < 4; j++) {
498       v[j][i] = 0;
499     }
500     // Implementing in/outflow BCs
501     if (u_normal > 0) {  // outflow
502       v[4][i] = -(1 - strong_form) * wdetJb * E * u_normal;
503     } else {  // inflow
504       v[4][i] = -(1 - strong_form) * wdetJb * E_wind * u_normal;
505     }
506   }
507   return 0;
508 }
509 
510 CEED_QFUNCTION(Advection_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
511   Advection_InOutFlowGeneric(ctx, Q, in, out, 3);
512   return 0;
513 }
514 
515 CEED_QFUNCTION(Advection2d_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
516   Advection_InOutFlowGeneric(ctx, Q, in, out, 2);
517   return 0;
518 }
519