Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Classes | Enumerations | Functions | Variables
JAANET Namespace Reference

Extensions to AAnet data format. More...

Classes

struct  JEvtEvaluator
 Auxiliary class to determine value of Evt objects. More...
 
struct  start_run
 Start of run record. More...
 
struct  String
 General purpose string class. More...
 
struct  detector
 Detector file. More...
 
struct  muon_desc_file
 Muon descriptor file. More...
 
struct  target
 Target. More...
 
struct  XSecFile
 Neutrino cross section file. More...
 
struct  drawing
 Drawing. More...
 
struct  cut
 General purpose class of phase space generation. More...
 
struct  cut_primary
 Phase space of primary particle. More...
 
struct  cut_seamuon
 Phase space of atmospheric muon generation. More...
 
struct  cut_in
 Phase space of incoming particle. More...
 
struct  cut_nu
 Phase space of incident neutrino. More...
 
struct  generator
 Description of Monte Carlo event generation applications. More...
 
struct  physics
 Generator for neutrino interaction. More...
 
struct  simul
 Generator for simulation. More...
 
struct  spectrum
 Neutrino energy spectrum. More...
 
struct  can
 The cylinder used for photon tracking. More...
 
struct  fixedcan
 The fixe cylinder used for photon tracking. More...
 
struct  genvol
 Neutrino vertex volume. More...
 
struct  coord_origin
 Coordinate origin. More...
 
struct  genhencut
 Phase space for incident neutrino. More...
 
struct  norma
 Normlisation of CORSIKA events. More...
 
struct  livetime
 Normalisation of MUPAGE events. More...
 
struct  flux
 Neutrino flux. More...
 
struct  seabottom
 The bottom of the sea. More...
 
struct  DAQ
 Livetime of DAQ data. More...
 
struct  tgen
 Scaling time of neutrino interaction generators. More...
 
struct  primary
 Primary particle. More...
 
struct  end_event
 End of event record. More...
 
struct  JHead
 Monte Carlo run header. More...
 
struct  getMUPAGEHeader
 Match header for MUPAGE. More...
 
struct  getGenhenHeader
 Match header for genhen. More...
 
struct  getGenieHeader
 Match header for genie. More...
 
struct  getGSeaGenHeader
 Match header for gseagen. More...
 
struct  getSireneHeader
 Match header for JSirene.cc. More...
 
struct  getKM3Header
 Match header for km3. More...
 
struct  getKM3SimHeader
 Match header for km3. More...
 
struct  getCorsikaHeader
 Match header for Corsika. More...
 
struct  getDAQHeader
 Match header for DAQ. More...
 
struct  JHit_t
 Auxiliary class to set-up Hit. More...
 
struct  JHits_t
 Auxiliary data structure for list of hits with hit merging capability. More...
 
struct  JParticle
 Auxiliary class to handle particle name, codes and mass. More...
 
struct  JPDB
 Collection of particles. More...
 
struct  JTrk_t
 Auxiliary class to set-up Trk. More...
 
struct  JWeightMupage
 Implementation of event weighing for MUPAGE data. More...
 
struct  JWeightEvent
 Low-level interface for event weighing. More...
 
struct  JWeightEventHelper
 Helper class for event weighing. More...
 
struct  JWeightFileScanner
 Template file scanner with event weight. More...
 

Enumerations

enum  JHitType_t {
  HIT_TYPE_MUON_DIRECT = +5, HIT_TYPE_MUON_SCATTERED = -5, HIT_TYPE_DELTARAYS_DIRECT = +4, HIT_TYPE_DELTARAYS_SCATTERED = -4,
  HIT_TYPE_BREMSSTRAHLUNG_DIRECT = +3, HIT_TYPE_BREMSSTRAHLUNG_SCATTERED = -3, HIT_TYPE_SHOWER_DIRECT = +99, HIT_TYPE_SHOWER_SCATTERED = -99,
  HIT_TYPE_NOISE = -1, HIT_TYPE_UNKNOWN = 0
}
 Enumeration of hit types based on km3 codes. More...
 
enum  JGeant4Type_t {
  GEANT4_TYPE_PHOTON = 1, GEANT4_TYPE_ANTIELECTRON = 2, GEANT4_TYPE_ELECTRON = 3, GEANT4_TYPE_NEUTRINO = 4,
  GEANT4_TYPE_ANTIMUON = 5, GEANT4_TYPE_MUON = 6, GEANT4_TYPE_NEUTRAL_PION = 7, GEANT4_TYPE_PION_PLUS = 8,
  GEANT4_TYPE_PION_MINUS = 9, GEANT4_TYPE_KAON_LONG = 10, GEANT4_TYPE_KAON_PLUS = 11, GEANT4_TYPE_KAON_MINUS = 12,
  GEANT4_TYPE_NEUTRON = 13, GEANT4_TYPE_PROTON = 14, GEANT4_TYPE_ANTIPROTON = 15, GEANT4_TYPE_KAON_SHORT = 16,
  GEANT4_TYPE_ETA = 17, GEANT4_TYPE_LAMBDA = 18, GEANT4_TYPE_SIGMA_PLUS = 19, GEANT4_TYPE_NEUTRAL_SIGMA = 20,
  GEANT4_TYPE_SIGMA_MINUS = 21, GEANT4_TYPE_NEUTRAL_XI = 22, GEANT4_TYPE_XI_MINUS = 23, GEANT4_TYPE_OMEGA_MINUS = 24,
  GEANT4_TYPE_ANTINEUTRON = 25, GEANT4_TYPE_ANTILAMBDA = 26, GEANT4_TYPE_ANTISIGMA_MINUS = 27, GEANT4_TYPE_NEUTRAL_ANTISIGMA = 28,
  GEANT4_TYPE_ANTISIGMA_PLUS = 29, GEANT4_TYPE_NEUTRAL_ANTIXI = 30, GEANT4_TYPE_ANTIXI_PLUS = 31, GEANT4_TYPE_ANTIOMEGA_PLUS = 32,
  GEANT4_TYPE_DEUTERON = 45, GEANT4_TYPE_TRITON = 46, GEANT4_TYPE_ALPHA = 47, GEANT4_TYPE_GEANTINO = 48,
  GEANT4_TYPE_HE3 = 49, GEANT4_TYPE_ANTITAU = 33, GEANT4_TYPE_TAU = 34
}
 Enumeration of hit types based on Geant4 codes, for compatbility with KM3Sim. More...
 
enum  JTrackType_t {
  TRACK_TYPE_ELECTRON = 11, TRACK_TYPE_NUE = 12, TRACK_TYPE_MUON = 13, TRACK_TYPE_NUMU = 14,
  TRACK_TYPE_TAU = 15, TRACK_TYPE_NUTAU = 16, TRACK_TYPE_PHOTON = 22, TRACK_TYPE_NEUTRAL_PION = 111,
  TRACK_TYPE_CHARGED_PION_PLUS = 211, TRACK_TYPE_PION_PLUS = 211, TRACK_TYPE_K_LONG = 130, TRACK_TYPE_K_SHORT = 310,
  TRACK_TYPE_K_PLUS = 321, TRACK_TYPE_PROTON = 2212, TRACK_TYPE_NEUTRON = 2112, TRACK_TYPE_LAMBDA = 3122,
  TRACK_TYPE_SIGMA_PLUS = 3222, TRACK_TYPE_NEUTRAL_SIGMA = 3212, TRACK_TYPE_SIGMA_MINUS = 3112, TRACK_TYPE_NEUTRAL_XI = 3322,
  TRACK_TYPE_XI_MINUS = 3312, TRACK_TYPE_OMEGA_MINUS = 3334, TRACK_TYPE_ANTIELECTRON = -11, TRACK_TYPE_ANTINUE = -12,
  TRACK_TYPE_ANTIMUON = -13, TRACK_TYPE_ANTINUMU = -14, TRACK_TYPE_ANTITAU = -15, TRACK_TYPE_ANTINUTAU = -16,
  TRACK_TYPE_NEUTRAL_ANTIPION = -111, TRACK_TYPE_CHARGED_PION_MINUS = -211, TRACK_TYPE_PION_MINUS = -211, TRACK_TYPE_ANTIK_LONG = -130,
  TRACK_TYPE_ANTIK_SHORT = -310, TRACK_TYPE_K_MINUS = -321, TRACK_TYPE_ANTIPROTON = -2212, TRACK_TYPE_ANTINEUTRON = -2112,
  TRACK_TYPE_ANTILAMBDA = -3122, TRACK_TYPE_ANTISIGMA_PLUS = -3222, TRACK_TYPE_ANTINEUTRAL_SIGMA = -3212, TRACK_TYPE_ANTISIGMA_MINUS = -3112,
  TRACK_TYPE_ANTINEUTRAL_XI = -3322, TRACK_TYPE_ANTIXI_MINUS = -3312, TRACK_TYPE_ANTIOMEGA_MINUS = -3334
}
 Enumeration of track types based on PDG codes. More...
 

Functions

double getTime (const Hit &hit)
 Get true time of hit. More...
 
double getNPE (const Hit &hit)
 Get true charge of hit. More...
 
bool is_noise (const Hit &hit)
 Verify hit origin. More...
 
JTimeRange getTimeRange (const Evt &event)
 Get time range (i.e. time between earliest and latest hit) of Monte Carlo event. More...
 
JTimeRange getTimeRange (const Evt &event, const JTimeRange &T_ns)
 Get time range (i.e. time between earliest and latest hit) of Monte Carlo event. More...
 
JPosition3D getPosition (const Vec &v)
 Get position. More...
 
JPosition3D getPosition (const Trk &track)
 Get position. More...
 
JDirection3D getDirection (const Vec &v)
 Get direction. More...
 
JDirection3D getDirection (const Trk &track)
 Get direction. More...
 
JAxis3D getAxis (const Trk &track)
 Get axis. More...
 
JTrack3E getTrack (const Trk &track)
 Get track. More...
 
JVertex3D getVertex (const Trk &track)
 Get vertex. More...
 
JTransformation3D getTransformation (const Trk &track)
 Get transformation. More...
 
JHitL0 getHit (const Hit &hit)
 Get transformation. More...
 
JHitL0 getHit (const Hit &hit, const JModuleRouter &router)
 Get transformation. More...
 
JHitL0 getHit (const Hit &hit, const JPMTRouter &router)
 Get transformation. More...
 
bool is_photon (const Trk &track)
 Test whether given track is a photon or neutral pion. More...
 
bool is_neutrino (const Trk &track)
 Test whether given track is a neutrino. More...
 
bool is_electron (const Trk &track)
 Test whether given track is a (anti-)electron. More...
 
bool is_muon (const Trk &track)
 Test whether given track is a (anti-)muon. More...
 
bool is_tau (const Trk &track)
 Test whether given track is a (anti-)tau. More...
 
bool is_pion (const Trk &track)
 Test whether given track is a charged pion. More...
 
bool is_proton (const Trk &track)
 Test whether given track is a (anti-)proton. More...
 
bool is_neutron (const Trk &track)
 Test whether given track is a (anti-)neutron. More...
 
bool is_hadron (const Trk &track)
 Test whether given track is a hadron. More...
 
bool has_particleID (const Trk &track, const int type)
 Test whether given track corresponds to given particle type. More...
 
bool has_neutrino (const Evt &evt)
 Test whether given event has an incoming neutrino. More...
 
const Trkget_neutrino (const Evt &evt)
 Get incoming neutrino. More...
 
int count_electrons (const Evt &evt)
 Count the number of electrons in a given event. More...
 
bool has_electron (const Evt &evt)
 Test whether given event has an electron. More...
 
const Trkget_electron (const Evt &evt)
 Get first electron from the event tracklist. More...
 
bool from_electron (const Hit &hit)
 Test whether given hit was produced by an electron. More...
 
int count_muons (const Evt &evt)
 Count the number of muons in a given event. More...
 
bool has_muon (const Evt &evt)
 Test whether given event has a muon. More...
 
const Trkget_muon (const Evt &evt)
 Get first muon from the event tracklist. More...
 
bool from_muon (const Hit &hit)
 Test whether given hit was produced by a muon. More...
 
int count_taus (const Evt &evt)
 Count the number of taus in a given event. More...
 
bool has_tau (const Evt &evt)
 Test whether given event has a tau. More...
 
const Trkget_tau (const Evt &evt)
 Get first tau from the event tracklist. More...
 
bool from_tau (const Hit &hit)
 Test whether given hit was produced by a tau. More...
 
int count_hadrons (const Evt &evt)
 Count the number of hadrons in a given event (not including neutral pions) More...
 
bool from_hadron (const Hit &hit)
 Test whether given hit was produced by a hadronic shower. More...
 
Trk get_hadronic_cascade (const Evt &evt)
 Auxiliary function to get average direction and total energy of a hadronic cascade. More...
 
void add_time_offset (Evt &evt, const double tOff)
 Add time offset to mc event; according to current definition, the absolute time of the event is defined by the track "t" attribute; this could change in the future if the global attribute mc_t is assigned to this purpose. More...
 
void copy (const Head &from, JHead &to)
 Copy header from from to to. More...
 
void copy (const JHead &from, Head &to)
 Copy header from from to to. More...
 
std::string getTag (const std::string &tag)
 Get tag without aanet extension "_<counter>" for identical tags. More...
 
bool operator== (const Head &first, const Head &second)
 Equal operator. More...
 
bool operator< (const Head &first, const Head &second)
 Less than operator. More...
 
bool is_mupage (const JHead &header)
 Check for generator. More...
 
bool is_genhen (const JHead &header)
 Check for generator. More...
 
bool is_genie (const JHead &header)
 Check for generator. More...
 
bool is_sirene (const JHead &header)
 Check for detector simulation. More...
 
bool is_km3 (const JHead &header)
 Check for detector simulation. More...
 
bool is_km3sim (const JHead &header)
 Check for detector simulation. More...
 
bool is_daq (const JHead &header)
 Check for real data. More...
 
template<class T >
T get (const JHead &header)
 Get object from header. More...
 
template<>
Vec get (const JHead &header)
 Get position offset of detector due to generator. More...
 
template<>
JPosition3D get (const JHead &header)
 Get position offset of detector due to generator. More...
 
template<>
JCylinder3D get (const JHead &header)
 Get cylinder corresponding to can. More...
 

Variables

static const int AASHOWER_RECONSTRUCTION_TYPE = 101
 AAshower reconstruction type for AAnet. More...
 
static const JEvtEvaluator getEvtValue
 Function object for evaluation of DAQ objects. More...
 
static const char AANET_TAG_SEPARATOR = '_'
 Separator for AAnet tag extension for multiple tags ("_<counter>"). More...
 

Detailed Description

Extensions to AAnet data format.

Author
mlincetto
mdejong

Enumeration Type Documentation

Enumeration of hit types based on km3 codes.

Enumerator
HIT_TYPE_MUON_DIRECT 

Direct light from muon.

HIT_TYPE_MUON_SCATTERED 

Scattered light from muon.

HIT_TYPE_DELTARAYS_DIRECT 

Direct light from delta-rays.

HIT_TYPE_DELTARAYS_SCATTERED 

Scattered light from delta-rays.

HIT_TYPE_BREMSSTRAHLUNG_DIRECT 

Direct light from Bremsstrahlung.

HIT_TYPE_BREMSSTRAHLUNG_SCATTERED 

Scattered light from Bremsstrahlung.

HIT_TYPE_SHOWER_DIRECT 

Direct light from primary shower.

HIT_TYPE_SHOWER_SCATTERED 

Scattered light from primary shower.

HIT_TYPE_NOISE 

Random noise.

HIT_TYPE_UNKNOWN 

Unknown source.

Definition at line 68 of file JAAnetToolkit.hh.

68  {
69  HIT_TYPE_MUON_DIRECT = +5, //!< Direct light from muon
70  HIT_TYPE_MUON_SCATTERED = -5, //!< Scattered light from muon
71  HIT_TYPE_DELTARAYS_DIRECT = +4, //!< Direct light from delta-rays
72  HIT_TYPE_DELTARAYS_SCATTERED = -4, //!< Scattered light from delta-rays
73  HIT_TYPE_BREMSSTRAHLUNG_DIRECT = +3, //!< Direct light from Bremsstrahlung
74  HIT_TYPE_BREMSSTRAHLUNG_SCATTERED = -3, //!< Scattered light from Bremsstrahlung
75  HIT_TYPE_SHOWER_DIRECT = +99, //!< Direct light from primary shower
76  HIT_TYPE_SHOWER_SCATTERED = -99, //!< Scattered light from primary shower
77  HIT_TYPE_NOISE = -1, //!< Random noise
78  HIT_TYPE_UNKNOWN = 0 //!< Unknown source
79  };
Scattered light from primary shower.
Direct light from delta-rays.
Scattered light from muon.
Scattered light from Bremsstrahlung.
Direct light from Bremsstrahlung.
Direct light from primary shower.
Direct light from muon.
Scattered light from delta-rays.

Enumeration of hit types based on Geant4 codes, for compatbility with KM3Sim.

Enumerator
GEANT4_TYPE_PHOTON 
GEANT4_TYPE_ANTIELECTRON 
GEANT4_TYPE_ELECTRON 
GEANT4_TYPE_NEUTRINO 
GEANT4_TYPE_ANTIMUON 
GEANT4_TYPE_MUON 
GEANT4_TYPE_NEUTRAL_PION 
GEANT4_TYPE_PION_PLUS 
GEANT4_TYPE_PION_MINUS 
GEANT4_TYPE_KAON_LONG 
GEANT4_TYPE_KAON_PLUS 
GEANT4_TYPE_KAON_MINUS 
GEANT4_TYPE_NEUTRON 
GEANT4_TYPE_PROTON 
GEANT4_TYPE_ANTIPROTON 
GEANT4_TYPE_KAON_SHORT 
GEANT4_TYPE_ETA 
GEANT4_TYPE_LAMBDA 
GEANT4_TYPE_SIGMA_PLUS 
GEANT4_TYPE_NEUTRAL_SIGMA 
GEANT4_TYPE_SIGMA_MINUS 
GEANT4_TYPE_NEUTRAL_XI 
GEANT4_TYPE_XI_MINUS 
GEANT4_TYPE_OMEGA_MINUS 
GEANT4_TYPE_ANTINEUTRON 
GEANT4_TYPE_ANTILAMBDA 
GEANT4_TYPE_ANTISIGMA_MINUS 
GEANT4_TYPE_NEUTRAL_ANTISIGMA 
GEANT4_TYPE_ANTISIGMA_PLUS 
GEANT4_TYPE_NEUTRAL_ANTIXI 
GEANT4_TYPE_ANTIXI_PLUS 
GEANT4_TYPE_ANTIOMEGA_PLUS 
GEANT4_TYPE_DEUTERON 
GEANT4_TYPE_TRITON 
GEANT4_TYPE_ALPHA 
GEANT4_TYPE_GEANTINO 
GEANT4_TYPE_HE3 
GEANT4_TYPE_ANTITAU 
GEANT4_TYPE_TAU 

Definition at line 18 of file JParticleTypes.hh.

18  { GEANT4_TYPE_PHOTON = 1,
23  GEANT4_TYPE_MUON = 6,
31  GEANT4_TYPE_PROTON = 14,
34  GEANT4_TYPE_ETA = 17,
35  GEANT4_TYPE_LAMBDA = 18,
51  GEANT4_TYPE_TRITON = 46,
52  GEANT4_TYPE_ALPHA = 47,
54  GEANT4_TYPE_HE3 = 49,
55  //
56  // KM3NeT specific codes
57  //
59  GEANT4_TYPE_TAU = 34 };

Enumeration of track types based on PDG codes.

Enumerator
TRACK_TYPE_ELECTRON 
TRACK_TYPE_NUE 
TRACK_TYPE_MUON 
TRACK_TYPE_NUMU 
TRACK_TYPE_TAU 
TRACK_TYPE_NUTAU 
TRACK_TYPE_PHOTON 
TRACK_TYPE_NEUTRAL_PION 
TRACK_TYPE_CHARGED_PION_PLUS 
TRACK_TYPE_PION_PLUS 
TRACK_TYPE_K_LONG 
TRACK_TYPE_K_SHORT 
TRACK_TYPE_K_PLUS 
TRACK_TYPE_PROTON 
TRACK_TYPE_NEUTRON 
TRACK_TYPE_LAMBDA 
TRACK_TYPE_SIGMA_PLUS 
TRACK_TYPE_NEUTRAL_SIGMA 
TRACK_TYPE_SIGMA_MINUS 
TRACK_TYPE_NEUTRAL_XI 
TRACK_TYPE_XI_MINUS 
TRACK_TYPE_OMEGA_MINUS 
TRACK_TYPE_ANTIELECTRON 
TRACK_TYPE_ANTINUE 
TRACK_TYPE_ANTIMUON 
TRACK_TYPE_ANTINUMU 
TRACK_TYPE_ANTITAU 
TRACK_TYPE_ANTINUTAU 
TRACK_TYPE_NEUTRAL_ANTIPION 
TRACK_TYPE_CHARGED_PION_MINUS 
TRACK_TYPE_PION_MINUS 
TRACK_TYPE_ANTIK_LONG 
TRACK_TYPE_ANTIK_SHORT 
TRACK_TYPE_K_MINUS 
TRACK_TYPE_ANTIPROTON 
TRACK_TYPE_ANTINEUTRON 
TRACK_TYPE_ANTILAMBDA 
TRACK_TYPE_ANTISIGMA_PLUS 
TRACK_TYPE_ANTINEUTRAL_SIGMA 
TRACK_TYPE_ANTISIGMA_MINUS 
TRACK_TYPE_ANTINEUTRAL_XI 
TRACK_TYPE_ANTIXI_MINUS 
TRACK_TYPE_ANTIOMEGA_MINUS 

Definition at line 64 of file JParticleTypes.hh.

64  { TRACK_TYPE_ELECTRON = 11,
65  TRACK_TYPE_NUE = 12,
66  TRACK_TYPE_MUON = 13,
67  TRACK_TYPE_NUMU = 14,
68  TRACK_TYPE_TAU = 15,
69  TRACK_TYPE_NUTAU = 16,
70  TRACK_TYPE_PHOTON = 22,
74  TRACK_TYPE_K_LONG = 130,
75  TRACK_TYPE_K_SHORT = 310,
76  TRACK_TYPE_K_PLUS = 321,
77  TRACK_TYPE_PROTON = 2212,
78  TRACK_TYPE_NEUTRON = 2112,
79  TRACK_TYPE_LAMBDA = 3122,
80  TRACK_TYPE_SIGMA_PLUS = 3222,
83  TRACK_TYPE_NEUTRAL_XI = 3322,
84  TRACK_TYPE_XI_MINUS = 3312,
86 
88  TRACK_TYPE_ANTINUE = -12,
89  TRACK_TYPE_ANTIMUON = -13,
90  TRACK_TYPE_ANTINUMU = -14,
91  TRACK_TYPE_ANTITAU = -15,
95  TRACK_TYPE_PION_MINUS = -211,
96  TRACK_TYPE_ANTIK_LONG = -130,
98  TRACK_TYPE_K_MINUS = -321,
99  TRACK_TYPE_ANTIPROTON = -2212,
100  TRACK_TYPE_ANTINEUTRON = -2112,
101  TRACK_TYPE_ANTILAMBDA = -3122,
106  TRACK_TYPE_ANTIXI_MINUS = -3312,
107  TRACK_TYPE_ANTIOMEGA_MINUS = -3334 };

Function Documentation

double JAANET::getTime ( const Hit hit)
inline

Get true time of hit.

Parameters
hithit
Returns
true time of hit

Definition at line 88 of file JAAnetToolkit.hh.

89  {
90  // if IGAIN = 0 in km3 program, then pure_xx values are not set.
91 
92  if (hit.pure_a == 0.0)
93  return hit.t;
94  else
95  return hit.pure_t;
96  }
double pure_t
photon time before pmt simultion (MC only)
Definition: Hit.hh:28
double pure_a
amptitude before pmt simution (MC only)
Definition: Hit.hh:29
double t
hit time (from calibration or MC truth)
Definition: Hit.hh:23
double JAANET::getNPE ( const Hit hit)
inline

Get true charge of hit.

Parameters
hithit
Returns
true charge of hit

Definition at line 105 of file JAAnetToolkit.hh.

106  {
107  // if IGAIN = 0 in km3 program, then pure_xx values are not set.
108 
109  if (hit.pure_a == 0.0)
110  return hit.a;
111  else
112  return hit.pure_a;
113  }
double a
hit amplitude (in p.e.)
Definition: Hit.hh:24
double pure_a
amptitude before pmt simution (MC only)
Definition: Hit.hh:29
bool JAANET::is_noise ( const Hit hit)
inline

Verify hit origin.

Parameters
hithit
Returns
true if noise; else false

Definition at line 122 of file JAAnetToolkit.hh.

123  {
124  return hit.type == HIT_TYPE_NOISE;
125  }
int type
particle type or parametrisation used for hit (mc only)
Definition: Hit.hh:30
JTimeRange JAANET::getTimeRange ( const Evt event)
inline

Get time range (i.e. time between earliest and latest hit) of Monte Carlo event.


Note that the global event time in not taken into account.

Parameters
eventevent
Returns
time range

Definition at line 135 of file JAAnetToolkit.hh.

136  {
137  JTimeRange time_range(JTimeRange::DEFAULT_RANGE);
138 
139  for (std::vector<Hit>::const_iterator hit = event.mc_hits.begin(); hit != event.mc_hits.end(); ++hit) {
140  if (!is_noise(*hit)) {
141  time_range.include(getTime(*hit));
142  }
143  }
144 
145  return time_range;
146  }
double getTime(const Hit &hit)
Get true time of hit.
bool is_noise(const Hit &hit)
Verify hit origin.
std::vector< Hit > mc_hits
MC: list of MC truth hits.
Definition: Evt.hh:44
JTimeRange JAANET::getTimeRange ( const Evt event,
const JTimeRange T_ns 
)
inline

Get time range (i.e. time between earliest and latest hit) of Monte Carlo event.


The resulting time range is larger than or equal to the given time window.
Note that the global event time in not taken into account.

Parameters
eventevent
T_nstime window [ns]
Returns
time range

Definition at line 158 of file JAAnetToolkit.hh.

159  {
160  JTimeRange time_range = getTimeRange(event);
161 
162  if (time_range.getLength() < T_ns.getLength()) {
163  time_range.setLength(T_ns.getLength());
164  }
165 
166  return time_range;
167  }
void setLength(argument_type length)
Set length (difference between upper and lower limit).
Definition: JRange.hh:313
JTimeRange getTimeRange(const Evt &event)
Get time range (i.e. time between earliest and latest hit) of Monte Carlo event.
T getLength() const
Get length (difference between upper and lower limit).
Definition: JRange.hh:302
JPosition3D JAANET::getPosition ( const Vec v)
inline

Get position.

Parameters
vvector
Returns
position

Definition at line 176 of file JAAnetToolkit.hh.

177  {
178  return JPosition3D(v.x, v.y, v.z);
179  }
double z
Definition: Vec.hh:14
double y
Definition: Vec.hh:14
double x
Definition: Vec.hh:14
Data structure for position in three dimensions.
Definition: JPosition3D.hh:35
JPosition3D JAANET::getPosition ( const Trk track)
inline

Get position.

Parameters
tracktrack
Returns
position

Definition at line 188 of file JAAnetToolkit.hh.

189  {
190  return getPosition(track.pos);
191  }
Vec pos
postion of the track at time t [m]
Definition: Trk.hh:15
JPosition3D getPosition(const Vec &v)
Get position.
JDirection3D JAANET::getDirection ( const Vec v)
inline

Get direction.

Parameters
vvector
Returns
direction

Definition at line 200 of file JAAnetToolkit.hh.

201  {
202  return JDirection3D(v.x, v.y, v.z);
203  }
Data structure for direction in three dimensions.
Definition: JDirection3D.hh:32
double z
Definition: Vec.hh:14
double y
Definition: Vec.hh:14
double x
Definition: Vec.hh:14
JDirection3D JAANET::getDirection ( const Trk track)
inline

Get direction.

Parameters
tracktrack
Returns
direction

Definition at line 211 of file JAAnetToolkit.hh.

212  {
213  return getDirection(track.dir);
214  }
Vec dir
track direction
Definition: Trk.hh:16
JDirection3D getDirection(const Vec &v)
Get direction.
JAxis3D JAANET::getAxis ( const Trk track)
inline

Get axis.

Parameters
tracktrack
Returns
axis

Definition at line 223 of file JAAnetToolkit.hh.

224  {
225  return JAxis3D(getPosition(track), getDirection(track));
226  }
Axis object.
Definition: JAxis3D.hh:38
JDirection3D getDirection(const Vec &v)
Get direction.
JPosition3D getPosition(const Vec &v)
Get position.
JTrack3E JAANET::getTrack ( const Trk track)
inline

Get track.

Parameters
tracktrack
Returns
track

Definition at line 235 of file JAAnetToolkit.hh.

236  {
237  return JTrack3E(JTrack3D(getAxis(track), track.t), track.E);
238  }
double t
track time [ns] (when the particle is at pos )
Definition: Trk.hh:17
double E
Energy [GeV] (either MC truth or reconstructed)
Definition: Trk.hh:18
3D track with energy.
Definition: JTrack3E.hh:29
JAxis3D getAxis(const Trk &track)
Get axis.
JVertex3D JAANET::getVertex ( const Trk track)
inline

Get vertex.

Parameters
tracktrack
Returns
vertex

Definition at line 247 of file JAAnetToolkit.hh.

248  {
249  return JVertex3D(getPosition(track), track.t);
250  }
double t
track time [ns] (when the particle is at pos )
Definition: Trk.hh:17
JPosition3D getPosition(const Vec &v)
Get position.
JTransformation3D JAANET::getTransformation ( const Trk track)
inline

Get transformation.

Parameters
tracktrack
Returns
transformation

Definition at line 259 of file JAAnetToolkit.hh.

260  {
261  return JTransformation3D(getPosition(track), getDirection(track));
262  }
JDirection3D getDirection(const Vec &v)
Get direction.
JPosition3D getPosition(const Vec &v)
Get position.
JHitL0 JAANET::getHit ( const Hit hit)
inline

Get transformation.

Parameters
hithit
Returns
hit

Definition at line 271 of file JAAnetToolkit.hh.

272  {
273  using namespace JPP;
274 
275  return JHitL0(JDAQPMTIdentifier(hit.dom_id, hit.channel_id), JAxis3D(getPosition(hit.pos), getDirection(hit.dir)), JHit(hit.t, hit.tot) );
276  }
Vec pos
hit position
Definition: Hit.hh:25
Acoustics hit.
Axis object.
Definition: JAxis3D.hh:38
Vec dir
hit direction; i.e. direction of the PMT
Definition: Hit.hh:26
double t
hit time (from calibration or MC truth)
Definition: Hit.hh:23
Data structure for L0 hit.
Definition: JHitL0.hh:27
int dom_id
module identifier from the data (unique in the detector).
Definition: Hit.hh:14
unsigned int tot
tot value as stored in raw data (int for pyroot)
Definition: Hit.hh:17
JDirection3D getDirection(const Vec &v)
Get direction.
unsigned int channel_id
PMT channel id {0,1, .., 31} local to moduke.
Definition: Hit.hh:15
JPosition3D getPosition(const Vec &v)
Get position.
JHitL0 JAANET::getHit ( const Hit hit,
const JModuleRouter router 
)
inline

Get transformation.

Parameters
hithit
routermodule router
Returns
hit

Definition at line 286 of file JAAnetToolkit.hh.

287  {
288  using namespace JPP;
289 
290  const JPMT& pmt = router.getModule(hit.dom_id).getPMT(hit.channel_id);
291 
292  return JHitL0(JDAQPMTIdentifier(hit.dom_id, hit.channel_id), pmt.getAxis(), JHit(getTime(hit.tdc, pmt.getCalibration()), getToT(hit.tot, pmt.getCalibration())));
293  }
const JModule & getModule(const JObjectID &id) const
Get module parameters.
unsigned int tdc
hit tdc (=time in ns)
Definition: Hit.hh:16
double getTime(const Hit &hit)
Get true time of hit.
Acoustics hit.
esac $JPP_DIR examples JDetector JTransitTime o $OUTPUT_FILE n N $NPE T $TTS_NS d $DEBUG for HISTOGRAM in tts tt2 pmt
Definition: JTransitTime.sh:36
double getToT(const T &tot, const JCalibration &cal)
Get calibrated time-over-threshold of hit.
Data structure for PMT geometry and calibration.
Definition: JPMT.hh:47
const JPMT & getPMT(const int index) const
Get PMT.
Definition: JModule.hh:174
Data structure for L0 hit.
Definition: JHitL0.hh:27
int dom_id
module identifier from the data (unique in the detector).
Definition: Hit.hh:14
unsigned int tot
tot value as stored in raw data (int for pyroot)
Definition: Hit.hh:17
unsigned int channel_id
PMT channel id {0,1, .., 31} local to moduke.
Definition: Hit.hh:15
JHitL0 JAANET::getHit ( const Hit hit,
const JPMTRouter router 
)
inline

Get transformation.

Parameters
hithit
routerpmt router
Returns
hit

Definition at line 303 of file JAAnetToolkit.hh.

304  {
305  using namespace JPP;
306 
307  const JPMT& pmt = router.getPMT(hit.pmt_id);
308 
309  return JHitL0(JDAQPMTIdentifier(hit.dom_id, hit.channel_id), pmt.getAxis(), JHit(getTime(hit.tdc, pmt.getCalibration()), getToT(hit.tot, pmt.getCalibration())));
310  }
int pmt_id
global PMT identifier as found in evt files
Definition: Hit.hh:20
const JPMT & getPMT(const JPMTAddress &address) const
Get PMT.
Definition: JPMTRouter.hh:90
unsigned int tdc
hit tdc (=time in ns)
Definition: Hit.hh:16
double getTime(const Hit &hit)
Get true time of hit.
Acoustics hit.
esac $JPP_DIR examples JDetector JTransitTime o $OUTPUT_FILE n N $NPE T $TTS_NS d $DEBUG for HISTOGRAM in tts tt2 pmt
Definition: JTransitTime.sh:36
double getToT(const T &tot, const JCalibration &cal)
Get calibrated time-over-threshold of hit.
Data structure for PMT geometry and calibration.
Definition: JPMT.hh:47
Data structure for L0 hit.
Definition: JHitL0.hh:27
int dom_id
module identifier from the data (unique in the detector).
Definition: Hit.hh:14
unsigned int tot
tot value as stored in raw data (int for pyroot)
Definition: Hit.hh:17
unsigned int channel_id
PMT channel id {0,1, .., 31} local to moduke.
Definition: Hit.hh:15
bool JAANET::is_photon ( const Trk track)
inline

Test whether given track is a photon or neutral pion.

Parameters
tracktrack
Returns
true if photon or neutral pion; else false

Definition at line 319 of file JAAnetToolkit.hh.

319  { return ( track.type == TRACK_TYPE_PHOTON ||
320  abs(track.type) == TRACK_TYPE_NEUTRAL_PION); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:22
bool JAANET::is_neutrino ( const Trk track)
inline

Test whether given track is a neutrino.

Parameters
tracktrack
Returns
true if neutrino; else false

Definition at line 328 of file JAAnetToolkit.hh.

328  { return (abs(track.type) == TRACK_TYPE_NUE ||
329  abs(track.type) == TRACK_TYPE_NUMU ||
330  abs(track.type) == TRACK_TYPE_NUTAU); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:22
bool JAANET::is_electron ( const Trk track)
inline

Test whether given track is a (anti-)electron.

Parameters
tracktrack
Returns
true if (anti-)electron; else false

Definition at line 338 of file JAAnetToolkit.hh.

338 { return (abs(track.type) == TRACK_TYPE_ELECTRON); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:22
bool JAANET::is_muon ( const Trk track)
inline

Test whether given track is a (anti-)muon.

Parameters
tracktrack
Returns
true if (anti-)muon; else false

Definition at line 346 of file JAAnetToolkit.hh.

346 { return (abs(track.type) == TRACK_TYPE_MUON); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:22
bool JAANET::is_tau ( const Trk track)
inline

Test whether given track is a (anti-)tau.

Parameters
tracktrack
Returns
true if (anti-)tau; else false

Definition at line 354 of file JAAnetToolkit.hh.

354 { return (abs(track.type) == TRACK_TYPE_TAU); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:22
bool JAANET::is_pion ( const Trk track)
inline

Test whether given track is a charged pion.

Parameters
tracktrack
Returns
true if charged pion; else false

Definition at line 362 of file JAAnetToolkit.hh.

362 { return (abs(track.type) == TRACK_TYPE_CHARGED_PION_PLUS); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:22
bool JAANET::is_proton ( const Trk track)
inline

Test whether given track is a (anti-)proton.

Parameters
tracktrack
Returns
true if (anti-)proton; else false

Definition at line 370 of file JAAnetToolkit.hh.

370 { return (abs(track.type) == TRACK_TYPE_PROTON); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:22
bool JAANET::is_neutron ( const Trk track)
inline

Test whether given track is a (anti-)neutron.

Parameters
tracktrack
Returns
true if (anti-)neutron; else false

Definition at line 378 of file JAAnetToolkit.hh.

378 { return (abs(track.type) == TRACK_TYPE_NEUTRON); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:22
bool JAANET::is_hadron ( const Trk track)
inline

Test whether given track is a hadron.

Parameters
tracktrack
Returns
true if hadron; else fails

Definition at line 386 of file JAAnetToolkit.hh.

386  { return (!is_neutrino(track) &&
387  !is_electron(track) &&
388  !is_muon (track) &&
389  !is_tau (track)); }
bool is_electron(const Trk &track)
Test whether given track is a (anti-)electron.
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
bool is_neutrino(const Trk &track)
Test whether given track is a neutrino.
bool is_tau(const Trk &track)
Test whether given track is a (anti-)tau.
bool JAANET::has_particleID ( const Trk track,
const int  type 
)
inline

Test whether given track corresponds to given particle type.

Parameters
tracktrack
typeparticle type
Returns
true if matches the corresponding particle; else false

Definition at line 398 of file JAAnetToolkit.hh.

398 { return (track.type == type); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:22
bool JAANET::has_neutrino ( const Evt evt)
inline

Test whether given event has an incoming neutrino.

Parameters
evtevent
Returns
true if neutrino; else false

Definition at line 406 of file JAAnetToolkit.hh.

407  {
408  return !evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0]);
409  }
bool is_neutrino(const Trk &track)
Test whether given track is a neutrino.
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:45
const Trk& JAANET::get_neutrino ( const Evt evt)
inline

Get incoming neutrino.

Parameters
evtevent
Returns
neutrino

Definition at line 417 of file JAAnetToolkit.hh.

418  {
419  if (has_neutrino(evt))
420  return evt.mc_trks[0];
421  else
422  THROW(JIndexOutOfRange, "This event has no neutrino.");
423  }
bool has_neutrino(const Evt &evt)
Test whether given event has an incoming neutrino.
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:670
Exception for accessing an index in a collection that is outside of its range.
Definition: JException.hh:90
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:45
int JAANET::count_electrons ( const Evt evt)
inline

Count the number of electrons in a given event.

Parameters
evtevent
Returns
number of electrons

Definition at line 431 of file JAAnetToolkit.hh.

432  {
433  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
434  }
bool is_electron(const Trk &track)
Test whether given track is a (anti-)electron.
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:45
bool JAANET::has_electron ( const Evt evt)
inline

Test whether given event has an electron.

Parameters
evtevent
Returns
true if event has electron; else false

Definition at line 442 of file JAAnetToolkit.hh.

443  {
444  return count_electrons(evt) != 0;
445  }
int count_electrons(const Evt &evt)
Count the number of electrons in a given event.
const Trk& JAANET::get_electron ( const Evt evt)
inline

Get first electron from the event tracklist.

Parameters
evtevent
Returns
electron

Definition at line 453 of file JAAnetToolkit.hh.

454  {
455  if (count_electrons(evt) > 0)
456  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
457  else
458  THROW(JIndexOutOfRange, "This event has no electron.");
459  }
bool is_electron(const Trk &track)
Test whether given track is a (anti-)electron.
int count_electrons(const Evt &evt)
Count the number of electrons in a given event.
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:670
Exception for accessing an index in a collection that is outside of its range.
Definition: JException.hh:90
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:45
bool JAANET::from_electron ( const Hit hit)
inline

Test whether given hit was produced by an electron.

Warning: This method will only work with the output of KM3Sim. For JSirene or KM3, you should check if track.id == hit.origin instead.

Parameters
hithit
Returns
true if hit was produced by an electron; else false

Definition at line 470 of file JAAnetToolkit.hh.

471  {
472  return hit.type == GEANT4_TYPE_ELECTRON || hit.type == GEANT4_TYPE_ANTIELECTRON;
473  }
int type
particle type or parametrisation used for hit (mc only)
Definition: Hit.hh:30
int JAANET::count_muons ( const Evt evt)
inline

Count the number of muons in a given event.

Parameters
evtevent
Returns
number of muons

Definition at line 481 of file JAAnetToolkit.hh.

482  {
483  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
484  }
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:45
bool JAANET::has_muon ( const Evt evt)
inline

Test whether given event has a muon.

Parameters
evtevent
Returns
true if event has muon; else false

Definition at line 492 of file JAAnetToolkit.hh.

493  {
494  return count_muons(evt) != 0;
495  }
int count_muons(const Evt &evt)
Count the number of muons in a given event.
const Trk& JAANET::get_muon ( const Evt evt)
inline

Get first muon from the event tracklist.

Parameters
evtevent
Returns
muon

Definition at line 503 of file JAAnetToolkit.hh.

504  {
505  if (count_muons(evt) > 0)
506  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
507  else
508  THROW(JIndexOutOfRange, "This event has no muon.");
509  }
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
int count_muons(const Evt &evt)
Count the number of muons in a given event.
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:670
Exception for accessing an index in a collection that is outside of its range.
Definition: JException.hh:90
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:45
bool JAANET::from_muon ( const Hit hit)
inline

Test whether given hit was produced by a muon.

Warning: This method will only work with the output of KM3Sim. For JSirene or KM3, you should check if track.id == hit.origin instead.

Parameters
hithit
Returns
true if hit was produced by a muon; else false

Definition at line 520 of file JAAnetToolkit.hh.

521  {
522  return hit.type == GEANT4_TYPE_MUON || hit.type == GEANT4_TYPE_ANTIMUON;
523  }
int type
particle type or parametrisation used for hit (mc only)
Definition: Hit.hh:30
int JAANET::count_taus ( const Evt evt)
inline

Count the number of taus in a given event.

Parameters
evtevent
Returns
number of taus

Definition at line 531 of file JAAnetToolkit.hh.

532  {
533  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
534  }
bool is_tau(const Trk &track)
Test whether given track is a (anti-)tau.
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:45
bool JAANET::has_tau ( const Evt evt)
inline

Test whether given event has a tau.

Parameters
evtevent
Returns
true if event has tau; else false

Definition at line 542 of file JAAnetToolkit.hh.

543  {
544  return count_taus(evt) != 0;
545  }
int count_taus(const Evt &evt)
Count the number of taus in a given event.
const Trk& JAANET::get_tau ( const Evt evt)
inline

Get first tau from the event tracklist.

Parameters
evtevent
Returns
tau

Definition at line 553 of file JAAnetToolkit.hh.

554  {
555  if (count_taus(evt) > 0)
556  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
557  else
558  THROW(JIndexOutOfRange, "This event has no tau.");
559  }
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:670
int count_taus(const Evt &evt)
Count the number of taus in a given event.
Exception for accessing an index in a collection that is outside of its range.
Definition: JException.hh:90
bool is_tau(const Trk &track)
Test whether given track is a (anti-)tau.
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:45
bool JAANET::from_tau ( const Hit hit)
inline

Test whether given hit was produced by a tau.

Warning: This method will only work with the output of KM3Sim. For JSirene or KM3, you should check if track.id == hit.origin instead.

Parameters
hithit
Returns
true if hit was produced by a tau; else false

Definition at line 570 of file JAAnetToolkit.hh.

571  {
572  return hit.type == GEANT4_TYPE_TAU || hit.type == GEANT4_TYPE_ANTITAU;
573  }
int type
particle type or parametrisation used for hit (mc only)
Definition: Hit.hh:30
int JAANET::count_hadrons ( const Evt evt)
inline

Count the number of hadrons in a given event (not including neutral pions)

Parameters
evtevent
Returns
number of hadrons

Definition at line 581 of file JAAnetToolkit.hh.

582  {
583  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_hadron);
584  }
bool is_hadron(const Trk &track)
Test whether given track is a hadron.
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:45
bool JAANET::from_hadron ( const Hit hit)
inline

Test whether given hit was produced by a hadronic shower.

Warning: This method will only work with the output of KM3Sim. For JSirene or KM3, you should check if track.id == hit.origin instead.

Parameters
hithit
Returns
true if hit was produced by a hadron; else false

Definition at line 595 of file JAAnetToolkit.hh.

596  {
597  return (!from_electron(hit) && !from_muon(hit) && !from_tau(hit));
598  }
bool from_tau(const Hit &hit)
Test whether given hit was produced by a tau.
bool from_electron(const Hit &hit)
Test whether given hit was produced by an electron.
bool from_muon(const Hit &hit)
Test whether given hit was produced by a muon.
Trk JAANET::get_hadronic_cascade ( const Evt evt)
inline

Auxiliary function to get average direction and total energy of a hadronic cascade.

Parameters
evtevent
Returns
cascade (default if there are no hadrons)

Definition at line 606 of file JAAnetToolkit.hh.

607  {
608  Trk cascade;
609 
610  for (std::vector<Trk>::const_iterator track = evt.mc_trks.begin(); track != evt.mc_trks.end(); ++track) {
611  if (is_hadron(*track)) {
612  cascade.E += track->E;
613  cascade.dir += track->dir * track->E;
614  cascade.t = track->t;
615  cascade.pos = track->pos;
616  }
617  }
618 
619  cascade.dir.normalize();
620 
621  return cascade;
622  }
double t
track time [ns] (when the particle is at pos )
Definition: Trk.hh:17
Vec dir
track direction
Definition: Trk.hh:16
double E
Energy [GeV] (either MC truth or reconstructed)
Definition: Trk.hh:18
Vec & normalize()
Normalise this vector.
Definition: Vec.hh:159
Vec pos
postion of the track at time t [m]
Definition: Trk.hh:15
bool is_hadron(const Trk &track)
Test whether given track is a hadron.
The Trk class represents a Monte Carlo (MC) particle as well as a reconstructed track/shower.
Definition: Trk.hh:12
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:45
void JAANET::add_time_offset ( Evt evt,
const double  tOff 
)
inline

Add time offset to mc event; according to current definition, the absolute time of the event is defined by the track "t" attribute; this could change in the future if the global attribute mc_t is assigned to this purpose.

Parameters
evtevent
tOfftime offset [ns]

Definition at line 632 of file JAAnetToolkit.hh.

633  {
634  for (std::vector<Trk>::iterator p = evt.mc_trks.begin(); p != evt.mc_trks.end(); p++) {
635  p->t += tOff;
636  }
637  }
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:45
void JAANET::copy ( const Head from,
JHead &  to 
)

Copy header from from to to.

Parameters
fromheader
toheader

Definition at line 153 of file JHead.cc.

154  {
155  using namespace std;
156  using namespace JPP;
157 
158  JRootReader reader(null, JHead::getEquationParameters(), JAAnetDictionary::getInstance());
159  JRootReadableClass cls(to);
160 
161  for (Head::const_iterator i = from.begin(); i != from.end(); ++i) {
162 
163  const JRootReadableClass& abc = cls.find(getTag(i->first).c_str());
164 
165  const string buffer = trim(i->second);
166 
167  if (abc.is_valid() && buffer != "") {
168 
169  JRedirectString redirect(reader, buffer);
170 
171  reader.getObject(abc);
172 
173  if (i->first == getTag(i->first)) {
174  to.insert(*i); // keep track of parsed tags
175  }
176 
177  } else {
178 
179  //to.insert(*i); // store data of unknown tags
180  }
181  }
182  }
const char * c_str() const
C-string.
Definition: JTag.hh:195
std::string trim(const std::string &buffer)
Trim string.
Definition: JLangToolkit.hh:79
JNET::JTag getTag(JLANG::JType< KM3NETDAQ::JDAQTimeslice >)
Definition: JDAQTags.hh:82
T & getInstance(const T &object)
Get static instance from temporary object.
Definition: JObject.hh:75
esac $JPP_BIN JLogger sh $LOGGER until pgrep JGetMessage</dev/null > dev null
void JAANET::copy ( const JHead &  from,
Head to 
)

Copy header from from to to.

Parameters
fromheader
toheader

Definition at line 191 of file JHead.cc.

192  {
193  using namespace std;
194 
195  to = Head();
196 
197  stringstream io;
198 
199  from.write(io);
200 
201  read(to, io);
202 
203  // copy pending data
204 
205  for (JHead::const_iterator i = from.begin(); i != from.end(); ++i) {
206  if (to.find(i->first) == to.end()) {
207  to.insert(*i);
208  }
209  }
210  }
bool read(Vec &v, std::istream &is)
Read a Vec(tor) from a stream.
Definition: io_ascii.hh:139
The Head class reflects the header of Monte-Carlo event files, which consists of keys (also referred ...
Definition: Head.hh:66
std::string JAANET::getTag ( const std::string &  tag)
inline

Get tag without aanet extension "_<counter>" for identical tags.

Parameters
tagtag
Returns
tag

Definition at line 61 of file JHead.hh.

62  {
63  using namespace std;
64 
65  const string::size_type pos = tag.find(AANET_TAG_SEPARATOR);
66 
67  if (pos != string::npos) {
68 
69  for (string::size_type i = pos + 1; i != tag.size(); ++i) {
70  if (!isdigit(tag[i])) {
71  return tag;
72  }
73  }
74 
75  return tag.substr(0, pos);
76  }
77 
78  return tag;
79  }
static const char AANET_TAG_SEPARATOR
Separator for AAnet tag extension for multiple tags (&quot;_&lt;counter&gt;&quot;).
Definition: JHead.hh:32
bool JAANET::operator== ( const Head first,
const Head second 
)
inline

Equal operator.

Note that this operator uses the JHead::match method.

Parameters
firstfirst header
secondsecond header
Returns
true if two headers are equal; else false

Definition at line 1249 of file JHead.hh.

1251  {
1252  return JHead(first).match(JHead(second));
1253  }
Monte Carlo run header.
Definition: JHead.hh:836
bool match(const JHead &header, const bool option=true) const
Test match of headers.
Definition: JHead.hh:991
bool JAANET::operator< ( const Head first,
const Head second 
)
inline

Less than operator.

Parameters
firstfirst header
secondsecond header
Returns
true if first header is less than second header; else false

Definition at line 1263 of file JHead.hh.

1265  {
1266  return JHead(first).less(JHead(second));
1267  }
bool less(const JHead &header) const
Comparison of headers.
Definition: JHead.hh:1018
Monte Carlo run header.
Definition: JHead.hh:836
bool JAANET::is_mupage ( const JHead &  header)
inline

Check for generator.

Parameters
headerheader
Returns
true if this Monte Carlo is produced by mupage; else false

Definition at line 199 of file JHeadToolkit.hh.

200  {
201  return header.match(getMUPAGEHeader(), false);
202  }
bool JAANET::is_genhen ( const JHead &  header)
inline

Check for generator.

Parameters
headerheader
Returns
true if this Monte Carlo is produced by genhen; else false

Definition at line 211 of file JHeadToolkit.hh.

212  {
213  return header.match(getGenhenHeader(), false);
214  }
bool JAANET::is_genie ( const JHead &  header)
inline

Check for generator.

Parameters
headerheader
Returns
true if this Monte Carlo is produced by genie; else false

Definition at line 223 of file JHeadToolkit.hh.

224  {
225  return (header.match(getGenieHeader(), false) ||
226  header.match(getGSeaGenHeader(), false));
227  }
bool JAANET::is_sirene ( const JHead &  header)
inline

Check for detector simulation.

Parameters
headerheader
Returns
true if this Monte Carlo is processed with JSirene; else false

Definition at line 236 of file JHeadToolkit.hh.

237  {
238  return header.match(getSireneHeader(), false);
239  }
bool JAANET::is_km3 ( const JHead &  header)
inline

Check for detector simulation.

Parameters
headerheader
Returns
true if this Monte Carlo is processed with km3; else false

Definition at line 248 of file JHeadToolkit.hh.

249  {
250  return header.match(getKM3Header(), false);
251  }
bool JAANET::is_km3sim ( const JHead &  header)
inline

Check for detector simulation.

Parameters
headerheader
Returns
true if this Monte Carlo is processed with KM3Sim; else false

Definition at line 260 of file JHeadToolkit.hh.

261  {
262  return header.match(getKM3SimHeader(), false);
263  }
bool JAANET::is_daq ( const JHead &  header)
inline

Check for real data.

Parameters
headerheader
Returns
true if this header corresponds to real data; else false

Definition at line 272 of file JHeadToolkit.hh.

273  {
274  return header.match(getDAQHeader(), false);
275  }
template<class T >
T JAANET::get ( const JHead &  header)
inline

Get object from header.

Parameters
headerheader
Returns
object

Get object from header.

Parameters
headerheader
Returns
position

Get object from header.

Parameters
headerheader
Returns
cylinder

Definition at line 295 of file JHeadToolkit.hh.

296  {
297  if (is_sirene(header) || is_km3(header)) {
298 
299  return header.coord_origin;
300 
301  } else {
302 
303  if (header.is_valid(&JHead::can))
304  return Vec(0.0, 0.0, -header.can.zmin);
305  else if (header.is_valid(&JHead::coord_origin))
306  return header.coord_origin;
307  else
308  return Vec(0.0, 0.0, 0.0);
309  }
310  }
The Vec class is a straightforward 3-d vector, which also works in pyroot.
Definition: Vec.hh:12
bool is_sirene(const JHead &header)
Check for detector simulation.
bool is_km3(const JHead &header)
Check for detector simulation.
template<>
Vec JAANET::get ( const JHead &  header)
inline

Get position offset of detector due to generator.

Get object from header.

Parameters
headerheader
Returns
position

Definition at line 295 of file JHeadToolkit.hh.

296  {
297  if (is_sirene(header) || is_km3(header)) {
298 
299  return header.coord_origin;
300 
301  } else {
302 
303  if (header.is_valid(&JHead::can))
304  return Vec(0.0, 0.0, -header.can.zmin);
305  else if (header.is_valid(&JHead::coord_origin))
306  return header.coord_origin;
307  else
308  return Vec(0.0, 0.0, 0.0);
309  }
310  }
The Vec class is a straightforward 3-d vector, which also works in pyroot.
Definition: Vec.hh:12
bool is_sirene(const JHead &header)
Check for detector simulation.
bool is_km3(const JHead &header)
Check for detector simulation.
template<>
JPosition3D JAANET::get ( const JHead &  header)
inline

Get position offset of detector due to generator.

Get object from header.

Parameters
headerheader
Returns
position

Definition at line 320 of file JHeadToolkit.hh.

321  {
322  const Vec pos = get<Vec>(header);
323 
324  return JPosition3D(pos.x, pos.y, pos.z);
325  }
double z
Definition: Vec.hh:14
double y
Definition: Vec.hh:14
double x
Definition: Vec.hh:14
The Vec class is a straightforward 3-d vector, which also works in pyroot.
Definition: Vec.hh:12
template<>
JCylinder3D JAANET::get ( const JHead &  header)
inline

Get cylinder corresponding to can.

Get object from header.

Parameters
headerheader
Returns
cylinder

Definition at line 335 of file JHeadToolkit.hh.

336  {
337  using namespace JGEOMETRY2D;
338 
339  return JCylinder3D(JCircle2D(JVector2D(), header.can.r),
340  header.can.zmin,
341  header.can.zmax);
342  }
Data structure for vector in two dimensions.
Definition: JVector2D.hh:31
Data structure for circle in two dimensions.
Definition: JCircle2D.hh:29

Variable Documentation

const int JAANET::AASHOWER_RECONSTRUCTION_TYPE = 101
static

AAshower reconstruction type for AAnet.

Definition at line 63 of file JAAnetToolkit.hh.

const JEvtEvaluator JAANET::getEvtValue
static

Function object for evaluation of DAQ objects.

Definition at line 48 of file JEvtEvaluator.hh.

const char JAANET::AANET_TAG_SEPARATOR = '_'
static

Separator for AAnet tag extension for multiple tags ("_<counter>").

Definition at line 32 of file JHead.hh.