MFEM  v3.4
Finite element discretization library
 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 -o 2 -hb
11 // ex4 -m ../data/fichera-q2.vtk
12 // ex4 -m ../data/fichera-q3.mesh -o 2 -sc
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 // ex4 -m ../data/amr-quad.mesh
18 // ex4 -m ../data/amr-hex.mesh
19 // ex4 -m ../data/amr-hex.mesh -o 2 -hb
20 // ex4 -m ../data/fichera-amr.mesh -o 2 -sc
21 // ex4 -m ../data/star-surf.mesh -o 1
22 //
23 // Description: This example code solves a simple 2D/3D H(div) diffusion
24 // problem corresponding to the second order definite equation
25 // -grad(alpha div F) + beta F = f with boundary condition F dot n
26 // = <given normal field>. Here, we use a given exact solution F
27 // and compute the corresponding r.h.s. f. We discretize with
28 // Raviart-Thomas finite elements.
29 //
30 // The example demonstrates the use of H(div) finite element
31 // spaces with the grad-div and H(div) vector finite element mass
32 // bilinear form, as well as the computation of discretization
33 // error when the exact solution is known. Bilinear form
34 // hybridization and static condensation are also illustrated.
35 //
36 // We recommend viewing examples 1-3 before viewing this example.
37 
38 #include "mfem.hpp"
39 #include <fstream>
40 #include <iostream>
41 
42 using namespace std;
43 using namespace mfem;
44 
45 // Exact solution, F, and r.h.s., f. See below for implementation.
46 void F_exact(const Vector &, Vector &);
47 void f_exact(const Vector &, Vector &);
48 double freq = 1.0, kappa;
49 
50 int main(int argc, char *argv[])
51 {
52  // 1. Parse command-line options.
53  const char *mesh_file = "../data/star.mesh";
54  int order = 1;
55  bool set_bc = true;
56  bool static_cond = false;
57  bool hybridization = false;
58  bool visualization = 1;
59 
60  OptionsParser args(argc, argv);
61  args.AddOption(&mesh_file, "-m", "--mesh",
62  "Mesh file to use.");
63  args.AddOption(&order, "-o", "--order",
64  "Finite element order (polynomial degree).");
65  args.AddOption(&set_bc, "-bc", "--impose-bc", "-no-bc", "--dont-impose-bc",
66  "Impose or not essential boundary conditions.");
67  args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
68  " solution.");
69  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
70  "--no-static-condensation", "Enable static condensation.");
71  args.AddOption(&hybridization, "-hb", "--hybridization", "-no-hb",
72  "--no-hybridization", "Enable hybridization.");
73  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
74  "--no-visualization",
75  "Enable or disable GLVis visualization.");
76  args.Parse();
77  if (!args.Good())
78  {
79  args.PrintUsage(cout);
80  return 1;
81  }
82  args.PrintOptions(cout);
83  kappa = freq * M_PI;
84 
85  // 2. Read the mesh from the given mesh file. We can handle triangular,
86  // quadrilateral, tetrahedral, hexahedral, surface and volume, as well as
87  // periodic meshes with the same code.
88  Mesh *mesh = new Mesh(mesh_file, 1, 1);
89  int dim = mesh->Dimension();
90  int sdim = mesh->SpaceDimension();
91 
92  // 3. Refine the mesh to increase the resolution. In this example we do
93  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
94  // largest number that gives a final mesh with no more than 25,000
95  // elements.
96  {
97  int ref_levels =
98  (int)floor(log(25000./mesh->GetNE())/log(2.)/dim);
99  for (int l = 0; l < ref_levels; l++)
100  {
101  mesh->UniformRefinement();
102  }
103  }
104 
105  // 4. Define a finite element space on the mesh. Here we use the
106  // Raviart-Thomas finite elements of the specified order.
107  FiniteElementCollection *fec = new RT_FECollection(order-1, dim);
108  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
109  cout << "Number of finite element unknowns: "
110  << fespace->GetTrueVSize() << endl;
111 
112  // 5. Determine the list of true (i.e. conforming) essential boundary dofs.
113  // In this example, the boundary conditions are defined by marking all
114  // the boundary attributes from the mesh as essential (Dirichlet) and
115  // converting them to a list of true dofs.
116  Array<int> ess_tdof_list;
117  if (mesh->bdr_attributes.Size())
118  {
119  Array<int> ess_bdr(mesh->bdr_attributes.Max());
120  ess_bdr = set_bc ? 1 : 0;
121  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
122  }
123 
124  // 6. Set up the linear form b(.) which corresponds to the right-hand side
125  // of the FEM linear system, which in this case is (f,phi_i) where f is
126  // given by the function f_exact and phi_i are the basis functions in the
127  // finite element fespace.
129  LinearForm *b = new LinearForm(fespace);
131  b->Assemble();
132 
133  // 7. Define the solution vector x as a finite element grid function
134  // corresponding to fespace. Initialize x by projecting the exact
135  // solution. Note that only values from the boundary faces will be used
136  // when eliminating the non-homogeneous boundary condition to modify the
137  // r.h.s. vector b.
138  GridFunction x(fespace);
140  x.ProjectCoefficient(F);
141 
142  // 8. Set up the bilinear form corresponding to the H(div) diffusion operator
143  // grad alpha div + beta I, by adding the div-div and the mass domain
144  // integrators.
146  Coefficient *beta = new ConstantCoefficient(1.0);
147  BilinearForm *a = new BilinearForm(fespace);
148  a->AddDomainIntegrator(new DivDivIntegrator(*alpha));
150 
151  // 9. Assemble the bilinear form and the corresponding linear system,
152  // applying any necessary transformations such as: eliminating boundary
153  // conditions, applying conforming constraints for non-conforming AMR,
154  // static condensation, hybridization, etc.
155  FiniteElementCollection *hfec = NULL;
156  FiniteElementSpace *hfes = NULL;
157  if (static_cond)
158  {
160  }
161  else if (hybridization)
162  {
163  hfec = new DG_Interface_FECollection(order-1, dim);
164  hfes = new FiniteElementSpace(mesh, hfec);
166  ess_tdof_list);
167  }
168  a->Assemble();
169 
170  SparseMatrix A;
171  Vector B, X;
172  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
173 
174  cout << "Size of linear system: " << A.Height() << endl;
175 
176 #ifndef MFEM_USE_SUITESPARSE
177  // 10. Define a simple symmetric Gauss-Seidel preconditioner and use it to
178  // solve the system A X = B with PCG.
179  GSSmoother M(A);
180  PCG(A, M, B, X, 1, 10000, 1e-20, 0.0);
181 #else
182  // 10. If compiled with SuiteSparse support, use UMFPACK to solve the system.
183  UMFPackSolver umf_solver;
184  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
185  umf_solver.SetOperator(A);
186  umf_solver.Mult(B, X);
187 #endif
188 
189  // 11. Recover the solution as a finite element grid function.
190  a->RecoverFEMSolution(X, *b, x);
191 
192  // 12. Compute and print the L^2 norm of the error.
193  cout << "\n|| F_h - F ||_{L^2} = " << x.ComputeL2Error(F) << '\n' << endl;
194 
195  // 13. Save the refined mesh and the solution. This output can be viewed
196  // later using GLVis: "glvis -m refined.mesh -g sol.gf".
197  {
198  ofstream mesh_ofs("refined.mesh");
199  mesh_ofs.precision(8);
200  mesh->Print(mesh_ofs);
201  ofstream sol_ofs("sol.gf");
202  sol_ofs.precision(8);
203  x.Save(sol_ofs);
204  }
205 
206  // 14. Send the solution by socket to a GLVis server.
207  if (visualization)
208  {
209  char vishost[] = "localhost";
210  int visport = 19916;
211  socketstream sol_sock(vishost, visport);
212  sol_sock.precision(8);
213  sol_sock << "solution\n" << *mesh << x << flush;
214  }
215 
216  // 15. Free the used memory.
217  delete hfes;
218  delete hfec;
219  delete a;
220  delete alpha;
221  delete beta;
222  delete b;
223  delete fespace;
224  delete fec;
225  delete mesh;
226 
227  return 0;
228 }
229 
230 
231 // The exact solution (for non-surface meshes)
232 void F_exact(const Vector &p, Vector &F)
233 {
234  int dim = p.Size();
235 
236  double x = p(0);
237  double y = p(1);
238  // double z = (dim == 3) ? p(2) : 0.0;
239 
240  F(0) = cos(kappa*x)*sin(kappa*y);
241  F(1) = cos(kappa*y)*sin(kappa*x);
242  if (dim == 3)
243  {
244  F(2) = 0.0;
245  }
246 }
247 
248 // The right hand side
249 void f_exact(const Vector &p, Vector &f)
250 {
251  int dim = p.Size();
252 
253  double x = p(0);
254  double y = p(1);
255  // double z = (dim == 3) ? p(2) : 0.0;
256 
257  double temp = 1 + 2*kappa*kappa;
258 
259  f(0) = temp*cos(kappa*x)*sin(kappa*y);
260  f(1) = temp*cos(kappa*y)*sin(kappa*x);
261  if (dim == 3)
262  {
263  f(2) = 0;
264  }
265 }
int Size() const
Logical size of the array.
Definition: array.hpp:133
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1034
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
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:263
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:51
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:120
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: fespace.cpp:365
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:618
int main(int argc, char *argv[])
Definition: ex1.cpp:45
(Q div u, div v) for RT elements
Data type for Gauss-Seidel smoother of sparse matrix.
void EnableHybridization(FiniteElementSpace *constr_space, BilinearFormIntegrator *constr_integ, const Array< int > &ess_tdof_list)
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:2349
int dim
Definition: ex3.cpp:47
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:36
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6741
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:109
void F_exact(const Vector &, Vector &)
Definition: ex4.cpp:232
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:457
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:256
int Dimension() const
Definition: mesh.hpp:645
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:172
int SpaceDimension() const
Definition: mesh.hpp:646
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:174
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:66
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
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1377
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:231
const double alpha
Definition: ex15.cpp:337
double kappa
Definition: ex3.cpp:46
Vector data type.
Definition: vector.hpp:48
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, SparseMatrix &A, Vector &X, Vector &B, int copy_interior=0)
Form a linear system, A X = B.
double freq
Definition: ex3.cpp:46
Integrator for (Q u, v) for VectorFiniteElements.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1785
void EnableStaticCondensation()
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:1688
bool Good() const
Definition: optparser.hpp:120