-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathboundary_jacobians.hpp
590 lines (475 loc) · 17.5 KB
/
boundary_jacobians.hpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
/* Copyright (C) 2017 Jerry Watkins
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program. If not, see <http://www.gnu.org/licenses/>.
*/
/* Slip-wall and Symmetry Boundary Condition */
template <size_t nVars, size_t nDims>
#ifdef _GPU
__device__ __forceinline__
#endif
void compute_Slip_Jac(double U[nVars], double Ub[nVars], double dUbdU[nVars][nVars],
double norm[nDims])
{
if (nDims == 2)
{
double nx = norm[0];
double ny = norm[1];
/* Primitive Variables */
double uL = U[1] / U[0];
double vL = U[2] / U[0];
double uR = Ub[1] / Ub[0];
double vR = Ub[2] / Ub[0];
/* Compute dUbdU */
dUbdU[0][0] = 1;
dUbdU[3][0] = 0.5 * (uL*uL + vL*vL - uR*uR - vR*vR);
dUbdU[1][1] = 1.0-nx*nx;
dUbdU[2][1] = -nx*ny;
dUbdU[3][1] = -uL + (1.0-nx*nx)*uR - nx*ny*vR;
dUbdU[1][2] = -nx*ny;
dUbdU[2][2] = 1.0-ny*ny;
dUbdU[3][2] = -vL - nx*ny*uR + (1.0-ny*ny)*vR;
dUbdU[3][3] = 1;
}
else if (nDims == 3)
{
double nx = norm[0];
double ny = norm[1];
double nz = norm[2];
/* Primitive Variables */
double uL = U[1] / U[0];
double vL = U[2] / U[0];
double wL = U[3] / U[0];
double uR = Ub[1] / Ub[0];
double vR = Ub[2] / Ub[0];
double wR = Ub[3] / Ub[0];
/* Compute dUbdU */
dUbdU[0][0] = 1;
dUbdU[4][0] = 0.5 * (uL*uL + vL*vL + wL*wL - uR*uR - vR*vR - wR*wR);
dUbdU[1][1] = 1.0-nx*nx;
dUbdU[2][1] = -nx*ny;
dUbdU[3][1] = -nx*nz;
dUbdU[4][1] = -uL + (1.0-nx*nx)*uR - nx*ny*vR - nx*nz*wR;
dUbdU[1][2] = -nx*ny;
dUbdU[2][2] = 1.0-ny*ny;
dUbdU[3][2] = -ny*nz;
dUbdU[4][2] = -vL - nx*ny*uR + (1.0-ny*ny)*vR - ny*nz*wR;
dUbdU[1][3] = -nx*nz;
dUbdU[2][3] = -ny*nz;
dUbdU[3][3] = 1.0-nz*nz;
dUbdU[4][3] = -wL - nx*nz*uR - ny*nz*vR + (1.0-nz*nz)*wR;
dUbdU[4][4] = 1;
}
}
/* No-Slip Adiabatic Wall */
template <size_t nVars, size_t nDims>
#ifdef _GPU
__device__ __forceinline__
#endif
void compute_NoSlip_Jac(double U[nVars], double Q[nVars][nDims], double dUbdU[nVars][nVars],
double dQbdU[nDims][nVars][nVars], double dQbdQ[nDims][nDims][nVars][nVars],
double norm[nDims])
{
/* Matrix Parameters */
double astar[nDims][nVars-1] = {0};
for (unsigned int dimi = 0; dimi < nDims; dimi++)
for (unsigned int dimj = 0; dimj < nDims; dimj++)
for (unsigned int var = 0; var < nVars-1; var++)
astar[dimi][var] += norm[dimi] * norm[dimj] * Q[var][dimj] / U[0];
double momF = 0.0;
for (unsigned int dim = 0; dim < nDims; dim++)
momF += U[dim+1] * U[dim+1];
momF /= U[0];
double estar1 = (U[nDims+1] - 2.0 * momF) / U[0];
double estar2 = (U[nDims+1] - momF) / U[0];
if (nDims == 2)
{
double nx = norm[0];
double ny = norm[1];
/* Primitive Variables */
double uL = U[1] / U[0];
double vL = U[2] / U[0];
/* Compute dUbdU */
dUbdU[0][0] = 1;
dUbdU[3][0] = 0.5 * (uL*uL + vL*vL);
dUbdU[3][1] = -uL;
dUbdU[3][2] = -vL;
dUbdU[3][3] = 1;
/* Compute dQbdU */
for (unsigned int dim = 0; dim < nDims; dim++)
{
dQbdU[dim][3][0] = estar1 * astar[dim][0] + uL * astar[dim][1] + vL * astar[dim][2];
dQbdU[dim][3][1] = astar[dim][1] - 2.0 * uL * astar[dim][0];
dQbdU[dim][3][2] = astar[dim][2] - 2.0 * vL * astar[dim][0];
dQbdU[dim][3][3] = astar[dim][0];
}
/* Compute dQxR/dQxL */
dQbdQ[0][0][0][0] = 1;
dQbdQ[0][0][3][0] = nx*nx * estar2;
dQbdQ[0][0][1][1] = 1;
dQbdQ[0][0][3][1] = nx*nx * uL;
dQbdQ[0][0][2][2] = 1;
dQbdQ[0][0][3][2] = nx*nx * vL;
dQbdQ[0][0][3][3] = 1.0 - nx*nx;
/* Compute dQyR/dQxL */
dQbdQ[1][0][3][0] = nx*ny * estar2;
dQbdQ[1][0][3][1] = nx*ny * uL;
dQbdQ[1][0][3][2] = nx*ny * vL;
dQbdQ[1][0][3][3] = -nx * ny;
/* Compute dQxR/dQyL */
dQbdQ[0][1][3][0] = nx*ny * estar2;
dQbdQ[0][1][3][1] = nx*ny * uL;
dQbdQ[0][1][3][2] = nx*ny * vL;
dQbdQ[0][1][3][3] = -nx * ny;
/* Compute dQyR/dQyL */
dQbdQ[1][1][0][0] = 1;
dQbdQ[1][1][3][0] = ny*ny * estar2;
dQbdQ[1][1][1][1] = 1;
dQbdQ[1][1][3][1] = ny*ny * uL;
dQbdQ[1][1][2][2] = 1;
dQbdQ[1][1][3][2] = ny*ny * vL;
dQbdQ[1][1][3][3] = 1.0 - ny*ny;
}
else if (nDims == 3)
{
double nx = norm[0];
double ny = norm[1];
double nz = norm[2];
/* Primitive Variables */
double uL = U[1] / U[0];
double vL = U[2] / U[0];
double wL = U[3] / U[0];
/* Compute dUbdU */
dUbdU[0][0] = 1;
dUbdU[4][0] = 0.5 * (uL*uL + vL*vL + wL*wL);
dUbdU[4][1] = -uL;
dUbdU[4][2] = -vL;
dUbdU[4][3] = -wL;
dUbdU[4][4] = 1;
/* Compute dQbdU */
for (unsigned int dim = 0; dim < nDims; dim++)
{
dQbdU[dim][4][0] = estar1 * astar[dim][0] + uL * astar[dim][1] +
vL * astar[dim][2] + wL * astar[dim][3];
dQbdU[dim][4][1] = astar[dim][1] - 2.0 * uL * astar[dim][0];
dQbdU[dim][4][2] = astar[dim][2] - 2.0 * vL * astar[dim][0];
dQbdU[dim][4][3] = astar[dim][3] - 2.0 * wL * astar[dim][0];
dQbdU[dim][4][4] = astar[dim][0];
}
/* Compute dQxR/dQxL */
dQbdQ[0][0][0][0] = 1;
dQbdQ[0][0][4][0] = nx*nx * estar2;
dQbdQ[0][0][1][1] = 1;
dQbdQ[0][0][4][1] = nx*nx * uL;
dQbdQ[0][0][2][2] = 1;
dQbdQ[0][0][4][2] = nx*nx * vL;
dQbdQ[0][0][3][3] = 1;
dQbdQ[0][0][4][3] = nx*nx * wL;
dQbdQ[0][0][4][4] = 1.0 - nx*nx;
/* Compute dQyR/dQxL */
dQbdQ[1][0][4][0] = nx*ny * estar2;
dQbdQ[1][0][4][1] = nx*ny * uL;
dQbdQ[1][0][4][2] = nx*ny * vL;
dQbdQ[1][0][4][3] = nx*ny * wL;
dQbdQ[1][0][4][4] = -nx * ny;
/* Compute dQzR/dQxL */
dQbdQ[2][0][4][0] = nx*nz * estar2;
dQbdQ[2][0][4][1] = nx*nz * uL;
dQbdQ[2][0][4][2] = nx*nz * vL;
dQbdQ[2][0][4][3] = nx*nz * wL;
dQbdQ[2][0][4][4] = -nx * nz;
/* Compute dQxR/dQyL */
dQbdQ[0][1][4][0] = nx*ny * estar2;
dQbdQ[0][1][4][1] = nx*ny * uL;
dQbdQ[0][1][4][2] = nx*ny * vL;
dQbdQ[0][1][4][3] = nx*ny * wL;
dQbdQ[0][1][4][4] = -nx * ny;
/* Compute dQyR/dQyL */
dQbdQ[1][1][0][0] = 1;
dQbdQ[1][1][4][0] = ny*ny * estar2;
dQbdQ[1][1][1][1] = 1;
dQbdQ[1][1][4][1] = ny*ny * uL;
dQbdQ[1][1][2][2] = 1;
dQbdQ[1][1][4][2] = ny*ny * vL;
dQbdQ[1][1][3][3] = 1;
dQbdQ[1][1][4][3] = ny*ny * wL;
dQbdQ[1][1][4][4] = 1.0 - ny*ny;
/* Compute dQzR/dQyL */
dQbdQ[2][1][4][0] = nz*ny * estar2;
dQbdQ[2][1][4][1] = nz*ny * uL;
dQbdQ[2][1][4][2] = nz*ny * vL;
dQbdQ[2][1][4][3] = nz*ny * wL;
dQbdQ[2][1][4][4] = -nz * ny;
/* Compute dQxR/dQzL */
dQbdQ[0][2][4][0] = nx*nz * estar2;
dQbdQ[0][2][4][1] = nx*nz * uL;
dQbdQ[0][2][4][2] = nx*nz * vL;
dQbdQ[0][2][4][3] = nx*nz * wL;
dQbdQ[0][2][4][4] = -nx * nz;
/* Compute dQyR/dQzL */
dQbdQ[1][2][4][0] = ny*nz * estar2;
dQbdQ[1][2][4][1] = ny*nz * uL;
dQbdQ[1][2][4][2] = ny*nz * vL;
dQbdQ[1][2][4][3] = ny*nz * wL;
dQbdQ[1][2][4][4] = -ny * nz;
/* Compute dQzR/dQzL */
dQbdQ[2][2][0][0] = 1;
dQbdQ[2][2][4][0] = nz*nz * estar2;
dQbdQ[2][2][1][1] = 1;
dQbdQ[2][2][4][1] = nz*nz * uL;
dQbdQ[2][2][2][2] = 1;
dQbdQ[2][2][4][2] = nz*nz * vL;
dQbdQ[2][2][3][3] = 1;
dQbdQ[2][2][4][3] = nz*nz * wL;
dQbdQ[2][2][4][4] = 1.0 - nz*nz;
}
}
/* Characteristic Riemann Invariant Far Field */
template <size_t nVars, size_t nDims>
#ifdef _GPU
__device__ __forceinline__
#endif
void compute_Char_Jac(double U[nVars], double Ub[nVars], double dUbdU[nVars][nVars],
double dQbdQ[nDims][nDims][nVars][nVars], double norm[nDims], double V_fs[nDims]
double rho_fs, double P_fs, double gamma, bool viscous)
{
/* Compute wall normal velocities */
double VnL = 0.0; double VnR = 0.0;
for (unsigned int dim = 0; dim < nDims; dim++)
{
VnL += U[dim+1] / U[0] * norm[dim];
VnR += V_fs[dim] * norm[dim];
}
/* Compute pressure */
double momF = 0.0;
for (unsigned int dim = 0; dim < nDims; dim++)
{
momF += U[dim + 1] * U[dim + 1];
}
momF /= U[0];
double PL = (gamma - 1.0) * (U[nDims + 1] - 0.5 * momF);
double PR = P_fs;
double cL = std::sqrt(gamma * PL / U[0]);
double cR = std::sqrt(gamma * PR / rho_fs);
/* Compute Riemann Invariants */
double RL = VnL + 2.0 / (gamma - 1) * cL;
double RB = VnR - 2.0 / (gamma - 1) * cR;
double cstar = 0.25 * (gamma - 1) * (RL - RB);
double ustarn = 0.5 * (RL + RB);
if (nDims == 2)
{
double nx = norm[0];
double ny = norm[1];
double gam = gamma;
/* Primitive Variables */
double rhoL = U[0];
double uL = U[1] / U[0];
double vL = U[2] / U[0];
double rhoR = Ub[0];
double uR = Ub[1] / Ub[0];
double vR = Ub[2] / Ub[0];
if (VnL < 0.0) /* Case 1: Inflow */
{
/* Matrix Parameters */
double a1 = 0.5 * rhoR / cstar;
double a2 = gam / (rhoL * cL);
double b1 = -VnL / rhoL - a2 / rhoL * (PL / (gam-1.0) - 0.5 * momF);
double b2 = nx / rhoL - a2 * uL;
double b3 = ny / rhoL - a2 * vL;
double b4 = a2 / cstar;
double c1 = cstar * cstar / ((gam-1.0) * gam) + 0.5 * (uR*uR + vR*vR);
double c2 = uR * nx + vR * ny + cstar / gam;
/* Compute dUbdU */
dUbdU[0][0] = a1 * b1;
dUbdU[1][0] = a1 * b1 * uR + 0.5 * rhoR * b1 * nx;
dUbdU[2][0] = a1 * b1 * vR + 0.5 * rhoR * b1 * ny;
dUbdU[3][0] = a1 * b1 * c1 + 0.5 * rhoR * b1 * c2;
dUbdU[0][1] = a1 * b2;
dUbdU[1][1] = a1 * b2 * uR + 0.5 * rhoR * b2 * nx;
dUbdU[2][1] = a1 * b2 * vR + 0.5 * rhoR * b2 * ny;
dUbdU[3][1] = a1 * b2 * c1 + 0.5 * rhoR * b2 * c2;
dUbdU[0][2] = a1 * b3;
dUbdU[1][2] = a1 * b3 * uR + 0.5 * rhoR * b3 * nx;
dUbdU[2][2] = a1 * b3 * vR + 0.5 * rhoR * b3 * ny;
dUbdU[3][2] = a1 * b3 * c1 + 0.5 * rhoR * b3 * c2;
dUbdU[0][3] = 0.5 * rhoR * b4;
dUbdU[1][3] = 0.5 * rhoR * (b4 * uR + a2 * nx);
dUbdU[2][3] = 0.5 * rhoR * (b4 * vR + a2 * ny);
dUbdU[3][3] = 0.5 * rhoR * (b4 * c1 + a2 * c2);
}
else /* Case 2: Outflow */
{
/* Matrix Parameters */
double a1 = gam * rhoR / (gam-1.0);
double a2 = gam / (rhoL * cL);
double a3 = (gam-1.0) / (gam * PL);
double a4 = (gam-1.0) / (2.0 * gam * cstar);
double a5 = rhoR * cstar * cstar / (gam-1.0) / (gam-1.0);
double a6 = rhoR * cstar / (2.0 * gam);
double b1 = -VnL / rhoL - a2 / rhoL * (PL / (gam-1.0) - 0.5 * momF);
double b2 = nx / rhoL - a2 * uL;
double b3 = ny / rhoL - a2 * vL;
double c1 = 0.5 * b1 * nx - (-VnL * nx + uL) / rhoL;
double c2 = 0.5 * b2 * nx + (1.0 - nx*nx) / rhoL;
double c3 = 0.5 * b3 * nx - nx * ny / rhoL;
double c4 = ustarn * nx + uL - VnL * nx;
double d1 = 0.5 * b1 * ny - (-VnL * ny + vL) / rhoL;
double d2 = 0.5 * b2 * ny - nx * ny / rhoL;
double d3 = 0.5 * b3 * ny + (1.0 - ny*ny) / rhoL;
double d4 = ustarn * ny + vL - VnL * ny;
double e1 = 1.0 / rhoL - 0.5 * a3 * momF / rhoL + a4 * b1;
double e2 = a3 * uL + a4 * b2;
double e3 = a3 * vL + a4 * b3;
double e4 = -a3 + a2 * a4;
double f1 = 0.5 * a1 * (c4*c4 + d4*d4) + a5;
/* Compute dUbdU */
dUbdU[0][0] = a1 * e1;
dUbdU[1][0] = a1 * e1 * c4 + rhoR * c1;
dUbdU[2][0] = a1 * e1 * d4 + rhoR * d1;
dUbdU[3][0] = rhoR * (c1*c4 + d1*d4) + e1 * f1 + a6 * b1;
dUbdU[0][1] = a1 * e2;
dUbdU[1][1] = a1 * e2 * c4 + rhoR * c2;
dUbdU[2][1] = a1 * e2 * d4 + rhoR * d2;
dUbdU[3][1] = rhoR * (c2*c4 + d2*d4) + e2 * f1 + a6 * b2;
dUbdU[0][2] = a1 * e3;
dUbdU[1][2] = a1 * e3 * c4 + rhoR * c3;
dUbdU[2][2] = a1 * e3 * d4 + rhoR * d3;
dUbdU[3][2] = rhoR * (c3*c4 + d3*d4) + e3 * f1 + a6 * b3;
dUbdU[0][3] = a1 * e4;
dUbdU[1][3] = a1 * e4 * c4 + 0.5 * rhoR * a2 * nx;
dUbdU[2][3] = a1 * e4 * d4 + 0.5 * rhoR * a2 * ny;
dUbdU[3][3] = 0.5 * rhoR * a2 * (c4*nx + d4*ny) + e4 * f1 + a2 * a6;
}
}
else if (nDims == 3)
{
double nx = norm[0];
double ny = norm[1];
double nz = norm[2];
double gam = gamma;
/* Primitive Variables */
double rhoL = U[0];
double uL = U[1] / U[0];
double vL = U[2] / U[0];
double wL = U[3] / U[0];
double rhoR = Ub[0];
double uR = Ub[1] / Ub[0];
double vR = Ub[2] / Ub[0];
double wR = Ub[3] / Ub[0];
if (VnL < 0.0) /* Case 1: Inflow */
{
/* Matrix Parameters */
double a1 = 0.5 * rhoR / cstar;
double a2 = gam / (rhoL * cL);
double b1 = -VnL / rhoL - a2 / rhoL * (PL / (gam-1.0) - 0.5 * momF);
double b2 = nx / rhoL - a2 * uL;
double b3 = ny / rhoL - a2 * vL;
double b4 = nz / rhoL - a2 * wL;
double b5 = a2 / cstar;
double c1 = cstar * cstar / ((gam-1.0) * gam) + 0.5 * (uR*uR + vR*vR + wR*wR);
double c2 = uR * nx + vR * ny + wR * nz + cstar / gam;
/* Compute dUbdU */
dUbdU[0][0] = a1 * b1;
dUbdU[1][0] = a1 * b1 * uR + 0.5 * rhoR * b1 * nx;
dUbdU[2][0] = a1 * b1 * vR + 0.5 * rhoR * b1 * ny;
dUbdU[3][0] = a1 * b1 * wR + 0.5 * rhoR * b1 * nz;
dUbdU[4][0] = a1 * b1 * c1 + 0.5 * rhoR * b1 * c2;
dUbdU[0][1] = a1 * b2;
dUbdU[1][1] = a1 * b2 * uR + 0.5 * rhoR * b2 * nx;
dUbdU[2][1] = a1 * b2 * vR + 0.5 * rhoR * b2 * ny;
dUbdU[3][1] = a1 * b2 * wR + 0.5 * rhoR * b2 * nz;
dUbdU[4][1] = a1 * b2 * c1 + 0.5 * rhoR * b2 * c2;
dUbdU[0][2] = a1 * b3;
dUbdU[1][2] = a1 * b3 * uR + 0.5 * rhoR * b3 * nx;
dUbdU[2][2] = a1 * b3 * vR + 0.5 * rhoR * b3 * ny;
dUbdU[3][2] = a1 * b3 * wR + 0.5 * rhoR * b3 * nz;
dUbdU[4][2] = a1 * b3 * c1 + 0.5 * rhoR * b3 * c2;
dUbdU[0][3] = a1 * b4;
dUbdU[1][3] = a1 * b4 * uR + 0.5 * rhoR * b4 * nx;
dUbdU[2][3] = a1 * b4 * vR + 0.5 * rhoR * b4 * ny;
dUbdU[3][3] = a1 * b4 * wR + 0.5 * rhoR * b4 * nz;
dUbdU[4][3] = a1 * b4 * c1 + 0.5 * rhoR * b4 * c2;
dUbdU[0][4] = 0.5 * rhoR * b5;
dUbdU[1][4] = 0.5 * rhoR * (b5 * uR + a2 * nx);
dUbdU[2][4] = 0.5 * rhoR * (b5 * vR + a2 * ny);
dUbdU[3][4] = 0.5 * rhoR * (b5 * wR + a2 * nz);
dUbdU[4][4] = 0.5 * rhoR * (b5 * c1 + a2 * c2);
}
else /* Case 2: Outflow */
{
/* Matrix Parameters */
double a1 = gam * rhoR / (gam-1.0);
double a2 = gam / (rhoL * cL);
double a3 = (gam-1.0) / (gam * PL);
double a4 = (gam-1.0) / (2.0 * gam * cstar);
double a5 = rhoR * cstar * cstar / (gam-1.0) / (gam-1.0);
double a6 = rhoR * cstar / (2.0 * gam);
double b1 = -VnL / rhoL - a2 / rhoL * (PL / (gam-1.0) - 0.5 * momF);
double b2 = nx / rhoL - a2 * uL;
double b3 = ny / rhoL - a2 * vL;
double b4 = nz / rhoL - a2 * wL;
double c1 = 0.5 * b1 * nx - (-VnL * nx + uL) / rhoL;
double c2 = 0.5 * b2 * nx + (1.0 - nx*nx) / rhoL;
double c3 = 0.5 * b3 * nx - nx * ny / rhoL;
double c4 = 0.5 * b4 * nx - nx * nz / rhoL;
double c5 = ustarn * nx + uL - VnL * nx;
double d1 = 0.5 * b1 * ny - (-VnL * ny + vL) / rhoL;
double d2 = 0.5 * b2 * ny - nx * ny / rhoL;
double d3 = 0.5 * b3 * ny + (1.0 - ny*ny) / rhoL;
double d4 = 0.5 * b4 * ny - ny * nz / rhoL;
double d5 = ustarn * ny + vL - VnL * ny;
double e1 = 0.5 * b1 * nz - (-VnL * nz + wL) / rhoL;
double e2 = 0.5 * b2 * nz - nx * nz / rhoL;
double e3 = 0.5 * b3 * nz - ny * nz / rhoL;
double e4 = 0.5 * b4 * nz + (1.0 - nz*nz) / rhoL;
double e5 = ustarn * nz + wL - VnL * nz;
double f1 = 1.0 / rhoL - 0.5 * a3 * momF / rhoL + a4 * b1;
double f2 = a3 * uL + a4 * b2;
double f3 = a3 * vL + a4 * b3;
double f4 = a3 * wL + a4 * b4;
double f5 = -a3 + a2 * a4;
double g1 = 0.5 * a1 * (c5*c5 + d5*d5 + e5*e5) + a5;
/* Compute dUbdU */
dUbdU[0][0] = a1 * f1;
dUbdU[1][0] = a1 * f1 * c5 + rhoR * c1;
dUbdU[2][0] = a1 * f1 * d5 + rhoR * d1;
dUbdU[3][0] = a1 * f1 * e5 + rhoR * e1;
dUbdU[4][0] = rhoR * (c1*c5 + d1*d5 + e1*e5) + f1 * g1 + a6 * b1;
dUbdU[0][1] = a1 * f2;
dUbdU[1][1] = a1 * f2 * c5 + rhoR * c2;
dUbdU[2][1] = a1 * f2 * d5 + rhoR * d2;
dUbdU[3][1] = a1 * f2 * e5 + rhoR * e2;
dUbdU[4][1] = rhoR * (c2*c5 + d2*d5 + e2*e5) + f2 * g1 + a6 * b2;
dUbdU[0][2] = a1 * f3;
dUbdU[1][2] = a1 * f3 * c5 + rhoR * c3;
dUbdU[2][2] = a1 * f3 * d5 + rhoR * d3;
dUbdU[3][2] = a1 * f3 * e5 + rhoR * e3;
dUbdU[4][2] = rhoR * (c3*c5 + d3*d5 + e3*e5) + f3 * g1 + a6 * b3;
dUbdU[0][3] = a1 * f4;
dUbdU[1][3] = a1 * f4 * c5 + rhoR * c4;
dUbdU[2][3] = a1 * f4 * d5 + rhoR * d4;
dUbdU[3][3] = a1 * f4 * e5 + rhoR * e4;
dUbdU[4][3] = rhoR * (c4*c5 + d4*d5 + e4*e5) + f4 * g1 + a6 * b4;
dUbdU[0][4] = a1 * f5;
dUbdU[1][4] = a1 * f5 * c5 + 0.5 * rhoR * a2 * nx;
dUbdU[2][4] = a1 * f5 * d5 + 0.5 * rhoR * a2 * ny;
dUbdU[3][4] = a1 * f5 * e5 + 0.5 * rhoR * a2 * nz;
dUbdU[4][4] = 0.5 * rhoR * a2 * (c5*nx + d5*ny + e5*nz) + f5 * g1 + a2 * a6;
}
}
/* Extrapolate gradients */
if (viscous)
for (unsigned int dim = 0; dim < nDims; dim++)
for (unsigned int var = 0; var < nVars; var++)
dQbdQ[dim][dim][var][var] = 1;
}