1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3a515125bSLeila Ghaffari // 4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5a515125bSLeila Ghaffari // 6727da7e7SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7a515125bSLeila Ghaffari 8a515125bSLeila Ghaffari /// @file 9a515125bSLeila Ghaffari /// Advection initial condition and operator for Navier-Stokes example using PETSc 10493642f1SJames Wright #include <ceed.h> 11d0cce58aSJeremy L Thompson #include <math.h> 12a515125bSLeila Ghaffari 13e88b842aSJames Wright #include "advection_types.h" 14ce192147SJames Wright #include "newtonian_state.h" 15ce192147SJames Wright #include "newtonian_types.h" 16e88b842aSJames Wright #include "stabilization_types.h" 171a74fa30SJames Wright #include "utils.h" 181a74fa30SJames Wright 19a515125bSLeila Ghaffari // ***************************************************************************** 209529d636SJames Wright // This QFunction sets the initial conditions and the boundary conditions 219529d636SJames Wright // for two test cases: ROTATION and TRANSLATION 229529d636SJames Wright // 239529d636SJames Wright // -- ROTATION (default) 249529d636SJames Wright // Initial Conditions: 259529d636SJames Wright // Mass Density: 269529d636SJames Wright // Constant mass density of 1.0 279529d636SJames Wright // Momentum Density: 289529d636SJames Wright // Rotational field in x,y 299529d636SJames Wright // Energy Density: 309529d636SJames Wright // Maximum of 1. x0 decreasing linearly to 0. as radial distance 319529d636SJames Wright // increases to (1.-r/rc), then 0. everywhere else 329529d636SJames Wright // 339529d636SJames Wright // Boundary Conditions: 349529d636SJames Wright // Mass Density: 359529d636SJames Wright // 0.0 flux 369529d636SJames Wright // Momentum Density: 379529d636SJames Wright // 0.0 389529d636SJames Wright // Energy Density: 399529d636SJames Wright // 0.0 flux 409529d636SJames Wright // 419529d636SJames Wright // -- TRANSLATION 429529d636SJames Wright // Initial Conditions: 439529d636SJames Wright // Mass Density: 449529d636SJames Wright // Constant mass density of 1.0 459529d636SJames Wright // Momentum Density: 469529d636SJames Wright // Constant rectilinear field in x,y 479529d636SJames Wright // Energy Density: 489529d636SJames Wright // Maximum of 1. x0 decreasing linearly to 0. as radial distance 499529d636SJames Wright // increases to (1.-r/rc), then 0. everywhere else 509529d636SJames Wright // 519529d636SJames Wright // Boundary Conditions: 529529d636SJames Wright // Mass Density: 539529d636SJames Wright // 0.0 flux 549529d636SJames Wright // Momentum Density: 559529d636SJames Wright // 0.0 569529d636SJames Wright // Energy Density: 579529d636SJames Wright // Inflow BCs: 589529d636SJames Wright // E = E_wind 599529d636SJames Wright // Outflow BCs: 609529d636SJames Wright // E = E(boundary) 619529d636SJames Wright // Both In/Outflow BCs for E are applied weakly in the 629529d636SJames Wright // QFunction "Advection2d_Sur" 639529d636SJames Wright // 649529d636SJames Wright // ***************************************************************************** 659529d636SJames Wright 669529d636SJames Wright // ***************************************************************************** 679529d636SJames Wright // This helper function provides the exact, time-dependent solution and IC formulation for 2D advection 689529d636SJames Wright // ***************************************************************************** 699529d636SJames Wright CEED_QFUNCTION_HELPER CeedInt Exact_AdvectionGeneric(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) { 709529d636SJames Wright const SetupContextAdv context = (SetupContextAdv)ctx; 719529d636SJames Wright const CeedScalar rc = context->rc; 729529d636SJames Wright const CeedScalar lx = context->lx; 739529d636SJames Wright const CeedScalar ly = context->ly; 749529d636SJames Wright const CeedScalar lz = dim == 2 ? 0. : context->lz; 759529d636SJames Wright const CeedScalar *wind = context->wind; 769529d636SJames Wright 779529d636SJames Wright const CeedScalar center[3] = {0.5 * lx, 0.5 * ly, 0.5 * lz}; 789529d636SJames Wright const CeedScalar theta = dim == 2 ? M_PI / 3 : M_PI; 799529d636SJames Wright const CeedScalar x0[3] = {center[0] + .25 * lx * cos(theta + time), center[1] + .25 * ly * sin(theta + time), 0.5 * lz}; 809529d636SJames Wright 819529d636SJames Wright const CeedScalar x = X[0], y = X[1], z = dim == 2 ? 0. : X[2]; 829529d636SJames Wright 839529d636SJames Wright CeedScalar r = 0.; 849529d636SJames Wright switch (context->initial_condition_type) { 859529d636SJames Wright case ADVECTIONIC_BUBBLE_SPHERE: 869529d636SJames Wright case ADVECTIONIC_BUBBLE_CYLINDER: 879529d636SJames Wright r = sqrt(Square(x - x0[0]) + Square(y - x0[1]) + Square(z - x0[2])); 889529d636SJames Wright break; 899529d636SJames Wright case ADVECTIONIC_COSINE_HILL: 909529d636SJames Wright r = sqrt(Square(x - center[0]) + Square(y - center[1])); 919529d636SJames Wright break; 929529d636SJames Wright case ADVECTIONIC_SKEW: 939529d636SJames Wright break; 949529d636SJames Wright } 959529d636SJames Wright 969529d636SJames Wright switch (context->wind_type) { 979529d636SJames Wright case WIND_ROTATION: 989529d636SJames Wright q[0] = 1.; 999529d636SJames Wright q[1] = -(y - center[1]); 1009529d636SJames Wright q[2] = (x - center[0]); 1019529d636SJames Wright q[3] = 0; 1029529d636SJames Wright break; 1039529d636SJames Wright case WIND_TRANSLATION: 1049529d636SJames Wright q[0] = 1.; 1059529d636SJames Wright q[1] = wind[0]; 1069529d636SJames Wright q[2] = wind[1]; 1079529d636SJames Wright q[3] = dim == 2 ? 0. : wind[2]; 1089529d636SJames Wright break; 1099529d636SJames Wright default: 1109529d636SJames Wright return 1; 1119529d636SJames Wright } 1129529d636SJames Wright 1139529d636SJames Wright switch (context->initial_condition_type) { 1149529d636SJames Wright case ADVECTIONIC_BUBBLE_SPHERE: 1159529d636SJames Wright case ADVECTIONIC_BUBBLE_CYLINDER: 1169529d636SJames Wright switch (context->bubble_continuity_type) { 1179529d636SJames Wright // original continuous, smooth shape 1189529d636SJames Wright case BUBBLE_CONTINUITY_SMOOTH: 1199529d636SJames Wright q[4] = r <= rc ? (1. - r / rc) : 0.; 1209529d636SJames Wright break; 1219529d636SJames Wright // discontinuous, sharp back half shape 1229529d636SJames Wright case BUBBLE_CONTINUITY_BACK_SHARP: 1239529d636SJames Wright q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) : 0.; 1249529d636SJames Wright break; 1259529d636SJames Wright // attempt to define a finite thickness that will get resolved under grid refinement 1269529d636SJames Wright case BUBBLE_CONTINUITY_THICK: 1279529d636SJames Wright q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) * fmin(1.0, (center[1] - y) / 1.25) : 0.; 1289529d636SJames Wright break; 1299529d636SJames Wright case BUBBLE_CONTINUITY_COSINE: 1309529d636SJames Wright q[4] = r <= rc ? .5 + .5 * cos(r * M_PI / rc) : 0; 1319529d636SJames Wright break; 1329529d636SJames Wright } 1339529d636SJames Wright break; 1349529d636SJames Wright case ADVECTIONIC_COSINE_HILL: { 1359529d636SJames Wright CeedScalar half_width = context->lx / 2; 1369529d636SJames Wright q[4] = r > half_width ? 0. : cos(2 * M_PI * r / half_width + M_PI) + 1.; 1379529d636SJames Wright } break; 1389529d636SJames Wright case ADVECTIONIC_SKEW: { 1399529d636SJames Wright CeedScalar skewed_barrier[3] = {wind[0], wind[1], 0}; 1409529d636SJames Wright CeedScalar inflow_to_point[3] = {x - context->lx / 2, y, 0}; 1419529d636SJames Wright CeedScalar cross_product[3] = {0}; 1429529d636SJames Wright const CeedScalar boundary_threshold = 20 * CEED_EPSILON; 1439529d636SJames Wright Cross3(skewed_barrier, inflow_to_point, cross_product); 1449529d636SJames Wright 1459529d636SJames Wright q[4] = cross_product[2] > boundary_threshold ? 0 : 1; 1469529d636SJames Wright if ((x < boundary_threshold && wind[0] < boundary_threshold) || // outflow at -x boundary 1479529d636SJames Wright (y < boundary_threshold && wind[1] < boundary_threshold) || // outflow at -y boundary 1489529d636SJames Wright (x > context->lx - boundary_threshold && wind[0] > boundary_threshold) || // outflow at +x boundary 1499529d636SJames Wright (y > context->ly - boundary_threshold && wind[1] > boundary_threshold) // outflow at +y boundary 1509529d636SJames Wright ) { 1519529d636SJames Wright q[4] = 0; 1529529d636SJames Wright } 1539529d636SJames Wright } break; 1549529d636SJames Wright } 1559529d636SJames Wright return 0; 1569529d636SJames Wright } 1579529d636SJames Wright 1589529d636SJames Wright // ***************************************************************************** 159a515125bSLeila Ghaffari // This QFunction sets the initial conditions for 3D advection 160a515125bSLeila Ghaffari // ***************************************************************************** 1612b916ea7SJeremy L Thompson CEED_QFUNCTION(ICsAdvection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 162a515125bSLeila Ghaffari const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 163a515125bSLeila Ghaffari CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 164a515125bSLeila Ghaffari 1653d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 166a515125bSLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 167139613f2SLeila Ghaffari CeedScalar q[5] = {0.}; 168a515125bSLeila Ghaffari 1690b3a1fabSJames Wright Exact_AdvectionGeneric(3, 0., x, 5, q, ctx); 170a515125bSLeila Ghaffari for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 1710b3a1fabSJames Wright } 172a515125bSLeila Ghaffari return 0; 173a515125bSLeila Ghaffari } 174a515125bSLeila Ghaffari 175a515125bSLeila Ghaffari // ***************************************************************************** 1769529d636SJames Wright // This QFunction sets the initial conditions for 2D advection 177a515125bSLeila Ghaffari // ***************************************************************************** 1789529d636SJames Wright CEED_QFUNCTION(ICsAdvection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 1799529d636SJames Wright const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 1809529d636SJames Wright CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 1819529d636SJames Wright const SetupContextAdv context = (SetupContextAdv)ctx; 1829529d636SJames Wright 1839529d636SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 1849529d636SJames Wright const CeedScalar x[] = {X[0][i], X[1][i]}; 1859529d636SJames Wright CeedScalar q[5] = {0.}; 1869529d636SJames Wright 1879529d636SJames Wright Exact_AdvectionGeneric(2, context->time, x, 5, q, ctx); 1889529d636SJames Wright for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 1899529d636SJames Wright } 190a515125bSLeila Ghaffari return 0; 191a515125bSLeila Ghaffari } 192a515125bSLeila Ghaffari 1939529d636SJames Wright CEED_QFUNCTION_HELPER void QdataUnpack_ND(CeedInt N, CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar *dXdx) { 19485efd435SJames Wright // Cannot directly use QdataUnpack* helper functions due to SYCL online compiler incompatabilities 1959529d636SJames Wright switch (N) { 1969529d636SJames Wright case 2: 19785efd435SJames Wright StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ); 19885efd435SJames Wright StoredValuesUnpack(Q, i, 1, 4, q_data, dXdx); 1999529d636SJames Wright break; 2009529d636SJames Wright case 3: 20185efd435SJames Wright StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ); 20285efd435SJames Wright StoredValuesUnpack(Q, i, 1, 9, q_data, dXdx); 2039529d636SJames Wright break; 2049529d636SJames Wright } 2059529d636SJames Wright } 2069529d636SJames Wright 2079529d636SJames Wright CEED_QFUNCTION_HELPER int QdataBoundaryUnpack_ND(CeedInt N, CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar *dXdx, 2089529d636SJames Wright CeedScalar *normal) { 20985efd435SJames Wright // Cannot directly use QdataBoundaryUnpack* helper functions due to SYCL online compiler incompatabilities 2109529d636SJames Wright switch (N) { 2119529d636SJames Wright case 2: 21285efd435SJames Wright if (wdetJ) StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ); 21385efd435SJames Wright if (normal) StoredValuesUnpack(Q, i, 1, 2, q_data, normal); 2149529d636SJames Wright break; 2159529d636SJames Wright case 3: 21685efd435SJames Wright if (wdetJ) StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ); 21785efd435SJames Wright if (normal) StoredValuesUnpack(Q, i, 1, 3, q_data, normal); 21885efd435SJames Wright if (dXdx) StoredValuesUnpack(Q, i, 4, 6, q_data, (CeedScalar *)dXdx); 2199529d636SJames Wright break; 2209529d636SJames Wright } 2219529d636SJames Wright return CEED_ERROR_SUCCESS; 2229529d636SJames Wright } 2239529d636SJames Wright 2249529d636SJames Wright CEED_QFUNCTION_HELPER void StatePhysicalGradientFromReference_ND(CeedInt N, CeedInt Q, CeedInt i, NewtonianIdealGasContext gas, State s, 2259529d636SJames Wright StateVariable state_var, const CeedScalar *grad_q, const CeedScalar *dXdx, 2269529d636SJames Wright State *grad_s) { 2279529d636SJames Wright switch (N) { 2289529d636SJames Wright case 2: { 2299529d636SJames Wright for (CeedInt k = 0; k < 2; k++) { 2309529d636SJames Wright CeedScalar dqi[5]; 2319529d636SJames Wright for (CeedInt j = 0; j < 5; j++) { 2329529d636SJames 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]; 2339529d636SJames Wright } 2349529d636SJames Wright grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var); 2359529d636SJames Wright } 2369529d636SJames Wright CeedScalar U[5] = {0.}; 2379529d636SJames Wright grad_s[2] = StateFromU(gas, U); 2389529d636SJames Wright } break; 2399529d636SJames Wright case 3: 24085efd435SJames Wright // Cannot directly use StatePhysicalGradientFromReference helper functions due to SYCL online compiler incompatabilities 24185efd435SJames Wright for (CeedInt k = 0; k < 3; k++) { 24285efd435SJames Wright CeedScalar dqi[5]; 24385efd435SJames Wright for (CeedInt j = 0; j < 5; j++) { 24485efd435SJames 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] + 24585efd435SJames Wright grad_q[(Q * 5) * 2 + Q * j + i] * dXdx[2 * N + k]; 24685efd435SJames Wright } 24785efd435SJames Wright grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var); 24885efd435SJames Wright } 2499529d636SJames Wright break; 2509529d636SJames Wright } 2519529d636SJames Wright } 2529529d636SJames Wright 253a78efa86SJames Wright // @brief Calculate the stabilization constant \tau 254a78efa86SJames Wright CEED_QFUNCTION_HELPER CeedScalar Tau(AdvectionContext context, const State s, const CeedScalar *dXdx, CeedInt dim) { 255a78efa86SJames Wright switch (context->stabilization_tau) { 256a78efa86SJames Wright case STAB_TAU_CTAU: { 257a78efa86SJames Wright CeedScalar uX[3] = {0.}; 258a78efa86SJames Wright 259a78efa86SJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 260a78efa86SJames Wright return context->CtauS / sqrt(DotN(uX, uX, dim)); 261a78efa86SJames Wright } break; 262a78efa86SJames Wright case STAB_TAU_ADVDIFF_SHAKIB: { 263a78efa86SJames Wright CeedScalar gijd_mat[9] = {0.}, gij_uj[3] = {0.}; 264a78efa86SJames Wright 265a78efa86SJames Wright MatMatN(dXdx, dXdx, dim, CEED_TRANSPOSE, CEED_NOTRANSPOSE, gijd_mat); 266a78efa86SJames Wright MatVecNM(gijd_mat, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, gij_uj); 267a78efa86SJames Wright return 1 / sqrt(Square(2 * context->Ctau_t / context->dt) + DotN(s.Y.velocity, gij_uj, dim) * context->Ctau_a); 268a78efa86SJames Wright } break; 269a78efa86SJames Wright default: 270a78efa86SJames Wright return 0.; 271a78efa86SJames Wright } 272a78efa86SJames Wright } 273a78efa86SJames Wright 2749529d636SJames Wright // ***************************************************************************** 2759529d636SJames Wright // This QFunction implements Advection for implicit time stepping method 2769529d636SJames Wright // ***************************************************************************** 2779529d636SJames Wright CEED_QFUNCTION_HELPER void IFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 2789529d636SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 2799529d636SJames Wright const CeedScalar(*grad_q) = in[1]; 2809529d636SJames Wright const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 2819529d636SJames Wright const CeedScalar(*q_data) = in[3]; 2829529d636SJames Wright 2839529d636SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 2849529d636SJames Wright CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 2859529d636SJames Wright CeedScalar *jac_data = out[2]; 2869529d636SJames Wright 2879529d636SJames Wright AdvectionContext context = (AdvectionContext)ctx; 2889529d636SJames Wright const CeedScalar zeros[14] = {0.}; 2899529d636SJames Wright NewtonianIdealGasContext gas; 2909529d636SJames Wright struct NewtonianIdealGasContext_ gas_struct = {0}; 2919529d636SJames Wright gas = &gas_struct; 2929529d636SJames Wright 2939529d636SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 2949529d636SJames Wright const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]}; 2959529d636SJames Wright const State s = StateFromU(gas, qi); 2969529d636SJames Wright 2979529d636SJames Wright CeedScalar wdetJ, dXdx[9]; 2989529d636SJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 2999529d636SJames Wright State grad_s[3]; 3009529d636SJames Wright StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s); 3019529d636SJames Wright 3029529d636SJames Wright const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total}; 3039529d636SJames Wright 3049529d636SJames Wright for (CeedInt f = 0; f < 4; f++) { 3059529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum 3069529d636SJames Wright v[f][i] = wdetJ * q_dot[f][i]; // K Mass/transient term 3079529d636SJames Wright } 3089529d636SJames Wright 3099529d636SJames Wright CeedScalar div_u = 0; 3109529d636SJames Wright for (CeedInt j = 0; j < dim; j++) { 3119529d636SJames Wright for (CeedInt k = 0; k < dim; k++) { 3129529d636SJames Wright div_u += grad_s[k].Y.velocity[j]; 3139529d636SJames Wright } 3149529d636SJames Wright } 3159529d636SJames Wright CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim); 3169529d636SJames Wright CeedScalar strong_res = q_dot[4][i] + strong_conv; 3179529d636SJames Wright 3189529d636SJames Wright v[4][i] = wdetJ * q_dot[4][i]; // transient part (ALWAYS) 3199529d636SJames Wright 3209529d636SJames Wright CeedScalar uX[3] = {0.}; 3219529d636SJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 3229529d636SJames Wright 3239529d636SJames Wright if (context->strong_form) { // Strong Galerkin convection term: v div(E u) 3249529d636SJames Wright v[4][i] += wdetJ * strong_conv; 3259529d636SJames Wright } else { // Weak Galerkin convection term: -dv \cdot (E u) 3269529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = -wdetJ * s.U.E_total * uX[j]; 3279529d636SJames Wright } 3289529d636SJames Wright 329*c8d249deSJames Wright { // Diffusion 330*c8d249deSJames Wright CeedScalar Fe[3], Fe_dXdx[3] = {0.}; 331*c8d249deSJames Wright 332*c8d249deSJames Wright for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total; 333*c8d249deSJames Wright MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx); 334*c8d249deSJames Wright for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] -= wdetJ * Fe_dXdx[k]; 335*c8d249deSJames Wright } 336*c8d249deSJames Wright 337a78efa86SJames Wright const CeedScalar TauS = Tau(context, s, dXdx, dim); 3389529d636SJames Wright for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) { 3399529d636SJames Wright case STAB_NONE: 3409529d636SJames Wright break; 3419529d636SJames Wright case STAB_SU: 3429529d636SJames Wright grad_v[j][4][i] += wdetJ * TauS * strong_conv * uX[j]; 3439529d636SJames Wright break; 3449529d636SJames Wright case STAB_SUPG: 3459529d636SJames Wright grad_v[j][4][i] += wdetJ * TauS * strong_res * uX[j]; 3469529d636SJames Wright break; 3479529d636SJames Wright } 3489529d636SJames Wright StoredValuesPack(Q, i, 0, 14, zeros, jac_data); 3499529d636SJames Wright } 3509529d636SJames Wright } 3519529d636SJames Wright 3522b916ea7SJeremy L Thompson CEED_QFUNCTION(IFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 353bd4b5413SJames Wright IFunction_AdvectionGeneric(ctx, Q, in, out, 3); 354a515125bSLeila Ghaffari return 0; 355a515125bSLeila Ghaffari } 356a515125bSLeila Ghaffari 3579529d636SJames Wright CEED_QFUNCTION(IFunction_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 3589529d636SJames Wright IFunction_AdvectionGeneric(ctx, Q, in, out, 2); 3599529d636SJames Wright return 0; 3609529d636SJames Wright } 3619529d636SJames Wright 362a78efa86SJames Wright CEED_QFUNCTION_HELPER void MassFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 363a78efa86SJames Wright const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 364a78efa86SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 365a78efa86SJames Wright const CeedScalar(*q_data) = in[2]; 366a78efa86SJames Wright 367a78efa86SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 368a78efa86SJames Wright CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 369a78efa86SJames Wright 370a78efa86SJames Wright AdvectionContext context = (AdvectionContext)ctx; 371a78efa86SJames Wright struct NewtonianIdealGasContext_ gas_struct = {0}; 372a78efa86SJames Wright NewtonianIdealGasContext gas = &gas_struct; 373a78efa86SJames Wright 374a78efa86SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 375a78efa86SJames Wright const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]}; 376a78efa86SJames Wright const State s = StateFromU(gas, qi); 377a78efa86SJames Wright CeedScalar wdetJ, dXdx[9]; 378a78efa86SJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 379a78efa86SJames Wright 380a78efa86SJames Wright for (CeedInt f = 0; f < 4; f++) { 381a78efa86SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum 382a78efa86SJames Wright v[f][i] = wdetJ * q_dot[f][i]; // K Mass/transient term 383a78efa86SJames Wright } 384a78efa86SJames Wright 385a78efa86SJames Wright // Unstabilized mass term 386a78efa86SJames Wright v[4][i] = wdetJ * q_dot[4][i]; 387a78efa86SJames Wright 388a78efa86SJames Wright // Stabilized mass term 389a78efa86SJames Wright CeedScalar uX[3] = {0.}; 390a78efa86SJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 391a78efa86SJames Wright const CeedScalar TauS = Tau(context, s, dXdx, dim); 392a78efa86SJames Wright for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) { 393a78efa86SJames Wright case STAB_NONE: 394a78efa86SJames Wright case STAB_SU: 395a78efa86SJames Wright grad_v[j][4][i] = 0; 396a78efa86SJames Wright break; // These should be run with the unstabilized mass matrix anyways 397a78efa86SJames Wright case STAB_SUPG: 398a78efa86SJames Wright grad_v[j][4][i] = wdetJ * TauS * q_dot[4][i] * uX[j]; 399a78efa86SJames Wright break; 400a78efa86SJames Wright } 401a78efa86SJames Wright } 402a78efa86SJames Wright } 403a78efa86SJames Wright 404a78efa86SJames Wright CEED_QFUNCTION(MassFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 405a78efa86SJames Wright MassFunction_AdvectionGeneric(ctx, Q, in, out, 3); 406a78efa86SJames Wright return 0; 407a78efa86SJames Wright } 408a78efa86SJames Wright 409a78efa86SJames Wright CEED_QFUNCTION(MassFunction_Advection2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 410a78efa86SJames Wright MassFunction_AdvectionGeneric(ctx, Q, in, out, 2); 411a78efa86SJames Wright return 0; 412a78efa86SJames Wright } 413a78efa86SJames Wright 4149529d636SJames Wright // ***************************************************************************** 4159529d636SJames Wright // This QFunction implements Advection for explicit time stepping method 4169529d636SJames Wright // ***************************************************************************** 4179529d636SJames Wright CEED_QFUNCTION_HELPER void RHSFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 4189529d636SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 4199529d636SJames Wright const CeedScalar(*grad_q) = in[1]; 4209529d636SJames Wright const CeedScalar(*q_data) = in[2]; 4219529d636SJames Wright 4229529d636SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 4239529d636SJames Wright CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 4249529d636SJames Wright 4259529d636SJames Wright AdvectionContext context = (AdvectionContext)ctx; 4269529d636SJames Wright struct NewtonianIdealGasContext_ gas_struct = {0}; 427a78efa86SJames Wright NewtonianIdealGasContext gas = &gas_struct; 4289529d636SJames Wright 4299529d636SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 4309529d636SJames Wright const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]}; 4319529d636SJames Wright const State s = StateFromU(gas, qi); 4329529d636SJames Wright 4339529d636SJames Wright CeedScalar wdetJ, dXdx[9]; 4349529d636SJames Wright QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 4359529d636SJames Wright State grad_s[3]; 4369529d636SJames Wright StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s); 4379529d636SJames Wright 4389529d636SJames Wright const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total}; 4399529d636SJames Wright 4409529d636SJames Wright for (CeedInt f = 0; f < 4; f++) { 4419529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum 4429529d636SJames Wright v[f][i] = 0.; 4439529d636SJames Wright } 4449529d636SJames Wright 4459529d636SJames Wright CeedScalar div_u = 0; 4469529d636SJames Wright for (CeedInt j = 0; j < dim; j++) { 4479529d636SJames Wright for (CeedInt k = 0; k < dim; k++) { 4489529d636SJames Wright div_u += grad_s[k].Y.velocity[j]; 4499529d636SJames Wright } 4509529d636SJames Wright } 4519529d636SJames Wright CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim); 4529529d636SJames Wright 4539529d636SJames Wright CeedScalar uX[3] = {0.}; 4549529d636SJames Wright MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 4559529d636SJames Wright 4569529d636SJames Wright if (context->strong_form) { // Strong Galerkin convection term: v div(E u) 4579529d636SJames Wright v[4][i] = -wdetJ * strong_conv; 4589529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = 0; 4599529d636SJames Wright } else { // Weak Galerkin convection term: -dv \cdot (E u) 4609529d636SJames Wright for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = wdetJ * s.U.E_total * uX[j]; 4619529d636SJames Wright v[4][i] = 0.; 4629529d636SJames Wright } 4639529d636SJames Wright 464*c8d249deSJames Wright { // Diffusion 465*c8d249deSJames Wright CeedScalar Fe[3], Fe_dXdx[3] = {0.}; 466*c8d249deSJames Wright 467*c8d249deSJames Wright for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total; 468*c8d249deSJames Wright MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx); 469*c8d249deSJames Wright for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] += wdetJ * Fe_dXdx[k]; 470*c8d249deSJames Wright } 471*c8d249deSJames Wright 472a78efa86SJames Wright const CeedScalar TauS = Tau(context, s, dXdx, dim); 4739529d636SJames Wright for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) { 4749529d636SJames Wright case STAB_NONE: 4759529d636SJames Wright break; 4769529d636SJames Wright case STAB_SU: 4779529d636SJames Wright case STAB_SUPG: 4789d860eefSJames Wright grad_v[j][4][i] -= wdetJ * TauS * strong_conv * uX[j]; 4799529d636SJames Wright break; 4809529d636SJames Wright } 4819529d636SJames Wright } 4829529d636SJames Wright } 4839529d636SJames Wright 4849529d636SJames Wright CEED_QFUNCTION(RHS_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 4859529d636SJames Wright RHSFunction_AdvectionGeneric(ctx, Q, in, out, 3); 4869529d636SJames Wright return 0; 4879529d636SJames Wright } 4889529d636SJames Wright 4899529d636SJames Wright CEED_QFUNCTION(RHS_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 4909529d636SJames Wright RHSFunction_AdvectionGeneric(ctx, Q, in, out, 2); 4919529d636SJames Wright return 0; 4929529d636SJames Wright } 4939529d636SJames Wright 4949529d636SJames Wright // ***************************************************************************** 4959529d636SJames Wright // This QFunction implements consistent outflow and inflow BCs 4969529d636SJames Wright // for advection 4979529d636SJames Wright // 4989529d636SJames Wright // Inflow and outflow faces are determined based on sign(dot(wind, normal)): 4999529d636SJames Wright // sign(dot(wind, normal)) > 0 : outflow BCs 5009529d636SJames Wright // sign(dot(wind, normal)) < 0 : inflow BCs 5019529d636SJames Wright // 5029529d636SJames Wright // Outflow BCs: 5039529d636SJames Wright // The validity of the weak form of the governing equations is extended to the outflow and the current values of E are applied. 5049529d636SJames Wright // 5059529d636SJames Wright // Inflow BCs: 5069529d636SJames Wright // A prescribed Total Energy (E_wind) is applied weakly. 5079529d636SJames Wright // ***************************************************************************** 5089529d636SJames Wright CEED_QFUNCTION(Advection_InOutFlowGeneric)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 5099529d636SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 5109529d636SJames Wright const CeedScalar(*q_data_sur) = in[2]; 5119529d636SJames Wright 5129529d636SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 5139529d636SJames Wright AdvectionContext context = (AdvectionContext)ctx; 5149529d636SJames Wright const CeedScalar E_wind = context->E_wind; 5159529d636SJames Wright const CeedScalar strong_form = context->strong_form; 5169529d636SJames Wright const bool is_implicit = context->implicit; 5179529d636SJames Wright 5189529d636SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 5199529d636SJames Wright const CeedScalar rho = q[0][i]; 5209529d636SJames Wright const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 5219529d636SJames Wright const CeedScalar E = q[4][i]; 5229529d636SJames Wright 5239529d636SJames Wright CeedScalar wdetJb, norm[3]; 5249529d636SJames Wright QdataBoundaryUnpack_ND(dim, Q, i, q_data_sur, &wdetJb, NULL, norm); 5259529d636SJames Wright wdetJb *= is_implicit ? -1. : 1.; 5269529d636SJames Wright 5279529d636SJames Wright const CeedScalar u_normal = DotN(norm, u, dim); 5289529d636SJames Wright 5299529d636SJames Wright // No Change in density or momentum 5309529d636SJames Wright for (CeedInt j = 0; j < 4; j++) { 5319529d636SJames Wright v[j][i] = 0; 5329529d636SJames Wright } 5339529d636SJames Wright // Implementing in/outflow BCs 5349529d636SJames Wright if (u_normal > 0) { // outflow 5359529d636SJames Wright v[4][i] = -(1 - strong_form) * wdetJb * E * u_normal; 5369529d636SJames Wright } else { // inflow 5379529d636SJames Wright v[4][i] = -(1 - strong_form) * wdetJb * E_wind * u_normal; 5389529d636SJames Wright } 5399529d636SJames Wright } 5409529d636SJames Wright return 0; 5419529d636SJames Wright } 5429529d636SJames Wright 5432b916ea7SJeremy L Thompson CEED_QFUNCTION(Advection_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 5448dba1efaSJames Wright Advection_InOutFlowGeneric(ctx, Q, in, out, 3); 545a515125bSLeila Ghaffari return 0; 546a515125bSLeila Ghaffari } 547a515125bSLeila Ghaffari 5489529d636SJames Wright CEED_QFUNCTION(Advection2d_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 5499529d636SJames Wright Advection_InOutFlowGeneric(ctx, Q, in, out, 2); 5509529d636SJames Wright return 0; 5519529d636SJames Wright } 552