MFEM  v4.0
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/escher.mesh -o 2
11 // ex3 -m ../data/fichera.mesh
12 // ex3 -m ../data/fichera-q2.vtk
13 // ex3 -m ../data/fichera-q3.mesh
14 // ex3 -m ../data/square-disc-nurbs.mesh
15 // ex3 -m ../data/beam-hex-nurbs.mesh
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 = new Mesh(mesh_file, 1, 1);
83  dim = mesh->Dimension();
84  int sdim = mesh->SpaceDimension();
85 
86  // 3. Refine the mesh to increase the resolution. In this example we do
87  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
88  // largest number that gives a final mesh with no more than 50,000
89  // elements.
90  {
91  int ref_levels =
92  (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
93  for (int l = 0; l < ref_levels; l++)
94  {
95  mesh->UniformRefinement();
96  }
97  }
98  mesh->ReorientTetMesh();
99 
100  // 4. Define a finite element space on the mesh. Here we use the Nedelec
101  // finite elements of the specified order.
102  FiniteElementCollection *fec = new ND_FECollection(order, dim);
103  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
104  cout << "Number of finite element unknowns: "
105  << fespace->GetTrueVSize() << endl;
106 
107  // 5. Determine the list of true (i.e. conforming) essential boundary dofs.
108  // In this example, the boundary conditions are defined by marking all
109  // the boundary attributes from the mesh as essential (Dirichlet) and
110  // converting them to a list of true dofs.
111  Array<int> ess_tdof_list;
112  if (mesh->bdr_attributes.Size())
113  {
114  Array<int> ess_bdr(mesh->bdr_attributes.Max());
115  ess_bdr = 1;
116  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
117  }
118 
119  // 6. Set up the linear form b(.) which corresponds to the right-hand side
120  // of the FEM linear system, which in this case is (f,phi_i) where f is
121  // given by the function f_exact and phi_i are the basis functions in the
122  // finite element fespace.
124  LinearForm *b = new LinearForm(fespace);
126  b->Assemble();
127 
128  // 7. Define the solution vector x as a finite element grid function
129  // corresponding to fespace. Initialize x by projecting the exact
130  // solution. Note that only values from the boundary edges will be used
131  // when eliminating the non-homogeneous boundary condition to modify the
132  // r.h.s. vector b.
133  GridFunction x(fespace);
135  x.ProjectCoefficient(E);
136 
137  // 8. Set up the bilinear form corresponding to the EM diffusion operator
138  // curl muinv curl + sigma I, by adding the curl-curl and the mass domain
139  // integrators.
140  Coefficient *muinv = new ConstantCoefficient(1.0);
142  BilinearForm *a = new BilinearForm(fespace);
143  a->AddDomainIntegrator(new CurlCurlIntegrator(*muinv));
145 
146  // 9. Assemble the bilinear form and the corresponding linear system,
147  // applying any necessary transformations such as: eliminating boundary
148  // conditions, applying conforming constraints for non-conforming AMR,
149  // static condensation, etc.
150  if (static_cond) { a->EnableStaticCondensation(); }
151  a->Assemble();
152 
153  SparseMatrix A;
154  Vector B, X;
155  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
156 
157  cout << "Size of linear system: " << A.Height() << endl;
158 
159 #ifndef MFEM_USE_SUITESPARSE
160  // 10. Define a simple symmetric Gauss-Seidel preconditioner and use it to
161  // solve the system Ax=b with PCG.
162  GSSmoother M(A);
163  PCG(A, M, B, X, 1, 500, 1e-12, 0.0);
164 #else
165  // 10. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
166  UMFPackSolver umf_solver;
167  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
168  umf_solver.SetOperator(A);
169  umf_solver.Mult(B, X);
170 #endif
171 
172  // 11. Recover the solution as a finite element grid function.
173  a->RecoverFEMSolution(X, *b, x);
174 
175  // 12. Compute and print the L^2 norm of the error.
176  cout << "\n|| E_h - E ||_{L^2} = " << x.ComputeL2Error(E) << '\n' << endl;
177 
178  // 13. Save the refined mesh and the solution. This output can be viewed
179  // later using GLVis: "glvis -m refined.mesh -g sol.gf".
180  {
181  ofstream mesh_ofs("refined.mesh");
182  mesh_ofs.precision(8);
183  mesh->Print(mesh_ofs);
184  ofstream sol_ofs("sol.gf");
185  sol_ofs.precision(8);
186  x.Save(sol_ofs);
187  }
188 
189  // 14. Send the solution by socket to a GLVis server.
190  if (visualization)
191  {
192  char vishost[] = "localhost";
193  int visport = 19916;
194  socketstream sol_sock(vishost, visport);
195  sol_sock.precision(8);
196  sol_sock << "solution\n" << *mesh << x << flush;
197  }
198 
199  // 15. Free the used memory.
200  delete a;
201  delete sigma;
202  delete muinv;
203  delete b;
204  delete fespace;
205  delete fec;
206  delete mesh;
207 
208  return 0;
209 }
210 
211 
212 void E_exact(const Vector &x, Vector &E)
213 {
214  if (dim == 3)
215  {
216  E(0) = sin(kappa * x(1));
217  E(1) = sin(kappa * x(2));
218  E(2) = sin(kappa * x(0));
219  }
220  else
221  {
222  E(0) = sin(kappa * x(1));
223  E(1) = sin(kappa * x(0));
224  if (x.Size() == 3) { E(2) = 0.0; }
225  }
226 }
227 
228 void f_exact(const Vector &x, Vector &f)
229 {
230  if (dim == 3)
231  {
232  f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
233  f(1) = (1. + kappa * kappa) * sin(kappa * x(2));
234  f(2) = (1. + kappa * kappa) * sin(kappa * x(0));
235  }
236  else
237  {
238  f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
239  f(1) = (1. + kappa * kappa) * sin(kappa * x(0));
240  if (x.Size() == 3) { f(2) = 0.0; }
241  }
242 }
int Size() const
Logical size of the array.
Definition: array.hpp:118
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1156
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:67
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:306
Integrator for (curl u, curl v) for Nedelec elements.
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Recover the solution of a linear system formed with FormLinearSystem().
int Size() const
Returns the size of the vector.
Definition: vector.hpp:150
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: fespace.cpp:366
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:676
int main(int argc, char *argv[])
Definition: ex1.cpp:58
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:344
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2627
int dim
Definition: ex3.cpp:48
Data type sparse matrix.
Definition: sparsemat.hpp:40
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:36
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7591
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
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:461
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:350
int Dimension() const
Definition: mesh.hpp:713
virtual void ReorientTetMesh()
Definition: mesh.cpp:4915
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
int SpaceDimension() const
Definition: mesh.hpp:714
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:179
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:355
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:85
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:76
void f_exact(const Vector &, Vector &)
Definition: ex3.cpp:228
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1582
void E_exact(const Vector &, Vector &)
Definition: ex3.cpp:212
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:231
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:240
double kappa
Definition: ex3.cpp:47
Vector data type.
Definition: vector.hpp:48
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.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1794
void EnableStaticCondensation()
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:1697
double sigma(const Vector &x)
Definition: maxwell.cpp:102
bool Good() const
Definition: optparser.hpp:122