1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
|
C Kerr metric in cartesian Boyer-Lindquist coordinates,
C as per MTW box 33.2.
C
C Author : D. Vulcanov (Timisoara, Romania)
C see ../../README for copyright & licensing info
C
C $Header$
#include "cctk.h"
#include "cctk_Parameters.h"
subroutine Exact__Kerr_BoyerLindquist(
$ x, y, z, t,
$ gdtt, gdtx, gdty, gdtz,
$ gdxx, gdyy, gdzz, gdxy, gdyz, gdzx,
$ gutt, gutx, guty, gutz,
$ guxx, guyy, guzz, guxy, guyz, guzx,
$ psi, Tmunu_flag)
implicit none
DECLARE_CCTK_PARAMETERS
c input arguments
CCTK_REAL x, y
CCTK_DECLARE(CCTK_REAL, z,)
CCTK_DECLARE(CCTK_REAL, t,)
c output arguments
CCTK_REAL gdtt, gdtx, gdty, gdtz,
$ gdxx, gdyy, gdzz, gdxy, gdyz, gdzx,
$ gutt, gutx, guty, gutz,
$ guxx, guyy, guzz, guxy, guyz, guzx
CCTK_DECLARE(CCTK_REAL, psi,)
LOGICAL Tmunu_flag
CCTK_REAL arad, marad
C This is a vacuum spacetime with no cosmological constant
Tmunu_flag = .false.
arad = Kerr_BoyerLindquist__spin
marad = Kerr_BoyerLindquist__mass
gdtt = -(y**2*arad**2+x**2-2*marad*x)/(x**2+y**2*arad**2)
gdtx = 2*(arad*marad*x*(y**2-1))/(x**2+y**2*arad**2)
gdty = 0.d0
gdtz = 0.d0
gdxx = -(x**4+x**2*arad**2+2*arad**2*marad*x+arad**2*y**2*x**2 - 2*arad**2*y**2*marad*x+arad**4*y**2)*(y**2-1)/(x**2+arad**2*y**2)
gdyy = (x**2+y**2*arad**2)/(x**2-2*marad*x+arad**2)
gdzz = -(x**2+y**2*arad**2)/(y**2-1)
gdxy = 0.d0
gdyz = 0.d0
gdzx = 0.d0
gutt = -(-x**4-x**2*arad**2-2*arad**2*marad*x-arad**2*y**2*x**2 +2*arad**2*y**2*marad*x-arad**4*y**2)/((x**2+arad**2*y**2)*(-x**2+2*marad*x-arad**2))
gutx = 2*(arad*marad*x)/((x**2+arad**2*y**2)*(-x**2+2*marad*x-arad**2))
guty = 0.d0
gutz = 0.d0
guxx = -(-arad**2*y**2-x**2+2*marad*x)/((x**2+arad**2*y**2)*(y**2-1)*(-x**2+2*marad*x-arad**2))
guyy = -(-x**2+2*marad*x-arad**2)/(x**2+arad**2*y**2)
guzz = -(y**2-1)/(x**2+arad**2*y**2)
guxy = 0.d0
guyz = 0.d0
guzx = 0.d0
return
end
|