{
 "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: InitSymBound.C\n",
    "\n",
    "## Authors: Leo Werneck & Zach Etienne\n",
    "\n",
    "<font color='red'>**This module is currently under development**</font>\n",
    "\n",
    "## This tutorial notebook explains the symmetry conditions used by the `IllinoisGRMHD` codes. \n",
    "\n",
    "### Note: This module will likely be absorbed by another one once 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='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__InitSymBound__C = os.path.join(IGM_src_dir_path,\"InitSymBound.C\")"
   ]
  },
  {
   "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='a_i_rhs_no_gauge_terms__c'></a>\n",
    "\n",
    "# Step 2: `A_i_rhs_no_gauge_terms.C` \\[Back to [top](#toc)\\]\n",
    "$$\\label{a_i_rhs_no_gauge_terms__c}$$"
   ]
  },
  {
   "cell_type": "code",
   "execution_count": 2,
   "metadata": {},
   "outputs": [
    {
     "name": "stdout",
     "output_type": "stream",
     "text": [
      "Writing ../src/InitSymBound.C\n"
     ]
    }
   ],
   "source": [
    "%%writefile $outfile_path__InitSymBound__C\n",
    "/*\n",
    "  Set the symmetries for the IllinoisGRMHD variables\n",
    "*/\n",
    "\n",
    "#include \"cctk.h\"\n",
    "#include <cstdio>\n",
    "#include <cstdlib>\n",
    "#include \"cctk_Arguments.h\"\n",
    "#include \"cctk_Parameters.h\"\n",
    "#include \"Symmetry.h\"\n",
    "#include \"IllinoisGRMHD_headers.h\"\n",
    "\n",
    "extern \"C\" void IllinoisGRMHD_InitSymBound(CCTK_ARGUMENTS)\n",
    "{\n",
    "  DECLARE_CCTK_ARGUMENTS;\n",
    "  DECLARE_CCTK_PARAMETERS;\n",
    "\n",
    "  if( ( CCTK_EQUALS(Matter_BC,\"frozen\") && !CCTK_EQUALS(EM_BC,\"frozen\") ) ||\n",
    "      ( !CCTK_EQUALS(Matter_BC,\"frozen\") && CCTK_EQUALS(EM_BC,\"frozen\") ) )\n",
    "    CCTK_VError(VERR_DEF_PARAMS,\"If Matter_BC or EM_BC is set to FROZEN, BOTH must be set to frozen!\");\n",
    "\n",
    "  if ((cctk_nghostzones[0]<3 || cctk_nghostzones[1]<3 || cctk_nghostzones[2]<3))\n",
    "    CCTK_VError(VERR_DEF_PARAMS,\"ERROR: The version of PPM in this thorn requires 3 ghostzones. You only have (%d,%d,%d) ghostzones!\",cctk_nghostzones[0],cctk_nghostzones[1],cctk_nghostzones[2]);\n",
    "\n",
    "  if(cctk_iteration==0) {\n",
    "    CCTK_VInfo(CCTK_THORNSTRING,\"Setting Symmetry = %s... at iteration = %d\",Symmetry,cctk_iteration);\n",
    "\n",
    "    int sym[3];\n",
    "\n",
    "    if(CCTK_EQUALS(Symmetry,\"none\")) {\n",
    "      /* FIRST SET NO SYMMETRY OPTION */\n",
    "      sym[0] = 1; sym[1] = 1; sym[2] = 1;\n",
    "      SetCartSymGN(cctkGH,sym,\"IllinoisGRMHD::grmhd_conservatives\");\n",
    "      SetCartSymGN(cctkGH,sym,\"IllinoisGRMHD::em_Ax\");\n",
    "      SetCartSymGN(cctkGH,sym,\"IllinoisGRMHD::em_Ay\");\n",
    "      SetCartSymGN(cctkGH,sym,\"IllinoisGRMHD::em_Az\");\n",
    "      SetCartSymGN(cctkGH,sym,\"IllinoisGRMHD::em_psi6phi\");\n",
    "      SetCartSymGN(cctkGH,sym,\"IllinoisGRMHD::grmhd_primitives_allbutBi\");\n",
    "    } else if(CCTK_EQUALS(Symmetry,\"equatorial\")) {\n",
    "      /* THEN SET EQUATORIAL SYMMETRY OPTION */\n",
    "      // Set default to no symmetry, which is correct for scalars and most vectors:\n",
    "      sym[0] = 1; sym[1] = 1; sym[2] = 1;\n",
    "      SetCartSymGN(cctkGH,sym,\"IllinoisGRMHD::grmhd_conservatives\");\n",
    "      // Don't worry about the wrong sym values since A_{\\mu} is staggered\n",
    "      // and we're going to impose the symmetry separately\n",
    "      SetCartSymGN(cctkGH,sym,\"IllinoisGRMHD::em_Ax\");\n",
    "      SetCartSymGN(cctkGH,sym,\"IllinoisGRMHD::em_Ay\");\n",
    "      SetCartSymGN(cctkGH,sym,\"IllinoisGRMHD::em_Az\");\n",
    "      SetCartSymGN(cctkGH,sym,\"IllinoisGRMHD::em_psi6phi\");\n",
    "\n",
    "      SetCartSymGN(cctkGH,sym,\"IllinoisGRMHD::grmhd_primitives_allbutBi\");\n",
    "\n",
    "      // Then set unstaggered B field variables\n",
    "      sym[2] = -Sym_Bz;\n",
    "      SetCartSymVN(cctkGH, sym,\"IllinoisGRMHD::Bx\");\n",
    "      SetCartSymVN(cctkGH, sym,\"IllinoisGRMHD::By\");\n",
    "      sym[2] = Sym_Bz;\n",
    "      SetCartSymVN(cctkGH, sym,\"IllinoisGRMHD::Bz\");\n",
    "\n",
    "      sym[2] = -1;\n",
    "      SetCartSymVN(cctkGH, sym,\"IllinoisGRMHD::mhd_st_z\");\n",
    "      SetCartSymVN(cctkGH, sym,\"IllinoisGRMHD::vz\");\n",
    "    } else {\n",
    "      CCTK_VError(VERR_DEF_PARAMS,\"IllinoisGRMHD_initsymbound: Should not be here; picked an impossible symmetry.\");\n",
    "    }\n",
    "  }\n",
    "}\n",
    "\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 InitSymBound.C: FAILED!\n",
      "Diff:\n",
      "18c18\n",
      "<   if( ( CCTK_EQUALS(Matter_BC,\"frozen\") && !CCTK_EQUALS(EM_BC,\"frozen\") ) || \n",
      "---\n",
      ">   if( ( CCTK_EQUALS(Matter_BC,\"frozen\") && !CCTK_EQUALS(EM_BC,\"frozen\") ) ||\n",
      "67a68\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/InitSymBound.C\"\n",
    "original_IGM_file_name = \"InitSymBound-original.C\"\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__InitSymBound__C  = !diff $original_IGM_file_path $outfile_path__InitSymBound__C\n",
    "\n",
    "if Validation__InitSymBound__C == []:\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 InitSymBound.C: PASSED!\")\n",
    "else:\n",
    "    # If the validation fails, we keep the original IGM source code file\n",
    "    print(\"Validation test for InitSymBound.C: 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__InitSymBound__C:\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__InitSymBound.pdf](Tutorial-IllinoisGRMHD__InitSymBound.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__InitSymBound.ipynb\n",
    "#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__InitSymBound.tex\n",
    "#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__InitSymBound.tex\n",
    "#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__InitSymBound.tex\n",
    "!rm -f Tut*.out Tut*.aux Tut*.log"
   ]
  }
 ],
 "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.1"
  }
 },
 "nbformat": 4,
 "nbformat_minor": 4
}