MFEM  v4.5.2
Finite element discretization library
ex4p.cpp
Go to the documentation of this file.
1 // MFEM Example 4 - Parallel Version
2 // PETSc Modification
3 //
4 // Compile with: make ex4p
5 //
6 // Sample runs:
7 // mpirun -np 4 ex4p -m ../../data/klein-bottle.mesh -o 2 --petscopts rc_ex4p
8 // mpirun -np 4 ex4p -m ../../data/klein-bottle.mesh -o 2 --petscopts rc_ex4p_bddc --nonoverlapping
9 //
10 // Description: This example code solves a simple 2D/3D H(div) diffusion
11 // problem corresponding to the second order definite equation
12 // -grad(alpha div F) + beta F = f with boundary condition F dot n
13 // = <given normal field>. Here, we use a given exact solution F
14 // and compute the corresponding r.h.s. f. We discretize with
15 // Raviart-Thomas finite elements.
16 //
17 // The example demonstrates the use of H(div) finite element
18 // spaces with the grad-div and H(div) vector finite element mass
19 // bilinear form, as well as the computation of discretization
20 // error when the exact solution is known. Bilinear form
21 // hybridization and static condensation are also illustrated.
22 //
23 // We recommend viewing examples 1-3 before viewing this example.
24 
25 #include "mfem.hpp"
26 #include <fstream>
27 #include <iostream>
28 
29 #ifndef MFEM_USE_PETSC
30 #error This example requires that MFEM is built with MFEM_USE_PETSC=YES
31 #endif
32 
33 using namespace std;
34 using namespace mfem;
35 
36 // Exact solution, F, and r.h.s., f. See below for implementation.
37 void F_exact(const Vector &, Vector &);
38 void f_exact(const Vector &, Vector &);
39 double freq = 1.0, kappa;
40 
41 int main(int argc, char *argv[])
42 {
43  // 1. Initialize MPI and HYPRE.
44  Mpi::Init(argc, argv);
45  int num_procs = Mpi::WorldSize();
46  int myid = Mpi::WorldRank();
47  Hypre::Init();
48 
49  // 2. Parse command-line options.
50  const char *mesh_file = "../../data/star.mesh";
51  int order = 1;
52  bool set_bc = true;
53  bool static_cond = false;
54  bool hybridization = false;
55  bool visualization = 1;
56  bool use_petsc = true;
57  const char *petscrc_file = "";
58  bool use_nonoverlapping = false;
59 
60  OptionsParser args(argc, argv);
61  args.AddOption(&mesh_file, "-m", "--mesh",
62  "Mesh file to use.");
63  args.AddOption(&order, "-o", "--order",
64  "Finite element order (polynomial degree).");
65  args.AddOption(&set_bc, "-bc", "--impose-bc", "-no-bc", "--dont-impose-bc",
66  "Impose or not essential boundary conditions.");
67  args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
68  " solution.");
69  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
70  "--no-static-condensation", "Enable static condensation.");
71  args.AddOption(&hybridization, "-hb", "--hybridization", "-no-hb",
72  "--no-hybridization", "Enable hybridization.");
73  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
74  "--no-visualization",
75  "Enable or disable GLVis visualization.");
76  args.AddOption(&use_petsc, "-usepetsc", "--usepetsc", "-no-petsc",
77  "--no-petsc",
78  "Use or not PETSc to solve the linear system.");
79  args.AddOption(&petscrc_file, "-petscopts", "--petscopts",
80  "PetscOptions file to use.");
81  args.AddOption(&use_nonoverlapping, "-nonoverlapping", "--nonoverlapping",
82  "-no-nonoverlapping", "--no-nonoverlapping",
83  "Use or not the block diagonal PETSc's matrix format "
84  "for non-overlapping domain decomposition.");
85  args.Parse();
86  if (!args.Good())
87  {
88  if (myid == 0)
89  {
90  args.PrintUsage(cout);
91  }
92  return 1;
93  }
94  if (myid == 0)
95  {
96  args.PrintOptions(cout);
97  }
98  // 2b. We initialize PETSc
99  if (use_petsc) { MFEMInitializePetsc(NULL,NULL,petscrc_file,NULL); }
100  kappa = freq * M_PI;
101 
102  // 3. Read the (serial) mesh from the given mesh file on all processors. We
103  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
104  // and volume, as well as periodic meshes with the same code.
105  Mesh *mesh = new Mesh(mesh_file, 1, 1);
106  int dim = mesh->Dimension();
107  int sdim = mesh->SpaceDimension();
108 
109  // 4. Refine the serial mesh on all processors to increase the resolution. In
110  // this example we do 'ref_levels' of uniform refinement. We choose
111  // 'ref_levels' to be the largest number that gives a final mesh with no
112  // more than 1,000 elements.
113  {
114  int ref_levels =
115  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
116  for (int l = 0; l < ref_levels; l++)
117  {
118  mesh->UniformRefinement();
119  }
120  }
121 
122  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
123  // this mesh further in parallel to increase the resolution. Once the
124  // parallel mesh is defined, the serial mesh can be deleted.
125  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
126  delete mesh;
127  {
128  int par_ref_levels = 2;
129  for (int l = 0; l < par_ref_levels; l++)
130  {
131  pmesh->UniformRefinement();
132  }
133  }
134 
135  // 6. Define a parallel finite element space on the parallel mesh. Here we
136  // use the Raviart-Thomas finite elements of the specified order.
137  FiniteElementCollection *fec = new RT_FECollection(order-1, dim);
138  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
139  HYPRE_BigInt size = fespace->GlobalTrueVSize();
140  if (myid == 0)
141  {
142  cout << "Number of finite element unknowns: " << size << endl;
143  }
144 
145  // 7. Determine the list of true (i.e. parallel conforming) essential
146  // boundary dofs. In this example, the boundary conditions are defined
147  // by marking all the boundary attributes from the mesh as essential
148  // (Dirichlet) and converting them to a list of true dofs.
149  Array<int> ess_tdof_list;
150  if (pmesh->bdr_attributes.Size())
151  {
152  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
153  ess_bdr = set_bc ? 1 : 0;
154  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
155  }
156 
157  // 8. Set up the parallel linear form b(.) which corresponds to the
158  // right-hand side of the FEM linear system, which in this case is
159  // (f,phi_i) where f is given by the function f_exact and phi_i are the
160  // basis functions in the finite element fespace.
162  ParLinearForm *b = new ParLinearForm(fespace);
163  b->AddDomainIntegrator(new VectorFEDomainLFIntegrator(f));
164  b->Assemble();
165 
166  // 9. Define the solution vector x as a parallel finite element grid function
167  // corresponding to fespace. Initialize x by projecting the exact
168  // solution. Note that only values from the boundary faces will be used
169  // when eliminating the non-homogeneous boundary condition to modify the
170  // r.h.s. vector b.
171  ParGridFunction x(fespace);
173  x.ProjectCoefficient(F);
174 
175  // 10. Set up the parallel bilinear form corresponding to the H(div)
176  // diffusion operator grad alpha div + beta I, by adding the div-div and
177  // the mass domain integrators.
179  Coefficient *beta = new ConstantCoefficient(1.0);
180  ParBilinearForm *a = new ParBilinearForm(fespace);
181  a->AddDomainIntegrator(new DivDivIntegrator(*alpha));
182  a->AddDomainIntegrator(new VectorFEMassIntegrator(*beta));
183 
184  // 11. Assemble the parallel bilinear form and the corresponding linear
185  // system, applying any necessary transformations such as: parallel
186  // assembly, eliminating boundary conditions, applying conforming
187  // constraints for non-conforming AMR, static condensation,
188  // hybridization, etc.
189  FiniteElementCollection *hfec = NULL;
190  ParFiniteElementSpace *hfes = NULL;
191  if (static_cond)
192  {
193  a->EnableStaticCondensation();
194  }
195  else if (hybridization)
196  {
197  hfec = new DG_Interface_FECollection(order-1, dim);
198  hfes = new ParFiniteElementSpace(pmesh, hfec);
199  a->EnableHybridization(hfes, new NormalTraceJumpIntegrator(),
200  ess_tdof_list);
201  }
202  a->Assemble();
203 
204  Vector B, X;
205  CGSolver *pcg = new CGSolver(MPI_COMM_WORLD);
206  pcg->SetRelTol(1e-12);
207  pcg->SetMaxIter(500);
208  pcg->SetPrintLevel(1);
209 
210  if (!use_petsc)
211  {
212  HypreParMatrix A;
213  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
214 
215  HYPRE_BigInt glob_size = A.GetGlobalNumRows();
216  if (myid == 0)
217  {
218  cout << "Size of linear system: " << glob_size << endl;
219  }
220 
221  // 12. Define and apply a parallel PCG solver for A X = B with the 2D AMS or
222  // the 3D ADS preconditioners from hypre. If using hybridization, the
223  // system is preconditioned with hypre's BoomerAMG.
224  HypreSolver *prec = NULL;
225  pcg->SetOperator(A);
226  if (hybridization) { prec = new HypreBoomerAMG(A); }
227  else
228  {
229  ParFiniteElementSpace *prec_fespace =
230  (a->StaticCondensationIsEnabled() ? a->SCParFESpace() : fespace);
231  if (dim == 2) { prec = new HypreAMS(A, prec_fespace); }
232  else { prec = new HypreADS(A, prec_fespace); }
233  }
234  pcg->SetPreconditioner(*prec);
235  pcg->Mult(B, X);
236  delete prec;
237  }
238  else
239  {
240  PetscParMatrix A;
241  PetscPreconditioner *prec = NULL;
242  a->SetOperatorType(use_nonoverlapping ?
243  Operator::PETSC_MATIS : Operator::PETSC_MATAIJ);
244  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
245 
246  if (myid == 0)
247  {
248  cout << "Size of linear system: " << A.M() << endl;
249  }
250 
251  pcg->SetOperator(A);
252  if (use_nonoverlapping)
253  {
254  ParFiniteElementSpace *prec_fespace =
255  (a->StaticCondensationIsEnabled() ? a->SCParFESpace() : fespace);
256 
257  // Auxiliary class for BDDC customization
259  // Inform the solver about the finite element space
260  opts.SetSpace(prec_fespace);
261  // Inform the solver about essential dofs
262  opts.SetEssBdrDofs(&ess_tdof_list);
263  // Create a BDDC solver with parameters
264  prec = new PetscBDDCSolver(A, opts);
265  }
266  else
267  {
268  // Create an empty preconditioner that can be customized at runtime.
269  prec = new PetscPreconditioner(A, "solver_");
270  }
271  pcg->SetPreconditioner(*prec);
272  pcg->Mult(B, X);
273  delete prec;
274  }
275  delete pcg;
276 
277  // 13. Recover the parallel grid function corresponding to X. This is the
278  // local finite element solution on each processor.
279  a->RecoverFEMSolution(X, *b, x);
280 
281  // 14. Compute and print the L^2 norm of the error.
282  {
283  double err = x.ComputeL2Error(F);
284  if (myid == 0)
285  {
286  cout << "\n|| F_h - F ||_{L^2} = " << err << '\n' << endl;
287  }
288  }
289 
290  // 15. Save the refined mesh and the solution in parallel. This output can
291  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
292  {
293  ostringstream mesh_name, sol_name;
294  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
295  sol_name << "sol." << setfill('0') << setw(6) << myid;
296 
297  ofstream mesh_ofs(mesh_name.str().c_str());
298  mesh_ofs.precision(8);
299  pmesh->Print(mesh_ofs);
300 
301  ofstream sol_ofs(sol_name.str().c_str());
302  sol_ofs.precision(8);
303  x.Save(sol_ofs);
304  }
305 
306  // 16. Send the solution by socket to a GLVis server.
307  if (visualization)
308  {
309  char vishost[] = "localhost";
310  int visport = 19916;
311  socketstream sol_sock(vishost, visport);
312  sol_sock << "parallel " << num_procs << " " << myid << "\n";
313  sol_sock.precision(8);
314  sol_sock << "solution\n" << *pmesh << x << flush;
315  }
316 
317  // 17. Free the used memory.
318  delete hfes;
319  delete hfec;
320  delete a;
321  delete alpha;
322  delete beta;
323  delete b;
324  delete fespace;
325  delete fec;
326  delete pmesh;
327 
328  // We finalize PETSc
329  if (use_petsc) { MFEMFinalizePetsc(); }
330 
331  return 0;
332 }
333 
334 
335 // The exact solution (for non-surface meshes)
336 void F_exact(const Vector &p, Vector &F)
337 {
338  int dim = p.Size();
339 
340  double x = p(0);
341  double y = p(1);
342  // double z = (dim == 3) ? p(2) : 0.0;
343 
344  F(0) = cos(kappa*x)*sin(kappa*y);
345  F(1) = cos(kappa*y)*sin(kappa*x);
346  if (dim == 3)
347  {
348  F(2) = 0.0;
349  }
350 }
351 
352 // The right hand side
353 void f_exact(const Vector &p, Vector &f)
354 {
355  int dim = p.Size();
356 
357  double x = p(0);
358  double y = p(1);
359  // double z = (dim == 3) ? p(2) : 0.0;
360 
361  double temp = 1 + 2*kappa*kappa;
362 
363  f(0) = temp*cos(kappa*x)*sin(kappa*y);
364  f(1) = temp*cos(kappa*y)*sin(kappa*x);
365  if (dim == 3)
366  {
367  f(2) = 0;
368  }
369 }
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1032
Conjugate gradient method.
Definition: solvers.hpp:493
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:1743
Abstract class for PETSc&#39;s preconditioners.
Definition: petsc.hpp:801
The Auxiliary-space Divergence Solver in hypre.
Definition: hypre.hpp:1820
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Definition: globals.hpp:71
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:324
int Dimension() const
Definition: mesh.hpp:1047
Wrapper for PETSc&#39;s matrix class.
Definition: petsc.hpp:315
void f_exact(const Vector &, Vector &)
Definition: ex4p.cpp:327
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
int main(int argc, char *argv[])
Definition: ex4p.cpp:59
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:712
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition: array.cpp:68
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150
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
STL namespace.
PetscInt M() const
Returns the global number of rows.
Definition: petsc.cpp:929
(Q div u, div v) for RT elements
void F_exact(const Vector &, Vector &)
Definition: ex4p.cpp:310
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1590
Class for parallel linear form.
Definition: plinearform.hpp:26
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
constexpr char vishost[]
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9878
constexpr int visport
Auxiliary class for BDDC customization.
Definition: petsc.hpp:827
void SetMaxIter(int max_it)
Definition: solvers.hpp:201
double freq
Definition: ex4p.cpp:57
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Definition: pgridfunc.hpp:281
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:285
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:373
A general vector function coefficient.
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
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:41
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 SpaceDimension() const
Definition: mesh.hpp:1048
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:873
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:936
double a
Definition: lissajous.cpp:41
void MFEMFinalizePetsc()
Definition: petsc.cpp:192
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
Definition: petsc.cpp:168
void SetSpace(ParFiniteElementSpace *fe)
Definition: petsc.hpp:842
int dim
Definition: ex24.cpp:53
double kappa
Definition: ex4p.cpp:57
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:1102
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:346
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:507
const double alpha
Definition: ex15.cpp:369
Vector data type.
Definition: vector.hpp:60
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:173
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4839
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
Class for parallel meshes.
Definition: pmesh.hpp:32
void SetEssBdrDofs(const Array< int > *essdofs, bool loc=false)
Specify dofs on the essential boundary.
Definition: petsc.hpp:847
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:635
double f(const Vector &p)