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
|
/*@@
@file ConfPhys.c
@date September 3rd 1999
@author Gabrielle Allen
@desc
Conversions between physical and conformal metrics.
Be very careful using these functions, note that
conformal_state is not changed, this is up to the
calling routine
@enddesc
@version $Header$
@@*/
#include "cctk.h"
static const char *rcsid = "$Header$";
CCTK_FILEVERSION(CactusEinstein_StaticConformal_ConfPhys_c)
/********************************************************************
********************* Local Data Types ***********************
********************************************************************/
/********************************************************************
********************* Local Routine Prototypes *********************
********************************************************************/
/********************************************************************
***************** Scheduled Routine Prototypes *********************
********************************************************************/
/********************************************************************
********************* Other Routine Prototypes *********************
********************************************************************/
/********************************************************************
********************* Local Data *****************************
********************************************************************/
/********************************************************************
********************* External Routines **********************
********************************************************************/
/*@@
@routine StaticConf_ConfToPhysInPlace
@date September 3rd 1999
@author Gabrielle Allen
@desc
Convert the conformal metric to the physical in place.
@enddesc
@calls
@calledby
@history
@endhistory
@var nx
@vdesc number of points in the x direction
@vtype int
@vio in
@var ny
@vdesc number of points in the y direction
@vtype int
@vio in
@var nz
@vdesc number of points in the z direction
@vtype int
@vio in
@var psi
@vdesc conformal factor
@vtype const CCTK_REAL *
@vio in
@var gxx
@vdesc xx componenent of metric
@vtype CCTK_REAL *
@vio inout
@var gxy
@vdesc xy componenent of metric
@vtype CCTK_REAL *
@vio inout
@var gxz
@vdesc xz componenent of metric
@vtype CCTK_REAL *
@vio inout
@var gyy
@vdesc yy componenent of metric
@vtype CCTK_REAL *
@vio inout
@var gyz
@vdesc yz componenent of metric
@vtype CCTK_REAL *
@vio inout
@var gzz
@vdesc zz componenent of metric
@vtype CCTK_REAL *
@vio inout
@@*/
void StaticConf_ConfToPhysInPlace (CCTK_INT const nx,
CCTK_INT const ny,
CCTK_INT const nz,
const CCTK_REAL * restrict const psi,
CCTK_REAL * restrict const gxx,
CCTK_REAL * restrict const gxy,
CCTK_REAL * restrict const gxz,
CCTK_REAL * restrict const gyy,
CCTK_REAL * restrict const gyz,
CCTK_REAL * restrict const gzz)
{
CCTK_REAL psi4;
int index;
CCTK_WARN (4, "Converting metric: conformal -> physical");
index = nx * ny * nz;
while (--index >= 0)
{
/* this should be faster than psi4 = pow (psi[index], 4) */
psi4 = psi[index] * psi[index];
psi4 = psi4 * psi4;
gxx[index] *= psi4;
gxy[index] *= psi4;
gxz[index] *= psi4;
gyy[index] *= psi4;
gyz[index] *= psi4;
gzz[index] *= psi4;
}
}
/*@@
@routine StaticConf_PhysToConfInPlace
@date September 3rd 1999
@author Gabrielle Allen
@desc
Convert the physical metric to the conformal in place.
@enddesc
@calls
@calledby
@history
@endhistory
@var nx
@vdesc number of points in the x direction
@vtype int
@vio in
@var ny
@vdesc number of points in the y direction
@vtype int
@vio in
@var nz
@vdesc number of points in the z direction
@vtype int
@vio in
@var psi
@vdesc conformal factor
@vtype const CCTK_REAL *
@vio in
@var gxx
@vdesc xx componenent of metric
@vtype CCTK_REAL *
@vio inout
@var gxy
@vdesc xy componenent of metric
@vtype CCTK_REAL *
@vio inout
@var gxz
@vdesc xz componenent of metric
@vtype CCTK_REAL *
@vio inout
@var gyy
@vdesc yy componenent of metric
@vtype CCTK_REAL *
@vio inout
@var gyz
@vdesc yz componenent of metric
@vtype CCTK_REAL *
@vio inout
@var gzz
@vdesc zz componenent of metric
@vtype CCTK_REAL *
@vio inout
@@*/
void StaticConf_PhysToConfInPlace (CCTK_INT const nx,
CCTK_INT const ny,
CCTK_INT const nz,
const CCTK_REAL * restrict const psi,
CCTK_REAL * restrict const gxx,
CCTK_REAL * restrict const gxy,
CCTK_REAL * restrict const gxz,
CCTK_REAL * restrict const gyy,
CCTK_REAL * restrict const gyz,
CCTK_REAL * restrict const gzz)
{
int index;
CCTK_REAL psi4;
CCTK_WARN (4, "Converting metric: physical -> conformal");
index = nx * ny * nz;
while (--index >= 0)
{
/* this should be faster than psi4 = pow (psi[index], 4) */
psi4 = psi[index] * psi[index];
psi4 = psi4 * psi4;
/* get the reciprocal for turning divisions into multiplications */
psi4 = 1.0 / psi4;
gxx[index] *= psi4;
gxy[index] *= psi4;
gxz[index] *= psi4;
gyy[index] *= psi4;
gyz[index] *= psi4;
gzz[index] *= psi4;
}
}
/********************************************************************
********************* Local Routines *************************
********************************************************************/
|