{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "\n", "# Equations of General Relativistic Magnetohydrodynamics (GRMHD)\n", "\n", "## Author: Zach Etienne\n", "\n", "## This notebook documents and constructs a number of quantities useful for building symbolic (SymPy) expressions for the equations of general relativistic magnetohydrodynamics (GRMHD), using the same (Valencia-like) formalism as `IllinoisGRMHD`.\n", "\n", "**Notebook Status:** Self-Validated; induction equation not yet implemented \n", "\n", "**Validation Notes:** This tutorial notebook has been confirmed to be self-consistent with its corresponding NRPy+ module, as documented [below](#code_validation). **Additional validation tests may have been performed, but are as yet, undocumented. (TODO)**\n", "\n", "## Introduction\n", "\n", "We write the equations of general relativistic magnetohydrodynamics in conservative form as follows (Eqs. 41-44 of [Duez *et al*](https://arxiv.org/pdf/astro-ph/0503420.pdf):\n", "\n", "\\begin{eqnarray}\n", "\\ \\partial_t \\rho_* &+& \\partial_j \\left(\\rho_* v^j\\right) = 0 \\\\\n", "\\partial_t \\tilde{\\tau} &+& \\partial_j \\left(\\alpha^2 \\sqrt{\\gamma} T^{0j} - \\rho_* v^j \\right) = s \\\\\n", "\\partial_t \\tilde{S}_i &+& \\partial_j \\left(\\alpha \\sqrt{\\gamma} T^j{}_i \\right) = \\frac{1}{2} \\alpha\\sqrt{\\gamma} T^{\\mu\\nu} g_{\\mu\\nu,i} \\\\\n", "\\partial_t \\tilde{B}^i &+& \\partial_j \\left(v^j \\tilde{B}^i - v^i \\tilde{B}^j\\right) = 0,\n", "\\end{eqnarray}\n", "where \n", "\n", "$$\n", "s = \\alpha \\sqrt{\\gamma}\\left[\\left(T^{00}\\beta^i\\beta^j + 2 T^{0i}\\beta^j + T^{ij} \\right)K_{ij}\n", "- \\left(T^{00}\\beta^i + T^{0i} \\right)\\partial_i\\alpha \\right].\n", "$$\n", "\n", "We represent $T^{\\mu\\nu}$ as the sum of the stress-energy tensor of a perfect fluid $T^{\\mu\\nu}_{\\rm GRHD}$, plus the stress-energy associated with the electromagnetic fields in the force-free electrodynamics approximation $T^{\\mu\\nu}_{\\rm GRFFE}$ (equivalently, $T^{\\mu\\nu}_{\\rm em}$ in Duez *et al*):\n", "\n", "$$\n", "T^{\\mu\\nu} = T^{\\mu\\nu}_{\\rm GRHD} + T^{\\mu\\nu}_{\\rm GRFFE},\n", "$$\n", "\n", "where \n", "\n", "* $T^{\\mu\\nu}_{\\rm GRHD}$ is constructed from rest-mass density $\\rho_0$, pressure $P$, internal energy $\\epsilon$, 4-velocity $u^\\mu$, and ADM metric quantities as described in the [NRPy+ GRHD equations tutorial notebook](Tutorial-GRHD_Equations-Cartesian.ipynb); and \n", "* $T^{\\mu\\nu}_{\\rm GRFFE}$ is constructed from the magnetic field vector $B^i$ and ADM metric quantities as described in the [NRPy+ GRFFE equations tutorial notebook](Tutorial-GRFFE_Equations-Cartesian.ipynb).\n", "\n", "All quantities can be written in terms of the full GRMHD stress-energy tensor $T^{\\mu\\nu}$ in precisely the same way they are defined in the GRHD equations. ***Therefore, we will not define special functions for generating these quantities, and instead refer the user to the appropriate functions in the [GRHD module](../edit/GRHD/equations.py)*** Namely,\n", "\n", "* The GRMHD conservative variables:\n", " * $\\rho_* = \\alpha\\sqrt{\\gamma} \\rho_0 u^0$, via `GRHD.compute_rho_star(alpha, sqrtgammaDET, rho_b,u4U)`\n", " * $\\tilde{\\tau} = \\alpha^2\\sqrt{\\gamma} T^{00} - \\rho_*$, via `GRHD.compute_tau_tilde(alpha, sqrtgammaDET, T4UU,rho_star)`\n", " * $\\tilde{S}_i = \\alpha \\sqrt{\\gamma} T^0{}_i$, via `GRHD.compute_S_tildeD(alpha, sqrtgammaDET, T4UD)`\n", "* The GRMHD fluxes:\n", " * $\\rho_*$ flux: $\\left(\\rho_* v^j\\right)$, via `GRHD.compute_rho_star_fluxU(vU, rho_star)`\n", " * $\\tilde{\\tau}$ flux: $\\left(\\alpha^2 \\sqrt{\\gamma} T^{0j} - \\rho_* v^j \\right)$, via `GRHD.compute_tau_tilde_fluxU(alpha, sqrtgammaDET, vU,T4UU,rho_star)`\n", " * $\\tilde{S}_i$ flux: $\\left(\\alpha \\sqrt{\\gamma} T^j{}_i \\right)$, via `GRHD.compute_S_tilde_fluxUD(alpha, sqrtgammaDET, T4UD)`\n", "* GRMHD source terms:\n", " * $\\tilde{\\tau}$ source term $s$: defined above, via `GRHD.compute_s_source_term(KDD,betaU,alpha, sqrtgammaDET,alpha_dD, T4UU)` \n", " * $\\tilde{S}_i$ source term: $\\frac{1}{2} \\alpha\\sqrt{\\gamma} T^{\\mu\\nu} g_{\\mu\\nu,i}$, via `GRHD.compute_S_tilde_source_termD(alpha, sqrtgammaDET,g4DD_zerotimederiv_dD, T4UU)`\n", "\n", "In summary, all terms in the GRMHD equations can be constructed once the full GRMHD stress-energy tensor $T^{\\mu\\nu} = T^{\\mu\\nu}_{\\rm GRHD} + T^{\\mu\\nu}_{\\rm GRFFE}$ is constructed. For completeness, the full set of input variables include:\n", "* Spacetime quantities:\n", " * ADM quantities $\\alpha$, $\\beta^i$, $\\gamma_{ij}$, $K_{ij}$\n", "* Hydrodynamical quantities:\n", " * Rest-mass density $\\rho_0$\n", " * Pressure $P$\n", " * Internal energy $\\epsilon$\n", " * 4-velocity $u^\\mu$\n", "* Electrodynamical quantities\n", " * Magnetic field $B^i= \\tilde{B}^i / \\gamma$\n", "\n", "### A Note on Notation\n", "\n", "As is standard in NRPy+, \n", "\n", "* Greek indices refer to four-dimensional quantities where the zeroth component indicates temporal (time) component.\n", "* Latin indices refer to three-dimensional quantities. This is somewhat counterintuitive since Python always indexes its lists starting from 0. As a result, the zeroth component of three-dimensional quantities will necessarily indicate the first *spatial* direction.\n", "\n", "For instance, in calculating the first term of $b^2 u^\\mu u^\\nu$, we use Greek indices:\n", "\n", "```python\n", "T4EMUU = ixp.zerorank2(DIM=4)\n", "for mu in range(4):\n", " for nu in range(4):\n", " # Term 1: b^2 u^{\\mu} u^{\\nu}\n", " T4EMUU[mu][nu] = smallb2*u4U[mu]*u4U[nu]\n", "```\n", "\n", "When we calculate $\\beta_i = \\gamma_{ij} \\beta^j$, we use Latin indices:\n", "```python\n", "betaD = ixp.zerorank1(DIM=3)\n", "for i in range(3):\n", " for j in range(3):\n", " betaD[i] += gammaDD[i][j] * betaU[j]\n", "```\n", "\n", "As a corollary, any expressions involving mixed Greek and Latin indices will need to offset one set of indices by one: A Latin index in a four-vector will be incremented and a Greek index in a three-vector will be decremented (however, the latter case does not occur in this tutorial notebook). This can be seen when we handle $\\frac{1}{2} \\alpha \\sqrt{\\gamma} T^{\\mu \\nu}_{\\rm EM} \\partial_i g_{\\mu \\nu}$:\n", "```python\n", "# \\alpha \\sqrt{\\gamma} T^{\\mu \\nu}_{\\rm EM} \\partial_i g_{\\mu \\nu} / 2\n", "for i in range(3):\n", " for mu in range(4):\n", " for nu in range(4):\n", " S_tilde_rhsD[i] += alpsqrtgam * T4EMUU[mu][nu] * g4DD_zerotimederiv_dD[mu][nu][i+1] / 2\n", "```" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Table of Contents\n", "$$\\label{toc}$$\n", "\n", "Each family of quantities is constructed within a given function (**boldfaced** below). This notebook is organized as follows\n", "\n", "\n", "1. [Step 1](#importmodules): Import needed NRPy+ & Python modules\n", "1. [Step 2](#stressenergy): Define the GRMHD stress-energy tensor $T^{\\mu\\nu}$ and $T^\\mu{}_\\nu$:\n", " * **compute_T4UU()**, **compute_T4UD()**: \n", "1. [Step 3](#declarevarsconstructgrhdeqs): Construct $T^{\\mu\\nu}$ from GRHD & GRFFE modules with ADM and GRMHD input variables, and construct GRMHD equations from the full GRMHD stress-energy tensor.\n", "1. [Step 4](#code_validation): Code Validation against `GRMHD.equations` NRPy+ module\n", "1. [Step 5](#latex_pdf_output): Output this notebook to $\\LaTeX$-formatted PDF file" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Step 1: Import needed NRPy+ & Python modules \\[Back to [top](#toc)\\]\n", "$$\\label{importmodules}$$" ] }, { "cell_type": "code", "execution_count": 1, "metadata": { "execution": { "iopub.execute_input": "2021-03-07T17:11:45.278752Z", "iopub.status.busy": "2021-03-07T17:11:45.277992Z", "iopub.status.idle": "2021-03-07T17:11:45.606239Z", "shell.execute_reply": "2021-03-07T17:11:45.605620Z" } }, "outputs": [], "source": [ "# Step 1: Import needed core NRPy+ modules\n", "import sympy as sp # SymPy: The Python computer algebra package upon which NRPy+ depends\n", "import indexedexp as ixp # NRPy+: Symbolic indexed expression (e.g., tensors, vectors, etc.) support" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Step 2: Define the GRMHD stress-energy tensor $T^{\\mu\\nu}$ and $T^\\mu{}_\\nu$ \\[Back to [top](#toc)\\]\n", "$$\\label{stressenergy}$$\n", "\n", "Recall from above that\n", "\n", "$$\n", "T^{\\mu\\nu} = T^{\\mu\\nu}_{\\rm GRHD} + T^{\\mu\\nu}_{\\rm GRFFE},\n", "$$\n", "where\n", "\n", "* $T^{\\mu\\nu}_{\\rm GRHD}$ is constructed from the `GRHD.compute_T4UU(gammaDD,betaU,alpha, rho_b,P,epsilon,u4U)` [GRHD](../edit/GRHD/equations.py) [(**tutorial**)](Tutorial-GRHD_Equations-Cartesian.ipynb) function, and\n", "* $T^{\\mu\\nu}_{\\rm GRFFE}$ is constructed from the `GRFFE.compute_TEM4UU(gammaDD,betaU,alpha, smallb4U, smallbsquared,u4U)` [GRFFE](../edit/GRFFE/equations.py) [(**tutorial**)](Tutorial-GRFFE_Equations-Cartesian.ipynb) function\n", "\n", "Since a lowering operation on a sum of tensors is equivalent to the lowering operation applied to the individual tensors in the sum,\n", "\n", "$$\n", "T^\\mu{}_{\\nu} = T^\\mu{}_{\\nu}{}^{\\rm GRHD} + T^\\mu{}_{\\nu}{}^{\\rm GRFFE},\n", "$$\n", "\n", "where\n", "\n", "* $T^\\mu{}_{\\nu}{}^{\\rm GRHD}$ is constructed from the `GRHD.compute_T4UD(gammaDD,betaU,alpha, T4UU)` [GRHD](../edit/GRHD/equations.py) [(**tutorial**)](Tutorial-GRHD_Equations-Cartesian.ipynb) function, and\n", "* $T^{\\mu\\nu}_{\\rm GRFFE}$ is constructed from the `GRFFE.compute_TEM4UD(gammaDD,betaU,alpha, TEM4UU)` [GRFFE](../edit/GRFFE/equations.py) [(**tutorial**)](Tutorial-GRFFE_Equations-Cartesian.ipynb) function." ] }, { "cell_type": "code", "execution_count": 2, "metadata": { "execution": { "iopub.execute_input": "2021-03-07T17:11:45.621186Z", "iopub.status.busy": "2021-03-07T17:11:45.620008Z", "iopub.status.idle": "2021-03-07T17:11:45.626575Z", "shell.execute_reply": "2021-03-07T17:11:45.627094Z" } }, "outputs": [], "source": [ "import GRHD.equations as GRHD\n", "import GRFFE.equations as GRFFE\n", "\n", "# Step 2.a: Define the GRMHD T^{mu nu} (a 4-dimensional tensor)\n", "def compute_GRMHD_T4UU(gammaDD,betaU,alpha, rho_b,P,epsilon,u4U, smallb4U, smallbsquared):\n", " global GRHDT4UU\n", " global GRFFET4UU\n", " global T4UU\n", "\n", " GRHD.compute_T4UU( gammaDD,betaU,alpha, rho_b,P,epsilon,u4U)\n", " GRFFE.compute_TEM4UU(gammaDD,betaU,alpha, smallb4U, smallbsquared,u4U)\n", "\n", " GRHDT4UU = ixp.zerorank2(DIM=4)\n", " GRFFET4UU = ixp.zerorank2(DIM=4)\n", " T4UU = ixp.zerorank2(DIM=4)\n", " for mu in range(4):\n", " for nu in range(4):\n", " GRHDT4UU[mu][nu] = GRHD.T4UU[mu][nu]\n", " GRFFET4UU[mu][nu] = GRFFE.TEM4UU[mu][nu]\n", " T4UU[mu][nu] = GRHD.T4UU[mu][nu] + GRFFE.TEM4UU[mu][nu]\n", "\n", "# Step 2.b: Define T^{mu}_{nu} (a 4-dimensional tensor)\n", "def compute_GRMHD_T4UD(gammaDD,betaU,alpha, GRHDT4UU,GRFFET4UU):\n", " global T4UD\n", "\n", " GRHD.compute_T4UD( gammaDD,betaU,alpha, GRHDT4UU)\n", " GRFFE.compute_TEM4UD(gammaDD,betaU,alpha, GRFFET4UU)\n", "\n", " T4UD = ixp.zerorank2(DIM=4)\n", " for mu in range(4):\n", " for nu in range(4):\n", " T4UD[mu][nu] = GRHD.T4UD[mu][nu] + GRFFE.TEM4UD[mu][nu]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Step 3: Declare ADM and hydrodynamical input variables, and construct all terms in GRMHD equations \\[Back to [top](#toc)\\]\n", "$$\\label{declarevarsconstructgrhdeqs}$$" ] }, { "cell_type": "code", "execution_count": 3, "metadata": { "execution": { "iopub.execute_input": "2021-03-07T17:11:45.644847Z", "iopub.status.busy": "2021-03-07T17:11:45.643516Z", "iopub.status.idle": "2021-03-07T17:11:45.973913Z", "shell.execute_reply": "2021-03-07T17:11:45.974465Z" } }, "outputs": [], "source": [ "# First define hydrodynamical quantities\n", "u4U = ixp.declarerank1(\"u4U\", DIM=4)\n", "rho_b,P,epsilon = sp.symbols('rho_b P epsilon',real=True)\n", "B_tildeU = ixp.declarerank1(\"B_tildeU\", DIM=3)\n", "\n", "# Then ADM quantities\n", "gammaDD = ixp.declarerank2(\"gammaDD\",\"sym01\",DIM=3)\n", "KDD = ixp.declarerank2(\"KDD\" ,\"sym01\",DIM=3)\n", "betaU = ixp.declarerank1(\"betaU\", DIM=3)\n", "alpha = sp.symbols('alpha', real=True)\n", "\n", "# Then numerical constant\n", "sqrt4pi = sp.symbols('sqrt4pi', real=True)\n", "\n", "# First compute smallb4U & smallbsquared from BtildeU, which are needed\n", "# for GRMHD stress-energy tensor T4UU and T4UD:\n", "GRHD.compute_sqrtgammaDET(gammaDD)\n", "GRFFE.compute_B_notildeU(GRHD.sqrtgammaDET, B_tildeU)\n", "GRFFE.compute_smallb4U( gammaDD,betaU,alpha, u4U,GRFFE.B_notildeU, sqrt4pi)\n", "GRFFE.compute_smallbsquared(gammaDD,betaU,alpha, GRFFE.smallb4U)\n", "\n", "# Then compute the GRMHD stress-energy tensor:\n", "compute_GRMHD_T4UU(gammaDD,betaU,alpha, rho_b,P,epsilon,u4U, GRFFE.smallb4U, GRFFE.smallbsquared)\n", "compute_GRMHD_T4UD(gammaDD,betaU,alpha, GRHDT4UU,GRFFET4UU)\n", "\n", "# Compute conservative variables in terms of primitive variables\n", "GRHD.compute_rho_star( alpha, GRHD.sqrtgammaDET, rho_b,u4U)\n", "GRHD.compute_tau_tilde(alpha, GRHD.sqrtgammaDET, T4UU,GRHD.rho_star)\n", "GRHD.compute_S_tildeD( alpha, GRHD.sqrtgammaDET, T4UD)\n", "\n", "# Then compute v^i from u^mu\n", "GRHD.compute_vU_from_u4U__no_speed_limit(u4U)\n", "\n", "# Next compute fluxes of conservative variables\n", "GRHD.compute_rho_star_fluxU( GRHD.vU, GRHD.rho_star)\n", "GRHD.compute_tau_tilde_fluxU(alpha, GRHD.sqrtgammaDET, GRHD.vU,T4UU,GRHD.rho_star)\n", "GRHD.compute_S_tilde_fluxUD( alpha, GRHD.sqrtgammaDET, T4UD)\n", "\n", "# Then declare derivatives & compute g4DD_zerotimederiv_dD\n", "gammaDD_dD = ixp.declarerank3(\"gammaDD_dD\",\"sym01\",DIM=3)\n", "betaU_dD = ixp.declarerank2(\"betaU_dD\" ,\"nosym\",DIM=3)\n", "alpha_dD = ixp.declarerank1(\"alpha_dD\" ,DIM=3)\n", "GRHD.compute_g4DD_zerotimederiv_dD(gammaDD,betaU,alpha, gammaDD_dD,betaU_dD,alpha_dD)\n", "\n", "# Then compute source terms on tau_tilde and S_tilde equations\n", "GRHD.compute_s_source_term(KDD,betaU,alpha, GRHD.sqrtgammaDET,alpha_dD, T4UU)\n", "GRHD.compute_S_tilde_source_termD( alpha, GRHD.sqrtgammaDET,GRHD.g4DD_zerotimederiv_dD, T4UU)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Step 4: Code Validation against `GRMHD.equations` NRPy+ module \\[Back to [top](#toc)\\]\n", "$$\\label{code_validation}$$\n", "\n", "As a code validation check, we verify agreement in the SymPy expressions for the GRHD equations generated in\n", "1. this tutorial versus\n", "2. the NRPy+ [GRMHD.equations](../edit/GRMHD/equations.py) module." ] }, { "cell_type": "code", "execution_count": 4, "metadata": { "execution": { "iopub.execute_input": "2021-03-07T17:11:45.980402Z", "iopub.status.busy": "2021-03-07T17:11:45.979414Z", "iopub.status.idle": "2021-03-07T17:11:45.986147Z", "shell.execute_reply": "2021-03-07T17:11:45.986660Z" } }, "outputs": [], "source": [ "import GRMHD.equations as GRMHD\n", "\n", "# Compute stress-energy tensor T4UU and T4UD:\n", "GRMHD.compute_GRMHD_T4UU(gammaDD,betaU,alpha, rho_b,P,epsilon,u4U, GRFFE.smallb4U, GRFFE.smallbsquared)\n", "GRMHD.compute_GRMHD_T4UD(gammaDD,betaU,alpha, GRMHD.GRHDT4UU,GRMHD.GRFFET4UU)" ] }, { "cell_type": "code", "execution_count": 5, "metadata": { "execution": { "iopub.execute_input": "2021-03-07T17:11:46.016877Z", "iopub.status.busy": "2021-03-07T17:11:46.015450Z", "iopub.status.idle": "2021-03-07T17:11:46.020673Z", "shell.execute_reply": "2021-03-07T17:11:46.020167Z" }, "scrolled": true }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "ALL TESTS PASSED!\n" ] } ], "source": [ "def comp_func(expr1,expr2,basename,prefixname2=\"Ge.\"):\n", " if str(expr1-expr2)!=\"0\":\n", " print(basename+\" - \"+prefixname2+basename+\" = \"+ str(expr1-expr2))\n", " return 1\n", " return 0\n", "\n", "def gfnm(basename,idx1,idx2=None,idx3=None):\n", " if idx2 is None:\n", " return basename+\"[\"+str(idx1)+\"]\"\n", " if idx3 is None:\n", " return basename+\"[\"+str(idx1)+\"][\"+str(idx2)+\"]\"\n", " return basename+\"[\"+str(idx1)+\"][\"+str(idx2)+\"][\"+str(idx3)+\"]\"\n", "\n", "expr_list = []\n", "exprcheck_list = []\n", "namecheck_list = []\n", "\n", "for mu in range(4):\n", " for nu in range(4):\n", " namecheck_list.extend([gfnm(\"GRMHD.GRHDT4UU\",mu,nu),gfnm(\"GRMHD.GRFFET4UU\",mu,nu),\n", " gfnm(\"GRMHD.T4UU\", mu,nu),gfnm(\"GRMHD.T4UD\", mu,nu)])\n", " exprcheck_list.extend([GRMHD.GRHDT4UU[mu][nu],GRMHD.GRFFET4UU[mu][nu],\n", " GRMHD.T4UU[mu][nu], GRMHD.T4UD[mu][nu]])\n", " expr_list.extend([GRHDT4UU[mu][nu],GRFFET4UU[mu][nu],\n", " T4UU[mu][nu], T4UD[mu][nu]])\n", "\n", "num_failures = 0\n", "for i in range(len(expr_list)):\n", " num_failures += comp_func(expr_list[i],exprcheck_list[i],namecheck_list[i])\n", "\n", "import sys\n", "if num_failures == 0:\n", " print(\"ALL TESTS PASSED!\")\n", "else:\n", " print(\"ERROR: \"+str(num_failures)+\" TESTS DID NOT PASS\")\n", " sys.exit(1)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Step 5: 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-GRMHD_Equations-Cartesian.pdf](Tutorial-GRMHD_Equations-Cartesian.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": 6, "metadata": { "execution": { "iopub.execute_input": "2021-03-07T17:11:46.025679Z", "iopub.status.busy": "2021-03-07T17:11:46.024658Z", "iopub.status.idle": "2021-03-07T17:11:49.521226Z", "shell.execute_reply": "2021-03-07T17:11:49.522012Z" } }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Created Tutorial-GRMHD_Equations-Cartesian.tex, and compiled LaTeX file to\n", " PDF file Tutorial-GRMHD_Equations-Cartesian.pdf\n" ] } ], "source": [ "import cmdline_helper as cmd # NRPy+: Multi-platform Python command-line interface\n", "cmd.output_Jupyter_notebook_to_LaTeXed_PDF(\"Tutorial-GRMHD_Equations-Cartesian\")" ] } ], "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.11.1" } }, "nbformat": 4, "nbformat_minor": 4 }