forked from mfem/mfem
-
Notifications
You must be signed in to change notification settings - Fork 0
/
ex3.cpp
272 lines (246 loc) · 9.66 KB
/
ex3.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
// MFEM Example 3
//
// Compile with: make ex3
//
// Sample runs: ex3 -m ../data/star.mesh
// ex3 -m ../data/beam-tri.mesh -o 2
// ex3 -m ../data/beam-tet.mesh
// ex3 -m ../data/beam-hex.mesh
// ex3 -m ../data/beam-hex.mesh -o 2 -pa
// ex3 -m ../data/escher.mesh
// ex3 -m ../data/escher.mesh -o 2
// ex3 -m ../data/fichera.mesh
// ex3 -m ../data/fichera-q2.vtk
// ex3 -m ../data/fichera-q3.mesh
// ex3 -m ../data/square-disc-nurbs.mesh
// ex3 -m ../data/beam-hex-nurbs.mesh
// ex3 -m ../data/amr-hex.mesh
// ex3 -m ../data/fichera-amr.mesh
// ex3 -m ../data/ref-prism.mesh -o 1
// ex3 -m ../data/octahedron.mesh -o 1
// ex3 -m ../data/star-surf.mesh -o 1
// ex3 -m ../data/mobius-strip.mesh -f 0.1
// ex3 -m ../data/klein-bottle.mesh -f 0.1
//
// Device sample runs:
// ex3 -m ../data/star.mesh -pa -d cuda
// ex3 -m ../data/star.mesh -pa -d raja-cuda
// ex3 -m ../data/star.mesh -pa -d raja-omp
// ex3 -m ../data/beam-hex.mesh -pa -d cuda
//
// Description: This example code solves a simple electromagnetic diffusion
// problem corresponding to the second order definite Maxwell
// equation curl curl E + E = f with boundary condition
// E x n = <given tangential field>. Here, we use a given exact
// solution E and compute the corresponding r.h.s. f.
// We discretize with Nedelec finite elements in 2D or 3D.
//
// The example demonstrates the use of H(curl) finite element
// spaces with the curl-curl and the (vector finite element) mass
// bilinear form, as well as the computation of discretization
// error when the exact solution is known. Static condensation is
// also illustrated.
//
// We recommend viewing examples 1-2 before viewing this example.
#include "mfem.hpp"
#include <fstream>
#include <iostream>
using namespace std;
using namespace mfem;
// Exact solution, E, and r.h.s., f. See below for implementation.
void E_exact(const Vector &, Vector &);
void f_exact(const Vector &, Vector &);
double freq = 1.0, kappa;
int dim;
int main(int argc, char *argv[])
{
// 1. Parse command-line options.
const char *mesh_file = "../data/beam-tet.mesh";
int order = 1;
bool static_cond = false;
bool pa = false;
const char *device_config = "cpu";
bool visualization = 1;
OptionsParser args(argc, argv);
args.AddOption(&mesh_file, "-m", "--mesh",
"Mesh file to use.");
args.AddOption(&order, "-o", "--order",
"Finite element order (polynomial degree).");
args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
" solution.");
args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
"--no-static-condensation", "Enable static condensation.");
args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
"--no-partial-assembly", "Enable Partial Assembly.");
args.AddOption(&device_config, "-d", "--device",
"Device configuration string, see Device::Configure().");
args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
"--no-visualization",
"Enable or disable GLVis visualization.");
args.Parse();
if (!args.Good())
{
args.PrintUsage(cout);
return 1;
}
args.PrintOptions(cout);
kappa = freq * M_PI;
// 2. Enable hardware devices such as GPUs, and programming models such as
// CUDA, OCCA, RAJA and OpenMP based on command line options.
Device device(device_config);
device.Print();
// 3. Read the mesh from the given mesh file. We can handle triangular,
// quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
// the same code.
Mesh *mesh = new Mesh(mesh_file, 1, 1);
dim = mesh->Dimension();
int sdim = mesh->SpaceDimension();
// 4. Refine the mesh to increase the resolution. In this example we do
// 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
// largest number that gives a final mesh with no more than 50,000
// elements.
{
int ref_levels =
(int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
for (int l = 0; l < ref_levels; l++)
{
mesh->UniformRefinement();
}
}
// 5. Define a finite element space on the mesh. Here we use the Nedelec
// finite elements of the specified order.
FiniteElementCollection *fec = new ND_FECollection(order, dim);
FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
cout << "Number of finite element unknowns: "
<< fespace->GetTrueVSize() << endl;
// 6. Determine the list of true (i.e. conforming) essential boundary dofs.
// In this example, the boundary conditions are defined by marking all
// the boundary attributes from the mesh as essential (Dirichlet) and
// converting them to a list of true dofs.
Array<int> ess_tdof_list;
if (mesh->bdr_attributes.Size())
{
Array<int> ess_bdr(mesh->bdr_attributes.Max());
ess_bdr = 1;
fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
}
// 7. Set up the linear form b(.) which corresponds to the right-hand side
// of the FEM linear system, which in this case is (f,phi_i) where f is
// given by the function f_exact and phi_i are the basis functions in the
// finite element fespace.
VectorFunctionCoefficient f(sdim, f_exact);
LinearForm *b = new LinearForm(fespace);
b->AddDomainIntegrator(new VectorFEDomainLFIntegrator(f));
b->Assemble();
// 8. Define the solution vector x as a finite element grid function
// corresponding to fespace. Initialize x by projecting the exact
// solution. Note that only values from the boundary edges will be used
// when eliminating the non-homogeneous boundary condition to modify the
// r.h.s. vector b.
GridFunction x(fespace);
VectorFunctionCoefficient E(sdim, E_exact);
x.ProjectCoefficient(E);
// 9. Set up the bilinear form corresponding to the EM diffusion operator
// curl muinv curl + sigma I, by adding the curl-curl and the mass domain
// integrators.
Coefficient *muinv = new ConstantCoefficient(1.0);
Coefficient *sigma = new ConstantCoefficient(1.0);
BilinearForm *a = new BilinearForm(fespace);
if (pa) { a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
a->AddDomainIntegrator(new CurlCurlIntegrator(*muinv));
a->AddDomainIntegrator(new VectorFEMassIntegrator(*sigma));
// 10. Assemble the bilinear form and the corresponding linear system,
// applying any necessary transformations such as: eliminating boundary
// conditions, applying conforming constraints for non-conforming AMR,
// static condensation, etc.
if (static_cond) { a->EnableStaticCondensation(); }
a->Assemble();
OperatorPtr A;
Vector B, X;
a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
cout << "Size of linear system: " << A->Height() << endl;
// 11. Solve the linear system A X = B.
if (pa) // Jacobi preconditioning in partial assembly mode
{
OperatorJacobiSmoother M(*a, ess_tdof_list);
PCG(*A, M, B, X, 1, 1000, 1e-12, 0.0);
}
else
{
#ifndef MFEM_USE_SUITESPARSE
// 11. Define a simple symmetric Gauss-Seidel preconditioner and use it to
// solve the system Ax=b with PCG.
GSSmoother M((SparseMatrix&)(*A));
PCG(*A, M, B, X, 1, 500, 1e-12, 0.0);
#else
// 11. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the
// system.
UMFPackSolver umf_solver;
umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
umf_solver.SetOperator(*A);
umf_solver.Mult(B, X);
#endif
}
// 12. Recover the solution as a finite element grid function.
a->RecoverFEMSolution(X, *b, x);
// 13. Compute and print the L^2 norm of the error.
cout << "\n|| E_h - E ||_{L^2} = " << x.ComputeL2Error(E) << '\n' << endl;
// 14. Save the refined mesh and the solution. This output can be viewed
// later using GLVis: "glvis -m refined.mesh -g sol.gf".
{
ofstream mesh_ofs("refined.mesh");
mesh_ofs.precision(8);
mesh->Print(mesh_ofs);
ofstream sol_ofs("sol.gf");
sol_ofs.precision(8);
x.Save(sol_ofs);
}
// 15. Send the solution by socket to a GLVis server.
if (visualization)
{
char vishost[] = "localhost";
int visport = 19916;
socketstream sol_sock(vishost, visport);
sol_sock.precision(8);
sol_sock << "solution\n" << *mesh << x << flush;
}
// 16. Free the used memory.
delete a;
delete sigma;
delete muinv;
delete b;
delete fespace;
delete fec;
delete mesh;
return 0;
}
void E_exact(const Vector &x, Vector &E)
{
if (dim == 3)
{
E(0) = sin(kappa * x(1));
E(1) = sin(kappa * x(2));
E(2) = sin(kappa * x(0));
}
else
{
E(0) = sin(kappa * x(1));
E(1) = sin(kappa * x(0));
if (x.Size() == 3) { E(2) = 0.0; }
}
}
void f_exact(const Vector &x, Vector &f)
{
if (dim == 3)
{
f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
f(1) = (1. + kappa * kappa) * sin(kappa * x(2));
f(2) = (1. + kappa * kappa) * sin(kappa * x(0));
}
else
{
f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
f(1) = (1. + kappa * kappa) * sin(kappa * x(0));
if (x.Size() == 3) { f(2) = 0.0; }
}
}