1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 /// @file 5 /// Advection initial condition and operator for Navier-Stokes example using PETSc 6 #include <ceed.h> 7 #include <math.h> 8 9 #include "advection_types.h" 10 #include "newtonian_state.h" 11 #include "newtonian_types.h" 12 #include "stabilization_types.h" 13 #include "utils.h" 14 15 // ***************************************************************************** 16 // This QFunction sets the initial conditions and the boundary conditions 17 // for two test cases: ROTATION and TRANSLATION 18 // 19 // -- ROTATION (default) 20 // Initial Conditions: 21 // Mass Density: 22 // Constant mass density of 1.0 23 // Momentum Density: 24 // Rotational field in x,y 25 // Energy Density: 26 // Maximum of 1. x0 decreasing linearly to 0. as radial distance 27 // increases to (1.-r/rc), then 0. everywhere else 28 // 29 // Boundary Conditions: 30 // Mass Density: 31 // 0.0 flux 32 // Momentum Density: 33 // 0.0 34 // Energy Density: 35 // 0.0 flux 36 // 37 // -- TRANSLATION 38 // Initial Conditions: 39 // Mass Density: 40 // Constant mass density of 1.0 41 // Momentum Density: 42 // Constant rectilinear field in x,y 43 // Energy Density: 44 // Maximum of 1. x0 decreasing linearly to 0. as radial distance 45 // increases to (1.-r/rc), then 0. everywhere else 46 // 47 // Boundary Conditions: 48 // Mass Density: 49 // 0.0 flux 50 // Momentum Density: 51 // 0.0 52 // Energy Density: 53 // Inflow BCs: 54 // E = E_wind 55 // Outflow BCs: 56 // E = E(boundary) 57 // Both In/Outflow BCs for E are applied weakly in the 58 // QFunction "Advection2d_Sur" 59 // 60 // ***************************************************************************** 61 62 // ***************************************************************************** 63 // This helper function provides the exact, time-dependent solution and IC formulation for 2D advection 64 // ***************************************************************************** 65 CEED_QFUNCTION_HELPER CeedInt Exact_AdvectionGeneric(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) { 66 const SetupContextAdv context = (SetupContextAdv)ctx; 67 const CeedScalar rc = context->rc; 68 const CeedScalar lx = context->lx; 69 const CeedScalar ly = context->ly; 70 const CeedScalar lz = dim == 2 ? 0. : context->lz; 71 const CeedScalar *wind = context->wind; 72 73 const CeedScalar center[3] = {0.5 * lx, 0.5 * ly, 0.5 * lz}; 74 const CeedScalar theta = dim == 2 ? M_PI / 3 : M_PI; 75 const CeedScalar x0[3] = {center[0] + .25 * lx * cos(theta + time), center[1] + .25 * ly * sin(theta + time), 0.5 * lz}; 76 77 const CeedScalar x = X[0], y = X[1], z = dim == 2 ? 0. : X[2]; 78 79 CeedScalar r = 0.; 80 switch (context->initial_condition_type) { 81 case ADVECTIONIC_BUBBLE_SPHERE: 82 case ADVECTIONIC_BUBBLE_CYLINDER: 83 r = sqrt(Square(x - x0[0]) + Square(y - x0[1]) + Square(z - x0[2])); 84 break; 85 case ADVECTIONIC_COSINE_HILL: 86 r = sqrt(Square(x - center[0]) + Square(y - center[1])); 87 break; 88 case ADVECTIONIC_SKEW: 89 break; 90 } 91 92 switch (context->wind_type) { 93 case WIND_ROTATION: 94 q[0] = 1.; 95 q[1] = -(y - center[1]); 96 q[2] = (x - center[0]); 97 q[3] = 0; 98 break; 99 case WIND_TRANSLATION: 100 q[0] = 1.; 101 q[1] = wind[0]; 102 q[2] = wind[1]; 103 q[3] = dim == 2 ? 0. : wind[2]; 104 break; 105 default: 106 return 1; 107 } 108 109 switch (context->initial_condition_type) { 110 case ADVECTIONIC_BUBBLE_SPHERE: 111 case ADVECTIONIC_BUBBLE_CYLINDER: 112 switch (context->bubble_continuity_type) { 113 // original continuous, smooth shape 114 case BUBBLE_CONTINUITY_SMOOTH: 115 q[4] = r <= rc ? (1. - r / rc) : 0.; 116 break; 117 // discontinuous, sharp back half shape 118 case BUBBLE_CONTINUITY_BACK_SHARP: 119 q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) : 0.; 120 break; 121 // attempt to define a finite thickness that will get resolved under grid refinement 122 case BUBBLE_CONTINUITY_THICK: 123 q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) * fmin(1.0, (center[1] - y) / 1.25) : 0.; 124 break; 125 case BUBBLE_CONTINUITY_COSINE: 126 q[4] = r <= rc ? .5 + .5 * cos(r * M_PI / rc) : 0; 127 break; 128 } 129 break; 130 case ADVECTIONIC_COSINE_HILL: { 131 CeedScalar half_width = context->lx / 2; 132 q[4] = r > half_width ? 0. : cos(2 * M_PI * r / half_width + M_PI) + 1.; 133 } break; 134 case ADVECTIONIC_SKEW: { 135 CeedScalar skewed_barrier[3] = {wind[0], wind[1], 0}; 136 CeedScalar inflow_to_point[3] = {x - context->lx / 2, y, 0}; 137 CeedScalar cross_product[3] = {0}; 138 const CeedScalar boundary_threshold = 20 * CEED_EPSILON; 139 Cross3(skewed_barrier, inflow_to_point, cross_product); 140 141 q[4] = cross_product[2] > boundary_threshold ? 0 : 1; 142 if ((x < boundary_threshold && wind[0] < boundary_threshold) || // outflow at -x boundary 143 (y < boundary_threshold && wind[1] < boundary_threshold) || // outflow at -y boundary 144 (x > context->lx - boundary_threshold && wind[0] > boundary_threshold) || // outflow at +x boundary 145 (y > context->ly - boundary_threshold && wind[1] > boundary_threshold) // outflow at +y boundary 146 ) { 147 q[4] = 0; 148 } 149 } break; 150 } 151 return 0; 152 } 153 154 // ***************************************************************************** 155 // This QFunction sets the initial conditions for 3D advection 156 // ***************************************************************************** 157 CEED_QFUNCTION(ICsAdvection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 158 const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 159 CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 160 161 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 162 const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 163 CeedScalar q[5] = {0.}; 164 165 Exact_AdvectionGeneric(3, 0., x, 5, q, ctx); 166 for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 167 } 168 return 0; 169 } 170 171 // ***************************************************************************** 172 // This QFunction sets the initial conditions for 2D advection 173 // ***************************************************************************** 174 CEED_QFUNCTION(ICsAdvection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 175 const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 176 CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 177 const SetupContextAdv context = (SetupContextAdv)ctx; 178 179 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 180 const CeedScalar x[] = {X[0][i], X[1][i]}; 181 CeedScalar q[5] = {0.}; 182 183 Exact_AdvectionGeneric(2, context->time, x, 5, q, ctx); 184 for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 185 } 186 return 0; 187 } 188 189 CEED_QFUNCTION_HELPER void StatePhysicalGradientFromReference_ND(CeedInt N, CeedInt Q, CeedInt i, NewtonianIdealGasContext gas, State s, 190 StateVariable state_var, const CeedScalar *grad_q, const CeedScalar *dXdx, 191 State *grad_s) { 192 switch (N) { 193 case 2: { 194 for (CeedInt k = 0; k < 2; k++) { 195 CeedScalar dqi[5]; 196 for (CeedInt j = 0; j < 5; j++) { 197 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]; 198 } 199 grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var); 200 } 201 CeedScalar U[5] = {0.}; 202 grad_s[2] = StateFromU(gas, U); 203 } break; 204 case 3: 205 // Cannot directly use StatePhysicalGradientFromReference helper functions due to SYCL online compiler incompatabilities 206 for (CeedInt k = 0; k < 3; k++) { 207 CeedScalar dqi[5]; 208 for (CeedInt j = 0; j < 5; j++) { 209 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] + 210 grad_q[(Q * 5) * 2 + Q * j + i] * dXdx[2 * N + k]; 211 } 212 grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var); 213 } 214 break; 215 } 216 } 217 218 // @brief Calculate the stabilization constant \tau 219 CEED_QFUNCTION_HELPER CeedScalar Tau(AdvectionContext context, const State s, const CeedScalar *dXdx, CeedInt dim) { 220 switch (context->stabilization_tau) { 221 case STAB_TAU_CTAU: { 222 CeedScalar uX[3] = {0.}; 223 224 MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 225 return context->CtauS / sqrt(DotN(uX, uX, dim)); 226 } break; 227 case STAB_TAU_ADVDIFF_SHAKIB: { 228 CeedScalar gijd_mat[9] = {0.}, gij_uj[3] = {0.}; 229 230 MatMatN(dXdx, dXdx, dim, CEED_TRANSPOSE, CEED_NOTRANSPOSE, gijd_mat); 231 MatVecNM(gijd_mat, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, gij_uj); 232 return 1 / sqrt(Square(2 * context->Ctau_t / context->dt) + DotN(s.Y.velocity, gij_uj, dim) * context->Ctau_a); 233 } break; 234 default: 235 return 0.; 236 } 237 } 238 239 // ***************************************************************************** 240 // This QFunction implements Advection for implicit time stepping method 241 // ***************************************************************************** 242 CEED_QFUNCTION_HELPER void IFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 243 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 244 const CeedScalar(*grad_q) = in[1]; 245 const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 246 const CeedScalar(*q_data) = in[3]; 247 248 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 249 CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 250 251 AdvectionContext context = (AdvectionContext)ctx; 252 NewtonianIdealGasContext gas; 253 struct NewtonianIdealGasContext_ gas_struct = {0}; 254 gas = &gas_struct; 255 256 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 257 const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]}; 258 const State s = StateFromU(gas, qi); 259 260 CeedScalar wdetJ, dXdx[9]; 261 QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 262 State grad_s[3]; 263 StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s); 264 265 const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total}; 266 267 for (CeedInt f = 0; f < 4; f++) { 268 for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum 269 v[f][i] = wdetJ * q_dot[f][i]; // K Mass/transient term 270 } 271 272 CeedScalar div_u = 0; 273 for (CeedInt j = 0; j < dim; j++) { 274 for (CeedInt k = 0; k < dim; k++) { 275 div_u += grad_s[k].Y.velocity[j]; 276 } 277 } 278 CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim); 279 CeedScalar strong_res = q_dot[4][i] + strong_conv; 280 281 v[4][i] = wdetJ * q_dot[4][i]; // transient part (ALWAYS) 282 283 CeedScalar uX[3] = {0.}; 284 MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 285 286 if (context->strong_form) { // Strong Galerkin convection term: v div(E u) 287 v[4][i] += wdetJ * strong_conv; 288 } else { // Weak Galerkin convection term: -dv \cdot (E u) 289 for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = -wdetJ * s.U.E_total * uX[j]; 290 } 291 292 { // Diffusion 293 CeedScalar Fe[3], Fe_dXdx[3] = {0.}; 294 295 for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total; 296 MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx); 297 for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] -= wdetJ * Fe_dXdx[k]; 298 } 299 300 const CeedScalar TauS = Tau(context, s, dXdx, dim); 301 for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) { 302 case STAB_NONE: 303 break; 304 case STAB_SU: 305 grad_v[j][4][i] += wdetJ * TauS * strong_conv * uX[j]; 306 break; 307 case STAB_SUPG: 308 grad_v[j][4][i] += wdetJ * TauS * strong_res * uX[j]; 309 break; 310 } 311 } 312 } 313 314 CEED_QFUNCTION(IFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 315 IFunction_AdvectionGeneric(ctx, Q, in, out, 3); 316 return 0; 317 } 318 319 CEED_QFUNCTION(IFunction_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 320 IFunction_AdvectionGeneric(ctx, Q, in, out, 2); 321 return 0; 322 } 323 324 CEED_QFUNCTION_HELPER void MassFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 325 const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 326 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 327 const CeedScalar(*q_data) = in[2]; 328 329 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 330 CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 331 332 AdvectionContext context = (AdvectionContext)ctx; 333 struct NewtonianIdealGasContext_ gas_struct = {0}; 334 NewtonianIdealGasContext gas = &gas_struct; 335 336 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 337 const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]}; 338 const State s = StateFromU(gas, qi); 339 CeedScalar wdetJ, dXdx[9]; 340 QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 341 342 for (CeedInt f = 0; f < 4; f++) { 343 for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum 344 v[f][i] = wdetJ * q_dot[f][i]; // K Mass/transient term 345 } 346 347 // Unstabilized mass term 348 v[4][i] = wdetJ * q_dot[4][i]; 349 350 // Stabilized mass term 351 CeedScalar uX[3] = {0.}; 352 MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 353 const CeedScalar TauS = Tau(context, s, dXdx, dim); 354 for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) { 355 case STAB_NONE: 356 case STAB_SU: 357 grad_v[j][4][i] = 0; 358 break; // These should be run with the unstabilized mass matrix anyways 359 case STAB_SUPG: 360 grad_v[j][4][i] = wdetJ * TauS * q_dot[4][i] * uX[j]; 361 break; 362 } 363 } 364 } 365 366 CEED_QFUNCTION(MassFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 367 MassFunction_AdvectionGeneric(ctx, Q, in, out, 3); 368 return 0; 369 } 370 371 CEED_QFUNCTION(MassFunction_Advection2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 372 MassFunction_AdvectionGeneric(ctx, Q, in, out, 2); 373 return 0; 374 } 375 376 // ***************************************************************************** 377 // This QFunction implements Advection for explicit time stepping method 378 // ***************************************************************************** 379 CEED_QFUNCTION_HELPER void RHSFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 380 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 381 const CeedScalar(*grad_q) = in[1]; 382 const CeedScalar(*q_data) = in[2]; 383 384 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 385 CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 386 387 AdvectionContext context = (AdvectionContext)ctx; 388 struct NewtonianIdealGasContext_ gas_struct = {0}; 389 NewtonianIdealGasContext gas = &gas_struct; 390 391 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 392 const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]}; 393 const State s = StateFromU(gas, qi); 394 395 CeedScalar wdetJ, dXdx[9]; 396 QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 397 State grad_s[3]; 398 StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s); 399 400 const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total}; 401 402 for (CeedInt f = 0; f < 4; f++) { 403 for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum 404 v[f][i] = 0.; 405 } 406 407 CeedScalar div_u = 0; 408 for (CeedInt j = 0; j < dim; j++) { 409 for (CeedInt k = 0; k < dim; k++) { 410 div_u += grad_s[k].Y.velocity[j]; 411 } 412 } 413 CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim); 414 415 CeedScalar uX[3] = {0.}; 416 MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 417 418 if (context->strong_form) { // Strong Galerkin convection term: v div(E u) 419 v[4][i] = -wdetJ * strong_conv; 420 for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = 0; 421 } else { // Weak Galerkin convection term: -dv \cdot (E u) 422 for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = wdetJ * s.U.E_total * uX[j]; 423 v[4][i] = 0.; 424 } 425 426 { // Diffusion 427 CeedScalar Fe[3], Fe_dXdx[3] = {0.}; 428 429 for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total; 430 MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx); 431 for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] += wdetJ * Fe_dXdx[k]; 432 } 433 434 const CeedScalar TauS = Tau(context, s, dXdx, dim); 435 for (CeedInt j = 0; j < dim; j++) switch (context->stabilization) { 436 case STAB_NONE: 437 break; 438 case STAB_SU: 439 case STAB_SUPG: 440 grad_v[j][4][i] -= wdetJ * TauS * strong_conv * uX[j]; 441 break; 442 } 443 } 444 } 445 446 CEED_QFUNCTION(RHS_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 447 RHSFunction_AdvectionGeneric(ctx, Q, in, out, 3); 448 return 0; 449 } 450 451 CEED_QFUNCTION(RHS_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 452 RHSFunction_AdvectionGeneric(ctx, Q, in, out, 2); 453 return 0; 454 } 455 456 // ***************************************************************************** 457 // This QFunction implements consistent outflow and inflow BCs 458 // for advection 459 // 460 // Inflow and outflow faces are determined based on sign(dot(wind, normal)): 461 // sign(dot(wind, normal)) > 0 : outflow BCs 462 // sign(dot(wind, normal)) < 0 : inflow BCs 463 // 464 // Outflow BCs: 465 // The validity of the weak form of the governing equations is extended to the outflow and the current values of E are applied. 466 // 467 // Inflow BCs: 468 // A prescribed Total Energy (E_wind) is applied weakly. 469 // ***************************************************************************** 470 CEED_QFUNCTION(Advection_InOutFlowGeneric)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 471 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 472 const CeedScalar(*q_data_sur) = in[2]; 473 474 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 475 AdvectionContext context = (AdvectionContext)ctx; 476 const CeedScalar E_wind = context->E_wind; 477 const CeedScalar strong_form = context->strong_form; 478 const bool is_implicit = context->implicit; 479 480 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 481 const CeedScalar rho = q[0][i]; 482 const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 483 const CeedScalar E = q[4][i]; 484 485 CeedScalar wdetJb, normal[3]; 486 QdataBoundaryUnpack_ND(dim, Q, i, q_data_sur, &wdetJb, NULL, normal); 487 wdetJb *= is_implicit ? -1. : 1.; 488 489 const CeedScalar u_normal = DotN(normal, u, dim); 490 491 // No Change in density or momentum 492 for (CeedInt j = 0; j < 4; j++) { 493 v[j][i] = 0; 494 } 495 // Implementing in/outflow BCs 496 if (u_normal > 0) { // outflow 497 v[4][i] = -(1 - strong_form) * wdetJb * E * u_normal; 498 } else { // inflow 499 v[4][i] = -(1 - strong_form) * wdetJb * E_wind * u_normal; 500 } 501 } 502 return 0; 503 } 504 505 CEED_QFUNCTION(Advection_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 506 Advection_InOutFlowGeneric(ctx, Q, in, out, 3); 507 return 0; 508 } 509 510 CEED_QFUNCTION(Advection2d_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 511 Advection_InOutFlowGeneric(ctx, Q, in, out, 2); 512 return 0; 513 } 514