MFEM  v3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
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/escher.mesh
10 // ex3 -m ../data/fichera.mesh
11 // ex3 -m ../data/fichera-q2.vtk
12 // ex3 -m ../data/fichera-q3.mesh
13 // ex3 -m ../data/square-disc-nurbs.mesh
14 // ex3 -m ../data/beam-hex-nurbs.mesh
15 // ex3 -m ../data/amr-hex.mesh
16 // ex3 -m ../data/fichera-amr.mesh
17 // ex3 -m ../data/star-surf.mesh -o 1
18 // ex3 -m ../data/mobius-strip.mesh -f 0.1
19 // ex3 -m ../data/klein-bottle.mesh -f 0.1
20 //
21 // Description: This example code solves a simple electromagnetic diffusion
22 // problem corresponding to the second order definite Maxwell
23 // equation curl curl E + E = f with boundary condition
24 // E x n = <given tangential field>. Here, we use a given exact
25 // solution E and compute the corresponding r.h.s. f.
26 // We discretize with Nedelec finite elements in 2D or 3D.
27 //
28 // The example demonstrates the use of H(curl) finite element
29 // spaces with the curl-curl and the (vector finite element) mass
30 // bilinear form, as well as the computation of discretization
31 // error when the exact solution is known. Static condensation is
32 // also illustrated.
33 //
34 // We recommend viewing examples 1-2 before viewing this example.
35 
36 #include "mfem.hpp"
37 #include <fstream>
38 #include <iostream>
39 
40 using namespace std;
41 using namespace mfem;
42 
43 // Exact solution, E, and r.h.s., f. See below for implementation.
44 void E_exact(const Vector &, Vector &);
45 void f_exact(const Vector &, Vector &);
46 double freq = 1.0, kappa;
47 int dim;
48 
49 int main(int argc, char *argv[])
50 {
51  // 1. Parse command-line options.
52  const char *mesh_file = "../data/beam-tet.mesh";
53  int order = 1;
54  bool static_cond = false;
55  bool visualization = 1;
56 
57  OptionsParser args(argc, argv);
58  args.AddOption(&mesh_file, "-m", "--mesh",
59  "Mesh file to use.");
60  args.AddOption(&order, "-o", "--order",
61  "Finite element order (polynomial degree).");
62  args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
63  " solution.");
64  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
65  "--no-static-condensation", "Enable static condensation.");
66  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
67  "--no-visualization",
68  "Enable or disable GLVis visualization.");
69  args.Parse();
70  if (!args.Good())
71  {
72  args.PrintUsage(cout);
73  return 1;
74  }
75  args.PrintOptions(cout);
76  kappa = freq * M_PI;
77 
78  // 2. Read the mesh from the given mesh file. We can handle triangular,
79  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
80  // the same code.
81  Mesh *mesh = new Mesh(mesh_file, 1, 1);
82  dim = mesh->Dimension();
83  int sdim = mesh->SpaceDimension();
84 
85  // 3. Refine the mesh to increase the resolution. In this example we do
86  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
87  // largest number that gives a final mesh with no more than 50,000
88  // elements.
89  {
90  int ref_levels =
91  (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
92  for (int l = 0; l < ref_levels; l++)
93  {
94  mesh->UniformRefinement();
95  }
96  }
97  mesh->ReorientTetMesh();
98 
99  // 4. Define a finite element space on the mesh. Here we use the Nedelec
100  // finite elements of the specified order.
101  FiniteElementCollection *fec = new ND_FECollection(order, dim);
102  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
103  cout << "Number of finite element unknowns: "
104  << fespace->GetTrueVSize() << endl;
105 
106  // 5. Determine the list of true (i.e. conforming) essential boundary dofs.
107  // In this example, the boundary conditions are defined by marking all
108  // the boundary attributes from the mesh as essential (Dirichlet) and
109  // converting them to a list of true dofs.
110  Array<int> ess_tdof_list;
111  if (mesh->bdr_attributes.Size())
112  {
113  Array<int> ess_bdr(mesh->bdr_attributes.Max());
114  ess_bdr = 1;
115  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
116  }
117 
118  // 6. Set up the linear form b(.) which corresponds to the right-hand side
119  // of the FEM linear system, which in this case is (f,phi_i) where f is
120  // given by the function f_exact and phi_i are the basis functions in the
121  // finite element fespace.
123  LinearForm *b = new LinearForm(fespace);
125  b->Assemble();
126 
127  // 7. Define the solution vector x as a finite element grid function
128  // corresponding to fespace. Initialize x by projecting the exact
129  // solution. Note that only values from the boundary edges will be used
130  // when eliminating the non-homogeneous boundary condition to modify the
131  // r.h.s. vector b.
132  GridFunction x(fespace);
134  x.ProjectCoefficient(E);
135 
136  // 8. Set up the bilinear form corresponding to the EM diffusion operator
137  // curl muinv curl + sigma I, by adding the curl-curl and the mass domain
138  // integrators.
139  Coefficient *muinv = new ConstantCoefficient(1.0);
140  Coefficient *sigma = new ConstantCoefficient(1.0);
141  BilinearForm *a = new BilinearForm(fespace);
142  a->AddDomainIntegrator(new CurlCurlIntegrator(*muinv));
144 
145  // 9. Assemble the bilinear form and the corresponding linear system,
146  // applying any necessary transformations such as: eliminating boundary
147  // conditions, applying conforming constraints for non-conforming AMR,
148  // static condensation, etc.
149  if (static_cond) { a->EnableStaticCondensation(); }
150  a->Assemble();
151 
152  SparseMatrix A;
153  Vector B, X;
154  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
155 
156  cout << "Size of linear system: " << A.Height() << endl;
157 
158 #ifndef MFEM_USE_SUITESPARSE
159  // 10. Define a simple symmetric Gauss-Seidel preconditioner and use it to
160  // solve the system Ax=b with PCG.
161  GSSmoother M(A);
162  PCG(A, M, B, X, 1, 500, 1e-12, 0.0);
163 #else
164  // 10. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
165  UMFPackSolver umf_solver;
166  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
167  umf_solver.SetOperator(A);
168  umf_solver.Mult(B, X);
169 #endif
170 
171  // 11. Recover the solution as a finite element grid function.
172  a->RecoverFEMSolution(X, *b, x);
173 
174  // 12. Compute and print the L^2 norm of the error.
175  cout << "\n|| E_h - E ||_{L^2} = " << x.ComputeL2Error(E) << '\n' << endl;
176 
177  // 13. Save the refined mesh and the solution. This output can be viewed
178  // later using GLVis: "glvis -m refined.mesh -g sol.gf".
179  {
180  ofstream mesh_ofs("refined.mesh");
181  mesh_ofs.precision(8);
182  mesh->Print(mesh_ofs);
183  ofstream sol_ofs("sol.gf");
184  sol_ofs.precision(8);
185  x.Save(sol_ofs);
186  }
187 
188  // 14. Send the solution by socket to a GLVis server.
189  if (visualization)
190  {
191  char vishost[] = "localhost";
192  int visport = 19916;
193  socketstream sol_sock(vishost, visport);
194  sol_sock.precision(8);
195  sol_sock << "solution\n" << *mesh << x << flush;
196  }
197 
198  // 15. Free the used memory.
199  delete a;
200  delete sigma;
201  delete muinv;
202  delete b;
203  delete fespace;
204  delete fec;
205  delete mesh;
206 
207  return 0;
208 }
209 
210 
211 void E_exact(const Vector &x, Vector &E)
212 {
213  if (dim == 3)
214  {
215  E(0) = sin(kappa * x(1));
216  E(1) = sin(kappa * x(2));
217  E(2) = sin(kappa * x(0));
218  }
219  else
220  {
221  E(0) = sin(kappa * x(1));
222  E(1) = sin(kappa * x(0));
223  if (x.Size() == 3) { E(2) = 0.0; }
224  }
225 }
226 
227 void f_exact(const Vector &x, Vector &f)
228 {
229  if (dim == 3)
230  {
231  f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
232  f(1) = (1. + kappa * kappa) * sin(kappa * x(2));
233  f(2) = (1. + kappa * kappa) * sin(kappa * x(0));
234  }
235  else
236  {
237  f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
238  f(1) = (1. + kappa * kappa) * sin(kappa * x(0));
239  if (x.Size() == 3) { f(2) = 0.0; }
240  }
241 }
int Size() const
Logical size of the array.
Definition: array.hpp:109
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
Subclass constant coefficient.
Definition: coefficient.hpp:57
Integrator for (curl u, curl v) for Nedelec elements.
Definition: bilininteg.hpp:408
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:34
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list)
Definition: fespace.cpp:274
void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
int Size() const
Returns the size of the vector.
Definition: vector.hpp:86
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:496
int main(int argc, char *argv[])
void FormLinearSystem(Array< int > &ess_tdof_list, Vector &x, Vector &b, SparseMatrix &A, Vector &X, Vector &B, int copy_interior=0)
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:333
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2116
int dim
Definition: ex3.cpp:47
double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:193
Data type sparse matrix.
Definition: sparsemat.hpp:38
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows.
Definition: operator.hpp:35
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6402
T Max() const
Definition: array.cpp:90
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:456
int Dimension() const
Definition: mesh.hpp:523
virtual void ReorientTetMesh()
Definition: mesh.cpp:4184
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
int SpaceDimension() const
Definition: mesh.hpp:524
virtual void Print(std::ostream &out=std::cout) const
Print the mesh to the given stream using the default MFEM mesh format.
Definition: mesh.cpp:6736
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Array< int > bdr_attributes
Definition: mesh.hpp:141
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:344
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
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)
Definition: optparser.hpp:74
void f_exact(const Vector &, Vector &)
Definition: ex3.cpp:227
virtual int GetTrueVSize()
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:164
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1164
void E_exact(const Vector &, Vector &)
Definition: ex3.cpp:211
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:171
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:231
double kappa
Definition: ex3.cpp:46
Vector data type.
Definition: vector.hpp:33
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
double freq
Definition: ex3.cpp:46
Integrator for (Q u, v) for VectorFiniteElements.
Definition: bilininteg.hpp:465
virtual void Mult(const Vector &b, Vector &x) const
Operator application.
Definition: solvers.cpp:1741
void EnableStaticCondensation()
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: solvers.cpp:1641
bool Good() const
Definition: optparser.hpp:120