xref: /libCEED/examples/fluids/qfunctions/advection.h (revision 93639ffbf2bbfa4b6788194663945beabeb0bf2d)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Advection initial condition and operator for Navier-Stokes example using PETSc
1077841947SLeila Ghaffari 
1177841947SLeila Ghaffari #ifndef advection_h
1277841947SLeila Ghaffari #define advection_h
1377841947SLeila Ghaffari 
14ba6664aeSJames Wright #include <ceed.h>
15c9c2c079SJeremy L Thompson #include <math.h>
1677841947SLeila Ghaffari 
17c44b1c7dSJames Wright #include "advection_types.h"
188f4d89c8SJames Wright #include "newtonian_state.h"
198f4d89c8SJames Wright #include "newtonian_types.h"
20c44b1c7dSJames Wright #include "stabilization_types.h"
218756a6ccSJames Wright #include "utils.h"
228756a6ccSJames Wright 
2377841947SLeila Ghaffari // *****************************************************************************
24*93639ffbSJames Wright // This QFunction sets the initial conditions and the boundary conditions
25*93639ffbSJames Wright //   for two test cases: ROTATION and TRANSLATION
26*93639ffbSJames Wright //
27*93639ffbSJames Wright // -- ROTATION (default)
28*93639ffbSJames Wright //      Initial Conditions:
29*93639ffbSJames Wright //        Mass Density:
30*93639ffbSJames Wright //          Constant mass density of 1.0
31*93639ffbSJames Wright //        Momentum Density:
32*93639ffbSJames Wright //          Rotational field in x,y
33*93639ffbSJames Wright //        Energy Density:
34*93639ffbSJames Wright //          Maximum of 1. x0 decreasing linearly to 0. as radial distance
35*93639ffbSJames Wright //            increases to (1.-r/rc), then 0. everywhere else
36*93639ffbSJames Wright //
37*93639ffbSJames Wright //      Boundary Conditions:
38*93639ffbSJames Wright //        Mass Density:
39*93639ffbSJames Wright //          0.0 flux
40*93639ffbSJames Wright //        Momentum Density:
41*93639ffbSJames Wright //          0.0
42*93639ffbSJames Wright //        Energy Density:
43*93639ffbSJames Wright //          0.0 flux
44*93639ffbSJames Wright //
45*93639ffbSJames Wright // -- TRANSLATION
46*93639ffbSJames Wright //      Initial Conditions:
47*93639ffbSJames Wright //        Mass Density:
48*93639ffbSJames Wright //          Constant mass density of 1.0
49*93639ffbSJames Wright //        Momentum Density:
50*93639ffbSJames Wright //           Constant rectilinear field in x,y
51*93639ffbSJames Wright //        Energy Density:
52*93639ffbSJames Wright //          Maximum of 1. x0 decreasing linearly to 0. as radial distance
53*93639ffbSJames Wright //            increases to (1.-r/rc), then 0. everywhere else
54*93639ffbSJames Wright //
55*93639ffbSJames Wright //      Boundary Conditions:
56*93639ffbSJames Wright //        Mass Density:
57*93639ffbSJames Wright //          0.0 flux
58*93639ffbSJames Wright //        Momentum Density:
59*93639ffbSJames Wright //          0.0
60*93639ffbSJames Wright //        Energy Density:
61*93639ffbSJames Wright //          Inflow BCs:
62*93639ffbSJames Wright //            E = E_wind
63*93639ffbSJames Wright //          Outflow BCs:
64*93639ffbSJames Wright //            E = E(boundary)
65*93639ffbSJames Wright //          Both In/Outflow BCs for E are applied weakly in the
66*93639ffbSJames Wright //            QFunction "Advection2d_Sur"
67*93639ffbSJames Wright //
68*93639ffbSJames Wright // *****************************************************************************
69*93639ffbSJames Wright 
70*93639ffbSJames Wright // *****************************************************************************
71*93639ffbSJames Wright // This helper function provides the exact, time-dependent solution and IC formulation for 2D advection
72*93639ffbSJames Wright // *****************************************************************************
73*93639ffbSJames Wright CEED_QFUNCTION_HELPER CeedInt Exact_AdvectionGeneric(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
74*93639ffbSJames Wright   const SetupContextAdv context = (SetupContextAdv)ctx;
75*93639ffbSJames Wright   const CeedScalar      rc      = context->rc;
76*93639ffbSJames Wright   const CeedScalar      lx      = context->lx;
77*93639ffbSJames Wright   const CeedScalar      ly      = context->ly;
78*93639ffbSJames Wright   const CeedScalar      lz      = dim == 2 ? 0. : context->lz;
79*93639ffbSJames Wright   const CeedScalar     *wind    = context->wind;
80*93639ffbSJames Wright 
81*93639ffbSJames Wright   const CeedScalar center[3] = {0.5 * lx, 0.5 * ly, 0.5 * lz};
82*93639ffbSJames Wright   const CeedScalar theta     = dim == 2 ? M_PI / 3 : M_PI;
83*93639ffbSJames Wright   const CeedScalar x0[3]     = {center[0] + .25 * lx * cos(theta + time), center[1] + .25 * ly * sin(theta + time), 0.5 * lz};
84*93639ffbSJames Wright 
85*93639ffbSJames Wright   const CeedScalar x = X[0], y = X[1], z = dim == 2 ? 0. : X[2];
86*93639ffbSJames Wright 
87*93639ffbSJames Wright   CeedScalar r = 0.;
88*93639ffbSJames Wright   switch (context->initial_condition_type) {
89*93639ffbSJames Wright     case ADVECTIONIC_BUBBLE_SPHERE:
90*93639ffbSJames Wright     case ADVECTIONIC_BUBBLE_CYLINDER:
91*93639ffbSJames Wright       r = sqrt(Square(x - x0[0]) + Square(y - x0[1]) + Square(z - x0[2]));
92*93639ffbSJames Wright       break;
93*93639ffbSJames Wright     case ADVECTIONIC_COSINE_HILL:
94*93639ffbSJames Wright       r = sqrt(Square(x - center[0]) + Square(y - center[1]));
95*93639ffbSJames Wright       break;
96*93639ffbSJames Wright     case ADVECTIONIC_SKEW:
97*93639ffbSJames Wright       break;
98*93639ffbSJames Wright   }
99*93639ffbSJames Wright 
100*93639ffbSJames Wright   switch (context->wind_type) {
101*93639ffbSJames Wright     case WIND_ROTATION:
102*93639ffbSJames Wright       q[0] = 1.;
103*93639ffbSJames Wright       q[1] = -(y - center[1]);
104*93639ffbSJames Wright       q[2] = (x - center[0]);
105*93639ffbSJames Wright       q[3] = 0;
106*93639ffbSJames Wright       break;
107*93639ffbSJames Wright     case WIND_TRANSLATION:
108*93639ffbSJames Wright       q[0] = 1.;
109*93639ffbSJames Wright       q[1] = wind[0];
110*93639ffbSJames Wright       q[2] = wind[1];
111*93639ffbSJames Wright       q[3] = dim == 2 ? 0. : wind[2];
112*93639ffbSJames Wright       break;
113*93639ffbSJames Wright     default:
114*93639ffbSJames Wright       return 1;
115*93639ffbSJames Wright   }
116*93639ffbSJames Wright 
117*93639ffbSJames Wright   switch (context->initial_condition_type) {
118*93639ffbSJames Wright     case ADVECTIONIC_BUBBLE_SPHERE:
119*93639ffbSJames Wright     case ADVECTIONIC_BUBBLE_CYLINDER:
120*93639ffbSJames Wright       switch (context->bubble_continuity_type) {
121*93639ffbSJames Wright         // original continuous, smooth shape
122*93639ffbSJames Wright         case BUBBLE_CONTINUITY_SMOOTH:
123*93639ffbSJames Wright           q[4] = r <= rc ? (1. - r / rc) : 0.;
124*93639ffbSJames Wright           break;
125*93639ffbSJames Wright         // discontinuous, sharp back half shape
126*93639ffbSJames Wright         case BUBBLE_CONTINUITY_BACK_SHARP:
127*93639ffbSJames Wright           q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) : 0.;
128*93639ffbSJames Wright           break;
129*93639ffbSJames Wright         // attempt to define a finite thickness that will get resolved under grid refinement
130*93639ffbSJames Wright         case BUBBLE_CONTINUITY_THICK:
131*93639ffbSJames Wright           q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) * fmin(1.0, (center[1] - y) / 1.25) : 0.;
132*93639ffbSJames Wright           break;
133*93639ffbSJames Wright         case BUBBLE_CONTINUITY_COSINE:
134*93639ffbSJames Wright           q[4] = r <= rc ? .5 + .5 * cos(r * M_PI / rc) : 0;
135*93639ffbSJames Wright           break;
136*93639ffbSJames Wright       }
137*93639ffbSJames Wright       break;
138*93639ffbSJames Wright     case ADVECTIONIC_COSINE_HILL: {
139*93639ffbSJames Wright       CeedScalar half_width = context->lx / 2;
140*93639ffbSJames Wright       q[4]                  = r > half_width ? 0. : cos(2 * M_PI * r / half_width + M_PI) + 1.;
141*93639ffbSJames Wright     } break;
142*93639ffbSJames Wright     case ADVECTIONIC_SKEW: {
143*93639ffbSJames Wright       CeedScalar       skewed_barrier[3]  = {wind[0], wind[1], 0};
144*93639ffbSJames Wright       CeedScalar       inflow_to_point[3] = {x - context->lx / 2, y, 0};
145*93639ffbSJames Wright       CeedScalar       cross_product[3]   = {0};
146*93639ffbSJames Wright       const CeedScalar boundary_threshold = 20 * CEED_EPSILON;
147*93639ffbSJames Wright       Cross3(skewed_barrier, inflow_to_point, cross_product);
148*93639ffbSJames Wright 
149*93639ffbSJames Wright       q[4] = cross_product[2] > boundary_threshold ? 0 : 1;
150*93639ffbSJames Wright       if ((x < boundary_threshold && wind[0] < boundary_threshold) ||                // outflow at -x boundary
151*93639ffbSJames Wright           (y < boundary_threshold && wind[1] < boundary_threshold) ||                // outflow at -y boundary
152*93639ffbSJames Wright           (x > context->lx - boundary_threshold && wind[0] > boundary_threshold) ||  // outflow at +x boundary
153*93639ffbSJames Wright           (y > context->ly - boundary_threshold && wind[1] > boundary_threshold)     // outflow at +y boundary
154*93639ffbSJames Wright       ) {
155*93639ffbSJames Wright         q[4] = 0;
156*93639ffbSJames Wright       }
157*93639ffbSJames Wright     } break;
158*93639ffbSJames Wright   }
159*93639ffbSJames Wright   return 0;
160*93639ffbSJames Wright }
161*93639ffbSJames Wright 
162*93639ffbSJames Wright // *****************************************************************************
16377841947SLeila Ghaffari // This QFunction sets the initial conditions for 3D advection
16477841947SLeila Ghaffari // *****************************************************************************
1652b730f8bSJeremy L Thompson CEED_QFUNCTION(ICsAdvection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
16677841947SLeila Ghaffari   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
16777841947SLeila Ghaffari   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
16877841947SLeila Ghaffari 
16946603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
17077841947SLeila Ghaffari     const CeedScalar x[]  = {X[0][i], X[1][i], X[2][i]};
171e6225c47SLeila Ghaffari     CeedScalar       q[5] = {0.};
17277841947SLeila Ghaffari 
17330e1b2c7SJames Wright     Exact_AdvectionGeneric(3, 0., x, 5, q, ctx);
17477841947SLeila Ghaffari     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
17530e1b2c7SJames Wright   }
17677841947SLeila Ghaffari   return 0;
17777841947SLeila Ghaffari }
17877841947SLeila Ghaffari 
17977841947SLeila Ghaffari // *****************************************************************************
180*93639ffbSJames Wright // This QFunction sets the initial conditions for 2D advection
18177841947SLeila Ghaffari // *****************************************************************************
182*93639ffbSJames Wright CEED_QFUNCTION(ICsAdvection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
183*93639ffbSJames Wright   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
184*93639ffbSJames Wright   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
185*93639ffbSJames Wright   const SetupContextAdv context    = (SetupContextAdv)ctx;
186*93639ffbSJames Wright 
187*93639ffbSJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
188*93639ffbSJames Wright     const CeedScalar x[]  = {X[0][i], X[1][i]};
189*93639ffbSJames Wright     CeedScalar       q[5] = {0.};
190*93639ffbSJames Wright 
191*93639ffbSJames Wright     Exact_AdvectionGeneric(2, context->time, x, 5, q, ctx);
192*93639ffbSJames Wright     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
193*93639ffbSJames Wright   }
19477841947SLeila Ghaffari   return 0;
19577841947SLeila Ghaffari }
19677841947SLeila Ghaffari 
197*93639ffbSJames Wright CEED_QFUNCTION_HELPER void QdataUnpack_ND(CeedInt N, CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar *dXdx) {
198*93639ffbSJames Wright   switch (N) {
199*93639ffbSJames Wright     case 2:
200*93639ffbSJames Wright       QdataUnpack_2D(Q, i, q_data, wdetJ, (CeedScalar(*)[2])dXdx);
201*93639ffbSJames Wright       break;
202*93639ffbSJames Wright     case 3:
203*93639ffbSJames Wright       QdataUnpack_3D(Q, i, q_data, wdetJ, (CeedScalar(*)[3])dXdx);
204*93639ffbSJames Wright       break;
205*93639ffbSJames Wright   }
206*93639ffbSJames Wright }
207*93639ffbSJames Wright 
208*93639ffbSJames Wright CEED_QFUNCTION_HELPER int QdataBoundaryUnpack_ND(CeedInt N, CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar *dXdx,
209*93639ffbSJames Wright                                                  CeedScalar *normal) {
210*93639ffbSJames Wright   switch (N) {
211*93639ffbSJames Wright     case 2:
212*93639ffbSJames Wright       QdataBoundaryUnpack_2D(Q, i, q_data, wdetJ, normal);
213*93639ffbSJames Wright       break;
214*93639ffbSJames Wright     case 3:
215*93639ffbSJames Wright       QdataBoundaryUnpack_3D(Q, i, q_data, wdetJ, (CeedScalar(*)[3])dXdx, normal);
216*93639ffbSJames Wright       break;
217*93639ffbSJames Wright   }
218*93639ffbSJames Wright   return CEED_ERROR_SUCCESS;
219*93639ffbSJames Wright }
220*93639ffbSJames Wright 
221*93639ffbSJames Wright CEED_QFUNCTION_HELPER void StatePhysicalGradientFromReference_ND(CeedInt N, CeedInt Q, CeedInt i, NewtonianIdealGasContext gas, State s,
222*93639ffbSJames Wright                                                                  StateVariable state_var, const CeedScalar *grad_q, const CeedScalar *dXdx,
223*93639ffbSJames Wright                                                                  State *grad_s) {
224*93639ffbSJames Wright   switch (N) {
225*93639ffbSJames Wright     case 2: {
226*93639ffbSJames Wright       for (CeedInt k = 0; k < 2; k++) {
227*93639ffbSJames Wright         CeedScalar dqi[5];
228*93639ffbSJames Wright         for (CeedInt j = 0; j < 5; j++) {
229*93639ffbSJames Wright           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*93639ffbSJames Wright         }
231*93639ffbSJames Wright         grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var);
232*93639ffbSJames Wright       }
233*93639ffbSJames Wright       CeedScalar U[5] = {0.};
234*93639ffbSJames Wright       grad_s[2]       = StateFromU(gas, U);
235*93639ffbSJames Wright     } break;
236*93639ffbSJames Wright     case 3:
237*93639ffbSJames Wright       StatePhysicalGradientFromReference(Q, i, gas, s, state_var, grad_q, (CeedScalar(*)[3])dXdx, grad_s);
238*93639ffbSJames Wright       break;
239*93639ffbSJames Wright   }
240*93639ffbSJames Wright }
241*93639ffbSJames Wright 
242*93639ffbSJames Wright // *****************************************************************************
243*93639ffbSJames Wright // This QFunction implements Advection for implicit time stepping method
244*93639ffbSJames Wright // *****************************************************************************
245*93639ffbSJames Wright CEED_QFUNCTION_HELPER void IFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
246*93639ffbSJames Wright   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
247*93639ffbSJames Wright   const CeedScalar(*grad_q)            = in[1];
248*93639ffbSJames Wright   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
249*93639ffbSJames Wright   const CeedScalar(*q_data)            = in[3];
250*93639ffbSJames Wright 
251*93639ffbSJames Wright   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
252*93639ffbSJames Wright   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
253*93639ffbSJames Wright   CeedScalar *jac_data               = out[2];
254*93639ffbSJames Wright 
255*93639ffbSJames Wright   AdvectionContext                 context   = (AdvectionContext)ctx;
256*93639ffbSJames Wright   const CeedScalar                 CtauS     = context->CtauS;
257*93639ffbSJames Wright   const CeedScalar                 zeros[14] = {0.};
258*93639ffbSJames Wright   NewtonianIdealGasContext         gas;
259*93639ffbSJames Wright   struct NewtonianIdealGasContext_ gas_struct = {0};
260*93639ffbSJames Wright   gas                                         = &gas_struct;
261*93639ffbSJames Wright 
262*93639ffbSJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
263*93639ffbSJames Wright     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
264*93639ffbSJames Wright     const State      s     = StateFromU(gas, qi);
265*93639ffbSJames Wright 
266*93639ffbSJames Wright     CeedScalar wdetJ, dXdx[9];
267*93639ffbSJames Wright     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
268*93639ffbSJames Wright     State grad_s[3];
269*93639ffbSJames Wright     StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
270*93639ffbSJames Wright 
271*93639ffbSJames Wright     const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
272*93639ffbSJames Wright 
273*93639ffbSJames Wright     for (CeedInt f = 0; f < 4; f++) {
274*93639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
275*93639ffbSJames Wright       v[f][i] = wdetJ * q_dot[f][i];                          // K Mass/transient term
276*93639ffbSJames Wright     }
277*93639ffbSJames Wright 
278*93639ffbSJames Wright     CeedScalar div_u = 0;
279*93639ffbSJames Wright     for (CeedInt j = 0; j < dim; j++) {
280*93639ffbSJames Wright       for (CeedInt k = 0; k < dim; k++) {
281*93639ffbSJames Wright         div_u += grad_s[k].Y.velocity[j];
282*93639ffbSJames Wright       }
283*93639ffbSJames Wright     }
284*93639ffbSJames Wright     CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
285*93639ffbSJames Wright     CeedScalar strong_res  = q_dot[4][i] + strong_conv;
286*93639ffbSJames Wright 
287*93639ffbSJames Wright     v[4][i] = wdetJ * q_dot[4][i];  // transient part (ALWAYS)
288*93639ffbSJames Wright 
289*93639ffbSJames Wright     CeedScalar uX[3] = {0.};
290*93639ffbSJames Wright     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
291*93639ffbSJames Wright 
292*93639ffbSJames Wright     if (context->strong_form) {  // Strong Galerkin convection term: v div(E u)
293*93639ffbSJames Wright       v[4][i] += wdetJ * strong_conv;
294*93639ffbSJames Wright     } else {  // Weak Galerkin convection term: -dv \cdot (E u)
295*93639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = -wdetJ * s.U.E_total * uX[j];
296*93639ffbSJames Wright     }
297*93639ffbSJames Wright 
298*93639ffbSJames Wright     const CeedScalar TauS = CtauS / sqrt(Dot3(uX, uX));
299*93639ffbSJames Wright     for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) {
300*93639ffbSJames Wright         case STAB_NONE:
301*93639ffbSJames Wright           break;
302*93639ffbSJames Wright         case STAB_SU:
303*93639ffbSJames Wright           grad_v[j][4][i] += wdetJ * TauS * strong_conv * uX[j];
304*93639ffbSJames Wright           break;
305*93639ffbSJames Wright         case STAB_SUPG:
306*93639ffbSJames Wright           grad_v[j][4][i] += wdetJ * TauS * strong_res * uX[j];
307*93639ffbSJames Wright           break;
308*93639ffbSJames Wright       }
309*93639ffbSJames Wright     StoredValuesPack(Q, i, 0, 14, zeros, jac_data);
310*93639ffbSJames Wright   }
311*93639ffbSJames Wright }
312*93639ffbSJames Wright 
3132b730f8bSJeremy L Thompson CEED_QFUNCTION(IFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
314372d1924SJames Wright   IFunction_AdvectionGeneric(ctx, Q, in, out, 3);
31577841947SLeila Ghaffari   return 0;
31677841947SLeila Ghaffari }
31777841947SLeila Ghaffari 
318*93639ffbSJames Wright CEED_QFUNCTION(IFunction_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
319*93639ffbSJames Wright   IFunction_AdvectionGeneric(ctx, Q, in, out, 2);
320*93639ffbSJames Wright   return 0;
321*93639ffbSJames Wright }
322*93639ffbSJames Wright 
323*93639ffbSJames Wright // *****************************************************************************
324*93639ffbSJames Wright // This QFunction implements Advection for explicit time stepping method
325*93639ffbSJames Wright // *****************************************************************************
326*93639ffbSJames Wright CEED_QFUNCTION_HELPER void RHSFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
327*93639ffbSJames Wright   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
328*93639ffbSJames Wright   const CeedScalar(*grad_q)        = in[1];
329*93639ffbSJames Wright   const CeedScalar(*q_data)        = in[2];
330*93639ffbSJames Wright 
331*93639ffbSJames Wright   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
332*93639ffbSJames Wright   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
333*93639ffbSJames Wright 
334*93639ffbSJames Wright   AdvectionContext                 context = (AdvectionContext)ctx;
335*93639ffbSJames Wright   const CeedScalar                 CtauS   = context->CtauS;
336*93639ffbSJames Wright   NewtonianIdealGasContext         gas;
337*93639ffbSJames Wright   struct NewtonianIdealGasContext_ gas_struct = {0};
338*93639ffbSJames Wright   gas                                         = &gas_struct;
339*93639ffbSJames Wright 
340*93639ffbSJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
341*93639ffbSJames Wright     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
342*93639ffbSJames Wright     const State      s     = StateFromU(gas, qi);
343*93639ffbSJames Wright 
344*93639ffbSJames Wright     CeedScalar wdetJ, dXdx[9];
345*93639ffbSJames Wright     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
346*93639ffbSJames Wright     State grad_s[3];
347*93639ffbSJames Wright     StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
348*93639ffbSJames Wright 
349*93639ffbSJames Wright     const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
350*93639ffbSJames Wright 
351*93639ffbSJames Wright     for (CeedInt f = 0; f < 4; f++) {
352*93639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
353*93639ffbSJames Wright       v[f][i] = 0.;
354*93639ffbSJames Wright     }
355*93639ffbSJames Wright 
356*93639ffbSJames Wright     CeedScalar div_u = 0;
357*93639ffbSJames Wright     for (CeedInt j = 0; j < dim; j++) {
358*93639ffbSJames Wright       for (CeedInt k = 0; k < dim; k++) {
359*93639ffbSJames Wright         div_u += grad_s[k].Y.velocity[j];
360*93639ffbSJames Wright       }
361*93639ffbSJames Wright     }
362*93639ffbSJames Wright     CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
363*93639ffbSJames Wright 
364*93639ffbSJames Wright     CeedScalar uX[3] = {0.};
365*93639ffbSJames Wright     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
366*93639ffbSJames Wright 
367*93639ffbSJames Wright     if (context->strong_form) {  // Strong Galerkin convection term: v div(E u)
368*93639ffbSJames Wright       v[4][i] = -wdetJ * strong_conv;
369*93639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = 0;
370*93639ffbSJames Wright     } else {  // Weak Galerkin convection term: -dv \cdot (E u)
371*93639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = wdetJ * s.U.E_total * uX[j];
372*93639ffbSJames Wright       v[4][i] = 0.;
373*93639ffbSJames Wright     }
374*93639ffbSJames Wright 
375*93639ffbSJames Wright     const CeedScalar TauS = CtauS / sqrt(Dot3(uX, uX));
376*93639ffbSJames Wright     for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) {
377*93639ffbSJames Wright         case STAB_NONE:
378*93639ffbSJames Wright           break;
379*93639ffbSJames Wright         case STAB_SU:
380*93639ffbSJames Wright         case STAB_SUPG:
381*93639ffbSJames Wright           grad_v[j][4][i] += wdetJ * TauS * strong_conv * uX[j];
382*93639ffbSJames Wright           break;
383*93639ffbSJames Wright       }
384*93639ffbSJames Wright   }
385*93639ffbSJames Wright }
386*93639ffbSJames Wright 
387*93639ffbSJames Wright CEED_QFUNCTION(RHS_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
388*93639ffbSJames Wright   RHSFunction_AdvectionGeneric(ctx, Q, in, out, 3);
389*93639ffbSJames Wright   return 0;
390*93639ffbSJames Wright }
391*93639ffbSJames Wright 
392*93639ffbSJames Wright CEED_QFUNCTION(RHS_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
393*93639ffbSJames Wright   RHSFunction_AdvectionGeneric(ctx, Q, in, out, 2);
394*93639ffbSJames Wright   return 0;
395*93639ffbSJames Wright }
396*93639ffbSJames Wright 
397*93639ffbSJames Wright // *****************************************************************************
398*93639ffbSJames Wright // This QFunction implements consistent outflow and inflow BCs
399*93639ffbSJames Wright //      for advection
400*93639ffbSJames Wright //
401*93639ffbSJames Wright //  Inflow and outflow faces are determined based on sign(dot(wind, normal)):
402*93639ffbSJames Wright //    sign(dot(wind, normal)) > 0 : outflow BCs
403*93639ffbSJames Wright //    sign(dot(wind, normal)) < 0 : inflow BCs
404*93639ffbSJames Wright //
405*93639ffbSJames Wright //  Outflow BCs:
406*93639ffbSJames Wright //    The validity of the weak form of the governing equations is extended to the outflow and the current values of E are applied.
407*93639ffbSJames Wright //
408*93639ffbSJames Wright //  Inflow BCs:
409*93639ffbSJames Wright //    A prescribed Total Energy (E_wind) is applied weakly.
410*93639ffbSJames Wright // *****************************************************************************
411*93639ffbSJames Wright CEED_QFUNCTION(Advection_InOutFlowGeneric)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
412*93639ffbSJames Wright   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
413*93639ffbSJames Wright   const CeedScalar(*q_data_sur)    = in[2];
414*93639ffbSJames Wright 
415*93639ffbSJames Wright   CeedScalar(*v)[CEED_Q_VLA]   = (CeedScalar(*)[CEED_Q_VLA])out[0];
416*93639ffbSJames Wright   AdvectionContext context     = (AdvectionContext)ctx;
417*93639ffbSJames Wright   const CeedScalar E_wind      = context->E_wind;
418*93639ffbSJames Wright   const CeedScalar strong_form = context->strong_form;
419*93639ffbSJames Wright   const bool       is_implicit = context->implicit;
420*93639ffbSJames Wright 
421*93639ffbSJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
422*93639ffbSJames Wright     const CeedScalar rho  = q[0][i];
423*93639ffbSJames Wright     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
424*93639ffbSJames Wright     const CeedScalar E    = q[4][i];
425*93639ffbSJames Wright 
426*93639ffbSJames Wright     CeedScalar wdetJb, norm[3];
427*93639ffbSJames Wright     QdataBoundaryUnpack_ND(dim, Q, i, q_data_sur, &wdetJb, NULL, norm);
428*93639ffbSJames Wright     wdetJb *= is_implicit ? -1. : 1.;
429*93639ffbSJames Wright 
430*93639ffbSJames Wright     const CeedScalar u_normal = DotN(norm, u, dim);
431*93639ffbSJames Wright 
432*93639ffbSJames Wright     // No Change in density or momentum
433*93639ffbSJames Wright     for (CeedInt j = 0; j < 4; j++) {
434*93639ffbSJames Wright       v[j][i] = 0;
435*93639ffbSJames Wright     }
436*93639ffbSJames Wright     // Implementing in/outflow BCs
437*93639ffbSJames Wright     if (u_normal > 0) {  // outflow
438*93639ffbSJames Wright       v[4][i] = -(1 - strong_form) * wdetJb * E * u_normal;
439*93639ffbSJames Wright     } else {  // inflow
440*93639ffbSJames Wright       v[4][i] = -(1 - strong_form) * wdetJb * E_wind * u_normal;
441*93639ffbSJames Wright     }
442*93639ffbSJames Wright   }
443*93639ffbSJames Wright   return 0;
444*93639ffbSJames Wright }
445*93639ffbSJames Wright 
4462b730f8bSJeremy L Thompson CEED_QFUNCTION(Advection_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
4474bdcf5bfSJames Wright   Advection_InOutFlowGeneric(ctx, Q, in, out, 3);
44877841947SLeila Ghaffari   return 0;
44977841947SLeila Ghaffari }
45077841947SLeila Ghaffari 
451*93639ffbSJames Wright CEED_QFUNCTION(Advection2d_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
452*93639ffbSJames Wright   Advection_InOutFlowGeneric(ctx, Q, in, out, 2);
453*93639ffbSJames Wright   return 0;
454*93639ffbSJames Wright }
455*93639ffbSJames Wright 
45677841947SLeila Ghaffari #endif  // advection_h
457