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
363
364
365
366
|
/*@@
@file GRHydro_Wrappers.F90
@date Sat Jun 29 18:57:07 PDT 2013
@author Roland Haas
@desc
Wrapper routines to provide a consistent interface to the aliased
functions even if the interface to the underlying implementation changes.
@enddesc
@@*/
#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"
#include "cctk_Functions.h"
subroutine Con2PrimPolyWrapper(handle, dens, sx, sy, sz, tau, rho, &
velx, vely, velz, epsilon, press, w_lorentz, uxx, uxy, uxz, uyy, &
uyz, uzz, detg, x, y, z, r, GRHydro_rho_min, GRHydro_reflevel, GRHydro_C2P_failed)
use Con2Prim_fortran_interfaces
implicit none
DECLARE_CCTK_FUNCTIONS
CCTK_REAL, INTENT(INOUT) :: dens, sx, sy, sz
CCTK_REAL, INTENT(OUT) :: tau
CCTK_REAL, INTENT(INOUT) :: rho
CCTK_REAL, INTENT(INOUT) :: velx, vely, velz, epsilon, press, w_lorentz
CCTK_REAL, INTENT(IN) :: uxx, uxy, uxz, uyy, uyz, uzz, detg
CCTK_REAL, INTENT(IN) :: x, y, z, r, GRHydro_rho_min
CCTK_INT, INTENT(IN) :: handle, GRHydro_reflevel
CCTK_REAL, INTENT(OUT) :: GRHydro_C2P_failed
CCTK_REAL :: sdetg
sdetg = sqrt(detg)
call Con2Prim_ptPolytype(handle, dens, sx, sy, sz, tau, rho, &
velx, vely, velz, epsilon, press, w_lorentz, uxx, uxy, uxz, uyy, &
uyz, uzz, sdetg, x, y, z, r, GRHydro_rho_min, GRHydro_reflevel, &
GRHydro_C2P_failed)
end subroutine
subroutine Con2PrimGenWrapper(handle, dens, sx, sy, sz, tau, rho, velx, vely, &
velz, epsilon, pressure, w_lorentz, uxx, uxy, &
uxz, uyy, uyz, uzz, det, x, y, z, r, epsnegative, &
GRHydro_rho_min, pmin, epsmin, GRHydro_reflevel, &
retval)
use Con2Prim_fortran_interfaces
implicit none
DECLARE_CCTK_FUNCTIONS
CCTK_INT, INTENT(IN) :: handle
CCTK_REAL, INTENT(INOUT) :: dens, sx, sy, sz, tau
CCTK_REAL, INTENT(INOUT) :: rho, velx, vely, velz, epsilon, pressure
CCTK_REAL, INTENT(INOUT) :: w_lorentz
CCTK_REAL, INTENT(IN) :: uxx, uxy, uxz, uyy, uyz, uzz, det
CCTK_REAL, INTENT(IN) :: x, y, z, r
CCTK_INT, INTENT(OUT) :: epsnegative
CCTK_REAL, INTENT(IN) :: GRHydro_rho_min, pmin, epsmin
CCTK_INT, INTENT(IN) :: GRHydro_reflevel
CCTK_REAL, INTENT(OUT) :: retval
CCTK_REAL :: sdetg
CCTK_INT, parameter :: cctk_iteration = -1, ii = -1, jj = -1, kk = -1
logical :: log_epsnegative
sdetg = sqrt(det)
call Con2Prim_pt(cctK_iteration, ii, jj, kk, &
handle, dens, sx, sy, sz, tau, rho, velx, vely, &
velz, epsilon, pressure, w_lorentz, uxx, uxy, &
uxz, uyy, uyz, uzz, sdetg, x, y, z, r, log_epsnegative, &
GRHydro_rho_min, pmin, epsmin, GRHydro_reflevel, &
retval)
if(log_epsnegative) then
epsnegative = 1
else
epsnegative = 0
end if
end subroutine
subroutine Con2PrimGenMWrapper(handle, keytemp, prec, gamma_eos, dens, sx, sy, sz, &
tau, Bconsx, Bconsy, Bconsz, y_e, temp, rho, velx, &
vely, velz, epsilon, pressure, Bvecx, Bvecy, Bvecz, &
Bvecsq, w_lorentz, gxx, gxy, gxz, gyy, gyz, gzz, uxx, &
uxy, uxz, uyy, uyz, uzz, det, epsnegative, retval)
use Con2PrimM_fortran_interfaces
implicit none
DECLARE_CCTK_FUNCTIONS
CCTK_INT, INTENT(IN) :: handle, keytemp
CCTK_REAL, INTENT(IN) :: prec, gamma_eos
CCTK_REAL, INTENT(INOUT) :: dens, sx, sy, sz, tau
CCTK_REAL, INTENT(IN) :: Bconsx, Bconsy, Bconsz
CCTK_REAL, INTENT(INOUT) :: y_e, temp
CCTK_REAL, INTENT(INOUT) :: rho, velx, vely, velz, epsilon, pressure
CCTK_REAL, INTENT(OUT) :: Bvecx, Bvecy, Bvecz, Bvecsq
CCTK_REAL, INTENT(INOUT) :: w_lorentz
CCTK_REAL, INTENT(IN) :: gxx, gxy, gxz, gyy, gyz, gzz
CCTK_REAL, INTENT(IN) :: uxx, uxy, uxz, uyy, uyz, uzz
CCTK_REAL, INTENT(IN) :: det
CCTK_INT, INTENT(OUT) :: epsnegative
CCTK_REAL, INTENT(OUT) :: retval
CCTK_REAL :: sdetg
sdetg = sqrt(det)
if(handle.eq.1 .or. handle.eq.2) then
call GRHydro_Con2PrimM_ptold(handle, gamma_eos, dens, sx, sy, sz, &
tau, Bconsx, Bconsy, Bconsz, rho, velx, &
vely, velz, epsilon, pressure, Bvecx, Bvecy, Bvecz, &
Bvecsq, w_lorentz, gxx, gxy, gxz, gyy, gyz, gzz, uxx, &
uxy, uxz, uyy, uyz, uzz, sdetg, epsnegative, retval)
else
call GRHydro_Con2PrimM_pt2(handle, keytemp, prec, prec, gamma_eos, dens, sx, sy, sz, &
tau, Bconsx, Bconsy, Bconsz, y_e, temp, rho, velx, &
vely, velz, epsilon, pressure, Bvecx, Bvecy, Bvecz, &
Bvecsq, w_lorentz, gxx, gxy, gxz, gyy, gyz, gzz, uxx, &
uxy, uxz, uyy, uyz, uzz, sdetg, epsnegative, retval)
end if
end subroutine
subroutine Con2PrimPolyMWrapper(handle, gamma_eos, dens, sx, sy, sz, sc,&
Bconsx, Bconsy, Bconsz, rho, velx, vely, velz,&
epsilon, pressure, Bvecx, Bvecy, Bvecz, Bvecsq,&
w_lorentz, gxx, gxy, gxz, gyy, gyz, gzz, uxx,&
uxy, uxz, uyy, uyz, uzz, det, epsnegative,&
retval)
use Con2PrimM_fortran_interfaces
implicit none
DECLARE_CCTK_FUNCTIONS
CCTK_INT, INTENT(IN) :: handle
CCTK_REAL, INTENT(IN) :: gamma_eos
CCTK_REAL, INTENT(INOUT) :: dens, sx, sy, sz, sc
CCTK_REAL, INTENT(INOUT) :: Bconsx, Bconsy, Bconsz
CCTK_REAL, INTENT(INOUT) :: rho, velx, vely, velz, epsilon, pressure
CCTK_REAL, INTENT(OUT) :: Bvecx, Bvecy, Bvecz
CCTK_REAL, INTENT(OUT) :: Bvecsq
CCTK_REAL, INTENT(INOUT) :: w_lorentz
CCTK_REAL, INTENT(IN) :: gxx, gxy, gxz, gyy, gyz, gzz
CCTK_REAL, INTENT(IN) :: uxx, uxy, uxz, uyy, uyz, uzz
CCTK_REAL, INTENT(IN) :: det
CCTK_INT, INTENT(OUT) :: epsnegative
CCTK_REAL, INTENT(OUT) :: retval
CCTK_REAL :: sdetg
sdetg = sqrt(det)
call GRHydro_Con2PrimM_Polytype_pt(handle, gamma_eos, dens, sx, sy, sz, sc, &
Bconsx, Bconsy, Bconsz, rho, velx, vely, velz, &
epsilon, pressure, Bvecx, Bvecy, Bvecz, Bvecsq, &
w_lorentz, gxx, gxy, gxz, gyy, gyz, gzz, uxx, &
uxy, uxz, uyy, uyz, uzz, sdetg, epsnegative, &
retval)
end subroutine
subroutine Prim2ConGenWrapper(handle, gxx, gxy, gxz, gyy, gyz, gzz, det, dens,&
sx, sy, sz, tau, rho, velx, vely, velz, epsilon,&
press, w_lorentz)
use Con2Prim_fortran_interfaces
implicit none
DECLARE_CCTK_FUNCTIONS
CCTK_INT, INTENT(IN) :: handle
CCTK_REAL, INTENT(IN) :: gxx, gxy, gxz, gyy, gyz, gzz
CCTK_REAL, INTENT(IN) :: det
CCTK_REAL, INTENT(OUT) :: dens, sx, sy, sz, tau
CCTK_REAL, INTENT(IN) :: rho, velx, vely, velz, epsilon
CCTK_REAL, INTENT(OUT) :: press, w_lorentz
CCTK_REAL :: sdetg
sdetg = sqrt(det)
call prim2con(handle, gxx, gxy, gxz, gyy, gyz, gzz, sdetg, dens,&
sx, sy, sz, tau, rho, velx, vely, velz, epsilon,&
press, w_lorentz)
end subroutine
subroutine Prim2ConPolyWrapper(handle, gxx, gxy, gxz, gyy, gyz, gzz, det, dens,&
sx, sy, sz, tau, rho, velx, vely, velz, epsilon,&
press, w_lorentz)
use Con2Prim_fortran_interfaces
implicit none
DECLARE_CCTK_FUNCTIONS
CCTK_INT, INTENT(IN) :: handle
CCTK_REAL, INTENT(IN) :: gxx, gxy, gxz, gyy, gyz, gzz
CCTK_REAL, INTENT(IN) :: det
CCTK_REAL, INTENT(OUT) :: dens, sx, sy, sz, tau
CCTK_REAL, INTENT(IN) :: rho, velx, vely, velz
CCTK_REAL, INTENT(OUT) :: epsilon, press, w_lorentz
CCTK_REAL :: sdetg
sdetg = sqrt(det)
call prim2conpolytype(handle, gxx, gxy, gxz, gyy, gyz, gzz, sdetg, dens,&
sx, sy, sz, tau, rho, velx, vely, velz, epsilon,&
press, w_lorentz)
end subroutine
subroutine Prim2ConGenMWrapper(handle, gxx, gxy, gxz, gyy, gyz, gzz, det, dens,&
sx, sy, sz, tau, Bconsx, Bconsy, Bconsz, rho,&
velx, vely, velz, epsilon, press, Bvecx, Bvecy,&
Bvecz, w_lorentz)
use Con2PrimM_fortran_interfaces
implicit none
DECLARE_CCTK_FUNCTIONS
CCTK_INT, INTENT(IN) :: handle
CCTK_REAL, INTENT(IN) :: gxx, gxy, gxz, gyy, gyz, gzz
CCTK_REAL, INTENT(IN) :: det
CCTK_REAL, INTENT(OUT) :: dens, sx, sy, sz, tau
CCTK_REAL, INTENT(OUT) :: Bconsx, Bconsy, Bconsz
CCTK_REAL, INTENT(IN) :: rho, velx, vely, velz, epsilon
CCTK_REAL, INTENT(OUT) :: press
CCTK_REAL, INTENT(IN) :: Bvecx, Bvecy, Bvecz
CCTK_REAL, INTENT(OUT) :: w_lorentz
CCTK_REAL :: sdetg
sdetg = sqrt(det)
call prim2conM(handle, gxx, gxy, gxz, gyy, gyz, gzz, sdetg, dens,&
sx, sy, sz, tau, Bconsx, Bconsy, Bconsz, rho,&
velx, vely, velz, epsilon, press, Bvecx, Bvecy,&
Bvecz, w_lorentz)
end subroutine
subroutine Prim2ConGenM_hotWrapper(handle, GRHydro_reflevel, i, j, k, x, y, z,&
gxx, gxy, gxz, gyy, gyz, gzz, det, dens, sx,&
sy, sz, tau, Bconsx, Bconsy, Bconsz, rho,&
velx, vely, velz, epsilon, press, Bvecx,&
Bvecy, Bvecz, w_lorentz, temperature, Y_e)
use Con2PrimM_fortran_interfaces
implicit none
DECLARE_CCTK_FUNCTIONS
CCTK_INT, INTENT(IN) :: handle, GRHydro_reflevel
CCTK_INT, INTENT(IN) :: i, j, k
CCTK_REAL, INTENT(IN) :: x, y, z
CCTK_REAL, INTENT(IN) :: gxx, gxy, gxz, gyy, gyz, gzz
CCTK_REAL, INTENT(IN) :: det
CCTK_REAL, INTENT(OUT) :: dens, sx, sy, sz, tau
CCTK_REAL, INTENT(OUT) :: Bconsx, Bconsy, Bconsz
CCTK_REAL, INTENT(IN) :: rho, velx, vely, velz, epsilon
CCTK_REAL, INTENT(OUT) :: press
CCTK_REAL, INTENT(IN) :: Bvecx, Bvecy, Bvecz
CCTK_REAL, INTENT(OUT) :: w_lorentz
CCTK_REAL, INTENT(INOUT) :: temperature
CCTK_REAL, INTENT(IN) :: Y_e
CCTK_REAL :: sdetg
sdetg = sqrt(det)
call prim2conM_hot(handle, GRHydro_reflevel, i, j, k, x, y, z, gxx,&
gxy, gxz, gyy, gyz, gzz, sdetg, dens, sx, sy, sz,&
tau, Bconsx, Bconsy, Bconsz, rho, velx, vely, velz,&
epsilon, press, Bvecx, Bvecy, Bvecz, w_lorentz,&
temperature, Y_e)
end subroutine
subroutine Prim2ConPolyMWrapper(handle, gxx, gxy, gxz, gyy, gyz, gzz, det,&
dens, sx, sy, sz, tau, Bconsx, Bconsy, Bconsz,&
rho, velx, vely, velz, epsilon, press, Bvecx,&
Bvecy, Bvecz, w_lorentz)
use Con2PrimM_fortran_interfaces
implicit none
DECLARE_CCTK_FUNCTIONS
CCTK_INT, INTENT(IN) :: handle
CCTK_REAL, INTENT(IN) :: gxx, gxy, gxz, gyy, gyz, gzz
CCTK_REAL, INTENT(IN) :: det
CCTK_REAL, INTENT(OUT) :: dens, sx, sy, sz, tau
CCTK_REAL, INTENT(OUT) :: Bconsx, Bconsy, Bconsz
CCTK_REAL, INTENT(IN) :: rho, velx, vely, velz
CCTK_REAL, INTENT(OUT) :: epsilon
CCTK_REAL, INTENT(OUT) :: press
CCTK_REAL, INTENT(IN) :: Bvecx, Bvecy, Bvecz
CCTK_REAL, INTENT(OUT) :: w_lorentz
CCTK_REAL :: sdetg
sdetg = sqrt(det)
call prim2conpolytypeM(handle, gxx, gxy, gxz, gyy, gyz, gzz, sdetg,&
dens, sx, sy, sz, tau, Bconsx, Bconsy, Bconsz,&
rho, velx, vely, velz, epsilon, press, Bvecx,&
Bvecy, Bvecz, w_lorentz)
end subroutine
subroutine Con2PrimGenMeeWrapper(handle, keytemp, prec, gamma_eos, dens, sx, sy,&
sz, tau, Bconsx, Bconsy, Bconsz, entropycons,&
y_e, temp, rho, velx, vely, velz, epsilon,&
pressure, Bvecx, Bvecy, Bvecz, Bvecsq,&
w_lorentz, gxx, gxy, gxz, gyy, gyz, gzz, uxx,&
uxy, uxz, uyy, uyz, uzz, det, epsnegative,&
retval)
use Con2Prim_fortran_interfaces
implicit none
DECLARE_CCTK_FUNCTIONS
CCTK_INT, INTENT(IN) :: handle, keytemp
CCTK_REAL, INTENT(IN) :: prec, gamma_eos
CCTK_REAL, INTENT(INOUT) :: dens, sx, sy, sz, tau
CCTK_REAL, INTENT(IN) :: Bconsx, Bconsy, Bconsz
CCTK_REAL, INTENT(INOUT) :: entropycons, y_e, temp
CCTK_REAL, INTENT(INOUT) :: rho, velx, vely, velz, epsilon, pressure
CCTK_REAL, INTENT(OUT) :: Bvecx, Bvecy, Bvecz
CCTK_REAL, INTENT(OUT) :: Bvecsq
CCTK_REAL, INTENT(INOUT) :: w_lorentz
CCTK_REAL, INTENT(IN) :: gxx, gxy, gxz, gyy, gyz, gzz
CCTK_REAL, INTENT(IN) :: uxx, uxy, uxz, uyy, uyz, uzz
CCTK_REAL, INTENT(IN) :: det
CCTK_INT, INTENT(OUT) :: epsnegative
CCTK_REAL, INTENT(OUT) :: retval
CCTK_REAL :: sdetg
sdetg = sqrt(det)
call GRHydro_Con2PrimM_ptee(handle, keytemp, prec, gamma_eos, dens, sx, sy, sz,&
tau, Bconsx, Bconsy, Bconsz, entropycons, y_e, temp,&
rho, velx, vely, velz, epsilon, pressure, Bvecx,&
Bvecy, Bvecz, Bvecsq, w_lorentz, gxx, gxy, gxz, gyy,&
gyz, gzz, uxx, uxy, uxz, uyy, uyz, uzz, sdetg,&
epsnegative, retval)
end subroutine
|