MFEM  v4.3.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/star-surf.mesh -o 2
20 // mpirun -np 4 ex3p -m ../data/mobius-strip.mesh -o 2 -f 0.1
21 // mpirun -np 4 ex3p -m ../data/klein-bottle.mesh -o 2 -f 0.1
22 //
23 // Device sample runs:
24 // mpirun -np 4 ex3p -m ../data/star.mesh -pa -d cuda
25 // mpirun -np 4 ex3p -m ../data/star.mesh -no-pa -d cuda
26 // mpirun -np 4 ex3p -m ../data/star.mesh -pa -d raja-cuda
27 // mpirun -np 4 ex3p -m ../data/star.mesh -pa -d raja-omp
28 // mpirun -np 4 ex3p -m ../data/beam-hex.mesh -pa -d cuda
29 //
30 // Description: This example code solves a simple electromagnetic diffusion
31 // problem corresponding to the second order definite Maxwell
32 // equation curl curl E + E = f with boundary condition
33 // E x n = <given tangential field>. Here, we use a given exact
34 // solution E and compute the corresponding r.h.s. f.
35 // We discretize with Nedelec finite elements in 2D or 3D.
36 //
37 // The example demonstrates the use of H(curl) finite element
38 // spaces with the curl-curl and the (vector finite element) mass
39 // bilinear form, as well as the computation of discretization
40 // error when the exact solution is known. Static condensation is
41 // also illustrated.
42 //
43 // We recommend viewing examples 1-2 before viewing this example.
44 
45 #include "mfem.hpp"
46 #include <fstream>
47 #include <iostream>
48 
49 using namespace std;
50 using namespace mfem;
51 
52 // Exact solution, E, and r.h.s., f. See below for implementation.
53 void E_exact(const Vector &, Vector &);
54 void f_exact(const Vector &, Vector &);
55 double freq = 1.0, kappa;
56 int dim;
57 
58 int main(int argc, char *argv[])
59 {
60  // 1. Initialize MPI.
61  int num_procs, myid;
62  MPI_Init(&argc, &argv);
63  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
64  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
65 
66  // 2. Parse command-line options.
67  const char *mesh_file = "../data/beam-tet.mesh";
68  int order = 1;
69  bool static_cond = false;
70  bool pa = false;
71  const char *device_config = "cpu";
72  bool visualization = true;
73 #ifdef MFEM_USE_AMGX
74  bool useAmgX = false;
75 #endif
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(&freq, "-f", "--frequency", "Set the frequency for the exact"
83  " solution.");
84  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
85  "--no-static-condensation", "Enable static condensation.");
86  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
87  "--no-partial-assembly", "Enable Partial Assembly.");
88  args.AddOption(&device_config, "-d", "--device",
89  "Device configuration string, see Device::Configure().");
90  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
91  "--no-visualization",
92  "Enable or disable GLVis visualization.");
93 #ifdef MFEM_USE_AMGX
94  args.AddOption(&useAmgX, "-amgx", "--useAmgX", "-no-amgx",
95  "--no-useAmgX",
96  "Enable or disable AmgX in MatrixFreeAMS.");
97 #endif
98 
99  args.Parse();
100  if (!args.Good())
101  {
102  if (myid == 0)
103  {
104  args.PrintUsage(cout);
105  }
106  // HYPRE_Finalize();
107  MPI_Finalize();
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. Tetrahedral
143  // meshes need to be reoriented before we can define high-order Nedelec
144  // spaces on them.
145  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
146  delete mesh;
147  {
148  int par_ref_levels = 2;
149  for (int l = 0; l < par_ref_levels; l++)
150  {
151  pmesh->UniformRefinement();
152  }
153  }
154  pmesh->ReorientTetMesh();
155 
156  // 7. Define a parallel finite element space on the parallel mesh. Here we
157  // use the Nedelec finite elements of the specified order.
158  FiniteElementCollection *fec = new ND_FECollection(order, dim);
159  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
160  HYPRE_BigInt size = fespace->GlobalTrueVSize();
161  if (myid == 0)
162  {
163  cout << "Number of finite element unknowns: " << size << endl;
164  }
165 
166  // 8. Determine the list of true (i.e. parallel conforming) essential
167  // boundary dofs. In this example, the boundary conditions are defined
168  // by marking all the boundary attributes from the mesh as essential
169  // (Dirichlet) and converting them to a list of true dofs.
170  Array<int> ess_tdof_list;
171  Array<int> ess_bdr;
172  if (pmesh->bdr_attributes.Size())
173  {
174  ess_bdr.SetSize(pmesh->bdr_attributes.Max());
175  ess_bdr = 1;
176  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
177  }
178 
179  // 9. Set up the parallel linear form b(.) which corresponds to the
180  // right-hand side of the FEM linear system, which in this case is
181  // (f,phi_i) where f is given by the function f_exact and phi_i are the
182  // basis functions in the finite element fespace.
184  ParLinearForm *b = new ParLinearForm(fespace);
186  b->Assemble();
187 
188  // 10. Define the solution vector x as a parallel finite element grid function
189  // corresponding to fespace. Initialize x by projecting the exact
190  // solution. Note that only values from the boundary edges will be used
191  // when eliminating the non-homogeneous boundary condition to modify the
192  // r.h.s. vector b.
193  ParGridFunction x(fespace);
195  x.ProjectCoefficient(E);
196 
197  // 11. Set up the parallel bilinear form corresponding to the EM diffusion
198  // operator curl muinv curl + sigma I, by adding the curl-curl and the
199  // mass domain integrators.
200  Coefficient *muinv = new ConstantCoefficient(1.0);
202  ParBilinearForm *a = new ParBilinearForm(fespace);
203  if (pa) { a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
204  a->AddDomainIntegrator(new CurlCurlIntegrator(*muinv));
206 
207  // 12. Assemble the parallel bilinear form and the corresponding linear
208  // system, applying any necessary transformations such as: parallel
209  // assembly, eliminating boundary conditions, applying conforming
210  // constraints for non-conforming AMR, static condensation, etc.
211  if (static_cond) { a->EnableStaticCondensation(); }
212  a->Assemble();
213 
214  OperatorPtr A;
215  Vector B, X;
216  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
217 
218  // 13. Solve the system AX=B using PCG with an AMS preconditioner.
219  if (pa)
220  {
221 #ifdef MFEM_USE_AMGX
222  MatrixFreeAMS ams(*a, *A, *fespace, muinv, sigma, NULL, ess_bdr, useAmgX);
223 #else
224  MatrixFreeAMS ams(*a, *A, *fespace, muinv, sigma, NULL, ess_bdr);
225 #endif
226  CGSolver cg(MPI_COMM_WORLD);
227  cg.SetRelTol(1e-12);
228  cg.SetMaxIter(1000);
229  cg.SetPrintLevel(1);
230  cg.SetOperator(*A);
231  cg.SetPreconditioner(ams);
232  cg.Mult(B, X);
233  }
234  else
235  {
236  if (myid == 0)
237  {
238  cout << "Size of linear system: "
239  << A.As<HypreParMatrix>()->GetGlobalNumRows() << endl;
240  }
241 
242  ParFiniteElementSpace *prec_fespace =
243  (a->StaticCondensationIsEnabled() ? a->SCParFESpace() : fespace);
244  HypreAMS ams(*A.As<HypreParMatrix>(), prec_fespace);
245  HyprePCG pcg(*A.As<HypreParMatrix>());
246  pcg.SetTol(1e-12);
247  pcg.SetMaxIter(500);
248  pcg.SetPrintLevel(2);
249  pcg.SetPreconditioner(ams);
250  pcg.Mult(B, X);
251  }
252 
253  // 14. Recover the parallel grid function corresponding to X. This is the
254  // local finite element solution on each processor.
255  a->RecoverFEMSolution(X, *b, x);
256 
257  // 15. Compute and print the L^2 norm of the error.
258  {
259  double err = x.ComputeL2Error(E);
260  if (myid == 0)
261  {
262  cout << "\n|| E_h - E ||_{L^2} = " << err << '\n' << endl;
263  }
264  }
265 
266  // 16. Save the refined mesh and the solution in parallel. This output can
267  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
268  {
269  ostringstream mesh_name, sol_name;
270  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
271  sol_name << "sol." << setfill('0') << setw(6) << myid;
272 
273  ofstream mesh_ofs(mesh_name.str().c_str());
274  mesh_ofs.precision(8);
275  pmesh->Print(mesh_ofs);
276 
277  ofstream sol_ofs(sol_name.str().c_str());
278  sol_ofs.precision(8);
279  x.Save(sol_ofs);
280  }
281 
282  // 17. Send the solution by socket to a GLVis server.
283  if (visualization)
284  {
285  char vishost[] = "localhost";
286  int visport = 19916;
287  socketstream sol_sock(vishost, visport);
288  sol_sock << "parallel " << num_procs << " " << myid << "\n";
289  sol_sock.precision(8);
290  sol_sock << "solution\n" << *pmesh << x << flush;
291  }
292 
293  // 18. Free the used memory.
294  delete a;
295  delete sigma;
296  delete muinv;
297  delete b;
298  delete fespace;
299  delete fec;
300  delete pmesh;
301 
302  MPI_Finalize();
303 
304  return 0;
305 }
306 
307 
308 void E_exact(const Vector &x, Vector &E)
309 {
310  if (dim == 3)
311  {
312  E(0) = sin(kappa * x(1));
313  E(1) = sin(kappa * x(2));
314  E(2) = sin(kappa * x(0));
315  }
316  else
317  {
318  E(0) = sin(kappa * x(1));
319  E(1) = sin(kappa * x(0));
320  if (x.Size() == 3) { E(2) = 0.0; }
321  }
322 }
323 
324 void f_exact(const Vector &x, Vector &f)
325 {
326  if (dim == 3)
327  {
328  f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
329  f(1) = (1. + kappa * kappa) * sin(kappa * x(2));
330  f(2) = (1. + kappa * kappa) * sin(kappa * x(0));
331  }
332  else
333  {
334  f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
335  f(1) = (1. + kappa * kappa) * sin(kappa * x(0));
336  if (x.Size() == 3) { f(2) = 0.0; }
337  }
338 }
void SetTol(double tol)
Definition: hypre.cpp:3569
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:790
Conjugate gradient method.
Definition: solvers.hpp:316
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:99
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:1540
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:616
virtual void ReorientTetMesh()
See the remarks for the serial version in mesh.hpp.
Definition: pmesh.cpp:2798
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:275
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:841
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:493
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3589
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
void SetPrintLevel(int print_lvl)
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:150
constexpr char vishost[]
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
constexpr int visport
void SetMaxIter(int max_it)
Definition: solvers.hpp:100
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.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4382
int Dimension() const
Definition: mesh.hpp:911
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3579
A general vector function coefficient.
int SpaceDimension() const
Definition: mesh.hpp:912
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:204
void SetRelTol(double rtol)
Definition: solvers.hpp:98
PCG solver in hypre.
Definition: hypre.hpp:1060
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 SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:674
void f_exact(const Vector &, Vector &)
Definition: ex3.cpp:257
HYPRE_Int HYPRE_BigInt
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
double a
Definition: lissajous.cpp:41
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:273
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:3594
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
void E_exact(const Vector &, Vector &)
Definition: ex3.cpp:241
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:330
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:416
Vector data type.
Definition: vector.hpp:60
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:92
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:277
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:3617
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