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