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