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 // 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  const char *petscrc_file = "";
68 
69  OptionsParser args(argc, argv);
70  args.AddOption(&mesh_file, "-m", "--mesh",
71  "Mesh file to use.");
72  args.AddOption(&order, "-o", "--order",
73  "Finite element order (polynomial degree).");
74  args.AddOption(&par_format, "-pf", "--parallel-format", "-sf",
75  "--serial-format",
76  "Format to use when saving the results for VisIt.");
77  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
78  "--no-visualization",
79  "Enable or disable GLVis visualization.");
80  args.AddOption(&use_petsc, "-usepetsc", "--usepetsc", "-no-petsc",
81  "--no-petsc",
82  "Use or not PETSc to solve the linear system.");
83  args.AddOption(&petscrc_file, "-petscopts", "--petscopts",
84  "PetscOptions file to use.");
85  args.AddOption(&use_nonoverlapping, "-nonoverlapping", "--nonoverlapping",
86  "-no-nonoverlapping", "--no-nonoverlapping",
87  "Use or not the block diagonal PETSc's matrix format "
88  "for non-overlapping domain decomposition.");
89  args.Parse();
90  if (!args.Good())
91  {
92  if (verbose)
93  {
94  args.PrintUsage(cout);
95  }
96  MPI_Finalize();
97  return 1;
98  }
99  if (verbose)
100  {
101  args.PrintOptions(cout);
102  }
103  // 2b. We initialize PETSc
104  if (use_petsc) { PetscInitialize(NULL,NULL,petscrc_file,NULL); }
105 
106  // 3. Read the (serial) mesh from the given mesh file on all processors. We
107  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
108  // and volume meshes with the same code.
109  Mesh *mesh = new Mesh(mesh_file, 1, 1);
110  int dim = mesh->Dimension();
111 
112  // 4. Refine the serial mesh on all processors to increase the resolution. In
113  // this example we do 'ref_levels' of uniform refinement. We choose
114  // 'ref_levels' to be the largest number that gives a final mesh with no
115  // more than 10,000 elements.
116  {
117  int ref_levels =
118  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
119  for (int l = 0; l < ref_levels; l++)
120  {
121  mesh->UniformRefinement();
122  }
123  }
124 
125  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
126  // this mesh further in parallel to increase the resolution. Once the
127  // parallel mesh is defined, the serial mesh can be deleted.
128  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
129  delete mesh;
130  {
131  int par_ref_levels = 2;
132  for (int l = 0; l < par_ref_levels; l++)
133  {
134  pmesh->UniformRefinement();
135  }
136  }
137 
138  // 6. Define a parallel finite element space on the parallel mesh. Here we
139  // use the Raviart-Thomas finite elements of the specified order.
140  FiniteElementCollection *hdiv_coll(new RT_FECollection(order, dim));
141  FiniteElementCollection *l2_coll(new L2_FECollection(order, dim));
142 
143  ParFiniteElementSpace *R_space = new ParFiniteElementSpace(pmesh, hdiv_coll);
144  ParFiniteElementSpace *W_space = new ParFiniteElementSpace(pmesh, l2_coll);
145 
146  HYPRE_Int dimR = R_space->GlobalTrueVSize();
147  HYPRE_Int dimW = W_space->GlobalTrueVSize();
148 
149  if (verbose)
150  {
151  std::cout << "***********************************************************\n";
152  std::cout << "dim(R) = " << dimR << "\n";
153  std::cout << "dim(W) = " << dimW << "\n";
154  std::cout << "dim(R+W) = " << dimR + dimW << "\n";
155  std::cout << "***********************************************************\n";
156  }
157 
158  // 7. Define the two BlockStructure of the problem. block_offsets is used
159  // for Vector based on dof (like ParGridFunction or ParLinearForm),
160  // block_trueOffstes is used for Vector based on trueDof (HypreParVector
161  // for the rhs and solution of the linear system). The offsets computed
162  // here are local to the processor.
163  Array<int> block_offsets(3); // number of variables + 1
164  block_offsets[0] = 0;
165  block_offsets[1] = R_space->GetVSize();
166  block_offsets[2] = W_space->GetVSize();
167  block_offsets.PartialSum();
168 
169  Array<int> block_trueOffsets(3); // number of variables + 1
170  block_trueOffsets[0] = 0;
171  block_trueOffsets[1] = R_space->TrueVSize();
172  block_trueOffsets[2] = W_space->TrueVSize();
173  block_trueOffsets.PartialSum();
174 
175  // 8. Define the coefficients, analytical solution, and rhs of the PDE.
176  ConstantCoefficient k(1.0);
177 
178  VectorFunctionCoefficient fcoeff(dim, fFun);
179  FunctionCoefficient fnatcoeff(f_natural);
180  FunctionCoefficient gcoeff(gFun);
181 
182  VectorFunctionCoefficient ucoeff(dim, uFun_ex);
184 
185  // 9. Define the parallel grid function and parallel linear forms, solution
186  // vector and rhs.
187  BlockVector x(block_offsets), rhs(block_offsets);
188  BlockVector trueX(block_trueOffsets), trueRhs(block_trueOffsets);
189 
190  ParLinearForm *fform(new ParLinearForm);
191  fform->Update(R_space, rhs.GetBlock(0), 0);
194  fform->Assemble();
195  fform->ParallelAssemble(trueRhs.GetBlock(0));
196 
197  ParLinearForm *gform(new ParLinearForm);
198  gform->Update(W_space, rhs.GetBlock(1), 0);
199  gform->AddDomainIntegrator(new DomainLFIntegrator(gcoeff));
200  gform->Assemble();
201  gform->ParallelAssemble(trueRhs.GetBlock(1));
202 
203  // 10. Assemble the finite element matrices for the Darcy operator
204  //
205  // D = [ M B^T ]
206  // [ B 0 ]
207  // where:
208  //
209  // M = \int_\Omega k u_h \cdot v_h d\Omega u_h, v_h \in R_h
210  // B = -\int_\Omega \div u_h q_h d\Omega u_h \in R_h, q_h \in W_h
211  ParBilinearForm *mVarf(new ParBilinearForm(R_space));
212  ParMixedBilinearForm *bVarf(new ParMixedBilinearForm(R_space, W_space));
213 
214  PetscParMatrix *pM = NULL, *pB = NULL, *pBT = NULL;
215  HypreParMatrix *M = NULL, *B = NULL, *BT = NULL;
216  Operator::Type tid =
217  !use_petsc ? Operator::Hypre_ParCSR :
218  (use_nonoverlapping ? Operator::PETSC_MATIS : Operator::PETSC_MATAIJ);
219  OperatorHandle Mh(tid), Bh(tid);
220 
222  mVarf->Assemble();
223  mVarf->Finalize();
224  mVarf->ParallelAssemble(Mh);
225  if (!use_petsc) { Mh.Get(M); }
226  else { Mh.Get(pM); }
227  Mh.SetOperatorOwner(false);
228 
230  bVarf->Assemble();
231  bVarf->Finalize();
232  if (!use_petsc)
233  {
234  B = bVarf->ParallelAssemble();
235  (*B) *= -1;
236  }
237  else
238  {
239  bVarf->ParallelAssemble(Bh);
240  Bh.Get(pB);
241  Bh.SetOperatorOwner(false);
242  (*pB) *= -1;
243  }
244 
245  if (!use_petsc) { BT = B->Transpose(); }
246  else { pBT = pB->Transpose(); };
247 
248  Operator *darcyOp = NULL;
249  if (!use_petsc)
250  {
251  BlockOperator *tdarcyOp = new BlockOperator(block_trueOffsets);
252  tdarcyOp->SetBlock(0,0,M);
253  tdarcyOp->SetBlock(0,1,BT);
254  tdarcyOp->SetBlock(1,0,B);
255  darcyOp = tdarcyOp;
256  }
257  else
258  {
259  // We construct the BlockOperator and we then convert it to a
260  // PetscParMatrix to avoid any conversion in the construction of the
261  // preconditioners.
262  BlockOperator *tdarcyOp = new BlockOperator(block_trueOffsets);
263  tdarcyOp->SetBlock(0,0,pM);
264  tdarcyOp->SetBlock(0,1,pBT);
265  tdarcyOp->SetBlock(1,0,pB);
266  darcyOp = new PetscParMatrix(pM->GetComm(),tdarcyOp,
267  use_nonoverlapping ? Operator::PETSC_MATIS :
268  Operator::PETSC_MATAIJ);
269  delete tdarcyOp;
270  }
271 
272  // 11. Construct the operators for preconditioner
273  //
274  // P = [ diag(M) 0 ]
275  // [ 0 B diag(M)^-1 B^T ]
276  //
277  // Here we use Symmetric Gauss-Seidel to approximate the inverse of the
278  // pressure Schur Complement.
279  PetscPreconditioner *pdarcyPr = NULL;
280  BlockDiagonalPreconditioner *darcyPr = NULL;
281  HypreSolver *invM = NULL, *invS = NULL;
282  HypreParMatrix *S = NULL;
283  HypreParMatrix *MinvBt = NULL;
284  HypreParVector *Md = NULL;
285  if (!use_petsc)
286  {
287  MinvBt = B->Transpose();
288  Md = new HypreParVector(MPI_COMM_WORLD, M->GetGlobalNumRows(),
289  M->GetRowStarts());
290  M->GetDiag(*Md);
291 
292  MinvBt->InvScaleRows(*Md);
293  S = ParMult(B, MinvBt);
294 
295  invM = new HypreDiagScale(*M);
296  invS = new HypreBoomerAMG(*S);
297 
298  invM->iterative_mode = false;
299  invS->iterative_mode = false;
300 
301  darcyPr = new BlockDiagonalPreconditioner(block_trueOffsets);
302  darcyPr->SetDiagonalBlock(0, invM);
303  darcyPr->SetDiagonalBlock(1, invS);
304  }
305  else
306  {
307  if (use_nonoverlapping)
308  {
309  // For saddle point problems, we need to provide BDDC the list of
310  // boundary dofs either essential or natural.
311  // Since R_space is the only space that may have boundary dofs and it
312  // is ordered first then W_space, we don't need any local offset when
313  // specifying the dofs.
314  Array<int> bdr_tdof_list;
315  bool local = false;
316  if (pmesh->bdr_attributes.Size())
317  {
318  Array<int> bdr(pmesh->bdr_attributes.Max());
319  bdr = 1;
320 
321  R_space->GetEssentialTrueDofs(bdr, bdr_tdof_list);
322  local = false;
323  // Alternatively, you can also provide the list of dofs in local
324  // ordering:
325  // R_space->GetEssentialVDofs(bdr, bdr_tdof_list);
326  // bdr_tdof_list.SetSize(R_space->GetVSize()+W_space->GetVSize(),0);
327  // local = true;
328  }
329  else
330  {
331  MFEM_ABORT("Need to know the boundary dofs");
332  }
333 
335  opts.SetNatBdrDofs(&bdr_tdof_list,local);
336  // See also command line options rc_ex5p_bddc
337  pdarcyPr = new PetscBDDCSolver(MPI_COMM_WORLD,*darcyOp,opts,"prec_");
338  }
339  else
340  {
341  // With PETSc, we can construct the (same) block-diagonal solver with
342  // command line options (see rc_ex5p_fieldsplit)
343  pdarcyPr = new PetscFieldSplitSolver(MPI_COMM_WORLD,*darcyOp,"prec_");
344  }
345  }
346 
347  // 12. Solve the linear system with MINRES.
348  // Check the norm of the unpreconditioned residual.
349 
350  int maxIter(500);
351  double rtol(1.e-6);
352  double atol(1.e-10);
353 
354  chrono.Clear();
355  chrono.Start();
356 
357  trueX = 0.0;
358  if (!use_petsc)
359  {
360  MINRESSolver solver(MPI_COMM_WORLD);
361  solver.SetAbsTol(atol);
362  solver.SetRelTol(rtol);
363  solver.SetMaxIter(maxIter);
364  solver.SetOperator(*darcyOp);
365  solver.SetPreconditioner(*darcyPr);
366  solver.SetPrintLevel(1);
367  solver.Mult(trueRhs, trueX);
368  if (verbose)
369  {
370  if (solver.GetConverged())
371  {
372  std::cout << "MINRES converged in " << solver.GetNumIterations()
373  << " iterations with a residual norm of "
374  << solver.GetFinalNorm() << ".\n";
375  }
376  else
377  {
378  std::cout << "MINRES did not converge in "
379  << solver.GetNumIterations()
380  << " iterations. Residual norm is "
381  << solver.GetFinalNorm() << ".\n";
382  }
383  std::cout << "MINRES solver took " << chrono.RealTime() << "s. \n";
384  }
385 
386  }
387  else
388  {
389  std::string solvertype;
390  PetscLinearSolver *solver;
391  if (use_nonoverlapping)
392  {
393  // We can use conjugate gradients to solve the problem
394  solver = new PetscPCGSolver(MPI_COMM_WORLD);
395  solvertype = "PCG";
396  }
397  else
398  {
399  solver = new PetscLinearSolver(MPI_COMM_WORLD);
400  solvertype = "MINRES";
401  }
402  solver->SetOperator(*darcyOp);
403  solver->SetPreconditioner(*pdarcyPr);
404  solver->SetAbsTol(atol);
405  solver->SetRelTol(rtol);
406  solver->SetMaxIter(maxIter);
407  solver->SetPrintLevel(2);
408  solver->Mult(trueRhs, trueX);
409  if (verbose)
410  {
411  if (solver->GetConverged())
412  {
413  std::cout << solvertype << " converged in "
414  << solver->GetNumIterations()
415  << " iterations with a residual norm of "
416  << solver->GetFinalNorm() << ".\n";
417  }
418  else
419  {
420  std::cout << solvertype << " did not converge in "
421  << solver->GetNumIterations()
422  << " iterations. Residual norm is "
423  << solver->GetFinalNorm() << ".\n";
424  }
425  std::cout << solvertype << " solver took "
426  << chrono.RealTime() << "s. \n";
427  }
428  delete solver;
429  }
430  chrono.Stop();
431 
432  // 13. Extract the parallel grid function corresponding to the finite element
433  // approximation X. This is the local solution on each processor. Compute
434  // L2 error norms.
437  u->MakeRef(R_space, x.GetBlock(0), 0);
438  p->MakeRef(W_space, x.GetBlock(1), 0);
439  u->Distribute(&(trueX.GetBlock(0)));
440  p->Distribute(&(trueX.GetBlock(1)));
441 
442  int order_quad = max(2, 2*order+1);
443  const IntegrationRule *irs[Geometry::NumGeom];
444  for (int i=0; i < Geometry::NumGeom; ++i)
445  {
446  irs[i] = &(IntRules.Get(i, order_quad));
447  }
448 
449  double err_u = u->ComputeL2Error(ucoeff, irs);
450  double norm_u = ComputeGlobalLpNorm(2, ucoeff, *pmesh, irs);
451  double err_p = p->ComputeL2Error(pcoeff, irs);
452  double norm_p = ComputeGlobalLpNorm(2, pcoeff, *pmesh, irs);
453 
454  if (verbose)
455  {
456  std::cout << "|| u_h - u_ex || / || u_ex || = " << err_u / norm_u << "\n";
457  std::cout << "|| p_h - p_ex || / || p_ex || = " << err_p / norm_p << "\n";
458  }
459 
460  // 14. Save the refined mesh and the solution in parallel. This output can be
461  // viewed later using GLVis: "glvis -np <np> -m mesh -g sol_*".
462  {
463  ostringstream mesh_name, u_name, p_name;
464  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
465  u_name << "sol_u." << setfill('0') << setw(6) << myid;
466  p_name << "sol_p." << setfill('0') << setw(6) << myid;
467 
468  ofstream mesh_ofs(mesh_name.str().c_str());
469  mesh_ofs.precision(8);
470  pmesh->Print(mesh_ofs);
471 
472  ofstream u_ofs(u_name.str().c_str());
473  u_ofs.precision(8);
474  u->Save(u_ofs);
475 
476  ofstream p_ofs(p_name.str().c_str());
477  p_ofs.precision(8);
478  p->Save(p_ofs);
479  }
480 
481  // 15. Save data in the VisIt format
482  VisItDataCollection visit_dc("Example5-Parallel", pmesh);
483  visit_dc.RegisterField("velocity", u);
484  visit_dc.RegisterField("pressure", p);
485  visit_dc.SetFormat(!par_format ?
486  DataCollection::SERIAL_FORMAT :
487  DataCollection::PARALLEL_FORMAT);
488  visit_dc.Save();
489 
490  // 16. Send the solution by socket to a GLVis server.
491  if (visualization)
492  {
493  char vishost[] = "localhost";
494  int visport = 19916;
495  socketstream u_sock(vishost, visport);
496  u_sock << "parallel " << num_procs << " " << myid << "\n";
497  u_sock.precision(8);
498  u_sock << "solution\n" << *pmesh << *u << "window_title 'Velocity'"
499  << endl;
500  // Make sure all ranks have sent their 'u' solution before initiating
501  // another set of GLVis connections (one from each rank):
502  MPI_Barrier(pmesh->GetComm());
503  socketstream p_sock(vishost, visport);
504  p_sock << "parallel " << num_procs << " " << myid << "\n";
505  p_sock.precision(8);
506  p_sock << "solution\n" << *pmesh << *p << "window_title 'Pressure'"
507  << endl;
508  }
509 
510  // 17. Free the used memory.
511  delete fform;
512  delete gform;
513  delete u;
514  delete p;
515  delete darcyOp;
516  delete darcyPr;
517  delete pdarcyPr;
518  delete invM;
519  delete invS;
520  delete S;
521  delete BT;
522  delete B;
523  delete M;
524  delete pBT;
525  delete pB;
526  delete pM;
527  delete MinvBt;
528  delete Md;
529  delete mVarf;
530  delete bVarf;
531  delete W_space;
532  delete R_space;
533  delete l2_coll;
534  delete hdiv_coll;
535  delete pmesh;
536 
537  // We finalize PETSc
538  if (use_petsc) { PetscFinalize(); }
539 
540  MPI_Finalize();
541 
542  return 0;
543 }
544 
545 
546 void uFun_ex(const Vector & x, Vector & u)
547 {
548  double xi(x(0));
549  double yi(x(1));
550  double zi(0.0);
551  if (x.Size() == 3)
552  {
553  zi = x(2);
554  }
555 
556  u(0) = - exp(xi)*sin(yi)*cos(zi);
557  u(1) = - exp(xi)*cos(yi)*cos(zi);
558 
559  if (x.Size() == 3)
560  {
561  u(2) = exp(xi)*sin(yi)*sin(zi);
562  }
563 }
564 
565 // Change if needed
566 double pFun_ex(const Vector & x)
567 {
568  double xi(x(0));
569  double yi(x(1));
570  double zi(0.0);
571 
572  if (x.Size() == 3)
573  {
574  zi = x(2);
575  }
576 
577  return exp(xi)*sin(yi)*cos(zi);
578 }
579 
580 void fFun(const Vector & x, Vector & f)
581 {
582  f = 0.0;
583 }
584 
585 double gFun(const Vector & x)
586 {
587  if (x.Size() == 3)
588  {
589  return -pFun_ex(x);
590  }
591  else
592  {
593  return 0;
594  }
595 }
596 
597 double f_natural(const Vector & x)
598 {
599  return (-pFun_ex(x));
600 }
int Size() const
Logical size of the array.
Definition: array.hpp:133
void SetPreconditioner(Solver &precond)
Sets the solver to perform preconditioning.
Definition: petsc.cpp:2003
void SetPrintLevel(int plev)
Definition: petsc.cpp:1407
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
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:626
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:83
int GetNumIterations() const
Definition: solvers.hpp:66
Abstract class for PETSc&#39;s preconditioners.
Definition: petsc.hpp:523
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
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: petsc.cpp:1839
Subclass constant coefficient.
Definition: coefficient.hpp:57
double GetFinalNorm()
Definition: petsc.cpp:1654
Wrapper for PETSc&#39;s matrix class.
Definition: petsc.hpp:133
Abstract class for PETSc&#39;s linear solvers.
Definition: petsc.hpp:476
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:51
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
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
void SetRelTol(double tol)
Definition: petsc.cpp:1330
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
void SetMaxIter(int max_iter)
Definition: petsc.cpp:1380
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
Auxiliary class for BDDC customization.
Definition: petsc.hpp:545
Timing object.
Definition: tic_toc.hpp:34
MPI_Comm GetComm() const
Get the associated MPI communicator.
Definition: petsc.hpp:244
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
Definition: solvers.hpp:63
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:109
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
int GetNumIterations()
Definition: petsc.cpp:1621
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
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:124
void SetAbsTol(double atol)
Definition: solvers.hpp:62
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:174
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.
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:106
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:215
int GetConverged()
Definition: petsc.cpp:1588
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 SetAbsTol(double tol)
Definition: petsc.cpp:1355
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
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
Definition: handle.hpp:100
Vector data type.
Definition: vector.hpp:48
HYPRE_Int * GetRowStarts() const
Definition: hypre.hpp:385
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:572
Abstract operator.
Definition: operator.hpp:21
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:175
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:2043
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