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