{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Kerr-Newman spacetime\n", "\n", "This notebook demonstrates a few capabilities of SageMath in computations regarding the Kerr-Newman spacetime. The corresponding tools have been developed within the [SageManifolds](https://sagemanifolds.obspm.fr) project." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "*NB:* a version of SageMath at least equal to 9.2 is required to run this notebook:" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "'SageMath version 10.1, Release Date: 2023-08-20'" ] }, "execution_count": 1, "metadata": {}, "output_type": "execute_result" } ], "source": [ "version()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "First we set up the notebook to display mathematical objects using LaTeX rendering:" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "%display latex" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "and we initialize a time counter for benchmarking:" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "import time\n", "comput_time0 = time.perf_counter()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Since some computations are quite heavy, we ask for running them in parallel on 8 threads:" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [], "source": [ "Parallelism().set(nproc=8)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Spacetime manifold\n", "\n", "We declare the Kerr-Newman spacetime (or more precisely the part of it covered by Boyer-Lindquist coordinates) as a 4-dimensional Lorentzian manifold $\\mathcal{M}$:" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [], "source": [ "M = Manifold(4, 'M', latex_name=r'\\mathcal{M}', structure='Lorentzian')" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We then introduce the standard Boyer-Lindquist coordinates as a chart `BL` (for *Boyer-Lindquist*) on $\\mathcal{M}$, via the method `chart()`, the argument of which is a string\n", "(delimited by `r\"...\"` because of the backslash symbols) expressing the coordinates names, their ranges (the default is $(-\\infty,+\\infty)$) and their LaTeX symbols:" ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Chart (M, (t, r, th, ph))\n" ] }, { "data": { "text/html": [ "\\(\\displaystyle \\left(\\mathcal{M},(t, r, {\\theta}, {\\phi})\\right)\\)" ], "text/latex": [ "$\\displaystyle \\left(\\mathcal{M},(t, r, {\\theta}, {\\phi})\\right)$" ], "text/plain": [ "Chart (M, (t, r, th, ph))" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "BL. = M.chart(r\"t r th:(0,pi):\\theta ph:(0,2*pi):\\phi\") \n", "print(BL); BL" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Metric tensor\n", "\n", "The 3 parameters $m$, $a$ and $q$ of the Kerr-Newman spacetime are declared as symbolic variables:" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\left(m, a, q\\right)\\)" ], "text/latex": [ "$\\displaystyle \\left(m, a, q\\right)$" ], "text/plain": [ "(m, a, q)" ] }, "execution_count": 7, "metadata": {}, "output_type": "execute_result" } ], "source": [ "var('m a q')" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We get the (yet undefined) spacetime metric:" ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [], "source": [ "g = M.metric()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

The metric is defined by its components in the coordinate frame associated with Boyer-Lindquist coordinates, which is the current manifold's default frame:

" ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle g = \\left( -\\frac{q^{2} - 2 \\, m r}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} - 1 \\right) \\mathrm{d} t\\otimes \\mathrm{d} t + \\left( \\frac{{\\left(q^{2} - 2 \\, m r\\right)} a \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\right) \\mathrm{d} t\\otimes \\mathrm{d} {\\phi} + \\left( \\frac{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}}{a^{2} + q^{2} - 2 \\, m r + r^{2}} \\right) \\mathrm{d} r\\otimes \\mathrm{d} r + \\left( a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2} \\right) \\mathrm{d} {\\theta}\\otimes \\mathrm{d} {\\theta} + \\left( \\frac{{\\left(q^{2} - 2 \\, m r\\right)} a \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\right) \\mathrm{d} {\\phi}\\otimes \\mathrm{d} t -{\\left(\\frac{{\\left(q^{2} - 2 \\, m r\\right)} a^{2} \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} - a^{2} - r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2} \\mathrm{d} {\\phi}\\otimes \\mathrm{d} {\\phi}\\)" ], "text/latex": [ "$\\displaystyle g = \\left( -\\frac{q^{2} - 2 \\, m r}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} - 1 \\right) \\mathrm{d} t\\otimes \\mathrm{d} t + \\left( \\frac{{\\left(q^{2} - 2 \\, m r\\right)} a \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\right) \\mathrm{d} t\\otimes \\mathrm{d} {\\phi} + \\left( \\frac{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}}{a^{2} + q^{2} - 2 \\, m r + r^{2}} \\right) \\mathrm{d} r\\otimes \\mathrm{d} r + \\left( a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2} \\right) \\mathrm{d} {\\theta}\\otimes \\mathrm{d} {\\theta} + \\left( \\frac{{\\left(q^{2} - 2 \\, m r\\right)} a \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\right) \\mathrm{d} {\\phi}\\otimes \\mathrm{d} t -{\\left(\\frac{{\\left(q^{2} - 2 \\, m r\\right)} a^{2} \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} - a^{2} - r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2} \\mathrm{d} {\\phi}\\otimes \\mathrm{d} {\\phi}$" ], "text/plain": [ "g = (-(q^2 - 2*m*r)/(a^2*cos(th)^2 + r^2) - 1) dt⊗dt + (q^2 - 2*m*r)*a*sin(th)^2/(a^2*cos(th)^2 + r^2) dt⊗dph + (a^2*cos(th)^2 + r^2)/(a^2 + q^2 - 2*m*r + r^2) dr⊗dr + (a^2*cos(th)^2 + r^2) dth⊗dth + (q^2 - 2*m*r)*a*sin(th)^2/(a^2*cos(th)^2 + r^2) dph⊗dt - ((q^2 - 2*m*r)*a^2*sin(th)^2/(a^2*cos(th)^2 + r^2) - a^2 - r^2)*sin(th)^2 dph⊗dph" ] }, "execution_count": 9, "metadata": {}, "output_type": "execute_result" } ], "source": [ "rho2 = r^2 + (a*cos(th))^2\n", "Delta = r^2 -2*m*r + a^2 + q^2\n", "g[0,0] = -1 + (2*m*r-q^2)/rho2\n", "g[0,3] = -a*sin(th)^2*(2*m*r-q^2)/rho2\n", "g[1,1], g[2,2] = rho2/Delta, rho2\n", "g[3,3] = (r^2 + a^2 + (2*m*r-q^2)*(a*sin(th))^2/rho2)*sin(th)^2\n", "g.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

The list of the non-vanishing components:

" ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\begin{array}{lcl} g_{ \\, t \\, t }^{ \\phantom{\\, t}\\phantom{\\, t} } & = & -\\frac{q^{2} - 2 \\, m r}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} - 1 \\\\ g_{ \\, t \\, {\\phi} }^{ \\phantom{\\, t}\\phantom{\\, {\\phi}} } & = & \\frac{{\\left(q^{2} - 2 \\, m r\\right)} a \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\\\ g_{ \\, r \\, r }^{ \\phantom{\\, r}\\phantom{\\, r} } & = & \\frac{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}}{a^{2} + q^{2} - 2 \\, m r + r^{2}} \\\\ g_{ \\, {\\theta} \\, {\\theta} }^{ \\phantom{\\, {\\theta}}\\phantom{\\, {\\theta}} } & = & a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2} \\\\ g_{ \\, {\\phi} \\, t }^{ \\phantom{\\, {\\phi}}\\phantom{\\, t} } & = & \\frac{{\\left(q^{2} - 2 \\, m r\\right)} a \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\\\ g_{ \\, {\\phi} \\, {\\phi} }^{ \\phantom{\\, {\\phi}}\\phantom{\\, {\\phi}} } & = & -{\\left(\\frac{{\\left(q^{2} - 2 \\, m r\\right)} a^{2} \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} - a^{2} - r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2} \\end{array}\\)" ], "text/latex": [ "$\\displaystyle \\begin{array}{lcl} g_{ \\, t \\, t }^{ \\phantom{\\, t}\\phantom{\\, t} } & = & -\\frac{q^{2} - 2 \\, m r}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} - 1 \\\\ g_{ \\, t \\, {\\phi} }^{ \\phantom{\\, t}\\phantom{\\, {\\phi}} } & = & \\frac{{\\left(q^{2} - 2 \\, m r\\right)} a \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\\\ g_{ \\, r \\, r }^{ \\phantom{\\, r}\\phantom{\\, r} } & = & \\frac{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}}{a^{2} + q^{2} - 2 \\, m r + r^{2}} \\\\ g_{ \\, {\\theta} \\, {\\theta} }^{ \\phantom{\\, {\\theta}}\\phantom{\\, {\\theta}} } & = & a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2} \\\\ g_{ \\, {\\phi} \\, t }^{ \\phantom{\\, {\\phi}}\\phantom{\\, t} } & = & \\frac{{\\left(q^{2} - 2 \\, m r\\right)} a \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\\\ g_{ \\, {\\phi} \\, {\\phi} }^{ \\phantom{\\, {\\phi}}\\phantom{\\, {\\phi}} } & = & -{\\left(\\frac{{\\left(q^{2} - 2 \\, m r\\right)} a^{2} \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} - a^{2} - r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2} \\end{array}$" ], "text/plain": [ "g_t,t = -(q^2 - 2*m*r)/(a^2*cos(th)^2 + r^2) - 1 \n", "g_t,ph = (q^2 - 2*m*r)*a*sin(th)^2/(a^2*cos(th)^2 + r^2) \n", "g_r,r = (a^2*cos(th)^2 + r^2)/(a^2 + q^2 - 2*m*r + r^2) \n", "g_th,th = a^2*cos(th)^2 + r^2 \n", "g_ph,t = (q^2 - 2*m*r)*a*sin(th)^2/(a^2*cos(th)^2 + r^2) \n", "g_ph,ph = -((q^2 - 2*m*r)*a^2*sin(th)^2/(a^2*cos(th)^2 + r^2) - a^2 - r^2)*sin(th)^2 " ] }, "execution_count": 10, "metadata": {}, "output_type": "execute_result" } ], "source": [ "g.display_comp()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

The component $g^{tt}$ of the inverse metric:

" ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\frac{a^{4} + 2 \\, a^{2} r^{2} + r^{4} - {\\left(a^{4} + a^{2} q^{2} - 2 \\, a^{2} m r + a^{2} r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}{2 \\, m r^{3} - r^{4} - {\\left(a^{2} + q^{2}\\right)} r^{2} - {\\left(a^{4} + a^{2} q^{2} - 2 \\, a^{2} m r + a^{2} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}}\\)" ], "text/latex": [ "$\\displaystyle \\frac{a^{4} + 2 \\, a^{2} r^{2} + r^{4} - {\\left(a^{4} + a^{2} q^{2} - 2 \\, a^{2} m r + a^{2} r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}{2 \\, m r^{3} - r^{4} - {\\left(a^{2} + q^{2}\\right)} r^{2} - {\\left(a^{4} + a^{2} q^{2} - 2 \\, a^{2} m r + a^{2} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}}$" ], "text/plain": [ "(a^4 + 2*a^2*r^2 + r^4 - (a^4 + a^2*q^2 - 2*a^2*m*r + a^2*r^2)*sin(th)^2)/(2*m*r^3 - r^4 - (a^2 + q^2)*r^2 - (a^4 + a^2*q^2 - 2*a^2*m*r + a^2*r^2)*cos(th)^2)" ] }, "execution_count": 11, "metadata": {}, "output_type": "execute_result" } ], "source": [ "g.inverse()[0,0]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

The lapse function:

" ] }, { "cell_type": "code", "execution_count": 12, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\text{Scalar field on the 4-dimensional Lorentzian manifold M}\\)" ], "text/latex": [ "$\\displaystyle \\text{Scalar field on the 4-dimensional Lorentzian manifold M}$" ], "text/plain": [ "Scalar field on the 4-dimensional Lorentzian manifold M" ] }, "execution_count": 12, "metadata": {}, "output_type": "execute_result" } ], "source": [ "N = 1/sqrt(-(g.inverse()[[0,0]])); N" ] }, { "cell_type": "code", "execution_count": 13, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\begin{array}{llcl} & \\mathcal{M} & \\longrightarrow & \\mathbb{R} \\\\ & \\left(t, r, {\\theta}, {\\phi}\\right) & \\longmapsto & \\frac{\\sqrt{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\sqrt{{\\left| a^{2} + q^{2} - 2 \\, m r + r^{2} \\right|}}}{\\sqrt{a^{4} + 2 \\, a^{2} r^{2} + r^{4} - {\\left(a^{4} + a^{2} q^{2} - 2 \\, a^{2} m r + a^{2} r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}} \\end{array}\\)" ], "text/latex": [ "$\\displaystyle \\begin{array}{llcl} & \\mathcal{M} & \\longrightarrow & \\mathbb{R} \\\\ & \\left(t, r, {\\theta}, {\\phi}\\right) & \\longmapsto & \\frac{\\sqrt{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\sqrt{{\\left| a^{2} + q^{2} - 2 \\, m r + r^{2} \\right|}}}{\\sqrt{a^{4} + 2 \\, a^{2} r^{2} + r^{4} - {\\left(a^{4} + a^{2} q^{2} - 2 \\, a^{2} m r + a^{2} r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}} \\end{array}$" ], "text/plain": [ "M → ℝ\n", "(t, r, th, ph) ↦ sqrt(a^2*cos(th)^2 + r^2)*sqrt(abs(a^2 + q^2 - 2*m*r + r^2))/sqrt(a^4 + 2*a^2*r^2 + r^4 - (a^4 + a^2*q^2 - 2*a^2*m*r + a^2*r^2)*sin(th)^2)" ] }, "execution_count": 13, "metadata": {}, "output_type": "execute_result" } ], "source": [ "N.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Electromagnetic field tensor\n", "\n", "Let us first get the 1-forms $(\\mathrm{d}t, \\mathrm{d}r, \\mathrm{d}\\theta, \\mathrm{d}\\phi)$ as those forming\n", "the coframe associated with Boyer-Lindquist coordinates:" ] }, { "cell_type": "code", "execution_count": 14, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\left(\\mathcal{M}, \\left(\\mathrm{d} t,\\mathrm{d} r,\\mathrm{d} {\\theta},\\mathrm{d} {\\phi}\\right)\\right)\\)" ], "text/latex": [ "$\\displaystyle \\left(\\mathcal{M}, \\left(\\mathrm{d} t,\\mathrm{d} r,\\mathrm{d} {\\theta},\\mathrm{d} {\\phi}\\right)\\right)$" ], "text/plain": [ "Coordinate coframe (M, (dt,dr,dth,dph))" ] }, "execution_count": 14, "metadata": {}, "output_type": "execute_result" } ], "source": [ "BL.coframe()" ] }, { "cell_type": "code", "execution_count": 15, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\left(\\mathrm{d} t, \\mathrm{d} {\\phi}\\right)\\)" ], "text/latex": [ "$\\displaystyle \\left(\\mathrm{d} t, \\mathrm{d} {\\phi}\\right)$" ], "text/plain": [ "(1-form dt on the 4-dimensional Lorentzian manifold M,\n", " 1-form dph on the 4-dimensional Lorentzian manifold M)" ] }, "execution_count": 15, "metadata": {}, "output_type": "execute_result" } ], "source": [ "dt, dr, dth, dph = BL.coframe()[:]\n", "dt, dph" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The electromagnetic 4-potential 1-form $A$ of the Kerr-Newman solution is" ] }, { "cell_type": "code", "execution_count": 16, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle A = \\left( -\\frac{q r}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\right) \\mathrm{d} t + \\left( \\frac{a q r \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\right) \\mathrm{d} {\\phi}\\)" ], "text/latex": [ "$\\displaystyle A = \\left( -\\frac{q r}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\right) \\mathrm{d} t + \\left( \\frac{a q r \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\right) \\mathrm{d} {\\phi}$" ], "text/plain": [ "A = -q*r/(a^2*cos(th)^2 + r^2) dt + a*q*r*sin(th)^2/(a^2*cos(th)^2 + r^2) dph" ] }, "execution_count": 16, "metadata": {}, "output_type": "execute_result" } ], "source": [ "A = - q*r/rho2 * (dt - a*sin(th)^2*dph)\n", "A.set_name('A')\n", "A.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The electromagnetic field tensor $F$ is computed as the exterior derivative of the 4-potential: $F = \\mathrm{d} A$. In Sage, the exterior derivative of a $p$-form is returned by the function `diff()`:" ] }, { "cell_type": "code", "execution_count": 17, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle F = \\left( \\frac{a^{2} q \\cos\\left({\\theta}\\right)^{2} - q r^{2}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} t\\wedge \\mathrm{d} r + \\left( \\frac{2 \\, a^{2} q r \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} t\\wedge \\mathrm{d} {\\theta} + \\left( \\frac{{\\left(a^{3} q \\cos\\left({\\theta}\\right)^{2} - a q r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} r\\wedge \\mathrm{d} {\\phi} + \\left( \\frac{2 \\, {\\left(a^{3} q r + a q r^{3}\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} {\\theta}\\wedge \\mathrm{d} {\\phi}\\)" ], "text/latex": [ "$\\displaystyle F = \\left( \\frac{a^{2} q \\cos\\left({\\theta}\\right)^{2} - q r^{2}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} t\\wedge \\mathrm{d} r + \\left( \\frac{2 \\, a^{2} q r \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} t\\wedge \\mathrm{d} {\\theta} + \\left( \\frac{{\\left(a^{3} q \\cos\\left({\\theta}\\right)^{2} - a q r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} r\\wedge \\mathrm{d} {\\phi} + \\left( \\frac{2 \\, {\\left(a^{3} q r + a q r^{3}\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} {\\theta}\\wedge \\mathrm{d} {\\phi}$" ], "text/plain": [ "F = (a^2*q*cos(th)^2 - q*r^2)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) dt∧dr + 2*a^2*q*r*cos(th)*sin(th)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) dt∧dth + (a^3*q*cos(th)^2 - a*q*r^2)*sin(th)^2/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) dr∧dph + 2*(a^3*q*r + a*q*r^3)*cos(th)*sin(th)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) dth∧dph" ] }, "execution_count": 17, "metadata": {}, "output_type": "execute_result" } ], "source": [ "F = diff(A)\n", "F.set_name('F')\n", "F.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "As a check, let us compare with Eq. (33.5) of [Misner, Thorne & Wheeler (1973)](https://press.princeton.edu/books/ebook/9781400889099/gravitation):" ] }, { "cell_type": "code", "execution_count": 18, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\mathrm{True}\\)" ], "text/latex": [ "$\\displaystyle \\mathrm{True}$" ], "text/plain": [ "True" ] }, "execution_count": 18, "metadata": {}, "output_type": "execute_result" } ], "source": [ "F == q/rho2^2 * (r^2-a^2*cos(th)^2)*dr.wedge( dt - a*sin(th)^2*dph ) \\\n", " + 2*q/rho2^2 * a*r*cos(th)*sin(th)*dth.wedge( (r^2+a^2)*dph - a*dt )" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We can get a short expression of $F$ by factoring the components:" ] }, { "cell_type": "code", "execution_count": 19, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle F = \\frac{{\\left(a \\cos\\left({\\theta}\\right) + r\\right)} {\\left(a \\cos\\left({\\theta}\\right) - r\\right)} q}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\mathrm{d} t\\wedge \\mathrm{d} r + \\frac{2 \\, a^{2} q r \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\mathrm{d} t\\wedge \\mathrm{d} {\\theta} + \\frac{{\\left(a \\cos\\left({\\theta}\\right) + r\\right)} {\\left(a \\cos\\left({\\theta}\\right) - r\\right)} a q \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\mathrm{d} r\\wedge \\mathrm{d} {\\phi} + \\frac{2 \\, {\\left(a^{2} + r^{2}\\right)} a q r \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\mathrm{d} {\\theta}\\wedge \\mathrm{d} {\\phi}\\)" ], "text/latex": [ "$\\displaystyle F = \\frac{{\\left(a \\cos\\left({\\theta}\\right) + r\\right)} {\\left(a \\cos\\left({\\theta}\\right) - r\\right)} q}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\mathrm{d} t\\wedge \\mathrm{d} r + \\frac{2 \\, a^{2} q r \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\mathrm{d} t\\wedge \\mathrm{d} {\\theta} + \\frac{{\\left(a \\cos\\left({\\theta}\\right) + r\\right)} {\\left(a \\cos\\left({\\theta}\\right) - r\\right)} a q \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\mathrm{d} r\\wedge \\mathrm{d} {\\phi} + \\frac{2 \\, {\\left(a^{2} + r^{2}\\right)} a q r \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\mathrm{d} {\\theta}\\wedge \\mathrm{d} {\\phi}$" ], "text/plain": [ "F = (a*cos(th) + r)*(a*cos(th) - r)*q/(a^2*cos(th)^2 + r^2)^2 dt∧dr + 2*a^2*q*r*cos(th)*sin(th)/(a^2*cos(th)^2 + r^2)^2 dt∧dth + (a*cos(th) + r)*(a*cos(th) - r)*a*q*sin(th)^2/(a^2*cos(th)^2 + r^2)^2 dr∧dph + 2*(a^2 + r^2)*a*q*r*cos(th)*sin(th)/(a^2*cos(th)^2 + r^2)^2 dth∧dph" ] }, "execution_count": 19, "metadata": {}, "output_type": "execute_result" } ], "source": [ "F.apply_map(factor)\n", "F.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The list of non-vanishing components of $F$:" ] }, { "cell_type": "code", "execution_count": 20, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\begin{array}{lcl} F_{ \\, t \\, r }^{ \\phantom{\\, t}\\phantom{\\, r} } & = & \\frac{{\\left(a \\cos\\left({\\theta}\\right) + r\\right)} {\\left(a \\cos\\left({\\theta}\\right) - r\\right)} q}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\\\ F_{ \\, t \\, {\\theta} }^{ \\phantom{\\, t}\\phantom{\\, {\\theta}} } & = & \\frac{2 \\, a^{2} q r \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\\\ F_{ \\, r \\, t }^{ \\phantom{\\, r}\\phantom{\\, t} } & = & -\\frac{a^{2} q \\cos\\left({\\theta}\\right)^{2} - q r^{2}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\\\ F_{ \\, r \\, {\\phi} }^{ \\phantom{\\, r}\\phantom{\\, {\\phi}} } & = & \\frac{{\\left(a \\cos\\left({\\theta}\\right) + r\\right)} {\\left(a \\cos\\left({\\theta}\\right) - r\\right)} a q \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\\\ F_{ \\, {\\theta} \\, t }^{ \\phantom{\\, {\\theta}}\\phantom{\\, t} } & = & -\\frac{2 \\, a^{2} q r \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\\\ F_{ \\, {\\theta} \\, {\\phi} }^{ \\phantom{\\, {\\theta}}\\phantom{\\, {\\phi}} } & = & \\frac{2 \\, {\\left(a^{2} + r^{2}\\right)} a q r \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\\\ F_{ \\, {\\phi} \\, r }^{ \\phantom{\\, {\\phi}}\\phantom{\\, r} } & = & -\\frac{{\\left(a^{3} q \\cos\\left({\\theta}\\right)^{2} - a q r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\\\ F_{ \\, {\\phi} \\, {\\theta} }^{ \\phantom{\\, {\\phi}}\\phantom{\\, {\\theta}} } & = & -\\frac{2 \\, {\\left(a^{3} q r + a q r^{3}\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\end{array}\\)" ], "text/latex": [ "$\\displaystyle \\begin{array}{lcl} F_{ \\, t \\, r }^{ \\phantom{\\, t}\\phantom{\\, r} } & = & \\frac{{\\left(a \\cos\\left({\\theta}\\right) + r\\right)} {\\left(a \\cos\\left({\\theta}\\right) - r\\right)} q}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\\\ F_{ \\, t \\, {\\theta} }^{ \\phantom{\\, t}\\phantom{\\, {\\theta}} } & = & \\frac{2 \\, a^{2} q r \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\\\ F_{ \\, r \\, t }^{ \\phantom{\\, r}\\phantom{\\, t} } & = & -\\frac{a^{2} q \\cos\\left({\\theta}\\right)^{2} - q r^{2}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\\\ F_{ \\, r \\, {\\phi} }^{ \\phantom{\\, r}\\phantom{\\, {\\phi}} } & = & \\frac{{\\left(a \\cos\\left({\\theta}\\right) + r\\right)} {\\left(a \\cos\\left({\\theta}\\right) - r\\right)} a q \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\\\ F_{ \\, {\\theta} \\, t }^{ \\phantom{\\, {\\theta}}\\phantom{\\, t} } & = & -\\frac{2 \\, a^{2} q r \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\\\ F_{ \\, {\\theta} \\, {\\phi} }^{ \\phantom{\\, {\\theta}}\\phantom{\\, {\\phi}} } & = & \\frac{2 \\, {\\left(a^{2} + r^{2}\\right)} a q r \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{2}} \\\\ F_{ \\, {\\phi} \\, r }^{ \\phantom{\\, {\\phi}}\\phantom{\\, r} } & = & -\\frac{{\\left(a^{3} q \\cos\\left({\\theta}\\right)^{2} - a q r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\\\ F_{ \\, {\\phi} \\, {\\theta} }^{ \\phantom{\\, {\\phi}}\\phantom{\\, {\\theta}} } & = & -\\frac{2 \\, {\\left(a^{3} q r + a q r^{3}\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\end{array}$" ], "text/plain": [ "F_t,r = (a*cos(th) + r)*(a*cos(th) - r)*q/(a^2*cos(th)^2 + r^2)^2 \n", "F_t,th = 2*a^2*q*r*cos(th)*sin(th)/(a^2*cos(th)^2 + r^2)^2 \n", "F_r,t = -(a^2*q*cos(th)^2 - q*r^2)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) \n", "F_r,ph = (a*cos(th) + r)*(a*cos(th) - r)*a*q*sin(th)^2/(a^2*cos(th)^2 + r^2)^2 \n", "F_th,t = -2*a^2*q*r*cos(th)*sin(th)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) \n", "F_th,ph = 2*(a^2 + r^2)*a*q*r*cos(th)*sin(th)/(a^2*cos(th)^2 + r^2)^2 \n", "F_ph,r = -(a^3*q*cos(th)^2 - a*q*r^2)*sin(th)^2/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) \n", "F_ph,th = -2*(a^3*q*r + a*q*r^3)*cos(th)*sin(th)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) " ] }, "execution_count": 20, "metadata": {}, "output_type": "execute_result" } ], "source": [ "F.display_comp()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

The Hodge dual of $F$:

" ] }, { "cell_type": "code", "execution_count": 21, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\star F = \\left( \\frac{2 \\, a q r \\sqrt{\\cos\\left({\\theta}\\right) + 1} \\sqrt{-\\cos\\left({\\theta}\\right) + 1} \\cos\\left({\\theta}\\right)}{{\\left(a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}\\right)} \\sin\\left({\\theta}\\right)} \\right) \\mathrm{d} t\\wedge \\mathrm{d} r + \\left( -\\frac{{\\left(a^{3} q \\cos\\left({\\theta}\\right)^{2} - a q r^{2}\\right)} \\sqrt{\\cos\\left({\\theta}\\right) + 1} \\sqrt{-\\cos\\left({\\theta}\\right) + 1}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} t\\wedge \\mathrm{d} {\\theta} + \\left( \\frac{2 \\, a^{2} q r \\sqrt{\\cos\\left({\\theta}\\right) + 1} \\sqrt{-\\cos\\left({\\theta}\\right) + 1} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} r\\wedge \\mathrm{d} {\\phi} + \\left( -\\frac{{\\left(a^{4} q - q r^{4} - {\\left(a^{4} q + a^{2} q r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}\\right)} \\sqrt{\\cos\\left({\\theta}\\right) + 1} \\sqrt{-\\cos\\left({\\theta}\\right) + 1}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} {\\theta}\\wedge \\mathrm{d} {\\phi}\\)" ], "text/latex": [ "$\\displaystyle \\star F = \\left( \\frac{2 \\, a q r \\sqrt{\\cos\\left({\\theta}\\right) + 1} \\sqrt{-\\cos\\left({\\theta}\\right) + 1} \\cos\\left({\\theta}\\right)}{{\\left(a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}\\right)} \\sin\\left({\\theta}\\right)} \\right) \\mathrm{d} t\\wedge \\mathrm{d} r + \\left( -\\frac{{\\left(a^{3} q \\cos\\left({\\theta}\\right)^{2} - a q r^{2}\\right)} \\sqrt{\\cos\\left({\\theta}\\right) + 1} \\sqrt{-\\cos\\left({\\theta}\\right) + 1}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} t\\wedge \\mathrm{d} {\\theta} + \\left( \\frac{2 \\, a^{2} q r \\sqrt{\\cos\\left({\\theta}\\right) + 1} \\sqrt{-\\cos\\left({\\theta}\\right) + 1} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} r\\wedge \\mathrm{d} {\\phi} + \\left( -\\frac{{\\left(a^{4} q - q r^{4} - {\\left(a^{4} q + a^{2} q r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}\\right)} \\sqrt{\\cos\\left({\\theta}\\right) + 1} \\sqrt{-\\cos\\left({\\theta}\\right) + 1}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} {\\theta}\\wedge \\mathrm{d} {\\phi}$" ], "text/plain": [ "*F = 2*a*q*r*sqrt(cos(th) + 1)*sqrt(-cos(th) + 1)*cos(th)/((a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4)*sin(th)) dt∧dr - (a^3*q*cos(th)^2 - a*q*r^2)*sqrt(cos(th) + 1)*sqrt(-cos(th) + 1)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) dt∧dth + 2*a^2*q*r*sqrt(cos(th) + 1)*sqrt(-cos(th) + 1)*cos(th)*sin(th)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) dr∧dph - (a^4*q - q*r^4 - (a^4*q + a^2*q*r^2)*sin(th)^2)*sqrt(cos(th) + 1)*sqrt(-cos(th) + 1)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) dth∧dph" ] }, "execution_count": 21, "metadata": {}, "output_type": "execute_result" } ], "source": [ "star_F = F.hodge_dual(g)\n", "star_F.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Maxwell equations

\n", "\n", "

Let us check that $F$ obeys the two (source-free) Maxwell equations:

" ] }, { "cell_type": "code", "execution_count": 22, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\mathrm{d}F = 0\\)" ], "text/latex": [ "$\\displaystyle \\mathrm{d}F = 0$" ], "text/plain": [ "dF = 0" ] }, "execution_count": 22, "metadata": {}, "output_type": "execute_result" } ], "source": [ "diff(F).display()" ] }, { "cell_type": "code", "execution_count": 23, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\mathrm{d}\\star F = 0\\)" ], "text/latex": [ "$\\displaystyle \\mathrm{d}\\star F = 0$" ], "text/plain": [ "d*F = 0" ] }, "execution_count": 23, "metadata": {}, "output_type": "execute_result" } ], "source": [ "diff(star_F).display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Levi-Civita Connection

\n", "\n", "

The Levi-Civita connection $\\nabla$ associated with $g$:

" ] }, { "cell_type": "code", "execution_count": 24, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Levi-Civita connection nabla_g associated with the Lorentzian metric g on the 4-dimensional Lorentzian manifold M\n" ] } ], "source": [ "nabla = g.connection()\n", "print(nabla)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Let us verify that the covariant derivative of $g$ with respect to $\\nabla$ vanishes identically:

" ] }, { "cell_type": "code", "execution_count": 25, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\mathrm{True}\\)" ], "text/latex": [ "$\\displaystyle \\mathrm{True}$" ], "text/plain": [ "True" ] }, "execution_count": 25, "metadata": {}, "output_type": "execute_result" } ], "source": [ "nabla(g) == 0" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Another view of the above property:

" ] }, { "cell_type": "code", "execution_count": 26, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\nabla_{g} g = 0\\)" ], "text/latex": [ "$\\displaystyle \\nabla_{g} g = 0$" ], "text/plain": [ "nabla_g(g) = 0" ] }, "execution_count": 26, "metadata": {}, "output_type": "execute_result" } ], "source": [ "nabla(g).display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

The nonzero Christoffel symbols (skipping those that can be deduced by symmetry of the last two indices):

" ] }, { "cell_type": "code", "execution_count": 27, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\begin{array}{lcl} \\Gamma_{ \\phantom{\\, t} \\, t \\, r }^{ \\, t \\phantom{\\, t} \\phantom{\\, r} } & = & \\frac{a^{4} m + a^{2} q^{2} r + q^{2} r^{3} - m r^{4} - {\\left(a^{4} m + a^{2} m r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}{2 \\, m r^{5} - r^{6} - {\\left(a^{2} + q^{2}\\right)} r^{4} - {\\left(a^{6} + a^{4} q^{2} - 2 \\, a^{4} m r + a^{4} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + 2 \\, {\\left(2 \\, a^{2} m r^{3} - a^{2} r^{4} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\\\ \\Gamma_{ \\phantom{\\, t} \\, t \\, {\\theta} }^{ \\, t \\phantom{\\, t} \\phantom{\\, {\\theta}} } & = & \\frac{{\\left(a^{2} q^{2} - 2 \\, a^{2} m r\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\\\ \\Gamma_{ \\phantom{\\, t} \\, r \\, {\\phi} }^{ \\, t \\phantom{\\, r} \\phantom{\\, {\\phi}} } & = & -\\frac{a^{3} q^{2} r - a^{3} m r^{2} + 2 \\, a q^{2} r^{3} - 3 \\, a m r^{4} - {\\left(a^{5} m + a^{3} q^{2} r - a^{3} m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + {\\left(a^{5} m - 2 \\, a q^{2} r^{3} + 3 \\, a m r^{4}\\right)} \\cos\\left({\\theta}\\right)^{2}}{2 \\, m r^{5} - r^{6} - {\\left(a^{2} + q^{2}\\right)} r^{4} - {\\left(a^{6} + a^{4} q^{2} - 2 \\, a^{4} m r + a^{4} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + 2 \\, {\\left(2 \\, a^{2} m r^{3} - a^{2} r^{4} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\\\ \\Gamma_{ \\phantom{\\, t} \\, {\\theta} \\, {\\phi} }^{ \\, t \\phantom{\\, {\\theta}} \\phantom{\\, {\\phi}} } & = & \\frac{{\\left(a^{5} q^{2} - 2 \\, a^{5} m r\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)^{5} - {\\left(a^{5} q^{2} - 2 \\, a^{5} m r + a^{3} q^{2} r^{2} - 2 \\, a^{3} m r^{3}\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)^{3}}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\\\ \\Gamma_{ \\phantom{\\, r} \\, t \\, t }^{ \\, r \\phantom{\\, t} \\phantom{\\, t} } & = & \\frac{m r^{4} - {\\left(2 \\, m^{2} + q^{2}\\right)} r^{3} + {\\left(a^{2} m + 3 \\, m q^{2}\\right)} r^{2} - {\\left(a^{4} m + a^{2} m q^{2} - 2 \\, a^{2} m^{2} r + a^{2} m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2} - {\\left(a^{2} q^{2} + q^{4}\\right)} r}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\\\ \\Gamma_{ \\phantom{\\, r} \\, t \\, {\\phi} }^{ \\, r \\phantom{\\, t} \\phantom{\\, {\\phi}} } & = & -\\frac{{\\left(a m r^{4} - {\\left(2 \\, a m^{2} + a q^{2}\\right)} r^{3} + {\\left(a^{3} m + 3 \\, a m q^{2}\\right)} r^{2} - {\\left(a^{5} m + a^{3} m q^{2} - 2 \\, a^{3} m^{2} r + a^{3} m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2} - {\\left(a^{3} q^{2} + a q^{4}\\right)} r\\right)} \\sin\\left({\\theta}\\right)^{2}}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\\\ \\Gamma_{ \\phantom{\\, r} \\, r \\, r }^{ \\, r \\phantom{\\, r} \\phantom{\\, r} } & = & -\\frac{a^{2} m + q^{2} r - m r^{2} - {\\left(a^{2} m - a^{2} r\\right)} \\sin\\left({\\theta}\\right)^{2}}{2 \\, m r^{3} - r^{4} - {\\left(a^{2} + q^{2}\\right)} r^{2} - {\\left(a^{4} + a^{2} q^{2} - 2 \\, a^{2} m r + a^{2} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\\\ \\Gamma_{ \\phantom{\\, r} \\, r \\, {\\theta} }^{ \\, r \\phantom{\\, r} \\phantom{\\, {\\theta}} } & = & -\\frac{a^{2} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\\\ \\Gamma_{ \\phantom{\\, r} \\, {\\theta} \\, {\\theta} }^{ \\, r \\phantom{\\, {\\theta}} \\phantom{\\, {\\theta}} } & = & \\frac{2 \\, m r^{2} - r^{3} - {\\left(a^{2} + q^{2}\\right)} r}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\\\ \\Gamma_{ \\phantom{\\, r} \\, {\\phi} \\, {\\phi} }^{ \\, r \\phantom{\\, {\\phi}} \\phantom{\\, {\\phi}} } & = & \\frac{{\\left(a^{2} m r^{4} - {\\left(2 \\, a^{2} m^{2} + a^{2} q^{2}\\right)} r^{3} + {\\left(a^{4} m + 3 \\, a^{2} m q^{2}\\right)} r^{2} - {\\left(a^{6} m + a^{4} m q^{2} - 2 \\, a^{4} m^{2} r + a^{4} m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2} - {\\left(a^{4} q^{2} + a^{2} q^{4}\\right)} r\\right)} \\sin\\left({\\theta}\\right)^{4} + {\\left(2 \\, m r^{6} - r^{7} - {\\left(a^{2} + q^{2}\\right)} r^{5} + {\\left(2 \\, a^{4} m r^{2} - a^{4} r^{3} - {\\left(a^{6} + a^{4} q^{2}\\right)} r\\right)} \\cos\\left({\\theta}\\right)^{4} + 2 \\, {\\left(2 \\, a^{2} m r^{4} - a^{2} r^{5} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{3}\\right)} \\cos\\left({\\theta}\\right)^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\\\ \\Gamma_{ \\phantom{\\, {\\theta}} \\, t \\, t }^{ \\, {\\theta} \\phantom{\\, t} \\phantom{\\, t} } & = & \\frac{{\\left(a^{2} q^{2} - 2 \\, a^{2} m r\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\\\ \\Gamma_{ \\phantom{\\, {\\theta}} \\, t \\, {\\phi} }^{ \\, {\\theta} \\phantom{\\, t} \\phantom{\\, {\\phi}} } & = & -\\frac{{\\left(a^{3} q^{2} - 2 \\, a^{3} m r + a q^{2} r^{2} - 2 \\, a m r^{3}\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\\\ \\Gamma_{ \\phantom{\\, {\\theta}} \\, r \\, r }^{ \\, {\\theta} \\phantom{\\, r} \\phantom{\\, r} } & = & -\\frac{a^{2} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{2 \\, m r^{3} - r^{4} - {\\left(a^{2} + q^{2}\\right)} r^{2} - {\\left(a^{4} + a^{2} q^{2} - 2 \\, a^{2} m r + a^{2} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\\\ \\Gamma_{ \\phantom{\\, {\\theta}} \\, r \\, {\\theta} }^{ \\, {\\theta} \\phantom{\\, r} \\phantom{\\, {\\theta}} } & = & \\frac{r}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\\\ \\Gamma_{ \\phantom{\\, {\\theta}} \\, {\\theta} \\, {\\theta} }^{ \\, {\\theta} \\phantom{\\, {\\theta}} \\phantom{\\, {\\theta}} } & = & -\\frac{a^{2} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\\\ \\Gamma_{ \\phantom{\\, {\\theta}} \\, {\\phi} \\, {\\phi} }^{ \\, {\\theta} \\phantom{\\, {\\phi}} \\phantom{\\, {\\phi}} } & = & -\\frac{{\\left({\\left(a^{6} + a^{4} q^{2} - 2 \\, a^{4} m r + a^{4} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{5} - 2 \\, {\\left(2 \\, a^{2} m r^{3} - a^{2} r^{4} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{3} - {\\left(a^{4} q^{2} - 2 \\, a^{4} m r + 2 \\, a^{2} q^{2} r^{2} - 4 \\, a^{2} m r^{3} - a^{2} r^{4} - r^{6}\\right)} \\cos\\left({\\theta}\\right)\\right)} \\sin\\left({\\theta}\\right)}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\\\ \\Gamma_{ \\phantom{\\, {\\phi}} \\, t \\, r }^{ \\, {\\phi} \\phantom{\\, t} \\phantom{\\, r} } & = & \\frac{a^{3} m \\cos\\left({\\theta}\\right)^{2} + a q^{2} r - a m r^{2}}{2 \\, m r^{5} - r^{6} - {\\left(a^{2} + q^{2}\\right)} r^{4} - {\\left(a^{6} + a^{4} q^{2} - 2 \\, a^{4} m r + a^{4} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + 2 \\, {\\left(2 \\, a^{2} m r^{3} - a^{2} r^{4} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\\\ \\Gamma_{ \\phantom{\\, {\\phi}} \\, t \\, {\\theta} }^{ \\, {\\phi} \\phantom{\\, t} \\phantom{\\, {\\theta}} } & = & \\frac{{\\left(a q^{2} - 2 \\, a m r\\right)} \\cos\\left({\\theta}\\right)}{{\\left(a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}\\right)} \\sin\\left({\\theta}\\right)} \\\\ \\Gamma_{ \\phantom{\\, {\\phi}} \\, r \\, {\\phi} }^{ \\, {\\phi} \\phantom{\\, r} \\phantom{\\, {\\phi}} } & = & -\\frac{a^{2} q^{2} r - a^{2} m r^{2} + q^{2} r^{3} - 2 \\, m r^{4} + r^{5} - {\\left(a^{4} m - a^{4} r\\right)} \\cos\\left({\\theta}\\right)^{4} + {\\left(a^{4} m - a^{2} m r^{2} + 2 \\, a^{2} r^{3}\\right)} \\cos\\left({\\theta}\\right)^{2}}{2 \\, m r^{5} - r^{6} - {\\left(a^{2} + q^{2}\\right)} r^{4} - {\\left(a^{6} + a^{4} q^{2} - 2 \\, a^{4} m r + a^{4} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + 2 \\, {\\left(2 \\, a^{2} m r^{3} - a^{2} r^{4} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\\\ \\Gamma_{ \\phantom{\\, {\\phi}} \\, {\\theta} \\, {\\phi} }^{ \\, {\\phi} \\phantom{\\, {\\theta}} \\phantom{\\, {\\phi}} } & = & \\frac{a^{4} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)^{4} - {\\left(2 \\, a^{4} + a^{2} q^{2} - 2 \\, a^{2} m r + 2 \\, a^{2} r^{2}\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)^{2} + {\\left(a^{4} + 2 \\, a^{2} r^{2} + r^{4}\\right)} \\cos\\left({\\theta}\\right)}{{\\left(a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}\\right)} \\sin\\left({\\theta}\\right)} \\end{array}\\)" ], "text/latex": [ "$\\displaystyle \\begin{array}{lcl} \\Gamma_{ \\phantom{\\, t} \\, t \\, r }^{ \\, t \\phantom{\\, t} \\phantom{\\, r} } & = & \\frac{a^{4} m + a^{2} q^{2} r + q^{2} r^{3} - m r^{4} - {\\left(a^{4} m + a^{2} m r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}{2 \\, m r^{5} - r^{6} - {\\left(a^{2} + q^{2}\\right)} r^{4} - {\\left(a^{6} + a^{4} q^{2} - 2 \\, a^{4} m r + a^{4} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + 2 \\, {\\left(2 \\, a^{2} m r^{3} - a^{2} r^{4} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\\\ \\Gamma_{ \\phantom{\\, t} \\, t \\, {\\theta} }^{ \\, t \\phantom{\\, t} \\phantom{\\, {\\theta}} } & = & \\frac{{\\left(a^{2} q^{2} - 2 \\, a^{2} m r\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\\\ \\Gamma_{ \\phantom{\\, t} \\, r \\, {\\phi} }^{ \\, t \\phantom{\\, r} \\phantom{\\, {\\phi}} } & = & -\\frac{a^{3} q^{2} r - a^{3} m r^{2} + 2 \\, a q^{2} r^{3} - 3 \\, a m r^{4} - {\\left(a^{5} m + a^{3} q^{2} r - a^{3} m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + {\\left(a^{5} m - 2 \\, a q^{2} r^{3} + 3 \\, a m r^{4}\\right)} \\cos\\left({\\theta}\\right)^{2}}{2 \\, m r^{5} - r^{6} - {\\left(a^{2} + q^{2}\\right)} r^{4} - {\\left(a^{6} + a^{4} q^{2} - 2 \\, a^{4} m r + a^{4} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + 2 \\, {\\left(2 \\, a^{2} m r^{3} - a^{2} r^{4} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\\\ \\Gamma_{ \\phantom{\\, t} \\, {\\theta} \\, {\\phi} }^{ \\, t \\phantom{\\, {\\theta}} \\phantom{\\, {\\phi}} } & = & \\frac{{\\left(a^{5} q^{2} - 2 \\, a^{5} m r\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)^{5} - {\\left(a^{5} q^{2} - 2 \\, a^{5} m r + a^{3} q^{2} r^{2} - 2 \\, a^{3} m r^{3}\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)^{3}}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\\\ \\Gamma_{ \\phantom{\\, r} \\, t \\, t }^{ \\, r \\phantom{\\, t} \\phantom{\\, t} } & = & \\frac{m r^{4} - {\\left(2 \\, m^{2} + q^{2}\\right)} r^{3} + {\\left(a^{2} m + 3 \\, m q^{2}\\right)} r^{2} - {\\left(a^{4} m + a^{2} m q^{2} - 2 \\, a^{2} m^{2} r + a^{2} m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2} - {\\left(a^{2} q^{2} + q^{4}\\right)} r}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\\\ \\Gamma_{ \\phantom{\\, r} \\, t \\, {\\phi} }^{ \\, r \\phantom{\\, t} \\phantom{\\, {\\phi}} } & = & -\\frac{{\\left(a m r^{4} - {\\left(2 \\, a m^{2} + a q^{2}\\right)} r^{3} + {\\left(a^{3} m + 3 \\, a m q^{2}\\right)} r^{2} - {\\left(a^{5} m + a^{3} m q^{2} - 2 \\, a^{3} m^{2} r + a^{3} m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2} - {\\left(a^{3} q^{2} + a q^{4}\\right)} r\\right)} \\sin\\left({\\theta}\\right)^{2}}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\\\ \\Gamma_{ \\phantom{\\, r} \\, r \\, r }^{ \\, r \\phantom{\\, r} \\phantom{\\, r} } & = & -\\frac{a^{2} m + q^{2} r - m r^{2} - {\\left(a^{2} m - a^{2} r\\right)} \\sin\\left({\\theta}\\right)^{2}}{2 \\, m r^{3} - r^{4} - {\\left(a^{2} + q^{2}\\right)} r^{2} - {\\left(a^{4} + a^{2} q^{2} - 2 \\, a^{2} m r + a^{2} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\\\ \\Gamma_{ \\phantom{\\, r} \\, r \\, {\\theta} }^{ \\, r \\phantom{\\, r} \\phantom{\\, {\\theta}} } & = & -\\frac{a^{2} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\\\ \\Gamma_{ \\phantom{\\, r} \\, {\\theta} \\, {\\theta} }^{ \\, r \\phantom{\\, {\\theta}} \\phantom{\\, {\\theta}} } & = & \\frac{2 \\, m r^{2} - r^{3} - {\\left(a^{2} + q^{2}\\right)} r}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\\\ \\Gamma_{ \\phantom{\\, r} \\, {\\phi} \\, {\\phi} }^{ \\, r \\phantom{\\, {\\phi}} \\phantom{\\, {\\phi}} } & = & \\frac{{\\left(a^{2} m r^{4} - {\\left(2 \\, a^{2} m^{2} + a^{2} q^{2}\\right)} r^{3} + {\\left(a^{4} m + 3 \\, a^{2} m q^{2}\\right)} r^{2} - {\\left(a^{6} m + a^{4} m q^{2} - 2 \\, a^{4} m^{2} r + a^{4} m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2} - {\\left(a^{4} q^{2} + a^{2} q^{4}\\right)} r\\right)} \\sin\\left({\\theta}\\right)^{4} + {\\left(2 \\, m r^{6} - r^{7} - {\\left(a^{2} + q^{2}\\right)} r^{5} + {\\left(2 \\, a^{4} m r^{2} - a^{4} r^{3} - {\\left(a^{6} + a^{4} q^{2}\\right)} r\\right)} \\cos\\left({\\theta}\\right)^{4} + 2 \\, {\\left(2 \\, a^{2} m r^{4} - a^{2} r^{5} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{3}\\right)} \\cos\\left({\\theta}\\right)^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\\\ \\Gamma_{ \\phantom{\\, {\\theta}} \\, t \\, t }^{ \\, {\\theta} \\phantom{\\, t} \\phantom{\\, t} } & = & \\frac{{\\left(a^{2} q^{2} - 2 \\, a^{2} m r\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\\\ \\Gamma_{ \\phantom{\\, {\\theta}} \\, t \\, {\\phi} }^{ \\, {\\theta} \\phantom{\\, t} \\phantom{\\, {\\phi}} } & = & -\\frac{{\\left(a^{3} q^{2} - 2 \\, a^{3} m r + a q^{2} r^{2} - 2 \\, a m r^{3}\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\\\ \\Gamma_{ \\phantom{\\, {\\theta}} \\, r \\, r }^{ \\, {\\theta} \\phantom{\\, r} \\phantom{\\, r} } & = & -\\frac{a^{2} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{2 \\, m r^{3} - r^{4} - {\\left(a^{2} + q^{2}\\right)} r^{2} - {\\left(a^{4} + a^{2} q^{2} - 2 \\, a^{2} m r + a^{2} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\\\ \\Gamma_{ \\phantom{\\, {\\theta}} \\, r \\, {\\theta} }^{ \\, {\\theta} \\phantom{\\, r} \\phantom{\\, {\\theta}} } & = & \\frac{r}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\\\ \\Gamma_{ \\phantom{\\, {\\theta}} \\, {\\theta} \\, {\\theta} }^{ \\, {\\theta} \\phantom{\\, {\\theta}} \\phantom{\\, {\\theta}} } & = & -\\frac{a^{2} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\\\ \\Gamma_{ \\phantom{\\, {\\theta}} \\, {\\phi} \\, {\\phi} }^{ \\, {\\theta} \\phantom{\\, {\\phi}} \\phantom{\\, {\\phi}} } & = & -\\frac{{\\left({\\left(a^{6} + a^{4} q^{2} - 2 \\, a^{4} m r + a^{4} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{5} - 2 \\, {\\left(2 \\, a^{2} m r^{3} - a^{2} r^{4} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{3} - {\\left(a^{4} q^{2} - 2 \\, a^{4} m r + 2 \\, a^{2} q^{2} r^{2} - 4 \\, a^{2} m r^{3} - a^{2} r^{4} - r^{6}\\right)} \\cos\\left({\\theta}\\right)\\right)} \\sin\\left({\\theta}\\right)}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\\\ \\Gamma_{ \\phantom{\\, {\\phi}} \\, t \\, r }^{ \\, {\\phi} \\phantom{\\, t} \\phantom{\\, r} } & = & \\frac{a^{3} m \\cos\\left({\\theta}\\right)^{2} + a q^{2} r - a m r^{2}}{2 \\, m r^{5} - r^{6} - {\\left(a^{2} + q^{2}\\right)} r^{4} - {\\left(a^{6} + a^{4} q^{2} - 2 \\, a^{4} m r + a^{4} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + 2 \\, {\\left(2 \\, a^{2} m r^{3} - a^{2} r^{4} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\\\ \\Gamma_{ \\phantom{\\, {\\phi}} \\, t \\, {\\theta} }^{ \\, {\\phi} \\phantom{\\, t} \\phantom{\\, {\\theta}} } & = & \\frac{{\\left(a q^{2} - 2 \\, a m r\\right)} \\cos\\left({\\theta}\\right)}{{\\left(a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}\\right)} \\sin\\left({\\theta}\\right)} \\\\ \\Gamma_{ \\phantom{\\, {\\phi}} \\, r \\, {\\phi} }^{ \\, {\\phi} \\phantom{\\, r} \\phantom{\\, {\\phi}} } & = & -\\frac{a^{2} q^{2} r - a^{2} m r^{2} + q^{2} r^{3} - 2 \\, m r^{4} + r^{5} - {\\left(a^{4} m - a^{4} r\\right)} \\cos\\left({\\theta}\\right)^{4} + {\\left(a^{4} m - a^{2} m r^{2} + 2 \\, a^{2} r^{3}\\right)} \\cos\\left({\\theta}\\right)^{2}}{2 \\, m r^{5} - r^{6} - {\\left(a^{2} + q^{2}\\right)} r^{4} - {\\left(a^{6} + a^{4} q^{2} - 2 \\, a^{4} m r + a^{4} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + 2 \\, {\\left(2 \\, a^{2} m r^{3} - a^{2} r^{4} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\\\ \\Gamma_{ \\phantom{\\, {\\phi}} \\, {\\theta} \\, {\\phi} }^{ \\, {\\phi} \\phantom{\\, {\\theta}} \\phantom{\\, {\\phi}} } & = & \\frac{a^{4} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)^{4} - {\\left(2 \\, a^{4} + a^{2} q^{2} - 2 \\, a^{2} m r + 2 \\, a^{2} r^{2}\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)^{2} + {\\left(a^{4} + 2 \\, a^{2} r^{2} + r^{4}\\right)} \\cos\\left({\\theta}\\right)}{{\\left(a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}\\right)} \\sin\\left({\\theta}\\right)} \\end{array}$" ], "text/plain": [ "Gam^t_t,r = (a^4*m + a^2*q^2*r + q^2*r^3 - m*r^4 - (a^4*m + a^2*m*r^2)*sin(th)^2)/(2*m*r^5 - r^6 - (a^2 + q^2)*r^4 - (a^6 + a^4*q^2 - 2*a^4*m*r + a^4*r^2)*cos(th)^4 + 2*(2*a^2*m*r^3 - a^2*r^4 - (a^4 + a^2*q^2)*r^2)*cos(th)^2) \n", "Gam^t_t,th = (a^2*q^2 - 2*a^2*m*r)*cos(th)*sin(th)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) \n", "Gam^t_r,ph = -(a^3*q^2*r - a^3*m*r^2 + 2*a*q^2*r^3 - 3*a*m*r^4 - (a^5*m + a^3*q^2*r - a^3*m*r^2)*cos(th)^4 + (a^5*m - 2*a*q^2*r^3 + 3*a*m*r^4)*cos(th)^2)/(2*m*r^5 - r^6 - (a^2 + q^2)*r^4 - (a^6 + a^4*q^2 - 2*a^4*m*r + a^4*r^2)*cos(th)^4 + 2*(2*a^2*m*r^3 - a^2*r^4 - (a^4 + a^2*q^2)*r^2)*cos(th)^2) \n", "Gam^t_th,ph = ((a^5*q^2 - 2*a^5*m*r)*cos(th)*sin(th)^5 - (a^5*q^2 - 2*a^5*m*r + a^3*q^2*r^2 - 2*a^3*m*r^3)*cos(th)*sin(th)^3)/(a^6*cos(th)^6 + 3*a^4*r^2*cos(th)^4 + 3*a^2*r^4*cos(th)^2 + r^6) \n", "Gam^r_t,t = (m*r^4 - (2*m^2 + q^2)*r^3 + (a^2*m + 3*m*q^2)*r^2 - (a^4*m + a^2*m*q^2 - 2*a^2*m^2*r + a^2*m*r^2)*cos(th)^2 - (a^2*q^2 + q^4)*r)/(a^6*cos(th)^6 + 3*a^4*r^2*cos(th)^4 + 3*a^2*r^4*cos(th)^2 + r^6) \n", "Gam^r_t,ph = -(a*m*r^4 - (2*a*m^2 + a*q^2)*r^3 + (a^3*m + 3*a*m*q^2)*r^2 - (a^5*m + a^3*m*q^2 - 2*a^3*m^2*r + a^3*m*r^2)*cos(th)^2 - (a^3*q^2 + a*q^4)*r)*sin(th)^2/(a^6*cos(th)^6 + 3*a^4*r^2*cos(th)^4 + 3*a^2*r^4*cos(th)^2 + r^6) \n", "Gam^r_r,r = -(a^2*m + q^2*r - m*r^2 - (a^2*m - a^2*r)*sin(th)^2)/(2*m*r^3 - r^4 - (a^2 + q^2)*r^2 - (a^4 + a^2*q^2 - 2*a^2*m*r + a^2*r^2)*cos(th)^2) \n", "Gam^r_r,th = -a^2*cos(th)*sin(th)/(a^2*cos(th)^2 + r^2) \n", "Gam^r_th,th = (2*m*r^2 - r^3 - (a^2 + q^2)*r)/(a^2*cos(th)^2 + r^2) \n", "Gam^r_ph,ph = ((a^2*m*r^4 - (2*a^2*m^2 + a^2*q^2)*r^3 + (a^4*m + 3*a^2*m*q^2)*r^2 - (a^6*m + a^4*m*q^2 - 2*a^4*m^2*r + a^4*m*r^2)*cos(th)^2 - (a^4*q^2 + a^2*q^4)*r)*sin(th)^4 + (2*m*r^6 - r^7 - (a^2 + q^2)*r^5 + (2*a^4*m*r^2 - a^4*r^3 - (a^6 + a^4*q^2)*r)*cos(th)^4 + 2*(2*a^2*m*r^4 - a^2*r^5 - (a^4 + a^2*q^2)*r^3)*cos(th)^2)*sin(th)^2)/(a^6*cos(th)^6 + 3*a^4*r^2*cos(th)^4 + 3*a^2*r^4*cos(th)^2 + r^6) \n", "Gam^th_t,t = (a^2*q^2 - 2*a^2*m*r)*cos(th)*sin(th)/(a^6*cos(th)^6 + 3*a^4*r^2*cos(th)^4 + 3*a^2*r^4*cos(th)^2 + r^6) \n", "Gam^th_t,ph = -(a^3*q^2 - 2*a^3*m*r + a*q^2*r^2 - 2*a*m*r^3)*cos(th)*sin(th)/(a^6*cos(th)^6 + 3*a^4*r^2*cos(th)^4 + 3*a^2*r^4*cos(th)^2 + r^6) \n", "Gam^th_r,r = -a^2*cos(th)*sin(th)/(2*m*r^3 - r^4 - (a^2 + q^2)*r^2 - (a^4 + a^2*q^2 - 2*a^2*m*r + a^2*r^2)*cos(th)^2) \n", "Gam^th_r,th = r/(a^2*cos(th)^2 + r^2) \n", "Gam^th_th,th = -a^2*cos(th)*sin(th)/(a^2*cos(th)^2 + r^2) \n", "Gam^th_ph,ph = -((a^6 + a^4*q^2 - 2*a^4*m*r + a^4*r^2)*cos(th)^5 - 2*(2*a^2*m*r^3 - a^2*r^4 - (a^4 + a^2*q^2)*r^2)*cos(th)^3 - (a^4*q^2 - 2*a^4*m*r + 2*a^2*q^2*r^2 - 4*a^2*m*r^3 - a^2*r^4 - r^6)*cos(th))*sin(th)/(a^6*cos(th)^6 + 3*a^4*r^2*cos(th)^4 + 3*a^2*r^4*cos(th)^2 + r^6) \n", "Gam^ph_t,r = (a^3*m*cos(th)^2 + a*q^2*r - a*m*r^2)/(2*m*r^5 - r^6 - (a^2 + q^2)*r^4 - (a^6 + a^4*q^2 - 2*a^4*m*r + a^4*r^2)*cos(th)^4 + 2*(2*a^2*m*r^3 - a^2*r^4 - (a^4 + a^2*q^2)*r^2)*cos(th)^2) \n", "Gam^ph_t,th = (a*q^2 - 2*a*m*r)*cos(th)/((a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4)*sin(th)) \n", "Gam^ph_r,ph = -(a^2*q^2*r - a^2*m*r^2 + q^2*r^3 - 2*m*r^4 + r^5 - (a^4*m - a^4*r)*cos(th)^4 + (a^4*m - a^2*m*r^2 + 2*a^2*r^3)*cos(th)^2)/(2*m*r^5 - r^6 - (a^2 + q^2)*r^4 - (a^6 + a^4*q^2 - 2*a^4*m*r + a^4*r^2)*cos(th)^4 + 2*(2*a^2*m*r^3 - a^2*r^4 - (a^4 + a^2*q^2)*r^2)*cos(th)^2) \n", "Gam^ph_th,ph = (a^4*cos(th)*sin(th)^4 - (2*a^4 + a^2*q^2 - 2*a^2*m*r + 2*a^2*r^2)*cos(th)*sin(th)^2 + (a^4 + 2*a^2*r^2 + r^4)*cos(th))/((a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4)*sin(th)) " ] }, "execution_count": 27, "metadata": {}, "output_type": "execute_result" } ], "source": [ "g.christoffel_symbols_display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Killing vectors

\n", "

The default vector frame on the spacetime manifold is the coordinate basis associated with Boyer-Lindquist coordinates:

" ] }, { "cell_type": "code", "execution_count": 28, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\mathrm{True}\\)" ], "text/latex": [ "$\\displaystyle \\mathrm{True}$" ], "text/plain": [ "True" ] }, "execution_count": 28, "metadata": {}, "output_type": "execute_result" } ], "source": [ "M.default_frame() is BL.frame()" ] }, { "cell_type": "code", "execution_count": 29, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\left(\\mathcal{M}, \\left(\\frac{\\partial}{\\partial t },\\frac{\\partial}{\\partial r },\\frac{\\partial}{\\partial {\\theta} },\\frac{\\partial}{\\partial {\\phi} }\\right)\\right)\\)" ], "text/latex": [ "$\\displaystyle \\left(\\mathcal{M}, \\left(\\frac{\\partial}{\\partial t },\\frac{\\partial}{\\partial r },\\frac{\\partial}{\\partial {\\theta} },\\frac{\\partial}{\\partial {\\phi} }\\right)\\right)$" ], "text/plain": [ "Coordinate frame (M, (∂/∂t,∂/∂r,∂/∂th,∂/∂ph))" ] }, "execution_count": 29, "metadata": {}, "output_type": "execute_result" } ], "source": [ "BL.frame()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Let us consider the first vector field of this frame:

" ] }, { "cell_type": "code", "execution_count": 30, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\frac{\\partial}{\\partial t }\\)" ], "text/latex": [ "$\\displaystyle \\frac{\\partial}{\\partial t }$" ], "text/plain": [ "Vector field ∂/∂t on the 4-dimensional Lorentzian manifold M" ] }, "execution_count": 30, "metadata": {}, "output_type": "execute_result" } ], "source": [ "xi = BL.frame()[0] ; xi" ] }, { "cell_type": "code", "execution_count": 31, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Vector field ∂/∂t on the 4-dimensional Lorentzian manifold M\n" ] } ], "source": [ "print(xi)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

The 1-form associated to it by metric duality is

" ] }, { "cell_type": "code", "execution_count": 32, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\left( -\\frac{q^{2} - 2 \\, m r}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} - 1 \\right) \\mathrm{d} t + \\left( \\frac{{\\left(q^{2} - 2 \\, m r\\right)} a \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\right) \\mathrm{d} {\\phi}\\)" ], "text/latex": [ "$\\displaystyle \\left( -\\frac{q^{2} - 2 \\, m r}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} - 1 \\right) \\mathrm{d} t + \\left( \\frac{{\\left(q^{2} - 2 \\, m r\\right)} a \\sin\\left({\\theta}\\right)^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\right) \\mathrm{d} {\\phi}$" ], "text/plain": [ "(-(q^2 - 2*m*r)/(a^2*cos(th)^2 + r^2) - 1) dt + (q^2 - 2*m*r)*a*sin(th)^2/(a^2*cos(th)^2 + r^2) dph" ] }, "execution_count": 32, "metadata": {}, "output_type": "execute_result" } ], "source": [ "xi_form = xi.down(g)\n", "xi_form.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Its covariant derivative is

" ] }, { "cell_type": "code", "execution_count": 33, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Tensor field of type (0,2) on the 4-dimensional Lorentzian manifold M\n" ] }, { "data": { "text/html": [ "\\(\\displaystyle \\left( \\frac{m r^{4} - {\\left(2 \\, m^{2} + q^{2}\\right)} r^{3} + {\\left(a^{2} m + 3 \\, m q^{2}\\right)} r^{2} - {\\left(a^{4} m + a^{2} m q^{2} - 2 \\, a^{2} m^{2} r + a^{2} m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2} - {\\left(a^{2} q^{2} + q^{4}\\right)} r}{2 \\, m r^{5} - r^{6} - {\\left(a^{2} + q^{2}\\right)} r^{4} - {\\left(a^{6} + a^{4} q^{2} - 2 \\, a^{4} m r + a^{4} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + 2 \\, {\\left(2 \\, a^{2} m r^{3} - a^{2} r^{4} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\right) \\mathrm{d} t\\otimes \\mathrm{d} r + \\left( -\\frac{{\\left(a^{2} q^{2} - 2 \\, a^{2} m r\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} t\\otimes \\mathrm{d} {\\theta} + \\left( -\\frac{m r^{4} - {\\left(2 \\, m^{2} + q^{2}\\right)} r^{3} + {\\left(a^{2} m + 3 \\, m q^{2}\\right)} r^{2} - {\\left(a^{4} m + a^{2} m q^{2} - 2 \\, a^{2} m^{2} r + a^{2} m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2} - {\\left(a^{2} q^{2} + q^{4}\\right)} r}{2 \\, m r^{5} - r^{6} - {\\left(a^{2} + q^{2}\\right)} r^{4} - {\\left(a^{6} + a^{4} q^{2} - 2 \\, a^{4} m r + a^{4} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + 2 \\, {\\left(2 \\, a^{2} m r^{3} - a^{2} r^{4} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\right) \\mathrm{d} r\\otimes \\mathrm{d} t + \\left( -\\frac{a^{3} m \\cos\\left({\\theta}\\right)^{4} - a q^{2} r + a m r^{2} - {\\left(a^{3} m - a q^{2} r + a m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} r\\otimes \\mathrm{d} {\\phi} + \\left( \\frac{{\\left(a^{2} q^{2} - 2 \\, a^{2} m r\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} {\\theta}\\otimes \\mathrm{d} t + \\left( -\\frac{{\\left(a^{3} q^{2} - 2 \\, a^{3} m r + a q^{2} r^{2} - 2 \\, a m r^{3}\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} {\\theta}\\otimes \\mathrm{d} {\\phi} + \\left( \\frac{a^{3} m \\cos\\left({\\theta}\\right)^{4} - a q^{2} r + a m r^{2} - {\\left(a^{3} m - a q^{2} r + a m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} {\\phi}\\otimes \\mathrm{d} r + \\left( \\frac{{\\left(a^{3} q^{2} - 2 \\, a^{3} m r + a q^{2} r^{2} - 2 \\, a m r^{3}\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} {\\phi}\\otimes \\mathrm{d} {\\theta}\\)" ], "text/latex": [ "$\\displaystyle \\left( \\frac{m r^{4} - {\\left(2 \\, m^{2} + q^{2}\\right)} r^{3} + {\\left(a^{2} m + 3 \\, m q^{2}\\right)} r^{2} - {\\left(a^{4} m + a^{2} m q^{2} - 2 \\, a^{2} m^{2} r + a^{2} m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2} - {\\left(a^{2} q^{2} + q^{4}\\right)} r}{2 \\, m r^{5} - r^{6} - {\\left(a^{2} + q^{2}\\right)} r^{4} - {\\left(a^{6} + a^{4} q^{2} - 2 \\, a^{4} m r + a^{4} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + 2 \\, {\\left(2 \\, a^{2} m r^{3} - a^{2} r^{4} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\right) \\mathrm{d} t\\otimes \\mathrm{d} r + \\left( -\\frac{{\\left(a^{2} q^{2} - 2 \\, a^{2} m r\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} t\\otimes \\mathrm{d} {\\theta} + \\left( -\\frac{m r^{4} - {\\left(2 \\, m^{2} + q^{2}\\right)} r^{3} + {\\left(a^{2} m + 3 \\, m q^{2}\\right)} r^{2} - {\\left(a^{4} m + a^{2} m q^{2} - 2 \\, a^{2} m^{2} r + a^{2} m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2} - {\\left(a^{2} q^{2} + q^{4}\\right)} r}{2 \\, m r^{5} - r^{6} - {\\left(a^{2} + q^{2}\\right)} r^{4} - {\\left(a^{6} + a^{4} q^{2} - 2 \\, a^{4} m r + a^{4} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + 2 \\, {\\left(2 \\, a^{2} m r^{3} - a^{2} r^{4} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\right) \\mathrm{d} r\\otimes \\mathrm{d} t + \\left( -\\frac{a^{3} m \\cos\\left({\\theta}\\right)^{4} - a q^{2} r + a m r^{2} - {\\left(a^{3} m - a q^{2} r + a m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} r\\otimes \\mathrm{d} {\\phi} + \\left( \\frac{{\\left(a^{2} q^{2} - 2 \\, a^{2} m r\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} {\\theta}\\otimes \\mathrm{d} t + \\left( -\\frac{{\\left(a^{3} q^{2} - 2 \\, a^{3} m r + a q^{2} r^{2} - 2 \\, a m r^{3}\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} {\\theta}\\otimes \\mathrm{d} {\\phi} + \\left( \\frac{a^{3} m \\cos\\left({\\theta}\\right)^{4} - a q^{2} r + a m r^{2} - {\\left(a^{3} m - a q^{2} r + a m r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} {\\phi}\\otimes \\mathrm{d} r + \\left( \\frac{{\\left(a^{3} q^{2} - 2 \\, a^{3} m r + a q^{2} r^{2} - 2 \\, a m r^{3}\\right)} \\cos\\left({\\theta}\\right) \\sin\\left({\\theta}\\right)}{a^{4} \\cos\\left({\\theta}\\right)^{4} + 2 \\, a^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + r^{4}} \\right) \\mathrm{d} {\\phi}\\otimes \\mathrm{d} {\\theta}$" ], "text/plain": [ "(m*r^4 - (2*m^2 + q^2)*r^3 + (a^2*m + 3*m*q^2)*r^2 - (a^4*m + a^2*m*q^2 - 2*a^2*m^2*r + a^2*m*r^2)*cos(th)^2 - (a^2*q^2 + q^4)*r)/(2*m*r^5 - r^6 - (a^2 + q^2)*r^4 - (a^6 + a^4*q^2 - 2*a^4*m*r + a^4*r^2)*cos(th)^4 + 2*(2*a^2*m*r^3 - a^2*r^4 - (a^4 + a^2*q^2)*r^2)*cos(th)^2) dt⊗dr - (a^2*q^2 - 2*a^2*m*r)*cos(th)*sin(th)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) dt⊗dth - (m*r^4 - (2*m^2 + q^2)*r^3 + (a^2*m + 3*m*q^2)*r^2 - (a^4*m + a^2*m*q^2 - 2*a^2*m^2*r + a^2*m*r^2)*cos(th)^2 - (a^2*q^2 + q^4)*r)/(2*m*r^5 - r^6 - (a^2 + q^2)*r^4 - (a^6 + a^4*q^2 - 2*a^4*m*r + a^4*r^2)*cos(th)^4 + 2*(2*a^2*m*r^3 - a^2*r^4 - (a^4 + a^2*q^2)*r^2)*cos(th)^2) dr⊗dt - (a^3*m*cos(th)^4 - a*q^2*r + a*m*r^2 - (a^3*m - a*q^2*r + a*m*r^2)*cos(th)^2)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) dr⊗dph + (a^2*q^2 - 2*a^2*m*r)*cos(th)*sin(th)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) dth⊗dt - (a^3*q^2 - 2*a^3*m*r + a*q^2*r^2 - 2*a*m*r^3)*cos(th)*sin(th)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) dth⊗dph + (a^3*m*cos(th)^4 - a*q^2*r + a*m*r^2 - (a^3*m - a*q^2*r + a*m*r^2)*cos(th)^2)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) dph⊗dr + (a^3*q^2 - 2*a^3*m*r + a*q^2*r^2 - 2*a*m*r^3)*cos(th)*sin(th)/(a^4*cos(th)^4 + 2*a^2*r^2*cos(th)^2 + r^4) dph⊗dth" ] }, "execution_count": 33, "metadata": {}, "output_type": "execute_result" } ], "source": [ "nab_xi = nabla(xi_form)\n", "print(nab_xi)\n", "nab_xi.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Let us check that the vector field $\\xi=\\frac{\\partial}{\\partial t}$ obeys Killing equation:

" ] }, { "cell_type": "code", "execution_count": 34, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\mathrm{True}\\)" ], "text/latex": [ "$\\displaystyle \\mathrm{True}$" ], "text/plain": [ "True" ] }, "execution_count": 34, "metadata": {}, "output_type": "execute_result" } ], "source": [ "nab_xi.symmetrize() == 0" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Similarly, let us check that $\\chi := \\frac{\\partial}{\\partial\\phi}$ is a Killing vector:

" ] }, { "cell_type": "code", "execution_count": 35, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\frac{\\partial}{\\partial {\\phi} }\\)" ], "text/latex": [ "$\\displaystyle \\frac{\\partial}{\\partial {\\phi} }$" ], "text/plain": [ "Vector field ∂/∂ph on the 4-dimensional Lorentzian manifold M" ] }, "execution_count": 35, "metadata": {}, "output_type": "execute_result" } ], "source": [ "chi = BL.frame()[3] ; chi" ] }, { "cell_type": "code", "execution_count": 36, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\mathrm{True}\\)" ], "text/latex": [ "$\\displaystyle \\mathrm{True}$" ], "text/plain": [ "True" ] }, "execution_count": 36, "metadata": {}, "output_type": "execute_result" } ], "source": [ "nabla(chi.down(g)).symmetrize() == 0" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Another way to check that $\\xi$ and $\\chi$ are Killing vectors is the vanishing of the Lie derivative of the metric tensor along them:

" ] }, { "cell_type": "code", "execution_count": 37, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\mathrm{True}\\)" ], "text/latex": [ "$\\displaystyle \\mathrm{True}$" ], "text/plain": [ "True" ] }, "execution_count": 37, "metadata": {}, "output_type": "execute_result" } ], "source": [ "g.lie_derivative(xi) == 0" ] }, { "cell_type": "code", "execution_count": 38, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\mathrm{True}\\)" ], "text/latex": [ "$\\displaystyle \\mathrm{True}$" ], "text/plain": [ "True" ] }, "execution_count": 38, "metadata": {}, "output_type": "execute_result" } ], "source": [ "g.lie_derivative(chi) == 0" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Curvature\n", "\n", "The Ricci tensor of $g$:" ] }, { "cell_type": "code", "execution_count": 39, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Field of symmetric bilinear forms Ric(g) on the 4-dimensional Lorentzian manifold M\n" ] } ], "source": [ "Ric = g.ricci()\n", "print(Ric)" ] }, { "cell_type": "code", "execution_count": 40, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\mathrm{Ric}\\left(g\\right) = \\left( -\\frac{a^{2} q^{2} \\cos\\left({\\theta}\\right)^{2} - 2 \\, a^{2} q^{2} - q^{4} + 2 \\, m q^{2} r - q^{2} r^{2}}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\right) \\mathrm{d} t\\otimes \\mathrm{d} t + \\left( \\frac{{\\left(2 \\, a^{5} q^{2} + a^{3} q^{4} - 2 \\, a^{3} m q^{2} r + 2 \\, a^{3} q^{2} r^{2}\\right)} \\sin\\left({\\theta}\\right)^{4} - {\\left(2 \\, a^{5} q^{2} + a^{3} q^{4} - 2 \\, a^{3} m q^{2} r - 2 \\, a m q^{2} r^{3} + 2 \\, a q^{2} r^{4} + {\\left(4 \\, a^{3} q^{2} + a q^{4}\\right)} r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}{a^{8} \\cos\\left({\\theta}\\right)^{8} + 4 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{6} + 6 \\, a^{4} r^{4} \\cos\\left({\\theta}\\right)^{4} + 4 \\, a^{2} r^{6} \\cos\\left({\\theta}\\right)^{2} + r^{8}} \\right) \\mathrm{d} t\\otimes \\mathrm{d} {\\phi} + \\left( \\frac{q^{2}}{2 \\, m r^{3} - r^{4} - {\\left(a^{2} + q^{2}\\right)} r^{2} - {\\left(a^{4} + a^{2} q^{2} - 2 \\, a^{2} m r + a^{2} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\right) \\mathrm{d} r\\otimes \\mathrm{d} r + \\left( \\frac{q^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\right) \\mathrm{d} {\\theta}\\otimes \\mathrm{d} {\\theta} + \\left( \\frac{{\\left(2 \\, a^{5} q^{2} + a^{3} q^{4} - 2 \\, a^{3} m q^{2} r + 2 \\, a^{3} q^{2} r^{2}\\right)} \\sin\\left({\\theta}\\right)^{4} - {\\left(2 \\, a^{5} q^{2} + a^{3} q^{4} - 2 \\, a^{3} m q^{2} r - 2 \\, a m q^{2} r^{3} + 2 \\, a q^{2} r^{4} + {\\left(4 \\, a^{3} q^{2} + a q^{4}\\right)} r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}{a^{8} \\cos\\left({\\theta}\\right)^{8} + 4 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{6} + 6 \\, a^{4} r^{4} \\cos\\left({\\theta}\\right)^{4} + 4 \\, a^{2} r^{6} \\cos\\left({\\theta}\\right)^{2} + r^{8}} \\right) \\mathrm{d} {\\phi}\\otimes \\mathrm{d} t + \\left( -\\frac{{\\left(a^{6} q^{2} + a^{4} q^{4} - 2 \\, a^{4} m q^{2} r + a^{4} q^{2} r^{2}\\right)} \\sin\\left({\\theta}\\right)^{6} - {\\left(a^{4} q^{4} - 2 \\, a^{4} m q^{2} r + a^{2} q^{4} r^{2} - 2 \\, a^{2} m q^{2} r^{3}\\right)} \\sin\\left({\\theta}\\right)^{4} - {\\left(a^{6} q^{2} + 3 \\, a^{4} q^{2} r^{2} + 3 \\, a^{2} q^{2} r^{4} + q^{2} r^{6}\\right)} \\sin\\left({\\theta}\\right)^{2}}{a^{8} \\cos\\left({\\theta}\\right)^{8} + 4 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{6} + 6 \\, a^{4} r^{4} \\cos\\left({\\theta}\\right)^{4} + 4 \\, a^{2} r^{6} \\cos\\left({\\theta}\\right)^{2} + r^{8}} \\right) \\mathrm{d} {\\phi}\\otimes \\mathrm{d} {\\phi}\\)" ], "text/latex": [ "$\\displaystyle \\mathrm{Ric}\\left(g\\right) = \\left( -\\frac{a^{2} q^{2} \\cos\\left({\\theta}\\right)^{2} - 2 \\, a^{2} q^{2} - q^{4} + 2 \\, m q^{2} r - q^{2} r^{2}}{a^{6} \\cos\\left({\\theta}\\right)^{6} + 3 \\, a^{4} r^{2} \\cos\\left({\\theta}\\right)^{4} + 3 \\, a^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} + r^{6}} \\right) \\mathrm{d} t\\otimes \\mathrm{d} t + \\left( \\frac{{\\left(2 \\, a^{5} q^{2} + a^{3} q^{4} - 2 \\, a^{3} m q^{2} r + 2 \\, a^{3} q^{2} r^{2}\\right)} \\sin\\left({\\theta}\\right)^{4} - {\\left(2 \\, a^{5} q^{2} + a^{3} q^{4} - 2 \\, a^{3} m q^{2} r - 2 \\, a m q^{2} r^{3} + 2 \\, a q^{2} r^{4} + {\\left(4 \\, a^{3} q^{2} + a q^{4}\\right)} r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}{a^{8} \\cos\\left({\\theta}\\right)^{8} + 4 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{6} + 6 \\, a^{4} r^{4} \\cos\\left({\\theta}\\right)^{4} + 4 \\, a^{2} r^{6} \\cos\\left({\\theta}\\right)^{2} + r^{8}} \\right) \\mathrm{d} t\\otimes \\mathrm{d} {\\phi} + \\left( \\frac{q^{2}}{2 \\, m r^{3} - r^{4} - {\\left(a^{2} + q^{2}\\right)} r^{2} - {\\left(a^{4} + a^{2} q^{2} - 2 \\, a^{2} m r + a^{2} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\right) \\mathrm{d} r\\otimes \\mathrm{d} r + \\left( \\frac{q^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\right) \\mathrm{d} {\\theta}\\otimes \\mathrm{d} {\\theta} + \\left( \\frac{{\\left(2 \\, a^{5} q^{2} + a^{3} q^{4} - 2 \\, a^{3} m q^{2} r + 2 \\, a^{3} q^{2} r^{2}\\right)} \\sin\\left({\\theta}\\right)^{4} - {\\left(2 \\, a^{5} q^{2} + a^{3} q^{4} - 2 \\, a^{3} m q^{2} r - 2 \\, a m q^{2} r^{3} + 2 \\, a q^{2} r^{4} + {\\left(4 \\, a^{3} q^{2} + a q^{4}\\right)} r^{2}\\right)} \\sin\\left({\\theta}\\right)^{2}}{a^{8} \\cos\\left({\\theta}\\right)^{8} + 4 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{6} + 6 \\, a^{4} r^{4} \\cos\\left({\\theta}\\right)^{4} + 4 \\, a^{2} r^{6} \\cos\\left({\\theta}\\right)^{2} + r^{8}} \\right) \\mathrm{d} {\\phi}\\otimes \\mathrm{d} t + \\left( -\\frac{{\\left(a^{6} q^{2} + a^{4} q^{4} - 2 \\, a^{4} m q^{2} r + a^{4} q^{2} r^{2}\\right)} \\sin\\left({\\theta}\\right)^{6} - {\\left(a^{4} q^{4} - 2 \\, a^{4} m q^{2} r + a^{2} q^{4} r^{2} - 2 \\, a^{2} m q^{2} r^{3}\\right)} \\sin\\left({\\theta}\\right)^{4} - {\\left(a^{6} q^{2} + 3 \\, a^{4} q^{2} r^{2} + 3 \\, a^{2} q^{2} r^{4} + q^{2} r^{6}\\right)} \\sin\\left({\\theta}\\right)^{2}}{a^{8} \\cos\\left({\\theta}\\right)^{8} + 4 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{6} + 6 \\, a^{4} r^{4} \\cos\\left({\\theta}\\right)^{4} + 4 \\, a^{2} r^{6} \\cos\\left({\\theta}\\right)^{2} + r^{8}} \\right) \\mathrm{d} {\\phi}\\otimes \\mathrm{d} {\\phi}$" ], "text/plain": [ "Ric(g) = -(a^2*q^2*cos(th)^2 - 2*a^2*q^2 - q^4 + 2*m*q^2*r - q^2*r^2)/(a^6*cos(th)^6 + 3*a^4*r^2*cos(th)^4 + 3*a^2*r^4*cos(th)^2 + r^6) dt⊗dt + ((2*a^5*q^2 + a^3*q^4 - 2*a^3*m*q^2*r + 2*a^3*q^2*r^2)*sin(th)^4 - (2*a^5*q^2 + a^3*q^4 - 2*a^3*m*q^2*r - 2*a*m*q^2*r^3 + 2*a*q^2*r^4 + (4*a^3*q^2 + a*q^4)*r^2)*sin(th)^2)/(a^8*cos(th)^8 + 4*a^6*r^2*cos(th)^6 + 6*a^4*r^4*cos(th)^4 + 4*a^2*r^6*cos(th)^2 + r^8) dt⊗dph + q^2/(2*m*r^3 - r^4 - (a^2 + q^2)*r^2 - (a^4 + a^2*q^2 - 2*a^2*m*r + a^2*r^2)*cos(th)^2) dr⊗dr + q^2/(a^2*cos(th)^2 + r^2) dth⊗dth + ((2*a^5*q^2 + a^3*q^4 - 2*a^3*m*q^2*r + 2*a^3*q^2*r^2)*sin(th)^4 - (2*a^5*q^2 + a^3*q^4 - 2*a^3*m*q^2*r - 2*a*m*q^2*r^3 + 2*a*q^2*r^4 + (4*a^3*q^2 + a*q^4)*r^2)*sin(th)^2)/(a^8*cos(th)^8 + 4*a^6*r^2*cos(th)^6 + 6*a^4*r^4*cos(th)^4 + 4*a^2*r^6*cos(th)^2 + r^8) dph⊗dt - ((a^6*q^2 + a^4*q^4 - 2*a^4*m*q^2*r + a^4*q^2*r^2)*sin(th)^6 - (a^4*q^4 - 2*a^4*m*q^2*r + a^2*q^4*r^2 - 2*a^2*m*q^2*r^3)*sin(th)^4 - (a^6*q^2 + 3*a^4*q^2*r^2 + 3*a^2*q^2*r^4 + q^2*r^6)*sin(th)^2)/(a^8*cos(th)^8 + 4*a^6*r^2*cos(th)^6 + 6*a^4*r^4*cos(th)^4 + 4*a^2*r^6*cos(th)^2 + r^8) dph⊗dph" ] }, "execution_count": 40, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Ric.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We can get a shorter expression by factorizing the components:" ] }, { "cell_type": "code", "execution_count": 41, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\mathrm{Ric}\\left(g\\right) = -\\frac{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} - 2 \\, a^{2} - q^{2} + 2 \\, m r - r^{2}\\right)} q^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{3}} \\mathrm{d} t\\otimes \\mathrm{d} t + \\frac{{\\left(a^{2} \\sin\\left({\\theta}\\right)^{2} - a^{2} - r^{2}\\right)} {\\left(2 \\, a^{2} + q^{2} - 2 \\, m r + 2 \\, r^{2}\\right)} a q^{2} \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{4}} \\mathrm{d} t\\otimes \\mathrm{d} {\\phi} -\\frac{q^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)} {\\left(a^{2} + q^{2} - 2 \\, m r + r^{2}\\right)}} \\mathrm{d} r\\otimes \\mathrm{d} r + \\left( \\frac{q^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\right) \\mathrm{d} {\\theta}\\otimes \\mathrm{d} {\\theta} + \\frac{{\\left(a^{2} \\sin\\left({\\theta}\\right)^{2} - a^{2} - r^{2}\\right)} {\\left(2 \\, a^{2} + q^{2} - 2 \\, m r + 2 \\, r^{2}\\right)} a q^{2} \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{4}} \\mathrm{d} {\\phi}\\otimes \\mathrm{d} t -\\frac{{\\left(a^{4} \\sin\\left({\\theta}\\right)^{2} + a^{2} q^{2} \\sin\\left({\\theta}\\right)^{2} - 2 \\, a^{2} m r \\sin\\left({\\theta}\\right)^{2} + a^{2} r^{2} \\sin\\left({\\theta}\\right)^{2} + a^{4} + 2 \\, a^{2} r^{2} + r^{4}\\right)} {\\left(a^{2} \\sin\\left({\\theta}\\right)^{2} - a^{2} - r^{2}\\right)} q^{2} \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{4}} \\mathrm{d} {\\phi}\\otimes \\mathrm{d} {\\phi}\\)" ], "text/latex": [ "$\\displaystyle \\mathrm{Ric}\\left(g\\right) = -\\frac{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} - 2 \\, a^{2} - q^{2} + 2 \\, m r - r^{2}\\right)} q^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{3}} \\mathrm{d} t\\otimes \\mathrm{d} t + \\frac{{\\left(a^{2} \\sin\\left({\\theta}\\right)^{2} - a^{2} - r^{2}\\right)} {\\left(2 \\, a^{2} + q^{2} - 2 \\, m r + 2 \\, r^{2}\\right)} a q^{2} \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{4}} \\mathrm{d} t\\otimes \\mathrm{d} {\\phi} -\\frac{q^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)} {\\left(a^{2} + q^{2} - 2 \\, m r + r^{2}\\right)}} \\mathrm{d} r\\otimes \\mathrm{d} r + \\left( \\frac{q^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} \\right) \\mathrm{d} {\\theta}\\otimes \\mathrm{d} {\\theta} + \\frac{{\\left(a^{2} \\sin\\left({\\theta}\\right)^{2} - a^{2} - r^{2}\\right)} {\\left(2 \\, a^{2} + q^{2} - 2 \\, m r + 2 \\, r^{2}\\right)} a q^{2} \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{4}} \\mathrm{d} {\\phi}\\otimes \\mathrm{d} t -\\frac{{\\left(a^{4} \\sin\\left({\\theta}\\right)^{2} + a^{2} q^{2} \\sin\\left({\\theta}\\right)^{2} - 2 \\, a^{2} m r \\sin\\left({\\theta}\\right)^{2} + a^{2} r^{2} \\sin\\left({\\theta}\\right)^{2} + a^{4} + 2 \\, a^{2} r^{2} + r^{4}\\right)} {\\left(a^{2} \\sin\\left({\\theta}\\right)^{2} - a^{2} - r^{2}\\right)} q^{2} \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{4}} \\mathrm{d} {\\phi}\\otimes \\mathrm{d} {\\phi}$" ], "text/plain": [ "Ric(g) = -(a^2*cos(th)^2 - 2*a^2 - q^2 + 2*m*r - r^2)*q^2/(a^2*cos(th)^2 + r^2)^3 dt⊗dt + (a^2*sin(th)^2 - a^2 - r^2)*(2*a^2 + q^2 - 2*m*r + 2*r^2)*a*q^2*sin(th)^2/(a^2*cos(th)^2 + r^2)^4 dt⊗dph - q^2/((a^2*cos(th)^2 + r^2)*(a^2 + q^2 - 2*m*r + r^2)) dr⊗dr + q^2/(a^2*cos(th)^2 + r^2) dth⊗dth + (a^2*sin(th)^2 - a^2 - r^2)*(2*a^2 + q^2 - 2*m*r + 2*r^2)*a*q^2*sin(th)^2/(a^2*cos(th)^2 + r^2)^4 dph⊗dt - (a^4*sin(th)^2 + a^2*q^2*sin(th)^2 - 2*a^2*m*r*sin(th)^2 + a^2*r^2*sin(th)^2 + a^4 + 2*a^2*r^2 + r^4)*(a^2*sin(th)^2 - a^2 - r^2)*q^2*sin(th)^2/(a^2*cos(th)^2 + r^2)^4 dph⊗dph" ] }, "execution_count": 41, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Ric.apply_map(factor)\n", "Ric.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "A matrix view of the components:" ] }, { "cell_type": "code", "execution_count": 42, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\left(\\begin{array}{rrrr}\n", "-\\frac{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} - 2 \\, a^{2} - q^{2} + 2 \\, m r - r^{2}\\right)} q^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{3}} & 0 & 0 & \\frac{{\\left(a^{2} \\sin\\left({\\theta}\\right)^{2} - a^{2} - r^{2}\\right)} {\\left(2 \\, a^{2} + q^{2} - 2 \\, m r + 2 \\, r^{2}\\right)} a q^{2} \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{4}} \\\\\n", "0 & -\\frac{q^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)} {\\left(a^{2} + q^{2} - 2 \\, m r + r^{2}\\right)}} & 0 & 0 \\\\\n", "0 & 0 & \\frac{q^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} & 0 \\\\\n", "\\frac{{\\left(a^{2} \\sin\\left({\\theta}\\right)^{2} - a^{2} - r^{2}\\right)} {\\left(2 \\, a^{2} + q^{2} - 2 \\, m r + 2 \\, r^{2}\\right)} a q^{2} \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{4}} & 0 & 0 & -\\frac{{\\left(a^{4} \\sin\\left({\\theta}\\right)^{2} + a^{2} q^{2} \\sin\\left({\\theta}\\right)^{2} - 2 \\, a^{2} m r \\sin\\left({\\theta}\\right)^{2} + a^{2} r^{2} \\sin\\left({\\theta}\\right)^{2} + a^{4} + 2 \\, a^{2} r^{2} + r^{4}\\right)} {\\left(a^{2} \\sin\\left({\\theta}\\right)^{2} - a^{2} - r^{2}\\right)} q^{2} \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{4}}\n", "\\end{array}\\right)\\)" ], "text/latex": [ "$\\displaystyle \\left(\\begin{array}{rrrr}\n", "-\\frac{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} - 2 \\, a^{2} - q^{2} + 2 \\, m r - r^{2}\\right)} q^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{3}} & 0 & 0 & \\frac{{\\left(a^{2} \\sin\\left({\\theta}\\right)^{2} - a^{2} - r^{2}\\right)} {\\left(2 \\, a^{2} + q^{2} - 2 \\, m r + 2 \\, r^{2}\\right)} a q^{2} \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{4}} \\\\\n", "0 & -\\frac{q^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)} {\\left(a^{2} + q^{2} - 2 \\, m r + r^{2}\\right)}} & 0 & 0 \\\\\n", "0 & 0 & \\frac{q^{2}}{a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}} & 0 \\\\\n", "\\frac{{\\left(a^{2} \\sin\\left({\\theta}\\right)^{2} - a^{2} - r^{2}\\right)} {\\left(2 \\, a^{2} + q^{2} - 2 \\, m r + 2 \\, r^{2}\\right)} a q^{2} \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{4}} & 0 & 0 & -\\frac{{\\left(a^{4} \\sin\\left({\\theta}\\right)^{2} + a^{2} q^{2} \\sin\\left({\\theta}\\right)^{2} - 2 \\, a^{2} m r \\sin\\left({\\theta}\\right)^{2} + a^{2} r^{2} \\sin\\left({\\theta}\\right)^{2} + a^{4} + 2 \\, a^{2} r^{2} + r^{4}\\right)} {\\left(a^{2} \\sin\\left({\\theta}\\right)^{2} - a^{2} - r^{2}\\right)} q^{2} \\sin\\left({\\theta}\\right)^{2}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{4}}\n", "\\end{array}\\right)$" ], "text/plain": [ "[ -(a^2*cos(th)^2 - 2*a^2 - q^2 + 2*m*r - r^2)*q^2/(a^2*cos(th)^2 + r^2)^3 0 0 (a^2*sin(th)^2 - a^2 - r^2)*(2*a^2 + q^2 - 2*m*r + 2*r^2)*a*q^2*sin(th)^2/(a^2*cos(th)^2 + r^2)^4]\n", "[ 0 -q^2/((a^2*cos(th)^2 + r^2)*(a^2 + q^2 - 2*m*r + r^2)) 0 0]\n", "[ 0 0 q^2/(a^2*cos(th)^2 + r^2) 0]\n", "[ (a^2*sin(th)^2 - a^2 - r^2)*(2*a^2 + q^2 - 2*m*r + 2*r^2)*a*q^2*sin(th)^2/(a^2*cos(th)^2 + r^2)^4 0 0 -(a^4*sin(th)^2 + a^2*q^2*sin(th)^2 - 2*a^2*m*r*sin(th)^2 + a^2*r^2*sin(th)^2 + a^4 + 2*a^2*r^2 + r^4)*(a^2*sin(th)^2 - a^2 - r^2)*q^2*sin(th)^2/(a^2*cos(th)^2 + r^2)^4]" ] }, "execution_count": 42, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Ric[:]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Let us check that in the Kerr case, i.e. when $q=0$, the Ricci tensor is zero:

" ] }, { "cell_type": "code", "execution_count": 43, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\left(\\begin{array}{rrrr}\n", "0 & 0 & 0 & 0 \\\\\n", "0 & 0 & 0 & 0 \\\\\n", "0 & 0 & 0 & 0 \\\\\n", "0 & 0 & 0 & 0\n", "\\end{array}\\right)\\)" ], "text/latex": [ "$\\displaystyle \\left(\\begin{array}{rrrr}\n", "0 & 0 & 0 & 0 \\\\\n", "0 & 0 & 0 & 0 \\\\\n", "0 & 0 & 0 & 0 \\\\\n", "0 & 0 & 0 & 0\n", "\\end{array}\\right)$" ], "text/plain": [ "[0 0 0 0]\n", "[0 0 0 0]\n", "[0 0 0 0]\n", "[0 0 0 0]" ] }, "execution_count": 43, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Ric_Kerr = Ric.copy()\n", "Ric_Kerr.apply_map(lambda f: f.subs({q: 0}))\n", "Ric_Kerr[:]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Riemann tensor\n", "\n", "The Riemann curvature tensor of $g$:" ] }, { "cell_type": "code", "execution_count": 44, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Tensor field Riem(g) of type (1,3) on the 4-dimensional Lorentzian manifold M\n" ] } ], "source": [ "R = g.riemann()\n", "R.apply_map(factor)\n", "print(R)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

The component $R^0_{\\ \\, 101}$ of the Riemann tensor is

" ] }, { "cell_type": "code", "execution_count": 45, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle -\\frac{a^{4} q^{2} \\cos\\left({\\theta}\\right)^{4} - 3 \\, a^{4} m r \\cos\\left({\\theta}\\right)^{4} - 2 \\, a^{4} q^{2} \\cos\\left({\\theta}\\right)^{2} + 9 \\, a^{4} m r \\cos\\left({\\theta}\\right)^{2} - 2 \\, a^{2} q^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + 7 \\, a^{2} m r^{3} \\cos\\left({\\theta}\\right)^{2} + 4 \\, a^{2} q^{2} r^{2} - 3 \\, a^{2} m r^{3} + 3 \\, q^{2} r^{4} - 2 \\, m r^{5}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{3} {\\left(a^{2} + q^{2} - 2 \\, m r + r^{2}\\right)}}\\)" ], "text/latex": [ "$\\displaystyle -\\frac{a^{4} q^{2} \\cos\\left({\\theta}\\right)^{4} - 3 \\, a^{4} m r \\cos\\left({\\theta}\\right)^{4} - 2 \\, a^{4} q^{2} \\cos\\left({\\theta}\\right)^{2} + 9 \\, a^{4} m r \\cos\\left({\\theta}\\right)^{2} - 2 \\, a^{2} q^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + 7 \\, a^{2} m r^{3} \\cos\\left({\\theta}\\right)^{2} + 4 \\, a^{2} q^{2} r^{2} - 3 \\, a^{2} m r^{3} + 3 \\, q^{2} r^{4} - 2 \\, m r^{5}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{3} {\\left(a^{2} + q^{2} - 2 \\, m r + r^{2}\\right)}}$" ], "text/plain": [ "-(a^4*q^2*cos(th)^4 - 3*a^4*m*r*cos(th)^4 - 2*a^4*q^2*cos(th)^2 + 9*a^4*m*r*cos(th)^2 - 2*a^2*q^2*r^2*cos(th)^2 + 7*a^2*m*r^3*cos(th)^2 + 4*a^2*q^2*r^2 - 3*a^2*m*r^3 + 3*q^2*r^4 - 2*m*r^5)/((a^2*cos(th)^2 + r^2)^3*(a^2 + q^2 - 2*m*r + r^2))" ] }, "execution_count": 45, "metadata": {}, "output_type": "execute_result" } ], "source": [ "R[0,1,0,1]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

The expression in the uncharged limit (Kerr spacetime) is

" ] }, { "cell_type": "code", "execution_count": 46, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\frac{{\\left(3 \\, a^{2} \\cos\\left({\\theta}\\right)^{2} - r^{2}\\right)} {\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} - 3 \\, a^{2} - 2 \\, r^{2}\\right)} m r}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{3} {\\left(a^{2} - 2 \\, m r + r^{2}\\right)}}\\)" ], "text/latex": [ "$\\displaystyle \\frac{{\\left(3 \\, a^{2} \\cos\\left({\\theta}\\right)^{2} - r^{2}\\right)} {\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} - 3 \\, a^{2} - 2 \\, r^{2}\\right)} m r}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{3} {\\left(a^{2} - 2 \\, m r + r^{2}\\right)}}$" ], "text/plain": [ "(3*a^2*cos(th)^2 - r^2)*(a^2*cos(th)^2 - 3*a^2 - 2*r^2)*m*r/((a^2*cos(th)^2 + r^2)^3*(a^2 - 2*m*r + r^2))" ] }, "execution_count": 46, "metadata": {}, "output_type": "execute_result" } ], "source": [ "R[0,1,0,1].expr().subs(q=0).factor()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

while in the non-rotating limit (Reissner-Nordström spacetime), it is

" ] }, { "cell_type": "code", "execution_count": 47, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle -\\frac{3 \\, q^{2} - 2 \\, m r}{{\\left(q^{2} - 2 \\, m r + r^{2}\\right)} r^{2}}\\)" ], "text/latex": [ "$\\displaystyle -\\frac{3 \\, q^{2} - 2 \\, m r}{{\\left(q^{2} - 2 \\, m r + r^{2}\\right)} r^{2}}$" ], "text/plain": [ "-(3*q^2 - 2*m*r)/((q^2 - 2*m*r + r^2)*r^2)" ] }, "execution_count": 47, "metadata": {}, "output_type": "execute_result" } ], "source": [ "R[0,1,0,1].expr().subs(a=0).factor()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

In the Schwarzschild limit, it reduces to

" ] }, { "cell_type": "code", "execution_count": 48, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle -\\frac{2 \\, m}{{\\left(2 \\, m r - r^{2}\\right)} r}\\)" ], "text/latex": [ "$\\displaystyle -\\frac{2 \\, m}{{\\left(2 \\, m r - r^{2}\\right)} r}$" ], "text/plain": [ "-2*m/((2*m*r - r^2)*r)" ] }, "execution_count": 48, "metadata": {}, "output_type": "execute_result" } ], "source": [ "R[0,1,0,1].expr().subs(a=0, q=0)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Obviously, it vanishes in the flat space limit:

" ] }, { "cell_type": "code", "execution_count": 49, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle 0\\)" ], "text/latex": [ "$\\displaystyle 0$" ], "text/plain": [ "0" ] }, "execution_count": 49, "metadata": {}, "output_type": "execute_result" } ], "source": [ "R[0,1,0,1].expr().subs(m=0, a=0, q=0)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Ricci scalar\n", "\n", "The Ricci scalar $R = g^{ab} R_{ab}$ of the Kerr-Newman spacetime vanishes identically:" ] }, { "cell_type": "code", "execution_count": 50, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\begin{array}{llcl} \\mathrm{r}\\left(g\\right):& \\mathcal{M} & \\longrightarrow & \\mathbb{R} \\\\ & \\left(t, r, {\\theta}, {\\phi}\\right) & \\longmapsto & 0 \\end{array}\\)" ], "text/latex": [ "$\\displaystyle \\begin{array}{llcl} \\mathrm{r}\\left(g\\right):& \\mathcal{M} & \\longrightarrow & \\mathbb{R} \\\\ & \\left(t, r, {\\theta}, {\\phi}\\right) & \\longmapsto & 0 \\end{array}$" ], "text/plain": [ "r(g): M → ℝ\n", " (t, r, th, ph) ↦ 0" ] }, "execution_count": 50, "metadata": {}, "output_type": "execute_result" } ], "source": [ "g.ricci_scalar().display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Einstein equation

\n", "

The Einstein tensor is

" ] }, { "cell_type": "code", "execution_count": 51, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Field of symmetric bilinear forms Ric(g)-unnamed metric on the 4-dimensional Lorentzian manifold M\n" ] } ], "source": [ "G = Ric - 1/2*g.ricci_scalar()*g\n", "print(G)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Since the Ricci scalar is zero, the Einstein tensor reduces to the Ricci tensor:

" ] }, { "cell_type": "code", "execution_count": 52, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\mathrm{True}\\)" ], "text/latex": [ "$\\displaystyle \\mathrm{True}$" ], "text/plain": [ "True" ] }, "execution_count": 52, "metadata": {}, "output_type": "execute_result" } ], "source": [ "G == Ric" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The invariant $F_{ab} F^{ab}$ of the electromagnetic field:" ] }, { "cell_type": "code", "execution_count": 53, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Scalar field on the 4-dimensional Lorentzian manifold M\n" ] } ], "source": [ "Fuu = F.up(g)\n", "F2 = F['_ab']*Fuu['^ab']\n", "print(F2)" ] }, { "cell_type": "code", "execution_count": 54, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\begin{array}{llcl} & \\mathcal{M} & \\longrightarrow & \\mathbb{R} \\\\ & \\left(t, r, {\\theta}, {\\phi}\\right) & \\longmapsto & -\\frac{2 \\, {\\left(a^{4} q^{2} \\cos\\left({\\theta}\\right)^{4} - 6 \\, a^{2} q^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + q^{2} r^{4}\\right)}}{a^{8} \\cos\\left({\\theta}\\right)^{8} + 4 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{6} + 6 \\, a^{4} r^{4} \\cos\\left({\\theta}\\right)^{4} + 4 \\, a^{2} r^{6} \\cos\\left({\\theta}\\right)^{2} + r^{8}} \\end{array}\\)" ], "text/latex": [ "$\\displaystyle \\begin{array}{llcl} & \\mathcal{M} & \\longrightarrow & \\mathbb{R} \\\\ & \\left(t, r, {\\theta}, {\\phi}\\right) & \\longmapsto & -\\frac{2 \\, {\\left(a^{4} q^{2} \\cos\\left({\\theta}\\right)^{4} - 6 \\, a^{2} q^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} + q^{2} r^{4}\\right)}}{a^{8} \\cos\\left({\\theta}\\right)^{8} + 4 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{6} + 6 \\, a^{4} r^{4} \\cos\\left({\\theta}\\right)^{4} + 4 \\, a^{2} r^{6} \\cos\\left({\\theta}\\right)^{2} + r^{8}} \\end{array}$" ], "text/plain": [ "M → ℝ\n", "(t, r, th, ph) ↦ -2*(a^4*q^2*cos(th)^4 - 6*a^2*q^2*r^2*cos(th)^2 + q^2*r^4)/(a^8*cos(th)^8 + 4*a^6*r^2*cos(th)^6 + 6*a^4*r^4*cos(th)^4 + 4*a^2*r^6*cos(th)^2 + r^8)" ] }, "execution_count": 54, "metadata": {}, "output_type": "execute_result" } ], "source": [ "F2.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

The energy-momentum tensor of the electromagnetic field:

" ] }, { "cell_type": "code", "execution_count": 55, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Tensor field of type (0,2) on the 4-dimensional Lorentzian manifold M\n" ] } ], "source": [ "Fud = F.up(g,0)\n", "T = 1/(4*pi)*( F['_k.']*Fud['^k_.'] - 1/4*F2 * g )\n", "T.apply_map(factor)\n", "print(T)" ] }, { "cell_type": "code", "execution_count": 56, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\left(\\begin{array}{rrrr}\n", "-\\frac{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} - 2 \\, a^{2} - q^{2} + 2 \\, m r - r^{2}\\right)} q^{2}}{8 \\, \\pi {\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{3}} & 0 & 0 & -\\frac{{\\left(2 \\, a^{2} + q^{2} - 2 \\, m r + 2 \\, r^{2}\\right)} a q^{2} \\sin\\left({\\theta}\\right)^{2}}{8 \\, \\pi {\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{3}} \\\\\n", "0 & -\\frac{q^{2}}{8 \\, \\pi {\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)} {\\left(a^{2} + q^{2} - 2 \\, m r + r^{2}\\right)}} & 0 & 0 \\\\\n", "0 & 0 & \\frac{q^{2}}{8 \\, \\pi {\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}} & 0 \\\\\n", "-\\frac{{\\left(2 \\, a^{2} + q^{2} - 2 \\, m r + 2 \\, r^{2}\\right)} a q^{2} \\sin\\left({\\theta}\\right)^{2}}{8 \\, \\pi {\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{3}} & 0 & 0 & \\frac{{\\left(a^{8} \\cos\\left({\\theta}\\right)^{6} + a^{6} r^{2} \\cos\\left({\\theta}\\right)^{6} + 2 \\, a^{8} \\cos\\left({\\theta}\\right)^{4} \\sin\\left({\\theta}\\right)^{2} + a^{6} q^{2} \\cos\\left({\\theta}\\right)^{4} \\sin\\left({\\theta}\\right)^{2} - 2 \\, a^{6} m r \\cos\\left({\\theta}\\right)^{4} \\sin\\left({\\theta}\\right)^{2} + 2 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{4} \\sin\\left({\\theta}\\right)^{2} - 5 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{4} - 5 \\, a^{4} r^{4} \\cos\\left({\\theta}\\right)^{4} - 4 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{2} \\sin\\left({\\theta}\\right)^{2} + 2 \\, a^{4} q^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} \\sin\\left({\\theta}\\right)^{2} - 4 \\, a^{4} m r^{3} \\cos\\left({\\theta}\\right)^{2} \\sin\\left({\\theta}\\right)^{2} - 4 \\, a^{4} r^{4} \\cos\\left({\\theta}\\right)^{2} \\sin\\left({\\theta}\\right)^{2} + 8 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{2} + 11 \\, a^{4} r^{4} \\cos\\left({\\theta}\\right)^{2} + 3 \\, a^{2} r^{6} \\cos\\left({\\theta}\\right)^{2} + 2 \\, a^{4} r^{4} \\sin\\left({\\theta}\\right)^{2} + a^{2} q^{2} r^{4} \\sin\\left({\\theta}\\right)^{2} - 2 \\, a^{2} m r^{5} \\sin\\left({\\theta}\\right)^{2} + 2 \\, a^{2} r^{6} \\sin\\left({\\theta}\\right)^{2} + a^{2} r^{6} + r^{8}\\right)} q^{2} \\sin\\left({\\theta}\\right)^{2}}{8 \\, \\pi {\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{5}}\n", "\\end{array}\\right)\\)" ], "text/latex": [ "$\\displaystyle \\left(\\begin{array}{rrrr}\n", "-\\frac{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} - 2 \\, a^{2} - q^{2} + 2 \\, m r - r^{2}\\right)} q^{2}}{8 \\, \\pi {\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{3}} & 0 & 0 & -\\frac{{\\left(2 \\, a^{2} + q^{2} - 2 \\, m r + 2 \\, r^{2}\\right)} a q^{2} \\sin\\left({\\theta}\\right)^{2}}{8 \\, \\pi {\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{3}} \\\\\n", "0 & -\\frac{q^{2}}{8 \\, \\pi {\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)} {\\left(a^{2} + q^{2} - 2 \\, m r + r^{2}\\right)}} & 0 & 0 \\\\\n", "0 & 0 & \\frac{q^{2}}{8 \\, \\pi {\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}} & 0 \\\\\n", "-\\frac{{\\left(2 \\, a^{2} + q^{2} - 2 \\, m r + 2 \\, r^{2}\\right)} a q^{2} \\sin\\left({\\theta}\\right)^{2}}{8 \\, \\pi {\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{3}} & 0 & 0 & \\frac{{\\left(a^{8} \\cos\\left({\\theta}\\right)^{6} + a^{6} r^{2} \\cos\\left({\\theta}\\right)^{6} + 2 \\, a^{8} \\cos\\left({\\theta}\\right)^{4} \\sin\\left({\\theta}\\right)^{2} + a^{6} q^{2} \\cos\\left({\\theta}\\right)^{4} \\sin\\left({\\theta}\\right)^{2} - 2 \\, a^{6} m r \\cos\\left({\\theta}\\right)^{4} \\sin\\left({\\theta}\\right)^{2} + 2 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{4} \\sin\\left({\\theta}\\right)^{2} - 5 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{4} - 5 \\, a^{4} r^{4} \\cos\\left({\\theta}\\right)^{4} - 4 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{2} \\sin\\left({\\theta}\\right)^{2} + 2 \\, a^{4} q^{2} r^{2} \\cos\\left({\\theta}\\right)^{2} \\sin\\left({\\theta}\\right)^{2} - 4 \\, a^{4} m r^{3} \\cos\\left({\\theta}\\right)^{2} \\sin\\left({\\theta}\\right)^{2} - 4 \\, a^{4} r^{4} \\cos\\left({\\theta}\\right)^{2} \\sin\\left({\\theta}\\right)^{2} + 8 \\, a^{6} r^{2} \\cos\\left({\\theta}\\right)^{2} + 11 \\, a^{4} r^{4} \\cos\\left({\\theta}\\right)^{2} + 3 \\, a^{2} r^{6} \\cos\\left({\\theta}\\right)^{2} + 2 \\, a^{4} r^{4} \\sin\\left({\\theta}\\right)^{2} + a^{2} q^{2} r^{4} \\sin\\left({\\theta}\\right)^{2} - 2 \\, a^{2} m r^{5} \\sin\\left({\\theta}\\right)^{2} + 2 \\, a^{2} r^{6} \\sin\\left({\\theta}\\right)^{2} + a^{2} r^{6} + r^{8}\\right)} q^{2} \\sin\\left({\\theta}\\right)^{2}}{8 \\, \\pi {\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{5}}\n", "\\end{array}\\right)$" ], "text/plain": [ "[ -1/8*(a^2*cos(th)^2 - 2*a^2 - q^2 + 2*m*r - r^2)*q^2/(pi*(a^2*cos(th)^2 + r^2)^3) 0 0 -1/8*(2*a^2 + q^2 - 2*m*r + 2*r^2)*a*q^2*sin(th)^2/(pi*(a^2*cos(th)^2 + r^2)^3)]\n", "[ 0 -1/8*q^2/(pi*(a^2*cos(th)^2 + r^2)*(a^2 + q^2 - 2*m*r + r^2)) 0 0]\n", "[ 0 0 1/8*q^2/(pi*(a^2*cos(th)^2 + r^2)) 0]\n", "[ -1/8*(2*a^2 + q^2 - 2*m*r + 2*r^2)*a*q^2*sin(th)^2/(pi*(a^2*cos(th)^2 + r^2)^3) 0 0 1/8*(a^8*cos(th)^6 + a^6*r^2*cos(th)^6 + 2*a^8*cos(th)^4*sin(th)^2 + a^6*q^2*cos(th)^4*sin(th)^2 - 2*a^6*m*r*cos(th)^4*sin(th)^2 + 2*a^6*r^2*cos(th)^4*sin(th)^2 - 5*a^6*r^2*cos(th)^4 - 5*a^4*r^4*cos(th)^4 - 4*a^6*r^2*cos(th)^2*sin(th)^2 + 2*a^4*q^2*r^2*cos(th)^2*sin(th)^2 - 4*a^4*m*r^3*cos(th)^2*sin(th)^2 - 4*a^4*r^4*cos(th)^2*sin(th)^2 + 8*a^6*r^2*cos(th)^2 + 11*a^4*r^4*cos(th)^2 + 3*a^2*r^6*cos(th)^2 + 2*a^4*r^4*sin(th)^2 + a^2*q^2*r^4*sin(th)^2 - 2*a^2*m*r^5*sin(th)^2 + 2*a^2*r^6*sin(th)^2 + a^2*r^6 + r^8)*q^2*sin(th)^2/(pi*(a^2*cos(th)^2 + r^2)^5)]" ] }, "execution_count": 56, "metadata": {}, "output_type": "execute_result" } ], "source": [ "T[:]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Check that the Einstein equation is satisfied:" ] }, { "cell_type": "code", "execution_count": 57, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\mathrm{True}\\)" ], "text/latex": [ "$\\displaystyle \\mathrm{True}$" ], "text/plain": [ "True" ] }, "execution_count": 57, "metadata": {}, "output_type": "execute_result" } ], "source": [ "G == 8*pi*T" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Bianchi identity

\n", "\n", "

Let us check the Bianchi identity $\\nabla_p R^i_{\\ \\, j kl} + \\nabla_k R^i_{\\ \\, jlp} + \\nabla_l R^i_{\\ \\, jpk} = 0$:

" ] }, { "cell_type": "code", "execution_count": 58, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Tensor field nabla_g(Riem(g)) of type (1,4) on the 4-dimensional Lorentzian manifold M\n" ] } ], "source": [ "DR = nabla(R) # long (takes a while)\n", "print(DR) " ] }, { "cell_type": "code", "execution_count": 59, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 " ] } ], "source": [ "for i in M.irange():\n", " for j in M.irange():\n", " for k in M.irange():\n", " for l in M.irange():\n", " for p in M.irange():\n", " print(DR[i,j,k,l,p] + DR[i,j,l,p,k] + DR[i,j,p,k,l], end=' ')" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

If the last sign in the Bianchi identity is changed to minus, the identity does no longer hold:

" ] }, { "cell_type": "code", "execution_count": 60, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle 0\\)" ], "text/latex": [ "$\\displaystyle 0$" ], "text/plain": [ "0" ] }, "execution_count": 60, "metadata": {}, "output_type": "execute_result" } ], "source": [ "DR[0,1,2,3,1] + DR[0,1,3,1,2] + DR[0,1,1,2,3] # should be zero (Bianchi identity)" ] }, { "cell_type": "code", "execution_count": 61, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle -\\frac{4 \\, {\\left({\\left(a^{5} q^{2} - 6 \\, a^{5} m r + a^{3} q^{2} r^{2} - 6 \\, a^{3} m r^{3}\\right)} \\cos\\left({\\theta}\\right)^{3} - {\\left(5 \\, a^{3} q^{2} r^{2} - 6 \\, a^{3} m r^{3} + 5 \\, a q^{2} r^{4} - 6 \\, a m r^{5}\\right)} \\cos\\left({\\theta}\\right)\\right)} \\sin\\left({\\theta}\\right)}{2 \\, m r^{7} - r^{8} - {\\left(a^{2} + q^{2}\\right)} r^{6} - {\\left(a^{8} + a^{6} q^{2} - 2 \\, a^{6} m r + a^{6} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{6} + 3 \\, {\\left(2 \\, a^{4} m r^{3} - a^{4} r^{4} - {\\left(a^{6} + a^{4} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + 3 \\, {\\left(2 \\, a^{2} m r^{5} - a^{2} r^{6} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{4}\\right)} \\cos\\left({\\theta}\\right)^{2}}\\)" ], "text/latex": [ "$\\displaystyle -\\frac{4 \\, {\\left({\\left(a^{5} q^{2} - 6 \\, a^{5} m r + a^{3} q^{2} r^{2} - 6 \\, a^{3} m r^{3}\\right)} \\cos\\left({\\theta}\\right)^{3} - {\\left(5 \\, a^{3} q^{2} r^{2} - 6 \\, a^{3} m r^{3} + 5 \\, a q^{2} r^{4} - 6 \\, a m r^{5}\\right)} \\cos\\left({\\theta}\\right)\\right)} \\sin\\left({\\theta}\\right)}{2 \\, m r^{7} - r^{8} - {\\left(a^{2} + q^{2}\\right)} r^{6} - {\\left(a^{8} + a^{6} q^{2} - 2 \\, a^{6} m r + a^{6} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{6} + 3 \\, {\\left(2 \\, a^{4} m r^{3} - a^{4} r^{4} - {\\left(a^{6} + a^{4} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{4} + 3 \\, {\\left(2 \\, a^{2} m r^{5} - a^{2} r^{6} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{4}\\right)} \\cos\\left({\\theta}\\right)^{2}}$" ], "text/plain": [ "-4*((a^5*q^2 - 6*a^5*m*r + a^3*q^2*r^2 - 6*a^3*m*r^3)*cos(th)^3 - (5*a^3*q^2*r^2 - 6*a^3*m*r^3 + 5*a*q^2*r^4 - 6*a*m*r^5)*cos(th))*sin(th)/(2*m*r^7 - r^8 - (a^2 + q^2)*r^6 - (a^8 + a^6*q^2 - 2*a^6*m*r + a^6*r^2)*cos(th)^6 + 3*(2*a^4*m*r^3 - a^4*r^4 - (a^6 + a^4*q^2)*r^2)*cos(th)^4 + 3*(2*a^2*m*r^5 - a^2*r^6 - (a^4 + a^2*q^2)*r^4)*cos(th)^2)" ] }, "execution_count": 61, "metadata": {}, "output_type": "execute_result" } ], "source": [ "DR[0,1,2,3,1] + DR[0,1,3,1,2] - DR[0,1,1,2,3] # note the change of the second + to -" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Kretschmann scalar\n", "\n", "The tensor $R^\\flat$, of components $R_{abcd} = g_{am} R^m_{\\ \\, bcd}$:" ] }, { "cell_type": "code", "execution_count": 62, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Tensor field of type (0,4) on the 4-dimensional Lorentzian manifold M\n" ] } ], "source": [ "dR = R.down(g)\n", "print(dR)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The tensor $R^\\sharp$, of components $R^{abcd} = g^{bp} g^{cq} g^{dr} R^a_{\\ \\, pqr}$:" ] }, { "cell_type": "code", "execution_count": 63, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Tensor field of type (4,0) on the 4-dimensional Lorentzian manifold M\n" ] } ], "source": [ "uR = R.up(g)\n", "print(uR)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The Kretschmann scalar $K := R^{abcd} R_{abcd}$:" ] }, { "cell_type": "code", "execution_count": 64, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\begin{array}{llcl} & \\mathcal{M} & \\longrightarrow & \\mathbb{R} \\\\ & \\left(t, r, {\\theta}, {\\phi}\\right) & \\longmapsto & -\\frac{8 \\, {\\left(6 \\, m^{2} r^{8} - 12 \\, {\\left(m^{3} + m q^{2}\\right)} r^{7} + {\\left(6 \\, a^{2} m^{2} + 30 \\, m^{2} q^{2} + 7 \\, q^{4}\\right)} r^{6} - 6 \\, {\\left(a^{8} m^{2} + a^{6} m^{2} q^{2} - 2 \\, a^{6} m^{3} r + a^{6} m^{2} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{6} - 2 \\, {\\left(6 \\, a^{2} m q^{2} + 13 \\, m q^{4}\\right)} r^{5} + 7 \\, {\\left(a^{2} q^{4} + q^{6}\\right)} r^{4} + {\\left(7 \\, a^{6} q^{4} + 7 \\, a^{4} q^{6} + 90 \\, a^{4} m^{2} r^{4} - 60 \\, {\\left(3 \\, a^{4} m^{3} + a^{4} m q^{2}\\right)} r^{3} + {\\left(90 \\, a^{6} m^{2} + 210 \\, a^{4} m^{2} q^{2} + 7 \\, a^{4} q^{4}\\right)} r^{2} - 2 \\, {\\left(30 \\, a^{6} m q^{2} + 37 \\, a^{4} m q^{4}\\right)} r\\right)} \\cos\\left({\\theta}\\right)^{4} - 2 \\, {\\left(45 \\, a^{2} m^{2} r^{6} - 30 \\, {\\left(3 \\, a^{2} m^{3} + 2 \\, a^{2} m q^{2}\\right)} r^{5} + {\\left(45 \\, a^{4} m^{2} + 165 \\, a^{2} m^{2} q^{2} + 17 \\, a^{2} q^{4}\\right)} r^{4} - 2 \\, {\\left(30 \\, a^{4} m q^{2} + 47 \\, a^{2} m q^{4}\\right)} r^{3} + 17 \\, {\\left(a^{4} q^{4} + a^{2} q^{6}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}\\right)}}{2 \\, m r^{13} - r^{14} - {\\left(a^{2} + q^{2}\\right)} r^{12} - {\\left(a^{14} + a^{12} q^{2} - 2 \\, a^{12} m r + a^{12} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{12} + 6 \\, {\\left(2 \\, a^{10} m r^{3} - a^{10} r^{4} - {\\left(a^{12} + a^{10} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{10} + 15 \\, {\\left(2 \\, a^{8} m r^{5} - a^{8} r^{6} - {\\left(a^{10} + a^{8} q^{2}\\right)} r^{4}\\right)} \\cos\\left({\\theta}\\right)^{8} + 20 \\, {\\left(2 \\, a^{6} m r^{7} - a^{6} r^{8} - {\\left(a^{8} + a^{6} q^{2}\\right)} r^{6}\\right)} \\cos\\left({\\theta}\\right)^{6} + 15 \\, {\\left(2 \\, a^{4} m r^{9} - a^{4} r^{10} - {\\left(a^{6} + a^{4} q^{2}\\right)} r^{8}\\right)} \\cos\\left({\\theta}\\right)^{4} + 6 \\, {\\left(2 \\, a^{2} m r^{11} - a^{2} r^{12} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{10}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\end{array}\\)" ], "text/latex": [ "$\\displaystyle \\begin{array}{llcl} & \\mathcal{M} & \\longrightarrow & \\mathbb{R} \\\\ & \\left(t, r, {\\theta}, {\\phi}\\right) & \\longmapsto & -\\frac{8 \\, {\\left(6 \\, m^{2} r^{8} - 12 \\, {\\left(m^{3} + m q^{2}\\right)} r^{7} + {\\left(6 \\, a^{2} m^{2} + 30 \\, m^{2} q^{2} + 7 \\, q^{4}\\right)} r^{6} - 6 \\, {\\left(a^{8} m^{2} + a^{6} m^{2} q^{2} - 2 \\, a^{6} m^{3} r + a^{6} m^{2} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{6} - 2 \\, {\\left(6 \\, a^{2} m q^{2} + 13 \\, m q^{4}\\right)} r^{5} + 7 \\, {\\left(a^{2} q^{4} + q^{6}\\right)} r^{4} + {\\left(7 \\, a^{6} q^{4} + 7 \\, a^{4} q^{6} + 90 \\, a^{4} m^{2} r^{4} - 60 \\, {\\left(3 \\, a^{4} m^{3} + a^{4} m q^{2}\\right)} r^{3} + {\\left(90 \\, a^{6} m^{2} + 210 \\, a^{4} m^{2} q^{2} + 7 \\, a^{4} q^{4}\\right)} r^{2} - 2 \\, {\\left(30 \\, a^{6} m q^{2} + 37 \\, a^{4} m q^{4}\\right)} r\\right)} \\cos\\left({\\theta}\\right)^{4} - 2 \\, {\\left(45 \\, a^{2} m^{2} r^{6} - 30 \\, {\\left(3 \\, a^{2} m^{3} + 2 \\, a^{2} m q^{2}\\right)} r^{5} + {\\left(45 \\, a^{4} m^{2} + 165 \\, a^{2} m^{2} q^{2} + 17 \\, a^{2} q^{4}\\right)} r^{4} - 2 \\, {\\left(30 \\, a^{4} m q^{2} + 47 \\, a^{2} m q^{4}\\right)} r^{3} + 17 \\, {\\left(a^{4} q^{4} + a^{2} q^{6}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{2}\\right)}}{2 \\, m r^{13} - r^{14} - {\\left(a^{2} + q^{2}\\right)} r^{12} - {\\left(a^{14} + a^{12} q^{2} - 2 \\, a^{12} m r + a^{12} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{12} + 6 \\, {\\left(2 \\, a^{10} m r^{3} - a^{10} r^{4} - {\\left(a^{12} + a^{10} q^{2}\\right)} r^{2}\\right)} \\cos\\left({\\theta}\\right)^{10} + 15 \\, {\\left(2 \\, a^{8} m r^{5} - a^{8} r^{6} - {\\left(a^{10} + a^{8} q^{2}\\right)} r^{4}\\right)} \\cos\\left({\\theta}\\right)^{8} + 20 \\, {\\left(2 \\, a^{6} m r^{7} - a^{6} r^{8} - {\\left(a^{8} + a^{6} q^{2}\\right)} r^{6}\\right)} \\cos\\left({\\theta}\\right)^{6} + 15 \\, {\\left(2 \\, a^{4} m r^{9} - a^{4} r^{10} - {\\left(a^{6} + a^{4} q^{2}\\right)} r^{8}\\right)} \\cos\\left({\\theta}\\right)^{4} + 6 \\, {\\left(2 \\, a^{2} m r^{11} - a^{2} r^{12} - {\\left(a^{4} + a^{2} q^{2}\\right)} r^{10}\\right)} \\cos\\left({\\theta}\\right)^{2}} \\end{array}$" ], "text/plain": [ "M → ℝ\n", "(t, r, th, ph) ↦ -8*(6*m^2*r^8 - 12*(m^3 + m*q^2)*r^7 + (6*a^2*m^2 + 30*m^2*q^2 + 7*q^4)*r^6 - 6*(a^8*m^2 + a^6*m^2*q^2 - 2*a^6*m^3*r + a^6*m^2*r^2)*cos(th)^6 - 2*(6*a^2*m*q^2 + 13*m*q^4)*r^5 + 7*(a^2*q^4 + q^6)*r^4 + (7*a^6*q^4 + 7*a^4*q^6 + 90*a^4*m^2*r^4 - 60*(3*a^4*m^3 + a^4*m*q^2)*r^3 + (90*a^6*m^2 + 210*a^4*m^2*q^2 + 7*a^4*q^4)*r^2 - 2*(30*a^6*m*q^2 + 37*a^4*m*q^4)*r)*cos(th)^4 - 2*(45*a^2*m^2*r^6 - 30*(3*a^2*m^3 + 2*a^2*m*q^2)*r^5 + (45*a^4*m^2 + 165*a^2*m^2*q^2 + 17*a^2*q^4)*r^4 - 2*(30*a^4*m*q^2 + 47*a^2*m*q^4)*r^3 + 17*(a^4*q^4 + a^2*q^6)*r^2)*cos(th)^2)/(2*m*r^13 - r^14 - (a^2 + q^2)*r^12 - (a^14 + a^12*q^2 - 2*a^12*m*r + a^12*r^2)*cos(th)^12 + 6*(2*a^10*m*r^3 - a^10*r^4 - (a^12 + a^10*q^2)*r^2)*cos(th)^10 + 15*(2*a^8*m*r^5 - a^8*r^6 - (a^10 + a^8*q^2)*r^4)*cos(th)^8 + 20*(2*a^6*m*r^7 - a^6*r^8 - (a^8 + a^6*q^2)*r^6)*cos(th)^6 + 15*(2*a^4*m*r^9 - a^4*r^10 - (a^6 + a^4*q^2)*r^8)*cos(th)^4 + 6*(2*a^2*m*r^11 - a^2*r^12 - (a^4 + a^2*q^2)*r^10)*cos(th)^2)" ] }, "execution_count": 64, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Kr_scalar = uR['^abcd']*dR['_abcd']\n", "Kr_scalar.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "A variant of this expression can be obtained by invoking the method `factor()` on the coordinate function representing the scalar field in the manifold's default chart:" ] }, { "cell_type": "code", "execution_count": 65, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle -\\frac{8 \\, {\\left(6 \\, a^{6} m^{2} \\cos\\left({\\theta}\\right)^{6} - 7 \\, a^{4} q^{4} \\cos\\left({\\theta}\\right)^{4} + 60 \\, a^{4} m q^{2} r \\cos\\left({\\theta}\\right)^{4} - 90 \\, a^{4} m^{2} r^{2} \\cos\\left({\\theta}\\right)^{4} + 34 \\, a^{2} q^{4} r^{2} \\cos\\left({\\theta}\\right)^{2} - 120 \\, a^{2} m q^{2} r^{3} \\cos\\left({\\theta}\\right)^{2} + 90 \\, a^{2} m^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} - 7 \\, q^{4} r^{4} + 12 \\, m q^{2} r^{5} - 6 \\, m^{2} r^{6}\\right)}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{6}}\\)" ], "text/latex": [ "$\\displaystyle -\\frac{8 \\, {\\left(6 \\, a^{6} m^{2} \\cos\\left({\\theta}\\right)^{6} - 7 \\, a^{4} q^{4} \\cos\\left({\\theta}\\right)^{4} + 60 \\, a^{4} m q^{2} r \\cos\\left({\\theta}\\right)^{4} - 90 \\, a^{4} m^{2} r^{2} \\cos\\left({\\theta}\\right)^{4} + 34 \\, a^{2} q^{4} r^{2} \\cos\\left({\\theta}\\right)^{2} - 120 \\, a^{2} m q^{2} r^{3} \\cos\\left({\\theta}\\right)^{2} + 90 \\, a^{2} m^{2} r^{4} \\cos\\left({\\theta}\\right)^{2} - 7 \\, q^{4} r^{4} + 12 \\, m q^{2} r^{5} - 6 \\, m^{2} r^{6}\\right)}}{{\\left(a^{2} \\cos\\left({\\theta}\\right)^{2} + r^{2}\\right)}^{6}}$" ], "text/plain": [ "-8*(6*a^6*m^2*cos(th)^6 - 7*a^4*q^4*cos(th)^4 + 60*a^4*m*q^2*r*cos(th)^4 - 90*a^4*m^2*r^2*cos(th)^4 + 34*a^2*q^4*r^2*cos(th)^2 - 120*a^2*m*q^2*r^3*cos(th)^2 + 90*a^2*m^2*r^4*cos(th)^2 - 7*q^4*r^4 + 12*m*q^2*r^5 - 6*m^2*r^6)/(a^2*cos(th)^2 + r^2)^6" ] }, "execution_count": 65, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Kr = Kr_scalar.coord_function()\n", "Kr.factor()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

As a check, we can compare Kr to the formula given by R. Conn Henry, Astrophys. J. 535, 350 (2000):

" ] }, { "cell_type": "code", "execution_count": 66, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\mathrm{True}\\)" ], "text/latex": [ "$\\displaystyle \\mathrm{True}$" ], "text/plain": [ "True" ] }, "execution_count": 66, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Kr == 8/(r^2+(a*cos(th))^2)^6 *( \n", " 6*m^2*(r^6 - 15*r^4*(a*cos(th))^2 + 15*r^2*(a*cos(th))^4 - (a*cos(th))^6) \n", " - 12*m*q^2*r*(r^4 - 10*(a*r*cos(th))^2 + 5*(a*cos(th))^4) \n", " + q^4*(7*r^4 - 34*(a*r*cos(th))^2 + 7*(a*cos(th))^4) )" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

The Schwarzschild value of the Kretschmann scalar is recovered by setting $a=0$ and $q=0$:

" ] }, { "cell_type": "code", "execution_count": 67, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\\(\\displaystyle \\frac{48 \\, m^{2}}{r^{6}}\\)" ], "text/latex": [ "$\\displaystyle \\frac{48 \\, m^{2}}{r^{6}}$" ], "text/plain": [ "48*m^2/r^6" ] }, "execution_count": 67, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Kr.expr().subs(a=0, q=0)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Let us plot the Kretschmann scalar for $m=1$, $a=0.9$ and $q=0.5$:

" ] }, { "cell_type": "code", "execution_count": 68, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\n", "\n" ], "text/plain": [ "Graphics3d Object" ] }, "execution_count": 68, "metadata": {}, "output_type": "execute_result" } ], "source": [ "K1 = Kr.expr().subs(m=1, a=0.9, q=0.5)\n", "plot3d(K1, (r,1,3), (th, 0, pi), axes_labels=['r', 'theta', 'Kr'])" ] }, { "cell_type": "code", "execution_count": 69, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Total elapsed time: 1108.1922762840022 s\n" ] } ], "source": [ "print(\"Total elapsed time: {} s\".format(time.perf_counter() - comput_time0))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "*NB:* most of the computational time is spent in checking the Bianchi identity." ] } ], "metadata": { "kernelspec": { "display_name": "SageMath 10.1", "language": "sage", "name": "sagemath" }, "language": "python", "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.12" } }, "nbformat": 4, "nbformat_minor": 4 }