Jpp  17.3.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Classes | Typedefs | Enumerations | Functions | Variables
JAANET Namespace Reference

Extensions to Evt data format. More...

Classes

struct  JDiffuseFlux
 Low-level interface for diffuse fluxes. More...
 
struct  JEvtEvaluator
 Auxiliary class to determine value of Evt objects. More...
 
struct  JEvtWeight
 Abstract base class for event weighing. More...
 
struct  JEvtWeightCorsika
 Implementation of event weighting for Corsika data. More...
 
struct  JEvtWeightDAQ
 Implementation of event weighing for DAQ data. More...
 
struct  JEvtWeightFactor
 Low-level interface for retrieving a specifiable multiplication factor
corresponding to a given event. More...
 
struct  JEvtWeightFactorFunction
 Implementation of event-weight factor interface. More...
 
struct  JEvtWeightFactorFunction< JDiffuseFluxFunction_t, JDiffuseFlux >
 Implementation of event-weight factor interface for diffuse flux objects. More...
 
struct  JEvtWeightFactorFunction< pEvtWeightFactor, JEvtWeightFactor_t >
 Implementation of C-style event-weight factor. More...
 
struct  JEvtWeightFactorFunction< pDiffuseFlux, JDiffuseFlux >
 Implementation of C-style diffuse flux event-weight factor. More...
 
struct  JEvtWeightFactorHelper
 Helper class for event-weight factor. More...
 
struct  JEvtWeightFactorHelper< JDiffuseFlux >
 Explicit emplate specialization of event-weight factor helper for diffuse flux objects. More...
 
struct  JEvtWeightFactorMultiParticle
 Implementation of event-weight factor for multiple particle types. More...
 
struct  JEvtWeightFactorMupage
 Implementation of reweighting factor for mupage events according to a specifiable ROOT TFormula. More...
 
struct  JEvtWeightGSeaGen
 Implementation of event weighting for GSeaGen data. More...
 
struct  JEvtWeightHelper
 Helper class for event weighing. More...
 
struct  JEvtWeightInterface
 Low-level interface for event weighing. More...
 
struct  JEvtWeightKM3BUU
 Implementation of event weighting for KM3BUU data. More...
 
struct  JEvtWeightMupage
 Implementation of event weighing for MUPAGE data. More...
 
struct  JEvtWeighter
 Look-up table for event weighters. More...
 
struct  JFlux
 Low-level interface for retrieving the flux corresponding to a given event. More...
 
struct  JFluxAtmosphericNeutrino
 Implementation of atmospheric neutrino flux using official KM3NeT atmospheric flux function. More...
 
struct  JRange_t
 Type definition of range. 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
 Physics information. 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 fixed 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  depth
 Depth. More...
 
struct  DAQ
 Livetime of DAQ data. More...
 
struct  tgen
 Time duration of event generation. More...
 
struct  time_interval
 UTC time interval for event generation. More...
 
struct  primary
 Primary particle. More...
 
struct  end_event
 End of event record. More...
 
struct  JHead
 Monte Carlo run header. More...
 
struct  JHeadHelper
 Auxiliary map of application to check method. More...
 
struct  JHeadWriter
 Implementation for Head output of JHead objects with ROOT dictionary. More...
 
struct  JOscFlux
 Implementation of oscillated neutrino flux. More...
 
struct  JParticle
 Auxiliary class to handle particle name, codes and mass. More...
 
struct  JPDB
 Collection of particles. More...
 
struct  JVolume
 Auxiliary class for histogramming of effective volume. More...
 

Typedefs

typedef double(* pEvtWeightFactor )(const Evt &)
 Type definition of event-weight factor pointer. More...
 
using pFlux = pEvtWeightFactor
 Type definition of flux function pointer. More...
 
typedef double(* pDiffuseFlux )(int, double, double)
 Type definition of pointer to diffuse flux function. More...
 
typedef JEvtWeightFactorHelper
< JFlux
JFluxHelper
 Type definition of event-weight factor helper for flux functions. More...
 
typedef JEvtWeightFactorHelper
< JDiffuseFlux
JDiffuseFluxHelper
 Type definition of event-weight factor helper for diffuse flux objects. More...
 
typedef
JEvtWeightFactorMultiParticle
< JFlux
JFluxMultiParticle
 Type-definition of multi-particle event-weight factor for fluxes. More...
 
typedef bool(* is_head )(const JHead &)
 Type definition of test function of header. 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_D_PLUS = 35, GEANT4_TYPE_D_MINUS = 36, GEANT4_TYPE_D_ZERO = 37, GEANT4_TYPE_ANTID_ZERO = 38,
  GEANT4_TYPE_D_PLUS_S = 39, GEANT4_TYPE_D_MINUS_S = 40, GEANT4_TYPE_LAMBDA_PLUS_C = 41, 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, GEANT4_TYPE_XI_PLUS_C = 0, GEANT4_TYPE_XI_MINUS_C = 0,
  GEANT4_TYPE_XI_ZERO_C = 0, GEANT4_TYPE_ANTIXI_ZERO_C = 0, GEANT4_TYPE_OMEGA_PLUS_B = 0, GEANT4_TYPE_ANTIOMEGA_ZERO_C = 0,
  GEANT4_TYPE_OMEGA_MINUS_B = 0, GEANT4_TYPE_OMEGA_ZERO_C = 0, GEANT4_TYPE_B_PLUS = 0, GEANT4_TYPE_B_MINUS = 0,
  GEANT4_TYPE_XI_MINUS_B = 0, GEANT4_TYPE_XI_PLUS_B = 0, GEANT4_TYPE_B_ZERO = 0, GEANT4_TYPE_ANTIB_ZERO = 0,
  GEANT4_TYPE_B_ZERO_S = 0, GEANT4_TYPE_ANTIB_ZERO_S = 0, GEANT4_TYPE_XI_ZERO_B = 0, GEANT4_TYPE_ANTIXI_ZERO_B = 0,
  GEANT4_TYPE_LAMBDA_B = 0, GEANT4_TYPE_ANTILAMBDA_B = 0, GEANT4_TYPE_B_PLUS_C = 0, GEANT4_TYPE_B_MINUS_C = 0
}
 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_D_PLUS = 411, TRACK_TYPE_D_ZERO = 421, TRACK_TYPE_D_PLUS_S = 431,
  TRACK_TYPE_B_ZERO = 511, TRACK_TYPE_B_PLUS = 521, TRACK_TYPE_B_ZERO_S = 531, TRACK_TYPE_B_PLUS_C = 541,
  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_LAMBDA_PLUS_C = 4122, TRACK_TYPE_XI_ZERO_C = 4132, TRACK_TYPE_XI_PLUS_C = 4232,
  TRACK_TYPE_OMEGA_ZERO_C = 4332, TRACK_TYPE_LAMBDA_B = 5122, TRACK_TYPE_XI_MINUS_B = 5132, TRACK_TYPE_XI_ZERO_B = 5232,
  TRACK_TYPE_OMEGA_MINUS_B = 5332, 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_D_MINUS = -411, TRACK_TYPE_D_MINUS_S = -431, TRACK_TYPE_ANTID_ZERO = -421,
  TRACK_TYPE_ANTIB_ZERO = -511, TRACK_TYPE_B_MINUS = -521, TRACK_TYPE_ANTIB_ZERO_S = -531, TRACK_TYPE_B_MINUS_C = -541,
  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, TRACK_TYPE_ANTIXI_ZERO_C = -4132, TRACK_TYPE_XI_MINUS_C = -4232, TRACK_TYPE_ANTIOMEGA_ZERO_C = -4332,
  TRACK_TYPE_ANTILAMBDA_B = -5122, TRACK_TYPE_XI_PLUS_B = -5132, TRACK_TYPE_ANTIXI_ZERO_B = -5232, TRACK_TYPE_OMEGA_PLUS_B = -5332
}
 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 &pos)
 Get position. More...
 
Vec getPosition (const JPosition3D &pos)
 Get position. More...
 
JPosition3D getPosition (const Trk &track)
 Get position. More...
 
JDirection3D getDirection (const Vec &dir)
 Get direction. More...
 
Vec getDirection (const JDirection3D &dir)
 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...
 
double getW (const Trk &track, const size_t index, const double value)
 Get track information. 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_lepton (const Trk &track)
 Test whether given track is a lepton. More...
 
bool is_charged_lepton (const Trk &track)
 Test whether given track is a charged lepton. More...
 
bool is_hadron (const Trk &track)
 Test whether given track is a hadron. More...
 
bool is_meson (const Trk &track)
 Test whether given track is a meson (is hadron and third digit of PDG code is zero) More...
 
bool is_baryon (const Trk &track)
 Test whether given track is a baryon (is hadron and third digit of PDG code is not zero) More...
 
bool has_particleID (const Trk &track, const int type)
 Test whether given track corresponds to given particle type. More...
 
bool is_initialstate (const Trk &track)
 Test whether given track corresponds to an initial state particle. More...
 
bool is_finalstate (const Trk &track)
 Test whether given track corresponds to a final state particle. 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...
 
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...
 
double getKineticEnergy (const double E, const double m)
 Get kinetic energy of particle with given mass. More...
 
double getKineticEnergy (const Trk &trk)
 Get track kinetic energy. More...
 
double getE0 (const Evt &evt)
 Get initial state energy of a neutrino interaction. More...
 
double getE1 (const Evt &evt)
 Get final state energy of a neutrino interaction. More...
 
Vec getP0 (const Evt &evt)
 Get momentum vector of the initial state of a neutrino interaction. More...
 
Vec getP1 (const Evt &evt)
 Get momentum vector of the final state of a neutrino interaction. More...
 
template<class JFunction_t , class JEvtWeightFactor_t = JEvtWeightFactor>
JEvtWeightFactorFunction
< JFunction_t,
JEvtWeightFactor_t > 
make_weightFactor (const JFunction_t &function)
 Auxiliary method for creating an interface to an event-weight factor. More...
 
template<class JEvtWeightFactor_t = JEvtWeightFactor>
JEvtWeightFactorFunction
< pEvtWeightFactor,
JEvtWeightFactor_t > 
make_weightFactor (pEvtWeightFactor function)
 Auxiliary method for creating an interface to an event-weight factor. More...
 
template<class JFunction_t >
JEvtWeightFactorFunction
< JFunction_t, JFlux
make_fluxFunction (const JFunction_t &flux)
 Auxiliary method for creating an interface to a flux function. More...
 
JEvtWeightFactorFunction
< pFlux, JFlux
make_fluxFunction (pFlux flux)
 Auxiliary method for creating an interface to a flux function. More...
 
template<class JFunction_t >
JEvtWeightFactorFunction
< JFunction_t, JDiffuseFlux
make_diffuseFluxFunction (const JFunction_t &flux)
 Auxiliary method for creating an interface to a diffuse flux function. More...
 
JEvtWeightFactorFunction
< pDiffuseFlux, JDiffuseFlux
make_diffuseFluxFunction (pDiffuseFlux flux)
 Auxiliary method for creating an interface to a diffuse flux function. 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...
 
std::string getTag (const std::string &tag, const int counter)
 Get tag with 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_genhen (const JHead &header)
 Check for generator. More...
 
bool is_gseagen (const JHead &header)
 Check for generator. More...
 
bool is_mupage (const JHead &header)
 Check for generator. More...
 
bool is_corsika (const JHead &header)
 Check for generator. More...
 
bool is_km3buu (const JHead &header)
 Check for generator. 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_sirene (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
 AAnet shower fit reconstruction type. More...
 
static const JEvtEvaluator getEvtValue
 Function object for evaluation of DAQ objects. More...
 
JEvtWeighter getEventWeighter
 Function object for mapping header to event weighter. More...
 
static const char AANET_TAG_SEPARATOR = '_'
 Separator for tag extension of multiple tags in Head ("_<counter>"). More...
 
static JHeadHelper get_is_head
 Function object to get check method for given application. More...
 

Detailed Description

Extensions to Evt data format.

Author
bjung
mlincetto
mdejong

Typedef Documentation

typedef double(* JAANET::pEvtWeightFactor)(const Evt &)

Type definition of event-weight factor pointer.

Definition at line 105 of file JEvtWeightFactorFunction.hh.

Type definition of flux function pointer.

Definition at line 109 of file JEvtWeightFactorFunction.hh.

typedef double(* JAANET::pDiffuseFlux)(int, double, double)

Type definition of pointer to diffuse flux function.

Definition at line 150 of file JEvtWeightFactorFunction.hh.

Type definition of event-weight factor helper for flux functions.

Definition at line 173 of file JEvtWeightFactorHelper.hh.

Type definition of event-weight factor helper for diffuse flux objects.

Definition at line 176 of file JEvtWeightFactorHelper.hh.

Type-definition of multi-particle event-weight factor for fluxes.

Definition at line 163 of file JEvtWeightFactorMultiParticle.hh.

typedef bool(* JAANET::is_head)(const JHead &)

Type definition of test function of header.

Definition at line 27 of file JHeadToolkit.hh.

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 73 of file JAAnetToolkit.hh.

73  {
74  HIT_TYPE_MUON_DIRECT = +5, //!< Direct light from muon
75  HIT_TYPE_MUON_SCATTERED = -5, //!< Scattered light from muon
76  HIT_TYPE_DELTARAYS_DIRECT = +4, //!< Direct light from delta-rays
77  HIT_TYPE_DELTARAYS_SCATTERED = -4, //!< Scattered light from delta-rays
78  HIT_TYPE_BREMSSTRAHLUNG_DIRECT = +3, //!< Direct light from Bremsstrahlung
79  HIT_TYPE_BREMSSTRAHLUNG_SCATTERED = -3, //!< Scattered light from Bremsstrahlung
80  HIT_TYPE_SHOWER_DIRECT = +99, //!< Direct light from primary shower
81  HIT_TYPE_SHOWER_SCATTERED = -99, //!< Scattered light from primary shower
82  HIT_TYPE_NOISE = -1, //!< Random noise
83  HIT_TYPE_UNKNOWN = 0 //!< Unknown source
84  };
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_D_PLUS 
GEANT4_TYPE_D_MINUS 
GEANT4_TYPE_D_ZERO 
GEANT4_TYPE_ANTID_ZERO 
GEANT4_TYPE_D_PLUS_S 
GEANT4_TYPE_D_MINUS_S 
GEANT4_TYPE_LAMBDA_PLUS_C 
GEANT4_TYPE_DEUTERON 
GEANT4_TYPE_TRITON 
GEANT4_TYPE_ALPHA 
GEANT4_TYPE_GEANTINO 
GEANT4_TYPE_HE3 
GEANT4_TYPE_ANTITAU 
GEANT4_TYPE_TAU 
GEANT4_TYPE_XI_PLUS_C 
GEANT4_TYPE_XI_MINUS_C 
GEANT4_TYPE_XI_ZERO_C 
GEANT4_TYPE_ANTIXI_ZERO_C 
GEANT4_TYPE_OMEGA_PLUS_B 
GEANT4_TYPE_ANTIOMEGA_ZERO_C 
GEANT4_TYPE_OMEGA_MINUS_B 
GEANT4_TYPE_OMEGA_ZERO_C 
GEANT4_TYPE_B_PLUS 
GEANT4_TYPE_B_MINUS 
GEANT4_TYPE_XI_MINUS_B 
GEANT4_TYPE_XI_PLUS_B 
GEANT4_TYPE_B_ZERO 
GEANT4_TYPE_ANTIB_ZERO 
GEANT4_TYPE_B_ZERO_S 
GEANT4_TYPE_ANTIB_ZERO_S 
GEANT4_TYPE_XI_ZERO_B 
GEANT4_TYPE_ANTIXI_ZERO_B 
GEANT4_TYPE_LAMBDA_B 
GEANT4_TYPE_ANTILAMBDA_B 
GEANT4_TYPE_B_PLUS_C 
GEANT4_TYPE_B_MINUS_C 

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,
50  GEANT4_TYPE_D_PLUS = 35,
52  GEANT4_TYPE_D_ZERO = 37,
58  GEANT4_TYPE_TRITON = 46,
59  GEANT4_TYPE_ALPHA = 47,
61  GEANT4_TYPE_HE3 = 49,
62  //
63  // %KM3NeT specific codes
64  //
66  GEANT4_TYPE_TAU = 34,
67  //
68  // %Particles with PDG code but not GEANT
69  //

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_D_PLUS 
TRACK_TYPE_D_ZERO 
TRACK_TYPE_D_PLUS_S 
TRACK_TYPE_B_ZERO 
TRACK_TYPE_B_PLUS 
TRACK_TYPE_B_ZERO_S 
TRACK_TYPE_B_PLUS_C 
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_LAMBDA_PLUS_C 
TRACK_TYPE_XI_ZERO_C 
TRACK_TYPE_XI_PLUS_C 
TRACK_TYPE_OMEGA_ZERO_C 
TRACK_TYPE_LAMBDA_B 
TRACK_TYPE_XI_MINUS_B 
TRACK_TYPE_XI_ZERO_B 
TRACK_TYPE_OMEGA_MINUS_B 
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_D_MINUS 
TRACK_TYPE_D_MINUS_S 
TRACK_TYPE_ANTID_ZERO 
TRACK_TYPE_ANTIB_ZERO 
TRACK_TYPE_B_MINUS 
TRACK_TYPE_ANTIB_ZERO_S 
TRACK_TYPE_B_MINUS_C 
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 
TRACK_TYPE_ANTIXI_ZERO_C 
TRACK_TYPE_XI_MINUS_C 
TRACK_TYPE_ANTIOMEGA_ZERO_C 
TRACK_TYPE_ANTILAMBDA_B 
TRACK_TYPE_XI_PLUS_B 
TRACK_TYPE_ANTIXI_ZERO_B 
TRACK_TYPE_OMEGA_PLUS_B 

Definition at line 96 of file JParticleTypes.hh.

96  { TRACK_TYPE_ELECTRON = 11,
97  TRACK_TYPE_NUE = 12,
98  TRACK_TYPE_MUON = 13,
99  TRACK_TYPE_NUMU = 14,
100  TRACK_TYPE_TAU = 15,
101  TRACK_TYPE_NUTAU = 16,
102  TRACK_TYPE_PHOTON = 22,
105  TRACK_TYPE_PION_PLUS = 211,
106  TRACK_TYPE_K_LONG = 130,
107  TRACK_TYPE_K_SHORT = 310,
108  TRACK_TYPE_K_PLUS = 321,
109  TRACK_TYPE_D_PLUS = 411,
110  TRACK_TYPE_D_ZERO = 421,
111  TRACK_TYPE_D_PLUS_S = 431,
112  TRACK_TYPE_B_ZERO = 511,
113  TRACK_TYPE_B_PLUS = 521,
114  TRACK_TYPE_B_ZERO_S = 531,
115  TRACK_TYPE_B_PLUS_C = 541,
116  TRACK_TYPE_PROTON = 2212,
117  TRACK_TYPE_NEUTRON = 2112,
118  TRACK_TYPE_LAMBDA = 3122,
119  TRACK_TYPE_SIGMA_PLUS = 3222,
121  TRACK_TYPE_SIGMA_MINUS = 3112,
122  TRACK_TYPE_NEUTRAL_XI = 3322,
123  TRACK_TYPE_XI_MINUS = 3312,
124  TRACK_TYPE_OMEGA_MINUS = 3334,
126  TRACK_TYPE_XI_ZERO_C = 4132,
127  TRACK_TYPE_XI_PLUS_C = 4232,
129  TRACK_TYPE_LAMBDA_B = 5122,
130  TRACK_TYPE_XI_MINUS_B = 5132,
131  TRACK_TYPE_XI_ZERO_B = 5232,
133 
135  TRACK_TYPE_ANTINUE = -12,
136  TRACK_TYPE_ANTIMUON = -13,
137  TRACK_TYPE_ANTINUMU = -14,
138  TRACK_TYPE_ANTITAU = -15,
139  TRACK_TYPE_ANTINUTAU = -16,
142  TRACK_TYPE_PION_MINUS = -211,
143  TRACK_TYPE_ANTIK_LONG = -130,
144  TRACK_TYPE_ANTIK_SHORT = -310,
145  TRACK_TYPE_K_MINUS = -321,
146  TRACK_TYPE_D_MINUS = -411,
147  TRACK_TYPE_D_MINUS_S = -431,
148  TRACK_TYPE_ANTID_ZERO = -421,
149  TRACK_TYPE_ANTIB_ZERO = -511,
150  TRACK_TYPE_B_MINUS = -521,
152  TRACK_TYPE_B_MINUS_C = -541,
153  TRACK_TYPE_ANTIPROTON = -2212,
154  TRACK_TYPE_ANTINEUTRON = -2112,
155  TRACK_TYPE_ANTILAMBDA = -3122,
160  TRACK_TYPE_ANTIXI_MINUS = -3312,
162  TRACK_TYPE_ANTIXI_ZERO_C = -4132,
163  TRACK_TYPE_XI_MINUS_C = -4232,
165  TRACK_TYPE_ANTILAMBDA_B = -5122,
166  TRACK_TYPE_XI_PLUS_B = -5132,
167  TRACK_TYPE_ANTIXI_ZERO_B = -5232,
168  TRACK_TYPE_OMEGA_PLUS_B = -5332 };

Function Documentation

double JAANET::getTime ( const Hit hit)
inline

Get true time of hit.

Parameters
hithit
Returns
true time of hit

Definition at line 93 of file JAAnetToolkit.hh.

94  {
95  return hit.t;
96  }
double t
hit time (from tdc+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  return hit.a;
108  }
double a
hit amplitude (in p.e.)
Definition: Hit.hh:24
bool JAANET::is_noise ( const Hit hit)
inline

Verify hit origin.

Parameters
hithit
Returns
true if noise; else false

Definition at line 117 of file JAAnetToolkit.hh.

118  {
119  return hit.type == HIT_TYPE_NOISE;
120  }
int type
particle type or parametrisation used for hit (mc only)
Definition: Hit.hh:28
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 130 of file JAAnetToolkit.hh.

131  {
132  JTimeRange time_range(JTimeRange::DEFAULT_RANGE);
133 
134  for (std::vector<Hit>::const_iterator hit = event.mc_hits.begin(); hit != event.mc_hits.end(); ++hit) {
135  if (!is_noise(*hit)) {
136  time_range.include(getTime(*hit));
137  }
138  }
139 
140  return time_range;
141  }
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:48
JTOOLS::JRange< double > JTimeRange
Type definition for time range (unit [s]).
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 153 of file JAAnetToolkit.hh.

154  {
155  JTimeRange time_range = getTimeRange(event);
156 
157  if (!time_range.is_valid()) {
158  time_range = T_ns;
159  }
160 
161  if (time_range.getLength() < T_ns.getLength()) {
162  time_range.setLength(T_ns.getLength());
163  }
164 
165  return time_range;
166  }
JTimeRange getTimeRange(const Evt &event)
Get time range (i.e. time between earliest and latest hit) of Monte Carlo event.
JTOOLS::JRange< double > JTimeRange
Type definition for time range (unit [s]).
JPosition3D JAANET::getPosition ( const Vec pos)
inline

Get position.

Parameters
posposition
Returns
position

Definition at line 175 of file JAAnetToolkit.hh.

176  {
177  return JPosition3D(pos.x, pos.y, pos.z);
178  }
double z
Definition: Vec.hh:14
double y
Definition: Vec.hh:14
double x
Definition: Vec.hh:14
Vec JAANET::getPosition ( const JPosition3D &  pos)
inline

Get position.

Parameters
posposition
Returns
position

Definition at line 187 of file JAAnetToolkit.hh.

188  {
189  return Vec(pos.getX(), pos.getY(), pos.getZ());
190  }
The Vec class is a straightforward 3-d vector, which also works in pyroot.
Definition: Vec.hh:12
JPosition3D JAANET::getPosition ( const Trk track)
inline

Get position.

Parameters
tracktrack
Returns
position

Definition at line 199 of file JAAnetToolkit.hh.

200  {
201  return getPosition(track.pos);
202  }
JPosition3D getPosition(const Vec &pos)
Get position.
Vec pos
postion [m] of the track at time t
Definition: Trk.hh:17
JDirection3D JAANET::getDirection ( const Vec dir)
inline

Get direction.

Parameters
dirdirection
Returns
direction

Definition at line 211 of file JAAnetToolkit.hh.

212  {
213  return JDirection3D(dir.x, dir.y, dir.z);
214  }
double z
Definition: Vec.hh:14
double y
Definition: Vec.hh:14
double x
Definition: Vec.hh:14
Vec JAANET::getDirection ( const JDirection3D &  dir)
inline

Get direction.

Parameters
dirdirection
Returns
direction

Definition at line 223 of file JAAnetToolkit.hh.

224  {
225  return Vec(dir.getDX(), dir.getDY(), dir.getDZ());
226  }
The Vec class is a straightforward 3-d vector, which also works in pyroot.
Definition: Vec.hh:12
JDirection3D JAANET::getDirection ( const Trk track)
inline

Get direction.

Parameters
tracktrack
Returns
direction

Definition at line 235 of file JAAnetToolkit.hh.

236  {
237  return getDirection(track.dir);
238  }
Vec dir
track direction
Definition: Trk.hh:18
JDirection3D getDirection(const Vec &dir)
Get direction.
JAxis3D JAANET::getAxis ( const Trk track)
inline

Get axis.

Parameters
tracktrack
Returns
axis

Definition at line 247 of file JAAnetToolkit.hh.

248  {
249  return JAxis3D(getPosition(track), getDirection(track));
250  }
JDirection3D getDirection(const Vec &dir)
Get direction.
JPosition3D getPosition(const Vec &pos)
Get position.
JTrack3E JAANET::getTrack ( const Trk track)
inline

Get track.

Parameters
tracktrack
Returns
track

Definition at line 259 of file JAAnetToolkit.hh.

260  {
261  return JTrack3E(JTrack3D(getAxis(track), track.t), track.E);
262  }
double t
track time [ns] (when the particle is at pos )
Definition: Trk.hh:19
double E
Energy [GeV] (either MC truth or reconstructed)
Definition: Trk.hh:20
JAxis3D getAxis(const Trk &track)
Get axis.
JVertex3D JAANET::getVertex ( const Trk track)
inline

Get vertex.

Parameters
tracktrack
Returns
vertex

Definition at line 271 of file JAAnetToolkit.hh.

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

Get transformation.

Parameters
tracktrack
Returns
transformation

Definition at line 283 of file JAAnetToolkit.hh.

284  {
285  return JTransformation3D(getPosition(track), getDirection(track));
286  }
JDirection3D getDirection(const Vec &dir)
Get direction.
JPosition3D getPosition(const Vec &pos)
Get position.
double JAANET::getW ( const Trk track,
const size_t  index,
const double  value 
)
inline

Get track information.

Parameters
tracktrack
indexindex
valuedefault value
Returns
actual value

Definition at line 297 of file JAAnetToolkit.hh.

298  {
299  if (index < track.fitinf.size())
300  return track.fitinf[index];
301  else
302  return value;
303  }
std::vector< double > fitinf
place to store additional fit info, see km3net-dataformat/definitions/fitparameters.csv
Definition: Trk.hh:32
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 312 of file JAAnetToolkit.hh.

312  { return ( track.type == TRACK_TYPE_PHOTON ||
313  abs(track.type) == TRACK_TYPE_NEUTRAL_PION); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
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 321 of file JAAnetToolkit.hh.

321  { return (abs(track.type) == TRACK_TYPE_NUE ||
322  abs(track.type) == TRACK_TYPE_NUMU ||
323  abs(track.type) == TRACK_TYPE_NUTAU); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
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 331 of file JAAnetToolkit.hh.

331 { return (abs(track.type) == TRACK_TYPE_ELECTRON); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
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 339 of file JAAnetToolkit.hh.

339 { return (abs(track.type) == TRACK_TYPE_MUON); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
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 347 of file JAAnetToolkit.hh.

347 { return (abs(track.type) == TRACK_TYPE_TAU); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
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 355 of file JAAnetToolkit.hh.

355 { return (abs(track.type) == TRACK_TYPE_CHARGED_PION_PLUS); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
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 363 of file JAAnetToolkit.hh.

363 { return (abs(track.type) == TRACK_TYPE_PROTON); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
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 371 of file JAAnetToolkit.hh.

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

Test whether given track is a lepton.

Parameters
tracktrack
Returns
true if lepton; else fails

Definition at line 379 of file JAAnetToolkit.hh.

379  { return (is_neutrino(track) ||
380  is_electron(track) ||
381  is_muon (track) ||
382  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::is_charged_lepton ( const Trk track)
inline

Test whether given track is a charged lepton.

Parameters
tracktrack
Returns
true if lepton; else fails

Definition at line 390 of file JAAnetToolkit.hh.

390  { return (is_electron(track) ||
391  is_muon (track) ||
392  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_tau(const Trk &track)
Test whether given track is a (anti-)tau.
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 400 of file JAAnetToolkit.hh.

400  { return (abs(track.type) != TRACK_TYPE_PHOTON &&
401  !is_lepton(track)); }
bool is_lepton(const Trk &track)
Test whether given track is a lepton.
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
bool JAANET::is_meson ( const Trk track)
inline

Test whether given track is a meson (is hadron and third digit of PDG code is zero)

Parameters
tracktrack
Returns
true if hadron; else fails

Definition at line 409 of file JAAnetToolkit.hh.

409  { return (is_hadron(track) &&
410  ((int)(track.type / 1000)) % 10 == 0); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
bool is_hadron(const Trk &track)
Test whether given track is a hadron.
bool JAANET::is_baryon ( const Trk track)
inline

Test whether given track is a baryon (is hadron and third digit of PDG code is not zero)

Parameters
tracktrack
Returns
true if hadron; else fails

Definition at line 418 of file JAAnetToolkit.hh.

418  { return (is_hadron(track) &&
419  ((int)(track.type / 1000)) % 10 != 0); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
bool is_hadron(const Trk &track)
Test whether given track is a hadron.
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 428 of file JAAnetToolkit.hh.

428 { return (track.type == type); }
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
bool JAANET::is_initialstate ( const Trk track)
inline

Test whether given track corresponds to an initial state particle.

Parameters
tracktrack
Returns
true if track corresponds to an initial state particle; else false

Definition at line 436 of file JAAnetToolkit.hh.

437  {
438  if (Evt::ROOT_IO_VERSION >= 14) {
439 
440  return (track.status == TRK_ST_PRIMARYNEUTRINO ||
441  track.status == TRK_ST_PRIMARYCOSMIC ||
442  track.status == TRK_ST_MUONBUNDLE ||
443  track.status == TRK_ST_ININUCLEI);
444 
445  } else if (Evt::ROOT_IO_VERSION >= 11) {
446 
447  return (track.mother_id == TRK_MOTHER_UNDEFINED ||
448  track.mother_id == TRK_MOTHER_NONE);
449 
450  } else {
451 
452  return is_neutrino(track) && track.id == 0;
453  }
454  }
static const int TRK_MOTHER_NONE
mother id of a particle if it has no parent
Definition: trkmembers.hh:13
bool is_neutrino(const Trk &track)
Test whether given track is a neutrino.
int mother_id
MC id of the parent particle.
Definition: Trk.hh:29
static const int TRK_ST_PRIMARYNEUTRINO
initial state neutrino (&#39;neutrino&#39; tag in evt files from gseagen and genhen).
Definition: trkmembers.hh:16
static const int TRK_ST_MUONBUNDLE
initial state muon bundle (mupage)
Definition: trkmembers.hh:18
static const int TRK_ST_ININUCLEI
Initial state nuclei (gseagen)
Definition: trkmembers.hh:19
static const int TRK_ST_PRIMARYCOSMIC
initial state cosmic ray (&#39;track_primary&#39; tag in evt files from corant).
Definition: trkmembers.hh:17
static const int TRK_MOTHER_UNDEFINED
KM3NeT Data Definitions v3.0.0-3-gef79250 https://git.km3net.de/common/km3net-dataformat.
Definition: trkmembers.hh:12
int id
track identifier
Definition: Trk.hh:16
int status
MC status code, see km3net-dataformat/definitions/trkmembers.csv for values.
Definition: Trk.hh:28
static int ROOT_IO_VERSION
Streamer version as obtained from ROOT file.
Definition: Evt.hh:274
bool JAANET::is_finalstate ( const Trk track)
inline

Test whether given track corresponds to a final state particle.

Parameters
tracktrack
Returns
true if track corresponds to a final state particle; else false

Definition at line 462 of file JAAnetToolkit.hh.

463  {
464  if (Evt::ROOT_IO_VERSION >= 16) {
465  return track.status == TRK_ST_FINALSTATE;
466  } else {
467  return !is_initialstate(track);
468  }
469  }
bool is_initialstate(const Trk &track)
Test whether given track corresponds to an initial state particle.
int status
MC status code, see km3net-dataformat/definitions/trkmembers.csv for values.
Definition: Trk.hh:28
static const int TRK_ST_FINALSTATE
for MC: the particle must be processed by detector simulation (&#39;track_in&#39; tag in evt files)...
Definition: trkmembers.hh:15
static int ROOT_IO_VERSION
Streamer version as obtained from ROOT file.
Definition: Evt.hh:274
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 477 of file JAAnetToolkit.hh.

478  {
479  if (Evt::ROOT_IO_VERSION >= 14) {
480 
481  std::vector<Trk>::const_iterator i = find_if(evt.mc_trks.cbegin(), evt.mc_trks.cend(),
483  return i != evt.mc_trks.cend();
484 
485  } else {
486 
487  return !evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0]);
488  }
489  }
JPredicate< JResult_t T::*, JComparison::eq > make_predicate(JResult_t T::*member, const JResult_t value)
Helper method to create predicate for data member.
Definition: JPredicate.hh:128
bool is_neutrino(const Trk &track)
Test whether given track is a neutrino.
static const int TRK_ST_PRIMARYNEUTRINO
initial state neutrino (&#39;neutrino&#39; tag in evt files from gseagen and genhen).
Definition: trkmembers.hh:16
int status
MC status code, see km3net-dataformat/definitions/trkmembers.csv for values.
Definition: Trk.hh:28
static int ROOT_IO_VERSION
Streamer version as obtained from ROOT file.
Definition: Evt.hh:274
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:49
const Trk& JAANET::get_neutrino ( const Evt evt)
inline

Get incoming neutrino.

Parameters
evtevent
Returns
neutrino

Definition at line 497 of file JAAnetToolkit.hh.

498  {
499  if (Evt::ROOT_IO_VERSION >= 14) {
500 
501  std::vector<Trk>::const_iterator i = find_if(evt.mc_trks.cbegin(), evt.mc_trks.cend(),
503 
504  if (i != evt.mc_trks.cend()) { return *i; }
505 
506  } else if (!evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0])) {
507 
508  return evt.mc_trks[0];
509  }
510 
511  THROW(JIndexOutOfRange, "JAANET::get_neutrino(): Event " << evt.id << " has no neutrino.");
512  }
JPredicate< JResult_t T::*, JComparison::eq > make_predicate(JResult_t T::*member, const JResult_t value)
Helper method to create predicate for data member.
Definition: JPredicate.hh:128
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:696
bool is_neutrino(const Trk &track)
Test whether given track is a neutrino.
static const int TRK_ST_PRIMARYNEUTRINO
initial state neutrino (&#39;neutrino&#39; tag in evt files from gseagen and genhen).
Definition: trkmembers.hh:16
int status
MC status code, see km3net-dataformat/definitions/trkmembers.csv for values.
Definition: Trk.hh:28
int id
offline event identifier
Definition: Evt.hh:22
static int ROOT_IO_VERSION
Streamer version as obtained from ROOT file.
Definition: Evt.hh:274
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:49
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 520 of file JAAnetToolkit.hh.

521  {
522  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
523  }
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:49
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 531 of file JAAnetToolkit.hh.

532  {
533  return count_electrons(evt) != 0;
534  }
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 542 of file JAAnetToolkit.hh.

543  {
544  if (count_electrons(evt) > 0)
545  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
546  else
547  THROW(JIndexOutOfRange, "This event has no electron.");
548  }
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:696
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:49
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 559 of file JAAnetToolkit.hh.

560  {
561  return hit.type == GEANT4_TYPE_ELECTRON || hit.type == GEANT4_TYPE_ANTIELECTRON;
562  }
int type
particle type or parametrisation used for hit (mc only)
Definition: Hit.hh:28
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 570 of file JAAnetToolkit.hh.

571  {
572  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
573  }
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:49
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 581 of file JAAnetToolkit.hh.

582  {
583  return count_muons(evt) != 0;
584  }
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 592 of file JAAnetToolkit.hh.

593  {
594  if (count_muons(evt) > 0)
595  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
596  else
597  THROW(JIndexOutOfRange, "This event has no muon.");
598  }
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:696
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:49
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 609 of file JAAnetToolkit.hh.

610  {
611  return hit.type == GEANT4_TYPE_MUON || hit.type == GEANT4_TYPE_ANTIMUON;
612  }
int type
particle type or parametrisation used for hit (mc only)
Definition: Hit.hh:28
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 620 of file JAAnetToolkit.hh.

621  {
622  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
623  }
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:49
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 631 of file JAAnetToolkit.hh.

632  {
633  return count_taus(evt) != 0;
634  }
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 642 of file JAAnetToolkit.hh.

643  {
644  if (count_taus(evt) > 0)
645  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
646  else
647  THROW(JIndexOutOfRange, "This event has no tau.");
648  }
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:696
int count_taus(const Evt &evt)
Count the number of taus in a given event.
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:49
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 659 of file JAAnetToolkit.hh.

660  {
661  return hit.type == GEANT4_TYPE_TAU || hit.type == GEANT4_TYPE_ANTITAU;
662  }
int type
particle type or parametrisation used for hit (mc only)
Definition: Hit.hh:28
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 670 of file JAAnetToolkit.hh.

671  {
672  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_hadron);
673  }
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:49
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 684 of file JAAnetToolkit.hh.

685  {
686  return (!from_electron(hit) && !from_muon(hit) && !from_tau(hit));
687  }
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.
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 698 of file JAAnetToolkit.hh.

699  {
700  for (std::vector<Trk>::iterator p = evt.mc_trks.begin(); p != evt.mc_trks.end(); p++) {
701  p->t += tOff;
702  }
703  }
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:49
double JAANET::getKineticEnergy ( const double  E,
const double  m 
)
inline

Get kinetic energy of particle with given mass.

Parameters
Eenergy [GeV]
mmass [GeV]
Returns
energy [GeV]

Definition at line 713 of file JAAnetToolkit.hh.

714  {
715  if (E > m) {
716  return sqrt((E - m) * (E + m));
717  } else {
718  return 0.0;
719  }
720  }
then usage E
Definition: JMuonPostfit.sh:35
double JAANET::getKineticEnergy ( const Trk trk)
inline

Get track kinetic energy.

Parameters
trktrack
Returns
kinetic energy [GeV]

Definition at line 729 of file JAAnetToolkit.hh.

730  {
731  const double energy = trk.E;
732  const double mass = JPDB::getInstance().getPDG(trk.type).mass;
733 
734  return getKineticEnergy(energy, mass);
735  }
double E
Energy [GeV] (either MC truth or reconstructed)
Definition: Trk.hh:20
double getKineticEnergy(const double E, const double m)
Get kinetic energy of particle with given mass.
T & getInstance(const T &object)
Get static instance from temporary object.
Definition: JObject.hh:75
int type
MC: particle type in PDG encoding.
Definition: Trk.hh:24
then for APP in event gandalf start energy
Definition: JMuonMCEvt.sh:44
double JAANET::getE0 ( const Evt evt)
inline

Get initial state energy of a neutrino interaction.


This method includes baryon number conservation.

Parameters
evtevent
Returns
energy [GeV]

Definition at line 745 of file JAAnetToolkit.hh.

746  {
747  const Trk& neutrino = get_neutrino(evt);
748 
749  double E0 = neutrino.E;
750 
751  for (std::vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.end(); ++track) {
752 
753  if (!is_finalstate(*track)) { continue; }
754 
755  if (track->type == TRACK_TYPE_NEUTRON ||
756  track->type == TRACK_TYPE_SIGMA_MINUS) {
757 
758  E0 += MASS_NEUTRON;
759 
760  } else if (track->type == TRACK_TYPE_PROTON ||
761  track->type == TRACK_TYPE_LAMBDA ||
762  track->type == TRACK_TYPE_SIGMA_PLUS) {
763 
764  E0 += MASS_PROTON;
765 
766  } else if (track->type == TRACK_TYPE_ANTINEUTRON) {
767 
768  E0 += MASS_NEUTRON;
769 
770  } else if (track->type == TRACK_TYPE_ANTIPROTON ||
771  track->type == TRACK_TYPE_ANTILAMBDA) {
772 
773  E0 += MASS_PROTON;
774  }
775  }
776 
777  return E0;
778  }
double E
Energy [GeV] (either MC truth or reconstructed)
Definition: Trk.hh:20
static const double MASS_NEUTRON
neutron mass [GeV]
bool is_finalstate(const Trk &track)
Test whether given track corresponds to a final state particle.
static const double MASS_PROTON
proton mass [GeV]
const Trk & get_neutrino(const Evt &evt)
Get incoming neutrino.
The Trk class represents a Monte Carlo (MC) particle as well as a reconstructed track/shower.
Definition: Trk.hh:14
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:49
double JAANET::getE1 ( const Evt evt)
inline

Get final state energy of a neutrino interaction.


This method includes muon energy loss.

Parameters
evtevent
Returns
energy [GeV]

Definition at line 788 of file JAAnetToolkit.hh.

789  {
790  double E1 = 0.0;
791 
792  const Trk& neutrino = get_neutrino(evt);
793 
794  for (std::vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.cend(); ++track) {
795 
796  if (!is_finalstate(*track)) { continue; }
797 
798  if (is_muon(*track)) {
799 
800  const Trk& mother = ( track->mother_id > TRK_MOTHER_UNDEFINED ?
801  evt.mc_trks[track->mother_id] : neutrino );
802 
803  const double distance = (track->pos - mother.pos).len();
804 
805  E1 += JPHYSICS::gWater(track->E, -distance);
806 
807  } else {
808 
809  E1 += track->E;
810  }
811  }
812 
813  return E1;
814  }
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
static const JGeaneWater gWater
Function object for energy loss of muon in sea water.
Definition: JGeane.hh:381
bool is_finalstate(const Trk &track)
Test whether given track corresponds to a final state particle.
static const int TRK_MOTHER_UNDEFINED
KM3NeT Data Definitions v3.0.0-3-gef79250 https://git.km3net.de/common/km3net-dataformat.
Definition: trkmembers.hh:12
Vec pos
postion [m] of the track at time t
Definition: Trk.hh:17
const Trk & get_neutrino(const Evt &evt)
Get incoming neutrino.
The Trk class represents a Monte Carlo (MC) particle as well as a reconstructed track/shower.
Definition: Trk.hh:14
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:49
Vec JAANET::getP0 ( const Evt evt)
inline

Get momentum vector of the initial state of a neutrino interaction.


This method assumes neutrino DIS on a stationary nucleus

Parameters
evtevent
Returns
energy [GeV]

Definition at line 824 of file JAAnetToolkit.hh.

825  {
826  const Trk& neutrino = get_neutrino(evt);
827 
828  return neutrino.E * neutrino.dir;
829  }
Vec dir
track direction
Definition: Trk.hh:18
double E
Energy [GeV] (either MC truth or reconstructed)
Definition: Trk.hh:20
const Trk & get_neutrino(const Evt &evt)
Get incoming neutrino.
The Trk class represents a Monte Carlo (MC) particle as well as a reconstructed track/shower.
Definition: Trk.hh:14
Vec JAANET::getP1 ( const Evt evt)
inline

Get momentum vector of the final state of a neutrino interaction.


This method includes muon energy losses.

Parameters
evtevent
Returns
final state momentum vector [GeV]

Definition at line 839 of file JAAnetToolkit.hh.

840  {
841  Vec P1(0,0,0);
842 
843  const Trk& neutrino = get_neutrino(evt);
844 
845  for (std::vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.cend(); ++track) {
846 
847  if (!is_finalstate(*track)) { continue; }
848 
849  double kineticEnergy = 0.0;
850 
851  if (is_muon(*track)) {
852 
853  const Trk& mother = ( track->mother_id > TRK_MOTHER_UNDEFINED ?
854  evt.mc_trks[track->mother_id] : neutrino );
855 
856  const double distance = (track->pos - mother.pos).len();
857  const double energy = JPHYSICS::gWater(track->E, -distance);
858 
859  kineticEnergy = getKineticEnergy(energy, JPHYSICS::MASS_MUON);
860 
861  } else {
862 
863  kineticEnergy = getKineticEnergy(*track);
864  }
865 
866  P1 += kineticEnergy * track->dir;
867  }
868 
869  return P1;
870  }
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
static const double MASS_MUON
muon mass [GeV]
static const JGeaneWater gWater
Function object for energy loss of muon in sea water.
Definition: JGeane.hh:381
bool is_finalstate(const Trk &track)
Test whether given track corresponds to a final state particle.
The Vec class is a straightforward 3-d vector, which also works in pyroot.
Definition: Vec.hh:12
double getKineticEnergy(const double E, const double m)
Get kinetic energy of particle with given mass.
static const int TRK_MOTHER_UNDEFINED
KM3NeT Data Definitions v3.0.0-3-gef79250 https://git.km3net.de/common/km3net-dataformat.
Definition: trkmembers.hh:12
then for APP in event gandalf start energy
Definition: JMuonMCEvt.sh:44
Vec pos
postion [m] of the track at time t
Definition: Trk.hh:17
then P1
const Trk & get_neutrino(const Evt &evt)
Get incoming neutrino.
The Trk class represents a Monte Carlo (MC) particle as well as a reconstructed track/shower.
Definition: Trk.hh:14
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:49
template<class JFunction_t , class JEvtWeightFactor_t = JEvtWeightFactor>
JEvtWeightFactorFunction<JFunction_t, JEvtWeightFactor_t> JAANET::make_weightFactor ( const JFunction_t &  function)
inline

Auxiliary method for creating an interface to an event-weight factor.

Parameters
functionfunction object
Returns
event-weight factor interface

Definition at line 199 of file JEvtWeightFactorFunction.hh.

199  {
200  return JEvtWeightFactorFunction<JFunction_t, JEvtWeightFactor_t>(function);
201  }
template<class JEvtWeightFactor_t = JEvtWeightFactor>
JEvtWeightFactorFunction<pEvtWeightFactor, JEvtWeightFactor_t> JAANET::make_weightFactor ( pEvtWeightFactor  function)
inline

Auxiliary method for creating an interface to an event-weight factor.

Parameters
functionfunction pointer
Returns
event-weight factor interface

Definition at line 211 of file JEvtWeightFactorFunction.hh.

211  {
212  return JEvtWeightFactorFunction<pEvtWeightFactor, JEvtWeightFactor_t>(function);
213  }
template<class JFunction_t >
JEvtWeightFactorFunction<JFunction_t, JFlux> JAANET::make_fluxFunction ( const JFunction_t &  flux)
inline

Auxiliary method for creating an interface to a flux function.

Parameters
fluxflux function object
Returns
flux function interface

Definition at line 223 of file JEvtWeightFactorFunction.hh.

223  {
224  return JEvtWeightFactorFunction<JFunction_t, JFlux>(flux);
225  }
JEvtWeightFactorFunction<pFlux, JFlux> JAANET::make_fluxFunction ( pFlux  flux)
inline

Auxiliary method for creating an interface to a flux function.

Parameters
fluxflux function pointer
Returns
flux function interface

Definition at line 234 of file JEvtWeightFactorFunction.hh.

234  {
235  return JEvtWeightFactorFunction<pFlux, JFlux>(flux);
236  }
template<class JFunction_t >
JEvtWeightFactorFunction<JFunction_t, JDiffuseFlux> JAANET::make_diffuseFluxFunction ( const JFunction_t &  flux)
inline

Auxiliary method for creating an interface to a diffuse flux function.

Parameters
fluxdiffuse flux function object
Returns
diffuse flux function interface

Definition at line 246 of file JEvtWeightFactorFunction.hh.

246  {
247  return JEvtWeightFactorFunction<JFunction_t, JDiffuseFlux>(flux);
248  }
JEvtWeightFactorFunction<pDiffuseFlux, JDiffuseFlux> JAANET::make_diffuseFluxFunction ( pDiffuseFlux  flux)
inline

Auxiliary method for creating an interface to a diffuse flux function.

Parameters
fluxdiffuse flux function pointer
Returns
diffuse flux function interface

Definition at line 257 of file JEvtWeightFactorFunction.hh.

257  {
258  return JEvtWeightFactorFunction<pDiffuseFlux, JDiffuseFlux>(flux);
259  }
void JAANET::copy ( const Head from,
JHead &  to 
)

Copy header from from to to.

Parameters
fromheader
toheader

Definition at line 162 of file JHead.cc.

163  {
164  using namespace std;
165  using namespace JPP;
166 
167  JRootReader reader(null, JHead::getEquationParameters(), JAAnetDictionary::getInstance());
168 
169  JRootReadableClass cls(to);
170 
171  for (Head::const_iterator i = from.begin(); i != from.end(); ++i) {
172 
173  const JRootReadableClass& abc = cls.find(getTag(i->first).c_str());
174 
175  const string buffer = trim(i->second);
176 
177  if (abc.is_valid() && buffer != "") {
178 
179  JRedirectString redirect(reader, buffer);
180 
181  reader.getObject(abc);
182 
183  if (i->first == getTag(i->first)) {
184  to.insert(*i); // keep track of parsed tags
185  }
186 
187  } else {
188 
189  to.insert(*i); // store data of unknown tags
190  }
191  }
192  }
const char * c_str() const
C-string.
Definition: JTag.hh:201
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 201 of file JHead.cc.

202  {
203  using namespace std;
204  using namespace JPP;
205 
206  to = Head();
207 
208  JHeadWriter writer(to, JHead::getEquationParameters(), JAAnetDictionary::getInstance());
209 
210  JRootWritableClass cls(from);
211 
212  unique_ptr<TIterator> i(cls.getClass()->GetListOfDataMembers()->MakeIterator());
213 
214  for (const TDataMember* p; (p = (const TDataMember*) i->Next()) != NULL; ) {
215  if (!JRootClass::is_static(*p)) {
216  if (from.find(p->GetName()) != from.end() ||
217  cls.get(*p) == JRootClass(&JHead::start_run) ||
218  cls.get(*p) == JRootClass(&JHead::end_event)) {
219  writer.put(p->GetName(), cls.get(*p), true);
220  }
221  }
222  }
223 
224  // copy pending data
225 
226  for (JHead::const_iterator i = from.begin(); i != from.end(); ++i) {
227  if (to.find(i->first) == to.end()) {
228  to.insert(*i);
229  }
230  }
231  }
T & getInstance(const T &object)
Get static instance from temporary object.
Definition: JObject.hh:75
The Head class reflects the header of Monte-Carlo event files, which consists of keys (also referred ...
Definition: Head.hh:65
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 94 of file JHead.hh.

95  {
96  using namespace std;
97 
98  const string::size_type pos = tag.find(AANET_TAG_SEPARATOR);
99 
100  if (pos != string::npos) {
101 
102  for (string::size_type i = pos + 1; i != tag.size(); ++i) {
103  if (!isdigit(tag[i])) {
104  return tag;
105  }
106  }
107 
108  return tag.substr(0, pos);
109  }
110 
111  return tag;
112  }
static const char AANET_TAG_SEPARATOR
Separator for tag extension of multiple tags in Head (&quot;_&lt;counter&gt;&quot;).
Definition: JHead.hh:65
std::string JAANET::getTag ( const std::string tag,
const int  counter 
)
inline

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

Parameters
tagtag
countercounter
Returns
tag

Definition at line 122 of file JHead.hh.

123  {
124  std::ostringstream os;
125 
126  os << tag << AANET_TAG_SEPARATOR << counter;
127 
128  return os.str();
129  }
static const char AANET_TAG_SEPARATOR
Separator for tag extension of multiple tags in Head (&quot;_&lt;counter&gt;&quot;).
Definition: JHead.hh:65
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 1725 of file JHead.hh.

1727  {
1728  return JHead(first).match(JHead(second));
1729  }
bool match(const JHead &header) const
Test match of headers.
Definition: JHead.hh:1403
Monte Carlo run header.
Definition: JHead.hh:1167
bool JAANET::operator< ( const Head first,
const Head second 
)
inline

Less than operator.

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

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

Definition at line 1741 of file JHead.hh.

1743  {
1744  return JHead(first).less(JHead(second));
1745  }
bool less(const JHead &header) const
Comparison of headers.
Definition: JHead.hh:1415
Monte Carlo run header.
Definition: JHead.hh:1167
bool JAANET::is_genhen ( const JHead &  header)
inline

Check for generator.

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

Definition at line 36 of file JHeadToolkit.hh.

37  {
38  for (const auto& i : header.simul) {
39  if (i.program == APPLICATION_GENHEN) {
40  return true;
41  }
42  }
43 
44  for (const auto& i : header.physics) { // legacy
45  if (i.buffer.find(APPLICATION_GENHEN) != std::string::npos) {
46  return true;
47  }
48  }
49 
50  return false;
51  }
static const char *const APPLICATION_GENHEN
KM3NeT Data Definitions v3.0.0-3-gef79250 https://git.km3net.de/common/km3net-dataformat.
Definition: applications.hh:12
bool JAANET::is_gseagen ( const JHead &  header)
inline

Check for generator.

Parameters
headerheader
Returns
true if this header is produced by gSeaGen; else false

Definition at line 60 of file JHeadToolkit.hh.

61  {
62  for (const auto& i : header.simul) {
63  if (i.program == APPLICATION_GSEAGEN) {
64  return true;
65  }
66  }
67 
68  for (const auto& i : header.physics) { // legacy
69  if (i.buffer.find(APPLICATION_GSEAGEN) != std::string::npos) {
70  return true;
71  }
72  }
73 
74  return false;
75  }
static const char *const APPLICATION_GSEAGEN
event generator
Definition: applications.hh:13
bool JAANET::is_mupage ( const JHead &  header)
inline

Check for generator.

Parameters
headerheader
Returns
true if this header is produced by MUPAGE; else false

Definition at line 84 of file JHeadToolkit.hh.

85  {
86  for (const auto& i : header.simul) {
87  if (i.program == APPLICATION_MUPAGE) {
88  return true;
89  }
90  }
91 
92  return false;
93  }
static const char *const APPLICATION_MUPAGE
event generator
Definition: applications.hh:14
bool JAANET::is_corsika ( const JHead &  header)
inline

Check for generator.

Parameters
headerheader
Returns
true if this header is produced by Corsika; else false

Definition at line 102 of file JHeadToolkit.hh.

103  {
104  for (const auto& i : header.simul) {
105  if (i.program == APPLICATION_CORSIKA) {
106  return true;
107  }
108  }
109 
110  return false;
111  }
static const char *const APPLICATION_CORSIKA
event generator
Definition: applications.hh:15
bool JAANET::is_km3buu ( const JHead &  header)
inline

Check for generator.

Parameters
headerheader
Returns
true if this header is produced by km3buu; else false

Definition at line 120 of file JHeadToolkit.hh.

121  {
122  for (const auto& i : header.simul) {
123  if (i.program == APPLICATION_KM3BUU) {
124  return true;
125  }
126  }
127 
128  return false;
129  }
static const char *const APPLICATION_KM3BUU
event generator
Definition: applications.hh:16
bool JAANET::is_km3 ( const JHead &  header)
inline

Check for detector simulation.

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

Definition at line 138 of file JHeadToolkit.hh.

139  {
140  for (const auto& i : header.simul) {
141  if (i.program == APPLICATION_KM3) {
142  return true;
143  }
144  }
145 
146  return false;
147  }
static const char *const APPLICATION_KM3
detector simulation
Definition: applications.hh:17
bool JAANET::is_km3sim ( const JHead &  header)
inline

Check for detector simulation.

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

Definition at line 156 of file JHeadToolkit.hh.

157  {
158  for (const auto& i : header.simul) {
159  if (i.program == APPLICATION_KM3SIM) {
160  return true;
161  }
162  }
163 
164  return false;
165  }
static const char *const APPLICATION_KM3SIM
detector simulation
Definition: applications.hh:18
bool JAANET::is_sirene ( const JHead &  header)
inline

Check for detector simulation.

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

Definition at line 174 of file JHeadToolkit.hh.

175  {
176  for (const auto& i : header.simul) {
177  if (i.program == APPLICATION_JSIRENE) {
178  return true;
179  }
180  }
181 
182  return false;
183  }
static const char *const APPLICATION_JSIRENE
detector simulation
Definition: applications.hh:19
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 192 of file JHeadToolkit.hh.

193  {
194  return header.DAQ.livetime_s > 0.0;
195  }
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 260 of file JHeadToolkit.hh.

261  {
262  if (is_sirene(header) || is_km3(header)) {
263 
264  return header.coord_origin;
265 
266  } else {
267 
268  if (header.is_valid(&JHead::can))
269  return Vec(0.0, 0.0, -header.can.zmin);
270  else if (header.is_valid(&JHead::coord_origin))
271  return header.coord_origin;
272  else
273  return Vec(0.0, 0.0, 0.0);
274  }
275  }
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 260 of file JHeadToolkit.hh.

261  {
262  if (is_sirene(header) || is_km3(header)) {
263 
264  return header.coord_origin;
265 
266  } else {
267 
268  if (header.is_valid(&JHead::can))
269  return Vec(0.0, 0.0, -header.can.zmin);
270  else if (header.is_valid(&JHead::coord_origin))
271  return header.coord_origin;
272  else
273  return Vec(0.0, 0.0, 0.0);
274  }
275  }
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 285 of file JHeadToolkit.hh.

286  {
287  const Vec pos = get<Vec>(header);
288 
289  return JPosition3D(pos.x, pos.y, pos.z);
290  }
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 300 of file JHeadToolkit.hh.

301  {
302  using namespace JGEOMETRY2D;
303 
304  return JCylinder3D(JCircle2D(JVector2D(), header.can.r),
305  header.can.zmin,
306  header.can.zmax);
307  }
Data structure for vector in two dimensions.
Definition: JVector2D.hh:32
Data structure for circle in two dimensions.
Definition: JCircle2D.hh:33

Variable Documentation

const int JAANET::AASHOWER_RECONSTRUCTION_TYPE = 101
static

AAnet shower fit reconstruction type.

Definition at line 68 of file JAAnetToolkit.hh.

const JEvtEvaluator JAANET::getEvtValue
static

Function object for evaluation of DAQ objects.

Definition at line 48 of file JEvtEvaluator.hh.

JEvtWeighter JAANET::getEventWeighter

Function object for mapping header to event weighter.

Definition at line 5 of file JEvtWeightToolkit.cc.

const char JAANET::AANET_TAG_SEPARATOR = '_'
static

Separator for tag extension of multiple tags in Head ("_<counter>").

Definition at line 65 of file JHead.hh.

JHeadHelper JAANET::get_is_head
static
Initial value:
= {
std::make_pair(APPLICATION_GENHEN, is_genhen),
std::make_pair(APPLICATION_GSEAGEN, is_gseagen),
std::make_pair(APPLICATION_MUPAGE, is_mupage),
std::make_pair(APPLICATION_CORSIKA, is_corsika),
std::make_pair(APPLICATION_KM3BUU, is_km3buu),
std::make_pair(APPLICATION_KM3, is_km3),
std::make_pair(APPLICATION_KM3SIM, is_km3sim),
std::make_pair(APPLICATION_JSIRENE, is_sirene),
std::make_pair("DAQ", is_daq)
}
bool is_mupage(const JHead &header)
Check for generator.
Definition: JHeadToolkit.hh:84
bool is_km3buu(const JHead &header)
Check for generator.
static const char *const APPLICATION_KM3
detector simulation
Definition: applications.hh:17
bool is_gseagen(const JHead &header)
Check for generator.
Definition: JHeadToolkit.hh:60
static const char *const APPLICATION_GSEAGEN
event generator
Definition: applications.hh:13
static const char *const APPLICATION_JSIRENE
detector simulation
Definition: applications.hh:19
bool is_corsika(const JHead &header)
Check for generator.
static const char *const APPLICATION_CORSIKA
event generator
Definition: applications.hh:15
bool is_daq(const JHead &header)
Check for real data.
static const char *const APPLICATION_KM3BUU
event generator
Definition: applications.hh:16
static const char *const APPLICATION_KM3SIM
detector simulation
Definition: applications.hh:18
static const char *const APPLICATION_MUPAGE
event generator
Definition: applications.hh:14
bool is_genhen(const JHead &header)
Check for generator.
Definition: JHeadToolkit.hh:36
static const char *const APPLICATION_GENHEN
KM3NeT Data Definitions v3.0.0-3-gef79250 https://git.km3net.de/common/km3net-dataformat.
Definition: applications.hh:12
bool is_sirene(const JHead &header)
Check for detector simulation.
bool is_km3sim(const JHead &header)
Check for detector simulation.
bool is_km3(const JHead &header)
Check for detector simulation.

Function object to get check method for given application.

Definition at line 230 of file JHeadToolkit.hh.