MFEM  v3.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex5p.cpp
Go to the documentation of this file.
1 // MFEM Example 5 - Parallel Version
2 //
3 // Compile with: make ex5p
4 //
5 // Sample runs: mpirun -np 4 ex5p -m ../data/square-disc.mesh
6 // mpirun -np 4 ex5p -m ../data/star.mesh
7 // mpirun -np 4 ex5p -m ../data/beam-tet.mesh
8 // mpirun -np 4 ex5p -m ../data/beam-hex.mesh
9 // mpirun -np 4 ex5p -m ../data/escher.mesh
10 // mpirun -np 4 ex5p -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 
32 using namespace std;
33 using namespace mfem;
34 
35 // Define the analytical solution and forcing terms / boundary conditions
36 void uFun_ex(const Vector & x, Vector & u);
37 double pFun_ex(const Vector & x);
38 void fFun(const Vector & x, Vector & f);
39 double gFun(const Vector & x);
40 double f_natural(const Vector & x);
41 
42 int main(int argc, char *argv[])
43 {
44  StopWatch chrono;
45 
46  // 1. Initialize MPI.
47  int num_procs, myid;
48  MPI_Init(&argc, &argv);
49  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
50  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
51  bool verbose = (myid == 0);
52 
53  // 2. Parse command-line options.
54  const char *mesh_file = "../data/star.mesh";
55  int order = 1;
56  bool visualization = 1;
57 
58  OptionsParser args(argc, argv);
59  args.AddOption(&mesh_file, "-m", "--mesh",
60  "Mesh file to use.");
61  args.AddOption(&order, "-o", "--order",
62  "Finite element order (polynomial degree).");
63  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
64  "--no-visualization",
65  "Enable or disable GLVis visualization.");
66  args.Parse();
67  if (!args.Good())
68  {
69  if (verbose)
70  {
71  args.PrintUsage(cout);
72  }
73  MPI_Finalize();
74  return 1;
75  }
76  if (verbose)
77  {
78  args.PrintOptions(cout);
79  }
80 
81  // 3. Read the (serial) mesh from the given mesh file on all processors. We
82  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
83  // and volume meshes with the same code.
84  Mesh *mesh;
85  ifstream imesh(mesh_file);
86  if (!imesh)
87  {
88  if (verbose)
89  {
90  cerr << "\nCan not open mesh file: " << mesh_file << '\n' << endl;
91  }
92  MPI_Finalize();
93  return 2;
94  }
95  mesh = new Mesh(imesh, 1, 1);
96  imesh.close();
97  int dim = mesh->Dimension();
98 
99  // 4. Refine the serial mesh on all processors to increase the resolution. In
100  // this example we do 'ref_levels' of uniform refinement. We choose
101  // 'ref_levels' to be the largest number that gives a final mesh with no
102  // more than 10,000 elements.
103  {
104  int ref_levels =
105  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
106  for (int l = 0; l < ref_levels; l++)
107  {
108  mesh->UniformRefinement();
109  }
110  }
111 
112  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
113  // this mesh further in parallel to increase the resolution. Once the
114  // parallel mesh is defined, the serial mesh can be deleted.
115  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
116  delete mesh;
117  {
118  int par_ref_levels = 2;
119  for (int l = 0; l < par_ref_levels; l++)
120  {
121  pmesh->UniformRefinement();
122  }
123  }
124 
125  // 6. Define a parallel finite element space on the parallel mesh. Here we
126  // use the Raviart-Thomas finite elements of the specified order.
127  FiniteElementCollection *hdiv_coll(new RT_FECollection(order, dim));
128  FiniteElementCollection *l2_coll(new L2_FECollection(order, dim));
129 
130  ParFiniteElementSpace *R_space = new ParFiniteElementSpace(pmesh, hdiv_coll);
131  ParFiniteElementSpace *W_space = new ParFiniteElementSpace(pmesh, l2_coll);
132 
133  HYPRE_Int dimR = R_space->GlobalTrueVSize();
134  HYPRE_Int dimW = W_space->GlobalTrueVSize();
135 
136  if (verbose)
137  {
138  std::cout << "***********************************************************\n";
139  std::cout << "dim(R) = " << dimR << "\n";
140  std::cout << "dim(W) = " << dimW << "\n";
141  std::cout << "dim(R+W) = " << dimR + dimW << "\n";
142  std::cout << "***********************************************************\n";
143  }
144 
145  // 7. Define the two BlockStructure of the problem. block_offsets is used
146  // for Vector based on dof (like ParGridFunction or ParLinearForm),
147  // block_trueOffstes is used for Vector based on trueDof (HypreParVector
148  // for the rhs and solution of the linear system). The offsets computed
149  // here are local to the processor.
150  Array<int> block_offsets(3); // number of variables + 1
151  block_offsets[0] = 0;
152  block_offsets[1] = R_space->GetVSize();
153  block_offsets[2] = W_space->GetVSize();
154  block_offsets.PartialSum();
155 
156  Array<int> block_trueOffsets(3); // number of variables + 1
157  block_trueOffsets[0] = 0;
158  block_trueOffsets[1] = R_space->TrueVSize();
159  block_trueOffsets[2] = W_space->TrueVSize();
160  block_trueOffsets.PartialSum();
161 
162  // 8. Define the coefficients, analytical solution, and rhs of the PDE.
163  ConstantCoefficient k(1.0);
164 
165  VectorFunctionCoefficient fcoeff(dim, fFun);
166  FunctionCoefficient fnatcoeff(f_natural);
167  FunctionCoefficient gcoeff(gFun);
168 
169  VectorFunctionCoefficient ucoeff(dim, uFun_ex);
171 
172  // 9. Define the parallel grid function and parallel linear forms, solution
173  // vector and rhs.
174  BlockVector x(block_offsets), rhs(block_offsets);
175  BlockVector trueX(block_trueOffsets), trueRhs(block_trueOffsets);
176 
177  ParLinearForm *fform(new ParLinearForm);
178  fform->Update(R_space, rhs.GetBlock(0), 0);
181  fform->Assemble();
182  fform->ParallelAssemble(trueRhs.GetBlock(0));
183 
184  ParLinearForm *gform(new ParLinearForm);
185  gform->Update(W_space, rhs.GetBlock(1), 0);
186  gform->AddDomainIntegrator(new DomainLFIntegrator(gcoeff));
187  gform->Assemble();
188  gform->ParallelAssemble(trueRhs.GetBlock(1));
189 
190  // 10. Assemble the finite element matrices for the Darcy operator
191  //
192  // D = [ M B^T ]
193  // [ B 0 ]
194  // where:
195  //
196  // M = \int_\Omega k u_h \cdot v_h d\Omega u_h, v_h \in R_h
197  // B = -\int_\Omega \div u_h q_h d\Omega u_h \in R_h, q_h \in W_h
198  ParBilinearForm *mVarf(new ParBilinearForm(R_space));
199  ParMixedBilinearForm *bVarf(new ParMixedBilinearForm(R_space, W_space));
200 
201  HypreParMatrix *M, *B;
202 
204  mVarf->Assemble();
205  mVarf->Finalize();
206  M = mVarf->ParallelAssemble();
207 
209  bVarf->Assemble();
210  bVarf->Finalize();
211  B = bVarf->ParallelAssemble();
212  (*B) *= -1;
213 
214  HypreParMatrix *BT = B->Transpose();
215 
216  BlockOperator *darcyOp = new BlockOperator(block_trueOffsets);
217  darcyOp->SetBlock(0,0, M);
218  darcyOp->SetBlock(0,1, BT);
219  darcyOp->SetBlock(1,0, B);
220 
221  // 11. Construct the operators for preconditioner
222  //
223  // P = [ diag(M) 0 ]
224  // [ 0 B diag(M)^-1 B^T ]
225  //
226  // Here we use Symmetric Gauss-Seidel to approximate the inverse of the
227  // pressure Schur Complement.
228  HypreParMatrix *MinvBt = B->Transpose();
229  HypreParVector *Md = new HypreParVector(MPI_COMM_WORLD, M->GetGlobalNumRows(),
230  M->GetRowStarts());
231  M->GetDiag(*Md);
232 
233  MinvBt->InvScaleRows(*Md);
234  HypreParMatrix *S = ParMult(B, MinvBt);
235 
236  HypreSolver *invM, *invS;
237  invM = new HypreDiagScale(*M);
238  invS = new HypreBoomerAMG(*S);
239 
240  invM->iterative_mode = false;
241  invS->iterative_mode = false;
242 
244  block_trueOffsets);
245  darcyPr->SetDiagonalBlock(0, invM);
246  darcyPr->SetDiagonalBlock(1, invS);
247 
248  // 12. Solve the linear system with MINRES.
249  // Check the norm of the unpreconditioned residual.
250 
251  int maxIter(500);
252  double rtol(1.e-6);
253  double atol(1.e-10);
254 
255  chrono.Clear();
256  chrono.Start();
257  MINRESSolver solver(MPI_COMM_WORLD);
258  solver.SetAbsTol(atol);
259  solver.SetRelTol(rtol);
260  solver.SetMaxIter(maxIter);
261  solver.SetOperator(*darcyOp);
262  solver.SetPreconditioner(*darcyPr);
263  solver.SetPrintLevel(verbose);
264  trueX = 0.0;
265  solver.Mult(trueRhs, trueX);
266  chrono.Stop();
267 
268  if (verbose)
269  {
270  if (solver.GetConverged())
271  std::cout << "MINRES converged in " << solver.GetNumIterations()
272  << " iterations with a residual norm of " << solver.GetFinalNorm() << ".\n";
273  else
274  std::cout << "MINRES did not converge in " << solver.GetNumIterations()
275  << " iterations. Residual norm is " << solver.GetFinalNorm() << ".\n";
276  std::cout << "MINRES solver took " << chrono.RealTime() << "s. \n";
277  }
278 
279  // 13. Extract the parallel grid function corresponding to the finite element
280  // approximation X. This is the local solution on each processor. Compute
281  // L2 error norms.
284  u->Update(R_space, x.GetBlock(0), 0);
285  p->Update(W_space, x.GetBlock(1), 0);
286  u->Distribute(&(trueX.GetBlock(0)));
287  p->Distribute(&(trueX.GetBlock(1)));
288 
289  int order_quad = max(2, 2*order+1);
290  const IntegrationRule *irs[Geometry::NumGeom];
291  for (int i=0; i < Geometry::NumGeom; ++i)
292  {
293  irs[i] = &(IntRules.Get(i, order_quad));
294  }
295 
296  double err_u = u->ComputeL2Error(ucoeff, irs);
297  double norm_u = ComputeGlobalLpNorm(2, ucoeff, *pmesh, irs);
298  double err_p = p->ComputeL2Error(pcoeff, irs);
299  double norm_p = ComputeGlobalLpNorm(2, pcoeff, *pmesh, irs);
300 
301  if (verbose)
302  {
303  std::cout << "|| u_h - u_ex || / || u_ex || = " << err_u / norm_u << "\n";
304  std::cout << "|| p_h - p_ex || / || p_ex || = " << err_p / norm_p << "\n";
305  }
306 
307  // 14. Save the refined mesh and the solution in parallel. This output can be
308  // viewed later using GLVis: "glvis -np <np> -m mesh -g sol_*".
309  {
310  ostringstream mesh_name, u_name, p_name;
311  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
312  u_name << "sol_u." << setfill('0') << setw(6) << myid;
313  p_name << "sol_p." << setfill('0') << setw(6) << myid;
314 
315  ofstream mesh_ofs(mesh_name.str().c_str());
316  mesh_ofs.precision(8);
317  pmesh->Print(mesh_ofs);
318 
319  ofstream u_ofs(u_name.str().c_str());
320  u_ofs.precision(8);
321  u->Save(u_ofs);
322 
323  ofstream p_ofs(p_name.str().c_str());
324  p_ofs.precision(8);
325  p->Save(p_ofs);
326  }
327 
328  // 15. Save data in the VisIt format
329  VisItDataCollection visit_dc("Example5-Parallel", pmesh);
330  visit_dc.RegisterField("velocity", u);
331  visit_dc.RegisterField("pressure", p);
332  visit_dc.Save();
333 
334  // 16. Send the solution by socket to a GLVis server.
335  if (visualization)
336  {
337  char vishost[] = "localhost";
338  int visport = 19916;
339  socketstream u_sock(vishost, visport);
340  u_sock << "parallel " << num_procs << " " << myid << "\n";
341  u_sock.precision(8);
342  u_sock << "solution\n" << *pmesh << *u << "window_title 'Velocity'"
343  << endl;
344  // Make sure all ranks have sent their 'u' solution before initiating
345  // another set of GLVis connections (one from each rank):
346  MPI_Barrier(pmesh->GetComm());
347  socketstream p_sock(vishost, visport);
348  p_sock << "parallel " << num_procs << " " << myid << "\n";
349  p_sock.precision(8);
350  p_sock << "solution\n" << *pmesh << *p << "window_title 'Pressure'"
351  << endl;
352  }
353 
354  // 17. Free the used memory.
355  delete fform;
356  delete gform;
357  delete u;
358  delete p;
359  delete darcyOp;
360  delete darcyPr;
361  delete invM;
362  delete invS;
363  delete S;
364  delete Md;
365  delete MinvBt;
366  delete BT;
367  delete B;
368  delete M;
369  delete mVarf;
370  delete bVarf;
371  delete W_space;
372  delete R_space;
373  delete l2_coll;
374  delete hdiv_coll;
375  delete pmesh;
376 
377  MPI_Finalize();
378 
379  return 0;
380 }
381 
382 
383 void uFun_ex(const Vector & x, Vector & u)
384 {
385  double xi(x(0));
386  double yi(x(1));
387  double zi(0.0);
388  if (x.Size() == 3)
389  {
390  zi = x(2);
391  }
392 
393  u(0) = - exp(xi)*sin(yi)*cos(zi);
394  u(1) = - exp(xi)*cos(yi)*cos(zi);
395 
396  if (x.Size() == 3)
397  {
398  u(2) = exp(xi)*sin(yi)*sin(zi);
399  }
400 }
401 
402 // Change if needed
403 double pFun_ex(const Vector & x)
404 {
405  double xi(x(0));
406  double yi(x(1));
407  double zi(0.0);
408 
409  if (x.Size() == 3)
410  {
411  zi = x(2);
412  }
413 
414  return exp(xi)*sin(yi)*cos(zi);
415 }
416 
417 void fFun(const Vector & x, Vector & f)
418 {
419  f = 0.0;
420 }
421 
422 double gFun(const Vector & x)
423 {
424  if (x.Size() == 3)
425  {
426  return -pFun_ex(x);
427  }
428  else
429  {
430  return 0;
431  }
432 }
433 
434 double f_natural(const Vector & x)
435 {
436  return (-pFun_ex(x));
437 }
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
int GetNumIterations() const
Definition: solvers.hpp:64
double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:148
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
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 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.
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:306
double RealTime()
Definition: tic_toc.cpp:317
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetBlock(int iRow, int iCol, Operator *op)
Add a block op in the block-entry (iblock, jblock).
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 Update(ParFiniteElementSpace *pf=NULL)
Definition: plinearform.cpp:21
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
double GetFinalNorm() const
Definition: solvers.hpp:66
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
The BoomerAMG solver in hypre.
Definition: hypre.hpp:698
virtual void Mult(const Vector &b, Vector &x) const
Operator application.
Definition: solvers.cpp:1007
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_Int GetGlobalNumRows() const
Definition: hypre.hpp:333
int dim
Definition: ex3.cpp:48
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
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
Jacobi preconditioner in hypre.
Definition: hypre.hpp:659
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 Assemble(int skip_zeros=1)
Assemble the local matrix.
double pFun_ex(const Vector &x)
Definition: ex5.cpp:328
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:828
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
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:58
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
MPI_Comm GetComm() const
Definition: pmesh.hpp:86
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
double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:79
HypreParMatrix * Transpose()
Returns the transpose of *this.
Definition: hypre.cpp:892
void uFun_ex(const Vector &x, Vector &u)
Definition: ex5.cpp:308
double gFun(const Vector &x)
Definition: ex5.cpp:347
Class for parallel bilinear form.
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
Definition: hypre.cpp:1081
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:522
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:171
class for C-function coefficient
Vector data type.
Definition: vector.hpp:33
HYPRE_Int * GetRowStarts() const
Definition: hypre.hpp:339
Class for parallel grid function.
Definition: pgridfunc.hpp:31
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:143
A class to handle Block systems in a matrix-free implementation.
Vector & GetBlock(int i)
Get the i-th vector in the block.
Class for parallel meshes.
Definition: pmesh.hpp:28
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:34
Integrator for (Q u, v) for VectorFiniteElements.
Definition: bilininteg.hpp:458
HypreParMatrix * ParMult(HypreParMatrix *A, HypreParMatrix *B)
Returns the matrix A * B.
Definition: hypre.cpp:1308
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
virtual void Print(std::ostream &out=std::cout) const
Definition: pmesh.cpp:2850
bool Good() const
Definition: optparser.hpp:120