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
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
|
/*@@
@file AHFinder_exp.F
@date April 1998
@author Miguel Alcubierre
@desc
Here I calculate the (normalized) expansion of the
horizon function. The expansion is defined by the
following expression:
__2 2 / ab __a__b \
exp = \/ f / u + d f d f / u | K - \/ \/ f / u | - trK
a b \ /
where:
/ mn \ 1/2
u = | g d f d f |
\ m n /
@enddesc
@version $Header$
@@*/
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
subroutine AHFinder_exp(CCTK_ARGUMENTS)
use AHFinder_dat
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
integer i,j,k
integer CCTK_Equals
CCTK_REAL det,idet
CCTK_REAL gdxx,gdyy,gdzz,gdxy,gdxz,gdyz
CCTK_REAL guxx,guyy,guzz,guxy,guxz,guyz
CCTK_REAL kdxx,kdyy,kdzz,kdxy,kdxz,kdyz
CCTK_REAL ddxxx,ddxyy,ddxzz,ddxxy,ddxxz,ddxyz
CCTK_REAL ddyxx,ddyyy,ddyzz,ddyxy,ddyxz,ddyyz
CCTK_REAL ddzxx,ddzyy,ddzzz,ddzxy,ddzxz,ddzyz
CCTK_REAL d2fxx,d2fyy,d2fzz,d2fxy,d2fxz,d2fyz
CCTK_REAL c2fxx,c2fyy,c2fzz,c2fxy,c2fxz,c2fyz
CCTK_REAL dfx,dfy,dfz,dfux,dfuy,dfuz
CCTK_REAL aux,T0,T1,T2,T3,T4
CCTK_REAL idx,idy,idz
CCTK_REAL idx2,idy2,idz2,idxy,idxz,idyz
CCTK_REAL zero,half,one,two
! Description of variables:
!
! i,j,k, counters
!
! gdij g
! ij
! ij
! guij g
!
!
! det det(g)
!
! idet 1/det
!
!
! kdij K
! ij
!
! ddijk d g
! i jk
!
! dfi d f
! i
! ij
! dfui g d f
! j
!
! d2fij d d f
! i j
!
! c2fij f
! ;ij
!
! aux,T* auxilliary variables
! *******************
! *** NUMBERS ***
! *******************
zero = 0.0d0
half = 0.5d0
one = 1.0d0
two = 2.0d0
! **************************
! *** FIND EXPANSION ***
! **************************
idx = half/dx
idy = half/dy
idz = half/dz
idx2 = one/dx**2
idy2 = one/dy**2
idz2 = one/dz**2
idxy = idx*idy
idxz = idx*idz
idyz = idy*idz
do k=2,nz-1
do j=2,ny-1
do i=2,nx-1
! Find spatial metric.
gdxx = gxx(i,j,k)
gdyy = gyy(i,j,k)
gdzz = gzz(i,j,k)
gdxy = gxy(i,j,k)
gdxz = gxz(i,j,k)
gdyz = gyz(i,j,k)
! Find extrinsic curvature.
kdxx = kxx(i,j,k)
kdyy = kyy(i,j,k)
kdzz = kzz(i,j,k)
kdxy = kxy(i,j,k)
kdxz = kxz(i,j,k)
kdyz = kyz(i,j,k)
! Find derivatives of metric using finite differences.
ddxxx = idx*(gxx(i+1,j,k) - gxx(i-1,j,k))
ddxyy = idx*(gyy(i+1,j,k) - gyy(i-1,j,k))
ddxzz = idx*(gzz(i+1,j,k) - gzz(i-1,j,k))
ddxxy = idx*(gxy(i+1,j,k) - gxy(i-1,j,k))
ddxxz = idx*(gxz(i+1,j,k) - gxz(i-1,j,k))
ddxyz = idx*(gyz(i+1,j,k) - gyz(i-1,j,k))
ddyxx = idy*(gxx(i,j+1,k) - gxx(i,j-1,k))
ddyyy = idy*(gyy(i,j+1,k) - gyy(i,j-1,k))
ddyzz = idy*(gzz(i,j+1,k) - gzz(i,j-1,k))
ddyxy = idy*(gxy(i,j+1,k) - gxy(i,j-1,k))
ddyxz = idy*(gxz(i,j+1,k) - gxz(i,j-1,k))
ddyyz = idy*(gyz(i,j+1,k) - gyz(i,j-1,k))
ddzxx = idz*(gxx(i,j,k+1) - gxx(i,j,k-1))
ddzyy = idz*(gyy(i,j,k+1) - gyy(i,j,k-1))
ddzzz = idz*(gzz(i,j,k+1) - gzz(i,j,k-1))
ddzxy = idz*(gxy(i,j,k+1) - gxy(i,j,k-1))
ddzxz = idz*(gxz(i,j,k+1) - gxz(i,j,k-1))
ddzyz = idz*(gyz(i,j,k+1) - gyz(i,j,k-1))
! Find determinant of spatial metric.
det = gdxx*gdyy*gdzz + two*gdxy*gdxz*gdyz
. - gdxx*gdyz**2 - gdyy*gdxz**2 - gdzz*gdxy**2
! If determinant is not zero proceed.
if (det.gt.zero) then
idet = one/det
! Find inverse spatial metric.
guxx = idet*(gdyy*gdzz - gdyz**2)
guyy = idet*(gdxx*gdzz - gdxz**2)
guzz = idet*(gdxx*gdyy - gdxy**2)
guxy = idet*(gdxz*gdyz - gdzz*gdxy)
guxz = idet*(gdxy*gdyz - gdyy*gdxz)
guyz = idet*(gdxy*gdxz - gdxx*gdyz)
! Find spatial derivatives of f.
T0 = two*ahfgrid(i,j,k)
T1 = ahfgrid(i+1,j,k)
T2 = ahfgrid(i-1,j,k)
dfx = (T1 - T2)*idx
d2fxx = (T1 - T0 + T2)*idx2
T1 = ahfgrid(i,j+1,k)
T2 = ahfgrid(i,j-1,k)
dfy = (T1 - T2)*idy
d2fyy = (T1 - T0 + T2)*idy2
T1 = ahfgrid(i,j,k+1)
T2 = ahfgrid(i,j,k-1)
dfz = (T1 - T2)*idz
d2fzz = (T1 - T0 + T2)*idz2
! Save gradient of horizon function and its norm
! (they will be needed later in the flow algorithm
! and in the gaussian curvature). Here I also try
! to avoid possible division by zero later by not
! allowing ahfgradn(i,j,k) = 0. This should only
! ever happen far from the horizon, so resetting this
! to 1 should have no important effects.
ahfgradx(i,j,k) = dfx
ahfgrady(i,j,k) = dfy
ahfgradz(i,j,k) = dfz
aux = guxx*dfx**2 + guyy*dfy**2 + guzz*dfz**2
. + two*(guxy*dfx*dfy + guxz*dfx*dfz + guyz*dfy*dfz)
ahfgradn(i,j,k) = sqrt(aux)
if (ahfgradn(i,j,k).eq.zero) ahfgradn(i,j,k) = one
! Find crossed derivatives.
d2fxy = (ahfgrid(i+1,j+1,k) + ahfgrid(i-1,j-1,k)
. - ahfgrid(i+1,j-1,k) - ahfgrid(i-1,j+1,k))*idxy
d2fxz = (ahfgrid(i+1,j,k+1) + ahfgrid(i-1,j,k-1)
. - ahfgrid(i+1,j,k-1) - ahfgrid(i-1,j,k+1))*idxz
d2fyz = (ahfgrid(i,j+1,k+1) + ahfgrid(i,j-1,k-1)
. - ahfgrid(i,j+1,k-1) - ahfgrid(i,j-1,k+1))*idyz
! Raise indices in first derivatives.
dfux = guxx*dfx + guxy*dfy + guxz*dfz
dfuy = guxy*dfx + guyy*dfy + guyz*dfz
dfuz = guxz*dfx + guyz*dfy + guzz*dfz
! Find second covariant derivatives of f.
c2fxx = d2fxx - half*(dfux*ddxxx
. + dfuy*(two*ddxxy - ddyxx)
. + dfuz*(two*ddxxz - ddzxx))
c2fyy = d2fyy - half*(dfuy*ddyyy
. + dfux*(two*ddyxy - ddxyy)
. + dfuz*(two*ddyyz - ddzyy))
c2fzz = d2fzz - half*(dfuz*ddzzz
. + dfux*(two*ddzxz - ddxzz)
. + dfuy*(two*ddzyz - ddyzz))
c2fxy = d2fxy - half*(dfux*ddyxx + dfuy*ddxyy
. + dfuz*(ddxyz + ddyxz - ddzxy))
c2fxz = d2fxz - half*(dfux*ddzxx + dfuz*ddxzz
. + dfuy*(ddxyz + ddzxy - ddyxz))
c2fyz = d2fyz - half*(dfuy*ddzyy + dfuz*ddyzz
. + dfux*(ddyxz + ddzxy - ddxyz))
! / m \ 1/2
! Find: u = | d f d f |
! \ m /
T0 = dfx*dfux + dfy*dfuy + dfz*dfuz
! if T0 is positive proceed.
if (T0.gt.zero) then
T0 = one/sqrt(T0)
! __2
! Find: \/ f / u
T1 = guxx*c2fxx + guyy*c2fyy + guzz*c2fzz
. + two*(guxy*c2fxy + guxz*c2fxz + guyz*c2fyz)
T1 = T1*T0
! a b 2
! Find: K d f d f / u
! ab
T2 = kdxx*dfux**2 + kdyy*dfuy**2 + kdzz*dfuz**2
. + two*(dfux*(kdxy*dfuy + kdxz*dfuz) + kdyz*dfuy*dfuz)
T2 = T2*T0**2
! __a__b 3
! Find: \/ \/ f d f d f / u
! a b
T3 = c2fxx*dfux**2 + c2fyy*dfuy**2 + c2fzz*dfuz**2
. + two*(dfux*(c2fxy*dfuy + c2fxz*dfuz)
. + c2fyz*dfuy*dfuz)
T3 = T3*T0**3
! Find: trK
T4 = guxx*kdxx + guyy*kdyy + guzz*kdzz
. + two*(guxy*kdxy + guxz*kdxz + guyz*kdyz)
! Find the expansion.
ahf_exp(i,j,k) = T1 + T2 - T3 - T4
! T0 was not positive.
else
ahf_exp(i,j,k) = zero
end if
! Determinant was zero.
else
ahf_exp(i,j,k) = zero
end if
end do
end do
end do
! Zero out the edges before doing the boundaries.
ahf_exp( 1, 1, :) = 0.0d0
ahf_exp( 1,ny, :) = 0.0d0
ahf_exp( 1, :, 1) = 0.0d0
ahf_exp( 1, :,nz) = 0.0d0
ahf_exp(nx, 1, :) = 0.0d0
ahf_exp(nx,ny, :) = 0.0d0
ahf_exp(nx, :, 1) = 0.0d0
ahf_exp(nx, :,nz) = 0.0d0
ahf_exp( :, 1, 1) = 0.0d0
ahf_exp( :,ny, 1) = 0.0d0
ahf_exp( :, 1,nz) = 0.0d0
ahf_exp( :,ny,nz) = 0.0d0
ahfgradx( 1, 1, :) = 0.0d0
ahfgradx( 1,ny, :) = 0.0d0
ahfgradx( 1, :, 1) = 0.0d0
ahfgradx( 1, :,nz) = 0.0d0
ahfgradx(nx, 1, :) = 0.0d0
ahfgradx(nx,ny, :) = 0.0d0
ahfgradx(nx, :, 1) = 0.0d0
ahfgradx(nx, :,nz) = 0.0d0
ahfgradx( :, 1, 1) = 0.0d0
ahfgradx( :,ny, 1) = 0.0d0
ahfgradx( :, 1,nz) = 0.0d0
ahfgradx( :,ny,nz) = 0.0d0
ahfgrady( 1, 1, :) = 0.0d0
ahfgrady( 1,ny, :) = 0.0d0
ahfgrady( 1, :, 1) = 0.0d0
ahfgrady( 1, :,nz) = 0.0d0
ahfgrady(nx, 1, :) = 0.0d0
ahfgrady(nx,ny, :) = 0.0d0
ahfgrady(nx, :, 1) = 0.0d0
ahfgrady(nx, :,nz) = 0.0d0
ahfgrady( :, 1, 1) = 0.0d0
ahfgrady( :,ny, 1) = 0.0d0
ahfgrady( :, 1,nz) = 0.0d0
ahfgrady( :,ny,nz) = 0.0d0
ahfgradz( 1, 1, :) = 0.0d0
ahfgradz( 1,ny, :) = 0.0d0
ahfgradz( 1, :, 1) = 0.0d0
ahfgradz( 1, :,nz) = 0.0d0
ahfgradz(nx, 1, :) = 0.0d0
ahfgradz(nx,ny, :) = 0.0d0
ahfgradz(nx, :, 1) = 0.0d0
ahfgradz(nx, :,nz) = 0.0d0
ahfgradz( :, 1, 1) = 0.0d0
ahfgradz( :,ny, 1) = 0.0d0
ahfgradz( :, 1,nz) = 0.0d0
ahfgradz( :,ny,nz) = 0.0d0
ahfgradn( 1, 1, :) = 0.0d0
ahfgradn( 1,ny, :) = 0.0d0
ahfgradn( 1, :, 1) = 0.0d0
ahfgradn( 1, :,nz) = 0.0d0
ahfgradn(nx, 1, :) = 0.0d0
ahfgradn(nx,ny, :) = 0.0d0
ahfgradn(nx, :, 1) = 0.0d0
ahfgradn(nx, :,nz) = 0.0d0
ahfgradn( :, 1, 1) = 0.0d0
ahfgradn( :,ny, 1) = 0.0d0
ahfgradn( :, 1,nz) = 0.0d0
ahfgradn( :,ny,nz) = 0.0d0
! Boundaries on x direction.
ahf_exp(1,2:ny-1,2:nz-1) = 2.0D0*ahf_exp(2,2:ny-1,2:nz-1)
. - ahf_exp(3,2:ny-1,2:nz-1)
ahfgradx(1,2:ny-1,2:nz-1) = 2.0D0*ahfgradx(2,2:ny-1,2:nz-1)
. - ahfgradx(3,2:ny-1,2:nz-1)
ahfgrady(1,2:ny-1,2:nz-1) = 2.0D0*ahfgrady(2,2:ny-1,2:nz-1)
. - ahfgrady(3,2:ny-1,2:nz-1)
ahfgradz(1,2:ny-1,2:nz-1) = 2.0D0*ahfgradz(2,2:ny-1,2:nz-1)
. - ahfgradz(3,2:ny-1,2:nz-1)
ahfgradn(1,2:ny-1,2:nz-1) = 2.0D0*ahfgradn(2,2:ny-1,2:nz-1)
. - ahfgradn(3,2:ny-1,2:nz-1)
ahf_exp(nx,2:ny-1,2:nz-1) = 2.0D0*ahf_exp(nx-1,2:ny-1,2:nz-1)
. - ahf_exp(nx-2,2:ny-1,2:nz-1)
ahfgradx(nx,2:ny-1,2:nz-1) = 2.0D0*ahfgradx(nx-1,2:ny-1,2:nz-1)
. - ahfgradx(nx-2,2:ny-1,2:nz-1)
ahfgrady(nx,2:ny-1,2:nz-1) = 2.0D0*ahfgrady(nx-1,2:ny-1,2:nz-1)
. - ahfgrady(nx-2,2:ny-1,2:nz-1)
ahfgradz(nx,2:ny-1,2:nz-1) = 2.0D0*ahfgradz(nx-1,2:ny-1,2:nz-1)
. - ahfgradz(nx-2,2:ny-1,2:nz-1)
ahfgradn(nx,2:ny-1,2:nz-1) = 2.0D0*ahfgradn(nx-1,2:ny-1,2:nz-1)
. - ahfgradn(nx-2,2:ny-1,2:nz-1)
! Boundaries on y direction.
ahf_exp(2:nx-1,1,2:nz-1) = 2.0D0*ahf_exp(2:nx-1,2,2:nz-1)
. - ahf_exp(2:nx-1,3,2:nz-1)
ahfgradx(2:nx-1,1,2:nz-1) = 2.0D0*ahfgradx(2:nx-1,2,2:nz-1)
. - ahfgradx(2:nx-1,3,2:nz-1)
ahfgrady(2:nx-1,1,2:nz-1) = 2.0D0*ahfgrady(2:nx-1,2,2:nz-1)
. - ahfgrady(2:nx-1,3,2:nz-1)
ahfgradz(2:nx-1,1,2:nz-1) = 2.0D0*ahfgradz(2:nx-1,2,2:nz-1)
. - ahfgradz(2:nx-1,3,2:nz-1)
ahfgradn(2:nx-1,1,2:nz-1) = 2.0D0*ahfgradn(2:nx-1,2,2:nz-1)
. - ahfgradn(2:nx-1,3,2:nz-1)
ahf_exp(2:nx-1,ny,2:nz-1) = 2.0D0*ahf_exp(2:nx-1,ny-1,2:nz-1)
. - ahf_exp(2:nx-1,ny-2,2:nz-1)
ahfgradx(2:nx-1,ny,2:nz-1) = 2.0D0*ahfgradx(2:nx-1,ny-1,2:nz-1)
. - ahfgradx(2:nx-1,ny-2,2:nz-1)
ahfgrady(2:nx-1,ny,2:nz-1) = 2.0D0*ahfgrady(2:nx-1,ny-1,2:nz-1)
. - ahfgrady(2:nx-1,ny-2,2:nz-1)
ahfgradz(2:nx-1,ny,2:nz-1) = 2.0D0*ahfgradz(2:nx-1,ny-1,2:nz-1)
. - ahfgradz(2:nx-1,ny-2,2:nz-1)
ahfgradn(2:nx-1,ny,2:nz-1) = 2.0D0*ahfgradn(2:nx-1,ny-1,2:nz-1)
. - ahfgradn(2:nx-1,ny-2,2:nz-1)
! Boundaries on z direction.
ahf_exp(2:nx-1,2:ny-1,1) = 2.0D0*ahf_exp(2:nx-1,2:ny-1,2)
. - ahf_exp(2:nx-1,2:ny-1,3)
ahfgradx(2:nx-1,2:ny-1,1) = 2.0D0*ahfgradx(2:nx-1,2:ny-1,2)
. - ahfgradx(2:nx-1,2:ny-1,3)
ahfgrady(2:nx-1,2:ny-1,1) = 2.0D0*ahfgrady(2:nx-1,2:ny-1,2)
. - ahfgrady(2:nx-1,2:ny-1,3)
ahfgradz(2:nx-1,2:ny-1,1) = 2.0D0*ahfgradz(2:nx-1,2:ny-1,2)
. - ahfgradz(2:nx-1,2:ny-1,3)
ahfgradn(2:nx-1,2:ny-1,1) = 2.0D0*ahfgradn(2:nx-1,2:ny-1,2)
. - ahfgradn(2:nx-1,2:ny-1,3)
ahf_exp(2:nx-1,2:ny-1,nz) = 2.0D0*ahf_exp(2:nx-1,2:ny-1,nz-1)
. - ahf_exp(2:nx-1,2:ny-1,nz-2)
ahfgradx(2:nx-1,2:ny-1,nz) = 2.0D0*ahfgradx(2:nx-1,2:ny-1,nz-1)
. - ahfgradx(2:nx-1,2:ny-1,nz-2)
ahfgrady(2:nx-1,2:ny-1,nz) = 2.0D0*ahfgrady(2:nx-1,2:ny-1,nz-1)
. - ahfgrady(2:nx-1,2:ny-1,nz-2)
ahfgradz(2:nx-1,2:ny-1,nz) = 2.0D0*ahfgradz(2:nx-1,2:ny-1,nz-1)
. - ahfgradz(2:nx-1,2:ny-1,nz-2)
ahfgradn(2:nx-1,2:ny-1,nz) = 2.0D0*ahfgradn(2:nx-1,2:ny-1,nz-1)
. - ahfgradn(2:nx-1,2:ny-1,nz-2)
! Synchronize.
call CCTK_SyncGroup(ierr,cctkGH,"ahfinder::ahfgradient")
call CCTK_SyncGroup(ierr,cctkGH,"ahfinder::ahfinderexp")
call CartSymGN(ierr,cctkGH,"ahfinder::ahfgradient")
call CartSymGN(ierr,cctkGH,"ahfinder::ahfinderexp")
! ***************
! *** END ***
! ***************
return
end
|