1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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
10c0b5abf0SJeremy L Thompson #include <ceed/types.h>
11c0b5abf0SJeremy L Thompson #ifndef CEED_RUNNING_JIT_PASS
12c9c2c079SJeremy L Thompson #include <math.h>
13c0b5abf0SJeremy L Thompson #include <stdbool.h>
14c0b5abf0SJeremy L Thompson #endif
1577841947SLeila Ghaffari
16c44b1c7dSJames Wright #include "advection_types.h"
178f4d89c8SJames Wright #include "newtonian_state.h"
188f4d89c8SJames Wright #include "newtonian_types.h"
19c44b1c7dSJames Wright #include "stabilization_types.h"
208756a6ccSJames Wright #include "utils.h"
218756a6ccSJames Wright
2277841947SLeila Ghaffari // *****************************************************************************
2393639ffbSJames Wright // This QFunction sets the initial conditions and the boundary conditions
2493639ffbSJames Wright // for two test cases: ROTATION and TRANSLATION
2593639ffbSJames Wright //
2693639ffbSJames Wright // -- ROTATION (default)
2793639ffbSJames Wright // Initial Conditions:
2893639ffbSJames Wright // Mass Density:
2993639ffbSJames Wright // Constant mass density of 1.0
3093639ffbSJames Wright // Momentum Density:
3193639ffbSJames Wright // Rotational field in x,y
3293639ffbSJames Wright // Energy Density:
3393639ffbSJames Wright // Maximum of 1. x0 decreasing linearly to 0. as radial distance
3493639ffbSJames Wright // increases to (1.-r/rc), then 0. everywhere else
3593639ffbSJames Wright //
3693639ffbSJames Wright // Boundary Conditions:
3793639ffbSJames Wright // Mass Density:
3893639ffbSJames Wright // 0.0 flux
3993639ffbSJames Wright // Momentum Density:
4093639ffbSJames Wright // 0.0
4193639ffbSJames Wright // Energy Density:
4293639ffbSJames Wright // 0.0 flux
4393639ffbSJames Wright //
4493639ffbSJames Wright // -- TRANSLATION
4593639ffbSJames Wright // Initial Conditions:
4693639ffbSJames Wright // Mass Density:
4793639ffbSJames Wright // Constant mass density of 1.0
4893639ffbSJames Wright // Momentum Density:
4993639ffbSJames Wright // Constant rectilinear field in x,y
5093639ffbSJames Wright // Energy Density:
5193639ffbSJames Wright // Maximum of 1. x0 decreasing linearly to 0. as radial distance
5293639ffbSJames Wright // increases to (1.-r/rc), then 0. everywhere else
5393639ffbSJames Wright //
5493639ffbSJames Wright // Boundary Conditions:
5593639ffbSJames Wright // Mass Density:
5693639ffbSJames Wright // 0.0 flux
5793639ffbSJames Wright // Momentum Density:
5893639ffbSJames Wright // 0.0
5993639ffbSJames Wright // Energy Density:
6093639ffbSJames Wright // Inflow BCs:
6193639ffbSJames Wright // E = E_wind
6293639ffbSJames Wright // Outflow BCs:
6393639ffbSJames Wright // E = E(boundary)
6493639ffbSJames Wright // Both In/Outflow BCs for E are applied weakly in the
6593639ffbSJames Wright // QFunction "Advection2d_Sur"
6693639ffbSJames Wright //
6793639ffbSJames Wright // *****************************************************************************
6893639ffbSJames Wright
6993639ffbSJames Wright // *****************************************************************************
7093639ffbSJames Wright // This helper function provides the exact, time-dependent solution and IC formulation for 2D advection
7193639ffbSJames Wright // *****************************************************************************
Exact_AdvectionGeneric(CeedInt dim,CeedScalar time,const CeedScalar X[],CeedInt Nf,CeedScalar q[],void * ctx)7293639ffbSJames Wright CEED_QFUNCTION_HELPER CeedInt Exact_AdvectionGeneric(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
7393639ffbSJames Wright const SetupContextAdv context = (SetupContextAdv)ctx;
7493639ffbSJames Wright const CeedScalar rc = context->rc;
7593639ffbSJames Wright const CeedScalar lx = context->lx;
7693639ffbSJames Wright const CeedScalar ly = context->ly;
7793639ffbSJames Wright const CeedScalar lz = dim == 2 ? 0. : context->lz;
7893639ffbSJames Wright const CeedScalar *wind = context->wind;
7993639ffbSJames Wright
8093639ffbSJames Wright const CeedScalar center[3] = {0.5 * lx, 0.5 * ly, 0.5 * lz};
8193639ffbSJames Wright const CeedScalar theta = dim == 2 ? M_PI / 3 : M_PI;
8293639ffbSJames Wright const CeedScalar x0[3] = {center[0] + .25 * lx * cos(theta + time), center[1] + .25 * ly * sin(theta + time), 0.5 * lz};
8393639ffbSJames Wright
8493639ffbSJames Wright const CeedScalar x = X[0], y = X[1], z = dim == 2 ? 0. : X[2];
8593639ffbSJames Wright
8693639ffbSJames Wright CeedScalar r = 0.;
8793639ffbSJames Wright switch (context->initial_condition_type) {
8893639ffbSJames Wright case ADVECTIONIC_BUBBLE_SPHERE:
8993639ffbSJames Wright case ADVECTIONIC_BUBBLE_CYLINDER:
9093639ffbSJames Wright r = sqrt(Square(x - x0[0]) + Square(y - x0[1]) + Square(z - x0[2]));
9193639ffbSJames Wright break;
9293639ffbSJames Wright case ADVECTIONIC_COSINE_HILL:
9393639ffbSJames Wright r = sqrt(Square(x - center[0]) + Square(y - center[1]));
9493639ffbSJames Wright break;
9593639ffbSJames Wright case ADVECTIONIC_SKEW:
9693639ffbSJames Wright break;
9793639ffbSJames Wright }
9893639ffbSJames Wright
9993639ffbSJames Wright switch (context->wind_type) {
10093639ffbSJames Wright case WIND_ROTATION:
10193639ffbSJames Wright q[0] = 1.;
10293639ffbSJames Wright q[1] = -(y - center[1]);
10393639ffbSJames Wright q[2] = (x - center[0]);
10493639ffbSJames Wright q[3] = 0;
10593639ffbSJames Wright break;
10693639ffbSJames Wright case WIND_TRANSLATION:
10793639ffbSJames Wright q[0] = 1.;
10893639ffbSJames Wright q[1] = wind[0];
10993639ffbSJames Wright q[2] = wind[1];
11093639ffbSJames Wright q[3] = dim == 2 ? 0. : wind[2];
11193639ffbSJames Wright break;
11293639ffbSJames Wright default:
11393639ffbSJames Wright return 1;
11493639ffbSJames Wright }
11593639ffbSJames Wright
11693639ffbSJames Wright switch (context->initial_condition_type) {
11793639ffbSJames Wright case ADVECTIONIC_BUBBLE_SPHERE:
11893639ffbSJames Wright case ADVECTIONIC_BUBBLE_CYLINDER:
11993639ffbSJames Wright switch (context->bubble_continuity_type) {
12093639ffbSJames Wright // original continuous, smooth shape
12193639ffbSJames Wright case BUBBLE_CONTINUITY_SMOOTH:
12293639ffbSJames Wright q[4] = r <= rc ? (1. - r / rc) : 0.;
12393639ffbSJames Wright break;
12493639ffbSJames Wright // discontinuous, sharp back half shape
12593639ffbSJames Wright case BUBBLE_CONTINUITY_BACK_SHARP:
12693639ffbSJames Wright q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) : 0.;
12793639ffbSJames Wright break;
12893639ffbSJames Wright // attempt to define a finite thickness that will get resolved under grid refinement
12993639ffbSJames Wright case BUBBLE_CONTINUITY_THICK:
13093639ffbSJames Wright q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) * fmin(1.0, (center[1] - y) / 1.25) : 0.;
13193639ffbSJames Wright break;
13293639ffbSJames Wright case BUBBLE_CONTINUITY_COSINE:
13393639ffbSJames Wright q[4] = r <= rc ? .5 + .5 * cos(r * M_PI / rc) : 0;
13493639ffbSJames Wright break;
13593639ffbSJames Wright }
13693639ffbSJames Wright break;
13793639ffbSJames Wright case ADVECTIONIC_COSINE_HILL: {
13893639ffbSJames Wright CeedScalar half_width = context->lx / 2;
13993639ffbSJames Wright q[4] = r > half_width ? 0. : cos(2 * M_PI * r / half_width + M_PI) + 1.;
14093639ffbSJames Wright } break;
14193639ffbSJames Wright case ADVECTIONIC_SKEW: {
14293639ffbSJames Wright CeedScalar skewed_barrier[3] = {wind[0], wind[1], 0};
14393639ffbSJames Wright CeedScalar inflow_to_point[3] = {x - context->lx / 2, y, 0};
14493639ffbSJames Wright CeedScalar cross_product[3] = {0};
14593639ffbSJames Wright const CeedScalar boundary_threshold = 20 * CEED_EPSILON;
14693639ffbSJames Wright Cross3(skewed_barrier, inflow_to_point, cross_product);
14793639ffbSJames Wright
14893639ffbSJames Wright q[4] = cross_product[2] > boundary_threshold ? 0 : 1;
14993639ffbSJames Wright if ((x < boundary_threshold && wind[0] < boundary_threshold) || // outflow at -x boundary
15093639ffbSJames Wright (y < boundary_threshold && wind[1] < boundary_threshold) || // outflow at -y boundary
15193639ffbSJames Wright (x > context->lx - boundary_threshold && wind[0] > boundary_threshold) || // outflow at +x boundary
15293639ffbSJames Wright (y > context->ly - boundary_threshold && wind[1] > boundary_threshold) // outflow at +y boundary
15393639ffbSJames Wright ) {
15493639ffbSJames Wright q[4] = 0;
15593639ffbSJames Wright }
15693639ffbSJames Wright } break;
15793639ffbSJames Wright }
15893639ffbSJames Wright return 0;
15993639ffbSJames Wright }
16093639ffbSJames Wright
16193639ffbSJames Wright // *****************************************************************************
16277841947SLeila Ghaffari // This QFunction sets the initial conditions for 3D advection
16377841947SLeila Ghaffari // *****************************************************************************
ICsAdvection(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)1642b730f8bSJeremy L Thompson CEED_QFUNCTION(ICsAdvection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
16577841947SLeila Ghaffari const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
16677841947SLeila Ghaffari CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
16777841947SLeila Ghaffari
16846603fc5SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
16977841947SLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
170e6225c47SLeila Ghaffari CeedScalar q[5] = {0.};
17177841947SLeila Ghaffari
17230e1b2c7SJames Wright Exact_AdvectionGeneric(3, 0., x, 5, q, ctx);
17377841947SLeila Ghaffari for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
17430e1b2c7SJames Wright }
17577841947SLeila Ghaffari return 0;
17677841947SLeila Ghaffari }
17777841947SLeila Ghaffari
17877841947SLeila Ghaffari // *****************************************************************************
17993639ffbSJames Wright // This QFunction sets the initial conditions for 2D advection
18077841947SLeila Ghaffari // *****************************************************************************
ICsAdvection2d(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)18193639ffbSJames Wright CEED_QFUNCTION(ICsAdvection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
18293639ffbSJames Wright const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
18393639ffbSJames Wright CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
18493639ffbSJames Wright const SetupContextAdv context = (SetupContextAdv)ctx;
18593639ffbSJames Wright
18693639ffbSJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
18793639ffbSJames Wright const CeedScalar x[] = {X[0][i], X[1][i]};
18893639ffbSJames Wright CeedScalar q[5] = {0.};
18993639ffbSJames Wright
19093639ffbSJames Wright Exact_AdvectionGeneric(2, context->time, x, 5, q, ctx);
19193639ffbSJames Wright for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
19293639ffbSJames Wright }
19377841947SLeila Ghaffari return 0;
19477841947SLeila Ghaffari }
19577841947SLeila Ghaffari
QdataUnpack_ND(CeedInt N,CeedInt Q,CeedInt i,const CeedScalar * q_data,CeedScalar * wdetJ,CeedScalar * dXdx)19693639ffbSJames Wright CEED_QFUNCTION_HELPER void QdataUnpack_ND(CeedInt N, CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar *dXdx) {
197b1559114SJames Wright // Cannot directly use QdataUnpack* helper functions due to SYCL online compiler incompatabilities
19893639ffbSJames Wright switch (N) {
19993639ffbSJames Wright case 2:
200b1559114SJames Wright StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ);
201b1559114SJames Wright StoredValuesUnpack(Q, i, 1, 4, q_data, dXdx);
20293639ffbSJames Wright break;
20393639ffbSJames Wright case 3:
204b1559114SJames Wright StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ);
205b1559114SJames Wright StoredValuesUnpack(Q, i, 1, 9, q_data, dXdx);
20693639ffbSJames Wright break;
20793639ffbSJames Wright }
20893639ffbSJames Wright }
20993639ffbSJames Wright
QdataBoundaryUnpack_ND(CeedInt N,CeedInt Q,CeedInt i,const CeedScalar * q_data,CeedScalar * wdetJ,CeedScalar * dXdx,CeedScalar * normal)21093639ffbSJames Wright CEED_QFUNCTION_HELPER int QdataBoundaryUnpack_ND(CeedInt N, CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar *dXdx,
21193639ffbSJames Wright CeedScalar *normal) {
212b1559114SJames Wright // Cannot directly use QdataBoundaryUnpack* helper functions due to SYCL online compiler incompatabilities
21393639ffbSJames Wright switch (N) {
21493639ffbSJames Wright case 2:
215b1559114SJames Wright if (wdetJ) StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ);
216b1559114SJames Wright if (normal) StoredValuesUnpack(Q, i, 1, 2, q_data, normal);
21793639ffbSJames Wright break;
21893639ffbSJames Wright case 3:
219b1559114SJames Wright if (wdetJ) StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ);
220b1559114SJames Wright if (normal) StoredValuesUnpack(Q, i, 1, 3, q_data, normal);
221b1559114SJames Wright if (dXdx) StoredValuesUnpack(Q, i, 4, 6, q_data, (CeedScalar *)dXdx);
22293639ffbSJames Wright break;
22393639ffbSJames Wright }
22493639ffbSJames Wright return CEED_ERROR_SUCCESS;
22593639ffbSJames Wright }
22693639ffbSJames Wright
StatePhysicalGradientFromReference_ND(CeedInt N,CeedInt Q,CeedInt i,NewtonianIdealGasContext gas,State s,StateVariable state_var,const CeedScalar * grad_q,const CeedScalar * dXdx,State * grad_s)22793639ffbSJames Wright CEED_QFUNCTION_HELPER void StatePhysicalGradientFromReference_ND(CeedInt N, CeedInt Q, CeedInt i, NewtonianIdealGasContext gas, State s,
22893639ffbSJames Wright StateVariable state_var, const CeedScalar *grad_q, const CeedScalar *dXdx,
22993639ffbSJames Wright State *grad_s) {
23093639ffbSJames Wright switch (N) {
23193639ffbSJames Wright case 2: {
23293639ffbSJames Wright for (CeedInt k = 0; k < 2; k++) {
23393639ffbSJames Wright CeedScalar dqi[5];
23493639ffbSJames Wright for (CeedInt j = 0; j < 5; j++) {
23593639ffbSJames 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];
23693639ffbSJames Wright }
23793639ffbSJames Wright grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var);
23893639ffbSJames Wright }
23993639ffbSJames Wright CeedScalar U[5] = {0.};
24093639ffbSJames Wright grad_s[2] = StateFromU(gas, U);
24193639ffbSJames Wright } break;
24293639ffbSJames Wright case 3:
243b1559114SJames Wright // Cannot directly use StatePhysicalGradientFromReference helper functions due to SYCL online compiler incompatabilities
244b1559114SJames Wright for (CeedInt k = 0; k < 3; k++) {
245b1559114SJames Wright CeedScalar dqi[5];
246b1559114SJames Wright for (CeedInt j = 0; j < 5; j++) {
247b1559114SJames 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] +
248b1559114SJames Wright grad_q[(Q * 5) * 2 + Q * j + i] * dXdx[2 * N + k];
249b1559114SJames Wright }
250b1559114SJames Wright grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var);
251b1559114SJames Wright }
25293639ffbSJames Wright break;
25393639ffbSJames Wright }
25493639ffbSJames Wright }
25593639ffbSJames Wright
256b4e9a8f8SJames Wright // @brief Calculate the stabilization constant \tau
Tau(AdvectionContext context,const State s,const CeedScalar * dXdx,CeedInt dim)257b4e9a8f8SJames Wright CEED_QFUNCTION_HELPER CeedScalar Tau(AdvectionContext context, const State s, const CeedScalar *dXdx, CeedInt dim) {
258b4e9a8f8SJames Wright switch (context->stabilization_tau) {
259b4e9a8f8SJames Wright case STAB_TAU_CTAU: {
260b4e9a8f8SJames Wright CeedScalar uX[3] = {0.};
261b4e9a8f8SJames Wright
262b4e9a8f8SJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
263b4e9a8f8SJames Wright return context->CtauS / sqrt(DotN(uX, uX, dim));
264b4e9a8f8SJames Wright } break;
265b4e9a8f8SJames Wright case STAB_TAU_ADVDIFF_SHAKIB: {
266b4e9a8f8SJames Wright CeedScalar gijd_mat[9] = {0.}, gij_uj[3] = {0.};
267b4e9a8f8SJames Wright
268b4e9a8f8SJames Wright MatMatN(dXdx, dXdx, dim, CEED_TRANSPOSE, CEED_NOTRANSPOSE, gijd_mat);
269b4e9a8f8SJames Wright MatVecNM(gijd_mat, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, gij_uj);
270b4e9a8f8SJames Wright return 1 / sqrt(Square(2 * context->Ctau_t / context->dt) + DotN(s.Y.velocity, gij_uj, dim) * context->Ctau_a);
271b4e9a8f8SJames Wright } break;
272b4e9a8f8SJames Wright default:
273b4e9a8f8SJames Wright return 0.;
274b4e9a8f8SJames Wright }
275b4e9a8f8SJames Wright }
276b4e9a8f8SJames Wright
27793639ffbSJames Wright // *****************************************************************************
27893639ffbSJames Wright // This QFunction implements Advection for implicit time stepping method
27993639ffbSJames Wright // *****************************************************************************
IFunction_AdvectionGeneric(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out,CeedInt dim)28093639ffbSJames Wright CEED_QFUNCTION_HELPER void IFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
28193639ffbSJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
28293639ffbSJames Wright const CeedScalar(*grad_q) = in[1];
28393639ffbSJames Wright const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
28493639ffbSJames Wright const CeedScalar(*q_data) = in[3];
28593639ffbSJames Wright
28693639ffbSJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
28793639ffbSJames Wright CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
28893639ffbSJames Wright CeedScalar *jac_data = out[2];
28993639ffbSJames Wright
29093639ffbSJames Wright AdvectionContext context = (AdvectionContext)ctx;
29193639ffbSJames Wright const CeedScalar zeros[14] = {0.};
29293639ffbSJames Wright NewtonianIdealGasContext gas;
29393639ffbSJames Wright struct NewtonianIdealGasContext_ gas_struct = {0};
29493639ffbSJames Wright gas = &gas_struct;
29593639ffbSJames Wright
29693639ffbSJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
29793639ffbSJames Wright const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
29893639ffbSJames Wright const State s = StateFromU(gas, qi);
29993639ffbSJames Wright
30093639ffbSJames Wright CeedScalar wdetJ, dXdx[9];
30193639ffbSJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
30293639ffbSJames Wright State grad_s[3];
30393639ffbSJames Wright StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
30493639ffbSJames Wright
30593639ffbSJames Wright const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
30693639ffbSJames Wright
30793639ffbSJames Wright for (CeedInt f = 0; f < 4; f++) {
30893639ffbSJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum
30993639ffbSJames Wright v[f][i] = wdetJ * q_dot[f][i]; // K Mass/transient term
31093639ffbSJames Wright }
31193639ffbSJames Wright
31293639ffbSJames Wright CeedScalar div_u = 0;
31393639ffbSJames Wright for (CeedInt j = 0; j < dim; j++) {
31493639ffbSJames Wright for (CeedInt k = 0; k < dim; k++) {
31593639ffbSJames Wright div_u += grad_s[k].Y.velocity[j];
31693639ffbSJames Wright }
31793639ffbSJames Wright }
31893639ffbSJames Wright CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
31993639ffbSJames Wright CeedScalar strong_res = q_dot[4][i] + strong_conv;
32093639ffbSJames Wright
32193639ffbSJames Wright v[4][i] = wdetJ * q_dot[4][i]; // transient part (ALWAYS)
32293639ffbSJames Wright
32393639ffbSJames Wright CeedScalar uX[3] = {0.};
32493639ffbSJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
32593639ffbSJames Wright
32693639ffbSJames Wright if (context->strong_form) { // Strong Galerkin convection term: v div(E u)
32793639ffbSJames Wright v[4][i] += wdetJ * strong_conv;
32893639ffbSJames Wright } else { // Weak Galerkin convection term: -dv \cdot (E u)
32993639ffbSJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = -wdetJ * s.U.E_total * uX[j];
33093639ffbSJames Wright }
33193639ffbSJames Wright
332d1d77723SJames Wright { // Diffusion
333d1d77723SJames Wright CeedScalar Fe[3], Fe_dXdx[3] = {0.};
334d1d77723SJames Wright
335d1d77723SJames Wright for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total;
336d1d77723SJames Wright MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx);
337d1d77723SJames Wright for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] -= wdetJ * Fe_dXdx[k];
338d1d77723SJames Wright }
339d1d77723SJames Wright
340b4e9a8f8SJames Wright const CeedScalar TauS = Tau(context, s, dXdx, dim);
34193639ffbSJames Wright for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) {
34293639ffbSJames Wright case STAB_NONE:
34393639ffbSJames Wright break;
34493639ffbSJames Wright case STAB_SU:
34593639ffbSJames Wright grad_v[j][4][i] += wdetJ * TauS * strong_conv * uX[j];
34693639ffbSJames Wright break;
34793639ffbSJames Wright case STAB_SUPG:
34893639ffbSJames Wright grad_v[j][4][i] += wdetJ * TauS * strong_res * uX[j];
34993639ffbSJames Wright break;
35093639ffbSJames Wright }
35193639ffbSJames Wright StoredValuesPack(Q, i, 0, 14, zeros, jac_data);
35293639ffbSJames Wright }
35393639ffbSJames Wright }
35493639ffbSJames Wright
IFunction_Advection(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)3552b730f8bSJeremy L Thompson CEED_QFUNCTION(IFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
356372d1924SJames Wright IFunction_AdvectionGeneric(ctx, Q, in, out, 3);
35777841947SLeila Ghaffari return 0;
35877841947SLeila Ghaffari }
35977841947SLeila Ghaffari
IFunction_Advection2d(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)36093639ffbSJames Wright CEED_QFUNCTION(IFunction_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
36193639ffbSJames Wright IFunction_AdvectionGeneric(ctx, Q, in, out, 2);
36293639ffbSJames Wright return 0;
36393639ffbSJames Wright }
36493639ffbSJames Wright
MassFunction_AdvectionGeneric(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out,CeedInt dim)365b4e9a8f8SJames Wright CEED_QFUNCTION_HELPER void MassFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
366b4e9a8f8SJames Wright const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
367b4e9a8f8SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1];
368b4e9a8f8SJames Wright const CeedScalar(*q_data) = in[2];
369b4e9a8f8SJames Wright
370b4e9a8f8SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
371b4e9a8f8SJames Wright CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
372b4e9a8f8SJames Wright
373b4e9a8f8SJames Wright AdvectionContext context = (AdvectionContext)ctx;
374b4e9a8f8SJames Wright struct NewtonianIdealGasContext_ gas_struct = {0};
375b4e9a8f8SJames Wright NewtonianIdealGasContext gas = &gas_struct;
376b4e9a8f8SJames Wright
377b4e9a8f8SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
378b4e9a8f8SJames Wright const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
379b4e9a8f8SJames Wright const State s = StateFromU(gas, qi);
380b4e9a8f8SJames Wright CeedScalar wdetJ, dXdx[9];
381b4e9a8f8SJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
382b4e9a8f8SJames Wright
383b4e9a8f8SJames Wright for (CeedInt f = 0; f < 4; f++) {
384b4e9a8f8SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum
385b4e9a8f8SJames Wright v[f][i] = wdetJ * q_dot[f][i]; // K Mass/transient term
386b4e9a8f8SJames Wright }
387b4e9a8f8SJames Wright
388b4e9a8f8SJames Wright // Unstabilized mass term
389b4e9a8f8SJames Wright v[4][i] = wdetJ * q_dot[4][i];
390b4e9a8f8SJames Wright
391b4e9a8f8SJames Wright // Stabilized mass term
392b4e9a8f8SJames Wright CeedScalar uX[3] = {0.};
393b4e9a8f8SJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
394b4e9a8f8SJames Wright const CeedScalar TauS = Tau(context, s, dXdx, dim);
395b4e9a8f8SJames Wright for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) {
396b4e9a8f8SJames Wright case STAB_NONE:
397b4e9a8f8SJames Wright case STAB_SU:
398b4e9a8f8SJames Wright grad_v[j][4][i] = 0;
399b4e9a8f8SJames Wright break; // These should be run with the unstabilized mass matrix anyways
400b4e9a8f8SJames Wright case STAB_SUPG:
401b4e9a8f8SJames Wright grad_v[j][4][i] = wdetJ * TauS * q_dot[4][i] * uX[j];
402b4e9a8f8SJames Wright break;
403b4e9a8f8SJames Wright }
404b4e9a8f8SJames Wright }
405b4e9a8f8SJames Wright }
406b4e9a8f8SJames Wright
MassFunction_Advection(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)407b4e9a8f8SJames Wright CEED_QFUNCTION(MassFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
408b4e9a8f8SJames Wright MassFunction_AdvectionGeneric(ctx, Q, in, out, 3);
409b4e9a8f8SJames Wright return 0;
410b4e9a8f8SJames Wright }
411b4e9a8f8SJames Wright
MassFunction_Advection2D(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)412b4e9a8f8SJames Wright CEED_QFUNCTION(MassFunction_Advection2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
413b4e9a8f8SJames Wright MassFunction_AdvectionGeneric(ctx, Q, in, out, 2);
414b4e9a8f8SJames Wright return 0;
415b4e9a8f8SJames Wright }
416b4e9a8f8SJames Wright
41793639ffbSJames Wright // *****************************************************************************
41893639ffbSJames Wright // This QFunction implements Advection for explicit time stepping method
41993639ffbSJames Wright // *****************************************************************************
RHSFunction_AdvectionGeneric(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out,CeedInt dim)42093639ffbSJames Wright CEED_QFUNCTION_HELPER void RHSFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
42193639ffbSJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
42293639ffbSJames Wright const CeedScalar(*grad_q) = in[1];
42393639ffbSJames Wright const CeedScalar(*q_data) = in[2];
42493639ffbSJames Wright
42593639ffbSJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
42693639ffbSJames Wright CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
42793639ffbSJames Wright
42893639ffbSJames Wright AdvectionContext context = (AdvectionContext)ctx;
42993639ffbSJames Wright struct NewtonianIdealGasContext_ gas_struct = {0};
430b4e9a8f8SJames Wright NewtonianIdealGasContext gas = &gas_struct;
43193639ffbSJames Wright
43293639ffbSJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
43393639ffbSJames Wright const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]};
43493639ffbSJames Wright const State s = StateFromU(gas, qi);
43593639ffbSJames Wright
43693639ffbSJames Wright CeedScalar wdetJ, dXdx[9];
43793639ffbSJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx);
43893639ffbSJames Wright State grad_s[3];
43993639ffbSJames Wright StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s);
44093639ffbSJames Wright
44193639ffbSJames Wright const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total};
44293639ffbSJames Wright
44393639ffbSJames Wright for (CeedInt f = 0; f < 4; f++) {
44493639ffbSJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum
44593639ffbSJames Wright v[f][i] = 0.;
44693639ffbSJames Wright }
44793639ffbSJames Wright
44893639ffbSJames Wright CeedScalar div_u = 0;
44993639ffbSJames Wright for (CeedInt j = 0; j < dim; j++) {
45093639ffbSJames Wright for (CeedInt k = 0; k < dim; k++) {
45193639ffbSJames Wright div_u += grad_s[k].Y.velocity[j];
45293639ffbSJames Wright }
45393639ffbSJames Wright }
45493639ffbSJames Wright CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim);
45593639ffbSJames Wright
45693639ffbSJames Wright CeedScalar uX[3] = {0.};
45793639ffbSJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX);
45893639ffbSJames Wright
45993639ffbSJames Wright if (context->strong_form) { // Strong Galerkin convection term: v div(E u)
46093639ffbSJames Wright v[4][i] = -wdetJ * strong_conv;
46193639ffbSJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = 0;
46293639ffbSJames Wright } else { // Weak Galerkin convection term: -dv \cdot (E u)
46393639ffbSJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = wdetJ * s.U.E_total * uX[j];
46493639ffbSJames Wright v[4][i] = 0.;
46593639ffbSJames Wright }
46693639ffbSJames Wright
467d1d77723SJames Wright { // Diffusion
468d1d77723SJames Wright CeedScalar Fe[3], Fe_dXdx[3] = {0.};
469d1d77723SJames Wright
470d1d77723SJames Wright for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total;
471d1d77723SJames Wright MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx);
472d1d77723SJames Wright for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] += wdetJ * Fe_dXdx[k];
473d1d77723SJames Wright }
474d1d77723SJames Wright
475b4e9a8f8SJames Wright const CeedScalar TauS = Tau(context, s, dXdx, dim);
47693639ffbSJames Wright for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) {
47793639ffbSJames Wright case STAB_NONE:
47893639ffbSJames Wright break;
47993639ffbSJames Wright case STAB_SU:
48093639ffbSJames Wright case STAB_SUPG:
481fa9a7b5dSJames Wright grad_v[j][4][i] -= wdetJ * TauS * strong_conv * uX[j];
48293639ffbSJames Wright break;
48393639ffbSJames Wright }
48493639ffbSJames Wright }
48593639ffbSJames Wright }
48693639ffbSJames Wright
RHS_Advection(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)48793639ffbSJames Wright CEED_QFUNCTION(RHS_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
48893639ffbSJames Wright RHSFunction_AdvectionGeneric(ctx, Q, in, out, 3);
48993639ffbSJames Wright return 0;
49093639ffbSJames Wright }
49193639ffbSJames Wright
RHS_Advection2d(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)49293639ffbSJames Wright CEED_QFUNCTION(RHS_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
49393639ffbSJames Wright RHSFunction_AdvectionGeneric(ctx, Q, in, out, 2);
49493639ffbSJames Wright return 0;
49593639ffbSJames Wright }
49693639ffbSJames Wright
49793639ffbSJames Wright // *****************************************************************************
49893639ffbSJames Wright // This QFunction implements consistent outflow and inflow BCs
49993639ffbSJames Wright // for advection
50093639ffbSJames Wright //
50193639ffbSJames Wright // Inflow and outflow faces are determined based on sign(dot(wind, normal)):
50293639ffbSJames Wright // sign(dot(wind, normal)) > 0 : outflow BCs
50393639ffbSJames Wright // sign(dot(wind, normal)) < 0 : inflow BCs
50493639ffbSJames Wright //
50593639ffbSJames Wright // Outflow BCs:
50693639ffbSJames Wright // The validity of the weak form of the governing equations is extended to the outflow and the current values of E are applied.
50793639ffbSJames Wright //
50893639ffbSJames Wright // Inflow BCs:
50993639ffbSJames Wright // A prescribed Total Energy (E_wind) is applied weakly.
51093639ffbSJames Wright // *****************************************************************************
Advection_InOutFlowGeneric(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out,CeedInt dim)51193639ffbSJames Wright CEED_QFUNCTION(Advection_InOutFlowGeneric)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) {
51293639ffbSJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
51393639ffbSJames Wright const CeedScalar(*q_data_sur) = in[2];
51493639ffbSJames Wright
51593639ffbSJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
51693639ffbSJames Wright AdvectionContext context = (AdvectionContext)ctx;
51793639ffbSJames Wright const CeedScalar E_wind = context->E_wind;
51893639ffbSJames Wright const CeedScalar strong_form = context->strong_form;
51993639ffbSJames Wright const bool is_implicit = context->implicit;
52093639ffbSJames Wright
52193639ffbSJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
52293639ffbSJames Wright const CeedScalar rho = q[0][i];
52393639ffbSJames Wright const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
52493639ffbSJames Wright const CeedScalar E = q[4][i];
52593639ffbSJames Wright
52693639ffbSJames Wright CeedScalar wdetJb, norm[3];
52793639ffbSJames Wright QdataBoundaryUnpack_ND(dim, Q, i, q_data_sur, &wdetJb, NULL, norm);
52893639ffbSJames Wright wdetJb *= is_implicit ? -1. : 1.;
52993639ffbSJames Wright
53093639ffbSJames Wright const CeedScalar u_normal = DotN(norm, u, dim);
53193639ffbSJames Wright
53293639ffbSJames Wright // No Change in density or momentum
53393639ffbSJames Wright for (CeedInt j = 0; j < 4; j++) {
53493639ffbSJames Wright v[j][i] = 0;
53593639ffbSJames Wright }
53693639ffbSJames Wright // Implementing in/outflow BCs
53793639ffbSJames Wright if (u_normal > 0) { // outflow
53893639ffbSJames Wright v[4][i] = -(1 - strong_form) * wdetJb * E * u_normal;
53993639ffbSJames Wright } else { // inflow
54093639ffbSJames Wright v[4][i] = -(1 - strong_form) * wdetJb * E_wind * u_normal;
54193639ffbSJames Wright }
54293639ffbSJames Wright }
54393639ffbSJames Wright return 0;
54493639ffbSJames Wright }
54593639ffbSJames Wright
Advection_InOutFlow(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)5462b730f8bSJeremy L Thompson CEED_QFUNCTION(Advection_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
5474bdcf5bfSJames Wright Advection_InOutFlowGeneric(ctx, Q, in, out, 3);
54877841947SLeila Ghaffari return 0;
54977841947SLeila Ghaffari }
55077841947SLeila Ghaffari
Advection2d_InOutFlow(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)55193639ffbSJames Wright CEED_QFUNCTION(Advection2d_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
55293639ffbSJames Wright Advection_InOutFlowGeneric(ctx, Q, in, out, 2);
55393639ffbSJames Wright return 0;
55493639ffbSJames Wright }
555