MFEM  v3.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex5.cpp
Go to the documentation of this file.
1 // MFEM Example 5
2 //
3 // Compile with: make ex5
4 //
5 // Sample runs: ex5 -m ../data/square-disc.mesh
6 // ex5 -m ../data/star.mesh
7 // ex5 -m ../data/beam-tet.mesh
8 // ex5 -m ../data/beam-hex.mesh
9 // ex5 -m ../data/escher.mesh
10 // ex5 -m ../data/fichera.mesh
11 //
12 // Description: This example code solves a simple 2D/3D mixed Darcy problem
13 // corresponding to the saddle point system
14 // k*u + grad p = f
15 // - div u = g
16 // with natural boundary condition -p = <given pressure>.
17 // Here, we use a given exact solution (u,p) and compute the
18 // corresponding r.h.s. (f,g). We discretize with Raviart-Thomas
19 // finite elements (velocity u) and piecewise discontinuous
20 // polynomials (pressure p).
21 //
22 // The example demonstrates the use of the BlockMatrix class, as
23 // well as the collective saving of several grid functions in a
24 // VisIt (visit.llnl.gov) visualization format.
25 //
26 // We recommend viewing examples 1-4 before viewing this example.
27 
28 #include "mfem.hpp"
29 #include <fstream>
30 #include <iostream>
31 #include <algorithm>
32 
33 using namespace std;
34 using namespace mfem;
35 
36 // Define the analytical solution and forcing terms / boundary conditions
37 void uFun_ex(const Vector & x, Vector & u);
38 double pFun_ex(const Vector & x);
39 void fFun(const Vector & x, Vector & f);
40 double gFun(const Vector & x);
41 double f_natural(const Vector & x);
42 
43 int main(int argc, char *argv[])
44 {
45  StopWatch chrono;
46 
47  // 1. Parse command-line options.
48  const char *mesh_file = "../data/star.mesh";
49  int order = 1;
50  bool visualization = 1;
51 
52  OptionsParser args(argc, argv);
53  args.AddOption(&mesh_file, "-m", "--mesh",
54  "Mesh file to use.");
55  args.AddOption(&order, "-o", "--order",
56  "Finite element order (polynomial degree).");
57  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
58  "--no-visualization",
59  "Enable or disable GLVis visualization.");
60  args.Parse();
61  if (!args.Good())
62  {
63  args.PrintUsage(cout);
64  return 1;
65  }
66  args.PrintOptions(cout);
67 
68  // 2. Read the mesh from the given mesh file. We can handle triangular,
69  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
70  // the same code.
71  Mesh *mesh;
72  ifstream imesh(mesh_file);
73  if (!imesh)
74  {
75  cerr << "\nCan not open mesh file: " << mesh_file << '\n' << endl;
76  return 2;
77  }
78  mesh = new Mesh(imesh, 1, 1);
79  imesh.close();
80  int dim = mesh->Dimension();
81 
82  // 3. Refine the mesh to increase the resolution. In this example we do
83  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
84  // largest number that gives a final mesh with no more than 10,000
85  // elements.
86  {
87  int ref_levels =
88  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
89  for (int l = 0; l < ref_levels; l++)
90  {
91  mesh->UniformRefinement();
92  }
93  }
94 
95  // 4. Define a finite element space on the mesh. Here we use the
96  // Raviart-Thomas finite elements of the specified order.
97  FiniteElementCollection *hdiv_coll(new RT_FECollection(order, dim));
98  FiniteElementCollection *l2_coll(new L2_FECollection(order, dim));
99 
100  FiniteElementSpace *R_space = new FiniteElementSpace(mesh, hdiv_coll);
101  FiniteElementSpace *W_space = new FiniteElementSpace(mesh, l2_coll);
102 
103  // 5. Define the BlockStructure of the problem, i.e. define the array of
104  // offsets for each variable. The last component of the Array is the sum
105  // of the dimensions of each block.
106  Array<int> block_offsets(3); // number of variables + 1
107  block_offsets[0] = 0;
108  block_offsets[1] = R_space->GetVSize();
109  block_offsets[2] = W_space->GetVSize();
110  block_offsets.PartialSum();
111 
112  std::cout << "***********************************************************\n";
113  std::cout << "dim(R) = " << block_offsets[1] - block_offsets[0] << "\n";
114  std::cout << "dim(W) = " << block_offsets[2] - block_offsets[1] << "\n";
115  std::cout << "dim(R+W) = " << block_offsets.Last() << "\n";
116  std::cout << "***********************************************************\n";
117 
118  // 6. Define the coefficients, analytical solution, and rhs of the PDE.
119  ConstantCoefficient k(1.0);
120 
121  VectorFunctionCoefficient fcoeff(dim, fFun);
122  FunctionCoefficient fnatcoeff(f_natural);
123  FunctionCoefficient gcoeff(gFun);
124 
125  VectorFunctionCoefficient ucoeff(dim, uFun_ex);
127 
128  // 7. Allocate memory (x, rhs) for the analytical solution and the right hand
129  // side. Define the GridFunction u,p for the finite element solution and
130  // linear forms fform and gform for the right hand side. The data
131  // allocated by x and rhs are passed as a reference to the grid functions
132  // (u,p) and the linear forms (fform, gform).
133  BlockVector x(block_offsets), rhs(block_offsets);
134 
135  LinearForm *fform(new LinearForm);
136  fform->Update(R_space, rhs.GetBlock(0), 0);
139  fform->Assemble();
140 
141  LinearForm *gform(new LinearForm);
142  gform->Update(W_space, rhs.GetBlock(1), 0);
143  gform->AddDomainIntegrator(new DomainLFIntegrator(gcoeff));
144  gform->Assemble();
145 
146  // 8. Assemble the finite element matrices for the Darcy operator
147  //
148  // D = [ M B^T ]
149  // [ B 0 ]
150  // where:
151  //
152  // M = \int_\Omega k u_h \cdot v_h d\Omega u_h, v_h \in R_h
153  // B = -\int_\Omega \div u_h q_h d\Omega u_h \in R_h, q_h \in W_h
154  BilinearForm *mVarf(new BilinearForm(R_space));
155  MixedBilinearForm *bVarf(new MixedBilinearForm(R_space, W_space));
156 
158  mVarf->Assemble();
159  mVarf->Finalize();
160  SparseMatrix &M(mVarf->SpMat());
161 
163  bVarf->Assemble();
164  bVarf->Finalize();
165  SparseMatrix & B(bVarf->SpMat());
166  B *= -1.;
167  SparseMatrix *BT = Transpose(B);
168 
169  BlockMatrix darcyMatrix(block_offsets);
170  darcyMatrix.SetBlock(0,0, &M);
171  darcyMatrix.SetBlock(0,1, BT);
172  darcyMatrix.SetBlock(1,0, &B);
173 
174  // 9. Construct the operators for preconditioner
175  //
176  // P = [ diag(M) 0 ]
177  // [ 0 B diag(M)^-1 B^T ]
178  //
179  // Here we use Symmetric Gauss-Seidel to approximate the inverse of the
180  // pressure Schur Complement
181  SparseMatrix *MinvBt = Transpose(B);
182  Vector Md(M.Height());
183  M.GetDiag(Md);
184  for (int i = 0; i < Md.Size(); i++)
185  {
186  MinvBt->ScaleRow(i, 1./Md(i));
187  }
188  SparseMatrix *S = Mult(B, *MinvBt);
189 
190  Solver *invM, *invS;
191  invM = new DSmoother(M);
192 #ifndef MFEM_USE_SUITESPARSE
193  invS = new GSSmoother(*S);
194 #else
195  invS = new UMFPackSolver(*S);
196 #endif
197 
198  invM->iterative_mode = false;
199  invS->iterative_mode = false;
200 
201  BlockDiagonalPreconditioner darcyPrec(block_offsets);
202  darcyPrec.SetDiagonalBlock(0, invM);
203  darcyPrec.SetDiagonalBlock(1, invS);
204 
205  // 10. Solve the linear system with MINRES.
206  // Check the norm of the unpreconditioned residual.
207  int maxIter(1000);
208  double rtol(1.e-6);
209  double atol(1.e-10);
210 
211  chrono.Clear();
212  chrono.Start();
213  MINRESSolver solver;
214  solver.SetAbsTol(atol);
215  solver.SetRelTol(rtol);
216  solver.SetMaxIter(maxIter);
217  solver.SetOperator(darcyMatrix);
218  solver.SetPreconditioner(darcyPrec);
219  solver.SetPrintLevel(1);
220  x = 0.0;
221  solver.Mult(rhs, x);
222  chrono.Stop();
223 
224  if (solver.GetConverged())
225  std::cout << "MINRES converged in " << solver.GetNumIterations()
226  << " iterations with a residual norm of " << solver.GetFinalNorm() << ".\n";
227  else
228  std::cout << "MINRES did not converge in " << solver.GetNumIterations()
229  << " iterations. Residual norm is " << solver.GetFinalNorm() << ".\n";
230  std::cout << "MINRES solver took " << chrono.RealTime() << "s. \n";
231 
232  // 11. Create the grid functions u and p. Compute the L2 error norms.
233  GridFunction u, p;
234  u.Update(R_space, x.GetBlock(0), 0);
235  p.Update(W_space, x.GetBlock(1), 0);
236 
237  int order_quad = max(2, 2*order+1);
238  const IntegrationRule *irs[Geometry::NumGeom];
239  for (int i=0; i < Geometry::NumGeom; ++i)
240  {
241  irs[i] = &(IntRules.Get(i, order_quad));
242  }
243 
244  double err_u = u.ComputeL2Error(ucoeff, irs);
245  double norm_u = ComputeLpNorm(2., ucoeff, *mesh, irs);
246  double err_p = p.ComputeL2Error(pcoeff, irs);
247  double norm_p = ComputeLpNorm(2., pcoeff, *mesh, irs);
248 
249  std::cout << "|| u_h - u_ex || / || u_ex || = " << err_u / norm_u << "\n";
250  std::cout << "|| p_h - p_ex || / || p_ex || = " << err_p / norm_p << "\n";
251 
252  // 12. Save the mesh and the solution. This output can be viewed later using
253  // GLVis: "glvis -m ex5.mesh -g sol_u.gf" or "glvis -m ex5.mesh -g
254  // sol_p.gf".
255  {
256  ofstream mesh_ofs("ex5.mesh");
257  mesh_ofs.precision(8);
258  mesh->Print(mesh_ofs);
259 
260  ofstream u_ofs("sol_u.gf");
261  u_ofs.precision(8);
262  u.Save(u_ofs);
263 
264  ofstream p_ofs("sol_p.gf");
265  p_ofs.precision(8);
266  p.Save(p_ofs);
267  }
268 
269  // 13. Save data in the VisIt format
270  VisItDataCollection visit_dc("Example5", mesh);
271  visit_dc.RegisterField("velocity", &u);
272  visit_dc.RegisterField("pressure", &p);
273  visit_dc.Save();
274 
275  // 14. Send the solution by socket to a GLVis server.
276  if (visualization)
277  {
278  char vishost[] = "localhost";
279  int visport = 19916;
280  socketstream u_sock(vishost, visport);
281  u_sock.precision(8);
282  u_sock << "solution\n" << *mesh << u << "window_title 'Velocity'" << endl;
283  socketstream p_sock(vishost, visport);
284  p_sock.precision(8);
285  p_sock << "solution\n" << *mesh << p << "window_title 'Pressure'" << endl;
286  }
287 
288  // 15. Free the used memory.
289  delete fform;
290  delete gform;
291  delete invM;
292  delete invS;
293  delete S;
294  delete MinvBt;
295  delete BT;
296  delete mVarf;
297  delete bVarf;
298  delete W_space;
299  delete R_space;
300  delete l2_coll;
301  delete hdiv_coll;
302  delete mesh;
303 
304  return 0;
305 }
306 
307 
308 void uFun_ex(const Vector & x, Vector & u)
309 {
310  double xi(x(0));
311  double yi(x(1));
312  double zi(0.0);
313  if (x.Size() == 3)
314  {
315  zi = x(2);
316  }
317 
318  u(0) = - exp(xi)*sin(yi)*cos(zi);
319  u(1) = - exp(xi)*cos(yi)*cos(zi);
320 
321  if (x.Size() == 3)
322  {
323  u(2) = exp(xi)*sin(yi)*sin(zi);
324  }
325 }
326 
327 // Change if needed
328 double pFun_ex(const Vector & x)
329 {
330  double xi(x(0));
331  double yi(x(1));
332  double zi(0.0);
333 
334  if (x.Size() == 3)
335  {
336  zi = x(2);
337  }
338 
339  return exp(xi)*sin(yi)*cos(zi);
340 }
341 
342 void fFun(const Vector & x, Vector & f)
343 {
344  f = 0.0;
345 }
346 
347 double gFun(const Vector & x)
348 {
349  if (x.Size() == 3)
350  {
351  return -pFun_ex(x);
352  }
353  else
354  {
355  return 0;
356  }
357 }
358 
359 double f_natural(const Vector & x)
360 {
361  return (-pFun_ex(x));
362 }
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:47
int GetVSize() const
Definition: fespace.hpp:164
Class for integration rule.
Definition: intrules.hpp:83
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
int GetNumIterations() const
Definition: solvers.hpp:64
Data type for scaled Jacobi-type smoother of sparse matrix.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:235
void fFun(const Vector &x, Vector &f)
Definition: ex5.cpp:342
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
virtual void RegisterField(const char *field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
Subclass constant coefficient.
Definition: coefficient.hpp:57
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:445
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:84
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:454
virtual void Save()
Save the collection and a VisIt root file.
double RealTime()
Definition: tic_toc.cpp:317
bool iterative_mode
If true, use the second argument of Mult as an initial guess.
Definition: operator.hpp:106
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
MINRES method.
Definition: solvers.hpp:218
void ScaleRow(const int row, const double scale)
Definition: sparsemat.cpp:2025
double GetFinalNorm() const
Definition: solvers.hpp:66
Data type for Gauss-Seidel smoother of sparse matrix.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
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
void SetBlock(int i, int j, SparseMatrix *mat)
Set A(i,j) = mat.
Definition: blockmatrix.cpp:59
virtual void Mult(const Vector &b, Vector &x) const
Operator application.
Definition: solvers.cpp:1007
int dim
Definition: ex3.cpp:48
double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
IntegrationRules IntRules(0)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:292
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:231
double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:190
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
Data type sparse matrix.
Definition: sparsemat.hpp:38
double f_natural(const Vector &x)
Definition: ex5.cpp:359
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7316
Timing object.
Definition: tic_toc.hpp:34
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
Definition: solvers.hpp:61
void Assemble(int skip_zeros=1)
void Transpose(const Table &A, Table &At, int _ncols_A)
Transpose a Table.
Definition: table.cpp:385
double pFun_ex(const Vector &x)
Definition: ex5.cpp:328
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:993
int GetConverged() const
Definition: solvers.hpp:65
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator.
Definition: linearform.cpp:24
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
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
void SetAbsTol(double atol)
Definition: solvers.hpp:60
int main(int argc, char *argv[])
Definition: ex1.cpp:45
void SetRelTol(double rtol)
Definition: solvers.hpp:59
Abstract finite element space.
Definition: fespace.hpp:62
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 PartialSum()
Partial Sum.
Definition: array.cpp:121
const SparseMatrix & SpMat() const
void uFun_ex(const Vector &x, Vector &u)
Definition: ex5.cpp:308
double gFun(const Vector &x)
Definition: ex5.cpp:347
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:304
T & Last()
Return the last element in the array.
Definition: array.hpp:403
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:171
class for C-function coefficient
Vector data type.
Definition: vector.hpp:33
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 matrix.
Base class for solvers.
Definition: operator.hpp:102
Vector & GetBlock(int i)
Get the i-th vector in the block.
Integrator for (Q u, v) for VectorFiniteElements.
Definition: bilininteg.hpp:458
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:95
bool Good() const
Definition: optparser.hpp:120