1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3a515125bSLeila Ghaffari 4a515125bSLeila Ghaffari /// @file 5a515125bSLeila Ghaffari /// Advection initial condition and operator for Navier-Stokes example using PETSc 6493642f1SJames Wright #include <ceed.h> 7d0cce58aSJeremy L Thompson #include <math.h> 8a515125bSLeila Ghaffari 9e88b842aSJames Wright #include "advection_types.h" 10ce192147SJames Wright #include "newtonian_state.h" 11ce192147SJames Wright #include "newtonian_types.h" 12e88b842aSJames Wright #include "stabilization_types.h" 131a74fa30SJames Wright #include "utils.h" 141a74fa30SJames Wright 15a515125bSLeila Ghaffari // ***************************************************************************** 169529d636SJames Wright // This QFunction sets the initial conditions and the boundary conditions 179529d636SJames Wright // for two test cases: ROTATION and TRANSLATION 189529d636SJames Wright // 199529d636SJames Wright // -- ROTATION (default) 209529d636SJames Wright // Initial Conditions: 219529d636SJames Wright // Mass Density: 229529d636SJames Wright // Constant mass density of 1.0 239529d636SJames Wright // Momentum Density: 249529d636SJames Wright // Rotational field in x,y 259529d636SJames Wright // Energy Density: 269529d636SJames Wright // Maximum of 1. x0 decreasing linearly to 0. as radial distance 279529d636SJames Wright // increases to (1.-r/rc), then 0. everywhere else 289529d636SJames Wright // 299529d636SJames Wright // Boundary Conditions: 309529d636SJames Wright // Mass Density: 319529d636SJames Wright // 0.0 flux 329529d636SJames Wright // Momentum Density: 339529d636SJames Wright // 0.0 349529d636SJames Wright // Energy Density: 359529d636SJames Wright // 0.0 flux 369529d636SJames Wright // 379529d636SJames Wright // -- TRANSLATION 389529d636SJames Wright // Initial Conditions: 399529d636SJames Wright // Mass Density: 409529d636SJames Wright // Constant mass density of 1.0 419529d636SJames Wright // Momentum Density: 429529d636SJames Wright // Constant rectilinear field in x,y 439529d636SJames Wright // Energy Density: 449529d636SJames Wright // Maximum of 1. x0 decreasing linearly to 0. as radial distance 459529d636SJames Wright // increases to (1.-r/rc), then 0. everywhere else 469529d636SJames Wright // 479529d636SJames Wright // Boundary Conditions: 489529d636SJames Wright // Mass Density: 499529d636SJames Wright // 0.0 flux 509529d636SJames Wright // Momentum Density: 519529d636SJames Wright // 0.0 529529d636SJames Wright // Energy Density: 539529d636SJames Wright // Inflow BCs: 549529d636SJames Wright // E = E_wind 559529d636SJames Wright // Outflow BCs: 569529d636SJames Wright // E = E(boundary) 579529d636SJames Wright // Both In/Outflow BCs for E are applied weakly in the 589529d636SJames Wright // QFunction "Advection2d_Sur" 599529d636SJames Wright // 609529d636SJames Wright // ***************************************************************************** 619529d636SJames Wright 629529d636SJames Wright // ***************************************************************************** 639529d636SJames Wright // This helper function provides the exact, time-dependent solution and IC formulation for 2D advection 649529d636SJames Wright // ***************************************************************************** 65*97cfd714SJames Wright CEED_QFUNCTION_HELPER int Exact_AdvectionGeneric(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) { 669529d636SJames Wright const SetupContextAdv context = (SetupContextAdv)ctx; 679529d636SJames Wright const CeedScalar rc = context->rc; 689529d636SJames Wright const CeedScalar lx = context->lx; 699529d636SJames Wright const CeedScalar ly = context->ly; 709529d636SJames Wright const CeedScalar lz = dim == 2 ? 0. : context->lz; 719529d636SJames Wright const CeedScalar *wind = context->wind; 729529d636SJames Wright 739529d636SJames Wright const CeedScalar center[3] = {0.5 * lx, 0.5 * ly, 0.5 * lz}; 749529d636SJames Wright const CeedScalar theta = dim == 2 ? M_PI / 3 : M_PI; 759529d636SJames Wright const CeedScalar x0[3] = {center[0] + .25 * lx * cos(theta + time), center[1] + .25 * ly * sin(theta + time), 0.5 * lz}; 769529d636SJames Wright 779529d636SJames Wright const CeedScalar x = X[0], y = X[1], z = dim == 2 ? 0. : X[2]; 789529d636SJames Wright 799529d636SJames Wright switch (context->wind_type) { 805f636aeaSJames Wright case ADVDIF_WIND_ROTATION: 819529d636SJames Wright q[0] = 1.; 829529d636SJames Wright q[1] = -(y - center[1]); 839529d636SJames Wright q[2] = (x - center[0]); 849529d636SJames Wright q[3] = 0; 859529d636SJames Wright break; 865f636aeaSJames Wright case ADVDIF_WIND_TRANSLATION: 879529d636SJames Wright q[0] = 1.; 889529d636SJames Wright q[1] = wind[0]; 899529d636SJames Wright q[2] = wind[1]; 909529d636SJames Wright q[3] = dim == 2 ? 0. : wind[2]; 919529d636SJames Wright break; 929529d636SJames Wright default: 939529d636SJames Wright return 1; 949529d636SJames Wright } 959529d636SJames Wright 969529d636SJames Wright switch (context->initial_condition_type) { 975f636aeaSJames Wright case ADVDIF_IC_BUBBLE_SPHERE: 985f636aeaSJames Wright case ADVDIF_IC_BUBBLE_CYLINDER: { 99a62be6baSJames Wright CeedScalar r = sqrt(Square(x - x0[0]) + Square(y - x0[1]) + Square(z - x0[2])); 100a62be6baSJames Wright 1019529d636SJames Wright switch (context->bubble_continuity_type) { 1029529d636SJames Wright // original continuous, smooth shape 1035f636aeaSJames Wright case ADVDIF_BUBBLE_CONTINUITY_SMOOTH: 1049529d636SJames Wright q[4] = r <= rc ? (1. - r / rc) : 0.; 1059529d636SJames Wright break; 1069529d636SJames Wright // discontinuous, sharp back half shape 1075f636aeaSJames Wright case ADVDIF_BUBBLE_CONTINUITY_BACK_SHARP: 1089529d636SJames Wright q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) : 0.; 1099529d636SJames Wright break; 1109529d636SJames Wright // attempt to define a finite thickness that will get resolved under grid refinement 1115f636aeaSJames Wright case ADVDIF_BUBBLE_CONTINUITY_THICK: 1129529d636SJames Wright q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) * fmin(1.0, (center[1] - y) / 1.25) : 0.; 1139529d636SJames Wright break; 1145f636aeaSJames Wright case ADVDIF_BUBBLE_CONTINUITY_COSINE: 1159529d636SJames Wright q[4] = r <= rc ? .5 + .5 * cos(r * M_PI / rc) : 0; 1169529d636SJames Wright break; 1179529d636SJames Wright } 1189529d636SJames Wright break; 119a62be6baSJames Wright } 120a62be6baSJames Wright 1215f636aeaSJames Wright case ADVDIF_IC_COSINE_HILL: { 122a62be6baSJames Wright CeedScalar r = sqrt(Square(x - center[0]) + Square(y - center[1])); 1239529d636SJames Wright CeedScalar half_width = context->lx / 2; 1249529d636SJames Wright q[4] = r > half_width ? 0. : cos(2 * M_PI * r / half_width + M_PI) + 1.; 1259529d636SJames Wright } break; 126a62be6baSJames Wright 1275f636aeaSJames Wright case ADVDIF_IC_SKEW: { 1289529d636SJames Wright CeedScalar skewed_barrier[3] = {wind[0], wind[1], 0}; 1299529d636SJames Wright CeedScalar inflow_to_point[3] = {x - context->lx / 2, y, 0}; 1309529d636SJames Wright CeedScalar cross_product[3] = {0}; 1319529d636SJames Wright const CeedScalar boundary_threshold = 20 * CEED_EPSILON; 1329529d636SJames Wright Cross3(skewed_barrier, inflow_to_point, cross_product); 1339529d636SJames Wright 1349529d636SJames Wright q[4] = cross_product[2] > boundary_threshold ? 0 : 1; 1359529d636SJames Wright if ((x < boundary_threshold && wind[0] < boundary_threshold) || // outflow at -x boundary 1369529d636SJames Wright (y < boundary_threshold && wind[1] < boundary_threshold) || // outflow at -y boundary 1379529d636SJames Wright (x > context->lx - boundary_threshold && wind[0] > boundary_threshold) || // outflow at +x boundary 1389529d636SJames Wright (y > context->ly - boundary_threshold && wind[1] > boundary_threshold) // outflow at +y boundary 1399529d636SJames Wright ) { 1409529d636SJames Wright q[4] = 0; 1419529d636SJames Wright } 1429529d636SJames Wright } break; 143a62be6baSJames Wright 1445f636aeaSJames Wright case ADVDIF_IC_WAVE: { 145a62be6baSJames Wright CeedScalar theta = context->wave_frequency * DotN(X, wind, dim) + context->wave_phase; 146a62be6baSJames Wright switch (context->wave_type) { 147a62be6baSJames Wright case ADVDIF_WAVE_SINE: 148a62be6baSJames Wright q[4] = sin(theta); 149a62be6baSJames Wright break; 150a62be6baSJames Wright case ADVDIF_WAVE_SQUARE: 151a62be6baSJames Wright q[4] = sin(theta) > 100 * CEED_EPSILON ? 1 : -1; 152a62be6baSJames Wright break; 153a62be6baSJames Wright } 154a62be6baSJames Wright } 1559529d636SJames Wright } 1569529d636SJames Wright return 0; 1579529d636SJames Wright } 1589529d636SJames Wright 1599529d636SJames Wright // ***************************************************************************** 160a515125bSLeila Ghaffari // This QFunction sets the initial conditions for 3D advection 161a515125bSLeila Ghaffari // ***************************************************************************** 1622b916ea7SJeremy L Thompson CEED_QFUNCTION(ICsAdvection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 163a515125bSLeila Ghaffari const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 164a515125bSLeila Ghaffari CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 165a515125bSLeila Ghaffari 1663d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 167a515125bSLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 168139613f2SLeila Ghaffari CeedScalar q[5] = {0.}; 169a515125bSLeila Ghaffari 1700b3a1fabSJames Wright Exact_AdvectionGeneric(3, 0., x, 5, q, ctx); 171a515125bSLeila Ghaffari for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 1720b3a1fabSJames Wright } 173a515125bSLeila Ghaffari return 0; 174a515125bSLeila Ghaffari } 175a515125bSLeila Ghaffari 176a515125bSLeila Ghaffari // ***************************************************************************** 1779529d636SJames Wright // This QFunction sets the initial conditions for 2D advection 178a515125bSLeila Ghaffari // ***************************************************************************** 1799529d636SJames Wright CEED_QFUNCTION(ICsAdvection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 1809529d636SJames Wright const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 1819529d636SJames Wright CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 1829529d636SJames Wright const SetupContextAdv context = (SetupContextAdv)ctx; 1839529d636SJames Wright 1849529d636SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 1859529d636SJames Wright const CeedScalar x[] = {X[0][i], X[1][i]}; 1869529d636SJames Wright CeedScalar q[5] = {0.}; 1879529d636SJames Wright 1889529d636SJames Wright Exact_AdvectionGeneric(2, context->time, x, 5, q, ctx); 1899529d636SJames Wright for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 1909529d636SJames Wright } 191a515125bSLeila Ghaffari return 0; 192a515125bSLeila Ghaffari } 193a515125bSLeila Ghaffari 1949529d636SJames Wright CEED_QFUNCTION_HELPER void StatePhysicalGradientFromReference_ND(CeedInt N, CeedInt Q, CeedInt i, NewtonianIdealGasContext gas, State s, 1959529d636SJames Wright StateVariable state_var, const CeedScalar *grad_q, const CeedScalar *dXdx, 1969529d636SJames Wright State *grad_s) { 1979529d636SJames Wright switch (N) { 1989529d636SJames Wright case 2: { 1999529d636SJames Wright for (CeedInt k = 0; k < 2; k++) { 2009529d636SJames Wright CeedScalar dqi[5]; 2019529d636SJames Wright for (CeedInt j = 0; j < 5; j++) { 2029529d636SJames 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]; 2039529d636SJames Wright } 2049529d636SJames Wright grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var); 2059529d636SJames Wright } 2069529d636SJames Wright CeedScalar U[5] = {0.}; 2079529d636SJames Wright grad_s[2] = StateFromU(gas, U); 2089529d636SJames Wright } break; 2099529d636SJames Wright case 3: 21085efd435SJames Wright // Cannot directly use StatePhysicalGradientFromReference helper functions due to SYCL online compiler incompatabilities 21185efd435SJames Wright for (CeedInt k = 0; k < 3; k++) { 21285efd435SJames Wright CeedScalar dqi[5]; 21385efd435SJames Wright for (CeedInt j = 0; j < 5; j++) { 21485efd435SJames 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] + 21585efd435SJames Wright grad_q[(Q * 5) * 2 + Q * j + i] * dXdx[2 * N + k]; 21685efd435SJames Wright } 21785efd435SJames Wright grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var); 21885efd435SJames Wright } 2199529d636SJames Wright break; 2209529d636SJames Wright } 2219529d636SJames Wright } 2229529d636SJames Wright 223a78efa86SJames Wright // @brief Calculate the stabilization constant \tau 224a78efa86SJames Wright CEED_QFUNCTION_HELPER CeedScalar Tau(AdvectionContext context, const State s, const CeedScalar *dXdx, CeedInt dim) { 225a78efa86SJames Wright switch (context->stabilization_tau) { 226a78efa86SJames Wright case STAB_TAU_CTAU: { 227a78efa86SJames Wright CeedScalar uX[3] = {0.}; 228a78efa86SJames Wright 229a78efa86SJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 230a78efa86SJames Wright return context->CtauS / sqrt(DotN(uX, uX, dim)); 231a78efa86SJames Wright } break; 232a78efa86SJames Wright case STAB_TAU_ADVDIFF_SHAKIB: { 233a78efa86SJames Wright CeedScalar gijd_mat[9] = {0.}, gij_uj[3] = {0.}; 234a78efa86SJames Wright 235a78efa86SJames Wright MatMatN(dXdx, dXdx, dim, CEED_TRANSPOSE, CEED_NOTRANSPOSE, gijd_mat); 236a78efa86SJames Wright MatVecNM(gijd_mat, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, gij_uj); 237fbabb365SJames Wright return 1 / sqrt(Square(2 * context->Ctau_t / context->dt) + DotN(s.Y.velocity, gij_uj, dim) * context->Ctau_a + 238fbabb365SJames Wright Square(context->diffusion_coeff) * DotN(gijd_mat, gijd_mat, dim * dim) * context->Ctau_d); 239a78efa86SJames Wright } break; 240a78efa86SJames Wright default: 241a78efa86SJames Wright return 0.; 242a78efa86SJames Wright } 243a78efa86SJames Wright } 244a78efa86SJames Wright 2459529d636SJames Wright // ***************************************************************************** 2469529d636SJames Wright // This QFunction implements Advection for implicit time stepping method 2479529d636SJames Wright // ***************************************************************************** 248*97cfd714SJames Wright CEED_QFUNCTION_HELPER int IFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 2494c5ab12fSJames Wright AdvectionContext context = (AdvectionContext)ctx; 2504c5ab12fSJames Wright 2519529d636SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 2529529d636SJames Wright const CeedScalar(*grad_q) = in[1]; 2539529d636SJames Wright const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 2549529d636SJames Wright const CeedScalar(*q_data) = in[3]; 2554c5ab12fSJames Wright const CeedScalar(*divFdiff) = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? in[5] : NULL; 2569529d636SJames Wright 2579529d636SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 2589529d636SJames Wright CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 2599529d636SJames Wright 2609529d636SJames Wright NewtonianIdealGasContext gas; 2619529d636SJames Wright struct NewtonianIdealGasContext_ gas_struct = {0}; 2629529d636SJames Wright gas = &gas_struct; 2639529d636SJames Wright 2649529d636SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 2659529d636SJames Wright const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]}; 2669529d636SJames Wright const State s = StateFromU(gas, qi); 2679529d636SJames Wright 2689529d636SJames Wright CeedScalar wdetJ, dXdx[9]; 2699529d636SJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 2709529d636SJames Wright State grad_s[3]; 2719529d636SJames Wright StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s); 2729529d636SJames Wright 2739529d636SJames Wright const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total}; 2749529d636SJames Wright 2759529d636SJames Wright for (CeedInt f = 0; f < 4; f++) { 2769529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum 2779529d636SJames Wright v[f][i] = wdetJ * q_dot[f][i]; // K Mass/transient term 2789529d636SJames Wright } 2799529d636SJames Wright 2809529d636SJames Wright CeedScalar div_u = 0; 2819529d636SJames Wright for (CeedInt j = 0; j < dim; j++) { 2829529d636SJames Wright for (CeedInt k = 0; k < dim; k++) { 2839529d636SJames Wright div_u += grad_s[k].Y.velocity[j]; 2849529d636SJames Wright } 2859529d636SJames Wright } 2869529d636SJames Wright CeedScalar uX[3] = {0.}; 2879529d636SJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 2884c5ab12fSJames Wright CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim); 2899529d636SJames Wright 2904c5ab12fSJames Wright v[4][i] = wdetJ * q_dot[4][i]; // transient part (ALWAYS) 2919529d636SJames Wright if (context->strong_form) { // Strong Galerkin convection term: v div(E u) 2929529d636SJames Wright v[4][i] += wdetJ * strong_conv; 2939529d636SJames Wright } else { // Weak Galerkin convection term: -dv \cdot (E u) 2949529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = -wdetJ * s.U.E_total * uX[j]; 2959529d636SJames Wright } 2969529d636SJames Wright 297c8d249deSJames Wright { // Diffusion 298c8d249deSJames Wright CeedScalar Fe[3], Fe_dXdx[3] = {0.}; 299c8d249deSJames Wright 300c8d249deSJames Wright for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total; 301c8d249deSJames Wright MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx); 302c8d249deSJames Wright for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] -= wdetJ * Fe_dXdx[k]; 303c8d249deSJames Wright } 304c8d249deSJames Wright 305a78efa86SJames Wright const CeedScalar TauS = Tau(context, s, dXdx, dim); 3064c5ab12fSJames Wright for (CeedInt j = 0; j < dim; j++) { 3074c5ab12fSJames Wright switch (context->stabilization) { 3089529d636SJames Wright case STAB_NONE: 3099529d636SJames Wright break; 3109529d636SJames Wright case STAB_SU: 3114c5ab12fSJames Wright grad_v[j][4][i] += wdetJ * TauS * uX[j] * strong_conv; 3129529d636SJames Wright break; 3134c5ab12fSJames Wright case STAB_SUPG: { 3144c5ab12fSJames Wright CeedScalar divFdiff_i = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? divFdiff[i] : 0.; 3154c5ab12fSJames Wright grad_v[j][4][i] += wdetJ * TauS * uX[j] * (q_dot[4][i] + strong_conv + divFdiff_i); 3164c5ab12fSJames Wright } break; 3174c5ab12fSJames Wright } 3189529d636SJames Wright } 3199529d636SJames Wright } 320*97cfd714SJames Wright return 0; 3219529d636SJames Wright } 3229529d636SJames Wright 3232b916ea7SJeremy L Thompson CEED_QFUNCTION(IFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 324*97cfd714SJames Wright return IFunction_AdvectionGeneric(ctx, Q, in, out, 3); 325a515125bSLeila Ghaffari } 326a515125bSLeila Ghaffari 3279529d636SJames Wright CEED_QFUNCTION(IFunction_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 328*97cfd714SJames Wright return IFunction_AdvectionGeneric(ctx, Q, in, out, 2); 3299529d636SJames Wright } 3309529d636SJames Wright 331*97cfd714SJames Wright CEED_QFUNCTION_HELPER int MassFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 332a78efa86SJames Wright const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 333a78efa86SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 334a78efa86SJames Wright const CeedScalar(*q_data) = in[2]; 335a78efa86SJames Wright 336a78efa86SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 337a78efa86SJames Wright CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 338a78efa86SJames Wright 339a78efa86SJames Wright AdvectionContext context = (AdvectionContext)ctx; 340a78efa86SJames Wright struct NewtonianIdealGasContext_ gas_struct = {0}; 341a78efa86SJames Wright NewtonianIdealGasContext gas = &gas_struct; 342a78efa86SJames Wright 343a78efa86SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 344a78efa86SJames Wright const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]}; 345a78efa86SJames Wright const State s = StateFromU(gas, qi); 346a78efa86SJames Wright CeedScalar wdetJ, dXdx[9]; 347a78efa86SJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 348a78efa86SJames Wright 349a78efa86SJames Wright for (CeedInt f = 0; f < 4; f++) { 350a78efa86SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum 351a78efa86SJames Wright v[f][i] = wdetJ * q_dot[f][i]; // K Mass/transient term 352a78efa86SJames Wright } 353a78efa86SJames Wright 354a78efa86SJames Wright // Unstabilized mass term 355a78efa86SJames Wright v[4][i] = wdetJ * q_dot[4][i]; 356a78efa86SJames Wright 357a78efa86SJames Wright // Stabilized mass term 358a78efa86SJames Wright CeedScalar uX[3] = {0.}; 359a78efa86SJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 360a78efa86SJames Wright const CeedScalar TauS = Tau(context, s, dXdx, dim); 361a78efa86SJames Wright for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) { 362a78efa86SJames Wright case STAB_NONE: 363a78efa86SJames Wright case STAB_SU: 364a78efa86SJames Wright grad_v[j][4][i] = 0; 365a78efa86SJames Wright break; // These should be run with the unstabilized mass matrix anyways 366a78efa86SJames Wright case STAB_SUPG: 367a78efa86SJames Wright grad_v[j][4][i] = wdetJ * TauS * q_dot[4][i] * uX[j]; 368a78efa86SJames Wright break; 369a78efa86SJames Wright } 370a78efa86SJames Wright } 371*97cfd714SJames Wright return 0; 372a78efa86SJames Wright } 373a78efa86SJames Wright 374a78efa86SJames Wright CEED_QFUNCTION(MassFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 375*97cfd714SJames Wright return MassFunction_AdvectionGeneric(ctx, Q, in, out, 3); 376a78efa86SJames Wright } 377a78efa86SJames Wright 378a78efa86SJames Wright CEED_QFUNCTION(MassFunction_Advection2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 379*97cfd714SJames Wright return MassFunction_AdvectionGeneric(ctx, Q, in, out, 2); 380a78efa86SJames Wright } 381a78efa86SJames Wright 3829529d636SJames Wright // ***************************************************************************** 3839529d636SJames Wright // This QFunction implements Advection for explicit time stepping method 3849529d636SJames Wright // ***************************************************************************** 385*97cfd714SJames Wright CEED_QFUNCTION_HELPER int RHSFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 3865f952e8dSJames Wright AdvectionContext context = (AdvectionContext)ctx; 3875f952e8dSJames Wright 3889529d636SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 3899529d636SJames Wright const CeedScalar(*grad_q) = in[1]; 3909529d636SJames Wright const CeedScalar(*q_data) = in[2]; 3915f952e8dSJames Wright const CeedScalar(*divFdiff) = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? in[4] : NULL; 3929529d636SJames Wright 3939529d636SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 3949529d636SJames Wright CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 3959529d636SJames Wright 3969529d636SJames Wright struct NewtonianIdealGasContext_ gas_struct = {0}; 397a78efa86SJames Wright NewtonianIdealGasContext gas = &gas_struct; 3989529d636SJames Wright 3999529d636SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 4009529d636SJames Wright const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]}; 4019529d636SJames Wright const State s = StateFromU(gas, qi); 4029529d636SJames Wright 4039529d636SJames Wright CeedScalar wdetJ, dXdx[9]; 4049529d636SJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 4059529d636SJames Wright State grad_s[3]; 4069529d636SJames Wright StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s); 4079529d636SJames Wright 4089529d636SJames Wright const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total}; 4099529d636SJames Wright 4109529d636SJames Wright for (CeedInt f = 0; f < 4; f++) { 4119529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum 4129529d636SJames Wright v[f][i] = 0.; 4139529d636SJames Wright } 4149529d636SJames Wright 4159529d636SJames Wright CeedScalar div_u = 0; 4169529d636SJames Wright for (CeedInt j = 0; j < dim; j++) { 4179529d636SJames Wright for (CeedInt k = 0; k < dim; k++) { 4189529d636SJames Wright div_u += grad_s[k].Y.velocity[j]; 4199529d636SJames Wright } 4209529d636SJames Wright } 4219529d636SJames Wright CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim); 4229529d636SJames Wright 4239529d636SJames Wright CeedScalar uX[3] = {0.}; 4249529d636SJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 4259529d636SJames Wright 4269529d636SJames Wright if (context->strong_form) { // Strong Galerkin convection term: v div(E u) 4279529d636SJames Wright v[4][i] = -wdetJ * strong_conv; 4289529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = 0; 4299529d636SJames Wright } else { // Weak Galerkin convection term: -dv \cdot (E u) 4309529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = wdetJ * s.U.E_total * uX[j]; 4319529d636SJames Wright v[4][i] = 0.; 4329529d636SJames Wright } 4339529d636SJames Wright 434c8d249deSJames Wright { // Diffusion 435c8d249deSJames Wright CeedScalar Fe[3], Fe_dXdx[3] = {0.}; 436c8d249deSJames Wright 437c8d249deSJames Wright for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total; 438c8d249deSJames Wright MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx); 439c8d249deSJames Wright for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] += wdetJ * Fe_dXdx[k]; 440c8d249deSJames Wright } 441c8d249deSJames Wright 442a78efa86SJames Wright const CeedScalar TauS = Tau(context, s, dXdx, dim); 4435f952e8dSJames Wright for (CeedInt j = 0; j < dim; j++) { 4445f952e8dSJames Wright switch (context->stabilization) { 4459529d636SJames Wright case STAB_NONE: 4469529d636SJames Wright break; 4479529d636SJames Wright case STAB_SU: 4485f952e8dSJames Wright case STAB_SUPG: { 4495f952e8dSJames Wright CeedScalar divFdiff_i = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? divFdiff[i] : 0.; 4505f952e8dSJames Wright grad_v[j][4][i] -= wdetJ * TauS * (strong_conv + divFdiff_i) * uX[j]; 4515f952e8dSJames Wright } break; 4525f952e8dSJames Wright } 4539529d636SJames Wright } 4549529d636SJames Wright } 455*97cfd714SJames Wright return 0; 4569529d636SJames Wright } 4579529d636SJames Wright 4589529d636SJames Wright CEED_QFUNCTION(RHS_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 459*97cfd714SJames Wright return RHSFunction_AdvectionGeneric(ctx, Q, in, out, 3); 4609529d636SJames Wright } 4619529d636SJames Wright 4629529d636SJames Wright CEED_QFUNCTION(RHS_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 463*97cfd714SJames Wright return RHSFunction_AdvectionGeneric(ctx, Q, in, out, 2); 4649529d636SJames Wright } 4659529d636SJames Wright 4669529d636SJames Wright // ***************************************************************************** 4679529d636SJames Wright // This QFunction implements consistent outflow and inflow BCs 4689529d636SJames Wright // for advection 4699529d636SJames Wright // 4709529d636SJames Wright // Inflow and outflow faces are determined based on sign(dot(wind, normal)): 4719529d636SJames Wright // sign(dot(wind, normal)) > 0 : outflow BCs 4729529d636SJames Wright // sign(dot(wind, normal)) < 0 : inflow BCs 4739529d636SJames Wright // 4749529d636SJames Wright // Outflow BCs: 4759529d636SJames Wright // The validity of the weak form of the governing equations is extended to the outflow and the current values of E are applied. 4769529d636SJames Wright // 4779529d636SJames Wright // Inflow BCs: 4789529d636SJames Wright // A prescribed Total Energy (E_wind) is applied weakly. 4799529d636SJames Wright // ***************************************************************************** 480*97cfd714SJames Wright CEED_QFUNCTION_HELPER int Advection_InOutFlowGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 4819529d636SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 4829529d636SJames Wright const CeedScalar(*q_data_sur) = in[2]; 4839529d636SJames Wright 4849529d636SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 4859529d636SJames Wright AdvectionContext context = (AdvectionContext)ctx; 4869529d636SJames Wright const CeedScalar E_wind = context->E_wind; 4879529d636SJames Wright const CeedScalar strong_form = context->strong_form; 4889529d636SJames Wright const bool is_implicit = context->implicit; 4899529d636SJames Wright 4909529d636SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 4919529d636SJames Wright const CeedScalar rho = q[0][i]; 4929529d636SJames Wright const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 4939529d636SJames Wright const CeedScalar E = q[4][i]; 4949529d636SJames Wright 49578e8b7daSJames Wright CeedScalar wdetJb, normal[3]; 49678e8b7daSJames Wright QdataBoundaryUnpack_ND(dim, Q, i, q_data_sur, &wdetJb, NULL, normal); 4979529d636SJames Wright wdetJb *= is_implicit ? -1. : 1.; 4989529d636SJames Wright 49978e8b7daSJames Wright const CeedScalar u_normal = DotN(normal, u, dim); 5009529d636SJames Wright 5019529d636SJames Wright // No Change in density or momentum 5029529d636SJames Wright for (CeedInt j = 0; j < 4; j++) { 5039529d636SJames Wright v[j][i] = 0; 5049529d636SJames Wright } 5059529d636SJames Wright // Implementing in/outflow BCs 5069529d636SJames Wright if (u_normal > 0) { // outflow 5079529d636SJames Wright v[4][i] = -(1 - strong_form) * wdetJb * E * u_normal; 5089529d636SJames Wright } else { // inflow 5099529d636SJames Wright v[4][i] = -(1 - strong_form) * wdetJb * E_wind * u_normal; 5109529d636SJames Wright } 5119529d636SJames Wright } 5129529d636SJames Wright return 0; 5139529d636SJames Wright } 5149529d636SJames Wright 5152b916ea7SJeremy L Thompson CEED_QFUNCTION(Advection_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 516*97cfd714SJames Wright return Advection_InOutFlowGeneric(ctx, Q, in, out, 3); 517a515125bSLeila Ghaffari } 518a515125bSLeila Ghaffari 5199529d636SJames Wright CEED_QFUNCTION(Advection2d_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 520*97cfd714SJames Wright return Advection_InOutFlowGeneric(ctx, Q, in, out, 2); 5219529d636SJames Wright } 522c2d90829SJames Wright 523c2d90829SJames Wright // @brief Volume integral for RHS of divergence of diffusive flux direct projection 524c2d90829SJames Wright CEED_QFUNCTION_HELPER int DivDiffusiveFluxVolumeRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, 525c2d90829SJames Wright const CeedInt dim) { 526c2d90829SJames Wright const CeedScalar(*Grad_q) = in[0]; 527c2d90829SJames Wright const CeedScalar(*q_data) = in[1]; 528c2d90829SJames Wright CeedScalar(*Grad_v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 529c2d90829SJames Wright 530c2d90829SJames Wright AdvectionContext context = (AdvectionContext)ctx; 531c2d90829SJames Wright 532c2d90829SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 533c2d90829SJames Wright CeedScalar wdetJ, dXdx[9], F_diff[3] = {0.}; 534c2d90829SJames Wright 535c2d90829SJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 536c2d90829SJames Wright { // Get physical diffusive flux 537c2d90829SJames Wright CeedScalar Grad_qn[15], grad_E_ref[3]; 538c2d90829SJames Wright 539c2d90829SJames Wright GradUnpackN(Q, i, 5, dim, Grad_q, Grad_qn); 540c2d90829SJames Wright CopyN(&Grad_qn[4 * dim], grad_E_ref, dim); 541c2d90829SJames Wright MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff); 542c2d90829SJames Wright ScaleN(F_diff, -context->diffusion_coeff, dim); 543c2d90829SJames Wright } 544c2d90829SJames Wright 545c2d90829SJames Wright CeedScalar F_diff_dXdx[3] = {0.}; 546c2d90829SJames Wright MatVecNM(dXdx, F_diff, dim, dim, CEED_NOTRANSPOSE, F_diff_dXdx); 547c2d90829SJames Wright for (CeedInt k = 0; k < dim; k++) Grad_v[k][i] = -wdetJ * F_diff_dXdx[k]; 548c2d90829SJames Wright } 549c2d90829SJames Wright return 0; 550c2d90829SJames Wright } 551c2d90829SJames Wright 552c2d90829SJames Wright CEED_QFUNCTION(DivDiffusiveFluxVolumeRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 553c2d90829SJames Wright return DivDiffusiveFluxVolumeRHS_AdvDif_Generic(ctx, Q, in, out, 2); 554c2d90829SJames Wright } 555c2d90829SJames Wright 556c2d90829SJames Wright CEED_QFUNCTION(DivDiffusiveFluxVolumeRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 557c2d90829SJames Wright return DivDiffusiveFluxVolumeRHS_AdvDif_Generic(ctx, Q, in, out, 3); 558c2d90829SJames Wright } 559c2d90829SJames Wright 560c2d90829SJames Wright // @brief Boundary integral for RHS of divergence of diffusive flux direct projection 561c2d90829SJames Wright CEED_QFUNCTION_HELPER int DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, 562c2d90829SJames Wright const CeedInt dim) { 563c2d90829SJames Wright const CeedScalar(*Grad_q) = in[0]; 564c2d90829SJames Wright const CeedScalar(*q_data) = in[1]; 565c2d90829SJames Wright CeedScalar(*v) = out[0]; 566c2d90829SJames Wright 567c2d90829SJames Wright AdvectionContext context = (AdvectionContext)ctx; 568c2d90829SJames Wright 569c2d90829SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 570c2d90829SJames Wright CeedScalar wdetJ, normal[3], dXdx[9], F_diff[3] = {0.}; 571c2d90829SJames Wright 572c2d90829SJames Wright QdataBoundaryGradientUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx, normal); 573c2d90829SJames Wright { // Get physical diffusive flux 574c2d90829SJames Wright CeedScalar Grad_qn[15], grad_E_ref[3]; 575c2d90829SJames Wright 576c2d90829SJames Wright GradUnpackN(Q, i, 5, dim, Grad_q, Grad_qn); 577c2d90829SJames Wright CopyN(&Grad_qn[4 * dim], grad_E_ref, dim); 578c2d90829SJames Wright MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff); 579c2d90829SJames Wright ScaleN(F_diff, -context->diffusion_coeff, dim); 580c2d90829SJames Wright } 581c2d90829SJames Wright 582c2d90829SJames Wright v[i] = wdetJ * DotN(F_diff, normal, dim); 583c2d90829SJames Wright } 584c2d90829SJames Wright return 0; 585c2d90829SJames Wright } 586c2d90829SJames Wright 587c2d90829SJames Wright CEED_QFUNCTION(DivDiffusiveFluxBoundaryRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 588c2d90829SJames Wright return DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(ctx, Q, in, out, 2); 589c2d90829SJames Wright } 590c2d90829SJames Wright 591c2d90829SJames Wright CEED_QFUNCTION(DivDiffusiveFluxBoundaryRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 592c2d90829SJames Wright return DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(ctx, Q, in, out, 3); 593c2d90829SJames Wright } 59440b78511SJames Wright 59540b78511SJames Wright // @brief Volume integral for RHS of diffusive flux indirect projection 59640b78511SJames Wright CEED_QFUNCTION_HELPER int DiffusiveFluxRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, 59740b78511SJames Wright const CeedInt dim) { 59840b78511SJames Wright const CeedScalar(*Grad_q) = in[0]; 59940b78511SJames Wright const CeedScalar(*q_data) = in[1]; 60040b78511SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 60140b78511SJames Wright 60240b78511SJames Wright AdvectionContext context = (AdvectionContext)ctx; 60340b78511SJames Wright 60440b78511SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 60540b78511SJames Wright CeedScalar wdetJ, dXdx[9], F_diff[3] = {0.}; 60640b78511SJames Wright 60740b78511SJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 60840b78511SJames Wright { // Get physical diffusive flux 60940b78511SJames Wright CeedScalar Grad_qn[15], grad_E_ref[3]; 61040b78511SJames Wright 61140b78511SJames Wright GradUnpackN(Q, i, 5, dim, Grad_q, Grad_qn); 61240b78511SJames Wright CopyN(&Grad_qn[4 * dim], grad_E_ref, dim); 61340b78511SJames Wright MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff); 61440b78511SJames Wright ScaleN(F_diff, -context->diffusion_coeff, dim); 61540b78511SJames Wright } 61640b78511SJames Wright for (CeedInt k = 0; k < dim; k++) v[k][i] = wdetJ * F_diff[k]; 61740b78511SJames Wright } 61840b78511SJames Wright return 0; 61940b78511SJames Wright } 62040b78511SJames Wright 62140b78511SJames Wright CEED_QFUNCTION(DiffusiveFluxRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 62240b78511SJames Wright return DiffusiveFluxRHS_AdvDif_Generic(ctx, Q, in, out, 2); 62340b78511SJames Wright } 62440b78511SJames Wright 62540b78511SJames Wright CEED_QFUNCTION(DiffusiveFluxRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 62640b78511SJames Wright return DiffusiveFluxRHS_AdvDif_Generic(ctx, Q, in, out, 3); 62740b78511SJames Wright } 628