linearFEM1d-simple.cc
Introduction
One dimensional FEM in Concepts. This code features first order Ansatz functions on a (possibly) uniformly refined mesh.
The equation solved is
on (-1,1) with homogeneous Dirichlet boundary conditions. The linear system is solved using CG.
The exact solution is u = -3 cosh(x)/cosh(1) + x^2 + 2.
The relative energy error of the FE solution is computed at the end of the program
Contents
Commented Program
System include files for input and output operations.
Concepts include files.
With this using directive, we do not need to prepend concepts:
: before Real
everytime.
Mesh
The mesh is taken from the meshes tutorial: Line
. Here, the relevant code is included:
Linear Form
The linear form is used to compute the element contributions to the load vector. The computations are done in the application operator of this class.
As we need numerical integration, the class has a member quad_
of type concepts::Quadrature. The template parameter 4 indicates the type of quadrature rule. It contains the quadrature points and weights for given number of quadrature points (one in our case).
The constructor takes as argument the number of quadrature points and initialises the quadrature object quad_
.
The application operator is called for every element elm
and is responsible for computing the contributions of this element to the global load vector. The local load vector is stored in em
.
Resize the element load vector to the correct size and set all entries to zero. This way, the computed values can be added into em
.
Check if the general element elm
really is an element of type linearFEM::Line which can be handled by our routine. This is achieved by trying to cast elm
to linearFEM::Line. If the cast succeeds (ie. elm
is really of type linearFEM::Line), then e
contains a pointer to this element. If the cast fails, e
is equal to 0 and the assertion fails too. If the assertions fails, an error message is displayed (indicating the name of the file, the line number in the file, the name of the routine where the error happened and the check which failed) and the execution of the program stops.
The element of type linearFEM::Line has a member cell
which returns a reference to a concepts::Edge1d and this class is able to compute the size of the cell (on which the element is based). The size of the cell is needed to compute the Jacobian in the numerical integration.
Next, we get pointers to the first elements of the arrays of quadrature weights and points. These data were computed in the constructor (which initialised the member quad_
).
Now, the computational loops of this routine start.
There are two nested loops here: the outer loop enumerates the shape funtions of this element and the inner loop enumerates the quadrature points.
m
is a reference to the entry of the local load vector which corresponds to the current shape function.
Now, all data is known to add to m:
w
[k] is the quadrature weight associated to the current quadrature pointq
[k].m += w[k]h
is the determinant of the Jacobian of the element map which maps the reference coordinates (0,1) to the physical coordinates of this element.1/2
is the Jacobian of the map from (-1,1) to (0,1). The integration takes place on (-1,1) (as described in concepts::Quadrature).* h / 2.0e->shapefct(i, q[k])
is the value of the shape functioni
at the pointq
[k]. As you can see in linearFEM::Line, it takes values in (-1,1).* e->shapefct(i, q[k])- The last line
pow(e->cell().chi( (1.0+q[k])/2.0 ), 2.0)
computes the value of the right hand side funtion f = x2.(1.0+q[k])/2.0
maps a point in (-1,1) to (0,1) and the functionpow
from the standard math library computes ab.
Main Program
Set up the mesh. The values -1, 0 and 1 give the coordinates of the left, middle and right points in the mesh.
The boundary conditions are given by assigning a boundary type to an attribute. The end points of the interval which are created above have the attribute 2.
The class Boundary
holds a type (Dirichlet or Neumann) and (if necessary) a function or value. These objects can be added to BoundaryConditions
which manages them.
Set up the space according to the current level. The space takes the mesh, the level of refinement and the boundary conditions as input. If there are Dirichlet boundary conditions somewhere, the respective degrees of freedom are not created.
Compute the load vector. This is done by creating an object of type RHS
which is a linear form. This linear form is then given to the vector. The vector fills its entries by evaluating the linear form on every element of the space and assembling the local load vectors according to the assembling information contained in the elements.
Compute the stiffness matrix and the mass matrix. As above, a bilinear form is created and given to a matrix which the assembles the mass and stiffness matrix respectively.
A linear combination L
is formed from the mass and stiffness matrix. The linear system is going to be solved iteratively. This linear combination is used to initialize a CG solver together with the maximal residual and the maximal number of iterations.
An empty vector to hold the solution of the linear system.
The CG solver is essentially a linear operator. Therefore, one can apply a vector to it and get the solution of the linear system which was used to set up the CG solver.
The CG solver throws a exception of type concepts::NoConvergence
which is caught here. Although the solver might not have converged, the solution might give some hints to the user. Therefore, the execution of the program is not stopped here and only a warning message is displayed.
Print some information about the solution process and the solution itself.
FE energy error.
Here, the exceptions are caught and some output is generated. Normally, this should not happen. We return a value different from 0 in case an exception is caught.
Sane return value.
Results
Output of the program:
@section complete Complete Source Code @author Philipp Frauenfelder, 2004