{
 "cells": [
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "# 1.4 Static Condensation\n",
    "\n",
    "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.\n",
    "\n",
    "In this tutorial, you will learn \n",
    "\n",
    "- how to turn on static condensation (using a keyword argument `condense`),\n",
    "- how to get the full solution from the condensed system using attributes called\n",
    "\n",
    "    - `harmonic_extension_trans`\n",
    "    - `harmonic_extension`\n",
    "    - `inner_solve`,\n",
    "- their relationship to a Schur complement factorization,  \n",
    "- types of degrees of freedom such as `COUPLING_TYPE.LOCAL_DOF`, and \n",
    "- an automatic utility for solving via static condensation."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "from ngsolve import *\n",
    "from ngsolve.webgui import Draw\n",
    "\n",
    "mesh = Mesh(unit_square.GenerateMesh(maxh=0.4))\n",
    "\n",
    "fes = H1(mesh, order=4, dirichlet='bottom|right')\n",
    "u, v = fes.TnT()"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "### Asking `BilinearForm` to `condense`\n",
    "\n",
    "Keyword arguments `condense`, or the deprecated keyword argument `eliminate_internal`, are synonymous. "
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "a = BilinearForm(grad(u) * grad(v) * dx, condense=True).Assemble()\n",
    "f = LinearForm(1 * v * dx).Assemble()\n",
    "u = GridFunction(fes)"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "* The assembled matrix $A=$ `a.mat` can be block partitioned into\n",
    "\n",
    "$$\n",
    "A =\n",
    "\\left(\n",
    "\\begin{array}{cc}\n",
    "A_{LL} & A_{LE}\n",
    "\\\\\n",
    "A_{EL} & A_{EE}\n",
    "\\end{array}\n",
    "\\right)\n",
    "$$\n",
    "where $L$ denotes the set of **local** (or internal)  degrees of freedom and $E$ denotes the set of **interface** degrees of freedom.\n",
    "\n",
    "* 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.)\n",
    "\n",
    "* The condensed system is known as the **Schur complement**:\n",
    "\n",
    "$$\n",
    "S = A_{EE} - A_{EL} A_{LL}^{-1} A_{LE}.\n",
    "$$\n",
    "\n",
    "* When `condense` is set to `True` in `a`, the statement `a.Assemble` actually assembles $S$. Thus \n",
    "    - `a.mat` $= \\left(\\begin{array}{cc}0 & 0 \\\\0 & S\\end{array}\\right),$ and\n",
    "\n",
    "    - `a.mat.Inverse(freedofs=fes.FreeDofs(coupling=True))` $= \\left(\\begin{array}{cc}\n",
    "0 & 0 \\\\\n",
    "0 & S^{-1}\n",
    "\\end{array}\\right).$\n",
    "Let us make this matrix:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "invS = a.mat.Inverse(freedofs=fes.FreeDofs(coupling=True))"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "Note that instead of the usual `fes.FreeDofs()`, we have used\n",
    "        `fes.FreeDofs(coupling=True)`\n",
    "        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.)"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "### A factorization of the inverse"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "* To solve\n",
    "    $$ \\left(\\begin{array}{cc}\n",
    "A_{LL} & A_{LE}\\\\\n",
    "A_{EL} & A_{EE}\n",
    "\\end{array}\\right)\n",
    "\\left(\\begin{array}{c}\n",
    "u_L \\\\ u_E\n",
    "\\end{array}\\right)= \n",
    "\\left(\\begin{array}{c}\n",
    "f_L \\\\ f_E\n",
    "\\end{array}\\right)\n",
    "    $$\n",
    "we use a factorization of $A^{-1}$ that uses $S^{-1}$. Namely, we use the following identity:\n",
    "\n",
    "$$\n",
    "\\left(\\begin{array}{c}\n",
    "u_L \\\\ u_E\n",
    "\\end{array}\\right)\n",
    "=\n",
    "\\left(\\begin{array}{cc}\n",
    "I & -A_{LL}^{-1} A_{LE} \\\\\n",
    "0 & I \n",
    "\\end{array}\\right)\n",
    "\\left(\\begin{array}{cc}\n",
    "A_{LL}^{-1} & 0 \\\\\n",
    "0 & S^{-1}\n",
    "\\end{array}\\right)\n",
    "\\left(\\begin{array}{cc}\n",
    "I & 0 \\\\\n",
    "-A_{EL} A_{LL}^{-1} & I \n",
    "\\end{array}\\right)\n",
    "\\left(\\begin{array}{c}\n",
    "f_L \\\\ f_E\n",
    "\\end{array}\\right)\n",
    "$$"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "* NGSolve provides \n",
    "\n",
    "    - `a.harmonic_extension_trans` $=\n",
    "\\left(\\begin{array}{cc}\n",
    "0 & 0 \\\\\n",
    "-A_{EL} A_{LL}^{-1} & 0 \n",
    "\\end{array}\\right)$\n",
    "   - `a.harmonic_extension` $=\n",
    "\\left(\\begin{array}{cc}\n",
    "0 & -A_{LL}^{-1} A_{LE}\n",
    "\\\\\n",
    "0 & 0\n",
    "\\end{array}\\right)$\n",
    "\n",
    "    - `a.inner_solve`  $=\\left(\\begin{array}{cc}\n",
    "A_{LL}^{-1} & 0 \n",
    "\\\\\n",
    "0 & 0\n",
    "\\end{array}\\right)$."
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "### Full solution from the condensed system "
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "Using these attributes of `a`, an implementation of the factorization of $A^{-1}$,\n",
    "$$\n",
    "\\left(\\begin{array}{c}\n",
    "u_L \\\\ u_E\n",
    "\\end{array}\\right)\n",
    "=\n",
    "\\left(\\begin{array}{cc}\n",
    "I & -A_{LL}^{-1} A_{LE} \\\\\n",
    "0 & I \n",
    "\\end{array}\\right)\n",
    "\\left(\\begin{array}{cc}\n",
    "A_{LL}^{-1} & 0 \\\\\n",
    "0 & S^{-1}\n",
    "\\end{array}\\right)\n",
    "\\left(\\begin{array}{cc}\n",
    "I & 0 \\\\\n",
    "-A_{EL} A_{LL}^{-1} & I \n",
    "\\end{array}\\right)\n",
    "\\left(\\begin{array}{c}\n",
    "f_L \\\\ f_E\n",
    "\\end{array}\\right)\n",
    "$$\n",
    "can be given in these few lines:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "invS = a.mat.Inverse(freedofs=fes.FreeDofs(coupling=True))\n",
    "ext = IdentityMatrix() + a.harmonic_extension\n",
    "extT = IdentityMatrix() + a.harmonic_extension_trans\n",
    "invA =  ext @ invS @ extT + a.inner_solve"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "u.vec.data = invA * f.vec\n",
    "Draw(u)"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "### Behind the scenes: `CouplingType`"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "How does NGSolve know what is in the index sets $L$ and $E$? \n",
    "\n",
    "* Look at `fes.CouplingType` to see a classification of degrees of freedom (dofs).\n"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "for i in range(fes.ndof):\n",
    "    print(fes.CouplingType(i))"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "dof_types = {}\n",
    "for i in range(fes.ndof):\n",
    "    ctype = fes.CouplingType(i)\n",
    "    if ctype in dof_types.keys():\n",
    "        dof_types[ctype] += 1\n",
    "    else:\n",
    "        dof_types[ctype] = 1\n",
    "dof_types"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "The `LOCAL_DOF` forms the set $L$ and the remainder forms the set $E$. All finite element spaces in NGSolve have such dof classification.\n",
    "\n",
    "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}$!)\n"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "### Inhomogeneous Dirichlet boundary conditions\n",
    "\n",
    "In case of inhomogeneous Dirichlet boundary conditions, we combine the technique of Dirichlet data extension in \n",
    "[Unit 1.3](../unit-1.3-dirichlet/dirichlet.ipynb) \n",
    "with the above static condensation principle in the following code."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "U = x*x*(1-y)*(1-y)          # U = manufactured solution\n",
    "DeltaU = 2*((1-y)*(1-y)+x*x) # Source: DeltaU = ∆U\n",
    "\n",
    "f = LinearForm(fes)\n",
    "f += -DeltaU * v * dx\n",
    "f.Assemble()\n",
    "\n",
    "u = GridFunction(fes)\n",
    "u.Set(U, BND)               \n",
    "u.vec.data += a.harmonic_extension * u.vec \n",
    "u.vec.data += invA * (f.vec - a.mat * u.vec)\n",
    "\n",
    "sqrt(Integrate((U-u)*(U-u),mesh))  # Compute error"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "### Automatic utility `BVP`\n",
    "\n",
    "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. "
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "u = GridFunction(fes)\n",
    "u.Set(U, BND) \n",
    "\n",
    "c = Preconditioner(a,\"direct\")\n",
    "c.Update()\n",
    "solvers.BVP(bf=a, lf=f, gf=u, pre=c)\n",
    "sqrt(Integrate((U-u)*(U-u),mesh))"
   ]
  }
 ],
 "metadata": {
  "kernelspec": {
   "display_name": "Python 3 (ipykernel)",
   "language": "python",
   "name": "python3"
  },
  "language_info": {
   "codemirror_mode": {
    "name": "ipython",
    "version": 3
   },
   "file_extension": ".py",
   "mimetype": "text/x-python",
   "name": "python",
   "nbconvert_exporter": "python",
   "pygments_lexer": "ipython3",
   "version": "3.11.3"
  }
 },
 "nbformat": 4,
 "nbformat_minor": 4
}
