xref: /honee/qfunctions/eulervortex.h (revision bb8a0c61f21224cefcdd60e71004bb99df1e9a58)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Euler traveling vortex initial condition and operator for Navier-Stokes
10a515125bSLeila Ghaffari /// example using PETSc
11a515125bSLeila Ghaffari 
12a515125bSLeila Ghaffari // Model from:
13a515125bSLeila Ghaffari //   On the Order of Accuracy and Numerical Performance of Two Classes of
14a515125bSLeila Ghaffari //   Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011).
15a515125bSLeila Ghaffari 
16a515125bSLeila Ghaffari #ifndef eulervortex_h
17a515125bSLeila Ghaffari #define eulervortex_h
18a515125bSLeila Ghaffari 
19a515125bSLeila Ghaffari #include <math.h>
203a8779fbSJames Wright #include <ceed.h>
21a515125bSLeila Ghaffari 
22a515125bSLeila Ghaffari #ifndef M_PI
23a515125bSLeila Ghaffari #define M_PI    3.14159265358979323846
24a515125bSLeila Ghaffari #endif
25a515125bSLeila Ghaffari 
26a515125bSLeila Ghaffari #ifndef euler_context_struct
27a515125bSLeila Ghaffari #define euler_context_struct
28a515125bSLeila Ghaffari typedef struct EulerContext_ *EulerContext;
29a515125bSLeila Ghaffari struct EulerContext_ {
30a515125bSLeila Ghaffari   CeedScalar center[3];
31a515125bSLeila Ghaffari   CeedScalar curr_time;
32a515125bSLeila Ghaffari   CeedScalar vortex_strength;
33d8a22b9eSJed Brown   CeedScalar c_tau;
34a515125bSLeila Ghaffari   CeedScalar mean_velocity[3];
35a515125bSLeila Ghaffari   bool implicit;
36139613f2SLeila Ghaffari   int euler_test;
37139613f2SLeila Ghaffari   int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG
38a515125bSLeila Ghaffari };
39a515125bSLeila Ghaffari #endif
40a515125bSLeila Ghaffari 
41a515125bSLeila Ghaffari // *****************************************************************************
42a515125bSLeila Ghaffari // This function sets the initial conditions
43a515125bSLeila Ghaffari //
44a515125bSLeila Ghaffari //   Temperature:
45a515125bSLeila Ghaffari //     T   = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2)
46a515125bSLeila Ghaffari //   Density:
47a515125bSLeila Ghaffari //     rho = (T/S_vortex)^(1 / (gamma - 1))
48a515125bSLeila Ghaffari //   Pressure:
49a515125bSLeila Ghaffari //     P   = rho * T
50a515125bSLeila Ghaffari //   Velocity:
51a515125bSLeila Ghaffari //     ui  = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi)
52a515125bSLeila Ghaffari //     r   = sqrt( (x - xc)**2 + (y - yc)**2 )
53a515125bSLeila Ghaffari //   Velocity/Momentum Density:
54a515125bSLeila Ghaffari //     Ui  = rho ui
55a515125bSLeila Ghaffari //   Total Energy:
56a515125bSLeila Ghaffari //     E   = P / (gamma - 1) + rho (u u)/2
57a515125bSLeila Ghaffari //
58a515125bSLeila Ghaffari // Constants:
59a515125bSLeila Ghaffari //   cv              ,  Specific heat, constant volume
60a515125bSLeila Ghaffari //   cp              ,  Specific heat, constant pressure
61a515125bSLeila Ghaffari //   vortex_strength ,  Strength of vortex
62a515125bSLeila Ghaffari //   center          ,  Location of bubble center
63a515125bSLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
64a515125bSLeila Ghaffari //
65a515125bSLeila Ghaffari // *****************************************************************************
66a515125bSLeila Ghaffari 
67a515125bSLeila Ghaffari // *****************************************************************************
68a515125bSLeila Ghaffari // This helper function provides support for the exact, time-dependent solution
69a515125bSLeila Ghaffari //   (currently not implemented) and IC formulation for Euler traveling vortex
70a515125bSLeila Ghaffari // *****************************************************************************
71a515125bSLeila Ghaffari CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time,
72a515125bSLeila Ghaffari                                       const CeedScalar X[], CeedInt Nf, CeedScalar q[],
73a515125bSLeila Ghaffari                                       void *ctx) {
74a515125bSLeila Ghaffari   // Context
75a515125bSLeila Ghaffari   const EulerContext context = (EulerContext)ctx;
76a515125bSLeila Ghaffari   const CeedScalar vortex_strength    = context->vortex_strength;
77a515125bSLeila Ghaffari   const CeedScalar *center            = context->center; // Center of the domain
78a515125bSLeila Ghaffari   const CeedScalar *mean_velocity = context->mean_velocity;
79a515125bSLeila Ghaffari 
80a515125bSLeila Ghaffari   // Setup
81a515125bSLeila Ghaffari   const CeedScalar gamma = 1.4;
82a515125bSLeila Ghaffari   const CeedScalar cv    = 2.5;
83a515125bSLeila Ghaffari   const CeedScalar R     = 1.;
84a515125bSLeila Ghaffari   const CeedScalar x     = X[0], y = X[1]; // Coordinates
85a515125bSLeila Ghaffari   // Vortex center
86a515125bSLeila Ghaffari   const CeedScalar xc = center[0] + mean_velocity[0] * time;
87a515125bSLeila Ghaffari   const CeedScalar yc = center[1] + mean_velocity[1] * time;
88a515125bSLeila Ghaffari 
89a515125bSLeila Ghaffari   const CeedScalar x0       = x - xc;
90a515125bSLeila Ghaffari   const CeedScalar y0       = y - yc;
91a515125bSLeila Ghaffari   const CeedScalar r        = sqrt( x0*x0 + y0*y0 );
92a515125bSLeila Ghaffari   const CeedScalar C        = vortex_strength * exp((1. - r*r)/2.) / (2. * M_PI);
93139613f2SLeila Ghaffari   const CeedScalar delta_T  = - (gamma - 1.) * vortex_strength * vortex_strength *
94139613f2SLeila Ghaffari                               exp(1 - r*r) / (8. * gamma * M_PI * M_PI);
95a515125bSLeila Ghaffari   const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma
96a515125bSLeila Ghaffari   const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength /
97a515125bSLeila Ghaffari                               (8.*gamma*M_PI*M_PI);
98a515125bSLeila Ghaffari   CeedScalar rho, P, T, E, u[3] = {0.};
99a515125bSLeila Ghaffari 
100a515125bSLeila Ghaffari   // Initial Conditions
101a515125bSLeila Ghaffari   switch (context->euler_test) {
102a515125bSLeila Ghaffari   case 0: // Traveling vortex
103a515125bSLeila Ghaffari     T = 1 + delta_T;
104a515125bSLeila Ghaffari     // P = rho * T
105a515125bSLeila Ghaffari     // P = S * rho^gamma
106a515125bSLeila Ghaffari     // Solve for rho, then substitute for P
107139613f2SLeila Ghaffari     rho  = pow(T/S_vortex, 1 / (gamma - 1.));
108a515125bSLeila Ghaffari     P    = rho * T;
109a515125bSLeila Ghaffari     u[0] = mean_velocity[0] - C*y0;
110a515125bSLeila Ghaffari     u[1] = mean_velocity[1] + C*x0;
111a515125bSLeila Ghaffari 
112a515125bSLeila Ghaffari     // Assign exact solution
113a515125bSLeila Ghaffari     q[0] = rho;
114a515125bSLeila Ghaffari     q[1] = rho * u[0];
115a515125bSLeila Ghaffari     q[2] = rho * u[1];
116a515125bSLeila Ghaffari     q[3] = rho * u[2];
117a515125bSLeila Ghaffari     q[4] = P / (gamma - 1.) + rho * (u[0]*u[0] + u[1]*u[1]) / 2.;
118a515125bSLeila Ghaffari     break;
119a515125bSLeila Ghaffari   case 1: // Constant zero velocity, density constant, total energy constant
120a515125bSLeila Ghaffari     rho  = 1.;
121a515125bSLeila Ghaffari     E    = 2.;
122a515125bSLeila Ghaffari 
123a515125bSLeila Ghaffari     // Assign exact solution
124a515125bSLeila Ghaffari     q[0] = rho;
125a515125bSLeila Ghaffari     q[1] = rho * u[0];
126a515125bSLeila Ghaffari     q[2] = rho * u[1];
127a515125bSLeila Ghaffari     q[3] = rho * u[2];
128a515125bSLeila Ghaffari     q[4] = E;
129a515125bSLeila Ghaffari     break;
130a515125bSLeila Ghaffari   case 2: // Constant nonzero velocity, density constant, total energy constant
131a515125bSLeila Ghaffari     rho  = 1.;
132a515125bSLeila Ghaffari     E    = 2.;
133a515125bSLeila Ghaffari     u[0] = mean_velocity[0];
134a515125bSLeila Ghaffari     u[1] = mean_velocity[1];
135a515125bSLeila Ghaffari 
136a515125bSLeila Ghaffari     // Assign exact solution
137a515125bSLeila Ghaffari     q[0] = rho;
138a515125bSLeila Ghaffari     q[1] = rho * u[0];
139a515125bSLeila Ghaffari     q[2] = rho * u[1];
140a515125bSLeila Ghaffari     q[3] = rho * u[2];
141a515125bSLeila Ghaffari     q[4] = E;
142a515125bSLeila Ghaffari     break;
143a515125bSLeila Ghaffari   case 3: // Velocity zero, pressure constant
144a515125bSLeila Ghaffari     // (so density and internal energy will be non-constant),
145a515125bSLeila Ghaffari     // but the velocity should stay zero and the bubble won't diffuse
146a515125bSLeila Ghaffari     // (for Euler, where there is no thermal conductivity)
147a515125bSLeila Ghaffari     P    = 1.;
148a515125bSLeila Ghaffari     T    = 1. - S_bubble * exp(1. - r*r);
149a515125bSLeila Ghaffari     rho  = P / (R*T);
150a515125bSLeila Ghaffari 
151a515125bSLeila Ghaffari     // Assign exact solution
152a515125bSLeila Ghaffari     q[0] = rho;
153a515125bSLeila Ghaffari     q[1] = rho * u[0];
154a515125bSLeila Ghaffari     q[2] = rho * u[1];
155a515125bSLeila Ghaffari     q[3] = rho * u[2];
156a515125bSLeila Ghaffari     q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.);
157a515125bSLeila Ghaffari     break;
158a515125bSLeila Ghaffari   case 4: // Constant nonzero velocity, pressure constant
159a515125bSLeila Ghaffari     // (so density and internal energy will be non-constant),
160a515125bSLeila Ghaffari     // it should be transported across the domain, but velocity stays constant
161a515125bSLeila Ghaffari     P    = 1.;
162a515125bSLeila Ghaffari     T    = 1. - S_bubble * exp(1. - r*r);
163a515125bSLeila Ghaffari     rho  = P / (R*T);
164a515125bSLeila Ghaffari     u[0] = mean_velocity[0];
165a515125bSLeila Ghaffari     u[1] = mean_velocity[1];
166a515125bSLeila Ghaffari 
167a515125bSLeila Ghaffari     // Assign exact solution
168a515125bSLeila Ghaffari     q[0] = rho;
169a515125bSLeila Ghaffari     q[1] = rho * u[0];
170a515125bSLeila Ghaffari     q[2] = rho * u[1];
171a515125bSLeila Ghaffari     q[3] = rho * u[2];
172a515125bSLeila Ghaffari     q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.);
173a515125bSLeila Ghaffari     break;
1740df2634dSLeila Ghaffari   case 5: // non-smooth thermal bubble - cylinder
1750df2634dSLeila Ghaffari     P    = 1.;
1760df2634dSLeila Ghaffari     T = 1. - (r < 1. ? S_bubble : 0.);
1770df2634dSLeila Ghaffari     rho  = P / (R*T);
1780df2634dSLeila Ghaffari     u[0] = mean_velocity[0];
1790df2634dSLeila Ghaffari     u[1] = mean_velocity[1];
1800df2634dSLeila Ghaffari 
1810df2634dSLeila Ghaffari     // Assign exact solution
1820df2634dSLeila Ghaffari     q[0] = rho;
1830df2634dSLeila Ghaffari     q[1] = rho * u[0];
1840df2634dSLeila Ghaffari     q[2] = rho * u[1];
1850df2634dSLeila Ghaffari     q[3] = rho * u[2];
1860df2634dSLeila Ghaffari     q[4] = rho * (cv * T + (u[0]*u[0] + u[1]*u[1])/2.);
1870df2634dSLeila Ghaffari     break;
188a515125bSLeila Ghaffari   }
189a515125bSLeila Ghaffari   // Return
190a515125bSLeila Ghaffari   return 0;
191a515125bSLeila Ghaffari }
192a515125bSLeila Ghaffari 
193a515125bSLeila Ghaffari // *****************************************************************************
194139613f2SLeila Ghaffari // Helper function for computing flux Jacobian
195139613f2SLeila Ghaffari // *****************************************************************************
196d8a22b9eSJed Brown CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5],
197139613f2SLeila Ghaffari     const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
198139613f2SLeila Ghaffari     const CeedScalar gamma) {
199139613f2SLeila Ghaffari   CeedScalar u_sq = u[0]*u[0] + u[1]*u[1] + u[2]*u[2]; // Velocity square
200139613f2SLeila Ghaffari   for (CeedInt i=0; i<3; i++) { // Jacobian matrices for 3 directions
201139613f2SLeila Ghaffari     for (CeedInt j=0; j<3; j++) { // Rows of each Jacobian matrix
202139613f2SLeila Ghaffari       dF[i][j+1][0] = ((i==j) ? ((gamma-1.)*(u_sq/2.)) : 0.) - u[i]*u[j];
203139613f2SLeila Ghaffari       for (CeedInt k=0; k<3; k++) { // Columns of each Jacobian matrix
204139613f2SLeila Ghaffari         dF[i][0][k+1]   = ((i==k) ? 1. : 0.);
205139613f2SLeila Ghaffari         dF[i][j+1][k+1] = ((j==k) ? u[i] : 0.) +
206139613f2SLeila Ghaffari                           ((i==k) ? u[j] : 0.) -
207139613f2SLeila Ghaffari                           ((i==j) ? u[k] : 0.) * (gamma-1.);
208139613f2SLeila Ghaffari         dF[i][4][k+1]   = ((i==k) ? (E*gamma/rho - (gamma-1.)*u_sq/2.) : 0.) -
209139613f2SLeila Ghaffari                           (gamma-1.)*u[i]*u[k];
210139613f2SLeila Ghaffari       }
211139613f2SLeila Ghaffari       dF[i][j+1][4] = ((i==j) ? (gamma-1.) : 0.);
212139613f2SLeila Ghaffari     }
213139613f2SLeila Ghaffari     dF[i][4][0] = u[i] * ((gamma-1.)*u_sq - E*gamma/rho);
214139613f2SLeila Ghaffari     dF[i][4][4] = u[i] * gamma;
215139613f2SLeila Ghaffari   }
216139613f2SLeila Ghaffari }
217139613f2SLeila Ghaffari 
218139613f2SLeila Ghaffari // *****************************************************************************
219d8a22b9eSJed Brown // Helper function for computing Tau elements (stabilization constant)
220d8a22b9eSJed Brown //   Model from:
221d8a22b9eSJed Brown //     Stabilized Methods for Compressible Flows, Hughes et al 2010
222d8a22b9eSJed Brown //
223d8a22b9eSJed Brown //   Spatial criterion #2 - Tau is a 3x3 diagonal matrix
224d8a22b9eSJed Brown //   Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
225d8a22b9eSJed Brown //
226d8a22b9eSJed Brown // Where
227d8a22b9eSJed Brown //   c_tau     = stabilization constant (0.5 is reported as "optimal")
228d8a22b9eSJed Brown //   h[i]      = 2 length(dxdX[i])
229d8a22b9eSJed Brown //   Pe        = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
230d8a22b9eSJed Brown //   Xi(Pe)    = coth Pe - 1. / Pe (1. at large local Peclet number )
231d8a22b9eSJed Brown //   rho(A[i]) = spectral radius of the convective flux Jacobian i,
232d8a22b9eSJed Brown //               wave speed in direction i
233d8a22b9eSJed Brown // *****************************************************************************
234d8a22b9eSJed Brown CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3],
235d8a22b9eSJed Brown                                        const CeedScalar dXdx[3][3], const CeedScalar u[3],
236d8a22b9eSJed Brown                                        const CeedScalar sound_speed, const CeedScalar c_tau) {
237d8a22b9eSJed Brown   for (int i=0; i<3; i++) {
238d8a22b9eSJed Brown     // length of element in direction i
239d8a22b9eSJed Brown     CeedScalar h = 2 / sqrt(dXdx[0][i]*dXdx[0][i] + dXdx[1][i]*dXdx[1][i] +
240d8a22b9eSJed Brown                             dXdx[2][i]*dXdx[2][i]);
241d8a22b9eSJed Brown     // fastest wave in direction i
242d8a22b9eSJed Brown     CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
243d8a22b9eSJed Brown     Tau_x[i] = c_tau * h / fastest_wave;
244d8a22b9eSJed Brown   }
245d8a22b9eSJed Brown }
246d8a22b9eSJed Brown 
247d8a22b9eSJed Brown // *****************************************************************************
248a515125bSLeila Ghaffari // This QFunction sets the initial conditions for Euler traveling vortex
249a515125bSLeila Ghaffari // *****************************************************************************
250a515125bSLeila Ghaffari CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q,
251a515125bSLeila Ghaffari                          const CeedScalar *const *in, CeedScalar *const *out) {
252a515125bSLeila Ghaffari   // Inputs
253a515125bSLeila Ghaffari   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
254a515125bSLeila Ghaffari 
255a515125bSLeila Ghaffari   // Outputs
256a515125bSLeila Ghaffari   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
257a515125bSLeila Ghaffari   const EulerContext context = (EulerContext)ctx;
258a515125bSLeila Ghaffari 
259a515125bSLeila Ghaffari   CeedPragmaSIMD
260a515125bSLeila Ghaffari   // Quadrature Point Loop
261a515125bSLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
262a515125bSLeila Ghaffari     const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
263139613f2SLeila Ghaffari     CeedScalar q[5] = {0.};
264a515125bSLeila Ghaffari 
265a515125bSLeila Ghaffari     Exact_Euler(3, context->curr_time, x, 5, q, ctx);
266a515125bSLeila Ghaffari 
267a515125bSLeila Ghaffari     for (CeedInt j=0; j<5; j++)
268a515125bSLeila Ghaffari       q0[j][i] = q[j];
269a515125bSLeila Ghaffari   } // End of Quadrature Point Loop
270a515125bSLeila Ghaffari 
271a515125bSLeila Ghaffari   // Return
272a515125bSLeila Ghaffari   return 0;
273a515125bSLeila Ghaffari }
274a515125bSLeila Ghaffari 
275a515125bSLeila Ghaffari // *****************************************************************************
276a515125bSLeila Ghaffari // This QFunction implements the following formulation of Euler equations
277a515125bSLeila Ghaffari //   with explicit time stepping method
278a515125bSLeila Ghaffari //
279a515125bSLeila Ghaffari // This is 3D Euler for compressible gas dynamics in conservation
280a515125bSLeila Ghaffari //   form with state variables of density, momentum density, and total
281a515125bSLeila Ghaffari //   energy density.
282a515125bSLeila Ghaffari //
283a515125bSLeila Ghaffari // State Variables: q = ( rho, U1, U2, U3, E )
284a515125bSLeila Ghaffari //   rho - Mass Density
285a515125bSLeila Ghaffari //   Ui  - Momentum Density,      Ui = rho ui
286a515125bSLeila Ghaffari //   E   - Total Energy Density,  E  = P / (gamma - 1) + rho (u u)/2
287a515125bSLeila Ghaffari //
288a515125bSLeila Ghaffari // Euler Equations:
289a515125bSLeila Ghaffari //   drho/dt + div( U )                   = 0
290a515125bSLeila Ghaffari //   dU/dt   + div( rho (u x u) + P I3 )  = 0
291a515125bSLeila Ghaffari //   dE/dt   + div( (E + P) u )           = 0
292a515125bSLeila Ghaffari //
293a515125bSLeila Ghaffari // Equation of State:
294a515125bSLeila Ghaffari //   P = (gamma - 1) (E - rho (u u) / 2)
295a515125bSLeila Ghaffari //
296a515125bSLeila Ghaffari // Constants:
297a515125bSLeila Ghaffari //   cv              ,  Specific heat, constant volume
298a515125bSLeila Ghaffari //   cp              ,  Specific heat, constant pressure
299a515125bSLeila Ghaffari //   g               ,  Gravity
300a515125bSLeila Ghaffari //   gamma  = cp / cv,  Specific heat ratio
301a515125bSLeila Ghaffari // *****************************************************************************
302a515125bSLeila Ghaffari CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q,
303a515125bSLeila Ghaffari                       const CeedScalar *const *in, CeedScalar *const *out) {
304a515125bSLeila Ghaffari   // *INDENT-OFF*
305a515125bSLeila Ghaffari   // Inputs
306a515125bSLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
307139613f2SLeila Ghaffari                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
308a515125bSLeila Ghaffari                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
309a515125bSLeila Ghaffari   // Outputs
310a515125bSLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
311a515125bSLeila Ghaffari              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
312a515125bSLeila Ghaffari 
313139613f2SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
314d8a22b9eSJed Brown   const CeedScalar c_tau = context->c_tau;
315a515125bSLeila Ghaffari   const CeedScalar gamma = 1.4;
316a515125bSLeila Ghaffari 
317a515125bSLeila Ghaffari   CeedPragmaSIMD
318a515125bSLeila Ghaffari   // Quadrature Point Loop
319a515125bSLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
320a515125bSLeila Ghaffari     // *INDENT-OFF*
321a515125bSLeila Ghaffari     // Setup
322a515125bSLeila Ghaffari     // -- Interp in
323a515125bSLeila Ghaffari     const CeedScalar rho        =   q[0][i];
324a515125bSLeila Ghaffari     const CeedScalar u[3]       =  {q[1][i] / rho,
325a515125bSLeila Ghaffari                                     q[2][i] / rho,
326a515125bSLeila Ghaffari                                     q[3][i] / rho
327a515125bSLeila Ghaffari                                    };
328a515125bSLeila Ghaffari     const CeedScalar E          =   q[4][i];
329139613f2SLeila Ghaffari     const CeedScalar drho[3]    =  {dq[0][0][i],
330139613f2SLeila Ghaffari                                     dq[1][0][i],
331139613f2SLeila Ghaffari                                     dq[2][0][i]
332139613f2SLeila Ghaffari                                    };
333139613f2SLeila Ghaffari     const CeedScalar dU[3][3]   = {{dq[0][1][i],
334139613f2SLeila Ghaffari                                     dq[1][1][i],
335139613f2SLeila Ghaffari                                     dq[2][1][i]},
336139613f2SLeila Ghaffari                                    {dq[0][2][i],
337139613f2SLeila Ghaffari                                     dq[1][2][i],
338139613f2SLeila Ghaffari                                     dq[2][2][i]},
339139613f2SLeila Ghaffari                                    {dq[0][3][i],
340139613f2SLeila Ghaffari                                     dq[1][3][i],
341139613f2SLeila Ghaffari                                     dq[2][3][i]}
342139613f2SLeila Ghaffari                                   };
343139613f2SLeila Ghaffari     const CeedScalar dE[3]      =  {dq[0][4][i],
344139613f2SLeila Ghaffari                                     dq[1][4][i],
345139613f2SLeila Ghaffari                                     dq[2][4][i]
346139613f2SLeila Ghaffari                                    };
347a515125bSLeila Ghaffari     // -- Interp-to-Interp q_data
348a515125bSLeila Ghaffari     const CeedScalar wdetJ      =   q_data[0][i];
349a515125bSLeila Ghaffari     // -- Interp-to-Grad q_data
350a515125bSLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
351a515125bSLeila Ghaffari     // *INDENT-OFF*
352a515125bSLeila Ghaffari     const CeedScalar dXdx[3][3] = {{q_data[1][i],
353a515125bSLeila Ghaffari                                     q_data[2][i],
354a515125bSLeila Ghaffari                                     q_data[3][i]},
355a515125bSLeila Ghaffari                                    {q_data[4][i],
356a515125bSLeila Ghaffari                                     q_data[5][i],
357a515125bSLeila Ghaffari                                     q_data[6][i]},
358a515125bSLeila Ghaffari                                    {q_data[7][i],
359a515125bSLeila Ghaffari                                     q_data[8][i],
360a515125bSLeila Ghaffari                                     q_data[9][i]}
361a515125bSLeila Ghaffari                                   };
362a515125bSLeila Ghaffari     // *INDENT-ON*
363139613f2SLeila Ghaffari     // dU/dx
364139613f2SLeila Ghaffari     CeedScalar drhodx[3] = {0.};
365139613f2SLeila Ghaffari     CeedScalar dEdx[3] = {0.};
366139613f2SLeila Ghaffari     CeedScalar dUdx[3][3] = {{0.}};
367139613f2SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
368139613f2SLeila Ghaffari     for (int j=0; j<3; j++) {
369139613f2SLeila Ghaffari       for (int k=0; k<3; k++) {
370139613f2SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
371139613f2SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
372139613f2SLeila Ghaffari         for (int l=0; l<3; l++) {
373139613f2SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
374139613f2SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l];  //dXdx_j,k * dXdx_k,j
375139613f2SLeila Ghaffari         }
376139613f2SLeila Ghaffari       }
377139613f2SLeila Ghaffari     }
378139613f2SLeila Ghaffari     // Pressure
379a515125bSLeila Ghaffari     const CeedScalar
380a515125bSLeila Ghaffari     E_kinetic  = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]),
381a515125bSLeila Ghaffari     E_internal = E - E_kinetic,
382139613f2SLeila Ghaffari     P          = E_internal * (gamma - 1.); // P = pressure
383a515125bSLeila Ghaffari 
384a515125bSLeila Ghaffari     // The Physics
385a515125bSLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
386a515125bSLeila Ghaffari     for (int j=0; j<5; j++) {
387139613f2SLeila Ghaffari       v[j][i] = 0.;
388a515125bSLeila Ghaffari       for (int k=0; k<3; k++)
389139613f2SLeila Ghaffari         dv[k][j][i] = 0.;
390a515125bSLeila Ghaffari     }
391a515125bSLeila Ghaffari 
392a515125bSLeila Ghaffari     // -- Density
393a515125bSLeila Ghaffari     // ---- u rho
394a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
395a515125bSLeila Ghaffari       dv[j][0][i]  += wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] +
396a515125bSLeila Ghaffari                              rho*u[2]*dXdx[j][2]);
397a515125bSLeila Ghaffari     // -- Momentum
398a515125bSLeila Ghaffari     // ---- rho (u x u) + P I3
399a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
400a515125bSLeila Ghaffari       for (int k=0; k<3; k++)
401139613f2SLeila Ghaffari         dv[k][j+1][i]  += wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] +
402139613f2SLeila Ghaffari                                  (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] +
403139613f2SLeila Ghaffari                                  (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]);
404a515125bSLeila Ghaffari     // -- Total Energy Density
405a515125bSLeila Ghaffari     // ---- (E + P) u
406a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
407a515125bSLeila Ghaffari       dv[j][4][i]  += wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] +
408a515125bSLeila Ghaffari                                          u[2]*dXdx[j][2]);
409139613f2SLeila Ghaffari 
410139613f2SLeila Ghaffari     // --Stabilization terms
411139613f2SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
412139613f2SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
413d8a22b9eSJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
414139613f2SLeila Ghaffari 
415139613f2SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
416139613f2SLeila Ghaffari     CeedScalar dqdx[5][3];
417139613f2SLeila Ghaffari     for (int j=0; j<3; j++) {
418139613f2SLeila Ghaffari       dqdx[0][j] = drhodx[j];
419139613f2SLeila Ghaffari       dqdx[4][j] = dEdx[j];
420139613f2SLeila Ghaffari       for (int k=0; k<3; k++)
421139613f2SLeila Ghaffari         dqdx[k+1][j] = dUdx[k][j];
422139613f2SLeila Ghaffari     }
423139613f2SLeila Ghaffari 
424139613f2SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
425139613f2SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
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           strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
430139613f2SLeila Ghaffari 
431d8a22b9eSJed Brown     // Stabilization
432d8a22b9eSJed Brown     // -- Tau elements
433d8a22b9eSJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
434d8a22b9eSJed Brown     CeedScalar Tau_x[3] = {0.};
435d8a22b9eSJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
436139613f2SLeila Ghaffari 
437d8a22b9eSJed Brown     // -- Stabilization method: none or SU
438*bb8a0c61SJames Wright     CeedScalar stab[5][3] = {{0.}};
439139613f2SLeila Ghaffari     switch (context->stabilization) {
440139613f2SLeila Ghaffari     case 0:        // Galerkin
441139613f2SLeila Ghaffari       break;
442139613f2SLeila Ghaffari     case 1:        // SU
443139613f2SLeila Ghaffari       for (int j=0; j<3; j++)
444139613f2SLeila Ghaffari         for (int k=0; k<5; k++)
445139613f2SLeila Ghaffari           for (int l=0; l<5; l++)
446*bb8a0c61SJames Wright             stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
447139613f2SLeila Ghaffari 
448139613f2SLeila Ghaffari       for (int j=0; j<5; j++)
449139613f2SLeila Ghaffari         for (int k=0; k<3; k++)
450139613f2SLeila Ghaffari           dv[k][j][i] -= wdetJ*(stab[j][0] * dXdx[k][0] +
451139613f2SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
452139613f2SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
453139613f2SLeila Ghaffari       break;
454139613f2SLeila Ghaffari     case 2:        // SUPG is not implemented for explicit scheme
455139613f2SLeila Ghaffari       break;
456139613f2SLeila Ghaffari     }
457139613f2SLeila Ghaffari 
458a515125bSLeila Ghaffari   } // End Quadrature Point Loop
459a515125bSLeila Ghaffari 
460a515125bSLeila Ghaffari   // Return
461a515125bSLeila Ghaffari   return 0;
462a515125bSLeila Ghaffari }
463a515125bSLeila Ghaffari 
464a515125bSLeila Ghaffari // *****************************************************************************
465a515125bSLeila Ghaffari // This QFunction implements the Euler equations with (mentioned above)
466a515125bSLeila Ghaffari //   with implicit time stepping method
467a515125bSLeila Ghaffari //
468a515125bSLeila Ghaffari // *****************************************************************************
469a515125bSLeila Ghaffari CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q,
470a515125bSLeila Ghaffari                                 const CeedScalar *const *in, CeedScalar *const *out) {
471a515125bSLeila Ghaffari   // *INDENT-OFF*
472a515125bSLeila Ghaffari   // Inputs
473a515125bSLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
474139613f2SLeila Ghaffari                    (*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1],
475a515125bSLeila Ghaffari                    (*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2],
476a515125bSLeila Ghaffari                    (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
477a515125bSLeila Ghaffari   // Outputs
478a515125bSLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0],
479a515125bSLeila Ghaffari              (*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
480a515125bSLeila Ghaffari 
481139613f2SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
482d8a22b9eSJed Brown   const CeedScalar c_tau = context->c_tau;
483a515125bSLeila Ghaffari   const CeedScalar gamma = 1.4;
484a515125bSLeila Ghaffari 
485a515125bSLeila Ghaffari   CeedPragmaSIMD
486a515125bSLeila Ghaffari   // Quadrature Point Loop
487a515125bSLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
488a515125bSLeila Ghaffari     // *INDENT-OFF*
489a515125bSLeila Ghaffari     // Setup
490a515125bSLeila Ghaffari     // -- Interp in
491a515125bSLeila Ghaffari     const CeedScalar rho        =   q[0][i];
492a515125bSLeila Ghaffari     const CeedScalar u[3]       =  {q[1][i] / rho,
493a515125bSLeila Ghaffari                                     q[2][i] / rho,
494a515125bSLeila Ghaffari                                     q[3][i] / rho
495a515125bSLeila Ghaffari                                    };
496a515125bSLeila Ghaffari     const CeedScalar E          =   q[4][i];
497139613f2SLeila Ghaffari     const CeedScalar drho[3]    =  {dq[0][0][i],
498139613f2SLeila Ghaffari                                     dq[1][0][i],
499139613f2SLeila Ghaffari                                     dq[2][0][i]
500139613f2SLeila Ghaffari                                    };
501139613f2SLeila Ghaffari     const CeedScalar dU[3][3]   = {{dq[0][1][i],
502139613f2SLeila Ghaffari                                     dq[1][1][i],
503139613f2SLeila Ghaffari                                     dq[2][1][i]},
504139613f2SLeila Ghaffari                                    {dq[0][2][i],
505139613f2SLeila Ghaffari                                     dq[1][2][i],
506139613f2SLeila Ghaffari                                     dq[2][2][i]},
507139613f2SLeila Ghaffari                                    {dq[0][3][i],
508139613f2SLeila Ghaffari                                     dq[1][3][i],
509139613f2SLeila Ghaffari                                     dq[2][3][i]}
510139613f2SLeila Ghaffari                                   };
511139613f2SLeila Ghaffari     const CeedScalar dE[3]      =  {dq[0][4][i],
512139613f2SLeila Ghaffari                                     dq[1][4][i],
513139613f2SLeila Ghaffari                                     dq[2][4][i]
514139613f2SLeila Ghaffari                                    };
515a515125bSLeila Ghaffari     // -- Interp-to-Interp q_data
516a515125bSLeila Ghaffari     const CeedScalar wdetJ      =   q_data[0][i];
517a515125bSLeila Ghaffari     // -- Interp-to-Grad q_data
518a515125bSLeila Ghaffari     // ---- Inverse of change of coordinate matrix: X_i,j
519a515125bSLeila Ghaffari     // *INDENT-OFF*
520a515125bSLeila Ghaffari     const CeedScalar dXdx[3][3] = {{q_data[1][i],
521a515125bSLeila Ghaffari                                     q_data[2][i],
522a515125bSLeila Ghaffari                                     q_data[3][i]},
523a515125bSLeila Ghaffari                                    {q_data[4][i],
524a515125bSLeila Ghaffari                                     q_data[5][i],
525a515125bSLeila Ghaffari                                     q_data[6][i]},
526a515125bSLeila Ghaffari                                    {q_data[7][i],
527a515125bSLeila Ghaffari                                     q_data[8][i],
528a515125bSLeila Ghaffari                                     q_data[9][i]}
529a515125bSLeila Ghaffari                                   };
530a515125bSLeila Ghaffari     // *INDENT-ON*
531139613f2SLeila Ghaffari     // dU/dx
532139613f2SLeila Ghaffari     CeedScalar drhodx[3] = {0.};
533139613f2SLeila Ghaffari     CeedScalar dEdx[3] = {0.};
534139613f2SLeila Ghaffari     CeedScalar dUdx[3][3] = {{0.}};
535139613f2SLeila Ghaffari     CeedScalar dXdxdXdxT[3][3] = {{0.}};
536139613f2SLeila Ghaffari     for (int j=0; j<3; j++) {
537139613f2SLeila Ghaffari       for (int k=0; k<3; k++) {
538139613f2SLeila Ghaffari         drhodx[j] += drho[k] * dXdx[k][j];
539139613f2SLeila Ghaffari         dEdx[j] += dE[k] * dXdx[k][j];
540139613f2SLeila Ghaffari         for (int l=0; l<3; l++) {
541139613f2SLeila Ghaffari           dUdx[j][k] += dU[j][l] * dXdx[l][k];
542139613f2SLeila Ghaffari           dXdxdXdxT[j][k] += dXdx[j][l]*dXdx[k][l];  //dXdx_j,k * dXdx_k,j
543139613f2SLeila Ghaffari         }
544139613f2SLeila Ghaffari       }
545139613f2SLeila Ghaffari     }
546a515125bSLeila Ghaffari     const CeedScalar
547a515125bSLeila Ghaffari     E_kinetic  = 0.5 * rho * (u[0]*u[0] + u[1]*u[1] + u[2]*u[2]),
548a515125bSLeila Ghaffari     E_internal = E - E_kinetic,
549139613f2SLeila Ghaffari     P          = E_internal * (gamma - 1.); // P = pressure
550a515125bSLeila Ghaffari 
551a515125bSLeila Ghaffari     // The Physics
552a515125bSLeila Ghaffari     // Zero v and dv so all future terms can safely sum into it
553a515125bSLeila Ghaffari     for (int j=0; j<5; j++) {
554139613f2SLeila Ghaffari       v[j][i] = 0.;
555a515125bSLeila Ghaffari       for (int k=0; k<3; k++)
556139613f2SLeila Ghaffari         dv[k][j][i] = 0.;
557a515125bSLeila Ghaffari     }
558a515125bSLeila Ghaffari     //-----mass matrix
559a515125bSLeila Ghaffari     for (int j=0; j<5; j++)
560a515125bSLeila Ghaffari       v[j][i] += wdetJ*q_dot[j][i];
561a515125bSLeila Ghaffari 
562a515125bSLeila Ghaffari     // -- Density
563a515125bSLeila Ghaffari     // ---- u rho
564a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
565a515125bSLeila Ghaffari       dv[j][0][i]  -= wdetJ*(rho*u[0]*dXdx[j][0] + rho*u[1]*dXdx[j][1] +
566a515125bSLeila Ghaffari                              rho*u[2]*dXdx[j][2]);
567a515125bSLeila Ghaffari     // -- Momentum
568a515125bSLeila Ghaffari     // ---- rho (u x u) + P I3
569a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
570a515125bSLeila Ghaffari       for (int k=0; k<3; k++)
571139613f2SLeila Ghaffari         dv[k][j+1][i]  -= wdetJ*((rho*u[j]*u[0] + (j==0?P:0.))*dXdx[k][0] +
572139613f2SLeila Ghaffari                                  (rho*u[j]*u[1] + (j==1?P:0.))*dXdx[k][1] +
573139613f2SLeila Ghaffari                                  (rho*u[j]*u[2] + (j==2?P:0.))*dXdx[k][2]);
574a515125bSLeila Ghaffari     // -- Total Energy Density
575a515125bSLeila Ghaffari     // ---- (E + P) u
576a515125bSLeila Ghaffari     for (int j=0; j<3; j++)
577a515125bSLeila Ghaffari       dv[j][4][i]  -= wdetJ * (E + P) * (u[0]*dXdx[j][0] + u[1]*dXdx[j][1] +
578a515125bSLeila Ghaffari                                          u[2]*dXdx[j][2]);
579139613f2SLeila Ghaffari 
580139613f2SLeila Ghaffari     // -- Stabilization terms
581139613f2SLeila Ghaffari     // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
582139613f2SLeila Ghaffari     CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
583d8a22b9eSJed Brown     ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
584139613f2SLeila Ghaffari 
585139613f2SLeila Ghaffari     // ---- dqdx collects drhodx, dUdx and dEdx in one vector
586139613f2SLeila Ghaffari     CeedScalar dqdx[5][3];
587139613f2SLeila Ghaffari     for (int j=0; j<3; j++) {
588139613f2SLeila Ghaffari       dqdx[0][j] = drhodx[j];
589139613f2SLeila Ghaffari       dqdx[4][j] = dEdx[j];
590139613f2SLeila Ghaffari       for (int k=0; k<3; k++)
591139613f2SLeila Ghaffari         dqdx[k+1][j] = dUdx[k][j];
592139613f2SLeila Ghaffari     }
593139613f2SLeila Ghaffari 
594139613f2SLeila Ghaffari     // ---- strong_conv = dF/dq * dq/dx    (Strong convection)
595139613f2SLeila Ghaffari     CeedScalar strong_conv[5] = {0.};
596139613f2SLeila Ghaffari     for (int j=0; j<3; j++)
597139613f2SLeila Ghaffari       for (int k=0; k<5; k++)
598139613f2SLeila Ghaffari         for (int l=0; l<5; l++)
599139613f2SLeila Ghaffari           strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
600139613f2SLeila Ghaffari 
601139613f2SLeila Ghaffari     // ---- Strong residual
602139613f2SLeila Ghaffari     CeedScalar strong_res[5];
603139613f2SLeila Ghaffari     for (int j=0; j<5; j++)
604139613f2SLeila Ghaffari       strong_res[j] = q_dot[j][i] + strong_conv[j];
605139613f2SLeila Ghaffari 
606d8a22b9eSJed Brown     // Stabilization
607d8a22b9eSJed Brown     // -- Tau elements
608d8a22b9eSJed Brown     const CeedScalar sound_speed = sqrt(gamma * P / rho);
609d8a22b9eSJed Brown     CeedScalar Tau_x[3] = {0.};
610d8a22b9eSJed Brown     Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
611139613f2SLeila Ghaffari 
612d8a22b9eSJed Brown     // -- Stabilization method: none, SU, or SUPG
613*bb8a0c61SJames Wright     CeedScalar stab[5][3] = {{0.}};
614139613f2SLeila Ghaffari     switch (context->stabilization) {
615139613f2SLeila Ghaffari     case 0:        // Galerkin
616139613f2SLeila Ghaffari       break;
617139613f2SLeila Ghaffari     case 1:        // SU
618139613f2SLeila Ghaffari       for (int j=0; j<3; j++)
619139613f2SLeila Ghaffari         for (int k=0; k<5; k++)
620139613f2SLeila Ghaffari           for (int l=0; l<5; l++)
621*bb8a0c61SJames Wright             stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
622139613f2SLeila Ghaffari 
623139613f2SLeila Ghaffari       for (int j=0; j<5; j++)
624139613f2SLeila Ghaffari         for (int k=0; k<3; k++)
625139613f2SLeila Ghaffari           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
626139613f2SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
627139613f2SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
628139613f2SLeila Ghaffari       break;
629139613f2SLeila Ghaffari     case 2:        // SUPG
630139613f2SLeila Ghaffari       for (int j=0; j<3; j++)
631139613f2SLeila Ghaffari         for (int k=0; k<5; k++)
632139613f2SLeila Ghaffari           for (int l=0; l<5; l++)
633*bb8a0c61SJames Wright             stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l];
634139613f2SLeila Ghaffari 
635139613f2SLeila Ghaffari       for (int j=0; j<5; j++)
636139613f2SLeila Ghaffari         for (int k=0; k<3; k++)
637139613f2SLeila Ghaffari           dv[k][j][i] += wdetJ*(stab[j][0] * dXdx[k][0] +
638139613f2SLeila Ghaffari                                 stab[j][1] * dXdx[k][1] +
639139613f2SLeila Ghaffari                                 stab[j][2] * dXdx[k][2]);
640139613f2SLeila Ghaffari       break;
641139613f2SLeila Ghaffari     }
642a515125bSLeila Ghaffari   } // End Quadrature Point Loop
643a515125bSLeila Ghaffari 
644a515125bSLeila Ghaffari   // Return
645a515125bSLeila Ghaffari   return 0;
646a515125bSLeila Ghaffari }
647a515125bSLeila Ghaffari // *****************************************************************************
648002797a3SLeila Ghaffari // This QFunction sets the inflow boundary conditions for
649002797a3SLeila Ghaffari //   the traveling vortex problem.
650a515125bSLeila Ghaffari //
651a515125bSLeila Ghaffari //  Prescribed T_inlet and P_inlet are converted to conservative variables
652a515125bSLeila Ghaffari //      and applied weakly.
653a515125bSLeila Ghaffari //
654a515125bSLeila Ghaffari // *****************************************************************************
655002797a3SLeila Ghaffari CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q,
656a515125bSLeila Ghaffari                                        const CeedScalar *const *in,
657a515125bSLeila Ghaffari                                        CeedScalar *const *out) {
658a515125bSLeila Ghaffari   // *INDENT-OFF*
659a515125bSLeila Ghaffari   // Inputs
660002797a3SLeila Ghaffari   const CeedScalar (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1];
661a515125bSLeila Ghaffari   // Outputs
662a515125bSLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
663a515125bSLeila Ghaffari   // *INDENT-ON*
664a515125bSLeila Ghaffari   EulerContext context = (EulerContext)ctx;
665a515125bSLeila Ghaffari   const int euler_test      = context->euler_test;
666a515125bSLeila Ghaffari   const bool implicit       = context->implicit;
667a515125bSLeila Ghaffari   CeedScalar *mean_velocity = context->mean_velocity;
668a515125bSLeila Ghaffari   const CeedScalar cv    = 2.5;
669a515125bSLeila Ghaffari   const CeedScalar R     = 1.;
670a515125bSLeila Ghaffari   CeedScalar T_inlet;
671a515125bSLeila Ghaffari   CeedScalar P_inlet;
672a515125bSLeila Ghaffari 
673a515125bSLeila Ghaffari   // For test cases 1 and 3 the background velocity is zero
674a515125bSLeila Ghaffari   if (euler_test == 1 || euler_test == 3)
675a515125bSLeila Ghaffari     for (CeedInt i=0; i<3; i++) mean_velocity[i] = 0.;
676a515125bSLeila Ghaffari 
677a515125bSLeila Ghaffari   // For test cases 1 and 2, T_inlet = T_inlet = 0.4
678a515125bSLeila Ghaffari   if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
679a515125bSLeila Ghaffari   else T_inlet = P_inlet = 1.;
680a515125bSLeila Ghaffari 
681a515125bSLeila Ghaffari   CeedPragmaSIMD
682a515125bSLeila Ghaffari   // Quadrature Point Loop
683a515125bSLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
684a515125bSLeila Ghaffari     // Setup
685a515125bSLeila Ghaffari     // -- Interp-to-Interp q_data
686a515125bSLeila Ghaffari     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
687a515125bSLeila Ghaffari     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
688a515125bSLeila Ghaffari     // We can effect this by swapping the sign on this weight
689a515125bSLeila Ghaffari     const CeedScalar wdetJb  = (implicit ? -1. : 1.) * q_data_sur[0][i];
690002797a3SLeila Ghaffari     // ---- Normal vect
691a515125bSLeila Ghaffari     const CeedScalar norm[3] = {q_data_sur[1][i],
692a515125bSLeila Ghaffari                                 q_data_sur[2][i],
693a515125bSLeila Ghaffari                                 q_data_sur[3][i]
694a515125bSLeila Ghaffari                                };
695a515125bSLeila Ghaffari 
696a515125bSLeila Ghaffari     // face_normal = Normal vector of the face
697a515125bSLeila Ghaffari     const CeedScalar face_normal = norm[0]*mean_velocity[0] +
698a515125bSLeila Ghaffari                                    norm[1]*mean_velocity[1] +
699a515125bSLeila Ghaffari                                    norm[2]*mean_velocity[2];
700a515125bSLeila Ghaffari     // The Physics
701a515125bSLeila Ghaffari     // Zero v so all future terms can safely sum into it
702139613f2SLeila Ghaffari     for (int j=0; j<5; j++) v[j][i] = 0.;
703a515125bSLeila Ghaffari 
704a515125bSLeila Ghaffari     // Implementing in/outflow BCs
705002797a3SLeila Ghaffari     if (face_normal > 0) {
706a515125bSLeila Ghaffari     } else { // inflow
707a515125bSLeila Ghaffari       const CeedScalar rho_inlet = P_inlet/(R*T_inlet);
708a515125bSLeila Ghaffari       const CeedScalar E_kinetic_inlet = (mean_velocity[0]*mean_velocity[0] +
709a515125bSLeila Ghaffari                                           mean_velocity[1]*mean_velocity[1]) / 2.;
710a515125bSLeila Ghaffari       // incoming total energy
711a515125bSLeila Ghaffari       const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
712a515125bSLeila Ghaffari 
713a515125bSLeila Ghaffari       // The Physics
714a515125bSLeila Ghaffari       // -- Density
715a515125bSLeila Ghaffari       v[0][i] -= wdetJb * rho_inlet * face_normal;
716a515125bSLeila Ghaffari 
717a515125bSLeila Ghaffari       // -- Momentum
718a515125bSLeila Ghaffari       for (int j=0; j<3; j++)
719a515125bSLeila Ghaffari         v[j+1][i] -= wdetJb *(rho_inlet * face_normal * mean_velocity[j] +
720a515125bSLeila Ghaffari                               norm[j] * P_inlet);
721a515125bSLeila Ghaffari 
722a515125bSLeila Ghaffari       // -- Total Energy Density
723a515125bSLeila Ghaffari       v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
724a515125bSLeila Ghaffari     }
725a515125bSLeila Ghaffari 
726a515125bSLeila Ghaffari   } // End Quadrature Point Loop
727a515125bSLeila Ghaffari   return 0;
728a515125bSLeila Ghaffari }
729a515125bSLeila Ghaffari 
730a515125bSLeila Ghaffari // *****************************************************************************
73168ef3d20SLeila Ghaffari // This QFunction sets the outflow boundary conditions for
73268ef3d20SLeila Ghaffari //   the Euler solver.
73368ef3d20SLeila Ghaffari //
73468ef3d20SLeila Ghaffari //  Outflow BCs:
73568ef3d20SLeila Ghaffari //    The validity of the weak form of the governing equations is
73668ef3d20SLeila Ghaffari //      extended to the outflow.
73768ef3d20SLeila Ghaffari //
73868ef3d20SLeila Ghaffari // *****************************************************************************
73968ef3d20SLeila Ghaffari CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q,
74068ef3d20SLeila Ghaffari                               const CeedScalar *const *in,
74168ef3d20SLeila Ghaffari                               CeedScalar *const *out) {
74268ef3d20SLeila Ghaffari   // *INDENT-OFF*
74368ef3d20SLeila Ghaffari   // Inputs
74468ef3d20SLeila Ghaffari   const CeedScalar (*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0],
74568ef3d20SLeila Ghaffari                    (*q_data_sur)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1];
74668ef3d20SLeila Ghaffari   // Outputs
74768ef3d20SLeila Ghaffari   CeedScalar (*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
74868ef3d20SLeila Ghaffari   // *INDENT-ON*
74968ef3d20SLeila Ghaffari   EulerContext context = (EulerContext)ctx;
75068ef3d20SLeila Ghaffari   const bool implicit       = context->implicit;
75168ef3d20SLeila Ghaffari   CeedScalar *mean_velocity = context->mean_velocity;
75268ef3d20SLeila Ghaffari 
75368ef3d20SLeila Ghaffari   const CeedScalar gamma = 1.4;
75468ef3d20SLeila Ghaffari 
75568ef3d20SLeila Ghaffari   CeedPragmaSIMD
75668ef3d20SLeila Ghaffari   // Quadrature Point Loop
75768ef3d20SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
75868ef3d20SLeila Ghaffari     // Setup
75968ef3d20SLeila Ghaffari     // -- Interp in
76068ef3d20SLeila Ghaffari     const CeedScalar rho      =  q[0][i];
76168ef3d20SLeila Ghaffari     const CeedScalar u[3]     = {q[1][i] / rho,
76268ef3d20SLeila Ghaffari                                  q[2][i] / rho,
76368ef3d20SLeila Ghaffari                                  q[3][i] / rho
76468ef3d20SLeila Ghaffari                                 };
76568ef3d20SLeila Ghaffari     const CeedScalar E        =  q[4][i];
76668ef3d20SLeila Ghaffari 
76768ef3d20SLeila Ghaffari     // -- Interp-to-Interp q_data
76868ef3d20SLeila Ghaffari     // For explicit mode, the surface integral is on the RHS of ODE q_dot = f(q).
76968ef3d20SLeila Ghaffari     // For implicit mode, it gets pulled to the LHS of implicit ODE/DAE g(q_dot, q).
77068ef3d20SLeila Ghaffari     // We can effect this by swapping the sign on this weight
77168ef3d20SLeila Ghaffari     const CeedScalar wdetJb     =   (implicit ? -1. : 1.) * q_data_sur[0][i];
77268ef3d20SLeila Ghaffari     // ---- Normal vectors
77368ef3d20SLeila Ghaffari     const CeedScalar norm[3]    =   {q_data_sur[1][i],
77468ef3d20SLeila Ghaffari                                      q_data_sur[2][i],
77568ef3d20SLeila Ghaffari                                      q_data_sur[3][i]
77668ef3d20SLeila Ghaffari                                     };
77768ef3d20SLeila Ghaffari 
77868ef3d20SLeila Ghaffari     // face_normal = Normal vector of the face
77968ef3d20SLeila Ghaffari     const CeedScalar face_normal = norm[0]*mean_velocity[0] +
78068ef3d20SLeila Ghaffari                                    norm[1]*mean_velocity[1] +
78168ef3d20SLeila Ghaffari                                    norm[2]*mean_velocity[2];
78268ef3d20SLeila Ghaffari     // The Physics
78368ef3d20SLeila Ghaffari     // Zero v so all future terms can safely sum into it
78468ef3d20SLeila Ghaffari     for (int j=0; j<5; j++) v[j][i] = 0;
78568ef3d20SLeila Ghaffari 
78668ef3d20SLeila Ghaffari     // Implementing in/outflow BCs
78768ef3d20SLeila Ghaffari     if (face_normal > 0) { // outflow
78868ef3d20SLeila Ghaffari       const CeedScalar E_kinetic = (u[0]*u[0] + u[1]*u[1]) / 2.;
78968ef3d20SLeila Ghaffari       const CeedScalar P         = (E - E_kinetic * rho) * (gamma - 1.); // pressure
79068ef3d20SLeila Ghaffari       const CeedScalar u_normal  = norm[0]*u[0] + norm[1]*u[1] +
79168ef3d20SLeila Ghaffari                                    norm[2]*u[2]; // Normal velocity
79268ef3d20SLeila Ghaffari       // The Physics
79368ef3d20SLeila Ghaffari       // -- Density
79468ef3d20SLeila Ghaffari       v[0][i] -= wdetJb * rho * u_normal;
79568ef3d20SLeila Ghaffari 
79668ef3d20SLeila Ghaffari       // -- Momentum
79768ef3d20SLeila Ghaffari       for (int j=0; j<3; j++)
79868ef3d20SLeila Ghaffari         v[j+1][i] -= wdetJb *(rho * u_normal * u[j] + norm[j] * P);
79968ef3d20SLeila Ghaffari 
80068ef3d20SLeila Ghaffari       // -- Total Energy Density
80168ef3d20SLeila Ghaffari       v[4][i] -= wdetJb * u_normal * (E + P);
80268ef3d20SLeila Ghaffari     }
80368ef3d20SLeila Ghaffari   } // End Quadrature Point Loop
80468ef3d20SLeila Ghaffari   return 0;
80568ef3d20SLeila Ghaffari }
80668ef3d20SLeila Ghaffari 
80768ef3d20SLeila Ghaffari // *****************************************************************************
808a515125bSLeila Ghaffari 
809a515125bSLeila Ghaffari #endif // eulervortex_h
810