-
Notifications
You must be signed in to change notification settings - Fork 2
/
AR-CHAIN.f
2835 lines (2660 loc) · 87.6 KB
/
AR-CHAIN.f
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
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
c THIS CODE USES THE ALGORITHMIC REGULARIZATION i.e. no coordinate transformation, only time transformation + leapfrog + extrapolation method
! NearFinalVersionOf_ARC.f
subroutine ARCparams_Dot_CH_is_here ! This subroutine is not used, but the contents is simply what is supposed to be in the file ARCparams.CH.
! If you do not have that file, create it and copy the content of this,
! without the 'end' statement, into it! Then you can compile this code.
! ARCparams.CH =filename
IMPLICIT REAL*8 (A-H,M,O-Z)
PARAMETER (NMX=200,NMX2=2*NMX,NMX3=3*NMX,NMX4=4*NMX,
& NMX8=8*NMX,NMXm=NMX*(NMX-1)/2)
COMMON/DataForRoutines1/X(NMX3),V(NMX3),WTTL,M(NMX),
& XC(NMX3),WC(NMX3),MC(NMX)
& ,XI(NMX3),VI(NMX3),MASS,RINV(NMXm),RSUM,INAME(NMX),N
COMMON/DataForChainRoutinesTwo/MMIJ,CMX(3),CMV(3)
& ,ENERGY,Energr,CHTIME
common/softening/ee,cmethod(3),Clight,NofBH
common/TIMECOMMON/Taika,timecomparison,steppi
common/spincommon/spin(3)! the relative spin of M(1) !Spin=spin*G*M^2/c
common/tolerancecommon/EPS
end
! The following publications contain info needed for understanding this code
! \bibitem[Hellstr{\"o}m and Mikkola(2010)]{2010CeMDA.106..143H} Hellstr{\"o}m, C., Mikkola, S.\ 2010.\ Explicit algorithmic regularization in the few-body problem for velocity-dependent perturbations.\ Celestial Mechanics and Dynamical Astronomy 106, 143-156.
! \bibitem[Mikkola and Merritt(2008)]{2008AJ....135.2398M} Mikkola, S., Merritt, D.\ 2008.\ Implementing Few-Body Algorithmic Regularization with Post-Newtonian Terms.\ The Astronomical Journal 135, 2398-2405.
! \bibitem[Mikkola and Merritt(2006)]{2006MNRAS.372..219M} Mikkola, S., Merritt, D.\ 2006.\ Algorithmic regularization with velocity-dependent forces.\ Monthly Notices of the Royal Astronomical Society 372, 219-223.
! \bibitem[Mikkola and Aarseth(2002)]{2002CeMDA..84..343M} Mikkola, S., Aarseth, S.\ 2002.\ A Time-Transformed Leapfrog Scheme.\ Celestial Mechanics and Dynamical Astronomy 84, 343-354.
! \bibitem[Mikkola and Aarseth(1996)]{1996CeMDA..64..197M} Mikkola, S., Aarseth, S.~J.\ 1996.\ A Slow-down Treatment for Close Binaries.\ Celestial Mechanics and Dynamical Astronomy 64, 197-208.
! \bibitem[Mikkola and Aarseth(1993)]{1993CeMDA..57..439M} Mikkola, S., Aarseth, S.~J.\ 1993.\ An implementation of N-body chain regularization.\ Celestial Mechanics and Dynamical Astronomy 57, 439-459.
! \bibitem[Mikkola and Tanikawa(2013)]{2013NewA...20...38M} Mikkola, S., Tanikawa, K.\ 2013.\ Implementation of an efficient logarithmic-Hamiltonian three-body code.\ New Astronomy 20, 38-41.
! \bibitem[Mikkola and Tanikawa(2013)]{2013MNRAS.430.2822M} Mikkola, S., Tanikawa, K.\ 2013.\ Regularizing dynamical problems with the symplectic logarithmic Hamiltonian leapfrog.\ Monthly Notices of the Royal Astronomical Society 430, 2822-2827.
! \bibitem[Mikkola and Tanikawa(1999)]{1999MNRAS.310..745M} Mikkola, S., Tanikawa, K.\ 1999.\ Algorithmic regularization of the few-body problem.\ Monthly Notices of the Royal Astronomical Society 310, 745-749.
Program ARCcode ! This is the main program, not subroutine.
! IMPLICIT REAL*8 (A-H,M,O-Z) ! this is in 'ARCparams.CH'
include 'ARCparams.CH'
COMMON/DIAGNOSTICS/GAMMA,H,IWR
common/justforfun/Tkin,Upot,dSkin,dSpot ! used in diagno
common/outputindex/index4output(200),N_ini
!REAL*8 G0(3),G(3),xw(3),vw(3)
REAL*8 cmet(3),xwr(NMX3)!,dum(3) !& ,ai(NMX),ei(NMX),inc(NMX),Omi(NMX),ooi(NMX)
& ,cmxx(3),cmvx(3)
LOGICAL NEWREG
CHARACTER*22 outfile
common/vindex/ivelocity
common/collision/icollision,ione,itwo,iwarning
! UNITS: I use G=1 (all other things are up to the user)
! BUT: I always use M_sun=1, lenght unit= 1AU.
! IN SUCH A SYSTEM: time is such that 1 year = 2*Pi
! C (vel. of light) is approximately = 10,000. (c=299792458 m / s, AU=149597870.7km ==> c=10065.1 au/t1)
666 CONTINUE ! jump here to start a new simulation in the same run.
icollision=0
iopen=0
TIME=0 ! initialization
SP0=0
DELTAT=1 ! something (but read below)
IWR=1 ! write some info ( set -1 to not to get that..), read below
! Initial values, UNITS: G=1
READ(5,*,err=999)IWR,N,DELTAT,TMAX,stepr,soft,cmet,Clight,!
& outfile,Ixc,Nbh ,spin,EPS,ivelocity,KSMX
! stepr is now obsolete
! Look the example cdr* file. It contains some info about the data read here.
! (look to the end part of the data file!)
if(N.lt.2)STOP ! if N<2 this code is not needed
ee=soft**2 ! square of softening lenght (often zero, but can be used)
open(66,file=outfile) ! output file
open(67,file='merge.dat') ! output for merger info
MASS=0
DO I=1,N
L=3*(I-1)
READ(5,*)M(I),(X(L+K),K=1,3),(V(L+K),K=1,3) ! Read masses, coordinates anf velocities
MASS=MASS+M(I) ! determine total mass
index4output(i)=i ! initialize output index (to be modified in case of mergers)
test=M(I)+cdot(x(L+1),X(L+1))+cdot(V(L+1),V(L+1))
if(test.eq.0.0)go to 2
END DO
2 N=I-1
!write(6,*)I,N,' = I N'
N_ini=N
call Reduce2cm(x,m,N,cmxx) ! x->cm-system, in this version output is in cm anyway
call Reduce2cm(v,m,N,cmvx) ! v->cm-system
do k=1,3
! cmxx(k)=0 ! cm coords ->0
! cmvx(k)=0 ! cm vels ->0
end do
ENER0=0
NEWREG=.TRUE. ! Now we are starting a new regularized integration
!KSMX=10 000 ! only this many steps without return, KSMX is actually read, normally big value recommended
goto 200 ! Go to write initial quantities (i.e. at time=0), after which code returns to statement 100
c
100 CONTINUE
call ARC !
& (N,X,V,M,TIME,DELTAT,EPS,NEWREG,KSMX,soft,cmet,clight,Ixc,NBH,
& spin,CMXX,CMVX) ! Here you get cm-coords (=Xa) aind cm-vels (=V). CMXX and CMCX are position and velocity of the centre-of-mass
! At this point U can do your own analysis/output etc....
200 CONTINUE
call Diagnostic output(time,xwr,cmxx) ! write errors and coords to the given outputfile (name read to variable 'outfile')
! some more coords are written in the routine MERGE_I1_I2 in case whén a BH swollows an other body
if(iwr.ge.0) call Write Elements(time,iopen) ! with respect to M(1), to files aexes.dat, eccs.dat, incs.dat, Omes.dat,omes.dat,
if(iwr.gt.1)call FIND BINARIES(time) ! this is usually unimportant, but may give info about binary formation
234 format(1x,f18.6,1p,600g13.5)
! rs=2*m(1)/clight**2 ! Radius of event horizon of M1
IF(TIME.LT.TMAX)then ! continue if not yet at maximum time
GOTO 100
else
goto 666 ! go read data for the next experiment (to stop set a line of zeros into the data file)
end if
999 END
subroutine Diagnostic output(time,xwr,cmxx)
include 'ARCparams.CH'
COMMON/DIAGNOSTICS/GAMMA,H,IWR
common/justforfun/Tkin,Upot,dSkin,dSpot ! used in diagno
common/outputindex/index4output(200),N_ini
REAL*8 G0(3),G(3),xwr(NMX3),CMXX(3)!,CMVX(3)
common/vindex/ivelocity
common/collision/icollision,ione,itwo,iwarning
save
if(time.eq.0.0)then
iENER0=0
goto 200
END IF
c Diagnostics !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
call CONSTANTS OF MOTION(ENER1,G,AL)
IF(ENER0.eq.0.0)THEN
c This CONTANTS OF MOTION works only after first cal of ARC ! Therefore this is here
ENER0=ENER1
g0(1)=g(1)
g0(2)=g(2)
g0(3)=g(3)
cangmo=mass**2.5d0*sqrt(Al)/abs(Ener0)
end if
am_error=sqrt(square(g,g0))/cangmo ! in case of perturbations this is not necessarily error but just a measure of change
WRITE(6,123)TIME!
& ,log((Tkin-ENERGY)/Upot),log(dSkin/dSpot),am_error
& ,ENERGY,ENERGR!logH = the primary constant (=0!)
! ENERGY and EnerGR (=Relatistic change of ENERGY) are integrated (due to possible external perturbations)
& ,N ! print time, logH, N (of bodies left)
call flush(6)
123 FORMAT(1x,'T: ',1p,g20.6,' log((T-E)/U)=',1p,g10.2,
& ' log(dSdot/sdotU)=',1p,g10.2,
& ' d(RxV)/Am=',1p,3g15.7,
& ' Nb=',0p,1x,i3)
200 continue
if(iwr.gt.-2)then
do i=1,3*N_ini
xwr(i)=1.d9 ! put all outside
end do
do i=1,N
j=index4output(i) ! take still existing particles 2 correct indicies
j0=3*j-3
i0=3*i-3
do k=1,3
xwr(j0+k)=x(i0+k)+cmxx(k) ! add centre-of-mass (remove cmxx if U need cm coords)
! or replace cmxx(k) by -x(k) if U want M1 to be
! origin (se also subroutine MERGE_I1_I2)
end do ! k
end do !i
write(66,234)time,(xwr(k),k=1,3*n_ini)
c & (xwr(k)-xwr(1),xwr(k+1)-xwr(2),xwr(k+2)-xwr(3),k=1,3*n_ini,3) ! Write coords to 66,
! in case of a collision write more in Merge_i1_i2
! this is just for figs and/or movies
call flush(66)
end if ! iwr.gt.-2
234 format(1x,f18.6,1p,600g13.5)
return
end
subroutine Write Elements(time,iopen) ! with respect to M(1)
! IMPLICIT REAL*8 (A-H,M,O-Z) ! this is in 'ARCparams.CH'
include 'ARCparams.CH'
COMMON/DIAGNOSTICS/GAMMA,H,IWR
common/justforfun/Tkin,Upot,dSkin,dSpot ! used in diagno
common/outputindex/index4output(200),N_ini
REAL*8 !G0(3),G(3), cmet(3),
& xw(3),vw(3)!,xwr(NMX3)!,dum(3)
& ,ai(NMX),ei(NMX),inc(NMX),Omi(NMX),ooi(NMX)!,cmxx(3),cmvx(3)
!LOGICAL NEWREG
!CHARACTER*22 outfile
common/vindex/ivelocity
common/collision/icollision,ione,itwo,iwarning
save
do j=2,N_ini! Needed if N_ini>N
ai(j)=0
ei(j)=0
inc(j)=0
Omi(j)=0
ooi(j)=0
end do
do i=2,N
i0=3*i-3
do k=1,3
xw(k)=x(i0+k)-x(k)
vw(k)=v(i0+k)-v(k)
end do
mw=m(1)+m(i)
! Orbital elements with respect to the central body evaluated here.
j=index4output(i)
call elmnts
& (xw,vw,mw,ai(j),ei(j),moi,inc(j),Omi(j),ooi(j),alfai,qi,tqi)
end do ! i=2,N
! in case of a collision write more in Merge_i1_i2
! this is just for figs and/or movies
if(iopen.eq.0 )then
iopen=1
open(71,file='axes.dat')
open(72,file='eccs.dat')
open(73,file='incs.dat')
open(74,file='Omes.dat')
open(75,file='omes.dat')
open(76,file='spin.dat')
end if
!
!Write Keplerian orbital elements (with respect to M1)
write(71,171)time,(ai(k),k=2,N_ini) ! a
write(72,171)time,(ei(k),k=2,N_ini) ! e
write(73,171)time,(inc(k),k=2,N_ini) ! i
write(74,171)time,(Omi(k),k=2,N_ini) ! \Omega
write(75,171)time,(ooi(k),k=2,N_ini) ! \omega
spa=sqrt(cdot(spin,spin))
if(sp0.eq.0)sp0=spa
dsp=spa-sp0
write(76,*)time,spin,dsp ! spin(k), k=1,3 of M1 (|spin|<1),
! dsp is error in the length of the spin vector
171 format(1x,f12.3,201g18.10)
call flush(71)
call flush(72)
call flush(73)
call flush(74)
call flush(75)
call flush(76)
return
end
Subroutine MERGE_I1_I2(time)!Merge the colliding body with the BH, ('time' is here only for write(66....)
include 'ARCparams.CH'
REAL*8 SM(NMX),XR(NMX3),XDR(NMX3),xwr(nmx3),ywr(nmx3)
COMMON/collision/icollision,Ione,Itwo,iwarning
common/outputindex/index4output(200),N_ini
COMMON/DIAGNOSTICS/GAMMA,H,IWR
save
! --------------------------------------------4 MOVIE -------------------------
if(iwr.gt.-2)then ! Set IWR=-2 if you do not want this output
do i=1,3*N_ini
xwr(i)=1.d9
end do
do i=1,N
j=index4output(i)
j0=3*j-3
i0=3*i-3
do k=1,3
xwr(j0+k)=x(i0+k)+cmx(k)
ywr(j0+k)=xwr(j0+k)
if(i.gt.3 .and. (i.eq.ione .or. i.eq.itwo))ywr(j0+k)=1.d9 ! this moves the body out of movie screen(?)
end do ! k
end do !i
ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
do i_delay=1,33!Make the movie FLASH in case of a MERGER. (thefore 33 outputs, each other time out of he movie frame)
write(66,234)time,
& (xwr(k)-xwr(1),xwr(k+1)-xwr(2),xwr(k+2)-xwr(3),k=1,3*n_ini,3)
write(66,234)time,
& (ywr(k)-ywr(1),ywr(k+1)-ywr(2),ywr(k+2)-ywr(3),k=1,3*n_ini,3)
234 format(1x,f18.6,1p,600g13.5)
end do
call flush(66)
cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
end if ! iwr.gt.-2
!--------------------------------------------------------------------------------
DO I=1,ione-1
SM(I)=M(I)
DO K=1,3
XR(3*I-3+K)=X(3*I-3+K)
XDR(3*I-3+K)=V(3*I-3+K)
end do
end do
Myks=M(ione)
Mkax=M(itwo)
SM(Ione)=M(Ione)+M(Itwo)
DO 6 K=1,3
XR(3*Ione-3+K)=(M(Ione)*X((Ione-1)*3+K)
& +M(Itwo)*X((Itwo-1)*3+K))/SM(Ione)
XDR(3*Ione-3+K)=(M(Ione)*V((Ione-1)*3+K)
& +M(Itwo)*V((Itwo-1)*3+K))/SM(Ione)
6 CONTINUE
DO I=Ione+1,Itwo-1
sm(i)=m(i)
do k=1,3
XR(3*I-3+K)=X(3*I-3+k)
XDR(3*I-3+K)=V(3*I-3+k)
end do
end do
do i=Itwo,N-1
index4output(i)=index4output(i+1)
end do
DO I=Itwo+1,N
sm(i-1)=m(i)
do k=1,3
XR(3*I-6+K)=X(3*I-3+k)
XDR(3*I-6+K)=V(3*I-3+k)
end do
end do
! MOVE THE REDUCED SYSTEM TO M,X,V
! New value of the number of bodies.
N=N-1
if(Itwo.le.NofBH)NofBH=NofBH-1 ! # of BH's reduced!
DO 8 I=1,N
M(I)=SM(I)
DO 7 K=1,3
X(3*i-3+k)=XR(3*i-3+k)
V(3*i-3+k)=XDR(3*i-3+k)
7 CONTINUE
8 CONTINUE
icollision=0
i1wr=index4output(ione)
i2wr=index4output(itwo) !?? wrong ?? because already changed (above)
write(6,*)' Merge:',i1wr,i2wr-1,' (N, NBH=) (',N,' , ',NofBH, ') '
& ,' masses ',(M(k),k=1,n),Myks,Mkax
write(67,*)' At time ',time,' bodies number ',ione,
&' and ',itwo ,' with masses ',Myks,Mkax,' merge '
call flush(67)
ione=0
itwo=0
if(N.eq.1)then
write(6,*)' Only one body left!'
STOP
end if
RETURN
END
! * * *
SUBROUTINE ARC ! This is the original Alforithmic Regularization Chain code.
& (NN,XX,VX,MX,TIME,DELTAT,TOL,NEWREG,KSMX,soft,cmet,cl,Ixc,NBH,
& spini,CMXX,CMVX)
! BETTER TO USE CM-coords & vels for XX & VX and CMXX CMVX
! FOR CM-position (needed in the Perturbations routine).
!-----------------------------------------------------------------
! NOTE: some variables (eg. Energy and EnerGR are only in the
! common. The internal NB-energy = ENERGY (should be)
! Energy= integrated E-value (excluding grav.radiation)
! EnerGr= Energy radiated away (grav.radiation if Clight.ne.0.0)
! CHAIN INTEGRATION. Perturbations & CM-motion included (in principle).
! NN=# of bodies; XX=(cm)coords, VX=(cm)vels, MX=masses,
!c CMXX=coords of CM, CMVX=vels of CM ! removed
! TIME=time, dettaT=time interval
! STEP=stepsize (set=0 initially)
! NEWREG=.true. iff chain membership has changed
! KSMX=max # of steps without return (use some large # )
! soft =optional softening( U=1/sqrt(r**2+soft**2) )
! cmet= 3-d vector that determines the method:
! (1,0,0) =logH, (0,1,0)=TTL,(0,0,1)=DIFSY2 without t-tranformation
!
! cl=speed of light
! NOTE: cl=0 => no relativistic terms !!!
! Ixc = 2 =>iteration to excat time, 1 => approximate exact time, =0 no exact time but return after CHTIME>DELTAT
! often IXC =0 is fastest, IXC=1 second fastest, IXC=2 can be slow (but usually accurate output time).
INCLUDE 'ARCparams.CH'
COMMON/DerOfTime/GTIME
COMMON/DIAGNOSTICS/GAMMA,H,IWR
common/omegacoefficients/OMEC(NMX,NMX)
common/collision/icollision,ione,itwo,iwarning
common/itemaxcommon/aitemax,itemax,itemax_used
common/turhia/rw,fr,frm,akiih(3)
REAL*8 G0(3),XX(*),VX(*),MX(*),cmet(3),spini(3),CMXX(3),CMVX(3)
REAL*8 Y(1500),SY(1500),Yold(1500)
LOGICAL MUSTSWITCH,NEWREG
data ntrue,nfalse,nwritten/3*0/
save
! Initial constants of motion
10 CONTINUE !
tnext0=time+deltat
knx=tnext0/deltat+0.1d0
tnext=knx*deltat
tstep=tnext-time
if(newreg)then
ntrue=ntrue+1
else
nfalse=nfalse+1
end if
if(ntrue.gt.nfalse+10 .and. nwritten.eq.0)then
nwritten=1
write(6,*)char(7),char(7)
write(6,*)' NEWREG should be set .TRUE. only'
write(6,*)' in the very beginning of a new simulation'
write(6,*)' NOT at every step!! (May reduce accuracy!!)'
write(6,*)' even if it may look like the contrary.'
end if
if(NN.gt.NMX)then
write(6,*)' THIS CODE CAN HANDLE ONLY ',NMX,' BODIES '
write(6,*)' Yuo are trying to use N=',NN
write(6,*)' Try increasing NMX in ARCparams.CH '
write(6,*)' and increase some (large) dimensions elsewhere'
write(6,*)' in the same proportion. STOPPING'
STOP
end if
! if(cmet(1).eq.0.0 .and. cmet(2).ne.0.0)then
! write(6,*)' In this version cmethod(1) should not be zero'
! write(6,*)' if cmethod(2).ne.0.0 '
! write(6,*)cmet,' = cmethod(k),k=1,3 '
! write(6,*)' STOPPING '
! STOP
! end if
if(deltat.eq.0 .and. Ixc .eq.1)then
write(6,*)' You cannot use DELTA=0 and Ixc=1 '
write(6,*)' since then every output will be at time=0 '
write(6,*)' STOPPING '
STOP
end if
if(cmet(1)+cmet(2)+cmet(3).eq.0)then
write(6,*)' You have not defined the time-transformation'
write(6,*)cmet,' = cmethod(k),k=1,3 '
write(6,*)' STOPPING '
STOP
end if
CHTIME=0.D0
icollision=0
Taika=TIME ! to common
NofBH=NBH ! - " -
IF(NEWREG)THEN
step=0
iwarning=0
itemax=12
itemax_used=0
ee=soft**2 ! to common
do k=1,3
spin(k)=spini(k) ! 2 common
cmethod(k)=cmet(k) ! 2 common
end do
clight=cl ! -"-
N=NN
mass=0
DO I=1,N
M(I)=MX(I)
mass=mass+m(i)
END DO
MMIJ=0.D0
DO I=1,N-1
DO J=I+1,N
MMIJ=MMIJ+M(I)*M(J)
END DO
END DO
MMIJ=MMIJ/(N*(N-1)/2.d0) !mean mass product
IF(MMIJ.eq.0.0.and.cmethod(2).ne.0.0)MMIJ=1
DO I=1,3*N
X(I)=XX(I)
V(I)=VX(I)
END DO
IF(MMIJ.eq.0)THEN
write(6,*)'You have at most one non-zero mass => t''=1/0 and'
write(6,*)'this does not work'
STOP
END IF
call FIND CHAIN INDICES ! necessary 4 forming the chain
IF(IWR.GT.0)WRITE(6,1232)time,(INAME(KW),KW=1,N)
call INITIALIZE XC and WC
call CONSTANTS OF MOTION(ENERGY,G0,ALAG)
EnerGr=0 ! Energy change due to relativisstic PN-terms
gtime=1/ALAG
do K=1,3
CMX(K)=CMXX(K)
CMV(K)=CMVX(K)
end do
call omegacoef
STIME=0.D0
NEWREG=.FALSE.
WTTL=Wfunction()
mmss=0
do i=1,n-1
do j=i+1,n
mmss=mmss+m(i)*m(j)
end do
end do
call Take Y from XC WC (Y,Nvar)
do i=1,Nvar
SY(i)=0.1 ! arbitrary assumption for the order of magnitude of quantities
end do
if(step.eq.0)then
step=1.e-3
!call Initial Stepsize(X,V,M,N,ee,step) ! New initial step determination
call Estimate Stepsize(tstep,step) !
end if
EPS=TOL
END IF ! NEWREG
KSTEPS=0
stimex=0
777 KSTEPS=KSTEPS+1
call Take Y from XC WC (Y,Nvar)
call Obtain Order of Y(SY,Y)
stime=0
f1=chtime-tstep!deltaT ! for exact time iteration
d1=gtime
call take y from XC WC(Yold,Nvar)
steppi=step ! to common
call DIFSYAB(Nvar,EPS,SY,step,stime,Y)
I_switch=1
call Put Y to XC WC (Y,Nvar) ! move parameters in Y to CHAIN variables
if(step.eq.0)stop
call CHECK SWITCHING CONDITIONS(MUST SWITCH) ! is it necessary to find new CHAIN?
IF(MUST SWITCH)THEN
I_switch=0
call Chain Transformation ! Here the new CHAIN is formed
WTTL=Wfunction() ! this may not be necessary, but probably OK.
call Take Y from XC WC(Y,Nvar) ! CHAIN variables to Y vector
!IF(IWR.GT.0) WRITE(6,1232)time+chtime,(INAME(KW),KW=1,N) !for checking what is going on
1232 FORMAT(1X,g12.4,' I-CHAIN',20I3)
END IF ! MUST SWITCH
f2=chtime-tstep!deltaT ! for exact time iteration
d2=gtime
x1=-stime
x2=0
IF(CHTIME.LT.tstep.and.(KSTEPS.lt.KSMX)
& .and.(icollision.eq.0))goto 777
c------------------------------------------------------------------------------------------
IF (icollision.NE.0) THEN ! handle a collison
nmerger = nmerger + 1
CALL Merge_i1_i2(time) ! MERGE the two particles
newreg=.TRUE. ! chain has changed=> restart needed
NN=N ! copy new chain
DO i=1,NN
MX(i)=M(i)
DO k=1,3
xx(3*i-3+k)=x(3*i-3+k)
vx(3*i-3+k)=v(3*i-3+k)
ENDDO
ENDDO ! done copying new chain
tstep=tnext-time ! set time step to continue
! chain integration
c write(16,*)time, tstep,' GOTO 10 '
IF ((abs(tstep).GT.1.d-6*deltat).AND.(NN.GT.1)) GOTO 10 ! continue if necessary (and if there are at least 2 bodies)
ENDIF
c..........................................................................................
if(KSTEPS.lt.KSMX .and.Ixc.gt.0.and.icollision.eq.0)then
! Integrate TO approximate EXACT OUTPUTTIME
IF(Ixc.eq.1)then ! approx outputtime with Stumpff-Weiss-priciple
if(abs(f1).lt.abs(f2)*I_switch)then ! I_switch prevents use of f1 if just SWITCHed
call put y to xc wc (yold,nvar)
call obtain order of y(sy,y)
1001 call Estimate Stepsize(-f1,step2) !lengt of s-step to get nearly exact output time
!write(6,*)' STEPPI2 ',step2
call DIFSYAB(Nvar,EPS,SY,step2,stime,Yold)
call Put Y to XC WC (Yold,Nvar)
else
call Estimate Stepsize(-f2,step2)
call obtain order of y (sy,y)
!write(6,*)' Steppi2 ',step2
call DIFSYAB(Nvar,EPS,SY,step2,stime,Y)
call Put Y to XC WC (Y,Nvar)
end if
stimex=stimex+stime! 4 estimating max next step
elseif(Ixc.eq.2)then ! Iteration to exact time
call Iterate2ExactTime(Y,Nvar,deltaT,f1,d1,f2,d2,x1,x2)
end if
end if
if(stimex.eq.0)stimex=step
call update x and v
DO I=1,3*N
XX(I)=X(I)
VX(I)=V(I)
END DO
DO I=1,3
spini(I)=spin(I)
CMXX(I)=CMX(I)
CMVX(I)=CMV(I)
END DO
TIME=TIME+CHTIME
if(chtime.lt.0.D0)write(6,*)time,chtime, ' t cht <0!'
!write(16,*)time,tstep,chtime,deltat,' t tstep chtime deltat '
RETURN
END
subroutine Iterate2ExactTime(Y0,Nvar,deltaT,f1,d1,f2,d2,x1,x2)
INCLUDE 'ARCparams.CH'
COMMON/DerOfTime/GTIME
COMMON/collision/icollision,Ione,Itwo,iwarning
REAL*8 Y(1500),SY(1500),Y0(*)
data tiny/1.d-6/
save
! This a quite slow iteration to exact time. Can also fail but that is rear.
iskeleita=0
it=0
hs=abs(x1-x2)
1111 CONTINUE
it=it+1
do i=1,nvar
y(i)=y0(i)
end do
stime=0
dx1=-f1/d1
dx2=-f2/d2
if(abs(dx1).lt.abs(dx2))then
xnew=x1+dx1
else
xnew=x2+dx2
end if
!
test=(x1-xnew)*(xnew-x2)
if(test.lt.(-tiny*hs).or.(it+1).eq.(it+1)/5*5)then
xnew=(x1+x2)/2 ! bisect if out of interval
end if
sfinal=xnew
call Put Y to XC WC (Y,Nvar)
!--------------------------------------------------------------------------
call Obtain Order of Y(SY,y)
do k=1,5
step=sfinal-stime
if(abs(step).gt.1.d-6*abs(hs).or.k.eq.1)then !!!!
call DIFSYAB(Nvar,EPS,SY,step,stime,Y)
iskeleita=iskeleita+1
! it=it+1
else
goto 222
end if
end do
222 continue
call Put Y to XC WC (Y,Nvar)
call UPDATE X AND V
fnew=chtime-deltaT
dfnew=gtime
! keep it bracketed
if(f1*fnew.le.0.D0)then
f2=fnew
d2=dfnew
x2=xnew
else
f1=fnew
d1=dfnew
x1=xnew
end if
if((abs(deltaT-chtime).gt.1.d-6*deltat).and.(it.lt.100))
& goto 1111
! ONE FINAL STEP SHOULD BE HERE (if above not-so-accurate test)
!--------------------------------------------------------------------
do i=1,Nvar
y0(i)=y(i)
end do
call Put Y to XC WC (Y,Nvar)
call UPDATE X AND V
return
end
subroutine LEAPFROG(STEP,Leaps,stime) !Leapfroging steps for the Bulirsch-Stoer extrapolation
implicit real*8 (a-h,M,o-z)
save
CALL PUT V 2 W
hs=step
h2=hs/2
call XCmotion(h2)
stime=stime+h2
do k=1,Leaps-1
call WCmotion(hs)
call XCmotion(hs)
stime=stime+hs
end do
call WCmotion(hs)
call XCmotion(h2)
stime=stime+h2
return
end
function Wfunction() ! evaluate the TTL-function for time tranformation ! eioo subroutine
INCLUDE 'ARCparams.CH'
common/omegacoefficients/OMEC(NMX,NMX)
save
OMEGA=0.0d0
DO I=1,N-1
DO J=I+1,N
if(omec(i,j).ne.0.D0)then
RIJ=SQRT(SQUARE(X(3*I-2),X(3*J-2)))
OMEGA=OMEGA+omec(i,j)/RIJ
end if
END DO
END DO
Wfunction=OMEGA
RETURN
END
subroutine omegacoef ! coefficients in \Omega=sum omec(i,j)/r_{ij}
INCLUDE 'ARCparams.CH'
common/omegacoefficients/OMEC(NMX,NMX)
SAVE
icount=0
do i=1,N-1
do j=i+1,N
! if both masses are=0, then they do not interact (and omec can be set =0)
if(m(i)+m(j).gt.0.D0 .and. cmethod(2).ne.0.D0)then
OMEC(I,J)=mmij
OMEC(J,I)=mmij
icount=icount+1
else
OMEC(I,J)=0
OMEC(J,I)=0
end if
end do
end do
if(icount.eq.0)cmethod(2)=0 ! all terms zero anyway (for some different settings of omec)
return
end
SUBROUTINE XCMOTION(hs) ! movement of coordinates
INCLUDE 'ARCparams.CH'
COMMON/IncrementCommon/WTTLinc,XCinc(NMX3),WCinc(NMX3),
& CMXinc(3),CMVinc(3),ENERGYinc,Energrinc,CHTIMEinc,spin inc(3)
COMMON/DerOfTime/G
COMMON/DIAGNOSTICS/GAMMA,H,IWR
save
Te=-ENERGY
if(cmethod(1).ne.0.0d0)then
call EVALUATE V(V,WC)
do I=1,N
I0=3*I-3
Te=Te+M(I)*(V(I0+1)**2+V(I0+2)**2+V(I0+3)**2)/2 ! Kinetic energyi-Energy
end do
end if ! cmethod(1).ne.0.0d0
G=1/(Te*cmethod(1)+WTTL*cmethod(2)+cmethod(3)) ! time transformation = t'
if(G.lt.0.D0.and.iwr.gt.0)then
write(6,*)1/G,' tdot <0 ! '
return ! this was seriously wrong, but may work because this step (likely) gets rejected
end if
dT= hs*G
DO I=1,N-1
L=3*(I-1)
DO K=1,3
XCinc(L+K)=XCinc(L+k)+WC(L+K)*dT ! incerements of coords
XC(L+K)=XC(L+K)+WC(L+K)*dT ! the new coords are needed anyway
END DO
END DO
CHTIMEinc=CHTIMEinc+dT ! time increment (time is a coordinate)
CHTIME=CHTIME+dT ! new time value
do k=1,3
CMXinc(k)=CMXinc(k)+dt*cmv(k) ! center of mass increment
cmx(k)=cmx(k)+dt*cmv(k) ! new cm-coordinates
end do
RETURN
END
subroutine PUT V 2 W ! these parameters are needed if there are v-ependent forces
include 'ARCparams.CH'
common/vwcommon/Ww(nmx3),WTTLw,cmvw(3),spinw(3)
save
do i=1,3*(N-1)
Ww(i)=WC(I)
end do
WTTLw=WTTL
do k=1,3
spinw(k)=spin(k)
cmvw(k)=cmv(k)
end do
return
end
subroutine Velocity Dependent Perturbations ! This name tells what this is!
& (dT,Va,spina,acc,dcmv,df,dfGR,dspin)
INCLUDE 'ARCparams.CH'
common/vindex/ivelocity
real*8 df(*),Va(*),dcmv(3),dfGR(*),dfR(nmx3),acc(nmx3)
real*8 dspin(3),spina(3)
save
! v-dependent perturbations are evaluated here
do i=1,3*n
dfr(i)=0
dfgr(i)=0
end do
do k=1,3
dspin(k)=0
end do
if(Clight.ne.0)then ! include only if Clight set >0
call Relativistic ACCELERATIONS(dfr,dfGR,Va,spina,dspin)
end if
if(ivelocity.gt.0)then ! USED ONLY IF ivelcity.gt.0
call Non relativistic v_dependent perturbations(dfr)! add v-dependent to dfr(), (e.g. friction)
end if
do i=1,3*n
df(i)=acc(i)+dfr(i)
end do
call reduce 2 cm(df,m,n,dcmv)
return
end
SUBROUTINE CHECK SWITCHING CONDITIONS(MUSTSWITCH)! Is it necessary to construct a new CHAIN?
INCLUDE 'ARCparams.CH'
LOGICAL MUSTSWITCH
DATA Ncall,NSWITCH/0,200000000/
save
MUST SWITCH=.FALSE.
Ncall=Ncall+1
! Switch anyway after every NSWITCHth step.
IF(Ncall.GE.NSWITCH)THEN
Ncall=0
MUST SWITCH=.TRUE.
RETURN
END IF
! Inspect the structure of the chain.
! NOTE: Inverse values 1/r are used instead of r itself.
ADISTI=0.5*(N-1)/RSUM
LRI=N-1
DO I=1,N-2
DO J=I+2,N
LRI=LRI+1
! Do not inspect if 1/r is small.
IF(RINV(LRI).GT.ADISTI)THEN
IF(J-I.GT.2)THEN
! Check for a dangerous long loop.
! RINVMX=MAX(RINV(I-1),RINV(I),RINV(J-1),RINV(J))
IF(I.GT.1)THEN
RINVMX=MAX(RINV(I-1),RINV(I))
ELSE
RINVMX=RINV(1)
END IF
RINVMX=MAX(RINVMX,RINV(J-1))
IF(J.LT.N)RINVMX=MAX(RINVMX,RINV(J))
IF(RINV(LRI).GT.RINVMX)THEN ! 0.7*RINVMX may be more careful
MUST SWITCH=.TRUE.
Ncall=0
RETURN
END IF
ELSE
! Is this a triangle with smallest size not regularised?
IF( RINV(LRI).GT.MAX(RINV(I),RINV(I+1)))THEN
MUST SWITCH=.TRUE.
Ncall=0
RETURN
END IF
END IF
END IF
END DO
END DO
RETURN
END
SUBROUTINE FIND CHAIN INDICES ! Here the new indecies along the new CHAIN are obtained
INCLUDE 'ARCparams.CH'
REAL*8 RIJ2(NMXM)
INTEGER IC(NMX2),IJ(NMXM,2),IND(NMXM)
LOGICAL USED(NMXM),SUC,LOOP
save
L=0
DO I=1,N-1
DO J=I+1,N
L=L+1
RIJ2(L)=SQUARE(X(3*I-2),X(3*J-2))
IJ(L,1)=I
IJ(L,2)=J
USED(L)=.FALSE.
END DO
END DO
call ARRANGE(L,RIJ2,IND)
LMIN=1+NMX
LMAX=2+NMX
IC(LMIN)=IJ(IND(1),1)
IC(LMAX)=IJ(IND(1),2)
USED(IND(1))=.TRUE.
1 DO I=2,L
LI=IND(I)
IF( .NOT.USED(LI))THEN
call CHECK CONNECTION(IC,LMIN,LMAX,IJ,LI,SUC,LOOP)
IF(SUC)THEN
USED(LI)=.TRUE.
GOTO 2
ELSE
USED(LI)=LOOP
END IF
END IF
END DO
2 IF(LMAX-LMIN+1.LT.N)GO TO 1
L=0
DO I=LMIN,LMAX
L=L+1
INAME(L)=IC(I)
END DO
RETURN
END
SUBROUTINE CHECK CONNECTION(IC,LMIN,LMAX,IJ,LI,SUC,LOOP) ! Loops in CHAIN are not allowed
INCLUDE 'ARCparams.CH'
INTEGER IC(*),ICC(2),IJ(NMXM,2)
LOGICAL SUC,LOOP
save
SUC=.FALSE.
LOOP=.FALSE.
ICC(1)=IC(LMIN)
ICC(2)=IC(LMAX)
DO I=1,2
DO J=1,2
IF(ICC(I).EQ.IJ(LI,J))THEN
JC=3-J
LOOP=.TRUE.
DO L=LMIN,LMAX
IF(IC(L).EQ.IJ(LI,JC))RETURN
END DO
SUC=.TRUE.
LOOP=.FALSE.
IF(I.EQ.1)THEN
LMIN=LMIN-1
IC(LMIN)=IJ(LI,JC)
RETURN
ELSE
LMAX=LMAX+1
IC(LMAX)=IJ(LI,JC)
RETURN
END IF
END IF
END DO
END DO
RETURN
END
SUBROUTINE ARRANGE(N,Array,Indx) ! this is a sorting routine (needed in CHAIN construction)
implicit real*8 (a-h,o-z)
dimension Array(*),Indx(*)
save
do 11 j=1,N
Indx(j)=J
11 continue
if(N.lt.2)RETURN
l=N/2+1
ir=N
10 CONTINUE
IF(l.gt.1)THEN
l=l-1
Indxt=Indx(l)
q=Array(Indxt)
ELSE
Indxt=Indx(ir)
q=Array(Indxt)
Indx(ir)=Indx(1)
ir=ir-1
IF(ir.eq.1)THEN
Indx(1)=Indxt
RETURN
END IF
END IF
i=l
j=l+l
20 IF(j.le.ir)THEN
IF(j.lt.ir)THEN
IF(Array(Indx(j)).lt.Array(Indx(j+1)))j=j+1
END IF
IF(q.lt.Array(Indx(j)))THEN
Indx(i)=Indx(j)
i=j
j=j+j
ELSE
j=ir+1
END IF
GOTO 20
END IF
Indx(i)=Indxt
GO TO 10
END