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