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