MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
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
51using namespace std;
52using namespace mfem;
53
54int 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();
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 int slu_npdep = 1;
71
72 OptionsParser args(argc, argv);
73 args.AddOption(&mesh_file, "-m", "--mesh",
74 "Mesh file to use.");
75 args.AddOption(&order, "-o", "--order",
76 "Finite element order (polynomial degree) or -1 for"
77 " isoparametric space.");
78 args.AddOption(&device_config, "-d", "--device",
79 "Device configuration string, see Device::Configure().");
80 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
81 "--no-visualization",
82 "Enable or disable GLVis visualization.");
83 args.AddOption(&slu_colperm, "-cp", "--colperm",
84 "SuperLU Column Permutation Method: 0-NATURAL, 1-MMD-ATA "
85 "2-MMD_AT_PLUS_A, 3-COLAMD, 4-METIS_AT_PLUS_A, 5-PARMETIS "
86 "6-ZOLTAN");
87 args.AddOption(&slu_rowperm, "-rp", "--rowperm",
88 "SuperLU Row Permutation Method: 0-NOROWPERM, 1-LargeDiag");
89 args.AddOption(&slu_iterref, "-ir", "--iterref",
90 "SuperLU Iterative Refinement: 0-NOREFINE, 1-Single, "
91 "2-Double, 3-Extra");
92 args.AddOption(&slu_npdep, "-npdep", "--npdepth",
93 "Depth of 3D parition for SuperLU (>= 7.2.0)");
94
95 args.Parse();
96 if (!args.Good())
97 {
98 if (myid == 0)
99 {
100 args.PrintUsage(cout);
101 }
102 return 1;
103 }
104 if (myid == 0)
105 {
106 args.PrintOptions(cout);
107 }
108
109 // 3. Enable hardware devices such as GPUs, and programming models such as
110 // CUDA, OCCA, RAJA and OpenMP based on command line options.
111 Device device(device_config);
112 if (myid == 0) { device.Print(); }
113
114 // 4. Read the (serial) mesh from the given mesh file on all processors. We
115 // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
116 // and volume meshes with the same code.
117 Mesh mesh(mesh_file, 1, 1);
118 int dim = mesh.Dimension();
119
120 // 5. Refine the serial mesh on all processors to increase the resolution. In
121 // this example we do 'ref_levels' of uniform refinement. We choose
122 // 'ref_levels' to be the largest number that gives a final mesh with no
123 // more than 1,000 elements.
124 {
125 int ref_levels =
126 (int)floor(log(1000./mesh.GetNE())/log(2.)/dim);
127 for (int l = 0; l < ref_levels; l++)
128 {
129 mesh.UniformRefinement();
130 }
131 }
132
133 // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
134 // this mesh further in parallel to increase the resolution. Once the
135 // parallel mesh is defined, the serial mesh can be deleted.
136 ParMesh pmesh(MPI_COMM_WORLD, mesh);
137 mesh.Clear();
138 {
139 int par_ref_levels = 2;
140 for (int l = 0; l < par_ref_levels; l++)
141 {
142 pmesh.UniformRefinement();
143 }
144 }
145
146 // 7. Define a parallel finite element space on the parallel mesh. Here we
147 // use continuous Lagrange finite elements of the specified order. If
148 // order < 1, we instead use an isoparametric/isogeometric space.
150 bool delete_fec;
151 if (order > 0)
152 {
153 fec = new H1_FECollection(order, dim);
154 delete_fec = true;
155 }
156 else if (pmesh.GetNodes())
157 {
158 fec = pmesh.GetNodes()->OwnFEC();
159 delete_fec = false;
160 if (myid == 0)
161 {
162 cout << "Using isoparametric FEs: " << fec->Name() << endl;
163 }
164 }
165 else
166 {
167 fec = new H1_FECollection(order = 1, dim);
168 delete_fec = true;
169 }
170 ParFiniteElementSpace fespace(&pmesh, fec);
171 HYPRE_BigInt size = fespace.GlobalTrueVSize();
172 if (myid == 0)
173 {
174 cout << "Number of finite element unknowns: " << size << endl;
175 }
176
177 // 8. Determine the list of true (i.e. parallel conforming) essential
178 // boundary dofs. In this example, the boundary conditions are defined
179 // by marking all the boundary attributes from the mesh as essential
180 // (Dirichlet) and converting them to a list of true dofs.
181 Array<int> ess_tdof_list;
182 if (pmesh.bdr_attributes.Size())
183 {
184 Array<int> ess_bdr(pmesh.bdr_attributes.Max());
185 ess_bdr = 1;
186 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
187 }
188
189 // 9. Set up the parallel linear form b(.) which corresponds to the
190 // right-hand side of the FEM linear system, which in this case is
191 // (1,phi_i) where phi_i are the basis functions in fespace.
192 ParLinearForm b(&fespace);
193 ConstantCoefficient one(1.0);
194 b.AddDomainIntegrator(new DomainLFIntegrator(one));
195 b.Assemble();
196
197 // 10. Define the solution vector x as a parallel finite element grid function
198 // corresponding to fespace. Initialize x with initial guess of zero,
199 // which satisfies the boundary conditions.
200 ParGridFunction x(&fespace);
201 x = 0.0;
202
203 // 11. Set up the parallel bilinear form a(.,.) on the finite element space
204 // corresponding to the Laplacian operator -Delta, by adding the Diffusion
205 // domain integrator.
206 ParBilinearForm a(&fespace);
207 a.AddDomainIntegrator(new DiffusionIntegrator(one));
208
209 // 12. Assemble the parallel bilinear form and the corresponding linear
210 // system, applying any necessary transformations such as: parallel
211 // assembly, eliminating boundary conditions, applying conforming
212 // constraints for non-conforming AMR, static condensation, etc.
213 a.Assemble();
214
215 OperatorPtr A;
216 Vector B, X;
217 a.FormLinearSystem(ess_tdof_list, x, b, A, X, B);
218
219 // 13. Solve the linear system A X = B utilizing SuperLU.
220 SuperLUSolver *superlu = new SuperLUSolver(MPI_COMM_WORLD, slu_npdep);
221 Operator *SLU_A = new SuperLURowLocMatrix(*A.As<HypreParMatrix>());
222 superlu->SetPrintStatistics(true);
223 superlu->SetSymmetricPattern(false);
224
225 if (slu_colperm == 0)
226 {
228 }
229 else if (slu_colperm == 1)
230 {
232 }
233 else if (slu_colperm == 2)
234 {
236 }
237 else if (slu_colperm == 3)
238 {
240 }
241 else if (slu_colperm == 4)
242 {
244 }
245 else if (slu_colperm == 5)
246 {
248 }
249 else if (slu_colperm == 6)
250 {
252 }
253
254 if (slu_rowperm == 0)
255 {
257 }
258 else if (slu_rowperm == 1)
259 {
260#ifdef MFEM_USE_SUPERLU5
262#else
264#endif
265 }
266
267 if (slu_iterref == 0)
268 {
270 }
271 else if (slu_iterref == 1)
272 {
274 }
275 else if (slu_iterref == 2)
276 {
278 }
279 else if (slu_iterref == 3)
280 {
282 }
283
284 superlu->SetOperator(*SLU_A);
285 superlu->SetPrintStatistics(true);
286 superlu->Mult(B, X);
287
288 delete superlu;
289 delete SLU_A;
290
291 // 14. Recover the parallel grid function corresponding to X. This is the
292 // local finite element solution on each processor.
293 a.RecoverFEMSolution(X, b, x);
294
295 // 15. Save the refined mesh and the solution in parallel. This output can
296 // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
297 {
298 ostringstream mesh_name, sol_name;
299 mesh_name << "mesh." << setfill('0') << setw(6) << myid;
300 sol_name << "sol." << setfill('0') << setw(6) << myid;
301
302 ofstream mesh_ofs(mesh_name.str().c_str());
303 mesh_ofs.precision(8);
304 pmesh.Print(mesh_ofs);
305
306 ofstream sol_ofs(sol_name.str().c_str());
307 sol_ofs.precision(8);
308 x.Save(sol_ofs);
309 }
310
311 // 16. Send the solution by socket to a GLVis server.
312 if (visualization)
313 {
314 char vishost[] = "localhost";
315 int visport = 19916;
316 socketstream sol_sock(vishost, visport);
317 sol_sock << "parallel " << num_procs << " " << myid << "\n";
318 sol_sock.precision(8);
319 sol_sock << "solution\n" << pmesh << x << flush;
320 }
321
322 // 17. Free the used memory.
323 if (delete_fec)
324 {
325 delete fec;
326 }
327
328 return 0;
329}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
int Size() const
Return the logical size of the array.
Definition array.hpp:144
A coefficient that is constant across space and time.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition device.hpp:123
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition device.cpp:286
Class for domain integration .
Definition lininteg.hpp:109
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
virtual const char * Name() const
Definition fe_coll.hpp:79
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.hpp:74
Mesh data type.
Definition mesh.hpp:56
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:282
void Clear()
Clear the contents of the Mesh.
Definition mesh.hpp:730
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:8973
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
Pointer to an Operator of a specified type.
Definition handle.hpp:34
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get().
Definition handle.hpp:104
Abstract operator.
Definition operator.hpp:25
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear sys...
Definition operator.cpp:148
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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 'var' to receive the value. Enable/disable tags are used to set the bool...
Definition optparser.hpp:82
bool Good() const
Return true if the command line options were parsed successfully.
Class for parallel bilinear form.
Abstract parallel finite element space.
Definition pfespace.hpp:29
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
HYPRE_BigInt GlobalTrueVSize() const
Definition pfespace.hpp:285
Class for parallel grid function.
Definition pgridfunc.hpp:33
void Save(std::ostream &out) const override
Class for parallel linear form.
Class for parallel meshes.
Definition pmesh.hpp:34
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Definition pmesh.cpp:4801
void SetColumnPermutation(superlu::ColPerm col_perm)
Specify how to permute the columns of the matrix.
Definition superlu.cpp:399
void SetRowPermutation(superlu::RowPerm row_perm)
Specify how to permute the rows of the matrix.
Definition superlu.cpp:415
void Mult(const Vector &x, Vector &y) const
Factor and solve the linear system .
Definition superlu.cpp:587
void SetSymmetricPattern(bool sym)
Specify whether the matrix has a symmetric pattern to avoid extra work (default false)
Definition superlu.cpp:454
void SetOperator(const Operator &op)
Set the operator/matrix.
Definition superlu.cpp:475
void SetIterativeRefine(superlu::IterRefine iter_ref)
Specify how to handle iterative refinement.
Definition superlu.cpp:427
void SetPrintStatistics(bool print_stat)
Specify whether to print the solver statistics (default true)
Definition superlu.cpp:385
Vector data type.
Definition vector.hpp:80
int dim
Definition ex24.cpp:53
int main()
HYPRE_Int HYPRE_BigInt
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
@ NATURAL
Natural ordering.
Definition superlu.hpp:60
@ MMD_ATA
Minimum degree ordering on structure of .
Definition superlu.hpp:62
@ METIS_AT_PLUS_A
Sequential ordering on structure of using the METIS package.
Definition superlu.hpp:68
@ MMD_AT_PLUS_A
Minimum degree ordering on structure of .
Definition superlu.hpp:64
@ PARMETIS
Sequential ordering on structure of using the PARMETIS package.
Definition superlu.hpp:71
@ ZOLTAN
Use the Zoltan library from Sandia to define the column ordering.
Definition superlu.hpp:73
@ COLAMD
Approximate minimum degree column ordering.
Definition superlu.hpp:66
@ SLU_SINGLE
Iterative refinement accumulating residuals in a float.
Definition superlu.hpp:84
@ SLU_DOUBLE
Iterative refinement accumulating residuals in a double.
Definition superlu.hpp:86
@ SLU_EXTRA
Iterative refinement accumulating residuals in a higher precision variable.
Definition superlu.hpp:88
@ NOREFINE
No interative refinement.
Definition superlu.hpp:82
@ NOROWPERM
No row permutation.
Definition superlu.hpp:34
@ LargeDiag_MC64
Duff/Koster algorithm to make the diagonals large compared to the off-diagonals. Use LargeDiag for Su...
Definition superlu.hpp:46
const int visport
const char vishost[]