This module is currently under development
IllinoisGRMHD
codes. This module will likely be absorbed by another one once we finish documenting the code.¶We will now use the cmdline_helper.py NRPy+ module to create the source directory within the IllinoisGRMHD
NRPy+ directory, if it does not exist yet.
# 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")
%%writefile $outfile_path__InitSymBound__C
/*
Set the symmetries for the IllinoisGRMHD variables
*/
#include "cctk.h"
#include <cstdio>
#include <cstdlib>
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"
#include "Symmetry.h"
#include "IllinoisGRMHD_headers.h"
extern "C" void IllinoisGRMHD_InitSymBound(CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS;
DECLARE_CCTK_PARAMETERS;
if( ( CCTK_EQUALS(Matter_BC,"frozen") && !CCTK_EQUALS(EM_BC,"frozen") ) ||
( !CCTK_EQUALS(Matter_BC,"frozen") && CCTK_EQUALS(EM_BC,"frozen") ) )
CCTK_VError(VERR_DEF_PARAMS,"If Matter_BC or EM_BC is set to FROZEN, BOTH must be set to frozen!");
if ((cctk_nghostzones[0]<3 || cctk_nghostzones[1]<3 || cctk_nghostzones[2]<3))
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]);
if(cctk_iteration==0) {
CCTK_VInfo(CCTK_THORNSTRING,"Setting Symmetry = %s... at iteration = %d",Symmetry,cctk_iteration);
int sym[3];
if(CCTK_EQUALS(Symmetry,"none")) {
/* FIRST SET NO SYMMETRY OPTION */
sym[0] = 1; sym[1] = 1; sym[2] = 1;
SetCartSymGN(cctkGH,sym,"IllinoisGRMHD::grmhd_conservatives");
SetCartSymGN(cctkGH,sym,"IllinoisGRMHD::em_Ax");
SetCartSymGN(cctkGH,sym,"IllinoisGRMHD::em_Ay");
SetCartSymGN(cctkGH,sym,"IllinoisGRMHD::em_Az");
SetCartSymGN(cctkGH,sym,"IllinoisGRMHD::em_psi6phi");
SetCartSymGN(cctkGH,sym,"IllinoisGRMHD::grmhd_primitives_allbutBi");
} else if(CCTK_EQUALS(Symmetry,"equatorial")) {
/* THEN SET EQUATORIAL SYMMETRY OPTION */
// Set default to no symmetry, which is correct for scalars and most vectors:
sym[0] = 1; sym[1] = 1; sym[2] = 1;
SetCartSymGN(cctkGH,sym,"IllinoisGRMHD::grmhd_conservatives");
// Don't worry about the wrong sym values since A_{\mu} is staggered
// and we're going to impose the symmetry separately
SetCartSymGN(cctkGH,sym,"IllinoisGRMHD::em_Ax");
SetCartSymGN(cctkGH,sym,"IllinoisGRMHD::em_Ay");
SetCartSymGN(cctkGH,sym,"IllinoisGRMHD::em_Az");
SetCartSymGN(cctkGH,sym,"IllinoisGRMHD::em_psi6phi");
SetCartSymGN(cctkGH,sym,"IllinoisGRMHD::grmhd_primitives_allbutBi");
// Then set unstaggered B field variables
sym[2] = -Sym_Bz;
SetCartSymVN(cctkGH, sym,"IllinoisGRMHD::Bx");
SetCartSymVN(cctkGH, sym,"IllinoisGRMHD::By");
sym[2] = Sym_Bz;
SetCartSymVN(cctkGH, sym,"IllinoisGRMHD::Bz");
sym[2] = -1;
SetCartSymVN(cctkGH, sym,"IllinoisGRMHD::mhd_st_z");
SetCartSymVN(cctkGH, sym,"IllinoisGRMHD::vz");
} else {
CCTK_VError(VERR_DEF_PARAMS,"IllinoisGRMHD_initsymbound: Should not be here; picked an impossible symmetry.");
}
}
}
Writing ../src/InitSymBound.C
# 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
!wget -O $original_IGM_file_path $original_IGM_file_url
# Perform validation
Validation__InitSymBound__C = !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
!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)
Validation test for InitSymBound.C: FAILED! Diff: 18c18 < if( ( CCTK_EQUALS(Matter_BC,"frozen") && !CCTK_EQUALS(EM_BC,"frozen") ) || --- > if( ( CCTK_EQUALS(Matter_BC,"frozen") && !CCTK_EQUALS(EM_BC,"frozen") ) || 67a68 >
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 (Note that clicking on this link may not work; you may need to open the PDF file through another means).
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
!rm -f Tut*.out Tut*.aux Tut*.log