MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex14p.cpp
Go to the documentation of this file.
1// MFEM Example 14 - Parallel Version
2//
3// Compile with: make ex14p
4//
5// Sample runs: mpirun -np 4 ex14p -m ../data/inline-quad.mesh -o 0
6// mpirun -np 4 ex14p -m ../data/star.mesh -o 2
7// mpirun -np 4 ex14p -m ../data/star-mixed.mesh -o 2
8// mpirun -np 4 ex14p -m ../data/star-mixed.mesh -o 2 -k 0 -e 1
9// mpirun -np 4 ex14p -m ../data/escher.mesh -s 1
10// mpirun -np 4 ex14p -m ../data/fichera.mesh -s 1 -k 1
11// mpirun -np 4 ex14p -m ../data/fichera-mixed.mesh -s 1 -k 1
12// mpirun -np 4 ex14p -m ../data/square-disc-p2.vtk -o 2
13// mpirun -np 4 ex14p -m ../data/square-disc-p3.mesh -o 3
14// mpirun -np 4 ex14p -m ../data/square-disc-nurbs.mesh -o 1
15// mpirun -np 4 ex14p -m ../data/disc-nurbs.mesh -rs 4 -o 2 -s 1 -k 0
16// mpirun -np 4 ex14p -m ../data/pipe-nurbs.mesh -o 1
17// mpirun -np 4 ex14p -m ../data/inline-segment.mesh -rs 5
18// mpirun -np 4 ex14p -m ../data/amr-quad.mesh -rs 3
19// mpirun -np 4 ex14p -m ../data/amr-hex.mesh
20// mpirun -np 4 ex14p -pa -rs 1 -rp 0 -o 3
21// mpirun -np 4 ex14p -pa -rs 1 -rp 0 -m ../data/fichera.mesh -o 3
22//
23// Device sample runs:
24// mpirun -np 4 ex14p -pa -rs 2 -rp 0 -d cuda -o 3
25// mpirun -np 4 ex14p -pa -rs 2 -rp 0 -d cuda -m ../data/fichera.mesh -o 3
26//
27// Description: This example code demonstrates the use of MFEM to define a
28// discontinuous Galerkin (DG) finite element discretization of
29// the Laplace problem -Delta u = 1 with homogeneous Dirichlet
30// boundary conditions. Finite element spaces of any order,
31// including zero on regular grids, are supported. The example
32// highlights the use of discontinuous spaces and DG-specific face
33// integrators.
34//
35// We recommend viewing examples 1 and 9 before viewing this
36// example.
37
38#include "mfem.hpp"
39#include <fstream>
40#include <iostream>
41
42using namespace std;
43using namespace mfem;
44
45class CustomSolverMonitor : public IterativeSolverMonitor
46{
47private:
48 const ParMesh &pmesh;
49 ParGridFunction &pgf;
50public:
51 CustomSolverMonitor(const ParMesh &pmesh_,
52 ParGridFunction &pgf_) :
53 pmesh(pmesh_),
54 pgf(pgf_) {}
55
56 void MonitorSolution(int i, real_t norm, const Vector &x, bool final)
57 {
58 char vishost[] = "localhost";
59 int visport = 19916;
60 int num_procs, myid;
61
62 MPI_Comm_size(pmesh.GetComm(), &num_procs);
63 MPI_Comm_rank(pmesh.GetComm(), &myid);
64
65 pgf.SetFromTrueDofs(x);
66
67 socketstream sol_sock(vishost, visport);
68 sol_sock << "parallel " << num_procs << " " << myid << "\n";
69 sol_sock.precision(8);
70 sol_sock << "solution\n" << pmesh << pgf
71 << "window_title 'Iteration no " << i << "'"
72 << "keys rRjlc\n" << flush;
73 }
74};
75
76int main(int argc, char *argv[])
77{
78 // 1. Initialize MPI and HYPRE.
79 Mpi::Init(argc, argv);
81
82 // 2. Parse command-line options.
83 const char *mesh_file = "../data/star.mesh";
84 int ser_ref_levels = -1;
85 int par_ref_levels = 2;
86 int order = 1;
87 real_t sigma = -1.0;
88 real_t kappa = -1.0;
89 real_t eta = 0.0;
90 bool pa = false;
91 bool visualization = 1;
92 const char *device_config = "cpu";
93
94 OptionsParser args(argc, argv);
95 args.AddOption(&mesh_file, "-m", "--mesh",
96 "Mesh file to use.");
97 args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
98 "Number of times to refine the mesh uniformly in serial,"
99 " -1 for auto.");
100 args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
101 "Number of times to refine the mesh uniformly in parallel.");
102 args.AddOption(&order, "-o", "--order",
103 "Finite element order (polynomial degree) >= 0.");
104 args.AddOption(&sigma, "-s", "--sigma",
105 "One of the three DG penalty parameters, typically +1/-1."
106 " See the documentation of class DGDiffusionIntegrator.");
107 args.AddOption(&kappa, "-k", "--kappa",
108 "One of the three DG penalty parameters, should be positive."
109 " Negative values are replaced with (order+1)^2.");
110 args.AddOption(&eta, "-e", "--eta", "BR2 penalty parameter.");
111 args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
112 "--no-partial-assembly", "Enable Partial Assembly.");
113 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
114 "--no-visualization",
115 "Enable or disable GLVis visualization.");
116 args.AddOption(&device_config, "-d", "--device",
117 "Device configuration string, see Device::Configure().");
118 args.Parse();
119 if (!args.Good())
120 {
121 if (Mpi::Root())
122 {
123 args.PrintUsage(cout);
124 }
125 return 1;
126 }
127 if (kappa < 0)
128 {
129 kappa = (order+1)*(order+1);
130 }
131 if (Mpi::Root())
132 {
133 args.PrintOptions(cout);
134 }
135
136 Device device(device_config);
137 if (Mpi::Root()) { device.Print(); }
138
139 // 3. Read the (serial) mesh from the given mesh file on all processors. We
140 // can handle triangular, quadrilateral, tetrahedral and hexahedral meshes
141 // with the same code. NURBS meshes are projected to second order meshes.
142 Mesh mesh(mesh_file);
143 int dim = mesh.Dimension();
144
145 // 4. Refine the serial mesh on all processors to increase the resolution. In
146 // this example we do 'ser_ref_levels' of uniform refinement. By default,
147 // or if ser_ref_levels < 0, we choose it to be the largest number that
148 // gives a final mesh with no more than 50,000 elements.
149 {
150 if (ser_ref_levels < 0)
151 {
152 ser_ref_levels = (int)floor(log(10000./mesh.GetNE())/log(2.)/dim);
153 }
154 for (int l = 0; l < ser_ref_levels; l++)
155 {
156 mesh.UniformRefinement();
157 }
158 }
159 if (mesh.NURBSext)
160 {
161 mesh.SetCurvature(max(order, 1));
162 }
163
164 // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
165 // this mesh further in parallel to increase the resolution. Once the
166 // parallel mesh is defined, the serial mesh can be deleted.
167 ParMesh pmesh(MPI_COMM_WORLD, mesh);
168 mesh.Clear();
169 {
170 for (int l = 0; l < par_ref_levels; l++)
171 {
172 pmesh.UniformRefinement();
173 }
174 }
175
176 // 6. Define a parallel finite element space on the parallel mesh. Here we
177 // use discontinuous finite elements of the specified order >= 0.
179 DG_FECollection fec(order, dim, bt);
180 ParFiniteElementSpace fespace(&pmesh, &fec);
181 HYPRE_BigInt size = fespace.GlobalTrueVSize();
182 if (Mpi::Root())
183 {
184 cout << "Number of unknowns: " << size << endl;
185 }
186
187 // 7. Set up the parallel linear form b(.) which corresponds to the
188 // right-hand side of the FEM linear system.
189 ParLinearForm b(&fespace);
190 ConstantCoefficient one(1.0);
191 ConstantCoefficient zero(0.0);
192 b.AddDomainIntegrator(new DomainLFIntegrator(one));
193 b.AddBdrFaceIntegrator(
194 new DGDirichletLFIntegrator(zero, one, sigma, kappa));
195 b.Assemble();
196
197 // 8. Define the solution vector x as a parallel finite element grid function
198 // corresponding to fespace. Initialize x with initial guess of zero.
199 ParGridFunction x(&fespace);
200 x = 0.0;
201
202 // 9. Set up the bilinear form a(.,.) on the finite element space
203 // corresponding to the Laplacian operator -Delta, by adding the Diffusion
204 // domain integrator and the interior and boundary DG face integrators.
205 // Note that boundary conditions are imposed weakly in the form, so there
206 // is no need for dof elimination. After serial and parallel assembly we
207 // extract the corresponding parallel matrix A.
208 ParBilinearForm a(&fespace);
209 a.AddDomainIntegrator(new DiffusionIntegrator(one));
210 a.AddInteriorFaceIntegrator(new DGDiffusionIntegrator(one, sigma, kappa));
211 a.AddBdrFaceIntegrator(new DGDiffusionIntegrator(one, sigma, kappa));
212 if (eta > 0)
213 {
214 MFEM_VERIFY(!pa, "BR2 not yet compatible with partial assembly.");
215 a.AddInteriorFaceIntegrator(new DGDiffusionBR2Integrator(fespace, eta));
216 a.AddBdrFaceIntegrator(new DGDiffusionBR2Integrator(fespace, eta));
217 }
218 if (pa) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
219 a.Assemble();
220 a.Finalize();
221
222 // 10. Define the parallel (hypre) matrix and vectors representing a(.,.),
223 // b(.) and the finite element approximation.
225
226 std::unique_ptr<HypreBoomerAMG> amg;
227 if (pa)
228 {
229 A.Reset(&a, false);
230 }
231 else
232 {
234 a.ParallelAssemble(A);
235 amg.reset(new HypreBoomerAMG(*A.As<HypreParMatrix>()));
236 }
237
238 // 11. Depending on the symmetry of A, define and apply a parallel PCG or
239 // GMRES solver for AX=B using the BoomerAMG preconditioner from hypre.
240 if (sigma == -1.0)
241 {
242 CGSolver cg(MPI_COMM_WORLD);
243 cg.SetRelTol(1e-12);
244 cg.SetMaxIter(500);
245 cg.SetPrintLevel(1);
246 cg.SetOperator(*A);
247 if (amg) { cg.SetPreconditioner(*amg); }
248 cg.Mult(b, x);
249 }
250 else
251 {
252 CustomSolverMonitor monitor(pmesh, x);
253 GMRESSolver gmres(MPI_COMM_WORLD);
254 gmres.SetAbsTol(0.0);
255 gmres.SetRelTol(1e-12);
256 gmres.SetMaxIter(500);
257 gmres.SetKDim(10);
258 gmres.SetPrintLevel(1);
259 gmres.SetOperator(*A);
260 if (amg) { gmres.SetPreconditioner(*amg); }
261 gmres.SetMonitor(monitor);
262 gmres.Mult(b, x);
263 }
264
265 // 12. Save the refined mesh and the solution in parallel. This output can
266 // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
267 {
268 ostringstream mesh_name, sol_name;
269 mesh_name << "mesh." << setfill('0') << setw(6) << Mpi::WorldRank();
270 sol_name << "sol." << setfill('0') << setw(6) << Mpi::WorldRank();
271
272 ofstream mesh_ofs(mesh_name.str().c_str());
273 mesh_ofs.precision(8);
274 pmesh.Print(mesh_ofs);
275
276 ofstream sol_ofs(sol_name.str().c_str());
277 sol_ofs.precision(8);
278 x.Save(sol_ofs);
279 }
280
281 // 13. Send the solution by socket to a GLVis server.
282 if (visualization)
283 {
284 char vishost[] = "localhost";
285 int visport = 19916;
286 socketstream sol_sock(vishost, visport);
287 sol_sock << "parallel " << Mpi::WorldSize() << " " << Mpi::WorldRank() << "\n";
288 sol_sock.precision(8);
289 sol_sock << "solution\n" << pmesh << x << flush;
290 }
291
292 return 0;
293}
@ GaussLobatto
Closed type.
Definition fe_base.hpp:32
@ GaussLegendre
Open type.
Definition fe_base.hpp:31
Conjugate gradient method.
Definition solvers.hpp:513
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition solvers.hpp:526
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the Conjugate Gradient method.
Definition solvers.cpp:718
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
GMRES method.
Definition solvers.hpp:547
void SetKDim(int dim)
Set the number of iteration to perform between restarts, default is 50.
Definition solvers.hpp:559
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the GMRES method.
Definition solvers.cpp:980
The BoomerAMG solver in hypre.
Definition hypre.hpp:1691
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
Abstract base class for an iterative solver monitor.
Definition solvers.hpp:37
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
Definition solvers.cpp:179
void SetRelTol(real_t rtol)
Definition solvers.hpp:209
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition solvers.cpp:173
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition solvers.cpp:71
void SetMaxIter(int max_it)
Definition solvers.hpp:211
void SetMonitor(IterativeSolverMonitor &m)
Set the iterative solver monitor.
Definition solvers.hpp:296
void SetAbsTol(real_t atol)
Definition solvers.hpp:210
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
Mesh data type.
Definition mesh.hpp:56
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition mesh.hpp:290
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
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
Definition mesh.cpp:6211
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
static bool Root()
Return true if the rank in MPI_COMM_WORLD is zero.
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
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition handle.hpp:132
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition handle.hpp:145
@ Hypre_ParCSR
ID for class HypreParMatrix.
Definition operator.hpp:287
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
HYPRE_BigInt GlobalTrueVSize() const
Definition pfespace.hpp:285
Class for parallel grid function.
Definition pgridfunc.hpp:33
void Save(std::ostream &out) const override
void SetFromTrueDofs(const Vector &tv) override
Set the GridFunction from the given true-dof vector.
Class for parallel linear form.
Class for parallel meshes.
Definition pmesh.hpp:34
MPI_Comm GetComm() const
Definition pmesh.hpp:402
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Definition pmesh.cpp:4801
Vector data type.
Definition vector.hpp:80
real_t sigma(const Vector &x)
Definition maxwell.cpp:102
real_t kappa
Definition ex24.cpp:54
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
const int visport
float real_t
Definition config.hpp:43
const char vishost[]