forked from mfem/mfem
-
Notifications
You must be signed in to change notification settings - Fork 0
/
ex1.cpp
283 lines (262 loc) · 10.4 KB
/
ex1.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
// MFEM Example 1
//
// Compile with: make ex1
//
// Sample runs: ex1 -m ../data/square-disc.mesh
// ex1 -m ../data/star.mesh
// ex1 -m ../data/star-mixed.mesh
// ex1 -m ../data/escher.mesh
// ex1 -m ../data/fichera.mesh
// ex1 -m ../data/fichera-mixed.mesh
// ex1 -m ../data/toroid-wedge.mesh
// ex1 -m ../data/octahedron.mesh -o 1
// ex1 -m ../data/periodic-annulus-sector.msh
// ex1 -m ../data/periodic-torus-sector.msh
// ex1 -m ../data/square-disc-p2.vtk -o 2
// ex1 -m ../data/square-disc-p3.mesh -o 3
// ex1 -m ../data/square-disc-nurbs.mesh -o -1
// ex1 -m ../data/star-mixed-p2.mesh -o 2
// ex1 -m ../data/disc-nurbs.mesh -o -1
// ex1 -m ../data/pipe-nurbs.mesh -o -1
// ex1 -m ../data/fichera-mixed-p2.mesh -o 2
// ex1 -m ../data/star-surf.mesh
// ex1 -m ../data/square-disc-surf.mesh
// ex1 -m ../data/inline-segment.mesh
// ex1 -m ../data/amr-quad.mesh
// ex1 -m ../data/amr-hex.mesh
// ex1 -m ../data/fichera-amr.mesh
// ex1 -m ../data/mobius-strip.mesh
// ex1 -m ../data/mobius-strip.mesh -o -1 -sc
//
// Device sample runs:
// ex1 -pa -d cuda
// ex1 -fa -d cuda
// ex1 -pa -d raja-cuda
// * ex1 -pa -d raja-hip
// ex1 -pa -d occa-cuda
// ex1 -pa -d raja-omp
// ex1 -pa -d occa-omp
// ex1 -pa -d ceed-cpu
// ex1 -pa -d ceed-cpu -o 4 -a
// ex1 -pa -d ceed-cpu -m ../data/square-mixed.mesh
// ex1 -pa -d ceed-cpu -m ../data/fichera-mixed.mesh
// * ex1 -pa -d ceed-cuda
// * ex1 -pa -d ceed-hip
// ex1 -pa -d ceed-cuda:/gpu/cuda/shared
// ex1 -pa -d ceed-cuda:/gpu/cuda/shared -m ../data/square-mixed.mesh
// ex1 -pa -d ceed-cuda:/gpu/cuda/shared -m ../data/fichera-mixed.mesh
// ex1 -m ../data/beam-hex.mesh -pa -d cuda
// ex1 -m ../data/beam-tet.mesh -pa -d ceed-cpu
// ex1 -m ../data/beam-tet.mesh -pa -d ceed-cuda:/gpu/cuda/ref
//
// Description: This example code demonstrates the use of MFEM to define a
// simple finite element discretization of the Laplace problem
// -Delta u = 1 with homogeneous Dirichlet boundary conditions.
// Specifically, we discretize using a FE space of the specified
// order, or if order < 1 using an isoparametric/isogeometric
// space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
// NURBS mesh, etc.)
//
// The example highlights the use of mesh refinement, finite
// element grid functions, as well as linear and bilinear forms
// corresponding to the left-hand side and right-hand side of the
// discrete linear system. We also cover the explicit elimination
// of essential boundary conditions, static condensation, and the
// optional connection to the GLVis tool for visualization.
#include "mfem.hpp"
#include <fstream>
#include <iostream>
using namespace std;
using namespace mfem;
int main(int argc, char *argv[])
{
// 1. Parse command-line options.
const char *mesh_file = "../data/star.mesh";
int order = 1;
bool static_cond = false;
bool pa = false;
bool fa = false;
const char *device_config = "cpu";
bool visualization = true;
bool algebraic_ceed = 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(&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(&fa, "-fa", "--full-assembly", "-no-fa",
"--no-full-assembly", "Enable Full Assembly.");
args.AddOption(&device_config, "-d", "--device",
"Device configuration string, see Device::Configure().");
#ifdef MFEM_USE_CEED
args.AddOption(&algebraic_ceed, "-a", "--algebraic", "-no-a", "--no-algebraic",
"Use algebraic Ceed solver");
#endif
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);
// 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(mesh_file, 1, 1);
int dim = mesh.Dimension();
// 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 continuous
// Lagrange finite elements of the specified order. If order < 1, we
// instead use an isoparametric/isogeometric space.
FiniteElementCollection *fec;
bool delete_fec;
if (order > 0)
{
fec = new H1_FECollection(order, dim);
delete_fec = true;
}
else if (mesh.GetNodes())
{
fec = mesh.GetNodes()->OwnFEC();
delete_fec = false;
cout << "Using isoparametric FEs: " << fec->Name() << endl;
}
else
{
fec = new H1_FECollection(order = 1, dim);
delete_fec = true;
}
FiniteElementSpace fespace(&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 (1,phi_i) where phi_i are
// the basis functions in the finite element fespace.
LinearForm b(&fespace);
ConstantCoefficient one(1.0);
b.AddDomainIntegrator(new DomainLFIntegrator(one));
b.Assemble();
// 8. Define the solution vector x as a finite element grid function
// corresponding to fespace. Initialize x with initial guess of zero,
// which satisfies the boundary conditions.
GridFunction x(&fespace);
x = 0.0;
// 9. Set up the bilinear form a(.,.) on the finite element space
// corresponding to the Laplacian operator -Delta, by adding the Diffusion
// domain integrator.
BilinearForm a(&fespace);
if (pa) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
if (fa)
{
a.SetAssemblyLevel(AssemblyLevel::FULL);
// Sort the matrix column indices when running on GPU or with OpenMP (i.e.
// when Device::IsEnabled() returns true). This makes the results
// bit-for-bit deterministic at the cost of somewhat longer run time.
a.EnableSparseMatrixSorting(Device::IsEnabled());
}
a.AddDomainIntegrator(new DiffusionIntegrator(one));
// 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)
{
#ifndef MFEM_USE_SUITESPARSE
// Use a simple symmetric Gauss-Seidel preconditioner with PCG.
GSSmoother M((SparseMatrix&)(*A));
PCG(*A, M, B, X, 1, 200, 1e-12, 0.0);
#else
// 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
}
else
{
if (UsesTensorBasis(fespace))
{
if (algebraic_ceed)
{
ceed::AlgebraicSolver M(a, ess_tdof_list);
PCG(*A, M, B, X, 1, 400, 1e-12, 0.0);
}
else
{
OperatorJacobiSmoother M(a, ess_tdof_list);
PCG(*A, M, B, X, 1, 400, 1e-12, 0.0);
}
}
else
{
CG(*A, B, X, 1, 400, 1e-12, 0.0);
}
}
// 12. Recover the solution as a finite element grid function.
a.RecoverFEMSolution(X, b, x);
// 13. 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);
// 14. 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;
}
// 15. Free the used memory.
if (delete_fec)
{
delete fec;
}
return 0;
}