xref: /honee/qfunctions/advection.h (revision b4fd18dfeb7fe20bc2ce09e18a422e556e44809a)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari 
4a515125bSLeila Ghaffari /// @file
5a515125bSLeila Ghaffari /// Advection initial condition and operator for Navier-Stokes example using PETSc
6493642f1SJames Wright #include <ceed.h>
7d0cce58aSJeremy L Thompson #include <math.h>
8a515125bSLeila Ghaffari 
9e88b842aSJames Wright #include "advection_types.h"
10ce192147SJames Wright #include "newtonian_state.h"
11ce192147SJames Wright #include "newtonian_types.h"
12e88b842aSJames Wright #include "stabilization_types.h"
131a74fa30SJames Wright #include "utils.h"
141a74fa30SJames Wright 
15a515125bSLeila Ghaffari // *****************************************************************************
169529d636SJames Wright // This QFunction sets the initial conditions and the boundary conditions
179529d636SJames Wright //   for two test cases: ROTATION and TRANSLATION
189529d636SJames Wright //
199529d636SJames Wright // -- ROTATION (default)
209529d636SJames Wright //      Initial Conditions:
219529d636SJames Wright //        Mass Density:
229529d636SJames Wright //          Constant mass density of 1.0
239529d636SJames Wright //        Momentum Density:
249529d636SJames Wright //          Rotational field in x,y
259529d636SJames Wright //        Energy Density:
269529d636SJames Wright //          Maximum of 1. x0 decreasing linearly to 0. as radial distance
279529d636SJames Wright //            increases to (1.-r/rc), then 0. everywhere else
289529d636SJames Wright //
299529d636SJames Wright //      Boundary Conditions:
309529d636SJames Wright //        Mass Density:
319529d636SJames Wright //          0.0 flux
329529d636SJames Wright //        Momentum Density:
339529d636SJames Wright //          0.0
349529d636SJames Wright //        Energy Density:
359529d636SJames Wright //          0.0 flux
369529d636SJames Wright //
379529d636SJames Wright // -- TRANSLATION
389529d636SJames Wright //      Initial Conditions:
399529d636SJames Wright //        Mass Density:
409529d636SJames Wright //          Constant mass density of 1.0
419529d636SJames Wright //        Momentum Density:
429529d636SJames Wright //           Constant rectilinear field in x,y
439529d636SJames Wright //        Energy Density:
449529d636SJames Wright //          Maximum of 1. x0 decreasing linearly to 0. as radial distance
459529d636SJames Wright //            increases to (1.-r/rc), then 0. everywhere else
469529d636SJames Wright //
479529d636SJames Wright //      Boundary Conditions:
489529d636SJames Wright //        Mass Density:
499529d636SJames Wright //          0.0 flux
509529d636SJames Wright //        Momentum Density:
519529d636SJames Wright //          0.0
529529d636SJames Wright //        Energy Density:
539529d636SJames Wright //          Inflow BCs:
549529d636SJames Wright //            E = E_wind
559529d636SJames Wright //          Outflow BCs:
569529d636SJames Wright //            E = E(boundary)
579529d636SJames Wright //          Both In/Outflow BCs for E are applied weakly in the
589529d636SJames Wright //            QFunction "Advection2d_Sur"
599529d636SJames Wright //
609529d636SJames Wright // *****************************************************************************
619529d636SJames Wright 
629529d636SJames Wright // *****************************************************************************
639529d636SJames Wright // This helper function provides the exact, time-dependent solution and IC formulation for 2D advection
649529d636SJames Wright // *****************************************************************************
6597cfd714SJames Wright CEED_QFUNCTION_HELPER int Exact_AdvectionGeneric(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
669529d636SJames Wright   const SetupContextAdv context = (SetupContextAdv)ctx;
679529d636SJames Wright   const CeedScalar      rc      = context->rc;
689529d636SJames Wright   const CeedScalar      lx      = context->lx;
699529d636SJames Wright   const CeedScalar      ly      = context->ly;
709529d636SJames Wright   const CeedScalar      lz      = dim == 2 ? 0. : context->lz;
719529d636SJames Wright   const CeedScalar     *wind    = context->wind;
729529d636SJames Wright 
739529d636SJames Wright   const CeedScalar center[3] = {0.5 * lx, 0.5 * ly, 0.5 * lz};
749529d636SJames Wright   const CeedScalar theta     = dim == 2 ? M_PI / 3 : M_PI;
759529d636SJames Wright   const CeedScalar x0[3]     = {center[0] + .25 * lx * cos(theta + time), center[1] + .25 * ly * sin(theta + time), 0.5 * lz};
769529d636SJames Wright 
779529d636SJames Wright   const CeedScalar x = X[0], y = X[1], z = dim == 2 ? 0. : X[2];
789529d636SJames Wright 
799529d636SJames Wright   switch (context->wind_type) {
805f636aeaSJames Wright     case ADVDIF_WIND_ROTATION:
819529d636SJames Wright       q[0] = 1.;
829529d636SJames Wright       q[1] = -(y - center[1]);
839529d636SJames Wright       q[2] = (x - center[0]);
849529d636SJames Wright       q[3] = 0;
859529d636SJames Wright       break;
865f636aeaSJames Wright     case ADVDIF_WIND_TRANSLATION:
879529d636SJames Wright       q[0] = 1.;
889529d636SJames Wright       q[1] = wind[0];
899529d636SJames Wright       q[2] = wind[1];
909529d636SJames Wright       q[3] = dim == 2 ? 0. : wind[2];
919529d636SJames Wright       break;
923d1afcc1SJames Wright     case ADVDIF_WIND_BOUNDARY_LAYER:
933d1afcc1SJames Wright       q[0] = 1.;
943d1afcc1SJames Wright       q[1] = y / ly;
953d1afcc1SJames Wright       q[2] = 0.;
963d1afcc1SJames Wright       q[3] = 0.;
973d1afcc1SJames Wright       break;
989529d636SJames Wright   }
999529d636SJames Wright 
1009529d636SJames Wright   switch (context->initial_condition_type) {
1015f636aeaSJames Wright     case ADVDIF_IC_BUBBLE_SPHERE:
1025f636aeaSJames Wright     case ADVDIF_IC_BUBBLE_CYLINDER: {
103a62be6baSJames Wright       CeedScalar r = sqrt(Square(x - x0[0]) + Square(y - x0[1]) + Square(z - x0[2]));
104a62be6baSJames Wright 
1059529d636SJames Wright       switch (context->bubble_continuity_type) {
1069529d636SJames Wright         // original continuous, smooth shape
1075f636aeaSJames Wright         case ADVDIF_BUBBLE_CONTINUITY_SMOOTH:
1089529d636SJames Wright           q[4] = r <= rc ? (1. - r / rc) : 0.;
1099529d636SJames Wright           break;
1109529d636SJames Wright         // discontinuous, sharp back half shape
1115f636aeaSJames Wright         case ADVDIF_BUBBLE_CONTINUITY_BACK_SHARP:
1129529d636SJames Wright           q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) : 0.;
1139529d636SJames Wright           break;
1149529d636SJames Wright         // attempt to define a finite thickness that will get resolved under grid refinement
1155f636aeaSJames Wright         case ADVDIF_BUBBLE_CONTINUITY_THICK:
1169529d636SJames Wright           q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) * fmin(1.0, (center[1] - y) / 1.25) : 0.;
1179529d636SJames Wright           break;
1185f636aeaSJames Wright         case ADVDIF_BUBBLE_CONTINUITY_COSINE:
1199529d636SJames Wright           q[4] = r <= rc ? .5 + .5 * cos(r * M_PI / rc) : 0;
1209529d636SJames Wright           break;
1219529d636SJames Wright       }
1229529d636SJames Wright       break;
123a62be6baSJames Wright     }
124a62be6baSJames Wright 
1255f636aeaSJames Wright     case ADVDIF_IC_COSINE_HILL: {
126a62be6baSJames Wright       CeedScalar r          = sqrt(Square(x - center[0]) + Square(y - center[1]));
1279529d636SJames Wright       CeedScalar half_width = context->lx / 2;
1289529d636SJames Wright       q[4]                  = r > half_width ? 0. : cos(2 * M_PI * r / half_width + M_PI) + 1.;
1299529d636SJames Wright     } break;
130a62be6baSJames Wright 
1315f636aeaSJames Wright     case ADVDIF_IC_SKEW: {
1329529d636SJames Wright       CeedScalar       skewed_barrier[3]  = {wind[0], wind[1], 0};
1339529d636SJames Wright       CeedScalar       inflow_to_point[3] = {x - context->lx / 2, y, 0};
1349529d636SJames Wright       CeedScalar       cross_product[3]   = {0};
1359529d636SJames Wright       const CeedScalar boundary_threshold = 20 * CEED_EPSILON;
1369529d636SJames Wright       Cross3(skewed_barrier, inflow_to_point, cross_product);
1379529d636SJames Wright 
1389529d636SJames Wright       q[4] = cross_product[2] > boundary_threshold ? 0 : 1;
1399529d636SJames Wright       if ((x < boundary_threshold && wind[0] < boundary_threshold) ||                // outflow at -x boundary
1409529d636SJames Wright           (y < boundary_threshold && wind[1] < boundary_threshold) ||                // outflow at -y boundary
1419529d636SJames Wright           (x > context->lx - boundary_threshold && wind[0] > boundary_threshold) ||  // outflow at +x boundary
1429529d636SJames Wright           (y > context->ly - boundary_threshold && wind[1] > boundary_threshold)     // outflow at +y boundary
1439529d636SJames Wright       ) {
1449529d636SJames Wright         q[4] = 0;
1459529d636SJames Wright       }
1469529d636SJames Wright     } break;
147a62be6baSJames Wright 
1485f636aeaSJames Wright     case ADVDIF_IC_WAVE: {
149a62be6baSJames Wright       CeedScalar theta = context->wave_frequency * DotN(X, wind, dim) + context->wave_phase;
150a62be6baSJames Wright       switch (context->wave_type) {
151a62be6baSJames Wright         case ADVDIF_WAVE_SINE:
152a62be6baSJames Wright           q[4] = sin(theta);
153a62be6baSJames Wright           break;
154a62be6baSJames Wright         case ADVDIF_WAVE_SQUARE:
155a62be6baSJames Wright           q[4] = sin(theta) > 100 * CEED_EPSILON ? 1 : -1;
156a62be6baSJames Wright           break;
157a62be6baSJames Wright       }
1583d1afcc1SJames Wright     } break;
1593d1afcc1SJames Wright     case ADVDIF_IC_BOUNDARY_LAYER: {
1603d1afcc1SJames Wright       const CeedScalar boundary_threshold = 20 * CEED_EPSILON;
1613d1afcc1SJames Wright 
1623d1afcc1SJames Wright       if ((x < boundary_threshold) || (y > ly - boundary_threshold)) {
1633d1afcc1SJames Wright         q[4] = 1;  // inflow and top boundary
1643d1afcc1SJames Wright       } else if (y < boundary_threshold) {
1653d1afcc1SJames Wright         q[4] = 0;  // lower wall
166*b4fd18dfSJames Wright       } else {     // interior and outflow boundary
167*b4fd18dfSJames Wright         CeedScalar bl_height = ly * context->bl_height_factor;
168*b4fd18dfSJames Wright         if (y > bl_height) q[4] = 1;
169*b4fd18dfSJames Wright         else q[4] = y / bl_height;
170a62be6baSJames Wright       }
1713d1afcc1SJames Wright     } break;
1729529d636SJames Wright   }
1739529d636SJames Wright   return 0;
1749529d636SJames Wright }
1759529d636SJames Wright 
1769529d636SJames Wright // *****************************************************************************
177a515125bSLeila Ghaffari // This QFunction sets the initial conditions for 3D advection
178a515125bSLeila Ghaffari // *****************************************************************************
1792b916ea7SJeremy L Thompson CEED_QFUNCTION(ICsAdvection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
180a515125bSLeila Ghaffari   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
181a515125bSLeila Ghaffari   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
182a515125bSLeila Ghaffari 
1833d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
184a515125bSLeila Ghaffari     const CeedScalar x[]  = {X[0][i], X[1][i], X[2][i]};
185139613f2SLeila Ghaffari     CeedScalar       q[5] = {0.};
186a515125bSLeila Ghaffari 
1870b3a1fabSJames Wright     Exact_AdvectionGeneric(3, 0., x, 5, q, ctx);
188a515125bSLeila Ghaffari     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
1890b3a1fabSJames Wright   }
190a515125bSLeila Ghaffari   return 0;
191a515125bSLeila Ghaffari }
192a515125bSLeila Ghaffari 
193a515125bSLeila Ghaffari // *****************************************************************************
1949529d636SJames Wright // This QFunction sets the initial conditions for 2D advection
195a515125bSLeila Ghaffari // *****************************************************************************
1969529d636SJames Wright CEED_QFUNCTION(ICsAdvection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
1979529d636SJames Wright   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
1989529d636SJames Wright   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
1999529d636SJames Wright   const SetupContextAdv context    = (SetupContextAdv)ctx;
2009529d636SJames Wright 
2019529d636SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
2029529d636SJames Wright     const CeedScalar x[]  = {X[0][i], X[1][i]};
2039529d636SJames Wright     CeedScalar       q[5] = {0.};
2049529d636SJames Wright 
2059529d636SJames Wright     Exact_AdvectionGeneric(2, context->time, x, 5, q, ctx);
2069529d636SJames Wright     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
2079529d636SJames Wright   }
208a515125bSLeila Ghaffari   return 0;
209a515125bSLeila Ghaffari }
210a515125bSLeila Ghaffari 
2119529d636SJames Wright CEED_QFUNCTION_HELPER void StatePhysicalGradientFromReference_ND(CeedInt N, CeedInt Q, CeedInt i, NewtonianIdealGasContext gas, State s,
2129529d636SJames Wright                                                                  StateVariable state_var, const CeedScalar *grad_q, const CeedScalar *dXdx,
2139529d636SJames Wright                                                                  State *grad_s) {
2149529d636SJames Wright   switch (N) {
2159529d636SJames Wright     case 2: {
2169529d636SJames Wright       for (CeedInt k = 0; k < 2; k++) {
2179529d636SJames Wright         CeedScalar dqi[5];
2189529d636SJames Wright         for (CeedInt j = 0; j < 5; j++) {
2199529d636SJames 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];
2209529d636SJames Wright         }
2219529d636SJames Wright         grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var);
2229529d636SJames Wright       }
2239529d636SJames Wright       CeedScalar U[5] = {0.};
2249529d636SJames Wright       grad_s[2]       = StateFromU(gas, U);
2259529d636SJames Wright     } break;
2269529d636SJames Wright     case 3:
22785efd435SJames Wright       // Cannot directly use StatePhysicalGradientFromReference helper functions due to SYCL online compiler incompatabilities
22885efd435SJames Wright       for (CeedInt k = 0; k < 3; k++) {
22985efd435SJames Wright         CeedScalar dqi[5];
23085efd435SJames Wright         for (CeedInt j = 0; j < 5; j++) {
23185efd435SJames 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] +
23285efd435SJames Wright                    grad_q[(Q * 5) * 2 + Q * j + i] * dXdx[2 * N + k];
23385efd435SJames Wright         }
23485efd435SJames Wright         grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var);
23585efd435SJames Wright       }
2369529d636SJames Wright       break;
2379529d636SJames Wright   }
2389529d636SJames Wright }
2399529d636SJames Wright 
240a78efa86SJames Wright // @brief Calculate the stabilization constant \tau
241a78efa86SJames Wright CEED_QFUNCTION_HELPER CeedScalar Tau(AdvectionContext context, const State s, const CeedScalar *dXdx, CeedInt dim) {
242a78efa86SJames Wright   switch (context->stabilization_tau) {
243a78efa86SJames Wright     case STAB_TAU_CTAU: {
244a78efa86SJames Wright       CeedScalar uX[3] = {0.};
245a78efa86SJames Wright 
246a78efa86SJames Wright       MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
247a78efa86SJames Wright       return context->CtauS / sqrt(DotN(uX, uX, dim));
248a78efa86SJames Wright     } break;
249a78efa86SJames Wright     case STAB_TAU_ADVDIFF_SHAKIB: {
250a78efa86SJames Wright       CeedScalar gijd_mat[9] = {0.}, gij_uj[3] = {0.};
251a78efa86SJames Wright 
252a78efa86SJames Wright       MatMatN(dXdx, dXdx, dim, CEED_TRANSPOSE, CEED_NOTRANSPOSE, gijd_mat);
2534ca5135bSJames Wright       // (1/2)^2 to account for reference element size; for length 1 square/cube element, gij should be identity matrix
2544ca5135bSJames Wright       ScaleN(gijd_mat, 0.25, Square(dim));
255a78efa86SJames Wright       MatVecNM(gijd_mat, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, gij_uj);
2564ca5135bSJames Wright       return 1 / sqrt(Square(2 * context->Ctau_t / context->dt) + DotN(s.Y.velocity, gij_uj, dim) * Square(context->Ctau_a) +
2574ca5135bSJames Wright                       Square(context->diffusion_coeff) * DotN(gijd_mat, gijd_mat, dim * dim) * Square(context->Ctau_d));
258a78efa86SJames Wright     } break;
259a78efa86SJames Wright     default:
260a78efa86SJames Wright       return 0.;
261a78efa86SJames Wright   }
262a78efa86SJames Wright }
263a78efa86SJames Wright 
2649529d636SJames Wright // *****************************************************************************
2659529d636SJames Wright // This QFunction implements Advection for implicit time stepping method
2669529d636SJames Wright // *****************************************************************************
26797cfd714SJames Wright CEED_QFUNCTION_HELPER int IFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
2684c5ab12fSJames Wright   AdvectionContext context = (AdvectionContext)ctx;
2694c5ab12fSJames Wright 
2709529d636SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
2719529d636SJames Wright   const CeedScalar(*grad_q)            = in[1];
2729529d636SJames Wright   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
2739529d636SJames Wright   const CeedScalar(*q_data)            = in[3];
2744c5ab12fSJames Wright   const CeedScalar(*divFdiff)          = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? in[5] : NULL;
2759529d636SJames Wright 
2769529d636SJames Wright   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
2779529d636SJames Wright   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
2789529d636SJames Wright 
2799529d636SJames Wright   NewtonianIdealGasContext         gas;
2809529d636SJames Wright   struct NewtonianIdealGasContext_ gas_struct = {0};
2819529d636SJames Wright   gas                                         = &gas_struct;
2829529d636SJames Wright 
2839529d636SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
2849529d636SJames Wright     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
2859529d636SJames Wright     const State      s     = StateFromU(gas, qi);
2869529d636SJames Wright 
2879529d636SJames Wright     CeedScalar wdetJ, dXdx[9];
2889529d636SJames Wright     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
2899529d636SJames Wright     State grad_s[3];
2909529d636SJames Wright     StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
2919529d636SJames Wright 
2929529d636SJames Wright     const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
2939529d636SJames Wright 
2949529d636SJames Wright     for (CeedInt f = 0; f < 4; f++) {
2959529d636SJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
2969529d636SJames Wright       v[f][i] = wdetJ * q_dot[f][i];                          // K Mass/transient term
2979529d636SJames Wright     }
2989529d636SJames Wright 
2999529d636SJames Wright     CeedScalar div_u = 0;
3009529d636SJames Wright     for (CeedInt j = 0; j < dim; j++) {
3019529d636SJames Wright       for (CeedInt k = 0; k < dim; k++) {
3029529d636SJames Wright         div_u += grad_s[k].Y.velocity[j];
3039529d636SJames Wright       }
3049529d636SJames Wright     }
3059529d636SJames Wright     CeedScalar uX[3] = {0.};
3069529d636SJames Wright     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
3074c5ab12fSJames Wright     CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
3089529d636SJames Wright 
3094c5ab12fSJames Wright     v[4][i] = wdetJ * q_dot[4][i];  // transient part (ALWAYS)
3109529d636SJames Wright     if (context->strong_form) {     // Strong Galerkin convection term: v div(E u)
3119529d636SJames Wright       v[4][i] += wdetJ * strong_conv;
3129529d636SJames Wright     } else {  // Weak Galerkin convection term: -dv \cdot (E u)
3139529d636SJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = -wdetJ * s.U.E_total * uX[j];
3149529d636SJames Wright     }
3159529d636SJames Wright 
316c8d249deSJames Wright     {  // Diffusion
317c8d249deSJames Wright       CeedScalar Fe[3], Fe_dXdx[3] = {0.};
318c8d249deSJames Wright 
319c8d249deSJames Wright       for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total;
320c8d249deSJames Wright       MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx);
321c8d249deSJames Wright       for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] -= wdetJ * Fe_dXdx[k];
322c8d249deSJames Wright     }
323c8d249deSJames Wright 
324a78efa86SJames Wright     const CeedScalar TauS = Tau(context, s, dXdx, dim);
3254c5ab12fSJames Wright     for (CeedInt j = 0; j < dim; j++) {
3264c5ab12fSJames Wright       switch (context->stabilization) {
3279529d636SJames Wright         case STAB_NONE:
3289529d636SJames Wright           break;
3299529d636SJames Wright         case STAB_SU:
3304c5ab12fSJames Wright           grad_v[j][4][i] += wdetJ * TauS * uX[j] * strong_conv;
3319529d636SJames Wright           break;
3324c5ab12fSJames Wright         case STAB_SUPG: {
3334c5ab12fSJames Wright           CeedScalar divFdiff_i = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? divFdiff[i] : 0.;
3344c5ab12fSJames Wright           grad_v[j][4][i] += wdetJ * TauS * uX[j] * (q_dot[4][i] + strong_conv + divFdiff_i);
3354c5ab12fSJames Wright         } break;
3364c5ab12fSJames Wright       }
3379529d636SJames Wright     }
3389529d636SJames Wright   }
33997cfd714SJames Wright   return 0;
3409529d636SJames Wright }
3419529d636SJames Wright 
3422b916ea7SJeremy L Thompson CEED_QFUNCTION(IFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
34397cfd714SJames Wright   return IFunction_AdvectionGeneric(ctx, Q, in, out, 3);
344a515125bSLeila Ghaffari }
345a515125bSLeila Ghaffari 
3469529d636SJames Wright CEED_QFUNCTION(IFunction_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
34797cfd714SJames Wright   return IFunction_AdvectionGeneric(ctx, Q, in, out, 2);
3489529d636SJames Wright }
3499529d636SJames Wright 
35097cfd714SJames Wright CEED_QFUNCTION_HELPER int MassFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
351a78efa86SJames Wright   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
352a78efa86SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[1];
353a78efa86SJames Wright   const CeedScalar(*q_data)            = in[2];
354a78efa86SJames Wright 
355a78efa86SJames Wright   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
356a78efa86SJames Wright   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
357a78efa86SJames Wright 
358a78efa86SJames Wright   AdvectionContext                 context    = (AdvectionContext)ctx;
359a78efa86SJames Wright   struct NewtonianIdealGasContext_ gas_struct = {0};
360a78efa86SJames Wright   NewtonianIdealGasContext         gas        = &gas_struct;
361a78efa86SJames Wright 
362a78efa86SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
363a78efa86SJames Wright     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
364a78efa86SJames Wright     const State      s     = StateFromU(gas, qi);
365a78efa86SJames Wright     CeedScalar       wdetJ, dXdx[9];
366a78efa86SJames Wright     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
367a78efa86SJames Wright 
368a78efa86SJames Wright     for (CeedInt f = 0; f < 4; f++) {
369a78efa86SJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
370a78efa86SJames Wright       v[f][i] = wdetJ * q_dot[f][i];                          // K Mass/transient term
371a78efa86SJames Wright     }
372a78efa86SJames Wright 
373a78efa86SJames Wright     // Unstabilized mass term
374a78efa86SJames Wright     v[4][i] = wdetJ * q_dot[4][i];
375a78efa86SJames Wright 
376a78efa86SJames Wright     // Stabilized mass term
377a78efa86SJames Wright     CeedScalar uX[3] = {0.};
378a78efa86SJames Wright     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
379a78efa86SJames Wright     const CeedScalar TauS = Tau(context, s, dXdx, dim);
38071acc5eeSJames Wright     for (CeedInt j = 0; j < dim; j++) {
38171acc5eeSJames Wright       switch (context->stabilization) {
382a78efa86SJames Wright         case STAB_NONE:
383a78efa86SJames Wright         case STAB_SU:
384a78efa86SJames Wright           grad_v[j][4][i] = 0;
385a78efa86SJames Wright           break;  // These should be run with the unstabilized mass matrix anyways
386a78efa86SJames Wright         case STAB_SUPG:
387a78efa86SJames Wright           grad_v[j][4][i] = wdetJ * TauS * q_dot[4][i] * uX[j];
388a78efa86SJames Wright           break;
389a78efa86SJames Wright       }
390a78efa86SJames Wright     }
39171acc5eeSJames Wright   }
39297cfd714SJames Wright   return 0;
393a78efa86SJames Wright }
394a78efa86SJames Wright 
395a78efa86SJames Wright CEED_QFUNCTION(MassFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
39697cfd714SJames Wright   return MassFunction_AdvectionGeneric(ctx, Q, in, out, 3);
397a78efa86SJames Wright }
398a78efa86SJames Wright 
399a78efa86SJames Wright CEED_QFUNCTION(MassFunction_Advection2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
40097cfd714SJames Wright   return MassFunction_AdvectionGeneric(ctx, Q, in, out, 2);
401a78efa86SJames Wright }
402a78efa86SJames Wright 
4039529d636SJames Wright // *****************************************************************************
4049529d636SJames Wright // This QFunction implements Advection for explicit time stepping method
4059529d636SJames Wright // *****************************************************************************
40697cfd714SJames Wright CEED_QFUNCTION_HELPER int RHSFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
4075f952e8dSJames Wright   AdvectionContext context = (AdvectionContext)ctx;
4085f952e8dSJames Wright 
4099529d636SJames Wright   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
4109529d636SJames Wright   const CeedScalar(*grad_q)        = in[1];
4119529d636SJames Wright   const CeedScalar(*q_data)        = in[2];
4125f952e8dSJames Wright   const CeedScalar(*divFdiff)      = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? in[4] : NULL;
4139529d636SJames Wright 
4149529d636SJames Wright   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
4159529d636SJames Wright   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
4169529d636SJames Wright 
4179529d636SJames Wright   struct NewtonianIdealGasContext_ gas_struct = {0};
418a78efa86SJames Wright   NewtonianIdealGasContext         gas        = &gas_struct;
4199529d636SJames Wright 
4209529d636SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
4219529d636SJames Wright     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
4229529d636SJames Wright     const State      s     = StateFromU(gas, qi);
4239529d636SJames Wright 
4249529d636SJames Wright     CeedScalar wdetJ, dXdx[9];
4259529d636SJames Wright     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
4269529d636SJames Wright     State grad_s[3];
4279529d636SJames Wright     StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
4289529d636SJames Wright 
4299529d636SJames Wright     const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
4309529d636SJames Wright 
4319529d636SJames Wright     for (CeedInt f = 0; f < 4; f++) {
4329529d636SJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
4339529d636SJames Wright       v[f][i] = 0.;
4349529d636SJames Wright     }
4359529d636SJames Wright 
4369529d636SJames Wright     CeedScalar div_u = 0;
4379529d636SJames Wright     for (CeedInt j = 0; j < dim; j++) {
4389529d636SJames Wright       for (CeedInt k = 0; k < dim; k++) {
4399529d636SJames Wright         div_u += grad_s[k].Y.velocity[j];
4409529d636SJames Wright       }
4419529d636SJames Wright     }
4429529d636SJames Wright     CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
4439529d636SJames Wright 
4449529d636SJames Wright     CeedScalar uX[3] = {0.};
4459529d636SJames Wright     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
4469529d636SJames Wright 
4479529d636SJames Wright     if (context->strong_form) {  // Strong Galerkin convection term: v div(E u)
4489529d636SJames Wright       v[4][i] = -wdetJ * strong_conv;
4499529d636SJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = 0;
4509529d636SJames Wright     } else {  // Weak Galerkin convection term: -dv \cdot (E u)
4519529d636SJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = wdetJ * s.U.E_total * uX[j];
4529529d636SJames Wright       v[4][i] = 0.;
4539529d636SJames Wright     }
4549529d636SJames Wright 
455c8d249deSJames Wright     {  // Diffusion
456c8d249deSJames Wright       CeedScalar Fe[3], Fe_dXdx[3] = {0.};
457c8d249deSJames Wright 
458c8d249deSJames Wright       for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total;
459c8d249deSJames Wright       MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx);
460c8d249deSJames Wright       for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] += wdetJ * Fe_dXdx[k];
461c8d249deSJames Wright     }
462c8d249deSJames Wright 
463a78efa86SJames Wright     const CeedScalar TauS = Tau(context, s, dXdx, dim);
4645f952e8dSJames Wright     for (CeedInt j = 0; j < dim; j++) {
4655f952e8dSJames Wright       switch (context->stabilization) {
4669529d636SJames Wright         case STAB_NONE:
4679529d636SJames Wright           break;
4689529d636SJames Wright         case STAB_SU:
4695f952e8dSJames Wright         case STAB_SUPG: {
4705f952e8dSJames Wright           CeedScalar divFdiff_i = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? divFdiff[i] : 0.;
4715f952e8dSJames Wright           grad_v[j][4][i] -= wdetJ * TauS * (strong_conv + divFdiff_i) * uX[j];
4725f952e8dSJames Wright         } break;
4735f952e8dSJames Wright       }
4749529d636SJames Wright     }
4759529d636SJames Wright   }
47697cfd714SJames Wright   return 0;
4779529d636SJames Wright }
4789529d636SJames Wright 
4799529d636SJames Wright CEED_QFUNCTION(RHS_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
48097cfd714SJames Wright   return RHSFunction_AdvectionGeneric(ctx, Q, in, out, 3);
4819529d636SJames Wright }
4829529d636SJames Wright 
4839529d636SJames Wright CEED_QFUNCTION(RHS_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
48497cfd714SJames Wright   return RHSFunction_AdvectionGeneric(ctx, Q, in, out, 2);
4859529d636SJames Wright }
4869529d636SJames Wright 
4879529d636SJames Wright // *****************************************************************************
4889529d636SJames Wright // This QFunction implements consistent outflow and inflow BCs
4899529d636SJames Wright //      for advection
4909529d636SJames Wright //
4919529d636SJames Wright //  Inflow and outflow faces are determined based on sign(dot(wind, normal)):
4929529d636SJames Wright //    sign(dot(wind, normal)) > 0 : outflow BCs
4939529d636SJames Wright //    sign(dot(wind, normal)) < 0 : inflow BCs
4949529d636SJames Wright //
4959529d636SJames Wright //  Outflow BCs:
4969529d636SJames Wright //    The validity of the weak form of the governing equations is extended to the outflow and the current values of E are applied.
4979529d636SJames Wright //
4989529d636SJames Wright //  Inflow BCs:
4999529d636SJames Wright //    A prescribed Total Energy (E_wind) is applied weakly.
5009529d636SJames Wright // *****************************************************************************
50197cfd714SJames Wright CEED_QFUNCTION_HELPER int Advection_InOutFlowGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
5029529d636SJames Wright   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
5039529d636SJames Wright   const CeedScalar(*q_data_sur)    = in[2];
5049529d636SJames Wright 
5059529d636SJames Wright   CeedScalar(*v)[CEED_Q_VLA]   = (CeedScalar(*)[CEED_Q_VLA])out[0];
5069529d636SJames Wright   AdvectionContext context     = (AdvectionContext)ctx;
5079529d636SJames Wright   const CeedScalar E_wind      = context->E_wind;
5089529d636SJames Wright   const CeedScalar strong_form = context->strong_form;
5099529d636SJames Wright   const bool       is_implicit = context->implicit;
5109529d636SJames Wright 
5119529d636SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
5129529d636SJames Wright     const CeedScalar rho  = q[0][i];
5139529d636SJames Wright     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
5149529d636SJames Wright     const CeedScalar E    = q[4][i];
5159529d636SJames Wright 
51678e8b7daSJames Wright     CeedScalar wdetJb, normal[3];
51778e8b7daSJames Wright     QdataBoundaryUnpack_ND(dim, Q, i, q_data_sur, &wdetJb, NULL, normal);
5189529d636SJames Wright     wdetJb *= is_implicit ? -1. : 1.;
5199529d636SJames Wright 
52078e8b7daSJames Wright     const CeedScalar u_normal = DotN(normal, u, dim);
5219529d636SJames Wright 
5229529d636SJames Wright     // No Change in density or momentum
5239529d636SJames Wright     for (CeedInt j = 0; j < 4; j++) {
5249529d636SJames Wright       v[j][i] = 0;
5259529d636SJames Wright     }
5269529d636SJames Wright     // Implementing in/outflow BCs
5279529d636SJames Wright     if (u_normal > 0) {  // outflow
5289529d636SJames Wright       v[4][i] = -(1 - strong_form) * wdetJb * E * u_normal;
5299529d636SJames Wright     } else {  // inflow
5309529d636SJames Wright       v[4][i] = -(1 - strong_form) * wdetJb * E_wind * u_normal;
5319529d636SJames Wright     }
5329529d636SJames Wright   }
5339529d636SJames Wright   return 0;
5349529d636SJames Wright }
5359529d636SJames Wright 
5362b916ea7SJeremy L Thompson CEED_QFUNCTION(Advection_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
53797cfd714SJames Wright   return Advection_InOutFlowGeneric(ctx, Q, in, out, 3);
538a515125bSLeila Ghaffari }
539a515125bSLeila Ghaffari 
5409529d636SJames Wright CEED_QFUNCTION(Advection2d_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
54197cfd714SJames Wright   return Advection_InOutFlowGeneric(ctx, Q, in, out, 2);
5429529d636SJames Wright }
543c2d90829SJames Wright 
544c2d90829SJames Wright // @brief Volume integral for RHS of divergence of diffusive flux direct projection
545c2d90829SJames Wright CEED_QFUNCTION_HELPER int DivDiffusiveFluxVolumeRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out,
546c2d90829SJames Wright                                                                    const CeedInt dim) {
547c2d90829SJames Wright   const CeedScalar(*Grad_q)       = in[0];
548c2d90829SJames Wright   const CeedScalar(*q_data)       = in[1];
549c2d90829SJames Wright   CeedScalar(*Grad_v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
550c2d90829SJames Wright 
551c2d90829SJames Wright   AdvectionContext context = (AdvectionContext)ctx;
552c2d90829SJames Wright 
553c2d90829SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
554c2d90829SJames Wright     CeedScalar wdetJ, dXdx[9], F_diff[3] = {0.};
555c2d90829SJames Wright 
556c2d90829SJames Wright     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
557c2d90829SJames Wright     {  // Get physical diffusive flux
558c2d90829SJames Wright       CeedScalar Grad_qn[15], grad_E_ref[3];
559c2d90829SJames Wright 
560c2d90829SJames Wright       GradUnpackN(Q, i, 5, dim, Grad_q, Grad_qn);
561c2d90829SJames Wright       CopyN(&Grad_qn[4 * dim], grad_E_ref, dim);
562c2d90829SJames Wright       MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff);
563c2d90829SJames Wright       ScaleN(F_diff, -context->diffusion_coeff, dim);
564c2d90829SJames Wright     }
565c2d90829SJames Wright 
566c2d90829SJames Wright     CeedScalar F_diff_dXdx[3] = {0.};
567c2d90829SJames Wright     MatVecNM(dXdx, F_diff, dim, dim, CEED_NOTRANSPOSE, F_diff_dXdx);
568c2d90829SJames Wright     for (CeedInt k = 0; k < dim; k++) Grad_v[k][i] = -wdetJ * F_diff_dXdx[k];
569c2d90829SJames Wright   }
570c2d90829SJames Wright   return 0;
571c2d90829SJames Wright }
572c2d90829SJames Wright 
573c2d90829SJames Wright CEED_QFUNCTION(DivDiffusiveFluxVolumeRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
574c2d90829SJames Wright   return DivDiffusiveFluxVolumeRHS_AdvDif_Generic(ctx, Q, in, out, 2);
575c2d90829SJames Wright }
576c2d90829SJames Wright 
577c2d90829SJames Wright CEED_QFUNCTION(DivDiffusiveFluxVolumeRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
578c2d90829SJames Wright   return DivDiffusiveFluxVolumeRHS_AdvDif_Generic(ctx, Q, in, out, 3);
579c2d90829SJames Wright }
580c2d90829SJames Wright 
581c2d90829SJames Wright // @brief Boundary integral for RHS of divergence of diffusive flux direct projection
582c2d90829SJames Wright CEED_QFUNCTION_HELPER int DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out,
583c2d90829SJames Wright                                                                      const CeedInt dim) {
584c2d90829SJames Wright   const CeedScalar(*Grad_q) = in[0];
585c2d90829SJames Wright   const CeedScalar(*q_data) = in[1];
586c2d90829SJames Wright   CeedScalar(*v)            = out[0];
587c2d90829SJames Wright 
588c2d90829SJames Wright   AdvectionContext context = (AdvectionContext)ctx;
589c2d90829SJames Wright 
590c2d90829SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
591c2d90829SJames Wright     CeedScalar wdetJ, normal[3], dXdx[9], F_diff[3] = {0.};
592c2d90829SJames Wright 
593c2d90829SJames Wright     QdataBoundaryGradientUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx, normal);
594c2d90829SJames Wright     {  // Get physical diffusive flux
595c2d90829SJames Wright       CeedScalar Grad_qn[15], grad_E_ref[3];
596c2d90829SJames Wright 
597c2d90829SJames Wright       GradUnpackN(Q, i, 5, dim, Grad_q, Grad_qn);
598c2d90829SJames Wright       CopyN(&Grad_qn[4 * dim], grad_E_ref, dim);
599c2d90829SJames Wright       MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff);
600c2d90829SJames Wright       ScaleN(F_diff, -context->diffusion_coeff, dim);
601c2d90829SJames Wright     }
602c2d90829SJames Wright 
603c2d90829SJames Wright     v[i] = wdetJ * DotN(F_diff, normal, dim);
604c2d90829SJames Wright   }
605c2d90829SJames Wright   return 0;
606c2d90829SJames Wright }
607c2d90829SJames Wright 
608c2d90829SJames Wright CEED_QFUNCTION(DivDiffusiveFluxBoundaryRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
609c2d90829SJames Wright   return DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(ctx, Q, in, out, 2);
610c2d90829SJames Wright }
611c2d90829SJames Wright 
612c2d90829SJames Wright CEED_QFUNCTION(DivDiffusiveFluxBoundaryRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
613c2d90829SJames Wright   return DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(ctx, Q, in, out, 3);
614c2d90829SJames Wright }
61540b78511SJames Wright 
61640b78511SJames Wright // @brief Volume integral for RHS of diffusive flux indirect projection
61740b78511SJames Wright CEED_QFUNCTION_HELPER int DiffusiveFluxRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out,
61840b78511SJames Wright                                                           const CeedInt dim) {
61940b78511SJames Wright   const CeedScalar(*Grad_q)  = in[0];
62040b78511SJames Wright   const CeedScalar(*q_data)  = in[1];
62140b78511SJames Wright   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
62240b78511SJames Wright 
62340b78511SJames Wright   AdvectionContext context = (AdvectionContext)ctx;
62440b78511SJames Wright 
62540b78511SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
62640b78511SJames Wright     CeedScalar wdetJ, dXdx[9], F_diff[3] = {0.};
62740b78511SJames Wright 
62840b78511SJames Wright     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
62940b78511SJames Wright     {  // Get physical diffusive flux
63040b78511SJames Wright       CeedScalar Grad_qn[15], grad_E_ref[3];
63140b78511SJames Wright 
63240b78511SJames Wright       GradUnpackN(Q, i, 5, dim, Grad_q, Grad_qn);
63340b78511SJames Wright       CopyN(&Grad_qn[4 * dim], grad_E_ref, dim);
63440b78511SJames Wright       MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff);
63540b78511SJames Wright       ScaleN(F_diff, -context->diffusion_coeff, dim);
63640b78511SJames Wright     }
63740b78511SJames Wright     for (CeedInt k = 0; k < dim; k++) v[k][i] = wdetJ * F_diff[k];
63840b78511SJames Wright   }
63940b78511SJames Wright   return 0;
64040b78511SJames Wright }
64140b78511SJames Wright 
64240b78511SJames Wright CEED_QFUNCTION(DiffusiveFluxRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
64340b78511SJames Wright   return DiffusiveFluxRHS_AdvDif_Generic(ctx, Q, in, out, 2);
64440b78511SJames Wright }
64540b78511SJames Wright 
64640b78511SJames Wright CEED_QFUNCTION(DiffusiveFluxRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
64740b78511SJames Wright   return DiffusiveFluxRHS_AdvDif_Generic(ctx, Q, in, out, 3);
64840b78511SJames Wright }
649