G4Tree.h
Go to the documentation of this file.
1 /*
2  * G4Tree.h
3  *
4  * Created on: Feb 18, 2021
5  * Author: chilgenb
6  */
7 
8 #ifndef GARANA_G4TREE_H_
9 #define GARANA_G4TREE_H_
10 
11 #include "garana/Base/TreeReader.h"
12 
13 using std::vector;
14 
15 namespace garana {
16 
17  class G4Tree : public TreeReader {
18 
19  public:
20 
21  // default const'or
22  virtual ~G4Tree() {};
23 
24  // pure virtual functions for structured or flat g4Tree
25  virtual const UInt_t NSim() const = 0; ///< number of particles
26  virtual const UInt_t NPoints(const UInt_t& iparticle) const = 0; ///< number of G4 steps (i.e. trajectory points)
27  virtual const UInt_t NRegions(const UInt_t& iparticle) const = 0; ///< number of regions traversed by particle
28  virtual const Int_t Region(const UInt_t& iparticle, const UInt_t& iregion) const = 0; ///< region number
29  virtual const vector<const TLorentzVector*>* SimMomEnter(const UInt_t& iparticle)const = 0; ///< particle 4-momentum at entry point, all regions
30  virtual const vector<const TLorentzVector*>* SimMomExit(const UInt_t& iparticle) const = 0; ///< particle 4-momentum at exit point, all regions
31  virtual const vector<const TLorentzVector*>* SimPosEnter(const UInt_t& iparticle)const = 0; ///< particle 4-position at entry point, all regions
32  virtual const vector<const TLorentzVector*>* SimPosExit(const UInt_t& iparticle) const = 0; ///< particle 4-position at exit point, all regions
33  virtual const TLorentzVector* SimMomEnter(const UInt_t& iparticle, const UInt_t& iregion)const = 0; ///< particle 4-momentum at entry point in region
34  virtual const TLorentzVector* SimMomExit(const UInt_t& iparticle, const UInt_t& iregion) const = 0; ///< particle 4-momentum at exit point in region
35  virtual const TLorentzVector* SimPosEnter(const UInt_t& iparticle, const UInt_t& iregion)const = 0; ///< particle 4-position at entry point in region
36  virtual const TLorentzVector* SimPosExit(const UInt_t& iparticle, const UInt_t& iregion) const = 0; ///< particle 4-position at exit point in region
37  virtual const bool IsPrimary(const UInt_t& iparticle) const = 0; ///< did particle come from generator?
38  virtual const Int_t PDG(const UInt_t& iparticle) const = 0; ///< particle PDG code
39  virtual const int ParentPDG(const UInt_t& iparticle) const = 0; ///< parent particle's PDG code
40  virtual const int ProgenitorPDG(const UInt_t& iparticle) const = 0; ///< PDG of primary that led this one
41  virtual const int TrackID(const UInt_t& iparticle) const = 0; ///< G4 track ID (can be <0 if it fell below trking threshold)
42  virtual const int ParentTrackID(const UInt_t& iparticle) const = 0; ///< G4 track ID of parent particle
43  virtual const int ProgenitorTrackID(const UInt_t& iparticle) const = 0; ///< G4 track ID of primary that led this one
44  virtual const Int_t ProcessI(const UInt_t& iparticle) const = 0; ///< code for process that created this one
45  virtual const Int_t ProcessF(const UInt_t& iparticle) const = 0; ///< code for process that killed this one
46 
47  // truth matching
48  UInt_t const GetTruthIndex(UInt_t iparticle) const; ///< index in gen tree subentry to truth match to this
49 
50  const UInt_t NPrimary() const;
51 
52  bool HasPassedTPC(const UInt_t& iparticle) const; ///< did the G4Particle pass through any TPC drift volume(s)?
53  bool HasPassedCalo(const UInt_t& iparticle) const; ///< did the G4Particle pass through any active ECal volume(s)?
54  bool IsStoppedTPC(const UInt_t& iparticle) const; ///< did the G4Particle stop/decay in any TPC drift volume(s)?
55  bool IsStoppedCalo(const UInt_t& iparticle) const; ///< did the G4Particle stop/decay in any active ECal volume(s)?
56  bool IsContainedTPC(const UInt_t& iparticle) const; ///< if the G4Particle was produced in any TPC drift volume, does it remain in either drift volume?
57  bool IsContainedCalo(const UInt_t& iparticle) const; ///< if the G4Particle was produced in any active ECal volume, does it remain there?
58  bool IsCathodeCrosser(const UInt_t& iparticle) const; ///< did the G4Particle cross the TPC central cathode?
59  bool IsContainedTPCEvent() const; ///< do all particles produced in any TPC drift volume in this event remain in either volume?
60  bool IsContainedTPCPrimaries() const; ///< do all primaries produced in any TPC drift volume in this event remain in either volume?
61  bool IsContainedCaloEvent() const; ///< do all particles produced in any active ECal volume in this event remain there?
62  bool IsContainedCaloPrimaries() const; ///< do all primaries produced in any active ECal volume in this event remain there?
63 
64  const TLorentzVector* SimMomBegin(const UInt_t& iparticle) const;
65  const TLorentzVector* SimMomEnd(const UInt_t& iparticle) const;
66  const TLorentzVector* SimPosBegin(const UInt_t& iparticle) const;
67  const TLorentzVector* SimPosEnd(const UInt_t& iparticle) const;
68 
69  protected:
70 
71  vector<UInt_t>* fG4TruthIndex = nullptr;
72  vector<UInt_t>* fG4FSIndex = nullptr;
73  TBranch* b_G4TruthIndex = nullptr;
74  TBranch* b_G4FSIndex = nullptr;
75  };//class
76 }//namespace
77 
78 #endif /* GARANA_G4TREE_H_ */
bool IsStoppedCalo(const UInt_t &iparticle) const
did the G4Particle stop/decay in any active ECal volume(s)?
Definition: G4Tree.cxx:66
virtual const int TrackID(const UInt_t &iparticle) const =0
G4 track ID (can be <0 if it fell below trking threshold)
virtual const int ProgenitorPDG(const UInt_t &iparticle) const =0
PDG of primary that led this one.
virtual const UInt_t NSim() const =0
number of particles
virtual const UInt_t NRegions(const UInt_t &iparticle) const =0
number of regions traversed by particle
virtual const Int_t PDG(const UInt_t &iparticle) const =0
particle PDG code
bool IsCathodeCrosser(const UInt_t &iparticle) const
did the G4Particle cross the TPC central cathode?
Definition: G4Tree.cxx:107
bool IsContainedCaloEvent() const
do all particles produced in any active ECal volume in this event remain there?
Definition: G4Tree.cxx:145
bool IsContainedTPCPrimaries() const
do all primaries produced in any TPC drift volume in this event remain in either volume?
Definition: G4Tree.cxx:133
const TLorentzVector * SimPosEnd(const UInt_t &iparticle) const
Definition: G4Tree.cxx:181
virtual const vector< const TLorentzVector * > * SimPosEnter(const UInt_t &iparticle) const =0
particle 4-position at entry point, all regions
struct vector vector
virtual const vector< const TLorentzVector * > * SimMomEnter(const UInt_t &iparticle) const =0
particle 4-momentum at entry point, all regions
virtual const int ParentPDG(const UInt_t &iparticle) const =0
parent particle&#39;s PDG code
bool IsStoppedTPC(const UInt_t &iparticle) const
did the G4Particle stop/decay in any TPC drift volume(s)?
Definition: G4Tree.cxx:58
const TLorentzVector * SimMomBegin(const UInt_t &iparticle) const
Definition: G4Tree.cxx:169
virtual const vector< const TLorentzVector * > * SimPosExit(const UInt_t &iparticle) const =0
particle 4-position at exit point, all regions
bool HasPassedCalo(const UInt_t &iparticle) const
did the G4Particle pass through any active ECal volume(s)?
Definition: G4Tree.cxx:47
const TLorentzVector * SimPosBegin(const UInt_t &iparticle) const
Definition: G4Tree.cxx:177
vector< UInt_t > * fG4TruthIndex
Definition: G4Tree.h:71
TBranch * b_G4TruthIndex
Definition: G4Tree.h:73
virtual const bool IsPrimary(const UInt_t &iparticle) const =0
did particle come from generator?
virtual const int ProgenitorTrackID(const UInt_t &iparticle) const =0
G4 track ID of primary that led this one.
bool IsContainedCaloPrimaries() const
do all primaries produced in any active ECal volume in this event remain there?
Definition: G4Tree.cxx:156
virtual const Int_t ProcessI(const UInt_t &iparticle) const =0
code for process that created this one
virtual ~G4Tree()
Definition: G4Tree.h:22
bool IsContainedCalo(const UInt_t &iparticle) const
if the G4Particle was produced in any active ECal volume, does it remain there?
Definition: G4Tree.cxx:90
bool HasPassedTPC(const UInt_t &iparticle) const
did the G4Particle pass through any TPC drift volume(s)?
Definition: G4Tree.cxx:36
vector< UInt_t > * fG4FSIndex
Definition: G4Tree.h:72
bool IsContainedTPCEvent() const
do all particles produced in any TPC drift volume in this event remain in either volume?
Definition: G4Tree.cxx:122
const TLorentzVector * SimMomEnd(const UInt_t &iparticle) const
Definition: G4Tree.cxx:173
virtual const Int_t ProcessF(const UInt_t &iparticle) const =0
code for process that killed this one
UInt_t const GetTruthIndex(UInt_t iparticle) const
index in gen tree subentry to truth match to this
Definition: G4Tree.cxx:12
virtual const vector< const TLorentzVector * > * SimMomExit(const UInt_t &iparticle) const =0
particle 4-momentum at exit point, all regions
virtual const int ParentTrackID(const UInt_t &iparticle) const =0
G4 track ID of parent particle.
bool IsContainedTPC(const UInt_t &iparticle) const
if the G4Particle was produced in any TPC drift volume, does it remain in either drift volume...
Definition: G4Tree.cxx:74
virtual const UInt_t NPoints(const UInt_t &iparticle) const =0
number of G4 steps (i.e. trajectory points)
virtual const Int_t Region(const UInt_t &iparticle, const UInt_t &iregion) const =0
region number
const UInt_t NPrimary() const
Definition: G4Tree.cxx:27
TBranch * b_G4FSIndex
Definition: G4Tree.h:74