#!/usr/bin/env python # coding: utf-8 # # # # # Tutorial-IllinoisGRMHD: postpostinitial__set_symmetries__copy_timelevels.C # # ## Authors: Leo Werneck & Zach Etienne # # **This module is currently under development** # # ## This tutorial module explains tasks that are performed after the set up of initial data. # # ### 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)). # # # # Table of Contents # $$\label{toc}$$ # # This module is organized as follows # # 0. [Step 0](#src_dir): **Source directory creation** # 1. [Step 1](#introduction): **Introduction** # 1. [Step 2](#postpostinitial__set_symmetries__copy_timelevels__c): **`postpostinitial__set_symmetries__copy_timelevels.C`** # 1. [Step n-1](#code_validation): **Code validation** # 1. [Step n](#latex_pdf_output): **Output this notebook to $\LaTeX$-formatted PDF file** # # # # 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__postpostinitial__set_symmetries__copy_timelevels__C = os.path.join(IGM_src_dir_path,"postpostinitial__set_symmetries__copy_timelevels.C") # # # # Step 1: Introduction \[Back to [top](#toc)\] # $$\label{introduction}$$ # # # # Step 2: `postpostinitial__set_symmetries__copy_timelevels.C` \[Back to [top](#toc)\] # $$\label{postpostinitial__set_symmetries__copy_timelevels__c}$$ # In[2]: get_ipython().run_cell_magic('writefile', '$outfile_path__postpostinitial__set_symmetries__copy_timelevels__C', '//-------------------------------------------------\n// Stuff to run right after initial data is set up\n//-------------------------------------------------\n\n#include "cctk.h"\n#include \n#include \n#include "cctk_Arguments.h"\n#include "cctk_Functions.h"\n#include "cctk_Parameters.h"\n#include "Symmetry.h"\n#include "IllinoisGRMHD_headers.h"\n#include "IllinoisGRMHD_EoS_lowlevel_functs.C"\n\nextern "C" void IllinoisGRMHD_PostPostInitial_Set_Symmetries__Copy_Timelevels(CCTK_ARGUMENTS) {\n DECLARE_CCTK_ARGUMENTS;\n DECLARE_CCTK_PARAMETERS;\n\n /**********************************\n * Piecewise Polytropic EOS Patch *\n * Printing the EOS table *\n **********************************/\n /*\n * The short piece of code below takes care\n * of initializing the EOS parameters.\n * Please refer to the "inlined_functions.C"\n * source file for the documentation on the\n * function.\n */\n eos_struct eos;\n initialize_EOS_struct_from_input(eos);\n\n /* For diagnostic and user convenience purposes, we print\n * out the EOS parameters (rho_ppoly_tab, K_ppoly_tab,\n * and Gamma_ppoly_tab) at t=0.\n */\n if(cctk_iteration==0 && (int)GetRefinementLevel(cctkGH)==0) { print_EOS_table(eos); }\n\n if(Gamma_th<0)\n CCTK_VError(VERR_DEF_PARAMS,"ERROR. Default Gamma_th (=-1) detected. You must set Gamma_th to the appropriate value in your initial data thorn, or your .par file!\\n");\n\n //For emfields, we assume that you\'ve set Bx, By, Bz (the UN-tilded B^i\'s)\n // or Ax, Ay, Az (if using constrained transport scheme of Del Zanna)\n\n if(CCTK_EQUALS(Symmetry,"equatorial")) {\n // SET SYMMETRY GHOSTZONES ON ALL CONSERVATIVE AND PRIMIIVE VARIABLES!\n int ierr;\n ierr=CartSymGN(cctkGH,"IllinoisGRMHD::grmhd_conservatives"); if(ierr!=0) CCTK_VError(VERR_DEF_PARAMS,"Microsoft error code #1874109358120048. Grep it in the source code");\n ierr=CartSymGN(cctkGH,"IllinoisGRMHD::grmhd_primitives_allbutBi"); if(ierr!=0) CCTK_VError(VERR_DEF_PARAMS,"Microsoft error code #1874109358120049. Grep it in the source code");\n\n // Finish up by setting symmetry ghostzones on Bx, By, Bz, and their staggered variants.\n CCTK_REAL gridfunc_syms_Bx[3] = {-1, 1,-Sym_Bz};\n IllinoisGRMHD_set_symmetry_gzs_staggered(cctkGH,cctk_lsh,x,y,z, Bx , gridfunc_syms_Bx,0,0,0);\n IllinoisGRMHD_set_symmetry_gzs_staggered(cctkGH,cctk_lsh,x,y,z, Bx_stagger, gridfunc_syms_Bx,1,0,0);\n CCTK_REAL gridfunc_syms_By[3] = { 1,-1,-Sym_Bz};\n IllinoisGRMHD_set_symmetry_gzs_staggered(cctkGH,cctk_lsh,x,y,z, By , gridfunc_syms_Bx,0,0,0);\n IllinoisGRMHD_set_symmetry_gzs_staggered(cctkGH,cctk_lsh,x,y,z, By_stagger, gridfunc_syms_By,0,1,0);\n CCTK_REAL gridfunc_syms_Bz[3] = { 1, 1, Sym_Bz};\n IllinoisGRMHD_set_symmetry_gzs_staggered(cctkGH,cctk_lsh,x,y,z, Bz , gridfunc_syms_Bz,0,0,0);\n IllinoisGRMHD_set_symmetry_gzs_staggered(cctkGH,cctk_lsh,x,y,z, Bz_stagger, gridfunc_syms_Bz,0,0,1);\n\n CCTK_REAL gridfunc_syms_psi6phi[3] = { 1, 1, 1};\n IllinoisGRMHD_set_symmetry_gzs_staggered(cctkGH,cctk_lsh,x,y,z,psi6phi , gridfunc_syms_psi6phi,1,1,1);\n CCTK_REAL gridfunc_syms_Ax[3] = {-1, 1, Sym_Bz};\n IllinoisGRMHD_set_symmetry_gzs_staggered(cctkGH,cctk_lsh,x,y,z, Ax , gridfunc_syms_Ax,0,1,1);\n CCTK_REAL gridfunc_syms_Ay[3] = { 1,-1, Sym_Bz};\n IllinoisGRMHD_set_symmetry_gzs_staggered(cctkGH,cctk_lsh,x,y,z, Ay , gridfunc_syms_Ay,1,0,1);\n CCTK_REAL gridfunc_syms_Az[3] = { 1, 1,-Sym_Bz};\n IllinoisGRMHD_set_symmetry_gzs_staggered(cctkGH,cctk_lsh,x,y,z, Az , gridfunc_syms_Az,1,1,0);\n }\n\n\n //------------------------------------------------------------------\n // FILL _p AND _p_p TIMELEVELS. Probably don\'t need to do this if\n // Carpet::init_fill_timelevels=yes and\n // MoL::initial_data_is_crap = yes\n // NOTE: We don\'t fill metric data here.\n // FIXME: Do we really need this?\n\n#pragma omp parallel for\n for(int k=0;k # # # 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/postpostinitial__set_symmetries__copy_timelevels.C" original_IGM_file_name = "postpostinitial__set_symmetries__copy_timelevels-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__postpostinitial__set_symmetries__copy_timelevels__C = get_ipython().getoutput('diff $original_IGM_file_path $outfile_path__postpostinitial__set_symmetries__copy_timelevels__C') if Validation__postpostinitial__set_symmetries__copy_timelevels__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 postpostinitial__set_symmetries__copy_timelevels.C: PASSED!") else: # If the validation fails, we keep the original IGM source code file print("Validation test for postpostinitial__set_symmetries__copy_timelevels.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__postpostinitial__set_symmetries__copy_timelevels__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__postpostinitial__set_symmetries__copy_timelevels.pdf](Tutorial-IllinoisGRMHD__postpostinitial__set_symmetries__copy_timelevels.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__postpostinitial__set_symmetries__copy_timelevels.ipynb #!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__postpostinitial__set_symmetries__copy_timelevels.tex #!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__postpostinitial__set_symmetries__copy_timelevels.tex #!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__postpostinitial__set_symmetries__copy_timelevels.tex get_ipython().system('rm -f Tut*.out Tut*.aux Tut*.log')