MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex14p.cpp
Go to the documentation of this file.
1 // MFEM Example 14 - Parallel Version
2 //
3 // Compile with: make ex14p
4 //
5 // Sample runs: mpirun -np 4 ex14p -m ../data/inline-quad.mesh -o 0
6 // mpirun -np 4 ex14p -m ../data/star.mesh -o 2
7 // mpirun -np 4 ex14p -m ../data/star-mixed.mesh -o 2
8 // mpirun -np 4 ex14p -m ../data/star-mixed.mesh -o 2 -k 0 -e 1
9 // mpirun -np 4 ex14p -m ../data/escher.mesh -s 1
10 // mpirun -np 4 ex14p -m ../data/fichera.mesh -s 1 -k 1
11 // mpirun -np 4 ex14p -m ../data/fichera-mixed.mesh -s 1 -k 1
12 // mpirun -np 4 ex14p -m ../data/square-disc-p2.vtk -o 2
13 // mpirun -np 4 ex14p -m ../data/square-disc-p3.mesh -o 3
14 // mpirun -np 4 ex14p -m ../data/square-disc-nurbs.mesh -o 1
15 // mpirun -np 4 ex14p -m ../data/disc-nurbs.mesh -rs 4 -o 2 -s 1 -k 0
16 // mpirun -np 4 ex14p -m ../data/pipe-nurbs.mesh -o 1
17 // mpirun -np 4 ex14p -m ../data/inline-segment.mesh -rs 5
18 // mpirun -np 4 ex14p -m ../data/amr-quad.mesh -rs 3
19 // mpirun -np 4 ex14p -m ../data/amr-hex.mesh
20 //
21 // Description: This example code demonstrates the use of MFEM to define a
22 // discontinuous Galerkin (DG) finite element discretization of
23 // the Laplace problem -Delta u = 1 with homogeneous Dirichlet
24 // boundary conditions. Finite element spaces of any order,
25 // including zero on regular grids, are supported. The example
26 // highlights the use of discontinuous spaces and DG-specific face
27 // integrators.
28 //
29 // We recommend viewing examples 1 and 9 before viewing this
30 // example.
31 
32 #include "mfem.hpp"
33 #include <fstream>
34 #include <iostream>
35 
36 using namespace std;
37 using namespace mfem;
38 
39 class CustomSolverMonitor : public IterativeSolverMonitor
40 {
41 public:
42  CustomSolverMonitor(const ParMesh *m,
43  ParGridFunction *f) :
44  pmesh(m),
45  pgf(f) {}
46 
47  void MonitorSolution(int i, double norm, const Vector &x, bool final)
48  {
49  char vishost[] = "localhost";
50  int visport = 19916;
51  int num_procs, myid;
52 
53  MPI_Comm_size(pmesh->GetComm(),&num_procs);
54  MPI_Comm_rank(pmesh->GetComm(),&myid);
55 
56  pgf->SetFromTrueDofs(x);
57 
58  socketstream sol_sock(vishost, visport);
59  sol_sock << "parallel " << num_procs << " " << myid << "\n";
60  sol_sock.precision(8);
61  sol_sock << "solution\n" << *pmesh << *pgf
62  << "window_title 'Iteration no " << i << "'"
63  << "keys rRjlc\n" << flush;
64  }
65 
66 private:
67  const ParMesh *pmesh;
68  ParGridFunction *pgf;
69 };
70 
71 int main(int argc, char *argv[])
72 {
73  // 1. Initialize MPI.
74  int num_procs, myid;
75  MPI_Init(&argc, &argv);
76  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
77  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
78 
79  // 2. Parse command-line options.
80  const char *mesh_file = "../data/star.mesh";
81  int ser_ref_levels = -1;
82  int par_ref_levels = 2;
83  int order = 1;
84  double sigma = -1.0;
85  double kappa = -1.0;
86  double eta = 0.0;
87  bool visualization = 1;
88 
89  OptionsParser args(argc, argv);
90  args.AddOption(&mesh_file, "-m", "--mesh",
91  "Mesh file to use.");
92  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
93  "Number of times to refine the mesh uniformly in serial,"
94  " -1 for auto.");
95  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
96  "Number of times to refine the mesh uniformly in parallel.");
97  args.AddOption(&order, "-o", "--order",
98  "Finite element order (polynomial degree) >= 0.");
99  args.AddOption(&sigma, "-s", "--sigma",
100  "One of the three DG penalty parameters, typically +1/-1."
101  " See the documentation of class DGDiffusionIntegrator.");
102  args.AddOption(&kappa, "-k", "--kappa",
103  "One of the three DG penalty parameters, should be positive."
104  " Negative values are replaced with (order+1)^2.");
105  args.AddOption(&eta, "-e", "--eta", "BR2 penalty parameter.");
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 (kappa < 0)
120  {
121  kappa = (order+1)*(order+1);
122  }
123  if (myid == 0)
124  {
125  args.PrintOptions(cout);
126  }
127 
128  // 3. Read the (serial) mesh from the given mesh file on all processors. We
129  // can handle triangular, quadrilateral, tetrahedral and hexahedral meshes
130  // with the same code. NURBS meshes are projected to second order meshes.
131  Mesh *mesh = new Mesh(mesh_file, 1, 1);
132  int dim = mesh->Dimension();
133 
134  // 4. Refine the serial mesh on all processors to increase the resolution. In
135  // this example we do 'ser_ref_levels' of uniform refinement. By default,
136  // or if ser_ref_levels < 0, we choose it to be the largest number that
137  // gives a final mesh with no more than 50,000 elements.
138  {
139  if (ser_ref_levels < 0)
140  {
141  ser_ref_levels = (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
142  }
143  for (int l = 0; l < ser_ref_levels; l++)
144  {
145  mesh->UniformRefinement();
146  }
147  }
148  if (mesh->NURBSext)
149  {
150  mesh->SetCurvature(max(order, 1));
151  }
152 
153  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
154  // this mesh further in parallel to increase the resolution. Once the
155  // parallel mesh is defined, the serial mesh can be deleted.
156  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
157  delete mesh;
158  {
159  for (int l = 0; l < par_ref_levels; l++)
160  {
161  pmesh->UniformRefinement();
162  }
163  }
164 
165  // 6. Define a parallel finite element space on the parallel mesh. Here we
166  // use discontinuous finite elements of the specified order >= 0.
167  FiniteElementCollection *fec = new DG_FECollection(order, dim);
168  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
169  HYPRE_BigInt size = fespace->GlobalTrueVSize();
170  if (myid == 0)
171  {
172  cout << "Number of unknowns: " << size << endl;
173  }
174 
175  // 7. Set up the parallel linear form b(.) which corresponds to the
176  // right-hand side of the FEM linear system.
177  ParLinearForm *b = new ParLinearForm(fespace);
178  ConstantCoefficient one(1.0);
179  ConstantCoefficient zero(0.0);
182  new DGDirichletLFIntegrator(zero, one, sigma, kappa));
183  b->Assemble();
184 
185  // 8. Define the solution vector x as a parallel finite element grid function
186  // corresponding to fespace. Initialize x with initial guess of zero.
187  ParGridFunction x(fespace);
188  x = 0.0;
189 
190  // 9. Set up the bilinear form a(.,.) on the finite element space
191  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
192  // domain integrator and the interior and boundary DG face integrators.
193  // Note that boundary conditions are imposed weakly in the form, so there
194  // is no need for dof elimination. After serial and parallel assembly we
195  // extract the corresponding parallel matrix A.
196  ParBilinearForm *a = new ParBilinearForm(fespace);
198  a->AddInteriorFaceIntegrator(new DGDiffusionIntegrator(one, sigma, kappa));
199  a->AddBdrFaceIntegrator(new DGDiffusionIntegrator(one, sigma, kappa));
200  if (eta > 0)
201  {
203  a->AddBdrFaceIntegrator(new DGDiffusionBR2Integrator(fespace, eta));
204  }
205  a->Assemble();
206  a->Finalize();
207 
208  // 10. Define the parallel (hypre) matrix and vectors representing a(.,.),
209  // b(.) and the finite element approximation.
213 
214  delete a;
215  delete b;
216 
217  // 11. Depending on the symmetry of A, define and apply a parallel PCG or
218  // GMRES solver for AX=B using the BoomerAMG preconditioner from hypre.
219  HypreSolver *amg = new HypreBoomerAMG(*A);
220  if (sigma == -1.0)
221  {
222  HyprePCG pcg(*A);
223  pcg.SetTol(1e-12);
224  pcg.SetMaxIter(200);
225  pcg.SetPrintLevel(2);
226  pcg.SetPreconditioner(*amg);
227  pcg.Mult(*B, *X);
228  }
229  else
230  {
231  CustomSolverMonitor monitor(pmesh, &x);
232  GMRESSolver gmres(MPI_COMM_WORLD);
233  gmres.SetAbsTol(0.0);
234  gmres.SetRelTol(1e-12);
235  gmres.SetMaxIter(200);
236  gmres.SetKDim(10);
237  gmres.SetPrintLevel(1);
238  gmres.SetOperator(*A);
239  gmres.SetPreconditioner(*amg);
240  gmres.SetMonitor(monitor);
241  gmres.Mult(*B, *X);
242  }
243  delete amg;
244 
245  // 12. Extract the parallel grid function corresponding to the finite element
246  // approximation X. This is the local solution on each processor.
247  x = *X;
248 
249  // 13. Save the refined mesh and the solution in parallel. This output can
250  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
251  {
252  ostringstream mesh_name, sol_name;
253  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
254  sol_name << "sol." << setfill('0') << setw(6) << myid;
255 
256  ofstream mesh_ofs(mesh_name.str().c_str());
257  mesh_ofs.precision(8);
258  pmesh->Print(mesh_ofs);
259 
260  ofstream sol_ofs(sol_name.str().c_str());
261  sol_ofs.precision(8);
262  x.Save(sol_ofs);
263  }
264 
265  // 14. Send the solution by socket to a GLVis server.
266  if (visualization)
267  {
268  char vishost[] = "localhost";
269  int visport = 19916;
270  socketstream sol_sock(vishost, visport);
271  sol_sock << "parallel " << num_procs << " " << myid << "\n";
272  sol_sock.precision(8);
273  sol_sock << "solution\n" << *pmesh << x << flush;
274  }
275 
276  // 15. Free the used memory.
277  delete X;
278  delete B;
279  delete A;
280 
281  delete fespace;
282  delete fec;
283  delete pmesh;
284 
285  MPI_Finalize();
286 
287  return 0;
288 }
void SetTol(double tol)
Definition: hypre.cpp:3569
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:97
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:275
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:841
double kappa
Definition: ex24.cpp:54
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3589
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:884
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1387
Class for parallel linear form.
Definition: plinearform.hpp:26
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:98
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[]
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
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:4882
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void SetKDim(int dim)
Set the number of iteration to perform between restarts, default is 50.
Definition: solvers.hpp:361
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4382
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:83
int Dimension() const
Definition: mesh.hpp:911
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3579
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
Definition: fe_coll.hpp:338
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:99
Abstract base class for an iterative solver monitor.
Definition: solvers.hpp:36
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
void SetAbsTol(double atol)
Definition: solvers.hpp:99
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
GMRES method.
Definition: solvers.hpp:348
void SetMonitor(IterativeSolverMonitor &m)
Set the iterative solver monitor.
Definition: solvers.hpp:114
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:206
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
int dim
Definition: ex24.cpp:53
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator. Assumes ownership of bfi.
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
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:970
Vector data type.
Definition: vector.hpp:60
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:92
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:277
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
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:87
int main()
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
Definition: pgridfunc.cpp:171
double f(const Vector &p)
double sigma(const Vector &x)
Definition: maxwell.cpp:102
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150