![]() |
Dune-Functions 2.11
|
The following example explains how to solve the Poisson problem with Lagrange finite elements using the Python interface of dune-functions. The full example file poisson-pq2.py can be found in the examples/ subdirectory.
Be warned that writing a quadrature loop for a stiffness matrix with dune-functions and Python is currently quite slow. This approach is therefore mainly useful for small problems, for educational purposes, or if you need complete control. For faster and more efficient ways to assemble finite element stiffness matrices, have a look at the dune-fufem and dune-fem modules.
Let
be a bounded open subset of
or
. The Poisson problem asks for a scalar function
that solves
for a given function
. To make the solution unique, it is additionally required that the solution
satisfies boundary conditions
The finite element method replaces the partial differential equation by a so-called weak formulation, and solves it in a space
of continuous, piecewise polynomial functions: The goal is now to find a function
still satisfying the boundary conditions, such that
where
is the space of finite element functions that are zero on the boundary. Writing this weak form in terms of a basis
of
, one obtains a sparse linear system of equations
which can be solved using standard algorithms from numerical linear algebra.
As every Python program, this one starts by importing required external modules. In this case, we need code from numpy and scipy for the sparse linear algebra, and from dune for grids, finite element bases, and quadrature rules.
The main program starts by creating a grid. For simplicity, the domain here is the unit square in two space dimensions, and the grid is a 32 by 32 quadrilateral grid:
The next line creates the basis of the finite element space. We use second-order Lagrange finite elements here, i.e., the space
consists of continuous functions that are quadratic in each variable on each element. The class functions.Lagrange encodes the nodal basis of that space:
Then, the Poisson problem is assembled. The source term
is specified in a lambda expression f, which here encodes the constant function
. The actual assembly then happens in the method assembleLaplaceMatrix, which we discuss below.
The assembleLaplaceMatrix method returns two objects, the stiffness matrix A and the assembled source term vector b.
The next code block manages the handling of the Dirichlet boundary conditions. In this example, the entire boundary is the Dirichet boundary. To find all degrees of freedom on the grid boundary, the method forEachBoundaryDOF visits each degree of freedom (DOF) associated to the grid boundary and calls the callback function markDOF with the global index of that degree of freedom. This callback then sets the corresponding entry in an array called isDirichlet to True.
Then, the function
from the problem formulation is defined. In this example the Dirichlet values are given by the closed-form expression
This closed form expression is evaluated at the boundary Lagrange nodes by the method basis.interpolate.
After this, the information about the boundary conditions is integrated into the linear system of equations. The code loops over all nonzero entries of the matrix. If an entry belongs to a matrix row that corresponds to a boundary degree of freedom, it is set to 1 if it is the diagonal entry, and to 0 otherwise. This pins the degree of freedom to the corresponding boundary value. This part of the code is completely independent from Dune.
At this point, the linear system of equations has been completely assembled, and can now be solved. For this, the code simply calls a sparse direct solver from the Python scipy package. This is not part of Dune, either:
The solver leaves the algebraic solution vector in the vector variable x. The solution is then written to a file, in the VTK file format. To do so, the solution vector must be reinterpreted as a finite element functions again. The asFunction method does this:
This will write a file called poisson-pq2-result.vtu, which contains the grid and the solution function
. It can be opened, for example, with the ParaView visualization program.
The assembler is implemented in two methods above the main method of the program.
The central function is assembleLaplaceMatrix, which computes the stiffness matrix
with entries
and the load vector
with entries
As described in many text books, this assembly is implemented as a loop over the grid elements. On each element, all relevant integrals are computed and stored in a matrix localA, the element stiffness matrix. Since the basis functions
are zero on most elements by construction, this element stiffness matrix is much smaller than the global stiffness matrix A, but it is dense. Computation of localA happens by calling yet another method localAssembler (discussed below). Afterwards, the entries of the element stiffness matrix are added to the appropriate places of the global stiffness matrix A. The basis object provides the correspondence between the entries of the element stiffness matrices and the entries of the global stiffness matrix.
Note that A is sparse (i.e., it contains mainly zeros). This property has to be exploited. Therefore, a lil_matrix object from the Python scipy package is used to build it. After assembly has been completed, the matrix is converted to Compressed Sparse Row (CSR) format, because that is more efficient for the solver.
The element stiffness matrices are computed by a free function called localAssembler. For a given element
, it computes the element stiffness matrix
and the element load vector
with entries
Both are dense objects. In this example, where the domain is two-dimensional, the grid consists of quadrilaterals, and the finite element space consists of second-order Lagrange elements, each
is a
matrix. The element and the set of basis functions on the element are passed to the method via a localView object. See the Dune::Functions::DefaultLocalView page for documentation of the corresponding C++ interface.
The computation of the integrals happens via numerical quadrature, and the quadrature rule is defined on a reference element
. The integral therefore needs to be transformed to the reference element, using a mapping
:
The integrals over
are then computed by a quadrature rule.