Hall A ROOT/C++ Analyzer (podd)
Loading...
Searching...
No Matches
THaTrack.h
Go to the documentation of this file.
1#ifndef Podd_THaTrack_h_
2#define Podd_THaTrack_h_
3
5//
6// THaTrack
7//
9
10#include "TObject.h"
11#include "TVector3.h"
12#include "THaPIDinfo.h"
13#include "DataType.h" // for kBig
14#include <cstring> // for memset
15
16class TClonesArray;
18class THaCluster;
19class THaTrackID;
20
21class THaTrack : public TObject {
22
23public:
24
25 // Bits for fType
26 enum {
27 kHasDet = BIT(0), // Detector coordinates set
28 kHasFP = BIT(1), // Focal plane coordinates set
29 kHasRot = BIT(2), // Rotating TRANSPORT coordinates set
30 kHasTarget = BIT(3), // Target coordinates reconstructed
31 kHasVertex = BIT(4) // Vertex reconstructed
32 };
33
34 // Default constructor
50
51 // Constructor with fp coordinates
52 // FIXME: this really should be setting detector coordinates
54 THaTrackingDetector* creator=nullptr, THaTrackID* id=nullptr,
55 THaPIDinfo* pid=nullptr )
56 : TObject(),
57 fX(x), fY(y), fTheta(theta), fPhi(phi), fP(kBig),
64 fNclusters(0), fClusters{}, fPIDinfo(pid), fCreator(creator), fIndex(-1),
65 fTrkNum(0), fID(id), fFlag(0), fType(kHasFP), fChi2(kBig), fNDoF(0),
69 {
70 if(pid) pid->Clear();
71 }
72
73 virtual ~THaTrack();
74
76 void Clear( Option_t* opt="" );
78 Int_t GetNclusters() const { return fNclusters; }
79 Int_t GetIndex() const { return fIndex; }
80 THaCluster* GetCluster( Int_t i ) { return fClusters[i]; }
81 UInt_t GetFlag() const { return fFlag; }
82 UInt_t GetType() const { return fType; }
83 THaTrackID* GetID() const { return fID; }
84 Int_t GetTrkNum() const { return fTrkNum; }
85
86 Double_t GetP() const { return fP; }
87 Double_t GetPhi() const { return fPhi; }
88 THaPIDinfo* GetPIDinfo() const { return fPIDinfo; }
89 Double_t GetTheta() const { return fTheta; }
90 Double_t GetX() const { return fX; }
91 Double_t GetY() const { return fY; }
92 Double_t GetX( Double_t z ) const { return fX + z*fTheta; }
93 Double_t GetY( Double_t z ) const { return fY + z*fPhi; }
94
95 Double_t GetChi2() const { return fChi2; }
96 Int_t GetNDoF() const { return fNDoF; }
97
98 Double_t GetDX() const { return fDX; }
99 Double_t GetDY() const { return fDY; }
100 Double_t GetDTheta() const { return fDTheta; }
101 Double_t GetDPhi() const { return fDPhi; }
102 Double_t GetRX() const { return fRX; }
103 Double_t GetRY() const { return fRY; }
104 Double_t GetRTheta() const { return fRTheta; }
105 Double_t GetRPhi() const { return fRPhi; }
106 Double_t GetTX() const { return fTX; }
107 Double_t GetTY() const { return fTY; }
108 Double_t GetTTheta() const { return fTTheta; }
109 Double_t GetTPhi() const { return fTPhi; }
110 Double_t GetDp() const { return fDp; }
111 Double_t GetLabPx() const { return fPvect.X(); }
112 Double_t GetLabPy() const { return fPvect.Y(); }
113 Double_t GetLabPz() const { return fPvect.Z(); }
114 Double_t GetVertexX() const { return fVertex.X(); }
115 Double_t GetVertexY() const { return fVertex.Y(); }
116 Double_t GetVertexZ() const { return fVertex.Z(); }
117 Double_t GetPathLen() const { return fPathl; }
118
119 TVector3& GetPvect() { return fPvect; }
122
123 Double_t GetTime() const { return fTime; } // at refplane (s)
124 Double_t GetdTime() const { return fdTime; } // (s)
125 Double_t GetBeta() const { return fBeta; } // from scint.
126 Double_t GetdBeta() const { return fdBeta; }
127 Double_t GetDedx() const { return fDedx; }
128 Double_t GetEnergy() const { return fEnergy; }
129 Int_t GetNPMT() const { return fNPMT; }
130 Double_t GetBetaChi2() const { return fBetaChi2; }
131 Double_t GetFPTime() const { return fFPTime; }
132 Int_t GetGoodPlane3() const { return fGoodPlane3; }
133 Int_t GetGoodPlane4() const { return fGoodPlane4; }
134
135 bool HasDet() const { return (fType&kHasDet); }
136 bool HasFP() const { return (fType&kHasFP); }
137 bool HasRot() const { return (fType&kHasRot); }
138 bool HasTarget() const { return (fType&kHasTarget); }
139 bool HasVertex() const { return (fType&kHasVertex); }
140
141 void Print( Option_t* opt="" ) const;
142
143 void Set( Double_t x, Double_t y, Double_t theta, Double_t phi )
144 { fX = x; fY = y; fTheta = theta; fPhi = phi; fType |= kHasFP; }
145 void SetR( Double_t x, Double_t y,
146 Double_t theta, Double_t phi );
147 void SetD( Double_t x, Double_t y,
148 Double_t theta, Double_t phi );
149 void SetTarget( Double_t x, Double_t y,
150 Double_t theta, Double_t phi );
151
152 void SetPathLen( Double_t pathl ) { fPathl = pathl; /* meters */ }
153 void SetTime( Double_t time ) { fTime = time; /* seconds */ }
154 void SetdTime( Double_t dt ) { fdTime = dt; /* seconds */ }
155 void SetBeta( Double_t beta ) { fBeta = beta; }
156 void SetdBeta( Double_t db ) { fdBeta = db; }
157 void SetDedx(Double_t dedx) { fDedx = dedx; }
158 void SetEnergy(Double_t energy) { fEnergy = energy; }
159 void SetNPMT(Int_t npmt) { fNPMT = npmt; }
160 void SetBetaChi2(Double_t betachi2) { fBetaChi2 = betachi2; }
161 void SetFPTime(Double_t fptime) { fFPTime = fptime; }
162 void SetGoodPlane3(Int_t gdplane3) { fGoodPlane3 = gdplane3; }
163 void SetGoodPlane4(Int_t gdplane4) { fGoodPlane4 = gdplane4; }
164
165 void SetChi2( Double_t chi2, Int_t ndof ) { fChi2=chi2; fNDoF=ndof; }
166
167 void SetID( THaTrackID* id ) { fID = id; }
168 void SetFlag( UInt_t flag ) { fFlag = flag; }
169 void SetType( UInt_t flag ) { fType = flag; }
170 void SetMomentum( Double_t p ) { fP = p; }
171 void SetDp( Double_t dp ) { fDp = dp; }
172 void SetTrkNum( Int_t n ) { fTrkNum = n; }
174 void SetIndex( Int_t idx ) { fIndex = idx; }
175 void SetPIDinfo( THaPIDinfo* pid ) { fPIDinfo = pid; }
176 void SetPvect( const TVector3& pvect ) { fPvect = pvect; }
177 void SetVertex( const TVector3& vert )
178 { fVertex = vert; fType |= kHasVertex; }
180 { fVertex.SetXYZ( x, y, z ); fType |= kHasVertex; }
181 void SetVertexError( const TVector3& err )
182 { fVertexError = err; }
185
186 virtual Bool_t IsSortable() const { return true; }
187 virtual Int_t Compare(const TObject* obj) const;
188
189protected:
190
191 enum { kMAXCL = 4 };
192
193 // Focal plane coordinates (TRANSPORT system projected to z=0)
194 Double_t fX; // x position in TRANSPORT plane (m)
195 Double_t fY; // y position in TRANSPORT plane (m)
196 Double_t fTheta; // Tangent of TRANSPORT Theta (x')
197 Double_t fPhi; // Tangent of TRANSPORT Phi (y')
198 Double_t fP; // Track momentum (GeV)
199
200 // coordinates in the detector system
201 Double_t fDX; // x position in DCS
202 Double_t fDY; // y position in DCS
203 Double_t fDTheta; // Tangent of DCS Theta
204 Double_t fDPhi; // Tangent of DCS Phi
205
206 // coordinates in the rotated TRANSPORT system
207 Double_t fRX; // x position in focal plane (m)
208 Double_t fRY; // y position in focal plane (m)
209 Double_t fRTheta; // Tangent of TRANSPORT Theta (x')
210 Double_t fRPhi; // Tangent of TRANSPORT Phi (y')
211
212 // reconstructed coordinates at the target in the target TRANSPORT system (z=0)
213 Double_t fTX; // x position at target (m)
214 Double_t fTY; // y position at target (m)
215 Double_t fTTheta; // Tangent of TRANSPORT Theta (out-of-plane angle) at target
216 Double_t fTPhi; // Tangent of TRANSPORT Phi (in-plane angle) at target
217 Double_t fDp; // dp/p_center -- fractional change in momentum
218
219 TVector3 fPvect; // Momentum vector at target in lab system (GeV)
220 TVector3 fVertex; // Vertex location in lab (m) valid if fHasVertex
221 TVector3 fVertexError; // Uncertainties in fVertex coordinates.
222
223 Double_t fPathl; // pathlength from target (z=0) (meters)
224
225 Double_t fTime; // time of track at focal plane (s)
226 Double_t fdTime; // uncertainty in fTime
227 Double_t fBeta; // beta of track
228 Double_t fdBeta; // uncertainty in fBeta
229
230 // Status variables and objects related to this track
235 Int_t fIndex; // Track index (-1 = none, 0 = first, etc.)
236 Int_t fTrkNum; // Track number (0 = unassigned)
237
239 UInt_t fFlag; // General status flag (for use by tracking det.)
240 UInt_t fType; // Flag indicating which vectors reconstructed
241
242 Double_t fChi2; // goodness of track fit
243 Int_t fNDoF; // number of hits on the track contributing to chi2
244
245 Double_t fDedx; // dEdX from hodoscopes
246 Double_t fEnergy; // Energy from calorimeter
247 // Needed for "prune" select best track method
248 Int_t fNPMT; // Number of PMTs hit in track
249 Double_t fBetaChi2; // (reduced) Chisq of fit on Beta
250 Double_t fFPTime; // Focal Plane time (same as fTime?)
251 Int_t fGoodPlane3; // Track hit a plane 3 paddle
252 Int_t fGoodPlane4; // Track hit a plane 4 paddle
253
254 ClassDef(THaTrack,5) // A generic particle track
255};
256
257//__________________ inlines __________________________________________________
258inline
260{
261 // set the coordinates in the rotated focal-plane frame
262 fDX = x;
263 fDY = y;
264 fDTheta = theta;
265 fDPhi = phi;
266 fType |= kHasDet;
267}
268
269//_____________________________________________________________________________
270inline
272{
273 // set the coordinates in the rotated focal-plane frame
274 fRX = x;
275 fRY = y;
276 fRTheta = theta;
277 fRPhi = phi;
278 fType |= kHasRot;
279}
280
281//_____________________________________________________________________________
282inline
284{
285 // set the coordinates in the target frame
286 fTX = x;
287 fTY = y;
288 fTTheta = theta;
289 fTPhi = phi;
290 fType |= kHasTarget;
291}
292
293#endif
int Int_t
unsigned int UInt_t
const Data_t kBig
Definition DataType.h:15
uint32_t time
#define d(i)
bool Bool_t
double Double_t
const char Option_t
#define ClassDef(name, id)
#define BIT(n)
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize id
Double_t GetBetaChi2() const
Definition THaTrack.h:130
bool HasVertex() const
Definition THaTrack.h:139
Double_t fDedx
Definition THaTrack.h:245
Double_t fTime
Definition THaTrack.h:225
TVector3 & GetVertex()
Definition THaTrack.h:120
Double_t fPathl
Definition THaTrack.h:223
Double_t GetLabPx() const
Definition THaTrack.h:111
Int_t GetNPMT() const
Definition THaTrack.h:129
void SetPvect(const TVector3 &pvect)
Definition THaTrack.h:176
void SetID(THaTrackID *id)
Definition THaTrack.h:167
TVector3 fVertexError
Definition THaTrack.h:221
UInt_t fType
Definition THaTrack.h:240
Int_t fGoodPlane4
Definition THaTrack.h:252
Double_t GetDPhi() const
Definition THaTrack.h:101
THaTrackingDetector * GetCreator() const
Definition THaTrack.h:77
Double_t GetX(Double_t z) const
Definition THaTrack.h:92
Double_t fBetaChi2
Definition THaTrack.h:249
Double_t GetPathLen() const
Definition THaTrack.h:117
Int_t fGoodPlane3
Definition THaTrack.h:251
Double_t fdTime
Definition THaTrack.h:226
Int_t fNPMT
Definition THaTrack.h:248
Double_t GetChi2() const
Definition THaTrack.h:95
void SetFPTime(Double_t fptime)
Definition THaTrack.h:161
void Set(Double_t x, Double_t y, Double_t theta, Double_t phi)
Definition THaTrack.h:143
Double_t GetY(Double_t z) const
Definition THaTrack.h:93
Double_t fX
Definition THaTrack.h:194
TVector3 & GetVertexError()
Definition THaTrack.h:121
Double_t GetDX() const
Definition THaTrack.h:98
Double_t GetDTheta() const
Definition THaTrack.h:100
Double_t fRTheta
Definition THaTrack.h:209
bool HasFP() const
Definition THaTrack.h:136
virtual Bool_t IsSortable() const
Definition THaTrack.h:186
bool HasRot() const
Definition THaTrack.h:137
Double_t fTheta
Definition THaTrack.h:196
Double_t fdBeta
Definition THaTrack.h:228
void Clear(Option_t *opt="")
Definition THaTrack.cxx:27
Double_t fRPhi
Definition THaTrack.h:210
Double_t GetTY() const
Definition THaTrack.h:107
Int_t GetNDoF() const
Definition THaTrack.h:96
Double_t fY
Definition THaTrack.h:195
void SetDedx(Double_t dedx)
Definition THaTrack.h:157
Double_t GetX() const
Definition THaTrack.h:90
Double_t GetLabPy() const
Definition THaTrack.h:112
Double_t GetTPhi() const
Definition THaTrack.h:109
void Print(Option_t *opt="") const
Definition THaTrack.cxx:68
Double_t fRY
Definition THaTrack.h:208
Double_t GetRTheta() const
Definition THaTrack.h:104
THaCluster * GetCluster(Int_t i)
Definition THaTrack.h:80
Double_t GetRY() const
Definition THaTrack.h:103
Double_t fTTheta
Definition THaTrack.h:215
TVector3 fPvect
Definition THaTrack.h:219
Int_t GetTrkNum() const
Definition THaTrack.h:84
Double_t GetVertexZ() const
Definition THaTrack.h:116
void SetVertex(const TVector3 &vert)
Definition THaTrack.h:177
Double_t GetRX() const
Definition THaTrack.h:102
Double_t GetPhi() const
Definition THaTrack.h:87
Double_t fFPTime
Definition THaTrack.h:250
Int_t AddCluster(THaCluster *c)
Definition THaTrack.cxx:55
void SetGoodPlane4(Int_t gdplane4)
Definition THaTrack.h:163
Int_t fTrkNum
Definition THaTrack.h:236
Double_t GetVertexY() const
Definition THaTrack.h:115
void SetType(UInt_t flag)
Definition THaTrack.h:169
void SetGoodPlane3(Int_t gdplane3)
Definition THaTrack.h:162
void SetdTime(Double_t dt)
Definition THaTrack.h:154
Double_t fTX
Definition THaTrack.h:213
virtual Int_t Compare(const TObject *obj) const
Definition THaTrack.cxx:92
Double_t fEnergy
Definition THaTrack.h:246
Double_t GetEnergy() const
Definition THaTrack.h:128
Int_t fNDoF
Definition THaTrack.h:243
Double_t GetDp() const
Definition THaTrack.h:110
UInt_t GetType() const
Definition THaTrack.h:82
void SetTime(Double_t time)
Definition THaTrack.h:153
virtual ~THaTrack()
Definition THaTrack.cxx:19
void SetdBeta(Double_t db)
Definition THaTrack.h:156
void SetNPMT(Int_t npmt)
Definition THaTrack.h:159
void SetVertex(Double_t x, Double_t y, Double_t z)
Definition THaTrack.h:179
Double_t GetTheta() const
Definition THaTrack.h:89
Int_t fNclusters
Definition THaTrack.h:231
Double_t fBeta
Definition THaTrack.h:227
void SetR(Double_t x, Double_t y, Double_t theta, Double_t phi)
Definition THaTrack.h:271
void SetVertexError(Double_t x, Double_t y, Double_t z)
Definition THaTrack.h:183
void SetVertexError(const TVector3 &err)
Definition THaTrack.h:181
TVector3 & GetPvect()
Definition THaTrack.h:119
Double_t fTPhi
Definition THaTrack.h:216
void SetMomentum(Double_t p)
Definition THaTrack.h:170
Double_t fDX
Definition THaTrack.h:201
THaTrackID * GetID() const
Definition THaTrack.h:83
Double_t fTY
Definition THaTrack.h:214
@ kHasRot
Definition THaTrack.h:29
@ kHasVertex
Definition THaTrack.h:31
@ kHasDet
Definition THaTrack.h:27
@ kHasTarget
Definition THaTrack.h:30
Double_t GetBeta() const
Definition THaTrack.h:125
Double_t GetY() const
Definition THaTrack.h:91
void SetDp(Double_t dp)
Definition THaTrack.h:171
bool HasTarget() const
Definition THaTrack.h:138
Int_t GetIndex() const
Definition THaTrack.h:79
Int_t GetNclusters() const
Definition THaTrack.h:78
TVector3 fVertex
Definition THaTrack.h:220
void SetTrkNum(Int_t n)
Definition THaTrack.h:172
Double_t fPhi
Definition THaTrack.h:197
Int_t fIndex
Detector creating this track.
Definition THaTrack.h:235
void SetFlag(UInt_t flag)
Definition THaTrack.h:168
Double_t GetVertexX() const
Definition THaTrack.h:114
void SetTarget(Double_t x, Double_t y, Double_t theta, Double_t phi)
Definition THaTrack.h:283
bool HasDet() const
Definition THaTrack.h:135
void SetBeta(Double_t beta)
Definition THaTrack.h:155
Double_t GetdBeta() const
Definition THaTrack.h:126
THaPIDinfo * fPIDinfo
Clusters of this track.
Definition THaTrack.h:233
Double_t GetDY() const
Definition THaTrack.h:99
Double_t fDTheta
Definition THaTrack.h:203
void SetIndex(Int_t idx)
Definition THaTrack.h:174
Int_t GetGoodPlane3() const
Definition THaTrack.h:132
Double_t fChi2
Definition THaTrack.h:242
Double_t GetTime() const
Definition THaTrack.h:123
Double_t fDp
Definition THaTrack.h:217
Double_t GetP() const
Definition THaTrack.h:86
UInt_t fFlag
Track identifier.
Definition THaTrack.h:239
void SetBetaChi2(Double_t betachi2)
Definition THaTrack.h:160
Double_t GetTTheta() const
Definition THaTrack.h:108
THaTrack()
Definition THaTrack.h:35
Double_t fRX
Definition THaTrack.h:207
void SetEnergy(Double_t energy)
Definition THaTrack.h:158
Double_t GetdTime() const
Definition THaTrack.h:124
Double_t GetFPTime() const
Definition THaTrack.h:131
THaCluster * fClusters[kMAXCL]
Number of clusters.
Definition THaTrack.h:232
Double_t GetTX() const
Definition THaTrack.h:106
THaTrackingDetector * fCreator
Particle ID information for this track.
Definition THaTrack.h:234
void SetD(Double_t x, Double_t y, Double_t theta, Double_t phi)
Definition THaTrack.h:259
Double_t fDY
Definition THaTrack.h:202
Double_t GetDedx() const
Definition THaTrack.h:127
void SetCreator(THaTrackingDetector *d)
Definition THaTrack.h:173
Int_t GetGoodPlane4() const
Definition THaTrack.h:133
Double_t GetLabPz() const
Definition THaTrack.h:113
THaTrackID * fID
Definition THaTrack.h:238
void SetPIDinfo(THaPIDinfo *pid)
Definition THaTrack.h:175
void SetPathLen(Double_t pathl)
Definition THaTrack.h:152
UInt_t GetFlag() const
Definition THaTrack.h:81
THaPIDinfo * GetPIDinfo() const
Definition THaTrack.h:88
void SetChi2(Double_t chi2, Int_t ndof)
Definition THaTrack.h:165
Double_t GetRPhi() const
Definition THaTrack.h:105
THaTrack(Double_t x, Double_t y, Double_t theta, Double_t phi, THaTrackingDetector *creator=nullptr, THaTrackID *id=nullptr, THaPIDinfo *pid=nullptr)
Definition THaTrack.h:53
Double_t fP
Definition THaTrack.h:198
Double_t fDPhi
Definition THaTrack.h:204
virtual void Clear(Option_t *="")
Double_t Z() const
void SetXYZ(Double_t x, Double_t y, Double_t z)
Double_t Y() const
Double_t X() const
double beta(double x, double y)
Double_t y[n]
Double_t x[n]
const Int_t n