MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex0.cpp
Go to the documentation of this file.
1 // MFEM Example 0
2 //
3 // Compile with: make ex0
4 //
5 // Sample runs: ex0
6 // ex0 -m ../data/fichera.mesh
7 // ex0 -m ../data/square-disc.mesh -o 2
8 //
9 // Description: This example code demonstrates the most basic usage of MFEM to
10 // define a simple finite element discretization of the Laplace
11 // problem -Delta u = 1 with zero Dirichlet boundary conditions.
12 // General 2D/3D mesh files and finite element polynomial degrees
13 // can be specified by command line options.
14 
15 #include "mfem.hpp"
16 #include <fstream>
17 #include <iostream>
18 
19 using namespace std;
20 using namespace mfem;
21 
22 int main(int argc, char *argv[])
23 {
24  // 1. Parse command line options.
25  const char *mesh_file = "../data/star.mesh";
26  int order = 1;
27 
28  OptionsParser args(argc, argv);
29  args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use.");
30  args.AddOption(&order, "-o", "--order", "Finite element polynomial degree");
31  args.ParseCheck();
32 
33  // 2. Read the mesh from the given mesh file, and refine once uniformly.
34  Mesh mesh(mesh_file);
35  mesh.UniformRefinement();
36 
37  // 3. Define a finite element space on the mesh. Here we use H1 continuous
38  // high-order Lagrange finite elements of the given order.
39  H1_FECollection fec(order, mesh.Dimension());
40  FiniteElementSpace fespace(&mesh, &fec);
41  cout << "Number of unknowns: " << fespace.GetTrueVSize() << endl;
42 
43  // 4. Extract the list of all the boundary DOFs. These will be marked as
44  // Dirichlet in order to enforce zero boundary conditions.
45  Array<int> boundary_dofs;
46  fespace.GetBoundaryTrueDofs(boundary_dofs);
47 
48  // 5. Define the solution x as a finite element grid function in fespace. Set
49  // the initial guess to zero, which also sets the boundary conditions.
50  GridFunction x(&fespace);
51  x = 0.0;
52 
53  // 6. Set up the linear form b(.) corresponding to the right-hand side.
54  ConstantCoefficient one(1.0);
55  LinearForm b(&fespace);
57  b.Assemble();
58 
59  // 7. Set up the bilinear form a(.,.) corresponding to the -Delta operator.
60  BilinearForm a(&fespace);
62  a.Assemble();
63 
64  // 8. Form the linear system A X = B. This includes eliminating boundary
65  // conditions, applying AMR constraints, and other transformations.
66  SparseMatrix A;
67  Vector B, X;
68  a.FormLinearSystem(boundary_dofs, x, b, A, X, B);
69 
70  // 9. Solve the system using PCG with symmetric Gauss-Seidel preconditioner.
71  GSSmoother M(A);
72  PCG(A, M, B, X, 1, 200, 1e-12, 0.0);
73 
74  // 10. Recover the solution x as a grid function and save to file. The output
75  // can be viewed using GLVis as follows: "glvis -m mesh.mesh -g sol.gf"
76  a.RecoverFEMSolution(X, b, x);
77  x.Save("sol.gf");
78  mesh.Save("mesh.mesh");
79 
80  return 0;
81 }
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:108
virtual void Save(const char *fname, int precision=16) const
Definition: mesh.cpp:10289
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
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:168
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().
Data type for Gauss-Seidel smoother of sparse matrix.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3676
Data type sparse matrix.
Definition: sparsemat.hpp:50
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9816
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
int Dimension() const
Definition: mesh.hpp:1006
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:41
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
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
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...
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
Vector data type.
Definition: vector.hpp:60
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:220
Vector with associated FE space and LinearFormIntegrators.
Definition: linearform.hpp:24
int main()
void ParseCheck(std::ostream &out=mfem::out)
Definition: optparser.cpp:252