-
Notifications
You must be signed in to change notification settings - Fork 26
/
evolv2.f
2340 lines (2340 loc) · 74 KB
/
evolv2.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
***
SUBROUTINE evolv2(kstar,mass0,mass,rad,lumin,massc,radc,
& menv,renv,ospin,epoch,tms,
& tphys,tphysf,dtp,z,zpars,tb,ecc,vkick)
implicit none
***
*
* B I N A R Y
* ***********
*
* Roche lobe overflow.
* --------------------
*
* Developed by Jarrod Hurley, IOA, Cambridge.
* .........................................................
*
* Advice by Christopher Tout, Onno Pols & Sverre Aarseth.
* ++++++++++++++++++++++++++++++++++++++++++++++++++
*
* Adapted from Aarseth's code 21st September 1996.
* Fully revised on 27th November 1996 to remove vestiges of N-body code and
* incorporate corrections.
* Fully revised on 1st April 1998 to include new stellar evolution formulae
* and associated binary evolution changes.
* Fully revised on 4th July 1998 to include eccentricity, tidal
* circularization, wind accretion, velocity kicks for supernovae and all
* associated orbital momentum changes.
* Revised on 31st October 2000 to upgrade restrictions imposed on the
* timestep owing to magnetic braking and orbital angular momentum changes.
*
***
*
* See Tout et al., 1997, MNRAS, 291, 732 for a description of many of the
* processes in this code as well as the relevant references mentioned
* within the code.
*
* Reference for the stellar evolution formulae is Hurley, Pols & Tout,
* 2000, MNRAS, 315, 543 (SSE paper).
* Reference for the binary evolution algorithm is Hurley, Tout & Pols,
* 2002, MNRAS, 329, 897 (BSE paper).
*
***
*
* March 2001 *
* Changes since version 3, i.e. since production of Paper3:
*
* 1) The Eddington limit flag (on/off) has been replaced by an
* Eddington limit multiplicative factor (eddfac). So if you
* want to neglect the Eddington limit you would set eddfac
* to a large value.
*
* 2) To determine whether material transferred during RLOF forms
* an accretion disk around the secondary or hits the secondary
* in a direct stream we calculate a minimum radial distance, rmin,
* of the mass stream from the secondary. This is taken from eq.(1)
* of Ulrich & Burger (1976, ApJ, 206, 509) which they fitted to
* the calculations of Lubow & Shu (1974, ApJ, 198, 383).
* If rmin is less than the radius of the secondary then an
* accretion disk is not formed.
* Note that the formula for rmin given by Ulrich & Burger is valid
* for all q whereas that given by Nelemans et al. (2001, A&A,
* submitted) in their eq.(6) is only valid for q < 1 where
* they define q = Mdonor/Maccretor, i.e. DD systems.
*
* 3) The changes to orbital and spin angular momentum owing to
* RLOF mass transfer have been improved, and an new input option
* now exists.
* When mass is lost from the system during RLOF there are now
* three choices as to how the orbital angular momentum is
* affected: a) the lost material carries with it a fraction
* gamma of the orbital angular momentum, i.e.
* dJorb = gamma*dm*a^2*omega_orb; b) the material carries with it
* the specific angular momentum of the primary, i.e.
* dJorb = dm*a_1^2*omega_orb; or c) assume the material is lost
* from the system as if a wind from the secondary, i.e.
* dJorb = dm*a_2^2*omega_orb.
* The parameter gamma is an input option.
* Choice c) is used if the mass transfer is super-Eddington
* or the system is experiencing novae eruptions.
* In all other cases choice a) is used if gamma > 0.0, b) if
* gamma = -1.0 and c) is used if gamma = -2.0.
* The primary spin angular momentum is reduced by an amount
* dm1*r_1^2*omega_1 when an amount of mass dm1 is transferred
* from the primary.
* If the secondary accretes through a disk then its spin
* angular momentum is altered by assuming that the material
* falls onto the star from the inner edge of a Keplerian
* disk and that the system is in a steady state, i.e.
* an amount dm2*SQRT(G*m_2*r_2).
* If there is no accretion disk then we calculate the angular
* momentum of the transferred material by using the radius at
* at which the disk would have formed (rdisk = 1.7*rmin, see
* Ulrich & Burger 1976) if allowed, i.e. the angular momentum
* of the inner Lagrangian point, and add this directly to
* the secondary, i.e. an amount dm2*SQRT(G*m_2*rdisk).
* Total angular momentum is conserved in this model.
*
* 4) Now using q_crit = 3.0 for MS-MS Roche systems (previously we
* had nothing). This corresponds roughly to R proportional to M^5
* which should be true for the majority of the MS (varies from
* (M^17 -> M^2). If q > q_crit then contact occurs.
* For CHeB primaries we also take q_crit = 3.0 and allow
* common-envelope to occur if this is exceeded.
*
* 5) The value of lambda used in calculations of the envelope binding
* energy for giants in common-envelope is now variable (see function
* in zfuncs). The lambda function has been fitted by Onno to detailed
* models ... he will write about this soon!
*
* 6) Note that eq.42 in the paper is missing a SQRT around the
* MR^2/a^5 part. This needs to be corrected in any code update
* paper with a thanks to Jeremy Sepinsky (student at NorthWestern).
* It is ok in the code.
*
* March 2003 *
* New input options added:
*
* ifflag - for the mass of a WD you can choose to use the mass that
* results from the evolution algorithm (basically a competition
* between core-mass growth and envelope mass-loss) or use the IFMR
* proposed by Han, Podsiadlowski & Eggleton, 1995, MNRAS, 272, 800
* [>0 activates HPE IFMR].
*
* wdflag - for the cooling of WDs you can choose to use either the standard
* Mestel cooling law (see SSE paper) or a modified-Mestel law that
* is better matched to detailed models (provided by Brad Hansen
* ... see Hurley & Shara, 2003, ApJ, May 20, in press)
* [>0 activates modified-Mestel].
*
* bhflag - choose whether or not black holes should get velocity kicks
* at formation
* [0= no kick; >0 kick].
*
* nsflag - for the mass of neutron stars and black holes you can use either
* the SSE prescription or the prescription presented by
* Belczynski et al. 2002, ApJ, 572, 407 who found that SSE was
* underestimating the masses of these stars. In either case you also
* need to set the maximum NS mass (mxns) for the prescription
* [0= SSE, mxns=1.8; >0 Belczynski, mxns=3.0].
*
* Sept 2004 *
* Input options added/changed:
*
* ceflag - set to 3 this uses de Kool (or Podsiadlowski) CE prescription,
* other options, such as Yungelson, could be added as well.
*
* hewind - factor to control the amount of He star mass-loss, i.e.
* 1.0e-13*hewind*L^(2/3) gives He star mass-loss.
*
*
* ++++++++++++++++++++++++++++++++++++++++++++++++++
***
*
INTEGER loop,iter,intpol,k,ip,jp,j1,j2
PARAMETER(loop=20000)
INTEGER kstar(2),kw,kst,kw1,kw2,kmin,kmax
INTEGER ktype(0:14,0:14)
COMMON /TYPES/ ktype
INTEGER ceflag,tflag,ifflag,nsflag,wdflag
COMMON /FLAGS/ ceflag,tflag,ifflag,nsflag,wdflag
*
REAL*8 km,km0,tphys,tphys0,dtm0,tphys00
REAL*8 tphysf,dtp,tsave
REAL*8 aj(2),aj0(2),epoch(2),tms(2),tbgb(2),tkh(2),dtmi(2)
REAL*8 mass0(2),mass(2),massc(2),menv(2),mass00(2),mcxx(2)
REAL*8 rad(2),rol(2),rol0(2),rdot(2),radc(2),renv(2),radx(2)
REAL*8 lumin(2),k2str(2),q(2),dms(2),dmr(2),dmt(2),vkick(2)
REAL*8 dml,vorb2,vwind2,omv2,ivsqm,lacc,vs(3)
REAL*8 sep,dr,tb,dme,tdyn,taum,dm1,dm2,dmchk,qc,dt,pd,rlperi
REAL*8 m1ce,m2ce,mch,tmsnew,dm22,mew
PARAMETER(mch=1.44d0)
REAL*8 yeardy,aursun
PARAMETER(yeardy=365.24d0,aursun=214.95d0)
REAL*8 acc1,tiny
PARAMETER(acc1=3.920659d+08,tiny=1.0d-14)
REAL*8 ecc,ecc1,tc,tcirc,ttid,ecc2,omecc2,sqome2,sqome3,sqome5
REAL*8 f1,f2,f3,f4,f5,f,raa2,raa6,eqspin,rg2,tcqr
REAL*8 k3,mr23yr,twopi
PARAMETER(k3=0.21d0,mr23yr=0.4311d0)
REAL*8 jspin(2),ospin(2),jorb,oorb,jspbru,ospbru
REAL*8 delet,delet1,dspint(2),djspint(2),djtx(2)
REAL*8 dtj,djorb,djgr,djmb,djt,djtt,rmin,rdisk
REAL*8 neta,bwind,hewind,mxns
COMMON /VALUE1/ neta,bwind,hewind,mxns
REAL*8 beta,xi,acc2,epsnov,eddfac,gamma
COMMON /VALUE5/ beta,xi,acc2,epsnov,eddfac,gamma
*
REAL*8 z,tm,tn,m0,mt,rm,lum,mc,rc,me,re,k2,age,dtm,dtr
REAL*8 tscls(20),lums(10),GB(10),zpars(20)
REAL*8 zero,ngtv,ngtv2,mt2,rrl1,rrl2,mcx,teff1,teff2
REAL*8 mass1i,mass2i,tbi,ecci
LOGICAL coel,com,prec,inttry,change,snova,sgl,bsymb,esymb,bss
LOGICAL supedd,novae,disk
LOGICAL isave,iplot
REAL*8 rl,mlwind,vrotf,corerd
EXTERNAL rl,mlwind,vrotf,corerd
REAL bcm(50000,34),bpp(80,10)
COMMON /BINARY/ bcm,bpp
*
* Save the initial state.
*
mass1i = mass0(1)
mass2i = mass0(2)
tbi = tb
ecci = ecc
*
zero = 0.d0
ngtv = -1.d0
ngtv2 = -2.d0
twopi = 2.d0*ACOS(-1.d0)
*
* Initialize the parameters.
*
kmin = 1
kmax = 2
sgl = .false.
mt2 = MIN(mass(1),mass(2))
kst = 0
*
if(mt2.lt.tiny.or.tb.le.0.d0)then
sgl = .true.
if(mt2.lt.tiny)then
mt2 = 0.d0
if(mass(1).lt.tiny)then
if(tphys.lt.tiny)then
mass0(1) = 0.01d0
mass(1) = mass0(1)
kst = 1
else
kmin = 2
lumin(1) = 1.0d-10
rad(1) = 1.0d-10
massc(1) = 0.d0
dmt(1) = 0.d0
dmr(1) = 0.d0
endif
ospin(1) = 1.0d-10
jspin(1) = 1.0d-10
else
if(tphys.lt.tiny)then
mass0(2) = 0.01d0
mass(2) = mass0(2)
kst = 2
else
kmax = 1
lumin(2) = 1.0d-10
rad(2) = 1.0d-10
massc(2) = 0.d0
dmt(2) = 0.d0
dmr(2) = 0.d0
endif
ospin(2) = 1.0d-10
jspin(2) = 1.0d-10
endif
endif
ecc = -1.d0
tb = 0.d0
sep = 1.0d+10
oorb = 0.d0
jorb = 0.d0
if(ospin(1).lt.0.0) ospin(1) = 1.0d-10
if(ospin(2).lt.0.0) ospin(2) = 1.0d-10
q(1) = 1.0d+10
q(2) = 1.0d+10
rol(1) = 1.0d+10
rol(2) = 1.0d+10
else
tb = tb/yeardy
sep = aursun*(tb*tb*(mass(1) + mass(2)))**(1.d0/3.d0)
oorb = twopi/tb
jorb = mass(1)*mass(2)/(mass(1)+mass(2))
& *SQRT(1.d0-ecc*ecc)*sep*sep*oorb
if(ospin(1).lt.0.d0) ospin(1) = oorb
if(ospin(2).lt.0.d0) ospin(2) = oorb
endif
*
do 500 , k = kmin,kmax
age = tphys - epoch(k)
CALL star(kstar(k),mass0(k),mass(k),tm,tn,tscls,lums,GB,zpars)
CALL hrdiag(mass0(k),age,mass(k),tm,tn,tscls,lums,GB,zpars,
& rm,lum,kstar(k),mc,rc,me,re,k2)
aj(k) = age
epoch(k) = tphys - age
rad(k) = rm
lumin(k) = lum
massc(k) = mc
radc(k) = rc
menv(k) = me
renv(k) = re
k2str(k) = k2
tms(k) = tm
tbgb(k) = tscls(1)
*
if(tphys.lt.tiny.and.ospin(k).le.0.001d0)then
ospin(k) = 45.35d0*vrotf(mass(k))/rm
endif
jspin(k) = ospin(k)*(k2*rm*rm*(mass(k)-mc)+k3*rc*rc*mc)
if(.not.sgl)then
q(k) = mass(k)/mass(3-k)
rol(k) = rl(q(k))*sep
endif
rol0(k) = rol(k)
dmr(k) = 0.d0
dmt(k) = 0.d0
djspint(k) = 0.d0
dtmi(k) = 1.0d+06
*
500 continue
*
if(mt2.lt.tiny)then
sep = 0.d0
if(kst.gt.0)then
mass0(kst) = 0.d0
mass(kst) = 0.d0
kmin = 3 - kst
kmax = kmin
endif
endif
*
* On the first entry the previous timestep is zero to prevent mass loss.
*
dtm = 0.d0
delet = 0.d0
djorb = 0.d0
bss = .false.
*
* Setup variables which control the output (if it is required).
*
ip = 0
jp = 0
tsave = tphys
isave = .true.
iplot = .false.
if(dtp.le.0.d0)then
iplot = .true.
isave = .false.
tsave = tphysf
elseif(dtp.gt.tphysf)then
isave = .false.
tsave = tphysf
endif
if(tphys.ge.tphysf) goto 140
*
4 iter = 0
intpol = 0
inttry = .false.
change = .false.
prec = .false.
snova = .false.
coel = .false.
com = .false.
bsymb = .false.
esymb = .false.
tphys0 = tphys
ecc1 = ecc
j1 = 1
j2 = 2
if(kstar(1).ge.10.and.kstar(1).le.14) dtmi(1) = 0.01d0
if(kstar(2).ge.10.and.kstar(2).le.14) dtmi(2) = 0.01d0
dm1 = 0.d0
dm2 = 0.d0
*
5 kw1 = kstar(1)
kw2 = kstar(2)
*
dt = 1.0d+06*dtm
eqspin = 0.d0
djtt = 0.d0
*
if(intpol.eq.0.and.ABS(dtm).gt.tiny.and..not.sgl)then
vorb2 = acc1*(mass(1)+mass(2))/sep
ivsqm = 1.d0/SQRT(1.d0-ecc*ecc)
do 501 , k = 1,2
*
* Calculate wind mass loss from the previous timestep.
*
if(neta.gt.tiny)then
rlperi = rol(k)*(1.d0-ecc)
dmr(k) = mlwind(kstar(k),lumin(k),rad(k),mass(k),
& massc(k),rlperi,z)
*
* Calculate how much of wind mass loss from companion will be
* accreted (Boffin & Jorissen, A&A 1988, 205, 155).
*
vwind2 = 2.d0*beta*acc1*mass(k)/rad(k)
omv2 = (1.d0 + vorb2/vwind2)**(3.d0/2.d0)
dmt(3-k) = ivsqm*acc2*dmr(k)*((acc1*mass(3-k)/vwind2)**2)
& /(2.d0*sep*sep*omv2)
dmt(3-k) = MIN(dmt(3-k),0.8d0*dmr(k))
else
dmr(k) = 0.d0
dmt(3-k) = 0.d0
endif
501 continue
*
* Diagnostic for Symbiotic-type stars.
*
if(neta.gt.tiny.and..not.esymb)then
lacc = 3.14d+07*mass(j2)*dmt(j2)/rad(j2)
lacc = lacc/lumin(j1)
if((lacc.gt.0.01d0.and..not.bsymb).or.
& (lacc.lt.0.01d0.and.bsymb))then
jp = MIN(80,jp + 1)
bpp(jp,1) = tphys
bpp(jp,2) = mass(1)
bpp(jp,3) = mass(2)
bpp(jp,4) = float(kstar(1))
bpp(jp,5) = float(kstar(2))
bpp(jp,6) = sep
bpp(jp,7) = ecc
bpp(jp,8) = rad(1)/rol(1)
bpp(jp,9) = rad(2)/rol(2)
if(bsymb)then
bpp(jp,10) = 13.0
esymb = .true.
else
bpp(jp,10) = 12.0
bsymb = .true.
endif
endif
endif
*
* Calculate orbital angular momentum change due to wind mass loss.
*
ecc2 = ecc*ecc
omecc2 = 1.d0 - ecc2
sqome2 = SQRT(omecc2)
*
djorb = ((dmr(1)+q(1)*dmt(1))*mass(2)*mass(2) +
& (dmr(2)+q(2)*dmt(2))*mass(1)*mass(1))*
& sep*sep*sqome2*oorb/(mass(1)+mass(2))**2
delet = ecc*(dmt(1)*(0.5d0/mass(1) + 1.d0/(mass(1)+mass(2))) +
& dmt(2)*(0.5d0/mass(2) + 1.d0/(mass(1)+mass(2))))
*
* For very close systems include angular momentum loss owing to
* gravitational radiation.
*
if(sep.le.10.d0)then
djgr = 8.315d-10*mass(1)*mass(2)*(mass(1)+mass(2))/
& (sep*sep*sep*sep)
f1 = (19.d0/6.d0) + (121.d0/96.d0)*ecc2
sqome5 = sqome2**5
delet1 = djgr*ecc*f1/sqome5
djgr = djgr*jorb*(1.d0+0.875d0*ecc2)/sqome5
djorb = djorb + djgr
delet = delet + delet1
endif
*
do 502 , k = 1,2
*
* Calculate change in the intrinsic spin of the star.
*
djtx(k) = (2.d0/3.d0)*xi*dmt(k)*rad(3-k)*rad(3-k)*ospin(3-k)
djspint(k) = (2.d0/3.d0)*(dmr(k)*rad(k)*rad(k)*ospin(k)) -
& djtx(k)
*
* Include magnetic braking for stars that have appreciable convective
* envelopes. This includes MS stars with M < 1.25, HG stars near the GB
* and giants. MB is not allowed for fully convective MS stars.
*
if(mass(k).gt.0.35d0.and.kstar(k).lt.10)then
djmb = 5.83d-16*menv(k)*(rad(k)*ospin(k))**3/mass(k)
djspint(k) = djspint(k) + djmb
*
* Limit to a 3% angular momentum change for the star owing to MB.
* This is found to work best with the maximum iteration of 20000,
* i.e. does not create an excessive number of iterations, while not
* affecting the evolution outcome when compared with a 2% restriction.
*
if(djmb.gt.tiny)then
dtj = 0.03d0*jspin(k)/ABS(djmb)
dt = MIN(dt,dtj)
endif
endif
*
* Calculate circularization, orbital shrinkage and spin up.
*
dspint(k) = 0.d0
if(((kstar(k).le.9.and.rad(k).ge.0.01d0*rol(k)).or.
& (kstar(k).ge.10.and.k.eq.j1)).and.tflag.gt.0)then
*
raa2 = (rad(k)/sep)**2
raa6 = raa2**3
*
* Hut's polynomials.
*
f5 = 1.d0+ecc2*(3.d0+ecc2*0.375d0)
f4 = 1.d0+ecc2*(1.5d0+ecc2*0.125d0)
f3 = 1.d0+ecc2*(3.75d0+ecc2*(1.875d0+ecc2*7.8125d-02))
f2 = 1.d0+ecc2*(7.5d0+ecc2*(5.625d0+ecc2*0.3125d0))
f1 = 1.d0+ecc2*(15.5d0+ecc2*(31.875d0+ecc2*(11.5625d0
& +ecc2*0.390625d0)))
*
if((kstar(k).eq.1.and.mass(k).ge.1.25d0).or.
& kstar(k).eq.4.or.kstar(k).eq.7)then
*
* Radiative damping (Zahn, 1977, A&A, 57, 383 and 1975, A&A, 41, 329).
*
tc = 1.592d-09*(mass(k)**2.84d0)
f = 1.9782d+04*SQRT((mass(k)*rad(k)*rad(k))/sep**5)*
& tc*(1.d0+q(3-k))**(5.d0/6.d0)
tcqr = f*q(3-k)*raa6
rg2 = k2str(k)
elseif(kstar(k).le.9)then
*
* Convective damping (Hut, 1981, A&A, 99, 126).
*
tc = mr23yr*(menv(k)*renv(k)*(rad(k)-0.5d0*renv(k))/
& (3.d0*lumin(k)))**(1.d0/3.d0)
ttid = twopi/(1.0d-10 + ABS(oorb - ospin(k)))
f = MIN(1.d0,(ttid/(2.d0*tc)**2))
tcqr = 2.d0*f*q(3-k)*raa6*menv(k)/(21.d0*tc*mass(k))
rg2 = (k2str(k)*(mass(k)-massc(k)))/mass(k)
else
*
* Degenerate damping (Campbell, 1984, MNRAS, 207, 433)
*
f = 7.33d-09*(lumin(k)/mass(k))**(5.d0/7.d0)
tcqr = f*q(3-k)*q(3-k)*raa2*raa2/(1.d0+q(3-k))
rg2 = k3
endif
*
* Circularization.
*
sqome3 = sqome2**3
delet1 = 27.d0*tcqr*(1.d0+q(3-k))*raa2*(ecc/sqome2**13)*
& (f3 - (11.d0/18.d0)*sqome3*f4*ospin(k)/oorb)
tcirc = ecc/(ABS(delet1) + 1.0d-20)
delet = delet + delet1
*
* Spin up of star.
*
dspint(k) = (3.d0*q(3-k)*tcqr/(rg2*omecc2**6))*
& (f2*oorb - sqome3*f5*ospin(k))
*
* Calculate the equilibrium spin at which no angular momentum
* can be transferred.
*
eqspin = oorb*f2/(sqome3*f5)
*
* Calculate angular momentum change for the star owing to tides.
*
djt = (k2str(k)*(mass(k)-massc(k))*rad(k)*rad(k) +
& k3*massc(k)*radc(k)*radc(k))*dspint(k)
if(kstar(k).le.6.or.ABS(djt)/jspin(k).gt.0.1d0)then
djtt = djtt + djt
endif
endif
502 continue
*
* Limit to 2% orbital angular momentum change.
*
djtt = djtt + djorb
if(ABS(djtt).gt.tiny)then
dtj = 0.02d0*jorb/ABS(djtt)
dt = MIN(dt,dtj)
endif
dtm = dt/1.0d+06
*
elseif(ABS(dtm).gt.tiny.and.sgl)then
do 503 , k = kmin,kmax
if(neta.gt.tiny)then
rlperi = 0.d0
dmr(k) = mlwind(kstar(k),lumin(k),rad(k),mass(k),
& massc(k),rlperi,z)
else
dmr(k) = 0.d0
endif
dmt(k) = 0.d0
djspint(k) = (2.d0/3.d0)*dmr(k)*rad(k)*rad(k)*ospin(k)
if(mass(k).gt.0.35d0.and.kstar(k).lt.10)then
djmb = 5.83d-16*menv(k)*(rad(k)*ospin(k))**3/mass(k)
djspint(k) = djspint(k) + djmb
if(djmb.gt.tiny)then
dtj = 0.03d0*jspin(k)/ABS(djmb)
dt = MIN(dt,dtj)
endif
endif
503 continue
dtm = dt/1.0d+06
endif
*
do 504 , k = kmin,kmax
*
dms(k) = (dmr(k) - dmt(k))*dt
if(kstar(k).lt.10)then
dml = mass(k) - massc(k)
if(dml.lt.dms(k))then
dml = MAX(dml,2.d0*tiny)
dtm = (dml/dms(k))*dtm
if(k.eq.2) dms(1) = dms(1)*dml/dms(2)
dms(k) = dml
dt = 1.0d+06*dtm
endif
*
* Limit to 1% mass loss.
*
if(dms(k).gt.0.01d0*mass(k))then
dtm = 0.01d0*mass(k)*dtm/dms(k)
if(k.eq.2) dms(1) = dms(1)*0.01d0*mass(2)/dms(2)
dms(k) = 0.01d0*mass(k)
dt = 1.0d+06*dtm
endif
endif
*
504 continue
*
* Update mass and intrinsic spin (checking that the star is not spun
* past the equilibrium) and reset epoch for a MS (and possibly a HG) star.
*
do 505 , k = kmin,kmax
*
if(eqspin.gt.0.d0.and.ABS(dspint(k)).gt.tiny)then
if(intpol.eq.0)then
if(dspint(k).ge.0.d0)then
dspint(k) = MIN(dspint(k),(eqspin-ospin(k))/dt)
else
dspint(k) = MAX(dspint(k),(eqspin-ospin(k))/dt)
endif
djt = (k2str(k)*(mass(k)-massc(k))*rad(k)*rad(k) +
& k3*massc(k)*radc(k)*radc(k))*dspint(k)
djorb = djorb + djt
djspint(k) = djspint(k) - djt
endif
endif
*
jspin(k) = MAX(1.0d-10,jspin(k) - djspint(k)*dt)
*
* Ensure that the star does not spin up beyond break-up.
*
ospbru = twopi*SQRT(mass(k)*aursun**3/rad(k)**3)
jspbru = (k2str(k)*(mass(k)-massc(k))*rad(k)*rad(k) +
& k3*massc(k)*radc(k)*radc(k))*ospbru
if(jspin(k).gt.jspbru.and.ABS(dtm).gt.tiny)then
mew = 1.d0
if(djtx(k).gt.0.d0)then
mew = MIN(mew,(jspin(k) - jspbru)/djtx(k))
endif
jspin(k) = jspbru
* If excess material should not be accreted, activate next line.
* dms(k) = dms(k) + (1.d0 - mew)*dmt(k)*dt
endif
*
if(ABS(dms(k)).gt.tiny)then
mass(k) = mass(k) - dms(k)
if(kstar(k).le.2.or.kstar(k).eq.7)then
m0 = mass0(k)
mass0(k) = mass(k)
CALL star(kstar(k),mass0(k),mass(k),tm,tn,tscls,
& lums,GB,zpars)
if(kstar(k).eq.2)then
if(GB(9).lt.massc(k).or.m0.gt.zpars(3))then
mass0(k) = m0
else
epoch(k) = tm + (tscls(1) - tm)*(aj(k)-tms(k))/
& (tbgb(k) - tms(k))
epoch(k) = tphys - epoch(k)
endif
else
epoch(k) = tphys - aj(k)*tm/tms(k)
endif
endif
endif
*
505 continue
*
if(.not.sgl)then
*
ecc1 = ecc1 - delet*dt
ecc = MAX(ecc1,0.d0)
if(ecc.lt.1.0d-10) ecc = 0.d0
*
if(ecc.ge.1.d0) goto 135
*
jorb = jorb - djorb*dt
sep = (mass(1) + mass(2))*jorb*jorb/
& ((mass(1)*mass(2)*twopi)**2*aursun**3*(1.d0-ecc*ecc))
tb = (sep/aursun)*SQRT(sep/(aursun*(mass(1)+mass(2))))
oorb = twopi/tb
endif
*
* Advance the time.
*
if(intpol.eq.0)then
tphys0 = tphys
dtm0 = dtm
endif
tphys = tphys + dtm
*
do 6 , k = kmin,kmax
*
* Acquire stellar parameters (M, R, L, Mc & K*) at apparent evolution age.
*
age = tphys - epoch(k)
aj0(k) = age
kw = kstar(k)
m0 = mass0(k)
mt = mass(k)
mc = massc(k)
if(intpol.eq.0) mcxx(k) = mc
if(intpol.gt.0) mc = mcxx(k)
mass00(k) = m0
*
* Masses over 100Msun should probably not be trusted in the
* evolution formulae.
*
if(mt.gt.100.d0)then
WRITE(99,*)' MASS EXCEEDED ',mass1i,mass2i,tbi,ecci,mt
goto 140
endif
*
CALL star(kw,m0,mt,tm,tn,tscls,lums,GB,zpars)
CALL hrdiag(m0,age,mt,tm,tn,tscls,lums,GB,zpars,
& rm,lum,kw,mc,rc,me,re,k2)
*
if(kw.ne.15)then
ospin(k) = jspin(k)/(k2*(mt-mc)*rm*rm+k3*mc*rc*rc)
endif
*
* At this point there may have been a supernova.
*
if(kw.ne.kstar(k).and.kstar(k).le.12.and.
& (kw.eq.13.or.kw.eq.14))then
if(sgl)then
CALL kick(kw,mass(k),mt,0.d0,0.d0,-1.d0,0.d0,vs)
vkick(k) = dsqrt(vs(1)*vs(1)+vs(2)*vs(2)+vs(3)*vs(3))
else
CALL kick(kw,mass(k),mt,mass(3-k),ecc,sep,jorb,vs)
vkick(k) = dsqrt(vs(1)*vs(1)+vs(2)*vs(2)+vs(3)*vs(3))
if(ecc.gt.1.d0)then
kstar(k) = kw
mass(k) = mt
epoch(k) = tphys - age
goto 135
endif
tb = (sep/aursun)*SQRT(sep/(aursun*(mt+mass(3-k))))
oorb = twopi/tb
endif
snova = .true.
endif
*
if(kw.ne.kstar(k))then
change = .true.
mass(k) = mt
dtmi(k) = 0.01d0
if(kw.eq.15)then
kstar(k) = kw
goto 135
endif
mass0(k) = m0
epoch(k) = tphys - age
if(kw.gt.6.and.kstar(k).le.6)then
bsymb = .false.
esymb = .false.
endif
endif
*
* Force new NS or BH to have a second period.
*
if(kstar(k).eq.13.or.kstar(k).eq.14)then
if(tphys-epoch(k).lt.tiny)then
ospin(k) = 2.0d+08
jspin(k) = k3*rc*rc*mc*ospin(k)
endif
endif
*
* Set radius derivative for later interpolation.
*
if(ABS(dtm).gt.tiny)then
rdot(k) = ABS(rm - rad(k))/dtm
else
rdot(k) = 0.d0
endif
*
* Base new time scale for changes in radius & mass on stellar type.
*
dt = dtmi(k)
CALL deltat(kw,age,tm,tn,tscls,dt,dtr)
*
* Choose minimum of time-scale and remaining interval.
*
dtmi(k) = MIN(dt,dtr)
*
* Save relevent solar quantities.
*
aj(k) = age
kstar(k) = kw
rad(k) = rm
lumin(k) = lum
massc(k) = mc
radc(k) = rc
menv(k) = me
renv(k) = re
k2str(k) = k2
tms(k) = tm
tbgb(k) = tscls(1)
*
* Check for blue straggler formation.
*
if(kw.le.1.and.tm.lt.tphys.and..not.bss)then
bss = .true.
jp = MIN(80,jp + 1)
bpp(jp,1) = tphys
bpp(jp,2) = mass(1)
bpp(jp,3) = mass(2)
bpp(jp,4) = float(kstar(1))
bpp(jp,5) = float(kstar(2))
bpp(jp,6) = sep
bpp(jp,7) = ecc
bpp(jp,8) = rad(1)/rol(1)
bpp(jp,9) = rad(2)/rol(2)
bpp(jp,10) = 14.0
endif
*
6 continue
*
if(.not.sgl)then
*
* Determine the mass ratios.
*
do 506 , k = 1,2
q(k) = mass(k)/mass(3-k)
506 continue
*
* Determine the Roche lobe radii and adjust the radius derivative.
*
do 507 , k = 1,2
rol(k) = rl(q(k))*sep
if(ABS(dtm).gt.tiny)then
rdot(k) = rdot(k) + (rol(k) - rol0(k))/dtm
rol0(k) = rol(k)
endif
507 continue
else
do 508 , k = kmin,kmax
rol(k) = 10000.d0*rad(k)
508 continue
endif
*
if((tphys.lt.tiny.and.ABS(dtm).lt.tiny.and.
& (mass2i.lt.0.1d0.or..not.sgl)).or.snova)then
jp = MIN(80,jp + 1)
bpp(jp,1) = tphys
bpp(jp,2) = mass(1)
bpp(jp,3) = mass(2)
bpp(jp,4) = float(kstar(1))
bpp(jp,5) = float(kstar(2))
bpp(jp,6) = sep
bpp(jp,7) = ecc
bpp(jp,8) = rad(1)/rol(1)
bpp(jp,9) = rad(2)/rol(2)
bpp(jp,10) = 1.0
if(snova)then
bpp(jp,10) = 2.0
dtm = 0.d0
goto 4
endif
endif
*
if((isave.and.tphys.ge.tsave).or.iplot)then
if(sgl.or.(rad(1).lt.rol(1).and.rad(2).lt.rol(2)).
& or.tphys.lt.tiny)then
ip = ip + 1
bcm(ip,1) = tphys
bcm(ip,2) = float(kstar(1))
bcm(ip,3) = mass0(1)
bcm(ip,4) = mass(1)
bcm(ip,5) = log10(lumin(1))
bcm(ip,6) = log10(rad(1))
teff1 = 1000.d0*((1130.d0*lumin(1)/
& (rad(1)**2.d0))**(1.d0/4.d0))
bcm(ip,7) = log10(teff1)
bcm(ip,8) = massc(1)
bcm(ip,9) = radc(1)
bcm(ip,10) = menv(1)
bcm(ip,11) = renv(1)
bcm(ip,12) = epoch(1)
bcm(ip,13) = ospin(1)
bcm(ip,14) = dmt(1) - dmr(1)
bcm(ip,15) = rad(1)/rol(1)
bcm(ip,16) = float(kstar(2))
bcm(ip,17) = mass0(2)
bcm(ip,18) = mass(2)
bcm(ip,19) = log10(lumin(2))
bcm(ip,20) = log10(rad(2))
teff2 = 1000.d0*((1130.d0*lumin(2)/
& (rad(2)**2.d0))**(1.d0/4.d0))
bcm(ip,21) = log10(teff2)
bcm(ip,22) = massc(2)
bcm(ip,23) = radc(2)
bcm(ip,24) = menv(2)
bcm(ip,25) = renv(2)
bcm(ip,26) = epoch(2)
bcm(ip,27) = ospin(2)
bcm(ip,28) = dmt(2) - dmr(2)
bcm(ip,29) = rad(2)/rol(2)
bcm(ip,30) = tb
bcm(ip,31) = sep
bcm(ip,32) = ecc
if(isave) tsave = tsave + dtp
endif
endif
*
* If not interpolating set the next timestep.
*
if(intpol.eq.0)then
dtm = MAX(1.0d-07*tphys,MIN(dtmi(1),dtmi(2)))
dtm = MIN(dtm,tsave-tphys)
if(iter.eq.0) dtm0 = dtm
endif
if(sgl) goto 98
*
* Set j1 to the donor - the primary
* and j2 to the accretor - the secondary.
*
if(intpol.eq.0)then
if(rad(1)/rol(1).ge.rad(2)/rol(2))then
j1 = 1
j2 = 2
else
j1 = 2
j2 = 1
endif
endif
*
* Test whether Roche lobe overflow has begun.
*
if(rad(j1).gt.rol(j1))then
*
* Interpolate back until the primary is just filling its Roche lobe.
*
if(rad(j1).ge.1.002d0*rol(j1))then
if(intpol.eq.0) tphys00 = tphys
intpol = intpol + 1
if(iter.eq.0) goto 7
if(inttry) goto 7
if(intpol.ge.100)then
WRITE(99,*)' INTPOL EXCEEDED ',mass1i,mass2i,tbi,ecci
goto 140
endif
dr = rad(j1) - 1.001d0*rol(j1)
if(ABS(rdot(j1)).lt.tiny.or.prec)then
goto 7
endif
dtm = -dr/ABS(rdot(j1))
if(ABS(tphys0-tphys).gt.tiny) dtm = MAX(dtm,tphys0-tphys)
if(kstar(1).ne.kw1)then
kstar(1) = kw1
mass0(1) = mass00(1)
epoch(1) = tphys - aj0(1)
endif
if(kstar(2).ne.kw2)then
kstar(2) = kw2
mass0(2) = mass00(2)
epoch(2) = tphys - aj0(2)
endif
change = .false.
else
*
* Enter Roche lobe overflow
*
if(tphys.ge.tphysf) goto 140
goto 7
endif
else
*
* Check if already interpolating.
*
if(intpol.gt.0)then
intpol = intpol + 1
if(intpol.ge.80)then
inttry = .true.
endif
if(ABS(rdot(j1)).lt.tiny)then
prec = .true.
dtm = 1.0d-07*tphys
else
dr = rad(j1) - 1.001d0*rol(j1)
dtm = -dr/ABS(rdot(j1))
endif
if((tphys+dtm).ge.tphys00)then
*
* If this occurs then most likely the star is a high mass type 4
* where the radius can change very sharply or possibly there is a
* discontinuity in the radius as a function of time and HRDIAG
* needs to be checked!
*
dtm = 0.5d0*(tphys00 - tphys0)
dtm = MAX(dtm,1.0d-10)
prec = .true.
endif
tphys0 = tphys
endif
endif
*
* Check for collision at periastron.
*
pd = sep*(1.d0 - ecc)
if(pd.lt.(rad(1)+rad(2)).and.intpol.eq.0) goto 130
*