{
 "cells": [
  {
   "attachments": {},
   "cell_type": "markdown",
   "id": "0",
   "metadata": {},
   "source": [
    "# VTK - Output"
   ]
  },
  {
   "attachments": {},
   "cell_type": "markdown",
   "id": "1",
   "metadata": {},
   "source": [
    "Netgen/NGSolve comes with several convenient ways to visualize your PDE solutions: the webgui, the legacy netgen GUI, ... For some use cases visualization toolboxes like [ParaView](https://www.paraview.org/) may still be able to obtain more insight or more complex visualizations from you computational results.\n",
    "\n",
    "To put your results in a corresponding format you can use the `VTKOutput` of NGSolve. \n",
    "\n",
    "A `VTKOutput` object obtains a list of coefficient functions (scalar or vector) and a list of labels. Together with the mesh information these information will be put into a vtk-formatted file on every `Do`-call of the object. As VTK essentially only draws piecewise linears you may want to describe a `subdivision` argument larger `0` to obtain higher resolution of a higher order FE function.\n",
    "\n",
    "Let's explore a few use cases by examples:"
   ]
  },
  {
   "attachments": {},
   "cell_type": "markdown",
   "id": "2",
   "metadata": {},
   "source": [
    "## Example 1: single solution output\n",
    "Let's recap the Poisson solution from [unit 1.1](../unit-1.1-poisson/poisson.ipynb):"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "3",
   "metadata": {},
   "outputs": [],
   "source": [
    "from ngsolve import *\n",
    "from netgen.geom2d import unit_square\n",
    "mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))\n",
    "fes = H1(mesh, order=2, dirichlet=\"bottom|right\")\n",
    "gfu = GridFunction(fes) \n",
    "a = BilinearForm(fes, symmetric=True)\n",
    "a += grad(fes.TrialFunction())*grad(fes.TestFunction())*dx\n",
    "a.Assemble()\n",
    "f = LinearForm(fes)\n",
    "f += x*fes.TestFunction()*dx\n",
    "f.Assemble()\n",
    "gfu.vec.data = a.mat.Inverse(freedofs=fes.FreeDofs()) * f.vec"
   ]
  },
  {
   "attachments": {},
   "cell_type": "markdown",
   "id": "4",
   "metadata": {},
   "source": [
    "Now, instead of using the webgui, let us export the result data to a vtk file:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "5",
   "metadata": {},
   "outputs": [],
   "source": [
    "vtk = VTKOutput(mesh,coefs=[gfu],names=[\"sol\"],filename=\"vtk_example1\",subdivision=2)\n",
    "vtk.Do()"
   ]
  },
  {
   "attachments": {},
   "cell_type": "markdown",
   "id": "6",
   "metadata": {},
   "source": [
    "Note, only the `Do()`-call writes the output. Let's see what we got:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "7",
   "metadata": {},
   "outputs": [],
   "source": [
    "!ls -lh vtk_example1*"
   ]
  },
  {
   "attachments": {},
   "cell_type": "markdown",
   "id": "8",
   "metadata": {},
   "source": [
    "We obtained the file `vtk_example1.vtu`. This file can be opened with ParaView, e.g. with `paraview vtk_example1.vtu` where you can select your solution field and use several \"filters\" to process the data to get nice visualizations as the following:\n",
    "<div>\n",
    "<img src=\"vtk_example1.png\" width=\"600\"/>\n",
    "</div>\n",
    "You can now play around with tools like paraview to do many different complex operations."
   ]
  },
  {
   "attachments": {},
   "cell_type": "markdown",
   "id": "9",
   "metadata": {},
   "source": [
    "## Example 2: Series of data\n",
    "\n",
    "In some cases you may want to analyse a sequence of simulation, e.g. for generating a video. In this case, you simply call `Do()` several times. You can also provide an additional time stamp with e.g. `Do(time=0.3)`.\n",
    "\n",
    "Let's take a time-dependent coefficient function and put this to a sequence of VTK outputs:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "10",
   "metadata": {},
   "outputs": [],
   "source": [
    "t = Parameter(0)\n",
    "coef = sin(4*(x+y+t))\n",
    "vtkout = VTKOutput(mesh,coefs=[coef],names=[\"sol\"],filename=\"vtk_example2\",subdivision=2)\n",
    "vtkout.Do(time = 0.0)\n",
    "for i in range(20):\n",
    "    t.Set((i+1)/20)\n",
    "    vtkout.Do(time = (i+1)/20.0)"
   ]
  },
  {
   "attachments": {},
   "cell_type": "markdown",
   "id": "11",
   "metadata": {},
   "source": [
    "Let's see what we got this time:"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "12",
   "metadata": {},
   "outputs": [],
   "source": [
    "ls -lh vtk_example2*"
   ]
  },
  {
   "attachments": {},
   "cell_type": "markdown",
   "id": "13",
   "metadata": {},
   "source": [
    "\n",
    "* We obtained 21 `.vtu`-files corresponding to the 21 `Do`-calls.\n",
    "* Only the first `Do`-call does not have an additional suffix, while later calls append the `_step...`-suffix. \n",
    "* Furthermore there is a `.pvd` file which bundles the individual files and associates them to the provided time stamps (if provided):"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "14",
   "metadata": {},
   "outputs": [],
   "source": [
    "!cat vtk_example2.pvd"
   ]
  },
  {
   "attachments": {},
   "cell_type": "markdown",
   "id": "15",
   "metadata": {},
   "source": [
    "Now you can open this \"meta file\" with paraview and visualize the time series with many different filter. Finally, you can export stills or videos from that. It will look like this:\n",
    "<div>\n",
    "<img src=\"vtk_example2.png\" width=\"600\"/>\n",
    "</div>"
   ]
  },
  {
   "attachments": {},
   "cell_type": "markdown",
   "id": "16",
   "metadata": {},
   "source": [
    "## Further options\n",
    "\n",
    "To learn more about possible options you show call `help(VTKOutput)`. Nevertheless, let us comment on a few of theme:\n",
    "\n",
    "* `VTKOutput(..., subdivision=i)` refines the mesh (locally and only temporarily) as often as `i`. Note that when opened with ParaView the visualized mesh no longer coincides with your computational mesh.\n",
    "* `VTKOutput(..., only_element=i)` puts only the data of one selected element out. This can be interesting for basis function visualization or alike.\n",
    "* `VTKOutput(..., floatsize=fs)` where `fs` is either `\"single\"` or `\"double\"` decides on the precision of the output. Default is `\"double\"`. Using `\"single\"` can significantly reduce the used disk space.\n",
    "\n",
    "there are further options for the `Do(...)`-call that we briefly discuss:\n",
    "\n",
    "* `....Do(time=t)` adds the time stamp `t` to the meta-file if a series of outputs is used.\n",
    "* `....Do(vb=vb)` allows you to switch from volume mesh output (`vb=VOL`, default) to surface mesh output (`vb=BND`)\n",
    "* `....Do(drawelems=els)` allows you to only write a submesh of the underlying mesh where the corresponding set of elements is prescribed by the `BitArray` `els`.\n"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": null,
   "id": "17",
   "metadata": {},
   "outputs": [],
   "source": [
    "help(VTKOutput)"
   ]
  }
 ],
 "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.10.10"
  }
 },
 "nbformat": 4,
 "nbformat_minor": 5
}
