1 #ifndef __JAANET__JAANETTOOLKIT__
2 #define __JAANET__JAANETTOOLKIT__
41 namespace JPP {
using namespace JAANET; }
135 time_range.include(
getTime(*hit));
156 if (!time_range.is_valid()) {
160 if (time_range.getLength() < T_ns.getLength()) {
161 time_range.setLength(T_ns.getLength());
296 inline double getW(
const Trk& track,
const size_t index,
const double value)
298 if (index < track.
fitinf.size())
299 return track.
fitinf[index];
409 ((
int)(track.
type / 1000)) % 10 == 0); }
418 ((
int)(track.
type / 1000)) % 10 != 0); }
663 return sqrt((E - m) * (E + m));
705 for (
size_t i = 1; i != evt.
mc_trks.size(); ++i) {
750 for (
size_t i = 1; i != evt.
mc_trks.size(); ++i) {
786 P0 = neutrino.
E * neutrino.
dir;
811 for (
size_t i = 1; i != evt.
mc_trks.size(); ++i) {
815 double kineticEnergy = 0.0;
829 P1 += kineticEnergy * track.
dir;
bool from_tau(const Hit &hit)
Test whether given hit was produced by a tau.
Vec getP1(const Evt &evt)
Get momentum vector of the final state of a neutrino interaction.
JVertex3D getVertex(const Trk &track)
Get vertex.
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.
JHitType_t
Enumeration of hit types based on km3 codes.
Scattered light from primary shower.
JPredicate< JResult_t T::*, JComparison::eq > make_predicate(JResult_t T::*member, const JResult_t value)
Helper method to create predicate for data member.
bool is_proton(const Trk &track)
Test whether given track is a (anti-)proton.
const Trk & get_muon(const Evt &evt)
Get first muon from the event tracklist.
Data structure for direction in three dimensions.
JTrack3E getTrack(const Trk &track)
Get track.
Direct light from delta-rays.
double t
track time [ns] (when the particle is at pos )
Scattered light from muon.
JTransformation3D getTransformation(const Trk &track)
Get transformation.
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
double getE0(const Evt &evt)
Get initial state energy of a neutrino interaction.
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 defin...
bool has_neutrino(const Evt &evt)
Test whether given event has an incoming neutrino.
bool is_neutrino() const
Check if this is a netrino.
bool is_muon(const Trk &track)
Test whether given track is a (anti-)muon.
bool has_electron(const Evt &evt)
Test whether given event has an electron.
JTOOLS::JRange< double > JTimeRange
Type definition for time range (unit [ns]).
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.
double getE1(const Evt &evt)
Get final state energy of a neutrino interaction.
static const double MASS_MUON
muon mass [GeV]
bool has_particleID(const Trk &track, const int type)
Test whether given track corresponds to given particle type.
bool is_charged_lepton(const Trk &track)
Test whether given track is a charged lepton.
Scattered light from Bremsstrahlung.
static const JGeaneWater gWater
Function object for energy loss of muon in sea water.
const JParticle & getPDG(const int pdg) const
Get particle corresponding to given PDG code.
Direct light from Bremsstrahlung.
double getTime(const Hit &hit)
Get true time of hit.
bool is_baryon(const Trk &track)
Test whether given track is a baryon (is hadron and third digit of PDG code is not zero) ...
double E
Energy [GeV] (either MC truth or reconstructed)
Direct light from primary shower.
static const int AASHOWER_RECONSTRUCTION_TYPE
AAnet shower fit reconstruction type.
static const JPDB & getInstance()
Get particle data book.
bool is_noise(const Hit &hit)
Verify hit origin.
JTimeRange getTimeRange(const Evt &event)
Get time range (i.e. time between earliest and latest hit) of Monte Carlo event.
static const double MASS_NEUTRON
neutron mass [GeV]
bool is_neutrino(const Trk &track)
Test whether given track is a neutrino.
bool is_lepton(const Trk &track)
Test whether given track is a lepton.
double getW(const Trk &track, const size_t index, const double value)
Get track information.
double a
hit amplitude (in p.e.)
The Vec class is a straightforward 3-d vector, which also works in pyroot.
double getKineticEnergy(const double E, const double m)
Get kinetic energy of particle with given mass.
JDirection3D getDirection(const Vec &dir)
Get direction.
const Trk & get_electron(const Evt &evt)
Get first electron from the event tracklist.
double getDY() const
Get y direction.
double getDX() const
Get x direction.
JAxis3D getAxis(const Trk &track)
Get axis.
JPosition3D getPosition(const Vec &pos)
Get position.
std::vector< double > fitinf
place to store additional fit info, for jgandalf, see JFitParameters.hh
static const int TRK_ST_PRIMARYNEUTRINO
initial state neutrino ('neutrino' tag in evt files from gseagen and genhen).
double getY() const
Get y position.
bool is_photon(const Trk &track)
Test whether given track is a photon or neutral pion.
bool from_hadron(const Hit &hit)
Test whether given hit was produced by a hadronic shower.
int type
MC: particle type in PDG encoding.
int status
MC status code, see km3net-dataformat/definitions/trkmembers.csv for values.
then for APP in event gandalf start energy
bool is_pion(const Trk &track)
Test whether given track is a charged pion.
Vec pos
postion of the track at time t [m]
double mass
mass of particle [GeV]
int count_taus(const Evt &evt)
Count the number of taus in a given event.
const Trk & get_tau(const Evt &evt)
Get first tau from the event tracklist.
bool has_muon(const Evt &evt)
Test whether given event has a muon.
std::vector< Hit > mc_hits
MC: list of MC truth hits.
double t
hit time (from tdc+calibration or MC truth)
Scattered light from delta-rays.
static const double MASS_PROTON
proton mass [GeV]
Definition of particle types.
double getX() const
Get x position.
bool has_tau(const Evt &evt)
Test whether given event has a tau.
int id
offline event identifier
double getNPE(const Hit &hit)
Get true charge of hit.
Data structure for position in three dimensions.
bool is_neutron(const Trk &track)
Test whether given track is a (anti-)neutron.
Exception for accessing a value in a collection that is outside of its range.
Exception for accessing an index in a collection that is outside of its range.
Vec getP0(const Evt &evt)
Get momentum vector of the initial state of a neutrino interaction.
const Trk & get_neutrino(const Evt &evt)
Get incoming neutrino.
bool is_hadron(const Trk &track)
Test whether given track is a hadron.
bool is_tau(const Trk &track)
Test whether given track is a (anti-)tau.
bool from_electron(const Hit &hit)
Test whether given hit was produced by an electron.
The Trk class represents a Monte Carlo (MC) particle as well as a reconstructed track/shower.
int type
particle type or parametrisation used for hit (mc only)
double getZ() const
Get z position.
double getDZ() const
Get z direction.
std::vector< Trk > mc_trks
MC: list of MC truth tracks.
bool is_meson(const Trk &track)
Test whether given track is a meson (is hadron and third digit of PDG code is zero) ...
The Evt class respresent a Monte Carlo (MC) event as well as an offline event.
int count_hadrons(const Evt &evt)
Count the number of hadrons in a given event (not including neutral pions)
bool from_muon(const Hit &hit)
Test whether given hit was produced by a muon.