{
 "cells": [
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "<script async src=\"https://www.googletagmanager.com/gtag/js?id=UA-59152712-8\"></script>\n",
    "<script>\n",
    "  window.dataLayer = window.dataLayer || [];\n",
    "  function gtag(){dataLayer.push(arguments);}\n",
    "  gtag('js', new Date());\n",
    "\n",
    "  gtag('config', 'UA-59152712-8');\n",
    "</script>\n",
    "\n",
    "# Tutorial-IllinoisGRMHD: IllinoisGRMHD_headers.h\n",
    "\n",
    "## Authors: Leo Werneck & Zach Etienne\n",
    "\n",
    "<font color='red'>**This module is currently under development**</font>\n",
    "\n",
    "## In this tutorial module we explain the main `IllinoisGRMHD` header file. This module will likely be absorbed by another one by the time we finish documenting the code\n",
    "\n",
    "### Required and recommended citations:\n",
    "\n",
    "* **(Required)** Etienne, Z. B., Paschalidis, V., Haas R., Mösta P., and Shapiro, S. L. IllinoisGRMHD: an open-source, user-friendly GRMHD code for dynamical spacetimes. Class. Quantum Grav. 32 (2015) 175009. ([arxiv:1501.07276](http://arxiv.org/abs/1501.07276)).\n",
    "* **(Required)** Noble, S. C., Gammie, C. F., McKinney, J. C., Del Zanna, L. Primitive Variable Solvers for Conservative General Relativistic Magnetohydrodynamics. Astrophysical Journal, 641, 626 (2006) ([astro-ph/0512420](https://arxiv.org/abs/astro-ph/0512420)).\n",
    "* **(Recommended)** Del Zanna, L., Bucciantini N., Londrillo, P. An efficient shock-capturing central-type scheme for multidimensional relativistic flows - II. Magnetohydrodynamics. A&A 400 (2) 397-413 (2003). DOI: 10.1051/0004-6361:20021641 ([astro-ph/0210618](https://arxiv.org/abs/astro-ph/0210618))."
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "<a id='toc'></a>\n",
    "\n",
    "# Table of Contents\n",
    "$$\\label{toc}$$\n",
    "\n",
    "This module is organized as follows\n",
    "\n",
    "0. [Step 0](#src_dir): **Source directory creation**\n",
    "1. [Step 1](#introduction): **Introduction**\n",
    "1. [Step 2](#igm_headers__h): **`IllinoisGRMHD_headers.h`**\n",
    "1. [Step n-1](#code_validation): **Code validation**\n",
    "1. [Step n](#latex_pdf_output): **Output this notebook to $\\LaTeX$-formatted PDF file**"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "<a id='src_dir'></a>\n",
    "\n",
    "# Step 0: Source directory creation \\[Back to [top](#toc)\\]\n",
    "$$\\label{src_dir}$$\n",
    "\n",
    "We will now use the [cmdline_helper.py NRPy+ module](Tutorial-Tutorial-cmdline_helper.ipynb) to create the source directory within the `IllinoisGRMHD` NRPy+ directory, if it does not exist yet."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 1,
   "metadata": {},
   "outputs": [],
   "source": [
    "# Step 0: Creation of the IllinoisGRMHD source directory\n",
    "# Step 0a: Add NRPy's directory to the path\n",
    "# https://stackoverflow.com/questions/16780014/import-file-from-parent-directory\n",
    "import os,sys\n",
    "nrpy_dir_path = os.path.join(\"..\",\"..\")\n",
    "if nrpy_dir_path not in sys.path:\n",
    "    sys.path.append(nrpy_dir_path)\n",
    "\n",
    "# Step 0b: Load up cmdline_helper and create the directory\n",
    "import cmdline_helper as cmd\n",
    "IGM_src_dir_path = os.path.join(\"..\",\"src\")\n",
    "cmd.mkdir(IGM_src_dir_path)\n",
    "\n",
    "# Step 0c: Create the output file path\n",
    "outfile_path__IllinoisGRMHD_headers__h = os.path.join(IGM_src_dir_path,\"IllinoisGRMHD_headers.h\")"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "<a id='introduction'></a>\n",
    "\n",
    "# Step 1: Introduction \\[Back to [top](#toc)\\]\n",
    "$$\\label{introduction}$$"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "<a id='igm_headers__h'></a>\n",
    "\n",
    "# Step 2: `IllinoisGRMHD_headers.h` \\[Back to [top](#toc)\\]\n",
    "$$\\label{igm_headers__h}$$"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 2,
   "metadata": {},
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "Writing ../src/IllinoisGRMHD_headers.h\n"
     ]
    }
   ],
   "source": [
    "%%writefile $outfile_path__IllinoisGRMHD_headers__h\n",
    "// To safeguard against double-including this header file:\n",
    "#ifndef ILLINOISGRMHD_HEADERS_H_\n",
    "#define ILLINOISGRMHD_HEADERS_H_\n",
    "\n",
    "#define MIN(a,b) ( ((a) < (b)) ? (a) : (b) )\n",
    "#define MAX(a,b) ( ((a) > (b)) ? (a) : (b) )\n",
    "#define SQR(x) ((x) * (x))\n",
    "#define ONE_OVER_SQRT_4PI 0.282094791773878143474039725780\n",
    "\n",
    "#define VERR_DEF_PARAMS __LINE__, __FILE__, CCTK_THORNSTRING\n",
    "\n",
    "// The order here MATTERS, as we assume that GUPXX+1=GUPYY, etc.\n",
    "static const int PHI=0,PSI=1,GXX=2,GXY=3,GXZ=4,GYY=5,GYZ=6,GZZ=7,\n",
    "  LAPM1=8,SHIFTX=9,SHIFTY=10,SHIFTZ=11,GUPXX=12,GUPYY=13,GUPZZ=14,\n",
    "  NUMVARS_FOR_METRIC_FACEVALS=15; //<-- Be _sure_ to set this correctly, or you'll have memory access bugs!\n",
    "\n",
    "// These are not used for facevals in the reconstruction step, but boy are they useful anyway.\n",
    "static const int GUPXY=15,GUPXZ=16,GUPYZ=17,\n",
    "  NUMVARS_FOR_METRIC=18; //<-- Be _sure_ to set this correctly, or you'll have memory access bugs!\n",
    "\n",
    "// The order here MATTERS, and must be consistent with the order in the in_prims[] array in driver_evaluate_MHD_rhs.C.\n",
    "static const int RHOB=0,PRESSURE=1,VX=2,VY=3,VZ=4,\n",
    "  BX_CENTER=5,BY_CENTER=6,BZ_CENTER=7,BX_STAGGER=8,BY_STAGGER=9,BZ_STAGGER=10,\n",
    "  VXR=11,VYR=12,VZR=13,VXL=14,VYL=15,VZL=16,MAXNUMVARS=17;  //<-- Be _sure_ to define MAXNUMVARS appropriately!\n",
    "\n",
    "static const int UT=0,UX=1,UY=2,UZ=3;\n",
    "\n",
    "// The \"I\" suffix denotes interpolation. In other words, these\n",
    "//    definitions are used for interpolation ONLY. The order here\n",
    "//    matters as well!\n",
    "static const int SHIFTXI=0,SHIFTYI=1,SHIFTZI=2,GUPXXI=3,GUPXYI=4,GUPXZI=5,GUPYYI=6,GUPYZI=7,GUPZZI=8,\n",
    "  PSII=9,LAPM1I=10,A_XI=11,A_YI=12,A_ZI=13,LAPSE_PSI2I=14,LAPSE_OVER_PSI6I=15,MAXNUMINTERP=16;\n",
    "\n",
    "// Again, the order here MATTERS, since we assume in the code that, e.g., smallb[0]=b^t, smallb[3]=b^z, etc.\n",
    "static const int SMALLBT=0,SMALLBX=1,SMALLBY=2,SMALLBZ=3,SMALLB2=4,NUMVARS_SMALLB=5;\n",
    "\n",
    "// Again, the order here MATTERS, since we assume in the code that, CONSERV[STILDEX+1] = \\tilde{S}_y\n",
    "static const int RHOSTAR=0,STILDEX=1,STILDEY=2,STILDEZ=3,TAUENERGY=4,NUM_CONSERVS=5;\n",
    "\n",
    "static const int LAPSE=0,PSI2=1,PSI4=2,PSI6=3,PSIM4=4,LAPSEINV=5,NUMVARS_METRIC_AUX=6;\n",
    "#define SET_LAPSE_PSI4(array_name,METRIC)   {                   \\\n",
    "      array_name[LAPSE] = METRIC[LAPM1]+1.0;                    \\\n",
    "      array_name[PSI2]  = exp(2.0*METRIC[PHI]);                 \\\n",
    "      array_name[PSI4]  = SQR(array_name[PSI2]);                \\\n",
    "      array_name[PSI6]  = array_name[PSI4]*array_name[PSI2];    \\\n",
    "      array_name[PSIM4]  = 1.0/array_name[PSI4];                \\\n",
    "      array_name[LAPSEINV]  = 1.0/array_name[LAPSE];            \\\n",
    "  }\n",
    "\n",
    "// Keeping track of ghostzones between routines is a nightmare, so\n",
    "//   we instead attach ghostzone info to each gridfunction and set\n",
    "//   the ghostzone information correctly within each routine.\n",
    "struct gf_and_gz_struct {\n",
    "  CCTK_REAL *gf;\n",
    "  int gz_lo[4],gz_hi[4];\n",
    "};\n",
    "\n",
    "#define MAX_EOS_PARAMS 10\n",
    "struct eos_struct {\n",
    "  int neos;\n",
    "  CCTK_REAL rho_ppoly_tab[MAX_EOS_PARAMS-1];\n",
    "  CCTK_REAL eps_integ_const[MAX_EOS_PARAMS],K_ppoly_tab[MAX_EOS_PARAMS],Gamma_ppoly_tab[MAX_EOS_PARAMS];\n",
    "};\n",
    "\n",
    "struct output_stats {\n",
    "  int font_fixed,vel_limited,failure_checker;\n",
    "  long n_iter;\n",
    "};\n",
    "\n",
    "\n",
    "// FIXME: For cosmetic purposes, we might want to make everything either zero-offset or one-offset, instead of a mixture.\n",
    "const int kronecker_delta[4][3] = { { 0,0,0 },\n",
    "                                    { 1,0,0 },\n",
    "                                    { 0,1,0 },\n",
    "                                    { 0,0,1 } };\n",
    "\n",
    "/* PUBLIC FUNCTIONS, USED OUTSIDE IllinoisGRMHD AS WELL */\n",
    "void IllinoisGRMHD_enforce_limits_on_primitives_and_recompute_conservs(const int already_computed_physical_metric_and_inverse,CCTK_REAL *U,struct output_stats &stats,eos_struct &eos,\n",
    "                                                                       CCTK_REAL *METRIC,CCTK_REAL g4dn[4][4],CCTK_REAL g4up[4][4], CCTK_REAL *TUPMUNU,CCTK_REAL *TDNMUNU,CCTK_REAL *CONSERVS);\n",
    "\n",
    "void IllinoisGRMHD_convert_ADM_to_BSSN__enforce_detgtij_eq_1__and_compute_gtupij\n",
    "(const cGH *cctkGH,const int *cctk_lsh,\n",
    " CCTK_REAL *gxx,CCTK_REAL *gxy,CCTK_REAL *gxz,CCTK_REAL *gyy,CCTK_REAL *gyz,CCTK_REAL *gzz,CCTK_REAL *alp,\n",
    " CCTK_REAL *gtxx,CCTK_REAL *gtxy,CCTK_REAL *gtxz,CCTK_REAL *gtyy,CCTK_REAL *gtyz,CCTK_REAL *gtzz,\n",
    " CCTK_REAL *gtupxx,CCTK_REAL *gtupxy,CCTK_REAL *gtupxz,CCTK_REAL *gtupyy,CCTK_REAL *gtupyz,CCTK_REAL *gtupzz,\n",
    " CCTK_REAL *phi,CCTK_REAL *psi,CCTK_REAL *lapm1);\n",
    "\n",
    "void IllinoisGRMHD_set_symmetry_gzs_staggered(const cGH *cctkGH, const int *cctk_lsh,CCTK_REAL *X,CCTK_REAL *Y,CCTK_REAL *Z,  CCTK_REAL *gridfunc,\n",
    "                                              CCTK_REAL *gridfunc_syms,int stagger_x,int stagger_y,int stagger_z);\n",
    "\n",
    "#include \"IllinoisGRMHD_EoS_lowlevel_functs.C\"\n",
    "#endif // ILLINOISGRMHD_HEADERS_H\n",
    "\n"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "<a id='code_validation'></a>\n",
    "\n",
    "# Step n-1: Code validation \\[Back to [top](#toc)\\]\n",
    "$$\\label{code_validation}$$\n",
    "\n",
    "First we download the original `IllinoisGRMHD` source code and then compare it to the source code generated by this tutorial notebook."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 3,
   "metadata": {},
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "Validation test for IllinoisGRMHD_headers.h: FAILED!\n",
      "Diff:\n",
      "17c17\n",
      "< // These are not used for facevals in the reconstruction step, but boy are they useful anyway. \n",
      "---\n",
      "> // These are not used for facevals in the reconstruction step, but boy are they useful anyway.\n",
      "61,64c61,63\n",
      "<   CCTK_REAL rho_tab[MAX_EOS_PARAMS],P_tab[MAX_EOS_PARAMS],eps_tab[MAX_EOS_PARAMS],k_tab[MAX_EOS_PARAMS],gamma_tab[MAX_EOS_PARAMS];\n",
      "<   CCTK_REAL gamma_th;\n",
      "<   CCTK_REAL K_poly;\n",
      "< };  \n",
      "---\n",
      ">   CCTK_REAL rho_ppoly_tab[MAX_EOS_PARAMS-1];\n",
      ">   CCTK_REAL eps_integ_const[MAX_EOS_PARAMS],K_ppoly_tab[MAX_EOS_PARAMS],Gamma_ppoly_tab[MAX_EOS_PARAMS];\n",
      "> };\n",
      "91a91\n",
      "> #include \"IllinoisGRMHD_EoS_lowlevel_functs.C\"\n",
      "92a93\n",
      "> \n"
     ]
    }
   ],
   "source": [
    "# Verify if the code generated by this tutorial module\n",
    "# matches the original IllinoisGRMHD source code\n",
    "\n",
    "# First download the original IllinoisGRMHD source code\n",
    "import urllib\n",
    "from os import path\n",
    "\n",
    "original_IGM_file_url  = \"https://bitbucket.org/zach_etienne/wvuthorns/raw/5611b2f0b17135538c9d9d17c7da062abe0401b6/IllinoisGRMHD/src/IllinoisGRMHD_headers.h\"\n",
    "original_IGM_file_name = \"IllinoisGRMHD_headers-original.h\"\n",
    "original_IGM_file_path = os.path.join(IGM_src_dir_path,original_IGM_file_name)\n",
    "\n",
    "# Then download the original IllinoisGRMHD source code\n",
    "# We try it here in a couple of ways in an attempt to keep\n",
    "# the code more portable\n",
    "try:\n",
    "    original_IGM_file_code = urllib.request.urlopen(original_IGM_file_url).read().decode(\"utf-8\")\n",
    "    # Write down the file the original IllinoisGRMHD source code\n",
    "    with open(original_IGM_file_path,\"w\") as file:\n",
    "        file.write(original_IGM_file_code)\n",
    "except:\n",
    "    try:\n",
    "        original_IGM_file_code = urllib.urlopen(original_IGM_file_url).read().decode(\"utf-8\")\n",
    "        # Write down the file the original IllinoisGRMHD source code\n",
    "        with open(original_IGM_file_path,\"w\") as file:\n",
    "            file.write(original_IGM_file_code)\n",
    "    except:\n",
    "        # If all else fails, hope wget does the job\n",
    "        !wget -O $original_IGM_file_path $original_IGM_file_url\n",
    "\n",
    "# Perform validation\n",
    "Validation__IllinoisGRMHD_headers__h  = !diff $original_IGM_file_path $outfile_path__IllinoisGRMHD_headers__h\n",
    "\n",
    "if Validation__IllinoisGRMHD_headers__h == []:\n",
    "    # If the validation passes, we do not need to store the original IGM source code file\n",
    "    !rm $original_IGM_file_path\n",
    "    print(\"Validation test for IllinoisGRMHD_headers.h: PASSED!\")\n",
    "else:\n",
    "    # If the validation fails, we keep the original IGM source code file\n",
    "    print(\"Validation test for IllinoisGRMHD_headers.h: FAILED!\")\n",
    "    # We also print out the difference between the code generated\n",
    "    # in this tutorial module and the original IGM source code\n",
    "    print(\"Diff:\")\n",
    "    for diff_line in Validation__IllinoisGRMHD_headers__h:\n",
    "        print(diff_line)"
   ]
  },
  {
   "cell_type": "markdown",
   "metadata": {},
   "source": [
    "<a id='latex_pdf_output'></a>\n",
    "\n",
    "# Step n: Output this notebook to $\\LaTeX$-formatted PDF file \\[Back to [top](#toc)\\]\n",
    "$$\\label{latex_pdf_output}$$\n",
    "\n",
    "The following code cell converts this Jupyter notebook into a proper, clickable $\\LaTeX$-formatted PDF file. After the cell is successfully run, the generated PDF may be found in the root NRPy+ tutorial directory, with filename\n",
    "[Tutorial-IllinoisGRMHD__IllinoisGRMHD_headers.pdf](Tutorial-IllinoisGRMHD__IllinoisGRMHD_headers.pdf) (Note that clicking on this link may not work; you may need to open the PDF file through another means)."
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 4,
   "metadata": {},
   "outputs": [],
   "source": [
    "latex_nrpy_style_path = os.path.join(nrpy_dir_path,\"latex_nrpy_style.tplx\")\n",
    "#!jupyter nbconvert --to latex --template $latex_nrpy_style_path --log-level='WARN' Tutorial-IllinoisGRMHD__IllinoisGRMHD_headers.ipynb\n",
    "#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__IllinoisGRMHD_headers.tex\n",
    "#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__IllinoisGRMHD_headers.tex\n",
    "#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__IllinoisGRMHD_headers.tex\n",
    "!rm -f Tut*.out Tut*.aux Tut*.log"
   ]
  }
 ],
 "metadata": {
  "kernelspec": {
   "display_name": "Python 3",
   "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.8.3"
  }
 },
 "nbformat": 4,
 "nbformat_minor": 2
}