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