aboutsummaryrefslogtreecommitdiff
path: root/src/metric.F
blob: f28c948e1695626bdd945a135ce0488d02915c73 (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
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
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
c This subroutine calculates the 4-metric and its inverse at an event,
c for a given model, by decoding  decoded_exact_model  and calling the
c appropriate subroutine for that model.
C $Header$

#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
#include "cctk_Functions.h"

#include "param_defs.inc"

      subroutine Exact__metric_for_model(
     $     decoded_exact_model,
     $     x, y, z, t,
     $     gdtt, gdtx, gdty, gdtz,
     $     gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $     gutt, gutx, guty, gutz,
     $     guxx, guyy, guzz, guxy, guyz, guxz,
     $     psi, Tmunu_flag)

      implicit none
      DECLARE_CCTK_FUNCTIONS

c arguments
      CCTK_INT decoded_exact_model
      CCTK_REAL x, y, z, t,
     $     gdtt, gdtx, gdty, gdtz,
     $     gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $     gutt, gutx, guty, gutz,
     $     guxx, guyy, guzz, guxy, guyz, guxz
      CCTK_REAL psi
      LOGICAL   Tmunu_flag

c local variables
      character*100 warn_buffer

cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc

c
c Minkowski spacetime
c

      if     (decoded_exact_model .eq. EXACT__Minkowski) then
         call Exact__Minkowski(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Minkowski_shift) then
         call Exact__Minkowski_shift(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Minkowski_funny) then
         call Exact__Minkowski_funny(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Minkowski_gauge_wave) then
         call Exact__Minkowski_gauge_wave(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Minkowski_shifted_gauge_wave) then
         call Exact__Minkowski_shifted_gauge_wave(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Minkowski_conf_wave) then
         call Exact__Minkowski_conf_wave(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc

c
c black hole spacetimes
c

      elseif (decoded_exact_model .eq. EXACT__Schwarzschild_EF) then
         call Exact__Schwarzschild_EF(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Schwarzschild_PG) then
         call Exact__Schwarzschild_PG(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Schwarzschild_BL) then
         call Exact__Schwarzschild_BL(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Schwarzschild_Novikov) then
         call Exact__Schwarzschild_Novikov(x,y,z,t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Kerr_BoyerLindquist) then
         call Exact__Kerr_BoyerLindquist(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Kerr_KerrSchild) then
         call Exact__Kerr_KerrSchild(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Kerr_KerrSchild_spherical) then
         call Exact__Kerr_KerrSchild_spherical(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Schwarzschild_Lemaitre) then
         call Exact__Schwarzschild_Lemaitre(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__multi_BH) then
         call Exact__multi_BH(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

c
c not fully implemented yet -- see Nina Jansen for details
c
c     elseif (decoded_exact_model .eq. EXACT__Alvi) then
c        call Exact__Alvi(
c    $        x, y, z, t,
c    $        gdtt, gdtx, gdty, gdtz,
c    $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
c    $        gutt, gutx, guty, gutz,
c    $        guxx, guyy, guzz, guxy, guyz, guxz,
c    $        psi, Tmunu_flag)
c

      elseif (decoded_exact_model .eq. EXACT__Thorne_fakebinary) then
         call Exact__Thorne_fakebinary(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc

c
c cosmological spacetimes
c

      elseif (decoded_exact_model .eq. EXACT__Lemaitre) then
         call Exact__Lemaitre(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__de_Sitter) then
         call Exact__de_Sitter(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__de_Sitter_Lambda) then
         call Exact__de_Sitter_Lambda(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__anti_de_Sitter_Lambda) then
         call Exact__anti_de_Sitter_Lambda(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Bianchi_I) then
         call Exact__Bianchi_I(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Goedel) then
         call Exact__Goedel(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Bertotti) then
         call Exact__Bertotti(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Kasner_like) then
         call Exact__Kasner_like(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Kasner_axisymmetric) then
         call Exact__Kasner_axisymmetric(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Kasner_generalized) then
         call Exact__Kasner_generalized(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Gowdy_wave) then
         call Exact__Gowdy_wave(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__Milne) then
         call Exact__Milne(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc

c
c miscellaneous spacetimes
c

      elseif (decoded_exact_model .eq. EXACT__boost_rotation_symmetric) then
         call Exact__boost_rotation_symmetric(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__bowl) then
         call Exact__bowl(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

      elseif (decoded_exact_model .eq. EXACT__constant_density_star) then
         call Exact__constant_density_star(
     $        x, y, z, t,
     $        gdtt, gdtx, gdty, gdtz,
     $        gdxx, gdyy, gdzz, gdxy, gdyz, gdxz,
     $        gutt, gutx, guty, gutz,
     $        guxx, guyy, guzz, guxy, guyz, guxz,
     $        psi, Tmunu_flag)

cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc

      else
         write (warn_buffer, '(a,i8)')
     $         'Unknown decoded_exact_model = ', decoded_exact_model
         call CCTK_WARN(0, warn_buffer)
      endif

      return
      end