MFEM  v3.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex4.cpp
Go to the documentation of this file.
1 // MFEM Example 4
2 //
3 // Compile with: make ex4
4 //
5 // Sample runs: ex4 -m ../data/square-disc.mesh
6 // ex4 -m ../data/star.mesh
7 // ex4 -m ../data/beam-tet.mesh
8 // ex4 -m ../data/beam-hex.mesh
9 // ex4 -m ../data/escher.mesh
10 // ex4 -m ../data/fichera.mesh
11 // ex4 -m ../data/fichera-q2.vtk
12 // ex4 -m ../data/fichera-q3.mesh
13 // ex4 -m ../data/square-disc-nurbs.mesh
14 // ex4 -m ../data/beam-hex-nurbs.mesh
15 // ex4 -m ../data/periodic-square.mesh -no-bc
16 // ex4 -m ../data/periodic-cube.mesh -no-bc
17 //
18 // Description: This example code solves a simple 2D/3D H(div) diffusion
19 // problem corresponding to the second order definite equation
20 // -grad(alpha div F) + beta F = f with boundary condition F dot n
21 // = <given normal field>. Here, we use a given exact solution F
22 // and compute the corresponding r.h.s. f. We discretize with the
23 // Raviart-Thomas finite elements.
24 //
25 // The example demonstrates the use of H(div) finite element
26 // spaces with the grad-div and H(div) vector finite element mass
27 // bilinear form, as well as the computation of discretization
28 // error when the exact solution is known.
29 //
30 // We recommend viewing examples 1-3 before viewing this example.
31 
32 #include "mfem.hpp"
33 #include <fstream>
34 #include <iostream>
35 
36 using namespace std;
37 using namespace mfem;
38 
39 // Exact solution, F, and r.h.s., f. See below for implementation.
40 void F_exact(const Vector &, Vector &);
41 void f_exact(const Vector &, Vector &);
42 
43 int main(int argc, char *argv[])
44 {
45  // 1. Parse command-line options.
46  const char *mesh_file = "../data/star.mesh";
47  int order = 1;
48  bool set_bc = true;
49  bool visualization = 1;
50 
51  OptionsParser args(argc, argv);
52  args.AddOption(&mesh_file, "-m", "--mesh",
53  "Mesh file to use.");
54  args.AddOption(&order, "-o", "--order",
55  "Finite element order (polynomial degree).");
56  args.AddOption(&set_bc, "-bc", "--impose-bc", "-no-bc", "--dont-impose-bc",
57  "Impose or not essential boundary conditions.");
58  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
59  "--no-visualization",
60  "Enable or disable GLVis visualization.");
61  args.Parse();
62  if (!args.Good())
63  {
64  args.PrintUsage(cout);
65  return 1;
66  }
67  args.PrintOptions(cout);
68 
69  // 2. Read the mesh from the given mesh file. We can handle triangular,
70  // quadrilateral, tetrahedral, hexahedral, surface and volume, as well as
71  // periodic meshes with the same code.
72  Mesh *mesh;
73  ifstream imesh(mesh_file);
74  if (!imesh)
75  {
76  cerr << "\nCan not open mesh file: " << mesh_file << '\n' << endl;
77  return 2;
78  }
79  mesh = new Mesh(imesh, 1, 1);
80  imesh.close();
81  int dim = mesh->Dimension();
82 
83  // 3. Refine the mesh to increase the resolution. In this example we do
84  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
85  // largest number that gives a final mesh with no more than 25,000
86  // elements.
87  {
88  int ref_levels =
89  (int)floor(log(25000./mesh->GetNE())/log(2.)/dim);
90  for (int l = 0; l < ref_levels; l++)
91  mesh->UniformRefinement();
92  }
93 
94  // 4. Define a finite element space on the mesh. Here we use the lowest order
95  // Raviart-Thomas finite elements, but we can easily switch to
96  // higher-order spaces by changing the value of p.
97  FiniteElementCollection *fec = new RT_FECollection(order-1, dim);
98  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
99  cout << "Number of unknowns: " << fespace->GetVSize() << endl;
100 
101  // 5. Set up the linear form b(.) which corresponds to the right-hand side
102  // of the FEM linear system, which in this case is (f,phi_i) where f is
103  // given by the function f_exact and phi_i are the basis functions in the
104  // finite element fespace.
106  LinearForm *b = new LinearForm(fespace);
108  b->Assemble();
109 
110  // 6. Define the solution vector x as a finite element grid function
111  // corresponding to fespace. Initialize x by projecting the exact
112  // solution. Note that only values from the boundary faces will be used
113  // when eliminating the non-homogeneous boundary condition to modify the
114  // r.h.s. vector b.
115  GridFunction x(fespace);
117  x.ProjectCoefficient(F);
118 
119  // 7. Set up the bilinear form corresponding to the H(div) diffusion operator
120  // grad alpha div + beta I, by adding the div-div and the mass domain
121  // integrators and finally imposing the non-homogeneous Dirichlet boundary
122  // conditions. The boundary conditions are implemented by marking all the
123  // boundary attributes from the mesh as essential (Dirichlet). After
124  // assembly and finalizing we extract the corresponding sparse matrix A.
125  Coefficient *alpha = new ConstantCoefficient(1.0);
126  Coefficient *beta = new ConstantCoefficient(1.0);
127  BilinearForm *a = new BilinearForm(fespace);
128  a->AddDomainIntegrator(new DivDivIntegrator(*alpha));
130  a->Assemble();
131  if (set_bc && mesh->bdr_attributes.Size())
132  {
133  Array<int> ess_bdr(mesh->bdr_attributes.Max());
134  ess_bdr = 1;
135  a->EliminateEssentialBC(ess_bdr, x, *b);
136  }
137  a->Finalize();
138  const SparseMatrix &A = a->SpMat();
139 
140 #ifndef MFEM_USE_SUITESPARSE
141  // 8. Define a simple symmetric Gauss-Seidel preconditioner and use it to
142  // solve the system Ax=b with PCG.
143  GSSmoother M(A);
144  x = 0.0;
145  PCG(A, M, *b, x, 1, 10000, 1e-20, 0.0);
146 #else
147  // 8. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
148  UMFPackSolver umf_solver;
149  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
150  umf_solver.SetOperator(A);
151  umf_solver.Mult(*b, x);
152 #endif
153 
154  // 9. Compute and print the L^2 norm of the error.
155  cout << "\n|| F_h - F ||_{L^2} = " << x.ComputeL2Error(F) << '\n' << endl;
156 
157  // 10. Save the refined mesh and the solution. This output can be viewed
158  // later using GLVis: "glvis -m refined.mesh -g sol.gf".
159  {
160  ofstream mesh_ofs("refined.mesh");
161  mesh_ofs.precision(8);
162  mesh->Print(mesh_ofs);
163  ofstream sol_ofs("sol.gf");
164  sol_ofs.precision(8);
165  x.Save(sol_ofs);
166  }
167 
168  // 11. Send the solution by socket to a GLVis server.
169  if (visualization)
170  {
171  char vishost[] = "localhost";
172  int visport = 19916;
173  socketstream sol_sock(vishost, visport);
174  sol_sock.precision(8);
175  sol_sock << "solution\n" << *mesh << x << flush;
176  }
177 
178  // 12. Free the used memory.
179  delete a;
180  delete alpha;
181  delete beta;
182  delete b;
183  delete fespace;
184  delete fec;
185  delete mesh;
186 
187  return 0;
188 }
189 
190 
191 // The exact solution
192 void F_exact(const Vector &p, Vector &F)
193 {
194  int dim = p.Size();
195 
196  double x = p(0);
197  double y = p(1);
198  // double z = (dim == 3) ? p(2) : 0.0;
199 
200  F(0) = cos(M_PI*x)*sin(M_PI*y);
201  F(1) = cos(M_PI*y)*sin(M_PI*x);
202  if (dim == 3)
203  F(2) = 0.0;
204 }
205 
206 // The right hand side
207 void f_exact(const Vector &p, Vector &f)
208 {
209  int dim = p.Size();
210 
211  double x = p(0);
212  double y = p(1);
213  // double z = (dim == 3) ? p(2) : 0.0;
214 
215  double temp = 1 + 2*M_PI*M_PI;
216 
217  f(0) = temp*cos(M_PI*x)*sin(M_PI*y);
218  f(1) = temp*cos(M_PI*y)*sin(M_PI*x);
219  if (dim == 3)
220  f(2) = 0;
221 }
int Size() const
Logical size of the array.
Definition: array.hpp:108
int GetVSize() const
Definition: fespace.hpp:151
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:26
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
Subclass constant coefficient.
Definition: coefficient.hpp:57
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:34
int Size() const
Returns the size of the vector.
Definition: vector.hpp:76
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:396
(Q div u, div v) for RT elements
Definition: bilininteg.hpp:494
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:329
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:1797
double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:166
Data type sparse matrix.
Definition: sparsemat.hpp:38
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6225
T Max() const
Definition: array.cpp:78
void F_exact(const Vector &, Vector &)
Definition: ex4.cpp:192
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:424
int Dimension() const
Definition: mesh.hpp:417
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:385
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:119
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:7136
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Array< int > bdr_attributes
Definition: mesh.hpp:305
int main(int argc, char *argv[])
Definition: ex1.cpp:39
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:340
Abstract finite element space.
Definition: fespace.hpp:61
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:195
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:266
void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:967
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:171
Vector data type.
Definition: vector.hpp:29
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
const SparseMatrix & SpMat() const
Returns a reference to the sparse martix.
void EliminateEssentialBC(Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs, int d=0)
Integrator for (Q u, v) for VectorFiniteElements.
Definition: bilininteg.hpp:434
virtual void Mult(const Vector &b, Vector &x) const
Operator application.
Definition: solvers.cpp:1597
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: solvers.cpp:1505
bool Good() const
Definition: optparser.hpp:120