This page was generated from unit-1.3-dirichlet/dirichlet.ipynb.

1.3 Dirichlet boundary conditions

This tutorial goes in depth into the mechanisms required to solve the Dirichlet problem

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

with a nonzero Dirichlet boundary condition

\[u|_{\Gamma_D} = g \quad \text{ on a boundary part } \Gamma_D \subset \partial\Omega.\]

The same mechanisms are used in solving boundary value problems involving operators other than the Laplacian.

You will see how to perform these tasks in NGSolve: - extend Dirichlet data from boundary parts, - convert boundary data into a volume source, - reduce inhomogeneous Dirichlet case to the homogeneous case, and - perform all these tasks automatically within a utility.

If you are just interested in the automatic utility, then please skip to the last part of the tutorial.

Spaces with Dirichlet conditions on part of the boundary

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

mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
mesh.GetBoundaries()
[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:

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

Compare this space with the one without the dirichlet flag:

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

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 if and only if dof is a free degree of freedom.

[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: 11111111111111111111111111111
free dofs of fes:
 0: 00001111000011110000111111111111111111101010110111
50: 11111111101101101111111111111110101110111111111111
100: 11111111111111111111111111111
  • 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.

Extension of boundary data

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\). These are the issues to consider in this approach:

  • How to define an extension \(u_D\) in the finite element space?

  • How to form and solve the system for \(u_0\)?

Let us address the first in the following example.

Suppose we are given that

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

We interpolate \(g\) on the boundary of the domain and extend it to zero on the elements not having an intersection with \(\Gamma_D\).

[6]:
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.

Thus, gfu now contains the extension \(u_D\). Next, we turn to set up the system for \(u_0\).

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.

[7]:
u, v = fes.TnT()

a = BilinearForm(fes, symmetric=True)
a += grad(u)*grad(v)*dx
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 may not be the case in practice) for ease of presentation. The first row gives

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

Since we have already constructed \(u_D\), we need to perform these next steps:

  • 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\).

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.\]
[8]:
f = LinearForm(fes)
f += 1*v*dx
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:

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

The automatic utility BVP

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

[10]:
gfu.Set(g, BND)
c = Preconditioner(a,"local")   #<- Jacobi preconditioner
#c = Preconditioner(a,"direct") #<- sparse direct solver
c.Update()
solvers.BVP(bf=a, lf=f, gf=gfu, pre=c)
Redraw()
iteration 0 error = 1.324197630852001
iteration 1 error = 0.7482753330831278
iteration 2 error = 0.6935365444529648
iteration 3 error = 0.5530815561584947
iteration 4 error = 0.453678306354409
iteration 5 error = 0.2731551456979582
iteration 6 error = 0.1951543065747879
iteration 7 error = 0.11089771808731903
iteration 8 error = 0.05456942143195438
iteration 9 error = 0.038487308353212846
iteration 10 error = 0.025822078243051938
iteration 11 error = 0.01821275053415462
iteration 12 error = 0.010781043321490713
iteration 13 error = 0.006592401333924963
iteration 14 error = 0.0038730411130962522
iteration 15 error = 0.0019102998892333533
iteration 16 error = 0.0008755946467775109
iteration 17 error = 0.0004089159296759819
iteration 18 error = 0.0002564404510203854
iteration 19 error = 0.00012480800222510315
iteration 20 error = 6.769009877447335e-05
iteration 21 error = 3.306871209493643e-05
iteration 22 error = 1.975978745569913e-05
iteration 23 error = 1.2224226967564665e-05
iteration 24 error = 6.288984116832294e-06
iteration 25 error = 3.337757839771981e-06
iteration 26 error = 1.570161093800567e-06
iteration 27 error = 8.680106877075956e-07
iteration 28 error = 4.271339256525377e-07
iteration 29 error = 1.9287597213516924e-07
iteration 30 error = 9.76863469480573e-08
iteration 31 error = 4.61169599076755e-08
iteration 32 error = 2.491303943496592e-08
iteration 33 error = 1.2560326821667961e-08