MFEM  v4.5.1
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 // SuperLU Modification
3 //
4 // Compile with: make ex1p
5 //
6 // Sample runs: mpirun -np 4 ex1p -m ../../data/square-disc.mesh
7 // mpirun -np 4 ex1p -m ../../data/star.mesh
8 // mpirun -np 4 ex1p -m ../../data/star-mixed.mesh
9 // mpirun -np 4 ex1p -m ../../data/escher.mesh
10 // mpirun -np 4 ex1p -m ../../data/fichera.mesh
11 // mpirun -np 4 ex1p -m ../../data/fichera-mixed.mesh
12 // mpirun -np 4 ex1p -m ../../data/toroid-wedge.mesh
13 // mpirun -np 4 ex1p -m ../../data/periodic-annulus-sector.msh
14 // mpirun -np 4 ex1p -m ../../data/periodic-torus-sector.msh
15 // mpirun -np 4 ex1p -m ../../data/square-disc-p2.vtk -o 2
16 // mpirun -np 4 ex1p -m ../../data/square-disc-nurbs.mesh -o -1
17 // mpirun -np 4 ex1p -m ../../data/star-mixed-p2.mesh -o 2
18 // mpirun -np 4 ex1p -m ../../data/disc-nurbs.mesh -o -1
19 // mpirun -np 4 ex1p -m ../../data/pipe-nurbs.mesh -o -1
20 // mpirun -np 4 ex1p -m ../../data/ball-nurbs.mesh -o 2
21 // mpirun -np 4 ex1p -m ../../data/star-surf.mesh
22 // mpirun -np 4 ex1p -m ../../data/square-disc-surf.mesh
23 // mpirun -np 4 ex1p -m ../../data/inline-segment.mesh
24 // mpirun -np 4 ex1p -m ../../data/amr-quad.mesh
25 // mpirun -np 4 ex1p -m ../../data/amr-hex.mesh
26 // mpirun -np 4 ex1p -m ../../data/mobius-strip.mesh
27 //
28 // Description: This example code demonstrates the use of MFEM to define a
29 // simple finite element discretization of the Laplace problem
30 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
31 // Specifically, we discretize using a FE space of the specified
32 // order, or if order < 1 using an isoparametric/isogeometric
33 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
34 // NURBS mesh, etc.)
35 //
36 // The example highlights the use of mesh refinement, finite
37 // element grid functions, as well as linear and bilinear forms
38 // corresponding to the left-hand side and right-hand side of the
39 // discrete linear system. We also cover the explicit elimination
40 // of essential boundary conditions, static condensation, and the
41 // optional connection to the GLVis tool for visualization.
42 
43 #include "mfem.hpp"
44 #include <fstream>
45 #include <iostream>
46 
47 #ifndef MFEM_USE_SUPERLU
48 #error This example requires that MFEM is built with MFEM_USE_SUPERLU=YES
49 #endif
50 
51 using namespace std;
52 using namespace mfem;
53 
54 int main(int argc, char *argv[])
55 {
56  // 1. Initialize MPI and HYPRE.
57  Mpi::Init(argc, argv);
58  int num_procs = Mpi::WorldSize();
59  int myid = Mpi::WorldRank();
60  Hypre::Init();
61 
62  // 2. Parse command-line options.
63  const char *mesh_file = "../../data/star.mesh";
64  int order = 1;
65  const char *device_config = "cpu";
66  bool visualization = true;
67  int slu_colperm = 4;
68  int slu_rowperm = 1;
69  int slu_iterref = 2;
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(&device_config, "-d", "--device",
78  "Device configuration string, see Device::Configure().");
79  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
80  "--no-visualization",
81  "Enable or disable GLVis visualization.");
82  args.AddOption(&slu_colperm, "-cp", "--colperm",
83  "SuperLU Column Permutation Method: 0-NATURAL, 1-MMD-ATA "
84  "2-MMD_AT_PLUS_A, 3-COLAMD, 4-METIS_AT_PLUS_A, 5-PARMETIS "
85  "6-ZOLTAN");
86  args.AddOption(&slu_rowperm, "-rp", "--rowperm",
87  "SuperLU Row Permutation Method: 0-NOROWPERM, 1-LargeDiag");
88  args.AddOption(&slu_iterref, "-rp", "--rowperm",
89  "SuperLU Iterative Refinement: 0-NOREFINE, 1-Single, "
90  "2-Double, 3-Extra");
91 
92  args.Parse();
93  if (!args.Good())
94  {
95  if (myid == 0)
96  {
97  args.PrintUsage(cout);
98  }
99  return 1;
100  }
101  if (myid == 0)
102  {
103  args.PrintOptions(cout);
104  }
105 
106  // 3. Enable hardware devices such as GPUs, and programming models such as
107  // CUDA, OCCA, RAJA and OpenMP based on command line options.
108  Device device(device_config);
109  if (myid == 0) { device.Print(); }
110 
111  // 4. Read the (serial) mesh from the given mesh file on all processors. We
112  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
113  // and volume meshes with the same code.
114  Mesh mesh(mesh_file, 1, 1);
115  int dim = mesh.Dimension();
116 
117  // 5. Refine the serial mesh on all processors to increase the resolution. In
118  // this example we do 'ref_levels' of uniform refinement. We choose
119  // 'ref_levels' to be the largest number that gives a final mesh with no
120  // more than 1,000 elements.
121  {
122  int ref_levels =
123  (int)floor(log(1000./mesh.GetNE())/log(2.)/dim);
124  for (int l = 0; l < ref_levels; l++)
125  {
126  mesh.UniformRefinement();
127  }
128  }
129 
130  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
131  // this mesh further in parallel to increase the resolution. Once the
132  // parallel mesh is defined, the serial mesh can be deleted.
133  ParMesh pmesh(MPI_COMM_WORLD, mesh);
134  mesh.Clear();
135  {
136  int par_ref_levels = 2;
137  for (int l = 0; l < par_ref_levels; l++)
138  {
139  pmesh.UniformRefinement();
140  }
141  }
142 
143  // 7. Define a parallel finite element space on the parallel mesh. Here we
144  // use continuous Lagrange finite elements of the specified order. If
145  // order < 1, we instead use an isoparametric/isogeometric space.
147  bool delete_fec;
148  if (order > 0)
149  {
150  fec = new H1_FECollection(order, dim);
151  delete_fec = true;
152  }
153  else if (pmesh.GetNodes())
154  {
155  fec = pmesh.GetNodes()->OwnFEC();
156  delete_fec = false;
157  if (myid == 0)
158  {
159  cout << "Using isoparametric FEs: " << fec->Name() << endl;
160  }
161  }
162  else
163  {
164  fec = new H1_FECollection(order = 1, dim);
165  delete_fec = true;
166  }
167  ParFiniteElementSpace fespace(&pmesh, fec);
168  HYPRE_BigInt size = fespace.GlobalTrueVSize();
169  if (myid == 0)
170  {
171  cout << "Number of finite element unknowns: " << size << endl;
172  }
173 
174  // 8. Determine the list of true (i.e. parallel conforming) essential
175  // boundary dofs. In this example, the boundary conditions are defined
176  // by marking all the boundary attributes from the mesh as essential
177  // (Dirichlet) and converting them to a list of true dofs.
178  Array<int> ess_tdof_list;
179  if (pmesh.bdr_attributes.Size())
180  {
181  Array<int> ess_bdr(pmesh.bdr_attributes.Max());
182  ess_bdr = 1;
183  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
184  }
185 
186  // 9. Set up the parallel linear form b(.) which corresponds to the
187  // right-hand side of the FEM linear system, which in this case is
188  // (1,phi_i) where phi_i are the basis functions in fespace.
189  ParLinearForm b(&fespace);
190  ConstantCoefficient one(1.0);
192  b.Assemble();
193 
194  // 10. Define the solution vector x as a parallel finite element grid function
195  // corresponding to fespace. Initialize x with initial guess of zero,
196  // which satisfies the boundary conditions.
197  ParGridFunction x(&fespace);
198  x = 0.0;
199 
200  // 11. Set up the parallel bilinear form a(.,.) on the finite element space
201  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
202  // domain integrator.
203  ParBilinearForm a(&fespace);
205 
206  // 12. Assemble the parallel bilinear form and the corresponding linear
207  // system, applying any necessary transformations such as: parallel
208  // assembly, eliminating boundary conditions, applying conforming
209  // constraints for non-conforming AMR, static condensation, etc.
210  a.Assemble();
211 
212  OperatorPtr A;
213  Vector B, X;
214  a.FormLinearSystem(ess_tdof_list, x, b, A, X, B);
215 
216  // 13. Solve the linear system A X = B utilizing SuperLU.
217  SuperLUSolver *superlu = new SuperLUSolver(MPI_COMM_WORLD);
218  Operator *SLU_A = new SuperLURowLocMatrix(*A.As<HypreParMatrix>());
219  superlu->SetPrintStatistics(true);
220  superlu->SetSymmetricPattern(false);
221 
222  if (slu_colperm == 0)
223  {
225  }
226  else if (slu_colperm == 1)
227  {
229  }
230  else if (slu_colperm == 2)
231  {
233  }
234  else if (slu_colperm == 3)
235  {
237  }
238  else if (slu_colperm == 4)
239  {
241  }
242  else if (slu_colperm == 5)
243  {
245  }
246  else if (slu_colperm == 6)
247  {
249  }
250 
251  if (slu_rowperm == 0)
252  {
254  }
255  else if (slu_rowperm == 1)
256  {
257 #ifdef MFEM_USE_SUPERLU5
259 #else
261 #endif
262  }
263 
264  if (slu_iterref == 0)
265  {
267  }
268  else if (slu_iterref == 1)
269  {
271  }
272  else if (slu_iterref == 2)
273  {
275  }
276  else if (slu_iterref == 3)
277  {
279  }
280 
281  superlu->SetOperator(*SLU_A);
282  superlu->SetPrintStatistics(true);
283  superlu->Mult(B, X);
284  superlu->DismantleGrid();
285 
286  delete SLU_A;
287  delete superlu;
288 
289  // 14. Recover the parallel grid function corresponding to X. This is the
290  // local finite element solution on each processor.
291  a.RecoverFEMSolution(X, b, x);
292 
293  // 15. Save the refined mesh and the solution in parallel. This output can
294  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
295  {
296  ostringstream mesh_name, sol_name;
297  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
298  sol_name << "sol." << setfill('0') << setw(6) << myid;
299 
300  ofstream mesh_ofs(mesh_name.str().c_str());
301  mesh_ofs.precision(8);
302  pmesh.Print(mesh_ofs);
303 
304  ofstream sol_ofs(sol_name.str().c_str());
305  sol_ofs.precision(8);
306  x.Save(sol_ofs);
307  }
308 
309  // 16. Send the solution by socket to a GLVis server.
310  if (visualization)
311  {
312  char vishost[] = "localhost";
313  int visport = 19916;
314  socketstream sol_sock(vishost, visport);
315  sol_sock << "parallel " << num_procs << " " << myid << "\n";
316  sol_sock.precision(8);
317  sol_sock << "solution\n" << pmesh << x << flush;
318  }
319 
320  // 17. Free the used memory.
321  if (delete_fec)
322  {
323  delete fec;
324  }
325 
326  return 0;
327 }
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:108
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1032
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:104
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: superlu.cpp:485
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
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(...
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:285
void SetRowPermutation(superlu::RowPerm row_perm, Array< int > *perm=NULL)
Definition: superlu.cpp:354
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:923
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:873
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:279
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetSymmetricPattern(bool sym)
Definition: superlu.cpp:425
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: superlu.cpp:591
Class for parallel linear form.
Definition: plinearform.hpp:26
void SetIterativeRefine(superlu::IterRefine iter_ref)
Definition: superlu.cpp:391
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
constexpr char vishost[]
void SetColumnPermutation(superlu::ColPerm col_perm)
Definition: superlu.cpp:345
void Assemble()
Assembles the ParLinearForm i.e. sums over all domain/bdr integrators.
Definition: plinearform.cpp:46
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9816
constexpr int visport
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void SetPrintStatistics(bool print_stat)
Definition: superlu.cpp:327
int Dimension() const
Definition: mesh.hpp:1006
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:41
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:270
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
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)
Add a boolean option and set &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
HYPRE_Int HYPRE_BigInt
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:65
int dim
Definition: ex24.cpp:53
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:324
Class for parallel bilinear form.
void Clear()
Clear the contents of the Mesh.
Definition: mesh.hpp:899
Vector data type.
Definition: vector.hpp:60
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7908
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:220
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4770
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:121
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
Class for parallel meshes.
Definition: pmesh.hpp:32
int main()
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150