1.3 Dirichlet boundary conditions

This section shows how to solve the Dirichlet problem

\[-\Delta u = f \quad \text{ in } \Omega\]

with a nonhomogeneous Dirichlet boundary condition \(u|_{\Gamma_D} = g\) on a boundary part \(\Gamma_D\).

The extension technique

We use the standard technique of reducing a problem with essential non-homogeneous boundary conditions to one with homogeneous boundary condition using an extension. The solution \(u\) in \(H^1\) satisfies

\[u|_{\Gamma_D} = g\]

and

\[\int_\Omega \nabla u \cdot \nabla v_0 = \int_\Omega f v_0\]

for all \(v_0\) in \(\in H_{0,D}^1 = \{ v \in H^1: v|_{\Gamma_D} = 0\}\). Split the solution

\[u = u_0 + u_D\]

where \(u_D\) is an extension of \(g\) into \(\Omega\). Then we only need to find \(u_0\) in \(H^1_{0,D}\) satisfying the homogeneous Dirichlet problem

\[\int_\Omega \nabla u_0 \cdot \nabla v_0 = \int_\Omega f v_0 - \int_\Omega \nabla u_D \cdot \nabla v_0\]

for all \(v_0\) in \(H_{0,D}^1\).

Issues to consider

  • How to define an extension \(u_D\) in the finite element space?
  • How to form and solve the system for \(u_0\)?

Finite element spaces with Dirichlet conditions

In [1]:
import netgen.gui
%gui tk
from ngsolve import *
from netgen.geom2d import unit_square

mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
mesh.GetBoundaries()
Out[1]:
('bottom', 'right', 'top', 'left')

The unit_square has its boundaries marked as left, right, top and bottom. Suppose we want non-homogeneous Dirichlet boundary conditions on

\[\Gamma_D = \Gamma_{left} \cup \Gamma_{right}.\]

Then, we set the space as follows:

In [2]:
fes = H1(mesh, order=2, dirichlet="left|right")

Compare this space with the one without the dirichlet flag:

In [3]:
fs2 = H1(mesh, order=2)
fes.ndof, fs2.ndof    # total number of unknowns
Out[3]:
(133, 133)

Thus, the dirichlet flag did not change ndof. In NGSolve the unknowns are split into two groups: * dirichlet dofs (or constrained dofs), * free dofs.

The facility FreeDofs gives a BitArray such that FreeDofs[dof] is True iff dof is a free degree of freedom.

In [4]:
print("free dofs of fs2 without \"dirichlet\" flag:\n",
      fs2.FreeDofs())
print("free dofs of fes:\n", fes.FreeDofs())
free dofs of fs2 without "dirichlet" flag:
 0: 11111111111111111111111111111111111111111111111111
50: 11111111111111111111111111111111111111111111111111
100: 111111111111111111111111111111111
free dofs of fes:
 0: 00001111000011110000111111111111111111110101011011
50: 11111111110110110111111111111111011011011111111111
100: 111111111111111111111111111111111
  • The space fs2 without dirichlet flag has only free dofs (no dirichlet dofs).
  • The other space fes has a few dofs that are marked as not free. These are the dofs that are located on the boundary regions we marked as dirichlet.

Forms and assembly

In NGSolve, bilinear and linear forms are defined independently of the dirichlet flags. Matrices and vectors are set up with respect to all unknowns (free or constrained) so they may be restricted to any group of unknowns later.

In [5]:
u = fes.TrialFunction()
v = fes.TestFunction()

a = BilinearForm(fes, symmetric=True)
a += SymbolicBFI(grad(u)*grad(v))
a.Assemble()

If \(A=\) a.mat is the matrix just assembled, then we want to solve for

\[A (u_0 + u_D) = f \quad \Rightarrow \quad A u_0 = f - A u_D\]

or

\[\begin{split}\left( \begin{array}{cc} A_{FF} & A_{FD} \\ A_{DF} & A_{DD} \end{array} \right) \left( \begin{array}{c} u_{0,F} \\ 0 \end{array} \right) = \left( \begin{array}{c} {f}_F \\ {f}_D \end{array} \right) - \left( \begin{array}{cc} A_{FF} & A_{FD} \\ A_{DF} & A_{DD} \end{array} \right) \left( \begin{array}{c} u_{D,F} \\ u_{D,D} \end{array} \right)\end{split}\]

where we have block partitioned using free dofs (\(F\)) and dirichlet dofs (\(D\)) as if they were numbered consecutively (which is typically not true). The first row gives

\[A_{FF} u_{0,F} = f_F - [A u_D]_F.\]

Hence, we need to:

  • Construct \(u_D\) from \(g\).
  • Set up the right hand side from \(f\) and \(u_D\).
  • Solve a linear system which involves only \(A_{FF}\).
  • Add solution: \(u = u_0 + u_D\).

Extending the boundary values

Suppose we are given that

\[g = \sin(y) \qquad \text{ on } \; \Gamma_D.\]
In [6]:
g = sin(y)

We interpolate \(g\) on the boundary of the domain and extend trivially to \(\Omega\) as follows:

In [7]:
gfu = GridFunction(fes)
gfu.Set(g, BND)
Draw(gfu)

The keyword BND tells Set that g need only be interpolated on those parts of the boundary that are marked dirichlet.

Solve for the free dofs

We need to assemble the right hand side of \(A_{FF} u_{0,F} = f_F - [A u_D]_F\), namely

\[r = f - A u_D.\]
In [8]:
f = LinearForm(fes)
f += SymbolicLFI(1*v)
f.Assemble()

r = f.vec.CreateVector()
r.data = f.vec - a.mat * gfu.vec

The implementation of

\[\begin{split}u = u_D + \left( \begin{array}{cc} A_{FF}^{-1} & 0 \\ 0 & 0 \end{array} \right) r\end{split}\]

by sparse solvers is achieved by the following:

In [9]:
gfu.vec.data += a.mat.Inverse(freedofs=fes.FreeDofs()) * r
Redraw()

Using BVP

NGSolve also provides a BVP facility, within which the above steps are performed automatically. You provide \(A\), \(f\), a grid function gfu with the boundary condition, and a preconditioner. Then BVP solves the problem with non-homogeneous Dirichlet boundary condition and overwrites gfu with the solution.

In [10]:
c = Preconditioner(a,"local")   #<- Jacobi preconditioner
#c = Preconditioner(a,"direct") #<- sparse direct solver
c.Update()
BVP(bf=a,lf=f,gf=gfu,pre=c).Do()
Redraw()