/
photon_simulator.py
1567 lines (1346 loc) · 63.4 KB
/
photon_simulator.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
"""
Classes for generating lists of photons and detected events
The SciPy Proceeding that describes this module in detail may be found at:
http://conference.scipy.org/proceedings/scipy2014/zuhone.html
The algorithms used here are based off of the method used by the
PHOX code (http://www.mpa-garching.mpg.de/~kdolag/Phox/),
developed by Veronica Biffi and Klaus Dolag. References for
PHOX may be found at:
Biffi, V., Dolag, K., Bohringer, H., & Lemson, G. 2012, MNRAS, 420, 3545
http://adsabs.harvard.edu/abs/2012MNRAS.420.3545B
Biffi, V., Dolag, K., Bohringer, H. 2013, MNRAS, 428, 1395
http://adsabs.harvard.edu/abs/2013MNRAS.428.1395B
"""
#-----------------------------------------------------------------------------
# Copyright (c) 2013, yt Development Team.
#
# Distributed under the terms of the Modified BSD License.
#
# The full license is in the file COPYING.txt, distributed with this software.
#-----------------------------------------------------------------------------
from yt.extern.six import string_types
from collections import defaultdict
import numpy as np
from yt.funcs import mylog, get_pbar, iterable, ensure_list
from yt.utilities.physical_constants import clight
from yt.utilities.cosmology import Cosmology
from yt.utilities.orientation import Orientation
from yt.visualization.fits_image import assert_same_wcs
from yt.utilities.parallel_tools.parallel_analysis_interface import \
communication_system, parallel_root_only, get_mpi_type, \
parallel_capable
from yt.units.yt_array import YTQuantity, YTArray, uconcatenate
from yt.utilities.on_demand_imports import _h5py as h5py
from yt.utilities.on_demand_imports import _astropy
import warnings
import os
comm = communication_system.communicators[-1]
axes_lookup = {"x":("y","z"),
"y":("z","x"),
"z":("x","y")}
def force_unicode(value):
if hasattr(value, 'decode'):
return value.decode('utf8')
else:
return value
def parse_value(value, default_units):
if isinstance(value, YTQuantity):
return value.in_units(default_units)
elif iterable(value):
return YTQuantity(value[0], value[1]).in_units(default_units)
else:
return YTQuantity(value, default_units)
def validate_parameters(first, second, skip=[]):
keys1 = list(first.keys())
keys2 = list(first.keys())
keys1.sort()
keys2.sort()
if keys1 != keys2:
raise RuntimeError("The two inputs do not have the same parameters!")
for k1, k2 in zip(keys1, keys2):
if k1 not in skip:
v1 = first[k1]
v2 = second[k2]
if isinstance(v1, string_types) or isinstance(v2, string_types):
check_equal = v1 == v2
else:
check_equal = np.allclose(v1, v2, rtol=0.0, atol=1.0e-10)
if not check_equal:
raise RuntimeError("The values for the parameter '%s' in the two inputs" % k1 +
" are not identical (%s vs. %s)!" % (v1, v2))
class PhotonList(object):
def __init__(self, photons, parameters, cosmo, p_bins):
self.photons = photons
self.parameters = parameters
self.cosmo = cosmo
self.p_bins = p_bins
self.num_cells = len(photons["x"])
def keys(self):
return self.photons.keys()
def items(self):
ret = []
for k, v in self.photons.items():
if k == "Energy":
ret.append((k, self[k]))
else:
ret.append((k,v))
return ret
def values(self):
ret = []
for k, v in self.photons.items():
if k == "Energy":
ret.append(self[k])
else:
ret.append(v)
return ret
def __getitem__(self, key):
if key == "Energy":
return [self.photons["Energy"][self.p_bins[i]:self.p_bins[i+1]]
for i in range(self.num_cells)]
else:
return self.photons[key]
def __contains__(self, key):
return key in self.photons
def __repr__(self):
return self.photons.__repr__()
@classmethod
def from_file(cls, filename):
r"""
Initialize a PhotonList from the HDF5 file *filename*.
"""
photons = {}
parameters = {}
f = h5py.File(filename, "r")
p = f["/parameters"]
parameters["FiducialExposureTime"] = YTQuantity(p["fid_exp_time"].value, "s")
parameters["FiducialArea"] = YTQuantity(p["fid_area"].value, "cm**2")
parameters["FiducialRedshift"] = p["fid_redshift"].value
parameters["FiducialAngularDiameterDistance"] = YTQuantity(p["fid_d_a"].value, "Mpc")
parameters["Dimension"] = p["dimension"].value
parameters["Width"] = YTQuantity(p["width"].value, "kpc")
parameters["HubbleConstant"] = p["hubble"].value
parameters["OmegaMatter"] = p["omega_matter"].value
parameters["OmegaLambda"] = p["omega_lambda"].value
d = f["/data"]
num_cells = d["x"][:].shape[0]
start_c = comm.rank*num_cells//comm.size
end_c = (comm.rank+1)*num_cells//comm.size
photons["x"] = YTArray(d["x"][start_c:end_c], "kpc")
photons["y"] = YTArray(d["y"][start_c:end_c], "kpc")
photons["z"] = YTArray(d["z"][start_c:end_c], "kpc")
photons["dx"] = YTArray(d["dx"][start_c:end_c], "kpc")
photons["vx"] = YTArray(d["vx"][start_c:end_c], "km/s")
photons["vy"] = YTArray(d["vy"][start_c:end_c], "km/s")
photons["vz"] = YTArray(d["vz"][start_c:end_c], "km/s")
n_ph = d["num_photons"][:]
if comm.rank == 0:
start_e = np.uint64(0)
else:
start_e = n_ph[:start_c].sum()
end_e = start_e + np.uint64(n_ph[start_c:end_c].sum())
photons["NumberOfPhotons"] = n_ph[start_c:end_c]
p_bins = np.cumsum(photons["NumberOfPhotons"])
p_bins = np.insert(p_bins, 0, [np.uint64(0)])
photons["Energy"] = YTArray(d["energy"][start_e:end_e], "keV")
f.close()
cosmo = Cosmology(hubble_constant=parameters["HubbleConstant"],
omega_matter=parameters["OmegaMatter"],
omega_lambda=parameters["OmegaLambda"])
return cls(photons, parameters, cosmo, p_bins)
@classmethod
def from_scratch(cls, data_source, redshift, area,
exp_time, photon_model, parameters=None,
center=None, dist=None, cosmology=None):
r"""
Initialize a PhotonList from a photon model. The redshift, collecting area,
exposure time, and cosmology are stored in the *parameters* dictionary which
is passed to the *photon_model* function.
Parameters
----------
data_source : `yt.data_objects.data_containers.YTSelectionContainer`
The data source from which the photons will be generated.
redshift : float
The cosmological redshift for the photons.
area : float
The collecting area to determine the number of photons in cm^2.
exp_time : float
The exposure time to determine the number of photons in seconds.
photon_model : function
A function that takes the *data_source* and the *parameters*
dictionary and returns a *photons* dictionary. Must be of the
form: photon_model(data_source, parameters)
parameters : dict, optional
A dictionary of parameters to be passed to the user function.
center : string or array_like, optional
The origin of the photons. Accepts "c", "max", or a coordinate.
dist : tuple, optional
The angular diameter distance in the form (value, unit), used
mainly for nearby sources. This may be optionally supplied
instead of it being determined from the *redshift* and given *cosmology*.
cosmology : `yt.utilities.cosmology.Cosmology`, optional
Cosmological information. If not supplied, we try to get
the cosmology from the dataset. Otherwise, \LambdaCDM with
the default yt parameters is assumed.
Examples
--------
This is the simplest possible example, where we call the built-in thermal model:
>>> thermal_model = ThermalPhotonModel(apec_model, Zmet=0.3)
>>> redshift = 0.05
>>> area = 6000.0
>>> time = 2.0e5
>>> sp = ds.sphere("c", (500., "kpc"))
>>> my_photons = PhotonList.from_user_model(sp, redshift, area,
... time, thermal_model)
If you wish to make your own photon model function, it must take as its
arguments the *data_source* and the *parameters* dictionary. However you
determine them, the *photons* dict needs to have the following items, corresponding
to cells which have photons:
"x" : the x-position of the cell relative to the source center in kpc, YTArray
"y" : the y-position of the cell relative to the source center in kpc, YTArray
"z" : the z-position of the cell relative to the source center in kpc, YTArray
"vx" : the x-velocity of the cell in km/s, YTArray
"vy" : the y-velocity of the cell in km/s, YTArray
"vz" : the z-velocity of the cell in km/s, YTArray
"dx" : the width of the cell in kpc, YTArray
"NumberOfPhotons" : the number of photons in the cell, NumPy array of unsigned 64-bit integers
"Energy" : the source rest-frame energies of the photons, YTArray
The last array is not the same size as the others because it contains the energies in all of
the cells in a single 1-D array. The first photons["NumberOfPhotons"][0] elements are
for the first cell, the next photons["NumberOfPhotons"][1] are for the second cell, and so on.
The following is a simple example where a point source with a single line emission
spectrum of photons is created. More complicated examples which actually
create photons based on the fields in the dataset could be created.
>>> import numpy as np
>>> import yt
>>> from yt.analysis_modules.photon_simulator.api import PhotonList
>>> def line_func(source, parameters):
...
... ds = source.ds
...
... num_photons = parameters["num_photons"]
... E0 = parameters["line_energy"] # Energies are in keV
... sigE = parameters["line_sigma"]
... src_ctr = parameters["center"]
...
... energies = norm.rvs(loc=E0, scale=sigE, size=num_photons)
...
... # Place everything in the center cell
... for i, ax in enumerate("xyz"):
... photons[ax] = (ds.domain_center[0]-src_ctr[0]).in_units("kpc")
... photons["vx"] = ds.arr([0], "km/s")
... photons["vy"] = ds.arr([0], "km/s")
... photons["vz"] = ds.arr([100.0], "km/s")
... photons["dx"] = ds.find_field_values_at_point("dx", ds.domain_center).in_units("kpc")
... photons["NumberOfPhotons"] = np.array(num_photons*np.ones(1), dtype="uint64")
... photons["Energy"] = ds.arr(energies, "keV")
>>>
>>> redshift = 0.05
>>> area = 6000.0
>>> time = 2.0e5
>>> parameters = {"num_photons" : 10000, "line_energy" : 5.0,
... "line_sigma" : 0.1}
>>> ddims = (128,128,128)
>>> random_data = {"density":(np.random.random(ddims),"g/cm**3")}
>>> ds = yt.load_uniform_grid(random_data, ddims)
>>> dd = ds.all_data()
>>> my_photons = PhotonList.from_user_model(dd, redshift, area,
... time, line_func,
... parameters=parameters)
"""
ds = data_source.ds
if parameters is None:
parameters = {}
if cosmology is None:
hubble = getattr(ds, "hubble_constant", None)
omega_m = getattr(ds, "omega_matter", None)
omega_l = getattr(ds, "omega_lambda", None)
if hubble == 0: hubble = None
if hubble is not None and \
omega_m is not None and \
omega_l is not None:
cosmo = Cosmology(hubble_constant=hubble,
omega_matter=omega_m,
omega_lambda=omega_l)
else:
cosmo = Cosmology()
else:
cosmo = cosmology
mylog.info("Cosmology: h = %g, omega_matter = %g, omega_lambda = %g" %
(cosmo.hubble_constant, cosmo.omega_matter, cosmo.omega_lambda))
if dist is None:
D_A = cosmo.angular_diameter_distance(0.0,redshift).in_units("Mpc")
else:
D_A = parse_value(dist, "Mpc")
redshift = 0.0
if center in ("center", "c"):
parameters["center"] = ds.domain_center
elif center in ("max", "m"):
parameters["center"] = ds.find_max("density")[-1]
elif iterable(center):
if isinstance(center, YTArray):
parameters["center"] = center.in_units("code_length")
elif isinstance(center, tuple):
if center[0] == "min":
parameters["center"] = ds.find_min(center[1])[-1]
elif center[0] == "max":
parameters["center"] = ds.find_max(center[1])[-1]
else:
raise RuntimeError
else:
parameters["center"] = ds.arr(center, "code_length")
elif center is None:
parameters["center"] = data_source.get_field_parameter("center")
parameters["FiducialExposureTime"] = parse_value(exp_time, "s")
parameters["FiducialArea"] = parse_value(area, "cm**2")
parameters["FiducialRedshift"] = redshift
parameters["FiducialAngularDiameterDistance"] = D_A
parameters["HubbleConstant"] = cosmo.hubble_constant
parameters["OmegaMatter"] = cosmo.omega_matter
parameters["OmegaLambda"] = cosmo.omega_lambda
dimension = 0
width = 0.0
for i, ax in enumerate("xyz"):
le, re = data_source.quantities.extrema(ax)
delta_min, delta_max = data_source.quantities.extrema("d%s"%ax)
le -= 0.5*delta_max
re += 0.5*delta_max
width = max(width, re-parameters["center"][i], parameters["center"][i]-le)
dimension = max(dimension, int(width/delta_min))
parameters["Dimension"] = 2*dimension
parameters["Width"] = 2.*width.in_units("kpc")
photons = photon_model(data_source, parameters)
mylog.info("Finished generating photons.")
p_bins = np.cumsum(photons["NumberOfPhotons"])
p_bins = np.insert(p_bins, 0, [np.uint64(0)])
return cls(photons, parameters, cosmo, p_bins)
def write_h5_file(self, photonfile):
"""
Write the photons to the HDF5 file *photonfile*.
"""
if parallel_capable:
mpi_long = get_mpi_type("int64")
mpi_double = get_mpi_type("float64")
local_num_cells = len(self.photons["x"])
sizes_c = comm.comm.gather(local_num_cells, root=0)
local_num_photons = np.sum(self.photons["NumberOfPhotons"])
sizes_p = comm.comm.gather(local_num_photons, root=0)
if comm.rank == 0:
num_cells = sum(sizes_c)
num_photons = sum(sizes_p)
disps_c = [sum(sizes_c[:i]) for i in range(len(sizes_c))]
disps_p = [sum(sizes_p[:i]) for i in range(len(sizes_p))]
x = np.zeros(num_cells)
y = np.zeros(num_cells)
z = np.zeros(num_cells)
vx = np.zeros(num_cells)
vy = np.zeros(num_cells)
vz = np.zeros(num_cells)
dx = np.zeros(num_cells)
n_ph = np.zeros(num_cells, dtype="uint64")
e = np.zeros(num_photons)
else:
sizes_c = []
sizes_p = []
disps_c = []
disps_p = []
x = np.empty([])
y = np.empty([])
z = np.empty([])
vx = np.empty([])
vy = np.empty([])
vz = np.empty([])
dx = np.empty([])
n_ph = np.empty([])
e = np.empty([])
comm.comm.Gatherv([self.photons["x"].d, local_num_cells, mpi_double],
[x, (sizes_c, disps_c), mpi_double], root=0)
comm.comm.Gatherv([self.photons["y"].d, local_num_cells, mpi_double],
[y, (sizes_c, disps_c), mpi_double], root=0)
comm.comm.Gatherv([self.photons["z"].d, local_num_cells, mpi_double],
[z, (sizes_c, disps_c), mpi_double], root=0)
comm.comm.Gatherv([self.photons["vx"].d, local_num_cells, mpi_double],
[vx, (sizes_c, disps_c), mpi_double], root=0)
comm.comm.Gatherv([self.photons["vy"].d, local_num_cells, mpi_double],
[vy, (sizes_c, disps_c), mpi_double], root=0)
comm.comm.Gatherv([self.photons["vz"].d, local_num_cells, mpi_double],
[vz, (sizes_c, disps_c), mpi_double], root=0)
comm.comm.Gatherv([self.photons["dx"].d, local_num_cells, mpi_double],
[dx, (sizes_c, disps_c), mpi_double], root=0)
comm.comm.Gatherv([self.photons["NumberOfPhotons"], local_num_cells, mpi_long],
[n_ph, (sizes_c, disps_c), mpi_long], root=0)
comm.comm.Gatherv([self.photons["Energy"].d, local_num_photons, mpi_double],
[e, (sizes_p, disps_p), mpi_double], root=0)
else:
x = self.photons["x"].d
y = self.photons["y"].d
z = self.photons["z"].d
vx = self.photons["vx"].d
vy = self.photons["vy"].d
vz = self.photons["vz"].d
dx = self.photons["dx"].d
n_ph = self.photons["NumberOfPhotons"]
e = self.photons["Energy"].d
if comm.rank == 0:
f = h5py.File(photonfile, "w")
# Parameters
p = f.create_group("parameters")
p.create_dataset("fid_area", data=float(self.parameters["FiducialArea"]))
p.create_dataset("fid_exp_time", data=float(self.parameters["FiducialExposureTime"]))
p.create_dataset("fid_redshift", data=self.parameters["FiducialRedshift"])
p.create_dataset("hubble", data=self.parameters["HubbleConstant"])
p.create_dataset("omega_matter", data=self.parameters["OmegaMatter"])
p.create_dataset("omega_lambda", data=self.parameters["OmegaLambda"])
p.create_dataset("fid_d_a", data=float(self.parameters["FiducialAngularDiameterDistance"]))
p.create_dataset("dimension", data=self.parameters["Dimension"])
p.create_dataset("width", data=float(self.parameters["Width"]))
# Data
d = f.create_group("data")
d.create_dataset("x", data=x)
d.create_dataset("y", data=y)
d.create_dataset("z", data=z)
d.create_dataset("vx", data=vx)
d.create_dataset("vy", data=vy)
d.create_dataset("vz", data=vz)
d.create_dataset("dx", data=dx)
d.create_dataset("num_photons", data=n_ph)
d.create_dataset("energy", data=e)
f.close()
comm.barrier()
def project_photons(self, normal, area_new=None, exp_time_new=None,
redshift_new=None, dist_new=None,
absorb_model=None, psf_sigma=None,
sky_center=None, responses=None,
convolve_energies=False, no_shifting=False,
north_vector=None, prng=np.random):
r"""
Projects photons onto an image plane given a line of sight.
Parameters
----------
normal : character or array_like
Normal vector to the plane of projection. If "x", "y", or "z", will
assume to be along that axis (and will probably be faster). Otherwise,
should be an off-axis normal vector, e.g [1.0,2.0,-3.0]
area_new : float, optional
New value for the effective area of the detector. If *responses*
are specified the value of this keyword is ignored.
exp_time_new : float, optional
The new value for the exposure time.
redshift_new : float, optional
The new value for the cosmological redshift.
dist_new : tuple, optional
The new value for the angular diameter distance in the form
(value, unit), used mainly for nearby sources. This may be optionally supplied
instead of it being determined from the cosmology.
absorb_model : 'yt.analysis_modules.photon_simulator.PhotonModel`, optional
A model for galactic absorption.
psf_sigma : float, optional
Quick-and-dirty psf simulation using Gaussian smoothing with
standard deviation *psf_sigma* in degrees.
sky_center : array_like, optional
Center RA, Dec of the events in degrees.
responses : list of strings, optional
The names of the ARF and/or RMF files to convolve the photons with.
convolve_energies : boolean, optional
If this is set, the photon energies will be convolved with the RMF.
no_shifting : boolean, optional
If set, the photon energies will not be Doppler shifted.
north_vector : a sequence of floats
A vector defining the 'up' direction. This option sets the orientation of
the plane of projection. If not set, an arbitrary grid-aligned north_vector
is chosen. Ignored in the case where a particular axis (e.g., "x", "y", or
"z") is explicitly specified.
prng : NumPy `RandomState` object or numpy.random
A pseudo-random number generator. Typically will only be specified if you
have a reason to generate the same set of random numbers, such as for a
test. Default is the numpy.random module.
Examples
--------
>>> L = np.array([0.1,-0.2,0.3])
>>> events = my_photons.project_photons(L, area_new="sim_arf.fits",
... redshift_new=0.05,
... psf_sigma=0.01)
"""
if redshift_new is not None and dist_new is not None:
mylog.error("You may specify a new redshift or distance, "+
"but not both!")
if sky_center is None:
sky_center = YTArray([30.,45.], "degree")
else:
sky_center = YTArray(sky_center, "degree")
dx = self.photons["dx"].d
nx = self.parameters["Dimension"]
if psf_sigma is not None:
psf_sigma = parse_value(psf_sigma, "degree")
if not isinstance(normal, string_types):
L = np.array(normal)
orient = Orientation(L, north_vector=north_vector)
x_hat = orient.unit_vectors[0]
y_hat = orient.unit_vectors[1]
z_hat = orient.unit_vectors[2]
n_ph = self.photons["NumberOfPhotons"]
n_ph_tot = n_ph.sum()
eff_area = None
parameters = {}
if responses is not None:
responses = ensure_list(responses)
parameters["ARF"] = responses[0]
if len(responses) == 2:
parameters["RMF"] = responses[1]
area_new = parameters["ARF"]
zobs0 = self.parameters["FiducialRedshift"]
D_A0 = self.parameters["FiducialAngularDiameterDistance"]
scale_factor = 1.0
# If we use an RMF, figure out where the response matrix actually is.
if "RMF" in parameters:
rmf = _astropy.pyfits.open(parameters["RMF"])
if "MATRIX" in rmf:
mat_key = "MATRIX"
elif "SPECRESP MATRIX" in rmf:
mat_key = "SPECRESP MATRIX"
else:
raise RuntimeError("Cannot find the response matrix in the RMF "
"file %s! " % parameters["RMF"]+"It should "
"be named \"MATRIX\" or \"SPECRESP MATRIX\".")
rmf.close()
else:
mat_key = None
if (exp_time_new is None and area_new is None and
redshift_new is None and dist_new is None):
my_n_obs = n_ph_tot
zobs = zobs0
D_A = D_A0
else:
if exp_time_new is None:
Tratio = 1.
else:
Tratio = parse_value(exp_time_new, "s")/self.parameters["FiducialExposureTime"]
if area_new is None:
Aratio = 1.
elif isinstance(area_new, string_types):
if comm.rank == 0:
mylog.info("Using energy-dependent effective area: %s" % (parameters["ARF"]))
f = _astropy.pyfits.open(area_new)
earf = 0.5*(f["SPECRESP"].data.field("ENERG_LO")+f["SPECRESP"].data.field("ENERG_HI"))
eff_area = np.nan_to_num(f["SPECRESP"].data.field("SPECRESP"))
if "RMF" in parameters:
weights = self._normalize_arf(parameters["RMF"], mat_key)
eff_area *= weights
else:
mylog.warning("You specified an ARF but not an RMF. This is ok if the "+
"responses are normalized properly. If not, you may "+
"get inconsistent results.")
f.close()
Aratio = eff_area.max()/self.parameters["FiducialArea"].v
else:
mylog.info("Using constant effective area.")
Aratio = parse_value(area_new, "cm**2")/self.parameters["FiducialArea"]
if redshift_new is None and dist_new is None:
Dratio = 1.
zobs = zobs0
D_A = D_A0
else:
if redshift_new is None:
zobs = 0.0
D_A = parse_value(dist_new, "Mpc")
else:
zobs = redshift_new
D_A = self.cosmo.angular_diameter_distance(0.0,zobs).in_units("Mpc")
scale_factor = (1.+zobs0)/(1.+zobs)
Dratio = D_A0*D_A0*(1.+zobs0)**3 / \
(D_A*D_A*(1.+zobs)**3)
fak = Aratio*Tratio*Dratio
if fak > 1:
raise ValueError("This combination of requested parameters results in "
"%g%% more photons collected than are " % (100.*(fak-1.)) +
"available in the sample. Please reduce the collecting "
"area, exposure time, or increase the distance/redshift "
"of the object. Alternatively, generate a larger sample "
"of photons.")
my_n_obs = np.uint64(n_ph_tot*fak)
n_obs_all = comm.mpi_allreduce(my_n_obs)
if comm.rank == 0:
mylog.info("Total number of photons to use: %d" % (n_obs_all))
if my_n_obs == n_ph_tot:
idxs = np.arange(my_n_obs,dtype='uint64')
else:
idxs = prng.permutation(n_ph_tot)[:my_n_obs].astype("uint64")
obs_cells = np.searchsorted(self.p_bins, idxs, side='right')-1
delta = dx[obs_cells]
if isinstance(normal, string_types):
xsky = prng.uniform(low=-0.5,high=0.5,size=my_n_obs)
ysky = prng.uniform(low=-0.5,high=0.5,size=my_n_obs)
xsky *= delta
ysky *= delta
xsky += self.photons[axes_lookup[normal][0]][obs_cells]
ysky += self.photons[axes_lookup[normal][1]][obs_cells]
if not no_shifting:
vz = self.photons["v%s" % normal]
else:
x = prng.uniform(low=-0.5,high=0.5,size=my_n_obs)
y = prng.uniform(low=-0.5,high=0.5,size=my_n_obs)
z = prng.uniform(low=-0.5,high=0.5,size=my_n_obs)
if not no_shifting:
vz = self.photons["vx"]*z_hat[0] + \
self.photons["vy"]*z_hat[1] + \
self.photons["vz"]*z_hat[2]
x *= delta
y *= delta
z *= delta
x += self.photons["x"][obs_cells].d
y += self.photons["y"][obs_cells].d
z += self.photons["z"][obs_cells].d
xsky = x*x_hat[0] + y*x_hat[1] + z*x_hat[2]
ysky = x*y_hat[0] + y*y_hat[1] + z*y_hat[2]
if no_shifting:
eobs = self.photons["Energy"][idxs]
else:
shift = -vz.in_cgs()/clight
shift = np.sqrt((1.-shift)/(1.+shift))
eobs = self.photons["Energy"][idxs]*shift[obs_cells]
eobs *= scale_factor
if absorb_model is None:
not_abs = np.ones(eobs.shape, dtype='bool')
else:
mylog.info("Absorbing.")
absorb_model.prepare_spectrum()
emid = absorb_model.emid
aspec = absorb_model.get_spectrum()
absorb = np.interp(eobs, emid, aspec, left=0.0, right=0.0)
randvec = aspec.max()*prng.uniform(size=eobs.shape)
not_abs = randvec < absorb
absorb_model.cleanup_spectrum()
if eff_area is None:
detected = np.ones(eobs.shape, dtype='bool')
else:
mylog.info("Applying energy-dependent effective area.")
earea = np.interp(eobs, earf, eff_area, left=0.0, right=0.0)
randvec = eff_area.max()*prng.uniform(size=eobs.shape)
detected = randvec < earea
detected = np.logical_and(not_abs, detected)
events = {}
dx_min = self.parameters["Width"]/self.parameters["Dimension"]
dtheta = YTQuantity(np.rad2deg(dx_min/D_A), "degree")
events["xpix"] = xsky[detected]/dx_min.v + 0.5*(nx+1)
events["ypix"] = ysky[detected]/dx_min.v + 0.5*(nx+1)
events["eobs"] = eobs[detected]
events = comm.par_combine_object(events, datatype="dict", op="cat")
if psf_sigma is not None:
events["xpix"] += prng.normal(sigma=psf_sigma/dtheta)
events["ypix"] += prng.normal(sigma=psf_sigma/dtheta)
num_events = len(events["xpix"])
if comm.rank == 0:
mylog.info("Total number of observed photons: %d" % num_events)
if "RMF" in parameters and convolve_energies:
events, info = self._convolve_with_rmf(parameters["RMF"], events,
mat_key, prng)
for k, v in info.items():
parameters[k] = v
if exp_time_new is None:
parameters["ExposureTime"] = self.parameters["FiducialExposureTime"]
else:
parameters["ExposureTime"] = exp_time_new
if area_new is None:
parameters["Area"] = self.parameters["FiducialArea"]
else:
parameters["Area"] = area_new
parameters["Redshift"] = zobs
parameters["AngularDiameterDistance"] = D_A.in_units("Mpc")
parameters["sky_center"] = sky_center
parameters["pix_center"] = np.array([0.5*(nx+1)]*2)
parameters["dtheta"] = dtheta
return EventList(events, parameters)
def _normalize_arf(self, respfile, mat_key):
rmf = _astropy.pyfits.open(respfile)
table = rmf[mat_key]
weights = np.array([w.sum() for w in table.data["MATRIX"]])
rmf.close()
return weights
def _convolve_with_rmf(self, respfile, events, mat_key, prng):
"""
Convolve the events with a RMF file.
"""
mylog.info("Reading response matrix file (RMF): %s" % (respfile))
hdulist = _astropy.pyfits.open(respfile)
tblhdu = hdulist[mat_key]
n_de = len(tblhdu.data["ENERG_LO"])
mylog.info("Number of energy bins in RMF: %d" % (n_de))
mylog.info("Energy limits: %g %g" % (min(tblhdu.data["ENERG_LO"]),
max(tblhdu.data["ENERG_HI"])))
tblhdu2 = hdulist["EBOUNDS"]
n_ch = len(tblhdu2.data["CHANNEL"])
mylog.info("Number of channels in RMF: %d" % (n_ch))
eidxs = np.argsort(events["eobs"])
phEE = events["eobs"][eidxs].d
detectedChannels = []
# run through all photon energies and find which bin they go in
k = 0
fcurr = 0
last = len(phEE)
pbar = get_pbar("Scattering energies with RMF:", last)
for low,high in zip(tblhdu.data["ENERG_LO"],tblhdu.data["ENERG_HI"]):
# weight function for probabilities from RMF
weights = np.nan_to_num(np.float64(tblhdu.data[k]["MATRIX"][:]))
weights /= weights.sum()
# build channel number list associated to array value,
# there are groups of channels in rmfs with nonzero probabilities
trueChannel = []
f_chan = np.nan_to_num(tblhdu.data["F_CHAN"][k])
n_chan = np.nan_to_num(tblhdu.data["N_CHAN"][k])
n_grp = np.nan_to_num(tblhdu.data["N_CHAN"][k])
if not iterable(f_chan):
f_chan = [f_chan]
n_chan = [n_chan]
n_grp = [n_grp]
for start,nchan in zip(f_chan, n_chan):
end = start + nchan
if start == end:
trueChannel.append(start)
else:
for j in range(start,end):
trueChannel.append(j)
if len(trueChannel) > 0:
for q in range(fcurr,last):
if phEE[q] >= low and phEE[q] < high:
channelInd = prng.choice(len(weights), p=weights)
fcurr += 1
detectedChannels.append(trueChannel[channelInd])
if phEE[q] >= high:
break
pbar.update(fcurr)
k += 1
pbar.finish()
dchannel = np.array(detectedChannels)
events["xpix"] = events["xpix"][eidxs]
events["ypix"] = events["ypix"][eidxs]
events["eobs"] = YTArray(phEE, "keV")
events[tblhdu.header["CHANTYPE"]] = dchannel.astype(int)
info = {"ChannelType" : tblhdu.header["CHANTYPE"],
"Telescope" : tblhdu.header["TELESCOP"],
"Instrument" : tblhdu.header["INSTRUME"]}
info["Mission"] = tblhdu.header.get("MISSION","")
return events, info
class EventList(object):
def __init__(self, events, parameters):
self.events = events
self.parameters = parameters
self.num_events = events["xpix"].shape[0]
self.wcs = _astropy.pywcs.WCS(naxis=2)
self.wcs.wcs.crpix = parameters["pix_center"]
self.wcs.wcs.crval = parameters["sky_center"].d
self.wcs.wcs.cdelt = [-parameters["dtheta"].value, parameters["dtheta"].value]
self.wcs.wcs.ctype = ["RA---TAN","DEC--TAN"]
self.wcs.wcs.cunit = ["deg"]*2
def keys(self):
return self.events.keys()
def has_key(self, key):
return key in self.keys()
def items(self):
return self.events.items()
def values(self):
return self.events.values()
def __getitem__(self,key):
if key not in self.events:
if key == "xsky" or key == "ysky":
x,y = self.wcs.wcs_pix2world(self.events["xpix"], self.events["ypix"], 1)
self.events["xsky"] = YTArray(x, "degree")
self.events["ysky"] = YTArray(y, "degree")
return self.events[key]
def __repr__(self):
return self.events.__repr__()
def __contains__(self, key):
return key in self.events
def __add__(self, other):
assert_same_wcs(self.wcs, other.wcs)
validate_parameters(self.parameters, other.parameters)
events = {}
for item1, item2 in zip(self.items(), other.items()):
k1, v1 = item1
k2, v2 = item2
events[k1] = uconcatenate([v1,v2])
return EventList(events, self.parameters)
def filter_events(self, region):
"""
Filter events using a ds9 *region*. Requires the pyregion package.
Returns a new EventList.
"""
import pyregion
import os
if os.path.exists(region):
reg = pyregion.open(region)
else:
reg = pyregion.parse(region)
r = reg.as_imagecoord(header=self.wcs.to_header())
f = r.get_filter()
idxs = f.inside_x_y(self["xpix"], self["ypix"])
if idxs.sum() == 0:
raise RuntimeError("No events are inside this region!")
new_events = {}
for k, v in self.items():
new_events[k] = v[idxs]
return EventList(new_events, self.parameters)
@classmethod
def from_h5_file(cls, h5file):
"""
Initialize an EventList from a HDF5 file with filename *h5file*.
"""
events = {}
parameters = {}
f = h5py.File(h5file, "r")
p = f["/parameters"]
parameters["ExposureTime"] = YTQuantity(p["exp_time"].value, "s")
area = force_unicode(p['area'].value)
if isinstance(area, string_types):
parameters["Area"] = area
else:
parameters["Area"] = YTQuantity(area, "cm**2")
parameters["Redshift"] = p["redshift"].value
parameters["AngularDiameterDistance"] = YTQuantity(p["d_a"].value, "Mpc")
parameters["sky_center"] = YTArray(p["sky_center"][:], "deg")
parameters["dtheta"] = YTQuantity(p["dtheta"].value, "deg")
parameters["pix_center"] = p["pix_center"][:]
if "rmf" in p:
parameters["RMF"] = force_unicode(p["rmf"].value)
if "arf" in p:
parameters["ARF"] = force_unicode(p["arf"].value)
if "channel_type" in p:
parameters["ChannelType"] = force_unicode(p["channel_type"].value)
if "mission" in p:
parameters["Mission"] = force_unicode(p["mission"].value)
if "telescope" in p:
parameters["Telescope"] = force_unicode(p["telescope"].value)
if "instrument" in p:
parameters["Instrument"] = force_unicode(p["instrument"].value)
d = f["/data"]
events["xpix"] = d["xpix"][:]
events["ypix"] = d["ypix"][:]
events["eobs"] = YTArray(d["eobs"][:], "keV")
if "pi" in d:
events["PI"] = d["pi"][:]
if "pha" in d:
events["PHA"] = d["pha"][:]
f.close()
return cls(events, parameters)
@classmethod
def from_fits_file(cls, fitsfile):
"""
Initialize an EventList from a FITS file with filename *fitsfile*.
"""
hdulist = _astropy.pyfits.open(fitsfile)
tblhdu = hdulist["EVENTS"]
events = {}
parameters = {}
parameters["ExposureTime"] = YTQuantity(tblhdu.header["EXPOSURE"], "s")
if isinstance(tblhdu.header["AREA"], (string_types, bytes)):
parameters["Area"] = tblhdu.header["AREA"]
else:
parameters["Area"] = YTQuantity(tblhdu.header["AREA"], "cm**2")
parameters["Redshift"] = tblhdu.header["REDSHIFT"]
parameters["AngularDiameterDistance"] = YTQuantity(tblhdu.header["D_A"], "Mpc")
if "RMF" in tblhdu.header:
parameters["RMF"] = tblhdu["RMF"]
if "ARF" in tblhdu.header:
parameters["ARF"] = tblhdu["ARF"]
if "CHANTYPE" in tblhdu.header:
parameters["ChannelType"] = tblhdu["CHANTYPE"]
if "MISSION" in tblhdu.header:
parameters["Mission"] = tblhdu["MISSION"]
if "TELESCOP" in tblhdu.header:
parameters["Telescope"] = tblhdu["TELESCOP"]
if "INSTRUME" in tblhdu.header:
parameters["Instrument"] = tblhdu["INSTRUME"]
parameters["sky_center"] = YTArray([tblhdu["TCRVL2"],tblhdu["TCRVL3"]], "deg")
parameters["pix_center"] = np.array([tblhdu["TCRVL2"],tblhdu["TCRVL3"]])
parameters["dtheta"] = YTQuantity(tblhdu["TCRVL3"], "deg")
events["xpix"] = tblhdu.data.field("X")
events["ypix"] = tblhdu.data.field("Y")
events["eobs"] = YTArray(tblhdu.data.field("ENERGY")/1000., "keV")