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
5ea615d4cSJames Wright /// Advection initial condition and operator for HONEE
63e17a7a1SJames Wright #include <ceed/types.h>
7a515125bSLeila Ghaffari
8e88b842aSJames Wright #include "advection_types.h"
9ce192147SJames Wright #include "newtonian_state.h"
10ce192147SJames Wright #include "newtonian_types.h"
11e88b842aSJames Wright #include "stabilization_types.h"
121a74fa30SJames Wright #include "utils.h"
131a74fa30SJames Wright
14a515125bSLeila Ghaffari // *****************************************************************************
159529d636SJames Wright // This QFunction sets the initial conditions and the boundary conditions
169529d636SJames Wright // for two test cases: ROTATION and TRANSLATION
179529d636SJames Wright //
189529d636SJames Wright // -- ROTATION (default)
199529d636SJames Wright // Initial Conditions:
209529d636SJames Wright // Mass Density:
219529d636SJames Wright // Constant mass density of 1.0
229529d636SJames Wright // Momentum Density:
239529d636SJames Wright // Rotational field in x,y
249529d636SJames Wright // Energy Density:
259529d636SJames Wright // Maximum of 1. x0 decreasing linearly to 0. as radial distance
269529d636SJames Wright // increases to (1.-r/rc), then 0. everywhere else
279529d636SJames Wright //
289529d636SJames Wright // Boundary Conditions:
299529d636SJames Wright // Mass Density:
309529d636SJames Wright // 0.0 flux
319529d636SJames Wright // Momentum Density:
329529d636SJames Wright // 0.0
339529d636SJames Wright // Energy Density:
349529d636SJames Wright // 0.0 flux
359529d636SJames Wright //
369529d636SJames Wright // -- TRANSLATION
379529d636SJames Wright // Initial Conditions:
389529d636SJames Wright // Mass Density:
399529d636SJames Wright // Constant mass density of 1.0
409529d636SJames Wright // Momentum Density:
419529d636SJames Wright // Constant rectilinear field in x,y
429529d636SJames Wright // Energy Density:
439529d636SJames Wright // Maximum of 1. x0 decreasing linearly to 0. as radial distance
449529d636SJames Wright // increases to (1.-r/rc), then 0. everywhere else
459529d636SJames Wright //
469529d636SJames Wright // Boundary Conditions:
479529d636SJames Wright // Mass Density:
489529d636SJames Wright // 0.0 flux
499529d636SJames Wright // Momentum Density:
509529d636SJames Wright // 0.0
519529d636SJames Wright // Energy Density:
529529d636SJames Wright // Inflow BCs:
539529d636SJames Wright // E = E_wind
549529d636SJames Wright // Outflow BCs:
559529d636SJames Wright // E = E(boundary)
569529d636SJames Wright // Both In/Outflow BCs for E are applied weakly in the
579529d636SJames Wright // QFunction "Advection2d_Sur"
589529d636SJames Wright //
599529d636SJames Wright // *****************************************************************************
609529d636SJames Wright
619529d636SJames Wright // *****************************************************************************
629529d636SJames Wright // This helper function provides the exact, time-dependent solution and IC formulation for 2D advection
639529d636SJames Wright // *****************************************************************************
Exact_AdvectionGeneric(CeedInt dim,CeedScalar time,const CeedScalar X[],CeedInt Nf,CeedScalar q[],void * ctx)6497cfd714SJames Wright CEED_QFUNCTION_HELPER int Exact_AdvectionGeneric(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
659529d636SJames Wright const SetupContextAdv context = (SetupContextAdv)ctx;
669529d636SJames Wright const CeedScalar rc = context->rc;
679529d636SJames Wright const CeedScalar lx = context->lx;
689529d636SJames Wright const CeedScalar ly = context->ly;
699529d636SJames Wright const CeedScalar lz = dim == 2 ? 0. : context->lz;
709529d636SJames Wright const CeedScalar *wind = context->wind;
719529d636SJames Wright
729529d636SJames Wright const CeedScalar center[3] = {0.5 * lx, 0.5 * ly, 0.5 * lz};
739529d636SJames Wright const CeedScalar theta = dim == 2 ? M_PI / 3 : M_PI;
749529d636SJames Wright const CeedScalar x0[3] = {center[0] + .25 * lx * cos(theta + time), center[1] + .25 * ly * sin(theta + time), 0.5 * lz};
759529d636SJames Wright
769529d636SJames Wright const CeedScalar x = X[0], y = X[1], z = dim == 2 ? 0. : X[2];
779529d636SJames Wright
789529d636SJames Wright switch (context->wind_type) {
795f636aeaSJames Wright case ADVDIF_WIND_ROTATION:
809529d636SJames Wright q[0] = 1.;
819529d636SJames Wright q[1] = -(y - center[1]);
829529d636SJames Wright q[2] = (x - center[0]);
839529d636SJames Wright q[3] = 0;
849529d636SJames Wright break;
855f636aeaSJames Wright case ADVDIF_WIND_TRANSLATION:
869529d636SJames Wright q[0] = 1.;
879529d636SJames Wright q[1] = wind[0];
889529d636SJames Wright q[2] = wind[1];
899529d636SJames Wright q[3] = dim == 2 ? 0. : wind[2];
909529d636SJames Wright break;
913d1afcc1SJames Wright case ADVDIF_WIND_BOUNDARY_LAYER:
923d1afcc1SJames Wright q[0] = 1.;
933d1afcc1SJames Wright q[1] = y / ly;
943d1afcc1SJames Wright q[2] = 0.;
953d1afcc1SJames Wright q[3] = 0.;
963d1afcc1SJames Wright break;
979529d636SJames Wright }
989529d636SJames Wright
999529d636SJames Wright switch (context->initial_condition_type) {
1005f636aeaSJames Wright case ADVDIF_IC_BUBBLE_SPHERE:
1015f636aeaSJames Wright case ADVDIF_IC_BUBBLE_CYLINDER: {
102a62be6baSJames Wright CeedScalar r = sqrt(Square(x - x0[0]) + Square(y - x0[1]) + Square(z - x0[2]));
1039529d636SJames Wright switch (context->bubble_continuity_type) {
1049529d636SJames Wright // original continuous, smooth shape
1055f636aeaSJames Wright case ADVDIF_BUBBLE_CONTINUITY_SMOOTH:
1069529d636SJames Wright q[4] = r <= rc ? (1. - r / rc) : 0.;
1079529d636SJames Wright break;
1089529d636SJames Wright // discontinuous, sharp back half shape
1095f636aeaSJames Wright case ADVDIF_BUBBLE_CONTINUITY_BACK_SHARP:
1109529d636SJames Wright q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) : 0.;
1119529d636SJames Wright break;
1129529d636SJames Wright // attempt to define a finite thickness that will get resolved under grid refinement
1135f636aeaSJames Wright case ADVDIF_BUBBLE_CONTINUITY_THICK:
1149529d636SJames Wright q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) * fmin(1.0, (center[1] - y) / 1.25) : 0.;
1159529d636SJames Wright break;
1165f636aeaSJames Wright case ADVDIF_BUBBLE_CONTINUITY_COSINE:
1179529d636SJames Wright q[4] = r <= rc ? .5 + .5 * cos(r * M_PI / rc) : 0;
1189529d636SJames Wright break;
1199529d636SJames Wright }
1209529d636SJames Wright break;
121a62be6baSJames Wright }
122a62be6baSJames Wright
1235f636aeaSJames Wright case ADVDIF_IC_COSINE_HILL: {
124a62be6baSJames Wright CeedScalar r = sqrt(Square(x - center[0]) + Square(y - center[1]));
1259529d636SJames Wright CeedScalar half_width = context->lx / 2;
1269529d636SJames Wright q[4] = r > half_width ? 0. : cos(2 * M_PI * r / half_width + M_PI) + 1.;
1279529d636SJames Wright } break;
128a62be6baSJames Wright
1295f636aeaSJames Wright case ADVDIF_IC_SKEW: {
1309529d636SJames Wright CeedScalar skewed_barrier[3] = {wind[0], wind[1], 0};
1319529d636SJames Wright CeedScalar inflow_to_point[3] = {x - context->lx / 2, y, 0};
1329529d636SJames Wright CeedScalar cross_product[3] = {0};
1339529d636SJames Wright const CeedScalar boundary_threshold = 20 * CEED_EPSILON;
1349529d636SJames Wright Cross3(skewed_barrier, inflow_to_point, cross_product);
1359529d636SJames Wright
1369529d636SJames Wright q[4] = cross_product[2] > boundary_threshold ? 0 : 1;
1379529d636SJames Wright if ((x < boundary_threshold && wind[0] < boundary_threshold) || // outflow at -x boundary
1389529d636SJames Wright (y < boundary_threshold && wind[1] < boundary_threshold) || // outflow at -y boundary
1399529d636SJames Wright (x > context->lx - boundary_threshold && wind[0] > boundary_threshold) || // outflow at +x boundary
1409529d636SJames Wright (y > context->ly - boundary_threshold && wind[1] > boundary_threshold) // outflow at +y boundary
1419529d636SJames Wright ) {
1429529d636SJames Wright q[4] = 0;
1439529d636SJames Wright }
1449529d636SJames Wright } break;
145a62be6baSJames Wright
1465f636aeaSJames Wright case ADVDIF_IC_WAVE: {
147a62be6baSJames Wright CeedScalar theta = context->wave_frequency * DotN(X, wind, dim) + context->wave_phase;
148a62be6baSJames Wright switch (context->wave_type) {
149a62be6baSJames Wright case ADVDIF_WAVE_SINE:
150a62be6baSJames Wright q[4] = sin(theta);
151a62be6baSJames Wright break;
152a62be6baSJames Wright case ADVDIF_WAVE_SQUARE:
153a62be6baSJames Wright q[4] = sin(theta) > 100 * CEED_EPSILON ? 1 : -1;
154a62be6baSJames Wright break;
155a62be6baSJames Wright }
1563d1afcc1SJames Wright } break;
1570fb1909eSJames Wright
1583d1afcc1SJames Wright case ADVDIF_IC_BOUNDARY_LAYER: {
1593d1afcc1SJames Wright const CeedScalar boundary_threshold = 20 * CEED_EPSILON;
1603d1afcc1SJames Wright
1613d1afcc1SJames Wright if ((x < boundary_threshold) || (y > ly - boundary_threshold)) {
1623d1afcc1SJames Wright q[4] = 1; // inflow and top boundary
1633d1afcc1SJames Wright } else if (y < boundary_threshold) {
1643d1afcc1SJames Wright q[4] = 0; // lower wall
165b4fd18dfSJames Wright } else { // interior and outflow boundary
166b4fd18dfSJames Wright CeedScalar bl_height = ly * context->bl_height_factor;
167b4fd18dfSJames Wright if (y > bl_height) q[4] = 1;
168b4fd18dfSJames Wright else q[4] = y / bl_height;
169a62be6baSJames Wright }
1703d1afcc1SJames Wright } break;
1719529d636SJames Wright }
1729529d636SJames Wright return 0;
1739529d636SJames Wright }
1749529d636SJames Wright
1759529d636SJames Wright // *****************************************************************************
176a515125bSLeila Ghaffari // This QFunction sets the initial conditions for 3D advection
177a515125bSLeila Ghaffari // *****************************************************************************
ICsAdvection(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)1782b916ea7SJeremy L Thompson CEED_QFUNCTION(ICsAdvection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
179a515125bSLeila Ghaffari const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
180a515125bSLeila Ghaffari CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
181a515125bSLeila Ghaffari
1823d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
183a515125bSLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
184139613f2SLeila Ghaffari CeedScalar q[5] = {0.};
185a515125bSLeila Ghaffari
1860b3a1fabSJames Wright Exact_AdvectionGeneric(3, 0., x, 5, q, ctx);
187a515125bSLeila Ghaffari for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
1880b3a1fabSJames Wright }
189a515125bSLeila Ghaffari return 0;
190a515125bSLeila Ghaffari }
191a515125bSLeila Ghaffari
192a515125bSLeila Ghaffari // *****************************************************************************
1939529d636SJames Wright // This QFunction sets the initial conditions for 2D advection
194a515125bSLeila Ghaffari // *****************************************************************************
ICsAdvection2d(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)1959529d636SJames Wright CEED_QFUNCTION(ICsAdvection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
1969529d636SJames Wright const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
1979529d636SJames Wright CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
1989529d636SJames Wright const SetupContextAdv context = (SetupContextAdv)ctx;
1999529d636SJames Wright
2009529d636SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
2019529d636SJames Wright const CeedScalar x[] = {X[0][i], X[1][i]};
2029529d636SJames Wright CeedScalar q[5] = {0.};
2039529d636SJames Wright
2049529d636SJames Wright Exact_AdvectionGeneric(2, context->time, x, 5, q, ctx);
2059529d636SJames Wright for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
2069529d636SJames Wright }
207a515125bSLeila Ghaffari return 0;
208a515125bSLeila Ghaffari }
209a515125bSLeila Ghaffari
StatePhysicalGradientFromReference_ND(CeedInt N,CeedInt Q,CeedInt i,NewtonianIGProperties gas,State s,StateVariable state_var,const CeedScalar * grad_q,const CeedScalar * dXdx,State * grad_s)210cde3d787SJames Wright CEED_QFUNCTION_HELPER void StatePhysicalGradientFromReference_ND(CeedInt N, CeedInt Q, CeedInt i, NewtonianIGProperties gas, State s,
2119529d636SJames Wright StateVariable state_var, const CeedScalar *grad_q, const CeedScalar *dXdx,
2129529d636SJames Wright State *grad_s) {
2139529d636SJames Wright switch (N) {
2149529d636SJames Wright case 2: {
2159529d636SJames Wright for (CeedInt k = 0; k < 2; k++) {
2169529d636SJames Wright CeedScalar dqi[5];
2179529d636SJames Wright for (CeedInt j = 0; j < 5; j++) {
2189529d636SJames 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];
2199529d636SJames Wright }
2209529d636SJames Wright grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var);
2219529d636SJames Wright }
2229529d636SJames Wright CeedScalar U[5] = {0.};
2239529d636SJames Wright grad_s[2] = StateFromU(gas, U);
2249529d636SJames Wright } break;
2259529d636SJames Wright case 3:
22685efd435SJames Wright // Cannot directly use StatePhysicalGradientFromReference helper functions due to SYCL online compiler incompatabilities
22785efd435SJames Wright for (CeedInt k = 0; k < 3; k++) {
22885efd435SJames Wright CeedScalar dqi[5];
22985efd435SJames Wright for (CeedInt j = 0; j < 5; j++) {
23085efd435SJames 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] +
23185efd435SJames Wright grad_q[(Q * 5) * 2 + Q * j + i] * dXdx[2 * N + k];
23285efd435SJames Wright }
23385efd435SJames Wright grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var);
23485efd435SJames Wright }
2359529d636SJames Wright break;
2369529d636SJames Wright }
2379529d636SJames Wright }
2389529d636SJames Wright
239a78efa86SJames Wright // @brief Calculate the stabilization constant \tau
Tau(AdvectionContext context,const State s,const CeedScalar * dXdx,CeedInt dim)240a78efa86SJames Wright CEED_QFUNCTION_HELPER CeedScalar Tau(AdvectionContext context, const State s, const CeedScalar *dXdx, CeedInt dim) {
241a78efa86SJames Wright switch (context->stabilization_tau) {
242a78efa86SJames Wright case STAB_TAU_CTAU: {
243a78efa86SJames Wright CeedScalar uX[3] = {0.};
244a78efa86SJames Wright
245a78efa86SJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
246a78efa86SJames Wright return context->CtauS / sqrt(DotN(uX, uX, dim));
247a78efa86SJames Wright } break;
248a78efa86SJames Wright case STAB_TAU_ADVDIFF_SHAKIB: {
249a78efa86SJames Wright CeedScalar gijd_mat[9] = {0.}, gij_uj[3] = {0.};
250a78efa86SJames Wright
251a78efa86SJames Wright MatMatN(dXdx, dXdx, dim, CEED_TRANSPOSE, CEED_NOTRANSPOSE, gijd_mat);
2524ca5135bSJames Wright // (1/2)^2 to account for reference element size; for length 1 square/cube element, gij should be identity matrix
2534ca5135bSJames Wright ScaleN(gijd_mat, 0.25, Square(dim));
254a78efa86SJames Wright MatVecNM(gijd_mat, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, gij_uj);
2554ca5135bSJames Wright return 1 / sqrt(Square(2 * context->Ctau_t / context->dt) + DotN(s.Y.velocity, gij_uj, dim) * Square(context->Ctau_a) +
2564ca5135bSJames Wright Square(context->diffusion_coeff) * DotN(gijd_mat, gijd_mat, dim * dim) * Square(context->Ctau_d));
257a78efa86SJames Wright } break;
258a78efa86SJames Wright default:
259a78efa86SJames Wright return 0.;
260a78efa86SJames Wright }
261a78efa86SJames Wright }
262a78efa86SJames Wright
2639529d636SJames Wright // *****************************************************************************
2649529d636SJames Wright // This QFunction implements Advection for implicit time stepping method
2659529d636SJames Wright // *****************************************************************************
IFunction_AdvectionGeneric(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out,CeedInt dim)26697cfd714SJames Wright CEED_QFUNCTION_HELPER int IFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
2674c5ab12fSJames Wright AdvectionContext context = (AdvectionContext)ctx;
2684c5ab12fSJames Wright
2699529d636SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
2709529d636SJames Wright const CeedScalar(*grad_q) = in[1];
2719529d636SJames Wright const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
2729529d636SJames Wright const CeedScalar(*q_data) = in[3];
2734c5ab12fSJames Wright const CeedScalar(*divFdiff) = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? in[5] : NULL;
2749529d636SJames Wright
2759529d636SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
2769529d636SJames Wright CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
2779529d636SJames Wright
278cde3d787SJames Wright NewtonianIGProperties gas = {0};
2799529d636SJames Wright
2809529d636SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
2819529d636SJames Wright const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
2829529d636SJames Wright const State s = StateFromU(gas, qi);
2839529d636SJames Wright
2849529d636SJames Wright CeedScalar wdetJ, dXdx[9];
2859529d636SJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
2869529d636SJames Wright State grad_s[3];
2879529d636SJames Wright StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
2889529d636SJames Wright
2899529d636SJames Wright const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
2909529d636SJames Wright
2919529d636SJames Wright for (CeedInt f = 0; f < 4; f++) {
2929529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum
2939529d636SJames Wright v[f][i] = wdetJ * q_dot[f][i]; // K Mass/transient term
2949529d636SJames Wright }
2959529d636SJames Wright
2969529d636SJames Wright CeedScalar div_u = 0;
2979529d636SJames Wright for (CeedInt j = 0; j < dim; j++) {
2989529d636SJames Wright for (CeedInt k = 0; k < dim; k++) {
2999529d636SJames Wright div_u += grad_s[k].Y.velocity[j];
3009529d636SJames Wright }
3019529d636SJames Wright }
3029529d636SJames Wright CeedScalar uX[3] = {0.};
3039529d636SJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
3044c5ab12fSJames Wright CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
3059529d636SJames Wright
3064c5ab12fSJames Wright v[4][i] = wdetJ * q_dot[4][i]; // transient part (ALWAYS)
3079529d636SJames Wright if (context->strong_form) { // Strong Galerkin convection term: v div(E u)
3089529d636SJames Wright v[4][i] += wdetJ * strong_conv;
3099529d636SJames Wright } else { // Weak Galerkin convection term: -dv \cdot (E u)
3109529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = -wdetJ * s.U.E_total * uX[j];
3119529d636SJames Wright }
3129529d636SJames Wright
313c8d249deSJames Wright { // Diffusion
314c8d249deSJames Wright CeedScalar Fe[3], Fe_dXdx[3] = {0.};
315c8d249deSJames Wright
316c8d249deSJames Wright for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total;
317c8d249deSJames Wright MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx);
318c8d249deSJames Wright for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] -= wdetJ * Fe_dXdx[k];
319c8d249deSJames Wright }
320c8d249deSJames Wright
321a78efa86SJames Wright const CeedScalar TauS = Tau(context, s, dXdx, dim);
3224c5ab12fSJames Wright for (CeedInt j = 0; j < dim; j++) {
3234c5ab12fSJames Wright switch (context->stabilization) {
3249529d636SJames Wright case STAB_NONE:
3259529d636SJames Wright break;
3269529d636SJames Wright case STAB_SU:
3274c5ab12fSJames Wright grad_v[j][4][i] += wdetJ * TauS * uX[j] * strong_conv;
3289529d636SJames Wright break;
3294c5ab12fSJames Wright case STAB_SUPG: {
3304c5ab12fSJames Wright CeedScalar divFdiff_i = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? divFdiff[i] : 0.;
3314c5ab12fSJames Wright grad_v[j][4][i] += wdetJ * TauS * uX[j] * (q_dot[4][i] + strong_conv + divFdiff_i);
3324c5ab12fSJames Wright } break;
3334c5ab12fSJames Wright }
3349529d636SJames Wright }
3359529d636SJames Wright }
33697cfd714SJames Wright return 0;
3379529d636SJames Wright }
3389529d636SJames Wright
IFunction_Advection(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)3392b916ea7SJeremy L Thompson CEED_QFUNCTION(IFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
34097cfd714SJames Wright return IFunction_AdvectionGeneric(ctx, Q, in, out, 3);
341a515125bSLeila Ghaffari }
342a515125bSLeila Ghaffari
IFunction_Advection2d(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)3439529d636SJames Wright CEED_QFUNCTION(IFunction_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
34497cfd714SJames Wright return IFunction_AdvectionGeneric(ctx, Q, in, out, 2);
3459529d636SJames Wright }
3469529d636SJames Wright
MassFunction_AdvectionGeneric(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out,CeedInt dim)34797cfd714SJames Wright CEED_QFUNCTION_HELPER int MassFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
348a78efa86SJames Wright const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
349a78efa86SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1];
350a78efa86SJames Wright const CeedScalar(*q_data) = in[2];
351a78efa86SJames Wright
352a78efa86SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
353a78efa86SJames Wright CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
354a78efa86SJames Wright
355a78efa86SJames Wright AdvectionContext context = (AdvectionContext)ctx;
356cde3d787SJames Wright NewtonianIGProperties gas = {0};
357a78efa86SJames Wright
358a78efa86SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
359a78efa86SJames Wright const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
360a78efa86SJames Wright const State s = StateFromU(gas, qi);
361a78efa86SJames Wright CeedScalar wdetJ, dXdx[9];
362a78efa86SJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
363a78efa86SJames Wright
364a78efa86SJames Wright for (CeedInt f = 0; f < 4; f++) {
365a78efa86SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum
366a78efa86SJames Wright v[f][i] = wdetJ * q_dot[f][i]; // K Mass/transient term
367a78efa86SJames Wright }
368a78efa86SJames Wright
369a78efa86SJames Wright // Unstabilized mass term
370a78efa86SJames Wright v[4][i] = wdetJ * q_dot[4][i];
371a78efa86SJames Wright
372a78efa86SJames Wright // Stabilized mass term
373a78efa86SJames Wright CeedScalar uX[3] = {0.};
374a78efa86SJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
375a78efa86SJames Wright const CeedScalar TauS = Tau(context, s, dXdx, dim);
37671acc5eeSJames Wright for (CeedInt j = 0; j < dim; j++) {
37771acc5eeSJames Wright switch (context->stabilization) {
378a78efa86SJames Wright case STAB_NONE:
379a78efa86SJames Wright case STAB_SU:
380a78efa86SJames Wright grad_v[j][4][i] = 0;
381a78efa86SJames Wright break; // These should be run with the unstabilized mass matrix anyways
382a78efa86SJames Wright case STAB_SUPG:
383a78efa86SJames Wright grad_v[j][4][i] = wdetJ * TauS * q_dot[4][i] * uX[j];
384a78efa86SJames Wright break;
385a78efa86SJames Wright }
386a78efa86SJames Wright }
38771acc5eeSJames Wright }
38897cfd714SJames Wright return 0;
389a78efa86SJames Wright }
390a78efa86SJames Wright
MassFunction_Advection(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)391a78efa86SJames Wright CEED_QFUNCTION(MassFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
39297cfd714SJames Wright return MassFunction_AdvectionGeneric(ctx, Q, in, out, 3);
393a78efa86SJames Wright }
394a78efa86SJames Wright
MassFunction_Advection2D(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)395a78efa86SJames Wright CEED_QFUNCTION(MassFunction_Advection2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
39697cfd714SJames Wright return MassFunction_AdvectionGeneric(ctx, Q, in, out, 2);
397a78efa86SJames Wright }
398a78efa86SJames Wright
3999529d636SJames Wright // *****************************************************************************
4009529d636SJames Wright // This QFunction implements Advection for explicit time stepping method
4019529d636SJames Wright // *****************************************************************************
RHSFunction_AdvectionGeneric(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out,CeedInt dim)40297cfd714SJames Wright CEED_QFUNCTION_HELPER int RHSFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
4035f952e8dSJames Wright AdvectionContext context = (AdvectionContext)ctx;
4045f952e8dSJames Wright
4059529d636SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
4069529d636SJames Wright const CeedScalar(*grad_q) = in[1];
4079529d636SJames Wright const CeedScalar(*q_data) = in[2];
4085f952e8dSJames Wright const CeedScalar(*divFdiff) = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? in[4] : NULL;
4099529d636SJames Wright
4109529d636SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
4119529d636SJames Wright CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
4129529d636SJames Wright
413cde3d787SJames Wright NewtonianIGProperties gas = {0};
4149529d636SJames Wright
4159529d636SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
4169529d636SJames Wright const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
4179529d636SJames Wright const State s = StateFromU(gas, qi);
4189529d636SJames Wright
4199529d636SJames Wright CeedScalar wdetJ, dXdx[9];
4209529d636SJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
4219529d636SJames Wright State grad_s[3];
4229529d636SJames Wright StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
4239529d636SJames Wright
4249529d636SJames Wright const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
4259529d636SJames Wright
4269529d636SJames Wright for (CeedInt f = 0; f < 4; f++) {
4279529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum
4289529d636SJames Wright v[f][i] = 0.;
4299529d636SJames Wright }
4309529d636SJames Wright
4319529d636SJames Wright CeedScalar div_u = 0;
4329529d636SJames Wright for (CeedInt j = 0; j < dim; j++) {
4339529d636SJames Wright for (CeedInt k = 0; k < dim; k++) {
4349529d636SJames Wright div_u += grad_s[k].Y.velocity[j];
4359529d636SJames Wright }
4369529d636SJames Wright }
4379529d636SJames Wright CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
4389529d636SJames Wright
4399529d636SJames Wright CeedScalar uX[3] = {0.};
4409529d636SJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
4419529d636SJames Wright
4429529d636SJames Wright if (context->strong_form) { // Strong Galerkin convection term: v div(E u)
4439529d636SJames Wright v[4][i] = -wdetJ * strong_conv;
4449529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = 0;
4459529d636SJames Wright } else { // Weak Galerkin convection term: -dv \cdot (E u)
4469529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = wdetJ * s.U.E_total * uX[j];
4479529d636SJames Wright v[4][i] = 0.;
4489529d636SJames Wright }
4499529d636SJames Wright
450c8d249deSJames Wright { // Diffusion
451c8d249deSJames Wright CeedScalar Fe[3], Fe_dXdx[3] = {0.};
452c8d249deSJames Wright
453c8d249deSJames Wright for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total;
454c8d249deSJames Wright MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx);
455c8d249deSJames Wright for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] += wdetJ * Fe_dXdx[k];
456c8d249deSJames Wright }
457c8d249deSJames Wright
458a78efa86SJames Wright const CeedScalar TauS = Tau(context, s, dXdx, dim);
4595f952e8dSJames Wright for (CeedInt j = 0; j < dim; j++) {
4605f952e8dSJames Wright switch (context->stabilization) {
4619529d636SJames Wright case STAB_NONE:
4629529d636SJames Wright break;
4639529d636SJames Wright case STAB_SU:
4645f952e8dSJames Wright case STAB_SUPG: {
4655f952e8dSJames Wright CeedScalar divFdiff_i = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? divFdiff[i] : 0.;
4665f952e8dSJames Wright grad_v[j][4][i] -= wdetJ * TauS * (strong_conv + divFdiff_i) * uX[j];
4675f952e8dSJames Wright } break;
4685f952e8dSJames Wright }
4699529d636SJames Wright }
4709529d636SJames Wright }
47197cfd714SJames Wright return 0;
4729529d636SJames Wright }
4739529d636SJames Wright
RHS_Advection(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)4749529d636SJames Wright CEED_QFUNCTION(RHS_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
47597cfd714SJames Wright return RHSFunction_AdvectionGeneric(ctx, Q, in, out, 3);
4769529d636SJames Wright }
4779529d636SJames Wright
RHS_Advection2d(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)4789529d636SJames Wright CEED_QFUNCTION(RHS_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
47997cfd714SJames Wright return RHSFunction_AdvectionGeneric(ctx, Q, in, out, 2);
4809529d636SJames Wright }
4819529d636SJames Wright
4829529d636SJames Wright // *****************************************************************************
4839529d636SJames Wright // This QFunction implements consistent outflow and inflow BCs
4849529d636SJames Wright // for advection
4859529d636SJames Wright //
4869529d636SJames Wright // Inflow and outflow faces are determined based on sign(dot(wind, normal)):
4879529d636SJames Wright // sign(dot(wind, normal)) > 0 : outflow BCs
4889529d636SJames Wright // sign(dot(wind, normal)) < 0 : inflow BCs
4899529d636SJames Wright //
4909529d636SJames Wright // Outflow BCs:
4919529d636SJames Wright // The validity of the weak form of the governing equations is extended to the outflow and the current values of E are applied.
4929529d636SJames Wright //
4939529d636SJames Wright // Inflow BCs:
4949529d636SJames Wright // A prescribed Total Energy (E_wind) is applied weakly.
4959529d636SJames Wright // *****************************************************************************
Advection_InOutFlowGeneric(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out,CeedInt dim)49697cfd714SJames Wright CEED_QFUNCTION_HELPER int Advection_InOutFlowGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
4979529d636SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
4989529d636SJames Wright const CeedScalar(*q_data_sur) = in[2];
4999529d636SJames Wright
5009529d636SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
5019529d636SJames Wright AdvectionContext context = (AdvectionContext)ctx;
5029529d636SJames Wright const CeedScalar E_wind = context->E_wind;
5039529d636SJames Wright const CeedScalar strong_form = context->strong_form;
5049529d636SJames Wright const bool is_implicit = context->implicit;
5059529d636SJames Wright
5069529d636SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
5079529d636SJames Wright const CeedScalar rho = q[0][i];
5089529d636SJames Wright const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
5099529d636SJames Wright const CeedScalar E = q[4][i];
5109529d636SJames Wright
51178e8b7daSJames Wright CeedScalar wdetJb, normal[3];
51278e8b7daSJames Wright QdataBoundaryUnpack_ND(dim, Q, i, q_data_sur, &wdetJb, NULL, normal);
5139529d636SJames Wright wdetJb *= is_implicit ? -1. : 1.;
5149529d636SJames Wright
51578e8b7daSJames Wright const CeedScalar u_normal = DotN(normal, u, dim);
5169529d636SJames Wright
5179529d636SJames Wright // No Change in density or momentum
5189529d636SJames Wright for (CeedInt j = 0; j < 4; j++) {
5199529d636SJames Wright v[j][i] = 0;
5209529d636SJames Wright }
5219529d636SJames Wright // Implementing in/outflow BCs
5229529d636SJames Wright if (u_normal > 0) { // outflow
5239529d636SJames Wright v[4][i] = -(1 - strong_form) * wdetJb * E * u_normal;
5249529d636SJames Wright } else { // inflow
5259529d636SJames Wright v[4][i] = -(1 - strong_form) * wdetJb * E_wind * u_normal;
5269529d636SJames Wright }
5279529d636SJames Wright }
5289529d636SJames Wright return 0;
5299529d636SJames Wright }
5309529d636SJames Wright
Advection_InOutFlow(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)5312b916ea7SJeremy L Thompson CEED_QFUNCTION(Advection_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
53297cfd714SJames Wright return Advection_InOutFlowGeneric(ctx, Q, in, out, 3);
533a515125bSLeila Ghaffari }
534a515125bSLeila Ghaffari
Advection2d_InOutFlow(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)5359529d636SJames Wright CEED_QFUNCTION(Advection2d_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
53697cfd714SJames Wright return Advection_InOutFlowGeneric(ctx, Q, in, out, 2);
5379529d636SJames Wright }
538c2d90829SJames Wright
539c2d90829SJames Wright // @brief Volume integral for RHS of divergence of diffusive flux direct projection
DivDiffusiveFluxVolumeRHS_AdvDif_Generic(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out,const CeedInt dim)540c2d90829SJames Wright CEED_QFUNCTION_HELPER int DivDiffusiveFluxVolumeRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out,
541c2d90829SJames Wright const CeedInt dim) {
542c2d90829SJames Wright const CeedScalar(*Grad_q) = in[0];
543c2d90829SJames Wright const CeedScalar(*q_data) = in[1];
544c2d90829SJames Wright CeedScalar(*Grad_v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
545c2d90829SJames Wright
546c2d90829SJames Wright AdvectionContext context = (AdvectionContext)ctx;
547c2d90829SJames Wright
548c2d90829SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
549c2d90829SJames Wright CeedScalar wdetJ, dXdx[9], F_diff[3] = {0.};
550c2d90829SJames Wright
551c2d90829SJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
552c2d90829SJames Wright { // Get physical diffusive flux
553c2d90829SJames Wright CeedScalar Grad_qn[15], grad_E_ref[3];
554c2d90829SJames Wright
555*22440147SJames Wright GradUnpackND(Q, i, 5, dim, Grad_q, Grad_qn);
556c2d90829SJames Wright CopyN(&Grad_qn[4 * dim], grad_E_ref, dim);
557c2d90829SJames Wright MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff);
558c2d90829SJames Wright ScaleN(F_diff, -context->diffusion_coeff, dim);
559c2d90829SJames Wright }
560c2d90829SJames Wright
561c2d90829SJames Wright CeedScalar F_diff_dXdx[3] = {0.};
562c2d90829SJames Wright MatVecNM(dXdx, F_diff, dim, dim, CEED_NOTRANSPOSE, F_diff_dXdx);
563c2d90829SJames Wright for (CeedInt k = 0; k < dim; k++) Grad_v[k][i] = -wdetJ * F_diff_dXdx[k];
564c2d90829SJames Wright }
565c2d90829SJames Wright return 0;
566c2d90829SJames Wright }
567c2d90829SJames Wright
DivDiffusiveFluxVolumeRHS_AdvDif_2D(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)568c2d90829SJames Wright CEED_QFUNCTION(DivDiffusiveFluxVolumeRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
569c2d90829SJames Wright return DivDiffusiveFluxVolumeRHS_AdvDif_Generic(ctx, Q, in, out, 2);
570c2d90829SJames Wright }
571c2d90829SJames Wright
DivDiffusiveFluxVolumeRHS_AdvDif_3D(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)572c2d90829SJames Wright CEED_QFUNCTION(DivDiffusiveFluxVolumeRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
573c2d90829SJames Wright return DivDiffusiveFluxVolumeRHS_AdvDif_Generic(ctx, Q, in, out, 3);
574c2d90829SJames Wright }
575c2d90829SJames Wright
576c2d90829SJames Wright // @brief Boundary integral for RHS of divergence of diffusive flux direct projection
DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out,const CeedInt dim)577c2d90829SJames Wright CEED_QFUNCTION_HELPER int DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out,
578c2d90829SJames Wright const CeedInt dim) {
579c2d90829SJames Wright const CeedScalar(*Grad_q) = in[0];
580c2d90829SJames Wright const CeedScalar(*q_data) = in[1];
581c2d90829SJames Wright CeedScalar(*v) = out[0];
582c2d90829SJames Wright
583c2d90829SJames Wright AdvectionContext context = (AdvectionContext)ctx;
584c2d90829SJames Wright
585c2d90829SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
586c2d90829SJames Wright CeedScalar wdetJ, normal[3], dXdx[9], F_diff[3] = {0.};
587c2d90829SJames Wright
588c2d90829SJames Wright QdataBoundaryGradientUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx, normal);
589c2d90829SJames Wright { // Get physical diffusive flux
590c2d90829SJames Wright CeedScalar Grad_qn[15], grad_E_ref[3];
591c2d90829SJames Wright
592*22440147SJames Wright GradUnpackND(Q, i, 5, dim, Grad_q, Grad_qn);
593c2d90829SJames Wright CopyN(&Grad_qn[4 * dim], grad_E_ref, dim);
594c2d90829SJames Wright MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff);
595c2d90829SJames Wright ScaleN(F_diff, -context->diffusion_coeff, dim);
596c2d90829SJames Wright }
597c2d90829SJames Wright
598c2d90829SJames Wright v[i] = wdetJ * DotN(F_diff, normal, dim);
599c2d90829SJames Wright }
600c2d90829SJames Wright return 0;
601c2d90829SJames Wright }
602c2d90829SJames Wright
DivDiffusiveFluxBoundaryRHS_AdvDif_2D(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)603c2d90829SJames Wright CEED_QFUNCTION(DivDiffusiveFluxBoundaryRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
604c2d90829SJames Wright return DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(ctx, Q, in, out, 2);
605c2d90829SJames Wright }
606c2d90829SJames Wright
DivDiffusiveFluxBoundaryRHS_AdvDif_3D(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)607c2d90829SJames Wright CEED_QFUNCTION(DivDiffusiveFluxBoundaryRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
608c2d90829SJames Wright return DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(ctx, Q, in, out, 3);
609c2d90829SJames Wright }
61040b78511SJames Wright
61140b78511SJames Wright // @brief Volume integral for RHS of diffusive flux indirect projection
DiffusiveFluxRHS_AdvDif_Generic(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out,const CeedInt dim)61240b78511SJames Wright CEED_QFUNCTION_HELPER int DiffusiveFluxRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out,
61340b78511SJames Wright const CeedInt dim) {
61440b78511SJames Wright const CeedScalar(*Grad_q) = in[0];
61540b78511SJames Wright const CeedScalar(*q_data) = in[1];
61640b78511SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
61740b78511SJames Wright
61840b78511SJames Wright AdvectionContext context = (AdvectionContext)ctx;
61940b78511SJames Wright
62040b78511SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
62140b78511SJames Wright CeedScalar wdetJ, dXdx[9], F_diff[3] = {0.};
62240b78511SJames Wright
62340b78511SJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
62440b78511SJames Wright { // Get physical diffusive flux
62540b78511SJames Wright CeedScalar Grad_qn[15], grad_E_ref[3];
62640b78511SJames Wright
627*22440147SJames Wright GradUnpackND(Q, i, 5, dim, Grad_q, Grad_qn);
62840b78511SJames Wright CopyN(&Grad_qn[4 * dim], grad_E_ref, dim);
62940b78511SJames Wright MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff);
63040b78511SJames Wright ScaleN(F_diff, -context->diffusion_coeff, dim);
63140b78511SJames Wright }
63240b78511SJames Wright for (CeedInt k = 0; k < dim; k++) v[k][i] = wdetJ * F_diff[k];
63340b78511SJames Wright }
63440b78511SJames Wright return 0;
63540b78511SJames Wright }
63640b78511SJames Wright
DiffusiveFluxRHS_AdvDif_2D(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)63740b78511SJames Wright CEED_QFUNCTION(DiffusiveFluxRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
63840b78511SJames Wright return DiffusiveFluxRHS_AdvDif_Generic(ctx, Q, in, out, 2);
63940b78511SJames Wright }
64040b78511SJames Wright
DiffusiveFluxRHS_AdvDif_3D(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)64140b78511SJames Wright CEED_QFUNCTION(DiffusiveFluxRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
64240b78511SJames Wright return DiffusiveFluxRHS_AdvDif_Generic(ctx, Q, in, out, 3);
64340b78511SJames Wright }
644