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