{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "\n", "# Tutorial-IllinoisGRMHD: MoL_registration.C\n", "\n", "## Authors: Leo Werneck & Zach Etienne\n", "\n", "**This module is currently under development**\n", "\n", "## This tutorial module explains the Method of Lines (MoL) options 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": [ "\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](#mol_registration__c): **`MoL_registration.C`**\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": [ "\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__MoL_registration__C = os.path.join(IGM_src_dir_path,\"MoL_registration.C\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Step 1: Introduction \\[Back to [top](#toc)\\]\n", "$$\\label{introduction}$$" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Step 2: `MoL_registration.C` \\[Back to [top](#toc)\\]\n", "$$\\label{mol_registration__c}$$" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Writing ../src/MoL_registration.C\n" ] } ], "source": [ "%%writefile $outfile_path__MoL_registration__C\n", "//--------------------------------------------------------------------------\n", "// Register with the time stepper\n", "// (MoL thorn, found in arrangements/CactusBase/MoL)\n", "// To understand this, read documentation in arrangements/CactusBase/MoL/doc\n", "//--------------------------------------------------------------------------\n", "\n", "#include \"cctk.h\"\n", "#include \n", "#include \n", "#include \n", "#include \"cctk_Parameters.h\"\n", "#include \"cctk_Arguments.h\"\n", "\n", "#include \"Symmetry.h\"\n", "\n", "extern \"C\" void IllinoisGRMHD_RegisterVars(CCTK_ARGUMENTS)\n", "{\n", " DECLARE_CCTK_ARGUMENTS;\n", " DECLARE_CCTK_PARAMETERS;\n", "\n", " CCTK_INT ierr = 0, group, rhs;\n", "\n", " // Register evolution & RHS gridfunction groups\n", "\n", " /* Ax and Ax_rhs */\n", " group = CCTK_GroupIndex(\"IllinoisGRMHD::em_Ax\");\n", " rhs = CCTK_GroupIndex(\"IllinoisGRMHD::em_Ax_rhs\");\n", " ierr += MoLRegisterEvolvedGroup(group, rhs);\n", "\n", " /* Ay and Ay_rhs */\n", " group = CCTK_GroupIndex(\"IllinoisGRMHD::em_Ay\");\n", " rhs = CCTK_GroupIndex(\"IllinoisGRMHD::em_Ay_rhs\");\n", " ierr += MoLRegisterEvolvedGroup(group, rhs);\n", "\n", " /* Az and Az_rhs */\n", " group = CCTK_GroupIndex(\"IllinoisGRMHD::em_Az\");\n", " rhs = CCTK_GroupIndex(\"IllinoisGRMHD::em_Az_rhs\");\n", " ierr += MoLRegisterEvolvedGroup(group, rhs);\n", "\n", " /* psi6phi and psi6phi_rhs */\n", " group = CCTK_GroupIndex(\"IllinoisGRMHD::em_psi6phi\");\n", " rhs = CCTK_GroupIndex(\"IllinoisGRMHD::em_psi6phi_rhs\");\n", " ierr += MoLRegisterEvolvedGroup(group, rhs);\n", "\n", " /* ALL OTHER EVOLVED VARIABLES (rho_star,tau,mhd_st_x,mhd_st_y,mhd_st_z) */\n", " group = CCTK_GroupIndex(\"IllinoisGRMHD::grmhd_conservatives\");\n", " rhs = CCTK_GroupIndex(\"IllinoisGRMHD::grmhd_conservatives_rhs\");\n", " ierr += MoLRegisterEvolvedGroup(group, rhs);\n", "\n", " if (ierr) CCTK_ERROR(\"Problems registering with MoL\");\n", " //***********************************************\n", "\n", " //***********************************************\n", " // Next register ADMBase variables needed by\n", " // IllinoisGRMHD as SaveAndRestore, so that\n", " // they are not set to NaN at the start of\n", " // each timestep (requiring that they be\n", " // e.g., recomputed from BSSN variables\n", " // in the BSSN solver, like Baikal or\n", " // ML_BSSN)\n", " ierr += MoLRegisterSaveAndRestoreGroup(CCTK_GroupIndex(\"admbase::lapse\"));\n", " ierr += MoLRegisterSaveAndRestoreGroup(CCTK_GroupIndex(\"admbase::shift\"));\n", " ierr += MoLRegisterSaveAndRestoreGroup(CCTK_GroupIndex(\"admbase::metric\"));\n", " ierr += MoLRegisterSaveAndRestoreGroup(CCTK_GroupIndex(\"admbase::curv\"));\n", " if (ierr) CCTK_ERROR(\"Problems registering with MoLRegisterSaveAndRestoreGroup\");\n", "}" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\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 MoL_registration.C: FAILED!\n", "Diff:\n", "2c2\n", "< // Register with the time stepper \n", "---\n", "> // Register with the time stepper\n", "20c20\n", "< \n", "---\n", "> \n", "52a53\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/MoL_registration.C\"\n", "original_IGM_file_name = \"MoL_registration-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__MoL_registration__C = !diff $original_IGM_file_path $outfile_path__MoL_registration__C\n", "\n", "if Validation__MoL_registration__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 MoL_registration.C: PASSED!\")\n", "else:\n", " # If the validation fails, we keep the original IGM source code file\n", " print(\"Validation test for MoL_registration.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__MoL_registration__C:\n", " print(diff_line)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\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__MoL_registration.pdf](Tutorial-IllinoisGRMHD__MoL_registration.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__MoL_registration.ipynb\n", "#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__MoL_registration.tex\n", "#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__MoL_registration.tex\n", "#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__MoL_registration.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 }