MFEM  v3.0
 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
10 // mpirun -np 4 ex4p -m ../data/fichera.mesh
11 // mpirun -np 4 ex4p -m ../data/fichera-q2.vtk
12 // mpirun -np 4 ex4p -m ../data/fichera-q3.mesh
13 // mpirun -np 4 ex4p -m ../data/square-disc-nurbs.mesh
14 // mpirun -np 4 ex4p -m ../data/beam-hex-nurbs.mesh
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 //
18 // Description: This example code solves a simple 2D/3D H(div) diffusion
19 // problem corresponding to the second order definite equation
20 // -grad(alpha div F) + beta F = f with boundary condition F dot n
21 // = <given normal field>. Here, we use a given exact solution F
22 // and compute the corresponding r.h.s. f. We discretize with
23 // Raviart-Thomas finite elements.
24 //
25 // The example demonstrates the use of H(div) finite element
26 // spaces with the grad-div and H(div) vector finite element mass
27 // bilinear form, as well as the computation of discretization
28 // error when the exact solution is known.
29 //
30 // We recommend viewing examples 1-3 before viewing this example.
31 
32 #include "mfem.hpp"
33 #include <fstream>
34 #include <iostream>
35 
36 using namespace std;
37 using namespace mfem;
38 
39 // Exact solution, F, and r.h.s., f. See below for implementation.
40 void F_exact(const Vector &, Vector &);
41 void f_exact(const Vector &, Vector &);
42 
43 int main(int argc, char *argv[])
44 {
45  // 1. Initialize MPI.
46  int num_procs, myid;
47  MPI_Init(&argc, &argv);
48  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
49  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
50 
51  // 2. Parse command-line options.
52  const char *mesh_file = "../data/star.mesh";
53  int order = 1;
54  bool set_bc = true;
55  bool visualization = 1;
56 
57  OptionsParser args(argc, argv);
58  args.AddOption(&mesh_file, "-m", "--mesh",
59  "Mesh file to use.");
60  args.AddOption(&order, "-o", "--order",
61  "Finite element order (polynomial degree).");
62  args.AddOption(&set_bc, "-bc", "--impose-bc", "-no-bc", "--dont-impose-bc",
63  "Impose or not essential boundary conditions.");
64  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
65  "--no-visualization",
66  "Enable or disable GLVis visualization.");
67  args.Parse();
68  if (!args.Good())
69  {
70  if (myid == 0)
71  args.PrintUsage(cout);
72  MPI_Finalize();
73  return 1;
74  }
75  if (myid == 0)
76  args.PrintOptions(cout);
77 
78  // 3. Read the (serial) mesh from the given mesh file on all processors. We
79  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
80  // and volume, as well as periodic meshes with the same code.
81  Mesh *mesh;
82  ifstream imesh(mesh_file);
83  if (!imesh)
84  {
85  if (myid == 0)
86  cerr << "\nCan not open mesh file: " << mesh_file << '\n' << endl;
87  MPI_Finalize();
88  return 2;
89  }
90  mesh = new Mesh(imesh, 1, 1);
91  imesh.close();
92  int dim = mesh->Dimension();
93 
94  // 4. Refine the serial mesh on all processors to increase the resolution. In
95  // this example we do 'ref_levels' of uniform refinement. We choose
96  // 'ref_levels' to be the largest number that gives a final mesh with no
97  // more than 1,000 elements.
98  {
99  int ref_levels =
100  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
101  for (int l = 0; l < ref_levels; l++)
102  mesh->UniformRefinement();
103  }
104 
105  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
106  // this mesh further in parallel to increase the resolution. Once the
107  // parallel mesh is defined, the serial mesh can be deleted. Tetrahedral
108  // meshes need to be reoriented before we can define high-order Nedelec
109  // spaces on them (this is needed in the ADS solver below).
110  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
111  delete mesh;
112  {
113  int par_ref_levels = 2;
114  for (int l = 0; l < par_ref_levels; l++)
115  pmesh->UniformRefinement();
116  }
117  pmesh->ReorientTetMesh();
118 
119  // 6. Define a parallel finite element space on the parallel mesh. Here we
120  // use the lowest order Raviart-Thomas finite elements, but we can easily
121  // switch to higher-order spaces by changing the value of p.
122  FiniteElementCollection *fec = new RT_FECollection(order-1, dim);
123  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
124  int size = fespace->GlobalTrueVSize();
125  if (myid == 0)
126  cout << "Number of unknowns: " << size << endl;
127 
128  // 7. Set up the parallel linear form b(.) which corresponds to the
129  // right-hand side of the FEM linear system, which in this case is
130  // (f,phi_i) where f is given by the function f_exact and phi_i are the
131  // basis functions in the finite element fespace.
133  ParLinearForm *b = new ParLinearForm(fespace);
135  b->Assemble();
136 
137  // 8. Define the solution vector x as a parallel finite element grid function
138  // corresponding to fespace. Initialize x by projecting the exact
139  // solution. Note that only values from the boundary faces will be used
140  // when eliminating the non-homogeneous boundary condition to modify the
141  // r.h.s. vector b.
142  ParGridFunction x(fespace);
144  x.ProjectCoefficient(F);
145 
146  // 9. Set up the parallel bilinear form corresponding to the H(div) diffusion
147  // operator grad alpha div + beta I, by adding the div-div and the
148  // mass domain integrators and finally imposing non-homogeneous Dirichlet
149  // boundary conditions. The boundary conditions are implemented by
150  // marking all the boundary attributes from the mesh as essential
151  // (Dirichlet). After serial and parallel assembly we extract the
152  // parallel matrix A.
153  Coefficient *alpha = new ConstantCoefficient(1.0);
154  Coefficient *beta = new ConstantCoefficient(1.0);
155  ParBilinearForm *a = new ParBilinearForm(fespace);
156  a->AddDomainIntegrator(new DivDivIntegrator(*alpha));
158  a->Assemble();
159  if (set_bc && pmesh->bdr_attributes.Size())
160  {
161  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
162  ess_bdr = 1;
163  a->EliminateEssentialBC(ess_bdr, x, *b);
164  }
165  a->Finalize();
166 
167  // 10. Define the parallel (hypre) matrix and vectors representing a(.,.),
168  // b(.) and the finite element approximation.
172  *X = 0.0;
173 
174  delete a;
175  delete alpha;
176  delete beta;
177  delete b;
178 
179  // 11. Define and apply a parallel PCG solver for AX=B with the 2D AMS or the
180  // 3D ADS preconditioners from hypre.
181  HypreSolver *prec;
182  if (dim == 2)
183  prec = new HypreAMS(*A, fespace);
184  else
185  prec = new HypreADS(*A, fespace);
186  HyprePCG *pcg = new HyprePCG(*A);
187  pcg->SetTol(1e-10);
188  pcg->SetMaxIter(500);
189  pcg->SetPrintLevel(2);
190  pcg->SetPreconditioner(*prec);
191  pcg->Mult(*B, *X);
192 
193  // 12. Extract the parallel grid function corresponding to the finite element
194  // approximation X. This is the local solution on each processor.
195  x = *X;
196 
197  // 13. Compute and print the L^2 norm of the error.
198  {
199  double err = x.ComputeL2Error(F);
200  if (myid == 0)
201  cout << "\n|| F_h - F ||_{L^2} = " << err << '\n' << endl;
202  }
203 
204  // 14. Save the refined mesh and the solution in parallel. This output can
205  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
206  {
207  ostringstream mesh_name, sol_name;
208  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
209  sol_name << "sol." << setfill('0') << setw(6) << myid;
210 
211  ofstream mesh_ofs(mesh_name.str().c_str());
212  mesh_ofs.precision(8);
213  pmesh->Print(mesh_ofs);
214 
215  ofstream sol_ofs(sol_name.str().c_str());
216  sol_ofs.precision(8);
217  x.Save(sol_ofs);
218  }
219 
220  // 15. Send the solution by socket to a GLVis server.
221  if (visualization)
222  {
223  char vishost[] = "localhost";
224  int visport = 19916;
225  socketstream sol_sock(vishost, visport);
226  sol_sock << "parallel " << num_procs << " " << myid << "\n";
227  sol_sock.precision(8);
228  sol_sock << "solution\n" << *pmesh << x << flush;
229  }
230 
231  // 16. Free the used memory.
232  delete pcg;
233  delete prec;
234  delete X;
235  delete B;
236  delete A;
237  delete fespace;
238  delete fec;
239  delete pmesh;
240 
241  MPI_Finalize();
242 
243  return 0;
244 }
245 
246 
247 // The exact solution
248 void F_exact(const Vector &p, Vector &F)
249 {
250  int dim = p.Size();
251 
252  double x = p(0);
253  double y = p(1);
254  // double z = (dim == 3) ? p(2) : 0.0;
255 
256  F(0) = cos(M_PI*x)*sin(M_PI*y);
257  F(1) = cos(M_PI*y)*sin(M_PI*x);
258  if (dim == 3)
259  F(2) = 0.0;
260 }
261 
262 // The right hand side
263 void f_exact(const Vector &p, Vector &f)
264 {
265  int dim = p.Size();
266 
267  double x = p(0);
268  double y = p(1);
269  // double z = (dim == 3) ? p(2) : 0.0;
270 
271  double temp = 1 + 2*M_PI*M_PI;
272 
273  f(0) = temp*cos(M_PI*x)*sin(M_PI*y);
274  f(1) = temp*cos(M_PI*y)*sin(M_PI*x);
275  if (dim == 3)
276  f(2) = 0;
277 }
void SetTol(double tol)
Definition: hypre.cpp:1339
int Size() const
Logical size of the array.
Definition: array.hpp:108
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:554
double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:128
The Auxiliary-space Divergence Solver in hypre.
Definition: hypre.hpp:581
Subclass constant coefficient.
Definition: coefficient.hpp:57
virtual void ReorientTetMesh()
See the remarks for the serial version in mesh.hpp.
Definition: pmesh.cpp:1200
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:76
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:396
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:250
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:229
(Q div u, div v) for RT elements
Definition: bilininteg.hpp:494
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:1354
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
Class for parallel linear form.
Definition: plinearform.hpp:26
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6225
T Max() const
Definition: array.cpp:78
void F_exact(const Vector &, Vector &)
Definition: ex4.cpp:192
void Assemble(int skip_zeros=1)
Assemble the local matrix.
int Dimension() const
Definition: mesh.hpp:417
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:385
void SetMaxIter(int max_iter)
Definition: hypre.cpp:1344
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:119
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:38
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Array< int > bdr_attributes
Definition: mesh.hpp:305
int main(int argc, char *argv[])
Definition: ex1.cpp:39
PCG solver in hypre.
Definition: hypre.hpp:381
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:195
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:1360
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:266
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:345
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:171
Vector data type.
Definition: vector.hpp:29
Class for parallel grid function.
Definition: pgridfunc.hpp:31
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:103
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:1377
Class for parallel meshes.
Definition: pmesh.hpp:27
void ParallelAverage(Vector &tv) const
Returns the vector averaged on the true dofs.
Definition: pgridfunc.cpp:110
void EliminateEssentialBC(Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs, int d=0)
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:34
Integrator for (Q u, v) for VectorFiniteElements.
Definition: bilininteg.hpp:434
virtual void Print(std::ostream &out=std::cout) const
Definition: pmesh.cpp:2444
bool Good() const
Definition: optparser.hpp:120