{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "\n", "# Tutorial-IllinoisGRMHD: Configuration files\n", "\n", "## Authors: Leo Werneck & Zach Etienne\n", "\n", "**This module is currently under development**\n", "\n", "## In this tutorial module we explain the files necessary to configure `IllinoisGRMHD` so that it can be properly used by the Einstein Toolkit.\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](#make_code_defn): **`make.code.defn`**\n", "1. [Step 3](#configuration__ccl): **`configuration.ccl`**\n", "1. [Step 4](#interface__ccl): **`interface.ccl`**\n", "1. [Step 5](#param__ccl): **`param.ccl`**\n", "1. [Step 6](#schedule__ccl): **`schedule.ccl`**\n", "1. [Step 7](#code_validation__txt): **`code_validation.txt`**\n", "1. [Step 8](#code_validation): **Code validation**\n", " 1. [Step 8.a](#code_validation__make_code_defn): *`make.code.defn`*\n", " 1. [Step 8.b](#code_validation__configuration__ccl): *`configuration.ccl`*\n", " 1. [Step 8.c](#code_validation__interface__ccl): *`interface.ccl`*\n", " 1. [Step 8.d](#code_validation__param__ccl): *`param.ccl`*\n", " 1. [Step 8.e](#code_validation__schedule__ccl): *`schedule.ccl`*\n", " 1. [Step 8.f](#code_validation__code_validation__txt): *`code_validation.txt`*\n", "1. [Step 9](#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: For this tutorial module we will also need IllinoisGRMHD's main directory path\n", "IGM_main_dir_path = \"..\"\n", "\n", "# Step 0d: Create the output file path\n", "outfile_path__make_code_defn = os.path.join(IGM_src_dir_path, \"make.code.defn\")\n", "outfile_path__configuration__ccl = os.path.join(IGM_main_dir_path,\"configuration.ccl\")\n", "outfile_path__interface__ccl = os.path.join(IGM_main_dir_path,\"interface.ccl\")\n", "outfile_path__param__ccl = os.path.join(IGM_main_dir_path,\"param.ccl\")\n", "outfile_path__schedule__ccl = os.path.join(IGM_main_dir_path,\"schedule.ccl\")\n", "outfile_path__code_validation__txt = \"code_validation.txt\"" ] }, { "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: `make.code.defn` \\[Back to [top](#toc)\\]\n", "$$\\label{make_code_defn}$$" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Writing ../src/make.code.defn\n" ] } ], "source": [ "%%writefile $outfile_path__make_code_defn\n", "# Main make.code.defn file for thorn IllinoisGRMHD\n", "\n", "# Source files in this directory\n", "SRCS = InitSymBound.C MoL_registration.C \\\n", " \\\n", " postpostinitial__set_symmetries__copy_timelevels.C \\\n", " \\\n", " convert_ADM_to_BSSN__enforce_detgtij_eq_1__and_compute_gtupij.C \\\n", " \\\n", " driver_evaluate_MHD_rhs.C \\\n", " compute_B_and_Bstagger_from_A.C \\\n", " driver_conserv_to_prims.C \\\n", " symmetry__set_gzs_staggered_gfs.C \\\n", " \\\n", " outer_boundaries.C\n", "\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Step 3: `configuration.ccl` \\[Back to [top](#toc)\\]\n", "$$\\label{configuration__ccl}$$" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Writing ../configuration.ccl\n" ] } ], "source": [ "%%writefile $outfile_path__configuration__ccl\n", "# Configuration definition for thorn IllinoisGRMHD\n", "REQUIRES EOS_Omni Boundary CartGrid3D SpaceMask\n", "\n", "PROVIDES IllinoisGRMHD\n", "{\n", "}\n", "\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Step 4: `interface.ccl` \\[Back to [top](#toc)\\]\n", "$$\\label{interface__ccl}$$" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Writing ../interface.ccl\n" ] } ], "source": [ "%%writefile $outfile_path__interface__ccl\n", "# Interface definition for thorn IllinoisGRMHD\n", "\n", "implements: IllinoisGRMHD\n", "inherits: ADMBase, Boundary, SpaceMask, Tmunubase, HydroBase, grid\n", "\n", "includes header: IllinoisGRMHD_headers.h in IllinoisGRMHD_headers.h\n", "\n", "USES INCLUDE: Symmetry.h\n", "\n", "public:\n", "\n", "#vvvvvvvv EVOLVED VARIABLES vvvvvvvv#\n", "#cctk_real grmhd_conservatives type = GF TAGS='prolongation=\"Lag3\"' Timelevels=3\n", "cctk_real grmhd_conservatives type = GF Timelevels=3\n", "{\n", " rho_star,tau,mhd_st_x,mhd_st_y,mhd_st_z # Note that st = Stilde, as mhd_st_i = \\tilde{S}_i.\n", "} \"Evolved mhd variables\"\n", "\n", "# These variables are semi-staggered:\n", "# Ax is defined on the semi-staggered grid (i,j+1/2,k+1/2)\n", "# WARNING: WILL NOT WORK PROPERLY WITHOUT SEMI-STAGGERED PROLONGATION/RESTRICTION:\n", "cctk_real em_Ax type = GF Timelevels=3 tags='Prolongation=\"STAGGER011\"'\n", "{\n", " Ax\n", "} \"x-component of the vector potential, evolved when constrained_transport_scheme==3\"\n", "\n", "# Ay is defined on the semi-staggered grid (i+1/2,j,k+1/2)\n", "# WARNING: WILL NOT WORK PROPERLY WITHOUT SEMI-STAGGERED PROLONGATION/RESTRICTION:\n", "cctk_real em_Ay type = GF Timelevels=3 tags='Prolongation=\"STAGGER101\"'\n", "{\n", " Ay\n", "} \"y-component of the vector potential, evolved when constrained_transport_scheme==3\"\n", "# WARNING: WILL NOT WORK PROPERLY WITHOUT SEMI-STAGGERED PROLONGATION/RESTRICTION:\n", "\n", "# Az is defined on the semi-staggered grid (i+1/2,j+1/2,k)\n", "cctk_real em_Az type = GF Timelevels=3 tags='Prolongation=\"STAGGER110\"'\n", "{\n", " Az\n", "} \"z-component of the vector potential, evolved when constrained_transport_scheme==3\"\n", "\n", "# psi6phi is defined on the staggered grid (i+1/2,j+1/2,k+1/2)\n", "# WARNING: WILL NOT WORK PROPERLY WITHOUT FULLY-STAGGERED PROLONGATION/RESTRICTION:\n", "#\n", "cctk_real em_psi6phi type = GF Timelevels=3 tags='Prolongation=\"STAGGER111\"'\n", "{\n", " psi6phi\n", "} \"sqrt{gamma} Phi, where Phi is the em scalar potential\"\n", "#^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^#\n", "\n", "\n", "#vvvvvvv PRIMITIVE VARIABLES vvvvvvv#\n", "# TODO: split into groups with well defined symmetry properties: (rho_b, P, u0), (vx,vy,vz)\n", "cctk_real grmhd_primitives_allbutBi type = GF TAGS='InterpNumTimelevels=1 prolongation=\"none\"'\n", "{\n", " rho_b,P,vx,vy,vz\n", "} \"Primitive variables density, pressure, and components of three velocity v^i. Note that v^i is defined in terms of 4-velocity as: v^i = u^i/u^0. Note that this definition differs from the Valencia formalism.\"\n", "# It is useful to split Bi from Bi_stagger, since we're generally only interested in outputting Bi for diagnostics\n", "cctk_real grmhd_primitives_Bi type = GF TAGS='InterpNumTimelevels=1 prolongation=\"none\"'\n", "{\n", " Bx,By,Bz\n", "} \"B-field components defined at vertices.\"\n", "cctk_real grmhd_primitives_Bi_stagger type = GF TAGS='InterpNumTimelevels=1 prolongation=\"none\"'\n", "{\n", " Bx_stagger,By_stagger,Bz_stagger\n", "} \"B-field components defined at staggered points [Bx_stagger at (i+1/2,j,k),By_stagger at (i,j+1/2,k),Bz_stagger at (i,j,k+1/2)].\"\n", "#^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^#\n", "\n", "#vvvvvvv BSSN-based quantities, computed from ADM quantities.v vvvvvvv#\n", "cctk_real BSSN_quantities type = GF TAGS='prolongation=\"none\" Checkpoint=\"no\"'\n", "{\n", " gtxx,gtxy,gtxz,gtyy,gtyz,gtzz,gtupxx,gtupxy,gtupxz,gtupyy,gtupyz,gtupzz,phi_bssn,psi_bssn,lapm1\n", "} \"BSSN quantities, computed from ADM quantities\"\n", "#^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^#\n", "\n", "private:\n", "\n", "#vvvvvvv DIAGNOSTIC GRIDFUNCTIONS vvvvvvv#\n", "cctk_real diagnostic_gfs type = GF TAGS='prolongation=\"none\" Checkpoint=\"no\"'\n", "{\n", " failure_checker\n", "} \"Gridfunction to track conservative-to-primitives solver fixes. Beware that this gridfunction is overwritten at each RK substep.\"\n", "#^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^#\n", "\n", "#vvvvvvv TEMPORARY VARIABLES FOR RECONSTRUCTION vvvvvvv#\n", "cctk_real grmhd_primitives_reconstructed_temps type = GF TAGS='prolongation=\"none\" Checkpoint=\"no\"'\n", "{\n", " ftilde_gf,temporary,\n", " rho_br,Pr,vxr,vyr,vzr,Bxr,Byr,Bzr,Bx_staggerr,By_staggerr,Bz_staggerr,\n", " rho_bl,Pl,vxl,vyl,vzl,Bxl,Byl,Bzl,Bx_staggerl,By_staggerl,Bz_staggerl,\n", " vxrr,vxrl,vyrr,vyrl,vzrr,vzrl,vxlr,vxll,vylr,vyll,vzlr,vzll\n", "} \"Temporary variables used for primitives reconstruction\"\n", "#^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^#\n", "\n", "#vvvvvvv RHS VARIABLES vvvvvvv#\n", "cctk_real grmhd_conservatives_rhs type = GF TAGS='prolongation=\"none\" Checkpoint=\"no\"'\n", "{\n", " rho_star_rhs,tau_rhs,st_x_rhs,st_y_rhs,st_z_rhs\n", "} \"Storage for the right-hand side of the partial_t rho_star, partial_t tau, and partial_t tilde{S}_i equations. Needed for MoL timestepping.\"\n", "\n", "cctk_real em_Ax_rhs type = GF TAGS='prolongation=\"none\" Checkpoint=\"no\"'\n", "{\n", " Ax_rhs\n", "} \"Storage for the right-hand side of the partial_t A_x equation. Needed for MoL timestepping.\"\n", "cctk_real em_Ay_rhs type = GF TAGS='prolongation=\"none\" Checkpoint=\"no\"'\n", "{\n", " Ay_rhs\n", "} \"Storage for the right-hand side of the partial_t A_y equation. Needed for MoL timestepping.\"\n", "cctk_real em_Az_rhs type = GF TAGS='prolongation=\"none\" Checkpoint=\"no\"'\n", "{\n", " Az_rhs\n", "} \"Storage for the right-hand side of the partial_t A_z equation. Needed for MoL timestepping.\"\n", "cctk_real em_psi6phi_rhs type = GF TAGS='prolongation=\"none\" Checkpoint=\"no\"'\n", "{\n", " psi6phi_rhs\n", "} \"Storage for the right-hand side of the partial_t (psi^6 Phi) equation. Needed for MoL timestepping.\"\n", "#^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^#\n", "\n", "#vvvvvvv TEMPORARY VARIABLES USEFUL FOR A-FIELD EVOLUTION vvvvvvv#\n", "cctk_real grmhd_cmin_cmax_temps type = GF TAGS='prolongation=\"none\" Checkpoint=\"no\"'\n", "{\n", " cmin_x,cmax_x,\n", " cmin_y,cmax_y,\n", " cmin_z,cmax_z\n", "} \"Store min and max characteristic speeds in all three directions.\"\n", "#^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^#\n", "\n", "#vvvvvvv TEMPORARY VARIABLES USEFUL FOR FLUX COMPUTATION vvvvvvv#\n", "cctk_real grmhd_flux_temps type = GF TAGS='prolongation=\"none\" Checkpoint=\"no\"'\n", "{\n", " rho_star_flux,tau_flux,st_x_flux,st_y_flux,st_z_flux\n", "} \"Temporary variables for storing the flux terms of tilde{S}_i.\"\n", "#^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^#\n", "\n", "#vvvvvvv T^{\\mu \\nu}, stored to avoid expensive recomputation vvvvvvv#\n", "cctk_real TUPmunu type = GF TAGS='prolongation=\"none\" Checkpoint=\"no\"'\n", "{\n", " TUPtt,TUPtx,TUPty,TUPtz,TUPxx,TUPxy,TUPxz,TUPyy,TUPyz,TUPzz\n", "} \"T^{mu nu}, stored to avoid expensive recomputation\"\n", "#^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^#\n", "\n", "###########################################################################\n", "####################################################\n", "### Functions provided by MoL for registration ###\n", "####################################################\n", "\n", "CCTK_INT FUNCTION MoLRegisterEvolved(CCTK_INT IN EvolvedIndex, \\\n", " CCTK_INT IN RHSIndex)\n", "CCTK_INT FUNCTION MoLRegisterEvolvedSlow(CCTK_INT IN EvolvedIndex, \\\n", " CCTK_INT IN RHSIndexSlow)\n", "CCTK_INT FUNCTION MoLRegisterConstrained(CCTK_INT IN ConstrainedIndex)\n", "CCTK_INT FUNCTION MoLRegisterEvolvedGroup(CCTK_INT IN EvolvedIndex, \\\n", " CCTK_INT IN RHSIndex)\n", "CCTK_INT FUNCTION MoLRegisterEvolvedGroupSlow(CCTK_INT IN EvolvedIndex, \\\n", " CCTK_INT IN RHSIndexSlow)\n", "CCTK_INT FUNCTION MoLRegisterConstrainedGroup(CCTK_INT IN ConstrainedIndex)\n", "CCTK_INT FUNCTION MoLRegisterSaveAndRestoreGroup(CCTK_INT IN SandRIndex)\n", "\n", "USES FUNCTION MoLRegisterEvolved\n", "USES FUNCTION MoLRegisterEvolvedSlow\n", "USES FUNCTION MoLRegisterConstrained\n", "USES FUNCTION MoLRegisterEvolvedGroup\n", "USES FUNCTION MoLRegisterEvolvedGroupSlow\n", "USES FUNCTION MoLRegisterConstrainedGroup\n", "USES FUNCTION MoLRegisterSaveAndRestoreGroup\n", "\n", "#########################################\n", "### Aliased functions from Boundary ###\n", "#########################################\n", "\n", "CCTK_INT FUNCTION Boundary_SelectVarForBC(CCTK_POINTER_TO_CONST IN GH, \\\n", " CCTK_INT IN faces, CCTK_INT IN boundary_width, CCTK_INT IN table_handle, \\\n", " CCTK_STRING IN var_name, CCTK_STRING IN bc_name)\n", "CCTK_INT FUNCTION Boundary_SelectGroupForBC(CCTK_POINTER_TO_CONST IN GH, \\\n", " CCTK_INT IN faces, CCTK_INT IN boundary_width, CCTK_INT IN table_handle, \\\n", " CCTK_STRING IN group_name, CCTK_STRING IN bc_name)\n", "\n", "USES FUNCTION Boundary_SelectVarForBC\n", "USES FUNCTION Boundary_SelectGroupForBC\n", "###########################################################################\n", "\n", "#########################################\n", "### Aliased functions from Carpet ###\n", "#########################################\n", "\n", "CCTK_INT FUNCTION \\\n", " GetRefinementLevel \\\n", " (CCTK_POINTER_TO_CONST IN cctkGH)\n", "USES FUNCTION GetRefinementLevel\n", "\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Step 5: `param.ccl` \\[Back to [top](#toc)\\]\n", "$$\\label{param__ccl}$$" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Writing ../param.ccl\n" ] } ], "source": [ "%%writefile $outfile_path__param__ccl\n", "# Parameter definitions for thorn IllinoisGRMHD\n", "\n", "shares: ADMBase\n", "USES CCTK_INT lapse_timelevels\n", "USES CCTK_INT shift_timelevels\n", "USES CCTK_INT metric_timelevels\n", "\n", "#########################################################\n", "restricted:\n", "#########################################################\n", "# Verbosity Level\n", "KEYWORD verbose \"Determines how much evolution information is output\" STEERABLE=ALWAYS\n", "{\n", " \"no\" :: \"Complete silence\"\n", " \"essential\" :: \"Essential health monitoring of the GRMHD evolution: Information about conservative-to-primitive fixes, etc.\"\n", " \"essential+iteration output\" :: \"Outputs health monitoring information, plus a record of which RK iteration. Very useful for backtracing a crash.\"\n", "} \"essential+iteration output\"\n", "#########################################################\n", "\n", "\n", "#########################################################\n", "# SPEED LIMIT: Set maximum relativistic gamma factor\n", "#\n", "REAL GAMMA_SPEED_LIMIT \"Maximum relativistic gamma factor.\"\n", "{\n", " 1:* :: \"Positive > 1, though you'll likely have troubles far above 10.\"\n", "} 10.0\n", "#########################################################\n", "\n", "#########################################################\n", "# CONSERV TO PRIMS PARAMETERS\n", "\n", "# FIXME: Enable this parameter! IllinoisGRMHD is currently hard-coded to tau_stildefix_enable=2.\n", "#INT tau_stildefix_enable \"tau<0 fix in primitive_vars_hybrid2 to reduce number of Font fixes, especially in puncture+matter evolutions\" STEERABLE=ALWAYS\n", "#{\n", "# 0:3 :: \"zero (disable), one (enable everywhere), or two (enable only where Psi6 > Psi6threshold [i.e., inside the horizon, where B's are set to zero], or three (kludge: set B=0 if tau<0 inside horizon))\"\n", "#} 0\n", "\n", "BOOLEAN update_Tmunu \"Update Tmunu, for RHS of Einstein's equations?\" STEERABLE=ALWAYS\n", "{\n", "} \"yes\"\n", "\n", "############################\n", "# Limiters on tau and rho_b:\n", "REAL tau_atm \"Floor value on the energy variable tau (cf. tau_stildefix_enable). Given the variety of systems this code may encounter, there *is no reasonable default*. Effectively the current (enormous) value should disable the tau_atm floor. Please set this in your initial data thorn, and reset at will during evolutions.\" STEERABLE=ALWAYS\n", "{\n", " 0:* :: \"Positive\"\n", "} 1e100\n", "\n", "REAL rho_b_atm \"Floor value on the baryonic rest mass density rho_b (atmosphere). Given the variety of systems this code may encounter, there *is no reasonable default*. Your run will die unless you override this default value in your initial data thorn.\" STEERABLE=ALWAYS\n", "{\n", " *:* :: \"Allow for negative values. This enables us to debug the code and verify if rho_b_atm is properly set.\"\n", "} 1e200\n", "\n", "REAL rho_b_max \"Ceiling value on the baryonic rest mass density rho_b. The enormous value effectively disables this ceiling by default. It can be quite useful after a black hole has accreted a lot of mass, leading to enormous densities inside the BH. To enable this trick, set rho_b_max in your initial data thorn! You are welcome to change this parameter mid-run (after restarting from a checkpoint).\" STEERABLE=ALWAYS\n", "{\n", " 0:* :: \"Note that you will have problems unless rho_b_atm Psi6threshold, we assume we're inside the horizon in the primitives solver, and certain limits are relaxed or imposed\" STEERABLE=ALWAYS\n", "{\n", " *:* :: \"Can set to anything\"\n", "} 1e100\n", "#########################################################\n", "\n", "#########################################################\n", "# EQUATION OF STATE PARAMS, LOOK FOR MORE IN interface.ccl!\n", "INT neos \"number of parameters in EOS table. If you want to increase from the default max value, you MUST also set eos_params_arrays1 and eos_params_arrays2 in interface.ccl to be consistent!\"\n", "{\n", " 1:10 :: \"Any integer between 1 and 10\"\n", "} 1\n", "\n", "REAL Gamma_th \"thermal gamma parameter\"\n", "{\n", " 0:* :: \"Physical values\"\n", "-1 :: \"forbidden value to make sure it is explicitly set in the parfile\"\n", "} -1\n", "\n", "REAL K_ppoly_tab0 \"Also known as k_tab[0], this is the polytropic constant for the lowest density piece of the piecewise polytrope. All other k_tab EOS array elements are set from user-defined rho_tab EOS array elements and by enforcing continuity in the equation of state.\"\n", "{\n", " 0:* :: \"Physical values\"\n", "-1 :: \"forbidden value to make sure it is explicitly set in the parfile\"\n", "} -1\n", "\n", "REAL rho_ppoly_tab_in[10] \"Set polytropic rho parameters\"\n", "{\n", " 0.0:* :: \"after this time (inclusively)\"\n", "-1.0 :: \"forbidden value to make sure it is explicitly set in the parfile\"\n", "} -1.0\n", "\n", "REAL Gamma_ppoly_tab_in[11] \"Set polytropic rho parameters\"\n", "{\n", " 0.0:* :: \"after this time (inclusively)\"\n", "-1.0 :: \"forbidden value to make sure it is explicitly set in the parfile\"\n", "} -1.0\n", "\n", "#########################################################\n", "\n", "#########################################################\n", "# OUTER BOUNDARY CONDITION CHOICE\n", "KEYWORD Matter_BC \"Chosen Matter boundary condition\"\n", "{\n", " \"outflow\" :: \"Outflow boundary conditions\"\n", " \"frozen\" :: \"Frozen boundaries\"\n", "} \"outflow\"\n", "\n", "KEYWORD EM_BC \"EM field boundary condition\"\n", "{\n", " \"copy\" :: \"Copy data from nearest boundary point\"\n", " \"frozen\" :: \"Frozen boundaries\"\n", "} \"copy\"\n", "#########################################################\n", "\n", "\n", "#########################################################\n", "# SYMMETRY BOUNDARY PARAMS. Needed for handling staggered gridfunctions.\n", "KEYWORD Symmetry \"Currently only no symmetry supported, though work has begun in adding equatorial-symmetry support. FIXME: Extend ET symmetry interface to support symmetries on staggered gridfunctions\"\n", "{\n", " \"none\" :: \"no symmetry, full 3d domain\"\n", "} \"none\"\n", "\n", "REAL Sym_Bz \"In-progress equatorial symmetry support: Symmetry parameter across z axis for magnetic fields = +/- 1\"\n", "{\n", " -1.0:1.0 :: \"Set to +1 or -1.\"\n", "} 1.0\n", "#########################################################\n", "\n", "###############################################################################################\n", "private:\n", "\n", "#########################################################\n", "# EVOLUTION PARAMS\n", "REAL damp_lorenz \"Damping factor for the generalized Lorenz gauge. Has units of 1/length = 1/M. Typically set this parameter to 1.5/(maximum Delta t on AMR grids).\" STEERABLE=ALWAYS\n", "{\n", " *:* :: \"any real\"\n", "} 0.0\n", "#########################################################\n", "\n", "\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Step 6: `schedule.ccl` \\[Back to [top](#toc)\\]\n", "$$\\label{schedule__ccl}$$" ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Writing ../schedule.ccl\n" ] } ], "source": [ "%%writefile $outfile_path__schedule__ccl\n", "# Scheduler setup for IllinoisGRMHD\n", "\n", "STORAGE: ADMBase::metric[metric_timelevels], ADMBase::curv[metric_timelevels], ADMBase::lapse[lapse_timelevels], ADMBase::shift[shift_timelevels]\n", "STORAGE: IllinoisGRMHD::BSSN_quantities\n", "\n", "STORAGE: grmhd_conservatives[3],em_Ax[3],em_Ay[3],em_Az[3],em_psi6phi[3]\n", "STORAGE: grmhd_primitives_allbutBi,grmhd_primitives_Bi,grmhd_primitives_Bi_stagger,grmhd_primitives_reconstructed_temps,grmhd_conservatives_rhs,em_Ax_rhs,em_Ay_rhs,em_Az_rhs,em_psi6phi_rhs,grmhd_cmin_cmax_temps,grmhd_flux_temps,TUPmunu,diagnostic_gfs\n", "\n", "####################\n", "# RUN INITIALLY ONLY\n", "schedule IllinoisGRMHD_RegisterVars in MoL_Register after BSSN_RegisterVars after lapse_RegisterVars\n", "{\n", " LANG: C\n", " OPTIONS: META\n", "} \"Register evolved, rhs variables in IllinoisGRMHD for MoL\"\n", "\n", "# Tells the symmetry thorn how to apply symmetries on each gridfunction:\n", "schedule IllinoisGRMHD_InitSymBound at BASEGRID after Lapse_InitSymBound\n", "{\n", " LANG: C\n", "} \"Schedule symmetries\"\n", "####################\n", "\n", "####################\n", "# POSTPOSTINITIAL\n", "# HydroBase_Con2Prim in CCTK_POSTPOSTINITIAL sets conserv to prim then\n", "# outer boundaries (OBs, which are technically disabled). The post OB\n", "# SYNCs actually reprolongate the conservative variables, making cons\n", "# and prims INCONSISTENT. So here we redo the con2prim, avoiding the\n", "# SYNC afterward, then copy the result to other timelevels\"\n", "schedule GROUP IllinoisGRMHD_PostPostInitial at CCTK_POSTPOSTINITIAL before MoL_PostStep after HydroBase_Con2Prim\n", "{\n", "} \"HydroBase_Con2Prim in CCTK_POSTPOSTINITIAL sets conserv to prim then outer boundaries (OBs, which are technically disabled). The post OB SYNCs actually reprolongate the conservative variables, making cons and prims INCONSISTENT. So here we redo the con2prim, avoiding the SYNC afterward, then copy the result to other timelevels\"\n", "\n", "schedule IllinoisGRMHD_InitSymBound in IllinoisGRMHD_PostPostInitial as postid before compute_b\n", "{\n", " SYNC: grmhd_conservatives,em_Ax,em_Ay,em_Az,em_psi6phi\n", " LANG: C\n", "} \"Schedule symmetries -- Actually just a placeholder function to ensure prolongations / processor syncs are done BEFORE outer boundaries are updated.\"\n", "\n", "# Easiest primitives to solve for: B^i\n", "schedule IllinoisGRMHD_compute_B_and_Bstagger_from_A in IllinoisGRMHD_PostPostInitial as compute_b after postid after empostid after lapsepostid\n", "{\n", " # This is strictly a processor sync, as prolongation is disabled for all primitives & B^i's.\n", " SYNC: grmhd_primitives_Bi,grmhd_primitives_Bi_stagger # FIXME: Are both SYNC's necessary?\n", " LANG: C\n", "} \"Compute B and B_stagger from A SYNC: grmhd_primitives_Bi,grmhd_primitives_Bi_stagger\"\n", "\n", "# Nontrivial primitives solve, for P,rho_b,vx,vy,vz:\n", "schedule IllinoisGRMHD_conserv_to_prims in IllinoisGRMHD_PostPostInitial after compute_b\n", "{\n", " LANG: C\n", "} \"Compute primitive variables from conservatives. This is non-trivial, requiring a Newton-Raphson root-finder.\"\n", "\n", "# Copy data to other timelevels.\n", "schedule IllinoisGRMHD_PostPostInitial_Set_Symmetries__Copy_Timelevels in IllinoisGRMHD_PostPostInitial as mhdpostid after compute_b after p2c\n", "{\n", " LANG: C\n", "} \"Compute post-initialdata quantities\"\n", "####################\n", "\n", "####################\n", "# RHS EVALUATION\n", "schedule IllinoisGRMHD_driver_evaluate_MHD_rhs in MoL_CalcRHS as IllinoisGRMHD_RHS_eval after bssn_rhs after shift_rhs\n", "{\n", " LANG: C\n", "} \"Evaluate RHSs of GR Hydro & GRMHD equations\"\n", "####################\n", "\n", "\n", "############################################################\n", "# COMPUTE B FROM A & RE-SOLVE FOR PRIMITIVES\n", "# After a full timestep, there are two types of boundaries that need filling:\n", "# (A) Outer boundaries (on coarsest level)\n", "# (B) AMR grid refinement boundaries\n", "\n", "# (A) OUTER BOUNDARY STEPS:\n", "# ( 0) Synchronize (prolongate/restrict) all evolved variables\n", "# ( 1) Apply outer boundary conditions (BCs) on A_{\\mu}\n", "# ( 2) Compute B^i from A_i everywhere, synchronize (processor sync) B^i\n", "# ( 3) Call con2prim to get consistent primitives {P,rho_b,vx,vy,vz} and conservatives at all points (if no restriction, really only need interior)\n", "# ( 4) Apply outer BCs on {P,rho_b,vx,vy,vz}, recompute conservatives.\n", "\n", "# (B) AMR GRID REFINEMENT BOUNDARY STEPS:\n", "# Same as steps 0,2,3 above. Just need if() statements in steps 1,4 to prevent \"outer boundaries\" being updated\n", "# Problem: all the sync's in outer boundary updates might just overwrite prolongated values.\n", "############################################################\n", "\n", "schedule IllinoisGRMHD_InitSymBound in HydroBase_Boundaries as Boundary_SYNCs before compute_B_postrestrict\n", "{\n", " SYNC: grmhd_conservatives,em_Ax,em_Ay,em_Az,em_psi6phi\n", " LANG: C\n", "} \"Schedule symmetries -- Actually just a placeholder function to ensure prolongations / processor syncs are done BEFORE outer boundaries are updated.\"\n", "\n", "schedule IllinoisGRMHD_outer_boundaries_on_A_mu in HydroBase_Boundaries after Boundary_SYNCs before mhd_conserv2prims_postrestrict\n", "{\n", " LANG: C\n", "} \"Apply linear extrapolation BCs on A_{mu}, so that BCs are flat on B^i.\"\n", "\n", "# Easiest primitives to solve for: B^i.\n", "# Note however that B^i depends on derivatives of A_{\\mu}, so a SYNC is necessary on B^i.\n", "schedule IllinoisGRMHD_compute_B_and_Bstagger_from_A in HydroBase_Boundaries after outer_boundaries_on_A_mu\n", "{\n", " # This is strictly a processor sync, as prolongation is disabled for all B^i's.\n", " SYNC: grmhd_primitives_Bi,grmhd_primitives_Bi_stagger # FIXME: Are both SYNC's necessary?\n", " LANG: C\n", "} \"Compute B and B_stagger from A, SYNC: grmhd_primitives_Bi,grmhd_primitives_Bi_stagger\"\n", "\n", "# Nontrivial primitives solve, for P,rho_b,vx,vy,vz.\n", "schedule IllinoisGRMHD_conserv_to_prims in AddToTmunu after compute_B_and_Bstagger_from_A\n", "{\n", " LANG: C\n", "} \"Compute primitive variables from conservatives. This is non-trivial, requiring a Newton-Raphson root-finder.\"\n", "\n", "schedule IllinoisGRMHD_outer_boundaries_on_P_rho_b_vx_vy_vz in AddToTmunu after IllinoisGRMHD_conserv_to_prims\n", "{\n", "# We must sync {P,rho_b,vx,vy,vz} here.\n", " SYNC: grmhd_primitives_allbutBi\n", " LANG: C\n", "} \"Apply outflow-only, flat BCs on {P,rho_b,vx,vy,vz}. Outflow only == velocities pointed inward from outer boundary are set to zero.\"\n", "##########################################################\n", "\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Step 7: `code_validation.txt` \\[Back to [top](#toc)\\]\n", "$$\\label{code_validation__txt}$$" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Writing code_validation.txt\n" ] } ], "source": [ "%%writefile $outfile_path__code_validation__txt\n", "0 0 0.129285232345409 0\n", "6 0.375 0.127243949890016 0\n", "12 0.75 0.12663194218958 9.99201e-16\n", "18 1.125 0.126538236778999 0\n", "24 1.5 0.126693696091085 0\n", "30 1.875 0.126354095699745 0\n", "36 2.25 0.125536381948334 0\n", "42 2.625 0.124535850791511 0\n", "48 3 0.123592701331336 -1.9984e-15\n", "54 3.375 0.122798174115381 6.00908e-15\n", "60 3.75 0.122119289857302 5.20001e-14\n", "66 4.125 0.121490429374158 -7.219e-12\n", "72 4.5 0.120874029334284 9.062e-12\n", "78 4.875 0.120259717466328 8.073e-12\n", "84 5.25 0.11967641900387 6.773e-12\n", "90 5.625 0.119123915533179 3.9644e-11\n", "96 6 0.118609688443716 8.7028e-11\n", "102 6.375 0.118142085593448 1.0356e-10\n", "108 6.75 0.117724554127494 -2.1806e-11\n", "114 7.125 0.117368613539 -2.98919e-10\n", "120 7.5 0.117103537221293 -6.55126e-10\n", "126 7.875 0.116916242747612 -8.58055e-10\n", "132 8.25 0.116788537512057 -6.40583e-10\n", "138 8.625 0.116760604471516 1.2299e-10\n", "144 9 0.116800459472318 1.17986e-09\n", "150 9.375 0.116876502360487 2.30931e-09\n", "156 9.75 0.117001584159343 3.22793e-09\n", "162 10.125 0.117152245129722 3.70731e-09\n", "168 10.5 0.117321427488885 3.80249e-09\n", "174 10.875 0.117504914930401 3.18643e-09\n", "180 11.25 0.117677337706038 1.53467e-09\n", "186 11.625 0.11783535130534 -1.0388e-09\n", "192 12 0.11797817411705 -4.09727e-09\n", "198 12.375 0.118099716148388 -7.36175e-09\n", "204 12.75 0.118197132134204 -1.11993e-08\n", "210 13.125 0.118268160681997 -1.54798e-08\n", "216 13.5 0.118308938021455 -1.879e-08\n", "222 13.875 0.118311817073909 -1.9184e-08\n", "228 14.25 0.118266030948925 -1.56794e-08\n", "234 14.625 0.118178894932027 -8.96657e-09\n", "240 15 0.118048178792271 -8.64983e-10\n", "\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Step 8: Code validation \\[Back to [top](#toc)\\]\n", "$$\\label{code_validation}$$\n", "\n", "\n", "\n", "## Step 8.a: `make.code.defn` \\[Back to [top](#toc)\\]\n", "$$\\label{code_validation__make_code_defn}$$\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": 8, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Validation test for make.code.defn: FAILED!\n", "Diff:\n", "15a16\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/make.code.defn\"\n", "original_IGM_file_name = \"make.code.defn-original\"\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__make_code_defn = !diff $original_IGM_file_path $outfile_path__make_code_defn\n", "\n", "if Validation__make_code_defn == []:\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 make.code.defn: PASSED!\")\n", "else:\n", " # If the validation fails, we keep the original IGM source code file\n", " print(\"Validation test for make.code.defn: 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__make_code_defn:\n", " print(diff_line)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "## Step 8.b: `configuration.ccl` \\[Back to [top](#toc)\\]\n", "$$\\label{code_validation__configuration__ccl}$$\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": 9, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Validation test for configuration.ccl: FAILED!\n", "Diff:\n", "6a7\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/configuration.ccl\"\n", "original_IGM_file_name = \"configuration-original.ccl\"\n", "original_IGM_file_path = os.path.join(IGM_main_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__configuration__ccl = !diff $original_IGM_file_path $outfile_path__configuration__ccl\n", "\n", "if Validation__configuration__ccl == []:\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 configuration.ccl: PASSED!\")\n", "else:\n", " # If the validation fails, we keep the original IGM source code file\n", " print(\"Validation test for configuration.ccl: 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__configuration__ccl:\n", " print(diff_line)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "## Step 8.c: `interface.ccl` \\[Back to [top](#toc)\\]\n", "$$\\label{code_validation__interface__ccl}$$\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": 10, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Validation test for interface.ccl: FAILED!\n", "Diff:\n", "46c46\n", "< psi6phi \n", "---\n", "> psi6phi\n", "68,84d67\n", "< #vvvvvvv EQUATION OF STATE VARIABLES vvvvvvv#\n", "< REAL eos_params_arrays1 TYPE=ARRAY DISTRIB=CONSTANT DIM=1 SIZE=10\n", "< {\n", "< rho_tab,P_tab,eps_tab\n", "< } \"Equation of state (EOS) parameters\"\n", "< \n", "< REAL eos_params_arrays2 TYPE=ARRAY DISTRIB=CONSTANT DIM=1 SIZE=11\n", "< {\n", "< k_tab,gamma_tab\n", "< } \"Equation of state (EOS) parameters\"\n", "< \n", "< REAL eos_params_scalar TYPE = SCALAR \n", "< {\n", "< n_poly\n", "< } \"polytropic index\"\n", "< #^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^#\n", "< \n", "118,119c101,102\n", "< { \n", "< Ax_rhs \n", "---\n", "> {\n", "> Ax_rhs\n", "123c106\n", "< Ay_rhs \n", "---\n", "> Ay_rhs\n", "126,127c109,110\n", "< { \n", "< Az_rhs \n", "---\n", "> {\n", "> Az_rhs\n", "131c114\n", "< psi6phi_rhs \n", "---\n", "> psi6phi_rhs\n", "205a189\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/interface.ccl\"\n", "original_IGM_file_name = \"interface-original.ccl\"\n", "original_IGM_file_path = os.path.join(IGM_main_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__interface__ccl = !diff $original_IGM_file_path $outfile_path__interface__ccl\n", "\n", "if Validation__interface__ccl == []:\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 interface.ccl: PASSED!\")\n", "else:\n", " # If the validation fails, we keep the original IGM source code file\n", " print(\"Validation test for interface.ccl: 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__interface__ccl:\n", " print(diff_line)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "## Step 8.d: `param.ccl` \\[Back to [top](#toc)\\]\n", "$$\\label{code_validation__param__ccl}$$\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": 11, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Validation test for param.ccl: FAILED!\n", "Diff:\n", "23c23\n", "< # \n", "---\n", "> #\n", "79c79\n", "< REAL gamma_th \"thermal gamma parameter\"\n", "---\n", "> REAL Gamma_th \"thermal gamma parameter\"\n", "85c85\n", "< REAL K_poly \"initial polytropic constant\"\n", "---\n", "> REAL K_ppoly_tab0 \"Also known as k_tab[0], this is the polytropic constant for the lowest density piece of the piecewise polytrope. All other k_tab EOS array elements are set from user-defined rho_tab EOS array elements and by enforcing continuity in the equation of state.\"\n", "87,88c87,102\n", "< 0:* :: \"Positive\"\n", "< } 1.0\n", "---\n", "> 0:* :: \"Physical values\"\n", "> -1 :: \"forbidden value to make sure it is explicitly set in the parfile\"\n", "> } -1\n", "> \n", "> REAL rho_ppoly_tab_in[10] \"Set polytropic rho parameters\"\n", "> {\n", "> 0.0:* :: \"after this time (inclusively)\"\n", "> -1.0 :: \"forbidden value to make sure it is explicitly set in the parfile\"\n", "> } -1.0\n", "> \n", "> REAL Gamma_ppoly_tab_in[11] \"Set polytropic rho parameters\"\n", "> {\n", "> 0.0:* :: \"after this time (inclusively)\"\n", "> -1.0 :: \"forbidden value to make sure it is explicitly set in the parfile\"\n", "> } -1.0\n", "> \n", "129a144,145\n", "> \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/param.ccl\"\n", "original_IGM_file_name = \"param-original.ccl\"\n", "original_IGM_file_path = os.path.join(IGM_main_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__param__ccl = !diff $original_IGM_file_path $outfile_path__param__ccl\n", "\n", "if Validation__param__ccl == []:\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 param.ccl: PASSED!\")\n", "else:\n", " # If the validation fails, we keep the original IGM source code file\n", " print(\"Validation test for param.ccl: 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__param__ccl:\n", " print(diff_line)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "## Step 8.e: `schedule.ccl` \\[Back to [top](#toc)\\]\n", "$$\\label{code_validation__schedule__ccl}$$\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": 12, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Validation test for schedule.ccl: FAILED!\n", "Diff:\n", "6c6\n", "< STORAGE: grmhd_conservatives[3],em_Ax[3],em_Ay[3],em_Az[3],em_psi6phi[3] \n", "---\n", "> STORAGE: grmhd_conservatives[3],em_Ax[3],em_Ay[3],em_Az[3],em_psi6phi[3]\n", "8d7\n", "< STORAGE: eos_params_arrays1,eos_params_arrays2,eos_params_scalar\n", "28c27\n", "< # outer boundaries (OBs, which are technically disabled). The post OB \n", "---\n", "> # outer boundaries (OBs, which are technically disabled). The post OB\n", "30c29\n", "< # and prims INCONSISTENT. So here we redo the con2prim, avoiding the \n", "---\n", "> # and prims INCONSISTENT. So here we redo the con2prim, avoiding the\n", "122a122\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/schedule.ccl\"\n", "original_IGM_file_name = \"schedule-original.ccl\"\n", "original_IGM_file_path = os.path.join(IGM_main_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__schedule__ccl = !diff $original_IGM_file_path $outfile_path__schedule__ccl\n", "\n", "if Validation__schedule__ccl == []:\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 schedule.ccl: PASSED!\")\n", "else:\n", " # If the validation fails, we keep the original IGM source code file\n", " print(\"Validation test for schedule.ccl: 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__schedule__ccl:\n", " print(diff_line)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "## Step 8.f: `code_validation.txt` \\[Back to [top](#toc)\\]\n", "$$\\label{code_validation__code_validation__txt}$$\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": 13, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Validation test for code_validation.txt: FAILED!\n", "Diff:\n", "41a42\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/code_validation.txt\"\n", "original_IGM_file_name = \"code_validation-original.txt\"\n", "original_IGM_file_path = os.path.join(IGM_main_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__code_validation__txt = !diff $original_IGM_file_path $outfile_path__code_validation__txt\n", "\n", "if Validation__code_validation__txt == []:\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 code_validation.txt: PASSED!\")\n", "else:\n", " # If the validation fails, we keep the original IGM source code file\n", " print(\"Validation test for code_validation.txt: 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__code_validation__txt:\n", " print(diff_line)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Step 9: 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__Configuration_files.pdf](Tutorial-IllinoisGRMHD__Configuration_files.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": 14, "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__Configuration_files.ipynb\n", "#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__Configuration_files.tex\n", "#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__Configuration_files.tex\n", "#!pdflatex -interaction=batchmode Tutorial-IllinoisGRMHD__Configuration_files.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.10.11" } }, "nbformat": 4, "nbformat_minor": 4 }