MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex1p.cpp
Go to the documentation of this file.
1 // MFEM Example 1 - Parallel Version
2 //
3 // Compile with: make ex1p
4 //
5 // Sample runs: mpirun -np 4 ex1p -m ../data/square-disc.mesh
6 // mpirun -np 4 ex1p -m ../data/star.mesh
7 // mpirun -np 4 ex1p -m ../data/star-mixed.mesh
8 // mpirun -np 4 ex1p -m ../data/escher.mesh
9 // mpirun -np 4 ex1p -m ../data/fichera.mesh
10 // mpirun -np 4 ex1p -m ../data/fichera-mixed.mesh
11 // mpirun -np 4 ex1p -m ../data/toroid-wedge.mesh
12 // mpirun -np 4 ex1p -m ../data/square-disc-p2.vtk -o 2
13 // mpirun -np 4 ex1p -m ../data/square-disc-p3.mesh -o 3
14 // mpirun -np 4 ex1p -m ../data/square-disc-nurbs.mesh -o -1
15 // mpirun -np 4 ex1p -m ../data/star-mixed-p2.mesh -o 2
16 // mpirun -np 4 ex1p -m ../data/disc-nurbs.mesh -o -1
17 // mpirun -np 4 ex1p -m ../data/pipe-nurbs.mesh -o -1
18 // mpirun -np 4 ex1p -m ../data/ball-nurbs.mesh -o 2
19 // mpirun -np 4 ex1p -m ../data/fichera-mixed-p2.mesh -o 2
20 // mpirun -np 4 ex1p -m ../data/star-surf.mesh
21 // mpirun -np 4 ex1p -m ../data/square-disc-surf.mesh
22 // mpirun -np 4 ex1p -m ../data/inline-segment.mesh
23 // mpirun -np 4 ex1p -m ../data/amr-quad.mesh
24 // mpirun -np 4 ex1p -m ../data/amr-hex.mesh
25 // mpirun -np 4 ex1p -m ../data/mobius-strip.mesh
26 // mpirun -np 4 ex1p -m ../data/mobius-strip.mesh -o -1 -sc
27 //
28 // Device sample runs:
29 // mpirun -np 4 ex1p -pa -d cuda
30 // mpirun -np 4 ex1p -pa -d occa-cuda
31 // mpirun -np 4 ex1p -pa -d raja-omp
32 //
33 // Description: This example code demonstrates the use of MFEM to define a
34 // simple finite element discretization of the Laplace problem
35 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
36 // Specifically, we discretize using a FE space of the specified
37 // order, or if order < 1 using an isoparametric/isogeometric
38 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
39 // NURBS mesh, etc.)
40 //
41 // The example highlights the use of mesh refinement, finite
42 // element grid functions, as well as linear and bilinear forms
43 // corresponding to the left-hand side and right-hand side of the
44 // discrete linear system. We also cover the explicit elimination
45 // of essential boundary conditions, static condensation, and the
46 // optional connection to the GLVis tool for visualization.
47 
48 #include "mfem.hpp"
49 #include <fstream>
50 #include <iostream>
51 
52 using namespace std;
53 using namespace mfem;
54 
55 int main(int argc, char *argv[])
56 {
57  // 1. Initialize MPI.
58  int num_procs, myid;
59  MPI_Init(&argc, &argv);
60  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
61  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
62 
63  // 2. Parse command-line options.
64  const char *mesh_file = "../data/star.mesh";
65  int order = 1;
66  bool static_cond = false;
67  bool pa = false;
68  const char *device_config = "cpu";
69  bool visualization = true;
70 
71  OptionsParser args(argc, argv);
72  args.AddOption(&mesh_file, "-m", "--mesh",
73  "Mesh file to use.");
74  args.AddOption(&order, "-o", "--order",
75  "Finite element order (polynomial degree) or -1 for"
76  " isoparametric space.");
77  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
78  "--no-static-condensation", "Enable static condensation.");
79  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
80  "--no-partial-assembly", "Enable Partial Assembly.");
81  args.AddOption(&device_config, "-d", "--device",
82  "Device configuration string, see Device::Configure().");
83  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
84  "--no-visualization",
85  "Enable or disable GLVis visualization.");
86  args.Parse();
87  if (!args.Good())
88  {
89  if (myid == 0)
90  {
91  args.PrintUsage(cout);
92  }
93  MPI_Finalize();
94  return 1;
95  }
96  if (myid == 0)
97  {
98  args.PrintOptions(cout);
99  }
100 
101  // 3. Enable hardware devices such as GPUs, and programming models such as
102  // CUDA, OCCA, RAJA and OpenMP based on command line options.
103  Device device(device_config);
104  if (myid == 0) { device.Print(); }
105 
106  // 4. Read the (serial) mesh from the given mesh file on all processors. We
107  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
108  // and volume meshes with the same code.
109  Mesh *mesh = new Mesh(mesh_file, 1, 1);
110  int dim = mesh->Dimension();
111 
112  // 5. Refine the serial mesh on all processors to increase the resolution. In
113  // this example we do 'ref_levels' of uniform refinement. We choose
114  // 'ref_levels' to be the largest number that gives a final mesh with no
115  // more than 10,000 elements.
116  {
117  int ref_levels =
118  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
119  for (int l = 0; l < ref_levels; l++)
120  {
121  mesh->UniformRefinement();
122  }
123  }
124 
125  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
126  // this mesh further in parallel to increase the resolution. Once the
127  // parallel mesh is defined, the serial mesh can be deleted.
128  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
129  delete mesh;
130  {
131  int par_ref_levels = 2;
132  for (int l = 0; l < par_ref_levels; l++)
133  {
134  pmesh->UniformRefinement();
135  }
136  }
137 
138  // 7. Define a parallel finite element space on the parallel mesh. Here we
139  // use continuous Lagrange finite elements of the specified order. If
140  // order < 1, we instead use an isoparametric/isogeometric space.
142  if (order > 0)
143  {
144  fec = new H1_FECollection(order, dim);
145  }
146  else if (pmesh->GetNodes())
147  {
148  fec = pmesh->GetNodes()->OwnFEC();
149  if (myid == 0)
150  {
151  cout << "Using isoparametric FEs: " << fec->Name() << endl;
152  }
153  }
154  else
155  {
156  fec = new H1_FECollection(order = 1, dim);
157  }
158  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
159  HYPRE_Int size = fespace->GlobalTrueVSize();
160  if (myid == 0)
161  {
162  cout << "Number of finite element unknowns: " << size << endl;
163  }
164 
165  // 8. Determine the list of true (i.e. parallel conforming) essential
166  // boundary dofs. In this example, the boundary conditions are defined
167  // by marking all the boundary attributes from the mesh as essential
168  // (Dirichlet) and converting them to a list of true dofs.
169  Array<int> ess_tdof_list;
170  if (pmesh->bdr_attributes.Size())
171  {
172  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
173  ess_bdr = 1;
174  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
175  }
176 
177  // 9. Set up the parallel linear form b(.) which corresponds to the
178  // right-hand side of the FEM linear system, which in this case is
179  // (1,phi_i) where phi_i are the basis functions in fespace.
180  ParLinearForm *b = new ParLinearForm(fespace);
181  ConstantCoefficient one(1.0);
183  b->Assemble();
184 
185  // 10. Define the solution vector x as a parallel finite element grid function
186  // corresponding to fespace. Initialize x with initial guess of zero,
187  // which satisfies the boundary conditions.
188  ParGridFunction x(fespace);
189  x = 0.0;
190 
191  // 11. Set up the parallel bilinear form a(.,.) on the finite element space
192  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
193  // domain integrator.
194  ParBilinearForm *a = new ParBilinearForm(fespace);
195  if (pa) { a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
197 
198  // 12. Assemble the parallel bilinear form and the corresponding linear
199  // system, applying any necessary transformations such as: parallel
200  // assembly, eliminating boundary conditions, applying conforming
201  // constraints for non-conforming AMR, static condensation, etc.
202  if (static_cond) { a->EnableStaticCondensation(); }
203  a->Assemble();
204 
205  OperatorPtr A;
206  Vector B, X;
207  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
208 
209  // 13. Solve the linear system A X = B.
210  // * With full assembly, use the BoomerAMG preconditioner from hypre.
211  // * With partial assembly, use no preconditioner, for now.
212  Solver *prec = NULL;
213  if (!pa) { prec = new HypreBoomerAMG; }
214  CGSolver cg(MPI_COMM_WORLD);
215  cg.SetRelTol(1e-12);
216  cg.SetMaxIter(2000);
217  cg.SetPrintLevel(1);
218  if (prec) { cg.SetPreconditioner(*prec); }
219  cg.SetOperator(*A);
220  cg.Mult(B, X);
221  delete prec;
222 
223  // 14. Recover the parallel grid function corresponding to X. This is the
224  // local finite element solution on each processor.
225  a->RecoverFEMSolution(X, *b, x);
226 
227  // 15. Save the refined mesh and the solution in parallel. This output can
228  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
229  {
230  ostringstream mesh_name, sol_name;
231  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
232  sol_name << "sol." << setfill('0') << setw(6) << myid;
233 
234  ofstream mesh_ofs(mesh_name.str().c_str());
235  mesh_ofs.precision(8);
236  pmesh->Print(mesh_ofs);
237 
238  ofstream sol_ofs(sol_name.str().c_str());
239  sol_ofs.precision(8);
240  x.Save(sol_ofs);
241  }
242 
243  // 16. Send the solution by socket to a GLVis server.
244  if (visualization)
245  {
246  char vishost[] = "localhost";
247  int visport = 19916;
248  socketstream sol_sock(vishost, visport);
249  sol_sock << "parallel " << num_procs << " " << myid << "\n";
250  sol_sock.precision(8);
251  sol_sock << "solution\n" << *pmesh << x << flush;
252  }
253 
254  // 17. Free the used memory.
255  delete a;
256  delete b;
257  delete fespace;
258  if (order > 0) { delete fec; }
259  delete pmesh;
260 
261  MPI_Finalize();
262 
263  return 0;
264 }
int Size() const
Logical size of the array.
Definition: array.hpp:118
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:737
Conjugate gradient method.
Definition: solvers.hpp:111
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 Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:290
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:676
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:488
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:98
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::FULL.
Abstract parallel finite element space.
Definition: pfespace.hpp:28
int main(int argc, char *argv[])
Definition: ex1.cpp:58
The BoomerAMG solver in hypre.
Definition: hypre.hpp:873
Class for parallel linear form.
Definition: plinearform.hpp:26
int dim
Definition: ex3.cpp:48
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:67
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7591
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: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 AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:179
void SetRelTol(double rtol)
Definition: solvers.hpp:61
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
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
virtual const char * Name() const
Definition: fe_coll.hpp:53
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
Class for parallel bilinear form.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:125
Vector data type.
Definition: vector.hpp:48
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5837
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:88
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:83
Base class for solvers.
Definition: operator.hpp:279
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:96
Class for parallel meshes.
Definition: pmesh.hpp:32
void EnableStaticCondensation()
bool Good() const
Definition: optparser.hpp:122