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 1977841947SLeila Ghaffari #include <math.h> 2088b783a1SJames Wright #include <ceed.h> 21*13fa47b2SJames Wright #include "utils.h" 2277841947SLeila Ghaffari 2377841947SLeila Ghaffari typedef struct EulerContext_ *EulerContext; 2477841947SLeila Ghaffari struct EulerContext_ { 2577841947SLeila Ghaffari CeedScalar center[3]; 2677841947SLeila Ghaffari CeedScalar curr_time; 2777841947SLeila Ghaffari CeedScalar vortex_strength; 28932417b3SJed Brown CeedScalar c_tau; 2977841947SLeila Ghaffari CeedScalar mean_velocity[3]; 3077841947SLeila Ghaffari bool implicit; 31e6225c47SLeila Ghaffari int euler_test; 32e6225c47SLeila Ghaffari int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG 3377841947SLeila Ghaffari }; 3477841947SLeila Ghaffari 3577841947SLeila Ghaffari // ***************************************************************************** 3677841947SLeila Ghaffari // This function sets the initial conditions 3777841947SLeila Ghaffari // 3877841947SLeila Ghaffari // Temperature: 3977841947SLeila Ghaffari // T = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2) 4077841947SLeila Ghaffari // Density: 4177841947SLeila Ghaffari // rho = (T/S_vortex)^(1 / (gamma - 1)) 4277841947SLeila Ghaffari // Pressure: 4377841947SLeila Ghaffari // P = rho * T 4477841947SLeila Ghaffari // Velocity: 4577841947SLeila Ghaffari // ui = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi) 4677841947SLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 ) 4777841947SLeila Ghaffari // Velocity/Momentum Density: 4877841947SLeila Ghaffari // Ui = rho ui 4977841947SLeila Ghaffari // Total Energy: 5077841947SLeila Ghaffari // E = P / (gamma - 1) + rho (u u)/2 5177841947SLeila Ghaffari // 5277841947SLeila Ghaffari // Constants: 5377841947SLeila Ghaffari // cv , Specific heat, constant volume 5477841947SLeila Ghaffari // cp , Specific heat, constant pressure 5577841947SLeila Ghaffari // vortex_strength , Strength of vortex 5677841947SLeila Ghaffari // center , Location of bubble center 5777841947SLeila Ghaffari // gamma = cp / cv, Specific heat ratio 5877841947SLeila Ghaffari // 5977841947SLeila Ghaffari // ***************************************************************************** 6077841947SLeila Ghaffari 6177841947SLeila Ghaffari // ***************************************************************************** 6277841947SLeila Ghaffari // This helper function provides support for the exact, time-dependent solution 6377841947SLeila Ghaffari // (currently not implemented) and IC formulation for Euler traveling vortex 6477841947SLeila Ghaffari // ***************************************************************************** 6577841947SLeila Ghaffari CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, 6677841947SLeila Ghaffari const CeedScalar X[], CeedInt Nf, CeedScalar q[], 6777841947SLeila Ghaffari void *ctx) { 6877841947SLeila Ghaffari // Context 6977841947SLeila Ghaffari const EulerContext context = (EulerContext)ctx; 7077841947SLeila Ghaffari const CeedScalar vortex_strength = context->vortex_strength; 7177841947SLeila Ghaffari const CeedScalar *center = context->center; // Center of the domain 7277841947SLeila Ghaffari const CeedScalar *mean_velocity = context->mean_velocity; 7377841947SLeila Ghaffari 7477841947SLeila Ghaffari // Setup 7577841947SLeila Ghaffari const CeedScalar gamma = 1.4; 7677841947SLeila Ghaffari const CeedScalar cv = 2.5; 7777841947SLeila Ghaffari const CeedScalar R = 1.; 7877841947SLeila Ghaffari const CeedScalar x = X[0], y = X[1]; // Coordinates 7977841947SLeila Ghaffari // Vortex center 8077841947SLeila Ghaffari const CeedScalar xc = center[0] + mean_velocity[0] * time; 8177841947SLeila Ghaffari const CeedScalar yc = center[1] + mean_velocity[1] * time; 8277841947SLeila Ghaffari 8377841947SLeila Ghaffari const CeedScalar x0 = x - xc; 8477841947SLeila Ghaffari const CeedScalar y0 = y - yc; 8577841947SLeila Ghaffari const CeedScalar r = sqrt( x0*x0 + y0*y0 ); 8677841947SLeila Ghaffari const CeedScalar C = vortex_strength * exp((1. - r*r)/2.) / (2. * M_PI); 87e6225c47SLeila Ghaffari const CeedScalar delta_T = - (gamma - 1.) * vortex_strength * vortex_strength * 88e6225c47SLeila Ghaffari exp(1 - r*r) / (8. * gamma * M_PI * M_PI); 8977841947SLeila Ghaffari const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma 9077841947SLeila Ghaffari const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / 9177841947SLeila Ghaffari (8.*gamma*M_PI*M_PI); 9277841947SLeila Ghaffari CeedScalar rho, P, T, E, u[3] = {0.}; 9377841947SLeila Ghaffari 9477841947SLeila Ghaffari // Initial Conditions 9577841947SLeila Ghaffari switch (context->euler_test) { 9677841947SLeila Ghaffari case 0: // Traveling vortex 9777841947SLeila Ghaffari T = 1 + delta_T; 9877841947SLeila Ghaffari // P = rho * T 9977841947SLeila Ghaffari // P = S * rho^gamma 10077841947SLeila Ghaffari // Solve for rho, then substitute for P 101e6225c47SLeila Ghaffari rho = pow(T/S_vortex, 1 / (gamma - 1.)); 10277841947SLeila Ghaffari P = rho * T; 10377841947SLeila Ghaffari u[0] = mean_velocity[0] - C*y0; 10477841947SLeila Ghaffari u[1] = mean_velocity[1] + C*x0; 10577841947SLeila Ghaffari 10677841947SLeila Ghaffari // Assign exact solution 10777841947SLeila Ghaffari q[0] = rho; 10877841947SLeila Ghaffari q[1] = rho * u[0]; 10977841947SLeila Ghaffari q[2] = rho * u[1]; 11077841947SLeila Ghaffari q[3] = rho * u[2]; 11177841947SLeila Ghaffari q[4] = P / (gamma - 1.) + rho * (u[0]*u[0] + u[1]*u[1]) / 2.; 11277841947SLeila Ghaffari break; 11377841947SLeila Ghaffari case 1: // Constant zero velocity, density constant, total energy constant 11477841947SLeila Ghaffari rho = 1.; 11577841947SLeila Ghaffari E = 2.; 11677841947SLeila Ghaffari 11777841947SLeila Ghaffari // Assign exact solution 11877841947SLeila Ghaffari q[0] = rho; 11977841947SLeila Ghaffari q[1] = rho * u[0]; 12077841947SLeila Ghaffari q[2] = rho * u[1]; 12177841947SLeila Ghaffari q[3] = rho * u[2]; 12277841947SLeila Ghaffari q[4] = E; 12377841947SLeila Ghaffari break; 12477841947SLeila Ghaffari case 2: // Constant nonzero velocity, density constant, total energy constant 12577841947SLeila Ghaffari rho = 1.; 12677841947SLeila Ghaffari E = 2.; 12777841947SLeila Ghaffari u[0] = mean_velocity[0]; 12877841947SLeila Ghaffari u[1] = mean_velocity[1]; 12977841947SLeila Ghaffari 13077841947SLeila Ghaffari // Assign exact solution 13177841947SLeila Ghaffari q[0] = rho; 13277841947SLeila Ghaffari q[1] = rho * u[0]; 13377841947SLeila Ghaffari q[2] = rho * u[1]; 13477841947SLeila Ghaffari q[3] = rho * u[2]; 13577841947SLeila Ghaffari q[4] = E; 13677841947SLeila Ghaffari break; 13777841947SLeila Ghaffari case 3: // Velocity zero, pressure constant 13877841947SLeila Ghaffari // (so density and internal energy will be non-constant), 13977841947SLeila Ghaffari // but the velocity should stay zero and the bubble won't diffuse 14077841947SLeila Ghaffari // (for Euler, where there is no thermal conductivity) 14177841947SLeila Ghaffari P = 1.; 14277841947SLeila Ghaffari T = 1. - S_bubble * exp(1. - r*r); 14377841947SLeila Ghaffari rho = P / (R*T); 14477841947SLeila Ghaffari 14577841947SLeila Ghaffari // Assign exact solution 14677841947SLeila Ghaffari q[0] = rho; 14777841947SLeila Ghaffari q[1] = rho * u[0]; 14877841947SLeila Ghaffari q[2] = rho * u[1]; 14977841947SLeila Ghaffari q[3] = rho * u[2]; 15077841947SLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 15177841947SLeila Ghaffari break; 15277841947SLeila Ghaffari case 4: // Constant nonzero velocity, pressure constant 15377841947SLeila Ghaffari // (so density and internal energy will be non-constant), 15477841947SLeila Ghaffari // it should be transported across the domain, but velocity stays constant 15577841947SLeila Ghaffari P = 1.; 15677841947SLeila Ghaffari T = 1. - S_bubble * exp(1. - r*r); 15777841947SLeila Ghaffari rho = P / (R*T); 15877841947SLeila Ghaffari u[0] = mean_velocity[0]; 15977841947SLeila Ghaffari u[1] = mean_velocity[1]; 16077841947SLeila Ghaffari 16177841947SLeila Ghaffari // Assign exact solution 16277841947SLeila Ghaffari q[0] = rho; 16377841947SLeila Ghaffari q[1] = rho * u[0]; 16477841947SLeila Ghaffari q[2] = rho * u[1]; 16577841947SLeila Ghaffari q[3] = rho * u[2]; 16677841947SLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 16777841947SLeila Ghaffari break; 16832f166c6SLeila Ghaffari case 5: // non-smooth thermal bubble - cylinder 16932f166c6SLeila Ghaffari P = 1.; 17032f166c6SLeila Ghaffari T = 1. - (r < 1. ? S_bubble : 0.); 17132f166c6SLeila Ghaffari rho = P / (R*T); 17232f166c6SLeila Ghaffari u[0] = mean_velocity[0]; 17332f166c6SLeila Ghaffari u[1] = mean_velocity[1]; 17432f166c6SLeila Ghaffari 17532f166c6SLeila Ghaffari // Assign exact solution 17632f166c6SLeila Ghaffari q[0] = rho; 17732f166c6SLeila Ghaffari q[1] = rho * u[0]; 17832f166c6SLeila Ghaffari q[2] = rho * u[1]; 17932f166c6SLeila Ghaffari q[3] = rho * u[2]; 18032f166c6SLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 18132f166c6SLeila Ghaffari break; 18277841947SLeila Ghaffari } 18377841947SLeila Ghaffari // Return 18477841947SLeila Ghaffari return 0; 18577841947SLeila Ghaffari } 18677841947SLeila Ghaffari 18777841947SLeila Ghaffari // ***************************************************************************** 188e6225c47SLeila Ghaffari // Helper function for computing flux Jacobian 189e6225c47SLeila Ghaffari // ***************************************************************************** 190932417b3SJed Brown CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], 191e6225c47SLeila Ghaffari const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 192e6225c47SLeila Ghaffari const CeedScalar gamma) { 193e6225c47SLeila Ghaffari CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 194e6225c47SLeila Ghaffari for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 195e6225c47SLeila Ghaffari for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix 196e6225c47SLeila Ghaffari dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2.)) : 0.) - u[i]*u[j]; 197e6225c47SLeila Ghaffari for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix 198e6225c47SLeila Ghaffari dF[i][0][k+1] = ((i==k) ? 1. : 0.); 199e6225c47SLeila Ghaffari dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) + 200e6225c47SLeila Ghaffari ((i==k) ? u[j] : 0.) - 201e6225c47SLeila Ghaffari ((i==j) ? u[k] : 0.) * (gamma-1.); 202e6225c47SLeila Ghaffari dF[i][4][k+1] = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) - 203e6225c47SLeila Ghaffari (gamma-1.)*u[i]*u[k]; 204e6225c47SLeila Ghaffari } 205e6225c47SLeila Ghaffari dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.); 206e6225c47SLeila Ghaffari } 207e6225c47SLeila Ghaffari dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho); 208e6225c47SLeila Ghaffari dF[i][4][4] = u[i] * gamma; 209e6225c47SLeila Ghaffari } 210e6225c47SLeila Ghaffari } 211e6225c47SLeila Ghaffari 212e6225c47SLeila Ghaffari // ***************************************************************************** 213932417b3SJed Brown // Helper function for computing Tau elements (stabilization constant) 214932417b3SJed Brown // Model from: 215932417b3SJed Brown // Stabilized Methods for Compressible Flows, Hughes et al 2010 216932417b3SJed Brown // 217932417b3SJed Brown // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 218932417b3SJed Brown // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 219932417b3SJed Brown // 220932417b3SJed Brown // Where 221932417b3SJed Brown // c_tau = stabilization constant (0.5 is reported as "optimal") 222932417b3SJed Brown // h[i] = 2 length(dxdX[i]) 223932417b3SJed Brown // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 224932417b3SJed Brown // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 225932417b3SJed Brown // rho(A[i]) = spectral radius of the convective flux Jacobian i, 226932417b3SJed Brown // wave speed in direction i 227932417b3SJed Brown // ***************************************************************************** 228932417b3SJed Brown CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], 229932417b3SJed Brown const CeedScalar dXdx[3][3], const CeedScalar u[3], 230932417b3SJed Brown const CeedScalar sound_speed, const CeedScalar c_tau) { 231ba6664aeSJames Wright for (CeedInt i=0; i<3; i++) { 232932417b3SJed Brown // length of element in direction i 233932417b3SJed Brown CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] + 234932417b3SJed Brown dXdx[2][i]*dXdx[2][i]); 235932417b3SJed Brown // fastest wave in direction i 236932417b3SJed Brown CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 237932417b3SJed Brown Tau_x[i] = c_tau * h / fastest_wave; 238932417b3SJed Brown } 239932417b3SJed Brown } 240932417b3SJed Brown 241932417b3SJed Brown // ***************************************************************************** 24277841947SLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex 24377841947SLeila Ghaffari // ***************************************************************************** 24477841947SLeila Ghaffari CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, 24577841947SLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 24677841947SLeila Ghaffari // Inputs 24777841947SLeila Ghaffari const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 24877841947SLeila Ghaffari 24977841947SLeila Ghaffari // Outputs 25077841947SLeila Ghaffari CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 25177841947SLeila Ghaffari const EulerContext context = (EulerContext)ctx; 25277841947SLeila Ghaffari 25377841947SLeila Ghaffari CeedPragmaSIMD 25477841947SLeila Ghaffari // Quadrature Point Loop 25577841947SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 25677841947SLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 257e6225c47SLeila Ghaffari CeedScalar q[5] = {0.}; 25877841947SLeila Ghaffari 25977841947SLeila Ghaffari Exact_Euler(3, context->curr_time, x, 5, q, ctx); 26077841947SLeila Ghaffari 26177841947SLeila Ghaffari for (CeedInt j=0; j<5; j++) 26277841947SLeila Ghaffari q0[j][i] = q[j]; 26377841947SLeila Ghaffari } // End of Quadrature Point Loop 26477841947SLeila Ghaffari 26577841947SLeila Ghaffari // Return 26677841947SLeila Ghaffari return 0; 26777841947SLeila Ghaffari } 26877841947SLeila Ghaffari 26977841947SLeila Ghaffari // ***************************************************************************** 27077841947SLeila Ghaffari // This QFunction implements the following formulation of Euler equations 27177841947SLeila Ghaffari // with explicit time stepping method 27277841947SLeila Ghaffari // 27377841947SLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation 27477841947SLeila Ghaffari // form with state variables of density, momentum density, and total 27577841947SLeila Ghaffari // energy density. 27677841947SLeila Ghaffari // 27777841947SLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E ) 27877841947SLeila Ghaffari // rho - Mass Density 27977841947SLeila Ghaffari // Ui - Momentum Density, Ui = rho ui 28077841947SLeila Ghaffari // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2 28177841947SLeila Ghaffari // 28277841947SLeila Ghaffari // Euler Equations: 28377841947SLeila Ghaffari // drho/dt + div( U ) = 0 28477841947SLeila Ghaffari // dU/dt + div( rho (u x u) + P I3 ) = 0 28577841947SLeila Ghaffari // dE/dt + div( (E + P) u ) = 0 28677841947SLeila Ghaffari // 28777841947SLeila Ghaffari // Equation of State: 28877841947SLeila Ghaffari // P = (gamma - 1) (E - rho (u u) / 2) 28977841947SLeila Ghaffari // 29077841947SLeila Ghaffari // Constants: 29177841947SLeila Ghaffari // cv , Specific heat, constant volume 29277841947SLeila Ghaffari // cp , Specific heat, constant pressure 29377841947SLeila Ghaffari // g , Gravity 29477841947SLeila Ghaffari // gamma = cp / cv, Specific heat ratio 29577841947SLeila Ghaffari // ***************************************************************************** 29677841947SLeila Ghaffari CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, 29777841947SLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 29877841947SLeila Ghaffari // *INDENT-OFF* 29977841947SLeila Ghaffari // Inputs 30077841947SLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 301e6225c47SLeila Ghaffari (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 30277841947SLeila Ghaffari (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 30377841947SLeila Ghaffari // Outputs 30477841947SLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 30577841947SLeila Ghaffari (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 30677841947SLeila Ghaffari 307e6225c47SLeila Ghaffari EulerContext context = (EulerContext)ctx; 308932417b3SJed Brown const CeedScalar c_tau = context->c_tau; 30977841947SLeila Ghaffari const CeedScalar gamma = 1.4; 31077841947SLeila Ghaffari 31177841947SLeila Ghaffari CeedPragmaSIMD 31277841947SLeila Ghaffari // Quadrature Point Loop 31377841947SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 31477841947SLeila Ghaffari // *INDENT-OFF* 31577841947SLeila Ghaffari // Setup 31677841947SLeila Ghaffari // -- Interp in 31777841947SLeila Ghaffari const CeedScalar rho = q[0][i]; 31877841947SLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 31977841947SLeila Ghaffari q[2][i] / rho, 32077841947SLeila Ghaffari q[3][i] / rho 32177841947SLeila Ghaffari }; 32277841947SLeila Ghaffari const CeedScalar E = q[4][i]; 323e6225c47SLeila Ghaffari const CeedScalar drho[3] = {dq[0][0][i], 324e6225c47SLeila Ghaffari dq[1][0][i], 325e6225c47SLeila Ghaffari dq[2][0][i] 326e6225c47SLeila Ghaffari }; 327e6225c47SLeila Ghaffari const CeedScalar dU[3][3] = {{dq[0][1][i], 328e6225c47SLeila Ghaffari dq[1][1][i], 329e6225c47SLeila Ghaffari dq[2][1][i]}, 330e6225c47SLeila Ghaffari {dq[0][2][i], 331e6225c47SLeila Ghaffari dq[1][2][i], 332e6225c47SLeila Ghaffari dq[2][2][i]}, 333e6225c47SLeila Ghaffari {dq[0][3][i], 334e6225c47SLeila Ghaffari dq[1][3][i], 335e6225c47SLeila Ghaffari dq[2][3][i]} 336e6225c47SLeila Ghaffari }; 337e6225c47SLeila Ghaffari const CeedScalar dE[3] = {dq[0][4][i], 338e6225c47SLeila Ghaffari dq[1][4][i], 339e6225c47SLeila Ghaffari dq[2][4][i] 340e6225c47SLeila Ghaffari }; 34177841947SLeila Ghaffari // -- Interp-to-Interp q_data 34277841947SLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 34377841947SLeila Ghaffari // -- Interp-to-Grad q_data 34477841947SLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 34577841947SLeila Ghaffari // *INDENT-OFF* 34677841947SLeila Ghaffari const CeedScalar dXdx[3][3] = {{q_data[1][i], 34777841947SLeila Ghaffari q_data[2][i], 34877841947SLeila Ghaffari q_data[3][i]}, 34977841947SLeila Ghaffari {q_data[4][i], 35077841947SLeila Ghaffari q_data[5][i], 35177841947SLeila Ghaffari q_data[6][i]}, 35277841947SLeila Ghaffari {q_data[7][i], 35377841947SLeila Ghaffari q_data[8][i], 35477841947SLeila Ghaffari q_data[9][i]} 35577841947SLeila Ghaffari }; 35677841947SLeila Ghaffari // *INDENT-ON* 357e6225c47SLeila Ghaffari // dU/dx 358e6225c47SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 359e6225c47SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 360e6225c47SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 361e6225c47SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 362ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) { 363ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) { 364e6225c47SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 365e6225c47SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 366ba6664aeSJames Wright for (CeedInt l=0; l<3; l++) { 367e6225c47SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 368e6225c47SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 369e6225c47SLeila Ghaffari } 370e6225c47SLeila Ghaffari } 371e6225c47SLeila Ghaffari } 372e6225c47SLeila Ghaffari // Pressure 37377841947SLeila Ghaffari const CeedScalar 37477841947SLeila Ghaffari E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 37577841947SLeila Ghaffari E_internal = E - E_kinetic, 376e6225c47SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 37777841947SLeila Ghaffari 37877841947SLeila Ghaffari // The Physics 37977841947SLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 380ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) { 381e6225c47SLeila Ghaffari v[j][i] = 0.; 382ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 383e6225c47SLeila Ghaffari dv[k][j][i] = 0.; 38477841947SLeila Ghaffari } 38577841947SLeila Ghaffari 38677841947SLeila Ghaffari // -- Density 38777841947SLeila Ghaffari // ---- u rho 388ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 38977841947SLeila Ghaffari dv[j][0][i] += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 39077841947SLeila Ghaffari rho*u[2]*dXdx[j][2]); 39177841947SLeila Ghaffari // -- Momentum 39277841947SLeila Ghaffari // ---- rho (u x u) + P I3 393ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 394ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 395e6225c47SLeila Ghaffari dv[k][j+1][i] += wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] + 396e6225c47SLeila Ghaffari (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] + 397e6225c47SLeila Ghaffari (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]); 39877841947SLeila Ghaffari // -- Total Energy Density 39977841947SLeila Ghaffari // ---- (E + P) u 400ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 40177841947SLeila Ghaffari dv[j][4][i] += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 40277841947SLeila Ghaffari u[2]*dXdx[j][2]); 403e6225c47SLeila Ghaffari 404e6225c47SLeila Ghaffari // --Stabilization terms 405e6225c47SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 406e6225c47SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 407932417b3SJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 408e6225c47SLeila Ghaffari 409e6225c47SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 410e6225c47SLeila Ghaffari CeedScalar dqdx[5][3]; 411ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) { 412e6225c47SLeila Ghaffari dqdx[0][j] = drhodx[j]; 413e6225c47SLeila Ghaffari dqdx[4][j] = dEdx[j]; 414ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 415e6225c47SLeila Ghaffari dqdx[k+1][j] = dUdx[k][j]; 416e6225c47SLeila Ghaffari } 417e6225c47SLeila Ghaffari 418e6225c47SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 419e6225c47SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 420ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 421ba6664aeSJames Wright for (CeedInt k=0; k<5; k++) 422ba6664aeSJames Wright for (CeedInt l=0; l<5; l++) 423e6225c47SLeila Ghaffari strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 424e6225c47SLeila Ghaffari 425932417b3SJed Brown // Stabilization 426932417b3SJed Brown // -- Tau elements 427932417b3SJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 428932417b3SJed Brown CeedScalar Tau_x[3] = {0.}; 429932417b3SJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 430e6225c47SLeila Ghaffari 431932417b3SJed Brown // -- Stabilization method: none or SU 43288626eedSJames Wright CeedScalar stab[5][3] = {{0.}}; 433e6225c47SLeila Ghaffari switch (context->stabilization) { 434e6225c47SLeila Ghaffari case 0: // Galerkin 435e6225c47SLeila Ghaffari break; 436e6225c47SLeila Ghaffari case 1: // SU 437ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 438ba6664aeSJames Wright for (CeedInt k=0; k<5; k++) 439ba6664aeSJames Wright for (CeedInt l=0; l<5; l++) 44088626eedSJames Wright stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 441e6225c47SLeila Ghaffari 442ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) 443ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 444e6225c47SLeila Ghaffari dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] + 445e6225c47SLeila Ghaffari stab[j][1] * dXdx[k][1] + 446e6225c47SLeila Ghaffari stab[j][2] * dXdx[k][2]); 447e6225c47SLeila Ghaffari break; 448e6225c47SLeila Ghaffari case 2: // SUPG is not implemented for explicit scheme 449e6225c47SLeila Ghaffari break; 450e6225c47SLeila Ghaffari } 451e6225c47SLeila Ghaffari 45277841947SLeila Ghaffari } // End Quadrature Point Loop 45377841947SLeila Ghaffari 45477841947SLeila Ghaffari // Return 45577841947SLeila Ghaffari return 0; 45677841947SLeila Ghaffari } 45777841947SLeila Ghaffari 45877841947SLeila Ghaffari // ***************************************************************************** 45977841947SLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above) 46077841947SLeila Ghaffari // with implicit time stepping method 46177841947SLeila Ghaffari // 46277841947SLeila Ghaffari // ***************************************************************************** 46377841947SLeila Ghaffari CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, 46477841947SLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 46577841947SLeila Ghaffari // *INDENT-OFF* 46677841947SLeila Ghaffari // Inputs 46777841947SLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 468e6225c47SLeila Ghaffari (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 46977841947SLeila Ghaffari (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 47077841947SLeila Ghaffari (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 47177841947SLeila Ghaffari // Outputs 47277841947SLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 47377841947SLeila Ghaffari (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 47477841947SLeila Ghaffari 475e6225c47SLeila Ghaffari EulerContext context = (EulerContext)ctx; 476932417b3SJed Brown const CeedScalar c_tau = context->c_tau; 47777841947SLeila Ghaffari const CeedScalar gamma = 1.4; 47877841947SLeila Ghaffari 47977841947SLeila Ghaffari CeedPragmaSIMD 48077841947SLeila Ghaffari // Quadrature Point Loop 48177841947SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 48277841947SLeila Ghaffari // *INDENT-OFF* 48377841947SLeila Ghaffari // Setup 48477841947SLeila Ghaffari // -- Interp in 48577841947SLeila Ghaffari const CeedScalar rho = q[0][i]; 48677841947SLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 48777841947SLeila Ghaffari q[2][i] / rho, 48877841947SLeila Ghaffari q[3][i] / rho 48977841947SLeila Ghaffari }; 49077841947SLeila Ghaffari const CeedScalar E = q[4][i]; 491e6225c47SLeila Ghaffari const CeedScalar drho[3] = {dq[0][0][i], 492e6225c47SLeila Ghaffari dq[1][0][i], 493e6225c47SLeila Ghaffari dq[2][0][i] 494e6225c47SLeila Ghaffari }; 495e6225c47SLeila Ghaffari const CeedScalar dU[3][3] = {{dq[0][1][i], 496e6225c47SLeila Ghaffari dq[1][1][i], 497e6225c47SLeila Ghaffari dq[2][1][i]}, 498e6225c47SLeila Ghaffari {dq[0][2][i], 499e6225c47SLeila Ghaffari dq[1][2][i], 500e6225c47SLeila Ghaffari dq[2][2][i]}, 501e6225c47SLeila Ghaffari {dq[0][3][i], 502e6225c47SLeila Ghaffari dq[1][3][i], 503e6225c47SLeila Ghaffari dq[2][3][i]} 504e6225c47SLeila Ghaffari }; 505e6225c47SLeila Ghaffari const CeedScalar dE[3] = {dq[0][4][i], 506e6225c47SLeila Ghaffari dq[1][4][i], 507e6225c47SLeila Ghaffari dq[2][4][i] 508e6225c47SLeila Ghaffari }; 50977841947SLeila Ghaffari // -- Interp-to-Interp q_data 51077841947SLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 51177841947SLeila Ghaffari // -- Interp-to-Grad q_data 51277841947SLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 51377841947SLeila Ghaffari // *INDENT-OFF* 51477841947SLeila Ghaffari const CeedScalar dXdx[3][3] = {{q_data[1][i], 51577841947SLeila Ghaffari q_data[2][i], 51677841947SLeila Ghaffari q_data[3][i]}, 51777841947SLeila Ghaffari {q_data[4][i], 51877841947SLeila Ghaffari q_data[5][i], 51977841947SLeila Ghaffari q_data[6][i]}, 52077841947SLeila Ghaffari {q_data[7][i], 52177841947SLeila Ghaffari q_data[8][i], 52277841947SLeila Ghaffari q_data[9][i]} 52377841947SLeila Ghaffari }; 52477841947SLeila Ghaffari // *INDENT-ON* 525e6225c47SLeila Ghaffari // dU/dx 526e6225c47SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 527e6225c47SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 528e6225c47SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 529e6225c47SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 530ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) { 531ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) { 532e6225c47SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 533e6225c47SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 534ba6664aeSJames Wright for (CeedInt l=0; l<3; l++) { 535e6225c47SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 536e6225c47SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 537e6225c47SLeila Ghaffari } 538e6225c47SLeila Ghaffari } 539e6225c47SLeila Ghaffari } 54077841947SLeila Ghaffari const CeedScalar 54177841947SLeila Ghaffari E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 54277841947SLeila Ghaffari E_internal = E - E_kinetic, 543e6225c47SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 54477841947SLeila Ghaffari 54577841947SLeila Ghaffari // The Physics 54677841947SLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 547ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) { 548e6225c47SLeila Ghaffari v[j][i] = 0.; 549ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 550e6225c47SLeila Ghaffari dv[k][j][i] = 0.; 55177841947SLeila Ghaffari } 55277841947SLeila Ghaffari //-----mass matrix 553ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) 55477841947SLeila Ghaffari v[j][i] += wdetJ*q_dot[j][i]; 55577841947SLeila Ghaffari 55677841947SLeila Ghaffari // -- Density 55777841947SLeila Ghaffari // ---- u rho 558ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 55977841947SLeila Ghaffari dv[j][0][i] -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 56077841947SLeila Ghaffari rho*u[2]*dXdx[j][2]); 56177841947SLeila Ghaffari // -- Momentum 56277841947SLeila Ghaffari // ---- rho (u x u) + P I3 563ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 564ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 565e6225c47SLeila Ghaffari dv[k][j+1][i] -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] + 566e6225c47SLeila Ghaffari (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] + 567e6225c47SLeila Ghaffari (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]); 56877841947SLeila Ghaffari // -- Total Energy Density 56977841947SLeila Ghaffari // ---- (E + P) u 570ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 57177841947SLeila Ghaffari dv[j][4][i] -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 57277841947SLeila Ghaffari u[2]*dXdx[j][2]); 573e6225c47SLeila Ghaffari 574e6225c47SLeila Ghaffari // -- Stabilization terms 575e6225c47SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 576e6225c47SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 577932417b3SJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 578e6225c47SLeila Ghaffari 579e6225c47SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 580e6225c47SLeila Ghaffari CeedScalar dqdx[5][3]; 581ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) { 582e6225c47SLeila Ghaffari dqdx[0][j] = drhodx[j]; 583e6225c47SLeila Ghaffari dqdx[4][j] = dEdx[j]; 584ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 585e6225c47SLeila Ghaffari dqdx[k+1][j] = dUdx[k][j]; 586e6225c47SLeila Ghaffari } 587e6225c47SLeila Ghaffari 588e6225c47SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 589e6225c47SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 590ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 591ba6664aeSJames Wright for (CeedInt k=0; k<5; k++) 592ba6664aeSJames Wright for (CeedInt l=0; l<5; l++) 593e6225c47SLeila Ghaffari strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 594e6225c47SLeila Ghaffari 595e6225c47SLeila Ghaffari // ---- Strong residual 596e6225c47SLeila Ghaffari CeedScalar strong_res[5]; 597ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) 598e6225c47SLeila Ghaffari strong_res[j] = q_dot[j][i] + strong_conv[j]; 599e6225c47SLeila Ghaffari 600932417b3SJed Brown // Stabilization 601932417b3SJed Brown // -- Tau elements 602932417b3SJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 603932417b3SJed Brown CeedScalar Tau_x[3] = {0.}; 604932417b3SJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 605e6225c47SLeila Ghaffari 606932417b3SJed Brown // -- Stabilization method: none, SU, or SUPG 60788626eedSJames Wright CeedScalar stab[5][3] = {{0.}}; 608e6225c47SLeila Ghaffari switch (context->stabilization) { 609e6225c47SLeila Ghaffari case 0: // Galerkin 610e6225c47SLeila Ghaffari break; 611e6225c47SLeila Ghaffari case 1: // SU 612ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 613ba6664aeSJames Wright for (CeedInt k=0; k<5; k++) 614ba6664aeSJames Wright for (CeedInt l=0; l<5; l++) 61588626eedSJames Wright stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l]; 616e6225c47SLeila Ghaffari 617ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) 618ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 619e6225c47SLeila Ghaffari dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 620e6225c47SLeila Ghaffari stab[j][1] * dXdx[k][1] + 621e6225c47SLeila Ghaffari stab[j][2] * dXdx[k][2]); 622e6225c47SLeila Ghaffari break; 623e6225c47SLeila Ghaffari case 2: // SUPG 624ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 625ba6664aeSJames Wright for (CeedInt k=0; k<5; k++) 626ba6664aeSJames Wright for (CeedInt l=0; l<5; l++) 62788626eedSJames Wright stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l]; 628e6225c47SLeila Ghaffari 629ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) 630ba6664aeSJames Wright for (CeedInt k=0; k<3; k++) 631e6225c47SLeila Ghaffari dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 632e6225c47SLeila Ghaffari stab[j][1] * dXdx[k][1] + 633e6225c47SLeila Ghaffari stab[j][2] * dXdx[k][2]); 634e6225c47SLeila Ghaffari break; 635e6225c47SLeila Ghaffari } 63677841947SLeila Ghaffari } // End Quadrature Point Loop 63777841947SLeila Ghaffari 63877841947SLeila Ghaffari // Return 63977841947SLeila Ghaffari return 0; 64077841947SLeila Ghaffari } 64177841947SLeila Ghaffari // ***************************************************************************** 6422fe7aee7SLeila Ghaffari // This QFunction sets the inflow boundary conditions for 6432fe7aee7SLeila Ghaffari // the traveling vortex problem. 64477841947SLeila Ghaffari // 64577841947SLeila Ghaffari // Prescribed T_inlet and P_inlet are converted to conservative variables 64677841947SLeila Ghaffari // and applied weakly. 64777841947SLeila Ghaffari // 64877841947SLeila Ghaffari // ***************************************************************************** 6492fe7aee7SLeila Ghaffari CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, 65077841947SLeila Ghaffari const CeedScalar *const *in, 65177841947SLeila Ghaffari CeedScalar *const *out) { 65277841947SLeila Ghaffari // *INDENT-OFF* 65377841947SLeila Ghaffari // Inputs 654e8b03feeSJames Wright const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 65577841947SLeila Ghaffari // Outputs 65677841947SLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 65777841947SLeila Ghaffari // *INDENT-ON* 65877841947SLeila Ghaffari EulerContext context = (EulerContext)ctx; 65977841947SLeila Ghaffari const int euler_test = context->euler_test; 66077841947SLeila Ghaffari const bool implicit = context->implicit; 66177841947SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 66277841947SLeila Ghaffari const CeedScalar cv = 2.5; 66377841947SLeila Ghaffari const CeedScalar R = 1.; 66477841947SLeila Ghaffari CeedScalar T_inlet; 66577841947SLeila Ghaffari CeedScalar P_inlet; 66677841947SLeila Ghaffari 66777841947SLeila Ghaffari // For test cases 1 and 3 the background velocity is zero 66877841947SLeila Ghaffari if (euler_test == 1 || euler_test == 3) 66977841947SLeila Ghaffari for (CeedInt i=0; i<3; i++) mean_velocity[i] = 0.; 67077841947SLeila Ghaffari 67177841947SLeila Ghaffari // For test cases 1 and 2, T_inlet = T_inlet = 0.4 67277841947SLeila Ghaffari if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4; 67377841947SLeila Ghaffari else T_inlet = P_inlet = 1.; 67477841947SLeila Ghaffari 67577841947SLeila Ghaffari CeedPragmaSIMD 67677841947SLeila Ghaffari // Quadrature Point Loop 67777841947SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 67877841947SLeila Ghaffari // Setup 67977841947SLeila Ghaffari // -- Interp-to-Interp q_data 68077841947SLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 68177841947SLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 68277841947SLeila Ghaffari // We can effect this by swapping the sign on this weight 68377841947SLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 6842fe7aee7SLeila Ghaffari // ---- Normal vect 68577841947SLeila Ghaffari const CeedScalar norm[3] = {q_data_sur[1][i], 68677841947SLeila Ghaffari q_data_sur[2][i], 68777841947SLeila Ghaffari q_data_sur[3][i] 68877841947SLeila Ghaffari }; 68977841947SLeila Ghaffari 69077841947SLeila Ghaffari // face_normal = Normal vector of the face 69177841947SLeila Ghaffari const CeedScalar face_normal = norm[0]*mean_velocity[0] + 69277841947SLeila Ghaffari norm[1]*mean_velocity[1] + 69377841947SLeila Ghaffari norm[2]*mean_velocity[2]; 69477841947SLeila Ghaffari // The Physics 69577841947SLeila Ghaffari // Zero v so all future terms can safely sum into it 696ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) v[j][i] = 0.; 69777841947SLeila Ghaffari 69877841947SLeila Ghaffari // Implementing in/outflow BCs 6992fe7aee7SLeila Ghaffari if (face_normal > 0) { 70077841947SLeila Ghaffari } else { // inflow 70177841947SLeila Ghaffari const CeedScalar rho_inlet = P_inlet/(R*T_inlet); 70277841947SLeila Ghaffari const CeedScalar E_kinetic_inlet = (mean_velocity[0]*mean_velocity[0] + 70377841947SLeila Ghaffari mean_velocity[1]*mean_velocity[1]) / 2.; 70477841947SLeila Ghaffari // incoming total energy 70577841947SLeila Ghaffari const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet); 70677841947SLeila Ghaffari 70777841947SLeila Ghaffari // The Physics 70877841947SLeila Ghaffari // -- Density 70977841947SLeila Ghaffari v[0][i] -= wdetJb * rho_inlet * face_normal; 71077841947SLeila Ghaffari 71177841947SLeila Ghaffari // -- Momentum 712ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 71377841947SLeila Ghaffari v[j+1][i] -= wdetJb *(rho_inlet * face_normal * mean_velocity[j] + 71477841947SLeila Ghaffari norm[j] * P_inlet); 71577841947SLeila Ghaffari 71677841947SLeila Ghaffari // -- Total Energy Density 71777841947SLeila Ghaffari v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet); 71877841947SLeila Ghaffari } 71977841947SLeila Ghaffari 72077841947SLeila Ghaffari } // End Quadrature Point Loop 72177841947SLeila Ghaffari return 0; 72277841947SLeila Ghaffari } 72377841947SLeila Ghaffari 72477841947SLeila Ghaffari // ***************************************************************************** 72555e76554SLeila Ghaffari // This QFunction sets the outflow boundary conditions for 72655e76554SLeila Ghaffari // the Euler solver. 72755e76554SLeila Ghaffari // 72855e76554SLeila Ghaffari // Outflow BCs: 72955e76554SLeila Ghaffari // The validity of the weak form of the governing equations is 73055e76554SLeila Ghaffari // extended to the outflow. 73155e76554SLeila Ghaffari // 73255e76554SLeila Ghaffari // ***************************************************************************** 73355e76554SLeila Ghaffari CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, 73455e76554SLeila Ghaffari const CeedScalar *const *in, 73555e76554SLeila Ghaffari CeedScalar *const *out) { 73655e76554SLeila Ghaffari // *INDENT-OFF* 73755e76554SLeila Ghaffari // Inputs 73855e76554SLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 739e8b03feeSJames Wright (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 74055e76554SLeila Ghaffari // Outputs 74155e76554SLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 74255e76554SLeila Ghaffari // *INDENT-ON* 74355e76554SLeila Ghaffari EulerContext context = (EulerContext)ctx; 74455e76554SLeila Ghaffari const bool implicit = context->implicit; 74555e76554SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 74655e76554SLeila Ghaffari 74755e76554SLeila Ghaffari const CeedScalar gamma = 1.4; 74855e76554SLeila Ghaffari 74955e76554SLeila Ghaffari CeedPragmaSIMD 75055e76554SLeila Ghaffari // Quadrature Point Loop 75155e76554SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 75255e76554SLeila Ghaffari // Setup 75355e76554SLeila Ghaffari // -- Interp in 75455e76554SLeila Ghaffari const CeedScalar rho = q[0][i]; 75555e76554SLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 75655e76554SLeila Ghaffari q[2][i] / rho, 75755e76554SLeila Ghaffari q[3][i] / rho 75855e76554SLeila Ghaffari }; 75955e76554SLeila Ghaffari const CeedScalar E = q[4][i]; 76055e76554SLeila Ghaffari 76155e76554SLeila Ghaffari // -- Interp-to-Interp q_data 76255e76554SLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 76355e76554SLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 76455e76554SLeila Ghaffari // We can effect this by swapping the sign on this weight 76555e76554SLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 76655e76554SLeila Ghaffari // ---- Normal vectors 76755e76554SLeila Ghaffari const CeedScalar norm[3] = {q_data_sur[1][i], 76855e76554SLeila Ghaffari q_data_sur[2][i], 76955e76554SLeila Ghaffari q_data_sur[3][i] 77055e76554SLeila Ghaffari }; 77155e76554SLeila Ghaffari 77255e76554SLeila Ghaffari // face_normal = Normal vector of the face 77355e76554SLeila Ghaffari const CeedScalar face_normal = norm[0]*mean_velocity[0] + 77455e76554SLeila Ghaffari norm[1]*mean_velocity[1] + 77555e76554SLeila Ghaffari norm[2]*mean_velocity[2]; 77655e76554SLeila Ghaffari // The Physics 77755e76554SLeila Ghaffari // Zero v so all future terms can safely sum into it 778ba6664aeSJames Wright for (CeedInt j=0; j<5; j++) v[j][i] = 0; 77955e76554SLeila Ghaffari 78055e76554SLeila Ghaffari // Implementing in/outflow BCs 78155e76554SLeila Ghaffari if (face_normal > 0) { // outflow 78255e76554SLeila Ghaffari const CeedScalar E_kinetic = (u[0]*u[0] + u[1]*u[1]) / 2.; 78355e76554SLeila Ghaffari const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure 78455e76554SLeila Ghaffari const CeedScalar u_normal = norm[0]*u[0] + norm[1]*u[1] + 78555e76554SLeila Ghaffari norm[2]*u[2]; // Normal velocity 78655e76554SLeila Ghaffari // The Physics 78755e76554SLeila Ghaffari // -- Density 78855e76554SLeila Ghaffari v[0][i] -= wdetJb * rho * u_normal; 78955e76554SLeila Ghaffari 79055e76554SLeila Ghaffari // -- Momentum 791ba6664aeSJames Wright for (CeedInt j=0; j<3; j++) 79255e76554SLeila Ghaffari v[j+1][i] -= wdetJb *(rho * u_normal * u[j] + norm[j] * P); 79355e76554SLeila Ghaffari 79455e76554SLeila Ghaffari // -- Total Energy Density 79555e76554SLeila Ghaffari v[4][i] -= wdetJb * u_normal * (E + P); 79655e76554SLeila Ghaffari } 79755e76554SLeila Ghaffari } // End Quadrature Point Loop 79855e76554SLeila Ghaffari return 0; 79955e76554SLeila Ghaffari } 80055e76554SLeila Ghaffari 80155e76554SLeila Ghaffari // ***************************************************************************** 80277841947SLeila Ghaffari 80377841947SLeila Ghaffari #endif // eulervortex_h 804