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
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
|
C Robertson-Walker universe
C See: J.N. Islam, An Introduction to
C Mathematical Cosmology, Cambridge, 1992
C
C Author : D. Vulcanov (Timisoara, Romania)
C see ../../README for copyright & licensing info
C
C $Header$
C
C FIXME:
C This metric doesn't work. The argument rama is the R(t) in the
C Robertson-Walker metric, and if this is passed correctly then this
C subroutine computes the correct metric. But the rest of this thorn
C doesn't know to pass this value. :( :( See Mitica's Cosmo thorn
C for a better way to get this metric.
#include "cctk.h"
#include "cctk_Parameters.h"
subroutine Exact__Robertson_Walker(
$ 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,
$ rama)
implicit none
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_FUNCTIONS
c input arguments
CCTK_REAL x, y, z, 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_REAL psi
LOGICAL Tmunu_flag
CCTK_REAL rama
c local static variables
logical firstcall
CCTK_REAL kapa
data firstcall /.true./
save firstcall, kapa
c local variables
CCTK_REAL r2,x2,y2,z2,am,ag
C This model may set the stress-energy tensor
Tmunu_flag = .true.
if (firstcall) then
kapa = Robertson_Walker__k
firstcall = .false.
end if
x2=x*x
y2=y*y
z2=z*z
r2 = x2+y2+z2
ag = kapa/(1.0D0-kapa*r2)
am = rama*rama
gdtt = -1.0D0
gdtx = 0.0D0
gdty = 0.0D0
gdtz = 0.0D0
gdxx = am*(1.0D0 + ag*x2)
gdyy = am*(1.0D0 + ag*y2)
gdzz = am*(1.0D0 + ag*z2)
gdxy = am*ag*x*y
gdyz = am*ag*y*z
gdzx = am*ag*z*x
gutt = -1.0D0
gutx = 0.0D0
guty = 0.0D0
gutz = 0.0D0
guxx = (1.0D0-kapa*x2)/am
guyy = (1.0D0-kapa*y2)/am
guzz = (1.0D0-kapa*z2)/am
guxy = -kapa*x*y/am
guyz = -kapa*y*z/am
guzx = -kapa*z*x/am
return
end
|