MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex4.cpp
Go to the documentation of this file.
1 // MFEM Example 4
2 //
3 // Compile with: make ex4
4 //
5 // Sample runs: ex4 -m ../data/square-disc.mesh
6 // ex4 -m ../data/star.mesh
7 // ex4 -m ../data/beam-tet.mesh
8 // ex4 -m ../data/beam-hex.mesh
9 // ex4 -m ../data/beam-hex.mesh -o 2 -pa
10 // ex4 -m ../data/escher.mesh
11 // ex4 -m ../data/fichera.mesh -o 2 -hb
12 // ex4 -m ../data/fichera-q2.vtk
13 // ex4 -m ../data/fichera-q3.mesh -o 2 -sc
14 // ex4 -m ../data/square-disc-nurbs.mesh
15 // ex4 -m ../data/beam-hex-nurbs.mesh
16 // ex4 -m ../data/periodic-square.mesh -no-bc
17 // ex4 -m ../data/periodic-cube.mesh -no-bc
18 // ex4 -m ../data/amr-quad.mesh
19 // ex4 -m ../data/amr-hex.mesh
20 // ex4 -m ../data/amr-hex.mesh -o 2 -hb
21 // ex4 -m ../data/fichera-amr.mesh -o 2 -sc
22 // ex4 -m ../data/star-surf.mesh -o 1
23 //
24 // Device sample runs:
25 // ex4 -m ../data/star.mesh -pa -d cuda
26 // ex4 -m ../data/star.mesh -pa -d raja-cuda
27 // ex4 -m ../data/star.mesh -pa -d raja-omp
28 // ex4 -m ../data/beam-hex.mesh -pa -d cuda
29 //
30 // Description: This example code solves a simple 2D/3D H(div) diffusion
31 // problem corresponding to the second order definite equation
32 // -grad(alpha div F) + beta F = f with boundary condition F dot n
33 // = <given normal field>. Here, we use a given exact solution F
34 // and compute the corresponding r.h.s. f. We discretize with
35 // Raviart-Thomas finite elements.
36 //
37 // The example demonstrates the use of H(div) finite element
38 // spaces with the grad-div and H(div) vector finite element mass
39 // bilinear form, as well as the computation of discretization
40 // error when the exact solution is known. Bilinear form
41 // hybridization and static condensation are also illustrated.
42 //
43 // We recommend viewing examples 1-3 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, F, and r.h.s., f. See below for implementation.
53 void F_exact(const Vector &, Vector &);
54 void f_exact(const Vector &, Vector &);
55 double freq = 1.0, kappa;
56 
57 int main(int argc, char *argv[])
58 {
59  // 1. Parse command-line options.
60  const char *mesh_file = "../data/star.mesh";
61  int order = 1;
62  bool set_bc = true;
63  bool static_cond = false;
64  bool hybridization = false;
65  bool pa = false;
66  const char *device_config = "cpu";
67  bool visualization = 1;
68 
69  OptionsParser args(argc, argv);
70  args.AddOption(&mesh_file, "-m", "--mesh",
71  "Mesh file to use.");
72  args.AddOption(&order, "-o", "--order",
73  "Finite element order (polynomial degree).");
74  args.AddOption(&set_bc, "-bc", "--impose-bc", "-no-bc", "--dont-impose-bc",
75  "Impose or not essential boundary conditions.");
76  args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
77  " solution.");
78  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
79  "--no-static-condensation", "Enable static condensation.");
80  args.AddOption(&hybridization, "-hb", "--hybridization", "-no-hb",
81  "--no-hybridization", "Enable hybridization.");
82  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
83  "--no-partial-assembly", "Enable Partial Assembly.");
84  args.AddOption(&device_config, "-d", "--device",
85  "Device configuration string, see Device::Configure().");
86  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
87  "--no-visualization",
88  "Enable or disable GLVis visualization.");
89  args.Parse();
90  if (!args.Good())
91  {
92  args.PrintUsage(cout);
93  return 1;
94  }
95  args.PrintOptions(cout);
96  kappa = freq * M_PI;
97 
98  // 2. Enable hardware devices such as GPUs, and programming models such as
99  // CUDA, OCCA, RAJA and OpenMP based on command line options.
100  Device device(device_config);
101  device.Print();
102 
103  // 3. Read the mesh from the given mesh file. We can handle triangular,
104  // quadrilateral, tetrahedral, hexahedral, surface and volume, as well as
105  // periodic meshes with the same code.
106  Mesh *mesh = new Mesh(mesh_file, 1, 1);
107  int dim = mesh->Dimension();
108  int sdim = mesh->SpaceDimension();
109 
110  // 4. Refine the mesh to increase the resolution. In this example we do
111  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
112  // largest number that gives a final mesh with no more than 25,000
113  // elements.
114  {
115  int ref_levels =
116  (int)floor(log(25000./mesh->GetNE())/log(2.)/dim);
117  for (int l = 0; l < ref_levels; l++)
118  {
119  mesh->UniformRefinement();
120  }
121  }
122 
123  // 5. Define a finite element space on the mesh. Here we use the
124  // Raviart-Thomas finite elements of the specified order.
125  FiniteElementCollection *fec = new RT_FECollection(order-1, dim);
126  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
127  cout << "Number of finite element unknowns: "
128  << fespace->GetTrueVSize() << endl;
129 
130  // 6. Determine the list of true (i.e. conforming) essential boundary dofs.
131  // In this example, the boundary conditions are defined by marking all
132  // the boundary attributes from the mesh as essential (Dirichlet) and
133  // converting them to a list of true dofs.
134  Array<int> ess_tdof_list;
135  if (mesh->bdr_attributes.Size())
136  {
137  Array<int> ess_bdr(mesh->bdr_attributes.Max());
138  ess_bdr = set_bc ? 1 : 0;
139  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
140  }
141 
142  // 7. Set up the linear form b(.) which corresponds to the right-hand side
143  // of the FEM linear system, which in this case is (f,phi_i) where f is
144  // given by the function f_exact and phi_i are the basis functions in the
145  // finite element fespace.
147  LinearForm *b = new LinearForm(fespace);
149  b->Assemble();
150 
151  // 8. Define the solution vector x as a finite element grid function
152  // corresponding to fespace. Initialize x by projecting the exact
153  // solution. Note that only values from the boundary faces will be used
154  // when eliminating the non-homogeneous boundary condition to modify the
155  // r.h.s. vector b.
156  GridFunction x(fespace);
158  x.ProjectCoefficient(F);
159 
160  // 9. Set up the bilinear form corresponding to the H(div) diffusion operator
161  // grad alpha div + beta I, by adding the div-div and the mass domain
162  // integrators.
164  Coefficient *beta = new ConstantCoefficient(1.0);
165  BilinearForm *a = new BilinearForm(fespace);
166  if (pa) { a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
167  a->AddDomainIntegrator(new DivDivIntegrator(*alpha));
169 
170  // 10. Assemble the bilinear form and the corresponding linear system,
171  // applying any necessary transformations such as: eliminating boundary
172  // conditions, applying conforming constraints for non-conforming AMR,
173  // static condensation, hybridization, etc.
174  FiniteElementCollection *hfec = NULL;
175  FiniteElementSpace *hfes = NULL;
176  if (static_cond)
177  {
179  }
180  else if (hybridization)
181  {
182  hfec = new DG_Interface_FECollection(order-1, dim);
183  hfes = new FiniteElementSpace(mesh, hfec);
185  ess_tdof_list);
186  }
187  a->Assemble();
188 
189  OperatorPtr A;
190  Vector B, X;
191  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
192 
193  cout << "Size of linear system: " << A->Height() << endl;
194 
195  // 11. Solve the linear system A X = B.
196  if (!pa)
197  {
198 #ifndef MFEM_USE_SUITESPARSE
199  // Use a simple symmetric Gauss-Seidel preconditioner with PCG.
200  GSSmoother M((SparseMatrix&)(*A));
201  PCG(*A, M, B, X, 1, 10000, 1e-20, 0.0);
202 #else
203  // If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
204  UMFPackSolver umf_solver;
205  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
206  umf_solver.SetOperator(*A);
207  umf_solver.Mult(B, X);
208 #endif
209  }
210  else // Jacobi preconditioning in partial assembly mode
211  {
212  if (UsesTensorBasis(*fespace))
213  {
214  OperatorJacobiSmoother M(*a, ess_tdof_list);
215  PCG(*A, M, B, X, 1, 10000, 1e-20, 0.0);
216  }
217  else
218  {
219  CG(*A, B, X, 1, 10000, 1e-20, 0.0);
220  }
221  }
222 
223  // 12. Recover the solution as a finite element grid function.
224  a->RecoverFEMSolution(X, *b, x);
225 
226  // 13. Compute and print the L^2 norm of the error.
227  cout << "\n|| F_h - F ||_{L^2} = " << x.ComputeL2Error(F) << '\n' << endl;
228 
229  // 14. Save the refined mesh and the solution. This output can be viewed
230  // later using GLVis: "glvis -m refined.mesh -g sol.gf".
231  {
232  ofstream mesh_ofs("refined.mesh");
233  mesh_ofs.precision(8);
234  mesh->Print(mesh_ofs);
235  ofstream sol_ofs("sol.gf");
236  sol_ofs.precision(8);
237  x.Save(sol_ofs);
238  }
239 
240  // 15. Send the solution by socket to a GLVis server.
241  if (visualization)
242  {
243  char vishost[] = "localhost";
244  int visport = 19916;
245  socketstream sol_sock(vishost, visport);
246  sol_sock.precision(8);
247  sol_sock << "solution\n" << *mesh << x << flush;
248  }
249 
250  // 16. Free the used memory.
251  delete hfes;
252  delete hfec;
253  delete a;
254  delete alpha;
255  delete beta;
256  delete b;
257  delete fespace;
258  delete fec;
259  delete mesh;
260 
261  return 0;
262 }
263 
264 
265 // The exact solution (for non-surface meshes)
266 void F_exact(const Vector &p, Vector &F)
267 {
268  int dim = p.Size();
269 
270  double x = p(0);
271  double y = p(1);
272  // double z = (dim == 3) ? p(2) : 0.0; // Uncomment if F is changed to depend on z
273 
274  F(0) = cos(kappa*x)*sin(kappa*y);
275  F(1) = cos(kappa*y)*sin(kappa*x);
276  if (dim == 3)
277  {
278  F(2) = 0.0;
279  }
280 }
281 
282 // The right hand side
283 void f_exact(const Vector &p, Vector &f)
284 {
285  int dim = p.Size();
286 
287  double x = p(0);
288  double y = p(1);
289  // double z = (dim == 3) ? p(2) : 0.0; // Uncomment if f is changed to depend on z
290 
291  double temp = 1 + 2*kappa*kappa;
292 
293  f(0) = temp*cos(kappa*x)*sin(kappa*y);
294  f(1) = temp*cos(kappa*y)*sin(kappa*x);
295  if (dim == 3)
296  {
297  f(2) = 0;
298  }
299 }
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1234
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
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:431
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
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 RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Recover the solution of a linear system formed with FormLinearSystem().
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition: fespace.cpp:415
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:261
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
double kappa
Definition: ex24.cpp:54
bool UsesTensorBasis(const FiniteElementSpace &fes)
Definition: fespace.hpp:983
int main(int argc, char *argv[])
Definition: ex1.cpp:66
(Q div u, div v) for RT elements
Data type for Gauss-Seidel smoother of sparse matrix.
void EnableHybridization(FiniteElementSpace *constr_space, BilinearFormIntegrator *constr_integ, const Array< int > &ess_tdof_list)
Enable hybridization.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:731
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3417
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:46
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:65
constexpr char vishost[]
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:128
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:8382
constexpr int visport
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
void F_exact(const Vector &, Vector &)
Definition: ex4.cpp:266
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:720
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:733
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:403
int Dimension() const
Definition: mesh.hpp:788
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:272
A general vector function coefficient.
int SpaceDimension() const
Definition: mesh.hpp:789
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:201
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:742
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
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:257
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...
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:304
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:2252
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:260
const double alpha
Definition: ex15.cpp:336
Vector data type.
Definition: vector.hpp:51
Vector with associated FE space and LinearFormIntegrators.
Definition: linearform.hpp:23
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:118
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:2818
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
double f(const Vector &p)
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:2721
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145