1 // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED: http://github.com/ceed
7
8 /// @file
9 /// Euler traveling vortex initial condition and operator for Navier-Stokes
10 /// example using PETSc
11
12 // Model from:
13 // On the Order of Accuracy and Numerical Performance of Two Classes of Finite Volume WENO Schemes, Zhang, Zhang, and Shu (2011).
14 #include <ceed/types.h>
15 #ifndef CEED_RUNNING_JIT_PASS
16 #include <math.h>
17 #include <stdbool.h>
18 #endif
19
20 #include "utils.h"
21
22 typedef struct EulerContext_ *EulerContext;
23 struct EulerContext_ {
24 CeedScalar center[3];
25 CeedScalar curr_time;
26 CeedScalar vortex_strength;
27 CeedScalar c_tau;
28 CeedScalar mean_velocity[3];
29 bool implicit;
30 int euler_test;
31 int stabilization; // See StabilizationType: 0=none, 1=SU, 2=SUPG
32 };
33
34 // *****************************************************************************
35 // This function sets the initial conditions
36 //
37 // Temperature:
38 // T = 1 - (gamma - 1) vortex_strength**2 exp(1 - r**2) / (8 gamma pi**2)
39 // Density:
40 // rho = (T/S_vortex)^(1 / (gamma - 1))
41 // Pressure:
42 // P = rho * T
43 // Velocity:
44 // ui = 1 + vortex_strength exp((1 - r**2)/2.) [yc - y, x - xc] / (2 pi)
45 // r = sqrt( (x - xc)**2 + (y - yc)**2 )
46 // Velocity/Momentum Density:
47 // Ui = rho ui
48 // Total Energy:
49 // E = P / (gamma - 1) + rho (u u)/2
50 //
51 // Constants:
52 // cv , Specific heat, constant volume
53 // cp , Specific heat, constant pressure
54 // vortex_strength , Strength of vortex
55 // center , Location of bubble center
56 // gamma = cp / cv, Specific heat ratio
57 //
58 // *****************************************************************************
59
60 // *****************************************************************************
61 // This helper function provides support for the exact, time-dependent solution (currently not implemented) and IC formulation for Euler traveling
62 // vortex
63 // *****************************************************************************
Exact_Euler(CeedInt dim,CeedScalar time,const CeedScalar X[],CeedInt Nf,CeedScalar q[],void * ctx)64 CEED_QFUNCTION_HELPER int Exact_Euler(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, CeedScalar q[], void *ctx) {
65 // Context
66 const EulerContext context = (EulerContext)ctx;
67 const CeedScalar vortex_strength = context->vortex_strength;
68 const CeedScalar *center = context->center; // Center of the domain
69 const CeedScalar *mean_velocity = context->mean_velocity;
70
71 // Setup
72 const CeedScalar gamma = 1.4;
73 const CeedScalar cv = 2.5;
74 const CeedScalar R = 1.;
75 const CeedScalar x = X[0], y = X[1]; // Coordinates
76 // Vortex center
77 const CeedScalar xc = center[0] + mean_velocity[0] * time;
78 const CeedScalar yc = center[1] + mean_velocity[1] * time;
79
80 const CeedScalar x0 = x - xc;
81 const CeedScalar y0 = y - yc;
82 const CeedScalar r = sqrt(x0 * x0 + y0 * y0);
83 const CeedScalar C = vortex_strength * exp((1. - r * r) / 2.) / (2. * M_PI);
84 const CeedScalar delta_T = -(gamma - 1.) * vortex_strength * vortex_strength * exp(1 - r * r) / (8. * gamma * M_PI * M_PI);
85 const CeedScalar S_vortex = 1; // no perturbation in the entropy P / rho^gamma
86 const CeedScalar S_bubble = (gamma - 1.) * vortex_strength * vortex_strength / (8. * gamma * M_PI * M_PI);
87 CeedScalar rho, P, T, E, u[3] = {0.};
88
89 // Initial Conditions
90 switch (context->euler_test) {
91 case 0: // Traveling vortex
92 T = 1 + delta_T;
93 // P = rho * T
94 // P = S * rho^gamma
95 // Solve for rho, then substitute for P
96 rho = pow(T / S_vortex, 1 / (gamma - 1.));
97 P = rho * T;
98 u[0] = mean_velocity[0] - C * y0;
99 u[1] = mean_velocity[1] + C * x0;
100
101 // Assign exact solution
102 q[0] = rho;
103 q[1] = rho * u[0];
104 q[2] = rho * u[1];
105 q[3] = rho * u[2];
106 q[4] = P / (gamma - 1.) + rho * (u[0] * u[0] + u[1] * u[1]) / 2.;
107 break;
108 case 1: // Constant zero velocity, density constant, total energy constant
109 rho = 1.;
110 E = 2.;
111
112 // Assign exact solution
113 q[0] = rho;
114 q[1] = rho * u[0];
115 q[2] = rho * u[1];
116 q[3] = rho * u[2];
117 q[4] = E;
118 break;
119 case 2: // Constant nonzero velocity, density constant, total energy constant
120 rho = 1.;
121 E = 2.;
122 u[0] = mean_velocity[0];
123 u[1] = mean_velocity[1];
124
125 // Assign exact solution
126 q[0] = rho;
127 q[1] = rho * u[0];
128 q[2] = rho * u[1];
129 q[3] = rho * u[2];
130 q[4] = E;
131 break;
132 case 3: // Velocity zero, pressure constant (so density and internal energy will be non-constant), but the velocity should stay zero and the
133 // bubble won't diffuse
134 // (for Euler, where there is no thermal conductivity)
135 P = 1.;
136 T = 1. - S_bubble * exp(1. - r * r);
137 rho = P / (R * T);
138
139 // Assign exact solution
140 q[0] = rho;
141 q[1] = rho * u[0];
142 q[2] = rho * u[1];
143 q[3] = rho * u[2];
144 q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
145 break;
146 case 4: // Constant nonzero velocity, pressure constant (so density and internal energy will be non-constant),
147 // It should be transported across the domain, but velocity stays constant
148 P = 1.;
149 T = 1. - S_bubble * exp(1. - r * r);
150 rho = P / (R * T);
151 u[0] = mean_velocity[0];
152 u[1] = mean_velocity[1];
153
154 // Assign exact solution
155 q[0] = rho;
156 q[1] = rho * u[0];
157 q[2] = rho * u[1];
158 q[3] = rho * u[2];
159 q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
160 break;
161 case 5: // non-smooth thermal bubble - cylinder
162 P = 1.;
163 T = 1. - (r < 1. ? S_bubble : 0.);
164 rho = P / (R * T);
165 u[0] = mean_velocity[0];
166 u[1] = mean_velocity[1];
167
168 // Assign exact solution
169 q[0] = rho;
170 q[1] = rho * u[0];
171 q[2] = rho * u[1];
172 q[3] = rho * u[2];
173 q[4] = rho * (cv * T + (u[0] * u[0] + u[1] * u[1]) / 2.);
174 break;
175 }
176 return 0;
177 }
178
179 // *****************************************************************************
180 // Helper function for computing flux Jacobian
181 // *****************************************************************************
ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5],const CeedScalar rho,const CeedScalar u[3],const CeedScalar E,const CeedScalar gamma)182 CEED_QFUNCTION_HELPER void ConvectiveFluxJacobian_Euler(CeedScalar dF[3][5][5], const CeedScalar rho, const CeedScalar u[3], const CeedScalar E,
183 const CeedScalar gamma) {
184 CeedScalar u_sq = u[0] * u[0] + u[1] * u[1] + u[2] * u[2]; // Velocity square
185 for (CeedInt i = 0; i < 3; i++) { // Jacobian matrices for 3 directions
186 for (CeedInt j = 0; j < 3; j++) { // Rows of each Jacobian matrix
187 dF[i][j + 1][0] = ((i == j) ? ((gamma - 1.) * (u_sq / 2.)) : 0.) - u[i] * u[j];
188 for (CeedInt k = 0; k < 3; k++) { // Columns of each Jacobian matrix
189 dF[i][0][k + 1] = ((i == k) ? 1. : 0.);
190 dF[i][j + 1][k + 1] = ((j == k) ? u[i] : 0.) + ((i == k) ? u[j] : 0.) - ((i == j) ? u[k] : 0.) * (gamma - 1.);
191 dF[i][4][k + 1] = ((i == k) ? (E * gamma / rho - (gamma - 1.) * u_sq / 2.) : 0.) - (gamma - 1.) * u[i] * u[k];
192 }
193 dF[i][j + 1][4] = ((i == j) ? (gamma - 1.) : 0.);
194 }
195 dF[i][4][0] = u[i] * ((gamma - 1.) * u_sq - E * gamma / rho);
196 dF[i][4][4] = u[i] * gamma;
197 }
198 }
199
200 // *****************************************************************************
201 // Helper function for computing Tau elements (stabilization constant)
202 // Model from:
203 // Stabilized Methods for Compressible Flows, Hughes et al 2010
204 //
205 // Spatial criterion #2 - Tau is a 3x3 diagonal matrix
206 // Tau[i] = c_tau h[i] Xi(Pe) / rho(A[i]) (no sum)
207 //
208 // Where
209 // c_tau = stabilization constant (0.5 is reported as "optimal")
210 // h[i] = 2 length(dxdX[i])
211 // Pe = Peclet number ( Pe = sqrt(u u) / dot(dXdx,u) diffusivity )
212 // Xi(Pe) = coth Pe - 1. / Pe (1. at large local Peclet number )
213 // rho(A[i]) = spectral radius of the convective flux Jacobian i, wave speed in direction i
214 // *****************************************************************************
Tau_spatial(CeedScalar Tau_x[3],const CeedScalar dXdx[3][3],const CeedScalar u[3],const CeedScalar sound_speed,const CeedScalar c_tau)215 CEED_QFUNCTION_HELPER void Tau_spatial(CeedScalar Tau_x[3], const CeedScalar dXdx[3][3], const CeedScalar u[3], const CeedScalar sound_speed,
216 const CeedScalar c_tau) {
217 for (CeedInt i = 0; i < 3; i++) {
218 // length of element in direction i
219 CeedScalar h = 2 / sqrt(dXdx[0][i] * dXdx[0][i] + dXdx[1][i] * dXdx[1][i] + dXdx[2][i] * dXdx[2][i]);
220 // fastest wave in direction i
221 CeedScalar fastest_wave = fabs(u[i]) + sound_speed;
222 Tau_x[i] = c_tau * h / fastest_wave;
223 }
224 }
225
226 // *****************************************************************************
227 // This QFunction sets the initial conditions for Euler traveling vortex
228 // *****************************************************************************
ICsEuler(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)229 CEED_QFUNCTION(ICsEuler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
230 const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
231 CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
232
233 const EulerContext context = (EulerContext)ctx;
234
235 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
236 const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
237 CeedScalar q[5] = {0.};
238
239 Exact_Euler(3, context->curr_time, x, 5, q, ctx);
240 for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
241 }
242 return 0;
243 }
244
245 // *****************************************************************************
246 // This QFunction implements the following formulation of Euler equations with explicit time stepping method
247 //
248 // This is 3D Euler for compressible gas dynamics in conservation form with state variables of density, momentum density, and total energy density.
249 //
250 // State Variables: q = ( rho, U1, U2, U3, E )
251 // rho - Mass Density
252 // Ui - Momentum Density, Ui = rho ui
253 // E - Total Energy Density, E = P / (gamma - 1) + rho (u u)/2
254 //
255 // Euler Equations:
256 // drho/dt + div( U ) = 0
257 // dU/dt + div( rho (u x u) + P I3 ) = 0
258 // dE/dt + div( (E + P) u ) = 0
259 //
260 // Equation of State:
261 // P = (gamma - 1) (E - rho (u u) / 2)
262 //
263 // Constants:
264 // cv , Specific heat, constant volume
265 // cp , Specific heat, constant pressure
266 // g , Gravity
267 // gamma = cp / cv, Specific heat ratio
268 // *****************************************************************************
Euler(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)269 CEED_QFUNCTION(Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
270 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
271 const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
272 const CeedScalar(*q_data) = in[2];
273 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
274 CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
275
276 EulerContext context = (EulerContext)ctx;
277 const CeedScalar c_tau = context->c_tau;
278 const CeedScalar gamma = 1.4;
279
280 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
281 // Setup
282 // -- Interp in
283 const CeedScalar rho = q[0][i];
284 const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
285 const CeedScalar E = q[4][i];
286 const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
287 const CeedScalar dU[3][3] = {
288 {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
289 {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
290 {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
291 };
292 const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
293 CeedScalar wdetJ, dXdx[3][3];
294 QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx);
295 // dU/dx
296 CeedScalar drhodx[3] = {0.};
297 CeedScalar dEdx[3] = {0.};
298 CeedScalar dUdx[3][3] = {{0.}};
299 CeedScalar dXdxdXdxT[3][3] = {{0.}};
300 for (CeedInt j = 0; j < 3; j++) {
301 for (CeedInt k = 0; k < 3; k++) {
302 drhodx[j] += drho[k] * dXdx[k][j];
303 dEdx[j] += dE[k] * dXdx[k][j];
304 for (CeedInt l = 0; l < 3; l++) {
305 dUdx[j][k] += dU[j][l] * dXdx[l][k];
306 dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j
307 }
308 }
309 }
310 // Pressure
311 const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic,
312 P = E_internal * (gamma - 1.); // P = pressure
313
314 // The Physics
315 // Zero v and dv so all future terms can safely sum into it
316 for (CeedInt j = 0; j < 5; j++) {
317 v[j][i] = 0.;
318 for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
319 }
320
321 // -- Density
322 // ---- u rho
323 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]);
324 // -- Momentum
325 // ---- rho (u x u) + P I3
326 for (CeedInt j = 0; j < 3; j++) {
327 for (CeedInt k = 0; k < 3; k++) {
328 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] +
329 (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
330 }
331 }
332 // -- Total Energy Density
333 // ---- (E + P) u
334 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]);
335
336 // --Stabilization terms
337 // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
338 CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
339 ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
340
341 // ---- dqdx collects drhodx, dUdx and dEdx in one vector
342 CeedScalar dqdx[5][3];
343 for (CeedInt j = 0; j < 3; j++) {
344 dqdx[0][j] = drhodx[j];
345 dqdx[4][j] = dEdx[j];
346 for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
347 }
348
349 // ---- strong_conv = dF/dq * dq/dx (Strong convection)
350 CeedScalar strong_conv[5] = {0.};
351 for (CeedInt j = 0; j < 3; j++) {
352 for (CeedInt k = 0; k < 5; k++) {
353 for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
354 }
355 }
356
357 // Stabilization
358 // -- Tau elements
359 const CeedScalar sound_speed = sqrt(gamma * P / rho);
360 CeedScalar Tau_x[3] = {0.};
361 Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
362
363 // -- Stabilization method: none or SU
364 CeedScalar stab[5][3] = {{0.}};
365 switch (context->stabilization) {
366 case 0: // Galerkin
367 break;
368 case 1: // SU
369 for (CeedInt j = 0; j < 3; j++) {
370 for (CeedInt k = 0; k < 5; k++) {
371 for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
372 }
373 }
374
375 for (CeedInt j = 0; j < 5; j++) {
376 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]);
377 }
378 break;
379 case 2: // SUPG is not implemented for explicit scheme
380 break;
381 }
382 }
383 return 0;
384 }
385
386 // *****************************************************************************
387 // This QFunction implements the Euler equations with (mentioned above) with implicit time stepping method
388 // *****************************************************************************
IFunction_Euler(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)389 CEED_QFUNCTION(IFunction_Euler)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
390 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
391 const CeedScalar(*dq)[5][CEED_Q_VLA] = (const CeedScalar(*)[5][CEED_Q_VLA])in[1];
392 const CeedScalar(*q_dot)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
393 const CeedScalar(*q_data) = in[3];
394 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
395 CeedScalar(*dv)[5][CEED_Q_VLA] = (CeedScalar(*)[5][CEED_Q_VLA])out[1];
396 CeedScalar *jac_data = out[2];
397
398 EulerContext context = (EulerContext)ctx;
399 const CeedScalar c_tau = context->c_tau;
400 const CeedScalar gamma = 1.4;
401 const CeedScalar zeros[14] = {0.};
402
403 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
404 // Setup
405 // -- Interp in
406 const CeedScalar rho = q[0][i];
407 const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
408 const CeedScalar E = q[4][i];
409 const CeedScalar drho[3] = {dq[0][0][i], dq[1][0][i], dq[2][0][i]};
410 const CeedScalar dU[3][3] = {
411 {dq[0][1][i], dq[1][1][i], dq[2][1][i]},
412 {dq[0][2][i], dq[1][2][i], dq[2][2][i]},
413 {dq[0][3][i], dq[1][3][i], dq[2][3][i]}
414 };
415 const CeedScalar dE[3] = {dq[0][4][i], dq[1][4][i], dq[2][4][i]};
416 CeedScalar wdetJ, dXdx[3][3];
417 QdataUnpack_3D(Q, i, q_data, &wdetJ, dXdx);
418 // dU/dx
419 CeedScalar drhodx[3] = {0.};
420 CeedScalar dEdx[3] = {0.};
421 CeedScalar dUdx[3][3] = {{0.}};
422 CeedScalar dXdxdXdxT[3][3] = {{0.}};
423 for (CeedInt j = 0; j < 3; j++) {
424 for (CeedInt k = 0; k < 3; k++) {
425 drhodx[j] += drho[k] * dXdx[k][j];
426 dEdx[j] += dE[k] * dXdx[k][j];
427 for (CeedInt l = 0; l < 3; l++) {
428 dUdx[j][k] += dU[j][l] * dXdx[l][k];
429 dXdxdXdxT[j][k] += dXdx[j][l] * dXdx[k][l]; // dXdx_j,k * dXdx_k,j
430 }
431 }
432 }
433 const CeedScalar E_kinetic = 0.5 * rho * (u[0] * u[0] + u[1] * u[1] + u[2] * u[2]), E_internal = E - E_kinetic,
434 P = E_internal * (gamma - 1.); // P = pressure
435
436 // The Physics
437 // Zero v and dv so all future terms can safely sum into it
438 for (CeedInt j = 0; j < 5; j++) {
439 v[j][i] = 0.;
440 for (CeedInt k = 0; k < 3; k++) dv[k][j][i] = 0.;
441 }
442 //-----mass matrix
443 for (CeedInt j = 0; j < 5; j++) v[j][i] += wdetJ * q_dot[j][i];
444
445 // -- Density
446 // ---- u rho
447 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]);
448 // -- Momentum
449 // ---- rho (u x u) + P I3
450 for (CeedInt j = 0; j < 3; j++) {
451 for (CeedInt k = 0; k < 3; k++) {
452 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] +
453 (rho * u[j] * u[2] + (j == 2 ? P : 0.)) * dXdx[k][2]);
454 }
455 }
456 // -- Total Energy Density
457 // ---- (E + P) u
458 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]);
459
460 // -- Stabilization terms
461 // ---- jacob_F_conv[3][5][5] = dF(convective)/dq at each direction
462 CeedScalar jacob_F_conv[3][5][5] = {{{0.}}};
463 ConvectiveFluxJacobian_Euler(jacob_F_conv, rho, u, E, gamma);
464
465 // ---- dqdx collects drhodx, dUdx and dEdx in one vector
466 CeedScalar dqdx[5][3];
467 for (CeedInt j = 0; j < 3; j++) {
468 dqdx[0][j] = drhodx[j];
469 dqdx[4][j] = dEdx[j];
470 for (CeedInt k = 0; k < 3; k++) dqdx[k + 1][j] = dUdx[k][j];
471 }
472
473 // ---- strong_conv = dF/dq * dq/dx (Strong convection)
474 CeedScalar strong_conv[5] = {0.};
475 for (CeedInt j = 0; j < 3; j++) {
476 for (CeedInt k = 0; k < 5; k++) {
477 for (CeedInt l = 0; l < 5; l++) strong_conv[k] += jacob_F_conv[j][k][l] * dqdx[l][j];
478 }
479 }
480
481 // ---- Strong residual
482 CeedScalar strong_res[5];
483 for (CeedInt j = 0; j < 5; j++) strong_res[j] = q_dot[j][i] + strong_conv[j];
484
485 // Stabilization
486 // -- Tau elements
487 const CeedScalar sound_speed = sqrt(gamma * P / rho);
488 CeedScalar Tau_x[3] = {0.};
489 Tau_spatial(Tau_x, dXdx, u, sound_speed, c_tau);
490
491 // -- Stabilization method: none, SU, or SUPG
492 CeedScalar stab[5][3] = {{0.}};
493 switch (context->stabilization) {
494 case 0: // Galerkin
495 break;
496 case 1: // SU
497 for (CeedInt j = 0; j < 3; j++) {
498 for (CeedInt k = 0; k < 5; k++) {
499 for (CeedInt l = 0; l < 5; l++) stab[k][j] += jacob_F_conv[j][k][l] * Tau_x[j] * strong_conv[l];
500 }
501 }
502
503 for (CeedInt j = 0; j < 5; j++) {
504 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]);
505 }
506 break;
507 case 2: // SUPG
508 for (CeedInt j = 0; j < 3; j++) {
509 for (CeedInt k = 0; k < 5; k++) {
510 for (CeedInt l = 0; l < 5; l++) stab[k][j] = jacob_F_conv[j][k][l] * Tau_x[j] * strong_res[l];
511 }
512 }
513
514 for (CeedInt j = 0; j < 5; j++) {
515 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]);
516 }
517 break;
518 }
519 StoredValuesPack(Q, i, 0, 14, zeros, jac_data);
520 }
521 return 0;
522 }
523 // *****************************************************************************
524 // This QFunction sets the inflow boundary conditions for the traveling vortex problem.
525 //
526 // Prescribed T_inlet and P_inlet are converted to conservative variables and applied weakly.
527 // *****************************************************************************
TravelingVortex_Inflow(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)528 CEED_QFUNCTION(TravelingVortex_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
529 const CeedScalar(*q_data_sur) = in[2];
530 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
531
532 EulerContext context = (EulerContext)ctx;
533 const int euler_test = context->euler_test;
534 const bool is_implicit = context->implicit;
535 CeedScalar *mean_velocity = context->mean_velocity;
536 const CeedScalar cv = 2.5;
537 const CeedScalar R = 1.;
538 CeedScalar T_inlet;
539 CeedScalar P_inlet;
540
541 // For test cases 1 and 3 the background velocity is zero
542 if (euler_test == 1 || euler_test == 3) {
543 for (CeedInt i = 0; i < 3; i++) mean_velocity[i] = 0.;
544 }
545
546 // For test cases 1 and 2, T_inlet = T_inlet = 0.4
547 if (euler_test == 1 || euler_test == 2) T_inlet = P_inlet = .4;
548 else T_inlet = P_inlet = 1.;
549
550 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
551 CeedScalar wdetJb, norm[3];
552 QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm);
553 wdetJb *= is_implicit ? -1. : 1.;
554
555 // face_normal = Normal vector of the face
556 const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
557 // The Physics
558 // Zero v so all future terms can safely sum into it
559 for (CeedInt j = 0; j < 5; j++) v[j][i] = 0.;
560
561 // Implementing in/outflow BCs
562 if (face_normal > 0) {
563 } else { // inflow
564 const CeedScalar rho_inlet = P_inlet / (R * T_inlet);
565 const CeedScalar E_kinetic_inlet = (mean_velocity[0] * mean_velocity[0] + mean_velocity[1] * mean_velocity[1]) / 2.;
566 // incoming total energy
567 const CeedScalar E_inlet = rho_inlet * (cv * T_inlet + E_kinetic_inlet);
568
569 // The Physics
570 // -- Density
571 v[0][i] -= wdetJb * rho_inlet * face_normal;
572
573 // -- Momentum
574 for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho_inlet * face_normal * mean_velocity[j] + norm[j] * P_inlet);
575
576 // -- Total Energy Density
577 v[4][i] -= wdetJb * face_normal * (E_inlet + P_inlet);
578 }
579 }
580 return 0;
581 }
582
583 // *****************************************************************************
584 // This QFunction sets the outflow boundary conditions for the Euler solver.
585 //
586 // Outflow BCs:
587 // The validity of the weak form of the governing equations is extended to the outflow.
588 // *****************************************************************************
Euler_Outflow(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)589 CEED_QFUNCTION(Euler_Outflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
590 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
591 const CeedScalar(*q_data_sur) = in[2];
592 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
593
594 EulerContext context = (EulerContext)ctx;
595 const bool is_implicit = context->implicit;
596 CeedScalar *mean_velocity = context->mean_velocity;
597
598 const CeedScalar gamma = 1.4;
599
600 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
601 // Setup
602 // -- Interp in
603 const CeedScalar rho = q[0][i];
604 const CeedScalar u[3] = {q[1][i] / rho, q[2][i] / rho, q[3][i] / rho};
605 const CeedScalar E = q[4][i];
606
607 CeedScalar wdetJb, norm[3];
608 QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm);
609 wdetJb *= is_implicit ? -1. : 1.;
610
611 // face_normal = Normal vector of the face
612 const CeedScalar face_normal = norm[0] * mean_velocity[0] + norm[1] * mean_velocity[1] + norm[2] * mean_velocity[2];
613 // The Physics
614 // Zero v so all future terms can safely sum into it
615 for (CeedInt j = 0; j < 5; j++) v[j][i] = 0;
616
617 // Implementing in/outflow BCs
618 if (face_normal > 0) { // outflow
619 const CeedScalar E_kinetic = (u[0] * u[0] + u[1] * u[1]) / 2.;
620 const CeedScalar P = (E - E_kinetic * rho) * (gamma - 1.); // pressure
621 const CeedScalar u_normal = norm[0] * u[0] + norm[1] * u[1] + norm[2] * u[2]; // Normal velocity
622 // The Physics
623 // -- Density
624 v[0][i] -= wdetJb * rho * u_normal;
625
626 // -- Momentum
627 for (CeedInt j = 0; j < 3; j++) v[j + 1][i] -= wdetJb * (rho * u_normal * u[j] + norm[j] * P);
628
629 // -- Total Energy Density
630 v[4][i] -= wdetJb * u_normal * (E + P);
631 }
632 }
633 return 0;
634 }
635