def make_point_sources(area, exp_time, positions, sky_center, spectra, prng=None): r""" Create a new :class:`~pyxsim.event_list.EventList` which contains point sources. Parameters ---------- area : float, (value, unit) tuple, :class:`~yt.units.yt_array.YTQuantity`, or :class:`~astropy.units.Quantity` The collecting area to determine the number of events. If units are not specified, it is assumed to be in cm^2. exp_time : float, (value, unit) tuple, :class:`~yt.units.yt_array.YTQuantity`, or :class:`~astropy.units.Quantity` The exposure time to determine the number of events. If units are not specified, it is assumed to be in seconds. positions : array of source positions, shape 2xN The positions of the point sources in RA, Dec, where N is the number of point sources. Coordinates should be in degrees. sky_center : array-like Center RA, Dec of the events in degrees. spectra : list (size N) of :class:`~soxs.spectra.Spectrum` objects The spectra for the point sources, where N is the number of point sources. Assumed to be in the observer frame. prng : integer or :class:`~numpy.random.RandomState` object 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 to use the :mod:`numpy.random` module. """ prng = parse_prng(prng) spectra = ensure_list(spectra) positions = ensure_list(positions) area = parse_value(area, "cm**2") exp_time = parse_value(exp_time, "s") t_exp = exp_time.value/comm.size x = [] y = [] e = [] for pos, spectrum in zip(positions, spectra): eobs = spectrum.generate_energies(t_exp, area.value, prng=prng) ne = eobs.size x.append(YTArray([pos[0]] * ne, "degree")) y.append(YTArray([pos[1]] * ne, "degree")) e.append(YTArray.from_astropy(eobs)) parameters = {"sky_center": YTArray(sky_center, "degree"), "exp_time": exp_time, "area": area} events = {} events["eobs"] = uconcatenate(e) events["xsky"] = uconcatenate(x) events["ysky"] = uconcatenate(y) return EventList(events, parameters)
def test_background(): tmpdir = tempfile.mkdtemp() curdir = os.getcwd() os.chdir(tmpdir) kT_sim = 1.0 Z_sim = 0.0 norm_sim = 4.0e-2 nH_sim = 0.04 redshift = 0.01 exp_time = (200., "ks") area = (1000., "cm**2") wcs = create_dummy_wcs() abs_model = WabsModel(nH_sim) events = EventList.create_empty_list(exp_time, area, wcs) spec_model = TableApecModel(0.05, 12.0, 5000, thermal_broad=False) spec = spec_model.return_spectrum(kT_sim, Z_sim, redshift, norm_sim) new_events = events.add_background(spec_model.ebins, spec, prng=prng, absorb_model=abs_model) new_events = ACIS_I(new_events, rebin=False, convolve_psf=False, prng=prng) new_events.write_spectrum("background_evt.pi", clobber=True) os.system("cp %s ." % new_events.parameters["ARF"]) os.system("cp %s ." % new_events.parameters["RMF"]) load_user_model(mymodel, "wapec") add_user_pars("wapec", ["nH", "kT", "metallicity", "redshift", "norm"], [0.01, 4.0, 0.2, redshift, norm_sim*0.8], parmins=[0.0, 0.1, 0.0, -20.0, 0.0], parmaxs=[10.0, 20.0, 10.0, 20.0, 1.0e9], parfrozen=[False, False, False, True, False]) load_pha("background_evt.pi") set_stat("cstat") set_method("simplex") ignore(":0.5, 8.0:") set_model("wapec") fit() set_covar_opt("sigma", 1.6) covar() res = get_covar_results() assert np.abs(res.parvals[0]-nH_sim) < res.parmaxes[0] assert np.abs(res.parvals[1]-kT_sim) < res.parmaxes[1] assert np.abs(res.parvals[2]-Z_sim) < res.parmaxes[2] assert np.abs(res.parvals[3]-norm_sim) < res.parmaxes[3] os.chdir(curdir) shutil.rmtree(tmpdir)
def make_background(area, exp_time, fov, sky_center, spectrum, prng=None): r""" Create a new :class:`~pyxsim.event_list.EventList` which is filled uniformly with background events. Parameters ---------- area : float, (value, unit) tuple, :class:`~yt.units.yt_array.YTQuantity`, or :class:`~astropy.units.Quantity` The collecting area to determine the number of events. If units are not specified, it is assumed to be in cm^2. exp_time : float, (value, unit) tuple, :class:`~yt.units.yt_array.YTQuantity`, or :class:`~astropy.units.Quantity` The exposure time to determine the number of events. If units are not specified, it is assumed to be in seconds. fov : float, (value, unit) tuple, :class:`~yt.units.yt_array.YTQuantity`, or :class:`~astropy.units.Quantity` The field of view of the event file. If units are not provided, they are assumed to be in arcminutes. sky_center : array-like Center RA, Dec of the events in degrees. spectrum : :class:`~soxs.spectra.Spectrum` The spectrum for the background. prng : integer or :class:`~numpy.random.RandomState` object 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 to use the :mod:`numpy.random` module. """ prng = parse_prng(prng) fov = parse_value(fov, "arcmin") exp_time = parse_value(exp_time, "s") area = parse_value(area, "cm**2") t_exp = exp_time.value / comm.size e = spectrum.generate_energies(t_exp, area.value, prng=prng) fov_model = FillFOVModel(sky_center[0], sky_center[1], fov.value) ra, dec = fov_model.generate_coords(e.size, prng=prng) parameters = { "sky_center": YTArray(sky_center, "degree"), "exp_time": exp_time, "area": area } events = {} events["xsky"] = YTArray(ra.value, "degree") events["ysky"] = YTArray(dec.value, "degree") events["eobs"] = YTArray(e.value, "keV") return EventList(events, parameters)
def __call__(self, events, rebin=True, convolve_psf=True, convolve_arf=True, convolve_rmf=True, prng=None): new_events = EventList(deepcopy(events.events), events.parameters.copy(), events.wcs.copy()) if prng is None: prng = np.random if rebin: self.rebin(new_events) if convolve_psf: self.convolve_with_psf(new_events, prng) if convolve_arf: new_events["xsky"] new_events["ysky"] self.apply_effective_area(new_events, prng) if convolve_rmf: self.convolve_energies(new_events, prng) return new_events
def test_point_source(): tmpdir = tempfile.mkdtemp() curdir = os.getcwd() os.chdir(tmpdir) nH_sim = 0.02 norm_sim = 1.0e-4 alpha_sim = 0.95 redshift = 0.02 exp_time = (100., "ks") area = (3000., "cm**2") wcs = create_dummy_wcs() ebins = np.linspace(0.1, 11.5, 2001) emid = 0.5*(ebins[1:]+ebins[:-1]) spec = norm_sim*(emid*(1.0+redshift))**(-alpha_sim) de = np.diff(ebins)[0] abs_model = TBabsModel(nH_sim) events = EventList.create_empty_list(exp_time, area, wcs) positions = [(30.01, 45.0)] new_events = events.add_point_sources(positions, ebins, spec, prng=prng, absorb_model=abs_model) new_events = ACIS_S(new_events, prng=prng) scalex = float(np.std(new_events['xpix'])*sigma_to_fwhm*new_events.parameters["dtheta"]) scaley = float(np.std(new_events['ypix'])*sigma_to_fwhm*new_events.parameters["dtheta"]) psf_scale = ACIS_S.psf_scale assert (scalex - psf_scale)/psf_scale < 0.01 assert (scaley - psf_scale)/psf_scale < 0.01 new_events.write_spectrum("point_source_evt.pi", clobber=True) os.system("cp %s ." % new_events.parameters["ARF"]) os.system("cp %s ." % new_events.parameters["RMF"]) load_user_model(mymodel, "tplaw") add_user_pars("tplaw", ["nH", "norm", "redshift", "alpha"], [0.01, norm_sim*0.8, redshift, 0.9], parmins=[0.0, 0.0, 0.0, 0.1], parmaxs=[10.0, 1.0e9, 10.0, 10.0], parfrozen=[False, False, True, False]) load_pha("point_source_evt.pi") set_stat("cstat") set_method("simplex") ignore(":0.5, 9.0:") set_model("tplaw") fit() set_covar_opt("sigma", 1.6) covar() res = get_covar_results() assert np.abs(res.parvals[0]-nH_sim) < res.parmaxes[0] assert np.abs(res.parvals[1]-norm_sim/de) < res.parmaxes[1] assert np.abs(res.parvals[2]-alpha_sim) < res.parmaxes[2] os.chdir(curdir) shutil.rmtree(tmpdir)
def generate_events(self, area, exp_time, angular_width, source_model, sky_center, parameters=None, velocity_fields=None, absorb_model=None, nH=None, no_shifting=False, sigma_pos=None, prng=None): """ Generate projected events from a light cone simulation. Parameters ---------- area : float, (value, unit) tuple, or :class:`~yt.units.yt_array.YTQuantity` The collecting area to determine the number of events. If units are not specified, it is assumed to be in cm^2. exp_time : float, (value, unit) tuple, or :class:`~yt.units.yt_array.YTQuantity` The exposure time to determine the number of events. If units are not specified, it is assumed to be in seconds. angular_width : float, (value, unit) tuple, or :class:`~yt.units.yt_array.YTQuantity` The angular width of the light cone simulation. If units are not specified, it is assumed to be in degrees. source_model : :class:`~pyxsim.source_models.SourceModel` A source model used to generate the events. sky_center : array-like Center RA, Dec of the events in degrees. parameters : dict, optional A dictionary of parameters to be passed for the source model to use, if necessary. velocity_fields : list of fields The yt fields to use for the velocity. If not specified, the following will be assumed: ['velocity_x', 'velocity_y', 'velocity_z'] for grid datasets ['particle_velocity_x', 'particle_velocity_y', 'particle_velocity_z'] for particle datasets absorb_model : string or :class:`~pyxsim.spectral_models.AbsorptionModel` A model for foreground galactic absorption, to simulate the absorption of events before being detected. This cannot be applied here if you already did this step previously in the creation of the :class:`~pyxsim.photon_list.PhotonList` instance. Known options for strings are "wabs" and "tbabs". nH : float, optional The foreground column density in units of 10^22 cm^{-2}. Only used if absorption is applied. no_shifting : boolean, optional If set, the photon energies will not be Doppler shifted. sigma_pos : float, optional Apply a gaussian smoothing operation to the sky positions of the events. This may be useful when the binned events appear blocky due to their uniform distribution within simulation cells. However, this will move the events away from their originating position on the sky, and so may distort surface brightness profiles and/or spectra. Should probably only be used for visualization purposes. Supply a float here to smooth with a standard deviation with this fraction of the cell size. Default: None prng : integer or :class:`~numpy.random.RandomState` object 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 to use the :mod:`numpy.random` module. """ prng = parse_prng(prng) area = parse_value(area, "cm**2") exp_time = parse_value(exp_time, "s") aw = parse_value(angular_width, "deg") tot_events = defaultdict(list) for output in self.light_cone_solution: ds = load(output["filename"]) ax = output["projection_axis"] c = output[ "projection_center"] * ds.domain_width + ds.domain_left_edge le = c.copy() re = c.copy() width = ds.quan(aw * output["box_width_per_angle"], "unitary").to("code_length") depth = ds.domain_width[ax].in_units( "code_length") * output["box_depth_fraction"] le[ax] -= 0.5 * depth re[ax] += 0.5 * depth for off_ax in axes_lookup[ax]: le[off_ax] -= 0.5 * width re[off_ax] += 0.5 * width reg = ds.box(le, re) photons = PhotonList.from_data_source( reg, output['redshift'], area, exp_time, source_model, parameters=parameters, center=c, velocity_fields=velocity_fields, cosmology=ds.cosmology) if sum(photons["num_photons"]) > 0: events = photons.project_photons("xyz"[ax], sky_center, absorb_model=absorb_model, nH=nH, no_shifting=no_shifting, sigma_pos=sigma_pos, prng=prng) if events.num_events > 0: tot_events["xsky"].append(events["xsky"]) tot_events["ysky"].append(events["ysky"]) tot_events["eobs"].append(events["eobs"]) del events del photons parameters = { "exp_time": exp_time, "area": area, "sky_center": YTArray(sky_center, "deg") } for key in tot_events: tot_events[key] = uconcatenate(tot_events[key]) return EventList(tot_events, parameters)
def project_photons(self, normal, sky_center, absorb_model=None, nH=None, no_shifting=False, north_vector=None, sigma_pos=None, kernel="top_hat", prng=None, **kwargs): r""" Projects photons onto an image plane given a line of sight. Returns a new :class:`~pyxsim.event_list.EventList`. 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] sky_center : array-like Center RA, Dec of the events in degrees. absorb_model : string or :class:`~pyxsim.spectral_models.AbsorptionModel` A model for foreground galactic absorption, to simulate the absorption of events before being detected. This cannot be applied here if you already did this step previously in the creation of the :class:`~pyxsim.photon_list.PhotonList` instance. Known options for strings are "wabs" and "tbabs". nH : float, optional The foreground column density in units of 10^22 cm^{-2}. Only used if absorption is applied. 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. sigma_pos : float, optional Apply a gaussian smoothing operation to the sky positions of the events. This may be useful when the binned events appear blocky due to their uniform distribution within simulation cells. However, this will move the events away from their originating position on the sky, and so may distort surface brightness profiles and/or spectra. Should probably only be used for visualization purposes. Supply a float here to smooth with a standard deviation with this fraction of the cell size. Default: None kernel : string, optional The kernel used when smoothing positions of X-rays originating from SPH particles, "gaussian" or "top_hat". Default: "top_hat". prng : integer or :class:`~numpy.random.RandomState` object 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 to use the :mod:`numpy.random` module. Examples -------- >>> L = np.array([0.1,-0.2,0.3]) >>> events = my_photons.project_photons(L, [30., 45.]) """ prng = parse_prng(prng) scale_shift = -1.0 / clight.to("km/s") if "smooth_positions" in kwargs: issue_deprecation_warning( "'smooth_positions' has been renamed to " "'sigma_pos' and the former is deprecated!") sigma_pos = kwargs["smooth_positions"] if "redshift_new" in kwargs or "area_new" in kwargs or \ "exp_time_new" in kwargs or "dist_new" in kwargs: issue_deprecation_warning( "Changing the redshift, distance, area, or " "exposure time has been deprecated in " "project_photons!") if sigma_pos is not None and self.parameters[ "data_type"] == "particles": raise RuntimeError( "The 'smooth_positions' argument should not be used with " "particle-based datasets!") if isinstance(absorb_model, string_types): if absorb_model not in absorb_models: raise KeyError("%s is not a known absorption model!" % absorb_model) absorb_model = absorb_models[absorb_model] if absorb_model is not None: if nH is None: raise RuntimeError( "You specified an absorption model, but didn't " "specify a value for nH!") absorb_model = absorb_model(nH) sky_center = YTArray(sky_center, "degree") n_ph = self.photons["num_photons"] 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] else: x_hat = np.zeros(3) y_hat = np.zeros(3) z_hat = np.zeros(3) parameters = {} D_A = self.parameters["fid_d_a"] events = {} eobs = self.photons["energy"].v if not no_shifting: if comm.rank == 0: mylog.info("Doppler-shifting photon energies.") if isinstance(normal, string_types): shift = self.photons["vel"][:, "xyz".index(normal)] * scale_shift else: shift = np.dot(self.photons["vel"], z_hat) * scale_shift doppler_shift(shift, n_ph, eobs) if absorb_model is None: det = np.ones(eobs.size, dtype='bool') num_det = eobs.size else: if comm.rank == 0: mylog.info("Foreground galactic absorption: using " "the %s model and nH = %g." % (absorb_model._name, nH)) det = absorb_model.absorb_photons(eobs, prng=prng) num_det = det.sum() events["eobs"] = YTArray(eobs[det], "keV") num_events = comm.mpi_allreduce(num_det) if comm.rank == 0: mylog.info("%d events have been detected." % num_events) if num_det > 0: if comm.rank == 0: mylog.info("Assigning positions to events.") if isinstance(normal, string_types): norm = "xyz".index(normal) else: norm = normal xsky, ysky = scatter_events(norm, prng, kernel, self.parameters["data_type"], num_det, det, self.photons["num_photons"], self.photons["pos"].d, self.photons["dx"].d, x_hat, y_hat) if self.parameters[ "data_type"] == "cells" and sigma_pos is not None: if comm.rank == 0: mylog.info("Optionally smoothing sky positions.") sigma = sigma_pos * np.repeat(self.photons["dx"].d, n_ph)[det] xsky += sigma * prng.normal(loc=0.0, scale=1.0, size=num_det) ysky += sigma * prng.normal(loc=0.0, scale=1.0, size=num_det) d_a = D_A.to("kpc").v xsky /= d_a ysky /= d_a if comm.rank == 0: mylog.info("Converting pixel to sky coordinates.") pixel_to_cel(xsky, ysky, sky_center) else: xsky = [] ysky = [] events["xsky"] = YTArray(xsky, "degree") events["ysky"] = YTArray(ysky, "degree") parameters["exp_time"] = self.parameters["fid_exp_time"] parameters["area"] = self.parameters["fid_area"] parameters["sky_center"] = sky_center return EventList(events, parameters)
def project_photons(self, normal, area_new=None, exp_time_new=None, redshift_new=None, dist_new=None, absorb_model=None, sky_center=None, no_shifting=False, north_vector=None, prng=None): r""" Projects photons onto an image plane given a line of sight. Returns a new :class:`~pyxsim.event_list.EventList`. 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, (value, unit) tuple, or :class:`~yt.units.yt_array.YTQuantity`, optional New value for the (constant) collecting area of the detector. If units are not specified, is assumed to be in cm**2. exp_time_new : float, (value, unit) tuple, or :class:`~yt.units.yt_array.YTQuantity`, optional The new value for the exposure time. If units are not specified it is assumed to be in seconds. redshift_new : float, optional The new value for the cosmological redshift. dist_new : float, (value, unit) tuple, or :class:`~yt.units.yt_array.YTQuantity`, optional The new value for the angular diameter distance, used for nearby sources. This may be optionally supplied instead of it being determined from the cosmology. If units are not specified, it is assumed to be in Mpc. To use this, the redshift must be zero. absorb_model : :class:`~pyxsim.spectral_models.AbsorptionModel` A model for foreground galactic absorption. sky_center : array-like, optional Center RA, Dec of the events in degrees. 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 : :class:`~numpy.random.RandomState` object or :mod:`~numpy.random`, optional 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 :mod:`numpy.random` module. Examples -------- >>> L = np.array([0.1,-0.2,0.3]) >>> events = my_photons.project_photons(L, area_new=10000., ... redshift_new=0.05) """ if prng is None: prng = np.random 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 if isinstance(normal, string_types): # if on-axis, just use the maximum width of the plane perpendicular # to that axis w = self.parameters["Width"].copy() w["xyz".index(normal)] = 0.0 ax_idx = np.argmax(w) else: # if off-axis, just use the largest width to make sure we get everything ax_idx = np.argmax(self.parameters["Width"]) nx = self.parameters["Dimension"][ax_idx] dx_min = (self.parameters["Width"] / self.parameters["Dimension"])[ax_idx] 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() parameters = {} zobs0 = self.parameters["FiducialRedshift"] D_A0 = self.parameters["FiducialAngularDiameterDistance"] scale_factor = 1.0 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: exp_time_new = parse_value(exp_time_new, "s") Tratio = exp_time_new / self.parameters["FiducialExposureTime"] if area_new is None: Aratio = 1. else: area_new = parse_value(area_new, "cm**2") Aratio = area_new / self.parameters["FiducialArea"] if redshift_new is None and dist_new is None: Dratio = 1. zobs = zobs0 D_A = D_A0 else: if dist_new is not None: if redshift_new is not None and redshift_new > 0.0: mylog.warning( "Redshift must be zero for nearby sources. Resetting redshift to 0.0." ) 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.int64(n_ph_tot * fak) Nn = 4294967294 if my_n_obs == n_ph_tot: if my_n_obs <= Nn: idxs = np.arange(my_n_obs, dtype='uint32') else: idxs = np.arange(my_n_obs, dtype='uint64') else: if n_ph_tot <= Nn: idxs = np.arange(n_ph_tot, dtype='uint32') prng.shuffle(idxs) idxs = idxs[:my_n_obs] else: Nc = np.int32(n_ph_tot / Nn) idxs = np.zeros(my_n_obs, dtype=np.uint64) Nup = np.uint32(my_n_obs / Nc) for i in range(Nc + 1): if (i + 1) * Nc < n_ph_tot: idtm = np.arange(i * Nc, (i + 1) * Nc, dtype='uint64') Nupt = Nup else: idtm = np.arange(i * Nc, n_ph_tot, dtype='uint64') Nupt = my_n_obs - i * Nup prng.shuffle(idtm) idxs[i * Nup, i * Nup + Nupt] = idtm[:Nupt] del (idtm) # idxs = prng.permutation(n_ph_tot)[:my_n_obs].astype("int64") obs_cells = np.searchsorted(self.p_bins, idxs, side='right') - 1 delta = dx[obs_cells] if isinstance(normal, string_types): if self.parameters["DataType"] == "cells": 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) elif self.parameters["DataType"] == "particles": xsky = prng.normal(loc=0.0, scale=1.0, size=my_n_obs) ysky = prng.normal(loc=0.0, scale=1.0, size=my_n_obs) xsky *= delta ysky *= delta xsky += self.photons[axes_lookup[normal][0]].d[obs_cells] ysky += self.photons[axes_lookup[normal][1]].d[obs_cells] if not no_shifting: vz = self.photons["v%s" % normal] else: if self.parameters["DataType"] == "cells": 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) elif self.parameters["DataType"] == "particles": x = prng.normal(loc=0.0, scale=1.0, size=my_n_obs) y = prng.normal(loc=0.0, scale=1.0, size=my_n_obs) z = prng.normal(loc=0.0, scale=1.0, 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"].d[obs_cells] y += self.photons["y"].d[obs_cells] z += self.photons["z"].d[obs_cells] 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] del (delta) 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] shift = -vz[obs_cells].in_cgs() / clight shift = np.sqrt((1. - shift) / (1. + shift)) eobs = self.photons["Energy"][idxs] eobs *= shift del (shift) eobs *= scale_factor if absorb_model is None: detected = np.ones(eobs.shape, dtype='bool') else: detected = absorb_model.absorb_photons(eobs, prng=prng) events = {} 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") num_events = len(events["xpix"]) if comm.rank == 0: mylog.info("Total number of observed photons: %d" % num_events) 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)