aboutsummaryrefslogtreecommitdiff
path: root/src/metrics/Schwarzschild_EF.F77
blob: 66a9db47b1cc6e063418d8525079f0c8ec4b55a7 (plain)
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
c Schwarzschild metric in Eddington-Finkelstein coordinates,
c as per MTW box 31.2
C
C Author: unknown
C Copyright/License: unknown
C
c $Header$

#include "cctk.h"
#include "cctk_Parameters.h"

      subroutine Exact__Schwarzschild_EF(
     $     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, 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

c local static variables
      CCTK_REAL eps, m

c local variables
      CCTK_REAL r

C This is a vacuum spacetime with no cosmological constant
      Tmunu_flag = .false.

C     Get parameters of the exact solution.

      eps = Schwarzschild_EF__epsilon
      m   = Schwarzschild_EF__mass

      r = max(sqrt(x**2 + y**2 + z**2), eps)

      gdtt = - (1.d0 - 2.d0 * m / r)
      gdtx = 2.d0 * m * x / r**2
      gdty = 2.d0 * m * y / r**2
      gdtz = 2.d0 * m * z / r**2
      gdxx = 1.d0 + 2.d0 * m * x**2 / r**3
      gdyy = 1.d0 + 2.d0 * m * y**2 / r**3
      gdzz = 1.d0 + 2.d0 * m * z**2 / r**3
      gdxy = 2.d0 * m * x * y / r**3
      gdyz = 2.d0 * m * y * z / r**3
      gdzx = 2.d0 * m * z * x / r**3

      gutt = - (1.d0 + 2.d0 * m / r)
      gutx = 2.d0 * m * x / r**2
      guty = 2.d0 * m * y / r**2
      gutz = 2.d0 * m * z / r**2
      guxx = 1.d0 - 2.d0 * m * x**2 / r**3
      guyy = 1.d0 - 2.d0 * m * y**2 / r**3
      guzz = 1.d0 - 2.d0 * m * z**2 / r**3
      guxy = - 2.d0 * m * x * y / r**3
      guyz = - 2.d0 * m * y * z / r**3
      guzx = - 2.d0 * m * z * x / r**3

      return
      end