MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex3p.cpp
Go to the documentation of this file.
1 // MFEM Example 3 - Parallel Version
2 //
3 // Compile with: make ex3p
4 //
5 // Sample runs: mpirun -np 4 ex3p -m ../data/star.mesh
6 // mpirun -np 4 ex3p -m ../data/square-disc.mesh -o 2
7 // mpirun -np 4 ex3p -m ../data/beam-tet.mesh
8 // mpirun -np 4 ex3p -m ../data/beam-hex.mesh
9 // mpirun -np 4 ex3p -m ../data/escher.mesh
10 // mpirun -np 4 ex3p -m ../data/escher.mesh -o 2
11 // mpirun -np 4 ex3p -m ../data/fichera.mesh
12 // mpirun -np 4 ex3p -m ../data/fichera-q2.vtk
13 // mpirun -np 4 ex3p -m ../data/fichera-q3.mesh
14 // mpirun -np 4 ex3p -m ../data/square-disc-nurbs.mesh
15 // mpirun -np 4 ex3p -m ../data/beam-hex-nurbs.mesh
16 // mpirun -np 4 ex3p -m ../data/amr-quad.mesh -o 2
17 // mpirun -np 4 ex3p -m ../data/amr-hex.mesh
18 // mpirun -np 4 ex3p -m ../data/star-surf.mesh -o 2
19 // mpirun -np 4 ex3p -m ../data/mobius-strip.mesh -o 2 -f 0.1
20 // mpirun -np 4 ex3p -m ../data/klein-bottle.mesh -o 2 -f 0.1
21 //
22 // Description: This example code solves a simple electromagnetic diffusion
23 // problem corresponding to the second order definite Maxwell
24 // equation curl curl E + E = f with boundary condition
25 // E x n = <given tangential field>. Here, we use a given exact
26 // solution E and compute the corresponding r.h.s. f.
27 // We discretize with Nedelec finite elements in 2D or 3D.
28 //
29 // The example demonstrates the use of H(curl) finite element
30 // spaces with the curl-curl and the (vector finite element) mass
31 // bilinear form, as well as the computation of discretization
32 // error when the exact solution is known. Static condensation is
33 // also illustrated.
34 //
35 // We recommend viewing examples 1-2 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, E, and r.h.s., f. See below for implementation.
45 void E_exact(const Vector &, Vector &);
46 void f_exact(const Vector &, Vector &);
47 double freq = 1.0, kappa;
48 int dim;
49 
50 int main(int argc, char *argv[])
51 {
52  // 1. Initialize MPI.
53  int num_procs, myid;
54  MPI_Init(&argc, &argv);
55  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
56  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
57 
58  // 2. Parse command-line options.
59  const char *mesh_file = "../data/beam-tet.mesh";
60  int order = 1;
61  bool static_cond = 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(&freq, "-f", "--frequency", "Set the frequency for the exact"
70  " solution.");
71  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
72  "--no-static-condensation", "Enable static condensation.");
73  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
74  "--no-visualization",
75  "Enable or disable GLVis visualization.");
76 
77  args.Parse();
78  if (!args.Good())
79  {
80  if (myid == 0)
81  {
82  args.PrintUsage(cout);
83  }
84  MPI_Finalize();
85  return 1;
86  }
87  if (myid == 0)
88  {
89  args.PrintOptions(cout);
90  }
91  kappa = freq * M_PI;
92 
93  // 3. Read the (serial) mesh from the given mesh file on all processors. We
94  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
95  // and volume meshes with the same code.
96  Mesh *mesh = new Mesh(mesh_file, 1, 1);
97  dim = mesh->Dimension();
98  int sdim = mesh->SpaceDimension();
99 
100  // 4. Refine the serial mesh on all processors to increase the resolution. In
101  // this example we do 'ref_levels' of uniform refinement. We choose
102  // 'ref_levels' to be the largest number that gives a final mesh with no
103  // more than 1,000 elements.
104  {
105  int ref_levels =
106  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
107  for (int l = 0; l < ref_levels; l++)
108  {
109  mesh->UniformRefinement();
110  }
111  }
112 
113  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
114  // this mesh further in parallel to increase the resolution. Once the
115  // parallel mesh is defined, the serial mesh can be deleted. Tetrahedral
116  // meshes need to be reoriented before we can define high-order Nedelec
117  // spaces on them.
118  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
119  delete mesh;
120  {
121  int par_ref_levels = 2;
122  for (int l = 0; l < par_ref_levels; l++)
123  {
124  pmesh->UniformRefinement();
125  }
126  }
127  pmesh->ReorientTetMesh();
128 
129  // 6. Define a parallel finite element space on the parallel mesh. Here we
130  // use the Nedelec finite elements of the specified order.
131  FiniteElementCollection *fec = new ND_FECollection(order, dim);
132  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
133  HYPRE_Int size = fespace->GlobalTrueVSize();
134  if (myid == 0)
135  {
136  cout << "Number of finite element unknowns: " << size << endl;
137  }
138 
139  // 7. Determine the list of true (i.e. parallel conforming) essential
140  // boundary dofs. In this example, the boundary conditions are defined
141  // by marking all the boundary attributes from the mesh as essential
142  // (Dirichlet) and converting them to a list of true dofs.
143  Array<int> ess_tdof_list;
144  if (pmesh->bdr_attributes.Size())
145  {
146  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
147  ess_bdr = 1;
148  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
149  }
150 
151  // 8. Set up the parallel linear form b(.) which corresponds to the
152  // right-hand side of the FEM linear system, which in this case is
153  // (f,phi_i) where f is given by the function f_exact and phi_i are the
154  // basis functions in the finite element fespace.
156  ParLinearForm *b = new ParLinearForm(fespace);
158  b->Assemble();
159 
160  // 9. Define the solution vector x as a parallel finite element grid function
161  // corresponding to fespace. Initialize x by projecting the exact
162  // solution. Note that only values from the boundary edges will be used
163  // when eliminating the non-homogeneous boundary condition to modify the
164  // r.h.s. vector b.
165  ParGridFunction x(fespace);
167  x.ProjectCoefficient(E);
168 
169  // 10. Set up the parallel bilinear form corresponding to the EM diffusion
170  // operator curl muinv curl + sigma I, by adding the curl-curl and the
171  // mass domain integrators.
172  Coefficient *muinv = new ConstantCoefficient(1.0);
174  ParBilinearForm *a = new ParBilinearForm(fespace);
175  a->AddDomainIntegrator(new CurlCurlIntegrator(*muinv));
177 
178  // 11. Assemble the parallel bilinear form and the corresponding linear
179  // system, applying any necessary transformations such as: parallel
180  // assembly, eliminating boundary conditions, applying conforming
181  // constraints for non-conforming AMR, static condensation, etc.
182  if (static_cond) { a->EnableStaticCondensation(); }
183  a->Assemble();
184 
185  HypreParMatrix A;
186  Vector B, X;
187  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
188 
189  if (myid == 0)
190  {
191  cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
192  }
193 
194  // 12. Define and apply a parallel PCG solver for AX=B with the AMS
195  // preconditioner from hypre.
196  ParFiniteElementSpace *prec_fespace =
197  (a->StaticCondensationIsEnabled() ? a->SCParFESpace() : fespace);
198  HypreSolver *ams = new HypreAMS(A, prec_fespace);
199  HyprePCG *pcg = new HyprePCG(A);
200  pcg->SetTol(1e-12);
201  pcg->SetMaxIter(500);
202  pcg->SetPrintLevel(2);
203  pcg->SetPreconditioner(*ams);
204  pcg->Mult(B, X);
205 
206  // 13. Recover the parallel grid function corresponding to X. This is the
207  // local finite element solution on each processor.
208  a->RecoverFEMSolution(X, *b, x);
209 
210  // 14. Compute and print the L^2 norm of the error.
211  {
212  double err = x.ComputeL2Error(E);
213  if (myid == 0)
214  {
215  cout << "\n|| E_h - E ||_{L^2} = " << err << '\n' << endl;
216  }
217  }
218 
219  // 15. Save the refined mesh and the solution in parallel. This output can
220  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
221  {
222  ostringstream mesh_name, sol_name;
223  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
224  sol_name << "sol." << setfill('0') << setw(6) << myid;
225 
226  ofstream mesh_ofs(mesh_name.str().c_str());
227  mesh_ofs.precision(8);
228  pmesh->Print(mesh_ofs);
229 
230  ofstream sol_ofs(sol_name.str().c_str());
231  sol_ofs.precision(8);
232  x.Save(sol_ofs);
233  }
234 
235  // 16. Send the solution by socket to a GLVis server.
236  if (visualization)
237  {
238  char vishost[] = "localhost";
239  int visport = 19916;
240  socketstream sol_sock(vishost, visport);
241  sol_sock << "parallel " << num_procs << " " << myid << "\n";
242  sol_sock.precision(8);
243  sol_sock << "solution\n" << *pmesh << x << flush;
244  }
245 
246  // 17. Free the used memory.
247  delete pcg;
248  delete ams;
249  delete a;
250  delete sigma;
251  delete muinv;
252  delete b;
253  delete fespace;
254  delete fec;
255  delete pmesh;
256 
257  MPI_Finalize();
258 
259  return 0;
260 }
261 
262 
263 void E_exact(const Vector &x, Vector &E)
264 {
265  if (dim == 3)
266  {
267  E(0) = sin(kappa * x(1));
268  E(1) = sin(kappa * x(2));
269  E(2) = sin(kappa * x(0));
270  }
271  else
272  {
273  E(0) = sin(kappa * x(1));
274  E(1) = sin(kappa * x(0));
275  if (x.Size() == 3) { E(2) = 0.0; }
276  }
277 }
278 
279 void f_exact(const Vector &x, Vector &f)
280 {
281  if (dim == 3)
282  {
283  f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
284  f(1) = (1. + kappa * kappa) * sin(kappa * x(2));
285  f(2) = (1. + kappa * kappa) * sin(kappa * x(0));
286  }
287  else
288  {
289  f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
290  f(1) = (1. + kappa * kappa) * sin(kappa * x(0));
291  if (x.Size() == 3) { f(2) = 0.0; }
292  }
293 }
void SetTol(double tol)
Definition: hypre.cpp:2182
int Size() const
Logical size of the array.
Definition: array.hpp:118
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:737
ParFiniteElementSpace * SCParFESpace() const
Return the parallel trace FE space associated with static condensation.
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:936
Subclass constant coefficient.
Definition: coefficient.hpp:67
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 ReorientTetMesh()
See the remarks for the serial version in mesh.hpp.
Definition: pmesh.cpp:2402
Integrator for (curl u, curl v) for Nedelec elements.
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
int Size() const
Returns the size of the vector.
Definition: vector.hpp:150
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:676
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:488
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:292
int main(int argc, char *argv[])
Definition: ex1.cpp:58
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2197
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:417
int dim
Definition: ex3.cpp:48
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7591
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:251
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:3863
int Dimension() const
Definition: mesh.hpp:713
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2187
int SpaceDimension() const
Definition: mesh.hpp:714
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
bool StaticCondensationIsEnabled() const
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:179
PCG solver in hypre.
Definition: hypre.hpp:706
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:76
void f_exact(const Vector &, Vector &)
Definition: ex3.cpp:228
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:248
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2202
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
void E_exact(const Vector &, Vector &)
Definition: ex3.cpp:212
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:670
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:231
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:240
double kappa
Definition: ex3.cpp:47
Vector data type.
Definition: vector.hpp:48
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2223
double freq
Definition: ex3.cpp:47
Class for parallel meshes.
Definition: pmesh.hpp:32
Integrator for (Q u, v) for VectorFiniteElements.
void EnableStaticCondensation()
double sigma(const Vector &x)
Definition: maxwell.cpp:102
bool Good() const
Definition: optparser.hpp:122