Chapter 2
Getting Started

To illustrate with an example, let us explain how FreeFem++ solves the Poisson’s equation: for a given function f(x,y), find a function u(x,y) satisfying

- Δu (x,y) =   f(x,y)  for all (x,y) ∈ Ω,                     (2.1)
    u(x,y) =   0   fo r all (x,y) on ∂Ω, .                   (2.2)
Here Ω is the boundary of the bounded open set Ω ℝ2 and Δu = ∂2u
∂x2- + ∂2u
∂y2-.

The following is a FreeFem++ program which computes u when f(x,y) = xy and Ω is the unit disk. The boundary C = Ω is

C  = {(x,y)|x = co s(t),y = sin(t),0 ≤ t ≤ 2 π}

Note that in FreeFem++ the domain Ω is assumed to described by its boundary that is on the left side of its boundary oriented by the parameter. As illustrated in Fig. 2.2, we can see the isovalue of u by using plot (see line 13 below).


PIC

Figure 2.1: mesh Th by build(C(50))

PIC

Figure 2.2: isovalue by plot(u)


Example 2.1  


     // defining the boundary
 1: border C(t=0,2⋆pi){x=cos(t); y=sin(t);}
     // the triangulated domain Th is on the left side of its boundary
 2: mesh Th = buildmesh (C(50));
     // the finite element space defined over Th is called here Vh
 3; fespace Vh(Th,P1);
 4: Vh u,v;   // defines u and v as piecewise-P1 continuous functions
 5: func f= x⋆y;   // definition of a called f function
 6: real cpu=clock();  // get the clock in second
 7: solve Poisson(u,v,solver=LU) =   // defines the PDE
 8:    int2d(Th)(dx(u)⋆dx(v) + dy(u)⋆dy(v))    // bilinear part
 9:    - int2d(Th)( f⋆v)           // right hand side
 10:    + on(C,u=0)  ;   // Dirichlet boundary condition
 11: plot(u);
 12: cout << " CPU time = " << clock()-cpu << endl;

Note that the qualifier solver=LU is not required and by default a multi-frontal LU would have been used. Note also that the lines containing clock are equally not required. Finally note how close to the mathematics FreeFem++ input language is. Line 8 and 9 correspond to the mathematical variational equation

∫                        ∫
   (∂u-∂v-+ ∂u-∂v)dxdy =     fvdxdy
 Th ∂x ∂x   ∂y ∂y         Th

for all v which are in the finite element space Vh and zero on the boundary C.

Exercise : Change P1 into P2 and run the program.

2.0.1 FEM by FreeFem++ : how does it work?

This first example shows how FreeFem++ executes with no effort all the usual steps required by the finite element method (FEM). Let us go through them one by one.

1st line: the boundary Γ is described analytically by a parametric equation for x and for y. When Γ = j=0JΓ j then each curve Γj, must be specified and crossings of Γj are not allowed except at end points .

The keyword “label” can be added to define a group of boundaries for later use (boundary conditions for instance). Hence the circle could also have been described as two half circle with the same label:


border Gamma1(t=0,pi)   {x=cos(t); y=sin(t); label=C}
border Gamma2(t=pi,2⋆pi){x=cos(t); y=sin(t); label=C}
 

Boundaries can be referred to either by name ( Gamma1 for example) or by label ( C here) or even by its internal number here 1 for the first half circle and 2 for the second (more examples are in Section 5.8).

2nd line: the triangulation Th of Ω is automatically generated by buildmesh(C(50)) using 50 points on C as in Fig. 2.1.

The domain is assumed to be on the left side of the boundary which is implicitly oriented by the parametrization. So an elliptic hole can be added by


 border C(t=2⋆pi,0){x=0.1+0.3⋆cos(t); y=0.5⋆sin(t);}
 

If by mistake one had written


 border C(t=0,2⋆pi){x=0.1+0.3⋆cos(t); y=0.5⋆sin(t);}
 

then the inside of the ellipse would be triangulated as well as the outside.

Automatic mesh generation is based on the Delaunay-Voronoi algorithm. Refinement of the mesh are done by increasing the number of points on Γ, for example, buildmesh(C(100)), because inner vertices are determined by the density of points on the boundary. Mesh adaptation can be performed also against a given function f by calling adaptmesh(Th,f).

Now the name Th (Th in FreeFem++ ) refers to the family {Tk}k=1,⋅⋅⋅,nt of triangles shown in figure 2.1. Traditionally h refers to the mesh size, nt to the number of triangles in Th and nv to the number of vertices, but it is seldom that we will have to use them explicitly. If Ω is not a polygonal domain, a “skin” remains between the exact domain Ω and its approximation Ωh = k=1ntT k. However, we notice that all corners of Γh = Ωh are on Γ.

3rd line: A finite element space is, usually, a space of polynomial functions on elements, triangles here only, with certain matching properties at edges, vertices etc. Here fespace Vh(Th,P1) defines Vh to be the space of continuous functions which are affine in x,y on each triangle of Th. As it is a linear vector space of finite dimension, basis can be found. The canonical basis is made of functions, called the hat functions ϕk which are continuous piecewise affine and are equal to 1 on one vertex and 0 on all others. A typical hat function is shown on figure 2.4 1. Then

            (       |                                        )
            ||{       ||         M∑                              ||}
Vh(Th, P1) = || w(x,y)|||w(x,y) =   wkϕk(x,y),wk are real num bers||
            (       |         k=1                            )
(2.3)

where M is the dimension of Vh, i.e. the number of vertices. The wk are called the degree of freedom of w and M the number of the degree of freedom.


PIC

Figure 2.3: mesh Th
PIC
Figure 2.4: Graph of ϕ1 (left hand side) and ϕ6


It is said also that the nodes of this finite element method are the vertices.
Currently FreeFem++ implements the following elements , (see section 6 for the full description)

P0 piecewise constant,

P1 continuous piecewise linear,

P2 continuous piecewise quadratic,

RT0 Raviart-Thomas piecewise constant,

P1nc piecewise linear non-conforming,

P1dc piecewise linear discontinuous,

P2dc piecewise quadratic discontinuous,

P1b piecewise linear continuous plus bubble,

P2b piecewise quadratic continuous plus bubble.

...

To get the full list do in a unix terminal, in directory examples++-tutorial do


FreeFem++ dumptable.edp
grep TypeOfFE lestables

The user can add other elements fairly easily if required.

Step3: Setting the problem
4th line: Vh u,v declares that u and v are approximated as above, namely

                 M -1
                 ∑
u(x, y) ≃ uh(x,y) =   ukϕk(x,y)
                  k=0
(2.4)

5th line: the right hand side f is defined analytically using the keyword func.

7th–9th lines: defines the bilinear form of equation (2.1) and its Dirichlet boundary conditions (2.2).
This variational formulation is derived by multiplying (2.1) by v(x,y) and integrating the result over Ω:

  ∫            ∫
-    vΔudx dy =    vfdxdy
   Ω             Ω
Then, by Green’s formula, the problem is converted into finding u such that
a(u,v) - ℓ(f,v)∫ = 0   ∀v satisfying v = 0∫on ∂Ω.                (2.5)

w ith a(u,v) =   ∇u ⋅ ∇v dxdy,  ℓ(f ,v) =   fvdx dy              (2.6)
               Ω                         Ω
In FreeFem++ the problem Poisson can be declared only (see below) for future use or declared and solved at the same time in which case

Vh u,v; solve Poisson(u,v) =

and (2.5) is written with dx(u) = ∂u∕∂x, dy(u) = ∂u∕∂y and

∫
   ∇u  ⋅ ∇vdx dy -→ in t2d(T h) ( dx (u) ⋆dx( v) +  dy(u )⋆dy (v)  )
∫ Ω

   f vdxdy -→  int2d (Th )(  f⋆ v )     (notice here, u is unused )
  Ω
In FreeFem++ there is no need to distinguish the bilinear form a from the linear form , as long as the terms are inside different integrals, FreeFem++ find out which one is the bilinear form by checking where both terms u and v are present.

The other way is to define the problem and then we solve it, write :


  Vh u,v; problem Poisson(u,v) =
   ...
  Poisson;  // the problem now is solved here

Step4: Solution and visualization

6th line: The current time in seconds is stored into the real-valued variable cpu.

7th line The problem is solved.

11th line: The visualization is done as illustrated in Fig. 2.2 (see Section 7.1 for zoom, postscript and other commands).

12th line: The computing time (not counting graphics) is written on the console Notice the C++-like syntax; the user needs not study C++ for using FreeFem++ , but it helps to guess what is allowed in the language.

Access to matrices and vectors
Internally FreeFem++ will solve a linear system of the type

M-1                                          ∫
∑
    Aijuj - Fi = 0,  i = 0,⋅⋅⋅,M - 1;    Fi =    fϕidxdy               (2.7)
j=0                                           Ω
which is found by using (2.4) and replacing v by ϕi in (2.5). And the Dirichlet conditions are implemented by penalty, namely by setting Aii = 1030 and Fi = 1030 *0 if i is a boundary degree of freedom. Note, that the number 1030 is called tgv (très grande valeur) and it is generally possible to change this value , see the index item solve!tgv=.

The matrix A = (Aij) is called stiffness matrix .

If the user wants to access A directly he can do so by using (see section 6.10 page 282 for details)


varf a(u,v) = int2d(Th)( dx(u)⋆dx(v) + dy(u)⋆dy(v))
               + on(C,u=0) ;
matrix A=a(Vh,Vh);   // stiffness matrix,

The vector F in (2.7) can also be constructed manually


varf l(unused,v) = int2d(Th)(f⋆v)+on(C,u=0);
Vh F;  F[] = l(0,Vh);  // F[] is the vector associated to the function F

The problem can then be solved at the level of algebra by


 u[]=A^-1⋆F[];  // u[] is the vector associated to the function u

Note 2.1 Here uand Fare finite element function, and u[]and F[]give the array of value associated ( u[](ui)i=0,,M-1 and F[](Fi)i=0,,M-1). So we have

         M∑-1                          M∑-1
u (x,y) =    u [][i]ϕi(x,y),    F(x,y) =    F [][i]ϕi(x,y)
         i=0                          i=0
where ϕi,i = 0...,,M -1 are the basis functions of Vhlike in equation (2.3), and M = Vh.ndof is the number of degree of freedom (i.e. the dimension of the space Vh).

The linear system (2.7) is solved by UMFPACK unless another option is mentioned specifically as in

Vh u,v; problem Poisson(u,v,solver=CG) = int2d(...

meaning that Poisson is declared only here and when it is called (by simply writing Poisson; ) then (2.7) will be solved by the Conjugate Gradient method.

2.0.2 Some Features of FreeFem++

The language of FreeFem++ is typed, polymorphic and reentrant with macro generation (see 9.12). Every variable must be typed and declared in a statement each statement separated from the next by a semicolon ”;”. The syntax is that of C++ by default augmented with something that is more akin to TEX. For the specialist, one key guideline is that FreeFem++ rarely generates an internal finite element array; this was adopted for speed and consequently FreeFem++ could be hard to beat in terms of execution speed, except for the time lost in the interpretation of the language (which can be reduced by a systematic usage of varf and matrices instead of problem.

2.1 The Development Cycle: Edit–Run/Visualize–Revise

An integrated environment is provided with FreeFem++ written by A. Le Hyaric; Many examples and tutorials are also given along with this documentation and it is best to study them and learn by example. Explanations for some of these examples are given in this book in the next chapter. If you are a FEM beginner, you also may want to read a book on variational formulations.

The development cycle will have the following steps:

Modeling:
From strong forms of PDE to weak forms, one must know the variational formulation to use FreeFem++ ; one should also have an eye on the reusability of the variational formulation so as to keep the same internal matrices; a typical example is the time dependent heat equation with an implicit time scheme: the internal matrix can be factorized only once and FreeFem++ can be taught to do so.
Programming:
Write the code in FreeFem++ language using a text editor such as the one provided in the integrated environment.
Run:
Run the code (e.g. written in file mycode.edp). If not from the integrated environment it can be done at the console level by

% FreeFem++ mycode.edp

Note, the name of the command FreeFem++ may depend on your installation.

Visualization:
Use the keyword plot to display functions while FreeFem++ is running. Use the plot-parameter wait=1 to stop the program to have time to see the plot. Use the plot-parameter ps="toto.eps" to generate a postscript file to archive the results.
Debugging:
A global variable ”debug” (for example) can help as in wait=true to wait=false.


bool debug = true;
border a(t=0,2⋆pi){ x=cos(t); y=sin(t);label=1;}
border b(t=0,2⋆pi){ x=0.8+0.3⋆cos(t); y=0.3⋆sin(t);label=2;}
plot(a(50)+b(-30),wait=debug);  // plot the borders to see the intersection
// (so change (0.8 in 0.3 in b) then needs a mouse click
mesh Th = buildmesh(a(50)+b(-30));
plot(Th,wait=debug);   // plot Th then needs a mouse click
fespace Vh(Th,P2);
Vh f = sin(pi⋆x)⋆cos(pi⋆y);
plot(f,wait=debug);   // plot the function f
Vh g = sin(pi⋆x + cos(pi⋆y));
plot(g,wait=debug);   // plot the function g

Changing debug to false will make the plots flow continuously; drinking coffee and watching the flow of graph on the screen can then become a pleasant experience.

Error messages are displayed in the console window. They are not always very explicit because of the template structure of the C++ code, sorry! Nevertheless they are displayed at the right place. For example, if you forget parenthesis as in


bool debug = true;
mesh Th = square(10,10;
plot(Th);

then you will get the following message from FreeFem++,


    2 : mesh Th = square(10,10;
 Error line number 2, in file bb.edp, before  token ;
parse error
  current line = 2
Compile error : parse error
        line number :2, ;
error Compile error : parse error
        line number :2, ;
 code = 1

If you use the same symbol twice as in


real aaa =1;
real aaa;

then you will get the message


    2 : real aaa; The identifier aaa exists
          the existing type is <Pd>
          the new  type is <Pd>

If you find that the program isn’t doing what you want you may also use cout to display in text format on the console window the value of variables.

The following example works:


...;
fespace Vh...; Vh u;...
cout<<u;...
matrix A=a(Vh,Vh);...
cout<<A;

Another trick is to comment in and out by using the“ //” as in C++. For example


real aaa =1;
// real aaa;