#!/usr/bin/env python # coding: utf-8 # # # # # Tutorial-IllinoisGRMHD: InitSymBound.C # # ## Authors: Leo Werneck & Zach Etienne # # **This module is currently under development** # # ## This tutorial notebook explains the symmetry conditions used by the `IllinoisGRMHD` codes. # # ### Note: This module will likely be absorbed by another one once we finish documenting the code. # # ### Required and recommended citations: # # * **(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)). # * **(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)). # * **(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)). # # # # Step 0: Source directory creation \[Back to [top](#toc)\] # $$\label{src_dir}$$ # # 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. # In[1]: # Step 0: Creation of the IllinoisGRMHD source directory # Step 0a: Add NRPy's directory to the path # https://stackoverflow.com/questions/16780014/import-file-from-parent-directory import os,sys nrpy_dir_path = os.path.join("..","..") if nrpy_dir_path not in sys.path: sys.path.append(nrpy_dir_path) # Step 0b: Load up cmdline_helper and create the directory import cmdline_helper as cmd IGM_src_dir_path = os.path.join("..","src") cmd.mkdir(IGM_src_dir_path) # Step 0c: Create the output file path outfile_path__InitSymBound__C = os.path.join(IGM_src_dir_path,"InitSymBound.C") # # # # Step 1: Introduction \[Back to [top](#toc)\] # $$\label{introduction}$$ # # # # Step 2: `A_i_rhs_no_gauge_terms.C` \[Back to [top](#toc)\] # $$\label{a_i_rhs_no_gauge_terms__c}$$ # In[2]: get_ipython().run_cell_magic('writefile', '$outfile_path__InitSymBound__C', '/*\n Set the symmetries for the IllinoisGRMHD variables\n*/\n\n#include "cctk.h"\n#include \n#include \n#include "cctk_Arguments.h"\n#include "cctk_Parameters.h"\n#include "Symmetry.h"\n#include "IllinoisGRMHD_headers.h"\n\nextern "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') # # # # Step n-1: Code validation \[Back to [top](#toc)\] # $$\label{code_validation}$$ # # First we download the original `IllinoisGRMHD` source code and then compare it to the source code generated by this tutorial notebook. # In[3]: # Verify if the code generated by this tutorial module # matches the original IllinoisGRMHD source code # First download the original IllinoisGRMHD source code import urllib from os import path original_IGM_file_url = "https://bitbucket.org/zach_etienne/wvuthorns/raw/5611b2f0b17135538c9d9d17c7da062abe0401b6/IllinoisGRMHD/src/InitSymBound.C" original_IGM_file_name = "InitSymBound-original.C" original_IGM_file_path = os.path.join(IGM_src_dir_path,original_IGM_file_name) # Then download the original IllinoisGRMHD source code # We try it here in a couple of ways in an attempt to keep # the code more portable try: original_IGM_file_code = urllib.request.urlopen(original_IGM_file_url).read().decode("utf-8") # Write down the file the original IllinoisGRMHD source code with open(original_IGM_file_path,"w") as file: file.write(original_IGM_file_code) except: try: original_IGM_file_code = urllib.urlopen(original_IGM_file_url).read().decode("utf-8") # Write down the file the original IllinoisGRMHD source code with open(original_IGM_file_path,"w") as file: file.write(original_IGM_file_code) except: # If all else fails, hope wget does the job get_ipython().system('wget -O $original_IGM_file_path $original_IGM_file_url') # Perform validation Validation__InitSymBound__C = get_ipython().getoutput('diff $original_IGM_file_path $outfile_path__InitSymBound__C') if Validation__InitSymBound__C == []: # If the validation passes, we do not need to store the original IGM source code file get_ipython().system('rm $original_IGM_file_path') print("Validation test for InitSymBound.C: PASSED!") else: # If the validation fails, we keep the original IGM source code file print("Validation test for InitSymBound.C: FAILED!") # We also print out the difference between the code generated # in this tutorial module and the original IGM source code print("Diff:") for diff_line in Validation__InitSymBound__C: print(diff_line) # # # # Step n: Output this notebook to $\LaTeX$-formatted PDF file \[Back to [top](#toc)\] # $$\label{latex_pdf_output}$$ # # 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 # [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). # In[4]: latex_nrpy_style_path = os.path.join(nrpy_dir_path,"latex_nrpy_style.tplx") #!jupyter nbconvert --to latex --template $latex_nrpy_style_path --log-level='WARN' Tutorial-IllinoisGRMHD__InitSymBound.ipynb #!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__InitSymBound.tex #!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__InitSymBound.tex #!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__InitSymBound.tex get_ipython().system('rm -f Tut*.out Tut*.aux Tut*.log')