MFEM  v4.4.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 //
3 // Compile with: make ex4p
4 //
5 // Sample runs: mpirun -np 4 ex4p -m ../data/square-disc.mesh
6 // mpirun -np 4 ex4p -m ../data/star.mesh
7 // mpirun -np 4 ex4p -m ../data/beam-tet.mesh
8 // mpirun -np 4 ex4p -m ../data/beam-hex.mesh
9 // mpirun -np 4 ex4p -m ../data/beam-hex.mesh -o 2 -pa
10 // mpirun -np 4 ex4p -m ../data/escher.mesh -o 2 -sc
11 // mpirun -np 4 ex4p -m ../data/fichera.mesh -o 2 -hb
12 // mpirun -np 4 ex4p -m ../data/fichera-q2.vtk
13 // mpirun -np 4 ex4p -m ../data/fichera-q3.mesh -o 2 -sc
14 // mpirun -np 4 ex4p -m ../data/square-disc-nurbs.mesh -o 3
15 // mpirun -np 4 ex4p -m ../data/beam-hex-nurbs.mesh -o 3
16 // mpirun -np 4 ex4p -m ../data/periodic-square.mesh -no-bc
17 // mpirun -np 4 ex4p -m ../data/periodic-cube.mesh -no-bc
18 // mpirun -np 4 ex4p -m ../data/amr-quad.mesh
19 // mpirun -np 3 ex4p -m ../data/amr-quad.mesh -o 2 -hb
20 // mpirun -np 4 ex4p -m ../data/amr-hex.mesh -o 2 -sc
21 // mpirun -np 4 ex4p -m ../data/amr-hex.mesh -o 2 -hb
22 // mpirun -np 4 ex4p -m ../data/ref-prism.mesh -o 1
23 // mpirun -np 4 ex4p -m ../data/octahedron.mesh -o 1
24 // mpirun -np 4 ex4p -m ../data/star-surf.mesh -o 3 -hb
25 //
26 // Device sample runs:
27 // mpirun -np 4 ex4p -m ../data/star.mesh -pa -d cuda
28 // mpirun -np 4 ex4p -m ../data/star.mesh -pa -d raja-cuda
29 // mpirun -np 4 ex4p -m ../data/star.mesh -pa -d raja-omp
30 // mpirun -np 4 ex4p -m ../data/beam-hex.mesh -pa -d cuda
31 //
32 // Description: This example code solves a simple 2D/3D H(div) diffusion
33 // problem corresponding to the second order definite equation
34 // -grad(alpha div F) + beta F = f with boundary condition F dot n
35 // = <given normal field>. Here, we use a given exact solution F
36 // and compute the corresponding r.h.s. f. We discretize with
37 // Raviart-Thomas finite elements.
38 //
39 // The example demonstrates the use of H(div) finite element
40 // spaces with the grad-div and H(div) vector finite element mass
41 // bilinear form, as well as the computation of discretization
42 // error when the exact solution is known. Bilinear form
43 // hybridization and static condensation are also illustrated.
44 //
45 // We recommend viewing examples 1-3 before viewing this example.
46 
47 #include "mfem.hpp"
48 #include <fstream>
49 #include <iostream>
50 
51 using namespace std;
52 using namespace mfem;
53 
54 // Exact solution, F, and r.h.s., f. See below for implementation.
55 void F_exact(const Vector &, Vector &);
56 void f_exact(const Vector &, Vector &);
57 double freq = 1.0, kappa;
58 
59 int main(int argc, char *argv[])
60 {
61  // 1. Initialize MPI and HYPRE.
62  Mpi::Init(argc, argv);
63  int num_procs = Mpi::WorldSize();
64  int myid = Mpi::WorldRank();
65  Hypre::Init();
66 
67  // 2. Parse command-line options.
68  const char *mesh_file = "../data/star.mesh";
69  int order = 1;
70  bool set_bc = true;
71  bool static_cond = false;
72  bool hybridization = false;
73  bool pa = false;
74  const char *device_config = "cpu";
75  bool visualization = 1;
76 
77  OptionsParser args(argc, argv);
78  args.AddOption(&mesh_file, "-m", "--mesh",
79  "Mesh file to use.");
80  args.AddOption(&order, "-o", "--order",
81  "Finite element order (polynomial degree).");
82  args.AddOption(&set_bc, "-bc", "--impose-bc", "-no-bc", "--dont-impose-bc",
83  "Impose or not essential boundary conditions.");
84  args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
85  " solution.");
86  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
87  "--no-static-condensation", "Enable static condensation.");
88  args.AddOption(&hybridization, "-hb", "--hybridization", "-no-hb",
89  "--no-hybridization", "Enable hybridization.");
90  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
91  "--no-partial-assembly", "Enable Partial Assembly.");
92  args.AddOption(&device_config, "-d", "--device",
93  "Device configuration string, see Device::Configure().");
94  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
95  "--no-visualization",
96  "Enable or disable GLVis visualization.");
97  args.Parse();
98  if (!args.Good())
99  {
100  if (myid == 0)
101  {
102  args.PrintUsage(cout);
103  }
104  return 1;
105  }
106  if (myid == 0)
107  {
108  args.PrintOptions(cout);
109  }
110  kappa = freq * M_PI;
111 
112  // 3. Enable hardware devices such as GPUs, and programming models such as
113  // CUDA, OCCA, RAJA and OpenMP based on command line options.
114  Device device(device_config);
115  if (myid == 0) { device.Print(); }
116 
117  // 4. 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  // 5. 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  int ref_levels =
130  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
131  for (int l = 0; l < ref_levels; l++)
132  {
133  mesh->UniformRefinement();
134  }
135  }
136 
137  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
138  // this mesh further in parallel to increase the resolution. Once the
139  // parallel mesh is defined, the serial mesh can be deleted.
140  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
141  delete mesh;
142  {
143  int par_ref_levels = 2;
144  for (int l = 0; l < par_ref_levels; l++)
145  {
146  pmesh->UniformRefinement();
147  }
148  }
149 
150  // 7. Define a parallel finite element space on the parallel mesh. Here we
151  // use the Raviart-Thomas finite elements of the specified order.
152  FiniteElementCollection *fec = new RT_FECollection(order-1, dim);
153  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
154  HYPRE_BigInt size = fespace->GlobalTrueVSize();
155  if (myid == 0)
156  {
157  cout << "Number of finite element unknowns: " << size << endl;
158  }
159 
160  // 8. Determine the list of true (i.e. parallel conforming) essential
161  // boundary dofs. In this example, the boundary conditions are defined
162  // by marking all the boundary attributes from the mesh as essential
163  // (Dirichlet) and converting them to a list of true dofs.
164  Array<int> ess_tdof_list;
165  if (pmesh->bdr_attributes.Size())
166  {
167  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
168  ess_bdr = set_bc ? 1 : 0;
169  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
170  }
171 
172  // 9. Set up the parallel linear form b(.) which corresponds to the
173  // right-hand side of the FEM linear system, which in this case is
174  // (f,phi_i) where f is given by the function f_exact and phi_i are the
175  // basis functions in the finite element fespace.
177  ParLinearForm *b = new ParLinearForm(fespace);
179  b->Assemble();
180 
181  // 10. Define the solution vector x as a parallel finite element grid function
182  // corresponding to fespace. Initialize x by projecting the exact
183  // solution. Note that only values from the boundary faces will be used
184  // when eliminating the non-homogeneous boundary condition to modify the
185  // r.h.s. vector b.
186  ParGridFunction x(fespace);
188  x.ProjectCoefficient(F);
189 
190  // 11. Set up the parallel bilinear form corresponding to the H(div)
191  // diffusion operator grad alpha div + beta I, by adding the div-div and
192  // the mass domain integrators.
194  Coefficient *beta = new ConstantCoefficient(1.0);
195  ParBilinearForm *a = new ParBilinearForm(fespace);
196  if (pa) { a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
197  a->AddDomainIntegrator(new DivDivIntegrator(*alpha));
199 
200  // 12. 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  {
210  }
211  else if (hybridization)
212  {
213  hfec = new DG_Interface_FECollection(order-1, dim);
214  hfes = new ParFiniteElementSpace(pmesh, hfec);
216  ess_tdof_list);
217  }
218  a->Assemble();
219 
220  OperatorPtr A;
221  Vector B, X;
222  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
223 
224  if (myid == 0 && !pa)
225  {
226  cout << "Size of linear system: "
227  << A.As<HypreParMatrix>()->GetGlobalNumRows() << endl;
228  }
229 
230  // 13. Define and apply a parallel PCG solver for A X = B with the 2D AMS or
231  // the 3D ADS preconditioners from hypre. If using hybridization, the
232  // system is preconditioned with hypre's BoomerAMG. In the partial
233  // assembly case, use Jacobi preconditioning.
234  Solver *prec = NULL;
235  CGSolver *pcg = new CGSolver(MPI_COMM_WORLD);
236  pcg->SetOperator(*A);
237  pcg->SetRelTol(1e-12);
238  pcg->SetMaxIter(2000);
239  pcg->SetPrintLevel(1);
240  if (hybridization) { prec = new HypreBoomerAMG(*A.As<HypreParMatrix>()); }
241  else if (pa) { prec = new OperatorJacobiSmoother(*a, ess_tdof_list); }
242  else
243  {
244  ParFiniteElementSpace *prec_fespace =
245  (a->StaticCondensationIsEnabled() ? a->SCParFESpace() : fespace);
246  if (dim == 2) { prec = new HypreAMS(*A.As<HypreParMatrix>(), prec_fespace); }
247  else { prec = new HypreADS(*A.As<HypreParMatrix>(), prec_fespace); }
248  }
249  pcg->SetPreconditioner(*prec);
250  pcg->Mult(B, X);
251 
252  // 14. Recover the parallel grid function corresponding to X. This is the
253  // local finite element solution on each processor.
254  a->RecoverFEMSolution(X, *b, x);
255 
256  // 15. Compute and print the L^2 norm of the error.
257  {
258  double error = x.ComputeL2Error(F);
259  if (myid == 0)
260  {
261  cout << "\n|| F_h - F ||_{L^2} = " << error << '\n' << endl;
262  }
263  }
264 
265  // 16. 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) << myid;
270  sol_name << "sol." << setfill('0') << setw(6) << myid;
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  // 17. 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 " << num_procs << " " << myid << "\n";
288  sol_sock.precision(8);
289  sol_sock << "solution\n" << *pmesh << x << flush;
290  }
291 
292  // 18. Free the used memory.
293  delete pcg;
294  delete prec;
295  delete hfes;
296  delete hfec;
297  delete a;
298  delete alpha;
299  delete beta;
300  delete b;
301  delete fespace;
302  delete fec;
303  delete pmesh;
304 
305  return 0;
306 }
307 
308 
309 // The exact solution (for non-surface meshes)
310 void F_exact(const Vector &p, Vector &F)
311 {
312  int dim = p.Size();
313 
314  double x = p(0);
315  double y = p(1);
316  // double z = (dim == 3) ? p(2) : 0.0; // Uncomment if F is changed to depend on z
317 
318  F(0) = cos(kappa*x)*sin(kappa*y);
319  F(1) = cos(kappa*y)*sin(kappa*x);
320  if (dim == 3)
321  {
322  F(2) = 0.0;
323  }
324 }
325 
326 // The right hand side
327 void f_exact(const Vector &p, Vector &f)
328 {
329  int dim = p.Size();
330 
331  double x = p(0);
332  double y = p(1);
333  // double z = (dim == 3) ? p(2) : 0.0; // Uncomment if f is changed to depend on z
334 
335  double temp = 1 + 2*kappa*kappa;
336 
337  f(0) = temp*cos(kappa*x)*sin(kappa*y);
338  f(1) = temp*cos(kappa*y)*sin(kappa*x);
339  if (dim == 3)
340  {
341  f(2) = 0;
342  }
343 }
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1032
Conjugate gradient method.
Definition: solvers.hpp:465
ParFiniteElementSpace * SCParFESpace() const
Return the parallel trace FE space associated with static condensation.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:104
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:1596
The Auxiliary-space Divergence Solver in hypre.
Definition: hypre.hpp:1635
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
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:711
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:285
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:928
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
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
double kappa
Definition: ex24.cpp:54
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:525
(Q div u, div v) for RT elements
void EnableHybridization(FiniteElementSpace *constr_space, BilinearFormIntegrator *constr_integ, const Array< int > &ess_tdof_list)
Enable hybridization.
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1443
Class for parallel linear form.
Definition: plinearform.hpp:26
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
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[]
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:274
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9737
constexpr int visport
void SetMaxIter(int max_it)
Definition: solvers.hpp:200
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
void F_exact(const Vector &, Vector &)
Definition: ex4.cpp:268
void Assemble(int skip_zeros=1)
Assemble the local matrix.
int Dimension() const
Definition: mesh.hpp:999
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
FDualNumber< tbase > cos(const FDualNumber< tbase > &f)
cos([dual number])
Definition: fdual.hpp:471
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:337
A general vector function coefficient.
int SpaceDimension() const
Definition: mesh.hpp:1000
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
bool StaticCondensationIsEnabled() const
Check if static condensation was actually enabled by a previous call to EnableStaticCondensation().
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:270
void SetRelTol(double rtol)
Definition: solvers.hpp:198
FDualNumber< tbase > sin(const FDualNumber< tbase > &f)
sin([dual number])
Definition: fdual.hpp:572
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
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
void f_exact(const Vector &, Vector &)
Definition: ex3.cpp:258
HYPRE_Int HYPRE_BigInt
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
FDualNumber< tbase > log(const FDualNumber< tbase > &f)
log([dual number])
Definition: fdual.hpp:515
double a
Definition: lissajous.cpp:41
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:281
int dim
Definition: ex24.cpp:53
double freq
Definition: ex24.cpp:54
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.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:266
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:479
const double alpha
Definition: ex15.cpp:369
Vector data type.
Definition: vector.hpp:60
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:173
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4684
Base class for solvers.
Definition: operator.hpp:651
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
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:337
Class for parallel meshes.
Definition: pmesh.hpp:32
int main()
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150