xref: /libCEED/examples/fluids/qfunctions/advection.h (revision b18328c48bac463522bc9df72304dbd8dd7879bd)
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 // *****************************************************************************
2493639ffbSJames Wright // This QFunction sets the initial conditions and the boundary conditions
2593639ffbSJames Wright //   for two test cases: ROTATION and TRANSLATION
2693639ffbSJames Wright //
2793639ffbSJames Wright // -- ROTATION (default)
2893639ffbSJames Wright //      Initial Conditions:
2993639ffbSJames Wright //        Mass Density:
3093639ffbSJames Wright //          Constant mass density of 1.0
3193639ffbSJames Wright //        Momentum Density:
3293639ffbSJames Wright //          Rotational field in x,y
3393639ffbSJames Wright //        Energy Density:
3493639ffbSJames Wright //          Maximum of 1. x0 decreasing linearly to 0. as radial distance
3593639ffbSJames Wright //            increases to (1.-r/rc), then 0. everywhere else
3693639ffbSJames Wright //
3793639ffbSJames Wright //      Boundary Conditions:
3893639ffbSJames Wright //        Mass Density:
3993639ffbSJames Wright //          0.0 flux
4093639ffbSJames Wright //        Momentum Density:
4193639ffbSJames Wright //          0.0
4293639ffbSJames Wright //        Energy Density:
4393639ffbSJames Wright //          0.0 flux
4493639ffbSJames Wright //
4593639ffbSJames Wright // -- TRANSLATION
4693639ffbSJames Wright //      Initial Conditions:
4793639ffbSJames Wright //        Mass Density:
4893639ffbSJames Wright //          Constant mass density of 1.0
4993639ffbSJames Wright //        Momentum Density:
5093639ffbSJames Wright //           Constant rectilinear field in x,y
5193639ffbSJames Wright //        Energy Density:
5293639ffbSJames Wright //          Maximum of 1. x0 decreasing linearly to 0. as radial distance
5393639ffbSJames Wright //            increases to (1.-r/rc), then 0. everywhere else
5493639ffbSJames Wright //
5593639ffbSJames Wright //      Boundary Conditions:
5693639ffbSJames Wright //        Mass Density:
5793639ffbSJames Wright //          0.0 flux
5893639ffbSJames Wright //        Momentum Density:
5993639ffbSJames Wright //          0.0
6093639ffbSJames Wright //        Energy Density:
6193639ffbSJames Wright //          Inflow BCs:
6293639ffbSJames Wright //            E = E_wind
6393639ffbSJames Wright //          Outflow BCs:
6493639ffbSJames Wright //            E = E(boundary)
6593639ffbSJames Wright //          Both In/Outflow BCs for E are applied weakly in the
6693639ffbSJames Wright //            QFunction "Advection2d_Sur"
6793639ffbSJames Wright //
6893639ffbSJames Wright // *****************************************************************************
6993639ffbSJames Wright 
7093639ffbSJames Wright // *****************************************************************************
7193639ffbSJames Wright // This helper function provides the exact, time-dependent solution and IC formulation for 2D advection
7293639ffbSJames Wright // *****************************************************************************
7393639ffbSJames Wright CEED_QFUNCTION_HELPER CeedInt Exact_AdvectionGeneric(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
7493639ffbSJames Wright   const SetupContextAdv context = (SetupContextAdv)ctx;
7593639ffbSJames Wright   const CeedScalar      rc      = context->rc;
7693639ffbSJames Wright   const CeedScalar      lx      = context->lx;
7793639ffbSJames Wright   const CeedScalar      ly      = context->ly;
7893639ffbSJames Wright   const CeedScalar      lz      = dim == 2 ? 0. : context->lz;
7993639ffbSJames Wright   const CeedScalar     *wind    = context->wind;
8093639ffbSJames Wright 
8193639ffbSJames Wright   const CeedScalar center[3] = {0.5 * lx, 0.5 * ly, 0.5 * lz};
8293639ffbSJames Wright   const CeedScalar theta     = dim == 2 ? M_PI / 3 : M_PI;
8393639ffbSJames Wright   const CeedScalar x0[3]     = {center[0] + .25 * lx * cos(theta + time), center[1] + .25 * ly * sin(theta + time), 0.5 * lz};
8493639ffbSJames Wright 
8593639ffbSJames Wright   const CeedScalar x = X[0], y = X[1], z = dim == 2 ? 0. : X[2];
8693639ffbSJames Wright 
8793639ffbSJames Wright   CeedScalar r = 0.;
8893639ffbSJames Wright   switch (context->initial_condition_type) {
8993639ffbSJames Wright     case ADVECTIONIC_BUBBLE_SPHERE:
9093639ffbSJames Wright     case ADVECTIONIC_BUBBLE_CYLINDER:
9193639ffbSJames Wright       r = sqrt(Square(x - x0[0]) + Square(y - x0[1]) + Square(z - x0[2]));
9293639ffbSJames Wright       break;
9393639ffbSJames Wright     case ADVECTIONIC_COSINE_HILL:
9493639ffbSJames Wright       r = sqrt(Square(x - center[0]) + Square(y - center[1]));
9593639ffbSJames Wright       break;
9693639ffbSJames Wright     case ADVECTIONIC_SKEW:
9793639ffbSJames Wright       break;
9893639ffbSJames Wright   }
9993639ffbSJames Wright 
10093639ffbSJames Wright   switch (context->wind_type) {
10193639ffbSJames Wright     case WIND_ROTATION:
10293639ffbSJames Wright       q[0] = 1.;
10393639ffbSJames Wright       q[1] = -(y - center[1]);
10493639ffbSJames Wright       q[2] = (x - center[0]);
10593639ffbSJames Wright       q[3] = 0;
10693639ffbSJames Wright       break;
10793639ffbSJames Wright     case WIND_TRANSLATION:
10893639ffbSJames Wright       q[0] = 1.;
10993639ffbSJames Wright       q[1] = wind[0];
11093639ffbSJames Wright       q[2] = wind[1];
11193639ffbSJames Wright       q[3] = dim == 2 ? 0. : wind[2];
11293639ffbSJames Wright       break;
11393639ffbSJames Wright     default:
11493639ffbSJames Wright       return 1;
11593639ffbSJames Wright   }
11693639ffbSJames Wright 
11793639ffbSJames Wright   switch (context->initial_condition_type) {
11893639ffbSJames Wright     case ADVECTIONIC_BUBBLE_SPHERE:
11993639ffbSJames Wright     case ADVECTIONIC_BUBBLE_CYLINDER:
12093639ffbSJames Wright       switch (context->bubble_continuity_type) {
12193639ffbSJames Wright         // original continuous, smooth shape
12293639ffbSJames Wright         case BUBBLE_CONTINUITY_SMOOTH:
12393639ffbSJames Wright           q[4] = r <= rc ? (1. - r / rc) : 0.;
12493639ffbSJames Wright           break;
12593639ffbSJames Wright         // discontinuous, sharp back half shape
12693639ffbSJames Wright         case BUBBLE_CONTINUITY_BACK_SHARP:
12793639ffbSJames Wright           q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) : 0.;
12893639ffbSJames Wright           break;
12993639ffbSJames Wright         // attempt to define a finite thickness that will get resolved under grid refinement
13093639ffbSJames Wright         case BUBBLE_CONTINUITY_THICK:
13193639ffbSJames Wright           q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) * fmin(1.0, (center[1] - y) / 1.25) : 0.;
13293639ffbSJames Wright           break;
13393639ffbSJames Wright         case BUBBLE_CONTINUITY_COSINE:
13493639ffbSJames Wright           q[4] = r <= rc ? .5 + .5 * cos(r * M_PI / rc) : 0;
13593639ffbSJames Wright           break;
13693639ffbSJames Wright       }
13793639ffbSJames Wright       break;
13893639ffbSJames Wright     case ADVECTIONIC_COSINE_HILL: {
13993639ffbSJames Wright       CeedScalar half_width = context->lx / 2;
14093639ffbSJames Wright       q[4]                  = r > half_width ? 0. : cos(2 * M_PI * r / half_width + M_PI) + 1.;
14193639ffbSJames Wright     } break;
14293639ffbSJames Wright     case ADVECTIONIC_SKEW: {
14393639ffbSJames Wright       CeedScalar       skewed_barrier[3]  = {wind[0], wind[1], 0};
14493639ffbSJames Wright       CeedScalar       inflow_to_point[3] = {x - context->lx / 2, y, 0};
14593639ffbSJames Wright       CeedScalar       cross_product[3]   = {0};
14693639ffbSJames Wright       const CeedScalar boundary_threshold = 20 * CEED_EPSILON;
14793639ffbSJames Wright       Cross3(skewed_barrier, inflow_to_point, cross_product);
14893639ffbSJames Wright 
14993639ffbSJames Wright       q[4] = cross_product[2] > boundary_threshold ? 0 : 1;
15093639ffbSJames Wright       if ((x < boundary_threshold && wind[0] < boundary_threshold) ||                // outflow at -x boundary
15193639ffbSJames Wright           (y < boundary_threshold && wind[1] < boundary_threshold) ||                // outflow at -y boundary
15293639ffbSJames Wright           (x > context->lx - boundary_threshold && wind[0] > boundary_threshold) ||  // outflow at +x boundary
15393639ffbSJames Wright           (y > context->ly - boundary_threshold && wind[1] > boundary_threshold)     // outflow at +y boundary
15493639ffbSJames Wright       ) {
15593639ffbSJames Wright         q[4] = 0;
15693639ffbSJames Wright       }
15793639ffbSJames Wright     } break;
15893639ffbSJames Wright   }
15993639ffbSJames Wright   return 0;
16093639ffbSJames Wright }
16193639ffbSJames Wright 
16293639ffbSJames 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 // *****************************************************************************
18093639ffbSJames Wright // This QFunction sets the initial conditions for 2D advection
18177841947SLeila Ghaffari // *****************************************************************************
18293639ffbSJames Wright CEED_QFUNCTION(ICsAdvection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
18393639ffbSJames Wright   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
18493639ffbSJames Wright   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
18593639ffbSJames Wright   const SetupContextAdv context    = (SetupContextAdv)ctx;
18693639ffbSJames Wright 
18793639ffbSJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
18893639ffbSJames Wright     const CeedScalar x[]  = {X[0][i], X[1][i]};
18993639ffbSJames Wright     CeedScalar       q[5] = {0.};
19093639ffbSJames Wright 
19193639ffbSJames Wright     Exact_AdvectionGeneric(2, context->time, x, 5, q, ctx);
19293639ffbSJames Wright     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
19393639ffbSJames Wright   }
19477841947SLeila Ghaffari   return 0;
19577841947SLeila Ghaffari }
19677841947SLeila Ghaffari 
19793639ffbSJames Wright CEED_QFUNCTION_HELPER void QdataUnpack_ND(CeedInt N, CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar *dXdx) {
19893639ffbSJames Wright   switch (N) {
19993639ffbSJames Wright     case 2:
20093639ffbSJames Wright       QdataUnpack_2D(Q, i, q_data, wdetJ, (CeedScalar(*)[2])dXdx);
20193639ffbSJames Wright       break;
20293639ffbSJames Wright     case 3:
20393639ffbSJames Wright       QdataUnpack_3D(Q, i, q_data, wdetJ, (CeedScalar(*)[3])dXdx);
20493639ffbSJames Wright       break;
20593639ffbSJames Wright   }
20693639ffbSJames Wright }
20793639ffbSJames Wright 
20893639ffbSJames Wright CEED_QFUNCTION_HELPER int QdataBoundaryUnpack_ND(CeedInt N, CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar *dXdx,
20993639ffbSJames Wright                                                  CeedScalar *normal) {
21093639ffbSJames Wright   switch (N) {
21193639ffbSJames Wright     case 2:
21293639ffbSJames Wright       QdataBoundaryUnpack_2D(Q, i, q_data, wdetJ, normal);
21393639ffbSJames Wright       break;
21493639ffbSJames Wright     case 3:
21593639ffbSJames Wright       QdataBoundaryUnpack_3D(Q, i, q_data, wdetJ, (CeedScalar(*)[3])dXdx, normal);
21693639ffbSJames Wright       break;
21793639ffbSJames Wright   }
21893639ffbSJames Wright   return CEED_ERROR_SUCCESS;
21993639ffbSJames Wright }
22093639ffbSJames Wright 
22193639ffbSJames Wright CEED_QFUNCTION_HELPER void StatePhysicalGradientFromReference_ND(CeedInt N, CeedInt Q, CeedInt i, NewtonianIdealGasContext gas, State s,
22293639ffbSJames Wright                                                                  StateVariable state_var, const CeedScalar *grad_q, const CeedScalar *dXdx,
22393639ffbSJames Wright                                                                  State *grad_s) {
22493639ffbSJames Wright   switch (N) {
22593639ffbSJames Wright     case 2: {
22693639ffbSJames Wright       for (CeedInt k = 0; k < 2; k++) {
22793639ffbSJames Wright         CeedScalar dqi[5];
22893639ffbSJames Wright         for (CeedInt j = 0; j < 5; j++) {
22993639ffbSJames 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];
23093639ffbSJames Wright         }
23193639ffbSJames Wright         grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var);
23293639ffbSJames Wright       }
23393639ffbSJames Wright       CeedScalar U[5] = {0.};
23493639ffbSJames Wright       grad_s[2]       = StateFromU(gas, U);
23593639ffbSJames Wright     } break;
23693639ffbSJames Wright     case 3:
23793639ffbSJames Wright       StatePhysicalGradientFromReference(Q, i, gas, s, state_var, grad_q, (CeedScalar(*)[3])dXdx, grad_s);
23893639ffbSJames Wright       break;
23993639ffbSJames Wright   }
24093639ffbSJames Wright }
24193639ffbSJames Wright 
24293639ffbSJames Wright // *****************************************************************************
24393639ffbSJames Wright // This QFunction implements Advection for implicit time stepping method
24493639ffbSJames Wright // *****************************************************************************
24593639ffbSJames Wright CEED_QFUNCTION_HELPER void IFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
24693639ffbSJames Wright   const CeedScalar(*q)[CEED_Q_VLA]     = (const CeedScalar(*)[CEED_Q_VLA])in[0];
24793639ffbSJames Wright   const CeedScalar(*grad_q)            = in[1];
24893639ffbSJames Wright   const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
24993639ffbSJames Wright   const CeedScalar(*q_data)            = in[3];
25093639ffbSJames Wright 
25193639ffbSJames Wright   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
25293639ffbSJames Wright   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
25393639ffbSJames Wright   CeedScalar *jac_data               = out[2];
25493639ffbSJames Wright 
25593639ffbSJames Wright   AdvectionContext                 context   = (AdvectionContext)ctx;
25693639ffbSJames Wright   const CeedScalar                 CtauS     = context->CtauS;
25793639ffbSJames Wright   const CeedScalar                 zeros[14] = {0.};
25893639ffbSJames Wright   NewtonianIdealGasContext         gas;
25993639ffbSJames Wright   struct NewtonianIdealGasContext_ gas_struct = {0};
26093639ffbSJames Wright   gas                                         = &gas_struct;
26193639ffbSJames Wright 
26293639ffbSJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
26393639ffbSJames Wright     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
26493639ffbSJames Wright     const State      s     = StateFromU(gas, qi);
26593639ffbSJames Wright 
26693639ffbSJames Wright     CeedScalar wdetJ, dXdx[9];
26793639ffbSJames Wright     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
26893639ffbSJames Wright     State grad_s[3];
26993639ffbSJames Wright     StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
27093639ffbSJames Wright 
27193639ffbSJames Wright     const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
27293639ffbSJames Wright 
27393639ffbSJames Wright     for (CeedInt f = 0; f < 4; f++) {
27493639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
27593639ffbSJames Wright       v[f][i] = wdetJ * q_dot[f][i];                          // K Mass/transient term
27693639ffbSJames Wright     }
27793639ffbSJames Wright 
27893639ffbSJames Wright     CeedScalar div_u = 0;
27993639ffbSJames Wright     for (CeedInt j = 0; j < dim; j++) {
28093639ffbSJames Wright       for (CeedInt k = 0; k < dim; k++) {
28193639ffbSJames Wright         div_u += grad_s[k].Y.velocity[j];
28293639ffbSJames Wright       }
28393639ffbSJames Wright     }
28493639ffbSJames Wright     CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
28593639ffbSJames Wright     CeedScalar strong_res  = q_dot[4][i] + strong_conv;
28693639ffbSJames Wright 
28793639ffbSJames Wright     v[4][i] = wdetJ * q_dot[4][i];  // transient part (ALWAYS)
28893639ffbSJames Wright 
28993639ffbSJames Wright     CeedScalar uX[3] = {0.};
29093639ffbSJames Wright     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
29193639ffbSJames Wright 
29293639ffbSJames Wright     if (context->strong_form) {  // Strong Galerkin convection term: v div(E u)
29393639ffbSJames Wright       v[4][i] += wdetJ * strong_conv;
29493639ffbSJames Wright     } else {  // Weak Galerkin convection term: -dv \cdot (E u)
29593639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = -wdetJ * s.U.E_total * uX[j];
29693639ffbSJames Wright     }
29793639ffbSJames Wright 
298*b18328c4SJames Wright     CeedScalar TauS = 0;
299*b18328c4SJames Wright     switch (context->stabilization_tau) {
300*b18328c4SJames Wright       case STAB_TAU_CTAU:
301*b18328c4SJames Wright         TauS = CtauS / sqrt(Dot3(uX, uX));
302*b18328c4SJames Wright         break;
303*b18328c4SJames Wright       case STAB_TAU_ADVDIFF_SHAKIB: {
304*b18328c4SJames Wright         CeedScalar gijd_mat[9] = {0.}, gij_uj[3] = {0.};
305*b18328c4SJames Wright         MatMatN(dXdx, dXdx, dim, CEED_TRANSPOSE, CEED_NOTRANSPOSE, gijd_mat);
306*b18328c4SJames Wright 
307*b18328c4SJames Wright         MatVecNM(gijd_mat, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, gij_uj);
308*b18328c4SJames Wright         TauS = 1 / sqrt(Square(2 * context->Ctau_t / context->dt) + DotN(s.Y.velocity, gij_uj, dim) * context->Ctau_a);
309*b18328c4SJames Wright       } break;
310*b18328c4SJames Wright     }
311*b18328c4SJames Wright 
31293639ffbSJames Wright     for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) {
31393639ffbSJames Wright         case STAB_NONE:
31493639ffbSJames Wright           break;
31593639ffbSJames Wright         case STAB_SU:
31693639ffbSJames Wright           grad_v[j][4][i] += wdetJ * TauS * strong_conv * uX[j];
31793639ffbSJames Wright           break;
31893639ffbSJames Wright         case STAB_SUPG:
31993639ffbSJames Wright           grad_v[j][4][i] += wdetJ * TauS * strong_res * uX[j];
32093639ffbSJames Wright           break;
32193639ffbSJames Wright       }
32293639ffbSJames Wright     StoredValuesPack(Q, i, 0, 14, zeros, jac_data);
32393639ffbSJames Wright   }
32493639ffbSJames Wright }
32593639ffbSJames Wright 
3262b730f8bSJeremy L Thompson CEED_QFUNCTION(IFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
327372d1924SJames Wright   IFunction_AdvectionGeneric(ctx, Q, in, out, 3);
32877841947SLeila Ghaffari   return 0;
32977841947SLeila Ghaffari }
33077841947SLeila Ghaffari 
33193639ffbSJames Wright CEED_QFUNCTION(IFunction_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
33293639ffbSJames Wright   IFunction_AdvectionGeneric(ctx, Q, in, out, 2);
33393639ffbSJames Wright   return 0;
33493639ffbSJames Wright }
33593639ffbSJames Wright 
33693639ffbSJames Wright // *****************************************************************************
33793639ffbSJames Wright // This QFunction implements Advection for explicit time stepping method
33893639ffbSJames Wright // *****************************************************************************
33993639ffbSJames Wright CEED_QFUNCTION_HELPER void RHSFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
34093639ffbSJames Wright   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
34193639ffbSJames Wright   const CeedScalar(*grad_q)        = in[1];
34293639ffbSJames Wright   const CeedScalar(*q_data)        = in[2];
34393639ffbSJames Wright 
34493639ffbSJames Wright   CeedScalar(*v)[CEED_Q_VLA]         = (CeedScalar(*)[CEED_Q_VLA])out[0];
34593639ffbSJames Wright   CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
34693639ffbSJames Wright 
34793639ffbSJames Wright   AdvectionContext                 context = (AdvectionContext)ctx;
34893639ffbSJames Wright   const CeedScalar                 CtauS   = context->CtauS;
34993639ffbSJames Wright   NewtonianIdealGasContext         gas;
35093639ffbSJames Wright   struct NewtonianIdealGasContext_ gas_struct = {0};
35193639ffbSJames Wright   gas                                         = &gas_struct;
35293639ffbSJames Wright 
35393639ffbSJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
35493639ffbSJames Wright     const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
35593639ffbSJames Wright     const State      s     = StateFromU(gas, qi);
35693639ffbSJames Wright 
35793639ffbSJames Wright     CeedScalar wdetJ, dXdx[9];
35893639ffbSJames Wright     QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
35993639ffbSJames Wright     State grad_s[3];
36093639ffbSJames Wright     StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
36193639ffbSJames Wright 
36293639ffbSJames Wright     const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
36393639ffbSJames Wright 
36493639ffbSJames Wright     for (CeedInt f = 0; f < 4; f++) {
36593639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0;  // No Change in density or momentum
36693639ffbSJames Wright       v[f][i] = 0.;
36793639ffbSJames Wright     }
36893639ffbSJames Wright 
36993639ffbSJames Wright     CeedScalar div_u = 0;
37093639ffbSJames Wright     for (CeedInt j = 0; j < dim; j++) {
37193639ffbSJames Wright       for (CeedInt k = 0; k < dim; k++) {
37293639ffbSJames Wright         div_u += grad_s[k].Y.velocity[j];
37393639ffbSJames Wright       }
37493639ffbSJames Wright     }
37593639ffbSJames Wright     CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
37693639ffbSJames Wright 
37793639ffbSJames Wright     CeedScalar uX[3] = {0.};
37893639ffbSJames Wright     MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
37993639ffbSJames Wright 
38093639ffbSJames Wright     if (context->strong_form) {  // Strong Galerkin convection term: v div(E u)
38193639ffbSJames Wright       v[4][i] = -wdetJ * strong_conv;
38293639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = 0;
38393639ffbSJames Wright     } else {  // Weak Galerkin convection term: -dv \cdot (E u)
38493639ffbSJames Wright       for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = wdetJ * s.U.E_total * uX[j];
38593639ffbSJames Wright       v[4][i] = 0.;
38693639ffbSJames Wright     }
38793639ffbSJames Wright 
38893639ffbSJames Wright     const CeedScalar TauS = CtauS / sqrt(Dot3(uX, uX));
38993639ffbSJames Wright     for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) {
39093639ffbSJames Wright         case STAB_NONE:
39193639ffbSJames Wright           break;
39293639ffbSJames Wright         case STAB_SU:
39393639ffbSJames Wright         case STAB_SUPG:
39493639ffbSJames Wright           grad_v[j][4][i] += wdetJ * TauS * strong_conv * uX[j];
39593639ffbSJames Wright           break;
39693639ffbSJames Wright       }
39793639ffbSJames Wright   }
39893639ffbSJames Wright }
39993639ffbSJames Wright 
40093639ffbSJames Wright CEED_QFUNCTION(RHS_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
40193639ffbSJames Wright   RHSFunction_AdvectionGeneric(ctx, Q, in, out, 3);
40293639ffbSJames Wright   return 0;
40393639ffbSJames Wright }
40493639ffbSJames Wright 
40593639ffbSJames Wright CEED_QFUNCTION(RHS_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
40693639ffbSJames Wright   RHSFunction_AdvectionGeneric(ctx, Q, in, out, 2);
40793639ffbSJames Wright   return 0;
40893639ffbSJames Wright }
40993639ffbSJames Wright 
41093639ffbSJames Wright // *****************************************************************************
41193639ffbSJames Wright // This QFunction implements consistent outflow and inflow BCs
41293639ffbSJames Wright //      for advection
41393639ffbSJames Wright //
41493639ffbSJames Wright //  Inflow and outflow faces are determined based on sign(dot(wind, normal)):
41593639ffbSJames Wright //    sign(dot(wind, normal)) > 0 : outflow BCs
41693639ffbSJames Wright //    sign(dot(wind, normal)) < 0 : inflow BCs
41793639ffbSJames Wright //
41893639ffbSJames Wright //  Outflow BCs:
41993639ffbSJames Wright //    The validity of the weak form of the governing equations is extended to the outflow and the current values of E are applied.
42093639ffbSJames Wright //
42193639ffbSJames Wright //  Inflow BCs:
42293639ffbSJames Wright //    A prescribed Total Energy (E_wind) is applied weakly.
42393639ffbSJames Wright // *****************************************************************************
42493639ffbSJames Wright CEED_QFUNCTION(Advection_InOutFlowGeneric)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
42593639ffbSJames Wright   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
42693639ffbSJames Wright   const CeedScalar(*q_data_sur)    = in[2];
42793639ffbSJames Wright 
42893639ffbSJames Wright   CeedScalar(*v)[CEED_Q_VLA]   = (CeedScalar(*)[CEED_Q_VLA])out[0];
42993639ffbSJames Wright   AdvectionContext context     = (AdvectionContext)ctx;
43093639ffbSJames Wright   const CeedScalar E_wind      = context->E_wind;
43193639ffbSJames Wright   const CeedScalar strong_form = context->strong_form;
43293639ffbSJames Wright   const bool       is_implicit = context->implicit;
43393639ffbSJames Wright 
43493639ffbSJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
43593639ffbSJames Wright     const CeedScalar rho  = q[0][i];
43693639ffbSJames Wright     const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
43793639ffbSJames Wright     const CeedScalar E    = q[4][i];
43893639ffbSJames Wright 
43993639ffbSJames Wright     CeedScalar wdetJb, norm[3];
44093639ffbSJames Wright     QdataBoundaryUnpack_ND(dim, Q, i, q_data_sur, &wdetJb, NULL, norm);
44193639ffbSJames Wright     wdetJb *= is_implicit ? -1. : 1.;
44293639ffbSJames Wright 
44393639ffbSJames Wright     const CeedScalar u_normal = DotN(norm, u, dim);
44493639ffbSJames Wright 
44593639ffbSJames Wright     // No Change in density or momentum
44693639ffbSJames Wright     for (CeedInt j = 0; j < 4; j++) {
44793639ffbSJames Wright       v[j][i] = 0;
44893639ffbSJames Wright     }
44993639ffbSJames Wright     // Implementing in/outflow BCs
45093639ffbSJames Wright     if (u_normal > 0) {  // outflow
45193639ffbSJames Wright       v[4][i] = -(1 - strong_form) * wdetJb * E * u_normal;
45293639ffbSJames Wright     } else {  // inflow
45393639ffbSJames Wright       v[4][i] = -(1 - strong_form) * wdetJb * E_wind * u_normal;
45493639ffbSJames Wright     }
45593639ffbSJames Wright   }
45693639ffbSJames Wright   return 0;
45793639ffbSJames Wright }
45893639ffbSJames Wright 
4592b730f8bSJeremy L Thompson CEED_QFUNCTION(Advection_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
4604bdcf5bfSJames Wright   Advection_InOutFlowGeneric(ctx, Q, in, out, 3);
46177841947SLeila Ghaffari   return 0;
46277841947SLeila Ghaffari }
46377841947SLeila Ghaffari 
46493639ffbSJames Wright CEED_QFUNCTION(Advection2d_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
46593639ffbSJames Wright   Advection_InOutFlowGeneric(ctx, Q, in, out, 2);
46693639ffbSJames Wright   return 0;
46793639ffbSJames Wright }
46893639ffbSJames Wright 
46977841947SLeila Ghaffari #endif  // advection_h
470