1a515125bSLeila Ghaffari // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2a515125bSLeila Ghaffari // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3a515125bSLeila Ghaffari // reserved. See files LICENSE and NOTICE for details. 4a515125bSLeila Ghaffari // 5a515125bSLeila Ghaffari // This file is part of CEED, a collection of benchmarks, miniapps, software 6a515125bSLeila Ghaffari // libraries and APIs for efficient high-order finite element and spectral 7a515125bSLeila Ghaffari // element discretizations for exascale applications. For more information and 8a515125bSLeila Ghaffari // source code availability see http://github.com/ceed. 9a515125bSLeila Ghaffari // 10a515125bSLeila Ghaffari // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11a515125bSLeila Ghaffari // a collaborative effort of two U.S. Department of Energy organizations (Office 12a515125bSLeila Ghaffari // of Science and the National Nuclear Security Administration) responsible for 13a515125bSLeila Ghaffari // the planning and preparation of a capable exascale ecosystem, including 14a515125bSLeila Ghaffari // software, applications, hardware, advanced system engineering and early 15a515125bSLeila Ghaffari // testbed platforms, in support of the nation's exascale computing imperative. 16a515125bSLeila Ghaffari 17a515125bSLeila Ghaffari /// @file 18a515125bSLeila Ghaffari /// Euler traveling vortex initial condition and operator for Navier-Stokes 19a515125bSLeila Ghaffari /// example using PETSc 20a515125bSLeila Ghaffari 21a515125bSLeila Ghaffari // Model from: 22a515125bSLeila Ghaffari // On the Order of Accuracy and Numerical Performance of Two Classes of 23a515125bSLeila Ghaffari // Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011). 24a515125bSLeila Ghaffari 25a515125bSLeila Ghaffari #ifndef eulervortex_h 26a515125bSLeila Ghaffari #define eulervortex_h 27a515125bSLeila Ghaffari 28a515125bSLeila Ghaffari #include <math.h> 29*3a8779fbSJames Wright #include <ceed.h> 30a515125bSLeila Ghaffari 31a515125bSLeila Ghaffari #ifndef M_PI 32a515125bSLeila Ghaffari #define M_PI 3.14159265358979323846 33a515125bSLeila Ghaffari #endif 34a515125bSLeila Ghaffari 35a515125bSLeila Ghaffari #ifndef euler_context_struct 36a515125bSLeila Ghaffari #define euler_context_struct 37a515125bSLeila Ghaffari typedef struct EulerContext_ *EulerContext; 38a515125bSLeila Ghaffari struct EulerContext_ { 39a515125bSLeila Ghaffari CeedScalar center[3]; 40a515125bSLeila Ghaffari CeedScalar curr_time; 41a515125bSLeila Ghaffari CeedScalar vortex_strength; 42d8a22b9eSJed Brown CeedScalar c_tau; 43a515125bSLeila Ghaffari CeedScalar mean_velocity[3]; 44a515125bSLeila Ghaffari bool implicit; 45139613f2SLeila Ghaffari int euler_test; 46139613f2SLeila Ghaffari int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG 47a515125bSLeila Ghaffari }; 48a515125bSLeila Ghaffari #endif 49a515125bSLeila Ghaffari 50a515125bSLeila Ghaffari // ***************************************************************************** 51a515125bSLeila Ghaffari // This function sets the initial conditions 52a515125bSLeila Ghaffari // 53a515125bSLeila Ghaffari // Temperature: 54a515125bSLeila Ghaffari // T = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2) 55a515125bSLeila Ghaffari // Density: 56a515125bSLeila Ghaffari // rho = (T/S_vortex)^(1 / (gamma - 1)) 57a515125bSLeila Ghaffari // Pressure: 58a515125bSLeila Ghaffari // P = rho * T 59a515125bSLeila Ghaffari // Velocity: 60a515125bSLeila Ghaffari // ui = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi) 61a515125bSLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 ) 62a515125bSLeila Ghaffari // Velocity/Momentum Density: 63a515125bSLeila Ghaffari // Ui = rho ui 64a515125bSLeila Ghaffari // Total Energy: 65a515125bSLeila Ghaffari // E = P / (gamma - 1) + rho (u u)/2 66a515125bSLeila Ghaffari // 67a515125bSLeila Ghaffari // Constants: 68a515125bSLeila Ghaffari // cv , Specific heat, constant volume 69a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 70a515125bSLeila Ghaffari // vortex_strength , Strength of vortex 71a515125bSLeila Ghaffari // center , Location of bubble center 72a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio 73a515125bSLeila Ghaffari // 74a515125bSLeila Ghaffari // ***************************************************************************** 75a515125bSLeila Ghaffari 76a515125bSLeila Ghaffari // ***************************************************************************** 77a515125bSLeila Ghaffari // This helper function provides support for the exact, time-dependent solution 78a515125bSLeila Ghaffari // (currently not implemented) and IC formulation for Euler traveling vortex 79a515125bSLeila Ghaffari // ***************************************************************************** 80a515125bSLeila Ghaffari CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, 81a515125bSLeila Ghaffari const CeedScalar X[], CeedInt Nf, CeedScalar q[], 82a515125bSLeila Ghaffari void *ctx) { 83a515125bSLeila Ghaffari // Context 84a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx; 85a515125bSLeila Ghaffari const CeedScalar vortex_strength = context->vortex_strength; 86a515125bSLeila Ghaffari const CeedScalar *center = context->center; // Center of the domain 87a515125bSLeila Ghaffari const CeedScalar *mean_velocity = context->mean_velocity; 88a515125bSLeila Ghaffari 89a515125bSLeila Ghaffari // Setup 90a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 91a515125bSLeila Ghaffari const CeedScalar cv = 2.5; 92a515125bSLeila Ghaffari const CeedScalar R = 1.; 93a515125bSLeila Ghaffari const CeedScalar x = X[0], y = X[1]; // Coordinates 94a515125bSLeila Ghaffari // Vortex center 95a515125bSLeila Ghaffari const CeedScalar xc = center[0] + mean_velocity[0] * time; 96a515125bSLeila Ghaffari const CeedScalar yc = center[1] + mean_velocity[1] * time; 97a515125bSLeila Ghaffari 98a515125bSLeila Ghaffari const CeedScalar x0 = x - xc; 99a515125bSLeila Ghaffari const CeedScalar y0 = y - yc; 100a515125bSLeila Ghaffari const CeedScalar r = sqrt( x0*x0 + y0*y0 ); 101a515125bSLeila Ghaffari const CeedScalar C = vortex_strength * exp((1. - r*r)/2.) / (2. * M_PI); 102139613f2SLeila Ghaffari const CeedScalar delta_T = - (gamma - 1.) * vortex_strength * vortex_strength * 103139613f2SLeila Ghaffari exp(1 - r*r) / (8. * gamma * M_PI * M_PI); 104a515125bSLeila Ghaffari const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma 105a515125bSLeila Ghaffari const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / 106a515125bSLeila Ghaffari (8.*gamma*M_PI*M_PI); 107a515125bSLeila Ghaffari CeedScalar rho, P, T, E, u[3] = {0.}; 108a515125bSLeila Ghaffari 109a515125bSLeila Ghaffari // Initial Conditions 110a515125bSLeila Ghaffari switch (context->euler_test) { 111a515125bSLeila Ghaffari case 0: // Traveling vortex 112a515125bSLeila Ghaffari T = 1 + delta_T; 113a515125bSLeila Ghaffari // P = rho * T 114a515125bSLeila Ghaffari // P = S * rho^gamma 115a515125bSLeila Ghaffari // Solve for rho, then substitute for P 116139613f2SLeila Ghaffari rho = pow(T/S_vortex, 1 / (gamma - 1.)); 117a515125bSLeila Ghaffari P = rho * T; 118a515125bSLeila Ghaffari u[0] = mean_velocity[0] - C*y0; 119a515125bSLeila Ghaffari u[1] = mean_velocity[1] + C*x0; 120a515125bSLeila Ghaffari 121a515125bSLeila Ghaffari // Assign exact solution 122a515125bSLeila Ghaffari q[0] = rho; 123a515125bSLeila Ghaffari q[1] = rho * u[0]; 124a515125bSLeila Ghaffari q[2] = rho * u[1]; 125a515125bSLeila Ghaffari q[3] = rho * u[2]; 126a515125bSLeila Ghaffari q[4] = P / (gamma - 1.) + rho * (u[0]*u[0] + u[1]*u[1]) / 2.; 127a515125bSLeila Ghaffari break; 128a515125bSLeila Ghaffari case 1: // Constant zero velocity, density constant, total energy constant 129a515125bSLeila Ghaffari rho = 1.; 130a515125bSLeila Ghaffari E = 2.; 131a515125bSLeila Ghaffari 132a515125bSLeila Ghaffari // Assign exact solution 133a515125bSLeila Ghaffari q[0] = rho; 134a515125bSLeila Ghaffari q[1] = rho * u[0]; 135a515125bSLeila Ghaffari q[2] = rho * u[1]; 136a515125bSLeila Ghaffari q[3] = rho * u[2]; 137a515125bSLeila Ghaffari q[4] = E; 138a515125bSLeila Ghaffari break; 139a515125bSLeila Ghaffari case 2: // Constant nonzero velocity, density constant, total energy constant 140a515125bSLeila Ghaffari rho = 1.; 141a515125bSLeila Ghaffari E = 2.; 142a515125bSLeila Ghaffari u[0] = mean_velocity[0]; 143a515125bSLeila Ghaffari u[1] = mean_velocity[1]; 144a515125bSLeila Ghaffari 145a515125bSLeila Ghaffari // Assign exact solution 146a515125bSLeila Ghaffari q[0] = rho; 147a515125bSLeila Ghaffari q[1] = rho * u[0]; 148a515125bSLeila Ghaffari q[2] = rho * u[1]; 149a515125bSLeila Ghaffari q[3] = rho * u[2]; 150a515125bSLeila Ghaffari q[4] = E; 151a515125bSLeila Ghaffari break; 152a515125bSLeila Ghaffari case 3: // Velocity zero, pressure constant 153a515125bSLeila Ghaffari // (so density and internal energy will be non-constant), 154a515125bSLeila Ghaffari // but the velocity should stay zero and the bubble won't diffuse 155a515125bSLeila Ghaffari // (for Euler, where there is no thermal conductivity) 156a515125bSLeila Ghaffari P = 1.; 157a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r*r); 158a515125bSLeila Ghaffari rho = P / (R*T); 159a515125bSLeila Ghaffari 160a515125bSLeila Ghaffari // Assign exact solution 161a515125bSLeila Ghaffari q[0] = rho; 162a515125bSLeila Ghaffari q[1] = rho * u[0]; 163a515125bSLeila Ghaffari q[2] = rho * u[1]; 164a515125bSLeila Ghaffari q[3] = rho * u[2]; 165a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 166a515125bSLeila Ghaffari break; 167a515125bSLeila Ghaffari case 4: // Constant nonzero velocity, pressure constant 168a515125bSLeila Ghaffari // (so density and internal energy will be non-constant), 169a515125bSLeila Ghaffari // it should be transported across the domain, but velocity stays constant 170a515125bSLeila Ghaffari P = 1.; 171a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r*r); 172a515125bSLeila Ghaffari rho = P / (R*T); 173a515125bSLeila Ghaffari u[0] = mean_velocity[0]; 174a515125bSLeila Ghaffari u[1] = mean_velocity[1]; 175a515125bSLeila Ghaffari 176a515125bSLeila Ghaffari // Assign exact solution 177a515125bSLeila Ghaffari q[0] = rho; 178a515125bSLeila Ghaffari q[1] = rho * u[0]; 179a515125bSLeila Ghaffari q[2] = rho * u[1]; 180a515125bSLeila Ghaffari q[3] = rho * u[2]; 181a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 182a515125bSLeila Ghaffari break; 1830df2634dSLeila Ghaffari case 5: // non-smooth thermal bubble - cylinder 1840df2634dSLeila Ghaffari P = 1.; 1850df2634dSLeila Ghaffari T = 1. - (r < 1. ? S_bubble : 0.); 1860df2634dSLeila Ghaffari rho = P / (R*T); 1870df2634dSLeila Ghaffari u[0] = mean_velocity[0]; 1880df2634dSLeila Ghaffari u[1] = mean_velocity[1]; 1890df2634dSLeila Ghaffari 1900df2634dSLeila Ghaffari // Assign exact solution 1910df2634dSLeila Ghaffari q[0] = rho; 1920df2634dSLeila Ghaffari q[1] = rho * u[0]; 1930df2634dSLeila Ghaffari q[2] = rho * u[1]; 1940df2634dSLeila Ghaffari q[3] = rho * u[2]; 1950df2634dSLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 1960df2634dSLeila Ghaffari break; 197a515125bSLeila Ghaffari } 198a515125bSLeila Ghaffari // Return 199a515125bSLeila Ghaffari return 0; 200a515125bSLeila Ghaffari } 201a515125bSLeila Ghaffari 202a515125bSLeila Ghaffari // ***************************************************************************** 203139613f2SLeila Ghaffari // Helper function for computing flux Jacobian 204139613f2SLeila Ghaffari // ***************************************************************************** 205d8a22b9eSJed Brown CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], 206139613f2SLeila Ghaffari const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 207139613f2SLeila Ghaffari const CeedScalar gamma) { 208139613f2SLeila Ghaffari CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 209139613f2SLeila Ghaffari for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 210139613f2SLeila Ghaffari for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix 211139613f2SLeila Ghaffari dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2.)) : 0.) - u[i]*u[j]; 212139613f2SLeila Ghaffari for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix 213139613f2SLeila Ghaffari dF[i][0][k+1] = ((i==k) ? 1. : 0.); 214139613f2SLeila Ghaffari dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) + 215139613f2SLeila Ghaffari ((i==k) ? u[j] : 0.) - 216139613f2SLeila Ghaffari ((i==j) ? u[k] : 0.) * (gamma-1.); 217139613f2SLeila Ghaffari dF[i][4][k+1] = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) - 218139613f2SLeila Ghaffari (gamma-1.)*u[i]*u[k]; 219139613f2SLeila Ghaffari } 220139613f2SLeila Ghaffari dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.); 221139613f2SLeila Ghaffari } 222139613f2SLeila Ghaffari dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho); 223139613f2SLeila Ghaffari dF[i][4][4] = u[i] * gamma; 224139613f2SLeila Ghaffari } 225139613f2SLeila Ghaffari } 226139613f2SLeila Ghaffari 227139613f2SLeila Ghaffari // ***************************************************************************** 228d8a22b9eSJed Brown // Helper function for computing Tau elements (stabilization constant) 229d8a22b9eSJed Brown // Model from: 230d8a22b9eSJed Brown // Stabilized Methods for Compressible Flows, Hughes et al 2010 231d8a22b9eSJed Brown // 232d8a22b9eSJed Brown // Spatial criterion #2 - Tau is a 3x3 diagonal matrix 233d8a22b9eSJed Brown // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum) 234d8a22b9eSJed Brown // 235d8a22b9eSJed Brown // Where 236d8a22b9eSJed Brown // c_tau = stabilization constant (0.5 is reported as "optimal") 237d8a22b9eSJed Brown // h[i] = 2 length(dxdX[i]) 238d8a22b9eSJed Brown // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity ) 239d8a22b9eSJed Brown // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number ) 240d8a22b9eSJed Brown // rho(A[i]) = spectral radius of the convective flux Jacobian i, 241d8a22b9eSJed Brown // wave speed in direction i 242d8a22b9eSJed Brown // ***************************************************************************** 243d8a22b9eSJed Brown CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], 244d8a22b9eSJed Brown const CeedScalar dXdx[3][3], const CeedScalar u[3], 245d8a22b9eSJed Brown const CeedScalar sound_speed, const CeedScalar c_tau) { 246d8a22b9eSJed Brown for (int i=0; i<3; i++) { 247d8a22b9eSJed Brown // length of element in direction i 248d8a22b9eSJed Brown CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] + 249d8a22b9eSJed Brown dXdx[2][i]*dXdx[2][i]); 250d8a22b9eSJed Brown // fastest wave in direction i 251d8a22b9eSJed Brown CeedScalar fastest_wave = fabs(u[i]) + sound_speed; 252d8a22b9eSJed Brown Tau_x[i] = c_tau * h / fastest_wave; 253d8a22b9eSJed Brown } 254d8a22b9eSJed Brown } 255d8a22b9eSJed Brown 256d8a22b9eSJed Brown // ***************************************************************************** 257a515125bSLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex 258a515125bSLeila Ghaffari // ***************************************************************************** 259a515125bSLeila Ghaffari CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, 260a515125bSLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 261a515125bSLeila Ghaffari // Inputs 262a515125bSLeila Ghaffari const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 263a515125bSLeila Ghaffari 264a515125bSLeila Ghaffari // Outputs 265a515125bSLeila Ghaffari CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 266a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx; 267a515125bSLeila Ghaffari 268a515125bSLeila Ghaffari CeedPragmaSIMD 269a515125bSLeila Ghaffari // Quadrature Point Loop 270a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 271a515125bSLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 272139613f2SLeila Ghaffari CeedScalar q[5] = {0.}; 273a515125bSLeila Ghaffari 274a515125bSLeila Ghaffari Exact_Euler(3, context->curr_time, x, 5, q, ctx); 275a515125bSLeila Ghaffari 276a515125bSLeila Ghaffari for (CeedInt j=0; j<5; j++) 277a515125bSLeila Ghaffari q0[j][i] = q[j]; 278a515125bSLeila Ghaffari } // End of Quadrature Point Loop 279a515125bSLeila Ghaffari 280a515125bSLeila Ghaffari // Return 281a515125bSLeila Ghaffari return 0; 282a515125bSLeila Ghaffari } 283a515125bSLeila Ghaffari 284a515125bSLeila Ghaffari // ***************************************************************************** 285a515125bSLeila Ghaffari // This QFunction implements the following formulation of Euler equations 286a515125bSLeila Ghaffari // with explicit time stepping method 287a515125bSLeila Ghaffari // 288a515125bSLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation 289a515125bSLeila Ghaffari // form with state variables of density, momentum density, and total 290a515125bSLeila Ghaffari // energy density. 291a515125bSLeila Ghaffari // 292a515125bSLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E ) 293a515125bSLeila Ghaffari // rho - Mass Density 294a515125bSLeila Ghaffari // Ui - Momentum Density, Ui = rho ui 295a515125bSLeila Ghaffari // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2 296a515125bSLeila Ghaffari // 297a515125bSLeila Ghaffari // Euler Equations: 298a515125bSLeila Ghaffari // drho/dt + div( U ) = 0 299a515125bSLeila Ghaffari // dU/dt + div( rho (u x u) + P I3 ) = 0 300a515125bSLeila Ghaffari // dE/dt + div( (E + P) u ) = 0 301a515125bSLeila Ghaffari // 302a515125bSLeila Ghaffari // Equation of State: 303a515125bSLeila Ghaffari // P = (gamma - 1) (E - rho (u u) / 2) 304a515125bSLeila Ghaffari // 305a515125bSLeila Ghaffari // Constants: 306a515125bSLeila Ghaffari // cv , Specific heat, constant volume 307a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 308a515125bSLeila Ghaffari // g , Gravity 309a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio 310a515125bSLeila Ghaffari // ***************************************************************************** 311a515125bSLeila Ghaffari CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, 312a515125bSLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 313a515125bSLeila Ghaffari // *INDENT-OFF* 314a515125bSLeila Ghaffari // Inputs 315a515125bSLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 316139613f2SLeila Ghaffari (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 317a515125bSLeila Ghaffari (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 318a515125bSLeila Ghaffari // Outputs 319a515125bSLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 320a515125bSLeila Ghaffari (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 321a515125bSLeila Ghaffari 322139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx; 323d8a22b9eSJed Brown const CeedScalar c_tau = context->c_tau; 324a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 325a515125bSLeila Ghaffari 326a515125bSLeila Ghaffari CeedPragmaSIMD 327a515125bSLeila Ghaffari // Quadrature Point Loop 328a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 329a515125bSLeila Ghaffari // *INDENT-OFF* 330a515125bSLeila Ghaffari // Setup 331a515125bSLeila Ghaffari // -- Interp in 332a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 333a515125bSLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 334a515125bSLeila Ghaffari q[2][i] / rho, 335a515125bSLeila Ghaffari q[3][i] / rho 336a515125bSLeila Ghaffari }; 337a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 338139613f2SLeila Ghaffari const CeedScalar drho[3] = {dq[0][0][i], 339139613f2SLeila Ghaffari dq[1][0][i], 340139613f2SLeila Ghaffari dq[2][0][i] 341139613f2SLeila Ghaffari }; 342139613f2SLeila Ghaffari const CeedScalar dU[3][3] = {{dq[0][1][i], 343139613f2SLeila Ghaffari dq[1][1][i], 344139613f2SLeila Ghaffari dq[2][1][i]}, 345139613f2SLeila Ghaffari {dq[0][2][i], 346139613f2SLeila Ghaffari dq[1][2][i], 347139613f2SLeila Ghaffari dq[2][2][i]}, 348139613f2SLeila Ghaffari {dq[0][3][i], 349139613f2SLeila Ghaffari dq[1][3][i], 350139613f2SLeila Ghaffari dq[2][3][i]} 351139613f2SLeila Ghaffari }; 352139613f2SLeila Ghaffari const CeedScalar dE[3] = {dq[0][4][i], 353139613f2SLeila Ghaffari dq[1][4][i], 354139613f2SLeila Ghaffari dq[2][4][i] 355139613f2SLeila Ghaffari }; 356a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 357a515125bSLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 358a515125bSLeila Ghaffari // -- Interp-to-Grad q_data 359a515125bSLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 360a515125bSLeila Ghaffari // *INDENT-OFF* 361a515125bSLeila Ghaffari const CeedScalar dXdx[3][3] = {{q_data[1][i], 362a515125bSLeila Ghaffari q_data[2][i], 363a515125bSLeila Ghaffari q_data[3][i]}, 364a515125bSLeila Ghaffari {q_data[4][i], 365a515125bSLeila Ghaffari q_data[5][i], 366a515125bSLeila Ghaffari q_data[6][i]}, 367a515125bSLeila Ghaffari {q_data[7][i], 368a515125bSLeila Ghaffari q_data[8][i], 369a515125bSLeila Ghaffari q_data[9][i]} 370a515125bSLeila Ghaffari }; 371a515125bSLeila Ghaffari // *INDENT-ON* 372139613f2SLeila Ghaffari // dU/dx 373139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 374139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 375139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 376139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 377139613f2SLeila Ghaffari for (int j=0; j<3; j++) { 378139613f2SLeila Ghaffari for (int k=0; k<3; k++) { 379139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 380139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 381139613f2SLeila Ghaffari for (int l=0; l<3; l++) { 382139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 383139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 384139613f2SLeila Ghaffari } 385139613f2SLeila Ghaffari } 386139613f2SLeila Ghaffari } 387139613f2SLeila Ghaffari // Pressure 388a515125bSLeila Ghaffari const CeedScalar 389a515125bSLeila Ghaffari E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 390a515125bSLeila Ghaffari E_internal = E - E_kinetic, 391139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 392a515125bSLeila Ghaffari 393a515125bSLeila Ghaffari // The Physics 394a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 395a515125bSLeila Ghaffari for (int j=0; j<5; j++) { 396139613f2SLeila Ghaffari v[j][i] = 0.; 397a515125bSLeila Ghaffari for (int k=0; k<3; k++) 398139613f2SLeila Ghaffari dv[k][j][i] = 0.; 399a515125bSLeila Ghaffari } 400a515125bSLeila Ghaffari 401a515125bSLeila Ghaffari // -- Density 402a515125bSLeila Ghaffari // ---- u rho 403a515125bSLeila Ghaffari for (int j=0; j<3; j++) 404a515125bSLeila Ghaffari dv[j][0][i] += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 405a515125bSLeila Ghaffari rho*u[2]*dXdx[j][2]); 406a515125bSLeila Ghaffari // -- Momentum 407a515125bSLeila Ghaffari // ---- rho (u x u) + P I3 408a515125bSLeila Ghaffari for (int j=0; j<3; j++) 409a515125bSLeila Ghaffari for (int k=0; k<3; k++) 410139613f2SLeila Ghaffari dv[k][j+1][i] += wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] + 411139613f2SLeila Ghaffari (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] + 412139613f2SLeila Ghaffari (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]); 413a515125bSLeila Ghaffari // -- Total Energy Density 414a515125bSLeila Ghaffari // ---- (E + P) u 415a515125bSLeila Ghaffari for (int j=0; j<3; j++) 416a515125bSLeila Ghaffari dv[j][4][i] += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 417a515125bSLeila Ghaffari u[2]*dXdx[j][2]); 418139613f2SLeila Ghaffari 419139613f2SLeila Ghaffari // --Stabilization terms 420139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 421139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 422d8a22b9eSJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 423139613f2SLeila Ghaffari 424139613f2SLeila Ghaffari // ---- Transpose of the Jacobian 425139613f2SLeila Ghaffari CeedScalar jacob_F_conv_T[3][5][5]; 426139613f2SLeila Ghaffari for (int j=0; j<3; j++) 427139613f2SLeila Ghaffari for (int k=0; k<5; k++) 428139613f2SLeila Ghaffari for (int l=0; l<5; l++) 429139613f2SLeila Ghaffari jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k]; 430139613f2SLeila Ghaffari 431139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 432139613f2SLeila Ghaffari CeedScalar dqdx[5][3]; 433139613f2SLeila Ghaffari for (int j=0; j<3; j++) { 434139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j]; 435139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j]; 436139613f2SLeila Ghaffari for (int k=0; k<3; k++) 437139613f2SLeila Ghaffari dqdx[k+1][j] = dUdx[k][j]; 438139613f2SLeila Ghaffari } 439139613f2SLeila Ghaffari 440139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 441139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 442139613f2SLeila Ghaffari for (int j=0; j<3; j++) 443139613f2SLeila Ghaffari for (int k=0; k<5; k++) 444139613f2SLeila Ghaffari for (int l=0; l<5; l++) 445139613f2SLeila Ghaffari strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 446139613f2SLeila Ghaffari 447d8a22b9eSJed Brown // Stabilization 448d8a22b9eSJed Brown // -- Tau elements 449d8a22b9eSJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 450d8a22b9eSJed Brown CeedScalar Tau_x[3] = {0.}; 451d8a22b9eSJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 452139613f2SLeila Ghaffari 453d8a22b9eSJed Brown // -- Stabilization method: none or SU 454139613f2SLeila Ghaffari CeedScalar stab[5][3]; 455139613f2SLeila Ghaffari switch (context->stabilization) { 456139613f2SLeila Ghaffari case 0: // Galerkin 457139613f2SLeila Ghaffari break; 458139613f2SLeila Ghaffari case 1: // SU 459139613f2SLeila Ghaffari for (int j=0; j<3; j++) 460139613f2SLeila Ghaffari for (int k=0; k<5; k++) 461139613f2SLeila Ghaffari for (int l=0; l<5; l++) 462d8a22b9eSJed Brown stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_conv[l]; 463139613f2SLeila Ghaffari 464139613f2SLeila Ghaffari for (int j=0; j<5; j++) 465139613f2SLeila Ghaffari for (int k=0; k<3; k++) 466139613f2SLeila Ghaffari dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] + 467139613f2SLeila Ghaffari stab[j][1] * dXdx[k][1] + 468139613f2SLeila Ghaffari stab[j][2] * dXdx[k][2]); 469139613f2SLeila Ghaffari break; 470139613f2SLeila Ghaffari case 2: // SUPG is not implemented for explicit scheme 471139613f2SLeila Ghaffari break; 472139613f2SLeila Ghaffari } 473139613f2SLeila Ghaffari 474a515125bSLeila Ghaffari } // End Quadrature Point Loop 475a515125bSLeila Ghaffari 476a515125bSLeila Ghaffari // Return 477a515125bSLeila Ghaffari return 0; 478a515125bSLeila Ghaffari } 479a515125bSLeila Ghaffari 480a515125bSLeila Ghaffari // ***************************************************************************** 481a515125bSLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above) 482a515125bSLeila Ghaffari // with implicit time stepping method 483a515125bSLeila Ghaffari // 484a515125bSLeila Ghaffari // ***************************************************************************** 485a515125bSLeila Ghaffari CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, 486a515125bSLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 487a515125bSLeila Ghaffari // *INDENT-OFF* 488a515125bSLeila Ghaffari // Inputs 489a515125bSLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 490139613f2SLeila Ghaffari (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 491a515125bSLeila Ghaffari (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 492a515125bSLeila Ghaffari (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 493a515125bSLeila Ghaffari // Outputs 494a515125bSLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 495a515125bSLeila Ghaffari (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 496a515125bSLeila Ghaffari 497139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx; 498d8a22b9eSJed Brown const CeedScalar c_tau = context->c_tau; 499a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 500a515125bSLeila Ghaffari 501a515125bSLeila Ghaffari CeedPragmaSIMD 502a515125bSLeila Ghaffari // Quadrature Point Loop 503a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 504a515125bSLeila Ghaffari // *INDENT-OFF* 505a515125bSLeila Ghaffari // Setup 506a515125bSLeila Ghaffari // -- Interp in 507a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 508a515125bSLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 509a515125bSLeila Ghaffari q[2][i] / rho, 510a515125bSLeila Ghaffari q[3][i] / rho 511a515125bSLeila Ghaffari }; 512a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 513139613f2SLeila Ghaffari const CeedScalar drho[3] = {dq[0][0][i], 514139613f2SLeila Ghaffari dq[1][0][i], 515139613f2SLeila Ghaffari dq[2][0][i] 516139613f2SLeila Ghaffari }; 517139613f2SLeila Ghaffari const CeedScalar dU[3][3] = {{dq[0][1][i], 518139613f2SLeila Ghaffari dq[1][1][i], 519139613f2SLeila Ghaffari dq[2][1][i]}, 520139613f2SLeila Ghaffari {dq[0][2][i], 521139613f2SLeila Ghaffari dq[1][2][i], 522139613f2SLeila Ghaffari dq[2][2][i]}, 523139613f2SLeila Ghaffari {dq[0][3][i], 524139613f2SLeila Ghaffari dq[1][3][i], 525139613f2SLeila Ghaffari dq[2][3][i]} 526139613f2SLeila Ghaffari }; 527139613f2SLeila Ghaffari const CeedScalar dE[3] = {dq[0][4][i], 528139613f2SLeila Ghaffari dq[1][4][i], 529139613f2SLeila Ghaffari dq[2][4][i] 530139613f2SLeila Ghaffari }; 531a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 532a515125bSLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 533a515125bSLeila Ghaffari // -- Interp-to-Grad q_data 534a515125bSLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 535a515125bSLeila Ghaffari // *INDENT-OFF* 536a515125bSLeila Ghaffari const CeedScalar dXdx[3][3] = {{q_data[1][i], 537a515125bSLeila Ghaffari q_data[2][i], 538a515125bSLeila Ghaffari q_data[3][i]}, 539a515125bSLeila Ghaffari {q_data[4][i], 540a515125bSLeila Ghaffari q_data[5][i], 541a515125bSLeila Ghaffari q_data[6][i]}, 542a515125bSLeila Ghaffari {q_data[7][i], 543a515125bSLeila Ghaffari q_data[8][i], 544a515125bSLeila Ghaffari q_data[9][i]} 545a515125bSLeila Ghaffari }; 546a515125bSLeila Ghaffari // *INDENT-ON* 547139613f2SLeila Ghaffari // dU/dx 548139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 549139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 550139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 551139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 552139613f2SLeila Ghaffari for (int j=0; j<3; j++) { 553139613f2SLeila Ghaffari for (int k=0; k<3; k++) { 554139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 555139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 556139613f2SLeila Ghaffari for (int l=0; l<3; l++) { 557139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 558139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 559139613f2SLeila Ghaffari } 560139613f2SLeila Ghaffari } 561139613f2SLeila Ghaffari } 562a515125bSLeila Ghaffari const CeedScalar 563a515125bSLeila Ghaffari E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 564a515125bSLeila Ghaffari E_internal = E - E_kinetic, 565139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 566a515125bSLeila Ghaffari 567a515125bSLeila Ghaffari // The Physics 568a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 569a515125bSLeila Ghaffari for (int j=0; j<5; j++) { 570139613f2SLeila Ghaffari v[j][i] = 0.; 571a515125bSLeila Ghaffari for (int k=0; k<3; k++) 572139613f2SLeila Ghaffari dv[k][j][i] = 0.; 573a515125bSLeila Ghaffari } 574a515125bSLeila Ghaffari //-----mass matrix 575a515125bSLeila Ghaffari for (int j=0; j<5; j++) 576a515125bSLeila Ghaffari v[j][i] += wdetJ*q_dot[j][i]; 577a515125bSLeila Ghaffari 578a515125bSLeila Ghaffari // -- Density 579a515125bSLeila Ghaffari // ---- u rho 580a515125bSLeila Ghaffari for (int j=0; j<3; j++) 581a515125bSLeila Ghaffari dv[j][0][i] -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 582a515125bSLeila Ghaffari rho*u[2]*dXdx[j][2]); 583a515125bSLeila Ghaffari // -- Momentum 584a515125bSLeila Ghaffari // ---- rho (u x u) + P I3 585a515125bSLeila Ghaffari for (int j=0; j<3; j++) 586a515125bSLeila Ghaffari for (int k=0; k<3; k++) 587139613f2SLeila Ghaffari dv[k][j+1][i] -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] + 588139613f2SLeila Ghaffari (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] + 589139613f2SLeila Ghaffari (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]); 590a515125bSLeila Ghaffari // -- Total Energy Density 591a515125bSLeila Ghaffari // ---- (E + P) u 592a515125bSLeila Ghaffari for (int j=0; j<3; j++) 593a515125bSLeila Ghaffari dv[j][4][i] -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 594a515125bSLeila Ghaffari u[2]*dXdx[j][2]); 595139613f2SLeila Ghaffari 596139613f2SLeila Ghaffari // -- Stabilization terms 597139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 598139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 599d8a22b9eSJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 600139613f2SLeila Ghaffari 601139613f2SLeila Ghaffari // ---- Transpose of the Jacobian 602139613f2SLeila Ghaffari CeedScalar jacob_F_conv_T[3][5][5]; 603139613f2SLeila Ghaffari for (int j=0; j<3; j++) 604139613f2SLeila Ghaffari for (int k=0; k<5; k++) 605139613f2SLeila Ghaffari for (int l=0; l<5; l++) 606139613f2SLeila Ghaffari jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k]; 607139613f2SLeila Ghaffari 608139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 609139613f2SLeila Ghaffari CeedScalar dqdx[5][3]; 610139613f2SLeila Ghaffari for (int j=0; j<3; j++) { 611139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j]; 612139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j]; 613139613f2SLeila Ghaffari for (int k=0; k<3; k++) 614139613f2SLeila Ghaffari dqdx[k+1][j] = dUdx[k][j]; 615139613f2SLeila Ghaffari } 616139613f2SLeila Ghaffari 617139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 618139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 619139613f2SLeila Ghaffari for (int j=0; j<3; j++) 620139613f2SLeila Ghaffari for (int k=0; k<5; k++) 621139613f2SLeila Ghaffari for (int l=0; l<5; l++) 622139613f2SLeila Ghaffari strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 623139613f2SLeila Ghaffari 624139613f2SLeila Ghaffari // ---- Strong residual 625139613f2SLeila Ghaffari CeedScalar strong_res[5]; 626139613f2SLeila Ghaffari for (int j=0; j<5; j++) 627139613f2SLeila Ghaffari strong_res[j] = q_dot[j][i] + strong_conv[j]; 628139613f2SLeila Ghaffari 629d8a22b9eSJed Brown // Stabilization 630d8a22b9eSJed Brown // -- Tau elements 631d8a22b9eSJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho); 632d8a22b9eSJed Brown CeedScalar Tau_x[3] = {0.}; 633d8a22b9eSJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau); 634139613f2SLeila Ghaffari 635d8a22b9eSJed Brown // -- Stabilization method: none, SU, or SUPG 636139613f2SLeila Ghaffari CeedScalar stab[5][3]; 637139613f2SLeila Ghaffari switch (context->stabilization) { 638139613f2SLeila Ghaffari case 0: // Galerkin 639139613f2SLeila Ghaffari break; 640139613f2SLeila Ghaffari case 1: // SU 641139613f2SLeila Ghaffari for (int j=0; j<3; j++) 642139613f2SLeila Ghaffari for (int k=0; k<5; k++) 643139613f2SLeila Ghaffari for (int l=0; l<5; l++) 644d8a22b9eSJed Brown stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_conv[l]; 645139613f2SLeila Ghaffari 646139613f2SLeila Ghaffari for (int j=0; j<5; j++) 647139613f2SLeila Ghaffari for (int k=0; k<3; k++) 648139613f2SLeila Ghaffari dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 649139613f2SLeila Ghaffari stab[j][1] * dXdx[k][1] + 650139613f2SLeila Ghaffari stab[j][2] * dXdx[k][2]); 651139613f2SLeila Ghaffari break; 652139613f2SLeila Ghaffari case 2: // SUPG 653139613f2SLeila Ghaffari for (int j=0; j<3; j++) 654139613f2SLeila Ghaffari for (int k=0; k<5; k++) 655139613f2SLeila Ghaffari for (int l=0; l<5; l++) 656d8a22b9eSJed Brown stab[k][j] = jacob_F_conv_T[j][k][l] * Tau_x[j] * strong_res[l]; 657139613f2SLeila Ghaffari 658139613f2SLeila Ghaffari for (int j=0; j<5; j++) 659139613f2SLeila Ghaffari for (int k=0; k<3; k++) 660139613f2SLeila Ghaffari dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 661139613f2SLeila Ghaffari stab[j][1] * dXdx[k][1] + 662139613f2SLeila Ghaffari stab[j][2] * dXdx[k][2]); 663139613f2SLeila Ghaffari break; 664139613f2SLeila Ghaffari } 665a515125bSLeila Ghaffari } // End Quadrature Point Loop 666a515125bSLeila Ghaffari 667a515125bSLeila Ghaffari // Return 668a515125bSLeila Ghaffari return 0; 669a515125bSLeila Ghaffari } 670a515125bSLeila Ghaffari // ***************************************************************************** 671002797a3SLeila Ghaffari // This QFunction sets the inflow boundary conditions for 672002797a3SLeila Ghaffari // the traveling vortex problem. 673a515125bSLeila Ghaffari // 674a515125bSLeila Ghaffari // Prescribed T_inlet and P_inlet are converted to conservative variables 675a515125bSLeila Ghaffari // and applied weakly. 676a515125bSLeila Ghaffari // 677a515125bSLeila Ghaffari // ***************************************************************************** 678002797a3SLeila Ghaffari CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, 679a515125bSLeila Ghaffari const CeedScalar *const *in, 680a515125bSLeila Ghaffari CeedScalar *const *out) { 681a515125bSLeila Ghaffari // *INDENT-OFF* 682a515125bSLeila Ghaffari // Inputs 683002797a3SLeila Ghaffari const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 684a515125bSLeila Ghaffari // Outputs 685a515125bSLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 686a515125bSLeila Ghaffari // *INDENT-ON* 687a515125bSLeila Ghaffari EulerContext context = (EulerContext)ctx; 688a515125bSLeila Ghaffari const int euler_test = context->euler_test; 689a515125bSLeila Ghaffari const bool implicit = context->implicit; 690a515125bSLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 691a515125bSLeila Ghaffari const CeedScalar cv = 2.5; 692a515125bSLeila Ghaffari const CeedScalar R = 1.; 693a515125bSLeila Ghaffari CeedScalar T_inlet; 694a515125bSLeila Ghaffari CeedScalar P_inlet; 695a515125bSLeila Ghaffari 696a515125bSLeila Ghaffari // For test cases 1 and 3 the background velocity is zero 697a515125bSLeila Ghaffari if (euler_test == 1 || euler_test == 3) 698a515125bSLeila Ghaffari for (CeedInt i=0; i<3; i++) mean_velocity[i] = 0.; 699a515125bSLeila Ghaffari 700a515125bSLeila Ghaffari // For test cases 1 and 2, T_inlet = T_inlet = 0.4 701a515125bSLeila Ghaffari if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4; 702a515125bSLeila Ghaffari else T_inlet = P_inlet = 1.; 703a515125bSLeila Ghaffari 704a515125bSLeila Ghaffari CeedPragmaSIMD 705a515125bSLeila Ghaffari // Quadrature Point Loop 706a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 707a515125bSLeila Ghaffari // Setup 708a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 709a515125bSLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 710a515125bSLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 711a515125bSLeila Ghaffari // We can effect this by swapping the sign on this weight 712a515125bSLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 713002797a3SLeila Ghaffari // ---- Normal vect 714a515125bSLeila Ghaffari const CeedScalar norm[3] = {q_data_sur[1][i], 715a515125bSLeila Ghaffari q_data_sur[2][i], 716a515125bSLeila Ghaffari q_data_sur[3][i] 717a515125bSLeila Ghaffari }; 718a515125bSLeila Ghaffari 719a515125bSLeila Ghaffari // face_normal = Normal vector of the face 720a515125bSLeila Ghaffari const CeedScalar face_normal = norm[0]*mean_velocity[0] + 721a515125bSLeila Ghaffari norm[1]*mean_velocity[1] + 722a515125bSLeila Ghaffari norm[2]*mean_velocity[2]; 723a515125bSLeila Ghaffari // The Physics 724a515125bSLeila Ghaffari // Zero v so all future terms can safely sum into it 725139613f2SLeila Ghaffari for (int j=0; j<5; j++) v[j][i] = 0.; 726a515125bSLeila Ghaffari 727a515125bSLeila Ghaffari // Implementing in/outflow BCs 728002797a3SLeila Ghaffari if (face_normal > 0) { 729a515125bSLeila Ghaffari } else { // inflow 730a515125bSLeila Ghaffari const CeedScalar rho_inlet = P_inlet/(R*T_inlet); 731a515125bSLeila Ghaffari const CeedScalar E_kinetic_inlet = (mean_velocity[0]*mean_velocity[0] + 732a515125bSLeila Ghaffari mean_velocity[1]*mean_velocity[1]) / 2.; 733a515125bSLeila Ghaffari // incoming total energy 734a515125bSLeila Ghaffari const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet); 735a515125bSLeila Ghaffari 736a515125bSLeila Ghaffari // The Physics 737a515125bSLeila Ghaffari // -- Density 738a515125bSLeila Ghaffari v[0][i] -= wdetJb * rho_inlet * face_normal; 739a515125bSLeila Ghaffari 740a515125bSLeila Ghaffari // -- Momentum 741a515125bSLeila Ghaffari for (int j=0; j<3; j++) 742a515125bSLeila Ghaffari v[j+1][i] -= wdetJb *(rho_inlet * face_normal * mean_velocity[j] + 743a515125bSLeila Ghaffari norm[j] * P_inlet); 744a515125bSLeila Ghaffari 745a515125bSLeila Ghaffari // -- Total Energy Density 746a515125bSLeila Ghaffari v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet); 747a515125bSLeila Ghaffari } 748a515125bSLeila Ghaffari 749a515125bSLeila Ghaffari } // End Quadrature Point Loop 750a515125bSLeila Ghaffari return 0; 751a515125bSLeila Ghaffari } 752a515125bSLeila Ghaffari 753a515125bSLeila Ghaffari // ***************************************************************************** 75468ef3d20SLeila Ghaffari // This QFunction sets the outflow boundary conditions for 75568ef3d20SLeila Ghaffari // the Euler solver. 75668ef3d20SLeila Ghaffari // 75768ef3d20SLeila Ghaffari // Outflow BCs: 75868ef3d20SLeila Ghaffari // The validity of the weak form of the governing equations is 75968ef3d20SLeila Ghaffari // extended to the outflow. 76068ef3d20SLeila Ghaffari // 76168ef3d20SLeila Ghaffari // ***************************************************************************** 76268ef3d20SLeila Ghaffari CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, 76368ef3d20SLeila Ghaffari const CeedScalar *const *in, 76468ef3d20SLeila Ghaffari CeedScalar *const *out) { 76568ef3d20SLeila Ghaffari // *INDENT-OFF* 76668ef3d20SLeila Ghaffari // Inputs 76768ef3d20SLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 76868ef3d20SLeila Ghaffari (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 76968ef3d20SLeila Ghaffari // Outputs 77068ef3d20SLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 77168ef3d20SLeila Ghaffari // *INDENT-ON* 77268ef3d20SLeila Ghaffari EulerContext context = (EulerContext)ctx; 77368ef3d20SLeila Ghaffari const bool implicit = context->implicit; 77468ef3d20SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 77568ef3d20SLeila Ghaffari 77668ef3d20SLeila Ghaffari const CeedScalar gamma = 1.4; 77768ef3d20SLeila Ghaffari 77868ef3d20SLeila Ghaffari CeedPragmaSIMD 77968ef3d20SLeila Ghaffari // Quadrature Point Loop 78068ef3d20SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 78168ef3d20SLeila Ghaffari // Setup 78268ef3d20SLeila Ghaffari // -- Interp in 78368ef3d20SLeila Ghaffari const CeedScalar rho = q[0][i]; 78468ef3d20SLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 78568ef3d20SLeila Ghaffari q[2][i] / rho, 78668ef3d20SLeila Ghaffari q[3][i] / rho 78768ef3d20SLeila Ghaffari }; 78868ef3d20SLeila Ghaffari const CeedScalar E = q[4][i]; 78968ef3d20SLeila Ghaffari 79068ef3d20SLeila Ghaffari // -- Interp-to-Interp q_data 79168ef3d20SLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 79268ef3d20SLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 79368ef3d20SLeila Ghaffari // We can effect this by swapping the sign on this weight 79468ef3d20SLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 79568ef3d20SLeila Ghaffari // ---- Normal vectors 79668ef3d20SLeila Ghaffari const CeedScalar norm[3] = {q_data_sur[1][i], 79768ef3d20SLeila Ghaffari q_data_sur[2][i], 79868ef3d20SLeila Ghaffari q_data_sur[3][i] 79968ef3d20SLeila Ghaffari }; 80068ef3d20SLeila Ghaffari 80168ef3d20SLeila Ghaffari // face_normal = Normal vector of the face 80268ef3d20SLeila Ghaffari const CeedScalar face_normal = norm[0]*mean_velocity[0] + 80368ef3d20SLeila Ghaffari norm[1]*mean_velocity[1] + 80468ef3d20SLeila Ghaffari norm[2]*mean_velocity[2]; 80568ef3d20SLeila Ghaffari // The Physics 80668ef3d20SLeila Ghaffari // Zero v so all future terms can safely sum into it 80768ef3d20SLeila Ghaffari for (int j=0; j<5; j++) v[j][i] = 0; 80868ef3d20SLeila Ghaffari 80968ef3d20SLeila Ghaffari // Implementing in/outflow BCs 81068ef3d20SLeila Ghaffari if (face_normal > 0) { // outflow 81168ef3d20SLeila Ghaffari const CeedScalar E_kinetic = (u[0]*u[0] + u[1]*u[1]) / 2.; 81268ef3d20SLeila Ghaffari const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure 81368ef3d20SLeila Ghaffari const CeedScalar u_normal = norm[0]*u[0] + norm[1]*u[1] + 81468ef3d20SLeila Ghaffari norm[2]*u[2]; // Normal velocity 81568ef3d20SLeila Ghaffari // The Physics 81668ef3d20SLeila Ghaffari // -- Density 81768ef3d20SLeila Ghaffari v[0][i] -= wdetJb * rho * u_normal; 81868ef3d20SLeila Ghaffari 81968ef3d20SLeila Ghaffari // -- Momentum 82068ef3d20SLeila Ghaffari for (int j=0; j<3; j++) 82168ef3d20SLeila Ghaffari v[j+1][i] -= wdetJb *(rho * u_normal * u[j] + norm[j] * P); 82268ef3d20SLeila Ghaffari 82368ef3d20SLeila Ghaffari // -- Total Energy Density 82468ef3d20SLeila Ghaffari v[4][i] -= wdetJb * u_normal * (E + P); 82568ef3d20SLeila Ghaffari } 82668ef3d20SLeila Ghaffari } // End Quadrature Point Loop 82768ef3d20SLeila Ghaffari return 0; 82868ef3d20SLeila Ghaffari } 82968ef3d20SLeila Ghaffari 83068ef3d20SLeila Ghaffari // ***************************************************************************** 831a515125bSLeila Ghaffari 832a515125bSLeila Ghaffari #endif // eulervortex_h 833