-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathex18.hpp
490 lines (403 loc) · 13.3 KB
/
ex18.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
// MFEM Example 18 - Serial/Parallel Shared Code
#include "mfem.hpp"
using namespace std;
using namespace mfem;
// Problem definition
extern int problem;
// Maximum characteristic speed (updated by integrators)
extern double max_char_speed;
extern const int num_equation;
extern const double specific_heat_ratio;
extern const double gas_constant;
// Time-dependent operator for the right-hand side of the ODE representing the
// DG weak form.
class FE_Evolution : public TimeDependentOperator
{
private:
const int dim;
FiniteElementSpace &vfes;
Operator &A;
SparseMatrix &Aflux;
DenseTensor Me_inv;
mutable Vector state;
mutable DenseMatrix f;
mutable DenseTensor flux;
mutable Vector z;
void GetFlux(const DenseMatrix &state_, DenseTensor &flux_) const;
public:
FE_Evolution(FiniteElementSpace &vfes_,
Operator &A_, SparseMatrix &Aflux_);
virtual void Mult(const Vector &x, Vector &y) const;
virtual ~FE_Evolution() { }
};
// Implements a simple Rusanov flux
class RiemannSolver
{
private:
Vector flux1;
Vector flux2;
public:
RiemannSolver();
double Eval(const Vector &state1, const Vector &state2,
const Vector &nor, Vector &flux);
};
// Interior face term: <F.n(u),[w]>
class FaceIntegrator : public NonlinearFormIntegrator
{
private:
RiemannSolver rsolver;
Vector shape1;
Vector shape2;
Vector funval1;
Vector funval2;
Vector nor;
Vector fluxN;
public:
FaceIntegrator(RiemannSolver &rsolver_, const int dim);
virtual void AssembleFaceVector(const FiniteElement &el1,
const FiniteElement &el2,
FaceElementTransformations &Tr,
const Vector &elfun, Vector &elvect);
};
// Implementation of class FE_Evolution
FE_Evolution::FE_Evolution(FiniteElementSpace &vfes_,
Operator &A_, SparseMatrix &Aflux_)
: TimeDependentOperator(A_.Height()),
dim(vfes_.GetFE(0)->GetDim()),
vfes(vfes_),
A(A_),
Aflux(Aflux_),
Me_inv(vfes.GetFE(0)->GetDof(), vfes.GetFE(0)->GetDof(), vfes.GetNE()),
state(num_equation),
f(num_equation, dim),
flux(vfes.GetNDofs(), dim, num_equation),
z(A.Height())
{
// Standard local assembly and inversion for energy mass matrices.
const int dof = vfes.GetFE(0)->GetDof();
DenseMatrix Me(dof);
DenseMatrixInverse inv(&Me);
MassIntegrator mi;
for (int i = 0; i < vfes.GetNE(); i++)
{
mi.AssembleElementMatrix(*vfes.GetFE(i), *vfes.GetElementTransformation(i), Me);
inv.Factor();
inv.GetInverseMatrix(Me_inv(i));
}
}
void FE_Evolution::Mult(const Vector &x, Vector &y) const
{
// 0. Reset wavespeed computation before operator application.
max_char_speed = 0.;
// 1. Create the vector z with the face terms -<F.n(u), [w]>.
A.Mult(x, z);
// 2. Add the element terms.
// i. computing the flux approximately as a grid function by interpolating
// at the solution nodes.
// ii. multiplying this grid function by a (constant) mixed bilinear form for
// each of the num_equation, computing (F(u), grad(w)) for each equation.
DenseMatrix xmat(x.GetData(), vfes.GetNDofs(), num_equation);
GetFlux(xmat, flux);
for (int k = 0; k < num_equation; k++)
{
Vector fk(flux(k).GetData(), dim * vfes.GetNDofs());
Vector zk(z.GetData() + k * vfes.GetNDofs(), vfes.GetNDofs());
Aflux.AddMult(fk, zk);
}
// 3. Multiply element-wise by the inverse mass matrices.
Vector zval;
Array<int> vdofs;
const int dof = vfes.GetFE(0)->GetDof();
DenseMatrix zmat, ymat(dof, num_equation);
for (int i = 0; i < vfes.GetNE(); i++)
{
// Return the vdofs ordered byNODES
vfes.GetElementVDofs(i, vdofs);
z.GetSubVector(vdofs, zval);
zmat.UseExternalData(zval.GetData(), dof, num_equation);
mfem::Mult(Me_inv(i), zmat, ymat);
y.SetSubVector(vdofs, ymat.GetData());
}
}
// Physicality check (at end)
bool StateIsPhysical(const Vector &state, const int dim);
// Pressure (EOS) computation
inline double ComputePressure(const Vector &state, int dim)
{
const double den = state(0);
const Vector den_vel(state.GetData() + 1, dim);
const double den_energy = state(1 + dim);
double den_vel2 = 0;
for (int d = 0; d < dim; d++) { den_vel2 += den_vel(d) * den_vel(d); }
den_vel2 /= den;
return (specific_heat_ratio - 1.0) * (den_energy - 0.5 * den_vel2);
}
// Compute the vector flux F(u)
void ComputeFlux(const Vector &state, int dim, DenseMatrix &flux)
{
const double den = state(0);
const Vector den_vel(state.GetData() + 1, dim);
const double den_energy = state(1 + dim);
MFEM_ASSERT(StateIsPhysical(state, dim), "");
const double pres = ComputePressure(state, dim);
for (int d = 0; d < dim; d++)
{
flux(0, d) = den_vel(d);
for (int i = 0; i < dim; i++)
{
flux(1+i, d) = den_vel(i) * den_vel(d) / den;
}
flux(1+d, d) += pres;
}
const double H = (den_energy + pres) / den;
for (int d = 0; d < dim; d++)
{
flux(1+dim, d) = den_vel(d) * H;
}
}
// Compute the scalar F(u).n
void ComputeFluxDotN(const Vector &state, const Vector &nor,
Vector &fluxN)
{
// NOTE: nor in general is not a unit normal
const int dim = nor.Size();
const double den = state(0);
const Vector den_vel(state.GetData() + 1, dim);
const double den_energy = state(1 + dim);
MFEM_ASSERT(StateIsPhysical(state, dim), "");
const double pres = ComputePressure(state, dim);
double den_velN = 0;
for (int d = 0; d < dim; d++) { den_velN += den_vel(d) * nor(d); }
fluxN(0) = den_velN;
for (int d = 0; d < dim; d++)
{
fluxN(1+d) = den_velN * den_vel(d) / den + pres * nor(d);
}
const double H = (den_energy + pres) / den;
fluxN(1 + dim) = den_velN * H;
}
// Compute the maximum characteristic speed.
inline double ComputeMaxCharSpeed(const Vector &state, const int dim)
{
const double den = state(0);
const Vector den_vel(state.GetData() + 1, dim);
double den_vel2 = 0;
for (int d = 0; d < dim; d++) { den_vel2 += den_vel(d) * den_vel(d); }
den_vel2 /= den;
const double pres = ComputePressure(state, dim);
const double sound = sqrt(specific_heat_ratio * pres / den);
const double vel = sqrt(den_vel2 / den);
return vel + sound;
}
// Compute the flux at solution nodes.
void FE_Evolution::GetFlux(const DenseMatrix &x_, DenseTensor &flux_) const
{
const int flux_dof = flux_.SizeI();
const int flux_dim = flux_.SizeJ();
for (int i = 0; i < flux_dof; i++)
{
for (int k = 0; k < num_equation; k++) { state(k) = x_(i, k); }
ComputeFlux(state, flux_dim, f);
for (int d = 0; d < flux_dim; d++)
{
for (int k = 0; k < num_equation; k++)
{
flux_(i, d, k) = f(k, d);
}
}
// Update max char speed
const double mcs = ComputeMaxCharSpeed(state, flux_dim);
if (mcs > max_char_speed) { max_char_speed = mcs; }
}
}
// Implementation of class RiemannSolver
RiemannSolver::RiemannSolver() :
flux1(num_equation),
flux2(num_equation) { }
double RiemannSolver::Eval(const Vector &state1, const Vector &state2,
const Vector &nor, Vector &flux)
{
// NOTE: nor in general is not a unit normal
const int dim = nor.Size();
MFEM_ASSERT(StateIsPhysical(state1, dim), "");
MFEM_ASSERT(StateIsPhysical(state2, dim), "");
const double maxE1 = ComputeMaxCharSpeed(state1, dim);
const double maxE2 = ComputeMaxCharSpeed(state2, dim);
const double maxE = max(maxE1, maxE2);
ComputeFluxDotN(state1, nor, flux1);
ComputeFluxDotN(state2, nor, flux2);
double normag = 0;
for (int i = 0; i < dim; i++)
{
normag += nor(i) * nor(i);
}
normag = sqrt(normag);
for (int i = 0; i < num_equation; i++)
{
flux(i) = 0.5 * (flux1(i) + flux2(i))
- 0.5 * maxE * (state2(i) - state1(i)) * normag;
}
return maxE;
}
// Implementation of class FaceIntegrator
FaceIntegrator::FaceIntegrator(RiemannSolver &rsolver_, const int dim) :
rsolver(rsolver_),
funval1(num_equation),
funval2(num_equation),
nor(dim),
fluxN(num_equation) { }
void FaceIntegrator::AssembleFaceVector(const FiniteElement &el1,
const FiniteElement &el2,
FaceElementTransformations &Tr,
const Vector &elfun, Vector &elvect)
{
// Compute the term <F.n(u),[w]> on the interior faces.
const int dof1 = el1.GetDof();
const int dof2 = el2.GetDof();
shape1.SetSize(dof1);
shape2.SetSize(dof2);
elvect.SetSize((dof1 + dof2) * num_equation);
elvect = 0.0;
DenseMatrix elfun1_mat(elfun.GetData(), dof1, num_equation);
DenseMatrix elfun2_mat(elfun.GetData() + dof1 * num_equation, dof2,
num_equation);
DenseMatrix elvect1_mat(elvect.GetData(), dof1, num_equation);
DenseMatrix elvect2_mat(elvect.GetData() + dof1 * num_equation, dof2,
num_equation);
// Integration order calculation from DGTraceIntegrator
int intorder;
if (Tr.Elem2No >= 0)
intorder = (min(Tr.Elem1->OrderW(), Tr.Elem2->OrderW()) +
2*max(el1.GetOrder(), el2.GetOrder()));
else
{
intorder = Tr.Elem1->OrderW() + 2*el1.GetOrder();
}
if (el1.Space() == FunctionSpace::Pk)
{
intorder++;
}
const IntegrationRule *ir = &IntRules.Get(Tr.GetGeometryType(), intorder);
for (int i = 0; i < ir->GetNPoints(); i++)
{
const IntegrationPoint &ip = ir->IntPoint(i);
Tr.SetAllIntPoints(&ip); // set face and element int. points
// Calculate basis functions on both elements at the face
el1.CalcShape(Tr.GetElement1IntPoint(), shape1);
el2.CalcShape(Tr.GetElement2IntPoint(), shape2);
// Interpolate elfun at the point
elfun1_mat.MultTranspose(shape1, funval1);
elfun2_mat.MultTranspose(shape2, funval2);
// Get the normal vector and the flux on the face
CalcOrtho(Tr.Jacobian(), nor);
const double mcs = rsolver.Eval(funval1, funval2, nor, fluxN);
// Update max char speed
if (mcs > max_char_speed) { max_char_speed = mcs; }
fluxN *= ip.weight;
for (int k = 0; k < num_equation; k++)
{
for (int s = 0; s < dof1; s++)
{
elvect1_mat(s, k) -= fluxN(k) * shape1(s);
}
for (int s = 0; s < dof2; s++)
{
elvect2_mat(s, k) += fluxN(k) * shape2(s);
}
}
}
}
// Check that the state is physical - enabled in debug mode
bool StateIsPhysical(const Vector &state, const int dim)
{
const double den = state(0);
const Vector den_vel(state.GetData() + 1, dim);
const double den_energy = state(1 + dim);
if (den < 0)
{
cout << "Negative density: ";
for (int i = 0; i < state.Size(); i++)
{
cout << state(i) << " ";
}
cout << endl;
return false;
}
if (den_energy <= 0)
{
cout << "Negative energy: ";
for (int i = 0; i < state.Size(); i++)
{
cout << state(i) << " ";
}
cout << endl;
return false;
}
double den_vel2 = 0;
for (int i = 0; i < dim; i++) { den_vel2 += den_vel(i) * den_vel(i); }
den_vel2 /= den;
const double pres = (specific_heat_ratio - 1.0) * (den_energy - 0.5 * den_vel2);
if (pres <= 0)
{
cout << "Negative pressure: " << pres << ", state: ";
for (int i = 0; i < state.Size(); i++)
{
cout << state(i) << " ";
}
cout << endl;
return false;
}
return true;
}
// Initial condition
void InitialCondition(const Vector &x, Vector &y)
{
MFEM_ASSERT(x.Size() == 2, "");
double radius = 0, Minf = 0, beta = 0;
if (problem == 1)
{
// "Fast vortex"
radius = 0.2;
Minf = 0.5;
beta = 1. / 5.;
}
else if (problem == 2)
{
// "Slow vortex"
radius = 0.2;
Minf = 0.05;
beta = 1. / 50.;
}
else
{
mfem_error("Cannot recognize problem."
"Options are: 1 - fast vortex, 2 - slow vortex");
}
const double xc = 0.0, yc = 0.0;
// Nice units
const double vel_inf = 1.;
const double den_inf = 1.;
// Derive remainder of background state from this and Minf
const double pres_inf = (den_inf / specific_heat_ratio) * (vel_inf / Minf) *
(vel_inf / Minf);
const double temp_inf = pres_inf / (den_inf * gas_constant);
double r2rad = 0.0;
r2rad += (x(0) - xc) * (x(0) - xc);
r2rad += (x(1) - yc) * (x(1) - yc);
r2rad /= (radius * radius);
const double shrinv1 = 1.0 / (specific_heat_ratio - 1.);
const double velX = vel_inf * (1 - beta * (x(1) - yc) / radius * exp(
-0.5 * r2rad));
const double velY = vel_inf * beta * (x(0) - xc) / radius * exp(-0.5 * r2rad);
const double vel2 = velX * velX + velY * velY;
const double specific_heat = gas_constant * specific_heat_ratio * shrinv1;
const double temp = temp_inf - 0.5 * (vel_inf * beta) *
(vel_inf * beta) / specific_heat * exp(-r2rad);
const double den = den_inf * pow(temp/temp_inf, shrinv1);
const double pres = den * gas_constant * temp;
const double energy = shrinv1 * pres / den + 0.5 * vel2;
y(0) = den;
y(1) = den * velX;
y(2) = den * velY;
y(3) = den * energy;
}