wimpsim-reader C++ API Reference

struct BookKeeper
#include <I3WimpSimReader.h>

a Bookkepping class to encapsulate different cases

Public Functions

inline BookKeeper()

constructor

Public Members

double nu_

total/collected neutrino weight/events in the read file

double lep_

total/collected lepton weight/events in the read file

double e_

total/collected electron weight/events in the read file

double mu_

total/collected muon weight/events in the read file

double tau_

total/collected tau weight/events in the read file

double nc_

total/collected neutral current weight/events in the read file

double had_

total/collected hadronic weight/events in the read file

double lowe_

total/collected low-energetic weight/events in the read file

struct GenBox

a struct to contain a about the GeneratedBox and the ‘Generated Volume’.

Public Members

I3Position particle_pos

position of particle in real IceCube coordinate system

I3Position box_pos

relative position of particle inside the box eigensystem

I3Direction box_vector

direction of the long-side of the box

double box_length

length of the long-side of the box

double box_horizontal

length of the horizontal-side of the box

double box_vertical

length of the vertical-side of the box

double vgen

generated volume enclosed by the box: V_eff_mu

double aproj

projected faceside area for the neutrino: A_eff_nu

class I3WimpSimReader : public I3Module
#include <I3WimpSimReader.h>

A Generator Module for that reads WimpSimFiles

Public Functions

I3WimpSimReader(const I3Context &context)

Constructor with arguments passed.

~I3WimpSimReader()

std Destructor

void Configure()

std Configure (configure from context, check correct params)

void Process()

std Process (the action goes down here)

void Finish()

std Finish (close files, report)

Protected Types

enum exit_codes

exit codes

Values:

enumerator ERROR
enumerator SUCCESS
enumerator FAILED_ENERGYCUTS
enumerator FAILED_LEPTONCUTS
enumerator FAILED_TIMECUTS
enumerator FAILED_ZENITHCUTS
enumerator END_EVENTSTREAM
enumerator ERROR_EVENTSTREAM

Protected Attributes

std::string randomServiceName_

OPTION: the name of the random-number generator.

std::vector<std::string> filenames_

OPTION: the filename mask with wildcards and stuff.

uint over_

OPTION: N oversampling.

uint number_events_

OPTION: number of events to process.

std::vector<double> position_limits_

OPTION: list of position limits.

double radius_

OPTION: radius of effVol.

double zenith_min_

OPTION: optional lower zenith cut.

double zenith_max_

OPTION: optional upper zenith cut.

bool electrons_

OPTION: Read and simulate electron vertices.

bool muons_

OPTION: Read and simulate muon vertices.

bool taus_

OPTION: Read and simulate tau vertices.

bool ncs_

OPTION: Read and simulate tau vertices.

double sensHeight_

OPTION: Muon box activated height.

double sensRadius_

OPTION: Muon box activated Radius.

double startmjd_

OPTION: MJD from which to start simulating.

double endmjd_

OPTION: MJD towards which to cease simulating.

bool writeInfoOpt_

OPTION: write Info about weights into an IceTray.Info frame.

std::string infoFileName_

OPTION: Write Information to this text-file.

bool writeDrivingTimeOpt_

OPTION: write DrivingTime to the frame (needed for correct processing 2011)

bool pFrameOpt_

OPTION: write P-Frames instead of Q-Frames, this is for backwards-compatibility.

bool flatMapOpt_

OPTION: write a flat I3StringDoubleMap for the WIMP.

double fixGenLength_

OPTION: fall back value if a dynamic length can be calculated.

int file_index_

the current index in file_vector_

boost::iostreams::filtering_istream wimpfile_

the currently opened file

uint issued_files_

number of files that have been successfully read

uint issued_events_

number of events that have been generated

uint processed_events_

number of processed events

BookKeeper col_weight_

collected weights

BookKeeper col_events_

collected events

BookKeeper tot_weight_

total weights

BookKeeper tot_events_

total events

Private Functions

SET_LOGGER ("I3WimpSimReader")

turn your logging on: ‘I3WimpSimReader

bool OpenNextFile()

Move through the vector of file names and open the next file in line.

Returns:

true upon success, false if not

bool ReadFileHeader(boost::shared_ptr<WimpSimTools::WimpSimBaseHeader> &WimpHeaderPtr)

Reads Headers from the WimpSim-file.

Parameters:

WimpHeaderPtr – pointer to the object to which to write the FileHeader

Returns:

true if header was successfully read, false if not

I3Position RandomPosition() const

Calculate a position inside position_limits (with ranomization)

Returns:

the randomized position.

double VolumeBoxOrCylinder() const

Calculate the constant Injection Volume as a Box or a Cylinder.

Returns:

a static return of the volume of either a box or a cylinder

double RandomPositionVolume() const

Calculate the Volume of the Injection Volume.

Returns:

the volume in [m3].

GenBox GenerateBox(const I3Particle &neutrino, const I3Particle &lepton) const

Calculate sensitive box around muon track and effective Volume.

Parameters:
  • neutrino – the neutrino to calculate the effective area for

  • lepton – the lepton to calculate the muonbox for

Returns:

MuonBox containing the calculated variables

double Length_mu(const double energy) const

Calculate Length of mu track in ice from second order parametrization.

Parameters:

energy – : muon energy

Returns:

length of muon track (in meter)

double Length_elmag(const double energy) const

Calculate Length of electron cascade in ice.

Parameters:

energy – : electron energy

Returns:

length of electron cascade (in meter)

double Length_tau(const double energy) const

Calculate Length of tau track in ice.

Parameters:

energy – : tau energy

Returns:

length of tau track (contained) (in meter)

double Length_had(const double energy) const

Calculate Length of hadronic cascades in ice.

Parameters:

energy – : hadronic energy

Returns:

length of hadronic cascades (contained) (in meter)

I3FramePtr WriteInfoFrame(const boost::shared_ptr<WimpSimTools::WimpSimBaseHeader> WimpHeaderPtr) const

Writes Info-Frames containing basic steering info from the .dat file.

Parameters:

WimpHeaderPtr – pointer to the container with the information from the file header

Returns:

pointer to the heads info frame

I3FramePtr WriteFinishFrame() const

A function that report on the processed events and weights collected and writes them to a Info-frame and pushes it to the Outbox.

Returns:

pointer to the last info frame

void Report() const

A function that report on the processed events and weights collected.

double RandomMJD() const

Compute a random MJD between starmjd_ and endmjd_.

Returns:

MJD randomized

I3Time RandomMJDTime() const

Compute a random MJD between starmjd_ and endmjd_.

Returns:

MJD randomized

I3Position RotateZenAzi(const I3Position &pos, const I3Direction &dir) const

rotate the box which was previously defined in the Cartesian-coordinate eigenvectors to the eigenvectors of the particle

Parameters:
  • pos – : the position

  • dir – : the rotation expressed in 3 space

Returns:

position in IceCube

I3FramePtr WriteEvent(const boost::shared_ptr<WimpSimTools::WimpSimBaseHeader> WimpHeaderPtr, const boost::shared_ptr<WimpSimTools::WimpSimBaseEvent> WimpEventPtr, const double mjd) const

write the event to a frame

Parameters:
  • WimpHeaderPtr – pointer to the FileHeader Information Container

  • WimpEventPtr – pointer to the Event Information container

  • mjd – time which to write to I3EventHeader and DrivingTime

Returns:

the frame to which the WimpEvent has been written

exit_codes ReadSunEvent(boost::shared_ptr<WimpSimTools::WimpSimSunEvent> &SunEventPtr)

write the Sun type event to the container

Parameters:

SunEventPtr – pointer to the object to which to write the event

Returns:

WimpSimReader internal exit_codes

exit_codes ReadEarthEvent(boost::shared_ptr<WimpSimTools::WimpSimEarthEvent> &EarthEventPtr)

write the Earth type event to the container

Parameters:

EarthEventPtr – pointer to the object to which to write the event

Returns:

WimpSimReader internal exit_codes

exit_codes CutsNWeights(const boost::shared_ptr<WimpSimTools::WimpSimBaseEvent> WimpEventPtr, const double mjd)

Evaluate Cuts and Weights of this event.

Parameters:
  • WimpEventPtr – pointer to the object from which the decision and values should be taken

  • mjd – the modjulian day that should be used

Returns:

WimpSimReader internal exit_codes

std::string HeaderString(const boost::shared_ptr<WimpSimTools::WimpSimBaseHeader> WimpHeaderPtr) const

write file header info to a string

std::string ReportString() const

write module state information

void WriteInfoToFile(const std::string &filename) const

write all info to a file

Private Members

I3RandomServicePtr randomService_

The RandomService installed in the tray: can also be passed as an OPTION.

class WimpSimBaseEvent
#include <WimpSimTools.h>

encapsule all base information WimpSim provides for every Events

Subclassed by WimpSimTools::WimpSimEarthEvent, WimpSimTools::WimpSimSunEvent

Public Functions

virtual ~WimpSimBaseEvent()
virtual double GetIceCubeNu_Azimuth() const = 0

nadir of primary neutrino (in deg)

Implementation to Get Azimuth of the Primary Neutrino

virtual double GetIceCubeNu_Zenith() const = 0

Implementation to Get Zenith of the Primary Neutrino.

virtual double GetIceCubeLep_Azimuth() const = 0

Implementation to Get Azimuth of the Secondary Lepton.

virtual double GetIceCubeLep_Zenith() const = 0

Implementation to Get Zenith of the Secondary Lepton.

virtual double GetIceCubeHad_Azimuth() const = 0

Implementation to Get Azimuth of the Secondary Hadron.

virtual double GetIceCubeHad_Zenith() const = 0

Implementation to Get Zenith of the Secondary Hadron.

virtual double GetIceCubeOrigin_Azimuth() const = 0

Implementation to Get Zenith of the Origin of this event.

virtual double GetIceCubeOrigin_Zenith() const = 0

Implementation to Get Zenith of the Origin of this event.

virtual double GetMJD() const = 0

Implementation to Get the MJD time of this event.

Public Members

uint eventnbr_

number of the event

std::string nu_type_

type of the primary neutrino (text)

double nu_energy_

total energy of the primary neutrino (in GeV)

double nu_weight_

angle of the primary neutrino (in deg)

std::string lep_type_

type of the secondary lepton inice (text)

double lep_energy_

total energy of the secondary lepton inice (in GeV)

double lep_weight_

weight of the secondary lepton inice

std::string had_type_

type of the secondary hadron inice (text)

double had_energy_

total energy of the secondary hadron inice (in GeV)

double had_weight_

angle of the secondary hadron inice (in deg)

class WimpSimBaseHeader
#include <WimpSimTools.h>

encapsulate all information WimpSim provides in the FileHeader

Public Members

double mass_

WIMP mass.

uint channel_

Annihilation channel.

std::string source_

Source of the WIMPs (default is Sun, Earth)

uint number_content_

Number of annihilations simulated this file.

uint number_total_

Number of annihilations simulated total.

uint files_index_

index of this file

uint files_total_

Number of files this run.

int seed_

initial seed

double theta_12_

weak mixing angle 12

double theta_13_

weak mixing angle 13

double theta_23_

weak mixing angle 23

double delta_

delta

double delta_m2_21_

neutrino mass squared difference 21

double delta_m2_31_

neutrino mass squared difference 31

double latitude_

latitude of the detector position in [deg]

double longitude_

longitude of the detector position in [deg]

double mjd1_

start MJD of the simulation

double mjd2_

end MJD of the simulation

class WimpSimEarthEvent : public WimpSimTools::WimpSimBaseEvent
#include <WimpSimTools.h>

encapsulate all information WimpSim Provides for Earth Events

Public Functions

virtual ~WimpSimEarthEvent()
virtual double GetIceCubeNu_Azimuth() const

Get Azimuth of the Primary Neutrino.

virtual double GetIceCubeNu_Zenith() const

Get Zenith of the Primary Neutrino.

virtual double GetIceCubeLep_Azimuth() const

Get Azimuth of the Secondary Lepton.

virtual double GetIceCubeLep_Zenith() const

Get Zenith of the Secondary Lepton.

virtual double GetIceCubeHad_Azimuth() const

Get Azimuth of the Secondary Hadron.

virtual double GetIceCubeHad_Zenith() const

Get Zenith of the Secondary Hadron.

virtual double GetIceCubeOrigin_Azimuth() const

Get Azimuth of the Origin of this event.

virtual double GetIceCubeOrigin_Zenith() const

Get Zenith of the Origin of this event.

virtual double GetMJD() const

Get the MJD time of this event.

Public Members

double nu_angle_

angle of the primary neutrino (in deg)

double lep_angle_

angle of the secondary lepton inice (in deg)

double had_angle_

angle of the secondary hadron inice (in deg)

double nu_nadir_

weight of the secondary hadron inice

double nu_az_

azimuth of primary neutrino (in deg)

double lep_nadir_

nadir of secondary lepton inice (in deg)

double lep_az_

azimuth of secondary lepton inice (in deg)

double had_nadir_

nadir of secondary hadron inice (in deg)

double had_az_

azimuth of secondary hadron inice (in deg)

double nu_el_

elevation of primary neutrino (in deg)

double lep_el_

elevation of secondary lepton inice (in deg)

double had_el_

elevation of secondary hadron inice (in deg)

class WimpSimSunEvent : public WimpSimTools::WimpSimBaseEvent
#include <WimpSimTools.h>

encapsulate all information WimpSim Provides for Sun Events

Public Functions

virtual ~WimpSimSunEvent()
virtual double GetIceCubeNu_Azimuth() const

Get Azimuth of the Primary Neutrino.

virtual double GetIceCubeNu_Zenith() const

Get Zenith of the Primary Neutrino.

virtual double GetIceCubeLep_Azimuth() const

Get Azimuth of the Secondary Lepton.

virtual double GetIceCubeLep_Zenith() const

Get Zenith of the Secondary Lepton.

virtual double GetIceCubeHad_Azimuth() const

Get Azimuth of the Secondary Hadron.

virtual double GetIceCubeHad_Zenith() const

Get Zenith of the Secondary Hadron.

virtual double GetIceCubeOrigin_Azimuth() const

Get Azimuth of the Origin of this event.

virtual double GetIceCubeOrigin_Zenith() const

Get Zenith of the Origin of this event.

virtual double GetMJD() const

Get the MJD time of this event.

Public Members

double mjd_

mjd of the event

double sun_nadir_

nadir of sun (in deg)

double sun_ra_

sun right asscention (in deg)

double sun_dec_

declination of the sun (in deg)

double sun_az_

azimuth of the sun (in deg)

double sun_el_

elevation (zenith) of the sun (in deg)

double nu_angle_

angle of the primary neutrino (in deg)

double lep_angle_

weight of the primary neutrino

double had_angle_

angle of the secondary hadron inice (in deg)

double nu_nadir_

nadir of primary neutrino (in deg)

double nu_az_

azimuth of primary neutrino (in deg)

double lep_nadir_

nadir of secondary lepton inice (in deg)

double lep_az_

azimuth of secondary lepton inice (in deg)

double had_nadir_

nadir of secondary hadron inice (in deg)

double had_az_

azimuth of secondary hadron inice (in deg)

double nu_el_

elevation of primary neutrino (in deg)

double lep_el_

elevation of secondary lepton inice (in deg)

double had_el_

elevation of secondary hadron inice (in deg)

namespace std

STL namespace.

namespace WimpSim
namespace WimpSimTools

Functions

double IceCubeAziRadFromWimpSimAziRad(const double wimpsim_azi)

Convert the azimuth given by WimpSim to IceCube Definitions of Azimuth [in rad].

Parameters:

wimpsim_azi – : from WimpSim [0…2*pi]

Returns:

Azimuth in IC coordinates

double IceCubeZenRadFromWimpSimElRad(const double wimpsim_el)

Convert the elevation given by WimpSim to IceCube Definitions of Zenith [in rad].

Parameters:

wimpsim_el – : from WimpSim [0…pi]

Returns:

Zenith in IC coordinates

double IceCubeZenRadFromWimpSimNadirRad(const double wimpsim_nadir)

Convert the Nadir given by WimpSim to IceCube Definitions of Zenith [in rad].

Parameters:

wimpsim_nadir – : from WimpSim [0…pi]

Returns:

Zenith in IC coordinates

double GetParticleMass(const I3Particle::ParticleType partType)

Returns the mass of a particle.

Parameters:

partType – type of that particle

Returns:

the mass of that particle

double TotEnergyToKinEnergy(const double tot_energy, const I3Particle::ParticleType partType)

Convert the energy from total energy for this particle type to the kinetic energy.

Parameters:
  • tot_energy – total energy of that particle

  • partType – type of that particle

Returns:

kinematic energy for that particle

file I3WimpSimReader.cxx
#include “icetray/I3Units.h”
#include “dataclasses/I3Constants.h”
#include “icetray/I3Frame.h”
#include “icetray/I3Module.h”
#include “phys-services/I3RandomService.h”
#include “icetray/I3Bool.h”
#include “icetray/I3Int.h”
#include “dataclasses/I3Double.h”
#include “dataclasses/I3String.h”
#include “dataclasses/I3Time.h”
#include “dataclasses/I3Map.h”
#include “dataclasses/I3Position.h”
#include “dataclasses/I3Direction.h”
#include “dataclasses/physics/I3EventHeader.h”
#include “dataclasses/physics/I3Particle.h”
#include “dataclasses/physics/I3MCTreeUtils.h”
#include <boost/assign/list_of.hpp>
#include <boost/assign/std/vector.hpp>
#include <boost/assert.hpp>
#include <boost/algorithm/string.hpp>
#include <boost/iostreams/filter/gzip.hpp>
#include <boost/iostreams/device/file.hpp>
#include <list>
#include <fstream>
#include <sstream>
#include <math.h>
#include <iostream>

Reads WimpSim data for Sun and Earth Events of the New WimpSim-v3.01 format This Module takes a WIMPSim Files, reads them event by event and puts Varibales into the frame: I3EventHeader, I3MCTree, WIMP_params.

copyright (C) 2012 the IceCube Collaboration

Rcs

Date

Rcs

2012-11-18

Author

mzoll marcel.zoll@fysik.su.se

file I3WimpSimReader.h
#include “simclasses/I3WimpParams.h”
#include <string>
#include <vector>
#include <fstream>
#include “icetray/I3Module.h”
#include “icetray/I3Frame.h”
#include “dataclasses/physics/I3Particle.h”
#include “dataclasses/I3Position.h”
#include “dataclasses/I3Direction.h”
#include “dataclasses/I3Time.h”
#include “phys-services/I3RandomService.h”
#include “boost/shared_ptr.hpp”
#include <boost/iostreams/filtering_stream.hpp>

copyright (C) 2012 the IceCube Collaboration

Rcs

Date

Rcs

2012-12-20

Author

mzoll marcel.zoll@fysik.su.se

Functions

I3_MODULE(I3WimpSimReader)
file index.dox
file WimpSimTools.cxx
#include <math.h>
#include “dataclasses/I3Constants.h”

copyright (C) 2012 the icecube collaboration

Version

Rcs

Date

Rcs

2012-11-18

Author

mzoll marcel.zoll@fysik.su.se this is a library that provides definitions and functions important to read WimpSim

file WimpSimTools.h
#include <string>
#include “dataclasses/physics/I3Particle.h”
#include “simclasses/I3WimpParams.h”

copyright (C) 2012 the icecube collaboration

Version

Rcs

WimpSimTools.h xxxxx 2012-10-10 17:01:33Z mzoll

Date

Rcs

2012-12-20

Author

mzoll marcel.zoll@fysik.su.se this is a library that provides definitions and functions important to read WimpSim

dir icetray
dir private
dir public
dir wimpsim-reader
dir wimpsim-reader
dir wimpsim-reader
page index

Author

mzoll mzoll@fysik.su.se

General Procedure

The way how wimpsim-reader works is the following way and order:

  • Configure()

    |- Configure itself from the parameters passed to the Module.

    |- Check these parameters for compatibility

  • Process() (in a continues loop until a frame is delivered)

    |- Check if a WimpSim-file is open, if not open a new one

    |- Check if a newly opened file has the correct EventHeader information, and store that information

    |- Push one Info-frame in front with that information

    |- Determine and call routine: ReadSunEvent() or ReadEarthEvent() and etore read event data to WimpEventContainer

    |- Apply preliminary Cuts and sum up weights (bookkeeping) in CutNWeights()

    |- Evaluate Exit-Conditions |- If everyting went fine, call the internal Algorithms to create I3Objects and compute Volumes: WriteEventFrame()

  • ReadSunEvent() ReadEarthEvent()

    |- The file can contain events in the format ‘EV I O O EE’ and ‘EV I EE’ where latter is the case for low energetic events which do not have daughters

    |- Read the event information linewise and store information to the WimpEventContainer

    |- Hand back to Process()

    -CutsNWeights()

    |- Convert variable units to IC units/conventions if necessary (zenith, azimuth, particletypes)

    |- Store all weights in ‘Total weight map’

    |- Apply pre-cuts (zenith, timing) and hand back to Process()

    |- Store all weights in ‘Delivered weight map’

    |- Hand back to Process()

  • WriteEventFrame()

    |- Prepare and store information to I3Objects

    |- Put according information into I3Objects

    |- Randomize vertex position according to Box-position, Injection Radius

    |- If this interaction (CC: leptype, NC) is selected:

    |- Use GenerateBox() to inject particle anew into sensitive Volume and compute Vgen

    |- Put all I3Objects to the frame

    |- Hand back to Process()

  • Finish()

    |- Create a new Info-frame

    |- Put “Collected weight map” and “Total weight map” to frame

    |- Push Info-frame

Objects in WimpSim files

Events are read by the following format from file

Event format for Sun events:

EV number MJD Sun-El (deg) Sun-Az (deg) Sun-Ra (deg) sun-dec (deg)

Tag particle-name energy (GeV) Sun-angle (deg) weight Az (deg) El (deg)

EE (end of event)

Event format for Earth Events:

EV number

Tag particle-name energy (GeV) nadir-angle (deg) azi-angle (deg) weight

EE (end of event)

There are two different formats of events:

Complete Events:

EV, I, O, O, EE

(EventHeader, Primary, Daughter, Daughter, EndEvent)

Incomplete Events due to too low energy:

EV, I, EE

(EventHeader, Primary, EndEvent)