{
 "cells": [
  {
   "cell_type": "markdown",
   "id": "0",
   "metadata": {},
   "source": [
    "# 9.1 Implementation of Finite Elements\n",
    "\n",
    "This tutorial shows how to implement our own finite elements in C++,\n",
    "and how to use them within the NGSolve language. We implement first order and second order triangular finite elements.\n",
    "\n",
    "* Finite elements implement the basis functions:\n",
    "[myElement.hpp](myElement.hpp) \n",
    "[myElement.cpp](myElement.cpp)\n",
    "\n",
    "* Differential operators define the mapping from coefficients to function values:\n",
    "[myDiffOp.hpp](myDiffOp.hpp) \n",
    "\n",
    "* Finite element spaces implement the enumeration of degrees of freedom, and creation of elements:\n",
    "[myFESpace.hpp](myFESpace.hpp) \n",
    "[myFESpace.cpp](myFESpace.cpp)\n",
    "\n",
    "* Everything is included from the main file [mymodule.cpp](mymodule.cpp)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "1",
   "metadata": {},
   "source": [
    "We read the main cpp-File to a string, and let NGSolve call the compiler, and load the new library as a Python module:"
   ]
  },
  {
   "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.MyFESpace(mesh, secondorder=True, 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*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(1*v*dx).Assemble()\n",
    "gfu.vec.data = a.mat.Inverse(fes.FreeDofs())*f.vec\n",
    "Draw (gfu);"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "11",
   "metadata": {},
   "source": [
    "Draw basis functions:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "12",
   "metadata": {},
   "outputs": [],
   "source": [
    "gfu.vec[:] = 0\n",
    "gfu.vec[mesh.nv-3] = 1\n",
    "gfu.vec[fes.ndof-1] = 1\n",
    "Draw (gfu, order=2);"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "13",
   "metadata": {},
   "source": [
    "Documentation provided in the DocInfo structure is available in Python - help. Look for 'secondorder'."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "14",
   "metadata": {},
   "outputs": [],
   "source": [
    "help (m.MyFESpace)"
   ]
  },
  {
   "cell_type": "markdown",
   "id": "15",
   "metadata": {},
   "source": [
    "**Exercises:**\n",
    "\n",
    "Extend MyFESpace by the following 1st and 2nd order elements:\n",
    "\n",
    "- 1D finite elements (ET_SEGM), as needed for boundary conditions\n",
    "- quadrilateral elements (ET_QUAD), use geom.GenerateMesh(quad_dominated=True)\n",
    "- tetrahedral elements (ET_TET), test it for 3D domains\n",
    "\n",
    "Next, implement 3rd order elements"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "16",
   "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
}
