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