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 int 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 switch (context->wind_type) { 80 case ADVDIF_WIND_ROTATION: 81 q[0] = 1.; 82 q[1] = -(y - center[1]); 83 q[2] = (x - center[0]); 84 q[3] = 0; 85 break; 86 case ADVDIF_WIND_TRANSLATION: 87 q[0] = 1.; 88 q[1] = wind[0]; 89 q[2] = wind[1]; 90 q[3] = dim == 2 ? 0. : wind[2]; 91 break; 92 case ADVDIF_WIND_BOUNDARY_LAYER: 93 q[0] = 1.; 94 q[1] = y / ly; 95 q[2] = 0.; 96 q[3] = 0.; 97 break; 98 } 99 100 switch (context->initial_condition_type) { 101 case ADVDIF_IC_BUBBLE_SPHERE: 102 case ADVDIF_IC_BUBBLE_CYLINDER: { 103 CeedScalar r = sqrt(Square(x - x0[0]) + Square(y - x0[1]) + Square(z - x0[2])); 104 105 switch (context->bubble_continuity_type) { 106 // original continuous, smooth shape 107 case ADVDIF_BUBBLE_CONTINUITY_SMOOTH: 108 q[4] = r <= rc ? (1. - r / rc) : 0.; 109 break; 110 // discontinuous, sharp back half shape 111 case ADVDIF_BUBBLE_CONTINUITY_BACK_SHARP: 112 q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) : 0.; 113 break; 114 // attempt to define a finite thickness that will get resolved under grid refinement 115 case ADVDIF_BUBBLE_CONTINUITY_THICK: 116 q[4] = ((r <= rc) && (y < center[1])) ? (1. - r / rc) * fmin(1.0, (center[1] - y) / 1.25) : 0.; 117 break; 118 case ADVDIF_BUBBLE_CONTINUITY_COSINE: 119 q[4] = r <= rc ? .5 + .5 * cos(r * M_PI / rc) : 0; 120 break; 121 } 122 break; 123 } 124 125 case ADVDIF_IC_COSINE_HILL: { 126 CeedScalar r = sqrt(Square(x - center[0]) + Square(y - center[1])); 127 CeedScalar half_width = context->lx / 2; 128 q[4] = r > half_width ? 0. : cos(2 * M_PI * r / half_width + M_PI) + 1.; 129 } break; 130 131 case ADVDIF_IC_SKEW: { 132 CeedScalar skewed_barrier[3] = {wind[0], wind[1], 0}; 133 CeedScalar inflow_to_point[3] = {x - context->lx / 2, y, 0}; 134 CeedScalar cross_product[3] = {0}; 135 const CeedScalar boundary_threshold = 20 * CEED_EPSILON; 136 Cross3(skewed_barrier, inflow_to_point, cross_product); 137 138 q[4] = cross_product[2] > boundary_threshold ? 0 : 1; 139 if ((x < boundary_threshold && wind[0] < boundary_threshold) || // outflow at -x boundary 140 (y < boundary_threshold && wind[1] < boundary_threshold) || // outflow at -y boundary 141 (x > context->lx - boundary_threshold && wind[0] > boundary_threshold) || // outflow at +x boundary 142 (y > context->ly - boundary_threshold && wind[1] > boundary_threshold) // outflow at +y boundary 143 ) { 144 q[4] = 0; 145 } 146 } break; 147 148 case ADVDIF_IC_WAVE: { 149 CeedScalar theta = context->wave_frequency * DotN(X, wind, dim) + context->wave_phase; 150 switch (context->wave_type) { 151 case ADVDIF_WAVE_SINE: 152 q[4] = sin(theta); 153 break; 154 case ADVDIF_WAVE_SQUARE: 155 q[4] = sin(theta) > 100 * CEED_EPSILON ? 1 : -1; 156 break; 157 } 158 } break; 159 case ADVDIF_IC_BOUNDARY_LAYER: { 160 const CeedScalar boundary_threshold = 20 * CEED_EPSILON; 161 162 if ((x < boundary_threshold) || (y > ly - boundary_threshold)) { 163 q[4] = 1; // inflow and top boundary 164 } else if (y < boundary_threshold) { 165 q[4] = 0; // lower wall 166 } else { 167 q[4] = y / ly; // interior and outflow boundary 168 } 169 } break; 170 } 171 return 0; 172 } 173 174 // ***************************************************************************** 175 // This QFunction sets the initial conditions for 3D advection 176 // ***************************************************************************** 177 CEED_QFUNCTION(ICsAdvection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 178 const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 179 CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 180 181 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 182 const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 183 CeedScalar q[5] = {0.}; 184 185 Exact_AdvectionGeneric(3, 0., x, 5, q, ctx); 186 for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 187 } 188 return 0; 189 } 190 191 // ***************************************************************************** 192 // This QFunction sets the initial conditions for 2D advection 193 // ***************************************************************************** 194 CEED_QFUNCTION(ICsAdvection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 195 const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 196 CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 197 const SetupContextAdv context = (SetupContextAdv)ctx; 198 199 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 200 const CeedScalar x[] = {X[0][i], X[1][i]}; 201 CeedScalar q[5] = {0.}; 202 203 Exact_AdvectionGeneric(2, context->time, x, 5, q, ctx); 204 for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 205 } 206 return 0; 207 } 208 209 CEED_QFUNCTION_HELPER void StatePhysicalGradientFromReference_ND(CeedInt N, CeedInt Q, CeedInt i, NewtonianIdealGasContext gas, State s, 210 StateVariable state_var, const CeedScalar *grad_q, const CeedScalar *dXdx, 211 State *grad_s) { 212 switch (N) { 213 case 2: { 214 for (CeedInt k = 0; k < 2; k++) { 215 CeedScalar dqi[5]; 216 for (CeedInt j = 0; j < 5; j++) { 217 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]; 218 } 219 grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var); 220 } 221 CeedScalar U[5] = {0.}; 222 grad_s[2] = StateFromU(gas, U); 223 } break; 224 case 3: 225 // Cannot directly use StatePhysicalGradientFromReference helper functions due to SYCL online compiler incompatabilities 226 for (CeedInt k = 0; k < 3; k++) { 227 CeedScalar dqi[5]; 228 for (CeedInt j = 0; j < 5; j++) { 229 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] + 230 grad_q[(Q * 5) * 2 + Q * j + i] * dXdx[2 * N + k]; 231 } 232 grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var); 233 } 234 break; 235 } 236 } 237 238 // @brief Calculate the stabilization constant \tau 239 CEED_QFUNCTION_HELPER CeedScalar Tau(AdvectionContext context, const State s, const CeedScalar *dXdx, CeedInt dim) { 240 switch (context->stabilization_tau) { 241 case STAB_TAU_CTAU: { 242 CeedScalar uX[3] = {0.}; 243 244 MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 245 return context->CtauS / sqrt(DotN(uX, uX, dim)); 246 } break; 247 case STAB_TAU_ADVDIFF_SHAKIB: { 248 CeedScalar gijd_mat[9] = {0.}, gij_uj[3] = {0.}; 249 250 MatMatN(dXdx, dXdx, dim, CEED_TRANSPOSE, CEED_NOTRANSPOSE, gijd_mat); 251 // (1/2)^2 to account for reference element size; for length 1 square/cube element, gij should be identity matrix 252 ScaleN(gijd_mat, 0.25, Square(dim)); 253 MatVecNM(gijd_mat, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, gij_uj); 254 return 1 / sqrt(Square(2 * context->Ctau_t / context->dt) + DotN(s.Y.velocity, gij_uj, dim) * Square(context->Ctau_a) + 255 Square(context->diffusion_coeff) * DotN(gijd_mat, gijd_mat, dim * dim) * Square(context->Ctau_d)); 256 } break; 257 default: 258 return 0.; 259 } 260 } 261 262 // ***************************************************************************** 263 // This QFunction implements Advection for implicit time stepping method 264 // ***************************************************************************** 265 CEED_QFUNCTION_HELPER int IFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 266 AdvectionContext context = (AdvectionContext)ctx; 267 268 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 269 const CeedScalar(*grad_q) = in[1]; 270 const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 271 const CeedScalar(*q_data) = in[3]; 272 const CeedScalar(*divFdiff) = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? in[5] : NULL; 273 274 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 275 CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 276 277 NewtonianIdealGasContext gas; 278 struct NewtonianIdealGasContext_ gas_struct = {0}; 279 gas = &gas_struct; 280 281 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 282 const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]}; 283 const State s = StateFromU(gas, qi); 284 285 CeedScalar wdetJ, dXdx[9]; 286 QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 287 State grad_s[3]; 288 StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s); 289 290 const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total}; 291 292 for (CeedInt f = 0; f < 4; f++) { 293 for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum 294 v[f][i] = wdetJ * q_dot[f][i]; // K Mass/transient term 295 } 296 297 CeedScalar div_u = 0; 298 for (CeedInt j = 0; j < dim; j++) { 299 for (CeedInt k = 0; k < dim; k++) { 300 div_u += grad_s[k].Y.velocity[j]; 301 } 302 } 303 CeedScalar uX[3] = {0.}; 304 MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 305 CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim); 306 307 v[4][i] = wdetJ * q_dot[4][i]; // transient part (ALWAYS) 308 if (context->strong_form) { // Strong Galerkin convection term: v div(E u) 309 v[4][i] += wdetJ * strong_conv; 310 } else { // Weak Galerkin convection term: -dv \cdot (E u) 311 for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = -wdetJ * s.U.E_total * uX[j]; 312 } 313 314 { // Diffusion 315 CeedScalar Fe[3], Fe_dXdx[3] = {0.}; 316 317 for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total; 318 MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx); 319 for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] -= wdetJ * Fe_dXdx[k]; 320 } 321 322 const CeedScalar TauS = Tau(context, s, dXdx, dim); 323 for (CeedInt j = 0; j < dim; j++) { 324 switch (context->stabilization) { 325 case STAB_NONE: 326 break; 327 case STAB_SU: 328 grad_v[j][4][i] += wdetJ * TauS * uX[j] * strong_conv; 329 break; 330 case STAB_SUPG: { 331 CeedScalar divFdiff_i = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? divFdiff[i] : 0.; 332 grad_v[j][4][i] += wdetJ * TauS * uX[j] * (q_dot[4][i] + strong_conv + divFdiff_i); 333 } break; 334 } 335 } 336 } 337 return 0; 338 } 339 340 CEED_QFUNCTION(IFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 341 return IFunction_AdvectionGeneric(ctx, Q, in, out, 3); 342 } 343 344 CEED_QFUNCTION(IFunction_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 345 return IFunction_AdvectionGeneric(ctx, Q, in, out, 2); 346 } 347 348 CEED_QFUNCTION_HELPER int MassFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 349 const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 350 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 351 const CeedScalar(*q_data) = in[2]; 352 353 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 354 CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 355 356 AdvectionContext context = (AdvectionContext)ctx; 357 struct NewtonianIdealGasContext_ gas_struct = {0}; 358 NewtonianIdealGasContext gas = &gas_struct; 359 360 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 361 const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]}; 362 const State s = StateFromU(gas, qi); 363 CeedScalar wdetJ, dXdx[9]; 364 QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 365 366 for (CeedInt f = 0; f < 4; f++) { 367 for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum 368 v[f][i] = wdetJ * q_dot[f][i]; // K Mass/transient term 369 } 370 371 // Unstabilized mass term 372 v[4][i] = wdetJ * q_dot[4][i]; 373 374 // Stabilized mass term 375 CeedScalar uX[3] = {0.}; 376 MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 377 const CeedScalar TauS = Tau(context, s, dXdx, dim); 378 for (CeedInt j = 0; j < dim; j++) { 379 switch (context->stabilization) { 380 case STAB_NONE: 381 case STAB_SU: 382 grad_v[j][4][i] = 0; 383 break; // These should be run with the unstabilized mass matrix anyways 384 case STAB_SUPG: 385 grad_v[j][4][i] = wdetJ * TauS * q_dot[4][i] * uX[j]; 386 break; 387 } 388 } 389 } 390 return 0; 391 } 392 393 CEED_QFUNCTION(MassFunction_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 394 return MassFunction_AdvectionGeneric(ctx, Q, in, out, 3); 395 } 396 397 CEED_QFUNCTION(MassFunction_Advection2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 398 return MassFunction_AdvectionGeneric(ctx, Q, in, out, 2); 399 } 400 401 // ***************************************************************************** 402 // This QFunction implements Advection for explicit time stepping method 403 // ***************************************************************************** 404 CEED_QFUNCTION_HELPER int RHSFunction_AdvectionGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 405 AdvectionContext context = (AdvectionContext)ctx; 406 407 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 408 const CeedScalar(*grad_q) = in[1]; 409 const CeedScalar(*q_data) = in[2]; 410 const CeedScalar(*divFdiff) = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? in[4] : NULL; 411 412 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 413 CeedScalar(*grad_v)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 414 415 struct NewtonianIdealGasContext_ gas_struct = {0}; 416 NewtonianIdealGasContext gas = &gas_struct; 417 418 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 419 const CeedScalar qi[5] = {q[0][i], q[1][i], q[2][i], q[3][i], q[4][i]}; 420 const State s = StateFromU(gas, qi); 421 422 CeedScalar wdetJ, dXdx[9]; 423 QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 424 State grad_s[3]; 425 StatePhysicalGradientFromReference_ND(dim, Q, i, gas, s, STATEVAR_CONSERVATIVE, grad_q, dXdx, grad_s); 426 427 const CeedScalar Grad_E[3] = {grad_s[0].U.E_total, grad_s[1].U.E_total, grad_s[2].U.E_total}; 428 429 for (CeedInt f = 0; f < 4; f++) { 430 for (CeedInt j = 0; j < dim; j++) grad_v[j][f][i] = 0; // No Change in density or momentum 431 v[f][i] = 0.; 432 } 433 434 CeedScalar div_u = 0; 435 for (CeedInt j = 0; j < dim; j++) { 436 for (CeedInt k = 0; k < dim; k++) { 437 div_u += grad_s[k].Y.velocity[j]; 438 } 439 } 440 CeedScalar strong_conv = s.U.E_total * div_u + DotN(s.Y.velocity, Grad_E, dim); 441 442 CeedScalar uX[3] = {0.}; 443 MatVecNM(dXdx, s.Y.velocity, dim, dim, CEED_NOTRANSPOSE, uX); 444 445 if (context->strong_form) { // Strong Galerkin convection term: v div(E u) 446 v[4][i] = -wdetJ * strong_conv; 447 for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = 0; 448 } else { // Weak Galerkin convection term: -dv \cdot (E u) 449 for (CeedInt j = 0; j < dim; j++) grad_v[j][4][i] = wdetJ * s.U.E_total * uX[j]; 450 v[4][i] = 0.; 451 } 452 453 { // Diffusion 454 CeedScalar Fe[3], Fe_dXdx[3] = {0.}; 455 456 for (CeedInt i = 0; i < dim; i++) Fe[i] = -context->diffusion_coeff * grad_s[i].U.E_total; 457 MatVecNM(dXdx, Fe, dim, dim, CEED_NOTRANSPOSE, Fe_dXdx); 458 for (CeedInt k = 0; k < dim; k++) grad_v[k][4][i] += wdetJ * Fe_dXdx[k]; 459 } 460 461 const CeedScalar TauS = Tau(context, s, dXdx, dim); 462 for (CeedInt j = 0; j < dim; j++) { 463 switch (context->stabilization) { 464 case STAB_NONE: 465 break; 466 case STAB_SU: 467 case STAB_SUPG: { 468 CeedScalar divFdiff_i = context->divFdiff_method != DIV_DIFF_FLUX_PROJ_NONE ? divFdiff[i] : 0.; 469 grad_v[j][4][i] -= wdetJ * TauS * (strong_conv + divFdiff_i) * uX[j]; 470 } break; 471 } 472 } 473 } 474 return 0; 475 } 476 477 CEED_QFUNCTION(RHS_Advection)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 478 return RHSFunction_AdvectionGeneric(ctx, Q, in, out, 3); 479 } 480 481 CEED_QFUNCTION(RHS_Advection2d)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 482 return RHSFunction_AdvectionGeneric(ctx, Q, in, out, 2); 483 } 484 485 // ***************************************************************************** 486 // This QFunction implements consistent outflow and inflow BCs 487 // for advection 488 // 489 // Inflow and outflow faces are determined based on sign(dot(wind, normal)): 490 // sign(dot(wind, normal)) > 0 : outflow BCs 491 // sign(dot(wind, normal)) < 0 : inflow BCs 492 // 493 // Outflow BCs: 494 // The validity of the weak form of the governing equations is extended to the outflow and the current values of E are applied. 495 // 496 // Inflow BCs: 497 // A prescribed Total Energy (E_wind) is applied weakly. 498 // ***************************************************************************** 499 CEED_QFUNCTION_HELPER int Advection_InOutFlowGeneric(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, CeedInt dim) { 500 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 501 const CeedScalar(*q_data_sur) = in[2]; 502 503 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 504 AdvectionContext context = (AdvectionContext)ctx; 505 const CeedScalar E_wind = context->E_wind; 506 const CeedScalar strong_form = context->strong_form; 507 const bool is_implicit = context->implicit; 508 509 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 510 const CeedScalar rho = q[0][i]; 511 const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 512 const CeedScalar E = q[4][i]; 513 514 CeedScalar wdetJb, normal[3]; 515 QdataBoundaryUnpack_ND(dim, Q, i, q_data_sur, &wdetJb, NULL, normal); 516 wdetJb *= is_implicit ? -1. : 1.; 517 518 const CeedScalar u_normal = DotN(normal, u, dim); 519 520 // No Change in density or momentum 521 for (CeedInt j = 0; j < 4; j++) { 522 v[j][i] = 0; 523 } 524 // Implementing in/outflow BCs 525 if (u_normal > 0) { // outflow 526 v[4][i] = -(1 - strong_form) * wdetJb * E * u_normal; 527 } else { // inflow 528 v[4][i] = -(1 - strong_form) * wdetJb * E_wind * u_normal; 529 } 530 } 531 return 0; 532 } 533 534 CEED_QFUNCTION(Advection_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 535 return Advection_InOutFlowGeneric(ctx, Q, in, out, 3); 536 } 537 538 CEED_QFUNCTION(Advection2d_InOutFlow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 539 return Advection_InOutFlowGeneric(ctx, Q, in, out, 2); 540 } 541 542 // @brief Volume integral for RHS of divergence of diffusive flux direct projection 543 CEED_QFUNCTION_HELPER int DivDiffusiveFluxVolumeRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, 544 const CeedInt dim) { 545 const CeedScalar(*Grad_q) = in[0]; 546 const CeedScalar(*q_data) = in[1]; 547 CeedScalar(*Grad_v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 548 549 AdvectionContext context = (AdvectionContext)ctx; 550 551 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 552 CeedScalar wdetJ, dXdx[9], F_diff[3] = {0.}; 553 554 QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 555 { // Get physical diffusive flux 556 CeedScalar Grad_qn[15], grad_E_ref[3]; 557 558 GradUnpackN(Q, i, 5, dim, Grad_q, Grad_qn); 559 CopyN(&Grad_qn[4 * dim], grad_E_ref, dim); 560 MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff); 561 ScaleN(F_diff, -context->diffusion_coeff, dim); 562 } 563 564 CeedScalar F_diff_dXdx[3] = {0.}; 565 MatVecNM(dXdx, F_diff, dim, dim, CEED_NOTRANSPOSE, F_diff_dXdx); 566 for (CeedInt k = 0; k < dim; k++) Grad_v[k][i] = -wdetJ * F_diff_dXdx[k]; 567 } 568 return 0; 569 } 570 571 CEED_QFUNCTION(DivDiffusiveFluxVolumeRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 572 return DivDiffusiveFluxVolumeRHS_AdvDif_Generic(ctx, Q, in, out, 2); 573 } 574 575 CEED_QFUNCTION(DivDiffusiveFluxVolumeRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 576 return DivDiffusiveFluxVolumeRHS_AdvDif_Generic(ctx, Q, in, out, 3); 577 } 578 579 // @brief Boundary integral for RHS of divergence of diffusive flux direct projection 580 CEED_QFUNCTION_HELPER int DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, 581 const CeedInt dim) { 582 const CeedScalar(*Grad_q) = in[0]; 583 const CeedScalar(*q_data) = in[1]; 584 CeedScalar(*v) = out[0]; 585 586 AdvectionContext context = (AdvectionContext)ctx; 587 588 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 589 CeedScalar wdetJ, normal[3], dXdx[9], F_diff[3] = {0.}; 590 591 QdataBoundaryGradientUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx, normal); 592 { // Get physical diffusive flux 593 CeedScalar Grad_qn[15], grad_E_ref[3]; 594 595 GradUnpackN(Q, i, 5, dim, Grad_q, Grad_qn); 596 CopyN(&Grad_qn[4 * dim], grad_E_ref, dim); 597 MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff); 598 ScaleN(F_diff, -context->diffusion_coeff, dim); 599 } 600 601 v[i] = wdetJ * DotN(F_diff, normal, dim); 602 } 603 return 0; 604 } 605 606 CEED_QFUNCTION(DivDiffusiveFluxBoundaryRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 607 return DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(ctx, Q, in, out, 2); 608 } 609 610 CEED_QFUNCTION(DivDiffusiveFluxBoundaryRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 611 return DivDiffusiveFluxBoundaryRHS_AdvDif_Generic(ctx, Q, in, out, 3); 612 } 613 614 // @brief Volume integral for RHS of diffusive flux indirect projection 615 CEED_QFUNCTION_HELPER int DiffusiveFluxRHS_AdvDif_Generic(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out, 616 const CeedInt dim) { 617 const CeedScalar(*Grad_q) = in[0]; 618 const CeedScalar(*q_data) = in[1]; 619 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 620 621 AdvectionContext context = (AdvectionContext)ctx; 622 623 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 624 CeedScalar wdetJ, dXdx[9], F_diff[3] = {0.}; 625 626 QdataUnpack_ND(dim, Q, i, q_data, &wdetJ, dXdx); 627 { // Get physical diffusive flux 628 CeedScalar Grad_qn[15], grad_E_ref[3]; 629 630 GradUnpackN(Q, i, 5, dim, Grad_q, Grad_qn); 631 CopyN(&Grad_qn[4 * dim], grad_E_ref, dim); 632 MatVecNM(dXdx, grad_E_ref, dim, dim, CEED_NOTRANSPOSE, F_diff); 633 ScaleN(F_diff, -context->diffusion_coeff, dim); 634 } 635 for (CeedInt k = 0; k < dim; k++) v[k][i] = wdetJ * F_diff[k]; 636 } 637 return 0; 638 } 639 640 CEED_QFUNCTION(DiffusiveFluxRHS_AdvDif_2D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 641 return DiffusiveFluxRHS_AdvDif_Generic(ctx, Q, in, out, 2); 642 } 643 644 CEED_QFUNCTION(DiffusiveFluxRHS_AdvDif_3D)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 645 return DiffusiveFluxRHS_AdvDif_Generic(ctx, Q, in, out, 3); 646 } 647