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
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
|
c/*@@
c @file Extract.F
c @date 28th October 1997
c @author Gab Allen
c @desc Waveform extraction
c
c @enddesc
c@@*/
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
c/*@@
c @routine Extract
c @date 28th October 1997
c @author Gab Allen
c @desc Entry point for waveform extraction routines
c
c @enddesc
c @calls D3_extract
c @calledby No idea
c @history
c
c @endhistory
c@@*/
SUBROUTINE Extract(CCTK_ARGUMENTS)
USE D3_extract_int
IMPLICIT NONE
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_FUNCTIONS
c Non-Cactus input variables for D3_extract
INTEGER ::
& igrid,lmode,mmode,Psi_power=1
INTEGER ::
& do_ADMmass(2)
CCTK_REAL ::
& orig(3),radius
CCTK_REAL ::
& x_1d(cctk_lsh(1)),y_1d(cctk_lsh(2)),z_1d(cctk_lsh(3))
c Output variables from D3_extract
CCTK_REAL ::
& dtaudt,ADMmass(2),mass,rsch,momentum(3),spin(3)
CCTK_REAL,ALLOCATABLE ::
& Qodd(:,:,:), Qeven(:,:,:)
c Local variables
CCTK_INT :: nchar
INTEGER ::
& ix,iy,iz,fn1,fn2,lmin,lmax,mmin,mmax,lstep,mstep,
& il,im,it,ndet,idet,ioutput
INTEGER,SAVE :: openfile = 1
INTEGER,PARAMETER ::
& max_detectors = 9,number_timecoords = 2
CCTK_REAL,PARAMETER ::
& two = 2.0D0
CCTK_REAL ::
& time,r1,r2,dr,r_max,xmin,xmax,ymin,ymax,zmin,zmax
CCTK_REAL ::
& detector(max_detectors)
CHARACTER*100 ::
& out1,out2,massfile,rschfile,ADMmassfile1,ADMmassfile2,
& momentumfile1,momentumfile2,momentumfile3,spinfile1,spinfile2,spinfile3
CHARACTER*3 ::
& timestring
LOGICAL :: use_proper,use_coordinate,do_output
CCTK_REAL,DIMENSION(10),SAVE :: proper_time
CCTK_REAL,SAVE :: last_time
INTEGER, save :: open_file_level(10) = 0
INTEGER myproc,ierr
CCTK_REAL dx,dy,dz
character*200 filestr
character*80 infoline
c ------------------------------------------------------------------
c Get the output directory sorted out
call CCTK_FortranString(nchar,outdir,filestr)
c ------------------------------------------------------------------
myproc = CCTK_MyProc(cctkGH)
dx = cctk_delta_space(1)
dy = cctk_delta_space(2)
dz = cctk_delta_space(3)
c ------------------------------------------------------------------
c
c 0. Check to see if Extract should have been called
c
c ------------------------------------------------------------------
IF ((itout) .LE. 0) RETURN
IF(MOD(int(cctk_iteration),int(itout)) .NE. 0) RETURN
if (verbose == 1) then
write(infoline,'(A24,G12.7)')
& 'Calling Extract at time ',cctk_time
call CCTK_INFO(infoline)
end if
c ------------------------------------------------------------------
c
c 1. Initial Stuff
c
c ------------------------------------------------------------------
c Get the number of polar divisions
IF (MOD(int(Nt),2) == 1 .OR. MOD(int(Np),2) == 1 .OR.
& Nt < 0 .OR. Np < 0) THEN
call CCTK_WARN(0,"Error in Nt or Np in Extract")
END IF
c Get the origin of spherical symmetry
orig(1) = origin_x
orig(2) = origin_y
orig(3) = origin_z
c Set the value of igrid
IF (CCTK_EQUALS(domain,"octant")) THEN
igrid = 1
ELSEIF (CCTK_EQUALS(domain,"bitant")) THEN
igrid = 3
ELSE
igrid = 0
END IF
#ifdef THORN_CARTOON_2D
IF (contains("cartoon_active","yes") .NE. 0) THEN
igrid = 2
Np = 1
END IF
#endif
c Create 1D coordinate arrays
do ix = 1, cctk_lsh(1)
x_1d(ix) = x(ix,1,1)
enddo
do iy = 1, cctk_lsh(2)
y_1d(iy) = y(1,iy,1)
enddo
do iz = 1, cctk_lsh(3)
z_1d(iz) = z(1,1,iz)
enddo
c See if the ADM mass should be calculated and output
IF (doADMmass == 1) THEN
do_ADMmass(1) = 1
IF (use_conformal == 1) THEN
do_ADMmass(2) = 1
ELSE
do_ADMmass(2) = 0
ENDIF
ELSE
do_ADMmass(:) = 0
END IF
c Set array containing gtt ... for the moment I am lazy and
c just make it the lapse and do not include the shift!
g00 = alp
c What kind of timecoordinate do I want to use
if (CCTK_EQUALS(timecoord,"both")) then
use_proper = .true.
use_coordinate = .true.
if (cctk_iteration == 0) then
proper_time = 0
last_time = cctk_time
end if
elseif (CCTK_EQUALS(timecoord,"proper")) then
use_proper = .true.
use_coordinate = .false.
if (cctk_iteration == 0) then
proper_time = 0
last_time = cctk_time
end if
elseif (CCTK_EQUALS(timecoord,"coord")) then
use_proper = .false.
use_coordinate = .true.
else
use_proper = .false.
use_coordinate = .true.
endif
c ------------------------------------------------------------------
c
c 2. Sort out which modes to extract
c
c ------------------------------------------------------------------
lmode = l_mode
IF (all_modes == 1) THEN
lmin = 2
lmax = lmode
c Get m_mode in this case too, or the T3E has a junk value (PW)
mmode = 0
ALLOCATE(Qodd(2,2:lmode,0:lmode),Qeven(2,2:lmode,0:lmode))
ELSE
lmin = lmode
lmax = lmode
mmode = m_mode
ALLOCATE(Qodd(2,2:lmode,0:mmode),Qeven(2,2:lmode,0:mmode))
ENDIF
c Check do not have bad values for the modes
IF (lmode < 2 .OR. mmode > lmode .OR. mmode < 0) THEN
call CCTK_WARN(0,"Error in lmode,mmode in Extract")
END IF
IF (lmode > 9) THEN
call CCTK_WARN(4,"Filenames will probably be crazy in Extract")
END IF
c If using an octant, do not do odd modes
IF (igrid == 1 .OR. igrid == 2) THEN
lstep = 2 ; mstep = 2
ELSE
lstep = 1 ; mstep = 1
END IF
c If we have cartoon, then jump past all the m-modes in loop
IF (igrid == 2) mstep = lmode+1
c ------------------------------------------------------------------
c
c 3. Find maximum radius of extraction on grid
c
c ------------------------------------------------------------------
call CCTK_CoordRange(ierr,cctkGH,xmin,xmax,-1,"x","cart3d")
call CCTK_CoordRange(ierr,cctkGH,ymin,ymax,-1,"y","cart3d")
call CCTK_CoordRange(ierr,cctkGH,zmin,zmax,-1,"z","cart3d")
IF (igrid == 2) THEN
r_max = MIN(xmax-two*Dx-orig(1),
& zmax-two*Dz-orig(3))
ELSE IF (igrid == 1) THEN
r_max = MIN(zmax-two*Dx-orig(1),
& ymax-two*Dy-orig(2),
& zmax-two*Dz-orig(3))
ELSE IF (igrid == 3) THEN
r_max = MIN(xmax-two*Dx-orig(1),
& ymax-two*Dy-orig(2),
& zmax-two*Dz-orig(3),
& DABS(xmin+two*Dx-orig(1)),
& DABS(ymin+two*Dy-orig(2)))
ELSE
r_max = MIN(xmax-two*Dx-orig(1),
& ymax-two*Dy-orig(2),
& zmax-two*Dz-orig(3),
& DABS(xmin+two*Dx-orig(1)),
& DABS(ymin+two*Dy-orig(2)),
& DABS(zmin+two*Dz-orig(3)))
END IF
c ------------------------------------------------------------------
c
c 4. Extract Cauchy initial data for a linear wave equation
c
c ------------------------------------------------------------------
test_Cauchy: IF (Cauchy == 1) THEN
it = Cauchy_timestep
c check_timestep: IF (cctk_iteration == it) THEN
check_timestep: IF (open_file_level(1) .eq. 0) THEN
c Open output files
write (*,*) "Open output file on level ", 1
open_file_level(1)=1;
test_myproc1: IF (CCTK_MyProc(cctkGH) == 0 ) THEN
out1 = filestr(1:nchar)//"/rsch_ini.rl"
call Extract_open(cctkGH,0,out1,77)
out1 = filestr(1:nchar)//"/mass_ini.rl"
call Extract_open(cctkGH,0,out1,78)
IF (do_ADMmass(1) == 1) THEN
out1 = filestr(1:nchar)//"/ADMmass_ini.rl"
call Extract_open(cctkGH,0,out1,79)
ENDIF
IF (do_ADMmass(2) == 1) THEN
out1 = filestr(1:nchar)//"/ADMmassc_ini.rl"
call Extract_open(cctkGH,0,out1,80)
ENDIF
IF (do_momentum == 1) THEN
out1 = filestr(1:nchar)//"/momentum_x_ini.rl"
call Extract_open(cctkGH,0,out1,781)
out1 = filestr(1:nchar)//"/momentum_y_ini.rl"
call Extract_open(cctkGH,0,out1,782)
out1 = filestr(1:nchar)//"/momentum_z_ini.rl"
call Extract_open(cctkGH,0,out1,783)
ENDIF
IF (do_spin == 1) THEN
out1 = filestr(1:nchar)//"/spin_x_ini.rl"
call Extract_open(cctkGH,0,out1,784)
out1 = filestr(1:nchar)//"/spin_y_ini.rl"
call Extract_open(cctkGH,0,out1,785)
out1 = filestr(1:nchar)//"/spin_z_ini.rl"
call Extract_open(cctkGH,0,out1,786)
ENDIF
loop_l1: DO il = lmin,lmax,lstep
IF (all_modes == 0) THEN
mmin = mmode ; mmax = mmode
ELSE
mmin = 0 ; mmax = il
END IF
loop_m1: DO im = mmin,mmax,lstep
fn1 = il*10+im
out1 = filestr(1:nchar)//"/Qeven_ini_"//
& CHAR(il+48)//CHAR(im+48)//".rl"
call Extract_open(cctkGH,0,out1,fn1)
c Only print odd-parity if full grid
IF (igrid == 0) THEN
fn2 = 100+il*10+im
out2 = filestr(1:nchar)//"/Qodd_ini_"//
& CHAR(il+48)//CHAR(im+48)//".rl"
call Extract_open(cctkGH,0,out2,fn2)
END IF
END DO loop_m1
END DO loop_l1
END IF test_myproc1
c Find range of extraction radii
r1 = Cauchy_r1
r2 = r_max
dr = Cauchy_dr
c Do extraction at each radius
IF (verbose == 1) THEN
WRITE(*,*) "Extracting Cauchy initial data"
WRITE(*,*) " r = ",r1," to ",r2," step ",dr
END IF
radius = r1
extract_at_each_radius: DO WHILE (radius < r2)
CALL D3_extract(cctkGH,do_ADMmass,do_momentum,do_spin,
& igrid,orig,myproc,Nt,Np,all_modes,lmode,
& mmode,x_1d,y_1d,z_1d,Dx,Dy,Dz,Psi_power,Psi,g00,
& gxx,gxy,gxz,gyy,gyz,gzz,kxx,kxy,kxz,kyy,kyz,kzz,
& radius,ADMmass,momentum,spin,mass,rsch,
& Qodd,Qeven,temp3d,dtaudt)
IF (verbose == 1) THEN
WRITE(*,*) "Extracted at r =",radius
WRITE(*,*) " Sch radius/mass =",rsch,mass
ENDIF
c Write to file at each radius
test_myproc2: IF (CCTK_MyProc(cctkGH) == 0) THEN
WRITE( 77,*) radius,rsch
WRITE( 78,*) radius,mass
IF (do_ADMmass(1) == 1) WRITE( 79,*) radius,ADMmass(1)
IF (do_ADMmass(2) == 1) WRITE( 80,*) radius,ADMmass(2)
IF (do_momentum == 1) then
WRITE(781,*) radius,momentum(1)
WRITE(782,*) radius,momentum(2)
WRITE(783,*) radius,momentum(3)
end if
IF (do_spin == 1) THEN
WRITE(784,*) radius,spin(1)
WRITE(785,*) radius,spin(2)
WRITE(786,*) radius,spin(3)
END IF
loop_l2: DO il = lmin,lmax,lstep
IF (all_modes == 0) THEN
mmin = mmode ; mmax = mmode
ELSE
mmin = 0 ; mmax = il
END IF
loop_m2: DO im = mmin,mmax,mstep
fn1 = il*10+im
WRITE(fn1,*) rsch,Qeven(:,il,im)
c Only print odd-parity if full grid
IF (igrid == 0) THEN
fn2 = 100+il*10+im
WRITE(fn2,*) rsch,Qeven(:,il,im)
END IF
END DO loop_m2
END DO loop_l2
END IF test_myproc2
radius = radius+dr
END DO extract_at_each_radius
c Now close all the files
test_myproc3: IF (myproc == 0) THEN
CLOSE( 77) ! Schwarzschild radius
CLOSE( 78) ! Schwarzschild mass
IF (do_ADMmass(1) == 1) CLOSE( 79) ! ADM mass
IF (do_ADMmass(2) == 1) CLOSE( 80) ! ADM mass
IF (do_momentum == 1) THEN
CLOSE(781) ! momentum
CLOSE(782) ! momentum
CLOSE(783) ! momentum
end if
IF (do_spin == 1) THEN
CLOSE(784) ! spin
CLOSE(785) ! spin
CLOSE(786) ! spin
END IF
loop_l3: DO il = lmin,lmax,lstep
IF (all_modes == 0) THEN
mmin = mmode ; mmax = mmode
ELSE
mmin = 0 ; mmax = il
END IF
loop_m3: DO im = mmin,mmax,mstep
fn1 = il*10+im
CLOSE(fn1) ! Qeven_ini_lm.dat
IF (igrid == 0) THEN
fn2 = 100+il*10+im
CLOSE(fn2) ! Qodd_ini_lm.dat
END IF
END DO loop_m3
END DO loop_l3
END IF test_myproc3
END IF check_timestep
END IF test_Cauchy
c ------------------------------------------------------------------
c
c 5. Extract waveforms at detectors
c
c ------------------------------------------------------------------
c Cannot use the conformal equation for ADM mass now
do_ADMmass(2) = 0
ndet = num_detectors
test_detectors: IF (ndet > 0) THEN
IF (ndet > 9) THEN
call CCTK_WARN(0,"Too many detectors in Extract")
END IF
detector(1) = detector1
detector(2) = detector2
detector(3) = detector3
detector(4) = detector4
detector(5) = detector5
detector(6) = detector6
detector(7) = detector7
detector(8) = detector8
detector(9) = detector9
DO idet = ndet+1, max_detectors
detector(idet) = 0.0D0
END DO
DO idet = 1, ndet
IF (detector(idet) > r_max) THEN
IF (openfile == 1 .OR. cctk_iteration == it) THEN
call CCTK_WARN(8,"Removing detectors outside coordinate range")
END IF
ndet =idet-1
GOTO 404
ENDIF
END DO
404 CONTINUE
IF (verbose == 1) THEN
IF (openfile == 1) THEN
WRITE(*,*) "Extracting at ",ndet," detectors "
DO idet=1,ndet
WRITE(*,*) " r = ",detector(idet)
ENDDO
END IF
END IF
detector_loop: DO idet = 1,ndet
radius = detector(idet)
IF (verbose == 1) THEN
WRITE(*,*) "Calling extract at radius ... ",radius
END IF
CALL D3_extract(cctkGH,do_ADMmass,do_momentum,do_spin,
& igrid,orig,myproc,Nt,Np,all_modes,lmode,mmode,x_1d,y_1d,z_1d,
& Dx,Dy,Dz,Psi_power,Psi,g00,gxx,gxy,gxz,gyy,gyz,gzz,kxx,kxy,kxz,
& kyy,kyz,kzz,radius,ADMmass,momentum,spin,mass,rsch,
& Qodd,Qeven,temp3d,dtaudt)
IF (verbose == 1) THEN
WRITE(*,*) "Detector ",idet," ..."
WRITE(*,*) " Sch radius/mass =",rsch,mass
END IF
c Output to files
test_myproc4: IF (myproc == 0) THEN
DO ioutput=1,number_timecoords
do_output = .false.
timestring = "ERR"
time = 0D0
IF (ioutput == 1 .AND. use_coordinate) THEN
timestring = ".tl"
do_output = .true.
time = cctk_time
ELSEIF (ioutput == 2 .AND. use_proper) THEN
timestring = ".ul"
do_output = .true.
time = proper_time(idet) + (cctk_time-last_time)*dtaudt
proper_time(idet) = time
if (idet == ndet) last_time = cctk_time
ENDIF
IF (do_output) THEN
c Output extracted radius and mass
rschfile = filestr(1:nchar)//"/rsch_"//"R"//
& CHAR(idet+48)//timestring
massfile = filestr(1:nchar)//"/mass_"//"R"//
& CHAR(idet+48)//timestring
ADMmassfile1 = filestr(1:nchar)//"/ADMmass_"//"R"//
& CHAR(idet+48)//timestring
ADMmassfile2 = filestr(1:nchar)//"/ADMmassc_"//"R"//
& CHAR(idet+48)//timestring
momentumfile1 = filestr(1:nchar)//"/momentum_x_"//"R"//
& CHAR(idet+48)//timestring
momentumfile2 = filestr(1:nchar)//"/momentum_y_"//"R"//
& CHAR(idet+48)//timestring
momentumfile3 = filestr(1:nchar)//"/momentum_z_"//"R"//
& CHAR(idet+48)//timestring
spinfile1 = filestr(1:nchar)//"/spin_x_"//"R"//
& CHAR(idet+48)//timestring
spinfile2 = filestr(1:nchar)//"/spin_y_"//"R"//
& CHAR(idet+48)//timestring
spinfile3 = filestr(1:nchar)//"/spin_z_"//"R"//
& CHAR(idet+48)//timestring
call Extract_write(cctkGH,openfile,rschfile,time,rsch)
call Extract_write(cctkGH,openfile,massfile,time,mass)
IF (do_ADMmass(1) == 1) THEN
call Extract_write(cctkGH,openfile,ADMmassfile1,time,ADMmass(1))
END IF
IF (do_ADMmass(2) == 1) THEN
call Extract_write(cctkGH,openfile,ADMmassfile2,time,ADMmass(2))
END IF
IF (do_momentum == 1) THEN
call Extract_write(cctkGH,openfile,momentumfile1,time,momentum(1))
call Extract_write(cctkGH,openfile,momentumfile2,time,momentum(2))
call Extract_write(cctkGH,openfile,momentumfile3,time,momentum(3))
END IF
IF (do_spin == 1) THEN
call Extract_write(cctkGH,openfile,spinfile1,cctk_time,spin(1))
call Extract_write(cctkGH,openfile,spinfile2,cctk_time,spin(2))
call Extract_write(cctkGH,openfile,spinfile3,cctk_time,spin(3))
END IF
c Output gauge invariant variables
loop_l5: DO il = lmin,lmax,lstep
IF (all_modes == 0) THEN
mmin = mmode ; mmax = mmode
ELSE
mmin = 0 ; mmax = il
END IF
loop_m5: DO im = mmin,mmax,mstep
out1 = filestr(1:nchar)//"/Qeven_"//"R"
& //CHAR(idet+48)//"_"//CHAR(il+48)//
& CHAR(im+48)//timestring
out2 = filestr(1:nchar)//"/Qodd_"//"R"
& //CHAR(idet+48)//"_"//CHAR(il+48)//
& CHAR(im+48)//timestring
c Write even parity waveforms
call Extract_write2(cctkGH,openfile,out1,time,Qeven(:,il,im))
c Only write odd parity waveforms if full grid
IF (igrid == 0) THEN
call Extract_write2(cctkGH,openfile,out2,time,Qodd(:,il,im))
END IF
END DO loop_m5
END DO loop_l5
ENDIF
END DO
END IF test_myproc4
END DO detector_loop
END IF test_detectors
DEALLOCATE(Qodd,Qeven)
openfile = 0
20 format(e24.15,e24.15,e24.15)
21 format(e24.15,e24.15)
END SUBROUTINE Extract
SUBROUTINE Extract_open(cctkGH,openfile,filename,filehandle)
implicit none
DECLARE_CCTK_PARAMETERS
integer filehandle,openfile
character*(*) filename
CCTK_POINTER cctkGH
if (openfile==1) then
OPEN(UNIT= filehandle,FILE=filename,STATUS="unknown")
call Extract_Advertise(cctkGH,filename)
else
OPEN(UNIT=filehandle,FILE=filename,STATUS="old",POSITION="append")
end if
END SUBROUTINE Extract_open
SUBROUTINE Extract_write(cctkGH, openfile,filename,value1,value2)
implicit none
DECLARE_CCTK_PARAMETERS
CCTK_POINTER cctkGH
CCTK_REAL value1,value2
integer openfile
character*(*) filename
if (openfile==1) then
OPEN(UNIT= 444,FILE=filename,STATUS="unknown")
call Extract_Advertise(cctkGH,filename)
else
OPEN(UNIT=444,FILE=filename,STATUS="old",POSITION="append")
end if
write(444,*) value1,value2
close(444)
END SUBROUTINE Extract_write
SUBROUTINE Extract_write2(cctkGH, openfile,filename,value1,value2)
implicit none
DECLARE_CCTK_PARAMETERS
CCTK_POINTER cctkGH
CCTK_REAL value1
CCTK_REAL value2(2)
integer openfile
character*(*) filename
if (openfile==1) then
OPEN(UNIT= 444,FILE=filename,STATUS="unknown")
call Extract_Advertise(cctkGH,filename)
else
OPEN(UNIT=444,FILE=filename,STATUS="old",POSITION="append")
end if
write(444,*) value1,value2(1),value2(2)
close(444)
END SUBROUTINE Extract_write2
|