G4Tree.cxx
Go to the documentation of this file.
1 /*
2  * G4Tree.cxx
3  *
4  * Created on: Feb 18, 2021
5  * Author: chilgenb
6  */
7 
8 #include "garana/Base/G4Tree.h"
9 
10 using namespace garana;
11 
12 const UInt_t G4Tree::GetTruthIndex(UInt_t iparticle) const {
13  try {
14  if(iparticle < fG4TruthIndex->size())
15  return fG4TruthIndex->at(iparticle);
16 
17  throw(iparticle);
18  }
19  catch(UInt_t index){
20  std::cerr << "G4TruthIndex out of range ("
21  << iparticle << " vs. " << fG4TruthIndex->size()
22  << std::endl;
23  return UINT_MAX;
24  }
25 }
26 
27 const UInt_t G4Tree::NPrimary() const {
28 
29  UInt_t nprimary=0;
30  while(nprimary<this->NSim() && this->IsPrimary(nprimary))
31  nprimary++;
32 
33  return nprimary;
34 }
35 
36 bool G4Tree::HasPassedTPC(const UInt_t& iparticle) const {
37 
38  for(size_t ireg=0; ireg<this->NRegions(iparticle); ireg++) {
39  if(this->Region(iparticle,ireg)==0 || this->Region(iparticle,ireg)==1) {
40  return true;
41  }
42  }
43 
44  return false;
45 }
46 
47 bool G4Tree::HasPassedCalo(const UInt_t& iparticle) const {
48 
49  for(size_t ireg=0; ireg<this->NRegions(iparticle); ireg++) {
50  if(this->Region(iparticle,ireg)==3 || this->Region(iparticle,ireg)==4) {
51  return true;
52  }
53  }
54 
55  return false;
56 }
57 
58 bool G4Tree::IsStoppedTPC(const UInt_t& iparticle) const {
59  UInt_t ireg = this->NRegions(iparticle)-1;
60  if(this->Region(iparticle,ireg)==0 || this->Region(iparticle,ireg)==1)
61  return true;
62  else
63  return false;
64 }
65 
66 bool G4Tree::IsStoppedCalo(const UInt_t& iparticle) const {
67  UInt_t ireg = this->NRegions(iparticle)-1;
68  if(this->Region(iparticle,ireg)==3 || this->Region(iparticle,ireg)==4)
69  return true;
70  else
71  return false;
72 }
73 
74 bool G4Tree::IsContainedTPC(const UInt_t& iparticle) const {
75 
76  bool isAV = false;
77  for(size_t ireg=0; ireg<this->NRegions(iparticle); ireg++) {
78  if(this->Region(iparticle,ireg)!=0 && this->Region(iparticle,ireg)!=1)
79  return false;
80  else
81  isAV = true;
82  }
83 
84  if(isAV)
85  return true;
86  else
87  return false;
88 }
89 
90 bool G4Tree::IsContainedCalo(const UInt_t& iparticle) const {
91 
92  bool isCalo = false;
93  for(size_t ireg=0; ireg<this->NRegions(iparticle); ireg++) {
94  if(this->Region(iparticle,ireg)!=3 && this->Region(iparticle,ireg)!=4)
95  return false;
96  else
97  isCalo = true;
98  }
99 
100  if(isCalo)
101  return true;
102  else
103  return false;
104 
105 }
106 
107 bool G4Tree::IsCathodeCrosser(const UInt_t& iparticle) const {
108 
109  if(this->NRegions(iparticle)<2)
110  return false;
111 
112  for(size_t ireg=0; ireg<this->NRegions(iparticle)-1; ireg++) {
113 
114  if( (this->Region(iparticle,ireg)==0 && this->Region(iparticle,ireg+1)==1) ||
115  (this->Region(iparticle,ireg)==1 && this->Region(iparticle,ireg+1)==0) )
116  return true;
117  }
118 
119  return false;
120 }
121 
123 
124  for(size_t isim=0; isim<this->NSim(); isim++){
125  if(this->HasPassedTPC(isim))
126  if(!this->IsContainedTPC(isim))
127  return false;
128  }
129 
130  return true;
131 }
132 
134 
135  for(size_t isim=0; isim<this->NSim(); isim++){
136  if(this->HasPassedTPC(isim))
137  if(this->IsPrimary(isim))
138  if(!this->IsContainedTPC(isim))
139  return false;
140  }
141 
142  return true;
143 }
144 
146 
147  for(size_t isim=0; isim<this->NSim(); isim++){
148  if(this->HasPassedCalo(isim))
149  if(!this->IsContainedCalo(isim))
150  return false;
151  }
152 
153  return true;
154 }
155 
157 
158  for(size_t isim=0; isim<this->NSim(); isim++){
159  if(this->HasPassedCalo(isim))
160  if(this->IsPrimary(isim))
161  if(!this->IsContainedCalo(isim))
162  return false;
163  }
164 
165  return true;
166 
167 }
168 
169 const TLorentzVector* G4Tree::SimMomBegin(const UInt_t& iparticle) const {
170  return this->SimMomEnter(iparticle,0);
171 }
172 
173 const TLorentzVector* G4Tree::SimMomEnd(const UInt_t& iparticle) const {
174  return this->SimMomExit(iparticle,this->NRegions(iparticle)-1);
175 }
176 
177 const TLorentzVector* G4Tree::SimPosBegin(const UInt_t& iparticle) const {
178  return this->SimPosEnter(iparticle,0);
179 }
180 
181 const TLorentzVector* G4Tree::SimPosEnd(const UInt_t& iparticle) const {
182  return this->SimPosExit(iparticle,this->NRegions(iparticle)-1);
183 }
bool IsStoppedCalo(const UInt_t &iparticle) const
did the G4Particle stop/decay in any active ECal volume(s)?
Definition: G4Tree.cxx:66
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
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
virtual const vector< const TLorentzVector * > * SimMomEnter(const UInt_t &iparticle) const =0
particle 4-momentum at entry point, all regions
bool IsStoppedTPC(const UInt_t &iparticle) const
did the G4Particle stop/decay in any TPC drift volume(s)?
Definition: G4Tree.cxx:58
decltype(auto) constexpr size(T &&obj)
ADL-aware version of std::size.
Definition: StdUtils.h:92
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
virtual const bool IsPrimary(const UInt_t &iparticle) const =0
did particle come from generator?
bool IsContainedCaloPrimaries() const
do all primaries produced in any active ECal volume in this event remain there?
Definition: G4Tree.cxx:156
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
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
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
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 Int_t Region(const UInt_t &iparticle, const UInt_t &iregion) const =0
region number
const UInt_t NPrimary() const
Definition: G4Tree.cxx:27
QTextStream & endl(QTextStream &s)