1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari
4a515125bSLeila Ghaffari /// @file
5*ea615d4cSJames Wright /// Euler traveling vortex initial condition and operator for HONEE
6a515125bSLeila Ghaffari
7a515125bSLeila Ghaffari // Model from:
804e40bb6SJeremy L Thompson // On the Order of Accuracy and Numerical Performance of Two Classes of Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011).
93e17a7a1SJames Wright #include <ceed/types.h>
103e17a7a1SJames Wright #ifndef CEED_RUNNING_JIT_PASS
113e17a7a1SJames Wright #include <stdbool.h>
123e17a7a1SJames Wright #endif
132b916ea7SJeremy L Thompson
14704b8bbeSJames Wright #include "utils.h"
15a515125bSLeila Ghaffari
16a515125bSLeila Ghaffari typedef struct EulerContext_ *EulerContext;
17a515125bSLeila Ghaffari struct EulerContext_ {
18a515125bSLeila Ghaffari CeedScalar center[3];
19a515125bSLeila Ghaffari CeedScalar curr_time;
20a515125bSLeila Ghaffari CeedScalar vortex_strength;
21d8a22b9eSJed Brown CeedScalar c_tau;
22a515125bSLeila Ghaffari CeedScalar mean_velocity[3];
23a515125bSLeila Ghaffari bool implicit;
24139613f2SLeila Ghaffari int euler_test;
25139613f2SLeila Ghaffari int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG
26a515125bSLeila Ghaffari };
27a515125bSLeila Ghaffari
28a515125bSLeila Ghaffari // *****************************************************************************
29a515125bSLeila Ghaffari // This function sets the initial conditions
30a515125bSLeila Ghaffari //
31a515125bSLeila Ghaffari // Temperature:
32a515125bSLeila Ghaffari // T = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2)
33a515125bSLeila Ghaffari // Density:
34a515125bSLeila Ghaffari // rho = (T/S_vortex)^(1 / (gamma - 1))
35a515125bSLeila Ghaffari // Pressure:
36a515125bSLeila Ghaffari // P = rho * T
37a515125bSLeila Ghaffari // Velocity:
38a515125bSLeila Ghaffari // ui = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi)
39a515125bSLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 )
40a515125bSLeila Ghaffari // Velocity/Momentum Density:
41a515125bSLeila Ghaffari // Ui = rho ui
42a515125bSLeila Ghaffari // Total Energy:
43a515125bSLeila Ghaffari // E = P / (gamma - 1) + rho (u u)/2
44a515125bSLeila Ghaffari //
45a515125bSLeila Ghaffari // Constants:
46a515125bSLeila Ghaffari // cv , Specific heat, constant volume
47a515125bSLeila Ghaffari // cp , Specific heat, constant pressure
48a515125bSLeila Ghaffari // vortex_strength , Strength of vortex
49a515125bSLeila Ghaffari // center , Location of bubble center
50a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio
51a515125bSLeila Ghaffari //
52a515125bSLeila Ghaffari // *****************************************************************************
53a515125bSLeila Ghaffari
54a515125bSLeila Ghaffari // *****************************************************************************
5504e40bb6SJeremy L Thompson // This helper function provides support for the exact, time-dependent solution (currently not implemented) and IC formulation for Euler traveling
5604e40bb6SJeremy L Thompson // vortex
57a515125bSLeila Ghaffari // *****************************************************************************
Exact_Euler(CeedInt dim,CeedScalar time,const CeedScalar X[],CeedInt Nf,CeedScalar q[],void * ctx)582b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
59a515125bSLeila Ghaffari // Context
60a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx;
61a515125bSLeila Ghaffari const CeedScalar vortex_strength = context->vortex_strength;
62a515125bSLeila Ghaffari const CeedScalar *center = context->center; // Center of the domain
63a515125bSLeila Ghaffari const CeedScalar *mean_velocity = context->mean_velocity;
64a515125bSLeila Ghaffari
65a515125bSLeila Ghaffari // Setup
66a515125bSLeila Ghaffari const CeedScalar gamma = 1.4;
67a515125bSLeila Ghaffari const CeedScalar cv = 2.5;
68a515125bSLeila Ghaffari const CeedScalar R = 1.;
69a515125bSLeila Ghaffari const CeedScalar x = X[0], y = X[1]; // Coordinates
70a515125bSLeila Ghaffari // Vortex center
71a515125bSLeila Ghaffari const CeedScalar xc = center[0] + mean_velocity[0] * time;
72a515125bSLeila Ghaffari const CeedScalar yc = center[1] + mean_velocity[1] * time;
73a515125bSLeila Ghaffari
74a515125bSLeila Ghaffari const CeedScalar x0 = x - xc;
75a515125bSLeila Ghaffari const CeedScalar y0 = y - yc;
7674960ff5SJames Wright const CeedScalar r = sqrt(Square(x0) + Square(y0));
77a515125bSLeila Ghaffari const CeedScalar C = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI);
782b916ea7SJeremy L Thompson const CeedScalar delta_T = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI);
79a515125bSLeila Ghaffari const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma
802b916ea7SJeremy L Thompson const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI);
81a515125bSLeila Ghaffari CeedScalar rho, P, T, E, u[3] = {0.};
82a515125bSLeila Ghaffari
83a515125bSLeila Ghaffari // Initial Conditions
84a515125bSLeila Ghaffari switch (context->euler_test) {
85a515125bSLeila Ghaffari case 0: // Traveling vortex
86a515125bSLeila Ghaffari T = 1 + delta_T;
87a515125bSLeila Ghaffari // P = rho * T
88a515125bSLeila Ghaffari // P = S * rho^gamma
89a515125bSLeila Ghaffari // Solve for rho, then substitute for P
90139613f2SLeila Ghaffari rho = pow(T / S_vortex, 1 / (gamma - 1.));
91a515125bSLeila Ghaffari P = rho * T;
92a515125bSLeila Ghaffari u[0] = mean_velocity[0] - C * y0;
93a515125bSLeila Ghaffari u[1] = mean_velocity[1] + C * x0;
94a515125bSLeila Ghaffari
95a515125bSLeila Ghaffari // Assign exact solution
96a515125bSLeila Ghaffari q[0] = rho;
97a515125bSLeila Ghaffari q[1] = rho * u[0];
98a515125bSLeila Ghaffari q[2] = rho * u[1];
99a515125bSLeila Ghaffari q[3] = rho * u[2];
100a515125bSLeila Ghaffari q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.;
101a515125bSLeila Ghaffari break;
102a515125bSLeila Ghaffari case 1: // Constant zero velocity, density constant, total energy constant
103a515125bSLeila Ghaffari rho = 1.;
104a515125bSLeila Ghaffari E = 2.;
105a515125bSLeila Ghaffari
106a515125bSLeila Ghaffari // Assign exact solution
107a515125bSLeila Ghaffari q[0] = rho;
108a515125bSLeila Ghaffari q[1] = rho * u[0];
109a515125bSLeila Ghaffari q[2] = rho * u[1];
110a515125bSLeila Ghaffari q[3] = rho * u[2];
111a515125bSLeila Ghaffari q[4] = E;
112a515125bSLeila Ghaffari break;
113a515125bSLeila Ghaffari case 2: // Constant nonzero velocity, density constant, total energy constant
114a515125bSLeila Ghaffari rho = 1.;
115a515125bSLeila Ghaffari E = 2.;
116a515125bSLeila Ghaffari u[0] = mean_velocity[0];
117a515125bSLeila Ghaffari u[1] = mean_velocity[1];
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] = E;
125a515125bSLeila Ghaffari break;
12604e40bb6SJeremy L Thompson case 3: // Velocity zero, pressure constant (so density and internal energy will be non-constant), but the velocity should stay zero and the
12704e40bb6SJeremy L Thompson // bubble won't diffuse
128a515125bSLeila Ghaffari // (for Euler, where there is no thermal conductivity)
129a515125bSLeila Ghaffari P = 1.;
130a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r * r);
131a515125bSLeila Ghaffari rho = P / (R * T);
132a515125bSLeila Ghaffari
133a515125bSLeila Ghaffari // Assign exact solution
134a515125bSLeila Ghaffari q[0] = rho;
135a515125bSLeila Ghaffari q[1] = rho * u[0];
136a515125bSLeila Ghaffari q[2] = rho * u[1];
137a515125bSLeila Ghaffari q[3] = rho * u[2];
138a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
139a515125bSLeila Ghaffari break;
14004e40bb6SJeremy L Thompson case 4: // Constant nonzero velocity, pressure constant (so density and internal energy will be non-constant),
14104e40bb6SJeremy L Thompson // It should be transported across the domain, but velocity stays constant
142a515125bSLeila Ghaffari P = 1.;
143a515125bSLeila Ghaffari T = 1. - S_bubble * exp(1. - r * r);
144a515125bSLeila Ghaffari rho = P / (R * T);
145a515125bSLeila Ghaffari u[0] = mean_velocity[0];
146a515125bSLeila Ghaffari u[1] = mean_velocity[1];
147a515125bSLeila Ghaffari
148a515125bSLeila Ghaffari // Assign exact solution
149a515125bSLeila Ghaffari q[0] = rho;
150a515125bSLeila Ghaffari q[1] = rho * u[0];
151a515125bSLeila Ghaffari q[2] = rho * u[1];
152a515125bSLeila Ghaffari q[3] = rho * u[2];
153a515125bSLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
154a515125bSLeila Ghaffari break;
1550df2634dSLeila Ghaffari case 5: // non-smooth thermal bubble - cylinder
1560df2634dSLeila Ghaffari P = 1.;
1570df2634dSLeila Ghaffari T = 1. - (r < 1. ? S_bubble : 0.);
1580df2634dSLeila Ghaffari rho = P / (R * T);
1590df2634dSLeila Ghaffari u[0] = mean_velocity[0];
1600df2634dSLeila Ghaffari u[1] = mean_velocity[1];
1610df2634dSLeila Ghaffari
1620df2634dSLeila Ghaffari // Assign exact solution
1630df2634dSLeila Ghaffari q[0] = rho;
1640df2634dSLeila Ghaffari q[1] = rho * u[0];
1650df2634dSLeila Ghaffari q[2] = rho * u[1];
1660df2634dSLeila Ghaffari q[3] = rho * u[2];
1670df2634dSLeila Ghaffari q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
1680df2634dSLeila Ghaffari break;
169a515125bSLeila Ghaffari }
170a515125bSLeila Ghaffari return 0;
171a515125bSLeila Ghaffari }
172a515125bSLeila Ghaffari
173a515125bSLeila Ghaffari // *****************************************************************************
174139613f2SLeila Ghaffari // Helper function for computing flux Jacobian
175139613f2SLeila Ghaffari // *****************************************************************************
ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5],const CeedScalar rho,const CeedScalar u[3],const CeedScalar E,const CeedScalar gamma)1762b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
177139613f2SLeila Ghaffari const CeedScalar gamma) {
178139613f2SLeila Ghaffari CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2]; // Velocity square
179139613f2SLeila Ghaffari for (CeedInt i = 0; i < 3; i++) { // Jacobian matrices for 3 directions
180139613f2SLeila Ghaffari for (CeedInt j = 0; j < 3; j++) { // Rows of each Jacobian matrix
181139613f2SLeila Ghaffari dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j];
182139613f2SLeila Ghaffari for (CeedInt k = 0; k < 3; k++) { // Columns of each Jacobian matrix
183139613f2SLeila Ghaffari dF[i][0][k + 1] = ((i == k) ? 1. : 0.);
1842b916ea7SJeremy L Thompson dF[i][j + 1][k + 1] = ((j == k) ? u[i] : 0.) + ((i == k) ? u[j] : 0.) - ((i == j) ? u[k] : 0.) * (gamma - 1.);
1852b916ea7SJeremy L Thompson dF[i][4][k + 1] = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k];
186139613f2SLeila Ghaffari }
187139613f2SLeila Ghaffari dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.);
188139613f2SLeila Ghaffari }
189139613f2SLeila Ghaffari dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho);
190139613f2SLeila Ghaffari dF[i][4][4] = u[i] * gamma;
191139613f2SLeila Ghaffari }
192139613f2SLeila Ghaffari }
193139613f2SLeila Ghaffari
194139613f2SLeila Ghaffari // *****************************************************************************
195d8a22b9eSJed Brown // Helper function for computing Tau elements (stabilization constant)
196d8a22b9eSJed Brown // Model from:
197d8a22b9eSJed Brown // Stabilized Methods for Compressible Flows, Hughes et al 2010
198d8a22b9eSJed Brown //
199d8a22b9eSJed Brown // Spatial criterion #2 - Tau is a 3x3 diagonal matrix
200d8a22b9eSJed Brown // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
201d8a22b9eSJed Brown //
202d8a22b9eSJed Brown // Where
203d8a22b9eSJed Brown // c_tau = stabilization constant (0.5 is reported as "optimal")
204d8a22b9eSJed Brown // h[i] = 2 length(dxdX[i])
205d8a22b9eSJed Brown // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
206d8a22b9eSJed Brown // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number )
20704e40bb6SJeremy L Thompson // rho(A[i]) = spectral radius of the convective flux Jacobian i, wave speed in direction i
208d8a22b9eSJed Brown // *****************************************************************************
Tau_spatial(CeedScalar Tau_x[3],const CeedScalar dXdx[3][3],const CeedScalar u[3],const CeedScalar sound_speed,const CeedScalar c_tau)2092b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], const CeedScalar dXdx[3][3], const CeedScalar u[3], const CeedScalar sound_speed,
2102b916ea7SJeremy L Thompson const CeedScalar c_tau) {
211493642f1SJames Wright for (CeedInt i = 0; i < 3; i++) {
212d8a22b9eSJed Brown // length of element in direction i
21374960ff5SJames Wright CeedScalar h = 2 / sqrt(Square(dXdx[0][i]) + Square(dXdx[1][i]) + Square(dXdx[2][i]));
214d8a22b9eSJed Brown // fastest wave in direction i
215d8a22b9eSJed Brown CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
216d8a22b9eSJed Brown Tau_x[i] = c_tau * h / fastest_wave;
217d8a22b9eSJed Brown }
218d8a22b9eSJed Brown }
219d8a22b9eSJed Brown
220d8a22b9eSJed Brown // *****************************************************************************
221a515125bSLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex
222a515125bSLeila Ghaffari // *****************************************************************************
ICsEuler(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)2232b916ea7SJeremy L Thompson CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
224a515125bSLeila Ghaffari const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
225a515125bSLeila Ghaffari CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
226b193fadcSJames Wright
227a515125bSLeila Ghaffari const EulerContext context = (EulerContext)ctx;
228a515125bSLeila Ghaffari
2293d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
230a515125bSLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
231139613f2SLeila Ghaffari CeedScalar q[5] = {0.};
232a515125bSLeila Ghaffari
233a515125bSLeila Ghaffari Exact_Euler(3, context->curr_time, x, 5, q, ctx);
2342b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
235b193fadcSJames Wright }
236a515125bSLeila Ghaffari return 0;
237a515125bSLeila Ghaffari }
238a515125bSLeila Ghaffari
239a515125bSLeila Ghaffari // *****************************************************************************
24004e40bb6SJeremy L Thompson // This QFunction implements the following formulation of Euler equations with explicit time stepping method
241a515125bSLeila Ghaffari //
24204e40bb6SJeremy L Thompson // This is 3D Euler for compressible gas dynamics in conservation form with state variables of density, momentum density, and total energy density.
243a515125bSLeila Ghaffari //
244a515125bSLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E )
245a515125bSLeila Ghaffari // rho - Mass Density
246a515125bSLeila Ghaffari // Ui - Momentum Density, Ui = rho ui
247a515125bSLeila Ghaffari // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2
248a515125bSLeila Ghaffari //
249a515125bSLeila Ghaffari // Euler Equations:
250a515125bSLeila Ghaffari // drho/dt + div( U ) = 0
251a515125bSLeila Ghaffari // dU/dt + div( rho (u x u) + P I3 ) = 0
252a515125bSLeila Ghaffari // dE/dt + div( (E + P) u ) = 0
253a515125bSLeila Ghaffari //
254a515125bSLeila Ghaffari // Equation of State:
255a515125bSLeila Ghaffari // P = (gamma - 1) (E - rho (u u) / 2)
256a515125bSLeila Ghaffari //
257a515125bSLeila Ghaffari // Constants:
258a515125bSLeila Ghaffari // cv , Specific heat, constant volume
259a515125bSLeila Ghaffari // cp , Specific heat, constant pressure
260a515125bSLeila Ghaffari // g , Gravity
261a515125bSLeila Ghaffari // gamma = cp / cv, Specific heat ratio
262a515125bSLeila Ghaffari // *****************************************************************************
Euler(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)2632b916ea7SJeremy L Thompson CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
2643d65b166SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
2653d65b166SJames Wright const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
266ade49511SJames Wright const CeedScalar(*q_data) = in[2];
2673d65b166SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
2683d65b166SJames Wright CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
269a515125bSLeila Ghaffari
270139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx;
271d8a22b9eSJed Brown const CeedScalar c_tau = context->c_tau;
272a515125bSLeila Ghaffari const CeedScalar gamma = 1.4;
273a515125bSLeila Ghaffari
2743d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
275a515125bSLeila Ghaffari // Setup
276a515125bSLeila Ghaffari // -- Interp in
277a515125bSLeila Ghaffari const CeedScalar rho = q[0][i];
2782b916ea7SJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
279a515125bSLeila Ghaffari const CeedScalar E = q[4][i];
2802b916ea7SJeremy L Thompson const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
2812b916ea7SJeremy L Thompson const CeedScalar dU[3][3] = {
2822b916ea7SJeremy L Thompson {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
2832b916ea7SJeremy L Thompson {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
2842b916ea7SJeremy L Thompson {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
285139613f2SLeila Ghaffari };
2862b916ea7SJeremy L Thompson const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
287ade49511SJames Wright CeedScalar wdetJ, dXdx[3][3];
288ade49511SJames Wright QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx);
289139613f2SLeila Ghaffari // dU/dx
290139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.};
291139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.};
292139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}};
293139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}};
294493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) {
295493642f1SJames Wright for (CeedInt k = 0; k < 3; k++) {
296139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j];
297139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j];
298493642f1SJames Wright for (CeedInt l = 0; l < 3; l++) {
299139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k];
300139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j
301139613f2SLeila Ghaffari }
302139613f2SLeila Ghaffari }
303139613f2SLeila Ghaffari }
304139613f2SLeila Ghaffari // Pressure
3052b916ea7SJeremy L Thompson const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic,
306139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure
307a515125bSLeila Ghaffari
308a515125bSLeila Ghaffari // The Physics
309a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it
310493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) {
311139613f2SLeila Ghaffari v[j][i] = 0.;
3122b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
313a515125bSLeila Ghaffari }
314a515125bSLeila Ghaffari
315a515125bSLeila Ghaffari // -- Density
316a515125bSLeila Ghaffari // ---- u rho
3172b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) dv[j][0][i] += wdetJ * (rho * u[0] * dXdx[j][0] + rho * u[1] * dXdx[j][1] + rho * u[2] * dXdx[j][2]);
318a515125bSLeila Ghaffari // -- Momentum
319a515125bSLeila Ghaffari // ---- rho (u x u) + P I3
3202b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) {
3212b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) {
3222b916ea7SJeremy L Thompson dv[k][j + 1][i] += wdetJ * ((rho * u[j] * u[0] + (j == 0 ? P : 0.)) * dXdx[k][0] + (rho * u[j] * u[1] + (j == 1 ? P : 0.)) * dXdx[k][1] +
323139613f2SLeila Ghaffari (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
3242b916ea7SJeremy L Thompson }
3252b916ea7SJeremy L Thompson }
326a515125bSLeila Ghaffari // -- Total Energy Density
327a515125bSLeila Ghaffari // ---- (E + P) u
3282b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) dv[j][4][i] += wdetJ * (E + P) * (u[0] * dXdx[j][0] + u[1] * dXdx[j][1] + u[2] * dXdx[j][2]);
329139613f2SLeila Ghaffari
330139613f2SLeila Ghaffari // --Stabilization terms
331139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
332139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
333d8a22b9eSJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
334139613f2SLeila Ghaffari
335139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector
336139613f2SLeila Ghaffari CeedScalar dqdx[5][3];
337493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) {
338139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j];
339139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j];
3402b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
341139613f2SLeila Ghaffari }
342139613f2SLeila Ghaffari
343139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection)
344139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.};
3452b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) {
3462b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) {
3472b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
3482b916ea7SJeremy L Thompson }
3492b916ea7SJeremy L Thompson }
350139613f2SLeila Ghaffari
351d8a22b9eSJed Brown // Stabilization
352d8a22b9eSJed Brown // -- Tau elements
353d8a22b9eSJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho);
354d8a22b9eSJed Brown CeedScalar Tau_x[3] = {0.};
355d8a22b9eSJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
356139613f2SLeila Ghaffari
357d8a22b9eSJed Brown // -- Stabilization method: none or SU
358bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}};
359139613f2SLeila Ghaffari switch (context->stabilization) {
360139613f2SLeila Ghaffari case 0: // Galerkin
361139613f2SLeila Ghaffari break;
362139613f2SLeila Ghaffari case 1: // SU
3632b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) {
3642b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) {
3652b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
3662b916ea7SJeremy L Thompson }
3672b916ea7SJeremy L Thompson }
368139613f2SLeila Ghaffari
3692b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) {
3702b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] -= wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]);
3712b916ea7SJeremy L Thompson }
372139613f2SLeila Ghaffari break;
373139613f2SLeila Ghaffari case 2: // SUPG is not implemented for explicit scheme
374139613f2SLeila Ghaffari break;
375139613f2SLeila Ghaffari }
376b193fadcSJames Wright }
377a515125bSLeila Ghaffari return 0;
378a515125bSLeila Ghaffari }
379a515125bSLeila Ghaffari
380a515125bSLeila Ghaffari // *****************************************************************************
38104e40bb6SJeremy L Thompson // This QFunction implements the Euler equations with (mentioned above) with implicit time stepping method
382a515125bSLeila Ghaffari // *****************************************************************************
IFunction_Euler(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)3832b916ea7SJeremy L Thompson CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
3843d65b166SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
3853d65b166SJames Wright const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
3863d65b166SJames Wright const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
387ade49511SJames Wright const CeedScalar(*q_data) = in[3];
3883d65b166SJames Wright CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
3893d65b166SJames Wright CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
390a515125bSLeila Ghaffari
391139613f2SLeila Ghaffari EulerContext context = (EulerContext)ctx;
392d8a22b9eSJed Brown const CeedScalar c_tau = context->c_tau;
393a515125bSLeila Ghaffari const CeedScalar gamma = 1.4;
394a515125bSLeila Ghaffari
3953d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
396a515125bSLeila Ghaffari // Setup
397a515125bSLeila Ghaffari // -- Interp in
398a515125bSLeila Ghaffari const CeedScalar rho = q[0][i];
3992b916ea7SJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
400a515125bSLeila Ghaffari const CeedScalar E = q[4][i];
4012b916ea7SJeremy L Thompson const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
4022b916ea7SJeremy L Thompson const CeedScalar dU[3][3] = {
4032b916ea7SJeremy L Thompson {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
4042b916ea7SJeremy L Thompson {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
4052b916ea7SJeremy L Thompson {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
406139613f2SLeila Ghaffari };
4072b916ea7SJeremy L Thompson const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
408ade49511SJames Wright CeedScalar wdetJ, dXdx[3][3];
409ade49511SJames Wright QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx);
410139613f2SLeila Ghaffari // dU/dx
411139613f2SLeila Ghaffari CeedScalar drhodx[3] = {0.};
412139613f2SLeila Ghaffari CeedScalar dEdx[3] = {0.};
413139613f2SLeila Ghaffari CeedScalar dUdx[3][3] = {{0.}};
414139613f2SLeila Ghaffari CeedScalar dXdxdXdxT[3][3] = {{0.}};
415493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) {
416493642f1SJames Wright for (CeedInt k = 0; k < 3; k++) {
417139613f2SLeila Ghaffari drhodx[j] += drho[k] * dXdx[k][j];
418139613f2SLeila Ghaffari dEdx[j] += dE[k] * dXdx[k][j];
419493642f1SJames Wright for (CeedInt l = 0; l < 3; l++) {
420139613f2SLeila Ghaffari dUdx[j][k] += dU[j][l] * dXdx[l][k];
421139613f2SLeila Ghaffari dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j
422139613f2SLeila Ghaffari }
423139613f2SLeila Ghaffari }
424139613f2SLeila Ghaffari }
4252b916ea7SJeremy L Thompson const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic,
426139613f2SLeila Ghaffari P = E_internal * (gamma - 1.); // P = pressure
427a515125bSLeila Ghaffari
428a515125bSLeila Ghaffari // The Physics
429a515125bSLeila Ghaffari // Zero v and dv so all future terms can safely sum into it
430493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) {
431139613f2SLeila Ghaffari v[j][i] = 0.;
4322b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
433a515125bSLeila Ghaffari }
434a515125bSLeila Ghaffari //-----mass matrix
4352b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i];
436a515125bSLeila Ghaffari
437a515125bSLeila Ghaffari // -- Density
438a515125bSLeila Ghaffari // ---- u rho
4392b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) dv[j][0][i] -= wdetJ * (rho * u[0] * dXdx[j][0] + rho * u[1] * dXdx[j][1] + rho * u[2] * dXdx[j][2]);
440a515125bSLeila Ghaffari // -- Momentum
441a515125bSLeila Ghaffari // ---- rho (u x u) + P I3
4422b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) {
4432b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) {
4442b916ea7SJeremy L Thompson dv[k][j + 1][i] -= wdetJ * ((rho * u[j] * u[0] + (j == 0 ? P : 0.)) * dXdx[k][0] + (rho * u[j] * u[1] + (j == 1 ? P : 0.)) * dXdx[k][1] +
445139613f2SLeila Ghaffari (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
4462b916ea7SJeremy L Thompson }
4472b916ea7SJeremy L Thompson }
448a515125bSLeila Ghaffari // -- Total Energy Density
449a515125bSLeila Ghaffari // ---- (E + P) u
4502b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) dv[j][4][i] -= wdetJ * (E + P) * (u[0] * dXdx[j][0] + u[1] * dXdx[j][1] + u[2] * dXdx[j][2]);
451139613f2SLeila Ghaffari
452139613f2SLeila Ghaffari // -- Stabilization terms
453139613f2SLeila Ghaffari // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
454139613f2SLeila Ghaffari CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
455d8a22b9eSJed Brown ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
456139613f2SLeila Ghaffari
457139613f2SLeila Ghaffari // ---- dqdx collects drhodx, dUdx and dEdx in one vector
458139613f2SLeila Ghaffari CeedScalar dqdx[5][3];
459493642f1SJames Wright for (CeedInt j = 0; j < 3; j++) {
460139613f2SLeila Ghaffari dqdx[0][j] = drhodx[j];
461139613f2SLeila Ghaffari dqdx[4][j] = dEdx[j];
4622b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
463139613f2SLeila Ghaffari }
464139613f2SLeila Ghaffari
465139613f2SLeila Ghaffari // ---- strong_conv = dF/dq * dq/dx (Strong convection)
466139613f2SLeila Ghaffari CeedScalar strong_conv[5] = {0.};
4672b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) {
4682b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) {
4692b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
4702b916ea7SJeremy L Thompson }
4712b916ea7SJeremy L Thompson }
472139613f2SLeila Ghaffari
473139613f2SLeila Ghaffari // ---- Strong residual
474139613f2SLeila Ghaffari CeedScalar strong_res[5];
4752b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j];
476139613f2SLeila Ghaffari
477d8a22b9eSJed Brown // Stabilization
478d8a22b9eSJed Brown // -- Tau elements
479d8a22b9eSJed Brown const CeedScalar sound_speed = sqrt(gamma * P / rho);
480d8a22b9eSJed Brown CeedScalar Tau_x[3] = {0.};
481d8a22b9eSJed Brown Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
482139613f2SLeila Ghaffari
483d8a22b9eSJed Brown // -- Stabilization method: none, SU, or SUPG
484bb8a0c61SJames Wright CeedScalar stab[5][3] = {{0.}};
485139613f2SLeila Ghaffari switch (context->stabilization) {
486139613f2SLeila Ghaffari case 0: // Galerkin
487139613f2SLeila Ghaffari break;
488139613f2SLeila Ghaffari case 1: // SU
4892b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) {
4902b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) {
4912b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
4922b916ea7SJeremy L Thompson }
4932b916ea7SJeremy L Thompson }
494139613f2SLeila Ghaffari
4952b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) {
4962b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] += wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]);
4972b916ea7SJeremy L Thompson }
498139613f2SLeila Ghaffari break;
499139613f2SLeila Ghaffari case 2: // SUPG
5002b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 3; j++) {
5012b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 5; k++) {
5022b916ea7SJeremy L Thompson for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l];
5032b916ea7SJeremy L Thompson }
5042b916ea7SJeremy L Thompson }
505139613f2SLeila Ghaffari
5062b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) {
5072b916ea7SJeremy L Thompson for (CeedInt k = 0; k < 3; k++) dv[k][j][i] += wdetJ * (stab[j][0] * dXdx[k][0] + stab[j][1] * dXdx[k][1] + stab[j][2] * dXdx[k][2]);
5082b916ea7SJeremy L Thompson }
509139613f2SLeila Ghaffari break;
510139613f2SLeila Ghaffari }
511b193fadcSJames Wright }
512a515125bSLeila Ghaffari return 0;
513a515125bSLeila Ghaffari }
514a515125bSLeila Ghaffari // *****************************************************************************
51504e40bb6SJeremy L Thompson // This QFunction sets the inflow boundary conditions for the traveling vortex problem.
516a515125bSLeila Ghaffari //
51704e40bb6SJeremy L Thompson // Prescribed T_inlet and P_inlet are converted to conservative variables and applied weakly.
518a515125bSLeila Ghaffari // *****************************************************************************
TravelingVortex_Inflow(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)5192b916ea7SJeremy L Thompson CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
520ade49511SJames Wright const CeedScalar(*q_data_sur) = in[2];
521a515125bSLeila Ghaffari CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
522b193fadcSJames Wright
523a515125bSLeila Ghaffari EulerContext context = (EulerContext)ctx;
524a515125bSLeila Ghaffari const int euler_test = context->euler_test;
525ade49511SJames Wright const bool is_implicit = context->implicit;
526a515125bSLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity;
527a515125bSLeila Ghaffari const CeedScalar cv = 2.5;
528a515125bSLeila Ghaffari const CeedScalar R = 1.;
529a515125bSLeila Ghaffari CeedScalar T_inlet;
530a515125bSLeila Ghaffari CeedScalar P_inlet;
531a515125bSLeila Ghaffari
532a515125bSLeila Ghaffari // For test cases 1 and 3 the background velocity is zero
5332b916ea7SJeremy L Thompson if (euler_test == 1 || euler_test == 3) {
534a515125bSLeila Ghaffari for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.;
5352b916ea7SJeremy L Thompson }
536a515125bSLeila Ghaffari
537a515125bSLeila Ghaffari // For test cases 1 and 2, T_inlet = T_inlet = 0.4
538a515125bSLeila Ghaffari if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
539a515125bSLeila Ghaffari else T_inlet = P_inlet = 1.;
540a515125bSLeila Ghaffari
5413d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
54278e8b7daSJames Wright CeedScalar wdetJb, normal[3];
54378e8b7daSJames Wright QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, normal);
544ade49511SJames Wright wdetJb *= is_implicit ? -1. : 1.;
545a515125bSLeila Ghaffari
546a515125bSLeila Ghaffari // face_normal = Normal vector of the face
54778e8b7daSJames Wright const CeedScalar face_normal = Dot3(normal, mean_velocity);
548a515125bSLeila Ghaffari // The Physics
549a515125bSLeila Ghaffari // Zero v so all future terms can safely sum into it
550493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.;
551a515125bSLeila Ghaffari
552a515125bSLeila Ghaffari // Implementing in/outflow BCs
553002797a3SLeila Ghaffari if (face_normal > 0) {
554a515125bSLeila Ghaffari } else { // inflow
555a515125bSLeila Ghaffari const CeedScalar rho_inlet = P_inlet / (R * T_inlet);
5562b916ea7SJeremy L Thompson const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.;
557a515125bSLeila Ghaffari // incoming total energy
558a515125bSLeila Ghaffari const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
559a515125bSLeila Ghaffari
560a515125bSLeila Ghaffari // The Physics
561a515125bSLeila Ghaffari // -- Density
562a515125bSLeila Ghaffari v[0][i] -= wdetJb * rho_inlet * face_normal;
563a515125bSLeila Ghaffari
564a515125bSLeila Ghaffari // -- Momentum
56578e8b7daSJames Wright for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho_inlet * face_normal * mean_velocity[j] + normal[j] * P_inlet);
566a515125bSLeila Ghaffari
567a515125bSLeila Ghaffari // -- Total Energy Density
568a515125bSLeila Ghaffari v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
569a515125bSLeila Ghaffari }
570b193fadcSJames Wright }
571a515125bSLeila Ghaffari return 0;
572a515125bSLeila Ghaffari }
573a515125bSLeila Ghaffari
574a515125bSLeila Ghaffari // *****************************************************************************
57504e40bb6SJeremy L Thompson // This QFunction sets the outflow boundary conditions for the Euler solver.
57668ef3d20SLeila Ghaffari //
57768ef3d20SLeila Ghaffari // Outflow BCs:
57804e40bb6SJeremy L Thompson // The validity of the weak form of the governing equations is extended to the outflow.
57968ef3d20SLeila Ghaffari // *****************************************************************************
Euler_Outflow(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)5802b916ea7SJeremy L Thompson CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
5813d65b166SJames Wright const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
582ade49511SJames Wright const CeedScalar(*q_data_sur) = in[2];
58368ef3d20SLeila Ghaffari CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
584b193fadcSJames Wright
58568ef3d20SLeila Ghaffari EulerContext context = (EulerContext)ctx;
586ade49511SJames Wright const bool is_implicit = context->implicit;
58768ef3d20SLeila Ghaffari CeedScalar *mean_velocity = context->mean_velocity;
58868ef3d20SLeila Ghaffari
58968ef3d20SLeila Ghaffari const CeedScalar gamma = 1.4;
59068ef3d20SLeila Ghaffari
5913d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
59268ef3d20SLeila Ghaffari // Setup
59368ef3d20SLeila Ghaffari // -- Interp in
59468ef3d20SLeila Ghaffari const CeedScalar rho = q[0][i];
5952b916ea7SJeremy L Thompson const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
59668ef3d20SLeila Ghaffari const CeedScalar E = q[4][i];
59768ef3d20SLeila Ghaffari
59878e8b7daSJames Wright CeedScalar wdetJb, normal[3];
59978e8b7daSJames Wright QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, normal);
600ade49511SJames Wright wdetJb *= is_implicit ? -1. : 1.;
60168ef3d20SLeila Ghaffari
60268ef3d20SLeila Ghaffari // face_normal = Normal vector of the face
60378e8b7daSJames Wright const CeedScalar face_normal = Dot3(normal, mean_velocity);
60468ef3d20SLeila Ghaffari // The Physics
60568ef3d20SLeila Ghaffari // Zero v so all future terms can safely sum into it
606493642f1SJames Wright for (CeedInt j = 0; j < 5; j++) v[j][i] = 0;
60768ef3d20SLeila Ghaffari
60868ef3d20SLeila Ghaffari // Implementing in/outflow BCs
60968ef3d20SLeila Ghaffari if (face_normal > 0) { // outflow
61068ef3d20SLeila Ghaffari const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.;
61168ef3d20SLeila Ghaffari const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure
61278e8b7daSJames Wright const CeedScalar u_normal = Dot3(normal, u); // Normal velocity
61368ef3d20SLeila Ghaffari // The Physics
61468ef3d20SLeila Ghaffari // -- Density
61568ef3d20SLeila Ghaffari v[0][i] -= wdetJb * rho * u_normal;
61668ef3d20SLeila Ghaffari
61768ef3d20SLeila Ghaffari // -- Momentum
61878e8b7daSJames Wright for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + normal[j] * P);
61968ef3d20SLeila Ghaffari
62068ef3d20SLeila Ghaffari // -- Total Energy Density
62168ef3d20SLeila Ghaffari v[4][i] -= wdetJb * u_normal * (E + P);
62268ef3d20SLeila Ghaffari }
623b193fadcSJames Wright }
62468ef3d20SLeila Ghaffari return 0;
62568ef3d20SLeila Ghaffari }
626