forked from mfem/mfem
-
Notifications
You must be signed in to change notification settings - Fork 0
/
ex6.cpp
279 lines (248 loc) · 10.4 KB
/
ex6.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
// MFEM Example 6
//
// Compile with: make ex6
//
// Sample runs: ex6 -m ../data/square-disc.mesh -o 1
// ex6 -m ../data/square-disc.mesh -o 2
// ex6 -m ../data/square-disc-nurbs.mesh -o 2
// ex6 -m ../data/star.mesh -o 3
// ex6 -m ../data/escher.mesh -o 2
// ex6 -m ../data/fichera.mesh -o 2
// ex6 -m ../data/disc-nurbs.mesh -o 2
// ex6 -m ../data/ball-nurbs.mesh
// ex6 -m ../data/pipe-nurbs.mesh
// ex6 -m ../data/star-surf.mesh -o 2
// ex6 -m ../data/square-disc-surf.mesh -o 2
// ex6 -m ../data/amr-quad.mesh
// ex6 -m ../data/inline-segment.mesh -o 1 -md 100
//
// Device sample runs:
// ex6 -pa -d cuda
// ex6 -pa -d occa-cuda
// ex6 -pa -d raja-omp
// ex6 -pa -d ceed-cpu
// * ex6 -pa -d ceed-cuda
// ex6 -pa -d ceed-cuda:/gpu/cuda/shared
//
// Description: This is a version of Example 1 with a simple adaptive mesh
// refinement loop. The problem being solved is again the Laplace
// equation -Delta u = 1 with homogeneous Dirichlet boundary
// conditions. The problem is solved on a sequence of meshes which
// are locally refined in a conforming (triangles, tetrahedrons)
// or non-conforming (quadrilaterals, hexahedra) manner according
// to a simple ZZ error estimator.
//
// The example demonstrates MFEM's capability to work with both
// conforming and nonconforming refinements, in 2D and 3D, on
// linear, curved and surface meshes. Interpolation of functions
// from coarse to fine meshes, as well as persistent GLVis
// visualization are also illustrated.
//
// We recommend viewing Example 1 before viewing this example.
#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 pa = false;
const char *device_config = "cpu";
int max_dofs = 50000;
bool LSZZ = false;
bool visualization = true;
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(&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(&max_dofs, "-md", "--max-dofs",
"Stop after reaching this many degrees of freedom.");
args.AddOption(&LSZZ, "-ls", "--ls-zz", "-no-ls",
"--no-ls-zz",
"Switch to least-squares ZZ estimator.");
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();
int sdim = mesh.SpaceDimension();
// 4. Since a NURBS mesh can currently only be refined uniformly, we need to
// convert it to a piecewise-polynomial curved mesh. First we refine the
// NURBS mesh a bit more and then project the curvature to quadratic Nodes.
if (mesh.NURBSext)
{
for (int i = 0; i < 2; i++)
{
mesh.UniformRefinement();
}
mesh.SetCurvature(2);
}
// 5. Define a finite element space on the mesh. The polynomial order is
// one (linear) by default, but this can be changed on the command line.
H1_FECollection fec(order, dim);
FiniteElementSpace fespace(&mesh, &fec);
// 6. As in Example 1, we set up bilinear and linear forms corresponding to
// the Laplace problem -\Delta u = 1. We don't assemble the discrete
// problem yet, this will be done in the main loop.
BilinearForm a(&fespace);
if (pa)
{
a.SetAssemblyLevel(AssemblyLevel::PARTIAL);
a.SetDiagonalPolicy(Operator::DIAG_ONE);
}
LinearForm b(&fespace);
ConstantCoefficient one(1.0);
ConstantCoefficient zero(0.0);
BilinearFormIntegrator *integ = new DiffusionIntegrator(one);
a.AddDomainIntegrator(integ);
b.AddDomainIntegrator(new DomainLFIntegrator(one));
// 7. The solution vector x and the associated finite element grid function
// will be maintained over the AMR iterations. We initialize it to zero.
GridFunction x(&fespace);
x = 0.0;
// 8. All boundary attributes will be used for essential (Dirichlet) BC.
MFEM_VERIFY(mesh.bdr_attributes.Size() > 0,
"Boundary attributes required in the mesh.");
Array<int> ess_bdr(mesh.bdr_attributes.Max());
ess_bdr = 1;
// 9. Connect to GLVis.
char vishost[] = "localhost";
int visport = 19916;
socketstream sol_sock;
if (visualization)
{
sol_sock.open(vishost, visport);
}
// 10. Set up an error estimator. Here we use the Zienkiewicz-Zhu estimator
// that uses the ComputeElementFlux method of the DiffusionIntegrator to
// recover a smoothed flux (gradient) that is subtracted from the element
// flux to get an error indicator. We need to supply the space for the
// smoothed flux: an (H1)^sdim (i.e., vector-valued) space is used here.
ErrorEstimator *estimator{nullptr};
if (LSZZ)
{
estimator = new LSZienkiewiczZhuEstimator(*integ, x);
if (dim == 3 && mesh.GetElementType(0) != Element::HEXAHEDRON)
{
dynamic_cast<LSZienkiewiczZhuEstimator *>
(estimator)->SetTichonovRegularization();
}
}
else
{
auto flux_fes = new FiniteElementSpace(&mesh, &fec, sdim);
estimator = new ZienkiewiczZhuEstimator(*integ, x, flux_fes);
dynamic_cast<ZienkiewiczZhuEstimator *>(estimator)->SetAnisotropic();
}
// 11. A refiner selects and refines elements based on a refinement strategy.
// The strategy here is to refine elements with errors larger than a
// fraction of the maximum element error. Other strategies are possible.
// The refiner will call the given error estimator.
ThresholdRefiner refiner(*estimator);
refiner.SetTotalErrorFraction(0.7);
// 12. The main AMR loop. In each iteration we solve the problem on the
// current mesh, visualize the solution, and refine the mesh.
for (int it = 0; ; it++)
{
int cdofs = fespace.GetTrueVSize();
cout << "\nAMR iteration " << it << endl;
cout << "Number of unknowns: " << cdofs << endl;
// 13. Assemble the right-hand side.
b.Assemble();
// 14. Set Dirichlet boundary values in the GridFunction x.
// Determine the list of Dirichlet true DOFs in the linear system.
Array<int> ess_tdof_list;
x.ProjectBdrCoefficient(zero, ess_bdr);
fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
// 15. Assemble the stiffness matrix.
a.Assemble();
// 16. Create the linear system: eliminate boundary conditions, constrain
// hanging nodes and possibly apply other transformations. The system
// will be solved for true (unconstrained) DOFs only.
OperatorPtr A;
Vector B, X;
const int copy_interior = 1;
a.FormLinearSystem(ess_tdof_list, x, b, A, X, B, copy_interior);
// 17. 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, 3, 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 // Diagonal preconditioning in partial assembly mode.
{
OperatorJacobiSmoother M(a, ess_tdof_list);
PCG(*A, M, B, X, 3, 2000, 1e-12, 0.0);
}
// 18. After solving the linear system, reconstruct the solution as a
// finite element GridFunction. Constrained nodes are interpolated
// from true DOFs (it may therefore happen that x.Size() >= X.Size()).
a.RecoverFEMSolution(X, b, x);
// 19. Send solution by socket to the GLVis server.
if (visualization && sol_sock.good())
{
sol_sock.precision(8);
sol_sock << "solution\n" << mesh << x << flush;
}
if (cdofs > max_dofs)
{
cout << "Reached the maximum number of dofs. Stop." << endl;
break;
}
// 20. Call the refiner to modify the mesh. The refiner calls the error
// estimator to obtain element errors, then it selects elements to be
// refined and finally it modifies the mesh. The Stop() method can be
// used to determine if a stopping criterion was met.
refiner.Apply(mesh);
if (refiner.Stop())
{
cout << "Stopping criterion satisfied. Stop." << endl;
break;
}
// 21. Update the space to reflect the new state of the mesh. Also,
// interpolate the solution x so that it lies in the new space but
// represents the same function. This saves solver iterations later
// since we'll have a good initial guess of x in the next step.
// Internally, FiniteElementSpace::Update() calculates an
// interpolation matrix which is then used by GridFunction::Update().
fespace.Update();
x.Update();
// 22. Inform also the bilinear and linear forms that the space has
// changed.
a.Update();
b.Update();
}
delete estimator;
return 0;
}