MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex14.cpp
Go to the documentation of this file.
1// MFEM Example 14
2//
3// Compile with: make ex14
4//
5// Sample runs: ex14 -m ../data/inline-quad.mesh -o 0
6// ex14 -m ../data/star.mesh -r 4 -o 2
7// ex14 -m ../data/star-mixed.mesh -r 4 -o 2
8// ex14 -m ../data/star-mixed.mesh -r 2 -o 2 -k 0 -e 1
9// ex14 -m ../data/escher.mesh -s 1
10// ex14 -m ../data/fichera.mesh -s 1 -k 1
11// ex14 -m ../data/fichera-mixed.mesh -s 1 -k 1
12// ex14 -m ../data/square-disc-p2.vtk -r 3 -o 2
13// ex14 -m ../data/square-disc-p3.mesh -r 2 -o 3
14// ex14 -m ../data/square-disc-nurbs.mesh -o 1
15// ex14 -m ../data/disc-nurbs.mesh -r 3 -o 2 -s 1 -k 0
16// ex14 -m ../data/pipe-nurbs.mesh -o 1
17// ex14 -m ../data/inline-segment.mesh -r 5
18// ex14 -m ../data/amr-quad.mesh -r 3
19// ex14 -m ../data/amr-hex.mesh
20// ex14 -m ../data/fichera-amr.mesh
21// ex14 -pa -r 1 -o 3
22// ex14 -pa -r 1 -o 3 -m ../data/fichera.mesh
23//
24// Device sample runs:
25// ex14 -pa -r 2 -d cuda -o 3
26// ex14 -pa -r 2 -d cuda -o 3 -m ../data/fichera.mesh
27//
28// Description: This example code demonstrates the use of MFEM to define a
29// discontinuous Galerkin (DG) finite element discretization of
30// the Laplace problem -Delta u = 1 with homogeneous Dirichlet
31// boundary conditions. Finite element spaces of any order,
32// including zero on regular grids, are supported. The example
33// highlights the use of discontinuous spaces and DG-specific face
34// integrators.
35//
36// We recommend viewing examples 1 and 9 before viewing this
37// example.
38
39#include "mfem.hpp"
40#include <fstream>
41#include <iostream>
42
43using namespace std;
44using namespace mfem;
45
46int main(int argc, char *argv[])
47{
48 // 1. Parse command-line options.
49 const char *mesh_file = "../data/star.mesh";
50 int ref_levels = -1;
51 int order = 1;
52 real_t sigma = -1.0;
53 real_t kappa = -1.0;
54 real_t eta = 0.0;
55 bool pa = false;
56 bool visualization = 1;
57 const char *device_config = "cpu";
58
59 OptionsParser args(argc, argv);
60 args.AddOption(&mesh_file, "-m", "--mesh",
61 "Mesh file to use.");
62 args.AddOption(&ref_levels, "-r", "--refine",
63 "Number of times to refine the mesh uniformly, -1 for auto.");
64 args.AddOption(&order, "-o", "--order",
65 "Finite element order (polynomial degree) >= 0.");
66 args.AddOption(&sigma, "-s", "--sigma",
67 "One of the three DG penalty parameters, typically +1/-1."
68 " See the documentation of class DGDiffusionIntegrator.");
69 args.AddOption(&kappa, "-k", "--kappa",
70 "One of the three DG penalty parameters, should be positive."
71 " Negative values are replaced with (order+1)^2.");
72 args.AddOption(&eta, "-e", "--eta", "BR2 penalty parameter.");
73 args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
74 "--no-partial-assembly", "Enable Partial Assembly.");
75 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
76 "--no-visualization",
77 "Enable or disable GLVis visualization.");
78 args.AddOption(&device_config, "-d", "--device",
79 "Device configuration string, see Device::Configure().");
80 args.Parse();
81 if (!args.Good())
82 {
83 args.PrintUsage(cout);
84 return 1;
85 }
86 if (kappa < 0)
87 {
88 kappa = (order+1)*(order+1);
89 }
90 args.PrintOptions(cout);
91
92 // 2. Enable hardware devices such as GPUs, and programming models such as
93 // CUDA, OCCA, RAJA and OpenMP based on command line options.
94 Device device(device_config);
95 device.Print();
96
97 // 3. Read the mesh from the given mesh file. We can handle triangular,
98 // quadrilateral, tetrahedral and hexahedral meshes with the same code.
99 // NURBS meshes are projected to second order meshes.
100 Mesh mesh(mesh_file);
101 const int dim = mesh.Dimension();
102
103 // 4. Refine the mesh to increase the resolution. In this example we do
104 // 'ref_levels' of uniform refinement. By default, or if ref_levels < 0,
105 // we choose it to be the largest number that gives a final mesh with no
106 // more than 50,000 elements.
107 {
108 if (ref_levels < 0)
109 {
110 ref_levels = (int)floor(log(50000./mesh.GetNE())/log(2.)/dim);
111 }
112 for (int l = 0; l < ref_levels; l++)
113 {
114 mesh.UniformRefinement();
115 }
116 }
117 if (mesh.NURBSext)
118 {
119 mesh.SetCurvature(max(order, 1));
120 }
121
122 // 5. Define a finite element space on the mesh. Here we use discontinuous
123 // finite elements of the specified order >= 0.
125 DG_FECollection fec(order, dim, bt);
126 FiniteElementSpace fespace(&mesh, &fec);
127 cout << "Number of unknowns: " << fespace.GetVSize() << endl;
128
129 // 6. Set up the linear form b(.) which corresponds to the right-hand side of
130 // the FEM linear system.
131 LinearForm b(&fespace);
132 ConstantCoefficient one(1.0);
133 ConstantCoefficient zero(0.0);
134 b.AddDomainIntegrator(new DomainLFIntegrator(one));
135 b.AddBdrFaceIntegrator(
136 new DGDirichletLFIntegrator(zero, one, sigma, kappa));
137 b.Assemble();
138
139 // 7. Define the solution vector x as a finite element grid function
140 // corresponding to fespace. Initialize x with initial guess of zero.
141 GridFunction x(&fespace);
142 x = 0.0;
143
144 // 8. Set up the bilinear form a(.,.) on the finite element space
145 // corresponding to the Laplacian operator -Delta, by adding the Diffusion
146 // domain integrator and the interior and boundary DG face integrators.
147 // Note that boundary conditions are imposed weakly in the form, so there
148 // is no need for dof elimination. After assembly and finalizing we
149 // extract the corresponding sparse matrix A.
150 BilinearForm a(&fespace);
151 a.AddDomainIntegrator(new DiffusionIntegrator(one));
152 a.AddInteriorFaceIntegrator(new DGDiffusionIntegrator(one, sigma, kappa));
153 a.AddBdrFaceIntegrator(new DGDiffusionIntegrator(one, sigma, kappa));
154 if (eta > 0)
155 {
156 MFEM_VERIFY(!pa, "BR2 not yet compatible with partial assembly.");
157 a.AddInteriorFaceIntegrator(new DGDiffusionBR2Integrator(fespace, eta));
158 a.AddBdrFaceIntegrator(new DGDiffusionBR2Integrator(fespace, eta));
159 }
160 if (pa) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
161 a.Assemble();
162 a.Finalize();
163
164 // 9. Define a simple symmetric Gauss-Seidel preconditioner and use it to
165 // solve the system Ax=b with PCG in the symmetric case, and GMRES in the
166 // non-symmetric one. (Note that tolerances are squared: 1e-12 corresponds
167 // to a relative tolerance of 1e-6).
168 //
169 // If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
170 if (pa)
171 {
172 MFEM_VERIFY(sigma == -1.0,
173 "The case of PA with sigma != -1 is not yet supported.");
174 CG(a, b, x, 1, 500, 1e-12, 0.0);
175 }
176 else
177 {
178 const SparseMatrix &A = a.SpMat();
179#ifndef MFEM_USE_SUITESPARSE
180 GSSmoother M(A);
181 if (sigma == -1.0)
182 {
183 PCG(A, M, b, x, 1, 500, 1e-12, 0.0);
184 }
185 else
186 {
187 GMRES(A, M, b, x, 1, 500, 10, 1e-12, 0.0);
188 }
189#else
190 UMFPackSolver umf_solver;
191 umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
192 umf_solver.SetOperator(A);
193 umf_solver.Mult(b, x);
194#endif
195 }
196
197 // 10. Save the refined mesh and the solution. This output can be viewed
198 // later using GLVis: "glvis -m refined.mesh -g sol.gf".
199 ofstream mesh_ofs("refined.mesh");
200 mesh_ofs.precision(8);
201 mesh.Print(mesh_ofs);
202 ofstream sol_ofs("sol.gf");
203 sol_ofs.precision(8);
204 x.Save(sol_ofs);
205
206 // 11. Send the solution by socket to a GLVis server.
207 if (visualization)
208 {
209 char vishost[] = "localhost";
210 int visport = 19916;
211 socketstream sol_sock(vishost, visport);
212 sol_sock.precision(8);
213 sol_sock << "solution\n" << mesh << x << flush;
214 }
215
216 return 0;
217}
@ GaussLobatto
Closed type.
Definition fe_base.hpp:32
@ GaussLegendre
Open type.
Definition fe_base.hpp:31
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
A coefficient that is constant across space and time.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition device.hpp:123
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition device.cpp:286
Class for domain integration .
Definition lininteg.hpp:109
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:713
Data type for Gauss-Seidel smoother of sparse matrix.
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
Vector with associated FE space and LinearFormIntegrators.
Mesh data type.
Definition mesh.hpp:56
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition mesh.hpp:290
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Definition mesh.hpp:2288
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
Definition mesh.cpp:6211
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
Definition optparser.hpp:82
bool Good() const
Return true if the command line options were parsed successfully.
Data type sparse matrix.
Definition sparsemat.hpp:51
Direct sparse solver using UMFPACK.
Definition solvers.hpp:1096
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition solvers.cpp:3102
real_t Control[UMFPACK_CONTROL]
Definition solvers.hpp:1106
virtual void Mult(const Vector &b, Vector &x) const
Direct solution of the linear system using UMFPACK.
Definition solvers.cpp:3197
real_t sigma(const Vector &x)
Definition maxwell.cpp:102
real_t kappa
Definition ex24.cpp:54
int dim
Definition ex24.cpp:53
int main()
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
const int visport
int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, int m, real_t &tol, real_t atol, int printit)
GMRES method. (tolerances are squared)
Definition solvers.cpp:1342
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
Definition solvers.cpp:913
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
Definition solvers.cpp:898
float real_t
Definition config.hpp:43
const char vishost[]