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