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> 29a515125bSLeila Ghaffari 30a515125bSLeila Ghaffari #ifndef M_PI 31a515125bSLeila Ghaffari #define M_PI 3.14159265358979323846 32a515125bSLeila Ghaffari #endif 33a515125bSLeila Ghaffari 34a515125bSLeila Ghaffari #ifndef euler_context_struct 35a515125bSLeila Ghaffari #define euler_context_struct 36a515125bSLeila Ghaffari typedef struct EulerContext_ *EulerContext; 37a515125bSLeila Ghaffari struct EulerContext_ { 38a515125bSLeila Ghaffari CeedScalar center[3]; 39a515125bSLeila Ghaffari CeedScalar curr_time; 40a515125bSLeila Ghaffari CeedScalar vortex_strength; 41a515125bSLeila Ghaffari CeedScalar mean_velocity[3]; 42a515125bSLeila Ghaffari bool implicit; 43*139613f2SLeila Ghaffari int euler_test; 44*139613f2SLeila Ghaffari int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG 45a515125bSLeila Ghaffari }; 46a515125bSLeila Ghaffari #endif 47a515125bSLeila Ghaffari 48a515125bSLeila Ghaffari // ***************************************************************************** 49a515125bSLeila Ghaffari // This function sets the initial conditions 50a515125bSLeila Ghaffari // 51a515125bSLeila Ghaffari // Temperature: 52a515125bSLeila Ghaffari // T = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2) 53a515125bSLeila Ghaffari // Density: 54a515125bSLeila Ghaffari // rho = (T/S_vortex)^(1 / (gamma - 1)) 55a515125bSLeila Ghaffari // Pressure: 56a515125bSLeila Ghaffari // P = rho * T 57a515125bSLeila Ghaffari // Velocity: 58a515125bSLeila Ghaffari // ui = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi) 59a515125bSLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 ) 60a515125bSLeila Ghaffari // Velocity/Momentum Density: 61a515125bSLeila Ghaffari // Ui = rho ui 62a515125bSLeila Ghaffari // Total Energy: 63a515125bSLeila Ghaffari // E = P / (gamma - 1) + rho (u u)/2 64a515125bSLeila Ghaffari // 65a515125bSLeila Ghaffari // Constants: 66a515125bSLeila Ghaffari // cv , Specific heat, constant volume 67a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 68a515125bSLeila Ghaffari // vortex_strength , Strength of vortex 69a515125bSLeila Ghaffari // center , Location of bubble center 70a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio 71a515125bSLeila Ghaffari // 72a515125bSLeila Ghaffari // ***************************************************************************** 73a515125bSLeila Ghaffari 74a515125bSLeila Ghaffari // ***************************************************************************** 75a515125bSLeila Ghaffari // This helper function provides support for the exact, time-dependent solution 76a515125bSLeila Ghaffari // (currently not implemented) and IC formulation for Euler traveling vortex 77a515125bSLeila Ghaffari // ***************************************************************************** 78a515125bSLeila Ghaffari CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, 79a515125bSLeila Ghaffari const CeedScalar X[], CeedInt Nf, CeedScalar q[], 80a515125bSLeila Ghaffari void *ctx) { 81a515125bSLeila Ghaffari // Context 82a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx; 83a515125bSLeila Ghaffari const CeedScalar vortex_strength = context->vortex_strength; 84a515125bSLeila Ghaffari const CeedScalar *center = context->center; // Center of the domain 85a515125bSLeila Ghaffari const CeedScalar *mean_velocity = context->mean_velocity; 86a515125bSLeila Ghaffari 87a515125bSLeila Ghaffari // Setup 88a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 89a515125bSLeila Ghaffari const CeedScalar cv = 2.5; 90a515125bSLeila Ghaffari const CeedScalar R = 1.; 91a515125bSLeila Ghaffari const CeedScalar x = X[0], y = X[1]; // Coordinates 92a515125bSLeila Ghaffari // Vortex center 93a515125bSLeila Ghaffari const CeedScalar xc = center[0] + mean_velocity[0] * time; 94a515125bSLeila Ghaffari const CeedScalar yc = center[1] + mean_velocity[1] * time; 95a515125bSLeila Ghaffari 96a515125bSLeila Ghaffari const CeedScalar x0 = x - xc; 97a515125bSLeila Ghaffari const CeedScalar y0 = y - yc; 98a515125bSLeila Ghaffari const CeedScalar r = sqrt( x0*x0 + y0*y0 ); 99a515125bSLeila Ghaffari const CeedScalar C = vortex_strength * exp((1. - r*r)/2.) / (2. * M_PI); 100*139613f2SLeila Ghaffari const CeedScalar delta_T = - (gamma - 1.) * vortex_strength * vortex_strength * 101*139613f2SLeila Ghaffari exp(1 - r*r) / (8. * gamma * M_PI * M_PI); 102a515125bSLeila Ghaffari const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma 103a515125bSLeila Ghaffari const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / 104a515125bSLeila Ghaffari (8.*gamma*M_PI*M_PI); 105a515125bSLeila Ghaffari CeedScalar rho, P, T, E, u[3] = {0.}; 106a515125bSLeila Ghaffari 107a515125bSLeila Ghaffari // Initial Conditions 108a515125bSLeila Ghaffari switch (context->euler_test) { 109a515125bSLeila Ghaffari case 0: // Traveling vortex 110a515125bSLeila Ghaffari T = 1 + delta_T; 111a515125bSLeila Ghaffari // P = rho * T 112a515125bSLeila Ghaffari // P = S * rho^gamma 113a515125bSLeila Ghaffari // Solve for rho, then substitute for P 114*139613f2SLeila Ghaffari rho = pow(T/S_vortex, 1 / (gamma - 1.)); 115a515125bSLeila Ghaffari P = rho * T; 116a515125bSLeila Ghaffari u[0] = mean_velocity[0] - C*y0; 117a515125bSLeila Ghaffari u[1] = mean_velocity[1] + C*x0; 118a515125bSLeila Ghaffari 119a515125bSLeila Ghaffari // Assign exact solution 120a515125bSLeila Ghaffari q[0] = rho; 121a515125bSLeila Ghaffari q[1] = rho * u[0]; 122a515125bSLeila Ghaffari q[2] = rho * u[1]; 123a515125bSLeila Ghaffari q[3] = rho * u[2]; 124a515125bSLeila Ghaffari q[4] = P / (gamma - 1.) + rho * (u[0]*u[0] + u[1]*u[1]) / 2.; 125a515125bSLeila Ghaffari break; 126a515125bSLeila Ghaffari case 1: // Constant zero velocity, density constant, total energy constant 127a515125bSLeila Ghaffari rho = 1.; 128a515125bSLeila Ghaffari E = 2.; 129a515125bSLeila Ghaffari 130a515125bSLeila Ghaffari // Assign exact solution 131a515125bSLeila Ghaffari q[0] = rho; 132a515125bSLeila Ghaffari q[1] = rho * u[0]; 133a515125bSLeila Ghaffari q[2] = rho * u[1]; 134a515125bSLeila Ghaffari q[3] = rho * u[2]; 135a515125bSLeila Ghaffari q[4] = E; 136a515125bSLeila Ghaffari break; 137a515125bSLeila Ghaffari case 2: // Constant nonzero velocity, density constant, total energy constant 138a515125bSLeila Ghaffari rho = 1.; 139a515125bSLeila Ghaffari E = 2.; 140a515125bSLeila Ghaffari u[0] = mean_velocity[0]; 141a515125bSLeila Ghaffari u[1] = mean_velocity[1]; 142a515125bSLeila Ghaffari 143a515125bSLeila Ghaffari // Assign exact solution 144a515125bSLeila Ghaffari q[0] = rho; 145a515125bSLeila Ghaffari q[1] = rho * u[0]; 146a515125bSLeila Ghaffari q[2] = rho * u[1]; 147a515125bSLeila Ghaffari q[3] = rho * u[2]; 148a515125bSLeila Ghaffari q[4] = E; 149a515125bSLeila Ghaffari break; 150a515125bSLeila Ghaffari case 3: // Velocity zero, pressure constant 151a515125bSLeila Ghaffari // (so density and internal energy will be non-constant), 152a515125bSLeila Ghaffari // but the velocity should stay zero and the bubble won't diffuse 153a515125bSLeila Ghaffari // (for Euler, where there is no thermal conductivity) 154a515125bSLeila Ghaffari P = 1.; 155a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r*r); 156a515125bSLeila Ghaffari rho = P / (R*T); 157a515125bSLeila Ghaffari 158a515125bSLeila Ghaffari // Assign exact solution 159a515125bSLeila Ghaffari q[0] = rho; 160a515125bSLeila Ghaffari q[1] = rho * u[0]; 161a515125bSLeila Ghaffari q[2] = rho * u[1]; 162a515125bSLeila Ghaffari q[3] = rho * u[2]; 163a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 164a515125bSLeila Ghaffari break; 165a515125bSLeila Ghaffari case 4: // Constant nonzero velocity, pressure constant 166a515125bSLeila Ghaffari // (so density and internal energy will be non-constant), 167a515125bSLeila Ghaffari // it should be transported across the domain, but velocity stays constant 168a515125bSLeila Ghaffari P = 1.; 169a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r*r); 170a515125bSLeila Ghaffari rho = P / (R*T); 171a515125bSLeila Ghaffari u[0] = mean_velocity[0]; 172a515125bSLeila Ghaffari u[1] = mean_velocity[1]; 173a515125bSLeila Ghaffari 174a515125bSLeila Ghaffari // Assign exact solution 175a515125bSLeila Ghaffari q[0] = rho; 176a515125bSLeila Ghaffari q[1] = rho * u[0]; 177a515125bSLeila Ghaffari q[2] = rho * u[1]; 178a515125bSLeila Ghaffari q[3] = rho * u[2]; 179a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.); 180a515125bSLeila Ghaffari break; 181a515125bSLeila Ghaffari } 182a515125bSLeila Ghaffari // Return 183a515125bSLeila Ghaffari return 0; 184a515125bSLeila Ghaffari } 185a515125bSLeila Ghaffari 186a515125bSLeila Ghaffari // ***************************************************************************** 187*139613f2SLeila Ghaffari // Helper function for computing flux Jacobian 188*139613f2SLeila Ghaffari // ***************************************************************************** 189*139613f2SLeila Ghaffari CEED_QFUNCTION_HELPER void computeFluxJacobian_Euler(CeedScalar dF[3][5][5], 190*139613f2SLeila Ghaffari const CeedScalar rho, const CeedScalar u[3], const CeedScalar E, 191*139613f2SLeila Ghaffari const CeedScalar gamma) { 192*139613f2SLeila Ghaffari CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square 193*139613f2SLeila Ghaffari for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions 194*139613f2SLeila Ghaffari for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix 195*139613f2SLeila Ghaffari dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2.)) : 0.) - u[i]*u[j]; 196*139613f2SLeila Ghaffari for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix 197*139613f2SLeila Ghaffari dF[i][0][k+1] = ((i==k) ? 1. : 0.); 198*139613f2SLeila Ghaffari dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) + 199*139613f2SLeila Ghaffari ((i==k) ? u[j] : 0.) - 200*139613f2SLeila Ghaffari ((i==j) ? u[k] : 0.) * (gamma-1.); 201*139613f2SLeila Ghaffari dF[i][4][k+1] = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) - 202*139613f2SLeila Ghaffari (gamma-1.)*u[i]*u[k]; 203*139613f2SLeila Ghaffari } 204*139613f2SLeila Ghaffari dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.); 205*139613f2SLeila Ghaffari } 206*139613f2SLeila Ghaffari dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho); 207*139613f2SLeila Ghaffari dF[i][4][4] = u[i] * gamma; 208*139613f2SLeila Ghaffari } 209*139613f2SLeila Ghaffari } 210*139613f2SLeila Ghaffari 211*139613f2SLeila Ghaffari // ***************************************************************************** 212a515125bSLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex 213a515125bSLeila Ghaffari // ***************************************************************************** 214a515125bSLeila Ghaffari CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, 215a515125bSLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 216a515125bSLeila Ghaffari // Inputs 217a515125bSLeila Ghaffari const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 218a515125bSLeila Ghaffari 219a515125bSLeila Ghaffari // Outputs 220a515125bSLeila Ghaffari CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 221a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx; 222a515125bSLeila Ghaffari 223a515125bSLeila Ghaffari CeedPragmaSIMD 224a515125bSLeila Ghaffari // Quadrature Point Loop 225a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 226a515125bSLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 227*139613f2SLeila Ghaffari CeedScalar q[5] = {0.}; 228a515125bSLeila Ghaffari 229a515125bSLeila Ghaffari Exact_Euler(3, context->curr_time, x, 5, q, ctx); 230a515125bSLeila Ghaffari 231a515125bSLeila Ghaffari for (CeedInt j=0; j<5; j++) 232a515125bSLeila Ghaffari q0[j][i] = q[j]; 233a515125bSLeila Ghaffari } // End of Quadrature Point Loop 234a515125bSLeila Ghaffari 235a515125bSLeila Ghaffari // Return 236a515125bSLeila Ghaffari return 0; 237a515125bSLeila Ghaffari } 238a515125bSLeila Ghaffari 239a515125bSLeila Ghaffari // ***************************************************************************** 240a515125bSLeila Ghaffari // This QFunction implements the following formulation of Euler equations 241a515125bSLeila Ghaffari // with explicit time stepping method 242a515125bSLeila Ghaffari // 243a515125bSLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation 244a515125bSLeila Ghaffari // form with state variables of density, momentum density, and total 245a515125bSLeila Ghaffari // energy density. 246a515125bSLeila Ghaffari // 247a515125bSLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E ) 248a515125bSLeila Ghaffari // rho - Mass Density 249a515125bSLeila Ghaffari // Ui - Momentum Density, Ui = rho ui 250a515125bSLeila Ghaffari // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2 251a515125bSLeila Ghaffari // 252a515125bSLeila Ghaffari // Euler Equations: 253a515125bSLeila Ghaffari // drho/dt + div( U ) = 0 254a515125bSLeila Ghaffari // dU/dt + div( rho (u x u) + P I3 ) = 0 255a515125bSLeila Ghaffari // dE/dt + div( (E + P) u ) = 0 256a515125bSLeila Ghaffari // 257a515125bSLeila Ghaffari // Equation of State: 258a515125bSLeila Ghaffari // P = (gamma - 1) (E - rho (u u) / 2) 259a515125bSLeila Ghaffari // 260a515125bSLeila Ghaffari // Constants: 261a515125bSLeila Ghaffari // cv , Specific heat, constant volume 262a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 263a515125bSLeila Ghaffari // g , Gravity 264a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio 265a515125bSLeila Ghaffari // ***************************************************************************** 266a515125bSLeila Ghaffari CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, 267a515125bSLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 268a515125bSLeila Ghaffari // *INDENT-OFF* 269a515125bSLeila Ghaffari // Inputs 270a515125bSLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 271*139613f2SLeila Ghaffari (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 272a515125bSLeila Ghaffari (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 273a515125bSLeila Ghaffari // Outputs 274a515125bSLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 275a515125bSLeila Ghaffari (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 276a515125bSLeila Ghaffari 277*139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx; 278a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 279*139613f2SLeila Ghaffari const CeedScalar cv = 2.5; 280a515125bSLeila Ghaffari 281a515125bSLeila Ghaffari CeedPragmaSIMD 282a515125bSLeila Ghaffari // Quadrature Point Loop 283a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 284a515125bSLeila Ghaffari // *INDENT-OFF* 285a515125bSLeila Ghaffari // Setup 286a515125bSLeila Ghaffari // -- Interp in 287a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 288a515125bSLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 289a515125bSLeila Ghaffari q[2][i] / rho, 290a515125bSLeila Ghaffari q[3][i] / rho 291a515125bSLeila Ghaffari }; 292a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 293*139613f2SLeila Ghaffari const CeedScalar drho[3] = {dq[0][0][i], 294*139613f2SLeila Ghaffari dq[1][0][i], 295*139613f2SLeila Ghaffari dq[2][0][i] 296*139613f2SLeila Ghaffari }; 297*139613f2SLeila Ghaffari const CeedScalar dU[3][3] = {{dq[0][1][i], 298*139613f2SLeila Ghaffari dq[1][1][i], 299*139613f2SLeila Ghaffari dq[2][1][i]}, 300*139613f2SLeila Ghaffari {dq[0][2][i], 301*139613f2SLeila Ghaffari dq[1][2][i], 302*139613f2SLeila Ghaffari dq[2][2][i]}, 303*139613f2SLeila Ghaffari {dq[0][3][i], 304*139613f2SLeila Ghaffari dq[1][3][i], 305*139613f2SLeila Ghaffari dq[2][3][i]} 306*139613f2SLeila Ghaffari }; 307*139613f2SLeila Ghaffari const CeedScalar dE[3] = {dq[0][4][i], 308*139613f2SLeila Ghaffari dq[1][4][i], 309*139613f2SLeila Ghaffari dq[2][4][i] 310*139613f2SLeila Ghaffari }; 311a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 312a515125bSLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 313a515125bSLeila Ghaffari // -- Interp-to-Grad q_data 314a515125bSLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 315a515125bSLeila Ghaffari // *INDENT-OFF* 316a515125bSLeila Ghaffari const CeedScalar dXdx[3][3] = {{q_data[1][i], 317a515125bSLeila Ghaffari q_data[2][i], 318a515125bSLeila Ghaffari q_data[3][i]}, 319a515125bSLeila Ghaffari {q_data[4][i], 320a515125bSLeila Ghaffari q_data[5][i], 321a515125bSLeila Ghaffari q_data[6][i]}, 322a515125bSLeila Ghaffari {q_data[7][i], 323a515125bSLeila Ghaffari q_data[8][i], 324a515125bSLeila Ghaffari q_data[9][i]} 325a515125bSLeila Ghaffari }; 326a515125bSLeila Ghaffari // *INDENT-ON* 327*139613f2SLeila Ghaffari // dU/dx 328*139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 329*139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 330*139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 331*139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 332*139613f2SLeila Ghaffari for (int j=0; j<3; j++) { 333*139613f2SLeila Ghaffari for (int k=0; k<3; k++) { 334*139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 335*139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 336*139613f2SLeila Ghaffari for (int l=0; l<3; l++) { 337*139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 338*139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 339*139613f2SLeila Ghaffari } 340*139613f2SLeila Ghaffari } 341*139613f2SLeila Ghaffari } 342*139613f2SLeila Ghaffari // Pressure 343a515125bSLeila Ghaffari const CeedScalar 344a515125bSLeila Ghaffari E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 345a515125bSLeila Ghaffari E_internal = E - E_kinetic, 346*139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 347a515125bSLeila Ghaffari 348a515125bSLeila Ghaffari // The Physics 349a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 350a515125bSLeila Ghaffari for (int j=0; j<5; j++) { 351*139613f2SLeila Ghaffari v[j][i] = 0.; 352a515125bSLeila Ghaffari for (int k=0; k<3; k++) 353*139613f2SLeila Ghaffari dv[k][j][i] = 0.; 354a515125bSLeila Ghaffari } 355a515125bSLeila Ghaffari 356a515125bSLeila Ghaffari // -- Density 357a515125bSLeila Ghaffari // ---- u rho 358a515125bSLeila Ghaffari for (int j=0; j<3; j++) 359a515125bSLeila Ghaffari dv[j][0][i] += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 360a515125bSLeila Ghaffari rho*u[2]*dXdx[j][2]); 361a515125bSLeila Ghaffari // -- Momentum 362a515125bSLeila Ghaffari // ---- rho (u x u) + P I3 363a515125bSLeila Ghaffari for (int j=0; j<3; j++) 364a515125bSLeila Ghaffari for (int k=0; k<3; k++) 365*139613f2SLeila Ghaffari dv[k][j+1][i] += wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] + 366*139613f2SLeila Ghaffari (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] + 367*139613f2SLeila Ghaffari (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]); 368a515125bSLeila Ghaffari // -- Total Energy Density 369a515125bSLeila Ghaffari // ---- (E + P) u 370a515125bSLeila Ghaffari for (int j=0; j<3; j++) 371a515125bSLeila Ghaffari dv[j][4][i] += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 372a515125bSLeila Ghaffari u[2]*dXdx[j][2]); 373*139613f2SLeila Ghaffari 374*139613f2SLeila Ghaffari // --Stabilization terms 375*139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 376*139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 377*139613f2SLeila Ghaffari computeFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 378*139613f2SLeila Ghaffari 379*139613f2SLeila Ghaffari // ---- Transpose of the Jacobian 380*139613f2SLeila Ghaffari CeedScalar jacob_F_conv_T[3][5][5]; 381*139613f2SLeila Ghaffari for (int j=0; j<3; j++) 382*139613f2SLeila Ghaffari for (int k=0; k<5; k++) 383*139613f2SLeila Ghaffari for (int l=0; l<5; l++) 384*139613f2SLeila Ghaffari jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k]; 385*139613f2SLeila Ghaffari 386*139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 387*139613f2SLeila Ghaffari CeedScalar dqdx[5][3]; 388*139613f2SLeila Ghaffari for (int j=0; j<3; j++) { 389*139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j]; 390*139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j]; 391*139613f2SLeila Ghaffari for (int k=0; k<3; k++) 392*139613f2SLeila Ghaffari dqdx[k+1][j] = dUdx[k][j]; 393*139613f2SLeila Ghaffari } 394*139613f2SLeila Ghaffari 395*139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 396*139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 397*139613f2SLeila Ghaffari for (int j=0; j<3; j++) 398*139613f2SLeila Ghaffari for (int k=0; k<5; k++) 399*139613f2SLeila Ghaffari for (int l=0; l<5; l++) 400*139613f2SLeila Ghaffari strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 401*139613f2SLeila Ghaffari 402*139613f2SLeila Ghaffari // ---- Tau elements 403*139613f2SLeila Ghaffari CeedScalar uX[3]; 404*139613f2SLeila Ghaffari for (int j=0; j<3; j++) 405*139613f2SLeila Ghaffari uX[j] = dXdx[j][0]*u[0] + dXdx[j][1]*u[1] + dXdx[j][2]*u[2]; 406*139613f2SLeila Ghaffari const CeedScalar uiujgij = uX[0]*uX[0] + uX[1]*uX[1] + uX[2]*uX[2]; 407*139613f2SLeila Ghaffari const CeedScalar Cc = 1.; 408*139613f2SLeila Ghaffari const CeedScalar Ce = 1.; 409*139613f2SLeila Ghaffari const CeedScalar f1 = rho * sqrt(uiujgij); 410*139613f2SLeila Ghaffari const CeedScalar TauC = (Cc * f1) / 411*139613f2SLeila Ghaffari (8. * (dXdxdXdxT[0][0] + dXdxdXdxT[1][1] + dXdxdXdxT[2][2])); 412*139613f2SLeila Ghaffari const CeedScalar TauM = 1. / (f1>1. ? f1 : 1.); 413*139613f2SLeila Ghaffari const CeedScalar TauE = TauM / (Ce * cv); 414*139613f2SLeila Ghaffari const CeedScalar Tau[5] = {TauC, TauM, TauM, TauM, TauE}; 415*139613f2SLeila Ghaffari 416*139613f2SLeila Ghaffari // ---- Stabilization 417*139613f2SLeila Ghaffari CeedScalar stab[5][3]; 418*139613f2SLeila Ghaffari switch (context->stabilization) { 419*139613f2SLeila Ghaffari case 0: // Galerkin 420*139613f2SLeila Ghaffari break; 421*139613f2SLeila Ghaffari case 1: // SU 422*139613f2SLeila Ghaffari for (int j=0; j<3; j++) 423*139613f2SLeila Ghaffari for (int k=0; k<5; k++) 424*139613f2SLeila Ghaffari for (int l=0; l<5; l++) 425*139613f2SLeila Ghaffari stab[k][j] = jacob_F_conv_T[j][k][l] * Tau[l] * strong_conv[l]; 426*139613f2SLeila Ghaffari 427*139613f2SLeila Ghaffari for (int j=0; j<5; j++) 428*139613f2SLeila Ghaffari for (int k=0; k<3; k++) 429*139613f2SLeila Ghaffari dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] + 430*139613f2SLeila Ghaffari stab[j][1] * dXdx[k][1] + 431*139613f2SLeila Ghaffari stab[j][2] * dXdx[k][2]); 432*139613f2SLeila Ghaffari break; 433*139613f2SLeila Ghaffari case 2: // SUPG is not implemented for explicit scheme 434*139613f2SLeila Ghaffari break; 435*139613f2SLeila Ghaffari } 436*139613f2SLeila Ghaffari 437a515125bSLeila Ghaffari } // End Quadrature Point Loop 438a515125bSLeila Ghaffari 439a515125bSLeila Ghaffari // Return 440a515125bSLeila Ghaffari return 0; 441a515125bSLeila Ghaffari } 442a515125bSLeila Ghaffari 443a515125bSLeila Ghaffari // ***************************************************************************** 444a515125bSLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above) 445a515125bSLeila Ghaffari // with implicit time stepping method 446a515125bSLeila Ghaffari // 447a515125bSLeila Ghaffari // ***************************************************************************** 448a515125bSLeila Ghaffari CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, 449a515125bSLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 450a515125bSLeila Ghaffari // *INDENT-OFF* 451a515125bSLeila Ghaffari // Inputs 452a515125bSLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 453*139613f2SLeila Ghaffari (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1], 454a515125bSLeila Ghaffari (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2], 455a515125bSLeila Ghaffari (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 456a515125bSLeila Ghaffari // Outputs 457a515125bSLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0], 458a515125bSLeila Ghaffari (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1]; 459a515125bSLeila Ghaffari 460*139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx; 461a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 462*139613f2SLeila Ghaffari const CeedScalar cv = 2.5; 463a515125bSLeila Ghaffari 464a515125bSLeila Ghaffari CeedPragmaSIMD 465a515125bSLeila Ghaffari // Quadrature Point Loop 466a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 467a515125bSLeila Ghaffari // *INDENT-OFF* 468a515125bSLeila Ghaffari // Setup 469a515125bSLeila Ghaffari // -- Interp in 470a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 471a515125bSLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 472a515125bSLeila Ghaffari q[2][i] / rho, 473a515125bSLeila Ghaffari q[3][i] / rho 474a515125bSLeila Ghaffari }; 475a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 476*139613f2SLeila Ghaffari const CeedScalar drho[3] = {dq[0][0][i], 477*139613f2SLeila Ghaffari dq[1][0][i], 478*139613f2SLeila Ghaffari dq[2][0][i] 479*139613f2SLeila Ghaffari }; 480*139613f2SLeila Ghaffari const CeedScalar dU[3][3] = {{dq[0][1][i], 481*139613f2SLeila Ghaffari dq[1][1][i], 482*139613f2SLeila Ghaffari dq[2][1][i]}, 483*139613f2SLeila Ghaffari {dq[0][2][i], 484*139613f2SLeila Ghaffari dq[1][2][i], 485*139613f2SLeila Ghaffari dq[2][2][i]}, 486*139613f2SLeila Ghaffari {dq[0][3][i], 487*139613f2SLeila Ghaffari dq[1][3][i], 488*139613f2SLeila Ghaffari dq[2][3][i]} 489*139613f2SLeila Ghaffari }; 490*139613f2SLeila Ghaffari const CeedScalar dE[3] = {dq[0][4][i], 491*139613f2SLeila Ghaffari dq[1][4][i], 492*139613f2SLeila Ghaffari dq[2][4][i] 493*139613f2SLeila Ghaffari }; 494a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 495a515125bSLeila Ghaffari const CeedScalar wdetJ = q_data[0][i]; 496a515125bSLeila Ghaffari // -- Interp-to-Grad q_data 497a515125bSLeila Ghaffari // ---- Inverse of change of coordinate matrix: X_i,j 498a515125bSLeila Ghaffari // *INDENT-OFF* 499a515125bSLeila Ghaffari const CeedScalar dXdx[3][3] = {{q_data[1][i], 500a515125bSLeila Ghaffari q_data[2][i], 501a515125bSLeila Ghaffari q_data[3][i]}, 502a515125bSLeila Ghaffari {q_data[4][i], 503a515125bSLeila Ghaffari q_data[5][i], 504a515125bSLeila Ghaffari q_data[6][i]}, 505a515125bSLeila Ghaffari {q_data[7][i], 506a515125bSLeila Ghaffari q_data[8][i], 507a515125bSLeila Ghaffari q_data[9][i]} 508a515125bSLeila Ghaffari }; 509a515125bSLeila Ghaffari // *INDENT-ON* 510*139613f2SLeila Ghaffari // dU/dx 511*139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.}; 512*139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.}; 513*139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}}; 514*139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}}; 515*139613f2SLeila Ghaffari for (int j=0; j<3; j++) { 516*139613f2SLeila Ghaffari for (int k=0; k<3; k++) { 517*139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j]; 518*139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j]; 519*139613f2SLeila Ghaffari for (int l=0; l<3; l++) { 520*139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k]; 521*139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l]; //dXdx_j,k * dXdx_k,j 522*139613f2SLeila Ghaffari } 523*139613f2SLeila Ghaffari } 524*139613f2SLeila Ghaffari } 525a515125bSLeila Ghaffari const CeedScalar 526a515125bSLeila Ghaffari E_kinetic = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]), 527a515125bSLeila Ghaffari E_internal = E - E_kinetic, 528*139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure 529a515125bSLeila Ghaffari 530a515125bSLeila Ghaffari // The Physics 531a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it 532a515125bSLeila Ghaffari for (int j=0; j<5; j++) { 533*139613f2SLeila Ghaffari v[j][i] = 0.; 534a515125bSLeila Ghaffari for (int k=0; k<3; k++) 535*139613f2SLeila Ghaffari dv[k][j][i] = 0.; 536a515125bSLeila Ghaffari } 537a515125bSLeila Ghaffari //-----mass matrix 538a515125bSLeila Ghaffari for (int j=0; j<5; j++) 539a515125bSLeila Ghaffari v[j][i] += wdetJ*q_dot[j][i]; 540a515125bSLeila Ghaffari 541a515125bSLeila Ghaffari // -- Density 542a515125bSLeila Ghaffari // ---- u rho 543a515125bSLeila Ghaffari for (int j=0; j<3; j++) 544a515125bSLeila Ghaffari dv[j][0][i] -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] + 545a515125bSLeila Ghaffari rho*u[2]*dXdx[j][2]); 546a515125bSLeila Ghaffari // -- Momentum 547a515125bSLeila Ghaffari // ---- rho (u x u) + P I3 548a515125bSLeila Ghaffari for (int j=0; j<3; j++) 549a515125bSLeila Ghaffari for (int k=0; k<3; k++) 550*139613f2SLeila Ghaffari dv[k][j+1][i] -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] + 551*139613f2SLeila Ghaffari (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] + 552*139613f2SLeila Ghaffari (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]); 553a515125bSLeila Ghaffari // -- Total Energy Density 554a515125bSLeila Ghaffari // ---- (E + P) u 555a515125bSLeila Ghaffari for (int j=0; j<3; j++) 556a515125bSLeila Ghaffari dv[j][4][i] -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] + 557a515125bSLeila Ghaffari u[2]*dXdx[j][2]); 558*139613f2SLeila Ghaffari 559*139613f2SLeila Ghaffari // -- Stabilization terms 560*139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction 561*139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}}; 562*139613f2SLeila Ghaffari computeFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma); 563*139613f2SLeila Ghaffari 564*139613f2SLeila Ghaffari // ---- Transpose of the Jacobian 565*139613f2SLeila Ghaffari CeedScalar jacob_F_conv_T[3][5][5]; 566*139613f2SLeila Ghaffari for (int j=0; j<3; j++) 567*139613f2SLeila Ghaffari for (int k=0; k<5; k++) 568*139613f2SLeila Ghaffari for (int l=0; l<5; l++) 569*139613f2SLeila Ghaffari jacob_F_conv_T[j][k][l] = jacob_F_conv[j][l][k]; 570*139613f2SLeila Ghaffari 571*139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector 572*139613f2SLeila Ghaffari CeedScalar dqdx[5][3]; 573*139613f2SLeila Ghaffari for (int j=0; j<3; j++) { 574*139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j]; 575*139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j]; 576*139613f2SLeila Ghaffari for (int k=0; k<3; k++) 577*139613f2SLeila Ghaffari dqdx[k+1][j] = dUdx[k][j]; 578*139613f2SLeila Ghaffari } 579*139613f2SLeila Ghaffari 580*139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection) 581*139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.}; 582*139613f2SLeila Ghaffari for (int j=0; j<3; j++) 583*139613f2SLeila Ghaffari for (int k=0; k<5; k++) 584*139613f2SLeila Ghaffari for (int l=0; l<5; l++) 585*139613f2SLeila Ghaffari strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j]; 586*139613f2SLeila Ghaffari 587*139613f2SLeila Ghaffari // ---- Strong residual 588*139613f2SLeila Ghaffari CeedScalar strong_res[5]; 589*139613f2SLeila Ghaffari for (int j=0; j<5; j++) 590*139613f2SLeila Ghaffari strong_res[j] = q_dot[j][i] + strong_conv[j]; 591*139613f2SLeila Ghaffari 592*139613f2SLeila Ghaffari // ---- Tau elements 593*139613f2SLeila Ghaffari CeedScalar uX[3]; 594*139613f2SLeila Ghaffari for (int j=0; j<3; j++) 595*139613f2SLeila Ghaffari uX[j] = dXdx[j][0]*u[0] + dXdx[j][1]*u[1] + dXdx[j][2]*u[2]; 596*139613f2SLeila Ghaffari const CeedScalar uiujgij = uX[0]*uX[0] + uX[1]*uX[1] + uX[2]*uX[2]; 597*139613f2SLeila Ghaffari const CeedScalar Cc = 1.; 598*139613f2SLeila Ghaffari const CeedScalar Ce = 1.; 599*139613f2SLeila Ghaffari const CeedScalar f1 = rho * sqrt(uiujgij); 600*139613f2SLeila Ghaffari const CeedScalar TauC = (Cc * f1) / 601*139613f2SLeila Ghaffari (8. * (dXdxdXdxT[0][0] + dXdxdXdxT[1][1] + dXdxdXdxT[2][2])); 602*139613f2SLeila Ghaffari const CeedScalar TauM = 1. / (f1>1. ? f1 : 1.); 603*139613f2SLeila Ghaffari const CeedScalar TauE = TauM / (Ce * cv); 604*139613f2SLeila Ghaffari const CeedScalar Tau[5] = {TauC, TauM, TauM, TauM, TauE}; 605*139613f2SLeila Ghaffari 606*139613f2SLeila Ghaffari // ---- Stabilization 607*139613f2SLeila Ghaffari CeedScalar stab[5][3]; 608*139613f2SLeila Ghaffari switch (context->stabilization) { 609*139613f2SLeila Ghaffari case 0: // Galerkin 610*139613f2SLeila Ghaffari break; 611*139613f2SLeila Ghaffari case 1: // SU 612*139613f2SLeila Ghaffari for (int j=0; j<3; j++) 613*139613f2SLeila Ghaffari for (int k=0; k<5; k++) 614*139613f2SLeila Ghaffari for (int l=0; l<5; l++) 615*139613f2SLeila Ghaffari stab[k][j] = jacob_F_conv_T[j][k][l] * Tau[l] * strong_conv[l]; 616*139613f2SLeila Ghaffari 617*139613f2SLeila Ghaffari for (int j=0; j<5; j++) 618*139613f2SLeila Ghaffari for (int k=0; k<3; k++) 619*139613f2SLeila Ghaffari dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 620*139613f2SLeila Ghaffari stab[j][1] * dXdx[k][1] + 621*139613f2SLeila Ghaffari stab[j][2] * dXdx[k][2]); 622*139613f2SLeila Ghaffari break; 623*139613f2SLeila Ghaffari case 2: // SUPG 624*139613f2SLeila Ghaffari for (int j=0; j<3; j++) 625*139613f2SLeila Ghaffari for (int k=0; k<5; k++) 626*139613f2SLeila Ghaffari for (int l=0; l<5; l++) 627*139613f2SLeila Ghaffari stab[k][j] = jacob_F_conv_T[j][k][l] * Tau[l] * strong_res[l]; 628*139613f2SLeila Ghaffari 629*139613f2SLeila Ghaffari for (int j=0; j<5; j++) 630*139613f2SLeila Ghaffari for (int k=0; k<3; k++) 631*139613f2SLeila Ghaffari dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] + 632*139613f2SLeila Ghaffari stab[j][1] * dXdx[k][1] + 633*139613f2SLeila Ghaffari stab[j][2] * dXdx[k][2]); 634*139613f2SLeila Ghaffari break; 635*139613f2SLeila Ghaffari } 636a515125bSLeila Ghaffari } // End Quadrature Point Loop 637a515125bSLeila Ghaffari 638a515125bSLeila Ghaffari // Return 639a515125bSLeila Ghaffari return 0; 640a515125bSLeila Ghaffari } 641a515125bSLeila Ghaffari // ***************************************************************************** 642a515125bSLeila Ghaffari // This QFunction sets the boundary conditions 643a515125bSLeila Ghaffari // In this problem, only in/outflow BCs are implemented 644a515125bSLeila Ghaffari // 645a515125bSLeila Ghaffari // Inflow and outflow faces are determined based on 646a515125bSLeila Ghaffari // sign(dot(mean_velocity, normal)): 647a515125bSLeila Ghaffari // sign(dot(mean_velocity, normal)) > 0 : outflow BCs 648a515125bSLeila Ghaffari // sign(dot(mean_velocity, normal)) < 0 : inflow BCs 649a515125bSLeila Ghaffari // 650a515125bSLeila Ghaffari // Outflow BCs: 651a515125bSLeila Ghaffari // The validity of the weak form of the governing equations is 652a515125bSLeila Ghaffari // extended to the outflow. 653a515125bSLeila Ghaffari // 654a515125bSLeila Ghaffari // Inflow BCs: 655a515125bSLeila Ghaffari // Prescribed T_inlet and P_inlet are converted to conservative variables 656a515125bSLeila Ghaffari // and applied weakly. 657a515125bSLeila Ghaffari // 658a515125bSLeila Ghaffari // ***************************************************************************** 659a515125bSLeila Ghaffari CEED_QFUNCTION(Euler_Sur)(void *ctx, CeedInt Q, 660a515125bSLeila Ghaffari const CeedScalar *const *in, 661a515125bSLeila Ghaffari CeedScalar *const *out) { 662a515125bSLeila Ghaffari // *INDENT-OFF* 663a515125bSLeila Ghaffari // Inputs 664a515125bSLeila Ghaffari const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 665a515125bSLeila Ghaffari (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 666a515125bSLeila Ghaffari // Outputs 667a515125bSLeila Ghaffari CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 668a515125bSLeila Ghaffari // *INDENT-ON* 669a515125bSLeila Ghaffari EulerContext context = (EulerContext)ctx; 670a515125bSLeila Ghaffari const int euler_test = context->euler_test; 671a515125bSLeila Ghaffari const bool implicit = context->implicit; 672a515125bSLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity; 673a515125bSLeila Ghaffari 674a515125bSLeila Ghaffari const CeedScalar gamma = 1.4; 675a515125bSLeila Ghaffari const CeedScalar cv = 2.5; 676a515125bSLeila Ghaffari const CeedScalar R = 1.; 677a515125bSLeila Ghaffari CeedScalar T_inlet; 678a515125bSLeila Ghaffari CeedScalar P_inlet; 679a515125bSLeila Ghaffari 680a515125bSLeila Ghaffari // For test cases 1 and 3 the background velocity is zero 681a515125bSLeila Ghaffari if (euler_test == 1 || euler_test == 3) 682a515125bSLeila Ghaffari for (CeedInt i=0; i<3; i++) mean_velocity[i] = 0.; 683a515125bSLeila Ghaffari 684a515125bSLeila Ghaffari // For test cases 1 and 2, T_inlet = T_inlet = 0.4 685a515125bSLeila Ghaffari if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4; 686a515125bSLeila Ghaffari else T_inlet = P_inlet = 1.; 687a515125bSLeila Ghaffari 688a515125bSLeila Ghaffari CeedPragmaSIMD 689a515125bSLeila Ghaffari // Quadrature Point Loop 690a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 691a515125bSLeila Ghaffari // Setup 692a515125bSLeila Ghaffari // -- Interp in 693a515125bSLeila Ghaffari const CeedScalar rho = q[0][i]; 694a515125bSLeila Ghaffari const CeedScalar u[3] = {q[1][i] / rho, 695a515125bSLeila Ghaffari q[2][i] / rho, 696a515125bSLeila Ghaffari q[3][i] / rho 697a515125bSLeila Ghaffari }; 698a515125bSLeila Ghaffari const CeedScalar E = q[4][i]; 699a515125bSLeila Ghaffari 700a515125bSLeila Ghaffari // -- Interp-to-Interp q_data 701a515125bSLeila Ghaffari // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q). 702a515125bSLeila Ghaffari // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q). 703a515125bSLeila Ghaffari // We can effect this by swapping the sign on this weight 704a515125bSLeila Ghaffari const CeedScalar wdetJb = (implicit ? -1. : 1.) * q_data_sur[0][i]; 705a515125bSLeila Ghaffari // ---- Normal vectors 706a515125bSLeila Ghaffari const CeedScalar norm[3] = {q_data_sur[1][i], 707a515125bSLeila Ghaffari q_data_sur[2][i], 708a515125bSLeila Ghaffari q_data_sur[3][i] 709a515125bSLeila Ghaffari }; 710a515125bSLeila Ghaffari 711a515125bSLeila Ghaffari // face_normal = Normal vector of the face 712a515125bSLeila Ghaffari const CeedScalar face_normal = norm[0]*mean_velocity[0] + 713a515125bSLeila Ghaffari norm[1]*mean_velocity[1] + 714a515125bSLeila Ghaffari norm[2]*mean_velocity[2]; 715a515125bSLeila Ghaffari // The Physics 716a515125bSLeila Ghaffari // Zero v so all future terms can safely sum into it 717*139613f2SLeila Ghaffari for (int j=0; j<5; j++) v[j][i] = 0.; 718a515125bSLeila Ghaffari 719a515125bSLeila Ghaffari // Implementing in/outflow BCs 720a515125bSLeila Ghaffari if (face_normal > 0) { // outflow 721a515125bSLeila Ghaffari const CeedScalar E_kinetic = (u[0]*u[0] + u[1]*u[1]) / 2.; 722a515125bSLeila Ghaffari const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure 723a515125bSLeila Ghaffari const CeedScalar u_normal = norm[0]*u[0] + norm[1]*u[1] + 724a515125bSLeila Ghaffari norm[2]*u[2]; // Normal velocity 725a515125bSLeila Ghaffari // The Physics 726a515125bSLeila Ghaffari // -- Density 727a515125bSLeila Ghaffari v[0][i] -= wdetJb * rho * u_normal; 728a515125bSLeila Ghaffari 729a515125bSLeila Ghaffari // -- Momentum 730a515125bSLeila Ghaffari for (int j=0; j<3; j++) 731a515125bSLeila Ghaffari v[j+1][i] -= wdetJb *(rho * u_normal * u[j] + norm[j] * P); 732a515125bSLeila Ghaffari 733a515125bSLeila Ghaffari // -- Total Energy Density 734a515125bSLeila Ghaffari v[4][i] -= wdetJb * u_normal * (E + P); 735a515125bSLeila Ghaffari 736a515125bSLeila Ghaffari } else { // inflow 737a515125bSLeila Ghaffari const CeedScalar rho_inlet = P_inlet/(R*T_inlet); 738a515125bSLeila Ghaffari const CeedScalar E_kinetic_inlet = (mean_velocity[0]*mean_velocity[0] + 739a515125bSLeila Ghaffari mean_velocity[1]*mean_velocity[1]) / 2.; 740a515125bSLeila Ghaffari // incoming total energy 741a515125bSLeila Ghaffari const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet); 742a515125bSLeila Ghaffari 743a515125bSLeila Ghaffari // The Physics 744a515125bSLeila Ghaffari // -- Density 745a515125bSLeila Ghaffari v[0][i] -= wdetJb * rho_inlet * face_normal; 746a515125bSLeila Ghaffari 747a515125bSLeila Ghaffari // -- Momentum 748a515125bSLeila Ghaffari for (int j=0; j<3; j++) 749a515125bSLeila Ghaffari v[j+1][i] -= wdetJb *(rho_inlet * face_normal * mean_velocity[j] + 750a515125bSLeila Ghaffari norm[j] * P_inlet); 751a515125bSLeila Ghaffari 752a515125bSLeila Ghaffari // -- Total Energy Density 753a515125bSLeila Ghaffari v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet); 754a515125bSLeila Ghaffari } 755a515125bSLeila Ghaffari 756a515125bSLeila Ghaffari } // End Quadrature Point Loop 757a515125bSLeila Ghaffari return 0; 758a515125bSLeila Ghaffari } 759a515125bSLeila Ghaffari 760a515125bSLeila Ghaffari // ***************************************************************************** 761a515125bSLeila Ghaffari 762a515125bSLeila Ghaffari #endif // eulervortex_h 763