MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex24p.cpp
Go to the documentation of this file.
1 // MFEM Example 24 - Parallel Version
2 //
3 // Compile with: make ex24p
4 //
5 // Sample runs: mpirun -np 4 ex24p -m ../data/star.mesh
6 // mpirun -np 4 ex24p -m ../data/square-disc.mesh -o 2
7 // mpirun -np 4 ex24p -m ../data/beam-tet.mesh
8 // mpirun -np 4 ex24p -m ../data/beam-hex.mesh -o 2 -pa
9 // mpirun -np 4 ex24p -m ../data/beam-hex.mesh -o 2 -pa -p 1
10 // mpirun -np 4 ex24p -m ../data/beam-hex.mesh -o 2 -pa -p 2
11 // mpirun -np 4 ex24p -m ../data/escher.mesh
12 // mpirun -np 4 ex24p -m ../data/escher.mesh -o 2
13 // mpirun -np 4 ex24p -m ../data/fichera.mesh
14 // mpirun -np 4 ex24p -m ../data/fichera-q2.vtk
15 // mpirun -np 4 ex24p -m ../data/fichera-q3.mesh
16 // mpirun -np 4 ex24p -m ../data/square-disc-nurbs.mesh
17 // mpirun -np 4 ex24p -m ../data/beam-hex-nurbs.mesh
18 // mpirun -np 4 ex24p -m ../data/amr-quad.mesh -o 2
19 // mpirun -np 4 ex24p -m ../data/amr-hex.mesh
20 //
21 // Device sample runs:
22 // mpirun -np 4 ex24p -m ../data/star.mesh -pa -d cuda
23 // mpirun -np 4 ex24p -m ../data/star.mesh -pa -d raja-cuda
24 // mpirun -np 4 ex24p -m ../data/star.mesh -pa -d raja-omp
25 // mpirun -np 4 ex24p -m ../data/beam-hex.mesh -pa -d cuda
26 //
27 // Description: This example code illustrates usage of mixed finite element
28 // spaces, with three variants:
29 //
30 // 1) (grad p, u) for p in H^1 tested against u in H(curl)
31 // 2) (curl v, u) for v in H(curl) tested against u in H(div), 3D
32 // 3) (div v, q) for v in H(div) tested against q in L_2
33 //
34 // Using different approaches, we project the gradient, curl, or
35 // divergence to the appropriate space.
36 //
37 // We recommend viewing examples 1, 3, and 5 before viewing this
38 // example.
39 
40 #include "mfem.hpp"
41 #include <fstream>
42 #include <iostream>
43 
44 using namespace std;
45 using namespace mfem;
46 
47 double p_exact(const Vector &x);
48 void gradp_exact(const Vector &, Vector &);
49 double div_gradp_exact(const Vector &x);
50 void v_exact(const Vector &x, Vector &v);
51 void curlv_exact(const Vector &x, Vector &cv);
52 
53 int dim;
54 double freq = 1.0, kappa;
55 
56 int main(int argc, char *argv[])
57 {
58  // 1. Initialize MPI.
59  int num_procs, myid;
60  MPI_Init(&argc, &argv);
61  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
62  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
63 
64  // 2. Parse command-line options.
65  const char *mesh_file = "../data/beam-hex.mesh";
66  int order = 1;
67  int prob = 0;
68  bool static_cond = false;
69  bool pa = false;
70  const char *device_config = "cpu";
71  bool visualization = 1;
72 
73  OptionsParser args(argc, argv);
74  args.AddOption(&mesh_file, "-m", "--mesh",
75  "Mesh file to use.");
76  args.AddOption(&order, "-o", "--order",
77  "Finite element order (polynomial degree).");
78  args.AddOption(&prob, "-p", "--problem-type",
79  "Choose between 0: grad, 1: curl, 2: div");
80  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
81  "--no-static-condensation", "Enable static condensation.");
82  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
83  "--no-partial-assembly", "Enable Partial Assembly.");
84  args.AddOption(&device_config, "-d", "--device",
85  "Device configuration string, see Device::Configure().");
86  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
87  "--no-visualization",
88  "Enable or disable GLVis visualization.");
89 
90  args.Parse();
91  if (!args.Good())
92  {
93  if (myid == 0)
94  {
95  args.PrintUsage(cout);
96  }
97  MPI_Finalize();
98  return 1;
99  }
100  if (myid == 0)
101  {
102  args.PrintOptions(cout);
103  }
104  kappa = freq * M_PI;
105 
106  // 3. Enable hardware devices such as GPUs, and programming models such as
107  // CUDA, OCCA, RAJA and OpenMP based on command line options.
108  Device device(device_config);
109  if (myid == 0) { device.Print(); }
110 
111  // 4. Read the (serial) mesh from the given mesh file on all processors. We
112  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
113  // and volume meshes with the same code.
114  Mesh *mesh = new Mesh(mesh_file, 1, 1);
115  dim = mesh->Dimension();
116  int sdim = mesh->SpaceDimension();
117 
118  // 5. Refine the serial mesh on all processors to increase the resolution. In
119  // this example we do 'ref_levels' of uniform refinement. We choose
120  // 'ref_levels' to be the largest number that gives a final mesh with no
121  // more than 1,000 elements.
122  {
123  int ref_levels = (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
124  for (int l = 0; l < ref_levels; l++)
125  {
126  mesh->UniformRefinement();
127  }
128  }
129 
130  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
131  // this mesh further in parallel to increase the resolution. Once the
132  // parallel mesh is defined, the serial mesh can be deleted. Tetrahedral
133  // meshes need to be reoriented before we can define high-order Nedelec
134  // spaces on them.
135  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
136  delete mesh;
137  {
138  int par_ref_levels = 1;
139  for (int l = 0; l < par_ref_levels; l++)
140  {
141  pmesh->UniformRefinement();
142  }
143  }
144  pmesh->ReorientTetMesh();
145 
146  // 7. Define a parallel finite element space on the parallel mesh. Here we
147  // use Nedelec or Raviart-Thomas finite elements of the specified order.
148  FiniteElementCollection *trial_fec = NULL;
149  FiniteElementCollection *test_fec = NULL;
150 
151  if (prob == 0)
152  {
153  trial_fec = new H1_FECollection(order, dim);
154  test_fec = new ND_FECollection(order, dim);
155  }
156  else if (prob == 1)
157  {
158  trial_fec = new ND_FECollection(order, dim);
159  test_fec = new RT_FECollection(order-1, dim);
160  }
161  else
162  {
163  trial_fec = new RT_FECollection(order-1, dim);
164  test_fec = new L2_FECollection(order-1, dim);
165  }
166 
167  ParFiniteElementSpace trial_fes(pmesh, trial_fec);
168  ParFiniteElementSpace test_fes(pmesh, test_fec);
169 
170  HYPRE_BigInt trial_size = trial_fes.GlobalTrueVSize();
171  HYPRE_BigInt test_size = test_fes.GlobalTrueVSize();
172 
173  if (myid == 0)
174  {
175  if (prob == 0)
176  {
177  cout << "Number of Nedelec finite element unknowns: " << test_size << endl;
178  cout << "Number of H1 finite element unknowns: " << trial_size << endl;
179  }
180  else if (prob == 1)
181  {
182  cout << "Number of Nedelec finite element unknowns: " << trial_size << endl;
183  cout << "Number of Raviart-Thomas finite element unknowns: " << test_size <<
184  endl;
185  }
186  else
187  {
188  cout << "Number of Raviart-Thomas finite element unknowns: "
189  << trial_size << endl;
190  cout << "Number of L2 finite element unknowns: " << test_size << endl;
191  }
192  }
193 
194  // 8. Define the solution vector as a parallel finite element grid function
195  // corresponding to the trial fespace.
196  ParGridFunction gftest(&test_fes);
197  ParGridFunction gftrial(&trial_fes);
198  ParGridFunction x(&test_fes);
200  VectorFunctionCoefficient gradp_coef(sdim, gradp_exact);
201  VectorFunctionCoefficient v_coef(sdim, v_exact);
202  VectorFunctionCoefficient curlv_coef(sdim, curlv_exact);
203  FunctionCoefficient divgradp_coef(div_gradp_exact);
204 
205  if (prob == 0)
206  {
207  gftrial.ProjectCoefficient(p_coef);
208  }
209  else if (prob == 1)
210  {
211  gftrial.ProjectCoefficient(v_coef);
212  }
213  else
214  {
215  gftrial.ProjectCoefficient(gradp_coef);
216  }
217 
218  gftrial.SetTrueVector();
219  gftrial.SetFromTrueVector();
220 
221  // 9. Set up the parallel bilinear forms for L2 projection.
222  ConstantCoefficient one(1.0);
223  ParBilinearForm a(&test_fes);
224  ParMixedBilinearForm a_mixed(&trial_fes, &test_fes);
225  if (pa)
226  {
227  a.SetAssemblyLevel(AssemblyLevel::PARTIAL);
228  a_mixed.SetAssemblyLevel(AssemblyLevel::PARTIAL);
229  }
230 
231  if (prob == 0)
232  {
235  }
236  else if (prob == 1)
237  {
240  }
241  else
242  {
245  }
246 
247  // 10. Assemble the parallel bilinear form and the corresponding linear
248  // system, applying any necessary transformations such as: parallel
249  // assembly, eliminating boundary conditions, applying conforming
250  // constraints for non-conforming AMR, static condensation, etc.
251  if (static_cond) { a.EnableStaticCondensation(); }
252 
253  a.Assemble();
254  if (!pa) { a.Finalize(); }
255 
256  a_mixed.Assemble();
257  if (!pa) { a_mixed.Finalize(); }
258 
259  Vector B(test_fes.GetTrueVSize());
260  Vector X(test_fes.GetTrueVSize());
261 
262  if (pa)
263  {
264  ParLinearForm b(&test_fes); // used as a vector
265  a_mixed.Mult(gftrial, b); // process-local multiplication
266  b.ParallelAssemble(B);
267  }
268  else
269  {
270  HypreParMatrix *mixed = a_mixed.ParallelAssemble();
271 
272  Vector P(trial_fes.GetTrueVSize());
273  gftrial.GetTrueDofs(P);
274 
275  mixed->Mult(P,B);
276 
277  delete mixed;
278  }
279 
280  // 11. Define and apply a parallel PCG solver for AX=B with Jacobi
281  // preconditioner.
282  if (pa)
283  {
284  Array<int> ess_tdof_list; // empty
285 
286  OperatorPtr A;
287  a.FormSystemMatrix(ess_tdof_list, A);
288 
289  OperatorJacobiSmoother Jacobi(a, ess_tdof_list);
290 
291  CGSolver cg(MPI_COMM_WORLD);
292  cg.SetRelTol(1e-12);
293  cg.SetMaxIter(1000);
294  cg.SetPrintLevel(1);
295  cg.SetOperator(*A);
296  cg.SetPreconditioner(Jacobi);
297  X = 0.0;
298  cg.Mult(B, X);
299  }
300  else
301  {
302  HypreParMatrix *Amat = a.ParallelAssemble();
303  HypreDiagScale Jacobi(*Amat);
304  HyprePCG pcg(*Amat);
305  pcg.SetTol(1e-12);
306  pcg.SetMaxIter(1000);
307  pcg.SetPrintLevel(2);
308  pcg.SetPreconditioner(Jacobi);
309  X = 0.0;
310  pcg.Mult(B, X);
311 
312  delete Amat;
313  }
314 
315  x.SetFromTrueDofs(X);
316 
317  // 12. Compute the same field by applying a DiscreteInterpolator.
318  ParGridFunction discreteInterpolant(&test_fes);
319  ParDiscreteLinearOperator dlo(&trial_fes, &test_fes);
320  if (prob == 0)
321  {
323  }
324  else if (prob == 1)
325  {
327  }
328  else
329  {
331  }
332 
333  dlo.Assemble();
334  dlo.Mult(gftrial, discreteInterpolant);
335 
336  // 13. Compute the projection of the exact field.
337  ParGridFunction exact_proj(&test_fes);
338  if (prob == 0)
339  {
340  exact_proj.ProjectCoefficient(gradp_coef);
341  }
342  else if (prob == 1)
343  {
344  exact_proj.ProjectCoefficient(curlv_coef);
345  }
346  else
347  {
348  exact_proj.ProjectCoefficient(divgradp_coef);
349  }
350 
351  exact_proj.SetTrueVector();
352  exact_proj.SetFromTrueVector();
353 
354  // 14. Compute and print the L_2 norm of the error.
355  if (prob == 0)
356  {
357  double errSol = x.ComputeL2Error(gradp_coef);
358  double errInterp = discreteInterpolant.ComputeL2Error(gradp_coef);
359  double errProj = exact_proj.ComputeL2Error(gradp_coef);
360 
361  if (myid == 0)
362  {
363  cout << "\n Solution of (E_h,v) = (grad p_h,v) for E_h and v in H(curl)"
364  ": || E_h - grad p ||_{L_2} = " << errSol << '\n' << endl;
365  cout << " Gradient interpolant E_h = grad p_h in H(curl): || E_h - grad"
366  " p ||_{L_2} = " << errInterp << '\n' << endl;
367  cout << " Projection E_h of exact grad p in H(curl): || E_h - grad p "
368  "||_{L_2} = " << errProj << '\n' << endl;
369  }
370  }
371  else if (prob == 1)
372  {
373  double errSol = x.ComputeL2Error(curlv_coef);
374  double errInterp = discreteInterpolant.ComputeL2Error(curlv_coef);
375  double errProj = exact_proj.ComputeL2Error(curlv_coef);
376 
377  if (myid == 0)
378  {
379  cout << "\n Solution of (E_h,w) = (curl v_h,w) for E_h and w in "
380  "H(div): || E_h - curl v ||_{L_2} = " << errSol << '\n' << endl;
381  cout << " Curl interpolant E_h = curl v_h in H(div): || E_h - curl v "
382  "||_{L_2} = " << errInterp << '\n' << endl;
383  cout << " Projection E_h of exact curl v in H(div): || E_h - curl v "
384  "||_{L_2} = " << errProj << '\n' << endl;
385  }
386  }
387  else
388  {
389  int order_quad = max(2, 2*order+1);
390  const IntegrationRule *irs[Geometry::NumGeom];
391  for (int i=0; i < Geometry::NumGeom; ++i)
392  {
393  irs[i] = &(IntRules.Get(i, order_quad));
394  }
395 
396  double errSol = x.ComputeL2Error(divgradp_coef, irs);
397  double errInterp = discreteInterpolant.ComputeL2Error(divgradp_coef, irs);
398  double errProj = exact_proj.ComputeL2Error(divgradp_coef, irs);
399 
400  if (myid == 0)
401  {
402  cout << "\n Solution of (f_h,q) = (div v_h,q) for f_h and q in L_2: "
403  "|| f_h - div v ||_{L_2} = " << errSol << '\n' << endl;
404  cout << " Divergence interpolant f_h = div v_h in L_2: || f_h - div v"
405  " ||_{L_2} = " << errInterp << '\n' << endl;
406  cout << " Projection f_h of exact div v in L_2: || f_h - div v "
407  "||_{L_2} = " << errProj << '\n' << endl;
408  }
409  }
410 
411  // 15. Save the refined mesh and the solution in parallel. This output can
412  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
413  {
414  ostringstream mesh_name, sol_name;
415  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
416  sol_name << "sol." << setfill('0') << setw(6) << myid;
417 
418  ofstream mesh_ofs(mesh_name.str().c_str());
419  mesh_ofs.precision(8);
420  pmesh->Print(mesh_ofs);
421 
422  ofstream sol_ofs(sol_name.str().c_str());
423  sol_ofs.precision(8);
424  x.Save(sol_ofs);
425  }
426 
427  // 16. Send the solution by socket to a GLVis server.
428  if (visualization)
429  {
430  char vishost[] = "localhost";
431  int visport = 19916;
432  socketstream sol_sock(vishost, visport);
433  sol_sock << "parallel " << num_procs << " " << myid << "\n";
434  sol_sock.precision(8);
435  sol_sock << "solution\n" << *pmesh << x << flush;
436  }
437 
438  // 17. Free the used memory.
439  delete trial_fec;
440  delete test_fec;
441  delete pmesh;
442 
443  MPI_Finalize();
444 
445  return 0;
446 }
447 
448 double p_exact(const Vector &x)
449 {
450  if (dim == 3)
451  {
452  return sin(x(0)) * sin(x(1)) * sin(x(2));
453  }
454  else if (dim == 2)
455  {
456  return sin(x(0)) * sin(x(1));
457  }
458 
459  return 0.0;
460 }
461 
462 void gradp_exact(const Vector &x, Vector &f)
463 {
464  if (dim == 3)
465  {
466  f(0) = cos(x(0)) * sin(x(1)) * sin(x(2));
467  f(1) = sin(x(0)) * cos(x(1)) * sin(x(2));
468  f(2) = sin(x(0)) * sin(x(1)) * cos(x(2));
469  }
470  else
471  {
472  f(0) = cos(x(0)) * sin(x(1));
473  f(1) = sin(x(0)) * cos(x(1));
474  if (x.Size() == 3) { f(2) = 0.0; }
475  }
476 }
477 
478 double div_gradp_exact(const Vector &x)
479 {
480  if (dim == 3)
481  {
482  return -3.0 * sin(x(0)) * sin(x(1)) * sin(x(2));
483  }
484  else if (dim == 2)
485  {
486  return -2.0 * sin(x(0)) * sin(x(1));
487  }
488 
489  return 0.0;
490 }
491 
492 void v_exact(const Vector &x, Vector &v)
493 {
494  if (dim == 3)
495  {
496  v(0) = sin(kappa * x(1));
497  v(1) = sin(kappa * x(2));
498  v(2) = sin(kappa * x(0));
499  }
500  else
501  {
502  v(0) = sin(kappa * x(1));
503  v(1) = sin(kappa * x(0));
504  if (x.Size() == 3) { v(2) = 0.0; }
505  }
506 }
507 
508 void curlv_exact(const Vector &x, Vector &cv)
509 {
510  if (dim == 3)
511  {
512  cv(0) = -kappa * cos(kappa * x(2));
513  cv(1) = -kappa * cos(kappa * x(0));
514  cv(2) = -kappa * cos(kappa * x(1));
515  }
516  else
517  {
518  cv = 0.0;
519  }
520 }
void SetTol(double tol)
Definition: hypre.cpp:3569
Conjugate gradient method.
Definition: solvers.hpp:316
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:920
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:143
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:616
virtual void ReorientTetMesh()
See the remarks for the serial version in mesh.hpp.
Definition: pmesh.cpp:2798
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:275
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:841
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:279
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
double kappa
Definition: ex24.cpp:54
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:493
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3589
void AddDomainInterpolator(DiscreteInterpolator *di)
Adds a domain interpolator. Assumes ownership of di.
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:279
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.
Class for parallel linear form.
Definition: plinearform.hpp:26
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:137
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: pgridfunc.hpp:164
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
void curlv_exact(const Vector &x, Vector &cv)
Definition: ex24.cpp:440
constexpr char vishost[]
Jacobi preconditioner in hypre.
Definition: hypre.hpp:1234
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:130
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
constexpr int visport
void SetMaxIter(int max_it)
Definition: solvers.hpp:100
void Assemble(int skip_zeros=1)
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4382
void v_exact(const Vector &x, Vector &v)
Definition: ex24.cpp:424
int Dimension() const
Definition: mesh.hpp:911
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
prob_type prob
Definition: ex25.cpp:153
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3579
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:341
double p_exact(const Vector &x)
Definition: ex24.cpp:380
A general vector function coefficient.
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACY.
int SpaceDimension() const
Definition: mesh.hpp:912
void SetRelTol(double rtol)
Definition: solvers.hpp:98
PCG solver in hypre.
Definition: hypre.hpp:1060
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
HYPRE_Int HYPRE_BigInt
double a
Definition: lissajous.cpp:41
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
Class for parallel bilinear form using different test and trial FE spaces.
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:273
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
int dim
Definition: ex24.cpp:53
double freq
Definition: ex24.cpp:54
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:3594
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:143
Class for parallel bilinear form.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:330
A general function coefficient.
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:416
Vector data type.
Definition: vector.hpp:60
double div_gradp_exact(const Vector &x)
Definition: ex24.cpp:410
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0) const
Computes y = alpha * A * x + beta * y.
Definition: hypre.cpp:1529
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:92
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
Class for parallel grid function.
Definition: pgridfunc.hpp:32
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:121
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:277
virtual void Assemble(int skip_zeros=1)
Construct the internal matrix representation of the discrete linear operator.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:3617
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:377
virtual void Mult(const Vector &x, Vector &y) const
Matrix multiplication: .
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:87
int main()
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:285
double f(const Vector &p)
void gradp_exact(const Vector &, Vector &)
Definition: ex24.cpp:394
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150