MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex24p.cpp
Go to the documentation of this file.
1 // MFEM Example 24 - Parallel Version
2 //
3 // Compile with: make ex24p
4 //
5 // Sample runs: mpirun -np 4 ex24p -m ../data/star.mesh
6 // mpirun -np 4 ex24p -m ../data/square-disc.mesh -o 2
7 // mpirun -np 4 ex24p -m ../data/beam-tet.mesh
8 // mpirun -np 4 ex24p -m ../data/beam-hex.mesh -o 2 -pa
9 // mpirun -np 4 ex24p -m ../data/escher.mesh
10 // mpirun -np 4 ex24p -m ../data/escher.mesh -o 2
11 // mpirun -np 4 ex24p -m ../data/fichera.mesh
12 // mpirun -np 4 ex24p -m ../data/fichera-q2.vtk
13 // mpirun -np 4 ex24p -m ../data/fichera-q3.mesh
14 // mpirun -np 4 ex24p -m ../data/square-disc-nurbs.mesh
15 // mpirun -np 4 ex24p -m ../data/beam-hex-nurbs.mesh
16 // mpirun -np 4 ex24p -m ../data/amr-quad.mesh -o 2
17 // mpirun -np 4 ex24p -m ../data/amr-hex.mesh
18 //
19 // Device sample runs:
20 // mpirun -np 4 ex24p -m ../data/star.mesh -pa -d cuda
21 // mpirun -np 4 ex24p -m ../data/star.mesh -pa -d raja-cuda
22 // mpirun -np 4 ex24p -m ../data/star.mesh -pa -d raja-omp
23 // mpirun -np 4 ex24p -m ../data/beam-hex.mesh -pa -d cuda
24 //
25 // Description: This example code illustrates usage of mixed finite element
26 // spaces. Using two different approaches, we project a gradient
27 // of a function in H^1 to H(curl). Other spaces and example
28 // computations are to be added in the future.
29 //
30 // We recommend viewing examples 1 and 3 before viewing this
31 // example.
32 
33 #include "mfem.hpp"
34 #include <fstream>
35 #include <iostream>
36 
37 using namespace std;
38 using namespace mfem;
39 
40 double p_exact(const Vector &x);
41 void gradp_exact(const Vector &, Vector &);
42 
43 int dim;
44 
45 int main(int argc, char *argv[])
46 {
47  // 1. Initialize MPI.
48  int num_procs, myid;
49  MPI_Init(&argc, &argv);
50  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
51  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
52 
53  // 2. Parse command-line options.
54  const char *mesh_file = "../data/beam-hex.mesh";
55  int order = 1;
56  bool static_cond = false;
57  bool pa = false;
58  const char *device_config = "cpu";
59  bool visualization = 1;
60 
61  OptionsParser args(argc, argv);
62  args.AddOption(&mesh_file, "-m", "--mesh",
63  "Mesh file to use.");
64  args.AddOption(&order, "-o", "--order",
65  "Finite element order (polynomial degree).");
66  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
67  "--no-static-condensation", "Enable static condensation.");
68  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
69  "--no-partial-assembly", "Enable Partial Assembly.");
70  args.AddOption(&device_config, "-d", "--device",
71  "Device configuration string, see Device::Configure().");
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 
91  // 3. Enable hardware devices such as GPUs, and programming models such as
92  // CUDA, OCCA, RAJA and OpenMP based on command line options.
93  Device device(device_config);
94  if (myid == 0) { device.Print(); }
95 
96  // 4. Read the (serial) mesh from the given mesh file on all processors. We
97  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
98  // and volume meshes with the same code.
99  Mesh *mesh = new Mesh(mesh_file, 1, 1);
100  dim = mesh->Dimension();
101  int sdim = mesh->SpaceDimension();
102 
103  // 5. Refine the serial mesh on all processors to increase the resolution. In
104  // this example we do 'ref_levels' of uniform refinement. We choose
105  // 'ref_levels' to be the largest number that gives a final mesh with no
106  // more than 1,000 elements.
107  {
108  int ref_levels = (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
109  for (int l = 0; l < ref_levels; l++)
110  {
111  mesh->UniformRefinement();
112  }
113  }
114 
115  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
116  // this mesh further in parallel to increase the resolution. Once the
117  // parallel mesh is defined, the serial mesh can be deleted. Tetrahedral
118  // meshes need to be reoriented before we can define high-order Nedelec
119  // spaces on them.
120  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
121  delete mesh;
122  {
123  int par_ref_levels = 1;
124  for (int l = 0; l < par_ref_levels; l++)
125  {
126  pmesh->UniformRefinement();
127  }
128  }
129  pmesh->ReorientTetMesh();
130 
131  // 7. Define a parallel finite element space on the parallel mesh. Here we
132  // use the Nedelec finite elements of the specified order.
133  FiniteElementCollection *fec = new ND_FECollection(order, dim);
134  FiniteElementCollection *H1fec = new H1_FECollection(order, dim);
135  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
136  ParFiniteElementSpace *H1fespace = new ParFiniteElementSpace(pmesh, H1fec);
137  HYPRE_Int size = fespace->GlobalTrueVSize();
138  HYPRE_Int H1size = H1fespace->GlobalTrueVSize();
139  if (myid == 0)
140  {
141  cout << "Number of Nedelec finite element unknowns: " << size << endl;
142  cout << "Number of H1 finite element unknowns: " << H1size << endl;
143  }
144 
145  // 8. Define the solution vector x as a parallel finite element grid function
146  // corresponding to fespace. Initialize x by projecting the exact
147  // solution. Note that only values from the boundary edges will be used
148  // when eliminating the non-homogeneous boundary condition to modify the
149  // r.h.s. vector b.
150  ParGridFunction x(fespace);
152  ParGridFunction p(H1fespace);
153  p.ProjectCoefficient(p_coef);
154  p.SetTrueVector();
155  p.SetFromTrueVector();
156 
157  VectorFunctionCoefficient gradp_coef(sdim, gradp_exact);
158 
159  // 9. Set up the parallel bilinear forms.
160  Coefficient *muinv = new ConstantCoefficient(1.0);
162  ParBilinearForm *a = new ParBilinearForm(fespace);
163  ParMixedBilinearForm *a_NDH1 = new ParMixedBilinearForm(H1fespace, fespace);
164  if (pa)
165  {
166  a->SetAssemblyLevel(AssemblyLevel::PARTIAL);
167  a_NDH1->SetAssemblyLevel(AssemblyLevel::PARTIAL);
168  }
169 
170  // First approach: L2 projection
173 
174  // 10. Assemble the parallel bilinear form and the corresponding linear
175  // system, applying any necessary transformations such as: parallel
176  // assembly, eliminating boundary conditions, applying conforming
177  // constraints for non-conforming AMR, static condensation, etc.
178  if (static_cond) { a->EnableStaticCondensation(); }
179 
180  a->Assemble();
181  if (!pa) { a->Finalize(); }
182 
183  a_NDH1->Assemble();
184  if (!pa) { a_NDH1->Finalize(); }
185 
186  Vector B(fespace->GetTrueVSize());
187  Vector X(fespace->GetTrueVSize());
188 
189  if (pa)
190  {
191  ParLinearForm *b = new ParLinearForm(fespace); // used as a vector
192  a_NDH1->Mult(p, *b); // process-local multiplication
193  b->ParallelAssemble(B);
194  delete b;
195  }
196  else
197  {
198  HypreParMatrix *NDH1 = a_NDH1->ParallelAssemble();
199 
200  Vector P(H1fespace->GetTrueVSize());
201  p.GetTrueDofs(P);
202 
203  NDH1->Mult(P,B);
204 
205  delete NDH1;
206  }
207 
208  // 11. Define and apply a parallel PCG solver for AX=B with Jacobi
209  // preconditioner.
210  if (pa)
211  {
212  Array<int> ess_tdof_list; // empty
213 
214  OperatorPtr A;
215  a->FormSystemMatrix(ess_tdof_list, A);
216 
217  OperatorJacobiSmoother Jacobi(*a, ess_tdof_list);
218 
219  CGSolver cg(MPI_COMM_WORLD);
220  cg.SetRelTol(1e-12);
221  cg.SetMaxIter(1000);
222  cg.SetPrintLevel(1);
223  cg.SetOperator(*A);
224  cg.SetPreconditioner(Jacobi);
225  X = 0.0;
226  cg.Mult(B, X);
227  }
228  else
229  {
230  HypreParMatrix *Amat = a->ParallelAssemble();
231  HypreDiagScale Jacobi(*Amat);
232  HyprePCG pcg(*Amat);
233  pcg.SetTol(1e-12);
234  pcg.SetMaxIter(1000);
235  pcg.SetPrintLevel(2);
236  pcg.SetPreconditioner(Jacobi);
237  X = 0.0;
238  pcg.Mult(B, X);
239 
240  delete Amat;
241  }
242 
243  x.SetFromTrueDofs(X);
244 
245  // 12. Second approach: compute the same solution by applying
246  // GradientInterpolator in H(curl).
247  ParDiscreteLinearOperator grad(H1fespace, fespace);
249  grad.Assemble();
250 
251  ParGridFunction gradp(fespace);
252  grad.Mult(p, gradp);
253 
254  // 13. Compute the projection of the exact grad p.
255  ParGridFunction exact_gradp(fespace);
256  exact_gradp.ProjectCoefficient(gradp_coef);
257  exact_gradp.SetTrueVector();
258  exact_gradp.SetFromTrueVector();
259 
260  // 14. Compute and print the L^2 norm of the error.
261  {
262  double errSol = x.ComputeL2Error(gradp_coef);
263  double errInterp = gradp.ComputeL2Error(gradp_coef);
264  double errProj = exact_gradp.ComputeL2Error(gradp_coef);
265 
266  if (myid == 0)
267  {
268  cout << "\n Solution of (E_h,v) = (grad p_h,v) for E_h and v in "
269  "H(curl): || E_h - grad p ||_{L^2} = " << errSol << '\n' << endl;
270  cout << " Gradient interpolant E_h = grad p_h in H(curl): || E_h - "
271  "grad p ||_{L^2} = " << errInterp << '\n' << endl;
272  cout << " Projection E_h of exact grad p in H(curl): || E_h - grad p "
273  "||_{L^2} = " << errProj << '\n' << endl;
274  }
275  }
276 
277  // 15. Save the refined mesh and the solution in parallel. This output can
278  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
279  {
280  ostringstream mesh_name, sol_name;
281  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
282  sol_name << "sol." << setfill('0') << setw(6) << myid;
283 
284  ofstream mesh_ofs(mesh_name.str().c_str());
285  mesh_ofs.precision(8);
286  pmesh->Print(mesh_ofs);
287 
288  ofstream sol_ofs(sol_name.str().c_str());
289  sol_ofs.precision(8);
290  x.Save(sol_ofs);
291  }
292 
293  // 16. Send the solution by socket to a GLVis server.
294  if (visualization)
295  {
296  char vishost[] = "localhost";
297  int visport = 19916;
298  socketstream sol_sock(vishost, visport);
299  sol_sock << "parallel " << num_procs << " " << myid << "\n";
300  sol_sock.precision(8);
301  sol_sock << "solution\n" << *pmesh << x << flush;
302  }
303 
304  // 17. Free the used memory.
305  delete a;
306  delete a_NDH1;
307  delete sigma;
308  delete muinv;
309  delete fespace;
310  delete H1fespace;
311  delete fec;
312  delete H1fec;
313  delete pmesh;
314 
315  MPI_Finalize();
316 
317  return 0;
318 }
319 
320 double p_exact(const Vector &x)
321 {
322  if (dim == 3)
323  {
324  return sin(x(0)) * sin(x(1)) * sin(x(2));
325  }
326  else if (dim == 2)
327  {
328  return sin(x(0)) * sin(x(1));
329  }
330 
331  return 0.0;
332 }
333 
334 void gradp_exact(const Vector &x, Vector &f)
335 {
336  if (dim == 3)
337  {
338  f(0) = cos(x(0)) * sin(x(1)) * sin(x(2));
339  f(1) = sin(x(0)) * cos(x(1)) * sin(x(2));
340  f(2) = sin(x(0)) * sin(x(1)) * cos(x(2));
341  }
342  else
343  {
344  f(0) = cos(x(0)) * sin(x(1));
345  f(1) = sin(x(0)) * cos(x(1));
346  if (x.Size() == 3) { f(2) = 0.0; }
347  }
348 }
void SetTol(double tol)
Definition: hypre.cpp:2305
Conjugate gradient method.
Definition: solvers.hpp:152
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:139
Subclass constant coefficient.
Definition: coefficient.hpp:67
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:358
virtual void ReorientTetMesh()
See the remarks for the serial version in mesh.hpp.
Definition: pmesh.cpp:2508
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:693
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:485
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:236
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::FULL.
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:301
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
int main(int argc, char *argv[])
Definition: ex1.cpp:62
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A * x + beta * y.
Definition: hypre.cpp:1021
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2320
void AddDomainInterpolator(DiscreteInterpolator *di)
Adds a domain interpolator. Assumes ownership of di.
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:255
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P_test^t A P_trial.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds a domain integrator. Assumes ownership of bfi.
Class for parallel linear form.
Definition: plinearform.hpp:26
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:133
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: pgridfunc.hpp:159
Jacobi preconditioner in hypre.
Definition: hypre.hpp:859
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:84
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7982
void SetMaxIter(int max_it)
Definition: solvers.hpp:65
void Assemble(int skip_zeros=1)
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:4102
int Dimension() const
Definition: mesh.hpp:744
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2310
double p_exact(const Vector &x)
Definition: ex24.cpp:258
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::FULL.
int SpaceDimension() const
Definition: mesh.hpp:745
void SetRelTol(double rtol)
Definition: solvers.hpp:63
PCG solver in hypre.
Definition: hypre.hpp:743
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
double a
Definition: lissajous.cpp:41
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
Class for parallel bilinear form using different test and trial FE spaces.
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:250
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
int dim
Definition: ex24.cpp:43
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2325
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:143
Class for parallel bilinear form.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:166
class for C-function coefficient
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:251
Vector data type.
Definition: vector.hpp:48
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:83
Class for parallel grid function.
Definition: pgridfunc.hpp:32
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:114
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
virtual void Assemble(int skip_zeros=1)
Construct the internal matrix representation of the discrete linear operator.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2348
Class for parallel meshes.
Definition: pmesh.hpp:32
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:34
void EnableStaticCondensation()
void gradp_exact(const Vector &, Vector &)
Definition: ex24.cpp:272
double sigma(const Vector &x)
Definition: maxwell.cpp:102
bool Good() const
Definition: optparser.hpp:125