{
"cells": [
{
"cell_type": "markdown",
"metadata": {},
"source": [
"\n",
"\n",
"\n",
"# Tutorial-IllinoisGRMHD: compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.C\n",
"\n",
"## Authors: Leo Werneck & Zach Etienne\n",
"\n",
"**This module is currently under development**\n",
"\n",
"## In this tutorial module we explain how we evaluate $\\partial_{t}\\tilde{\\tau}$, terms that depend on the extrinsic curvature, $K_{ij}$, and the energy-momentum tensor $T^{\\mu\\nu}$\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](#compute_tau_rhs_extrinsic_curvature_terms_and_tupmunu__c): **`compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.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__compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu__C = os.path.join(IGM_src_dir_path,\"compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.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: `compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.C` \\[Back to [top](#toc)\\]\n",
"$$\\label{compute_tau_rhs_extrinsic_curvature_terms_and_tupmunu__c}$$"
]
},
{
"cell_type": "code",
"execution_count": 2,
"metadata": {},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"Writing ../src/compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.C\n"
]
}
],
"source": [
"%%writefile $outfile_path__compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu__C\n",
"static void compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu\n",
"(const cGH *cctkGH,const int *cctk_lsh,const int *cctk_nghostzones,CCTK_REAL *dX,CCTK_REAL **metric,gf_and_gz_struct *prims,\n",
" CCTK_REAL **TUPmunu,eos_struct &eos,CCTK_REAL Gamma_th,\n",
" CCTK_REAL *gupxy,CCTK_REAL *gupxz,CCTK_REAL *gupyz,\n",
" CCTK_REAL *kxx,CCTK_REAL *kxy,CCTK_REAL *kxz,CCTK_REAL *kyy,CCTK_REAL *kyz,CCTK_REAL *kzz,\n",
" CCTK_REAL *tau_rhs) {\n",
"\n",
" // These loop extents must be consistent with add_fluxes_and_source_terms_to_hydro_rhss(), since we use TUPmunu there as well.\n",
"#pragma omp parallel for\n",
" for(int k=cctk_nghostzones[2];k\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 compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.C: FAILED!\n",
"Diff:\n",
"3c3\n",
"< CCTK_REAL **TUPmunu,eos_struct &eos,\n",
"---\n",
"> CCTK_REAL **TUPmunu,eos_struct &eos,CCTK_REAL Gamma_th,\n",
"20c20\n",
"< // forces us to read in only once, storing in local variables \n",
"---\n",
"> // forces us to read in only once, storing in local variables\n",
"43c43\n",
"< \n",
"---\n",
"> \n",
"50c50\n",
"< compute_smallba_b2_and_u_i_over_u0_psi4(METRIC,METRIC_LAP_PSI4,U,u0L,ONE_OVER_LAPSE_SQRT_4PI, \n",
"---\n",
"> compute_smallba_b2_and_u_i_over_u0_psi4(METRIC,METRIC_LAP_PSI4,U,u0L,ONE_OVER_LAPSE_SQRT_4PI,\n",
"56c56\n",
"< CCTK_REAL h=0, P_cold,eps_cold,dPcold_drho,eps_th,gamma_cold; /* <- Note that in setting h, we need to define several \n",
"---\n",
"> CCTK_REAL h=0, P_cold,eps_cold,dPcold_drho,eps_th,Gamma_cold; /* <- Note that in setting h, we need to define several\n",
"59c59\n",
"< compute_P_cold__eps_cold__dPcold_drho__eps_th__h__gamma_cold(U,eos,P_cold,eps_cold,dPcold_drho,eps_th,h,gamma_cold);\n",
"---\n",
"> compute_P_cold__eps_cold__dPcold_drho__eps_th__h__Gamma_cold(U,eos,Gamma_th,P_cold,eps_cold,dPcold_drho,eps_th,h,Gamma_cold);\n",
"99c99\n",
"< \n",
"---\n",
"> \n",
"103c103\n",
"< TUP[0][0] * (SQR(METRIC[SHIFTX])*KxxL + SQR(METRIC[SHIFTY])*KyyL + SQR(METRIC[SHIFTZ])*KzzL + \n",
"---\n",
"> TUP[0][0] * (SQR(METRIC[SHIFTX])*KxxL + SQR(METRIC[SHIFTY])*KyyL + SQR(METRIC[SHIFTZ])*KzzL +\n",
"105c105\n",
"< 2.0*(TUP[0][1]*METRIC[SHIFTX]*KxxL + TUP[0][2]*METRIC[SHIFTY]*KyyL + TUP[0][3]*METRIC[SHIFTZ]*KzzL + \n",
"---\n",
"> 2.0*(TUP[0][1]*METRIC[SHIFTX]*KxxL + TUP[0][2]*METRIC[SHIFTY]*KyyL + TUP[0][3]*METRIC[SHIFTZ]*KzzL +\n",
"109c109\n",
"< TUP[1][1]*KxxL + TUP[2][2]*KyyL + TUP[3][3]*KzzL + \n",
"---\n",
"> TUP[1][1]*KxxL + TUP[2][2]*KyyL + TUP[3][3]*KzzL +\n",
"119a120\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/compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.C\"\n",
"original_IGM_file_name = \"compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu-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__compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu__C = !diff $original_IGM_file_path $outfile_path__compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu__C\n",
"\n",
"if Validation__compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu__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 compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.C: PASSED!\")\n",
"else:\n",
" # If the validation fails, we keep the original IGM source code file\n",
" print(\"Validation test for compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.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__compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu__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__compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.pdf](Tutorial-IllinoisGRMHD__compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.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__compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.ipynb\n",
"#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.tex\n",
"#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.tex\n",
"#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__compute_tau_rhs_extrinsic_curvature_terms_and_TUPmunu.tex\n",
"!rm -f Tut*.out Tut*.aux Tut*.log"
]
}
],
"metadata": {
"kernelspec": {
"display_name": "Python 3",
"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.8.3"
}
},
"nbformat": 4,
"nbformat_minor": 2
}