13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 377841947SLeila Ghaffari // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 577841947SLeila Ghaffari // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 777841947SLeila Ghaffari 877841947SLeila Ghaffari /// @file 977841947SLeila Ghaffari /// Euler traveling vortex initial condition and operator for Navier-Stokes 1077841947SLeila Ghaffari /// example using PETSc 1177841947SLeila Ghaffari 1277841947SLeila Ghaffari // Model from: 1377841947SLeila Ghaffari // On the Order of Accuracy and Numerical Performance of Two Classes of 1477841947SLeila Ghaffari // Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011). 1577841947SLeila Ghaffari 1677841947SLeila Ghaffari #ifndef eulervortex_h 1777841947SLeila Ghaffari #define eulervortex_h 1877841947SLeila Ghaffari 1988b783a1SJames Wright #include <ceed.h> 20c9c2c079SJeremy L Thompson #include <math.h> 212b730f8bSJeremy L Thompson 2213fa47b2SJames Wright #include "utils.h" 2377841947SLeila Ghaffari 2477841947SLeila Ghaffari typedef struct EulerContext_ *EulerContext; 2577841947SLeila Ghaffari struct EulerContext_ { 2677841947SLeila Ghaffari CeedScalar center[3]; 2777841947SLeila Ghaffari CeedScalar curr_time; 2877841947SLeila Ghaffari CeedScalar vortex_strength; 29932417b3SJed Brown CeedScalar c_tau; 3077841947SLeila Ghaffari CeedScalar mean_velocity[3]; 3177841947SLeila Ghaffari bool implicit; 32e6225c47SLeila Ghaffari int euler_test; 33e6225c47SLeila Ghaffari int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG 3477841947SLeila Ghaffari }; 3577841947SLeila Ghaffari 3677841947SLeila Ghaffari // ***************************************************************************** 3777841947SLeila Ghaffari // This function sets the initial conditions 3877841947SLeila Ghaffari // 3977841947SLeila Ghaffari // Temperature: 4077841947SLeila Ghaffari // T = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2) 4177841947SLeila Ghaffari // Density: 4277841947SLeila Ghaffari // rho = (T/S_vortex)^(1 / (gamma - 1)) 4377841947SLeila Ghaffari // Pressure: 4477841947SLeila Ghaffari // P = rho * T 4577841947SLeila Ghaffari // Velocity: 4677841947SLeila Ghaffari // ui = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi) 4777841947SLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 ) 4877841947SLeila Ghaffari // Velocity/Momentum Density: 4977841947SLeila Ghaffari // Ui = rho ui 5077841947SLeila Ghaffari // Total Energy: 5177841947SLeila Ghaffari // E = P / (gamma - 1) + rho (u u)/2 5277841947SLeila Ghaffari // 5377841947SLeila Ghaffari // Constants: 5477841947SLeila Ghaffari // cv , Specific heat, constant volume 5577841947SLeila Ghaffari // cp , Specific heat, constant pressure 5677841947SLeila Ghaffari // vortex_strength , Strength of vortex 5777841947SLeila Ghaffari // center , Location of bubble center 5877841947SLeila Ghaffari // gamma = cp / cv, Specific heat ratio 5977841947SLeila Ghaffari // 6077841947SLeila Ghaffari // ***************************************************************************** 6177841947SLeila Ghaffari 6277841947SLeila Ghaffari // ***************************************************************************** 6377841947SLeila Ghaffari // This helper function provides support for the exact, time-dependent solution 6477841947SLeila Ghaffari // (currently not implemented) and IC formulation for Euler traveling vortex 6577841947SLeila Ghaffari // ***************************************************************************** 662b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) { 6777841947SLeila Ghaffari // Context 6877841947SLeila Ghaffari const EulerContext context = (EulerContext)ctx; 6977841947SLeila Ghaffari const CeedScalar vortex_strength = context->vortex_strength; 7077841947SLeila Ghaffari const CeedScalar *center = context->center; // Center of the domain 7177841947SLeila Ghaffari const CeedScalar *mean_velocity = context->mean_velocity; 7277841947SLeila Ghaffari 7377841947SLeila Ghaffari // Setup 7477841947SLeila Ghaffari const CeedScalar gamma = 1.4; 7577841947SLeila Ghaffari const CeedScalar cv = 2.5; 7677841947SLeila Ghaffari const CeedScalar R = 1.; 7777841947SLeila Ghaffari const CeedScalar x = X[0], y = X[1]; // Coordinates 7877841947SLeila Ghaffari // Vortex center 7977841947SLeila Ghaffari const CeedScalar xc = center[0] + mean_velocity[0] * time; 8077841947SLeila Ghaffari const CeedScalar yc = center[1] + mean_velocity[1] * time; 8177841947SLeila Ghaffari 8277841947SLeila Ghaffari const CeedScalar x0 = x - xc; 8377841947SLeila Ghaffari const CeedScalar y0 = y - yc; 8477841947SLeila Ghaffari const CeedScalar r = sqrt(x0 * x0 + y0 * y0); 8577841947SLeila Ghaffari const CeedScalar C = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI); 862b730f8bSJeremy L Thompson const CeedScalar delta_T = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI); 8777841947SLeila Ghaffari const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma 882b730f8bSJeremy L Thompson const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI); 8977841947SLeila Ghaffari CeedScalar rho, P, T, E, u[3] = {0.}; 9077841947SLeila Ghaffari 9177841947SLeila Ghaffari // Initial Conditions 9277841947SLeila Ghaffari switch (context->euler_test) { 9377841947SLeila Ghaffari case 0: // Traveling vortex 9477841947SLeila Ghaffari T = 1 + delta_T; 9577841947SLeila Ghaffari // P = rho * T 9677841947SLeila Ghaffari // P = S * rho^gamma 9777841947SLeila Ghaffari // Solve for rho, then substitute for P 98e6225c47SLeila Ghaffari rho = pow(T / S_vortex, 1 / (gamma - 1.)); 9977841947SLeila Ghaffari P = rho * T; 10077841947SLeila Ghaffari u[0] = mean_velocity[0] - C * y0; 10177841947SLeila Ghaffari u[1] = mean_velocity[1] + C * x0; 10277841947SLeila Ghaffari 10377841947SLeila Ghaffari // Assign exact solution 10477841947SLeila Ghaffari q[0] = rho; 10577841947SLeila Ghaffari q[1] = rho * u[0]; 10677841947SLeila Ghaffari q[2] = rho * u[1]; 10777841947SLeila Ghaffari q[3] = rho * u[2]; 10877841947SLeila Ghaffari q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.; 10977841947SLeila Ghaffari break; 11077841947SLeila Ghaffari case 1: // Constant zero velocity, density constant, total energy constant 11177841947SLeila Ghaffari rho = 1.; 11277841947SLeila Ghaffari E = 2.; 11377841947SLeila Ghaffari 11477841947SLeila Ghaffari // Assign exact solution 11577841947SLeila Ghaffari q[0] = rho; 11677841947SLeila Ghaffari q[1] = rho * u[0]; 11777841947SLeila Ghaffari q[2] = rho * u[1]; 11877841947SLeila Ghaffari q[3] = rho * u[2]; 11977841947SLeila Ghaffari q[4] = E; 12077841947SLeila Ghaffari break; 12177841947SLeila Ghaffari case 2: // Constant nonzero velocity, density constant, total energy constant 12277841947SLeila Ghaffari rho = 1.; 12377841947SLeila Ghaffari E = 2.; 12477841947SLeila Ghaffari u[0] = mean_velocity[0]; 12577841947SLeila Ghaffari u[1] = mean_velocity[1]; 12677841947SLeila Ghaffari 12777841947SLeila Ghaffari // Assign exact solution 12877841947SLeila Ghaffari q[0] = rho; 12977841947SLeila Ghaffari q[1] = rho * u[0]; 13077841947SLeila Ghaffari q[2] = rho * u[1]; 13177841947SLeila Ghaffari q[3] = rho * u[2]; 13277841947SLeila Ghaffari q[4] = E; 13377841947SLeila Ghaffari break; 13477841947SLeila Ghaffari case 3: // Velocity zero, pressure constant 13577841947SLeila Ghaffari // (so density and internal energy will be non-constant), 13677841947SLeila Ghaffari // but the velocity should stay zero and the bubble won't diffuse 13777841947SLeila Ghaffari // (for Euler, where there is no thermal conductivity) 13877841947SLeila Ghaffari P = 1.; 13977841947SLeila Ghaffari T = 1. - S_bubble * exp(1. - r * r); 14077841947SLeila Ghaffari rho = P / (R * T); 14177841947SLeila Ghaffari 14277841947SLeila Ghaffari // Assign exact solution 14377841947SLeila Ghaffari q[0] = rho; 14477841947SLeila Ghaffari q[1] = rho * u[0]; 14577841947SLeila Ghaffari q[2] = rho * u[1]; 14677841947SLeila Ghaffari q[3] = rho * u[2]; 14777841947SLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 14877841947SLeila Ghaffari break; 14977841947SLeila Ghaffari case 4: // Constant nonzero velocity, pressure constant 15077841947SLeila Ghaffari // (so density and internal energy will be non-constant), 15177841947SLeila Ghaffari // it should be transported across the domain, but velocity stays constant 15277841947SLeila Ghaffari P = 1.; 15377841947SLeila Ghaffari T = 1. - S_bubble * exp(1. - r * r); 15477841947SLeila Ghaffari rho = P / (R * T); 15577841947SLeila Ghaffari u[0] = mean_velocity[0]; 15677841947SLeila Ghaffari u[1] = mean_velocity[1]; 15777841947SLeila Ghaffari 15877841947SLeila Ghaffari // Assign exact solution 15977841947SLeila Ghaffari q[0] = rho; 16077841947SLeila Ghaffari q[1] = rho * u[0]; 16177841947SLeila Ghaffari q[2] = rho * u[1]; 16277841947SLeila Ghaffari q[3] = rho * u[2]; 16377841947SLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 16477841947SLeila Ghaffari break; 16532f166c6SLeila Ghaffari case 5: // non-smooth thermal bubble - cylinder 16632f166c6SLeila Ghaffari P = 1.; 16732f166c6SLeila Ghaffari T = 1. - (r < 1. ? S_bubble : 0.); 16832f166c6SLeila Ghaffari rho = P / (R * T); 16932f166c6SLeila Ghaffari u[0] = mean_velocity[0]; 17032f166c6SLeila Ghaffari u[1] = mean_velocity[1]; 17132f166c6SLeila Ghaffari 17232f166c6SLeila Ghaffari // Assign exact solution 17332f166c6SLeila Ghaffari q[0] = rho; 17432f166c6SLeila Ghaffari q[1] = rho * u[0]; 17532f166c6SLeila Ghaffari q[2] = rho * u[1]; 17632f166c6SLeila Ghaffari q[3] = rho * u[2]; 17732f166c6SLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.); 17832f166c6SLeila Ghaffari break; 17977841947SLeila Ghaffari } 18077841947SLeila Ghaffari // Return 18177841947SLeila Ghaffari return 0; 18277841947SLeila Ghaffari } 18377841947SLeila Ghaffari 18477841947SLeila Ghaffari // ***************************************************************************** 185e6225c47SLeila Ghaffari // Helper function for computing flux Jacobian 186e6225c47SLeila Ghaffari // ***************************************************************************** 1872b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 188e6225c47SLeila Ghaffari const CeedScalar gamma) { 189e6225c47SLeila Ghaffari CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2]; // Velocity square 190e6225c47SLeila Ghaffari for (CeedInt i = 0; i < 3; i++) { // Jacobian matrices for 3 directions 191e6225c47SLeila Ghaffari for (CeedInt j = 0; j < 3; j++) { // Rows of each Jacobian matrix 192e6225c47SLeila Ghaffari dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j]; 193e6225c47SLeila Ghaffari for (CeedInt k = 0; k < 3; k++) { // Columns of each Jacobian matrix 194e6225c47SLeila Ghaffari dF[i][0][k + 1] = ((i == k) ? 1. : 0.); 1952b730f8bSJeremy 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.); 1962b730f8bSJeremy L Thompson dF[i][4][k + 1] = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k]; 197e6225c47SLeila Ghaffari } 198e6225c47SLeila Ghaffari dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.); 199e6225c47SLeila Ghaffari } 200e6225c47SLeila Ghaffari dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho); 201e6225c47SLeila Ghaffari dF[i][4][4] = u[i] * gamma; 202e6225c47SLeila Ghaffari } 203e6225c47SLeila Ghaffari } 204e6225c47SLeila Ghaffari 205e6225c47SLeila Ghaffari // ***************************************************************************** 206932417b3SJed Brown // Helper function for computing Tau elements (stabilization constant) 207932417b3SJed Brown // Model from: 208932417b3SJed Brown // Stabilized Methods for Compressible Flows, Hughes et al 2010 209932417b3SJed Brown // 210932417b3SJed Brown // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 211932417b3SJed Brown // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 212932417b3SJed Brown // 213932417b3SJed Brown // Where 214932417b3SJed Brown // c_tau = stabilization constant (0.5 is reported as "optimal") 215932417b3SJed Brown // h[i] = 2 length(dxdX[i]) 216932417b3SJed Brown // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 217932417b3SJed Brown // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 218932417b3SJed Brown // rho(A[i]) = spectral radius of the convective flux Jacobian i, 219932417b3SJed Brown // wave speed in direction i 220932417b3SJed Brown // ***************************************************************************** 2212b730f8bSJeremy 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, 2222b730f8bSJeremy L Thompson const CeedScalar c_tau) { 223ba6664aeSJames Wright for (CeedInt i = 0; i < 3; i++) { 224932417b3SJed Brown // length of element in direction i 2252b730f8bSJeremy L Thompson CeedScalar h = 2 / sqrt(dXdx[0][i] * dXdx[0][i] + dXdx[1][i] * dXdx[1][i] + dXdx[2][i] * dXdx[2][i]); 226932417b3SJed Brown // fastest wave in direction i 227932417b3SJed Brown CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 228932417b3SJed Brown Tau_x[i] = c_tau * h / fastest_wave; 229932417b3SJed Brown } 230932417b3SJed Brown } 231932417b3SJed Brown 232932417b3SJed Brown // ***************************************************************************** 23377841947SLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex 23477841947SLeila Ghaffari // ***************************************************************************** 2352b730f8bSJeremy L Thompson CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 23677841947SLeila Ghaffari // Inputs 23777841947SLeila Ghaffari const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 23877841947SLeila Ghaffari 23977841947SLeila Ghaffari // Outputs 24077841947SLeila Ghaffari CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 24177841947SLeila Ghaffari const EulerContext context = (EulerContext)ctx; 24277841947SLeila Ghaffari 24377841947SLeila Ghaffari // Quadrature Point Loop 244*46603fc5SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 24577841947SLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 246e6225c47SLeila Ghaffari CeedScalar q[5] = {0.}; 24777841947SLeila Ghaffari 24877841947SLeila Ghaffari Exact_Euler(3, context->curr_time, x, 5, q, ctx); 24977841947SLeila Ghaffari 2502b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 25177841947SLeila Ghaffari } // End of Quadrature Point Loop 25277841947SLeila Ghaffari 25377841947SLeila Ghaffari // Return 25477841947SLeila Ghaffari return 0; 25577841947SLeila Ghaffari } 25677841947SLeila Ghaffari 25777841947SLeila Ghaffari // ***************************************************************************** 25877841947SLeila Ghaffari // This QFunction implements the following formulation of Euler equations 25977841947SLeila Ghaffari // with explicit time stepping method 26077841947SLeila Ghaffari // 26177841947SLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation 26277841947SLeila Ghaffari // form with state variables of density, momentum density, and total 26377841947SLeila Ghaffari // energy density. 26477841947SLeila Ghaffari // 26577841947SLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E ) 26677841947SLeila Ghaffari // rho - Mass Density 26777841947SLeila Ghaffari // Ui - Momentum Density, Ui = rho ui 26877841947SLeila Ghaffari // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2 26977841947SLeila Ghaffari // 27077841947SLeila Ghaffari // Euler Equations: 27177841947SLeila Ghaffari // drho/dt + div( U ) = 0 27277841947SLeila Ghaffari // dU/dt + div( rho (u x u) + P I3 ) = 0 27377841947SLeila Ghaffari // dE/dt + div( (E + P) u ) = 0 27477841947SLeila Ghaffari // 27577841947SLeila Ghaffari // Equation of State: 27677841947SLeila Ghaffari // P = (gamma - 1) (E - rho (u u) / 2) 27777841947SLeila Ghaffari // 27877841947SLeila Ghaffari // Constants: 27977841947SLeila Ghaffari // cv , Specific heat, constant volume 28077841947SLeila Ghaffari // cp , Specific heat, constant pressure 28177841947SLeila Ghaffari // g , Gravity 28277841947SLeila Ghaffari // gamma = cp / cv, Specific heat ratio 28377841947SLeila Ghaffari // ***************************************************************************** 2842b730f8bSJeremy L Thompson CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 28577841947SLeila Ghaffari // Inputs 286*46603fc5SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 287*46603fc5SJames Wright const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1]; 288*46603fc5SJames Wright const CeedScalar(*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 289*46603fc5SJames Wright 29077841947SLeila Ghaffari // Outputs 291*46603fc5SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 292*46603fc5SJames Wright CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 29377841947SLeila Ghaffari 294e6225c47SLeila Ghaffari EulerContext context = (EulerContext)ctx; 295932417b3SJed Brown const CeedScalar c_tau = context->c_tau; 29677841947SLeila Ghaffari const CeedScalar gamma = 1.4; 29777841947SLeila Ghaffari 29877841947SLeila Ghaffari // Quadrature Point Loop 299*46603fc5SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 30077841947SLeila Ghaffari // Setup 30177841947SLeila Ghaffari // -- Interp in 30277841947SLeila Ghaffari const CeedScalar rho = q[0][i]; 3032b730f8bSJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 30477841947SLeila Ghaffari const CeedScalar E = q[4][i]; 3052b730f8bSJeremy L Thompson const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]}; 3062b730f8bSJeremy L Thompson const CeedScalar dU[3][3] = { 3072b730f8bSJeremy L Thompson {dq[0][1][i], dq[1][1][i], dq[2][1][i]}, 3082b730f8bSJeremy L Thompson {dq[0][2][i], dq[1][2][i], dq[2][2][i]}, 3092b730f8bSJeremy L Thompson {dq[0][3][i], dq[1][3][i], dq[2][3][i]} 310e6225c47SLeila Ghaffari }; 3112b730f8bSJeremy L Thompson const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]}; 31277841947SLeila Ghaffari // -- Interp-to-Interp q_data 31377841947SLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 31477841947SLeila Ghaffari // -- Interp-to-Grad q_data 31577841947SLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 3162b730f8bSJeremy L Thompson const CeedScalar dXdx[3][3] = { 3172b730f8bSJeremy L Thompson {q_data[1][i], q_data[2][i], q_data[3][i]}, 3182b730f8bSJeremy L Thompson {q_data[4][i], q_data[5][i], q_data[6][i]}, 3192b730f8bSJeremy L Thompson {q_data[7][i], q_data[8][i], q_data[9][i]} 32077841947SLeila Ghaffari }; 321e6225c47SLeila Ghaffari // dU/dx 322e6225c47SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 323e6225c47SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 324e6225c47SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 325e6225c47SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 326ba6664aeSJames Wright for (CeedInt j = 0; j < 3; j++) { 327ba6664aeSJames Wright for (CeedInt k = 0; k < 3; k++) { 328e6225c47SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 329e6225c47SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 330ba6664aeSJames Wright for (CeedInt l = 0; l < 3; l++) { 331e6225c47SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 332e6225c47SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j 333e6225c47SLeila Ghaffari } 334e6225c47SLeila Ghaffari } 335e6225c47SLeila Ghaffari } 336e6225c47SLeila Ghaffari // Pressure 3372b730f8bSJeremy 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, 338e6225c47SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 33977841947SLeila Ghaffari 34077841947SLeila Ghaffari // The Physics 34177841947SLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 342ba6664aeSJames Wright for (CeedInt j = 0; j < 5; j++) { 343e6225c47SLeila Ghaffari v[j][i] = 0.; 3442b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.; 34577841947SLeila Ghaffari } 34677841947SLeila Ghaffari 34777841947SLeila Ghaffari // -- Density 34877841947SLeila Ghaffari // ---- u rho 3492b730f8bSJeremy 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]); 35077841947SLeila Ghaffari // -- Momentum 35177841947SLeila Ghaffari // ---- rho (u x u) + P I3 3522b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3532b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { 3542b730f8bSJeremy 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] + 355e6225c47SLeila Ghaffari (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]); 3562b730f8bSJeremy L Thompson } 3572b730f8bSJeremy L Thompson } 35877841947SLeila Ghaffari // -- Total Energy Density 35977841947SLeila Ghaffari // ---- (E + P) u 3602b730f8bSJeremy 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]); 361e6225c47SLeila Ghaffari 362e6225c47SLeila Ghaffari // --Stabilization terms 363e6225c47SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 364e6225c47SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 365932417b3SJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 366e6225c47SLeila Ghaffari 367e6225c47SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 368e6225c47SLeila Ghaffari CeedScalar dqdx[5][3]; 369ba6664aeSJames Wright for (CeedInt j = 0; j < 3; j++) { 370e6225c47SLeila Ghaffari dqdx[0][j] = drhodx[j]; 371e6225c47SLeila Ghaffari dqdx[4][j] = dEdx[j]; 3722b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j]; 373e6225c47SLeila Ghaffari } 374e6225c47SLeila Ghaffari 375e6225c47SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 376e6225c47SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 3772b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3782b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 3792b730f8bSJeremy L Thompson for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 3802b730f8bSJeremy L Thompson } 3812b730f8bSJeremy L Thompson } 382e6225c47SLeila Ghaffari 383932417b3SJed Brown // Stabilization 384932417b3SJed Brown // -- Tau elements 385932417b3SJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 386932417b3SJed Brown CeedScalar Tau_x[3] = {0.}; 387932417b3SJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 388e6225c47SLeila Ghaffari 389932417b3SJed Brown // -- Stabilization method: none or SU 39088626eedSJames Wright CeedScalar stab[5][3] = {{0.}}; 391e6225c47SLeila Ghaffari switch (context->stabilization) { 392e6225c47SLeila Ghaffari case 0: // Galerkin 393e6225c47SLeila Ghaffari break; 394e6225c47SLeila Ghaffari case 1: // SU 3952b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3962b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 3972b730f8bSJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 3982b730f8bSJeremy L Thompson } 3992b730f8bSJeremy L Thompson } 400e6225c47SLeila Ghaffari 4012b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 4022b730f8bSJeremy 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]); 4032b730f8bSJeremy L Thompson } 404e6225c47SLeila Ghaffari break; 405e6225c47SLeila Ghaffari case 2: // SUPG is not implemented for explicit scheme 406e6225c47SLeila Ghaffari break; 407e6225c47SLeila Ghaffari } 408e6225c47SLeila Ghaffari 40977841947SLeila Ghaffari } // End Quadrature Point Loop 41077841947SLeila Ghaffari 41177841947SLeila Ghaffari // Return 41277841947SLeila Ghaffari return 0; 41377841947SLeila Ghaffari } 41477841947SLeila Ghaffari 41577841947SLeila Ghaffari // ***************************************************************************** 41677841947SLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above) 41777841947SLeila Ghaffari // with implicit time stepping method 41877841947SLeila Ghaffari // 41977841947SLeila Ghaffari // ***************************************************************************** 4202b730f8bSJeremy L Thompson CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 42177841947SLeila Ghaffari // Inputs 422*46603fc5SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 423*46603fc5SJames Wright const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1]; 424*46603fc5SJames Wright const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 425*46603fc5SJames Wright const CeedScalar(*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 426*46603fc5SJames Wright 42777841947SLeila Ghaffari // Outputs 428*46603fc5SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 429*46603fc5SJames Wright CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 43077841947SLeila Ghaffari 431e6225c47SLeila Ghaffari EulerContext context = (EulerContext)ctx; 432932417b3SJed Brown const CeedScalar c_tau = context->c_tau; 43377841947SLeila Ghaffari const CeedScalar gamma = 1.4; 43477841947SLeila Ghaffari 43577841947SLeila Ghaffari // Quadrature Point Loop 436*46603fc5SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 43777841947SLeila Ghaffari // Setup 43877841947SLeila Ghaffari // -- Interp in 43977841947SLeila Ghaffari const CeedScalar rho = q[0][i]; 4402b730f8bSJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 44177841947SLeila Ghaffari const CeedScalar E = q[4][i]; 4422b730f8bSJeremy L Thompson const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]}; 4432b730f8bSJeremy L Thompson const CeedScalar dU[3][3] = { 4442b730f8bSJeremy L Thompson {dq[0][1][i], dq[1][1][i], dq[2][1][i]}, 4452b730f8bSJeremy L Thompson {dq[0][2][i], dq[1][2][i], dq[2][2][i]}, 4462b730f8bSJeremy L Thompson {dq[0][3][i], dq[1][3][i], dq[2][3][i]} 447e6225c47SLeila Ghaffari }; 4482b730f8bSJeremy L Thompson const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]}; 44977841947SLeila Ghaffari // -- Interp-to-Interp q_data 45077841947SLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 45177841947SLeila Ghaffari // -- Interp-to-Grad q_data 45277841947SLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 4532b730f8bSJeremy L Thompson const CeedScalar dXdx[3][3] = { 4542b730f8bSJeremy L Thompson {q_data[1][i], q_data[2][i], q_data[3][i]}, 4552b730f8bSJeremy L Thompson {q_data[4][i], q_data[5][i], q_data[6][i]}, 4562b730f8bSJeremy L Thompson {q_data[7][i], q_data[8][i], q_data[9][i]} 45777841947SLeila Ghaffari }; 458e6225c47SLeila Ghaffari // dU/dx 459e6225c47SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 460e6225c47SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 461e6225c47SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 462e6225c47SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 463ba6664aeSJames Wright for (CeedInt j = 0; j < 3; j++) { 464ba6664aeSJames Wright for (CeedInt k = 0; k < 3; k++) { 465e6225c47SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 466e6225c47SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 467ba6664aeSJames Wright for (CeedInt l = 0; l < 3; l++) { 468e6225c47SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 469e6225c47SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j 470e6225c47SLeila Ghaffari } 471e6225c47SLeila Ghaffari } 472e6225c47SLeila Ghaffari } 4732b730f8bSJeremy 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, 474e6225c47SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 47577841947SLeila Ghaffari 47677841947SLeila Ghaffari // The Physics 47777841947SLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 478ba6664aeSJames Wright for (CeedInt j = 0; j < 5; j++) { 479e6225c47SLeila Ghaffari v[j][i] = 0.; 4802b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.; 48177841947SLeila Ghaffari } 48277841947SLeila Ghaffari //-----mass matrix 4832b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i]; 48477841947SLeila Ghaffari 48577841947SLeila Ghaffari // -- Density 48677841947SLeila Ghaffari // ---- u rho 4872b730f8bSJeremy 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]); 48877841947SLeila Ghaffari // -- Momentum 48977841947SLeila Ghaffari // ---- rho (u x u) + P I3 4902b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 4912b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { 4922b730f8bSJeremy 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] + 493e6225c47SLeila Ghaffari (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]); 4942b730f8bSJeremy L Thompson } 4952b730f8bSJeremy L Thompson } 49677841947SLeila Ghaffari // -- Total Energy Density 49777841947SLeila Ghaffari // ---- (E + P) u 4982b730f8bSJeremy 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]); 499e6225c47SLeila Ghaffari 500e6225c47SLeila Ghaffari // -- Stabilization terms 501e6225c47SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 502e6225c47SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 503932417b3SJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 504e6225c47SLeila Ghaffari 505e6225c47SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 506e6225c47SLeila Ghaffari CeedScalar dqdx[5][3]; 507ba6664aeSJames Wright for (CeedInt j = 0; j < 3; j++) { 508e6225c47SLeila Ghaffari dqdx[0][j] = drhodx[j]; 509e6225c47SLeila Ghaffari dqdx[4][j] = dEdx[j]; 5102b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j]; 511e6225c47SLeila Ghaffari } 512e6225c47SLeila Ghaffari 513e6225c47SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 514e6225c47SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 5152b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 5162b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 5172b730f8bSJeremy L Thompson for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 5182b730f8bSJeremy L Thompson } 5192b730f8bSJeremy L Thompson } 520e6225c47SLeila Ghaffari 521e6225c47SLeila Ghaffari // ---- Strong residual 522e6225c47SLeila Ghaffari CeedScalar strong_res[5]; 5232b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j]; 524e6225c47SLeila Ghaffari 525932417b3SJed Brown // Stabilization 526932417b3SJed Brown // -- Tau elements 527932417b3SJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 528932417b3SJed Brown CeedScalar Tau_x[3] = {0.}; 529932417b3SJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 530e6225c47SLeila Ghaffari 531932417b3SJed Brown // -- Stabilization method: none, SU, or SUPG 53288626eedSJames Wright CeedScalar stab[5][3] = {{0.}}; 533e6225c47SLeila Ghaffari switch (context->stabilization) { 534e6225c47SLeila Ghaffari case 0: // Galerkin 535e6225c47SLeila Ghaffari break; 536e6225c47SLeila Ghaffari case 1: // SU 5372b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 5382b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 5392b730f8bSJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 5402b730f8bSJeremy L Thompson } 5412b730f8bSJeremy L Thompson } 542e6225c47SLeila Ghaffari 5432b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 5442b730f8bSJeremy 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]); 5452b730f8bSJeremy L Thompson } 546e6225c47SLeila Ghaffari break; 547e6225c47SLeila Ghaffari case 2: // SUPG 5482b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 5492b730f8bSJeremy L Thompson for (CeedInt k = 0; k < 5; k++) { 5502b730f8bSJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l]; 5512b730f8bSJeremy L Thompson } 5522b730f8bSJeremy L Thompson } 553e6225c47SLeila Ghaffari 5542b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 5; j++) { 5552b730f8bSJeremy 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]); 5562b730f8bSJeremy L Thompson } 557e6225c47SLeila Ghaffari break; 558e6225c47SLeila Ghaffari } 55977841947SLeila Ghaffari } // End Quadrature Point Loop 56077841947SLeila Ghaffari 56177841947SLeila Ghaffari // Return 56277841947SLeila Ghaffari return 0; 56377841947SLeila Ghaffari } 56477841947SLeila Ghaffari // ***************************************************************************** 5652fe7aee7SLeila Ghaffari // This QFunction sets the inflow boundary conditions for 5662fe7aee7SLeila Ghaffari // the traveling vortex problem. 56777841947SLeila Ghaffari // 56877841947SLeila Ghaffari // Prescribed T_inlet and P_inlet are converted to conservative variables 56977841947SLeila Ghaffari // and applied weakly. 57077841947SLeila Ghaffari // 57177841947SLeila Ghaffari // ***************************************************************************** 5722b730f8bSJeremy L Thompson CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 57377841947SLeila Ghaffari // Inputs 574e8b03feeSJames Wright const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 57577841947SLeila Ghaffari // Outputs 57677841947SLeila Ghaffari CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 57777841947SLeila Ghaffari EulerContext context = (EulerContext)ctx; 57877841947SLeila Ghaffari const int euler_test = context->euler_test; 57977841947SLeila Ghaffari const bool implicit = context->implicit; 58077841947SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 58177841947SLeila Ghaffari const CeedScalar cv = 2.5; 58277841947SLeila Ghaffari const CeedScalar R = 1.; 58377841947SLeila Ghaffari CeedScalar T_inlet; 58477841947SLeila Ghaffari CeedScalar P_inlet; 58577841947SLeila Ghaffari 58677841947SLeila Ghaffari // For test cases 1 and 3 the background velocity is zero 5872b730f8bSJeremy L Thompson if (euler_test == 1 || euler_test == 3) { 58877841947SLeila Ghaffari for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.; 5892b730f8bSJeremy L Thompson } 59077841947SLeila Ghaffari 59177841947SLeila Ghaffari // For test cases 1 and 2, T_inlet = T_inlet = 0.4 59277841947SLeila Ghaffari if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4; 59377841947SLeila Ghaffari else T_inlet = P_inlet = 1.; 59477841947SLeila Ghaffari 59577841947SLeila Ghaffari // Quadrature Point Loop 596*46603fc5SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 59777841947SLeila Ghaffari // Setup 59877841947SLeila Ghaffari // -- Interp-to-Interp q_data 59977841947SLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 60077841947SLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 60177841947SLeila Ghaffari // We can effect this by swapping the sign on this weight 60277841947SLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 6032fe7aee7SLeila Ghaffari // ---- Normal vect 6042b730f8bSJeremy L Thompson const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]}; 60577841947SLeila Ghaffari 60677841947SLeila Ghaffari // face_normal = Normal vector of the face 6072b730f8bSJeremy L Thompson const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2]; 60877841947SLeila Ghaffari // The Physics 60977841947SLeila Ghaffari // Zero v so all future terms can safely sum into it 610ba6664aeSJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.; 61177841947SLeila Ghaffari 61277841947SLeila Ghaffari // Implementing in/outflow BCs 6132fe7aee7SLeila Ghaffari if (face_normal > 0) { 61477841947SLeila Ghaffari } else { // inflow 61577841947SLeila Ghaffari const CeedScalar rho_inlet = P_inlet / (R * T_inlet); 6162b730f8bSJeremy L Thompson const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.; 61777841947SLeila Ghaffari // incoming total energy 61877841947SLeila Ghaffari const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet); 61977841947SLeila Ghaffari 62077841947SLeila Ghaffari // The Physics 62177841947SLeila Ghaffari // -- Density 62277841947SLeila Ghaffari v[0][i] -= wdetJb * rho_inlet * face_normal; 62377841947SLeila Ghaffari 62477841947SLeila Ghaffari // -- Momentum 6252b730f8bSJeremy 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); 62677841947SLeila Ghaffari 62777841947SLeila Ghaffari // -- Total Energy Density 62877841947SLeila Ghaffari v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet); 62977841947SLeila Ghaffari } 63077841947SLeila Ghaffari 63177841947SLeila Ghaffari } // End Quadrature Point Loop 63277841947SLeila Ghaffari return 0; 63377841947SLeila Ghaffari } 63477841947SLeila Ghaffari 63577841947SLeila Ghaffari // ***************************************************************************** 63655e76554SLeila Ghaffari // This QFunction sets the outflow boundary conditions for 63755e76554SLeila Ghaffari // the Euler solver. 63855e76554SLeila Ghaffari // 63955e76554SLeila Ghaffari // Outflow BCs: 64055e76554SLeila Ghaffari // The validity of the weak form of the governing equations is 64155e76554SLeila Ghaffari // extended to the outflow. 64255e76554SLeila Ghaffari // 64355e76554SLeila Ghaffari // ***************************************************************************** 6442b730f8bSJeremy L Thompson CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 64555e76554SLeila Ghaffari // Inputs 646*46603fc5SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 647*46603fc5SJames Wright const CeedScalar(*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 648*46603fc5SJames Wright 64955e76554SLeila Ghaffari // Outputs 65055e76554SLeila Ghaffari CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 65155e76554SLeila Ghaffari EulerContext context = (EulerContext)ctx; 65255e76554SLeila Ghaffari const bool implicit = context->implicit; 65355e76554SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 65455e76554SLeila Ghaffari 65555e76554SLeila Ghaffari const CeedScalar gamma = 1.4; 65655e76554SLeila Ghaffari 65755e76554SLeila Ghaffari // Quadrature Point Loop 658*46603fc5SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 65955e76554SLeila Ghaffari // Setup 66055e76554SLeila Ghaffari // -- Interp in 66155e76554SLeila Ghaffari const CeedScalar rho = q[0][i]; 6622b730f8bSJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho}; 66355e76554SLeila Ghaffari const CeedScalar E = q[4][i]; 66455e76554SLeila Ghaffari 66555e76554SLeila Ghaffari // -- Interp-to-Interp q_data 66655e76554SLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 66755e76554SLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 66855e76554SLeila Ghaffari // We can effect this by swapping the sign on this weight 66955e76554SLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 67055e76554SLeila Ghaffari // ---- Normal vectors 6712b730f8bSJeremy L Thompson const CeedScalar norm[3] = {q_data_sur[1][i], q_data_sur[2][i], q_data_sur[3][i]}; 67255e76554SLeila Ghaffari 67355e76554SLeila Ghaffari // face_normal = Normal vector of the face 6742b730f8bSJeremy L Thompson const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2]; 67555e76554SLeila Ghaffari // The Physics 67655e76554SLeila Ghaffari // Zero v so all future terms can safely sum into it 677ba6664aeSJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0; 67855e76554SLeila Ghaffari 67955e76554SLeila Ghaffari // Implementing in/outflow BCs 68055e76554SLeila Ghaffari if (face_normal > 0) { // outflow 68155e76554SLeila Ghaffari const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.; 68255e76554SLeila Ghaffari const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure 6832b730f8bSJeremy L Thompson const CeedScalar u_normal = norm[0] * u[0] + norm[1] * u[1] + norm[2] * u[2]; // Normal velocity 68455e76554SLeila Ghaffari // The Physics 68555e76554SLeila Ghaffari // -- Density 68655e76554SLeila Ghaffari v[0][i] -= wdetJb * rho * u_normal; 68755e76554SLeila Ghaffari 68855e76554SLeila Ghaffari // -- Momentum 6892b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P); 69055e76554SLeila Ghaffari 69155e76554SLeila Ghaffari // -- Total Energy Density 69255e76554SLeila Ghaffari v[4][i] -= wdetJb * u_normal * (E + P); 69355e76554SLeila Ghaffari } 69455e76554SLeila Ghaffari } // End Quadrature Point Loop 69555e76554SLeila Ghaffari return 0; 69655e76554SLeila Ghaffari } 69755e76554SLeila Ghaffari 69855e76554SLeila Ghaffari // ***************************************************************************** 69977841947SLeila Ghaffari 70077841947SLeila Ghaffari #endif // eulervortex_h 701