MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex24.cpp
Go to the documentation of this file.
1 // MFEM Example 24
2 //
3 // Compile with: make ex24
4 //
5 // Sample runs: ex24 -m ../data/star.mesh
6 // ex24 -m ../data/square-disc.mesh -o 2
7 // ex24 -m ../data/beam-tet.mesh
8 // ex24 -m ../data/beam-hex.mesh -o 2 -pa
9 // ex24 -m ../data/beam-hex.mesh -o 2 -pa -p 1
10 // ex24 -m ../data/beam-hex.mesh -o 2 -pa -p 2
11 // ex24 -m ../data/escher.mesh
12 // ex24 -m ../data/escher.mesh -o 2
13 // ex24 -m ../data/fichera.mesh
14 // ex24 -m ../data/fichera-q2.vtk
15 // ex24 -m ../data/fichera-q3.mesh
16 // ex24 -m ../data/square-disc-nurbs.mesh
17 // ex24 -m ../data/beam-hex-nurbs.mesh
18 // ex24 -m ../data/amr-quad.mesh -o 2
19 // ex24 -m ../data/amr-hex.mesh
20 //
21 // Device sample runs:
22 // ex24 -m ../data/star.mesh -pa -d cuda
23 // ex24 -m ../data/star.mesh -pa -d raja-cuda
24 // ex24 -m ../data/star.mesh -pa -d raja-omp
25 // ex24 -m ../data/beam-hex.mesh -pa -d cuda
26 //
27 // Description: This example code illustrates usage of mixed finite element
28 // spaces, with three variants:
29 //
30 // 1) (grad p, u) for p in H^1 tested against u in H(curl)
31 // 2) (curl v, u) for v in H(curl) tested against u in H(div), 3D
32 // 3) (div v, q) for v in H(div) tested against q in L_2
33 //
34 // Using different approaches, we project the gradient, curl, or
35 // divergence to the appropriate space.
36 //
37 // We recommend viewing examples 1, 3, and 5 before viewing this
38 // example.
39 
40 #include "mfem.hpp"
41 #include <fstream>
42 #include <iostream>
43 
44 using namespace std;
45 using namespace mfem;
46 
47 double p_exact(const Vector &x);
48 void gradp_exact(const Vector &, Vector &);
49 double div_gradp_exact(const Vector &x);
50 void v_exact(const Vector &x, Vector &v);
51 void curlv_exact(const Vector &x, Vector &cv);
52 
53 int dim;
54 double freq = 1.0, kappa;
55 
56 int main(int argc, char *argv[])
57 {
58  // 1. Parse command-line options.
59  const char *mesh_file = "../data/beam-hex.mesh";
60  int order = 1;
61  int prob = 0;
62  bool static_cond = false;
63  bool pa = false;
64  const char *device_config = "cpu";
65  bool visualization = 1;
66 
67  OptionsParser args(argc, argv);
68  args.AddOption(&mesh_file, "-m", "--mesh",
69  "Mesh file to use.");
70  args.AddOption(&order, "-o", "--order",
71  "Finite element order (polynomial degree).");
72  args.AddOption(&prob, "-p", "--problem-type",
73  "Choose between 0: grad, 1: curl, 2: div");
74  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
75  "--no-static-condensation", "Enable static condensation.");
76  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
77  "--no-partial-assembly", "Enable Partial Assembly.");
78  args.AddOption(&device_config, "-d", "--device",
79  "Device configuration string, see Device::Configure().");
80  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
81  "--no-visualization",
82  "Enable or disable GLVis visualization.");
83 
84  args.Parse();
85  if (!args.Good())
86  {
87  args.PrintUsage(cout);
88  return 1;
89  }
90  args.PrintOptions(cout);
91  kappa = freq * M_PI;
92 
93  // 2. Enable hardware devices such as GPUs, and programming models such as
94  // CUDA, OCCA, RAJA and OpenMP based on command line options.
95  Device device(device_config);
96  device.Print();
97 
98  // 3. Read the mesh from the given mesh file. We can handle triangular,
99  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
100  // the same code.
101  Mesh *mesh = new Mesh(mesh_file, 1, 1);
102  dim = mesh->Dimension();
103  int sdim = mesh->SpaceDimension();
104 
105  // 4. Refine the mesh to increase the resolution. In this example we do
106  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
107  // largest number that gives a final mesh with no more than 50,000
108  // elements.
109  {
110  int ref_levels = (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
111  for (int l = 0; l < ref_levels; l++)
112  {
113  mesh->UniformRefinement();
114  }
115  }
116  mesh->ReorientTetMesh();
117 
118  // 5. Define a finite element space on the mesh. Here we use Nedelec or
119  // Raviart-Thomas finite elements of the specified order.
120  FiniteElementCollection *trial_fec = NULL;
121  FiniteElementCollection *test_fec = NULL;
122 
123  if (prob == 0)
124  {
125  trial_fec = new H1_FECollection(order, dim);
126  test_fec = new ND_FECollection(order, dim);
127  }
128  else if (prob == 1)
129  {
130  trial_fec = new ND_FECollection(order, dim);
131  test_fec = new RT_FECollection(order-1, dim);
132  }
133  else
134  {
135  trial_fec = new RT_FECollection(order-1, dim);
136  test_fec = new L2_FECollection(order-1, dim);
137  }
138 
139  FiniteElementSpace trial_fes(mesh, trial_fec);
140  FiniteElementSpace test_fes(mesh, test_fec);
141 
142  int trial_size = trial_fes.GetTrueVSize();
143  int test_size = test_fes.GetTrueVSize();
144 
145  if (prob == 0)
146  {
147  cout << "Number of Nedelec finite element unknowns: " << test_size << endl;
148  cout << "Number of H1 finite element unknowns: " << trial_size << endl;
149  }
150  else if (prob == 1)
151  {
152  cout << "Number of Nedelec finite element unknowns: " << trial_size << endl;
153  cout << "Number of Raviart-Thomas finite element unknowns: " << test_size <<
154  endl;
155  }
156  else
157  {
158  cout << "Number of Raviart-Thomas finite element unknowns: "
159  << trial_size << endl;
160  cout << "Number of L2 finite element unknowns: " << test_size << endl;
161  }
162 
163  // 6. Define the solution vector as a finite element grid function
164  // corresponding to the trial fespace.
165  GridFunction gftest(&test_fes);
166  GridFunction gftrial(&trial_fes);
167  GridFunction x(&test_fes);
169  VectorFunctionCoefficient gradp_coef(sdim, gradp_exact);
170  VectorFunctionCoefficient v_coef(sdim, v_exact);
171  VectorFunctionCoefficient curlv_coef(sdim, curlv_exact);
172  FunctionCoefficient divgradp_coef(div_gradp_exact);
173 
174  if (prob == 0)
175  {
176  gftrial.ProjectCoefficient(p_coef);
177  }
178  else if (prob == 1)
179  {
180  gftrial.ProjectCoefficient(v_coef);
181  }
182  else
183  {
184  gftrial.ProjectCoefficient(gradp_coef);
185  }
186 
187  gftrial.SetTrueVector();
188  gftrial.SetFromTrueVector();
189 
190  // 7. Set up the bilinear forms for L2 projection.
191  ConstantCoefficient one(1.0);
192  BilinearForm a(&test_fes);
193  MixedBilinearForm a_mixed(&trial_fes, &test_fes);
194  if (pa)
195  {
196  a.SetAssemblyLevel(AssemblyLevel::PARTIAL);
197  a_mixed.SetAssemblyLevel(AssemblyLevel::PARTIAL);
198  }
199 
200  if (prob == 0)
201  {
204  }
205  else if (prob == 1)
206  {
209  }
210  else
211  {
214  }
215 
216  // 8. Assemble the bilinear form and the corresponding linear system,
217  // applying any necessary transformations such as: eliminating boundary
218  // conditions, applying conforming constraints for non-conforming AMR,
219  // static condensation, etc.
220  if (static_cond) { a.EnableStaticCondensation(); }
221 
222  a.Assemble();
223  if (!pa) { a.Finalize(); }
224 
225  a_mixed.Assemble();
226  if (!pa) { a_mixed.Finalize(); }
227 
228  if (pa)
229  {
230  a_mixed.Mult(gftrial, x);
231  }
232  else
233  {
234  SparseMatrix& mixed = a_mixed.SpMat();
235  mixed.Mult(gftrial, x);
236  }
237 
238  // 9. Define and apply a PCG solver for Ax = b with Jacobi preconditioner.
239  {
240  GridFunction rhs(&test_fes);
241  rhs = x;
242  x = 0.0;
243 
244  CGSolver cg;
245  cg.SetRelTol(1e-12);
246  cg.SetMaxIter(1000);
247  cg.SetPrintLevel(1);
248  if (pa)
249  {
250  Array<int> ess_tdof_list; // empty
251  OperatorJacobiSmoother Jacobi(a, ess_tdof_list);
252 
253  cg.SetOperator(a);
254  cg.SetPreconditioner(Jacobi);
255  cg.Mult(rhs, x);
256  }
257  else
258  {
259  SparseMatrix& Amat = a.SpMat();
260  DSmoother Jacobi(Amat);
261 
262  cg.SetOperator(Amat);
263  cg.SetPreconditioner(Jacobi);
264  cg.Mult(rhs, x);
265  }
266  }
267 
268  // 10. Compute the same field by applying a DiscreteInterpolator.
269  GridFunction discreteInterpolant(&test_fes);
270  DiscreteLinearOperator dlo(&trial_fes, &test_fes);
271  if (prob == 0)
272  {
274  }
275  else if (prob == 1)
276  {
278  }
279  else
280  {
282  }
283 
284  dlo.Assemble();
285  dlo.Mult(gftrial, discreteInterpolant);
286 
287  // 11. Compute the projection of the exact field.
288  GridFunction exact_proj(&test_fes);
289  if (prob == 0)
290  {
291  exact_proj.ProjectCoefficient(gradp_coef);
292  }
293  else if (prob == 1)
294  {
295  exact_proj.ProjectCoefficient(curlv_coef);
296  }
297  else
298  {
299  exact_proj.ProjectCoefficient(divgradp_coef);
300  }
301 
302  exact_proj.SetTrueVector();
303  exact_proj.SetFromTrueVector();
304 
305  // 12. Compute and print the L_2 norm of the error.
306  if (prob == 0)
307  {
308  double errSol = x.ComputeL2Error(gradp_coef);
309  double errInterp = discreteInterpolant.ComputeL2Error(gradp_coef);
310  double errProj = exact_proj.ComputeL2Error(gradp_coef);
311 
312  cout << "\n Solution of (E_h,v) = (grad p_h,v) for E_h and v in H(curl): "
313  "|| E_h - grad p ||_{L_2} = " << errSol << '\n' << endl;
314  cout << " Gradient interpolant E_h = grad p_h in H(curl): || E_h - grad p"
315  " ||_{L_2} = " << errInterp << '\n' << endl;
316  cout << " Projection E_h of exact grad p in H(curl): || E_h - grad p "
317  "||_{L_2} = " << errProj << '\n' << endl;
318  }
319  else if (prob == 1)
320  {
321  double errSol = x.ComputeL2Error(curlv_coef);
322  double errInterp = discreteInterpolant.ComputeL2Error(curlv_coef);
323  double errProj = exact_proj.ComputeL2Error(curlv_coef);
324 
325  cout << "\n Solution of (E_h,w) = (curl v_h,w) for E_h and w in H(div): "
326  "|| E_h - curl v ||_{L_2} = " << errSol << '\n' << endl;
327  cout << " Curl interpolant E_h = curl v_h in H(div): || E_h - curl v "
328  "||_{L_2} = " << errInterp << '\n' << endl;
329  cout << " Projection E_h of exact curl v in H(div): || E_h - curl v "
330  "||_{L_2} = " << errProj << '\n' << endl;
331  }
332  else
333  {
334  int order_quad = max(2, 2*order+1);
335  const IntegrationRule *irs[Geometry::NumGeom];
336  for (int i=0; i < Geometry::NumGeom; ++i)
337  {
338  irs[i] = &(IntRules.Get(i, order_quad));
339  }
340 
341  double errSol = x.ComputeL2Error(divgradp_coef, irs);
342  double errInterp = discreteInterpolant.ComputeL2Error(divgradp_coef, irs);
343  double errProj = exact_proj.ComputeL2Error(divgradp_coef, irs);
344 
345  cout << "\n Solution of (f_h,q) = (div v_h,q) for f_h and q in L_2: "
346  "|| f_h - div v ||_{L_2} = " << errSol << '\n' << endl;
347  cout << " Divergence interpolant f_h = div v_h in L_2: || f_h - div v "
348  "||_{L_2} = " << errInterp << '\n' << endl;
349  cout << " Projection f_h of exact div v in L_2: || f_h - div v "
350  "||_{L_2} = " << errProj << '\n' << endl;
351  }
352 
353  // 13. Save the refined mesh and the solution. This output can be viewed
354  // later using GLVis: "glvis -m refined.mesh -g sol.gf".
355  ofstream mesh_ofs("refined.mesh");
356  mesh_ofs.precision(8);
357  mesh->Print(mesh_ofs);
358  ofstream sol_ofs("sol.gf");
359  sol_ofs.precision(8);
360  x.Save(sol_ofs);
361 
362  // 14. Send the solution by socket to a GLVis server.
363  if (visualization)
364  {
365  char vishost[] = "localhost";
366  int visport = 19916;
367  socketstream sol_sock(vishost, visport);
368  sol_sock.precision(8);
369  sol_sock << "solution\n" << *mesh << x << flush;
370  }
371 
372  // 15. Free the used memory.
373  delete trial_fec;
374  delete test_fec;
375  delete mesh;
376 
377  return 0;
378 }
379 
380 double p_exact(const Vector &x)
381 {
382  if (dim == 3)
383  {
384  return sin(x(0)) * sin(x(1)) * sin(x(2));
385  }
386  else if (dim == 2)
387  {
388  return sin(x(0)) * sin(x(1));
389  }
390 
391  return 0.0;
392 }
393 
394 void gradp_exact(const Vector &x, Vector &f)
395 {
396  if (dim == 3)
397  {
398  f(0) = cos(x(0)) * sin(x(1)) * sin(x(2));
399  f(1) = sin(x(0)) * cos(x(1)) * sin(x(2));
400  f(2) = sin(x(0)) * sin(x(1)) * cos(x(2));
401  }
402  else
403  {
404  f(0) = cos(x(0)) * sin(x(1));
405  f(1) = sin(x(0)) * cos(x(1));
406  if (x.Size() == 3) { f(2) = 0.0; }
407  }
408 }
409 
410 double div_gradp_exact(const Vector &x)
411 {
412  if (dim == 3)
413  {
414  return -3.0 * sin(x(0)) * sin(x(1)) * sin(x(2));
415  }
416  else if (dim == 2)
417  {
418  return -2.0 * sin(x(0)) * sin(x(1));
419  }
420 
421  return 0.0;
422 }
423 
424 void v_exact(const Vector &x, Vector &v)
425 {
426  if (dim == 3)
427  {
428  v(0) = sin(kappa * x(1));
429  v(1) = sin(kappa * x(2));
430  v(2) = sin(kappa * x(0));
431  }
432  else
433  {
434  v(0) = sin(kappa * x(1));
435  v(1) = sin(kappa * x(0));
436  if (x.Size() == 3) { v(2) = 0.0; }
437  }
438 }
439 
440 void curlv_exact(const Vector &x, Vector &cv)
441 {
442  if (dim == 3)
443  {
444  cv(0) = -kappa * cos(kappa * x(2));
445  cv(1) = -kappa * cos(kappa * x(0));
446  cv(2) = -kappa * cos(kappa * x(1));
447  }
448  else
449  {
450  cv = 0.0;
451  }
452 }
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1387
Conjugate gradient method.
Definition: solvers.hpp:316
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
Data type for scaled Jacobi-type smoother of sparse matrix.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:920
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:143
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:456
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:616
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
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
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void AddDomainInterpolator(DiscreteInterpolator *di)
Adds a domain interpolator. Assumes ownership of di.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds a domain integrator. Assumes ownership of bfi.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3484
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:137
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
Data type sparse matrix.
Definition: sparsemat.hpp:41
void curlv_exact(const Vector &x, Vector &cv)
Definition: ex24.cpp:440
constexpr char vishost[]
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:130
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
void Assemble(int skip_zeros=1)
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:543
void v_exact(const Vector &x, Vector &v)
Definition: ex24.cpp:424
int Dimension() const
Definition: mesh.hpp:911
virtual void ReorientTetMesh()
Definition: mesh.cpp:6376
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
prob_type prob
Definition: ex25.cpp:153
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:341
double p_exact(const Vector &x)
Definition: ex24.cpp:380
A general vector function coefficient.
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACY.
int SpaceDimension() const
Definition: mesh.hpp:912
void SetRelTol(double rtol)
Definition: solvers.hpp:98
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:613
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
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix: .
double a
Definition: lissajous.cpp:41
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
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:327
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2278
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:330
A general function coefficient.
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:416
Vector data type.
Definition: vector.hpp:60
double div_gradp_exact(const Vector &x)
Definition: ex24.cpp:410
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:92
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:121
virtual void Assemble(int skip_zeros=1)
Construct the internal matrix representation of the discrete linear operator.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:377
virtual void Mult(const Vector &x, Vector &y) const
Matrix multiplication: .
int main()
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:285
double f(const Vector &p)
void gradp_exact(const Vector &, Vector &)
Definition: ex24.cpp:394
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150