{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "## Black branes in Lifshitz-like spacetimes" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "This Jupyter/SageMath 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://github.com/sagemanifolds/SageManifolds/raw/master/Worksheets/v0.9/SM_black_brane.ipynb)." ] }, { "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": true }, "outputs": [], "source": [ "%display latex" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Space manifold\n", "\n", "Let us declare $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 a coordinate system on $M$:" ] }, { "cell_type": "code", "execution_count": 3, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "Chart (M, (t, x, y1, y2, r))" ] }, "execution_count": 3, "metadata": {}, "output_type": "execute_result" } ], "source": [ "X. = M.chart('t x y1:y_1 y2:y_2 r')\n", "X" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Next, we define the metric tensor, which depends on some real number $\\nu$ and some arbitary function $f$:" ] }, { "cell_type": "code", "execution_count": 4, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "g = -e^(2*nu*r)*f(r) dt*dt + e^(2*nu*r) dx*dx + e^(2*r) dy1*dy1 + e^(2*r) dy2*dy2 + 1/f(r) dr*dr" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" } ], "source": [ "g = M.lorentzian_metric('g')\n", "var('nu', latex_name=r'\\nu', domain='real')\n", "ff = function('f')(r)\n", "g[0,0] = -ff*exp(2*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[4,4] = 1/ff\n", "g.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "If $f(r)=1$, this is the metric of a Lifshitz spacetime; if, in addition $\\nu=1$, $(M,g)$ is a Poincaré patch of $\\mathrm{AdS}_5$." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Curvature\n", "\n", "The Riemann tensor is" ] }, { "cell_type": "code", "execution_count": 5, "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": "code", "execution_count": 6, "metadata": { "collapsed": false, "scrolled": true }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "Riem(g)^t_x,t,x = -nu^2*e^(2*nu*r)*f(r) - 1/2*nu*e^(2*nu*r)*d(f)/dr \n", "Riem(g)^t_y1,t,y1 = -nu*e^(2*r)*f(r) - 1/2*e^(2*r)*d(f)/dr \n", "Riem(g)^t_y2,t,y2 = -nu*e^(2*r)*f(r) - 1/2*e^(2*r)*d(f)/dr \n", "Riem(g)^t_r,t,r = -1/2*(2*nu^2*f(r) + 3*nu*d(f)/dr + d^2(f)/dr^2)/f(r) \n", "Riem(g)^x_t,t,x = -nu^2*e^(2*nu*r)*f(r)^2 - 1/2*nu*e^(2*nu*r)*f(r)*d(f)/dr \n", "Riem(g)^x_y1,x,y1 = -nu*e^(2*r)*f(r) \n", "Riem(g)^x_y2,x,y2 = -nu*e^(2*r)*f(r) \n", "Riem(g)^x_r,x,r = -1/2*(2*nu^2*f(r) + nu*d(f)/dr)/f(r) \n", "Riem(g)^y1_t,t,y1 = -nu*e^(2*nu*r)*f(r)^2 - 1/2*e^(2*nu*r)*f(r)*d(f)/dr \n", "Riem(g)^y1_x,x,y1 = nu*e^(2*nu*r)*f(r) \n", "Riem(g)^y1_y2,y1,y2 = -e^(2*r)*f(r) \n", "Riem(g)^y1_r,y1,r = -1/2*(2*f(r) + d(f)/dr)/f(r) \n", "Riem(g)^y2_t,t,y2 = -nu*e^(2*nu*r)*f(r)^2 - 1/2*e^(2*nu*r)*f(r)*d(f)/dr \n", "Riem(g)^y2_x,x,y2 = nu*e^(2*nu*r)*f(r) \n", "Riem(g)^y2_y1,y1,y2 = e^(2*r)*f(r) \n", "Riem(g)^y2_r,y2,r = -1/2*(2*f(r) + d(f)/dr)/f(r) \n", "Riem(g)^r_t,t,r = -nu^2*e^(2*nu*r)*f(r)^2 - 3/2*nu*e^(2*nu*r)*f(r)*d(f)/dr - 1/2*e^(2*nu*r)*f(r)*d^2(f)/dr^2 \n", "Riem(g)^r_x,x,r = nu^2*e^(2*nu*r)*f(r) + 1/2*nu*e^(2*nu*r)*d(f)/dr \n", "Riem(g)^r_y1,y1,r = e^(2*r)*f(r) + 1/2*e^(2*r)*d(f)/dr \n", "Riem(g)^r_y2,y2,r = e^(2*r)*f(r) + 1/2*e^(2*r)*d(f)/dr " ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Riem.display_comp(only_nonredundant=True)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The Ricci tensor:" ] }, { "cell_type": "code", "execution_count": 7, "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": 8, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "Ric(g) = (2*(nu^2 + nu)*e^(2*nu*r)*f(r)^2 + (2*nu + 1)*e^(2*nu*r)*f(r)*d(f)/dr + 1/2*e^(2*nu*r)*f(r)*d^2(f)/dr^2) dt*dt + (-2*(nu^2 + nu)*e^(2*nu*r)*f(r) - nu*e^(2*nu*r)*d(f)/dr) dx*dx + (-2*(nu + 1)*e^(2*r)*f(r) - e^(2*r)*d(f)/dr) dy1*dy1 + (-2*(nu + 1)*e^(2*r)*f(r) - e^(2*r)*d(f)/dr) dy2*dy2 - 1/2*(4*(nu^2 + 1)*f(r) + 2*(2*nu + 1)*d(f)/dr + d^2(f)/dr^2)/f(r) dr*dr" ] }, "execution_count": 8, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Ric.display()" ] }, { "cell_type": "code", "execution_count": 9, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "Ric(g)_t,t = 2*(nu^2 + nu)*e^(2*nu*r)*f(r)^2 + (2*nu + 1)*e^(2*nu*r)*f(r)*d(f)/dr + 1/2*e^(2*nu*r)*f(r)*d^2(f)/dr^2 \n", "Ric(g)_x,x = -2*(nu^2 + nu)*e^(2*nu*r)*f(r) - nu*e^(2*nu*r)*d(f)/dr \n", "Ric(g)_y1,y1 = -2*(nu + 1)*e^(2*r)*f(r) - e^(2*r)*d(f)/dr \n", "Ric(g)_y2,y2 = -2*(nu + 1)*e^(2*r)*f(r) - e^(2*r)*d(f)/dr \n", "Ric(g)_r,r = -1/2*(4*(nu^2 + 1)*f(r) + 2*(2*nu + 1)*d(f)/dr + d^2(f)/dr^2)/f(r) " ] }, "execution_count": 9, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Ric.display_comp()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The Ricci scalar:" ] }, { "cell_type": "code", "execution_count": 10, "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": 11, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "r(g): M --> R\n", " (t, x, y1, y2, r) |--> -2*(3*nu^2 + 4*nu + 3)*f(r) - (5*nu + 4)*d(f)/dr - d^2(f)/dr^2" ] }, "execution_count": 11, "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 cosmological constant $\\bar{\\Lambda} = -\\Lambda/2$ with $\\Lambda>0$, 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} \\right) \\sqrt{-g} \\, \\mathrm{d}^5 x \\qquad\\qquad \\mbox{(1)}$$\n", "\n", "where $R(g)$ is the Ricci scalar of metric $g$ and $\\lambda$ is the dilatonic coupling constant." ] }, { "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) \\iff e^{\\lambda\\phi} = \\mu e^{4r},$$\n", "where $\\mu$ is a constant. " ] }, { "cell_type": "code", "execution_count": 12, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "phi: M --> R\n", " (t, x, y1, y2, r) |--> (4*r + log(mu))/lamb" ] }, "execution_count": 12, "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": 13, "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": 14, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "dphi = 4/lamb dr" ] }, "execution_count": 14, "metadata": {}, "output_type": "execute_result" } ], "source": [ "dphi.display()" ] }, { "cell_type": "code", "execution_count": 15, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "[0, 0, 0, 0, 4/lamb]" ] }, "execution_count": 15, "metadata": {}, "output_type": "execute_result" } ], "source": [ "dphi[:] # all the components in the default frame" ] }, { "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", "Let us first get the 1-forms $\\mathrm{d}y_1$ and $\\mathrm{d}y_2$:" ] }, { "cell_type": "code", "execution_count": 16, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "Coordinate coframe (M, (dt,dx,dy1,dy2,dr))" ] }, "execution_count": 16, "metadata": {}, "output_type": "execute_result" } ], "source": [ "X.coframe()" ] }, { "cell_type": "code", "execution_count": 17, "metadata": { "collapsed": false }, "outputs": [], "source": [ "dy1 = X.coframe()[2]\n", "dy2 = X.coframe()[3]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Then we can form $F$ according to the above ansatz:" ] }, { "cell_type": "code", "execution_count": 18, "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": 18, "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": 19, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "3-form dF on the 5-dimensional differentiable manifold M\n" ] } ], "source": [ "print(F.exterior_der())" ] }, { "cell_type": "code", "execution_count": 20, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "dF = 0" ] }, "execution_count": 20, "metadata": {}, "output_type": "execute_result" } ], "source": [ "F.exterior_der().display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Let us evaluate the square $F_{mn} F^{mn}$ of $F$:" ] }, { "cell_type": "code", "execution_count": 21, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Tensor field of type (2,0) on the 5-dimensional differentiable manifold M\n" ] } ], "source": [ "Fu = F.up(g)\n", "print(Fu)" ] }, { "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", "(t, x, y1, y2, r) |--> 1/2*q^2*e^(-4*r)" ] }, "execution_count": 22, "metadata": {}, "output_type": "execute_result" } ], "source": [ "F2 = F['_{mn}']*Fu['^{mn}'] # using LaTeX notations for 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": [ "### Einstein equation\n", "\n", "Let us first introduce the cosmological constant:" ] }, { "cell_type": "code", "execution_count": 26, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "Lamb" ] }, "execution_count": 26, "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} = 0 $$\n", "We write it as\n", "\n", " EE == 0\n", "\n", "with `EE` defined by" ] }, { "cell_type": "code", "execution_count": 27, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Field of symmetric bilinear forms E 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\n", "EE.set_name('E')\n", "print(EE)" ] }, { "cell_type": "code", "execution_count": 28, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "E_t,t = 2*(nu^2 + nu)*e^(2*nu*r)*f(r)^2 + (2*nu + 1)*e^(2*nu*r)*f(r)*d(f)/dr - 1/24*(mu*q^2 + 8*Lamb)*e^(2*nu*r)*f(r) + 1/2*e^(2*nu*r)*f(r)*d^2(f)/dr^2 \n", "E_x,x = -2*(nu^2 + nu)*e^(2*nu*r)*f(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(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(r) - 1/12*(mu*q^2 - 4*Lamb)*e^(2*r) - e^(2*r)*d(f)/dr \n", "E_r,r = 1/24*(lamb^2*mu*q^2 + 8*Lamb*lamb^2 - 12*lamb^2*d^2(f)/dr^2 - 48*(lamb^2*nu^2 + lamb^2 + 4)*f(r) - 24*(2*lamb^2*nu + lamb^2)*d(f)/dr)/(lamb^2*f(r)) " ] }, "execution_count": 28, "metadata": {}, "output_type": "execute_result" } ], "source": [ "EE.display_comp(only_nonredundant=True)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We note that `EE==0` leads to only 4 independent equations:" ] }, { "cell_type": "code", "execution_count": 29, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "2*(nu^2 + nu)*f(r)^2 + (2*nu + 1)*f(r)*d(f)/dr - 1/24*(mu*q^2 + 8*Lamb)*f(r) + 1/2*f(r)*d^2(f)/dr^2" ] }, "execution_count": 29, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq0 = EE[0,0]/exp(2*nu*r)\n", "eq0" ] }, { "cell_type": "code", "execution_count": 30, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "1/24*mu*q^2 - 2*(nu^2 + nu)*f(r) - nu*d(f)/dr + 1/3*Lamb" ] }, "execution_count": 30, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq1 = EE[1,1]/exp(2*nu*r)\n", "eq1" ] }, { "cell_type": "code", "execution_count": 31, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "-1/12*mu*q^2 - 2*(nu + 1)*f(r) + 1/3*Lamb - d(f)/dr" ] }, "execution_count": 31, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq2 = EE[2,2]/exp(2*r)\n", "eq2" ] }, { "cell_type": "code", "execution_count": 32, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "1/24*lamb^2*mu*q^2 + 1/3*Lamb*lamb^2 - 1/2*lamb^2*d^2(f)/dr^2 - 2*(lamb^2*nu^2 + lamb^2 + 4)*f(r) - (2*lamb^2*nu + lamb^2)*d(f)/dr" ] }, "execution_count": 32, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq3 = EE[4,4]*lamb^2*f(r)\n", "eq3" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Dilaton field equation\n", "\n", "First we evaluate $\\nabla_m \\nabla^m \\phi$:" ] }, { "cell_type": "code", "execution_count": 33, "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" ] } ], "source": [ "nab = g.connection()\n", "print(nab)" ] }, { "cell_type": "code", "execution_count": 34, "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", "(t, x, y1, y2, r) |--> 4*(2*(nu + 1)*f(r) + d(f)/dr)/lamb" ] }, "execution_count": 34, "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": 35, "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", "(t, x, y1, y2, r) |--> -1/8*(lamb^2*mu*q^2 - 64*(nu + 1)*f(r) - 32*d(f)/dr)/lamb" ] }, "execution_count": 35, "metadata": {}, "output_type": "execute_result" } ], "source": [ "DE = box_phi - lamb/4*exp(lamb*phi) * F2\n", "print(DE)\n", "DE.display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Hence the dilaton field equation provides a fourth equation:" ] }, { "cell_type": "code", "execution_count": 36, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "-1/8*lamb^2*mu*q^2 + 8*(nu + 1)*f(r) + 4*d(f)/dr" ] }, "execution_count": 36, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq4 = DE.coord_function()*lamb\n", "eq4" ] }, { "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": 37, "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": 37, "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; indeed the tensor $\\nabla_p (e^{\\lambda\\phi} F^{mn})$ has a vanishing trace, as we can check:" ] }, { "cell_type": "code", "execution_count": 38, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "mu*q d/dy1*d/dy2*dr - 1/2*mu*q*e^(2*r)*f(r) d/dy1*d/dr*dy2 - mu*q d/dy2*d/dy1*dr + 1/2*mu*q*e^(2*r)*f(r) d/dy2*d/dr*dy1 + 1/2*mu*q*e^(2*r)*f(r) d/dr*d/dy1*dy2 - 1/2*mu*q*e^(2*r)*f(r) d/dr*d/dy2*dy1" ] }, "execution_count": 38, "metadata": {}, "output_type": "execute_result" } ], "source": [ "nab(exp(lamb*phi)*Fu).display()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Solving the field equations\n", "\n", "The system to solve is" ] }, { "cell_type": "code", "execution_count": 39, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "2*(nu^2 + nu)*f(r)^2 + (2*nu + 1)*f(r)*d(f)/dr - 1/24*(mu*q^2 + 8*Lamb)*f(r) + 1/2*f(r)*d^2(f)/dr^2 ' = 0'" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/html": [ "" ], "text/plain": [ "1/24*mu*q^2 - 2*(nu^2 + nu)*f(r) - nu*d(f)/dr + 1/3*Lamb ' = 0'" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/html": [ "" ], "text/plain": [ "-1/12*mu*q^2 - 2*(nu + 1)*f(r) + 1/3*Lamb - d(f)/dr ' = 0'" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/html": [ "" ], "text/plain": [ "1/24*lamb^2*mu*q^2 + 1/3*Lamb*lamb^2 - 1/2*lamb^2*d^2(f)/dr^2 - 2*(lamb^2*nu^2 + lamb^2 + 4)*f(r) - (2*lamb^2*nu + lamb^2)*d(f)/dr ' = 0'" ] }, "metadata": {}, "output_type": "display_data" }, { "data": { "text/html": [ "" ], "text/plain": [ "-1/8*lamb^2*mu*q^2 + 8*(nu + 1)*f(r) + 4*d(f)/dr ' = 0'" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "eqs = [eq0, eq1, eq2, eq3, eq4]\n", "for eq in eqs:\n", " pretty_print(eq, ' = 0')" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Let us solve `eq1` for $f(r)$:" ] }, { "cell_type": "code", "execution_count": 40, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "_C*e^(-2*(nu + 1)*r) + 1/48*mu*q^2/((nu + 1)*nu) + 1/6*Lamb/((nu + 1)*nu)" ] }, "execution_count": 40, "metadata": {}, "output_type": "execute_result" } ], "source": [ "sol_f = desolve(eq1.expr() == 0, f(r), ivar=r)\n", "sol_f.expand()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Hence, up to some rescaling the solution is of the type\n", "\n", "$$ f(r) = 1 - m e^{-(2\\nu +2)r}, $$\n", "\n", "where $m$ is a constant. Hence we declare" ] }, { "cell_type": "code", "execution_count": 41, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "r |--> -m*e^(-2*(nu + 1)*r) + 1" ] }, "execution_count": 41, "metadata": {}, "output_type": "execute_result" } ], "source": [ "var('m', domain='real')\n", "fm(r) = 1 - m*exp(-(2*nu+2)*r)\n", "fm" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "and substitute this function for $f(r)$ in all the equations:" ] }, { "cell_type": "code", "execution_count": 42, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "1/24*(m*mu*q^2 - 48*m*nu^2 + 8*Lamb*m - 48*m*nu - (mu*q^2 - 48*nu^2 + 8*Lamb - 48*nu)*e^(2*nu*r + 2*r))*e^(-2*nu*r - 2*r)" ] }, "execution_count": 42, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq0m = eq0.expr().substitute_function(f, fm).simplify_full()\n", "eq0m" ] }, { "cell_type": "code", "execution_count": 43, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "1/24*m*mu*q^2 - 2*m*nu^2 + 1/3*Lamb*m - 2*m*nu - 1/24*(mu*q^2 - 48*nu^2 + 8*Lamb - 48*nu)*e^(2*nu*r + 2*r)" ] }, "execution_count": 43, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq0m = (eq0m * exp(2*nu*r+2*r)).simplify_full()\n", "eq0m" ] }, { "cell_type": "code", "execution_count": 44, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "1/24*mu*q^2 - 2*nu^2 + 1/3*Lamb - 2*nu" ] }, "execution_count": 44, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq1m = eq1.expr().substitute_function(f, fm).simplify_full()\n", "eq1m" ] }, { "cell_type": "code", "execution_count": 45, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "-1/12*mu*q^2 + 1/3*Lamb - 2*nu - 2" ] }, "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*(48*lamb^2*m*nu - 48*(lamb^2 + 4)*m - (lamb^2*mu*q^2 - 48*lamb^2*nu^2 + 8*(Lamb - 6)*lamb^2 - 192)*e^(2*nu*r + 2*r))*e^(-2*nu*r - 2*r)" ] }, "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": [ "-2*lamb^2*m*nu + 2*(lamb^2 + 4)*m + 1/24*(lamb^2*mu*q^2 - 48*lamb^2*nu^2 + 8*(Lamb - 6)*lamb^2 - 192)*e^(2*nu*r + 2*r)" ] }, "execution_count": 47, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq3m = (eq3m * exp(2*nu*r+2*r)).simplify_full()\n", "eq3m" ] }, { "cell_type": "code", "execution_count": 48, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "-1/8*lamb^2*mu*q^2 + 8*nu + 8" ] }, "execution_count": 48, "metadata": {}, "output_type": "execute_result" } ], "source": [ "eq4m = eq4.expr().substitute_function(f, fm).simplify_full()\n", "eq4m" ] }, { "cell_type": "code", "execution_count": 49, "metadata": { "collapsed": false }, "outputs": [], "source": [ "eqs = [eq0m, eq1m, eq2m, eq3m, eq4m]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Solution for $\\nu = 2$" ] }, { "cell_type": "code", "execution_count": 50, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "[1/24*m*mu*q^2 + 1/3*(Lamb - 36)*m - 1/24*(mu*q^2 + 8*Lamb - 288)*e^(6*r) == 0,\n", " 1/24*mu*q^2 + 1/3*Lamb - 12 == 0,\n", " -1/12*mu*q^2 + 1/3*Lamb - 6 == 0,\n", " -2*(lamb^2 - 4)*m + 1/24*(lamb^2*mu*q^2 + 8*(Lamb - 30)*lamb^2 - 192)*e^(6*r) == 0,\n", " -1/8*lamb^2*mu*q^2 + 24 == 0]" ] }, "execution_count": 50, "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": "code", "execution_count": 51, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "[[lamb == 2, mu == 48/r1^2, Lamb == 30, q == r1, m == r2, r == r3], [lamb == -2, mu == 48/r4^2, Lamb == 30, q == r4, m == r5, r == r6]]" ] }, "execution_count": 51, "metadata": {}, "output_type": "execute_result" } ], "source": [ "solve([eq == 0 for eq in neqs], lamb, mu, Lamb, q, m, r)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "In the above solutions, $r_i$, with $i$ an integer, stands for an arbitrary parameter. In particular, we notice that $\\mu$ and $q$ are related by $\\mu q^2 = 48$ and that the value of $m$ can be chosen arbitrarily." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Solution for $\\nu=4$" ] }, { "cell_type": "code", "execution_count": 52, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "[1/24*m*mu*q^2 + 1/3*(Lamb - 120)*m - 1/24*(mu*q^2 + 8*Lamb - 960)*e^(10*r) == 0,\n", " 1/24*mu*q^2 + 1/3*Lamb - 40 == 0,\n", " -1/12*mu*q^2 + 1/3*Lamb - 10 == 0,\n", " -2*(3*lamb^2 - 4)*m + 1/24*(lamb^2*mu*q^2 + 8*(Lamb - 102)*lamb^2 - 192)*e^(10*r) == 0,\n", " -1/8*lamb^2*mu*q^2 + 40 == 0]" ] }, "execution_count": 52, "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": "code", "execution_count": 53, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/html": [ "" ], "text/plain": [ "[[lamb == 2/3*sqrt(3), mu == 240/r7^2, Lamb == 90, q == r7, m == r8, r == r9], [lamb == -2/3*sqrt(3), mu == 240/r10^2, Lamb == 90, q == r10, m == r11, r == r12]]" ] }, "execution_count": 53, "metadata": {}, "output_type": "execute_result" } ], "source": [ "solve([eq == 0 for eq in neqs], lamb, mu, Lamb, q, m, r)" ] }, { "cell_type": "markdown", "metadata": { "collapsed": true }, "source": [ "As above, $r_i$, with $i$ an integer, stands for an arbitrary parameter. The constants $\\mu$ and $q$ are now related by $\\mu q^2 = 240$ and the value of $m$ is still arbitrary." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [] } ], "metadata": { "kernelspec": { "display_name": "Sage 6.9", "language": "", "name": "sage_6_9" }, "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.9" }, "name": "Lifshitz_black_brane.ipynb" }, "nbformat": 4, "nbformat_minor": 0 }