MFEM  v4.5.1
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 and HYPRE.
59  Mpi::Init(argc, argv);
60  int num_procs = Mpi::WorldSize();
61  int myid = Mpi::WorldRank();
62  Hypre::Init();
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  return 1;
98  }
99  if (myid == 0)
100  {
101  args.PrintOptions(cout);
102  }
103  kappa = freq * M_PI;
104 
105  // 3. Enable hardware devices such as GPUs, and programming models such as
106  // CUDA, OCCA, RAJA and OpenMP based on command line options.
107  Device device(device_config);
108  if (myid == 0) { device.Print(); }
109 
110  // 4. Read the (serial) mesh from the given mesh file on all processors. We
111  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
112  // and volume meshes with the same code.
113  Mesh *mesh = new Mesh(mesh_file, 1, 1);
114  dim = mesh->Dimension();
115  int sdim = mesh->SpaceDimension();
116 
117  // 5. Refine the serial mesh on all processors to increase the resolution. In
118  // this example we do 'ref_levels' of uniform refinement. We choose
119  // 'ref_levels' to be the largest number that gives a final mesh with no
120  // more than 1,000 elements.
121  {
122  int ref_levels = (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
123  for (int l = 0; l < ref_levels; l++)
124  {
125  mesh->UniformRefinement();
126  }
127  }
128 
129  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
130  // this mesh further in parallel to increase the resolution. Once the
131  // parallel mesh is defined, the serial mesh can be deleted. Tetrahedral
132  // meshes need to be reoriented before we can define high-order Nedelec
133  // spaces on them.
134  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
135  delete mesh;
136  {
137  int par_ref_levels = 1;
138  for (int l = 0; l < par_ref_levels; l++)
139  {
140  pmesh->UniformRefinement();
141  }
142  }
143 
144  // 7. Define a parallel finite element space on the parallel mesh. Here we
145  // use Nedelec or Raviart-Thomas finite elements of the specified order.
146  FiniteElementCollection *trial_fec = NULL;
147  FiniteElementCollection *test_fec = NULL;
148 
149  if (prob == 0)
150  {
151  trial_fec = new H1_FECollection(order, dim);
152  test_fec = new ND_FECollection(order, dim);
153  }
154  else if (prob == 1)
155  {
156  trial_fec = new ND_FECollection(order, dim);
157  test_fec = new RT_FECollection(order-1, dim);
158  }
159  else
160  {
161  trial_fec = new RT_FECollection(order-1, dim);
162  test_fec = new L2_FECollection(order-1, dim);
163  }
164 
165  ParFiniteElementSpace trial_fes(pmesh, trial_fec);
166  ParFiniteElementSpace test_fes(pmesh, test_fec);
167 
168  HYPRE_BigInt trial_size = trial_fes.GlobalTrueVSize();
169  HYPRE_BigInt test_size = test_fes.GlobalTrueVSize();
170 
171  if (myid == 0)
172  {
173  if (prob == 0)
174  {
175  cout << "Number of Nedelec finite element unknowns: " << test_size << endl;
176  cout << "Number of H1 finite element unknowns: " << trial_size << endl;
177  }
178  else if (prob == 1)
179  {
180  cout << "Number of Nedelec finite element unknowns: " << trial_size << endl;
181  cout << "Number of Raviart-Thomas finite element unknowns: " << test_size <<
182  endl;
183  }
184  else
185  {
186  cout << "Number of Raviart-Thomas finite element unknowns: "
187  << trial_size << endl;
188  cout << "Number of L2 finite element unknowns: " << test_size << endl;
189  }
190  }
191 
192  // 8. Define the solution vector as a parallel finite element grid function
193  // corresponding to the trial fespace.
194  ParGridFunction gftest(&test_fes);
195  ParGridFunction gftrial(&trial_fes);
196  ParGridFunction x(&test_fes);
198  VectorFunctionCoefficient gradp_coef(sdim, gradp_exact);
199  VectorFunctionCoefficient v_coef(sdim, v_exact);
200  VectorFunctionCoefficient curlv_coef(sdim, curlv_exact);
201  FunctionCoefficient divgradp_coef(div_gradp_exact);
202 
203  if (prob == 0)
204  {
205  gftrial.ProjectCoefficient(p_coef);
206  }
207  else if (prob == 1)
208  {
209  gftrial.ProjectCoefficient(v_coef);
210  }
211  else
212  {
213  gftrial.ProjectCoefficient(gradp_coef);
214  }
215 
216  gftrial.SetTrueVector();
217  gftrial.SetFromTrueVector();
218 
219  // 9. Set up the parallel bilinear forms for L2 projection.
220  ConstantCoefficient one(1.0);
221  ParBilinearForm a(&test_fes);
222  ParMixedBilinearForm a_mixed(&trial_fes, &test_fes);
223  if (pa)
224  {
225  a.SetAssemblyLevel(AssemblyLevel::PARTIAL);
226  a_mixed.SetAssemblyLevel(AssemblyLevel::PARTIAL);
227  }
228 
229  if (prob == 0)
230  {
233  }
234  else if (prob == 1)
235  {
238  }
239  else
240  {
243  }
244 
245  // 10. Assemble the parallel bilinear form and the corresponding linear
246  // system, applying any necessary transformations such as: parallel
247  // assembly, eliminating boundary conditions, applying conforming
248  // constraints for non-conforming AMR, static condensation, etc.
249  if (static_cond) { a.EnableStaticCondensation(); }
250 
251  a.Assemble();
252  if (!pa) { a.Finalize(); }
253 
254  a_mixed.Assemble();
255  if (!pa) { a_mixed.Finalize(); }
256 
257  Vector B(test_fes.GetTrueVSize());
258  Vector X(test_fes.GetTrueVSize());
259 
260  if (pa)
261  {
262  ParLinearForm b(&test_fes); // used as a vector
263  a_mixed.Mult(gftrial, b); // process-local multiplication
264  b.ParallelAssemble(B);
265  }
266  else
267  {
268  HypreParMatrix *mixed = a_mixed.ParallelAssemble();
269 
270  Vector P(trial_fes.GetTrueVSize());
271  gftrial.GetTrueDofs(P);
272 
273  mixed->Mult(P,B);
274 
275  delete mixed;
276  }
277 
278  // 11. Define and apply a parallel PCG solver for AX=B with Jacobi
279  // preconditioner.
280  if (pa)
281  {
282  Array<int> ess_tdof_list; // empty
283 
284  OperatorPtr A;
285  a.FormSystemMatrix(ess_tdof_list, A);
286 
287  OperatorJacobiSmoother Jacobi(a, ess_tdof_list);
288 
289  CGSolver cg(MPI_COMM_WORLD);
290  cg.SetRelTol(1e-12);
291  cg.SetMaxIter(1000);
292  cg.SetPrintLevel(1);
293  cg.SetOperator(*A);
294  cg.SetPreconditioner(Jacobi);
295  X = 0.0;
296  cg.Mult(B, X);
297  }
298  else
299  {
300  HypreParMatrix *Amat = a.ParallelAssemble();
301  HypreDiagScale Jacobi(*Amat);
302  HyprePCG pcg(*Amat);
303  pcg.SetTol(1e-12);
304  pcg.SetMaxIter(1000);
305  pcg.SetPrintLevel(2);
306  pcg.SetPreconditioner(Jacobi);
307  X = 0.0;
308  pcg.Mult(B, X);
309 
310  delete Amat;
311  }
312 
313  x.SetFromTrueDofs(X);
314 
315  // 12. Compute the same field by applying a DiscreteInterpolator.
316  ParGridFunction discreteInterpolant(&test_fes);
317  ParDiscreteLinearOperator dlo(&trial_fes, &test_fes);
318  if (prob == 0)
319  {
321  }
322  else if (prob == 1)
323  {
325  }
326  else
327  {
329  }
330 
331  dlo.Assemble();
332  dlo.Mult(gftrial, discreteInterpolant);
333 
334  // 13. Compute the projection of the exact field.
335  ParGridFunction exact_proj(&test_fes);
336  if (prob == 0)
337  {
338  exact_proj.ProjectCoefficient(gradp_coef);
339  }
340  else if (prob == 1)
341  {
342  exact_proj.ProjectCoefficient(curlv_coef);
343  }
344  else
345  {
346  exact_proj.ProjectCoefficient(divgradp_coef);
347  }
348 
349  exact_proj.SetTrueVector();
350  exact_proj.SetFromTrueVector();
351 
352  // 14. Compute and print the L_2 norm of the error.
353  if (prob == 0)
354  {
355  double errSol = x.ComputeL2Error(gradp_coef);
356  double errInterp = discreteInterpolant.ComputeL2Error(gradp_coef);
357  double errProj = exact_proj.ComputeL2Error(gradp_coef);
358 
359  if (myid == 0)
360  {
361  cout << "\n Solution of (E_h,v) = (grad p_h,v) for E_h and v in H(curl)"
362  ": || E_h - grad p ||_{L_2} = " << errSol << '\n' << endl;
363  cout << " Gradient interpolant E_h = grad p_h in H(curl): || E_h - grad"
364  " p ||_{L_2} = " << errInterp << '\n' << endl;
365  cout << " Projection E_h of exact grad p in H(curl): || E_h - grad p "
366  "||_{L_2} = " << errProj << '\n' << endl;
367  }
368  }
369  else if (prob == 1)
370  {
371  double errSol = x.ComputeL2Error(curlv_coef);
372  double errInterp = discreteInterpolant.ComputeL2Error(curlv_coef);
373  double errProj = exact_proj.ComputeL2Error(curlv_coef);
374 
375  if (myid == 0)
376  {
377  cout << "\n Solution of (E_h,w) = (curl v_h,w) for E_h and w in "
378  "H(div): || E_h - curl v ||_{L_2} = " << errSol << '\n' << endl;
379  cout << " Curl interpolant E_h = curl v_h in H(div): || E_h - curl v "
380  "||_{L_2} = " << errInterp << '\n' << endl;
381  cout << " Projection E_h of exact curl v in H(div): || E_h - curl v "
382  "||_{L_2} = " << errProj << '\n' << endl;
383  }
384  }
385  else
386  {
387  int order_quad = max(2, 2*order+1);
388  const IntegrationRule *irs[Geometry::NumGeom];
389  for (int i=0; i < Geometry::NumGeom; ++i)
390  {
391  irs[i] = &(IntRules.Get(i, order_quad));
392  }
393 
394  double errSol = x.ComputeL2Error(divgradp_coef, irs);
395  double errInterp = discreteInterpolant.ComputeL2Error(divgradp_coef, irs);
396  double errProj = exact_proj.ComputeL2Error(divgradp_coef, irs);
397 
398  if (myid == 0)
399  {
400  cout << "\n Solution of (f_h,q) = (div v_h,q) for f_h and q in L_2: "
401  "|| f_h - div v ||_{L_2} = " << errSol << '\n' << endl;
402  cout << " Divergence interpolant f_h = div v_h in L_2: || f_h - div v"
403  " ||_{L_2} = " << errInterp << '\n' << endl;
404  cout << " Projection f_h of exact div v in L_2: || f_h - div v "
405  "||_{L_2} = " << errProj << '\n' << endl;
406  }
407  }
408 
409  // 15. Save the refined mesh and the solution in parallel. This output can
410  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
411  {
412  ostringstream mesh_name, sol_name;
413  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
414  sol_name << "sol." << setfill('0') << setw(6) << myid;
415 
416  ofstream mesh_ofs(mesh_name.str().c_str());
417  mesh_ofs.precision(8);
418  pmesh->Print(mesh_ofs);
419 
420  ofstream sol_ofs(sol_name.str().c_str());
421  sol_ofs.precision(8);
422  x.Save(sol_ofs);
423  }
424 
425  // 16. Send the solution by socket to a GLVis server.
426  if (visualization)
427  {
428  char vishost[] = "localhost";
429  int visport = 19916;
430  socketstream sol_sock(vishost, visport);
431  sol_sock << "parallel " << num_procs << " " << myid << "\n";
432  sol_sock.precision(8);
433  sol_sock << "solution\n" << *pmesh << x << flush;
434  }
435 
436  // 17. Free the used memory.
437  delete trial_fec;
438  delete test_fec;
439  delete pmesh;
440 
441  return 0;
442 }
443 
444 double p_exact(const Vector &x)
445 {
446  if (dim == 3)
447  {
448  return sin(x(0)) * sin(x(1)) * sin(x(2));
449  }
450  else if (dim == 2)
451  {
452  return sin(x(0)) * sin(x(1));
453  }
454 
455  return 0.0;
456 }
457 
458 void gradp_exact(const Vector &x, Vector &f)
459 {
460  if (dim == 3)
461  {
462  f(0) = cos(x(0)) * sin(x(1)) * sin(x(2));
463  f(1) = sin(x(0)) * cos(x(1)) * sin(x(2));
464  f(2) = sin(x(0)) * sin(x(1)) * cos(x(2));
465  }
466  else
467  {
468  f(0) = cos(x(0)) * sin(x(1));
469  f(1) = sin(x(0)) * cos(x(1));
470  if (x.Size() == 3) { f(2) = 0.0; }
471  }
472 }
473 
474 double div_gradp_exact(const Vector &x)
475 {
476  if (dim == 3)
477  {
478  return -3.0 * sin(x(0)) * sin(x(1)) * sin(x(2));
479  }
480  else if (dim == 2)
481  {
482  return -2.0 * sin(x(0)) * sin(x(1));
483  }
484 
485  return 0.0;
486 }
487 
488 void v_exact(const Vector &x, Vector &v)
489 {
490  if (dim == 3)
491  {
492  v(0) = sin(kappa * x(1));
493  v(1) = sin(kappa * x(2));
494  v(2) = sin(kappa * x(0));
495  }
496  else
497  {
498  v(0) = sin(kappa * x(1));
499  v(1) = sin(kappa * x(0));
500  if (x.Size() == 3) { v(2) = 0.0; }
501  }
502 }
503 
504 void curlv_exact(const Vector &x, Vector &cv)
505 {
506  if (dim == 3)
507  {
508  cv(0) = -kappa * cos(kappa * x(2));
509  cv(1) = -kappa * cos(kappa * x(0));
510  cv(2) = -kappa * cos(kappa * x(1));
511  }
512  else
513  {
514  cv = 0.0;
515  }
516 }
void SetTol(double tol)
Definition: hypre.cpp:3995
Conjugate gradient method.
Definition: solvers.hpp:465
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:923
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:152
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:711
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:285
int Size() const
Returns the size of the vector.
Definition: vector.hpp:200
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:923
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:873
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:525
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:4015
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:289
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:146
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: pgridfunc.hpp:169
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
void curlv_exact(const Vector &x, Vector &cv)
Definition: ex24.cpp:439
constexpr char vishost[]
Jacobi preconditioner in hypre.
Definition: hypre.hpp:1410
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:274
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9816
constexpr int visport
void SetMaxIter(int max_it)
Definition: solvers.hpp:200
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Definition: pgridfunc.hpp:281
void Assemble(int skip_zeros=1)
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void v_exact(const Vector &x, Vector &v)
Definition: ex24.cpp:423
int Dimension() const
Definition: mesh.hpp:1006
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
prob_type prob
Definition: ex25.cpp:153
void SetMaxIter(int max_iter)
Definition: hypre.cpp:4005
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:341
double p_exact(const Vector &x)
Definition: ex24.cpp:379
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:1007
void SetRelTol(double rtol)
Definition: solvers.hpp:198
PCG solver in hypre.
Definition: hypre.hpp:1204
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 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:4020
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:324
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:479
A general function coefficient.
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:415
Vector data type.
Definition: vector.hpp:60
double div_gradp_exact(const Vector &x)
Definition: ex24.cpp:409
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:1732
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:173
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:220
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4770
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:343
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:4043
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:379
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:95
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:288
double f(const Vector &p)
void gradp_exact(const Vector &, Vector &)
Definition: ex24.cpp:393
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150