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 convergence of the relative energy error and the maximal nodal error are shown together with plots of the solution.
System include files. The first is for input and output operations, the second one declares unique_ptr
and the third is for stringstream
.
Concepts include files.
With this using directive, we do not need to prepend concepts:
: before Real
everytime.
The mesh is taken from the meshes tutorial: Line
. Here, the relevant code is included:
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 (gauss_p
).
The constructor takes as argument the number of quadrature points gauss_p
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 point q
[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). e->shapefct(i, q[k])
is the value of the shape function i
at the point q
[k]. As you can see in linearFEM::Line, it takes values in (-1,1). 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 function pow
from the standard math library computes ab. The only parameters of the routine fem
are the maximal level of refinement and the number of quadrature points. The additional argument output
is used to for the gathering the output.
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.
Start the loop over the different levels.
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.
Gather the results of the computation. Using output
, all data can be displayed and/or saved at the end.
The main program just does some comand line parameter handling calls the main function fem
. It also features the main try
...catch brace to catch exceptions thrown by the program.
The inputParser
takes care of all input parameters and the results which are gathered (see input output tutorial). At the end, they can be printed.
Default values for the variable which could be set by a command line argument: the maximal level of refinement and the number of points for the Gauss quadrature.
Parse the command line options using the standard library method getopt
. Have a look at getopt(3)
to find out more about this function.
Prepare the results part to gather the results and compute the exact FE energy and the element size.
The values for the level and the quadrature points are fixed now, call fem
.
Print the results.
In order to be able to plot the numbers computed in the fem
routine, it is convenient to have them in table suitable for Gnuplot.
Every array from output
which should be present in the table has to be added to table
.
Then, a file is created with high numerical output precision and the table is printed to the file. The file is closed automatically. See the results section on how to process this file.
Here, the exceptions are caught and some output is generated. Normally, this should not happen.
Sane return value.
Output of the program:
Error diagrams: <img src="linearFEM1d-conv.png" width="640" height="480" alt="error diagram"> The diagram above was created using the following commands: @code
set logscale x set logscale y set grid set grid mxtics set grid mytics set key top left set xlabel 'h' set ylabel 'error' plot 'linearFEM1d-1.gnuplot' using 3:2 title 'relative energy error, gauss_p = 1', \ '' using 3:4 title 'max nodal error, gauss_p = 1', \ 'linearFEM1d-2.gnuplot' using 3:2 title 'relative energy error, gauss_p = 2', \ '' using 3:4 title 'max nodal error, gauss_p = 2'
Exact and numerical solution: The diagram above was created using the following commands:
@section complete Complete Source Code @author Philipp Frauenfelder, 2004