MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex3p.cpp
Go to the documentation of this file.
1 // MFEM Example 3 - Parallel Version
2 // PETSc Modification
3 //
4 // Compile with: make ex3p
5 //
6 // Sample runs:
7 // mpirun -np 4 ex3p -m ../../data/klein-bottle.mesh -o 2 -f 0.1 --petscopts rc_ex3p
8 // mpirun -np 4 ex3p -m ../../data/klein-bottle.mesh -o 2 -f 0.1 --petscopts rc_ex3p_bddc --nonoverlapping
9 //
10 // Description: This example code solves a simple electromagnetic diffusion
11 // problem corresponding to the second order definite Maxwell
12 // equation curl curl E + E = f with boundary condition
13 // E x n = <given tangential field>. Here, we use a given exact
14 // solution E and compute the corresponding r.h.s. f.
15 // We discretize with Nedelec finite elements in 2D or 3D.
16 //
17 // The example demonstrates the use of H(curl) finite element
18 // spaces with the curl-curl and the (vector finite element) mass
19 // bilinear form, as well as the computation of discretization
20 // error when the exact solution is known. Static condensation is
21 // also illustrated.
22 //
23 // The example also show how to use the non-overlapping feature of
24 // the ParBilinearForm class to obtain the linear operator in
25 // a format suitable for the BDDC preconditioner in PETSc.
26 //
27 // We recommend viewing examples 1-2 before viewing this example.
28 
29 #include "mfem.hpp"
30 #include <fstream>
31 #include <iostream>
32 
33 #ifndef MFEM_USE_PETSC
34 #error This example requires that MFEM is built with MFEM_USE_PETSC=YES
35 #endif
36 
37 using namespace std;
38 using namespace mfem;
39 
40 // Exact solution, E, and r.h.s., f. See below for implementation.
41 void E_exact(const Vector &, Vector &);
42 void f_exact(const Vector &, Vector &);
43 double freq = 1.0, kappa;
44 int dim;
45 
46 int main(int argc, char *argv[])
47 {
48  // 1. Initialize MPI and HYPRE.
49  Mpi::Init(argc, argv);
50  int num_procs = Mpi::WorldSize();
51  int myid = Mpi::WorldRank();
52  Hypre::Init();
53 
54  // 2. Parse command-line options.
55  const char *mesh_file = "../../data/beam-tet.mesh";
56  int order = 1;
57  bool static_cond = false;
58  bool visualization = 1;
59  bool use_petsc = true;
60  const char *petscrc_file = "";
61  bool use_nonoverlapping = false;
62 
63  OptionsParser args(argc, argv);
64  args.AddOption(&mesh_file, "-m", "--mesh",
65  "Mesh file to use.");
66  args.AddOption(&order, "-o", "--order",
67  "Finite element order (polynomial degree).");
68  args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
69  " solution.");
70  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
71  "--no-static-condensation", "Enable static condensation.");
72  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
73  "--no-visualization",
74  "Enable or disable GLVis visualization.");
75  args.AddOption(&use_petsc, "-usepetsc", "--usepetsc", "-no-petsc",
76  "--no-petsc",
77  "Use or not PETSc to solve the linear system.");
78  args.AddOption(&petscrc_file, "-petscopts", "--petscopts",
79  "PetscOptions file to use.");
80  args.AddOption(&use_nonoverlapping, "-nonoverlapping", "--nonoverlapping",
81  "-no-nonoverlapping", "--no-nonoverlapping",
82  "Use or not the block diagonal PETSc's matrix format "
83  "for non-overlapping domain decomposition.");
84  args.Parse();
85  if (!args.Good())
86  {
87  if (myid == 0)
88  {
89  args.PrintUsage(cout);
90  }
91  return 1;
92  }
93  if (myid == 0)
94  {
95  args.PrintOptions(cout);
96  }
97  // 2b. We initialize PETSc
98  if (use_petsc) { MFEMInitializePetsc(NULL,NULL,petscrc_file,NULL); }
99  kappa = freq * M_PI;
100 
101  // 3. Read the (serial) mesh from the given mesh file on all processors. We
102  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
103  // and volume meshes with the same code.
104  Mesh *mesh = new Mesh(mesh_file, 1, 1);
105  dim = mesh->Dimension();
106  int sdim = mesh->SpaceDimension();
107 
108  // 4. Refine the serial mesh on all processors to increase the resolution. In
109  // this example we do 'ref_levels' of uniform refinement. We choose
110  // 'ref_levels' to be the largest number that gives a final mesh with no
111  // more than 1,000 elements.
112  {
113  int ref_levels =
114  (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
115  for (int l = 0; l < ref_levels; l++)
116  {
117  mesh->UniformRefinement();
118  }
119  }
120 
121  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
122  // this mesh further in parallel to increase the resolution. Once the
123  // parallel mesh is defined, the serial mesh can be deleted.
124  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
125  delete mesh;
126  {
127  int par_ref_levels = 2;
128  for (int l = 0; l < par_ref_levels; l++)
129  {
130  pmesh->UniformRefinement();
131  }
132  }
133 
134  // 6. Define a parallel finite element space on the parallel mesh. Here we
135  // use the Nedelec finite elements of the specified order.
136  FiniteElementCollection *fec = new ND_FECollection(order, dim);
137  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
138  HYPRE_BigInt size = fespace->GlobalTrueVSize();
139  if (myid == 0)
140  {
141  cout << "Number of finite element unknowns: " << size << endl;
142  }
143 
144  // 7. Determine the list of true (i.e. parallel conforming) essential
145  // boundary dofs. In this example, the boundary conditions are defined
146  // by marking all the boundary attributes from the mesh as essential
147  // (Dirichlet) and converting them to a list of true dofs.
148  Array<int> ess_tdof_list;
149  if (pmesh->bdr_attributes.Size())
150  {
151  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
152  ess_bdr = 1;
153  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
154  }
155 
156  // 8. Set up the parallel linear form b(.) which corresponds to the
157  // right-hand side of the FEM linear system, which in this case is
158  // (f,phi_i) where f is given by the function f_exact and phi_i are the
159  // basis functions in the finite element fespace.
161  ParLinearForm *b = new ParLinearForm(fespace);
163  b->Assemble();
164 
165  // 9. Define the solution vector x as a parallel finite element grid function
166  // corresponding to fespace. Initialize x by projecting the exact
167  // solution. Note that only values from the boundary edges will be used
168  // when eliminating the non-homogeneous boundary condition to modify the
169  // r.h.s. vector b.
170  ParGridFunction x(fespace);
172  x.ProjectCoefficient(E);
173 
174  // 10. Set up the parallel bilinear form corresponding to the EM diffusion
175  // operator curl muinv curl + sigma I, by adding the curl-curl and the
176  // mass domain integrators.
177  Coefficient *muinv = new ConstantCoefficient(1.0);
179  ParBilinearForm *a = new ParBilinearForm(fespace);
180  a->AddDomainIntegrator(new CurlCurlIntegrator(*muinv));
182 
183  // 11. Assemble the parallel bilinear form and the corresponding linear
184  // system, applying any necessary transformations such as: parallel
185  // assembly, eliminating boundary conditions, applying conforming
186  // constraints for non-conforming AMR, static condensation, etc.
187  if (static_cond) { a->EnableStaticCondensation(); }
188  a->Assemble();
189 
190  Vector B, X;
191  if (!use_petsc)
192  {
193  HypreParMatrix A;
194  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
195 
196  if (myid == 0)
197  {
198  cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
199  }
200 
201  // 12. Define and apply a parallel PCG solver for AX=B with the AMS
202  // preconditioner from hypre.
203  ParFiniteElementSpace *prec_fespace =
204  (a->StaticCondensationIsEnabled() ? a->SCParFESpace() : fespace);
205  HypreSolver *ams = new HypreAMS(A, prec_fespace);
206  HyprePCG *pcg = new HyprePCG(A);
207  pcg->SetTol(1e-10);
208  pcg->SetMaxIter(500);
209  pcg->SetPrintLevel(2);
210  pcg->SetPreconditioner(*ams);
211  pcg->Mult(B, X);
212  delete pcg;
213  delete ams;
214  }
215  else
216  {
217  PetscParMatrix A;
218  a->SetOperatorType(use_nonoverlapping ?
219  Operator::PETSC_MATIS : Operator::PETSC_MATAIJ);
220  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
221 
222  if (myid == 0)
223  {
224  cout << "Size of linear system: " << A.M() << endl;
225  }
226 
227  // 12. Define and apply a parallel PCG solver.
228  ParFiniteElementSpace *prec_fespace =
229  (a->StaticCondensationIsEnabled() ? a->SCParFESpace() : fespace);
230  PetscPCGSolver *pcg = new PetscPCGSolver(A);
231  PetscPreconditioner *prec = NULL;
232  pcg->SetTol(1e-10);
233  pcg->SetMaxIter(500);
234  pcg->SetPrintLevel(2);
235  if (use_nonoverlapping)
236  {
237  // Auxiliary class for BDDC customization
239  // Inform the solver about the finite element space
240  opts.SetSpace(prec_fespace);
241  // Inform the solver about essential dofs
242  opts.SetEssBdrDofs(&ess_tdof_list);
243  // Create a BDDC solver with parameters
244  prec = new PetscBDDCSolver(A,opts);
245  }
246  else
247  {
248  // Create an empty preconditioner object that can
249  // be customized at runtime
250  prec = new PetscPreconditioner(A,"solver_");
251  }
252  pcg->SetPreconditioner(*prec);
253  pcg->Mult(B, X);
254  delete pcg;
255  delete prec;
256  }
257 
258  // 13. Recover the parallel grid function corresponding to X. This is the
259  // local finite element solution on each processor.
260  a->RecoverFEMSolution(X, *b, x);
261 
262  // 14. Compute and print the L^2 norm of the error.
263  {
264  double err = x.ComputeL2Error(E);
265  if (myid == 0)
266  {
267  cout << "\n|| E_h - E ||_{L^2} = " << err << '\n' << endl;
268  }
269  }
270 
271  // 15. Save the refined mesh and the solution in parallel. This output can
272  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
273  {
274  ostringstream mesh_name, sol_name;
275  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
276  sol_name << "sol." << setfill('0') << setw(6) << myid;
277 
278  ofstream mesh_ofs(mesh_name.str().c_str());
279  mesh_ofs.precision(8);
280  pmesh->Print(mesh_ofs);
281 
282  ofstream sol_ofs(sol_name.str().c_str());
283  sol_ofs.precision(8);
284  x.Save(sol_ofs);
285  }
286 
287  // 16. Send the solution by socket to a GLVis server.
288  if (visualization)
289  {
290  char vishost[] = "localhost";
291  int visport = 19916;
292  socketstream sol_sock(vishost, visport);
293  sol_sock << "parallel " << num_procs << " " << myid << "\n";
294  sol_sock.precision(8);
295  sol_sock << "solution\n" << *pmesh << x << flush;
296  }
297 
298  // 17. Free the used memory.
299  delete a;
300  delete sigma;
301  delete muinv;
302  delete b;
303  delete fespace;
304  delete fec;
305  delete pmesh;
306 
307  // We finalize PETSc
308  if (use_petsc) { MFEMFinalizePetsc(); }
309 
310  return 0;
311 }
312 
313 
314 void E_exact(const Vector &x, Vector &E)
315 {
316  if (dim == 3)
317  {
318  E(0) = sin(kappa * x(1));
319  E(1) = sin(kappa * x(2));
320  E(2) = sin(kappa * x(0));
321  }
322  else
323  {
324  E(0) = sin(kappa * x(1));
325  E(1) = sin(kappa * x(0));
326  if (x.Size() == 3) { E(2) = 0.0; }
327  }
328 }
329 
330 void f_exact(const Vector &x, Vector &f)
331 {
332  if (dim == 3)
333  {
334  f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
335  f(1) = (1. + kappa * kappa) * sin(kappa * x(2));
336  f(2) = (1. + kappa * kappa) * sin(kappa * x(0));
337  }
338  else
339  {
340  f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
341  f(1) = (1. + kappa * kappa) * sin(kappa * x(0));
342  if (x.Size() == 3) { f(2) = 0.0; }
343  }
344 }
void SetTol(double tol)
Definition: hypre.cpp:3775
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
void SetPreconditioner(Solver &precond)
Definition: petsc.cpp:3021
void SetPrintLevel(int plev)
Definition: petsc.cpp:2370
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1032
ParFiniteElementSpace * SCParFESpace() const
Return the parallel trace FE space associated with static condensation.
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:1596
Abstract class for PETSc&#39;s preconditioners.
Definition: petsc.hpp:775
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(...
Wrapper for PETSc&#39;s matrix class.
Definition: petsc.hpp:307
Integrator for (curl u, curl v) for Nedelec elements.
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
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
void SetMaxIter(int max_iter)
Definition: petsc.cpp:2343
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3795
PetscInt M() const
Returns the global number of rows.
Definition: petsc.cpp:924
void SetOperatorType(Operator::Type tid)
Set the operator type id for the parallel matrix/operator when using AssemblyLevel::LEGACY.
Class for parallel linear form.
Definition: plinearform.hpp:26
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[]
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9737
constexpr int visport
Auxiliary class for BDDC customization.
Definition: petsc.hpp:801
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.
int Dimension() const
Definition: mesh.hpp:999
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3785
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
FDualNumber< tbase > sin(const FDualNumber< tbase > &f)
sin([dual number])
Definition: fdual.hpp:572
PCG solver in hypre.
Definition: hypre.hpp:1116
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
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
void MFEMFinalizePetsc()
Definition: petsc.cpp:192
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:281
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:629
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
Definition: petsc.cpp:168
void SetSpace(ParFiniteElementSpace *fe)
Definition: petsc.hpp:816
int dim
Definition: ex24.cpp:53
double freq
Definition: ex24.cpp:54
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:3800
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
void SetTol(double tol)
Definition: petsc.cpp:2288
void E_exact(const Vector &, Vector &)
Definition: ex3.cpp:242
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:1026
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:266
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:411
Vector data type.
Definition: vector.hpp:60
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4684
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:337
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:3823
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
Definition: petsc.cpp:3101
Class for parallel meshes.
Definition: pmesh.hpp:32
void SetEssBdrDofs(const Array< int > *essdofs, bool loc=false)
Specify dofs on the essential boundary.
Definition: petsc.hpp:821
int main()
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
double f(const Vector &p)
double sigma(const Vector &x)
Definition: maxwell.cpp:102
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150