Jpp  18.2.0-rc.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  JAtmosphericNeutrinoFlux
 Implementation of atmospheric neutrino flux using official KM3NeT atmospheric flux function. More...
 
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  JEvtWeightFactorGSeaGen
 Implementation of reweighting factor for simulated neutrino interactions according to a specifiable ROOT TFormula. 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  JEvtWeightMiscellaneous
 Implementation of event weighing for miscellaneous data
such as a merged offline file containing neutrinos and atmospheric muons. More...
 
struct  JEvtWeightMupage
 Implementation of event weighing for MUPAGE data. More...
 
struct  JEvtWeightNormalisation
 Auxiliary data structure for storing pairs of header UUIDs and event-weight normalisations. More...
 
struct  JEvtWeighter
 Look-up table for event weighters. More...
 
struct  JNeutrinoTypeCollection
 Auxiliary class for parsing a vector of neutrino PDG identifiers. More...
 
struct  JFluxMapParser
 Auxiliary class for parsing multiparticle fluxes. More...
 
struct  JFlatFlux
 Function object for constant fluxes. More...
 
struct  JFlux
 Low-level interface for retrieving the flux corresponding to a given event. 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  calibration
 Calibration. 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  K40
 Livetime of noise 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...
 
class  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  JLorentzBoost
 Auxiliary class for performing Lorentz boosts. More...
 
struct  JMultiHead
 Auxiliary data structure to store multiple headers
and bookkeep event-weight normalisations. 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  JPowerLawFlux
 Example function object for computing power-law flux. More...
 
struct  JVolume
 Auxiliary class for histogramming of effective volume. More...
 
struct  JEventSelector
 Event selector. 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_NEUTRAL_ANTISIGMA = -3212, TRACK_TYPE_ANTISIGMA_MINUS = -3112, TRACK_TYPE_NEUTRAL_ANTIXI = -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...
 
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...
 
JVertex3D getVertex (const Trk &track)
 Get vertex. More...
 
JVertex3D getVertex (const Evt &event)
 Get event vertex. 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...
 
double getInvariantMass (const Evt &event)
 Get final state invariant mass. More...
 
template<class JOscProbInterpolator_t = JOscProbInterpolator<>>
JEvtWeightFactorFunction
< JAtmosphericNeutrinoFlux,
JFlux
make_atmosphericNeutrinoFluxFunction (const std::string &oscProbTableFile, const JOscParameters &oscParameters)
 Auxiliary method for creating an interface to an atmospheric neutrino flux function
using an oscillation probability interpolation table. 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...
 
JEvtWeightFactorFunction
< JFlatFlux, JFlux
make_flatFluxFunction (const double flux)
 Auxiliary method for creating an interface to a flat 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...
 
JLorentzBoost getBoostToCOM (const Evt &event)
 Get Lorentz boost to the Center of Mass (COM) frame for a given neutrino interaction. More...
 
void boostToCOM (Evt &event)
 Boost event to the Center of Mass (COM) frame. More...
 
JEvtWeightFactorFunction
< JPowerLawFlux, JFlux
make_powerLawFluxFunction (const double normalisation, const double spectralIndex)
 Auxiliary method for creating an interface to a power-law flux function. More...
 
bool select (const Trk &trk, const Evt &evt)
 Event selection. More...
 
bool hasW (const Trk &trk, const int i)
 Check availability of value. More...
 
double getW (const Trk &trk, const int i)
 Get associated value. More...
 
double getW (const Trk &trk, const int i, const double value)
 Get associated value. More...
 
void setW (Trk &trk, const int i, const double value)
 Set associated value. 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, bjjung
mdejong
bjjung

Typedef Documentation

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

Type definition of event-weight factor pointer.

Definition at line 108 of file JEvtWeightFactorFunction.hh.

Type definition of flux function pointer.

Definition at line 112 of file JEvtWeightFactorFunction.hh.

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

Type definition of pointer to diffuse flux function.

Definition at line 153 of file JEvtWeightFactorFunction.hh.

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

Definition at line 185 of file JEvtWeightFactorHelper.hh.

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

Definition at line 188 of file JEvtWeightFactorHelper.hh.

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

Definition at line 169 of file JEvtWeightFactorMultiParticle.hh.

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

Type definition of test function of header.

Definition at line 28 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_NEUTRAL_ANTISIGMA 
TRACK_TYPE_ANTISIGMA_MINUS 
TRACK_TYPE_NEUTRAL_ANTIXI 
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.
JTransformation3D JAANET::getTransformation ( const Trk track)
inline

Get transformation.

Parameters
tracktrack
Returns
transformation

Definition at line 271 of file JAAnetToolkit.hh.

272  {
273  return JTransformation3D(getPosition(track), getDirection(track));
274  }
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 285 of file JAAnetToolkit.hh.

286  {
287  if (index < track.fitinf.size())
288  return track.fitinf[index];
289  else
290  return value;
291  }
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 300 of file JAAnetToolkit.hh.

300  { return ( track.type == TRACK_TYPE_PHOTON ||
301  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 309 of file JAAnetToolkit.hh.

309  { return (abs(track.type) == TRACK_TYPE_NUE ||
310  abs(track.type) == TRACK_TYPE_NUMU ||
311  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 319 of file JAAnetToolkit.hh.

319 { 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 327 of file JAAnetToolkit.hh.

327 { 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 335 of file JAAnetToolkit.hh.

335 { 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 343 of file JAAnetToolkit.hh.

343 { 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 351 of file JAAnetToolkit.hh.

351 { 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 359 of file JAAnetToolkit.hh.

359 { 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 367 of file JAAnetToolkit.hh.

367  { return (is_neutrino(track) ||
368  is_electron(track) ||
369  is_muon (track) ||
370  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 378 of file JAAnetToolkit.hh.

378  { return (is_electron(track) ||
379  is_muon (track) ||
380  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 388 of file JAAnetToolkit.hh.

388  { return (abs(track.type) != TRACK_TYPE_PHOTON &&
389  !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 397 of file JAAnetToolkit.hh.

397  { return (is_hadron(track) &&
398  ((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 406 of file JAAnetToolkit.hh.

406  { return (is_hadron(track) &&
407  ((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 416 of file JAAnetToolkit.hh.

416 { 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 424 of file JAAnetToolkit.hh.

425  {
426  if (Evt::ROOT_IO_VERSION >= 14) {
427 
428  return (track.status == TRK_ST_PRIMARYNEUTRINO ||
429  track.status == TRK_ST_PRIMARYCOSMIC ||
430  track.status == TRK_ST_MUONBUNDLE ||
431  track.status == TRK_ST_ININUCLEI);
432 
433  } else if (Evt::ROOT_IO_VERSION >= 11) {
434 
435  return (track.mother_id == TRK_MOTHER_UNDEFINED ||
436  track.mother_id == TRK_MOTHER_NONE);
437 
438  } else {
439 
440  return is_neutrino(track) && track.id == 0;
441  }
442  }
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.1.0-13-g917945e 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 450 of file JAAnetToolkit.hh.

451  {
452  if (Evt::ROOT_IO_VERSION >= 16) {
453  return track.status == TRK_ST_FINALSTATE;
454  } else {
455  return !is_initialstate(track);
456  }
457  }
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 465 of file JAAnetToolkit.hh.

466  {
467  if (Evt::ROOT_IO_VERSION >= 14) {
468 
469  std::vector<Trk>::const_iterator i = find_if(evt.mc_trks.cbegin(), evt.mc_trks.cend(),
471  return i != evt.mc_trks.cend();
472 
473  } else {
474 
475  return !evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0]);
476  }
477  }
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 485 of file JAAnetToolkit.hh.

486  {
487  if (Evt::ROOT_IO_VERSION >= 14) {
488 
489  std::vector<Trk>::const_iterator i = find_if(evt.mc_trks.cbegin(), evt.mc_trks.cend(),
491 
492  if (i != evt.mc_trks.cend()) { return *i; }
493 
494  } else if (!evt.mc_trks.empty() && is_neutrino(evt.mc_trks[0])) {
495 
496  return evt.mc_trks[0];
497  }
498 
499  THROW(JIndexOutOfRange, "JAANET::get_neutrino(): Event " << evt.id << " has no neutrino.");
500  }
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:712
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
JVertex3D JAANET::getVertex ( const Trk track)
inline

Get vertex.

Parameters
tracktrack
Returns
vertex

Definition at line 508 of file JAAnetToolkit.hh.

509  {
510  return JVertex3D(getPosition(track), track.t);
511  }
double t
track time [ns] (when the particle is at pos )
Definition: Trk.hh:19
JPosition3D getPosition(const Vec &pos)
Get position.
JVertex3D JAANET::getVertex ( const Evt event)
inline

Get event vertex.

Parameters
eventevent
Returns
event vertex

Definition at line 519 of file JAAnetToolkit.hh.

520  {
521  using namespace std;
522  using namespace JPP;
523 
524  if (has_neutrino(event)) {
525 
526  const Trk& neutrino = get_neutrino(event);
527 
528  return getVertex(neutrino);
529 
530  } else if (!event.mc_trks.empty()) {
531 
532  vector<Trk>::const_iterator i = find_if(event.mc_trks.begin(), event.mc_trks.end(), &is_initialstate);
533 
534  if (i != event.mc_trks.end()) {
535  return getVertex(*i);
536  } else {
537  THROW(JValueOutOfRange, "getVertex(): No initial state track found.");
538  }
539 
540  } else if (!event.trks.empty()) { // For reconstructed DAQ events
541 
542  const Trk& track = get_best_reconstructed_track<JPP_RECONSTRUCTION_TYPE>(event);
543 
544  return getVertex(track);
545 
546  } else {
547 
548  THROW(JValueOutOfRange, "getVertex(): Could not find a valid event vertex.");
549  }
550  }
JVertex3D getVertex(const Trk &track)
Get vertex.
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:712
bool is_initialstate(const Trk &track)
Test whether given track corresponds to an initial state particle.
std::vector< Trk > trks
list of reconstructed tracks (can be several because of prefits,showers, etc).
Definition: Evt.hh:39
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
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 558 of file JAAnetToolkit.hh.

559  {
560  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
561  }
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 569 of file JAAnetToolkit.hh.

570  {
571  return count_electrons(evt) != 0;
572  }
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 580 of file JAAnetToolkit.hh.

581  {
582  if (count_electrons(evt) > 0)
583  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_electron);
584  else
585  THROW(JIndexOutOfRange, "This event has no electron.");
586  }
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:712
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 597 of file JAAnetToolkit.hh.

598  {
599  return hit.type == GEANT4_TYPE_ELECTRON || hit.type == GEANT4_TYPE_ANTIELECTRON;
600  }
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 608 of file JAAnetToolkit.hh.

609  {
610  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
611  }
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 619 of file JAAnetToolkit.hh.

620  {
621  return count_muons(evt) != 0;
622  }
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 630 of file JAAnetToolkit.hh.

631  {
632  if (count_muons(evt) > 0)
633  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_muon);
634  else
635  THROW(JIndexOutOfRange, "This event has no muon.");
636  }
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:712
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 647 of file JAAnetToolkit.hh.

648  {
649  return hit.type == GEANT4_TYPE_MUON || hit.type == GEANT4_TYPE_ANTIMUON;
650  }
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 658 of file JAAnetToolkit.hh.

659  {
660  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
661  }
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 669 of file JAAnetToolkit.hh.

670  {
671  return count_taus(evt) != 0;
672  }
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 680 of file JAAnetToolkit.hh.

681  {
682  if (count_taus(evt) > 0)
683  return *std::find_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_tau);
684  else
685  THROW(JIndexOutOfRange, "This event has no tau.");
686  }
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:712
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 697 of file JAAnetToolkit.hh.

698  {
699  return hit.type == GEANT4_TYPE_TAU || hit.type == GEANT4_TYPE_ANTITAU;
700  }
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 708 of file JAAnetToolkit.hh.

709  {
710  return std::count_if(evt.mc_trks.begin(), evt.mc_trks.end(), is_hadron);
711  }
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 722 of file JAAnetToolkit.hh.

723  {
724  return (!from_electron(hit) && !from_muon(hit) && !from_tau(hit));
725  }
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 736 of file JAAnetToolkit.hh.

737  {
738  for (std::vector<Trk>::iterator p = evt.mc_trks.begin(); p != evt.mc_trks.end(); p++) {
739  p->t += tOff;
740  }
741  }
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 751 of file JAAnetToolkit.hh.

752  {
753  if (E > m) {
754  return sqrt((E - m) * (E + m));
755  } else {
756  return 0.0;
757  }
758  }
then usage $script< input file >[option[primary[working directory]]] nWhere option can be E
Definition: JMuonPostfit.sh:40
double JAANET::getKineticEnergy ( const Trk trk)
inline

Get track kinetic energy.

Parameters
trktrack
Returns
kinetic energy [GeV]

Definition at line 767 of file JAAnetToolkit.hh.

768  {
769  const double energy = trk.E;
770  const double mass = JPDB::getInstance().getPDG(trk.type).mass;
771 
772  return getKineticEnergy(energy, mass);
773  }
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 783 of file JAAnetToolkit.hh.

784  {
785  using namespace std;
786 
787  const Trk& neutrino = get_neutrino(evt);
788 
789  double E0 = neutrino.E;
790 
791  for (vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.end(); ++track) {
792 
793  if (!is_finalstate(*track)) { continue; }
794 
795  if (track->type == TRACK_TYPE_NEUTRON ||
796  track->type == TRACK_TYPE_SIGMA_MINUS ||
797  track->type == TRACK_TYPE_NEUTRAL_SIGMA) {
798 
799  E0 += MASS_NEUTRON;
800 
801  } else if (track->type == TRACK_TYPE_PROTON ||
802  track->type == TRACK_TYPE_LAMBDA ||
803  track->type == TRACK_TYPE_SIGMA_PLUS) {
804 
805  E0 += MASS_PROTON;
806 
807  } else if (track->type == TRACK_TYPE_ANTINEUTRON) {
808 
809  E0 -= MASS_NEUTRON;
810 
811  } else if (track->type == TRACK_TYPE_ANTIPROTON ||
812  track->type == TRACK_TYPE_ANTILAMBDA) {
813 
814  E0 -= MASS_PROTON;
815  }
816  }
817 
818  return E0;
819  }
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 829 of file JAAnetToolkit.hh.

830  {
831  using namespace std;
832  using namespace JPP;
833 
834  double E1 = 0.0;
835 
836  const Trk& neutrino = get_neutrino(evt);
837 
838  for (vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.cend(); ++track) {
839 
840  if (!is_finalstate(*track)) { continue; }
841 
842  if (is_muon(*track)) {
843 
844  const Trk& mother = ( track->mother_id > TRK_MOTHER_UNDEFINED ?
845  evt.mc_trks[track->mother_id] : neutrino );
846 
847  const double distance = (track->pos - mother.pos).len();
848 
849  E1 += gWater(track->E, -distance);
850 
851  } else {
852 
853  E1 += track->E;
854  }
855  }
856 
857  return E1;
858  }
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.1.0-13-g917945e 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 868 of file JAAnetToolkit.hh.

869  {
870  const Trk& neutrino = get_neutrino(evt);
871 
872  return neutrino.E * neutrino.dir;
873  }
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 883 of file JAAnetToolkit.hh.

884  {
885  using namespace std;
886  using namespace JPP;
887 
888  Vec P1(0,0,0);
889 
890  const Trk& neutrino = get_neutrino(evt);
891 
892  for (vector<Trk>::const_iterator track = evt.mc_trks.cbegin(); track != evt.mc_trks.cend(); ++track) {
893 
894  if (!is_finalstate(*track)) { continue; }
895 
896  double kineticEnergy = 0.0;
897 
898  if (is_muon(*track)) {
899 
900  const Trk& mother = ( track->mother_id > TRK_MOTHER_UNDEFINED ?
901  evt.mc_trks[track->mother_id] : neutrino );
902 
903  const double distance = (track->pos - mother.pos).len();
904  const double energy = gWater(track->E, -distance);
905 
906  kineticEnergy = getKineticEnergy(energy, MASS_MUON);
907 
908  } else {
909 
910  kineticEnergy = getKineticEnergy(*track);
911  }
912 
913  P1 += kineticEnergy * track->dir;
914  }
915 
916  return P1;
917  }
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.1.0-13-g917945e 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
double JAANET::getInvariantMass ( const Evt event)
inline

Get final state invariant mass.

Parameters
eventevent
Returns
invariant mass [GeV]

Definition at line 926 of file JAAnetToolkit.hh.

927  {
928  using namespace std;
929  using namespace JPP;
930 
931  double Minv = 0.0;
932 
933  for (vector<Trk>::const_iterator track = event.mc_trks.cbegin(); track != event.mc_trks.cend(); ++track) {
934 
935  if (!is_finalstate(*track)) { continue; }
936 
937  const double mass = JPDB::getInstance().getPDG(track->type).mass;
938 
939  Minv += mass;
940  }
941 
942  return Minv;
943  }
bool is_finalstate(const Trk &track)
Test whether given track corresponds to a final state particle.
T & getInstance(const T &object)
Get static instance from temporary object.
Definition: JObject.hh:75
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
Definition: Evt.hh:49
template<class JOscProbInterpolator_t = JOscProbInterpolator<>>
JEvtWeightFactorFunction<JAtmosphericNeutrinoFlux, JFlux> JAANET::make_atmosphericNeutrinoFluxFunction ( const std::string oscProbTableFile,
const JOscParameters &  oscParameters 
)
inline

Auxiliary method for creating an interface to an atmospheric neutrino flux function
using an oscillation probability interpolation table.

Parameters
oscProbTableFileoscillation probability interpolation table
oscParametersoscillation parameters

Definition at line 88 of file JAtmosphericNeutrinoFlux.hh.

89  {
90 
91  const JOscProbInterpolator_t interpolator(oscProbTableFile.c_str(), oscParameters);
92 
93  const JAtmosphericNeutrinoFlux flux(interpolator);
94 
95  return make_fluxFunction(flux);
96  }
JEvtWeightFactorFunction< JFunction_t, JFlux > make_fluxFunction(const JFunction_t &flux)
Auxiliary method for creating an interface to a flux function.
Implementation of atmospheric neutrino flux using official KM3NeT atmospheric flux function...
Neutrino flux.
Definition: JHead.hh:906
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 202 of file JEvtWeightFactorFunction.hh.

202  {
204  }
Implementation of event-weight factor interface.
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 214 of file JEvtWeightFactorFunction.hh.

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 226 of file JEvtWeightFactorFunction.hh.

226  {
228  }
Implementation of event-weight factor interface.
Neutrino flux.
Definition: JHead.hh:906
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 237 of file JEvtWeightFactorFunction.hh.

237  {
239  }
Implementation of event-weight factor interface.
Neutrino flux.
Definition: JHead.hh:906
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 249 of file JEvtWeightFactorFunction.hh.

249  {
251  }
Implementation of event-weight factor interface.
Neutrino flux.
Definition: JHead.hh:906
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 260 of file JEvtWeightFactorFunction.hh.

260  {
262  }
Implementation of C-style diffuse flux event-weight factor.
Neutrino flux.
Definition: JHead.hh:906
JEvtWeightFactorFunction<JFlatFlux, JFlux> JAANET::make_flatFluxFunction ( const double  flux)
inline

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

Parameters
fluxflux [GeV * m^-2 * sr^-1 * s^-1]

Definition at line 92 of file JFlatFlux.hh.

92  {
93 
94  const JFlatFlux flatFlux(flux);
95 
96  return make_fluxFunction(flatFlux);
97  }
JEvtWeightFactorFunction< JFunction_t, JFlux > make_fluxFunction(const JFunction_t &flux)
Auxiliary method for creating an interface to a flux function.
Neutrino flux.
Definition: JHead.hh:906
Function object for constant fluxes.
Definition: JFlatFlux.hh:20
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:88
T & getInstance(const T &object)
Get static instance from temporary object.
Definition: JObject.hh:75
static JNullStream null
Null I/O stream.
Definition: JNullStream.hh:51
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 1798 of file JHead.hh.

1800  {
1801  return JHead(first).match(JHead(second));
1802  }
Monte Carlo run header.
Definition: JHead.hh:1234
bool match(const JHead &header) const
Test match of headers.
Definition: JHead.hh:1472
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 1814 of file JHead.hh.

1816  {
1817  return JHead(first).less(JHead(second));
1818  }
bool less(const JHead &header) const
Comparison of headers.
Definition: JHead.hh:1484
Monte Carlo run header.
Definition: JHead.hh:1234
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 37 of file JHeadToolkit.hh.

38  {
39  for (const auto& i : header.simul) {
40  if (i.program == APPLICATION_GENHEN) {
41  return true;
42  }
43  }
44 
45  for (const auto& i : header.physics) { // legacy
46  if (i.buffer.find(APPLICATION_GENHEN) != std::string::npos) {
47  return true;
48  }
49  }
50 
51  return false;
52  }
static const char *const APPLICATION_GENHEN
KM3NeT Data Definitions v3.1.0-13-g917945e 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 61 of file JHeadToolkit.hh.

62  {
63  for (const auto& i : header.simul) {
64  if (i.program == APPLICATION_GSEAGEN) {
65  return true;
66  }
67  }
68 
69  for (const auto& i : header.physics) { // legacy
70  if (i.buffer.find(APPLICATION_GSEAGEN) != std::string::npos) {
71  return true;
72  }
73  }
74 
75  return false;
76  }
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 85 of file JHeadToolkit.hh.

86  {
87  for (const auto& i : header.simul) {
88  if (i.program == APPLICATION_MUPAGE) {
89  return true;
90  }
91  }
92 
93  return false;
94  }
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 103 of file JHeadToolkit.hh.

104  {
105  for (const auto& i : header.simul) {
106  if (i.program == APPLICATION_CORSIKA) {
107  return true;
108  }
109  }
110 
111  return false;
112  }
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 121 of file JHeadToolkit.hh.

122  {
123  for (const auto& i : header.simul) {
124  if (i.program == APPLICATION_KM3BUU) {
125  return true;
126  }
127  }
128 
129  return false;
130  }
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 139 of file JHeadToolkit.hh.

140  {
141  for (const auto& i : header.simul) {
142  if (i.program == APPLICATION_KM3) {
143  return true;
144  }
145  }
146 
147  return false;
148  }
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 157 of file JHeadToolkit.hh.

158  {
159  for (const auto& i : header.simul) {
160  if (i.program == APPLICATION_KM3SIM) {
161  return true;
162  }
163  }
164 
165  return false;
166  }
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 175 of file JHeadToolkit.hh.

176  {
177  for (const auto& i : header.simul) {
178  if (i.program == APPLICATION_JSIRENE) {
179  return true;
180  }
181  }
182 
183  return false;
184  }
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 193 of file JHeadToolkit.hh.

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

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

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

291  {
292  const Vec pos = get<Vec>(header);
293 
294  return JPosition3D(pos.x, pos.y, pos.z);
295  }
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 305 of file JHeadToolkit.hh.

306  {
307  using namespace JGEOMETRY2D;
308 
309  if (header.is_valid(&JHead::fixedcan)) {
310 
311  return JCylinder3D(JCircle2D(JVector2D(header.fixedcan.xcenter,
312  header.fixedcan.ycenter),
313  header.fixedcan.radius),
314  header.fixedcan.zmin,
315  header.fixedcan.zmax);
316 
317  } else if (header.is_valid(&JHead::can)) {
318 
319  return JCylinder3D(JCircle2D(JVector2D(), header.can.r),
320  header.can.zmin,
321  header.can.zmax);
322  } else {
323 
324  return JCylinder3D();
325  }
326  }
Data structure for vector in two dimensions.
Definition: JVector2D.hh:32
Data structure for circle in two dimensions.
Definition: JCircle2D.hh:33
JLorentzBoost JAANET::getBoostToCOM ( const Evt event)

Get Lorentz boost to the Center of Mass (COM) frame for a given neutrino interaction.

Parameters
eventevent

Definition at line 247 of file JLorentzBoost.hh.

248  {
249  using namespace JPP;
250 
251  if (has_neutrino(event)) {
252 
253  const Trk& nu = get_neutrino(event);
254 
255  const Vec beta = -nu.E / getE0(event) * nu.dir;
256 
257  return Boost(beta.x, beta.y, beta.z);
258 
259  } else {
260  THROW(JValueOutOfRange, "getBoostToCOM(): Given event does not correspond to a neutrino interaction.");
261  }
262  }
double z
Definition: Vec.hh:14
double getE0(const Evt &evt)
Get initial state energy of a neutrino interaction.
bool has_neutrino(const Evt &evt)
Test whether given event has an incoming neutrino.
Vec dir
track direction
Definition: Trk.hh:18
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:712
double E
Energy [GeV] (either MC truth or reconstructed)
Definition: Trk.hh:20
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
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
void JAANET::boostToCOM ( Evt event)

Boost event to the Center of Mass (COM) frame.

Parameters
eventevent

Definition at line 270 of file JLorentzBoost.hh.

271  {
272  const JLorentzBoost boost = getBoostToCOM(event);
273  boost(event);
274  }
JLorentzBoost getBoostToCOM(const Evt &event)
Get Lorentz boost to the Center of Mass (COM) frame for a given neutrino interaction.
Auxiliary class for performing Lorentz boosts.
JEvtWeightFactorFunction<JPowerLawFlux, JFlux> JAANET::make_powerLawFluxFunction ( const double  normalisation,
const double  spectralIndex 
)
inline

Auxiliary method for creating an interface to a power-law flux function.

Parameters
normalisationnormalisation [GeV * m^-2 * sr^-1 * s^-1]
spectralIndexspectral index

Definition at line 112 of file JPowerLawFlux.hh.

113  {
114 
115  const JPowerLawFlux flux(normalisation, spectralIndex);
116 
117  return make_fluxFunction(flux);
118  }
JEvtWeightFactorFunction< JFunction_t, JFlux > make_fluxFunction(const JFunction_t &flux)
Auxiliary method for creating an interface to a flux function.
Neutrino flux.
Definition: JHead.hh:906
Example function object for computing power-law flux.
bool JAANET::select ( const Trk trk,
const Evt evt 
)

Event selection.

Parameters
trktrack
evtevent
Returns
true if selected; else false

Definition at line 16 of file JAAnet/event_selector.cc.

17  {
18  using namespace std;
19  using namespace JPP;
20  /*
21  if (has_muon(evt)) {
22 
23  Trk muon;
24 
25  for (const auto& t1 : evt.mc_trks) {
26  if (is_muon(t1)) {
27  if (t1.E > muon.E) {
28  muon = t1;
29  }
30  }
31  }
32 
33  return (is_muon(muon) && getAngle(getDirection(trk), getDirection(muon)) > 5.0);
34  }
35  */
36  return true;
37  }
bool JAANET::hasW ( const Trk trk,
const int  i 
)
inline

Check availability of value.

Parameters
trktrack
iindex
Returns
true if available; else false

Definition at line 100 of file JAAnet/JEvD.cc.

101  {
102  return (i >= 0 && i < (int) trk.fitinf.size());
103  }
std::vector< double > fitinf
place to store additional fit info, see km3net-dataformat/definitions/fitparameters.csv
Definition: Trk.hh:32
double JAANET::getW ( const Trk trk,
const int  i 
)
inline

Get associated value.

Parameters
trktrack
iindex
Returns
value

Definition at line 113 of file JAAnet/JEvD.cc.

114  {
115  return trk.fitinf.at(i);
116  }
std::vector< double > fitinf
place to store additional fit info, see km3net-dataformat/definitions/fitparameters.csv
Definition: Trk.hh:32
double JAANET::getW ( const Trk trk,
const int  i,
const double  value 
)
inline

Get associated value.

Parameters
trktrack
iindex
valuedefault value
Returns
value

Definition at line 127 of file JAAnet/JEvD.cc.

128  {
129  if (hasW(trk,i))
130  return trk.fitinf.at(i);
131  else
132  return value;
133  }
std::vector< double > fitinf
place to store additional fit info, see km3net-dataformat/definitions/fitparameters.csv
Definition: Trk.hh:32
bool hasW(const Trk &trk, const int i)
Check availability of value.
Definition: JAAnet/JEvD.cc:100
void JAANET::setW ( Trk trk,
const int  i,
const double  value 
)

Set associated value.

Parameters
trktrack
iindex
valuevalue

Definition at line 143 of file JAAnet/JEvD.cc.

144  {
145  if (i >= (int) trk.fitinf.size()) {
146  trk.fitinf.resize(i + 1, 0.0);
147  }
148 
149  trk.fitinf[i] = value;
150  }
std::vector< double > fitinf
place to store additional fit info, see km3net-dataformat/definitions/fitparameters.csv
Definition: Trk.hh:32

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 12 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:85
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:61
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:37
static const char *const APPLICATION_GENHEN
KM3NeT Data Definitions v3.1.0-13-g917945e 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 231 of file JHeadToolkit.hh.