-
Notifications
You must be signed in to change notification settings - Fork 0
/
ex33p.cpp
446 lines (407 loc) · 14.9 KB
/
ex33p.cpp
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
// MFEM Example 33 - Parallel Version
//
// Compile with: make ex33p
//
// Sample runs: mpirun -np 4 ex33p -m ../data/square-disc.mesh -alpha 0.33 -o 2
// mpirun -np 4 ex33p -m ../data/square-disc.mesh -alpha 4.5 -o 3
// mpirun -np 4 ex33p -m ../data/star.mesh -alpha 1.4 -o 3
// mpirun -np 4 ex33p -m ../data/star.mesh -alpha 0.99 -o 3
// mpirun -np 4 ex33p -m ../data/inline-quad.mesh -alpha 0.5 -o 3
// mpirun -np 4 ex33p -m ../data/amr-quad.mesh -alpha 1.5 -o 3
// mpirun -np 4 ex33p -m ../data/disc-nurbs.mesh -alpha 0.33 -o 3 -r 2
// mpirun -np 4 ex33p -m ../data/disc-nurbs.mesh -alpha 2.4 -o 3 -r 4
// mpirun -np 4 ex33p -m ../data/l-shape.mesh -alpha 0.33 -o 3 -r 4
// mpirun -np 4 ex33p -m ../data/l-shape.mesh -alpha 1.7 -o 3 -r 5
//
// Verification runs:
// mpirun -np 4 ex33p -m ../data/inline-segment.mesh -ver -alpha 1.7 -o 2 -r 2
// mpirun -np 4 ex33p -m ../data/inline-quad.mesh -ver -alpha 1.2 -o 2 -r 2
// mpirun -np 4 ex33p -m ../data/amr-quad.mesh -ver -alpha 2.6 -o 2 -r 2
// mpirun -np 4 ex33p -m ../data/inline-hex.mesh -ver -alpha 0.3 -o 2 -r 1
// Note: the analytic solution to this problem is u = ∏_{i=0}^{dim-1} sin(π x_i)
// for all alpha.
//
// Description:
//
// In this example we solve the following fractional PDE with MFEM:
//
// ( - Δ )^α u = f in Ω, u = 0 on ∂Ω, 0 < α,
//
// To solve this FPDE, we apply the operator ( - Δ )^(-N), where the integer
// N is given by floor(α). By doing so, we obtain
//
// ( - Δ )^(α-N) u = ( - Δ )^(-N) f in Ω, u = 0 on ∂Ω, 0 < α.
//
// We first compute the right hand side by solving the integer order PDE
//
// ( - Δ )^N g = f in Ω, g = ( - Δ )^k g = 0 on ∂Ω, k = 1,..,N-1
//
// The remaining FPDE is then given by
//
// ( - Δ )^(α-N) u = g in Ω, u = 0 on ∂Ω.
//
// We rely on a rational approximation [2] of the normal linear operator
// A^{-α + N}, where A = - Δ (with associated homogeneous boundary conditions)
// and (a-N) in (0,1). We approximate the operator
//
// A^{-α+N} ≈ Σ_{i=0}^M c_i (A + d_i I)^{-1}, d_0 = 0, d_i > 0,
//
// where I is the L2-identity operator and the coefficients c_i and d_i
// are generated offline to a prescribed accuracy in a pre-processing step.
// We use the triple-A algorithm [1] to generate the rational approximation
// that this partial fractional expansion derives from. We then solve M+1
// independent integer-order PDEs,
//
// A u_i + d_i u_i = c_i g in Ω, u_i = 0 on ∂Ω, i=0,...,M,
//
// using MFEM and sum u_i to arrive at an approximate solution of the FPDE
//
// u ≈ Σ_{i=0}^M u_i.
//
// (If alpha is an integer, we stop after the first PDE was solved.)
//
// References:
//
// [1] Nakatsukasa, Y., Sète, O., & Trefethen, L. N. (2018). The AAA algorithm
// for rational approximation. SIAM Journal on Scientific Computing, 40(3),
// A1494-A1522.
//
// [2] Harizanov, S., Lazarov, R., Margenov, S., Marinov, P., & Pasciak, J.
// (2020). Analysis of numerical methods for spectral fractional elliptic
// equations based on the best uniform rational approximation. Journal of
// Computational Physics, 408, 109285.
//
#include "mfem.hpp"
#include <fstream>
#include <iostream>
#include <math.h>
#include <string>
#include "ex33.hpp"
using namespace std;
using namespace mfem;
int main(int argc, char *argv[])
{
#ifdef MFEM_USE_SINGLE
cout << "This example is not supported in single precision.\n\n";
return MFEM_SKIP_RETURN_VALUE;
#endif
// 0. Initialize MPI.
Mpi::Init(argc, argv);
int num_procs = Mpi::WorldSize();
int myid = Mpi::WorldRank();
Hypre::Init();
// 1. Parse command-line options.
const char *mesh_file = "../data/star.mesh";
int order = 1;
int num_refs = 3;
real_t alpha = 0.5;
bool visualization = true;
bool verification = false;
OptionsParser args(argc, argv);
args.AddOption(&mesh_file, "-m", "--mesh",
"Mesh file to use.");
args.AddOption(&order, "-o", "--order",
"Finite element order (polynomial degree) or -1 for"
" isoparametric space.");
args.AddOption(&num_refs, "-r", "--refs",
"Number of uniform refinements");
args.AddOption(&alpha, "-alpha", "--alpha",
"Fractional exponent");
args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
"--no-visualization",
"Enable or disable GLVis visualization.");
args.AddOption(&verification, "-ver", "--verification", "-no-ver",
"--no-verification",
"Use sinusoidal function (f) for analytic comparison.");
args.Parse();
if (!args.Good())
{
args.PrintUsage(cout);
return 1;
}
if (Mpi::Root())
{
args.PrintOptions(cout);
}
Array<real_t> coeffs, poles;
int progress_steps = 1;
// 2. Compute the rational expansion coefficients that define the
// integer-order PDEs.
const int power_of_laplace = floor(alpha);
real_t exponent_to_approximate = alpha - power_of_laplace;
bool integer_order = false;
// Check if alpha is an integer or not.
if (abs(exponent_to_approximate) > 1e-12)
{
if (Mpi::Root())
{
mfem::out << "Approximating the fractional exponent "
<< exponent_to_approximate
<< endl;
}
ComputePartialFractionApproximation(exponent_to_approximate, coeffs,
poles);
// If the example is build without LAPACK, the exponent_to_approximate
// might be modified by the function call above.
alpha = exponent_to_approximate + power_of_laplace;
}
else
{
integer_order = true;
if (Mpi::Root())
{
mfem::out << "Treating integer order PDE." << endl;
}
}
// 3. Read the mesh from the given mesh file.
Mesh mesh(mesh_file, 1, 1);
int dim = mesh.Dimension();
// 4. Refine the mesh to increase the resolution.
for (int i = 0; i < num_refs; i++)
{
mesh.UniformRefinement();
}
ParMesh pmesh(MPI_COMM_WORLD, mesh);
mesh.Clear();
// 5. Define a finite element space on the mesh.
H1_FECollection fec(order, dim);
ParFiniteElementSpace fespace(&pmesh, &fec);
if (Mpi::Root())
{
cout << "Number of finite element unknowns: "
<< fespace.GetTrueVSize() << endl;
}
// 6. Determine the list of true (i.e. conforming) essential boundary dofs.
Array<int> ess_tdof_list;
if (pmesh.bdr_attributes.Size())
{
Array<int> ess_bdr(pmesh.bdr_attributes.Max());
ess_bdr = 1;
fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
}
// 7. Define diffusion coefficient, load, and solution GridFunction.
auto func = [&alpha](const Vector &x)
{
real_t val = 1.0;
for (int i=0; i<x.Size(); i++)
{
val *= sin(M_PI*x(i));
}
return pow(x.Size()*pow(M_PI,2), alpha) * val;
};
FunctionCoefficient f(func);
ConstantCoefficient one(1.0);
ParGridFunction u(&fespace);
ParGridFunction x(&fespace);
ParGridFunction g(&fespace);
u = 0.0;
x = 0.0;
g = 0.0;
// 8. Prepare for visualization.
char vishost[] = "localhost";
int visport = 19916;
// 9. Set up the linear form b(.) for integer-order PDE solves.
ParLinearForm b(&fespace);
if (verification)
{
// This statement is only relevant for the verification of the code. It
// uses a different f such that an analytic solution is known and easy
// to compare with the numerical one. The FPDE becomes:
// (-Δ)^α u = (2\pi ^2)^α sin(\pi x) sin(\pi y) on [0,1]^2
// -> u(x,y) = sin(\pi x) sin(\pi y)
b.AddDomainIntegrator(new DomainLFIntegrator(f));
}
else
{
b.AddDomainIntegrator(new DomainLFIntegrator(one));
}
b.Assemble();
// ------------------------------------------------------------------------
// 10. Solve the PDE (-Δ)^N g = f, i.e. compute g = (-Δ)^{-1}^N f.
// ------------------------------------------------------------------------
if (power_of_laplace > 0)
{
// 10.1 Compute Stiffnes Matrix
ParBilinearForm k(&fespace);
k.AddDomainIntegrator(new DiffusionIntegrator(one));
k.Assemble();
// 10.2 Compute Mass Matrix
ParBilinearForm m(&fespace);
m.AddDomainIntegrator(new MassIntegrator(one));
m.Assemble();
HypreParMatrix mass;
Array<int> empty;
m.FormSystemMatrix(empty, mass);
// 10.3 Form the system of equations
Vector B, X;
OperatorPtr Op;
k.FormLinearSystem(ess_tdof_list, g, b, Op, X, B);
HypreBoomerAMG prec;
prec.SetPrintLevel(-1);
CGSolver cg(MPI_COMM_WORLD);
cg.SetRelTol(1e-12);
cg.SetMaxIter(2000);
cg.SetPrintLevel(3);
cg.SetPreconditioner(prec);
cg.SetOperator(*Op);
if (Mpi::Root())
{
mfem::out << "\nComputing (-Δ) ^ -" << power_of_laplace
<< " ( f ) " << endl;
}
for (int i = 0; i < power_of_laplace; i++)
{
// 10.4 Solve the linear system Op X = B (N times).
cg.Mult(B, X);
// 10.5 Visualize the solution g of -Δ ^ N g = f in the last step
if (i == power_of_laplace - 1)
{
// Needed for visualization and solution verification.
k.RecoverFEMSolution(X, b, g);
if (integer_order && verification)
{
// For an integer order PDE, g is also our solution u.
u+=g;
}
if (visualization)
{
socketstream fout;
ostringstream oss_f;
fout.open(vishost, visport);
fout.precision(8);
oss_f.str(""); oss_f.clear();
oss_f << "Step " << progress_steps++ << ": Solution of PDE -Δ ^ "
<< power_of_laplace
<< " g = f";
fout << "parallel " << num_procs << " " << myid << "\n"
<< "solution\n" << pmesh << g
<< "window_title '" << oss_f.str() << "'" << flush;
}
}
// 10.6 Prepare for next iteration (primal / dual space)
mass.Mult(X, B);
X.SetSubVectorComplement(ess_tdof_list,0.0);
}
// 10.7 Extract solution for the next step. The b now corresponds to the
// function g in the PDE.
const SparseMatrix* rm = fespace.GetRestrictionMatrix();
rm->MultTranspose(B, b);
}
// ------------------------------------------------------------------------
// 11. Solve the fractional PDE by solving M integer order PDEs and adding
// up the solutions.
// ------------------------------------------------------------------------
if (!integer_order)
{
// Setup visualization.
socketstream xout, uout;
ostringstream oss_x, oss_u;
if (visualization)
{
xout.open(vishost, visport);
xout.precision(8);
uout.open(vishost, visport);
uout.precision(8);
}
// Iterate over all expansion coefficient that contribute to the
// solution.
for (int i = 0; i < coeffs.Size(); i++)
{
if (Mpi::Root())
{
mfem::out << "\nSolving PDE -Δ u + " << -poles[i]
<< " u = " << coeffs[i] << " g " << endl;
}
// 11.1 Reset GridFunction for integer-order PDE solve.
x = 0.0;
// 11.2 Set up the bilinear form a(.,.) for integer-order PDE solve.
ParBilinearForm a(&fespace);
a.AddDomainIntegrator(new DiffusionIntegrator(one));
ConstantCoefficient d_i(-poles[i]);
a.AddDomainIntegrator(new MassIntegrator(d_i));
a.Assemble();
// 11.3 Assemble the bilinear form and the corresponding linear system.
OperatorPtr A;
Vector B, X;
a.FormLinearSystem(ess_tdof_list, x, b, A, X, B);
// 11.4 Solve the linear system A X = B.
HypreBoomerAMG prec;
prec.SetPrintLevel(-1);
CGSolver cg(MPI_COMM_WORLD);
cg.SetRelTol(1e-12);
cg.SetMaxIter(2000);
cg.SetPrintLevel(3);
cg.SetPreconditioner(prec);
cg.SetOperator(*A);
cg.Mult(B, X);
// 11.5 Recover the solution as a finite element grid function.
a.RecoverFEMSolution(X, b, x);
// 11.6 Accumulate integer-order PDE solutions.
x *= coeffs[i];
u += x;
// 11.7 Send fractional PDE solution to a GLVis server.
if (visualization)
{
oss_x.str(""); oss_x.clear();
oss_x << "Step " << progress_steps
<< ": Solution of PDE -Δ u + " << -poles[i]
<< " u = " << coeffs[i] << " g";
xout << "parallel " << num_procs << " " << myid << "\n"
<< "solution\n" << pmesh << x
<< "window_title '" << oss_x.str() << "'" << flush;
oss_u.str(""); oss_u.clear();
oss_u << "Step " << progress_steps + 1
<< ": Solution of fractional PDE (-Δ)^" << alpha
<< " u = f";
uout << "parallel " << num_procs << " " << myid << "\n"
<< "solution\n" << pmesh << u
<< "window_title '" << oss_u.str() << "'"
<< flush;
}
}
}
// ------------------------------------------------------------------------
// 12. (optional) Verify the solution.
// ------------------------------------------------------------------------
if (verification)
{
auto solution = [] (const Vector &x)
{
real_t val = 1.0;
for (int i=0; i<x.Size(); i++)
{
val *= sin(M_PI*x(i));
}
return val;
};
FunctionCoefficient sol(solution);
real_t l2_error = u.ComputeL2Error(sol);
if (Mpi::Root())
{
string analytic_solution,expected_mesh;
switch (dim)
{
case 1:
analytic_solution = "sin(π x)";
expected_mesh = "inline_segment.mesh";
break;
case 2:
analytic_solution = "sin(π x) sin(π y)";
expected_mesh = "inline_quad.mesh";
break;
default:
analytic_solution = "sin(π x) sin(π y) sin(π z)";
expected_mesh = "inline_hex.mesh";
break;
}
mfem::out << "\n" << string(80,'=')
<< "\n\nSolution Verification in "<< dim << "D \n\n"
<< "Analytic solution : " << analytic_solution << "\n"
<< "Expected mesh : " << expected_mesh <<"\n"
<< "Your mesh : " << mesh_file << "\n"
<< "L2 error : " << l2_error << "\n\n"
<< string(80,'=') << endl;
}
}
return 0;
}