MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
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  if (myid == 0)
148  {
149  cout << "\nMFEM SIMD width: " << MFEM_SIMD_BYTES/sizeof(double)
150  << " doubles\n" << endl;
151  }
152 
153  // See class BasisType in fem/fe_coll.hpp for available basis types
154  int basis = BasisType::GetType(basis_type[0]);
155  if (myid == 0)
156  {
157  cout << "Using " << BasisType::Name(basis) << " basis ..." << endl;
158  }
159 
160  // 3. Read the (serial) mesh from the given mesh file on all processors. We
161  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
162  // and volume meshes with the same code.
163  Mesh *mesh = new Mesh(mesh_file, 1, 1);
164  int dim = mesh->Dimension();
165 
166  // 4. Check if the optimized version matches the given mesh
167  if (perf)
168  {
169  if (myid == 0)
170  {
171  cout << "High-performance version using integration rule with "
172  << int_rule_t::qpts << " points ..." << endl;
173  }
174  if (!mesh_t::MatchesGeometry(*mesh))
175  {
176  if (myid == 0)
177  {
178  cout << "The given mesh does not match the optimized 'geom' parameter.\n"
179  << "Recompile with suitable 'geom' value." << endl;
180  }
181  delete mesh;
182  MPI_Finalize();
183  return 4;
184  }
185  else if (!mesh_t::MatchesNodes(*mesh))
186  {
187  if (myid == 0)
188  {
189  cout << "Switching the mesh curvature to match the "
190  << "optimized value (order " << mesh_p << ") ..." << endl;
191  }
192  mesh->SetCurvature(mesh_p, false, -1, Ordering::byNODES);
193  }
194  }
195 
196  // 5. Refine the serial mesh on all processors to increase the resolution. In
197  // this example we do 'ref_levels' of uniform refinement. We choose
198  // 'ref_levels' to be the largest number that gives a final mesh with no
199  // more than 10,000 elements, or as specified on the command line with the
200  // option '--refine-serial'.
201  {
202  int ref_levels = (ser_ref_levels != -1) ? ser_ref_levels :
203  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
204  for (int l = 0; l < ref_levels; l++)
205  {
206  mesh->UniformRefinement();
207  }
208  }
209  if (!perf && mesh->NURBSext)
210  {
211  const int new_mesh_p = std::min(sol_p, mesh_p);
212  if (myid == 0)
213  {
214  cout << "NURBS mesh: switching the mesh curvature to be "
215  << "min(sol_p, mesh_p) = " << new_mesh_p << " ..." << endl;
216  }
217  mesh->SetCurvature(new_mesh_p, false, -1, Ordering::byNODES);
218  }
219 
220  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
221  // this mesh further in parallel to increase the resolution. Once the
222  // parallel mesh is defined, the serial mesh can be deleted.
223  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
224  delete mesh;
225  {
226  for (int l = 0; l < par_ref_levels; l++)
227  {
228  pmesh->UniformRefinement();
229  }
230  }
231  if (pmesh->MeshGenerator() & 1) // simplex mesh
232  {
233  MFEM_VERIFY(pc_choice != LOR, "triangle and tet meshes do not support"
234  " the LOR preconditioner yet");
235  }
236 
237  // 7. Define a parallel finite element space on the parallel mesh. Here we
238  // use continuous Lagrange finite elements of the specified order. If
239  // order < 1, we instead use an isoparametric/isogeometric space.
241  if (order > 0)
242  {
243  fec = new H1_FECollection(order, dim, basis);
244  }
245  else if (pmesh->GetNodes())
246  {
247  fec = pmesh->GetNodes()->OwnFEC();
248  if (myid == 0)
249  {
250  cout << "Using isoparametric FEs: " << fec->Name() << endl;
251  }
252  }
253  else
254  {
255  fec = new H1_FECollection(order = 1, dim, basis);
256  }
257  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
258  HYPRE_Int size = fespace->GlobalTrueVSize();
259  if (myid == 0)
260  {
261  cout << "Number of finite element unknowns: " << size << endl;
262  }
263 
264  ParMesh *pmesh_lor = NULL;
265  FiniteElementCollection *fec_lor = NULL;
266  ParFiniteElementSpace *fespace_lor = NULL;
267  if (pc_choice == LOR)
268  {
269  int basis_lor = basis;
270  if (basis == BasisType::Positive) { basis_lor=BasisType::ClosedUniform; }
271  pmesh_lor = new ParMesh(pmesh, order, basis_lor);
272  fec_lor = new H1_FECollection(1, dim);
273  fespace_lor = new ParFiniteElementSpace(pmesh_lor, fec_lor);
274  }
275 
276  // 8. Check if the optimized version matches the given space
277  if (perf && !sol_fes_t::Matches(*fespace))
278  {
279  if (myid == 0)
280  {
281  cout << "The given order does not match the optimized parameter.\n"
282  << "Recompile with suitable 'sol_p' value." << endl;
283  }
284  delete fespace;
285  delete fec;
286  delete mesh;
287  MPI_Finalize();
288  return 5;
289  }
290 
291  // 9. Determine the list of true (i.e. parallel conforming) essential
292  // boundary dofs. In this example, the boundary conditions are defined
293  // by marking all the boundary attributes from the mesh as essential
294  // (Dirichlet) and converting them to a list of true dofs.
295  Array<int> ess_tdof_list;
296  if (pmesh->bdr_attributes.Size())
297  {
298  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
299  ess_bdr = 1;
300  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
301  }
302 
303  // 10. Set up the parallel linear form b(.) which corresponds to the
304  // right-hand side of the FEM linear system, which in this case is
305  // (1,phi_i) where phi_i are the basis functions in fespace.
306  ParLinearForm *b = new ParLinearForm(fespace);
307  ConstantCoefficient one(1.0);
309  b->Assemble();
310 
311  // 11. Define the solution vector x as a parallel finite element grid
312  // function corresponding to fespace. Initialize x with initial guess of
313  // zero, which satisfies the boundary conditions.
314  ParGridFunction x(fespace);
315  x = 0.0;
316 
317  // 12. Set up the parallel bilinear form a(.,.) on the finite element space
318  // that will hold the matrix corresponding to the Laplacian operator.
319  ParBilinearForm *a = new ParBilinearForm(fespace);
320  ParBilinearForm *a_pc = NULL;
321  if (pc_choice == LOR) { a_pc = new ParBilinearForm(fespace_lor); }
322  if (pc_choice == HO) { a_pc = new ParBilinearForm(fespace); }
323 
324  // 13. Assemble the parallel bilinear form and the corresponding linear
325  // system, applying any necessary transformations such as: parallel
326  // assembly, eliminating boundary conditions, applying conforming
327  // constraints for non-conforming AMR, static condensation, etc.
328  if (static_cond)
329  {
331  MFEM_VERIFY(pc_choice != LOR,
332  "cannot use LOR preconditioner with static condensation");
333  }
334 
335  if (myid == 0)
336  {
337  cout << "Assembling the matrix ..." << flush;
338  }
339  tic_toc.Clear();
340  tic_toc.Start();
341  // Pre-allocate sparsity assuming dense element matrices
343 
344  HPCBilinearForm *a_hpc = NULL;
345  Operator *a_oper = NULL;
346 
347  if (!perf)
348  {
349  // Standard assembly using a diffusion domain integrator
351  a->Assemble();
352  }
353  else
354  {
355  // High-performance assembly/evaluation using the templated operator type
356  a_hpc = new HPCBilinearForm(integ_t(coeff_t(1.0)), *fespace);
357  if (matrix_free)
358  {
359  a_hpc->Assemble(); // partial assembly
360  }
361  else
362  {
363  a_hpc->AssembleBilinearForm(*a); // full matrix assembly
364  }
365  }
366  tic_toc.Stop();
367  if (myid == 0)
368  {
369  cout << " done, " << tic_toc.RealTime() << "s." << endl;
370  }
371 
372  // 14. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
373  // preconditioner from hypre.
374 
375  // Setup the operator matrix (if applicable)
376  HypreParMatrix A;
377  Vector B, X;
378  if (perf && matrix_free)
379  {
380  a_hpc->FormLinearSystem(ess_tdof_list, x, *b, a_oper, X, B);
381  HYPRE_Int glob_size = fespace->GlobalTrueVSize();
382  if (myid == 0)
383  {
384  cout << "Size of linear system: " << glob_size << endl;
385  }
386  }
387  else
388  {
389  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
390  HYPRE_Int glob_size = A.GetGlobalNumRows();
391  if (myid == 0)
392  {
393  cout << "Size of linear system: " << glob_size << endl;
394  }
395  a_oper = &A;
396  }
397 
398  // Setup the matrix used for preconditioning
399  if (myid == 0)
400  {
401  cout << "Assembling the preconditioning matrix ..." << flush;
402  }
403  tic_toc.Clear();
404  tic_toc.Start();
405 
406  HypreParMatrix A_pc;
407  if (pc_choice == LOR)
408  {
409  // TODO: assemble the LOR matrix using the performance code
411  a_pc->UsePrecomputedSparsity();
412  a_pc->Assemble();
413  a_pc->FormSystemMatrix(ess_tdof_list, A_pc);
414  }
415  else if (pc_choice == HO)
416  {
417  if (!matrix_free)
418  {
419  A_pc.MakeRef(A); // matrix already assembled, reuse it
420  }
421  else
422  {
423  a_pc->UsePrecomputedSparsity();
424  a_hpc->AssembleBilinearForm(*a_pc);
425  a_pc->FormSystemMatrix(ess_tdof_list, A_pc);
426  }
427  }
428  tic_toc.Stop();
429  if (myid == 0)
430  {
431  cout << " done, " << tic_toc.RealTime() << "s." << endl;
432  }
433 
434  // Solve with CG or PCG, depending if the matrix A_pc is available
435  CGSolver *pcg;
436  pcg = new CGSolver(MPI_COMM_WORLD);
437  pcg->SetRelTol(1e-6);
438  pcg->SetMaxIter(500);
439  pcg->SetPrintLevel(1);
440 
441  HypreSolver *amg = NULL;
442 
443  pcg->SetOperator(*a_oper);
444  if (pc_choice != NONE)
445  {
446  amg = new HypreBoomerAMG(A_pc);
447  pcg->SetPreconditioner(*amg);
448  }
449 
450  tic_toc.Clear();
451  tic_toc.Start();
452 
453  pcg->Mult(B, X);
454 
455  tic_toc.Stop();
456  delete amg;
457 
458  if (myid == 0)
459  {
460  // Note: In the pcg algorithm, the number of operator Mult() calls is
461  // N_iter and the number of preconditioner Mult() calls is N_iter+1.
462  cout << "Time per CG step: "
463  << tic_toc.RealTime() / pcg->GetNumIterations() << "s." << endl;
464  }
465 
466  // 15. Recover the parallel grid function corresponding to X. This is the
467  // local finite element solution on each processor.
468  if (perf && matrix_free)
469  {
470  a_hpc->RecoverFEMSolution(X, *b, x);
471  }
472  else
473  {
474  a->RecoverFEMSolution(X, *b, x);
475  }
476 
477  // 16. Save the refined mesh and the solution in parallel. This output can
478  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
479  {
480  ostringstream mesh_name, sol_name;
481  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
482  sol_name << "sol." << setfill('0') << setw(6) << myid;
483 
484  ofstream mesh_ofs(mesh_name.str().c_str());
485  mesh_ofs.precision(8);
486  pmesh->Print(mesh_ofs);
487 
488  ofstream sol_ofs(sol_name.str().c_str());
489  sol_ofs.precision(8);
490  x.Save(sol_ofs);
491  }
492 
493  // 17. Send the solution by socket to a GLVis server.
494  if (visualization)
495  {
496  char vishost[] = "localhost";
497  int visport = 19916;
498  socketstream sol_sock(vishost, visport);
499  sol_sock << "parallel " << num_procs << " " << myid << "\n";
500  sol_sock.precision(8);
501  sol_sock << "solution\n" << *pmesh << x << flush;
502  }
503 
504  // 18. Free the used memory.
505  delete a;
506  delete a_hpc;
507  if (a_oper != &A) { delete a_oper; }
508  delete a_pc;
509  delete b;
510  delete fespace;
511  delete fespace_lor;
512  delete fec_lor;
513  delete pmesh_lor;
514  if (order > 0) { delete fec; }
515  delete pmesh;
516  delete pcg;
517 
518  MPI_Finalize();
519 
520  return 0;
521 }
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
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:775
Templated bilinear form class, cf. bilinearform.?pp.
Conjugate gradient method.
Definition: solvers.hpp:258
int GetNumIterations() const
Definition: solvers.hpp:101
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:535
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to &#39;master&#39;.
Definition: hypre.cpp:785
const Geometry::Type geom
Definition: ex1.cpp:40
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
const int rdim
Definition: ex1.cpp:43
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:819
StopWatch tic_toc
Definition: tic_toc.cpp:447
double RealTime()
Definition: tic_toc.cpp:426
Abstract parallel finite element space.
Definition: pfespace.hpp:28
int main(int argc, char *argv[])
Definition: ex1.cpp:66
The Integrator class combines a kernel and a coefficient.
Definition: tbilininteg.hpp:26
H1_FiniteElement< geom, sol_p > sol_fe_t
Definition: ex1.cpp:52
void Stop()
Stop the stopwatch.
Definition: tic_toc.cpp:416
TConstantCoefficient coeff_t
Definition: ex1.cpp:57
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1079
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:415
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
constexpr char vishost[]
void AssembleBilinearForm(BilinearForm &a) const
Assemble element matrices and add them to the bilinear form.
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:153
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:8382
constexpr int visport
const int sol_p
Definition: ex1.cpp:42
void SetMaxIter(int max_it)
Definition: solvers.hpp:98
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:4165
void Assemble()
Partial assembly of quadrature point data.
int MeshGenerator()
Get the mesh generator/type.
Definition: mesh.hpp:730
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
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:4169
int Dimension() const
Definition: mesh.hpp:788
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
void Start()
Clear the elapsed time and start the stopwatch.
Definition: tic_toc.cpp:411
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:201
H1_FiniteElementSpace< mesh_fe_t > mesh_fes_t
Definition: ex1.cpp:48
void SetRelTol(double rtol)
Definition: solvers.hpp:96
void UsePrecomputedSparsity(int ps=1)
For scalar FE spaces, precompute the sparsity pattern of the matrix (assuming dense element matrices)...
H1_FiniteElementSpace< sol_fe_t > sol_fes_t
Definition: ex1.cpp:53
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
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)
Add a boolean option and set &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
H1_FiniteElement< geom, mesh_p > mesh_fe_t
Definition: ex1.cpp:47
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:51
TMesh< mesh_fes_t > mesh_t
Definition: ex1.cpp:49
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:85
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
TBilinearForm< mesh_t, sol_fes_t, int_rule_t, integ_t > HPCBilinearForm
Definition: ex1.cpp:61
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:203
virtual 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:61
const int mesh_p
Definition: ex1.cpp:41
int dim
Definition: ex24.cpp:53
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:304
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:700
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:272
Vector data type.
Definition: vector.hpp:51
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:6603
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
TIntegrationRule< geom, ir_order > int_rule_t
Definition: ex1.cpp:56
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:159
TIntegrator< coeff_t, TDiffusionKernel > integ_t
Definition: ex1.cpp:58
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
Class for parallel meshes.
Definition: pmesh.hpp:32
void Clear()
Clear the elapsed time on the stopwatch and restart it if it&#39;s running.
Definition: tic_toc.cpp:406
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
const int ir_order
Definition: ex1.cpp:44
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145