MFEM  v3.1
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 //
3 // Compile with: make ex4p
4 //
5 // Sample runs: mpirun -np 4 ex4p -m ../data/square-disc.mesh
6 // mpirun -np 4 ex4p -m ../data/star.mesh
7 // mpirun -np 4 ex4p -m ../data/beam-tet.mesh
8 // mpirun -np 4 ex4p -m ../data/beam-hex.mesh
9 // mpirun -np 4 ex4p -m ../data/escher.mesh -o 2 -sc
10 // mpirun -np 4 ex4p -m ../data/fichera.mesh -o 2 -hb
11 // mpirun -np 4 ex4p -m ../data/fichera-q2.vtk
12 // mpirun -np 4 ex4p -m ../data/fichera-q3.mesh -o 2 -sc
13 // mpirun -np 4 ex4p -m ../data/square-disc-nurbs.mesh -o 3
14 // mpirun -np 4 ex4p -m ../data/beam-hex-nurbs.mesh -o 3
15 // mpirun -np 4 ex4p -m ../data/periodic-square.mesh -no-bc
16 // mpirun -np 4 ex4p -m ../data/periodic-cube.mesh -no-bc
17 // mpirun -np 4 ex4p -m ../data/amr-quad.mesh
18 // mpirun -np 4 ex4p -m ../data/amr-hex.mesh -o 2 -sc
19 // mpirun -np 4 ex4p -m ../data/star-surf.mesh -o 3 -hb
20 //
21 // Description: This example code solves a simple 2D/3D H(div) diffusion
22 // problem corresponding to the second order definite equation
23 // -grad(alpha div F) + beta F = f with boundary condition F dot n
24 // = <given normal field>. Here, we use a given exact solution F
25 // and compute the corresponding r.h.s. f. We discretize with
26 // Raviart-Thomas finite elements.
27 //
28 // The example demonstrates the use of H(div) finite element
29 // spaces with the grad-div and H(div) vector finite element mass
30 // bilinear form, as well as the computation of discretization
31 // error when the exact solution is known. Bilinear form
32 // hybridization and static condensation are also illustrated.
33 //
34 // We recommend viewing examples 1-3 before viewing this example.
35 
36 #include "mfem.hpp"
37 #include <fstream>
38 #include <iostream>
39 
40 using namespace std;
41 using namespace mfem;
42 
43 // Exact solution, F, and r.h.s., f. See below for implementation.
44 void F_exact(const Vector &, Vector &);
45 void f_exact(const Vector &, Vector &);
46 double freq = 1.0, kappa;
47 
48 int main(int argc, char *argv[])
49 {
50  // 1. Initialize MPI.
51  int num_procs, myid;
52  MPI_Init(&argc, &argv);
53  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
54  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
55 
56  // 2. Parse command-line options.
57  const char *mesh_file = "../data/star.mesh";
58  int order = 1;
59  bool set_bc = true;
60  bool static_cond = false;
61  bool hybridization = false;
62  bool visualization = 1;
63 
64  OptionsParser args(argc, argv);
65  args.AddOption(&mesh_file, "-m", "--mesh",
66  "Mesh file to use.");
67  args.AddOption(&order, "-o", "--order",
68  "Finite element order (polynomial degree).");
69  args.AddOption(&set_bc, "-bc", "--impose-bc", "-no-bc", "--dont-impose-bc",
70  "Impose or not essential boundary conditions.");
71  args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
72  " solution.");
73  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
74  "--no-static-condensation", "Enable static condensation.");
75  args.AddOption(&hybridization, "-hb", "--hybridization", "-no-hb",
76  "--no-hybridization", "Enable hybridization.");
77  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
78  "--no-visualization",
79  "Enable or disable GLVis visualization.");
80  args.Parse();
81  if (!args.Good())
82  {
83  if (myid == 0)
84  {
85  args.PrintUsage(cout);
86  }
87  MPI_Finalize();
88  return 1;
89  }
90  if (myid == 0)
91  {
92  args.PrintOptions(cout);
93  }
94  kappa = freq * M_PI;
95 
96  // 3. Read the (serial) mesh from the given mesh file on all processors. We
97  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
98  // and volume, as well as periodic meshes with the same code.
99  Mesh *mesh;
100  ifstream imesh(mesh_file);
101  if (!imesh)
102  {
103  if (myid == 0)
104  {
105  cerr << "\nCan not open mesh file: " << mesh_file << '\n' << endl;
106  }
107  MPI_Finalize();
108  return 2;
109  }
110  mesh = new Mesh(imesh, 1, 1);
111  imesh.close();
112  int dim = mesh->Dimension();
113  int sdim = mesh->SpaceDimension();
114 
115  // 4. Refine the serial mesh on all processors to increase the resolution. In
116  // this example we do 'ref_levels' of uniform refinement. We choose
117  // 'ref_levels' to be the largest number that gives a final mesh with no
118  // more than 1,000 elements.
119  {
120  int ref_levels =
121  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
122  for (int l = 0; l < ref_levels; l++)
123  {
124  mesh->UniformRefinement();
125  }
126  }
127 
128  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
129  // this mesh further in parallel to increase the resolution. Once the
130  // parallel mesh is defined, the serial mesh can be deleted. Tetrahedral
131  // meshes need to be reoriented before we can define high-order Nedelec
132  // spaces on them (this is needed in the ADS solver below).
133  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
134  delete mesh;
135  {
136  int par_ref_levels = 2;
137  for (int l = 0; l < par_ref_levels; l++)
138  {
139  pmesh->UniformRefinement();
140  }
141  }
142  pmesh->ReorientTetMesh();
143 
144  // 6. Define a parallel finite element space on the parallel mesh. Here we
145  // use the Raviart-Thomas finite elements of the specified order.
146  FiniteElementCollection *fec = new RT_FECollection(order-1, dim);
147  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
148  HYPRE_Int size = fespace->GlobalTrueVSize();
149  if (myid == 0)
150  {
151  cout << "Number of finite element unknowns: " << size << endl;
152  }
153 
154  // 7. Determine the list of true (i.e. parallel conforming) essential
155  // boundary dofs. In this example, the boundary conditions are defined
156  // by marking all the boundary attributes from the mesh as essential
157  // (Dirichlet) and converting them to a list of true dofs.
158  Array<int> ess_tdof_list;
159  if (pmesh->bdr_attributes.Size())
160  {
161  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
162  ess_bdr = set_bc ? 1 : 0;
163  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
164  }
165 
166  // 8. Set up the parallel linear form b(.) which corresponds to the
167  // right-hand side of the FEM linear system, which in this case is
168  // (f,phi_i) where f is given by the function f_exact and phi_i are the
169  // basis functions in the finite element fespace.
171  ParLinearForm *b = new ParLinearForm(fespace);
173  b->Assemble();
174 
175  // 9. Define the solution vector x as a parallel finite element grid function
176  // corresponding to fespace. Initialize x by projecting the exact
177  // solution. Note that only values from the boundary faces will be used
178  // when eliminating the non-homogeneous boundary condition to modify the
179  // r.h.s. vector b.
180  ParGridFunction x(fespace);
182  x.ProjectCoefficient(F);
183 
184  // 10. Set up the parallel bilinear form corresponding to the H(div)
185  // diffusion operator grad alpha div + beta I, by adding the div-div and
186  // the mass domain integrators.
187  Coefficient *alpha = new ConstantCoefficient(1.0);
188  Coefficient *beta = new ConstantCoefficient(1.0);
189  ParBilinearForm *a = new ParBilinearForm(fespace);
190  a->AddDomainIntegrator(new DivDivIntegrator(*alpha));
192 
193  // 11. Assemble the parallel bilinear form and the corresponding linear
194  // system, applying any necessary transformations such as: parallel
195  // assembly, eliminating boundary conditions, applying conforming
196  // constraints for non-conforming AMR, static condensation,
197  // hybridization, etc.
198  FiniteElementCollection *hfec = NULL;
199  ParFiniteElementSpace *hfes = NULL;
200  if (static_cond)
201  {
203  }
204  else if (hybridization)
205  {
206  hfec = new DG_Interface_FECollection(order-1, dim);
207  hfes = new ParFiniteElementSpace(pmesh, hfec);
209  ess_tdof_list);
210  }
211  a->Assemble();
212 
213  HypreParMatrix A;
214  Vector B, X;
215  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
216 
217  HYPRE_Int glob_size = A.GetGlobalNumRows();
218  if (myid == 0)
219  {
220  cout << "Size of linear system: " << glob_size << endl;
221  }
222 
223  // 12. Define and apply a parallel PCG solver for A X = B with the 2D AMS or
224  // the 3D ADS preconditioners from hypre. If using hybridization, the
225  // system is preconditioned with hypre's BoomerAMG.
226  HypreSolver *prec = NULL;
227  CGSolver *pcg = new CGSolver(A.GetComm());
228  pcg->SetOperator(A);
229  pcg->SetRelTol(1e-12);
230  pcg->SetMaxIter(500);
231  pcg->SetPrintLevel(1);
232  if (hybridization) { prec = new HypreBoomerAMG(A); }
233  else
234  {
235  ParFiniteElementSpace *prec_fespace =
236  (a->StaticCondensationIsEnabled() ? a->SCParFESpace() : fespace);
237  if (dim == 2) { prec = new HypreAMS(A, prec_fespace); }
238  else { prec = new HypreADS(A, prec_fespace); }
239  }
240  pcg->SetPreconditioner(*prec);
241  pcg->Mult(B, X);
242 
243  // 13. Recover the parallel grid function corresponding to X. This is the
244  // local finite element solution on each processor.
245  a->RecoverFEMSolution(X, *b, x);
246 
247  // 14. Compute and print the L^2 norm of the error.
248  {
249  double err = x.ComputeL2Error(F);
250  if (myid == 0)
251  {
252  cout << "\n|| F_h - F ||_{L^2} = " << err << '\n' << endl;
253  }
254  }
255 
256  // 15. Save the refined mesh and the solution in parallel. This output can
257  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
258  {
259  ostringstream mesh_name, sol_name;
260  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
261  sol_name << "sol." << setfill('0') << setw(6) << myid;
262 
263  ofstream mesh_ofs(mesh_name.str().c_str());
264  mesh_ofs.precision(8);
265  pmesh->Print(mesh_ofs);
266 
267  ofstream sol_ofs(sol_name.str().c_str());
268  sol_ofs.precision(8);
269  x.Save(sol_ofs);
270  }
271 
272  // 16. Send the solution by socket to a GLVis server.
273  if (visualization)
274  {
275  char vishost[] = "localhost";
276  int visport = 19916;
277  socketstream sol_sock(vishost, visport);
278  sol_sock << "parallel " << num_procs << " " << myid << "\n";
279  sol_sock.precision(8);
280  sol_sock << "solution\n" << *pmesh << x << flush;
281  }
282 
283  // 17. Free the used memory.
284  delete pcg;
285  delete prec;
286  delete hfes;
287  delete hfec;
288  delete a;
289  delete alpha;
290  delete beta;
291  delete b;
292  delete fespace;
293  delete fec;
294  delete pmesh;
295 
296  MPI_Finalize();
297 
298  return 0;
299 }
300 
301 
302 // The exact solution (for non-surface meshes)
303 void F_exact(const Vector &p, Vector &F)
304 {
305  int dim = p.Size();
306 
307  double x = p(0);
308  double y = p(1);
309  // double z = (dim == 3) ? p(2) : 0.0;
310 
311  F(0) = cos(kappa*x)*sin(kappa*y);
312  F(1) = cos(kappa*y)*sin(kappa*x);
313  if (dim == 3)
314  {
315  F(2) = 0.0;
316  }
317 }
318 
319 // The right hand side
320 void f_exact(const Vector &p, Vector &f)
321 {
322  int dim = p.Size();
323 
324  double x = p(0);
325  double y = p(1);
326  // double z = (dim == 3) ? p(2) : 0.0;
327 
328  double temp = 1 + 2*kappa*kappa;
329 
330  f(0) = temp*cos(kappa*x)*sin(kappa*y);
331  f(1) = temp*cos(kappa*y)*sin(kappa*x);
332  if (dim == 3)
333  {
334  f(2) = 0;
335  }
336 }
int Size() const
Logical size of the array.
Definition: array.hpp:109
Conjugate gradient method.
Definition: solvers.hpp:109
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:261
ParFiniteElementSpace * SCParFESpace() const
Return the parallel trace FE space associated with static condensation.
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:761
double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:148
The Auxiliary-space Divergence Solver in hypre.
Definition: hypre.hpp:793
Subclass constant coefficient.
Definition: coefficient.hpp:57
virtual void ReorientTetMesh()
See the remarks for the serial version in mesh.hpp.
Definition: pmesh.cpp:1432
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:34
int Size() const
Returns the size of the vector.
Definition: vector.hpp:84
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:454
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:306
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:241
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
Definition: solvers.cpp:294
(Q div u, div v) for RT elements
Definition: bilininteg.hpp:519
void EnableHybridization(FiniteElementSpace *constr_space, BilinearFormIntegrator *constr_integ, const Array< int > &ess_tdof_list)
The BoomerAMG solver in hypre.
Definition: hypre.hpp:698
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_Int GetGlobalNumRows() const
Definition: hypre.hpp:333
int dim
Definition: ex3.cpp:48
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7316
void SetMaxIter(int max_it)
Definition: solvers.hpp:61
T Max() const
Definition: array.cpp:90
void F_exact(const Vector &, Vector &)
Definition: ex4.cpp:239
void Assemble(int skip_zeros=1)
Assemble the local matrix.
int Dimension() const
Definition: mesh.hpp:475
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:133
int SpaceDimension() const
Definition: mesh.hpp:476
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
bool StaticCondensationIsEnabled() const
Array< int > bdr_attributes
Definition: mesh.hpp:140
int main(int argc, char *argv[])
Definition: ex1.cpp:45
void SetRelTol(double rtol)
Definition: solvers.hpp:59
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list)
Definition: pfespace.cpp:545
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)
Definition: optparser.hpp:74
void f_exact(const Vector &, Vector &)
Definition: ex3.cpp:236
void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
void FormLinearSystem(Array< int > &ess_tdof_list, Vector &x, Vector &b, HypreParMatrix &A, Vector &X, Vector &B, int copy_interior=0)
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:522
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:171
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:123
double kappa
Definition: ex3.cpp:47
Vector data type.
Definition: vector.hpp:33
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:92
Class for parallel grid function.
Definition: pgridfunc.hpp:31
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:143
double freq
Definition: ex3.cpp:47
Class for parallel meshes.
Definition: pmesh.hpp:28
Integrator for (Q u, v) for VectorFiniteElements.
Definition: bilininteg.hpp:458
void EnableStaticCondensation()
virtual void Print(std::ostream &out=std::cout) const
Definition: pmesh.cpp:2850
bool Good() const
Definition: optparser.hpp:120