{ "cells": [ { "cell_type": "markdown", "metadata": { "collapsed": false }, "source": [ "# Vaidya-Lifshitz solution" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "This Jupyter/SageMath \n", "worksheet implements some computations of the article\n", "- I. Ya. Aref'eva, A. A. Golubtsova & E. Gourgoulhon: *Analytic black branes in \n", " Lifshitz-like backgrounds and thermalization*,\n", " [arXiv:1601.06046](http://arxiv.org/abs/1601.06046)\n", "\n", "These computations are based on [SageManifolds](http://sagemanifolds.obspm.fr) (v0.9).\n", "\n", "The worksheet file (ipynb format) can be downloaded from [here](https://cloud.sagemath.com/projects/3edbca82-97d6-41b3-9b6f-d83ea06fc1e9/files/Vaidya-Lifshitz.html)." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "First we set up the notebook to display mathematical objects using LaTeX formatting:" ] }, { "cell_type": "code", "execution_count": 1, "metadata": { "collapsed": false }, "outputs": [], "source": [ "%display latex" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Spacetime manifold and coordinates\n", "\n", "Let us declare the spacetime $M$ as a 5-dimensional manifold:" ] }, { "cell_type": "code", "execution_count": 2, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "5-dimensional differentiable manifold M\n" ] } ], "source": [ "M = Manifold(5, 'M')\n", "print M" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We introduce coordinates of Eddington-Finkelstein type:" ] }, { "cell_type": "code", "execution_count": 3, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "Chart (M, (v, x, y1, y2, r))" ] }, "execution_count": 3, "metadata": {}, "output_type": "execute_result" } ], "source": [ "X. = M.chart('v x y1:y_1 y2:y_2 r')\n", "X" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The metric tensor:" ] }, { "cell_type": "code", "execution_count": 4, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "g = -e^(2*nu*r)*f(v, r) dv*dv + e^(nu*r) dv*dr + e^(2*nu*r) dx*dx + e^(2*r) dy1*dy1 + e^(2*r) dy2*dy2 + e^(nu*r) dr*dv" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" } ], "source": [ "g = M.lorentzian_metric('g')\n", "nu = var('nu', latex_name=r'\\nu', domain='real')\n", "ff = function('f')(v, r)\n", "g[0,0] = -exp(2*nu*r)*ff\n", "g[0,4] = exp(nu*r)\n", "g[1,1] = exp(2*nu*r)\n", "g[2,2] = exp(2*r)\n", "g[3,3] = exp(2*r)\n", "g.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The non-vanishing components of $g$:" ] }, { "cell_type": "code", "execution_count": 5, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "g_v,v = -e^(2*nu*r)*f(v, r) \n", "g_v,r = e^(nu*r) \n", "g_x,x = e^(2*nu*r) \n", "g_y1,y1 = e^(2*r) \n", "g_y2,y2 = e^(2*r) \n", "g_r,v = e^(nu*r) " ] }, "execution_count": 5, "metadata": {}, "output_type": "execute_result" } ], "source": [ "g.display_comp()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "A matrix view of the components:" ] }, { "cell_type": "code", "execution_count": 6, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "[-e^(2*nu*r)*f(v, r) 0 0 0 e^(nu*r)]\n", "[ 0 e^(2*nu*r) 0 0 0]\n", "[ 0 0 e^(2*r) 0 0]\n", "[ 0 0 0 e^(2*r) 0]\n", "[ e^(nu*r) 0 0 0 0]" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "g[:]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Curvature\n", "\n", "The Riemann tensor is" ] }, { "cell_type": "code", "execution_count": 7, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Tensor field Riem(g) of type (1,3) on the 5-dimensional differentiable manifold M\n" ] } ], "source": [ "Riem = g.riemann()\n", "print Riem" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Some component, e.g. $\\mathrm{Riem}^0_{\\ \\, 004}$:" ] }, { "cell_type": "code", "execution_count": 8, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "-nu^2*e^(nu*r)*f(v, r) - 3/2*nu*e^(nu*r)*d(f)/dr - 1/2*e^(nu*r)*d^2(f)/dr^2" ] }, "execution_count": 8, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Riem[0,0,0,4]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The Ricci tensor:" ] }, { "cell_type": "code", "execution_count": 9, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Field of symmetric bilinear forms Ric(g) on the 5-dimensional differentiable manifold M\n" ] } ], "source": [ "Ric = g.ricci()\n", "print Ric" ] }, { "cell_type": "code", "execution_count": 10, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "Ric(g) = (2*(nu^2 + nu)*e^(2*nu*r)*f(v, r)^2 + (2*nu + 1)*e^(2*nu*r)*f(v, r)*d(f)/dr - 1/2*(nu + 2)*e^(nu*r)*d(f)/dv + 1/2*e^(2*nu*r)*f(v, r)*d^2(f)/dr^2) dv*dv + (-2*(nu^2 + nu)*e^(nu*r)*f(v, r) - (2*nu + 1)*e^(nu*r)*d(f)/dr - 1/2*e^(nu*r)*d^2(f)/dr^2) dv*dr + (-2*(nu^2 + nu)*e^(2*nu*r)*f(v, r) - nu*e^(2*nu*r)*d(f)/dr) dx*dx + (-2*(nu + 1)*e^(2*r)*f(v, r) - e^(2*r)*d(f)/dr) dy1*dy1 + (-2*(nu + 1)*e^(2*r)*f(v, r) - e^(2*r)*d(f)/dr) dy2*dy2 + (-2*(nu^2 + nu)*e^(nu*r)*f(v, r) - (2*nu + 1)*e^(nu*r)*d(f)/dr - 1/2*e^(nu*r)*d^2(f)/dr^2) dr*dv + (2*nu - 2) dr*dr" ] }, "execution_count": 10, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Ric.display()" ] }, { "cell_type": "code", "execution_count": 11, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "Ric(g)_v,v = 2*(nu^2 + nu)*e^(2*nu*r)*f(v, r)^2 + (2*nu + 1)*e^(2*nu*r)*f(v, r)*d(f)/dr - 1/2*(nu + 2)*e^(nu*r)*d(f)/dv + 1/2*e^(2*nu*r)*f(v, r)*d^2(f)/dr^2 \n", "Ric(g)_v,r = -2*(nu^2 + nu)*e^(nu*r)*f(v, r) - (2*nu + 1)*e^(nu*r)*d(f)/dr - 1/2*e^(nu*r)*d^2(f)/dr^2 \n", "Ric(g)_x,x = -2*(nu^2 + nu)*e^(2*nu*r)*f(v, r) - nu*e^(2*nu*r)*d(f)/dr \n", "Ric(g)_y1,y1 = -2*(nu + 1)*e^(2*r)*f(v, r) - e^(2*r)*d(f)/dr \n", "Ric(g)_y2,y2 = -2*(nu + 1)*e^(2*r)*f(v, r) - e^(2*r)*d(f)/dr \n", "Ric(g)_r,v = -2*(nu^2 + nu)*e^(nu*r)*f(v, r) - (2*nu + 1)*e^(nu*r)*d(f)/dr - 1/2*e^(nu*r)*d^2(f)/dr^2 \n", "Ric(g)_r,r = 2*nu - 2 " ] }, "execution_count": 11, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Ric.display_comp()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The Ricci scalar:" ] }, { "cell_type": "code", "execution_count": 12, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Scalar field r(g) on the 5-dimensional differentiable manifold M\n" ] } ], "source": [ "Rscal = g.ricci_scalar()\n", "print Rscal" ] }, { "cell_type": "code", "execution_count": 13, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "r(g): M --> R\n", " (v, x, y1, y2, r) |--> -2*(3*nu^2 + 4*nu + 3)*f(v, r) - (5*nu + 4)*d(f)/dr - d^2(f)/dr^2" ] }, "execution_count": 13, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Rscal.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Source model\n", "Let us consider a model based on the following action, involving a dilaton scalar field $\\phi$ and a Maxwell 2-form $F$:\n", "\n", "$$ S = \\int \\left( R(g) + \\Lambda - \\frac{1}{2} \\nabla_m \\phi \\nabla^m \\phi - \\frac{1}{4} e^{\\lambda\\phi} F_{mn} F^{mn} + 8\\pi \\mathcal{L}_{\\rm shell} \\right) \\sqrt{-g} \\, \\mathrm{d}^5 x \\qquad\\qquad \\mbox{(1)}$$\n", "\n", "where $R(g)$ is the Ricci scalar of metric $g$, $\\Lambda$ is the cosmological constant, $\\lambda$ is the dilatonic coupling constant and $\\mathcal{L}_{\\rm shell}$ is the Lagrangian of some infalling shell of massless matter." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### The dilaton scalar field\n", "\n", "We consider the following ansatz for the dilaton scalar field $\\phi$:\n", "$$ \\phi = \\frac{1}{\\lambda} \\left( 4 r + \\ln\\mu \\right),$$\n", "where $\\mu$ is a constant. " ] }, { "cell_type": "code", "execution_count": 14, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "phi: M --> R\n", " (v, x, y1, y2, r) |--> (4*r + log(mu))/lamb" ] }, "execution_count": 14, "metadata": {}, "output_type": "execute_result" } ], "source": [ "var('mu', latex_name=r'\\mu', domain='real')\n", "var('lamb', latex_name=r'\\lambda', domain='real')\n", "phi = M.scalar_field({X: (4*r + ln(mu))/lamb}, \n", " name='phi', latex_name=r'\\phi')\n", "phi.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The 1-form $\\mathrm{d}\\phi$ is" ] }, { "cell_type": "code", "execution_count": 15, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "1-form dphi on the 5-dimensional differentiable manifold M\n" ] } ], "source": [ "dphi = phi.differential()\n", "print dphi" ] }, { "cell_type": "code", "execution_count": 16, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "dphi = 4/lamb dr" ] }, "execution_count": 16, "metadata": {}, "output_type": "execute_result" } ], "source": [ "dphi.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The components of $\\mathrm{d}\\phi$ is the default frame of $M$ are" ] }, { "cell_type": "code", "execution_count": 17, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "[0, 0, 0, 0, 4/lamb]" ] }, "execution_count": 17, "metadata": {}, "output_type": "execute_result" } ], "source": [ "dphi[:]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### The 2-form field\n", "\n", "We consider the following ansatz for $F$:\n", "$$ F = \\frac{1}{2} q \\, \\mathrm{d}y_1\\wedge \\mathrm{d}y_2, $$\n", "where $q$ is a constant. \n", "\n", "Let us first get the 1-forms $\\mathrm{d}y_1$ and $\\mathrm{d}y_2$:" ] }, { "cell_type": "code", "execution_count": 18, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "Coordinate coframe (M, (dv,dx,dy1,dy2,dr))" ] }, "execution_count": 18, "metadata": {}, "output_type": "execute_result" } ], "source": [ "X.coframe()" ] }, { "cell_type": "code", "execution_count": 19, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "1-form dy1 on the 5-dimensional differentiable manifold M\n", "1-form dy2 on the 5-dimensional differentiable manifold M\n" ] }, { "data": { "text/html": [ "" ], "text/plain": [ "(1-form dy1 on the 5-dimensional differentiable manifold M,\n", " 1-form dy2 on the 5-dimensional differentiable manifold M)" ] }, "execution_count": 19, "metadata": {}, "output_type": "execute_result" } ], "source": [ "dy1 = X.coframe()[2]\n", "dy2 = X.coframe()[3]\n", "print dy1\n", "print dy2\n", "dy1, dy2" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Then we can form $F$ according to the above ansatz:" ] }, { "cell_type": "code", "execution_count": 20, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "2-form F on the 5-dimensional differentiable manifold M\n" ] }, { "data": { "text/html": [ "" ], "text/plain": [ "F = 1/2*q dy1/\\dy2" ] }, "execution_count": 20, "metadata": {}, "output_type": "execute_result" } ], "source": [ "var('q', domain='real')\n", "F = q/2 * dy1.wedge(dy2)\n", "F.set_name('F')\n", "print F\n", "F.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "By construction, the 2-form $F$ is closed (since $q$ is constant):" ] }, { "cell_type": "code", "execution_count": 21, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "dF = 0" ] }, "execution_count": 21, "metadata": {}, "output_type": "execute_result" } ], "source": [ "xder(F).display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Let us evaluate the square $F_{mn} F^{mn}$ of $F$:" ] }, { "cell_type": "code", "execution_count": 22, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Scalar field on the 5-dimensional differentiable manifold M\n" ] }, { "data": { "text/html": [ "" ], "text/plain": [ "M --> R\n", "(v, x, y1, y2, r) |--> 1/2*q^2*e^(-4*r)" ] }, "execution_count": 22, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Fu = F.up(g)\n", "F2 = F['_{mn}']*Fu['^{mn}'] # using LaTeX notations to denote contraction\n", "print F2\n", "F2.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We shall also need the tensor $\\mathcal{F}_{mn} := F_{mp} F_n^{\\ \\, p}$:" ] }, { "cell_type": "code", "execution_count": 23, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Tensor field of type (0,2) on the 5-dimensional differentiable manifold M\n" ] }, { "data": { "text/html": [ "" ], "text/plain": [ "1/4*q^2*e^(-2*r) dy1*dy1 + 1/4*q^2*e^(-2*r) dy2*dy2" ] }, "execution_count": 23, "metadata": {}, "output_type": "execute_result" } ], "source": [ "FF = F['_mp'] * F.up(g,1)['^p_n']\n", "print FF\n", "FF.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The tensor field $\\mathcal{F}$ is symmetric:" ] }, { "cell_type": "code", "execution_count": 24, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "True" ] }, "execution_count": 24, "metadata": {}, "output_type": "execute_result" } ], "source": [ "FF == FF.symmetrize()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Therefore, from now on, we set" ] }, { "cell_type": "code", "execution_count": 25, "metadata": { "collapsed": true }, "outputs": [], "source": [ "FF = FF.symmetrize()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### The infalling shell of massless particles\n", "\n", "Energy-momentum tensor of the infalling shell:" ] }, { "cell_type": "code", "execution_count": 26, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "T = T00(v, r) dv*dv" ] }, "execution_count": 26, "metadata": {}, "output_type": "execute_result" } ], "source": [ "dv = X.coframe()[0]\n", "T = function('T00', latex_name='T_{00}')(v,r)*dv*dv\n", "T.set_name('T')\n", "T.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Einstein equation\n", "\n", "Let us first introduce (minus twice) the cosmological constant:" ] }, { "cell_type": "code", "execution_count": 27, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "Lamb" ] }, "execution_count": 27, "metadata": {}, "output_type": "execute_result" } ], "source": [ "var('Lamb', latex_name=r'\\Lambda', domain='real')" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "From the action (1), the field equation for the metric $g$ is\n", "$$ R_{mn} + \\frac{\\Lambda}{3} \\, g - \\frac{1}{2}\\partial_m\\phi \\partial_n\\phi -\\frac{1}{2} e^{\\lambda\\phi} F_{mp} F^{\\ \\, p}_n + \\frac{1}{12} e^{\\lambda\\phi} F_{rs} F^{rs} \\, g_{mn} - 8\\pi T_{mn} = 0 \n", "$$\n", "We write it as\n", "\n", " EE == 0\n", "\n", "with `EE` defined by" ] }, { "cell_type": "code", "execution_count": 28, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Tensor field E of type (0,2) on the 5-dimensional differentiable manifold M\n" ] } ], "source": [ "EE = Ric + Lamb/3*g - 1/2* (dphi*dphi) - 1/2*exp(lamb*phi)*FF \\\n", " + 1/12*exp(lamb*phi)*F2*g -8*pi*T\n", "EE.set_name('E')\n", "print EE" ] }, { "cell_type": "code", "execution_count": 29, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "E_v,v = 2*(nu^2 + nu)*e^(2*nu*r)*f(v, r)^2 + (2*nu + 1)*e^(2*nu*r)*f(v, r)*d(f)/dr - 1/24*(mu*q^2 + 8*Lamb)*e^(2*nu*r)*f(v, r) - 1/2*(nu + 2)*e^(nu*r)*d(f)/dv + 1/2*e^(2*nu*r)*f(v, r)*d^2(f)/dr^2 - 8*pi*T00(v, r) \n", "E_v,r = -2*(nu^2 + nu)*e^(nu*r)*f(v, r) - (2*nu + 1)*e^(nu*r)*d(f)/dr + 1/24*(mu*q^2 + 8*Lamb)*e^(nu*r) - 1/2*e^(nu*r)*d^2(f)/dr^2 \n", "E_x,x = -2*(nu^2 + nu)*e^(2*nu*r)*f(v, r) - nu*e^(2*nu*r)*d(f)/dr + 1/24*(mu*q^2 + 8*Lamb)*e^(2*nu*r) \n", "E_y1,y1 = -2*(nu + 1)*e^(2*r)*f(v, r) - 1/12*(mu*q^2 - 4*Lamb)*e^(2*r) - e^(2*r)*d(f)/dr \n", "E_y2,y2 = -2*(nu + 1)*e^(2*r)*f(v, r) - 1/12*(mu*q^2 - 4*Lamb)*e^(2*r) - e^(2*r)*d(f)/dr \n", "E_r,v = -2*(nu^2 + nu)*e^(nu*r)*f(v, r) - (2*nu + 1)*e^(nu*r)*d(f)/dr + 1/24*(mu*q^2 + 8*Lamb)*e^(nu*r) - 1/2*e^(nu*r)*d^2(f)/dr^2 \n", "E_r,r = 2*(lamb^2*nu - lamb^2 - 4)/lamb^2 " ] }, "execution_count": 29, "metadata": {}, "output_type": "execute_result" } ], "source": [ "EE.display_comp(only_nonredundant=True)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We note that `EE==0` leads to 5 independent equations:" ] }, { "cell_type": "code", "execution_count": 30, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "2*(nu^2 + nu)*e^(2*nu*r)*f(v, r)^2 + (2*nu + 1)*e^(2*nu*r)*f(v, r)*d(f)/dr - 1/24*(mu*q^2 + 8*Lamb)*e^(2*nu*r)*f(v, r) - 1/2*(nu + 2)*e^(nu*r)*d(f)/dv + 1/2*e^(2*nu*r)*f(v, r)*d^2(f)/dr^2 - 8*pi*T00(v, r)" ] }, "execution_count": 30, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq1 = EE[0,0]\n", "eq1" ] }, { "cell_type": "code", "execution_count": 31, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "1/24*mu*q^2 - 2*(nu^2 + nu)*f(v, r) - (2*nu + 1)*d(f)/dr + 1/3*Lamb - 1/2*d^2(f)/dr^2" ] }, "execution_count": 31, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq2 = EE[0,4]/exp(nu*r)\n", "eq2" ] }, { "cell_type": "code", "execution_count": 32, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "1/24*mu*q^2 - 2*(nu^2 + nu)*f(v, r) - nu*d(f)/dr + 1/3*Lamb" ] }, "execution_count": 32, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq3 = EE[1,1]/exp(2*nu*r)\n", "eq3" ] }, { "cell_type": "code", "execution_count": 33, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "-1/12*mu*q^2 - 2*(nu + 1)*f(v, r) + 1/3*Lamb - d(f)/dr" ] }, "execution_count": 33, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq4 = EE[2,2]/exp(2*r)\n", "eq4" ] }, { "cell_type": "code", "execution_count": 34, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "lamb^2*nu - lamb^2 - 4" ] }, "execution_count": 34, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq5 = EE[4,4]*lamb^2/2\n", "eq5" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Dilaton field equation\n", "\n", "First we evaluate $\\nabla_m \\nabla^m \\phi$:" ] }, { "cell_type": "code", "execution_count": 35, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Levi-Civita connection nabla_g associated with the Lorentzian metric g on the 5-dimensional differentiable manifold M\n" ] }, { "data": { "text/html": [ "" ], "text/plain": [ "Levi-Civita connection nabla_g associated with the Lorentzian metric g on the 5-dimensional differentiable manifold M" ] }, "execution_count": 35, "metadata": {}, "output_type": "execute_result" } ], "source": [ "nab = g.connection()\n", "print nab\n", "nab" ] }, { "cell_type": "code", "execution_count": 36, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Scalar field on the 5-dimensional differentiable manifold M\n" ] }, { "data": { "text/html": [ "" ], "text/plain": [ "M --> R\n", "(v, x, y1, y2, r) |--> 4*(2*(nu + 1)*f(v, r) + d(f)/dr)/lamb" ] }, "execution_count": 36, "metadata": {}, "output_type": "execute_result" } ], "source": [ "box_phi = nab(nab(phi).up(g)).trace()\n", "print box_phi\n", "box_phi.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "From the action (1), the field equation for $\\phi$ is \n", "$$ \\nabla_m \\nabla^m \\phi = \\frac{\\lambda}{4} e^{\\lambda\\phi} F_{mn} F^{mn}$$\n", "We write it as\n", "\n", " DE == 0\n", " \n", "with `DE` defined by" ] }, { "cell_type": "code", "execution_count": 37, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Scalar field on the 5-dimensional differentiable manifold M\n" ] } ], "source": [ "DE = box_phi - lamb/4*exp(lamb*phi) * F2\n", "print DE" ] }, { "cell_type": "code", "execution_count": 38, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "M --> R\n", "(v, x, y1, y2, r) |--> -1/8*(lamb^2*mu*q^2 - 64*(nu + 1)*f(v, r) - 32*d(f)/dr)/lamb" ] }, "execution_count": 38, "metadata": {}, "output_type": "execute_result" } ], "source": [ "DE.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Hence the dilaton field equation provides a 6th equation:" ] }, { "cell_type": "code", "execution_count": 39, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "-lamb^2*mu*q^2 + 64*(nu + 1)*f(v, r) + 32*d(f)/dr" ] }, "execution_count": 39, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq6 = DE.coord_function()*8*lamb\n", "eq6" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Maxwell equation\n", "\n", "From the action (1), the field equation for $F$ is \n", "$$ \\nabla_m \\left( e^{\\lambda\\phi} F^{mn} \\right)= 0 $$\n", "We write it as\n", "\n", " ME == 0\n", " \n", "with `ME` defined by" ] }, { "cell_type": "code", "execution_count": 40, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Vector field on the 5-dimensional differentiable manifold M\n" ] }, { "data": { "text/html": [ "" ], "text/plain": [ "0" ] }, "execution_count": 40, "metadata": {}, "output_type": "execute_result" } ], "source": [ "ME = nab(exp(lamb*phi)*Fu).trace(0,2)\n", "print ME\n", "ME.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We get identically zero. Hence the Maxwell equation does not bring another equation to be solved." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Summary\n", "\n", "We have 6 equations to solve:" ] }, { "cell_type": "code", "execution_count": 41, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "2*(nu^2 + nu)*e^(2*nu*r)*f(v, r)^2 + (2*nu + 1)*e^(2*nu*r)*f(v, r)*d(f)/dr - 1/24*(mu*q^2 + 8*Lamb)*e^(2*nu*r)*f(v, r) - 1/2*(nu + 2)*e^(nu*r)*d(f)/dv + 1/2*e^(2*nu*r)*f(v, r)*d^2(f)/dr^2 - 8*pi*T00(v, r)" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/html": [ "" ], "text/plain": [ "1/24*mu*q^2 - 2*(nu^2 + nu)*f(v, r) - (2*nu + 1)*d(f)/dr + 1/3*Lamb - 1/2*d^2(f)/dr^2" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/html": [ "" ], "text/plain": [ "1/24*mu*q^2 - 2*(nu^2 + nu)*f(v, r) - nu*d(f)/dr + 1/3*Lamb" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/html": [ "" ], "text/plain": [ "-1/12*mu*q^2 - 2*(nu + 1)*f(v, r) + 1/3*Lamb - d(f)/dr" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/html": [ "" ], "text/plain": [ "lamb^2*nu - lamb^2 - 4" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/html": [ "" ], "text/plain": [ "-lamb^2*mu*q^2 + 64*(nu + 1)*f(v, r) + 32*d(f)/dr" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "eqs = [eq1, eq2, eq3, eq4, eq5, eq6]\n", "for eq in eqs:\n", " pretty_print(eq)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Solutions\n", "\n", "Let us show that solutions exists for $\\nu=2$ and $\\nu=4$ with the following specific form of the blackening function:\n", "\n", "$$ f(v,r) = 1 - m(v) e^{-(2\\nu +2)r} $$\n", " \n", "To this aim, we declare" ] }, { "cell_type": "code", "execution_count": 42, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "(v, r) |--> -e^(-2*(nu + 1)*r)*m(v) + 1" ] }, "execution_count": 42, "metadata": {}, "output_type": "execute_result" } ], "source": [ "fm(v,r) = 1 - function('m')(v)*exp(-(2*nu+2)*r)\n", "fm" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "and substitute this function for $f(v,r)$ in all the equations:" ] }, { "cell_type": "code", "execution_count": 43, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "-1/24*(192*pi*T00(v, r)*e^(nu*r + 2*r) - (mu*q^2 - 48*nu^2 + 8*Lamb - 48*nu)*e^(nu*r)*m(v) + (mu*q^2 - 48*nu^2 + 8*Lamb - 48*nu)*e^(3*nu*r + 2*r) - 12*(nu + 2)*D[0](m)(v))*e^(-nu*r - 2*r)" ] }, "execution_count": 43, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq1m = eq1.expr().substitute_function(f, fm).simplify_full()\n", "eq1m" ] }, { "cell_type": "code", "execution_count": 44, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "-8*pi*T00(v, r)*e^(nu*r + 2*r) + 1/24*(mu*q^2 - 48*nu^2 + 8*Lamb - 48*nu)*e^(nu*r)*m(v) - 1/24*(mu*q^2 - 48*nu^2 + 8*Lamb - 48*nu)*e^(3*nu*r + 2*r) + 1/2*(nu + 2)*D[0](m)(v)" ] }, "execution_count": 44, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq1m = (eq1m * exp(nu*r+2*r)).simplify_full()\n", "eq1m" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Note that in this expression $D[0](m)(v)$ stands for $\\mathrm{d}m/\\mathrm{d}v$." ] }, { "cell_type": "code", "execution_count": 45, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "1/24*mu*q^2 - 2*nu^2 + 1/3*Lamb - 2*nu" ] }, "execution_count": 45, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq2m = eq2.expr().substitute_function(f, fm).simplify_full()\n", "eq2m" ] }, { "cell_type": "code", "execution_count": 46, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "1/24*mu*q^2 - 2*nu^2 + 1/3*Lamb - 2*nu" ] }, "execution_count": 46, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq3m = eq3.expr().substitute_function(f, fm).simplify_full()\n", "eq3m" ] }, { "cell_type": "code", "execution_count": 47, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "-1/12*mu*q^2 + 1/3*Lamb - 2*nu - 2" ] }, "execution_count": 47, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq4m = eq4.expr().substitute_function(f, fm).simplify_full()\n", "eq4m" ] }, { "cell_type": "code", "execution_count": 48, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "lamb^2*nu - lamb^2 - 4" ] }, "execution_count": 48, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq5m = eq5.expr().substitute_function(f, fm).simplify_full()\n", "eq5m" ] }, { "cell_type": "code", "execution_count": 49, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "-lamb^2*mu*q^2 + 64*nu + 64" ] }, "execution_count": 49, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq6m = eq6.expr().substitute_function(f, fm).simplify_full()\n", "eq6m" ] }, { "cell_type": "code", "execution_count": 50, "metadata": { "collapsed": true }, "outputs": [], "source": [ "eqs = [eq1m, eq2m, eq3m, eq4m, eq5m, eq6m]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Solution for $\\nu = 2$" ] }, { "cell_type": "code", "execution_count": 51, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "[-8*pi*T00(v, r)*e^(4*r) + 1/24*(mu*q^2 + 8*Lamb - 288)*e^(2*r)*m(v) - 1/24*(mu*q^2 + 8*Lamb - 288)*e^(8*r) + 2*D[0](m)(v) == 0,\n", " 1/24*mu*q^2 + 1/3*Lamb - 12 == 0,\n", " 1/24*mu*q^2 + 1/3*Lamb - 12 == 0,\n", " -1/12*mu*q^2 + 1/3*Lamb - 6 == 0,\n", " lamb^2 - 4 == 0,\n", " -lamb^2*mu*q^2 + 192 == 0]" ] }, "execution_count": 51, "metadata": {}, "output_type": "execute_result" } ], "source": [ "neqs = [eq.subs(nu=2).simplify_full() for eq in eqs]\n", "[eq == 0 for eq in neqs]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Let us first search for a solution of all the equations but the first one:" ] }, { "cell_type": "code", "execution_count": 52, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "[1/24*mu*q^2 + 1/3*Lamb - 12 == 0,\n", " 1/24*mu*q^2 + 1/3*Lamb - 12 == 0,\n", " -1/12*mu*q^2 + 1/3*Lamb - 6 == 0,\n", " lamb^2 - 4 == 0,\n", " -lamb^2*mu*q^2 + 192 == 0]" ] }, "execution_count": 52, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eqs2 = eqs[1:] # we skip eq1m\n", "neqs = [eq.subs(nu=2).simplify_full() for eq in eqs2]\n", "[eq == 0 for eq in neqs]" ] }, { "cell_type": "code", "execution_count": 53, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "[[lamb == -2, mu == 48/r1^2, Lamb == 30, q == r1], [lamb == 2, mu == 48/r2^2, Lamb == 30, q == r2]]" ] }, "execution_count": 53, "metadata": {}, "output_type": "execute_result" } ], "source": [ "solve([eq == 0 for eq in neqs], lamb, mu, Lamb, q)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Then, we substitute the found solution in the first equation:" ] }, { "cell_type": "code", "execution_count": 54, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "-8*pi*T00(v, r)*e^(4*r) + 2*D[0](m)(v)" ] }, "execution_count": 54, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq = eq1m.subs({nu:2, Lamb:30, lamb:2, mu:48/q^2})\n", "eq" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Hence, for any choice of $m(v)$, we get a solution with \n", "$$\n", "\\Lambda=30,\\quad \\lambda = \\pm 2, \\quad\\mu = \\frac{48}{q^2} \n", "$$\n", "and\n", "$$\n", " T_{00}(\\nu,r) = \\frac{1}{4\\pi e^{4 r}} \\frac{\\mathrm{d}m}{\\mathrm{d}v}\n", "$$" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Solution for $\\nu=4$" ] }, { "cell_type": "code", "execution_count": 55, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "[-8*pi*T00(v, r)*e^(6*r) + 1/24*(mu*q^2 + 8*Lamb - 960)*e^(4*r)*m(v) - 1/24*(mu*q^2 + 8*Lamb - 960)*e^(14*r) + 3*D[0](m)(v) == 0,\n", " 1/24*mu*q^2 + 1/3*Lamb - 40 == 0,\n", " 1/24*mu*q^2 + 1/3*Lamb - 40 == 0,\n", " -1/12*mu*q^2 + 1/3*Lamb - 10 == 0,\n", " 3*lamb^2 - 4 == 0,\n", " -lamb^2*mu*q^2 + 320 == 0]" ] }, "execution_count": 55, "metadata": {}, "output_type": "execute_result" } ], "source": [ "neqs = [eq.subs(nu=4).simplify_full() for eq in eqs]\n", "[eq == 0 for eq in neqs]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Let us first search for a solution of all the equations but the first one:" ] }, { "cell_type": "code", "execution_count": 56, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "[1/24*mu*q^2 + 1/3*Lamb - 40 == 0,\n", " 1/24*mu*q^2 + 1/3*Lamb - 40 == 0,\n", " -1/12*mu*q^2 + 1/3*Lamb - 10 == 0,\n", " 3*lamb^2 - 4 == 0,\n", " -lamb^2*mu*q^2 + 320 == 0]" ] }, "execution_count": 56, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eqs2 = eqs[1:] # we skip eq1m\n", "neqs = [eq.subs(nu=4).simplify_full() for eq in eqs2]\n", "[eq == 0 for eq in neqs]" ] }, { "cell_type": "code", "execution_count": 57, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "[[lamb == 2/3*sqrt(3), mu == 240/r3^2, Lamb == 90, q == r3], [lamb == -2/3*sqrt(3), mu == 240/r4^2, Lamb == 90, q == r4]]" ] }, "execution_count": 57, "metadata": {}, "output_type": "execute_result" } ], "source": [ "solve([eq == 0 for eq in neqs], lamb, mu, Lamb, q)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Then, we substitute the found solution in the first equation:" ] }, { "cell_type": "code", "execution_count": 58, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "-8*pi*T00(v, r)*e^(6*r) + 3*D[0](m)(v)" ] }, "execution_count": 58, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq = eq1m.subs({nu:4, Lamb:90, lamb:2/sqrt(3), mu:240/q^2})\n", "eq" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Hence, for any choice of $m(v)$, we get a solution with \n", "$$\n", "\\Lambda=90,\\quad \\lambda = \\pm \\frac{2}{\\sqrt{3}}, \\quad\\mu = \\frac{240}{q^2} \n", "$$\n", "and\n", "$$\n", " T_{00}(\\nu,r) = \\frac{3}{8\\pi e^{6 r}} \\frac{\\mathrm{d}m}{\\mathrm{d}v}\n", "$$" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [] } ], "metadata": { "kernelspec": { "display_name": "SageMath 6.10", "language": "", "name": "sagemath" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 2 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython2", "version": "2.7.10" }, "name": "Vaidya-Lifshitz.ipynb" }, "nbformat": 4, "nbformat_minor": 0 }