MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
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.
44  int num_procs, myid;
45  MPI_Init(&argc, &argv);
46  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
47  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
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  MPI_Finalize();
93  return 1;
94  }
95  if (myid == 0)
96  {
97  args.PrintOptions(cout);
98  }
99  // 2b. We initialize PETSc
100  if (use_petsc) { MFEMInitializePetsc(NULL,NULL,petscrc_file,NULL); }
101  kappa = freq * M_PI;
102 
103  // 3. Read the (serial) mesh from the given mesh file on all processors. We
104  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
105  // and volume, as well as periodic meshes with the same code.
106  Mesh *mesh = new Mesh(mesh_file, 1, 1);
107  int dim = mesh->Dimension();
108  int sdim = mesh->SpaceDimension();
109 
110  // 4. Refine the serial mesh on all processors to increase the resolution. In
111  // this example we do 'ref_levels' of uniform refinement. We choose
112  // 'ref_levels' to be the largest number that gives a final mesh with no
113  // more than 1,000 elements.
114  {
115  int ref_levels =
116  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
117  for (int l = 0; l < ref_levels; l++)
118  {
119  mesh->UniformRefinement();
120  }
121  }
122 
123  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
124  // this mesh further in parallel to increase the resolution. Once the
125  // parallel mesh is defined, the serial mesh can be deleted. Tetrahedral
126  // meshes need to be reoriented before we can define high-order Nedelec
127  // spaces on them (this is needed in the ADS solver below).
128  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
129  delete mesh;
130  {
131  int par_ref_levels = 2;
132  for (int l = 0; l < par_ref_levels; l++)
133  {
134  pmesh->UniformRefinement();
135  }
136  }
137  pmesh->ReorientTetMesh();
138 
139  // 6. Define a parallel finite element space on the parallel mesh. Here we
140  // use the Raviart-Thomas finite elements of the specified order.
141  FiniteElementCollection *fec = new RT_FECollection(order-1, dim);
142  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
143  HYPRE_BigInt size = fespace->GlobalTrueVSize();
144  if (myid == 0)
145  {
146  cout << "Number of finite element unknowns: " << size << endl;
147  }
148 
149  // 7. Determine the list of true (i.e. parallel conforming) essential
150  // boundary dofs. In this example, the boundary conditions are defined
151  // by marking all the boundary attributes from the mesh as essential
152  // (Dirichlet) and converting them to a list of true dofs.
153  Array<int> ess_tdof_list;
154  if (pmesh->bdr_attributes.Size())
155  {
156  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
157  ess_bdr = set_bc ? 1 : 0;
158  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
159  }
160 
161  // 8. Set up the parallel linear form b(.) which corresponds to the
162  // right-hand side of the FEM linear system, which in this case is
163  // (f,phi_i) where f is given by the function f_exact and phi_i are the
164  // basis functions in the finite element fespace.
166  ParLinearForm *b = new ParLinearForm(fespace);
168  b->Assemble();
169 
170  // 9. Define the solution vector x as a parallel finite element grid function
171  // corresponding to fespace. Initialize x by projecting the exact
172  // solution. Note that only values from the boundary faces will be used
173  // when eliminating the non-homogeneous boundary condition to modify the
174  // r.h.s. vector b.
175  ParGridFunction x(fespace);
177  x.ProjectCoefficient(F);
178 
179  // 10. Set up the parallel bilinear form corresponding to the H(div)
180  // diffusion operator grad alpha div + beta I, by adding the div-div and
181  // the mass domain integrators.
183  Coefficient *beta = new ConstantCoefficient(1.0);
184  ParBilinearForm *a = new ParBilinearForm(fespace);
185  a->AddDomainIntegrator(new DivDivIntegrator(*alpha));
187 
188  // 11. Assemble the parallel bilinear form and the corresponding linear
189  // system, applying any necessary transformations such as: parallel
190  // assembly, eliminating boundary conditions, applying conforming
191  // constraints for non-conforming AMR, static condensation,
192  // hybridization, etc.
193  FiniteElementCollection *hfec = NULL;
194  ParFiniteElementSpace *hfes = NULL;
195  if (static_cond)
196  {
198  }
199  else if (hybridization)
200  {
201  hfec = new DG_Interface_FECollection(order-1, dim);
202  hfes = new ParFiniteElementSpace(pmesh, hfec);
204  ess_tdof_list);
205  }
206  a->Assemble();
207 
208  Vector B, X;
209  CGSolver *pcg = new CGSolver(MPI_COMM_WORLD);
210  pcg->SetRelTol(1e-12);
211  pcg->SetMaxIter(500);
212  pcg->SetPrintLevel(1);
213 
214  if (!use_petsc)
215  {
216  HypreParMatrix A;
217  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
218 
219  HYPRE_BigInt glob_size = A.GetGlobalNumRows();
220  if (myid == 0)
221  {
222  cout << "Size of linear system: " << glob_size << endl;
223  }
224 
225  // 12. Define and apply a parallel PCG solver for A X = B with the 2D AMS or
226  // the 3D ADS preconditioners from hypre. If using hybridization, the
227  // system is preconditioned with hypre's BoomerAMG.
228  HypreSolver *prec = NULL;
229  pcg->SetOperator(A);
230  if (hybridization) { prec = new HypreBoomerAMG(A); }
231  else
232  {
233  ParFiniteElementSpace *prec_fespace =
234  (a->StaticCondensationIsEnabled() ? a->SCParFESpace() : fespace);
235  if (dim == 2) { prec = new HypreAMS(A, prec_fespace); }
236  else { prec = new HypreADS(A, prec_fespace); }
237  }
238  pcg->SetPreconditioner(*prec);
239  pcg->Mult(B, X);
240  delete prec;
241  }
242  else
243  {
244  PetscParMatrix A;
245  PetscPreconditioner *prec = NULL;
246  a->SetOperatorType(use_nonoverlapping ?
247  Operator::PETSC_MATIS : Operator::PETSC_MATAIJ);
248  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
249 
250  if (myid == 0)
251  {
252  cout << "Size of linear system: " << A.M() << endl;
253  }
254 
255  pcg->SetOperator(A);
256  if (use_nonoverlapping)
257  {
258  ParFiniteElementSpace *prec_fespace =
259  (a->StaticCondensationIsEnabled() ? a->SCParFESpace() : fespace);
260 
261  // Auxiliary class for BDDC customization
263  // Inform the solver about the finite element space
264  opts.SetSpace(prec_fespace);
265  // Inform the solver about essential dofs
266  opts.SetEssBdrDofs(&ess_tdof_list);
267  // Create a BDDC solver with parameters
268  prec = new PetscBDDCSolver(A, opts);
269  }
270  else
271  {
272  // Create an empty preconditioner that can be customized at runtime.
273  prec = new PetscPreconditioner(A, "solver_");
274  }
275  pcg->SetPreconditioner(*prec);
276  pcg->Mult(B, X);
277  delete prec;
278  }
279  delete pcg;
280 
281  // 13. Recover the parallel grid function corresponding to X. This is the
282  // local finite element solution on each processor.
283  a->RecoverFEMSolution(X, *b, x);
284 
285  // 14. Compute and print the L^2 norm of the error.
286  {
287  double err = x.ComputeL2Error(F);
288  if (myid == 0)
289  {
290  cout << "\n|| F_h - F ||_{L^2} = " << err << '\n' << endl;
291  }
292  }
293 
294  // 15. Save the refined mesh and the solution in parallel. This output can
295  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
296  {
297  ostringstream mesh_name, sol_name;
298  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
299  sol_name << "sol." << setfill('0') << setw(6) << myid;
300 
301  ofstream mesh_ofs(mesh_name.str().c_str());
302  mesh_ofs.precision(8);
303  pmesh->Print(mesh_ofs);
304 
305  ofstream sol_ofs(sol_name.str().c_str());
306  sol_ofs.precision(8);
307  x.Save(sol_ofs);
308  }
309 
310  // 16. Send the solution by socket to a GLVis server.
311  if (visualization)
312  {
313  char vishost[] = "localhost";
314  int visport = 19916;
315  socketstream sol_sock(vishost, visport);
316  sol_sock << "parallel " << num_procs << " " << myid << "\n";
317  sol_sock.precision(8);
318  sol_sock << "solution\n" << *pmesh << x << flush;
319  }
320 
321  // 17. Free the used memory.
322  delete hfes;
323  delete hfec;
324  delete a;
325  delete alpha;
326  delete beta;
327  delete b;
328  delete fespace;
329  delete fec;
330  delete pmesh;
331 
332  // We finalize PETSc
333  if (use_petsc) { MFEMFinalizePetsc(); }
334 
335  MPI_Finalize();
336 
337  return 0;
338 }
339 
340 
341 // The exact solution (for non-surface meshes)
342 void F_exact(const Vector &p, Vector &F)
343 {
344  int dim = p.Size();
345 
346  double x = p(0);
347  double y = p(1);
348  // double z = (dim == 3) ? p(2) : 0.0;
349 
350  F(0) = cos(kappa*x)*sin(kappa*y);
351  F(1) = cos(kappa*y)*sin(kappa*x);
352  if (dim == 3)
353  {
354  F(2) = 0.0;
355  }
356 }
357 
358 // The right hand side
359 void f_exact(const Vector &p, Vector &f)
360 {
361  int dim = p.Size();
362 
363  double x = p(0);
364  double y = p(1);
365  // double z = (dim == 3) ? p(2) : 0.0;
366 
367  double temp = 1 + 2*kappa*kappa;
368 
369  f(0) = temp*cos(kappa*x)*sin(kappa*y);
370  f(1) = temp*cos(kappa*y)*sin(kappa*x);
371  if (dim == 3)
372  {
373  f(2) = 0;
374  }
375 }
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:790
Conjugate gradient method.
Definition: solvers.hpp:316
ParFiniteElementSpace * SCParFESpace() const
Return the parallel trace FE space associated with static condensation.
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:1540
Abstract class for PETSc&#39;s preconditioners.
Definition: petsc.hpp:775
The Auxiliary-space Divergence Solver in hypre.
Definition: hypre.hpp:1579
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:616
Wrapper for PETSc&#39;s matrix class.
Definition: petsc.hpp:307
virtual void ReorientTetMesh()
See the remarks for the serial version in mesh.hpp.
Definition: pmesh.cpp:2798
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:275
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
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
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:493
(Q div u, div v) for RT elements
PetscInt M() const
Returns the global number of rows.
Definition: petsc.cpp:928
void EnableHybridization(FiniteElementSpace *constr_space, BilinearFormIntegrator *constr_integ, const Array< int > &ess_tdof_list)
Enable hybridization.
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1387
void SetOperatorType(Operator::Type tid)
Set the operator type id for the parallel matrix/operator when using AssemblyLevel::LEGACY.
Class for parallel linear form.
Definition: plinearform.hpp:26
void SetPrintLevel(int print_lvl)
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: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
Auxiliary class for BDDC customization.
Definition: petsc.hpp:801
void SetMaxIter(int max_it)
Definition: solvers.hpp:100
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
void F_exact(const Vector &, Vector &)
Definition: ex4.cpp:266
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4382
int Dimension() const
Definition: mesh.hpp:911
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:341
A general vector function coefficient.
int SpaceDimension() const
Definition: mesh.hpp:912
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
bool StaticCondensationIsEnabled() const
Check if static condensation was actually enabled by a previous call to EnableStaticCondensation().
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:204
void SetRelTol(double rtol)
Definition: solvers.hpp:98
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
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
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
void f_exact(const Vector &, Vector &)
Definition: ex3.cpp:257
HYPRE_Int HYPRE_BigInt
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
double a
Definition: lissajous.cpp:41
void MFEMFinalizePetsc()
Definition: petsc.cpp:192
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:273
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:573
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
Definition: petsc.cpp:168
void SetSpace(ParFiniteElementSpace *fe)
Definition: petsc.hpp:816
int dim
Definition: ex24.cpp:53
double freq
Definition: ex24.cpp:54
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
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:266
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:330
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:92
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:277
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:821
int main()
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150