NuRadioReco.modules.io.coreas.coreas module
- NuRadioReco.modules.io.coreas.coreas.make_sim_shower(*args, **kwargs)[source]
DEPRECATED: This function has been moved to create_sim_shower_from_hdf5(), however its functionality has been modified heavily. You will probably want to move to create_sim_shower() instead, which has an easier interface. Please refer to the documentation of create_sim_shower() for more information.
- NuRadioReco.modules.io.coreas.coreas.make_sim_station(*args, **kwargs)[source]
DEPRECATED: This function has been moved to create_sim_station(), however its functionality has been modified. Please refer to the documentation of create_sim_station() for more information.
- NuRadioReco.modules.io.coreas.coreas.get_angles(corsika, declination)[source]
Converting angles in corsika coordinates to local coordinates.
Corsika positive x-axis points to the magnetic north, NRR coordinates positive x-axis points to the geographic east. Corsika positive y-axis points to the west, NRR coordinates positive y-axis points to the geographic north. Corsika z-axis points upwards, NuRadio z-axis points upwards.
Corsika’s zenith angle of a particle trajectory is defined between the particle momentum vector and the negative z-axis, meaning that the particle is described in the direction where it is going to. The azimuthal angle is described between the positive x-axis and the horizontal component of the particle momentum vector (i.e. with respect to the magnetic north) proceeding counterclockwise.
NRR describes the particle zenith and azimuthal angle in the direction where the particle is coming from. Therefore, the zenith angle is the same, but the azimuthal angle has to be shifted by 180 + 90 degrees. The north has to be shifted by 90 degrees plus difference between geomagnetic and magnetic north.
- Parameters:
- corsikahdf5 file object
the open hdf5 file object of the corsika hdf5 file
- declinationfloat
declination of the magnetic field, in internal units
- Returns:
- zenithfloat
zenith angle
- azimuthfloat
azimuth angle
- magnetic_field_vectornp.ndarray
magnetic field vector
Examples
The declinations can be obtained using the functions in the radiotools helper package, if you have the magnetic field for the site you are interested in.
>>> magnet = hp.get_magnetic_field_vector('mooresbay') >>> dec = hp.get_declination(magnet) >>> evt = h5py.File('NuRadioReco/examples/example_data/example_data.hdf5', 'r') >>> get_angles(corsika, dec)[2] / units.gauss array([ 0.05646405, -0.08733734, 0.614 ]) >>> magnet array([ 0.058457, -0.09042 , 0.61439 ])
- NuRadioReco.modules.io.coreas.coreas.get_geomagnetic_angle(zenith, azimuth, magnetic_field_vector)[source]
Calculates the angle between the geomagnetic field and the shower axis defined by zenith and azimuth.
- Parameters:
- zenithfloat
zenith angle (in internal units)
- azimuthfloat
azimuth angle (in internal units)
- magnetic_field_vectornp.ndarray
The magnetic field vector in the NRR coordinate system (x points East, y points North, z points up)
- Returns:
- geomagnetic_anglefloat
geomagnetic angle
- NuRadioReco.modules.io.coreas.coreas.convert_obs_to_nuradio_efield(observer, zenith, azimuth, magnetic_field_vector)[source]
Converts the electric field from one CoREAS observer to NuRadio units and the on-sky coordinate system.
The on-sky CS in NRR has basis vectors eR, eTheta, ePhi. Before going to the on-sky CS, we account for the magnetic field which does not point strictly North. To get the zenith, azimuth and magnetic field vector, one can use get_angles(). The observer array should have the shape (n_samples, 4) with the columns (time, Ey, -Ex, Ez), where (x, y, z) is the NuRadio CS.
- Parameters:
- observernp.ndarray
The observer as in the HDF5 file, e.g. list(corsika[‘CoREAS’][‘observers’].values())[i].
- zenithfloat
zenith angle (in internal units)
- azimuthfloat
azimuth angle (in internal units)
- magnetic_field_vectornp.ndarray
magnetic field vector
- Returns:
- efield: np.array (3, n_samples)
Electric field in the on-sky CS (r, theta, phi)
- efield_times: np.array (n_samples)
The time values corresponding to the electric field samples
- NuRadioReco.modules.io.coreas.coreas.convert_obs_positions_to_nuradio_on_ground(observer_pos, declination=0)[source]
Convert observer positions from the CORSIKA CS to the NRR ground CS.
First, the observer position is converted to the NRR coordinate conventions (i.e. x pointing East, y pointing North, z pointing up). Then, the observer position is corrected for magnetic north (as CORSIKA only has two components to its magnetic field vector) and put in the geographic CS. To get the zenith, azimuth and magnetic field vector, one can use the get_angles() function. If multiple observers are to be converted, the observer array should have the shape (n_observers, 3).
- Parameters:
- observer_posnp.ndarray
The observer’s position as extracted from the HDF5 file, e.g. corsika[‘CoREAS’][‘my_observer’].attrs[‘position’]
- declinationfloat (default: 0)
Declination of the magnetic field.
- Returns:
- obs_positions_geo: np.ndarray
observer positions in geographic coordinates, shaped as (n_observers, 3).
- NuRadioReco.modules.io.coreas.coreas.read_CORSIKA7(input_file, declination=None, site=None)[source]
This function reads in a CORSIKA/CoREAS HDF5 file and returns an Event object with all relevant information from the file. This Event object can then be used to create an interpolator, for example.
The structure of the Event is as follows. It contains one station with ID equal to 0, which hosts a SimStation. The SimStation stores the simulated ElectricField objects in on-sky coordinates (ie eR, eTheta, ePhi). They are equipped with a position attribute that contains the position of the simulated antenna (in the NRR ground coordinate system) and are associated to a channel with an ID equal to the index of the channel as it was read from the HDF5 file.
Next to the (Sim)Station, the Event also contains a SimShower object, which stores the CORSIKA input parameters. For a list of stored parameters, see the create_sim_shower_from_hdf5() function.
Note that the function assumes the energy has been fixed to a single value, as is typical with a CoREAS simulation.
- Parameters:
- input_file: str
Path to the CORSIKA HDF5 file
- declination: float, default=0
The declination to use for the magnetic field, in internal units
- site: str, default=None
Instead of declination, a site name can be given to retrieve the declination using the magnetic field as obtained from the radiotools.helper.get_magnetic_field_vector() function
- Returns:
- evt: NuRadioReco.framework.event.Event
Event object containing the CORSIKA information
- NuRadioReco.modules.io.coreas.coreas.create_sim_shower_from_hdf5(corsika, declination=0)[source]
Creates an NuRadioReco RadioShower from a CoREAS HDF5 file, which contains the simulation inputs shower parameters. These include
the primary particle type
the observation level
the zenith and azimuth angles
the magnetic field vector
the energy of the primary particle
The following parameters are retrieved from the REAS file:
the core position
the depth of the shower maximum (in g/cm2, converted to internal units)
the distance of the shower maximum (in cm, converted to internal units)
the refractive index at ground level
the declination of the magnetic field
Lastly, these parameters are also saved IF they are available in the HDF5 file:
the atmospheric model used for the simulation
the electromagnetic energy of the shower (only present in high-level quantities are present)
This function is used in the read_CORSIKA7() function to create the SimShower object. In order to copy a SimShower object from an Event object, use the create_sim_shower() method.
- Parameters:
- corsikahdf5 file object
the open hdf5 file object of the corsika hdf5 file
- declinationfloat
- Returns:
- sim_shower: RadioShower
simulated shower object
- NuRadioReco.modules.io.coreas.coreas.create_sim_shower(evt, core_shift=None)[source]
Create an NuRadioReco SimShower from an Event object created with e.g. read_CORSIKA7(), and apply a core shift if desired. If no core shift is given, the returned SimShower will have the same core as used in the CoREAS simulation.
- Parameters:
- evt: Event
event object containing the CoREAS output, e.g. created with read_CORSIKA7()
- core_shift: np.ndarray, default=None
The core shift to apply to the core position, in internal units. Must be 3D array.
- Returns:
- sim_shower: RadioShower
simulated shower object
- NuRadioReco.modules.io.coreas.coreas.create_sim_station(station_id, evt, weight=None)[source]
Creates an NuRadioReco SimStation with the information from an Event object created with e.g. read_CORSIKA7().
Optionally, station can be assigned a weight. Note that the station is empty. To add an electric field the function add_electric_field_to_sim_station() has to be used.
- Parameters:
- station_idint
The id to assign to the new station
- evtEvent
event object containing the CoREAS output
- weightfloat
weight corresponds to area covered by station
- Returns:
- sim_station: SimStation
simulated station object
- NuRadioReco.modules.io.coreas.coreas.add_electric_field_to_sim_station(sim_station, channel_ids, efield, efield_start_time, zenith, azimuth, sampling_rate, efield_position=None)[source]
Adds an electric field trace to an existing SimStation, with the provided attributes.
All variables should be provided in internal units.
- Parameters:
- sim_stationSimStation
The simulated station object, e.g. from make_empty_sim_station()
- channel_idslist of int
Channel IDs to associate to the electric field
- efieldnp.ndarray
Electric field trace, shaped as (3, n_samples)
- efield_start_timefloat
Start time of the electric field trace
- zenithfloat
Zenith angle
- azimuthfloat
Azimuth angle
- sampling_ratefloat
Sampling rate of the trace
- efield_positionnp.ndarray or list of float
Position to associate to the electric field
- NuRadioReco.modules.io.coreas.coreas.calculate_simulation_weights(positions, zenith, azimuth, site='summit', debug=False)[source]
Calculate weights according to the area that one observer position in a starshape pattern represents. Weights are therefore given in units of area. Note: The volume of a 2d convex hull is the area.
- Parameters:
- positionslist
station position with [x, y, z] on ground
- zenithfloat
zenith angle of the shower
- azimuthfloat
azimuth angle of the shower
- sitestr
site of the simulation, default is ‘summit’
- debugbool
if true, plots are created for debugging
- Returns:
- weightsnp.array
weights of the observer position
- NuRadioReco.modules.io.coreas.coreas.set_fluence_of_efields(function, sim_station, quantity=<electricFieldParameters.signal_energy_fluence: 4>)[source]
This helper function is used to set the fluence quantity of all electric fields in a SimStation. Use this to calculate the fluences to use for interpolation.
One option to use as function is trace_utilities.get_electric_field_energy_fluence().
- Parameters:
- function: callable
The function to apply to the traces in order to calculate the fluence. Should take in a (3, n_samples) shaped array and return a float (or an array with 3 elements if you want the fluence per polarisation).
- sim_station: SimStation
The simulated station object
- quantity: electric field parameter, default=efp.signal_energy_fluence
The parameter where to store the result of the fluence calculation