MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex3.cpp
Go to the documentation of this file.
1 // MFEM Example 3
2 //
3 // Compile with: make ex3
4 //
5 // Sample runs: ex3 -m ../data/star.mesh
6 // ex3 -m ../data/beam-tri.mesh -o 2
7 // ex3 -m ../data/beam-tet.mesh
8 // ex3 -m ../data/beam-hex.mesh
9 // ex3 -m ../data/beam-hex.mesh -o 2 -pa
10 // ex3 -m ../data/escher.mesh
11 // ex3 -m ../data/escher.mesh -o 2
12 // ex3 -m ../data/fichera.mesh
13 // ex3 -m ../data/fichera-q2.vtk
14 // ex3 -m ../data/fichera-q3.mesh
15 // ex3 -m ../data/square-disc-nurbs.mesh
16 // ex3 -m ../data/beam-hex-nurbs.mesh
17 // ex3 -m ../data/amr-hex.mesh
18 // ex3 -m ../data/fichera-amr.mesh
19 // ex3 -m ../data/ref-prism.mesh -o 1
20 // ex3 -m ../data/octahedron.mesh -o 1
21 // ex3 -m ../data/star-surf.mesh -o 1
22 // ex3 -m ../data/mobius-strip.mesh -f 0.1
23 // ex3 -m ../data/klein-bottle.mesh -f 0.1
24 //
25 // Device sample runs:
26 // ex3 -m ../data/star.mesh -pa -d cuda
27 // ex3 -m ../data/star.mesh -pa -d raja-cuda
28 // ex3 -m ../data/star.mesh -pa -d raja-omp
29 // ex3 -m ../data/beam-hex.mesh -pa -d cuda
30 //
31 // Description: This example code solves a simple electromagnetic diffusion
32 // problem corresponding to the second order definite Maxwell
33 // equation curl curl E + E = f with boundary condition
34 // E x n = <given tangential field>. Here, we use a given exact
35 // solution E and compute the corresponding r.h.s. f.
36 // We discretize with Nedelec finite elements in 2D or 3D.
37 //
38 // The example demonstrates the use of H(curl) finite element
39 // spaces with the curl-curl and the (vector finite element) mass
40 // bilinear form, as well as the computation of discretization
41 // error when the exact solution is known. Static condensation is
42 // also illustrated.
43 //
44 // We recommend viewing examples 1-2 before viewing this example.
45 
46 #include "mfem.hpp"
47 #include <fstream>
48 #include <iostream>
49 
50 using namespace std;
51 using namespace mfem;
52 
53 // Exact solution, E, and r.h.s., f. See below for implementation.
54 void E_exact(const Vector &, Vector &);
55 void f_exact(const Vector &, Vector &);
56 double freq = 1.0, kappa;
57 int dim;
58 
59 int main(int argc, char *argv[])
60 {
61  // 1. Parse command-line options.
62  const char *mesh_file = "../data/beam-tet.mesh";
63  int order = 1;
64  bool static_cond = 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(&freq, "-f", "--frequency", "Set the frequency for the exact"
75  " solution.");
76  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
77  "--no-static-condensation", "Enable static condensation.");
78  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
79  "--no-partial-assembly", "Enable Partial Assembly.");
80  args.AddOption(&device_config, "-d", "--device",
81  "Device configuration string, see Device::Configure().");
82  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
83  "--no-visualization",
84  "Enable or disable GLVis visualization.");
85  args.Parse();
86  if (!args.Good())
87  {
88  args.PrintUsage(cout);
89  return 1;
90  }
91  args.PrintOptions(cout);
92  kappa = freq * M_PI;
93 
94  // 2. Enable hardware devices such as GPUs, and programming models such as
95  // CUDA, OCCA, RAJA and OpenMP based on command line options.
96  Device device(device_config);
97  device.Print();
98 
99  // 3. Read the mesh from the given mesh file. We can handle triangular,
100  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
101  // the same code.
102  Mesh *mesh = new Mesh(mesh_file, 1, 1);
103  dim = mesh->Dimension();
104  int sdim = mesh->SpaceDimension();
105 
106  // 4. Refine the mesh to increase the resolution. In this example we do
107  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
108  // largest number that gives a final mesh with no more than 50,000
109  // elements.
110  {
111  int ref_levels =
112  (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
113  for (int l = 0; l < ref_levels; l++)
114  {
115  mesh->UniformRefinement();
116  }
117  }
118 
119  // 5. Define a finite element space on the mesh. Here we use the Nedelec
120  // finite elements of the specified order.
121  FiniteElementCollection *fec = new ND_FECollection(order, dim);
122  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
123  cout << "Number of finite element unknowns: "
124  << fespace->GetTrueVSize() << endl;
125 
126  // 6. Determine the list of true (i.e. conforming) essential boundary dofs.
127  // In this example, the boundary conditions are defined by marking all
128  // the boundary attributes from the mesh as essential (Dirichlet) and
129  // converting them to a list of true dofs.
130  Array<int> ess_tdof_list;
131  if (mesh->bdr_attributes.Size())
132  {
133  Array<int> ess_bdr(mesh->bdr_attributes.Max());
134  ess_bdr = 1;
135  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
136  }
137 
138  // 7. Set up the linear form b(.) which corresponds to the right-hand side
139  // of the FEM linear system, which in this case is (f,phi_i) where f is
140  // given by the function f_exact and phi_i are the basis functions in the
141  // finite element fespace.
143  LinearForm *b = new LinearForm(fespace);
145  b->Assemble();
146 
147  // 8. Define the solution vector x as a finite element grid function
148  // corresponding to fespace. Initialize x by projecting the exact
149  // solution. Note that only values from the boundary edges will be used
150  // when eliminating the non-homogeneous boundary condition to modify the
151  // r.h.s. vector b.
152  GridFunction x(fespace);
154  x.ProjectCoefficient(E);
155 
156  // 9. Set up the bilinear form corresponding to the EM diffusion operator
157  // curl muinv curl + sigma I, by adding the curl-curl and the mass domain
158  // integrators.
159  Coefficient *muinv = new ConstantCoefficient(1.0);
161  BilinearForm *a = new BilinearForm(fespace);
162  if (pa) { a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
163  a->AddDomainIntegrator(new CurlCurlIntegrator(*muinv));
165 
166  // 10. Assemble the bilinear form and the corresponding linear system,
167  // applying any necessary transformations such as: eliminating boundary
168  // conditions, applying conforming constraints for non-conforming AMR,
169  // static condensation, etc.
170  if (static_cond) { a->EnableStaticCondensation(); }
171  a->Assemble();
172 
173  OperatorPtr A;
174  Vector B, X;
175  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
176 
177  cout << "Size of linear system: " << A->Height() << endl;
178 
179  // 11. Solve the linear system A X = B.
180  if (pa) // Jacobi preconditioning in partial assembly mode
181  {
182  OperatorJacobiSmoother M(*a, ess_tdof_list);
183  PCG(*A, M, B, X, 1, 1000, 1e-12, 0.0);
184  }
185  else
186  {
187 #ifndef MFEM_USE_SUITESPARSE
188  // 11. Define a simple symmetric Gauss-Seidel preconditioner and use it to
189  // solve the system Ax=b with PCG.
190  GSSmoother M((SparseMatrix&)(*A));
191  PCG(*A, M, B, X, 1, 500, 1e-12, 0.0);
192 #else
193  // 11. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the
194  // system.
195  UMFPackSolver umf_solver;
196  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
197  umf_solver.SetOperator(*A);
198  umf_solver.Mult(B, X);
199 #endif
200  }
201 
202  // 12. Recover the solution as a finite element grid function.
203  a->RecoverFEMSolution(X, *b, x);
204 
205  // 13. Compute and print the L^2 norm of the error.
206  cout << "\n|| E_h - E ||_{L^2} = " << x.ComputeL2Error(E) << '\n' << endl;
207 
208  // 14. Save the refined mesh and the solution. This output can be viewed
209  // later using GLVis: "glvis -m refined.mesh -g sol.gf".
210  {
211  ofstream mesh_ofs("refined.mesh");
212  mesh_ofs.precision(8);
213  mesh->Print(mesh_ofs);
214  ofstream sol_ofs("sol.gf");
215  sol_ofs.precision(8);
216  x.Save(sol_ofs);
217  }
218 
219  // 15. Send the solution by socket to a GLVis server.
220  if (visualization)
221  {
222  char vishost[] = "localhost";
223  int visport = 19916;
224  socketstream sol_sock(vishost, visport);
225  sol_sock.precision(8);
226  sol_sock << "solution\n" << *mesh << x << flush;
227  }
228 
229  // 16. Free the used memory.
230  delete a;
231  delete sigma;
232  delete muinv;
233  delete b;
234  delete fespace;
235  delete fec;
236  delete mesh;
237 
238  return 0;
239 }
240 
241 
242 void E_exact(const Vector &x, Vector &E)
243 {
244  if (dim == 3)
245  {
246  E(0) = sin(kappa * x(1));
247  E(1) = sin(kappa * x(2));
248  E(2) = sin(kappa * x(0));
249  }
250  else
251  {
252  E(0) = sin(kappa * x(1));
253  E(1) = sin(kappa * x(0));
254  if (x.Size() == 3) { E(2) = 0.0; }
255  }
256 }
257 
258 void f_exact(const Vector &x, Vector &f)
259 {
260  if (dim == 3)
261  {
262  f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
263  f(1) = (1. + kappa * kappa) * sin(kappa * x(2));
264  f(2) = (1. + kappa * kappa) * sin(kappa * x(0));
265  }
266  else
267  {
268  f(0) = (1. + kappa * kappa) * sin(kappa * x(1));
269  f(1) = (1. + kappa * kappa) * sin(kappa * x(0));
270  if (x.Size() == 3) { f(2) = 0.0; }
271  }
272 }
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:84
Integrator for (curl u, curl v) for Nedelec elements.
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:168
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:200
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:923
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
Data type for Gauss-Seidel smoother of sparse matrix.
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:3676
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:50
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:67
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:9816
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 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:590
int Dimension() const
Definition: mesh.hpp:1006
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
A general vector function coefficient.
int SpaceDimension() const
Definition: mesh.hpp:1007
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:41
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
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:41
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
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Definition: gridfunc.cpp:2766
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:1671
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:2399
void E_exact(const Vector &, Vector &)
Definition: ex3.cpp:242
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:346
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:415
Vector data type.
Definition: vector.hpp:60
Vector with associated FE space and LinearFormIntegrators.
Definition: linearform.hpp:24
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:3189
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:3094
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