1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, 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 /// Euler traveling vortex initial condition and operator for Navier-Stokes 10a515125bSLeila Ghaffari /// example using PETSc 11a515125bSLeila Ghaffari 12a515125bSLeila Ghaffari // Model from: 13*04e40bb6SJeremy L Thompson // On the Order of Accuracy and Numerical Performance of Two Classes of Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011). 14a515125bSLeila Ghaffari 15a515125bSLeila Ghaffari #ifndef eulervortex_h 16a515125bSLeila Ghaffari #define eulervortex_h 17a515125bSLeila Ghaffari 183a8779fbSJames Wright #include <ceed.h> 19d0cce58aSJeremy L Thompson #include <math.h> 202b916ea7SJeremy L Thompson 21704b8bbeSJames Wright #include "utils.h" 22a515125bSLeila Ghaffari 23a515125bSLeila Ghaffari typedef struct EulerContext_ *EulerContext; 24a515125bSLeila Ghaffari struct EulerContext_ { 25a515125bSLeila Ghaffari CeedScalar center[3]; 26a515125bSLeila Ghaffari CeedScalar curr_time; 27a515125bSLeila Ghaffari CeedScalar vortex_strength; 28d8a22b9eSJed Brown CeedScalar c_tau; 29a515125bSLeila Ghaffari CeedScalar mean_velocity[3]; 30a515125bSLeila Ghaffari bool implicit; 31139613f2SLeila Ghaffari int euler_test; 32139613f2SLeila Ghaffari int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG 33a515125bSLeila Ghaffari }; 34a515125bSLeila Ghaffari 35a515125bSLeila Ghaffari // ***************************************************************************** 36a515125bSLeila Ghaffari // This function sets the initial conditions 37a515125bSLeila Ghaffari // 38a515125bSLeila Ghaffari // Temperature: 39a515125bSLeila Ghaffari // T = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2) 40a515125bSLeila Ghaffari // Density: 41a515125bSLeila Ghaffari // rho = (T/S_vortex)^(1 / (gamma - 1)) 42a515125bSLeila Ghaffari // Pressure: 43a515125bSLeila Ghaffari // P = rho * T 44a515125bSLeila Ghaffari // Velocity: 45a515125bSLeila Ghaffari // ui = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi) 46a515125bSLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 ) 47a515125bSLeila Ghaffari // Velocity/Momentum Density: 48a515125bSLeila Ghaffari // Ui = rho ui 49a515125bSLeila Ghaffari // Total Energy: 50a515125bSLeila Ghaffari // E = P / (gamma - 1) + rho (u u)/2 51a515125bSLeila Ghaffari // 52a515125bSLeila Ghaffari // Constants: 53a515125bSLeila Ghaffari // cv , Specific heat, constant volume 54a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 55a515125bSLeila Ghaffari // vortex_strength , Strength of vortex 56a515125bSLeila Ghaffari // center , Location of bubble center 57a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio 58a515125bSLeila Ghaffari // 59a515125bSLeila Ghaffari // ***************************************************************************** 60a515125bSLeila Ghaffari 61a515125bSLeila Ghaffari // ***************************************************************************** 62*04e40bb6SJeremy L Thompson // This helper function provides support for the exact, time-dependent solution (currently not implemented) and IC formulation for Euler traveling 63*04e40bb6SJeremy L Thompson // vortex 64a515125bSLeila Ghaffari // ***************************************************************************** 652b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) { 66a515125bSLeila Ghaffari // Context 67a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx; 68a515125bSLeila Ghaffari const CeedScalar vortex_strength = context->vortex_strength; 69a515125bSLeila Ghaffari const CeedScalar *center = context->center; // Center of the domain 70a515125bSLeila Ghaffari const CeedScalar *mean_velocity = context->mean_velocity; 71a515125bSLeila Ghaffari 72a515125bSLeila Ghaffari // Setup 73a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 74a515125bSLeila Ghaffari const CeedScalar cv = 2.5; 75a515125bSLeila Ghaffari const CeedScalar R = 1.; 76a515125bSLeila Ghaffari const CeedScalar x = X[0], y = X[1]; // Coordinates 77a515125bSLeila Ghaffari // Vortex center 78a515125bSLeila Ghaffari const CeedScalar xc = center[0] + mean_velocity[0] * time; 79a515125bSLeila Ghaffari const CeedScalar yc = center[1] + mean_velocity[1] * time; 80a515125bSLeila Ghaffari 81a515125bSLeila Ghaffari const CeedScalar x0 = x - xc; 82a515125bSLeila Ghaffari const CeedScalar y0 = y - yc; 83a515125bSLeila Ghaffari const CeedScalar r = sqrt(x0 * x0 + y0 * y0); 84a515125bSLeila Ghaffari const CeedScalar C = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI); 852b916ea7SJeremy L Thompson const CeedScalar delta_T = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI); 86a515125bSLeila Ghaffari const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma 872b916ea7SJeremy L Thompson const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI); 88a515125bSLeila Ghaffari CeedScalar rho, P, T, E, u[3] = {0.}; 89a515125bSLeila Ghaffari 90a515125bSLeila Ghaffari // Initial Conditions 91a515125bSLeila Ghaffari switch (context->euler_test) { 92a515125bSLeila Ghaffari case 0: // Traveling vortex 93a515125bSLeila Ghaffari T = 1 + delta_T; 94a515125bSLeila Ghaffari // P = rho * T 95a515125bSLeila Ghaffari // P = S * rho^gamma 96a515125bSLeila Ghaffari // Solve for rho, then substitute for P 97139613f2SLeila Ghaffari rho = pow(T / S_vortex, 1 / (gamma - 1.)); 98a515125bSLeila Ghaffari P = rho * T; 99a515125bSLeila Ghaffari u[0] = mean_velocity[0] - C * y0; 100a515125bSLeila Ghaffari u[1] = mean_velocity[1] + C * x0; 101a515125bSLeila Ghaffari 102a515125bSLeila Ghaffari // Assign exact solution 103a515125bSLeila Ghaffari q[0] = rho; 104a515125bSLeila Ghaffari q[1] = rho * u[0]; 105a515125bSLeila Ghaffari q[2] = rho * u[1]; 106a515125bSLeila Ghaffari q[3] = rho * u[2]; 107a515125bSLeila Ghaffari q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.; 108a515125bSLeila Ghaffari break; 109a515125bSLeila Ghaffari case 1: // Constant zero velocity, density constant, total energy constant 110a515125bSLeila Ghaffari rho = 1.; 111a515125bSLeila Ghaffari E = 2.; 112a515125bSLeila Ghaffari 113a515125bSLeila Ghaffari // Assign exact solution 114a515125bSLeila Ghaffari q[0] = rho; 115a515125bSLeila Ghaffari q[1] = rho * u[0]; 116a515125bSLeila Ghaffari q[2] = rho * u[1]; 117a515125bSLeila Ghaffari q[3] = rho * u[2]; 118a515125bSLeila Ghaffari q[4] = E; 119a515125bSLeila Ghaffari break; 120a515125bSLeila Ghaffari case 2: // Constant nonzero velocity, density constant, total energy constant 121a515125bSLeila Ghaffari rho = 1.; 122a515125bSLeila Ghaffari E = 2.; 123a515125bSLeila Ghaffari u[0] = mean_velocity[0]; 124a515125bSLeila Ghaffari u[1] = mean_velocity[1]; 125a515125bSLeila Ghaffari 126a515125bSLeila Ghaffari // Assign exact solution 127a515125bSLeila Ghaffari q[0] = rho; 128a515125bSLeila Ghaffari q[1] = rho * u[0]; 129a515125bSLeila Ghaffari q[2] = rho * u[1]; 130a515125bSLeila Ghaffari q[3] = rho * u[2]; 131a515125bSLeila Ghaffari q[4] = E; 132a515125bSLeila Ghaffari break; 133*04e40bb6SJeremy L Thompson case 3: // Velocity zero, pressure constant (so density and internal energy will be non-constant), but the velocity should stay zero and the 134*04e40bb6SJeremy L Thompson // bubble won't diffuse 135a515125bSLeila Ghaffari // (for Euler, where there is no thermal conductivity) 136a515125bSLeila Ghaffari P = 1.; 137a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r * r); 138a515125bSLeila Ghaffari rho = P / (R * T); 139a515125bSLeila Ghaffari 140a515125bSLeila Ghaffari // Assign exact solution 141a515125bSLeila Ghaffari q[0] = rho; 142a515125bSLeila Ghaffari q[1] = rho * u[0]; 143a515125bSLeila Ghaffari q[2] = rho * u[1]; 144a515125bSLeila Ghaffari q[3] = rho * u[2]; 145a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 146a515125bSLeila Ghaffari break; 147*04e40bb6SJeremy L Thompson case 4: // Constant nonzero velocity, pressure constant (so density and internal energy will be non-constant), 148*04e40bb6SJeremy L Thompson // It should be transported across the domain, but velocity stays constant 149a515125bSLeila Ghaffari P = 1.; 150a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r * r); 151a515125bSLeila Ghaffari rho = P / (R * T); 152a515125bSLeila Ghaffari u[0] = mean_velocity[0]; 153a515125bSLeila Ghaffari u[1] = mean_velocity[1]; 154a515125bSLeila Ghaffari 155a515125bSLeila Ghaffari // Assign exact solution 156a515125bSLeila Ghaffari q[0] = rho; 157a515125bSLeila Ghaffari q[1] = rho * u[0]; 158a515125bSLeila Ghaffari q[2] = rho * u[1]; 159a515125bSLeila Ghaffari q[3] = rho * u[2]; 160a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 161a515125bSLeila Ghaffari break; 1620df2634dSLeila Ghaffari case 5: // non-smooth thermal bubble - cylinder 1630df2634dSLeila Ghaffari P = 1.; 1640df2634dSLeila Ghaffari T = 1. - (r < 1. ? S_bubble : 0.); 1650df2634dSLeila Ghaffari rho = P / (R * T); 1660df2634dSLeila Ghaffari u[0] = mean_velocity[0]; 1670df2634dSLeila Ghaffari u[1] = mean_velocity[1]; 1680df2634dSLeila Ghaffari 1690df2634dSLeila Ghaffari // Assign exact solution 1700df2634dSLeila Ghaffari q[0] = rho; 1710df2634dSLeila Ghaffari q[1] = rho * u[0]; 1720df2634dSLeila Ghaffari q[2] = rho * u[1]; 1730df2634dSLeila Ghaffari q[3] = rho * u[2]; 1740df2634dSLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 1750df2634dSLeila Ghaffari break; 176a515125bSLeila Ghaffari } 177a515125bSLeila Ghaffari // Return 178a515125bSLeila Ghaffari return 0; 179a515125bSLeila Ghaffari } 180a515125bSLeila Ghaffari 181a515125bSLeila Ghaffari // ***************************************************************************** 182139613f2SLeila Ghaffari // Helper function for computing flux Jacobian 183139613f2SLeila Ghaffari // ***************************************************************************** 1842b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 185139613f2SLeila Ghaffari const CeedScalar gamma) { 186139613f2SLeila Ghaffari CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2]; // Velocity square 187139613f2SLeila Ghaffari for (CeedInt i = 0; i < 3; i++) { // Jacobian matrices for 3 directions 188139613f2SLeila Ghaffari for (CeedInt j = 0; j < 3; j++) { // Rows of each Jacobian matrix 189139613f2SLeila Ghaffari dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j]; 190139613f2SLeila Ghaffari for (CeedInt k = 0; k < 3; k++) { // Columns of each Jacobian matrix 191139613f2SLeila Ghaffari dF[i][0][k + 1] = ((i == k) ? 1. : 0.); 1922b916ea7SJeremy L Thompson dF[i][j + 1][k + 1] = ((j == k) ? u[i] : 0.) + ((i == k) ? u[j] : 0.) - ((i == j) ? u[k] : 0.) * (gamma - 1.); 1932b916ea7SJeremy L Thompson dF[i][4][k + 1] = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k]; 194139613f2SLeila Ghaffari } 195139613f2SLeila Ghaffari dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.); 196139613f2SLeila Ghaffari } 197139613f2SLeila Ghaffari dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho); 198139613f2SLeila Ghaffari dF[i][4][4] = u[i] * gamma; 199139613f2SLeila Ghaffari } 200139613f2SLeila Ghaffari } 201139613f2SLeila Ghaffari 202139613f2SLeila Ghaffari // ***************************************************************************** 203d8a22b9eSJed Brown // Helper function for computing Tau elements (stabilization constant) 204d8a22b9eSJed Brown // Model from: 205d8a22b9eSJed Brown // Stabilized Methods for Compressible Flows, Hughes et al 2010 206d8a22b9eSJed Brown // 207d8a22b9eSJed Brown // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 208d8a22b9eSJed Brown // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 209d8a22b9eSJed Brown // 210d8a22b9eSJed Brown // Where 211d8a22b9eSJed Brown // c_tau = stabilization constant (0.5 is reported as "optimal") 212d8a22b9eSJed Brown // h[i] = 2 length(dxdX[i]) 213d8a22b9eSJed Brown // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 214d8a22b9eSJed Brown // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 215*04e40bb6SJeremy L Thompson // rho(A[i]) = spectral radius of the convective flux Jacobian i, wave speed in direction i 216d8a22b9eSJed Brown // ***************************************************************************** 2172b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], const CeedScalar dXdx[3][3], const CeedScalar u[3], const CeedScalar sound_speed, 2182b916ea7SJeremy L Thompson const CeedScalar c_tau) { 219493642f1SJames Wright for (CeedInt i = 0; i < 3; i++) { 220d8a22b9eSJed Brown // length of element in direction i 2212b916ea7SJeremy L Thompson CeedScalar h = 2 / sqrt(dXdx[0][i] * dXdx[0][i] + dXdx[1][i] * dXdx[1][i] + dXdx[2][i] * dXdx[2][i]); 222d8a22b9eSJed Brown // fastest wave in direction i 223d8a22b9eSJed Brown CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 224d8a22b9eSJed Brown Tau_x[i] = c_tau * h / fastest_wave; 225d8a22b9eSJed Brown } 226d8a22b9eSJed Brown } 227d8a22b9eSJed Brown 228d8a22b9eSJed Brown // ***************************************************************************** 229a515125bSLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex 230a515125bSLeila Ghaffari // ***************************************************************************** 2312b916ea7SJeremy L Thompson CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 232a515125bSLeila Ghaffari // Inputs 233a515125bSLeila Ghaffari const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 234a515125bSLeila Ghaffari 235a515125bSLeila Ghaffari // Outputs 236a515125bSLeila Ghaffari CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 237a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx; 238a515125bSLeila Ghaffari 239a515125bSLeila Ghaffari // Quadrature Point Loop 2403d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 241a515125bSLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 242139613f2SLeila Ghaffari CeedScalar q[5] = {0.}; 243a515125bSLeila Ghaffari 244a515125bSLeila Ghaffari Exact_Euler(3, context->curr_time, x, 5, q, ctx); 245a515125bSLeila Ghaffari 2462b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 247a515125bSLeila Ghaffari } // End of Quadrature Point Loop 248a515125bSLeila Ghaffari 249a515125bSLeila Ghaffari // Return 250a515125bSLeila Ghaffari return 0; 251a515125bSLeila Ghaffari } 252a515125bSLeila Ghaffari 253a515125bSLeila Ghaffari // ***************************************************************************** 254*04e40bb6SJeremy L Thompson // This QFunction implements the following formulation of Euler equations with explicit time stepping method 255a515125bSLeila Ghaffari // 256*04e40bb6SJeremy L Thompson // This is 3D Euler for compressible gas dynamics in conservation form with state variables of density, momentum density, and total energy density. 257a515125bSLeila Ghaffari // 258a515125bSLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E ) 259a515125bSLeila Ghaffari // rho - Mass Density 260a515125bSLeila Ghaffari // Ui - Momentum Density, Ui = rho ui 261a515125bSLeila Ghaffari // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2 262a515125bSLeila Ghaffari // 263a515125bSLeila Ghaffari // Euler Equations: 264a515125bSLeila Ghaffari // drho/dt + div( U ) = 0 265a515125bSLeila Ghaffari // dU/dt + div( rho (u x u) + P I3 ) = 0 266a515125bSLeila Ghaffari // dE/dt + div( (E + P) u ) = 0 267a515125bSLeila Ghaffari // 268a515125bSLeila Ghaffari // Equation of State: 269a515125bSLeila Ghaffari // P = (gamma - 1) (E - rho (u u) / 2) 270a515125bSLeila Ghaffari // 271a515125bSLeila Ghaffari // Constants: 272a515125bSLeila Ghaffari // cv , Specific heat, constant volume 273a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 274a515125bSLeila Ghaffari // g , Gravity 275a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio 276a515125bSLeila Ghaffari // ***************************************************************************** 2772b916ea7SJeremy L Thompson CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 278a515125bSLeila Ghaffari // Inputs 2793d65b166SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 2803d65b166SJames Wright const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1]; 2813d65b166SJames Wright const CeedScalar(*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 2823d65b166SJames Wright 283a515125bSLeila Ghaffari // Outputs 2843d65b166SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 2853d65b166SJames Wright CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 286a515125bSLeila Ghaffari 287139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx; 288d8a22b9eSJed Brown const CeedScalar c_tau = context->c_tau; 289a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 290a515125bSLeila Ghaffari 291a515125bSLeila Ghaffari // Quadrature Point Loop 2923d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 293a515125bSLeila Ghaffari // Setup 294a515125bSLeila Ghaffari // -- Interp in 295a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 2962b916ea7SJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 297a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 2982b916ea7SJeremy L Thompson const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]}; 2992b916ea7SJeremy L Thompson const CeedScalar dU[3][3] = { 3002b916ea7SJeremy L Thompson {dq[0][1][i], dq[1][1][i], dq[2][1][i]}, 3012b916ea7SJeremy L Thompson {dq[0][2][i], dq[1][2][i], dq[2][2][i]}, 3022b916ea7SJeremy L Thompson {dq[0][3][i], dq[1][3][i], dq[2][3][i]} 303139613f2SLeila Ghaffari }; 3042b916ea7SJeremy L Thompson const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]}; 305a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 306a515125bSLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 307a515125bSLeila Ghaffari // -- Interp-to-Grad q_data 308a515125bSLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 3092b916ea7SJeremy L Thompson const CeedScalar dXdx[3][3] = { 3102b916ea7SJeremy L Thompson {q_data[1][i], q_data[2][i], q_data[3][i]}, 3112b916ea7SJeremy L Thompson {q_data[4][i], q_data[5][i], q_data[6][i]}, 3122b916ea7SJeremy L Thompson {q_data[7][i], q_data[8][i], q_data[9][i]} 313a515125bSLeila Ghaffari }; 314139613f2SLeila Ghaffari // dU/dx 315139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 316139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 317139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 318139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 319493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) { 320493642f1SJames Wright for (CeedInt k = 0; k < 3; k++) { 321139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 322139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 323493642f1SJames Wright for (CeedInt l = 0; l < 3; l++) { 324139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 325139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j 326139613f2SLeila Ghaffari } 327139613f2SLeila Ghaffari } 328139613f2SLeila Ghaffari } 329139613f2SLeila Ghaffari // Pressure 3302b916ea7SJeremy L Thompson const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic, 331139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 332a515125bSLeila Ghaffari 333a515125bSLeila Ghaffari // The Physics 334a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 335493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) { 336139613f2SLeila Ghaffari v[j][i] = 0.; 3372b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.; 338a515125bSLeila Ghaffari } 339a515125bSLeila Ghaffari 340a515125bSLeila Ghaffari // -- Density 341a515125bSLeila Ghaffari // ---- u rho 3422b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) dv[j][0][i] += wdetJ * (rho * u[0] * dXdx[j][0] + rho * u[1] * dXdx[j][1] + rho * u[2] * dXdx[j][2]); 343a515125bSLeila Ghaffari // -- Momentum 344a515125bSLeila Ghaffari // ---- rho (u x u) + P I3 3452b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3462b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { 3472b916ea7SJeremy L Thompson dv[k][j + 1][i] += wdetJ * ((rho * u[j] * u[0] + (j == 0 ? P : 0.)) * dXdx[k][0] + (rho * u[j] * u[1] + (j == 1 ? P : 0.)) * dXdx[k][1] + 348139613f2SLeila Ghaffari (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]); 3492b916ea7SJeremy L Thompson } 3502b916ea7SJeremy L Thompson } 351a515125bSLeila Ghaffari // -- Total Energy Density 352a515125bSLeila Ghaffari // ---- (E + P) u 3532b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) dv[j][4][i] += wdetJ * (E + P) * (u[0] * dXdx[j][0] + u[1] * dXdx[j][1] + u[2] * dXdx[j][2]); 354139613f2SLeila Ghaffari 355139613f2SLeila Ghaffari // --Stabilization terms 356139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 357139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 358d8a22b9eSJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 359139613f2SLeila Ghaffari 360139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 361139613f2SLeila Ghaffari CeedScalar dqdx[5][3]; 362493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) { 363139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j]; 364139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j]; 3652b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j]; 366139613f2SLeila Ghaffari } 367139613f2SLeila Ghaffari 368139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 369139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 3702b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3712b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 3722b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 3732b916ea7SJeremy L Thompson } 3742b916ea7SJeremy L Thompson } 375139613f2SLeila Ghaffari 376d8a22b9eSJed Brown // Stabilization 377d8a22b9eSJed Brown // -- Tau elements 378d8a22b9eSJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 379d8a22b9eSJed Brown CeedScalar Tau_x[3] = {0.}; 380d8a22b9eSJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 381139613f2SLeila Ghaffari 382d8a22b9eSJed Brown // -- Stabilization method: none or SU 383bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}}; 384139613f2SLeila Ghaffari switch (context->stabilization) { 385139613f2SLeila Ghaffari case 0: // Galerkin 386139613f2SLeila Ghaffari break; 387139613f2SLeila Ghaffari case 1: // SU 3882b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3892b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 3902b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 3912b916ea7SJeremy L Thompson } 3922b916ea7SJeremy L Thompson } 393139613f2SLeila Ghaffari 3942b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 3952b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] -= wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]); 3962b916ea7SJeremy L Thompson } 397139613f2SLeila Ghaffari break; 398139613f2SLeila Ghaffari case 2: // SUPG is not implemented for explicit scheme 399139613f2SLeila Ghaffari break; 400139613f2SLeila Ghaffari } 401139613f2SLeila Ghaffari 402a515125bSLeila Ghaffari } // End Quadrature Point Loop 403a515125bSLeila Ghaffari 404a515125bSLeila Ghaffari // Return 405a515125bSLeila Ghaffari return 0; 406a515125bSLeila Ghaffari } 407a515125bSLeila Ghaffari 408a515125bSLeila Ghaffari // ***************************************************************************** 409*04e40bb6SJeremy L Thompson // This QFunction implements the Euler equations with (mentioned above) with implicit time stepping method 410a515125bSLeila Ghaffari // ***************************************************************************** 4112b916ea7SJeremy L Thompson CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 412a515125bSLeila Ghaffari // Inputs 4133d65b166SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 4143d65b166SJames Wright const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1]; 4153d65b166SJames Wright const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 4163d65b166SJames Wright const CeedScalar(*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 4173d65b166SJames Wright 418a515125bSLeila Ghaffari // Outputs 4193d65b166SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 4203d65b166SJames Wright CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 421a515125bSLeila Ghaffari 422139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx; 423d8a22b9eSJed Brown const CeedScalar c_tau = context->c_tau; 424a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 425a515125bSLeila Ghaffari 426a515125bSLeila Ghaffari // Quadrature Point Loop 4273d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 428a515125bSLeila Ghaffari // Setup 429a515125bSLeila Ghaffari // -- Interp in 430a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 4312b916ea7SJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 432a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 4332b916ea7SJeremy L Thompson const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]}; 4342b916ea7SJeremy L Thompson const CeedScalar dU[3][3] = { 4352b916ea7SJeremy L Thompson {dq[0][1][i], dq[1][1][i], dq[2][1][i]}, 4362b916ea7SJeremy L Thompson {dq[0][2][i], dq[1][2][i], dq[2][2][i]}, 4372b916ea7SJeremy L Thompson {dq[0][3][i], dq[1][3][i], dq[2][3][i]} 438139613f2SLeila Ghaffari }; 4392b916ea7SJeremy L Thompson const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]}; 440a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 441a515125bSLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 442a515125bSLeila Ghaffari // -- Interp-to-Grad q_data 443a515125bSLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 4442b916ea7SJeremy L Thompson const CeedScalar dXdx[3][3] = { 4452b916ea7SJeremy L Thompson {q_data[1][i], q_data[2][i], q_data[3][i]}, 4462b916ea7SJeremy L Thompson {q_data[4][i], q_data[5][i], q_data[6][i]}, 4472b916ea7SJeremy L Thompson {q_data[7][i], q_data[8][i], q_data[9][i]} 448a515125bSLeila Ghaffari }; 449139613f2SLeila Ghaffari // dU/dx 450139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 451139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 452139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 453139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 454493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) { 455493642f1SJames Wright for (CeedInt k = 0; k < 3; k++) { 456139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 457139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 458493642f1SJames Wright for (CeedInt l = 0; l < 3; l++) { 459139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 460139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j 461139613f2SLeila Ghaffari } 462139613f2SLeila Ghaffari } 463139613f2SLeila Ghaffari } 4642b916ea7SJeremy L Thompson const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic, 465139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 466a515125bSLeila Ghaffari 467a515125bSLeila Ghaffari // The Physics 468a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 469493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) { 470139613f2SLeila Ghaffari v[j][i] = 0.; 4712b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.; 472a515125bSLeila Ghaffari } 473a515125bSLeila Ghaffari //-----mass matrix 4742b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i]; 475a515125bSLeila Ghaffari 476a515125bSLeila Ghaffari // -- Density 477a515125bSLeila Ghaffari // ---- u rho 4782b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) dv[j][0][i] -= wdetJ * (rho * u[0] * dXdx[j][0] + rho * u[1] * dXdx[j][1] + rho * u[2] * dXdx[j][2]); 479a515125bSLeila Ghaffari // -- Momentum 480a515125bSLeila Ghaffari // ---- rho (u x u) + P I3 4812b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 4822b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { 4832b916ea7SJeremy L Thompson dv[k][j + 1][i] -= wdetJ * ((rho * u[j] * u[0] + (j == 0 ? P : 0.)) * dXdx[k][0] + (rho * u[j] * u[1] + (j == 1 ? P : 0.)) * dXdx[k][1] + 484139613f2SLeila Ghaffari (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]); 4852b916ea7SJeremy L Thompson } 4862b916ea7SJeremy L Thompson } 487a515125bSLeila Ghaffari // -- Total Energy Density 488a515125bSLeila Ghaffari // ---- (E + P) u 4892b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) dv[j][4][i] -= wdetJ * (E + P) * (u[0] * dXdx[j][0] + u[1] * dXdx[j][1] + u[2] * dXdx[j][2]); 490139613f2SLeila Ghaffari 491139613f2SLeila Ghaffari // -- Stabilization terms 492139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 493139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 494d8a22b9eSJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 495139613f2SLeila Ghaffari 496139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 497139613f2SLeila Ghaffari CeedScalar dqdx[5][3]; 498493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) { 499139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j]; 500139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j]; 5012b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j]; 502139613f2SLeila Ghaffari } 503139613f2SLeila Ghaffari 504139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 505139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 5062b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 5072b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 5082b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 5092b916ea7SJeremy L Thompson } 5102b916ea7SJeremy L Thompson } 511139613f2SLeila Ghaffari 512139613f2SLeila Ghaffari // ---- Strong residual 513139613f2SLeila Ghaffari CeedScalar strong_res[5]; 5142b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j]; 515139613f2SLeila Ghaffari 516d8a22b9eSJed Brown // Stabilization 517d8a22b9eSJed Brown // -- Tau elements 518d8a22b9eSJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 519d8a22b9eSJed Brown CeedScalar Tau_x[3] = {0.}; 520d8a22b9eSJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 521139613f2SLeila Ghaffari 522d8a22b9eSJed Brown // -- Stabilization method: none, SU, or SUPG 523bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}}; 524139613f2SLeila Ghaffari switch (context->stabilization) { 525139613f2SLeila Ghaffari case 0: // Galerkin 526139613f2SLeila Ghaffari break; 527139613f2SLeila Ghaffari case 1: // SU 5282b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 5292b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 5302b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 5312b916ea7SJeremy L Thompson } 5322b916ea7SJeremy L Thompson } 533139613f2SLeila Ghaffari 5342b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 5352b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] += wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]); 5362b916ea7SJeremy L Thompson } 537139613f2SLeila Ghaffari break; 538139613f2SLeila Ghaffari case 2: // SUPG 5392b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 5402b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 5412b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l]; 5422b916ea7SJeremy L Thompson } 5432b916ea7SJeremy L Thompson } 544139613f2SLeila Ghaffari 5452b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 5462b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] += wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]); 5472b916ea7SJeremy L Thompson } 548139613f2SLeila Ghaffari break; 549139613f2SLeila Ghaffari } 550a515125bSLeila Ghaffari } // End Quadrature Point Loop 551a515125bSLeila Ghaffari 552a515125bSLeila Ghaffari // Return 553a515125bSLeila Ghaffari return 0; 554a515125bSLeila Ghaffari } 555a515125bSLeila Ghaffari // ***************************************************************************** 556*04e40bb6SJeremy L Thompson // This QFunction sets the inflow boundary conditions for the traveling vortex problem. 557a515125bSLeila Ghaffari // 558*04e40bb6SJeremy L Thompson // Prescribed T_inlet and P_inlet are converted to conservative variables and applied weakly. 559a515125bSLeila Ghaffari // ***************************************************************************** 5602b916ea7SJeremy L Thompson CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 561a515125bSLeila Ghaffari // Inputs 562dd64951cSJames Wright const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 563a515125bSLeila Ghaffari // Outputs 564a515125bSLeila Ghaffari CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 565a515125bSLeila Ghaffari EulerContext context = (EulerContext)ctx; 566a515125bSLeila Ghaffari const int euler_test = context->euler_test; 567a515125bSLeila Ghaffari const bool implicit = context->implicit; 568a515125bSLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 569a515125bSLeila Ghaffari const CeedScalar cv = 2.5; 570a515125bSLeila Ghaffari const CeedScalar R = 1.; 571a515125bSLeila Ghaffari CeedScalar T_inlet; 572a515125bSLeila Ghaffari CeedScalar P_inlet; 573a515125bSLeila Ghaffari 574a515125bSLeila Ghaffari // For test cases 1 and 3 the background velocity is zero 5752b916ea7SJeremy L Thompson if (euler_test == 1 || euler_test == 3) { 576a515125bSLeila Ghaffari for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.; 5772b916ea7SJeremy L Thompson } 578a515125bSLeila Ghaffari 579a515125bSLeila Ghaffari // For test cases 1 and 2, T_inlet = T_inlet = 0.4 580a515125bSLeila Ghaffari if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4; 581a515125bSLeila Ghaffari else T_inlet = P_inlet = 1.; 582a515125bSLeila Ghaffari 583a515125bSLeila Ghaffari // Quadrature Point Loop 5843d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 585a515125bSLeila Ghaffari // Setup 586a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 587a515125bSLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 588a515125bSLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 589a515125bSLeila Ghaffari // We can effect this by swapping the sign on this weight 590a515125bSLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 591002797a3SLeila Ghaffari // ---- Normal vect 5922b916ea7SJeremy L Thompson const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]}; 593a515125bSLeila Ghaffari 594a515125bSLeila Ghaffari // face_normal = Normal vector of the face 5952b916ea7SJeremy L Thompson const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2]; 596a515125bSLeila Ghaffari // The Physics 597a515125bSLeila Ghaffari // Zero v so all future terms can safely sum into it 598493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.; 599a515125bSLeila Ghaffari 600a515125bSLeila Ghaffari // Implementing in/outflow BCs 601002797a3SLeila Ghaffari if (face_normal > 0) { 602a515125bSLeila Ghaffari } else { // inflow 603a515125bSLeila Ghaffari const CeedScalar rho_inlet = P_inlet / (R * T_inlet); 6042b916ea7SJeremy L Thompson const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.; 605a515125bSLeila Ghaffari // incoming total energy 606a515125bSLeila Ghaffari const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet); 607a515125bSLeila Ghaffari 608a515125bSLeila Ghaffari // The Physics 609a515125bSLeila Ghaffari // -- Density 610a515125bSLeila Ghaffari v[0][i] -= wdetJb * rho_inlet * face_normal; 611a515125bSLeila Ghaffari 612a515125bSLeila Ghaffari // -- Momentum 6132b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho_inlet * face_normal * mean_velocity[j] + norm[j] * P_inlet); 614a515125bSLeila Ghaffari 615a515125bSLeila Ghaffari // -- Total Energy Density 616a515125bSLeila Ghaffari v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet); 617a515125bSLeila Ghaffari } 618a515125bSLeila Ghaffari 619a515125bSLeila Ghaffari } // End Quadrature Point Loop 620a515125bSLeila Ghaffari return 0; 621a515125bSLeila Ghaffari } 622a515125bSLeila Ghaffari 623a515125bSLeila Ghaffari // ***************************************************************************** 624*04e40bb6SJeremy L Thompson // This QFunction sets the outflow boundary conditions for the Euler solver. 62568ef3d20SLeila Ghaffari // 62668ef3d20SLeila Ghaffari // Outflow BCs: 627*04e40bb6SJeremy L Thompson // The validity of the weak form of the governing equations is extended to the outflow. 62868ef3d20SLeila Ghaffari // ***************************************************************************** 6292b916ea7SJeremy L Thompson CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 63068ef3d20SLeila Ghaffari // Inputs 6313d65b166SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 6323d65b166SJames Wright const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 6333d65b166SJames Wright 63468ef3d20SLeila Ghaffari // Outputs 63568ef3d20SLeila Ghaffari CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 63668ef3d20SLeila Ghaffari EulerContext context = (EulerContext)ctx; 63768ef3d20SLeila Ghaffari const bool implicit = context->implicit; 63868ef3d20SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 63968ef3d20SLeila Ghaffari 64068ef3d20SLeila Ghaffari const CeedScalar gamma = 1.4; 64168ef3d20SLeila Ghaffari 64268ef3d20SLeila Ghaffari // Quadrature Point Loop 6433d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 64468ef3d20SLeila Ghaffari // Setup 64568ef3d20SLeila Ghaffari // -- Interp in 64668ef3d20SLeila Ghaffari const CeedScalar rho = q[0][i]; 6472b916ea7SJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 64868ef3d20SLeila Ghaffari const CeedScalar E = q[4][i]; 64968ef3d20SLeila Ghaffari 65068ef3d20SLeila Ghaffari // -- Interp-to-Interp q_data 65168ef3d20SLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 65268ef3d20SLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 65368ef3d20SLeila Ghaffari // We can effect this by swapping the sign on this weight 65468ef3d20SLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 65568ef3d20SLeila Ghaffari // ---- Normal vectors 6562b916ea7SJeremy L Thompson const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]}; 65768ef3d20SLeila Ghaffari 65868ef3d20SLeila Ghaffari // face_normal = Normal vector of the face 6592b916ea7SJeremy L Thompson const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2]; 66068ef3d20SLeila Ghaffari // The Physics 66168ef3d20SLeila Ghaffari // Zero v so all future terms can safely sum into it 662493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0; 66368ef3d20SLeila Ghaffari 66468ef3d20SLeila Ghaffari // Implementing in/outflow BCs 66568ef3d20SLeila Ghaffari if (face_normal > 0) { // outflow 66668ef3d20SLeila Ghaffari const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.; 66768ef3d20SLeila Ghaffari const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure 6682b916ea7SJeremy L Thompson const CeedScalar u_normal = norm[0] * u[0] + norm[1] * u[1] + norm[2] * u[2]; // Normal velocity 66968ef3d20SLeila Ghaffari // The Physics 67068ef3d20SLeila Ghaffari // -- Density 67168ef3d20SLeila Ghaffari v[0][i] -= wdetJb * rho * u_normal; 67268ef3d20SLeila Ghaffari 67368ef3d20SLeila Ghaffari // -- Momentum 6742b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P); 67568ef3d20SLeila Ghaffari 67668ef3d20SLeila Ghaffari // -- Total Energy Density 67768ef3d20SLeila Ghaffari v[4][i] -= wdetJb * u_normal * (E + P); 67868ef3d20SLeila Ghaffari } 67968ef3d20SLeila Ghaffari } // End Quadrature Point Loop 68068ef3d20SLeila Ghaffari return 0; 68168ef3d20SLeila Ghaffari } 68268ef3d20SLeila Ghaffari 68368ef3d20SLeila Ghaffari // ***************************************************************************** 684a515125bSLeila Ghaffari 685a515125bSLeila Ghaffari #endif // eulervortex_h 686