MFEM  v3.3.2
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 = new Mesh(mesh_file, 1, 1);
72  int dim = mesh->Dimension();
73 
74  // 3. Refine the mesh to increase the resolution. In this example we do
75  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
76  // largest number that gives a final mesh with no more than 10,000
77  // elements.
78  {
79  int ref_levels =
80  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
81  for (int l = 0; l < ref_levels; l++)
82  {
83  mesh->UniformRefinement();
84  }
85  }
86 
87  // 4. Define a finite element space on the mesh. Here we use the
88  // Raviart-Thomas finite elements of the specified order.
89  FiniteElementCollection *hdiv_coll(new RT_FECollection(order, dim));
90  FiniteElementCollection *l2_coll(new L2_FECollection(order, dim));
91 
92  FiniteElementSpace *R_space = new FiniteElementSpace(mesh, hdiv_coll);
93  FiniteElementSpace *W_space = new FiniteElementSpace(mesh, l2_coll);
94 
95  // 5. Define the BlockStructure of the problem, i.e. define the array of
96  // offsets for each variable. The last component of the Array is the sum
97  // of the dimensions of each block.
98  Array<int> block_offsets(3); // number of variables + 1
99  block_offsets[0] = 0;
100  block_offsets[1] = R_space->GetVSize();
101  block_offsets[2] = W_space->GetVSize();
102  block_offsets.PartialSum();
103 
104  std::cout << "***********************************************************\n";
105  std::cout << "dim(R) = " << block_offsets[1] - block_offsets[0] << "\n";
106  std::cout << "dim(W) = " << block_offsets[2] - block_offsets[1] << "\n";
107  std::cout << "dim(R+W) = " << block_offsets.Last() << "\n";
108  std::cout << "***********************************************************\n";
109 
110  // 6. Define the coefficients, analytical solution, and rhs of the PDE.
111  ConstantCoefficient k(1.0);
112 
113  VectorFunctionCoefficient fcoeff(dim, fFun);
114  FunctionCoefficient fnatcoeff(f_natural);
115  FunctionCoefficient gcoeff(gFun);
116 
117  VectorFunctionCoefficient ucoeff(dim, uFun_ex);
119 
120  // 7. Allocate memory (x, rhs) for the analytical solution and the right hand
121  // side. Define the GridFunction u,p for the finite element solution and
122  // linear forms fform and gform for the right hand side. The data
123  // allocated by x and rhs are passed as a reference to the grid functions
124  // (u,p) and the linear forms (fform, gform).
125  BlockVector x(block_offsets), rhs(block_offsets);
126 
127  LinearForm *fform(new LinearForm);
128  fform->Update(R_space, rhs.GetBlock(0), 0);
131  fform->Assemble();
132 
133  LinearForm *gform(new LinearForm);
134  gform->Update(W_space, rhs.GetBlock(1), 0);
135  gform->AddDomainIntegrator(new DomainLFIntegrator(gcoeff));
136  gform->Assemble();
137 
138  // 8. Assemble the finite element matrices for the Darcy operator
139  //
140  // D = [ M B^T ]
141  // [ B 0 ]
142  // where:
143  //
144  // M = \int_\Omega k u_h \cdot v_h d\Omega u_h, v_h \in R_h
145  // B = -\int_\Omega \div u_h q_h d\Omega u_h \in R_h, q_h \in W_h
146  BilinearForm *mVarf(new BilinearForm(R_space));
147  MixedBilinearForm *bVarf(new MixedBilinearForm(R_space, W_space));
148 
150  mVarf->Assemble();
151  mVarf->Finalize();
152  SparseMatrix &M(mVarf->SpMat());
153 
155  bVarf->Assemble();
156  bVarf->Finalize();
157  SparseMatrix & B(bVarf->SpMat());
158  B *= -1.;
159  SparseMatrix *BT = Transpose(B);
160 
161  BlockMatrix darcyMatrix(block_offsets);
162  darcyMatrix.SetBlock(0,0, &M);
163  darcyMatrix.SetBlock(0,1, BT);
164  darcyMatrix.SetBlock(1,0, &B);
165 
166  // 9. Construct the operators for preconditioner
167  //
168  // P = [ diag(M) 0 ]
169  // [ 0 B diag(M)^-1 B^T ]
170  //
171  // Here we use Symmetric Gauss-Seidel to approximate the inverse of the
172  // pressure Schur Complement
173  SparseMatrix *MinvBt = Transpose(B);
174  Vector Md(M.Height());
175  M.GetDiag(Md);
176  for (int i = 0; i < Md.Size(); i++)
177  {
178  MinvBt->ScaleRow(i, 1./Md(i));
179  }
180  SparseMatrix *S = Mult(B, *MinvBt);
181 
182  Solver *invM, *invS;
183  invM = new DSmoother(M);
184 #ifndef MFEM_USE_SUITESPARSE
185  invS = new GSSmoother(*S);
186 #else
187  invS = new UMFPackSolver(*S);
188 #endif
189 
190  invM->iterative_mode = false;
191  invS->iterative_mode = false;
192 
193  BlockDiagonalPreconditioner darcyPrec(block_offsets);
194  darcyPrec.SetDiagonalBlock(0, invM);
195  darcyPrec.SetDiagonalBlock(1, invS);
196 
197  // 10. Solve the linear system with MINRES.
198  // Check the norm of the unpreconditioned residual.
199  int maxIter(1000);
200  double rtol(1.e-6);
201  double atol(1.e-10);
202 
203  chrono.Clear();
204  chrono.Start();
205  MINRESSolver solver;
206  solver.SetAbsTol(atol);
207  solver.SetRelTol(rtol);
208  solver.SetMaxIter(maxIter);
209  solver.SetOperator(darcyMatrix);
210  solver.SetPreconditioner(darcyPrec);
211  solver.SetPrintLevel(1);
212  x = 0.0;
213  solver.Mult(rhs, x);
214  chrono.Stop();
215 
216  if (solver.GetConverged())
217  std::cout << "MINRES converged in " << solver.GetNumIterations()
218  << " iterations with a residual norm of " << solver.GetFinalNorm() << ".\n";
219  else
220  std::cout << "MINRES did not converge in " << solver.GetNumIterations()
221  << " iterations. Residual norm is " << solver.GetFinalNorm() << ".\n";
222  std::cout << "MINRES solver took " << chrono.RealTime() << "s. \n";
223 
224  // 11. Create the grid functions u and p. Compute the L2 error norms.
225  GridFunction u, p;
226  u.MakeRef(R_space, x.GetBlock(0), 0);
227  p.MakeRef(W_space, x.GetBlock(1), 0);
228 
229  int order_quad = max(2, 2*order+1);
230  const IntegrationRule *irs[Geometry::NumGeom];
231  for (int i=0; i < Geometry::NumGeom; ++i)
232  {
233  irs[i] = &(IntRules.Get(i, order_quad));
234  }
235 
236  double err_u = u.ComputeL2Error(ucoeff, irs);
237  double norm_u = ComputeLpNorm(2., ucoeff, *mesh, irs);
238  double err_p = p.ComputeL2Error(pcoeff, irs);
239  double norm_p = ComputeLpNorm(2., pcoeff, *mesh, irs);
240 
241  std::cout << "|| u_h - u_ex || / || u_ex || = " << err_u / norm_u << "\n";
242  std::cout << "|| p_h - p_ex || / || p_ex || = " << err_p / norm_p << "\n";
243 
244  // 12. Save the mesh and the solution. This output can be viewed later using
245  // GLVis: "glvis -m ex5.mesh -g sol_u.gf" or "glvis -m ex5.mesh -g
246  // sol_p.gf".
247  {
248  ofstream mesh_ofs("ex5.mesh");
249  mesh_ofs.precision(8);
250  mesh->Print(mesh_ofs);
251 
252  ofstream u_ofs("sol_u.gf");
253  u_ofs.precision(8);
254  u.Save(u_ofs);
255 
256  ofstream p_ofs("sol_p.gf");
257  p_ofs.precision(8);
258  p.Save(p_ofs);
259  }
260 
261  // 13. Save data in the VisIt format
262  VisItDataCollection visit_dc("Example5", mesh);
263  visit_dc.RegisterField("velocity", &u);
264  visit_dc.RegisterField("pressure", &p);
265  visit_dc.Save();
266 
267  // 14. Send the solution by socket to a GLVis server.
268  if (visualization)
269  {
270  char vishost[] = "localhost";
271  int visport = 19916;
272  socketstream u_sock(vishost, visport);
273  u_sock.precision(8);
274  u_sock << "solution\n" << *mesh << u << "window_title 'Velocity'" << endl;
275  socketstream p_sock(vishost, visport);
276  p_sock.precision(8);
277  p_sock << "solution\n" << *mesh << p << "window_title 'Pressure'" << endl;
278  }
279 
280  // 15. Free the used memory.
281  delete fform;
282  delete gform;
283  delete invM;
284  delete invS;
285  delete S;
286  delete MinvBt;
287  delete BT;
288  delete mVarf;
289  delete bVarf;
290  delete W_space;
291  delete R_space;
292  delete l2_coll;
293  delete hdiv_coll;
294  delete mesh;
295 
296  return 0;
297 }
298 
299 
300 void uFun_ex(const Vector & x, Vector & u)
301 {
302  double xi(x(0));
303  double yi(x(1));
304  double zi(0.0);
305  if (x.Size() == 3)
306  {
307  zi = x(2);
308  }
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  {
315  u(2) = exp(xi)*sin(yi)*sin(zi);
316  }
317 }
318 
319 // Change if needed
320 double pFun_ex(const Vector & x)
321 {
322  double xi(x(0));
323  double yi(x(1));
324  double zi(0.0);
325 
326  if (x.Size() == 3)
327  {
328  zi = x(2);
329  }
330 
331  return exp(xi)*sin(yi)*cos(zi);
332 }
333 
334 void fFun(const Vector & x, Vector & f)
335 {
336  f = 0.0;
337 }
338 
339 double gFun(const Vector & x)
340 {
341  if (x.Size() == 3)
342  {
343  return -pFun_ex(x);
344  }
345  else
346  {
347  return 0;
348  }
349 }
350 
351 double f_natural(const Vector & x)
352 {
353  return (-pFun_ex(x));
354 }
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
int GetVSize() const
Definition: fespace.hpp:163
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1019
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:83
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
int GetNumIterations() const
Definition: solvers.hpp:66
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:861
void fFun(const Vector &x, Vector &f)
Definition: ex5.cpp:334
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:233
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:468
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:51
int Size() const
Returns the size of the vector.
Definition: vector.hpp:113
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:614
virtual void Save()
Save the collection and a VisIt root file.
double RealTime()
Definition: tic_toc.cpp:426
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:263
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
MINRES method.
Definition: solvers.hpp:220
void ScaleRow(const int row, const double scale)
Definition: sparsemat.cpp:2162
double GetFinalNorm() const
Definition: solvers.hpp:68
Data type for Gauss-Seidel smoother of sparse matrix.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:343
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2259
void SetBlock(int i, int j, SparseMatrix *mat)
Set A(i,j) = mat.
Definition: blockmatrix.cpp:60
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1046
int dim
Definition: ex3.cpp:47
double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:72
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:233
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:351
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6718
Timing object.
Definition: tic_toc.hpp:34
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
Definition: solvers.hpp:63
void Assemble(int skip_zeros=1)
void Transpose(const Table &A, Table &At, int _ncols_A)
Transpose a Table.
Definition: table.cpp:408
double pFun_ex(const Vector &x)
Definition: ex5.cpp:320
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:181
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1032
int GetConverged() const
Definition: solvers.hpp:67
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator.
Definition: linearform.cpp:33
int Dimension() const
Definition: mesh.hpp:641
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:239
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
int main()
void SetAbsTol(double atol)
Definition: solvers.hpp:62
void SetRelTol(double rtol)
Definition: solvers.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:139
const SparseMatrix & SpMat() const
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
void uFun_ex(const Vector &x, Vector &u)
Definition: ex5.cpp:300
double gFun(const Vector &x)
Definition: ex5.cpp:339
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:581
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:231
class for C-function coefficient
Vector data type.
Definition: vector.hpp:41
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:259
Vector & GetBlock(int i)
Get the i-th vector in the block.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:343
Integrator for (Q u, v) for VectorFiniteElements.
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:195
bool Good() const
Definition: optparser.hpp:120