MFEM  v3.3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
miniapps/performance/ex1p.cpp
Go to the documentation of this file.
1 // MFEM Example 1 - Parallel High-Performance Version
2 //
3 // Compile with: make ex1p
4 //
5 // Sample runs: mpirun -np 4 ex1p -m ../../data/fichera.mesh -perf -mf -pc lor
6 // mpirun -np 4 ex1p -m ../../data/fichera.mesh -perf -asm -pc ho
7 // mpirun -np 4 ex1p -m ../../data/fichera.mesh -perf -asm -pc ho -sc
8 // mpirun -np 4 ex1p -m ../../data/fichera.mesh -std -asm -pc ho
9 // mpirun -np 4 ex1p -m ../../data/fichera.mesh -std -asm -pc ho -sc
10 // mpirun -np 4 ex1p -m ../../data/amr-hex.mesh -perf -asm -pc ho -sc
11 // mpirun -np 4 ex1p -m ../../data/amr-hex.mesh -std -asm -pc ho -sc
12 // mpirun -np 4 ex1p -m ../../data/ball-nurbs.mesh -perf -asm -pc ho -sc
13 // mpirun -np 4 ex1p -m ../../data/ball-nurbs.mesh -std -asm -pc ho -sc
14 // mpirun -np 4 ex1p -m ../../data/pipe-nurbs.mesh -perf -mf -pc lor
15 // mpirun -np 4 ex1p -m ../../data/pipe-nurbs.mesh -std -asm -pc ho -sc
16 //
17 // Description: This example code demonstrates the use of MFEM to define a
18 // simple finite element discretization of the Laplace problem
19 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
20 // Specifically, we discretize using a FE space of the specified
21 // order, or if order < 1 using an isoparametric/isogeometric
22 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
23 // NURBS mesh, etc.)
24 //
25 // The example highlights the use of mesh refinement, finite
26 // element grid functions, as well as linear and bilinear forms
27 // corresponding to the left-hand side and right-hand side of the
28 // discrete linear system. We also cover the explicit elimination
29 // of essential boundary conditions, static condensation, and the
30 // optional connection to the GLVis tool for visualization.
31 
32 #include "mfem-performance.hpp"
33 #include <fstream>
34 #include <iostream>
35 
36 using namespace std;
37 using namespace mfem;
38 
39 // Define template parameters for optimized build.
40 const Geometry::Type geom = Geometry::CUBE; // mesh elements (default: hex)
41 const int mesh_p = 3; // mesh curvature (default: 3)
42 const int sol_p = 3; // solution order (default: 3)
44 const int ir_order = 2*sol_p+rdim-1;
45 
46 // Static mesh type
50 
51 // Static solution finite element space type
54 
55 // Static quadrature, coefficient and integrator types
59 
60 // Static bilinear form type, combining the above types
62 
63 int main(int argc, char *argv[])
64 {
65  // 1. Initialize MPI.
66  int num_procs, myid;
67  MPI_Init(&argc, &argv);
68  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
69  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
70 
71  // 2. Parse command-line options.
72  const char *mesh_file = "../../data/fichera.mesh";
73  int ser_ref_levels = -1;
74  int par_ref_levels = 1;
75  int order = sol_p;
76  const char *basis_type = "G"; // Gauss-Lobatto
77  bool static_cond = false;
78  const char *pc = "lor";
79  bool perf = true;
80  bool matrix_free = true;
81  bool visualization = 1;
82 
83  OptionsParser args(argc, argv);
84  args.AddOption(&mesh_file, "-m", "--mesh",
85  "Mesh file to use.");
86  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
87  "Number of times to refine the mesh uniformly in serial;"
88  " -1 = auto: <= 10,000 elements.");
89  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
90  "Number of times to refine the mesh uniformly in parallel.");
91  args.AddOption(&order, "-o", "--order",
92  "Finite element order (polynomial degree) or -1 for"
93  " isoparametric space.");
94  args.AddOption(&basis_type, "-b", "--basis-type",
95  "Basis: G - Gauss-Lobatto, P - Positive, U - Uniform");
96  args.AddOption(&perf, "-perf", "--hpc-version", "-std", "--standard-version",
97  "Enable high-performance, tensor-based, assembly/evaluation.");
98  args.AddOption(&matrix_free, "-mf", "--matrix-free", "-asm", "--assembly",
99  "Use matrix-free evaluation or efficient matrix assembly in "
100  "the high-performance version.");
101  args.AddOption(&pc, "-pc", "--preconditioner",
102  "Preconditioner: lor - low-order-refined (matrix-free) AMG, "
103  "ho - high-order (assembled) AMG, none.");
104  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
105  "--no-static-condensation", "Enable static condensation.");
106  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
107  "--no-visualization",
108  "Enable or disable GLVis visualization.");
109  args.Parse();
110  if (!args.Good())
111  {
112  if (myid == 0)
113  {
114  args.PrintUsage(cout);
115  }
116  MPI_Finalize();
117  return 1;
118  }
119  if (static_cond && perf && matrix_free)
120  {
121  if (myid == 0)
122  {
123  cout << "\nStatic condensation can not be used with matrix-free"
124  " evaluation!\n" << endl;
125  }
126  MPI_Finalize();
127  return 2;
128  }
129  MFEM_VERIFY(perf || !matrix_free,
130  "--standard-version is not compatible with --matrix-free");
131  if (myid == 0)
132  {
133  args.PrintOptions(cout);
134  }
135 
136  enum PCType { NONE, LOR, HO };
137  PCType pc_choice;
138  if (!strcmp(pc, "ho")) { pc_choice = HO; }
139  else if (!strcmp(pc, "lor")) { pc_choice = LOR; }
140  else if (!strcmp(pc, "none")) { pc_choice = NONE; }
141  else
142  {
143  mfem_error("Invalid Preconditioner specified");
144  return 3;
145  }
146 
147  // See class BasisType in fem/fe_coll.hpp for available basis types
148  int basis = BasisType::GetType(basis_type[0]);
149  if (myid == 0)
150  {
151  cout << "Using " << BasisType::Name(basis) << " basis ..." << endl;
152  }
153 
154  // 3. Read the (serial) mesh from the given mesh file on all processors. We
155  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
156  // and volume meshes with the same code.
157  Mesh *mesh = new Mesh(mesh_file, 1, 1);
158  int dim = mesh->Dimension();
159 
160  // 4. Check if the optimized version matches the given mesh
161  if (perf)
162  {
163  if (myid == 0)
164  {
165  cout << "High-performance version using integration rule with "
166  << int_rule_t::qpts << " points ..." << endl;
167  }
168  if (!mesh_t::MatchesGeometry(*mesh))
169  {
170  if (myid == 0)
171  {
172  cout << "The given mesh does not match the optimized 'geom' parameter.\n"
173  << "Recompile with suitable 'geom' value." << endl;
174  }
175  delete mesh;
176  MPI_Finalize();
177  return 4;
178  }
179  else if (!mesh_t::MatchesNodes(*mesh))
180  {
181  if (myid == 0)
182  {
183  cout << "Switching the mesh curvature to match the "
184  << "optimized value (order " << mesh_p << ") ..." << endl;
185  }
186  mesh->SetCurvature(mesh_p, false, -1, Ordering::byNODES);
187  }
188  }
189 
190  // 5. Refine the serial mesh on all processors to increase the resolution. In
191  // this example we do 'ref_levels' of uniform refinement. We choose
192  // 'ref_levels' to be the largest number that gives a final mesh with no
193  // more than 10,000 elements, or as specified on the command line with the
194  // option '--refine-serial'.
195  {
196  int ref_levels = (ser_ref_levels != -1) ? ser_ref_levels :
197  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
198  for (int l = 0; l < ref_levels; l++)
199  {
200  mesh->UniformRefinement();
201  }
202  }
203  if (!perf && mesh->NURBSext)
204  {
205  const int new_mesh_p = std::min(sol_p, mesh_p);
206  if (myid == 0)
207  {
208  cout << "NURBS mesh: switching the mesh curvature to be "
209  << "min(sol_p, mesh_p) = " << new_mesh_p << " ..." << endl;
210  }
211  mesh->SetCurvature(new_mesh_p, false, -1, Ordering::byNODES);
212  }
213 
214  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
215  // this mesh further in parallel to increase the resolution. Once the
216  // parallel mesh is defined, the serial mesh can be deleted.
217  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
218  delete mesh;
219  {
220  for (int l = 0; l < par_ref_levels; l++)
221  {
222  pmesh->UniformRefinement();
223  }
224  }
225  if (pmesh->MeshGenerator() & 1) // simplex mesh
226  {
227  MFEM_VERIFY(pc_choice != LOR, "triangle and tet meshes do not support"
228  " the LOR preconditioner yet");
229  }
230 
231  // 7. Define a parallel finite element space on the parallel mesh. Here we
232  // use continuous Lagrange finite elements of the specified order. If
233  // order < 1, we instead use an isoparametric/isogeometric space.
235  if (order > 0)
236  {
237  fec = new H1_FECollection(order, dim, basis);
238  }
239  else if (pmesh->GetNodes())
240  {
241  fec = pmesh->GetNodes()->OwnFEC();
242  if (myid == 0)
243  {
244  cout << "Using isoparametric FEs: " << fec->Name() << endl;
245  }
246  }
247  else
248  {
249  fec = new H1_FECollection(order = 1, dim, basis);
250  }
251  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
252  HYPRE_Int size = fespace->GlobalTrueVSize();
253  if (myid == 0)
254  {
255  cout << "Number of finite element unknowns: " << size << endl;
256  }
257 
258  ParMesh *pmesh_lor = NULL;
259  FiniteElementCollection *fec_lor = NULL;
260  ParFiniteElementSpace *fespace_lor = NULL;
261  if (pc_choice == LOR)
262  {
263  int basis_lor = basis;
264  if (basis == BasisType::Positive) { basis_lor=BasisType::ClosedUniform; }
265  pmesh_lor = new ParMesh(pmesh, order, basis_lor);
266  fec_lor = new H1_FECollection(1, dim);
267  fespace_lor = new ParFiniteElementSpace(pmesh_lor, fec_lor);
268  }
269 
270  // 8. Check if the optimized version matches the given space
271  if (perf && !sol_fes_t::Matches(*fespace))
272  {
273  if (myid == 0)
274  {
275  cout << "The given order does not match the optimized parameter.\n"
276  << "Recompile with suitable 'sol_p' value." << endl;
277  }
278  delete fespace;
279  delete fec;
280  delete mesh;
281  MPI_Finalize();
282  return 5;
283  }
284 
285  // 9. Determine the list of true (i.e. parallel conforming) essential
286  // boundary dofs. In this example, the boundary conditions are defined
287  // by marking all the boundary attributes from the mesh as essential
288  // (Dirichlet) and converting them to a list of true dofs.
289  Array<int> ess_tdof_list;
290  if (pmesh->bdr_attributes.Size())
291  {
292  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
293  ess_bdr = 1;
294  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
295  }
296 
297  // 10. Set up the parallel linear form b(.) which corresponds to the
298  // right-hand side of the FEM linear system, which in this case is
299  // (1,phi_i) where phi_i are the basis functions in fespace.
300  ParLinearForm *b = new ParLinearForm(fespace);
301  ConstantCoefficient one(1.0);
303  b->Assemble();
304 
305  // 11. Define the solution vector x as a parallel finite element grid
306  // function corresponding to fespace. Initialize x with initial guess of
307  // zero, which satisfies the boundary conditions.
308  ParGridFunction x(fespace);
309  x = 0.0;
310 
311  // 12. Set up the parallel bilinear form a(.,.) on the finite element space
312  // that will hold the matrix corresponding to the Laplacian operator.
313  ParBilinearForm *a = new ParBilinearForm(fespace);
314  ParBilinearForm *a_pc = NULL;
315  if (pc_choice == LOR) { a_pc = new ParBilinearForm(fespace_lor); }
316  if (pc_choice == HO) { a_pc = new ParBilinearForm(fespace); }
317 
318  // 13. Assemble the parallel bilinear form and the corresponding linear
319  // system, applying any necessary transformations such as: parallel
320  // assembly, eliminating boundary conditions, applying conforming
321  // constraints for non-conforming AMR, static condensation, etc.
322  if (static_cond)
323  {
325  MFEM_VERIFY(pc_choice != LOR,
326  "cannot use LOR preconditioner with static condensation");
327  }
328 
329  if (myid == 0)
330  {
331  cout << "Assembling the matrix ..." << flush;
332  }
333  tic_toc.Clear();
334  tic_toc.Start();
335  // Pre-allocate sparsity assuming dense element matrices
337 
338  HPCBilinearForm *a_hpc = NULL;
339  Operator *a_oper = NULL;
340 
341  if (!perf)
342  {
343  // Standard assembly using a diffusion domain integrator
345  a->Assemble();
346  }
347  else
348  {
349  // High-performance assembly/evaluation using the templated operator type
350  a_hpc = new HPCBilinearForm(integ_t(coeff_t(1.0)), *fespace);
351  if (matrix_free)
352  {
353  a_hpc->Assemble(); // partial assembly
354  }
355  else
356  {
357  a_hpc->AssembleBilinearForm(*a); // full matrix assembly
358  }
359  }
360  tic_toc.Stop();
361  if (myid == 0)
362  {
363  cout << " done, " << tic_toc.RealTime() << "s." << endl;
364  }
365 
366  // 14. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
367  // preconditioner from hypre.
368 
369  // Setup the operator matrix (if applicable)
370  HypreParMatrix A;
371  Vector B, X;
372  if (perf && matrix_free)
373  {
374  a_hpc->FormLinearSystem(ess_tdof_list, x, *b, a_oper, X, B);
375  HYPRE_Int glob_size = fespace->GlobalTrueVSize();
376  if (myid == 0)
377  {
378  cout << "Size of linear system: " << glob_size << endl;
379  }
380  }
381  else
382  {
383  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
384  HYPRE_Int glob_size = A.GetGlobalNumRows();
385  if (myid == 0)
386  {
387  cout << "Size of linear system: " << glob_size << endl;
388  }
389  a_oper = &A;
390  }
391 
392  // Setup the matrix used for preconditioning
393  if (myid == 0)
394  {
395  cout << "Assembling the preconditioning matrix ..." << flush;
396  }
397  tic_toc.Clear();
398  tic_toc.Start();
399 
400  HypreParMatrix A_pc;
401  if (pc_choice == LOR)
402  {
403  // TODO: assemble the LOR matrix using the performance code
405  a_pc->UsePrecomputedSparsity();
406  a_pc->Assemble();
407  a_pc->FormSystemMatrix(ess_tdof_list, A_pc);
408  }
409  else if (pc_choice == HO)
410  {
411  if (!matrix_free)
412  {
413  A_pc.MakeRef(A); // matrix already assembled, reuse it
414  }
415  else
416  {
417  a_pc->UsePrecomputedSparsity();
418  a_hpc->AssembleBilinearForm(*a_pc);
419  a_pc->FormSystemMatrix(ess_tdof_list, A_pc);
420  }
421  }
422  tic_toc.Stop();
423  if (myid == 0)
424  {
425  cout << " done, " << tic_toc.RealTime() << "s." << endl;
426  }
427 
428  // Solve with CG or PCG, depending if the matrix A_pc is available
429  CGSolver *pcg;
430  pcg = new CGSolver(MPI_COMM_WORLD);
431  pcg->SetRelTol(1e-6);
432  pcg->SetMaxIter(500);
433  pcg->SetPrintLevel(1);
434 
435  HypreSolver *amg = NULL;
436 
437  pcg->SetOperator(*a_oper);
438  if (pc_choice != NONE)
439  {
440  amg = new HypreBoomerAMG(A_pc);
441  pcg->SetPreconditioner(*amg);
442  }
443 
444  tic_toc.Clear();
445  tic_toc.Start();
446 
447  pcg->Mult(B, X);
448 
449  tic_toc.Stop();
450  delete amg;
451 
452  if (myid == 0)
453  {
454  // Note: In the pcg algorithm, the number of operator Mult() calls is
455  // N_iter and the number of preconditioner Mult() calls is N_iter+1.
456  cout << "Time per CG step: "
457  << tic_toc.RealTime() / pcg->GetNumIterations() << "s." << endl;
458  }
459 
460  // 15. Recover the parallel grid function corresponding to X. This is the
461  // local finite element solution on each processor.
462  if (perf && matrix_free)
463  {
464  a_hpc->RecoverFEMSolution(X, *b, x);
465  }
466  else
467  {
468  a->RecoverFEMSolution(X, *b, x);
469  }
470 
471  // 16. Save the refined mesh and the solution in parallel. This output can
472  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
473  {
474  ostringstream mesh_name, sol_name;
475  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
476  sol_name << "sol." << setfill('0') << setw(6) << myid;
477 
478  ofstream mesh_ofs(mesh_name.str().c_str());
479  mesh_ofs.precision(8);
480  pmesh->Print(mesh_ofs);
481 
482  ofstream sol_ofs(sol_name.str().c_str());
483  sol_ofs.precision(8);
484  x.Save(sol_ofs);
485  }
486 
487  // 17. Send the solution by socket to a GLVis server.
488  if (visualization)
489  {
490  char vishost[] = "localhost";
491  int visport = 19916;
492  socketstream sol_sock(vishost, visport);
493  sol_sock << "parallel " << num_procs << " " << myid << "\n";
494  sol_sock.precision(8);
495  sol_sock << "solution\n" << *pmesh << x << flush;
496  }
497 
498  // 18. Free the used memory.
499  delete a;
500  delete a_hpc;
501  if (a_oper != &A) { delete a_oper; }
502  delete a_pc;
503  delete b;
504  delete fespace;
505  delete fespace_lor;
506  delete fec_lor;
507  delete pmesh_lor;
508  if (order > 0) { delete fec; }
509  delete pmesh;
510  delete pcg;
511 
512  MPI_Finalize();
513 
514  return 0;
515 }
int Size() const
Logical size of the array.
Definition: array.hpp:110
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:529
Conjugate gradient method.
Definition: solvers.hpp:111
int GetNumIterations() const
Definition: solvers.hpp:66
Subclass constant coefficient.
Definition: coefficient.hpp:57
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:295
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to &#39;master&#39;.
Definition: hypre.cpp:770
const Geometry::Type geom
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:51
void AssembleBilinearForm(BilinearForm &a) const
const int rdim
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:614
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:361
StopWatch tic_toc
Definition: tic_toc.cpp:447
double RealTime()
Definition: tic_toc.cpp:426
Abstract parallel finite element space.
Definition: pfespace.hpp:31
H1_FiniteElement< geom, sol_p > sol_fe_t
TConstantCoefficient coeff_t
The BoomerAMG solver in hypre.
Definition: hypre.hpp:783
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_Int GetGlobalNumRows() const
Definition: hypre.hpp:374
int dim
Definition: ex3.cpp:47
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:72
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6718
const int sol_p
void SetMaxIter(int max_it)
Definition: solvers.hpp:63
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3315
int MeshGenerator()
Get the mesh generator/type.
Definition: mesh.hpp:607
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:108
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:178
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:3399
int Dimension() const
Definition: mesh.hpp:641
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
int main()
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:172
H1_FiniteElementSpace< mesh_fe_t > mesh_fes_t
void SetRelTol(double rtol)
Definition: solvers.hpp:61
void UsePrecomputedSparsity(int ps=1)
H1_FiniteElementSpace< sol_fe_t > sol_fes_t
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 mfem_error(const char *msg)
Definition: error.cpp:107
H1_FiniteElement< geom, mesh_p > mesh_fe_t
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
Form a constrained linear system using a matrix-free approach.
Definition: operator.cpp:21
TMesh< mesh_fes_t > mesh_t
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear sys...
Definition: operator.cpp:57
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
TBilinearForm< mesh_t, sol_fes_t, int_rule_t, integ_t > HPCBilinearForm
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:174
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
virtual const char * Name() const
Definition: fe_coll.hpp:117
const int mesh_p
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:607
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:125
Vector data type.
Definition: vector.hpp:41
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5401
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:93
TIntegrationRule< geom, ir_order > int_rule_t
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:146
TIntegrator< coeff_t, TDiffusionKernel > integ_t
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Abstract operator.
Definition: operator.hpp:21
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:175
Class for parallel meshes.
Definition: pmesh.hpp:29
void EnableStaticCondensation()
const int ir_order
bool Good() const
Definition: optparser.hpp:120