aboutsummaryrefslogtreecommitdiff
path: root/src/metrics/Schwarzschild_BL.F
blob: ed6e220abe3728083f3c5f6530749660d6a7e2ea (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
c Schwarzschild spacetime in Brill-Lindquist coordinates.
C
c $Header$

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

      subroutine Exact__Schwarzschild_BL(
     $     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
      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

c local variables
      CCTK_REAL eps, m

c local variables
      CCTK_REAL r, psi4

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

C     Get parameters of the exact solution.

      eps = Schwarzschild_BL__epsilon
      m   = Schwarzschild_BL__mass

      r = ((x**2 + y**2 + z**2)**2 + eps**4) ** 0.25d0
      psi4 = (1 + m / (2 * r)) ** 4

      gdtt = -1
      gdtx = 0
      gdty = 0
      gdtz = 0
      gdxx = psi4
      gdyy = psi4
      gdzz = psi4
      gdxy = 0
      gdyz = 0
      gdzx = 0

      gutt = -1
      gutx = 0
      guty = 0
      gutz = 0
      guxx = 1 / psi4
      guyy = 1 / psi4
      guzz = 1 / psi4
      guxy = 0
      guyz = 0
      guzx = 0

      return
      end