xref: /libCEED/examples/fluids/qfunctions/advection.h (revision b4e9a8f894a0481205b699a80d40c97cb033c4b2)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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
10ba6664aeSJames Wright #include <ceed.h>
11c9c2c079SJeremy L Thompson #include <math.h>
1277841947SLeila Ghaffari 
13c44b1c7dSJames Wright #include "advection_types.h"
148f4d89c8SJames Wright #include "newtonian_state.h"
158f4d89c8SJames Wright #include "newtonian_types.h"
16c44b1c7dSJames Wright #include "stabilization_types.h"
178756a6ccSJames Wright #include "utils.h"
188756a6ccSJames Wright 
1977841947SLeila Ghaffari // *****************************************************************************
2093639ffbSJames Wright // This QFunction sets the initial conditions and the boundary conditions
2193639ffbSJames Wright //   for two test cases: ROTATION and TRANSLATION
2293639ffbSJames Wright //
2393639ffbSJames Wright // -- ROTATION (default)
2493639ffbSJames Wright //      Initial Conditions:
2593639ffbSJames Wright //        Mass Density:
2693639ffbSJames Wright //          Constant mass density of 1.0
2793639ffbSJames Wright //        Momentum Density:
2893639ffbSJames Wright //          Rotational field in x,y
2993639ffbSJames Wright //        Energy Density:
3093639ffbSJames Wright //          Maximum of 1. x0 decreasing linearly to 0. as radial distance
3193639ffbSJames Wright //            increases to (1.-r/rc), then 0. everywhere else
3293639ffbSJames Wright //
3393639ffbSJames Wright //      Boundary Conditions:
3493639ffbSJames Wright //        Mass Density:
3593639ffbSJames Wright //          0.0 flux
3693639ffbSJames Wright //        Momentum Density:
3793639ffbSJames Wright //          0.0
3893639ffbSJames Wright //        Energy Density:
3993639ffbSJames Wright //          0.0 flux
4093639ffbSJames Wright //
4193639ffbSJames Wright // -- TRANSLATION
4293639ffbSJames Wright //      Initial Conditions:
4393639ffbSJames Wright //        Mass Density:
4493639ffbSJames Wright //          Constant mass density of 1.0
4593639ffbSJames Wright //        Momentum Density:
4693639ffbSJames Wright //           Constant rectilinear field in x,y
4793639ffbSJames Wright //        Energy Density:
4893639ffbSJames Wright //          Maximum of 1. x0 decreasing linearly to 0. as radial distance
4993639ffbSJames Wright //            increases to (1.-r/rc), then 0. everywhere else
5093639ffbSJames Wright //
5193639ffbSJames Wright //      Boundary Conditions:
5293639ffbSJames Wright //        Mass Density:
5393639ffbSJames Wright //          0.0 flux
5493639ffbSJames Wright //        Momentum Density:
5593639ffbSJames Wright //          0.0
5693639ffbSJames Wright //        Energy Density:
5793639ffbSJames Wright //          Inflow BCs:
5893639ffbSJames Wright //            E = E_wind
5993639ffbSJames Wright //          Outflow BCs:
6093639ffbSJames Wright //            E = E(boundary)
6193639ffbSJames Wright //          Both In/Outflow BCs for E are applied weakly in the
6293639ffbSJames Wright //            QFunction "Advection2d_Sur"
6393639ffbSJames Wright //
6493639ffbSJames Wright // *****************************************************************************
6593639ffbSJames Wright 
6693639ffbSJames Wright // *****************************************************************************
6793639ffbSJames Wright // This helper function provides the exact, time-dependent solution and IC formulation for 2D advection
6893639ffbSJames Wright // *****************************************************************************
6993639ffbSJames Wright CEED_QFUNCTION_HELPER CeedInt Exact_AdvectionGeneric(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
7093639ffbSJames Wright   const SetupContextAdv context = (SetupContextAdv)ctx;
7193639ffbSJames Wright   const CeedScalar      rc      = context->rc;
7293639ffbSJames Wright   const CeedScalar      lx      = context->lx;
7393639ffbSJames Wright   const CeedScalar      ly      = context->ly;
7493639ffbSJames Wright   const CeedScalar      lz      = dim == 2 ? 0. : context->lz;
7593639ffbSJames Wright   const CeedScalar     *wind    = context->wind;
7693639ffbSJames Wright 
7793639ffbSJames Wright   const CeedScalar center[3] = {0.5 * lx, 0.5 * ly, 0.5 * lz};
7893639ffbSJames Wright   const CeedScalar theta     = dim == 2 ? M_PI / 3 : M_PI;
7993639ffbSJames Wright   const CeedScalar x0[3]     = {center[0] + .25 * lx * cos(theta + time), center[1] + .25 * ly * sin(theta + time), 0.5 * lz};
8093639ffbSJames Wright 
8193639ffbSJames Wright   const CeedScalar x = X[0], y = X[1], z = dim == 2 ? 0. : X[2];
8293639ffbSJames Wright 
8393639ffbSJames Wright   CeedScalar r = 0.;
8493639ffbSJames Wright   switch (context->initial_condition_type) {
8593639ffbSJames Wright     case ADVECTIONIC_BUBBLE_SPHERE:
8693639ffbSJames Wright     case ADVECTIONIC_BUBBLE_CYLINDER:
8793639ffbSJames Wright       r = sqrt(Square(x - x0[0]) + Square(y - x0[1]) + Square(z - x0[2]));
8893639ffbSJames Wright       break;
8993639ffbSJames Wright     case ADVECTIONIC_COSINE_HILL:
9093639ffbSJames Wright       r = sqrt(Square(x - center[0]) + Square(y - center[1]));
9193639ffbSJames Wright       break;
9293639ffbSJames Wright     case ADVECTIONIC_SKEW:
9393639ffbSJames Wright       break;
9493639ffbSJames Wright   }
9593639ffbSJames Wright 
9693639ffbSJames Wright   switch (context->wind_type) {
9793639ffbSJames Wright     case WIND_ROTATION:
9893639ffbSJames Wright       q[0] = 1.;
9993639ffbSJames Wright       q[1] = -(y - center[1]);
10093639ffbSJames Wright       q[2] = (x - center[0]);
10193639ffbSJames Wright       q[3] = 0;
10293639ffbSJames Wright       break;
10393639ffbSJames Wright     case WIND_TRANSLATION:
10493639ffbSJames Wright       q[0] = 1.;
10593639ffbSJames Wright       q[1] = wind[0];
10693639ffbSJames Wright       q[2] = wind[1];
10793639ffbSJames Wright       q[3] = dim == 2 ? 0. : wind[2];
10893639ffbSJames Wright       break;
10993639ffbSJames Wright     default:
11093639ffbSJames Wright       return 1;
11193639ffbSJames Wright   }
11293639ffbSJames Wright 
11393639ffbSJames Wright   switch (context->initial_condition_type) {
11493639ffbSJames Wright     case ADVECTIONIC_BUBBLE_SPHERE:
11593639ffbSJames Wright     case ADVECTIONIC_BUBBLE_CYLINDER:
11693639ffbSJames Wright       switch (context->bubble_continuity_type) {
11793639ffbSJames Wright         // original continuous, smooth shape
11893639ffbSJames Wright         case BUBBLE_CONTINUITY_SMOOTH:
11993639ffbSJames Wright           q[4] = r <= rc ? (1. - r / rc) : 0.;
12093639ffbSJames Wright           break;
12193639ffbSJames Wright         // discontinuous, sharp back half shape
12293639ffbSJames Wright         case BUBBLE_CONTINUITY_BACK_SHARP:
12393639ffbSJames Wright           q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) : 0.;
12493639ffbSJames Wright           break;
12593639ffbSJames Wright         // attempt to define a finite thickness that will get resolved under grid refinement
12693639ffbSJames Wright         case BUBBLE_CONTINUITY_THICK:
12793639ffbSJames Wright           q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) * fmin(1.0, (center[1] - y) / 1.25) : 0.;
12893639ffbSJames Wright           break;
12993639ffbSJames Wright         case BUBBLE_CONTINUITY_COSINE:
13093639ffbSJames Wright           q[4] = r <= rc ? .5 + .5 * cos(r * M_PI / rc) : 0;
13193639ffbSJames Wright           break;
13293639ffbSJames Wright       }
13393639ffbSJames Wright       break;
13493639ffbSJames Wright     case ADVECTIONIC_COSINE_HILL: {
13593639ffbSJames Wright       CeedScalar half_width = context->lx / 2;
13693639ffbSJames Wright       q[4]                  = r > half_width ? 0. : cos(2 * M_PI * r / half_width + M_PI) + 1.;
13793639ffbSJames Wright     } break;
13893639ffbSJames Wright     case ADVECTIONIC_SKEW: {
13993639ffbSJames Wright       CeedScalar       skewed_barrier[3]  = {wind[0], wind[1], 0};
14093639ffbSJames Wright       CeedScalar       inflow_to_point[3] = {x - context->lx / 2, y, 0};
14193639ffbSJames Wright       CeedScalar       cross_product[3]   = {0};
14293639ffbSJames Wright       const CeedScalar boundary_threshold = 20 * CEED_EPSILON;
14393639ffbSJames Wright       Cross3(skewed_barrier, inflow_to_point, cross_product);
14493639ffbSJames Wright 
14593639ffbSJames Wright       q[4] = cross_product[2] > boundary_threshold ? 0 : 1;
14693639ffbSJames Wright       if ((x < boundary_threshold && wind[0] < boundary_threshold) ||                // outflow at -x boundary
14793639ffbSJames Wright           (y < boundary_threshold && wind[1] < boundary_threshold) ||                // outflow at -y boundary
14893639ffbSJames Wright           (x > context->lx - boundary_threshold && wind[0] > boundary_threshold) ||  // outflow at +x boundary
14993639ffbSJames Wright           (y > context->ly - boundary_threshold && wind[1] > boundary_threshold)     // outflow at +y boundary
15093639ffbSJames Wright       ) {
15193639ffbSJames Wright         q[4] = 0;
15293639ffbSJames Wright       }
15393639ffbSJames Wright     } break;
15493639ffbSJames Wright   }
15593639ffbSJames Wright   return 0;
15693639ffbSJames Wright }
15793639ffbSJames Wright 
15893639ffbSJames Wright // *****************************************************************************
15977841947SLeila Ghaffari // This QFunction sets the initial conditions for 3D advection
16077841947SLeila Ghaffari // *****************************************************************************
1612b730f8bSJeremy L Thompson CEED_QFUNCTION(ICsAdvection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
16277841947SLeila Ghaffari   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
16377841947SLeila Ghaffari   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
16477841947SLeila Ghaffari 
16546603fc5SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
16677841947SLeila Ghaffari     const CeedScalar x[]  = {X[0][i], X[1][i], X[2][i]};
167e6225c47SLeila Ghaffari     CeedScalar       q[5] = {0.};
16877841947SLeila Ghaffari 
16930e1b2c7SJames Wright     Exact_AdvectionGeneric(3, 0., x, 5, q, ctx);
17077841947SLeila Ghaffari     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
17130e1b2c7SJames Wright   }
17277841947SLeila Ghaffari   return 0;
17377841947SLeila Ghaffari }
17477841947SLeila Ghaffari 
17577841947SLeila Ghaffari // *****************************************************************************
17693639ffbSJames Wright // This QFunction sets the initial conditions for 2D advection
17777841947SLeila Ghaffari // *****************************************************************************
17893639ffbSJames Wright CEED_QFUNCTION(ICsAdvection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
17993639ffbSJames Wright   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
18093639ffbSJames Wright   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
18193639ffbSJames Wright   const SetupContextAdv context    = (SetupContextAdv)ctx;
18293639ffbSJames Wright 
18393639ffbSJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
18493639ffbSJames Wright     const CeedScalar x[]  = {X[0][i], X[1][i]};
18593639ffbSJames Wright     CeedScalar       q[5] = {0.};
18693639ffbSJames Wright 
18793639ffbSJames Wright     Exact_AdvectionGeneric(2, context->time, x, 5, q, ctx);
18893639ffbSJames Wright     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
18993639ffbSJames Wright   }
19077841947SLeila Ghaffari   return 0;
19177841947SLeila Ghaffari }
19277841947SLeila Ghaffari 
19393639ffbSJames Wright CEED_QFUNCTION_HELPER void QdataUnpack_ND(CeedInt N, CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar *dXdx) {
194b1559114SJames Wright   // Cannot directly use QdataUnpack* helper functions due to SYCL online compiler incompatabilities
19593639ffbSJames Wright   switch (N) {
19693639ffbSJames Wright     case 2:
197b1559114SJames Wright       StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ);
198b1559114SJames Wright       StoredValuesUnpack(Q, i, 1, 4, q_data, dXdx);
19993639ffbSJames Wright       break;
20093639ffbSJames Wright     case 3:
201b1559114SJames Wright       StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ);
202b1559114SJames Wright       StoredValuesUnpack(Q, i, 1, 9, q_data, dXdx);
20393639ffbSJames Wright       break;
20493639ffbSJames Wright   }
20593639ffbSJames Wright }
20693639ffbSJames Wright 
20793639ffbSJames Wright CEED_QFUNCTION_HELPER int QdataBoundaryUnpack_ND(CeedInt N, CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar *dXdx,
20893639ffbSJames Wright                                                  CeedScalar *normal) {
209b1559114SJames Wright   // Cannot directly use QdataBoundaryUnpack* helper functions due to SYCL online compiler incompatabilities
21093639ffbSJames Wright   switch (N) {
21193639ffbSJames Wright     case 2:
212b1559114SJames Wright       if (wdetJ) StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ);
213b1559114SJames Wright       if (normal) StoredValuesUnpack(Q, i, 1, 2, q_data, normal);
21493639ffbSJames Wright       break;
21593639ffbSJames Wright     case 3:
216b1559114SJames Wright       if (wdetJ) StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ);
217b1559114SJames Wright       if (normal) StoredValuesUnpack(Q, i, 1, 3, q_data, normal);
218b1559114SJames Wright       if (dXdx) StoredValuesUnpack(Q, i, 4, 6, q_data, (CeedScalar *)dXdx);
21993639ffbSJames Wright       break;
22093639ffbSJames Wright   }
22193639ffbSJames Wright   return CEED_ERROR_SUCCESS;
22293639ffbSJames Wright }
22393639ffbSJames Wright 
22493639ffbSJames Wright CEED_QFUNCTION_HELPER void StatePhysicalGradientFromReference_ND(CeedInt N, CeedInt Q, CeedInt i, NewtonianIdealGasContext gas, State s,
22593639ffbSJames Wright                                                                  StateVariable state_var, const CeedScalar *grad_q, const CeedScalar *dXdx,
22693639ffbSJames Wright                                                                  State *grad_s) {
22793639ffbSJames Wright   switch (N) {
22893639ffbSJames Wright     case 2: {
22993639ffbSJames Wright       for (CeedInt k = 0; k < 2; k++) {
23093639ffbSJames Wright         CeedScalar dqi[5];
23193639ffbSJames Wright         for (CeedInt j = 0; j < 5; j++) {
23293639ffbSJames 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];
23393639ffbSJames Wright         }
23493639ffbSJames Wright         grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var);
23593639ffbSJames Wright       }
23693639ffbSJames Wright       CeedScalar U[5] = {0.};
23793639ffbSJames Wright       grad_s[2]       = StateFromU(gas, U);
23893639ffbSJames Wright     } break;
23993639ffbSJames Wright     case 3:
240b1559114SJames Wright       // Cannot directly use StatePhysicalGradientFromReference helper functions due to SYCL online compiler incompatabilities
241b1559114SJames Wright       for (CeedInt k = 0; k < 3; k++) {
242b1559114SJames Wright         CeedScalar dqi[5];
243b1559114SJames Wright         for (CeedInt j = 0; j < 5; j++) {
244b1559114SJames 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] +
245b1559114SJames Wright                    grad_q[(Q * 5) * 2 + Q * j + i] * dXdx[2 * N + k];
246b1559114SJames Wright         }
247b1559114SJames Wright         grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var);
248b1559114SJames Wright       }
24993639ffbSJames Wright       break;
25093639ffbSJames Wright   }
25193639ffbSJames Wright }
25293639ffbSJames Wright 
253*b4e9a8f8SJames Wright // @brief Calculate the stabilization constant \tau
254*b4e9a8f8SJames Wright CEED_QFUNCTION_HELPER CeedScalar Tau(AdvectionContext context, const State s, const CeedScalar *dXdx, CeedInt dim) {
255*b4e9a8f8SJames Wright   switch (context->stabilization_tau) {
256*b4e9a8f8SJames Wright     case STAB_TAU_CTAU: {
257*b4e9a8f8SJames Wright       CeedScalar uX[3] = {0.};
258*b4e9a8f8SJames Wright 
259*b4e9a8f8SJames Wright       MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
260*b4e9a8f8SJames Wright       return context->CtauS / sqrt(DotN(uX, uX, dim));
261*b4e9a8f8SJames Wright     } break;
262*b4e9a8f8SJames Wright     case STAB_TAU_ADVDIFF_SHAKIB: {
263*b4e9a8f8SJames Wright       CeedScalar gijd_mat[9] = {0.}, gij_uj[3] = {0.};
264*b4e9a8f8SJames Wright 
265*b4e9a8f8SJames Wright       MatMatN(dXdx, dXdx, dim, CEED_TRANSPOSE, CEED_NOTRANSPOSE, gijd_mat);
266*b4e9a8f8SJames Wright       MatVecNM(gijd_mat, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, gij_uj);
267*b4e9a8f8SJames Wright       return 1 / sqrt(Square(2 * context->Ctau_t / context->dt) + DotN(s.Y.velocity, gij_uj, dim) * context->Ctau_a);
268*b4e9a8f8SJames Wright     } break;
269*b4e9a8f8SJames Wright     default:
270*b4e9a8f8SJames Wright       return 0.;
271*b4e9a8f8SJames Wright   }
272*b4e9a8f8SJames Wright }
273*b4e9a8f8SJames Wright 
27493639ffbSJames Wright // *****************************************************************************
27593639ffbSJames Wright // This QFunction implements Advection for implicit time stepping method
27693639ffbSJames Wright // *****************************************************************************
27793639ffbSJames Wright CEED_QFUNCTION_HELPER void IFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
27893639ffbSJames Wright   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
27993639ffbSJames Wright   const CeedScalar(*grad_q)            = in[1];
28093639ffbSJames Wright   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
28193639ffbSJames Wright   const CeedScalar(*q_data)            = in[3];
28293639ffbSJames Wright 
28393639ffbSJames Wright   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
28493639ffbSJames Wright   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
28593639ffbSJames Wright   CeedScalar *jac_data               = out[2];
28693639ffbSJames Wright 
28793639ffbSJames Wright   AdvectionContext                 context   = (AdvectionContext)ctx;
28893639ffbSJames Wright   const CeedScalar                 zeros[14] = {0.};
28993639ffbSJames Wright   NewtonianIdealGasContext         gas;
29093639ffbSJames Wright   struct NewtonianIdealGasContext_ gas_struct = {0};
29193639ffbSJames Wright   gas                                         = &gas_struct;
29293639ffbSJames Wright 
29393639ffbSJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
29493639ffbSJames Wright     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
29593639ffbSJames Wright     const State      s     = StateFromU(gas, qi);
29693639ffbSJames Wright 
29793639ffbSJames Wright     CeedScalar wdetJ, dXdx[9];
29893639ffbSJames Wright     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
29993639ffbSJames Wright     State grad_s[3];
30093639ffbSJames Wright     StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
30193639ffbSJames Wright 
30293639ffbSJames Wright     const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
30393639ffbSJames Wright 
30493639ffbSJames Wright     for (CeedInt f = 0; f < 4; f++) {
30593639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
30693639ffbSJames Wright       v[f][i] = wdetJ * q_dot[f][i];                          // K Mass/transient term
30793639ffbSJames Wright     }
30893639ffbSJames Wright 
30993639ffbSJames Wright     CeedScalar div_u = 0;
31093639ffbSJames Wright     for (CeedInt j = 0; j < dim; j++) {
31193639ffbSJames Wright       for (CeedInt k = 0; k < dim; k++) {
31293639ffbSJames Wright         div_u += grad_s[k].Y.velocity[j];
31393639ffbSJames Wright       }
31493639ffbSJames Wright     }
31593639ffbSJames Wright     CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
31693639ffbSJames Wright     CeedScalar strong_res  = q_dot[4][i] + strong_conv;
31793639ffbSJames Wright 
31893639ffbSJames Wright     v[4][i] = wdetJ * q_dot[4][i];  // transient part (ALWAYS)
31993639ffbSJames Wright 
32093639ffbSJames Wright     CeedScalar uX[3] = {0.};
32193639ffbSJames Wright     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
32293639ffbSJames Wright 
32393639ffbSJames Wright     if (context->strong_form) {  // Strong Galerkin convection term: v div(E u)
32493639ffbSJames Wright       v[4][i] += wdetJ * strong_conv;
32593639ffbSJames Wright     } else {  // Weak Galerkin convection term: -dv \cdot (E u)
32693639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = -wdetJ * s.U.E_total * uX[j];
32793639ffbSJames Wright     }
32893639ffbSJames Wright 
329*b4e9a8f8SJames Wright     const CeedScalar TauS = Tau(context, s, dXdx, dim);
33093639ffbSJames Wright     for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) {
33193639ffbSJames Wright         case STAB_NONE:
33293639ffbSJames Wright           break;
33393639ffbSJames Wright         case STAB_SU:
33493639ffbSJames Wright           grad_v[j][4][i] += wdetJ * TauS * strong_conv * uX[j];
33593639ffbSJames Wright           break;
33693639ffbSJames Wright         case STAB_SUPG:
33793639ffbSJames Wright           grad_v[j][4][i] += wdetJ * TauS * strong_res * uX[j];
33893639ffbSJames Wright           break;
33993639ffbSJames Wright       }
34093639ffbSJames Wright     StoredValuesPack(Q, i, 0, 14, zeros, jac_data);
34193639ffbSJames Wright   }
34293639ffbSJames Wright }
34393639ffbSJames Wright 
3442b730f8bSJeremy L Thompson CEED_QFUNCTION(IFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
345372d1924SJames Wright   IFunction_AdvectionGeneric(ctx, Q, in, out, 3);
34677841947SLeila Ghaffari   return 0;
34777841947SLeila Ghaffari }
34877841947SLeila Ghaffari 
34993639ffbSJames Wright CEED_QFUNCTION(IFunction_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
35093639ffbSJames Wright   IFunction_AdvectionGeneric(ctx, Q, in, out, 2);
35193639ffbSJames Wright   return 0;
35293639ffbSJames Wright }
35393639ffbSJames Wright 
354*b4e9a8f8SJames Wright CEED_QFUNCTION_HELPER void MassFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
355*b4e9a8f8SJames Wright   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
356*b4e9a8f8SJames Wright   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[1];
357*b4e9a8f8SJames Wright   const CeedScalar(*q_data)            = in[2];
358*b4e9a8f8SJames Wright 
359*b4e9a8f8SJames Wright   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
360*b4e9a8f8SJames Wright   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
361*b4e9a8f8SJames Wright 
362*b4e9a8f8SJames Wright   AdvectionContext                 context    = (AdvectionContext)ctx;
363*b4e9a8f8SJames Wright   struct NewtonianIdealGasContext_ gas_struct = {0};
364*b4e9a8f8SJames Wright   NewtonianIdealGasContext         gas        = &gas_struct;
365*b4e9a8f8SJames Wright 
366*b4e9a8f8SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
367*b4e9a8f8SJames Wright     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
368*b4e9a8f8SJames Wright     const State      s     = StateFromU(gas, qi);
369*b4e9a8f8SJames Wright     CeedScalar       wdetJ, dXdx[9];
370*b4e9a8f8SJames Wright     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
371*b4e9a8f8SJames Wright 
372*b4e9a8f8SJames Wright     for (CeedInt f = 0; f < 4; f++) {
373*b4e9a8f8SJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
374*b4e9a8f8SJames Wright       v[f][i] = wdetJ * q_dot[f][i];                          // K Mass/transient term
375*b4e9a8f8SJames Wright     }
376*b4e9a8f8SJames Wright 
377*b4e9a8f8SJames Wright     // Unstabilized mass term
378*b4e9a8f8SJames Wright     v[4][i] = wdetJ * q_dot[4][i];
379*b4e9a8f8SJames Wright 
380*b4e9a8f8SJames Wright     // Stabilized mass term
381*b4e9a8f8SJames Wright     CeedScalar uX[3] = {0.};
382*b4e9a8f8SJames Wright     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
383*b4e9a8f8SJames Wright     const CeedScalar TauS = Tau(context, s, dXdx, dim);
384*b4e9a8f8SJames Wright     for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) {
385*b4e9a8f8SJames Wright         case STAB_NONE:
386*b4e9a8f8SJames Wright         case STAB_SU:
387*b4e9a8f8SJames Wright           grad_v[j][4][i] = 0;
388*b4e9a8f8SJames Wright           break;  // These should be run with the unstabilized mass matrix anyways
389*b4e9a8f8SJames Wright         case STAB_SUPG:
390*b4e9a8f8SJames Wright           grad_v[j][4][i] = wdetJ * TauS * q_dot[4][i] * uX[j];
391*b4e9a8f8SJames Wright           break;
392*b4e9a8f8SJames Wright       }
393*b4e9a8f8SJames Wright   }
394*b4e9a8f8SJames Wright }
395*b4e9a8f8SJames Wright 
396*b4e9a8f8SJames Wright CEED_QFUNCTION(MassFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
397*b4e9a8f8SJames Wright   MassFunction_AdvectionGeneric(ctx, Q, in, out, 3);
398*b4e9a8f8SJames Wright   return 0;
399*b4e9a8f8SJames Wright }
400*b4e9a8f8SJames Wright 
401*b4e9a8f8SJames Wright CEED_QFUNCTION(MassFunction_Advection2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
402*b4e9a8f8SJames Wright   MassFunction_AdvectionGeneric(ctx, Q, in, out, 2);
403*b4e9a8f8SJames Wright   return 0;
404*b4e9a8f8SJames Wright }
405*b4e9a8f8SJames Wright 
40693639ffbSJames Wright // *****************************************************************************
40793639ffbSJames Wright // This QFunction implements Advection for explicit time stepping method
40893639ffbSJames Wright // *****************************************************************************
40993639ffbSJames Wright CEED_QFUNCTION_HELPER void RHSFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
41093639ffbSJames Wright   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
41193639ffbSJames Wright   const CeedScalar(*grad_q)        = in[1];
41293639ffbSJames Wright   const CeedScalar(*q_data)        = in[2];
41393639ffbSJames Wright 
41493639ffbSJames Wright   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
41593639ffbSJames Wright   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
41693639ffbSJames Wright 
41793639ffbSJames Wright   AdvectionContext                 context    = (AdvectionContext)ctx;
41893639ffbSJames Wright   struct NewtonianIdealGasContext_ gas_struct = {0};
419*b4e9a8f8SJames Wright   NewtonianIdealGasContext         gas        = &gas_struct;
42093639ffbSJames Wright 
42193639ffbSJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
42293639ffbSJames Wright     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
42393639ffbSJames Wright     const State      s     = StateFromU(gas, qi);
42493639ffbSJames Wright 
42593639ffbSJames Wright     CeedScalar wdetJ, dXdx[9];
42693639ffbSJames Wright     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
42793639ffbSJames Wright     State grad_s[3];
42893639ffbSJames Wright     StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
42993639ffbSJames Wright 
43093639ffbSJames Wright     const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
43193639ffbSJames Wright 
43293639ffbSJames Wright     for (CeedInt f = 0; f < 4; f++) {
43393639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
43493639ffbSJames Wright       v[f][i] = 0.;
43593639ffbSJames Wright     }
43693639ffbSJames Wright 
43793639ffbSJames Wright     CeedScalar div_u = 0;
43893639ffbSJames Wright     for (CeedInt j = 0; j < dim; j++) {
43993639ffbSJames Wright       for (CeedInt k = 0; k < dim; k++) {
44093639ffbSJames Wright         div_u += grad_s[k].Y.velocity[j];
44193639ffbSJames Wright       }
44293639ffbSJames Wright     }
44393639ffbSJames Wright     CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
44493639ffbSJames Wright 
44593639ffbSJames Wright     CeedScalar uX[3] = {0.};
44693639ffbSJames Wright     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
44793639ffbSJames Wright 
44893639ffbSJames Wright     if (context->strong_form) {  // Strong Galerkin convection term: v div(E u)
44993639ffbSJames Wright       v[4][i] = -wdetJ * strong_conv;
45093639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = 0;
45193639ffbSJames Wright     } else {  // Weak Galerkin convection term: -dv \cdot (E u)
45293639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = wdetJ * s.U.E_total * uX[j];
45393639ffbSJames Wright       v[4][i] = 0.;
45493639ffbSJames Wright     }
45593639ffbSJames Wright 
456*b4e9a8f8SJames Wright     const CeedScalar TauS = Tau(context, s, dXdx, dim);
45793639ffbSJames Wright     for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) {
45893639ffbSJames Wright         case STAB_NONE:
45993639ffbSJames Wright           break;
46093639ffbSJames Wright         case STAB_SU:
46193639ffbSJames Wright         case STAB_SUPG:
462fa9a7b5dSJames Wright           grad_v[j][4][i] -= wdetJ * TauS * strong_conv * uX[j];
46393639ffbSJames Wright           break;
46493639ffbSJames Wright       }
46593639ffbSJames Wright   }
46693639ffbSJames Wright }
46793639ffbSJames Wright 
46893639ffbSJames Wright CEED_QFUNCTION(RHS_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
46993639ffbSJames Wright   RHSFunction_AdvectionGeneric(ctx, Q, in, out, 3);
47093639ffbSJames Wright   return 0;
47193639ffbSJames Wright }
47293639ffbSJames Wright 
47393639ffbSJames Wright CEED_QFUNCTION(RHS_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
47493639ffbSJames Wright   RHSFunction_AdvectionGeneric(ctx, Q, in, out, 2);
47593639ffbSJames Wright   return 0;
47693639ffbSJames Wright }
47793639ffbSJames Wright 
47893639ffbSJames Wright // *****************************************************************************
47993639ffbSJames Wright // This QFunction implements consistent outflow and inflow BCs
48093639ffbSJames Wright //      for advection
48193639ffbSJames Wright //
48293639ffbSJames Wright //  Inflow and outflow faces are determined based on sign(dot(wind, normal)):
48393639ffbSJames Wright //    sign(dot(wind, normal)) > 0 : outflow BCs
48493639ffbSJames Wright //    sign(dot(wind, normal)) < 0 : inflow BCs
48593639ffbSJames Wright //
48693639ffbSJames Wright //  Outflow BCs:
48793639ffbSJames Wright //    The validity of the weak form of the governing equations is extended to the outflow and the current values of E are applied.
48893639ffbSJames Wright //
48993639ffbSJames Wright //  Inflow BCs:
49093639ffbSJames Wright //    A prescribed Total Energy (E_wind) is applied weakly.
49193639ffbSJames Wright // *****************************************************************************
49293639ffbSJames Wright CEED_QFUNCTION(Advection_InOutFlowGeneric)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
49393639ffbSJames Wright   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
49493639ffbSJames Wright   const CeedScalar(*q_data_sur)    = in[2];
49593639ffbSJames Wright 
49693639ffbSJames Wright   CeedScalar(*v)[CEED_Q_VLA]   = (CeedScalar(*)[CEED_Q_VLA])out[0];
49793639ffbSJames Wright   AdvectionContext context     = (AdvectionContext)ctx;
49893639ffbSJames Wright   const CeedScalar E_wind      = context->E_wind;
49993639ffbSJames Wright   const CeedScalar strong_form = context->strong_form;
50093639ffbSJames Wright   const bool       is_implicit = context->implicit;
50193639ffbSJames Wright 
50293639ffbSJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
50393639ffbSJames Wright     const CeedScalar rho  = q[0][i];
50493639ffbSJames Wright     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
50593639ffbSJames Wright     const CeedScalar E    = q[4][i];
50693639ffbSJames Wright 
50793639ffbSJames Wright     CeedScalar wdetJb, norm[3];
50893639ffbSJames Wright     QdataBoundaryUnpack_ND(dim, Q, i, q_data_sur, &wdetJb, NULL, norm);
50993639ffbSJames Wright     wdetJb *= is_implicit ? -1. : 1.;
51093639ffbSJames Wright 
51193639ffbSJames Wright     const CeedScalar u_normal = DotN(norm, u, dim);
51293639ffbSJames Wright 
51393639ffbSJames Wright     // No Change in density or momentum
51493639ffbSJames Wright     for (CeedInt j = 0; j < 4; j++) {
51593639ffbSJames Wright       v[j][i] = 0;
51693639ffbSJames Wright     }
51793639ffbSJames Wright     // Implementing in/outflow BCs
51893639ffbSJames Wright     if (u_normal > 0) {  // outflow
51993639ffbSJames Wright       v[4][i] = -(1 - strong_form) * wdetJb * E * u_normal;
52093639ffbSJames Wright     } else {  // inflow
52193639ffbSJames Wright       v[4][i] = -(1 - strong_form) * wdetJb * E_wind * u_normal;
52293639ffbSJames Wright     }
52393639ffbSJames Wright   }
52493639ffbSJames Wright   return 0;
52593639ffbSJames Wright }
52693639ffbSJames Wright 
5272b730f8bSJeremy L Thompson CEED_QFUNCTION(Advection_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
5284bdcf5bfSJames Wright   Advection_InOutFlowGeneric(ctx, Q, in, out, 3);
52977841947SLeila Ghaffari   return 0;
53077841947SLeila Ghaffari }
53177841947SLeila Ghaffari 
53293639ffbSJames Wright CEED_QFUNCTION(Advection2d_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
53393639ffbSJames Wright   Advection_InOutFlowGeneric(ctx, Q, in, out, 2);
53493639ffbSJames Wright   return 0;
53593639ffbSJames Wright }
536