{
 "cells": [
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "# 1.1 First NGSolve example\n",
    "\n",
    "Let us solve the Poisson problem of finding $u$ satisfying \n",
    "\n",
    "$$\n",
    "\\begin{aligned}\n",
    "-\\Delta u & = f && \\text { in  the unit square},\n",
    "\\\\\n",
    "u & = 0 && \\text{ on the bottom and right parts of the boundary},\n",
    "\\\\\n",
    "\\frac{\\partial u }{\\partial n } & = 0 \n",
    "&& \\text{ on the remaining  boundary parts}.\n",
    "\\end{aligned}\n",
    "$$"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "## Quick steps to solution:"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "#### 1. Import NGSolve and Netgen Python modules:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "from ngsolve import *\n",
    "from ngsolve.webgui import Draw"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "#### 2. Generate an unstructured mesh"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))\n",
    "mesh.nv, mesh.ne   # number of vertices & elements "
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "Here we prescribed a maximal mesh-size of 0.2 using the `maxh` flag. "
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "Draw(mesh);"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "#### 3. Declare a finite element space:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "fes = H1(mesh, order=2, dirichlet=\"bottom|right\")\n",
    "fes.ndof  # number of unknowns in this space"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "Python's help system displays further documentation."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "help(fes)"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "#### 4. Declare test function, trial function, and grid function \n",
    "\n",
    "* Test and trial function are symbolic objects - called `ProxyFunctions` -  that help you construct bilinear forms (and have no space to hold solutions). \n",
    "\n",
    "* `GridFunctions`, on the other hand, represent functions in the finite element space and contains memory to hold coefficient vectors."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "u = fes.TrialFunction()  # symbolic object\n",
    "v = fes.TestFunction()   # symbolic object\n",
    "gfu = GridFunction(fes)  # solution "
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "Alternately, you can get both the trial and test variables at once:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "u, v = fes.TnT()"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "#### 5. Define and assemble linear and bilinear forms:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "a = BilinearForm(fes)\n",
    "a += grad(u)*grad(v)*dx\n",
    "a.Assemble()\n",
    "\n",
    "f = LinearForm(fes)\n",
    "f += x*v*dx\n",
    "f.Assemble();"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "Alternately, we can do one-liners: "
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "a = BilinearForm(grad(u)*grad(v)*dx).Assemble()\n",
    "f = LinearForm(x*v*dx).Assemble()"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "You can examine the linear system in more detail:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "print(f.vec)"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "print(a.mat)"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "#### 6. Solve the system:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "gfu.vec.data = \\\n",
    "    a.mat.Inverse(freedofs=fes.FreeDofs()) * f.vec\n",
    "Draw(gfu);"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "The Dirichlet boundary condition constrains some degrees of freedom. The argument `fes.FreeDofs()` indicates that only the remaining \"free\" degrees of freedom should participate in the linear solve."
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "You can examine the coefficient vector of solution if needed:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "metadata": {},
   "outputs": [],
   "source": [
    "print(gfu.vec)"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "You can see the zeros coming from the zero boundary conditions."
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "## Ways to interact with NGSolve"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "* A jupyter notebook (like this one) gives you one way to interact with NGSolve. When you have a complex sequence of tasks to perform, the notebook may not be adequate.\n",
    "\n",
    "\n",
    "* You can write an entire python module in a text editor and call python on the command line. (A script of the above is provided in `poisson.py`.)\n",
    "    ```\n",
    "    python3 poisson.py\n",
    "    ```\n",
    "  \n",
    "* If you want the Netgen GUI, then use `netgen` on the command line:\n",
    "    ```\n",
    "    netgen poisson.py\n",
    "    ```\n",
    "  You can then ask for a python shell from the GUI's menu options (`Solve -> Python shell`).\n",
    "  "
   ]
  }
 ],
 "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
}
