MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex5p.cpp
Go to the documentation of this file.
1 // MFEM Example 5 - Parallel Version
2 //
3 // Compile with: make ex5p
4 //
5 // Sample runs: mpirun -np 4 ex5p -m ../data/square-disc.mesh
6 // mpirun -np 4 ex5p -m ../data/star.mesh
7 // mpirun -np 4 ex5p -m ../data/star.mesh -r 2 -pa
8 // mpirun -np 4 ex5p -m ../data/beam-tet.mesh
9 // mpirun -np 4 ex5p -m ../data/beam-hex.mesh
10 // mpirun -np 4 ex5p -m ../data/beam-hex.mesh -pa
11 // mpirun -np 4 ex5p -m ../data/escher.mesh
12 // mpirun -np 4 ex5p -m ../data/fichera.mesh
13 //
14 // Device sample runs:
15 // mpirun -np 4 ex5p -m ../data/star.mesh -r 2 -pa -d cuda
16 // mpirun -np 4 ex5p -m ../data/star.mesh -r 2 -pa -d raja-cuda
17 // mpirun -np 4 ex5p -m ../data/star.mesh -r 2 -pa -d raja-omp
18 // mpirun -np 4 ex5p -m ../data/beam-hex.mesh -pa -d cuda
19 //
20 // Description: This example code solves a simple 2D/3D mixed Darcy problem
21 // corresponding to the saddle point system
22 //
23 // k*u + grad p = f
24 // - div u = g
25 //
26 // with natural boundary condition -p = <given pressure>.
27 // Here, we use a given exact solution (u,p) and compute the
28 // corresponding r.h.s. (f,g). We discretize with Raviart-Thomas
29 // finite elements (velocity u) and piecewise discontinuous
30 // polynomials (pressure p).
31 //
32 // The example demonstrates the use of the BlockMatrix class, as
33 // well as the collective saving of several grid functions in
34 // VisIt (visit.llnl.gov) and ParaView (paraview.org) formats.
35 // Optional saving with ADIOS2 (adios2.readthedocs.io) streams is
36 // also illustrated.
37 //
38 // We recommend viewing examples 1-4 before viewing this example.
39 
40 #include "mfem.hpp"
41 #include <fstream>
42 #include <iostream>
43 
44 using namespace std;
45 using namespace mfem;
46 
47 // Define the analytical solution and forcing terms / boundary conditions
48 void uFun_ex(const Vector & x, Vector & u);
49 double pFun_ex(const Vector & x);
50 void fFun(const Vector & x, Vector & f);
51 double gFun(const Vector & x);
52 double f_natural(const Vector & x);
53 
54 int main(int argc, char *argv[])
55 {
56  StopWatch chrono;
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  bool verbose = (myid == 0);
64 
65  // 2. Parse command-line options.
66  const char *mesh_file = "../data/star.mesh";
67  int ref_levels = -1;
68  int order = 1;
69  bool par_format = false;
70  bool pa = false;
71  const char *device_config = "cpu";
72  bool visualization = 1;
73  bool adios2 = false;
74 
75  OptionsParser args(argc, argv);
76  args.AddOption(&mesh_file, "-m", "--mesh",
77  "Mesh file to use.");
78  args.AddOption(&ref_levels, "-r", "--refine",
79  "Number of times to refine the mesh uniformly.");
80  args.AddOption(&order, "-o", "--order",
81  "Finite element order (polynomial degree).");
82  args.AddOption(&par_format, "-pf", "--parallel-format", "-sf",
83  "--serial-format",
84  "Format to use when saving the results for VisIt.");
85  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
86  "--no-partial-assembly", "Enable Partial Assembly.");
87  args.AddOption(&device_config, "-d", "--device",
88  "Device configuration string, see Device::Configure().");
89  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
90  "--no-visualization",
91  "Enable or disable GLVis visualization.");
92  args.AddOption(&adios2, "-adios2", "--adios2-streams", "-no-adios2",
93  "--no-adios2-streams",
94  "Save data using adios2 streams.");
95  args.Parse();
96  if (!args.Good())
97  {
98  if (verbose)
99  {
100  args.PrintUsage(cout);
101  }
102  MPI_Finalize();
103  return 1;
104  }
105  if (verbose)
106  {
107  args.PrintOptions(cout);
108  }
109 
110  // 3. Enable hardware devices such as GPUs, and programming models such as
111  // CUDA, OCCA, RAJA and OpenMP based on command line options.
112  Device device(device_config);
113  if (myid == 0) { device.Print(); }
114 
115  // 4. Read the (serial) mesh from the given mesh file on all processors. We
116  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
117  // and volume meshes with the same code.
118  Mesh *mesh = new Mesh(mesh_file, 1, 1);
119  int dim = mesh->Dimension();
120 
121  // 5. Refine the serial mesh on all processors to increase the resolution. In
122  // this example we do 'ref_levels' of uniform refinement. We choose
123  // 'ref_levels' to be the largest number that gives a final mesh with no
124  // more than 10,000 elements, unless the user specifies it as input.
125  {
126  if (ref_levels == -1)
127  {
128  ref_levels = (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
129  }
130 
131  for (int l = 0; l < ref_levels; l++)
132  {
133  mesh->UniformRefinement();
134  }
135  }
136 
137  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
138  // this mesh further in parallel to increase the resolution. Once the
139  // parallel mesh is defined, the serial mesh can be deleted.
140  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
141  delete mesh;
142  {
143  int par_ref_levels = 2;
144  for (int l = 0; l < par_ref_levels; l++)
145  {
146  pmesh->UniformRefinement();
147  }
148  }
149 
150  // 7. Define a parallel finite element space on the parallel mesh. Here we
151  // use the Raviart-Thomas finite elements of the specified order.
152  FiniteElementCollection *hdiv_coll(new RT_FECollection(order, dim));
153  FiniteElementCollection *l2_coll(new L2_FECollection(order, dim));
154 
155  ParFiniteElementSpace *R_space = new ParFiniteElementSpace(pmesh, hdiv_coll);
156  ParFiniteElementSpace *W_space = new ParFiniteElementSpace(pmesh, l2_coll);
157 
158  HYPRE_Int dimR = R_space->GlobalTrueVSize();
159  HYPRE_Int dimW = W_space->GlobalTrueVSize();
160 
161  if (verbose)
162  {
163  std::cout << "***********************************************************\n";
164  std::cout << "dim(R) = " << dimR << "\n";
165  std::cout << "dim(W) = " << dimW << "\n";
166  std::cout << "dim(R+W) = " << dimR + dimW << "\n";
167  std::cout << "***********************************************************\n";
168  }
169 
170  // 8. Define the two BlockStructure of the problem. block_offsets is used
171  // for Vector based on dof (like ParGridFunction or ParLinearForm),
172  // block_trueOffstes is used for Vector based on trueDof (HypreParVector
173  // for the rhs and solution of the linear system). The offsets computed
174  // here are local to the processor.
175  Array<int> block_offsets(3); // number of variables + 1
176  block_offsets[0] = 0;
177  block_offsets[1] = R_space->GetVSize();
178  block_offsets[2] = W_space->GetVSize();
179  block_offsets.PartialSum();
180 
181  Array<int> block_trueOffsets(3); // number of variables + 1
182  block_trueOffsets[0] = 0;
183  block_trueOffsets[1] = R_space->TrueVSize();
184  block_trueOffsets[2] = W_space->TrueVSize();
185  block_trueOffsets.PartialSum();
186 
187  // 9. Define the coefficients, analytical solution, and rhs of the PDE.
188  ConstantCoefficient k(1.0);
189 
190  VectorFunctionCoefficient fcoeff(dim, fFun);
191  FunctionCoefficient fnatcoeff(f_natural);
192  FunctionCoefficient gcoeff(gFun);
193 
194  VectorFunctionCoefficient ucoeff(dim, uFun_ex);
196 
197  // 10. Define the parallel grid function and parallel linear forms, solution
198  // vector and rhs.
199  MemoryType mt = device.GetMemoryType();
200  BlockVector x(block_offsets, mt), rhs(block_offsets, mt);
201  BlockVector trueX(block_trueOffsets, mt), trueRhs(block_trueOffsets, mt);
202 
203  ParLinearForm *fform(new ParLinearForm);
204  fform->Update(R_space, rhs.GetBlock(0), 0);
207  fform->Assemble();
208  fform->SyncAliasMemory(rhs);
209  fform->ParallelAssemble(trueRhs.GetBlock(0));
210  trueRhs.GetBlock(0).SyncAliasMemory(trueRhs);
211 
212  ParLinearForm *gform(new ParLinearForm);
213  gform->Update(W_space, rhs.GetBlock(1), 0);
214  gform->AddDomainIntegrator(new DomainLFIntegrator(gcoeff));
215  gform->Assemble();
216  gform->SyncAliasMemory(rhs);
217  gform->ParallelAssemble(trueRhs.GetBlock(1));
218  trueRhs.GetBlock(1).SyncAliasMemory(trueRhs);
219 
220  // 11. Assemble the finite element matrices for the Darcy operator
221  //
222  // D = [ M B^T ]
223  // [ B 0 ]
224  // where:
225  //
226  // M = \int_\Omega k u_h \cdot v_h d\Omega u_h, v_h \in R_h
227  // B = -\int_\Omega \div u_h q_h d\Omega u_h \in R_h, q_h \in W_h
228  ParBilinearForm *mVarf(new ParBilinearForm(R_space));
229  ParMixedBilinearForm *bVarf(new ParMixedBilinearForm(R_space, W_space));
230 
231  HypreParMatrix *M = NULL;
232  HypreParMatrix *B = NULL;
233 
234  if (pa) { mVarf->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
236  mVarf->Assemble();
237  if (!pa) { mVarf->Finalize(); }
238 
239  if (pa) { bVarf->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
241  bVarf->Assemble();
242  if (!pa) { bVarf->Finalize(); }
243 
244  BlockOperator *darcyOp = new BlockOperator(block_trueOffsets);
245 
246  Array<int> empty_tdof_list; // empty
247  OperatorPtr opM, opB;
248 
249  TransposeOperator *Bt = NULL;
250 
251  if (pa)
252  {
253  mVarf->FormSystemMatrix(empty_tdof_list, opM);
254  bVarf->FormRectangularSystemMatrix(empty_tdof_list, empty_tdof_list, opB);
255  Bt = new TransposeOperator(opB.Ptr());
256 
257  darcyOp->SetBlock(0,0, opM.Ptr());
258  darcyOp->SetBlock(0,1, Bt, -1.0);
259  darcyOp->SetBlock(1,0, opB.Ptr(), -1.0);
260  }
261  else
262  {
263  M = mVarf->ParallelAssemble();
264  B = bVarf->ParallelAssemble();
265  (*B) *= -1;
266  Bt = new TransposeOperator(B);
267 
268  darcyOp->SetBlock(0,0, M);
269  darcyOp->SetBlock(0,1, Bt);
270  darcyOp->SetBlock(1,0, B);
271  }
272 
273  // 12. Construct the operators for preconditioner
274  //
275  // P = [ diag(M) 0 ]
276  // [ 0 B diag(M)^-1 B^T ]
277  //
278  // Here we use Symmetric Gauss-Seidel to approximate the inverse of the
279  // pressure Schur Complement.
280  HypreParMatrix *MinvBt = NULL;
281  HypreParVector *Md = NULL;
282  HypreParMatrix *S = NULL;
283  Vector Md_PA;
284  Solver *invM, *invS;
285 
286  if (pa)
287  {
288  Md_PA.SetSize(R_space->GetTrueVSize());
289  mVarf->AssembleDiagonal(Md_PA);
290  auto Md_host = Md_PA.HostRead();
291  Vector invMd(Md_PA.Size());
292  for (int i=0; i<Md_PA.Size(); ++i)
293  {
294  invMd(i) = 1.0 / Md_host[i];
295  }
296 
297  Vector BMBt_diag(W_space->GetTrueVSize());
298  bVarf->AssembleDiagonal_ADAt(invMd, BMBt_diag);
299 
300  Array<int> ess_tdof_list; // empty
301 
302  invM = new OperatorJacobiSmoother(Md_PA, ess_tdof_list);
303  invS = new OperatorJacobiSmoother(BMBt_diag, ess_tdof_list);
304  }
305  else
306  {
307  Md = new HypreParVector(MPI_COMM_WORLD, M->GetGlobalNumRows(),
308  M->GetRowStarts());
309  M->GetDiag(*Md);
310 
311  MinvBt = B->Transpose();
312  MinvBt->InvScaleRows(*Md);
313  S = ParMult(B, MinvBt);
314 
315  invM = new HypreDiagScale(*M);
316  invS = new HypreBoomerAMG(*S);
317  }
318 
319  invM->iterative_mode = false;
320  invS->iterative_mode = false;
321 
323  block_trueOffsets);
324  darcyPr->SetDiagonalBlock(0, invM);
325  darcyPr->SetDiagonalBlock(1, invS);
326 
327  // 13. Solve the linear system with MINRES.
328  // Check the norm of the unpreconditioned residual.
329  int maxIter(pa ? 1000 : 500);
330  double rtol(1.e-6);
331  double atol(1.e-10);
332 
333  chrono.Clear();
334  chrono.Start();
335  MINRESSolver solver(MPI_COMM_WORLD);
336  solver.SetAbsTol(atol);
337  solver.SetRelTol(rtol);
338  solver.SetMaxIter(maxIter);
339  solver.SetOperator(*darcyOp);
340  solver.SetPreconditioner(*darcyPr);
341  solver.SetPrintLevel(verbose);
342  trueX = 0.0;
343  solver.Mult(trueRhs, trueX);
344  if (device.IsEnabled()) { trueX.HostRead(); }
345  chrono.Stop();
346 
347  if (verbose)
348  {
349  if (solver.GetConverged())
350  std::cout << "MINRES converged in " << solver.GetNumIterations()
351  << " iterations with a residual norm of " << solver.GetFinalNorm() << ".\n";
352  else
353  std::cout << "MINRES did not converge in " << solver.GetNumIterations()
354  << " iterations. Residual norm is " << solver.GetFinalNorm() << ".\n";
355  std::cout << "MINRES solver took " << chrono.RealTime() << "s. \n";
356  }
357 
358  // 14. Extract the parallel grid function corresponding to the finite element
359  // approximation X. This is the local solution on each processor. Compute
360  // L2 error norms.
363  u->MakeRef(R_space, x.GetBlock(0), 0);
364  p->MakeRef(W_space, x.GetBlock(1), 0);
365  u->Distribute(&(trueX.GetBlock(0)));
366  p->Distribute(&(trueX.GetBlock(1)));
367 
368  int order_quad = max(2, 2*order+1);
369  const IntegrationRule *irs[Geometry::NumGeom];
370  for (int i=0; i < Geometry::NumGeom; ++i)
371  {
372  irs[i] = &(IntRules.Get(i, order_quad));
373  }
374 
375  double err_u = u->ComputeL2Error(ucoeff, irs);
376  double norm_u = ComputeGlobalLpNorm(2, ucoeff, *pmesh, irs);
377  double err_p = p->ComputeL2Error(pcoeff, irs);
378  double norm_p = ComputeGlobalLpNorm(2, pcoeff, *pmesh, irs);
379 
380  if (verbose)
381  {
382  std::cout << "|| u_h - u_ex || / || u_ex || = " << err_u / norm_u << "\n";
383  std::cout << "|| p_h - p_ex || / || p_ex || = " << err_p / norm_p << "\n";
384  }
385 
386  // 15. Save the refined mesh and the solution in parallel. This output can be
387  // viewed later using GLVis: "glvis -np <np> -m mesh -g sol_*".
388  {
389  ostringstream mesh_name, u_name, p_name;
390  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
391  u_name << "sol_u." << setfill('0') << setw(6) << myid;
392  p_name << "sol_p." << setfill('0') << setw(6) << myid;
393 
394  ofstream mesh_ofs(mesh_name.str().c_str());
395  mesh_ofs.precision(8);
396  pmesh->Print(mesh_ofs);
397 
398  ofstream u_ofs(u_name.str().c_str());
399  u_ofs.precision(8);
400  u->Save(u_ofs);
401 
402  ofstream p_ofs(p_name.str().c_str());
403  p_ofs.precision(8);
404  p->Save(p_ofs);
405  }
406 
407  // 16. Save data in the VisIt format
408  VisItDataCollection visit_dc("Example5-Parallel", pmesh);
409  visit_dc.RegisterField("velocity", u);
410  visit_dc.RegisterField("pressure", p);
411  visit_dc.SetFormat(!par_format ?
412  DataCollection::SERIAL_FORMAT :
413  DataCollection::PARALLEL_FORMAT);
414  visit_dc.Save();
415 
416  // 17. Save data in the ParaView format
417  ParaViewDataCollection paraview_dc("Example5P", pmesh);
418  paraview_dc.SetPrefixPath("ParaView");
419  paraview_dc.SetLevelsOfDetail(order);
420  paraview_dc.SetDataFormat(VTKFormat::BINARY);
421  paraview_dc.SetHighOrderOutput(true);
422  paraview_dc.SetCycle(0);
423  paraview_dc.SetTime(0.0);
424  paraview_dc.RegisterField("velocity",u);
425  paraview_dc.RegisterField("pressure",p);
426  paraview_dc.Save();
427 
428  // 18. Optionally output a BP (binary pack) file using ADIOS2. This can be
429  // visualized with the ParaView VTX reader.
430 #ifdef MFEM_USE_ADIOS2
431  if (adios2)
432  {
433  std::string postfix(mesh_file);
434  postfix.erase(0, std::string("../data/").size() );
435  postfix += "_o" + std::to_string(order);
436  const std::string collection_name = "ex5-p_" + postfix + ".bp";
437 
438  ADIOS2DataCollection adios2_dc(MPI_COMM_WORLD, collection_name, pmesh);
439  adios2_dc.SetLevelsOfDetail(1);
440  adios2_dc.SetCycle(1);
441  adios2_dc.SetTime(0.0);
442  adios2_dc.RegisterField("velocity",u);
443  adios2_dc.RegisterField("pressure",p);
444  adios2_dc.Save();
445  }
446 #endif
447 
448  // 19. Send the solution by socket to a GLVis server.
449  if (visualization)
450  {
451  char vishost[] = "localhost";
452  int visport = 19916;
453  socketstream u_sock(vishost, visport);
454  u_sock << "parallel " << num_procs << " " << myid << "\n";
455  u_sock.precision(8);
456  u_sock << "solution\n" << *pmesh << *u << "window_title 'Velocity'"
457  << endl;
458  // Make sure all ranks have sent their 'u' solution before initiating
459  // another set of GLVis connections (one from each rank):
460  MPI_Barrier(pmesh->GetComm());
461  socketstream p_sock(vishost, visport);
462  p_sock << "parallel " << num_procs << " " << myid << "\n";
463  p_sock.precision(8);
464  p_sock << "solution\n" << *pmesh << *p << "window_title 'Pressure'"
465  << endl;
466  }
467 
468  // 20. Free the used memory.
469  delete fform;
470  delete gform;
471  delete u;
472  delete p;
473  delete darcyOp;
474  delete darcyPr;
475  delete invM;
476  delete invS;
477  delete S;
478  delete Md;
479  delete MinvBt;
480  delete Bt;
481  delete B;
482  delete M;
483  delete mVarf;
484  delete bVarf;
485  delete W_space;
486  delete R_space;
487  delete l2_coll;
488  delete hdiv_coll;
489  delete pmesh;
490 
491  MPI_Finalize();
492 
493  return 0;
494 }
495 
496 
497 void uFun_ex(const Vector & x, Vector & u)
498 {
499  double xi(x(0));
500  double yi(x(1));
501  double zi(0.0);
502  if (x.Size() == 3)
503  {
504  zi = x(2);
505  }
506 
507  u(0) = - exp(xi)*sin(yi)*cos(zi);
508  u(1) = - exp(xi)*cos(yi)*cos(zi);
509 
510  if (x.Size() == 3)
511  {
512  u(2) = exp(xi)*sin(yi)*sin(zi);
513  }
514 }
515 
516 // Change if needed
517 double pFun_ex(const Vector & x)
518 {
519  double xi(x(0));
520  double yi(x(1));
521  double zi(0.0);
522 
523  if (x.Size() == 3)
524  {
525  zi = x(2);
526  }
527 
528  return exp(xi)*sin(yi)*cos(zi);
529 }
530 
531 void fFun(const Vector & x, Vector & f)
532 {
533  f = 0.0;
534 }
535 
536 double gFun(const Vector & x)
537 {
538  if (x.Size() == 3)
539  {
540  return -pFun_ex(x);
541  }
542  else
543  {
544  return 0;
545  }
546 }
547 
548 double f_natural(const Vector & x)
549 {
550  return (-pFun_ex(x));
551 }
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:400
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
int GetNumIterations() const
Definition: solvers.hpp:101
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:915
void fFun(const Vector &x, Vector &f)
Definition: ex5.cpp:419
A class to handle Vectors in a block fashion.
Definition: blockvector.hpp:30
void SetDataFormat(VTKFormat fmt)
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:459
std::string to_string(int i)
Convert an integer to an std::string.
Definition: text.hpp:51
Helper class for ParaView visualization data.
int TrueVSize() const
Obsolete, kept for backward compatibility.
Definition: pfespace.hpp:382
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
virtual void Save()
Save the collection and a VisIt root file.
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:819
void AssembleDiagonal(Vector &diag) const
Assemble the diagonal of the bilinear form into diag.
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:261
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
double RealTime()
Definition: tic_toc.cpp:426
Abstract parallel finite element space.
Definition: pfespace.hpp:28
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:638
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
MINRES method.
Definition: solvers.hpp:368
static bool IsEnabled()
Return true if any backend other than Backend::CPU is enabled.
Definition: device.hpp:237
int main(int argc, char *argv[])
Definition: ex1.cpp:66
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
Definition: hypre.cpp:1623
void Update(ParFiniteElementSpace *pf=NULL)
Update the object according to the given new FE space *pf.
Definition: plinearform.cpp:21
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
double GetFinalNorm() const
Definition: solvers.hpp:103
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:255
void Stop()
Stop the stopwatch.
Definition: tic_toc.cpp:416
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.
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Definition: pgridfunc.cpp:102
const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:376
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1079
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1365
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:415
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:381
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
constexpr char vishost[]
Jacobi preconditioner in hypre.
Definition: hypre.hpp:931
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:128
double f_natural(const Vector &x)
Definition: ex5.cpp:436
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:8382
constexpr int visport
Timing object.
Definition: tic_toc.hpp:34
static MemoryType GetMemoryType()
(DEPRECATED) Equivalent to GetDeviceMemoryType().
Definition: device.hpp:265
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
Definition: solvers.hpp:98
void Assemble(int skip_zeros=1)
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
virtual void SetFormat(int fmt)
Set the desired output mesh and data format.
double pFun_ex(const Vector &x)
Definition: ex5.cpp:405
virtual void FormRectangularSystemMatrix(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)
Return in A a parallel (on truedofs) version of this operator.
HypreParMatrix * Transpose() const
Returns the transpose of *this.
Definition: hypre.cpp:951
void SetHighOrderOutput(bool high_order_output_)
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:887
void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const
Assemble the diagonal of ADA^T into diag, where A is this mixed bilinear form and D is a diagonal...
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1351
int GetConverged() const
Definition: solvers.hpp:102
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:53
int Dimension() const
Definition: mesh.hpp:788
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
void SetTime(double t)
Set physical time (for time-dependent simulations)
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
void Start()
Clear the elapsed time and start the stopwatch.
Definition: tic_toc.cpp:411
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:272
A general vector function coefficient.
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACYFULL.
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:70
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
void SetAbsTol(double atol)
Definition: solvers.hpp:97
void SetRelTol(double rtol)
Definition: solvers.hpp:96
MPI_Comm GetComm() const
Definition: pmesh.hpp:236
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:28
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
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
Definition: array.cpp:103
double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:132
The transpose of a given operator. Switches the roles of the methods Mult() and MultTranspose().
Definition: operator.hpp:697
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
void uFun_ex(const Vector &x, Vector &u)
Definition: ex5.cpp:385
double gFun(const Vector &x)
Definition: ex5.cpp:424
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:267
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
Definition: hypre.cpp:1213
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.
void SetLevelsOfDetail(int levels_of_detail_)
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:260
A general function coefficient.
Vector data type.
Definition: vector.hpp:51
void SetLevelsOfDetail(const int levels_of_detail) noexcept
virtual void Save() override
HYPRE_Int * GetRowStarts() const
Return the parallel row partitioning array.
Definition: hypre.hpp:425
Base class for solvers.
Definition: operator.hpp:634
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:118
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
A class to handle Block systems in a matrix-free implementation.
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:378
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:46
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
void Clear()
Clear the elapsed time on the stopwatch and restart it if it&#39;s running.
Definition: tic_toc.cpp:406
void SetBlock(int iRow, int iCol, Operator *op, double c=1.0)
Add a block op in the block-entry (iblock, jblock).
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
Vector & GetBlock(int i)
Get the i-th vector in the block.
Definition: blockvector.hpp:87
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:221
double f(const Vector &p)
void SyncAliasMemory(const Vector &v)
Update the alias memory location of the vector to match v.
Definition: vector.hpp:193
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145