MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex3.cpp
Go to the documentation of this file.
1// MFEM Example 3
2//
3// Compile with: make ex3
4//
5// Sample runs: ex3 -m ../data/star.mesh
6// ex3 -m ../data/beam-tri.mesh -o 2
7// ex3 -m ../data/beam-tet.mesh
8// ex3 -m ../data/beam-hex.mesh
9// ex3 -m ../data/beam-hex.mesh -o 2 -pa
10// ex3 -m ../data/escher.mesh
11// ex3 -m ../data/escher.mesh -o 2
12// ex3 -m ../data/fichera.mesh
13// ex3 -m ../data/fichera-q2.vtk
14// ex3 -m ../data/fichera-q3.mesh
15// ex3 -m ../data/square-disc-nurbs.mesh
16// ex3 -m ../data/beam-hex-nurbs.mesh
17// ex3 -m ../data/amr-hex.mesh
18// ex3 -m ../data/fichera-amr.mesh
19// ex3 -m ../data/ref-prism.mesh -o 1
20// ex3 -m ../data/octahedron.mesh -o 1
21// ex3 -m ../data/star-surf.mesh -o 1
22// ex3 -m ../data/mobius-strip.mesh -f 0.1
23// ex3 -m ../data/klein-bottle.mesh -f 0.1
24//
25// Device sample runs:
26// ex3 -m ../data/star.mesh -pa -d cuda
27// ex3 -m ../data/star.mesh -pa -d raja-cuda
28// ex3 -m ../data/star.mesh -pa -d raja-omp
29// ex3 -m ../data/beam-hex.mesh -pa -d cuda
30//
31// Description: This example code solves a simple electromagnetic diffusion
32// problem corresponding to the second order definite Maxwell
33// equation curl curl E + E = f with boundary condition
34// E x n = <given tangential field>. Here, we use a given exact
35// solution E and compute the corresponding r.h.s. f.
36// We discretize with Nedelec finite elements in 2D or 3D.
37//
38// The example demonstrates the use of H(curl) finite element
39// spaces with the curl-curl and the (vector finite element) mass
40// bilinear form, as well as the computation of discretization
41// error when the exact solution is known. Static condensation is
42// also illustrated.
43//
44// We recommend viewing examples 1-2 before viewing this example.
45
46#include "mfem.hpp"
47#include <fstream>
48#include <iostream>
49
50using namespace std;
51using namespace mfem;
52
53// Exact solution, E, and r.h.s., f. See below for implementation.
54void E_exact(const Vector &, Vector &);
55void f_exact(const Vector &, Vector &);
57int dim;
58
59int main(int argc, char *argv[])
60{
61 // 1. Parse command-line options.
62 const char *mesh_file = "../data/beam-tet.mesh";
63 int order = 1;
64 bool static_cond = false;
65 bool pa = false;
66 bool nc = false;
67 const char *device_config = "cpu";
68 bool visualization = 1;
69
70 OptionsParser args(argc, argv);
71 args.AddOption(&mesh_file, "-m", "--mesh",
72 "Mesh file to use.");
73 args.AddOption(&order, "-o", "--order",
74 "Finite element order (polynomial degree).");
75 args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
76 " solution.");
77 args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
78 "--no-static-condensation", "Enable static condensation.");
79 args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
80 "--no-partial-assembly", "Enable Partial Assembly.");
81 args.AddOption(&nc, "-nc", "--non-conforming", "-c",
82 "--conforming",
83 "Mark the mesh as nonconforming before partitioning.");
84 args.AddOption(&device_config, "-d", "--device",
85 "Device configuration string, see Device::Configure().");
86 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
87 "--no-visualization",
88 "Enable or disable GLVis visualization.");
89 args.Parse();
90 if (!args.Good())
91 {
92 args.PrintUsage(cout);
93 return 1;
94 }
95 args.PrintOptions(cout);
96 kappa = freq * M_PI;
97
98 // 2. Enable hardware devices such as GPUs, and programming models such as
99 // CUDA, OCCA, RAJA and OpenMP based on command line options.
100 Device device(device_config);
101 device.Print();
102
103 // 3. Read the mesh from the given mesh file. We can handle triangular,
104 // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
105 // the same code.
106 Mesh *mesh = new Mesh(mesh_file, 1, 1);
107 dim = mesh->Dimension();
108 int sdim = mesh->SpaceDimension();
109 if (nc)
110 {
111 // Can set to false to use conformal refinement for simplices.
112 mesh->EnsureNCMesh(true);
113 }
114
115 // 4. Refine the mesh to increase the resolution. In this example we do
116 // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
117 // largest number that gives a final mesh with no more than 50,000
118 // elements.
119 {
120 int ref_levels =
121 (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
122 for (int l = 0; l < ref_levels; l++)
123 {
124 mesh->UniformRefinement();
125 }
126 }
127
128 // 5. Define a finite element space on the mesh. Here we use the Nedelec
129 // finite elements of the specified order.
131 FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
132 cout << "Number of finite element unknowns: "
133 << fespace->GetTrueVSize() << endl;
134
135 // 6. Determine the list of true (i.e. conforming) essential boundary dofs.
136 // In this example, the boundary conditions are defined by marking all
137 // the boundary attributes from the mesh as essential (Dirichlet) and
138 // converting them to a list of true dofs.
139 Array<int> ess_tdof_list;
140 if (mesh->bdr_attributes.Size())
141 {
142 Array<int> ess_bdr(mesh->bdr_attributes.Max());
143 ess_bdr = 1;
144 fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
145 }
146
147 // 7. Set up the linear form b(.) which corresponds to the right-hand side
148 // of the FEM linear system, which in this case is (f,phi_i) where f is
149 // given by the function f_exact and phi_i are the basis functions in the
150 // finite element fespace.
152 LinearForm *b = new LinearForm(fespace);
153 b->AddDomainIntegrator(new VectorFEDomainLFIntegrator(f));
154 b->Assemble();
155
156 // 8. Define the solution vector x as a finite element grid function
157 // corresponding to fespace. Initialize x by projecting the exact
158 // solution. Note that only values from the boundary edges will be used
159 // when eliminating the non-homogeneous boundary condition to modify the
160 // r.h.s. vector b.
161 GridFunction x(fespace);
164
165 // 9. Set up the bilinear form corresponding to the EM diffusion operator
166 // curl muinv curl + sigma I, by adding the curl-curl and the mass domain
167 // integrators.
168 Coefficient *muinv = new ConstantCoefficient(1.0);
170 BilinearForm *a = new BilinearForm(fespace);
171 if (pa) { a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
172 a->AddDomainIntegrator(new CurlCurlIntegrator(*muinv));
173 a->AddDomainIntegrator(new VectorFEMassIntegrator(*sigma));
174
175 // 10. Assemble the bilinear form and the corresponding linear system,
176 // applying any necessary transformations such as: eliminating boundary
177 // conditions, applying conforming constraints for non-conforming AMR,
178 // static condensation, etc.
179 if (static_cond) { a->EnableStaticCondensation(); }
180 a->Assemble();
181
182 OperatorPtr A;
183 Vector B, X;
184 a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
185
186 cout << "Size of linear system: " << A->Height() << endl;
187
188 // 11. Solve the linear system A X = B.
189 if (pa) // Jacobi preconditioning in partial assembly mode
190 {
191 OperatorJacobiSmoother M(*a, ess_tdof_list);
192 PCG(*A, M, B, X, 1, 1000, 1e-12, 0.0);
193 }
194 else
195 {
196#ifndef MFEM_USE_SUITESPARSE
197 // 11. Define a simple symmetric Gauss-Seidel preconditioner and use it to
198 // solve the system Ax=b with PCG.
199 GSSmoother M((SparseMatrix&)(*A));
200 PCG(*A, M, B, X, 1, 500, 1e-12, 0.0);
201#else
202 // 11. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the
203 // system.
204 UMFPackSolver umf_solver;
205 umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
206 umf_solver.SetOperator(*A);
207 umf_solver.Mult(B, X);
208#endif
209 }
210
211 // 12. Recover the solution as a finite element grid function.
212 a->RecoverFEMSolution(X, *b, x);
213
214 // 13. Compute and print the L^2 norm of the error.
215 cout << "\n|| E_h - E ||_{L^2} = " << x.ComputeL2Error(E) << '\n' << endl;
216
217 // 14. Save the refined mesh and the solution. This output can be viewed
218 // later using GLVis: "glvis -m refined.mesh -g sol.gf".
219 {
220 ofstream mesh_ofs("refined.mesh");
221 mesh_ofs.precision(8);
222 mesh->Print(mesh_ofs);
223 ofstream sol_ofs("sol.gf");
224 sol_ofs.precision(8);
225 x.Save(sol_ofs);
226 }
227
228 // 15. Send the solution by socket to a GLVis server.
229 if (visualization)
230 {
231 char vishost[] = "localhost";
232 int visport = 19916;
233 socketstream sol_sock(vishost, visport);
234 sol_sock.precision(8);
235 sol_sock << "solution\n" << *mesh << x << flush;
236 }
237
238 // 16. Free the used memory.
239 delete a;
240 delete sigma;
241 delete muinv;
242 delete b;
243 delete fespace;
244 delete fec;
245 delete mesh;
246
247 return 0;
248}
249
250
251void E_exact(const Vector &x, Vector &E)
252{
253 if (dim == 3)
254 {
255 E(0) = sin(kappa * x(1));
256 E(1) = sin(kappa * x(2));
257 E(2) = sin(kappa * x(0));
258 }
259 else
260 {
261 E(0) = sin(kappa * x(1));
262 E(1) = sin(kappa * x(0));
263 if (x.Size() == 3) { E(2) = 0.0; }
264 }
265}
266
267void f_exact(const Vector &x, Vector &f)
268{
269 if (dim == 3)
270 {
271 f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
272 f(1) = (1. + kappa * kappa) * sin(kappa * x(2));
273 f(2) = (1. + kappa * kappa) * sin(kappa * x(0));
274 }
275 else
276 {
277 f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
278 f(1) = (1. + kappa * kappa) * sin(kappa * x(0));
279 if (x.Size() == 3) { f(2) = 0.0; }
280 }
281}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
int Size() const
Return the logical size of the array.
Definition array.hpp:144
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
A coefficient that is constant across space and time.
Integrator for for Nedelec elements.
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
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:716
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition fespace.cpp:588
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.
virtual real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Vector with associated FE space and LinearFormIntegrators.
Mesh data type.
Definition mesh.hpp:56
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:282
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
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1163
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition mesh.cpp:10626
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition fe_coll.hpp:465
Pointer to an Operator of a specified type.
Definition handle.hpp:34
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition solvers.hpp:313
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
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
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition lininteg.hpp:347
A general vector function coefficient.
Vector data type.
Definition vector.hpp:80
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
real_t sigma(const Vector &x)
Definition maxwell.cpp:102
real_t kappa
Definition ex3.cpp:56
int dim
Definition ex3.cpp:57
void E_exact(const Vector &, Vector &)
Definition ex3.cpp:251
void f_exact(const Vector &, Vector &)
Definition ex3.cpp:267
real_t freq
Definition ex3.cpp:56
int main()
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
const int visport
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
float real_t
Definition config.hpp:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
const char vishost[]