-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathpyHALMA.py
1854 lines (1592 loc) · 95.6 KB
/
pyHALMA.py
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
# Created by Óscar Monllor-Berbegal
# FIRST REALESE OF pyHALMA 15/02/2024
# Description: This is the main file of the pyHALMA code. It reads MASCLET files and calculates the properties of the haloes.
# WARNINGS:
# 1. When RPS_FLAG or POT_ENERGY_FLAG are activated, they will significantly slow down the code.
# 2. When ESCAPE_CLEANING is activated, it will significantly slow down the code.
# 3. If not FULL_DOMAIN_FLAG, the code will only read the region defined by X1, X2, Y1, Y2, Z1, Z2.
# hence the indices will be referred to the region defined by X1, X2, Y1, Y2, Z1, Z2. Not the whole domain.
# 4. LL calculated in the code as a function of stellar particle mass is EXPERIMENTAL.
# 5. Distance variables are in COMOVING kpc
# 6. RMAX is calculated with respect to the center of mass of the halo, not the density peak or the most bound particle.
# 7. Code should not be executed in Windows, since multiprocessing by Fork is not available (just Spawn).
import time
import numpy as np
import sys
import os
import json
import warnings
warnings.filterwarnings('ignore')
import multiprocessing
import numba
from tqdm import tqdm
from math import acos
from scipy.interpolate import RegularGridInterpolator
from scipy.spatial import KDTree
import gc
########## ########## ########## ########## ##########
########## INPUT PARAMETERS FROM pyHALMA.dat #########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
with open('pyHALMA.dat', 'r') as f:
f.readline() # MAIN BLOCK
f.readline()
FIRST,LAST,STEP = np.array(f.readline().split()[0].split(','), dtype = np.int32)
f.readline()
NX, NY, NZ = np.array(f.readline().split()[0].split(','), dtype = np.int32)
f.readline()
NAMRX, NAMRY, NAMRZ, NLEVELS = np.array(f.readline().split()[0].split(','), dtype = np.int32)
f.readline()
NPALEV = int(f.readline())
f.readline()
ACHE, OMEGA0, T0 = np.array(f.readline().split()[0].split(','), dtype = np.float64)
f.readline()
Z0, L = np.array(f.readline().split()[0].split(','), dtype = np.float64)
f.readline()
FULL_DOMAIN_FLAG = bool(int(f.readline()))
f.readline() #INPUT DOMAIN
X1, X2, Y1, Y2, Z1, Z2 = np.array(f.readline().split()[0].split(','), dtype = np.float64)
f.readline()
MASS_PART = float(f.readline())
f.readline()
f.readline()
f.readline()
LL, = np.array(f.readline().split()[0].split(','), dtype = np.float64)
LL /= 1e3 #to Mpc
f.readline()
CALCULATE_LL_FLAG = bool(int(f.readline()))
f.readline()
PYFOF_FLAG = bool(int(f.readline()))
f.readline()
PARALLEL_FOF = bool(int(f.readline()))
f.readline()
MINP, = np.array(f.readline().split()[0].split(','), dtype = np.int32)
f.readline()
PSC_FLAG = bool(int(f.readline()))
f.readline()
SIG_FIL, Q_FIL = np.array(f.readline().split()[0].split(','), dtype = np.float64)
f.readline()
ESCAPE_CLEANING, FACTOR_V = np.array(f.readline().split()[0].split(','), dtype = np.float64)
ESCAPE_CLEANING = bool(int(ESCAPE_CLEANING))
f.readline()
RPS_FLAG = bool(int(f.readline()))
f.readline()
POT_ENERGY_FLAG = bool(int(f.readline()))
f.readline()
FACTOR_R12_POT = float(f.readline())
f.readline()
BRUTE_FORCE_LIM = int(f.readline())
f.readline()
f.readline()
f.readline()
NCORE = int(f.readline())
f.readline()
PATH_MASCLET_FRAMEWORK = f.readline()[:-1]
f.readline()
f.readline() #SUBSTRUCTURE BLOCK
f.readline()
SUBSTRUCTURE_FLAG = bool(int(f.readline()))
f.readline()
MINP_SUB, DEC_FACTOR, MIN_LL, DIST_FACTOR = np.array(f.readline().split()[0].split(','), dtype = np.float64)
MINP_SUB = int(MINP_SUB)
MIN_LL /= 1e3 #to Mpc
f.readline()
f.readline() #CALIPSO BLOCK
f.readline()
CALIPSO_FLAG = bool(int(f.readline()))
f.readline()
CLIGHT = float(f.readline())
f.readline()
USUN, GSUN, RSUN, ISUN = np.array(f.readline().split()[0].split(','), dtype = np.float64)
f.readline()
METSUN = float(f.readline())
f.readline()
L_START, L_END = np.array(f.readline().split()[0].split(','), dtype = np.float64)
f.readline()
f.readline() #ASOHF BLOCK
f.readline()
DM_FLAG = bool(int(f.readline()))
f.readline()
ASOHF_FLAG = bool(int(f.readline()))
f.readline()
FACTOR_R12_DM = float(f.readline())
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## FUNCTIONS CALLED BY MAIN ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
#functions to search for asohf correspondences
if ASOHF_FLAG:
@numba.njit(fastmath = True)
def match_finder(R1, x1, y1, z1, R2, x2, y2, z2, asohf_matches, thres = 1.):
dim2 = len(R2)
for j in range(dim2):
if (np.sqrt((x1-x2[j])**2 + (y1-y2[j])**2 + (z1-z2[j])**2) < thres*(R1 + R2[j]) + LL):
asohf_matches[j] = 1
return asohf_matches
@numba.njit(fastmath = True)
def DM_halo_finder(R1, x1, y1, z1, R2, x2, y2, z2,
asohf_dm_matches, matches_distance, thres = 1., fact_Rvir = 0.5):
dim2 = len(x2)
for j in range(dim2):
dist = np.sqrt((x1-x2[j])**2 + (y1-y2[j])**2 + (z1-z2[j])**2)
if dist < (thres*R1 + fact_Rvir*R2[j]):
asohf_dm_matches[j] = 1
matches_distance[j] = dist
return asohf_dm_matches, matches_distance
def good_groups(groups, minp):
return np.array([len(groups[ig]) >= minp for ig in range(len(groups))], dtype = bool)
def write_catalogue(haloes, iteration_data, PYHALMA_OUTPUT):
catalogue = open(PYHALMA_OUTPUT+'/stellar_haloes'+string_it, 'w')
catalogue.write('**************************************************************' + '\n')
catalogue.write(' '.join(str(x) for x in [*iteration_data.values()]) + '\n')
catalogue.write('**************************************************************' + '\n')
gap = ''
first_strings = ['Halo','n','Mass', 'xcm','ycm','zcm', 'xpeak','ypeak','zpeak', 'xbound','ybound','zbound','id','Mass','frac', 'm_rps','m_rps','m_SFR','m_in',
' R ','R_05','R_05','R_05','R_05','R_05', 'sigma','sigma','sig_x',
'sig_y','sig_z','J','Jx','Jy','Jz', 'V_x','V_y','V_z','Pro.','Pro.',
'n','type', 'age','age', 'Z','Z', 'V/Sigma', 'lambda', 'k_co', 'v_TF', 'a', 'b', 'c', 'sersic',
'lum_u','lum_g','lum_r','lum_i','sb_u','sb_g','sb_r','sb_i','ur_color','gr_color','Mass',
'ASOHF', 'Mass', 'R_vir', 'Mass']
second_strings = ['ID',' part', ' * ', 'kpc','kpc','kpc', 'kpc','kpc','kpc', 'kpc','kpc','kpc', 'bound', 'gas','g_cold', 'cold','hot',' * ',' * ', 'max','3D',
'1D','1D_x','1D_y','1D_z', '05_3D','05_1D','05_1D','05_1D','05_1D',
' ',' ',' ',' ', 'km/s','km/s','km/s',
'(1)','(2)','merg','merg','m_weig','mean', 'm_weig','mean', ' ', ' ', ' ',
'km/s', 'kpc', 'kpc', 'kpc', ' ',
' ',' ',' ',' ',' ',' ',' ',' ',' ',' ','BH',
'ID', 'ASOHF', ' ', 'DM']
first_line = ''
for element in first_strings:
first_line += f'{element:15s}'
second_line = ''
for element in second_strings:
second_line += f'{element:15s}'
catalogue.write('------------------------------------------------------------------------------------------------------------')
catalogue.write('------------------------------------------------------------------------------------------------------------')
catalogue.write('------------------------------------------------------------------------------------------------------------')
catalogue.write('------------------------------------------------------------------------------------------------------------')
catalogue.write('------------------------------------------------------------------------------------------------------------')
catalogue.write('------------------------------------------------------------------------------------------------------------')
catalogue.write('--------------------------')
catalogue.write('\n')
catalogue.write(' '+first_line)
catalogue.write('\n')
catalogue.write(' '+second_line)
catalogue.write('\n')
catalogue.write('------------------------------------------------------------------------------------------------------------')
catalogue.write('------------------------------------------------------------------------------------------------------------')
catalogue.write('------------------------------------------------------------------------------------------------------------')
catalogue.write('------------------------------------------------------------------------------------------------------------')
catalogue.write('------------------------------------------------------------------------------------------------------------')
catalogue.write('------------------------------------------------------------------------------------------------------------')
catalogue.write('-------------------------')
catalogue.write('\n')
nhal = iteration_data['nhal']
for ih in range(nhal):
ih_values = [*haloes[ih].values()]
catalogue_line = ''
for element in ih_values:
if type(element) == int or type(element) == np.int32 or type(element) == np.int64:
catalogue_line += f'{element:15d}'
elif type(element) == float or type(element) == np.float32 or type(element) == np.float64:
if np.log10(element) >= 5. or np.log10(element) <= -2:
catalogue_line += f'{element:15.2e}'
else:
catalogue_line += f'{element:15.2f}'
catalogue.write(catalogue_line)
catalogue.write('\n')
catalogue.close()
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
# M A I N ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
#########################################################
# SETTING UP masclet_framework and import scripts/modules
#########################################################
sys.path.append(PATH_MASCLET_FRAMEWORK)
from masclet_framework import read_masclet, units, read_asohf
from python_scripts import halo_properties, halo_gas, pycalipso, fof
#########################################################
##############################################
# SETTING UP FOLDERS
##############################################
ASOHF_OUTPUT = 'asohf_results'
SIMU_MASCLET = 'simu_masclet'
##############################################
##############################################
# CREATE PARAMETERS.JSON FILE FOR READ_MASCLET
##############################################
PARAMETERS_FILE = 'masclet_parameters.json'
parameters = {'NMAX': int(NX), 'NMAY': int(NY), 'NMAZ': int(NZ),
'NPALEV': int(NPALEV), 'NLEVELS': int(NLEVELS),
'NAMRX': int(NAMRX), 'NAMRY': int(NAMRY), 'NAMRZ': int(NAMRZ),
'SIZE': float(L)}
with open(os.path.join(SIMU_MASCLET, PARAMETERS_FILE), 'w') as json_file:
json.dump(parameters, json_file)
##############################################
##############################################
# IMPORT PYFOF IF CHOSEN
##############################################
if PYFOF_FLAG:
import pyfof
##############################################
##################################
# SETTING UP NUMBER OF CORES
##################################
if __name__ == '__main__':
numba.set_num_threads(NCORE)
##################################
print('*************************************')
print('*************************************')
print('----> Available CPU cores:', numba.config.NUMBA_NUM_THREADS)
print('----> Using', numba.get_num_threads(), 'cores')
##################################
#############################################################################################
# CREATE pyHALMA_output FOLDER INSIDE SIMU_MASCLET if it does not exist
#############################################################################################
if not os.path.exists(SIMU_MASCLET+'/pyHALMA_output'):
os.makedirs(SIMU_MASCLET+'/pyHALMA_output')
PYHALMA_OUTPUT = SIMU_MASCLET+'/pyHALMA_output'
#############################################################################################
##############################################
# SETTING UP BASIC COSMOLOGY
##############################################
#Scale factor at z = 0, which is 1 Mpc, that is, 1/10.98 u.l.
RE0 = 1.0/10.98
#Background density at z = 0 in Msun/Mpc^3
RODO = 3 * OMEGA0 * (ACHE*3.66e-3)**2 * units.mass_to_sun / units.length_to_mpc**3
##############################################
########## ########## ########## ########## ##########
### CALIPSO BASICS
########## ########## ########## ########## ##########
DIRSSP = 'E-MILES'
if CALIPSO_FLAG:
from astropy import cosmology
print('*************************************')
print('Before FoF, reading calipso files..')
WAVELENGHTS, SSP, AGE_SPAN, Z_SPAN, MH_SPAN, N_AGES, N_Z, N_W = pycalipso.readSSPfiles(DIRSSP, L_START, L_END)
N_F, N_LINES_FILTERS, W_FILTERS, RESPONSE_FILTERS = pycalipso.readFilters()
print('.. done')
dlum = 1.e-5 #10 pc --> absolute mag
magssp = np.zeros((N_AGES,N_Z,N_F))
fluxssp = np.zeros((N_AGES,N_Z,N_F))
LUMU = np.zeros((N_AGES,N_Z))
LUMG = np.zeros((N_AGES,N_Z))
LUMR = np.zeros((N_AGES,N_Z))
LUMI = np.zeros((N_AGES,N_Z))
# HERE WE CALCULATE ABSOLUTE MAGNITUDES (SINCE DLUM = 10 PC) AND THUS LUMINOSITIES
# Calculating luminosity in u,g,r filters of each SSP
for iage in range(N_AGES):
for iZ in range(N_Z):
pycalipso.mag_v1_0(WAVELENGHTS, SSP[iage, iZ, :], N_W, magssp[iage, iZ, :],
fluxssp[iage, iZ, :], N_F, N_LINES_FILTERS, W_FILTERS, RESPONSE_FILTERS, dlum, zeta=0.0)
LUMU[iage,iZ]=10.**(-0.4*(magssp[iage,iZ, 0]-USUN)) #luminosity in u filter (U band more or less)
LUMG[iage,iZ]=10.**(-0.4*(magssp[iage,iZ, 1]-GSUN)) #luminosity in g filter (B band more or less)
LUMR[iage,iZ]=10.**(-0.4*(magssp[iage,iZ, 2]-RSUN)) #luminosity in r filter (R band more or less)
LUMI[iage,iZ]=10.**(-0.4*(magssp[iage,iZ, 3]-ISUN)) #luminosity in i filter (I band more or less)
I_START = np.argmin(np.abs(L_START - WAVELENGHTS))
I_END = np.argmin(np.abs(L_END - WAVELENGHTS))
print('Spectrum will be written in the range:', WAVELENGHTS[I_START], WAVELENGHTS[I_END], '(A)')
DISP=WAVELENGHTS[2]-WAVELENGHTS[1] #resolución espectral
print('*************************************')
print()
print()
print()
#############################################################################################
# CREATE calipso OUTPUT FOLDERS INSIDE SIMU_MASCLET IF THEY DON'T EXIST
#############################################################################################
if not os.path.exists(SIMU_MASCLET+'/images_x'):
os.makedirs(SIMU_MASCLET+'/images_x')
if not os.path.exists(SIMU_MASCLET+'/images_y'):
os.makedirs(SIMU_MASCLET+'/images_y')
if not os.path.exists(SIMU_MASCLET+'/images_z'):
os.makedirs(SIMU_MASCLET+'/images_z')
if not os.path.exists(SIMU_MASCLET+'/spectra_x'):
os.makedirs(SIMU_MASCLET+'/spectra_x')
if not os.path.exists(SIMU_MASCLET+'/spectra_y'):
os.makedirs(SIMU_MASCLET+'/spectra_y')
if not os.path.exists(SIMU_MASCLET+'/spectra_z'):
os.makedirs(SIMU_MASCLET+'/spectra_z')
#############################################################################################
########## ########## ########## ########## ##########
########## ########## ########## ########## ##########
##############################################
### CALCULATING LINKING LENGHT IF CHOSEN
##############################################
# THIS IS MADE SUCH THAT IF MASS_PART = 10^6 Msun, then LL = 2.4 ckpc,
# which is the most suitable configuration found
if CALCULATE_LL_FLAG:
gamma = np.log10(MASS_PART)
power = (gamma - 7.5)/3.
LL = 7.5 * 10**power #in Kpc
#round to 2 decimals
LL = np.round(LL, 2)
#to Mpc
LL /= 1e3
del gamma, power
##############################################
##############################################
print('****************************************************')
print('******************** pyHALMA ***********************')
print('****************************************************')
print('!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!')
print('CHECK!!!! --> LINKING LENGHT (kpc)', LL*1e3)
print('!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!')
print()
#############################################################################################
#############################################################################################
#############################################################################################
#############################################################################################
# LOOP OVER MASCLET SNAPSHOTS
#############################################################################################
# Variables for accreted/in-situ mass calculation
part_insitu_before2 = np.array([]) #For every particle, if it was insitu, two iterations before
part_ih_before2 = np.array([]) #For every particle, the halo it belonged to two iterations before
part_oripas_before2 = np.array([]) #For every particle, the oripas iteration two iterations before
pro1_before2 = np.array([]) #For every halo, the main progenitor two iterations before
part_insitu_before = np.array([]) #For every particle, if it was insitu the previous iteration
part_ih_before = np.array([]) #For every particle, the halo it belonged to in the previous iteration
part_oripas_before = np.array([]) #For every particle, the oripas iteration before
####################################
oripas_before = np.array([]) #For every halo, the oripas iteration before
omm = np.array([]) #Masses of the haloes in the iteration before
cosmo_time_before = 0. #Time of the iteration before
for it_count, iteration in enumerate(range(FIRST, LAST+STEP, STEP)):
print()
print('!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!')
print(' Iteration', iteration)
print('!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!')
print()
print('Opening MASCLET files')
print()
#String to save outputs
string_it = f'{iteration:05d}'
#OPEN MASCLET files
#READ GRID
grid_data = read_masclet.read_grids(iteration, path=SIMU_MASCLET, parameters_path=SIMU_MASCLET, digits=5,
read_general=True, read_patchnum=True, read_dmpartnum=False,
read_patchcellextension=True, read_patchcellposition=True, read_patchposition=True,
read_patchparent=False)
cosmo_time = grid_data[1]
mass_dm_part = grid_data[3]*units.mass_to_sun
zeta = grid_data[4]
rete = RE0 / (1.0+zeta) # Scale factor at this iteration
rho_B = RODO / rete**3 # Background density of the universe at this iteration
dt = 0. #time step between iterations
if it_count>0:
dt = (cosmo_time - cosmo_time_before)*units.time_to_yr/1e9
else:
try:
grid_data_bef = read_masclet.read_grids(iteration-STEP, path=SIMU_MASCLET, parameters_path=SIMU_MASCLET, digits=5, read_general=True)
cosmo_time_before = grid_data_bef[1]
dt = (cosmo_time - cosmo_time_before)*units.time_to_yr/1e9
except:
pass
cosmo_time = grid_data[1]
print(f'Cosmo time (Gyr): {cosmo_time*units.time_to_yr/1e9:.2f}')
print(f'Redshift (z): {zeta:.2f}')
if FULL_DOMAIN_FLAG:
masclet_st_data = read_masclet.read_clst(iteration, path = SIMU_MASCLET, parameters_path=SIMU_MASCLET,
digits=5, max_refined_level=1000,
output_deltastar=False, verbose=False, output_position=True,
output_velocity=True, output_mass=True, output_time=True,
output_metalicity=True, output_id=True, are_BH = True,
output_BH=True)
else:
if PARALLEL_FOF:
print(' WARNING: PARALLEL FoF is not available for non FULL_DOMAIN')
print(' PARALLEL FoF will be deactivated')
print()
PARALLEL_FOF = False
read_region = ("box",X1,X2,Y1,Y2,Z1,Z2)
masclet_st_data = read_masclet.lowlevel_read_clst(iteration, path = SIMU_MASCLET, parameters_path=SIMU_MASCLET,
digits=5, max_refined_level=1000,
output_deltastar=False, verbose=False, output_position=True,
output_velocity=True, output_mass=True, output_time=True,
output_metalicity=True, output_id=True, are_BH = True,
output_BH=True, read_region=read_region)
st_x = masclet_st_data[0]
st_y = masclet_st_data[1]
st_z = masclet_st_data[2]
st_vx = masclet_st_data[3]*CLIGHT #in km/s
st_vy = masclet_st_data[4]*CLIGHT #in km/s
st_vz = masclet_st_data[5]*CLIGHT #in km/s
st_mass = masclet_st_data[6]*units.mass_to_sun #in Msun
st_age = masclet_st_data[7]*units.time_to_yr/1e9 #in Gyr
st_met = masclet_st_data[8]
st_oripa = masclet_st_data[9] #necessary for mergers
#insitu formed stars
st_insitu = np.zeros(len(st_x), dtype = np.int32)
bh_x = masclet_st_data[10]
bh_y = masclet_st_data[11]
bh_z = masclet_st_data[12]
bh_mass = masclet_st_data[16]*units.mass_to_sun #in Msun
if st_x.shape[0] > 0: #THERE ARE STAR PARTICLES
##############################################
# BUILD KD-TREE FOR STAR CALCULATIONS
##############################################
data = np.array((st_x, st_y, st_z)).T
data = data.astype(np.float64)
st_kdtree = KDTree(data)
##############################################
print()
print(f'Number of star particles: {len(st_x):.2e}', )
#APPLY FOF
print()
print('----------> FoF begins <--------')
print()
if PYFOF_FLAG:
print(' PYFOF ACTIVATED')
print()
if PARALLEL_FOF:
print(' PARALLEL FoF ACTIVATED')
print()
t0 = time.time()
groups = fof.pyfof_friends_of_friends_parallel(st_x, st_y, st_z,
LL, L, MINP, st_mass)
groups = np.array(groups, dtype=object)
print(' time in FoF:', time.time()-t0)
else:
t0 = time.time()
groups = fof.pyfof_friends_of_friends_serial(st_x, st_y, st_z, LL)
groups = np.array(groups, dtype=object)
print(' time in FoF:', time.time()-t0)
else:
print(' INTERNAL FoF ACTIVATED')
print()
if PARALLEL_FOF:
print(' PARALLEL FoF ACTIVATED')
print()
t0 = time.time()
groups = fof.friends_of_friends_parallel(st_x, st_y, st_z, LL, L, MINP, st_mass)
groups = np.array(groups, dtype=object)
print(' time in FoF:', time.time()-t0)
else:
t0 = time.time()
groups = fof.friends_of_friends_serial(st_x, st_y, st_z, LL)
groups = np.array(groups, dtype=object)
print(' time in FoF:', time.time()-t0)
#CLEAN THOSE HALOES with npart < minp
groups = groups[good_groups(groups, MINP)]
print()
print(len(groups),'haloes found')
print(np.sum([len(groups[ihal]) for ihal in range(len(groups))]), 'particles in haloes')
print()
if len(groups) > 0: #THERE ARE HALOES
#in order to check if dark matter data is needed
masclet_dm_data = None
#READ GAS IF RPS or ESCAPE VELOCITY CLEANING
if RPS_FLAG or ESCAPE_CLEANING:
print(' Opening grid, clus and DM files')
#READ CLUS
print(' Reading clus file')
if FULL_DOMAIN_FLAG:
gas_data = read_masclet.read_clus(iteration, path=SIMU_MASCLET, parameters_path=SIMU_MASCLET, digits=5, max_refined_level=1000, output_delta=True,
output_v=True, output_pres=False, output_pot=False, output_opot=False, output_temp=True, output_metalicity=False,
output_cr0amr=True, output_solapst=True, is_mascletB=False, output_B=False, is_cooling=True, verbose=False)
else:
read_region = ("box",X1,X2,Y1,Y2,Z1,Z2)
gas_data = read_masclet.lowlevel_read_clus(iteration, path=SIMU_MASCLET, parameters_path=SIMU_MASCLET, digits=5, max_refined_level=1000, output_delta=True,
output_v=True, output_pres=False, output_pot=False, output_opot=False, output_temp=True, output_metalicity=False,
output_cr0amr=True, output_solapst=True, is_mascletB=False, output_B=False, is_cooling=True, verbose=False,
read_region=read_region)
#READ DM
print(' Reading DM file')
if FULL_DOMAIN_FLAG:
masclet_dm_data = read_masclet.read_cldm(iteration, path = SIMU_MASCLET, parameters_path=SIMU_MASCLET,
digits=5, max_refined_level=1000, output_deltadm = False,
output_position=True, output_velocity=False, output_mass=True)
else:
read_region = ("box",X1,X2,Y1,Y2,Z1,Z2)
masclet_dm_data = read_masclet.lowlevel_read_cldm(iteration, path = SIMU_MASCLET, parameters_path=SIMU_MASCLET,
digits=5, max_refined_level=1000, output_deltadm = False,
output_position=True, output_velocity=False, output_mass=True, read_region=read_region)
print(' Done')
print()
#READ DM IF DM_FLAG and not (RPS_FLAG or ESCAPE_CLEANING)
if DM_FLAG and not (RPS_FLAG or ESCAPE_CLEANING):
#READ DM
print(' Reading DM file')
if FULL_DOMAIN_FLAG:
masclet_dm_data = read_masclet.read_cldm(iteration, path = SIMU_MASCLET, parameters_path=SIMU_MASCLET,
digits=5, max_refined_level=1000, output_deltadm = False,
output_position=True, output_velocity=False, output_mass=True)
else:
read_region = ("box",X1,X2,Y1,Y2,Z1,Z2)
masclet_dm_data = read_masclet.lowlevel_read_cldm(iteration, path = SIMU_MASCLET, parameters_path=SIMU_MASCLET,
digits=5, max_refined_level=1000, output_deltadm = False,
output_position=True, output_velocity=False, output_mass=True, read_region=read_region)
print(' Done')
print()
#READ ASOHF STELLAR AND DARK MATTER CATALOGUES
if ASOHF_FLAG:
print(' Reading ASOHF files')
if FULL_DOMAIN_FLAG:
asohf_dm_data = read_asohf.read_families(iteration, path=ASOHF_OUTPUT, output_format='arrays')
asohf_st_data = read_asohf.read_stellar_haloes(iteration, path=ASOHF_OUTPUT, output_format='arrays')
else:
read_region = ("box",X1,X2,Y1,Y2,Z1,Z2)
asohf_dm_data = read_asohf.read_families(iteration, path=ASOHF_OUTPUT, output_format='arrays', read_region=read_region)
asohf_st_data = read_asohf.read_stellar_haloes(iteration, path=ASOHF_OUTPUT, output_format='arrays', read_region=read_region)
asohf_st_num = len(asohf_st_data['id'])
asohf_dm_num = len(asohf_dm_data['id'])
print(' Done')
print()
if not (RPS_FLAG or ESCAPE_CLEANING): #if not opened before, open dm file
print(' Reading DM file')
if FULL_DOMAIN_FLAG:
masclet_dm_data = read_masclet.read_cldm(iteration, path = SIMU_MASCLET, parameters_path=SIMU_MASCLET,
digits=5, max_refined_level=1000, output_deltadm = False,
output_position=True, output_velocity=False, output_mass=True)
else:
read_region = ("box",X1,X2,Y1,Y2,Z1,Z2)
masclet_dm_data = read_masclet.lowlevel_read_cldm(iteration, path = SIMU_MASCLET, parameters_path=SIMU_MASCLET,
digits=5, max_refined_level=1000, output_deltadm = False,
output_position=True, output_velocity=False, output_mass=True, read_region=read_region)
print(' Done')
# BUILD DARK MATTER KDTREE FOR POT. ENERGY PURPOSES
if DM_FLAG or RPS_FLAG or POT_ENERGY_FLAG:
##############################################
print(' Building DM KDtree')
t0 = time.time()
data_dm = np.array((masclet_dm_data[0], masclet_dm_data[1], masclet_dm_data[2])).T
dm_kdtree = KDTree(data_dm)
print(' Time to build DM kdtree:', time.time()-t0, 's')
##############################################
print()
###########################################
###########################################
#### SUBSTRUCTURE IDENTIFICATION BLOCK ####
if SUBSTRUCTURE_FLAG:
print('---> Identifying substructure (splitting bridges between separated FoF groups) <---')
counter_sub = 0
groups_before = np.copy(groups)
groups = groups.tolist()
for ihal in tqdm(range(len(groups_before))):
part_list = np.array(groups_before[ihal])
starting_npart = len(part_list)
new_groups = []
sub_LL = LL/DEC_FACTOR
while sub_LL > MIN_LL:
#Center of mass
(cx, cy, cz, mass) = halo_properties.center_of_mass(part_list, st_x, st_y, st_z, st_mass)
#Half mass radius
r12 = halo_properties.half_mass_radius(cx, cy, cz, mass, part_list, st_x, st_y, st_z, st_mass)
#Sub FoF
if PYFOF_FLAG:
sub_groups = fof.pyfof_friends_of_friends_serial(st_x[part_list], st_y[part_list], st_z[part_list], sub_LL)
sub_groups = np.array(sub_groups, dtype=object)
sub_groups = sub_groups[good_groups(sub_groups, MINP_SUB)]
else:
sub_groups = fof.friends_of_friends_serial(st_x[part_list], st_y[part_list], st_z[part_list], sub_LL)
sub_groups = np.array(sub_groups, dtype=object)
sub_groups = sub_groups[good_groups(sub_groups, MINP_SUB)]
for group in sub_groups:
part_fraction = len(group)/len(part_list)
if 0.6 > part_fraction > 0.01:
#Center of mass sub
(cx_sub, cy_sub, cz_sub, mass_sub) = halo_properties.center_of_mass(
part_list[group], st_x, st_y, st_z, st_mass
)
#Distance between centers of mass
dist_sub = np.sqrt((cx-cx_sub)**2 + (cy-cy_sub)**2 + (cz-cz_sub)**2)
if dist_sub > DIST_FACTOR * r12:
new_groups.append(part_list[group].tolist())
counter_sub += 1
#eliminate particles in group from part_list
for new_group in new_groups:
part_list = np.setdiff1d(part_list, new_group)
sub_LL = sub_LL/DEC_FACTOR
#If more than 70% of the particles were eliminated, do not save the halo
#This can happen if a MAJOR MERGER is being connected by a bridge
if len(part_list)/starting_npart < 0.3:
groups[ihal] = []
print('WARNING! Empty halo due to substructure identification!')
else:
groups[ihal] = part_list.tolist()
#Append new groups:
if len(new_groups) > 0:
groups.extend(new_groups)
groups = np.array(groups, dtype=object)
#Remove empty groups
groups = groups[good_groups(groups, MINP)]
print('Substructures / Cut FoF bridges:', counter_sub)
print()
###########################################
###########################################
###########################################
#CALCULATE CM, CM_vel AND phase-space cleaning
print('---> Phase-space cleaning begins <---')
center_x = []
center_y = []
center_z = []
masses = []
rmax = []
num_particles = []
num_cells = []
new_groups = []
which_cell_list_x = []
which_cell_list_y = []
which_cell_list_z = []
for ihal in tqdm(range(len(groups))):
# Particle indices, center of mass and mass
part_list = np.array(groups[ihal], dtype=np.int32)
(cx, cy, cz, mass) = halo_properties.center_of_mass(part_list, st_x, st_y, st_z, st_mass)
(vx, vy, vz) = halo_properties.CM_velocity(mass, part_list, st_vx, st_vy, st_vz, st_mass)
most_distant_r = halo_properties.furthest_particle(cx, cy, cz, part_list, st_x, st_y, st_z)
#PHASE-SPACE CLEANING
#Create 3D grid with cellsize 2*LL for cleaning. Calculate most distant particle to center of mass
grid = np.arange( - ( most_distant_r + LL ), most_distant_r + LL, 2*LL)
n_cell = len(grid)
vcm_cell = np.zeros((n_cell, n_cell, n_cell, 3))
mass_cell = np.zeros((n_cell, n_cell, n_cell))
quantas_cell = np.zeros((n_cell, n_cell, n_cell))
sig3D_cell = np.zeros((n_cell, n_cell, n_cell))
#CALCULATE cell quantities
(vcm_cell, mass_cell, quantas_cell, sig3D_cell,
which_cell_x, which_cell_y, which_cell_z) = halo_properties.calc_cell(cx, cy, cz, grid, part_list,
st_x, st_y, st_z, st_vx, st_vy, st_vz,
st_mass, vcm_cell, mass_cell, quantas_cell, sig3D_cell)
#DO THE CLEANING: Right CM, Mass and Furthest particle
(cx, cy, cz, mass,
most_distant_r, control,
which_cell_x, which_cell_y, which_cell_z) = halo_properties.clean_cell(cx, cy, cz, mass, most_distant_r, grid, part_list, st_x, st_y,
st_z, st_vx, st_vy, st_vz, st_mass, vcm_cell,
mass_cell, quantas_cell, sig3D_cell, LL, Q_FIL, SIG_FIL,
which_cell_x, which_cell_y, which_cell_z)
if not PSC_FLAG: #IF NOT PHASE-SPACE CLEANING, ALL PARTICLES ARE GOOD
control = np.ones(len(part_list)).astype(bool)
if len(part_list) < 500: #if low number of particles, do not clean, as it may not converge
control = np.ones(len(part_list)).astype(bool)
#CLEANING DONE --> CALCULATE AGAIN CENTER OF MASS AND MASS
control = control.astype(bool)
npart = len(part_list[control])
cx, cy, cz, mass = halo_properties.center_of_mass(part_list[control], st_x, st_y, st_z, st_mass)
most_distant_r = halo_properties.furthest_particle(cx, cy, cz, part_list[control], st_x, st_y, st_z)
#If after cleaning there are more than MINP particles, save the halo
if npart >= MINP:
center_x.append(cx)
center_y.append(cy)
center_z.append(cz)
masses.append(mass)
rmax.append(most_distant_r)
num_particles.append(npart)
num_cells.append(n_cell)
new_groups.append(part_list[control])
which_cell_list_x.append(which_cell_x)
which_cell_list_y.append(which_cell_y)
which_cell_list_z.append(which_cell_z)
#From list to arrays for convenience
center_x = np.array(center_x)
center_y = np.array(center_y)
center_z = np.array(center_z)
masses = np.array(masses)
rmax = np.array(rmax)
num_particles = np.array(num_particles).astype(np.int32)
num_halos = len(new_groups)
n_part_in_halos = np.sum(num_particles)
num_cells = np.array(num_cells).astype(np.int32)
which_cell_list_x = np.array(which_cell_list_x, dtype=object)
which_cell_list_y = np.array(which_cell_list_y, dtype=object)
which_cell_list_z = np.array(which_cell_list_z, dtype=object)
print('Number of haloes after phase-space cleaning:', num_halos)
print('Number of particles in haloes after cleaning:', n_part_in_halos)
##########################################################################
##########################################################################
##########################################################################
#CALCULATE HALO PROPERTIES
print()
print('Calculating properties')
rad05 = np.zeros(num_halos)
rad05_x = np.zeros(num_halos)
rad05_y = np.zeros(num_halos)
rad05_z = np.zeros(num_halos)
velocities_x = np.zeros(num_halos)
velocities_y = np.zeros(num_halos)
velocities_z = np.zeros(num_halos)
specific_angular_momentum_x = np.zeros(num_halos)
specific_angular_momentum_y = np.zeros(num_halos)
specific_angular_momentum_z = np.zeros(num_halos)
specific_angular_momentum = np.zeros(num_halos)
density_peak_x = np.zeros(num_halos)
density_peak_y = np.zeros(num_halos)
density_peak_z = np.zeros(num_halos)
star_formation_masses = np.zeros(num_halos)
insitu_masses = np.zeros(num_halos)
sig_3D = np.zeros(num_halos)
sig_1D_x = np.zeros(num_halos)
sig_1D_y = np.zeros(num_halos)
sig_1D_z = np.zeros(num_halos)
stellar_ages = np.zeros(num_halos)
stellar_age_mass = np.zeros(num_halos)
metallicities = np.zeros(num_halos)
metallicities_mass = np.zeros(num_halos)
vsigma = np.zeros(num_halos)
kinematic_morphologies = np.zeros(num_halos)
lambda_ensellem = np.zeros(num_halos)
tully_fisher_velocities = np.zeros(num_halos)
s_axis_major = np.zeros(num_halos)
s_axis_intermediate = np.zeros(num_halos)
s_axis_minor = np.zeros(num_halos)
sersic_indices = np.zeros(num_halos)
gas_masses = np.zeros(num_halos)
cold_bound_gas_fractions = np.zeros(num_halos)
cold_unbound_gas_masses = np.zeros(num_halos)
hot_unbound_gas_masses = np.zeros(num_halos)
SMBH_masses = np.zeros(num_halos)
asohf_IDs = np.zeros(num_halos, dtype=np.int32)
asohf_mass = np.zeros(num_halos)
asohf_Rvir = np.zeros(num_halos)
darkmatter_mass = np.zeros(num_halos)
most_bound_x = np.zeros(num_halos)
most_bound_y = np.zeros(num_halos)
most_bound_z = np.zeros(num_halos)
most_bound_id = np.zeros(num_halos, dtype=np.int32)
for ihal in tqdm(range(num_halos)):
#HALO PARTICLE INDICES
part_list = new_groups[ihal]
#DENSITY PEAK
density_peak_x[ihal], density_peak_y[ihal], density_peak_z[ihal] = halo_properties.density_peak(part_list, st_x, st_y, st_z, st_mass, LL)
#ASSUME SOME CENTERING, CENTER OF MASS OR DENSITY PEAK
#CARE HERE, RMAX IS CALCULATED WITH RESPECT TO THE CENTER OF MASS, NOT THE DENSITY PEAK
cx = center_x[ihal]
cy = center_y[ihal]
cz = center_z[ihal]
mass = masses[ihal]
#EFFECTIVE RADIUS
rad05[ihal] = halo_properties.half_mass_radius(cx, cy, cz, mass, part_list, st_x, st_y, st_z, st_mass)
(rad05_x[ihal],
rad05_y[ihal],
rad05_z[ihal]) = halo_properties.half_mass_radius_proj(cx, cy, cz, mass, part_list, st_x, st_y, st_z, st_mass)
#BULK VELOCITY
(velocities_x[ihal],
velocities_y[ihal],
velocities_z[ihal]) = halo_properties.CM_velocity(mass, part_list, st_vx, st_vy, st_vz, st_mass)
#ANGULAR MOMENTUM
(specific_angular_momentum_x[ihal],
specific_angular_momentum_y[ihal],
specific_angular_momentum_z[ihal]) = halo_properties.angular_momentum(mass, part_list, st_x, st_y, st_z, st_vx, st_vy, st_vz, st_mass,
cx, cy, cz, velocities_x[ihal], velocities_y[ihal], velocities_z[ihal])
specific_angular_momentum[ihal] = ( specific_angular_momentum_x[ihal]**2 +
specific_angular_momentum_y[ihal]**2 +
specific_angular_momentum_z[ihal]**2 )**0.5
#TULLY-FISHER
tully_fisher_velocities[ihal] = halo_properties.tully_fisher_velocity(part_list, cx, cy, cz, st_x, st_y, st_z,
velocities_x[ihal], velocities_y[ihal], velocities_z[ihal],
specific_angular_momentum_x[ihal],
specific_angular_momentum_y[ihal],
specific_angular_momentum_z[ihal],
st_vx, st_vy, st_vz, st_mass, rad05[ihal])
#KINEMATIC MORPHOLOGY (Correa et al. 2017 [EAGLE])
kinematic_morphologies[ihal] = halo_properties.kinematic_morphology(part_list, st_x, st_y, st_z, st_vx, st_vy, st_vz, st_mass,
cx, cy, cz, velocities_x[ihal], velocities_y[ihal], velocities_z[ihal],
specific_angular_momentum_x[ihal], specific_angular_momentum_y[ihal],
specific_angular_momentum_z[ihal])
#SHAPE TENSOR EIGENVALUES (ELLIPSOID semi-axes)
(s_axis_major[ihal],
s_axis_intermediate[ihal],
s_axis_minor[ihal]) = halo_properties.halo_shape_fortran(part_list, st_x, st_y, st_z, st_mass, cx, cy, cz, rad05[ihal])
#SÉRSIC INDEX
sersic_indices[ihal] = halo_properties.simple_sersic_index(
part_list, st_x, st_y, st_z, center_x[ihal],
center_y[ihal], center_z[ihal],
rad05[ihal], LL, num_particles[ihal],
True, specific_angular_momentum_x[ihal],
specific_angular_momentum_y[ihal],
specific_angular_momentum_z[ihal],
)
star_formation_masses[ihal] = halo_properties.star_formation(part_list, st_mass, st_age, cosmo_time*units.time_to_yr/1e9, dt)
#SIGMA
sig_3D[ihal] = halo_properties.sigma_effective(part_list, rad05[ihal], st_x, st_y, st_z, st_vx, st_vy, st_vz, cx, cy, cz,
velocities_x[ihal], velocities_y[ihal], velocities_z[ihal])
# FAST / SLOW ROTATOR
grid = np.arange(-(rmax[ihal]+LL), rmax[ihal]+LL, 2*LL) #centers of the cells
n_cell = len(grid)
(sig_1D_x[ihal],
sig_1D_y[ihal],
sig_1D_z[ihal],
vsigma[ihal],
lambda_ensellem[ihal]) = halo_properties.sigma_projections_fortran(grid, n_cell, part_list, st_x, st_y, st_z,
st_vx, st_vy, st_vz,
velocities_x[ihal], velocities_y[ihal], velocities_z[ihal],
st_mass, cx, cy, cz,
rad05_x[ihal], rad05_y[ihal], rad05_z[ihal],
LL)
# METALLICITY AND AGE
(stellar_ages[ihal],
stellar_age_mass[ihal],
metallicities[ihal],
metallicities_mass[ihal]) = halo_properties.avg_age_metallicity(part_list, st_age, st_met, st_mass, cosmo_time*units.time_to_yr/1e9)
# GAS, DM AND ST PARTICLES INSIDE RADIUS FOR POTENTIAL ENERGY
POT_RADIUS = FACTOR_R12_POT*rad05[ihal]
if RPS_FLAG or POT_ENERGY_FLAG:
(gas_x_in, gas_y_in, gas_z_in,
gas_vx_in, gas_vy_in, gas_vz_in,
gas_mass_in, gas_temp_in,
dm_x_in, dm_y_in, dm_z_in, dm_mass_in,
st_x_in, st_y_in, st_z_in, st_mass_in,
st_oripa_in) = halo_gas.st_gas_dm_particles_inside(
rete, L, NX, grid_data,
gas_data, masclet_dm_data, masclet_st_data,
st_kdtree, dm_kdtree,
cx, cy, cz, POT_RADIUS, rho_B
)
# RPS EFFECTS
if RPS_FLAG:
(
gas_masses[ihal],
cold_bound_gas_fractions[ihal],
cold_unbound_gas_masses[ihal],
hot_unbound_gas_masses[ihal] ) = halo_gas.RPS(
gas_x_in, gas_y_in, gas_z_in,
gas_vx_in, gas_vy_in, gas_vz_in,
gas_mass_in, gas_temp_in,
dm_x_in, dm_y_in, dm_z_in, dm_mass_in,
st_x_in, st_y_in, st_z_in, st_mass_in,
velocities_x[ihal], velocities_y[ihal],
velocities_z[ihal], BRUTE_FORCE_LIM,
mass_dm_part