MFEM  v4.1.0
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 // PETSc Modification
3 //
4 // Compile with: make ex5p
5 //
6 // Sample runs:
7 // mpirun -np 4 ex5p -m ../../data/beam-tet.mesh --petscopts rc_ex5p_fieldsplit
8 // mpirun -np 4 ex5p -m ../../data/star.mesh --petscopts rc_ex5p_bddc --nonoverlapping
9 //
10 // Description: This example code solves a simple 2D/3D mixed Darcy problem
11 // corresponding to the saddle point system
12 // k*u + grad p = f
13 // - div u = g
14 // with natural boundary condition -p = <given pressure>.
15 // Here, we use a given exact solution (u,p) and compute the
16 // corresponding r.h.s. (f,g). We discretize with Raviart-Thomas
17 // finite elements (velocity u) and piecewise discontinuous
18 // polynomials (pressure p).
19 //
20 // The example demonstrates the use of the BlockMatrix class, as
21 // well as the collective saving of several grid functions in a
22 // VisIt (visit.llnl.gov) visualization format.
23 //
24 // Two types of PETSc solvers can be used: BDDC or fieldsplit.
25 // When using BDDC, the nonoverlapping assembly feature should be
26 // used. This specific example needs PETSc compiled with support
27 // for SuiteSparse and/or MUMPS for using BDDC.
28 //
29 // We recommend viewing examples 1-4 before viewing this example.
30 
31 #include "mfem.hpp"
32 #include <fstream>
33 #include <iostream>
34 
35 #ifndef MFEM_USE_PETSC
36 #error This example requires that MFEM is built with MFEM_USE_PETSC=YES
37 #endif
38 
39 using namespace std;
40 using namespace mfem;
41 
42 // Define the analytical solution and forcing terms / boundary conditions
43 void uFun_ex(const Vector & x, Vector & u);
44 double pFun_ex(const Vector & x);
45 void fFun(const Vector & x, Vector & f);
46 double gFun(const Vector & x);
47 double f_natural(const Vector & x);
48 
49 int main(int argc, char *argv[])
50 {
51  StopWatch chrono;
52 
53  // 1. Initialize MPI.
54  int num_procs, myid;
55  MPI_Init(&argc, &argv);
56  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
57  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
58  bool verbose = (myid == 0);
59 
60  // 2. Parse command-line options.
61  const char *mesh_file = "../../data/star.mesh";
62  int order = 1;
63  bool par_format = false;
64  bool visualization = 1;
65  bool use_petsc = true;
66  bool use_nonoverlapping = false;
67  bool local_bdr_spec = false;
68  const char *petscrc_file = "";
69 
70  OptionsParser args(argc, argv);
71  args.AddOption(&mesh_file, "-m", "--mesh",
72  "Mesh file to use.");
73  args.AddOption(&order, "-o", "--order",
74  "Finite element order (polynomial degree).");
75  args.AddOption(&par_format, "-pf", "--parallel-format", "-sf",
76  "--serial-format",
77  "Format to use when saving the results for VisIt.");
78  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
79  "--no-visualization",
80  "Enable or disable GLVis visualization.");
81  args.AddOption(&use_petsc, "-usepetsc", "--usepetsc", "-no-petsc",
82  "--no-petsc",
83  "Use or not PETSc to solve the linear system.");
84  args.AddOption(&petscrc_file, "-petscopts", "--petscopts",
85  "PetscOptions file to use.");
86  args.AddOption(&use_nonoverlapping, "-nonoverlapping", "--nonoverlapping",
87  "-no-nonoverlapping", "--no-nonoverlapping",
88  "Use or not the block diagonal PETSc's matrix format "
89  "for non-overlapping domain decomposition.");
90  args.AddOption(&local_bdr_spec, "-local-bdr", "--local-bdr", "-no-local-bdr",
91  "--no-local-bdr",
92  "Specify boundary dofs in local (Vdofs) ordering.");
93  args.Parse();
94  if (!args.Good())
95  {
96  if (verbose)
97  {
98  args.PrintUsage(cout);
99  }
100  MPI_Finalize();
101  return 1;
102  }
103  if (verbose)
104  {
105  args.PrintOptions(cout);
106  }
107  // 2b. We initialize PETSc
108  if (use_petsc) { MFEMInitializePetsc(NULL,NULL,petscrc_file,NULL); }
109 
110  // 3. Read the (serial) mesh from the given mesh file on all processors. We
111  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
112  // and volume meshes with the same code.
113  Mesh *mesh = new Mesh(mesh_file, 1, 1);
114  int dim = mesh->Dimension();
115 
116  // 4. Refine the serial mesh on all processors to increase the resolution. In
117  // this example we do 'ref_levels' of uniform refinement. We choose
118  // 'ref_levels' to be the largest number that gives a final mesh with no
119  // more than 10,000 elements.
120  {
121  int ref_levels =
122  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
123  for (int l = 0; l < ref_levels; l++)
124  {
125  mesh->UniformRefinement();
126  }
127  }
128 
129  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
130  // this mesh further in parallel to increase the resolution. Once the
131  // parallel mesh is defined, the serial mesh can be deleted.
132  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
133  delete mesh;
134  {
135  int par_ref_levels = 2;
136  for (int l = 0; l < par_ref_levels; l++)
137  {
138  pmesh->UniformRefinement();
139  }
140  }
141 
142  // 6. Define a parallel finite element space on the parallel mesh. Here we
143  // use the Raviart-Thomas finite elements of the specified order.
144  FiniteElementCollection *hdiv_coll(new RT_FECollection(order, dim));
145  FiniteElementCollection *l2_coll(new L2_FECollection(order, dim));
146 
147  ParFiniteElementSpace *R_space = new ParFiniteElementSpace(pmesh, hdiv_coll);
148  ParFiniteElementSpace *W_space = new ParFiniteElementSpace(pmesh, l2_coll);
149 
150  HYPRE_Int dimR = R_space->GlobalTrueVSize();
151  HYPRE_Int dimW = W_space->GlobalTrueVSize();
152 
153  if (verbose)
154  {
155  std::cout << "***********************************************************\n";
156  std::cout << "dim(R) = " << dimR << "\n";
157  std::cout << "dim(W) = " << dimW << "\n";
158  std::cout << "dim(R+W) = " << dimR + dimW << "\n";
159  std::cout << "***********************************************************\n";
160  }
161 
162  // 7. Define the two BlockStructure of the problem. block_offsets is used
163  // for Vector based on dof (like ParGridFunction or ParLinearForm),
164  // block_trueOffstes is used for Vector based on trueDof (HypreParVector
165  // for the rhs and solution of the linear system). The offsets computed
166  // here are local to the processor.
167  Array<int> block_offsets(3); // number of variables + 1
168  block_offsets[0] = 0;
169  block_offsets[1] = R_space->GetVSize();
170  block_offsets[2] = W_space->GetVSize();
171  block_offsets.PartialSum();
172 
173  Array<int> block_trueOffsets(3); // number of variables + 1
174  block_trueOffsets[0] = 0;
175  block_trueOffsets[1] = R_space->TrueVSize();
176  block_trueOffsets[2] = W_space->TrueVSize();
177  block_trueOffsets.PartialSum();
178 
179  // 8. Define the coefficients, analytical solution, and rhs of the PDE.
180  ConstantCoefficient k(1.0);
181 
182  VectorFunctionCoefficient fcoeff(dim, fFun);
183  FunctionCoefficient fnatcoeff(f_natural);
184  FunctionCoefficient gcoeff(gFun);
185 
186  VectorFunctionCoefficient ucoeff(dim, uFun_ex);
188 
189  // 9. Define the parallel grid function and parallel linear forms, solution
190  // vector and rhs.
191  BlockVector x(block_offsets), rhs(block_offsets);
192  BlockVector trueX(block_trueOffsets), trueRhs(block_trueOffsets);
193 
194  ParLinearForm *fform(new ParLinearForm);
195  fform->Update(R_space, rhs.GetBlock(0), 0);
198  fform->Assemble();
199  fform->ParallelAssemble(trueRhs.GetBlock(0));
200 
201  ParLinearForm *gform(new ParLinearForm);
202  gform->Update(W_space, rhs.GetBlock(1), 0);
203  gform->AddDomainIntegrator(new DomainLFIntegrator(gcoeff));
204  gform->Assemble();
205  gform->ParallelAssemble(trueRhs.GetBlock(1));
206 
207  // 10. Assemble the finite element matrices for the Darcy operator
208  //
209  // D = [ M B^T ]
210  // [ B 0 ]
211  // where:
212  //
213  // M = \int_\Omega k u_h \cdot v_h d\Omega u_h, v_h \in R_h
214  // B = -\int_\Omega \div u_h q_h d\Omega u_h \in R_h, q_h \in W_h
215  ParBilinearForm *mVarf(new ParBilinearForm(R_space));
216  ParMixedBilinearForm *bVarf(new ParMixedBilinearForm(R_space, W_space));
217 
218  PetscParMatrix *pM = NULL, *pB = NULL, *pBT = NULL;
219  HypreParMatrix *M = NULL, *B = NULL, *BT = NULL;
220  Operator::Type tid =
221  !use_petsc ? Operator::Hypre_ParCSR :
222  (use_nonoverlapping ? Operator::PETSC_MATIS : Operator::PETSC_MATAIJ);
223  OperatorHandle Mh(tid), Bh(tid);
224 
226  mVarf->Assemble();
227  mVarf->Finalize();
228  mVarf->ParallelAssemble(Mh);
229  if (!use_petsc) { Mh.Get(M); }
230  else { Mh.Get(pM); }
231  Mh.SetOperatorOwner(false);
232 
234  bVarf->Assemble();
235  bVarf->Finalize();
236  if (!use_petsc)
237  {
238  B = bVarf->ParallelAssemble();
239  (*B) *= -1;
240  }
241  else
242  {
243  bVarf->ParallelAssemble(Bh);
244  Bh.Get(pB);
245  Bh.SetOperatorOwner(false);
246  (*pB) *= -1;
247  }
248 
249  if (!use_petsc) { BT = B->Transpose(); }
250  else { pBT = pB->Transpose(); };
251 
252  Operator *darcyOp = NULL;
253  if (!use_petsc)
254  {
255  BlockOperator *tdarcyOp = new BlockOperator(block_trueOffsets);
256  tdarcyOp->SetBlock(0,0,M);
257  tdarcyOp->SetBlock(0,1,BT);
258  tdarcyOp->SetBlock(1,0,B);
259  darcyOp = tdarcyOp;
260  }
261  else
262  {
263  // We construct the BlockOperator and we then convert it to a
264  // PetscParMatrix to avoid any conversion in the construction of the
265  // preconditioners.
266  BlockOperator *tdarcyOp = new BlockOperator(block_trueOffsets);
267  tdarcyOp->SetBlock(0,0,pM);
268  tdarcyOp->SetBlock(0,1,pBT);
269  tdarcyOp->SetBlock(1,0,pB);
270  darcyOp = new PetscParMatrix(pM->GetComm(),tdarcyOp,
271  use_nonoverlapping ? Operator::PETSC_MATIS :
272  Operator::PETSC_MATAIJ);
273  delete tdarcyOp;
274  }
275 
276  // 11. Construct the operators for preconditioner
277  //
278  // P = [ diag(M) 0 ]
279  // [ 0 B diag(M)^-1 B^T ]
280  //
281  // Here we use Symmetric Gauss-Seidel to approximate the inverse of the
282  // pressure Schur Complement.
283  PetscPreconditioner *pdarcyPr = NULL;
284  BlockDiagonalPreconditioner *darcyPr = NULL;
285  HypreSolver *invM = NULL, *invS = NULL;
286  HypreParMatrix *S = NULL;
287  HypreParMatrix *MinvBt = NULL;
288  HypreParVector *Md = NULL;
289  if (!use_petsc)
290  {
291  MinvBt = B->Transpose();
292  Md = new HypreParVector(MPI_COMM_WORLD, M->GetGlobalNumRows(),
293  M->GetRowStarts());
294  M->GetDiag(*Md);
295 
296  MinvBt->InvScaleRows(*Md);
297  S = ParMult(B, MinvBt);
298 
299  invM = new HypreDiagScale(*M);
300  invS = new HypreBoomerAMG(*S);
301 
302  invM->iterative_mode = false;
303  invS->iterative_mode = false;
304 
305  darcyPr = new BlockDiagonalPreconditioner(block_trueOffsets);
306  darcyPr->SetDiagonalBlock(0, invM);
307  darcyPr->SetDiagonalBlock(1, invS);
308  }
309  else
310  {
311  if (use_nonoverlapping)
312  {
314 
315  // For saddle point problems, we need to provide BDDC the list of
316  // boundary dofs either essential or natural.
317  // Since R_space is the only space that may have boundary dofs and it
318  // is ordered first then W_space, we don't need any local offset when
319  // specifying the dofs.
320  Array<int> bdr_tdof_list;
321  if (pmesh->bdr_attributes.Size())
322  {
323  Array<int> bdr(pmesh->bdr_attributes.Max());
324  bdr = 1;
325 
326  if (!local_bdr_spec)
327  {
328  // Essential dofs in global ordering
329  R_space->GetEssentialTrueDofs(bdr, bdr_tdof_list);
330  }
331  else
332  {
333  // Alternatively, you can also provide the list of dofs in local
334  // ordering
335  R_space->GetEssentialVDofs(bdr, bdr_tdof_list);
336  bdr_tdof_list.SetSize(R_space->GetVSize()+W_space->GetVSize(),0);
337  }
338  opts.SetNatBdrDofs(&bdr_tdof_list,local_bdr_spec);
339  }
340  else
341  {
342  MFEM_WARNING("Missing boundary dofs. This may cause solver failures.");
343  }
344 
345  // See also command line options rc_ex5p_bddc
346  pdarcyPr = new PetscBDDCSolver(MPI_COMM_WORLD,*darcyOp,opts,"prec_");
347  }
348  else
349  {
350  // With PETSc, we can construct the (same) block-diagonal solver with
351  // command line options (see rc_ex5p_fieldsplit)
352  pdarcyPr = new PetscFieldSplitSolver(MPI_COMM_WORLD,*darcyOp,"prec_");
353  }
354  }
355 
356  // 12. Solve the linear system with MINRES.
357  // Check the norm of the unpreconditioned residual.
358 
359  int maxIter(500);
360  double rtol(1.e-6);
361  double atol(1.e-10);
362 
363  chrono.Clear();
364  chrono.Start();
365 
366  trueX = 0.0;
367  if (!use_petsc)
368  {
369  MINRESSolver solver(MPI_COMM_WORLD);
370  solver.SetAbsTol(atol);
371  solver.SetRelTol(rtol);
372  solver.SetMaxIter(maxIter);
373  solver.SetOperator(*darcyOp);
374  solver.SetPreconditioner(*darcyPr);
375  solver.SetPrintLevel(1);
376  solver.Mult(trueRhs, trueX);
377  if (verbose)
378  {
379  if (solver.GetConverged())
380  {
381  std::cout << "MINRES converged in " << solver.GetNumIterations()
382  << " iterations with a residual norm of "
383  << solver.GetFinalNorm() << ".\n";
384  }
385  else
386  {
387  std::cout << "MINRES did not converge in "
388  << solver.GetNumIterations()
389  << " iterations. Residual norm is "
390  << solver.GetFinalNorm() << ".\n";
391  }
392  std::cout << "MINRES solver took " << chrono.RealTime() << "s. \n";
393  }
394 
395  }
396  else
397  {
398  std::string solvertype;
399  PetscLinearSolver *solver;
400  if (use_nonoverlapping)
401  {
402  // We can use conjugate gradients to solve the problem
403  solver = new PetscPCGSolver(MPI_COMM_WORLD);
404  solvertype = "PCG";
405  }
406  else
407  {
408  solver = new PetscLinearSolver(MPI_COMM_WORLD);
409  solvertype = "MINRES";
410  }
411  solver->SetOperator(*darcyOp);
412  solver->SetPreconditioner(*pdarcyPr);
413  solver->SetAbsTol(atol);
414  solver->SetRelTol(rtol);
415  solver->SetMaxIter(maxIter);
416  solver->SetPrintLevel(2);
417  solver->Mult(trueRhs, trueX);
418  if (verbose)
419  {
420  if (solver->GetConverged())
421  {
422  std::cout << solvertype << " converged in "
423  << solver->GetNumIterations()
424  << " iterations with a residual norm of "
425  << solver->GetFinalNorm() << ".\n";
426  }
427  else
428  {
429  std::cout << solvertype << " did not converge in "
430  << solver->GetNumIterations()
431  << " iterations. Residual norm is "
432  << solver->GetFinalNorm() << ".\n";
433  }
434  std::cout << solvertype << " solver took "
435  << chrono.RealTime() << "s. \n";
436  }
437  delete solver;
438  }
439  chrono.Stop();
440 
441  // 13. Extract the parallel grid function corresponding to the finite element
442  // approximation X. This is the local solution on each processor. Compute
443  // L2 error norms.
446  u->MakeRef(R_space, x.GetBlock(0), 0);
447  p->MakeRef(W_space, x.GetBlock(1), 0);
448  u->Distribute(&(trueX.GetBlock(0)));
449  p->Distribute(&(trueX.GetBlock(1)));
450 
451  int order_quad = max(2, 2*order+1);
452  const IntegrationRule *irs[Geometry::NumGeom];
453  for (int i=0; i < Geometry::NumGeom; ++i)
454  {
455  irs[i] = &(IntRules.Get(i, order_quad));
456  }
457 
458  double err_u = u->ComputeL2Error(ucoeff, irs);
459  double norm_u = ComputeGlobalLpNorm(2, ucoeff, *pmesh, irs);
460  double err_p = p->ComputeL2Error(pcoeff, irs);
461  double norm_p = ComputeGlobalLpNorm(2, pcoeff, *pmesh, irs);
462 
463  if (verbose)
464  {
465  std::cout << "|| u_h - u_ex || / || u_ex || = " << err_u / norm_u << "\n";
466  std::cout << "|| p_h - p_ex || / || p_ex || = " << err_p / norm_p << "\n";
467  }
468 
469  // 14. Save the refined mesh and the solution in parallel. This output can be
470  // viewed later using GLVis: "glvis -np <np> -m mesh -g sol_*".
471  {
472  ostringstream mesh_name, u_name, p_name;
473  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
474  u_name << "sol_u." << setfill('0') << setw(6) << myid;
475  p_name << "sol_p." << setfill('0') << setw(6) << myid;
476 
477  ofstream mesh_ofs(mesh_name.str().c_str());
478  mesh_ofs.precision(8);
479  pmesh->Print(mesh_ofs);
480 
481  ofstream u_ofs(u_name.str().c_str());
482  u_ofs.precision(8);
483  u->Save(u_ofs);
484 
485  ofstream p_ofs(p_name.str().c_str());
486  p_ofs.precision(8);
487  p->Save(p_ofs);
488  }
489 
490  // 15. Save data in the VisIt format
491  VisItDataCollection visit_dc("Example5-Parallel", pmesh);
492  visit_dc.RegisterField("velocity", u);
493  visit_dc.RegisterField("pressure", p);
494  visit_dc.SetFormat(!par_format ?
495  DataCollection::SERIAL_FORMAT :
496  DataCollection::PARALLEL_FORMAT);
497  visit_dc.Save();
498 
499  // 16. Send the solution by socket to a GLVis server.
500  if (visualization)
501  {
502  char vishost[] = "localhost";
503  int visport = 19916;
504  socketstream u_sock(vishost, visport);
505  u_sock << "parallel " << num_procs << " " << myid << "\n";
506  u_sock.precision(8);
507  u_sock << "solution\n" << *pmesh << *u << "window_title 'Velocity'"
508  << endl;
509  // Make sure all ranks have sent their 'u' solution before initiating
510  // another set of GLVis connections (one from each rank):
511  MPI_Barrier(pmesh->GetComm());
512  socketstream p_sock(vishost, visport);
513  p_sock << "parallel " << num_procs << " " << myid << "\n";
514  p_sock.precision(8);
515  p_sock << "solution\n" << *pmesh << *p << "window_title 'Pressure'"
516  << endl;
517  }
518 
519  // 17. Free the used memory.
520  delete fform;
521  delete gform;
522  delete u;
523  delete p;
524  delete darcyOp;
525  delete darcyPr;
526  delete pdarcyPr;
527  delete invM;
528  delete invS;
529  delete S;
530  delete BT;
531  delete B;
532  delete M;
533  delete pBT;
534  delete pB;
535  delete pM;
536  delete MinvBt;
537  delete Md;
538  delete mVarf;
539  delete bVarf;
540  delete W_space;
541  delete R_space;
542  delete l2_coll;
543  delete hdiv_coll;
544  delete pmesh;
545 
546  // We finalize PETSc
547  if (use_petsc) { MFEMFinalizePetsc(); }
548 
549  MPI_Finalize();
550 
551  return 0;
552 }
553 
554 
555 void uFun_ex(const Vector & x, Vector & u)
556 {
557  double xi(x(0));
558  double yi(x(1));
559  double zi(0.0);
560  if (x.Size() == 3)
561  {
562  zi = x(2);
563  }
564 
565  u(0) = - exp(xi)*sin(yi)*cos(zi);
566  u(1) = - exp(xi)*cos(yi)*cos(zi);
567 
568  if (x.Size() == 3)
569  {
570  u(2) = exp(xi)*sin(yi)*sin(zi);
571  }
572 }
573 
574 // Change if needed
575 double pFun_ex(const Vector & x)
576 {
577  double xi(x(0));
578  double yi(x(1));
579  double zi(0.0);
580 
581  if (x.Size() == 3)
582  {
583  zi = x(2);
584  }
585 
586  return exp(xi)*sin(yi)*cos(zi);
587 }
588 
589 void fFun(const Vector & x, Vector & f)
590 {
591  f = 0.0;
592 }
593 
594 double gFun(const Vector & x)
595 {
596  if (x.Size() == 3)
597  {
598  return -pFun_ex(x);
599  }
600  else
601  {
602  return 0;
603  }
604 }
605 
606 double f_natural(const Vector & x)
607 {
608  return (-pFun_ex(x));
609 }
int Size() const
Logical size of the array.
Definition: array.hpp:124
void SetPreconditioner(Solver &precond)
Definition: petsc.cpp:2500
void SetPrintLevel(int plev)
Definition: petsc.cpp:1856
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:381
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:770
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
int GetNumIterations() const
Definition: solvers.hpp:68
Abstract class for PETSc&#39;s preconditioners.
Definition: petsc.hpp:663
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:890
void fFun(const Vector &x, Vector &f)
Definition: ex5.cpp:346
A class to handle Vectors in a block fashion.
Definition: blockvector.hpp:30
virtual void SetOperator(const Operator &op)
Definition: petsc.cpp:2354
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs, int component=-1) const
Determine the boundary degrees of freedom.
Definition: pfespace.cpp:756
Subclass constant coefficient.
Definition: coefficient.hpp:67
double GetFinalNorm()
Definition: petsc.cpp:2122
Wrapper for PETSc&#39;s matrix class.
Definition: petsc.hpp:197
Abstract class for PETSc&#39;s linear solvers.
Definition: petsc.hpp:610
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:693
virtual void Save()
Save the collection and a VisIt root file.
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:485
double RealTime()
Definition: tic_toc.cpp:426
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetRelTol(double tol)
Definition: petsc.cpp:1779
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:504
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
MINRES method.
Definition: solvers.hpp:262
void SetMaxIter(int max_iter)
Definition: petsc.cpp:1829
int main(int argc, char *argv[])
Definition: ex1.cpp:62
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
Definition: hypre.cpp:1610
void Update(ParFiniteElementSpace *pf=NULL)
Update the object according to the given new FE space *pf.
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:70
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P_test^t A P_trial.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds a domain integrator. Assumes ownership of bfi.
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Definition: pgridfunc.cpp:102
The BoomerAMG solver in hypre.
Definition: hypre.hpp:951
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1144
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:421
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:275
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
Jacobi preconditioner in hypre.
Definition: hypre.hpp:859
double f_natural(const Vector &x)
Definition: ex5.cpp:363
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7982
Auxiliary class for BDDC customization.
Definition: petsc.hpp:689
Timing object.
Definition: tic_toc.hpp:34
MPI_Comm GetComm() const
Get the associated MPI communicator.
Definition: petsc.cpp:823
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
Definition: solvers.hpp:65
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
void Assemble(int skip_zeros=1)
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:251
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4102
virtual void SetFormat(int fmt)
Set the desired output mesh and data format.
double pFun_ex(const Vector &x)
Definition: ex5.cpp:332
HypreParMatrix * Transpose() const
Returns the transpose of *this.
Definition: hypre.cpp:1003
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:939
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1130
int GetConverged() const
Definition: solvers.hpp:69
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:53
int Dimension() const
Definition: mesh.hpp:744
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
int GetNumIterations()
Definition: petsc.cpp:2089
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:191
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:73
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:229
void SetAbsTol(double atol)
Definition: solvers.hpp:64
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:189
void SetRelTol(double rtol)
Definition: solvers.hpp:63
MPI_Comm GetComm() const
Definition: pmesh.hpp:230
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:76
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:635
void PartialSum()
Partial Sum.
Definition: array.cpp:103
double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:132
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
void MFEMFinalizePetsc()
Definition: petsc.cpp:182
void uFun_ex(const Vector &x, Vector &u)
Definition: ex5.cpp:312
double gFun(const Vector &x)
Definition: ex5.cpp:351
Class for parallel bilinear form using different test and trial FE spaces.
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:112
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:250
int GetConverged()
Definition: petsc.cpp:2056
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:1235
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
Definition: petsc.cpp:165
int dim
Definition: ex24.cpp:43
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void SetAbsTol(double tol)
Definition: petsc.cpp:1804
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:684
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:231
class for C-function coefficient
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
Definition: handle.hpp:106
Vector data type.
Definition: vector.hpp:48
HYPRE_Int * GetRowStarts() const
Return the parallel row partitioning array.
Definition: hypre.hpp:431
Class for parallel grid function.
Definition: pgridfunc.hpp:32
void SetNatBdrDofs(const Array< int > *natdofs, bool loc=false)
Specify dofs on the natural boundary.
Definition: petsc.hpp:717
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
A class to handle Block systems in a matrix-free implementation.
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
Definition: petsc.cpp:2580
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:376
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:34
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:84
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:143
bool Good() const
Definition: optparser.hpp:125