def get_efield_antenna_factor(station, frequencies, channels, detector, zenith, azimuth, antenna_pattern_provider): """ Returns the antenna response to a radio signal coming from a specific direction Parameters --------------- station: Station frequencies: array of complex frequencies of the radio signal for which the antenna response is needed channels: array of int IDs of the channels detector: Detector zenith, azimuth: float, float incoming direction of the signal. Note that refraction and reflection at the ice/air boundary are taken into account antenna_pattern_provider: AntennaPatternProvider """ n_ice = ice.get_refractive_index(-0.01, detector.get_site(station.get_id())) efield_antenna_factor = np.zeros((len(channels), 2, len(frequencies)), dtype=np.complex) # from antenna model in e_theta, e_phi for iCh, channel_id in enumerate(channels): zenith_antenna = zenith t_theta = 1. t_phi = 1. # first check case if signal comes from above if zenith <= 0.5 * np.pi and station.is_cosmic_ray(): # is antenna below surface? position = detector.get_relative_position(station.get_id(), channel_id) if position[2] <= 0: zenith_antenna = geo_utl.get_fresnel_angle(zenith, n_ice, 1) t_theta = geo_utl.get_fresnel_t_p(zenith, n_ice, 1) t_phi = geo_utl.get_fresnel_t_s(zenith, n_ice, 1) logger.info("channel {:d}: electric field is refracted into the firn. theta {:.0f} -> {:.0f}. Transmission coefficient p (eTheta) {:.2f} s (ePhi) {:.2f}".format(iCh, zenith / units.deg, zenith_antenna / units.deg, t_theta, t_phi)) else: # now the signal is coming from below, do we have an antenna above the surface? position = detector.get_relative_position(station.get_id(), channel_id) if(position[2] > 0): zenith_antenna = geo_utl.get_fresnel_angle(zenith, 1., n_ice) if(zenith_antenna is None): logger.warning("fresnel reflection at air-firn boundary leads to unphysical results, no reconstruction possible") return None logger.debug("angles: zenith {0:.0f}, zenith antenna {1:.0f}, azimuth {2:.0f}".format(np.rad2deg(zenith), np.rad2deg(zenith_antenna), np.rad2deg(azimuth))) antenna_model = detector.get_antenna_model(station.get_id(), channel_id, zenith_antenna) antenna_pattern = antenna_pattern_provider.load_antenna_pattern(antenna_model) ori = detector.get_antenna_orientation(station.get_id(), channel_id) VEL = antenna_pattern.get_antenna_response_vectorized(frequencies, zenith_antenna, azimuth, *ori) efield_antenna_factor[iCh] = np.array([VEL['theta'] * t_theta, VEL['phi'] * t_phi]) return efield_antenna_factor
def run(self, evt, station, det, channels_to_use=None, cosmic_ray=False): """ Fits the direction using templates Parameters ---------- evt: event station: station det: detector channels_to_use: list (default: [0, 1, 2, 3] antenna to use for fit cosmic_ray: bool type to set correlation template """ if channels_to_use is None: channels_to_use = [0, 1, 2, 3] if (cosmic_ray): type_str = 'cr' xcorrelations = chp.cr_xcorrelations else: type_str = 'nu' xcorrelations = chp.nu_xcorrelations station_id = station.get_id() channels = station.iter_channels(channels_to_use) times = [] positions = [] for iCh, channel in enumerate(channels): channel_id = channel.get_id() times.append( channel[xcorrelations]['{}_ref_xcorr_time'.format(type_str)] + channel.get_trace_start_time()) positions.append(det.get_relative_position(station_id, channel_id)) times = np.array(times) positions = np.array(positions) site = det.get_site(station_id) n_ice = ice.get_refractive_index(-0.01, site) from scipy import optimize as opt def obj_plane(params, positions, t_measured): zenith, azimuth = params if cosmic_ray: if ((zenith < 0) or (zenith > 0.5 * np.pi)): return np.inf else: if ((zenith < 0.5 * np.pi) or (zenith > np.pi)): return np.inf v = hp.spherical_to_cartesian(zenith, azimuth) c = constants.c * units.m / units.s if not cosmic_ray: c = c / n_ice logger.debug("using speed of light = {:.4g}".format(c)) t_expected = -(np.dot(v, positions.T) / c) sigma = 1 * units.ns chi2 = np.sum(((t_expected - t_expected.mean()) - (t_measured - t_measured.mean()))**2 / sigma**2) logger.debug("texp = {texp}, tm = {tmeas}, {chi2}".format( texp=t_expected, tmeas=t_measured, chi2=chi2)) return chi2 method = "Nelder-Mead" options = {'maxiter': 1000, 'disp': False} zenith_start = 135 * units.deg if cosmic_ray: zenith_start = 45 * units.deg starting_chi2 = {} for starting_az in np.array([0, 90, 180, 270]) * units.degree: starting_chi2[starting_az] = obj_plane((zenith_start, starting_az), positions, times) azimuth_start = min(starting_chi2, key=starting_chi2.get) res = opt.minimize(obj_plane, x0=[zenith_start, azimuth_start], args=(positions, times), method=method, options=options) output_str = "reconstucted angles theta = {:.1f}, phi = {:.1f}".format( res.x[0] / units.deg, hp.get_normalized_angle(res.x[1]) / units.deg) if station.has_sim_station(): sim_zen = station.get_sim_station()[stnp.zenith] sim_az = station.get_sim_station()[stnp.azimuth] dOmega = hp.get_angle( hp.spherical_to_cartesian(sim_zen, sim_az), hp.spherical_to_cartesian(res.x[0], res.x[1])) output_str += " MC theta = {:.1f}, phi = {:.1f}, dOmega = {:.2f}".format( sim_zen / units.deg, sim_az / units.deg, dOmega / units.deg) logger.info(output_str) station[stnp.zenith] = res.x[0] station[stnp.azimuth] = hp.get_normalized_angle(res.x[1]) if (cosmic_ray): station[stnp.cr_zenith] = res.x[0] station[stnp.cr_azimuth] = hp.get_normalized_angle(res.x[1]) else: station[stnp.nu_zenith] = res.x[0] station[stnp.nu_azimuth] = hp.get_normalized_angle(res.x[1])
def run( self, event, station, detector, passband=None ): """ Adds noise resulting from galactic radio emission to the channel traces Parameters -------------- event: Event object The event containing the station to whose channels noise shall be added station: Station object The station whose channels noise shall be added to detector: Detector object The detector description passband: list of float Lower and upper bound of the frequency range in which noise shall be added """ if passband is None: passband = [10 * units.MHz, 1000 * units.MHz] self.__gdsm = pygdsm.pygsm.GlobalSkyModel() site_latitude, site_longitude = detector.get_site_coordinates(station.get_id()) site_location = astropy.coordinates.EarthLocation(lat=site_latitude * astropy.units.deg, lon=site_longitude * astropy.units.deg) station_time = station.get_station_time() station_time.format = 'iso' noise_temperatures = np.zeros((len(self.__interpolaiton_frequencies), healpy.pixelfunc.nside2npix(self.__n_side))) local_cs = astropy.coordinates.AltAz(location=site_location, obstime=station_time) solid_angle = healpy.pixelfunc.nside2pixarea(self.__n_side, degrees=False) pixel_longitudes, pixel_latitudes = healpy.pixelfunc.pix2ang(self.__n_side, range(healpy.pixelfunc.nside2npix(self.__n_side)), lonlat=True) pixel_longitudes *= units.deg pixel_latitudes *= units.deg galactic_coordinates = astropy.coordinates.Galactic(l=pixel_longitudes * astropy.units.rad, b=pixel_latitudes * astropy.units.rad) local_coordinates = galactic_coordinates.transform_to(local_cs) n_ice = ice.get_refractive_index(-0.01, detector.get_site(station.get_id())) # save noise temperatures for all directions and frequencies for i_freq, noise_freq in enumerate(self.__interpolaiton_frequencies): radio_sky = self.__gdsm.generate(noise_freq / units.MHz) radio_sky = healpy.pixelfunc.ud_grade(radio_sky, self.__n_side) noise_temperatures[i_freq] = radio_sky for channel in station.iter_channels(): antenna_pattern = self.__antenna_pattern_provider.load_antenna_pattern(detector.get_antenna_model(station.get_id(), channel.get_id())) freqs = channel.get_frequencies() d_f = freqs[2] - freqs[1] sampling_rate = channel.get_sampling_rate() channel_spectrum = channel.get_frequency_spectrum() passband_filter = (freqs > passband[0]) & (freqs < passband[1]) noise_spec_sum = np.zeros_like(channel.get_frequency_spectrum()) flux_sum = np.zeros(freqs[passband_filter].shape) efield_sum = np.zeros((3, freqs.shape[0]), dtype=np.complex) if self.__debug: plt.close('all') fig = plt.figure(figsize=(12, 8)) ax_1 = fig.add_subplot(221) ax_2 = fig.add_subplot(222) ax_3 = fig.add_subplot(223) ax_4 = fig.add_subplot(224) ax_1.grid() ax_1.set_yscale('log') ax_2.grid() ax_2.set_yscale('log') ax_3.grid() ax_3.set_yscale('log') ax_4.grid() ax_4.set_yscale('log') ax_1.set_xlabel('f [MHz]') ax_2.set_xlabel('f [MHz]') ax_3.set_xlabel('f [MHz]') ax_4.set_xlabel('f [MHz]') ax_1.set_ylabel('T [K]') ax_2.set_ylabel('S [W/m²/MHz]') ax_3.set_ylabel('E [V/m]') ax_4.set_ylabel('U [V]') ax_4.plot(channel.get_frequencies() / units.MHz, np.abs(channel.get_frequency_spectrum()), c='C0') ax_4.set_ylim([1.e-8, None]) fig1 = plt.figure() ax1 = fig1.add_subplot(211) ax2 = fig1.add_subplot(212) ax1.plot(channel.get_times(), channel.get_trace(), label='original trace') ax2.plot(channel.get_frequencies() / units.MHz, np.abs(channel.get_frequency_spectrum()), label='original spectrum') ax1.grid() ax2.grid() for i_pixel in range(healpy.pixelfunc.nside2npix(self.__n_side)): azimuth = local_coordinates[i_pixel].az.rad zenith = np.pi / 2. - local_coordinates[i_pixel].alt.rad if zenith > 90. * units.deg: continue temperature_interpolator = scipy.interpolate.interp1d(self.__interpolaiton_frequencies, np.log10(noise_temperatures[:, i_pixel]), kind='quadratic') noise_temperature = np.power(10, temperature_interpolator(freqs[passband_filter])) # calculate spectral radiance of radio signal using rayleigh-jeans law S = (2. * (scipy.constants.Boltzmann * units.joule / units.kelvin) * freqs[passband_filter]**2 / (scipy.constants.c * units.m / units.s)**2 * (noise_temperature) * solid_angle) S[np.isnan(S)] = 0 # calculate radiance per energy bin S_per_bin = S * d_f flux_sum += S_per_bin # calculate electric field per energy bin from the radiance per bin E = np.sqrt(S_per_bin / (scipy.constants.c * units.m / units.s * scipy.constants.epsilon_0 * (units.coulomb / units.V / units.m))) / (d_f) if self.__debug: ax_1.scatter(self.__interpolaiton_frequencies / units.MHz, noise_temperatures[:, i_pixel] / units.kelvin, c='k', alpha=.01) ax_1.plot(freqs[passband_filter] / units.MHz, noise_temperature, c='k', alpha=.02) ax_2.plot(freqs[passband_filter] / units.MHz, S_per_bin / d_f / (units.watt / units.m**2 / units.MHz), c='k', alpha=.02) ax_3.plot(freqs[passband_filter] / units.MHz, E / (units.V / units.m), c='k', alpha=.02) # assign random phases and polarizations to electric field noise_spectrum = np.zeros((3, freqs.shape[0]), dtype=np.complex) phases = np.random.uniform(0, 2. * np.pi, len(S)) polarizations = np.random.uniform(0, 2. * np.pi, len(S)) noise_spectrum[1][passband_filter] = np.exp(1j * phases) * E * np.cos(polarizations) noise_spectrum[2][passband_filter] = np.exp(1j * phases) * E * np.sin(polarizations) efield_sum += noise_spectrum antenna_orientation = detector.get_antenna_orientation(station.get_id(), channel.get_id()) # consider signal reflection at ice surface if detector.get_relative_position(station.get_id(), channel.get_id())[2] < 0: t_theta = geometryUtilities.get_fresnel_t_p(zenith, n_ice, 1) t_phi = geometryUtilities.get_fresnel_t_s(zenith, n_ice, 1) fresnel_zenith = geometryUtilities.get_fresnel_angle(zenith, 1., n_ice) if fresnel_zenith is None: continue else: t_theta = 1 t_phi = 1 fresnel_zenith = zenith # fold electric field with antenna response antenna_response = antenna_pattern.get_antenna_response_vectorized(freqs, fresnel_zenith, azimuth, *antenna_orientation) channel_noise_spectrum = antenna_response['theta'] * noise_spectrum[1] * t_theta + antenna_response['phi'] * noise_spectrum[2] * t_phi if self.__debug: ax_4.plot(freqs / units.MHz, np.abs(channel_noise_spectrum) / units.V, c='k', alpha=.01) channel_spectrum += channel_noise_spectrum noise_spec_sum += channel_noise_spectrum channel.set_frequency_spectrum(channel_spectrum, sampling_rate) if self.__debug: ax_2.plot(freqs[passband_filter] / units.MHz, flux_sum / d_f / (units.watt / units.m**2 / units.MHz), c='C0', label='total flux') ax_3.plot(freqs[passband_filter] / units.MHz, np.sqrt(np.abs(efield_sum[1])**2 + np.abs(efield_sum[2])**2)[passband_filter] / (units.V / units.m), c='k', linestyle='-', label='sum of E-fields') ax_3.plot(freqs[passband_filter] / units.MHz, np.sqrt(flux_sum / (scipy.constants.c * (units.m / units.s)) / (scipy.constants.epsilon_0 * (units.farad / units.m))) / d_f / (units.V / units.m), c='C2', label='E-field from total flux') ax_4.plot(channel.get_frequencies() / units.MHz, np.abs(noise_spec_sum), c='k', linestyle=':', label='total noise') ax_2.legend() ax_3.legend() ax_4.legend() fig.tight_layout() ax1.plot(channel.get_times(), channel.get_trace(), label='new trace') ax1.plot(channel.get_times(), fft.freq2time(noise_spec_sum, channel.get_sampling_rate()), label='noise') ax2.plot(channel.get_frequencies() / units.MHz, np.abs(channel.get_frequency_spectrum()), label='new spectrum') ax2.plot(channel.get_frequencies() / units.MHz, np.abs(noise_spec_sum), label='noise') ax1.legend() ax2.legend() plt.show()
def get_array_of_channels(station, det, zenith, azimuth, polarization): """ Returns an array of the channel traces that is cut to the physical overlapping time Parameters ---------- station: Station Station from which to take the channels det: Detector The detector description zenith: float Arrival zenith angle at antenna azimuth: float Arrival azimuth angle at antenna polarization: int 0: eTheta 1: ePhi """ time_shifts = np.zeros(8) t_geos = np.zeros(8) sampling_rate = station.get_channel(0).get_sampling_rate() station_id = station.get_id() site = det.get_site(station_id) for iCh, channel in enumerate(station.get_electric_fields()): channel_id = channel.get_channel_ids()[0] antenna_position = det.get_relative_position(station_id, channel_id) # determine refractive index of signal propagation speed between antennas refractive_index = ice.get_refractive_index( 1, site) # if signal comes from above, in-air propagation speed if (zenith > 0.5 * np.pi): refractive_index = ice.get_refractive_index( antenna_position[2], site ) # if signal comes from below, use refractivity at antenna position refractive_index = 1.353 time_shift = -geo_utl.get_time_delay_from_direction( zenith, azimuth, antenna_position, n=refractive_index) t_geos[iCh] = time_shift time_shift += channel.get_trace_start_time() time_shifts[iCh] = time_shift delta_t = time_shifts.max() - time_shifts.min() tmin = time_shifts.min() tmax = time_shifts.max() trace_length = station.get_electric_fields()[0].get_times( )[-1] - station.get_electric_fields()[0].get_times()[0] traces = [] n_samples = None for iCh, channel in enumerate(station.get_electric_fields()): tstart = delta_t - (time_shifts[iCh] - tmin) tstop = tmax - time_shifts[iCh] - delta_t + trace_length iStart = int(round(tstart * sampling_rate)) iStop = int(round(tstop * sampling_rate)) if (n_samples is None): n_samples = iStop - iStart if (n_samples % 2): n_samples -= 1 trace = copy.copy(channel.get_trace() [polarization]) # copy to not modify data structure trace = trace[iStart:(iStart + n_samples)] base_trace = NuRadioReco.framework.base_trace.BaseTrace( ) # create base trace class to do the fft with correct normalization etc. base_trace.set_trace(trace, sampling_rate) traces.append(base_trace) return traces
def get_array_of_channels(station, use_channels, det, zenith, azimuth, antenna_pattern_provider, time_domain=False): time_shifts = np.zeros(len(use_channels)) t_cables = np.zeros(len(use_channels)) t_geos = np.zeros(len(use_channels)) station_id = station.get_id() site = det.get_site(station_id) for iCh, channel in enumerate(station.iter_channels(use_channels)): channel_id = channel.get_id() antenna_position = det.get_relative_position(station_id, channel_id) # determine refractive index of signal propagation speed between antennas refractive_index = ice.get_refractive_index( 1, site) # if signal comes from above, in-air propagation speed if station.is_cosmic_ray(): if (zenith > 0.5 * np.pi): refractive_index = ice.get_refractive_index( antenna_position[2], site ) # if signal comes from below, use refractivity at antenna position if station.is_neutrino(): refractive_index = ice.get_refractive_index( antenna_position[2], site) time_shift = -geo_utl.get_time_delay_from_direction( zenith, azimuth, antenna_position, n=refractive_index) t_geos[iCh] = time_shift t_cables[iCh] = channel.get_trace_start_time() logger.debug( "time shift channel {}: {:.2f}ns (signal prop), {:.2f}ns (trace start time)" .format(channel.get_id(), time_shift, channel.get_trace_start_time())) time_shift += channel.get_trace_start_time() time_shifts[iCh] = time_shift delta_t = time_shifts.max() - time_shifts.min() tmin = time_shifts.min() tmax = time_shifts.max() logger.debug("adding relative station time = {:.0f}ns".format( (t_cables.min() + t_geos.max()) / units.ns)) logger.debug("delta t is {:.2f}".format(delta_t / units.ns)) trace_length = station.get_channel( 0).get_times()[-1] - station.get_channel(0).get_times()[0] debug_cut = 0 if (debug_cut): fig, ax = plt.subplots(len(use_channels), 1) traces = [] n_samples = None for iCh, channel in enumerate(station.iter_channels(use_channels)): tstart = delta_t - (time_shifts[iCh] - tmin) tstop = tmax - time_shifts[iCh] - delta_t + trace_length iStart = int(round(tstart * channel.get_sampling_rate())) iStop = int(round(tstop * channel.get_sampling_rate())) if (n_samples is None): n_samples = iStop - iStart if (n_samples % 2): n_samples -= 1 trace = copy.copy( channel.get_trace()) # copy to not modify data structure trace = trace[iStart:(iStart + n_samples)] if (debug_cut): ax[iCh].plot(trace) base_trace = NuRadioReco.framework.base_trace.BaseTrace( ) # create base trace class to do the fft with correct normalization etc. base_trace.set_trace(trace, channel.get_sampling_rate()) traces.append(base_trace) times = traces[0].get_times( ) # assumes that all channels have the same sampling rate if (time_domain): # save time domain traces first to avoid extra fft V_timedomain = np.zeros((len(use_channels), len(times))) for iCh, trace in enumerate(traces): V_timedomain[iCh] = trace.get_trace() frequencies = traces[0].get_frequencies( ) # assumes that all channels have the same sampling rate V = np.zeros((len(use_channels), len(frequencies)), dtype=np.complex) for iCh, trace in enumerate(traces): V[iCh] = trace.get_frequency_spectrum() efield_antenna_factor = trace_utilities.get_efield_antenna_factor( station, frequencies, use_channels, det, zenith, azimuth, antenna_pattern_provider) if (debug_cut): plt.show() if (time_domain): return efield_antenna_factor, V, V_timedomain return efield_antenna_factor, V
def run(self, evt, station, det, use_channels=None, use_MC_direction=False, force_Polarization=''): """ run method. This function is executed for each event Parameters --------- evt station det use_channels: array of ints (default: [0, 1, 2, 3]) the channel ids to use for the electric field reconstruction use_MC_direction: bool if True uses zenith and azimuth direction from simulated station if False uses reconstructed direction from station parameters. force_Polarization: str if eTheta or ePhi, then only reconstructs chosen polarization of electric field, assuming the other is 0. Otherwise, reconstructs electric field for both eTheta and ePhi """ if use_channels is None: use_channels = [0, 1, 2, 3] station_id = station.get_id() if use_MC_direction: zenith = station.get_sim_station()[stnp.zenith] azimuth = station.get_sim_station()[stnp.azimuth] else: logger.info( "Using reconstructed (or starting) angles as no signal arrival angles are present" ) zenith = station[stnp.zenith] azimuth = station[stnp.azimuth] efield_antenna_factor, V = get_array_of_channels( station, use_channels, det, zenith, azimuth, self.antenna_provider) n_frequencies = len(V[0]) denom = (efield_antenna_factor[0][0] * efield_antenna_factor[1][1] - efield_antenna_factor[0][1] * efield_antenna_factor[1][0]) mask = np.abs(denom) != 0 # solving for electric field using just two orthorgonal antennas E1 = np.zeros_like(V[0]) E2 = np.zeros_like(V[0]) E1[mask] = (V[0] * efield_antenna_factor[1][1] - V[1] * efield_antenna_factor[0][1])[mask] / denom[mask] E2[mask] = (V[1] - efield_antenna_factor[1][0] * E1)[mask] / efield_antenna_factor[1][1][mask] denom = (efield_antenna_factor[0][0] * efield_antenna_factor[-1][1] - efield_antenna_factor[0][1] * efield_antenna_factor[-1][0]) mask = np.abs(denom) != 0 E1[mask] = (V[0] * efield_antenna_factor[-1][1] - V[-1] * efield_antenna_factor[0][1])[mask] / denom[mask] E2[mask] = (V[-1] - efield_antenna_factor[-1][0] * E1)[mask] / efield_antenna_factor[-1][1][mask] # solve it in a vectorized way efield3_f = np.zeros((2, n_frequencies), dtype=np.complex) if force_Polarization == 'eTheta': efield3_f[:1, mask] = np.moveaxis( stacked_lstsq( np.moveaxis(efield_antenna_factor[:, 0, mask], 1, 0)[:, :, np.newaxis], np.moveaxis(V[:, mask], 1, 0)), 0, 1) elif force_Polarization == 'ePhi': efield3_f[1:, mask] = np.moveaxis( stacked_lstsq( np.moveaxis(efield_antenna_factor[:, 1, mask], 1, 0)[:, :, np.newaxis], np.moveaxis(V[:, mask], 1, 0)), 0, 1) else: efield3_f[:, mask] = np.moveaxis( stacked_lstsq( np.moveaxis(efield_antenna_factor[:, :, mask], 2, 0), np.moveaxis(V[:, mask], 1, 0)), 0, 1) # add eR direction efield3_f = np.array([ np.zeros_like(efield3_f[0], dtype=np.complex), efield3_f[0], efield3_f[1] ]) electric_field = NuRadioReco.framework.electric_field.ElectricField( use_channels, [0, 0, 0]) electric_field.set_frequency_spectrum( efield3_f, station.get_channel(0).get_sampling_rate()) electric_field.set_parameter(efp.zenith, zenith) electric_field.set_parameter(efp.azimuth, azimuth) # figure out the timing of the E-field t_shifts = np.zeros(V.shape[0]) site = det.get_site(station_id) if (zenith > 0.5 * np.pi): logger.warning( "Module has not been optimized for neutrino reconstruction yet. Results may be nonsense." ) refractive_index = ice.get_refractive_index( -1, site ) # if signal comes from below, use refractivity at antenna position else: refractive_index = ice.get_refractive_index( 1, site) # if signal comes from above, in-air propagation speed for i_ch, channel_id in enumerate(use_channels): antenna_position = det.get_relative_position( station.get_id(), channel_id) t_shifts[i_ch] = station.get_channel( channel_id).get_trace_start_time( ) - geo_utl.get_time_delay_from_direction( zenith, azimuth, antenna_position, n=refractive_index) electric_field.set_trace_start_time(t_shifts.max()) station.add_electric_field(electric_field)
def run(self, evt, station, det): t = time.time() # access simulated efield and high level parameters sim_station = station.get_sim_station() if (len(sim_station.get_electric_fields()) == 0): raise LookupError(f"station {station.get_id()} has no efields") sim_station_id = sim_station.get_id() # first we determine the trace start time of all channels and correct # for different cable delays times_min = [] times_max = [] for iCh in det.get_channel_ids(sim_station_id): for electric_field in sim_station.get_electric_fields_for_channels( [iCh]): time_resolution = 1. / electric_field.get_sampling_rate() cab_delay = det.get_cable_delay(sim_station_id, iCh) t0 = electric_field.get_trace_start_time() + cab_delay # if we have a cosmic ray event, the different signal travel time to the antennas has to be taken into account if sim_station.is_cosmic_ray(): site = det.get_site(sim_station_id) antenna_position = det.get_relative_position( sim_station_id, iCh) - electric_field.get_position() if sim_station.get_parameter( stnp.zenith ) > 90 * units.deg: # signal is coming from below, so we take IOR of ice index_of_refraction = ice.get_refractive_index( antenna_position[2], site) else: # signal is coming from above, so we take IOR of air index_of_refraction = ice.get_refractive_index(1, site) # For cosmic ray events, we only have one electric field for all channels, so we have to account # for the difference in signal travel between channels. IMPORTANT: This is only accurate # if all channels have the same z coordinate travel_time_shift = geo_utl.get_time_delay_from_direction( sim_station.get_parameter(stnp.zenith), sim_station.get_parameter(stnp.azimuth), antenna_position, index_of_refraction) t0 += travel_time_shift if ( not np.isnan(t0) ): # trace start time is None if no ray tracing solution was found and channel contains only zeros times_min.append(t0) times_max.append(t0 + electric_field.get_number_of_samples() / electric_field.get_sampling_rate()) self.logger.debug( "trace start time {}, cab_delty {}, tracelength {}". format( electric_field.get_trace_start_time(), cab_delay, electric_field.get_number_of_samples() / electric_field.get_sampling_rate())) times_min = np.array(times_min) times_max = np.array(times_max) if times_min.min() < 0: times_min -= times_min.min() times_max -= times_min.min() times_min = np.array(times_min) - self.__pre_pulse_time times_max = np.array(times_max) + self.__post_pulse_time trace_length = times_max.max() - times_min.min() trace_length_samples = int(round(trace_length / time_resolution)) if trace_length_samples % 2 != 0: trace_length_samples += 1 self.logger.debug( "smallest trace start time {:.1f}, largest trace time {:.1f} -> n_samples = {:d} {:.0f}ns)" .format(times_min.min(), times_max.max(), trace_length_samples, trace_length / units.ns)) # loop over all channels for channel_id in det.get_channel_ids(station.get_id()): # one channel might contain multiple channels to store the signals from multiple ray paths, # so we loop over all simulated channels with the same id, # convolve each trace with the antenna response for the given angles # and everything up in the time domain self.logger.debug('channel id {}'.format(channel_id)) channel = NuRadioReco.framework.channel.Channel(channel_id) channel_spectrum = None trace_object = None if (self.__debug): from matplotlib import pyplot as plt fig, axes = plt.subplots(2, 1) for electric_field in sim_station.get_electric_fields_for_channels( [channel_id]): # all simulated channels have a different trace start time # in a measurement, all channels have the same physical start time # so we need to create one long trace that can hold all the different channel times # to achieve a good time resolution, we upsample the trace first. new_efield = NuRadioReco.framework.base_trace.BaseTrace( ) # create new data structure with new efield length new_efield.set_trace(copy.copy(electric_field.get_trace()), electric_field.get_sampling_rate()) new_trace = np.zeros((3, trace_length_samples)) # calculate the start bin if (not np.isnan(electric_field.get_trace_start_time())): cab_delay = det.get_cable_delay(sim_station_id, channel_id) if sim_station.is_cosmic_ray(): site = det.get_site(sim_station_id) antenna_position = det.get_relative_position( sim_station_id, channel_id) - electric_field.get_position() if sim_station.get_parameter( stnp.zenith ) > 90 * units.deg: # signal is coming from below, so we take IOR of ice index_of_refraction = ice.get_refractive_index( antenna_position[2], site) else: # signal is coming from above, so we take IOR of air index_of_refraction = ice.get_refractive_index( 1, site) travel_time_shift = geo_utl.get_time_delay_from_direction( sim_station.get_parameter(stnp.zenith), sim_station.get_parameter(stnp.azimuth), antenna_position, index_of_refraction) start_time = electric_field.get_trace_start_time( ) + cab_delay - times_min.min() + travel_time_shift start_bin = int(round(start_time / time_resolution)) time_remainder = start_time - start_bin * time_resolution else: start_time = electric_field.get_trace_start_time( ) + cab_delay - times_min.min() start_bin = int(round(start_time / time_resolution)) time_remainder = start_time - start_bin * time_resolution self.logger.debug( 'channel {}, start time {:.1f} = bin {:d}, ray solution {}' .format( channel_id, electric_field.get_trace_start_time() + cab_delay, start_bin, electric_field[efp.ray_path_type])) new_efield.apply_time_shift(time_remainder) new_trace[:, start_bin:(start_bin + new_efield.get_number_of_samples() )] = new_efield.get_trace() trace_object = NuRadioReco.framework.base_trace.BaseTrace() trace_object.set_trace(new_trace, 1. / time_resolution) trace_object.set_trace_start_time( np.min(times_min) - cab_delay) if (self.__debug): axes[0].plot(trace_object.get_times(), new_trace[1], label="eTheta {}".format( electric_field[efp.ray_path_type]), c='C0') axes[0].plot(trace_object.get_times(), new_trace[2], label="ePhi {}".format( electric_field[efp.ray_path_type]), c='C0', linestyle=':') axes[0].plot(electric_field.get_times(), electric_field.get_trace()[1], c='C1', linestyle='-', alpha=.5) axes[0].plot(electric_field.get_times(), electric_field.get_trace()[2], c='C1', linestyle=':', alpha=.5) ff = trace_object.get_frequencies() efield_fft = trace_object.get_frequency_spectrum() zenith = electric_field[efp.zenith] azimuth = electric_field[efp.azimuth] # get antenna pattern for current channel VEL = trace_utilities.get_efield_antenna_factor( sim_station, ff, [channel_id], det, zenith, azimuth, self.antenna_provider) if VEL is None: # this can happen if there is not signal path to the antenna voltage_fft = np.zeros_like( efield_fft[1]) # set voltage trace to zeros else: # Apply antenna response to electric field VEL = VEL[ 0] # we only requested the VEL for one channel, so selecting it voltage_fft = np.sum( VEL * np.array([efield_fft[1], efield_fft[2]]), axis=0) # Remove DC offset voltage_fft[np.where(ff < 5 * units.MHz)] = 0. if (self.__debug): axes[1].plot(trace_object.get_times(), fft.freq2time( voltage_fft, electric_field.get_sampling_rate()), label="{}, zen = {:.0f}deg".format( electric_field[efp.ray_path_type], zenith / units.deg)) if ('amp' in self.__uncertainty): voltage_fft *= np.random.normal( 1, self.__uncertainty['amp'][channel_id]) if ('sys_amp' in self.__uncertainty): voltage_fft *= self.__uncertainty['sys_amp'][channel_id] if (channel_spectrum is None): channel_spectrum = voltage_fft else: channel_spectrum += voltage_fft if (self.__debug): axes[0].legend(loc='upper left') axes[1].legend(loc='upper left') plt.show() if trace_object is None: # this happens if don't have any efield for this channel # set the trace to zeros channel.set_trace(np.zeros(trace_length_samples), 1. / time_resolution) else: channel.set_frequency_spectrum( channel_spectrum, trace_object.get_sampling_rate()) channel.set_trace_start_time(times_min.min()) station.add_channel(channel) self.__t += time.time() - t
def run(self, evt, station, det, channels_to_use=None, cosmic_ray=False): """ Parameters ---------------- evt: Event The event to run the module on station: Station The station to run the module on det: Detector The detector description channels_to_use: list of int (default: [0, 1, 2, 3]) List with the IDs of channels to use for reconstruction cosmic_ray: Bool (default: False) Flag to mark event as cosmic ray """ if channels_to_use is None: channels_to_use = [0, 1, 2, 3] station_id = station.get_id() times = [] times_error = [] positions = [] for iCh, efield in enumerate(station.get_electric_fields()): if (len(efield.get_channel_ids()) > 1): # FIXME: this can be changed later if each efield has a position and absolute time raise AttributeError( "found efield that is valid for more than one channel. Position can't be determined." ) channel_id = efield.get_channel_ids()[0] if (channel_id not in channels_to_use): continue times.append(efield[efp.signal_time]) if (efield.has_parameter_error(efp.signal_time)): times_error.append( (efield.get_parameter_error(efp.signal_time)**2 + self.__time_uncertainty**2)**0.5) else: times_error.append(self.__time_uncertainty) positions.append(det.get_relative_position(station_id, channel_id)) times = np.array(times) times_error = np.array(times_error) positions = np.array(positions) site = det.get_site(station_id) n_ice = ice.get_refractive_index(-0.01, site) from scipy import optimize as opt def get_expected_times(params, channel_positions): zenith, azimuth = params if cosmic_ray: if ((zenith < 0) or (zenith > 0.5 * np.pi)): return np.ones(len(channel_positions)) * np.inf else: if ((zenith < 0.5 * np.pi) or (zenith > np.pi)): return np.ones(len(channel_positions)) * np.inf v = hp.spherical_to_cartesian(zenith, azimuth) c = constants.c * units.m / units.s if not cosmic_ray: c = c / n_ice logger.debug("using speed of light = {:.4g}".format(c)) t_expected = -(np.dot(v, channel_positions.T) / c) return t_expected def obj_plane(params, pos, t_measured): t_expected = get_expected_times(params, pos) chi_squared = np.sum( ((t_expected - t_expected.mean()) - (t_measured - t_measured.mean()))**2 / times_error**2) logger.debug("texp = {texp}, tm = {tmeas}, {chi2}".format( texp=t_expected, tmeas=t_measured, chi2=chi_squared)) return chi_squared method = "Nelder-Mead" options = {'maxiter': 1000, 'disp': False} zenith_start = 135 * units.deg if cosmic_ray: zenith_start = 45 * units.deg starting_chi2 = {} for starting_az in np.array([0, 90, 180, 270]) * units.degree: starting_chi2[starting_az] = obj_plane((zenith_start, starting_az), positions, times) azimuth_start = min(starting_chi2, key=starting_chi2.get) res = opt.minimize(obj_plane, x0=[zenith_start, azimuth_start], args=(positions, times), method=method, options=options) chi2 = res.fun df = len(channels_to_use) - 3 if (df == 0): chi2ndf = chi2 chi2prob = stats.chi2.sf(chi2, 1) else: chi2ndf = chi2 / df chi2prob = stats.chi2.sf(chi2, df) output_str = "reconstucted angles theta = {:.1f}, phi = {:.1f}, chi2/ndf = {:.2g}/{:d} = {:.2g}, chi2prob = {:.3g}".format( res.x[0] / units.deg, hp.get_normalized_angle(res.x[1]) / units.deg, res.fun, df, chi2ndf, chi2prob) logger.info(output_str) station[stnp.zenith] = res.x[0] station[stnp.azimuth] = hp.get_normalized_angle(res.x[1]) station[stnp.chi2_efield_time_direction_fit] = chi2 station[stnp.ndf_efield_time_direction_fit] = df if (cosmic_ray): station[stnp.cr_zenith] = res.x[0] station[stnp.cr_azimuth] = hp.get_normalized_angle(res.x[1]) if (self.__debug): # calculate residuals t_exp = get_expected_times(res.x, positions) from matplotlib import pyplot as plt fig, ax = plt.subplots(1, 1) ax.errorbar(channels_to_use, ((times - times.mean()) - (t_exp - t_exp.mean())) / units.ns, fmt='o', yerr=times_error / units.ns) ax.set_xlabel("channel id") ax.set_ylabel(r"$t_\mathrm{meas} - t_\mathrm{exp}$ [ns]") pass
def run(self, evt, station, det): t = time.time() # access simulated efield and high level parameters sim_station = station.get_sim_station() if(len(sim_station.get_electric_fields()) == 0): raise LookupError(f"station {station.get_id()} has no efields") for channel_id in det.get_channel_ids(station.get_id()): # one channel might contain multiple channels to store the signals from multiple ray paths and showers, # so we loop over all simulated channels with the same id, self.logger.debug('channel id {}'.format(channel_id)) for electric_field in sim_station.get_electric_fields_for_channels([channel_id]): sim_channel = NuRadioReco.framework.sim_channel.SimChannel(channel_id, shower_id=electric_field.get_shower_id(), ray_tracing_id=electric_field.get_ray_tracing_solution_id()) ff = electric_field.get_frequencies() efield_fft = electric_field.get_frequency_spectrum() zenith = electric_field[efp.zenith] azimuth = electric_field[efp.azimuth] # get antenna pattern for current channel VEL = trace_utilities.get_efield_antenna_factor(sim_station, ff, [channel_id], det, zenith, azimuth, self.antenna_provider) if VEL is None: # this can happen if there is not signal path to the antenna voltage_fft = np.zeros_like(efield_fft[1]) # set voltage trace to zeros else: # Apply antenna response to electric field VEL = VEL[0] # we only requested the VEL for one channel, so selecting it voltage_fft = np.sum(VEL * np.array([efield_fft[1], efield_fft[2]]), axis=0) # Remove DC offset voltage_fft[np.where(ff < 5 * units.MHz)] = 0. if sim_station.is_cosmic_ray(): site = det.get_site(station.get_id()) antenna_position = det.get_relative_position(station.get_id(), channel_id) - electric_field.get_position() if zenith > 90 * units.deg: # signal is coming from below, so we take IOR of ice index_of_refraction = ice.get_refractive_index(antenna_position[2], site) else: # signal is coming from above, so we take IOR of air index_of_refraction = ice.get_refractive_index(1, site) # For cosmic ray events, we only have one electric field for all channels, so we have to account # for the difference in signal travel between channels. IMPORTANT: This is only accurate # if all channels have the same z coordinate travel_time_shift = geometryUtilities.get_time_delay_from_direction( zenith, azimuth, antenna_position, index_of_refraction ) else: travel_time_shift = 0 # set the trace to zeros sim_channel.set_frequency_spectrum(voltage_fft, electric_field.get_sampling_rate()) sim_channel.set_trace_start_time(electric_field.get_trace_start_time() + travel_time_shift) sim_station.add_channel(sim_channel) self.__t += time.time() - t