MFEM  v4.1.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/escher.mesh
10 // ex24 -m ../data/escher.mesh -o 2
11 // ex24 -m ../data/fichera.mesh
12 // ex24 -m ../data/fichera-q2.vtk
13 // ex24 -m ../data/fichera-q3.mesh
14 // ex24 -m ../data/square-disc-nurbs.mesh
15 // ex24 -m ../data/beam-hex-nurbs.mesh
16 // ex24 -m ../data/amr-quad.mesh -o 2
17 // ex24 -m ../data/amr-hex.mesh
18 //
19 // Device sample runs:
20 // ex24 -m ../data/star.mesh -pa -d cuda
21 // ex24 -m ../data/star.mesh -pa -d raja-cuda
22 // ex24 -m ../data/star.mesh -pa -d raja-omp
23 // ex24 -m ../data/beam-hex.mesh -pa -d cuda
24 //
25 // Description: This example code illustrates usage of mixed finite element
26 // spaces. Using two different approaches, we project a gradient
27 // of a function in H^1 to H(curl). Other spaces and example
28 // computations are to be added in the future.
29 //
30 // We recommend viewing examples 1 and 3 before viewing this
31 // example.
32 
33 #include "mfem.hpp"
34 #include <fstream>
35 #include <iostream>
36 
37 using namespace std;
38 using namespace mfem;
39 
40 double p_exact(const Vector &x);
41 void gradp_exact(const Vector &, Vector &);
42 
43 int dim;
44 
45 int main(int argc, char *argv[])
46 {
47  // 1. Parse command-line options.
48  const char *mesh_file = "../data/beam-hex.mesh";
49  int order = 1;
50  bool static_cond = false;
51  bool pa = false;
52  const char *device_config = "cpu";
53  bool visualization = 1;
54 
55  OptionsParser args(argc, argv);
56  args.AddOption(&mesh_file, "-m", "--mesh",
57  "Mesh file to use.");
58  args.AddOption(&order, "-o", "--order",
59  "Finite element order (polynomial degree).");
60  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
61  "--no-static-condensation", "Enable static condensation.");
62  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
63  "--no-partial-assembly", "Enable Partial Assembly.");
64  args.AddOption(&device_config, "-d", "--device",
65  "Device configuration string, see Device::Configure().");
66  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
67  "--no-visualization",
68  "Enable or disable GLVis visualization.");
69 
70  args.Parse();
71  if (!args.Good())
72  {
73  args.PrintUsage(cout);
74  return 1;
75  }
76  args.PrintOptions(cout);
77 
78  // 2. Enable hardware devices such as GPUs, and programming models such as
79  // CUDA, OCCA, RAJA and OpenMP based on command line options.
80  Device device(device_config);
81  device.Print();
82 
83  // 3. Read the mesh from the given mesh file. We can handle triangular,
84  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
85  // the same code.
86  Mesh *mesh = new Mesh(mesh_file, 1, 1);
87  dim = mesh->Dimension();
88  int sdim = mesh->SpaceDimension();
89 
90  // 4. Refine the mesh to increase the resolution. In this example we do
91  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
92  // largest number that gives a final mesh with no more than 50,000
93  // elements.
94  {
95  int ref_levels = (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
96  for (int l = 0; l < ref_levels; l++)
97  {
98  mesh->UniformRefinement();
99  }
100  }
101  mesh->ReorientTetMesh();
102 
103  // 5. Define a parallel finite element space on the parallel mesh. Here we
104  // use the Nedelec finite elements of the specified order.
105  FiniteElementCollection *fec = new ND_FECollection(order, dim);
106  FiniteElementCollection *H1fec = new H1_FECollection(order, dim);
107  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
108  FiniteElementSpace *H1fespace = new FiniteElementSpace(mesh, H1fec);
109 
110  int size = fespace->GetTrueVSize();
111  int H1size = H1fespace->GetTrueVSize();
112  cout << "Number of Nedelec finite element unknowns: " << size << endl;
113  cout << "Number of H1 finite element unknowns: " << H1size << endl;
114 
115  // 6. Define the solution vector x as a parallel finite element grid function
116  // corresponding to fespace. Initialize x by projecting the exact
117  // solution. Note that only values from the boundary edges will be used
118  // when eliminating the non-homogeneous boundary condition to modify the
119  // r.h.s. vector b.
120  GridFunction x(fespace);
122  GridFunction p(H1fespace);
123  p.ProjectCoefficient(p_coef);
124  p.SetTrueVector();
125  p.SetFromTrueVector();
126 
127  VectorFunctionCoefficient gradp_coef(sdim, gradp_exact);
128 
129  // 7. Set up the bilinear forms.
130  Coefficient *muinv = new ConstantCoefficient(1.0);
132  BilinearForm *a = new BilinearForm(fespace);
133  MixedBilinearForm *a_NDH1 = new MixedBilinearForm(H1fespace, fespace);
134  if (pa)
135  {
136  a->SetAssemblyLevel(AssemblyLevel::PARTIAL);
137  a_NDH1->SetAssemblyLevel(AssemblyLevel::PARTIAL);
138  }
139 
140  // First approach: L2 projection
143 
144  // 8. Assemble the parallel bilinear form and the corresponding linear
145  // system, applying any necessary transformations such as: parallel
146  // assembly, eliminating boundary conditions, applying conforming
147  // constraints for non-conforming AMR, static condensation, etc.
148  if (static_cond) { a->EnableStaticCondensation(); }
149 
150  a->Assemble();
151  if (!pa) { a->Finalize(); }
152 
153  a_NDH1->Assemble();
154  if (!pa) { a_NDH1->Finalize(); }
155 
156  if (pa)
157  {
158  a_NDH1->Mult(p, x);
159  }
160  else
161  {
162  SparseMatrix& NDH1 = a_NDH1->SpMat();
163  NDH1.Mult(p, x);
164  }
165 
166  // 9. Define and apply a PCG solver for Ax = b with Jacobi preconditioner.
167  {
168  GridFunction rhs(fespace);
169  rhs = x;
170  x = 0.0;
171 
172  CGSolver cg;
173  cg.SetRelTol(1e-12);
174  cg.SetMaxIter(1000);
175  cg.SetPrintLevel(1);
176  if (pa)
177  {
178  Array<int> ess_tdof_list; // empty
179  OperatorJacobiSmoother Jacobi(*a, ess_tdof_list);
180 
181  cg.SetOperator(*a);
182  cg.SetPreconditioner(Jacobi);
183  cg.Mult(rhs, x);
184  }
185  else
186  {
187  SparseMatrix& Amat = a->SpMat();
188  DSmoother Jacobi(Amat);
189 
190  cg.SetOperator(Amat);
191  cg.SetPreconditioner(Jacobi);
192  cg.Mult(rhs, x);
193  }
194  }
195 
196  // 10. Second approach: compute the same solution by applying
197  // GradientInterpolator in H(curl).
198  DiscreteLinearOperator grad(H1fespace, fespace);
200  grad.Assemble();
201 
202  GridFunction gradp(fespace);
203  grad.Mult(p, gradp);
204 
205  // 11. Compute the projection of the exact grad p.
206  GridFunction exact_gradp(fespace);
207  exact_gradp.ProjectCoefficient(gradp_coef);
208  exact_gradp.SetTrueVector();
209  exact_gradp.SetFromTrueVector();
210 
211  // 12. Compute and print the L^2 norm of the error.
212  {
213  double errSol = x.ComputeL2Error(gradp_coef);
214  double errInterp = gradp.ComputeL2Error(gradp_coef);
215  double errProj = exact_gradp.ComputeL2Error(gradp_coef);
216 
217  cout << "\n Solution of (E_h,v) = (grad p_h,v) for E_h and v in H(curl): "
218  "|| E_h - grad p ||_{L^2} = " << errSol << '\n' << endl;
219  cout << " Gradient interpolant E_h = grad p_h in H(curl): || E_h - grad p"
220  "||_{L^2} = " << errInterp << '\n' << endl;
221  cout << " Projection E_h of exact grad p in H(curl): || E_h - grad p "
222  "||_{L^2} = " << errProj << '\n' << endl;
223  }
224 
225  // 13. Save the refined mesh and the solution. This output can be viewed
226  // later using GLVis: "glvis -m refined.mesh -g sol.gf".
227  ofstream mesh_ofs("refined.mesh");
228  mesh_ofs.precision(8);
229  mesh->Print(mesh_ofs);
230  ofstream sol_ofs("sol.gf");
231  sol_ofs.precision(8);
232  x.Save(sol_ofs);
233 
234  // 14. Send the solution by socket to a GLVis server.
235  if (visualization)
236  {
237  char vishost[] = "localhost";
238  int visport = 19916;
239  socketstream sol_sock(vishost, visport);
240  sol_sock.precision(8);
241  sol_sock << "solution\n" << *mesh << x << flush;
242  }
243 
244  // 15. Free the used memory.
245  delete a;
246  delete a_NDH1;
247  delete sigma;
248  delete muinv;
249  delete fespace;
250  delete H1fespace;
251  delete fec;
252  delete H1fec;
253  delete mesh;
254 
255  return 0;
256 }
257 
258 double p_exact(const Vector &x)
259 {
260  if (dim == 3)
261  {
262  return sin(x(0)) * sin(x(1)) * sin(x(2));
263  }
264  else if (dim == 2)
265  {
266  return sin(x(0)) * sin(x(1));
267  }
268 
269  return 0.0;
270 }
271 
272 void gradp_exact(const Vector &x, Vector &f)
273 {
274  if (dim == 3)
275  {
276  f(0) = cos(x(0)) * sin(x(1)) * sin(x(2));
277  f(1) = sin(x(0)) * cos(x(1)) * sin(x(2));
278  f(2) = sin(x(0)) * sin(x(1)) * cos(x(2));
279  }
280  else
281  {
282  f(0) = cos(x(0)) * sin(x(1));
283  f(1) = sin(x(0)) * cos(x(1));
284  if (x.Size() == 3) { f(2) = 0.0; }
285  }
286 }
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1188
Conjugate gradient method.
Definition: solvers.hpp:152
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
Data type for scaled Jacobi-type smoother of sparse matrix.
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:139
Subclass constant coefficient.
Definition: coefficient.hpp:67
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:318
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:358
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:693
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:236
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::FULL.
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
int main(int argc, char *argv[])
Definition: ex1.cpp:62
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:2731
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:133
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
Data type sparse matrix.
Definition: sparsemat.hpp:40
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:84
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7982
void SetMaxIter(int max_it)
Definition: solvers.hpp:65
void Assemble(int skip_zeros=1)
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:384
int Dimension() const
Definition: mesh.hpp:744
virtual void ReorientTetMesh()
Definition: mesh.cpp:5250
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
double p_exact(const Vector &x)
Definition: ex24.cpp:258
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::FULL.
int SpaceDimension() const
Definition: mesh.hpp:745
void SetRelTol(double rtol)
Definition: solvers.hpp:63
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:548
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
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)
Definition: optparser.hpp:76
const SparseMatrix & SpMat() const
double a
Definition: lissajous.cpp:41
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
int dim
Definition: ex24.cpp:43
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1685
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:166
class for C-function coefficient
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:251
Vector data type.
Definition: vector.hpp:48
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:83
const SparseMatrix & SpMat() const
Returns a 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:114
virtual void Assemble(int skip_zeros=1)
Construct the internal matrix representation of the discrete linear operator.
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void EnableStaticCondensation()
void gradp_exact(const Vector &, Vector &)
Definition: ex24.cpp:272
double sigma(const Vector &x)
Definition: maxwell.cpp:102
bool Good() const
Definition: optparser.hpp:125