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 //
3 // Compile with: make ex3p
4 //
5 // Sample runs: mpirun -np 4 ex3p -m ../data/star.mesh
6 // mpirun -np 4 ex3p -m ../data/square-disc.mesh -o 2
7 // mpirun -np 4 ex3p -m ../data/beam-tet.mesh
8 // mpirun -np 4 ex3p -m ../data/beam-hex.mesh
9 // mpirun -np 4 ex3p -m ../data/beam-hex.mesh -o 2 -pa
10 // mpirun -np 4 ex3p -m ../data/escher.mesh
11 // mpirun -np 4 ex3p -m ../data/escher.mesh -o 2
12 // mpirun -np 4 ex3p -m ../data/fichera.mesh
13 // mpirun -np 4 ex3p -m ../data/fichera-q2.vtk
14 // mpirun -np 4 ex3p -m ../data/fichera-q3.mesh
15 // mpirun -np 4 ex3p -m ../data/square-disc-nurbs.mesh
16 // mpirun -np 4 ex3p -m ../data/beam-hex-nurbs.mesh
17 // mpirun -np 4 ex3p -m ../data/amr-quad.mesh -o 2
18 // mpirun -np 4 ex3p -m ../data/amr-hex.mesh
19 // mpirun -np 4 ex3p -m ../data/ref-prism.mesh -o 1
20 // mpirun -np 4 ex3p -m ../data/octahedron.mesh -o 1
21 // mpirun -np 4 ex3p -m ../data/star-surf.mesh -o 2
22 // mpirun -np 4 ex3p -m ../data/mobius-strip.mesh -o 2 -f 0.1
23 // mpirun -np 4 ex3p -m ../data/klein-bottle.mesh -o 2 -f 0.1
24 //
25 // Device sample runs:
26 // mpirun -np 4 ex3p -m ../data/star.mesh -pa -d cuda
27 // mpirun -np 4 ex3p -m ../data/star.mesh -no-pa -d cuda
28 // mpirun -np 4 ex3p -m ../data/star.mesh -pa -d raja-cuda
29 // mpirun -np 4 ex3p -m ../data/star.mesh -pa -d raja-omp
30 // mpirun -np 4 ex3p -m ../data/beam-hex.mesh -pa -d cuda
31 //
32 // Description: This example code solves a simple electromagnetic diffusion
33 // problem corresponding to the second order definite Maxwell
34 // equation curl curl E + E = f with boundary condition
35 // E x n = <given tangential field>. Here, we use a given exact
36 // solution E and compute the corresponding r.h.s. f.
37 // We discretize with Nedelec finite elements in 2D or 3D.
38 //
39 // The example demonstrates the use of H(curl) finite element
40 // spaces with the curl-curl and the (vector finite element) mass
41 // bilinear form, as well as the computation of discretization
42 // error when the exact solution is known. Static condensation is
43 // also illustrated.
44 //
45 // We recommend viewing examples 1-2 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, E, and r.h.s., f. See below for implementation.
55 void E_exact(const Vector &, Vector &);
56 void f_exact(const Vector &, Vector &);
57 double freq = 1.0, kappa;
58 int dim;
59 
60 int main(int argc, char *argv[])
61 {
62  // 1. Initialize MPI and HYPRE.
63  Mpi::Init(argc, argv);
64  int num_procs = Mpi::WorldSize();
65  int myid = Mpi::WorldRank();
66  Hypre::Init();
67 
68  // 2. Parse command-line options.
69  const char *mesh_file = "../data/beam-tet.mesh";
70  int order = 1;
71  bool static_cond = false;
72  bool pa = false;
73  const char *device_config = "cpu";
74  bool visualization = true;
75 #ifdef MFEM_USE_AMGX
76  bool useAmgX = false;
77 #endif
78 
79  OptionsParser args(argc, argv);
80  args.AddOption(&mesh_file, "-m", "--mesh",
81  "Mesh file to use.");
82  args.AddOption(&order, "-o", "--order",
83  "Finite element order (polynomial degree).");
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(&pa, "-pa", "--partial-assembly", "-no-pa",
89  "--no-partial-assembly", "Enable Partial Assembly.");
90  args.AddOption(&device_config, "-d", "--device",
91  "Device configuration string, see Device::Configure().");
92  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
93  "--no-visualization",
94  "Enable or disable GLVis visualization.");
95 #ifdef MFEM_USE_AMGX
96  args.AddOption(&useAmgX, "-amgx", "--useAmgX", "-no-amgx",
97  "--no-useAmgX",
98  "Enable or disable AmgX in MatrixFreeAMS.");
99 #endif
100 
101  args.Parse();
102  if (!args.Good())
103  {
104  if (myid == 0)
105  {
106  args.PrintUsage(cout);
107  }
108  return 1;
109  }
110  if (myid == 0)
111  {
112  args.PrintOptions(cout);
113  }
114  kappa = freq * M_PI;
115 
116  // 3. Enable hardware devices such as GPUs, and programming models such as
117  // CUDA, OCCA, RAJA and OpenMP based on command line options.
118  Device device(device_config);
119  if (myid == 0) { device.Print(); }
120 
121  // 4. Read the (serial) mesh from the given mesh file on all processors. We
122  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
123  // and volume meshes with the same code.
124  Mesh *mesh = new Mesh(mesh_file, 1, 1);
125  dim = mesh->Dimension();
126  int sdim = mesh->SpaceDimension();
127 
128  // 5. Refine the serial mesh on all processors to increase the resolution. In
129  // this example we do 'ref_levels' of uniform refinement. We choose
130  // 'ref_levels' to be the largest number that gives a final mesh with no
131  // more than 1,000 elements.
132  {
133  int ref_levels = (int)floor(log(1000./mesh->GetNE())/log(2.)/dim);
134  for (int l = 0; l < ref_levels; l++)
135  {
136  mesh->UniformRefinement();
137  }
138  }
139 
140  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
141  // this mesh further in parallel to increase the resolution. Once the
142  // parallel mesh is defined, the serial mesh can be deleted.
143  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
144  delete mesh;
145  {
146  int par_ref_levels = 2;
147  for (int l = 0; l < par_ref_levels; l++)
148  {
149  pmesh->UniformRefinement();
150  }
151  }
152 
153  // 7. Define a parallel finite element space on the parallel mesh. Here we
154  // use the Nedelec finite elements of the specified order.
155  FiniteElementCollection *fec = new ND_FECollection(order, dim);
156  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
157  HYPRE_BigInt size = fespace->GlobalTrueVSize();
158  if (myid == 0)
159  {
160  cout << "Number of finite element unknowns: " << size << endl;
161  }
162 
163  // 8. Determine the list of true (i.e. parallel conforming) essential
164  // boundary dofs. In this example, the boundary conditions are defined
165  // by marking all the boundary attributes from the mesh as essential
166  // (Dirichlet) and converting them to a list of true dofs.
167  Array<int> ess_tdof_list;
168  Array<int> ess_bdr;
169  if (pmesh->bdr_attributes.Size())
170  {
171  ess_bdr.SetSize(pmesh->bdr_attributes.Max());
172  ess_bdr = 1;
173  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
174  }
175 
176  // 9. Set up the parallel linear form b(.) which corresponds to the
177  // right-hand side of the FEM linear system, which in this case is
178  // (f,phi_i) where f is given by the function f_exact and phi_i are the
179  // basis functions in the finite element fespace.
181  ParLinearForm *b = new ParLinearForm(fespace);
183  b->Assemble();
184 
185  // 10. Define the solution vector x as a parallel finite element grid function
186  // corresponding to fespace. Initialize x by projecting the exact
187  // solution. Note that only values from the boundary edges will be used
188  // when eliminating the non-homogeneous boundary condition to modify the
189  // r.h.s. vector b.
190  ParGridFunction x(fespace);
192  x.ProjectCoefficient(E);
193 
194  // 11. Set up the parallel bilinear form corresponding to the EM diffusion
195  // operator curl muinv curl + sigma I, by adding the curl-curl and the
196  // mass domain integrators.
197  Coefficient *muinv = new ConstantCoefficient(1.0);
199  ParBilinearForm *a = new ParBilinearForm(fespace);
200  if (pa) { a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
201  a->AddDomainIntegrator(new CurlCurlIntegrator(*muinv));
203 
204  // 12. Assemble the parallel bilinear form and the corresponding linear
205  // system, applying any necessary transformations such as: parallel
206  // assembly, eliminating boundary conditions, applying conforming
207  // constraints for non-conforming AMR, static condensation, etc.
208  if (static_cond) { a->EnableStaticCondensation(); }
209  a->Assemble();
210 
211  OperatorPtr A;
212  Vector B, X;
213  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
214 
215  // 13. Solve the system AX=B using PCG with an AMS preconditioner.
216  if (pa)
217  {
218 #ifdef MFEM_USE_AMGX
219  MatrixFreeAMS ams(*a, *A, *fespace, muinv, sigma, NULL, ess_bdr, useAmgX);
220 #else
221  MatrixFreeAMS ams(*a, *A, *fespace, muinv, sigma, NULL, ess_bdr);
222 #endif
223  CGSolver cg(MPI_COMM_WORLD);
224  cg.SetRelTol(1e-12);
225  cg.SetMaxIter(1000);
226  cg.SetPrintLevel(1);
227  cg.SetOperator(*A);
228  cg.SetPreconditioner(ams);
229  cg.Mult(B, X);
230  }
231  else
232  {
233  if (myid == 0)
234  {
235  cout << "Size of linear system: "
236  << A.As<HypreParMatrix>()->GetGlobalNumRows() << endl;
237  }
238 
239  ParFiniteElementSpace *prec_fespace =
240  (a->StaticCondensationIsEnabled() ? a->SCParFESpace() : fespace);
241  HypreAMS ams(*A.As<HypreParMatrix>(), prec_fespace);
242  HyprePCG pcg(*A.As<HypreParMatrix>());
243  pcg.SetTol(1e-12);
244  pcg.SetMaxIter(500);
245  pcg.SetPrintLevel(2);
246  pcg.SetPreconditioner(ams);
247  pcg.Mult(B, X);
248  }
249 
250  // 14. Recover the parallel grid function corresponding to X. This is the
251  // local finite element solution on each processor.
252  a->RecoverFEMSolution(X, *b, x);
253 
254  // 15. Compute and print the L^2 norm of the error.
255  {
256  double error = x.ComputeL2Error(E);
257  if (myid == 0)
258  {
259  cout << "\n|| E_h - E ||_{L^2} = " << error << '\n' << endl;
260  }
261  }
262 
263  // 16. Save the refined mesh and the solution in parallel. This output can
264  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
265  {
266  ostringstream mesh_name, sol_name;
267  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
268  sol_name << "sol." << setfill('0') << setw(6) << myid;
269 
270  ofstream mesh_ofs(mesh_name.str().c_str());
271  mesh_ofs.precision(8);
272  pmesh->Print(mesh_ofs);
273 
274  ofstream sol_ofs(sol_name.str().c_str());
275  sol_ofs.precision(8);
276  x.Save(sol_ofs);
277  }
278 
279  // 17. Send the solution by socket to a GLVis server.
280  if (visualization)
281  {
282  char vishost[] = "localhost";
283  int visport = 19916;
284  socketstream sol_sock(vishost, visport);
285  sol_sock << "parallel " << num_procs << " " << myid << "\n";
286  sol_sock.precision(8);
287  sol_sock << "solution\n" << *pmesh << x << flush;
288  }
289 
290  // 18. Free the used memory.
291  delete a;
292  delete sigma;
293  delete muinv;
294  delete b;
295  delete fespace;
296  delete fec;
297  delete pmesh;
298 
299  return 0;
300 }
301 
302 
303 void E_exact(const Vector &x, Vector &E)
304 {
305  if (dim == 3)
306  {
307  E(0) = sin(kappa * x(1));
308  E(1) = sin(kappa * x(2));
309  E(2) = sin(kappa * x(0));
310  }
311  else
312  {
313  E(0) = sin(kappa * x(1));
314  E(1) = sin(kappa * x(0));
315  if (x.Size() == 3) { E(2) = 0.0; }
316  }
317 }
318 
319 void f_exact(const Vector &x, Vector &f)
320 {
321  if (dim == 3)
322  {
323  f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
324  f(1) = (1. + kappa * kappa) * sin(kappa * x(2));
325  f(2) = (1. + kappa * kappa) * sin(kappa * x(0));
326  }
327  else
328  {
329  f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
330  f(1) = (1. + kappa * kappa) * sin(kappa * x(0));
331  if (x.Size() == 3) { f(2) = 0.0; }
332  }
333 }
void SetTol(double tol)
Definition: hypre.cpp:3775
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
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
Integrator for (curl u, curl v) for Nedelec elements.
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
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3795
Class for parallel linear form.
Definition: plinearform.hpp:26
An auxiliary Maxwell solver for a high-order curl-curl system without high-order assembly.
Definition: auxiliary.hpp:168
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[]
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 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
void SetRelTol(double rtol)
Definition: solvers.hpp:198
FDualNumber< tbase > sin(const FDualNumber< tbase > &f)
sin([dual number])
Definition: fdual.hpp:572
PCG solver in hypre.
Definition: hypre.hpp:1116
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 SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
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 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 E_exact(const Vector &, Vector &)
Definition: ex3.cpp:242
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
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:411
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
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
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:3823
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)
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