MFEM v4.8.0
Finite element discretization library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex4p.cpp
Go to the documentation of this file.
1// MFEM Example 4 - Parallel Version
2// PETSc Modification
3//
4// Compile with: make ex4p
5//
6// Sample runs:
7// mpirun -np 4 ex4p -m ../../data/klein-bottle.mesh -o 2 --petscopts rc_ex4p
8// mpirun -np 4 ex4p -m ../../data/klein-bottle.mesh -o 2 --petscopts rc_ex4p_bddc --nonoverlapping
9//
10// Description: This example code solves a simple 2D/3D H(div) diffusion
11// problem corresponding to the second order definite equation
12// -grad(alpha div F) + beta F = f with boundary condition F dot n
13// = <given normal field>. Here, we use a given exact solution F
14// and compute the corresponding r.h.s. f. We discretize with
15// Raviart-Thomas finite elements.
16//
17// The example demonstrates the use of H(div) finite element
18// spaces with the grad-div and H(div) vector finite element mass
19// bilinear form, as well as the computation of discretization
20// error when the exact solution is known. Bilinear form
21// hybridization and static condensation are also illustrated.
22//
23// We recommend viewing examples 1-3 before viewing this example.
24
25#include "mfem.hpp"
26#include <fstream>
27#include <iostream>
28
29#ifndef MFEM_USE_PETSC
30#error This example requires that MFEM is built with MFEM_USE_PETSC=YES
31#endif
32
33using namespace std;
34using namespace mfem;
35
36// Exact solution, F, and r.h.s., f. See below for implementation.
37void F_exact(const Vector &, Vector &);
38void f_exact(const Vector &, Vector &);
40
41int main(int argc, char *argv[])
42{
43 // 1. Initialize MPI and HYPRE.
44 Mpi::Init(argc, argv);
45 int num_procs = Mpi::WorldSize();
46 int myid = Mpi::WorldRank();
48
49 // 2. Parse command-line options.
50 const char *mesh_file = "../../data/star.mesh";
51 int ser_ref_levels = -1;
52 int par_ref_levels = 2;
53 int order = 1;
54 bool set_bc = true;
55 bool static_cond = false;
56 bool hybridization = false;
57 bool visualization = 1;
58 bool use_petsc = true;
59 const char *petscrc_file = "";
60 bool use_nonoverlapping = false;
61 const char *device_config = "cpu";
62
63 OptionsParser args(argc, argv);
64 args.AddOption(&mesh_file, "-m", "--mesh",
65 "Mesh file to use.");
66 args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
67 "Number of times to refine the mesh uniformly in serial.");
68 args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
69 "Number of times to refine the mesh uniformly in parallel.");
70 args.AddOption(&order, "-o", "--order",
71 "Finite element order (polynomial degree).");
72 args.AddOption(&set_bc, "-bc", "--impose-bc", "-no-bc", "--dont-impose-bc",
73 "Impose or not essential boundary conditions.");
74 args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
75 " solution.");
76 args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
77 "--no-static-condensation", "Enable static condensation.");
78 args.AddOption(&hybridization, "-hb", "--hybridization", "-no-hb",
79 "--no-hybridization", "Enable hybridization.");
80 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
81 "--no-visualization",
82 "Enable or disable GLVis visualization.");
83 args.AddOption(&use_petsc, "-usepetsc", "--usepetsc", "-no-petsc",
84 "--no-petsc",
85 "Use or not PETSc to solve the linear system.");
86 args.AddOption(&petscrc_file, "-petscopts", "--petscopts",
87 "PetscOptions file to use.");
88 args.AddOption(&use_nonoverlapping, "-nonoverlapping", "--nonoverlapping",
89 "-no-nonoverlapping", "--no-nonoverlapping",
90 "Use or not the block diagonal PETSc's matrix format "
91 "for non-overlapping domain decomposition.");
92 args.AddOption(&device_config, "-d", "--device",
93 "Device configuration string, see Device::Configure().");
94 args.Parse();
95 if (!args.Good())
96 {
97 if (myid == 0)
98 {
99 args.PrintUsage(cout);
100 }
101 return 1;
102 }
103 if (myid == 0)
104 {
105 args.PrintOptions(cout);
106 }
107 kappa = freq * M_PI;
108
109 // 2b. 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 // 2c. We initialize PETSc
115 if (use_petsc) { MFEMInitializePetsc(NULL,NULL,petscrc_file,NULL); }
116
117 // 3. Read the (serial) mesh from the given mesh file on all processors. We
118 // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
119 // and volume, as well as periodic meshes with the same code.
120 Mesh *mesh = new Mesh(mesh_file, 1, 1);
121 int dim = mesh->Dimension();
122 int sdim = mesh->SpaceDimension();
123
124 // 4. Refine the serial mesh on all processors to increase the resolution. In
125 // this example we do 'ref_levels' of uniform refinement. We choose
126 // 'ref_levels' to be the largest number that gives a final mesh with no
127 // more than 1,000 elements.
128 {
129 if (ser_ref_levels < 0)
130 {
131 ser_ref_levels = (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
132 }
133 for (int l = 0; l < ser_ref_levels; l++)
134 {
135 mesh->UniformRefinement();
136 }
137 }
138
139 // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
140 // this mesh further in parallel to increase the resolution. Once the
141 // parallel mesh is defined, the serial mesh can be deleted.
142 ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
143 delete mesh;
144 {
145 for (int l = 0; l < par_ref_levels; l++)
146 {
147 pmesh->UniformRefinement();
148 }
149 }
150
151 // 6. Define a parallel finite element space on the parallel mesh. Here we
152 // use the Raviart-Thomas finite elements of the specified order.
153 FiniteElementCollection *fec = new RT_FECollection(order-1, dim);
154 ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
155 HYPRE_BigInt size = fespace->GlobalTrueVSize();
156 if (myid == 0)
157 {
158 cout << "Number of finite element unknowns: " << size << endl;
159 }
160
161 // 7. Determine the list of true (i.e. parallel conforming) essential
162 // boundary dofs. In this example, the boundary conditions are defined
163 // by marking all the boundary attributes from the mesh as essential
164 // (Dirichlet) and converting them to a list of true dofs.
165 Array<int> ess_tdof_list;
166 if (pmesh->bdr_attributes.Size())
167 {
168 Array<int> ess_bdr(pmesh->bdr_attributes.Max());
169 ess_bdr = set_bc ? 1 : 0;
170 fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
171 }
172
173 // 8. Set up the parallel linear form b(.) which corresponds to the
174 // right-hand side of the FEM linear system, which in this case is
175 // (f,phi_i) where f is given by the function f_exact and phi_i are the
176 // basis functions in the finite element fespace.
178 ParLinearForm *b = new ParLinearForm(fespace);
179 b->AddDomainIntegrator(new VectorFEDomainLFIntegrator(f));
180 b->Assemble();
181
182 // 9. Define the solution vector x as a parallel finite element grid function
183 // corresponding to fespace. Initialize x by projecting the exact
184 // solution. Note that only values from the boundary faces will be used
185 // when eliminating the non-homogeneous boundary condition to modify the
186 // r.h.s. vector b.
187 ParGridFunction x(fespace);
190
191 // 10. Set up the parallel bilinear form corresponding to the H(div)
192 // diffusion operator grad alpha div + beta I, by adding the div-div and
193 // the mass domain integrators.
196 ParBilinearForm *a = new ParBilinearForm(fespace);
197 a->AddDomainIntegrator(new DivDivIntegrator(*alpha));
198 a->AddDomainIntegrator(new VectorFEMassIntegrator(*beta));
199
200 // 11. Assemble the parallel bilinear form and the corresponding linear
201 // system, applying any necessary transformations such as: parallel
202 // assembly, eliminating boundary conditions, applying conforming
203 // constraints for non-conforming AMR, static condensation,
204 // hybridization, etc.
205 FiniteElementCollection *hfec = NULL;
206 ParFiniteElementSpace *hfes = NULL;
207 if (static_cond)
208 {
209 a->EnableStaticCondensation();
210 }
211 else if (hybridization)
212 {
213 hfec = new DG_Interface_FECollection(order-1, dim);
214 hfes = new ParFiniteElementSpace(pmesh, hfec);
215 a->EnableHybridization(hfes, new NormalTraceJumpIntegrator(),
216 ess_tdof_list);
217 }
218 a->Assemble();
219
220 Vector B, X;
221 CGSolver *pcg = new CGSolver(MPI_COMM_WORLD);
222 pcg->SetRelTol(1e-12);
223 pcg->SetMaxIter(500);
224 pcg->SetPrintLevel(1);
225
226 if (!use_petsc)
227 {
229 a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
230
231 HYPRE_BigInt glob_size = A.GetGlobalNumRows();
232 if (myid == 0)
233 {
234 cout << "Size of linear system: " << glob_size << endl;
235 }
236
237 // 12. Define and apply a parallel PCG solver for A X = B with the 2D AMS or
238 // the 3D ADS preconditioners from hypre. If using hybridization, the
239 // system is preconditioned with hypre's BoomerAMG.
240 HypreSolver *prec = NULL;
241 pcg->SetOperator(A);
242 if (hybridization) { prec = new HypreBoomerAMG(A); }
243 else
244 {
245 ParFiniteElementSpace *prec_fespace =
246 (a->StaticCondensationIsEnabled() ? a->SCParFESpace() : fespace);
247 if (dim == 2) { prec = new HypreAMS(A, prec_fespace); }
248 else { prec = new HypreADS(A, prec_fespace); }
249 }
250 pcg->SetPreconditioner(*prec);
251 pcg->Mult(B, X);
252 delete prec;
253 }
254 else
255 {
257 PetscPreconditioner *prec = NULL;
258 a->SetOperatorType(use_nonoverlapping ?
260 a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
261
262 if (myid == 0)
263 {
264 cout << "Size of linear system: " << A.M() << endl;
265 }
266
267 pcg->SetOperator(A);
268 if (use_nonoverlapping)
269 {
270 ParFiniteElementSpace *prec_fespace =
271 (a->StaticCondensationIsEnabled() ? a->SCParFESpace() :
272 (hfes ? NULL : fespace));
273
274 // Auxiliary class for BDDC customization
276 // Inform the solver about the finite element space
277 opts.SetSpace(prec_fespace);
278 // Inform the solver about essential dofs
279 opts.SetEssBdrDofs(&ess_tdof_list);
280 // Create a BDDC solver with parameters
281 prec = new PetscBDDCSolver(A, opts);
282 }
283 else
284 {
285 // Create an empty preconditioner that can be customized at runtime.
286 prec = new PetscPreconditioner(A, "solver_");
287 }
288 pcg->SetPreconditioner(*prec);
289 pcg->Mult(B, X);
290 delete prec;
291 }
292 delete pcg;
293
294 // 13. Recover the parallel grid function corresponding to X. This is the
295 // local finite element solution on each processor.
296 a->RecoverFEMSolution(X, *b, x);
297
298 // 14. Compute and print the L^2 norm of the error.
299 {
301 if (myid == 0)
302 {
303 cout << "\n|| F_h - F ||_{L^2} = " << err << '\n' << endl;
304 }
305 }
306
307 // 15. Save the refined mesh and the solution in parallel. This output can
308 // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
309 {
310 ostringstream mesh_name, sol_name;
311 mesh_name << "mesh." << setfill('0') << setw(6) << myid;
312 sol_name << "sol." << setfill('0') << setw(6) << myid;
313
314 ofstream mesh_ofs(mesh_name.str().c_str());
315 mesh_ofs.precision(8);
316 pmesh->Print(mesh_ofs);
317
318 ofstream sol_ofs(sol_name.str().c_str());
319 sol_ofs.precision(8);
320 x.Save(sol_ofs);
321 }
322
323 // 16. Send the solution by socket to a GLVis server.
324 if (visualization)
325 {
326 char vishost[] = "localhost";
327 int visport = 19916;
328 socketstream sol_sock(vishost, visport);
329 sol_sock << "parallel " << num_procs << " " << myid << "\n";
330 sol_sock.precision(8);
331 sol_sock << "solution\n" << *pmesh << x << flush;
332 }
333
334 // 17. Free the used memory.
335 delete hfes;
336 delete hfec;
337 delete a;
338 delete alpha;
339 delete beta;
340 delete b;
341 delete fespace;
342 delete fec;
343 delete pmesh;
344
345 // We finalize PETSc
346 if (use_petsc) { MFEMFinalizePetsc(); }
347
348 return 0;
349}
350
351
352// The exact solution (for non-surface meshes)
353void F_exact(const Vector &p, Vector &F)
354{
355 int dim = p.Size();
356
357 real_t x = p(0);
358 real_t y = p(1);
359 // real_t z = (dim == 3) ? p(2) : 0.0;
360
361 F(0) = cos(kappa*x)*sin(kappa*y);
362 F(1) = cos(kappa*y)*sin(kappa*x);
363 if (dim == 3)
364 {
365 F(2) = 0.0;
366 }
367}
368
369// The right hand side
370void f_exact(const Vector &p, Vector &f)
371{
372 int dim = p.Size();
373
374 real_t x = p(0);
375 real_t y = p(1);
376 // real_t z = (dim == 3) ? p(2) : 0.0;
377
378 real_t temp = 1 + 2*kappa*kappa;
379
380 f(0) = temp*cos(kappa*x)*sin(kappa*y);
381 f(1) = temp*cos(kappa*y)*sin(kappa*x);
382 if (dim == 3)
383 {
384 f(2) = 0;
385 }
386}
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:147
Conjugate gradient method.
Definition solvers.hpp:538
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
Definition solvers.cpp:751
void SetOperator(const Operator &op) override
Set/update the solver for the given operator.
Definition solvers.hpp:551
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
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:297
for Raviart-Thomas elements
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
The Auxiliary-space Divergence Solver in hypre.
Definition hypre.hpp:1948
The Auxiliary-space Maxwell Solver in hypre.
Definition hypre.hpp:1871
The BoomerAMG solver in hypre.
Definition hypre.hpp:1717
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:408
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition hypre.hpp:699
Abstract class for hypre's solvers and preconditioners.
Definition hypre.hpp:1188
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.cpp:33
void SetRelTol(real_t rtol)
Definition solvers.hpp:229
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition solvers.cpp:174
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition solvers.cpp:72
void SetMaxIter(int max_it)
Definition solvers.hpp:231
Mesh data type.
Definition mesh.hpp:64
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:290
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1282
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1216
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1219
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:11295
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).
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
Form a constrained linear system using a matrix-free approach.
Definition operator.cpp:114
@ PETSC_MATIS
ID for class PetscParMatrix, MATIS format.
Definition operator.hpp:289
@ PETSC_MATAIJ
ID for class PetscParMatrix, MATAIJ format.
Definition operator.hpp:288
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:346
Class for parallel grid function.
Definition pgridfunc.hpp:50
void Save(std::ostream &out) const override
real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
Returns ||u_ex - u_h||_L2 in parallel for H1 or L2 elements.
void ProjectCoefficient(Coefficient &coeff) override
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
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:4800
Auxiliary class for BDDC customization.
Definition petsc.hpp:831
void SetEssBdrDofs(const Array< int > *essdofs, bool loc=false)
Specify dofs on the essential boundary.
Definition petsc.hpp:850
void SetSpace(ParFiniteElementSpace *fe)
Definition petsc.hpp:845
Wrapper for PETSc's matrix class.
Definition petsc.hpp:319
PetscInt M() const
Returns the global number of rows.
Definition petsc.cpp:1010
Abstract class for PETSc's preconditioners.
Definition petsc.hpp:805
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition fe_coll.hpp:403
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition lininteg.hpp:344
A general vector function coefficient.
Vector data type.
Definition vector.hpp:82
int Size() const
Returns the size of the vector.
Definition vector.hpp:226
Vector beta
const real_t alpha
Definition ex15.cpp:369
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
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
Definition petsc.cpp:206
void MFEMFinalizePetsc()
Definition petsc.cpp:246
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
Definition globals.hpp:71
float real_t
Definition config.hpp:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
const char vishost[]
real_t p(const Vector &x, real_t t)
real_t kappa
Definition ex4p.cpp:39
void f_exact(const Vector &, Vector &)
Definition ex4p.cpp:370
real_t freq
Definition ex4p.cpp:39
void F_exact(const Vector &, Vector &)
Definition ex4p.cpp:353