aboutsummaryrefslogtreecommitdiff
path: root/src/metrics/Bertotti.F77
blob: 003a2829ab78c82c370a9adafefd5f86e3a5438b (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 Bertotti spacetime with cosmological constant
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__Bertotti(
     $     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 variables
      CCTK_REAL baza
      CCTK_REAL unu, doi
 
C this model has a cosmological constant
C ==> it sets the stress-energy tensor in the "CalcTmunu" code
      Tmunu_flag = .true.

      baza = Bertotti__Lambda

      unu=exp(2.0D0*sqrt(-baza)*x)
      doi=exp(2.0D0*sqrt(-baza)*z)   

       
      gdtt = -unu
      gdtx = 0.0D0 
      gdty = 0.0D0
      gdtz = 0.0D0
      gdxx = 1.0D0
      gdyy = doi
      gdzz = 1.0D0
      gdxy = 0.0D0
      gdyz = 0.0D0
      gdzx = 0.0D0

      gutt = -1.0D0/unu
      gutx = 0.0D0
      guty = 0.0D0
      gutz = 0.0D0
      guxx = 1.0D0
      guyy = 1.0D0/doi
      guzz = 1.0D0
      guxy = 0.0D0
      guyz = 0.0D0
      guzx = 0.0D0


      return
      end