-
Notifications
You must be signed in to change notification settings - Fork 25
/
main.c
5332 lines (4559 loc) · 166 KB
/
main.c
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
/***************************************************************************
* Copyright (C) 2009 by Andreas H.W. Kuepper *
* *
* This program is free software; you can redistribute it and/or modify *
* it under the terms of the GNU General Public License as published by *
* the Free Software Foundation; either version 2 of the License, or *
* (at your option) any later version. *
* *
* This program is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
* GNU General Public License for more details. *
* *
* You should have received a copy of the GNU General Public License *
* along with this program; if not, write to the *
* Free Software Foundation, Inc., *
* 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. *
***************************************************************************/
/***************************************************************************
* Compile using the command: cc -lm -o mcluster mcluster.c *
* or use the Makefile, i.e. type: make mcluster or make mcluster_sse *
***************************************************************************/
/***************************************************************************
* Main array with stellar parameters: *
* star[][0] = mass(t = epoch) *
* star[][1-3] = {x, y, z} *
* star[][4-6] = {vx, vy, vz} *
* star[][7] = mass(t = 0) *
* star[][8] = kstar, stellar type in (SSE/BSE) *
* star[][9] = epoch1, age within evolutionary phase (SSE/BSE) *
* star[][10] = ospin, spin of star (SSE/BSE) *
* star[][11] = rstar, radius of star (SSE/BSE) *
* star[][12] = lstar, luminosity (SSE/BSE) *
* star[][13] = epochstar, age of star (SSE/BSE) *
* star[][14] = Zstar, metallicity of star *
***************************************************************************/
#include<stdio.h>
#include<fenv.h>
#include<stdlib.h>
#include<math.h>
#include<time.h>
#include<string.h>
#include<sys/stat.h>
#include<getopt.h>
#ifdef GPU
#include <cuda.h>
#include <cuda_runtime.h>
#endif
#include "main.h"
//use OpenMP if not specified otherwise
#ifndef NOOMP
#include<omp.h>
#endif
int main (int argv, char **argc) {
/*******************
* Input variables *
*******************/
//Basic physical parameters
int N = 0; //Number of stars, Mcl will be set to 0 if specified!
double Mcl = 1000.0; //Total mass of the cluster, only used when N is set to 0, necessary for usage of maximum stellar mass relation of Weidner & Kroupa (2006)
int profile = 0; //Density profile; =0 Plummer model, =1 King model (based on king0.f by D.C. Heggie), =2 mass segregated Subr model (Subr et al. 2007), =3 2-dimensional EFF template (Elson, Fall & Freeman 1987) or Nuker template (Lauer et al. 1995), =-1 no density gradient
double W0 = 5.0; //King's W0 paramter [0.3-12.0]
double S = 0.0; //Fraction of mass segregation for profile; =0.0 unsegregated, =1.0 completely segregated (take maximally S=0.99 for profile=2)
double D = 3.0; //Fractal dimension; =3.0 unfractal, =2.6 2/8 fractal, =2.0 4/8 fractal, =1.6 6/8 fractal, (2^D children per parent, following Goodwin & Whitworth 2004);
double Q = 0.5; //Initial virial ratio; =0.5 virial equilibrium, <0.5 collapsing, >0.5 expanding
double Rh = 0.8; //Half-mass radius [pc], ignored if profile = 3, set =-1 for using Marks & Kroupa (2012) Mcl-Rh relation
double gamma[3] = {2.0, 0.0, 2.0}; //Power-law slopes of EFF/Nuker templates (outer slope, inner slope, transition); set gamma[1] = 0.0 and gamma[2] = 2.0 for EFF (profile = 3)
double a = 1.0; //Scale radius of EFF/Nuker template (profile = 3) [pc]
double Rmax = 100.0; //Cut-off radius for EFF/Nuker template (profile = 3) [pc]
double tcrit = 100.0; //Simulation time [Myr]
int tf = 3; //Tidal field: =0 no tidal field, =1 Near-field approximation, =2 point-mass galaxy, =3 Allen & Santillan (1991) MW potential (or Sverre's version of it)
double RG[3] = {8500.0,0.0,0.0}; //Initial Galactic coordinates of the cluster [pc]
double VG[3] = {0.0,220.0,0.0}; //Initial velocity of the cluster [km/s]
//Mass function parameters
int mfunc = 1; //0 = single mass stars; 1 = use Kroupa (2001) mass function; 2 = use multi power law (based on mufu.c by L.Subr)
double single_mass = 1.0; //Stellar mass in case of single-mass cluster
double mlow = 0.08; //Lower mass limit for mfunc = 1 & mfunc = 4
double mup = 100.0; //Upper mass limit for mfunc = 1 & mfunc = 4
double alpha[MAX_AN] = {-1.35, -2.35, -2.7, 0.0, 0.0}; //alpha slopes for mfunc = 2
double mlim[MAX_MN] = {0.08, 0.5, 4.0, 100.0, 0.0, 0.0}; //mass limits for mfunc = 2
double alpha_L3 = 2.3; //alpha slope for mfunc = 4 (L3 mass function, Maschberger 2012)
double beta_L3 = 1.4; //beta slope for mfunc = 4
double mu_L3 = 0.2; //mu parameter for mfunc = 4
int weidner = 0; //Usage of Weidner & Kroupa (2006) relation for most massive star; =0 off, =1 on
int mloss = 3; //Stellar evolution; 0 = off, 3 = Eggleton, Tout & Hurley [KZ19]
int remnant = 1; //Use random kick velocity and present-day escape velocity to determine retention of compact remnants & evolved binary components (only for SSE/BSE version); =0 off, =1 on
double epoch = 0.0; //Age of the cluster, i.e. star burst has been ... Myr before [e.g. 1000.0, default = 0.0] [needs special compiling and SSE by Hurley, Pols & Tout]
double Z = 0.02; //Metallicity [0.0001-0.03, 0.02 = solar]
double FeH = -1.41; //Metallicity [Fe/H], only used when Z is set to 0
int prantzos = 0; //Usage of Prantzos (2007) relation for the life-times of stars. Set upper mass limit to Lifetime(mup) >= epoch
//Binary parameters
int nbin = 0; //Number of primordial binary systems
double fbin = 0.0; //Primordial binary fraction, number of binary systems = 0.5*N*fbin, only used when nbin is set to 0
int pairing = 3; //Pairing of binary components; 0= random pairing, 1= ordered pairing for components with masses M>msort, 2= random but separate pairing for components with masses m>Msort; 3= Use period distribution for M>msort from Sana et al. (2011,2012) and Oh et al. (2015).
double msort = 5.0; //Stars with masses > msort will be sorted and preferentially paired into binaries if pairing = 1
int adis = 3; //Semi-major axis distribution; 0= flat ranging from amin to amax, 1= based on Kroupa (1995) period distribution, 2= based on Duquennoy & Mayor (1991) period distribution, 3= based on Kroupa (1995) period distribution
int OBperiods = 2; //1: Use period distribution for massive binaries with M_primary > msort from Sana & Evans (2011); 2: Uniform distribution of mass ratio (0.1<q<1.0) for m>Msort and random pairing for remaining (Kiminki & Kobulnicky 2012; Sana et al., 2012; Kobulnicky et al. 2014; Oh, S., Kroupa, P., & Pflamm-Altenburg, J. (2015) period distribution for M>Msort (implemented by Long Wang)
double amin = 0.0001; //Minimum semi-major axis for adis = 0 [pc]
double amax = 0.01; //Maximum semi-major axis for adis = 0 [pc]
#ifdef SSE
int eigen = 0; //Use Kroupa (1995) eigenevolution for pre-main sequence short-period binaries; =0 off, =1 on [use either eigenevolution or BSE; BSE recommended when using SSE]
int BSE = 1; //Apply binary star evolution using BSE (Hurley, Tout & Pols 2002) =0 off, =1 on [use either eigenevolution or BSE; BSE recommended when using SSE]
#else
int eigen = 2; // 1. Use Kroupa (1995) eigenevolution for pre-main sequence short-period binaries; 2. Use modified eigenevolution from Belloni et al. (2017); 0. Off
int BSE = 0; //Apply binary star evolution using BSE (Hurley, Tout & Pols 2002) [needs special compiling and BSE]; =0 off, =1 on [use either eigenevolution or BSE; BSE recommended when using SSE]
#endif
//Gas parameters (only used for Nbody6 input)
double extmass = 0.0; //external Plummer (gas) sphere mass [Msun]
double extrad = 0.0; //external Plummer (gas) sphere scale factor [pc]
double extdecay = 0.0; //decay time for gas expulsion (set 0 for no decay) [Myr]
double extstart = 0.0; //delay time for start of gas expulsion [Myr]
//Code parameters
int code = 3; //Nbody version: =0 Nbody6, =1 Nbody4, =2 Nbody6 custom, =3 only create output list of stars, =4 Nbody7 (not yet fully functional), =5 Nbody6++GPU
unsigned int seed = 0; //Number seed for random number generator; =0 for randomization by local time
char *output = "test"; //Name of output files
double dtadj = 1.0; //DTADJ [N-body units (Myr in Nbody6 custom)], energy-check time step
double dtout = 1.0; //DELTAT [N-body units (Myr in Nbody6 custom)], output interval, must be multiple of DTADJ
double dtplot = 1.0; //DTPLOT [N-body units (Myr in Nbody6 custom)], output of HRdiagnostics, should be multiple of DTOUT, set to zero if output not desired
int gpu = 0; //Use of GPU, 0= off, 1= on
int regupdate = 0; //Update of regularization parameters during computation; 0 = off, 0 > on
int etaupdate = 0; //Update of ETAI & ETAR during computation; 0 = off, 0 > on
int esc = 2; //Removal of escapers; 0 = no removal, 1 = regular removal at 2*R_tide; 2 = removal and output in ESC
int units = 1; //Units of McLuster output; 0= Nbody-Units, 1= astrophysical units
//McLuster internal parameters
int match = 1; //Make cluster half-mass radius exactly match the expected/desired value; =0 off, =1 on (recommended)
int symmetry = 1; //Force spherical symmetry for fractal clusters; =0 off, =1 on (recommended)
int check = 0; //Make energy check at end of McLuster; =0 off, =1 on
int create_radial_profile = 1; //Creates a radial density profile and prints it to the screen; =0 off, =1 on
int create_cumulative_profile = 1; //Creates a radial cumulative profile and prints it to the screen; =0 off, =1 on
double Rgal = 10000.0; //Distance of cluster from sun for artificial CMD with observational errors [pc]
double Zsun = 0.02; //Solar metallicity
int NNBMAX_NBODY6 = 500; //Maximum number of neighbours allowed in NBODY6
double upper_IMF_limit = 150.0; //Maximum stellar mass allowed in McLuster [Msun]
int an = 0; //Counter for number of alpha slopes for mfunc = 2
int mn = 0; //Counter for number of mass limits for mfunc = 1, 2 & 4
int xn = 0; //Counter for components of galactocentric radius vector
int vn = 0; //Counter for components of cluster velocity vector
int xx = 0; //Counter for external potential input parameters
int gn = 0; //Counter for EFF/Nuker profile parameters
double extgas[4]; //Input array for external potential parameters
feenableexcept(FE_DIVBYZERO | FE_INVALID | FE_OVERFLOW);
//SSE internal parameters (see Hurley, Pols & Tout 2000)
value1_.neta = 0.5; //Reimers mass-loss coefficent (neta*4x10^-13; 0.5 normally)
value1_.bwind = 0.0; //Binary enhanced mass loss parameter (inactive for single)
value1_.hewind = 1.0; //Helium star mass loss factor (1.0 normally)
value1_.mxns = 3.0; //Maximum NS mass (1.8, nsflag=0; 3.0, nsflag=1)
points_.pts1 = 0.05; //Time-step parameter in evolution phase: MS (0.05)
points_.pts2 = 0.01; //Time-step parameter in evolution phase: GB, CHeB, AGB, HeGB (0.01)
points_.pts3 = 0.02; //Time-step parameter in evolution phase: HG, HeMS (0.02)
value4_.sigma = 190.0; //Kick velocities
value4_.bhflag = 1; //bhflag > 0 allows velocity kick at BH formation
//BSE internal parameters (see Hurley, Pols & Tout 2002)
flags_.ceflag = 0; //ceflag > 0 activates spin-energy correction in common-envelope (0) #defunct#
flags_.tflag = 1; //tflag > 0 activates tidal circularisation (1)
flags_.ifflag = 0; //ifflag > 0 uses WD IFMR of HPE, 1995, MNRAS, 272, 800 (0)
flags_.nsflag = 1; //nsflag > 0 takes NS/BH mass from Belczynski et al. 2002, ApJ, 572, 407 (1)
flags_.wdflag = 1; //wdflag > 0 uses modified-Mestel cooling for WDs (0)
value5_.beta = 1.0/8.0; //beta is wind velocity factor: proportional to vwind**2 (1/8)
value5_.xi = 1.0; //xi is the wind accretion efficiency factor (1.0)
value5_.acc2 = 3.0/2.0; //acc2 is the Bondi-Hoyle wind accretion factor (3/2)
value5_.epsnov = 0.001; //epsnov is the fraction of accreted matter retained in nova eruption (0.001)
value5_.eddfac = 1.0; //eddfac is Eddington limit factor for mass transfer (1.0)
value5_.gamma = -1.0; //gamma is the angular momentum factor for mass lost during Roche (-1.0)
value2_.alpha1 = 1.0; //alpha1 is the common-envelope efficiency parameter (1.0)
value2_.lambda = 0.5; //lambda is the binding energy factor for common envelope evolution (0.5)
//Command line input
int option;
while ((option = getopt(argv, argc, "N:M:P:W:R:r:c:g:S:D:T:Q:C:A:O:G:o:f:a:m:B:b:p:s:t:e:Z:X:x:V:u:h:?")) != -1) switch (option)
{
case 'N': N = atoi(optarg); Mcl = 0.0; break;
case 'M': Mcl = atof(optarg); N = 0; break;
case 'P':
profile = atoi(optarg);
if ((profile >= 0 && profile <= 3) || profile == -1)
break;
else {
printf("\nError: Density profile (P) needs to "
"be 0, 1, 2, 3 or -1, %d was given\n", profile);
exit(1);
}
case 'W':
W0 = atof(optarg);
if (W0 >= 0.3 && W0 <= 12.0)
break;
else {
printf("\nError: King's parameter (W) needs to be between 0.3 "
"and 12.0, %f was given\n", W0);
exit(1);
}
case 'R':
Rh = atof(optarg);
break;
case 'r':
a = atof(optarg);
break;
case 'c':
Rmax = atof(optarg);
break;
case 'g':
if (gn < 3) {
gamma[gn] = atof(optarg);
gn++;
break;
}
case 'S':
S = atof(optarg);
if (S >= 0 && S <= 0.99)
break;
else {
printf("\nError: Fraction of mass segretation (S) needs to "
"between 0 and 0.99, %f was given\n", S);
exit(1);
}
case 'D':
D = atof(optarg);
if (D == 3.0 || D == 2.6 || D == 2.0 || D == 1.6)
break;
else {
printf("\nError: Fractal dimension (D) needs to "
"be 1.6, 2.0, 2.6 or 3.0, %f was given\n", D);
exit(1);
}
case 'T':
tcrit = atof(optarg);
break;
case 'Q':
Q = atof(optarg);
break;
case 'C':
code = atoi(optarg);
if (code >= 0 && code <= 5)
break;
else {
printf("\nError: Nbody version (C) needs to "
"be between 0 and 5, %d was given\n", code);
exit(1);
}
case 'A':
dtadj = atof(optarg);
break;
case 'O':
dtout = atof(optarg);
break;
case 'G':
gpu = atoi(optarg);
if (gpu == 0 || gpu == 1)
break;
else {
printf("\nError: Use of GPU (G) needs to "
"be between 0 or 1, %d was given\n", gpu);
exit(1);
}
case 'o':
output = optarg;
break;
case 'f':
mfunc = atoi(optarg);
if (mfunc >= 0 && mfunc <= 2)
break;
else {
printf("\nError: Mass function (f) needs to "
"be between 0 and 2, %d was given\n", mfunc);
exit(1);
}
case 'a' :
if (an < MAX_AN) {
alpha[an] = atof(optarg);
if (an == 0) alpha_L3 = atof(optarg);
if (an == 1) beta_L3 = atof(optarg);
if (an == 2) mu_L3 = atof(optarg);
an++;
break;
} else { printf("\nError: Number of alphas exceeded maximum limit of %d\n", MAX_AN); return 1; }
case 'm' :
if (mn < MAX_MN) {
mlim[mn] = atof(optarg);
if (mn == 0) mlow = atof(optarg);
if (mn == MAX_MN-1) mup = atof(optarg);
mn++;
break;
} else { printf("\nError: Number of mass params exceded maximum limit of %d\n", MAX_MN); return 1; }
case 'B': nbin = atoi(optarg); break;
case 'b': fbin = atof(optarg); break;
case 'p':
pairing = atoi(optarg);
if (pairing >= 0 && pairing <= 3)
break;
else {
printf("\nError: Pairing of binary components (p) needs to "
"be between 0 and 3, %d was given\n", pairing);
exit(1);
}
case 's': seed = atoi(optarg); break;
case 't':
tf = atoi(optarg);
if (tf >= 0 && tf <= 3)
break;
else {
printf("\nError: Tidal field (t) needs to "
"be between 0 and 3, %d was given\n", tf);
exit(1);
}
case 'e': epoch = atof(optarg); break;
case 'Z':
Z = atof(optarg);
if (Z >= 0.0001 && Z <= 0.03)
break;
else {
printf("\nError: Metallicity (Z) needs to "
"be between 0.0001 and 0.03, %f was given\n", Z);
exit(1);
}
case 'X' :
if (xn < 3) { RG[xn] = atof(optarg); xn++; break; }
case 'V' :
if (vn < 3) { VG[vn] = atof(optarg); vn++; break; }
case 'x' :
if (xx < 4) { extgas[xx] = atof(optarg); xx++; break; }
case 'u': units = atoi(optarg); break;
case ':':
case 'h': help(msort); return 1;
case '?': help(msort); return 1;
};
// some parameters ought to be checked for validity
if (S<0 || S>1) {
printf("Bad value of S=%g\n",S);
exit(1);
}
if (mn-1 > 0) mup = mlim[mn-1];
if (seed) srand48(seed); //initialize random number generator by seed
else {
seed = (unsigned) time(NULL); //initialize random number generator by local time
seed %= 100000;
srand48(seed);
}
//print summary of input parameters to .info file
info(output, N, Mcl, profile, W0, S, D, Q, Rh, gamma, a, Rmax, tcrit, tf, RG, VG, mfunc, single_mass, mlow, mup, alpha, mlim, alpha_L3, beta_L3, mu_L3, weidner, mloss, remnant, epoch, Z, prantzos, nbin, fbin, pairing, msort, adis, amin, amax, eigen, BSE, extmass, extrad, extdecay, extstart, code, seed, dtadj, dtout, dtplot, gpu, regupdate, etaupdate, esc, units, match, symmetry, OBperiods);
/*********
* Start *
*********/
printf("\n-----START---- \n");
#ifdef NOOMP
clock_t t1, t2;
t1 = clock(); //start stop-watch
#else
double t1, t2;
#pragma omp parallel
{
if (omp_get_thread_num() == 0)
if (omp_get_num_threads() > 1) printf("\n\nUsing OpenMP with %d threads\n", omp_get_num_threads());
t1 = omp_get_wtime(); //start stop-watch
}
#endif
printf ("\n\nRandom seed = %i\n\n\n", seed);
if (seed) value3_.idum = seed; //idum is the random number seed used in the kick routine.
else value3_.idum = 10000000.0*drand48();
int i,j;
double M; //Total mass [M_sun]
double mmean; //Mean stellar mass [M_sun]
int NNBMAX; //Maximum neighbour number (Nbody6 only)
double RS0; //Initial radius of neighbour sphere [pc], Nbody6 only
double rtide; //Tidal radius [pc]
double omega; //Angular velocity of cluster around the galaxy
double rvir; //Virial radius [pc]
double cmr[7]; //For CoM correction
double rking, rplummer; //King-, Plummer radius
double MMAX; //most massive star
double tscale; //time-scale factor
double ekin = 0.0; //kinetic energy
double epot = 0.0; //potential energy
double sigma = 0.0; //velocity dispersion
int bin; //KZ(22) parameter (binaries yes/no)
int sse; //(evolved stellar population yes/no)
double submass[MAX_AN], subcount[MAX_AN], norm[MAX_AN], N_tmp, M_tmp; //mass function parameters for mfunc = 2
double Rh2D, Rh3D; //actual 2D/3D half-mass radius of the model
if (profile == 3) Rh = a; //set half-mass radius temporarily to scale radius for computation of escape velocity
if ((Mcl) && (N)) {
printf("\nWARNING: specify either Mcl (-M) or N (-N)!\n\n");
exit (1);
} else if ((!Mcl) && (mfunc == 3)) {
printf("\nWARNING: specify Mcl (-M) when using optimal sampling (-f 3)!\n\n");
exit (1);
}
if (xx > 0) {
if ((xx == 4) && (extgas[2])) {
extmass = extgas[0];
extrad = extgas[1];
extdecay = extgas[2];
extstart = extgas[3];
} else {
printf("\nWARNING: Insufficient or incorrect parameters specified for external Plummer potential!\n\n");
exit (1);
}
}
/***********************
* Generate star array *
***********************/
int columns = 15;
const int NMAX = N*2>1500000?N*2:1500000; //Maximum number of stars & orbits allowed in McLuster
double **star;
star = (double **)calloc(NMAX,sizeof(double *));
for (j=0;j<NMAX;j++){
star[j] = (double *)calloc(columns,sizeof(double));
if (star[j] == NULL) {
printf("\nMemory allocation failed!\n");
return 0;
}
}
/*******************************************
* Evaluate Z from [Fe/H] if Z is set to 0 *
*******************************************/
if (!Z) {
Z = pow(10.0, 0.977*FeH)*Zsun; //Bertelli, Bressan, Chiosi, Fagotto, Nasi, 1994, A&AS, 106, 275
printf("\nUsing Bertelli et al. (1994) relation to convert FeH = %.3f into Z = %.3f\n", FeH, Z);
}
/**********************************
* Calculate maximum stellar mass *
**********************************/
if (mfunc == 3) weidner = 1; //always use Weidner relation when using optimal sampling!
if (!N && weidner && mfunc) {
mup = upper_IMF_limit;
printf("\nUsing maximum stellar mass-cluster mass relation for upper stellar mass limit\n");
printf("\n(Weidner & Kroupa 2006, Pflamm-Altenburg & Kroupa 2007)\n");
//analytic fit to the observational data from Pflamm-Altenburg & Kroupa (2007), implementation by M. Kruckow
MMAX = pow(10,2.56*log10(Mcl)*pow(pow(3.82,9.17)+pow(log10(Mcl),9.17),-1/9.17)-0.38);
/*
//three-part power-law fit to the observational data
if (Mcl < 1000.0) {
MMAX = (log10(Mcl)*0.540563175) - 0.14120167;
MMAX = pow(10.0,MMAX);
} else if (Mcl < 3300.0) {
MMAX = (log10(Mcl)*0.19186051) + 0.9058611;
MMAX = pow(10.0,MMAX);
} else {
MMAX = (log10(Mcl)*0.360268003) + 0.313342031;
MMAX = pow(10.0,MMAX);
}*/
} else {
MMAX = mup;
}
if (mfunc && epoch && prantzos) {
printf("\nUsing Prantzos (2007) relation to reduce upper mass limit to Lifetime(mup) > epoch\n");
while (Lifetime(MMAX) < sqrt(pow(epoch,2))) {
MMAX -= 0.01;
}
}
/*******************
* Generate masses *
*******************/
printf("\n\n-----GENERATE MASSES----- \n");
//This loop has to be called multiple times with different N (or Mcl), Z and epoch in order to create multiple stellar populations
//make sure that epoch and Z are stored together with the stars, then concatenate all arrays
//for now, all Z and epochs are the same, i.e. a single stellar population
if (mfunc == 1) {
printf("\nMaximum stellar mass set to: %.2f\n", MMAX);
generate_m1(&N, star, mlow, mup, &M, &mmean, MMAX, Mcl, epoch, Z, Rh, remnant);
} else if (mfunc == 2) {
if (mn) {
for (i = mn+1; i < MAX_MN; i++) mlim[i] = 0.0;
} else {
for (i=0; i<MAX_MN; i++) {
if (mlim[i]) mn++;
}
}
if (an) {
for (i = an+1; i < MAX_AN; i++) alpha[i] = 0.0;
} else {
for (i=0; i<MAX_AN; i++) {
if (alpha[i]) an++;
}
}
if (an >= mn) an = mn - 1;
mn = an + 1;
if (!mn){
printf("\nError: at least one mass limit has to be specified\n");
return 1;
} else if (mn == 1) {
single_mass = mlim[0];
printf("\nSetting stellar masses to %g solar mass\n",single_mass);
if (!N) N = Mcl/single_mass;
for (j=0;j<N;j++) star[j][0] = 1.0/N;
mmean = single_mass;
M = N*mmean;
printf("\nM = %g\n", M);
} else {
printf("\nMaximum stellar mass set to: %.2f\n",MMAX);
norm[an-1] = 1.; //normalization factor of integral
N_tmp = subcount[an-1] = subint(mlim[an-1], mlim[an], alpha[an-1] + 1.); //integrated number of stars in interval [mlim[an-1]:mlim[an]]
M_tmp = submass[an-1] = subint(mlim[an-1], mlim[an], alpha[an-1] + 2.); //integrated mass of stars in interval [mlim[an-1]:mlim[an]]
for (i = an - 2; i >= 0; i--) {
norm[i] = norm[i+1] * pow(mlim[i+1], alpha[i+1] - alpha[i]);
subcount[i] = norm[i] * subint(mlim[i], mlim[i+1], alpha[i] + 1.);
N_tmp += subcount[i];
submass[i] = norm[i] * subint(mlim[i], mlim[i+1], alpha[i] + 2.);
M_tmp += submass[i];
}
generate_m2(an, mlim, alpha, Mcl, M_tmp, subcount, &N, &mmean, &M, star, MMAX, epoch, Z, Rh, remnant);
}
} else if (mfunc == 3) {
printf("\nMaximum stellar mass set to: %.2f\n",MMAX);
generate_m3(&N, star, mlow, mup, &M, &mmean, MMAX, Mcl);
randomize(star, N);
} else if (mfunc == 4) {
printf("\nMaximum stellar mass set to: %.2f\n", MMAX);
printf("\nUsing L3 IMF (Maschberger 2012)\n");
generate_m4(&N, star, alpha_L3, beta_L3, mu_L3, mlow, mup, &M, &mmean, MMAX, Mcl, epoch, Z, Rh, remnant);
} else {
printf("\nSetting stellar masses to %.1f solar mass\n", single_mass);
if (!N) N = Mcl/single_mass;
for (j=0;j<N;j++) {
star[j][0] = single_mass;
star[j][7] = single_mass;
star[j][8] = 0;
star[j][9] = 0.0;
star[j][10] = 0.0;
star[j][11] = 0.0;
star[j][12] = 0.0;
star[j][13] = 0.0;
star[j][14] = 0.0;
}
mmean = single_mass;
M = N*mmean;
printf("\nM = %g\n", M);
mloss = 0;
}
//set all stars to the same metallicity and age for now
double epochstar, zstar;
epochstar = 0.0; //age compared to the oldest stars in the cluster [Myr]
zstar = Z;
for (j=0;j<N;j++) {
star[j][13] = epochstar;
star[j][14] = zstar;
}
//Pair binary masses and convert to centre-of-mass particles
int Nstars;
if (!nbin) nbin = 0.5*N*fbin;
double **mbin; //component mass & stellar evol parameter array
mbin = (double **)calloc(nbin,sizeof(double *));
for (j=0;j<nbin;j++){
mbin[j] = (double *)calloc(20,sizeof(double));
if (mbin[j] == NULL) {
printf("\nMemory allocation failed!\n");
return 0;
}
}
Nstars = N;
if (nbin) {
printf("\nPreparing binary components.\n");
//Specify component pairing
if (pairing) {
if (pairing == 1) printf("\nApplying ordered pairing for stars with masses > %.1f Msun.\n",msort);
else if (pairing == 2) printf("\nApplying random pairing for stars with masses > %.1f Msun.\n",msort);
else if (pairing == 3) printf("\nApplying Sana et al. (2012) mass atio distribution for stars with masses > %.1f Msun.\n",msort);
order(star, N, M, msort, pairing);
} else {
randomize(star, N);
printf("\nApplying random pairing.\n");
}
N -= nbin;
for (j=0;j<nbin;j++) {
mbin[j][0] = star[2*j][0]+star[2*j+1][0]; //system mass
mbin[j][1] = star[2*j][0]; //primary mass
mbin[j][2] = star[2*j+1][0]; //secondary mass
mbin[j][3] = star[2*j][7]; //primary m0
mbin[j][4] = star[2*j][8]; //primary kw
mbin[j][5] = star[2*j][9]; //primary epoch1
mbin[j][6] = star[2*j][10]; //primary spin
mbin[j][7] = star[2*j][11]; //primary r
mbin[j][8] = star[2*j][12]; //primary lum
mbin[j][9] = star[2*j+1][7]; //secondary m0
mbin[j][10] = star[2*j+1][8]; //secondary kw
mbin[j][11] = star[2*j+1][9]; //secondary epoch1
mbin[j][12] = star[2*j+1][10]; //secondary spin
mbin[j][13] = star[2*j+1][11]; //secondary r
mbin[j][14] = star[2*j+1][12]; //secondary lum
mbin[j][15] = 1000+j; //identifier
mbin[j][16] = star[2*j][13]; //primary epochstar
mbin[j][17] = star[2*j+1][13]; //secondary epochstar
mbin[j][18] = star[2*j][14]; //primary zstar
mbin[j][19] = star[2*j+1][14]; //secondary zstar
star[2*j][0] += star[2*j+1][0]; //system mass
star[2*j+1][0] = 0.0;
star[2*j][7] = 1000+j; //identifier
star[2*j+1][7] = 0.0; //identifier
star[2*j][12] += star[2*j+1][12]; //system luminosity
star[2*j+1][12] = 0.0;
}
order(star, Nstars, M, 0.0, 0);
randomize(star, N);
}
//prepare mass segregation
double mlowest = MMAX; //search lowest mass star
double mhighest = 0; //search highest mass star
double mmeancom = 0.0;
for (i=0;i<N;i++) {
star[i][0] /= M; //scale masses to Nbody units
mmeancom += star[i][0];
if (star[i][0] < mlowest) mlowest = star[i][0];
if (star[i][0] > mhighest) mhighest = star[i][0];
}
mmeancom /= N;
int Nseg = ceil(N*mmeancom/mlowest); //number of necessary pos & vel pairs for Baumgardt et al. (2008) mass segregation routine
int Nunseg = N;
double *Mcum;
Mcum = (double *)calloc(N,sizeof(double));
if ((S) && !(profile == 2)) {//sort masses when mass segregation parameter > 0
printf("\nApplying mass segregation with S = %f\n",S);
segregate(star, N, S);
for (i=0;i<N;i++) {//calculate cumulative mass function Mcum
Mcum[i] = 0.0;
for (j=0;j<=i;j++) Mcum[i] = Mcum[i] + star[j][0];
}
N = Nseg;
}
/*************************************
* Generate positions and velocities *
*************************************/
printf("\n\n-----GENERATE POSITIONS & VELOCITIES----- \n");
//calculate half-mass radius according to Marks & Kroupa 2012 if Rh is set to -1
if ((Rh == -1) && (Mcl)) {
Rh = 0.1*pow(Mcl,0.13); //Marks & Kroupa (2012), implementation by M. Kruckow
printf("\nUsing Marks & Kroupa (2012) relation to derive half-mass radius from cluster mass: %g (pc)\n", Rh);
}
//evaluate approximate tidal radius assuming circular orbit
if (tf == 3) {
//in the case of Allen & Santillan potential, assume kappa = 1.4omega (eq. 9 in Kuepper et al. 2010)
omega = sqrt(VG[0]*VG[0]+VG[1]*VG[1]+VG[2]*VG[2])/sqrt(RG[0]*RG[0]+RG[1]*RG[1]+RG[2]*RG[2]);
rtide = pow(G*M/(2.0*omega*omega),1.0/3.0);
} else if (!tf) {
rtide = 1.0E5;
} else if ((tf == 1) && (code == 0 || code == 4 || code == 5)) {
//in case of Sverre's Nbody6 standard tidal field
rtide = pow(1.0*M/(3.0*M1pointmass),1.0/3.0)*8000.0;
} else {
//in the case of a point mass potential or near field approximation
rtide = pow(1.0*M/(3.0*M1pointmass),1.0/3.0)*sqrt(RG[0]*RG[0]+RG[1]*RG[1]+RG[2]*RG[2]);
}
printf("\nApproximate tidal radius: %g (pc)\n", rtide);
//generate scaled pos & vel, postpone scaling for Plummer and King in case of mass segregation
double rhtemp, rvirtemp;
if (profile == 1) {
printf("\nGenerating King model with parameters: N = %i\t W0 = %g\t Rh = %.3f\t D = %.2f\n",N, W0, Rh, D);
generate_king(N, W0, star, &rvirtemp, &rhtemp, &rking, D, symmetry);
} else if (profile == 2) {
N = Nunseg;
printf("\nGenerating segregated Subr model with parameters: N = %i\t S = %g\t Rh = %.3f\n",N, S, Rh);
rvir = Rh/0.76857063065978; //(value provided by L. Subr)
generate_subr(N, S, star, rtide, rvir);
printf ("\nrvir = %.5f\t rh = %.5f\t rtide = %.5f (pc)\n", rvir, Rh, rtide);
} else if (profile == 3) {
if (gamma[1] == 0.0 && gamma[2] == 2.0) printf("\nGenerating EFF model with parameters: N = %i\t a = %.1f\t gamma = %.3f\t Rmax = %.1f\t D = %.2f\n",N, a, gamma[0], Rmax, D);
else printf("\nGenerating Nuker model with parameters: N = %i\t a = %.1f\t gamma (outer) = %.3f\t gamma (inner) = %.3f\t transition = %.3f\t Rmax = %.1f\t D = %.2f\n",N, a, gamma[0],gamma[1],gamma[2], Rmax, D);
double p[6];
p[1] = 1.0; //rho0, will be scaled according to Rmax and Mtot
p[2] = a;//scale radius
p[3] = gamma[0];//outer power-law slope (>0.5)
p[4] = gamma[1];//inner power-law slope (0.0 for EFF template)
p[5] = gamma[2];//transition parameter (2.0 for EFF template)
generate_profile(N, star, Rmax, M, p, &Rh, D, symmetry);
printf("\nRh = %.1f pc\n", Rh);
} else if (profile == -1) {
printf("\nGenerating fractal distribution with parameters: N = %i\t Rh = %.3f\t D = %.2f\n", N, Rh, D);
fractalize(D, N, star, 0, symmetry);
rvir = Rh;
} else {
printf("\nGenerating Plummer model with parameters: N = %i\t Rh = %.3f\t D = %.2f\n", N, Rh, D);
rvir = Rh/0.772764;
rplummer = Rh/1.305;
generate_plummer(N, star, rtide, rvir, D, symmetry);
}
//Apply Baumgardt et al. (2008) mass segregation
if (!(profile == 2) && (S)) {
double **m_temp;
m_temp = (double **)calloc(Nunseg,sizeof(double *));
for (j=0;j<Nunseg;j++){
m_temp[j] = (double *)calloc(9,sizeof(double));
if (m_temp[j] == NULL) {
printf("\nMemory allocation failed!\n");
return 0;
}
}
for (i=0;i<Nunseg;i++) {//temporarily store the masses & stellar evol parameters
m_temp[i][0] = star[i][0];
m_temp[i][1] = star[i][7];
m_temp[i][2] = star[i][8];
m_temp[i][3] = star[i][9];
m_temp[i][4] = star[i][10];
m_temp[i][5] = star[i][11];
m_temp[i][6] = star[i][12];
m_temp[i][7] = star[i][13];
m_temp[i][8] = star[i][14];
}
printf("\nOrdering orbits by energy.\n");
energy_order(star, N, Nstars);
int nlow, nhigh, nrandom;
for (i=0;i<Nunseg;i++) {
nhigh = Nseg*Mcum[i];
if (i) {
nlow = Nseg*Mcum[i-1];
} else {
nlow = 0;
}
nrandom = (nhigh-nlow)*drand48()+nlow;
star[i][0] = m_temp[i][0];
star[i][1] = star[nrandom][1];
star[i][2] = star[nrandom][2];
star[i][3] = star[nrandom][3];
star[i][4] = star[nrandom][4];
star[i][5] = star[nrandom][5];
star[i][6] = star[nrandom][6];
star[i][7] = m_temp[i][1];
star[i][8] = m_temp[i][2];
star[i][9] = m_temp[i][3];
star[i][10] = m_temp[i][4];
star[i][11] = m_temp[i][5];
star[i][12] = m_temp[i][6];
star[i][13] = m_temp[i][7];
star[i][14] = m_temp[i][8];
}
for (j=0;j<Nunseg;j++) free (m_temp[j]);
free(m_temp);
N = Nunseg;
}
//CoM correction
printf("\nApplying centre-of-mass correction.\n");
for (j=0; j<7; j++) cmr[j] = 0.0;
for (j=0; j<N; j++) {
for (i=1;i<7;i++)
cmr[i] += star[j][0]*star[j][i];
}
for (j=0; j<N; j++) {
for (i=1;i<7;i++)
star[j][i] -= cmr[i];
}
//apply scaling to Nbody-units
if (profile == 0) {
printf("\nRe-scaling of orbits (dt ~ N^2!)\n");
double ke = 0.0;
double pe = 0.0;
double sx, sv, r2;
#ifdef GPU
gpupot(N,star,&pe);
#ifndef NOOMP
#pragma omp parallel shared(N, star) private(i)
{
#pragma omp for reduction(+: ke) schedule(dynamic)
#endif
for (i=0;i<N;i++) {
ke += star[i][0]*(pow(star[i][4],2)+pow(star[i][5],2)+pow(star[i][6],2));
}
#ifndef NOOMP
}
#endif
// printf("PE_GPU %lg ke %lg\n",pe,ke);
#else
#ifndef NOOMP
#pragma omp parallel shared(N, star) private(i, j, r2)
{
#pragma omp for reduction(+: pe, ke) schedule(dynamic)
#endif
for (i=0;i<N;i++) {
if (i) {
for (j=0;j<i-1;j++) {
r2 = (star[i][1]-star[j][1])*(star[i][1]-star[j][1]) + (star[i][2]-star[j][2])*(star[i][2]-star[j][2]) +(star[i][3]-star[j][3])*(star[i][3]-star[j][3]) ;
pe -= star[i][0]*star[j][0]/sqrt(r2);
}
}
ke += star[i][0]*(pow(star[i][4],2)+pow(star[i][5],2)+pow(star[i][6],2));
}
#ifndef NOOMP
}
#endif
// printf("PE_HOST %lg ke %lg\n",pe,ke);
#endif
rvir = -GNBODY*pow(MNBODY,2)/(2.0*pe);
sx = 1.0/rvir;
ke *= 0.5;
sv = sqrt(4.0*ke);
for (i=0;i<N;i++) {
star[i][1] *= sx;
star[i][2] *= sx;
star[i][3] *= sx;
star[i][4] /= sv;
star[i][5] /= sv;
star[i][6] /= sv;
}
ke = 0;
for (i=0;i<N;i++) {
ke += M*star[i][0]*(pow(star[i][4],2)+pow(star[i][5],2)+pow(star[i][6],2));
}
ke /= 0.5*N;
printf("Dynamical temperature of centre-of-mass particles kT = %lf\n\n",ke);
//make half-mass radius of the system match the desired one
radial_profile(star, N, rvir, M, 0, 0, code, &NNBMAX, &RS0, &Rh2D, &Rh3D, NNBMAX_NBODY6);
if (match) {
printf("\nmeasuring half-mass radius: %.7f \t %.7f (should/is)\nand correcting for this factor\n",Rh, Rh3D);
rvir = rvir *Rh/Rh3D;
}
printf ("\nrvir = %.5f\t rh = %.5f\t rplummer = %.5f\t rtide = %.5f (pc)\n", rvir, Rh, rplummer, rtide);
} else if ((profile == 1) || (profile == 3) || (profile == -1)) {
printf("\nRe-scaling of orbits (dt ~ N^2!)\n");
double pe = 0.0;
double ke = 0.0;
double r2, vscale;
#ifdef GPU
gpupot(N,star,&pe);
#ifndef NOOMP
#pragma omp parallel shared(N, star) private(i)
{
#pragma omp for reduction(+: ke) schedule(dynamic)
#endif
for (i=0;i<N;i++) {
ke += star[i][0]*(pow(star[i][4],2)+pow(star[i][5],2)+pow(star[i][6],2));
}
#ifndef NOOMP
}
#endif
// printf("PE_GPU %lg ke %lg\n",pe,ke);
#else
#ifndef NOOMP
#pragma omp parallel shared(N, star) private(i, j, r2)
{
#pragma omp for reduction(+: pe, ke) schedule(dynamic)
#endif
for (i=0;i<N;i++) {
for (j=0;j<i-1;j++) {
r2 = (star[i][1]-star[j][1])*(star[i][1]-star[j][1]) + (star[i][2]-star[j][2])*(star[i][2]-star[j][2]) +(star[i][3]-star[j][3])*(star[i][3]-star[j][3]) ;
pe -= star[i][0]*star[j][0]/sqrt(r2);
}
ke += star[i][0]*(pow(star[i][4],2)+pow(star[i][5],2)+pow(star[i][6],2));
}
#ifndef NOOMP
}
#endif
// printf("PE_HOST %lg ke %lg\n",pe,ke);
#endif
ke *= 0.5;
rvir = -GNBODY*pow(MNBODY,2)/(2.0*pe);
vscale = sqrt(4.0*ke);
ke = 0;
for (i=0;i<N;i++) {
for (j=0;j<3;j++) {
star[i][j+1] /= rvir;
star[i][j+4] /= vscale;
}
ke += M*star[i][0]*(pow(star[i][4],2)+pow(star[i][5],2)+pow(star[i][6],2));
}
ke /= 0.5*N;
printf("Dynamical temperature of centre-of-mass particles kT = %lf\n\n",ke);
//make half-mass radius of the system match the desired one
radial_profile(star, N, rvir, M, 0, 0, code, &NNBMAX, &RS0, &Rh2D, &Rh3D, NNBMAX_NBODY6);
if (match) {
printf("\nmeasuring half-mass radius: %.7f \t %.7f (should/is)\nand correcting for this factor\n",Rh, Rh3D);
rvir = rvir *Rh/Rh3D;
}
//printf ("\nrvir = %.5f\t rh = %.5f\t rtide = %.5f (pc)\n", rvir, Rh, rtide);
//printf("Edge radius (King units) = %g\t(Nbody units) = %g\n", rking, rking/rvirtemp);
//printf("Core radius (King units) = %g\t(Nbody units) = %g\n\n", 1.0, 1.0/rvirtemp);
//printf("Concentration = %g\n", log10(rking));
} else if (profile == 2) {
double ke = 0;
#ifndef NOOMP
#pragma omp parallel shared(N, star) private(i)
{
#pragma omp for reduction(+: ke) schedule(dynamic)
#endif
for (i=0;i<N;i++) {
ke += M*star[i][0]*(pow(star[i][4],2)+pow(star[i][5],2)+pow(star[i][6],2));
}
#ifndef NOOMP
}
#endif
ke /= 0.5*N;
printf("Dynamical temperature of centre-of-mass particles kT = %lf\n\n",ke);
//make half-mass radius of the system match the desired one
radial_profile(star, N, rvir, M, 0, 0, code, &NNBMAX, &RS0, &Rh2D, &Rh3D, NNBMAX_NBODY6);
printf("\nmeasuring half-mass radius: %.7f \t %.7f (should/is)\nand correcting for this factor\n",Rh, Rh3D);
rvir = rvir *Rh/Rh3D;
}
//Calculate radial density profile, estimate NNBMAX and RS0 (important for Nbody6 only)
radial_profile(star, N, rvir, M, create_radial_profile, create_cumulative_profile, code, &NNBMAX, &RS0, &Rh2D, &Rh3D, NNBMAX_NBODY6);
printf("\nActual half-mass radius of the cluster = (%.4f / %.4f) pc (3D / 2D)\n", Rh3D, Rh2D);
//scale RS0 to nbody units for Nbody6
RS0 /= 1.0*rvir;
/*********************
* Generate Binaries *
*********************/
printf("\n\n-----GENERATE BINARIES----- \n");
if ((!fbin) && (!nbin)) {
printf("\nNo primordial binaries!\n");
} else {
//re-create original array with Nstars (original N) entries
columns = 15;
double **star_temp;
star_temp = (double **)calloc(Nstars,sizeof(double *));
for (j=0;j<Nstars;j++){
star_temp[j] = (double *)calloc(columns,sizeof(double));
if (star_temp[j] == NULL) {
printf("\nMemory allocation failed!\n");
return 0;
}