MFEM  v3.0
 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(Vector & x);
39 void fFun(const Vector & x, Vector & f);
40 double gFun(Vector & x);
41 double f_natural(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  mesh->UniformRefinement();
91  }
92 
93  // 4. Define a finite element space on the mesh. Here we use the
94  // Raviart-Thomas finite elements of the specified order.
95  FiniteElementCollection *hdiv_coll(new RT_FECollection(order, dim));
96  FiniteElementCollection *l2_coll(new L2_FECollection(order, dim));
97 
98  FiniteElementSpace *R_space = new FiniteElementSpace(mesh, hdiv_coll);
99  FiniteElementSpace *W_space = new FiniteElementSpace(mesh, l2_coll);
100 
101  // 5. Define the BlockStructure of the problem, i.e. define the array of
102  // offsets for each variable. The last component of the Array is the sum
103  // of the dimensions of each block.
104  Array<int> block_offsets(3); // number of variables + 1
105  block_offsets[0] = 0;
106  block_offsets[1] = R_space->GetVSize();
107  block_offsets[2] = W_space->GetVSize();
108  block_offsets.PartialSum();
109 
110  std::cout << "***********************************************************\n";
111  std::cout << "dim(R) = " << block_offsets[1] - block_offsets[0] << "\n";
112  std::cout << "dim(W) = " << block_offsets[2] - block_offsets[1] << "\n";
113  std::cout << "dim(R+W) = " << block_offsets.Last() << "\n";
114  std::cout << "***********************************************************\n";
115 
116  // 6. Define the coefficients, analytical solution, and rhs of the PDE.
117  ConstantCoefficient k(1.0);
118 
119  VectorFunctionCoefficient fcoeff(dim, fFun);
120  FunctionCoefficient fnatcoeff(f_natural);
121  FunctionCoefficient gcoeff(gFun);
122 
123  VectorFunctionCoefficient ucoeff(dim, uFun_ex);
125 
126  // 7. Allocate memory (x, rhs) for the analytical solution and the right hand
127  // side. Define the GridFunction u,p for the finite element solution and
128  // linear forms fform and gform for the right hand side. The data
129  // allocated by x and rhs are passed as a reference to the grid fuctions
130  // (u,p) and the linear forms (fform, gform).
131  BlockVector x(block_offsets), rhs(block_offsets);
132 
133  LinearForm *fform(new LinearForm);
134  fform->Update(R_space, rhs.GetBlock(0), 0);
137  fform->Assemble();
138 
139  LinearForm *gform(new LinearForm);
140  gform->Update(W_space, rhs.GetBlock(1), 0);
141  gform->AddDomainIntegrator(new DomainLFIntegrator(gcoeff));
142  gform->Assemble();
143 
144  // 8. Assemble the finite element matrices for the Darcy operator
145  //
146  // D = [ M B^T ]
147  // [ B 0 ]
148  // where:
149  //
150  // M = \int_\Omega k u_h \cdot v_h d\Omega u_h, v_h \in R_h
151  // B = -\int_\Omega \div u_h q_h d\Omega u_h \in R_h, q_h \in W_h
152  BilinearForm *mVarf(new BilinearForm(R_space));
153  MixedBilinearForm *bVarf(new MixedBilinearForm(R_space, W_space));
154 
156  mVarf->Assemble();
157  mVarf->Finalize();
158  SparseMatrix &M(mVarf->SpMat());
159 
161  bVarf->Assemble();
162  bVarf->Finalize();
163  SparseMatrix & B(bVarf->SpMat());
164  B *= -1.;
165  SparseMatrix *BT = Transpose(B);
166 
167  BlockMatrix darcyMatrix(block_offsets);
168  darcyMatrix.SetBlock(0,0, &M);
169  darcyMatrix.SetBlock(0,1, BT);
170  darcyMatrix.SetBlock(1,0, &B);
171 
172  // 9. Construct the operators for preconditioner
173  //
174  // P = [ diag(M) 0 ]
175  // [ 0 B diag(M)^-1 B^T ]
176  //
177  // Here we use Symmetric Gauss-Seidel to approximate the inverse of the
178  // pressure Schur Complement
179  SparseMatrix *MinvBt = Transpose(B);
180  Vector Md(M.Height());
181  M.GetDiag(Md);
182  for (int i = 0; i < Md.Size(); i++)
183  MinvBt->ScaleRow(i, 1./Md(i));
184  SparseMatrix *S = Mult(B, *MinvBt);
185 
186  Solver *invM, *invS;
187  invM = new DSmoother(M);
188 #ifndef MFEM_USE_SUITESPARSE
189  invS = new GSSmoother(*S);
190 #else
191  invS = new UMFPackSolver(*S);
192 #endif
193 
194  invM->iterative_mode = false;
195  invS->iterative_mode = false;
196 
197  BlockDiagonalPreconditioner darcyPrec(block_offsets);
198  darcyPrec.SetDiagonalBlock(0, invM);
199  darcyPrec.SetDiagonalBlock(1, invS);
200 
201  // 10. Solve the linear system with MINRES.
202  // Check the norm of the unpreconditioned residual.
203  int maxIter(1000);
204  double rtol(1.e-6);
205  double atol(1.e-10);
206 
207  chrono.Clear();
208  chrono.Start();
209  MINRESSolver solver;
210  solver.SetAbsTol(atol);
211  solver.SetRelTol(rtol);
212  solver.SetMaxIter(maxIter);
213  solver.SetOperator(darcyMatrix);
214  solver.SetPreconditioner(darcyPrec);
215  solver.SetPrintLevel(1);
216  x = 0.0;
217  solver.Mult(rhs, x);
218  chrono.Stop();
219 
220  if (solver.GetConverged())
221  std::cout << "MINRES converged in " << solver.GetNumIterations()
222  << " iterations with a residual norm of " << solver.GetFinalNorm() << ".\n";
223  else
224  std::cout << "MINRES did not converge in " << solver.GetNumIterations()
225  << " iterations. Residual norm is " << solver.GetFinalNorm() << ".\n";
226  std::cout << "MINRES solver took " << chrono.RealTime() << "s. \n";
227 
228  // 11. Create the grid functions u and p. Compute the L2 error norms.
229  GridFunction u, p;
230  u.Update(R_space, x.GetBlock(0), 0);
231  p.Update(W_space, x.GetBlock(1), 0);
232 
233  int order_quad = max(2, 2*order+1);
234  const IntegrationRule *irs[Geometry::NumGeom];
235  for (int i=0; i < Geometry::NumGeom; ++i)
236  irs[i] = &(IntRules.Get(i, order_quad));
237 
238  double err_u = u.ComputeL2Error(ucoeff, irs);
239  double norm_u = ComputeLpNorm(2., ucoeff, *mesh, irs);
240  double err_p = p.ComputeL2Error(pcoeff, irs);
241  double norm_p = ComputeLpNorm(2., pcoeff, *mesh, irs);
242 
243  std::cout << "|| u_h - u_ex || / || u_ex || = " << err_u / norm_u << "\n";
244  std::cout << "|| p_h - p_ex || / || p_ex || = " << err_p / norm_p << "\n";
245 
246  // 12. Save the mesh and the solution. This output can be viewed later using
247  // GLVis: "glvis -m ex5.mesh -g sol_u.gf" or "glvis -m ex5.mesh -g
248  // sol_p.gf".
249  {
250  ofstream mesh_ofs("ex5.mesh");
251  mesh_ofs.precision(8);
252  mesh->Print(mesh_ofs);
253 
254  ofstream u_ofs("sol_u.gf");
255  u_ofs.precision(8);
256  u.Save(u_ofs);
257 
258  ofstream p_ofs("sol_p.gf");
259  p_ofs.precision(8);
260  p.Save(p_ofs);
261  }
262 
263  // 13. Save data in the VisIt format
264  VisItDataCollection visit_dc("Example5", mesh);
265  visit_dc.RegisterField("velocity", &u);
266  visit_dc.RegisterField("pressure", &p);
267  visit_dc.Save();
268 
269  // 14. Send the solution by socket to a GLVis server.
270  if (visualization)
271  {
272  char vishost[] = "localhost";
273  int visport = 19916;
274  socketstream u_sock(vishost, visport);
275  u_sock.precision(8);
276  u_sock << "solution\n" << *mesh << u << "window_title 'Velocity'" << endl;
277  socketstream p_sock(vishost, visport);
278  p_sock.precision(8);
279  p_sock << "solution\n" << *mesh << p << "window_title 'Pressure'" << endl;
280  }
281 
282  // 15. Free the used memory.
283  delete fform;
284  delete gform;
285  delete invM;
286  delete invS;
287  delete S;
288  delete MinvBt;
289  delete BT;
290  delete mVarf;
291  delete bVarf;
292  delete W_space;
293  delete R_space;
294  delete l2_coll;
295  delete hdiv_coll;
296  delete mesh;
297 
298  return 0;
299 }
300 
301 
302 void uFun_ex(const Vector & x, Vector & u)
303 {
304  double xi(x(0));
305  double yi(x(1));
306  double zi(0.0);
307  if (x.Size() == 3)
308  zi = x(2);
309 
310  u(0) = - exp(xi)*sin(yi)*cos(zi);
311  u(1) = - exp(xi)*cos(yi)*cos(zi);
312 
313  if (x.Size() == 3)
314  u(2) = exp(xi)*sin(yi)*sin(zi);
315 }
316 
317 // Change if needed
318 double pFun_ex(Vector & x)
319 {
320  double xi(x(0));
321  double yi(x(1));
322  double zi(0.0);
323 
324  if (x.Size() == 3)
325  zi = x(2);
326 
327  return exp(xi)*sin(yi)*cos(zi);
328 }
329 
330 void fFun(const Vector & x, Vector & f)
331 {
332  f = 0.0;
333 }
334 
335 double gFun(Vector & x)
336 {
337  if (x.Size() == 3)
338  return -pFun_ex(x);
339  else
340  return 0;
341 }
342 
343 double f_natural(Vector & x)
344 {
345  return (-pFun_ex(x));
346 }
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:47
int GetVSize() const
Definition: fespace.hpp:151
Class for integration rule.
Definition: intrules.hpp:63
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:26
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:194
void fFun(const Vector &x, Vector &f)
Definition: ex5.cpp:330
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:351
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
virtual void Save()
Save the collection and a VisIt root file.
double RealTime()
Definition: tic_toc.cpp:309
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:1760
Data type for Gauss-Seidel smoother of sparse matrix.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
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
void SetBlock(int i, int j, SparseMatrix *mat)
Set A(i,j) = mat.
Definition: blockmatrix.cpp:54
virtual void Mult(const Vector &b, Vector &x) const
Operator application.
Definition: solvers.cpp:946
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:264
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:166
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
Data type sparse matrix.
Definition: sparsemat.hpp:38
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6225
double f_natural(Vector &x)
Definition: ex5.cpp:343
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:305
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:934
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator.
Definition: linearform.cpp:24
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
void SetAbsTol(double atol)
Definition: solvers.hpp:60
int main(int argc, char *argv[])
Definition: ex1.cpp:39
void SetRelTol(double rtol)
Definition: solvers.hpp:59
double gFun(Vector &x)
Definition: ex5.cpp:335
Abstract finite element space.
Definition: fespace.hpp:61
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:120
const SparseMatrix & SpMat() const
void uFun_ex(const Vector &x, Vector &u)
Definition: ex5.cpp:302
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
T & Last()
Return the last element in the array.
Definition: array.hpp:361
double pFun_ex(Vector &x)
Definition: ex5.cpp:318
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:171
class for C-function coefficient
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.
Base class for solvers.
Definition: operator.hpp:102
double GetFinalNorm()
Definition: solvers.hpp:66
Vector & GetBlock(int i)
Get the i-th vector in the block.
Integrator for (Q u, v) for VectorFiniteElements.
Definition: bilininteg.hpp:434
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:83
bool Good() const
Definition: optparser.hpp:120