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
|
/*@@
@file AHFinder_output.F
@date Fri 2 Aug 2001
@author Thomas Radke
@desc
Does the output of AHFinder data into files.
@enddesc
@@*/
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
subroutine AHFinder_Output(CCTK_ARGUMENTS, report, status, horizon, mtype, intarea_h)
use AHFinder_dat
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
integer mtype
logical report, status, horizon
CCTK_REAL intarea_h
integer l, m, ierror
character*3 ext
character*200 almf,areaf,massf,radf,circeqf,merip1f,merip2f
character*200 asymxf,asymyf,asymzf
character*200 options
! ****************************
! *** BUILD FILE NAMES ***
! ****************************
c some compilers cannot trim an empty string, so we add at least one char
if (find3) then
ext = "_" // char(ichar('0') + mfind - 1) // "."
else
ext = "."
end if
almf = filestr(1:nfile) // "/ahf_coeff" // trim(ext) // "alm"
areaf = filestr(1:nfile) // "/ahf_area" // trim(ext) // "tl"
massf = filestr(1:nfile) // "/ahf_mass" // trim(ext) // "tl"
radf = filestr(1:nfile) // "/ahf_rad" // trim(ext) // "tl"
circeqf = filestr(1:nfile) // "/ahf_circ_eq" // trim(ext) // "tl"
merip1f = filestr(1:nfile) // "/ahf_meri_p1" // trim(ext) // "tl"
merip2f = filestr(1:nfile) // "/ahf_meri_p2" // trim(ext) // "tl"
asymxf = filestr(1:nfile) // "/ahf_asymx" // trim(ext) // "tl"
asymyf = filestr(1:nfile) // "/ahf_asymy" // trim(ext) // "tl"
asymzf = filestr(1:nfile) // "/ahf_asymz" // trim(ext) // "tl"
! *********************************
! *** OUTPUT GRID FUNCTIONS ***
! *********************************
if (ahf_2Doutput.ne.0) then
if (.not.find3) then
call CCTK_OutputVarByMethod(ierror,cctkGH,
. "ahfinder::ahfgrid","IOFlexIO_2D")
call CCTK_OutputVarByMethod(ierror,cctkGH,
. "ahfinder::ahf_exp","IOFlexIO_2D")
else if (mfind.eq.3) then
call CCTK_OutputVarByMethod(ierror,cctkGH,
. "ahfinder::ahfgrid3","IOFlexIO_2D")
call CCTK_OutputVarByMethod(ierror,cctkGH,
. "ahfinder::ahf_exp3","IOFlexIO_2D")
end if
end if
if (ahf_3Doutput.ne.0) then
if (.not.find3) then
call CCTK_OutputVarByMethod(ierror,cctkGH,
. "ahfinder::ahfgrid","IOFlexIO_3D")
call CCTK_OutputVarByMethod(ierror,cctkGH,
. "ahfinder::ahf_exp","IOFlexIO_3D")
else if (mfind.eq.3) then
call CCTK_OutputVarByMethod(ierror,cctkGH,
. "ahfinder::ahfgrid","IOFlexIO_3D")
call CCTK_OutputVarByMethod(ierror,cctkGH,
. "ahfinder::ahf_exp","IOFlexIO_3D")
end if
end if
! **************************************************
! *** OUTPUT ASCII DATA ONLY ON PROCESSOR 0 ***
! **************************************************
if (myproc.eq.0) then
! Mass.
if (ahf_ncall.eq.1) then
open(1,file=massf,form='formatted',status='replace')
write(1,"(A14)") '" Horizon mass'
else
open(1,file=massf,form='formatted',status='old',position='append')
end if
if (status.and.report) then
out_mass = 0.141047396D0 * sqrt(intarea_h)
else
out_mass = 0.0
end if
write(1,"(2ES14.6)") cctk_time, out_mass
close(1)
! Radius.
if (ahf_ncall.eq.1) then
open(1,file=radf,form='formatted',status='replace')
write(1,"(A16)") '" Horizon radius'
else
open(1,file=radf,form='formatted',status='old',position='append')
end if
if (status.and.report) then
out_radius = c0(0)
else
out_radius = 0.0
end if
write(1,"(2ES14.6)") cctk_time, out_radius
close(1)
! Equatorial circumference.
if (ahf_ncall.eq.1) then
open(1,file=circeqf,form='formatted',status='replace')
write(1,"(A26)") '" Equatorial circumference'
else
open(1,file=circeqf,form='formatted',status='old',position='append')
end if
if (status.and.report) then
out_perimeter = circ_eq
else
out_perimeter = 0.0
end if
write(1,"(2ES14.6)") cctk_time, out_perimeter
close(1)
! Meridians.
if (ahf_ncall.eq.1) then
open(1,file=merip1f,form='formatted',status='replace')
write(1,"(A27)") '" Length of meridian, phi=0'
else
open(1,file=merip1f,form='formatted',status='old',position='append')
end if
if (status.and.report) then
out_meridian1 = meri_p1
else
out_meridian1 = 0.0
end if
write(1,"(2ES14.6)") cctk_time, out_meridian1
close(1)
if (ahf_ncall.eq.1) then
open(1,file=merip2f,form='formatted',status='replace')
write(1,"(A30)") '" Length of meridian, phi=pi/2'
else
open(1,file=merip2f,form='formatted',status='old',position='append')
end if
if (status.and.report) then
out_meridian2 = meri_p2
else
out_meridian2 = 0.0
end if
write(1,"(2ES14.6)") cctk_time, out_meridian2
close(1)
! Area.
if (ahf_ncall.eq.1) then
open(1,file=areaf,form='formatted',status='replace')
write(1,"(A14)") '" Horizon area'
else
open(1,file=areaf,form='formatted',status='old',position='append')
end if
if (status.and.report) then
out_area = intarea_h
else
out_area = 0.0
end if
write(1,"(2ES14.6)") cctk_time, out_area
close(1)
! Asymmetries in x. Coefficients cc with odd m, and
! coefficients cs with even m indicate asymmetries
! on refelction x<->-x.
if (ahf_ncall.eq.1) then
open(1,file=asymxf,form='formatted',status='replace')
write(1,"(A29)") '" Asymmetries on x reflection'
else
open(1,file=asymxf,form='formatted',status='old',position='append')
end if
out_asymx = 0.0D0
if (status.and.report.and..not.refx) then
do l=1,lmax
do m=1,l
if (mod(m,2).ne.0) then
out_asymx = out_asymx + abs(cc(l,m))
else
out_asymx = out_asymx + abs(cs(l,m))
end if
end do
end do
out_asymx = out_asymx / c0(0)
end if
write(1,"(2ES14.6)") cctk_time, out_asymx
close(1)
! Asymmetries in y. Any cs coefficient different from
! zero indicates asymmetries on reflection y<->-y.
if (ahf_ncall.eq.1) then
open(1,file=asymyf,form='formatted',status='replace')
write(1,"(A29)") '" Asymmetries on y reflection'
else
open(1,file=asymyf,form='formatted',status='old',position='append')
end if
out_asymy = 0.0D0
if (status.and.report.and..not.refy) then
do l=1,lmax
do m=1,l
out_asymy = out_asymy + abs(cs(l,m))
end do
end do
out_asymy = out_asymy / c0(0)
end if
write(1,"(2ES14.6)") cctk_time, out_asymy
close(1)
! Asymmetries in z. Any coefficient c0 with odd l,
! or any coefficient cc and cs with odd (l-m) indicate
! asymmetries on reflection z<->-z.
if (ahf_ncall.eq.1) then
open(1,file=asymzf,form='formatted',status='replace')
write(1,"(A29)") '" Asymmetries on z reflection'
else
open(1,file=asymzf,form='formatted',status='old',position='append')
end if
out_asymz = 0.0D0
if (status.and.report.and..not.refz) then
do l=1,lmax
if (mod(l,2).ne.0) out_asymz = out_asymz + abs(c0(l))
do m=1,l
if (mod(l-m,2).ne.0) then
out_asymz = out_asymz + abs(cc(l,m)) + abs(cs(l,m))
end if
end do
end do
out_asymz = out_asymz / c0(0)
end if
write(1,"(2ES14.6)") cctk_time, out_asymz
close(1)
! Expansion coefficients.
if (ahf_ncall.eq.1) then
open(1,file=almf,form='formatted',status='replace')
write(1,"(A21)") '# Radial coefficients'
write(1,"(A1)") '#'
else
open(1,file=almf,form='formatted',status='old',position='append')
write(1,*)
end if
if (status.and.report) then
out_centerx = xc
out_centery = yc
out_centerz = zc
write(1,"(A15,3ES14.6)") '# centered on: ', out_centerx, out_centery, out_centerz
endif
write(1,"(A12,I4)") '# Time step ',cctk_iteration
write(1,"(A7,ES14.6)") '# Time ',cctk_time
write(1,"(A7,I4)") '# Call ',ahf_ncall
if (status.and.report) then
if (mtype.eq.1) then
if (horizon) then
write(1,"(A30,I4)") '# Surface found: Outer horizon'
else
write(1,"(A31,I4)") '# Surface found: Outer horizon?'
end if
else if (mtype.eq.2) then
if (horizon) then
write(1,"(A30,I4)") '# Surface found: Inner horizon'
else
write(1,"(A31,I4)") '# Surface found: Inner horizon?'
end if
else if (mtype.eq.3) then
if (horizon) then
write(1,"(A34,I4)") '# Surface found: Marginal horizon?'
else
write(1,"(A32,I4)") '# Surface found: Trapped surface'
end if
else if (mtype.eq.4) then
if (horizon) then
write(1,"(A34,I4)") '# Surface found: Marginal horizon?'
else
write(1,"(A30,I4)") '# Surface found: Not a horizon'
end if
end if
write(1,"(A1)") '#'
write(1,"(A22)") '# a_lm l m'
write(1,"(A1)") '#'
write(1,"(ES14.6,2I4)") c0(0),0,0
if (lmax .gt. 0) then
out_c0 = c0(1:)
out_cs = cs
out_cc = cc
do l=1,lmax
do m=l,1,-1
write(1,"(ES14.6,2I4)") out_cs(l,m),l,-m
end do
write(1,"(ES14.6,2I4)") out_c0(l),l,0
do m=1,l
write(1,"(ES14.6,2I4)") out_cc(l,m),l,m
end do
end do
end if
else
write(1,"(A18,I4)") '# No surface found'
end if
close(1)
end if
! ***********************************************
! *** OUTPUT AHFINDER DATA IN HDF5 FORMAT ***
! ***********************************************
options = "mfind=" // char(48 + mfind) // " mtype=" // char(48 + mtype)
if (status) then
options = options(1:15) // " status=1"
else
options = options(1:15) // " status=0"
end if
options = options // " "
if (horizon) then
options = options(1:24) // " horizon=1"
else
options = options(1:24) // " horizon=0"
end if
if (ahf_HDF5output.ne.0) then
call CCTK_OutputVarByMethod (ierror, cctkGH, options, "IOAHFinderHDF5")
if (ierror .ne. 0) then
call CCTK_WARN (1, 'Error calling I/O method "IOAHFinderHDF5" (I/O method not activated ?)')
end if
end if
end subroutine AHFinder_Output
|