Name

A step by step basic example — This example shows the basic usage of getfem.

Description

This example shows the basic usage of getfem, on the über-canonical problem above all others: solving the Laplacian, Δu+f=0 on a square, with the Dirichlet condition u=g(x) on the domain boundary.

The first step is to create a mesh. Since getfem++ does not come with its own mesher, one has to rely on an external mesher (see gf_mesh('import')), or use very simple meshes. For this example, we just consider a regular mesh whose nodes are {xi=0...10,j=0..10=(i/10,j/10)}.

 
// creation of a simple cartesian mesh
--> m = gf_mesh('cartesian',[0:.1:1],[0:.1:1])
m =
     id: 0
    cid: 0
 

If you try to look at the value of m, you'll notice that it appears to be a structure containing two integers. The first one is its identifier, the second one is its class-id, i.e. an identifier of its type. This small structure is just an "handle" or "descriptor" to the real object, which is stored in the getfem++ memory and cannot be represented via Matlab data structures. Anyway, you can still inspect the getfem++ objects via the command gf_workspace('stats').

Now we can try to have a look at the mesh, with its vertices numbering and the convexes numbering:

 
// we enable vertices and convexes labels
--> gf_plot_mesh(m, 'vertices', 'on', 'convexes', 'on');
 

As you can see, the mesh is regular, and the numbering of its nodes and convexes is also regular (this is guaranteed for cartesian meshes, but do not hope a similar numbering for the degrees of freedom).

The next step is to create a mesh_fem object. This one links a mesh with a set of FEM.

 
--> mf = gf_mesh_fem(m,1);    // create a mesh_fem of for a field of dimension 1 (i.e. a scalar field)
--> gf_mesh_fem_set(mf,'fem',gf_fem('FEM_QK(2,2)'));
 

The first instruction builds a new mesh_fem object, the second argument specifies that this object will be used to interpolate scalar fields (since the unknown is a scalar field). The second instruction assigns the Q2 FEM to every convex (each basis function is a polynomial of degree 4, remember that Pk => polynomials of degree k, while Qk => polynomials of degree 2k). As Q2 is a polynomial FEM, you can view the expression of its basis functions on the reference convex:

 
--> gf_fem_get(gf_fem('FEM_QK(2,2)'), 'poly_str')
ans =
    '1 - 3*x - 3*y + 2*x^2 + 9*x*y + 2*y^2 - 6*x^2*y - 6*x*y^2 + 4*x^2*y^2'
    '4*x - 4*x^2 - 12*x*y + 12*x^2*y + 8*x*y^2 - 8*x^2*y^2'
    '-x + 2*x^2 + 3*x*y - 6*x^2*y - 2*x*y^2 + 4*x^2*y^2'
    '4*y - 12*x*y - 4*y^2 + 8*x^2*y + 12*x*y^2 - 8*x^2*y^2'
    '16*x*y - 16*x^2*y - 16*x*y^2 + 16*x^2*y^2'
    '-4*x*y + 8*x^2*y + 4*x*y^2 - 8*x^2*y^2'
    '-y + 3*x*y + 2*y^2 - 2*x^2*y - 6*x*y^2 + 4*x^2*y^2'
    '-4*x*y + 4*x^2*y + 8*x*y^2 - 8*x^2*y^2'
    'x*y - 2*x^2*y - 2*x*y^2 + 4*x^2*y^2'
 

Now, in order to perform numerical integrations on mf, we need to build a mesh_im object:

 
// assign the same integration method on all convexes 
--> mim=gf_mesh_im(m, gf_integ('IM_EXACT_PARALLELEPIPED(2)'));
  

The integration method will be used to compute the various integrals on each element: here we choose to perform exact computations (no quadrature formula), which is possible since the geometric transformation of these convexes from the reference convex is linear (this is true for all simplices, and this is also true for the parallelepipeds of our regular mesh, but it is not true for general quadrangles), and the chosen FEM is polynomial. Hence it is possible to analytically integrate every basis function/product of basis functions/gradients/etc. There are many alternative FEM methods and integration methods (see the description of finite element and integration methods).

Note however that in the general case, approximate integration methods are a better choice than exact integration methods.

Now we have to find the "boundary" of the domain, in order to set a condition. A mesh object has the ability to store some sets of convexes and convex faces. These sets (called "regions") are accessed via an integer #id:

 
--> border = gf_mesh_get(m,'outer faces');
--> gf_mesh_set(m, 'region', 42, border); // create the region #42
--> gf_plot_mesh(m, 'regions', [42]); // the boundary edges appears in red
  

Here we find the faces of the convexes which are on the boundary of the mesh (i.e. the faces which are not shared by two convexes). remark: we could have used gf_mesh_get(m, 'OuTEr_faCes') , as the interface is case-insensitive, and whitespaces can be replaced by underscores. The array border has two rows, on the first row is a convex number, on the second row is a face number (which is local to the convex, there is no global numbering of faces). Then this set of faces is assigned to the region number 42.

At this point, we just have to stack some model bricks and run the solver to get the solution! The "model bricks" are created with the gf_mdbrick constructor. A model brick is basically an object which modifies a global linear system (tangent matrix for non-linear problems) and its associated right hand side. Typical modifications are insertion of the stiffness matrix for the problem considered (linear elasticity, laplacian, etc), handling of a set of contraints, Dirichlet condition, addition of a source term to the right hand side etc. The global tangent matrix and its right hand side are stored in a "model state" structure, created with the gf_mdstate constructor.

Let us build a problem with an easy solution: , then we have (the FEM won't be able to catch the exact solution since we use a Q2 method).

We start with a "generic elliptic" brick, which handles problems, where A can be a scalar field, a matrix field, or an order 4 tensor field. By default, A=1.

 
--> b0=gf_mdbrick('generic elliptic',mim,mf)
  

Each brick embeds a number of parameter fields. In the case of the generic elliptic brick, there is only one parameter field, the A(x) coefficient in . It is possible to view the list of parameters of the brick with

 
--> gf_mdbrick_get(b0, 'param list')
ans =

    'A'
--> gf_mdbrick_get(b0, 'param', 'A')

ans =

     1
 

Next we add a Dirichlet condition on the domain boundary:

 
--> b1=gf_mdbrick('dirichlet',b0,42,mf,'penalized')
  

Here the number 42 is the region number to which the dirichlet condition is applied. The 'penalized' says that the Dirichlet condition should be imposed via a penalization technique. Other ways are possible (augmented system, direct elimination). A mesh fem argument is also required, as the Dirichlet condition u=r is imposed in a weak form

where v is taken in the space of multipliers given by here by mf.

By default, the Dirichlet brick imposes u=0 on the specified boundary. We change this to :

 
--> R=gf_mesh_fem_get(mf, 'eval', list('(x-.5).^2 + (y-.5).^2 + x/5 - y/3'));
--> gf_mdbrick_set(b1, 'param', 'R', mf, R); 
  

Remark: the polynomial expression was interpolated on mf. It is possible only if mf is of Lagrange type. In this first example we use the same mesh fem for the unknown and for the data such as R, but in the general case, mf won't be Lagrangian and another (Lagrangian) mesh fem will be used for the description of Dirichlet conditions, source terms etc.

A "model state" variable is created, and the solver is launched:

 
--> mds=gf_mdstate('real')
--> gf_mdbrick_get(b1, 'solve', mds)
 

The model state now contains the solution (as well as other things, such as the linear system which was solved). It is extracted, a display into a matlab figure.

 
--> U=gf_mdstate_get(mds, 'state');
--> gf_plot(mf, U, 'mesh','on');
 

See Also

gf_workspace, gf_mesh, gf_fem, gf_plot

Authors

Y. Collette