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$.
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$.
import netgen.gui
%gui tk
from ngsolve import *
from netgen.geom2d import unit_square
mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
mesh.GetBoundaries()
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:
fes = H1(mesh, order=2, dirichlet="left|right")
Compare this space with the one without the dirichlet
flag:
fs2 = H1(mesh, order=2)
fes.ndof, fs2.ndof # total number of unknowns
Thus, the dirichlet
flag did not change ndof
. In NGSolve the unknowns are split into two groups:
The facility FreeDofs
gives a BitArray
such that FreeDofs[dof] is True iff dof is a free degree of freedom.
print("free dofs of fs2 without \"dirichlet\" flag:\n",
fs2.FreeDofs())
print("free dofs of fes:\n", fes.FreeDofs())
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
.
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.
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
$$
\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)
% \begin{array}{l} \leftarrow \text{ solved for } u_{0,F} \\ \leftarrow \text{ replaced with } u_{0,D}
% \quad \begin{array}{l} \leftarrow \text{ solved for } u_{0,F} \\ \leftarrow \text{ replaced with } u_{0,D} = g_D\end{array}
$$
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:
Suppose we are given that $$ g = \sin(y) \qquad \text{ on } \; \Gamma_D. $$
g = sin(y)
We interpolate $g$ on the boundary of the domain and extend trivially to $\Omega$ as follows:
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
.
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. $$
f = LinearForm(fes)
f += SymbolicLFI(1*v)
f.Assemble()
r = f.vec.CreateVector()
r.data = f.vec - a.mat * gfu.vec
The implementation of $$ u = u_D + \left( \begin{array}{cc} A_{FF}^{-1} & 0 \\ 0 & 0 \end{array} \right) r $$ by sparse solvers is achieved by the following:
gfu.vec.data += \
a.mat.Inverse(freedofs=fes.FreeDofs()) * r
Redraw()
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.
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()