MFEM  v4.6.0
Finite element discretization library
ex26p.cpp
Go to the documentation of this file.
1 // MFEM Example 26 - Parallel Version
2 //
3 // Compile with: make ex26p
4 //
5 // Sample runs: mpirun -np 4 ex26p -m ../data/star.mesh
6 // mpirun -np 4 ex26p -m ../data/fichera.mesh
7 // mpirun -np 4 ex26p -m ../data/beam-hex.mesh
8 //
9 // Device sample runs:
10 // mpirun -np 4 ex26p -d cuda
11 // mpirun -np 4 ex26p -d occa-cuda
12 // mpirun -np 4 ex26p -d raja-omp
13 // mpirun -np 4 ex26p -d ceed-cpu
14 // mpirun -np 4 ex26p -d ceed-cuda
15 //
16 // Description: This example code demonstrates the use of MFEM to define a
17 // simple finite element discretization of the Laplace problem
18 // -Delta u = 1 with homogeneous Dirichlet boundary conditions
19 // as in Example 1.
20 //
21 // It highlights on the creation of a hierarchy of discretization
22 // spaces with partial assembly and the construction of an
23 // efficient multigrid preconditioner for the iterative solver.
24 //
25 // We recommend viewing Example 1 before viewing this example.
26 
27 #include "mfem.hpp"
28 #include <fstream>
29 #include <iostream>
30 
31 using namespace std;
32 using namespace mfem;
33 
34 // Class for constructing a multigrid preconditioner for the diffusion operator.
35 // This example multigrid preconditioner class demonstrates the creation of the
36 // parallel diffusion bilinear forms and operators using partial assembly for
37 // all spaces except the coarsest one in the ParFiniteElementSpaceHierarchy.
38 // The multigrid uses a PCG solver preconditioned with AMG on the coarsest level
39 // and second order Chebyshev accelerated smoothers on the other levels.
40 class DiffusionMultigrid : public GeometricMultigrid
41 {
42 private:
44  HypreBoomerAMG* amg;
45 
46 public:
47  // Constructs a diffusion multigrid for the ParFiniteElementSpaceHierarchy
48  // and the array of essential boundaries
49  DiffusionMultigrid(ParFiniteElementSpaceHierarchy& fespaces,
50  Array<int>& ess_bdr)
51  : GeometricMultigrid(fespaces), one(1.0)
52  {
53  ConstructCoarseOperatorAndSolver(fespaces.GetFESpaceAtLevel(0), ess_bdr);
54 
55  for (int level = 1; level < fespaces.GetNumLevels(); ++level)
56  {
57  ConstructOperatorAndSmoother(fespaces.GetFESpaceAtLevel(level), ess_bdr);
58  }
59  }
60 
61  virtual ~DiffusionMultigrid()
62  {
63  delete amg;
64  }
65 
66 private:
67  void ConstructBilinearForm(ParFiniteElementSpace& fespace, Array<int>& ess_bdr,
68  bool partial_assembly)
69  {
70  ParBilinearForm* form = new ParBilinearForm(&fespace);
71  if (partial_assembly)
72  {
73  form->SetAssemblyLevel(AssemblyLevel::PARTIAL);
74  }
76  form->Assemble();
77  bfs.Append(form);
78 
79  essentialTrueDofs.Append(new Array<int>());
80  fespace.GetEssentialTrueDofs(ess_bdr, *essentialTrueDofs.Last());
81  }
82 
83  void ConstructCoarseOperatorAndSolver(ParFiniteElementSpace& coarse_fespace,
84  Array<int>& ess_bdr)
85  {
86  ConstructBilinearForm(coarse_fespace, ess_bdr, false);
87 
88  HypreParMatrix* hypreCoarseMat = new HypreParMatrix();
89  bfs.Last()->FormSystemMatrix(*essentialTrueDofs.Last(), *hypreCoarseMat);
90 
91  amg = new HypreBoomerAMG(*hypreCoarseMat);
92  amg->SetPrintLevel(-1);
93 
94  CGSolver* pcg = new CGSolver(MPI_COMM_WORLD);
95  pcg->SetPrintLevel(-1);
96  pcg->SetMaxIter(10);
97  pcg->SetRelTol(sqrt(1e-4));
98  pcg->SetAbsTol(0.0);
99  pcg->SetOperator(*hypreCoarseMat);
100  pcg->SetPreconditioner(*amg);
101 
102  AddLevel(hypreCoarseMat, pcg, true, true);
103  }
104 
105  void ConstructOperatorAndSmoother(ParFiniteElementSpace& fespace,
106  Array<int>& ess_bdr)
107  {
108  ConstructBilinearForm(fespace, ess_bdr, true);
109 
110  OperatorPtr opr;
111  opr.SetType(Operator::ANY_TYPE);
112  bfs.Last()->FormSystemMatrix(*essentialTrueDofs.Last(), opr);
113  opr.SetOperatorOwner(false);
114 
115  Vector diag(fespace.GetTrueVSize());
116  bfs.Last()->AssembleDiagonal(diag);
117 
118  Solver* smoother = new OperatorChebyshevSmoother(*opr, diag,
119  *essentialTrueDofs.Last(), 2, fespace.GetParMesh()->GetComm());
120 
121  AddLevel(opr.Ptr(), smoother, true, true);
122  }
123 };
124 
125 
126 int main(int argc, char *argv[])
127 {
128  // 1. Initialize MPI and HYPRE.
129  Mpi::Init(argc, argv);
130  int num_procs = Mpi::WorldSize();
131  int myid = Mpi::WorldRank();
132  Hypre::Init();
133 
134  // 2. Parse command-line options.
135  const char *mesh_file = "../data/star.mesh";
136  int geometric_refinements = 0;
137  int order_refinements = 2;
138  const char *device_config = "cpu";
139  bool visualization = true;
140 
141  OptionsParser args(argc, argv);
142  args.AddOption(&mesh_file, "-m", "--mesh",
143  "Mesh file to use.");
144  args.AddOption(&geometric_refinements, "-gr", "--geometric-refinements",
145  "Number of geometric refinements done prior to order refinements.");
146  args.AddOption(&order_refinements, "-or", "--order-refinements",
147  "Number of order refinements. Finest level in the hierarchy has order 2^{or}.");
148  args.AddOption(&device_config, "-d", "--device",
149  "Device configuration string, see Device::Configure().");
150  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
151  "--no-visualization",
152  "Enable or disable GLVis visualization.");
153  args.Parse();
154  if (!args.Good())
155  {
156  if (myid == 0)
157  {
158  args.PrintUsage(cout);
159  }
160  return 1;
161  }
162  if (myid == 0)
163  {
164  args.PrintOptions(cout);
165  }
166 
167  // 3. Enable hardware devices such as GPUs, and programming models such as
168  // CUDA, OCCA, RAJA and OpenMP based on command line options.
169  Device device(device_config);
170  if (myid == 0) { device.Print(); }
171 
172  // 4. Read the (serial) mesh from the given mesh file on all processors. We
173  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
174  // and volume meshes with the same code.
175  Mesh *mesh = new Mesh(mesh_file, 1, 1);
176  int dim = mesh->Dimension();
177 
178  // 5. Refine the serial mesh on all processors to increase the resolution. In
179  // this example we do 'ref_levels' of uniform refinement. We choose
180  // 'ref_levels' to be the largest number that gives a final mesh with no
181  // more than 1,000 elements.
182  {
183  int ref_levels =
184  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
185  for (int l = 0; l < ref_levels; l++)
186  {
187  mesh->UniformRefinement();
188  }
189  }
190 
191  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
192  // this mesh further in parallel to increase the resolution. Once the
193  // parallel mesh is defined, the serial mesh can be deleted.
194  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
195  delete mesh;
196  {
197  int par_ref_levels = 2;
198  for (int l = 0; l < par_ref_levels; l++)
199  {
200  pmesh->UniformRefinement();
201  }
202  }
203 
204  // 7. Define a parallel finite element space hierarchy on the parallel mesh.
205  // Here we use continuous Lagrange finite elements. We start with order 1
206  // on the coarse level and geometrically refine the spaces by the specified
207  // amount. Afterwards, we increase the order of the finite elements by a
208  // factor of 2 for each additional level.
210  ParFiniteElementSpace *coarse_fespace = new ParFiniteElementSpace(pmesh, fec);
211 
213  collections.Append(fec);
215  pmesh, coarse_fespace, true, true);
216  for (int level = 0; level < geometric_refinements; ++level)
217  {
218  fespaces->AddUniformlyRefinedLevel();
219  }
220  for (int level = 0; level < order_refinements; ++level)
221  {
222  collections.Append(new H1_FECollection((int)std::pow(2, level+1), dim));
223  fespaces->AddOrderRefinedLevel(collections.Last());
224  }
225 
226  HYPRE_BigInt size = fespaces->GetFinestFESpace().GlobalTrueVSize();
227  if (myid == 0)
228  {
229  cout << "Number of finite element unknowns: " << size << endl;
230  }
231 
232  // 8. Set up the parallel linear form b(.) which corresponds to the
233  // right-hand side of the FEM linear system, which in this case is
234  // (1,phi_i) where phi_i are the basis functions in fespace.
235  ParLinearForm *b = new ParLinearForm(&fespaces->GetFinestFESpace());
236  ConstantCoefficient one(1.0);
237  b->AddDomainIntegrator(new DomainLFIntegrator(one));
238  b->Assemble();
239 
240  // 9. Define the solution vector x as a parallel finite element grid function
241  // corresponding to fespace. Initialize x with initial guess of zero,
242  // which satisfies the boundary conditions.
243  ParGridFunction x(&fespaces->GetFinestFESpace());
244  x = 0.0;
245 
246  // 10. Create the multigrid operator using the previously created parallel
247  // FiniteElementSpaceHierarchy and additional boundary information. This
248  // operator is then used to create the MultigridSolver as preconditioner
249  // in the iterative solver.
250  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
251  if (pmesh->bdr_attributes.Size())
252  {
253  ess_bdr = 1;
254  }
255 
256  DiffusionMultigrid* M = new DiffusionMultigrid(*fespaces, ess_bdr);
257  M->SetCycleType(Multigrid::CycleType::VCYCLE, 1, 1);
258 
259  OperatorPtr A;
260  Vector X, B;
261  M->FormFineLinearSystem(x, *b, A, X, B);
262 
263  // 11. Solve the linear system A X = B.
264  CGSolver cg(MPI_COMM_WORLD);
265  cg.SetRelTol(1e-12);
266  cg.SetMaxIter(2000);
267  cg.SetPrintLevel(1);
268  cg.SetOperator(*A);
269  cg.SetPreconditioner(*M);
270  cg.Mult(B, X);
271 
272  // 12. Recover the parallel grid function corresponding to X. This is the
273  // local finite element solution on each processor.
274  M->RecoverFineFEMSolution(X, *b, x);
275 
276  // 13. Save the refined mesh and the solution in parallel. This output can be
277  // viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
278  {
279  ostringstream mesh_name, sol_name;
280  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
281  sol_name << "sol." << setfill('0') << setw(6) << myid;
282 
283  ofstream mesh_ofs(mesh_name.str().c_str());
284  mesh_ofs.precision(8);
285  fespaces->GetFinestFESpace().GetParMesh()->Print(mesh_ofs);
286 
287  ofstream sol_ofs(sol_name.str().c_str());
288  sol_ofs.precision(8);
289  x.Save(sol_ofs);
290  }
291 
292  // 14. Send the solution by socket to a GLVis server.
293  if (visualization)
294  {
295  char vishost[] = "localhost";
296  int visport = 19916;
297  socketstream sol_sock(vishost, visport);
298  sol_sock << "parallel " << num_procs << " " << myid << "\n";
299  sol_sock.precision(8);
300  sol_sock << "solution\n" << *fespaces->GetFinestFESpace().GetParMesh()
301  << x << flush;
302  }
303 
304  // 15. Free the used memory.
305  delete M;
306  delete b;
307  delete fespaces;
308  for (int level = 0; level < collections.Size(); ++level)
309  {
310  delete collections[level];
311  }
312 
313  return 0;
314 }
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:108
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1031
int visport
Conjugate gradient method.
Definition: solvers.hpp:493
Chebyshev accelerated smoothing with given vector, no matrix necessary.
Definition: solvers.hpp:379
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:331
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:462
Geometric multigrid associated with a hierarchy of finite element spaces.
Definition: multigrid.hpp:165
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:718
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition: array.cpp:68
const ParFiniteElementSpace & GetFinestFESpace() const override
Returns the finite element space at the finest level.
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.
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:159
Abstract parallel finite element space.
Definition: pfespace.hpp:28
STL namespace.
int main(int argc, char *argv[])
Definition: ex26p.cpp:126
const ParFiniteElementSpace & GetFESpaceAtLevel(int level) const override
Returns the finite element space at the given level.
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1590
int GetNumLevels() const
Returns the number of levels in the hierarchy.
Class for parallel linear form.
Definition: plinearform.hpp:26
void AddOrderRefinedLevel(FiniteElementCollection *fec, int dim=1, int ordering=Ordering::byVDIM) override
Adds one level to the hierarchy by using a different finite element order defined through FiniteEleme...
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
void AddUniformlyRefinedLevel(int dim=1, int ordering=Ordering::byVDIM) override
Adds one level to the hierarchy by uniformly refining the mesh on the previous level.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
char vishost[]
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:759
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:10232
void SetPrintLevel(int print_level)
Definition: hypre.hpp:1670
ParMesh * GetParMesh() const
Definition: pfespace.hpp:273
void SetMaxIter(int max_it)
Definition: solvers.hpp:201
void Assemble(int skip_zeros=1)
Assemble the local matrix.
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:281
MPI_Comm GetComm() const
Definition: pmesh.hpp:351
void SetAbsTol(double atol)
Definition: solvers.hpp:200
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:275
void SetRelTol(double rtol)
Definition: solvers.hpp:199
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
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:1086
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:285
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:120
int dim
Definition: ex24.cpp:53
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
Class for parallel bilinear form.
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:87
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
T & Last()
Return the last element in the array.
Definition: array.hpp:792
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:507
Vector data type.
Definition: vector.hpp:58
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:259
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4825
Base class for solvers.
Definition: operator.hpp:682
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
Class for parallel meshes.
Definition: pmesh.hpp:32
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition: handle.hpp:132