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