MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
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
19using namespace std;
20using namespace mfem;
21
22int main(int argc, char *argv[])
23{
24 // 1. Parse command line options.
25 string 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);
56 b.AddDomainIntegrator(new DomainLFIntegrator(one));
57 b.Assemble();
58
59 // 7. Set up the bilinear form a(.,.) corresponding to the -Delta operator.
60 BilinearForm a(&fespace);
61 a.AddDomainIntegrator(new DiffusionIntegrator);
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.
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}
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
A coefficient that is constant across space and time.
Class for domain integration .
Definition lininteg.hpp:109
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:716
void GetBoundaryTrueDofs(Array< int > &boundary_dofs, int component=-1)
Get a list of all boundary true dofs, boundary_dofs. For spaces with 'vdim' > 1, the 'component' para...
Definition fespace.cpp:632
Data type for Gauss-Seidel smoother of sparse matrix.
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
Vector with associated FE space and LinearFormIntegrators.
Mesh data type.
Definition mesh.hpp:56
virtual void Save(const std::string &fname, int precision=16) const
Definition mesh.cpp:11481
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
void ParseCheck(std::ostream &out=mfem::out)
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 'var' to receive the value. Enable/disable tags are used to set the bool...
Definition optparser.hpp:82
Data type sparse matrix.
Definition sparsemat.hpp:51
Vector data type.
Definition vector.hpp:80
int main()
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
Definition solvers.cpp:913