{
 "cells": [
  {
   "cell_type": "markdown",
   "id": "0",
   "metadata": {},
   "source": [
    "# 9.3 High Order Finite Elements\n",
    "\n"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "1",
   "metadata": {},
   "source": [
    "* Finite elements implement the basis functions:\n",
    "[myHOElement.hpp](myHOElement.hpp) \n",
    "[myHOElement.cpp](myHOElement.cpp)\n",
    "\n",
    "* Finite element spaces implement the enumeration of degrees of freedom, and creation of elements:\n",
    "[myHOFESpace.hpp](myHOFESpace.hpp) \n",
    "[myHOFESpace.cpp](myHOFESpace.cpp)\n",
    "\n"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "2",
   "metadata": {},
   "outputs": [],
   "source": [
    "from ngsolve.fem import CompilePythonModule\n",
    "from pathlib import Path\n",
    "\n",
    "m = CompilePythonModule(Path('mymodule.cpp'), init_function_name='mymodule')"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "3",
   "metadata": {},
   "outputs": [],
   "source": [
    "from netgen.occ import unit_square\n",
    "from ngsolve import *\n",
    "from ngsolve.webgui import Draw\n",
    "\n",
    "mesh = Mesh(unit_square.GenerateMesh(maxh=0.2, quad_dominated=False))"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "4",
   "metadata": {},
   "source": [
    "We can now create an instance of our own finite element space:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "5",
   "metadata": {},
   "outputs": [],
   "source": [
    "fes = m.MyHighOrderFESpace(mesh, order=4, dirichlet=\"left|bottom|top\")"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "6",
   "metadata": {},
   "source": [
    "and use it within NGSolve such as the builtin finite element spaces:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "7",
   "metadata": {},
   "outputs": [],
   "source": [
    "print (\"ndof = \", fes.ndof)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "8",
   "metadata": {},
   "outputs": [],
   "source": [
    "gfu = GridFunction(fes)\n",
    "gfu.Set(x*x*y*y)\n",
    "\n",
    "Draw (gfu)\n",
    "Draw (grad(gfu)[0], mesh);"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "9",
   "metadata": {},
   "source": [
    "and solve the standard problem:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "10",
   "metadata": {},
   "outputs": [],
   "source": [
    "u,v = fes.TnT()\n",
    "a = BilinearForm(grad(u)*grad(v)*dx).Assemble()\n",
    "f = LinearForm(10*v*dx).Assemble()\n",
    "gfu.vec.data = a.mat.Inverse(fes.FreeDofs())*f.vec\n",
    "Draw (gfu, order=3);"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "11",
   "metadata": {},
   "outputs": [],
   "source": [
    "errlist = []\n",
    "for p in range(1,13):\n",
    "    fes = m.MyHighOrderFESpace(mesh, order=p)\n",
    "    func = sin(pi*x)*sin(pi*y)\n",
    "    gfu = GridFunction(fes)\n",
    "    gfu.Set(func)\n",
    "    err = sqrt(Integrate( (func-gfu)**2, mesh, order=5+2*p))\n",
    "    errlist.append((p,err))\n",
    "print (errlist)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "12",
   "metadata": {},
   "outputs": [],
   "source": [
    "import matplotlib.pyplot as plt\n",
    "n,err = zip(*errlist)\n",
    "plt.yscale('log')\n",
    "plt.plot(n,err);"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "13",
   "metadata": {},
   "source": [
    "**Exercises:**\n",
    "\n",
    "Extend MyHighOrderFESpace by high order quadrilateral elements.\n",
    "\n",
    "http://www.numa.uni-linz.ac.at/Teaching/PhD/Finished/zaglmayr-diss.pdf, \n",
    "page 68 ff"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "14",
   "metadata": {},
   "outputs": [],
   "source": []
  }
 ],
 "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.12.1"
  }
 },
 "nbformat": 4,
 "nbformat_minor": 5
}
