{
"cells": [
{
"cell_type": "markdown",
"metadata": {},
"source": [
"# Solving Einstein equation to get Kottler solution\n",
"\n",
"This Jupyter/SageMath worksheet is relative to the lectures\n",
"[Geometry and physics of black holes](http://luth.obspm.fr/~luthier/gourgoulhon/bh16/)\n",
" \n",
"These computations are based on [SageManifolds](http://sagemanifolds.obspm.fr) (v0.9)\n",
"\n",
"Click [here](https://raw.githubusercontent.com/egourgoulhon/BHLectures/master/sage/Kottler_solution.ipynb) to download the worksheet file (ipynb format). To run it, you must start SageMath with the Jupyter notebook, with the command `sage -n jupyter`"
]
},
{
"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": [
"## Spacetime\n",
"\n",
"We declare the spacetime manifold $M$:"
]
},
{
"cell_type": "code",
"execution_count": 2,
"metadata": {
"collapsed": false
},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"4-dimensional differentiable manifold M\n"
]
}
],
"source": [
"M = Manifold(4, 'M')\n",
"print(M)"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"We declare the chart of spherical coordinates $(t,r,\\theta,\\phi)$:"
]
},
{
"cell_type": "code",
"execution_count": 3,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"Chart (M, (t, r, th, ph))"
]
},
"execution_count": 3,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"X. = M.chart(r't r:(0,+oo) th:(0,pi):\\theta ph:(0,2*pi):\\phi')\n",
"X"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"The static and spherically symmetric metric ansatz, with the unknown functions $A(r)$ and $B(r)$:"
]
},
{
"cell_type": "code",
"execution_count": 4,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"g = -A(r) dt*dt + B(r) dr*dr + r^2 dth*dth + r^2*sin(th)^2 dph*dph"
]
},
"execution_count": 4,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"g = M.lorentzian_metric('g')\n",
"A = function('A')\n",
"B = function('B')\n",
"g[0,0] = -A(r)\n",
"g[1,1] = B(r)\n",
"g[2,2] = r^2\n",
"g[3,3] = (r*sin(th))^2\n",
"g.display()"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"The Christoffel symbols of $g$, with respect to the default chart:"
]
},
{
"cell_type": "code",
"execution_count": 5,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"Gam^t_t,r = 1/2*d(A)/dr/A(r) \n",
"Gam^r_t,t = 1/2*d(A)/dr/B(r) \n",
"Gam^r_r,r = 1/2*d(B)/dr/B(r) \n",
"Gam^r_th,th = -r/B(r) \n",
"Gam^r_ph,ph = -r*sin(th)^2/B(r) \n",
"Gam^th_r,th = 1/r \n",
"Gam^th_ph,ph = -cos(th)*sin(th) \n",
"Gam^ph_r,ph = 1/r \n",
"Gam^ph_th,ph = cos(th)/sin(th) "
]
},
"execution_count": 5,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"g.christoffel_symbols_display()"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## Einstein equation"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"The cosmological constant:"
]
},
{
"cell_type": "code",
"execution_count": 6,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"Lamb"
]
},
"execution_count": 6,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"var('Lamb', latex_name='\\Lambda')"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"The Einstein equation:"
]
},
{
"cell_type": "code",
"execution_count": 7,
"metadata": {
"collapsed": false
},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"Field of symmetric bilinear forms E on the 4-dimensional differentiable manifold M\n"
]
}
],
"source": [
"EE = g.ricci() - 1/2*g.ricci_scalar()*g + Lamb*g\n",
"EE.set_name('E')\n",
"print(EE)"
]
},
{
"cell_type": "code",
"execution_count": 8,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"E_t,t = -((Lamb*r^2 - 1)*A(r)*B(r)^2 - r*A(r)*d(B)/dr + A(r)*B(r))/(r^2*B(r)^2) \n",
"E_r,r = ((Lamb*r^2 - 1)*A(r)*B(r) + r*d(A)/dr + A(r))/(r^2*A(r)) \n",
"E_th,th = 1/4*(4*Lamb*r^2*A(r)^2*B(r)^2 - r^2*B(r)*(d(A)/dr)^2 + 2*r^2*A(r)*B(r)*d^2(A)/dr^2 + 2*r*A(r)*B(r)*d(A)/dr - (r^2*A(r)*d(A)/dr + 2*r*A(r)^2)*d(B)/dr)/(A(r)^2*B(r)^2) \n",
"E_ph,ph = 1/4*(4*Lamb*r^2*A(r)^2*B(r)^2 - r^2*B(r)*(d(A)/dr)^2 + 2*r^2*A(r)*B(r)*d^2(A)/dr^2 + 2*r*A(r)*B(r)*d(A)/dr - (r^2*A(r)*d(A)/dr + 2*r*A(r)^2)*d(B)/dr)*sin(th)^2/(A(r)^2*B(r)^2) "
]
},
"execution_count": 8,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"EE.display_comp()"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"### Simplifying and rearranging the equations"
]
},
{
"cell_type": "code",
"execution_count": 9,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"-(Lamb*r^2 - 1)*B(r)^2 + r*d(B)/dr - B(r)"
]
},
"execution_count": 9,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"eq0 = EE[0,0]*r^2*B(r)^2/A(r); eq0"
]
},
{
"cell_type": "code",
"execution_count": 10,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"(Lamb*r^2 - 1)*A(r)*B(r) + r*d(A)/dr + A(r)"
]
},
"execution_count": 10,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"eq1 = EE[1,1]*r^2*A(r); eq1"
]
},
{
"cell_type": "code",
"execution_count": 11,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"4*Lamb*r^2*A(r)^2*B(r)^2 - r^2*B(r)*(d(A)/dr)^2 + 2*r^2*A(r)*B(r)*d^2(A)/dr^2 + 2*r*A(r)*B(r)*d(A)/dr - (r^2*A(r)*d(A)/dr + 2*r*A(r)^2)*d(B)/dr"
]
},
"execution_count": 11,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"eq2 = EE[2,2]*4*A(r)^2*B(r)^2; eq2"
]
},
{
"cell_type": "code",
"execution_count": 12,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"4*Lamb*r^2*A(r)^2*B(r)^2 - r^2*B(r)*(d(A)/dr)^2 + 2*r^2*A(r)*B(r)*d^2(A)/dr^2 + 2*r*A(r)*B(r)*d(A)/dr - (r^2*A(r)*d(A)/dr + 2*r*A(r)^2)*d(B)/dr"
]
},
"execution_count": 12,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"eq3 = EE[3,3]*4*A(r)^2*B(r)^2/sin(th)^2; eq3"
]
},
{
"cell_type": "code",
"execution_count": 13,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"True"
]
},
"execution_count": 13,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"eq3 == eq2"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"## Solving Einstein equation\n",
"\n",
"The following combination of eq0 and eq1 is particularly simple:"
]
},
{
"cell_type": "code",
"execution_count": 14,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"B(r)*d(A)/dr + A(r)*d(B)/dr"
]
},
"execution_count": 14,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"eq4 = (eq0*A(r) + eq1*B(r))/r; eq4"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"The solution is $A(r)B(r)=C$, where $C$ is a constant:"
]
},
{
"cell_type": "code",
"execution_count": 15,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"_C/A(r)"
]
},
"execution_count": 15,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"s = desolve(eq4.expr() == 0, B(r), ivar=r)\n",
"s"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"Let us rename the constant to $\\alpha$:"
]
},
{
"cell_type": "code",
"execution_count": 16,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"r |--> alpha/A(r)"
]
},
"execution_count": 16,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"var('alpha', latex_name=r'\\alpha')\n",
"B_sol(r) = s.subs(_C=alpha); B_sol"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"We replace $B(r)$ by the above value in the remaining equations:"
]
},
{
"cell_type": "code",
"execution_count": 17,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"(Lamb*r^2 - 1)*alpha + r*d(A)/dr + A(r)"
]
},
"execution_count": 17,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"eq5 = X.function(eq1.expr().substitute_function(B, B_sol)); eq5"
]
},
{
"cell_type": "code",
"execution_count": 18,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"4*Lamb*alpha^2*r^2 - alpha*r^2*(d(A)/dr)^2/A(r) + 2*alpha*r^2*d^2(A)/dr^2 + 2*alpha*r*d(A)/dr + (r^2*A(r)*d(A)/dr + 2*r*A(r)^2)*alpha*d(A)/dr/A(r)^2"
]
},
"execution_count": 18,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"eq6 = X.function(eq2.expr().substitute_function(B, B_sol)); eq6"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"Let us solve eq5 for $A(r)$. Note that we are using `eq5.expr()` to get a symbolic expression, as expected by the function `desolve`, while `eq5` is a coordinate function. "
]
},
{
"cell_type": "code",
"execution_count": 19,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"-1/3*Lamb*alpha*r^2 + alpha + _C/r"
]
},
"execution_count": 19,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"s = desolve(eq5.expr() == 0, A(r), ivar=r)\n",
"s.expand()"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"We rename the constant $C$ to $-2m$ and set the value of constant $\\alpha$ to $1$:"
]
},
{
"cell_type": "code",
"execution_count": 20,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"r |--> -1/3*Lamb*r^2 - 2*m/r + 1"
]
},
"execution_count": 20,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"var('m')\n",
"A_sol(r) = s.subs(_C=-2*m, alpha=1).expand()\n",
"A_sol"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"Let us check that `eq6` is fulfilled by the found value of $A(r)$:"
]
},
{
"cell_type": "code",
"execution_count": 21,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"0"
]
},
"execution_count": 21,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"eq6.expr().substitute_function(A, A_sol).subs(alpha=1).simplify_full()"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"### Final expression of the metric\n",
"\n",
"We have got the Kottler metric:"
]
},
{
"cell_type": "code",
"execution_count": 22,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"g = 1/3*(Lamb*r^3 + 6*m - 3*r)/r dt*dt - 3*r/(Lamb*r^3 + 6*m - 3*r) dr*dr + r^2 dth*dth + r^2*sin(th)^2 dph*dph"
]
},
"execution_count": 22,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"g[0,0] = -A_sol(r)\n",
"g[1,1] = 1/A_sol(r)\n",
"g.display()"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"which reduces to Schwarzschild metric as soons as the cosmological constant vanishes."
]
},
{
"cell_type": "code",
"execution_count": 23,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"Gam^t_t,r = (Lamb*r^3 - 3*m)/(Lamb*r^4 + 6*m*r - 3*r^2) \n",
"Gam^r_t,t = 1/9*(Lamb^2*r^6 + 3*Lamb*m*r^3 - 3*Lamb*r^4 - 18*m^2 + 9*m*r)/r^3 \n",
"Gam^r_r,r = -(Lamb*r^3 - 3*m)/(Lamb*r^4 + 6*m*r - 3*r^2) \n",
"Gam^r_th,th = 1/3*Lamb*r^3 + 2*m - r \n",
"Gam^r_ph,ph = 1/3*(Lamb*r^3 + 6*m - 3*r)*sin(th)^2 \n",
"Gam^th_r,th = 1/r \n",
"Gam^th_ph,ph = -cos(th)*sin(th) \n",
"Gam^ph_r,ph = 1/r \n",
"Gam^ph_th,ph = cos(th)/sin(th) "
]
},
"execution_count": 23,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"g.christoffel_symbols_display()"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"Let us check that Einstein equation is satisfied by the above metric:"
]
},
{
"cell_type": "code",
"execution_count": 24,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"E = 0"
]
},
"execution_count": 24,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"EE = g.ricci() - 1/2*g.ricci_scalar()*g + Lamb*g\n",
"EE.set_name('E')\n",
"EE.display()"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"The Ricci scalar is constant for this solution:"
]
},
{
"cell_type": "code",
"execution_count": 25,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"r(g): M --> R\n",
" (t, r, th, ph) |--> 4*Lamb"
]
},
"execution_count": 25,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"g.ricci_scalar().display()"
]
},
{
"cell_type": "markdown",
"metadata": {},
"source": [
"The Ricci tensor is proportional to the metric tensor:"
]
},
{
"cell_type": "code",
"execution_count": 26,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"Ric(g) = 1/3*(Lamb^2*r^3 + 6*Lamb*m - 3*Lamb*r)/r dt*dt - 3*Lamb*r/(Lamb*r^3 + 6*m - 3*r) dr*dr + Lamb*r^2 dth*dth + Lamb*r^2*sin(th)^2 dph*dph"
]
},
"execution_count": 26,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"g.ricci().display()"
]
},
{
"cell_type": "code",
"execution_count": 27,
"metadata": {
"collapsed": false
},
"outputs": [
{
"data": {
"text/html": [
""
],
"text/plain": [
"True"
]
},
"execution_count": 27,
"metadata": {},
"output_type": "execute_result"
}
],
"source": [
"g.ricci() == Lamb * g"
]
},
{
"cell_type": "code",
"execution_count": null,
"metadata": {
"collapsed": true
},
"outputs": [],
"source": []
}
],
"metadata": {
"kernelspec": {
"display_name": "SageMath 7.1",
"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"
}
},
"nbformat": 4,
"nbformat_minor": 0
}