This page was generated from unit-1.4-staticcond/staticcond.ipynb.

1.4 Static Condensation

Static condensation refers to the process of eliminating unknowns that are internal to elements from the global linear system. They are useful in standard methods and critical in methods like the HDG method. NGSolve automates static condensation across a variety of methods via a classification of degrees of freedom.

In this tutorial, you will learn

  • how to turn on static condensation (using a keyword argument condense),

  • how to get the full solution from the condensed system using attributes called

    • harmonic_extension_trans

    • harmonic_extension

    • inner_solve,

  • their relationship to a Schur complement factorization,
  • types of degrees of freedom such as COUPLING_TYPE.LOCAL_DOF, and

  • an automatic utility for solving via static condensation.

[1]:
from ngsolve import *
from ngsolve.webgui import Draw

mesh = Mesh(unit_square.GenerateMesh(maxh=0.4))

fes = H1(mesh, order=4, dirichlet='bottom|right')
u, v = fes.TnT()

Asking BilinearForm to condense

Keyword arguments condense, or the deprecated keyword argument eliminate_internal, are synonymous.

[2]:
a = BilinearForm(grad(u) * grad(v) * dx, condense=True).Assemble()
f = LinearForm(1 * v * dx).Assemble()
u = GridFunction(fes)
  • The assembled matrix \(A=\) a.mat can be block partitioned into

\[\begin{split}A = \left( \begin{array}{cc} A_{LL} & A_{LE} \\ A_{EL} & A_{EE} \end{array} \right)\end{split}\]

where \(L\) denotes the set of local (or internal) degrees of freedom and \(E\) denotes the set of interface degrees of freedom.

  • In our current example \(E\) consists of edge and vertex dofs while \(L\) consists of triangle dofs. (Note that in practice, \(L\) and \(E\) may not be ordered contiguously and \(L\) need not appear before \(E\), but such details are immaterial for our discussion here.)

  • The condensed system is known as the Schur complement:

\[S = A_{EE} - A_{EL} A_{LL}^{-1} A_{LE}.\]
  • When condense is set to True in a, the statement a.Assemble actually assembles \(S\). Thus

    • a.mat \(= \left(\begin{array}{cc}0 & 0 \\0 & S\end{array}\right),\) and

    • a.mat.Inverse(freedofs=fes.FreeDofs(coupling=True)) \(= \left(\begin{array}{cc} 0 & 0 \\ 0 & S^{-1} \end{array}\right).\) Let us make this matrix:

[3]:
invS = a.mat.Inverse(freedofs=fes.FreeDofs(coupling=True))

Note that instead of the usual fes.FreeDofs(), we have used fes.FreeDofs(coupling=True) or simply fes.FreeDofs(True) to specify that only the degrees of freedom that are not local and not Dirichlet should participate in the inverse computations. (The underlying assumption is that Dirichlet dofs cannot be local dofs.)

A factorization of the inverse

  • To solve

    \[\begin{split} \left(\begin{array}{cc} A_{LL} & A_{LE}\\ A_{EL} & A_{EE} \end{array}\right) \left(\begin{array}{c} u_L \\ u_E \end{array}\right)= \left(\begin{array}{c} f_L \\ f_E \end{array}\right)\end{split}\]

    we use a factorization of \(A^{-1}\) that uses \(S^{-1}\). Namely, we use the following identity:

\[\begin{split}\left(\begin{array}{c} u_L \\ u_E \end{array}\right) = \left(\begin{array}{cc} I & -A_{LL}^{-1} A_{LE} \\ 0 & I \end{array}\right) \left(\begin{array}{cc} A_{LL}^{-1} & 0 \\ 0 & S^{-1} \end{array}\right) \left(\begin{array}{cc} I & 0 \\ -A_{EL} A_{LL}^{-1} & I \end{array}\right) \left(\begin{array}{c} f_L \\ f_E \end{array}\right)\end{split}\]
  • NGSolve provides

    • a.harmonic_extension_trans \(= \left(\begin{array}{cc} 0 & 0 \\ -A_{EL} A_{LL}^{-1} & 0 \end{array}\right)\)

    • a.harmonic_extension \(= \left(\begin{array}{cc} 0 & -A_{LL}^{-1} A_{LE} \\ 0 & 0 \end{array}\right)\)

    • a.inner_solve \(=\left(\begin{array}{cc} A_{LL}^{-1} & 0 \\ 0 & 0 \end{array}\right)\).

Full solution from the condensed system

Using these attributes of a, an implementation of the factorization of \(A^{-1}\),

\[\begin{split}\left(\begin{array}{c} u_L \\ u_E \end{array}\right) = \left(\begin{array}{cc} I & -A_{LL}^{-1} A_{LE} \\ 0 & I \end{array}\right) \left(\begin{array}{cc} A_{LL}^{-1} & 0 \\ 0 & S^{-1} \end{array}\right) \left(\begin{array}{cc} I & 0 \\ -A_{EL} A_{LL}^{-1} & I \end{array}\right) \left(\begin{array}{c} f_L \\ f_E \end{array}\right)\end{split}\]

can be given in these few lines:

[4]:
invS = a.mat.Inverse(freedofs=fes.FreeDofs(coupling=True))
ext = IdentityMatrix() + a.harmonic_extension
extT = IdentityMatrix() + a.harmonic_extension_trans
invA =  ext @ invS @ extT + a.inner_solve
[5]:
u.vec.data = invA * f.vec
Draw(u)
[5]:
BaseWebGuiScene

Behind the scenes: CouplingType

How does NGSolve know what is in the index sets \(L\) and \(E\)?

  • Look at fes.CouplingType to see a classification of degrees of freedom (dofs).

[6]:
for i in range(fes.ndof):
    print(fes.CouplingType(i))
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.WIREBASKET_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.INTERFACE_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
COUPLING_TYPE.LOCAL_DOF
[7]:
dof_types = {}
for i in range(fes.ndof):
    ctype = fes.CouplingType(i)
    if ctype in dof_types.keys():
        dof_types[ctype] += 1
    else:
        dof_types[ctype] = 1
dof_types
[7]:
{<COUPLING_TYPE.WIREBASKET_DOF: 8>: 16,
 <COUPLING_TYPE.INTERFACE_DOF: 4>: 99,
 <COUPLING_TYPE.LOCAL_DOF: 2>: 54}

The LOCAL_DOF forms the set \(L\) and the remainder forms the set \(E\). All finite element spaces in NGSolve have such dof classification.

Through this classification a bilinear form is able to automatically compute the Schur complement and the accompanying extension operators. Users need only specify the condense option. (Of course, users should also make sure their method has an invertible \(A_{LL}\)!)

Inhomogeneous Dirichlet boundary conditions

In case of inhomogeneous Dirichlet boundary conditions, we combine the technique of Dirichlet data extension in Unit 1.3 with the above static condensation principle in the following code.

[8]:
U = x*x*(1-y)*(1-y)          # U = manufactured solution
DeltaU = 2*((1-y)*(1-y)+x*x) # Source: DeltaU = ∆U

f = LinearForm(fes)
f += -DeltaU * v * dx
f.Assemble()

u = GridFunction(fes)
u.Set(U, BND)
u.vec.data += a.harmonic_extension * u.vec
u.vec.data += invA * (f.vec - a.mat * u.vec)

sqrt(Integrate((U-u)*(U-u),mesh))  # Compute error
[8]:
7.594386074433171e-16

Automatic utility BVP

The automatic utility solvers.BVP, which we have seen previously, will perform the above steps under the hood, if condense is turned on in its input bilinear form.

[9]:
u = GridFunction(fes)
u.Set(U, BND)

c = Preconditioner(a,"direct")
c.Update()
solvers.BVP(bf=a, lf=f, gf=u, pre=c)
sqrt(Integrate((U-u)*(U-u),mesh))
CG iteration 1, residual = 0.3057548260912156
CG iteration 2, residual = 4.65936468236554e-16
[9]:
7.564237816969331e-16